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