Project
Loading...
Searching...
No Matches
TPCTrackStudy.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>
30#include "GPUO2InterfaceRefit.h"
34
35namespace o2::trackstudy
36{
37
38using namespace o2::framework;
41
47
49
50class TPCTrackStudySpec : public Task
51{
52 public:
53 TPCTrackStudySpec(std::shared_ptr<DataRequest> dr, std::shared_ptr<o2::base::GRPGeomRequest> gr, const o2::tpc::CorrectionMapsLoaderGloOpts& sclOpts, GTrackID::mask_t src, bool useMC)
54 : mDataRequest(dr), mGGCCDBRequest(gr), mTracksSrc(src), mUseMC(useMC)
55 {
56 mTPCCorrMapsLoader.setLumiScaleType(sclOpts.lumiType);
57 mTPCCorrMapsLoader.setLumiScaleMode(sclOpts.lumiMode);
58 }
59 ~TPCTrackStudySpec() final = default;
60 void init(InitContext& ic) final;
61 void run(ProcessingContext& pc) final;
62 void endOfStream(EndOfStreamContext& ec) final;
63 void finaliseCCDB(ConcreteDataMatcher& matcher, void* obj) final;
64 void process(o2::globaltracking::RecoContainer& recoData);
65
66 private:
67 void updateTimeDependentParams(ProcessingContext& pc);
68 std::shared_ptr<DataRequest> mDataRequest;
69 std::shared_ptr<o2::base::GRPGeomRequest> mGGCCDBRequest;
70 o2::tpc::VDriftHelper mTPCVDriftHelper{};
71 o2::tpc::CorrectionMapsLoader mTPCCorrMapsLoader{};
72 bool mUseMC{false};
73 bool mUseGPUModel{false};
74 float mXRef = 0.;
75 int mNMoves = 6;
76 int mTFStart = 0;
77 int mTFEnd = 999999999;
78 int mTFCount = -1;
79 bool mUseR = false;
80 bool mRepRef = false;
81 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOut;
82 std::unique_ptr<o2::utils::TreeStreamRedirector> mDBGOutCl;
83 float mITSROFrameLengthMUS = 0.;
84 GTrackID::mask_t mTracksSrc{};
85 o2::steer::MCKinematicsReader mcReader; // reader of MC information
86 //
87 // TPC data
88 gsl::span<const o2::tpc::TPCClRefElem> mTPCTrackClusIdx;
89 gsl::span<const o2::tpc::TrackTPC> mTPCTracksArray;
90 gsl::span<const unsigned char> mTPCRefitterShMap;
91 gsl::span<const unsigned int> mTPCRefitterOccMap;
92 const o2::tpc::ClusterNativeAccess* mTPCClusterIdxStruct = nullptr;
93 gsl::span<const o2::MCCompLabel> mTPCTrkLabels;
94 std::unique_ptr<o2::gpu::GPUO2InterfaceRefit> mTPCRefitter;
95};
96
98{
100 mXRef = ic.options().get<float>("target-x");
101 mNMoves = std::max(2, ic.options().get<int>("n-moves"));
102 mUseR = ic.options().get<bool>("use-r-as-x");
103 mRepRef = ic.options().get<bool>("repeat-ini-ref");
104 mUseGPUModel = ic.options().get<bool>("use-gpu-fitter");
105 mTFStart = ic.options().get<int>("tf-start");
106 mTFEnd = ic.options().get<int>("tf-end");
107 if (mXRef < 0.) {
108 mXRef = 0.;
109 }
110 mTPCCorrMapsLoader.init(ic);
111 mDBGOut = std::make_unique<o2::utils::TreeStreamRedirector>("tpc-trackStudy.root", "recreate");
112 if (ic.options().get<bool>("dump-clusters")) {
113 mDBGOutCl = std::make_unique<o2::utils::TreeStreamRedirector>("tpc-trackStudy-cl.root", "recreate");
114 }
115}
116
118{
119 mTFCount++;
120 if (mTFCount < mTFStart || mTFCount > mTFEnd) {
121 LOGP(info, "Skipping TF {}", mTFCount);
122 return;
123 }
124
126 recoData.collectData(pc, *mDataRequest.get()); // select tracks of needed type, with minimal cuts, the real selected will be done in the vertexer
127 updateTimeDependentParams(pc); // Make sure this is called after recoData.collectData, which may load some conditions
128 process(recoData);
129
130 if (mTFCount > mTFEnd) {
131 LOGP(info, "Stopping processing after TF {}", mTFCount);
133 return;
134 }
135}
136
137void TPCTrackStudySpec::updateTimeDependentParams(ProcessingContext& pc)
138{
140 mTPCVDriftHelper.extractCCDBInputs(pc);
141 mTPCCorrMapsLoader.extractCCDBInputs(pc);
142 static bool initOnceDone = false;
143 if (!initOnceDone) { // this params need to be queried only once
144 initOnceDone = true;
145 // none at the moment
146 }
147 // we may have other params which need to be queried regularly
148 bool updateMaps = false;
149 if (mTPCCorrMapsLoader.isUpdated()) {
150 mTPCCorrMapsLoader.acknowledgeUpdate();
151 updateMaps = true;
152 }
153 if (mTPCVDriftHelper.isUpdated()) {
154 LOGP(info, "Updating TPC fast transform map with new VDrift factor of {} wrt reference {} and DriftTimeOffset correction {} wrt {} from source {}",
155 mTPCVDriftHelper.getVDriftObject().corrFact, mTPCVDriftHelper.getVDriftObject().refVDrift,
156 mTPCVDriftHelper.getVDriftObject().timeOffsetCorr, mTPCVDriftHelper.getVDriftObject().refTimeOffset,
157 mTPCVDriftHelper.getSourceName());
158 mTPCVDriftHelper.acknowledgeUpdate();
159 updateMaps = true;
160 }
161 if (updateMaps) {
162 mTPCCorrMapsLoader.updateVDrift(mTPCVDriftHelper.getVDriftObject().corrFact, mTPCVDriftHelper.getVDriftObject().refVDrift, mTPCVDriftHelper.getVDriftObject().getTimeOffset());
163 }
164}
165
167{
168 static long counter = -1;
169 auto prop = o2::base::Propagator::Instance();
170
171 mTPCTracksArray = recoData.getTPCTracks();
172 mTPCTrackClusIdx = recoData.getTPCTracksClusterRefs();
173 mTPCClusterIdxStruct = &recoData.inputsTPCclusters->clusterIndex;
174 mTPCRefitterShMap = recoData.clusterShMapTPC;
175 mTPCRefitterOccMap = recoData.occupancyMapTPC;
176
177 std::vector<o2::InteractionTimeRecord> intRecs;
178 if (mUseMC) { // extract MC tracks
179 const o2::steer::DigitizationContext* digCont = nullptr;
180 if (!mcReader.initFromDigitContext("collisioncontext.root")) {
181 throw std::invalid_argument("initialization of MCKinematicsReader failed");
182 }
183 digCont = mcReader.getDigitizationContext();
184 intRecs = digCont->getEventRecords();
185 mTPCTrkLabels = recoData.getTPCTracksMCLabels();
186 }
187 if (mTPCTracksArray.size()) {
188 LOGP(info, "Found {} TPC tracks", mTPCTracksArray.size());
189 mTPCRefitter = std::make_unique<o2::gpu::GPUO2InterfaceRefit>(mTPCClusterIdxStruct, &mTPCCorrMapsLoader, prop->getNominalBz(), mTPCTrackClusIdx.data(), 0, mTPCRefitterShMap.data(), mTPCRefitterOccMap.data(), mTPCRefitterOccMap.size(), nullptr, o2::base::Propagator::Instance());
190 mTPCRefitter->setTrackReferenceX(900); // disable propagation after refit by setting reference to value > 500
191 }
192 float vdriftTB = mTPCVDriftHelper.getVDriftObject().getVDrift() * o2::tpc::ParameterElectronics::Instance().ZbinWidth; // VDrift expressed in cm/TimeBin
193 float tpcTBBias = mTPCVDriftHelper.getVDriftObject().getTimeOffset() / (8 * o2::constants::lhc::LHCBunchSpacingMUS);
194 std::vector<short> clSector, clRow;
195 std::vector<float> clX, clY, clZ;
196
197 auto dumpClusters = [this] {
198 static int tf = 0;
199 const auto* corrMap = this->mTPCCorrMapsLoader.getCorrMap();
200 for (int sector = 0; sector < 36; sector++) {
201 float alp = ((sector % 18) * 20 + 10) * TMath::DegToRad();
202 float sn = TMath::Sin(alp), cs = TMath::Cos(alp);
203 for (int row = 0; row < 152; row++) {
204 for (int ic = 0; ic < this->mTPCClusterIdxStruct->nClusters[sector][row]; ic++) {
205 const auto cl = this->mTPCClusterIdxStruct->clusters[sector][row][ic];
206 float x, y, z, xG, yG;
207 corrMap->TransformIdeal(sector, row, cl.getPad(), cl.getTime(), x, y, z, 0);
208 o2::math_utils::detail::rotateZ(x, y, xG, yG, sn, cs);
209 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);
210 (*mDBGOutCl) << "tpccl"
211 << "tf=" << tf << "sect=" << sector << "row=" << row << "pad=" << cl.getPad() << "time=" << cl.getTime() << "qmax=" << cl.getQmax() << "qtot=" << cl.getQtot()
212 << "sigT=" << cl.getSigmaTime() << "sigP=" << cl.getSigmaPad()
213 << "flags=" << cl.getFlags()
214 << "x=" << x << "y=" << y << "z=" << z << "xg=" << xG << "yg=" << yG
215 << "\n";
216 }
217 }
218 }
219 tf++;
220 };
221
222 if (mDBGOutCl) {
223 dumpClusters();
224 }
225
226 for (size_t itr = 0; itr < mTPCTracksArray.size(); itr++) {
227 auto tr = mTPCTracksArray[itr]; // create track copy
228 int side = 0;
229 if (tr.hasBothSidesClusters()) {
230 continue;
231 }
232 side = tr.hasASideClustersOnly() ? 1 : -1;
233 //=========================================================================
234 // create refitted copy
235 auto trackRefit = [itr, this](o2::track::TrackParCov& trc, float t) -> bool {
236 float chi2Out = 0;
237 int retVal = mUseGPUModel ? this->mTPCRefitter->RefitTrackAsGPU(trc, this->mTPCTracksArray[itr].getClusterRef(), t, &chi2Out, false, true) : this->mTPCRefitter->RefitTrackAsTrackParCov(trc, this->mTPCTracksArray[itr].getClusterRef(), t, &chi2Out, false, true);
238 if (retVal < 0) {
239 LOGP(warn, "Refit failed ({}) with time={}: track#{}[{}]", retVal, t, counter, trc.asString());
240 return false;
241 }
242 return true;
243 };
244
245 auto trackProp = [&tr, itr, prop, this](o2::track::TrackParCov& trc) -> bool {
246 if (!trc.rotate(tr.getAlpha())) {
247 LOGP(warn, "Rotation to original track alpha {} failed, track#{}[{}]", tr.getAlpha(), counter, trc.asString());
248 return false;
249 }
250 float xtgt = this->mXRef;
251 if (mUseR && !trc.getXatLabR(this->mXRef, xtgt, prop->getNominalBz(), o2::track::DirInward)) {
252 xtgt = 0;
253 return false;
254 }
255 if (!prop->PropagateToXBxByBz(trc, xtgt)) {
256 LOGP(warn, "Propagation to X={} failed, track#{}[{}]", xtgt, counter, trc.asString());
257 return false;
258 }
259 return true;
260 };
261
262 auto prepClus = [this, &tr, &clSector, &clRow, &clX, &clY, &clZ](float t) { // extract cluster info
263 clSector.clear();
264 clRow.clear();
265 clX.clear();
266 clY.clear();
267 clZ.clear();
268 int count = tr.getNClusters();
269 const o2::tpc::ClusterNative* cl = nullptr;
270 for (int ic = count; ic--;) {
271 uint8_t sector, row;
272 cl = &tr.getCluster(this->mTPCTrackClusIdx, ic, *this->mTPCClusterIdxStruct, sector, row);
273 clSector.push_back(sector);
274 clRow.push_back(row);
275 float x, y, z;
276 mTPCCorrMapsLoader.Transform(sector, row, cl->getPad(), cl->getTime(), x, y, z, t); // nominal time of the track
277 clX.push_back(x);
278 clY.push_back(y);
279 clZ.push_back(z);
280 }
281 };
282
283 //=========================================================================
284
285 auto trf = tr.getOuterParam(); // we refit inward original track
286 if (!trackRefit(trf, tr.getTime0()) || !trackProp(trf)) {
287 continue;
288 }
289
290 // propagate original track
291 if (!trackProp(tr)) {
292 continue;
293 }
294
295 prepClus(tr.getTime0()); // original clusters
296 counter++;
297 // store results
298 (*mDBGOut) << "tpcIni"
299 << "counter=" << counter
300 << "iniTrack=" << tr
301 << "iniTrackRef=" << trf
302 << "side=" << side
303 << "time=" << tr.getTime0()
304 << "clSector=" << clSector
305 << "clRow=" << clRow
306 << "clX=" << clX
307 << "clY=" << clY
308 << "clZ=" << clZ
309 << "\n";
310
311 float dz = 0;
312
313 while (mUseMC) { // impose MC time in TPC timebin and refit inward after resetted covariance
314 // extract MC truth
315 const o2::MCTrack* mcTrack = nullptr;
316 auto lbl = mTPCTrkLabels[itr];
317 if (!lbl.isValid() || !(mcTrack = mcReader.getTrack(lbl))) {
318 break;
319 }
320 long bc = intRecs[lbl.getEventID()].toLong(); // bunch crossing of the interaction
321 float bcTB = bc / 8. + tpcTBBias; // the same in TPC timebins, accounting for the TPC time bias
322 // create MC truth track in O2 format
323 std::array<float, 3> xyz{(float)mcTrack->GetStartVertexCoordinatesX(), (float)mcTrack->GetStartVertexCoordinatesY(), (float)mcTrack->GetStartVertexCoordinatesZ()},
324 pxyz{(float)mcTrack->GetStartVertexMomentumX(), (float)mcTrack->GetStartVertexMomentumY(), (float)mcTrack->GetStartVertexMomentumZ()};
325 TParticlePDG* pPDG = TDatabasePDG::Instance()->GetParticle(mcTrack->GetPdgCode());
326 if (!pPDG) {
327 break;
328 }
329 o2::track::TrackPar mctrO2(xyz, pxyz, TMath::Nint(pPDG->Charge() / 3), false);
330 //
331 // propagate it to the alpha/X of the reconstructed track
332 if (!mctrO2.rotate(tr.getAlpha()) || !prop->PropagateToXBxByBz(mctrO2, tr.getX())) {
333 break;
334 }
335 // now create a properly refitted track with correct time and distortions correction
336 {
337 auto trfm = tr.getOuterParam(); // we refit inward
338 // impose MC time in TPC timebin and refit inward after resetted covariance
339 if (!trackRefit(trfm, bcTB) || !trfm.rotate(tr.getAlpha()) || !prop->PropagateToXBxByBz(trfm, tr.getX())) {
340 LOGP(warn, "Failed to propagate MC-time refitted track#{} [{}] to X/alpha of original track [{}]", counter, trfm.asString(), tr.asString());
341 break;
342 }
343 // estimate Z shift in case of no-distortions
344 dz = (tr.getTime0() - bcTB) * vdriftTB;
345 if (tr.hasCSideClustersOnly()) {
346 dz = -dz;
347 }
348 //
349 prepClus(bcTB); // clusters for MC time
350 (*mDBGOut) << "tpcMC"
351 << "counter=" << counter
352 << "movTrackRef=" << trfm
353 << "mcTrack=" << mctrO2
354 << "side=" << side
355 << "imposedTB=" << bcTB
356 << "dz=" << dz
357 << "clX=" << clX
358 << "clY=" << clY
359 << "clZ=" << clZ
360 << "\n";
361 }
362 break;
363 }
364 // refit and store the same track for a few compatible times
365 float tmin = tr.getTime0() - tr.getDeltaTBwd();
366 float tmax = tr.getTime0() + tr.getDeltaTFwd();
367 for (int it = 0; it < mNMoves; it++) {
368 float tb = tmin + it * (tmax - tmin) / (mNMoves - 1);
369 auto trfm = tr.getOuterParam(); // we refit inward
370 // impose time in TPC timebin and refit inward after resetted covariance
371 if (!trackRefit(trfm, tb) || !trfm.rotate(tr.getAlpha()) || !prop->PropagateToXBxByBz(trfm, tr.getX())) {
372 LOGP(warn, "Failed to propagate time={} refitted track#{} [{}] to X/alpha of original track [{}]", tb, counter, trfm.asString(), tr.asString());
373 continue;
374 }
375 // estimate Z shift in case of no-distortions
376 dz = (tr.getTime0() - tb) * vdriftTB;
377 if (tr.hasCSideClustersOnly()) {
378 dz = -dz;
379 }
380 //
381 int mnm = mNMoves - 1;
382 prepClus(tb); // clusters for MC time
383 (*mDBGOut) << "tpcMov"
384 << "counter=" << counter
385 << "copy=" << it
386 << "maxCopy=" << mnm
387 << "movTrackRef=" << trfm;
388 if (mRepRef) {
389 (*mDBGOut) << "tpcMov"
390 << "iniTrackRef=" << trf << "time=" << tr.getTime0();
391 }
392 (*mDBGOut) << "tpcMov"
393 << "side=" << side
394 << "imposedTB=" << tb
395 << "dz=" << dz
396 << "clX=" << clX
397 << "clY=" << clY
398 << "clZ=" << clZ
399 << "\n";
400 }
401 }
402}
403
405{
406 mDBGOut.reset();
407 mDBGOutCl.reset();
408}
409
411{
413 return;
414 }
415 if (mTPCVDriftHelper.accountCCDBInputs(matcher, obj)) {
416 return;
417 }
418 if (mTPCCorrMapsLoader.accountCCDBInputs(matcher, obj)) {
419 return;
420 }
421}
422
424{
425 std::vector<OutputSpec> outputs;
426 Options opts{
427 {"target-x", VariantType::Float, 70.f, {"Try to propagate to this radius"}},
428 {"n-moves", VariantType::Int, 6, {"Number of moves in allow range"}},
429 {"dump-clusters", VariantType::Bool, false, {"dump clusters"}},
430 {"tf-start", VariantType::Int, 0, {"1st TF to process"}},
431 {"tf-end", VariantType::Int, 999999999, {"last TF to process"}},
432 {"use-gpu-fitter", VariantType::Bool, false, {"use GPU track model for refit instead of TrackParCov"}},
433 {"repeat-ini-ref", VariantType::Bool, false, {"store ini-refit param with every moved track"}},
434 {"use-r-as-x", VariantType::Bool, false, {"Use radius instead of target sector X"}}};
435 auto dataRequest = std::make_shared<DataRequest>();
436
437 dataRequest->requestTracks(srcTracks, useMC);
438 dataRequest->requestClusters(srcClusters, useMC);
439 auto ggRequest = std::make_shared<o2::base::GRPGeomRequest>(false, // orbitResetTime
440 true, // GRPECS=true
441 false, // GRPLHCIF
442 true, // GRPMagField
443 true, // askMatLUT
445 dataRequest->inputs,
446 true);
447 o2::tpc::VDriftHelper::requestCCDBInputs(dataRequest->inputs);
448 o2::tpc::CorrectionMapsLoader::requestCCDBInputs(dataRequest->inputs, opts, sclOpts);
449
450 return DataProcessorSpec{
451 "tpc-track-study",
452 dataRequest->inputs,
453 outputs,
454 AlgorithmSpec{adaptFromTask<TPCTrackStudySpec>(dataRequest, ggRequest, sclOpts, srcTracks, useMC)},
455 opts};
456}
457
458} // namespace o2::trackstudy
Helper class to access load maps from CCDB.
Wrapper container for different reconstructed object types.
Definition of the GeometryManager class.
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...
Utility functions for MC particles.
Definition of the Names Generator class.
Definition of the parameter class for the detector electronics.
uint32_t side
Definition RawData.h:0
Wrapper container for different reconstructed object types.
Helper class to extract VDrift from different sources.
std::ostringstream debug
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:143
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)
recalculate inverse correction
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
void endOfStream(EndOfStreamContext &ec) final
This is invoked whenever we have an EndOfStream event.
TPCTrackStudySpec(std::shared_ptr< DataRequest > dr, std::shared_ptr< o2::base::GRPGeomRequest > gr, const o2::tpc::CorrectionMapsLoaderGloOpts &sclOpts, GTrackID::mask_t src, bool useMC)
void init(InitContext &ic) final
void run(ProcessingContext &pc) final
void process(o2::globaltracking::RecoContainer &recoData)
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
constexpr double LHCBunchSpacingMUS
Defining PrimaryVertex explicitly as messageable.
Definition TFIDInfo.h:20
std::vector< ConfigParamSpec > Options
std::tuple< T, T > rotateZ(T xL, T yL, T snAlp, T csAlp)
detail::Bracket< float > Bracketf_t
Definition Primitive2D.h:40
TrackParCovF TrackParCov
Definition Track.h:33
o2::math_utils::Bracketf_t TBracket
o2::framework::DataProcessorSpec getTPCTrackStudySpec(o2::dataformats::GlobalTrackID::mask_t srcTracks, o2::dataformats::GlobalTrackID::mask_t srcClus, bool useMC, const o2::tpc::CorrectionMapsLoaderGloOpts &sclOpts)
create a processor spec
o2::dataformats::VtxTrackRef V2TRef
o2::dataformats::VtxTrackIndex VTIndex
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
Defining DataPointCompositeObject explicitly as copiable.
std::unique_ptr< GPUReconstructionTimeframe > tf
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
unsigned int nClusters[constants::MAXSECTOR][constants::MAXGLOBALPADROW]
const ClusterNative * clusters[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::vector< int > row