CbmRoot
Loading...
Searching...
No Matches
CbmCaTrackFitQa.cxx
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#include "CbmCaTrackFitQa.h"
11
12#include "CaMcPoint.h"
13#include "CbmQaCanvas.h"
14#include "CbmQaUtil.h"
15#include "TF1.h"
16#include "TFormula.h"
17#include "TH1.h"
18#include "TProfile.h"
19#include "TString.h"
20
21#include <algorithm>
22
23using namespace cbm::algo::ca;
25
26// *******************************************************
27// ** Implementation of cbm::ca::TrackFitQa functions **
28// *******************************************************
29
32
33using namespace cbm::algo; // TODO: SZh 16.05.2023: Remove this line
34
35// ---------------------------------------------------------------------------------------------------------------------
36//
37TrackFitQa::TrackFitQa(const char* pointTag, const char* prefix, std::shared_ptr<ObjList_t> pObjList)
38 : CbmQaIO(Form("%s_%s", prefix, pointTag), pObjList)
39{
42}
43
44// ---------------------------------------------------------------------------------------------------------------------
45//
47{
48 auto* pCanv = MakeQaObject<CbmQaCanvas>("canv_residuals", " residuals", kCXSIZEPX * 4, kCYSIZEPX * 2);
49 pCanv->Divide2D(7);
50
51
52 for (int iType = static_cast<int>(ETrackParType::BEGIN); iType < static_cast<int>(ETrackParType::END); ++iType) {
53 ETrackParType type = static_cast<ETrackParType>(iType);
54 if (fvbParIgnored[type]) {
55 continue;
56 }
57 pCanv->cd(iType + 1);
58 fvphResiduals[type]->Draw();
59 }
60
61 return pCanv;
62}
63
64// ---------------------------------------------------------------------------------------------------------------------
65//
67{
68 auto* pCanv = MakeQaObject<CbmQaCanvas>(fsPrefix + "_canv_pulls", " pulls", kCXSIZEPX * 4, kCYSIZEPX * 2);
69 pCanv->Divide2D(7);
70
71 for (int iType = static_cast<int>(ETrackParType::BEGIN); iType < static_cast<int>(ETrackParType::END); ++iType) {
72 ETrackParType type = static_cast<ETrackParType>(iType);
73 if (fvbParIgnored[type]) {
74 continue;
75 }
76 auto fit = TF1("fitpull", "[0] * TMath::Exp(TMath::ASinH(-0.5*[3]*((x-[1])/[2])**2)/[3])", -10., 10.);
77 fit.SetParameters(100, 0., 1., .3);
78 fit.SetParLimits(3, 0., 2.);
79 pCanv->cd(iType + 1);
80 fvphPulls[type]->Draw();
82 }
83
84 return pCanv;
85}
86
87// ---------------------------------------------------------------------------------------------------------------------
88void TrackFitQa::FillResAndPull(ETrackParType type, double recoVal, double recoErr, double trueVal)
89{
90 if (fvbParIgnored[type]) {
91 return;
92 }
93 double res = recoVal - trueVal;
94 double pull = res / recoErr;
95 if (type == ETrackParType::kQP) { // for the q/p parameter, the residual is calculated for the momentum
96 res = (recoVal / trueVal - 1.);
97 }
98 fvphResiduals[type]->Fill(res);
99 fvphPulls[type]->Fill(pull);
100}
101
102// ---------------------------------------------------------------------------------------------------------------------
103//
105{
106 auto CreateResidualHisto = [&](ETrackParType t, const char* name, const char* title) {
107 if (fvbParIgnored[t]) {
108 return;
109 }
110 TString sPrefix = (fsTitle.Length() > 0) ? fsTitle + " point residual for " : "residual for ";
111 fvphResiduals[t] = MakeQaObject<TH1F>(name, sPrefix + title, fvRBins[t], fvRLo[t], fvRUp[t]);
112 };
113
114 auto CreatePullHisto = [&](ETrackParType t, const char* name, const char* title) {
115 if (fvbParIgnored[t]) {
116 return;
117 }
118 TString sPrefix = (fsTitle.Length() > 0) ? fsTitle + " point pull for " : "pull for ";
119 fvphPulls[t] = MakeQaObject<TH1F>(name, sPrefix + title, fvPBins[t], fvPLo[t], fvPUp[t]);
120 };
121
122 CreateResidualHisto(ETrackParType::kX, "res_x", "x-coordinate;x^{reco} - x^{MC} [cm]");
123 CreateResidualHisto(ETrackParType::kY, "res_y", "y-coordinate;y^{reco} - y^{MC} [cm]");
124 CreateResidualHisto(ETrackParType::kTX, "res_tx", "slope along x-axis;T_{x}^{reco} - T_{x}^{MC}");
125 CreateResidualHisto(ETrackParType::kTY, "res_ty", "slope along y-axis;T_{y}^{reco} - T_{y}^{MC}");
126 CreateResidualHisto(ETrackParType::kQP, "res_P", "momentum; (p^{reco} - p^{MC})/p^{MC} []");
127 CreateResidualHisto(ETrackParType::kTIME, "res_t", "time;t^{reco} - t^{MC} [ns]");
128 CreateResidualHisto(ETrackParType::kVI, "res_vi", "inverse speed;v^{-1}_{reco} - v^{-1}_{MC} [c^{-1}]");
129
130 CreatePullHisto(ETrackParType::kX, "pull_x", "x-coordinate;(x^{reco} - x^{MC})/#sigma_{x}");
131 CreatePullHisto(ETrackParType::kY, "pull_y", "y-coordinate;(y^{reco} - y^{MC})/#sigma_{y}");
132 CreatePullHisto(ETrackParType::kTX, "pull_tx", "slope along x-axis;(T_{x}^{reco} - T_{x}^{MC})/#sigma_{T_{x}}");
133 CreatePullHisto(ETrackParType::kTY, "pull_ty", "slope along y-axis;(T_{y}^{reco} - T_{y}^{MC})/#sigma_{T_{y}}");
134 CreatePullHisto(ETrackParType::kQP, "pull_qp", "charge over mom.;((q/p)^{reco} - (q/p)^{MC})/#sigma_{q/p}");
135 CreatePullHisto(ETrackParType::kTIME, "pull_t", "time;(t^{reco} - t^{MC})/#sigma_{t}");
136 CreatePullHisto(ETrackParType::kVI, "pull_vi", "inverse speed;(v^{-1}_{reco} - v^{-1}_{MC})/#sigma_{v^{-1}}");
137
138 fph_chiSqNdf = MakeQaObject<TH1F>("chiSqNdf", "ChiSq / Ndf", 100, 0., 10.);
139 fph_prob = MakeQaObject<TH1F>("prob", "probability of the fit", 100, 0., 1.);
140
141 // FIXME: Replace hardcoded parameters with variables
142 fph_res_p_pMC = MakeQaObject<TProfile>("res_p_vs_pMC", "", 100, 0.0, 10.0, -2., 2.);
143 fph_res_phi_phiMC = MakeQaObject<TProfile>("res_phi_vs_phiMC", "", 100, -3.2, 3.2, -2., 2.);
144 fph_res_theta_thetaMC = MakeQaObject<TProfile>("res_theta_vs_phiMC", "", 100, 0., 3.2, -2., 2.);
145
146 // Set histogram titles
147 TString sPrefix = (fsTitle.Length() > 0) ? fsTitle + " point " : "";
148 fph_res_p_pMC->SetTitle(sPrefix + " resolution of momentum;p^{MC} [GeV/c];#delta p [GeV/c]");
149 fph_res_phi_phiMC->SetTitle(sPrefix + " resolution of polar angle;#phi^{MC} [rad];#delta#phi [rad]");
150 fph_res_theta_thetaMC->SetTitle(sPrefix + " resolution of polar angle;#theta^{MC} [rad];#delta#theta [rad]");
151}
152
153// ---------------------------------------------------------------------------------------------------------------------
154//
155void TrackFitQa::Fill(const TrackParamD& trPar, const McPoint& mcPoint, bool bTimeMeasured, double /*weight*/)
156{
157 if (fabs(trPar.GetZ()) - fabs(mcPoint.GetZ()) > 1.e-3) {
158 LOG(error) << "TrackFitQa::Fill: Track z-coordinate " << trPar.GetZ() << " does not match MC point z-coordinate "
159 << mcPoint.GetZ() << ". Skipping track.";
160 return;
161 }
162
163 // ** Time-independent measurements **
164 FillResAndPull(ETrackParType::kX, trPar.GetX(), trPar.GetXError(), mcPoint.GetX());
165 FillResAndPull(ETrackParType::kY, trPar.GetY(), trPar.GetYError(), mcPoint.GetY());
166 FillResAndPull(ETrackParType::kTX, trPar.GetTx(), trPar.GetTxError(), mcPoint.GetTx());
167 FillResAndPull(ETrackParType::kTY, trPar.GetTy(), trPar.GetTyError(), mcPoint.GetTy());
168 FillResAndPull(ETrackParType::kQP, trPar.GetQp(), trPar.GetQpError(), mcPoint.GetQp());
169
170 // Fill chi^2/ndf
171 if (trPar.GetNdf() > 0) {
172 fph_chiSqNdf->Fill(trPar.GetChiSq() / trPar.GetNdf());
173 }
174
175 fph_prob->Fill(TMath::Prob(trPar.GetChiSq(), trPar.GetNdf()));
176
177 // Momentum resolution
178 double recoP = std::fabs(mcPoint.GetCharge() / trPar.GetQp()); // reco mom. (with MC charge)
179 double resP = recoP - mcPoint.GetP(); // residual of total momentum
180
181 // Phi resolution
182 double mcPhi = mcPoint.GetPhi();
183 double recoTx = trPar.Tx();
184 double recoTy = trPar.Ty();
185 double resPhi = atan2(recoTx * cos(mcPhi) + recoTy * sin(mcPhi), -recoTx * cos(mcPhi) + recoTy * sin(mcPhi));
186
187 // Theta resolution
188 double resTheta = trPar.GetTheta() - mcPoint.GetTheta(); // residual of polar angle
189
190 resPhi = std::atan2(std::sin(resPhi), std::cos(resPhi)); // overflow over (-pi, pi] protection
191
192 fph_res_p_pMC->Fill(mcPoint.GetP(), resP);
193 fph_res_phi_phiMC->Fill(mcPoint.GetPhi(), resPhi);
194 fph_res_theta_thetaMC->Fill(mcPoint.GetTheta(), resTheta);
195
196 // ** Time-dependent measurements **
197 if (bTimeMeasured) {
198 FillResAndPull(ETrackParType::kTIME, trPar.GetTime(), trPar.GetTimeError(), mcPoint.GetTime());
199 FillResAndPull(ETrackParType::kVI, trPar.GetVi(), trPar.GetViError(), mcPoint.GetInvSpeed());
200 }
201}
202
203// ---------------------------------------------------------------------------------------------------------------------
204//
206{
207 // ** Residual distribution parameters **
209 fvRLo[ETrackParType::kX] = -4.e-3;
210 fvRUp[ETrackParType::kX] = +4.e-3;
212 fvRLo[ETrackParType::kY] = -4.e-2;
213 fvRUp[ETrackParType::kY] = +4.e-2;
215 fvRLo[ETrackParType::kTX] = -4.e-3;
216 fvRUp[ETrackParType::kTX] = +4.e-3;
218 fvRLo[ETrackParType::kTY] = -4.e-3;
219 fvRUp[ETrackParType::kTY] = +4.e-3;
229
230 // ** Pulls distribution parameters **
232 fvPLo[ETrackParType::kX] = -10.;
233 fvPUp[ETrackParType::kX] = +10.;
235 fvPLo[ETrackParType::kY] = -10.;
236 fvPUp[ETrackParType::kY] = +10.;
238 fvPLo[ETrackParType::kTX] = -10.;
239 fvPUp[ETrackParType::kTX] = +10.;
241 fvPLo[ETrackParType::kTY] = -10.;
242 fvPUp[ETrackParType::kTY] = +10.;
244 fvPLo[ETrackParType::kQP] = -10.;
245 fvPUp[ETrackParType::kQP] = +10.;
250 fvPLo[ETrackParType::kVI] = -10.;
251 fvPUp[ETrackParType::kVI] = +10.;
252}
253
254// ---------------------------------------------------------------------------------------------------------------------
255//
256void TrackFitQa::SetResidualHistoProperties(ETrackParType type, int nBins, double lo, double up)
257{
258 fvRBins[type] = nBins;
259 fvRLo[type] = lo;
260 fvRUp[type] = up;
261}
262
263// ---------------------------------------------------------------------------------------------------------------------
264//
265void TrackFitQa::SetPullHistoProperties(ETrackParType type, int nBins, double lo, double up)
266{
267 fvPBins[type] = nBins;
268 fvPLo[type] = lo;
269 fvPUp[type] = up;
270}
Internal class describing a MC point for CA tracking QA and performance (header)
QA submodule for track fit results (residuals and puls) at selected z-coordinate (header)
Definition of the CbmQaCanvas class.
Useful utilities for CBM QA tasks.
friend fvec cos(const fvec &a)
friend fvec sin(const fvec &a)
T * MakeQaObject(TString sName, TString sTitle, Args... args)
Definition CbmQaIO.h:260
TString fsPrefix
Unique prefix for all writeable root.
Definition CbmQaIO.h:158
CbmQaIO(TString prefixName, std::shared_ptr< ObjList_t > pObjList=nullptr)
Constructor.
Definition CbmQaIO.cxx:19
@ kSAMEDIR
Objects of different type will be stored to root directory.
Definition CbmQaIO.h:49
EStoringMode fStoringMode
Objects storing mode.
Definition CbmQaIO.h:160
TrackFitQa(const char *pointTag, const char *prefix, std::shared_ptr< ObjList_t > pObjList)
Constructor.
Class describes a unified MC-point, used in CA tracking QA analysis.
Definition CaMcPoint.h:33
int GetInvSpeed() const
Gets inverse speed at reference z of station [ns/cm].
Definition CaMcPoint.h:90
double GetPhi() const
Gets track azimuthal angle at reference z of station [rad].
Definition CaMcPoint.h:126
double GetTime() const
Gets time [ns].
Definition CaMcPoint.h:207
double GetP() const
Gets track momentum absolute value at reference z of station [GeV/c].
Definition CaMcPoint.h:120
double GetY() const
Gets y coordinate at reference z of station [cm].
Definition CaMcPoint.h:240
double GetQp() const
Gets track charge over momentum at reference z of station [ec/GeV].
Definition CaMcPoint.h:189
double GetCharge() const
Gets charge of the particle [e].
Definition CaMcPoint.h:64
double GetX() const
Gets x coordinate at reference z of station [cm].
Definition CaMcPoint.h:231
double GetTx() const
Gets slope along x-axis at reference z of station.
Definition CaMcPoint.h:213
double GetTheta() const
Gets polar angle at reference z of station [rad].
Definition CaMcPoint.h:198
double GetZ() const
Gets z coordinate at reference z of station [cm].
Definition CaMcPoint.h:249
double GetTy() const
Gets slope along x-axis at reference z of station.
Definition CaMcPoint.h:222
T GetTxError() const
Gets error of slope along x-axis.
T GetVi() const
Gets inverse velocity [ns/cm] in downstream direction.
T GetXError() const
Gets x position error [cm].
T GetYError() const
Gets y position error [cm].
T GetTy() const
Gets slope along y-axis.
T GetTimeError() const
Gets time error [ns].
T GetZ() const
Gets z position [cm].
T Tx() const
Gets slope along x-axis.
T GetTime() const
Gets time [ns].
T GetTyError() const
Gets error of slope along y-axis.
T GetTheta() const
Gets polar angle [rad].
T GetTx() const
Gets slope along x-axis.
T GetViError() const
Gets inverse velocity error [ns/cm].
T GetQp() const
Gets charge over momentum [ec/GeV].
T GetNdf() const
Gets NDF of track fit model.
T GetQpError() const
Gets error of charge over momentum [ec/GeV].
T GetY() const
Gets y position [cm].
T GetChiSq() const
Gets Chi-square of track fit model.
T Ty() const
Gets slope along y-axis.
T GetX() const
Gets x position [cm].
Set of histograms to monitor track parameters.
static constexpr int kCXSIZEPX
Canvas size along x-axis [px].
TrackParArray_t< double > fvRUp
Upper boundary, residuals.
TrackParArray_t< TH1F * > fvphResiduals
Residuals for different track parameters.
TrackParArray_t< double > fvRLo
Lower boundary, residuals.
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.
TrackParArray_t< int > fvPBins
Number of bins, pulls.
TProfile * fph_res_theta_thetaMC
Resolution of polar angle [rad].
TrackParArray_t< TH1F * > fvphPulls
Pulls for different track parameters.
TString fsTitle
Title of the point.
void Init()
Initializes histograms.
static constexpr int kCYSIZEPX
Canvas size along y-axis [px].
void SetDefaultProperties()
Sets default histogram and track parameter properties.
TProfile * fph_res_phi_phiMC
Resolution of azimuthal angle [rad].
TrackParArray_t< bool > fvbParIgnored
Flag: true - parameter is ignored.
TrackParArray_t< int > fvRBins
Number of bins, residuals.
void FillResAndPull(ETrackParType type, double recoVal, double recoErr, double trueVal)
Fills residual and pull for a given track parameter.
CbmQaCanvas * CreatePullPlot()
Creates pulls plot.
void SetPullHistoProperties(ETrackParType type, int nBins, double lo, double up)
Sets properties for a pull histogram.
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
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
@ END
end of enumeration
@ kTY
slope along y-axis
@ kTX
slope along x-axis
std::tuple< double, double > FitKaniadakisGaussian(TH1 *pHist)
Fit a histogram with Kaniadakis Gaussian Distribution.
Definition CbmQaUtil.cxx:88