CbmRoot
Loading...
Searching...
No Matches
CbmCaTrackTypeQa.h
Go to the documentation of this file.
1/* Copyright (C) 2023 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Sergei Zharko [committer] */
4
9
10
11#pragma once
12
13#include "CaParameters.h"
14#include "CbmCaTrackFitQa.h"
15#include "CbmL1DetectorID.h"
16#include "CbmL1Hit.h"
17#include "CbmQaIO.h"
18#include "KfFieldRegion.h"
19#include "KfTrackKalmanFilter.h"
20
21#include <map>
22#include <string>
23
24// Forward declarations
25namespace cbm::ca::tools
26{
27 class MCData;
28 class MCTrack;
29} // namespace cbm::ca::tools
30class CbmL1Track;
32
33class TH1F;
34class TH2F;
35class TProfile;
36class TProfile2D;
37
38namespace cbm::ca
39{
44 class TrackTypeQa : virtual public CbmQaIO {
45 public:
52 TrackTypeQa(const char* typeName, const char* prefixName, bool bUseMC, std::shared_ptr<ObjList_t> pObjList);
53
55 ~TrackTypeQa() = default;
56
58 TrackTypeQa(const TrackTypeQa&) = delete;
59
62
65
68
71
74 {
75 assert(fph_stations_hit);
76 return fph_stations_hit->GetBinContent(1);
77 }
78
81 {
82 assert(fph_stations_point);
83 return fph_stations_point->GetBinContent(1);
84 }
85
88
90 double GetIntegratedEff() const
91 {
92 assert(fph_rate_reco);
93 return fph_rate_reco->GetBinContent(1);
94 }
95
97 double GetClonesRate() const
98 {
99 assert(fph_rate_clones);
100 return fph_rate_clones->GetBinContent(1);
101 }
102
104 double GetKilledRate() const
105 {
106 assert(fbUseMC);
107 return fph_rate_killed->GetBinContent(1);
108 }
109
111 double GetNofMCTracks() const { return fCounterMC; }
112
115 {
116 //return fCounterMC * GetIntegratedEff();
117 return fCounterRecoTotal;
118 }
119
122
124 double GetNofClones() const { return GetClonesRate() * fCounterMC; }
125
127 const char* GetTitle() const { return fsTitle.Data(); }
128
130 void Init();
131
135 bool IsMCUsed() const { return fbUseMC; }
136
141 void FillRecoTrack(int iTrkReco, double weight = 1);
142
147 void FillMCTrack(int iTrkMC, double weight = 1);
148
152
156
159 void RegisterMCData(tools::MCData& mcData) { fpMCData = &mcData; }
160
163 void RegisterParameters(std::shared_ptr<ca::Parameters<float>>& pParameters) { fpParameters = pParameters; }
164
170 void SetDrawAtt(Color_t markerCol = 1, Style_t markerSty = 20, Color_t lineCol = -1, Style_t lineSty = 1);
171
174 void SetTitle(const char* title) { fsTitle = title; }
175
180 void SetMCTrackCut(std::function<bool(const tools::MCTrack&)> cut) { fMCTrackCut = cut; }
181
186 void SetRecoTrackCut(std::function<bool(const CbmL1Track&)> cut) { fRecoTrackCut = cut; }
187
190
191
192 // ************************
193 // ** List of histograms **
194 // ************************
195
196 // NOTE: Naming convention:
197 // "fph_" + <"reco"/"mc"> + <quantity> + <""/"MC">. Here <"reco"/"mc"> stands for reconstructed track or true MC
198 // track, <quantity> is a name of a quantity, versus which one wants to build the dependency, <""/"MC"> stands
199 // for an origin of the quantity - reconstructed or true Monte-Carlo. If there are several quantities, they are
200 // to be separated with "_" (example: "p_yMC" -- reconstructed momentum vs MC rapidity)
201
202 // ** Histograms from reconstructed information only **
203 TH1F* fph_reco_p = nullptr;
204 TH1F* fph_reco_pt = nullptr;
205 TH1F* fph_reco_phi = nullptr;
206 TH1F* fph_reco_theta = nullptr;
207 TH2F* fph_reco_theta_phi = nullptr;
208 TH1F* fph_reco_tx = nullptr;
209 TH1F* fph_reco_ty = nullptr;
210 TH2F* fph_reco_ty_tx = nullptr;
211 TH1F* fph_reco_eta = nullptr;
212 TH1F* fph_reco_nhits = nullptr;
213 TH1F* fph_reco_fhitR = nullptr;
214 TH1F* fph_reco_fsta = nullptr;
215 TH1F* fph_reco_lsta = nullptr;
216 TH1F* fph_reco_chi2_ndf = nullptr;
217 TH1F* fph_reco_chi2_ndf_time = nullptr;
218
219 // ** Reconstructed track distributions, requiring MC information **
220 TH1F* fph_reco_pMC = nullptr;
221 TH1F* fph_reco_ptMC = nullptr;
222 TH1F* fph_reco_yMC = nullptr;
223 TH1F* fph_reco_etaMC = nullptr;
224 TH2F* fph_reco_ptMC_yMC = nullptr;
225 TH1F* fph_reco_phiMC = nullptr;
226 TH1F* fph_reco_thetaMC = nullptr;
227 TH2F* fph_reco_thetaMC_phiMC = nullptr;
228 TH1F* fph_reco_txMC = nullptr;
229 TH1F* fph_reco_tyMC = nullptr;
230
231 // ** MC track distributions **
232 TH1F* fph_mc_pMC = nullptr;
233 TH1F* fph_mc_yMC = nullptr;
234 TH1F* fph_mc_etaMC = nullptr;
235 TH2F* fph_mc_ptMC_yMC = nullptr;
236 TH1F* fph_mc_ptMC = nullptr;
237 TH1F* fph_mc_phiMC = nullptr;
238 TH1F* fph_mc_thetaMC = nullptr;
239 TH2F* fph_mc_thetaMC_phiMC = nullptr;
240 TH1F* fph_mc_txMC = nullptr;
241 TH1F* fph_mc_tyMC = nullptr;
242 TH2F* fph_mc_tyMC_txMC = nullptr;
243
244 // ** Efficiencies **
245 TProfile* fph_eff_int = nullptr;
246 TProfile* fph_eff_pMC = nullptr;
247 TProfile* fph_eff_yMC = nullptr;
248 TProfile* fph_eff_ptMC = nullptr;
249 TProfile* fph_eff_thetaMC = nullptr;
250 TProfile* fph_eff_etaMC = nullptr;
251 TProfile* fph_eff_phiMC = nullptr;
252 TProfile* fph_eff_nhitsMC = nullptr;
253 TProfile* fph_eff_txMC = nullptr;
254 TProfile* fph_eff_tyMC = nullptr;
255
256 TProfile2D* fph_eff_thetaMC_phiMC = nullptr;
257 TProfile2D* fph_eff_ptMC_yMC = nullptr;
258 TProfile2D* fph_eff_tyMC_txMC = nullptr;
259
260 // ** Clone rate **
261
262 TProfile* fph_clone_pMC = nullptr;
263 TProfile* fph_clone_yMC = nullptr;
264 TProfile* fph_clone_ptMC = nullptr;
265 TProfile* fph_clone_thetaMC = nullptr;
266 TProfile* fph_clone_etaMC = nullptr;
267 TProfile* fph_clone_phiMC = nullptr;
268 TProfile* fph_clone_nhitsMC = nullptr;
269 TProfile* fph_clone_txMC = nullptr;
270 TProfile* fph_clone_tyMC = nullptr;
271
272 TProfile2D* fph_clone_thetaMC_phiMC = nullptr;
273 TProfile2D* fph_clone_ptMC_yMC = nullptr;
274 TProfile2D* fph_clone_tyMC_txMC = nullptr;
275
276 // ** Fit QA **
277 std::unique_ptr<TrackFitQa> fpFitQaFirstMCpoint = nullptr;
278 std::unique_ptr<TrackFitQa> fpFitQaFirstHit = nullptr;
279 std::unique_ptr<TrackFitQa> fpFitQaLastHit = nullptr;
280 std::unique_ptr<TrackFitQa> fpFitQaVertex = nullptr;
281
282 // TODO: Provide integrated efficiencies
283
284 private:
287 virtual void SetTH1Properties(TH1* pHist) const override;
288
289 // ** Technical profiles and histograms (counters) **
290 TProfile* fph_rate_reco = nullptr;
291 TProfile* fph_rate_killed = nullptr;
292 TProfile* fph_rate_clones = nullptr;
293 TProfile* fph_stations_point = nullptr;
294 TProfile* fph_stations_hit = nullptr;
295
298
299 int fCounterMC = 0;
302 double fRecoLength = 0.;
303 double fMCLength = 0.;
304 double fFakeLength = 0.;
305
306
307 bool fbUseMC = false;
308 TString fsTitle = "";
309
310 // TODO: SZh 20.03.2024: Maybe replace CbmL1Track with CaTrack?
314 std::shared_ptr<ca::Parameters<float>> fpParameters = nullptr;
315
316
317 // ** Cuts on tracks for a given track class **
318
320 std::function<bool(const tools::MCTrack&)> fMCTrackCut = [](const tools::MCTrack&) { return true; };
321
323 std::function<bool(const CbmL1Track&)> fRecoTrackCut = [](const CbmL1Track&) { return true; };
324
325
326 // **************************
327 // ** Histogram properties **
328 // **************************
329
330 // ** Binning **
331 // TODO: SZh 20.03.2024: Make yaml-configurable
332 static constexpr int kBinsP = 120;
333 static constexpr double kLoP = 0.;
334 static constexpr double kUpP = 12.;
335 static constexpr int kBinsPT = 100;
336 static constexpr double kLoPT = 0.;
337 static constexpr double kUpPT = 4.;
338 static constexpr int kBinsETA = 40;
339 static constexpr double kLoETA = 0.;
340 static constexpr double kUpETA = 4.;
341 static constexpr int kBinsY = 40;
342 static constexpr double kLoY = 0.;
343 static constexpr double kUpY = 4.;
344 static constexpr int kBinsPHI = 68;
345 static constexpr double kLoPHI = -3.2;
346 static constexpr double kUpPHI = +3.2;
347 static constexpr int kBinsTHETA = 68;
348 static constexpr double kLoTHETA = 0.;
349 static constexpr double kUpTHETA = 3.2;
350 static constexpr int kBinsTX = 80;
351 static constexpr double kLoTX = -2.;
352 static constexpr double kUpTX = +2.;
353 static constexpr int kBinsTY = 80;
354 static constexpr double kLoTY = -2.;
355 static constexpr double kUpTY = +2.;
356 static constexpr int kBinsNHITS = 15;
357 static constexpr double kLoNHITS = -0.5;
358 static constexpr double kUpNHITS = 14.5;
359 static constexpr int kBinsFHITR = 50;
360 static constexpr double kLoFHITR = 0.;
361 static constexpr double kUpFHITR = 150.;
362 static constexpr int kBinsNSTA = 12;
363 static constexpr double kLoNSTA = -0.5;
364 static constexpr double kUpNSTA = 11.5;
365 static constexpr int kBinsCHI2NDF = 200;
366 static constexpr double kLoCHI2NDF = 0.;
367 static constexpr double kUpCHI2NDF = 20.;
368
369 // ** Drawing attributes **
370 Color_t fMarkerColor = 1;
371 Style_t fMarkerStyle = 20;
372 Color_t fLineColor = 1;
373 Style_t fLineStyle = 1;
374
376 };
377
378} // namespace cbm::ca
QA submodule for track fit results (residuals and puls) at selected z-coordinate (header)
Implementation of L1DetectorID enum class for CBM.
Module for ROOT objects IO interface (header)
Magnetic flux density interpolation along the track vs. z-coordinate (header)
Track fit utilities for the CA tracking based on the Kalman filter.
ROOT object IO interface for QA.
Definition CbmQaIO.h:44
A container for all external parameters of the CA tracking algorithm.
Magnetic field region, corresponding to a hit triplet.
Unified QA for a group of tracks.
static constexpr double kLoTY
Lower boundary, slope along y.
static constexpr double kLoCHI2NDF
Lower boundary, chi2 over NDF.
int fCounterMC
Counter of MC tracks.
double GetNofMCTracks() const
Gets number of MC tracks.
static constexpr double kUpCHI2NDF
Upper boundary, chi2 over NDF.
tools::MCData * fpMCData
Pointer to MC data object.
bool fbUseMC
Flag: true - MC information is used.
double GetAverageFakeLength() const
Gets average fake length.
TH1F * fph_mc_ptMC
Transverse momentum over charge of MC tracks.
TH1F * fph_reco_nhits
Hit number of reconstructed tracks.
std::unique_ptr< TrackFitQa > fpFitQaLastHit
void FillRecoTrack(int iTrkReco, double weight=1)
Fills histograms with various track information.
static constexpr double kLoY
Lower boundary, rapidity.
TH1F * fph_reco_etaMC
MC pseudo-rapidity of reconstructed tracks.
void SetDrawAtt(Color_t markerCol=1, Style_t markerSty=20, Color_t lineCol=-1, Style_t lineSty=1)
Sets drawing attributes for histograms.
TH2F * fph_reco_ptMC_yMC
MC transverse momentum vs MC rapidity of reconstructed tracks.
TH1F * fph_reco_fhitR
Distance of the first hit from z-axis for reconstructed tracks.
void SetMCTrackCut(std::function< bool(const tools::MCTrack &)> cut)
Sets selection cuts on MC tracks.
double GetAverageNofStationsWithHit() const
Gets average reconstructed track length.
TProfile * fph_clone_ptMC
Efficiency vs. MC transverse momentum.
TProfile * fph_clone_nhitsMC
Efficiency vs. MC number of hits (total number of stations with a)
TProfile * fph_stations_hit
Average number of stations with hit.
TProfile2D * fph_eff_ptMC_yMC
Efficiency vs. MC transverse momentum and MC rapidity.
static constexpr double kLoP
Lower boundary, total momentum [GeV/c].
static constexpr double kLoETA
Lower boundary, pseudo-rapidity.
Style_t fLineStyle
Line style.
static constexpr double kLoNSTA
Lower boundary, number of stations.
TrackTypeQa(const TrackTypeQa &)=delete
Copy constructor.
TH2F * fph_mc_thetaMC_phiMC
Polar angle vs. azimuthal angle of MC tracks.
static constexpr int kBinsTHETA
Number of bins, polar angle.
TH1F * fph_reco_tyMC
MC Slope along y-axis of reconstructed tracks.
static constexpr int kBinsCHI2NDF
Number of bins, chi2 over NDF.
TH1F * fph_reco_fsta
First station index of reconstructed tracks.
TH1F * fph_reco_phiMC
MC Azimuthal angle of reconstructed tracks.
TProfile2D * fph_clone_tyMC_txMC
Efficiency vs. MC slopes.
double GetClonesRate() const
Gets clones rate.
int fCounterRecoTotal
Counter of reco tracks (total = reco + ghost + clones)
std::unique_ptr< TrackFitQa > fpFitQaFirstMCpoint
TH1F * fph_reco_txMC
MC Slope along x-axis of reconstructed tracks.
TH1F * fph_mc_phiMC
Azimuthal angle of MC tracks.
ca::Vector< CbmL1Track > * fpvRecoTracks
Pointer to vector of reconstructed tracks.
static constexpr int kBinsTX
Number of bins, slope along x.
static constexpr int kBinsFHITR
Number of bins, transverse dist. of the 1st hit from z-axis.
TProfile * fph_eff_pMC
Efficiency vs. MC momentum.
TrackTypeQa & operator=(const TrackTypeQa &)=delete
Copy assignment operator.
void SetRecoTrackCut(std::function< bool(const CbmL1Track &)> cut)
Sets selection cuts on reconstructed tracks.
ca::Vector< CbmL1HitDebugInfo > * fpvHits
Pointer to vector of reconstructed hits.
double fMCLength
Total length of MC tracks.
TH1F * fph_mc_pMC
MC total momentum over charge of MC tracks.
TProfile * fph_eff_yMC
Efficiency vs. MC rapidity.
TProfile * fph_eff_thetaMC
Efficiency vs. MC polar angle.
double GetAverageNofStationsWithPoint() const
Gets average MC track length.
TH1F * fph_reco_lsta
Last station index of reconstructed tracks.
bool IsMCUsed() const
Flag on MC usage by the track type.
TH1F * fph_mc_txMC
Slope along x-axis of MC tracks.
TProfile * fph_clone_etaMC
Efficiency vs. MC pseudorapidity.
static constexpr int kBinsPT
Number of bins, transverse momentum.
TProfile * fph_rate_clones
Rate of clone tracks / mc.
TH1F * fph_reco_tx
Slope along x-axis of reconstructed tracks.
TProfile * fph_eff_ptMC
Efficiency vs. MC transverse momentum.
TH2F * fph_mc_ptMC_yMC
MC transverse momentum vs. MC rapidity of MC tracks.
static constexpr double kUpTX
Upper boundary, slope along x.
void Init()
Initializes histograms.
static constexpr double kLoPHI
Lower boundary, azimuthal angle [rad].
static constexpr double kUpTHETA
Upper boundary, polar angle [rad].
static constexpr double kUpNSTA
Upper boundary, number of stations.
static constexpr double kUpY
Upper boundary, rapidity.
void RegisterRecoTracks(ca::Vector< CbmL1Track > &vTracks)
Registers the reconstructed tracks container.
Color_t fMarkerColor
Marker color.
TH1F * fph_reco_eta
Pseudo-rapidity of reconstructed tracks.
std::unique_ptr< TrackFitQa > fpFitQaFirstHit
std::function< bool(const CbmL1Track &)> fRecoTrackCut
Cut function on reconstructed tracks.
double GetNofRecoTracksTotal() const
Gets total number of reconstructed tracks (reco + ghost + clones)
TH1F * fph_reco_phi
Azimuthal angle of reconstructed tracks.
TString fsTitle
Title of the track category.
TH1F * fph_reco_ptMC
MC transverse momentum over charge of reconstructed tracks.
static constexpr int kBinsNHITS
Number of bins, number of hits.
TProfile2D * fph_clone_ptMC_yMC
Efficiency vs. MC transverse momentum and MC rapidity.
TH1F * fph_reco_chi2_ndf_time
Fit chi2 over NDF of reconstructed tracks for time.
static constexpr double kUpFHITR
Upper boundary, transverse dist. of the 1st hit from z-axis [cm].
TProfile * fph_clone_tyMC
Efficiency vs. MC slope along y-axis.
TH2F * fph_reco_ty_tx
Slope along x-axis vs y-axis of reconstructed tracks.
TProfile * fph_clone_phiMC
Efficiency vs. MC azimuthal angle.
double fFakeLength
Total length of fake tracks.
static constexpr int kBinsPHI
Number of bins, azimuthal angle.
TH2F * fph_reco_theta_phi
Polar angle vs. azimuthal angle of reconstructed tracks.
TH1F * fph_reco_chi2_ndf
Fit chi2 over NDF of reconstructed tracks.
static constexpr int kBinsP
Number of bins, total momentum.
TH1F * fph_mc_yMC
MC rapidity of MC tracks.
TH1F * fph_reco_ty
Slope along y-axis of reconstructed tracks.
TH1F * fph_mc_thetaMC
Polar angle of MC tracks.
static constexpr double kUpPT
Upper boundary, transverse momentum [GeV/c].
TProfile * fph_clone_pMC
Clone rate vs. MC momentum.
int fCounterClones
Counter of clone tracks.
TProfile * fph_clone_txMC
Efficiency vs. MC slope along x-axis.
~TrackTypeQa()=default
Destructor.
static constexpr double kUpNHITS
Upper boundary, number of hits.
const char * GetTitle() const
Gets title of track group.
void RegisterMCData(tools::MCData &mcData)
Registers the MC data.
TrackTypeQa(TrackTypeQa &&)=delete
Move constructor.
TProfile * fph_rate_killed
Rate of killed tracks / mc.
TProfile * fph_rate_reco
Rate of reconstructed tracks / mc.
TH2F * fph_mc_tyMC_txMC
Slope along x-axis vs y-axis of MC tracks.
double GetAverageRecoLength() const
Gets average reconstructed track length.
TH1F * fph_reco_pMC
MC total momentum over charge of reconstructed tracks.
double GetNofRecoTracksMatched() const
Gets number of total reco tracks.
TProfile * fph_eff_txMC
Efficiency vs. MC slope along x-axis.
TH1F * fph_reco_yMC
MC rapidity of reconstructed tracks.
void SetTitle(const char *title)
Sets title, which is to be reflected on legends and titles.
std::shared_ptr< ca::Parameters< float > > fpParameters
Pointer to parameters object.
TProfile * fph_eff_nhitsMC
Efficiency vs. MC number of hits (total number of stations with a)
static constexpr double kUpPHI
Upper boundary, azimuthal angle [rad].
Color_t fLineColor
Line color.
void RegisterParameters(std::shared_ptr< ca::Parameters< float > > &pParameters)
Registers CA parameters object.
TProfile2D * fph_clone_thetaMC_phiMC
Efficiency vs. MC theta and MC phi.
TProfile * fph_eff_phiMC
Efficiency vs. MC azimuthal angle.
void RegisterRecoHits(ca::Vector< CbmL1HitDebugInfo > &vHits)
Register the reconstructed hits container.
cbm::algo::kf::TrackKalmanFilter< double > fTrackFit
Track fitter.
TProfile2D * fph_eff_thetaMC_phiMC
Efficiency vs. MC theta and MC phi.
TH1F * fph_mc_etaMC
MC pseudo-rapidity of MC tracks.
TH1F * fph_reco_pt
Transverse momentum over charge of reconstructed tracks.
TH1F * fph_reco_p
Total momentum over charge of reconstructed tracks.
std::unique_ptr< TrackFitQa > fpFitQaVertex
TrackTypeQa & operator=(TrackTypeQa &&)=delete
Move assignment operator.
static constexpr double kLoTHETA
Lower boundary, polar angle [rad].
ClassDefNV(TrackTypeQa, 0)
static constexpr double kLoNHITS
Lower boundary, number of hits.
static constexpr double kUpTY
Upper boundary, slope along y.
double GetIntegratedEff() const
Gets integrated efficiency.
TH1F * fph_reco_thetaMC
MC Polar angle of reconstructed tracks.
double fRecoLength
Total length of reconstructed tracks.
static constexpr double kLoPT
Lower boundary, transverse momentum [GeV/c].
TH2F * fph_reco_thetaMC_phiMC
MC Polar angle vs. azimuthal angle of reconstructed tracks.
cbm::algo::kf::FieldRegion< double > fFieldRegion
Magnetic field.
TProfile * fph_eff_tyMC
Efficiency vs. MC slope along y-axis.
TProfile * fph_eff_int
Integrated efficiency.
double GetKilledRate() const
Gets killed rate.
std::function< bool(const tools::MCTrack &)> fMCTrackCut
Cut function on MC tracks.
static constexpr int kBinsTY
Number of bins, slope along y.
double GetNofClones() const
Gets number of clones.
static constexpr double kLoTX
Lower boundary, slope along x.
TProfile * fph_eff_etaMC
Efficiency vs. MC pseudorapidity.
void FillPerEvent()
Fills counters per event.
virtual void SetTH1Properties(TH1 *pHist) const override
Overrided virtual function of the CbmQaIO class, defines properties of the histograms.
TH1F * fph_reco_theta
Polar angle of reconstructed tracks.
static constexpr int kBinsNSTA
Number of bins, number of stations.
TProfile * fph_clone_thetaMC
Efficiency vs. MC polar angle.
TProfile2D * fph_eff_tyMC_txMC
Efficiency vs. MC slopes.
TProfile * fph_clone_yMC
Efficiency vs. MC rapidity.
void FillMCTrack(int iTrkMC, double weight=1)
Fills histograms with mc track information.
static constexpr int kBinsY
Number of bins, rapidity.
static constexpr double kLoFHITR
Lower boundary, transverse dist. of the 1st hit from z-axis [cm].
static constexpr int kBinsETA
Number of bins, pseudo-rapidity.
TH1F * fph_mc_tyMC
Slope along y-axis of MC tracks.
static constexpr double kUpETA
Upper boundary, pseudo-rapidity.
Style_t fMarkerStyle
Marker style.
static constexpr double kUpP
Upper boundary, total momentum [GeV/c].
TrackTypeQa(const char *typeName, const char *prefixName, bool bUseMC, std::shared_ptr< ObjList_t > pObjList)
Constructor.
TProfile * fph_stations_point
Average number of stations with MC point.
This class represents a package for tracking-related data.