Project
Loading...
Searching...
No Matches
MatchITSTPCQC.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
15#include "Framework/InputSpec.h"
19#include <algorithm>
20#include "TGraphAsymmErrors.h"
23
24using namespace o2::globaltracking;
25using namespace o2::mcutils;
27
33
34//_______________________________________________________
35
37{
38
39 for (int i = 0; i < matchType::SIZE; ++i) {
40 // Pt
41 delete mPtNum[i];
42 delete mPtDen[i];
43 delete mFractionITSTPCmatch[i];
44 delete mPtNum_noEta0[i];
45 delete mPtDen_noEta0[i];
46 delete mFractionITSTPCmatch_noEta0[i];
47 delete mPtPhysPrimNum[i];
48 delete mPtPhysPrimDen[i];
49 delete mFractionITSTPCmatchPhysPrim[i];
50
51 // Phi
52 delete mPhiNum[i];
53 delete mPhiDen[i];
54 delete mFractionITSTPCmatchPhi[i];
55 delete mPhiPhysPrimNum[i];
56 delete mPhiPhysPrimDen[i];
57 delete mFractionITSTPCmatchPhiPhysPrim[i];
58 delete mPhiVsPtNum[i];
59 delete mPhiVsPtDen[i];
60 delete mFractionITSTPCmatchPhiVsPt[i];
61
62 // Eta
63 delete mEtaNum[i];
64 delete mEtaDen[i];
65 delete mFractionITSTPCmatchEta[i];
66 delete mEtaPhysPrimNum[i];
67 delete mEtaPhysPrimDen[i];
68 delete mFractionITSTPCmatchEtaPhysPrim[i];
69 delete mEtaVsPtNum[i];
70 delete mEtaVsPtDen[i];
71 delete mFractionITSTPCmatchEtaVsPt[i];
72
73 // Clusters
74 delete mClsVsPtNum[i];
75 delete mClsVsPtDen[i];
76 delete mFractionITSTPCmatchClsVsPt[i];
77
78 // Chi2
79 delete mChi2VsPtNum[i];
80 delete mChi2VsPtDen[i];
81 delete mFractionITSTPCmatchChi2VsPt[i];
82
83 if (mUseTrkPID) { // Vs Tracking PID hypothesis
84 for (int j = 0; j < o2::track::PID::NIDs; ++j) {
85 // Pt
86 delete mPtNumVsTrkPID[i][j];
87 delete mPtDenVsTrkPID[i][j];
88 delete mFractionITSTPCmatchPtVsTrkPID[i][j];
89 // Phi
90 delete mPhiNumVsTrkPID[i][j];
91 delete mPhiDenVsTrkPID[i][j];
92 delete mFractionITSTPCmatchPhiVsTrkPID[i][j];
93 // Eta
94 delete mEtaNumVsTrkPID[i][j];
95 delete mEtaDenVsTrkPID[i][j];
96 delete mFractionITSTPCmatchEtaVsTrkPID[i][j];
97 }
98 }
99
100 // 1/Pt
101 delete m1OverPtNum[i];
102 delete m1OverPtDen[i];
103 delete mFractionITSTPCmatch1OverPt[i];
104 delete m1OverPtPhysPrimNum[i];
105 delete m1OverPtPhysPrimDen[i];
106 delete mFractionITSTPCmatchPhysPrim1OverPt[i];
107 }
108
109 // Residuals
110 delete mResidualPt;
111 delete mResidualPhi;
112 delete mResidualEta;
113 // Others
114 delete mChi2Matching;
115 delete mChi2Refit;
116 delete mTimeResVsPt;
117 delete mDCAr;
118 delete mDCArVsPtNum;
119 delete mDCArVsPtDen;
120 delete mFractionITSTPCmatchDCArVsPt;
121}
122
123//__________________________________________________________
124
126{
127 for (int i = 0; i < matchType::SIZE; ++i) {
128 // Pt
129 mPtNum[i]->Reset();
130 mPtDen[i]->Reset();
131 mPtNum_noEta0[i]->Reset();
132 mPtDen_noEta0[i]->Reset();
133
134 // Phi
135 mPhiNum[i]->Reset();
136 mPhiDen[i]->Reset();
137 mPhiVsPtNum[i]->Reset();
138 mPhiVsPtDen[i]->Reset();
139
140 // Eta
141 mEtaNum[i]->Reset();
142 mEtaDen[i]->Reset();
143 mEtaVsPtNum[i]->Reset();
144 mEtaVsPtDen[i]->Reset();
145
146 // Clusters
147 mClsVsPtNum[i]->Reset();
148 mClsVsPtDen[i]->Reset();
149
150 // Chi2
151 mChi2VsPtNum[i]->Reset();
152 mChi2VsPtDen[i]->Reset();
153
154 // 1/Pt
155 m1OverPtNum[i]->Reset();
156 m1OverPtDen[i]->Reset();
157
158 if (mUseTrkPID) { // Vs Tracking PID hypothesis
159 for (int j = 0; j < o2::track::PID::NIDs; ++j) {
160 // Pt
161 mPtNumVsTrkPID[i][j]->Reset();
162 mPtDenVsTrkPID[i][j]->Reset();
163 // Phi
164 mPhiNumVsTrkPID[i][j]->Reset();
165 mPhiDenVsTrkPID[i][j]->Reset();
166 // Eta
167 mEtaNumVsTrkPID[i][j]->Reset();
168 mEtaDenVsTrkPID[i][j]->Reset();
169 }
170 }
171
172 if (mUseMC) {
173 mPtPhysPrimNum[i]->Reset();
174 mPtPhysPrimDen[i]->Reset();
175
176 mPhiPhysPrimNum[i]->Reset();
177 mPhiPhysPrimDen[i]->Reset();
178
179 mEtaPhysPrimNum[i]->Reset();
180 mEtaPhysPrimDen[i]->Reset();
181
182 m1OverPtPhysPrimNum[i]->Reset();
183 m1OverPtPhysPrimDen[i]->Reset();
184 }
185 }
186
187 // Residuals
188 mResidualPt->Reset();
189 mResidualPhi->Reset();
190 mResidualEta->Reset();
191 // Others
192 mChi2Matching->Reset();
193 mChi2Refit->Reset();
194 mTimeResVsPt->Reset();
195 mDCAr->Reset();
196 mDCArVsPtNum->Reset();
197 mDCArVsPtDen->Reset();
198}
199
200//__________________________________________________________
202{
203 LOGP(debug, "Creating Variable Binning");
204 std::array<std::string, 2> title{"TPC", "ITS"};
205 std::array<std::string, 2> etaSel{"", ", |eta| < 0.9"};
206 std::array<int, 2> maxNCls{156, 7};
207 // log binning for pT
208 const Int_t nbinsPt = 100;
209 const Double_t xminPt = 0.01;
210 const Double_t xmaxPt = 20;
211 Double_t* xbinsPt = new Double_t[nbinsPt + 1];
212 Double_t xlogminPt = TMath::Log10(xminPt);
213 Double_t xlogmaxPt = TMath::Log10(xmaxPt);
214 Double_t dlogxPt = (xlogmaxPt - xlogminPt) / nbinsPt;
215 for (int i = 0; i <= nbinsPt; i++) {
216 Double_t xlogPt = xlogminPt + i * dlogxPt;
217 xbinsPt[i] = TMath::Exp(TMath::Log(10) * xlogPt);
218 }
219
220 LOGP(debug, "Creating Histograms");
221 // Data and MC
222 for (int i = 0; i < matchType::SIZE; ++i) {
223 // Pt
224 mPtNum[i] = new TH1D(Form("mPtNum_%s", title[i].c_str()), Form("Pt distribution of ITSTPC matched tracks, wrt %s tracks %s; Pt [GeV/c]; dNdPt", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 20.f);
225 mPtNum[i]->Sumw2();
226 mPtNum[i]->SetOption("logy");
227 mPtNum[i]->GetYaxis()->SetTitleOffset(1.4);
228 mPtDen[i] = new TH1D(Form("mPtDen_%s", title[i].c_str()), Form("Pt distribution of %s tracks %s; Pt [GeV/c]; dNdPt", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 20.f);
229 mPtDen[i]->Sumw2();
230 mPtDen[i]->SetOption("logy");
231 mPtDen[i]->GetYaxis()->SetTitleOffset(1.4);
232 mFractionITSTPCmatch[i] = new TEfficiency(Form("mFractionITSTPCmatch_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks wrt %s tracks vs Pt %s; Pt [GeV/c]; Eff", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 20.f);
233 mPtNum_noEta0[i] = new TH1D(Form("mPtNum_noEta0_%s", title[i].c_str()), Form("Pt distribution of ITSTPC matched tracks without |eta| < 0.05, wrt %s tracks %s; Pt [GeV/c]; dNdPt", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 20.f);
234 mPtNum_noEta0[i]->Sumw2();
235 mPtNum_noEta0[i]->SetOption("logy");
236 mPtNum_noEta0[i]->GetYaxis()->SetTitleOffset(1.4);
237 mPtDen_noEta0[i] = new TH1D(Form("mPtDen_noEta0_%s", title[i].c_str()), Form("Pt distribution of %s tracks without |eta| < 0.05 %s; Pt [GeV/c]; dNdPt", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 20.f);
238 mPtDen_noEta0[i]->Sumw2();
239 mPtDen_noEta0[i]->SetOption("logy");
240 mPtDen_noEta0[i]->GetYaxis()->SetTitleOffset(1.4);
241 mFractionITSTPCmatch_noEta0[i] = new TEfficiency(Form("mFractionITSTPCmatch_noEta0_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks wrt %s tracks vs Pt without |eta| < 0.05 %s; Pt [GeV/c]; Eff", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 20.f);
242
243 // Phi
244 mPhiNum[i] = new TH1F(Form("mPhiNum_%s", title[i].c_str()), Form("Phi distribution of ITSTPC matched tracks, wrt %s tracks %s; Phi [rad]; dNdPhi", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 2 * TMath::Pi());
245 mPhiNum[i]->Sumw2();
246 mPhiDen[i] = new TH1F(Form("mPhiDen_%s", title[i].c_str()), Form("Phi distribution of %s tracks %s; Phi [rad]; dNdPhi", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 2 * TMath::Pi());
247 mPhiDen[i]->Sumw2();
248 mFractionITSTPCmatchPhi[i] = new TEfficiency(Form("mFractionITSTPCmatchPhi_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks vs Phi wrt %s tracks %s; Phi [rad]; Eff", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 2 * TMath::Pi());
249 mPhiVsPtNum[i] = new TH2F(Form("mPhiVsPtNum_%s", title[i].c_str()), Form("Phi vs Pt distribution of ITSTPC matched tracks wrt %s %s; #it{p}_{T} [GeV#it{c}]; Phi [rad]; dNdPhi", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 20.f, 100, 0.f, 2 * TMath::Pi());
250 mPhiVsPtNum[i]->Sumw2();
251 mPhiVsPtDen[i] = new TH2F(Form("mPhiVsPtDen_%s", title[i].c_str()), Form("Phi vs Pt distribution of %s tracks %s; #it{p}_{T} [GeV#it{c}]; Phi [rad]; dNdPhi", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 20.f, 100, 0.f, 2 * TMath::Pi());
252 mPhiVsPtDen[i]->Sumw2();
253 mFractionITSTPCmatchPhiVsPt[i] = new TEfficiency(Form("mFractionITSTPCmatchPhiVsPt_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks wrt %s tracks %s, Phi vs Pt; #it{p}_{T} [GeV#it{c}]; Phi [rad]; Eff", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 20.f, 100, 0.f, 2 * TMath::Pi());
254
255 // Eta
256 mEtaNum[i] = new TH1F(Form("mEtaNum_%s", title[i].c_str()), Form("Eta distribution of ITSTPC matched tracks, wrt %s tracks; Eta; dNdEta", title[i].c_str()), 100, -2.f, 2.f);
257 mEtaNum[i]->Sumw2();
258 mEtaNum[i]->GetYaxis()->SetTitleOffset(1.4);
259 mEtaDen[i] = new TH1F(Form("mEtaDen_%s", title[i].c_str()), Form("Eta distribution of %s tracks; Eta; dNdEta", title[i].c_str()), 100, -2.f, 2.f);
260 mEtaDen[i]->Sumw2();
261 mEtaDen[i]->GetYaxis()->SetTitleOffset(1.4);
262 mFractionITSTPCmatchEta[i] = new TEfficiency(Form("mFractionITSTPCmatchEta_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks , wrt %s tracks, vs Eta; Eta; Eff", title[i].c_str()), 100, -2.f, 2.f);
263 mEtaVsPtNum[i] = new TH2F(Form("mEtaVsPtNum_%s", title[i].c_str()), Form("Eta vs Pt distribution of ITSTPC matched tracks, wrt %s tracks; #it{p}_{T} [GeV#it{c}]; Eta", title[i].c_str()), 100, 0.f, 20.f, 100, -2.f, 2.f);
264 mEtaVsPtNum[i]->Sumw2();
265 mEtaVsPtDen[i] = new TH2F(Form("mEtaVsPtDen_%s", title[i].c_str()), Form("Eta vs Pt distribution of %s tracks; #it{p}_{T} [GeV#it{c}]; Eta", title[i].c_str()), 100, 0.f, 20.f, 100, -2.f, 2.f);
266 mEtaVsPtDen[i]->Sumw2();
267 mFractionITSTPCmatchEtaVsPt[i] = new TEfficiency(Form("mFractionITSTPCmatchEtaVsPt_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks, wrt %s tracks, Eta vs Pt; #it{p}_{T} [GeV#it{c}]; Eta; Eff", title[i].c_str()), 100, 0.f, 20.f, 100, -2.f, 2.f);
268
269 // Clusters
270 mClsVsPtNum[i] = new TH2F(Form("mClsVsPtNum_%s", title[i].c_str()), Form("#Clusters vs Pt distribution of ITSTPC matched tracks, wrt %s tracks; #it{p}_{T} [GeV#it{c}]; #Clusters", title[i].c_str()), 100, 0.f, 20.f, maxNCls[i], 0, maxNCls[i]);
271 mClsVsPtNum[i]->Sumw2();
272 mClsVsPtDen[i] = new TH2F(Form("mClsVsPtDen_%s", title[i].c_str()), Form("#Clusters vs Pt distribution of %s tracks; #it{p}_{T} [GeV#it{c}]; #Clusters", title[i].c_str()), 100, 0.f, 20.f, maxNCls[i], 0, maxNCls[i]);
273 mClsVsPtDen[i]->Sumw2();
274 mFractionITSTPCmatchClsVsPt[i] = new TEfficiency(Form("mFractionITSTPCmatchClsVsPt_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks, wrt %s tracks, #Clusters vs Pt; #it{p}_{T} [GeV#it{c}]; #Clusters; Eff", title[i].c_str()), 100, 0.f, 20.f, maxNCls[i], 0, maxNCls[i]);
275
276 // Chi2
277 mChi2VsPtNum[i] = new TH2F(Form("mChi2VsPtNum_%s", title[i].c_str()), Form("Chi2 vs Pt distribution of ITSTPC matched tracks, wrt %s tracks; #it{p}_{T} [GeV#it{c}]; Chi2", title[i].c_str()), 100, 0.f, 20.f, 200, 0, 300);
278 mChi2VsPtNum[i]->Sumw2();
279 mChi2VsPtDen[i] = new TH2F(Form("mChi2VsPtDen_%s", title[i].c_str()), Form("Chi2 vs Pt distribution of %s tracks; #it{p}_{T} [GeV#it{c}]; Chi2", title[i].c_str()), 100, 0.f, 20.f, 200, 0, 300);
280 mChi2VsPtDen[i]->Sumw2();
281 mFractionITSTPCmatchChi2VsPt[i] = new TEfficiency(Form("mFractionITSTPCmatchChi2VsPt_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks, wrt %s tracks, Chi2 vs Pt; #it{p}_{T} [GeV#it{c}]; Chi2; Eff", title[i].c_str()), 100, 0.f, 20.f, 200, 0, 300);
282
283 // 1/pt
284 m1OverPtNum[i] = new TH1D(Form("m1OverPtNum_%s", title[i].c_str()), Form("1/Pt distribution of matched tracks, wrt %s tracks %s; 1/Pt [c/GeV]; dNdPt", title[i].c_str(), etaSel[i].c_str()), 100, -20.f, 20.f);
285 m1OverPtNum[i]->Sumw2();
286 m1OverPtDen[i] = new TH1D(Form("m1OverPtDen_%s", title[i].c_str()), Form("1/Pt distribution of %s tracks %s; 1/Pt [c/GeV]; dNdPt", title[i].c_str(), etaSel[i].c_str()), 100, -20.f, 20.f);
287 m1OverPtDen[i]->Sumw2();
288 mFractionITSTPCmatch1OverPt[i] = new TEfficiency(Form("mFractionITSTPCmatch1OverPt_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks vs 1/Pt, wrt %s tracks %s; 1/Pt [c/GeV]; Eff", title[i].c_str(), etaSel[i].c_str()), 100, -20.f, 20.f);
289
290 if (mUseTrkPID) { // Vs Tracking PID hypothesis
291 for (int j = 0; j < o2::track::PID::NIDs; ++j) {
292 // Pt
293 mPtNumVsTrkPID[i][j] = new TH1D(Form("mPtNumVsTrkPID_%s_PID%i", title[i].c_str(), j), Form("Pt distribution of ITSTPC matched tracks, wrt %s tracks %s, TrkPID %i; Pt [GeV/c]; dNdPt", title[i].c_str(), etaSel[i].c_str(), j), 100, 0.f, 20.f);
294 mPtNumVsTrkPID[i][j]->Sumw2();
295 mPtDenVsTrkPID[i][j] = new TH1D(Form("mPtDenVsTrkPID_%s_PID%i", title[i].c_str(), j), Form("Pt distribution of %s tracks %s, TrkPID %i; Pt [GeV/c]; dNdPt", title[i].c_str(), etaSel[i].c_str(), j), 100, 0.f, 20.f);
296 mPtDenVsTrkPID[i][j]->Sumw2();
297 mFractionITSTPCmatchPtVsTrkPID[i][j] = new TEfficiency(Form("mFractionITSTPCmatchPtVsTrkPID_%s_PID%i", title[i].c_str(), j), Form("Fraction of ITSTPC matched tracks wrt %s tracks vs Pt %s, TrkPID %i; Pt [GeV/c]; Eff", title[i].c_str(), etaSel[i].c_str(), j), 100, 0.f, 20.f);
298
299 // Phi
300 mPhiNumVsTrkPID[i][j] = new TH1D(Form("mPhiNumVsTrkPID_%s_PID%i", title[i].c_str(), j), Form("Phi distribution of ITSTPC matched tracks, wrt %s tracks %s, TrkPID %i; Phi [rad]; dNdPhi", title[i].c_str(), etaSel[i].c_str(), j), 100, 0.f, 2 * TMath::Pi());
301 mPhiNumVsTrkPID[i][j]->Sumw2();
302 mPhiDenVsTrkPID[i][j] = new TH1D(Form("mPhiDenVsTrkPID_%s_PID%i", title[i].c_str(), j), Form("Phi distribution of %s tracks %s, TrkPID %i; Phi [rad]; dNdPhi", title[i].c_str(), etaSel[i].c_str(), j), 100, 0.f, 2 * TMath::Pi());
303 mPhiDenVsTrkPID[i][j]->Sumw2();
304 mFractionITSTPCmatchPhiVsTrkPID[i][j] = new TEfficiency(Form("mFractionITSTPCmatchPhiVsTrkPID_%s_PID%i", title[i].c_str(), j), Form("Fraction of ITSTPC matched tracks wrt %s tracks vs Phi %s, TrkPID %i; Phi [rad]; Eff", title[i].c_str(), etaSel[i].c_str(), j), 100, 0.f, 2 * TMath::Pi());
305
306 // Eta
307 mEtaNumVsTrkPID[i][j] = new TH1D(Form("mEtaNumVsTrkPID_%s_PID%i", title[i].c_str(), j), Form("Eta distribution of ITSTPC matched tracks, wrt %s tracks %s, TrkPID %i; Eta; dNdEta", title[i].c_str(), etaSel[i].c_str(), j), 100, -2.f, 2.f);
308 mEtaNumVsTrkPID[i][j]->Sumw2();
309 mEtaDenVsTrkPID[i][j] = new TH1D(Form("mEtaDenVsTrkPID_%s_PID%i", title[i].c_str(), j), Form("Eta distribution of %s tracks %s, TrkPID %i; Eta; dNdEta", title[i].c_str(), etaSel[i].c_str(), j), 100, -2.f, 2.f);
310 mEtaDenVsTrkPID[i][j]->Sumw2();
311 mFractionITSTPCmatchEtaVsTrkPID[i][j] = new TEfficiency(Form("mFractionITSTPCmatchEtaVsTrkPID_%s_PID%i", title[i].c_str(), j), Form("Fraction of ITSTPC matched tracks wrt %s tracks vs Eta %s, TrkPID %i; Eta; Eff", title[i].c_str(), etaSel[i].c_str(), j), 100, -2.f, 2.f);
312 }
313 }
314 }
315
316 mResidualPt = new TH2F("mResidualPt", "Residuals of ITS-TPC matching in #it{p}_{T}; #it{p}_{T}^{ITS-TPC} [GeV/c]; #it{p}_{T}^{ITS-TPC} - #it{p}_{T}^{TPC} [GeV/c]", 100, 0.f, 20.f, 100, -1.f, 1.f);
317 mResidualPhi = new TH2F("mResidualPhi", "Residuals of ITS-TPC matching in #it{#phi}; #it{#phi}^{ITS-TPC} [rad]; #it{#phi}^{ITS-TPC} - #it{#phi}^{TPC} [rad]", 100, 0.f, 2 * TMath::Pi(), 100, -1.f, 1.f);
318 mResidualEta = new TH2F("mResidualEta", "Residuals of ITS-TPC matching in #it{#eta}; #it{#eta}^{ITS-TPC}; #it{#eta}^{ITS-TPC} - #it{#eta}^{TPC}", 100, -2.f, 2.f, 100, -1.f, 1.f);
319 mChi2Matching = new TH1F("mChi2Matching", "Chi2 of matching; chi2", 200, 0, 300);
320 mChi2Matching->SetOption("logy");
321 mChi2Matching->GetYaxis()->SetTitleOffset(1.4);
322 mChi2Refit = new TH1F("mChi2Refit", "Chi2 of refit; chi2", 200, 0, 300);
323 mChi2Refit->SetOption("logy");
324 mChi2Refit->GetYaxis()->SetTitleOffset(1.4);
325 mDCAr = new TH1F("mDCAr", "DCA of TPC tracks; DCAr", 200, -100, 100);
326 mDCArVsPtNum = new TH2F("mDCArVsPtNum", "DCA of TPC tracks Vs Pt Num; #it{p}_{T} [GeV/c]; DCAr", 100, 0, 20., 200, -30, 30);
327 mDCArVsPtNum->Sumw2();
328 mDCArVsPtDen = new TH2F("mDCArVsPtDen", "DCA of TPC tracks Vs Pt Den; #it{p}_{T} [GeV/c]; DCAr", 100, 0, 20., 200, -30, 30);
329 mDCArVsPtDen->Sumw2();
330 mFractionITSTPCmatchDCArVsPt = new TEfficiency("mFractionITSTPCmatchDCArVsPt", "Fraction of ITSTPC matched tracks wrt TPC vs DCAr; #it{p}_{T} [GeV#it{c}]; DCAr; Eff", 100, 0, 20., 200, -30, 30);
331
332 mTimeResVsPt = new TH2F("mTimeResVsPt", "Time resolution vs Pt; Pt [GeV/c]; time res [us]", nbinsPt, xbinsPt, 100, 0.f, 2.f);
333 mTimeResVsPt->SetOption("colz logz logy logx");
334 mTimeResVsPt->GetYaxis()->SetTitleOffset(1.4);
335
336 if (mUseMC) {
337 mcReader.initFromDigitContext("collisioncontext.root");
338
339 for (int i = 0; i < matchType::SIZE; ++i) {
340 mPtPhysPrimNum[i] = new TH1F(Form("mPtPhysPrimNum_%s", title[i].c_str()), Form("Pt distribution of matched tracks (physical primary), wrt %s tracks %s; Pt [GeV/c]; dNdPt", title[i].c_str(), etaSel[i].c_str()), nbinsPt, xbinsPt);
341 mPtPhysPrimNum[i]->Sumw2();
342 mPtPhysPrimDen[i] = new TH1F(Form("mPtPhysPrimDen_%s", title[i].c_str()), Form("Pt distribution of %s tracks (physical primary) %s; Pt [GeV/c]; dNdPt", title[i].c_str(), etaSel[i].c_str()), nbinsPt, xbinsPt);
343 mPtPhysPrimDen[i]->Sumw2();
344 mFractionITSTPCmatchPhiPhysPrim[i] = new TEfficiency(Form("mFractionITSTPCmatchPhiPhysPrim_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks vs Phi (physical primary), wrt %s tracks %s; Phi [rad]; Eff", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 2 * TMath::Pi());
345
346 mEtaPhysPrimNum[i] = new TH1F(Form("mEtaPhysPrimNum_%s", title[i].c_str()), Form("Eta distribution of matched tracks (physical primary), wrt %s tracks; Eta; dNdEta", title[i].c_str()), 100, -2.f, 2.f);
347 mEtaPhysPrimNum[i]->Sumw2();
348 mEtaPhysPrimDen[i] = new TH1F(Form("mEtaPhysPrimDen_%s", title[i].c_str()), Form("Eta distribution of %s tracks (physical primary); Eta; dNdEta", title[i].c_str()), 100, -2.f, 2.f);
349 mEtaPhysPrimDen[i]->Sumw2();
350 mFractionITSTPCmatchEtaPhysPrim[i] = new TEfficiency(Form("mFractionITSTPCmatchEtaPhysPrim_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks vs Eta (physical primary), wrt %s tracks; Eta; Eff", title[i].c_str()), 100, -2.f, 2.f);
351
352 mPhiPhysPrimNum[i] = new TH1F(Form("mPhiPhysPrimNum_%s", title[i].c_str()), Form("Phi distribution of matched tracks (physical primary), wrt %s tracks %s; Phi [rad]; dNdPhi", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 2 * TMath::Pi());
353 mPhiPhysPrimNum[i]->Sumw2();
354 mPhiPhysPrimDen[i] = new TH1F(Form("mPhiPhysPrimDen_%s", title[i].c_str()), Form("Phi distribution of %s tracks (physical primary) %s; Phi [rad]; dNdPhi", title[i].c_str(), etaSel[i].c_str()), 100, 0.f, 2 * TMath::Pi());
355 mPhiPhysPrimDen[i]->Sumw2();
356 mFractionITSTPCmatchPhysPrim[i] = new TEfficiency(Form("mFractionITSTPCmatchPhysPrim_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks vs Pt (physical primary), wrt %s tracks %s; Pt [GeV/c]; Eff", title[i].c_str(), etaSel[i].c_str()), nbinsPt, xbinsPt);
357
358 m1OverPtPhysPrimNum[i] = new TH1D(Form("m1OverPtPhysPrimNum_%s", title[i].c_str()), Form("1/Pt distribution of matched tracks (physical primary), wrt %s tracks %s; 1/Pt [c/GeV]; dNd1/Pt", title[i].c_str(), etaSel[i].c_str()), 100, -20.f, 20.f);
359 m1OverPtPhysPrimNum[i]->Sumw2();
360 m1OverPtPhysPrimDen[i] = new TH1D(Form("m1OverPtPhysPrimDen_%s", title[i].c_str()), Form("1/PtPt distribution of %s tracks (physical primary) %s; 1/Pt [c/GeV]; dNd1/Pt", title[i].c_str(), etaSel[i].c_str()), 100, -20.f, 20.f);
361 m1OverPtPhysPrimDen[i]->Sumw2();
362 mFractionITSTPCmatchPhysPrim1OverPt[i] = new TEfficiency(Form("mFractionITSTPCmatchPhysPrim1OverPt_%s", title[i].c_str()), Form("Fraction of ITSTPC matched tracks vs 1/Pt (physical primary), wrt %s tracks %s; 1/Pt [c/GeV]; Eff", title[i].c_str(), etaSel[i].c_str()), 100, -20.f, 20.f);
363 }
364 }
365
366 return true;
367}
368
369//__________________________________________________________
370
372{
373
374 // initialize data request, if it was not already done
375
376 mSrc &= mAllowedSources;
377
378 if (mSrc[GID::Source::ITSTPC] == 0 || mSrc[GID::Source::TPC] == 0 || mSrc[GID::Source::ITS] == 0) {
379 LOG(fatal) << "We cannot do ITSTPC QC, some sources are missing, check sources in " << mSrc;
380 }
381
382 mDataRequest = std::make_shared<o2::globaltracking::DataRequest>();
383 mDataRequest->requestTracks(mSrc, mUseMC);
384}
385
386//__________________________________________________________
387
389{
390
391 // Getting the B field
392 mBz = o2::base::Propagator::Instance()->getNominalBz();
393
394 static int evCount = 0;
395 mRecoCont.collectData(ctx, *mDataRequest.get());
396 mTPCTracks = mRecoCont.getTPCTracks();
397 mITSTracks = mRecoCont.getITSTracks();
398 mITSTPCTracks = mRecoCont.getTPCITSTracks();
399
400 LOG(debug) << "****** Number of found ITSTPC tracks = " << mITSTPCTracks.size();
401 LOG(debug) << "****** Number of found TPC tracks = " << mTPCTracks.size();
402 LOG(debug) << "****** Number of found ITS tracks = " << mITSTracks.size();
403
404 // cache selection for TPC and ITS tracks
405 std::vector<bool> isTPCTrackSelectedEntry(mTPCTracks.size(), false);
406 std::vector<bool> isITSTrackSelectedEntry(mITSTracks.size(), false);
407 TrackCuts cuts;
408 // ITS track
409 cuts.setMinPtITSCut(mPtITSCut);
410 cuts.setEtaITSCut(mEtaITSCut);
411 cuts.setMinNClustersITS(mMinNClustersITS);
412 cuts.setMaxChi2PerClusterITS(mMaxChi2PerClusterITS);
413 for (auto it = mRequiredITSHits.begin(); it != mRequiredITSHits.end(); it++) {
414 cuts.setRequireHitsInITSLayers((*it).first, (*it).second);
415 }
416 // TPC track
417 cuts.setMinPtTPCCut(mPtTPCCut);
418 cuts.setEtaTPCCut(mEtaTPCCut);
419 cuts.setMinNTPCClustersCut(mNTPCClustersCut);
420 cuts.setMaxDCATPCCut(mDCATPCCut);
421 cuts.setMaxDCATPCCutY(mDCATPCCutY);
422 // ITS-TPC track kinematics
423 cuts.setMinPtCut(mPtCut);
424 cuts.setMaxPtCut(mPtMaxCut);
425 cuts.setEtaCut(-mEtaCut, mEtaCut);
426
427 for (size_t itrk = 0; itrk < mTPCTracks.size(); ++itrk) {
428 auto const& trkTpc = mTPCTracks[itrk];
430 if (cuts.isSelected(id, mRecoCont)) {
431 // NB: same cuts for numerator and denominator tracks of ITS-TPC matching
432 // To change cuts only for numerator, something like o2::dataformats::GlobalTrackID id(itrk, GID::ITSTPC) is necessary
433 isTPCTrackSelectedEntry[itrk] = true;
434 }
435 }
436
437 for (size_t itrk = 0; itrk < mITSTracks.size(); ++itrk) {
438 auto const& trkIts = mITSTracks[itrk];
440 if (cuts.isSelected(id, mRecoCont)) {
441 // NB: same cuts for numerator and denominator tracks of ITS-TPC matching
442 // To change cuts only for numerator, something like o2::dataformats::GlobalTrackID id(itrk, GID::ITSTPC) is necessary
443 isITSTrackSelectedEntry[itrk] = true;
444 }
445 }
446
447 // numerator + eta, chi2...
448 if (mUseMC) {
449 for (int i = 0; i < matchType::SIZE; ++i) {
450 mMapLabels[i].clear();
451 }
452 for (int itrk = 0; itrk < static_cast<int>(mITSTPCTracks.size()); ++itrk) {
453 auto const& trk = mITSTPCTracks[itrk];
454 auto idxTrkTpc = trk.getRefTPC().getIndex();
455 if (trk.getRefITS().getSource() != GID::ITS) {
456 continue;
457 }
458 if (isTPCTrackSelectedEntry[idxTrkTpc] == true) {
459 auto lbl = mRecoCont.getTrackMCLabel({(unsigned int)(itrk), GID::Source::ITSTPC});
460 if (!lbl.isValid()) {
461 continue;
462 }
463 if (mMapLabels[matchType::TPC].find(lbl) == mMapLabels[matchType::TPC].end()) {
464 int source = lbl.getSourceID();
465 int event = lbl.getEventID();
466 const std::vector<o2::MCTrack>& pcontainer = mcReader.getTracks(source, event);
467 const o2::MCTrack& p = pcontainer[lbl.getTrackID()];
468 if (MCTrackNavigator::isPhysicalPrimary(p, pcontainer)) {
469 mMapLabels[matchType::TPC].insert({lbl, {itrk, true}});
470 } else {
471 mMapLabels[matchType::TPC].insert({lbl, {itrk, false}});
472 }
473 } else {
474 // winner (if more tracks have the same label) has the highest pt
475 if (mITSTPCTracks[mMapLabels[matchType::TPC].at(lbl).mIdx].getPt() < trk.getPt()) {
476 mMapLabels[matchType::TPC].at(lbl).mIdx = itrk;
477 }
478 }
479 }
480 auto idxTrkIts = trk.getRefITS().getIndex();
481 if (isITSTrackSelectedEntry[idxTrkIts] == true) {
482 auto lbl = mRecoCont.getTrackMCLabel({(unsigned int)(itrk), GID::Source::ITSTPC});
483 if (!lbl.isValid()) {
484 continue;
485 }
486 if (mMapLabels[matchType::ITS].find(lbl) == mMapLabels[matchType::ITS].end()) {
487 int source = lbl.getSourceID();
488 int event = lbl.getEventID();
489 const std::vector<o2::MCTrack>& pcontainer = mcReader.getTracks(source, event);
490 const o2::MCTrack& p = pcontainer[lbl.getTrackID()];
491 if (MCTrackNavigator::isPhysicalPrimary(p, pcontainer)) {
492 mMapLabels[matchType::ITS].insert({lbl, {itrk, true}});
493 } else {
494 mMapLabels[matchType::ITS].insert({lbl, {itrk, false}});
495 }
496 } else {
497 // winner (if more tracks have the same label) has the highest pt
498 if (mITSTPCTracks[mMapLabels[matchType::ITS].at(lbl).mIdx].getPt() < trk.getPt()) {
499 mMapLabels[matchType::ITS].at(lbl).mIdx = itrk;
500 }
501 }
502 }
503 }
504 LOG(info) << "number of entries in map for nominator (without duplicates) = " << mMapLabels.size();
505 // now we use only the tracks in the map to fill the histograms (--> tracks have passed the
506 // track selection and there are no duplicated tracks wrt the same MC label)
507 for (int i = 0; i < matchType::SIZE; ++i) {
508 for (auto const& el : mMapLabels[i]) {
509 auto const& trk = mITSTPCTracks[el.second.mIdx];
511 bool isEtaITSOk = true;
512 if (i == matchType::TPC) {
513 trkDen = mTPCTracks[trk.getRefTPC()];
514 } else {
515 trkDen = mITSTracks[trk.getRefITS()];
516 if (std::abs(trkDen.getEta()) > 0.9) {
517 // ITS track outside |eta | < 0.9, we don't fill pt, nor phi , nor phi vs pt histos
518 isEtaITSOk = false;
519 }
520 }
521 if (isEtaITSOk) {
522 mPtNum[i]->Fill(trkDen.getPt());
523 if (std::abs(trkDen.getEta()) > 0.05) {
524 mPtNum_noEta0[i]->Fill(trkDen.getPt());
525 }
526 mPhiNum[i]->Fill(trkDen.getPhi());
527 mPhiVsPtNum[i]->Fill(trkDen.getPt(), trkDen.getPhi());
528 m1OverPtNum[i]->Fill(trkDen.getSign() * trkDen.getPtInv());
529 // we fill also the denominator
530 mPtDen[i]->Fill(trkDen.getPt());
531 if (std::abs(trkDen.getEta()) > 0.05) {
532 mPtDen_noEta0[i]->Fill(trkDen.getPt());
533 }
534 mPhiDen[i]->Fill(trkDen.getPhi());
535 mPhiVsPtDen[i]->Fill(trkDen.getPt(), trkDen.getPhi());
536 m1OverPtDen[i]->Fill(trkDen.getSign() * trkDen.getPtInv());
537 if (mUseTrkPID) { // Vs Tracking PID hypothesis
538 mPtNumVsTrkPID[i][trkDen.getPID()]->Fill(trkDen.getPt());
539 mPhiNumVsTrkPID[i][trkDen.getPID()]->Fill(trkDen.getPhi());
540 // we fill also the denominator
541 mPtDenVsTrkPID[i][trkDen.getPID()]->Fill(trkDen.getPt());
542 mPhiDenVsTrkPID[i][trkDen.getPID()]->Fill(trkDen.getPhi());
543 }
544 }
545 mEtaNum[i]->Fill(trkDen.getEta());
546 mEtaVsPtNum[i]->Fill(trkDen.getPt(), trkDen.getEta());
547 // we fill also the denominator
548 mEtaDen[i]->Fill(trkDen.getEta());
549 mEtaVsPtDen[i]->Fill(trkDen.getPt(), trkDen.getEta());
550 if (i == matchType::TPC) {
551 auto tpcTrk = mTPCTracks[trk.getRefTPC()];
552 mClsVsPtNum[i]->Fill(tpcTrk.getPt(), tpcTrk.getNClusters());
553 mChi2VsPtNum[i]->Fill(tpcTrk.getPt(), tpcTrk.getChi2());
554 mClsVsPtDen[i]->Fill(tpcTrk.getPt(), tpcTrk.getNClusters());
555 mChi2VsPtDen[i]->Fill(tpcTrk.getPt(), tpcTrk.getChi2());
557 std::array<float, 2> dca{};
558 if (tpcTrk.propagateParamToDCA(v, mBz, &dca)) {
559 mDCArVsPtNum->Fill(tpcTrk.getPt(), dca[0]);
560 mDCArVsPtDen->Fill(tpcTrk.getPt(), dca[0]);
561 }
562 } else {
563 const auto& itsTrk = mITSTracks[trk.getRefITS()];
564 mClsVsPtNum[i]->Fill(itsTrk.getPt(), itsTrk.getNClusters());
565 mChi2VsPtNum[i]->Fill(itsTrk.getPt(), itsTrk.getChi2());
566 mClsVsPtDen[i]->Fill(itsTrk.getPt(), itsTrk.getNClusters());
567 mChi2VsPtDen[i]->Fill(itsTrk.getPt(), itsTrk.getChi2());
568 }
569 if (mUseTrkPID) { // Vs Tracking PID hypothesis
570 mEtaNumVsTrkPID[i][trkDen.getPID()]->Fill(trkDen.getEta());
571 // we fill also the denominator
572 mEtaDenVsTrkPID[i][trkDen.getPID()]->Fill(trkDen.getEta());
573 }
574 if (el.second.mIsPhysicalPrimary) {
575 if (isEtaITSOk) {
576 mPtPhysPrimNum[i]->Fill(trkDen.getPt());
577 mPhiPhysPrimNum[i]->Fill(trkDen.getPhi());
578 m1OverPtPhysPrimNum[i]->Fill(trkDen.getSign() * trkDen.getPtInv());
579 // we fill also the denominator
580 mPtPhysPrimDen[i]->Fill(trkDen.getPt());
581 mPhiPhysPrimDen[i]->Fill(trkDen.getPhi());
582 m1OverPtPhysPrimDen[i]->Fill(trkDen.getSign() * trkDen.getPtInv());
583 }
584 mEtaPhysPrimNum[i]->Fill(trkDen.getEta());
585 // we fill also the denominator
586 mEtaPhysPrimDen[i]->Fill(trkDen.getEta());
587 }
588 ++mNITSTPCSelectedTracks[i];
589 }
590 }
591 }
592 int iITSTPC = 0;
593 for (auto const& trk : mITSTPCTracks) {
594 if (trk.getRefTPC().getIndex() >= mTPCTracks.size()) {
595 LOG(fatal) << "******************** ATTENTION! for TPC track associated to matched track: idx = " << trk.getRefTPC().getIndex() << ", size of container = " << mTPCTracks.size() << " in TF " << evCount;
596 }
597 std::array<std::string, 2> title{"TPC", "ITS"};
598 for (int i = 0; i < matchType::SIZE; ++i) {
600 int idxTrkRef;
601 bool fillHisto = false;
602 bool isEtaITSOk = true;
603 if (i == matchType::TPC) {
604 trkRef = mTPCTracks[trk.getRefTPC()];
605 idxTrkRef = trk.getRefTPC().getIndex();
606 if (isTPCTrackSelectedEntry[idxTrkRef] == true) {
607 fillHisto = true;
608 ++mNITSTPCSelectedTracks[i];
609 }
610 } else {
611 idxTrkRef = trk.getRefITS().getIndex();
612 if (trk.getRefITS().getSource() == GID::ITSAB) {
613 // do not use afterburner tracks
614 LOG(debug) << "Track (ITS) with id " << idxTrkRef << " for ITSTPC track " << iITSTPC << " is from afterburner";
615 continue;
616 }
617 if (idxTrkRef >= mITSTracks.size()) {
618 LOG(fatal) << "******************** ATTENTION! for ITS track associated to matched track (NOT from AB): idx = " << trk.getRefITS().getIndex() << ", size of container = " << mITSTracks.size() << " in TF " << evCount;
619 }
620 trkRef = mITSTracks[trk.getRefITS()];
621 LOG(debug) << "Checking track (ITS) with id " << idxTrkRef << " for ITSTPC track " << iITSTPC << " and pt = " << trkRef.getPt();
622 if (isITSTrackSelectedEntry[idxTrkRef] == true) {
623 LOG(debug) << "Track was selected (ITS), with id " << idxTrkRef << " for ITSTPC track " << iITSTPC << " , we keep it in the numerator, pt = " << trkRef.getPt();
624 fillHisto = true;
625 ++mNITSTPCSelectedTracks[i];
626 } else {
627 LOG(debug) << "Track was not selected (ITS), with id " << idxTrkRef << " for ITSTPC track " << iITSTPC << " , we don't keep it in the numerator, pt = " << trkRef.getPt();
628 }
629 if (std::abs(trkRef.getEta()) > 0.9) {
630 // ITS track outside |eta | < 0.9, we don't fill pt, nor phi , nor phi vs pt histos
631 isEtaITSOk = false;
632 LOG(debug) << "Track (ITS), with id " << idxTrkRef << " for ITSTPC track " << iITSTPC << " will be discarded when filling pt of phi related histograms, since eta = " << trkRef.getEta() << " , we don't keep it in the numerator, pt = " << trkRef.getPt();
633 }
634 }
635 if (fillHisto == true) {
636 if (!mUseMC) {
637 LOG(debug) << "Filling num (" << title[i] << ") with track with id " << idxTrkRef << " for ITSTPC track " << iITSTPC << " with pt = " << trkRef.getPt();
638 if (isEtaITSOk) {
639 mPtNum[i]->Fill(trkRef.getPt());
640 if (std::abs(trkRef.getEta()) > 0.05) {
641 mPtNum_noEta0[i]->Fill(trkRef.getPt());
642 }
643 mPhiNum[i]->Fill(trkRef.getPhi());
644 if (mUseTrkPID) { // Vs Tracking PID hypothesis
645 mPtNumVsTrkPID[i][trkRef.getPID()]->Fill(trkRef.getPt());
646 mPhiNumVsTrkPID[i][trkRef.getPID()]->Fill(trkRef.getPhi());
647 }
648 mPhiVsPtNum[i]->Fill(trkRef.getPt(), trkRef.getPhi());
649 m1OverPtNum[i]->Fill(trkRef.getSign() * trkRef.getPtInv());
650 }
651 mEtaNum[i]->Fill(trkRef.getEta());
652 if (mUseTrkPID) { // Vs Tracking PID hypothesis
653 mEtaNumVsTrkPID[i][trkRef.getPID()]->Fill(trkRef.getEta());
654 }
655 mEtaVsPtNum[i]->Fill(trkRef.getPt(), trkRef.getEta());
656 if (i == matchType::TPC) {
657 const auto& tpcTrk = mTPCTracks[trk.getRefTPC()];
658 mClsVsPtNum[i]->Fill(tpcTrk.getPt(), tpcTrk.getNClusters());
659 mChi2VsPtNum[i]->Fill(tpcTrk.getPt(), tpcTrk.getChi2());
660 } else {
661 const auto& itsTrk = mITSTracks[trk.getRefITS()];
662 mClsVsPtNum[i]->Fill(itsTrk.getPt(), itsTrk.getNClusters());
663 mChi2VsPtNum[i]->Fill(itsTrk.getPt(), itsTrk.getChi2());
664 }
665 }
666 if (i == matchType::TPC) {
667 mResidualPt->Fill(trk.getPt(), trk.getPt() - trkRef.getPt());
668 mResidualPhi->Fill(trk.getPhi(), trk.getPhi() - trkRef.getPhi());
669 mResidualEta->Fill(trk.getEta(), trk.getEta() - trkRef.getEta());
670 mChi2Matching->Fill(trk.getChi2Match());
671 mChi2Refit->Fill(trk.getChi2Refit());
672 mTimeResVsPt->Fill(trkRef.getPt(), trk.getTimeMUS().getTimeStampError());
674 std::array<float, 2> dca;
675 if (trkRef.propagateParamToDCA(v, mBz, &dca)) {
676 mDCAr->Fill(dca[0]);
677 if (!mUseMC) {
678 mDCArVsPtNum->Fill(trkRef.getPt(), dca[0]);
679 }
680 }
681 LOG(debug) << "*** chi2Matching = " << trk.getChi2Match() << ", chi2refit = " << trk.getChi2Refit() << ", timeResolution = " << trk.getTimeMUS().getTimeStampError();
682 }
683 } else {
684 LOG(debug) << "Not filling num (" << title[i] << ") for ITSTPC track " << iITSTPC << " for track with pt " << trkRef.getPt();
685 }
686 }
687 ++iITSTPC;
688 }
689
690 // now filling the denominator for the efficiency calculation
691 if (mUseMC) {
692 for (int i = 0; i < matchType::SIZE; ++i) {
693 mMapRefLabels[i].clear();
694 }
695 // filling the map where we store for each MC label, the track id of the reconstructed
696 // track with the highest number of TPC clusters
697 for (int itrk = 0; itrk < static_cast<int>(mTPCTracks.size()); ++itrk) {
698 auto const& trk = mTPCTracks[itrk];
699 if (isTPCTrackSelectedEntry[itrk] == true) {
700 auto lbl = mRecoCont.getTrackMCLabel({(unsigned int)(itrk), GID::Source::TPC});
701 if (!lbl.isValid()) {
702 continue;
703 }
704 if (mMapLabels[matchType::TPC].find(lbl) != mMapLabels[matchType::TPC].end()) {
705 // the track was already added to the denominator
706 continue;
707 }
708 if (mMapRefLabels[matchType::TPC].find(lbl) == mMapRefLabels[matchType::TPC].end()) {
709 int source = lbl.getSourceID();
710 int event = lbl.getEventID();
711 const std::vector<o2::MCTrack>& pcontainer = mcReader.getTracks(source, event);
712 const o2::MCTrack& p = pcontainer[lbl.getTrackID()];
713 if (MCTrackNavigator::isPhysicalPrimary(p, pcontainer)) {
714 mMapRefLabels[matchType::TPC].insert({lbl, {itrk, true}});
715 } else {
716 mMapRefLabels[matchType::TPC].insert({lbl, {itrk, false}});
717 }
718 } else {
719 // winner (if more tracks have the same label) has the highest number of TPC clusters
720 if (mTPCTracks[mMapRefLabels[matchType::TPC].at(lbl).mIdx].getNClusters() < trk.getNClusters()) {
721 mMapRefLabels[matchType::TPC].at(lbl).mIdx = itrk;
722 }
723 }
724 }
725 }
726 // same for ITS
727 // filling the map where we store for each MC label, the track id of the reconstructed
728 // track with the highest number of ITS clusters
729 for (int itrk = 0; itrk < static_cast<int>(mITSTracks.size()); ++itrk) {
730 auto const& trk = mITSTracks[itrk];
731 if (isITSTrackSelectedEntry[itrk] == true) {
732 auto lbl = mRecoCont.getTrackMCLabel({(unsigned int)(itrk), GID::Source::ITS});
733 if (!lbl.isValid()) {
734 continue;
735 }
736 if (mMapLabels[matchType::ITS].find(lbl) != mMapLabels[matchType::ITS].end()) {
737 // the track was already added to the denominator
738 continue;
739 }
740 if (mMapRefLabels[matchType::ITS].find(lbl) == mMapRefLabels[matchType::ITS].end()) {
741 int source = lbl.getSourceID();
742 int event = lbl.getEventID();
743 const std::vector<o2::MCTrack>& pcontainer = mcReader.getTracks(source, event);
744 const o2::MCTrack& p = pcontainer[lbl.getTrackID()];
745 if (MCTrackNavigator::isPhysicalPrimary(p, pcontainer)) {
746 mMapRefLabels[matchType::ITS].insert({lbl, {itrk, true}});
747 } else {
748 mMapRefLabels[matchType::ITS].insert({lbl, {itrk, false}});
749 }
750 } else {
751 // winner (if more tracks have the same label) has the highest number of ITS clusters
752 if (mITSTracks[mMapRefLabels[matchType::ITS].at(lbl).mIdx].getNClusters() < trk.getNClusters()) {
753 mMapRefLabels[matchType::ITS].at(lbl).mIdx = itrk;
754 }
755 }
756 }
757 }
758 LOG(info) << "number of entries in map for denominator of TPC tracks (without duplicates) = " << mMapRefLabels[matchType::TPC].size() + mMapLabels[matchType::TPC].size();
759 LOG(info) << "number of entries in map for denominator of ITS tracks (without duplicates) = " << mMapRefLabels[matchType::ITS].size() + mMapLabels[matchType::ITS].size();
760 // now we use only the tracks in the map to fill the histograms (--> tracks have passed the
761 // track selection and there are no duplicated tracks wrt the same MC label)
762 for (auto const& el : mMapRefLabels[matchType::TPC]) {
763 auto const& trk = mTPCTracks[el.second.mIdx];
764 mPtDen[matchType::TPC]->Fill(trk.getPt());
765 if (std::abs(trk.getEta()) > 0.05) {
766 mPtDen_noEta0[matchType::TPC]->Fill(trk.getPt());
767 }
768 mPhiDen[matchType::TPC]->Fill(trk.getPhi());
769 mPhiVsPtDen[matchType::TPC]->Fill(trk.getPt(), trk.getPhi());
770 mEtaDen[matchType::TPC]->Fill(trk.getEta());
771 mEtaVsPtDen[matchType::TPC]->Fill(trk.getPt(), trk.getEta());
772 m1OverPtDen[matchType::TPC]->Fill(trk.getSign() * trk.getPtInv());
773 mClsVsPtDen[matchType::TPC]->Fill(trk.getPt(), trk.getNClusters());
774 mChi2VsPtDen[matchType::TPC]->Fill(trk.getPt(), trk.getChi2());
776 std::array<float, 2> dca{};
777 if (auto trc = trk; trc.propagateParamToDCA(v, mBz, &dca)) {
778 mDCArVsPtDen->Fill(trc.getPt(), dca[0]);
779 }
780 if (el.second.mIsPhysicalPrimary) {
781 mPtPhysPrimDen[matchType::TPC]->Fill(trk.getPt());
782 mPhiPhysPrimDen[matchType::TPC]->Fill(trk.getPhi());
783 mEtaPhysPrimDen[matchType::TPC]->Fill(trk.getEta());
784 m1OverPtPhysPrimDen[matchType::TPC]->Fill(trk.getSign() * trk.getPtInv());
785 }
786 ++mNTPCSelectedTracks;
787 }
788 for (auto const& el : mMapRefLabels[matchType::ITS]) {
789 auto const& trk = mITSTracks[el.second.mIdx];
790 if (std::abs(trk.getEta()) < 0.9) {
791 mPtDen[matchType::ITS]->Fill(trk.getPt());
792 if (std::abs(trk.getEta()) > 0.05) {
793 mPtDen_noEta0[matchType::ITS]->Fill(trk.getPt());
794 }
795 mPhiDen[matchType::ITS]->Fill(trk.getPhi());
796 mPhiVsPtDen[matchType::ITS]->Fill(trk.getPt(), trk.getPhi());
797 m1OverPtDen[matchType::ITS]->Fill(trk.getSign() * trk.getPtInv());
798 }
799 mEtaDen[matchType::ITS]->Fill(trk.getEta());
800 mEtaVsPtDen[matchType::ITS]->Fill(trk.getPt(), trk.getEta());
801 mClsVsPtDen[matchType::ITS]->Fill(trk.getPt(), trk.getNClusters());
802 mChi2VsPtDen[matchType::ITS]->Fill(trk.getPt(), trk.getChi2());
803 if (el.second.mIsPhysicalPrimary) {
804 if (std::abs(trk.getEta()) < 0.9) {
805 mPtPhysPrimDen[matchType::ITS]->Fill(trk.getPt());
806 mPhiPhysPrimDen[matchType::ITS]->Fill(trk.getPhi());
807 m1OverPtPhysPrimDen[matchType::ITS]->Fill(trk.getSign() * trk.getPtInv());
808 }
809 mEtaPhysPrimDen[matchType::ITS]->Fill(trk.getEta());
810 }
811 ++mNITSSelectedTracks;
812 }
813 } else {
814 // if we are in data, we loop over all tracks (no check on the label)
815 for (size_t itrk = 0; itrk < mTPCTracks.size(); ++itrk) {
816 auto const& trk = mTPCTracks[itrk];
817 if (isTPCTrackSelectedEntry[itrk] == true) {
818 LOG(debug) << "Filling den (TPC) with track with pt = " << trk.getPt();
819 mPtDen[matchType::TPC]->Fill(trk.getPt());
820 if (std::abs(trk.getEta()) > 0.05) {
821 mPtDen_noEta0[matchType::TPC]->Fill(trk.getPt());
822 } else {
823 LOG(debug) << "Track (ITS) " << itrk << " with pt = " << trk.getPt() << " and eta = " << trk.getEta() << " not used for den pt, phi, phi vs pt, 1.pt histos";
824 }
825 mPhiDen[matchType::TPC]->Fill(trk.getPhi());
826 mPhiVsPtDen[matchType::TPC]->Fill(trk.getPt(), trk.getPhi());
827 mEtaDen[matchType::TPC]->Fill(trk.getEta());
828 mEtaVsPtDen[matchType::TPC]->Fill(trk.getPt(), trk.getEta());
829 m1OverPtDen[matchType::TPC]->Fill(trk.getSign() * trk.getPtInv());
830 mClsVsPtDen[matchType::TPC]->Fill(trk.getPt(), trk.getNClusters());
831 mChi2VsPtDen[matchType::TPC]->Fill(trk.getPt(), trk.getChi2());
833 std::array<float, 2> dca{};
834 if (auto trc = trk; trc.propagateParamToDCA(v, mBz, &dca)) {
835 mDCArVsPtDen->Fill(trc.getPt(), dca[0]);
836 }
837 ++mNTPCSelectedTracks;
838 }
839 }
840 for (size_t itrk = 0; itrk < mITSTracks.size(); ++itrk) {
841 auto const& trk = mITSTracks[itrk];
842 LOG(debug) << "Checking den for track (ITS) " << itrk << " with pt " << trk.getPt() << " and eta = " << trk.getEta();
843 if (isITSTrackSelectedEntry[itrk] == true) {
844 if (std::abs(trk.getEta()) < 0.9) {
845 LOG(debug) << "Filling den for track (ITS) " << itrk << " with pt = " << trk.getPt() << " and eta = " << trk.getEta();
846 mPtDen[matchType::ITS]->Fill(trk.getPt());
847 if (std::abs(trk.getEta()) > 0.05) {
848 mPtDen_noEta0[matchType::ITS]->Fill(trk.getPt());
849 }
850 mPhiDen[matchType::ITS]->Fill(trk.getPhi());
851 mPhiVsPtDen[matchType::ITS]->Fill(trk.getPt(), trk.getPhi());
852 m1OverPtDen[matchType::ITS]->Fill(trk.getSign() * trk.getPtInv());
853 } else {
854 LOG(debug) << "Track (ITS) " << itrk << " with pt = " << trk.getPt() << " and eta = " << trk.getEta() << " not used for num pt, phi, phi vs pt, 1.pt histos";
855 }
856 mEtaDen[matchType::ITS]->Fill(trk.getEta());
857 mEtaVsPtDen[matchType::ITS]->Fill(trk.getPt(), trk.getEta());
858 mClsVsPtDen[matchType::ITS]->Fill(trk.getPt(), trk.getNClusters());
859 mChi2VsPtDen[matchType::ITS]->Fill(trk.getPt(), trk.getChi2());
860 ++mNITSSelectedTracks;
861 } else {
862 LOG(debug) << "Not filling for this track (ITS) " << itrk << " with pt = " << trk.getPt();
863 }
864 }
865 }
866 evCount++;
867}
868
869//__________________________________________________________
871{
872
873 std::array<std::string, 2> title{"TPC", "ITS"};
874
875 // first we use denominators and nominators to set the TEfficiency; later they are scaled
876
877 // some checks
878 for (int ti = 0; ti < matchType::SIZE; ++ti) {
879 for (int i = 0; i < mPtDen[ti]->GetNbinsX(); ++i) {
880 if (mPtDen[ti]->GetBinContent(i + 1) < mPtNum[ti]->GetBinContent(i + 1)) {
881 LOG(error) << title[ti] << ": bin " << i + 1 << " in [" << mPtNum[ti]->GetBinLowEdge(i + 1) << " , " << mPtNum[ti]->GetBinLowEdge(i + 1) + mPtNum[ti]->GetBinWidth(i + 1) << "]: mPtDen[i] = " << mPtDen[ti]->GetBinContent(i + 1) << ", mPtNum[i] = " << mPtNum[ti]->GetBinContent(i + 1);
882 }
883 }
884 for (int i = 0; i < mPtDen_noEta0[ti]->GetNbinsX(); ++i) {
885 if (mPtDen_noEta0[ti]->GetBinContent(i + 1) < mPtNum_noEta0[ti]->GetBinContent(i + 1)) {
886 LOG(error) << title[ti] << ": bin " << i + 1 << " in [" << mPtNum_noEta0[ti]->GetBinLowEdge(i + 1) << " , " << mPtNum_noEta0[ti]->GetBinLowEdge(i + 1) + mPtNum_noEta0[ti]->GetBinWidth(i + 1) << "]: mPtDen_noEta0[i] = " << mPtDen_noEta0[ti]->GetBinContent(i + 1) << ", mPtNum_noEta0[i] = " << mPtNum_noEta0[ti]->GetBinContent(i + 1);
887 }
888 }
889 for (int i = 0; i < mPhiDen[ti]->GetNbinsX(); ++i) {
890 if (mPhiDen[ti]->GetBinContent(i + 1) < mPhiNum[ti]->GetBinContent(i + 1)) {
891 LOG(error) << title[ti] << ": bin " << i + 1 << " in [" << mPhiNum[ti]->GetBinLowEdge(i + 1) << " , " << mPhiNum[ti]->GetBinLowEdge(i + 1) + mPhiNum[ti]->GetBinWidth(i + 1) << "]: mPhiDen[i] = " << mPhiDen[ti]->GetBinContent(i + 1) << ", mPhiNum[i] = " << mPhiNum[ti]->GetBinContent(i + 1);
892 }
893 }
894 for (int i = 0; i < mEtaDen[ti]->GetNbinsX(); ++i) {
895 if (mEtaDen[ti]->GetBinContent(i + 1) < mEtaNum[ti]->GetBinContent(i + 1)) {
896 LOG(error) << title[ti] << ": bin " << i + 1 << " in [" << mEtaNum[ti]->GetBinLowEdge(i + 1) << " , " << mEtaNum[ti]->GetBinLowEdge(i + 1) + mEtaNum[ti]->GetBinWidth(i + 1) << "]: mEtaDen[i] = " << mEtaDen[ti]->GetBinContent(i + 1) << ", mEtaNum[i] = " << mEtaNum[ti]->GetBinContent(i + 1);
897 }
898 }
899
900 // filling the efficiency
901 setEfficiency(mFractionITSTPCmatch[ti], mPtNum[ti], mPtDen[ti]);
902 setEfficiency(mFractionITSTPCmatch_noEta0[ti], mPtNum_noEta0[ti], mPtDen_noEta0[ti]);
903 setEfficiency(mFractionITSTPCmatchPhi[ti], mPhiNum[ti], mPhiDen[ti]);
904 setEfficiency(mFractionITSTPCmatchEta[ti], mEtaNum[ti], mEtaDen[ti]);
905 setEfficiency(mFractionITSTPCmatchPhiVsPt[ti], mPhiVsPtNum[ti], mPhiVsPtDen[ti], true);
906 setEfficiency(mFractionITSTPCmatchEtaVsPt[ti], mEtaVsPtNum[ti], mEtaVsPtDen[ti], true);
907 setEfficiency(mFractionITSTPCmatch1OverPt[ti], m1OverPtNum[ti], m1OverPtDen[ti]);
908 setEfficiency(mFractionITSTPCmatchClsVsPt[ti], mClsVsPtNum[ti], mClsVsPtDen[ti], true);
909 setEfficiency(mFractionITSTPCmatchChi2VsPt[ti], mChi2VsPtNum[ti], mChi2VsPtDen[ti], true);
910 if (mUseTrkPID) { // Vs Tracking PID hypothesis
911 for (int j = 0; j < o2::track::PID::NIDs; ++j) {
912 setEfficiency(mFractionITSTPCmatchPtVsTrkPID[ti][j], mPtNumVsTrkPID[ti][j], mPtDenVsTrkPID[ti][j]);
913 setEfficiency(mFractionITSTPCmatchPhiVsTrkPID[ti][j], mPhiNumVsTrkPID[ti][j], mPhiDenVsTrkPID[ti][j]);
914 setEfficiency(mFractionITSTPCmatchEtaVsTrkPID[ti][j], mEtaNumVsTrkPID[ti][j], mEtaDenVsTrkPID[ti][j]);
915 }
916 }
917 if (mUseMC) {
918 setEfficiency(mFractionITSTPCmatchPhysPrim[ti], mPtPhysPrimNum[ti], mPtPhysPrimDen[ti]);
919 setEfficiency(mFractionITSTPCmatchPhiPhysPrim[ti], mPhiPhysPrimNum[ti], mPhiPhysPrimDen[ti]);
920 setEfficiency(mFractionITSTPCmatchEtaPhysPrim[ti], mEtaPhysPrimNum[ti], mEtaPhysPrimDen[ti]);
921 setEfficiency(mFractionITSTPCmatchPhysPrim1OverPt[ti], m1OverPtPhysPrimNum[ti], m1OverPtPhysPrimDen[ti]);
922 }
923 }
924 setEfficiency(mFractionITSTPCmatchDCArVsPt, mDCArVsPtNum, mDCArVsPtDen, true);
925 /*
926 mPtTPC->Scale(scaleFactTPC);
927 mPt->Scale(scaleFactITSTPC);
928 mPhiTPC->Scale(scaleFactTPC);
929 mPhi->Scale(scaleFactITSTPC);
930 if (mUseMC) {
931 mPtTPCPhysPrim->Scale(scaleFactTPC);
932 mPtPhysPrim->Scale(scaleFactITSTPC);
933 mPhiTPCPhysPrim->Scale(scaleFactTPC);
934 mPhiPhysPrim->Scale(scaleFactITSTPC);
935 }
936 mEta->Scale(scaleFactITSTPC);
937 mChi2Matching->Scale(scaleFactITSTPC);
938 mChi2Refit->Scale(scaleFactITSTPC);
939 //mTimeResVsPt->Scale(scaleFactITSTPC); // if to few entries, one sees nothing after normalization --> let's not normalize
940 */
941}
942
943//__________________________________________________________
944
945void MatchITSTPCQC::setEfficiency(TEfficiency* eff, TH1* hnum, TH1* hden, bool is2D)
946{
947 if (eff == nullptr) {
948 LOG(fatal) << "Cannot get TEfficiency object ";
949 }
950 if (hnum == nullptr) {
951 LOG(fatal) << "Cannot get numerator histogram for TEfficiency object " << eff->GetName();
952 }
953 if (hden == nullptr) {
954 LOG(fatal) << "Cannot get denominator histogram for TEfficiency object " << eff->GetName();
955 }
956
957 // we need to force to replace the total histogram, otherwise it will compare it to the previous passed one, and it might get an error of inconsistency in the bin contents
958 if constexpr (false) { // checking
959 bool bad{false};
960 LOG(info) << "Setting efficiency " << eff->GetName() << " from " << hnum->GetName() << " and " << hden->GetName();
961 LOG(info) << "Num " << hnum->GetName() << " " << hnum->GetNbinsX() << " " << hnum->GetNbinsY() << " with " << hnum->GetEntries() << " entries";
962 LOG(info) << "Den " << hden->GetName() << " " << hden->GetNbinsX() << " " << hden->GetNbinsY() << " with " << hden->GetEntries() << " entries";
963 if (hnum->GetDimension() != hden->GetDimension()) {
964 LOGP(warning, "Histograms have different dimensions (num={} to den={})", hnum->GetDimension(), hden->GetDimension());
965 bad = true;
966 }
967 if (!TEfficiency::CheckBinning(*hnum, *hden)) {
968 LOGP(warning, "Histograms do not have a compatible binning");
969 bad = true;
970 }
971 if (!is2D) {
972 for (int i = 1; i <= hden->GetNbinsX(); i++) {
973 if (hden->GetBinContent(i) < hnum->GetBinContent(i)) {
974 LOG(warning) << "bin " << i << " den: " << hden->GetBinContent(i) << " < num: " << hnum->GetBinContent(i) << " should be the opposite";
975 bad = true;
976 }
977 }
978 } else {
979 for (int i = 1; i <= hden->GetNbinsX(); i++) {
980 for (int j = 1; j <= hden->GetNbinsY(); j++) {
981 if (hden->GetBinContent(i, j) < hnum->GetBinContent(i, j)) {
982 LOGP(warning, "bin {}/{} -> den: {} < num: {}", i, j, hden->GetBinContent(i, j), hnum->GetBinContent(i, j));
983 bad = true;
984 }
985 }
986 }
987 }
988 if (bad) {
989 return;
990 }
991 }
992 // we need to force to replace the total histogram, otherwise it will compare it to the previous passed one, and it might get an error of inconsistency in the bin contents
993 if (!eff->SetTotalHistogram(*hden, "f")) {
994 LOG(fatal) << "Something went wrong when defining the efficiency denominator " << eff->GetName() << " from " << hnum->GetName();
995 }
996 if (!eff->SetPassedHistogram(*hnum, "")) {
997 LOG(fatal) << "Something went wrong when defining the efficiency numerator " << eff->GetName() << " from " << hnum->GetName();
998 }
999 if (is2D) {
1000 eff->SetTitle(Form("%s;%s;%s;%s", eff->GetTitle(), hnum->GetXaxis()->GetTitle(), hnum->GetYaxis()->GetTitle(), "Efficiency"));
1001 } else {
1002 eff->SetTitle(Form("%s;%s;%s", eff->GetTitle(), hnum->GetXaxis()->GetTitle(), "Efficiency"));
1003 }
1004}
1005
1006//__________________________________________________________
1007
1008void MatchITSTPCQC::getHistos(TObjArray& objar)
1009{
1010
1011 for (int i = 0; i < matchType::SIZE; ++i) {
1012 objar.Add(mPtNum[i]);
1013 objar.Add(mPtDen[i]);
1014 objar.Add(mFractionITSTPCmatch[i]);
1015
1016 objar.Add(mPtNum_noEta0[i]);
1017 objar.Add(mPtDen_noEta0[i]);
1018 objar.Add(mFractionITSTPCmatch_noEta0[i]);
1019
1020 objar.Add(mPtPhysPrimNum[i]);
1021 objar.Add(mPtPhysPrimDen[i]);
1022 objar.Add(mFractionITSTPCmatchPhysPrim[i]);
1023
1024 objar.Add(mPhiNum[i]);
1025 objar.Add(mPhiDen[i]);
1026 objar.Add(mFractionITSTPCmatchPhi[i]);
1027
1028 if (mUseTrkPID) { // Vs Tracking PID hypothesis
1029 for (int j = 0; j < o2::track::PID::NIDs; ++j) {
1030 // Pt
1031 objar.Add(mPtNumVsTrkPID[i][j]);
1032 objar.Add(mPtDenVsTrkPID[i][j]);
1033 objar.Add(mFractionITSTPCmatchPtVsTrkPID[i][j]);
1034 // Phi
1035 objar.Add(mPhiNumVsTrkPID[i][j]);
1036 objar.Add(mPhiDenVsTrkPID[i][j]);
1037 objar.Add(mFractionITSTPCmatchPhiVsTrkPID[i][j]);
1038 // Eta
1039 objar.Add(mEtaNumVsTrkPID[i][j]);
1040 objar.Add(mEtaDenVsTrkPID[i][j]);
1041 objar.Add(mFractionITSTPCmatchEtaVsTrkPID[i][j]);
1042 }
1043 }
1044
1045 objar.Add(mPhiPhysPrimNum[i]);
1046 objar.Add(mPhiPhysPrimDen[i]);
1047 objar.Add(mFractionITSTPCmatchPhiPhysPrim[i]);
1048
1049 objar.Add(mPhiVsPtNum[i]);
1050 objar.Add(mPhiVsPtDen[i]);
1051 objar.Add(mFractionITSTPCmatchPhiVsPt[i]);
1052
1053 objar.Add(mEtaNum[i]);
1054 objar.Add(mEtaDen[i]);
1055 objar.Add(mFractionITSTPCmatchEta[i]);
1056
1057 objar.Add(mEtaPhysPrimNum[i]);
1058 objar.Add(mEtaPhysPrimDen[i]);
1059 objar.Add(mFractionITSTPCmatchEtaPhysPrim[i]);
1060
1061 objar.Add(mEtaVsPtNum[i]);
1062 objar.Add(mEtaVsPtDen[i]);
1063 objar.Add(mFractionITSTPCmatchEtaVsPt[i]);
1064
1065 objar.Add(mClsVsPtNum[i]);
1066 objar.Add(mClsVsPtDen[i]);
1067 objar.Add(mFractionITSTPCmatchClsVsPt[i]);
1068
1069 objar.Add(mChi2VsPtNum[i]);
1070 objar.Add(mChi2VsPtDen[i]);
1071 objar.Add(mFractionITSTPCmatchChi2VsPt[i]);
1072
1073 objar.Add(m1OverPtNum[i]);
1074 objar.Add(m1OverPtDen[i]);
1075 objar.Add(mFractionITSTPCmatch1OverPt[i]);
1076
1077 objar.Add(m1OverPtPhysPrimNum[i]);
1078 objar.Add(m1OverPtPhysPrimDen[i]);
1079 objar.Add(mFractionITSTPCmatchPhysPrim1OverPt[i]);
1080 }
1081 objar.Add(mChi2Matching);
1082 objar.Add(mChi2Refit);
1083 objar.Add(mTimeResVsPt);
1084 objar.Add(mResidualPt);
1085 objar.Add(mResidualPhi);
1086 objar.Add(mResidualEta);
1087 objar.Add(mDCAr);
1088 objar.Add(mDCArVsPtNum);
1089 objar.Add(mDCArVsPtDen);
1090 objar.Add(mFractionITSTPCmatchDCArVsPt);
1091}
Class to perform track cuts.
int32_t i
Helper for geometry and GRP related CCDB requests.
Utility functions for MC particles.
uint32_t j
Definition RawData.h:0
Result of refitting TPC-ITS matched track.
std::ostringstream debug
int getSourceID() const
void setMinPtCut(float value)
ITS+TPC kinematics.
Definition TrackCuts.h:73
void setMinNClustersITS(int32_t value)
Definition TrackCuts.h:59
void setMinNTPCClustersCut(int32_t value)
Definition TrackCuts.h:69
bool isSelected(GID trackIndex, o2::globaltracking::RecoContainer &data)
Definition TrackCuts.h:82
void setMaxPtCut(float value)
Definition TrackCuts.h:74
void setEtaITSCut(float value)
Definition TrackCuts.h:58
void setMaxChi2PerClusterITS(float value)
Definition TrackCuts.h:60
void setEtaTPCCut(float value)
Definition TrackCuts.h:68
void setMaxDCATPCCut(float value)
Definition TrackCuts.h:70
void setMinPtITSCut(float value)
Definition TrackCuts.h:57
void setMinPtTPCCut(float value)
TPC.
Definition TrackCuts.h:67
void setEtaCut(float valueMin, float valueMax)
Definition TrackCuts.h:75
void setMaxDCATPCCutY(float value)
Definition TrackCuts.h:71
void setRequireHitsInITSLayers(int8_t minNRequiredHits, std::set< uint8_t > requiredLayers)
Definition TrackCuts.h:61
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:143
void run(o2::framework::ProcessingContext &ctx)
static bool isPhysicalPrimary(o2::MCTrack const &p, std::vector< o2::MCTrack > const &pcontainer)
Definition MCUtils.cxx:73
bool initFromDigitContext(std::string_view filename)
std::vector< MCTrack > const & getTracks(int source, int event) const
variant returning all tracks for source and event at once
static constexpr ID NIDs
number of defined IDs
Definition PID.h:106
struct _cl_event * event
Definition glcorearb.h:2982
const GLdouble * v
Definition glcorearb.h:832
GLsizei GLsizei GLchar * source
Definition glcorearb.h:798
GLuint id
Definition glcorearb.h:650
TrackParCovF TrackParCov
Definition Track.h:33
o2::MCCompLabel getTrackMCLabel(GTrackID id) const
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"