Project
Loading...
Searching...
No Matches
TrackingStudy.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>
45#include "GPUO2InterfaceRefit.h"
46#include "GPUO2ExternalUser.h" // Needed for propper settings in GPUParam.h
47#include "GPUParam.h"
48#include "GPUParam.inc"
49#include "GPUTPCGeometry.h"
51#include "MathUtils/fit.h"
52#include "TPCFastTransformPOD.h"
53#include <TF1.h>
54
55namespace o2::trackstudy
56{
57
58using namespace o2::framework;
61
67
69
70class TrackingStudySpec final : public Task
71{
72 public:
73 TrackingStudySpec(std::shared_ptr<DataRequest> dr, std::shared_ptr<o2::base::GRPGeomRequest> gr, GTrackID::mask_t src, bool useMC)
74 : mDataRequest(dr), mGGCCDBRequest(gr), mTracksSrc(src), mUseMC(useMC) {}
75 ~TrackingStudySpec() final = default;
76 void init(InitContext& ic) final;
77 void run(ProcessingContext& pc) final;
78 void endOfStream(EndOfStreamContext& ec) final;
79 void finaliseCCDB(ConcreteDataMatcher& matcher, void* obj) final;
80 void process(o2::globaltracking::RecoContainer& recoData);
81
82 private:
83 void updateTimeDependentParams(ProcessingContext& pc);
84 float getDCAYCut(float pt) const;
85 float getDCAZCut(float pt) const;
86 std::shared_ptr<DataRequest> mDataRequest;
87 std::shared_ptr<o2::base::GRPGeomRequest> mGGCCDBRequest;
88 o2::tpc::VDriftHelper mTPCVDriftHelper{};
89 const o2::gpu::TPCFastTransformPOD* mTPCCorrMaps{nullptr};
90 bool mUseMC{false};
91 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOut;
92 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOutVtx;
93 std::unique_ptr<o2::gpu::GPUO2InterfaceRefit> mTPCRefitter;
94 std::vector<float> mMltHistTB, mTBinClOccAft, mTBinClOccBef, mTBinClOccWgh;
95 std::unique_ptr<TF1> mOccWghFun;
96 float mITSROFrameLengthMUS = 0.f;
97 float mTPCTBinMUS = 0.f; // TPC bin in microseconds
98 float mTPCTBinMUSInv = 0.f;
99 int mMaxNeighbours = 3;
100 float mMaxVTTimeDiff = 80.; // \mus
101 float mTPCDCAYCut = 2.;
102 float mTPCDCAZCut = 2.;
103 float mMinX = 46.;
104 float mMaxEta = 0.8;
105 float mMinPt = 0.1;
106 int mNOccBinsDrift = 10;
107 int mMinTPCClusters = 60;
108 int mNTPCOccBinLength = 0;
109 int mNHBPerTF = 0;
110 float mNTPCOccBinLengthInv;
111 bool mStoreWithITSOnly = false;
112 bool mDoPairsCorr = false;
113 std::string mDCAYFormula = "0.0105 + 0.0350 / pow(x, 1.1)";
114 std::string mDCAZFormula = "0.0105 + 0.0350 / pow(x, 1.1)";
115 GTrackID::mask_t mTracksSrc{};
117 o2::steer::MCKinematicsReader mcReader; // reader of MC information
118};
119
121{
123 int lane = ic.services().get<const o2::framework::DeviceSpec>().inputTimesliceId;
124 int maxLanes = ic.services().get<const o2::framework::DeviceSpec>().maxInputTimeslices;
125 std::string dbgnm = maxLanes == 1 ? "trackStudy.root" : fmt::format("trackStudy_{}.root", lane);
126 mDBGOut = std::make_unique<o2::utils::TreeStreamRedirector>(dbgnm.c_str(), "recreate");
127 dbgnm = maxLanes == 1 ? "trackStudyVtx.root" : fmt::format("trackStudyVtx_{}.root", lane);
128 mDBGOutVtx = std::make_unique<o2::utils::TreeStreamRedirector>(dbgnm.c_str(), "recreate");
129 mStoreWithITSOnly = ic.options().get<bool>("with-its-only");
130 mMaxVTTimeDiff = ic.options().get<float>("max-vtx-timediff");
131 mMaxNeighbours = ic.options().get<int>("max-vtx-neighbours");
132 mTPCDCAYCut = ic.options().get<float>("max-tpc-dcay");
133 mTPCDCAZCut = ic.options().get<float>("max-tpc-dcaz");
134 mMinX = ic.options().get<float>("min-x-prop");
135 mMaxEta = ic.options().get<float>("max-eta");
136 mMinPt = ic.options().get<float>("min-pt");
137 mMinTPCClusters = ic.options().get<int>("min-tpc-clusters");
138 mDCAYFormula = ic.options().get<std::string>("dcay-vs-pt");
139 mDCAZFormula = ic.options().get<std::string>("dcaz-vs-pt");
140 mDoPairsCorr = ic.options().get<bool>("pair-correlations");
141 mNOccBinsDrift = ic.options().get<int>("noccbins");
142 if (mNOccBinsDrift < 3) {
143 mNOccBinsDrift = 3;
144 }
145 auto str = ic.options().get<std::string>("occ-weight-fun");
146 if (!str.empty()) {
147 mOccWghFun = std::make_unique<TF1>("occFun", str.c_str(), -100., 100.);
148 }
149}
150
152{
154 recoData.collectData(pc, *mDataRequest.get()); // select tracks of needed type, with minimal cuts, the real selected will be done in the vertexer
155 updateTimeDependentParams(pc); // Make sure this is called after recoData.collectData, which may load some conditions
156 if (recoData.inputsTPCclusters) {
157 mTPCRefitter = std::make_unique<o2::gpu::GPUO2InterfaceRefit>(&recoData.inputsTPCclusters->clusterIndex, mTPCCorrMaps, o2::base::Propagator::Instance()->getNominalBz(),
158 recoData.getTPCTracksClusterRefs().data(), 0, recoData.clusterShMapTPC.data(), recoData.occupancyMapTPC.data(),
159 recoData.occupancyMapTPC.size(), nullptr, o2::base::Propagator::Instance());
160 mNTPCOccBinLength = mTPCRefitter->getParam()->rec.tpc.occupancyMapTimeBins;
161 mTBinClOccBef.clear();
162 mTBinClOccAft.clear();
163 mTBinClOccWgh.clear();
164 }
165
166 // prepare TPC occupancy data
167 if (mNTPCOccBinLength > 1 && recoData.occupancyMapTPC.size()) {
168 mNTPCOccBinLengthInv = 1. / mNTPCOccBinLength;
169 int nTPCBins = mNHBPerTF * o2::constants::lhc::LHCMaxBunches / 8, ninteg = 0;
170 int nTPCOccBins = nTPCBins * mNTPCOccBinLengthInv, sumBins = std::max(1, int(o2::constants::lhc::LHCMaxBunches / 8 * mNTPCOccBinLengthInv));
171 mTBinClOccAft.resize(nTPCOccBins);
172 mTBinClOccBef.resize(nTPCOccBins);
173 float sm = 0., tb = 0.5 * mNTPCOccBinLength;
174 mMltHistTB.resize(nTPCOccBins);
175 for (int i = 0; i < nTPCOccBins; i++) {
176 mMltHistTB[i] = mTPCRefitter->getParam()->GetUnscaledMult(tb);
177 tb += mNTPCOccBinLength;
178 }
179 for (int i = nTPCOccBins; i--;) {
180 sm += mMltHistTB[i];
181 if (i + sumBins < nTPCOccBins) {
182 sm -= mMltHistTB[i + sumBins];
183 }
184 mTBinClOccAft[i] = sm;
185 }
186 sm = 0;
187 for (int i = 0; i < nTPCOccBins; i++) {
188 sm += mMltHistTB[i];
189 if (i - sumBins > 0) {
190 sm -= mMltHistTB[i - sumBins];
191 }
192 mTBinClOccBef[i] = sm;
193 }
194 } else {
195 mTBinClOccBef.resize(1);
196 mTBinClOccAft.resize(1);
197 }
198
199 process(recoData);
200}
201
202void TrackingStudySpec::updateTimeDependentParams(ProcessingContext& pc)
203{
205 mTPCVDriftHelper.extractCCDBInputs(pc);
206 auto const& raw = pc.inputs().get<const char*>("corrMap");
207 mTPCCorrMaps = &o2::gpu::TPCFastTransformPOD::get(raw);
208 static bool initOnceDone = false;
209 if (!initOnceDone) { // this params need to be queried only once
210 initOnceDone = true;
211 // Note: reading of the ITS AlpideParam needed for ITS timing is done by the RecoContainer
214 if (!grp->isDetContinuousReadOut(DetID::ITS)) {
215 mITSROFrameLengthMUS = alpParams.roFrameLengthTrig / 1.e3; // ITS ROFrame duration in \mus
216 } else {
217 mITSROFrameLengthMUS = alpParams.roFrameLengthInBC * o2::constants::lhc::LHCBunchSpacingNS * 1e-3; // ITS ROFrame duration in \mus
218 }
220 mNHBPerTF = o2::base::GRPGeomHelper::instance().getGRPECS()->getNHBFPerTF();
222 mTPCTBinMUS = elParam.ZbinWidth; // TPC bin in microseconds
223 mTPCTBinMUSInv = 1. / mTPCTBinMUS;
224 }
225}
226
228{
229 auto pvvec = recoData.getPrimaryVertices();
230 auto trackIndex = recoData.getPrimaryVertexMatchedTracks(); // Global ID's for associated tracks
231 auto vtxRefs = recoData.getPrimaryVertexMatchedTrackRefs(); // references from vertex to these track IDs
232 auto prop = o2::base::Propagator::Instance();
233 auto FITInfo = recoData.getFT0RecPoints();
234 static int TFCount = 0;
235 int nv = vtxRefs.size();
237 o2::dataformats::PrimaryVertexExt vtxDummy(mMeanVtx.getPos(), {}, {}, 0);
238 std::vector<o2::dataformats::PrimaryVertexExt> pveVec(nv);
239 std::vector<float> tpcOccAftV, tpcOccBefV;
240 pveVec.back() = vtxDummy;
242 float tBiasITS = alpParams.roFrameBiasInBC * o2::constants::lhc::LHCBunchSpacingMUS;
244 std::vector<o2::dataformats::TrackInfoExt> trcExtVec;
245 std::vector<o2::trackstudy::TrackPairInfo> trcPairsVec;
246 auto vdrift = mTPCVDriftHelper.getVDriftObject().getVDrift();
247 float maxDriftTB = 250.f / vdrift / (o2::constants::lhc::LHCBunchSpacingMUS * 8);
248 int groupOcc = std::ceil(maxDriftTB / mNOccBinsDrift / mNTPCOccBinLength);
249
250 bool tpcTrackOK = recoData.isTrackSourceLoaded(GTrackID::TPC);
251
252 auto fillTPCClInfo = [&recoData, this](const o2::tpc::TrackTPC& trc, o2::dataformats::TrackInfoExt& trExt, float timestampTB = -1e9) {
253 const auto clRefs = recoData.getTPCTracksClusterRefs();
254 const auto tpcClusAcc = recoData.getTPCClusters();
255 const auto shMap = recoData.clusterShMapTPC;
256
257 if (recoData.inputsTPCclusters) {
258 uint8_t clSect = 0, clRow = 0, lowestR = -1;
259 uint32_t clIdx = 0;
260 for (int ic = 0; ic < trc.getNClusterReferences(); ic++) { // outside -> inside ordering, but on the sector boundaries backward jumps are possible
261 trc.getClusterReference(clRefs, ic, clSect, clRow, clIdx);
262 if (clRow < lowestR) {
263 trExt.rowCountTPC++;
264 lowestR = clRow;
265 }
266 unsigned int absoluteIndex = tpcClusAcc.clusterOffset[clSect][clRow] + clIdx;
267 if (shMap[absoluteIndex] & o2::gpu::GPUTPCGMMergedTrackHit::flagShared) {
268 trExt.nClTPCShared++;
269 }
270 }
271 trExt.rowMinTPC = lowestR;
272 const auto& clus = tpcClusAcc.clusters[clSect][clRow][clIdx];
273 trExt.padFromEdge = uint8_t(clus.getPad());
274 int npads = o2::gpu::GPUTPCGeometry::NPads(lowestR);
275 if (trExt.padFromEdge > npads / 2) {
276 trExt.padFromEdge = npads - 1 - trExt.padFromEdge;
277 }
278 this->mTPCCorrMaps->Transform(clSect, clRow, clus.getPad(), clus.getTime(), trExt.innerTPCPos0[0], trExt.innerTPCPos0[1], trExt.innerTPCPos0[2], trc.getTime0()); // nominal time of the track
279 if (timestampTB > -1e8) {
280 this->mTPCCorrMaps->Transform(clSect, clRow, clus.getPad(), clus.getTime(), trExt.innerTPCPos[0], trExt.innerTPCPos[1], trExt.innerTPCPos[2], timestampTB); // time assigned from the global track track
281 } else {
282 trExt.innerTPCPos = trExt.innerTPCPos0;
283 }
284 trc.getClusterReference(clRefs, 0, clSect, clRow, clIdx);
285 trExt.rowMaxTPC = clRow;
286 }
287 };
288
289 auto getTPCPairSharing = [&recoData, this](const o2::tpc::TrackTPC& trc0, const o2::tpc::TrackTPC& trc1) {
290 const auto clRefs = recoData.getTPCTracksClusterRefs();
291 uint8_t nsh = 0, nshRows = 0, lastSharedRow = -1;
292 if (recoData.inputsTPCclusters) {
293 uint8_t clSect0 = 0, clRow0 = 0, clSect1 = 0, clRow1 = 0;
294 uint32_t clIdx0 = 0, clIdx1 = 0;
295 int ic1Start = 0;
296 for (int ic0 = 0; ic0 < trc0.getNClusterReferences(); ic0++) { // outside -> inside, but on the sector boundaries backward jumps are possible
297 trc0.getClusterReference(clRefs, ic0, clSect0, clRow0, clIdx0);
298 for (int ic1 = ic1Start; ic1 < trc1.getNClusterReferences(); ic1++) { // outside -> inside, but on the sector boundaries backward jumps are possible
299 trc1.getClusterReference(clRefs, ic1, clSect1, clRow1, clIdx1);
300 if (clRow1 > clRow0) {
301 ic1Start = ic1 + 1;
302 continue; // catch up ic0
303 }
304 if (clRow1 == clRow0) {
305 if (clSect0 == clSect1 && clIdx0 == clIdx1) {
306 nsh++;
307 if (lastSharedRow != clRow0) {
308 lastSharedRow = clRow0;
309 nshRows++;
310 }
311 ic1Start = ic1 + 1;
312 break; // check next ic0
313 }
314 }
315 }
316 }
317 }
318 return std::make_pair(nsh, nshRows);
319 };
320
321 auto assignRecTrack = [&recoData, this](const o2::dataformats::TrackInfoExt& src, o2::trackstudy::RecTrack& dst) {
322 dst.track = src.track;
323 dst.gid = src.gid;
324 dst.ts.setTimeStamp(src.ttime);
325 dst.ts.setTimeStampError(src.ttimeE);
326 dst.nClITS = src.nClITS;
327 dst.nClTPC = src.nClTPC;
328 dst.pattITS = src.pattITS;
329 if (src.q2ptITS == 0. && dst.nClITS > 0) {
330 dst.pattITS |= 0x1 << 7;
331 }
332 dst.lowestPadRow = src.rowMinTPC;
333 if (this->mUseMC) {
334 auto gidSet = recoData.getSingleDetectorRefs(src.gid);
335 if (recoData.getTrackMCLabel(src.gid).isFake()) {
336 dst.flags |= RecTrack::FakeGLO;
337 }
338 auto msk = src.gid.getSourceDetectorsMask();
339 if (msk[DetID::ITS]) {
340 if (gidSet[GTrackID::ITS].isSourceSet()) { // has ITS track rather than AB tracklet
341 auto lblITS = recoData.getTrackMCLabel(gidSet[GTrackID::ITS]);
342 if (lblITS.isFake()) {
343 dst.flags |= RecTrack::FakeITS;
344 }
345 } else { // AB ITS tracklet
346 if (recoData.getTrackMCLabel(gidSet[GTrackID::ITSAB]).isFake()) {
347 dst.flags |= RecTrack::FakeITS;
348 }
349 }
350 if (msk[DetID::TPC]) { // has both ITS and TPC contribution
351 if (recoData.getTrackMCLabel(gidSet[GTrackID::ITSTPC]).isFake()) {
352 dst.flags |= RecTrack::FakeITSTPC;
353 }
354 }
355 }
356 if (msk[DetID::TPC]) {
357 if (recoData.getTrackMCLabel(gidSet[GTrackID::TPC]).isFake()) {
358 dst.flags |= RecTrack::FakeTPC;
359 }
360 }
361 }
362 };
363 tpcOccAftV.resize(mNOccBinsDrift);
364 tpcOccBefV.resize(mNOccBinsDrift);
365
366 for (int iv = 0; iv < nv; iv++) {
367 LOGP(debug, "processing PV {} of {}", iv, nv);
368 const auto& vtref = vtxRefs[iv];
369 if (iv != nv - 1) {
370 auto& pve = pveVec[iv];
371 static_cast<o2::dataformats::PrimaryVertex&>(pve) = pvvec[iv];
372 // find best matching FT0 signal
373 float bestTimeDiff = 1000, bestTime = -999;
374 int bestFTID = -1;
375 if (mTracksSrc[GTrackID::FT0]) {
376 for (int ift0 = vtref.getFirstEntryOfSource(GTrackID::FT0); ift0 < vtref.getFirstEntryOfSource(GTrackID::FT0) + vtref.getEntriesOfSource(GTrackID::FT0); ift0++) {
377 const auto& ft0 = FITInfo[trackIndex[ift0]];
378 if (ft0Params.isSelected(ft0)) {
379 auto fitTime = ft0.getInteractionRecord().differenceInBCMUS(recoData.startIR);
380 if (std::abs(fitTime - pve.getTimeStamp().getTimeStamp()) < bestTimeDiff) {
381 bestTimeDiff = fitTime - pve.getTimeStamp().getTimeStamp();
382 bestFTID = trackIndex[ift0];
383 }
384 }
385 }
386 } else {
387 LOGP(warn, "FT0 is not requested, cannot set complete vertex info");
388 }
389 if (bestFTID >= 0) {
390 pve.FT0A = FITInfo[bestFTID].getTrigger().getAmplA();
391 pve.FT0C = FITInfo[bestFTID].getTrigger().getAmplC();
392 pve.FT0Time = double(FITInfo[bestFTID].getInteractionRecord().differenceInBCMUS(recoData.startIR)) + FITInfo[bestFTID].getCollisionTimeMean() * 1e-6; // time in \mus
393 }
394 pve.VtxID = iv;
395 }
396 trcExtVec.clear();
397 trcPairsVec.clear();
398 float q2ptITS, q2ptTPC, q2ptITSTPC, q2ptITSTPCTRD;
399 for (int is = 0; is < GTrackID::NSources; is++) {
400 DetID::mask_t dm = GTrackID::getSourceDetectorsMask(is);
401 bool skipTracks = !mTracksSrc[is] || !recoData.isTrackSourceLoaded(is) || !(dm[DetID::ITS] || dm[DetID::TPC]);
402 int idMin = vtref.getFirstEntryOfSource(is), idMax = idMin + vtref.getEntriesOfSource(is);
403 for (int i = idMin; i < idMax; i++) {
404 auto vid = trackIndex[i];
405 bool pvCont = vid.isPVContributor();
406 if (pvCont) {
407 pveVec[iv].nSrc[is]++;
408 }
409 if (skipTracks) {
410 continue;
411 }
412 GTrackID tpcTrID;
413 const o2::tpc::TrackTPC* tpcTr = nullptr;
414 int nclTPC = 0;
415 if (dm[DetID::TPC] && tpcTrackOK) {
416 tpcTrID = recoData.getTPCContributorGID(vid);
417 tpcTr = &recoData.getTPCTrack(tpcTrID);
418 nclTPC = tpcTr->getNClusters();
419 if (nclTPC < mMinTPCClusters) {
420 continue;
421 }
422 }
423 bool ambig = vid.isAmbiguous();
424 auto trc = recoData.getTrackParam(vid);
425 if (fabs(trc.getEta()) > mMaxEta) {
426 continue;
427 }
428 if (iv < nv - 1 && is == GTrackID::TPC && tpcTr && !tpcTr->hasBothSidesClusters()) { // for unconstrained TPC tracks correct track Z
429 float corz = vdrift * (tpcTr->getTime0() * mTPCTBinMUS - pvvec[iv].getTimeStamp().getTimeStamp());
430 if (tpcTr->hasASideClustersOnly()) {
431 corz = -corz; // A-side
432 }
433 trc.setZ(trc.getZ() + corz);
434 }
435 float xmin = trc.getX();
437 if (!prop->propagateToDCA(iv == nv - 1 ? vtxDummy : pvvec[iv], trc, prop->getNominalBz(), 2., o2::base::PropagatorF::MatCorrType::USEMatCorrLUT, &dca)) {
438 continue;
439 }
440 bool hasITS = GTrackID::getSourceDetectorsMask(is)[GTrackID::ITS];
441 if (std::abs(dca.getY()) > (hasITS ? getDCAYCut(trc.getPt()) : mTPCDCAYCut) ||
442 std::abs(dca.getZ()) > (hasITS ? getDCAZCut(trc.getPt()) : mTPCDCAZCut)) {
443 continue;
444 }
445 if (trc.getPt() < mMinPt) {
446 continue;
447 }
448 if (iv != nv - 1) {
449 pveVec[iv].nSrcA[is]++;
450 if (ambig) {
451 pveVec[iv].nSrcAU[is]++;
452 }
453 }
454 if (!hasITS && mStoreWithITSOnly) {
455 continue;
456 }
457 {
458 auto& trcExt = trcExtVec.emplace_back();
459 recoData.getTrackTime(vid, trcExt.ttime, trcExt.ttimeE);
460 trcExt.track = trc;
461 trcExt.hashIU = trc.hash();
462 trcExt.dca = dca;
463 trcExt.gid = vid;
464 trcExt.xmin = xmin;
465 trcExt.dcaTPC.set(-999.f, -999.f);
466
467 if (tpcTr) {
468 float tsuse = trcExt.ttime / (8 * o2::constants::lhc::LHCBunchSpacingMUS);
469 if (tpcTr->hasASideClusters()) {
470 trcExt.setTPCA();
471 }
472 if (tpcTr->hasCSideClusters()) {
473 trcExt.setTPCC();
474 }
475 if (is == GTrackID::TPC) {
476 trcExt.dcaTPC = dca;
477 tsuse = -1e9;
478 } else {
479 o2::track::TrackParCov tmpTPC(*tpcTr);
480 if (iv < nv - 1 && is == GTrackID::TPC && tpcTr && !tpcTr->hasBothSidesClusters()) { // for unconstrained TPC tracks correct track Z
481 float corz = vdrift * (tpcTr->getTime0() * mTPCTBinMUS - pvvec[iv].getTimeStamp().getTimeStamp());
482 if (tpcTr->hasASideClustersOnly()) {
483 corz = -corz; // A-side
484 }
485 tmpTPC.setZ(tmpTPC.getZ() + corz);
486 }
487 if (!prop->propagateToDCA(iv == nv - 1 ? vtxDummy : pvvec[iv], tmpTPC, prop->getNominalBz(), 2., o2::base::PropagatorF::MatCorrType::USEMatCorrLUT, &trcExt.dcaTPC)) {
488 trcExt.dcaTPC.set(-999.f, -999.f);
489 }
490 }
491 fillTPCClInfo(*tpcTr, trcExt, tsuse);
492 trcExt.chi2TPC = tpcTr->getChi2();
493 }
494 auto gidRefs = recoData.getSingleDetectorRefs(vid);
495 if (gidRefs[GTrackID::ITS].isIndexSet()) {
496 const auto& itsTr = recoData.getITSTrack(gidRefs[GTrackID::ITS]);
497 trcExt.q2ptITS = itsTr.getQ2Pt();
498 trcExt.nClITS = itsTr.getNClusters();
499 for (int il = 0; il < 7; il++) {
500 if (itsTr.hasHitOnLayer(il)) {
501 trcExt.pattITS |= 0x1 << il;
502 }
503 }
504 } else if (gidRefs[GTrackID::ITSAB].isIndexSet()) {
505 const auto& itsTrf = recoData.getITSABRefs()[gidRefs[GTrackID::ITSAB]];
506 trcExt.nClITS = itsTrf.getNClusters();
507 for (int il = 0; il < 7; il++) {
508 if (itsTrf.hasHitOnLayer(il)) {
509 trcExt.pattITS |= 0x1 << il;
510 }
511 }
512 }
513 if (gidRefs[GTrackID::TPC].isIndexSet()) {
514 trcExt.q2ptTPC = recoData.getTrackParam(gidRefs[GTrackID::TPC]).getQ2Pt();
515 trcExt.nClTPC = nclTPC;
516 }
517 if (gidRefs[GTrackID::ITSTPC].isIndexSet()) {
518 const auto& trTPCITS = recoData.getTPCITSTrack(gidRefs[GTrackID::ITSTPC]);
519 trcExt.q2ptITSTPC = trTPCITS.getQ2Pt();
520 trcExt.chi2ITSTPC = trTPCITS.getChi2Match();
521 }
522 if (gidRefs[GTrackID::TRD].isIndexSet()) {
523 trcExt.q2ptITSTPCTRD = recoData.getTrackParam(gidRefs[GTrackID::TRD]).getQ2Pt();
524 }
525 if (gidRefs[GTrackID::TOF].isIndexSet()) {
526 trcExt.infoTOF = recoData.getTOFMatch(vid);
527 }
528 }
529 }
530 }
531 float tpcOccBef = 0., tpcOccAft = 0.;
532 if (iv != nv - 1) {
533 int tb = pveVec[iv].getTimeStamp().getTimeStamp() * mTPCTBinMUSInv * mNTPCOccBinLengthInv;
534 tpcOccBef = tb < 0 ? mTBinClOccBef[0] : (tb >= mTBinClOccBef.size() ? mTBinClOccBef.back() : mTBinClOccBef[tb]);
535 tpcOccAft = tb < 0 ? mTBinClOccAft[0] : (tb >= mTBinClOccAft.size() ? mTBinClOccAft.back() : mTBinClOccAft[tb]);
536 int tbc = pveVec[iv].getTimeStamp().getTimeStamp() * mTPCTBinMUSInv * mNTPCOccBinLengthInv - groupOcc / 2.;
537 for (int iob = 0; iob < mNOccBinsDrift; iob++) {
538 float sm = 0;
539 for (int ig = 0; ig < groupOcc; ig++) {
540 int ocb = tbc + ig + groupOcc * iob;
541 if (ocb < 0 || ocb >= (int)mMltHistTB.size()) {
542 sm = -1;
543 break;
544 }
545 sm += mMltHistTB[ocb];
546 }
547 tpcOccAftV[iob] = sm;
548 //
549 sm = 0;
550 for (int ig = 0; ig < groupOcc; ig++) {
551 int ocb = tbc + ig - groupOcc * iob;
552 if (ocb < 0 || ocb >= (int)mMltHistTB.size()) {
553 sm = -1;
554 break;
555 }
556 sm += mMltHistTB[ocb];
557 }
558 tpcOccBefV[iob] = sm;
559 }
560 }
561 (*mDBGOut) << "trpv"
562 << "orbit=" << recoData.startIR.orbit << "tfID=" << TFCount
563 << "tpcOccBef=" << tpcOccBef << "tpcOccAft=" << tpcOccAft
564 << "tpcOccBefV=" << tpcOccBefV << "tpcOccAftV=" << tpcOccAftV
565 << "pve=" << pveVec[iv] << "trc=" << trcExtVec << "\n";
566
567 if (mDoPairsCorr) {
568 for (int it0 = 0; it0 < (int)trcExtVec.size(); it0++) {
569 const auto& tr0 = trcExtVec[it0];
570 if (tr0.nClTPC < 1) {
571 continue;
572 }
573 for (int it1 = it0 + 1; it1 < (int)trcExtVec.size(); it1++) {
574 const auto& tr1 = trcExtVec[it1];
575 if (tr1.nClTPC < 1) {
576 continue;
577 }
578
579 if (std::abs(tr0.track.getTgl() - tr1.track.getTgl()) > 0.25) {
580 continue;
581 }
582 auto dphi = tr0.track.getPhi() - tr1.track.getPhi();
583 if (dphi < -o2::constants::math::PI) {
585 } else if (dphi > o2::constants::math::PI) {
587 }
588 if (std::abs(dphi) > 0.25) {
589 continue;
590 }
591 auto& pr = trcPairsVec.emplace_back();
592 assignRecTrack(tr0, pr.tr0);
593 assignRecTrack(tr1, pr.tr1);
594 auto shinfo = getTPCPairSharing(recoData.getTPCTrack(recoData.getTPCContributorGID(tr0.gid)), recoData.getTPCTrack(recoData.getTPCContributorGID(tr1.gid)));
595 pr.nshTPC = shinfo.first;
596 pr.nshTPCRow = shinfo.second;
597 }
598 }
599 (*mDBGOut) << "pairs" << "pr=" << trcPairsVec << "\n";
600 }
601 }
602
603 int nvtot = mMaxNeighbours < 0 ? -1 : (int)pveVec.size();
604
605 auto insSlot = [maxSlots = mMaxNeighbours](std::vector<float>& vc, float v, int slot, std::vector<int>& vid, int id) {
606 for (int i = maxSlots - 1; i > slot; i--) {
607 std::swap(vc[i], vc[i - 1]);
608 std::swap(vid[i], vid[i - 1]);
609 }
610 vc[slot] = v;
611 vid[slot] = id;
612 };
613
614 for (int cnt = 0; cnt < nvtot; cnt++) {
615 const auto& pve = pveVec[cnt];
616 float tv = pve.getTimeStamp().getTimeStamp();
617 std::vector<o2::dataformats::PrimaryVertexExt> pveT(mMaxNeighbours); // neighbours in time
618 std::vector<o2::dataformats::PrimaryVertexExt> pveZ(mMaxNeighbours); // neighbours in Z
619 std::vector<int> idT(mMaxNeighbours), idZ(mMaxNeighbours);
620 std::vector<float> dT(mMaxNeighbours), dZ(mMaxNeighbours);
621 for (int i = 0; i < mMaxNeighbours; i++) {
622 idT[i] = idZ[i] = -1;
623 dT[i] = mMaxVTTimeDiff;
624 dZ[i] = 1e9;
625 }
626 int cntM = cnt - 1, cntP = cnt + 1;
627 for (; cntM >= 0; cntM--) { // backward
628 const auto& vt = pveVec[cntM];
629 auto dtime = std::abs(tv - vt.getTimeStamp().getTimeStamp());
630 if (dtime > mMaxVTTimeDiff) {
631 continue;
632 }
633 for (int i = 0; i < mMaxNeighbours; i++) {
634 if (dT[i] > dtime) {
635 insSlot(dT, dtime, i, idT, cntM);
636 break;
637 }
638 }
639 auto dz = std::abs(pve.getZ() - vt.getZ());
640 for (int i = 0; i < mMaxNeighbours; i++) {
641 if (dZ[i] > dz) {
642 insSlot(dZ, dz, i, idZ, cntM);
643 break;
644 }
645 }
646 }
647 for (; cntP < nvtot; cntP++) { // forward
648 const auto& vt = pveVec[cntP];
649 auto dtime = std::abs(tv - vt.getTimeStamp().getTimeStamp());
650 if (dtime > mMaxVTTimeDiff) {
651 continue;
652 }
653 for (int i = 0; i < mMaxNeighbours; i++) {
654 if (dT[i] > dtime) {
655 insSlot(dT, dtime, i, idT, cntP);
656 break;
657 }
658 }
659 auto dz = std::abs(pve.getZ() - vt.getZ());
660 for (int i = 0; i < mMaxNeighbours; i++) {
661 if (dZ[i] > dz) {
662 insSlot(dZ, dz, i, idZ, cntP);
663 break;
664 }
665 }
666 }
667 for (int i = 0; i < mMaxNeighbours; i++) {
668 if (idT[i] != -1) {
669 pveT[i] = pveVec[idT[i]];
670 } else {
671 break;
672 }
673 }
674 for (int i = 0; i < mMaxNeighbours; i++) {
675 if (idZ[i] != -1) {
676 pveZ[i] = pveVec[idZ[i]];
677 } else {
678 break;
679 }
680 }
681 (*mDBGOutVtx) << "pvExt"
682 << "pve=" << pve
683 << "pveT=" << pveT
684 << "pveZ=" << pveZ
685 << "tfID=" << TFCount
686 << "\n";
687 }
688
689 TFCount++;
690}
691
693{
694 mDBGOut.reset();
695 mDBGOutVtx.reset();
696}
697
699{
701 return;
702 }
703 if (mTPCVDriftHelper.accountCCDBInputs(matcher, obj)) {
704 return;
705 }
706 if (matcher == ConcreteDataMatcher("GLO", "MEANVERTEX", 0)) {
707 LOG(info) << "Imposing new MeanVertex: " << ((const o2::dataformats::MeanVertexObject*)obj)->asString();
708 mMeanVtx = *(const o2::dataformats::MeanVertexObject*)obj;
709 return;
710 }
711}
712
713float TrackingStudySpec::getDCAYCut(float pt) const
714{
715 static TF1 fun("dcayvspt", mDCAYFormula.c_str(), 0, 20);
716 return fun.Eval(pt);
717}
718
719float TrackingStudySpec::getDCAZCut(float pt) const
720{
721 static TF1 fun("dcazvspt", mDCAZFormula.c_str(), 0, 20);
722 return fun.Eval(pt);
723}
724
726{
727 std::vector<OutputSpec> outputs;
728 auto dataRequest = std::make_shared<DataRequest>();
729
730 dataRequest->requestTracks(srcTracks, useMC);
731 dataRequest->requestClusters(srcClusters, useMC);
732 dataRequest->requestPrimaryVertices(useMC);
733 dataRequest->inputs.emplace_back("meanvtx", "GLO", "MEANVERTEX", 0, Lifetime::Condition, ccdbParamSpec("GLO/Calib/MeanVertex", {}, 1));
734 auto ggRequest = std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
735 true, // GRPECS=true
736 true, // GRPLHCIF
737 true, // GRPMagField
738 true, // askMatLUT
740 dataRequest->inputs,
741 true);
742
743 Options opts{
744 {"max-vtx-neighbours", VariantType::Int, 3, {"Max PV neighbours fill, no PV study if < 0"}},
745 {"max-vtx-timediff", VariantType::Float, 90.f, {"Max PV time difference to consider"}},
746 {"dcay-vs-pt", VariantType::String, "0.0105 + 0.0350 / pow(x, 1.1)", {"Formula for global tracks DCAy vs pT cut"}},
747 {"dcaz-vs-pt", VariantType::String, "0.0105 + 0.0350 / pow(x, 1.1)", {"Formula for global tracks DCAy vs pT cut"}},
748 {"min-tpc-clusters", VariantType::Int, 60, {"Cut on TPC clusters"}},
749 {"max-tpc-dcay", VariantType::Float, 5.f, {"Cut on TPC dcaY"}},
750 {"max-tpc-dcaz", VariantType::Float, 5.f, {"Cut on TPC dcaZ"}},
751 {"max-eta", VariantType::Float, 1.0f, {"Cut on track eta"}},
752 {"min-pt", VariantType::Float, 0.1f, {"Cut on track pT"}},
753 {"with-its-only", VariantType::Bool, false, {"Store tracks with ITS only"}},
754 {"pair-correlations", VariantType::Bool, false, {"Do pairs correlation"}},
755 {"occ-weight-fun", VariantType::String, "(x>=-40&&x<-5) ? (1./1225*pow(x+40,2)) : ((x>-5&&x<15) ? 1. : ((x>=15&&x<40) ? (-0.4/25*x+1.24 ) : ( (x>40&&x<100) ? -0.4/60*x+0.6+0.8/3 : 0)))", {"Occupancy weighting f-n vs time in musec"}},
756 {"noccbins", VariantType::Int, 10, {"Number of occupancy bins per full drift time"}},
757 {"min-x-prop", VariantType::Float, 100.f, {"track should be propagated to this X at least"}},
758 };
759 o2::tpc::VDriftHelper::requestCCDBInputs(dataRequest->inputs);
760 dataRequest->inputs.emplace_back("corrMap", o2::header::gDataOriginTPC, "TPCCORRMAP", 0, Lifetime::Timeframe);
761
762 return DataProcessorSpec{
763 "track-study",
764 dataRequest->inputs,
765 outputs,
766 AlgorithmSpec{adaptFromTask<TrackingStudySpec>(dataRequest, ggRequest, srcTracks, useMC)},
767 opts};
768}
769
770} // namespace o2::trackstudy
Wrapper container for different reconstructed object types.
Definition of the GeometryManager class.
std::ostringstream debug
Definition of the FIT RecPoints class.
int32_t i
o2::raw::RawFileWriter * raw
Helper for geometry and GRP related CCDB requests.
Global index for barrel track: provides provenance (detectors combination), index in respective array...
Utility functions for MC particles.
Definition of the Names Generator class.
Definition of the parameter class for the detector electronics.
Wrapper container for different reconstructed object types.
POD correction map.
Result of refitting TPC-ITS matched track.
Reference on ITS/MFT clusters set.
Helper class to extract VDrift from different sources.
Referenc on track indices contributing to the vertex, with possibility chose tracks from specific sou...
bool isFake() const
Definition MCCompLabel.h:78
void checkUpdates(o2::framework::ProcessingContext &pc)
static GRPGeomHelper & instance()
void setRequest(std::shared_ptr< GRPGeomRequest > req)
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:178
math_utils::Point3D< float > & getPos()
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
ServiceRegistryRef services()
Definition InitContext.h:34
ConfigParamRegistry const & options()
Definition InitContext.h:33
decltype(auto) get(R binding, int part=0) const
InputRecord & inputs()
The inputs associated with this processing context.
static 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
bool accountCCDBInputs(const o2::framework::ConcreteDataMatcher &matcher, void *obj)
void process(o2::globaltracking::RecoContainer &recoData)
TrackingStudySpec(std::shared_ptr< DataRequest > dr, std::shared_ptr< o2::base::GRPGeomRequest > gr, GTrackID::mask_t src, bool useMC)
void run(ProcessingContext &pc) final
void finaliseCCDB(ConcreteDataMatcher &matcher, void *obj) final
void endOfStream(EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
void init(InitContext &ic) final
GLenum src
Definition glcorearb.h:1767
const GLdouble * v
Definition glcorearb.h:832
GLenum GLenum dst
Definition glcorearb.h:1767
GLuint id
Definition glcorearb.h:650
constexpr o2::header::DataOrigin gDataOriginTPC
Definition DataHeader.h:576
constexpr double LHCBunchSpacingMUS
constexpr int LHCMaxBunches
constexpr double LHCBunchSpacingNS
constexpr float TwoPI
constexpr float PI
Defining ITS Vertex explicitly as messageable.
Definition Cartesian.h:288
std::vector< ConfigParamSpec > ccdbParamSpec(std::string const &path, int runDependent, std::vector< CCDBMetadata > metadata={}, int qrate=0)
std::vector< ConfigParamSpec > Options
detail::Bracket< float > Bracketf_t
Definition Primitive2D.h:40
TrackParCovF TrackParCov
Definition Track.h:33
o2::math_utils::Bracketf_t TBracket
o2::dataformats::VtxTrackRef V2TRef
o2::dataformats::VtxTrackIndex VTIndex
o2::framework::DataProcessorSpec getTrackingStudySpec(o2::dataformats::GlobalTrackID::mask_t srcTracks, o2::dataformats::GlobalTrackID::mask_t srcClus, bool useMC)
create a processor spec
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
uint32_t orbit
LHC orbit.
bool isSelected(const RecPoints &rp) 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
const o2::dataformats::TrackTPCITS & getTPCITSTrack(GTrackID gid) const
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
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
const o2::dataformats::MatchInfoTOF & getTOFMatch(GTrackID id) const
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
const std::string str