Project
Loading...
Searching...
No Matches
TrackMCStudy.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#include <vector>
13#include <TStopwatch.h>
21#include "ITStracking/IOUtils.h"
55#include "GPUO2InterfaceRefit.h"
56#include "GPUParam.h"
57#include "GPUParam.inc"
58#include "MathUtils/fit.h"
59#include "TPCFastTransformPOD.h"
60#include <TRandom.h>
61#include <map>
62#include <unordered_map>
63#include <array>
64#include <utility>
65#include <gsl/span>
66
67// workflow to study relation of reco tracks to MCTruth
68// o2-trackmc-study-workflow --device-verbosity 3 -b --run
69
70namespace o2::trackstudy
71{
72
73using namespace o2::framework;
76
80using VTIndexV = std::pair<int, o2::dataformats::VtxTrackIndex>;
83
85
86class TrackMCStudy final : public Task
87{
88 public:
89 TrackMCStudy(std::shared_ptr<DataRequest> dr, std::shared_ptr<o2::base::GRPGeomRequest> gr, GTrackID::mask_t src, bool checkSV)
90 : mDataRequest(dr), mGGCCDBRequest(gr), mTracksSrc(src), mCheckSV(checkSV) {}
91 ~TrackMCStudy() final = default;
92 void init(InitContext& ic) final;
93 void run(ProcessingContext& pc) final;
94 void endOfStream(EndOfStreamContext& ec) final;
95 void finaliseCCDB(ConcreteDataMatcher& matcher, void* obj) final;
96 void process(const o2::globaltracking::RecoContainer& recoData);
97
98 private:
99 void processTPCTrackRefs();
100 void processITSTracks(const o2::globaltracking::RecoContainer& recoData);
101 void loadTPCOccMap(const o2::globaltracking::RecoContainer& recoData);
102 void fillMCClusterInfo(const o2::globaltracking::RecoContainer& recoData);
103 void prepareITSData(const o2::globaltracking::RecoContainer& recoData);
104 bool processMCParticle(int src, int ev, int trid);
105 bool addMCParticle(const MCTrack& mctr, const o2::MCCompLabel& lb, TParticlePDG* pPDG = nullptr);
106 bool acceptMCCharged(const MCTrack& tr, const o2::MCCompLabel& lb, int followDec = -1);
107 bool propagateToRefX(o2::track::TrackParCov& trcTPC, o2::track::TrackParCov& trcITS);
108 bool refitV0(int i, o2::dataformats::V0& v0, const o2::globaltracking::RecoContainer& recoData);
109 void updateTimeDependentParams(ProcessingContext& pc);
110 float getDCAYCut(float pt) const;
111
112 const std::vector<o2::MCTrack>* mCurrMCTracks = nullptr;
113 TVector3 mCurrMCVertex;
114 o2::tpc::VDriftHelper mTPCVDriftHelper{};
115 const o2::gpu::TPCFastTransformPOD* mTPCCorrMaps{nullptr};
116 std::shared_ptr<DataRequest> mDataRequest;
117 std::shared_ptr<o2::base::GRPGeomRequest> mGGCCDBRequest;
118 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOut;
119 std::vector<float> mTBinClOcc;
120 std::vector<float> mTBinClOccHist; //< original occupancy
121 std::vector<long> mIntBC;
122 std::vector<float> mTPCOcc;
123 std::vector<int> mITSOcc; //< N ITS clusters in the ROF containing collision
124 std::vector<o2::BaseCluster<float>> mITSClustersArray;
125 const o2::itsmft::TopologyDictionary* mITSDict = nullptr;
126
127 bool mCheckSV = false; //< check SV binding (apart from prongs availability)
128 bool mRecProcStage = false; //< flag that the MC particle was added only at the stage of reco tracks processing
129 int mNTPCOccBinLength = 0;
130 float mNTPCOccBinLengthInv = -1.f;
131 int mVerbose = 0;
132 float mITSTimeBiasMUS = 0.f;
133 float mITSROFrameLengthMUS = 0.f;
134 float mTPCTBinMUS = 0.;
135
136 int mNCheckDecays = 0;
137
138 GTrackID::mask_t mTracksSrc{};
139 o2::steer::MCKinematicsReader mcReader; // reader of MC information
140 std::vector<int> mITSROF;
141 std::vector<TBracket> mITSROFBracket;
142 std::vector<o2::MCCompLabel> mDecProdLblPool; // labels of decay products to watch, added to MC map
143 std::vector<MCVertex> mMCVtVec{};
144
145 struct DecayRef {
146 o2::MCCompLabel mother{};
147 o2::track::TrackPar parent{};
148 int pdg = 0;
149 int daughterFirst = -1;
150 int daughterLast = -1;
151 int foundSVID = -1;
152 };
153 std::vector<std::vector<DecayRef>> mDecaysMaps; // for every parent particle to watch, store its label and entries of 1st/last decay product labels in mDecProdLblPool
154 std::unordered_map<o2::MCCompLabel, TrackFamily> mSelMCTracks;
155 std::unordered_map<o2::MCCompLabel, std::pair<int, int>> mSelTRefIdx;
156 std::vector<o2::track::TrackPar> mSelTRefs;
158 static constexpr float MaxSnp = 0.9; // max snp of ITS or TPC track at xRef to be matched
159};
160
162{
164 mcReader.initFromDigitContext("collisioncontext.root");
165
166 mDBGOut = std::make_unique<o2::utils::TreeStreamRedirector>("trackMCStudy.root", "recreate");
167 mVerbose = ic.options().get<int>("device-verbosity");
168
170 for (int id = 0; id < sizeof(params.decayPDG) / sizeof(int); id++) {
171 if (params.decayPDG[id] < 0) {
172 break;
173 }
174 mNCheckDecays++;
175 }
176 mDecaysMaps.resize(mNCheckDecays);
177}
178
180{
182 for (int i = 0; i < mNCheckDecays; i++) {
183 mDecaysMaps[i].clear();
184 }
185 mDecProdLblPool.clear();
186 mMCVtVec.clear();
187 mCurrMCTracks = nullptr;
188
189 recoData.collectData(pc, *mDataRequest.get()); // select tracks of needed type, with minimal cuts, the real selected will be done in the vertexer
190 updateTimeDependentParams(pc); // Make sure this is called after recoData.collectData, which may load some conditions
191 mRecProcStage = false;
192 process(recoData);
193}
194
195void TrackMCStudy::updateTimeDependentParams(ProcessingContext& pc)
196{
198 mTPCVDriftHelper.extractCCDBInputs(pc);
199 auto const& raw = pc.inputs().get<const char*>("corrMap");
200 mTPCCorrMaps = &o2::gpu::TPCFastTransformPOD::get(raw);
201 static bool initOnceDone = false;
202 if (!initOnceDone) { // this params need to be queried only once
203 initOnceDone = true;
205 mITSROFrameLengthMUS = o2::base::GRPGeomHelper::instance().getGRPECS()->isDetContinuousReadOut(o2::detectors::DetID::ITS) ? alpParamsITS.roFrameLengthInBC * o2::constants::lhc::LHCBunchSpacingMUS : alpParamsITS.roFrameLengthTrig * 1.e-3;
206 LOGP(info, "VertexTrackMatcher ITSROFrameLengthMUS:{}", mITSROFrameLengthMUS);
207
209 mTPCTBinMUS = elParam.ZbinWidth;
211 if (mCheckSV) {
212 const auto& svparam = o2::vertexing::SVertexerParams::Instance();
213 mFitterV0.setBz(o2::base::Propagator::Instance()->getNominalBz());
214 mFitterV0.setUseAbsDCA(svparam.useAbsDCA);
215 mFitterV0.setPropagateToPCA(false);
216 mFitterV0.setMaxR(svparam.maxRIni);
217 mFitterV0.setMinParamChange(svparam.minParamChange);
218 mFitterV0.setMinRelChi2Change(svparam.minRelChi2Change);
219 mFitterV0.setMaxDZIni(svparam.maxDZIni);
220 mFitterV0.setMaxDXYIni(svparam.maxDXYIni);
221 mFitterV0.setMaxChi2(svparam.maxChi2);
222 mFitterV0.setMatCorrType(o2::base::Propagator::MatCorrType(svparam.matCorr));
223 mFitterV0.setUsePropagator(svparam.usePropagator);
224 mFitterV0.setRefitWithMatCorr(svparam.refitWithMatCorr);
225 mFitterV0.setMaxStep(svparam.maxStep);
226 mFitterV0.setMaxSnp(svparam.maxSnp);
227 mFitterV0.setMinXSeed(svparam.minXSeed);
228 }
229 }
230}
231
233{
234 constexpr float SQRT12Inv = 0.288675f;
236 auto pvvec = recoData.getPrimaryVertices();
237 auto pvvecLbl = recoData.getPrimaryVertexMCLabels();
238 auto trackIndex = recoData.getPrimaryVertexMatchedTracks(); // Global ID's for associated tracks
239 auto vtxRefs = recoData.getPrimaryVertexMatchedTrackRefs(); // references from vertex to these track IDs
240 auto prop = o2::base::Propagator::Instance();
241 int nv = vtxRefs.size();
242 float vdriftTB = mTPCVDriftHelper.getVDriftObject().getVDrift() * o2::tpc::ParameterElectronics::Instance().ZbinWidth; // VDrift expressed in cm/TimeBin
243 float itsBias = 0.5 * mITSROFrameLengthMUS + o2::itsmft::DPLAlpideParam<o2::detectors::DetID::ITS>::Instance().roFrameBiasInBC * o2::constants::lhc::LHCBunchSpacingMUS; // ITS time is supplied in \mus as beginning of ROF
244
245 prepareITSData(recoData);
246 loadTPCOccMap(recoData);
247 auto getITSPatt = [&](GTrackID gid, uint8_t& ncl) {
248 int8_t patt = 0;
249 if (gid.getSource() == VTIndex::ITSAB) {
250 const auto& itsTrf = recoData.getITSABRefs()[gid];
251 ncl = itsTrf.getNClusters();
252 for (int il = 0; il < 7; il++) {
253 if (itsTrf.hasHitOnLayer(il)) {
254 patt |= 0x1 << il;
255 }
256 }
257 patt |= 0x1 << 7;
258 } else {
259 const auto& itsTr = recoData.getITSTrack(gid);
260 for (int il = 0; il < 7; il++) {
261 if (itsTr.hasHitOnLayer(il)) {
262 patt |= 0x1 << il;
263 ncl++;
264 }
265 }
266 }
267 return patt;
268 };
269
270 auto fillTPCClusterInfo = [&recoData](const o2::tpc::TrackTPC& trc, RecTrack& tref) {
271 if (recoData.inputsTPCclusters) {
272 uint8_t clSect = 0, clRow = 0, lowestR = -1;
273 uint32_t clIdx = 0;
274 const auto clRefs = recoData.getTPCTracksClusterRefs();
275 const auto tpcClusAcc = recoData.getTPCClusters();
276 const auto shMap = recoData.clusterShMapTPC;
277 for (int ic = 0; ic < trc.getNClusterReferences(); ic++) { // outside -> inside ordering, but on the sector boundaries backward jumps are possible
278 trc.getClusterReference(clRefs, ic, clSect, clRow, clIdx);
279 if (clRow < lowestR) {
280 tref.rowCountTPC++;
281 lowestR = clRow;
282 }
283 unsigned int absoluteIndex = tpcClusAcc.clusterOffset[clSect][clRow] + clIdx;
284 if (shMap[absoluteIndex] & o2::gpu::GPUTPCGMMergedTrackHit::flagShared) {
285 tref.nClTPCShared++;
286 }
287 }
288 tref.lowestPadRow = lowestR;
289 const auto& clus = tpcClusAcc.clusters[clSect][clRow][clIdx];
290 int padFromEdge = int(clus.getPad()), npads = o2::gpu::GPUTPCGeometry::NPads(clRow);
291 if (padFromEdge > npads / 2) {
292 padFromEdge = npads - 1 - padFromEdge;
293 }
294 tref.padFromEdge = uint8_t(padFromEdge);
295 trc.getClusterReference(clRefs, 0, clSect, clRow, clIdx);
296 tref.rowMaxTPC = clRow;
297 }
298 };
299
300 auto flagTPCClusters = [&recoData](const o2::tpc::TrackTPC& trc, o2::MCCompLabel lbTrc) {
301 if (recoData.inputsTPCclusters) {
302 const auto clRefs = recoData.getTPCTracksClusterRefs();
303 const auto* TPCClMClab = recoData.inputsTPCclusters->clusterIndex.clustersMCTruth;
304 const auto& TPCClusterIdxStruct = recoData.inputsTPCclusters->clusterIndex;
305 for (int ic = 0; ic < trc.getNClusterReferences(); ic++) {
306 uint8_t clSect = 0, clRow = 0;
307 uint32_t clIdx = 0;
308 trc.getClusterReference(clRefs, ic, clSect, clRow, clIdx);
309 auto labels = TPCClMClab->getLabels(clIdx + TPCClusterIdxStruct.clusterOffset[clSect][clRow]);
310 for (auto& lbl : labels) {
311 if (lbl == lbTrc) {
312 const_cast<o2::MCCompLabel&>(lbl).setFakeFlag(true); // actually, in this way we are flagging that this cluster was correctly attached
313 break;
314 }
315 }
316 }
317 }
318 };
319
320 {
321 const auto* digconst = mcReader.getDigitizationContext();
322 const auto& mcEvRecords = digconst->getEventRecords(false);
323 int ITSTimeBias = o2::itsmft::DPLAlpideParam<o2::detectors::DetID::ITS>::Instance().roFrameBiasInBC;
324 int ITSROFLen = o2::itsmft::DPLAlpideParam<o2::detectors::DetID::ITS>::Instance().roFrameLengthInBC;
325 unsigned int rofCount = 0;
326 const auto ITSClusROFRec = recoData.getITSClustersROFRecords();
327 for (const auto& mcIR : mcEvRecords) {
328 long tbc = mcIR.differenceInBC(recoData.startIR);
329 auto& mcVtx = mMCVtVec.emplace_back();
330 mcVtx.ts = tbc * o2::constants::lhc::LHCBunchSpacingMUS + mcIR.getTimeOffsetWrtBC() * 1e-3;
331 mcVtx.ID = mIntBC.size();
332 mIntBC.push_back(tbc);
333 int occBin = tbc / 8 * mNTPCOccBinLengthInv;
334 mTPCOcc.push_back(occBin < 0 ? mTBinClOcc[0] : (occBin >= mTBinClOcc.size() ? mTBinClOcc.back() : mTBinClOcc[occBin]));
335 // fill ITS occupancy
336 long gbc = mcIR.toLong();
337 while (rofCount < ITSClusROFRec.size()) {
338 long rofbcMin = ITSClusROFRec[rofCount].getBCData().toLong() + ITSTimeBias, rofbcMax = rofbcMin + ITSROFLen;
339 if (gbc < rofbcMin) { // IRs and ROFs are sorted, so this IR is prior of all ROFs
340 mITSOcc.push_back(0);
341 } else if (gbc < rofbcMax) {
342 mITSOcc.push_back(ITSClusROFRec[rofCount].getNEntries());
343 } else {
344 rofCount++; // test next ROF
345 continue;
346 }
347 break;
348 }
349 if (mNTPCOccBinLengthInv > 0.f) {
350 mcVtx.occTPCV.resize(params.nOccBinsDrift);
351 int grp = TMath::Max(1, TMath::Nint(params.nTBPerOccBin * mNTPCOccBinLengthInv));
352 for (int ib = 0; ib < params.nOccBinsDrift; ib++) {
353 float smb = 0;
354 int tbs = occBin + TMath::Nint(ib * params.nTBPerOccBin * mNTPCOccBinLengthInv);
355 for (int ig = 0; ig < grp; ig++) {
356 if (tbs >= 0 && tbs < int(mTBinClOccHist.size())) {
357 smb += mTBinClOccHist[tbs];
358 }
359 tbs++;
360 }
361 mcVtx.occTPCV[ib] = smb;
362 }
363 }
364 if (rofCount >= ITSClusROFRec.size()) {
365 mITSOcc.push_back(0); // IR after the last ROF
366 }
367 }
368 }
369 // collect interesting MC particle (tracks and parents)
370 int curSrcMC = 0, curEvMC = 0;
371 for (curSrcMC = 0; curSrcMC < (int)mcReader.getNSources(); curSrcMC++) {
372 if (mVerbose > 1) {
373 LOGP(info, "Source {}", curSrcMC);
374 }
375 int nev = mcReader.getNEvents(curSrcMC);
376 bool okAccVtx = true;
377 if (nev != (int)mMCVtVec.size()) {
378 LOGP(debug, "source {} has {} events while {} MC vertices were booked", curSrcMC, nev, mMCVtVec.size());
379 okAccVtx = false;
380 if (nev > (int)mMCVtVec.size()) { // QED
381 continue;
382 }
383 }
384 for (curEvMC = 0; curEvMC < nev; curEvMC++) {
385 if (mVerbose > 1) {
386 LOGP(info, "Event {}", curEvMC);
387 }
388 mCurrMCTracks = &mcReader.getTracks(curSrcMC, curEvMC);
389 const_cast<o2::dataformats::MCEventHeader&>(mcReader.getMCEventHeader(curSrcMC, curEvMC)).GetVertex(mCurrMCVertex);
390 if (okAccVtx) {
391 auto& pos = mMCVtVec[curEvMC].pos;
392 if (pos[2] < -999) {
393 pos[0] = mCurrMCVertex.X();
394 pos[1] = mCurrMCVertex.Y();
395 pos[2] = mCurrMCVertex.Z();
396 }
397 }
398 for (int itr = 0; itr < mCurrMCTracks->size(); itr++) {
399 processMCParticle(curSrcMC, curEvMC, itr);
400 }
401 }
402 }
403 if (mVerbose > 0) {
404 for (int id = 0; id < mNCheckDecays; id++) {
405 LOGP(info, "Decay PDG={} : {} entries", params.decayPDG[id], mDecaysMaps[id].size());
406 }
407 }
408
409 // add reconstruction info to MC particles. If MC particle was not selected before but was reconstrected, account MC info
410 mRecProcStage = true; // MC particles accepted only at this stage will be flagged
411 for (int iv = 0; iv < nv; iv++) {
412 if (mVerbose > 1) {
413 LOGP(info, "processing PV {} of {}", iv, nv);
414 }
415 o2::MCEventLabel pvLbl;
416 int pvID = -1;
417 if (iv < (int)pvvecLbl.size()) {
418 pvLbl = pvvecLbl[iv];
419 pvID = iv;
420 if (pvLbl.isSet() && pvLbl.getEventID() < mMCVtVec.size()) {
421 mMCVtVec[pvLbl.getEventID()].recVtx.emplace_back(RecPV{pvvec[iv], pvLbl});
422 }
423 }
424 const auto& vtref = vtxRefs[iv];
425 for (int is = GTrackID::NSources; is--;) {
426 DetID::mask_t dm = GTrackID::getSourceDetectorsMask(is);
427 if (!mTracksSrc[is] || !recoData.isTrackSourceLoaded(is) || !(dm[DetID::ITS] || dm[DetID::TPC])) {
428 continue;
429 }
430 int idMin = vtref.getFirstEntryOfSource(is), idMax = idMin + vtref.getEntriesOfSource(is);
431 for (int i = idMin; i < idMax; i++) {
432 auto vid = trackIndex[i];
433 const auto& trc = recoData.getTrackParam(vid);
434 if (trc.getPt() < params.minPt || std::abs(trc.getTgl()) > params.maxTgl) {
435 continue;
436 }
437 auto lbl = recoData.getTrackMCLabel(vid);
438 if (lbl.isValid()) {
439 lbl.setFakeFlag(false);
440 auto entry = mSelMCTracks.find(lbl);
441 if (entry == mSelMCTracks.end()) { // add the track which was not added during MC scan
442 if (lbl.getSourceID() != curSrcMC || lbl.getEventID() != curEvMC) {
443 curSrcMC = lbl.getSourceID();
444 curEvMC = lbl.getEventID();
445 mCurrMCTracks = &mcReader.getTracks(curSrcMC, curEvMC);
446 const_cast<o2::dataformats::MCEventHeader&>(mcReader.getMCEventHeader(curSrcMC, curEvMC)).GetVertex(mCurrMCVertex);
447 }
448 if (!acceptMCCharged((*mCurrMCTracks)[lbl.getTrackID()], lbl)) {
449 continue;
450 }
451 entry = mSelMCTracks.find(lbl);
452 }
453 auto& trackFamily = entry->second;
454 if (vid.isAmbiguous()) { // do not repeat ambiguous tracks
455 if (trackFamily.contains(vid)) {
456 continue;
457 }
458 }
459 auto& trf = trackFamily.recTracks.emplace_back();
460 trf.gid = vid; // account(iv, vid);
461 trf.pvID = pvID;
462 trf.pvLabel = pvLbl;
463 while (dm[DetID::ITS] && dm[DetID::TPC]) { // this track should have both ITS and TPC parts, if ITS was mismatched, fill it to its proper MC track slot
464 auto gidSet = recoData.getSingleDetectorRefs(vid);
465 if (!gidSet[GTrackID::ITS].isSourceSet()) {
466 break; // AB track, nothing to check
467 }
468 auto lblITS = recoData.getTrackMCLabel(gidSet[GTrackID::ITS]);
469 if (lblITS == trackFamily.mcTrackInfo.label) {
470 break; // correct match, no need for special treatment
471 }
472 const auto& trcITSF = recoData.getTrackParam(gidSet[GTrackID::ITS]);
473 if (trcITSF.getPt() < params.minPt || std::abs(trcITSF.getTgl()) > params.maxTgl) {
474 break; // ignore this track
475 }
476 auto entryOfFake = mSelMCTracks.find(lblITS);
477 if (entryOfFake == mSelMCTracks.end()) { // this MC track was not selected
478 break;
479 }
480 auto& trackFamilyOfFake = entryOfFake->second;
481 auto& trfOfFake = trackFamilyOfFake.recTracks.emplace_back();
482 trfOfFake.gid = gidSet[GTrackID::ITS]; // account(iv, vid);
483 break;
484 }
485 if (mVerbose > 1) {
486 LOGP(info, "Matched rec track {} to MC track {}", vid.asString(), entry->first.asString());
487 }
488 } else {
489 continue;
490 }
491 }
492 }
493 }
494
495 LOGP(info, "collected {} MC tracks", mSelMCTracks.size());
496 if (params.minTPCRefsToExtractClRes > 0 || params.storeTPCTrackRefs) { // prepare MC trackrefs for TPC
497 processTPCTrackRefs();
498 }
499
500 int mcnt = 0;
501 for (auto& entry : mSelMCTracks) {
502 auto& trackFam = entry.second;
503 auto& tracks = trackFam.recTracks;
504 mcnt++;
505 if (tracks.empty()) {
506 continue;
507 }
508 if (mVerbose > 1) {
509 LOGP(info, "Processing MC track#{} {} -> {} reconstructed tracks", mcnt - 1, entry.first.asString(), tracks.size());
510 }
511 // sort according to the gid complexity (in principle, should be already sorted due to the backwards loop over NSources above
512 std::sort(tracks.begin(), tracks.end(), [](const RecTrack& lhs, const RecTrack& rhs) {
513 const auto mskL = lhs.gid.getSourceDetectorsMask();
514 const auto mskR = rhs.gid.getSourceDetectorsMask();
515 bool itstpcL = mskL[DetID::ITS] && mskL[DetID::TPC], itstpcR = mskR[DetID::ITS] && mskR[DetID::TPC];
516 if (itstpcL && !itstpcR) { // to avoid TPC/TRD or TPC/TOF shadowing ITS/TPC
517 return true;
518 }
519 return lhs.gid.getSource() > rhs.gid.getSource();
520 });
521 if (params.storeTPCTrackRefs) {
522 auto rft = mSelTRefIdx.find(entry.first);
523 if (rft != mSelTRefIdx.end()) {
524 auto rfent = rft->second;
525 for (int irf = rfent.first; irf < rfent.second; irf++) {
526 trackFam.mcTrackInfo.trackRefsTPC.push_back(mSelTRefs[irf]);
527 }
528 }
529 }
530 // fill track params
531 int tcnt = 0;
532 for (auto& tref : tracks) {
533 if (tref.gid.isSourceSet()) {
534 auto gidSet = recoData.getSingleDetectorRefs(tref.gid);
535 tref.track = recoData.getTrackParam(tref.gid);
536 if (recoData.getTrackMCLabel(tref.gid).isFake()) {
537 tref.flags |= RecTrack::FakeGLO;
538 }
539 auto msk = tref.gid.getSourceDetectorsMask();
540 if (msk[DetID::ITS]) {
541 if (gidSet[GTrackID::ITS].isSourceSet()) { // has ITS track rather than AB tracklet
542 tref.pattITS = getITSPatt(gidSet[GTrackID::ITS], tref.nClITS);
543 if (trackFam.entITS < 0) {
544 trackFam.entITS = tcnt;
545 }
546 auto lblITS = recoData.getTrackMCLabel(gidSet[GTrackID::ITS]);
547 if (lblITS.isFake()) {
548 tref.flags |= RecTrack::FakeITS;
549 }
550 if (lblITS == trackFam.mcTrackInfo.label) {
551 trackFam.entITSFound = tcnt;
552 }
553 } else { // AB ITS tracklet
554 tref.pattITS = getITSPatt(gidSet[GTrackID::ITSAB], tref.nClITS);
555 if (recoData.getTrackMCLabel(gidSet[GTrackID::ITSAB]).isFake()) {
556 tref.flags |= RecTrack::FakeITS;
557 }
558 }
559 if (msk[DetID::TPC]) {
560 if (trackFam.entITSTPC < 0) { // has both ITS and TPC contribution
561 trackFam.entITSTPC = tcnt;
562 }
563 if (recoData.getTrackMCLabel(gidSet[GTrackID::ITSTPC]).isFake()) {
564 tref.flags |= RecTrack::FakeITSTPC;
565 }
566
567 if (msk[DetID::TRD]) {
568 if (recoData.getTrackMCLabel(gidSet[GTrackID::ITSTPCTRD]).isFake()) {
569 tref.flags |= RecTrack::FakeTRD;
570 }
571 if (msk[DetID::TOF]) {
572 if (recoData.getTrackMCLabel(gidSet[GTrackID::ITSTPCTRDTOF]).isFake()) {
573 tref.flags |= RecTrack::FakeTOF;
574 }
575 }
576 } else {
577 if (msk[DetID::TOF]) {
578 if (recoData.getTrackMCLabel(gidSet[GTrackID::ITSTPCTOF]).isFake()) {
579 tref.flags |= RecTrack::FakeTOF;
580 }
581 }
582 }
583 }
584 }
585 if (msk[DetID::TPC]) {
586 const auto& trtpc = recoData.getTPCTrack(gidSet[GTrackID::TPC]);
587 tref.nClTPC = trtpc.getNClusters();
588 if (trtpc.hasBothSidesClusters()) {
589 tref.flags |= RecTrack::HASACSides;
590 }
591 fillTPCClusterInfo(trtpc, tref);
592 flagTPCClusters(trtpc, entry.first);
593 if (trackFam.entTPC < 0) {
594 trackFam.entTPC = tcnt;
595 trackFam.tpcT0 = trtpc.getTime0();
596 }
597 if (recoData.getTrackMCLabel(gidSet[GTrackID::TPC]).isFake()) {
598 tref.flags |= RecTrack::FakeTPC;
599 }
600 if (!msk[DetID::ITS]) {
601 if (msk[DetID::TRD]) {
602 if (recoData.getTrackMCLabel(gidSet[GTrackID::TPCTRD]).isFake()) {
603 tref.flags |= RecTrack::FakeTRD;
604 }
605 if (msk[DetID::TOF]) {
606 if (recoData.getTrackMCLabel(gidSet[GTrackID::TPCTRDTOF]).isFake()) {
607 tref.flags |= RecTrack::FakeTOF;
608 }
609 }
610 } else {
611 if (msk[DetID::TOF]) {
612 if (recoData.getTrackMCLabel(gidSet[GTrackID::TPCTOF]).isFake()) {
613 tref.flags |= RecTrack::FakeTOF;
614 }
615 }
616 }
617 }
618 }
619 float ts = 0, terr = 0;
620 if (tref.gid.getSource() != GTrackID::ITS) {
621 recoData.getTrackTime(tref.gid, ts, terr);
622 tref.ts = timeEst{ts, terr};
623 } else {
624 const auto& itsBra = mITSROFBracket[mITSROF[tref.gid.getIndex()]];
625 tref.ts = timeEst{itsBra.mean(), itsBra.delta() * SQRT12Inv};
626 }
627 } else {
628 LOGP(info, "Invalid entry {} of {} getTrackMCLabel {}", tcnt, tracks.size(), tref.gid.asString());
629 }
630 tcnt++;
631 }
632 if (trackFam.entITS > -1 && trackFam.entTPC > -1) { // ITS and TPC were found but matching failed
633 auto vidITS = recoData.getITSContributorGID(tracks[trackFam.entITS].gid);
634 auto vidTPC = recoData.getTPCContributorGID(tracks[trackFam.entTPC].gid);
635 auto trcTPC = recoData.getTrackParam(vidTPC);
636 auto trcITS = recoData.getTrackParamOut(vidITS);
637 if (propagateToRefX(trcTPC, trcITS)) {
638 trackFam.trackITSProp = trcITS;
639 trackFam.trackTPCProp = trcTPC;
640 } else {
641 trackFam.trackITSProp.invalidate();
642 trackFam.trackTPCProp.invalidate();
643 }
644 } else {
645 trackFam.trackITSProp.invalidate();
646 trackFam.trackTPCProp.invalidate();
647 }
648 }
649
650 // SVertices (V0s)
651 if (mCheckSV) {
652 auto v0s = recoData.getV0sIdx();
653 auto prpr = [](o2::trackstudy::TrackFamily& f) {
654 std::string s;
655 s += fmt::format(" par {} Ntpccl={} Nitscl={} ", f.mcTrackInfo.pdgParent, f.mcTrackInfo.nTPCCl, f.mcTrackInfo.nITSCl);
656 for (auto& t : f.recTracks) {
657 s += t.gid.asString();
658 s += " ";
659 }
660 return s;
661 };
662 for (int svID; svID < (int)v0s.size(); svID++) {
663 const auto& v0idx = v0s[svID];
664 int nOKProngs = 0, realMCSVID = -1;
665 int8_t decTypeID = -1;
666 for (int ipr = 0; ipr < v0idx.getNProngs(); ipr++) {
667 auto mcl = recoData.getTrackMCLabel(v0idx.getProngID(ipr)); // was this MC particle selected?
668 auto itl = mSelMCTracks.find(mcl);
669 if (itl == mSelMCTracks.end()) {
670 nOKProngs = -1; // was not selected as interesting one, ignore
671 break;
672 }
673 auto& trackFamily = itl->second;
674 int decayParentIndex = trackFamily.mcTrackInfo.parentEntry;
675 if (decayParentIndex < 0) { // does not come from decay
676 break;
677 }
678 if (ipr == 0) {
679 realMCSVID = decayParentIndex;
680 decTypeID = trackFamily.mcTrackInfo.parentDecID;
681 nOKProngs = 1;
682 LOGP(debug, "Prong{} {} comes from {}/{}", ipr, prpr(trackFamily), decTypeID, realMCSVID);
683 continue;
684 }
685 if (realMCSVID != decayParentIndex || decTypeID != trackFamily.mcTrackInfo.parentDecID) {
686 break;
687 }
688 LOGP(debug, "Prong{} {} comes from {}/{}", ipr, prpr(trackFamily), decTypeID, realMCSVID);
689 nOKProngs++;
690 }
691 if (nOKProngs == v0idx.getNProngs()) { // all prongs are from the decay of MC parent which deemed to be interesting, flag it
692 LOGP(debug, "Decay {}/{} was found", decTypeID, realMCSVID);
693 mDecaysMaps[decTypeID][realMCSVID].foundSVID = svID;
694 }
695 }
696 }
697
698 // collect ITS/TPC cluster info for selected MC particles
699 fillMCClusterInfo(recoData);
700
701 // single tracks
702 for (auto& entry : mSelMCTracks) {
703 auto& trackFam = entry.second;
704 (*mDBGOut) << "tracks" << "tr=" << trackFam << "\n";
705 }
706
707 // decays
708 std::vector<TrackFamily> decFam;
709 for (int id = 0; id < mNCheckDecays; id++) {
710 std::string decTreeName = fmt::format("dec{}", params.decayPDG[id]);
711 for (const auto& dec : mDecaysMaps[id]) {
712 decFam.clear();
713 bool skip = false;
714 for (int idd = dec.daughterFirst; idd <= dec.daughterLast; idd++) {
715 auto dtLbl = mDecProdLblPool[idd]; // daughter MC label
716 const auto& dtFamily = mSelMCTracks[dtLbl];
717 if (dtFamily.mcTrackInfo.pdgParent != dec.pdg) {
718 LOGP(error, "{}-th decay (pdg={}): {} in {}:{} range refers to MC track with pdgParent = {}", id, params.decayPDG[id], idd, dec.daughterFirst, dec.daughterLast, dtFamily.mcTrackInfo.pdgParent);
719 skip = true;
720 break;
721 }
722 decFam.push_back(dtFamily);
723 }
724 if (!skip) {
726 if (dec.foundSVID >= 0 && !refitV0(dec.foundSVID, v0, recoData)) {
727 v0.invalidate();
728 }
729 (*mDBGOut) << decTreeName.c_str() << "pdgPar=" << dec.pdg << "trPar=" << dec.parent << "prod=" << decFam << "found=" << dec.foundSVID << "sv=" << v0 << "\n";
730 }
731 }
732 }
733
734 for (auto& mcVtx : mMCVtVec) { // sort rec.vertices in mult. order
735 std::sort(mcVtx.recVtx.begin(), mcVtx.recVtx.end(), [](const RecPV& lhs, const RecPV& rhs) {
736 return lhs.pv.getNContributors() > rhs.pv.getNContributors();
737 });
738 (*mDBGOut) << "mcVtxTree" << "mcVtx=" << mcVtx << "\n";
739 }
740
741 if (params.storeITSInfo) {
742 processITSTracks(recoData);
743 }
744}
745
746void TrackMCStudy::processTPCTrackRefs()
747{
748 constexpr float alpsec[18] = {0.174533, 0.523599, 0.872665, 1.221730, 1.570796, 1.919862, 2.268928, 2.617994, 2.967060, 3.316126, 3.665191, 4.014257, 4.363323, 4.712389, 5.061455, 5.410521, 5.759587, 6.108652};
749 constexpr float sinAlp[18] = {0.173648, 0.500000, 0.766044, 0.939693, 1.000000, 0.939693, 0.766044, 0.500000, 0.173648, -0.173648, -0.500000, -0.766044, -0.939693, -1.000000, -0.939693, -0.766044, -0.500000, -0.173648};
750 constexpr float cosAlp[18] = {0.984808, 0.866025, 0.642788, 0.342020, 0.000000, -0.342020, -0.642788, -0.866025, -0.984808, -0.984808, -0.866025, -0.642788, -0.342020, -0.000000, 0.342020, 0.642788, 0.866025, 0.984808};
752 for (auto& entry : mSelMCTracks) {
753 auto lb = entry.first;
754 auto trspan = mcReader.getTrackRefs(lb.getSourceID(), lb.getEventID(), lb.getTrackID());
755 int q = entry.second.mcTrackInfo.track.getCharge();
756 if (q * q != 1) {
757 continue;
758 }
759 int ref0entry = mSelTRefs.size(), nrefsSel = 0;
760 for (const auto& trf : trspan) {
761 if (trf.getDetectorId() != 1) { // process TPC only
762 continue;
763 }
764 float pT = std::sqrt(trf.Px() * trf.Px() + trf.Py() * trf.Py());
765 if (pT < 0.05) {
766 continue;
767 }
768 float secX, secY, phi = std::atan2(trf.Y(), trf.X());
769 int sector = o2::math_utils::angle2Sector(phi);
770 o2::math_utils::rotateZInv(trf.X(), trf.Y(), secX, secY, sinAlp[sector], cosAlp[sector]); // sector coordinates
771 float phiPt = std::atan2(trf.Py(), trf.Px());
772 o2::math_utils::bringTo02Pi(phiPt);
773 auto dphiPt = phiPt - alpsec[sector];
774 if (dphiPt > o2::constants::math::PI) { // account for wraps
776 } else if (dphiPt < -o2::constants::math::PI) {
778 } else if (std::abs(dphiPt) > o2::constants::math::PIHalf * 0.8) {
779 continue; // ignore backward going or parallel to padrows tracks
780 }
781 float tgL = trf.Pz() / pT;
782 std::array<float, 5> pars = {secY, trf.Z(), std::sin(dphiPt), tgL, q / pT};
783 auto& refTrack = mSelTRefs.emplace_back(secX, alpsec[sector], pars);
784 refTrack.setUserField(uint16_t(sector));
785 nrefsSel++;
786 }
787 if (nrefsSel < params.minTPCRefsToExtractClRes) {
788 mSelTRefs.resize(ref0entry); // discard unused tracks
789 continue;
790 } else {
791 mSelTRefIdx[lb] = std::make_pair(ref0entry, ref0entry + nrefsSel);
792 }
793 }
794}
795
796void TrackMCStudy::fillMCClusterInfo(const o2::globaltracking::RecoContainer& recoData)
797{
798 // TPC clusters info
799 const auto& TPCClusterIdxStruct = recoData.inputsTPCclusters->clusterIndex;
800 const auto* TPCClMClab = recoData.inputsTPCclusters->clusterIndex.clustersMCTruth;
802
803 ClResTPC clRes{};
804 for (uint8_t row = 0; row < 152; row++) { // we need to go in increasing row, so this should be the outer loop
805 for (uint8_t sector = 0; sector < 36; sector++) {
806 unsigned int offs = TPCClusterIdxStruct.clusterOffset[sector][row];
807 for (unsigned int icl0 = 0; icl0 < TPCClusterIdxStruct.nClusters[sector][row]; icl0++) {
808 const auto labels = TPCClMClab->getLabels(icl0 + offs);
809 int ncontLb = 0; // number of real contrubutors to this label (w/o noise)
810 for (const auto& lbl : labels) {
811 if (!lbl.isValid()) {
812 continue;
813 }
814 ncontLb++;
815 }
816 const auto& clus = TPCClusterIdxStruct.clusters[sector][row][icl0];
817 int tbinH = int(clus.getTime() * mNTPCOccBinLengthInv); // time bin converted to slot of the occ. histo
818 clRes.contTracks.clear();
819 bool doClusRes = (params.minTPCRefsToExtractClRes > 0) && (params.rejectClustersResStat <= 0. || gRandom->Rndm() < params.rejectClustersResStat);
820 for (auto lbl : labels) {
821 bool corrAttach = lbl.isFake(); // was this flagged in the flagTPCClusters called from process ?
822 lbl.setFakeFlag(false);
823 auto entry = mSelMCTracks.find(lbl);
824 if (entry == mSelMCTracks.end()) { // not selected
825 continue;
826 }
827 auto& mctr = entry->second.mcTrackInfo;
828 mctr.nTPCCl++;
829 if (row > mctr.maxTPCRow) {
830 mctr.maxTPCRow = row;
831 mctr.maxTPCRowSect = sector;
832 mctr.nUsedPadRows++;
833 } else if (row == 0 && mctr.nUsedPadRows == 0) {
834 mctr.nUsedPadRows++;
835 }
836 if (row < mctr.minTPCRow) {
837 mctr.minTPCRow = row;
838 mctr.minTPCRowSect = sector;
839 }
840 if (mctr.minTPCRowSect == sector && row > mctr.maxTPCRowInner) {
841 mctr.maxTPCRowInner = row;
842 }
843 if (ncontLb > 1) {
844 mctr.nTPCClShared++;
845 }
846 // try to extract ideal track position
847 if (doClusRes) {
848 auto entTRefIDsIt = mSelTRefIdx.find(lbl);
849 if (entTRefIDsIt == mSelTRefIdx.end()) {
850 continue;
851 }
852 float xc, yc, zc;
853 mTPCCorrMaps->Transform(sector, row, clus.getPad(), clus.getTime(), xc, yc, zc, mctr.bcInTF / 8.); // nominal time of the track
854
855 const auto& entTRefIDs = entTRefIDsIt->second;
856 // find bracketing TRef params
857 int entIDBelow = -1, entIDAbove = -1;
858 float xBelow = -1e6, xAbove = 1e6;
859
860 for (int entID = entTRefIDs.first; entID < entTRefIDs.second; entID++) {
861 const auto& refTr = mSelTRefs[entID];
862 if (refTr.getUserField() != sector % 18) {
863 continue;
864 }
865 if ((refTr.getX() < xc) && (refTr.getX() > xBelow) && (refTr.getX() > xc - params.maxTPCRefExtrap)) {
866 xBelow = refTr.getX();
867 entIDBelow = entID;
868 }
869 if ((refTr.getX() > xc) && (refTr.getX() < xAbove) && (refTr.getX() < xc + params.maxTPCRefExtrap)) {
870 xAbove = refTr.getX();
871 entIDAbove = entID;
872 }
873 }
874 if ((entIDBelow < 0 && entIDAbove < 0) || (params.requireTopBottomRefs && (entIDBelow < 0 || entIDAbove < 0))) {
875 continue;
876 }
877 auto prop = o2::base::Propagator::Instance();
878 o2::track::TrackPar tparAbove, tparBelow;
879 bool okBelow = entIDBelow >= 0 && prop->PropagateToXBxByBz((tparBelow = mSelTRefs[entIDBelow]), xc, 0.99, 2.);
880 bool okAbove = entIDAbove >= 0 && prop->PropagateToXBxByBz((tparAbove = mSelTRefs[entIDAbove]), xc, 0.99, 2.);
881 if ((!okBelow && !okAbove) || (params.requireTopBottomRefs && (!okBelow || !okAbove))) {
882 continue;
883 }
884
885 int nmeas = 0;
886 auto& clCont = clRes.contTracks.emplace_back();
887 clCont.corrAttach = corrAttach;
888 if (okBelow) {
889 clCont.below = {mSelTRefs[entIDBelow].getX(), tparBelow.getY(), tparBelow.getZ()};
890 clCont.snp += tparBelow.getSnp();
891 clCont.tgl += tparBelow.getTgl();
892 clCont.q2pt += tparBelow.getQ2Pt();
893 nmeas++;
894 }
895 if (okAbove) {
896 clCont.above = {mSelTRefs[entIDAbove].getX(), tparAbove.getY(), tparAbove.getZ()};
897 clCont.snp += tparAbove.getSnp();
898 clCont.tgl += tparAbove.getTgl();
899 clCont.q2pt += tparAbove.getQ2Pt();
900 nmeas++;
901 }
902 if (nmeas) {
903 if (clRes.contTracks.size() == 1) {
904 int occBin = mctr.bcInTF / 8 * mNTPCOccBinLengthInv;
905 clRes.occ = occBin < 0 ? mTBinClOcc[0] : (occBin >= mTBinClOcc.size() ? mTBinClOcc.back() : mTBinClOcc[occBin]);
906 }
907 clCont.xyz = {xc, yc, zc};
908 if (nmeas > 1) {
909 clCont.snp *= 0.5;
910 clCont.tgl *= 0.5;
911 clCont.q2pt *= 0.5;
912 }
913 } else {
914 clRes.contTracks.pop_back();
915 }
916 }
917 }
918 if (clRes.getNCont()) {
919 clRes.sect = sector;
920 clRes.row = row;
921 clRes.qtot = clus.getQtot();
922 clRes.qmax = clus.getQmax();
923 clRes.flags = clus.getFlags();
924 clRes.sigmaTimePacked = clus.sigmaTimePacked;
925 clRes.sigmaPadPacked = clus.sigmaPadPacked;
926 clRes.ncont = ncontLb;
927 clRes.sortCont();
928
929 if (tbinH < 0) {
930 tbinH = 0;
931 } else if (tbinH >= int(mTBinClOccHist.size())) {
932 tbinH = (int)mTBinClOccHist.size() - 1;
933 }
934 clRes.occBin = mTBinClOccHist[tbinH];
935
936 (*mDBGOut) << "clres" << "clr=" << clRes << "\n";
937 }
938 }
939 }
940 }
941 // fill ITS cluster info
942 const auto* mcITSClusters = recoData.getITSClustersMCLabels();
943 const auto& ITSClusters = recoData.getITSClusters();
944 for (unsigned int icl = 0; icl < ITSClusters.size(); icl++) {
945 const auto labels = mcITSClusters->getLabels(icl);
946 for (const auto& lbl : labels) {
947 auto entry = mSelMCTracks.find(lbl);
948 if (entry == mSelMCTracks.end()) { // not selected
949 continue;
950 }
951 auto& mctr = entry->second.mcTrackInfo;
952 mctr.nITSCl++;
953 mctr.pattITSCl |= 0x1 << o2::itsmft::ChipMappingITS::getLayer(ITSClusters[icl].getChipID());
954 }
955 }
956}
957
958bool TrackMCStudy::propagateToRefX(o2::track::TrackParCov& trcTPC, o2::track::TrackParCov& trcITS)
959{
960 bool refReached = false;
961 constexpr float TgHalfSector = 0.17632698f;
963 int trialsLeft = 2;
964 while (o2::base::Propagator::Instance()->PropagateToXBxByBz(trcTPC, par.XMatchingRef, MaxSnp, 2., par.matCorr)) {
965 if (refReached) {
966 break;
967 }
968 // make sure the track is indeed within the sector defined by alpha
969 if (fabs(trcTPC.getY()) < par.XMatchingRef * TgHalfSector) {
970 refReached = true;
971 break; // ok, within
972 }
973 if (!trialsLeft--) {
974 break;
975 }
976 auto alphaNew = o2::math_utils::angle2Alpha(trcTPC.getPhiPos());
977 if (!trcTPC.rotate(alphaNew) != 0) {
978 break; // failed (RS: check effect on matching tracks to neighbouring sector)
979 }
980 }
981 if (!refReached) {
982 return false;
983 }
984 refReached = false;
985 float alp = trcTPC.getAlpha();
986 if (!trcITS.rotate(alp) != 0 || !o2::base::Propagator::Instance()->PropagateToXBxByBz(trcITS, par.XMatchingRef, MaxSnp, 2., par.matCorr)) {
987 return false;
988 }
989 return true;
990}
991
992void TrackMCStudy::endOfStream(EndOfStreamContext& ec)
993{
994 mDBGOut.reset();
995}
996
997void TrackMCStudy::finaliseCCDB(ConcreteDataMatcher& matcher, void* obj)
998{
999 if (o2::base::GRPGeomHelper::instance().finaliseCCDB(matcher, obj)) {
1000 return;
1001 }
1002 if (mTPCVDriftHelper.accountCCDBInputs(matcher, obj)) {
1003 return;
1004 }
1005 if (matcher == ConcreteDataMatcher("ITS", "ALPIDEPARAM", 0)) {
1006 LOG(info) << "ITS Alpide param updated";
1008 par.printKeyValues();
1009 mITSTimeBiasMUS = par.roFrameBiasInBC * o2::constants::lhc::LHCBunchSpacingNS * 1e-3;
1010 mITSROFrameLengthMUS = par.roFrameLengthInBC * o2::constants::lhc::LHCBunchSpacingNS * 1e-3;
1011 return;
1012 }
1013 if (matcher == ConcreteDataMatcher("ITS", "CLUSDICT", 0)) {
1014 LOG(info) << "cluster dictionary updated";
1015 mITSDict = (const o2::itsmft::TopologyDictionary*)obj;
1016 return;
1017 }
1018}
1019
1020//_____________________________________________________
1021void TrackMCStudy::prepareITSData(const o2::globaltracking::RecoContainer& recoData)
1022{
1023 const auto ITSTracksArray = recoData.getITSTracks();
1024 const auto ITSTrackROFRec = recoData.getITSTracksROFRecords();
1025 int nROFs = ITSTrackROFRec.size();
1026 mITSROF.clear();
1027 mITSROFBracket.clear();
1028 mITSROF.reserve(ITSTracksArray.size());
1029 mITSROFBracket.reserve(ITSTracksArray.size());
1030 for (int irof = 0; irof < nROFs; irof++) {
1031 const auto& rofRec = ITSTrackROFRec[irof];
1032 long nBC = rofRec.getBCData().differenceInBC(recoData.startIR);
1033 float tMin = nBC * o2::constants::lhc::LHCBunchSpacingMUS + mITSTimeBiasMUS;
1034 float tMax = tMin + mITSROFrameLengthMUS;
1035 mITSROFBracket.emplace_back(tMin, tMax);
1036 for (int it = 0; it < rofRec.getNEntries(); it++) {
1037 mITSROF.push_back(irof);
1038 }
1039 }
1040}
1041/*
1042float TrackMCStudy::getDCAYCut(float pt) const
1043{
1044 static TF1 fun("dcayvspt", mDCAYFormula.c_str(), 0, 20);
1045 return fun.Eval(pt);
1046}
1047*/
1048
1049bool TrackMCStudy::processMCParticle(int src, int ev, int trid)
1050{
1051 const auto& mcPart = (*mCurrMCTracks)[trid];
1052 int pdg = mcPart.GetPdgCode();
1053 bool res = false;
1054 while (true) {
1055 auto lbl = o2::MCCompLabel(trid, ev, src);
1056 int decay = -1; // is this decay to watch?
1058 if (mcPart.T() < params.decayMotherMaxT) {
1059 for (int id = 0; id < mNCheckDecays; id++) {
1060 if (params.decayPDG[id] == std::abs(pdg)) {
1061 decay = id;
1062 break;
1063 }
1064 }
1065 if (decay >= 0) { // check if decay and kinematics is acceptable
1066 auto& decayPool = mDecaysMaps[decay];
1067 int idd0 = mcPart.getFirstDaughterTrackId(), idd1 = mcPart.getLastDaughterTrackId(); // we want only charged and trackable daughters
1068 int dtStart = mDecProdLblPool.size(), dtEnd = -1;
1069 if (idd0 < 0) {
1070 break;
1071 }
1072 for (int idd = idd0; idd <= idd1; idd++) {
1073 const auto& product = (*mCurrMCTracks)[idd];
1074 auto lbld = o2::MCCompLabel(idd, ev, src);
1075 if (!acceptMCCharged(product, lbld, decay)) {
1076 decay = -1; // discard decay
1077 mDecProdLblPool.resize(dtStart);
1078 break;
1079 }
1080 mDecProdLblPool.push_back(lbld); // register prong entry and label
1081 }
1082 if (decay >= 0) {
1083 // account decay
1084 dtEnd = mDecProdLblPool.size();
1085 for (int dtid = dtStart; dtid < dtEnd; dtid++) { // flag selected decay parent entry in the prongs MCs
1086 mSelMCTracks[mDecProdLblPool[dtid]].mcTrackInfo.parentEntry = decayPool.size();
1087 mSelMCTracks[mDecProdLblPool[dtid]].mcTrackInfo.parentDecID = int8_t(decay);
1088 }
1089 dtEnd--;
1090 std::array<float, 3> xyz{(float)mcPart.GetStartVertexCoordinatesX(), (float)mcPart.GetStartVertexCoordinatesY(), (float)mcPart.GetStartVertexCoordinatesZ()};
1091 std::array<float, 3> pxyz{(float)mcPart.GetStartVertexMomentumX(), (float)mcPart.GetStartVertexMomentumY(), (float)mcPart.GetStartVertexMomentumZ()};
1092 decayPool.emplace_back(DecayRef{lbl,
1093 o2::track::TrackPar(xyz, pxyz, TMath::Nint(O2DatabasePDG::Instance()->GetParticle(mcPart.GetPdgCode())->Charge() / 3), false),
1094 mcPart.GetPdgCode(), dtStart, dtEnd});
1095 if (mVerbose > 1) {
1096 LOGP(info, "Adding MC parent pdg={} {}, with prongs in {}:{} range", pdg, lbl.asString(), dtStart, dtEnd);
1097 }
1098 res = true; // Accept!
1099 }
1100 break;
1101 }
1102 }
1103 // check if this is a charged which should be processed but was not accounted as a decay product
1104 if (mSelMCTracks.find(lbl) == mSelMCTracks.end()) {
1105 res = acceptMCCharged(mcPart, lbl);
1106 }
1107 break;
1108 }
1109 return res;
1110}
1111
1112bool TrackMCStudy::acceptMCCharged(const MCTrack& tr, const o2::MCCompLabel& lb, int followDecay)
1113{
1115 if (tr.GetPt() < params.minPtMC ||
1116 std::abs(tr.GetTgl()) > params.maxTglMC ||
1117 tr.R2() > params.maxRMC * params.maxRMC) {
1118 if (mVerbose > 1 && followDecay > -1) {
1119 LOGP(info, "rejecting decay {} prong : pdg={}, pT={}, tgL={}, r={}", followDecay, tr.GetPdgCode(), tr.GetPt(), tr.GetTgl(), std::sqrt(tr.R2()));
1120 }
1121 return false;
1122 }
1123 float dx = tr.GetStartVertexCoordinatesX() - mCurrMCVertex.X(), dy = tr.GetStartVertexCoordinatesY() - mCurrMCVertex.Y(), dz = tr.GetStartVertexCoordinatesZ() - mCurrMCVertex.Z();
1124 float r2 = dx * dx + dy * dy;
1125 float posTgl2 = r2 > 1 && std::abs(dz) < 20 ? dz * dz / r2 : 0;
1126 if (posTgl2 > params.maxPosTglMC * params.maxPosTglMC) {
1127 if (mVerbose > 1 && followDecay > -1) {
1128 LOGP(info, "rejecting decay {} prong : pdg={}, pT={}, tgL={}, dr={}, dz={} r={}, z={}, posTgl={}", followDecay, tr.GetPdgCode(), tr.GetPt(), tr.GetTgl(), std::sqrt(r2), dz, std::sqrt(tr.R2()), tr.GetStartVertexCoordinatesZ(), std::sqrt(posTgl2));
1129 }
1130 return false;
1131 }
1132 if (params.requireITSorTPCTrackRefs) {
1133 auto trspan = mcReader.getTrackRefs(lb.getSourceID(), lb.getEventID(), lb.getTrackID());
1134 bool ok = false;
1135 for (const auto& trf : trspan) {
1136 if (trf.getDetectorId() == DetID::ITS || trf.getDetectorId() == DetID::TPC) {
1137 ok = true;
1138 break;
1139 }
1140 }
1141 if (!ok) {
1142 return false;
1143 }
1144 }
1145 TParticlePDG* pPDG = O2DatabasePDG::Instance()->GetParticle(tr.GetPdgCode());
1146 if (!pPDG) {
1147 LOGP(debug, "Unknown particle {}", tr.GetPdgCode());
1148 return false;
1149 }
1150 if (pPDG->Charge() == 0.) {
1151 return false;
1152 }
1153 return addMCParticle(tr, lb, pPDG);
1154}
1155
1156bool TrackMCStudy::addMCParticle(const MCTrack& mcPart, const o2::MCCompLabel& lb, TParticlePDG* pPDG)
1157{
1158 std::array<float, 3> xyz{(float)mcPart.GetStartVertexCoordinatesX(), (float)mcPart.GetStartVertexCoordinatesY(), (float)mcPart.GetStartVertexCoordinatesZ()};
1159 std::array<float, 3> pxyz{(float)mcPart.GetStartVertexMomentumX(), (float)mcPart.GetStartVertexMomentumY(), (float)mcPart.GetStartVertexMomentumZ()};
1160 if (!pPDG && !(pPDG = O2DatabasePDG::Instance()->GetParticle(mcPart.GetPdgCode()))) {
1161 LOGP(debug, "Unknown particle {}", mcPart.GetPdgCode());
1162 return false;
1163 }
1164 auto& mcEntry = mSelMCTracks[lb];
1165 mcEntry.mcTrackInfo.pdg = mcPart.GetPdgCode();
1166 mcEntry.mcTrackInfo.track = o2::track::TrackPar(xyz, pxyz, TMath::Nint(pPDG->Charge() / 3), true);
1167 mcEntry.mcTrackInfo.label = lb;
1168 mcEntry.mcTrackInfo.bcInTF = mIntBC[lb.getEventID()];
1169 mcEntry.mcTrackInfo.occTPC = mTPCOcc[lb.getEventID()];
1170 mcEntry.mcTrackInfo.occITS = mITSOcc[lb.getEventID()];
1171 mcEntry.mcTrackInfo.occTPCV = mMCVtVec[lb.getEventID()].occTPCV;
1172 if (mRecProcStage) {
1173 mcEntry.mcTrackInfo.setAddedAtRecStage();
1174 }
1175 if (o2::mcutils::MCTrackNavigator::isPhysicalPrimary(mcPart, *mCurrMCTracks)) {
1176 mcEntry.mcTrackInfo.setPrimary();
1177 }
1178 int moth = -1;
1179 o2::MCCompLabel mclbPar;
1180 if ((moth = mcPart.getMotherTrackId()) >= 0) {
1181 const auto& mcPartPar = (*mCurrMCTracks)[moth];
1182 mcEntry.mcTrackInfo.pdgParent = mcPartPar.GetPdgCode();
1183 }
1184 if (mcPart.isPrimary() && mcReader.getNEvents(lb.getSourceID()) == mMCVtVec.size()) {
1185 mMCVtVec[lb.getEventID()].nTrackSel++;
1186 }
1187 if (mVerbose > 1) {
1188 LOGP(info, "Adding charged MC pdg={} {} ", mcPart.GetPdgCode(), lb.asString());
1189 }
1190 return true;
1191}
1192
1193bool TrackMCStudy::refitV0(int i, o2::dataformats::V0& v0, const o2::globaltracking::RecoContainer& recoData)
1194{
1195 const auto& id = recoData.getV0sIdx()[i];
1196 auto seedP = recoData.getTrackParam(id.getProngID(0));
1197 auto seedN = recoData.getTrackParam(id.getProngID(1));
1198 bool isTPConly = (id.getProngID(0).getSource() == GTrackID::TPC) || (id.getProngID(1).getSource() == GTrackID::TPC);
1199 const auto& svparam = o2::vertexing::SVertexerParams::Instance();
1200 if (svparam.mTPCTrackPhotonTune && isTPConly) {
1201 mFitterV0.setMaxDZIni(svparam.mTPCTrackMaxDZIni);
1202 mFitterV0.setMaxDXYIni(svparam.mTPCTrackMaxDXYIni);
1203 mFitterV0.setMaxChi2(svparam.mTPCTrackMaxChi2);
1204 mFitterV0.setCollinear(true);
1205 }
1206 int nCand = mFitterV0.process(seedP, seedN);
1207 if (svparam.mTPCTrackPhotonTune && isTPConly) { // restore
1208 // Reset immediately to the defaults
1209 mFitterV0.setMaxDZIni(svparam.maxDZIni);
1210 mFitterV0.setMaxDXYIni(svparam.maxDXYIni);
1211 mFitterV0.setMaxChi2(svparam.maxChi2);
1212 mFitterV0.setCollinear(false);
1213 }
1214 if (nCand == 0) { // discard this pair
1215 return false;
1216 }
1217 const int cand = 0;
1218 if (!mFitterV0.isPropagateTracksToVertexDone(cand) && !mFitterV0.propagateTracksToVertex(cand)) {
1219 return false;
1220 }
1221 const auto& trPProp = mFitterV0.getTrack(0, cand);
1222 const auto& trNProp = mFitterV0.getTrack(1, cand);
1223 std::array<float, 3> pP{}, pN{};
1224 trPProp.getPxPyPzGlo(pP);
1225 trNProp.getPxPyPzGlo(pN);
1226 std::array<float, 3> pV0 = {pP[0] + pN[0], pP[1] + pN[1], pP[2] + pN[2]};
1227 auto p2V0 = pV0[0] * pV0[0] + pV0[1] * pV0[1] + pV0[2] * pV0[2];
1228 const auto& pv = recoData.getPrimaryVertex(id.getVertexID());
1229 const auto v0XYZ = mFitterV0.getPCACandidatePos(cand);
1230 float dx = v0XYZ[0] - pv.getX(), dy = v0XYZ[1] - pv.getY(), dz = v0XYZ[2] - pv.getZ(), prodXYZv0 = dx * pV0[0] + dy * pV0[1] + dz * pV0[2];
1231 float cosPA = prodXYZv0 / std::sqrt((dx * dx + dy * dy + dz * dz) * p2V0);
1232 new (&v0) o2::dataformats::V0(v0XYZ, pV0, mFitterV0.calcPCACovMatrixFlat(cand), trPProp, trNProp);
1233 v0.setDCA(mFitterV0.getChi2AtPCACandidate(cand));
1234 v0.setCosPA(cosPA);
1235 return true;
1236}
1237
1238void TrackMCStudy::loadTPCOccMap(const o2::globaltracking::RecoContainer& recoData)
1239{
1240 auto NHBPerTF = o2::base::GRPGeomHelper::instance().getGRPECS()->getNHBFPerTF();
1241 const auto& TPCOccMap = recoData.occupancyMapTPC;
1242 auto prop = o2::base::Propagator::Instance();
1243 auto TPCRefitter = std::make_unique<o2::gpu::GPUO2InterfaceRefit>(&recoData.inputsTPCclusters->clusterIndex, mTPCCorrMaps, prop->getNominalBz(),
1244 recoData.getTPCTracksClusterRefs().data(), 0, recoData.clusterShMapTPC.data(), TPCOccMap.data(), TPCOccMap.size(), nullptr, prop);
1245 mNTPCOccBinLength = TPCRefitter->getParam()->rec.tpc.occupancyMapTimeBins;
1246 mTBinClOcc.clear();
1247 if (mNTPCOccBinLength > 1 && TPCOccMap.size()) {
1248 mNTPCOccBinLengthInv = 1. / mNTPCOccBinLength;
1249 int nTPCBins = NHBPerTF * o2::constants::lhc::LHCMaxBunches / 8, ninteg = 0;
1250 int nTPCOccBins = nTPCBins * mNTPCOccBinLengthInv, sumBins = std::max(1, int(o2::constants::lhc::LHCMaxBunches / 8 * mNTPCOccBinLengthInv));
1251 mTBinClOcc.resize(nTPCOccBins);
1252 mTBinClOccHist.resize(nTPCOccBins);
1253 float sm = 0., tb = 0.5 * mNTPCOccBinLength;
1254 for (int i = 0; i < nTPCOccBins; i++) {
1255 mTBinClOccHist[i] = TPCRefitter->getParam()->GetUnscaledMult(tb);
1256 tb += mNTPCOccBinLength;
1257 }
1258 for (int i = nTPCOccBins; i--;) {
1259 sm += mTBinClOccHist[i];
1260 if (i + sumBins < nTPCOccBins) {
1261 sm -= mTBinClOccHist[i + sumBins];
1262 }
1263 mTBinClOcc[i] = sm;
1264 }
1265 } else {
1266 mTBinClOcc.resize(1);
1267 mTBinClOccHist.resize(1);
1268 }
1269}
1270
1271void TrackMCStudy::processITSTracks(const o2::globaltracking::RecoContainer& recoData)
1272{
1273 if (!mITSDict) {
1274 LOGP(warn, "ITS data is not loaded");
1275 return;
1276 }
1277 const auto itsTracks = recoData.getITSTracks();
1278 const auto itsLbls = recoData.getITSTracksMCLabels();
1279 const auto itsClRefs = recoData.getITSTracksClusterRefs();
1280 const auto clusITS = recoData.getITSClusters();
1281 const auto patterns = recoData.getITSClustersPatterns();
1283 auto pattIt = patterns.begin();
1284 mITSClustersArray.clear();
1285 mITSClustersArray.reserve(clusITS.size());
1286
1287 o2::its::ioutils::convertCompactClusters(clusITS, pattIt, mITSClustersArray, mITSDict);
1288 auto geom = o2::its::GeometryTGeo::Instance();
1289 int ntr = itsLbls.size();
1290 LOGP(info, "We have {} ITS clusters and the number of patterns is {}, ITSdict:{} NMCLabels: {}", clusITS.size(), patterns.size(), mITSDict != nullptr, itsLbls.size());
1291
1292 std::vector<int> evord(ntr);
1293 std::iota(evord.begin(), evord.end(), 0);
1294 std::sort(evord.begin(), evord.end(), [&](int i, int j) { return itsLbls[i] < itsLbls[j]; });
1295 std::vector<ITSHitInfo> outHitInfo;
1296 std::array<int, 7> cl2arr{};
1297
1298 for (int itr0 = 0; itr0 < ntr; itr0++) {
1299 auto itr = evord[itr0];
1300 const auto& itsTr = itsTracks[itr];
1301 const auto& itsLb = itsLbls[itr];
1302 // LOGP(info,"proc {} {} {}",itr0, itr, itsLb.asString());
1303 int nCl = itsTr.getNClusters();
1304 if (itsLb.isFake() || nCl < params.minITSClForITSoutput) {
1305 continue;
1306 }
1307 auto entrySel = mSelMCTracks.find(itsLb);
1308 if (entrySel == mSelMCTracks.end()) {
1309 continue;
1310 }
1311 outHitInfo.clear();
1312 cl2arr.fill(-1);
1313 auto clEntry = itsTr.getFirstClusterEntry();
1314 for (int iCl = nCl; iCl--;) { // clusters are stored from outer to inner layers
1315 const auto& cls = mITSClustersArray[itsClRefs[clEntry + iCl]];
1316 int hpos = outHitInfo.size();
1317 auto& hinf = outHitInfo.emplace_back();
1318 hinf.clus = cls;
1319 hinf.clus.setCount(geom->getLayer(cls.getSensorID()));
1320 geom->getSensorXAlphaRefPlane(cls.getSensorID(), hinf.chipX, hinf.chipAlpha);
1321 cl2arr[hinf.clus.getCount()] = hpos; // to facilitate finding the cluster of the layer
1322 }
1323 auto trspan = mcReader.getTrackRefs(itsLb.getSourceID(), itsLb.getEventID(), itsLb.getTrackID());
1324 int ilrc = -1, nrefAcc = 0;
1325 for (const auto& trf : trspan) {
1326 if (trf.getDetectorId() != 0) { // process ITS only
1327 continue;
1328 }
1329 int lrt = trf.getUserId(); // layer of the reference, but there might be multiple hits on the same layer
1330 int clEnt = cl2arr[lrt];
1331 if (clEnt < 0) {
1332 continue;
1333 }
1334 auto& hinf = outHitInfo[clEnt];
1335 float traX, traY;
1336 o2::math_utils::rotateZInv(trf.X(), trf.Y(), traX, traY, std::sin(hinf.chipAlpha), std::cos(hinf.chipAlpha)); // tracking coordinates of the reference
1337 if (hinf.trefXT < 1 || std::abs(traX - hinf.chipX) < std::abs(hinf.trefXT - hinf.chipX)) {
1338 if (hinf.trefXT < 1) {
1339 nrefAcc++;
1340 }
1341 hinf.tref = trf;
1342 hinf.trefXT = traX;
1343 hinf.trefYT = traY;
1344 }
1345 }
1346 (*mDBGOut) << "itsTree" << "hits=" << outHitInfo << "trIn=" << ((o2::track::TrackParCov&)itsTr) << "trOut=" << itsTr.getParamOut() << "mcTr=" << entrySel->second.mcTrackInfo.track << "mcPDG=" << entrySel->second.mcTrackInfo.pdg << "nTrefs=" << nrefAcc << "\n";
1347 }
1348}
1349
1351{
1352 std::vector<OutputSpec> outputs;
1353 Options opts{
1354 {"device-verbosity", VariantType::Int, 0, {"Verbosity level"}},
1355 {"dcay-vs-pt", VariantType::String, "0.0105 + 0.0350 / pow(x, 1.1)", {"Formula for global tracks DCAy vs pT cut"}},
1356 {"min-tpc-clusters", VariantType::Int, 60, {"Cut on TPC clusters"}},
1357 {"max-tpc-dcay", VariantType::Float, 2.f, {"Cut on TPC dcaY"}},
1358 {"max-tpc-dcaz", VariantType::Float, 2.f, {"Cut on TPC dcaZ"}},
1359 {"min-x-prop", VariantType::Float, 6.f, {"track should be propagated to this X at least"}}};
1360 auto dataRequest = std::make_shared<DataRequest>();
1361 bool useMC = true;
1362 dataRequest->requestTracks(srcTracks, useMC);
1363 dataRequest->requestClusters(srcClusters, useMC);
1364 dataRequest->requestPrimaryVertices(useMC);
1365 if (checkSV) {
1366 dataRequest->requestSecondaryVertices(useMC);
1367 }
1368 o2::tpc::VDriftHelper::requestCCDBInputs(dataRequest->inputs);
1369 dataRequest->inputs.emplace_back("corrMap", o2::header::gDataOriginTPC, "TPCCORRMAP", 0, Lifetime::Timeframe);
1370 auto ggRequest = std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
1371 true, // GRPECS=true
1372 true, // GRPLHCIF
1373 true, // GRPMagField
1374 true, // askMatLUT
1376 dataRequest->inputs,
1377 true);
1378
1379 return DataProcessorSpec{
1380 "track-mc-study",
1381 dataRequest->inputs,
1382 outputs,
1383 AlgorithmSpec{adaptFromTask<TrackMCStudy>(dataRequest, ggRequest, srcTracks, checkSV)},
1384 opts};
1385}
1386
1387} // namespace o2::trackstudy
std::vector< std::string > labels
Defintions for N-prongs secondary vertex fit.
Wrapper container for different reconstructed object types.
Definition of the GeometryManager class.
std::ostringstream debug
Definition of the FIT RecPoints class.
int32_t i
o2::raw::RawFileWriter * raw
Helper for geometry and GRP related CCDB requests.
Global index for barrel track: provides provenance (detectors combination), index in respective array...
Definition of the GeometryTGeo class.
Utility functions for MC particles.
Configurable params for TPC ITS matching.
Definition of the Names Generator class.
Definition of the parameter class for the detector electronics.
uint16_t pos
Definition RawData.h:3
uint32_t j
Definition RawData.h:0
uint32_t res
Definition RawData.h:0
Wrapper container for different reconstructed object types.
o2::track::TrackParCov TrackParCov
Definition Recon.h:39
Configurable params for secondary vertexer.
POD correction map.
Result of refitting TPC-ITS matched track.
Reference on ITS/MFT clusters set.
Helper class to extract VDrift from different sources.
Referenc on track indices contributing to the vertex, with possibility chose tracks from specific sou...
bool isFake() const
Definition MCCompLabel.h:78
void setFakeFlag(bool v=true)
int getTrackID() const
int getSourceID() const
int getEventID() const
std::string asString() const
bool isSet() const
int getEventID() const
Double_t GetStartVertexMomentumZ() const
Definition MCTrack.h:81
Double_t GetStartVertexMomentumX() const
Definition MCTrack.h:79
bool isPrimary() const
Definition MCTrack.h:75
Double_t GetStartVertexCoordinatesY() const
Definition MCTrack.h:83
Double_t GetPt() const
Definition MCTrack.h:109
Double_t GetStartVertexCoordinatesZ() const
Definition MCTrack.h:84
Double_t R2() const
production radius squared
Definition MCTrack.h:88
Double_t GetStartVertexMomentumY() const
Definition MCTrack.h:80
Double_t GetTgl() const
Definition MCTrack.h:142
Double_t GetStartVertexCoordinatesX() const
Definition MCTrack.h:82
Int_t GetPdgCode() const
Accessors.
Definition MCTrack.h:72
Int_t getMotherTrackId() const
Definition MCTrack.h:73
static TDatabasePDG * Instance()
void checkUpdates(o2::framework::ProcessingContext &pc)
static GRPGeomHelper & instance()
void setRequest(std::shared_ptr< GRPGeomRequest > req)
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:178
void setDCA(float d)
Definition V0.h:47
Static class with identifiers, bitmasks and names for ALICE detectors.
Definition DetID.h:58
static constexpr ID ITS
Definition DetID.h:63
static constexpr ID TRD
Definition DetID.h:65
static constexpr ID TPC
Definition DetID.h:64
static constexpr ID TOF
Definition DetID.h:66
ConfigParamRegistry const & options()
Definition InitContext.h:33
decltype(auto) get(R binding, int part=0) const
InputRecord & inputs()
The inputs associated with this processing context.
static GeometryTGeo * Instance()
void fillMatrixCache(int mask) override
static constexpr int getLayer(int chipSW)
static bool isPhysicalPrimary(o2::MCTrack const &p, std::vector< o2::MCTrack > const &pcontainer)
Definition MCUtils.cxx:73
std::vector< o2::InteractionTimeRecord > & getEventRecords(bool withQED=false)
bool initFromDigitContext(std::string_view filename)
DigitizationContext const * getDigitizationContext() const
size_t getNEvents(int source) const
Get number of events.
o2::dataformats::MCEventHeader const & getMCEventHeader(int source, int event) const
retrieves the MCEventHeader for a given eventID and sourceID
size_t getNSources() const
Get number of sources.
std::vector< MCTrack > const & getTracks(int source, int event) const
variant returning all tracks for source and event at once
static void requestCCDBInputs(std::vector< o2::framework::InputSpec > &inputs, bool laser=true, bool itstpcTgl=true)
void extractCCDBInputs(o2::framework::ProcessingContext &pc, bool laser=true, bool itstpcTgl=true)
const VDriftCorrFact & getVDriftObject() const
void endOfStream(EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
void init(InitContext &ic) final
void finaliseCCDB(ConcreteDataMatcher &matcher, void *obj) final
void run(ProcessingContext &pc) final
void process(const o2::globaltracking::RecoContainer &recoData)
TrackMCStudy(std::shared_ptr< DataRequest > dr, std::shared_ptr< o2::base::GRPGeomRequest > gr, GTrackID::mask_t src, bool checkSV)
~TrackMCStudy() final=default
GLenum src
Definition glcorearb.h:1767
GLuint entry
Definition glcorearb.h:5735
GLdouble f
Definition glcorearb.h:310
GLenum const GLfloat * params
Definition glcorearb.h:272
GLfloat v0
Definition glcorearb.h:811
GLuint id
Definition glcorearb.h:650
@ ITSClusters
constexpr o2::header::DataOrigin gDataOriginTPC
Definition DataHeader.h:576
constexpr double LHCBunchSpacingMUS
constexpr int LHCMaxBunches
constexpr double LHCBunchSpacingNS
constexpr float TwoPI
constexpr float PI
constexpr float PIHalf
Node par(int index)
Parameters.
Defining ITS Vertex explicitly as messageable.
Definition Cartesian.h:288
std::vector< ConfigParamSpec > Options
void convertCompactClusters(gsl::span< const itsmft::CompClusterExt > clusters, gsl::span< const unsigned char >::iterator &pattIt, std::vector< o2::BaseCluster< float > > &output, const itsmft::TopologyDictionary *dict)
convert compact clusters to 3D spacepoints
Definition IOUtils.cxx:35
int int int float float float int nCl
detail::Bracket< float > Bracketf_t
Definition Primitive2D.h:40
float angle2Alpha(float phi)
Definition Utils.h:203
int angle2Sector(float phi)
Definition Utils.h:183
std::tuple< float, float > rotateZInv(float xG, float yG, float snAlp, float csAlp)
Definition Utils.h:142
TrackParCovF TrackParCov
Definition Track.h:33
TrackParF TrackPar
Definition Track.h:29
std::pair< int, o2::dataformats::VtxTrackIndex > VTIndexV
o2::framework::DataProcessorSpec getTrackMCStudySpec(o2::dataformats::GlobalTrackID::mask_t srcTracks, o2::dataformats::GlobalTrackID::mask_t srcClus, bool checkSV)
create a processor spec
o2::dataformats::VtxTrackRef V2TRef
o2::dataformats::VtxTrackIndex VTIndex
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
GTrackID getITSContributorGID(GTrackID source) const
GlobalIDSet getSingleDetectorRefs(GTrackID gidx) const
const o2::tpc::TrackTPC & getTPCTrack(GTrackID id) const
const o2::tpc::ClusterNativeAccess & getTPCClusters() const
o2::MCCompLabel getTrackMCLabel(GTrackID id) const
GTrackID getTPCContributorGID(GTrackID source) const
const o2::track::TrackParCov & getTrackParam(GTrackID gidx) const
gsl::span< const unsigned char > clusterShMapTPC
externally set TPC clusters sharing map
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
const o2::track::TrackParCov & getTrackParamOut(GTrackID gidx) const
const o2::dataformats::PrimaryVertex & getPrimaryVertex(int i) const
const o2::its::TrackITS & getITSTrack(GTrackID gid) const
std::unique_ptr< o2::tpc::internal::getWorkflowTPCInput_ret > inputsTPCclusters
void getTrackTime(GTrackID gid, float &t, float &tErr) const
gsl::span< const unsigned int > occupancyMapTPC
externally set TPC clusters occupancy map
static constexpr int T2L
Definition Cartesian.h:55
static constexpr int T2GRot
Definition Cartesian.h:57
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
std::vector< int > row