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 "CaDebugger.h"
14#include "CaMcData.h"
15#include "CaParameters.h"
16#include "CbmCaTrackFitQa.h"
17#include "CbmL1DetectorID.h"
18#include "CbmQaIO.h"
19
20#include <map>
21#include <string>
22
23namespace cbm::algo::ca
24{
25 class McHitInfo;
26} // namespace cbm::algo::ca
27
28class CbmL1Track;
29
30class TH1F;
31class TH2F;
32class TProfile;
33class TProfile2D;
34
35namespace cbm::ca
36{
41 class TrackTypeQa : virtual public CbmQaIO {
42 public:
49 TrackTypeQa(const char* typeName, const char* prefixName, bool bUseMC, std::shared_ptr<ObjList_t> pObjList);
50
52 ~TrackTypeQa() = default;
53
55 TrackTypeQa(const TrackTypeQa&) = delete;
56
59
62
65
66 void SetDebugger(cbm::algo::ca::utils::Debugger* pDebugger) { fpDebugger = pDebugger; }
67
73
81
84
87 {
88 assert(fph_stations_hit);
89 return fph_stations_hit->GetBinContent(1);
90 }
91
94 {
95 assert(fph_stations_point);
96 return fph_stations_point->GetBinContent(1);
97 }
98
101
103 double GetIntegratedEff() const
104 {
105 assert(fph_rate_reco);
106 return fph_rate_reco->GetBinContent(1);
107 }
108
110 double GetClonesRate() const
111 {
112 assert(fph_rate_clones);
113 return fph_rate_clones->GetBinContent(1);
114 }
115
117 double GetKilledRate() const
118 {
119 assert(fbUseMC);
120 return fph_rate_killed->GetBinContent(1);
121 }
122
124 double GetNofMCTracks() const { return fCounterMC; }
125
128 {
129 //return fCounterMC * GetIntegratedEff();
130 return fCounterRecoTotal;
131 }
132
135
137 double GetNofClones() const { return GetClonesRate() * fCounterMC; }
138
140 const char* GetTitle() const { return fsTitle.Data(); }
141
143 void Init();
144
148 bool IsMCUsed() const { return fbUseMC; }
149
154 void FillRecoTrack(int iTrkReco, const FitQaPointSet& fitPoints, double weight = 1);
155
160 void FillMCTrack(int iTrkMC, double weight = 1);
161
165
169
172 void RegisterMCData(cbm::algo::ca::McData& mcData) { fpMCData = &mcData; }
173
176 void RegisterParameters(std::shared_ptr<const ca::Parameters<double>>& pParameters) { fpParameters = pParameters; }
177
183 void SetDrawAtt(Color_t markerCol = 1, Style_t markerSty = 20, Color_t lineCol = -1, Style_t lineSty = 1);
184
187 void SetTitle(const char* title) { fsTitle = title; }
188
193 void SetMCTrackCut(std::function<bool(const cbm::algo::ca::McTrack&)> cut) { fMCTrackCut = cut; }
194
199 void SetRecoTrackCut(std::function<bool(const CbmL1Track&)> cut) { fRecoTrackCut = cut; }
200
203
206
207 // ************************
208 // ** List of histograms **
209 // ************************
210
211 // NOTE: Naming convention:
212 // "fph_" + <"reco"/"mc"> + <quantity> + <""/"MC">. Here <"reco"/"mc"> stands for reconstructed track or true MC
213 // track, <quantity> is a name of a quantity, versus which one wants to build the dependency, <""/"MC"> stands
214 // for an origin of the quantity - reconstructed or true Monte-Carlo. If there are several quantities, they are
215 // to be separated with "_" (example: "p_yMC" -- reconstructed momentum vs MC rapidity)
216
217 // ** Histograms from reconstructed information only **
218 TH1F* fph_reco_p = nullptr;
219 TH1F* fph_reco_pt = nullptr;
220 TH1F* fph_reco_phi = nullptr;
221 TH1F* fph_reco_theta = nullptr;
222 TH2F* fph_reco_theta_phi = nullptr;
223 TH1F* fph_reco_tx = nullptr;
224 TH1F* fph_reco_ty = nullptr;
225 TH2F* fph_reco_ty_tx = nullptr;
226 TH1F* fph_reco_eta = nullptr;
227 TH1F* fph_reco_nhits = nullptr;
228 TH1F* fph_reco_fhitR = nullptr;
229 TH1F* fph_reco_fsta = nullptr;
230 TH1F* fph_reco_lsta = nullptr;
231 TH1F* fph_reco_chi2_ndf = nullptr;
232 TH1F* fph_reco_chi2_ndf_time = nullptr;
233
234 // ** Reconstructed track distributions, requiring MC information **
235 TH1F* fph_reco_pMC = nullptr;
236 TH1F* fph_reco_ptMC = nullptr;
237 TH1F* fph_reco_yMC = nullptr;
238 TH1F* fph_reco_etaMC = nullptr;
239 TH2F* fph_reco_ptMC_yMC = nullptr;
240 TH1F* fph_reco_phiMC = nullptr;
241 TH1F* fph_reco_thetaMC = nullptr;
242 TH2F* fph_reco_thetaMC_phiMC = nullptr;
243 TH1F* fph_reco_txMC = nullptr;
244 TH1F* fph_reco_tyMC = nullptr;
245
246 // ** MC track distributions **
247 TH1F* fph_mc_pMC = nullptr;
248 TH1F* fph_mc_yMC = nullptr;
249 TH1F* fph_mc_etaMC = nullptr;
250 TH2F* fph_mc_ptMC_yMC = nullptr;
251 TH1F* fph_mc_ptMC = nullptr;
252 TH1F* fph_mc_phiMC = nullptr;
253 TH1F* fph_mc_thetaMC = nullptr;
254 TH2F* fph_mc_thetaMC_phiMC = nullptr;
255 TH1F* fph_mc_txMC = nullptr;
256 TH1F* fph_mc_tyMC = nullptr;
257 TH2F* fph_mc_tyMC_txMC = nullptr;
258
259 // ** Efficiencies **
260 TProfile* fph_eff_int = nullptr;
261 TProfile* fph_eff_pMC = nullptr;
262 TProfile* fph_eff_yMC = nullptr;
263 TProfile* fph_eff_ptMC = nullptr;
264 TProfile* fph_eff_thetaMC = nullptr;
265 TProfile* fph_eff_etaMC = nullptr;
266 TProfile* fph_eff_phiMC = nullptr;
267 TProfile* fph_eff_nhitsMC = nullptr;
268 TProfile* fph_eff_txMC = nullptr;
269 TProfile* fph_eff_tyMC = nullptr;
270
271 TProfile2D* fph_eff_thetaMC_phiMC = nullptr;
272 TProfile2D* fph_eff_ptMC_yMC = nullptr;
273 TProfile2D* fph_eff_tyMC_txMC = nullptr;
274
275 // ** Clone rate **
276
277 TProfile* fph_clone_pMC = nullptr;
278 TProfile* fph_clone_yMC = nullptr;
279 TProfile* fph_clone_ptMC = nullptr;
280 TProfile* fph_clone_thetaMC = nullptr;
281 TProfile* fph_clone_etaMC = nullptr;
282 TProfile* fph_clone_phiMC = nullptr;
283 TProfile* fph_clone_nhitsMC = nullptr;
284 TProfile* fph_clone_txMC = nullptr;
285 TProfile* fph_clone_tyMC = nullptr;
286
287 TProfile2D* fph_clone_thetaMC_phiMC = nullptr;
288 TProfile2D* fph_clone_ptMC_yMC = nullptr;
289 TProfile2D* fph_clone_tyMC_txMC = nullptr;
290
291 // ** Fit QA **
292 std::unique_ptr<TrackFitQa> fpFitQaFirstHit = nullptr;
293 std::unique_ptr<TrackFitQa> fpFitQaLastHit = nullptr;
294 std::unique_ptr<TrackFitQa> fpFitQaVertex = nullptr;
295 std::unique_ptr<TrackFitQa> fpFitQaTrd = nullptr;
296
297 // TODO: Provide integrated efficiencies
298
299 private:
302 virtual void SetTH1Properties(TH1* pHist) const override;
303
304 // ** Technical profiles and histograms (counters) **
305 TProfile* fph_rate_reco = nullptr;
306 TProfile* fph_rate_killed = nullptr;
307 TProfile* fph_rate_clones = nullptr;
308 TProfile* fph_stations_point = nullptr;
309 TProfile* fph_stations_hit = nullptr;
310
311 int fCounterMC = 0;
314 double fRecoLength = 0.;
315 double fMCLength = 0.;
316 double fFakeLength = 0.;
317
318
319 bool fbUseMC = false;
320 TString fsTitle = "";
321
322 // TODO: SZh 20.03.2024: Maybe replace CbmL1Track with CaTrack?
326 std::shared_ptr<const ca::Parameters<double>> fpParameters = nullptr;
327
328
329 // ** Cuts on tracks for a given track class **
330
332 std::function<bool(const cbm::algo::ca::McTrack&)> fMCTrackCut = [](const cbm::algo::ca::McTrack&) { return true; };
333
335 std::function<bool(const CbmL1Track&)> fRecoTrackCut = [](const CbmL1Track&) { return true; };
336
337
338 // **************************
339 // ** Histogram properties **
340 // **************************
341
342 // ** Binning **
343 // TODO: SZh 20.03.2024: Make yaml-configurable
344 static constexpr int kBinsP = 120;
345 static constexpr double kLoP = 0.;
346 static constexpr double kUpP = 12.;
347 static constexpr int kBinsPT = 100;
348 static constexpr double kLoPT = 0.;
349 static constexpr double kUpPT = 4.;
350 static constexpr int kBinsETA = 40;
351 static constexpr double kLoETA = 0.;
352 static constexpr double kUpETA = 4.;
353 static constexpr int kBinsY = 40;
354 static constexpr double kLoY = 0.;
355 static constexpr double kUpY = 4.;
356 static constexpr int kBinsPHI = 68;
357 static constexpr double kLoPHI = -3.2;
358 static constexpr double kUpPHI = +3.2;
359 static constexpr int kBinsTHETA = 68;
360 static constexpr double kLoTHETA = 0.;
361 static constexpr double kUpTHETA = 3.2;
362 static constexpr int kBinsTX = 80;
363 static constexpr double kLoTX = -2.;
364 static constexpr double kUpTX = +2.;
365 static constexpr int kBinsTY = 80;
366 static constexpr double kLoTY = -2.;
367 static constexpr double kUpTY = +2.;
368 static constexpr int kBinsNHITS = 15;
369 static constexpr double kLoNHITS = -0.5;
370 static constexpr double kUpNHITS = 14.5;
371 static constexpr int kBinsFHITR = 50;
372 static constexpr double kLoFHITR = 0.;
373 static constexpr double kUpFHITR = 150.;
374 static constexpr int kBinsNSTA = 12;
375 static constexpr double kLoNSTA = -0.5;
376 static constexpr double kUpNSTA = 11.5;
377 static constexpr int kBinsCHI2NDF = 200;
378 static constexpr double kLoCHI2NDF = 0.;
379 static constexpr double kUpCHI2NDF = 20.;
380
381 // ** Drawing attributes **
382 Color_t fMarkerColor = 1;
383 Style_t fMarkerStyle = 20;
384 Color_t fLineColor = 1;
385 Style_t fLineStyle = 1;
386
388 }; // namespace cbm::ca
389
390} // namespace cbm::ca
Tracking Debugger class (implementation)
Data structure for internal tracking MC-information (header)
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)
CbmQaIO(TString prefixName, std::shared_ptr< ObjList_t > pObjList=nullptr)
Constructor.
Definition CbmQaIO.cxx:19
This class represents a package for tracking-related data.
Definition CaMcData.h:33
Class describes a unified MC-point, used in CA tracking QA analysis.
Definition CaMcPoint.h:33
A container for all external parameters of the CA tracking algorithm.
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.
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
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.
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.
void RegisterMCData(cbm::algo::ca::McData &mcData)
Registers the MC data.
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)
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.
void FillRecoTrack(int iTrkReco, const FitQaPointSet &fitPoints, double weight=1)
Fills histograms with various track information.
double fMCLength
Total length of MC tracks.
cbm::algo::ca::utils::Debugger * fpDebugger
debugger
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.
std::shared_ptr< const ca::Parameters< double > > fpParameters
Pointer to parameters object.
void RegisterParameters(std::shared_ptr< const ca::Parameters< double > > &pParameters)
Registers CA parameters object.
static constexpr int kBinsPT
Number of bins, transverse momentum.
TProfile * fph_rate_clones
Rate of clone tracks / mc.
void SetDebugger(cbm::algo::ca::utils::Debugger *pDebugger)
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].
std::function< bool(const cbm::algo::ca::McTrack &)> fMCTrackCut
Cut function on MC tracks.
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.
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.
cbm::algo::ca::McData * fpMCData
Pointer to MC data object.
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.
void SetMCTrackCut(std::function< bool(const cbm::algo::ca::McTrack &)> cut)
Sets selection cuts on MC tracks.
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.
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.
TProfile2D * fph_clone_thetaMC_phiMC
Efficiency vs. MC theta and MC phi.
TProfile * fph_eff_phiMC
Efficiency vs. MC azimuthal angle.
TProfile2D * fph_eff_thetaMC_phiMC
Efficiency vs. MC theta and MC phi.
void RegisterRecoHits(ca::Vector< cbm::algo::ca::McHitInfo > &vHits)
Register the reconstructed hits container.
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.
ca::Vector< cbm::algo::ca::McHitInfo > * fpvHits
Pointer to vector of reconstructed hits.
TProfile * fph_eff_tyMC
Efficiency vs. MC slope along y-axis.
TProfile * fph_eff_int
Integrated efficiency.
double GetKilledRate() const
Gets killed rate.
std::unique_ptr< TrackFitQa > fpFitQaTrd
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.
TODO: SZh 8.11.2022: add selection of parameterisation.
Definition CaBranch.h:14
TrackParam< double > TrackParamD
FitQaPoint firstHit
First hit point.
bool isTimeFitted
Flag: true - time is fitted, false - time is not fitted.
FitQaPoint lastHit
Last hit point.
cbm::algo::kf::TrackParamD trackParam