Project
Loading...
Searching...
No Matches
TPCRefitter.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 <random>
13#include <vector>
14#include <TStopwatch.h>
25#include "Framework/Task.h"
26#include "MathUtils/Tsallis.h"
33#include "GPUO2InterfaceRefit.h"
38
39namespace o2::trackstudy
40{
41
42using namespace o2::framework;
45
51
53
54class TPCRefitterSpec final : public Task
55{
56 public:
57 enum StudyType {
58 TPC = 0x1,
59 ITSTPC = 0x2,
60 Cosmics = 0x4,
61 };
63 Streamer = 0x1,
64 TFVectors = 0x2,
65 };
66 TPCRefitterSpec(std::shared_ptr<DataRequest> dr, std::shared_ptr<o2::base::GRPGeomRequest> gr, const o2::tpc::CorrectionMapsLoaderGloOpts& sclOpts, GTrackID::mask_t src, bool useMC)
67 : mDataRequest(dr), mGGCCDBRequest(gr), mTracksSrc(src), mUseMC(useMC)
68 {
69 mTPCCorrMapsLoader.setLumiScaleType(sclOpts.lumiType);
70 mTPCCorrMapsLoader.setLumiScaleMode(sclOpts.lumiMode);
71 mTPCCorrMapsLoader.setCheckCTPIDCConsistency(sclOpts.checkCTPIDCconsistency);
72 }
73 ~TPCRefitterSpec() final = default;
74 void init(InitContext& ic) final;
75 void run(ProcessingContext& pc) final;
76 void endOfStream(EndOfStreamContext& ec) final;
77 void finaliseCCDB(ConcreteDataMatcher& matcher, void* obj) final;
78 void process(o2::globaltracking::RecoContainer& recoData);
79 bool getDCAs(const o2::track::TrackPar& track, float& dcar, float& dcaz);
80
81 private:
82 void updateTimeDependentParams(ProcessingContext& pc);
83 std::shared_ptr<DataRequest> mDataRequest;
84 std::shared_ptr<o2::base::GRPGeomRequest> mGGCCDBRequest;
85 o2::tpc::VDriftHelper mTPCVDriftHelper{};
86 o2::tpc::CorrectionMapsLoader mTPCCorrMapsLoader{};
87 bool mUseMC{false};
88 bool mUseGPUModel{false};
89 float mXRef = 83.;
90 float mDCAMinPt = 1.;
91 int mTFStart = 0;
92 int mTFEnd = 999999999;
93 int mTFCount = -1;
94 int mDCAMinNCl = 80;
95 int mStudyType = 0;
96 int mWriterType = 0;
97 float mSqrt{13600};
98 float mSamplingFactor{0.1};
99 bool mUseR = false;
100 bool mEnableDCA = false;
101 int mWriteTrackClusters = 0;
102 bool mDoSampling{false};
103 bool mDoRefit{true};
104 std::vector<size_t> mClusterOccupancy;
105 std::vector<size_t> mITSTPCTrackOccupanyTPCTime;
106 std::vector<size_t> mITSTPCTrackOccupanyCombinedTime;
107 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOutTPC;
108 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOutITSTPC;
109 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOutTPCTF;
110 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOutITSTPCTF;
111 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOutCosmics;
112 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOutCl;
113 float mITSROFrameLengthMUS = 0.;
114 GTrackID::mask_t mTracksSrc{};
116 std::mt19937 mGenerator;
117 float mVdriftTB;
118 float mTPCTBBias;
119 uint32_t mTimeBinsPerTF{};
120 uint32_t mOccupancyBinsPerTF{};
121 uint32_t mTimeBinsPerDrift{500};
122 //
123 // Input data
124 //
125 gsl::span<const o2::tpc::TPCClRefElem> mTPCTrackClusIdx;
126 gsl::span<const o2::tpc::TrackTPC> mTPCTracksArray;
127 gsl::span<const o2::dataformats::TrackTPCITS> mITSTPCTracksArray;
128 gsl::span<const o2::its::TrackITS> mITSTracksArray;
129 gsl::span<const o2::dataformats::TrackCosmics> mCosmics;
130 gsl::span<const unsigned char> mTPCRefitterShMap;
131 gsl::span<const unsigned int> mTPCRefitterOccMap;
132 const o2::tpc::ClusterNativeAccess* mTPCClusterIdxStruct = nullptr;
133 gsl::span<const o2::MCCompLabel> mTPCTrkLabels;
134 std::unique_ptr<o2::gpu::GPUO2InterfaceRefit> mTPCRefitter;
135 std::vector<o2::InteractionTimeRecord> mIntRecs;
136
137 void fillOccupancyVectors(o2::globaltracking::RecoContainer& recoData);
138 bool processTPCTrack(o2::tpc::TrackTPC tr, o2::MCCompLabel lbl, o2::utils::TreeStreamRedirector* streamer, const o2::its::TrackITS* its = nullptr, const o2::dataformats::TrackTPCITS* itstpc = nullptr, bool outward = false, float time0custom = -1);
139 void processCosmics(o2::globaltracking::RecoContainer& recoData);
140};
141
143{
145 mXRef = ic.options().get<float>("target-x");
146 mUseR = ic.options().get<bool>("use-r-as-x");
147 mEnableDCA = ic.options().get<bool>("enable-dcas");
148 mUseGPUModel = ic.options().get<bool>("use-gpu-fitter");
149 mTFStart = ic.options().get<int>("tf-start");
150 mTFEnd = ic.options().get<int>("tf-end");
151 mDCAMinPt = ic.options().get<float>("dcaMinPt");
152 mDCAMinNCl = ic.options().get<float>("dcaMinNCl");
153 mSqrt = ic.options().get<float>("sqrts");
154 mSamplingFactor = ic.options().get<float>("sampling-factor");
155 mDoSampling = ic.options().get<bool>("do-sampling");
156 mDoRefit = ic.options().get<bool>("do-refit");
157 mStudyType = ic.options().get<int>("study-type");
158 mWriterType = ic.options().get<int>("writer-type");
159 mWriteTrackClusters = ic.options().get<int>("write-track-clusters");
160 const auto occBinsPerDrift = ic.options().get<uint32_t>("occupancy-bins-per-drift");
161 mTimeBinsPerTF = (o2::raw::HBFUtils::Instance().nHBFPerTF * o2::constants::lhc::LHCMaxBunches) / 8 + 2 * mTimeBinsPerDrift; // add one drift before and after the TF
162 mOccupancyBinsPerTF = static_cast<uint32_t>(std::ceil(float(mTimeBinsPerTF * occBinsPerDrift) / mTimeBinsPerDrift));
163 mClusterOccupancy.resize(mOccupancyBinsPerTF);
164 mITSTPCTrackOccupanyTPCTime.resize(mOccupancyBinsPerTF);
165 mITSTPCTrackOccupanyCombinedTime.resize(mOccupancyBinsPerTF);
166 LOGP(info, "Using {} bins for the occupancy per TF", mOccupancyBinsPerTF);
167
168 if ((mWriterType & WriterType::Streamer) == WriterType::Streamer) {
169 if ((mStudyType & StudyType::TPC) == StudyType::TPC) {
170 mDBGOutTPC = std::make_unique<o2::utils::TreeStreamRedirector>("tpctracks-study-streamer.root", "recreate");
171 }
172 if ((mStudyType & StudyType::ITSTPC) == StudyType::ITSTPC) {
173 mDBGOutITSTPC = std::make_unique<o2::utils::TreeStreamRedirector>("itstpctracks-study-streamer.root", "recreate");
174 }
175 if ((mStudyType & StudyType::Cosmics) == StudyType::Cosmics) {
176 mDBGOutCosmics = std::make_unique<o2::utils::TreeStreamRedirector>("cosmics-study-streamer.root", "recreate");
177 }
178 }
179 if (ic.options().get<bool>("dump-clusters")) {
180 mDBGOutCl = std::make_unique<o2::utils::TreeStreamRedirector>("tpc-trackStudy-cl.root", "recreate");
181 }
182
183 if (mXRef < 0.) {
184 mXRef = 0.;
185 }
186 mGenerator = std::mt19937(std::random_device{}());
187 mTPCCorrMapsLoader.init(ic);
188}
189
191{
192 ++mTFCount;
193 if (mTFCount < mTFStart || mTFCount > mTFEnd) {
194 LOGP(info, "Skipping TF {}", mTFCount);
195 return;
196 }
197
199 recoData.collectData(pc, *mDataRequest.get()); // select tracks of needed type, with minimal cuts, the real selected will be done in the vertexer
200 updateTimeDependentParams(pc); // Make sure this is called after recoData.collectData, which may load some conditions
201 fillOccupancyVectors(recoData);
202 process(recoData);
203
204 if (mTFCount >= mTFEnd) {
205 LOGP(info, "Stopping processing after TF {}", mTFCount);
207 return;
208 }
209}
210
211void TPCRefitterSpec::updateTimeDependentParams(ProcessingContext& pc)
212{
214 mTPCVDriftHelper.extractCCDBInputs(pc);
215 mTPCCorrMapsLoader.extractCCDBInputs(pc);
216 static bool initOnceDone = false;
217 if (!initOnceDone) { // this params need to be queried only once
218 initOnceDone = true;
219 // none at the moment
220 }
221 // we may have other params which need to be queried regularly
222 bool updateMaps = false;
223 if (mTPCCorrMapsLoader.isUpdated()) {
224 mTPCCorrMapsLoader.acknowledgeUpdate();
225 updateMaps = true;
226 }
227 if (mTPCVDriftHelper.isUpdated()) {
228 LOGP(info, "Updating TPC fast transform map with new VDrift factor of {} wrt reference {} and DriftTimeOffset correction {} wrt {} from source {}",
229 mTPCVDriftHelper.getVDriftObject().corrFact, mTPCVDriftHelper.getVDriftObject().refVDrift,
230 mTPCVDriftHelper.getVDriftObject().timeOffsetCorr, mTPCVDriftHelper.getVDriftObject().refTimeOffset,
231 mTPCVDriftHelper.getSourceName());
232 mTPCVDriftHelper.acknowledgeUpdate();
233 updateMaps = true;
234 }
235 if (updateMaps) {
236 mTPCCorrMapsLoader.updateVDrift(mTPCVDriftHelper.getVDriftObject().corrFact, mTPCVDriftHelper.getVDriftObject().refVDrift, mTPCVDriftHelper.getVDriftObject().getTimeOffset());
237 }
238}
239
240void TPCRefitterSpec::fillOccupancyVectors(o2::globaltracking::RecoContainer& recoData)
241{
242 // reset counters
243 std::fill(mClusterOccupancy.begin(), mClusterOccupancy.end(), 0u);
244 std::fill(mITSTPCTrackOccupanyTPCTime.begin(), mITSTPCTrackOccupanyTPCTime.end(), 0u);
245 std::fill(mITSTPCTrackOccupanyCombinedTime.begin(), mITSTPCTrackOccupanyCombinedTime.end(), 0u);
246
247 // fill cluster occupancy
248 const auto& clusterIndex = recoData.inputsTPCclusters->clusterIndex;
249 using namespace o2::tpc::constants;
250 for (int sector = 0; sector < MAXSECTOR; ++sector) {
251 for (int padrow = 0; padrow < MAXGLOBALPADROW; ++padrow) {
252 for (size_t icl = 0; icl < clusterIndex.nClusters[sector][padrow]; ++icl) {
253 const auto& cl = clusterIndex.clusters[sector][padrow][icl];
254 // shift by one TPC drift to allow seeing pile-up
255 const auto tpcTime = cl.getTime() + mTimeBinsPerDrift;
256 const uint32_t clOccPos = static_cast<uint32_t>(tpcTime * mOccupancyBinsPerTF / mTimeBinsPerTF);
257 if (clOccPos >= mOccupancyBinsPerTF) {
258 LOGP(error, "cluster with time {} outside TPC acceptanc", cl.getTime());
259 } else {
260 ++mClusterOccupancy[clOccPos];
261 }
262 }
263 }
264 }
265
266 // fill track occupancy for its-tpc matched tracks
267 auto tpcTracks = recoData.getTPCTracks();
268 auto itstpcTracks = recoData.getTPCITSTracks();
269 const auto& paramEle = o2::tpc::ParameterElectronics::Instance();
270
271 for (const auto& tpcitsTrack : itstpcTracks) {
272 const auto idxTPC = tpcitsTrack.getRefTPC().getIndex();
273 if (idxTPC >= tpcTracks.size()) {
274 LOGP(fatal, "TPC index {} out of array size {}", idxTPC, tpcTracks.size());
275 }
276 const auto& tpcTrack = tpcTracks[idxTPC];
277 // shift by one TPC drift to allow seeing pile-up
278 const auto tpcTime = tpcTrack.getTime0() + mTimeBinsPerDrift;
279 if (tpcTime >= 0) {
280 const uint32_t clOccPosTPC = static_cast<uint32_t>(tpcTime * mOccupancyBinsPerTF / mTimeBinsPerTF);
281 if (clOccPosTPC < mITSTPCTrackOccupanyTPCTime.size()) {
282 ++mITSTPCTrackOccupanyTPCTime[clOccPosTPC];
283 } else {
284 LOGP(warn, "TF {}: TPC occupancy index {} out of range {}", mTFCount, clOccPosTPC, mITSTPCTrackOccupanyTPCTime.size());
285 }
286 }
287 // convert mus to time bins
288 // shift by one TPC drift to allow seeing pile-up
289 const auto tpcitsTime = tpcitsTrack.getTimeMUS().getTimeStamp() / paramEle.ZbinWidth + mTimeBinsPerDrift;
290 if (tpcitsTime > 0) {
291 const uint32_t clOccPosITSTPC = static_cast<uint32_t>(tpcitsTime * mOccupancyBinsPerTF / mTimeBinsPerTF);
292 if (clOccPosITSTPC < mITSTPCTrackOccupanyCombinedTime.size()) {
293 ++mITSTPCTrackOccupanyCombinedTime[clOccPosITSTPC];
294 }
295 }
296 }
297
298 auto fillDebug = [this](o2::utils::TreeStreamRedirector* streamer) {
299 if (streamer) {
300 *streamer << "occupancy"
301 << "tfCounter=" << mTFCount
302 << "clusterOcc=" << mClusterOccupancy
303 << "tpcTrackTimeOcc=" << mITSTPCTrackOccupanyTPCTime
304 << "itstpcTrackTimeOcc=" << mITSTPCTrackOccupanyCombinedTime
305 << "\n";
306 }
307 };
308
309 fillDebug(mDBGOutTPC.get());
310 fillDebug(mDBGOutITSTPC.get());
311}
312
314{
315 auto prop = o2::base::Propagator::Instance();
316
317 mITSTracksArray = recoData.getITSTracks();
318 mTPCTracksArray = recoData.getTPCTracks();
319 mITSTPCTracksArray = recoData.getTPCITSTracks();
320 mCosmics = recoData.getCosmicTracks();
321
322 mTPCTrackClusIdx = recoData.getTPCTracksClusterRefs();
323 mTPCClusterIdxStruct = &recoData.inputsTPCclusters->clusterIndex;
324 mTPCRefitterShMap = recoData.clusterShMapTPC;
325 mTPCRefitterOccMap = recoData.occupancyMapTPC;
326
327 LOGP(info, "Processing TF {} with {} its, {} tpc, {} its-tpc tracks and {} comsmics", mTFCount, mITSTracksArray.size(), mTPCTracksArray.size(), mITSTPCTracksArray.size(), mCosmics.size());
328 if (mUseMC) { // extract MC tracks
329 const o2::steer::DigitizationContext* digCont = nullptr;
330 if (!mcReader.initFromDigitContext("collisioncontext.root")) {
331 throw std::invalid_argument("initialization of MCKinematicsReader failed");
332 }
333 digCont = mcReader.getDigitizationContext();
334 mIntRecs = digCont->getEventRecords();
335 mTPCTrkLabels = recoData.getTPCTracksMCLabels();
336 }
337
338 mTPCRefitter = std::make_unique<o2::gpu::GPUO2InterfaceRefit>(mTPCClusterIdxStruct, &mTPCCorrMapsLoader, prop->getNominalBz(), mTPCTrackClusIdx.data(), 0, mTPCRefitterShMap.data(), mTPCRefitterOccMap.data(), mTPCRefitterOccMap.size(), nullptr, prop);
339 mTPCRefitter->setTrackReferenceX(900); // disable propagation after refit by setting reference to value > 500
340
341 mVdriftTB = mTPCVDriftHelper.getVDriftObject().getVDrift() * o2::tpc::ParameterElectronics::Instance().ZbinWidth; // VDrift expressed in cm/TimeBin
342 mTPCTBBias = mTPCVDriftHelper.getVDriftObject().getTimeOffset() / (8 * o2::constants::lhc::LHCBunchSpacingMUS);
343
344 auto dumpClusters = [this] {
345 static int tf = 0;
346 const auto* corrMap = this->mTPCCorrMapsLoader.getCorrMap();
347 for (int sector = 0; sector < 36; sector++) {
348 float alp = ((sector % 18) * 20 + 10) * TMath::DegToRad();
349 float sn = TMath::Sin(alp), cs = TMath::Cos(alp);
350 for (int row = 0; row < 152; row++) {
351 for (int ic = 0; ic < this->mTPCClusterIdxStruct->nClusters[sector][row]; ic++) {
352 const auto cl = this->mTPCClusterIdxStruct->clusters[sector][row][ic];
353 float x, y, z, xG, yG;
354 corrMap->TransformIdeal(sector, row, cl.getPad(), cl.getTime(), x, y, z, 0);
355 o2::math_utils::detail::rotateZ(x, y, xG, yG, sn, cs);
356 LOGP(debug, "tf:{} s:{} r:{} p:{} t:{} qm:{} qt:{} f:{} x:{} y:{} z:{}", tf, sector, row, cl.getPad(), cl.getTime(), cl.getQmax(), cl.getQtot(), cl.getFlags(), x, y, z);
357 (*mDBGOutCl) << "tpccl"
358 << "tf=" << tf << "sect=" << sector << "row=" << row << "pad=" << cl.getPad() << "time=" << cl.getTime() << "qmax=" << cl.getQmax() << "qtot=" << cl.getQtot()
359 << "sigT=" << cl.getSigmaTime() << "sigP=" << cl.getSigmaPad()
360 << "flags=" << cl.getFlags()
361 << "x=" << x << "y=" << y << "z=" << z << "xg=" << xG << "yg=" << yG
362 << "\n";
363 }
364 }
365 }
366 tf++;
367 };
368
369 if (mDBGOutCl) {
370 dumpClusters();
371 }
372
373 if ((mStudyType & StudyType::TPC) == StudyType::TPC) {
374 for (size_t itr = 0; itr < mTPCTracksArray.size(); itr++) {
375 processTPCTrack(mTPCTracksArray[itr], mUseMC ? mTPCTrkLabels[itr] : o2::MCCompLabel{}, mDBGOutTPC.get());
376 }
377 }
378
379 if ((mStudyType & StudyType::ITSTPC) == StudyType::ITSTPC) {
380 for (const auto& tpcitsTrack : mITSTPCTracksArray) {
381 const auto idxTPC = tpcitsTrack.getRefTPC().getIndex();
382 const auto idxITS = tpcitsTrack.getRefITS().getIndex();
383 if (idxITS >= mITSTracksArray.size()) {
384 LOGP(fatal, "ITS index {} out of array size {}", idxITS, mITSTracksArray.size());
385 }
386 if (idxTPC >= mTPCTracksArray.size()) {
387 LOGP(fatal, "TPC index {} out of array size {}", idxTPC, mTPCTracksArray.size());
388 }
389 processTPCTrack(mTPCTracksArray[idxTPC], mUseMC ? mTPCTrkLabels[idxTPC] : o2::MCCompLabel{}, mDBGOutITSTPC.get(), &mITSTracksArray[idxITS], &tpcitsTrack);
390 }
391 }
392
393 if (mCosmics.size() > 0) {
394 LOGP(info, "Procssing {} cosmics", mCosmics.size());
395 processCosmics(recoData);
396 }
397}
398
400{
401 mDBGOutTPC.reset();
402 mDBGOutITSTPC.reset();
403 mDBGOutTPCTF.reset();
404 mDBGOutITSTPCTF.reset();
405 mDBGOutCosmics.reset();
406 mDBGOutCl.reset();
407}
408
410{
412 return;
413 }
414 if (mTPCVDriftHelper.accountCCDBInputs(matcher, obj)) {
415 return;
416 }
417 if (mTPCCorrMapsLoader.accountCCDBInputs(matcher, obj)) {
418 return;
419 }
420}
421
422bool TPCRefitterSpec::getDCAs(const o2::track::TrackPar& track, float& dcar, float& dcaz)
423{
424 auto propagator = o2::base::Propagator::Instance();
425 std::array<float, 2> dca;
426 const o2::math_utils::Point3D<float> refPoint{0, 0, 0};
427 o2::track::TrackPar propTrack(track);
428 const auto ok = propagator->propagateToDCABxByBz(refPoint, propTrack, 2., o2::base::Propagator::MatCorrType::USEMatCorrLUT, &dca);
429 dcar = dca[0];
430 dcaz = dca[1];
431 if (!ok) {
432 dcar = 9998.;
433 dcaz = 9998.;
434 }
435 return ok;
436}
437
438bool TPCRefitterSpec::processTPCTrack(o2::tpc::TrackTPC tr, o2::MCCompLabel lbl, o2::utils::TreeStreamRedirector* streamer, const o2::its::TrackITS* its, const o2::dataformats::TrackTPCITS* itstpc, bool outward, float time0custom)
439{
440 auto prop = o2::base::Propagator::Instance();
441 static long counter = -1;
442
443 struct ClusterData {
444 std::vector<int> occCl;
445 std::vector<short> clSector, clRow;
446 std::vector<float> clX, clY, clZ, clXI, clYI, clZI; // *I are the uncorrected cluster positions
447 std::vector<tpc::ClusterNative> clNative;
448 } clData;
449 float dcar, dcaz, dcarRef, dcazRef;
450
451 // auto tr = mTPCTracksArray[itr]; // create track copy
452 if (tr.hasBothSidesClusters()) {
453 return false;
454 }
455
456 bool sampleTsallis = false;
457 bool sampleMB = false;
458 float tsallisWeight = 0;
459 if (mDoSampling) {
460 std::uniform_real_distribution<> distr(0., 1.);
461 if (o2::math_utils::Tsallis::downsampleTsallisCharged(tr.getPt(), mSamplingFactor, mSqrt, tsallisWeight, distr(mGenerator))) {
462 sampleTsallis = true;
463 }
464 if (distr(mGenerator) < mSamplingFactor) {
465 sampleMB = true;
466 }
467
468 if (!sampleMB && !sampleTsallis) {
469 return false;
470 }
471 }
472 //=========================================================================
473 // create refitted copy
474 auto trackRefit = [&tr, this](o2::track::TrackParCov& trc, float t, float chi2refit, bool outward = false) -> bool {
475 int retVal = mUseGPUModel ? this->mTPCRefitter->RefitTrackAsGPU(trc, tr.getClusterRef(), t, &chi2refit, outward, true)
476 : this->mTPCRefitter->RefitTrackAsTrackParCov(trc, tr.getClusterRef(), t, &chi2refit, outward, true);
477 if (retVal < 0) {
478 LOGP(warn, "Refit failed ({}) with time={}: track#{}[{}]", retVal, t, counter, trc.asString());
479 return false;
480 }
481 return true;
482 };
483
484 auto trackProp = [&tr, prop, this](o2::track::TrackParCov& trc) -> bool {
485 if (!trc.rotate(tr.getAlpha())) {
486 LOGP(warn, "Rotation to original track alpha {} failed, track#{}[{}]", tr.getAlpha(), counter, trc.asString());
487 return false;
488 }
489 float xtgt = this->mXRef;
490 if (mUseR && !trc.getXatLabR(this->mXRef, xtgt, prop->getNominalBz(), o2::track::DirInward)) {
491 xtgt = 0;
492 return false;
493 }
494 if (!prop->PropagateToXBxByBz(trc, xtgt)) {
495 LOGP(warn, "Propagation to X={} failed, track#{}[{}]", xtgt, counter, trc.asString());
496 return false;
497 }
498 return true;
499 };
500
501 // auto prepClus = [this, &tr, &clSector, &clRow, &clX, &clY, &clZ, &clXI, &clYI, &clZI, &clNative](float t) { // extract cluster info
502 auto prepClus = [this, &tr, &clData](float t) { // extract cluster info
503 int count = tr.getNClusters();
504 const auto* corrMap = this->mTPCCorrMapsLoader.getCorrMap();
505 const o2::tpc::ClusterNative* cl = nullptr;
506 for (int ic = count; ic--;) {
507 uint8_t sector, row;
508 uint32_t clusterIndex;
509 o2::tpc::TrackTPC::getClusterReference(mTPCTrackClusIdx, ic, sector, row, clusterIndex, tr.getClusterRef());
510 unsigned int absoluteIndex = mTPCClusterIdxStruct->clusterOffset[sector][row] + clusterIndex;
511 cl = &mTPCClusterIdxStruct->clusters[sector][row][clusterIndex];
512 uint8_t clflags = cl->getFlags();
513 if (mTPCRefitterShMap[absoluteIndex] & o2::gpu::GPUTPCGMMergedTrackHit::flagShared) {
514 clflags |= 0x10;
515 }
516 clData.clSector.emplace_back(sector);
517 clData.clRow.emplace_back(row);
518 auto& clCopy = clData.clNative.emplace_back(*cl);
519 clCopy.setFlags(clflags);
520
521 float x, y, z;
522 // ideal transformation without distortions
523 corrMap->TransformIdeal(sector, row, cl->getPad(), cl->getTime(), x, y, z, t); // nominal time of the track
524 clData.clXI.emplace_back(x);
525 clData.clYI.emplace_back(y);
526 clData.clZI.emplace_back(z);
527
528 // transformation without distortions
529 mTPCCorrMapsLoader.Transform(sector, row, cl->getPad(), cl->getTime(), x, y, z, t); // nominal time of the track
530 clData.clX.emplace_back(x);
531 clData.clY.emplace_back(y);
532 clData.clZ.emplace_back(z);
533
534 // occupancy estimator
535 const auto tpcTime = cl->getTime() + mTimeBinsPerDrift;
536 const uint32_t clOccPosTPC = static_cast<uint32_t>(tpcTime * mOccupancyBinsPerTF / mTimeBinsPerTF);
537 clData.occCl.emplace_back((clOccPosTPC < mClusterOccupancy.size()) ? mClusterOccupancy[clOccPosTPC] : -1);
538 }
539 };
540
541 //=========================================================================
542
543 auto trf = tr.getOuterParam(); // we refit inward original track
544 float chi2refit = 0;
545 float time0 = tr.getTime0();
546 if (time0custom > 0) {
547 time0 = time0custom;
548 }
549 if (mDoRefit) {
550 if (!trackRefit(trf, time0, chi2refit) || !trackProp(trf)) {
551 return false;
552 }
553 }
554
555 // propagate original track
556 if (!trackProp(tr)) {
557 return false;
558 }
559
560 if (mWriteTrackClusters) {
561 prepClus(time0); // original clusters
562 }
563
564 if (mEnableDCA) {
565 dcar = dcaz = dcarRef = dcazRef = 9999.f;
566 if ((trf.getPt() > mDCAMinPt) && (tr.getNClusters() > mDCAMinNCl)) {
567 getDCAs(trf, dcarRef, dcazRef);
568 getDCAs(tr, dcar, dcaz);
569 }
570 }
571
572 counter++;
573 // store results
574 if (streamer) {
575 (*streamer) << "tpc"
576 << "counter=" << counter
577 << "tfCounter=" << mTFCount
578 << "tpc=" << tr;
579
580 if (mDoRefit) {
581 (*streamer) << "tpc"
582 << "tpcRF=" << trf
583 << "time0=" << time0
584 << "chi2refit=" << chi2refit;
585 }
586
587 if (mDoSampling) {
588 (*streamer) << "tpc"
589 << "tsallisWeight=" << tsallisWeight
590 << "sampleTsallis=" << sampleTsallis
591 << "sampleMB=" << sampleMB;
592 }
593
594 if (mWriteTrackClusters) {
595 (*streamer) << "tpc"
596 << "clSector=" << clData.clSector
597 << "clRow=" << clData.clRow;
598
599 if ((mWriteTrackClusters & 0x1) == 0x1) {
600 (*streamer) << "tpc"
601 << "cl=" << clData.clNative;
602 }
603
604 if ((mWriteTrackClusters & 0x2) == 0x2) {
605 (*streamer) << "tpc"
606 << "clX=" << clData.clX
607 << "clY=" << clData.clY
608 << "clZ=" << clData.clZ;
609 }
610
611 if ((mWriteTrackClusters & 0x4) == 0x4) {
612 (*streamer) << "tpc"
613 << "clXI=" << clData.clXI // ideal (uncorrected) cluster positions
614 << "clYI=" << clData.clYI // ideal (uncorrected) cluster positions
615 << "clZI=" << clData.clZI; // ideal (uncorrected) cluster positions
616 }
617
618 if ((mWriteTrackClusters & 0x8) == 0x8) {
619 (*streamer) << "tpc"
620 << "clOcc=" << clData.occCl;
621 }
622 }
623
624 if (its) {
625 (*streamer) << "tpc"
626 << "its=" << *its;
627 }
628 if (itstpc) {
629 (*streamer) << "tpc"
630 << "itstpc=" << *itstpc;
631 }
632
633 if (mEnableDCA) {
634 (*streamer) << "tpc"
635 << "dcar=" << dcar
636 << "dcaz=" << dcaz
637 << "dcarRef=" << dcarRef
638 << "dcazRef=" << dcazRef;
639 }
640
641 if (mUseMC) {
642 (*streamer) << "tpc"
643 << "mcLabel=" << lbl;
644 }
645
646 (*streamer) << "tpc"
647 << "\n";
648 }
649
650 float dz = 0;
651
652 if (mUseMC) { // impose MC time in TPC timebin and refit inward after resetted covariance
653 // extract MC truth
654 const o2::MCTrack* mcTrack = nullptr;
655 if (!lbl.isValid() || !(mcTrack = mcReader.getTrack(lbl))) {
656 return false;
657 }
658 long bc = mIntRecs[lbl.getEventID()].toLong(); // bunch crossing of the interaction
659 float bcTB = bc / 8. + mTPCTBBias; // the same in TPC timebins, accounting for the TPC time bias
660 // create MC truth track in O2 format
661 std::array<float, 3> xyz{(float)mcTrack->GetStartVertexCoordinatesX(), (float)mcTrack->GetStartVertexCoordinatesY(), (float)mcTrack->GetStartVertexCoordinatesZ()},
662 pxyz{(float)mcTrack->GetStartVertexMomentumX(), (float)mcTrack->GetStartVertexMomentumY(), (float)mcTrack->GetStartVertexMomentumZ()};
663 TParticlePDG* pPDG = TDatabasePDG::Instance()->GetParticle(mcTrack->GetPdgCode());
664 if (!pPDG) {
665 return false;
666 }
667 o2::track::TrackPar mctrO2(xyz, pxyz, TMath::Nint(pPDG->Charge() / 3), false);
668 //
669 // propagate it to the alpha/X of the reconstructed track
670 if (!mctrO2.rotate(tr.getAlpha()) || !prop->PropagateToXBxByBz(mctrO2, tr.getX())) {
671 return false;
672 }
673 // now create a properly refitted track with correct time and distortions correction
674 {
675 auto trfm = tr.getOuterParam(); // we refit inward
676 // impose MC time in TPC timebin and refit inward after resetted covariance
677 float chi2refit = 0;
678 if (!trackRefit(trfm, bcTB, chi2refit) || !trfm.rotate(tr.getAlpha()) || !prop->PropagateToXBxByBz(trfm, tr.getX())) {
679 LOGP(warn, "Failed to propagate MC-time refitted track#{} [{}] to X/alpha of original track [{}]", counter, trfm.asString(), tr.asString());
680 return false;
681 }
682 // estimate Z shift in case of no-distortions
683 dz = (tr.getTime0() - bcTB) * mVdriftTB;
684 if (tr.hasCSideClustersOnly()) {
685 dz = -dz;
686 }
687 //
688 prepClus(bcTB); // clusters for MC time
689 if (streamer) {
690 (*streamer) << "tpcMC"
691 << "counter=" << counter
692 << "movTrackRef=" << trfm
693 << "mcTrack=" << mctrO2
694 << "imposedTB=" << bcTB
695 << "chi2refit=" << chi2refit
696 << "dz=" << dz
697 << "clX=" << clData.clX
698 << "clY=" << clData.clY
699 << "clZ=" << clData.clZ
700 << "\n";
701 }
702 }
703 return false;
704 }
705
706 return true;
707}
708
709void TPCRefitterSpec::processCosmics(o2::globaltracking::RecoContainer& recoData)
710{
711 auto tof = recoData.getTOFClusters();
713 const auto invBinWidth = 1.f / par.ZbinWidth;
714
715 for (const auto& cosmic : mCosmics) {
716 //
717 const auto& gidtop = cosmic.getRefTop();
718 const auto& gidbot = cosmic.getRefBottom();
719
720 // LOGP(info, "Sources: {} - {}", o2::dataformats::GlobalTrackID::getSourceName(gidtop.getSource()), o2::dataformats::GlobalTrackID::getSourceName(gidbot.getSource()));
721
722 std::array<GTrackID, GTrackID::NSources> contributorsGID[2] = {recoData.getSingleDetectorRefs(cosmic.getRefTop()), recoData.getSingleDetectorRefs(cosmic.getRefBottom())};
723 const auto trackTime = cosmic.getTimeMUS().getTimeStamp() * invBinWidth;
724
725 // check if track has TPC & TOF for top and bottom part
726 // loop over both parts
727 for (const auto& comsmicInfo : contributorsGID) {
728 auto& tpcGlobal = comsmicInfo[GTrackID::TPC];
729 auto& tofGlobal = comsmicInfo[GTrackID::TOF];
730 if (tpcGlobal.isIndexSet() && tofGlobal.isIndexSet()) {
731 const auto itrTPC = tpcGlobal.getIndex();
732 const auto itrTOF = tofGlobal.getIndex();
733 const auto& tofCl = tof[itrTOF];
734 const auto tofTime = tofCl.getTime() * 1e-6 * invBinWidth; // ps -> us -> time bins
735 const auto tofTimeRaw = tofCl.getTimeRaw() * 1e-6 * invBinWidth; // ps -> us -> time bins
736 const auto& trackTPC = mTPCTracksArray[itrTPC];
737 // LOGP(info, "Cosmic time: {}, TOF time: {}, TOF time raw: {}, TPC time: {}", trackTime, tofTime, tofTimeRaw, trackTPC.getTime0());
738 processTPCTrack(trackTPC, mUseMC ? mTPCTrkLabels[itrTPC] : o2::MCCompLabel{}, mDBGOutCosmics.get(), nullptr, nullptr, false, tofTime);
739 }
740 }
741 }
742}
743
744DataProcessorSpec getTPCRefitterSpec(GTrackID::mask_t srcTracks, GTrackID::mask_t srcClusters, bool useMC, const o2::tpc::CorrectionMapsLoaderGloOpts& sclOpts, bool requestCosmics)
745{
746 std::vector<OutputSpec> outputs;
747 Options opts{
748 {"target-x", VariantType::Float, 83.f, {"Try to propagate to this radius"}},
749 {"dump-clusters", VariantType::Bool, false, {"dump all clusters"}},
750 {"write-track-clusters", VariantType::Int, 3, {"Bitmask write clusters associated to the track, full native cluster (0x1), corrected (0x2) and uncorrected (0x4) positions, (0x8) occupancy info"}},
751 {"tf-start", VariantType::Int, 0, {"1st TF to process"}},
752 {"tf-end", VariantType::Int, 999999999, {"last TF to process"}},
753 {"use-gpu-fitter", VariantType::Bool, false, {"use GPU track model for refit instead of TrackParCov"}},
754 {"do-refit", VariantType::Bool, false, {"do refitting of TPC track"}},
755 {"use-r-as-x", VariantType::Bool, false, {"Use radius instead of target sector X"}},
756 {"enable-dcas", VariantType::Bool, false, {"Propagate to DCA and add it to the tree"}},
757 {"dcaMinPt", VariantType::Float, 1.f, {"Min pT of tracks propagated to DCA"}},
758 {"dcaMinNCl", VariantType::Int, 80, {"Min number of clusters for tracks propagated to DCA"}},
759 {"sqrts", VariantType::Float, 13600.f, {"Centre of mass energy used for downsampling"}},
760 {"do-sampling", VariantType::Bool, false, {"Perform sampling, min. bias and on Tsallis function, using 'sampling-factor'"}},
761 {"sampling-factor", VariantType::Float, 0.1f, {"Sampling factor in case sample-unbinned-tsallis is used"}},
762 {"study-type", VariantType::Int, 1, {"Bitmask of study type: 0x1 = TPC only, 0x2 = TPC + ITS, 0x4 = Cosmics"}},
763 {"writer-type", VariantType::Int, 1, {"Bitmask of writer type: 0x1 = per track streamer, 0x2 = per TF vectors"}},
764 {"occupancy-bins-per-drift", VariantType::UInt32, 31u, {"number of bin for occupancy histogram per drift time (500tb)"}},
765 };
766 auto dataRequest = std::make_shared<DataRequest>();
767
768 dataRequest->requestTracks(srcTracks, useMC);
769 dataRequest->requestClusters(srcClusters, useMC);
770 if (requestCosmics) {
771 dataRequest->requestCoscmicTracks(useMC);
772 }
773 auto ggRequest = std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
774 true, // GRPECS=true
775 false, // GRPLHCIF
776 true, // GRPMagField
777 true, // askMatLUT
779 dataRequest->inputs,
780 true);
781 o2::tpc::VDriftHelper::requestCCDBInputs(dataRequest->inputs);
782 o2::tpc::CorrectionMapsLoader::requestCCDBInputs(dataRequest->inputs, opts, sclOpts);
783
784 return DataProcessorSpec{
785 "tpc-refitter",
786 dataRequest->inputs,
787 outputs,
788 AlgorithmSpec{adaptFromTask<TPCRefitterSpec>(dataRequest, ggRequest, sclOpts, srcTracks, useMC)},
789 opts};
790}
791
792} // namespace o2::trackstudy
Helper class to access load maps from CCDB.
Wrapper container for different reconstructed object types.
uint64_t bc
Definition RawEventData.h:5
int32_t retVal
Helper for geometry and GRP related CCDB requests.
Global index for barrel track: provides provenance (detectors combination), index in respective array...
Header to collect LHC related constants.
Definition of the parameter class for the detector electronics.
uint32_t padrow
Definition RawData.h:5
Wrapper container for different reconstructed object types.
Result of top-bottom cosmic tracks leg matching.
Helper class to extract VDrift from different sources.
std::ostringstream debug
int getEventID() const
bool isValid() const
Definition MCCompLabel.h:75
Double_t GetStartVertexMomentumZ() const
Definition MCTrack.h:81
Double_t GetStartVertexMomentumX() const
Definition MCTrack.h:79
Double_t GetStartVertexCoordinatesY() const
Definition MCTrack.h:83
Double_t GetStartVertexCoordinatesZ() const
Definition MCTrack.h:84
Double_t GetStartVertexMomentumY() const
Definition MCTrack.h:80
Double_t GetStartVertexCoordinatesX() const
Definition MCTrack.h:82
Int_t GetPdgCode() const
Accessors.
Definition MCTrack.h:72
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:147
Static class with identifiers, bitmasks and names for ALICE detectors.
Definition DetID.h:58
ConfigParamRegistry const & options()
Definition InitContext.h:33
ServiceRegistryRef services()
The services registry associated with this processing context.
std::vector< o2::InteractionTimeRecord > & getEventRecords(bool withQED=false)
bool initFromDigitContext(std::string_view filename)
DigitizationContext const * getDigitizationContext() const
MCTrack const * getTrack(o2::MCCompLabel const &) const
void extractCCDBInputs(o2::framework::ProcessingContext &pc)
void updateVDrift(float vdriftCorr, float vdrifRef, float driftTimeOffset=0)
bool accountCCDBInputs(const o2::framework::ConcreteDataMatcher &matcher, void *obj)
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
bool accountCCDBInputs(const o2::framework::ConcreteDataMatcher &matcher, void *obj)
static std::string_view getSourceName(Source s)
bool isUpdated() const
TPCRefitterSpec(std::shared_ptr< DataRequest > dr, std::shared_ptr< o2::base::GRPGeomRequest > gr, const o2::tpc::CorrectionMapsLoaderGloOpts &sclOpts, GTrackID::mask_t src, bool useMC)
void process(o2::globaltracking::RecoContainer &recoData)
void run(ProcessingContext &pc) final
void init(InitContext &ic) final
void endOfStream(EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
@ ITSTPC
TPC + ITS matched tracks.
bool getDCAs(const o2::track::TrackPar &track, float &dcar, float &dcaz)
@ Streamer
Write per track streamer information.
@ TFVectors
Writer vectors per TF.
void finaliseCCDB(ConcreteDataMatcher &matcher, void *obj) final
GLint GLenum GLint x
Definition glcorearb.h:403
GLenum src
Definition glcorearb.h:1767
GLint GLsizei count
Definition glcorearb.h:399
GLint y
Definition glcorearb.h:270
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
GLuint counter
Definition glcorearb.h:3987
uint8_t itsSharedClusterMap uint8_t
constexpr double LHCBunchSpacingMUS
constexpr int LHCMaxBunches
Node par(int index)
Parameters.
Defining PrimaryVertex explicitly as messageable.
Definition TFIDInfo.h:20
std::vector< ConfigParamSpec > Options
auto get(const std::byte *buffer, size_t=0)
Definition DataHeader.h:454
std::tuple< T, T > rotateZ(T xL, T yL, T snAlp, T csAlp)
detail::Bracket< float > Bracketf_t
Definition Primitive2D.h:40
return * this
constexpr int MAXSECTOR
Definition Constants.h:28
constexpr int MAXGLOBALPADROW
Definition Constants.h:34
TrackParCovF TrackParCov
Definition Track.h:33
o2::math_utils::Bracketf_t TBracket
o2::dataformats::VtxTrackRef V2TRef
o2::dataformats::VtxTrackIndex VTIndex
o2::framework::DataProcessorSpec getTPCRefitterSpec(o2::dataformats::GlobalTrackID::mask_t srcTracks, o2::dataformats::GlobalTrackID::mask_t srcClus, bool useMC, const o2::tpc::CorrectionMapsLoaderGloOpts &sclOpts, bool requestCosmics=false)
create a processor spec
@ sampleTsallis
perform sampling on tsallis pdf
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
Defining DataPointCompositeObject explicitly as copiable.
std::unique_ptr< GPUReconstructionTimeframe > tf
GlobalIDSet getSingleDetectorRefs(GTrackID gidx) const
gsl::span< const unsigned char > clusterShMapTPC
externally set TPC clusters sharing map
void collectData(o2::framework::ProcessingContext &pc, const DataRequest &request)
std::unique_ptr< o2::tpc::internal::getWorkflowTPCInput_ret > inputsTPCclusters
gsl::span< const unsigned int > occupancyMapTPC
externally set TPC clusters occupancy map
static bool downsampleTsallisCharged(float pt, float factorPt, float sqrts, float &weight, float rnd, float mass=0.13957)
Definition Tsallis.cxx:31
int nHBFPerTF
number of orbits per BC
Definition HBFUtils.h:138
unsigned int nClusters[constants::MAXSECTOR][constants::MAXGLOBALPADROW]
const ClusterNative * clusters[constants::MAXSECTOR][constants::MAXGLOBALPADROW]
unsigned int clusterOffset[constants::MAXSECTOR][constants::MAXGLOBALPADROW]
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)
std::uniform_int_distribution< unsigned long long > distr
std::vector< int > row