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