Project
Loading...
Searching...
No Matches
TrackInterpolation.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
17
20#include "ITStracking/IOUtils.h"
23#include "TOFBase/Geo.h"
25#include "DataFormatsTPC/Defs.h"
30#include "MathUtils/Tsallis.h"
31#include "TRDBase/PadPlane.h"
32#include "TMath.h"
34#include "Framework/Logger.h"
36#include "GPUO2InterfaceUtils.h"
38#include "GPUO2InterfaceRefit.h"
39#include "GPUParam.h"
40#include "GPUParam.inc"
41#include <set>
42#include <algorithm>
43#include <random>
44
45using namespace o2::tpc;
48
49bool UnbinnedResid::gInitDone = false;
50
52{
53 if (!isITS()) {
55 }
56 // ITS alpha repends on the chip ID
59}
60
62{
63 if (isTPC()) {
64 return param::RowX[row];
65 }
67 if (isITS()) {
68 return o2::its::GeometryTGeo::Instance()->getSensorRefX(channel); // ITS X repends on the chip ID
69 }
70 if (isTRD()) {
71 auto geo = o2::trd::Geometry::instance();
72 ROOT::Math::Impl::Transform3D<double>::Point local{geo->cdrHght() - 0.5 - 0.279, 0., 0.}; // see TrackletTransformer::transformTracklet
73 return (geo->getMatrixT2L(channel) ^ local).X();
74 }
75 if (isTOF()) {
76 int det[5];
78 float pos[3] = {0.f, 0.f, 0.f};
80 float posl[3] = {pos[0], pos[1], pos[2]};
82 return pos[2]; // coordinates in sector frame: note that the rotation above puts z in pos[1], the radial coordinate in pos[2], and the tangent coordinate in pos[0] (this is to match the TOF residual system, where we don't use the radial component), so we swap their positions.
83 }
84 LOGP(fatal, "Did not recognize detector type: row:{}, sec:{}, channel:{}", row, sec, channel);
85 return 0.;
86}
87
89{
90 if (!gInitDone) {
91 LOGP(warn, "geometry initialization was not done, doing this for the current timestamp");
92 init();
93 if (!gInitDone) {
94 LOGP(fatal, "geometry initialization failed");
95 }
96 }
97}
98
99void UnbinnedResid::init(long timestamp)
100{
101 if (gInitDone) {
102 LOGP(warn, "Initialization was already done");
103 return;
104 }
105 if (!gGeoManager) {
106 o2::ccdb::BasicCCDBManager::instance().getSpecific<TGeoManager>("GLO/Config/GeometryAligned", timestamp);
107 }
108 auto geoTRD = o2::trd::Geometry::instance();
109 geoTRD->createPadPlaneArray();
110 geoTRD->createClusterMatrixArray();
111 gInitDone = true;
112}
113
115{
116 // perform initialization
117 LOG(info) << "Start initializing TrackInterpolation";
118 if (mInitDone) {
119 LOG(error) << "Initialization already performed.";
120 return;
121 }
122
123 const auto& elParam = ParameterElectronics::Instance();
124 mTPCTimeBinMUS = elParam.ZbinWidth;
125
126 mFastTransform = std::move(TPCFastTransformHelperO2::instance()->create(0));
127
128 mBz = o2::base::Propagator::Instance()->getNominalBz();
129 mRecoParam.init(mBz);
130 mGeoTRD = o2::trd::Geometry::instance();
132
133 mSourcesConfigured = src;
134 mSourcesConfiguredMap = srcMap;
135 mSingleSourcesConfigured = (mSourcesConfigured == mSourcesConfiguredMap);
136 mTrackTypes.insert({GTrackID::ITSTPC, 0});
137 mTrackTypes.insert({GTrackID::ITSTPCTRD, 1});
138 mTrackTypes.insert({GTrackID::ITSTPCTOF, 2});
139 mTrackTypes.insert({GTrackID::ITSTPCTRDTOF, 3});
140
142 geom->fillMatrixCache(o2::math_utils::bit2Mask(o2::math_utils::TransformType::T2L, o2::math_utils::TransformType::L2G));
143 mTPCParam = o2::gpu::GPUO2InterfaceUtils::getFullParamShared(0.f, mNHBPerTF);
144 mInitDone = true;
145 LOGP(info, "Done initializing TrackInterpolation. Configured track input: {}. Track input specifically for map: {}",
146 GTrackID::getSourcesNames(mSourcesConfigured), mSingleSourcesConfigured ? "identical" : GTrackID::getSourcesNames(mSourcesConfiguredMap));
147}
148
150{
151 LOGP(debug, "Check if input track {} is accepted", gid.asString());
152 bool hasOuterPoint = gidTable[GTrackID::TRD].isIndexSet() || gidTable[GTrackID::TOF].isIndexSet();
153 if (!hasOuterPoint && !mProcessITSTPConly) {
154 return false; // don't do ITS-only extrapolation through TPC
155 }
156 const auto itsTrk = mRecoCont->getITSTrack(gidTable[GTrackID::ITS]);
157 const auto tpcTrk = mRecoCont->getTPCTrack(gidTable[GTrackID::TPC]);
158
159 if (gidTable[GTrackID::TRD].isIndexSet()) {
160 // TRD specific cuts
161 const auto& trdTrk = mRecoCont->getITSTPCTRDTrack<o2::trd::TrackTRD>(gidTable[GTrackID::ITSTPCTRD]);
162 if (trdTrk.getNtracklets() < mParams->minTRDNTrklts) {
163 return false;
164 }
165 }
166 // reduced chi2 cut is the same for all track types
167 if (itsTrk.getChi2() / itsTrk.getNumberOfClusters() > mParams->maxITSChi2 || tpcTrk.getChi2() / tpcTrk.getNClusterReferences() > mParams->maxTPCChi2) {
168 return false;
169 }
170 if (!hasOuterPoint) {
171 // ITS-TPC track (does not have outer points in TRD or TOF)
172 if (itsTrk.getNumberOfClusters() < mParams->minITSNClsNoOuterPoint || tpcTrk.getNClusterReferences() < mParams->minTPCNClsNoOuterPoint) {
173 return false;
174 }
175 if (itsTrk.getPt() < mParams->minPtNoOuterPoint) {
176 return false;
177 }
178 } else {
179 if (itsTrk.getNumberOfClusters() < mParams->minITSNCls || tpcTrk.getNClusterReferences() < mParams->minTPCNCls) {
180 return false;
181 }
182 }
183
184 auto trc = mRecoCont->getTrackParam(gid);
186 auto prop = o2::base::Propagator::Instance();
187 if (!prop->propagateToDCA(pv, trc, prop->getNominalBz(), 2., o2::base::PropagatorF::MatCorrType::USEMatCorrLUT, &dca)) {
188 return false;
189 }
190 if (dca.getR2() > mParams->maxDCA * mParams->maxDCA) {
191 return false;
192 }
193 return true;
194}
195
197{
198 LOGP(debug, "Trying to find valid source for {} in {}", GTrackID::getSourceName(src), GTrackID::getSourcesNames(mask));
201 return GTrackID::ITSTPCTRD;
202 } else if (mask[GTrackID::ITSTPC]) {
203 return GTrackID::ITSTPC;
204 } else {
205 return GTrackID::NSources;
206 }
207 } else if (src == GTrackID::ITSTPCTRD || src == GTrackID::ITSTPCTOF) {
208 if (mask[GTrackID::ITSTPC]) {
209 return GTrackID::ITSTPC;
210 } else {
211 return GTrackID::NSources;
212 }
213 } else {
214 return GTrackID::NSources;
215 }
216}
217
219{
220 mRecoCont = &inp;
221 uint32_t nTrackSeeds = 0;
222 uint32_t countSeedCandidates[4] = {0};
223 auto pvvec = mRecoCont->getPrimaryVertices();
224 auto trackIndex = mRecoCont->getPrimaryVertexMatchedTracks(); // Global ID's for associated tracks
225 auto vtxRefs = mRecoCont->getPrimaryVertexMatchedTrackRefs(); // references from vertex to these track IDs
226 int nv = vtxRefs.size() - 1;
227 GTrackID::mask_t allowedSources = GTrackID::getSourcesMask("ITS-TPC,ITS-TPC-TRD,ITS-TPC-TOF,ITS-TPC-TRD-TOF") & mSourcesConfigured;
228 constexpr std::array<int, 3> SrcFast = {int(GTrackID::ITSTPCTRD), int(GTrackID::ITSTPCTOF), int(GTrackID::ITSTPCTRDTOF)};
229 if (mParams->refitITS) {
230 mITSRefitSeedID.resize(mRecoCont->getITSTracks().size(), -1);
231 }
232
233 for (int iv = 0; iv < nv; iv++) {
234 LOGP(debug, "processing PV {} of {}", iv, nv);
235
236 const auto& vtref = vtxRefs[iv];
237 auto pv = pvvec[iv];
238 if (mParams->minTOFTRDPVContributors > 0) { // we want only PVs constrained by fast detectors
239 int nfound = 0;
240 bool usePV = false;
241 for (uint32_t is = 0; is < SrcFast.size() && !usePV; is++) {
242 int src = SrcFast[is], idMin = vtref.getFirstEntryOfSource(src), idMax = idMin + vtref.getEntriesOfSource(src);
243 for (int i = idMin; i < idMax; i++) {
244 if (trackIndex[i].isPVContributor() && (++nfound == mParams->minTOFTRDPVContributors)) {
245 usePV = true;
246 break;
247 }
248 }
249 }
250 if (!usePV) {
251 continue;
252 }
253 }
254
255 for (int is = GTrackID::NSources; is >= 0; is--) {
256 if (!allowedSources[is]) {
257 continue;
258 }
259 LOGP(debug, "Checking source {}", is);
260 int idMin = vtref.getFirstEntryOfSource(is), idMax = idMin + vtref.getEntriesOfSource(is);
261 for (int i = idMin; i < idMax; i++) {
262 auto vid = trackIndex[i];
263 auto vidOrig = vid; // in case only ITS-TPC tracks are configured vid might be overwritten. We need to remember it for the PID
264 if (mParams->ignoreNonPVContrib && !vid.isPVContributor()) {
265 continue;
266 }
267 if (vid.isAmbiguous()) {
268 continue;
269 }
270 auto gidTable = mRecoCont->getSingleDetectorRefs(vid);
271 if (!mSourcesConfigured[is]) {
272 auto src = findValidSource(mSourcesConfigured, static_cast<GTrackID::Source>(vid.getSource()));
274 LOGP(debug, "prepareInputTrackSample: Found valid source {}", GTrackID::getSourceName(src));
275 vid = gidTable[src];
276 gidTable = mRecoCont->getSingleDetectorRefs(vid);
277 } else {
278 break; // no valid source for this vertex track source exists
279 }
280 }
281 ++countSeedCandidates[mTrackTypes[vid.getSource()]];
282 LOGP(debug, "Checking vid {}", vid.asString());
283 if (!isInputTrackAccepted(vid, gidTable, pv)) {
284 continue;
285 }
286 mSeeds.push_back(mRecoCont->getITSTrack(gidTable[GTrackID::ITS]).getParamOut());
287 mSeeds.back().setPID(mRecoCont->getTrackParam(vidOrig).getPID(), true);
288 mGIDs.push_back(vid);
289 mGIDtables.push_back(gidTable);
290 mTrackTimes.push_back(pv.getTimeStamp().getTimeStamp());
291 mTrackIndices[mTrackTypes[vid.getSource()]].push_back(nTrackSeeds++);
292 mTrackPVID.push_back(iv);
293 }
294 }
295 }
296
297 LOGP(info, "Created {} seeds. {} out of {} ITS-TPC-TRD-TOF, {} out of {} ITS-TPC-TRD, {} out of {} ITS-TPC-TOF, {} out of {} ITS-TPC",
298 nTrackSeeds,
299 mTrackIndices[mTrackTypes[GTrackID::ITSTPCTRDTOF]].size(), countSeedCandidates[mTrackTypes[GTrackID::ITSTPCTRDTOF]],
300 mTrackIndices[mTrackTypes[GTrackID::ITSTPCTRD]].size(), countSeedCandidates[mTrackTypes[GTrackID::ITSTPCTRD]],
301 mTrackIndices[mTrackTypes[GTrackID::ITSTPCTOF]].size(), countSeedCandidates[mTrackTypes[GTrackID::ITSTPCTOF]],
302 mTrackIndices[mTrackTypes[GTrackID::ITSTPC]].size(), countSeedCandidates[mTrackTypes[GTrackID::ITSTPC]]);
303}
304
306{
307 std::random_device rd;
308 std::mt19937 g(rd());
309 std::uniform_real_distribution<> distr(0., 1.);
310 float weight = 0;
312}
313
315{
316 // main processing function
317 if (!mInitDone) {
318 LOG(error) << "Initialization not yet done. Aborting...";
319 return;
320 }
321 // set the input containers
322 mTPCTracksClusIdx = mRecoCont->getTPCTracksClusterRefs();
323 mTPCClusterIdxStruct = &mRecoCont->getTPCClusters();
324 int nbOccTOT = o2::gpu::GPUO2InterfaceRefit::fillOccupancyMapGetSize(mNHBPerTF, mTPCParam.get());
325 o2::gpu::GPUO2InterfaceUtils::paramUseExternalOccupancyMap(mTPCParam.get(), mNHBPerTF, mRecoCont->occupancyMapTPC.data(), nbOccTOT);
326 mNTPCOccBinLength = mTPCParam->rec.tpc.occupancyMapTimeBins;
327 mNTPCOccBinLengthInv = 1.f / mNTPCOccBinLength;
328 {
329 if (!mITSDict) {
330 LOG(error) << "No ITS dictionary available";
331 return;
332 }
333 mITSTrackClusIdx = mRecoCont->getITSTracksClusterRefs();
334 const auto clusITS = mRecoCont->getITSClusters();
335 const auto patterns = mRecoCont->getITSClustersPatterns();
336 auto pattIt = patterns.begin();
337 mITSClustersArray.clear();
338 mITSClustersArray.reserve(clusITS.size());
339 LOGP(info, "We have {} ITS clusters and the number of patterns is {}", clusITS.size(), patterns.size());
340 o2::its::ioutils::convertCompactClusters(clusITS, pattIt, mITSClustersArray, mITSDict);
341 }
342
343 // In case we have more input tracks available than are required per TF
344 // we want to sample them. But we still prefer global ITS-TPC-TRD-TOF tracks
345 // over ITS-TPC-TRD tracks and so on. So we have to shuffle the indices
346 // in blocks.
347 std::random_device rd;
348 std::mt19937 g(rd());
349 std::vector<uint32_t> trackIndices; // here we keep the GIDs for all track types in a single vector to use in loop
350 std::shuffle(mTrackIndices[mTrackTypes[GTrackID::ITSTPCTRDTOF]].begin(), mTrackIndices[mTrackTypes[GTrackID::ITSTPCTRDTOF]].end(), g);
351 std::shuffle(mTrackIndices[mTrackTypes[GTrackID::ITSTPCTRD]].begin(), mTrackIndices[mTrackTypes[GTrackID::ITSTPCTRD]].end(), g);
352 std::shuffle(mTrackIndices[mTrackTypes[GTrackID::ITSTPCTOF]].begin(), mTrackIndices[mTrackTypes[GTrackID::ITSTPCTOF]].end(), g);
353 std::shuffle(mTrackIndices[mTrackTypes[GTrackID::ITSTPC]].begin(), mTrackIndices[mTrackTypes[GTrackID::ITSTPC]].end(), g);
354 trackIndices.insert(trackIndices.end(), mTrackIndices[mTrackTypes[GTrackID::ITSTPCTRDTOF]].begin(), mTrackIndices[mTrackTypes[GTrackID::ITSTPCTRDTOF]].end());
355 trackIndices.insert(trackIndices.end(), mTrackIndices[mTrackTypes[GTrackID::ITSTPCTRD]].begin(), mTrackIndices[mTrackTypes[GTrackID::ITSTPCTRD]].end());
356 trackIndices.insert(trackIndices.end(), mTrackIndices[mTrackTypes[GTrackID::ITSTPCTOF]].begin(), mTrackIndices[mTrackTypes[GTrackID::ITSTPCTOF]].end());
357 trackIndices.insert(trackIndices.end(), mTrackIndices[mTrackTypes[GTrackID::ITSTPC]].begin(), mTrackIndices[mTrackTypes[GTrackID::ITSTPC]].end());
358
359 int nSeeds = mSeeds.size(), lastChecked = 0;
360 mParentID.clear();
361 mParentID.resize(nSeeds, -1);
362
363 int maxOutputTracks = (mMaxTracksPerTF >= 0) ? mMaxTracksPerTF + mAddTracksForMapPerTF : nSeeds;
364 mTrackData.reserve(maxOutputTracks);
365 mClRes.reserve(maxOutputTracks * param::NPadRows);
366 mDetInfoRes.reserve(maxOutputTracks * param::NPadRows);
367 bool maxTracksReached = false;
368 for (int iSeed = 0; iSeed < nSeeds; ++iSeed) {
369 if (mMaxTracksPerTF >= 0 && mTrackDataCompact.size() >= mMaxTracksPerTF + mAddTracksForMapPerTF) {
370 LOG(info) << "Maximum number of tracks per TF reached. Skipping the remaining " << nSeeds - iSeed << " tracks.";
371 break;
372 }
373 int seedIndex = trackIndices[iSeed];
374 if (mParams->enableTrackDownsampling && !isTrackSelected(mSeeds[seedIndex])) {
375 continue;
376 }
377 auto addPart = [this, seedIndex](GTrackID::Source src) {
378 this->mGIDs.push_back(this->mGIDtables[seedIndex][src]);
379 this->mGIDtables.push_back(this->mRecoCont->getSingleDetectorRefs(this->mGIDs.back()));
380 this->mTrackTimes.push_back(this->mTrackTimes[seedIndex]);
381 this->mSeeds.push_back(this->mSeeds[seedIndex]);
382 this->mParentID.push_back(seedIndex); // store parent seed id
383 this->mTrackPVID.push_back(this->mTrackPVID[seedIndex]);
384 };
385
386 GTrackID::mask_t partsAdded;
387 if (!mSingleSourcesConfigured && !mSourcesConfiguredMap[mGIDs[seedIndex].getSource()]) {
388 auto src = findValidSource(mSourcesConfiguredMap, static_cast<GTrackID::Source>(mGIDs[seedIndex].getSource()));
390 LOGP(debug, "process {}: Found valid source {} for {} | nseeds:{} mSeeds:{} used: {}", iSeed, GTrackID::getSourceName(src), GTrackID::getSourceName(mGIDs[seedIndex].getSource()), nSeeds, mSeeds.size(), mTrackDataCompact.size());
391 addPart(src);
392 }
393 }
394 if (mMaxTracksPerTF >= 0 && mTrackDataCompact.size() >= mMaxTracksPerTF) {
395 if (!maxTracksReached) {
396 LOGP(info, "We already have reached mMaxTracksPerTF={}, but we continue to create seeds until mAddTracksForMapPerTF={} is also reached, iSeed: {} of {} inital seeds", mMaxTracksPerTF, mAddTracksForMapPerTF, iSeed, nSeeds);
397 }
398 maxTracksReached = true;
399 continue;
400 }
401 if (mGIDs[seedIndex].includesDet(DetID::TRD) || mGIDs[seedIndex].includesDet(DetID::TOF)) {
402 interpolateTrack(seedIndex);
403 LOGP(debug, "interpolateTrack {} {}, accepted: {}", iSeed, GTrackID::getSourceName(mGIDs[seedIndex].getSource()), mTrackDataCompact.size());
404 if (mProcessSeeds) {
405 if (mGIDs[seedIndex].includesDet(DetID::TRD) && mGIDs[seedIndex].includesDet(DetID::TOF) && !partsAdded[GTrackID::ITSTPCTRD]) {
406 addPart(GTrackID::ITSTPCTRD);
407 }
408 if (!partsAdded[GTrackID::ITSTPC]) {
409 addPart(GTrackID::ITSTPC);
410 }
411 }
412 } else {
413 extrapolateTrack(seedIndex);
414 LOGP(debug, "extrapolateTrack {} {}, accepted: {}", iSeed, GTrackID::getSourceName(mGIDs[seedIndex].getSource()), mTrackDataCompact.size());
415 }
416 lastChecked = iSeed;
417 }
418 std::vector<int> remSeeds;
419 if (mSeeds.size() > ++lastChecked) {
420 remSeeds.resize(mSeeds.size() - lastChecked);
421 std::iota(remSeeds.begin(), remSeeds.end(), lastChecked);
422 std::shuffle(remSeeds.begin(), remSeeds.end(), g);
423 LOGP(info, "Up to {} tracks out of {} additional seeds will be processed in random order, of which {} are stripped versions, accepted seeds: {}", mAddTracksForMapPerTF, remSeeds.size(), mSeeds.size() - nSeeds, mTrackDataCompact.size());
424 }
425 int extraChecked = 0;
426 for (int iSeed : remSeeds) {
427 if (mAddTracksForMapPerTF > 0 && mTrackDataCompact.size() >= mMaxTracksPerTF + mAddTracksForMapPerTF) {
428 LOGP(info, "Maximum number {} of additional tracks per TF reached. Skipping the remaining {} tracks", mAddTracksForMapPerTF, remSeeds.size() - extraChecked);
429 break;
430 }
431 extraChecked++;
432 if (mGIDs[iSeed].includesDet(DetID::TRD) || mGIDs[iSeed].includesDet(DetID::TOF)) {
433 interpolateTrack(iSeed);
434 LOGP(debug, "extra check {} of {}, seed {} interpolateTrack {}, used: {}", extraChecked, remSeeds.size(), iSeed, GTrackID::getSourceName(mGIDs[iSeed].getSource()), mTrackDataCompact.size());
435 } else {
436 LOGP(debug, "extra check {} of {}, seed {} extrapolateTrack {}, used: {}", extraChecked, remSeeds.size(), iSeed, GTrackID::getSourceName(mGIDs[iSeed].getSource()), mTrackDataCompact.size());
437 extrapolateTrack(iSeed);
438 }
439 }
440 LOG(info) << "Could process " << mTrackData.size() << " tracks successfully. " << mRejectedResiduals << " residuals were rejected. " << mClRes.size() << " residuals were accepted.";
441 mRejectedResiduals = 0;
442}
443
445{
446 LOGP(debug, "Starting track interpolation for GID {}", mGIDs[iSeed].asString());
447 TrackData trackData;
448 o2::trd::Tracklet64 trkl64;
450 std::unique_ptr<TrackDataExtended> trackDataExtended;
451 std::vector<TPCClusterResiduals> clusterResiduals;
452 auto propagator = o2::base::Propagator::Instance();
453 const auto& gidTable = mGIDtables[iSeed];
454 const auto& trkTPC = mRecoCont->getTPCTrack(gidTable[GTrackID::TPC]);
455 const auto& trkITS = mRecoCont->getITSTrack(gidTable[GTrackID::ITS]);
456 if (mDumpTrackPoints) {
457 trackDataExtended = std::make_unique<TrackDataExtended>();
458 (*trackDataExtended).gid = mGIDs[iSeed];
459 (*trackDataExtended).clIdx.setFirstEntry(mClRes.size());
460 (*trackDataExtended).trkITS = trkITS;
461 (*trackDataExtended).trkTPC = trkTPC;
462 auto nCl = trkITS.getNumberOfClusters();
463 auto clEntry = trkITS.getFirstClusterEntry();
464 for (int iCl = nCl - 1; iCl >= 0; iCl--) { // clusters are stored from outer to inner layers
465 const auto& clsITS = mITSClustersArray[mITSTrackClusIdx[clEntry + iCl]];
466 (*trackDataExtended).clsITS.push_back(clsITS);
467 }
468 }
469 if (mParams->refitITS && !refITSTrack(gidTable[GTrackID::ITS], iSeed)) {
470 return;
471 }
472 trackData.gid = mGIDs[iSeed];
473 trackData.par = mSeeds[iSeed];
474 auto trkWork = mSeeds[iSeed];
475 o2::track::TrackPar trkInner{trkWork};
476 // reset the cache array (sufficient to set cluster available to zero)
477 for (auto& elem : mCache) {
478 elem.clAvailable = 0;
479 }
480 trackData.clIdx.setFirstEntry(mClRes.size()); // reference the first cluster residual belonging to this track
481 float clusterTimeBinOffset = mTrackTimes[iSeed] / mTPCTimeBinMUS;
482
483 // store the TPC cluster positions in the cache, as well as dedx info
484 std::array<std::pair<uint16_t, uint16_t>, constants::MAXGLOBALPADROW> mCacheDEDX{};
485 std::array<short, constants::MAXGLOBALPADROW> multBins{};
486 for (int iCl = trkTPC.getNClusterReferences(); iCl--;) {
487 uint8_t sector, row;
488 uint32_t clusterIndexInRow;
489 const auto& clTPC = trkTPC.getCluster(mTPCTracksClusIdx, iCl, *mTPCClusterIdxStruct, sector, row);
490 float clTPCX;
491 std::array<float, 2> clTPCYZ;
492 mFastTransform->TransformIdeal(sector, row, clTPC.getPad(), clTPC.getTime(), clTPCX, clTPCYZ[0], clTPCYZ[1], clusterTimeBinOffset);
493 mCache[row].clSec = sector;
494 mCache[row].clAvailable = 1;
495 mCache[row].clY = clTPCYZ[0];
496 mCache[row].clZ = clTPCYZ[1];
497 mCache[row].clAngle = o2::math_utils::sector2Angle(sector);
498 mCacheDEDX[row].first = clTPC.getQtot();
499 mCacheDEDX[row].second = clTPC.getQmax();
500 int imb = int(clTPC.getTime() * mNTPCOccBinLengthInv);
501 if (imb < mTPCParam->occupancyMapSize) {
502 multBins[row] = 1 + std::max(0, imb);
503 }
504 }
505
506 // extrapolate seed through TPC and store track position at each pad row
507 for (int iRow = 0; iRow < param::NPadRows; ++iRow) {
508 if (!mCache[iRow].clAvailable) {
509 continue;
510 }
511 if (!trkWork.rotate(mCache[iRow].clAngle)) {
512 LOG(debug) << "Failed to rotate track during first extrapolation";
513 return;
514 }
515 if (!propagator->PropagateToXBxByBz(trkWork, param::RowX[iRow], mParams->maxSnp, mParams->maxStep, mMatCorr)) {
516 LOG(debug) << "Failed on first extrapolation";
517 return;
518 }
519 mCache[iRow].y[ExtOut] = trkWork.getY();
520 mCache[iRow].z[ExtOut] = trkWork.getZ();
521 mCache[iRow].sy2[ExtOut] = trkWork.getSigmaY2();
522 mCache[iRow].szy[ExtOut] = trkWork.getSigmaZY();
523 mCache[iRow].sz2[ExtOut] = trkWork.getSigmaZ2();
524 mCache[iRow].snp[ExtOut] = trkWork.getSnp();
525 // printf("Track alpha at row %i: %.2f, Y(%.2f), Z(%.2f)\n", iRow, trkWork.getAlpha(), trkWork.getY(), trkWork.getZ());
526 }
527
528 // start from outermost cluster with outer refit and back propagation
529 if (gidTable[GTrackID::TOF].isIndexSet()) {
530 LOG(debug) << "TOF point available";
531 const auto& clTOF = mRecoCont->getTOFClusters()[gidTable[GTrackID::TOF]];
532 if (mDumpTrackPoints) {
533 (*trackDataExtended).clsTOF = clTOF;
534 (*trackDataExtended).matchTOF = mRecoCont->getTOFMatch(mGIDs[iSeed]);
535 }
536 const int clTOFSec = clTOF.getCount();
537 const float clTOFAlpha = o2::math_utils::sector2Angle(clTOFSec);
538 if (!trkWork.rotate(clTOFAlpha)) {
539 LOG(debug) << "Failed to rotate into TOF cluster sector frame";
540 return;
541 }
542 float clTOFxyz[3] = {clTOF.getX(), clTOF.getY(), clTOF.getZ()};
543 if (!clTOF.isInNominalSector()) {
544 o2::tof::Geo::alignedToNominalSector(clTOFxyz, clTOFSec); // go from the aligned to nominal sector frame
545 }
546 std::array<float, 2> clTOFYZ{clTOFxyz[1], clTOFxyz[2]};
547 std::array<float, 3> clTOFCov{mParams->sigYZ2TOF, 0.f, mParams->sigYZ2TOF}; // assume no correlation between y and z and equal cluster error sigma^2 = (3cm)^2 / 12
548 if (!propagator->PropagateToXBxByBz(trkWork, clTOFxyz[0], mParams->maxSnp, mParams->maxStep, mMatCorr)) {
549 LOG(debug) << "Failed final propagation to TOF radius";
550 return;
551 }
552 // TODO: check if reset of covariance matrix is needed here (or, in case TOF point is not available at outermost TRD layer)
553 if (!trkWork.update(clTOFYZ, clTOFCov)) {
554 LOG(debug) << "Failed to update extrapolated ITS track with TOF cluster";
555 // LOGF(info, "trkWork.y=%f, cl.y=%f, trkWork.z=%f, cl.z=%f", trkWork.getY(), clTOFYZ[0], trkWork.getZ(), clTOFYZ[1]);
556 return;
557 }
558 }
559 if (gidTable[GTrackID::TRD].isIndexSet()) {
560 LOG(debug) << "TRD available";
561 const auto& trkTRD = mRecoCont->getITSTPCTRDTrack<o2::trd::TrackTRD>(gidTable[GTrackID::ITSTPCTRD]);
562 if (mDumpTrackPoints) {
563 (*trackDataExtended).trkTRD = trkTRD;
564 }
565 for (int iLayer = o2::trd::constants::NLAYER - 1; iLayer >= 0; --iLayer) {
566 std::array<float, 2> trkltTRDYZ{};
567 std::array<float, 3> trkltTRDCov{};
568 int res = processTRDLayer(trkTRD, iLayer, trkWork, &trkltTRDYZ, &trkltTRDCov);
569 if (res == -1) { // no TRD tracklet in this layer
570 continue;
571 }
572 if (res < -1) { // failed to reach this layer
573 return;
574 }
575 if (!trkWork.update(trkltTRDYZ, trkltTRDCov)) {
576 LOG(debug) << "Failed to update track at TRD layer " << iLayer;
577 return;
578 }
579 }
580 }
581
582 if (mDumpTrackPoints) {
583 (*trackDataExtended).trkOuter = trkWork;
584 }
585 auto trkOuter = trkWork; // outer param
586
587 // go back through the TPC and store updated track positions
588 bool outerParamStored = false;
589 for (int iRow = param::NPadRows; iRow--;) {
590 if (!mCache[iRow].clAvailable) {
591 continue;
592 }
593 if (mProcessSeeds && !outerParamStored) {
594 // for debug purposes we store the track parameters
595 // of the refitted ITS-(TRD)-(TOF) track at the
596 // outermose TPC cluster if we are processing all seeds
597 // i.e. if we in any case also process the ITS-TPC only
598 // part of the same track
599 trackData.par = trkWork;
600 outerParamStored = true;
601 }
602 if (!trkWork.rotate(mCache[iRow].clAngle)) {
603 LOG(debug) << "Failed to rotate track during back propagation";
604 return;
605 }
606 if (!propagator->PropagateToXBxByBz(trkWork, param::RowX[iRow], mParams->maxSnp, mParams->maxStep, mMatCorr)) {
607 LOG(debug) << "Failed on back propagation";
608 // printf("trkX(%.2f), clX(%.2f), clY(%.2f), clZ(%.2f), alphaTOF(%.2f)\n", trkWork.getX(), param::RowX[iRow], clTOFYZ[0], clTOFYZ[1], clTOFAlpha);
609 return;
610 }
611 mCache[iRow].y[ExtIn] = trkWork.getY();
612 mCache[iRow].z[ExtIn] = trkWork.getZ();
613 mCache[iRow].sy2[ExtIn] = trkWork.getSigmaY2();
614 mCache[iRow].szy[ExtIn] = trkWork.getSigmaZY();
615 mCache[iRow].sz2[ExtIn] = trkWork.getSigmaZ2();
616 mCache[iRow].snp[ExtIn] = trkWork.getSnp();
617 }
618
619 // calculate weighted mean at each pad row (assume for now y and z are uncorrelated) and store residuals to TPC clusters
620 unsigned short deltaRow = 0;
621 for (int iRow = 0; iRow < param::NPadRows; ++iRow) {
622 if (!mCache[iRow].clAvailable) {
623 ++deltaRow;
624 continue;
625 }
626 float wTotY = 1.f / mCache[iRow].sy2[ExtOut] + 1.f / mCache[iRow].sy2[ExtIn];
627 float wTotZ = 1.f / mCache[iRow].sz2[ExtOut] + 1.f / mCache[iRow].sz2[ExtIn];
628 mCache[iRow].y[Int] = (mCache[iRow].y[ExtOut] / mCache[iRow].sy2[ExtOut] + mCache[iRow].y[ExtIn] / mCache[iRow].sy2[ExtIn]) / wTotY;
629 mCache[iRow].z[Int] = (mCache[iRow].z[ExtOut] / mCache[iRow].sz2[ExtOut] + mCache[iRow].z[ExtIn] / mCache[iRow].sz2[ExtIn]) / wTotZ;
630
631 // simple average w/o weighting for angle
632 mCache[iRow].snp[Int] = (mCache[iRow].snp[ExtOut] + mCache[iRow].snp[ExtIn]) / 2.f;
633
634 const auto dY = mCache[iRow].clY - mCache[iRow].y[Int];
635 const auto dZ = mCache[iRow].clZ - mCache[iRow].z[Int];
636 const auto y = mCache[iRow].y[Int];
637 const auto z = mCache[iRow].z[Int];
638 const auto snp = mCache[iRow].snp[Int];
639 const auto sec = mCache[iRow].clSec;
640 clusterResiduals.emplace_back(dY, dZ, y, z, snp, sec, deltaRow);
641
642 deltaRow = 1;
643 }
644 trackData.chi2TRD = gidTable[GTrackID::TRD].isIndexSet() ? mRecoCont->getITSTPCTRDTrack<o2::trd::TrackTRD>(gidTable[GTrackID::ITSTPCTRD]).getChi2() : 0;
645 trackData.chi2TPC = trkTPC.getChi2();
646 trackData.chi2ITS = trkITS.getChi2();
647 trackData.nClsTPC = trkTPC.getNClusterReferences();
648 trackData.nClsITS = trkITS.getNumberOfClusters();
649 trackData.nTrkltsTRD = gidTable[GTrackID::TRD].isIndexSet() ? mRecoCont->getITSTPCTRDTrack<o2::trd::TrackTRD>(gidTable[GTrackID::ITSTPCTRD]).getNtracklets() : 0;
650
651 double t0forTOF = 0.; // to be set if TOF is matched
652 float t0forTOFwithinBC = 0.f;
653 float t0forTOFres = 9999.f;
654
655 if (gidTable[GTrackID::TOF].isIndexSet()) {
656 const auto& tofMatch = mRecoCont->getTOFMatch(mGIDs[iSeed]);
657 ULong64_t bclongtof = (tofMatch.getSignal() - 10000) * o2::tof::Geo::BC_TIME_INPS_INV;
658 t0forTOF = tofMatch.getFT0Best(); // setting t0 for TOF
659 t0forTOFwithinBC = t0forTOF - bclongtof * o2::tof::Geo::BC_TIME_INPS;
660 t0forTOFres = tofMatch.getFT0BestRes();
661 trackData.deltaTOF = tofMatch.getSignal() - t0forTOF - tofMatch.getLTIntegralOut().getTOF(trkTPC.getPID().getID());
662 trackData.clAvailTOF = uint16_t(t0forTOFres);
663 } else {
664 trackData.clAvailTOF = 0;
665 }
666 trackData.dEdxTPC = trkTPC.getdEdx().dEdxTotTPC;
667
668 TrackParams params; // for refitted track parameters and flagging rejected clusters
669 if (mParams->skipOutlierFiltering || validateTrack(trackData, params, clusterResiduals)) {
670 // track is good
671 int nClValidated = 0;
672 int iRow = 0;
673 for (unsigned int iCl = 0; iCl < clusterResiduals.size(); ++iCl) {
674 iRow += clusterResiduals[iCl].dRow;
675 if (params.flagRej[iCl]) {
676 // skip masked cluster residual
677 continue;
678 }
679 const float tgPhi = clusterResiduals[iCl].snp / std::sqrt((1.f - clusterResiduals[iCl].snp) * (1.f + clusterResiduals[iCl].snp));
680 const auto dy = clusterResiduals[iCl].dy;
681 const auto dz = clusterResiduals[iCl].dz;
682 const auto y = clusterResiduals[iCl].y;
683 const auto z = clusterResiduals[iCl].z;
684 const auto sec = clusterResiduals[iCl].sec;
685 if ((std::abs(dy) < param::MaxResid) && (std::abs(dz) < param::MaxResid) && (std::abs(y) < param::MaxY) && (std::abs(z) < param::MaxZ) && (std::abs(tgPhi) < param::MaxTgSlp)) {
686 mClRes.emplace_back(dy, dz, tgPhi, y, z, iRow, sec);
687 mDetInfoRes.emplace_back().setTPC(mCacheDEDX[iRow].first, mCacheDEDX[iRow].second); // qtot, qmax
688 ++nClValidated;
689 } else {
690 ++mRejectedResiduals;
691 }
692 }
693 trackData.clIdx.setEntries(nClValidated);
694
695 // store multiplicity info
696 for (int ist = 0; ist < NSTACKS; ist++) {
697 int mltBinMin = 0x7ffff, mltBinMax = -1, prevBin = -1;
698 for (int ir = STACKROWS[ist]; ir < STACKROWS[ist + 1]; ir++) {
699 if (multBins[ir] != prevBin && multBins[ir] > 0) { // there is a cluster different from previous one
700 prevBin = multBins[ir];
701 if (multBins[ir] > mltBinMax) {
702 mltBinMax = multBins[ir];
703 }
704 if (multBins[ir] < mltBinMin) {
705 mltBinMin = multBins[ir];
706 }
707 }
708 }
709 if (--mltBinMin >= 0) { // we were offsetting bin IDs by 1!
710 float avMlt = 0;
711 for (int ib = mltBinMin; ib < mltBinMax; ib++) {
712 avMlt += mTPCParam->occupancyMap[ib];
713 }
714 avMlt /= (mltBinMax - mltBinMin);
715 trackData.setMultStack(avMlt, ist);
716 }
717 }
718
719 bool stopPropagation = !mExtDetResid;
720 if (!stopPropagation) {
721 // do we have TRD residuals to add?
722 trkWork = trkOuter;
723 if (gidTable[GTrackID::TRD].isIndexSet()) {
724 const auto& trkTRD = mRecoCont->getITSTPCTRDTrack<o2::trd::TrackTRD>(gidTable[GTrackID::ITSTPCTRD]);
725 for (int iLayer = 0; iLayer < o2::trd::constants::NLAYER; iLayer++) {
726 std::array<float, 2> trkltTRDYZ{};
727 int res = processTRDLayer(trkTRD, iLayer, trkWork, &trkltTRDYZ, nullptr, &trackData, &trkl64, &trklCalib);
728 if (res == -1) { // no traklet on this layer
729 continue;
730 }
731 if (res < -1) { // failed to reach this layer
732 stopPropagation = true;
733 break;
734 }
735
736 float tgPhi = trkWork.getSnp() / std::sqrt((1.f - trkWork.getSnp()) * (1.f + trkWork.getSnp()));
737 auto dy = trkltTRDYZ[0] - trkWork.getY();
738 auto dz = trkltTRDYZ[1] - trkWork.getZ();
739 if ((std::abs(dy) < param::MaxResid) && (std::abs(dz) < param::MaxResid) && (std::abs(trkWork.getY()) < param::MaxY) && (std::abs(trkWork.getZ()) < param::MaxZ) && (std::abs(tgPhi) < param::MaxTgSlp)) {
740 mClRes.emplace_back(dy, dz, tgPhi, trkWork.getY(), trkWork.getZ(), 160 + iLayer, o2::math_utils::angle2Sector(trkWork.getAlpha()), (short)res);
741 mDetInfoRes.emplace_back().setTRD(trkl64.getQ0(), trkl64.getQ1(), trkl64.getQ2(), trklCalib.getDy()); // q0,q1,q2,slope
742 trackData.nExtDetResid++;
743 }
744 }
745 }
746
747 // do we have TOF residual to add?
748 while (gidTable[GTrackID::TOF].isIndexSet() && !stopPropagation) {
749 const auto& clTOF = mRecoCont->getTOFClusters()[gidTable[GTrackID::TOF]];
750 float clTOFxyz[3] = {clTOF.getX(), clTOF.getY(), clTOF.getZ()};
751 if (!clTOF.isInNominalSector()) {
752 o2::tof::Geo::alignedToNominalSector(clTOFxyz, clTOF.getCount()); // go from the aligned to nominal sector frame
753 }
754 const float clTOFAlpha = o2::math_utils::sector2Angle(clTOF.getCount());
755 if (trkWork.getAlpha() != clTOFAlpha && !trkWork.rotate(clTOFAlpha)) {
756 LOG(debug) << "Failed to rotate into TOF cluster sector frame";
757 stopPropagation = true;
758 break;
759 }
760 if (!propagator->PropagateToXBxByBz(trkWork, clTOFxyz[0], mParams->maxSnp, mParams->maxStep, mMatCorr)) {
761 LOG(debug) << "Failed final propagation to TOF radius";
762 break;
763 }
764
765 float tgPhi = trkWork.getSnp() / std::sqrt((1.f - trkWork.getSnp()) * (1.f + trkWork.getSnp()));
766 auto dy = clTOFxyz[1] - trkWork.getY();
767 auto dz = clTOFxyz[2] - trkWork.getZ();
768 // get seeding track time
769
770 if ((std::abs(dy) < param::MaxResid) && (std::abs(dz) < param::MaxResid) && (std::abs(trkWork.getY()) < param::MaxY) && (std::abs(trkWork.getZ()) < param::MaxZ) && (std::abs(tgPhi) < param::MaxTgSlp)) {
771 mClRes.emplace_back(dy, dz, tgPhi, trkWork.getY(), trkWork.getZ(), 170, clTOF.getCount(), clTOF.getPadInSector());
772 // get seeding track time
773 if (!gidTable[GTrackID::ITSTPC].isIndexSet()) {
774 LOGP(fatal, "ITS-TPC seed index is not set for TOF track");
775 }
776 float tdif = static_cast<float>(clTOF.getTime() - t0forTOF); // time in \mus wrt interaction time0
777 mDetInfoRes.emplace_back().setTOF(tdif * 1e-6);
778 trackData.nExtDetResid++;
779 }
780 break;
781 }
782
783 // add ITS residuals
784 while (!stopPropagation) {
785 auto& trkWorkITS = trkInner; // this is ITS outer param
786 auto nCl = trkITS.getNumberOfClusters();
787 auto clEntry = trkITS.getFirstClusterEntry();
789 for (int iCl = 0; iCl < nCl; iCl++) { // clusters are stored from outer to inner layers
790 const auto& cls = mITSClustersArray[mITSTrackClusIdx[clEntry + iCl]];
791 int chip = cls.getSensorID();
792 float chipX, chipAlpha;
793 geom->getSensorXAlphaRefPlane(cls.getSensorID(), chipX, chipAlpha);
794 if (!trkWorkITS.rotate(chipAlpha) || !propagator->PropagateToXBxByBz(trkWorkITS, chipX, mParams->maxSnp, mParams->maxStep, mMatCorr)) {
795 LOGP(debug, "Failed final propagation to ITS X={} alpha={}", chipX, chipAlpha);
796 stopPropagation = true;
797 break;
798 }
799 float tgPhi = trkWorkITS.getSnp() / std::sqrt((1.f - trkWorkITS.getSnp()) * (1.f + trkWorkITS.getSnp()));
800 auto dy = cls.getY() - trkWorkITS.getY();
801 auto dz = cls.getZ() - trkWorkITS.getZ();
802 if ((std::abs(dy) < param::MaxResid) && (std::abs(dz) < param::MaxResid) && (std::abs(trkWorkITS.getY()) < param::MaxY) && (std::abs(trkWorkITS.getZ()) < param::MaxZ) && (std::abs(tgPhi) < param::MaxTgSlp)) {
803 mClRes.emplace_back(dy, dz, tgPhi, trkWorkITS.getY(), trkWorkITS.getZ(), 180 + geom->getLayer(cls.getSensorID()), -1, cls.getSensorID());
804 mDetInfoRes.emplace_back(); // empty placeholder
805 trackData.nExtDetResid++;
806 }
807 }
808 if (!stopPropagation) { // add residual to PV
809 const auto& pv = mRecoCont->getPrimaryVertices()[mTrackPVID[iSeed]];
810 o2::math_utils::Point3D<float> vtx{pv.getX(), pv.getY(), pv.getZ()};
811 if (!propagator->propagateToDCA(vtx, trkWorkITS, mBz, mParams->maxStep, mMatCorr)) {
812 LOGP(debug, "Failed propagation to DCA to PV ({} {} {}), {}", pv.getX(), pv.getY(), pv.getZ(), trkWorkITS.asString());
813 stopPropagation = true;
814 break;
815 }
816 // rotate PV to the track frame
817 float sn, cs, alpha = trkWorkITS.getAlpha();
818 math_utils::detail::bringToPMPi(alpha);
819 math_utils::detail::sincos<float>(alpha, sn, cs);
820 float xv = vtx.X() * cs + vtx.Y() * sn, yv = -vtx.X() * sn + vtx.Y() * cs, zv = vtx.Z();
821 auto dy = yv - trkWorkITS.getY();
822 auto dz = zv - trkWorkITS.getZ();
823 if ((std::abs(dy) < param::MaxResid) && (std::abs(dz) < param::MaxResid) && (std::abs(trkWorkITS.getY()) < param::MaxY) && (std::abs(trkWorkITS.getZ()) < param::MaxZ) && abs(xv) < param::MaxVtxX) {
824 short compXV = static_cast<short>(xv * 0x7fff / param::MaxVtxX);
825 mClRes.emplace_back(dy, dz, alpha / TMath::Pi(), trkWorkITS.getY(), trkWorkITS.getZ(), 190, -1, compXV);
826 if (!gidTable[GTrackID::ITSTPC].isIndexSet()) {
827 LOGP(fatal, "ITS-TPC seed index is not set for TOF track");
828 }
829 float tdif = pv.getTimeStamp().getTimeStamp() - mRecoCont->getTPCITSTrack(gidTable[GTrackID::ITSTPC]).getTimeMUS().getTimeStamp();
830 mDetInfoRes.emplace_back().setPV(tdif); // time in \mus wrt seeding ITS-TPC track
831 trackData.nExtDetResid++;
832 }
833 }
834 break;
835 }
836 }
837
838 mGIDsSuccess.push_back(mGIDs[iSeed]);
839 mTrackDataCompact.emplace_back(trackData.clIdx.getFirstEntry(), trackData.multStack, nClValidated, mGIDs[iSeed].getSource(), trackData.nExtDetResid);
840 mTrackData.push_back(std::move(trackData));
841 if (mDumpTrackPoints) {
842 (*trackDataExtended).clIdx.setEntries(nClValidated);
843 (*trackDataExtended).nExtDetResid = trackData.nExtDetResid;
844 mTrackDataExtended.push_back(std::move(*trackDataExtended));
845 }
846 }
847 if (mParams->writeUnfiltered) {
848 TrackData trkDataTmp = trackData;
849 trkDataTmp.clIdx.setFirstEntry(mClResUnfiltered.size());
850 trkDataTmp.clIdx.setEntries(clusterResiduals.size());
851 mTrackDataUnfiltered.push_back(std::move(trkDataTmp));
852 mClResUnfiltered.insert(mClResUnfiltered.end(), clusterResiduals.begin(), clusterResiduals.end());
853 }
854}
855
857 std::array<float, 2>* trkltTRDYZ, std::array<float, 3>* trkltTRDCov, TrackData* trkData,
859{
860 // return chamber ID (0:539) in case of successful processing, -1 if there is no TRD tracklet at given layer, -2 if processing failed
861 int trkltIdx = trkTRD.getTrackletIndex(iLayer);
862 if (trkltIdx < 0) {
863 return -1; // no TRD tracklet in this layer
864 }
865 const auto& trdSP = mRecoCont->getTRDCalibratedTracklets()[trkltIdx];
866 const auto& trdTrklt = mRecoCont->getTRDTracklets()[trkltIdx];
867 auto trkltDet = trdTrklt.getDetector();
868 auto trkltSec = trkltDet / (o2::trd::constants::NLAYER * o2::trd::constants::NSTACK);
869 if (trkltSec != o2::math_utils::angle2Sector(trkWork.getAlpha())) {
870 if (!trkWork.rotate(o2::math_utils::sector2Angle(trkltSec))) {
871 LOG(debug) << "Track could not be rotated in TRD tracklet coordinate system in layer " << iLayer;
872 return -2;
873 }
874 }
875 if (!o2::base::Propagator::Instance()->PropagateToXBxByBz(trkWork, trdSP.getX(), mParams->maxSnp, mParams->maxStep, mMatCorr)) {
876 LOG(debug) << "Failed propagation to TRD layer " << iLayer;
877 return -2;
878 }
879 if (trkltTRDYZ) {
880 const auto* pad = mGeoTRD->getPadPlane(trkltDet);
881 float tilt = tan(TMath::DegToRad() * pad->getTiltingAngle()); // tilt is signed! and returned in degrees
882 float tiltCorrUp = tilt * (trdSP.getZ() - trkWork.getZ());
883 float zPosCorrUp = trdSP.getZ() + mRecoParam.getZCorrCoeffNRC() * trkWork.getTgl(); // maybe Z can be corrected on avarage already by the tracklet transformer?
884 float padLength = pad->getRowSize(trdTrklt.getPadRow());
885 if (!((trkWork.getSigmaZ2() < (padLength * padLength / 12.f)) && (std::fabs(trdSP.getZ() - trkWork.getZ()) < padLength))) {
886 tiltCorrUp = 0.f;
887 }
888 (*trkltTRDYZ)[0] = trdSP.getY() - tiltCorrUp;
889 (*trkltTRDYZ)[1] = zPosCorrUp;
890 if (trkltTRDCov) {
891 mRecoParam.recalcTrkltCov(tilt, trkWork.getSnp(), pad->getRowSize(trdTrklt.getPadRow()), *trkltTRDCov);
892 }
893 }
894 if (trkData) {
895 auto slope = trdSP.getDy();
896 if (std::abs(slope) < param::MaxTRDSlope) {
897 trkData->TRDTrkltSlope[iLayer] = slope * 0x7fff / param::MaxTRDSlope;
898 }
899 }
900 if (trk64) {
901 *trk64 = trdTrklt;
902 }
903 if (trkCalib) {
904 *trkCalib = trdSP;
905 }
906 return trkltDet;
907}
908
910{
911 // extrapolate ITS-only track through TPC and store residuals to TPC clusters in the output vectors
912 LOGP(debug, "Starting track extrapolation for GID {}", mGIDs[iSeed].asString());
913 const auto& gidTable = mGIDtables[iSeed];
914 TrackData trackData;
915 o2::trd::Tracklet64 trkl64;
917 std::unique_ptr<TrackDataExtended> trackDataExtended;
918 std::vector<TPCClusterResiduals> clusterResiduals;
919 trackData.clIdx.setFirstEntry(mClRes.size());
920 const auto& trkITS = mRecoCont->getITSTrack(gidTable[GTrackID::ITS]);
921 const auto& trkTPC = mRecoCont->getTPCTrack(gidTable[GTrackID::TPC]);
922 if (mDumpTrackPoints) {
923 trackDataExtended = std::make_unique<TrackDataExtended>();
924 (*trackDataExtended).gid = mGIDs[iSeed];
925 (*trackDataExtended).clIdx.setFirstEntry(mClRes.size());
926 (*trackDataExtended).trkITS = trkITS;
927 (*trackDataExtended).trkTPC = trkTPC;
928 auto nCl = trkITS.getNumberOfClusters();
929 auto clEntry = trkITS.getFirstClusterEntry();
930 for (int iCl = nCl - 1; iCl >= 0; iCl--) { // clusters are stored from outer to inner layers
931 const auto& clsITS = mITSClustersArray[mITSTrackClusIdx[clEntry + iCl]];
932 (*trackDataExtended).clsITS.push_back(clsITS);
933 }
934 }
935 if (mParams->refitITS && !refITSTrack(gidTable[GTrackID::ITS], iSeed)) {
936 return;
937 }
938 trackData.gid = mGIDs[iSeed];
939 trackData.par = mSeeds[iSeed];
940
941 auto trkWork = mSeeds[iSeed];
942 float clusterTimeBinOffset = mTrackTimes[iSeed] / mTPCTimeBinMUS;
943 auto propagator = o2::base::Propagator::Instance();
944 unsigned short rowPrev = 0; // used to calculate dRow of two consecutive cluster residuals
945 unsigned short nMeasurements = 0;
946 uint8_t clRowPrev = constants::MAXGLOBALPADROW; // used to identify and skip split clusters on the same pad row
947 std::array<std::pair<uint16_t, uint16_t>, constants::MAXGLOBALPADROW> mCacheDEDX{};
948 std::array<short, constants::MAXGLOBALPADROW> multBins{};
949 for (int iCl = trkTPC.getNClusterReferences(); iCl--;) {
950 uint8_t sector, row;
951 uint32_t clusterIndexInRow;
952 const auto& cl = trkTPC.getCluster(mTPCTracksClusIdx, iCl, *mTPCClusterIdxStruct, sector, row);
953 if (clRowPrev == row) {
954 // if there are split clusters we only take the first one on the pad row
955 continue;
956 } else if (clRowPrev < constants::MAXGLOBALPADROW && clRowPrev > row) {
957 // we seem to be looping, abort this track
958 LOGP(debug, "TPC track with pT={} GeV and {} clusters has cluster {} on row {} while the previous cluster was on row {}",
959 mSeeds[iSeed].getPt(), trkTPC.getNClusterReferences(), iCl, row, clRowPrev);
960 return;
961 } else {
962 // this is the first cluster we see on this pad row
963 clRowPrev = row;
964 }
965 float x = 0, y = 0, z = 0;
966 mFastTransform->TransformIdeal(sector, row, cl.getPad(), cl.getTime(), x, y, z, clusterTimeBinOffset);
967 if (!trkWork.rotate(o2::math_utils::sector2Angle(sector))) {
968 return;
969 }
970 if (!propagator->PropagateToXBxByBz(trkWork, x, mParams->maxSnp, mParams->maxStep, mMatCorr)) {
971 return;
972 }
973
974 const auto dY = y - trkWork.getY();
975 const auto dZ = z - trkWork.getZ();
976 const auto ty = trkWork.getY();
977 const auto tz = trkWork.getZ();
978 const auto snp = trkWork.getSnp();
979 const auto sec = sector;
980 clusterResiduals.emplace_back(dY, dZ, ty, tz, snp, sec, row - rowPrev);
981 mCacheDEDX[row].first = cl.getQtot();
982 mCacheDEDX[row].second = cl.getQmax();
983 rowPrev = row;
984 int imb = int(cl.getTime() * mNTPCOccBinLengthInv);
985 if (imb < mTPCParam->occupancyMapSize) {
986 multBins[row] = 1 + std::max(0, imb);
987 }
988 ++nMeasurements;
989 }
990
991 TrackParams params; // for refitted track parameters and flagging rejected clusters
992 if (clusterResiduals.size() > constants::MAXGLOBALPADROW) {
993 LOGP(warn, "Extrapolated ITS-TPC track and found more residuals than possible ({})", clusterResiduals.size());
994 return;
995 }
996
997 trackData.chi2TPC = trkTPC.getChi2();
998 trackData.chi2ITS = trkITS.getChi2();
999 trackData.nClsTPC = trkTPC.getNClusterReferences();
1000 trackData.nClsITS = trkITS.getNumberOfClusters();
1001 trackData.clIdx.setEntries(nMeasurements);
1002 trackData.dEdxTPC = trkTPC.getdEdx().dEdxTotTPC;
1003 if (mDumpTrackPoints) {
1004 (*trackDataExtended).trkOuter = trkWork;
1005 }
1006
1007 if (mParams->skipOutlierFiltering || validateTrack(trackData, params, clusterResiduals)) {
1008 // track is good, store TPC part
1009
1010 int nClValidated = 0, iRow = 0;
1011 unsigned int iCl = 0;
1012 for (iCl = 0; iCl < clusterResiduals.size(); ++iCl) {
1013 iRow += clusterResiduals[iCl].dRow;
1014 if (iRow < param::NPadRows && params.flagRej[iCl]) { // skip masked cluster residual
1015 continue;
1016 }
1017 const float tgPhi = clusterResiduals[iCl].snp / std::sqrt((1.f - clusterResiduals[iCl].snp) * (1.f + clusterResiduals[iCl].snp));
1018 const auto dy = clusterResiduals[iCl].dy;
1019 const auto dz = clusterResiduals[iCl].dz;
1020 const auto y = clusterResiduals[iCl].y;
1021 const auto z = clusterResiduals[iCl].z;
1022 if ((std::abs(dy) < param::MaxResid) && (std::abs(dz) < param::MaxResid) && (std::abs(y) < param::MaxY) && (std::abs(z) < param::MaxZ) && (std::abs(tgPhi) < param::MaxTgSlp)) {
1023 mClRes.emplace_back(dy, dz, tgPhi, y, z, iRow, clusterResiduals[iCl].sec);
1024 mDetInfoRes.emplace_back().setTPC(mCacheDEDX[iRow].first, mCacheDEDX[iRow].second); // qtot, qmax
1025 ++nClValidated;
1026 } else {
1027 ++mRejectedResiduals;
1028 }
1029 }
1030 trackData.clIdx.setEntries(nClValidated);
1031
1032 // store multiplicity info
1033 for (int ist = 0; ist < NSTACKS; ist++) {
1034 int mltBinMin = 0x7ffff, mltBinMax = -1, prevBin = -1;
1035 for (int ir = STACKROWS[ist]; ir < STACKROWS[ist + 1]; ir++) {
1036 if (multBins[ir] != prevBin && multBins[ir] > 0) { // there is a cluster
1037 prevBin = multBins[ir];
1038 if (multBins[ir] > mltBinMax) {
1039 mltBinMax = multBins[ir];
1040 }
1041 if (multBins[ir] < mltBinMin) {
1042 mltBinMin = multBins[ir];
1043 }
1044 }
1045 }
1046 if (--mltBinMin >= 0) { // we were offsetting bin IDs by 1!
1047 float avMlt = 0;
1048 for (int ib = mltBinMin; ib < mltBinMax; ib++) {
1049 avMlt += mTPCParam->occupancyMap[ib];
1050 }
1051 avMlt /= (mltBinMax - mltBinMin);
1052 trackData.setMultStack(avMlt, ist);
1053 }
1054 }
1055
1056 bool stopPropagation = !mExtDetResid;
1057 if (!stopPropagation) {
1058 // do we have TRD residuals to add?
1059 int iSeedFull = mParentID[iSeed] == -1 ? iSeed : mParentID[iSeed];
1060 auto gidFull = mGIDs[iSeedFull];
1061 const auto& gidTableFull = mGIDtables[iSeedFull];
1062 if (gidTableFull[GTrackID::TRD].isIndexSet()) {
1063 const auto& trkTRD = mRecoCont->getITSTPCTRDTrack<o2::trd::TrackTRD>(gidTableFull[GTrackID::ITSTPCTRD]);
1064 trackData.nTrkltsTRD = trkTRD.getNtracklets();
1065 for (int iLayer = 0; iLayer < o2::trd::constants::NLAYER; iLayer++) {
1066 std::array<float, 2> trkltTRDYZ{};
1067 int res = processTRDLayer(trkTRD, iLayer, trkWork, &trkltTRDYZ, nullptr, &trackData, &trkl64, &trklCalib);
1068 if (res == -1) { // no traklet on this layer
1069 continue;
1070 }
1071 if (res < -1) { // failed to reach this layer
1072 stopPropagation = true;
1073 break;
1074 }
1075
1076 float tgPhi = trkWork.getSnp() / std::sqrt((1.f - trkWork.getSnp()) * (1.f + trkWork.getSnp()));
1077 auto dy = trkltTRDYZ[0] - trkWork.getY();
1078 auto dz = trkltTRDYZ[1] - trkWork.getZ();
1079 const auto sec = clusterResiduals[iCl].sec;
1080 if ((std::abs(dy) < param::MaxResid) && (std::abs(dz) < param::MaxResid) && (std::abs(trkWork.getY()) < param::MaxY) && (std::abs(trkWork.getZ()) < param::MaxZ) && (std::abs(tgPhi) < param::MaxTgSlp)) {
1081 mClRes.emplace_back(dy, dz, tgPhi, trkWork.getY(), trkWork.getZ(), 160 + iLayer, o2::math_utils::angle2Sector(trkWork.getAlpha()), (short)res);
1082 mDetInfoRes.emplace_back().setTRD(trkl64.getQ0(), trkl64.getQ1(), trkl64.getQ2(), trklCalib.getDy()); // q0,q1,q2,slope
1083 trackData.nExtDetResid++;
1084 }
1085 }
1086 }
1087
1088 // do we have TOF residual to add?
1089 trackData.clAvailTOF = 0;
1090 while (gidTableFull[GTrackID::TOF].isIndexSet() && !stopPropagation) {
1091 const auto& tofMatch = mRecoCont->getTOFMatch(gidFull);
1092 ULong64_t bclongtof = (tofMatch.getSignal() - 10000) * o2::tof::Geo::BC_TIME_INPS_INV;
1093 double t0forTOF = tofMatch.getFT0Best(); // setting t0 for TOF
1094 float t0forTOFwithinBC = t0forTOF - bclongtof * o2::tof::Geo::BC_TIME_INPS;
1095 float t0forTOFres = tofMatch.getFT0BestRes();
1096 trackData.deltaTOF = tofMatch.getSignal() - t0forTOF - tofMatch.getLTIntegralOut().getTOF(trkTPC.getPID().getID());
1097 trackData.clAvailTOF = uint16_t(t0forTOFres);
1098 const auto& clTOF = mRecoCont->getTOFClusters()[gidTableFull[GTrackID::TOF]];
1099 const float clTOFAlpha = o2::math_utils::sector2Angle(clTOF.getCount());
1100 float clTOFxyz[3] = {clTOF.getX(), clTOF.getY(), clTOF.getZ()};
1101 if (!clTOF.isInNominalSector()) {
1102 o2::tof::Geo::alignedToNominalSector(clTOFxyz, clTOF.getCount()); // go from the aligned to nominal sector frame
1103 }
1104 if (trkWork.getAlpha() != clTOFAlpha && !trkWork.rotate(clTOFAlpha)) {
1105 LOG(debug) << "Failed to rotate into TOF cluster sector frame";
1106 stopPropagation = true;
1107 break;
1108 }
1109 if (!propagator->PropagateToXBxByBz(trkWork, clTOFxyz[0], mParams->maxSnp, mParams->maxStep, mMatCorr)) {
1110 LOG(debug) << "Failed final propagation to TOF radius";
1111 break;
1112 }
1113
1114 float tgPhi = trkWork.getSnp() / std::sqrt((1.f - trkWork.getSnp()) * (1.f + trkWork.getSnp()));
1115 auto dy = clTOFxyz[1] - trkWork.getY();
1116 auto dz = clTOFxyz[2] - trkWork.getZ();
1117 if ((std::abs(dy) < param::MaxResid) && (std::abs(dz) < param::MaxResid) && (std::abs(trkWork.getY()) < param::MaxY) && (std::abs(trkWork.getZ()) < param::MaxZ) && (std::abs(tgPhi) < param::MaxTgSlp)) {
1118 mClRes.emplace_back(dy, dz, tgPhi, trkWork.getY(), trkWork.getZ(), 170, clTOF.getCount(), clTOF.getPadInSector());
1119 // get seeding track time
1120 if (!gidTableFull[GTrackID::ITSTPC].isIndexSet()) {
1121 LOGP(fatal, "ITS-TPC seed index is not set for TOF track");
1122 }
1123
1124 float tdif = static_cast<float>(clTOF.getTime() - t0forTOF); // time in \mus wrt interaction time0
1125 mDetInfoRes.emplace_back().setTOF(tdif * 1e-6); // time in \mus wrt seeding ITS-TPC track
1126 trackData.nExtDetResid++;
1127 }
1128 break;
1129 }
1130
1131 // add ITS residuals
1132 while (!stopPropagation) {
1133 o2::track::TrackPar trkWorkITS{trackData.par}; // this is ITS outer param
1134 auto nCl = trkITS.getNumberOfClusters();
1135 auto clEntry = trkITS.getFirstClusterEntry();
1136 auto geom = o2::its::GeometryTGeo::Instance();
1137 for (int iCl = 0; iCl < nCl; iCl++) { // clusters are stored from outer to inner layers
1138 const auto& cls = mITSClustersArray[mITSTrackClusIdx[clEntry + iCl]];
1139 int chip = cls.getSensorID();
1140 float chipX, chipAlpha;
1141 geom->getSensorXAlphaRefPlane(cls.getSensorID(), chipX, chipAlpha);
1142 if (!trkWorkITS.rotate(chipAlpha) || !propagator->propagateToX(trkWorkITS, chipX, mBz, mParams->maxSnp, mParams->maxStep, mMatCorr)) {
1143 LOGP(debug, "Failed final propagation to ITS X={} alpha={}", chipX, chipAlpha);
1144 stopPropagation = true;
1145 break;
1146 }
1147 float tgPhi = trkWorkITS.getSnp() / std::sqrt((1.f - trkWorkITS.getSnp()) * (1.f + trkWorkITS.getSnp()));
1148 auto dy = cls.getY() - trkWorkITS.getY();
1149 auto dz = cls.getZ() - trkWorkITS.getZ();
1150 if ((std::abs(dy) < param::MaxResid) && (std::abs(dz) < param::MaxResid) && (std::abs(trkWorkITS.getY()) < param::MaxY) && (std::abs(trkWorkITS.getZ()) < param::MaxZ) && (std::abs(tgPhi) < param::MaxTgSlp)) {
1151 mClRes.emplace_back(dy, dz, tgPhi, trkWorkITS.getY(), trkWorkITS.getZ(), 180 + geom->getLayer(cls.getSensorID()), -1, cls.getSensorID());
1152 mDetInfoRes.emplace_back(); // empty placeholder
1153 trackData.nExtDetResid++;
1154 }
1155 }
1156 if (!stopPropagation) { // add residual to PV
1157 const auto& pv = mRecoCont->getPrimaryVertices()[mTrackPVID[iSeed]];
1158 o2::math_utils::Point3D<float> vtx{pv.getX(), pv.getY(), pv.getZ()};
1159 if (!propagator->propagateToDCA(vtx, trkWorkITS, mBz, mParams->maxStep, mMatCorr)) {
1160 LOGP(debug, "Failed propagation to DCA to PV ({} {} {}), {}", pv.getX(), pv.getY(), pv.getZ(), trkWorkITS.asString());
1161 stopPropagation = true;
1162 break;
1163 }
1164 // rotate PV to the track frame
1165 float sn, cs, alpha = trkWorkITS.getAlpha();
1166 math_utils::detail::bringToPMPi(alpha);
1167 math_utils::detail::sincos<float>(alpha, sn, cs);
1168 float xv = vtx.X() * cs + vtx.Y() * sn, yv = -vtx.X() * sn + vtx.Y() * cs, zv = vtx.Z();
1169 auto dy = yv - trkWorkITS.getY();
1170 auto dz = zv - trkWorkITS.getZ();
1171 if ((std::abs(dy) < param::MaxResid) && (std::abs(dz) < param::MaxResid) && (std::abs(trkWorkITS.getY()) < param::MaxY) && (std::abs(trkWorkITS.getZ()) < param::MaxZ) && abs(xv) < param::MaxVtxX) {
1172 short compXV = static_cast<short>(xv * 0x7fff / param::MaxVtxX);
1173 mClRes.emplace_back(dy, dz, alpha / TMath::Pi(), trkWorkITS.getY(), trkWorkITS.getZ(), 190, -1, compXV);
1174 if (!gidTableFull[GTrackID::ITSTPC].isIndexSet()) {
1175 LOGP(fatal, "ITS-TPC seed index is not set for TOF track");
1176 }
1177 float tdif = pv.getTimeStamp().getTimeStamp() - mRecoCont->getTPCITSTrack(gidTableFull[GTrackID::ITSTPC]).getTimeMUS().getTimeStamp();
1178 mDetInfoRes.emplace_back().setPV(tdif); // time in \mus wrt seeding ITS-TPC track
1179 trackData.nExtDetResid++;
1180 }
1181 }
1182 break;
1183 }
1184 }
1185 mTrackData.push_back(std::move(trackData));
1186 mGIDsSuccess.push_back(mGIDs[iSeed]);
1187 mTrackDataCompact.emplace_back(trackData.clIdx.getFirstEntry(), trackData.multStack, nClValidated, mGIDs[iSeed].getSource(), trackData.nExtDetResid);
1188 if (mDumpTrackPoints) {
1189 (*trackDataExtended).clIdx.setEntries(nClValidated);
1190 (*trackDataExtended).nExtDetResid = trackData.nExtDetResid;
1191 mTrackDataExtended.push_back(std::move(*trackDataExtended));
1192 }
1193 }
1194 if (mParams->writeUnfiltered) {
1195 TrackData trkDataTmp = trackData;
1196 trkDataTmp.clIdx.setFirstEntry(mClResUnfiltered.size());
1197 trkDataTmp.clIdx.setEntries(clusterResiduals.size());
1198 mTrackDataUnfiltered.push_back(std::move(trkDataTmp));
1199 mClResUnfiltered.insert(mClResUnfiltered.end(), clusterResiduals.begin(), clusterResiduals.end());
1200 }
1201}
1202
1203bool TrackInterpolation::validateTrack(const TrackData& trk, TrackParams& params, const std::vector<TPCClusterResiduals>& clsRes) const
1204{
1205 if (clsRes.size() < mParams->minNCl) {
1206 // no enough clusters for this track to be considered
1207 LOG(debug) << "Skipping track with too few clusters: " << clsRes.size();
1208 return false;
1209 }
1210
1211 bool resHelix = compareToHelix(trk, params, clsRes);
1212 if (!resHelix) {
1213 LOG(debug) << "Skipping track too far from helix approximation";
1214 return false;
1215 }
1216 if (fabsf(mBz) > 0.01 && fabsf(params.qpt) > mParams->maxQ2Pt) {
1217 LOG(debug) << "Skipping track with too high q/pT: " << params.qpt;
1218 return false;
1219 }
1220 if (!outlierFiltering(trk, params, clsRes)) {
1221 return false;
1222 }
1223 return true;
1224}
1225
1226bool TrackInterpolation::compareToHelix(const TrackData& trk, TrackParams& params, const std::vector<TPCClusterResiduals>& clsRes) const
1227{
1228 std::array<float, param::NPadRows> residHelixY;
1229 std::array<float, param::NPadRows> residHelixZ;
1230
1231 std::array<float, param::NPadRows> xLab;
1232 std::array<float, param::NPadRows> yLab;
1233 std::array<float, param::NPadRows> sPath;
1234
1235 float curvature = fabsf(trk.par.getQ2Pt() * mBz * o2::constants::physics::LightSpeedCm2S * 1e-14f);
1236 int secFirst = clsRes[0].sec;
1237 float phiSect = (secFirst + .5f) * o2::constants::math::SectorSpanRad;
1238 float snPhi = sin(phiSect);
1239 float csPhi = cos(phiSect);
1240 sPath[0] = 0.f;
1241
1242 int iRow = 0;
1243 int nCl = clsRes.size();
1244 for (unsigned int iP = 0; iP < nCl; ++iP) {
1245 iRow += clsRes[iP].dRow;
1246 float yTrk = clsRes[iP].y;
1247 // LOGF(info, "iRow(%i), yTrk(%f)", iRow, yTrk);
1248 xLab[iP] = param::RowX[iRow];
1249 if (clsRes[iP].sec != secFirst) {
1250 float phiSectCurrent = (clsRes[iP].sec + .5f) * o2::constants::math::SectorSpanRad;
1251 float cs = cos(phiSectCurrent - phiSect);
1252 float sn = sin(phiSectCurrent - phiSect);
1253 xLab[iP] = param::RowX[iRow] * cs - yTrk * sn;
1254 yLab[iP] = yTrk * cs + param::RowX[iRow] * sn;
1255 } else {
1256 xLab[iP] = param::RowX[iRow];
1257 yLab[iP] = yTrk;
1258 }
1259 // this is needed only later, but we retrieve it already now to save another loop
1260 params.zTrk[iP] = clsRes[iP].z;
1261 params.xTrk[iP] = param::RowX[iRow];
1262 params.dy[iP] = clsRes[iP].dy;
1263 params.dz[iP] = clsRes[iP].dz;
1264 // done retrieving values for later
1265 if (iP > 0) {
1266 float dx = xLab[iP] - xLab[iP - 1];
1267 float dy = yLab[iP] - yLab[iP - 1];
1268 float ds2 = dx * dx + dy * dy;
1269 float ds = sqrt(ds2); // circular path (linear approximation)
1270 // if the curvature of the track or the (approximated) chord length is too large the more exact formula is used:
1271 // chord length = 2r * asin(ds/(2r))
1272 // using the first two terms of the tailer expansion for asin(x) ~ x + x^3 / 6
1273 if (ds * curvature > 0.05) {
1274 ds *= (1.f + ds2 * curvature * curvature / 24.f);
1275 }
1276 sPath[iP] = sPath[iP - 1] + ds;
1277 }
1278 }
1279 if (fabsf(mBz) < 0.01) {
1280 // for B=0 we don't need to try a circular fit...
1281 return true;
1282 }
1283 float xcSec = 0.f;
1284 float ycSec = 0.f;
1285 float r = 0.f;
1286 TrackResiduals::fitCircle(nCl, xLab, yLab, xcSec, ycSec, r, residHelixY);
1287 // LOGF(info, "Done with circle fit. nCl(%i), xcSec(%f), ycSec(%f), r(%f).", nCl, xcSec, ycSec, r);
1288 /*
1289 for (int i=0; i<nCl; ++i) {
1290 LOGF(info, "i(%i), xLab(%f), yLab(%f).", i, xLab[i], yLab[i]);
1291 }
1292 */
1293 // determine curvature
1294 float phiI = TMath::ATan2(yLab[0], xLab[0]);
1295 float phiF = TMath::ATan2(yLab[nCl - 1], xLab[nCl - 1]);
1296 if (phiI < 0) {
1298 }
1299 if (phiF < 0) {
1301 }
1302 float dPhi = phiF - phiI;
1303 float curvSign = -1.f;
1304 if (dPhi > 0) {
1305 if (dPhi < o2::constants::math::PI) {
1306 curvSign = 1.f;
1307 }
1308 } else if (dPhi < -o2::constants::math::PI) {
1309 curvSign = 1.f;
1310 }
1311 params.qpt = std::copysign(1.f / (r * mBz * o2::constants::physics::LightSpeedCm2S * 1e-14f), curvSign);
1312
1313 // calculate circle coordinates in the lab frame
1314 float xc = xcSec * csPhi - ycSec * snPhi;
1315 float yc = xcSec * snPhi + ycSec * csPhi;
1316
1317 std::array<float, 2> pol1Z;
1318 TrackResiduals::fitPoly1(nCl, sPath, params.zTrk, pol1Z);
1319
1320 params.tgl = pol1Z[0];
1321
1322 // max deviations in both directions from helix fit in y and z
1323 float hMinY = 1e9f;
1324 float hMaxY = -1e9f;
1325 float hMinZ = 1e9f;
1326 float hMaxZ = -1e9f;
1327 // extract residuals in Z and fill track slopes in sector frame
1328 int secCurr = secFirst;
1329 iRow = 0;
1330 for (unsigned int iCl = 0; iCl < nCl; ++iCl) {
1331 iRow += clsRes[iCl].dRow;
1332 float resZ = params.zTrk[iCl] - (pol1Z[1] + sPath[iCl] * pol1Z[0]);
1333 residHelixZ[iCl] = resZ;
1334 if (resZ < hMinZ) {
1335 hMinZ = resZ;
1336 }
1337 if (resZ > hMaxZ) {
1338 hMaxZ = resZ;
1339 }
1340 if (residHelixY[iCl] < hMinY) {
1341 hMinY = residHelixY[iCl];
1342 }
1343 if (residHelixY[iCl] > hMaxY) {
1344 hMaxY = residHelixY[iCl];
1345 }
1346 int sec = clsRes[iCl].sec;
1347 if (sec != secCurr) {
1348 secCurr = sec;
1349 phiSect = (.5f + sec) * o2::constants::math::SectorSpanRad;
1350 snPhi = sin(phiSect);
1351 csPhi = cos(phiSect);
1352 xcSec = xc * csPhi + yc * snPhi; // recalculate circle center in the sector frame
1353 }
1354 float cstalp = (param::RowX[iRow] - xcSec) / r;
1355 if (fabsf(cstalp) > 1.f - sFloatEps) {
1356 // track cannot reach this pad row
1357 cstalp = std::copysign(1.f - sFloatEps, cstalp);
1358 }
1359 params.tglArr[iCl] = cstalp / sqrt((1 - cstalp) * (1 + cstalp)); // 1 / tan(acos(cstalp)) = cstalp / sqrt(1 - cstalp^2)
1360
1361 // In B+ the slope of q- should increase with x. Just look on q * B
1362 if (params.qpt * mBz > 0) {
1363 params.tglArr[iCl] *= -1.f;
1364 }
1365 }
1366 // LOGF(info, "CompareToHelix: hMaxY(%f), hMinY(%f), hMaxZ(%f), hMinZ(%f). Max deviation allowed: y(%.2f), z(%.2f)", hMaxY, hMinY, hMaxZ, hMinZ, mParams->maxDevHelixY, mParams->maxDevHelixZ);
1367 // LOGF(info, "New pt/Q (%f), old pt/Q (%f)", 1./params.qpt, 1./trk.qPt);
1368 return fabsf(hMaxY - hMinY) < mParams->maxDevHelixY && fabsf(hMaxZ - hMinZ) < mParams->maxDevHelixZ;
1369}
1370
1371bool TrackInterpolation::outlierFiltering(const TrackData& trk, TrackParams& params, const std::vector<TPCClusterResiduals>& clsRes) const
1372{
1373 if (clsRes.size() < mParams->nMALong) {
1374 LOG(debug) << "Skipping track with too few clusters for long moving average: " << clsRes.size();
1375 return false;
1376 }
1377 float rmsLong = checkResiduals(trk, params, clsRes);
1378 if (static_cast<float>(params.flagRej.count()) / clsRes.size() > mParams->maxRejFrac) {
1379 LOGP(debug, "Skipping track with too many clusters rejected: {} out of {}", params.flagRej.count(), clsRes.size());
1380 return false;
1381 }
1382 if (rmsLong > mParams->maxRMSLong) {
1383 LOG(debug) << "Skipping track with too large RMS: " << rmsLong;
1384 return false;
1385 }
1386 return true;
1387}
1388
1389float TrackInterpolation::checkResiduals(const TrackData& trk, TrackParams& params, const std::vector<TPCClusterResiduals>& clsRes) const
1390{
1391 float rmsLong = 0.f;
1392
1393 int nCl = clsRes.size();
1394 int iClFirst = 0;
1395 int iClLast = nCl - 1;
1396 int secStart = clsRes[0].sec;
1397
1398 // arrays with differences / abs(differences) of points to their neighbourhood, initialized to zero
1399 std::array<float, param::NPadRows> yDiffLL{};
1400 std::array<float, param::NPadRows> zDiffLL{};
1401 std::array<float, param::NPadRows> absDevY{};
1402 std::array<float, param::NPadRows> absDevZ{};
1403
1404 for (unsigned int iCl = 0; iCl < nCl; ++iCl) {
1405 if (iCl < iClLast && clsRes[iCl].sec == secStart) {
1406 continue;
1407 }
1408 // sector changed or last cluster reached
1409 // now run estimators for all points in the same sector
1410 int nClSec = iCl - iClFirst;
1411 if (iCl == iClLast) {
1412 ++nClSec;
1413 }
1414 diffToLocLine(nClSec, iClFirst, params.xTrk, params.dy, yDiffLL);
1415 diffToLocLine(nClSec, iClFirst, params.xTrk, params.dz, zDiffLL);
1416 iClFirst = iCl;
1417 secStart = clsRes[iCl].sec;
1418 }
1419 // store abs deviations
1420 int nAccY = 0;
1421 int nAccZ = 0;
1422 for (int iCl = nCl; iCl--;) {
1423 if (fabsf(yDiffLL[iCl]) > param::sEps) {
1424 absDevY[nAccY++] = fabsf(yDiffLL[iCl]);
1425 }
1426 if (fabsf(zDiffLL[iCl]) > param::sEps) {
1427 absDevZ[nAccZ++] = fabsf(zDiffLL[iCl]);
1428 }
1429 }
1430 if (nAccY < mParams->minNumberOfAcceptedResiduals || nAccZ < mParams->minNumberOfAcceptedResiduals) {
1431 // mask all clusters
1432 LOGP(debug, "Accepted {} clusters for dY {} clusters for dZ, but required at least {} for both", nAccY, nAccZ, mParams->minNumberOfAcceptedResiduals);
1433 params.flagRej.set();
1434 return 0.f;
1435 }
1436 // estimate rms on 90% of the smallest deviations
1437 int nKeepY = static_cast<int>(.9 * nAccY);
1438 int nKeepZ = static_cast<int>(.9 * nAccZ);
1439 std::nth_element(absDevY.begin(), absDevY.begin() + nKeepY, absDevY.begin() + nAccY);
1440 std::nth_element(absDevZ.begin(), absDevZ.begin() + nKeepZ, absDevZ.begin() + nAccZ);
1441 float rmsYkeep = 0.f;
1442 float rmsZkeep = 0.f;
1443 for (int i = nKeepY; i--;) {
1444 rmsYkeep += absDevY[i] * absDevY[i];
1445 }
1446 for (int i = nKeepZ; i--;) {
1447 rmsZkeep += absDevZ[i] * absDevZ[i];
1448 }
1449 rmsYkeep = std::sqrt(rmsYkeep / nKeepY);
1450 rmsZkeep = std::sqrt(rmsZkeep / nKeepZ);
1451 if (rmsYkeep < param::sEps || rmsZkeep < param::sEps) {
1452 LOG(warning) << "Too small RMS: " << rmsYkeep << "(y), " << rmsZkeep << "(z).";
1453 params.flagRej.set();
1454 return 0.f;
1455 }
1456 float rmsYkeepI = 1.f / rmsYkeep;
1457 float rmsZkeepI = 1.f / rmsZkeep;
1458 int nAcc = 0;
1459 std::array<float, param::NPadRows> yAcc;
1460 std::array<float, param::NPadRows> yDiffLong;
1461 for (int iCl = 0; iCl < nCl; ++iCl) {
1462 yDiffLL[iCl] *= rmsYkeepI;
1463 zDiffLL[iCl] *= rmsZkeepI;
1464 if (yDiffLL[iCl] * yDiffLL[iCl] + zDiffLL[iCl] * zDiffLL[iCl] > mParams->maxStdDevMA) {
1465 params.flagRej.set(iCl);
1466 } else {
1467 yAcc[nAcc++] = params.dy[iCl];
1468 }
1469 }
1470 if (nAcc > mParams->nMALong) {
1471 diffToMA(nAcc, yAcc, yDiffLong);
1472 float average = 0.f;
1473 float rms = 0.f;
1474 for (int i = 0; i < nAcc; ++i) {
1475 average += yDiffLong[i];
1476 rms += yDiffLong[i] * yDiffLong[i];
1477 }
1478 average /= nAcc;
1479 rmsLong = rms / nAcc - average * average;
1480 rmsLong = (rmsLong > 0) ? std::sqrt(rmsLong) : 0.f;
1481 }
1482 return rmsLong;
1483}
1484
1485void TrackInterpolation::diffToLocLine(const int np, int idxOffset, const std::array<float, param::NPadRows>& x, const std::array<float, param::NPadRows>& y, std::array<float, param::NPadRows>& diffY) const
1486{
1487 // Calculate the difference between the points and the linear extrapolations from the neighbourhood.
1488 // Nothing more than multiple 1-d fits at once. Instead of building 4 sums (x, x^2, y, xy), 4 * nPoints sums are calculated at once
1489 // compare to TrackResiduals::fitPoly1() method
1490
1491 // adding one entry to the vectors saves an additional if statement when calculating the cumulants
1492 std::vector<float> sumX1vec(np + 1);
1493 std::vector<float> sumX2vec(np + 1);
1494 std::vector<float> sumY1vec(np + 1);
1495 std::vector<float> sumXYvec(np + 1);
1496 auto sumX1 = &(sumX1vec[1]);
1497 auto sumX2 = &(sumX2vec[1]);
1498 auto sumY1 = &(sumY1vec[1]);
1499 auto sumXY = &(sumXYvec[1]);
1500
1501 // accumulate sums for all points
1502 for (int iCl = 0; iCl < np; ++iCl) {
1503 int idx = iCl + idxOffset;
1504 sumX1[iCl] = sumX1[iCl - 1] + x[idx];
1505 sumX2[iCl] = sumX2[iCl - 1] + x[idx] * x[idx];
1506 sumY1[iCl] = sumY1[iCl - 1] + y[idx];
1507 sumXY[iCl] = sumXY[iCl - 1] + x[idx] * y[idx];
1508 }
1509
1510 for (int iCl = 0; iCl < np; ++iCl) {
1511 int iClLeft = iCl - mParams->nMAShort;
1512 int iClRight = iCl + mParams->nMAShort;
1513 if (iClLeft < 0) {
1514 iClLeft = 0;
1515 }
1516 if (iClRight >= np) {
1517 iClRight = np - 1;
1518 }
1519 int nPoints = iClRight - iClLeft;
1520 if (nPoints < mParams->nMAShort) {
1521 continue;
1522 }
1523 float nPointsInv = 1.f / nPoints;
1524 int iClLeftP = iClLeft - 1;
1525 int iClCurrP = iCl - 1;
1526 // extract sum from iClLeft to iClRight from cumulants, excluding iCl from the fit
1527 float sX1 = sumX1[iClRight] - sumX1[iClLeftP] - (sumX1[iCl] - sumX1[iClCurrP]);
1528 float sX2 = sumX2[iClRight] - sumX2[iClLeftP] - (sumX2[iCl] - sumX2[iClCurrP]);
1529 float sY1 = sumY1[iClRight] - sumY1[iClLeftP] - (sumY1[iCl] - sumY1[iClCurrP]);
1530 float sXY = sumXY[iClRight] - sumXY[iClLeftP] - (sumXY[iCl] - sumXY[iClCurrP]);
1531 float det = sX2 - nPointsInv * sX1 * sX1;
1532 if (fabsf(det) < 1e-12f) {
1533 continue;
1534 }
1535 float slope = (sXY - nPointsInv * sX1 * sY1) / det;
1536 float offset = nPointsInv * sY1 - nPointsInv * slope * sX1;
1537 diffY[iCl + idxOffset] = y[iCl + idxOffset] - slope * x[iCl + idxOffset] - offset;
1538 }
1539}
1540
1541void TrackInterpolation::diffToMA(const int np, const std::array<float, param::NPadRows>& y, std::array<float, param::NPadRows>& diffMA) const
1542{
1543 // Calculate
1544 std::vector<float> sumVec(np + 1);
1545 auto sum = &(sumVec[1]);
1546 for (int i = 0; i < np; ++i) {
1547 sum[i] = sum[i - 1] + y[i];
1548 }
1549 for (int i = 0; i < np; ++i) {
1550 diffMA[i] = 0;
1551 int iLeft = i - mParams->nMALong;
1552 int iRight = i + mParams->nMALong;
1553 if (iLeft < 0) {
1554 iLeft = 0;
1555 }
1556 if (iRight >= np) {
1557 iRight = np - 1;
1558 }
1559 int nPoints = iRight - iLeft;
1560 if (nPoints < mParams->nMALong) {
1561 // this cannot happen, since at least mParams->nMALong points are required as neighbours for this function to be called
1562 continue;
1563 }
1564 float movingAverage = (sum[iRight] - sum[iLeft - 1] - (sum[i] - sum[i - 1])) / nPoints;
1565 diffMA[i] = y[i] - movingAverage;
1566 }
1567}
1568
1570{
1571 mTrackData.clear();
1572 mTrackDataCompact.clear();
1573 mTrackDataExtended.clear();
1574 mClRes.clear();
1575 mDetInfoRes.clear();
1576 mTrackDataUnfiltered.clear();
1577 mClResUnfiltered.clear();
1578 mGIDsSuccess.clear();
1579 for (auto& vec : mTrackIndices) {
1580 vec.clear();
1581 }
1582 mGIDs.clear();
1583 mGIDtables.clear();
1584 mTrackTimes.clear();
1585 mSeeds.clear();
1586 mITSRefitSeedID.clear();
1587 mTrackPVID.clear();
1588}
1589
1590//______________________________________________
1592{
1593 // Attention! For the refit we are using reference VDrift and TDriftOffest rather than high-rate calibrated, since we want to have fixed reference over the run
1594 if (v.refVDrift != mTPCVDriftRef) {
1595 mTPCVDriftRef = v.refVDrift;
1596 mTPCDriftTimeOffsetRef = v.refTimeOffset;
1597 LOGP(info, "Imposing reference VDrift={}/TDrift={} for TPC residuals extraction", mTPCVDriftRef, mTPCDriftTimeOffsetRef);
1598 o2::tpc::TPCFastTransformHelperO2::instance()->updateCalibration(*mFastTransform, 0, 1.0, mTPCVDriftRef, mTPCDriftTimeOffsetRef);
1599 }
1600}
1601
1602//______________________________________________
1604{
1605 // refit ITS track outwards taking PID (unless already refitted) from the seed and reassign to the seed
1606 auto& seed = mSeeds[seedID];
1607 int refitID = mITSRefitSeedID[gid.getIndex()];
1608 if (refitID >= 0) { // track was already refitted
1609 if (mSeeds[refitID].getPID() == seed.getPID()) {
1610 seed = mSeeds[refitID];
1611 }
1612 return true;
1613 }
1614 const auto& trkITS = mRecoCont->getITSTrack(gid);
1615 // fetch clusters
1616 auto nCl = trkITS.getNumberOfClusters();
1617 auto clEntry = trkITS.getFirstClusterEntry();
1618 o2::track::TrackParCov track(trkITS); // start from the inner param
1619 track.resetCovariance();
1620 track.setCov(track.getQ2Pt() * track.getQ2Pt() * track.getCov()[o2::track::CovLabels::kSigQ2Pt2], o2::track::CovLabels::kSigQ2Pt2);
1621 track.setPID(seed.getPID());
1622 o2::track::TrackPar refLin(track); // and use it also as linearization reference
1623 auto geom = o2::its::GeometryTGeo::Instance();
1624 auto prop = o2::base::Propagator::Instance();
1625 for (int iCl = nCl - 1; iCl >= 0; iCl--) { // clusters are stored from outer to inner layers
1626 const auto& cls = mITSClustersArray[mITSTrackClusIdx[clEntry + iCl]];
1627 int chip = cls.getSensorID();
1628 float chipX, chipAlpha;
1629 geom->getSensorXAlphaRefPlane(cls.getSensorID(), chipX, chipAlpha);
1630 if (!track.rotate(chipAlpha, refLin, mBz)) {
1631 LOGP(debug, "failed to rotate ITS tracks to alpha={} for the refit: {}", chipAlpha, track.asString());
1632 return false;
1633 }
1634 if (!prop->propagateToX(track, refLin, cls.getX(), mBz, o2::base::PropagatorImpl<float>::MAX_SIN_PHI, o2::base::PropagatorImpl<float>::MAX_STEP, o2::base::PropagatorF::MatCorrType::USEMatCorrLUT)) {
1635 LOGP(debug, "failed to propagate ITS tracks to X={}: {}", cls.getX(), track.asString());
1636 return false;
1637 }
1638 std::array<float, 2> posTF{cls.getY(), cls.getZ()};
1639 std::array<float, 3> covTF{cls.getSigmaY2(), cls.getSigmaYZ(), cls.getSigmaZ2()};
1640 if (!track.update(posTF, covTF)) {
1641 LOGP(debug, "failed to update ITS tracks by cluster ({},{})/({},{},{})", track.asString(), cls.getY(), cls.getZ(), cls.getSigmaY2(), cls.getSigmaYZ(), cls.getSigmaZ2());
1642 return false;
1643 }
1644 if (mParams->shiftRefToCluster) {
1645 refLin.setY(posTF[0]);
1646 refLin.setZ(posTF[1]);
1647 }
1648 }
1649 seed = track;
1650 // memorize that this ITS track was already refitted
1651 mITSRefitSeedID[gid.getIndex()] = seedID;
1652 return true;
1653}
std::string asString(TDataMember const &dm, char *pointer)
Global TRD definitions and constants.
std::ostringstream debug
int32_t i
Definition of the GeometryTGeo class.
Definition of the parameter class for the detector electronics.
uint16_t pos
Definition RawData.h:3
uint16_t slope
Definition RawData.h:1
uint32_t res
Definition RawData.h:0
Definition of the TrackInterpolation class.
Definition of the TrackResiduals class.
calibration data from laser track calibration
Referenc on track indices contributing to the vertex, with possibility chose tracks from specific sou...
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:178
static BasicCCDBManager & instance()
T * getSpecific(std::string const &path, long timestamp=-1, MD metaData=MD(), std::map< std::string, std::string > *headers=nullptr)
retrieve an object of type T from CCDB as stored under path, timestamp and metaData
static mask_t getSourcesMask(const std::string_view srcList)
static std::string getSourcesNames(mask_t srcm)
Static class with identifiers, bitmasks and names for ALICE detectors.
Definition DetID.h:58
static constexpr ID TRD
Definition DetID.h:65
static constexpr ID TOF
Definition DetID.h:66
static size_t fillOccupancyMapGetSize(uint32_t nHbfPerTf, const GPUParam *param=nullptr)
static void paramUseExternalOccupancyMap(GPUParam *param, uint32_t nHbfPerTf, const uint32_t *occupancymap, int32_t occupancyMapSize)
static std::shared_ptr< GPUParam > getFullParamShared(float solenoidBz, uint32_t nHbfPerTf=0, std::unique_ptr< GPUO2InterfaceConfiguration > *pConfiguration=nullptr, std::unique_ptr< GPUSettingsO2 > *pO2Settings=nullptr, bool *autoMaxTimeBin=nullptr)
void init(float bz, const GPUSettingsRec *rec=nullptr)
Load parameterization for given magnetic field.
float getSensorRefAlpha(int isn) const
static GeometryTGeo * Instance()
float getSensorRefX(int isn) const
int getFirstClusterEntry() const
Definition TrackITS.h:66
static constexpr Double_t BC_TIME_INPS
Definition Geo.h:105
static constexpr Double_t BC_TIME_INPS_INV
Definition Geo.h:106
static void rotateToSector(Float_t *xyz, Int_t isector)
Definition Geo.cxx:1029
static void getPos(Int_t *det, Float_t *pos)
Definition Geo.cxx:491
static constexpr Int_t NPADSXSECTOR
Definition Geo.h:120
static void getVolumeIndices(Int_t index, Int_t *detId)
Definition Geo.cxx:543
static void alignedToNominalSector(Float_t *xyz, Int_t isector)
Definition Geo.cxx:1049
int updateCalibration(TPCFastTransform &transform, Long_t TimeStamp, float vDriftFactor=1.f, float vDriftRef=0.f, float driftTimeOffset=0.f)
Updates the transformation with the new time stamp.
static TPCFastTransformHelperO2 * instance()
Singleton.
void prepareInputTrackSample(const o2::globaltracking::RecoContainer &inp)
Prepare input track sample (not relying on CreateTracksVariadic functionality)
int processTRDLayer(const o2::trd::TrackTRD &trkTRD, int iLayer, o2::track::TrackParCov &trkWork, std::array< float, 2 > *trkltTRDYZ=nullptr, std::array< float, 3 > *trkltTRDCov=nullptr, TrackData *trkData=nullptr, o2::trd::Tracklet64 *trk64=nullptr, o2::trd::CalibratedTracklet *trkCalib=nullptr)
bool isInputTrackAccepted(const o2::dataformats::GlobalTrackID &gid, const o2::globaltracking::RecoContainer::GlobalIDSet &gidTable, const o2::dataformats::PrimaryVertex &pv) const
Check if input track passes configured cuts.
void diffToMA(const int np, const std::array< float, param::NPadRows > &y, std::array< float, param::NPadRows > &diffMA) const
For a given set of points, calculate their deviation from the moving average (build from the neighbou...
void diffToLocLine(const int np, int idxOffset, const std::array< float, param::NPadRows > &x, const std::array< float, param::NPadRows > &y, std::array< float, param::NPadRows > &diffY) const
For a given set of points, calculate the differences from each point to the fitted lines from all oth...
o2::dataformats::GlobalTrackID::Source findValidSource(const o2::dataformats::GlobalTrackID::mask_t mask, const o2::dataformats::GlobalTrackID::Source src) const
For given vertex track source which is not in mSourcesConfigured find the seeding source which is ena...
bool refITSTrack(o2::dataformats::GlobalTrackID, int iSeed)
float checkResiduals(const TrackData &trk, TrackParams &params, const std::vector< TPCClusterResiduals > &clsRes) const
bool compareToHelix(const TrackData &trk, TrackParams &params, const std::vector< TPCClusterResiduals > &clsRes) const
void process()
Main processing function.
@ ExtIn
extrapolation inwards of TRD/TOF track
@ Int
interpolation (mean positions of both extrapolations)
@ ExtOut
extrapolation outwards of ITS track
void setTPCVDrift(const o2::tpc::VDriftCorrFact &v)
bool isTrackSelected(const o2::track::TrackParCov &trk) const
Given the defined downsampling factor tsalisThreshold check if track is selected.
void reset()
Reset cache and output vectors.
bool outlierFiltering(const TrackData &trk, TrackParams &params, const std::vector< TPCClusterResiduals > &clsRes) const
void init(o2::dataformats::GlobalTrackID::mask_t src, o2::dataformats::GlobalTrackID::mask_t srcMap)
Initialize everything, set the requested track sources.
bool validateTrack(const TrackData &trk, TrackParams &params, const std::vector< TPCClusterResiduals > &clsRes) const
static bool fitPoly1(int nCl, std::array< float, param::NPadRows > &x, std::array< float, param::NPadRows > &y, std::array< float, 2 > &res)
static void fitCircle(int nCl, std::array< float, param::NPadRows > &x, std::array< float, param::NPadRows > &y, float &xc, float &yc, float &r, std::array< float, param::NPadRows > &residHelixY)
static Geometry * instance()
Definition Geometry.h:33
GLfloat GLfloat GLfloat alpha
Definition glcorearb.h:279
GLint GLenum GLint x
Definition glcorearb.h:403
GLenum src
Definition glcorearb.h:1767
GLsizeiptr size
Definition glcorearb.h:659
GLuint GLuint end
Definition glcorearb.h:469
const GLdouble * v
Definition glcorearb.h:832
GLuint GLuint GLfloat weight
Definition glcorearb.h:5477
GLenum const GLfloat * params
Definition glcorearb.h:272
GLintptr offset
Definition glcorearb.h:660
GLboolean GLboolean g
Definition glcorearb.h:1233
GLboolean r
Definition glcorearb.h:1233
GLint GLuint mask
Definition glcorearb.h:291
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
constexpr float TwoPI
constexpr float SectorSpanRad
constexpr float PI
constexpr float LightSpeedCm2S
constexpr double MassPionCharged
void convertCompactClusters(gsl::span< const itsmft::CompClusterExt > clusters, gsl::span< const unsigned char >::iterator &pattIt, std::vector< o2::BaseCluster< float > > &output, const itsmft::TopologyDictionary *dict)
convert compact clusters to 3D spacepoints
Definition IOUtils.cxx:35
int angle2Sector(float phi)
Definition Utils.h:183
float sector2Angle(int sect)
Definition Utils.h:193
constexpr int MAXGLOBALPADROW
Definition Constants.h:34
Global TPC definitions and constants.
Definition SimTraits.h:168
Enum< T >::Iterator begin(Enum< T >)
Definition Defs.h:156
DataT sum(const Vector< DataT, N > &a)
Definition Vector.h:107
TrackParCovF TrackParCov
Definition Track.h:33
constexpr int NLAYER
the number of layers
Definition Constants.h:27
constexpr int NSTACK
the number of stacks per sector
Definition Constants.h:26
GlobalIDSet getSingleDetectorRefs(GTrackID gidx) const
const o2::tpc::TrackTPC & getTPCTrack(GTrackID id) const
gsl::span< const o2::trd::CalibratedTracklet > getTRDCalibratedTracklets() const
const o2::tpc::ClusterNativeAccess & getTPCClusters() const
const o2::track::TrackParCov & getTrackParam(GTrackID gidx) const
const o2::dataformats::TrackTPCITS & getTPCITSTrack(GTrackID gid) const
const o2::its::TrackITS & getITSTrack(GTrackID gid) const
auto getITSTPCTRDTrack(GTrackID id) const
gsl::span< const unsigned int > occupancyMapTPC
externally set TPC clusters occupancy map
std::array< GTrackID, GTrackID::NSources > GlobalIDSet
gsl::span< const o2::trd::Tracklet64 > getTRDTracklets() const
const o2::dataformats::MatchInfoTOF & getTOFMatch(GTrackID id) const
static constexpr int L2G
Definition Cartesian.h:54
static constexpr int T2L
Definition Cartesian.h:55
static bool downsampleTsallisCharged(float pt, float factorPt, float sqrts, float &weight, float rnd, float mass=0.13957)
Definition Tsallis.cxx:31
float maxRejFrac
if the fraction of rejected clusters of a track is higher, the full track is invalidated
float maxDevHelixZ
max deviation in Z for clusters wrt helix fit
int nMALong
number of points to be used for moving average (long range)
float minPtNoOuterPoint
minimum pt for ITS-TPC tracks to be considered for extrapolation
int minTRDNTrklts
min number of TRD space points
bool refitITS
refit ITS tracks with PID attached to the seed
int nMAShort
number of points to be used for estimation of distance from local line (short range)
float maxDevHelixY
max deviation in Y for clusters wrt helix fit
float maxStdDevMA
max cluster std. deviation (Y^2 + Z^2) wrt moving average to accept
bool skipOutlierFiltering
if set, the outlier filtering will not be applied at all
bool ignoreNonPVContrib
flag if tracks which did not contribute to the PV should be ignored or not
int minITSNClsNoOuterPoint
min number of ITS clusters if no hit in TRD or TOF exists
int minTPCNCls
min number of TPC clusters
float tsalisThreshold
in case the sampling functions returns a value smaller than this the track is discarded (1....
bool shiftRefToCluster
when reftting the ITS track, shift the lin.reference to cluster after every update (better material m...
float maxStep
maximum step for propagation
int minTOFTRDPVContributors
min contributors from TRD or TOF (fast detectors) to consider tracks of this PV
bool enableTrackDownsampling
flag if track sampling shall be enabled or not
float maxQ2Pt
max fitted q/pt for a track to be used for calibration
float sigYZ2TOF
for now assume cluster error for TOF equal for all clusters in both Y and Z
int minNumberOfAcceptedResiduals
min number of accepted residuals for
float maxRMSLong
maximum variance of the cluster residuals wrt moving avarage for a track to be considered
int minITSNCls
min number of ITS clusters
float maxSnp
max snp when propagating tracks
bool writeUnfiltered
if set, all residuals and track parameters will be aggregated and dumped additionally without outlier...
int minNCl
min number of clusters in a track to be used for calibration
int minTPCNClsNoOuterPoint
min number of TPC clusters if no hit in TRD or TOF exists
Structure filled for each track with track quality information and a vector with TPCClusterResiduals.
float chi2TRD
chi2 of TRD track
float dEdxTPC
TPC dEdx information.
unsigned short nTrkltsTRD
number of attached TRD tracklets
unsigned short clAvailTOF
whether or not track seed has a matched TOF cluster, if so, gives the resolution of the T0 in ps
std::array< uint8_t, 4 > multStack
short TRDTrkltSlope[6]
TRD tracklet slope 0x7fff / param::MaxTRDSlope.
unsigned short nClsITS
number of attached ITS clusters
o2::dataformats::RangeReference clIdx
index of first cluster residual and total number of TPC cluster residuals of this track
float chi2TPC
chi2 of TPC track
unsigned short nClsTPC
number of attached TPC clusters
void setMultStack(float v, int stack)
o2::dataformats::GlobalTrackID gid
global track ID for seeding track
float deltaTOF
TOFsignal - T0 - texp(PID), if T0 is available.
uint8_t nExtDetResid
number of external detectors (to TPC) residuals stored, on top of clIdx.getEntries
o2::track::TrackPar par
ITS track at inner TPC radius.
float chi2ITS
chi2 of ITS track
Structure for on-the-fly re-calculated track parameters at the validation stage.
unsigned char sec
TPC sector (0..35)
unsigned char row
TPC pad row.
static void init(long timestamp=-1)
short channel
extra channel info (ITS chip ID, TRD chamber, TOF main pad within the sector)
o2::mch::DsIndex ds
std::vector< o2::ctf::BufferType > vec
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
o2::InteractionRecord ir(0, 0)
std::uniform_int_distribution< unsigned long long > distr
std::random_device rd
std::vector< int > row