Project
Loading...
Searching...
No Matches
AvgClusSize.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
15
18
19#include "Framework/Task.h"
23#include "ITStracking/IOUtils.h"
37
38#include <numeric>
39#include <TH1F.h>
40#include <TH2F.h>
41#include <THStack.h>
42#include <TFile.h>
43#include <TCanvas.h>
44#include <TLine.h>
45#include <TStyle.h>
46#include <TNtuple.h>
47
48namespace o2
49{
50namespace its
51{
52namespace study
53{
54using namespace o2::framework;
55using namespace o2::globaltracking;
56
66
67class AvgClusSizeStudy : public Task
68{
69 public:
70 AvgClusSizeStudy(std::shared_ptr<DataRequest> dr,
71 std::shared_ptr<o2::base::GRPGeomRequest> gr,
72 bool isMC,
73 std::shared_ptr<o2::steer::MCKinematicsReader> kineReader) : mDataRequest{dr}, mGGCCDBRequest(gr), mUseMC(isMC), mKineReader(kineReader){};
74 ~AvgClusSizeStudy() final = default;
75 void init(InitContext& ic) final;
76 void run(ProcessingContext&) final;
78 void finaliseCCDB(ConcreteDataMatcher&, void*) final;
79 void setClusterDictionary(const o2::itsmft::TopologyDictionary* d) { mDict = d; }
80
81 private:
82 // Other functions
85
86 // Helper functions
87 void prepareOutput();
88 void setStyle();
89 void updateTimeDependentParams(ProcessingContext& pc);
90 float getAverageClusterSize(o2::its::TrackITS*);
91 float calcV0HypoMass(const V0&, PID, PID);
92 void calcAPVars(const V0&, float*, float*);
93 void getClusterSizes(std::vector<int>&, const gsl::span<const o2::itsmft::CompClusterExt>, gsl::span<const unsigned char>::iterator&, const o2::itsmft::TopologyDictionary*);
94 void saveHistograms();
95 void plotHistograms();
96 void fillEtaBin(float eta, float clusSize, int i);
97
98 // Running options
99 bool mUseMC;
100
101 // Data
102 std::shared_ptr<o2::base::GRPGeomRequest> mGGCCDBRequest;
103 std::shared_ptr<DataRequest> mDataRequest;
104 std::vector<int> mClusterSizes;
105 gsl::span<const int> mInputITSidxs;
106 std::vector<o2::MCTrack> mMCTracks;
107 const o2::itsmft::TopologyDictionary* mDict = nullptr;
108
109 // Output plots
110 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOut;
111 std::unique_ptr<TNtuple> mOutputNtupleAll;
112 std::unique_ptr<TNtuple> mOutputNtupleCut;
113
114 std::unique_ptr<THStack> mMCR{};
115 std::unique_ptr<THStack> mMCCosPA{};
116 std::unique_ptr<THStack> mMCPosACS{};
117 std::unique_ptr<THStack> mMCNegACS{};
118 std::unique_ptr<THStack> mMCDauDCA{};
119 std::unique_ptr<THStack> mMCPosPVDCA{};
120 std::unique_ptr<THStack> mMCNegPVDCA{};
121 std::unique_ptr<THStack> mMCV0PVDCA{};
122
123 std::unique_ptr<TH1F> mMCRisTg{};
124 std::unique_ptr<TH1F> mMCRisBg{};
125 std::unique_ptr<TH1F> mMCCosPAisTg{};
126 std::unique_ptr<TH1F> mMCCosPAisBg{};
127 std::unique_ptr<TH1F> mMCPosACSisTg{};
128 std::unique_ptr<TH1F> mMCPosACSisBg{};
129 std::unique_ptr<TH1F> mMCNegACSisTg{};
130 std::unique_ptr<TH1F> mMCNegACSisBg{};
131 std::unique_ptr<TH1F> mMCDauDCAisTg{};
132 std::unique_ptr<TH1F> mMCDauDCAisBg{};
133 std::unique_ptr<TH1F> mMCPosPVDCAisTg{};
134 std::unique_ptr<TH1F> mMCPosPVDCAisBg{};
135 std::unique_ptr<TH1F> mMCNegPVDCAisTg{};
136 std::unique_ptr<TH1F> mMCNegPVDCAisBg{};
137 std::unique_ptr<TH1F> mMCV0PVDCAisTg{};
138 std::unique_ptr<TH1F> mMCV0PVDCAisBg{};
139
140 std::unique_ptr<TH2F> mMCArmPodolisTg{};
141 std::unique_ptr<TH2F> mMCArmPodolisBg{};
142
143 std::unique_ptr<THStack> mAvgClusSizeCEta{};
144 std::vector<std::unique_ptr<TH1F>> mAvgClusSizeCEtaVec{};
145
146 std::vector<float> mEtaBinUL; // upper edges for eta bins
147
148 // Counters for target V0 identification
149 int nNotValid = 0;
150 int nNullptrs = 0;
151 int nMotherIDMismatch = 0;
152 int nEvIDMismatch = 0;
153 int nTargetV0 = 0;
154 int nNotTargetV0 = 0;
155 int nV0OutOfEtaRange = 0;
156
157 std::string mOutName;
158 std::shared_ptr<o2::steer::MCKinematicsReader> mKineReader;
159};
160
162{
164 LOGP(info, "Starting average cluster size study...");
165
166 prepareOutput();
167
168 if (mUseMC) { // for counting the missed K0shorts
169 mKineReader = std::make_unique<o2::steer::MCKinematicsReader>("collisioncontext.root");
170 for (int iEvent{0}; iEvent < mKineReader->getNEvents(0); iEvent++) {
171 auto mctrk = mKineReader->getTracks(0, iEvent);
172 mMCTracks.insert(mMCTracks.end(), mctrk.begin(), mctrk.end());
173 }
174 }
175
176 LOGP(info, "Cluster size study initialized.");
177}
178
179void AvgClusSizeStudy::prepareOutput()
180{
182 mOutName = params.outFileName;
183 mDBGOut = std::make_unique<o2::utils::TreeStreamRedirector>(mOutName.c_str(), "recreate");
184 mOutputNtupleAll = std::make_unique<TNtuple>("v0_data", "v0 data", "dPosACS:dNegACS:cosPA:V0R:eta:dauDCA:dPospvDCA:dNegpvDCA:v0pvDCA:alpha:pT:K0mass:lambdaMass:antilambdaMass:v0PDGcode");
185 mOutputNtupleCut = std::make_unique<TNtuple>("cut_v0_data", "v0 data (cut)", "dPosACS:dNegACS:cosPA:V0R:eta:dauDCA:dPospvDCA:dNegpvDCA:v0pvDCA:alpha:pT:K0mass:lambdaMass:antilambdaMass:v0PDGcode");
186
187 mMCR = std::make_unique<THStack>("R", "V0 decay length R;R (cm?)");
188 mMCCosPA = std::make_unique<THStack>("cosPA", "cos(#theta_{p})");
189 mMCPosACS = std::make_unique<THStack>("acsPos", "Average cluster size per track;pixels / cluster / track");
190 mMCNegACS = std::make_unique<THStack>("acsNeg", "Average cluster size per track;pixels / cluster / track");
191 mMCDauDCA = std::make_unique<THStack>("dauDCA", "Prong-prong DCA;cm?");
192 mMCPosPVDCA = std::make_unique<THStack>("posPVDCA", "Positive prong-primary vertex DCA;cm?");
193 mMCNegPVDCA = std::make_unique<THStack>("negPVDCA", "Negative prong-primary vertex DCA;cm?");
194 mMCV0PVDCA = std::make_unique<THStack>("v0PVDCA", "V0 reconstructed track-primary vertex DCA;cm?");
195
196 mMCRisTg = std::make_unique<TH1F>("mcRTg", "target V0", 40, 0, 2);
197 mMCRisBg = std::make_unique<TH1F>("mcRBg", "background", 100, 0, 2);
198 mMCCosPAisTg = std::make_unique<TH1F>("mcCosPATg", "target V0", 40, 0, 1);
199 mMCCosPAisBg = std::make_unique<TH1F>("mcCosPABg", "background", 100, 0, 1);
200 mMCPosACSisTg = std::make_unique<TH1F>("mcPosACSTg", "target V0", 60, 0, 30);
201 mMCPosACSisBg = std::make_unique<TH1F>("mcPosACSBg", "background", 150, 0, 30);
202 mMCNegACSisTg = std::make_unique<TH1F>("mcNegACSTg", "target V0", 60, 0, 30);
203 mMCNegACSisBg = std::make_unique<TH1F>("mcNegACSBg", "background", 150, 0, 30);
204 mMCDauDCAisTg = std::make_unique<TH1F>("mcDauDCATg", "target V0", 40, 0, 2);
205 mMCDauDCAisBg = std::make_unique<TH1F>("mcDauDCABg", "background", 100, 0, 2);
206 mMCPosPVDCAisTg = std::make_unique<TH1F>("mcPosPVDCATg", "target V0", 40, 0, 15);
207 mMCPosPVDCAisBg = std::make_unique<TH1F>("mcPosPVDCABg", "background", 100, 0, 15);
208 mMCNegPVDCAisTg = std::make_unique<TH1F>("mcNegPVDCATg", "target V0", 40, 0, 15);
209 mMCNegPVDCAisBg = std::make_unique<TH1F>("mcNegPVDCABg", "background", 100, 0, 15);
210 mMCV0PVDCAisTg = std::make_unique<TH1F>("mcV0PVDCATg", "target V0", 40, 0, 15);
211 mMCV0PVDCAisBg = std::make_unique<TH1F>("mcV0PVDCABg", "background", 100, 0, 15);
212
213 mMCArmPodolisTg = std::make_unique<TH2F>("apPlotTg", "Armenteros-Podolanski;#alpha_{Arm};p_{T}^{Arm} (GeV/c)", 150, -2, 2, 150, 0, 0.5);
214 mMCArmPodolisBg = std::make_unique<TH2F>("apPlotBg", "Armenteros-Podolanski;#alpha_{Arm};p_{T}^{Arm} (GeV/c)", 150, -2, 2, 150, 0, 0.5);
215
216 mAvgClusSizeCEta = std::make_unique<THStack>("avgclussizeeta", "Average cluster size per track;pixels / cluster / track"); // auto-set axis ranges
217 float binWidth = (params.etaMax - params.etaMin) / (float)params.etaNBins;
218 mEtaBinUL.reserve(params.etaNBins);
219 for (int i = 0; i < params.etaNBins; i++) {
220 mEtaBinUL.emplace_back(params.etaMin + (binWidth * (i + 1)));
221 mAvgClusSizeCEtaVec.push_back(std::make_unique<TH1F>(Form("avgclussize%i", i), Form("%.2f < #eta < %.2f", mEtaBinUL[i] - binWidth, mEtaBinUL[i]), params.sizeNBins, 0, params.sizeMax));
222 mAvgClusSizeCEtaVec[i]->SetDirectory(nullptr);
223 mAvgClusSizeCEta->Add(mAvgClusSizeCEtaVec[i].get());
224 }
225
226 mMCRisTg->SetDirectory(nullptr);
227 mMCRisBg->SetDirectory(nullptr);
228 mMCCosPAisTg->SetDirectory(nullptr);
229 mMCCosPAisBg->SetDirectory(nullptr);
230 mMCDauDCAisTg->SetDirectory(nullptr);
231 mMCDauDCAisBg->SetDirectory(nullptr);
232 mMCPosACSisTg->SetDirectory(nullptr);
233 mMCPosACSisBg->SetDirectory(nullptr);
234 mMCNegACSisTg->SetDirectory(nullptr);
235 mMCNegACSisBg->SetDirectory(nullptr);
236 mMCPosPVDCAisTg->SetDirectory(nullptr);
237 mMCPosPVDCAisBg->SetDirectory(nullptr);
238 mMCNegPVDCAisTg->SetDirectory(nullptr);
239 mMCNegPVDCAisBg->SetDirectory(nullptr);
240 mMCV0PVDCAisTg->SetDirectory(nullptr);
241 mMCV0PVDCAisBg->SetDirectory(nullptr);
242
243 mMCArmPodolisTg->SetDirectory(nullptr);
244 mMCArmPodolisBg->SetDirectory(nullptr);
245
246 mOutputNtupleAll->SetDirectory(nullptr);
247 mOutputNtupleCut->SetDirectory(nullptr);
248
249 mMCR->Add(mMCRisTg.get());
250 mMCR->Add(mMCRisBg.get());
251 mMCCosPA->Add(mMCCosPAisTg.get());
252 mMCCosPA->Add(mMCCosPAisBg.get());
253 mMCPosACS->Add(mMCPosACSisTg.get());
254 mMCPosACS->Add(mMCPosACSisBg.get());
255 mMCNegACS->Add(mMCNegACSisTg.get());
256 mMCNegACS->Add(mMCNegACSisBg.get());
257 mMCDauDCA->Add(mMCDauDCAisTg.get());
258 mMCDauDCA->Add(mMCDauDCAisBg.get());
259 mMCPosPVDCA->Add(mMCPosPVDCAisTg.get());
260 mMCPosPVDCA->Add(mMCPosPVDCAisBg.get());
261 mMCNegPVDCA->Add(mMCNegPVDCAisTg.get());
262 mMCNegPVDCA->Add(mMCNegPVDCAisBg.get());
263 mMCV0PVDCA->Add(mMCV0PVDCAisTg.get());
264 mMCV0PVDCA->Add(mMCV0PVDCAisBg.get());
265}
266
267void AvgClusSizeStudy::setStyle()
268{
269 gStyle->SetPalette(kRainBow);
270 std::vector<int> colors = {1, 2, 3, 4, 6, 7, 41, 47};
271 std::vector<int> markers = {2, 3, 4, 5, 25, 26, 27, 28, 32};
272 for (int i = 0; i < mAvgClusSizeCEtaVec.size(); i++) {
273 mAvgClusSizeCEtaVec[i]->SetMarkerStyle(markers[i]);
274 mAvgClusSizeCEtaVec[i]->SetMarkerColor(colors[i]);
275 mAvgClusSizeCEtaVec[i]->SetLineColor(colors[i]);
276 }
277
278 mMCRisTg->SetLineColor(kRed);
279 mMCCosPAisTg->SetLineColor(kRed);
280 mMCPosACSisTg->SetLineColor(kRed);
281 mMCNegACSisTg->SetLineColor(kRed);
282 mMCDauDCAisTg->SetLineColor(kRed);
283 mMCPosPVDCAisTg->SetLineColor(kRed);
284 mMCNegPVDCAisTg->SetLineColor(kRed);
285 mMCV0PVDCAisTg->SetLineColor(kRed);
286}
287
289{
290 // auto geom = o2::its::GeometryTGeo::Instance();
292 recoData.collectData(pc, *mDataRequest.get());
293 updateTimeDependentParams(pc); // Make sure this is called after recoData.collectData, which may load some conditions
294 process(recoData);
295}
296
297void AvgClusSizeStudy::getClusterSizes(std::vector<int>& clusSizeVec, const gsl::span<const o2::itsmft::CompClusterExt> ITSclus, gsl::span<const unsigned char>::iterator& pattIt, const o2::itsmft::TopologyDictionary* mdict)
298{
299 for (unsigned int iClus{0}; iClus < ITSclus.size(); ++iClus) {
300 auto& clus = ITSclus[iClus];
301 auto pattID = clus.getPatternID();
302 int npix;
304
305 if (pattID == o2::itsmft::CompCluster::InvalidPatternID || mdict->isGroup(pattID)) {
306 patt.acquirePattern(pattIt);
307 npix = patt.getNPixels();
308 } else {
309 npix = mdict->getNpixels(pattID);
310 patt = mdict->getPattern(pattID);
311 }
312 clusSizeVec[iClus] = npix;
313 }
314}
315
316void AvgClusSizeStudy::loadData(o2::globaltracking::RecoContainer& recoData)
317{
318 mInputITSidxs = recoData.getITSTracksClusterRefs();
319 auto compClus = recoData.getITSClusters();
320 auto clusPatt = recoData.getITSClustersPatterns();
321 mClusterSizes.resize(compClus.size());
322 auto pattIt = clusPatt.begin();
323 getClusterSizes(mClusterSizes, compClus, pattIt, mDict);
324}
325
326void AvgClusSizeStudy::process(o2::globaltracking::RecoContainer& recoData)
327{
329 float dPosACS, dNegACS, dauDCA, cosPA, v0R, eta, dPospvDCA, dNegpvDCA, v0pvDCA, tgV0HypoMass, bgV0HypoMass, alphaArm, pT; // ACS=average cluster size
330 bool isMCTarget = false;
331 int targetPDGCode = 310;
332 TrackITS dPosRecoTrk, dNegRecoTrk; // daughter ITS tracks
333 DCA dPosDCA, dNegDCA, v0DCA; // DCA object for prong to primary vertex (DCA = o2::dataformats::DCA)
334 PVertex pv;
335
336 PID targetV0 = PID("K0");
337 PID backgroundV0 = PID("Lambda");
338 PID tgPos, tgNeg, bgPos, bgNeg;
339 if (params.targetV0 == "K0") {
340 bgPos = PID("Proton");
341 LOGP(info, "V0 target set to K0-short.");
342 } else if (params.targetV0 == "Lambda") {
343 targetV0 = PID("Lambda");
344 tgPos = PID("Proton");
345 backgroundV0 = PID("K0");
346 targetPDGCode = 3122;
347 LOGP(info, "V0 target set to Lambda.");
348 } else {
349 LOGP(warning, "Given V0 target not recognized, defaulting to K0-short.");
350 bgPos = PID("Proton");
351 }
352
353 // Variables for MC analysis
354 gsl::span<const o2::MCCompLabel> mcLabels;
355 o2::MCCompLabel dPosLab, dNegLab;
356 const o2::MCTrack *V0mcTrk, *dPosMCTrk, *dNegMCTrk;
357 int V0PdgCode, mPosTrkId, mNegTrkId;
358
359 loadData(recoData);
360 auto V0s = recoData.getV0s();
361 auto V0sIdx = recoData.getV0sIdx();
362 size_t nV0s = V0sIdx.size();
363 if (nV0s && nV0s != V0s.size()) {
364 LOGP(fatal, "This data has not secondary vertices kinematics filled");
365 }
366
367 LOGP(info, "Found {} reconstructed V0s.", nV0s);
368 LOGP(info, "Found {} ITS tracks.", recoData.getITSTracks().size());
369 LOGP(info, "Found {} ROFs.", recoData.getITSTracksROFRecords().size());
370 if (mUseMC) {
371 mcLabels = recoData.getITSTracksMCLabels();
372 LOGP(info, "Found {} labels.", mcLabels.size());
373 }
374
375 for (size_t iv = 0; iv < nV0s; iv++) {
376 auto v0 = V0s[iv];
377 const auto& v0Idx = V0sIdx[iv];
378 dPosRecoTrk = recoData.getITSTrack(v0Idx.getProngID(0)); // cannot use v0.getProng() since it returns TrackParCov not TrackITS
379 dNegRecoTrk = recoData.getITSTrack(v0Idx.getProngID(1));
380
381 pv = recoData.getPrimaryVertex(v0Idx.getVertexID()); // extract primary vertex
382 dPosRecoTrk.propagateToDCA(pv, params.b, &dPosDCA); // calculate and store DCA objects for both prongs
383 dNegRecoTrk.propagateToDCA(pv, params.b, &dNegDCA);
384 v0.propagateToDCA(pv, params.b, &v0DCA);
385
386 dPospvDCA = std::sqrt(dPosDCA.getR2()); // extract DCA distance from DCA object
387 dNegpvDCA = std::sqrt(dNegDCA.getR2());
388 v0pvDCA = std::sqrt(v0DCA.getR2());
389
390 eta = v0.getEta();
391 dauDCA = v0.getDCA();
392 cosPA = v0.getCosPA();
393 v0R = std::sqrt(v0.calcR2()); // gives distance from pvertex to origin? in centimeters (?) NOTE: unsure if this is to the primary vertex or to origin
394
395 dPosACS = getAverageClusterSize(&dPosRecoTrk);
396 dNegACS = getAverageClusterSize(&dNegRecoTrk);
397
398 tgV0HypoMass = calcV0HypoMass(v0, tgPos, tgNeg);
399 bgV0HypoMass = calcV0HypoMass(v0, bgPos, bgNeg);
400 calcAPVars(v0, &alphaArm, &pT);
401
402 if (mUseMC) { // check whether queried V0 is the target V0 in MC, and fill the cut validation plots
403 V0PdgCode = 0;
404 isMCTarget = false;
405 dPosLab = mcLabels[v0Idx.getProngID(0)]; // extract MC labels for the prongs
406 dNegLab = mcLabels[v0Idx.getProngID(1)];
407 if (!dPosLab.isValid() || !dNegLab.isValid()) {
408 LOGP(debug, "Daughter MCCompLabel not valid: {}(+) and {}(-). Skipping.", dPosLab.isValid(), dNegLab.isValid());
409 nNotValid++;
410 } else {
411 dPosMCTrk = mKineReader->getTrack(dPosLab);
412 dNegMCTrk = mKineReader->getTrack(dNegLab);
413 if (dPosMCTrk == nullptr || dNegMCTrk == nullptr) {
414 LOGP(debug, "Nullptr found: {}(+) and {}(-). Skipping.", (void*)dPosMCTrk, (void*)dNegMCTrk);
415 nNullptrs++;
416 } else {
417 mPosTrkId = dPosMCTrk->getMotherTrackId();
418 mNegTrkId = dNegMCTrk->getMotherTrackId();
419 LOGP(debug, "Daughter PDG codes: {}(+) and {}(-)", dPosMCTrk->GetPdgCode(), dNegMCTrk->GetPdgCode());
420 if (mPosTrkId != mNegTrkId || mPosTrkId == -1 || mNegTrkId == -1) {
421 LOGP(debug, "Mother track ID mismatch or default -1: {}(+) and {}(-). Skipping.", mPosTrkId, mNegTrkId);
422 nMotherIDMismatch++;
423 } else {
424 if (dNegLab.getEventID() != dPosLab.getEventID()) {
425 LOGP(debug, "Daughter EvID mismatch: {}(+) and {}(-). Skipping.", dPosLab.getEventID(), dNegLab.getEventID());
426 nEvIDMismatch++;
427 } else {
428 V0mcTrk = mKineReader->getTrack(dNegLab.getEventID(), mPosTrkId); // assume daughter MCTracks are in same event as V0 MCTrack
429 V0PdgCode = V0mcTrk->GetPdgCode();
430 if (V0PdgCode == targetPDGCode) {
431 isMCTarget = true;
432 nTargetV0++;
433 } else {
434 nNotTargetV0++;
435 }
436 }
437 }
438 }
439 }
440 if (isMCTarget) {
441 mMCCosPAisTg->Fill(cosPA);
442 mMCDauDCAisTg->Fill(dauDCA);
443 mMCRisTg->Fill(v0R);
444 mMCPosPVDCAisTg->Fill(dPospvDCA);
445 mMCNegPVDCAisTg->Fill(dNegpvDCA);
446 mMCPosACSisTg->Fill(dPosACS);
447 mMCNegACSisTg->Fill(dNegACS);
448 mMCV0PVDCAisTg->Fill(v0pvDCA);
449 mMCArmPodolisTg->Fill(alphaArm, pT);
450 } else {
451 mMCCosPAisBg->Fill(cosPA);
452 mMCDauDCAisBg->Fill(dauDCA);
453 mMCRisBg->Fill(v0R);
454 mMCPosPVDCAisBg->Fill(dPospvDCA);
455 mMCNegPVDCAisBg->Fill(dNegpvDCA);
456 mMCPosACSisBg->Fill(dPosACS);
457 mMCNegACSisBg->Fill(dNegACS);
458 mMCV0PVDCAisBg->Fill(v0pvDCA);
459 mMCArmPodolisBg->Fill(alphaArm, pT);
460 }
461 }
462
463 mOutputNtupleAll->Fill(dPosACS, dNegACS, cosPA, v0R, eta, dauDCA, dPospvDCA, dNegpvDCA, v0pvDCA, alphaArm, pT, calcV0HypoMass(v0, PID::Pion, PID::Pion), calcV0HypoMass(v0, PID::Proton, PID::Pion), calcV0HypoMass(v0, PID::Pion, PID::Proton), (float)V0PdgCode);
464 if ((cosPA > params.cosPAmin || params.disableCosPA) && (v0R < params.Rmax || params.disableRmax) && (v0R > params.Rmin || params.disableRmin) && (dauDCA < params.prongDCAmax || params.disableProngDCAmax) && (dPospvDCA > params.dauPVDCAmin || params.disableDauPVDCAmin) && (dNegpvDCA > params.dauPVDCAmin || params.disableDauPVDCAmin) && (v0pvDCA < params.v0PVDCAmax || params.disableV0PVDCAmax) && (abs(bgV0HypoMass - backgroundV0.getMass()) > params.bgV0window || params.disableMassHypoth) && (abs(tgV0HypoMass - targetV0.getMass()) < params.tgV0window || params.disableMassHypoth)) {
465 mOutputNtupleCut->Fill(dPosACS, dNegACS, cosPA, v0R, eta, dauDCA, dPospvDCA, dNegpvDCA, v0pvDCA, alphaArm, pT, calcV0HypoMass(v0, PID::Pion, PID::Pion), calcV0HypoMass(v0, PID::Proton, PID::Pion), calcV0HypoMass(v0, PID::Pion, PID::Proton), (float)V0PdgCode);
466 if (eta > params.etaMin && eta < params.etaMax) {
467 fillEtaBin(eta, dPosACS, 0);
468 fillEtaBin(eta, dNegACS, 0);
469 } else {
470 nV0OutOfEtaRange++;
471 }
472 }
473 }
474
475 if (mUseMC) {
476 LOGP(info, "MONTE CARLO OVERALL STATISTICS: {} total V0s, {} nonvalid daughter labels, {} nullptrs, {} motherID mismatches, {} evID mismatches, {} matching target, {} not matching target", V0s.size(), nNotValid, nNullptrs, nMotherIDMismatch, nEvIDMismatch, nTargetV0, nNotTargetV0);
477 int nPrimaryTargetV0 = 0;
478 for (auto& mcTrk : mMCTracks) { // search through all MC tracks to find primary target V0s, whether reconstructed or not
479 if (mcTrk.GetPdgCode() == targetPDGCode && mcTrk.isPrimary()) {
480 nPrimaryTargetV0++;
481 }
482 }
483 LOGP(info, "MONTE CARLO OVERALL STATISTICS: {} MC target V0s (isPrimary) found in MC tracks out of {} total MC tracks", nPrimaryTargetV0, mMCTracks.size());
484 }
485 LOGP(info, "{} V0s out of eta range ({}, {})", nV0OutOfEtaRange, params.etaMin, params.etaMax);
486}
487
488float AvgClusSizeStudy::getAverageClusterSize(TrackITS* daughter)
489{
490 int totalSize{0};
491 auto firstClus = daughter->getFirstClusterEntry();
492 auto ncl = daughter->getNumberOfClusters();
493 for (int icl = 0; icl < ncl; icl++) {
494 totalSize += mClusterSizes[mInputITSidxs[firstClus + icl]];
495 }
496 return (float)totalSize / (float)ncl;
497}
498
499float AvgClusSizeStudy::calcV0HypoMass(const V0& v0, PID hypothPIDPos, PID hypothPIDNeg)
500{
501 // Mass hypothesis calculation; taken from o2::strangeness_tracking::StrangenessTracker::calcMotherMass()
502 std::array<float, 3> pPos, pNeg, pV0;
503 v0.getProng(0).getPxPyPzGlo(pPos);
504 v0.getProng(1).getPxPyPzGlo(pNeg);
505 v0.getPxPyPzGlo(pV0);
506 double m2Pos = PID::getMass2(hypothPIDPos);
507 double m2Neg = PID::getMass2(hypothPIDNeg);
508 double p2Pos = (pPos[0] * pPos[0]) + (pPos[1] * pPos[1]) + (pPos[2] * pPos[2]);
509 double p2Neg = (pNeg[0] * pNeg[0]) + (pNeg[1] * pNeg[1]) + (pNeg[2] * pNeg[2]);
510 double ePos = std::sqrt(p2Pos + m2Pos), eNeg = std::sqrt(p2Neg + m2Neg);
511 double e2V0 = (ePos + eNeg) * (ePos + eNeg);
512 double pxV0 = (pPos[0] + pNeg[0]);
513 double pyV0 = (pPos[1] + pNeg[1]);
514 double pzV0 = (pPos[2] + pNeg[2]);
515 double p2V0 = (pxV0 * pxV0) + (pyV0 * pyV0) + (pzV0 * pzV0);
516 return (float)std::sqrt(e2V0 - p2V0);
517}
518
519void AvgClusSizeStudy::calcAPVars(const V0& v0, float* alphaArm, float* pT)
520{
521 // Calculation of the Armenteros-Podolanski variables
522 std::array<float, 3> pV0, pPos, pNeg;
523 v0.getProng(0).getPxPyPzGlo(pPos);
524 v0.getProng(1).getPxPyPzGlo(pNeg);
525 v0.getPxPyPzGlo(pV0);
526 double p2V0 = pV0[0] * pV0[0] + pV0[1] * pV0[1] + pV0[2] * pV0[2];
527 double qNeg = pNeg[0] * pV0[0] + pNeg[1] * pV0[1] + pNeg[2] * pV0[2];
528 double qPos = pPos[0] * pV0[0] + pPos[1] * pV0[1] + pPos[2] * pV0[2];
529 *alphaArm = (float)(qPos - qNeg) / (qPos + qNeg);
530 double p2Pos = pPos[0] * pPos[0] + pPos[1] * pPos[1] + pPos[2] * pPos[2];
531 *pT = (float)std::sqrt(p2Pos - ((qPos * qPos) / p2V0));
532};
533
534void AvgClusSizeStudy::fillEtaBin(float eta, float clusSize, int i)
535{
536 if (eta < mEtaBinUL[i]) {
537 mAvgClusSizeCEtaVec[i]->Fill(clusSize);
538 } else {
539 fillEtaBin(eta, clusSize, i + 1);
540 }
541}
542
543void AvgClusSizeStudy::updateTimeDependentParams(ProcessingContext& pc)
544{
546 static bool initOnceDone = false;
547 if (!initOnceDone) { // this param need to be queried only once
548 initOnceDone = true;
551 }
552}
553
554void AvgClusSizeStudy::saveHistograms()
555{
556 mDBGOut.reset();
557 TFile fout(mOutName.c_str(), "RECREATE");
558
559 fout.WriteTObject(mOutputNtupleAll.get());
560 fout.WriteTObject(mOutputNtupleCut.get());
561 fout.WriteTObject(mMCR.get());
562 fout.WriteTObject(mMCCosPA.get());
563 fout.WriteTObject(mMCPosACS.get());
564 fout.WriteTObject(mMCNegACS.get());
565 fout.WriteTObject(mMCDauDCA.get());
566 fout.WriteTObject(mMCPosPVDCA.get());
567 fout.WriteTObject(mMCNegPVDCA.get());
568 fout.WriteTObject(mMCV0PVDCA.get());
569 fout.WriteTObject(mMCArmPodolisTg.get());
570 fout.WriteTObject(mMCArmPodolisBg.get());
571
572 fout.WriteTObject(mAvgClusSizeCEta.get());
573 fout.Close();
574
575 LOGP(info, "Stored histograms into {}", mOutName.c_str());
576}
577
578void AvgClusSizeStudy::plotHistograms()
579{
580 if (mUseMC) {
581 TCanvas* cMCR = new TCanvas();
582 mMCR->Draw("nostack");
583 cMCR->BuildLegend();
584 cMCR->Print("mcR.png");
585 TCanvas* cMCCosPA = new TCanvas();
586 mMCCosPA->Draw("nostack");
587 cMCCosPA->BuildLegend();
588 cMCCosPA->Print("mcCosPA.png");
589 TCanvas* cMCPosACS = new TCanvas();
590 mMCPosACS->Draw("nostack");
591 cMCPosACS->BuildLegend();
592 cMCPosACS->Print("mcPosACS.png");
593 TCanvas* cMCNegACS = new TCanvas();
594 mMCNegACS->Draw("nostack");
595 cMCNegACS->BuildLegend();
596 cMCNegACS->Print("mcNegACS.png");
597 TCanvas* cMCDauDCA = new TCanvas();
598 mMCDauDCA->Draw("nostack");
599 cMCDauDCA->BuildLegend();
600 cMCDauDCA->Print("mcDauDCA.png");
601 TCanvas* cMCPosPVDCA = new TCanvas();
602 mMCPosPVDCA->Draw("nostack");
603 cMCPosPVDCA->BuildLegend();
604 cMCPosPVDCA->Print("mcPosPVDCA.png");
605 TCanvas* cMCNegPVDCA = new TCanvas();
606 mMCNegPVDCA->Draw("nostack");
607 cMCNegPVDCA->BuildLegend();
608 cMCNegPVDCA->Print("mcNegPVDCA.png");
609 TCanvas* cMCV0PVDCA = new TCanvas();
610 mMCV0PVDCA->Draw("nostack");
611 cMCV0PVDCA->BuildLegend();
612 cMCV0PVDCA->Print("mcV0PVDCA.png");
613 TCanvas* cMCArmPodolTg = new TCanvas();
614 mMCArmPodolisTg->Draw("COLZ");
615 cMCArmPodolTg->Print("mcArmPodolTg.png");
616 TCanvas* cMCArmPodolBg = new TCanvas();
617 mMCArmPodolisBg->Draw("COLZ");
618 cMCArmPodolBg->Print("mcArmPodolBg.png");
619 }
620
621 TCanvas* c10 = new TCanvas();
622 mAvgClusSizeCEta->Draw("P NOSTACK");
623 c10->BuildLegend(0.6, 0.6, 0.8, 0.8);
624 c10->Print("clusSizeEta.png");
625}
626
628{
630 setStyle();
631 saveHistograms();
632 if (params.generatePlots) {
633 plotHistograms();
634 }
635}
636
638{
640 return;
641 }
642 // o2::base::GRPGeomHelper::instance().finaliseCCDB(matcher, obj);
643 if (matcher == ConcreteDataMatcher("ITS", "CLUSDICT", 0)) {
645 return;
646 }
647}
648
649DataProcessorSpec getAvgClusSizeStudy(mask_t srcTracksMask, mask_t srcClustersMask, bool useMC, std::shared_ptr<o2::steer::MCKinematicsReader> kineReader)
650{
651 std::vector<OutputSpec> outputs;
652 auto dataRequest = std::make_shared<DataRequest>();
653 dataRequest->requestTracks(srcTracksMask, useMC);
654 dataRequest->requestClusters(srcClustersMask, useMC);
655 dataRequest->requestSecondaryVertices(useMC);
656 dataRequest->requestPrimaryVertices(useMC); // NOTE: may be necessary to use requestPrimaryVerticesTMP()...
657 // dataRequest->requestPrimaryVerticesTMP(useMC);
658
659 auto ggRequest = std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
660 true, // GRPECS=true
661 false, // GRPLHCIF
662 false, // GRPMagField
663 false, // askMatLUT
665 dataRequest->inputs,
666 true);
667 return DataProcessorSpec{
668 "its-study-AvgClusSize",
669 dataRequest->inputs,
670 outputs,
671 AlgorithmSpec{adaptFromTask<AvgClusSizeStudy>(dataRequest, ggRequest, useMC, kineReader)},
672 Options{}};
673}
674} // namespace study
675} // namespace its
676} // namespace o2
Wrapper container for different reconstructed object types.
particle ids, masses, names class definition
Base track model for the Barrel, params only, w/o covariance.
std::ostringstream debug
int32_t i
Helper for geometry and GRP related CCDB requests.
Header of the General Run Parameters object.
Global index for barrel track: provides provenance (detectors combination), index in respective array...
Definition of the GeometryTGeo class.
Definition of the MCTrack class.
Definition of the ITS track.
int getEventID() const
bool isValid() const
Definition MCCompLabel.h:75
Int_t GetPdgCode() const
Accessors.
Definition MCTrack.h:72
Int_t getMotherTrackId() const
Definition MCTrack.h:73
void checkUpdates(o2::framework::ProcessingContext &pc)
static GRPGeomHelper & instance()
void setRequest(std::shared_ptr< GRPGeomRequest > req)
o2::gpu::gpustd::bitset< 32 > mask_t
static GeometryTGeo * Instance()
void fillMatrixCache(int mask) override
int getFirstClusterEntry() const
Definition TrackITS.h:66
AvgClusSizeStudy(std::shared_ptr< DataRequest > dr, std::shared_ptr< o2::base::GRPGeomRequest > gr, bool isMC, std::shared_ptr< o2::steer::MCKinematicsReader > kineReader)
void init(InitContext &ic) final
void setClusterDictionary(const o2::itsmft::TopologyDictionary *d)
void finaliseCCDB(ConcreteDataMatcher &, void *) final
void run(ProcessingContext &) final
void endOfStream(EndOfStreamContext &) final
This is invoked whenever we have an EndOfStream event.
void acquirePattern(iterator &pattIt)
int getNPixels() const
Returns the number of fired pixels.
static constexpr unsigned short InvalidPatternID
Definition CompCluster.h:46
const ClusterPattern & getPattern(int n) const
Returns the pattern of the topology.
int getNpixels(int n) const
Returns the number of fired pixels of the n_th element.
bool isGroup(int n) const
Returns true if the element corresponds to a group of rare topologies.
static constexpr ID Pion
Definition PID.h:96
static constexpr ID Proton
Definition PID.h:98
GLenum const GLfloat * params
Definition glcorearb.h:272
GLfloat v0
Definition glcorearb.h:811
Defining PrimaryVertex explicitly as messageable.
std::vector< ConfigParamSpec > Options
auto get(const std::byte *buffer, size_t=0)
Definition DataHeader.h:454
o2::dataformats::DCA DCA
o2::track::TrackParCov Track
o2::track::PID PID
o2::framework::DataProcessorSpec getAvgClusSizeStudy(mask_t srcTracksMask, mask_t srcClustersMask, bool useMC, std::shared_ptr< o2::steer::MCKinematicsReader > kineReader)
TrackParCovF TrackParCov
Definition Track.h:33
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
const o2::dataformats::PrimaryVertex & getPrimaryVertex(int i) const
const o2::its::TrackITS & getITSTrack(GTrackID gid) const
static constexpr int T2L
Definition Cartesian.h:55
static constexpr int T2GRot
Definition Cartesian.h:57
static constexpr int T2G
Definition Cartesian.h:56