CbmRoot
Loading...
Searching...
No Matches
CbmMvdMcQaTask.h
Go to the documentation of this file.
1/* Copyright (C) 2025 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Ajit Kumar [committer] */
4
9
10#ifndef CbmMvdMcQaTask_h
11#define CbmMvdMcQaTask_h 1
12
13#include "CbmQaTask.h"
14
15class TClonesArray;
16class CbmMCEventList;
17class CbmMCDataArray;
19class FairRootManager;
20namespace cbm::mvd
21{
23 class CbmMvdMcQaTask : public CbmQaTask {
24 public:
25 CbmMvdMcQaTask(int verbose, bool isMcUsed);
26
27 void Check() override;
28
29 void CreateSummary() override;
30
31 void ExecQa() override;
32
33 InitStatus InitQa() override;
34
36 void SetKeVRange(int maxKeV, int binKeV = 1)
37 {
38 fKeVMax = maxKeV;
39 fKeVBin = (binKeV > 0 ? binKeV : 1);
40 }
41
42 enum class Origin
43 {
44 All = 0,
49 };
50 int MapPdgToSpeciesBin_(int pdg) const;
51 int GetStationIdFromZ_(double z) const;
52 // static std::string PdgToName_(int pdg);
55 void SetNumberOfStations(int nstation) { fNStations = nstation; }
56
57 private:
58 void ProcessMcTracks(const int& fileId, const int& eventId);
59 void ProcessMvdPoints(const int& fileId, const int& eventId);
60
62 FairRootManager* fManager = nullptr;
64
69
71 std::vector<int> fPdgBuf;
72 std::vector<Origin> fOriginBuf;
73
75
77 Long64_t fAllNPrim = 0;
78 Long64_t fAllNSec = 0;
79 std::unordered_map<int, uint64_t> fAllCountsPrim;
80 std::unordered_map<int, uint64_t> fAllCountsSec;
81
82 static constexpr int kNumOrigins = 5;
83 std::array<std::unordered_map<int, TH1F*>, kNumOrigins> fHistsElossKeVByPDG;
84 std::array<TList*, kNumOrigins> fLists{};
85 static constexpr const char* kOriginName[kNumOrigins] = {"all", "primary", "secondary", "deltaSec", "other"};
86 TH1F* GetOrBookElossHist(int pdg, Origin origin);
87
89 int fKeVMax = 2000;
90 int fKeVBin = 1;
91
93 double fZMinStation0 = -38, fZMaxStation0 = -34;
94 double fZMinStation1 = -34, fZMaxStation1 = -30;
95 double fZMinStation2 = -30, fZMaxStation2 = -26;
96 double fZMinStation3 = -26, fZMaxStation3 = -22;
97
99 TH1D* fh_nTracksAll = nullptr;
100 TH1D* fh_nMvdPointsEvt = nullptr;
101 TH1D* fh_nHitTracksEvt = nullptr;
102 TH1D* fh_nDeltaEvt = nullptr;
103
105 TH1D* fh_mc_pdg_all = nullptr;
106 TH1D* fh_mc_pdg_prim = nullptr;
107 TH1D* fh_mc_pdg_sec = nullptr;
108 TH1D* fh_mc_p_prim = nullptr;
109 TH1D* fh_mc_p_sec = nullptr;
110 TH1D* fh_mc_startZ = nullptr;
111 TH2D* fh_mc_startX_vs_Z = nullptr;
112 TH2D* fh_mc_startY_vs_Z = nullptr;
113 TH2D* fh_mc_startY_vs_X = nullptr;
114
116 // per-event per-station
117 std::vector<TH1D*> fh_mvd_nPointsEvt_station;
118
119 // XY maps per station, split by origin
120 std::vector<TH2D*> fh_mvd_xy_prim_station;
121 std::vector<TH2D*> fh_mvd_xy_sec_station;
122
123 // species-labelled PDG histos at MVD point level
124 TH1D* fh_mvd_pdg_all = nullptr;
125 TH1D* fh_mvd_pdg_prim = nullptr;
126 TH1D* fh_mvd_pdg_sec = nullptr;
127
128 // momentum distributions (p from parent MCTrack)
129 TH1D* fh_mvd_p_prim = nullptr;
130 TH1D* fh_mvd_p_sec = nullptr;
131
132 // energy loss correlations
133 TH2D* fh_mvd_dEkeV_vs_p = nullptr;
134 TH2D* fh_mvd_dEkeV_vs_pdg = nullptr;
135
136 TH1D* fh_sec_count_bySpecies = nullptr;
137 TH2D* fh_dEkeV_vs_pidClass = nullptr;
138 TH2D* fh_xy_points_nonprim = nullptr;
139 std::vector<TH2D*> fh_xy_points;
140 std::vector<TH1D*> fhCountMvdPoint;
141
142 TH1D* fhZdistMvd = nullptr;
143 TH2D* fh_dEkeV_vs_p = nullptr;
144
146 void BookHists_();
147 static bool DecodeIonPDG_(int pdg, int& Z, int& A);
148
149 int PidClass_(int pdg) const;
150 Origin ClassifyTrack_(const int& fileId, const int& eventId, Int_t trackindex) const;
151 double PointELossGeV_(const void* vp) const;
152
154 };
155
156} // namespace cbm::mvd
157
158#endif // CbmMvdCbmMvdMcQaTask_h
A base class for CBM QA task logic.
int Int_t
Some class A.
Access to a MC data branch for time-based analysis.
Task class creating and managing CbmMCDataArray objects.
Container class for MC events with number, file and start time.
CbmQaTask(const char *name, int verbose, bool isMCUsed, ECbmRecoMode recoMode=ECbmRecoMode::Timeslice)
Constructor from parameters.
Definition CbmQaTask.cxx:32
TH1D * fh_nMvdPointsEvt
MVD points per event.
void Check() override
Function to check, if the QA results are acceptable.
void ExecQa() override
Method to fill histograms per event or time-slice.
Long64_t fAllNPrim
global counters (over all events)
TH1D * fh_sec_count_bySpecies
overall counts (non-primary tracks with hits) by coarse species class
Int_t fNStations
Number of stations in MVD: TODO: Get Number of station from data file itself.
std::array< TList *, kNumOrigins > fLists
TH2D * fh_dEkeV_vs_pidClass
dE per point vs species class (non-primaries with hits)
void SetKeVRange(int maxKeV, int binKeV=1)
dE histogram settings (keV)
TH2D * fh_xy_points_nonprim
x vs y map of points (non-primaries)
TH2D * fh_mc_startY_vs_Z
start Y vs Z
Origin
track origin tag (public to avoid any ROOT dictionary issues)
std::vector< Origin > fOriginBuf
origin by track index
double fZMinStation0
TODO: Remove when we have getter for number of station directly from CbmMvdPoint.
TH2D * fh_mvd_dEkeV_vs_p
dE (keV) vs p (GeV/c) at MVD point
ClassDefOverride(CbmMvdMcQaTask, 1)
TH1D * fh_mc_pdg_sec
species for secondary MCTracks
TH2D * fh_mc_startY_vs_X
start Y vs X
CbmMCDataArray * fMcTracks
CbmMCEventList * fMcEventList
containers
int PidClass_(int pdg) const
static bool DecodeIonPDG_(int pdg, int &Z, int &A)
TH2D * fh_mvd_dEkeV_vs_pdg
dE (keV) vs species bin
std::unordered_map< int, uint64_t > fAllCountsSec
PDG -> count (all secondaries)
static constexpr const char * kOriginName[kNumOrigins]
std::vector< int > fPdgBuf
per-event caches
TH1D * fh_mc_p_sec
momentum of secondaries
CbmMCDataManager * fMcManager
std::vector< TH1D * > fh_mvd_nPointsEvt_station
MvdPoint histograms.
TH2D * fh_mc_startX_vs_Z
start X vs Z
int MapPdgToSpeciesBin_(int pdg) const
std::vector< TH2D * > fh_mvd_xy_sec_station
X-Y for secondary MVD points per station.
std::array< std::unordered_map< int, TH1F * >, kNumOrigins > fHistsElossKeVByPDG
TH1D * fh_mc_pdg_all
MCTrack-level histograms.
TH1D * fh_mc_startZ
start-z of MCTracks
TH1F * GetOrBookElossHist(int pdg, Origin origin)
double PointELossGeV_(const void *vp) const
CbmMvdMcQaTask(int verbose, bool isMcUsed)
std::unordered_map< int, uint64_t > fAllCountsPrim
PDG -> count (all primaries)
std::vector< TH2D * > fh_xy_points
x vs y map of points vs. station ID (non-primaries)
TH1D * fh_mvd_pdg_sec
species for secondary MVD points
TH1D * fh_mvd_pdg_prim
species for primary MVD points
void ProcessMcTracks(const int &fileId, const int &eventId)
FairRootManager * fManager
setup
TH1D * fh_mc_p_prim
momentum of primaries
InitStatus InitQa() override
Initializes the task.
void ProcessMvdPoints(const int &fileId, const int &eventId)
Origin ClassifyTrack_(const int &fileId, const int &eventId, Int_t trackindex) const
TH1D * fh_mc_pdg_prim
species for primary MCTracks
int GetStationIdFromZ_(double z) const
TH1D * fh_nDeltaEvt
number of secondary e+- with hits per event
TH1D * fh_nHitTracksEvt
unique non-primary tracks with hits per event
std::vector< TH2D * > fh_mvd_xy_prim_station
X-Y for primary MVD points per station.
void CreateSummary() override
Initializes QA-task summary: canvases, tables etc.
static constexpr int kNumOrigins
CbmMCDataArray * fMvdPoints
TH1D * fh_mvd_pdg_all
species for all MVD points
TH2D * fh_dEkeV_vs_p
dE (keV) vs p (GeV/c) at point (non-primaries)
void SetNumberOfStations(int nstation)
std::vector< TH1D * > fhCountMvdPoint
TH1D * fh_nTracksAll
Histograms.
int fKeVBin
1 keV bins by default