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