CbmRoot
Loading...
Searching...
No Matches
CbmCaOutputQa.cxx
Go to the documentation of this file.
1/* Copyright (C) 2023-2024 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Sergei Zharko [committer] */
4
9
10#include "CbmCaOutputQa.h"
11
12#include "CbmCaMCModule.h"
14#include "FairRootManager.h"
15#include "Logger.h"
16#include "TAttLine.h"
17#include "THStack.h"
18#include "TMarker.h"
19#include "TPad.h"
20#include "TText.h"
21
28
29// ---------------------------------------------------------------------------------------------------------------------
30//
31OutputQa::OutputQa(int verbose, bool isMCUsed, ECbmRecoMode recoMode)
32 : CbmQaTask("CbmCaOutputQa", verbose, isMCUsed, recoMode)
33{
34 // Create TS reader
35 fpTSReader = std::make_unique<TimeSliceReader>();
36
37 // Turn on default track classes
44
57
70
71 //AddTrackType(ETrackType::kAllE);
72 //AddTrackType(ETrackType::kAllMU);
73 //AddTrackType(ETrackType::kAllPI);
74 //AddTrackType(ETrackType::kAllK);
75 //AddTrackType(ETrackType::kAllPPBAR);
76
77 // Init track type histograms drawing attributes
79}
80
81// ---------------------------------------------------------------------------------------------------------------------
82//
84{
85 // Create summary table
86 if (IsMCUsed()) {
87 fpMCModule->Finish();
88
89 // TODO: Add cuts on entries from fmSummaryTableEntries
90 std::vector<ETrackType> vTypesToPlot;
91 int nRows = std::count_if(fmSummaryTableEntries.begin(), fmSummaryTableEntries.end(),
92 [&](const auto& f) { return fvbTrackTypeOn[f] && fvpTrackHistograms[f]->IsMCUsed(); })
93 + 1;
94
95 CbmQaTable* aTable = MakeQaObject<CbmQaTable>("summary_table", "Tracking summary table", nRows + 1, 9);
96 int iRow = 0;
97 std::vector<std::string> colNames = {"Efficiency", "Killed", "Length", "Fakes", "Clones",
98 "Reco/Evt", "MC/Evt", "Nst(hit)", "Nst(point)"};
99 aTable->SetNamesOfCols(colNames);
100 aTable->SetColWidth(14);
101 double nEvents = static_cast<double>(GetEventNumber());
102 LOG(info) << "Number of events: " << GetEventNumber();
103 for (auto trType : fmSummaryTableEntries) {
104 if (!fvbTrackTypeOn[trType] || !fvpTrackHistograms[trType]->IsMCUsed()) {
105 continue;
106 }
107 aTable->SetRowName(iRow, fvpTrackHistograms[trType]->GetTitle());
108 aTable->SetCell(iRow, 0, fvpTrackHistograms[trType]->GetIntegratedEff());
109 aTable->SetCell(iRow, 1, fvpTrackHistograms[trType]->GetKilledRate());
110 aTable->SetCell(iRow, 2, fvpTrackHistograms[trType]->GetAverageRecoLength());
111 aTable->SetCell(iRow, 3, fvpTrackHistograms[trType]->GetAverageFakeLength());
112 aTable->SetCell(iRow, 4, fvpTrackHistograms[trType]->GetClonesRate());
113 aTable->SetCell(iRow, 5, fvpTrackHistograms[trType]->GetNofRecoTracksMatched() / nEvents);
114 aTable->SetCell(iRow, 6, fvpTrackHistograms[trType]->GetNofMCTracks() / nEvents);
115 aTable->SetCell(iRow, 7, fvpTrackHistograms[trType]->GetAverageNofStationsWithHit());
116 aTable->SetCell(iRow, 8, fvpTrackHistograms[trType]->GetAverageNofStationsWithPoint());
117 ++iRow;
118 }
119 double nGhosts = 0.;
121 nGhosts = fvpTrackHistograms[ETrackType::kGhost]->fph_reco_p->GetEntries();
122 aTable->SetRowName(iRow, "N ghosts");
123 aTable->SetCell(iRow, 0, nGhosts);
124 aTable->SetRowName(iRow + 1, "Ghost rate");
125 aTable->SetCell(iRow + 1, 0, nGhosts / fMonitor.GetCounterValue(EMonitorKey::kTrack));
126 }
127 LOG(info) << '\n' << aTable->ToString(3);
128 }
129
130 LOG(info) << '\n' << fMonitor.ToString();
131}
132
133// ---------------------------------------------------------------------------------------------------------------------
134//
136{
137 constexpr double kZmin = 0.;
138 constexpr double kZmax = 300.;
139 constexpr double kXmin = -100.;
140 constexpr double kXmax = +100.;
141 constexpr double kYmin = -100.;
142 constexpr double kYmax = +100.;
143 constexpr std::array<Color_t, 11> kColorMC = {205, 209, 213, 217, 225, 208, 213, 216, 219, 224, 227};
144 constexpr std::array<Color_t, 3> kColorGhost = {201, 202, 203};
145 constexpr int kCanvX = 1920;
146 constexpr int kCanvY = 1080;
147 constexpr double kLMargin = 0.05; // Left margin
148 constexpr double kVEMargin = 0.15; // Vertical margin (top + bottom)
149 constexpr double kRMarginDispl = 0.4;
150 constexpr double kVIMargin = 0.0001;
151 constexpr Marker_t kMarkerPointWHit = 25;
152 constexpr Marker_t kMarkerPointWOHit = 5;
153 constexpr Marker_t kMarkerHitWPoint = 24;
154 constexpr Marker_t kMarkerHitWOPoint = 28;
155 constexpr double kFontSize = 0.035;
156 constexpr Width_t kLineWidth = 1;
157 constexpr Style_t kLineMCTrackReconstructed = 9;
158 constexpr Style_t kLineMCTrackNotReconstructed = 10;
159 constexpr Style_t kLineRecoTrackGhost = 2;
160 constexpr Style_t kLineRecoTrackNotGhost = 1;
161
162 int iEvent = GetEventNumber();
163 auto* pCanv = MakeQaObject<CbmQaCanvas>(Form("events/event_%d", iEvent), Form("event_%d", iEvent), kCanvX, kCanvY);
164 pCanv->Divide(1, 2, kVIMargin, kVIMargin);
165 pCanv->cd(1);
166 gPad->SetMargin(kLMargin, kRMarginDispl, kVIMargin, kVEMargin);
167 gPad->SetGrid();
168 auto* pHistX = new TH2F(Form("hFrameX_%d", iEvent), ";z [cm];x [cm]", 2, kZmin, kZmax, 2, kXmin, kXmax);
169 pHistX->GetYaxis()->SetTitleOffset(0.6);
170 pHistX->GetXaxis()->SetLabelSize(kFontSize);
171 pHistX->GetYaxis()->SetLabelSize(kFontSize);
172 pHistX->GetXaxis()->SetTitleSize(kFontSize);
173 pHistX->GetYaxis()->SetTitleSize(kFontSize);
174 pHistX->SetStats(false);
175 pHistX->Draw();
176 pCanv->cd(2);
177 gPad->SetMargin(kLMargin, kRMarginDispl, kVEMargin, kVIMargin);
178 gPad->SetGrid();
179 auto* pHistY = new TH2F(Form("hFrameY_%d", iEvent), ";z [cm];y [cm]", 2, kZmin, kZmax, 2, kYmin, kYmax);
180 pHistY->GetYaxis()->SetTitleOffset(0.6);
181 pHistY->GetXaxis()->SetLabelSize(kFontSize);
182 pHistY->GetYaxis()->SetLabelSize(kFontSize);
183 pHistY->GetXaxis()->SetTitleSize(kFontSize);
184 pHistY->GetYaxis()->SetTitleSize(kFontSize);
185 pHistY->SetStats(false);
186 pHistY->Draw();
187
188 pCanv->cd(1);
189
190 auto* pHeader = new TLegend(kLMargin, 1 - kVEMargin + 0.01, 0.99, 0.99);
191 pHeader->SetNColumns(6);
192 pHeader->SetTextSize(kFontSize);
193 pHeader->SetMargin(0.1);
194 pHeader->AddEntry(static_cast<TObject*>(nullptr), Form("event #%d", iEvent), "");
195 pHeader->AddEntry(new TMarker(0, 0, kMarkerPointWHit), "point w/ hit", "p");
196 pHeader->AddEntry(new TMarker(0, 0, kMarkerPointWOHit), "point w/o hit", "p");
197 pHeader->AddEntry(new TMarker(0, 0, kMarkerHitWPoint), "hit w/ point", "p");
198 pHeader->AddEntry(new TMarker(0, 0, kMarkerHitWOPoint), "hit w/o point", "p");
199 pHeader->Draw("same");
200
201 auto* pLegendReco = new TLegend(1 - kRMarginDispl, kVIMargin, 0.99, 1 - kVEMargin, "Reco tracks");
202 pLegendReco->SetMargin(0.1);
203 pLegendReco->SetTextSize(kFontSize);
204 //pLegendReco->SetHeader("Reco Tracks", "C");
205 pCanv->cd(2);
206 auto* pLegendMC = new TLegend(1 - kRMarginDispl, kVEMargin, 0.99, 1 - kVIMargin, "MC tracks");
207 pLegendMC->SetMargin(0.1);
208 pLegendMC->SetTextSize(kFontSize);
209 //pLegendMC->SetHeader("MC Tracks", "C");
210
211 int iColorRec = 0;
212 int iColorGhost = 0;
213
214 // Draw MC tracks
215 std::map<int, Color_t> mMCtrkColors; // Trk ID vs. color
216 if (IsMCUsed()) {
217 // Draw MC tracks
219 for (int iTmc = 0; iTmc < nMCTracks; ++iTmc) {
220 const auto& trk = fMCData.GetTrack(iTmc);
221 int nPoints = trk.GetNofPoints();
222 if (nPoints == 0) {
223 continue;
224 }
225 std::vector<double> trkPointX(nPoints);
226 std::vector<double> trkPointY(nPoints);
227 std::vector<double> trkPointZ(nPoints);
228 for (int iPLoc = 0; iPLoc < nPoints; ++iPLoc) {
229 const auto& point = fMCData.GetPoint(trk.GetPointIndexes()[iPLoc]);
230 trkPointX[iPLoc] = point.GetX();
231 trkPointY[iPLoc] = point.GetY();
232 trkPointZ[iPLoc] = point.GetZ();
233 }
234 Color_t currColor = 1;
235 Style_t currStyle = trk.IsReconstructed() ? kLineMCTrackReconstructed : kLineMCTrackNotReconstructed;
236 currColor = kColorMC[iColorRec];
237 iColorRec = (iColorRec + 1) % static_cast<int>(kColorMC.size());
238 mMCtrkColors[trk.GetId()] = currColor;
239 {
240 pCanv->cd(1);
241 auto* gr = new TGraph(nPoints, trkPointZ.data(), trkPointX.data());
242 gr->SetMarkerStyle(1);
243 gr->SetMarkerColor(currColor);
244 gr->SetLineColor(currColor);
245 gr->SetLineStyle(currStyle);
246 gr->SetLineWidth(kLineWidth);
247 gr->Draw("LSAME");
248
249 std::stringstream msg;
250 msg << "ID=" << trk.GetId() << ", ";
251 msg << "PDG=" << trk.GetPdgCode() << ", ";
252 msg << "p=" << trk.GetP() << " GeV/c, ";
253 msg << "rec-able=" << trk.IsReconstructable() << ", ";
254 msg << "rec-ed=" << trk.IsReconstructed() << ", ";
255 if (trk.GetNofRecoTracks() > 0) {
256 msg << "reco_IDs=(";
257 for (int iTr : trk.GetRecoTrackIndexes()) {
258 msg << iTr << ' ';
259 }
260 msg << "), ";
261 }
262 if (trk.GetNofTouchTracks() > 0) {
263 msg << "touch_IDs=(";
264 for (int iTr : trk.GetTouchTrackIndexes()) {
265 msg << iTr << ' ';
266 }
267 msg << "), ";
268 }
269 pLegendMC->AddEntry(gr, msg.str().c_str(), "l");
270 }
271 {
272 pCanv->cd(2);
273 auto* gr = new TGraph(nPoints, trkPointZ.data(), trkPointY.data());
274 gr->SetMarkerStyle(1);
275 gr->SetMarkerColor(currColor);
276 gr->SetLineColor(currColor);
277 gr->SetLineStyle(currStyle);
278 gr->SetLineWidth(kLineWidth);
279 gr->Draw("LSAME");
280 }
281 }
282
283 // Draw MC points
284 int nPoints = fMCData.GetNofPoints();
285 for (int iP = 0; iP < nPoints; ++iP) {
286 const auto& point = fMCData.GetPoint(iP);
287 bool bHasHit = point.GetHitIndexes().size() > 0;
288 Marker_t style = bHasHit ? kMarkerPointWHit : kMarkerPointWOHit;
289 Color_t color = mMCtrkColors.at(point.GetTrackId());
290 {
291 pCanv->cd(1);
292 auto* marker = new TMarker(point.GetZ(), point.GetX(), style);
293 marker->SetMarkerColor(color);
294 marker->Draw("same");
295
296 auto* pText = new TText(point.GetZ() + 2., point.GetX() + 2., Form("%d", point.GetStationId()));
297 pText->SetTextColor(color);
298 pText->SetTextSize(kFontSize);
299 pText->Draw("same");
300 }
301 {
302 pCanv->cd(2);
303 auto* marker = new TMarker(point.GetZ(), point.GetY(), style);
304 marker->SetMarkerColor(color);
305 marker->Draw("same");
306
307 auto* pText = new TText(point.GetZ() + 2., point.GetY() + 2., Form("%d", point.GetStationId()));
308 pText->SetTextColor(color);
309 pText->SetTextSize(kFontSize);
310 pText->Draw("same");
311 }
312 }
313
314 // Draw reconstructed tracks
315 int nRecoTracks = fvRecoTracks.size();
316 std::vector<char> vbHitUsed(fvHits.size());
317 std::vector<Color_t> vRecoTrkColors(fvHits.size()); // Reco hit ID vs. color
318 if (nRecoTracks > 0) {
319 for (int iTr = 0; iTr < nRecoTracks; ++iTr) {
320 const auto& trk = fvRecoTracks[iTr];
321 Color_t currColor = 1;
322 Style_t currStyle = trk.IsGhost() ? kLineRecoTrackGhost : kLineRecoTrackNotGhost;
323 if (trk.IsGhost()) {
324 currColor = kColorGhost[iColorGhost];
325 iColorGhost = (iColorGhost + 1) % static_cast<int>(kColorGhost.size());
326 }
327 else {
328 int iTmc = trk.GetMatchedMCTrackIndex();
329 currColor = iTmc > -1 ? mMCtrkColors[iTmc] : 1;
330 }
331
332 int nHits = trk.GetNofHits();
333 std::vector<double> trkHitX(nHits);
334 std::vector<double> trkHitY(nHits);
335 std::vector<double> trkHitZ(nHits);
336 for (int iHLoc = 0; iHLoc < nHits; ++iHLoc) {
337 int iH = trk.GetHitIndexes()[iHLoc];
338 const auto& hit = fvHits[iH];
339 vbHitUsed[iH] = true;
340 trkHitX[iHLoc] = hit.GetX();
341 trkHitY[iHLoc] = hit.GetY();
342 trkHitZ[iHLoc] = hit.GetZ();
343 vRecoTrkColors[iH] = currColor;
344 }
345
346 {
347 pCanv->cd(1);
348 auto* gr = new TGraph(nHits, trkHitZ.data(), trkHitX.data());
349 gr->SetMarkerStyle(1);
350 gr->SetMarkerColor(currColor);
351 gr->SetLineColor(currColor);
352 gr->SetLineStyle(currStyle);
353 gr->SetLineWidth(kLineWidth + 2);
354 gr->Draw("LSAME");
355 }
356 {
357 pCanv->cd(2);
358 auto* gr = new TGraph(nHits, trkHitZ.data(), trkHitY.data());
359 gr->SetMarkerStyle(1);
360 gr->SetMarkerColor(currColor);
361 gr->SetLineColor(currColor);
362 gr->SetLineStyle(currStyle);
363 gr->SetLineWidth(kLineWidth + 2);
364 gr->Draw("LSAME");
365 std::stringstream msg;
366 msg << "ID=" << trk.index << ", ";
367 msg << "MC_IDs=(";
368 for (int iTmc : trk.GetMCTrackIndexes()) {
369 msg << iTmc << ' ';
370 }
371 msg << "), ";
372 msg << "purity=" << trk.GetMaxPurity();
373 pLegendReco->AddEntry(gr, msg.str().c_str(), "l");
374 }
375 }
376 }
377 // Draw hits
378 int nHits = fvHits.size();
379 if (nHits > 0) {
380 for (int iH = 0; iH < nHits; ++iH) {
381 const auto& hit = fvHits[iH];
382 bool bFake = hit.GetBestMcPointId() == -1;
383 bool bUsed = vbHitUsed[iH];
384
385 Marker_t style = bFake ? kMarkerHitWOPoint : kMarkerHitWPoint;
386 Color_t color = bUsed ? vRecoTrkColors[iH] : 1;
387 {
388 pCanv->cd(1);
389 auto* marker = new TMarker(hit.GetZ(), hit.GetX(), style);
390 marker->SetMarkerColor(color);
391 marker->Draw("same");
392 }
393 {
394 pCanv->cd(2);
395 auto* marker = new TMarker(hit.GetZ(), hit.GetY(), style);
396 marker->SetMarkerColor(color);
397 marker->Draw("same");
398 }
399 }
400 }
401 }
402 pCanv->cd(1);
403 pLegendReco->Draw();
404 pCanv->cd(2);
405 pLegendMC->Draw();
406}
407
408// ---------------------------------------------------------------------------------------------------------------------
409//
411{
413
414 // Read reconstructed input
415 fpTSReader->ReadEvent(this->GetCurrentEvent());
416 int nHits = fvHits.size();
417 int nRecoTracks = fvRecoTracks.size();
418
421
422 // Match tracks and points
423 // Read MC tracks and points
424 int nMCTracks = 0;
425 int nMCPoints = 0;
426 if (IsMCUsed()) {
427 fpMCModule->InitEvent(this->GetCurrentEvent());
428 nMCPoints = fMCData.GetNofPoints();
430 fpMCModule->MatchHits();
431 fpMCModule->MatchTracks();
434
435 if (fbDrawEvents && nMCPoints > std::max(0, fEvtDisplayMinNofPoints)) {
436 DrawEvent();
437 }
438 }
439 LOG_IF(info, fVerbose > 2) << fName << ": Data sample consists of " << nHits << " hits, " << nRecoTracks
440 << " reco tracks, " << nMCTracks << " MC tracks, " << nMCPoints << " MC points";
441
442
443 // ************************************************************************
444 // ** Fill distributions for reconstructed tracks (including ghost ones) **
445 // ************************************************************************
446
447 for (int iTrkReco = 0; iTrkReco < nRecoTracks; ++iTrkReco) {
448 const auto& recoTrk = fvRecoTracks[iTrkReco];
449
450 // Reject tracks, which do not contain hits
451 if (recoTrk.GetNofHits() < 1) {
452 continue;
453 }
454
456
457 if (IsMCUsed()) {
458 // NOTE: The ghost status of track is now defined by its purity, thus it can still contain MC information
459 if (recoTrk.IsGhost()) {
461 }
462
463 int iTrkMC = recoTrk.GetMatchedMCTrackIndex();
464 if (iTrkMC > -1) {
465 const auto& mcTrk = fMCData.GetTrack(iTrkMC);
466 int pdg = mcTrk.GetPdgCode();
467 bool isPrimary = mcTrk.IsPrimary();
468
469 // Cut tracks, which did not leave hits in tracker
470 if (mcTrk.GetNofHits() == 0) {
471 continue;
472 }
473
474 if (isPrimary) {
476 bool bFast = mcTrk.GetP() > CbmL1Constants::MinFastMom;
477 bool bLong = mcTrk.GetTotNofStationsWithHit() == fpParameters->GetNstationsActive();
478 if (bFast) {
480 }
481 if (bLong) {
483 }
484 if (bLong && bFast) {
486 }
487 }
488 else {
490 if (mcTrk.GetP() > CbmL1Constants::MinFastMom) {
492 }
493 }
494
495 // Track distributions for different particle species
496 switch (std::abs(pdg)) {
497 case 211: // pion
499 if (isPrimary) {
502 }
503 else {
506 }
507 break;
508 case 2212: // proton
510 if (isPrimary) {
513 }
514 else {
517 }
518 break;
519 case 321: // kaon
521 if (isPrimary) {
524 }
525 else {
527 FillRecoTrack((pdg > 0) ? ETrackType::kSecKP : ETrackType::kSecKM, iTrkReco);
528 }
529 break;
530 case 11: // electron
532 if (isPrimary) {
535 }
536 else {
538 FillRecoTrack((pdg < 0) ? ETrackType::kSecEP : ETrackType::kSecEM, iTrkReco);
539 }
540 break;
541 case 13: // muon
543 if (isPrimary) {
546 }
547 else {
550 }
551 break;
552 } // switch abs(pdg): end
553 }
554 }
555 } // loop over recoTrk: end
556
557
558 // *************************************
559 // ** Fill distributions of MC-tracks **
560 // *************************************
561 if (IsMCUsed()) {
562 for (int iTrkMC = 0; iTrkMC < fMCData.GetNofTracks(); ++iTrkMC) {
563 const auto& mcTrk = fMCData.GetTrack(iTrkMC);
564
565 // ----- CUTS ON MC TRACKS
566 // Cut tracks, which did not leave hits in tracker
567 if (mcTrk.GetNofHits() == 0) {
568 continue;
569 }
570
571 // Cut tracks, which cannot be reconstructed
572 if (!mcTrk.IsReconstructable()) {
573 continue;
574 }
575 int pdg = mcTrk.GetPdgCode();
576 bool isPrimary = mcTrk.IsPrimary();
577
578 // Fill different track categories
580 if (isPrimary) {
582 bool bFast = mcTrk.GetP() > CbmL1Constants::MinFastMom;
583 bool bLong = mcTrk.GetTotNofStationsWithHit() == fpParameters->GetNstationsActive();
584 if (bFast) {
586 }
587 if (bLong) {
589 }
590 if (bLong && bFast) {
592 }
593 }
594 else {
596 if (mcTrk.GetP() > CbmL1Constants::MinFastMom) {
598 }
599 }
600
601 // Track distributions for different particle species
602 switch (std::abs(pdg)) {
603 case 211: // pion
605 if (isPrimary) {
608 }
609 else {
612 }
613 break;
614 case 2212: // proton
616 if (isPrimary) {
619 }
620 else {
623 }
624 break;
625 case 321: // kaon
627 if (isPrimary) {
630 }
631 else {
634 }
635 break;
636 case 11: // electron
638 if (isPrimary) {
641 }
642 else {
645 }
646 break;
647 case 13: // muon
649 if (isPrimary) {
652 }
653 else {
656 }
657 break;
658 } // switch abs(pdg): end
659 } // iTrkMC
660 } // IsMCUsed()
661}
662
663// ---------------------------------------------------------------------------------------------------------------------
664//
666{
668 std::vector<ETrackType> vCmpTypesGeneral = {kAll, kPrim, kSec};
669 std::vector<ETrackType> vCmpTypesPrim = {kPrim, kPrimE, kPrimMU, kPrimPI, kPrimK, kPrimPPBAR};
670 std::vector<ETrackType> vCmpTypesSec = {kSec, kSecE, kSecMU, kSecPI, kSecK, kSecPPBAR};
671 std::vector<ETrackType> vCmpTypesPions = {kAllPI, kPrimPIP, kPrimPIM, kSecPIP, kSecPIM};
672 std::vector<ETrackType> vCmpTypesKaons = {kAllK, kPrimKP, kPrimKM, kSecKP, kSecKM};
673 std::vector<ETrackType> vCmpTypesProtons = {kAllPPBAR, kPrimP, kPrimPBAR, kSecP, kSecPBAR};
674
676 auto DrawTrackDistributions = [&](TCanvas* pCanv, std::function<TH1F*(ETrackType)> Hist) {
677 pCanv->Divide(3, 2);
678 pCanv->cd(1);
679 gPad->SetLogy();
680 DrawSetOf<TH1F>(vCmpTypesGeneral, Hist);
681 pCanv->cd(2);
682 gPad->SetLogy();
683 DrawSetOf<TH1F>(vCmpTypesPrim, Hist);
684 pCanv->cd(3);
685 gPad->SetLogy();
686 DrawSetOf<TH1F>(vCmpTypesSec, Hist);
687 pCanv->cd(4);
688 gPad->SetLogy();
689 DrawSetOf<TH1F>(vCmpTypesPions, Hist);
690 pCanv->cd(5);
691 gPad->SetLogy();
692 DrawSetOf<TH1F>(vCmpTypesKaons, Hist);
693 pCanv->cd(6);
694 gPad->SetLogy();
695 DrawSetOf<TH1F>(vCmpTypesProtons, Hist);
696 };
697
699 auto DrawTrackEfficiens = [&](TCanvas* pCanv, std::function<TProfile*(ETrackType)> Prof) {
700 pCanv->Divide(3, 1);
701 pCanv->cd(1);
702 DrawSetOf<TProfile>(vCmpTypesGeneral, Prof);
703 pCanv->cd(2);
704 DrawSetOf<TProfile>(vCmpTypesPrim, Prof);
705 pCanv->cd(3);
706 DrawSetOf<TProfile>(vCmpTypesSec, Prof);
707 };
708
709
710 if (IsMCUsed()) {
711 // ** Reconstructed track distributions **
712 // Reconstructed pseudorapidity
713 auto* pc_reco_eta =
714 MakeQaObject<TCanvas>("reco_eta", "Reconstructed track pseudorapidity", kCXSIZEPX * 3, kCYSIZEPX * 2);
715 DrawTrackDistributions(pc_reco_eta, [&](ETrackType t) -> TH1F* { return fvpTrackHistograms[t]->fph_reco_eta; });
716
717 // MC pseudorapidity
718 auto* pc_reco_etaMC =
719 MakeQaObject<TCanvas>("reco_etaMC", "Reconstructed track MC pseudorapidity", kCXSIZEPX * 3, kCYSIZEPX * 2);
720 DrawTrackDistributions(pc_reco_etaMC, [&](ETrackType t) -> TH1F* { return fvpTrackHistograms[t]->fph_reco_etaMC; });
721
722 // MC momentum
723 auto* pc_reco_pMC =
724 MakeQaObject<TCanvas>("reco_pMC", "Reconstructed track MC momentum", kCXSIZEPX * 3, kCYSIZEPX * 2);
725 DrawTrackDistributions(pc_reco_pMC, [&](ETrackType t) -> TH1F* { return fvpTrackHistograms[t]->fph_reco_pMC; });
726
727 // MC rapidity
728 auto* pc_reco_yMC =
729 MakeQaObject<TCanvas>("reco_yMC", "Reconstructed track MC rapidity", kCXSIZEPX * 3, kCYSIZEPX * 2);
730 DrawTrackDistributions(pc_reco_yMC, [&](ETrackType t) -> TH1F* { return fvpTrackHistograms[t]->fph_reco_yMC; });
731
732 // ** MC track distributions **
733
734 // MC momentum
735 auto* pc_mc_pMC =
736 MakeQaObject<TCanvas>("mc_pMC", "MC reconstructable track MC momentum", kCXSIZEPX * 3, kCYSIZEPX * 2);
737 DrawTrackDistributions(pc_mc_pMC, [&](ETrackType t) -> TH1F* { return fvpTrackHistograms[t]->fph_mc_pMC; });
738
739 // MC rapidity
740 auto* pc_mc_yMC =
741 MakeQaObject<TCanvas>("mc_yMC", "MC reconstructable track MC rapidity", kCXSIZEPX * 3, kCYSIZEPX * 2);
742 DrawTrackDistributions(pc_mc_yMC, [&](ETrackType t) -> TH1F* { return fvpTrackHistograms[t]->fph_mc_yMC; });
743
744 // MC rapidity vs. MC momentum
745 // auto* pc_mc_pMC_yMC =
746 MakeQaObject<TCanvas>("mc_ptMC_yMC", "MC track MC transverse mom. vs. rapidity ", kCXSIZEPX * 3, kCYSIZEPX * 2);
747 DrawSetOf<TH2F>(vCmpTypesGeneral, [&](ETrackType t) -> TH2F* { return fvpTrackHistograms[t]->fph_reco_ptMC_yMC; });
748
749 // ** Efficiencies **
750
751 // MC momentum
752 auto* pc_eff_pMC = MakeQaObject<TCanvas>("eff_pMC", "Tracking Eff. vs. MC momentum", kCXSIZEPX * 3, kCYSIZEPX);
753 DrawTrackEfficiens(pc_eff_pMC, [&](ETrackType t) -> TProfile* { return fvpTrackHistograms[t]->fph_eff_pMC; });
754
755 auto* pc_eff_yMC = MakeQaObject<TCanvas>("eff_yMC", "Tracking Eff. vs. MC rapidity", kCXSIZEPX * 3, kCYSIZEPX);
756 DrawTrackEfficiens(pc_eff_yMC, [&](ETrackType t) -> TProfile* { return fvpTrackHistograms[t]->fph_eff_yMC; });
757
758 auto* pc_eff_thetaMC =
759 MakeQaObject<TCanvas>("eff_thetaMC", "Tracking Eff. vs. MC polar angle", kCXSIZEPX * 3, kCYSIZEPX);
760 DrawTrackEfficiens(pc_eff_thetaMC,
761 [&](ETrackType t) -> TProfile* { return fvpTrackHistograms[t]->fph_eff_thetaMC; });
762
763 auto* pc_eff_phiMC =
764 MakeQaObject<TCanvas>("eff_phiMC", "Tracking Eff. vs. MC azimuthal angle", kCXSIZEPX * 3, kCYSIZEPX);
765 DrawTrackEfficiens(pc_eff_phiMC, [&](ETrackType t) -> TProfile* { return fvpTrackHistograms[t]->fph_eff_phiMC; });
766
767 auto* pc_eff_etaMC =
768 MakeQaObject<TCanvas>("eff_etaMC", "Tracking Eff. vs. MC pseudorapidity", kCXSIZEPX * 3, kCYSIZEPX);
769 DrawTrackEfficiens(pc_eff_etaMC, [&](ETrackType t) -> TProfile* { return fvpTrackHistograms[t]->fph_eff_etaMC; });
770
771
772 // ** Pulls and residuals **
773 // NOTE: stored in a subdirectory for a given track type and point type
774 for (int iType = 0; iType < ETrackType::END; ++iType) {
775 if (fvbTrackTypeOn[iType] && fvpTrackHistograms[iType]->IsMCUsed()) {
776 fvpTrackHistograms[iType]->fpFitQaFirstHit->CreateResidualPlot();
777 fvpTrackHistograms[iType]->fpFitQaFirstHit->CreatePullPlot();
778 fvpTrackHistograms[iType]->fpFitQaLastHit->CreateResidualPlot();
779 fvpTrackHistograms[iType]->fpFitQaLastHit->CreatePullPlot();
780 }
781 }
782 }
783}
784
785// ---------------------------------------------------------------------------------------------------------------------
786//
788try {
789
790 if (fsParametersFilename.empty()) {
791 std::stringstream errMsg;
792 errMsg << "CA parameters input filename is not set. Please, provide initializer or read parameters from binary "
793 << "via OutputQa::ReadParameters(filename) from the qa macro\n";
794 throw std::runtime_error(errMsg.str());
795 }
796
798 LOG(info) << fName << ": parameters instance, reference: " << fpParameters.use_count();
799
800
801 // Turn off detectors that are not used in the reconstruction setup
802
803 fbUseMvd = fbUseMvd && (fpParameters->GetNstationsActive(ca::EDetectorID::kMvd) > 0);
804 fbUseSts = fbUseSts && (fpParameters->GetNstationsActive(ca::EDetectorID::kSts) > 0);
805 fbUseMuch = fbUseMuch && (fpParameters->GetNstationsActive(ca::EDetectorID::kMuch) > 0);
806 fbUseTrd = fbUseTrd && (fpParameters->GetNstationsActive(ca::EDetectorID::kTrd) > 0);
807 fbUseTof = fbUseTof && (fpParameters->GetNstationsActive(ca::EDetectorID::kTof) > 0);
808
809 // Turn off detectors, which hits are not presented in input tree
810
811 auto* fairManager = FairRootManager::Instance();
812 fbUseMvd = fbUseMvd && fairManager->GetObject("MvdHit");
813 fbUseSts = fbUseSts && fairManager->GetObject("StsHit");
814 fbUseMuch = fbUseMuch && fairManager->GetObject("MuchPixelHit");
815 fbUseTrd = fbUseTrd && fairManager->GetObject("TrdHit");
816 fbUseTof = fbUseTof && fairManager->GetObject("TofHit");
817
818
819 LOG(info) << fName << ": Detector subsystems used:";
820 LOG(info) << "\tMVD: " << (fbUseMvd ? "ON" : "OFF");
821 LOG(info) << "\tSTS: " << (fbUseSts ? "ON" : "OFF");
822 LOG(info) << "\tMuCh: " << (fbUseMuch ? "ON" : "OFF");
823 LOG(info) << "\tTRD: " << (fbUseTrd ? "ON" : "OFF");
824 LOG(info) << "\tTOF: " << (fbUseTof ? "ON" : "OFF");
825
826 LOG(info) << fName << ": Initializing data branches";
827
828 // Initialize IO data manager
829 if (!fpDataManager.get()) {
830 fpDataManager = std::make_shared<ca::DataManager>();
831 }
832
833 // Initialize time slice reader instance
834 fpTSReader->SetTrackingMode(fTrackingMode);
835 fpTSReader->SetDetector(ca::EDetectorID::kMvd, fbUseMvd);
836 fpTSReader->SetDetector(ca::EDetectorID::kSts, fbUseSts);
837 fpTSReader->SetDetector(ca::EDetectorID::kMuch, fbUseMuch);
838 fpTSReader->SetDetector(ca::EDetectorID::kTrd, fbUseTrd);
839 fpTSReader->SetDetector(ca::EDetectorID::kTof, fbUseTof);
840
841 fpTSReader->RegisterParameters(fpParameters);
842 fpTSReader->RegisterTracksContainer(fvRecoTracks);
843 fpTSReader->RegisterQaHitContainer(fvHits);
844 fpTSReader->RegisterHitIndexContainer(fvHitIds);
845 fpTSReader->SetSortQaHits(true);
846 if (!fpTSReader->InitRun()) {
847 throw std::runtime_error("Initialization of the CbmCaTimesliceReader object failed");
848 }
849
850 // Initialize MC module
851 if (IsMCUsed()) {
852 fpMCModule = std::make_shared<MCModule>(fVerbose, fPerformanceMode);
853 fpMCModule->SetDetector(ca::EDetectorID::kMvd, fbUseMvd);
854 fpMCModule->SetDetector(ca::EDetectorID::kSts, fbUseSts);
855 fpMCModule->SetDetector(ca::EDetectorID::kMuch, fbUseMuch);
856 fpMCModule->SetDetector(ca::EDetectorID::kTrd, fbUseTrd);
857 fpMCModule->SetDetector(ca::EDetectorID::kTof, fbUseTof);
858
859 fpMCModule->RegisterMCData(fMCData);
860 fpMCModule->RegisterRecoTrackContainer(fvRecoTracks);
861 fpMCModule->RegisterHitIndexContainer(fvHitIds);
862 fpMCModule->RegisterQaHitContainer(fvHits);
863 fpMCModule->RegisterParameters(fpParameters);
864 fpMCModule->RegisterFirstHitIndexes(fpTSReader->GetHitFirstIndexDet());
865 if (!fpMCModule->InitRun()) {
866 throw std::runtime_error("Initialization of the CbmCaMCModule object failed");
867 }
868 }
869
870 // Initialize monitor
877
878 // ----- Histograms initialization
879 //
880 auto RegisterTrackQa = [&](const char* typeName, const char* title, ETrackType type, bool bSuppressMC = false) {
881 if (!fvbTrackTypeOn[type]) {
882 return;
883 }
884 bool bUseMC = IsMCUsed() && !bSuppressMC;
885 fvsTrackTypeName[type] = typeName;
886 fvpTrackHistograms[type] = std::make_unique<TrackTypeQa>(typeName, fsPrefix.Data(), bUseMC, fpvObjList);
887 fvpTrackHistograms[type]->SetRootFolderName(fsRootFolderName + "/" + typeName);
888 fvpTrackHistograms[type]->SetTitle(title);
889 fvpTrackHistograms[type]->RegisterParameters(fpParameters);
890 fvpTrackHistograms[type]->RegisterRecoHits(fvHits);
891 fvpTrackHistograms[type]->RegisterRecoTracks(fvRecoTracks);
892 fvpTrackHistograms[type]->RegisterMCData(fMCData);
893 fvpTrackHistograms[type]->SetDrawAtt(fvTrackDrawAtts[type].fColor, fvTrackDrawAtts[type].fMarker);
894 fvpTrackHistograms[type]->Init();
895 };
896
897 //for (int i = 0; i < ETrackType::END; ++i) {
898 //LOG(info) << i << ' ' << fvpTrackHistograms[i].get() << ' ' << fvbTrackTypeOn[i];
899 //}
900
901 // TODO: Replace these parameters into the AddTrackType method!!!
902 RegisterTrackQa("all", "all", ETrackType::kAll);
903 if (IsMCUsed()) {
904 RegisterTrackQa("prim_long_fast", "primary long fast", ETrackType::kPrimLongFast);
905 RegisterTrackQa("prim_long", "primary long", ETrackType::kPrimLong);
906 RegisterTrackQa("ghost", "ghost", ETrackType::kGhost, true);
907 RegisterTrackQa("prim", "primary", ETrackType::kPrim);
908 RegisterTrackQa("prim_fast", "primary fast", ETrackType::kPrimFast);
909 RegisterTrackQa("sec", "secondary", ETrackType::kSec);
910 RegisterTrackQa("sec_fast", "secondary fast", ETrackType::kSecFast);
911 RegisterTrackQa("all_pi", "all #pi^{#pm}", ETrackType::kAllPI);
912 RegisterTrackQa("prim_pi", "primary #pi^{#pm}", ETrackType::kPrimPI);
913 RegisterTrackQa("prim_pip", "primary #pi^{#plus}", ETrackType::kPrimPIP);
914 RegisterTrackQa("prim_pim", "primary #pi^{#minus}", ETrackType::kPrimPIM);
915 RegisterTrackQa("sec_pi", "secondary #pi^{#pm}", ETrackType::kSecPI);
916 RegisterTrackQa("sec_pip", "secondary #pi^{#plus}", ETrackType::kSecPIP);
917 RegisterTrackQa("sec_pim", "secondary #pi^{#minus}", ETrackType::kSecPIM);
918 RegisterTrackQa("all_e", "all e^{#pm}", ETrackType::kAllE);
919 RegisterTrackQa("prim_e", "primary e^{#pm}", ETrackType::kPrimE);
920 RegisterTrackQa("prim_ep", "primary e^{#plus}", ETrackType::kPrimEP);
921 RegisterTrackQa("prim_em", "primary e^{#minus}", ETrackType::kPrimEM);
922 RegisterTrackQa("sec_e", "secondary e^{#pm}", ETrackType::kSecE);
923 RegisterTrackQa("sec_ep", "secondary e^{#plus}", ETrackType::kSecEP);
924 RegisterTrackQa("sec_em", "secondary e^{#minus}", ETrackType::kSecEM);
925 RegisterTrackQa("all_mu", "all #mu^{#pm}", ETrackType::kAllMU);
926 RegisterTrackQa("prim_mu", "primary #mu^{#pm}", ETrackType::kPrimMU);
927 RegisterTrackQa("prim_mup", "primary #mu^{#plus}", ETrackType::kPrimMUP);
928 RegisterTrackQa("prim_mum", "primary #mu^{#minus}", ETrackType::kPrimMUM);
929 RegisterTrackQa("sec_mu", "secondary #mu^{#pm}", ETrackType::kSecMU);
930 RegisterTrackQa("sec_mup", "secondary #mu^{#plus}", ETrackType::kSecMUP);
931 RegisterTrackQa("sec_mum", "secondary #mu^{#minus}", ETrackType::kSecMUM);
932 RegisterTrackQa("all_k", "all K^{#pm}", ETrackType::kAllK);
933 RegisterTrackQa("prim_k", "primary K^{#pm}", ETrackType::kPrimK);
934 RegisterTrackQa("prim_kp", "primary K^{#plus}", ETrackType::kPrimKP);
935 RegisterTrackQa("prim_km", "primary K^{#minus}", ETrackType::kPrimKM);
936 RegisterTrackQa("sec_k", "secondary K^{#pm}", ETrackType::kSecK);
937 RegisterTrackQa("sec_kp", "secondary K^{#plus}", ETrackType::kSecKP);
938 RegisterTrackQa("sec_km", "secondary K^{#minus}", ETrackType::kSecKM);
939 RegisterTrackQa("all_ppbar", "all p/#bar{p}", ETrackType::kAllPPBAR);
940 RegisterTrackQa("prim_ppbar", "primary p/#bar{p}", ETrackType::kPrimPPBAR);
941 RegisterTrackQa("prim_p", "primary p", ETrackType::kPrimP);
942 RegisterTrackQa("prim_pbar", "primary #bar{p}", ETrackType::kPrimPBAR);
943 RegisterTrackQa("sec_ppbar", "secondary p/#bar{p}", ETrackType::kSecPPBAR);
944 RegisterTrackQa("sec_p", "secondary p", ETrackType::kSecP);
945 RegisterTrackQa("sec_pbar", "secondary #bar{p}", ETrackType::kSecPBAR);
946 }
947
948 // Init default track types for the summary table
949 if (!fmSummaryTableEntries.size()) {
950 // clang-format off
958 };
959 // clang-format on
960 }
961
962 return kSUCCESS;
963}
964catch (const std::exception& err) {
965 LOG(error) << fName << ": Initialization failed. Reason: " << err.what();
966 return kFATAL;
967}
968
969
970// ---------------------------------------------------------------------------------------------------------------------
971//
973{
975 fvTrackDrawAtts[ETrackType::kGhost] = {kGray, 20};
976 fvTrackDrawAtts[ETrackType::kPrim] = {kGray + 3, 21};
977 fvTrackDrawAtts[ETrackType::kSec] = {kGray + 2, 25};
978
979 fvTrackDrawAtts[ETrackType::kAllPI] = {kRed - 4, 20};
980 fvTrackDrawAtts[ETrackType::kPrimPI] = {kRed - 2, 21};
981 fvTrackDrawAtts[ETrackType::kPrimPIP] = {kRed - 1, 22};
982 fvTrackDrawAtts[ETrackType::kPrimPIM] = {kRed - 3, 23};
983 fvTrackDrawAtts[ETrackType::kSecPI] = {kRed - 8, 25};
984 fvTrackDrawAtts[ETrackType::kSecPIP] = {kRed - 6, 26};
985 fvTrackDrawAtts[ETrackType::kSecPIM] = {kRed - 10, 32};
986
987 fvTrackDrawAtts[ETrackType::kAllK] = {kBlue - 4, 20};
988 fvTrackDrawAtts[ETrackType::kPrimK] = {kBlue - 2, 21};
989 fvTrackDrawAtts[ETrackType::kPrimKP] = {kBlue - 1, 22};
990 fvTrackDrawAtts[ETrackType::kPrimKM] = {kBlue - 3, 23};
991 fvTrackDrawAtts[ETrackType::kSecK] = {kBlue - 8, 25};
992 fvTrackDrawAtts[ETrackType::kSecKP] = {kBlue - 6, 26};
993 fvTrackDrawAtts[ETrackType::kSecKM] = {kBlue - 10, 32};
994
995 fvTrackDrawAtts[ETrackType::kAllPPBAR] = {kGreen - 4, 20};
996 fvTrackDrawAtts[ETrackType::kPrimPPBAR] = {kGreen - 2, 21};
997 fvTrackDrawAtts[ETrackType::kPrimP] = {kGreen - 1, 22};
998 fvTrackDrawAtts[ETrackType::kPrimPBAR] = {kGreen - 3, 23};
999 fvTrackDrawAtts[ETrackType::kSecPPBAR] = {kGreen - 8, 25};
1000 fvTrackDrawAtts[ETrackType::kSecP] = {kGreen - 6, 26};
1001 fvTrackDrawAtts[ETrackType::kSecPBAR] = {kGreen - 10, 32};
1002
1003 fvTrackDrawAtts[ETrackType::kAllE] = {kCyan - 4, 20};
1004 fvTrackDrawAtts[ETrackType::kPrimE] = {kCyan - 2, 21};
1005 fvTrackDrawAtts[ETrackType::kPrimEP] = {kCyan - 1, 22};
1006 fvTrackDrawAtts[ETrackType::kPrimEM] = {kCyan - 3, 23};
1007 fvTrackDrawAtts[ETrackType::kSecE] = {kCyan - 8, 25};
1008 fvTrackDrawAtts[ETrackType::kSecEP] = {kCyan - 6, 26};
1009 fvTrackDrawAtts[ETrackType::kSecEM] = {kCyan - 10, 32};
1010
1011 fvTrackDrawAtts[ETrackType::kAllMU] = {kMagenta - 4, 20};
1012 fvTrackDrawAtts[ETrackType::kPrimMU] = {kMagenta - 2, 21};
1013 fvTrackDrawAtts[ETrackType::kPrimMUP] = {kMagenta - 1, 22};
1014 fvTrackDrawAtts[ETrackType::kPrimMUM] = {kMagenta - 3, 23};
1015 fvTrackDrawAtts[ETrackType::kSecMU] = {kMagenta - 8, 25};
1016 fvTrackDrawAtts[ETrackType::kSecMUP] = {kMagenta - 6, 26};
1017 fvTrackDrawAtts[ETrackType::kSecMUM] = {kMagenta - 10, 32};
1018}
CA Tracking performance interface for CBM (header)
Tracking output QA-task (header)
Handles an instance of the CA-parameters as a shared pointer (header)
ECbmRecoMode
Reconstruct the full time slice or event-by-event.
Definition CbmDefs.h:162
Char_t * gr
Int_t nMCTracks
T * MakeQaObject(TString sName, TString sTitle, Args... args)
Definition CbmQaIO.h:261
TString fsPrefix
Unique prefix for all writeable root.
Definition CbmQaIO.h:159
std::shared_ptr< ObjList_t > fpvObjList
List of registered ROOT objects.
Definition CbmQaIO.h:162
TString fsRootFolderName
Name of root folder.
Definition CbmQaIO.h:157
TODO: SZh, 30.01.2023: Override THistPainter::PaintText() to add zeroes in tables.
Definition CbmQaTable.h:24
void SetNamesOfCols(const std::vector< std::string > &names)
Sets the names of table columns.
void SetCell(Int_t iRow, Int_t iCol, Double_t content, Double_t error=0.)
void SetRowName(Int_t iRow, const char *name)
Set name of a row.
void SetColWidth(Int_t width)
Definition CbmQaTable.h:74
std::string ToString(int prec, bool useErr=false) const
CbmEvent * GetCurrentEvent()
Gets pointer to current event.
Definition CbmQaTask.h:216
int GetEventNumber() const
Get current event number.
Definition CbmQaTask.h:183
bool IsMCUsed() const
Returns flag, whether MC information is used or not in the task.
Definition CbmQaTask.h:126
A CA Parameters object initialization class.
std::string ToString() const
Prints counters summary to string.
Definition CaMonitor.h:167
void IncrementCounter(ECounterKey key)
Increments key counter by 1.
Definition CaMonitor.h:99
void SetCounterName(ECounterKey key, std::string_view name)
Sets name of counter.
Definition CaMonitor.h:115
int GetCounterValue(ECounterKey key) const
Gets counter value.
Definition CaMonitor.h:81
void SetRatioKeys(const std::vector< ECounterKey > &vKeys)
Sets keys of counters, which are used as denominators for ratios.
Definition CaMonitor.h:111
A container for all external parameters of the CA tracking algorithm.
Class CbmCaPerformance is an interface to communicate between.
QA-task for CA tracking output results.
void ExecQa() override
Fills histograms from time slice or event.
static constexpr int kCYSIZEPX
Canvas size along y-axis [px].
void DrawSetOf(const std::vector< ETrackType > &vTypes, std::function< TObj *(ETrackType)> GetObj)
Utility function to draw a generic comparison of histograms from different track types.
OutputQa(int verbose, bool isMCUsed, ECbmRecoMode recoMode=ECbmRecoMode::EventByEvent)
Constructor from parameters.
TTypeArr_t< bool > fvbTrackTypeOn
Usage flag for different track types.
std::unique_ptr< TimeSliceReader > fpTSReader
Reader of the time slice.
bool fbUseTrd
is TRD used
ECbmCaTrackingMode fTrackingMode
Tracking mode.
TTypeArr_t< std::string > fvsTrackTypeName
Array of track type unique names.
void AddTrackType(ETrackType type, bool flag=true)
Adds track type.
bool fbUseMvd
is MVD used
ca::Vector< CbmL1HitDebugInfo > fvHits
void CreateSummary() override
Creates summary cavases, tables etc.
void FillRecoTrack(ETrackType type, int iTrkReco)
Fills reconstructed track by its index.
bool fbUseMuch
is MuCh used
InitStatus InitQa() override
Initialises the QA-task.
int fPerformanceMode
Performance mode.
void InitDrawingAttributes()
Defines drawing attributes for histograms of different track types.
bool fbUseTof
is TOF used
std::set< ETrackType > fmSummaryTableEntries
Which track types should be listed in the summary table.
ca::Monitor< EMonitorKey > fMonitor
TTypeArr_t< DrawAtt > fvTrackDrawAtts
Drawing attributes for track types.
void FillMCTrack(ETrackType type, int iTrkMC)
Fills MC track by its index.
int fEvtDisplayMinNofPoints
minimum number of MC points in the event display
tools::MCData fMCData
Input MC data (points and tracks)
void DrawEvent()
Draws event.
TTypeArr_t< std::unique_ptr< TrackTypeQa > > fvpTrackHistograms
Histogrammers for different track types.
std::shared_ptr< ca::Parameters< float > > fpParameters
Tracking parameters object.
static constexpr int kCXSIZEPX
Canvas size along x-axis [px].
std::shared_ptr< MCModule > fpMCModule
MC module.
std::string fsParametersFilename
std::shared_ptr< ca::DataManager > fpDataManager
Data manager.
ca::Vector< CbmL1HitId > fvHitIds
bool fbDrawEvents
flag to draw events with the event display
void Check() override
Method to check, if the QA results are acceptable.
ca::Vector< CbmL1Track > fvRecoTracks
bool fbUseSts
is STS used
const ParametersPtr_t Get(const std::string &filename)
Returns an shared pointer to the parameters instance.
static ParametersHandler * Instance()
Instance access.
int GetNofPoints() const
Gets number of points in this event/TS.
const auto & GetTrack(int idx) const
Gets a reference to MC track by its internal index.
int GetNofTracks() const
Gets number of tracks in this event/TS.
const auto & GetPoint(int idx) const
Gets a reference to MC point by its index.
const double MinFastMom
ETrackType
Enumeration fors track category.
@ kPrimPIP
primary pi+
@ kSecMUP
secondary mu+
@ kSecMUM
secondary mu-
@ kAllE
all e-/e+
@ kGhost
ghost tracks (no MC is used)
@ kPrim
primary
@ kPrimLongFast
primary long tracks (all stations in set)
@ kPrimEP
primary e+
@ kSecFast
secondary fast
@ kPrimP
primary p
@ kPrimKM
primary K-
@ kPrimPBAR
primary p-bar
@ kSecK
secondary K+/K-
@ kSecPIP
secondary pi+
@ kPrimMUM
primary mu-
@ kSecE
secondary e-/e+
@ kPrimPI
primary pi+/pi-
@ kPrimPIM
primary pi-
@ kSecEM
secondary e-
@ kAllPI
all pi+/pi-
@ kPrimK
primary K+/K-
@ kSec
secondary
@ kPrimMU
primary mu+/mu-
@ kPrimPPBAR
primary p/p-bar
@ kSecPI
secondary pi+/pi-
@ kPrimE
primary e-/e+
@ kSecP
secondary p
@ kAll
all tracks
@ kAllPPBAR
all p/p-bar
@ kSecPPBAR
secondary p/p-bar
@ kPrimFast
primary fast
@ kPrimEM
primary e-
@ kSecKM
secondary K-
@ kAllK
all K+/K-
@ kSecEP
secondary e+
@ kAllMU
all mu+/mu-
@ kPrimMUP
primary mu+
@ kSecKP
secondary K+
@ kSecMU
secondary mu+/mu-
@ kPrimLong
long primary
@ kPrimKP
primary K+
@ kSecPBAR
secondary p-bar
@ kSecPIM
secondary pi-