CbmRoot
Loading...
Searching...
No Matches
CbmCaInputQaBase.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 "CbmCaInputQaBase.h"
11
12#include "CaDefs.h"
13#include "CbmAddress.h"
14#include "CbmL1DetectorID.h"
15#include "CbmMCDataArray.h"
16#include "CbmMCEventList.h"
17#include "CbmMCTrack.h"
18#include "CbmMatch.h"
19#include "CbmMuchPixelHit.h"
20#include "CbmMuchPoint.h"
21#include "CbmMvdHit.h"
22#include "CbmMvdPoint.h"
23#include "CbmQaTable.h"
24#include "CbmQaUtil.h"
25#include "CbmStsCluster.h"
26#include "CbmStsHit.h"
27#include "CbmStsPoint.h"
28#include "CbmTimeSlice.h"
29#include "CbmTofAddress.h"
30#include "CbmTofHit.h"
31#include "CbmTofPoint.h"
33#include "CbmTrdHit.h"
34#include "CbmTrdPoint.h"
35#include "FairMCPoint.h"
36#include "FairRootManager.h"
37#include "Logger.h"
38#include "TBox.h"
39#include "TCanvas.h"
40#include "TClonesArray.h"
41#include "TEllipse.h"
42#include "TF1.h"
43#include "TFormula.h"
44#include "TGraphAsymmErrors.h"
45#include "TH1F.h"
46#include "TH2F.h"
47#include "TMath.h"
48#include "TProfile.h"
49#include "TProfile2D.h"
50#include "TStyle.h"
51
52#include <algorithm>
53#include <fstream>
54#include <iomanip>
55#include <numeric>
56#include <tuple>
57
58using namespace cbm::algo::ca::constants;
59
60namespace
61{
62 namespace ca = cbm::algo::ca;
63}
64
65// ---------------------------------------------------------------------------------------------------------------------
66//
67template<ca::EDetectorID DetID>
68CbmCaInputQaBase<DetID>::CbmCaInputQaBase(const char* name, int verbose, bool isMCUsed)
69 : CbmQaTask(name, verbose, isMCUsed)
70{
72}
73
74// ---------------------------------------------------------------------------------------------------------------------
75//
76template<ca::EDetectorID DetID>
78{
79 LOG(info) << "\n\n ** Checking the " << fName << " **\n\n";
80
81 int nSt = fpDetInterface->GetNtrackingStations();
82
83 // **************************************************************
84 // ** Basic checks, available both for real and simulated data **
85 // **************************************************************
86
87 // Fill efficiency distributions
88
89 if (IsMCUsed()) {
90 for (int iSt = 0; iSt < nSt; ++iSt) {
91 TProfile2D* effxy = fvpe_reco_eff_vs_xy[iSt];
92 for (int i = 1; i < effxy->GetNbinsX() - 1; i++) {
93 for (int j = 1; j < effxy->GetNbinsY() - 1; j++) {
94 int bin = effxy->GetBin(i, j);
95 if (effxy->GetBinEntries(bin) >= 1) {
96 fvph_reco_eff[iSt]->Fill(effxy->GetBinContent(bin));
97 fvph_reco_eff[nSt]->Fill(effxy->GetBinContent(bin));
98 }
99 }
100 }
101 }
102 }
103
104 // ----- Checks for mismatches in the ordering of the stations
105 //
106 {
107 bool res = true;
108
109 std::vector<double> vStationPos(nSt, 0.);
110 for (int iSt = 0; iSt < nSt; ++iSt) {
111 vStationPos[iSt] = fpDetInterface->GetZref(iSt);
113
114 if (!std::is_sorted(vStationPos.cbegin(), vStationPos.cend(), [](int l, int r) { return l <= r; })) {
115 if (fVerbose > 0) {
116 LOG(error) << fName << ": stations are ordered improperly along the beam axis:";
117 for (auto z : vStationPos) {
118 LOG(error) << "\t- " << z;
119 }
120 }
121 res = false;
122 }
123 StoreCheckResult("station_position_ordering", res);
124 }
125 // ----- Checks for mismatch between station and hit z positions
126 // The purpose of this block is to be ensured, that hits belong to the correct tracking station. For each tracking
127 // station a unified position z-coordinate is defined, which generally differs from the corresponding positions of
128 // reconstructed hits. This is due to non-regular position along the beam axis for each detector sensor. Thus, the
129 // positions of the hits along z-axis are distributed relatively to the position of the station in some range.
130 // If hits goes out this range, it can signal about a mismatch between hits and geometry. For limits of the range
131 // one should select large enough values to cover the difference described above and in the same time small enough
132 // to avoid overlaps with the neighboring stations.
133 // For each station, a distribution of z_{hit} - z_{st} is integrated over the defined range and scaled by the
134 // total number of entries to the distribution. If this value is smaller then unity, then some of the hits belong to
135 // another station.
136 //
137 {
138 bool res = true;
139 std::stringstream msgs;
140 msgs.precision(3);
141 for (int iSt = 0; iSt < nSt; ++iSt) {
142 int nHits = fvph_hit_station_delta_z[iSt]->GetEntries();
143 if (!nHits) {
144 LOG_IF(error, fVerbose > 0) << fName << ": station " << iSt << " does not have hits";
145 res = false;
146 continue;
147 }
148 int iBinMin = fvph_hit_station_delta_z[iSt]->FindBin(-fConfig.fMaxDiffZStHit);
149 int iBinMax = fvph_hit_station_delta_z[iSt]->FindBin(+fConfig.fMaxDiffZStHit);
150
151 auto nHitsWithin = fvph_hit_station_delta_z[iSt]->Integral(iBinMin, iBinMax);
152 if (nHitsWithin < nHits) {
153 if (!msgs.str().empty()) {
154 msgs << ", ";
155 }
156 msgs << (static_cast<double>((nHits - nHitsWithin) * 100) / nHits) << "% (st. " << iSt << ")";
157 res = false;
158 }
159 }
160 std::string msg = msgs.str().empty() ? "" : Form("Out of range z = +-%.2f cm: ", fConfig.fMaxDiffZStHit);
161 if (!msg.empty()) {
162 msg += msgs.str();
163 }
164 StoreCheckResult("station_position_hit_delta_z", res, msg);
165 }
166
167 // *******************************************************
168 // ** Additional checks, if MC information is available **
169 // *******************************************************
170
171 if (IsMCUsed()) {
172 // ----- Check efficiencies
173 // Error is raised, if any station has integrated efficiency lower then a selected threshold.
174 // Integrated efficiency is estimated with fitting the efficiency curve for hit distances (r) with a uniform
175 // distribution in the range from kEffRangeMin to kEffRangeMax, where the efficiency is assigned to be constant
176 //
177 // Fit efficiency curves
178 {
179 LOG(info) << "-- Hit efficiency integrated over hit distance from station center";
180
181 auto* pEffTable = MakeQaObject<CbmQaTable>("vs Station/eff_table", "Efficiency table", nSt, 1);
182 pEffTable->SetNamesOfCols({"Efficiency"});
183 pEffTable->SetColWidth(20);
184
185 for (int iSt = 0; iSt < nSt; ++iSt) {
186 auto eff = fvph_reco_eff[iSt]->GetMean(); // FIXME (GetMean(>>>2<<<))
187 pEffTable->SetRowName(iSt, Form("station %d", iSt));
188 pEffTable->SetCell(iSt, 0, eff);
189 bool res = CheckRange("Hit finder efficiency in station " + std::to_string(iSt), eff, fConfig.fEffThrsh, 1.000);
190 std::string msg = (res ? "" : Form("efficiency = %f lower then threshold = %f", eff, fConfig.fEffThrsh));
191 StoreCheckResult(Form("hit_efficiency_station_%d", iSt), res, msg);
192 }
193 LOG(info) << '\n' << pEffTable->ToString(3);
194 }
195
196 // ----- Checks for residuals
197 // Check hits for potential biases from central values
198 {
199 auto* pResidualsTable =
200 MakeQaObject<CbmQaTable>("vs Station/residuals_mean", "Residual mean values in different stations", nSt, 3);
201 pResidualsTable->SetNamesOfCols({"Residual(x) [cm]", "Residual(y) [cm]", "Residual(t) [ns]"});
202 pResidualsTable->SetColWidth(20);
203
204 // Fit residuals
205 for (int iSt = 0; iSt <= nSt; ++iSt) {
206
207 cbm::qa::util::SetLargeStats(fvph_res_x[iSt]);
208 cbm::qa::util::SetLargeStats(fvph_res_y[iSt]);
209 cbm::qa::util::SetLargeStats(fvph_res_t[iSt]);
210
211 pResidualsTable->SetRowName(iSt, Form("station %d", iSt));
212 pResidualsTable->SetCell(iSt, 0, fvph_res_x[iSt]->GetStdDev());
213 pResidualsTable->SetCell(iSt, 1, fvph_res_y[iSt]->GetStdDev());
214 pResidualsTable->SetCell(iSt, 2, fvph_res_t[iSt]->GetStdDev());
215 }
216 LOG(info) << '\n' << pResidualsTable->ToString(8);
217 }
218 // ----- Checks for pulls
219 //
220 // For the hit pull is defined as a ration of the difference between hit and MC-point position or time component
221 // values and an error of the hit position (or time) component, calculated in the same z-positions. In ideal case,
222 // when the resolutions are well determined for detecting stations, the pull distributions should have a RMS equal
223 // to unity. Here we allow a RMS of the pull distributions to be varied in a predefined range. If the RMS runs out
224 // this range, QA task fails.
225 {
226 std::vector<CbmQaTable*> vpPullTables = {
227 // {x, y, t}
228 MakeQaObject<CbmQaTable>("vs Station/pulls_x", "x-coordinate pull quantities in diff. stations", nSt, 2),
229 MakeQaObject<CbmQaTable>("vs Station/pulls_y", "y-coordinate pull quantities in diff. stations", nSt, 2),
230 MakeQaObject<CbmQaTable>("vs Station/pulls_t", "Time pull quantities in diff. stations", nSt, 2)};
231
232 for (auto* pullTable : vpPullTables) {
233 pullTable->SetNamesOfCols({"mean", "std.dev."});
234 pullTable->SetColWidth(20);
235 }
236
237 // NOTE: The table format is assumed: mean | 3.5 * mean err. | sigm | 3.5 * sigm err. |
238 auto DefinePullTableRow = [&](const TH1* h, CbmQaTable* table, int iSt) {
239 table->SetCell(iSt, 0, h->GetMean(), 3.5 * h->GetMeanError());
240 table->SetCell(iSt, 1, h->GetStdDev(), 3.5 * h->GetStdDevError());
241 };
242
243 for (int iSt = 0; iSt < nSt + 1; ++iSt) {
244 // Fit pull distributions for nicer representation. Fit results are not used in further checks.
245 cbm::qa::util::SetLargeStats(fvph_pull_x[iSt]);
246 cbm::qa::util::SetLargeStats(fvph_pull_y[iSt]);
247 cbm::qa::util::SetLargeStats(fvph_pull_t[iSt]);
248
249 for (auto* pullTable : vpPullTables) {
250 pullTable->SetRowName(iSt, (iSt == nSt ? "all stations" : Form("station %u", iSt)));
251 }
252 // Check the pull quality
253 {
254 auto [msg, status] = CheckRangePull(fvph_pull_x[iSt]);
255 StoreCheckResult(Form("pull_x_station_%d", iSt), status, msg);
256 DefinePullTableRow(fvph_pull_x[iSt], vpPullTables[0], iSt);
257 }
258 {
259 auto [msg, status] = CheckRangePull(fvph_pull_y[iSt]);
260 StoreCheckResult(Form("pull_y_station_%d", iSt), status, msg);
261 DefinePullTableRow(fvph_pull_y[iSt], vpPullTables[1], iSt);
262 }
263 {
264 auto [msg, status] = CheckRangePull(fvph_pull_t[iSt]);
265 StoreCheckResult(Form("pull_t_station_%d", iSt), status, msg);
266 DefinePullTableRow(fvph_pull_t[iSt], vpPullTables[2], iSt);
267 }
268 }
269
270 for (auto* pullTable : vpPullTables) {
271 LOG(info) << '\n' << pullTable->ToString(3, true);
272 }
273 }
274 } // McUsed
275
276 // Print out monitor
277 LOG(info) << fMonitor.ToString();
278}
279
280// ---------------------------------------------------------------------------------------------------------------------
281//
282template<ca::EDetectorID DetID>
283std::pair<std::string, bool> CbmCaInputQaBase<DetID>::CheckRangePull(TH1* h) const
284{
285 constexpr double factor = 3.5;
286 // Function to check overflow and underflow
287 constexpr auto Check = [&](double val, double err, double min, double max) -> std::pair<std::string, bool> {
288 std::stringstream msg;
289 bool res;
290 if (val + factor * err < min) {
291 msg << "underflow: " << val << " < " << min << " - " << factor << " x " << err;
292 res = false;
293 }
294 if (val - factor * err > max) {
295 msg << "overflow: " << val << " > " << max << " + " << factor << " x " << err;
296 res = false;
297 }
298 return std::make_pair(msg.str(), res);
299 };
300
301 std::pair<std::string, bool> ret = {"", true};
302 if (h->GetEntries() >= 10) {
303 // Mean check
304 {
305 auto [msg, res] = Check(h->GetMean(), h->GetMeanError(), -fConfig.fPullMeanThrsh, +fConfig.fPullMeanThrsh);
306 if (!msg.empty()) {
307 ret.first += "mean ";
308 ret.first += msg;
309 ret.second = res;
310 }
311 }
312 // Sigma check
313 {
314 auto [msg, res] =
315 Check(h->GetStdDev(), h->GetStdDevError(), 1. - fConfig.fPullWidthThrsh, 1. + fConfig.fPullWidthThrsh);
316 if (!msg.empty()) {
317 if (!ret.first.empty()) {
318 ret.first += "; ";
319 }
320 ret.first += "std.dev. ";
321 ret.first += msg;
322 ret.second = res;
323 }
324 }
325 }
326 return ret;
327}
328
329
330// ---------------------------------------------------------------------------------------------------------------------
331//
332template<ca::EDetectorID DetID>
334{
335 // Vectors with pointers to histograms
336 fvph_hit_xy.clear();
337 fvph_hit_zx.clear();
338 fvph_hit_zy.clear();
339
340 fvph_hit_station_delta_z.clear();
341
342 fvph_hit_dx.clear();
343 fvph_hit_dy.clear();
344 fvph_hit_du.clear();
345 fvph_hit_dv.clear();
346 fvph_hit_kuv.clear();
347 fvph_hit_dt.clear();
348
349 fvph_n_points_per_hit.clear();
350
351 fvph_point_xy.clear();
352 fvph_point_zx.clear();
353 fvph_point_zy.clear();
354
355 fvph_point_hit_delta_z.clear();
356
357 fvph_res_x.clear();
358 fvph_res_y.clear();
359 fvph_res_u.clear();
360 fvph_res_v.clear();
361 fvph_res_t.clear();
362
363 fvph_pull_x.clear();
364 fvph_pull_y.clear();
365 fvph_pull_u.clear();
366 fvph_pull_v.clear();
367 fvph_pull_t.clear();
368
369 fvph_res_x_vs_x.clear();
370 fvph_res_y_vs_y.clear();
371 fvph_res_u_vs_u.clear();
372 fvph_res_v_vs_v.clear();
373 fvph_res_t_vs_t.clear();
374
375 fvph_pull_x_vs_x.clear();
376 fvph_pull_y_vs_y.clear();
377 fvph_pull_u_vs_u.clear();
378 fvph_pull_v_vs_v.clear();
379 fvph_pull_t_vs_t.clear();
380
381 fvpe_reco_eff_vs_xy.clear();
382 fvph_reco_eff.clear();
383}
384
385// ---------------------------------------------------------------------------------------------------------------------
386//
387template<ca::EDetectorID DetID>
389{
390 const int nSt = fpDetInterface->GetNtrackingStations();
391 const int nHits = fpHits->GetEntriesFast();
392 const int nMCevents = (IsMCUsed()) ? fpMCEventList->GetNofEvents() : -1;
393
394 // TODO: SZh 06.09.2023: Probably, this approach can fail, if there are several input files are used. Thus I propose
395 // to use unordered_map with a CbmLink key type.
396 std::vector<std::vector<std::vector<int>>>
397 vNofHitsPerMcTrack; // Number of hits per MC track per station in different MC events
398
399 if (IsMCUsed()) {
400 vNofHitsPerMcTrack.resize(nMCevents);
401 for (int iE = 0; iE < nMCevents; ++iE) {
402 int iFile = fpMCEventList->GetFileIdByIndex(iE);
403 int iEvent = fpMCEventList->GetEventIdByIndex(iE);
404 int nMcTracks = fpMCTracks->Size(iFile, iEvent);
405 vNofHitsPerMcTrack[iE].resize(nSt);
406 for (int iSt = 0; iSt < nSt; iSt++) {
407 vNofHitsPerMcTrack[iE][iSt].resize(nMcTracks, 0);
408 }
409 }
410 }
411
412 for (int iH = 0; iH < nHits; ++iH) {
413 fHitQaData.Reset();
414 fHitQaData.SetHitIndex(iH);
415
416 const auto* pHit = dynamic_cast<const CbmPixelHit*>(fpHits->At(iH));
417 if (!pHit) {
418 LOG(error) << fName << ": hit with iH = " << iH << " is not an CbmStsHit (dynamic cast failed)";
419 continue;
420 }
421
422 if constexpr (ca::EDetectorID::kTof == DetID) {
423 auto address = pHit->GetAddress();
424 if (5 == CbmTofAddress::GetSmType(address)) {
425 continue;
426 } // Skip Bmon hits
427 }
428
429 fMonitor.IncrementCounter(EMonitorKey::kHit);
430 fMonitor.IncrementCounter(EMonitorKey::kHitAccepted);
431
432 // *************************
433 // ** Reconstructed hit QA
434
435 // Hit station geometry info
436 int iSt = fpDetInterface->GetTrackingStationIndex(pHit);
437 if (iSt < 0) {
438 continue;
439 }
440
441 if (iSt >= nSt) {
442 LOG(error) << fName << ": index of station (" << iSt << ") is out of range for hit with id = " << iH;
443 continue;
444 }
445
446 auto [stPhiU, stPhiV] = fpDetInterface->GetStereoAnglesSensor(pHit->GetAddress());
447
448 fHitQaData.SetPhiU(stPhiU);
449 fHitQaData.SetPhiV(stPhiV);
450 fHitQaData.SetHitX(pHit->GetX());
451 fHitQaData.SetHitY(pHit->GetY());
452 fHitQaData.SetHitZ(pHit->GetZ());
453 fHitQaData.SetHitTime(pHit->GetTime());
454
455 fHitQaData.SetHitDx(pHit->GetDx());
456 fHitQaData.SetHitDy(pHit->GetDy());
457 fHitQaData.SetHitDxy(pHit->GetDxy());
458 fHitQaData.SetHitTimeError(pHit->GetTimeError());
459 fHitQaData.SetStationID(iSt);
460
461 // Per station distributions
462 fvph_hit_xy[iSt]->Fill(fHitQaData.GetHitX(), fHitQaData.GetHitY());
463 fvph_hit_zx[iSt]->Fill(fHitQaData.GetHitZ(), fHitQaData.GetHitX());
464 fvph_hit_zy[iSt]->Fill(fHitQaData.GetHitZ(), fHitQaData.GetHitY());
465
466 fvph_hit_station_delta_z[iSt]->Fill(fHitQaData.GetHitZ() - fpDetInterface->GetZref(iSt));
467
468 fvph_hit_dx[iSt]->Fill(fHitQaData.GetHitDx());
469 fvph_hit_dy[iSt]->Fill(fHitQaData.GetHitDy());
470 fvph_hit_du[iSt]->Fill(fHitQaData.GetHitDu());
471 fvph_hit_dv[iSt]->Fill(fHitQaData.GetHitDv());
472 fvph_hit_kuv[iSt]->Fill(fHitQaData.GetHitRuv());
473 fvph_hit_dt[iSt]->Fill(fHitQaData.GetHitTimeError());
474
475 // Sum over station distributions
476 fvph_hit_xy[nSt]->Fill(fHitQaData.GetHitX(), fHitQaData.GetHitY());
477 fvph_hit_zx[nSt]->Fill(fHitQaData.GetHitZ(), fHitQaData.GetHitX());
478 fvph_hit_zy[nSt]->Fill(fHitQaData.GetHitZ(), fHitQaData.GetHitY());
479
480 fvph_hit_station_delta_z[nSt]->Fill(fHitQaData.GetHitZ() - fpDetInterface->GetZref(iSt));
481
482 fvph_hit_dx[nSt]->Fill(fHitQaData.GetHitDx());
483 fvph_hit_dy[nSt]->Fill(fHitQaData.GetHitDy());
484 fvph_hit_du[nSt]->Fill(fHitQaData.GetHitDu());
485 fvph_hit_dv[nSt]->Fill(fHitQaData.GetHitDv());
486 fvph_hit_kuv[nSt]->Fill(fHitQaData.GetHitRuv());
487 fvph_hit_dt[nSt]->Fill(fHitQaData.GetHitTimeError());
488
489 // **********************
490 // ** MC information QA
491
492 if (IsMCUsed()) {
493 const auto* pHitMatch = dynamic_cast<CbmMatch*>(fpHitMatches->At(iH));
494 assert(pHitMatch);
495
496 // Evaluate number of hits per one MC point
497 int nMCpoints = 0; // Number of MC points for this hit
498 for (int iLink = 0; iLink < pHitMatch->GetNofLinks(); ++iLink) {
499 const auto& link = pHitMatch->GetLink(iLink);
500
501 int iP = link.GetIndex(); // Index of MC point
502
503 // Skip noisy links
504 if (iP < 0) {
505 continue;
506 }
507
508 ++nMCpoints;
509
510 int iE = fpMCEventList->GetEventIndex(link);
511 // TODO: debug
512 //if (iE < 0) continue;
513
514 if (iE < 0 || iE >= nMCevents) {
515 LOG(error) << fName << ": id of MC event is out of range (hit id = " << iH << ", link id = " << iLink
516 << ", event id = " << iE << ", mc point ID = " << iP << ')';
517 continue;
518 }
519
520 // matched point
521 const auto* pMCPoint = dynamic_cast<const Point_t*>(fpMCPoints->Get(link));
522 if (!pMCPoint) {
523 LOG(error) << fName << ": MC point object does not exist for hit " << iH;
524 continue;
525 }
526
527 int iTr = pMCPoint->GetTrackID();
528
529 if (iTr >= static_cast<int>(vNofHitsPerMcTrack[iE][iSt].size())) {
530 LOG(error) << fName << ": index of MC track is out of range (hit id = " << iH << ", link id = " << iLink
531 << ", event id = " << iE << ", track id = " << iTr << ')';
532 continue;
533 }
534
535 vNofHitsPerMcTrack[iE][iSt][iTr]++;
536 }
537
538 fvph_n_points_per_hit[iSt]->Fill(nMCpoints);
539 fvph_n_points_per_hit[nSt]->Fill(nMCpoints);
540
541 //
542 // Fill the following histograms exclusively for isolated hits with a one-to-one correspondence to the mc track
543 // The mc track must satisfy the cuts
544 //
545
546 if (pHitMatch->GetNofLinks() != 1) {
547 continue;
548 }
549
550 // The best link to in the match (probably, the cut on nMCpoints is meaningless)
551 assert(pHitMatch->GetNofLinks() > 0); // Should be always true due to the cut above
552 const auto& bestPointLink = pHitMatch->GetMatchedLink();
553
554 // Skip noisy links
555 if (bestPointLink.GetIndex() < 0) {
556 continue;
557 }
558
559 // Point matched by the best link
560 const auto* pMCPoint = dynamic_cast<const Point_t*>(fpMCPoints->Get(bestPointLink));
561 fHitQaData.SetPointID(bestPointLink.GetIndex(), bestPointLink.GetEntry(), bestPointLink.GetFile());
562 if (!pMCPoint) {
563 LOG(error) << fName << ": MC point object does not exist for hit " << iH;
564 continue;
565 }
566
567 // MC track for this point
568 CbmLink bestTrackLink = bestPointLink;
569 bestTrackLink.SetIndex(pMCPoint->GetTrackID());
570 const auto* pMCTrack = dynamic_cast<const CbmMCTrack*>(fpMCTracks->Get(bestTrackLink));
571 if (!pMCTrack) {
572 LOG(error) << fName << ": MC track object does not exist for hit " << iH << " and link: ";
573 continue;
574 }
575
576 double t0MC = fpMCEventList->GetEventTime(bestPointLink);
577 if (t0MC < 0) {
578 LOG(error) << fName << ": MC time zero is lower then 0 ns: " << t0MC;
579 continue;
580 }
581
582 // cut on the mc track quality
583 if (!IsTrackSelected(pMCTrack, pMCPoint)) {
584 continue;
585 }
586
587 { // skip the case when the mc point participates to several hits
588 int iE = fpMCEventList->GetEventIndex(bestTrackLink);
589 if (vNofHitsPerMcTrack[iE][iSt][pMCPoint->GetTrackID()] != 1) {
590 continue;
591 }
592 }
593
594 // ----- MC point properties
595 //
596 double mass = pMCTrack->GetMass();
597 //double charge = pMCTrack->GetCharge();
598 //double pdg = pMCTrack->GetPdgCode();
599
600 // Entrance position and time
601 // NOTE: SZh 04.09.2023: Methods GetX(), GetY(), GetZ() for MVD, STS, MUCH, TRD and TOF always return
602 // positions of track at entrance to the active volume.
603 double xMC = pMCPoint->FairMCPoint::GetX();
604 double yMC = pMCPoint->FairMCPoint::GetY();
605 double zMC = pMCPoint->FairMCPoint::GetZ();
606 double tMC = pMCPoint->GetTime() + t0MC;
607
608 // MC point entrance momenta
609 // NOTE: SZh 04.09.2023: Methods GetPx(), GetPy(), GetPz() for MVD, STS, MUCH, TRD and TOF always return
610 // the momentum components of track at entrance to the active volume.
611 double pxMC = pMCPoint->GetPx();
612 double pyMC = pMCPoint->GetPy();
613 double pzMC = pMCPoint->GetPz();
614 double pMC = sqrt(pxMC * pxMC + pyMC * pyMC + pzMC * pzMC);
615
616 // MC point exit momenta
617 // double pxMCo = pMCPoint->GetPxOut();
618 // double pyMCo = pMCPoint->GetPyOut();
619 // double pzMCo = pMCPoint->GetPzOut();
620 // double pMCo = sqrt(pxMCo * pxMCo + pyMCo * pyMCo + pzMCo * pzMCo);
621
622 // Position and time shifted to z-coordinate of the hit
623 // TODO: SZh 04.09.2023: Probably, the approximation
624 double shiftZ = fHitQaData.GetHitZ() - zMC; // Difference between hit and point z positions
625 double xMCs = xMC + shiftZ * pxMC / pzMC;
626 double yMCs = yMC + shiftZ * pyMC / pzMC;
627 double tMCs = tMC + shiftZ / (pzMC * phys::SpeedOfLight) * sqrt(mass * mass + pMC * pMC);
628
629 fHitQaData.SetPointTime(tMCs);
630 fHitQaData.SetPointX(xMCs);
631 fHitQaData.SetPointY(yMCs);
632 fHitQaData.SetPointZ(fHitQaData.GetHitZ());
633
634 double zRes = kNAN;
635 if constexpr (ca::EDetectorID::kTof == DetID) {
636 zRes = fHitQaData.GetHitZ() - pMCPoint->GetZ();
637 }
638 else {
639 zRes = fHitQaData.GetHitZ() - 0.5 * (pMCPoint->GetZ() + pMCPoint->GetZOut());
640 }
641 fvph_point_hit_delta_z[iSt]->Fill(zRes);
642
643 double xRes = fHitQaData.GetResidualX();
644 double yRes = fHitQaData.GetResidualY();
645 double uRes = fHitQaData.GetResidualU();
646 double vRes = fHitQaData.GetResidualV();
647 double tRes = fHitQaData.GetResidualTime();
648
649 double xPull = fHitQaData.GetPullX();
650 double yPull = fHitQaData.GetPullY();
651 double uPull = fHitQaData.GetPullU();
652 double vPull = fHitQaData.GetPullV();
653 double tPull = fHitQaData.GetPullTime();
654
655 fvph_res_x[iSt]->Fill(xRes);
656 fvph_res_y[iSt]->Fill(yRes);
657 fvph_res_u[iSt]->Fill(uRes);
658 fvph_res_v[iSt]->Fill(vRes);
659 fvph_res_t[iSt]->Fill(tRes);
660
661 fvph_pull_x[iSt]->Fill(xPull);
662 fvph_pull_y[iSt]->Fill(yPull);
663 fvph_pull_u[iSt]->Fill(uPull);
664 fvph_pull_v[iSt]->Fill(vPull);
665 fvph_pull_t[iSt]->Fill(tPull);
666
667 fvph_res_x_vs_x[iSt]->Fill(fHitQaData.GetPointX(), xRes);
668 fvph_res_y_vs_y[iSt]->Fill(fHitQaData.GetPointY(), yRes);
669 fvph_res_u_vs_u[iSt]->Fill(fHitQaData.GetPointU(), uRes);
670 fvph_res_v_vs_v[iSt]->Fill(fHitQaData.GetPointV(), vRes);
671 fvph_res_t_vs_t[iSt]->Fill(fHitQaData.GetPointTime(), tRes);
672
673 fvph_pull_x_vs_x[iSt]->Fill(fHitQaData.GetPointX(), xPull);
674 fvph_pull_y_vs_y[iSt]->Fill(fHitQaData.GetPointY(), yPull);
675 fvph_pull_u_vs_u[iSt]->Fill(fHitQaData.GetPointU(), uPull);
676 fvph_pull_v_vs_v[iSt]->Fill(fHitQaData.GetPointV(), vPull);
677 fvph_pull_t_vs_t[iSt]->Fill(fHitQaData.GetPointTime(), tPull);
678
679 // fill summary histos
680
681 fvph_point_hit_delta_z[nSt]->Fill(zRes);
682
683 fvph_res_x[nSt]->Fill(xRes);
684 fvph_res_y[nSt]->Fill(yRes);
685 fvph_res_u[nSt]->Fill(uRes);
686 fvph_res_v[nSt]->Fill(vRes);
687 fvph_res_t[nSt]->Fill(tRes);
688
689 fvph_pull_x[nSt]->Fill(xPull);
690 fvph_pull_y[nSt]->Fill(yPull);
691 fvph_pull_u[nSt]->Fill(uPull);
692 fvph_pull_v[nSt]->Fill(vPull);
693 fvph_pull_t[nSt]->Fill(tPull);
694
695 fvph_res_x_vs_x[nSt]->Fill(fHitQaData.GetPointX(), xRes);
696 fvph_res_y_vs_y[nSt]->Fill(fHitQaData.GetPointY(), yRes);
697 fvph_res_u_vs_u[nSt]->Fill(fHitQaData.GetPointU(), uRes);
698 fvph_res_v_vs_v[nSt]->Fill(fHitQaData.GetPointV(), vRes);
699 fvph_res_t_vs_t[nSt]->Fill(fHitQaData.GetPointTime(), tRes);
700
701 fvph_pull_x_vs_x[nSt]->Fill(fHitQaData.GetPointX(), xPull);
702 fvph_pull_y_vs_y[nSt]->Fill(fHitQaData.GetPointY(), yPull);
703 fvph_pull_u_vs_u[nSt]->Fill(fHitQaData.GetPointU(), uPull);
704 fvph_pull_v_vs_v[nSt]->Fill(fHitQaData.GetPointV(), vPull);
705 fvph_pull_t_vs_t[nSt]->Fill(fHitQaData.GetPointTime(), tPull);
706 }
707 FillHistogramsPerHit();
708 } // loop over hits: end
709
710 // Fill hit reconstruction efficiencies
711 if (IsMCUsed()) {
712 for (int iE = 0; iE < nMCevents; ++iE) {
713 int iFile = fpMCEventList->GetFileIdByIndex(iE);
714 int iEvent = fpMCEventList->GetEventIdByIndex(iE);
715 int nPoints = fpMCPoints->Size(iFile, iEvent);
716 int nTracks = fpMCTracks->Size(iFile, iEvent);
717
718 // If efficiency for the track at the station is already evaluated
719 std::vector<std::vector<bool>> vIsTrackProcessed(nSt);
720 for (int iSt = 0; iSt < nSt; iSt++) {
721 vIsTrackProcessed[iSt].resize(nTracks, 0);
722 }
723
724 for (int iP = 0; iP < nPoints; ++iP) {
725 fHitQaData.Reset();
726 fHitQaData.SetPointID(iP, iEvent, iFile);
727 fMonitor.IncrementCounter(EMonitorKey::kMcPoint);
728
729 const auto* pMCPoint = dynamic_cast<const Point_t*>(fpMCPoints->Get(iFile, iEvent, iP));
730 if (!pMCPoint) {
731 LOG(error) << fName << ": MC point does not exist for iFile = " << iFile << ", iEvent = " << iEvent
732 << ", iP = " << iP;
733 continue;
734 }
735
736 int address = pMCPoint->GetDetectorID();
737 int iSt = fpDetInterface->GetTrackingStationIndex(pMCPoint);
738 if (iSt < 0) {
739 continue; // NOTE: legit, such sensors must be ignored in tracking
740 }
741
742 if (iSt >= nSt) {
743 LOG(error) << fName << ": MC point for FEI = " << iFile << ", " << iEvent << ", " << iP << " and address "
744 << address << " has wrong station index: iSt = " << iSt;
745 fMonitor.IncrementCounter(EMonitorKey::kMcPointWrongStation);
746 continue;
747 }
748
749 // NOTE: SZh 04.09.2023: Methods GetX(), GetY(), GetZ() for MVD, STS, MUCH, TRD and TOF always return
750 // positions of track at entrance to the active volume.
751 fHitQaData.SetPointX(pMCPoint->FairMCPoint::GetX());
752 fHitQaData.SetPointY(pMCPoint->FairMCPoint::GetY());
753 fHitQaData.SetPointZ(pMCPoint->FairMCPoint::GetZ());
754
755 fvph_point_xy[iSt]->Fill(fHitQaData.GetPointX(), fHitQaData.GetPointY());
756 fvph_point_zx[iSt]->Fill(fHitQaData.GetPointZ(), fHitQaData.GetPointX());
757 fvph_point_zy[iSt]->Fill(fHitQaData.GetPointZ(), fHitQaData.GetPointY());
758
759 fvph_point_xy[nSt]->Fill(fHitQaData.GetPointX(), fHitQaData.GetPointY());
760 fvph_point_zx[nSt]->Fill(fHitQaData.GetPointZ(), fHitQaData.GetPointX());
761 fvph_point_zy[nSt]->Fill(fHitQaData.GetPointZ(), fHitQaData.GetPointY());
762
763 int iTr = pMCPoint->GetTrackID();
764
765 if (iTr >= nTracks) {
766 LOG(error) << fName << ": index of MC track is out of range (point id = " << iP << ", event id = " << iE
767 << ", track id = " << iTr << ')';
768 continue;
769 }
770 const auto* pMCTrack = dynamic_cast<const CbmMCTrack*>(fpMCTracks->Get(iFile, iEvent, iTr));
771
772 fHitQaData.SetIfTrackHasHits(vNofHitsPerMcTrack[iE][iSt][iTr] > 0);
773 fHitQaData.SetIfTrackSelected(IsTrackSelected(pMCTrack, pMCPoint));
774
775 if (!pMCTrack) {
776 LOG(error) << fName << ": null MC track pointer for file id = " << iFile << ", event id = " << iEvent
777 << ", track id = " << iTr;
778 continue;
779 }
780
781 // check efficiency only once per mc track per station
782 if (vIsTrackProcessed[iSt][iTr]) {
783 continue;
784 }
785
786 vIsTrackProcessed[iSt][iTr] = true;
787
788 // cut on the mc track quality
789 if (!fHitQaData.GetIfTrackSelected()) {
790 continue;
791 }
792
793 // Conditions under which point is accounted as reconstructed: point
794 bool ifTrackHasHits = fHitQaData.GetIfTrackHasHits();
795
796 fvpe_reco_eff_vs_xy[iSt]->Fill(fHitQaData.GetPointX(), fHitQaData.GetPointY(), ifTrackHasHits);
797 fvpe_reco_eff_vs_xy[nSt]->Fill(fHitQaData.GetPointX(), fHitQaData.GetPointY(), ifTrackHasHits);
798 } // loop over MC-points: end
799 } // loop over MC-events: end
800 } // MC is used: end
801}
802
803// ---------------------------------------------------------------------------------------------------------------------
804//
805template<ca::EDetectorID DetID>
807{
808 // ----- Specific configuration initialization
809 fConfig = ReadSpecificConfig<CbmCaInputQaBase<DetID>::Config>().value_or(Config{});
810
811
812 LOG_IF(fatal, !fpDetInterface) << "\033[1;31m" << fName << ": tracking detector interface is undefined\033[0m";
813
814 // FairRootManager
815 auto* pFairRootManager = FairRootManager::Instance();
816 LOG_IF(fatal, !pFairRootManager) << "\033[1;31m" << fName << ": FairRootManager instance is a null pointer\033[0m";
817
818 // Time-slice
819 fpTimeSlice = dynamic_cast<CbmTimeSlice*>(pFairRootManager->GetObject("TimeSlice."));
820 LOG_IF(fatal, !fpTimeSlice) << "\033[1;31m" << fName << ": time-slice branch is not found\033[0m";
821
822 // ----- Reconstructed data-branches initialization
823
824 // Hits container
825 if constexpr (DetID == EDetectorID::kTof) {
826 fpHits = dynamic_cast<TClonesArray*>(pFairRootManager->GetObject("TofCalHit"));
827 if (!fpHits) {
828 fpHits = dynamic_cast<TClonesArray*>(pFairRootManager->GetObject("TofHit"));
829 }
830 }
831 else {
832 fpHits = dynamic_cast<TClonesArray*>(pFairRootManager->GetObject(cbm::ca::kDetHitBrName[DetID]));
833 }
834
835 LOG_IF(fatal, !fpHits) << "\033[1;31m" << fName << ": container of reconstructed hits is not found\033[0m";
836
837 // ----- MC-information branches initialization
838 if (IsMCUsed()) {
839 // MC manager
840 fpMCDataManager = dynamic_cast<CbmMCDataManager*>(pFairRootManager->GetObject("MCDataManager"));
841 LOG_IF(fatal, !fpMCDataManager) << "\033[1;31m" << fName << ": MC data manager branch is not found\033[0m";
842
843 // MC event list
844 fpMCEventList = dynamic_cast<CbmMCEventList*>(pFairRootManager->GetObject("MCEventList."));
845 LOG_IF(fatal, !fpMCEventList) << "\033[1;31m" << fName << ": MC event list branch is not found\033[0m";
846
847 // MC tracks
848 fpMCTracks = fpMCDataManager->InitBranch("MCTrack");
849 LOG_IF(fatal, !fpMCTracks) << "\033[1;31m" << fName << ": MC track branch is not found\033[0m";
850
851 // MC points
852 fpMCPoints = fpMCDataManager->InitBranch(cbm::ca::kDetPointBrName[DetID]);
853 LOG_IF(fatal, !fpMCTracks) << "\033[1;31m" << fName << ": MC point branch is not found\033[0m";
854
855 // Hit matches
856 const char* hitMatchName = Form("%sMatch", cbm::ca::kDetHitBrName[DetID]);
857 fpHitMatches = dynamic_cast<TClonesArray*>(pFairRootManager->GetObject(hitMatchName));
858 LOG_IF(fatal, !fpHitMatches) << "\033[1;31m]" << fName << ": hit match branch is not found\033[0m";
859 }
860
861 // get hit distribution ranges from the geometry
862 int nSt = fpDetInterface->GetNtrackingStations();
863
864 frXmin.resize(nSt + 1, 0.);
865 frXmax.resize(nSt + 1, 0.);
866
867 frYmin.resize(nSt + 1, 0.);
868 frYmax.resize(nSt + 1, 0.);
869
870 frZmin.resize(nSt + 1, 0.);
871 frZmax.resize(nSt + 1, 0.);
872
873 for (int i = nSt; i >= 0; --i) {
874 { // read the ranges for station i
875 int j = (i == nSt ? 0 : i);
876
877 frXmin[i] = fpDetInterface->GetXmin(j);
878 frXmax[i] = fpDetInterface->GetXmax(j);
879
880 frYmin[i] = fpDetInterface->GetYmin(j);
881 frYmax[i] = fpDetInterface->GetYmax(j);
882
883 frZmin[i] = fpDetInterface->GetZmin(j);
884 frZmax[i] = fpDetInterface->GetZmax(j);
885 }
886
887 // update overall ranges
888
889 frXmin[nSt] = std::min(frXmin[nSt], frXmin[i]);
890 frXmax[nSt] = std::max(frXmax[nSt], frXmax[i]);
891
892 frYmin[nSt] = std::min(frYmin[nSt], frYmin[i]);
893 frYmax[nSt] = std::max(frYmax[nSt], frYmax[i]);
894
895 frZmin[nSt] = std::min(frZmin[nSt], frZmin[i]);
896 frZmax[nSt] = std::max(frZmax[nSt], frZmax[i]);
897 }
898
899 // add 5% margins for a better view
900
901 for (int i = 0; i <= nSt; ++i) {
902 double dx = 0.05 * fabs(frXmax[i] - frXmin[i]);
903 frXmin[i] -= dx;
904 frXmax[i] += dx;
905
906 double dy = 0.05 * fabs(frYmax[i] - frYmin[i]);
907 frYmin[i] -= dy;
908 frYmax[i] += dy;
909
910 if constexpr (ca::EDetectorID::kMuch == DetID) {
911 frZmin[i] -= 40;
912 frZmax[i] += 40;
913 }
914 else {
915 double dz = 0.05 * fabs(frZmax[i] - frZmin[i]);
916 frZmin[i] -= dz;
917 frZmax[i] += dz;
918 }
919 }
920
921 // ----- Monitor initialization
922 //
923 fMonitor.SetName(Form("Monitor for %s", fName.Data()));
924 fMonitor.SetCounterName(EMonitorKey::kEvent, "N events");
925 fMonitor.SetCounterName(EMonitorKey::kHit, "N hits total");
926 fMonitor.SetCounterName(EMonitorKey::kHitAccepted, "N hits accepted");
927 fMonitor.SetCounterName(EMonitorKey::kMcPoint, "N MC points total");
928 fMonitor.SetCounterName(EMonitorKey::kMcPointWrongStation, "N MC points total");
929 fMonitor.SetRatioKeys({EMonitorKey::kEvent});
930
931
932 // ----- Histograms initialization
933 //
934 //SetStoringMode(EStoringMode::kSAMEDIR);
935 MakeQaDirectory("Summary");
936 MakeQaDirectory("Summary/vs Station");
937 if constexpr (ca::EDetectorID::kSts == DetID) {
938 MakeQaDirectory("Summary/vs N digi");
939 }
940 MakeQaDirectory("All stations");
941
942 fvph_hit_xy.resize(nSt + 1, nullptr);
943 fvph_hit_zy.resize(nSt + 1, nullptr);
944 fvph_hit_zx.resize(nSt + 1, nullptr);
945
946 fvph_hit_station_delta_z.resize(nSt + 1, nullptr);
947
948 fvph_hit_dx.resize(nSt + 1, nullptr);
949 fvph_hit_dy.resize(nSt + 1, nullptr);
950 fvph_hit_dt.resize(nSt + 1, nullptr);
951 fvph_hit_dv.resize(nSt + 1, nullptr);
952 fvph_hit_kuv.resize(nSt + 1, nullptr);
953 fvph_hit_du.resize(nSt + 1, nullptr);
954
955 std::string detName = fpDetInterface->GetDetectorName();
956
957 for (int iSt = 0; iSt <= nSt; ++iSt) {
958 TString sF = ""; // Subfolder
959 sF += (iSt == nSt) ? "All stations/" : Form("Station %d/", iSt);
960 TString nsuff = (iSt == nSt) ? "" : Form("_st%d", iSt); // Histogram name suffix
961 TString tsuff = (iSt == nSt) ? "" : Form(" in %s station %d", detName.c_str(), iSt); // Histogram title suffix
962 TString sN = "";
963 TString sT = "";
964
965 // place directories in a prefered order
966 MakeQaDirectory(sF + "occup/");
967 if (IsMCUsed()) {
968 MakeQaDirectory(sF + "res/");
969 MakeQaDirectory(sF + "pull/");
970 MakeQaDirectory(sF + "eff/");
971 }
972 MakeQaDirectory(sF + "err/");
973
974 int nBinsZ = (iSt < nSt) ? fNbins : fNbinsZ;
975
976 // Hit occupancy
977 sN = (TString) "hit_xy" + nsuff;
978 sT = (TString) "Hit occupancy in xy-Plane" + tsuff + ";x_{hit} [cm];y_{hit} [cm]";
979 fvph_hit_xy[iSt] =
980 MakeQaObject<TH2F>(sF + "occup/" + sN, sT, fNbins, frXmin[iSt], frXmax[iSt], fNbins, frYmin[iSt], frYmax[iSt]);
981
982 sN = (TString) "hit_zx" + nsuff;
983 sT = (TString) "Hit occupancy in xz-Plane" + tsuff + ";z_{hit} [cm];x_{hit} [cm]";
984 fvph_hit_zx[iSt] =
985 MakeQaObject<TH2F>(sF + "occup/" + sN, sT, nBinsZ, frZmin[iSt], frZmax[iSt], fNbins, frXmin[iSt], frXmax[iSt]);
986
987 sN = (TString) "hit_zy" + nsuff;
988 sT = (TString) "Hit occupancy in yz-plane" + tsuff + ";z_{hit} [cm];y_{hit} [cm]";
989 fvph_hit_zy[iSt] =
990 MakeQaObject<TH2F>(sF + "occup/" + sN, sT, nBinsZ, frZmin[iSt], frZmax[iSt], fNbins, frYmin[iSt], frYmax[iSt]);
991
992 // Hit errors
993 sN = (TString) "hit_dx" + nsuff;
994 sT = (TString) "Hit position error along x-axis" + tsuff + ";dx_{hit} [cm]";
995 fvph_hit_dx[iSt] = MakeQaObject<TH1F>(sF + "err/" + sN, sT, fNbins, fRHitDx[0], fRHitDx[1]);
996
997 sN = (TString) "hit_dy" + nsuff;
998 sT = (TString) "Hit position error along y-axis" + tsuff + ";dy_{hit} [cm]";
999 fvph_hit_dy[iSt] = MakeQaObject<TH1F>(sF + "err/" + sN, sT, fNbins, fRHitDy[0], fRHitDy[1]);
1000
1001 sN = (TString) "hit_du" + nsuff;
1002 sT = (TString) "Hit position error along the major detector coordinate U" + tsuff + ";du_{hit} [cm]";
1003 fvph_hit_du[iSt] = MakeQaObject<TH1F>(sF + "err/" + sN, sT, fNbins, fRHitDu[0], fRHitDu[1]);
1004
1005 sN = (TString) "hit_dv" + nsuff;
1006 sT = (TString) "Hit position error along the minor detector coordinate V" + tsuff + ";dv_{hit} [cm]";
1007 fvph_hit_dv[iSt] = MakeQaObject<TH1F>(sF + "err/" + sN, sT, fNbins, fRHitDv[0], fRHitDv[1]);
1008
1009 sN = (TString) "hit_kuv" + nsuff;
1010 sT = (TString) "Hit error correlation between the major (U) and the minor detector coordinate Vs" + tsuff
1011 + ";kuv_{hit} [unitless]";
1012 fvph_hit_kuv[iSt] = MakeQaObject<TH1F>(sF + "err/" + sN, sT, fNbins / 2 * 2 + 1, -1.1, 1.1);
1013
1014 sN = (TString) "hit_dt" + nsuff;
1015 sT = (TString) "Hit time error" + tsuff + ";dt_{hit} [ns]";
1016 fvph_hit_dt[iSt] = MakeQaObject<TH1F>(sF + "err/" + sN, sT, fNbins, fRHitDt[0], fRHitDt[1]);
1017
1018 sN = (TString) "hit_station_delta_z" + nsuff;
1019 sT = (TString) "Different between hit and station z-positions" + tsuff + ";z_{hit} - z_{st} [cm]";
1020 fvph_hit_station_delta_z[iSt] = MakeQaObject<TH1F>(sF + sN, sT, fNbins, -5., 5.);
1021
1022 } // loop over station index: end
1023
1024 // ----- Initialize histograms, which are use MC-information
1025 if (IsMCUsed()) {
1026 // Resize histogram vectors
1027 fvph_n_points_per_hit.resize(nSt + 1, nullptr);
1028 fvph_point_xy.resize(nSt + 1, nullptr);
1029 fvph_point_zx.resize(nSt + 1, nullptr);
1030 fvph_point_zy.resize(nSt + 1, nullptr);
1031 fvph_point_hit_delta_z.resize(nSt + 1, nullptr);
1032 fvph_res_x.resize(nSt + 1, nullptr);
1033 fvph_res_y.resize(nSt + 1, nullptr);
1034 fvph_res_u.resize(nSt + 1, nullptr);
1035 fvph_res_v.resize(nSt + 1, nullptr);
1036 fvph_res_t.resize(nSt + 1, nullptr);
1037 fvph_pull_x.resize(nSt + 1, nullptr);
1038 fvph_pull_y.resize(nSt + 1, nullptr);
1039 fvph_pull_u.resize(nSt + 1, nullptr);
1040 fvph_pull_v.resize(nSt + 1, nullptr);
1041 fvph_pull_t.resize(nSt + 1, nullptr);
1042 fvph_res_x_vs_x.resize(nSt + 1, nullptr);
1043 fvph_res_y_vs_y.resize(nSt + 1, nullptr);
1044 fvph_res_u_vs_u.resize(nSt + 1, nullptr);
1045 fvph_res_v_vs_v.resize(nSt + 1, nullptr);
1046 fvph_res_t_vs_t.resize(nSt + 1, nullptr);
1047 fvph_pull_x_vs_x.resize(nSt + 1, nullptr);
1048 fvph_pull_y_vs_y.resize(nSt + 1, nullptr);
1049 fvph_pull_u_vs_u.resize(nSt + 1, nullptr);
1050 fvph_pull_v_vs_v.resize(nSt + 1, nullptr);
1051 fvph_pull_t_vs_t.resize(nSt + 1, nullptr);
1052 fvpe_reco_eff_vs_xy.resize(nSt + 1, nullptr);
1053 fvph_reco_eff.resize(nSt + 1, nullptr);
1054
1055 for (int iSt = 0; iSt <= nSt; ++iSt) {
1056 TString sF = ""; // Subfolder
1057 sF += (iSt == nSt) ? "All stations/" : Form("Station %d/", iSt);
1058 TString nsuff = (iSt == nSt) ? "" : Form("_st%d", iSt); // Histogram name suffix
1059 TString tsuff = (iSt == nSt) ? "" : Form(" in %s station %d", detName.c_str(), iSt); // Histogram title suffix
1060 TString sN = "";
1061 TString sT = "";
1062
1063 sN = (TString) "n_points_per_hit" + nsuff;
1064 sT = (TString) "Number of points per hit" + tsuff + ";N_{point}/hit";
1065 fvph_n_points_per_hit[iSt] = MakeQaObject<TH1F>(sF + sN, sT, 10, -0.5, 9.5);
1066
1067 // Point occupancy
1068 sN = (TString) "point_xy" + nsuff;
1069 sT = (TString) "Point occupancy in XY plane" + tsuff + ";x_{MC} [cm];y_{MC} [cm]";
1070 fvph_point_xy[iSt] =
1071 MakeQaObject<TH2F>(sF + "occup/" + sN, sT, fNbins, frXmin[iSt], frXmax[iSt], fNbins, frYmin[iSt], frYmax[iSt]);
1072
1073 sN = (TString) "point_zx" + nsuff;
1074 sT = (TString) "Point Occupancy in XZ plane" + tsuff + ";z_{MC} [cm];x_{MC} [cm]";
1075 fvph_point_zx[iSt] =
1076 MakeQaObject<TH2F>(sF + "occup/" + sN, sT, fNbinsZ, frZmin[iSt], frZmax[iSt], fNbins, frXmin[iSt], frXmax[iSt]);
1077
1078 sN = (TString) "point_zy" + nsuff;
1079 sT = (TString) "Point Occupancy in YZ Plane" + tsuff + ";z_{MC} [cm];y_{MC} [cm]";
1080 fvph_point_zy[iSt] =
1081 MakeQaObject<TH2F>(sF + "occup/" + sN, sT, fNbinsZ, frZmin[iSt], frZmax[iSt], fNbins, frYmin[iSt], frYmax[iSt]);
1082
1083 // Difference between MC input z and hit z coordinates
1084 sN = (TString) "point_hit_delta_z" + nsuff;
1085 sT = (TString) "Distance between " + detName + " point and hit along z axis" + tsuff + ";z_{reco} - z_{MC} [cm]";
1086 fvph_point_hit_delta_z[iSt] = MakeQaObject<TH1F>(sF + sN, sT, fNbins, fRangeDzHitPoint[0], fRangeDzHitPoint[1]);
1087
1088 sN = (TString) "res_x" + nsuff;
1089 sT = (TString) "Residuals for X" + tsuff + ";x_{reco} - x_{MC} [cm]";
1090 fvph_res_x[iSt] = MakeQaObject<TH1F>(sF + "res/" + sN, sT, fNbins, fRResX[0], fRResX[1]);
1091
1092 sN = (TString) "res_y" + nsuff;
1093 sT = (TString) "Residuals for Y" + tsuff + ";y_{reco} - y_{MC} [cm]";
1094 fvph_res_y[iSt] = MakeQaObject<TH1F>(sF + "res/" + sN, sT, fNbins, fRResY[0], fRResY[1]);
1095
1096 sN = (TString) "res_u" + nsuff;
1097 sT = (TString) "Residuals for the major detector coordinate U" + tsuff + ";u_{reco} - u_{MC} [cm]";
1098 fvph_res_u[iSt] = MakeQaObject<TH1F>(sF + "res/" + sN, sT, fNbins, fRResU[0], fRResU[1]);
1099
1100 sN = (TString) "res_v" + nsuff;
1101 sT = (TString) "Residuals for the minor detector coordinate V" + tsuff + ";v_{reco} - v_{MC} [cm]";
1102 fvph_res_v[iSt] = MakeQaObject<TH1F>(sF + "res/" + sN, sT, fNbins, fRResV[0], fRResV[1]);
1103
1104 sN = (TString) "res_t" + nsuff;
1105 sT = (TString) "Residuals for Time" + tsuff + ";t_{reco} - t_{MC} [ns]";
1106 fvph_res_t[iSt] = MakeQaObject<TH1F>(sF + "res/" + sN, sT, fNbins, fRResT[0], fRResT[1]);
1107
1108 sN = (TString) "pull_x" + nsuff;
1109 sT = (TString) "Pulls for X" + tsuff + ";(x_{reco} - x_{MC}) / #sigma_{x}^{reco}";
1110 fvph_pull_x[iSt] = MakeQaObject<TH1F>(sF + "pull/" + sN, sT, kNbinsPull, kRPull[0], kRPull[1]);
1111
1112 sN = (TString) "pull_y" + nsuff;
1113 sT = (TString) "Pulls for Y" + tsuff + ";(y_{reco} - y_{MC}) / #sigma_{y}^{reco}";
1114 fvph_pull_y[iSt] = MakeQaObject<TH1F>(sF + "pull/" + sN, sT, kNbinsPull, kRPull[0], kRPull[1]);
1115
1116 sN = (TString) "pull_u" + nsuff;
1117 sT = (TString) "Pulls for the major detector coordinate U" + tsuff + ";(u_{reco} - u_{MC}) / #sigma_{u}^{reco}";
1118 fvph_pull_u[iSt] = MakeQaObject<TH1F>(sF + "pull/" + sN, sT, kNbinsPull, kRPull[0], kRPull[1]);
1119
1120 sN = (TString) "pull_v" + nsuff;
1121 sT = (TString) "Pulls for the minor detector coordinate V" + tsuff + ";(v_{reco} - v_{MC}) / #sigma_{v}^{reco}";
1122 fvph_pull_v[iSt] = MakeQaObject<TH1F>(sF + "pull/" + sN, sT, kNbinsPull, kRPull[0], kRPull[1]);
1123
1124 sN = (TString) "pull_t" + nsuff;
1125 sT = (TString) "Pulls for Time" + tsuff + ";(t_{reco} - t_{MC}) / #sigma_{t}^{reco}";
1126 fvph_pull_t[iSt] = MakeQaObject<TH1F>(sF + "pull/" + sN, sT, kNbinsPull, kRPull[0], kRPull[1]);
1127
1128 sN = (TString) "res_x_vs_x" + nsuff;
1129 sT = (TString) "Residuals for X" + tsuff + ";x_{MC} [cm];x_{reco} - x_{MC} [cm]";
1130 fvph_res_x_vs_x[iSt] =
1131 MakeQaObject<TH2F>(sF + "res/" + sN, sT, fNbins, frXmin[iSt], frXmax[iSt], fNbins, fRResX[0], fRResX[1]);
1132
1133 sN = (TString) "res_y_vs_y" + nsuff;
1134 sT = (TString) "Residuals for Y" + tsuff + ";y_{MC} [cm];y_{reco} - y_{MC} [cm]";
1135 fvph_res_y_vs_y[iSt] =
1136 MakeQaObject<TH2F>(sF + "res/" + sN, sT, fNbins, frYmin[iSt], frYmax[iSt], fNbins, fRResY[0], fRResY[1]);
1137
1138 sN = (TString) "res_u_vs_u" + nsuff;
1139 sT = (TString) "Residuals for the major detector coordinate U" + tsuff + ";u_{MC} [cm];u_{reco} - u_{MC} [cm]";
1140 fvph_res_u_vs_u[iSt] =
1141 MakeQaObject<TH2F>(sF + "res/" + sN, sT, fNbins, frXmin[iSt], frXmax[iSt], fNbins, fRResU[0], fRResU[1]);
1142
1143 sN = (TString) "res_v_vs_v" + nsuff;
1144 sT = (TString) "Residuals for the minor detector coordinate V" + tsuff + ";v_{MC} [cm];v_{reco} - v_{MC} [cm]";
1145 fvph_res_v_vs_v[iSt] =
1146 MakeQaObject<TH2F>(sF + "res/" + sN, sT, fNbins, frXmin[iSt], frXmax[iSt], fNbins, fRResV[0], fRResV[1]);
1147
1148 sN = (TString) "res_t_vs_t" + nsuff;
1149 sT = (TString) "Residuals for Time" + tsuff + ";t_{MC} [ns];t_{reco} - t_{MC} [ns]";
1150 fvph_res_t_vs_t[iSt] = MakeQaObject<TH2F>(sF + "res/" + sN, sT, fNbins, 0, 0, fNbins, fRResT[0], fRResT[1]);
1151
1152 sN = (TString) "pull_x_vs_x" + nsuff;
1153 sT = (TString) "Pulls for X" + tsuff + ";x_{MC} [cm];(x_{reco} - x_{MC}) / #sigma_{x}^{reco}";
1154 fvph_pull_x_vs_x[iSt] =
1155 MakeQaObject<TH2F>(sF + "pull/" + sN, sT, fNbins, frXmin[iSt], frXmax[iSt], kNbinsPull, kRPull[0], kRPull[1]);
1156
1157 sN = (TString) "pull_y_vs_y" + nsuff;
1158 sT = (TString) "Pulls for Y" + tsuff + ";y_{MC} [cm];(y_{reco} - y_{MC}) / #sigma_{y}^{reco}";
1159 fvph_pull_y_vs_y[iSt] =
1160 MakeQaObject<TH2F>(sF + "pull/" + sN, sT, fNbins, frYmin[iSt], frYmax[iSt], kNbinsPull, kRPull[0], kRPull[1]);
1161
1162 sN = (TString) "pull_u_vs_u" + nsuff;
1163 sT = (TString) "Pulls for the major detector coordinate U" + tsuff
1164 + ";u_{MC} [cm];(u_{reco} - u_{MC}) / #sigma_{u}^{reco}";
1165 fvph_pull_u_vs_u[iSt] =
1166 MakeQaObject<TH2F>(sF + "pull/" + sN, sT, fNbins, frXmin[iSt], frXmax[iSt], kNbinsPull, kRPull[0], kRPull[1]);
1167
1168 sN = (TString) "pull_v_vs_v" + nsuff;
1169 sT = (TString) "Pulls for the minor detector coordinate V" + tsuff
1170 + ";v_{MC} [cm];(v_{reco} - v_{MC}) / #sigma_{v}^{reco}";
1171 fvph_pull_v_vs_v[iSt] =
1172 MakeQaObject<TH2F>(sF + "pull/" + sN, sT, fNbins, frXmin[iSt], frXmax[iSt], kNbinsPull, kRPull[0], kRPull[1]);
1173
1174 sN = (TString) "pull_t_vs_t" + nsuff;
1175 sT = (TString) "Pulls for Time" + tsuff + ";t_{MC} [ns];(t_{reco} - t_{MC}) / #sigma_{t}^{reco}";
1176 fvph_pull_t_vs_t[iSt] = MakeQaObject<TH2F>(sF + "pull/" + sN, sT, fNbins, 0, 0, fNbins, kRPull[0], kRPull[1]);
1177
1178 sN = (TString) "reco_eff_vs_xy" + nsuff;
1179 sT = (TString) "Hit rec. efficiency in XY" + tsuff + ";x_{MC} [cm];y_{MC} [cm];#epsilon";
1180 fvpe_reco_eff_vs_xy[iSt] = MakeQaObject<TProfile2D>(sF + "eff/" + sN, sT, fNbins, frXmin[iSt], frXmax[iSt],
1181 fNbins, frYmin[iSt], frYmax[iSt]);
1182
1183 sN = (TString) "reco_eff" + nsuff;
1184 sT = (TString) "Hit rec. efficiency in XY bins" + tsuff + ";eff";
1185 fvph_reco_eff[iSt] = MakeQaObject<TH1F>(sF + "eff/" + sN, sT, 130, -0.005, 1.30 - 0.005);
1186 } // stations
1187 }
1188 return kSUCCESS;
1189}
1190
1191// ---------------------------------------------------------------------------------------------------------------------
1192//
1193template<ca::EDetectorID DetID>
1195{
1197 gStyle->SetOptFit(1);
1198 int nSt = fpDetInterface->GetNtrackingStations();
1199
1200 // ********************
1201 // ** Drawing options
1202
1203 // Contours
1204 constexpr auto contColor = kOrange + 7;
1205 constexpr auto contWidth = 2; // Line width
1206 constexpr auto contStyle = 2; // Line style
1207 constexpr auto contFill = 0; // Fill style
1208
1209 // **********************
1210 // ** summary for all stations
1211
1212 {
1213 auto* canv = MakeQaObject<TCanvas>("occ_hit", "Hit Occupancy", 1600, 800);
1214 canv->Divide(1, 3);
1215 canv->cd(1);
1216 fvph_hit_xy[nSt]->DrawCopy("colz", "");
1217 canv->cd(2);
1218 fvph_hit_zx[nSt]->DrawCopy("colz", "");
1219 canv->cd(3);
1220 fvph_hit_zy[nSt]->DrawCopy("colz", "");
1221 }
1222
1223 if (IsMCUsed()) {
1224 auto* canv = MakeQaObject<TCanvas>("occ_point", "Point Occupancy", 1600, 800);
1225 canv->Divide(1, 3);
1226 canv->cd(1);
1227 fvph_point_xy[nSt]->DrawCopy("colz", "");
1228 canv->cd(2);
1229 fvph_point_zx[nSt]->DrawCopy("colz", "");
1230 canv->cd(3);
1231 fvph_point_zy[nSt]->DrawCopy("colz", "");
1232 }
1233
1234 if (IsMCUsed()) {
1235 auto* canv = MakeQaObject<TCanvas>("residual", "Hit Residuals", 1600, 800);
1236 canv->Divide(1, 5);
1237 canv->cd(1);
1238 fvph_res_x[nSt]->DrawCopy("colz", "");
1239 canv->cd(2);
1240 fvph_res_y[nSt]->DrawCopy("colz", "");
1241 canv->cd(3);
1242 fvph_res_t[nSt]->DrawCopy("colz", "");
1243 canv->cd(4);
1244 fvph_res_u[nSt]->DrawCopy("colz", "");
1245 canv->cd(5);
1246 fvph_res_v[nSt]->DrawCopy("colz", "");
1247 }
1248
1249 if (IsMCUsed()) {
1250 auto* canv = MakeQaObject<TCanvas>("pull", "Hit Pulls", 1600, 800);
1251 canv->Divide(1, 5);
1252 canv->cd(1);
1253 fvph_pull_x[nSt]->DrawCopy("colz", "");
1254 canv->cd(2);
1255 fvph_pull_y[nSt]->DrawCopy("colz", "");
1256 canv->cd(3);
1257 fvph_pull_t[nSt]->DrawCopy("colz", "");
1258 canv->cd(4);
1259 fvph_pull_u[nSt]->DrawCopy("colz", "");
1260 canv->cd(5);
1261 fvph_pull_v[nSt]->DrawCopy("colz", "");
1262 }
1263
1264 if (IsMCUsed()) {
1265 auto* canv = MakeQaObject<TCanvas>("eff", "Hit Reconstruction Efficiency", 1600, 800);
1266 canv->Divide(1, 2);
1267 canv->cd(1);
1268 fvpe_reco_eff_vs_xy[nSt]->DrawCopy("colz", "");
1269 canv->cd(2);
1270 fvph_reco_eff[nSt]->DrawCopy("colz", "");
1271 }
1272
1273 {
1274 auto* canv = MakeQaObject<TCanvas>("err", "Hit Errors", 1600, 800);
1275 canv->Divide(2, 3);
1276 canv->cd(1);
1277 fvph_hit_dx[nSt]->DrawCopy();
1278 canv->cd(2);
1279 fvph_hit_dy[nSt]->DrawCopy();
1280 canv->cd(3);
1281 fvph_hit_dt[nSt]->DrawCopy();
1282 canv->cd(4);
1283 fvph_hit_du[nSt]->DrawCopy();
1284 canv->cd(5);
1285 fvph_hit_dv[nSt]->DrawCopy();
1286 canv->cd(6);
1287 fvph_hit_kuv[nSt]->DrawCopy();
1288 }
1289
1290 // TODO: Split the histograms for MC and real data
1291 if (IsMCUsed()) {
1292 auto* canv = MakeQaObject<TCanvas>("other", "Other histograms", 1600, 800);
1293 canv->Divide(1, 3);
1294 canv->cd(1);
1295 fvph_n_points_per_hit[nSt]->DrawCopy("colz", "");
1296 canv->cd(2);
1297 fvph_hit_station_delta_z[nSt]->DrawCopy("colz", "");
1298 canv->cd(3);
1299 fvph_point_hit_delta_z[nSt]->DrawCopy("colz", "");
1300 }
1301
1302
1303 // *********************************
1304 // ** Hit occupancies
1305
1306 { // ** Occupancy in XY plane
1307 auto* canv = MakeQaObject<TCanvas>("vs Station/hit_xy", "", 1600, 800);
1308 canv->DivideSquare(nSt);
1309 for (int iSt = 0; iSt < nSt; ++iSt) {
1310 canv->cd(iSt + 1);
1311 fvph_hit_xy[iSt]->DrawCopy("colz", "");
1312
1313 // Square boarders of the active area of the station
1314 double stXmin = fpDetInterface->GetXmin(iSt);
1315 double stXmax = fpDetInterface->GetXmax(iSt);
1316 double stYmin = fpDetInterface->GetYmin(iSt);
1317 double stYmax = fpDetInterface->GetYmax(iSt);
1318
1319 auto* pBox = new TBox(stXmin, stYmin, stXmax, stYmax);
1320 pBox->SetLineWidth(contWidth);
1321 pBox->SetLineStyle(contStyle);
1322 pBox->SetLineColor(contColor);
1323 pBox->SetFillStyle(contFill);
1324 pBox->Draw("SAME");
1325 }
1326 }
1327
1328 { // ** Occupancy in XZ plane
1329 auto* canv = MakeQaObject<TCanvas>("vs Station/hit_zx", "", 1600, 800);
1330 canv->cd();
1331 fvph_hit_zx[nSt]->DrawCopy("colz", "");
1332 for (int iSt = 0; iSt < nSt; ++iSt) {
1333 // Station positions in detector IFS
1334 double stZmin = fpDetInterface->GetZmin(iSt);
1335 double stZmax = fpDetInterface->GetZmax(iSt);
1336 double stHmin = fpDetInterface->GetXmin(iSt);
1337 double stHmax = fpDetInterface->GetXmax(iSt);
1338
1339 auto* pBox = new TBox(stZmin, stHmin, stZmax, stHmax);
1340 pBox->SetLineWidth(contWidth);
1341 pBox->SetLineStyle(contStyle);
1342 pBox->SetLineColor(contColor);
1343 pBox->SetFillStyle(contFill);
1344 pBox->Draw("SAME");
1345 }
1346 }
1347
1348 { // ** Occupancy in YZ plane
1349 auto* canv = MakeQaObject<TCanvas>("vs Station/hit_zy", "", 1600, 800);
1350 canv->cd();
1351 fvph_hit_zy[nSt]->DrawCopy("colz", "");
1352 for (int iSt = 0; iSt < nSt; ++iSt) {
1353 // Station positions in detector IFS
1354 double stZmin = fpDetInterface->GetZmin(iSt);
1355 double stZmax = fpDetInterface->GetZmax(iSt);
1356 double stHmin = fpDetInterface->GetYmin(iSt);
1357 double stHmax = fpDetInterface->GetYmax(iSt);
1358
1359 auto* pBox = new TBox(stZmin, stHmin, stZmax, stHmax);
1360 pBox->SetLineWidth(contWidth);
1361 pBox->SetLineStyle(contStyle);
1362 pBox->SetLineColor(contColor);
1363 pBox->SetFillStyle(contFill);
1364 pBox->Draw("SAME");
1365 }
1366 }
1367
1368 // ************
1369 // ************
1370
1371 if (IsMCUsed()) {
1372
1373 // **********************
1374 // ** Point occupancies
1375
1376 { // ** Occupancy in XY plane
1377 auto* canv = MakeQaObject<TCanvas>("vs Station/point_xy", "", 1600, 800);
1378 canv->DivideSquare(nSt);
1379 for (int iSt = 0; iSt < nSt; ++iSt) {
1380 canv->cd(iSt + 1);
1381 fvph_point_xy[iSt]->DrawCopy("colz", "");
1382
1383 // Square boarders of the active area of the station
1384 double stXmin = fpDetInterface->GetXmin(iSt);
1385 double stXmax = fpDetInterface->GetXmax(iSt);
1386 double stYmin = fpDetInterface->GetYmin(iSt);
1387 double stYmax = fpDetInterface->GetYmax(iSt);
1388
1389 auto* pBox = new TBox(stXmin, stYmin, stXmax, stYmax);
1390 pBox->SetLineWidth(contWidth);
1391 pBox->SetLineStyle(contStyle);
1392 pBox->SetLineColor(contColor);
1393 pBox->SetFillStyle(contFill);
1394 pBox->Draw("SAME");
1395 }
1396 }
1397
1398 { // ** Occupancy in XZ plane
1399 auto* canv = MakeQaObject<TCanvas>("vs Station/point_zx", "", 1600, 800);
1400 canv->cd();
1401 fvph_point_zx[nSt]->DrawCopy("colz", "");
1402 for (int iSt = 0; iSt < nSt; ++iSt) {
1403 // Station positions in detector IFS
1404 double stZmin = fpDetInterface->GetZmin(iSt);
1405 double stZmax = fpDetInterface->GetZmax(iSt);
1406 double stHmin = fpDetInterface->GetXmin(iSt);
1407 double stHmax = fpDetInterface->GetXmax(iSt);
1408
1409 auto* pBox = new TBox(stZmin, stHmin, stZmax, stHmax);
1410 pBox->SetLineWidth(contWidth);
1411 pBox->SetLineStyle(contStyle);
1412 pBox->SetLineColor(contColor);
1413 pBox->SetFillStyle(contFill);
1414 pBox->Draw("SAME");
1415 }
1416 }
1417
1418 { // ** Occupancy in YZ plane
1419 auto* canv = MakeQaObject<TCanvas>("vs Station/point_zy", "", 1600, 800);
1420 canv->cd();
1421 fvph_point_zy[nSt]->DrawCopy("colz", "");
1422 for (int iSt = 0; iSt < nSt; ++iSt) {
1423 // Station positions in detector IFS
1424 double stZmin = fpDetInterface->GetZmin(iSt);
1425 double stZmax = fpDetInterface->GetZmax(iSt);
1426 double stHmin = fpDetInterface->GetYmin(iSt);
1427 double stHmax = fpDetInterface->GetYmax(iSt);
1428
1429 auto* pBox = new TBox(stZmin, stHmin, stZmax, stHmax);
1430 pBox->SetLineWidth(contWidth);
1431 pBox->SetLineStyle(contStyle);
1432 pBox->SetLineColor(contColor);
1433 pBox->SetFillStyle(contFill);
1434 pBox->Draw("SAME");
1435 }
1436 }
1437
1438 // **************
1439 // ** Residuals
1440
1441 { // x-coordinate
1442 auto* canv = MakeQaObject<TCanvas>("vs Station/res_x", "Residuals for x coordinate", 1600, 800);
1443 canv->DivideSquare(nSt);
1444 for (int iSt = 0; iSt < nSt; ++iSt) {
1445 canv->cd(iSt + 1);
1446 fvph_res_x[iSt]->DrawCopy("", "");
1447 }
1448 }
1449
1450 { // y-coordinate
1451 auto* canv = MakeQaObject<TCanvas>("vs Station/res_y", "Residuals for y coordinate", 1600, 800);
1452 canv->DivideSquare(nSt);
1453 for (int iSt = 0; iSt < nSt; ++iSt) {
1454 canv->cd(iSt + 1);
1455 fvph_res_y[iSt]->DrawCopy("", "");
1456 }
1457 }
1458
1459 { // u-coordinate
1460 auto* canv = MakeQaObject<TCanvas>("vs Station/res_u", "Residuals for u coordinate", 1600, 800);
1461 canv->DivideSquare(nSt);
1462 for (int iSt = 0; iSt < nSt; ++iSt) {
1463 canv->cd(iSt + 1);
1464 fvph_res_u[iSt]->DrawCopy("", "");
1465 }
1466 }
1467
1468 { // v-coordinate
1469 auto* canv = MakeQaObject<TCanvas>("vs Station/res_v", "Residuals for v coordinate", 1600, 800);
1470 canv->DivideSquare(nSt);
1471 for (int iSt = 0; iSt < nSt; ++iSt) {
1472 canv->cd(iSt + 1);
1473 fvph_res_v[iSt]->DrawCopy("", "");
1474 }
1475 }
1476
1477 { // t-coordinate
1478 auto* canv = MakeQaObject<TCanvas>("vs Station/res_t", "Residuals for time", 1600, 800);
1479 canv->DivideSquare(nSt);
1480 for (int iSt = 0; iSt < nSt; ++iSt) {
1481 canv->cd(iSt + 1);
1482 fvph_res_t[iSt]->DrawCopy("", "");
1483 }
1484 }
1485
1486
1487 // **********
1488 // ** Pulls
1489
1490 { // x-coordinate
1491 auto* canv = MakeQaObject<TCanvas>("vs Station/pull_x", "Pulls for x coordinate", 1600, 800);
1492 canv->DivideSquare(nSt);
1493 for (int iSt = 0; iSt < nSt; ++iSt) {
1494 canv->cd(iSt + 1);
1495 fvph_pull_x[iSt]->DrawCopy("", "");
1496 }
1497 }
1498
1499 { // y-coordinate
1500 auto* canv = MakeQaObject<TCanvas>("vs Station/pull_y", "Pulls for y coordinate", 1600, 800);
1501 canv->DivideSquare(nSt);
1502 for (int iSt = 0; iSt < nSt; ++iSt) {
1503 canv->cd(iSt + 1);
1504 fvph_pull_y[iSt]->DrawCopy("", "");
1505 }
1506 }
1507
1508 { // u-coordinate
1509 auto* canv = MakeQaObject<TCanvas>("vs Station/pull_u", "Pulls for u coordinate", 1600, 800);
1510 canv->DivideSquare(nSt);
1511 for (int iSt = 0; iSt < nSt; ++iSt) {
1512 canv->cd(iSt + 1);
1513 fvph_pull_u[iSt]->DrawCopy("", "");
1514 }
1515 }
1516
1517 { // v-coordinate
1518 auto* canv = MakeQaObject<TCanvas>("vs Station/pull_v", "Pulls for v coordinate", 1600, 800);
1519 canv->DivideSquare(nSt);
1520 for (int iSt = 0; iSt < nSt; ++iSt) {
1521 canv->cd(iSt + 1);
1522 fvph_pull_v[iSt]->DrawCopy("", "");
1523 }
1524 }
1525
1526 { // t-coordinate
1527 auto* canv = MakeQaObject<TCanvas>("vs Station/pull_t", "Pulls for time", 1600, 800);
1528 canv->DivideSquare(nSt);
1529 for (int iSt = 0; iSt < nSt; ++iSt) {
1530 canv->cd(iSt + 1);
1531 fvph_pull_t[iSt]->DrawCopy("", "");
1532 }
1533 }
1534
1535 // ************************************
1536 // ** Hit reconstruction efficiencies
1537
1538 {
1539 auto* canv = MakeQaObject<TCanvas>("vs Station/reco_eff", "Hit efficiencies in xy bins", 1600, 800);
1540 canv->DivideSquare(nSt);
1541 for (int iSt = 0; iSt < nSt; ++iSt) {
1542 canv->cd(iSt + 1);
1543 //fvph_reco_eff[iSt]->SetMarkerStyle(20);
1544 fvph_reco_eff[iSt]->DrawCopy();
1545 }
1546 }
1547
1548 {
1549 auto* canv = MakeQaObject<TCanvas>("vs Station/reco_eff_vs_xy", "Hit efficiencies wrt x and y", 1600, 800);
1550 canv->DivideSquare(nSt);
1551 for (int iSt = 0; iSt < nSt; ++iSt) {
1552 canv->cd(iSt + 1);
1553 fvpe_reco_eff_vs_xy[iSt]->DrawCopy("colz", "");
1554 }
1555 }
1556 }
1557}
1558
1559// ---------------------------------------------------------------------------------------------------------------------
1560//
1561template<ca::EDetectorID DetID>
1563{
1564 assert(point);
1565
1566 if (fConfig.fMcTrackCuts.fbPrimary && track->GetMotherId() >= 0) {
1567 return false;
1568 }
1569
1570 double px = std::numeric_limits<double>::signaling_NaN();
1571 double py = std::numeric_limits<double>::signaling_NaN();
1572 double pz = std::numeric_limits<double>::signaling_NaN();
1573 if constexpr (ca::EDetectorID::kTof == DetID) {
1574 px = point->GetPx();
1575 py = point->GetPy();
1576 pz = point->GetPz();
1577 }
1578 else {
1579 px = point->GetPxOut();
1580 py = point->GetPyOut();
1581 pz = point->GetPzOut();
1582 }
1583 double p = sqrt(px * px + py * py + pz * pz);
1584
1585 if (pz <= 0.) {
1586 return false;
1587 }
1588
1589 if (p < fConfig.fMcTrackCuts.fMinMom) {
1590 return false;
1591 }
1592
1593 if (TMath::ATan2(sqrt(px * px + py * py), pz) * TMath::RadToDeg() > fConfig.fMcTrackCuts.fMaxTheta) {
1594 return false;
1595 }
1596
1597 return true;
1598}
1599
Compile-time constants definition for the CA tracking algorithm.
Class providing basic functionality for CA input QA-tasks.
Implementation of L1DetectorID enum class for CBM.
Class for pixel hits in MUCH detector.
Definition of CbmQaTable class.
Useful utilities for CBM QA tasks.
Data class for STS clusters.
Data class for a reconstructed hit in the STS.
Base abstract class for tracking detector interface to L1 (implementation of Checker)
Class for hits in TRD detector.
friend fvec sqrt(const fvec &a)
friend fscal max(fscal x, fscal y)
friend fscal min(fscal x, fscal y)
Generates beam ions for transport simulation.
A QA-task class, which provides assurance of MuCh hits and geometry.
void CreateSummary() override
Creates summary cavases, tables etc.
void ExecQa() override
Fills histograms per event or time-slice.
std::pair< std::string, bool > CheckRangePull(TH1 *h) const
Checks ranges for mean and standard deviation.
cbm::ca::PointTypes_t::at< DetID > Point_t
Point type for a given detector ID.
void Check() override
Method to check, if the QA results are acceptable.
CbmCaInputQaBase(const char *name, int verbose, bool isMCUsed)
Constructor from parameters.
bool IsTrackSelected(const CbmMCTrack *track, const Point_t *point) const
checks if the track at mc point passes the track selection cuts
InitStatus InitQa() override
Initializes QA task.
void DeInit() override
De-initializes histograms.
Task class creating and managing CbmMCDataArray objects.
Container class for MC events with number, file and start time.
int32_t GetMotherId() const
Definition CbmMCTrack.h:69
double GetMass() const
Mass of the associated particle.
void SetStoringMode(EStoringMode mode)
Set storing mode.
Definition CbmQaIO.h:112
@ kSUBDIR
Objects of different type will be stored in different subdirectories like histograms/ canvases/.
TODO: SZh, 30.01.2023: Override THistPainter::PaintText() to add zeroes in tables.
Definition CbmQaTable.h:24
Data class with information on a STS local track.
Bookkeeping of time-slice content.
static int32_t GetSmType(uint32_t address)
constexpr fscal SpeedOfLight
Speed of light [cm/ns].
Definition CaDefs.h:88
TODO: SZh 8.11.2022: add selection of parameterisation.
Definition CaBranch.h:14
constexpr DetIdArr_t< const char * > kDetPointBrName
Name of point branches for each detector.
constexpr DetIdArr_t< const char * > kDetHitBrName
Name of hit branches for each detector.
void SetLargeStats(TH1 *pHist)
Set large stat. window.
Specific configuration for the CA input QA task.