Project
Loading...
Searching...
No Matches
GPUQA.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
14
15#define QA_DEBUG 0
16#define QA_TIMING 0
17
18#include "Rtypes.h" // Include ROOT header first, to use ROOT and disable replacements
19
20#include "TH1F.h"
21#include "TH2F.h"
22#include "TH1D.h"
23#include "TGraphAsymmErrors.h"
24#include "TCanvas.h"
25#include "TPad.h"
26#include "TLegend.h"
27#include "TColor.h"
28#include "TPaveText.h"
29#include "TF1.h"
30#include "TFile.h"
31#include "TTree.h"
32#include "TStyle.h"
33#include "TLatex.h"
34#include "TObjArray.h"
35#include <sys/stat.h>
36
37#include "GPUQA.h"
38#include "GPUTPCDef.h"
39#include "GPUTPCTrackingData.h"
40#include "GPUChainTracking.h"
41#include "GPUChainTrackingGetters.inc"
42#include "GPUTPCTrack.h"
43#include "GPUTPCTracker.h"
44#include "GPUTPCGMMergedTrack.h"
45#include "GPUTPCGMPropagator.h"
47#include "GPUTPCMCInfo.h"
48#include "GPUTPCClusterData.h"
49#include "GPUO2DataTypes.h"
50#include "GPUParam.inc"
52#include "GPUTPCConvertImpl.h"
53#include "TPCFastTransform.h"
55#include "GPUROOTDump.h"
58#include "GPUSettings.h"
59#ifdef GPUCA_O2_LIB
69#include "TPDGCode.h"
70#include "TParticlePDG.h"
71#include "TDatabasePDG.h"
72#endif
73#include "GPUQAHelper.h"
74#include <algorithm>
75#include <cstdio>
76#include <cinttypes>
77
78#include "utils/timer.h"
79
80#include <oneapi/tbb.h>
81
82using namespace o2::gpu;
83
84#ifdef GPUCA_MERGER_BY_MC_LABEL
85#define CHECK_CLUSTER_STATE_INIT_LEG_BY_MC() \
86 if (!unattached && mTrackMCLabels[id].isValid()) { \
87 int32_t mcLabel = mTrackMCLabels[id].getTrackID(); \
88 int32_t mcEvent = mTrackMCLabels[id].getEventID(); \
89 int32_t mcSource = mTrackMCLabels[id].getSourceID(); \
90 if (mTrackMCLabelsReverse[mMCEventOffset[mcSource] + mcEvent][mcLabel] != id) { \
91 attach &= (~gputpcgmmergertypes::attachGoodLeg); \
92 } \
93 }
94#else
95#define CHECK_CLUSTER_STATE_INIT_LEG_BY_MC()
96#endif
97
98#define CHECK_CLUSTER_STATE_INIT() \
99 bool unattached = attach == 0; \
100 float qpt = 0; \
101 bool lowPt = false; \
102 bool mev200 = false; \
103 bool mergedLooper = false; \
104 int32_t id = attach & gputpcgmmergertypes::attachTrackMask; \
105 if (!unattached) { \
106 qpt = fabsf(mTracking->mIOPtrs.mergedTracks[id].GetParam().GetQPt()); \
107 lowPt = qpt * mTracking->GetParam().qptB5Scaler > mTracking->GetParam().rec.tpc.rejectQPtB5; \
108 mev200 = qpt > 5; \
109 mergedLooper = mTracking->mIOPtrs.mergedTracks[id].MergedLooper(); \
110 } \
111 bool physics = false, protect = false; \
112 CHECK_CLUSTER_STATE_INIT_LEG_BY_MC();
113
114#define CHECK_CLUSTER_STATE() \
115 CHECK_CLUSTER_STATE_INIT() \
116 if (mev200) { \
117 mClusterCounts.n200MeV++; \
118 } \
119 if (lowPt) { \
120 mClusterCounts.nLowPt++; \
121 } else if (mergedLooper) { \
122 mClusterCounts.nMergedLooper++; \
123 } else { \
124 GPUTPCClusterRejection::GetProtectionStatus<true>(attach, physics, protect, &mClusterCounts, &mev200); \
125 }
126
127#define CHECK_CLUSTER_STATE_NOCOUNT() \
128 CHECK_CLUSTER_STATE_INIT() \
129 (void)mev200; /* silence unused variable warning*/ \
130 if (!lowPt && !mergedLooper) { \
131 GPUTPCClusterRejection::GetProtectionStatus<false>(attach, physics, protect); \
132 }
133
134static const GPUSettingsQA& GPUQA_GetConfig(GPUChainTracking* chain)
135{
136 static GPUSettingsQA defaultConfig;
137 if (chain && chain->mConfigQA) {
138 return *chain->mConfigQA;
139 } else {
140 return defaultConfig;
141 }
142}
143
144// static const constexpr bool PLOT_ROOT = 0;
145// static const constexpr bool FIX_SCALES = 0;
146static const constexpr bool PERF_FIGURE = 0;
147// static const constexpr float FIXED_SCALES_MIN[5] = {-0.05, -0.05, -0.2, -0.2, -0.5};
148// static const constexpr float FIXED_SCALES_MAX[5] = {0.4, 0.7, 5, 3, 6.5};
149static const constexpr float LOG_PT_MIN = -1.;
150
151static constexpr float Y_MAX = 40;
152static constexpr float Z_MAX = 100;
153static constexpr float PT_MIN = GPUCA_MIN_TRACK_PTB5_DEFAULT;
154// static constexpr float PT_MIN2 = 0.1;
155static constexpr float PT_MIN_PRIM = 0.1;
156static constexpr float PT_MIN_CLUST = GPUCA_MIN_TRACK_PTB5_DEFAULT;
157static constexpr float PT_MAX = 20;
158static constexpr float ETA_MAX = 1.5;
159static constexpr float ETA_MAX2 = 0.9;
160
161static constexpr bool CLUST_HIST_INT_SUM = false;
162
163static constexpr const int32_t COLORCOUNT = 12;
164
165static const constexpr char* EFF_TYPES[5] = {"Rec", "Clone", "Fake", "All", "RecAndClone"};
166static const constexpr char* FINDABLE_NAMES[2] = {"All", "Findable"};
167static const constexpr char* PRIM_NAMES[2] = {"Prim", "Sec"};
168static const constexpr char* PARAMETER_NAMES[5] = {"Y", "Z", "#Phi", "#lambda", "Relative #it{p}_{T}"};
169static const constexpr char* PARAMETER_NAMES_NATIVE[5] = {"Y", "Z", "sin(#Phi)", "tan(#lambda)", "q/#it{p}_{T} (curvature)"};
170static const constexpr char* VSPARAMETER_NAMES[6] = {"Y", "Z", "Phi", "Eta", "Pt", "Pt_log"};
171static const constexpr char* EFF_NAMES[3] = {"Efficiency", "Clone Rate", "Fake Rate"};
172static const constexpr char* EFFICIENCY_TITLES[4] = {"Efficiency (Primary Tracks, Findable)", "Efficiency (Secondary Tracks, Findable)", "Efficiency (Primary Tracks)", "Efficiency (Secondary Tracks)"};
173static const constexpr double SCALE[5] = {10., 10., 1000., 1000., 100.};
174static const constexpr double SCALE_NATIVE[5] = {10., 10., 1000., 1000., 1.};
175static const constexpr char* XAXIS_TITLES[5] = {"#it{y}_{mc} (cm)", "#it{z}_{mc} (cm)", "#Phi_{mc} (rad)", "#eta_{mc}", "#it{p}_{Tmc} (GeV/#it{c})"};
176static const constexpr char* AXIS_TITLES[5] = {"#it{y}-#it{y}_{mc} (mm) (Resolution)", "#it{z}-#it{z}_{mc} (mm) (Resolution)", "#phi-#phi_{mc} (mrad) (Resolution)", "#lambda-#lambda_{mc} (mrad) (Resolution)", "(#it{p}_{T} - #it{p}_{Tmc}) / #it{p}_{Tmc} (%) (Resolution)"};
177static const constexpr char* AXIS_TITLES_NATIVE[5] = {"#it{y}-#it{y}_{mc} (mm) (Resolution)", "#it{z}-#it{z}_{mc} (mm) (Resolution)", "sin(#phi)-sin(#phi_{mc}) (Resolution)", "tan(#lambda)-tan(#lambda_{mc}) (Resolution)", "q*(q/#it{p}_{T} - q/#it{p}_{Tmc}) (Resolution)"};
178static const constexpr char* AXIS_TITLES_PULL[5] = {"#it{y}-#it{y}_{mc}/#sigma_{y} (Pull)", "#it{z}-#it{z}_{mc}/#sigma_{z} (Pull)", "sin(#phi)-sin(#phi_{mc})/#sigma_{sin(#phi)} (Pull)", "tan(#lambda)-tan(#lambda_{mc})/#sigma_{tan(#lambda)} (Pull)",
179 "q*(q/#it{p}_{T} - q/#it{p}_{Tmc})/#sigma_{q/#it{p}_{T}} (Pull)"};
180static const constexpr char* CLUSTER_NAMES[GPUQA::N_CLS_HIST] = {"Correctly attached clusters", "Fake attached clusters", "Attached + adjacent clusters", "Fake adjacent clusters", "Clusters of reconstructed tracks", "Used in Physics", "Protected", "All clusters"};
181static const constexpr char* CLUSTER_TITLES[GPUQA::N_CLS_TYPE] = {"Clusters Pt Distribution / Attachment", "Clusters Pt Distribution / Attachment (relative to all clusters)", "Clusters Pt Distribution / Attachment (integrated)"};
182static const constexpr char* CLUSTER_NAMES_SHORT[GPUQA::N_CLS_HIST] = {"Attached", "Fake", "AttachAdjacent", "FakeAdjacent", "FoundTracks", "Physics", "Protected", "All"};
183static const constexpr char* CLUSTER_TYPES[GPUQA::N_CLS_TYPE] = {"", "Ratio", "Integral"};
184static const constexpr int32_t COLORS_HEX[COLORCOUNT] = {0xB03030, 0x00A000, 0x0000C0, 0x9400D3, 0x19BBBF, 0xF25900, 0x7F7F7F, 0xFFD700, 0x07F707, 0x07F7F7, 0xF08080, 0x000000};
185
186static const constexpr int32_t CONFIG_DASHED_MARKERS = 0;
187
188static const constexpr float AXES_MIN[5] = {-Y_MAX, -Z_MAX, 0.f, -ETA_MAX, PT_MIN};
189static const constexpr float AXES_MAX[5] = {Y_MAX, Z_MAX, 2.f * M_PI, ETA_MAX, PT_MAX};
190static const constexpr int32_t AXIS_BINS[5] = {51, 51, 144, 31, 50};
191static const constexpr int32_t RES_AXIS_BINS[] = {1017, 113}; // Consecutive bin sizes, histograms are binned down until the maximum entry is 50, each bin size should evenly divide its predecessor.
192static const constexpr float RES_AXES[5] = {1., 1., 0.03, 0.03, 1.0};
193static const constexpr float RES_AXES_NATIVE[5] = {1., 1., 0.1, 0.1, 5.0};
194static const constexpr float PULL_AXIS = 10.f;
195
196std::vector<TColor*> GPUQA::mColors;
197int32_t GPUQA::initColors()
198{
199 mColors.reserve(COLORCOUNT);
200 for (int32_t i = 0; i < COLORCOUNT; i++) {
201 float f1 = (float)((COLORS_HEX[i] >> 16) & 0xFF) / (float)0xFF;
202 float f2 = (float)((COLORS_HEX[i] >> 8) & 0xFF) / (float)0xFF;
203 float f3 = (float)((COLORS_HEX[i] >> 0) & 0xFF) / (float)0xFF;
204 mColors.emplace_back(new TColor(10000 + i, f1, f2, f3));
205 }
206 return 0;
207}
208static constexpr Color_t defaultColorNums[COLORCOUNT] = {kRed, kBlue, kGreen, kMagenta, kOrange, kAzure, kBlack, kYellow, kGray, kTeal, kSpring, kPink};
209
210#define TRACK_EXPECTED_REFERENCE_X_DEFAULT 81
211#ifdef GPUCA_TPC_GEOMETRY_O2
212static inline int32_t GPUQA_O2_ConvertFakeLabel(int32_t label) { return label >= 0x7FFFFFFE ? -1 : label; }
213inline uint32_t GPUQA::GetNMCCollissions() const { return mMCInfosCol.size(); }
214inline uint32_t GPUQA::GetNMCTracks(int32_t iCol) const { return mMCInfosCol[iCol].num; }
215inline uint32_t GPUQA::GetNMCTracks(const mcLabelI_t& label) const { return mMCInfosCol[mMCEventOffset[label.getSourceID()] + label.getEventID()].num; }
216inline uint32_t GPUQA::GetNMCLabels() const { return mClNative->clustersMCTruth ? mClNative->clustersMCTruth->getIndexedSize() : 0; }
217inline const GPUQA::mcInfo_t& GPUQA::GetMCTrack(uint32_t iTrk, uint32_t iCol) { return mMCInfos[mMCInfosCol[iCol].first + iTrk]; }
218inline const GPUQA::mcInfo_t& GPUQA::GetMCTrack(const mcLabel_t& label) { return mMCInfos[mMCInfosCol[mMCEventOffset[label.getSourceID()] + label.getEventID()].first + label.getTrackID()]; }
219inline GPUQA::mcLabels_t GPUQA::GetMCLabel(uint32_t i) { return mClNative->clustersMCTruth->getLabels(i); }
220inline int32_t GPUQA::GetMCLabelNID(const mcLabels_t& label) { return label.size(); }
221inline int32_t GPUQA::GetMCLabelNID(uint32_t i) { return mClNative->clustersMCTruth->getLabels(i).size(); }
222inline GPUQA::mcLabel_t GPUQA::GetMCLabel(uint32_t i, uint32_t j) { return mClNative->clustersMCTruth->getLabels(i)[j]; }
223inline int32_t GPUQA::GetMCLabelID(uint32_t i, uint32_t j) { return GPUQA_O2_ConvertFakeLabel(mClNative->clustersMCTruth->getLabels(i)[j].getTrackID()); }
224inline int32_t GPUQA::GetMCLabelID(const mcLabels_t& label, uint32_t j) { return GPUQA_O2_ConvertFakeLabel(label[j].getTrackID()); }
225inline int32_t GPUQA::GetMCLabelID(const mcLabel_t& label) { return GPUQA_O2_ConvertFakeLabel(label.getTrackID()); }
226inline uint32_t GPUQA::GetMCLabelCol(uint32_t i, uint32_t j) { return mMCEventOffset[mClNative->clustersMCTruth->getLabels(i)[j].getSourceID()] + mClNative->clustersMCTruth->getLabels(i)[j].getEventID(); }
227inline const auto& GPUQA::GetClusterLabels() { return mClNative->clustersMCTruth; }
228inline float GPUQA::GetMCLabelWeight(uint32_t i, uint32_t j) { return 1; }
229inline float GPUQA::GetMCLabelWeight(const mcLabels_t& label, uint32_t j) { return 1; }
230inline float GPUQA::GetMCLabelWeight(const mcLabel_t& label) { return 1; }
231inline bool GPUQA::mcPresent() { return !mConfig.noMC && mTracking && mClNative && mClNative->clustersMCTruth && mMCInfos.size(); }
232uint32_t GPUQA::GetMCLabelCol(const mcLabel_t& label) const { return !label.isValid() ? 0 : (mMCEventOffset[label.getSourceID()] + label.getEventID()); }
233GPUQA::mcLabelI_t GPUQA::GetMCTrackLabel(uint32_t trackId) const { return trackId >= mTrackMCLabels.size() ? MCCompLabel() : mTrackMCLabels[trackId]; }
234bool GPUQA::CompareIgnoreFake(const mcLabelI_t& l1, const mcLabelI_t& l2) { return l1.compare(l2) >= 0; }
235#define TRACK_EXPECTED_REFERENCE_X 78
236#else
237inline GPUQA::mcLabelI_t::mcLabelI_t(const GPUQA::mcLabel_t& l) : track(l.fMCID) {}
238inline bool GPUQA::mcLabelI_t::operator==(const GPUQA::mcLabel_t& l) { return AbsLabelID(track) == l.fMCID; }
239inline uint32_t GPUQA::GetNMCCollissions() const { return 1; }
240inline uint32_t GPUQA::GetNMCTracks(int32_t iCol) const { return mTracking->mIOPtrs.nMCInfosTPC; }
241inline uint32_t GPUQA::GetNMCTracks(const mcLabelI_t& label) const { return mTracking->mIOPtrs.nMCInfosTPC; }
242inline uint32_t GPUQA::GetNMCLabels() const { return mTracking->mIOPtrs.nMCLabelsTPC; }
243inline const GPUQA::mcInfo_t& GPUQA::GetMCTrack(uint32_t iTrk, uint32_t iCol) { return mTracking->mIOPtrs.mcInfosTPC[AbsLabelID(iTrk)]; }
244inline const GPUQA::mcInfo_t& GPUQA::GetMCTrack(const mcLabel_t& label) { return GetMCTrack(label.fMCID, 0); }
245inline const GPUQA::mcInfo_t& GPUQA::GetMCTrack(const mcLabelI_t& label) { return GetMCTrack(label.track, 0); }
246inline const GPUQA::mcLabels_t& GPUQA::GetMCLabel(uint32_t i) { return mTracking->mIOPtrs.mcLabelsTPC[i]; }
247inline const GPUQA::mcLabel_t& GPUQA::GetMCLabel(uint32_t i, uint32_t j) { return mTracking->mIOPtrs.mcLabelsTPC[i].fClusterID[j]; }
248inline int32_t GPUQA::GetMCLabelNID(const mcLabels_t& label) { return 3; }
249inline int32_t GPUQA::GetMCLabelNID(uint32_t i) { return 3; }
250inline int32_t GPUQA::GetMCLabelID(uint32_t i, uint32_t j) { return mTracking->mIOPtrs.mcLabelsTPC[i].fClusterID[j].fMCID; }
251inline int32_t GPUQA::GetMCLabelID(const mcLabels_t& label, uint32_t j) { return label.fClusterID[j].fMCID; }
252inline int32_t GPUQA::GetMCLabelID(const mcLabel_t& label) { return label.fMCID; }
253inline uint32_t GPUQA::GetMCLabelCol(uint32_t i, uint32_t j) { return 0; }
254
255inline const auto& GPUQA::GetClusterLabels() { return mTracking->mIOPtrs.mcLabelsTPC; }
256inline float GPUQA::GetMCLabelWeight(uint32_t i, uint32_t j) { return mTracking->mIOPtrs.mcLabelsTPC[i].fClusterID[j].fWeight; }
257inline float GPUQA::GetMCLabelWeight(const mcLabels_t& label, uint32_t j) { return label.fClusterID[j].fWeight; }
258inline float GPUQA::GetMCLabelWeight(const mcLabel_t& label) { return label.fWeight; }
259inline int32_t GPUQA::FakeLabelID(int32_t id) { return id < 0 ? id : (-2 - id); }
260inline int32_t GPUQA::AbsLabelID(int32_t id) { return id >= 0 ? id : (-id - 2); }
261inline bool GPUQA::mcPresent() { return !mConfig.noMC && mTracking && GetNMCLabels() && GetNMCTracks(0); }
262uint32_t GPUQA::GetMCLabelCol(const mcLabel_t& label) const { return 0; }
263GPUQA::mcLabelI_t GPUQA::GetMCTrackLabel(uint32_t trackId) const { return trackId >= mTrackMCLabels.size() ? mcLabelI_t() : mTrackMCLabels[trackId]; }
264bool GPUQA::CompareIgnoreFake(const mcLabelI_t& l1, const mcLabelI_t& l2) { return AbsLabelID(l1) == AbsLabelID(l2); }
265#define TRACK_EXPECTED_REFERENCE_X TRACK_EXPECTED_REFERENCE_X_DEFAULT
266#endif
267template <class T>
268inline auto& GPUQA::GetMCTrackObj(T& obj, const GPUQA::mcLabelI_t& l)
269{
270 return obj[mMCEventOffset[l.getSourceID()] + l.getEventID()][l.getTrackID()];
271}
272
273template <>
274auto GPUQA::getHistArray<TH1F>()
275{
276 return std::make_pair(mHist1D, &mHist1D_pos);
277}
278template <>
279auto GPUQA::getHistArray<TH2F>()
280{
281 return std::make_pair(mHist2D, &mHist2D_pos);
282}
283template <>
284auto GPUQA::getHistArray<TH1D>()
285{
286 return std::make_pair(mHist1Dd, &mHist1Dd_pos);
287}
288template <>
289auto GPUQA::getHistArray<TGraphAsymmErrors>()
290{
291 return std::make_pair(mHistGraph, &mHistGraph_pos);
292}
293template <class T, typename... Args>
294void GPUQA::createHist(T*& h, const char* name, Args... args)
295{
296 const auto& p = getHistArray<T>();
297 if (mHaveExternalHists) {
298 if (p.first->size() <= p.second->size()) {
299 GPUError("Array sizes mismatch: Histograms %lu <= Positions %lu", p.first->size(), p.second->size());
300 throw std::runtime_error("Incoming histogram array incomplete");
301 }
302 if (strcmp((*p.first)[p.second->size()].GetName(), name)) {
303 GPUError("Histogram name mismatch: in array %s, trying to create %s", (*p.first)[p.second->size()].GetName(), name);
304 throw std::runtime_error("Incoming histogram has incorrect name");
305 }
306 } else {
307 if constexpr (std::is_same_v<T, TGraphAsymmErrors>) {
308 p.first->emplace_back();
309 p.first->back().SetName(name);
310 } else {
311 p.first->emplace_back(name, args...);
312 }
313 }
314 h = &((*p.first)[p.second->size()]);
315 p.second->emplace_back(&h);
316}
317
318namespace o2::gpu::internal
319{
321 std::tuple<std::vector<std::unique_ptr<TCanvas>>, std::vector<std::unique_ptr<TLegend>>, std::vector<std::unique_ptr<TPad>>, std::vector<std::unique_ptr<TLatex>>, std::vector<std::unique_ptr<TH1D>>> v;
322};
323} // namespace o2::gpu::internal
324
325template <class T, typename... Args>
326T* GPUQA::createGarbageCollected(Args... args)
327{
328 auto& v = std::get<std::vector<std::unique_ptr<T>>>(mGarbageCollector->v);
329 v.emplace_back(std::make_unique<T>(args...));
330 return v.back().get();
331}
332void GPUQA::clearGarbagageCollector()
333{
334 std::get<std::vector<std::unique_ptr<TPad>>>(mGarbageCollector->v).clear(); // Make sure to delete TPad first due to ROOT ownership (std::tuple has no defined order in its destructor)
335 std::apply([](auto&&... args) { ((args.clear()), ...); }, mGarbageCollector->v);
336}
337
338GPUQA::GPUQA(GPUChainTracking* chain, const GPUSettingsQA* config, const GPUParam* param) : mTracking(chain), mConfig(config ? *config : GPUQA_GetConfig(chain)), mParam(param ? param : &chain->GetParam()), mGarbageCollector(std::make_unique<internal::GPUQAGarbageCollection>())
339{
340 mMCEventOffset.resize(1, 0);
341}
342
344{
345 if (mQAInitialized && !mHaveExternalHists) {
346 delete mHist1D;
347 delete mHist2D;
348 delete mHist1Dd;
349 delete mHistGraph;
350 }
351 clearGarbagageCollector(); // Needed to guarantee correct order for ROOT ownership
352}
353
354bool GPUQA::clusterRemovable(int32_t attach, bool prot) const
355{
357 if (prot) {
358 return protect || physics;
359 }
360 return (!unattached && !physics && !protect);
361}
362
363template <class T>
364void GPUQA::SetAxisSize(T* e)
365{
366 e->GetYaxis()->SetTitleOffset(1.0);
367 e->GetYaxis()->SetTitleSize(0.045);
368 e->GetYaxis()->SetLabelSize(0.045);
369 e->GetXaxis()->SetTitleOffset(1.03);
370 e->GetXaxis()->SetTitleSize(0.045);
371 e->GetXaxis()->SetLabelOffset(-0.005);
372 e->GetXaxis()->SetLabelSize(0.045);
373}
374
375void GPUQA::SetLegend(TLegend* l)
376{
377 l->SetTextFont(72);
378 l->SetTextSize(0.016);
379 l->SetFillColor(0);
380}
381
382double* GPUQA::CreateLogAxis(int32_t nbins, float xmin, float xmax)
383{
384 float logxmin = std::log10(xmin);
385 float logxmax = std::log10(xmax);
386 float binwidth = (logxmax - logxmin) / nbins;
387
388 double* xbins = new double[nbins + 1];
389
390 xbins[0] = xmin;
391 for (int32_t i = 1; i <= nbins; i++) {
392 xbins[i] = std::pow(10, logxmin + i * binwidth);
393 }
394 return xbins;
395}
396
397void GPUQA::ChangePadTitleSize(TPad* p, float size)
398{
399 p->Update();
400 TPaveText* pt = (TPaveText*)(p->GetPrimitive("title"));
401 if (pt == nullptr) {
402 GPUError("Error changing title");
403 } else {
404 pt->SetTextSize(size);
405 p->Modified();
406 }
407}
408
409void GPUQA::DrawHisto(TH1* histo, char* filename, char* options)
410{
411 TCanvas tmp;
412 tmp.cd();
413 histo->Draw(options);
414 tmp.Print(filename);
415}
416
417void GPUQA::doPerfFigure(float x, float y, float size)
418{
419 if (!PERF_FIGURE) {
420 return;
421 }
422 TLatex* t = createGarbageCollected<TLatex>();
423 t->SetNDC(kTRUE);
424 t->SetTextColor(1);
425 t->SetTextSize(size);
426 t->DrawLatex(x, y, str_perf_figure_1);
427 t->DrawLatex(x, y - 0.01 - size, str_perf_figure_2);
428}
429
430void GPUQA::SetMCTrackRange(int32_t min, int32_t max)
431{
432 mMCTrackMin = min;
433 mMCTrackMax = max;
434}
435
436int32_t GPUQA::InitQACreateHistograms()
437{
438 char name[2048], fname[1024];
439 if (mQATasks & taskTrackingEff) {
440 // Create Efficiency Histograms
441 for (int32_t i = 0; i < 5; i++) {
442 for (int32_t j = 0; j < 2; j++) {
443 for (int32_t k = 0; k < 2; k++) {
444 for (int32_t l = 0; l < 5; l++) {
445 snprintf(name, 2048, "%s%s%s%sVs%s", "tracks", EFF_TYPES[i], FINDABLE_NAMES[j], PRIM_NAMES[k], VSPARAMETER_NAMES[l]);
446 if (l == 4) {
447 std::unique_ptr<double[]> binsPt{CreateLogAxis(AXIS_BINS[4], k == 0 ? PT_MIN_PRIM : AXES_MIN[4], AXES_MAX[4])};
448 createHist(mEff[i][j][k][l], name, name, AXIS_BINS[l], binsPt.get());
449 } else {
450 createHist(mEff[i][j][k][l], name, name, AXIS_BINS[l], AXES_MIN[l], AXES_MAX[l]);
451 }
452 if (!mHaveExternalHists) {
453 mEff[i][j][k][l]->Sumw2();
454 }
455 strcat(name, "_eff");
456 if (i < 4) {
457 createHist(mEffResult[i][j][k][l], name);
458 }
459 }
460 }
461 }
462 }
463 }
464
465 // Create Resolution Histograms
466 if (mQATasks & taskTrackingRes) {
467 for (int32_t i = 0; i < 5; i++) {
468 for (int32_t j = 0; j < 5; j++) {
469 snprintf(name, 2048, "rms_%s_vs_%s", VSPARAMETER_NAMES[i], VSPARAMETER_NAMES[j]);
470 snprintf(fname, 1024, "mean_%s_vs_%s", VSPARAMETER_NAMES[i], VSPARAMETER_NAMES[j]);
471 if (j == 4) {
472 std::unique_ptr<double[]> binsPt{CreateLogAxis(AXIS_BINS[4], mConfig.resPrimaries == 1 ? PT_MIN_PRIM : AXES_MIN[4], AXES_MAX[4])};
473 createHist(mRes[i][j][0], name, name, AXIS_BINS[j], binsPt.get());
474 createHist(mRes[i][j][1], fname, fname, AXIS_BINS[j], binsPt.get());
475 } else {
476 createHist(mRes[i][j][0], name, name, AXIS_BINS[j], AXES_MIN[j], AXES_MAX[j]);
477 createHist(mRes[i][j][1], fname, fname, AXIS_BINS[j], AXES_MIN[j], AXES_MAX[j]);
478 }
479 snprintf(name, 2048, "res_%s_vs_%s", VSPARAMETER_NAMES[i], VSPARAMETER_NAMES[j]);
480 const float* axis = mConfig.nativeFitResolutions ? RES_AXES_NATIVE : RES_AXES;
481 const int32_t nbins = i == 4 && mConfig.nativeFitResolutions ? (10 * RES_AXIS_BINS[0]) : RES_AXIS_BINS[0];
482 if (j == 4) {
483 std::unique_ptr<double[]> binsPt{CreateLogAxis(AXIS_BINS[4], mConfig.resPrimaries == 1 ? PT_MIN_PRIM : AXES_MIN[4], AXES_MAX[4])};
484 createHist(mRes2[i][j], name, name, nbins, -axis[i], axis[i], AXIS_BINS[j], binsPt.get());
485 } else {
486 createHist(mRes2[i][j], name, name, nbins, -axis[i], axis[i], AXIS_BINS[j], AXES_MIN[j], AXES_MAX[j]);
487 }
488 }
489 }
490 }
491
492 // Create Pull Histograms
493 if (mQATasks & taskTrackingResPull) {
494 for (int32_t i = 0; i < 5; i++) {
495 for (int32_t j = 0; j < 5; j++) {
496 snprintf(name, 2048, "pull_rms_%s_vs_%s", VSPARAMETER_NAMES[i], VSPARAMETER_NAMES[j]);
497 snprintf(fname, 1024, "pull_mean_%s_vs_%s", VSPARAMETER_NAMES[i], VSPARAMETER_NAMES[j]);
498 if (j == 4) {
499 std::unique_ptr<double[]> binsPt{CreateLogAxis(AXIS_BINS[4], AXES_MIN[4], AXES_MAX[4])};
500 createHist(mPull[i][j][0], name, name, AXIS_BINS[j], binsPt.get());
501 createHist(mPull[i][j][1], fname, fname, AXIS_BINS[j], binsPt.get());
502 } else {
503 createHist(mPull[i][j][0], name, name, AXIS_BINS[j], AXES_MIN[j], AXES_MAX[j]);
504 createHist(mPull[i][j][1], fname, fname, AXIS_BINS[j], AXES_MIN[j], AXES_MAX[j]);
505 }
506 snprintf(name, 2048, "pull_%s_vs_%s", VSPARAMETER_NAMES[i], VSPARAMETER_NAMES[j]);
507 if (j == 4) {
508 std::unique_ptr<double[]> binsPt{CreateLogAxis(AXIS_BINS[4], AXES_MIN[4], AXES_MAX[4])};
509 createHist(mPull2[i][j], name, name, RES_AXIS_BINS[0], -PULL_AXIS, PULL_AXIS, AXIS_BINS[j], binsPt.get());
510 } else {
511 createHist(mPull2[i][j], name, name, RES_AXIS_BINS[0], -PULL_AXIS, PULL_AXIS, AXIS_BINS[j], AXES_MIN[j], AXES_MAX[j]);
512 }
513 }
514 }
515 }
516
517 // Create Cluster Histograms
518 if (mQATasks & taskClusterAttach) {
519 for (int32_t i = 0; i < N_CLS_TYPE * N_CLS_HIST - 1; i++) {
520 int32_t ioffset = i >= (2 * N_CLS_HIST - 1) ? (2 * N_CLS_HIST - 1) : i >= N_CLS_HIST ? N_CLS_HIST : 0;
521 int32_t itype = i >= (2 * N_CLS_HIST - 1) ? 2 : i >= N_CLS_HIST ? 1 : 0;
522 snprintf(name, 2048, "clusters%s%s", CLUSTER_NAMES_SHORT[i - ioffset], CLUSTER_TYPES[itype]);
523 std::unique_ptr<double[]> binsPt{CreateLogAxis(AXIS_BINS[4], PT_MIN_CLUST, PT_MAX)};
524 createHist(mClusters[i], name, name, AXIS_BINS[4], binsPt.get());
525 }
526 }
527
528 if (mQATasks & taskTrackStatistics) {
529 // Create Tracks Histograms
530 for (int32_t i = 0; i < 2; i++) {
531 snprintf(name, 2048, i ? "nrows_with_cluster" : "nclusters");
532 createHist(mNCl[i], name, name, 160, 0, 159);
533 }
534 snprintf(name, 2048, "tracks");
535 std::unique_ptr<double[]> binsPt{CreateLogAxis(AXIS_BINS[4], PT_MIN_CLUST, PT_MAX)};
536 createHist(mTracks, name, name, AXIS_BINS[4], binsPt.get());
537 createHist(mClXY, "clXY", "clXY", 1000, -250, 250, 1000, -250, 250);
538 }
539
540 if ((mQATasks & taskClusterCounts) && mConfig.clusterRejectionHistograms) {
541 int32_t num = DoClusterCounts(nullptr, 2);
542 mHistClusterCount.resize(num);
543 DoClusterCounts(nullptr, 1);
544 }
545
546 for (uint32_t i = 0; i < mHist1D->size(); i++) {
547 *mHist1D_pos[i] = &(*mHist1D)[i];
548 }
549 for (uint32_t i = 0; i < mHist2D->size(); i++) {
550 *mHist2D_pos[i] = &(*mHist2D)[i];
551 }
552 for (uint32_t i = 0; i < mHist1Dd->size(); i++) {
553 *mHist1Dd_pos[i] = &(*mHist1Dd)[i];
554 }
555 for (uint32_t i = 0; i < mHistGraph->size(); i++) {
556 *mHistGraph_pos[i] = &(*mHistGraph)[i];
557 }
558
559 return 0;
560}
561
562int32_t GPUQA::loadHistograms(std::vector<TH1F>& i1, std::vector<TH2F>& i2, std::vector<TH1D>& i3, std::vector<TGraphAsymmErrors>& i4, int32_t tasks)
563{
564 if (tasks == -1) {
565 tasks = taskDefaultPostprocess;
566 }
567 if (mQAInitialized && (!mHaveExternalHists || tasks != mQATasks)) {
568 throw std::runtime_error("QA not initialized or initialized with different task array");
569 }
570 mHist1D = &i1;
571 mHist2D = &i2;
572 mHist1Dd = &i3;
573 mHistGraph = &i4;
574 mHist1D_pos.clear();
575 mHist2D_pos.clear();
576 mHist1Dd_pos.clear();
577 mHistGraph_pos.clear();
578 mHaveExternalHists = true;
579 if (mConfig.noMC) {
580 tasks &= tasksNoQC;
581 }
582 mQATasks = tasks;
583 if (InitQACreateHistograms()) {
584 return 1;
585 }
586 mQAInitialized = true;
587 return 0;
588}
589
590void GPUQA::DumpO2MCData(const char* filename) const
591{
592 FILE* fp = fopen(filename, "w+b");
593 if (fp == nullptr) {
594 return;
595 }
596 uint32_t n = mMCInfos.size();
597 fwrite(&n, sizeof(n), 1, fp);
598 fwrite(mMCInfos.data(), sizeof(mMCInfos[0]), n, fp);
599 n = mMCInfosCol.size();
600 fwrite(&n, sizeof(n), 1, fp);
601 fwrite(mMCInfosCol.data(), sizeof(mMCInfosCol[0]), n, fp);
602 n = mMCEventOffset.size();
603 fwrite(&n, sizeof(n), 1, fp);
604 fwrite(mMCEventOffset.data(), sizeof(mMCEventOffset[0]), n, fp);
605 fclose(fp);
606}
607
608int32_t GPUQA::ReadO2MCData(const char* filename)
609{
610 FILE* fp = fopen(filename, "rb");
611 if (fp == nullptr) {
612 return 1;
613 }
614 uint32_t n;
615 uint32_t x;
616 if ((x = fread(&n, sizeof(n), 1, fp)) != 1) {
617 fclose(fp);
618 return 1;
619 }
620 mMCInfos.resize(n);
621 if (fread(mMCInfos.data(), sizeof(mMCInfos[0]), n, fp) != n) {
622 fclose(fp);
623 return 1;
624 }
625 if ((x = fread(&n, sizeof(n), 1, fp)) != 1) {
626 fclose(fp);
627 return 1;
628 }
629 mMCInfosCol.resize(n);
630 if (fread(mMCInfosCol.data(), sizeof(mMCInfosCol[0]), n, fp) != n) {
631 fclose(fp);
632 return 1;
633 }
634 if ((x = fread(&n, sizeof(n), 1, fp)) != 1) {
635 fclose(fp);
636 return 1;
637 }
638 mMCEventOffset.resize(n);
639 if (fread(mMCEventOffset.data(), sizeof(mMCEventOffset[0]), n, fp) != n) {
640 fclose(fp);
641 return 1;
642 }
643 if (mTracking && mTracking->GetProcessingSettings().debugLevel >= 2) {
644 printf("Read %ld bytes MC Infos\n", ftell(fp));
645 }
646 fclose(fp);
647 if (mTracking) {
648 CopyO2MCtoIOPtr(&mTracking->mIOPtrs);
649 }
650 return 0;
651}
652
653void GPUQA::CopyO2MCtoIOPtr(GPUTrackingInOutPointers* ptr)
654{
655 ptr->mcInfosTPC = mMCInfos.data();
656 ptr->nMCInfosTPC = mMCInfos.size();
657 ptr->mcInfosTPCCol = mMCInfosCol.data();
658 ptr->nMCInfosTPCCol = mMCInfosCol.size();
659}
660
661void GPUQA::InitO2MCData(GPUTrackingInOutPointers* updateIOPtr)
662{
663#ifdef GPUCA_O2_LIB
664 if (!mO2MCDataLoaded) {
665 HighResTimer timer;
666 if (mTracking && mTracking->GetProcessingSettings().debugLevel) {
667 GPUInfo("Start reading O2 Track MC information");
668 timer.Start();
669 }
670 static constexpr float PRIM_MAX_T = 0.01f;
671
672 o2::steer::MCKinematicsReader mcReader("collisioncontext.root");
673 std::vector<int32_t> refId;
674
675 auto dc = o2::steer::DigitizationContext::loadFromFile("collisioncontext.root");
676 const auto& evrec = dc->getEventRecords();
677 const auto& evparts = dc->getEventParts();
678 std::vector<std::vector<float>> evTimeBins(mcReader.getNSources());
679 for (uint32_t i = 0; i < evTimeBins.size(); i++) {
680 evTimeBins[i].resize(mcReader.getNEvents(i), -100.f);
681 }
682 for (uint32_t i = 0; i < evrec.size(); i++) {
683 const auto& ir = evrec[i];
684 for (uint32_t j = 0; j < evparts[i].size(); j++) {
685 const int iSim = evparts[i][j].sourceID;
686 const int iEv = evparts[i][j].entryID;
687 if (iSim == o2::steer::QEDSOURCEID || ir.differenceInBC(o2::raw::HBFUtils::Instance().getFirstIR()) >= 0) {
690 if (evTimeBins[iSim][iEv] >= 0) {
691 throw std::runtime_error("Multiple time bins for same MC collision found");
692 }
693 evTimeBins[iSim][iEv] = timebin;
694 }
695 }
696 }
697
698 uint32_t nSimSources = mcReader.getNSources();
699 mMCEventOffset.resize(nSimSources);
700 uint32_t nSimTotalEvents = 0;
701 uint32_t nSimTotalTracks = 0;
702 for (uint32_t i = 0; i < nSimSources; i++) {
703 mMCEventOffset[i] = nSimTotalEvents;
704 nSimTotalEvents += mcReader.getNEvents(i);
705 }
706
707 mMCInfosCol.resize(nSimTotalEvents);
708 for (int32_t iSim = 0; iSim < mcReader.getNSources(); iSim++) {
709 for (int32_t i = 0; i < mcReader.getNEvents(iSim); i++) {
710 const float timebin = evTimeBins[iSim][i];
711
712 const std::vector<o2::MCTrack>& tracks = mcReader.getTracks(iSim, i);
713 const std::vector<o2::TrackReference>& trackRefs = mcReader.getTrackRefsByEvent(iSim, i);
714
715 refId.resize(tracks.size());
716 std::fill(refId.begin(), refId.end(), -1);
717 for (uint32_t j = 0; j < trackRefs.size(); j++) {
718 if (trackRefs[j].getDetectorId() == o2::detectors::DetID::TPC) {
719 int32_t trkId = trackRefs[j].getTrackID();
720 if (refId[trkId] == -1) {
721 refId[trkId] = j;
722 }
723 }
724 }
725 mMCInfosCol[mMCEventOffset[iSim] + i].first = mMCInfos.size();
726 mMCInfosCol[mMCEventOffset[iSim] + i].num = tracks.size();
727 mMCInfos.resize(mMCInfos.size() + tracks.size());
728 for (uint32_t j = 0; j < tracks.size(); j++) {
729 auto& info = mMCInfos[mMCInfosCol[mMCEventOffset[iSim] + i].first + j];
730 const auto& trk = tracks[j];
731 TParticlePDG* particle = TDatabasePDG::Instance()->GetParticle(trk.GetPdgCode());
732 Int_t pid = -1;
733 if (abs(trk.GetPdgCode()) == kElectron) {
734 pid = 0;
735 }
736 if (abs(trk.GetPdgCode()) == kMuonMinus) {
737 pid = 1;
738 }
739 if (abs(trk.GetPdgCode()) == kPiPlus) {
740 pid = 2;
741 }
742 if (abs(trk.GetPdgCode()) == kKPlus) {
743 pid = 3;
744 }
745 if (abs(trk.GetPdgCode()) == kProton) {
746 pid = 4;
747 }
748
749 info.charge = particle ? particle->Charge() : 0;
750 info.prim = trk.T() < PRIM_MAX_T;
751 info.primDaughters = 0;
752 if (trk.getFirstDaughterTrackId() != -1) {
753 for (int32_t k = trk.getFirstDaughterTrackId(); k <= trk.getLastDaughterTrackId(); k++) {
754 if (tracks[k].T() < PRIM_MAX_T) {
755 info.primDaughters = 1;
756 break;
757 }
758 }
759 }
760 info.pid = pid;
761 info.t0 = timebin;
762 if (refId[j] >= 0) {
763 const auto& trkRef = trackRefs[refId[j]];
764 info.x = trkRef.X();
765 info.y = trkRef.Y();
766 info.z = trkRef.Z();
767 info.pX = trkRef.Px();
768 info.pY = trkRef.Py();
769 info.pZ = trkRef.Pz();
770 info.genRadius = std::sqrt(trk.GetStartVertexCoordinatesX() * trk.GetStartVertexCoordinatesX() + trk.GetStartVertexCoordinatesY() * trk.GetStartVertexCoordinatesY() + trk.GetStartVertexCoordinatesZ() * trk.GetStartVertexCoordinatesZ());
771 } else {
772 info.x = info.y = info.z = info.pX = info.pY = info.pZ = 0;
773 info.genRadius = 0;
774 }
775 }
776 }
777 }
778 if (mTracking && mTracking->GetProcessingSettings().debugLevel) {
779 GPUInfo("Finished reading O2 Track MC information (%f seconds)", timer.GetCurrentElapsedTime());
780 }
781 mO2MCDataLoaded = true;
782 }
783 if (updateIOPtr) {
784 CopyO2MCtoIOPtr(updateIOPtr);
785 }
786#endif
787}
788
789int32_t GPUQA::InitQA(int32_t tasks)
790{
791 if (mQAInitialized) {
792 throw std::runtime_error("QA already initialized");
793 }
794 if (tasks == -1) {
795 tasks = taskDefault;
796 }
797
798 mHist1D = new std::vector<TH1F>;
799 mHist2D = new std::vector<TH2F>;
800 mHist1Dd = new std::vector<TH1D>;
801 mHistGraph = new std::vector<TGraphAsymmErrors>;
802 if (mConfig.noMC) {
803 tasks &= tasksNoQC;
804 }
805 mQATasks = tasks;
806
807 if (mTracking->GetProcessingSettings().qcRunFraction != 100.f && mQATasks != taskClusterCounts) {
808 throw std::runtime_error("QA with qcRunFraction only supported for taskClusterCounts");
809 }
810
811 if (mTracking) {
812 mClNative = mTracking->mIOPtrs.clustersNative;
813 }
814
815 if (InitQACreateHistograms()) {
816 return 1;
817 }
818
819 if (mConfig.enableLocalOutput) {
820 mkdir("plots", S_IRWXU | S_IRWXG | S_IROTH | S_IXOTH);
821 }
822
823#ifdef GPUCA_O2_LIB
824 if (!mConfig.noMC) {
825 InitO2MCData(mTracking ? &mTracking->mIOPtrs : nullptr);
826 }
827#endif
828
829 if (mConfig.matchMCLabels.size()) {
830 uint32_t nFiles = mConfig.matchMCLabels.size();
831 std::vector<std::unique_ptr<TFile>> files;
832 std::vector<std::vector<std::vector<int32_t>>*> labelsBuffer(nFiles);
833 std::vector<std::vector<std::vector<int32_t>>*> effBuffer(nFiles);
834 for (uint32_t i = 0; i < nFiles; i++) {
835 files.emplace_back(std::make_unique<TFile>(mConfig.matchMCLabels[i].c_str()));
836 labelsBuffer[i] = (std::vector<std::vector<int32_t>>*)files[i]->Get("mcLabelBuffer");
837 effBuffer[i] = (std::vector<std::vector<int32_t>>*)files[i]->Get("mcEffBuffer");
838 if (labelsBuffer[i] == nullptr || effBuffer[i] == nullptr) {
839 GPUError("Error opening / reading from labels file %u/%s: %p %p", i, mConfig.matchMCLabels[i].c_str(), (void*)labelsBuffer[i], (void*)effBuffer[i]);
840 exit(1);
841 }
842 }
843
844 mGoodTracks.resize(labelsBuffer[0]->size());
845 mGoodHits.resize(labelsBuffer[0]->size());
846 for (uint32_t iEvent = 0; iEvent < labelsBuffer[0]->size(); iEvent++) {
847 std::vector<bool> labelsOK((*effBuffer[0])[iEvent].size());
848 for (uint32_t k = 0; k < (*effBuffer[0])[iEvent].size(); k++) {
849 labelsOK[k] = false;
850 for (uint32_t l = 0; l < nFiles; l++) {
851 if ((*effBuffer[0])[iEvent][k] != (*effBuffer[l])[iEvent][k]) {
852 labelsOK[k] = true;
853 break;
854 }
855 }
856 }
857 mGoodTracks[iEvent].resize((*labelsBuffer[0])[iEvent].size());
858 for (uint32_t k = 0; k < (*labelsBuffer[0])[iEvent].size(); k++) {
859 if ((*labelsBuffer[0])[iEvent][k] == MC_LABEL_INVALID) {
860 continue;
861 }
862 mGoodTracks[iEvent][k] = labelsOK[abs((*labelsBuffer[0])[iEvent][k])];
863 }
864 }
865 }
866 mQAInitialized = true;
867 return 0;
868}
869
870void GPUQA::RunQA(bool matchOnly, const std::vector<o2::tpc::TrackTPC>* tracksExternal, const std::vector<o2::MCCompLabel>* tracksExtMC, const o2::tpc::ClusterNativeAccess* clNative)
871{
872 if (!mQAInitialized) {
873 throw std::runtime_error("QA not initialized");
874 }
875 if (mTracking && mTracking->GetProcessingSettings().debugLevel >= 2) {
876 GPUInfo("Running QA - Mask %d, Efficiency %d, Resolution %d, Pulls %d, Cluster Attachment %d, Track Statistics %d, Cluster Counts %d", mQATasks, (int32_t)(mQATasks & taskTrackingEff), (int32_t)(mQATasks & taskTrackingRes), (int32_t)(mQATasks & taskTrackingResPull), (int32_t)(mQATasks & taskClusterAttach), (int32_t)(mQATasks & taskTrackStatistics), (int32_t)(mQATasks & taskClusterCounts));
877 }
878 if (!clNative && mTracking) {
879 clNative = mTracking->mIOPtrs.clustersNative;
880 }
881 mClNative = clNative;
882
883#ifdef GPUCA_TPC_GEOMETRY_O2
884 uint32_t nSimEvents = GetNMCCollissions();
885 if (mTrackMCLabelsReverse.size() < nSimEvents) {
886 mTrackMCLabelsReverse.resize(nSimEvents);
887 }
888 if (mRecTracks.size() < nSimEvents) {
889 mRecTracks.resize(nSimEvents);
890 }
891 if (mFakeTracks.size() < nSimEvents) {
892 mFakeTracks.resize(nSimEvents);
893 }
894 if (mMCParam.size() < nSimEvents) {
895 mMCParam.resize(nSimEvents);
896 }
897#endif
898
899 // Initialize Arrays
900 uint32_t nReconstructedTracks = 0;
901 if (tracksExternal) {
902#ifdef GPUCA_O2_LIB
903 nReconstructedTracks = tracksExternal->size();
904#endif
905 } else {
906 nReconstructedTracks = mTracking->mIOPtrs.nMergedTracks;
907 }
908 mTrackMCLabels.resize(nReconstructedTracks);
909 for (uint32_t iCol = 0; iCol < GetNMCCollissions(); iCol++) {
910 mTrackMCLabelsReverse[iCol].resize(GetNMCTracks(iCol));
911 mRecTracks[iCol].resize(GetNMCTracks(iCol));
912 mFakeTracks[iCol].resize(GetNMCTracks(iCol));
913 mMCParam[iCol].resize(GetNMCTracks(iCol));
914 memset(mRecTracks[iCol].data(), 0, mRecTracks[iCol].size() * sizeof(mRecTracks[iCol][0]));
915 memset(mFakeTracks[iCol].data(), 0, mFakeTracks[iCol].size() * sizeof(mFakeTracks[iCol][0]));
916 for (size_t i = 0; i < mTrackMCLabelsReverse[iCol].size(); i++) {
917 mTrackMCLabelsReverse[iCol][i] = -1;
918 }
919 }
920 if (mQATasks & taskClusterAttach && GetNMCLabels()) {
921 mClusterParam.resize(GetNMCLabels());
922 memset(mClusterParam.data(), 0, mClusterParam.size() * sizeof(mClusterParam[0]));
923 }
924 HighResTimer timer;
925
926 mNEvents++;
927 if (mConfig.writeMCLabels) {
928 mcEffBuffer.resize(mNEvents);
929 mcLabelBuffer.resize(mNEvents);
930 mcEffBuffer[mNEvents - 1].resize(GetNMCTracks(0));
931 mcLabelBuffer[mNEvents - 1].resize(nReconstructedTracks);
932 }
933
934 bool mcAvail = mcPresent() || tracksExtMC;
935
936 if (mcAvail) {
937 // Assign Track MC Labels
938 timer.Start();
939 if (tracksExternal) {
940#ifdef GPUCA_O2_LIB
941 for (uint32_t i = 0; i < tracksExternal->size(); i++) {
942 mTrackMCLabels[i] = (*tracksExtMC)[i];
943 }
944#endif
945 } else {
946 tbb::parallel_for(tbb::blocked_range<uint32_t>(0, nReconstructedTracks, (QA_DEBUG == 0) ? 32 : nReconstructedTracks), [&](const tbb::blocked_range<uint32_t>& range) {
947 auto acc = GPUTPCTrkLbl<true, mcLabelI_t>(GetClusterLabels(), 1.f - mConfig.recThreshold);
948 for (auto i = range.begin(); i < range.end(); i++) {
949 acc.reset();
950 int32_t nClusters = 0;
951 const GPUTPCGMMergedTrack& track = mTracking->mIOPtrs.mergedTracks[i];
952 std::vector<mcLabel_t> labels;
953 for (uint32_t k = 0; k < track.NClusters(); k++) {
954 if (mTracking->mIOPtrs.mergedTrackHits[track.FirstClusterRef() + k].state & GPUTPCGMMergedTrackHit::flagReject) {
955 continue;
956 }
957 nClusters++;
958 uint32_t hitId = mTracking->mIOPtrs.mergedTrackHits[track.FirstClusterRef() + k].num;
959 if (hitId >= GetNMCLabels()) {
960 GPUError("Invalid hit id %u > %d (nClusters %d)", hitId, GetNMCLabels(), mTracking->mIOPtrs.clustersNative ? mTracking->mIOPtrs.clustersNative->nClustersTotal : 0);
961 throw std::runtime_error("qa error");
962 }
963 acc.addLabel(hitId);
964 for (int32_t j = 0; j < GetMCLabelNID(hitId); j++) {
965 if (GetMCLabelID(hitId, j) >= (int32_t)GetNMCTracks(GetMCLabelCol(hitId, j))) {
966 GPUError("Invalid label %d > %d (hit %d, label %d, col %d)", GetMCLabelID(hitId, j), GetNMCTracks(GetMCLabelCol(hitId, j)), hitId, j, (int32_t)GetMCLabelCol(hitId, j));
967 throw std::runtime_error("qa error");
968 }
969 if (GetMCLabelID(hitId, j) >= 0) {
970 if (QA_DEBUG >= 3 && track.OK()) {
971 GPUInfo("Track %d Cluster %u Label %d: %d (%f)", i, k, j, GetMCLabelID(hitId, j), GetMCLabelWeight(hitId, j));
972 }
973 }
974 }
975 }
976
977 float maxweight, sumweight;
978 int32_t maxcount;
979 auto maxLabel = acc.computeLabel(&maxweight, &sumweight, &maxcount);
980 mTrackMCLabels[i] = maxLabel;
981 if (QA_DEBUG && track.OK() && GetNMCTracks(maxLabel) > (uint32_t)maxLabel.getTrackID()) {
982 const mcInfo_t& mc = GetMCTrack(maxLabel);
983 GPUInfo("Track %d label %d (fake %d) weight %f clusters %d (fitted %d) (%f%% %f%%) Pt %f", i, maxLabel.getTrackID(), (int32_t)(maxLabel.isFake()), maxweight, nClusters, track.NClustersFitted(), 100.f * maxweight / sumweight, 100.f * (float)maxcount / (float)nClusters,
984 std::sqrt(mc.pX * mc.pX + mc.pY * mc.pY));
985 }
986 }
987 });
988 }
989 if (QA_TIMING || (mTracking && mTracking->GetProcessingSettings().debugLevel >= 3)) {
990 GPUInfo("QA Time: Assign Track Labels:\t\t%6.0f us", timer.GetCurrentElapsedTime(true) * 1e6);
991 }
992
993 for (uint32_t i = 0; i < nReconstructedTracks; i++) {
994 const GPUTPCGMMergedTrack* track = mTracking ? &mTracking->mIOPtrs.mergedTracks[i] : nullptr;
995 mcLabelI_t label = mTrackMCLabels[i];
996 if (mQATasks & taskClusterAttach) {
997 // fill cluster attachment status
998 if (!track->OK()) {
999 continue;
1000 }
1001 if (!mTrackMCLabels[i].isValid()) {
1002 for (uint32_t k = 0; k < track->NClusters(); k++) {
1003 if (mTracking->mIOPtrs.mergedTrackHits[track->FirstClusterRef() + k].state & GPUTPCGMMergedTrackHit::flagReject) {
1004 continue;
1005 }
1006 mClusterParam[mTracking->mIOPtrs.mergedTrackHits[track->FirstClusterRef() + k].num].fakeAttached++;
1007 }
1008 continue;
1009 }
1010 if (mMCTrackMin == -1 || (label.getTrackID() >= mMCTrackMin && label.getTrackID() < mMCTrackMax)) {
1011 for (uint32_t k = 0; k < track->NClusters(); k++) {
1012 if (mTracking->mIOPtrs.mergedTrackHits[track->FirstClusterRef() + k].state & GPUTPCGMMergedTrackHit::flagReject) {
1013 continue;
1014 }
1015 int32_t hitId = mTracking->mIOPtrs.mergedTrackHits[track->FirstClusterRef() + k].num;
1016 bool correct = false;
1017 for (int32_t j = 0; j < GetMCLabelNID(hitId); j++) {
1018 if (label == GetMCLabel(hitId, j)) {
1019 correct = true;
1020 break;
1021 }
1022 }
1023 if (correct) {
1024 mClusterParam[hitId].attached++;
1025 } else {
1026 mClusterParam[hitId].fakeAttached++;
1027 }
1028 }
1029 }
1030 }
1031
1032 if (mTrackMCLabels[i].isFake()) {
1033 (GetMCTrackObj(mFakeTracks, label))++;
1034 } else if (tracksExternal || !track->MergedLooper()) {
1035 GetMCTrackObj(mRecTracks, label)++;
1036 if (mMCTrackMin == -1 || (label.getTrackID() >= mMCTrackMin && label.getTrackID() < mMCTrackMax)) {
1037 int32_t& revLabel = GetMCTrackObj(mTrackMCLabelsReverse, label);
1038 if (tracksExternal) {
1039#ifdef GPUCA_O2_LIB
1040 if (revLabel == -1 || fabsf((*tracksExternal)[i].getZ()) < fabsf((*tracksExternal)[revLabel].getZ())) {
1041 revLabel = i;
1042 }
1043#endif
1044 } else {
1045 const auto* trks = mTracking->mIOPtrs.mergedTracks;
1046 bool comp;
1047 if (revLabel == -1) {
1048 comp = true;
1049 } else if (mTracking->GetParam().par.earlyTpcTransform) {
1050 comp = fabsf(trks[i].GetParam().GetZ() + trks[i].GetParam().GetTZOffset()) < fabsf(trks[revLabel].GetParam().GetZ() + trks[revLabel].GetParam().GetTZOffset());
1051 } else {
1052 float shift1 = mTracking->GetTPCTransformHelper()->getCorrMap()->convDeltaTimeToDeltaZinTimeFrame(trks[i].CSide() * GPUChainTracking::NSECTORS / 2, trks[i].GetParam().GetTZOffset());
1053 float shift2 = mTracking->GetTPCTransformHelper()->getCorrMap()->convDeltaTimeToDeltaZinTimeFrame(trks[revLabel].CSide() * GPUChainTracking::NSECTORS / 2, trks[revLabel].GetParam().GetTZOffset());
1054 comp = fabsf(trks[i].GetParam().GetZ() + shift1) < fabsf(trks[revLabel].GetParam().GetZ() + shift2);
1055 }
1056 if (revLabel == -1 || !trks[revLabel].OK() || (trks[i].OK() && comp)) {
1057 revLabel = i;
1058 }
1059 }
1060 }
1061 }
1062 }
1063 if ((mQATasks & taskClusterAttach) && mTracking->mIOPtrs.mergedTrackHitAttachment) {
1064 // fill cluster adjacent status
1065 for (uint32_t i = 0; i < GetNMCLabels(); i++) {
1066 if (mClusterParam[i].attached == 0 && mClusterParam[i].fakeAttached == 0) {
1067 int32_t attach = mTracking->mIOPtrs.mergedTrackHitAttachment[i];
1069 int32_t track = attach & gputpcgmmergertypes::attachTrackMask;
1070 mcLabelI_t trackL = mTrackMCLabels[track];
1071 bool fake = true;
1072 for (int32_t j = 0; j < GetMCLabelNID(i); j++) {
1073 // GPUInfo("Attach %x Track %d / %d:%d", attach, track, j, GetMCLabelID(i, j));
1074 if (trackL == GetMCLabel(i, j)) {
1075 fake = false;
1076 break;
1077 }
1078 }
1079 if (fake) {
1080 mClusterParam[i].fakeAdjacent++;
1081 } else {
1082 mClusterParam[i].adjacent++;
1083 }
1084 }
1085 }
1086 }
1087 }
1088
1089 if (mConfig.matchMCLabels.size()) {
1090 mGoodHits[mNEvents - 1].resize(GetNMCLabels());
1091 std::vector<bool> allowMCLabels(GetNMCTracks(0));
1092 for (uint32_t k = 0; k < GetNMCTracks(0); k++) {
1093 allowMCLabels[k] = false;
1094 }
1095 for (uint32_t i = 0; i < nReconstructedTracks; i++) {
1096 if (!mGoodTracks[mNEvents - 1][i]) {
1097 continue;
1098 }
1099 if (mConfig.matchDisplayMinPt > 0) {
1100 if (!mTrackMCLabels[i].isValid()) {
1101 continue;
1102 }
1103 const mcInfo_t& info = GetMCTrack(mTrackMCLabels[i]);
1104 if (info.pX * info.pX + info.pY * info.pY < mConfig.matchDisplayMinPt * mConfig.matchDisplayMinPt) {
1105 continue;
1106 }
1107 }
1108
1109 const GPUTPCGMMergedTrack& track = mTracking->mIOPtrs.mergedTracks[i];
1110 for (uint32_t j = 0; j < track.NClusters(); j++) {
1111 int32_t hitId = mTracking->mIOPtrs.mergedTrackHits[track.FirstClusterRef() + j].num;
1112 if (GetMCLabelNID(hitId)) {
1113 int32_t mcID = GetMCLabelID(hitId, 0);
1114 if (mcID >= 0) {
1115 allowMCLabels[mcID] = true;
1116 }
1117 }
1118 }
1119 }
1120 for (uint32_t i = 0; i < GetNMCLabels(); i++) {
1121 for (int32_t j = 0; j < GetMCLabelNID(i); j++) {
1122 int32_t mcID = GetMCLabelID(i, j);
1123 if (mcID >= 0 && allowMCLabels[mcID]) {
1124 mGoodHits[mNEvents - 1][i] = true;
1125 }
1126 }
1127 }
1128 }
1129 if (QA_TIMING || (mTracking && mTracking->GetProcessingSettings().debugLevel >= 3)) {
1130 GPUInfo("QA Time: Cluster attach status:\t\t%6.0f us", timer.GetCurrentElapsedTime(true) * 1e6);
1131 }
1132
1133 if (matchOnly) {
1134 return;
1135 }
1136
1137 // Recompute fNWeightCls (might have changed after merging events into timeframes)
1138 for (uint32_t iCol = 0; iCol < GetNMCCollissions(); iCol++) {
1139 for (uint32_t i = 0; i < GetNMCTracks(iCol); i++) {
1140 mMCParam[iCol][i].nWeightCls = 0.;
1141 }
1142 }
1143 for (uint32_t i = 0; i < GetNMCLabels(); i++) {
1144 float weightTotal = 0.f;
1145 for (int32_t j = 0; j < GetMCLabelNID(i); j++) {
1146 if (GetMCLabelID(i, j) >= 0) {
1147 weightTotal += GetMCLabelWeight(i, j);
1148 }
1149 }
1150 for (int32_t j = 0; j < GetMCLabelNID(i); j++) {
1151 if (GetMCLabelID(i, j) >= 0) {
1152 GetMCTrackObj(mMCParam, GetMCLabel(i, j)).nWeightCls += GetMCLabelWeight(i, j) / weightTotal;
1153 }
1154 }
1155 }
1156 if (QA_TIMING || (mTracking && mTracking->GetProcessingSettings().debugLevel >= 3)) {
1157 GPUInfo("QA Time: Compute cluster label weights:\t%6.0f us", timer.GetCurrentElapsedTime(true) * 1e6);
1158 }
1159
1160 // Compute MC Track Parameters for MC Tracks
1161 tbb::parallel_for<uint32_t>(0, GetNMCCollissions(), [&](auto iCol) {
1162 for (uint32_t i = 0; i < GetNMCTracks(iCol); i++) {
1163 const mcInfo_t& info = GetMCTrack(i, iCol);
1164 additionalMCParameters& mc2 = mMCParam[iCol][i];
1165 mc2.pt = std::sqrt(info.pX * info.pX + info.pY * info.pY);
1166 mc2.phi = M_PI + std::atan2(-info.pY, -info.pX);
1167 float p = info.pX * info.pX + info.pY * info.pY + info.pZ * info.pZ;
1168 if (p < 1e-18) {
1169 mc2.theta = mc2.eta = 0.f;
1170 } else {
1171 mc2.theta = info.pZ == 0 ? (M_PI / 2) : (std::acos(info.pZ / std::sqrt(p)));
1172 mc2.eta = -std::log(std::tan(0.5 * mc2.theta));
1173 }
1174 if (mConfig.writeMCLabels) {
1175 std::vector<int32_t>& effBuffer = mcEffBuffer[mNEvents - 1];
1176 effBuffer[i] = mRecTracks[iCol][i] * 1000 + mFakeTracks[iCol][i];
1177 }
1178 } // clang-format off
1179 }, tbb::simple_partitioner()); // clang-format on
1180 if (QA_TIMING || (mTracking && mTracking->GetProcessingSettings().debugLevel >= 3)) {
1181 GPUInfo("QA Time: Compute track mc parameters:\t%6.0f us", timer.GetCurrentElapsedTime(true) * 1e6);
1182 }
1183
1184 // Fill Efficiency Histograms
1185 if (mQATasks & taskTrackingEff) {
1186 for (uint32_t iCol = 0; iCol < GetNMCCollissions(); iCol++) {
1187 for (uint32_t i = 0; i < GetNMCTracks(iCol); i++) {
1188 if ((mMCTrackMin != -1 && (int32_t)i < mMCTrackMin) || (mMCTrackMax != -1 && (int32_t)i >= mMCTrackMax)) {
1189 continue;
1190 }
1191 const mcInfo_t& info = GetMCTrack(i, iCol);
1192 const additionalMCParameters& mc2 = mMCParam[iCol][i];
1193 if (mc2.nWeightCls == 0.f) {
1194 continue;
1195 }
1196 const float& mcpt = mc2.pt;
1197 const float& mcphi = mc2.phi;
1198 const float& mceta = mc2.eta;
1199
1200 if (info.primDaughters) {
1201 continue;
1202 }
1203 if (mc2.nWeightCls < mConfig.minNClEff) {
1204 continue;
1205 }
1206 int32_t findable = mc2.nWeightCls >= mConfig.minNClFindable;
1207 if (info.pid < 0) {
1208 continue;
1209 }
1210 if (info.charge == 0.f) {
1211 continue;
1212 }
1213 if (mConfig.filterCharge && info.charge * mConfig.filterCharge < 0) {
1214 continue;
1215 }
1216 if (mConfig.filterPID >= 0 && info.pid != mConfig.filterPID) {
1217 continue;
1218 }
1219
1220 if (fabsf(mceta) > ETA_MAX || mcpt < PT_MIN || mcpt > PT_MAX) {
1221 continue;
1222 }
1223
1224 float alpha = std::atan2(info.y, info.x);
1225 alpha /= M_PI / 9.f;
1226 alpha = std::floor(alpha);
1227 alpha *= M_PI / 9.f;
1228 alpha += M_PI / 18.f;
1229
1230 float c = std::cos(alpha);
1231 float s = std::sin(alpha);
1232 float localY = -info.x * s + info.y * c;
1233
1234 if (mConfig.dumpToROOT) {
1235 static auto effdump = GPUROOTDump<TNtuple>::getNew("eff", "alpha:x:y:z:mcphi:mceta:mcpt:rec:fake:findable:prim:ncls");
1236 float localX = info.x * c + info.y * s;
1237 effdump.Fill(alpha, localX, localY, info.z, mcphi, mceta, mcpt, mRecTracks[iCol][i], mFakeTracks[iCol][i], findable, info.prim, mc2.nWeightCls);
1238 }
1239
1240 for (int32_t j = 0; j < 4; j++) {
1241 for (int32_t k = 0; k < 2; k++) {
1242 if (k == 0 && findable == 0) {
1243 continue;
1244 }
1245
1246 int32_t val = (j == 0) ? (mRecTracks[iCol][i] ? 1 : 0) : (j == 1) ? (mRecTracks[iCol][i] ? mRecTracks[iCol][i] - 1 : 0) : (j == 2) ? mFakeTracks[iCol][i] : 1;
1247 if (val == 0) {
1248 continue;
1249 }
1250
1251 for (int32_t l = 0; l < 5; l++) {
1252 if (info.prim && mcpt < PT_MIN_PRIM) {
1253 continue;
1254 }
1255 if (l != 3 && fabsf(mceta) > ETA_MAX2) {
1256 continue;
1257 }
1258 if (l < 4 && mcpt < 1.f / mConfig.qpt) {
1259 continue;
1260 }
1261
1262 float pos = l == 0 ? localY : l == 1 ? info.z : l == 2 ? mcphi : l == 3 ? mceta : mcpt;
1263
1264 mEff[j][k][!info.prim][l]->Fill(pos, val);
1265 }
1266 }
1267 }
1268 }
1269 }
1270 if (QA_TIMING || (mTracking && mTracking->GetProcessingSettings().debugLevel >= 3)) {
1271 GPUInfo("QA Time: Fill efficiency histograms:\t%6.0f us", timer.GetCurrentElapsedTime(true) * 1e6);
1272 }
1273 }
1274
1275 // Fill Resolution Histograms
1276 if (mQATasks & (taskTrackingRes | taskTrackingResPull)) {
1277 GPUTPCGMPropagator prop;
1278 prop.SetMaxSinPhi(.999);
1279 prop.SetMaterialTPC();
1280 prop.SetPolynomialField(&mParam->polynomialField);
1281
1282 for (uint32_t i = 0; i < mTrackMCLabels.size(); i++) {
1283 if (mConfig.writeMCLabels) {
1284 std::vector<int32_t>& labelBuffer = mcLabelBuffer[mNEvents - 1];
1285 labelBuffer[i] = mTrackMCLabels[i].getTrackID();
1286 }
1287 if (mTrackMCLabels[i].isFake()) {
1288 continue;
1289 }
1290 const mcInfo_t& mc1 = GetMCTrack(mTrackMCLabels[i]);
1291 const additionalMCParameters& mc2 = GetMCTrackObj(mMCParam, mTrackMCLabels[i]);
1292
1293 if (mc1.primDaughters) {
1294 continue;
1295 }
1296 if (!tracksExternal) {
1297 if (!mTracking->mIOPtrs.mergedTracks[i].OK()) {
1298 continue;
1299 }
1300 if (mTracking->mIOPtrs.mergedTracks[i].MergedLooper()) {
1301 continue;
1302 }
1303 }
1304 if ((mMCTrackMin != -1 && mTrackMCLabels[i].getTrackID() < mMCTrackMin) || (mMCTrackMax != -1 && mTrackMCLabels[i].getTrackID() >= mMCTrackMax)) {
1305 continue;
1306 }
1307 if (fabsf(mc2.eta) > ETA_MAX || mc2.pt < PT_MIN || mc2.pt > PT_MAX) {
1308 continue;
1309 }
1310 if (mc1.charge == 0.f) {
1311 continue;
1312 }
1313 if (mc1.pid < 0) {
1314 continue;
1315 }
1316 if (mc1.t0 == -100.f) {
1317 continue;
1318 }
1319 if (mConfig.filterCharge && mc1.charge * mConfig.filterCharge < 0) {
1320 continue;
1321 }
1322 if (mConfig.filterPID >= 0 && mc1.pid != mConfig.filterPID) {
1323 continue;
1324 }
1325 if (mc2.nWeightCls < mConfig.minNClRes) {
1326 continue;
1327 }
1328 if (mConfig.resPrimaries == 1 && !mc1.prim) {
1329 continue;
1330 } else if (mConfig.resPrimaries == 2 && mc1.prim) {
1331 continue;
1332 }
1333 if (GetMCTrackObj(mTrackMCLabelsReverse, mTrackMCLabels[i]) != (int32_t)i) {
1334 continue;
1335 }
1336
1338 float alpha = 0.f;
1339 int32_t side;
1340 if (tracksExternal) {
1341#ifdef GPUCA_O2_LIB
1342 for (int32_t k = 0; k < 5; k++) {
1343 param.Par()[k] = (*tracksExternal)[i].getParams()[k];
1344 }
1345 for (int32_t k = 0; k < 15; k++) {
1346 param.Cov()[k] = (*tracksExternal)[i].getCov()[k];
1347 }
1348 param.X() = (*tracksExternal)[i].getX();
1349 param.TZOffset() = (*tracksExternal)[i].getTime0();
1350 alpha = (*tracksExternal)[i].getAlpha();
1351 side = (*tracksExternal)[i].hasBothSidesClusters() ? 2 : ((*tracksExternal)[i].hasCSideClusters() ? 1 : 0);
1352#endif
1353 } else {
1354 param = mTracking->mIOPtrs.mergedTracks[i].GetParam();
1355 alpha = mTracking->mIOPtrs.mergedTracks[i].GetAlpha();
1356 side = mTracking->mIOPtrs.mergedTracks[i].CCE() ? 2 : (mTracking->mIOPtrs.mergedTracks[i].CSide() ? 1 : 0);
1357 }
1358
1359 float mclocal[4]; // Rotated x,y,Px,Py mc-coordinates - the MC data should be rotated since the track is propagated best along x
1360 float c = std::cos(alpha);
1361 float s = std::sin(alpha);
1362 float x = mc1.x;
1363 float y = mc1.y;
1364 mclocal[0] = x * c + y * s;
1365 mclocal[1] = -x * s + y * c;
1366 float px = mc1.pX;
1367 float py = mc1.pY;
1368 mclocal[2] = px * c + py * s;
1369 mclocal[3] = -px * s + py * c;
1370
1371 if (mclocal[0] < TRACK_EXPECTED_REFERENCE_X - 3) {
1372 continue;
1373 }
1374 if (mclocal[0] > param.GetX() + 20) {
1375 continue;
1376 }
1377 if (param.GetX() > mConfig.maxResX) {
1378 continue;
1379 }
1380
1381 auto getdz = [this, &param, &mc1, &side, tracksExternal]() {
1382 if (tracksExternal) {
1383 return param.GetZ();
1384 }
1385 if (!mParam->continuousMaxTimeBin) {
1386 return param.GetZ() - mc1.z;
1387 }
1388#ifdef GPUCA_TPC_GEOMETRY_O2
1389 if (!mParam->par.earlyTpcTransform) {
1390 float shift = side == 2 ? 0 : mTracking->GetTPCTransformHelper()->getCorrMap()->convDeltaTimeToDeltaZinTimeFrame(side * GPUChainTracking::NSECTORS / 2, param.GetTZOffset() - mc1.t0);
1391 return param.GetZ() + shift - mc1.z;
1392 }
1393#endif
1394 return param.Z() + param.TZOffset() - mc1.z;
1395 };
1396
1397 prop.SetTrack(&param, alpha);
1398 bool inFlyDirection = 0;
1399 if (mConfig.strict) {
1400 const float dx = param.X() - std::max<float>(mclocal[0], TRACK_EXPECTED_REFERENCE_X_DEFAULT); // Limit distance check
1401 const float dy = param.Y() - mclocal[1];
1402 const float dz = getdz();
1403 if (dx * dx + dy * dy + dz * dz > 5.f * 5.f) {
1404 continue;
1405 }
1406 }
1407
1408 if (prop.PropagateToXAlpha(mclocal[0], alpha, inFlyDirection)) {
1409 continue;
1410 }
1411 if (fabsf(param.Y() - mclocal[1]) > (mConfig.strict ? 1.f : 4.f) || fabsf(getdz()) > (mConfig.strict ? 1.f : 4.f)) {
1412 continue;
1413 }
1414 float charge = mc1.charge > 0 ? 1.f : -1.f;
1415
1416 float deltaY = param.GetY() - mclocal[1];
1417 float deltaZ = getdz();
1418 float deltaPhiNative = param.GetSinPhi() - mclocal[3] / mc2.pt;
1419 float deltaPhi = std::asin(param.GetSinPhi()) - std::atan2(mclocal[3], mclocal[2]);
1420 float deltaLambdaNative = param.GetDzDs() - mc1.pZ / mc2.pt;
1421 float deltaLambda = std::atan(param.GetDzDs()) - std::atan2(mc1.pZ, mc2.pt);
1422 float deltaPtNative = (param.GetQPt() - charge / mc2.pt) * charge;
1423 float deltaPt = (fabsf(1.f / param.GetQPt()) - mc2.pt) / mc2.pt;
1424
1425 float paramval[5] = {mclocal[1], mc1.z, mc2.phi, mc2.eta, mc2.pt};
1426 float resval[5] = {deltaY, deltaZ, mConfig.nativeFitResolutions ? deltaPhiNative : deltaPhi, mConfig.nativeFitResolutions ? deltaLambdaNative : deltaLambda, mConfig.nativeFitResolutions ? deltaPtNative : deltaPt};
1427 float pullval[5] = {deltaY / std::sqrt(param.GetErr2Y()), deltaZ / std::sqrt(param.GetErr2Z()), deltaPhiNative / std::sqrt(param.GetErr2SinPhi()), deltaLambdaNative / std::sqrt(param.GetErr2DzDs()), deltaPtNative / std::sqrt(param.GetErr2QPt())};
1428
1429 for (int32_t j = 0; j < 5; j++) {
1430 for (int32_t k = 0; k < 5; k++) {
1431 if (k != 3 && fabsf(mc2.eta) > ETA_MAX2) {
1432 continue;
1433 }
1434 if (k < 4 && mc2.pt < 1.f / mConfig.qpt) {
1435 continue;
1436 }
1437 if (mQATasks & taskTrackingRes) {
1438 mRes2[j][k]->Fill(resval[j], paramval[k]);
1439 }
1440 if (mQATasks & taskTrackingResPull) {
1441 mPull2[j][k]->Fill(pullval[j], paramval[k]);
1442 }
1443 }
1444 }
1445 }
1446 if (QA_TIMING || (mTracking && mTracking->GetProcessingSettings().debugLevel >= 3)) {
1447 GPUInfo("QA Time: Fill resolution histograms:\t%6.0f us", timer.GetCurrentElapsedTime(true) * 1e6);
1448 }
1449 }
1450
1451 if (mQATasks & taskClusterAttach) {
1452 // Fill cluster histograms
1453 for (uint32_t iTrk = 0; iTrk < nReconstructedTracks; iTrk++) {
1454 const GPUTPCGMMergedTrack& track = mTracking->mIOPtrs.mergedTracks[iTrk];
1455 if (!track.OK()) {
1456 continue;
1457 }
1458 if (!mTrackMCLabels[iTrk].isValid()) {
1459 for (uint32_t k = 0; k < track.NClusters(); k++) {
1460 if (mTracking->mIOPtrs.mergedTrackHits[track.FirstClusterRef() + k].state & GPUTPCGMMergedTrackHit::flagReject) {
1461 continue;
1462 }
1463 int32_t hitId = mTracking->mIOPtrs.mergedTrackHits[track.FirstClusterRef() + k].num;
1464 float totalWeight = 0.;
1465 for (int32_t j = 0; j < GetMCLabelNID(hitId); j++) {
1466 if (GetMCLabelID(hitId, j) >= 0 && GetMCTrackObj(mMCParam, GetMCLabel(hitId, j)).pt > GPUCA_MIN_TRACK_PTB5_DEFAULT) {
1467 totalWeight += GetMCLabelWeight(hitId, j);
1468 }
1469 }
1470 int32_t attach = mTracking->mIOPtrs.mergedTrackHitAttachment[hitId];
1472 if (totalWeight > 0) {
1473 float weight = 1.f / (totalWeight * (mClusterParam[hitId].attached + mClusterParam[hitId].fakeAttached));
1474 for (int32_t j = 0; j < GetMCLabelNID(hitId); j++) {
1475 mcLabelI_t label = GetMCLabel(hitId, j);
1476 if (!label.isFake() && GetMCTrackObj(mMCParam, label).pt > GPUCA_MIN_TRACK_PTB5_DEFAULT) {
1477 float pt = GetMCTrackObj(mMCParam, label).pt;
1478 if (pt < PT_MIN_CLUST) {
1479 pt = PT_MIN_CLUST;
1480 }
1481 mClusters[CL_fake]->Fill(pt, GetMCLabelWeight(hitId, j) * weight);
1482 mClusters[CL_att_adj]->Fill(pt, GetMCLabelWeight(hitId, j) * weight);
1483 if (GetMCTrackObj(mRecTracks, label)) {
1484 mClusters[CL_tracks]->Fill(pt, GetMCLabelWeight(hitId, j) * weight);
1485 }
1486 mClusters[CL_all]->Fill(pt, GetMCLabelWeight(hitId, j) * weight);
1487 if (protect || physics) {
1488 mClusters[CL_prot]->Fill(pt, GetMCLabelWeight(hitId, j) * weight);
1489 }
1490 if (physics) {
1491 mClusters[CL_physics]->Fill(pt, GetMCLabelWeight(hitId, j) * weight);
1492 }
1493 }
1494 }
1495 } else {
1496 float weight = 1.f / (mClusterParam[hitId].attached + mClusterParam[hitId].fakeAttached);
1497 mClusters[CL_fake]->Fill(0.f, weight);
1498 mClusters[CL_att_adj]->Fill(0.f, weight);
1499 mClusters[CL_all]->Fill(0.f, weight);
1500 mClusterCounts.nUnaccessible += weight;
1501 if (protect || physics) {
1502 mClusters[CL_prot]->Fill(0.f, weight);
1503 }
1504 if (physics) {
1505 mClusters[CL_physics]->Fill(0.f, weight);
1506 }
1507 }
1508 }
1509 continue;
1510 }
1511 mcLabelI_t label = mTrackMCLabels[iTrk];
1512 if (mMCTrackMin != -1 && (label.getTrackID() < mMCTrackMin || label.getTrackID() >= mMCTrackMax)) {
1513 continue;
1514 }
1515 for (uint32_t k = 0; k < track.NClusters(); k++) {
1516 if (mTracking->mIOPtrs.mergedTrackHits[track.FirstClusterRef() + k].state & GPUTPCGMMergedTrackHit::flagReject) {
1517 continue;
1518 }
1519 int32_t hitId = mTracking->mIOPtrs.mergedTrackHits[track.FirstClusterRef() + k].num;
1520 float pt = GetMCTrackObj(mMCParam, label).pt;
1521 if (pt < PT_MIN_CLUST) {
1522 pt = PT_MIN_CLUST;
1523 }
1524 float weight = 1.f / (mClusterParam[hitId].attached + mClusterParam[hitId].fakeAttached);
1525 bool correct = false;
1526 for (int32_t j = 0; j < GetMCLabelNID(hitId); j++) {
1527 if (label == GetMCLabel(hitId, j)) {
1528 correct = true;
1529 break;
1530 }
1531 }
1532 if (correct) {
1533 mClusters[CL_attached]->Fill(pt, weight);
1534 mClusters[CL_tracks]->Fill(pt, weight);
1535 } else {
1536 mClusters[CL_fake]->Fill(pt, weight);
1537 }
1538 mClusters[CL_att_adj]->Fill(pt, weight);
1539 mClusters[CL_all]->Fill(pt, weight);
1540 int32_t attach = mTracking->mIOPtrs.mergedTrackHitAttachment[hitId];
1542 if (protect || physics) {
1543 mClusters[CL_prot]->Fill(pt, weight);
1544 }
1545 if (physics) {
1546 mClusters[CL_physics]->Fill(pt, weight);
1547 }
1548 }
1549 }
1550 for (uint32_t i = 0; i < GetNMCLabels(); i++) {
1551 if ((mMCTrackMin != -1 && GetMCLabelID(i, 0) < mMCTrackMin) || (mMCTrackMax != -1 && GetMCLabelID(i, 0) >= mMCTrackMax)) {
1552 continue;
1553 }
1554 if (mClusterParam[i].attached || mClusterParam[i].fakeAttached) {
1555 continue;
1556 }
1557 int32_t attach = mTracking->mIOPtrs.mergedTrackHitAttachment[i];
1559 if (mClusterParam[i].adjacent) {
1560 int32_t label = mTracking->mIOPtrs.mergedTrackHitAttachment[i] & gputpcgmmergertypes::attachTrackMask;
1561 if (!mTrackMCLabels[label].isValid()) {
1562 float totalWeight = 0.;
1563 for (int32_t j = 0; j < GetMCLabelNID(i); j++) {
1564 mcLabelI_t labelT = GetMCLabel(i, j);
1565 if (!labelT.isFake() && GetMCTrackObj(mMCParam, labelT).pt > GPUCA_MIN_TRACK_PTB5_DEFAULT) {
1566 totalWeight += GetMCLabelWeight(i, j);
1567 }
1568 }
1569 float weight = 1.f / totalWeight;
1570 if (totalWeight > 0) {
1571 for (int32_t j = 0; j < GetMCLabelNID(i); j++) {
1572 mcLabelI_t labelT = GetMCLabel(i, j);
1573 if (!labelT.isFake() && GetMCTrackObj(mMCParam, labelT).pt > GPUCA_MIN_TRACK_PTB5_DEFAULT) {
1574 float pt = GetMCTrackObj(mMCParam, labelT).pt;
1575 if (pt < PT_MIN_CLUST) {
1576 pt = PT_MIN_CLUST;
1577 }
1578 if (GetMCTrackObj(mRecTracks, labelT)) {
1579 mClusters[CL_tracks]->Fill(pt, GetMCLabelWeight(i, j) * weight);
1580 }
1581 mClusters[CL_att_adj]->Fill(pt, GetMCLabelWeight(i, j) * weight);
1582 mClusters[CL_fakeAdj]->Fill(pt, GetMCLabelWeight(i, j) * weight);
1583 mClusters[CL_all]->Fill(pt, GetMCLabelWeight(i, j) * weight);
1584 if (protect || physics) {
1585 mClusters[CL_prot]->Fill(pt, GetMCLabelWeight(i, j) * weight);
1586 }
1587 if (physics) {
1588 mClusters[CL_physics]->Fill(pt, GetMCLabelWeight(i, j) * weight);
1589 }
1590 }
1591 }
1592 } else {
1593 mClusters[CL_att_adj]->Fill(0.f, 1.f);
1594 mClusters[CL_fakeAdj]->Fill(0.f, 1.f);
1595 mClusters[CL_all]->Fill(0.f, 1.f);
1596 mClusterCounts.nUnaccessible++;
1597 if (protect || physics) {
1598 mClusters[CL_prot]->Fill(0.f, 1.f);
1599 }
1600 if (physics) {
1601 mClusters[CL_physics]->Fill(0.f, 1.f);
1602 }
1603 }
1604 } else {
1605 float pt = GetMCTrackObj(mMCParam, mTrackMCLabels[label]).pt;
1606 if (pt < PT_MIN_CLUST) {
1607 pt = PT_MIN_CLUST;
1608 }
1609 mClusters[CL_att_adj]->Fill(pt, 1.f);
1610 mClusters[CL_tracks]->Fill(pt, 1.f);
1611 mClusters[CL_all]->Fill(pt, 1.f);
1612 if (protect || physics) {
1613 mClusters[CL_prot]->Fill(pt, 1.f);
1614 }
1615 if (physics) {
1616 mClusters[CL_physics]->Fill(pt, 1.f);
1617 }
1618 }
1619 } else {
1620 float totalWeight = 0.;
1621 for (int32_t j = 0; j < GetMCLabelNID(i); j++) {
1622 mcLabelI_t labelT = GetMCLabel(i, j);
1623 if (!labelT.isFake() && GetMCTrackObj(mMCParam, labelT).pt > GPUCA_MIN_TRACK_PTB5_DEFAULT) {
1624 totalWeight += GetMCLabelWeight(i, j);
1625 }
1626 }
1627 if (totalWeight > 0) {
1628 for (int32_t j = 0; j < GetMCLabelNID(i); j++) {
1629 mcLabelI_t label = GetMCLabel(i, j);
1630 if (!label.isFake() && GetMCTrackObj(mMCParam, label).pt > GPUCA_MIN_TRACK_PTB5_DEFAULT) {
1631 float pt = GetMCTrackObj(mMCParam, label).pt;
1632 if (pt < PT_MIN_CLUST) {
1633 pt = PT_MIN_CLUST;
1634 }
1635 float weight = GetMCLabelWeight(i, j) / totalWeight;
1636 if (mClusterParam[i].fakeAdjacent) {
1637 mClusters[CL_fakeAdj]->Fill(pt, weight);
1638 }
1639 if (mClusterParam[i].fakeAdjacent) {
1640 mClusters[CL_att_adj]->Fill(pt, weight);
1641 }
1642 if (GetMCTrackObj(mRecTracks, label)) {
1643 mClusters[CL_tracks]->Fill(pt, weight);
1644 }
1645 mClusters[CL_all]->Fill(pt, weight);
1646 if (protect || physics) {
1647 mClusters[CL_prot]->Fill(pt, weight);
1648 }
1649 if (physics) {
1650 mClusters[CL_physics]->Fill(pt, weight);
1651 }
1652 }
1653 }
1654 } else {
1655 if (mClusterParam[i].fakeAdjacent) {
1656 mClusters[CL_fakeAdj]->Fill(0.f, 1.f);
1657 }
1658 if (mClusterParam[i].fakeAdjacent) {
1659 mClusters[CL_att_adj]->Fill(0.f, 1.f);
1660 }
1661 mClusters[CL_all]->Fill(0.f, 1.f);
1662 mClusterCounts.nUnaccessible++;
1663 if (protect || physics) {
1664 mClusters[CL_prot]->Fill(0.f, 1.f);
1665 }
1666 if (physics) {
1667 mClusters[CL_physics]->Fill(0.f, 1.f);
1668 }
1669 }
1670 }
1671 }
1672
1673 if (QA_TIMING || (mTracking && mTracking->GetProcessingSettings().debugLevel >= 3)) {
1674 GPUInfo("QA Time: Fill cluster histograms:\t%6.0f us", timer.GetCurrentElapsedTime(true) * 1e6);
1675 }
1676 }
1677 } else if (!mConfig.inputHistogramsOnly && !mConfig.noMC && (mQATasks & (taskTrackingEff | taskTrackingRes | taskTrackingResPull | taskClusterAttach))) {
1678 GPUWarning("No MC information available, only running partial TPC QA!");
1679 }
1680
1681 if (mQATasks & taskTrackStatistics) {
1682 // Fill track statistic histograms
1683 std::vector<std::array<float, 3>> clusterAttachCounts;
1684 if (mcAvail) {
1685 clusterAttachCounts.resize(GetNMCLabels(), {0.f, 0.f});
1686 }
1687 for (uint32_t i = 0; i < nReconstructedTracks; i++) {
1688 const GPUTPCGMMergedTrack& track = mTracking->mIOPtrs.mergedTracks[i];
1689 if (!track.OK()) {
1690 continue;
1691 }
1692 mTracks->Fill(1.f / fabsf(track.GetParam().GetQPt()));
1693 mNCl[0]->Fill(track.NClustersFitted());
1694 uint32_t nClCorrected = 0;
1695 const auto& trackClusters = mTracking->mIOPtrs.mergedTrackHits;
1696 uint32_t jNext = 0;
1697 for (uint32_t j = 0; j < track.NClusters(); j = jNext) {
1698 uint32_t rowClCount = !(trackClusters[track.FirstClusterRef() + j].state & GPUTPCGMMergedTrackHit::flagReject);
1699 for (jNext = j + 1; j < track.NClusters(); jNext++) {
1700 if (trackClusters[track.FirstClusterRef() + j].sector != trackClusters[track.FirstClusterRef() + jNext].sector || trackClusters[track.FirstClusterRef() + j].row != trackClusters[track.FirstClusterRef() + jNext].row) {
1701 break;
1702 }
1703 rowClCount += !(trackClusters[track.FirstClusterRef() + jNext].state & GPUTPCGMMergedTrackHit::flagReject);
1704 }
1705 if (trackClusters[track.FirstClusterRef() + j].leg == trackClusters[track.FirstClusterRef() + track.NClusters() - 1].leg && rowClCount) {
1706 nClCorrected++;
1707 }
1708 if (mcAvail && rowClCount) {
1709 for (uint32_t k = j; k < jNext; k++) {
1710 const auto& cl = trackClusters[track.FirstClusterRef() + k];
1711 if (cl.state & GPUTPCGMMergedTrackHit::flagReject) {
1712 continue;
1713 }
1714 bool labelOk = false, labelOkNonFake = false;
1715 const mcLabelI_t& trkLabel = mTrackMCLabels[i];
1716 if (trkLabel.isValid() && !trkLabel.isNoise()) {
1717 for (int32_t l = 0; l < GetMCLabelNID(cl.num); l++) {
1718 const mcLabelI_t& clLabel = GetMCLabel(cl.num, l);
1719 if (clLabel.isValid() && !clLabel.isNoise() && CompareIgnoreFake(trkLabel, clLabel)) {
1720 labelOk = true;
1721 if (!trkLabel.isFake()) {
1722 labelOkNonFake = true;
1723 }
1724 break;
1725 }
1726 }
1727 }
1728 clusterAttachCounts[cl.num][0] += 1.0f;
1729 clusterAttachCounts[cl.num][1] += (float)labelOk / rowClCount;
1730 clusterAttachCounts[cl.num][2] += (float)labelOkNonFake / rowClCount;
1731 }
1732 }
1733 }
1734 mNCl[1]->Fill(nClCorrected);
1735 }
1736 if (mClNative && mTracking && mTracking->GetTPCTransformHelper()) {
1737 for (uint32_t i = 0; i < GPUChainTracking::NSECTORS; i++) {
1738 for (uint32_t j = 0; j < GPUCA_ROW_COUNT; j++) {
1739 for (uint32_t k = 0; k < mClNative->nClusters[i][j]; k++) {
1740 const auto& cl = mClNative->clusters[i][j][k];
1741 float x, y, z;
1742 GPUTPCConvertImpl::convert(*mTracking->GetTPCTransformHelper()->getCorrMap(), mTracking->GetParam(), i, j, cl.getPad(), cl.getTime(), x, y, z);
1743 mTracking->GetParam().Sector2Global(i, x, y, z, &x, &y, &z);
1744 mClXY->Fill(x, y);
1745 }
1746 }
1747 }
1748 }
1749 if (mcAvail) {
1750 double clusterAttachNormalizedCount = 0, clusterAttachNormalizedCountNonFake = 0;
1751 for (uint32_t i = 0; i < clusterAttachCounts.size(); i++) {
1752 if (clusterAttachCounts[i][0]) {
1753 clusterAttachNormalizedCount += clusterAttachCounts[i][1] / clusterAttachCounts[i][0];
1754 clusterAttachNormalizedCountNonFake += clusterAttachCounts[i][2] / clusterAttachCounts[i][0];
1755 }
1756 }
1757 mClusterCounts.nCorrectlyAttachedNormalized = clusterAttachNormalizedCount;
1758 mClusterCounts.nCorrectlyAttachedNormalizedNonFake = clusterAttachNormalizedCountNonFake;
1759 clusterAttachCounts.clear();
1760 }
1761
1762 if (QA_TIMING || (mTracking && mTracking->GetProcessingSettings().debugLevel >= 3)) {
1763 GPUInfo("QA Time: Fill track statistics:\t%6.0f us", timer.GetCurrentElapsedTime(true) * 1e6);
1764 }
1765 }
1766
1767 uint32_t nCl = clNative ? clNative->nClustersTotal : mTracking->GetProcessors()->tpcMerger.NMaxClusters();
1768 mClusterCounts.nTotal += nCl;
1769 if (mQATasks & taskClusterCounts) {
1770 for (uint32_t i = 0; i < nCl; i++) {
1771 int32_t attach = mTracking->mIOPtrs.mergedTrackHitAttachment[i];
1773
1774 if (mcAvail) {
1775 float totalWeight = 0, weight400 = 0, weight40 = 0;
1776 for (int32_t j = 0; j < GetMCLabelNID(i); j++) {
1777 const auto& label = GetMCLabel(i, j);
1778 if (GetMCLabelID(label) >= 0) {
1779 totalWeight += GetMCLabelWeight(label);
1780 if (GetMCTrackObj(mMCParam, label).pt >= 0.4) {
1781 weight400 += GetMCLabelWeight(label);
1782 }
1783 if (GetMCTrackObj(mMCParam, label).pt <= 0.04) {
1784 weight40 += GetMCLabelWeight(label);
1785 }
1786 }
1787 }
1788 if (totalWeight > 0 && 10.f * weight400 >= totalWeight) {
1789 if (!unattached && !protect && !physics) {
1790 mClusterCounts.nFakeRemove400++;
1791 int32_t totalFake = weight400 < 0.9f * totalWeight;
1792 if (totalFake) {
1793 mClusterCounts.nFullFakeRemove400++;
1794 }
1795 /*printf("Fake removal (%d): Hit %7d, attached %d lowPt %d looper %d tube200 %d highIncl %d tube %d bad %d recPt %7.2f recLabel %6d", totalFake, i, (int32_t) (mClusterParam[i].attached || mClusterParam[i].fakeAttached),
1796 (int32_t) lowPt, (int32_t) ((attach & gputpcgmmergertypes::attachGoodLeg) == 0), (int32_t) ((attach & gputpcgmmergertypes::attachTube) && mev200),
1797 (int32_t) ((attach & gputpcgmmergertypes::attachHighIncl) != 0), (int32_t) ((attach & gputpcgmmergertypes::attachTube) != 0), (int32_t) ((attach & gputpcgmmergertypes::attachGood) == 0),
1798 fabsf(qpt) > 0 ? 1.f / qpt : 0.f, id);
1799 for (int32_t j = 0;j < GetMCLabelNID(i);j++)
1800 {
1801 //if (GetMCLabelID(i, j) < 0) break;
1802 printf(" - label%d %6d weight %5d", j, GetMCLabelID(i, j), (int32_t) GetMCLabelWeight(i, j));
1803 if (GetMCLabelID(i, j) >= 0) printf(" - pt %7.2f", mMCParam[GetMCLabelID(i, j)].pt);
1804 else printf(" ");
1805 }
1806 printf("\n");*/
1807 }
1808 mClusterCounts.nAbove400++;
1809 }
1810 if (totalWeight > 0 && weight40 >= 0.9 * totalWeight) {
1811 mClusterCounts.nBelow40++;
1812 if (protect || physics) {
1813 mClusterCounts.nFakeProtect40++;
1814 }
1815 }
1816 }
1817 if (physics) {
1818 mClusterCounts.nPhysics++;
1819 }
1820 if (physics || protect) {
1821 mClusterCounts.nProt++;
1822 }
1823 if (unattached) {
1824 mClusterCounts.nUnattached++;
1825 }
1826 }
1827 }
1828
1829 // Process cluster count statistics
1830 if ((mQATasks & taskClusterCounts) && mConfig.clusterRejectionHistograms) {
1831 DoClusterCounts(nullptr);
1832 mClusterCounts = counts_t();
1833 }
1834
1835 if (QA_TIMING || (mTracking && mTracking->GetProcessingSettings().debugLevel >= 3)) {
1836 GPUInfo("QA Time: Cluster Counts:\t%6.0f us", timer.GetCurrentElapsedTime(true) * 1e6);
1837 }
1838
1839 if (mConfig.dumpToROOT) {
1840 if (!clNative || !mTracking || !mTracking->mIOPtrs.mergedTrackHitAttachment || !mTracking->mIOPtrs.mergedTracks) {
1841 throw std::runtime_error("Cannot dump non o2::tpc::clusterNative clusters, need also hit attachmend and GPU tracks");
1842 }
1843 uint32_t clid = 0;
1844 for (uint32_t i = 0; i < GPUChainTracking::NSECTORS; i++) {
1845 for (uint32_t j = 0; j < GPUCA_ROW_COUNT; j++) {
1846 for (uint32_t k = 0; k < mClNative->nClusters[i][j]; k++) {
1847 const auto& cl = mClNative->clusters[i][j][k];
1848 uint32_t attach = mTracking->mIOPtrs.mergedTrackHitAttachment[clid];
1849 float x = 0, y = 0, z = 0;
1851 uint32_t track = attach & gputpcgmmergertypes::attachTrackMask;
1852 const auto& trk = mTracking->mIOPtrs.mergedTracks[track];
1853 mTracking->GetTPCTransformHelper()->Transform(i, j, cl.getPad(), cl.getTime(), x, y, z, trk.GetParam().GetTZOffset());
1854 mTracking->GetParam().Sector2Global(i, x, y, z, &x, &y, &z);
1855 }
1856 uint32_t extState = mTracking->mIOPtrs.mergedTrackHitStates ? mTracking->mIOPtrs.mergedTrackHitStates[clid] : 0;
1857
1858 if (mConfig.dumpToROOT >= 2) {
1861 memset((void*)&trk, 0, sizeof(trk));
1862 memset((void*)&trkHit, 0, sizeof(trkHit));
1864 uint32_t track = attach & gputpcgmmergertypes::attachTrackMask;
1865 trk = mTracking->mIOPtrs.mergedTracks[track];
1866 for (uint32_t l = 0; l < trk.NClusters(); l++) {
1867 const auto& tmp = mTracking->mIOPtrs.mergedTrackHits[trk.FirstClusterRef() + l];
1868 if (tmp.num == clid) {
1869 trkHit = tmp;
1870 break;
1871 }
1872 }
1873 }
1874 static auto cldump = GPUROOTDump<o2::tpc::ClusterNative, GPUTPCGMMergedTrack, GPUTPCGMMergedTrackHit, uint32_t, uint32_t, float, float, float, uint32_t, uint32_t, uint32_t>::getNew("cluster", "track", "trackHit", "attach", "extState", "x", "y", "z", "sector", "row", "nEv", "clusterTree");
1875 cldump.Fill(cl, trk, trkHit, attach, extState, x, y, z, i, j, mNEvents - 1);
1876 } else {
1877 static auto cldump = GPUROOTDump<o2::tpc::ClusterNative, uint32_t, uint32_t, float, float, float, uint32_t, uint32_t, uint32_t>::getNew("cluster", "attach", "extState", "x", "y", "z", "sector", "row", "nEv", "clusterTree");
1878 cldump.Fill(cl, attach, extState, x, y, z, i, j, mNEvents - 1);
1879 }
1880 clid++;
1881 }
1882 }
1883 }
1884
1885 static auto trkdump = GPUROOTDump<uint32_t, GPUTPCGMMergedTrack>::getNew("nEv", "track", "tracksTree");
1886 for (uint32_t i = 0; i < mTracking->mIOPtrs.nMergedTracks; i++) {
1887 if (mTracking->mIOPtrs.mergedTracks[i].OK()) {
1888 trkdump.Fill(mNEvents - 1, mTracking->mIOPtrs.mergedTracks[i]);
1889 }
1890 }
1891
1892 if (mTracking && mTracking->GetProcessingSettings().createO2Output) {
1893 static auto o2trkdump = GPUROOTDump<uint32_t, o2::tpc::TrackTPC>::getNew("nEv", "track", "tracksO2Tree");
1894 for (uint32_t i = 0; i < mTracking->mIOPtrs.nOutputTracksTPCO2; i++) {
1895 o2trkdump.Fill(mNEvents - 1, mTracking->mIOPtrs.outputTracksTPCO2[i]);
1896 }
1897 }
1898 }
1899 mTrackingScratchBuffer.clear();
1900 mTrackingScratchBuffer.shrink_to_fit();
1901}
1902
1903void GPUQA::GetName(char* fname, int32_t k)
1904{
1905 const int32_t nNewInput = mConfig.inputHistogramsOnly ? 0 : 1;
1906 if (k || mConfig.inputHistogramsOnly || mConfig.name.size()) {
1907 if (!(mConfig.inputHistogramsOnly || k)) {
1908 snprintf(fname, 1024, "%s - ", mConfig.name.c_str());
1909 } else if (mConfig.compareInputNames.size() > (unsigned)(k - nNewInput)) {
1910 snprintf(fname, 1024, "%s - ", mConfig.compareInputNames[k - nNewInput].c_str());
1911 } else {
1912 strcpy(fname, mConfig.compareInputs[k - nNewInput].c_str());
1913 if (strlen(fname) > 5 && strcmp(fname + strlen(fname) - 5, ".root") == 0) {
1914 fname[strlen(fname) - 5] = 0;
1915 }
1916 strcat(fname, " - ");
1917 }
1918 } else {
1919 fname[0] = 0;
1920 }
1921}
1922
1923template <class T>
1924T* GPUQA::GetHist(T*& ee, std::vector<std::unique_ptr<TFile>>& tin, int32_t k, int32_t nNewInput)
1925{
1926 T* e = ee;
1927 if ((mConfig.inputHistogramsOnly || k) && (e = dynamic_cast<T*>(tin[k - nNewInput]->Get(e->GetName()))) == nullptr) {
1928 GPUWarning("Missing histogram in input %s: %s", mConfig.compareInputs[k - nNewInput].c_str(), ee->GetName());
1929 return (nullptr);
1930 }
1931 ee = e;
1932 return (e);
1933}
1934
1935void GPUQA::DrawQAHistogramsCleanup()
1936{
1937 clearGarbagageCollector();
1938}
1939
1940void GPUQA::resetHists()
1941{
1942 if (!mQAInitialized) {
1943 throw std::runtime_error("QA not initialized");
1944 }
1945 if (mHaveExternalHists) {
1946 throw std::runtime_error("Cannot reset external hists");
1947 }
1948 for (auto& h : *mHist1D) {
1949 h.Reset();
1950 }
1951 for (auto& h : *mHist2D) {
1952 h.Reset();
1953 }
1954 for (auto& h : *mHist1Dd) {
1955 h.Reset();
1956 }
1957 for (auto& h : *mHistGraph) {
1958 h = TGraphAsymmErrors();
1959 }
1960 mClusterCounts = counts_t();
1961}
1962
1963int32_t GPUQA::DrawQAHistograms(TObjArray* qcout)
1964{
1965 const auto oldRootIgnoreLevel = gErrorIgnoreLevel;
1966 gErrorIgnoreLevel = kWarning;
1967 if (!mQAInitialized) {
1968 throw std::runtime_error("QA not initialized");
1969 }
1970
1971 if (mTracking && mTracking->GetProcessingSettings().debugLevel >= 2) {
1972 printf("Creating QA Histograms\n");
1973 }
1974
1975 std::vector<Color_t> colorNums(COLORCOUNT);
1976 if (!qcout) {
1977 static int32_t initColorsInitialized = initColors();
1978 (void)initColorsInitialized;
1979 }
1980 for (int32_t i = 0; i < COLORCOUNT; i++) {
1981 colorNums[i] = qcout ? defaultColorNums[i] : mColors[i]->GetNumber();
1982 }
1983
1984 bool mcAvail = mcPresent();
1985 char name[2048], fname[1024];
1986
1987 const int32_t nNewInput = mConfig.inputHistogramsOnly ? 0 : 1;
1988 const int32_t ConfigNumInputs = nNewInput + mConfig.compareInputs.size();
1989
1990 std::vector<std::unique_ptr<TFile>> tin;
1991 for (uint32_t i = 0; i < mConfig.compareInputs.size(); i++) {
1992 tin.emplace_back(std::make_unique<TFile>(mConfig.compareInputs[i].c_str()));
1993 }
1994 std::unique_ptr<TFile> tout = nullptr;
1995 if (mConfig.output.size()) {
1996 tout = std::make_unique<TFile>(mConfig.output.c_str(), "RECREATE");
1997 }
1998
1999 if (mConfig.enableLocalOutput || mConfig.shipToQCAsCanvas) {
2000 float legendSpacingString = 0.025;
2001 for (int32_t i = 0; i < ConfigNumInputs; i++) {
2002 GetName(fname, i);
2003 if (strlen(fname) * 0.006 > legendSpacingString) {
2004 legendSpacingString = strlen(fname) * 0.006;
2005 }
2006 }
2007
2008 // Create Canvas / Pads for Efficiency Histograms
2009 if (mQATasks & taskTrackingEff) {
2010 for (int32_t ii = 0; ii < 6; ii++) {
2011 int32_t i = ii == 5 ? 4 : ii;
2012 snprintf(fname, 1024, "eff_vs_%s_layout", VSPARAMETER_NAMES[ii]);
2013 snprintf(name, 2048, "Efficiency versus %s", VSPARAMETER_NAMES[i]);
2014 mCEff[ii] = createGarbageCollected<TCanvas>(fname, name, 0, 0, 700, 700. * 2. / 3.);
2015 mCEff[ii]->cd();
2016 float dy = 1. / 2.;
2017 mPEff[ii][0] = createGarbageCollected<TPad>("p0", "", 0.0, dy * 0, 0.5, dy * 1);
2018 mPEff[ii][0]->Draw();
2019 mPEff[ii][0]->SetRightMargin(0.04);
2020 mPEff[ii][1] = createGarbageCollected<TPad>("p1", "", 0.5, dy * 0, 1.0, dy * 1);
2021 mPEff[ii][1]->Draw();
2022 mPEff[ii][1]->SetRightMargin(0.04);
2023 mPEff[ii][2] = createGarbageCollected<TPad>("p2", "", 0.0, dy * 1, 0.5, dy * 2 - .001);
2024 mPEff[ii][2]->Draw();
2025 mPEff[ii][2]->SetRightMargin(0.04);
2026 mPEff[ii][3] = createGarbageCollected<TPad>("p3", "", 0.5, dy * 1, 1.0, dy * 2 - .001);
2027 mPEff[ii][3]->Draw();
2028 mPEff[ii][3]->SetRightMargin(0.04);
2029 mLEff[ii] = createGarbageCollected<TLegend>(0.92 - legendSpacingString * 1.45, 0.83 - (0.93 - 0.82) / 2. * (float)ConfigNumInputs, 0.98, 0.849);
2030 SetLegend(mLEff[ii]);
2031 }
2032 }
2033
2034 // Create Canvas / Pads for Resolution Histograms
2035 if (mQATasks & taskTrackingRes) {
2036 for (int32_t ii = 0; ii < 7; ii++) {
2037 int32_t i = ii == 5 ? 4 : ii;
2038 if (ii == 6) {
2039 snprintf(fname, 1024, "res_integral_layout");
2040 snprintf(name, 2048, "Integral Resolution");
2041 } else {
2042 snprintf(fname, 1024, "res_vs_%s_layout", VSPARAMETER_NAMES[ii]);
2043 snprintf(name, 2048, "Resolution versus %s", VSPARAMETER_NAMES[i]);
2044 }
2045 mCRes[ii] = createGarbageCollected<TCanvas>(fname, name, 0, 0, 700, 700. * 2. / 3.);
2046 mCRes[ii]->cd();
2047 gStyle->SetOptFit(1);
2048
2049 float dy = 1. / 2.;
2050 mPRes[ii][3] = createGarbageCollected<TPad>("p0", "", 0.0, dy * 0, 0.5, dy * 1);
2051 mPRes[ii][3]->Draw();
2052 mPRes[ii][3]->SetRightMargin(0.04);
2053 mPRes[ii][4] = createGarbageCollected<TPad>("p1", "", 0.5, dy * 0, 1.0, dy * 1);
2054 mPRes[ii][4]->Draw();
2055 mPRes[ii][4]->SetRightMargin(0.04);
2056 mPRes[ii][0] = createGarbageCollected<TPad>("p2", "", 0.0, dy * 1, 1. / 3., dy * 2 - .001);
2057 mPRes[ii][0]->Draw();
2058 mPRes[ii][0]->SetRightMargin(0.04);
2059 mPRes[ii][0]->SetLeftMargin(0.15);
2060 mPRes[ii][1] = createGarbageCollected<TPad>("p3", "", 1. / 3., dy * 1, 2. / 3., dy * 2 - .001);
2061 mPRes[ii][1]->Draw();
2062 mPRes[ii][1]->SetRightMargin(0.04);
2063 mPRes[ii][1]->SetLeftMargin(0.135);
2064 mPRes[ii][2] = createGarbageCollected<TPad>("p4", "", 2. / 3., dy * 1, 1.0, dy * 2 - .001);
2065 mPRes[ii][2]->Draw();
2066 mPRes[ii][2]->SetRightMargin(0.06);
2067 mPRes[ii][2]->SetLeftMargin(0.135);
2068 if (ii < 6) {
2069 mLRes[ii] = createGarbageCollected<TLegend>(0.9 - legendSpacingString * 1.45, 0.93 - (0.93 - 0.86) / 2. * (float)ConfigNumInputs, 0.98, 0.949);
2070 SetLegend(mLRes[ii]);
2071 }
2072 }
2073 }
2074
2075 // Create Canvas / Pads for Pull Histograms
2076 if (mQATasks & taskTrackingResPull) {
2077 for (int32_t ii = 0; ii < 7; ii++) {
2078 int32_t i = ii == 5 ? 4 : ii;
2079
2080 if (ii == 6) {
2081 snprintf(fname, 1024, "pull_integral_layout");
2082 snprintf(name, 2048, "Integral Pull");
2083 } else {
2084 snprintf(fname, 1024, "pull_vs_%s_layout", VSPARAMETER_NAMES[ii]);
2085 snprintf(name, 2048, "Pull versus %s", VSPARAMETER_NAMES[i]);
2086 }
2087 mCPull[ii] = createGarbageCollected<TCanvas>(fname, name, 0, 0, 700, 700. * 2. / 3.);
2088 mCPull[ii]->cd();
2089 gStyle->SetOptFit(1);
2090
2091 float dy = 1. / 2.;
2092 mPPull[ii][3] = createGarbageCollected<TPad>("p0", "", 0.0, dy * 0, 0.5, dy * 1);
2093 mPPull[ii][3]->Draw();
2094 mPPull[ii][3]->SetRightMargin(0.04);
2095 mPPull[ii][4] = createGarbageCollected<TPad>("p1", "", 0.5, dy * 0, 1.0, dy * 1);
2096 mPPull[ii][4]->Draw();
2097 mPPull[ii][4]->SetRightMargin(0.04);
2098 mPPull[ii][0] = createGarbageCollected<TPad>("p2", "", 0.0, dy * 1, 1. / 3., dy * 2 - .001);
2099 mPPull[ii][0]->Draw();
2100 mPPull[ii][0]->SetRightMargin(0.04);
2101 mPPull[ii][0]->SetLeftMargin(0.15);
2102 mPPull[ii][1] = createGarbageCollected<TPad>("p3", "", 1. / 3., dy * 1, 2. / 3., dy * 2 - .001);
2103 mPPull[ii][1]->Draw();
2104 mPPull[ii][1]->SetRightMargin(0.04);
2105 mPPull[ii][1]->SetLeftMargin(0.135);
2106 mPPull[ii][2] = createGarbageCollected<TPad>("p4", "", 2. / 3., dy * 1, 1.0, dy * 2 - .001);
2107 mPPull[ii][2]->Draw();
2108 mPPull[ii][2]->SetRightMargin(0.06);
2109 mPPull[ii][2]->SetLeftMargin(0.135);
2110 if (ii < 6) {
2111 mLPull[ii] = createGarbageCollected<TLegend>(0.9 - legendSpacingString * 1.45, 0.93 - (0.93 - 0.86) / 2. * (float)ConfigNumInputs, 0.98, 0.949);
2112 SetLegend(mLPull[ii]);
2113 }
2114 }
2115 }
2116
2117 // Create Canvas for Cluster Histos
2118 if (mQATasks & taskClusterAttach) {
2119 for (int32_t i = 0; i < 3; i++) {
2120 snprintf(fname, 1024, "clusters_%s_layout", CLUSTER_TYPES[i]);
2121 mCClust[i] = createGarbageCollected<TCanvas>(fname, CLUSTER_TITLES[i], 0, 0, 700, 700. * 2. / 3.);
2122 mCClust[i]->cd();
2123 mPClust[i] = createGarbageCollected<TPad>("p0", "", 0.0, 0.0, 1.0, 1.0);
2124 mPClust[i]->Draw();
2125 float y1 = i != 1 ? 0.77 : 0.27, y2 = i != 1 ? 0.9 : 0.42;
2126 mLClust[i] = createGarbageCollected<TLegend>(i == 2 ? 0.1 : (0.65 - legendSpacingString * 1.45), y2 - (y2 - y1) * (ConfigNumInputs + (i != 1) / 2.) + 0.005, i == 2 ? (0.3 + legendSpacingString * 1.45) : 0.9, y2);
2127 SetLegend(mLClust[i]);
2128 }
2129 }
2130
2131 // Create Canvas for track statistic histos
2132 if (mQATasks & taskTrackStatistics) {
2133 mCTracks = createGarbageCollected<TCanvas>("ctracks", "Track Pt", 0, 0, 700, 700. * 2. / 3.);
2134 mCTracks->cd();
2135 mPTracks = createGarbageCollected<TPad>("p0", "", 0.0, 0.0, 1.0, 1.0);
2136 mPTracks->Draw();
2137 mLTracks = createGarbageCollected<TLegend>(0.9 - legendSpacingString * 1.45, 0.93 - (0.93 - 0.86) / 2. * (float)ConfigNumInputs, 0.98, 0.949);
2138 SetLegend(mLTracks);
2139
2140 for (int32_t i = 0; i < 2; i++) {
2141 snprintf(name, 2048, "cncl%d Pull", i);
2142 mCNCl[i] = createGarbageCollected<TCanvas>(name, i ? "Number of clusters (corrected for multiple per row)" : "Number of clusters per track", 0, 0, 700, 700. * 2. / 3.);
2143 mCNCl[i]->cd();
2144 mPNCl[i] = createGarbageCollected<TPad>("p0", "", 0.0, 0.0, 1.0, 1.0);
2145 mPNCl[i]->Draw();
2146 mLNCl[i] = createGarbageCollected<TLegend>(0.9 - legendSpacingString * 1.45, 0.93 - (0.93 - 0.86) / 2. * (float)ConfigNumInputs, 0.98, 0.949);
2147 SetLegend(mLNCl[i]);
2148 }
2149
2150 mCClXY = createGarbageCollected<TCanvas>("clxy", "Number of clusters per X / Y", 0, 0, 700, 700. * 2. / 3.);
2151 mCClXY->cd();
2152 mPClXY = createGarbageCollected<TPad>("p0", "", 0.0, 0.0, 1.0, 1.0);
2153 mPClXY->Draw();
2154 }
2155 }
2156
2157 if (mConfig.enableLocalOutput && !mConfig.inputHistogramsOnly && (mQATasks & taskTrackingEff) && mcPresent()) {
2158 GPUInfo("QA Stats: Eff: Tracks Prim %d (Eta %d, Pt %d) %f%% (%f%%) Sec %d (Eta %d, Pt %d) %f%% (%f%%) - Res: Tracks %d (Eta %d, Pt %d)", (int32_t)mEff[3][1][0][0]->GetEntries(), (int32_t)mEff[3][1][0][3]->GetEntries(), (int32_t)mEff[3][1][0][4]->GetEntries(),
2159 mEff[0][0][0][0]->GetSumOfWeights() / std::max(1., mEff[3][0][0][0]->GetSumOfWeights()), mEff[0][1][0][0]->GetSumOfWeights() / std::max(1., mEff[3][1][0][0]->GetSumOfWeights()), (int32_t)mEff[3][1][1][0]->GetEntries(), (int32_t)mEff[3][1][1][3]->GetEntries(),
2160 (int32_t)mEff[3][1][1][4]->GetEntries(), mEff[0][0][1][0]->GetSumOfWeights() / std::max(1., mEff[3][0][1][0]->GetSumOfWeights()), mEff[0][1][1][0]->GetSumOfWeights() / std::max(1., mEff[3][1][1][0]->GetSumOfWeights()), (int32_t)mRes2[0][0]->GetEntries(),
2161 (int32_t)mRes2[0][3]->GetEntries(), (int32_t)mRes2[0][4]->GetEntries());
2162 }
2163
2164 int32_t flagShowVsPtLog = (mConfig.enableLocalOutput || mConfig.shipToQCAsCanvas) ? 1 : 0;
2165
2166 if (mQATasks & taskTrackingEff) {
2167 // Process / Draw Efficiency Histograms
2168 for (int32_t ii = 0; ii < 5 + flagShowVsPtLog; ii++) {
2169 int32_t i = ii == 5 ? 4 : ii;
2170 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2171 for (int32_t j = 0; j < 4; j++) {
2172 if (mConfig.enableLocalOutput || mConfig.shipToQCAsCanvas) {
2173 mPEff[ii][j]->cd();
2174 if (ii == 5) {
2175 mPEff[ii][j]->SetLogx();
2176 }
2177 }
2178 for (int32_t l = 0; l < 3; l++) {
2179 if (k == 0 && mConfig.inputHistogramsOnly == 0 && ii != 5) {
2180 if (l == 0) {
2181 // Divide eff, compute all for fake/clone
2182 auto oldLevel = gErrorIgnoreLevel;
2183 gErrorIgnoreLevel = kError;
2184 mEffResult[0][j / 2][j % 2][i]->Divide(mEff[l][j / 2][j % 2][i], mEff[3][j / 2][j % 2][i], "cl=0.683 b(1,1) mode");
2185 gErrorIgnoreLevel = oldLevel;
2186 mEff[3][j / 2][j % 2][i]->Reset(); // Sum up rec + clone + fake for fake rate
2187 mEff[3][j / 2][j % 2][i]->Add(mEff[0][j / 2][j % 2][i]);
2188 mEff[3][j / 2][j % 2][i]->Add(mEff[1][j / 2][j % 2][i]);
2189 mEff[3][j / 2][j % 2][i]->Add(mEff[2][j / 2][j % 2][i]);
2190 mEff[4][j / 2][j % 2][i]->Reset(); // Sum up rec + clone for clone rate
2191 mEff[4][j / 2][j % 2][i]->Add(mEff[0][j / 2][j % 2][i]);
2192 mEff[4][j / 2][j % 2][i]->Add(mEff[1][j / 2][j % 2][i]);
2193 } else {
2194 // Divide fake/clone
2195 auto oldLevel = gErrorIgnoreLevel;
2196 gErrorIgnoreLevel = kError;
2197 mEffResult[l][j / 2][j % 2][i]->Divide(mEff[l][j / 2][j % 2][i], mEff[l == 1 ? 4 : 3][j / 2][j % 2][i], "cl=0.683 b(1,1) mode");
2198 gErrorIgnoreLevel = oldLevel;
2199 }
2200 }
2201
2202 TGraphAsymmErrors* e = mEffResult[l][j / 2][j % 2][i];
2203
2204 if (!mConfig.inputHistogramsOnly && k == 0) {
2205 if (tout) {
2206 mEff[l][j / 2][j % 2][i]->Write();
2207 e->Write();
2208 if (l == 2) {
2209 mEff[3][j / 2][j % 2][i]->Write(); // Store also all histogram!
2210 mEff[4][j / 2][j % 2][i]->Write(); // Store also all histogram!
2211 }
2212 }
2213 } else if (GetHist(e, tin, k, nNewInput) == nullptr) {
2214 continue;
2215 }
2216 e->SetTitle(EFFICIENCY_TITLES[j]);
2217 e->GetYaxis()->SetTitle("(Efficiency)");
2218 e->GetXaxis()->SetTitle(XAXIS_TITLES[i]);
2219
2220 e->SetLineWidth(1);
2221 e->SetLineStyle(CONFIG_DASHED_MARKERS ? k + 1 : 1);
2222 SetAxisSize(e);
2223 if (qcout && !mConfig.shipToQCAsCanvas) {
2224 qcout->Add(e);
2225 }
2226 if (!mConfig.enableLocalOutput && !mConfig.shipToQCAsCanvas) {
2227 continue;
2228 }
2229 e->SetMarkerColor(kBlack);
2230 e->SetLineColor(colorNums[(l == 2 ? (ConfigNumInputs * 2 + k) : (k * 2 + l)) % COLORCOUNT]);
2231 e->GetHistogram()->GetYaxis()->SetRangeUser(-0.02, 1.02);
2232 e->Draw(k || l ? "same P" : "AP");
2233 if (j == 0) {
2234 GetName(fname, k);
2235 snprintf(name, 2048, "%s%s", fname, EFF_NAMES[l]);
2236 mLEff[ii]->AddEntry(e, name, "l");
2237 }
2238 }
2239 if (!mConfig.enableLocalOutput && !mConfig.shipToQCAsCanvas) {
2240 continue;
2241 }
2242 mCEff[ii]->cd();
2243 ChangePadTitleSize(mPEff[ii][j], 0.056);
2244 }
2245 }
2246 if (!mConfig.enableLocalOutput && !mConfig.shipToQCAsCanvas) {
2247 continue;
2248 }
2249
2250 mLEff[ii]->Draw();
2251
2252 if (qcout) {
2253 qcout->Add(mCEff[ii]);
2254 }
2255 if (!mConfig.enableLocalOutput) {
2256 continue;
2257 }
2258 doPerfFigure(0.2, 0.295, 0.025);
2259 mCEff[ii]->Print(Form("plots/eff_vs_%s.pdf", VSPARAMETER_NAMES[ii]));
2260 if (mConfig.writeRootFiles) {
2261 mCEff[ii]->Print(Form("plots/eff_vs_%s.root", VSPARAMETER_NAMES[ii]));
2262 }
2263 }
2264 }
2265
2266 if (mQATasks & (taskTrackingRes | taskTrackingResPull)) {
2267 // Process / Draw Resolution Histograms
2268 TH1D *resIntegral[5] = {}, *pullIntegral[5] = {};
2269 TCanvas* cfit = nullptr;
2270 std::unique_ptr<TF1> customGaus = std::make_unique<TF1>("G", "[0]*exp(-(x-[1])*(x-[1])/(2.*[2]*[2]))");
2271 for (int32_t p = 0; p < 2; p++) {
2272 if ((p == 0 && (mQATasks & taskTrackingRes) == 0) || (p == 1 && (mQATasks & taskTrackingResPull) == 0)) {
2273 continue;
2274 }
2275 for (int32_t ii = 0; ii < 5 + flagShowVsPtLog; ii++) {
2276 TCanvas* can = p ? mCPull[ii] : mCRes[ii];
2277 TLegend* leg = p ? mLPull[ii] : mLRes[ii];
2278 int32_t i = ii == 5 ? 4 : ii;
2279 for (int32_t j = 0; j < 5; j++) {
2280 TH2F* src = p ? mPull2[j][i] : mRes2[j][i];
2281 TH1F** dst = p ? mPull[j][i] : mRes[j][i];
2282 TH1D*& dstIntegral = p ? pullIntegral[j] : resIntegral[j];
2283 TPad* pad = p ? mPPull[ii][j] : mPRes[ii][j];
2284
2285 if (!mConfig.inputHistogramsOnly && ii != 5) {
2286 if (cfit == nullptr) {
2287 cfit = createGarbageCollected<TCanvas>();
2288 }
2289 cfit->cd();
2290
2291 TAxis* axis = src->GetYaxis();
2292 int32_t nBins = axis->GetNbins();
2293 int32_t integ = 1;
2294 for (int32_t bin = 1; bin <= nBins; bin++) {
2295 int32_t bin0 = std::max(bin - integ, 0);
2296 int32_t bin1 = std::min(bin + integ, nBins);
2297 std::unique_ptr<TH1D> proj{src->ProjectionX("proj", bin0, bin1)};
2298 proj->ClearUnderflowAndOverflow();
2299 if (proj->GetEntries()) {
2300 uint32_t rebin = 1;
2301 while (proj->GetMaximum() < 50 && rebin < sizeof(RES_AXIS_BINS) / sizeof(RES_AXIS_BINS[0])) {
2302 proj->Rebin(RES_AXIS_BINS[rebin - 1] / RES_AXIS_BINS[rebin]);
2303 rebin++;
2304 }
2305
2306 if (proj->GetEntries() < 20 || proj->GetRMS() < 0.00001) {
2307 dst[0]->SetBinContent(bin, proj->GetRMS());
2308 dst[0]->SetBinError(bin, std::sqrt(proj->GetRMS()));
2309 dst[1]->SetBinContent(bin, proj->GetMean());
2310 dst[1]->SetBinError(bin, std::sqrt(proj->GetRMS()));
2311 } else {
2312 proj->GetXaxis()->SetRange(0, 0);
2313 proj->GetXaxis()->SetRangeUser(std::max(proj->GetXaxis()->GetXmin(), proj->GetMean() - 3. * proj->GetRMS()), std::min(proj->GetXaxis()->GetXmax(), proj->GetMean() + 3. * proj->GetRMS()));
2314 bool forceLogLike = proj->GetMaximum() < 20;
2315 for (int32_t k = forceLogLike ? 2 : 0; k < 3; k++) {
2316 proj->Fit("gaus", forceLogLike || k == 2 ? "sQl" : k ? "sQww" : "sQ");
2317 TF1* fitFunc = proj->GetFunction("gaus");
2318
2319 if (k && !forceLogLike) {
2320 customGaus->SetParameters(fitFunc->GetParameter(0), fitFunc->GetParameter(1), fitFunc->GetParameter(2));
2321 proj->Fit(customGaus.get(), "sQ");
2322 fitFunc = customGaus.get();
2323 }
2324
2325 const float sigma = fabs(fitFunc->GetParameter(2));
2326 dst[0]->SetBinContent(bin, sigma);
2327 dst[1]->SetBinContent(bin, fitFunc->GetParameter(1));
2328 dst[0]->SetBinError(bin, fitFunc->GetParError(2));
2329 dst[1]->SetBinError(bin, fitFunc->GetParError(1));
2330
2331 const bool fail1 = sigma <= 0.f;
2332 const bool fail2 = fabs(proj->GetMean() - dst[1]->GetBinContent(bin)) > std::min<float>(p ? PULL_AXIS : mConfig.nativeFitResolutions ? RES_AXES_NATIVE[j] : RES_AXES[j], 3.f * proj->GetRMS());
2333 const bool fail3 = dst[0]->GetBinContent(bin) > 3.f * proj->GetRMS() || dst[0]->GetBinError(bin) > 1 || dst[1]->GetBinError(bin) > 1;
2334 const bool fail4 = fitFunc->GetParameter(0) < proj->GetMaximum() / 5.;
2335 const bool fail = fail1 || fail2 || fail3 || fail4;
2336 // if (p == 0 && ii == 4 && j == 2) DrawHisto(proj, Form("Hist_bin_%d-%d_vs_%d____%d_%d___%f-%f___%f-%f___%d.pdf", p, j, ii, bin, k, dst[0]->GetBinContent(bin), proj->GetRMS(), dst[1]->GetBinContent(bin), proj->GetMean(), (int32_t) fail), "");
2337
2338 if (!fail) {
2339 break;
2340 } else if (k >= 2) {
2341 dst[0]->SetBinContent(bin, proj->GetRMS());
2342 dst[0]->SetBinError(bin, std::sqrt(proj->GetRMS()));
2343 dst[1]->SetBinContent(bin, proj->GetMean());
2344 dst[1]->SetBinError(bin, std::sqrt(proj->GetRMS()));
2345 }
2346 }
2347 }
2348 } else {
2349 dst[0]->SetBinContent(bin, 0.f);
2350 dst[0]->SetBinError(bin, 0.f);
2351 dst[1]->SetBinContent(bin, 0.f);
2352 dst[1]->SetBinError(bin, 0.f);
2353 }
2354 }
2355 if (ii == 0) {
2356 dstIntegral = src->ProjectionX(mConfig.nativeFitResolutions ? PARAMETER_NAMES_NATIVE[j] : PARAMETER_NAMES[j], 0, nBins + 1);
2357 uint32_t rebin = 1;
2358 while (dstIntegral->GetMaximum() < 50 && rebin < sizeof(RES_AXIS_BINS) / sizeof(RES_AXIS_BINS[0])) {
2359 dstIntegral->Rebin(RES_AXIS_BINS[rebin - 1] / RES_AXIS_BINS[rebin]);
2360 rebin++;
2361 }
2362 }
2363 }
2364 if (ii == 0) {
2365 if (mConfig.inputHistogramsOnly) {
2366 dstIntegral = createGarbageCollected<TH1D>();
2367 }
2368 snprintf(fname, 1024, p ? "IntPull%s" : "IntRes%s", VSPARAMETER_NAMES[j]);
2369 snprintf(name, 2048, p ? "%s Pull" : "%s Resolution", p || mConfig.nativeFitResolutions ? PARAMETER_NAMES_NATIVE[j] : PARAMETER_NAMES[j]);
2370 dstIntegral->SetName(fname);
2371 dstIntegral->SetTitle(name);
2372 }
2373 if (mConfig.enableLocalOutput || mConfig.shipToQCAsCanvas) {
2374 pad->cd();
2375 }
2376 int32_t numColor = 0;
2377 float tmpMax = -1000.;
2378 float tmpMin = 1000.;
2379
2380 for (int32_t l = 0; l < 2; l++) {
2381 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2382 TH1F* e = dst[l];
2383 if (GetHist(e, tin, k, nNewInput) == nullptr) {
2384 continue;
2385 }
2386 if (nNewInput && k == 0 && ii != 5) {
2387 if (p == 0) {
2388 e->Scale(mConfig.nativeFitResolutions ? SCALE_NATIVE[j] : SCALE[j]);
2389 }
2390 }
2391 if (ii == 4) {
2392 e->GetXaxis()->SetRangeUser(0.2, PT_MAX);
2393 } else if (LOG_PT_MIN > 0 && ii == 5) {
2394 e->GetXaxis()->SetRangeUser(LOG_PT_MIN, PT_MAX);
2395 } else if (ii == 5) {
2396 e->GetXaxis()->SetRange(1, 0);
2397 }
2398 e->SetMinimum(-1111);
2399 e->SetMaximum(-1111);
2400
2401 if (e->GetMaximum() > tmpMax) {
2402 tmpMax = e->GetMaximum();
2403 }
2404 if (e->GetMinimum() < tmpMin) {
2405 tmpMin = e->GetMinimum();
2406 }
2407 }
2408 }
2409
2410 float tmpSpan;
2411 tmpSpan = tmpMax - tmpMin;
2412 tmpMax += tmpSpan * .02;
2413 tmpMin -= tmpSpan * .02;
2414 if (j == 2 && i < 3) {
2415 tmpMax += tmpSpan * 0.13 * ConfigNumInputs;
2416 }
2417
2418 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2419 for (int32_t l = 0; l < 2; l++) {
2420 TH1F* e = dst[l];
2421 if (!mConfig.inputHistogramsOnly && k == 0) {
2422 snprintf(name, 2048, p ? "%s Pull" : "%s Resolution", p || mConfig.nativeFitResolutions ? PARAMETER_NAMES_NATIVE[j] : PARAMETER_NAMES[j]);
2423 e->SetTitle(name);
2424 e->SetStats(kFALSE);
2425 if (tout) {
2426 if (l == 0) {
2427 mRes2[j][i]->SetOption("colz");
2428 mRes2[j][i]->Write();
2429 }
2430 e->Write();
2431 }
2432 } else if (GetHist(e, tin, k, nNewInput) == nullptr) {
2433 continue;
2434 }
2435 e->SetMaximum(tmpMax);
2436 e->SetMinimum(tmpMin);
2437 e->SetLineWidth(1);
2438 e->SetLineStyle(CONFIG_DASHED_MARKERS ? k + 1 : 1);
2439 SetAxisSize(e);
2440 e->GetYaxis()->SetTitle(p ? AXIS_TITLES_PULL[j] : mConfig.nativeFitResolutions ? AXIS_TITLES_NATIVE[j] : AXIS_TITLES[j]);
2441 e->GetXaxis()->SetTitle(XAXIS_TITLES[i]);
2442 if (LOG_PT_MIN > 0 && ii == 5) {
2443 e->GetXaxis()->SetRangeUser(LOG_PT_MIN, PT_MAX);
2444 }
2445
2446 if (j == 0) {
2447 e->GetYaxis()->SetTitleOffset(1.5);
2448 } else if (j < 3) {
2449 e->GetYaxis()->SetTitleOffset(1.4);
2450 }
2451 if (qcout && !mConfig.shipToQCAsCanvas) {
2452 qcout->Add(e);
2453 }
2454 if (!mConfig.enableLocalOutput && !mConfig.shipToQCAsCanvas) {
2455 continue;
2456 }
2457
2458 e->SetMarkerColor(kBlack);
2459 e->SetLineColor(colorNums[numColor++ % COLORCOUNT]);
2460 e->Draw(k || l ? "same" : "");
2461 if (j == 0) {
2462 GetName(fname, k);
2463 if (p) {
2464 snprintf(name, 2048, "%s%s", fname, l ? "Mean" : "Pull");
2465 } else {
2466 snprintf(name, 2048, "%s%s", fname, l ? "Mean" : "Resolution");
2467 }
2468 leg->AddEntry(e, name, "l");
2469 }
2470 }
2471 }
2472 if (!mConfig.enableLocalOutput && !mConfig.shipToQCAsCanvas) {
2473 continue;
2474 }
2475
2476 if (ii == 5) {
2477 pad->SetLogx();
2478 }
2479 can->cd();
2480 if (j == 4) {
2481 ChangePadTitleSize(pad, 0.056);
2482 }
2483 }
2484 if (!mConfig.enableLocalOutput && !mConfig.shipToQCAsCanvas) {
2485 continue;
2486 }
2487
2488 leg->Draw();
2489
2490 if (qcout) {
2491 qcout->Add(can);
2492 }
2493 if (!mConfig.enableLocalOutput) {
2494 continue;
2495 }
2496 doPerfFigure(0.2, 0.295, 0.025);
2497 can->Print(Form(p ? "plots/pull_vs_%s.pdf" : "plots/res_vs_%s.pdf", VSPARAMETER_NAMES[ii]));
2498 if (mConfig.writeRootFiles) {
2499 can->Print(Form(p ? "plots/pull_vs_%s.root" : "plots/res_vs_%s.root", VSPARAMETER_NAMES[ii]));
2500 }
2501 }
2502 }
2503
2504 // Process Integral Resolution Histogreams
2505 for (int32_t p = 0; p < 2; p++) {
2506 if ((p == 0 && (mQATasks & taskTrackingRes) == 0) || (p == 1 && (mQATasks & taskTrackingResPull) == 0)) {
2507 continue;
2508 }
2509 TCanvas* can = p ? mCPull[6] : mCRes[6];
2510 for (int32_t i = 0; i < 5; i++) {
2511 TPad* pad = p ? mPPull[6][i] : mPRes[6][i];
2512 TH1D* hist = p ? pullIntegral[i] : resIntegral[i];
2513 int32_t numColor = 0;
2514 if (mConfig.enableLocalOutput || mConfig.shipToQCAsCanvas) {
2515 pad->cd();
2516 }
2517 if (!mConfig.inputHistogramsOnly && mcAvail) {
2518 TH1D* e = hist;
2519 if (e && e->GetEntries()) {
2520 e->Fit("gaus", "sQ");
2521 }
2522 }
2523
2524 float tmpMax = 0;
2525 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2526 TH1D* e = hist;
2527 if (GetHist(e, tin, k, nNewInput) == nullptr) {
2528 continue;
2529 }
2530 e->SetMaximum(-1111);
2531 if (e->GetMaximum() > tmpMax) {
2532 tmpMax = e->GetMaximum();
2533 }
2534 }
2535
2536 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2537 TH1D* e = hist;
2538 if (GetHist(e, tin, k, nNewInput) == nullptr) {
2539 continue;
2540 }
2541 e->SetMaximum(tmpMax * 1.02);
2542 e->SetMinimum(tmpMax * -0.02);
2543 if (tout && !mConfig.inputHistogramsOnly && k == 0) {
2544 e->Write();
2545 }
2546 if (qcout && !mConfig.shipToQCAsCanvas) {
2547 qcout->Add(e);
2548 }
2549 if (!mConfig.enableLocalOutput && !mConfig.shipToQCAsCanvas) {
2550 continue;
2551 }
2552
2553 e->SetLineColor(colorNums[numColor++ % COLORCOUNT]);
2554 e->Draw(k == 0 ? "" : "same");
2555 }
2556 if (!mConfig.enableLocalOutput && !mConfig.shipToQCAsCanvas) {
2557 continue;
2558 }
2559 can->cd();
2560 }
2561 if (qcout) {
2562 qcout->Add(can);
2563 }
2564 if (!mConfig.enableLocalOutput) {
2565 continue;
2566 }
2567
2568 can->Print(p ? "plots/pull_integral.pdf" : "plots/res_integral.pdf");
2569 if (mConfig.writeRootFiles) {
2570 can->Print(p ? "plots/pull_integral.root" : "plots/res_integral.root");
2571 }
2572 }
2573 }
2574
2575 uint64_t attachClusterCounts[N_CLS_HIST];
2576 if (mQATasks & taskClusterAttach) {
2577 // Process Cluster Attachment Histograms
2578 if (mConfig.inputHistogramsOnly == 0) {
2579 for (int32_t i = N_CLS_HIST; i < N_CLS_TYPE * N_CLS_HIST - 1; i++) {
2580 mClusters[i]->Sumw2(true);
2581 }
2582 double totalVal = 0;
2583 if (!CLUST_HIST_INT_SUM) {
2584 for (int32_t j = 0; j < mClusters[N_CLS_HIST - 1]->GetXaxis()->GetNbins() + 2; j++) {
2585 totalVal += mClusters[N_CLS_HIST - 1]->GetBinContent(j);
2586 }
2587 }
2588 if (totalVal == 0.) {
2589 totalVal = 1.;
2590 }
2591 for (int32_t i = 0; i < N_CLS_HIST; i++) {
2592 double val = 0;
2593 for (int32_t j = 0; j < mClusters[i]->GetXaxis()->GetNbins() + 2; j++) {
2594 val += mClusters[i]->GetBinContent(j);
2595 mClusters[2 * N_CLS_HIST - 1 + i]->SetBinContent(j, val / totalVal);
2596 }
2597 attachClusterCounts[i] = val;
2598 }
2599
2600 if (!CLUST_HIST_INT_SUM) {
2601 for (int32_t i = 0; i < N_CLS_HIST; i++) {
2602 mClusters[2 * N_CLS_HIST - 1 + i]->SetMaximum(1.02);
2603 mClusters[2 * N_CLS_HIST - 1 + i]->SetMinimum(-0.02);
2604 }
2605 }
2606
2607 for (int32_t i = 0; i < N_CLS_HIST - 1; i++) {
2608 auto oldLevel = gErrorIgnoreLevel;
2609 gErrorIgnoreLevel = kError;
2610 mClusters[N_CLS_HIST + i]->Divide(mClusters[i], mClusters[N_CLS_HIST - 1], 1, 1, "B");
2611 gErrorIgnoreLevel = oldLevel;
2612 mClusters[N_CLS_HIST + i]->SetMinimum(-0.02);
2613 mClusters[N_CLS_HIST + i]->SetMaximum(1.02);
2614 }
2615 }
2616
2617 float tmpMax[2] = {0, 0}, tmpMin[2] = {0, 0};
2618 for (int32_t l = 0; l <= CLUST_HIST_INT_SUM; l++) {
2619 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2620 TH1* e = mClusters[l ? (N_CLS_TYPE * N_CLS_HIST - 2) : (N_CLS_HIST - 1)];
2621 if (GetHist(e, tin, k, nNewInput) == nullptr) {
2622 continue;
2623 }
2624 e->SetMinimum(-1111);
2625 e->SetMaximum(-1111);
2626 if (l == 0) {
2627 e->GetXaxis()->SetRange(2, AXIS_BINS[4]);
2628 }
2629 if (e->GetMaximum() > tmpMax[l]) {
2630 tmpMax[l] = e->GetMaximum();
2631 }
2632 if (e->GetMinimum() < tmpMin[l]) {
2633 tmpMin[l] = e->GetMinimum();
2634 }
2635 }
2636 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2637 for (int32_t i = 0; i < N_CLS_HIST; i++) {
2638 TH1* e = mClusters[l ? (2 * N_CLS_HIST - 1 + i) : i];
2639 if (GetHist(e, tin, k, nNewInput) == nullptr) {
2640 continue;
2641 }
2642 e->SetMaximum(tmpMax[l] * 1.02);
2643 e->SetMinimum(tmpMax[l] * -0.02);
2644 }
2645 }
2646 }
2647
2648 for (int32_t i = 0; i < N_CLS_TYPE; i++) {
2649 if (mConfig.enableLocalOutput || mConfig.shipToQCAsCanvas) {
2650 mPClust[i]->cd();
2651 mPClust[i]->SetLogx();
2652 }
2653 int32_t begin = i == 2 ? (2 * N_CLS_HIST - 1) : i == 1 ? N_CLS_HIST : 0;
2654 int32_t end = i == 2 ? (3 * N_CLS_HIST - 1) : i == 1 ? (2 * N_CLS_HIST - 1) : N_CLS_HIST;
2655 int32_t numColor = 0;
2656 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2657 for (int32_t j = end - 1; j >= begin; j--) {
2658 TH1* e = mClusters[j];
2659 if (GetHist(e, tin, k, nNewInput) == nullptr) {
2660 continue;
2661 }
2662
2663 e->SetTitle(CLUSTER_TITLES[i]);
2664 e->GetYaxis()->SetTitle(i == 0 ? "Number of TPC clusters" : i == 1 ? "Fraction of TPC clusters" : CLUST_HIST_INT_SUM ? "Total TPC clusters (integrated)" : "Fraction of TPC clusters (integrated)");
2665 e->GetXaxis()->SetTitle("#it{p}_{Tmc} (GeV/#it{c})");
2666 e->GetXaxis()->SetTitleOffset(1.1);
2667 e->GetXaxis()->SetLabelOffset(-0.005);
2668 if (tout && !mConfig.inputHistogramsOnly && k == 0) {
2669 e->Write();
2670 }
2671 e->SetStats(kFALSE);
2672 e->SetLineWidth(1);
2673 e->SetLineStyle(CONFIG_DASHED_MARKERS ? j + 1 : 1);
2674 if (i == 0) {
2675 e->GetXaxis()->SetRange(2, AXIS_BINS[4]);
2676 }
2677 if (qcout && !mConfig.shipToQCAsCanvas) {
2678 qcout->Add(e);
2679 }
2680 if (!mConfig.enableLocalOutput && !mConfig.shipToQCAsCanvas) {
2681 continue;
2682 }
2683
2684 e->SetMarkerColor(kBlack);
2685 e->SetLineColor(colorNums[numColor++ % COLORCOUNT]);
2686 e->Draw(j == end - 1 && k == 0 ? "" : "same");
2687 GetName(fname, k);
2688 snprintf(name, 2048, "%s%s", fname, CLUSTER_NAMES[j - begin]);
2689 mLClust[i]->AddEntry(e, name, "l");
2690 }
2691 }
2692 if (ConfigNumInputs == 1) {
2693 TH1* e = reinterpret_cast<TH1F*>(mClusters[begin + CL_att_adj]->Clone());
2694 e->Add(mClusters[begin + CL_prot], -1);
2695 if (qcout && !mConfig.shipToQCAsCanvas) {
2696 qcout->Add(e);
2697 }
2698 if (!mConfig.enableLocalOutput && !mConfig.shipToQCAsCanvas) {
2699 continue;
2700 }
2701
2702 e->SetLineColor(colorNums[numColor++ % COLORCOUNT]);
2703 e->Draw("same");
2704 mLClust[i]->AddEntry(e, "Removed (Strategy A)", "l");
2705 }
2706 if (!mConfig.enableLocalOutput && !mConfig.shipToQCAsCanvas) {
2707 continue;
2708 }
2709
2710 mLClust[i]->Draw();
2711
2712 if (qcout) {
2713 qcout->Add(mCClust[i]);
2714 }
2715 if (!mConfig.enableLocalOutput) {
2716 continue;
2717 }
2718 doPerfFigure(i != 2 ? 0.37 : 0.6, 0.295, 0.030);
2719 mCClust[i]->cd();
2720 mCClust[i]->Print(i == 2 ? "plots/clusters_integral.pdf" : i == 1 ? "plots/clusters_relative.pdf" : "plots/clusters.pdf");
2721 if (mConfig.writeRootFiles) {
2722 mCClust[i]->Print(i == 2 ? "plots/clusters_integral.root" : i == 1 ? "plots/clusters_relative.root" : "plots/clusters.root");
2723 }
2724 }
2725 }
2726
2727 // Process cluster count statistics
2728 if ((mQATasks & taskClusterCounts) && !mHaveExternalHists && !mConfig.clusterRejectionHistograms && !mConfig.inputHistogramsOnly) {
2729 DoClusterCounts(attachClusterCounts);
2730 }
2731 if ((qcout || tout) && (mQATasks & taskClusterCounts) && mConfig.clusterRejectionHistograms) {
2732 for (uint32_t i = 0; i < mHistClusterCount.size(); i++) {
2733 if (tout) {
2734 mHistClusterCount[i]->Write();
2735 }
2736 if (qcout) {
2737 qcout->Add(mHistClusterCount[i]);
2738 }
2739 }
2740 }
2741
2742 if (mQATasks & taskTrackStatistics) {
2743 // Process track statistic histograms
2744 float tmpMax = 0.;
2745 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2746 TH1F* e = mTracks;
2747 if (GetHist(e, tin, k, nNewInput) == nullptr) {
2748 continue;
2749 }
2750 e->SetMaximum(-1111);
2751 if (e->GetMaximum() > tmpMax) {
2752 tmpMax = e->GetMaximum();
2753 }
2754 }
2755 mPTracks->cd();
2756 mPTracks->SetLogx();
2757 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2758 TH1F* e = mTracks;
2759 if (GetHist(e, tin, k, nNewInput) == nullptr) {
2760 continue;
2761 }
2762 if (tout && !mConfig.inputHistogramsOnly && k == 0) {
2763 e->Write();
2764 }
2765 e->SetMaximum(tmpMax * 1.02);
2766 e->SetMinimum(tmpMax * -0.02);
2767 e->SetStats(kFALSE);
2768 e->SetLineWidth(1);
2769 e->GetYaxis()->SetTitle("a.u.");
2770 e->GetXaxis()->SetTitle("#it{p}_{Tmc} (GeV/#it{c})");
2771 if (qcout) {
2772 qcout->Add(e);
2773 }
2774 e->SetMarkerColor(kBlack);
2775 e->SetLineColor(colorNums[k % COLORCOUNT]);
2776 e->Draw(k == 0 ? "" : "same");
2777 GetName(fname, k);
2778 snprintf(name, 2048, "%sTrack Pt", fname);
2779 mLTracks->AddEntry(e, name, "l");
2780 }
2781 mLTracks->Draw();
2782 mCTracks->cd();
2783 mCTracks->Print("plots/tracks.pdf");
2784 if (mConfig.writeRootFiles) {
2785 mCTracks->Print("plots/tracks.root");
2786 }
2787
2788 for (int32_t i = 0; i < 2; i++) {
2789 tmpMax = 0.;
2790 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2791 TH1F* e = mNCl[i];
2792 if (GetHist(e, tin, k, nNewInput) == nullptr) {
2793 continue;
2794 }
2795 e->SetMaximum(-1111);
2796 if (e->GetMaximum() > tmpMax) {
2797 tmpMax = e->GetMaximum();
2798 }
2799 }
2800 mPNCl[i]->cd();
2801 for (int32_t k = 0; k < ConfigNumInputs; k++) {
2802 TH1F* e = mNCl[i];
2803 if (GetHist(e, tin, k, nNewInput) == nullptr) {
2804 continue;
2805 }
2806 if (tout && !mConfig.inputHistogramsOnly && k == 0) {
2807 e->Write();
2808 }
2809 e->SetMaximum(tmpMax * 1.02);
2810 e->SetMinimum(tmpMax * -0.02);
2811 e->SetStats(kFALSE);
2812 e->SetLineWidth(1);
2813 e->GetYaxis()->SetTitle("a.u.");
2814 e->GetXaxis()->SetTitle("NClusters");
2815 if (qcout) {
2816 qcout->Add(e);
2817 }
2818 e->SetMarkerColor(kBlack);
2819 e->SetLineColor(colorNums[k % COLORCOUNT]);
2820 e->Draw(k == 0 ? "" : "same");
2821 GetName(fname, k);
2822 snprintf(name, 2048, "%sNClusters%d", fname, i);
2823 mLNCl[i]->AddEntry(e, name, "l");
2824 }
2825 mLNCl[i]->Draw();
2826 mCNCl[i]->cd();
2827 snprintf(name, 2048, "plots/nClusters%s.pdf", i ? "_corrected" : "");
2828 mCNCl[i]->Print(name);
2829 if (mConfig.writeRootFiles) {
2830 snprintf(name, 2048, "plots/nClusters%s.root", i ? "_corrected" : "");
2831 mCNCl[i]->Print(name);
2832 }
2833 }
2834
2835 mPClXY->cd();
2836 mClXY->SetOption("colz");
2837 mClXY->Draw();
2838 mCClXY->cd();
2839 mCClXY->Print("plots/clustersXY.pdf");
2840 if (mConfig.writeRootFiles) {
2841 mCClXY->Print("plots/clustersXY.root");
2842 }
2843 }
2844
2845 if (tout && !mConfig.inputHistogramsOnly && mConfig.writeMCLabels) {
2846 gInterpreter->GenerateDictionary("vector<vector<int32_t>>", "");
2847 tout->WriteObject(&mcEffBuffer, "mcEffBuffer");
2848 tout->WriteObject(&mcLabelBuffer, "mcLabelBuffer");
2849 remove("AutoDict_vector_vector_int__.cxx");
2850 remove("AutoDict_vector_vector_int___cxx_ACLiC_dict_rdict.pcm");
2851 remove("AutoDict_vector_vector_int___cxx.d");
2852 remove("AutoDict_vector_vector_int___cxx.so");
2853 }
2854
2855 if (tout) {
2856 tout->Close();
2857 }
2858 for (uint32_t i = 0; i < mConfig.compareInputs.size(); i++) {
2859 tin[i]->Close();
2860 }
2861 if (!qcout) {
2862 clearGarbagageCollector();
2863 }
2864 GPUInfo("GPU TPC QA histograms have been written to %s files", mConfig.writeRootFiles ? ".pdf and .root" : ".pdf");
2865 gErrorIgnoreLevel = oldRootIgnoreLevel;
2866 return (0);
2867}
2868
2869void GPUQA::PrintClusterCount(int32_t mode, int32_t& num, const char* name, uint64_t n, uint64_t normalization)
2870{
2871 if (mode == 2) {
2872 // do nothing, just count num
2873 } else if (mode == 1) {
2874 char name2[128];
2875 snprintf(name2, 128, "clusterCount%d_", num);
2876 char* ptr = name2 + strlen(name2);
2877 for (uint32_t i = 0; i < strlen(name); i++) {
2878 if ((name[i] >= 'a' && name[i] <= 'z') || (name[i] >= 'A' && name[i] <= 'Z') || (name[i] >= '0' && name[i] <= '9')) {
2879 *(ptr++) = name[i];
2880 }
2881 }
2882 *ptr = 0;
2883 createHist(mHistClusterCount[num], name2, name, 1000, 0, mConfig.histMaxNClusters, 1000, 0, 100);
2884 } else if (mode == 0) {
2885 if (normalization && mConfig.enableLocalOutput) {
2886 printf("\t%40s: %'12" PRIu64 " (%6.2f%%)\n", name, n, 100.f * n / normalization);
2887 }
2888 if (mConfig.clusterRejectionHistograms) {
2889 float ratio = 100.f * n / std::max<uint64_t>(normalization, 1);
2890 mHistClusterCount[num]->Fill(normalization, ratio, 1);
2891 }
2892 }
2893 num++;
2894}
2895
2896int32_t GPUQA::DoClusterCounts(uint64_t* attachClusterCounts, int32_t mode)
2897{
2898 int32_t num = 0;
2899 if (mcPresent() && (mQATasks & taskClusterAttach) && attachClusterCounts) {
2900 for (int32_t i = 0; i < N_CLS_HIST; i++) {
2901 PrintClusterCount(mode, num, CLUSTER_NAMES[i], attachClusterCounts[i], mClusterCounts.nTotal);
2902 }
2903 PrintClusterCount(mode, num, "Unattached", attachClusterCounts[N_CLS_HIST - 1] - attachClusterCounts[CL_att_adj], mClusterCounts.nTotal);
2904 PrintClusterCount(mode, num, "Removed (Strategy A)", attachClusterCounts[CL_att_adj] - attachClusterCounts[CL_prot], mClusterCounts.nTotal); // Attached + Adjacent (also fake) - protected
2905 PrintClusterCount(mode, num, "Unaccessible", mClusterCounts.nUnaccessible, mClusterCounts.nTotal); // No contribution from track >= 10 MeV, unattached or fake-attached/adjacent
2906 } else {
2907 PrintClusterCount(mode, num, "All Clusters", mClusterCounts.nTotal, mClusterCounts.nTotal);
2908 PrintClusterCount(mode, num, "Used in Physics", mClusterCounts.nPhysics, mClusterCounts.nTotal);
2909 PrintClusterCount(mode, num, "Protected", mClusterCounts.nProt, mClusterCounts.nTotal);
2910 PrintClusterCount(mode, num, "Unattached", mClusterCounts.nUnattached, mClusterCounts.nTotal);
2911 PrintClusterCount(mode, num, "Removed (Strategy A)", mClusterCounts.nTotal - mClusterCounts.nUnattached - mClusterCounts.nProt, mClusterCounts.nTotal);
2912 PrintClusterCount(mode, num, "Removed (Strategy B)", mClusterCounts.nTotal - mClusterCounts.nProt, mClusterCounts.nTotal);
2913 }
2914
2915 PrintClusterCount(mode, num, "Merged Loopers (Afterburner)", mClusterCounts.nMergedLooper, mClusterCounts.nTotal);
2916 PrintClusterCount(mode, num, "High Inclination Angle", mClusterCounts.nHighIncl, mClusterCounts.nTotal);
2917 PrintClusterCount(mode, num, "Rejected", mClusterCounts.nRejected, mClusterCounts.nTotal);
2918 PrintClusterCount(mode, num, "Tube (> 200 MeV)", mClusterCounts.nTube, mClusterCounts.nTotal);
2919 PrintClusterCount(mode, num, "Tube (< 200 MeV)", mClusterCounts.nTube200, mClusterCounts.nTotal);
2920 PrintClusterCount(mode, num, "Looping Legs", mClusterCounts.nLoopers, mClusterCounts.nTotal);
2921 PrintClusterCount(mode, num, "Low Pt < 50 MeV", mClusterCounts.nLowPt, mClusterCounts.nTotal);
2922 PrintClusterCount(mode, num, "Low Pt < 200 MeV", mClusterCounts.n200MeV, mClusterCounts.nTotal);
2923
2924 if (mcPresent() && (mQATasks & taskClusterAttach)) {
2925 PrintClusterCount(mode, num, "Tracks > 400 MeV", mClusterCounts.nAbove400, mClusterCounts.nTotal);
2926 PrintClusterCount(mode, num, "Fake Removed (> 400 MeV)", mClusterCounts.nFakeRemove400, mClusterCounts.nAbove400);
2927 PrintClusterCount(mode, num, "Full Fake Removed (> 400 MeV)", mClusterCounts.nFullFakeRemove400, mClusterCounts.nAbove400);
2928 PrintClusterCount(mode, num, "Tracks < 40 MeV", mClusterCounts.nBelow40, mClusterCounts.nTotal);
2929 PrintClusterCount(mode, num, "Fake Protect (< 40 MeV)", mClusterCounts.nFakeProtect40, mClusterCounts.nBelow40);
2930 }
2931 if (mcPresent() && (mQATasks & taskTrackStatistics)) {
2932 PrintClusterCount(mode, num, "Correctly Attached all-trk normalized", mClusterCounts.nCorrectlyAttachedNormalized, mClusterCounts.nTotal);
2933 PrintClusterCount(mode, num, "Correctly Attached non-fake normalized", mClusterCounts.nCorrectlyAttachedNormalizedNonFake, mClusterCounts.nTotal);
2934 }
2935 return num;
2936}
2937
2938void* GPUQA::AllocateScratchBuffer(size_t nBytes)
2939{
2940 mTrackingScratchBuffer.resize((nBytes + sizeof(mTrackingScratchBuffer[0]) - 1) / sizeof(mTrackingScratchBuffer[0]));
2941 return mTrackingScratchBuffer.data();
2942}
A const (ready only) version of MCTruthContainer.
Helper class to access correction maps.
int16_t charge
Definition RawEventData.h:5
int32_t i
#define GPUCA_MIN_TRACK_PTB5_DEFAULT
#define TRACK_EXPECTED_REFERENCE_X_DEFAULT
Definition GPUQA.cxx:210
#define TRACK_EXPECTED_REFERENCE_X
Definition GPUQA.cxx:265
#define QA_DEBUG
Definition GPUQA.cxx:15
#define QA_TIMING
Definition GPUQA.cxx:16
#define CHECK_CLUSTER_STATE_NOCOUNT()
Definition GPUQA.cxx:127
#define CHECK_CLUSTER_STATE()
Definition GPUQA.cxx:114
int16_t Color_t
Definition GPUQA.h:32
GPUChain * chain
uint8_t leg
#define GPUCA_ROW_COUNT
Definition of the MCTrack class.
Definition of the Names Generator class.
uint16_t pos
Definition RawData.h:3
uint32_t j
Definition RawData.h:0
uint32_t side
Definition RawData.h:0
uint16_t pid
Definition RawData.h:2
uint32_t c
Definition RawData.h:2
Definition of TPCFastTransform class.
TBranch * ptr
int nClusters
double num
void Start()
Definition timer.cxx:57
double GetCurrentElapsedTime(bool reset=false)
Definition timer.cxx:110
Class for time synchronization of RawReader instances.
static constexpr ID TPC
Definition DetID.h:64
static constexpr int32_t NSECTORS
Definition GPUChain.h:58
int32_t ReadO2MCData(const char *filename)
Definition GPUQA.h:54
bool clusterRemovable(int32_t attach, bool prot) const
Definition GPUQA.h:52
~GPUQA()=default
Definition GPUQA.cxx:343
void * AllocateScratchBuffer(size_t nBytes)
Definition GPUQA.h:55
void SetMCTrackRange(int32_t min, int32_t max)
Definition GPUQA.h:47
mcLabelI_t GetMCTrackLabel(uint32_t trackId) const
Definition GPUQA.h:51
int32_t DrawQAHistograms()
Definition GPUQA.h:46
int32_t InitQA(int32_t tasks=0)
Definition GPUQA.h:44
void DumpO2MCData(const char *filename) const
Definition GPUQA.h:53
int32_t mcLabelI_t
Definition GPUQA.h:43
void RunQA(bool matchOnly=false)
Definition GPUQA.h:45
GPUQA(void *chain)
Definition GPUQA.h:41
static GPUROOTDump< T, Args... > getNew(const char *name1, Names... names)
Definition GPUROOTDump.h:64
static DigitizationContext * loadFromFile(std::string_view filename="")
GLdouble n
Definition glcorearb.h:1982
GLfloat GLfloat GLfloat alpha
Definition glcorearb.h:279
GLint GLenum GLint x
Definition glcorearb.h:403
GLenum mode
Definition glcorearb.h:266
GLenum src
Definition glcorearb.h:1767
GLuint GLfloat GLfloat GLfloat GLfloat y1
Definition glcorearb.h:5034
GLsizeiptr size
Definition glcorearb.h:659
GLuint GLuint end
Definition glcorearb.h:469
const GLdouble * v
Definition glcorearb.h:832
GLuint const GLchar * name
Definition glcorearb.h:781
GLuint GLuint GLfloat weight
Definition glcorearb.h:5477
GLenum GLint * range
Definition glcorearb.h:1899
GLint y
Definition glcorearb.h:270
GLenum GLenum dst
Definition glcorearb.h:1767
GLboolean * data
Definition glcorearb.h:298
GLuint GLsizei const GLchar * label
Definition glcorearb.h:2519
typedef void(APIENTRYP PFNGLCULLFACEPROC)(GLenum mode)
GLuint GLfloat * val
Definition glcorearb.h:1582
GLenum GLfloat param
Definition glcorearb.h:271
GLuint id
Definition glcorearb.h:650
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
float float float y2
Definition MathUtils.h:44
constexpr int LHCBCPERTIMEBIN
Definition Constants.h:38
Enum< T >::Iterator begin(Enum< T >)
Definition Defs.h:173
value_T f3
Definition TrackUtils.h:93
value_T f1
Definition TrackUtils.h:91
value_T f2
Definition TrackUtils.h:92
struct o2::upgrades_utils::@454 tracks
structure to keep trigger-related info
Defining DataPointCompositeObject explicitly as copiable.
std::string filename()
bool isValid(std::string alias)
int64_t differenceInBC(const InteractionRecord &other) const
std::tuple< std::vector< std::unique_ptr< TCanvas > >, std::vector< std::unique_ptr< TLegend > >, std::vector< std::unique_ptr< TPad > >, std::vector< std::unique_ptr< TLatex > >, std::vector< std::unique_ptr< TH1D > > > v
Definition GPUQA.cxx:321
IR getFirstIRofTF(const IR &rec) const
get 1st IR of TF corresponding to the 1st sampled orbit (in MC)
Definition HBFUtils.h:71
IR getFirstIR() const
Definition HBFUtils.h:47
constexpr size_t min
constexpr size_t max
o2::InteractionRecord ir(0, 0)
vec clear()
o2::InteractionRecord ir0(3, 5)