CbmRoot
Loading...
Searching...
No Matches
CbmCaTrackFitQa.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 "CaDefs.h"
14#include "CaEnumArray.h"
15#include "CbmL1DetectorID.h"
16#include "CbmQaIO.h"
17#include "KfFieldRegion.h"
18#include "KfTrackParam.h"
19
20#include <array>
21
22
23class CbmL1Track;
24class CbmL1HitDebugInfo;
25
26class TFolder;
27class TH1F;
28class TProfile;
29class CbmQaCanvas;
30
31namespace cbm::algo::ca
32{
33 class McPoint;
34}
35
36namespace cbm::ca
37{
38 namespace phys = cbm::algo::ca::constants::phys;
39 namespace constants = cbm::algo::ca::constants;
40 using cbm::algo::ca::McPoint;
41
56
59
62
63
65 template<typename T>
67
72 class TrackFitQa : public CbmQaIO {
73 public:
78 TrackFitQa(const char* pointTag, const char* prefix, std::shared_ptr<ObjList_t> pObjList);
79
81 ~TrackFitQa() = default;
82
84 TrackFitQa(const TrackFitQa&) = delete;
85
88
90 TrackFitQa& operator=(const TrackFitQa&) = delete;
91
94
96 const char* GetTitle() const { return fsTitle; }
97
99 void Init();
100
105 void Fill(const cbm::algo::kf::TrackParamD& trPar, const McPoint& mcPoint, bool bTimeMeasured, double weight = 1);
106
109 void SetParticleMass(double mass) { fMass = mass; }
110
113 void SetTitle(const char* title) { fsTitle = title; }
114
117
120
123
125 CbmQaCanvas* CreateResolutionPlot() { return nullptr; }
126
132 void SetResidualHistoProperties(ETrackParType type, int nBins, double lo, double up);
133
139 void SetPullHistoProperties(ETrackParType type, int nBins, double lo, double up);
140
143
144 // ************************
145 // ** List of histograms **
146 // ************************
147
148
149 TH1F* fph_res_x = nullptr;
150 TH1F* fph_res_y = nullptr;
151 TH1F* fph_res_tx = nullptr;
152 TH1F* fph_res_ty = nullptr;
153 TH1F* fph_res_qp = nullptr;
154 TH1F* fph_res_t = nullptr;
155 TH1F* fph_res_vi = nullptr;
156
157 // ** Pulls **
158 TH1F* fph_pull_x = nullptr;
159 TH1F* fph_pull_y = nullptr;
160 TH1F* fph_pull_tx = nullptr;
161 TH1F* fph_pull_ty = nullptr;
162 TH1F* fph_pull_qp = nullptr;
163 TH1F* fph_pull_t = nullptr;
164 TH1F* fph_pull_vi = nullptr;
165
166 // ** Chi^2 **
167 TH1F* fph_chiSqNdf = nullptr;
168 TH1F* fph_prob = nullptr;
169
170 // ** Resolution profiles **
171 TProfile* fph_res_p_pMC = nullptr;
172 TProfile* fph_res_phi_phiMC = nullptr;
173 TProfile* fph_res_theta_thetaMC = nullptr;
174
175 // **************************
176 // ** Histogram properties **
177 // **************************
178
179 // ** Binning **
180 static constexpr int kCXSIZEPX = 600;
181 static constexpr int kCYSIZEPX = 600;
182
183 private:
184 using FnVal_t = std::function<double()>;
190 void FillResAndPull(ETrackParType type, double recoVal, double recoErr, double trueVal);
191
194
196
197 // ** Distribution properties **
201
205
206 TString fsTitle = "";
207
208 double fMass = constants::phys::MuonMass;
209
211 };
212
213} // namespace cbm::ca
Compile-time constants definition for the CA tracking algorithm.
Implementation of cbm::algo::ca::EnumArray class.
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)
CbmQaIO(TString prefixName, std::shared_ptr< ObjList_t > pObjList=nullptr)
Constructor.
Definition CbmQaIO.cxx:19
Class of arrays, which can be accessed by an enum class entry as an index.
Class describes a unified MC-point, used in CA tracking QA analysis.
Definition CaMcPoint.h:33
static constexpr int kCXSIZEPX
Canvas size along x-axis [px].
~TrackFitQa()=default
Destructor.
TrackParArray_t< double > fvRUp
Upper boundary, residuals.
void SetParticleMass(double mass)
Sets particle mass, used for fitting a track.
TrackParArray_t< TH1F * > fvphResiduals
Residuals for different track parameters.
TrackParArray_t< double > fvRLo
Lower boundary, residuals.
std::function< double()> FnVal_t
TrackFitQa(const char *pointTag, const char *prefix, std::shared_ptr< ObjList_t > pObjList)
Constructor.
TH1F * fph_chiSqNdf
Chi^2/ndf of the track fit.
void Fill(const cbm::algo::kf::TrackParamD &trPar, const McPoint &mcPoint, bool bTimeMeasured, double weight=1)
Fills pull and residual histograms.
TH1F * fph_pull_vi
Pull of inverse speed.
TH1F * fph_res_qp
Residual of q/p [ec/GeV].
TrackParArray_t< int > fvPBins
Number of bins, pulls.
TrackFitQa(TrackFitQa &&)=delete
Move constructor.
TH1F * fph_res_x
Residual of x-coordinate [cm].
TProfile * fph_res_theta_thetaMC
Resolution of polar angle [rad].
TrackParArray_t< TH1F * > fvphPulls
Pulls for different track parameters.
void SetTitle(const char *title)
Sets title, which is to be reflected on legends and titles.
TH1F * fph_pull_ty
Pull of slope along y-axis.
void FitHistograms()
Fit histograms.
TString fsTitle
Title of the point.
void Init()
Initializes histograms.
TrackFitQa & operator=(TrackFitQa &&)=delete
Move assignment operator.
TH1F * fph_pull_t
Pull of time.
static constexpr int kCYSIZEPX
Canvas size along y-axis [px].
TH1F * fph_res_vi
Residual of inverse speed [1/c].
void SetDefaultProperties()
Sets default histogram and track parameter properties.
TH1F * fph_pull_qp
Pull of q/p.
TH1F * fph_pull_tx
Pull of slope along x-axis.
TProfile * fph_res_phi_phiMC
Resolution of azimuthal angle [rad].
CbmQaCanvas * CreateResolutionPlot()
Creates resolutionis plot.
TH1F * fph_pull_y
Pull of y-coordinate.
TH1F * fph_res_ty
Residual of slope along y-axis.
ClassDefNV(TrackFitQa, 0)
Mass of particle.
const char * GetTitle() const
Gets title of fit parameters.
TrackParArray_t< bool > fvbParIgnored
Flag: true - parameter is ignored.
TH1F * fph_res_tx
Residual of slope along x-axis.
TrackParArray_t< int > fvRBins
Number of bins, residuals.
TH1F * fph_res_y
Residual of y-coordinate [cm].
void FillResAndPull(ETrackParType type, double recoVal, double recoErr, double trueVal)
Fills residual and pull for a given track parameter.
TrackFitQa & operator=(const TrackFitQa &)=delete
Copy assignment operator.
CbmQaCanvas * CreatePullPlot()
Creates pulls plot.
TrackFitQa(const TrackFitQa &)=delete
Copy constructor.
TH1F * fph_res_t
Residual of time [ns].
void SetPullHistoProperties(ETrackParType type, int nBins, double lo, double up)
Sets properties for a pull histogram.
TH1F * fph_pull_x
Pull of x-coordinate.
TrackParArray_t< double > fvPUp
Upper boundary, pulls.
TrackParArray_t< double > fvPLo
Lower boundary, pulls.
void SetResidualHistoProperties(ETrackParType type, int nBins, double lo, double up)
Sets properties for a residual histogram.
TProfile * fph_res_p_pMC
Resolution of momentum [GeV/c].
CbmQaCanvas * CreateResidualPlot()
Creates residuals plot.
TH1F * fph_prob
prob() of the track fit
Physics constants.
Definition CaDefs.h:79
TODO: SZh 8.11.2022: add selection of parameterisation.
Definition CaBranch.h:14
TrackParam< double > TrackParamD
ETrackParType
Enumeration for track parameter type.
@ BEGIN
begin of enumeration
@ kQP
charge over total momentum
@ kTY
slope along y-axis
@ kTX
slope along x-axis
ETrackParType operator++(ETrackParType &type)
Prefix increment operator for ETrackParType.
ca::EnumArray< ETrackParType, T > TrackParArray_t
Alias to array, indexed with ETrackParType enumeration.