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