CbmRoot
Loading...
Searching...
No Matches
CbmCaInputQaSetup.cxx
Go to the documentation of this file.
1/* Copyright (C) 2023 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Sergei Zharko [committer] */
4
9
10#include "CbmCaInputQaSetup.h"
11
13#include "CbmL1DetectorID.h"
14#include "CbmMCDataManager.h"
15#include "CbmSetup.h"
16#include "FairRootManager.h"
17#include "TAxis.h"
18
19#include <Logger.h>
20
24
25// ---------------------------------------------------------------------------------------------------------------------
26//
27InputQaSetup::InputQaSetup(int verbose, bool isMCUsed) : CbmQaTask("CbmCaInputQaSetup", verbose, isMCUsed) {}
28
29// ---------------------------------------------------------------------------------------------------------------------
30//
32
33
34// ---------------------------------------------------------------------------------------------------------------------
35//
37{
38 if (IsMCUsed() && !fpMCEventHeader) {
39 throw std::logic_error("MC event header branch is unavailable");
40 }
41 for (int iD = 0; iD < static_cast<int>(ca::EDetectorID::END); ++iD) {
42 if (fvbUseDet[iD]) {
43 if (!fvpBrHits[iD]) {
44 throw std::logic_error(Form("Hit branch is unavailable for %s", kDetName[iD]));
45 }
46 if (IsMCUsed() && !fvpBrPoints[iD]) {
47 throw std::logic_error(Form("MC point branch is unavailable for %s", kDetName[iD]));
48 }
49 }
50 }
51}
52
53// ---------------------------------------------------------------------------------------------------------------------
54//
56{
57 auto CheckDetector = [&](ca::EDetectorID detID) {
58 if (CbmSetup::Instance()->IsActive(ToCbmModuleId(detID))) {
59 fvbUseDet[detID] = true;
60 }
61 L_(info) << fName << ": " << ToString(ToCbmModuleId(detID)) << " " << fvbUseDet[detID];
62 };
63 CheckDetector(ca::EDetectorID::kMvd);
64 CheckDetector(ca::EDetectorID::kSts);
65 CheckDetector(ca::EDetectorID::kMuch);
66 CheckDetector(ca::EDetectorID::kTrd);
67 CheckDetector(ca::EDetectorID::kTof);
68}
69
70// ---------------------------------------------------------------------------------------------------------------------
71//
73{
74 // Tracking setup by hits
75 {
76 auto* canv = MakeQaObject<TCanvas>("c_setup_hits", "Setup by hits", 1000, 1000);
77 canv->Divide(1, 2, 0.000001, 0.000001);
78 canv->cd(1);
79 gPad->SetLogz();
80 gPad->SetBottomMargin(0.000001);
81 gPad->SetRightMargin(0.05);
82 gPad->SetTitle("");
83 auto* ph_hit_xz = dynamic_cast<TH2F*>(fph_hit_xz.back()->Clone("h_hit_xz_clone"));
84 ph_hit_xz->SetTitle(";z [cm];x [cm]");
85 ph_hit_xz->Draw("col");
86 this->PutSetupNameOnPad(0.08, 0.83, 0.50, 0.89);
87 canv->cd(2);
88 gPad->SetLogz();
89 gPad->SetTopMargin(0.000001);
90 gPad->SetRightMargin(0.05);
91 gPad->SetTitle("");
92 auto* ph_hit_yz = dynamic_cast<TH2F*>(fph_hit_yz.back()->Clone("h_hit_yz_clone"));
93 ph_hit_yz->SetTitle(";z [cm];y [cm]");
94 ph_hit_yz->Draw("col");
95 this->PutSetupNameOnPad(0.08, 0.93, 0.50, 0.99);
96
97
98 auto LoBinEdge = [](TAxis* pAxis, double val) { return pAxis->GetBinLowEdge(pAxis->FindBin(val)); };
99 auto UpBinEdge = [](TAxis* pAxis, double val) { return pAxis->GetBinUpEdge(pAxis->FindBin(val)); };
100
101 for (int iDet = 0; iDet < static_cast<int>(fvpDetInterface.size()); ++iDet) {
102 if (!fvbUseDet[iDet]) {
103 continue;
104 }
105 int nSt = fvpDetInterface[iDet]->GetNtrackingStations();
106 for (int iSt = 0; iSt < nSt; ++iSt) {
107 int iStActive = fpParameters->GetStationIndexActive(iSt, static_cast<ca::EDetectorID>(iDet));
108 Color_t boxColor = (iStActive < 0) ? kGray + 1 : kOrange + 2;
109 const char* detName = Form("%s-%d", fvpDetInterface[iDet]->GetDetectorName().c_str(), iSt);
110 {
111 double zMin = LoBinEdge(fph_hit_xz.back()->GetXaxis(), fvZmin[iDet][iSt]);
112 double zMax = UpBinEdge(fph_hit_xz.back()->GetXaxis(), fvZmax[iDet][iSt]);
113 double xMin = LoBinEdge(fph_hit_xz.back()->GetYaxis(), fvXmin[iDet][iSt]);
114 double xMax = UpBinEdge(fph_hit_xz.back()->GetYaxis(), fvXmax[iDet][iSt]);
115 canv->cd(1);
116 auto* pBox = new TBox(zMin, xMin, zMax, xMax);
117 pBox->SetFillStyle(0);
118 pBox->SetLineWidth(2);
119 pBox->SetLineColor(boxColor);
120 pBox->Draw("same");
121 auto* pText = new TText(zMin, xMax, detName);
122 pText->SetTextColor(boxColor);
123 pText->SetTextSize(0.035);
124 pText->SetTextAngle(45);
125 pText->Draw("same");
126 }
127 {
128 double zMin = LoBinEdge(fph_hit_yz.back()->GetXaxis(), fvZmin[iDet][iSt]);
129 double zMax = UpBinEdge(fph_hit_yz.back()->GetXaxis(), fvZmax[iDet][iSt]);
130 double yMin = LoBinEdge(fph_hit_yz.back()->GetYaxis(), fvYmin[iDet][iSt]);
131 double yMax = UpBinEdge(fph_hit_yz.back()->GetYaxis(), fvYmax[iDet][iSt]);
132 canv->cd(2);
133 auto* pBox = new TBox(zMin, yMin, zMax, yMax);
134 pBox->SetFillStyle(0);
135 pBox->SetLineWidth(2);
136 pBox->SetLineColor(boxColor);
137 pBox->Draw("same");
138 auto* pText = new TText(zMin, yMax, detName);
139 pText->SetTextSize(0.035);
140 pText->SetTextColor(boxColor);
141 pText->SetTextAngle(45);
142 pText->Draw("same");
143 }
144 }
145 }
146 }
147}
148
149// -----------------------------------------------------------------------------------------------------------------
150//
151template<cbm::algo::ca::EDetectorID DetID>
153{
154 using Hit_t = HitTypes_t::at<DetID>;
155 using McPoint_t = PointTypes_t::at<DetID>;
156 int nHits = fvpBrHits[DetID]->GetEntriesFast();
157
158 for (int iH = 0; iH < nHits; ++iH) {
159 const auto* pHit = dynamic_cast<const Hit_t*>(fvpBrHits[DetID]->At(iH));
160 if (!pHit) {
161 LOG(warn) << fName << ": hit with iH = " << iH << " for detector " << kDetName[DetID] << " is not found";
162 }
163 auto address = pHit->GetAddress();
164
165 // skip Bmon hits (rule?)
166 if constexpr (ca::EDetectorID::kTof == DetID) {
167 if (5 == CbmTofAddress::GetSmType(address)) {
168 continue;
169 }
170 }
171
172 int iStLoc = fvpDetInterface[DetID]->GetTrackingStationIndex(address);
173 if (iStLoc < 0) {
174 continue;
175 }
176
177 int iStGeo = fpParameters->GetStationIndexGeometry(iStLoc, DetID);
178 auto xHit = pHit->GetX();
179 auto yHit = pHit->GetY();
180 auto zHit = pHit->GetZ();
181 if (fvXmin[DetID][iStLoc] > xHit) {
182 fvXmin[DetID][iStLoc] = xHit;
183 }
184 if (fvXmax[DetID][iStLoc] < xHit) {
185 fvXmax[DetID][iStLoc] = xHit;
186 }
187 if (fvYmin[DetID][iStLoc] > yHit) {
188 fvYmin[DetID][iStLoc] = yHit;
189 }
190 if (fvYmax[DetID][iStLoc] < yHit) {
191 fvYmax[DetID][iStLoc] = yHit;
192 }
193 if (fvZmin[DetID][iStLoc] > zHit) {
194 fvZmin[DetID][iStLoc] = zHit;
195 }
196 if (fvZmax[DetID][iStLoc] < zHit) {
197 fvZmax[DetID][iStLoc] = zHit;
198 }
199 if (iStGeo >= 0) {
200 fph_hit_xy[iStGeo]->Fill(xHit, yHit);
201 fph_hit_xz[iStGeo]->Fill(zHit, xHit);
202 fph_hit_yz[iStGeo]->Fill(zHit, yHit);
203 }
204 fph_hit_xz.back()->Fill(zHit, xHit);
205 fph_hit_yz.back()->Fill(zHit, yHit);
206 } // iH
207
208 if (IsMCUsed()) {
209 int nMCevents = fpMCEventList->GetNofEvents();
210 for (int iE = 0; iE < nMCevents; ++iE) {
211 int iFile = fpMCEventList->GetFileIdByIndex(iE);
212 int iEvent = fpMCEventList->GetEventIdByIndex(iE);
213 int nPoints = fvpBrPoints[DetID]->Size(iFile, iEvent);
214
215 for (int iP = 0; iP < nPoints; ++iP) {
216 const auto* pPoint = dynamic_cast<const McPoint_t*>(fvpBrPoints[DetID]->Get(iFile, iEvent, iP));
217 if (!pPoint) {
218 LOG(error) << fName << ": point with iFile=" << iFile << ", iEvent=" << iEvent << ", iP=" << iP
219 << " for detector " << kDetName[DetID] << " is not found";
220 continue;
221 }
222 auto address = pPoint->GetDetectorID();
223
224 int iStLoc = fvpDetInterface[DetID]->GetTrackingStationIndex(address);
225 if (iStLoc < 0) {
226 continue;
227 }
228
229 int iStGeo = fpParameters->GetStationIndexGeometry(iStLoc, DetID);
230 auto xPoint = pPoint->FairMCPoint::GetX();
231 auto yPoint = pPoint->FairMCPoint::GetY();
232 auto zPoint = pPoint->FairMCPoint::GetZ();
233
234 if (iStGeo >= 0) {
235 fph_point_xy[iStGeo]->Fill(xPoint, yPoint);
236 fph_point_xz[iStGeo]->Fill(zPoint, xPoint);
237 fph_point_yz[iStGeo]->Fill(zPoint, yPoint);
238 }
239 fph_point_xz.back()->Fill(zPoint, xPoint);
240 fph_point_yz.back()->Fill(zPoint, yPoint);
241 } // iP
242 } // iE
243 }
244}
245
246// ---------------------------------------------------------------------------------------------------------------------
247//
249{
250 if (fvbUseDet[ca::EDetectorID::kMvd]) {
252 }
253 if (fvbUseDet[ca::EDetectorID::kSts]) {
255 }
256 if (fvbUseDet[ca::EDetectorID::kMuch]) {
258 }
259 if (fvbUseDet[ca::EDetectorID::kTrd]) {
261 }
262 if (fvbUseDet[ca::EDetectorID::kTof]) {
264 }
265}
266
267
268// ---------------------------------------------------------------------------------------------------------------------
269//
271try {
272 LOG(info) << fName << ": initializing... ";
274
275 // Tracking parameters
277 LOG(info) << fName << ": parameters instance, reference: " << fpParameters.use_count();
278
279 auto pFairManager = FairRootManager::Instance();
280 assert(pFairManager);
281
282 auto pMcManager = IsMCUsed() ? dynamic_cast<CbmMCDataManager*>(pFairManager->GetObject("MCDataManager")) : nullptr;
283 if (IsMCUsed()) {
284 fpMCEventList = dynamic_cast<CbmMCEventList*>(pFairManager->GetObject("MCEventList."));
285 fpMCEventHeader = pMcManager->GetObject("MCEventHeader.");
286 }
287
288 fvpBrPoints.fill(nullptr);
289 fvpBrHits.fill(nullptr);
290
291 auto InitDetInterface = [&](CbmTrackingDetectorInterfaceBase* ifs, ca::EDetectorID detID) {
292 if (fvbUseDet[detID]) {
293 fvpDetInterface[detID] = ifs;
294 int nSt = ifs->GetNtrackingStations();
295 fvXmin[detID].resize(nSt, std::numeric_limits<double>::max());
296 fvXmax[detID].resize(nSt, std::numeric_limits<double>::lowest());
297 fvYmin[detID].resize(nSt, std::numeric_limits<double>::max());
298 fvYmax[detID].resize(nSt, std::numeric_limits<double>::lowest());
299 fvZmin[detID].resize(nSt, std::numeric_limits<double>::max());
300 fvZmax[detID].resize(nSt, std::numeric_limits<double>::lowest());
301 }
302 };
303
304 auto InitHitBranch = [&](const char* brName, ca::EDetectorID detID) {
305 if (fvbUseDet[detID]) {
306 fvpBrHits[detID] = dynamic_cast<TClonesArray*>(pFairManager->GetObject(brName));
307 if (!fvpBrHits[detID]) {
308 fvbUseDet[detID] = false; // Disable detectors without hits
309 }
310 }
311 return static_cast<bool>(fvpBrHits[detID]);
312 };
313
314 auto InitPointBranch = [&](const char* brName, ca::EDetectorID detID) {
315 if (IsMCUsed() && fvbUseDet[detID]) {
316 fvpBrPoints[detID] = pMcManager->InitBranch(brName);
317 }
318 };
319
320 InitDetInterface(CbmMvdTrackingInterface::Instance(), ca::EDetectorID::kMvd);
321 InitDetInterface(CbmStsTrackingInterface::Instance(), ca::EDetectorID::kSts);
322 InitDetInterface(CbmMuchTrackingInterface::Instance(), ca::EDetectorID::kMuch);
323 InitDetInterface(CbmTrdTrackingInterface::Instance(), ca::EDetectorID::kTrd);
324 InitDetInterface(CbmTofTrackingInterface::Instance(), ca::EDetectorID::kTof);
325
326 InitHitBranch("MvdHit", ca::EDetectorID::kMvd);
327 InitHitBranch("StsHit", ca::EDetectorID::kSts);
328 InitHitBranch("MuchPixelHit", ca::EDetectorID::kMuch);
329 InitHitBranch("TrdHit", ca::EDetectorID::kTrd);
330 if (!InitHitBranch("TofCalHit", ca::EDetectorID::kTof)) {
331 InitHitBranch("TofHit", ca::EDetectorID::kTof);
332 }
333
334 InitPointBranch("MvdPoint", ca::EDetectorID::kMvd);
335 InitPointBranch("StsPoint", ca::EDetectorID::kSts);
336 InitPointBranch("MuchPoint", ca::EDetectorID::kMuch);
337 InitPointBranch("TrdPoint", ca::EDetectorID::kTrd);
338 InitPointBranch("TofPoint", ca::EDetectorID::kTof);
339
340 this->CheckInit();
341
342
343 int nStGeo = fpParameters->GetNstationsGeometry();
344
345 MakeQaDirectory("hit_occupancy");
346 for (int iD = 0; iD < static_cast<int>(ca::EDetectorID::END); ++iD) {
347 if (fvbUseDet[iD]) {
348 MakeQaDirectory(Form("hit_occupancy/%s", kDetName[iD]));
349 }
350 }
351 MakeQaDirectory("point_occupancy");
352 for (int iD = 0; iD < static_cast<int>(ca::EDetectorID::END); ++iD) {
353 if (fvbUseDet[iD]) {
354 MakeQaDirectory(Form("point_occupancy/%s", kDetName[iD]));
355 }
356 }
357
358 fph_hit_xy.resize(nStGeo);
359 fph_hit_xz.resize(nStGeo + 1);
360 fph_hit_yz.resize(nStGeo + 1);
361 fph_point_xy.resize(nStGeo);
362 fph_point_xz.resize(nStGeo + 1);
363 fph_point_yz.resize(nStGeo + 1);
364
366 constexpr int kNbinsZ = 300;
367 constexpr int kNbinsX = 200;
368 constexpr int kNbinsY = 200;
369 constexpr float kMinZ = -5.;
370 constexpr float kMaxZ = 350.;
371 constexpr float kMinX = -100.;
372 constexpr float kMaxX = 100.;
373 constexpr float kMinY = -100.;
374 constexpr float kMaxY = 100.;
375
376 for (int iStGeo = 0; iStGeo <= nStGeo; ++iStGeo) {
377 auto [detID, iStLoc] = fpParameters->GetStationIndexLocal(iStGeo);
378 TString sF = "";
379 TString sN = "";
380 TString sT = "";
381 TString sNsuff = (iStGeo == nStGeo) ? "" : Form("_st_%s%d", kDetName[detID], iStLoc);
382 TString sTsuff = (iStGeo == nStGeo) ? "" : Form(" in %s station %d", kDetName[detID], iStLoc);
383
384 // Hit occupancy
385 sF = Form("hit_occupancy%s", ((iStGeo == nStGeo) ? "" : Form("/%s", kDetName[detID])));
386 sN = Form("hit_xz%s", sNsuff.Data());
387 sT = Form("hit occupancy in xz-plane%s;z_{hit} [cm];x_{hit} [cm]", sTsuff.Data());
388 fph_hit_xz[iStGeo] = MakeQaObject<TH2F>(sF + "/" + sN, sT, kNbinsZ, kMinZ, kMaxZ, kNbinsX, kMinX, kMaxX);
389
390 sN = Form("hit_yz%s", sNsuff.Data());
391 sT = Form("hit occupancy in yz-plane%s;z_{hit} [cm];y_{hit} [cm]", sTsuff.Data());
392 fph_hit_yz[iStGeo] = MakeQaObject<TH2F>(sF + "/" + sN, sT, kNbinsZ, kMinZ, kMaxZ, kNbinsY, kMinY, kMaxY);
393
394 if (iStGeo < nStGeo) {
395 sN = Form("hit_xy%s", sNsuff.Data());
396 sT = Form("hit occupancy in xy-plane%s;x_{hit} [cm];y_{hit} [cm]", sTsuff.Data());
397 fph_hit_xy[iStGeo] = MakeQaObject<TH2F>(sF + "/" + sN, sT, kNbinsX, kMinX, kMaxX, kNbinsY, kMinY, kMaxY);
398 }
399
400 if (IsMCUsed()) {
401 // Point occupancy
402 sF = Form("point_occupancy%s", ((iStGeo == nStGeo) ? "" : Form("/%s", kDetName[detID])));
403 sN = Form("point_xz%s", sNsuff.Data());
404 sT = Form("point occupancy in xz-plane%s;z_{point} [cm];x_{point} [cm]", sTsuff.Data());
405 fph_point_xz[iStGeo] = MakeQaObject<TH2F>(sF + "/" + sN, sT, kNbinsZ, kMinZ, kMaxZ, kNbinsX, kMinX, kMaxX);
406
407 sN = Form("point_yz%s", sNsuff.Data());
408 sT = Form("point occupancy in yz-plane%s;z_{point} [cm];y_{point} [cm]", sTsuff.Data());
409 fph_point_yz[iStGeo] = MakeQaObject<TH2F>(sF + "/" + sN, sT, kNbinsZ, kMinZ, kMaxZ, kNbinsY, kMinY, kMaxY);
410
411 if (iStGeo < nStGeo) {
412 sN = Form("point_xy%s", sNsuff.Data());
413 sT = Form("point occupancy in xy-plane%s;x_{point} [cm];y_{point} [cm]", sTsuff.Data());
414 fph_point_xy[iStGeo] = MakeQaObject<TH2F>(sF + "/" + sN, sT, kNbinsX, kMinX, kMaxX, kNbinsY, kMinY, kMaxY);
415 }
416 }
417 }
418
419
420 auto CreateMaterialBudgetHistograms = [&](const auto& kfSetup, const TString& dir) {
421 for (int iLayer = 0; iLayer < kfSetup.GetNofLayers(); ++iLayer) {
422 const auto& matMap{kfSetup.GetMaterial(iLayer)};
423 auto [detID, stationID] = kfSetup.GetIndexMap().GlobalToLocal(iLayer);
424 TString sN = Form("%s/mat_budget_%s_st%d", dir.Data(), kDetName[detID], stationID);
425 TString sT =
426 Form("Material budget map for %s station %d;x [cm];y [cm]; X/X_{0} [%%]", kDetName[detID], stationID);
427 auto nBins{matMap.GetNbins()};
428 auto xMin{-matMap.GetXYmax()};
429 auto xMax{+matMap.GetXYmax()};
430 auto* pHist = MakeQaObject<TH2F>(sN, sT, nBins, xMin, xMax, nBins, xMin, xMax);
431 for (int iX = 0; iX < nBins; ++iX) {
432 for (int iY = 0; iY < nBins; ++iY) {
433 pHist->SetBinContent(iX + 1, iY + 1, 100. * matMap.template GetBinThicknessX0<float>(iX, iY));
434 }
435 }
436 }
437 };
438
439 MakeQaDirectory("TrackingKFSetup");
440 MakeQaDirectory("TrackingKFSetup/geometry");
441 CreateMaterialBudgetHistograms(fpParameters->GetGeometrySetup(), "TrackingKFSetup/geometry");
442 MakeQaDirectory("TrackingKFSetup/active");
443 CreateMaterialBudgetHistograms(fpParameters->GetActiveSetup(), "TrackingKFSetup/active");
444
445 LOG(info) << fName << ": initializing... \033[1;32mDone\033[0m";
446 return kSUCCESS;
447}
448catch (const std::exception& e) {
449 LOG(error) << fName << ": initializing... \033[1;31mFailed\033[0m\nReason: " << e.what();
450 return kFATAL;
451}
#define L_(level)
QA task for tracking detector interfaces.
Handles an instance of the CA-parameters as a shared pointer (header)
std::string ToString(ECbmModuleId modId)
Definition CbmDefs.cxx:70
ECbmModuleId ToCbmModuleId(std::string modIdStr)
Definition CbmDefs.cxx:78
Implementation of L1DetectorID enum class for CBM.
Task class creating and managing CbmMCDataArray objects.
Container class for MC events with number, file and start time.
int32_t GetFileIdByIndex(uint32_t index)
File number by index @value File number for event at given index in list.
int32_t GetEventIdByIndex(uint32_t index)
Event number by index @value Event number for event at given index in list.
std::size_t GetNofEvents() const
Number of events in the list @value Number of events.
static CbmMuchTrackingInterface * Instance()
Gets pointer to the instance of the CbmMuchTrackingInterface.
static CbmMvdTrackingInterface * Instance()
Gets pointer to the instance of the CbmMvdTrackingInterface.
T * MakeQaObject(TString sName, TString sTitle, Args... args)
Definition CbmQaIO.h:261
void MakeQaDirectory(TString sName)
Definition CbmQaIO.cxx:132
bool IsMCUsed() const
Returns flag, whether MC information is used or not in the task.
Definition CbmQaTask.h:126
void PutSetupNameOnPad(double xMin, double yMin, double xMax, double yMax)
Puts setup title on the canvas.
static CbmSetup * Instance()
Definition CbmSetup.cxx:160
static CbmStsTrackingInterface * Instance()
Gets pointer to the instance of the CbmStsTrackingInterface class.
static int32_t GetSmType(uint32_t address)
static CbmTofTrackingInterface * Instance()
Gets pointer to the instance of the CbmTofTrackingInterface.
Abstract class, which should be inherited by every detecting subsystem tracking interface class.
int GetNtrackingStations() const
Gets actual number of stations, provided by the current geometry setup.
static CbmTrdTrackingInterface * Instance()
Gets pointer to the instance of the CbmTrdTrackingInterface.
A CA Parameters object initialization class.
A container for all external parameters of the CA tracking algorithm.
A QA task to analyze hit and MC point occupancy distributions in different tracking stations.
void FillHistogramsDet()
Fill histograms for a given detector type.
std::vector< TH2F * > fph_hit_xz
Hit occupancy in x-z plane vs. station.
std::shared_ptr< ca::Parameters< float > > fpParameters
Pointer to CA parameters object.
std::vector< TH2F * > fph_hit_xy
Hit occupancy in x-y plane vs. station.
DetIdArr_t< std::vector< double > > fvZmax
std::vector< TH2F * > fph_point_yz
MC point occupancy in y-z plane vs. station.
std::vector< TH2F * > fph_point_xz
MC point occupancy in x-z plane vs. station.
DetIdArr_t< std::vector< double > > fvXmin
DetIdArr_t< bool > fvbUseDet
Detector subsystem usage flag.
DetIdArr_t< std::vector< double > > fvZmin
std::string fsParametersFilename
Filename for the tracking parameters.
DetIdArr_t< std::vector< double > > fvXmax
std::vector< TH2F * > fph_hit_yz
Hit occupancy in y-z plane vs. station.
void CheckInit() const
Checks branches initialization.
void Check() override
Checks results of the QA and returns a success flag.
DetIdArr_t< TClonesArray * > fvpBrHits
Input branch for hits.
CbmMCDataObject * fpMCEventHeader
Pointer to MC event header.
InputQaSetup(int verbose, bool isMCUsed)
Constructor from parameters.
InitStatus InitQa() override
Initializes QA-task.
std::vector< TH2F * > fph_point_xy
MC point occupancy in x-z plane vs. station.
void ExecQa() override
Fills histograms.
DetIdArr_t< std::vector< double > > fvYmax
DetIdArr_t< CbmMCDataArray * > fvpBrPoints
Input branch for MC points.
DetIdArr_t< const CbmTrackingDetectorInterfaceBase * > fvpDetInterface
Pointers to the tracking detector interfaces for each subsystem.
CbmMCEventList * fpMCEventList
Pointer to MC event list.
void CheckoutDetectors()
Check out initialized detectors.
void CreateSummary() override
Creates summary cavases, tables etc.
DetIdArr_t< std::vector< double > > fvYmin
const ParametersPtr_t Get(const std::string &filename)
Returns an shared pointer to the parameters instance.
static ParametersHandler * Instance()
Instance access.
EDetectorID
Enumeration for the tracking detector subsystems in CBM-CA.
Definition CbmDefs.h:176
constexpr DetIdArr_t< const char * > kDetName
Names of detector subsystems.
std::tuple_element_t< static_cast< std::size_t >(DetID), std::tuple< Types... > > at