CbmRoot
Loading...
Searching...
No Matches
CbmTrdCalibTracker.cxx
Go to the documentation of this file.
1/* Copyright (C) 2022 Horia Hulubei National Institute of Physics and Nuclear Engineering, Bucharest
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Alexandru Bercuci[committer]*/
4
8
10
11#include "CbmDefs.h"
12#include "CbmDigiManager.h"
13#include "CbmLink.h"
14#include "CbmMCDataArray.h"
15#include "CbmMCDataManager.h"
16#include "CbmMCEventList.h"
17#include "CbmMCTrack.h"
18#include "CbmMatch.h"
19#include "CbmQaCanvas.h"
20//#include "CbmSetup.h"
21#include "CbmQaUtil.h"
22#include "CbmTimeSlice.h"
23#include "CbmTrdCluster.h"
24#include "CbmTrdDigi.h"
25#include "CbmTrdHit.h"
26#include "CbmTrdHitMC.h"
27#include "CbmTrdParModDigi.h" // for CbmTrdModule
28#include "CbmTrdParModGeo.h"
29#include "CbmTrdParSetDigi.h" // for CbmTrdParSetDigi
30#include "CbmTrdParSetGeo.h"
31#include "CbmTrdPoint.h"
32
33#include <FairRootManager.h>
34#include <FairRunAna.h>
35#include <FairRuntimeDb.h>
36#include <FairSink.h>
37#include <FairTask.h>
38#include <Logger.h>
39
40#include <TClonesArray.h>
41#include <TDatabasePDG.h>
42#include <TGenericClassInfo.h>
43#include <TGeoManager.h>
44#include <TGeoNode.h>
45#include <TMathBase.h>
46#include <TObjArray.h>
47#include <TParameter.h>
48#include <TParticlePDG.h>
49#include <TString.h>
50#include <TStyle.h>
51
52#include <algorithm>
53#include <cstdio>
54#include <iostream>
55#include <map>
56#include <vector>
57
59
60// -------------------------------------------------------------------------
61
62CbmTrdCalibTracker::CbmTrdCalibTracker(const char* name, Int_t verbose) : FairTask(name, verbose)
63{
64 // Constructor
65
66 // Create a list of histogramms
67
68 fHistList.clear();
69 fHistList.reserve(20);
70 fHistList.push_back(&fh1DresidualU);
71 fHistList.push_back(&fh1DresidualV);
72 fHistList.push_back(&fh1DresidualT);
73 fHistList.push_back(&fh2DresidualX);
74 fHistList.push_back(&fh2DresidualY);
75 fHistList.push_back(&fh2DresidualT);
76 fHistList.push_back(&fh1DpullU);
77 fHistList.push_back(&fh1DpullV);
78 fHistList.push_back(&fh1DpullT);
79 fHistList.push_back(&fh2DpullX);
80 fHistList.push_back(&fh2DpullY);
81 fHistList.push_back(&fh2DpullT);
82
83 // Keep the ownership on the histograms to avoid double deletion
84 for (unsigned int i = 0; i < fHistList.size(); i++) {
85 fHistList[i]->SetDirectory(0);
86 }
87}
88
89// -------------------------------------------------------------------------
94
95
96// -------------------------------------------------------------------------
98{
99 fIsTrdInSetup = 0;
100 fIsMcPresent = false;
102
103 fTimeSlice = nullptr;
104 fDigiManager = nullptr;
105
106 fMcManager = nullptr;
107 fMcTracks = nullptr;
108 fMcPoints = nullptr;
109
110 fClusters = nullptr;
111 fHits = nullptr;
112 fHitMatches = nullptr;
113
114 fOutFolder.Clear();
115
116 fHistFolder = nullptr;
117
118 fNevents.SetVal(0);
119
120 fhPointsPerHit.clear();
121 fhHitsPerPoint.clear();
122
123 if (fHitsMc) {
124 fHitsMc->Clear("C");
125 fHitsMc->Delete();
126 delete fHitsMc;
127 }
128}
129// -------------------------------------------------------------------------
130
132{
133 fTrdDigiPar = nullptr;
134 fTrdGeoPar = nullptr;
135
136 // Get run and runtime database
137 FairRunAna* ana = FairRunAna::Instance();
138 if (!ana) {
139 LOG(fatal) << "No FairRunAna instance";
140 }
141
142 FairRuntimeDb* rtdb = ana->GetRuntimeDb();
143 if (!rtdb) {
144 LOG(fatal) << "No FairRuntimeDb instance";
145 }
146
147 fTrdDigiPar = dynamic_cast<CbmTrdParSetDigi*>(rtdb->getContainer("CbmTrdParSetDigi"));
148 fTrdGeoPar = dynamic_cast<CbmTrdParSetGeo*>(rtdb->getContainer("CbmTrdParSetGeo"));
149}
150
151// -------------------------------------------------------------------------
153{
154 // Initialize and check the setup
155
156 DeInit();
157
158 // Check if TRD is present in the CBM setup
159 // A straightforward way to do so is CbmSetup::Instance()->IsActive(ECbmModuleId::kTrd),
160 // but unfortunatelly CbmSetup class is unaccesible.
161 // ( CbmSimSteer library requires libfreetype that has linking problems on MacOsX )
162 // Therefore let's simply check if any TRD material is present in the geometry.
163 //
164 fIsTrdInSetup = 0;
165 {
166 TObjArray* topNodes = gGeoManager->GetTopNode()->GetNodes();
167 for (Int_t iTopNode = 0; iTopNode < topNodes->GetEntriesFast(); iTopNode++) {
168 TGeoNode* topNode = static_cast<TGeoNode*>(topNodes->At(iTopNode));
169 if (TString(topNode->GetName()).Contains("trd")) {
170 fIsTrdInSetup = 1;
171 break;
172 }
173 }
174 }
175
176 if (!fIsTrdInSetup) {
177 LOG(info) << "TRD is not in the setup, do nothing";
178 return kSUCCESS;
179 }
180
181 // check parameter containers
182
183 if (!fTrdDigiPar) {
184 LOG(error) << "No CbmTrdParSetDigi container in FairRuntimeDb";
185 return kERROR;
186 }
187
188 if (!fTrdGeoPar) {
189 LOG(error) << "No CbmTrdParSetGeo container in FairRuntimeDb";
190 return kERROR;
191 }
192
193 FairRootManager* manager = FairRootManager::Instance();
195 fDigiManager->Init();
196
197 fTimeSlice = dynamic_cast<CbmTimeSlice*>(manager->GetObject("TimeSlice."));
198 if (!fTimeSlice) {
199 LOG(error) << "No time slice found";
200 return kERROR;
201 }
202
203 fHits = dynamic_cast<TClonesArray*>(manager->GetObject("TrdHit"));
204
205 if (!fHits) {
206 LOG(error) << "No TRD hit array found";
207 return kERROR;
208 }
209
210 fClusters = dynamic_cast<TClonesArray*>(manager->GetObject("TrdCluster"));
211
212 if (!fClusters) {
213 LOG(error) << "No TRD cluster array found";
214 return kERROR;
215 }
216
217 fMcManager = dynamic_cast<CbmMCDataManager*>(manager->GetObject("MCDataManager"));
218
219 fIsMcPresent = (fMcManager != nullptr);
220
221 if (fIsMcPresent) {
222
223 fMcEventList = dynamic_cast<CbmMCEventList*>(manager->GetObject("MCEventList."));
224 fMcTracks = fMcManager->InitBranch("MCTrack");
225 fMcPoints = fMcManager->InitBranch("TrdPoint");
226 fHitMatches = dynamic_cast<TClonesArray*>(manager->GetObject("TrdHitMatch"));
227
228 // we assume that when TRD is in the setup and the MC manager is present,
229 // the TRD MC data must be present too
230
231 if (!fMcEventList) {
232 LOG(error) << ": No MCEventList data!";
233 return kERROR;
234 }
235
236 if (!fMcTracks) {
237 LOG(error) << "No MC tracks found";
238 return kERROR;
239 }
240
241 if (!fMcPoints) {
242 LOG(error) << "No MC points found";
243 return kERROR;
244 }
245
246 if (!fHitMatches) {
247 LOG(error) << "No TRD hit matches found";
248 return kERROR;
249 }
250
251 if (!fDigiManager->IsMatchPresent(ECbmModuleId::kTrd)) {
252 LOG(error) << "No TRD digi matches found";
253 return kERROR;
254 }
255 }
256
257 // count TRD stations
258 // TODO: put the code below to some TRD geometry interface
259
261 {
262 Int_t layerCounter = 0;
263 TObjArray* topNodes = gGeoManager->GetTopNode()->GetNodes();
264 for (Int_t iTopNode = 0; iTopNode < topNodes->GetEntriesFast(); iTopNode++) {
265 TGeoNode* topNode = static_cast<TGeoNode*>(topNodes->At(iTopNode));
266 if (TString(topNode->GetName()).Contains("trd")) {
267 TObjArray* layers = topNode->GetNodes();
268 for (Int_t iLayer = 0; iLayer < layers->GetEntriesFast(); iLayer++) {
269 TGeoNode* layer = static_cast<TGeoNode*>(layers->At(iLayer));
270 if (TString(layer->GetName()).Contains("layer")) {
271 layerCounter++;
272 }
273 }
274 }
275 }
276 fNtrackingStations = layerCounter;
277 }
278
279 LOG(debug) << "Init(): fNtrackingStations = " << fNtrackingStations;
280
281 if (fNtrackingStations <= 0) {
282 LOG(error) << "can not count TRD tracking stations";
283 return kERROR;
284 }
285
286 // init output tree
287 fHitsMc = new TClonesArray("CbmTrdHitMC", 100);
288 manager->Register("TrdHitMC", "TRD", fHitsMc, IsOutputBranchPersistent("TrdHitMC"));
289
290 // initialise histograms
291 fOutFolder.SetOwner(false);
292 fHistFolder = fOutFolder.AddFolder("rawHist", "Raw Histograms");
293 gStyle->SetOptStat(0);
294
295 fNevents.SetVal(0);
296 fHistFolder->Add(&fNevents);
297
298 for (unsigned int i = 0; i < fHistList.size(); i++) {
299 fHistList[i]->Reset();
300 fHistFolder->Add(fHistList[i]);
301 }
302
303 fhPointsPerHit.clear();
304 fhHitsPerPoint.clear();
305 fhEfficiencyR.clear();
306 fhEfficiencyXY.clear();
307
312
313 for (Int_t i = 0; i < fNtrackingStations; i++) {
314 fhPointsPerHit.emplace_back(Form("hMcPointsPerHit%i", i),
315 Form("MC Points per Hit: Station %i;N mc points / hit", i), 10, -0.5, 9.5);
316 fhPointsPerHit[i].SetDirectory(0);
317 fhPointsPerHit[i].SetLineWidth(2);
318 fhPointsPerHit[i].SetOptStat(110);
319 fHistFolder->Add(&fhPointsPerHit[i]);
320 }
321
322 for (Int_t i = 0; i < fNtrackingStations; i++) {
323 fhHitsPerPoint.emplace_back(Form("hHitsPerMcPoint%i", i),
324 Form("Hits per MC Point: Station %i; N hits / mc point", i), 10, -0.5, 9.5);
325 fhHitsPerPoint[i].SetDirectory(0);
326 fhHitsPerPoint[i].SetLineWidth(2);
327 fhHitsPerPoint[i].SetOptStat(110);
328 fHistFolder->Add(&fhHitsPerPoint[i]);
329 }
330
331 for (Int_t i = 0; i < fNtrackingStations; i++) {
332
333 double dx = 350;
334 double dy = 350;
335 double dr = sqrt(dx * dx + dy * dy);
336
337 fhEfficiencyR.emplace_back(Form("hEfficiencyR%i", i), Form("Efficiency R: Station %i;R [cm]", i), 100, 0, dr);
338 fhEfficiencyR[i].SetDirectory(0);
339 fhEfficiencyR[i].SetOptStat(1110);
340 fHistFolder->Add(&fhEfficiencyR[i]);
341
342 fhEfficiencyXY.emplace_back(Form("hEfficiencyXY%i", i), Form("Efficiency XY: Station %i;X [cm];Y [cm]", i), 50, -dx,
343 dx, 50, -dy, dy);
344 fhEfficiencyXY[i].SetDirectory(0);
345 fhEfficiencyXY[i].SetOptStat(10);
346 fhEfficiencyXY[i].GetYaxis()->SetTitleOffset(1.4);
347 fHistFolder->Add(&fhEfficiencyXY[i]);
348 }
349
350 fCanvResidual.Clear();
351 fCanvResidual.Divide2D(6);
352
353 fCanvPull.Clear();
354 fCanvPull.Divide2D(6);
355
356 fCanvEfficiencyR.Clear();
358
359 fCanvEfficiencyXY.Clear();
361
362 fCanvPointsPerHit.Clear();
364
365 fCanvHitsPerPoint.Clear();
367
369 fOutFolder.Add(&fCanvPull);
374
375 // analyse the geometry setup
376 return GeometryQa();
377}
378
379
380// -------------------------------------------------------------------------
382{
383 if (!fIsTrdInSetup) {
384 return;
385 }
386
387 // update number of processed events
388 fNevents.SetVal(fNevents.GetVal() + 1);
389 LOG(debug) << "Event: " << fNevents.GetVal();
390 ResolutionQa();
391}
392
393
394// -------------------------------------------------------------------------
396{
397 FairSink* sink = FairRootManager::Instance()->GetSink();
398 if (sink) {
399 sink->WriteObject(&GetQa(), nullptr);
400 }
401}
402
403
404// -------------------------------------------------------------------------
406{
407 // check geometry of the TRD modules
408
409 double lastZ = -1;
410
411 // load alignment matrices
412 fTrdGeoPar->LoadAlignVolumes();
413 for (int iStation = 0; iStation < fNtrackingStations; iStation++) {
414
415 int ModuleId = fTrdDigiPar->GetModuleId(iStation);
416
417 const CbmTrdParModGeo* pGeo = dynamic_cast<const CbmTrdParModGeo*>(fTrdGeoPar->GetModulePar(ModuleId));
418 if (!pGeo) {
419 LOG(fatal) << " No Geo params for station " << iStation << " module " << ModuleId << " found ";
420 return kERROR;
421 }
422
423 double staZ = pGeo->GetZ();
424
425 // check that the stations are properly ordered in Z
426 if ((iStation > 0) && (staZ <= lastZ)) {
427 LOG(error) << "TRD trackig stations are not properly ordered in Z";
428 return kERROR;
429 }
430 lastZ = staZ;
431
432 //TODO: what are these 3 and 6 here in L1 code? Why are they hardcoed?
433 // if (num == 0 || num == 2 || num == 4) geo.push_back(3);
434 // if (num == 1 || num == 3) geo.push_back(6);
435 }
436
437 return kSUCCESS;
438}
439
440// -------------------------------------------------------------------------
442{
443 if (!fDigiManager->IsMatchPresent(ECbmModuleId::kTrd)) return;
444
445 if (!fIsMcPresent) return;
446
447 Int_t nHits = fHits->GetEntriesFast();
448 Int_t nClusters = fClusters->GetEntriesFast();
449 Int_t nDigis = fDigiManager->GetNofDigis(ECbmModuleId::kTrd);
450
451 int nMcEvents = fMcEventList->GetNofEvents();
452 int imc(0); // index of hit->MC QA objects
453 fHitsMc->Delete();
454
455 // Vector of integers parallel to mc points
456 std::vector<std::vector<int>> pointNhits;
457 pointNhits.resize(nMcEvents);
458 for (int iE = 0; iE < nMcEvents; iE++) {
459 int fileId = fMcEventList->GetFileIdByIndex(iE);
460 int eventId = fMcEventList->GetEventIdByIndex(iE);
461 int nPoints = fMcPoints->Size(fileId, eventId);
462 pointNhits[iE].resize(nPoints, 0);
463 }
464
465 for (Int_t iHit = 0; iHit < nHits; iHit++) {
466
467 const CbmTrdHit* hit = dynamic_cast<const CbmTrdHit*>(fHits->At(iHit));
468 if (!hit) {
469 LOG(error) << "TRD hit N " << iHit << " doesn't exist";
470 return;
471 }
472
473 if ((int) hit->GetClassType() > 1) {
474 // class type 0: trd-1D
475 // class type 1: trd-2D
476 // the return type is (currently) boolean, so this error should never happen
477 LOG(error) << "Unknown detector class type " << hit->GetClassType() << " for TRD hit N " << iHit;
478 return;
479 }
480
481 const int address = hit->GetAddress();
482 int StationIndex = CbmTrdAddress::GetLayerId(address);
483
484 if (StationIndex < 0 || StationIndex >= fNtrackingStations) {
485 LOG(fatal) << "TRD hit layer id " << StationIndex << " is out of range";
486 return;
487 }
488
489 Int_t clusterId = hit->GetRefId();
490 if (clusterId < 0 || clusterId >= nClusters) {
491 LOG(error) << "TRD cluster id " << clusterId << " is out of range";
492 return;
493 }
494
495 CbmTrdCluster* cluster = dynamic_cast<CbmTrdCluster*>(fClusters->At(clusterId));
496 if (!cluster) {
497 LOG(error) << "TRD cluster N " << clusterId << " doesn't exist";
498 return;
499 }
500
501 if (cluster->GetAddress() != address) {
502 LOG(error) << "TRD hit address " << address << " differs from its cluster address " << cluster->GetAddress();
503 return;
504 }
505
506 Int_t nClusterDigis = cluster->GetNofDigis();
507
508 if (nClusterDigis <= 0) {
509 LOG(error) << "hit match for TRD cluster " << clusterId << " has no digis ";
510 return;
511 }
512
513 // construct QA object
514 CbmTrdHitMC* hmc = new ((*fHitsMc)[imc++]) CbmTrdHitMC(*hit);
515 hmc->AddCluster(cluster);
516 uint64_t tdigi = 0;
517
518 // custom finder of the digi matches
519 CbmMatch myHitMatch;
520 for (Int_t iDigi = 0; iDigi < nClusterDigis; iDigi++) {
521 Int_t digiIdx = cluster->GetDigi(iDigi);
522 if (digiIdx < 0 || digiIdx >= nDigis) {
523 LOG(fatal) << "TRD cluster: digi index " << digiIdx << " out of range ";
524 return;
525 }
526 const CbmTrdDigi* digi = CbmDigiManager::Instance()->Get<CbmTrdDigi>(digiIdx);
527 if (!digi) {
528 LOG(fatal) << "TRD digi " << iDigi << " not found";
529 return;
530 }
531
532 if (digi->GetAddressModule() != address) {
533 std::stringstream ss;
534 ss << "TRD hit address " << address << " differs from its digi address " << digi->GetAddressModule();
535 hmc->SetErrorMsg(ss.str());
536 LOG(error) << ss.str();
537 return;
538 }
539 switch (digi->GetType()) {
541 if (!tdigi) tdigi = digi->GetTimeDAQ();
542 break;
544 if (!tdigi) tdigi = digi->GetTime();
545 break;
546 default: LOG(fatal) << "TRD digi type neither SPADIC or FASP"; return;
547 }
548 hmc->AddSignal(digi, tdigi);
549
550 const CbmMatch* match = dynamic_cast<const CbmMatch*>(fDigiManager->GetMatch(ECbmModuleId::kTrd, digiIdx));
551 if (!match) {
552 LOG(fatal) << "TRD digi match " << digiIdx << " not found";
553 return;
554 }
555 myHitMatch.AddLinks(*match);
556 }
557 hmc->PurgeSignals();
558
559 if (myHitMatch.GetNofLinks() == 0) {
560 continue;
561 }
562
563 const CbmLink& bestLink = myHitMatch.GetMatchedLink();
564
565 { // check if the hit match is correct
566 CbmMatch* hitMatch = dynamic_cast<CbmMatch*>(fHitMatches->At(iHit));
567 if (hitMatch == nullptr) {
568 std::stringstream ss;
569 ss << "hit match for TRD hit " << iHit << " doesn't exist";
570 hmc->SetErrorMsg(ss.str());
571 LOG(error) << ss.str();
572 }
573 else {
574 const CbmLink& link = hitMatch->GetMatchedLink();
575 if ((link != bestLink) && (link.GetWeight() != bestLink.GetWeight())) {
576 std::stringstream ss;
577 ss << "hit match for TRD hit " << iHit << " doesn't correspond to digi matches";
578 hmc->SetErrorMsg(ss.str());
579 LOG(error) << ss.str();
580 }
581 }
582 }
583
584 // how many MC points? fill N hits per mc point
585
586 int nHitPoints = 0;
587 for (int i = 0; i < myHitMatch.GetNofLinks(); i++) {
588 const CbmLink& link = myHitMatch.GetLink(i);
589 if (link.GetIndex() >= 0) { // not a noise digi
590 nHitPoints++;
591 int iE = fMcEventList->GetEventIndex(link);
592 if (iE < 0 || iE >= nMcEvents) {
593 std::stringstream ss;
594 ss << "link points to a non-existing MC event";
595 hmc->SetErrorMsg(ss.str());
596 LOG(error) << ss.str();
597 return;
598 }
599 if (link.GetIndex() >= (int) pointNhits[iE].size()) {
600 std::stringstream ss;
601 ss << "link points to a non-existing MC index";
602 hmc->SetErrorMsg(ss.str());
603 LOG(error) << ss.str();
604 return;
605 }
606 pointNhits[iE][link.GetIndex()]++;
607 }
608 }
609
610 fhPointsPerHit[StationIndex].Fill(nHitPoints);
611
612 // take corresponding MC point
613
614 // skip hits from the noise digis
615 if (bestLink.GetIndex() < 0) {
616 std::stringstream ss;
617 ss << "hit from noise [INFO]";
618 hmc->SetErrorMsg(ss.str());
619 continue;
620 }
621
622 CbmTrdPoint* p = dynamic_cast<CbmTrdPoint*>(fMcPoints->Get(bestLink));
623 if (p == nullptr) {
624 std::stringstream ss;
625 ss << "link points to a non-existing MC point";
626 hmc->SetErrorMsg(ss.str());
627 LOG(error) << ss.str();
628 return;
629 }
630
631 if (p->GetModuleAddress() != (int) CbmTrdAddress::GetModuleAddress(address)) {
632 std::stringstream ss;
633 ss << "mc point module address differs from the hit module address";
634 hmc->SetErrorMsg(ss.str());
635 LOG(error) << ss.str();
636 return;
637 }
638
639 // check that the nominal station Z is not far from the active material
640 {
641 // the cut of 1 cm is choosen arbitrary and can be changed
642
643 int ModuleId = fTrdDigiPar->GetModuleId(StationIndex);
644
645 const CbmTrdParModGeo* pGeo = dynamic_cast<const CbmTrdParModGeo*>(fTrdGeoPar->GetModulePar(ModuleId));
646 if (!pGeo) {
647 LOG(fatal) << " No Geo params for station " << StationIndex << " module " << ModuleId << " found ";
648 return;
649 }
650
651 double staZ = pGeo->GetZ(); // module->GetZ(); //+ 410;
652 if ((staZ < p->GetZIn() - 1.) || (staZ > p->GetZOut() + 1.)) {
653 std::stringstream ss;
654 ss << "TRD station " << StationIndex << ": active material Z[" << p->GetZIn() << " cm," << p->GetZOut()
655 << " cm] is too far from the nominal station Z " << staZ << " cm";
656 hmc->SetErrorMsg(ss.str());
657 LOG(error) << ss.str();
658 return;
659 }
660 // the cut of 1 cm is choosen arbitrary and can be changed
661 if (fabs(hit->GetZ() - staZ) > 1.) {
662 std::stringstream ss;
663 ss << "TRD station " << StationIndex << ": hit Z " << hit->GetZ()
664 << " cm, is too far from the nominal station Z " << staZ << " cm";
665 hmc->SetErrorMsg(ss.str());
666 LOG(error) << ss.str();
667 return;
668 }
669 }
670
671 // residual and pull
672
673 if (nHitPoints != 1) {
674 std::stringstream ss;
675 ss << "hit from mixed hit [INFO] nPoints=" << nHitPoints;
676 hmc->SetErrorMsg(ss.str());
677 continue; // only check residual for non-mixed hits
678 }
679
680 double t0 = fMcEventList->GetEventTime(bestLink);
681 if (t0 < 0) {
682 std::stringstream ss;
683 ss << "MC event time doesn't exist for a TRD link";
684 hmc->SetErrorMsg(ss.str());
685 LOG(error) << ss.str();
686 return;
687 }
688
689 double x = p->GetX(); // take the "In" point since the time is stored only for this point
690 double y = p->GetY();
691 double z = p->GetZ();
692 double t = t0 + p->GetTime();
693 double px = p->GetPx();
694 double py = p->GetPy();
695 double pz = p->GetPz();
696
697 // skip very slow particles
698 if (fabs(p->GetPzOut()) < 0.01) continue;
699
700 // transport the particle from the MC point to the hit Z
701 double dz = hit->GetZ() - z;
702 x += dz * px / pz;
703 y += dz * py / pz;
704
705 // get particle mass
706 Int_t pdg(0);
707 double mass = 0;
708 {
709 CbmLink mcTrackLink = bestLink;
710 mcTrackLink.SetIndex(p->GetTrackID());
711 CbmMCTrack* mcTrack = dynamic_cast<CbmMCTrack*>(fMcTracks->Get(mcTrackLink));
712 if (!mcTrack) {
713 std::stringstream ss;
714 ss << "MC track " << p->GetTrackID() << " doesn't exist";
715 hmc->SetErrorMsg(ss.str());
716 LOG(error) << ss.str();
717 return;
718 }
719 pdg = mcTrack->GetPdgCode();
720 if (pdg < 9999999 && ((TParticlePDG*) TDatabasePDG::Instance()->GetParticle(pdg))) {
721 mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
722 }
723 }
724 hmc->AddPoint(p, t0, pdg);
725 //std::cout << hmc->ToString() << "\n";
726 constexpr double speedOfLight = 29.979246; // cm/ns
727 TVector3 mom3;
728 p->Momentum(mom3);
729 t += dz / (pz * speedOfLight) * sqrt(mass * mass + mom3.Mag2());
730
731 double du = hit->GetX() - x;
732 double dv = hit->GetY() - y; //hmc->GetDy(); //hit->GetY() - y;
733 double dt = hit->GetTime() - t; // hmc->GetDt(); //hit->GetTime() - t;
734 double su = hmc->GetDx();
735 double sv = hit->GetDy();
736 double st = hit->GetTimeError();
737
738 // residuals
739 if (hit->GetClassType() == 0) {
740 if (StationIndex % 2) {
741 std::swap(du, dv);
742 std::swap(su, sv);
743 }
744 fh1DresidualU.Fill(du);
745 fh1DresidualV.Fill(dv);
746 fh1DresidualT.Fill(hit->GetTime() - t);
747 }
748 else {
749 fh2DresidualX.Fill(du);
750 fh2DresidualY.Fill(dv);
751 fh2DresidualT.Fill(hit->GetTime() - t);
752 }
753
754 // pulls
755 if ((su < 1.e-5) || (sv < 1.e-5) || (st < 1.e-5)) {
756 std::stringstream ss;
757 ss << "TRD hit errors are not properly set: errX " << hit->GetDx() << " errY " << hit->GetDy() << " errT ";
758 hmc->SetErrorMsg(ss.str());
759 LOG(error) << ss.str();
760 return;
761 }
762
763 if (hit->GetClassType() == 0) {
764 fh1DpullU.Fill(du / su);
765 fh1DpullV.Fill(dv / sv);
766 fh1DpullT.Fill(dt / st);
767 }
768 else {
769 fh2DpullX.Fill(du / su);
770 fh2DpullY.Fill(dv / sv);
771 fh2DpullT.Fill(dt / st);
772 }
773
774 } // iHit
775
776 // fill efficiency: N hits per MC point
777
778 for (int iE = 0; iE < nMcEvents; iE++) {
779 int fileId = fMcEventList->GetFileIdByIndex(iE);
780 int eventId = fMcEventList->GetEventIdByIndex(iE);
781 int nPoints = fMcPoints->Size(fileId, eventId);
782 for (int ip = 0; ip < nPoints; ip++) {
783 CbmTrdPoint* p = dynamic_cast<CbmTrdPoint*>(fMcPoints->Get(fileId, eventId, ip));
784 if (p == nullptr) {
785 LOG(error) << "MC point doesn't exist";
786 return;
787 }
788 auto address = p->GetModuleAddress();
789 int StationIndex = CbmTrdAddress::GetLayerId(address);
790 if (StationIndex < 0 || StationIndex >= fNtrackingStations) {
791 LOG(fatal) << "TRD hit layer id " << StationIndex << " for module " << address << " is out of range";
792 return;
793 }
794 fhHitsPerPoint[StationIndex].Fill(pointNhits[iE][ip]);
795 fhEfficiencyR[StationIndex].Fill(sqrt(p->GetX() * p->GetX() + p->GetY() * p->GetY()), (pointNhits[iE][ip] > 0));
796 fhEfficiencyXY[StationIndex].Fill(p->GetX(), p->GetY(), (pointNhits[iE][ip] > 0));
797 }
798 }
799}
800
801
802// -------------------------------------------------------------------------
804{
805 gStyle->SetPaperSize(20, 20);
806
807 for (Int_t i = 0; i < fNtrackingStations; i++) {
808 fCanvHitsPerPoint.cd(i + 1);
809 fhHitsPerPoint[i].DrawCopy("", "");
810 fCanvPointsPerHit.cd(i + 1);
811 fhPointsPerHit[i].DrawCopy("", "");
812
813 fCanvEfficiencyR.cd(i + 1);
814 fhEfficiencyR[i].DrawCopy("colz", "");
815
816 fCanvEfficiencyXY.cd(i + 1);
817 fhEfficiencyXY[i].DrawCopy("colz", "");
818 }
819
820 for (UInt_t i = 0; i < fHistList.size(); i++) {
822 }
823
824 for (UInt_t i = 0; i < 6; i++) {
825 fCanvResidual.cd(i + 1);
826 fHistList[i]->DrawCopy("", "");
827 fCanvPull.cd(i + 1);
828 fHistList[6 + i]->DrawCopy("", "");
829 }
830
831 return fOutFolder;
832}
@ kTrd
Transition Radiation Detector.
Definition CbmDefs.h:51
Definition of the CbmQaCanvas class.
Useful utilities for CBM QA tasks.
ClassImp(CbmTrdCalibTracker)
Data Container for TRD clusters.
Class for hits in TRD detector.
friend fvec sqrt(const fvec &a)
int Int_t
int32_t GetAddress() const
Definition CbmCluster.h:90
int32_t GetDigi(int32_t index) const
Get digi at position index.
Definition CbmCluster.h:76
int32_t GetNofDigis() const
Number of digis in cluster.
Definition CbmCluster.h:69
const Digi * Get(Int_t index) const
Get a digi object.
static CbmDigiManager * Instance()
Static instance.
double GetTimeError() const
Definition CbmHit.h:80
double GetTime() const
Definition CbmHit.h:79
int32_t GetAddress() const
Definition CbmHit.h:77
double GetZ() const
Definition CbmHit.h:74
int32_t GetRefId() const
Definition CbmHit.h:76
Task class creating and managing CbmMCDataArray objects.
Container class for MC events with number, file and start time.
int32_t GetPdgCode() const
Definition CbmMCTrack.h:67
void AddLinks(const CbmMatch &match)
Definition CbmMatch.cxx:39
const CbmLink & GetLink(int32_t i) const
Definition CbmMatch.h:39
int32_t GetNofLinks() const
Definition CbmMatch.h:42
const CbmLink & GetMatchedLink() const
Definition CbmMatch.h:41
double GetDy() const
Definition CbmPixelHit.h:76
double GetDx() const
Definition CbmPixelHit.h:75
double GetY() const
Definition CbmPixelHit.h:74
double GetX() const
Definition CbmPixelHit.h:73
Bookkeeping of time-slice content.
static uint32_t GetLayerId(uint32_t address)
Return layer ID from address.
static uint32_t GetModuleAddress(uint32_t address)
Return unique module ID from address.
std::vector< CbmQaHist< TH1D > * > fHistList
List of the above histogramms.
CbmQaHist< TH1D > fh1DpullV
CbmQaHist< TH1D > fh2DresidualY
TClonesArray * fClusters
Data.
CbmQaHist< TH1D > fh2DresidualX
CbmMCDataArray * fMcTracks
std::vector< CbmQaHist< TH1I > > fhHitsPerPoint
hits efficiency
CbmQaHist< TH1D > fh1DpullU
Histogram for PULL Distribution.
CbmTimeSlice * fTimeSlice
CbmQaHist< TH1D > fh1DresidualV
void Finish()
FairTask: Action at end of run. For this task and all of the subtasks.
CbmQaHist< TH1D > fh1DresidualU
number of processed events
TParameter< int > fNevents
subfolder for histograms
std::vector< CbmQaHist< TProfile > > fhEfficiencyR
CbmTrdCalibTracker(const char *name="TrdQaTracker", Int_t verbose=1)
Constructor.
void DeInit()
Free the memory etc.
CbmQaHist< TH1D > fh2DpullY
CbmQaCanvas fCanvResidual
Canvaces: collection of histogramms.
CbmMCDataManager * fMcManager
CbmTrdParSetDigi * fTrdDigiPar
CbmQaHist< TH1D > fh2DresidualT
CbmMCDataArray * fMcPoints
CbmDigiManager * fDigiManager
CbmQaHist< TH1D > fh2DpullT
CbmTrdParSetGeo * fTrdGeoPar
InitStatus ReInit()
FairTask: Reinitialisation.
void ResolutionQa()
Analysis of hit uncertainty (pull) distributions.
std::vector< CbmQaHist< TProfile2D > > fhEfficiencyXY
hits efficiency
void SetParContainers()
FairTask: Intialise parameter containers.
TClonesArray * fHitMatches
void Exec(Option_t *)
TTask: Process a timeslice.
InitStatus GeometryQa()
Check the geometry.
CbmQaHist< TH1D > fh2DpullX
std::vector< CbmQaHist< TH1I > > fhPointsPerHit
hits purity
TClonesArray * fHitsMc
Output.
CbmQaHist< TH1D > fh1DresidualT
TFolder * fHistFolder
output folder with histos and canvases
CbmQaHist< TH1D > fh1DpullT
~CbmTrdCalibTracker()
Destructor.
CbmMCEventList * fMcEventList
MC data.
Data Container for TRD clusters.
int32_t GetAddressModule() const
Getter module address in the experiment.
uint64_t GetTimeDAQ() const
Getter for global DAQ time [clk]. Differs for each ASIC. [In FASP case DAQ time is already stored in ...
Definition CbmTrdDigi.h:158
eCbmTrdAsicType GetType() const
Channel FEE SPADIC/FASP according to CbmTrdAsicType.
Definition CbmTrdDigi.h:170
double GetTime() const
Getter for physical time [ns]. Accounts for clock representation of each ASIC. In SPADIC case physica...
Definition CbmTrdDigi.h:155
TRD hit to MC point correlation class.
Definition CbmTrdHitMC.h:28
size_t AddPoint(const CbmTrdPoint *p, double t, int id)
Add MC points to the hit. The first time this function is called is for the best matched MC point.
void SetErrorMsg(std::string msg)
Store error message.
void AddCluster(const CbmTrdCluster *c)
Copy cluster details.
size_t PurgeSignals()
Applies to TRD2D and remove 0 charges from the boundaries of the cluster.
size_t AddSignal(const CbmTrdDigi *d, uint64_t t0)
Add signal values in the increasing order of pad index.
double GetDx() const
Calculate residuals in the bending plane.
data class for a reconstructed Energy-4D measurement in the TRD
Definition CbmTrdHit.h:40
bool GetClassType() const
Definition CbmTrdHit.h:80
Definition of geometry for one TRD module.
virtual Double_t GetZ() const
double GetPzOut() const
Definition CbmTrdPoint.h:73
double GetZIn() const
Definition CbmTrdPoint.h:67
int32_t GetModuleAddress() const
Definition CbmTrdPoint.h:78
double GetZOut() const
Definition CbmTrdPoint.h:70
std::tuple< double, double > FitKaniadakisGaussian(TH1 *pHist)
Fit a histogram with Kaniadakis Gaussian Distribution.
Definition CbmQaUtil.cxx:88