Project
Loading...
Searching...
No Matches
TrackerTraits.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.
15
16#include <algorithm>
17#include <iterator>
18#include <cmath>
19#include <type_traits>
20
21#include <oneapi/tbb/blocked_range.h>
22#include <oneapi/tbb/enumerable_thread_specific.h>
23
25#include "GPUCommonMath.h"
27#include "ITStracking/Cell.h"
35
36namespace o2::its
37{
38
39struct PassMode {
40 using OnePass = std::integral_constant<int, 0>;
41 using TwoPassCount = std::integral_constant<int, 1>;
42 using TwoPassInsert = std::integral_constant<int, 2>;
43};
44
45template <int NLayers>
46void TrackerTraits<NLayers>::computeLayerTracklets(const int iteration, int iVertex)
47{
48 const auto topology = mTimeFrame->getTrackingTopologyView();
49 for (int transitionId = 0; transitionId < topology.nTransitions; ++transitionId) {
50 mTimeFrame->getTracklets()[transitionId].clear();
51 mTimeFrame->getTrackletsLabel(transitionId).clear();
52 std::fill(mTimeFrame->getTrackletsLookupTable()[transitionId].begin(), mTimeFrame->getTrackletsLookupTable()[transitionId].end(), 0);
53 }
54
55 const Vertex diamondVert(mTrkParams[iteration].Diamond, mTrkParams[iteration].DiamondCov, 1, 1.f);
56 gsl::span<const Vertex> diamondSpan(&diamondVert, 1);
57
58 mTaskArena->execute([&] {
59 auto forTracklets = [&](auto Tag, int transitionId, int pivotROF, int base, int& offset) -> int {
60 const auto& transition = topology.getTransition(transitionId);
61 if (!mTimeFrame->getROFMaskView().isROFEnabled(transition.fromLayer, pivotROF)) {
62 return 0;
63 }
64 gsl::span<const Vertex> primaryVertices = mTrkParams[iteration].UseDiamond ? diamondSpan : mTimeFrame->getPrimaryVertices(transition.fromLayer, pivotROF);
65 if (primaryVertices.empty()) {
66 return 0;
67 }
68 const int startVtx = iVertex >= 0 ? iVertex : 0;
69 const int endVtx = iVertex >= 0 ? o2::gpu::CAMath::Min(iVertex + 1, int(primaryVertices.size())) : int(primaryVertices.size());
70 if (endVtx <= startVtx || (iVertex + 1) > primaryVertices.size()) {
71 return 0;
72 }
73
74 const auto& rofOverlap = mTimeFrame->getROFOverlapTableView().getOverlap(transition.fromLayer, transition.toLayer, pivotROF);
75 if (!rofOverlap.getEntries()) {
76 return 0;
77 }
78
79 int localCount = 0;
80 auto& tracklets = mTimeFrame->getTracklets()[transitionId];
81 auto layer0 = mTimeFrame->getClustersOnLayer(pivotROF, transition.fromLayer);
82 if (layer0.empty()) {
83 return 0;
84 }
85
86 const float meanDeltaR = mTrkParams[iteration].LayerRadii[transition.toLayer] - mTrkParams[iteration].LayerRadii[transition.fromLayer];
87 const float phiCut = mTimeFrame->getTransitionPhiCut(transitionId);
88 const float msAngle = mTimeFrame->getTransitionMSAngle(transitionId);
89
90 for (int iCluster = 0; iCluster < int(layer0.size()); ++iCluster) {
91 const Cluster& currentCluster = layer0[iCluster];
92 const int currentSortedIndex = mTimeFrame->getSortedIndex(pivotROF, transition.fromLayer, iCluster);
93 if (mTimeFrame->isClusterUsed(transition.fromLayer, currentCluster.clusterId)) {
94 continue;
95 }
96 const float inverseR0 = 1.f / currentCluster.radius;
97
98 for (int iV = startVtx; iV < endVtx; ++iV) {
99 const auto& pv = primaryVertices[iV];
100 if (!mTimeFrame->getROFVertexLookupTableView().isVertexCompatible(transition.fromLayer, pivotROF, pv)) {
101 continue;
102 }
103 if (pv.isFlagSet(Vertex::Flags::UPCMode) != mTrkParams[iteration].PassFlags[IterationStep::SelectUPCVertices]) {
104 continue;
105 }
106 const float resolution = o2::gpu::CAMath::Sqrt(math_utils::Sq(mTimeFrame->getPositionResolution(transition.fromLayer)) + math_utils::Sq(mTrkParams[iteration].PVres) / float(pv.getNContributors()));
107 const float tanLambda = (currentCluster.zCoordinate - pv.getZ()) * inverseR0;
108 const float zAtRmin = tanLambda * (mTimeFrame->getMinR(transition.toLayer) - currentCluster.radius) + currentCluster.zCoordinate;
109 const float zAtRmax = tanLambda * (mTimeFrame->getMaxR(transition.toLayer) - currentCluster.radius) + currentCluster.zCoordinate;
110 const float sqInvDeltaZ0 = 1.f / (math_utils::Sq(currentCluster.zCoordinate - pv.getZ()) + constants::Tolerance);
111 const float sigmaZ = o2::gpu::CAMath::Sqrt((math_utils::Sq(resolution) * math_utils::Sq(tanLambda) * ((math_utils::Sq(inverseR0) + sqInvDeltaZ0) * math_utils::Sq(meanDeltaR) + 1.f)) + math_utils::Sq(meanDeltaR * msAngle));
112 const auto bins = o2::its::getBinsRect(currentCluster, transition.toLayer, zAtRmin, zAtRmax,
113 sigmaZ * mTrkParams[iteration].NSigmaCut, phiCut,
114 mTimeFrame->getIndexTableUtils());
115 if (bins.x < 0) {
116 continue;
117 }
118 int phiBinsNum = bins.w - bins.y + 1;
119 if (phiBinsNum < 0) {
120 phiBinsNum += mTrkParams[iteration].PhiBins;
121 }
122
123 for (int targetROF = rofOverlap.getFirstEntry(); targetROF < rofOverlap.getEntriesBound(); ++targetROF) {
124 if (!mTimeFrame->getROFMaskView().isROFEnabled(transition.toLayer, targetROF)) {
125 continue;
126 }
127 auto layer1 = mTimeFrame->getClustersOnLayer(targetROF, transition.toLayer);
128 if (layer1.empty()) {
129 continue;
130 }
131 const auto ts = mTimeFrame->getROFOverlapTableView().getTimeStamp(transition.fromLayer, pivotROF, transition.toLayer, targetROF);
132 if (!ts.isCompatible(pv.getTimeStamp())) {
133 continue;
134 }
135 const auto& targetIndexTable = mTimeFrame->getIndexTable(targetROF, transition.toLayer);
136 const int zBinRange = (bins.z - bins.x) + 1;
137 for (int iPhi = 0; iPhi < phiBinsNum; ++iPhi) {
138 const int iPhiBin = (bins.y + iPhi) % mTrkParams[iteration].PhiBins;
139 const int firstBinIdx = mTimeFrame->getIndexTableUtils().getBinIndex(bins.x, iPhiBin);
140 const int maxBinIdx = firstBinIdx + zBinRange;
141 const int firstRow = targetIndexTable[firstBinIdx];
142 const int lastRow = targetIndexTable[maxBinIdx];
143 for (int iNext = firstRow; iNext < lastRow; ++iNext) {
144 if (iNext >= int(layer1.size())) {
145 break;
146 }
147 const Cluster& nextCluster = layer1[iNext];
148 if (mTimeFrame->isClusterUsed(transition.toLayer, nextCluster.clusterId)) {
149 continue;
150 }
151 const float deltaZ = o2::gpu::CAMath::Abs((tanLambda * (nextCluster.radius - currentCluster.radius)) + currentCluster.zCoordinate - nextCluster.zCoordinate);
152
153 if (deltaZ / sigmaZ < mTrkParams[iteration].NSigmaCut &&
154 math_utils::isPhiDifferenceBelow(currentCluster.phi, nextCluster.phi, phiCut)) {
155 const float phi{o2::gpu::CAMath::ATan2(currentCluster.yCoordinate - nextCluster.yCoordinate, currentCluster.xCoordinate - nextCluster.xCoordinate)};
156 const float tanL = (currentCluster.zCoordinate - nextCluster.zCoordinate) / (currentCluster.radius - nextCluster.radius);
157 if constexpr (decltype(Tag)::value == PassMode::OnePass::value) {
158 tracklets.emplace_back(currentSortedIndex, mTimeFrame->getSortedIndex(targetROF, transition.toLayer, iNext), tanL, phi, ts);
159 } else if constexpr (decltype(Tag)::value == PassMode::TwoPassCount::value) {
160 ++localCount;
161 } else if constexpr (decltype(Tag)::value == PassMode::TwoPassInsert::value) {
162 const int idx = base + offset++;
163 tracklets[idx] = Tracklet(currentSortedIndex, mTimeFrame->getSortedIndex(targetROF, transition.toLayer, iNext), tanL, phi, ts);
164 }
165 }
166 }
167 }
168 }
169 }
170 }
171 return localCount;
172 };
173
174 int dummy{0};
175 if (mTaskArena->max_concurrency() <= 1) {
176 for (int transitionId{0}; transitionId < topology.nTransitions; ++transitionId) {
177 const int fromLayer = topology.getTransition(transitionId).fromLayer;
178 const int startROF = 0, endROF = mTimeFrame->getROFOverlapTableView().getLayer(fromLayer).mNROFsTF;
179 for (int pivotROF{startROF}; pivotROF < endROF; ++pivotROF) {
180 forTracklets(PassMode::OnePass{}, transitionId, pivotROF, 0, dummy);
181 }
182 }
183 } else {
184 tbb::parallel_for(0, static_cast<int>(topology.nTransitions), [&](const int transitionId) {
185 const int fromLayer = topology.getTransition(transitionId).fromLayer;
186 const int startROF = 0, endROF = mTimeFrame->getROFOverlapTableView().getLayer(fromLayer).mNROFsTF;
187 bounded_vector<int> perROFCount((endROF - startROF) + 1, mMemoryPool.get());
188 tbb::parallel_for(startROF, endROF, [&](const int pivotROF) {
189 perROFCount[pivotROF - startROF] = forTracklets(PassMode::TwoPassCount{}, transitionId, pivotROF, 0, dummy);
190 });
191 std::exclusive_scan(perROFCount.begin(), perROFCount.end(), perROFCount.begin(), 0);
192 const int nTracklets = perROFCount.back();
193 mTimeFrame->getTracklets()[transitionId].resize(nTracklets);
194 if (nTracklets == 0) {
195 return;
196 }
197 tbb::parallel_for(startROF, endROF, [&](const int pivotROF) {
198 int baseIdx = perROFCount[pivotROF - startROF];
199 if (baseIdx == perROFCount[pivotROF + 1 - startROF]) {
200 return;
201 }
202 int localIdx = 0;
203 forTracklets(PassMode::TwoPassInsert{}, transitionId, pivotROF, baseIdx, localIdx);
204 });
205 });
206 }
207
208 tbb::parallel_for(0, static_cast<int>(topology.nTransitions), [&](const int transitionId) {
210 // duplicates can exist simply since we evaluate per vertex
211 auto& trkl{mTimeFrame->getTracklets()[transitionId]};
212 std::sort(trkl.begin(), trkl.end());
213 trkl.erase(std::unique(trkl.begin(), trkl.end()), trkl.end());
214 trkl.shrink_to_fit();
215 auto& lut{mTimeFrame->getTrackletsLookupTable()[transitionId]};
216 if (!trkl.empty()) {
217 for (const auto& tkl : trkl) {
218 lut[tkl.firstClusterIndex + 1]++;
219 }
220 std::inclusive_scan(lut.begin(), lut.end(), lut.begin());
221 }
222 });
223
225 if (mTimeFrame->hasMCinformation() && mTrkParams[iteration].CreateArtefactLabels) {
226 tbb::parallel_for(0, static_cast<int>(topology.nTransitions), [&](const int transitionId) {
227 const auto& transition = topology.getTransition(transitionId);
228 for (auto& trk : mTimeFrame->getTracklets()[transitionId]) {
229 MCCompLabel label;
230 int currentId{mTimeFrame->getClusters()[transition.fromLayer][trk.firstClusterIndex].clusterId};
231 int nextId{mTimeFrame->getClusters()[transition.toLayer][trk.secondClusterIndex].clusterId};
232 for (const auto& lab1 : mTimeFrame->getClusterLabels(transition.fromLayer, currentId)) {
233 for (const auto& lab2 : mTimeFrame->getClusterLabels(transition.toLayer, nextId)) {
234 if (lab1 == lab2 && lab1.isValid()) {
235 label = lab1;
236 break;
237 }
238 }
239 if (label.isValid()) {
240 break;
241 }
242 }
243 mTimeFrame->getTrackletsLabel(transitionId).emplace_back(label);
244 }
245 });
246 }
247 });
248}
249
250template <int NLayers>
252{
253 const auto topology = mTimeFrame->getTrackingTopologyView();
254 for (int cellTopologyId = 0; cellTopologyId < topology.nCells; ++cellTopologyId) {
255 deepVectorClear(mTimeFrame->getCells()[cellTopologyId]);
256 deepVectorClear(mTimeFrame->getCellsLookupTable()[cellTopologyId]);
257 if (mTimeFrame->hasMCinformation() && mTrkParams[iteration].CreateArtefactLabels) {
258 deepVectorClear(mTimeFrame->getCellsLabel(cellTopologyId));
259 }
260 }
261
262 mTaskArena->execute([&] {
263 auto forTrackletCells = [&](auto Tag, int cellTopologyId, bounded_vector<CellSeed>& layerCells, int iTracklet, int offset = 0) -> int {
264 const auto& cellTopology = topology.getCell(cellTopologyId);
265 const auto& firstTransition = topology.getTransition(cellTopology.firstTransition);
266 const auto& secondTransition = topology.getTransition(cellTopology.secondTransition);
267 const Tracklet& currentTracklet{mTimeFrame->getTracklets()[cellTopology.firstTransition][iTracklet]};
268 const int nextLayerClusterIndex{currentTracklet.secondClusterIndex};
269 const int nextLayerFirstTrackletIndex{mTimeFrame->getTrackletsLookupTable()[cellTopology.secondTransition][nextLayerClusterIndex]};
270 const int nextLayerLastTrackletIndex{mTimeFrame->getTrackletsLookupTable()[cellTopology.secondTransition][nextLayerClusterIndex + 1]};
271 int foundCells{0};
272 for (int iNextTracklet{nextLayerFirstTrackletIndex}; iNextTracklet < nextLayerLastTrackletIndex; ++iNextTracklet) {
273 const Tracklet& nextTracklet{mTimeFrame->getTracklets()[cellTopology.secondTransition][iNextTracklet]};
274 if (nextTracklet.firstClusterIndex != nextLayerClusterIndex) {
275 break;
276 }
277 if (!currentTracklet.getTimeStamp().isCompatible(nextTracklet.getTimeStamp())) {
278 continue;
279 }
280
281 const float deltaTanLambdaSigma = std::abs(currentTracklet.tanLambda - nextTracklet.tanLambda) / mTrkParams[iteration].CellDeltaTanLambdaSigma;
282 if (deltaTanLambdaSigma < mTrkParams[iteration].NSigmaCut) {
283
285 const int clusId[3]{
286 mTimeFrame->getClusters()[firstTransition.fromLayer][currentTracklet.firstClusterIndex].clusterId,
287 mTimeFrame->getClusters()[firstTransition.toLayer][nextTracklet.firstClusterIndex].clusterId,
288 mTimeFrame->getClusters()[secondTransition.toLayer][nextTracklet.secondClusterIndex].clusterId};
289 const int hitLayers[3]{firstTransition.fromLayer, firstTransition.toLayer, secondTransition.toLayer};
290 const auto& cluster1_glo = mTimeFrame->getUnsortedClusters()[firstTransition.fromLayer][clusId[0]];
291 const auto& cluster2_glo = mTimeFrame->getUnsortedClusters()[firstTransition.toLayer][clusId[1]];
292 const auto& cluster3_tf = mTimeFrame->getTrackingFrameInfoOnLayer(secondTransition.toLayer)[clusId[2]];
293 auto track{o2::its::track::buildTrackSeed(cluster1_glo, cluster2_glo, cluster3_tf, mBz)};
294
295 float chi2{0.f};
296 bool good{false};
297 for (int iC{2}; iC--;) {
298 const int hitLayer = hitLayers[iC];
299 const TrackingFrameInfo& trackingHit = mTimeFrame->getTrackingFrameInfoOnLayer(hitLayer)[clusId[iC]];
300
301 if (!track.rotate(trackingHit.alphaTrackingFrame)) {
302 break;
303 }
304
305 if (!track.propagateTo(trackingHit.xTrackingFrame, getBz())) {
306 break;
307 }
308
309 if (!track.correctForMaterial(mTrkParams[iteration].LayerxX0[hitLayer], mTrkParams[iteration].LayerxX0[hitLayer] * constants::Radl * constants::Rho, true)) {
310 break;
311 }
312
313 const auto predChi2{track.getPredictedChi2Quiet(trackingHit.positionTrackingFrame, trackingHit.covarianceTrackingFrame)};
314 if (!iC && predChi2 > mTrkParams[iteration].MaxChi2ClusterAttachment) {
315 break;
316 }
317
318 if (!track.o2::track::TrackParCov::update(trackingHit.positionTrackingFrame, trackingHit.covarianceTrackingFrame)) {
319 break;
320 }
321
322 good = !iC;
323 chi2 += predChi2;
324 }
325 if (good) {
326 TimeEstBC ts = currentTracklet.getTimeStamp();
327 ts += nextTracklet.getTimeStamp();
328 if constexpr (decltype(Tag)::value == PassMode::OnePass::value) {
329 layerCells.emplace_back(cellTopology.hitLayerMask, clusId[0], clusId[1], clusId[2], iTracklet, iNextTracklet, track, chi2, ts);
330 ++foundCells;
331 } else if constexpr (decltype(Tag)::value == PassMode::TwoPassCount::value) {
332 ++foundCells;
333 } else if constexpr (decltype(Tag)::value == PassMode::TwoPassInsert::value) {
334 layerCells[offset++] = CellSeed(cellTopology.hitLayerMask, clusId[0], clusId[1], clusId[2], iTracklet, iNextTracklet, track, chi2, ts);
335 ++foundCells;
336 } else {
337 static_assert(false, "Unknown mode!");
338 }
339 }
340 }
341 }
342 return foundCells;
343 };
344
345 for (int cellTopologyId = 0; cellTopologyId < topology.nCells; ++cellTopologyId) {
346 const auto& cellTopology = topology.getCell(cellTopologyId);
347 if (mTimeFrame->getTracklets()[cellTopology.firstTransition].empty() ||
348 mTimeFrame->getTracklets()[cellTopology.secondTransition].empty()) {
349 continue;
350 }
351
352 auto& layerCells = mTimeFrame->getCells()[cellTopologyId];
353 const int currentLayerTrackletsNum{static_cast<int>(mTimeFrame->getTracklets()[cellTopology.firstTransition].size())};
354 bounded_vector<int> perTrackletCount(currentLayerTrackletsNum + 1, 0, mMemoryPool.get());
355 if (mTaskArena->max_concurrency() <= 1) {
356 for (int iTracklet{0}; iTracklet < currentLayerTrackletsNum; ++iTracklet) {
357 perTrackletCount[iTracklet] = forTrackletCells(PassMode::OnePass{}, cellTopologyId, layerCells, iTracklet);
358 }
359 std::exclusive_scan(perTrackletCount.begin(), perTrackletCount.end(), perTrackletCount.begin(), 0);
360 } else {
361 tbb::parallel_for(0, currentLayerTrackletsNum, [&](const int iTracklet) {
362 perTrackletCount[iTracklet] = forTrackletCells(PassMode::TwoPassCount{}, cellTopologyId, layerCells, iTracklet);
363 });
364
365 std::exclusive_scan(perTrackletCount.begin(), perTrackletCount.end(), perTrackletCount.begin(), 0);
366 auto totalCells{perTrackletCount.back()};
367 if (totalCells == 0) {
368 auto& lut = mTimeFrame->getCellsLookupTable()[cellTopologyId];
369 lut.resize(currentLayerTrackletsNum + 1);
370 std::fill(lut.begin(), lut.end(), 0);
371 continue;
372 }
373 layerCells.resize(totalCells);
374
375 tbb::parallel_for(0, currentLayerTrackletsNum, [&](const int iTracklet) {
376 int offset = perTrackletCount[iTracklet];
377 if (offset == perTrackletCount[iTracklet + 1]) {
378 return;
379 }
380 forTrackletCells(PassMode::TwoPassInsert{}, cellTopologyId, layerCells, iTracklet, offset);
381 });
382 }
383
384 auto& lut = mTimeFrame->getCellsLookupTable()[cellTopologyId];
385 lut.resize(currentLayerTrackletsNum + 1);
386 std::copy_n(perTrackletCount.begin(), currentLayerTrackletsNum + 1, lut.begin());
387
388 if (mTimeFrame->hasMCinformation() && mTrkParams[iteration].CreateArtefactLabels) {
389 auto& labels = mTimeFrame->getCellsLabel(cellTopologyId);
390 labels.reserve(layerCells.size());
391 for (const auto& cell : layerCells) {
392 MCCompLabel currentLab{mTimeFrame->getTrackletsLabel(cellTopology.firstTransition)[cell.getFirstTrackletIndex()]};
393 MCCompLabel nextLab{mTimeFrame->getTrackletsLabel(cellTopology.secondTransition)[cell.getSecondTrackletIndex()]};
394 labels.emplace_back(currentLab == nextLab ? currentLab : MCCompLabel());
395 }
396 }
397 }
398 });
399
400 for (int transitionId = 0; transitionId < topology.nTransitions; ++transitionId) {
401 deepVectorClear(mTimeFrame->getTracklets()[transitionId]);
402 deepVectorClear(mTimeFrame->getTrackletsLabel(transitionId));
403 }
404}
405
406template <int NLayers>
408{
409 const auto topology = mTimeFrame->getTrackingTopologyView();
410 mTaskArena->execute([&] {
411 std::vector<bounded_vector<CellNeighbour>> cellsNeighboursByTarget;
412 cellsNeighboursByTarget.reserve(topology.nCells);
413 for (int cellTopologyId{0}; cellTopologyId < topology.nCells; ++cellTopologyId) {
414 deepVectorClear(mTimeFrame->getCellsNeighbours()[cellTopologyId]);
415 deepVectorClear(mTimeFrame->getCellsNeighboursTopology()[cellTopologyId]);
416 deepVectorClear(mTimeFrame->getCellsNeighboursLUT()[cellTopologyId]);
417 cellsNeighboursByTarget.emplace_back(mMemoryPool.get());
418 }
419
420 for (int outerLayer{0}; outerLayer < NLayers; ++outerLayer) {
421 for (int cellTopologyId{0}; cellTopologyId < topology.nCells; ++cellTopologyId) {
422 const auto& cellTopology = topology.getCell(cellTopologyId);
423 if (cellTopology.hitLayerMask.last() != outerLayer ||
424 mTimeFrame->getCells()[cellTopologyId].empty()) {
425 continue;
426 }
427 const auto successors = topology.getCellsStartingWithTransition(cellTopology.secondTransition);
428 if (!successors.getEntries()) {
429 continue;
430 }
431
432 tbb::enumerable_thread_specific<bounded_vector<CellNeighbour>> sourceNeighbours([&]() { return bounded_vector<CellNeighbour>{mMemoryPool.get()}; });
433 tbb::parallel_for(0, static_cast<int>(mTimeFrame->getCells()[cellTopologyId].size()), [&](const int iCell) {
434 auto& localNeighbours = sourceNeighbours.local();
435 const auto& currentCellSeed{mTimeFrame->getCells()[cellTopologyId][iCell]};
436 const int nextLayerTrackletIndex{currentCellSeed.getSecondTrackletIndex()};
437 for (int iSuccessor{0}; iSuccessor < successors.getEntries(); ++iSuccessor) {
438 const int nextCellTopologyId = topology.cellsByFirstTransition[successors.getFirstEntry() + iSuccessor];
439 if (mTimeFrame->getCells()[nextCellTopologyId].empty() ||
440 mTimeFrame->getCellsLookupTable()[nextCellTopologyId].empty()) {
441 continue;
442 }
443 const auto& nextCellLUT = mTimeFrame->getCellsLookupTable()[nextCellTopologyId];
444 if (nextLayerTrackletIndex + 1 >= static_cast<int>(nextCellLUT.size())) {
445 continue;
446 }
447 const int nextLayerFirstCellIndex{nextCellLUT[nextLayerTrackletIndex]};
448 const int nextLayerLastCellIndex{nextCellLUT[nextLayerTrackletIndex + 1]};
449 for (int iNextCell{nextLayerFirstCellIndex}; iNextCell < nextLayerLastCellIndex; ++iNextCell) {
450 const auto& nextCellSeedRef{mTimeFrame->getCells()[nextCellTopologyId][iNextCell]};
451 if (nextCellSeedRef.getFirstTrackletIndex() != nextLayerTrackletIndex || !currentCellSeed.getTimeStamp().isCompatible(nextCellSeedRef.getTimeStamp())) {
452 break;
453 }
454
455 auto nextCellSeed{mTimeFrame->getCells()[nextCellTopologyId][iNextCell]};
456 if (!nextCellSeed.rotate(currentCellSeed.getAlpha()) ||
457 !nextCellSeed.propagateTo(currentCellSeed.getX(), getBz())) {
458 continue;
459 }
460
461 float chi2 = currentCellSeed.getPredictedChi2(nextCellSeed);
462 if (chi2 > mTrkParams[iteration].MaxChi2ClusterAttachment) {
463 continue;
464 }
465
466 const int nextLevel = currentCellSeed.getLevel() + 1;
467 localNeighbours.emplace_back(cellTopologyId, iCell, nextCellTopologyId, iNextCell, nextLevel);
468 }
469 }
470 });
471
472 bounded_vector<size_t> count(topology.nCells, 0, mMemoryPool.get());
473 for (const auto& localNeighbours : sourceNeighbours) {
474 for (const auto& neigh : localNeighbours) {
475 ++count[neigh.nextCellTopology];
476 }
477 }
478 for (size_t i{0}; i < topology.nCells; ++i) {
479 cellsNeighboursByTarget[i].reserve(count[i]);
480 }
481 for (const auto& localNeighbours : sourceNeighbours) {
482 for (const auto& neigh : localNeighbours) {
483 cellsNeighboursByTarget[neigh.nextCellTopology].emplace_back(neigh);
484 if (neigh.level > mTimeFrame->getCells()[neigh.nextCellTopology][neigh.nextCell].getLevel()) {
485 mTimeFrame->getCells()[neigh.nextCellTopology][neigh.nextCell].setLevel(neigh.level);
486 }
487 }
488 }
489 }
490 }
491
492 for (int cellTopologyId{0}; cellTopologyId < topology.nCells; ++cellTopologyId) {
493 auto& cellsNeighbours = cellsNeighboursByTarget[cellTopologyId];
494 if (cellsNeighbours.empty()) {
495 continue;
496 }
497
498 std::sort(cellsNeighbours.begin(), cellsNeighbours.end(), [](const auto& a, const auto& b) {
499 return a.nextCell < b.nextCell;
500 });
501
502 auto& cellsNeighbourLUT = mTimeFrame->getCellsNeighboursLUT()[cellTopologyId];
503 cellsNeighbourLUT.assign(mTimeFrame->getCells()[cellTopologyId].size(), 0);
504 for (const auto& neigh : cellsNeighbours) {
505 ++cellsNeighbourLUT[neigh.nextCell];
506 }
507 std::inclusive_scan(cellsNeighbourLUT.begin(), cellsNeighbourLUT.end(), cellsNeighbourLUT.begin());
508
509 mTimeFrame->getCellsNeighbours()[cellTopologyId].reserve(cellsNeighbours.size());
510 mTimeFrame->getCellsNeighboursTopology()[cellTopologyId].reserve(cellsNeighbours.size());
511 std::ranges::transform(cellsNeighbours, std::back_inserter(mTimeFrame->getCellsNeighbours()[cellTopologyId]), [](const auto& neigh) { return neigh.cell; });
512 std::ranges::transform(cellsNeighbours, std::back_inserter(mTimeFrame->getCellsNeighboursTopology()[cellTopologyId]), [](const auto& neigh) { return neigh.cellTopology; });
513 }
514
515 // clean up LUTs
516 for (auto& cellLUT : mTimeFrame->getCellsLookupTable()) {
517 deepVectorClear(cellLUT);
518 }
519 });
520}
521
522template <int NLayers>
523template <typename InputSeed>
524void TrackerTraits<NLayers>::processNeighbours(int iteration, int defaultCellTopologyId, int iLevel, const bounded_vector<InputSeed>& currentCellSeed, const bounded_vector<int>& currentCellId, const bounded_vector<int>& currentCellTopologyId, bounded_vector<TrackSeedN>& updatedCellSeeds, bounded_vector<int>& updatedCellsIds, bounded_vector<int>& updatedCellsTopologyIds)
525{
526 auto propagator = o2::base::Propagator::Instance();
527
528 mTaskArena->execute([&] {
529 auto forCellNeighbours = [&](auto Tag, int iCell, int offset = 0) -> int {
530 const auto& currentCell{currentCellSeed[iCell]};
531 const int cellTopologyId = currentCellTopologyId.empty() ? defaultCellTopologyId : currentCellTopologyId[iCell];
532
533 if constexpr (decltype(Tag)::value != PassMode::TwoPassInsert::value) {
534 if (currentCell.getLevel() != iLevel) {
535 return 0;
536 }
537 if (currentCellId.empty()) {
538 for (int layer = 0; layer < NLayers; ++layer) {
539 const int clusterIndex = currentCell.getCluster(layer);
540 if (clusterIndex != constants::UnusedIndex && mTimeFrame->isClusterUsed(layer, clusterIndex)) {
541 return 0;
542 }
543 }
544 }
545 }
546
547 const int cellId = currentCellId.empty() ? iCell : currentCellId[iCell];
548 if (cellTopologyId < 0 || mTimeFrame->getCellsNeighboursLUT()[cellTopologyId].empty()) {
549 return 0;
550 }
551 const int startNeighbourId{cellId ? mTimeFrame->getCellsNeighboursLUT()[cellTopologyId][cellId - 1] : 0};
552 const int endNeighbourId{mTimeFrame->getCellsNeighboursLUT()[cellTopologyId][cellId]};
553 int foundSeeds{0};
554 for (int iNeighbourCell{startNeighbourId}; iNeighbourCell < endNeighbourId; ++iNeighbourCell) {
555 const int neighbourCellTopologyId = mTimeFrame->getCellsNeighboursTopology()[cellTopologyId][iNeighbourCell];
556 const int neighbourCellId = mTimeFrame->getCellsNeighbours()[cellTopologyId][iNeighbourCell];
557 const auto& neighbourCell = mTimeFrame->getCells()[neighbourCellTopologyId][neighbourCellId];
558 if (neighbourCell.getSecondTrackletIndex() != currentCell.getFirstTrackletIndex()) {
559 continue;
560 }
561 if (!currentCell.getTimeStamp().isCompatible(neighbourCell.getTimeStamp())) {
562 continue;
563 }
564 if (currentCell.getLevel() - 1 != neighbourCell.getLevel()) {
565 continue;
566 }
567 const int neighbourLayer = neighbourCell.getInnerLayer();
568 const int neighbourCluster = neighbourCell.getFirstClusterIndex();
569 if (mTimeFrame->isClusterUsed(neighbourLayer, neighbourCluster)) {
570 continue;
571 }
572
574 TrackSeedN seed{currentCell};
575 seed.getTimeStamp() = currentCell.getTimeStamp();
576 seed.getTimeStamp() += neighbourCell.getTimeStamp();
577 const auto& trHit = mTimeFrame->getTrackingFrameInfoOnLayer(neighbourLayer)[neighbourCluster];
578
579 if (!seed.rotate(trHit.alphaTrackingFrame)) {
580 continue;
581 }
582
583 if (!propagator->propagateToX(seed, trHit.xTrackingFrame, getBz(), o2::base::PropagatorImpl<float>::MAX_SIN_PHI, o2::base::PropagatorImpl<float>::MAX_STEP, mTrkParams[iteration].CorrType)) {
584 continue;
585 }
586
587 if (mTrkParams[iteration].CorrType == o2::base::PropagatorF::MatCorrType::USEMatCorrNONE) {
588 if (!seed.correctForMaterial(mTrkParams[iteration].LayerxX0[neighbourLayer], mTrkParams[iteration].LayerxX0[neighbourLayer] * constants::Radl * constants::Rho, true)) {
589 continue;
590 }
591 }
592
593 auto predChi2{seed.getPredictedChi2Quiet(trHit.positionTrackingFrame, trHit.covarianceTrackingFrame)};
594 if ((predChi2 > mTrkParams[iteration].MaxChi2ClusterAttachment) || predChi2 < 0.f) {
595 continue;
596 }
597 seed.setChi2(seed.getChi2() + predChi2);
598 if (!seed.o2::track::TrackParCov::update(trHit.positionTrackingFrame, trHit.covarianceTrackingFrame)) {
599 continue;
600 }
601
602 if constexpr (decltype(Tag)::value != PassMode::TwoPassCount::value) {
603 seed.getClusters()[neighbourLayer] = neighbourCluster;
604 auto mask = seed.getHitLayerMask();
605 mask.set(neighbourLayer);
606 seed.setHitLayerMask(mask);
607 seed.setLevel(neighbourCell.getLevel());
608 seed.setFirstTrackletIndex(neighbourCell.getFirstTrackletIndex());
609 seed.setSecondTrackletIndex(neighbourCell.getSecondTrackletIndex());
610 }
611
612 if constexpr (decltype(Tag)::value == PassMode::OnePass::value) {
613 updatedCellSeeds.push_back(seed);
614 updatedCellsIds.push_back(neighbourCellId);
615 updatedCellsTopologyIds.push_back(neighbourCellTopologyId);
616 } else if constexpr (decltype(Tag)::value == PassMode::TwoPassCount::value) {
617 ++foundSeeds;
618 } else if constexpr (decltype(Tag)::value == PassMode::TwoPassInsert::value) {
619 updatedCellSeeds[offset] = seed;
620 updatedCellsIds[offset] = neighbourCellId;
621 updatedCellsTopologyIds[offset++] = neighbourCellTopologyId;
622 } else {
623 static_assert(false, "Unknown mode!");
624 }
625 }
626 return foundSeeds;
627 };
628
629 const int nCells = static_cast<int>(currentCellSeed.size());
630 if (mTaskArena->max_concurrency() <= 1) {
631 for (int iCell{0}; iCell < nCells; ++iCell) {
632 forCellNeighbours(PassMode::OnePass{}, iCell);
633 }
634 } else {
635 bounded_vector<int> perCellCount(nCells + 1, 0, mMemoryPool.get());
636 tbb::parallel_for(0, nCells, [&](const int iCell) {
637 perCellCount[iCell] = forCellNeighbours(PassMode::TwoPassCount{}, iCell);
638 });
639
640 std::exclusive_scan(perCellCount.begin(), perCellCount.end(), perCellCount.begin(), 0);
641 auto totalNeighbours{perCellCount.back()};
642 if (totalNeighbours == 0) {
643 return;
644 }
645 updatedCellSeeds.resize(totalNeighbours);
646 updatedCellsIds.resize(totalNeighbours);
647 updatedCellsTopologyIds.resize(totalNeighbours);
648
649 tbb::parallel_for(0, nCells, [&](const int iCell) {
650 int offset = perCellCount[iCell];
651 if (offset == perCellCount[iCell + 1]) {
652 return;
653 }
654 forCellNeighbours(PassMode::TwoPassInsert{}, iCell, offset);
655 });
656 }
657 });
658}
659
660template <int NLayers>
661void TrackerTraits<NLayers>::findRoads(const int iteration)
662{
663 bounded_vector<bounded_vector<int>> firstClusters(mTrkParams[iteration].NLayers, bounded_vector<int>(mMemoryPool.get()), mMemoryPool.get());
664 firstClusters.resize(mTrkParams[iteration].NLayers);
665 const auto propagator = o2::base::Propagator::Instance();
666 const TrackingFrameInfo* tfInfos[NLayers]{};
667 const Cluster* unsortedClusters[NLayers]{};
668 for (int iLayer = 0; iLayer < NLayers; ++iLayer) {
669 tfInfos[iLayer] = mTimeFrame->getTrackingFrameInfoOnLayer(iLayer).data();
670 unsortedClusters[iLayer] = mTimeFrame->getUnsortedClusters()[iLayer].data();
671 }
672 const auto topology = mTimeFrame->getTrackingTopologyView();
673 for (int startLevel{mTrkParams[iteration].CellsPerRoad()}; startLevel >= mTrkParams[iteration].CellMinimumLevel(); --startLevel) {
674
675 auto seedFilter = [&](const auto& seed) {
676 return seed.getHitLayerMask().isAllowed(mTrkParams[iteration].MaxHoles, mTrkParams[iteration].HoleLayerMask) &&
677 seed.getHitLayerMask().length() >= mTrkParams[iteration].MinTrackLength &&
678 seed.getQ2Pt() <= 1.e3 && seed.getChi2() <= mTrkParams[iteration].MaxChi2NDF * ((startLevel + 2) * 2 - 5);
679 };
680
681 bounded_vector<TrackSeedN> trackSeeds(mMemoryPool.get());
682 for (int startCellTopologyId{0}; startCellTopologyId < topology.nCells; ++startCellTopologyId) {
683 const int startLayer = topology.getCell(startCellTopologyId).hitLayerMask.last();
684 if (!(mTrkParams[iteration].StartLayerMask.has(startLayer)) || mTimeFrame->getCells()[startCellTopologyId].empty()) {
685 continue;
686 }
687
688 bounded_vector<int> lastCellId(mMemoryPool.get()), updatedCellId(mMemoryPool.get());
689 bounded_vector<int> lastCellTopologyId(mMemoryPool.get()), updatedCellTopologyId(mMemoryPool.get());
690 bounded_vector<TrackSeedN> lastCellSeed(mMemoryPool.get()), updatedCellSeed(mMemoryPool.get());
691
692 processNeighbours(iteration, startCellTopologyId, startLevel, mTimeFrame->getCells()[startCellTopologyId], lastCellId, lastCellTopologyId, updatedCellSeed, updatedCellId, updatedCellTopologyId);
693
694 int level = startLevel;
695 while (level > 2 && !updatedCellSeed.empty()) {
696 lastCellSeed.swap(updatedCellSeed);
697 lastCellId.swap(updatedCellId);
698 lastCellTopologyId.swap(updatedCellTopologyId);
699 deepVectorClear(updatedCellSeed);
700 deepVectorClear(updatedCellId);
701 deepVectorClear(updatedCellTopologyId);
702 processNeighbours(iteration, constants::UnusedIndex, --level, lastCellSeed, lastCellId, lastCellTopologyId, updatedCellSeed, updatedCellId, updatedCellTopologyId);
703 }
704 deepVectorClear(lastCellId);
705 deepVectorClear(lastCellTopologyId);
706 deepVectorClear(lastCellSeed);
707
708 if (!updatedCellSeed.empty()) {
709 trackSeeds.reserve(trackSeeds.size() + std::count_if(updatedCellSeed.begin(), updatedCellSeed.end(), seedFilter));
710 std::copy_if(updatedCellSeed.begin(), updatedCellSeed.end(), std::back_inserter(trackSeeds), seedFilter);
711 }
712 }
713
714 if (trackSeeds.empty()) {
715 continue;
716 }
717
718 bounded_vector<TrackITSExt> tracks(mMemoryPool.get());
719 mTaskArena->execute([&] {
720 auto forSeed = [&](auto Tag, int iSeed, int offset = 0) {
721 TrackITSExt temporaryTrack;
722 bool refitSuccess = track::refitTrack<NLayers>(trackSeeds[iSeed],
723 temporaryTrack,
724 mTrkParams[iteration].MaxChi2ClusterAttachment,
725 mTrkParams[iteration].MaxChi2NDF,
726 mBz,
727 tfInfos,
728 unsortedClusters,
729 mTrkParams[iteration].LayerxX0.data(),
730 mTrkParams[iteration].LayerRadii.data(),
731 mTrkParams[iteration].MinPt.data(),
732 propagator,
733 mTrkParams[iteration].CorrType,
734 mTrkParams[iteration].ReseedIfShorter,
735 mTrkParams[iteration].ShiftRefToCluster,
736 mTrkParams[iteration].RepeatRefitOut);
737
738 if (refitSuccess) {
739 if constexpr (decltype(Tag)::value == PassMode::OnePass::value) {
740 tracks.push_back(temporaryTrack);
741 } else if constexpr (decltype(Tag)::value == PassMode::TwoPassCount::value) {
742 // nothing to do
743 } else if constexpr (decltype(Tag)::value == PassMode::TwoPassInsert::value) {
744 tracks[offset] = temporaryTrack;
745 } else {
746 static_assert(false, "Unknown mode!");
747 }
748 return 1;
749 }
750 return 0;
751 };
752
753 const int nSeeds = static_cast<int>(trackSeeds.size());
754 if (mTaskArena->max_concurrency() <= 1) {
755 for (int iSeed{0}; iSeed < nSeeds; ++iSeed) {
756 forSeed(PassMode::OnePass{}, iSeed);
757 }
758 } else {
759 // The double-pass allows us to avoid sizeable memory spikes
760 bounded_vector<int> perSeedCount(nSeeds + 1, 0, mMemoryPool.get());
761 tbb::parallel_for(0, nSeeds, [&](const int iSeed) {
762 perSeedCount[iSeed] = forSeed(PassMode::TwoPassCount{}, iSeed);
763 });
764
765 std::exclusive_scan(perSeedCount.begin(), perSeedCount.end(), perSeedCount.begin(), 0);
766 auto totalTracks{perSeedCount.back()};
767 if (totalTracks == 0) {
768 return;
769 }
770 tracks.resize(totalTracks);
771
772 tbb::parallel_for(0, nSeeds, [&](const int iSeed) {
773 if (perSeedCount[iSeed] == perSeedCount[iSeed + 1]) {
774 return;
775 }
776 forSeed(PassMode::TwoPassInsert{}, iSeed, perSeedCount[iSeed]);
777 });
778 }
779
780 deepVectorClear(trackSeeds);
781 });
782
783 std::sort(tracks.begin(), tracks.end(), [](const auto& a, const auto& b) {
784 return track::isBetter(a, b);
785 });
786
787 acceptTracks(iteration, tracks, firstClusters);
788 }
789 markTracks(iteration);
790}
791
792template <int NLayers>
794{
795 auto& trks = mTimeFrame->getTracks();
796 trks.reserve(trks.size() + tracks.size());
797 const float smallestROFHalf = mTimeFrame->getROFOverlapTableView().getClockLayer().mROFLength * 0.5f;
798 for (auto& track : tracks) {
799 int nShared = 0;
800 bool isFirstShared{false};
801 int firstLayer{-1}, firstCluster{-1};
802 for (int iLayer{0}; iLayer < mTrkParams[iteration].NLayers; ++iLayer) {
803 if (track.getClusterIndex(iLayer) == constants::UnusedIndex) {
804 continue;
805 }
806 bool isShared = mTimeFrame->isClusterUsed(iLayer, track.getClusterIndex(iLayer));
807 nShared += int(isShared);
808 if (firstLayer < 0) {
809 firstCluster = track.getClusterIndex(iLayer);
810 isFirstShared = isShared && mTrkParams[iteration].AllowSharingFirstCluster && std::find(firstClusters[iLayer].begin(), firstClusters[iLayer].end(), firstCluster) != firstClusters[iLayer].end();
811 firstLayer = iLayer;
812 }
813 }
814
816 if (nShared - int(isFirstShared && mTrkParams[iteration].AllowSharingFirstCluster) > mTrkParams[iteration].SharedMaxClusters) {
817 continue;
818 }
819
820 bool firstCls{true}, nominalCompatible{true};
821 TimeEstBC nominalTS, expandedTS;
822 for (int iLayer{0}; iLayer < mTrkParams[iteration].NLayers; ++iLayer) {
823 if (track.getClusterIndex(iLayer) == constants::UnusedIndex) {
824 continue;
825 }
826 mTimeFrame->markUsedCluster(iLayer, track.getClusterIndex(iLayer));
827 int currentROF = mTimeFrame->getClusterROF(iLayer, track.getClusterIndex(iLayer));
828 const auto nominalROFTS = mTimeFrame->getROFOverlapTableView().getLayer(iLayer).getROFTimeBounds(currentROF);
829 const auto expandedROFTS = mTimeFrame->getROFOverlapTableView().getLayer(iLayer).getROFTimeBounds(currentROF, true);
830 if (firstCls) {
831 firstCls = false;
832 nominalTS = nominalROFTS;
833 expandedTS = expandedROFTS;
834 } else {
835 if (nominalCompatible) {
836 if (nominalTS.isCompatible(nominalROFTS)) {
837 nominalTS += nominalROFTS;
838 } else {
839 nominalCompatible = false;
840 }
841 }
842 if (!expandedTS.isCompatible(expandedROFTS)) {
843 LOGP(fatal, "TS {}+/-{} are incompatible with {}+/-{}, this should not happen!", expandedROFTS.getTimeStamp(), expandedROFTS.getTimeStampError(), expandedTS.getTimeStamp(), expandedTS.getTimeStampError());
844 }
845 expandedTS += expandedROFTS;
846 }
847 }
848 track.getTimeStamp() = (nominalCompatible ? nominalTS : expandedTS).makeSymmetrical();
849 // this is a sanity clamp
850 // we cannot be worse than the clock so we clamp to this
851 if (track.getTimeStamp().getTimeStampError() > smallestROFHalf) {
852 track.getTimeStamp().setTimeStampError(smallestROFHalf);
853 }
854 track.setUserField(0);
855 track.getParamOut().setUserField(0);
856 trks.emplace_back(track);
857
858 if (mTrkParams[iteration].AllowSharingFirstCluster) {
859 firstClusters[firstLayer].push_back(firstCluster);
860 }
861 }
862}
863
864template <int NLayers>
866{
867 if (mTrkParams[iteration].AllowSharingFirstCluster) {
869 auto& tracks = mTimeFrame->getTracks();
870
871 bounded_vector<int> fclusSort(tracks.size(), mMemoryPool.get());
872 std::iota(fclusSort.begin(), fclusSort.end(), 0);
873 std::sort(fclusSort.begin(), fclusSort.end(), [&tracks](int a, int b) {
874 return tracks[a].getFirstLayerClusterIndex() < tracks[b].getFirstLayerClusterIndex();
875 });
876
877 auto areTracksSelected = [this, iteration](const TrackITSExt& t1, const TrackITSExt& t2) {
878 const auto t1FirstLayer{t1.getFirstClusterLayer()}, t2FirstLayer{t2.getFirstClusterLayer()};
879 if (t1FirstLayer != t2FirstLayer) {
880 return false;
881 }
882 if (mTimeFrame->getClusterROF(t1FirstLayer, t1.getClusterIndex(t1FirstLayer)) != mTimeFrame->getClusterROF(t2FirstLayer, t2.getClusterIndex(t2FirstLayer))) {
883 return false;
884 }
885 if (!math_utils::isPhiDifferenceBelow(t1.getPhi(), t2.getPhi(), mTrkParams[iteration].SharedClusterMaxDeltaPhi)) {
886 return false;
887 }
888 if (std::abs(t1.getEta() - t2.getEta()) > mTrkParams[iteration].SharedClusterMaxDeltaEta) {
889 return false;
890 }
891 if (mTrkParams[iteration].SharedClusterOppositeSign && t1.getSign() == t2.getSign()) {
892 return false;
893 }
894 return true;
895 };
896
897 for (int i{0}; i < static_cast<int>(fclusSort.size()); ++i) {
898 auto& track = tracks[fclusSort[i]];
899 for (int j{i + 1}; j < static_cast<int>(fclusSort.size()) && tracks[fclusSort[j]].getFirstLayerClusterIndex() == track.getFirstLayerClusterIndex(); ++j) {
900 auto& track2 = tracks[fclusSort[j]];
901 if (areTracksSelected(track, track2)) {
902 track.setSharedClusters();
903 track2.setSharedClusters();
904 }
905 }
906 }
907 }
908}
909
910template <int NLayers>
912{
913 mBz = bz;
914 mTimeFrame->setBz(bz);
915}
916
917template <int NLayers>
918void TrackerTraits<NLayers>::setNThreads(int n, std::shared_ptr<tbb::task_arena>& arena)
919{
920#if defined(OPTIMISATION_OUTPUT)
921 mTaskArena = std::make_shared<tbb::task_arena>(1);
922#else
923 if (arena == nullptr) {
924 mTaskArena = std::make_shared<tbb::task_arena>(std::abs(n));
925 LOGP(info, "Setting tracker with {} threads.", n);
926 } else {
927 mTaskArena = arena;
928 }
929#endif
930}
931
932template class TrackerTraits<7>;
935// ALICE3 upgrade
936#ifdef ENABLE_UPGRADES
937template class TrackerTraits<11>;
940#endif
941
942} // namespace o2::its
std::vector< std::string > labels
int32_t i
int32_t lastRow
bounded_vector< float > bins
uint32_t j
Definition RawData.h:0
Shared host/device helpers for ITS tracker trait implementations.
GPUd() value_type estimateLTFast(o2 static GPUd() float estimateLTIncrement(const o2 PropagatorImpl * Instance(bool uninitialized=false)
Definition Propagator.h:178
CellSeed: connections of three clusters.
Definition Cell.h:81
virtual void computeLayerTracklets(const int iteration, int iVertex)
GLdouble n
Definition glcorearb.h:1982
GLint GLsizei count
Definition glcorearb.h:399
GLuint GLuint end
Definition glcorearb.h:469
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLsizei const GLfloat * value
Definition glcorearb.h:819
GLintptr offset
Definition glcorearb.h:660
GLenum GLuint GLint GLint layer
Definition glcorearb.h:1310
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLint GLuint mask
Definition glcorearb.h:291
GLuint GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat t1
Definition glcorearb.h:5034
constexpr float Tolerance
Definition Constants.h:30
void deepVectorClear(std::vector< T > &vec)
std::pmr::vector< T > bounded_vector
float yCoordinate
Definition Cluster.h:48
float zCoordinate
Definition Cluster.h:49
float xCoordinate
Definition Cluster.h:47
std::integral_constant< int, 0 > OnePass
std::integral_constant< int, 1 > TwoPassCount
std::integral_constant< int, 2 > TwoPassInsert
std::array< float, 2 > positionTrackingFrame
Definition Cluster.h:75
std::array< float, 3 > covarianceTrackingFrame
Definition Cluster.h:76
std::vector< Tracklet64 > tracklets