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