CbmRoot
Loading...
Searching...
No Matches
CbmTrackingTrdQa.cxx
Go to the documentation of this file.
1/* Copyright (C) 2006-2021 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Sergey Gorbunov, Denis Bertini [committer], Volker Friese, Florian Uhlig */
4
5// -------------------------------------------------------------------------
6// ----- CbmTrackingTrdQa source file -----
7// ----- Created 11/01/06 by V. Friese -----
8// -------------------------------------------------------------------------
9
10// Includes class header
11#include "CbmTrackingTrdQa.h"
12
13// Includes from C++
14#include <cassert>
15#include <iomanip>
16
17// Includes from ROOT
18#include "TClonesArray.h"
19#include "TGeoManager.h"
20#include "TH1F.h"
21#include "TH2F.h"
22
23#include <TPDGCode.h>
24
25// Includes from FairRoot
26#include "FairEventHeader.h"
27#include "FairRun.h"
28
29#include <Logger.h>
30
31// Includes from CbmRoot
32#include "CbmGlobalTrack.h"
33#include "CbmMCDataArray.h"
34#include "CbmMCDataManager.h"
35#include "CbmMCTrack.h"
36#include "CbmRecoSetupManager.h"
37#include "CbmStsTrack.h"
38#include "CbmTimeSlice.h"
39#include "CbmTrackMatchNew.h"
40#include "CbmTrdHit.h"
41#include "CbmTrdPoint.h"
42#include "CbmTrdTrack.h"
43#include "FairRunAna.h"
44
45using std::fixed;
46using std::right;
47using std::setprecision;
48using std::setw;
49
50const char* CbmTrackingTrdQa::fgkIdxName[fgkNpdg] = {"e", "#mu", "#pi", "K", "#bf{p}", "any"};
51const char* CbmTrackingTrdQa::fgkIdxSymb[fgkNpdg] = {"e", "mu", "pi", "K", "p", "x"};
52
53// -------------------------------------------------------------------------
55{
56 int idx(5);
57 switch (pdg) {
58 case kElectron:
59 case kPositron: idx = 0; break;
60 case kMuonMinus:
61 case kMuonPlus: idx = 1; break;
62 case kPiMinus:
63 case kPiPlus: idx = 2; break;
64 case kKMinus:
65 case kKPlus: idx = 3; break;
66 case kProton:
67 case kProtonBar: idx = 4; break;
68 default: idx = 5; break;
69 }
70 return idx;
71}
73{
74 if (idx < 0 || idx >= 5) return 0;
75 int pdg = 0;
76 switch (idx) {
77 case 0: pdg = kElectron; break;
78 case 1: pdg = kMuonMinus; break;
79 case 2: pdg = kPiPlus; break;
80 case 3: pdg = kKPlus; break;
81 case 4: pdg = kProton; break;
82 }
83 return pdg;
84}
85const char* CbmTrackingTrdQa::Idx2Name(int idx)
86{
87 if (idx < 0 || idx >= 5) return "non";
88 return fgkIdxName[idx];
89}
90const char* CbmTrackingTrdQa::Idx2Symb(int idx)
91{
92 if (idx < 0 || idx >= 5) return "o";
93 return fgkIdxSymb[idx];
94}
95
96
97// ----- Default constructor -------------------------------------------
98CbmTrackingTrdQa::CbmTrackingTrdQa(Int_t iVerbose) : FairTask("CbmTrackingTrdQa", iVerbose) {}
99
100// -------------------------------------------------------------------------
101
102
103// ----- Standard constructor ------------------------------------------
104CbmTrackingTrdQa::CbmTrackingTrdQa(Int_t minStations, Double_t quota, Int_t iVerbose)
105 : FairTask("CbmTrackingTrdQa", iVerbose)
106 , fMinStations(minStations)
107 , fQuota(quota)
108{
109}
110// -------------------------------------------------------------------------
111
112
113// ----- Destructor ----------------------------------------------------
115{
116 fHistoList->Delete();
117 delete fHistoList;
118}
119// -------------------------------------------------------------------------
120
121
122// ----- Task execution ------------------------------------------------
123void CbmTrackingTrdQa::Exec(Option_t* /*opt*/)
124{
125
126 LOG(debug) << GetName() << ": Process event ";
127
128 // Timer
129 fTimer.Start();
130
131 // Eventwise counters
132 // Int_t nMCTracks = 0;
133 Int_t nTracks = 0;
134 Int_t nGhosts = 0;
135 Int_t nClones = 0;
136 Int_t nAll = 0;
137 Int_t nAcc = 0;
138 Int_t nRecAll = 0;
139 Int_t nPrim = 0;
140 Int_t nRecPrim = 0;
141 Int_t nFast = 0;
142 Int_t nRecFast = 0;
143 Int_t nFastLong = 0;
144 Int_t nRecFastLong = 0;
145 Int_t nSec = 0;
146 Int_t nRecSec = 0;
147 TVector3 momentum;
148 TVector3 vertex;
149
150 // check consistency
151 //assert(fGlobalTracks->GetEntriesFast() == fGlobalTrackMatches->GetEntriesFast());
152 assert(fTrdTracks->GetEntriesFast() == fTrdTrackMatches->GetEntriesFast());
153 assert(fStsTracks->GetEntriesFast() == fStsTrackMatches->GetEntriesFast());
154
155 std::vector<CbmLink> events = fTimeSlice->GetMatch().GetLinks();
156 std::sort(events.begin(), events.end());
157
158 {
159 fMcTrackInfoMap.clear();
160 McTrackInfo info;
161 info.fIsAccepted = false;
162 info.fNtimesReconstructed = 0;
163 info.fIsPrimary = false;
164 info.fIsFast = false;
165 info.fIsLong = false;
166
167 std::cout << "MC events in slice : " << events.size() << std::endl;
168
169 for (uint iLink = 0; iLink < events.size(); iLink++) {
170 CbmLink link = events[iLink];
171 Int_t nMCTracks = fMCTracks->Size(link);
172 std::cout << "MC event " << iLink << " n mc tracks " << nMCTracks << std::endl;
173 for (Int_t iTr = 0; iTr < nMCTracks; iTr++) {
174 link.SetIndex(iTr);
175 fMcTrackInfoMap.insert({link, info});
176 }
177 }
178 }
179
180 // Fill hit and track maps
181 FillHitMap();
182 FillTrackMatchMap(nTracks, nGhosts, nClones);
183
184 int nMcTracks = fMcTrackInfoMap.size();
185
186 const auto* pTrdIfs = cbm::RecoSetupManager::Instance()->GetSetup().GetTrd();
187
188 // Loop over MCTracks
189 int iMcTrack = 0;
190 for (auto itTrack = fMcTrackInfoMap.begin(); itTrack != fMcTrackInfoMap.end(); ++itTrack, ++iMcTrack) {
191 const CbmLink& link = itTrack->first;
192 McTrackInfo& info = itTrack->second;
193 CbmMCTrack* mcTrack = dynamic_cast<CbmMCTrack*>(fMCTracks->Get(link));
194 assert(mcTrack);
195
196 nAll++;
197
198 // Check origin of MCTrack
199 // TODO: Track origin should rather be compared to MC event vertex
200 // But that is not available from MCDataManager
201 mcTrack->GetStartVertex(vertex);
202 if (TMath::Abs(vertex.Z() - fTargetPos.Z()) < 1.) {
203 info.fIsPrimary = true;
204 }
205
206 // Get momentum
207 mcTrack->GetMomentum(momentum);
208 info.fP = momentum.Mag();
209 info.fPt = momentum.Pt();
210 info.fPdg = mcTrack->GetPdgCode();
211 info.fY = mcTrack->GetRapidity() - fYCM;
212 info.fIsFast = (info.fP > 0.1);
213
214 // Continue only for reconstructible tracks
215
216 Int_t nStations = info.fHitMap.size();
217
218 int nContStations = 0; // Number of continious stations
219 {
220 int istaprev = -100;
221 int len = 0;
222 for (auto itSta = info.fHitMap.begin(); itSta != info.fHitMap.end(); itSta++) {
223 if (itSta->first == istaprev + 1) {
224 len++;
225 }
226 else {
227 len = 1;
228 }
229 if (nContStations < len) {
230 nContStations = len;
231 }
232 istaprev = itSta->first;
233 }
234 }
235
236 info.fIsLong = (nContStations >= fTrdNstations);
237
238 if (nStations < fMinStations) continue; // Too few stations
239 if (nContStations < fMinStations) continue; // Too few stations
240
241 info.fIsAccepted = true;
242 nAcc++;
243
244 if (info.fIsPrimary) {
245 nPrim++;
246 }
247 else {
248 nSec++;
249 }
250
251 if (info.fIsFast) {
252 nFast++;
253 }
254
255 Bool_t isFastLong = (info.fIsFast && info.fIsLong);
256 if (isFastLong) {
257 nFastLong++;
258 }
259
260 // Fill histograms for reconstructible tracks
261
262 fhPtAccAll->Fill(info.fPt);
263 fhNpAccAll->Fill(Double_t(nStations));
264 if (info.fIsPrimary) {
265 fhPtAccPrim->Fill(info.fPt);
266 fhPidPtY["prmY"][Pdg2Idx(info.fPdg)]->Fill(info.fY, info.fPt);
267 fhNpAccPrim->Fill(Double_t(nStations));
268 }
269 else {
270 fhPtAccSec->Fill(info.fPt);
271 fhPidPtY["secY"][Pdg2Idx(info.fPdg)]->Fill(info.fY, info.fPt);
272 fhNpAccSec->Fill(Double_t(nStations));
273 fhZAccSec->Fill(vertex.Z());
274 }
275
276 // Get matched GlobalTrack
277 Int_t globalTrackId = info.fGlobalTrackMatch;
278 Int_t trdTrackId = info.fTrdTrackMatch;
279 Double_t quali = info.fQuali;
280 // Bool_t isRec = kFALSE;
281 if (globalTrackId < 0) continue;
282 if (trdTrackId < 0) continue;
283 CbmGlobalTrack* globalTrack = (CbmGlobalTrack*) fGlobalTracks->At(globalTrackId);
284
285 CbmTrdTrack* trdTrack = (CbmTrdTrack*) fTrdTracks->At(trdTrackId);
286 assert(trdTrack);
287 assert(quali >= fQuota);
288 CbmTrackMatchNew* match = (CbmTrackMatchNew*) fTrdTrackMatches->At(trdTrackId);
289 assert(match);
290 Int_t nTrue = match->GetNofTrueHits();
291 Int_t nWrong = match->GetNofWrongHits();
292 //Int_t nFake = match->GetNofFakeHits();
293 Int_t nFake = 0;
294 Int_t nAllHits = trdTrack->GetNofHits();
295 assert(nTrue + nWrong + nFake == nAllHits);
296
297 int nStsHits = 0;
298 {
299 Int_t stsTrackId = info.fStsTrackMatch;
300 if (stsTrackId >= 0) {
301 CbmStsTrack* stsTrack = (CbmStsTrack*) fStsTracks->At(stsTrackId);
302 assert(stsTrack);
303 nStsHits = stsTrack->GetTotalNofHits();
304 }
305 }
306
307 double qp = globalTrack->GetParamFirst()->GetQp();
308 //double q = (qp >= 1.) ? 1. : -1.;
309 double p = (fabs(qp) > 1. / 1000.) ? 1. / fabs(qp) : 1000.;
310 double tx = globalTrack->GetParamFirst()->GetTx();
311 double ty = globalTrack->GetParamFirst()->GetTy();
312 double pt = sqrt((tx * tx + ty * ty) / (1. + tx * tx + ty * ty)) * p;
313
314 // Verbose output
315 LOG(debug1) << GetName() << ": MCTrack " << iMcTrack << ", stations " << nStations << ", hits " << nAllHits
316 << ", true hits " << nTrue;
317
318 // Fill histograms for reconstructed tracks
319 nRecAll++;
320 fhPtRecAll->Fill(info.fPt);
321 fhNpRecAll->Fill(Double_t(nAllHits));
322 //fhPidXY[Pdg2Idx[info.fPdg]]->Fill(info.fY, info.fPt);
323 if (info.fIsPrimary) {
324 nRecPrim++;
325 fhPtRecPrim->Fill(info.fPt);
326 fhPidPtY["prmE"][Pdg2Idx(info.fPdg)]->Fill(info.fY, info.fPt);
327 fhNpRecPrim->Fill(Double_t(nAllHits));
328 if (info.fPt > 0.001) {
329 fhPtResPrim->Fill((pt / info.fPt - 1.));
330 }
331 if (info.fP > 0.001) {
332 double dp = p / info.fP - 1.;
333 fhPResPrim->Fill(dp);
334 if (nStsHits == 0) {
335 fhPResPrimSts0->Fill(dp);
336 }
337 if (nStsHits == 1) {
338 fhPResPrimSts1->Fill(dp);
339 }
340 if (nStsHits == 2) {
341 fhPResPrimSts2->Fill(dp);
342 }
343 if (nStsHits >= 3) {
344 fhPResPrimSts3->Fill(dp);
345 }
346 }
347 }
348 else {
349 nRecSec++;
350 fhPtRecSec->Fill(info.fPt);
351 fhPidPtY["secE"][Pdg2Idx(info.fPdg)]->Fill(info.fY, info.fPt);
352 fhNpRecSec->Fill(Double_t(nAllHits));
353 fhZRecSec->Fill(vertex.Z());
354 }
355 if (info.fIsFast) nRecFast++;
356 if (isFastLong) nRecFastLong++;
357
358 } // Loop over MCTracks
359
360 // Loop over MC points
361
362 for (uint iLink = 0; iLink < events.size(); iLink++) {
363 CbmLink link = events[iLink];
364 Int_t nMcPoints = fTrdPoints->Size(link);
365 std::cout << "MC event " << iLink << " n mc points " << nMcPoints << std::endl;
366 for (Int_t ip = 0; ip < nMcPoints; ip++) {
367 link.SetIndex(ip);
368 CbmTrdPoint* pt = (CbmTrdPoint*) fTrdPoints->Get(link);
369 assert(pt);
370 Int_t station = pTrdIfs->GetTrackingStationId(pt->GetDetectorID());
371 link.SetIndex(pt->GetTrackID());
372 McTrackInfo& info = getMcTrackInfo(link);
373 if (info.fIsAccepted) {
374 fhStationEffXY[station].Fill(pt->GetX(), pt->GetY(), (info.fNtimesReconstructed > 0) ? 100. : 0.);
375 }
376 }
377 }
378
379 // Calculate efficiencies
380 Double_t effAll = Double_t(nRecAll) / Double_t(nAcc);
381 Double_t effPrim = Double_t(nRecPrim) / Double_t(nPrim);
382 Double_t effFast = Double_t(nRecFast) / Double_t(nFast);
383 Double_t effFastLong = Double_t(nRecFastLong) / Double_t(nFastLong);
384 Double_t effSec = Double_t(nRecSec) / Double_t(nSec);
385
386 fTimer.Stop();
387
388
389 // Event summary
390 LOG(info) << "+ " << setw(20) << GetName() << ": Event " << setw(6) << right << fNEvents << ", real time " << fixed
391 << setprecision(6) << fTimer.RealTime() << " s, MC tracks: all " << nMcTracks << ", acc. " << nAcc
392 << ", rec. " << nRecAll << ", eff. " << setprecision(2) << 100. * effAll << " %";
393 if (fair::Logger::Logging(fair::Severity::debug)) {
394 LOG(debug) << "---------- CbmTrackingTrdQa : Event summary ------------";
395 LOG(debug) << "MCTracks : " << nAll << ", reconstructible: " << nAcc << ", reconstructed: " << nRecAll;
396 LOG(debug) << "Vertex : reconstructible: " << nPrim << ", reconstructed: " << nRecPrim << ", efficiency "
397 << effPrim * 100. << "%";
398 LOG(debug) << "Fast : reconstructible: " << nFast << ", reconstructed: " << nRecFast << ", efficiency "
399 << effFast * 100. << "%";
400 LOG(debug) << "Fast long : reconstructible: " << nFastLong << ", reconstructed: " << nRecFastLong << ", efficiency "
401 << effFastLong * 100. << "%";
402 LOG(debug) << "Non-vertex : reconstructible: " << nSec << ", reconstructed: " << nRecSec << ", efficiency "
403 << effSec * 100. << "%";
404 LOG(debug) << "TrdTracks " << nTracks << ", ghosts " << nGhosts << ", clones " << nClones;
405 LOG(debug) << "-----------------------------------------------------------\n";
406 }
407
408
409 // Increase counters
410 fNAll += nAll;
411 fNAccAll += nAcc;
412 fNAccPrim += nPrim;
413 fNAccFast += nFast;
414 fNAccFastLong += nFastLong;
415 fNAccSec += nSec;
416 fNRecAll += nRecAll;
417 fNRecPrim += nRecPrim;
418 fNRecFast += nRecFast;
419 fNRecFastLong += nRecFastLong;
420 fNRecSec += nRecSec;
421 fNGhosts += nGhosts;
422 fNClones += nClones;
423 fNEvents++;
424 fTime += fTimer.RealTime();
425}
426// -------------------------------------------------------------------------
427
428
429// ----- Public method SetParContainers --------------------------------
431// -------------------------------------------------------------------------
432
433
434// ----- Public method Init --------------------------------------------
436{
437
438 LOG(info) << "\n\n====================================================";
439 LOG(info) << GetName() << ": Initialising...";
440
441 fManager = FairRootManager::Instance();
442 assert(fManager);
443
444 fMcManager = dynamic_cast<CbmMCDataManager*>(fManager->GetObject("MCDataManager"));
445
446 assert(fMcManager);
447
448 fTimeSlice = static_cast<CbmTimeSlice*>(fManager->GetObject("TimeSlice."));
449
450 if (fTimeSlice == nullptr) {
451 LOG(fatal) << "CbmTrackingTrdQa: No time slice object";
452 }
453
454 if (fMcManager) {
455 fMCTracks = fMcManager->InitBranch("MCTrack");
456 fTrdPoints = fMcManager->InitBranch("TrdPoint");
457 }
458
459 assert(fMCTracks);
460 assert(fTrdPoints);
461
462 // Get the geometry
463 InitStatus geoStatus = GetGeometry();
464 if (geoStatus != kSUCCESS) {
465 LOG(error) << GetName() << "::Init: Error in reading geometry!";
466 return geoStatus;
467 }
468
469
470 // TRD
471
472 // Get TrdHit array
473 fTrdHits = (TClonesArray*) fManager->GetObject("TrdHit");
474 assert(fTrdHits);
475
476 // Get TrdHitMatch array
477 fTrdHitMatch = (TClonesArray*) fManager->GetObject("TrdHitMatch");
478 assert(fTrdHitMatch);
479
480 // Get GlobalTrack array
481 fGlobalTracks = (TClonesArray*) fManager->GetObject("GlobalTrack");
482 assert(fGlobalTracks);
483
484 // Get GlobalTrackMatch array
485 //fGlobalTrackMatches = (TClonesArray*) fManager->GetObject("GlobalTrackMatch");
486 //assert(fGlobalTrackMatches);
487
488 // Get TrdTrack array
489 fTrdTracks = (TClonesArray*) fManager->GetObject("TrdTrack");
490 assert(fTrdTracks);
491
492 // Get TrdTrackMatch array
493 fTrdTrackMatches = (TClonesArray*) fManager->GetObject("TrdTrackMatch");
494 assert(fTrdTrackMatches);
495
496 // Get StsTrack array
497 fStsTracks = (TClonesArray*) fManager->GetObject("StsTrack");
498 assert(fStsTracks);
499
500 // Get StsTrackMatch array
501 fStsTrackMatches = (TClonesArray*) fManager->GetObject("StsTrackMatch");
502 assert(fStsTrackMatches);
503
504 // Create histograms
505 CreateHistos();
506 Reset();
507
508 // Output
509 LOG(info) << " Number of Trd stations : " << fTrdNstations;
510 LOG(info) << " Target position ( " << fTargetPos.X() << ", " << fTargetPos.Y() << ", " << fTargetPos.Z() << ") cm";
511 LOG(info) << " Minimum number of Trd stations : " << fMinStations;
512 LOG(info) << " Matching quota : " << fQuota;
513 LOG(info) << "====================================================";
514
515 return geoStatus;
516}
517// -------------------------------------------------------------------------
518
519
520// ----- Public method ReInit ------------------------------------------
522{
523
524 LOG(info) << "\n\n====================================================";
525 LOG(info) << GetName() << ": Re-initialising...";
526
527 // Get the geometry of target and Trd
528 InitStatus geoStatus = GetGeometry();
529 if (geoStatus != kSUCCESS) {
530 LOG(error) << GetName() << "::Init: Error in reading geometry!";
531 return geoStatus;
532 }
533
534 // --- Screen log
535 LOG(info) << " Number of TRD stations : " << fTrdNstations;
536 LOG(info) << " Target position ( " << fTargetPos.X() << ", " << fTargetPos.Y() << ", " << fTargetPos.Z() << ") cm";
537 LOG(info) << " Minimum number of TRD stations : " << fMinStations;
538 LOG(info) << " Matching quota : " << fQuota;
539 LOG(info) << "====================================================";
540
541 return geoStatus;
542}
543// -------------------------------------------------------------------------
544
545
546// ----- Private method Finish -----------------------------------------
548{
549
550 // Divide histograms for efficiency calculation
558 for (int idx(0); idx < fgkNpdg; idx++) {
559 fhPidPtY["allY"][idx]->Add(fhPidPtY["prmY"][idx]);
560 fhPidPtY["allY"][idx]->Add(fhPidPtY["secY"][idx]);
561 fhPidPtY["allE"][idx]->Add(fhPidPtY["prmE"][idx]);
562 fhPidPtY["allE"][idx]->Add(fhPidPtY["secE"][idx]);
563 DivideHistos(fhPidPtY["prmE"][idx], fhPidPtY["prmY"][idx], fhPidPtY["prmE"][idx], "2D");
564 DivideHistos(fhPidPtY["secE"][idx], fhPidPtY["secY"][idx], fhPidPtY["secE"][idx], "2D");
565 DivideHistos(fhPidPtY["allE"][idx], fhPidPtY["allY"][idx], fhPidPtY["allE"][idx], "2D");
566 }
567 // Normalise histos for clones and ghosts to one event
568 if (fNEvents) {
569 fhNhClones->Scale(1. / Double_t(fNEvents));
570 fhNhGhosts->Scale(1. / Double_t(fNEvents));
571 }
572
573 // Calculate integrated efficiencies and rates
574 Double_t effAll = Double_t(fNRecAll) / Double_t(fNAccAll);
575 Double_t effPrim = Double_t(fNRecPrim) / Double_t(fNAccPrim);
576 Double_t effFast = Double_t(fNRecFast) / Double_t(fNAccFast);
577 Double_t effFastLong = Double_t(fNRecFastLong) / Double_t(fNAccFastLong);
578 Double_t effSec = Double_t(fNRecSec) / Double_t(fNAccSec);
579 Double_t rateGhosts = Double_t(fNGhosts) / Double_t(fNRecAll);
580 Double_t rateClones = Double_t(fNClones) / Double_t(fNRecAll);
581
582 // Run summary to screen
583 std::cout << std::endl;
584 LOG(info) << "=====================================";
585 LOG(info) << fName << ": Run summary ";
586 LOG(info) << "Events processed : " << fNEvents << setprecision(2);
587 LOG(info) << "Eff. all tracks : " << effAll * 100 << " % (" << fNRecAll << "/" << fNAccAll << ")";
588 LOG(info) << "Eff. vertex tracks : " << effPrim * 100 << " % (" << fNRecPrim << "/" << fNAccPrim << ")";
589 LOG(info) << "Eff. fast tracks : " << effFast * 100 << " % (" << fNRecFast << "/" << fNAccFast << ")";
590 LOG(info) << "Eff. fast long tracks : " << effFastLong * 100 << " % (" << fNRecFastLong << "/" << fNAccFastLong
591 << ")";
592 LOG(info) << "Eff. secondary tracks : " << effSec * 100 << " % (" << fNRecSec << "/" << fNAccSec << ")";
593 LOG(info) << "Ghost rate : " << rateGhosts * 100 << " % (" << fNGhosts << "/" << fNRecAll << ")";
594 LOG(info) << "Clone rate : " << rateClones * 100 << " % (" << fNClones << "/" << fNRecAll << ")";
595 LOG(info) << "mc tracks/event " << fNAll / fNEvents << " accepted " << fNRecAll / fNEvents;
596 LOG(info) << "Time per event : " << setprecision(6) << fTime / Double_t(fNEvents) << " s";
597
598
599 LOG(info) << "=====================================";
600
601 // Write histos to output
602 /*
603 gDirectory->mkdir("CbmTrackingTrdQa");
604 gDirectory->cd("CbmTrackingTrdQa");
605 TIter next(fHistoList);
606 while (TH1* histo = ((TH1*) next()))
607 histo->Write();
608 gDirectory->cd("..");
609*/
610 FairSink* sink = FairRootManager::Instance()->GetSink();
611 sink->WriteObject(&fOutFolder, nullptr);
612}
613// -------------------------------------------------------------------------
614
615
616// ----- Private method GetGeometry ------------------------------------
618{
619 // Get target geometry
621
622 // Check the reco setup manager
623 if (!cbm::RecoSetupManager::Instance()->IsInitialized()) {
624 return kFATAL;
625 }
626
627 const auto* pTrdIfs = cbm::RecoSetupManager::Instance()->GetSetup().GetTrd();
628 if (!pTrdIfs) {
629 LOG(error) << fName << ": TRD reco unit is not in the setup";
630 return kFATAL;
631 }
632
633 fTrdNstations = pTrdIfs->GetNofTrackingStations();
634
635 return kSUCCESS;
636}
637// -------------------------------------------------------------------------
638
639
640// ----- Get target node -----------------------------------------------
642{
643
644 TGeoNode* target = NULL;
645
646 gGeoManager->CdTop();
647 TGeoNode* cave = gGeoManager->GetCurrentNode();
648 for (Int_t iNode1 = 0; iNode1 < cave->GetNdaughters(); iNode1++) {
649 TString name = cave->GetDaughter(iNode1)->GetName();
650 if (name.Contains("pipe", TString::kIgnoreCase)) {
651 LOG(debug) << "Found pipe node " << name;
652 gGeoManager->CdDown(iNode1);
653 break;
654 }
655 }
656 for (Int_t iNode2 = 0; iNode2 < gGeoManager->GetCurrentNode()->GetNdaughters(); iNode2++) {
657 TString name = gGeoManager->GetCurrentNode()->GetDaughter(iNode2)->GetName();
658 if (name.Contains("pipevac1", TString::kIgnoreCase)) {
659 LOG(debug) << "Found vacuum node " << name;
660 gGeoManager->CdDown(iNode2);
661 break;
662 }
663 }
664 for (Int_t iNode3 = 0; iNode3 < gGeoManager->GetCurrentNode()->GetNdaughters(); iNode3++) {
665 TString name = gGeoManager->GetCurrentNode()->GetDaughter(iNode3)->GetName();
666 if (name.Contains("target", TString::kIgnoreCase)) {
667 LOG(debug) << "Found target node " << name;
668 gGeoManager->CdDown(iNode3);
669 target = gGeoManager->GetCurrentNode();
670 break;
671 }
672 }
673 if (!target) {
674 fTargetPos[0] = 0.;
675 fTargetPos[1] = 0.;
676 fTargetPos[2] = 0.;
677 }
678 else {
679 TGeoHMatrix* glbMatrix = gGeoManager->GetCurrentMatrix();
680 Double_t* pos = glbMatrix->GetTranslation();
681 fTargetPos[0] = pos[0];
682 fTargetPos[1] = pos[1];
683 fTargetPos[2] = pos[2];
684 }
685
686 gGeoManager->CdTop();
687}
688// -------------------------------------------------------------------------
689
690
691// ----- Private method CreateHistos -----------------------------------
693{
694
695 fOutFolder.Clear();
696
697 // Histogram list
698 fHistoList = new TList();
699
700 // Momentum distributions
701 Double_t minPt = 0.;
702 Double_t maxPt = 2.;
703 Int_t nBinsPt = 40;
704
705 fhStationEffXY.clear();
706
707 for (int i = 0; i < fTrdNstations; i++) {
708 //auto trdInterface = CbmTrdTrackingInterface::Instance();
709 double dx = 150; //trdInterface->GetXmax(i);
710 double dy = 150; //trdInterface->GetYmax(i);
711 fhStationEffXY.emplace_back(Form("fhStationEffXY%i", i), Form("Efficiency XY: Station %i;X [cm];Y [cm]", i), 300,
712 -dx, dx, 300, -dy, dy);
713 fhStationEffXY[i].SetDirectory(0);
714 fhStationEffXY[i].SetOptStat(10);
715 fhStationEffXY[i].GetYaxis()->SetTitleOffset(1.4);
716 }
717
718 for (int i = 0; i < fTrdNstations; i++) {
719 fHistoList->Add(&fhStationEffXY[i]);
720 }
721
722 fhPtAccAll = new TH1F("hPtAccAll", "all reconstructable tracks", nBinsPt, minPt, maxPt);
723 fhPtRecAll = new TH1F("hPtRecAll", "all reconstructed tracks", nBinsPt, minPt, maxPt);
724 fhPtEffAll = new TH1F("hPtEffAll", "efficiency all tracks", nBinsPt, minPt, maxPt);
725 fhPtAccPrim = new TH1F("hPtAccPrim", "reconstructable vertex tracks", nBinsPt, minPt, maxPt);
726 fhPtRecPrim = new TH1F("hPtRecPrim", "reconstructed vertex tracks", nBinsPt, minPt, maxPt);
727 fhPtEffPrim = new TH1F("hPtEffPrim", "efficiency vertex tracks", nBinsPt, minPt, maxPt);
728 fhPtAccSec = new TH1F("hPtAccSec", "reconstructable non-vertex tracks", nBinsPt, minPt, maxPt);
729 fhPtRecSec = new TH1F("hPtRecSec", "reconstructed non-vertex tracks", nBinsPt, minPt, maxPt);
730 fhPtEffSec = new TH1F("hPtEffSec", "efficiency non-vertex tracks", nBinsPt, minPt, maxPt);
740
741 const char* pltTitle[] = {"Yield", "Yield", "Yield", "Efficiency", "Efficiency", "Efficiency"};
742 const char* pltLab[] = {"yield", "yield", "yield", "#epsilon (%)", "#epsilon (%)", "#epsilon (%)"};
743 const char* vxTyp[] = {"allY", "prmY", "secY", "allE", "prmE", "secE"};
744 for (int ivx(0); ivx < 6; ivx++) {
745 for (int ipid(0); ipid < fgkNpdg; ipid++) {
746 fhPidPtY[vxTyp[ivx]][ipid] = new TH2F(
747 Form("hPtY_%s%s", vxTyp[ivx], Idx2Symb(ipid)),
748 Form("%s %s(%s); y - y_{cm}; p_{T} (GeV/c); %s", pltTitle[ivx], Idx2Name(ipid), vxTyp[ivx], pltLab[ivx]), 50,
749 -2.5, 2.5, nBinsPt, minPt, maxPt);
750 fHistoList->Add(fhPidPtY[vxTyp[ivx]][ipid]);
751 }
752 }
753 // Number-of-points distributions
754 Double_t minNp = -0.5;
755 Double_t maxNp = 15.5;
756 Int_t nBinsNp = 16;
757 fhNpAccAll = new TH1F("hNpAccAll", "all reconstructable tracks", nBinsNp, minNp, maxNp);
758 fhNpRecAll = new TH1F("hNpRecAll", "all reconstructed tracks", nBinsNp, minNp, maxNp);
759 fhNpEffAll = new TH1F("hNpEffAll", "efficiency all tracks", nBinsNp, minNp, maxNp);
760 fhNpAccPrim = new TH1F("hNpAccPrim", "reconstructable vertex tracks", nBinsNp, minNp, maxNp);
761 fhNpRecPrim = new TH1F("hNpRecPrim", "reconstructed vertex tracks", nBinsNp, minNp, maxNp);
762 fhNpEffPrim = new TH1F("hNpEffPrim", "efficiency vertex tracks", nBinsNp, minNp, maxNp);
763 fhNpAccSec = new TH1F("hNpAccSec", "reconstructable non-vertex tracks", nBinsNp, minNp, maxNp);
764 fhNpRecSec = new TH1F("hNpRecSec", "reconstructed non-vertex tracks", nBinsNp, minNp, maxNp);
765 fhNpEffSec = new TH1F("hNpEffSec", "efficiency non-vertex tracks", nBinsNp, minNp, maxNp);
775
776 // z(vertex) distributions
777 Double_t minZ = 0.;
778 Double_t maxZ = 50.;
779 Int_t nBinsZ = 50;
780 fhZAccSec = new TH1F("hZAccSec", "reconstructable non-vertex tracks", nBinsZ, minZ, maxZ);
781 fhZRecSec = new TH1F("hZRecSecl", "reconstructed non-vertex tracks", nBinsZ, minZ, maxZ);
782 fhZEffSec = new TH1F("hZEffRec", "efficiency non-vertex tracks", nBinsZ, minZ, maxZ);
783 fHistoList->Add(fhZAccSec);
784 fHistoList->Add(fhZRecSec);
785 fHistoList->Add(fhZEffSec);
786
787 // Number-of-hit distributions
788 fhNhClones = new TH1F("hNhClones", "number of hits for clones", nBinsNp, minNp, maxNp);
789 fhNhGhosts = new TH1F("hNhGhosts", "number of hits for ghosts", nBinsNp, minNp, maxNp);
792
793 fhPtResPrim = new TH1F("hPtPrim", "Resolution Pt Primaries [100%]", 100, -1., 1.);
795
796 fhPResPrim = new TH1F("hPPrim", "Resolution P Primaries [100%]", 100, -1., 1.);
798
799 fhPResPrimSts0 = new TH1F("hPPrimSts0", "Resolution P Primaries [100%], No Sts hits", 100, -1., 1.);
801
802 fhPResPrimSts1 = new TH1F("hPPrimSts1", "Resolution P Primaries [100%], 1 Sts hit", 100, -1., 1.);
804
805 fhPResPrimSts2 = new TH1F("hPPrimSts2", "Resolution P Primaries [100%], 2 Sts hits", 100, -1., 1.);
807
808 fhPResPrimSts3 = new TH1F("hPPrimSts3", "Resolution P Primaries [100%], >=3 Sts hits", 100, -1., 1.);
810
811 TIter next(fHistoList);
812 while (TH1* histo = ((TH1*) next())) {
813 fOutFolder.Add(histo);
814 }
815}
816// -------------------------------------------------------------------------
817
818
819// ----- Private method Reset ------------------------------------------
821{
822
823 TIter next(fHistoList);
824 while (TH1* histo = ((TH1*) next()))
825 histo->Reset();
826
829 fNGhosts = fNClones = fNEvents = 0;
830}
831// -------------------------------------------------------------------------
832
833
834// ----- Private method FillHitMap -------------------------------------
836{
837
838 // --- Fill hit map ( mcTrack -> ( station -> number of hits ) )
839
840 // pocess Trd hits
841 for (Int_t iHit = 0; iHit < fTrdHits->GetEntriesFast(); iHit++) {
842 CbmTrdHit* hit = (CbmTrdHit*) fTrdHits->At(iHit);
843 assert(hit);
844
845 if ((int) hit->GetClassType() != 1) {
846 // skip TRD-1D hit
847 continue;
848 }
849
850 const auto* pTrdIfs = cbm::RecoSetupManager::Instance()->GetSetup().GetTrd();
851 Int_t station = pTrdIfs->GetTrackingStationId(hit->GetAddress());
852
853 // carefully check the hit match by looking at the strips from both sides
854
855 const CbmMatch* match = dynamic_cast<const CbmMatch*>(fTrdHitMatch->At(iHit));
856 assert(match);
857 if (match->GetNofLinks() <= 0) continue;
858
859 CbmLink link = match->GetMatchedLink();
860
861 CbmTrdPoint* pt = (CbmTrdPoint*) fTrdPoints->Get(link);
862 assert(pt);
863 assert(pTrdIfs->GetTrackingStationId(pt->GetDetectorID()) == station);
864
865 link.SetIndex(pt->GetTrackID());
866 McTrackInfo& info = getMcTrackInfo(link);
867 info.fHitMap[station]++;
868 }
869 LOG(debug) << GetName() << ": Filled hit map from " << fTrdHits->GetEntriesFast() << " Trd hits";
870}
871// -------------------------------------------------------------------------
872
873
874// ------ Private method FillMatchMap ----------------------------------
876{
877 // Clear matching maps
878 for (auto it = fMcTrackInfoMap.begin(); it != fMcTrackInfoMap.end(); ++it) {
879 McTrackInfo& info = it->second;
880 info.fGlobalTrackMatch = -1;
881 info.fQuali = 0.;
882 info.fMatchedNHitsAll = 0;
883 info.fMatchedNHitsTrue = 0;
884 info.fNtimesReconstructed = 0;
885 }
886
887 // Loop over GlobalTracks. Check matched MCtrack and fill maps.
888 nGhosts = 0;
889 nClones = 0;
890 nRec = 0;
891
892 //assert(fGlobalTrackMatches);
893
894 for (Int_t iGlobalTrack = 0; iGlobalTrack < fGlobalTracks->GetEntriesFast(); iGlobalTrack++) {
895
896 // --- GlobalTrack
897 CbmGlobalTrack* globalTrack = (CbmGlobalTrack*) fGlobalTracks->At(iGlobalTrack);
898 assert(globalTrack);
899
900 // --- TrackMatch
901
902 //assert(iGlobalTrack >= 0 && iGlobalTrack < fGlobalTrackMatches->GetEntriesFast());
903 //CbmTrackMatchNew* globalMatch = (CbmTrackMatchNew*) fGlobalTrackMatches->At(iGlobalTrack);
904 //assert(globalMatch);
905
906 int iTrdTrack = globalTrack->GetTrdTrackIndex();
907
908 if (iTrdTrack < 0) continue;
909 CbmTrdTrack* trdTrack = (CbmTrdTrack*) fTrdTracks->At(iTrdTrack);
910 assert(trdTrack);
911 nRec++;
912
913 int iStsTrack = globalTrack->GetStsTrackIndex();
914
915 // --- TrackMatch
916
917 assert(iTrdTrack >= 0 && iTrdTrack < fTrdTrackMatches->GetEntriesFast());
918 CbmTrackMatchNew* trdMatch = (CbmTrackMatchNew*) fTrdTrackMatches->At(iTrdTrack);
919 assert(trdMatch);
920
921 Int_t nHits = trdTrack->GetNofHits();
922 Int_t nTrue = trdMatch->GetNofTrueHits();
923
924 // Check matching criterion (quota)
925 Double_t quali = Double_t(nTrue) / Double_t(nHits);
926
927 // Quality isn't good, it's a ghost
928
929 if (quali < fQuota) {
930 fhNhGhosts->Fill(nHits);
931 nGhosts++;
932 continue;
933 }
934
935 // Quality is good
936
937 // --- Matched MCTrack
938 assert(trdMatch->GetNofLinks() > 0);
939 CbmLink link = trdMatch->GetMatchedLink();
940 assert(link.GetIndex() >= 0);
941 link.SetFile(0); //??
942
943 McTrackInfo& info = getMcTrackInfo(link);
944
945 // previous match is better, this track is a clone
946 if ((quali < info.fQuali) || ((quali == info.fQuali) && (nTrue < info.fMatchedNHitsTrue))) {
947 if (info.fIsAccepted) {
948 fhNhClones->Fill(nHits);
949 nClones++;
950 }
951 continue;
952 }
953
954 // this track is better than the old one
955 if (info.fMatchedNHitsAll > 0) {
956 if (info.fIsAccepted) {
957 fhNhClones->Fill(info.fMatchedNHitsAll);
958 nClones++;
959 }
960 }
961 info.fGlobalTrackMatch = iGlobalTrack;
962 info.fTrdTrackMatch = iTrdTrack;
963 info.fStsTrackMatch = iStsTrack;
964 info.fQuali = quali;
965 info.fMatchedNHitsAll = nHits;
966 info.fMatchedNHitsTrue = nTrue;
968 } // Loop over GlobalTracks
969
970 LOG(debug) << GetName() << ": Filled match map for " << nRec << " Trd tracks. Ghosts " << nGhosts << " Clones "
971 << nClones;
972}
973// -------------------------------------------------------------------------
974
975
976// ----- Private method DivideHistos -----------------------------------
977void CbmTrackingTrdQa::DivideHistos(TH1* histo1, TH1* histo2, TH1* histo3, Option_t* opt)
978{
979
980 if (!histo1 || !histo2 || !histo3) {
981 LOG(fatal) << GetName() << "::DivideHistos: "
982 << "NULL histogram pointer";
983 }
984
985 Int_t nBins = histo1->GetNbinsX();
986 if (histo2->GetNbinsX() != nBins || histo3->GetNbinsX() != nBins) {
987 LOG(error) << GetName() << "::DivideHistos: "
988 << "Different bin numbers in histos";
989 LOG(error) << histo1->GetName() << " " << histo1->GetNbinsX();
990 LOG(error) << histo2->GetName() << " " << histo2->GetNbinsX();
991 LOG(error) << histo3->GetName() << " " << histo3->GetNbinsX();
992 return;
993 }
994 if (strcmp(opt, "2D") == 0) {
995 Int_t nBinsY = histo1->GetNbinsY();
996 if (histo2->GetNbinsY() != nBinsY || histo3->GetNbinsY() != nBinsY) {
997 LOG(error) << GetName() << "::DivideHistos: "
998 << "Different bin numbers in histos";
999 LOG(error) << histo1->GetName() << " " << histo1->GetNbinsY();
1000 LOG(error) << histo2->GetName() << " " << histo2->GetNbinsY();
1001 LOG(error) << histo3->GetName() << " " << histo3->GetNbinsY();
1002 return;
1003 }
1004 nBins *= nBinsY;
1005 }
1006
1007 Double_t c1, c2, c3, ce;
1008 for (Int_t iBin = 0; iBin < nBins; iBin++) {
1009 c1 = histo1->GetBinContent(iBin);
1010 c2 = histo2->GetBinContent(iBin);
1011 if (c2 != 0.) {
1012 c3 = c1 / c2;
1013 Double_t c4 = (c3 * (1. - c3) / c2);
1014 if (c4 >= 0.) {
1015 ce = TMath::Sqrt(c3 * (1. - c3) / c2);
1016 }
1017 else {
1018 ce = 0;
1019 }
1020 }
1021 else {
1022 c3 = 0.;
1023 ce = 0.;
1024 }
1025 histo3->SetBinContent(iBin, 1.e2 * c3);
1026 histo3->SetBinError(iBin, 1.e2 * ce);
1027 }
1028}
1029// -------------------------------------------------------------------------
1030
1031
ClassImp(CbmConverterManager)
Int_t nStsHits
Int_t nMCTracks
A manager for setup representation in CBM reconstruction.
Data class for STS tracks.
Class for hits in TRD detector.
friend fvec sqrt(const fvec &a)
int Int_t
bool Bool_t
int32_t GetStsTrackIndex() const
const FairTrackParam * GetParamFirst() const
int32_t GetTrdTrackIndex() const
int32_t GetAddress() const
Definition CbmHit.h:77
Task class creating and managing CbmMCDataArray objects.
void GetMomentum(TVector3 &momentum) const
Definition CbmMCTrack.h:169
int32_t GetPdgCode() const
Definition CbmMCTrack.h:67
double GetRapidity() const
void GetStartVertex(TVector3 &vertex) const
Definition CbmMCTrack.h:175
int32_t GetNofLinks() const
Definition CbmMatch.h:42
const CbmLink & GetMatchedLink() const
Definition CbmMatch.h:41
virtual int32_t GetTotalNofHits() const
Definition CbmStsTrack.h:81
Bookkeeping of time-slice content.
int32_t GetNofWrongHits() const
int32_t GetNofTrueHits() const
virtual int32_t GetNofHits() const
Definition CbmTrack.h:58
TClonesArray * fTrdHitMatch
TrdHits.
CbmMCDataArray * fMCTracks
MC tracks.
static const int fgkNpdg
static int Pdg2Idx(int pdg)
TClonesArray * fTrdHits
TrdPoints.
McTrackInfo & getMcTrackInfo(const CbmLink &link)
std::map< const char *, std::array< TH2F *, fgkNpdg > > fhPidPtY
FairRootManager * fManager
map track link -> track info
Int_t fTrdNstations
MCtrack.
TClonesArray * fStsTrackMatches
StsTrack.
TClonesArray * fTrdTracks
GlobalTrack.
std::map< CbmLink, McTrackInfo > fMcTrackInfoMap
static const char * Idx2Symb(int idx)
CbmMCDataArray * fTrdPoints
static const char * fgkIdxName[fgkNpdg]
virtual void Finish()
TClonesArray * fTrdTrackMatches
TrdTrack.
std::vector< CbmQaHist< TProfile2D > > fhStationEffXY
output folder with histos and canvases
virtual void Exec(Option_t *opt)
static int Idx2Pdg(int idx)
static const char * Idx2Name(int idx)
virtual InitStatus Init()
virtual void SetParContainers()
CbmMCDataManager * fMcManager
TClonesArray * fStsTracks
TrdTrackMatch.
TClonesArray * fGlobalTracks
TrdHitMatch.
void DivideHistos(TH1 *histo1, TH1 *histo2, TH1 *histo3, Option_t *opt="")
TVector3 fTargetPos
StsTrackMatch.
CbmTrackingTrdQa(Int_t iVerbose=1)
CbmTimeSlice * fTimeSlice
static const char * fgkIdxSymb[fgkNpdg]
virtual InitStatus ReInit()
void FillTrackMatchMap(Int_t &nRec, Int_t &nGhosts, Int_t &nClones)
data class for a reconstructed Energy-4D measurement in the TRD
Definition CbmTrdHit.h:40
bool GetClassType() const
Definition CbmTrdHit.h:80
const algo::RecoSetup & GetSetup() const
Setup accessor.
static RecoSetupManager * Instance()
Instance access.
const trd::RecoSetupUnit * GetTrd() const
Access to TRD reconstruction setup unit.
Definition RecoSetup.h:148
Int_t fStsTrackMatch
matched TrdTrack index
Int_t fMatchedNHitsTrue
number of matched hits
Double_t fQuali
matched StsTrack index
std::map< Int_t, Int_t > fHitMap
Int_t fMatchedNHitsAll
percentage of matched hits
Int_t fTrdTrackMatch
matched GlobalTrack index
Bool_t fIsAccepted
number of matched hits
Int_t fGlobalTrackMatch
Trd station -> number of attached hits.