Project
Loading...
Searching...
No Matches
MFTAssessment.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
13#include "Framework/InputSpec.h"
14#include "MFTBase/GeometryTGeo.h"
15#include "MathUtils/Utils.h"
17#include <TPaveText.h>
18#include <TLegend.h>
19#include <TStyle.h>
20#include <TFile.h>
21
22using namespace o2::mft;
25
26//__________________________________________________________
27void MFTAssessment::init(bool finalizeAnalysis)
28{
29 mLastTrackType = mUseMC ? kNumberOfTrackTypes : kReco + 1;
30 mFinalizeAnalysis = finalizeAnalysis;
32 mUnusedChips.fill(true);
33}
34
35//__________________________________________________________
37{
38
39 mTrackNumberOfClusters->Reset();
40 mCATrackNumberOfClusters->Reset();
41 mLTFTrackNumberOfClusters->Reset();
42 mTrackInvQPt->Reset();
43 mTrackChi2->Reset();
44 mTrackCharge->Reset();
45 mTrackPhi->Reset();
46 mPositiveTrackPhi->Reset();
47 mNegativeTrackPhi->Reset();
48 mTrackEta->Reset();
49 mTrackChi2pT->Reset();
50
51 mMFTClsZ->Reset();
52 mMFTClsOfTracksZ->Reset();
53
54 mUnusedChips.fill(true);
55 mNumberTFs = 0;
56
57 for (auto nMFTLayer = 0; nMFTLayer < 10; nMFTLayer++) {
58 mMFTClsXYinLayer[nMFTLayer]->Reset();
59 mMFTClsRinLayer[nMFTLayer]->Reset();
60 mMFTClsOfTracksXYinLayer[nMFTLayer]->Reset();
61 }
62
63 for (auto nMFTDisk = 0; nMFTDisk < 5; nMFTDisk++) {
64 mMFTClsXYRedundantInDisk[nMFTDisk]->Reset();
65 }
66
67 for (auto minNClusters : sMinNClustersList) {
68 auto nHisto = minNClusters - sMinNClustersList[0];
69 mTrackEtaNCls[nHisto]->Reset();
70 mTrackPhiNCls[nHisto]->Reset();
71 mTrackXYNCls[nHisto]->Reset();
72 mTrackEtaPhiNCls[nHisto]->Reset();
73 }
74 mCATrackEta->Reset();
75 mLTFTrackEta->Reset();
76 mTrackCotl->Reset();
77
78 mTrackROFNEntries->Reset();
79 mClusterROFNEntries->Reset();
80
81 mClusterSensorIndex->Reset();
82 mClusterPatternIndex->Reset();
83
84 for (int trackType = 0; trackType < mLastTrackType; trackType++) {
85 mHistPhiVsEta[trackType]->Reset();
86 mHistPtVsEta[trackType]->Reset();
87 mHistPhiVsPt[trackType]->Reset();
88 if (trackType != kGen || trackType != kTrackable || trackType != kRecoTrueMC) {
89 mHistTrackChi2[trackType]->Reset();
90 }
91 if (trackType != kReco || trackType != kRecoFake) {
92 mHistZvtxVsEta[trackType]->Reset();
93 }
94 if (trackType == kGen || trackType == kTrackable || trackType == kRecoTrueMC) {
95 mHistRVsZ[trackType]->Reset();
96 mHistIsPrimary[trackType]->Reset();
97 }
98 }
99
100 if (mUseMC) {
101 mMFTTrackables.clear();
102 mTrueTracksMap.clear();
103 mFakeTracksVec.clear();
104 mTrackableTracksMap.clear();
105
106 mHistPhiRecVsPhiGen->Reset();
107 mHistEtaRecVsEtaGen->Reset();
108
109 auto h = mChargeMatchEff->GetCopyTotalHisto();
110 h->Reset();
111 mChargeMatchEff->SetTotalHistogram(*h, "");
112 mChargeMatchEff->SetPassedHistogram(*h, "");
113
114 for (auto& h : mTH3Histos) {
115 h->Reset();
116 }
117 }
118}
119
120//__________________________________________________________
122{
123 auto MaxClusterROFSize = 25000;
124 auto MaxTrackROFSize = 5000;
125 auto ROFLengthInBC = 198;
126 auto ROFsPerOrbit = o2::constants::lhc::LHCMaxBunches / ROFLengthInBC;
127 auto MaxDuration = 60000.f;
128 auto TimeBinSize = 0.1f;
129 auto NofTimeBins = static_cast<int>(MaxDuration / TimeBinSize);
130
131 // Creating data-only histos
132 mTrackNumberOfClusters = std::make_unique<TH1F>("mMFTTrackNumberOfClusters",
133 "Number Of Clusters Per Track; # clusters; # entries", 10, 0.5, 10.5);
134
135 mCATrackNumberOfClusters = std::make_unique<TH1F>("mMFTCATrackNumberOfClusters",
136 "Number Of Clusters Per CA Track; # clusters; # tracks", 10, 0.5, 10.5);
137
138 mLTFTrackNumberOfClusters = std::make_unique<TH1F>("mMFTLTFTrackNumberOfClusters",
139 "Number Of Clusters Per LTF Track; # clusters; # entries", 10, 0.5, 10.5);
140
141 mTrackInvQPt = std::make_unique<TH1F>("mMFTTrackInvQPt", "Track q/p_{T}; q/p_{T} [1/GeV]; # entries", 200, -10, 10);
142
143 mTrackChi2 = std::make_unique<TH1F>("mMFTTrackChi2", "Track #chi^{2}/NDF; #chi^{2}/NDF; # entries", 210, -0.5, 20.5);
144
145 mTrackCharge = std::make_unique<TH1F>("mMFTTrackCharge", "Track Charge; q; # entries", 3, -1.5, 1.5);
146
147 mTrackPhi = std::make_unique<TH1F>("mMFTTrackPhi", "Track #phi; #phi; # entries", 100, -3.2, 3.2);
148
149 mPositiveTrackPhi = std::make_unique<TH1F>("mMFTPositiveTrackPhi", "Positive Track #phi; #phi; # entries", 100, -3.2, 3.2);
150
151 mNegativeTrackPhi = std::make_unique<TH1F>("mMFTNegativeTrackPhi", "Negative Track #phi; #phi; # entries", 100, -3.2, 3.2);
152
153 mTrackEta = std::make_unique<TH1F>("mMFTTrackEta", "Track #eta; #eta; # entries", 50, -4, -2);
154
155 mTrackChi2pT = std::make_unique<TH2F>("mMFTTrackChi2pT", "Track #chi^{2}/NDF vs p_{T}; #it{p}_{T} (GeV/c); #chi^{2}/NDF", 210, -0.5, 20.5, 210, -0.5, 20.5);
156
157 //----------------------------------------------------------------------------
158
159 mMFTClsZ = std::make_unique<TH1F>("mMFTClsZ", "Z of all clusters; Z (cm); # entries", 400, -80, -40);
160
161 mMFTClsOfTracksZ = std::make_unique<TH1F>("mMFTClsOfTracksZ", "Z of clusters belonging to MFT tracks; Z (cm); # entries", 400, -80, -40);
162
163 for (auto nMFTLayer = 0; nMFTLayer < 10; nMFTLayer++) {
164 mMFTClsXYinLayer[nMFTLayer] = std::make_unique<TH2F>(Form("mMFTClsXYinLayer%d", nMFTLayer), Form("Cluster Position in Layer %d; x (cm); y (cm)", nMFTLayer), 400, -20, 20, 400, -20, 20);
165 mMFTClsRinLayer[nMFTLayer] = std::make_unique<TH1F>(Form("mMFTClsRinLayer%d", nMFTLayer), Form("Cluster Radial Position in Layer %d; r (cm); # entries", nMFTLayer), 400, 0, 20);
166 mMFTClsOfTracksXYinLayer[nMFTLayer] = std::make_unique<TH2F>(Form("mMFTClsOfTracksXYinLayer%d", nMFTLayer), Form("Cluster (of MFT tracks) Position in Layer %d; x (cm); y (cm)", nMFTLayer), 400, -20, 20, 400, -20, 20);
167 }
168
169 for (auto nMFTDisk = 0; nMFTDisk < 5; nMFTDisk++) {
170 mMFTClsXYRedundantInDisk[nMFTDisk] = std::make_unique<TH2F>(Form("mMFTClsXYRedundantInDisk%d", nMFTDisk), Form("Redondant Cluster Position in disk %d; x (cm); y (cm)", nMFTDisk), 400, -20, 20, 400, -20, 20);
171 }
172
173 for (auto minNClusters : sMinNClustersList) {
174 auto nHisto = minNClusters - sMinNClustersList[0];
175 mTrackEtaNCls[nHisto] = std::make_unique<TH1F>(Form("mMFTTrackEta_%d_MinClusters", minNClusters), Form("Track #eta (NCls >= %d); #eta; # entries", minNClusters), 50, -4, -2);
176
177 mTrackPhiNCls[nHisto] = std::make_unique<TH1F>(Form("mMFTTrackPhi_%d_MinClusters", minNClusters), Form("Track #phi (NCls >= %d); #phi; # entries", minNClusters), 100, -3.2, 3.2);
178
179 mTrackXYNCls[nHisto] = std::make_unique<TH2F>(Form("mMFTTrackXY_%d_MinClusters", minNClusters), Form("Track Position (NCls >= %d); x; y", minNClusters), 320, -16, 16, 320, -16, 16);
180 mTrackXYNCls[nHisto]->SetOption("COLZ");
181
182 mTrackEtaPhiNCls[nHisto] = std::make_unique<TH2F>(Form("mMFTTrackEtaPhi_%d_MinClusters", minNClusters), Form("Track #eta , #phi (NCls >= %d); #eta; #phi", minNClusters), 50, -4, -2, 100, -3.2, 3.2);
183 mTrackEtaPhiNCls[nHisto]->SetOption("COLZ");
184 }
185
186 mCATrackEta = std::make_unique<TH1F>("mMFTCATrackEta", "CA Track #eta; #eta; # entries", 50, -4, -2);
187
188 mLTFTrackEta = std::make_unique<TH1F>("mMFTLTFTrackEta", "LTF Track #eta; #eta; # entries", 50, -4, -2);
189
190 mTrackCotl = std::make_unique<TH1F>("mMFTTrackCotl", "Track cot #lambda; cot #lambda; # entries", 100, -0.25, 0);
191
192 mClusterROFNEntries = std::make_unique<TH1F>("mMFTClustersROFSize", "MFT Cluster ROFs size; ROF Size; # entries", MaxClusterROFSize, 0, MaxClusterROFSize);
193
194 mTrackROFNEntries = std::make_unique<TH1F>("mMFTTrackROFSize", "MFT Track ROFs size; ROF Size; # entries", MaxTrackROFSize, 0, MaxTrackROFSize);
195
196 mTracksBC = std::make_unique<TH1F>("mMFTTracksBC", "Tracks per BC (sum over orbits); BCid; # entries", ROFsPerOrbit, 0, o2::constants::lhc::LHCMaxBunches);
197 mTracksBC->SetMinimum(0.1);
198
199 mNOfTracksTime = std::make_unique<TH1F>("mNOfTracksTime", "Number of tracks per time bin; time (s); # entries", NofTimeBins, 0, MaxDuration);
200 mNOfTracksTime->SetMinimum(0.1);
201
202 mNOfClustersTime = std::make_unique<TH1F>("mNOfClustersTime", "Number of clusters per time bin; time (s); # entries", NofTimeBins, 0, MaxDuration);
203 mNOfClustersTime->SetMinimum(0.1);
204
205 mClusterSensorIndex = std::make_unique<TH1F>("mMFTClusterSensorIndex", "Chip Cluster Occupancy;Chip ID;#Entries", 936, -0.5, 935.5);
206
207 mClusterPatternIndex = std::make_unique<TH1F>("mMFTClusterPatternIndex", "Cluster Pattern ID;Pattern ID;#Entries", 300, -0.5, 299.5);
208
209 for (int trackType = 0; trackType < mLastTrackType; trackType++) {
210 // mHistPhiVsEta
211 mHistPhiVsEta[trackType] = std::make_unique<TH2F>((std::string("mHistPhiVsEta") + mNameOfTrackTypes[trackType]).c_str(), (std::string("Phi Vs Eta of ") + mNameOfTrackTypes[trackType]).c_str(), 35, -4.5, -1, 24, 0, 2 * TMath::Pi());
212 mHistPhiVsEta[trackType]->SetXTitle((std::string("#eta of ") + mNameOfTrackTypes[trackType]).c_str());
213 mHistPhiVsEta[trackType]->SetYTitle((std::string("#phi of ") + mNameOfTrackTypes[trackType]).c_str());
214 mHistPhiVsEta[trackType]->Sumw2();
215 mHistPhiVsEta[trackType]->SetOption("COLZ");
216
217 // mHistPtVsEta
218 mHistPtVsEta[trackType] = std::make_unique<TH2F>((std::string("mHistPtVsEta") + mNameOfTrackTypes[trackType]).c_str(), (std::string("Pt Vs Eta of ") + mNameOfTrackTypes[trackType]).c_str(), 35, -4.5, -1, 40, 0., 10.);
219 mHistPtVsEta[trackType]->SetXTitle((std::string("#eta of ") + mNameOfTrackTypes[trackType]).c_str());
220 mHistPtVsEta[trackType]->SetYTitle((std::string("p_{T} (GeV/c) of ") + mNameOfTrackTypes[trackType]).c_str());
221 mHistPtVsEta[trackType]->Sumw2();
222 mHistPtVsEta[trackType]->SetOption("COLZ");
223
224 // mHistPhiVsPt
225 mHistPhiVsPt[trackType] = std::make_unique<TH2F>((std::string("mHistPhiVsPt") + mNameOfTrackTypes[trackType]).c_str(), (std::string("Phi Vs Pt of ") + mNameOfTrackTypes[trackType]).c_str(), 40, 0., 10., 24, 0, 2 * TMath::Pi());
226 mHistPhiVsPt[trackType]->SetXTitle((std::string("p_{T} (GeV/c) of ") + mNameOfTrackTypes[trackType]).c_str());
227 mHistPhiVsPt[trackType]->SetYTitle((std::string("#phi of ") + mNameOfTrackTypes[trackType]).c_str());
228 mHistPhiVsPt[trackType]->Sumw2();
229 mHistPhiVsPt[trackType]->SetOption("COLZ");
230
231 if (trackType != kGen || trackType != kTrackable || trackType != kRecoTrueMC) {
232 // mHistTrackChi2
233 mHistTrackChi2[trackType] = std::make_unique<TH1F>((std::string("mHistTrackChi2") + mNameOfTrackTypes[trackType]).c_str(), (std::string("Track #chi^{2}/NDF of ") + mNameOfTrackTypes[trackType]).c_str(), 210, -0.5, 20.5);
234 mHistTrackChi2[trackType]->SetXTitle("#chi^{2}/NDF");
235 mHistTrackChi2[trackType]->SetYTitle("Entries");
236 }
237
238 if (trackType == kGen || trackType == kTrackable || trackType == kRecoTrue || trackType == kRecoTrueMC) {
239 // mHistZvtxVsEta
240 mHistZvtxVsEta[trackType] = std::make_unique<TH2F>((std::string("mHistZvtxVsEta") + mNameOfTrackTypes[trackType]).c_str(), (std::string("Z_{vtx} Vs Eta of ") + mNameOfTrackTypes[trackType]).c_str(), 60, -15, 15, 70, -4.5, -1);
241 mHistZvtxVsEta[trackType]->SetXTitle((std::string("z_{vtx} (cm) of ") + mNameOfTrackTypes[trackType]).c_str());
242 mHistZvtxVsEta[trackType]->SetYTitle((std::string("#eta of ") + mNameOfTrackTypes[trackType]).c_str());
243 mHistZvtxVsEta[trackType]->Sumw2();
244 mHistZvtxVsEta[trackType]->SetOption("COLZ");
245 }
246 // mHistRVsZ
247 if (trackType != kReco || trackType != kRecoFake) {
248 mHistRVsZ[trackType] = std::make_unique<TH2F>((std::string("mHistRVsZ") + mNameOfTrackTypes[trackType]).c_str(), (std::string("R Vs Z of ") + mNameOfTrackTypes[trackType]).c_str(), 400, -80., 20., 400, 0., 80.);
249 mHistRVsZ[trackType]->SetXTitle((std::string("z (cm) origin of ") + mNameOfTrackTypes[trackType]).c_str());
250 mHistRVsZ[trackType]->SetYTitle((std::string("R (cm) radius of origin of ") + mNameOfTrackTypes[trackType]).c_str());
251 mHistRVsZ[trackType]->Sumw2();
252 mHistRVsZ[trackType]->SetOption("COLZ");
253
254 mHistIsPrimary[trackType] = std::make_unique<TH1F>((std::string("mHistIsPrimary") + mNameOfTrackTypes[trackType]).c_str(), (std::string("Is Primary - ") + mNameOfTrackTypes[trackType]).c_str(), 2, -.5, 1.5);
255 mHistIsPrimary[trackType]->SetXTitle((mNameOfTrackTypes[trackType] + std::string(" primaries")).c_str());
256 mHistIsPrimary[trackType]->SetYTitle("Entries");
257 mHistIsPrimary[trackType]->Sumw2();
258 mHistIsPrimary[trackType]->SetOption("COLZ");
259 }
260 }
261
262 // Creating MC-based histos
263 if (mUseMC) {
264
265 LOG(info) << "Initializing MC Reader";
266 if (!mMCReader.isInitialized()) {
267 if (!mMCReader.initFromDigitContext("collisioncontext.root")) {
268 throw std::invalid_argument("initialization of MCKinematicsReader failed");
269 }
270 }
271
272 mHistPhiRecVsPhiGen = std::make_unique<TH2F>("mHistPhiRecVsPhiGen", "Phi Rec Vs Phi Gen of true reco tracks ", 24, 0, 2 * TMath::Pi(), 24, 0, 2 * TMath::Pi());
273 mHistPhiRecVsPhiGen->SetXTitle((std::string("#phi of ") + mNameOfTrackTypes[kGen]).c_str());
274 mHistPhiRecVsPhiGen->SetYTitle((std::string("#phi of ") + mNameOfTrackTypes[kRecoTrue]).c_str());
275 mHistPhiRecVsPhiGen->Sumw2();
276 mHistPhiRecVsPhiGen->SetOption("COLZ");
277
278 mHistEtaRecVsEtaGen = std::make_unique<TH2F>("mHistEtaRecVsEtaGen", "Eta Rec Vs Eta Gen of true reco tracks ", 35, -4.5, -1.0, 35, -4.5, -1.0);
279 mHistEtaRecVsEtaGen->SetXTitle((std::string("#eta of ") + mNameOfTrackTypes[kGen]).c_str());
280 mHistEtaRecVsEtaGen->SetYTitle((std::string("#eta of ") + mNameOfTrackTypes[kRecoTrue]).c_str());
281 mHistEtaRecVsEtaGen->Sumw2();
282 mHistEtaRecVsEtaGen->SetOption("COLZ");
283
284 // Histos for Reconstruction assessment
285
286 mChargeMatchEff = std::make_unique<TEfficiency>("QMatchEff", "Charge Match;p_t [GeV];#epsilon", 50, 0, 20);
287
288 const int nTH3Histos = TH3Names.size();
289 auto n3Histo = 0;
290 for (auto& h : mTH3Histos) {
291 h = std::make_unique<TH3F>(TH3Names[n3Histo], TH3Titles[n3Histo],
292 (int)TH3Binning[n3Histo][0],
293 TH3Binning[n3Histo][1],
294 TH3Binning[n3Histo][2],
295 (int)TH3Binning[n3Histo][3],
296 TH3Binning[n3Histo][4],
297 TH3Binning[n3Histo][5],
298 (int)TH3Binning[n3Histo][6],
299 TH3Binning[n3Histo][7],
300 TH3Binning[n3Histo][8]);
301 h->GetXaxis()->SetTitle(TH3XaxisTitles[n3Histo]);
302 h->GetYaxis()->SetTitle(TH3YaxisTitles[n3Histo]);
303 h->GetZaxis()->SetTitle(TH3ZaxisTitles[n3Histo]);
304 ++n3Histo;
305 }
306 }
307}
308
309//__________________________________________________________
311{
312 mNumberTFs++; // TF Counter
313
314 // get tracks
315 mMFTTracks = ctx.inputs().get<gsl::span<o2::mft::TrackMFT>>("tracks");
316 mMFTTracksROF = ctx.inputs().get<gsl::span<o2::itsmft::ROFRecord>>("tracksrofs");
317 mMFTTrackClusIdx = ctx.inputs().get<gsl::span<int>>("trackClIdx");
318
319 // get clusters
320 mMFTClusters = ctx.inputs().get<gsl::span<o2::itsmft::CompClusterExt>>("compClusters");
321 mMFTClustersROF = ctx.inputs().get<gsl::span<o2::itsmft::ROFRecord>>("clustersrofs");
322 mMFTClusterPatterns = ctx.inputs().get<gsl::span<unsigned char>>("patterns");
323 pattIt = mMFTClusterPatterns.begin();
324 mMFTClustersGlobal.clear();
325 mMFTClustersGlobal.reserve(mMFTClusters.size());
326 o2::mft::ioutils::convertCompactClusters(mMFTClusters, pattIt, mMFTClustersGlobal, mDictionary);
327
328 if (mUseMC) {
329 // get labels
330 mMFTClusterLabels = ctx.inputs().get<const dataformats::MCTruthContainer<MCCompLabel>*>("clslabels");
331 mMFTTrackLabels = ctx.inputs().get<gsl::span<MCCompLabel>>("trklabels");
332 }
333
334 // Fill the clusters histograms
335 for (const auto& rof : mMFTClustersROF) {
336 mClusterROFNEntries->Fill(rof.getNEntries());
337 float seconds = orbitToSeconds(rof.getBCData().orbit, mRefOrbit) + rof.getBCData().bc * o2::constants::lhc::LHCBunchSpacingNS * 1e-9;
338 mNOfClustersTime->Fill(seconds, rof.getNEntries());
339 }
340
341 for (int icls = 0; icls < mMFTClusters.size(); ++icls) {
342 auto const oneCluster = mMFTClusters[icls];
343 auto const globalCluster = mMFTClustersGlobal[icls];
344
345 mClusterSensorIndex->Fill(oneCluster.getSensorID());
346 mClusterPatternIndex->Fill(oneCluster.getPatternID());
347
348 mMFTClsZ->Fill(globalCluster.getZ());
349
350 auto clsLayer = mMFTChipMapper.chip2Layer(oneCluster.getChipID());
351 mMFTClsXYinLayer[clsLayer]->Fill(globalCluster.getX(), globalCluster.getY());
352 mMFTClsRinLayer[clsLayer]->Fill(std::sqrt(std::pow(globalCluster.getX(), 2) + std::pow(globalCluster.getY(), 2)));
353 mUnusedChips[oneCluster.getChipID()] = false; // this chipID is used
354 }
355
356 // fill the tracks histogram
357
358 for (const auto& rof : mMFTTracksROF) {
359 mTrackROFNEntries->Fill(rof.getNEntries());
360 mTracksBC->Fill(rof.getBCData().bc, rof.getNEntries());
361 float seconds = orbitToSeconds(rof.getBCData().orbit, mRefOrbit) + rof.getBCData().bc * o2::constants::lhc::LHCBunchSpacingNS * 1e-9;
362 mNOfTracksTime->Fill(seconds, rof.getNEntries());
363 }
364
365 std::array<std::array<int, 2>, 5> clsEntriesForRedundancy;
366
367 for (auto& oneTrack : mMFTTracks) {
368 mTrackNumberOfClusters->Fill(oneTrack.getNumberOfPoints());
369 mTrackChi2->Fill(oneTrack.getChi2OverNDF());
370 mTrackCharge->Fill(oneTrack.getCharge());
371 mTrackPhi->Fill(oneTrack.getPhi());
372 mTrackEta->Fill(oneTrack.getEta());
373 mTrackCotl->Fill(1. / oneTrack.getTanl());
374 mTrackChi2pT->Fill(oneTrack.getPt(), oneTrack.getChi2OverNDF());
375
376 for (auto idisk = 0; idisk < 5; idisk++) {
377 clsEntriesForRedundancy[idisk] = {-1, -1};
378 }
379
380 auto ncls = oneTrack.getNumberOfPoints();
381 auto offset = oneTrack.getExternalClusterIndexOffset();
382
383 for (int icls = 0; icls < ncls; ++icls) // cluster loop
384 {
385
386 auto clsEntry = mMFTTrackClusIdx[offset + icls];
387 auto globalCluster = mMFTClustersGlobal[clsEntry];
388
389 mMFTClsOfTracksZ->Fill(globalCluster.getZ());
390
391 auto layer = mMFTMapping.ChipID2Layer[globalCluster.getSensorID()];
392
393 mMFTClsOfTracksXYinLayer[layer]->Fill(globalCluster.getX(), globalCluster.getY());
394
395 int clsMFTdiskID = layer / 2;
396
397 if (clsEntriesForRedundancy[clsMFTdiskID][0] != -1) {
398 clsEntriesForRedundancy[clsMFTdiskID][1] = clsEntry;
399 } else {
400 clsEntriesForRedundancy[clsMFTdiskID][0] = clsEntry;
401 }
402 }
403
404 for (auto idisk = 0; idisk < 5; idisk++) {
405 if ((clsEntriesForRedundancy[idisk][0] != -1) && (clsEntriesForRedundancy[idisk][1] != -1)) {
406 auto globalCluster1 = mMFTClustersGlobal[clsEntriesForRedundancy[idisk][0]];
407
408 mMFTClsXYRedundantInDisk[idisk]->Fill(globalCluster1.getX(), globalCluster1.getY());
409 }
410 }
411
412 for (auto minNClusters : sMinNClustersList) {
413 if (oneTrack.getNumberOfPoints() >= minNClusters) {
414 mTrackEtaNCls[minNClusters - sMinNClustersList[0]]->Fill(oneTrack.getEta());
415 mTrackPhiNCls[minNClusters - sMinNClustersList[0]]->Fill(oneTrack.getPhi());
416 mTrackXYNCls[minNClusters - sMinNClustersList[0]]->Fill(oneTrack.getX(), oneTrack.getY());
417 mTrackEtaPhiNCls[minNClusters - sMinNClustersList[0]]->Fill(oneTrack.getEta(), oneTrack.getPhi());
418 }
419 }
420
421 if (oneTrack.getCharge() == +1) {
422 mPositiveTrackPhi->Fill(oneTrack.getPhi());
423 mTrackInvQPt->Fill(1 / oneTrack.getPt());
424 }
425
426 if (oneTrack.getCharge() == -1) {
427 mNegativeTrackPhi->Fill(oneTrack.getPhi());
428 mTrackInvQPt->Fill(-1 / oneTrack.getPt());
429 }
430
431 if (oneTrack.isCA()) {
432 mCATrackNumberOfClusters->Fill(oneTrack.getNumberOfPoints());
433 mCATrackEta->Fill(oneTrack.getEta());
434 }
435 if (oneTrack.isLTF()) {
436 mLTFTrackNumberOfClusters->Fill(oneTrack.getNumberOfPoints());
437 mLTFTrackEta->Fill(oneTrack.getEta());
438 }
439 }
440}
441
442//__________________________________________________________
444{
445 for (auto src = 0; src < mMCReader.getNSources(); src++) {
446 for (Int_t event = 0; event < mMCReader.getNEvents(src); event++) {
447 auto evh = mMCReader.getMCEventHeader(src, event);
448 const auto& mcTracks = mMCReader.getTracks(src, event);
449 for (const auto& mcParticle : mcTracks) {
450 addMCParticletoHistos(&mcParticle, kGen, evh);
451 } // mcTracks
453 } // events
454 } // sources
455}
456
457//__________________________________________________________
459{
460 int trackID = 0, evnID = 0, srcID = 0;
461 bool fake = false;
462
463 std::unordered_map<o2::MCCompLabel, std::array<bool, 5>> mcTrackHasClusterInMFTDisks;
464
465 auto nClusters = mMFTClusters.size();
466 for (int icls = 0; icls < nClusters; ++icls) {
467 auto const cluster = mMFTClusters[icls];
468 auto labelSize = (mMFTClusterLabels->getLabels(icls)).size();
469 for (auto il = 0; il < labelSize; il++) {
470 auto clsLabel = (mMFTClusterLabels->getLabels(icls))[il];
471 clsLabel.get(trackID, evnID, srcID, fake);
472 if (fake) {
473 continue;
474 }
475 auto clsLayer = mMFTChipMapper.chip2Layer(cluster.getChipID());
476 int clsMFTdiskID = clsLayer / 2;
477 mcTrackHasClusterInMFTDisks[clsLabel][clsMFTdiskID] = true;
478 }
479 } // loop on clusters
480
481 // Identify trackables
482 mTrackableTracksMap.resize(mMCReader.getNSources());
483 auto src = 0;
484 for (auto& map : mTrackableTracksMap) {
485 map.resize(mMCReader.getNEvents(src++));
486 }
487
488 for (auto& trackClsInDisk : mcTrackHasClusterInMFTDisks) {
489 auto& clsdisk = trackClsInDisk.second;
490 auto nMFTDisks = 0;
491 for (auto disk : {0, 1, 2, 3, 4}) {
492 nMFTDisks += int(clsdisk[disk]);
493 }
494 if (nMFTDisks >= 4) {
495 mMFTTrackables[trackClsInDisk.first] = true;
496 mTrackableTracksMap[trackClsInDisk.first.getSourceID()][trackClsInDisk.first.getEventID()].push_back(trackClsInDisk.first);
497 }
498 }
499
500 // Process trackables
501 for (auto src = 0; src < mMCReader.getNSources(); src++) {
502 for (Int_t event = 0; event < mMCReader.getNEvents(src); event++) {
503 auto evH = mMCReader.getMCEventHeader(src, event);
504 for (auto& trackable : mTrackableTracksMap[src][event]) {
505 auto const* mcParticle = mMCReader.getTrack(trackable);
506 addMCParticletoHistos(mcParticle, kTrackable, evH);
507 }
509 } // events
510 } // sources
511}
512
513//__________________________________________________________
514void MFTAssessment::addMCParticletoHistos(const MCTrack* mcTr, const int trackType, const o2::dataformats::MCEventHeader& evH)
515{
516 auto zVtx = evH.GetZ();
517
518 auto pt = mcTr->GetPt();
519 auto eta = mcTr->GetEta();
520 float phi = TMath::ATan2(mcTr->Py(), mcTr->Px());
521 o2::math_utils::bringTo02Pi(phi);
522 auto z = mcTr->GetStartVertexCoordinatesZ();
523 auto R = sqrt(pow(mcTr->GetStartVertexCoordinatesX(), 2) + pow(mcTr->GetStartVertexCoordinatesY(), 2));
524
525 mHistPtVsEta[trackType]->Fill(eta, pt);
526 mHistPhiVsEta[trackType]->Fill(eta, phi);
527 mHistPhiVsPt[trackType]->Fill(pt, phi);
528 mHistZvtxVsEta[trackType]->Fill(z, eta);
529 mHistRVsZ[trackType]->Fill(z, R);
530 mHistIsPrimary[trackType]->Fill(mcTr->isPrimary());
531}
532
533//__________________________________________________________
535{
536 // For this moment this is used for MC-based assessment, but could be merged into runASyncQC(...)
537 for (auto mftTrack : mMFTTracks) {
538 const auto& pt_Rec = mftTrack.getPt();
539 const auto& eta_Rec = mftTrack.getEta();
540 float phi_Rec = mftTrack.getPhi();
541 o2::math_utils::bringTo02Pi(phi_Rec);
542 const auto& Chi2_Rec = mftTrack.getChi2OverNDF();
543
544 mHistPtVsEta[kReco]->Fill(eta_Rec, pt_Rec);
545 mHistPhiVsEta[kReco]->Fill(eta_Rec, phi_Rec);
546 mHistPhiVsPt[kReco]->Fill(pt_Rec, phi_Rec);
547 mHistTrackChi2[kReco]->Fill(Chi2_Rec);
548 }
549}
550
551//__________________________________________________________
553{
555
556 for (const auto& fakeMFTTrackID : mFakeTracksVec) {
557 auto mftTrack = mMFTTracks[fakeMFTTrackID];
558 const auto& pt_Rec = mftTrack.getPt();
559 const auto& eta_Rec = mftTrack.getEta();
560 float phi_Rec = mftTrack.getPhi();
561 o2::math_utils::bringTo02Pi(phi_Rec);
562 const auto& Chi2_Rec = mftTrack.getChi2OverNDF();
563
564 mHistPtVsEta[kRecoFake]->Fill(eta_Rec, pt_Rec);
565 mHistPhiVsEta[kRecoFake]->Fill(eta_Rec, phi_Rec);
566 mHistPhiVsPt[kRecoFake]->Fill(pt_Rec, phi_Rec);
567 mHistTrackChi2[kRecoFake]->Fill(Chi2_Rec);
568 }
569
570 for (auto src = 0; src < mMCReader.getNSources(); src++) {
571 for (Int_t event = 0; event < mMCReader.getNEvents(src); event++) {
572 auto evH = mMCReader.getMCEventHeader(src, event);
573 auto zVtx = evH.GetZ();
574
575 for (const auto& trueMFTTrackID : mTrueTracksMap[src][event]) {
576 auto mftTrack = mMFTTracks[trueMFTTrackID];
577 const auto& trackLabel = mMFTTrackLabels[trueMFTTrackID];
578 if (trackLabel.isCorrect()) {
579 auto const* mcParticle = mMCReader.getTrack(trackLabel);
580 auto pdgcode_MC = mcParticle->GetPdgCode();
581 int Q_Gen;
582 if (TDatabasePDG::Instance()->GetParticle(pdgcode_MC)) {
583 Q_Gen = TDatabasePDG::Instance()->GetParticle(pdgcode_MC)->Charge() / 3;
584 } else {
585 continue;
586 }
587
588 auto etaGen = mcParticle->GetEta();
589 float phiGen = TMath::ATan2(mcParticle->Py(), mcParticle->Px());
590 o2::math_utils::bringTo02Pi(phiGen);
591 auto ptGen = mcParticle->GetPt();
592 auto vxGen = mcParticle->GetStartVertexCoordinatesX();
593 auto vyGen = mcParticle->GetStartVertexCoordinatesY();
594 auto vzGen = mcParticle->GetStartVertexCoordinatesZ();
595 auto R = sqrt(pow(vxGen, 2) + pow(vyGen, 2));
596 auto tanlGen = mcParticle->Pz() / mcParticle->GetPt();
597 auto invQPtGen = 1.0 * Q_Gen / ptGen;
598
599 if (std::abs(mftTrack.getInvQPt()) > o2::constants::math::Almost0) {
600 mftTrack.propagateToZ(vzGen, mBz);
601 } else {
602 mftTrack.propagateToZlinear(vzGen);
603 }
604
605 const auto& pt_Rec = mftTrack.getPt();
606 const auto& invQPt_Rec = mftTrack.getInvQPt();
607 const auto& invQPt_Seed = mftTrack.getInvQPtSeed();
608 const auto& eta_Rec = mftTrack.getEta();
609 float phi_Rec = mftTrack.getPhi();
610 o2::math_utils::bringTo02Pi(phi_Rec);
611 const auto& nClusters = mftTrack.getNumberOfPoints();
612 const auto& Chi2_Rec = mftTrack.getChi2OverNDF();
613 int Q_Rec = mftTrack.getCharge();
614 // Residuals at vertex
615 auto x_res = mftTrack.getX() - vxGen;
616 auto y_res = mftTrack.getY() - vyGen;
617 auto eta_res = mftTrack.getEta() - etaGen;
618 auto phi_res = mftTrack.getPhi() - phiGen;
619 auto cotl_res = (1. / mftTrack.getTanl()) - (1. / tanlGen);
620 auto invQPt_res = invQPt_Rec - invQPtGen;
621 mHistPtVsEta[kRecoTrue]->Fill(eta_Rec, pt_Rec);
622 mHistPhiVsEta[kRecoTrue]->Fill(eta_Rec, phi_Rec);
623 mHistPhiVsPt[kRecoTrue]->Fill(pt_Rec, phi_Rec);
624 mHistTrackChi2[kRecoTrue]->Fill(Chi2_Rec);
625 mHistZvtxVsEta[kRecoTrue]->Fill(vzGen, eta_Rec);
626 mHistRVsZ[kRecoTrueMC]->Fill(vzGen, R);
627 mHistIsPrimary[kRecoTrueMC]->Fill(mcParticle->isPrimary());
628 mHistPtVsEta[kRecoTrueMC]->Fill(etaGen, ptGen);
629 mHistPhiVsEta[kRecoTrueMC]->Fill(etaGen, phiGen);
630 mHistPhiVsPt[kRecoTrueMC]->Fill(ptGen, phiGen);
631 mHistZvtxVsEta[kRecoTrueMC]->Fill(vzGen, etaGen);
632
633 mHistPhiRecVsPhiGen->Fill(phiGen, phi_Rec);
634 mHistEtaRecVsEtaGen->Fill(etaGen, eta_Rec);
636 auto d_Charge = Q_Rec - Q_Gen;
637 mChargeMatchEff->Fill(!d_Charge, ptGen);
638
639 mTH3Histos[kTH3TrackDeltaXVertexPtEta]->Fill(ptGen, etaGen, 1e4 * x_res);
640 mTH3Histos[kTH3TrackDeltaYVertexPtEta]->Fill(ptGen, etaGen, 1e4 * y_res);
641 mTH3Histos[kTH3TrackDeltaXDeltaYEta]->Fill(etaGen, 1e4 * x_res, 1e4 * y_res);
642 mTH3Histos[kTH3TrackDeltaXDeltaYPt]->Fill(ptGen, 1e4 * x_res, 1e4 * y_res);
643 mTH3Histos[kTH3TrackXPullPtEta]->Fill(ptGen, etaGen, x_res / sqrt(mftTrack.getCovariances()(0, 0)));
644 mTH3Histos[kTH3TrackYPullPtEta]->Fill(ptGen, etaGen, y_res / sqrt(mftTrack.getCovariances()(1, 1)));
645 mTH3Histos[kTH3TrackPhiPullPtEta]->Fill(ptGen, etaGen, phi_res / sqrt(mftTrack.getCovariances()(2, 2)));
646 mTH3Histos[kTH3TrackCotlPullPtEta]->Fill(ptGen, etaGen, cotl_res / sqrt(1. / mftTrack.getCovariances()(3, 3)));
647 mTH3Histos[kTH3TrackInvQPtPullPtEta]->Fill(ptGen, etaGen, invQPt_res / sqrt(mftTrack.getCovariances()(4, 4)));
648 mTH3Histos[kTH3TrackInvQPtResolutionPtEta]->Fill(ptGen, etaGen, (invQPt_Rec - invQPtGen) / invQPtGen);
649 mTH3Histos[kTH3TrackInvQPtResSeedPtEta]->Fill(ptGen, etaGen, (invQPt_Seed - invQPtGen) / invQPtGen);
650 mTH3Histos[kTH3TrackReducedChi2PtEta]->Fill(ptGen, etaGen, Chi2_Rec / (2 * nClusters - 5)); // 5: number of fitting parameters
651 }
652 }
654 } // events
655 } // sources
656}
657
658//__________________________________________________________
659void MFTAssessment::getHistos(TObjArray& objar)
660{
661 TH1F* mMFTDeadChipID = new TH1F("mMFTDeadChipID", "chipID of the dead chips; chipID; # entries", 936, -0.5, 935.5);
662 TH1F* mTFsCounter = new TH1F("mTFsCounter", "counter of TFs; count bin; # entries", 3, 0, 2);
663
664 auto chipID = 0;
665 for (auto chipState : mUnusedChips) {
666 mMFTDeadChipID->Fill(chipID, float(chipState));
667 chipID++;
668 }
669 mTFsCounter->Fill(1, mNumberTFs);
670
671 mMFTDeadChipID->Scale(mNumberTFs);
672
673 objar.Add(mTrackNumberOfClusters.get());
674 objar.Add(mCATrackNumberOfClusters.get());
675 objar.Add(mLTFTrackNumberOfClusters.get());
676 objar.Add(mTrackInvQPt.get());
677 objar.Add(mTrackChi2.get());
678 objar.Add(mTrackCharge.get());
679 objar.Add(mTrackPhi.get());
680 objar.Add(mPositiveTrackPhi.get());
681 objar.Add(mNegativeTrackPhi.get());
682 objar.Add(mTrackEta.get());
683 objar.Add(mTrackChi2pT.get());
684
685 //------
686 objar.Add(mMFTClsZ.get());
687 objar.Add(mMFTClsOfTracksZ.get());
688 objar.Add(mMFTDeadChipID);
689 objar.Add(mTFsCounter);
690 for (auto nMFTLayer = 0; nMFTLayer < 10; nMFTLayer++) {
691 objar.Add(mMFTClsXYinLayer[nMFTLayer].get());
692 objar.Add(mMFTClsRinLayer[nMFTLayer].get());
693 objar.Add(mMFTClsOfTracksXYinLayer[nMFTLayer].get());
694 }
695
696 for (auto nMFTDisk = 0; nMFTDisk < 5; nMFTDisk++) {
697 objar.Add(mMFTClsXYRedundantInDisk[nMFTDisk].get());
698 }
699
700 for (auto minNClusters : sMinNClustersList) {
701 auto nHisto = minNClusters - sMinNClustersList[0];
702 objar.Add(mTrackEtaNCls[nHisto].get());
703 objar.Add(mTrackPhiNCls[nHisto].get());
704 objar.Add(mTrackXYNCls[nHisto].get());
705 objar.Add(mTrackEtaPhiNCls[nHisto].get());
706 }
707 objar.Add(mCATrackEta.get());
708 objar.Add(mLTFTrackEta.get());
709 objar.Add(mTrackCotl.get());
710
711 objar.Add(mTrackROFNEntries.get());
712 objar.Add(mClusterROFNEntries.get());
713
714 objar.Add(mTracksBC.get());
715 objar.Add(mNOfTracksTime.get());
716 objar.Add(mNOfClustersTime.get());
717
718 objar.Add(mClusterSensorIndex.get());
719 objar.Add(mClusterPatternIndex.get());
720
721 for (int trackType = 0; trackType < mLastTrackType; trackType++) {
722 objar.Add(mHistPhiVsEta[trackType].get());
723 objar.Add(mHistPtVsEta[trackType].get());
724 objar.Add(mHistPhiVsPt[trackType].get());
725 if (trackType != kGen || trackType != kTrackable || trackType != kRecoTrueMC) {
726 objar.Add(mHistTrackChi2[trackType].get());
727 }
728 if (trackType != kReco || trackType != kRecoFake) {
729 objar.Add(mHistZvtxVsEta[trackType].get());
730 }
731 if (trackType == kGen || trackType == kTrackable || trackType == kRecoTrueMC) {
732 objar.Add(mHistRVsZ[trackType].get());
733 objar.Add(mHistIsPrimary[trackType].get());
734 }
735 }
736
737 if (mUseMC) {
738 objar.Add(mHistPhiRecVsPhiGen.get());
739 objar.Add(mHistEtaRecVsEtaGen.get());
740
741 // Histos for Reconstruction assessment
742
743 for (auto& h : mTH3Histos) {
744 objar.Add(h.get());
745 }
746
747 if (mFinalizeAnalysis) {
748 objar.Add(mHistVxtOffsetProjection.get());
749 }
750
751 objar.Add(mChargeMatchEff.get());
752
753 if (mFinalizeAnalysis) {
754 for (int slicedCanvas = 0; slicedCanvas < kNSlicedTH3; slicedCanvas++) {
755 objar.Add(mSlicedCanvas[slicedCanvas]);
756 }
757 }
758 }
759}
760
761//__________________________________________________________
762void MFTAssessment::TH3Slicer(TCanvas* canvas, std::unique_ptr<TH3F>& histo3D, std::vector<float> list, double window, int iPar, float marker_size)
763{
764 gStyle->SetOptTitle(kFALSE);
765 gStyle->SetOptStat(0); // Remove title of first histogram from canvas
766 gStyle->SetMarkerStyle(kFullCircle);
767 gStyle->SetMarkerSize(marker_size);
768 canvas->UseCurrentStyle();
769 canvas->cd();
770 std::string cname = canvas->GetName();
771 std::string ctitle = cname;
772 std::string option;
773 std::string option2 = "PLC PMC same";
774
775 TObjArray aSlices;
776 histo3D->GetYaxis()->SetRange(0, 0);
777 histo3D->GetXaxis()->SetRange(0, 0);
778 bool first = true;
779 if (cname.find("VsEta") < cname.length()) {
780 for (auto ptmin : list) {
781 auto ptmax = ptmin + window;
782 histo3D->GetXaxis()->SetRangeUser(ptmin, ptmax);
783
784 std::string ytitle = "\\sigma (";
785 ytitle += histo3D->GetZaxis()->GetTitle();
786 ytitle += ")";
787 auto title = Form("_%1.2f_%1.2f_yz", ptmin, ptmax);
788 auto aDBG = (TH2F*)histo3D->Project3D(title);
789 aDBG->GetXaxis()->SetRangeUser(0, 0);
790
791 aDBG->FitSlicesX(nullptr, 0, -1, 4, "QNR", &aSlices);
792 auto th1DBG = (TH1F*)aSlices[iPar];
793 th1DBG->SetTitle(Form("%1.2f < p_t < %1.2f", ptmin, ptmax));
794 th1DBG->SetStats(0);
795 th1DBG->SetYTitle(ytitle.c_str());
796 if (first) {
797 option = "PLC PMC";
798 } else {
799 option = "SAME PLC PMC";
800 }
801 first = false;
802 th1DBG->DrawClone(option.c_str());
803 }
804 } else if (cname.find("VsPt") < cname.length()) {
805 for (auto etamax : list) {
806 auto etamin = etamax + window;
807 histo3D->GetYaxis()->SetRangeUser(etamin, etamax);
808 std::string ytitle = "\\sigma (" + std::string(histo3D->GetZaxis()->GetTitle()) + ")";
809 auto title = Form("_%1.2f_%1.2f_xz", etamin, etamax);
810 auto aDBG = (TH2F*)histo3D->Project3D(title);
811 aDBG->FitSlicesX(nullptr, 0, -1, 4, "QNR", &aSlices);
812 auto th1DBG = (TH1F*)aSlices[iPar];
813 th1DBG->SetTitle(Form("%1.2f > \\eta > %1.2f", etamax, etamin));
814 th1DBG->SetStats(0);
815 th1DBG->SetYTitle(ytitle.c_str());
816 if (first) {
817 option = "PLC PMC";
818 } else {
819 option = "SAME PLC PMC";
820 }
821 first = false;
822 th1DBG->DrawClone(option.c_str());
823 }
824 } else {
825 exit(1);
826 }
827
828 histo3D->GetYaxis()->SetRange(0, 0);
829 histo3D->GetXaxis()->SetRange(0, 0);
830
831 TPaveText* t = new TPaveText(0.2223748, 0.9069355, 0.7776252, 0.965, "brNDC"); // left-up
832 t->SetBorderSize(0);
833 t->SetFillColor(gStyle->GetTitleFillColor());
834 t->AddText(ctitle.c_str());
835 t->Draw();
836
837 canvas->BuildLegend();
838 canvas->SetTicky();
839 canvas->SetGridy();
840 if (0) {
841 cname += ".png";
842 canvas->Print(cname.c_str());
843 }
844}
845
846//__________________________________________________________
848{
849 if (mFinalizeAnalysis) {
850 throw std::runtime_error("MFTAssessment error: data already loaded");
851 }
852 mFinalizeAnalysis = true;
853
854 TObjArray* objar;
855
856 TFile* f = new TFile(Form("MFTAssessment.root"));
857
858 mTrackNumberOfClusters = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTTrackNumberOfClusters"));
859
860 mCATrackNumberOfClusters = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTCATrackNumberOfClusters"));
861
862 mLTFTrackNumberOfClusters = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTLTFTrackNumberOfClusters"));
863
864 mTrackInvQPt = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTTrackInvQPt"));
865
866 mTrackChi2 = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTTrackChi2"));
867
868 mTrackCharge = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTTrackCharge"));
869
870 mTrackPhi = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTTrackPhi"));
871
872 mPositiveTrackPhi = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTPositiveTrackPhi"));
873
874 mNegativeTrackPhi = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTNegativeTrackPhi"));
875
876 mTrackEta = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTTrackEta"));
877
878 mTrackChi2pT = std::unique_ptr<TH2F>((TH2F*)f->Get("mMFTTrackChi2pT"));
879
880 //---------------------------------------------------------------------------
881
882 mMFTClsZ = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTClsZ"));
883
884 mMFTClsOfTracksZ = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTClsOfTracksZ"));
885
886 for (auto nMFTLayer = 0; nMFTLayer < 10; nMFTLayer++) {
887 mMFTClsXYinLayer[nMFTLayer] = std::unique_ptr<TH2F>((TH2F*)f->Get(Form("mMFTClsXYinLayer%d", nMFTLayer)));
888 mMFTClsRinLayer[nMFTLayer] = std::unique_ptr<TH1F>((TH1F*)f->Get(Form("mMFTClsRinLayer%d", nMFTLayer)));
889 mMFTClsOfTracksXYinLayer[nMFTLayer] = std::unique_ptr<TH2F>((TH2F*)f->Get(Form("mMFTClsOfTracksXYinLayer%d", nMFTLayer)));
890 }
891
892 for (auto nMFTDisk = 0; nMFTDisk < 5; nMFTDisk++) {
893 mMFTClsXYRedundantInDisk[nMFTDisk] = std::unique_ptr<TH2F>((TH2F*)f->Get(Form("mMFTClsXYRedundantInDisk%d", nMFTDisk)));
894 }
895
896 for (auto minNClusters : sMinNClustersList) {
897 auto nHisto = minNClusters - sMinNClustersList[0];
898 mTrackEtaNCls[nHisto] = std::unique_ptr<TH1F>((TH1F*)f->Get(Form("mMFTTrackEta_%d_MinClusters", minNClusters)));
899
900 mTrackPhiNCls[nHisto] = std::unique_ptr<TH1F>((TH1F*)f->Get(Form("mMFTTrackPhi_%d_MinClusters", minNClusters)));
901
902 mTrackXYNCls[nHisto] = std::unique_ptr<TH2F>((TH2F*)f->Get(Form("mMFTTrackXY_%d_MinClusters", minNClusters)));
903
904 mTrackEtaPhiNCls[nHisto] = std::unique_ptr<TH2F>((TH2F*)f->Get(Form("mMFTTrackEtaPhi_%d_MinClusters", minNClusters)));
905 }
906
907 mCATrackEta = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTCATrackEta"));
908
909 mLTFTrackEta = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTLTFTrackEta"));
910
911 mTrackCotl = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTTrackCotl"));
912
913 mClusterROFNEntries = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTClustersROFSize"));
914
915 mTrackROFNEntries = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTTrackROFSize"));
916
917 mTracksBC = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTTracksBC"));
918
919 mNOfTracksTime = std::unique_ptr<TH1F>((TH1F*)f->Get("mNOfTracksTime"));
920
921 mNOfClustersTime = std::unique_ptr<TH1F>((TH1F*)f->Get("mNOfClustersTime"));
922
923 mClusterSensorIndex = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTClusterSensorIndex"));
924
925 mClusterPatternIndex = std::unique_ptr<TH1F>((TH1F*)f->Get("mMFTClusterPatternIndex"));
926
927 for (int trackType = 0; trackType < mLastTrackType; trackType++) {
928 mHistPhiVsEta[trackType] = std::unique_ptr<TH2F>((TH2F*)f->Get((std::string("mHistPhiVsEta") + mNameOfTrackTypes[trackType]).c_str()));
929
930 mHistPtVsEta[trackType] = std::unique_ptr<TH2F>((TH2F*)f->Get((std::string("mHistPtVsEta") + mNameOfTrackTypes[trackType]).c_str()));
931
932 mHistPhiVsPt[trackType] = std::unique_ptr<TH2F>((TH2F*)f->Get((std::string("mHistPhiVsPt") + mNameOfTrackTypes[trackType]).c_str()));
933
934 if (trackType != kGen || trackType != kTrackable || trackType != kRecoTrueMC) {
935 mHistTrackChi2[trackType] = std::unique_ptr<TH1F>((TH1F*)f->Get((std::string("mHistTrackChi2") + mNameOfTrackTypes[trackType]).c_str()));
936 }
937
938 if (trackType != kReco || trackType != kRecoFake) {
939 mHistZvtxVsEta[trackType] = std::unique_ptr<TH2F>((TH2F*)f->Get((std::string("mHistZvtxVsEta") + mNameOfTrackTypes[trackType]).c_str()));
940 }
941 if (trackType == kGen || trackType == kTrackable || trackType == kRecoTrueMC) {
942 mHistRVsZ[trackType] = std::unique_ptr<TH2F>((TH2F*)f->Get((std::string("mHistRVsZ") + mNameOfTrackTypes[trackType]).c_str()));
943 mHistIsPrimary[trackType] = std::unique_ptr<TH1F>((TH1F*)f->Get((std::string("mHistIsPrimary") + mNameOfTrackTypes[trackType]).c_str()));
944 }
945 }
946
947 // Creating MC-based histos
948 if (mUseMC) {
949
950 mHistPhiRecVsPhiGen = std::unique_ptr<TH2F>((TH2F*)f->Get("mHistPhiRecVsPhiGen"));
951
952 mHistEtaRecVsEtaGen = std::unique_ptr<TH2F>((TH2F*)f->Get("mHistEtaRecVsEtaGen"));
953
954 // Histos for Reconstruction assessment
955 mChargeMatchEff = std::unique_ptr<TEfficiency>((TEfficiency*)f->Get("QMatchEff"));
956
957 const int nTH3Histos = TH3Names.size();
958 auto n3Histo = 0;
959 for (auto& h : mTH3Histos) {
960 h = std::unique_ptr<TH3F>((TH3F*)f->Get(TH3Names[n3Histo]));
961 ++n3Histo;
962 }
963 }
964 return true;
965}
966
967//__________________________________________________________
969{
970 mHistVxtOffsetProjection = std::unique_ptr<TH2F>((TH2F*)mTH3Histos[kTH3TrackDeltaXDeltaYEta]->Project3D("colz yz"));
971 mHistVxtOffsetProjection->SetNameTitle("Vertex_XY_OffSet", "Track-MC_Vertex Offset");
972 mHistVxtOffsetProjection->SetOption("COLZ");
973
974 if (mFinalizeAnalysis) {
975 std::vector<float> ptList({.5, 1.5, 5., 10., 15., 18.0});
976 float ptWindow = 0.4;
977 std::vector<float> etaList({-2.5, -2.8, -3.1});
978 float etaWindow = -0.2;
979
980 std::vector<float> sliceList;
981 float sliceWindow;
982
983 for (int nCanvas = 0; nCanvas < kNSlicedTH3; nCanvas++) {
984 if (nCanvas % 2) {
985 sliceList = etaList;
986 sliceWindow = etaWindow;
987 } else {
988 sliceList = ptList;
989 sliceWindow = ptWindow;
990 }
991 mSlicedCanvas[nCanvas] = new TCanvas(TH3SlicedNames[nCanvas], TH3SlicedNames[nCanvas], 1080, 1080);
992 TH3Slicer(mSlicedCanvas[nCanvas], mTH3Histos[TH3SlicedMap[nCanvas]], sliceList, sliceWindow, 2);
993 }
994 }
995}
General auxilliary methods.
o2::itsmft::ChipMappingMFT mMFTMapping
Class to perform assessment of MFT.
int nClusters
Class for time synchronization of RawReader instances.
bool isPrimary() const
Definition MCTrack.h:75
Double_t GetStartVertexCoordinatesY() const
Definition MCTrack.h:83
Double_t GetPt() const
Definition MCTrack.h:109
Double_t GetStartVertexCoordinatesZ() const
Definition MCTrack.h:84
Double_t Py() const
Definition MCTrack.h:102
Double_t Px() const
Definition MCTrack.h:101
Double_t GetEta() const
Definition MCTrack.h:131
Double_t GetStartVertexCoordinatesX() const
Definition MCTrack.h:82
Int_t GetPdgCode() const
Accessors.
Definition MCTrack.h:72
A container to hold and manage MC truth information/labels.
decltype(auto) get(R binding, int part=0) const
InputRecord & inputs()
The inputs associated with this processing context.
Int_t chip2Layer(Int_t chipID) const
impose user defined FEEId -> ruSW (staveID) conversion, to be used only for forced decoding of corrup...
static constexpr std::array< int, NChips > ChipID2Layer
void addMCParticletoHistos(const MCTrack *mcTr, const int TrackType, const o2::dataformats::MCEventHeader &evH)
void runASyncQC(o2::framework::ProcessingContext &ctx)
void init(bool finalizeAnalysis)
double orbitToSeconds(uint32_t orbit, uint32_t refOrbit)
void getHistos(TObjArray &objar)
bool initFromDigitContext(std::string_view filename)
MCTrack const * getTrack(o2::MCCompLabel const &) const
size_t getNEvents(int source) const
Get number of events.
o2::dataformats::MCEventHeader const & getMCEventHeader(int source, int event) const
retrieves the MCEventHeader for a given eventID and sourceID
void releaseTracksForSourceAndEvent(int source, int event)
API to ask releasing tracks (freeing memory) for source + event.
size_t getNSources() const
Get number of sources.
std::vector< MCTrack > const & getTracks(int source, int event) const
variant returning all tracks for source and event at once
struct _cl_event * event
Definition glcorearb.h:2982
GLenum src
Definition glcorearb.h:1767
GLsizeiptr size
Definition glcorearb.h:659
GLdouble f
Definition glcorearb.h:310
GLintptr offset
Definition glcorearb.h:660
GLenum GLuint GLint GLint layer
Definition glcorearb.h:1310
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
constexpr int LHCMaxBunches
constexpr double LHCBunchSpacingNS
constexpr float Almost0
auto get(const std::byte *buffer, size_t=0)
Definition DataHeader.h:454
void convertCompactClusters(gsl::span< const itsmft::CompClusterExt > clusters, gsl::span< const unsigned char >::iterator &pattIt, std::vector< o2::BaseCluster< float > > &output, const itsmft::TopologyDictionary *dict)
convert compact clusters to 3D spacepoints into std::vector<o2::BaseCluster<float>>
Definition IOUtils.cxx:111
@ kNumberOfTrackTypes
Definition list.h:40
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"