Project
Loading...
Searching...
No Matches
Tracks.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
12#define _USE_MATH_DEFINES
13
14#include <cmath>
15#include <memory>
16
17// root includes
18#include "TFile.h"
19#include "TRandom3.h"
20
21// o2 includes
24#include "GPUCommonArray.h"
26#include "TPCQC/Tracks.h"
27#include "TPCQC/Helpers.h"
28
30
31using namespace o2::tpc::qc;
32struct binning {
33 int bins;
34 double min;
35 double max;
36};
37// DCA histograms
38const std::vector<std::string_view> types{"A_Pos", "A_Neg", "C_Pos", "C_Neg"};
39const binning binsDCAr{200, -5., 5.};
40const binning binsDCArLargerRange{400, -10., 10.};
41const binning binsEta{200, -1., 1.};
42const binning binsClus{120, 60., 180.};
43const binning binsClusLargerRange{140, 60., 200.};
44//______________________________________________________________________________
45void Tracks::initializeHistograms()
46{
47
48 TH1::AddDirectory(false);
49 const auto logPtBinning = helpers::makeLogBinning(100, 0.05, 20);
50 // 1d hitograms
51 mMapHist["hNClustersBeforeCuts"] = std::make_unique<TH1F>("hNClustersBeforeCuts", "Number of clusters (before cuts);# TPC clusters", binsClus.max, 0, binsClus.max);
52 mMapHist["hNClustersAfterCuts"] = std::make_unique<TH1F>("hNClustersAfterCuts", "Number of clusters;# TPC clusters", binsClus.bins, binsClus.min, binsClus.max);
53 mMapHist["hEta"] = std::make_unique<TH1F>("hEta", "Pseudorapidity;#eta", binsEta.bins, binsEta.min, binsEta.max);
54 mMapHist["hPhiAside"] = std::make_unique<TH1F>("hPhiAside", "Azimuthal angle, A side;#phi", 360, 0., 2 * M_PI);
55 mMapHist["hPhiCside"] = std::make_unique<TH1F>("hPhiCside", "Azimuthal angle, C side;#phi", 360, 0., 2 * M_PI);
56 mMapHist["hPt"] = std::make_unique<TH1F>("hPt", "Transverse momentum;#it{p}_{T} (GeV/#it{c})", logPtBinning.size() - 1, logPtBinning.data());
57 mMapHist["hSign"] = std::make_unique<TH1F>("hSign", "Sign of electric charge;charge sign", 3, -1.5, 1.5);
58 if (!mTurnOffHistosForAsync) {
59 mMapHist["hEtaNeg"] = std::make_unique<TH1F>("hEtaNeg", "Pseudorapidity, neg. tracks;#eta", binsEta.bins, binsEta.min, binsEta.max);
60 mMapHist["hEtaPos"] = std::make_unique<TH1F>("hEtaPos", "Pseudorapidity, pos. tracks;#eta", binsEta.bins, binsEta.min, binsEta.max);
61 mMapHist["hPhiAsideNeg"] = std::make_unique<TH1F>("hPhiAsideNeg", "Azimuthal angle, A side, neg. tracks;#phi", 360, 0., 2 * M_PI);
62 mMapHist["hPhiAsidePos"] = std::make_unique<TH1F>("hPhiAsidePos", "Azimuthal angle, A side, pos. tracks;#phi", 360, 0., 2 * M_PI);
63 mMapHist["hPhiCsideNeg"] = std::make_unique<TH1F>("hPhiCsideNeg", "Azimuthal angle, C side, neg. tracks;#phi", 360, 0., 2 * M_PI);
64 mMapHist["hPhiCsidePos"] = std::make_unique<TH1F>("hPhiCsidePos", "Azimuthal angle, C side, pos. tracks;#phi", 360, 0., 2 * M_PI);
65 // 1d ratio histograms
66 mMapHist["hEtaRatio"] = std::make_unique<TH1F>("hEtaRatio", "Pseudorapidity, ratio neg./pos.;#eta", binsEta.bins, binsEta.min, binsEta.max);
67 mMapHist["hPhiAsideRatio"] = std::make_unique<TH1F>("hPhiAsideRatio", "Azimuthal angle, A side, ratio neg./pos.;#phi", 360, 0., 2 * M_PI);
68 mMapHist["hPhiCsideRatio"] = std::make_unique<TH1F>("hPhiCsideRatio", "Azimuthal angle, C side, ratio neg./pos.;#phi", 360, 0., 2 * M_PI);
69 mMapHist["hPtRatio"] = std::make_unique<TH1F>("hPtRatio", "Transverse momentum, ratio neg./pos. ;#it{p}_{T}", logPtBinning.size() - 1, logPtBinning.data());
70 mMapHist["hPhiBothSides"] = std::make_unique<TH1F>("hPhiBothSides", "Azimuthal angle, both sides clusters;#phi", 360, 0., 2 * M_PI);
71 mMapHist["hPtNeg"] = std::make_unique<TH1F>("hPtNeg", "Transverse momentum, neg. tracks;#it{p}_{T} (GeV/#it{c})", logPtBinning.size() - 1, logPtBinning.data());
72 mMapHist["hPtPos"] = std::make_unique<TH1F>("hPtPos", "Transverse momentum, pos. tracks;#it{p}_{T} (GeV/#it{c})", logPtBinning.size() - 1, logPtBinning.data());
73 }
74 mMapHist["hEtaBeforeCuts"] = std::make_unique<TH1F>("hEtaBeforeCuts", "Pseudorapidity (before cuts);#eta", 400, -2., 2.);
75 mMapHist["hPtBeforeCuts"] = std::make_unique<TH1F>("hPtBeforeCuts", "Transverse momentum (before cuts);#it{p}_{T} (GeV/#it{c})", logPtBinning.size() - 1, logPtBinning.data());
76 mMapHist["hQOverPt"] = std::make_unique<TH1F>("hQOverPt", "Charge over transverse momentum;q/#it{p}_{T}", 400, -20., 20.);
77 // 2d histograms
78 mMapHist["h2DNClustersEta"] = std::make_unique<TH2F>("h2DNClustersEta", "Number of clusters vs. #eta;#eta;# TPC clusters", binsEta.bins, binsEta.min, binsEta.max, binsClusLargerRange.bins, binsClusLargerRange.min, binsClusLargerRange.max);
79 mMapHist["h2DNClustersPhiAside"] = std::make_unique<TH2F>("h2DNClustersPhiAside", "Number of clusters vs. #phi, A side ;#phi;# TPC clusters", 360, 0., 2 * M_PI, binsClus.bins, binsClus.min, binsClus.max);
80 mMapHist["h2DNClustersPhiCside"] = std::make_unique<TH2F>("h2DNClustersPhiCside", "Number of clusters vs. #phi, C side ;#phi;# TPC clusters", 360, 0., 2 * M_PI, binsClus.bins, binsClus.min, binsClus.max);
81 mMapHist["h2DNClustersPt"] = std::make_unique<TH2F>("h2DNClustersPt", "Number of clusters vs. #it{p}_{T};#it{p}_{T} (GeV/#it{c});# TPC clusters", logPtBinning.size() - 1, logPtBinning.data(), binsClusLargerRange.bins, binsClusLargerRange.min, binsClusLargerRange.max);
82 mMapHist["h2DEtaPhi"] = std::make_unique<TH2F>("h2DEtaPhi", "Tracks in #eta vs. #phi;#phi;#eta", 360, 0., 2 * M_PI, binsEta.bins, binsEta.min, binsEta.max);
83 mMapHist["h2DNClustersEtaBeforeCuts"] = std::make_unique<TH2F>("h2DNClustersEtaBeforeCuts", "NClusters vs. #eta (before cuts);#eta;# TPC clusters", 400, -2., 2., 200, -0.5, 199.5);
84 mMapHist["h2DNClustersPtBeforeCuts"] = std::make_unique<TH2F>("h2DNClustersPtBeforeCuts", "NClusters vs. #it{p}_{T} (before cuts);#it{p}_{T} (GeV/#it{c});# TPC clusters", logPtBinning.size() - 1, logPtBinning.data(), 200, -0.5, 199.5);
85 mMapHist["h2DEtaPhiBeforeCuts"] = std::make_unique<TH2F>("h2DEtaPhiBeforeCuts", "Tracks in #eta vs. #phi (before cuts);#phi;#eta", 360, 0., 2 * M_PI, 400, -2., 2.);
86 if (!mTurnOffHistosForAsync) {
87 mMapHist["h2DQOverPtPhiAside"] = std::make_unique<TH2F>("h2DQOverPtPhiAside", "Charger over #it{p}_{T} vs. #phi, A side;#phi;q/#it{p}_{T}", 360, 0., 2 * M_PI, 400, -20., 20.);
88 mMapHist["h2DQOverPtPhiCside"] = std::make_unique<TH2F>("h2DQOverPtPhiCside", "Charger over #it{p}_{T} vs. #phi, C side;#phi;q/#it{p}_{T}", 360, 0., 2 * M_PI, 400, -20., 20.);
89 mMapHist["h2DEtaPhiNeg"] = std::make_unique<TH2F>("h2DEtaPhiNeg", "Negative tracks in #eta vs. #phi;#phi;#eta", 360, 0., 2 * M_PI, binsEta.bins, binsEta.min, binsEta.max);
90 mMapHist["h2DEtaPhiPos"] = std::make_unique<TH2F>("h2DEtaPhiPos", "Positive tracks in #eta vs. #phi;#phi;#eta", 360, 0., 2 * M_PI, binsEta.bins, binsEta.min, binsEta.max);
91 // eta vs pt and phi vs pt possitive and negative signs
92 mMapHist["hEtaVsPtPos"] = std::make_unique<TH2F>("hEtaVsPtPos", "#eta vs. #it{p}_{T} (Pos.);#it{p}_{T} (GeV/#it{c});#eta", logPtBinning.size() - 1, logPtBinning.data(), binsEta.bins, binsEta.min, binsEta.max);
93 mMapHist["hEtaVsPtNeg"] = std::make_unique<TH2F>("hEtaVsPtNeg", "#eta vs. #it{p}_{T} (Neg.);#it{p}_{T} (GeV/#it{c});#eta", logPtBinning.size() - 1, logPtBinning.data(), binsEta.bins, binsEta.min, binsEta.max);
94 mMapHist["hPhiVsPtPos"] = std::make_unique<TH2F>("hPhiVsPtPos", "#phi vs. #it{p}_{T} (Pos.);#it{p}_{T} (GeV/#it{c});#phi", logPtBinning.size() - 1, logPtBinning.data(), 360, 0., 2 * M_PI);
95 mMapHist["hPhiVsPtNeg"] = std::make_unique<TH2F>("hPhiVsPtNeg", "#phi vs. #it{p}_{T} (Neg.);#it{p}_{T} (GeV/#it{c});#phi", logPtBinning.size() - 1, logPtBinning.data(), 360, 0., 2 * M_PI);
96 }
97 // DCA Histograms
98 for (const auto type : types) {
99 mMapHist[fmt::format("hDCAr_{}", type).data()] = std::make_unique<TH2F>(fmt::format("hDCAr_{}", type).data(), fmt::format("DCAr {};#phi;DCAr (cm)", type).data(), 360, 0, o2::math_utils::twoPid(), binsDCAr.bins, binsDCAr.min, binsDCAr.max);
100 mMapHist[fmt::format("hDCAr_{}_pTmin", type).data()] = std::make_unique<TH2F>(fmt::format("hDCAr_{}_pTmin", type).data(), fmt::format("DCAr {} #it{{p}}_{{T}}^{{min}};#phi;DCAr (cm)", type).data(), 360, 0, o2::math_utils::twoPid(), binsDCAr.bins, binsDCAr.min, binsDCAr.max);
101 }
102 // DCA vs variables Histograms
103 mMapHist["hDCArVsPtPos"] = std::make_unique<TH2F>("hDCArVsPtPos", "DCAr Pos;#it{p}_{T} (GeV/#it{c});DCAr (cm)", logPtBinning.size() - 1, logPtBinning.data(), binsDCArLargerRange.bins, binsDCArLargerRange.min, binsDCArLargerRange.max);
104 mMapHist["hDCArVsEtaPos"] = std::make_unique<TH2F>("hDCArVsEtaPos", "DCAr Pos;#eta;DCAr (cm)", binsEta.bins, binsEta.min, binsEta.max, binsDCArLargerRange.bins, binsDCArLargerRange.min, binsDCArLargerRange.max);
105 mMapHist["hDCArVsNClsPos"] = std::make_unique<TH2F>("hDCArVsNClsPos", "DCAr Pos;# TPC clusters;DCAr (cm)", binsClus.bins, binsClus.min, binsClus.max, binsDCAr.bins, binsDCAr.min, binsDCAr.max);
106 mMapHist["hDCArVsPtNeg"] = std::make_unique<TH2F>("hDCArVsPtNeg", "DCAr Neg;#it{p}_{T} (GeV/#it{c});DCAr (cm)", logPtBinning.size() - 1, logPtBinning.data(), binsDCArLargerRange.bins, binsDCArLargerRange.min, binsDCArLargerRange.max);
107 mMapHist["hDCArVsEtaNeg"] = std::make_unique<TH2F>("hDCArVsEtaNeg", "DCAr Neg;#eta;DCAr (cm)", binsEta.bins, binsEta.min, binsEta.max, binsDCArLargerRange.bins, binsDCArLargerRange.min, binsDCArLargerRange.max);
108 mMapHist["hDCArVsNClsNeg"] = std::make_unique<TH2F>("hDCArVsNClsNeg", "DCAr Neg;# TPC clusters;DCAr (cm)", binsClus.bins, binsClus.min, binsClus.max, binsDCAr.bins, binsDCAr.min, binsDCAr.max);
109 // DCA vs variables Histogram with pT min selection
110 mMapHist["hDCArVsEtaPos_pTmin"] = std::make_unique<TH2F>("hDCArVsEtaPos_pTmin", "DCAr Pos #it{p}_{T}^{min};#eta;DCAr (cm)", binsEta.bins, binsEta.min, binsEta.max, binsDCAr.bins, binsDCAr.min, binsDCAr.max);
111 mMapHist["hDCArVsNClsPos_pTmin"] = std::make_unique<TH2F>("hDCArVsNClsPos_pTmin", "DCAr Pos #it{p}_{T}^{min};# TPC clusters;DCAr (cm)", binsClus.bins, binsClus.min, binsClus.max, binsDCAr.bins, binsDCAr.min, binsDCAr.max);
112 mMapHist["hDCArVsEtaNeg_pTmin"] = std::make_unique<TH2F>("hDCArVsEtaNeg_pTmin", "DCAr Neg #it{p}_{T}^{min};#eta;DCAr (cm)", binsEta.bins, binsEta.min, binsEta.max, binsDCAr.bins, binsDCAr.min, binsDCAr.max);
113 mMapHist["hDCArVsNClsNeg_pTmin"] = std::make_unique<TH2F>("hDCArVsNClsNeg_pTmin", "DCAr Neg #it{p}_{T}^{min};# TPC clusters;DCAr (cm)", binsClus.bins, binsClus.min, binsClus.max, binsDCAr.bins, binsDCAr.min, binsDCAr.max);
114}
115//______________________________________________________________________________
117{
118 for (const auto& pair : mMapHist) {
119 pair.second->Reset();
120 }
121}
122//______________________________________________________________________________
124{
125 // ===| variables required for cutting and filling |===
126 const auto eta = track.getEta();
127 const auto phi = track.getPhi();
128 const auto pt = track.getPt();
129 const auto sign = track.getSign();
130 const auto qOverPt = track.getQ2Pt();
131 const auto nCls = track.getNClusterReferences();
132 const auto dEdxTot = track.getdEdx().dEdxTotTPC;
133 const auto hasASideOnly = track.hasASideClustersOnly();
134 const auto hasCSideOnly = track.hasCSideClustersOnly();
135
136 const auto absEta = std::abs(eta);
137
138 // ===| histogram filling before cuts |===
139 mMapHist["hNClustersBeforeCuts"]->Fill(nCls);
140 mMapHist["hEtaBeforeCuts"]->Fill(eta);
141 mMapHist["hPtBeforeCuts"]->Fill(pt);
142 mMapHist["h2DNClustersEtaBeforeCuts"]->Fill(eta, nCls);
143 mMapHist["h2DNClustersPtBeforeCuts"]->Fill(pt, nCls);
144 mMapHist["h2DEtaPhiBeforeCuts"]->Fill(phi, eta);
145
146 // ===| histogram filling including cuts |===
147 if (absEta < mCutAbsEta && nCls > mCutMinnCls && dEdxTot > mCutMindEdxTot) {
148
149 // ===| 1D histogram filling |===
150 mMapHist["hNClustersAfterCuts"]->Fill(nCls);
151 mMapHist["hEta"]->Fill(eta);
152
153 //---| propagate to 0,0,0 |---
154 //
155 // propagator instance must be configured before (LUT, MagField)
156 auto propagator = o2::base::Propagator::Instance(true);
157 const int type = (track.getQ2Pt() < 0) + 2 * track.hasCSideClustersOnly();
158 auto dcaHist = mMapHist[fmt::format("hDCAr_{}", types[type]).data()].get();
159 auto dcaHist_pTmin = mMapHist[fmt::format("hDCAr_{}_pTmin", types[type]).data()].get();
160 const std::string signType((sign < 0) ? "Neg" : "Pos");
161 auto dcaHistPT = mMapHist["hDCArVsPt" + signType].get();
162 auto dcaHistEta = mMapHist["hDCArVsEta" + signType].get();
163 auto dcaHistNCluster = mMapHist["hDCArVsNCls" + signType].get();
164 auto dcaHistEta_pTmin = mMapHist["hDCArVsEta" + signType + "_pTmin"].get();
165 auto dcaHistNCluster_pTmin = mMapHist["hDCArVsNCls" + signType + "_pTmin"].get();
166
167 // set-up sampling for the DCA calculation
168 Double_t sampleProb = 2;
169
170 if (mSamplingFractionDCAr > 0) { // for now no SEED is given.
171 TRandom3 randomGenerator(0);
172 sampleProb = randomGenerator.Uniform(1);
173 }
174
175 bool isDcaCalculated = false;
176 float dcaValue = -999;
177
178 if (sampleProb > (Double_t)(1. - mSamplingFractionDCAr)) {
179
180 if (propagator->getMatLUT() && propagator->hasMagFieldSet()) {
181 // ---| fill DCA histos |---
183 const o2::math_utils::Point3D<float> refPoint{0, 0, 0};
184 o2::track::TrackPar propTrack(track);
185 if (propagator->propagateToDCABxByBz(refPoint, propTrack, 2.f, o2::base::Propagator::MatCorrType::USEMatCorrLUT, &dca)) {
186 const auto phi = o2::math_utils::to02PiGen(track.getPhi());
187 dcaHistPT->Fill(pt, dca[0]);
188 dcaHist->Fill(phi, dca[0]);
189 dcaHistEta->Fill(eta, dca[0]);
190 dcaHistNCluster->Fill(nCls, dca[0]);
191 if (pt > mCutMinPtDCAr) {
192 dcaHist_pTmin->Fill(phi, dca[0]);
193 dcaHistEta_pTmin->Fill(eta, dca[0]);
194 dcaHistNCluster_pTmin->Fill(nCls, dca[0]);
195 }
196 dcaValue = dca[0];
197 isDcaCalculated = true;
198 }
199 } else {
200 static bool reported = false;
201 if (!reported) {
202 LOGP(error, "o2::base::Propagator not properly initialized, MatLUT ({}) and / or Field ({}) missing, will not fill DCA histograms", (void*)propagator->getMatLUT(), propagator->hasMagFieldSet());
203 dcaHist->SetTitle(fmt::format("DCAr {} o2::base::Propagator not properly initialized", types[type]).data());
204 dcaHistPT->SetTitle(fmt::format("DCAr #it{{p}}_{{T}} {} o2::base::Propagator not properly initialized", signType).data());
205 dcaHist_pTmin->SetTitle(fmt::format("DCAr {} #it{{p}}_{{T}}^{{min}} o2::base::Propagator not properly initialized", types[type]).data());
206 dcaHistEta->SetTitle(fmt::format("DCAr #eta {} o2::base::Propagator not properly initialized", signType).data());
207 dcaHistNCluster->SetTitle(fmt::format("DCAr nClusters {} o2::base::Propagator not properly initialized", signType).data());
208 dcaHistEta_pTmin->SetTitle(fmt::format("DCAr #eta {} #it{{p}}_{{T}}^{{min}} o2::base::Propagator not properly initialized", signType).data());
209 dcaHistNCluster_pTmin->SetTitle(fmt::format("DCAr nClusters {} #it{{p}}_{{T}}^{{min}} o2::base::Propagator not properly initialized", signType).data());
210 reported = true;
211 }
212 }
213 }
214
215 if (mUseCutMaxAbsDCArOnHistos && !isDcaCalculated) {
216 // In this case the DCAr selection should be applied but is not available for the track, hence we simply return and report back
217 // For ease we just report back for every histogram that the propagator was not initialized
218 static bool reported = false;
219 if (!reported) {
220 LOGP(error, "o2::base::Propagator not properly initialized, MatLUT ({}) and / or Field ({}) missing, will not fill apply DCAr selection on histograms", (void*)propagator->getMatLUT(), propagator->hasMagFieldSet());
221 mMapHist["hPhiAside"]->SetTitle("hPhiAside o2::base::Propagator not properly initialized");
222 mMapHist["hPhiCside"]->SetTitle("hPhiCside o2::base::Propagator not properly initialized");
223 mMapHist["hPhiBothSides"]->SetTitle("hPhiBothSides o2::base::Propagator not properly initialized");
224 mMapHist["hPt"]->SetTitle("hPt o2::base::Propagator not properly initialized");
225 mMapHist["hSign"]->SetTitle("hSign o2::base::Propagator not properly initialized");
226 mMapHist["hQOverPt"]->SetTitle("hQOverPt o2::base::Propagator not properly initialized");
227 mMapHist["hEtaNeg"]->SetTitle("hEtaNeg o2::base::Propagator not properly initialized");
228 mMapHist["hPhiAsideNeg"]->SetTitle("hPhiAsideNeg o2::base::Propagator not properly initialized");
229 mMapHist["hPhiCsideNeg"]->SetTitle("hPhiCsideNeg o2::base::Propagator not properly initialized");
230 mMapHist["hEtaPos"]->SetTitle("hEtaPos o2::base::Propagator not properly initialized");
231 mMapHist["hPtPos"]->SetTitle("hPtPos o2::base::Propagator not properly initialized");
232 mMapHist["hPhiAsidePos"]->SetTitle("hPhiAsidePos o2::base::Propagator not properly initialized");
233 mMapHist["hPhiCsidePos"]->SetTitle("hPhiCsidePos o2::base::Propagator not properly initialized");
234 mMapHist["h2DNClustersEta"]->SetTitle("h2DNClustersEta o2::base::Propagator not properly initialized");
235 mMapHist["h2DNClustersPhiAside"]->SetTitle("h2DNClustersPhiAside o2::base::Propagator not properly initialized");
236 mMapHist["h2DQOverPtPhiAside"]->SetTitle("h2DQOverPtPhiAside o2::base::Propagator not properly initialized");
237 mMapHist["h2DNClustersPhiCside"]->SetTitle("h2DNClustersPhiCside o2::base::Propagator not properly initialized");
238 mMapHist["h2DQOverPtPhiCside"]->SetTitle("h2DQOverPtPhiCside o2::base::Propagator not properly initialized");
239 mMapHist["h2DNClustersPt"]->SetTitle("h2DNClustersPt o2::base::Propagator not properly initialized");
240 mMapHist["h2DEtaPhi"]->SetTitle("h2DEtaPhi o2::base::Propagator not properly initialized");
241 mMapHist["h2DEtaPhiNeg"]->SetTitle("h2DEtaPhiNeg o2::base::Propagator not properly initialized");
242 mMapHist["hEtaVsPtNeg"]->SetTitle("hEtaVsPtNeg o2::base::Propagator not properly initialized");
243 mMapHist["hPhiVsPtNeg"]->SetTitle("hPhiVsPtNeg o2::base::Propagator not properly initialized");
244 mMapHist["h2DEtaPhiPos"]->SetTitle("h2DEtaPhiPos o2::base::Propagator not properly initialized");
245 mMapHist["hEtaVsPtPos"]->SetTitle("hEtaVsPtPos o2::base::Propagator not properly initialized");
246 mMapHist["hPhiVsPtPos"]->SetTitle("hPhiVsPtPos o2::base::Propagator not properly initialized");
247 }
248 return true;
249 }
250
251 if (mUseCutMaxAbsDCArOnHistos && (std::abs(dcaValue) > mCutMaxAbsDCAr)) {
252 return true;
253 }
254
255 if (hasASideOnly == 1) {
256 mMapHist["hPhiAside"]->Fill(phi);
257 } else if (hasCSideOnly == 1) {
258 mMapHist["hPhiCside"]->Fill(phi);
259 } else {
260 if (!mTurnOffHistosForAsync) {
261 mMapHist["hPhiBothSides"]->Fill(phi);
262 }
263 }
264
265 mMapHist["hPt"]->Fill(pt);
266 mMapHist["hSign"]->Fill(sign);
267 mMapHist["hQOverPt"]->Fill(qOverPt);
268
269 if (sign < 0.) {
270 if (!mTurnOffHistosForAsync) {
271 mMapHist["hEtaNeg"]->Fill(eta);
272 mMapHist["hPtNeg"]->Fill(pt);
273 }
274 if (!mTurnOffHistosForAsync) {
275 if (hasASideOnly == 1) {
276 mMapHist["hPhiAsideNeg"]->Fill(phi);
277 } else if (hasCSideOnly == 1) {
278 mMapHist["hPhiCsideNeg"]->Fill(phi);
279 }
280 }
281 } else {
282 if (!mTurnOffHistosForAsync) {
283 mMapHist["hEtaPos"]->Fill(eta);
284 mMapHist["hPtPos"]->Fill(pt);
285 }
286 if (!mTurnOffHistosForAsync) {
287 if (hasASideOnly == 1) {
288 mMapHist["hPhiAsidePos"]->Fill(phi);
289 } else if (hasCSideOnly == 1) {
290 mMapHist["hPhiCsidePos"]->Fill(phi);
291 }
292 }
293 }
294
295 // ===| 2D histogram filling |===
296 mMapHist["h2DNClustersEta"]->Fill(eta, nCls);
297
298 if (hasASideOnly == 1) {
299 mMapHist["h2DNClustersPhiAside"]->Fill(phi, nCls);
300 if (!mTurnOffHistosForAsync) {
301 mMapHist["h2DQOverPtPhiAside"]->Fill(phi, qOverPt);
302 }
303 } else if (hasCSideOnly == 1) {
304 mMapHist["h2DNClustersPhiCside"]->Fill(phi, nCls);
305 if (!mTurnOffHistosForAsync) {
306 mMapHist["h2DQOverPtPhiCside"]->Fill(phi, qOverPt);
307 }
308 }
309
310 mMapHist["h2DNClustersPt"]->Fill(pt, nCls);
311 mMapHist["h2DEtaPhi"]->Fill(phi, eta);
312
313 if (!mTurnOffHistosForAsync) {
314 if (sign < 0.) {
315 mMapHist["h2DEtaPhiNeg"]->Fill(phi, eta);
316 mMapHist["hEtaVsPtNeg"]->Fill(pt, eta);
317 mMapHist["hPhiVsPtNeg"]->Fill(pt, phi);
318 } else {
319 mMapHist["h2DEtaPhiPos"]->Fill(phi, eta);
320 mMapHist["hEtaVsPtPos"]->Fill(pt, eta);
321 mMapHist["hPhiVsPtPos"]->Fill(pt, phi);
322 }
323 }
324 }
325
326 return true;
327}
328
329//______________________________________________________________________________
331{
332 if (!mTurnOffHistosForAsync) {
333 // ===| Dividing of 1D histograms -> Ratios |===
334 mMapHist["hEtaRatio"]->Divide(mMapHist["hEtaNeg"].get(), mMapHist["hEtaPos"].get());
335 mMapHist["hPhiAsideRatio"]->Divide(mMapHist["hPhiAsideNeg"].get(), mMapHist["hPhiAsidePos"].get());
336 mMapHist["hPhiCsideRatio"]->Divide(mMapHist["hPhiCsideNeg"].get(), mMapHist["hPhiCsidePos"].get());
337 mMapHist["hPtRatio"]->Divide(mMapHist["hPtNeg"].get(), mMapHist["hPtPos"].get());
338 }
339}
340
341//______________________________________________________________________________
342void Tracks::dumpToFile(std::string_view filename)
343{
344 auto f = std::unique_ptr<TFile>(TFile::Open(filename.data(), "recreate"));
345 for (const auto& [name, hist] : mMapHist) {
346 TObjArray arr;
347 arr.SetName(name.data());
348 arr.Add(hist.get());
349 arr.Write(arr.GetName(), TObject::kSingleKey);
350 }
351 f->Close();
352}
const binning binsEta
Definition Tracks.cxx:41
const binning binsDCArLargerRange
Definition Tracks.cxx:40
const binning binsClus
Definition Tracks.cxx:42
ClassImp(o2::tpc::qc::Tracks)
const binning binsClusLargerRange
Definition Tracks.cxx:43
const binning binsDCAr
Definition Tracks.cxx:39
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:143
Tracks quality control class.
Definition Tracks.h:49
void processEndOfCycle()
Function to be called at each endOfCycle.
Definition Tracks.cxx:330
void dumpToFile(std::string_view filename)
Dump results to a file.
Definition Tracks.cxx:342
bool processTrack(const o2::tpc::TrackTPC &track)
Definition Tracks.cxx:123
void resetHistograms()
Reset all histograms.
Definition Tracks.cxx:116
GLsizei GLenum GLenum * types
Definition glcorearb.h:2516
GLuint const GLchar * name
Definition glcorearb.h:781
GLdouble f
Definition glcorearb.h:310
GLint GLint GLsizei GLint GLenum GLenum type
Definition glcorearb.h:275
GLboolean * data
Definition glcorearb.h:298
std::array< T, N > array
float to02PiGen(float phi)
Definition Utils.h:70
std::vector< double > makeLogBinning(const int nbins, const double min, const double max)
get a vector containing binning info for constant sized bins on a log axis
Definition Helpers.cxx:22
std::string filename()
int bins
Definition PID.cxx:34
double max
Definition PID.cxx:36
double min
Definition PID.cxx:35