-
Notifications
You must be signed in to change notification settings - Fork 490
Expand file tree
/
Copy pathcheckTOFMatching.C
More file actions
528 lines (484 loc) · 21.2 KB
/
checkTOFMatching.C
File metadata and controls
528 lines (484 loc) · 21.2 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
// All rights not expressly granted are reserved.
//
// This software is distributed under the terms of the GNU General Public
// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
//
// In applying this license CERN does not waive the privileges and immunities
// granted to it by virtue of its status as an Intergovernmental Organization
// or submit itself to any jurisdiction.
#if !defined(__CLING__) || defined(__ROOTCLING__)
#include "TFile.h"
#include "TF1.h"
#include "TTree.h"
#include "TH1F.h"
#include "TH2F.h"
#include "TLine.h"
#include "TProfile.h"
#include "TMath.h"
#include "TCanvas.h"
#include "GlobalTracking/MatchTOF.h"
#include "ReconstructionDataFormats/TrackTPCITS.h"
#include "SimulationDataFormat/MCTruthContainer.h"
#include "SimulationDataFormat/MCCompLabel.h"
#include "DataFormatsTOF/Cluster.h"
#include <fairlogger/Logger.h>
#endif
//#define DEBUG
void checkTOFMatching(bool batchMode = true)
{
// macro to check the matching TOF-ITSTPC tracks
// getting TOF info
TFile* fmatchTOF = new TFile("o2match_tof_itstpc.root");
TTree* matchTOF = (TTree*)fmatchTOF->Get("matchTOF");
std::vector<o2::dataformats::MatchInfoTOF>* TOFMatchInfo;
TOFMatchInfo = new std::vector<o2::dataformats::MatchInfoTOF>;
matchTOF->SetBranchAddress("TOFMatchInfo", &TOFMatchInfo);
// getting the ITSTPCtracks
TFile* fmatchITSTPC = new TFile("o2match_itstpc.root");
TTree* matchTPCITS = (TTree*)fmatchITSTPC->Get("matchTPCITS");
std::vector<o2::dataformats::TrackTPCITS>* mTracksArrayInp = new std::vector<o2::dataformats::TrackTPCITS>;
matchTPCITS->SetBranchAddress("TPCITS", &mTracksArrayInp);
// getting the TPC tracks
TFile* ftracksTPC = new TFile("tpctracks.root");
TTree* tpcTree = (TTree*)ftracksTPC->Get("tpcrec");
std::vector<o2::tpc::TrackTPC>* mTPCTracksArrayInp = new std::vector<o2::tpc::TrackTPC>;
tpcTree->SetBranchAddress("TPCTracks", &mTPCTracksArrayInp);
std::vector<o2::MCCompLabel>* mcTPC = new std::vector<o2::MCCompLabel>();
tpcTree->SetBranchAddress("TPCTracksMCTruth", &mcTPC);
// getting the ITS tracks
TFile* ftracksITS = new TFile("o2trac_its.root");
TTree* itsTree = (TTree*)ftracksITS->Get("o2sim");
std::vector<o2::its::TrackITS>* mITSTracksArrayInp = new std::vector<o2::its::TrackITS>;
itsTree->SetBranchAddress("ITSTrack", &mITSTracksArrayInp);
std::vector<o2::MCCompLabel>* mcITS = new std::vector<o2::MCCompLabel>();
itsTree->SetBranchAddress("ITSTrackMCTruth", &mcITS);
// getting the TOF clusters
TFile* fclustersTOF = new TFile("tofclusters.root");
TTree* tofClTree = (TTree*)fclustersTOF->Get("o2sim");
std::vector<o2::tof::Cluster>* mTOFClustersArrayInp = new std::vector<o2::tof::Cluster>;
tofClTree->SetBranchAddress("TOFCluster", &mTOFClustersArrayInp);
o2::dataformats::MCTruthContainer<o2::MCCompLabel>* mcTOF = new o2::dataformats::MCTruthContainer<o2::MCCompLabel>();
tofClTree->SetBranchAddress("TOFClusterMCTruth", &mcTOF);
int nMatches = 0;
int nGoodMatches = 0;
int nBadMatches = 0;
TH1F* htrack = new TH1F("htrack", "tracks;p_{T} (GeV/c);N", 50, 0, 5);
TH1F* htof = new TH1F("htof", "tof matching |#eta|<0.9;p_{T} (GeV/c);#varepsilon", 50, 0, 5);
TH1F* htofGood = new TH1F("htofGood", "tof matching |#eta|<0.9;p_{T} (GeV/c);#varepsilon", 50, 0, 5);
TH1F* htofMism = new TH1F("htofMism", "tof matching |#eta|<0.9;p_{T} (GeV/c);#varepsilon", 50, 0, 5);
htofMism->SetLineColor(2);
TH1F* htrack_t = new TH1F("htrack_t", "tracks p_{T} > 0.5 GeV/c;Timestamp (#mus);N", 100, 0, 20000);
TH1F* htof_t = new TH1F("htof_t", "tof matching p_{T} > 0.5 GeV/c, |#eta|<0.9;Timestamp (#mus);#varepsilon", 100, 0, 20000);
TH1F* htof_tMism = new TH1F("htof_tMism", "tof matching p_{T} > 0.5 GeV/c, |#eta|<0.9;Timestamp (#mus);#varepsilon", 100, 0, 20000);
htof_tMism->SetLineColor(2);
TH1F* htrack_eta = new TH1F("htrack_eta", "tracks p_{T} > 0.5 GeV/c;#eta;N", 18, -0.9, 0.9);
TH1F* htof_eta = new TH1F("htof_eta", "tof matching p_{T} > 0.5 GeV/c;#eta;#varepsilon", 18, -0.9, 0.9);
TH1F* htof_etaMism = new TH1F("htof_etaMism", "tof matching p_{T} > 0.5 GeV/c;#eta;#varepsilon", 18, -0.9, 0.9);
htof_etaMism->SetLineColor(2);
TH1F* htof_res = new TH1F("htof_res", "tof residuals;residuals (cm)", 100, 0, 10);
TH1F* htof_resMism = new TH1F("htof_resMism", "tof residuals;residuals (cm)", 100, 0, 10);
htof_resMism->SetLineColor(2);
TH1F* hdeltatime = new TH1F("hdeltatime", "#Deltat TOF-Track;t_{TOF} - t_{track} (#mus);N", 100, -1, 1);
TH1F* hdeltatime_sigma = new TH1F("hdeltatime_sigma", ";#Deltat/#sigma;N", 100, -30, 30);
TH1F* hdeltatimeMism = new TH1F("hdeltatimeMism", "#Deltat TOF-Track;t_{TOF} - t_{track} (#mus);N", 100, -1, 1);
TH1F* hdeltatime_sigmaMism = new TH1F("hdeltatime_sigmaMism", ";#Deltat/#sigma;N", 100, -30, 30);
hdeltatimeMism->SetLineColor(2);
hdeltatime_sigmaMism->SetLineColor(2);
TH2F* hchi2 = new TH2F("hchi2", "#Sum of residuals distribution;p_{T} (GeV/c);residuals (cm)", 25, 0, 5, 100, 0, 10);
TH2F* hchi2sh = new TH2F("hchi2sh", "#Sum of residuals distribution, single hit;p_{T} (GeV/c);residuals (cm)", 25, 0, 5, 100, 0, 10);
TH2F* hchi2dh = new TH2F("hchi2dh", "#Sum of residuals distribution, double hits;p_{T} (GeV/c);residuals (cm)", 25, 0, 5, 100, 0, 10);
TH2F* hchi2th = new TH2F("hchi2th", "#Sum of residuals distribution, triple hits;p_{T} (GeV/c);residuals (cm)", 25, 0, 5, 100, 0, 10);
TH1F* hMaterial = new TH1F("hMaterial", ";Material budget (X0X2 units);", 100, 0, 20);
TProfile* hMismVsMaterial = new TProfile("hMismVsMaterial", "Mismatch fraction;Material budget (X0X2 units);Mismatch", 100, 0, 20);
hMismVsMaterial->SetLineColor(2);
TProfile* hResVsMaterial = new TProfile("hResVsMaterial", "Mismatch fraction;Material budget (X0X2 units);Mismatch", 100, 0, 20);
hResVsMaterial->SetLineColor(4);
// now looping over the entries in the matching tree
for (int ientry = 0; ientry < matchTOF->GetEntries(); ientry++) {
matchTOF->GetEvent(ientry);
matchTPCITS->GetEntry(ientry);
tpcTree->GetEntry(ientry);
tofClTree->GetEntry(ientry);
itsTree->GetEntry(ientry);
// loop over tracks
for (uint i = 0; i < mTracksArrayInp->size(); i++) {
o2::dataformats::TrackTPCITS trackITSTPC = mTracksArrayInp->at(i);
if (std::abs(trackITSTPC.getEta()) < 0.9) {
htrack->Fill(trackITSTPC.getPt());
if (trackITSTPC.getPt() > 0.5) {
htrack_t->Fill(trackITSTPC.getTimeMUS().getTimeStamp());
htrack_eta->Fill(trackITSTPC.getEta());
}
}
}
// now looping over the matched tracks
nMatches += TOFMatchInfo->size();
for (uint imatch = 0; imatch < TOFMatchInfo->size(); imatch++) {
// get ITS-TPC track
int indexITSTPCtrack = TOFMatchInfo->at(imatch).getTrackIndex();
o2::dataformats::MatchInfoTOF infoTOF = TOFMatchInfo->at(imatch);
int tofClIndex = infoTOF.getTOFClIndex();
float chi2 = infoTOF.getChi2();
#ifdef DEBUG
LOGF(info, "nentry in tree %d, matching %d, indexITSTPCtrack = %d, tofClIndex = %d, chi2 = %f", ientry, imatch, indexITSTPCtrack, tofClIndex, chi2);
#endif
float matBud = infoTOF.getLTIntegralOut().getX2X0();
// o2::MCCompLabel label = mcTOF->getElement(mcTOF->getMCTruthHeader(tofClIndex).index);
const auto& labelsTOF = mcTOF->getLabels(tofClIndex);
int trackIdTOF = -1;
int eventIdTOF = -1;
int sourceIdTOF = -1;
for (uint ilabel = 0; ilabel < labelsTOF.size(); ilabel++) {
#ifdef DEBUG
LOGF(info, "TOF label %d: trackID = %d, eventID = %d, sourceID = %d", ilabel, labelsTOF[ilabel].getTrackID(), labelsTOF[ilabel].getEventID(), labelsTOF[ilabel].getSourceID());
#endif
if (ilabel == 0) {
trackIdTOF = labelsTOF[ilabel].getTrackID();
eventIdTOF = labelsTOF[ilabel].getEventID();
sourceIdTOF = labelsTOF[ilabel].getSourceID();
}
}
o2::tof::Cluster tofCluster = mTOFClustersArrayInp->at(tofClIndex);
int nContributingChannels = tofCluster.getNumOfContributingChannels();
int mainContributingChannel = tofCluster.getMainContributingChannel();
#ifdef DEBUG
LOGF(info, "The TOF cluster has %d contributing channels, and the main one is %d", nContributingChannels, mainContributingChannel);
#endif
int indices[5];
o2::tof::Geo::getVolumeIndices(mainContributingChannel, indices);
#ifdef DEBUG
LOGF(info, "Indices of main contributing channel are %d, %d, %d, %d, %d", indices[0], indices[1], indices[2], indices[3], indices[4]);
#endif
bool isUpLeft = tofCluster.isAdditionalChannelSet(o2::tof::Cluster::kUpLeft);
bool isUp = tofCluster.isAdditionalChannelSet(o2::tof::Cluster::kUp);
bool isUpRight = tofCluster.isAdditionalChannelSet(o2::tof::Cluster::kUpRight);
bool isRight = tofCluster.isAdditionalChannelSet(o2::tof::Cluster::kRight);
bool isDownRight = tofCluster.isAdditionalChannelSet(o2::tof::Cluster::kDownRight);
bool isDown = tofCluster.isAdditionalChannelSet(o2::tof::Cluster::kDown);
bool isDownLeft = tofCluster.isAdditionalChannelSet(o2::tof::Cluster::kDownLeft);
bool isLeft = tofCluster.isAdditionalChannelSet(o2::tof::Cluster::kLeft);
#ifdef DEBUG
LOGF(info, "isUpLeft = %d, isUp = %d, isUpRight = %d, isRight = %d, isDownRight = %d, isDown = %d, isDownLeft = %d, isLeft = %d", isUpLeft, isUp, isUpRight, isRight, isDownRight, isDown, isDownLeft, isLeft);
#endif
int indexCont[5];
indexCont[0] = indices[0];
indexCont[1] = indices[1];
indexCont[2] = indices[2];
indexCont[3] = indices[3];
indexCont[4] = indices[4];
int numberOfSecondaryContributingChannels = 0;
int secondaryContributingChannel = -1;
if (isDown) {
indexCont[3]--;
numberOfSecondaryContributingChannels++;
secondaryContributingChannel = o2::tof::Geo::getIndex(indexCont);
#ifdef DEBUG
LOGF(info, "secondaryContributingChannel[down] = %d", secondaryContributingChannel);
#endif
indexCont[3] = indices[3];
}
if (isDownRight) {
indexCont[3]--;
indexCont[4]++;
numberOfSecondaryContributingChannels++;
secondaryContributingChannel = o2::tof::Geo::getIndex(indexCont);
#ifdef DEBUG
LOGF(info, "secondaryContributingChannel[downright] = %d", secondaryContributingChannel);
#endif
indexCont[3] = indices[3];
indexCont[4] = indices[4];
}
if (isDownLeft) {
indexCont[3]--;
indexCont[4]--;
numberOfSecondaryContributingChannels++;
secondaryContributingChannel = o2::tof::Geo::getIndex(indexCont);
#ifdef DEBUG
LOGF(info, "secondaryContributingChannel[downleft] = %d", secondaryContributingChannel);
#endif
indexCont[3] = indices[3];
indexCont[4] = indices[4];
}
if (isUp) {
indexCont[3]++;
numberOfSecondaryContributingChannels++;
secondaryContributingChannel = o2::tof::Geo::getIndex(indexCont);
#ifdef DEBUG
LOGF(info, "secondaryContributingChannel[up] = %d", secondaryContributingChannel);
#endif
indexCont[3] = indices[3];
}
if (isUpRight) {
indexCont[3]++;
indexCont[4]++;
numberOfSecondaryContributingChannels++;
secondaryContributingChannel = o2::tof::Geo::getIndex(indexCont);
#ifdef DEBUG
LOGF(info, "secondaryContributingChannel[upright] = %d", secondaryContributingChannel);
#endif
indexCont[3] = indices[3];
indexCont[4] = indices[4];
}
if (isUpLeft) { // increase padZ
indexCont[3]++;
indexCont[4]--;
numberOfSecondaryContributingChannels++;
secondaryContributingChannel = o2::tof::Geo::getIndex(indexCont);
#ifdef DEBUG
LOGF(info, "secondaryContributingChannel[upleft] = %d", secondaryContributingChannel);
#endif
indexCont[3] = indices[3];
indexCont[4] = indices[4];
}
if (isRight) { // increase padX
indexCont[4]++;
numberOfSecondaryContributingChannels++;
secondaryContributingChannel = o2::tof::Geo::getIndex(indexCont);
#ifdef DEBUG
LOGF(info, "secondaryContributingChannel[right] = %d", secondaryContributingChannel);
#endif
indexCont[4] = indices[4];
}
if (isLeft) { // decrease padX
indexCont[4]--;
numberOfSecondaryContributingChannels++;
secondaryContributingChannel = o2::tof::Geo::getIndex(indexCont);
#ifdef DEBUG
LOGF(info, "secondaryContributingChannel[left] = %d", secondaryContributingChannel);
#endif
indexCont[4] = indices[4];
}
#ifdef DEBUG
LOGF(info, "Total number of secondary channels= %d", numberOfSecondaryContributingChannels);
#endif
o2::dataformats::TrackTPCITS trackITSTPC = mTracksArrayInp->at(indexITSTPCtrack);
const auto evIdxTPC = trackITSTPC.getRefTPC();
#ifdef DEBUG
LOG(info) << "matched TPCtrack index:", evIdxTPC;
#endif
const auto evIdxITS = trackITSTPC.getRefITS();
#ifdef DEBUG
LOG(info) << "matched ITStrack index: ", evIdxITS;
#endif
// getting the TPC labels
const auto& labelsTPC = (*mcTPC)[evIdxTPC];
#ifdef DEBUG
LOGF(info, "TPC label: trackID = %d, eventID = %d, sourceID = %d", labelsTPC.getTrackID(), labelsTPC.getEventID(), labelsTPC.getSourceID());
#endif
// getting the ITS labels
const auto& labelsITS = (*mcITS)[evIdxITS];
#ifdef DEBUG
LOGF(info, "ITS label: trackID = %d, eventID = %d, sourceID = %d", labelsITS.getTrackID(), labelsITS.getEventID(), labelsITS.getSourceID());
#endif
bool bMatched = kFALSE;
for (uint ilabel = 0; ilabel < labelsTOF.size(); ilabel++) {
if ((std::abs(labelsTPC.getTrackID()) == labelsTOF[ilabel].getTrackID() && labelsTPC.getEventID() == labelsTOF[ilabel].getEventID() && labelsTPC.getSourceID() == labelsTOF[ilabel].getSourceID()) || (labelsITS.getTrackID() == labelsTOF[ilabel].getTrackID() && labelsITS.getEventID() == labelsTOF[ilabel].getEventID() && labelsITS.getSourceID() == labelsTOF[ilabel].getSourceID())) {
nGoodMatches++;
bMatched = kTRUE;
break;
}
}
if (!bMatched)
nBadMatches++;
hMaterial->Fill(matBud);
hMismVsMaterial->Fill(matBud, !bMatched);
if (bMatched)
hResVsMaterial->Fill(matBud, chi2);
if (bMatched)
htof_res->Fill(chi2);
else
htof_resMism->Fill(chi2);
float deltatime = tofCluster.getTime() * 1E-6 - trackITSTPC.getTimeMUS().getTimeStamp(); // in mus
if (bMatched) {
hdeltatime->Fill(deltatime);
hdeltatime_sigma->Fill(deltatime / trackITSTPC.getTimeMUS().getTimeStampError());
hchi2->Fill(trackITSTPC.getPt(), TOFMatchInfo->at(imatch).getChi2());
} else {
hdeltatimeMism->Fill(deltatime);
hdeltatime_sigmaMism->Fill(deltatime / trackITSTPC.getTimeMUS().getTimeStampError());
}
if (bMatched) {
if (nContributingChannels == 1)
hchi2sh->Fill(trackITSTPC.getPt(), TOFMatchInfo->at(imatch).getChi2());
else if (nContributingChannels == 2)
hchi2dh->Fill(trackITSTPC.getPt(), TOFMatchInfo->at(imatch).getChi2());
else if ((isRight || isLeft) && (isUp || isDown)) {
hchi2th->Fill(trackITSTPC.getPt(), TOFMatchInfo->at(imatch).getChi2());
}
}
if (std::abs(trackITSTPC.getEta()) < 0.9) {
htof->Fill(trackITSTPC.getPt());
if (bMatched)
htofGood->Fill(trackITSTPC.getPt());
else
htofMism->Fill(trackITSTPC.getPt());
if (trackITSTPC.getPt() > 0.5) {
if (bMatched)
htof_t->Fill(trackITSTPC.getTimeMUS().getTimeStamp());
else
htof_tMism->Fill(trackITSTPC.getTimeMUS().getTimeStamp());
if (bMatched)
htof_eta->Fill(trackITSTPC.getEta());
else
htof_etaMism->Fill(trackITSTPC.getEta());
}
}
bool TPCfound = false;
bool ITSfound = false;
for (uint i = 0; i < mTracksArrayInp->size(); i++) {
o2::dataformats::TrackTPCITS trackITSTPC = mTracksArrayInp->at(i);
const auto idxTPCcheck = trackITSTPC.getRefTPC();
const auto idxITScheck = trackITSTPC.getRefITS();
const auto& labelsTPCcheck = (*mcTPC)[idxTPCcheck.getIndex()];
if (std::abs(labelsTPCcheck.getTrackID()) == trackIdTOF && labelsTPCcheck.getEventID() == eventIdTOF && labelsTPCcheck.getSourceID() == sourceIdTOF) {
#ifdef DEBUG
LOGF(info, "The TPC track that should have been matched to TOF is number %d", i);
#endif
TPCfound = true;
}
const auto& labelsITScheck = (*mcITS)[idxITScheck.getIndex()];
if (labelsITScheck.getTrackID() == trackIdTOF && labelsITScheck.getEventID() == eventIdTOF && labelsITScheck.getSourceID() == sourceIdTOF) {
#ifdef DEBUG
LOGF(info, "The ITS track that should have been matched to TOF is number %d", i);
#endif
ITSfound = true;
}
}
#ifdef DEBUG
if (!TPCfound)
LOGF(info, "There is no TPC track found that should have corresponded to this TOF cluster!");
if (!ITSfound)
LOGF(info, "There is no ITS track found that should have corresponded to this TOF cluster!");
#endif
}
}
new TCanvas;
LOGF(info, "Number of matches = %d", nMatches);
LOGF(info, "Number of GOOD matches = %d (%.2f)", nGoodMatches, (float)nGoodMatches / nMatches);
LOGF(info, "Number of BAD matches = %d (%.2f)", nBadMatches, (float)nBadMatches / nMatches);
TFile* fout = nullptr;
if (batchMode)
fout = new TFile("tofmatching_qa.root", "RECREATE");
htofMism->Divide(htofMism, htof, 1, 1, "B");
htof->Divide(htofGood, htrack, 1, 1, "B");
htof->SetMarkerStyle(20);
if (batchMode) {
htof->Write();
htofMism->Write();
} else {
htof->Draw("P");
htofMism->Draw("SAME");
}
htof_eta->Divide(htof_eta, htrack_eta, 1, 1, "B");
htof_etaMism->Divide(htof_etaMism, htrack_eta, 1, 1, "B");
htof_eta->SetMarkerStyle(20);
if (batchMode) {
htof_eta->Write();
htof_etaMism->Write();
} else {
new TCanvas;
htof_eta->Draw("P");
htof_etaMism->Draw("same");
}
htof_t->Divide(htof_t, htrack_t, 1, 1, "B");
htof_tMism->Divide(htof_tMism, htrack_t, 1, 1, "B");
htof_t->SetMarkerStyle(20);
if (batchMode) {
htof_t->Write();
htof_tMism->Write();
hdeltatime->Write();
hdeltatimeMism->Write();
hdeltatime_sigma->Write();
hdeltatime_sigmaMism->Write();
} else {
new TCanvas;
htof_t->Draw("P");
htof_tMism->Draw("same");
new TCanvas;
hdeltatime->Draw();
hdeltatimeMism->Draw("SAME");
new TCanvas;
hdeltatime_sigma->Draw();
hdeltatime_sigmaMism->Draw("SAME");
}
if (batchMode) {
hchi2->Write();
hchi2sh->Write();
hchi2dh->Write();
hchi2th->Write();
} else {
TCanvas* cres = new TCanvas();
cres->Divide(2, 2);
cres->cd(1)->SetLogz();
hchi2->Draw("colz");
hchi2->ProfileX()->Draw("same");
TLine* l = new TLine(0, 0.983575, 5, 0.983575);
l->Draw("SAME");
l->SetLineStyle(2);
l->SetLineWidth(2);
l->SetLineColor(4);
cres->cd(2)->SetLogz();
hchi2sh->Draw("colz");
hchi2sh->ProfileX()->Draw("same");
TLine* l2 = new TLine(0, 1.044939, 5, 1.044939);
l2->Draw("SAME");
l2->SetLineStyle(2);
l2->SetLineWidth(2);
l2->SetLineColor(4);
cres->cd(3)->SetLogz();
hchi2dh->Draw("colz");
hchi2dh->ProfileX()->Draw("same");
TLine* l3 = new TLine(0, 0.73811975, 5, 0.73811975);
l3->Draw("SAME");
l3->SetLineStyle(2);
l3->SetLineWidth(2);
l3->SetLineColor(4);
cres->cd(4)->SetLogz();
hchi2th->Draw("colz");
hchi2th->ProfileX()->Draw("same");
TLine* l4 = new TLine(0, 0.3, 5, 0.3);
l4->Draw("SAME");
l4->SetLineStyle(2);
l4->SetLineWidth(2);
l4->SetLineColor(4);
}
if (batchMode) {
htof_res->Write();
htof_resMism->Write();
hMaterial->Write();
hMismVsMaterial->Write();
hResVsMaterial->Write();
} else {
TCanvas* cresiduals = new TCanvas("cresiduals", "cresiduals");
htof_res->Draw();
htof_resMism->Draw("SAME");
TCanvas* cmaterial = new TCanvas("cmaterial", "cmaterial");
hMaterial->DrawNormalized("", 10);
hMismVsMaterial->Draw("SAME");
hResVsMaterial->Draw("SAME");
}
float fraction = hchi2dh->GetEntries() * 1. / hchi2->GetEntries();
float fractionErr = TMath::Sqrt(fraction * (1 - fraction) / hchi2->GetEntries());
LOGF(info, "Fraction of multiple hits = (%.1f +/- %.1f)%c", fraction * 100, fractionErr * 100, '%');
htof->Fit("pol0", "", "", 1, 5);
float effMatch = 0;
if (htof->GetListOfFunctions()->At(0))
effMatch = ((TF1*)htof->GetListOfFunctions()->At(0))->GetParameter(0);
float effMatchErr = 0;
if (htof->GetListOfFunctions()->At(0))
effMatchErr = ((TF1*)htof->GetListOfFunctions()->At(0))->GetParError(0);
LOGF(info, "TOF matching eff (pt > 1) = %f +/- %f", effMatch, effMatchErr);
htofMism->Fit("pol0", "", "", 1, 5);
float mismMatch = 0;
if (htofMism->GetListOfFunctions()->At(0))
mismMatch = ((TF1*)htofMism->GetListOfFunctions()->At(0))->GetParameter(0);
float mismMatchErr = 0;
if (htofMism->GetListOfFunctions()->At(0))
mismMatchErr = ((TF1*)htofMism->GetListOfFunctions()->At(0))->GetParError(0);
LOGF(info, "TOF-track mismatch (pt > 1) = %f +/- %f", mismMatch, mismMatchErr);
if (fout)
fout->Close();
return;
}