Project
Loading...
Searching...
No Matches
VertexerTraits.cxx
Go to the documentation of this file.
1// Copyright 2019-2026 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.
12
13#include <algorithm>
14#include <memory>
15#include <ranges>
16#include <span>
17#include <unordered_map>
18
19#include <oneapi/tbb/blocked_range.h>
20#include <oneapi/tbb/parallel_for.h>
21#include <oneapi/tbb/combinable.h>
22
34
35namespace o2::its
36{
37namespace
38{
39
40template <TrackletMode Mode, bool EvalRun, int NLayers>
41void trackleterKernelHost(
42 const gsl::span<const Cluster>& clustersNextLayer, // 0 2
43 const gsl::span<const Cluster>& clustersCurrentLayer, // 1 1
44 const gsl::span<uint8_t>& usedClustersNextLayer, // 0 2
45 const int* indexTableNext,
46 const float phiCut,
47 bounded_vector<Tracklet>& tracklets,
48 gsl::span<int> foundTracklets,
49 const IndexTableUtils<NLayers>& utils,
50 const TimeEstBC& timErr,
51 gsl::span<int> rofFoundTrackletsOffsets,
52 const int globalOffsetNextLayer,
53 const int globalOffsetCurrentLayer,
54 const int maxTrackletsPerCluster)
55{
56 const int PhiBins{utils.getNphiBins()};
57 const int ZBins{utils.getNzBins()};
58 // loop on layer1 clusters
59 for (int iCurrentLayerClusterIndex = 0; iCurrentLayerClusterIndex < clustersCurrentLayer.size(); ++iCurrentLayerClusterIndex) {
60 int storedTracklets{0};
61 const Cluster& currentCluster{clustersCurrentLayer[iCurrentLayerClusterIndex]};
62 const int4 selectedBinsRect{o2::its::getBinsRect(currentCluster, (int)Mode + 1, 0.f, 0.f, 100.f, phiCut / 2, utils)};
63 if (selectedBinsRect.x >= 0) {
64 int phiBinsNum{selectedBinsRect.w - selectedBinsRect.y + 1};
65 if (phiBinsNum < 0) {
66 phiBinsNum += PhiBins;
67 }
68 // loop on phi bins next layer
69 for (int iPhiBin{selectedBinsRect.y}, iPhiCount{0}; iPhiCount < phiBinsNum && storedTracklets < maxTrackletsPerCluster; iPhiBin = ++iPhiBin == PhiBins ? 0 : iPhiBin, iPhiCount++) {
70 const int firstBinIndex{utils.getBinIndex(selectedBinsRect.x, iPhiBin)};
71 const int firstRowClusterIndex{indexTableNext[firstBinIndex]};
72 const int maxRowClusterIndex{indexTableNext[firstBinIndex + ZBins]};
73 // loop on clusters next layer
74 for (int iNextLayerClusterIndex{firstRowClusterIndex}; iNextLayerClusterIndex < maxRowClusterIndex && iNextLayerClusterIndex < static_cast<int>(clustersNextLayer.size()) && storedTracklets < maxTrackletsPerCluster; ++iNextLayerClusterIndex) {
75 if (usedClustersNextLayer[iNextLayerClusterIndex]) {
76 continue;
77 }
78 const Cluster& nextCluster{clustersNextLayer[iNextLayerClusterIndex]};
79 if (math_utils::isPhiDifferenceBelow(currentCluster.phi, nextCluster.phi, phiCut)) {
80 if (storedTracklets < maxTrackletsPerCluster) {
81 if constexpr (!EvalRun) {
82 if constexpr (Mode == TrackletMode::Layer0Layer1) {
83 tracklets[rofFoundTrackletsOffsets[iCurrentLayerClusterIndex] + storedTracklets] = Tracklet{globalOffsetNextLayer + iNextLayerClusterIndex, globalOffsetCurrentLayer + iCurrentLayerClusterIndex, nextCluster, currentCluster, timErr};
84 } else {
85 tracklets[rofFoundTrackletsOffsets[iCurrentLayerClusterIndex] + storedTracklets] = Tracklet{globalOffsetCurrentLayer + iCurrentLayerClusterIndex, globalOffsetNextLayer + iNextLayerClusterIndex, currentCluster, nextCluster, timErr};
86 }
87 }
88 ++storedTracklets;
89 }
90 }
91 }
92 }
93 }
94 if constexpr (EvalRun) {
95 foundTracklets[iCurrentLayerClusterIndex] += storedTracklets;
96 } else {
97 rofFoundTrackletsOffsets[iCurrentLayerClusterIndex] += storedTracklets;
98 }
99 }
100}
101
102void trackletSelectionKernelHost(
103 const Cluster* clusters0, // global layer 0 clusters
104 const Cluster* clusters1, // global layer 1 clusters
105 gsl::span<unsigned char> usedClusters0, // global layer 0 used clusters
106 gsl::span<unsigned char> usedClusters2, // global layer 2 used clusters
107 const gsl::span<const Tracklet>& tracklets01,
108 const gsl::span<const Tracklet>& tracklets12,
109 bounded_vector<uint8_t>& usedTracklets,
110 const gsl::span<int> foundTracklets01,
111 const gsl::span<int> foundTracklets12,
112 bounded_vector<Line>& lines,
113 const gsl::span<const o2::MCCompLabel>& trackletLabels,
114 bounded_vector<o2::MCCompLabel>& linesLabels,
115 const int nLayer1Clusters,
116 const float tanLambdaCut,
117 const float phiCut,
118 const int maxTracklets)
119{
120 int offset01{0}, offset12{0};
121 for (int iCurrentLayerClusterIndex{0}; iCurrentLayerClusterIndex < nLayer1Clusters; ++iCurrentLayerClusterIndex) {
122 int validTracklets{0};
123 const int endTracklet01 = offset01 + foundTracklets01[iCurrentLayerClusterIndex];
124 const int endTracklet12 = offset12 + foundTracklets12[iCurrentLayerClusterIndex];
125 for (int iTracklet12{offset12}; iTracklet12 < endTracklet12 && validTracklets != maxTracklets; ++iTracklet12) {
126 const auto& tracklet12{tracklets12[iTracklet12]};
127 for (int iTracklet01{offset01}; iTracklet01 < endTracklet01 && validTracklets != maxTracklets; ++iTracklet01) {
128 if (usedTracklets[iTracklet01]) {
129 continue;
130 }
131
132 const auto& tracklet01{tracklets01[iTracklet01]};
133 if (!tracklet01.getTimeStamp().isCompatible(tracklet12.getTimeStamp())) {
134 continue;
135 }
136
137 const float deltaTanLambda{o2::gpu::GPUCommonMath::Abs(tracklet01.tanLambda - tracklet12.tanLambda)};
138 if (deltaTanLambda >= tanLambdaCut) {
139 continue;
140 }
141 if (math_utils::isPhiDifferenceBelow(tracklet01.phi, tracklet12.phi, phiCut) && validTracklets != maxTracklets) {
142 usedClusters0[tracklet01.firstClusterIndex] = 1;
143 usedClusters2[tracklet12.secondClusterIndex] = 1;
144 usedTracklets[iTracklet01] = true;
145 lines.emplace_back(tracklet01, clusters0, clusters1);
146 if (!trackletLabels.empty()) {
147 linesLabels.emplace_back(trackletLabels[iTracklet01]);
148 }
149 ++validTracklets;
150 }
151 }
152 }
153 offset01 += foundTracklets01[iCurrentLayerClusterIndex];
154 offset12 += foundTracklets12[iCurrentLayerClusterIndex];
155 }
156}
157} // namespace
158
159template <int NLayers>
161{
162 mTimeFrame->initialise(trackingParams, 3);
163}
164
165template <int NLayers>
166void VertexerTraits<NLayers>::updateVertexingParameters(const std::vector<VertexingParameters>& vrtPar)
167{
168 mVrtParams = vrtPar;
169 mIndexTableUtils.setTrackingParameters(vrtPar[0]);
170 for (auto& par : mVrtParams) {
171 par.phiSpan = static_cast<int>(std::ceil(mIndexTableUtils.getNphiBins() * par.phiCut / o2::constants::math::TwoPI));
172 par.zSpan = static_cast<int>(std::ceil(par.zCut * mIndexTableUtils.getInverseZCoordinate(0)));
173 }
174}
175
176// Main functions
177template <int NLayers>
179{
180 mTaskArena->execute([&] {
181 tbb::parallel_for(0, mTimeFrame->getNrof(1), [&](const short pivotRofId) {
182 bool skip = skipROF(iteration, pivotRofId);
183 const auto& rofRange01 = mTimeFrame->getROFOverlapTableView().getOverlap(1, 0, pivotRofId);
184 for (auto targetRofId = rofRange01.getFirstEntry(); targetRofId < rofRange01.getEntriesBound(); ++targetRofId) {
185 const auto timeErr = mTimeFrame->getROFOverlapTableView().getTimeStamp(0, targetRofId, 1, pivotRofId);
186 trackleterKernelHost<TrackletMode::Layer0Layer1, true>(
187 !skip ? mTimeFrame->getClustersOnLayer(targetRofId, 0) : gsl::span<Cluster>(), // Clusters to be matched with the next layer in target rof
188 !skip ? mTimeFrame->getClustersOnLayer(pivotRofId, 1) : gsl::span<Cluster>(), // Clusters to be matched with the current layer in pivot rof
189 mTimeFrame->getUsedClustersROF(targetRofId, 0), // Span of the used clusters in the target rof
190 mTimeFrame->getIndexTable(targetRofId, 0).data(), // Index table to access the data on the next layer in target rof
191 mVrtParams[iteration].phiCut,
192 mTimeFrame->getTracklets()[0], // Flat tracklet buffer
193 mTimeFrame->getNTrackletsCluster(pivotRofId, 0), // Span of the number of tracklets per each cluster in pivot rof
194 mIndexTableUtils,
195 timeErr,
196 gsl::span<int>(), // Offset in the tracklet buffer
197 0,
198 0,
199 mVrtParams[iteration].maxTrackletsPerCluster);
200 }
201 const auto& rofRange12 = mTimeFrame->getROFOverlapTableView().getOverlap(1, 2, pivotRofId);
202 for (auto targetRofId = rofRange12.getFirstEntry(); targetRofId < rofRange12.getEntriesBound(); ++targetRofId) {
203 const auto timeErr = mTimeFrame->getROFOverlapTableView().getTimeStamp(2, targetRofId, 1, pivotRofId);
204 trackleterKernelHost<TrackletMode::Layer1Layer2, true>(
205 !skip ? mTimeFrame->getClustersOnLayer(targetRofId, 2) : gsl::span<Cluster>(),
206 !skip ? mTimeFrame->getClustersOnLayer(pivotRofId, 1) : gsl::span<Cluster>(),
207 mTimeFrame->getUsedClustersROF(targetRofId, 2),
208 mTimeFrame->getIndexTable(targetRofId, 2).data(),
209 mVrtParams[iteration].phiCut,
210 mTimeFrame->getTracklets()[1],
211 mTimeFrame->getNTrackletsCluster(pivotRofId, 1), // Span of the number of tracklets per each cluster in pivot rof
212 mIndexTableUtils,
213 timeErr,
214 gsl::span<int>(), // Offset in the tracklet buffer
215 0,
216 0,
217 mVrtParams[iteration].maxTrackletsPerCluster);
218 }
219 mTimeFrame->getNTrackletsROF(pivotRofId, 0) = std::accumulate(mTimeFrame->getNTrackletsCluster(pivotRofId, 0).begin(), mTimeFrame->getNTrackletsCluster(pivotRofId, 0).end(), 0);
220 mTimeFrame->getNTrackletsROF(pivotRofId, 1) = std::accumulate(mTimeFrame->getNTrackletsCluster(pivotRofId, 1).begin(), mTimeFrame->getNTrackletsCluster(pivotRofId, 1).end(), 0);
221 });
222
223 mTimeFrame->computeTrackletsPerROFScans();
224 if (auto tot0 = mTimeFrame->getTotalTrackletsTF(0), tot1 = mTimeFrame->getTotalTrackletsTF(1);
225 tot0 == 0 || tot1 == 0) {
226 return;
227 } else {
228 mTimeFrame->getTracklets()[0].resize(tot0);
229 mTimeFrame->getTracklets()[1].resize(tot1);
230 }
231
232 tbb::parallel_for(0, mTimeFrame->getNrof(1), [&](const short pivotRofId) {
233 bool skip = skipROF(iteration, pivotRofId);
234 const int globalOffsetPivot = mTimeFrame->getSortedStartIndex(pivotRofId, 1);
235 const auto& rofRange01 = mTimeFrame->getROFOverlapTableView().getOverlap(1, 0, pivotRofId);
236 for (auto targetRofId = rofRange01.getFirstEntry(); targetRofId < rofRange01.getEntriesBound(); ++targetRofId) {
237 const auto timeErr = mTimeFrame->getROFOverlapTableView().getTimeStamp(0, targetRofId, 1, pivotRofId);
238 trackleterKernelHost<TrackletMode::Layer0Layer1, false>(
239 !skip ? mTimeFrame->getClustersOnLayer(targetRofId, 0) : gsl::span<Cluster>(),
240 !skip ? mTimeFrame->getClustersOnLayer(pivotRofId, 1) : gsl::span<Cluster>(),
241 mTimeFrame->getUsedClustersROF(targetRofId, 0),
242 mTimeFrame->getIndexTable(targetRofId, 0).data(),
243 mVrtParams[iteration].phiCut,
244 mTimeFrame->getTracklets()[0],
245 mTimeFrame->getNTrackletsCluster(pivotRofId, 0),
246 mIndexTableUtils,
247 timeErr,
248 mTimeFrame->getExclusiveNTrackletsCluster(pivotRofId, 0),
249 mTimeFrame->getSortedStartIndex(targetRofId, 0),
250 globalOffsetPivot,
251 mVrtParams[iteration].maxTrackletsPerCluster);
252 }
253 const auto& rofRange12 = mTimeFrame->getROFOverlapTableView().getOverlap(1, 2, pivotRofId);
254 for (auto targetRofId = rofRange12.getFirstEntry(); targetRofId < rofRange12.getEntriesBound(); ++targetRofId) {
255 const auto timeErr = mTimeFrame->getROFOverlapTableView().getTimeStamp(2, targetRofId, 1, pivotRofId);
256 trackleterKernelHost<TrackletMode::Layer1Layer2, false>(
257 !skip ? mTimeFrame->getClustersOnLayer(targetRofId, 2) : gsl::span<Cluster>(),
258 !skip ? mTimeFrame->getClustersOnLayer(pivotRofId, 1) : gsl::span<Cluster>(),
259 mTimeFrame->getUsedClustersROF(targetRofId, 2),
260 mTimeFrame->getIndexTable(targetRofId, 2).data(),
261 mVrtParams[iteration].phiCut,
262 mTimeFrame->getTracklets()[1],
263 mTimeFrame->getNTrackletsCluster(pivotRofId, 1),
264 mIndexTableUtils,
265 timeErr,
266 mTimeFrame->getExclusiveNTrackletsCluster(pivotRofId, 1),
267 mTimeFrame->getSortedStartIndex(targetRofId, 2),
268 globalOffsetPivot,
269 mVrtParams[iteration].maxTrackletsPerCluster);
270 }
271 });
272 });
273
275 if (mTimeFrame->hasMCinformation()) {
276 for (const auto& trk : mTimeFrame->getTracklets()[0]) {
278 int sortedId0{trk.firstClusterIndex};
279 int sortedId1{trk.secondClusterIndex};
280 for (const auto& lab0 : mTimeFrame->getClusterLabels(0, mTimeFrame->getClusters()[0][sortedId0].clusterId)) {
281 for (const auto& lab1 : mTimeFrame->getClusterLabels(1, mTimeFrame->getClusters()[1][sortedId1].clusterId)) {
282 if (lab0 == lab1 && lab0.isValid()) {
283 label = lab0;
284 break;
285 }
286 }
287 if (label.isValid()) {
288 break;
289 }
290 }
291 mTimeFrame->getTrackletsLabel(0).emplace_back(label);
292 }
293 }
294}
295
296template <int NLayers>
298{
299 mTaskArena->execute([&] {
300 tbb::combinable<int> totalLines{0};
301 tbb::parallel_for(
302 tbb::blocked_range<short>(0, (short)mTimeFrame->getNrof(1)),
303 [&](const tbb::blocked_range<short>& Rofs) {
304 for (short pivotRofId = Rofs.begin(); pivotRofId < Rofs.end(); ++pivotRofId) {
305 if (mTimeFrame->getFoundTracklets(pivotRofId, 0).empty() || skipROF(iteration, pivotRofId)) {
306 continue;
307 }
308 mTimeFrame->getLines(pivotRofId).reserve(std::min(mTimeFrame->getFoundTracklets(pivotRofId, 0).size(), mTimeFrame->getNTrackletsCluster(pivotRofId, 0).size() * constants::MaxSelectedTrackletsPerCluster));
309 bounded_vector<uint8_t> usedTracklets(mTimeFrame->getFoundTracklets(pivotRofId, 0).size(), 0, mMemoryPool.get());
310 trackletSelectionKernelHost(
311 mTimeFrame->getClusters()[0].data(),
312 mTimeFrame->getClusters()[1].data(),
313 mTimeFrame->getUsedClusters(0),
314 mTimeFrame->getUsedClusters(2),
315 mTimeFrame->getFoundTracklets(pivotRofId, 0),
316 mTimeFrame->getFoundTracklets(pivotRofId, 1),
317 usedTracklets,
318 mTimeFrame->getNTrackletsCluster(pivotRofId, 0),
319 mTimeFrame->getNTrackletsCluster(pivotRofId, 1),
320 mTimeFrame->getLines(pivotRofId),
321 mTimeFrame->getLabelsFoundTracklets(pivotRofId, 0),
322 mTimeFrame->getLinesLabel(pivotRofId),
323 static_cast<int>(mTimeFrame->getClustersOnLayer(pivotRofId, 1).size()),
324 mVrtParams[iteration].tanLambdaCut,
325 mVrtParams[iteration].phiCut,
326 constants::MaxSelectedTrackletsPerCluster);
327 totalLines.local() += mTimeFrame->getLines(pivotRofId).size();
328 }
329 });
330 mTimeFrame->setNLinesTotal(totalLines.combine(std::plus<int>()));
331 });
332
333 // from here on we do not use tracklets anymore, so let's free them
334 deepVectorClear(mTimeFrame->getTracklets());
335}
336
337template <int NLayers>
339{
340 const int nRofs = mTimeFrame->getNrof(1);
341 std::vector<std::vector<Vertex>> rofVertices(nRofs);
342 std::vector<std::vector<VertexLabel>> rofLabels(nRofs);
343 const float pairCut2 = mVrtParams[iteration].pairCut * mVrtParams[iteration].pairCut;
344 const float duplicateZCut = mVrtParams[iteration].duplicateZCut > 0.f ? mVrtParams[iteration].duplicateZCut : std::max(4.f * mVrtParams[iteration].pairCut, 0.5f * mVrtParams[iteration].clusterCut);
345 const float duplicateDistance2Cut = mVrtParams[iteration].duplicateDistance2Cut > 0.f ? mVrtParams[iteration].duplicateDistance2Cut : std::max(16.f * pairCut2, 0.0625f * mVrtParams[iteration].clusterCut * mVrtParams[iteration].clusterCut);
347 settings.beamX = mTimeFrame->getBeamX();
348 settings.beamY = mTimeFrame->getBeamY();
349 settings.pairCut = mVrtParams[iteration].pairCut;
350 settings.pairCut2 = pairCut2;
351 settings.clusterCut = mVrtParams[iteration].clusterCut;
352 settings.coarseZWindow = mVrtParams[iteration].coarseZWindow;
353 settings.seedDedupZCut = mVrtParams[iteration].seedDedupZCut;
354 settings.refitDedupZCut = mVrtParams[iteration].refitDedupZCut;
355 settings.duplicateZCut = duplicateZCut;
356 settings.duplicateDistance2Cut = duplicateDistance2Cut;
357 settings.finalSelectionZCut = mVrtParams[iteration].finalSelectionZCut;
358 settings.maxZ = mVrtParams[iteration].maxZPositionAllowed;
359 settings.seedMemberRadiusTime = mVrtParams[iteration].seedMemberRadiusTime;
360 settings.seedMemberRadiusZ = mVrtParams[iteration].seedMemberRadiusZ;
361 settings.memoryPool = mMemoryPool;
362
363 const auto processROF = [&](const int rofId) {
364 if (skipROF(iteration, rofId)) {
365 return;
366 }
367 auto& lines = mTimeFrame->getLines(rofId);
368 auto clusters = line_vertexer::buildClusters(std::span<const Line>{lines.data(), lines.size()}, settings);
369 deepVectorClear(lines); // not needed after
370 auto clusterBeamDistance2 = [&](const ClusterLines& cluster) {
371 return (mTimeFrame->getBeamX() - cluster.getVertex()[0]) * (mTimeFrame->getBeamX() - cluster.getVertex()[0]) +
372 (mTimeFrame->getBeamY() - cluster.getVertex()[1]) * (mTimeFrame->getBeamY() - cluster.getVertex()[1]);
373 };
374 auto clusterBetter = [&](const ClusterLines& lhs, const ClusterLines& rhs) {
375 if (lhs.getSize() != rhs.getSize()) {
376 return lhs.getSize() > rhs.getSize();
377 }
378 if (o2::gpu::GPUCommonMath::Abs(lhs.getAvgDistance2() - rhs.getAvgDistance2()) > constants::Tolerance) {
379 return lhs.getAvgDistance2() < rhs.getAvgDistance2();
380 }
381 const auto lhsBeam = clusterBeamDistance2(lhs);
382 const auto rhsBeam = clusterBeamDistance2(rhs);
383 if (o2::gpu::GPUCommonMath::Abs(lhsBeam - rhsBeam) > constants::Tolerance) {
384 return lhsBeam < rhsBeam;
385 }
386 return lhs.getVertex()[2] < rhs.getVertex()[2];
387 };
388
389 // Cluster deduplication by local non-maximum suppression in time/space
390 std::sort(clusters.begin(), clusters.end(), clusterBetter);
391 float minClusterZ = std::numeric_limits<float>::max();
392 for (const auto& cluster : clusters) {
393 minClusterZ = std::min(minClusterZ, cluster.getVertex()[2]);
394 }
395 bounded_vector<ClusterLines> deduplicated(mMemoryPool.get());
396 deduplicated.reserve(clusters.size());
397 std::unordered_map<int, std::vector<int>> keptByZBin;
398 for (auto& candidate : clusters) {
399 bool duplicate = false;
400 const auto candidateZ = candidate.getVertex()[2];
401 const auto zBin = static_cast<int>(std::floor((candidateZ - minClusterZ) / settings.duplicateZCut));
402 for (int neighborBin = zBin - 1; neighborBin <= zBin + 1 && !duplicate; ++neighborBin) {
403 const auto found = keptByZBin.find(neighborBin);
404 if (found == keptByZBin.end()) {
405 continue;
406 }
407 for (const auto ownerId : found->second) {
408 const auto& owner = deduplicated[ownerId];
409 if (!candidate.getTimeStamp().isCompatible(owner.getTimeStamp())) {
410 continue;
411 }
412 if (o2::gpu::GPUCommonMath::Abs(candidate.getVertex()[2] - owner.getVertex()[2]) >= settings.duplicateZCut) {
413 continue;
414 }
415 const auto dx = candidate.getVertex()[0] - owner.getVertex()[0];
416 const auto dy = candidate.getVertex()[1] - owner.getVertex()[1];
417 const auto dz = candidate.getVertex()[2] - owner.getVertex()[2];
418 const auto distance2 = math_utils::SqSum(dx, dy, dz);
419 if (distance2 < settings.duplicateDistance2Cut) {
420 duplicate = true;
421 break;
422 }
423 }
424 }
425 if (duplicate) {
426 continue;
427 }
428
429 const auto ownerId = static_cast<int>(deduplicated.size());
430 keptByZBin[zBin].push_back(ownerId);
431 deduplicated.push_back(std::move(candidate));
432 }
433 clusters = std::move(deduplicated);
434 int nClusters = static_cast<int>(clusters.size());
435
436 // Vertex filtering with score-based local NMS
437 std::sort(clusters.begin(), clusters.end(), clusterBetter);
438 std::vector<int> candidateIndices;
439 candidateIndices.reserve(nClusters);
440 for (int iCluster{0}; iCluster < nClusters; ++iCluster) {
441 const bool zCompatible = o2::gpu::GPUCommonMath::Abs(clusters[iCluster].getVertex()[2]) < mVrtParams[iteration].maxZPositionAllowed;
442
443 if (zCompatible) {
444 candidateIndices.push_back(iCluster);
445 }
446 }
447
448 if (candidateIndices.empty()) {
449 return;
450 }
451
452 auto countSharedLabels = [](const ClusterLines& lhs, const ClusterLines& rhs) {
453 size_t shared = 0;
454 auto lhsIt = lhs.getLabels().begin();
455 auto rhsIt = rhs.getLabels().begin();
456 while (lhsIt != lhs.getLabels().end() && rhsIt != rhs.getLabels().end()) {
457 if (*lhsIt == *rhsIt) {
458 ++shared;
459 ++lhsIt;
460 ++rhsIt;
461 } else if (*lhsIt < *rhsIt) {
462 ++lhsIt;
463 } else {
464 ++rhsIt;
465 }
466 }
467 return shared;
468 };
469
470 float minCandidateZ = std::numeric_limits<float>::max();
471 for (const auto clusterId : candidateIndices) {
472 minCandidateZ = std::min(minCandidateZ, clusters[clusterId].getVertex()[2]);
473 }
474 std::unordered_map<int, std::vector<int>> selectedByZBin;
475 std::vector<int> selectedIndices;
476 selectedIndices.reserve(candidateIndices.size());
477 for (const auto clusterId : candidateIndices) {
478 const auto& candidate = clusters[clusterId];
479 const auto candidateZ = candidate.getVertex()[2];
480 const auto zBin = static_cast<int>((candidateZ - minCandidateZ) / settings.finalSelectionZCut);
481 bool suppressed = false;
482 for (int neighborBin = zBin - 1; neighborBin <= zBin + 1 && !suppressed; ++neighborBin) {
483 const auto found = selectedByZBin.find(neighborBin);
484 if (found == selectedByZBin.end()) {
485 continue;
486 }
487 for (const auto selectedId : found->second) {
488 const auto& selected = clusters[selectedId];
489 if (!candidate.getTimeStamp().isCompatible(selected.getTimeStamp())) {
490 continue;
491 }
492 const auto zDelta = o2::gpu::GPUCommonMath::Abs(candidateZ - selected.getVertex()[2]);
493 const auto sharedLabels = countSharedLabels(candidate, selected);
494 const auto minSize = std::min(candidate.getSize(), selected.getSize());
495 const bool overlapDuplicate = sharedLabels > 0 && sharedLabels * 4 >= minSize;
496 const bool strongZDuplicate = zDelta < settings.finalSelectionZCut;
497 const bool clearlyBetterMultiplicity = selected.getSize() >= candidate.getSize() + 3;
498 const bool clearlyBetterQuality = selected.getSize() > candidate.getSize() &&
499 selected.getAvgDistance2() + constants::Tolerance < 0.8f * candidate.getAvgDistance2();
500 const bool weakCandidate = clearlyBetterMultiplicity || clearlyBetterQuality;
501 if (overlapDuplicate || (strongZDuplicate && weakCandidate)) {
502 suppressed = true;
503 break;
504 }
505 }
506 }
507 if (suppressed) {
508 continue;
509 }
510 selectedByZBin[zBin].push_back(clusterId);
511 selectedIndices.push_back(clusterId);
512 }
513
514 // sort vertices by their multiplicity to opt. suppress lower mult. debris
515 std::vector<int> sortedIndices(selectedIndices.size());
516 std::iota(sortedIndices.begin(), sortedIndices.end(), 0);
517 std::sort(sortedIndices.begin(), sortedIndices.end(), [&selectedIndices, &clusters](int i, int j) {
518 return clusters[selectedIndices[i]].getSize() > clusters[selectedIndices[j]].getSize();
519 });
520 for (const auto sortedId : sortedIndices) {
521 const auto& cluster = clusters[selectedIndices[sortedId]];
522 const auto beamDistance2 = clusterBeamDistance2(cluster);
523 if (!(beamDistance2 < mVrtParams[iteration].NSigmaCut)) {
524 continue;
525 }
526 if (cluster.getSize() < mVrtParams[iteration].clusterContributorsCut) {
527 continue;
528 }
529 if (!rofVertices[rofId].empty() && cluster.getSize() < mVrtParams[iteration].suppressLowMultDebris) {
530 continue;
531 }
532
533 Vertex vertex{cluster.getVertex().data(),
534 cluster.getRMS2(),
535 (ushort)cluster.getSize(),
536 cluster.getAvgDistance2()};
537 if (mVrtParams[iteration].PassFlags[IterationStep::MarkVerticesAsUPC]) {
538 vertex.setFlags(Vertex::UPCMode);
539 }
540 vertex.setTimeStamp(cluster.getTimeStamp());
541 rofVertices[rofId].push_back(vertex);
542 if (mTimeFrame->hasMCinformation()) {
543 auto& lineLabels = mTimeFrame->getLinesLabel(rofId);
544 bounded_vector<o2::MCCompLabel> labels(mMemoryPool.get());
545 for (auto& index : cluster.getLabels()) {
546 labels.push_back(lineLabels[index]);
547 }
548 const auto mainLabel = computeMain(labels);
549 rofLabels[rofId].push_back(mainLabel);
550 }
551 }
552 };
553
554 if (mTaskArena->max_concurrency() <= 1) {
555 for (int rofId{0}; rofId < nRofs; ++rofId) {
556 processROF(rofId);
557 }
558 } else {
559 mTaskArena->execute([&] {
560 tbb::parallel_for(0, nRofs, [&](const int rofId) {
561 processROF(rofId);
562 });
563 });
564 }
565 // add vertices, these anyways get sorted afterward
566 for (int rofId{0}; rofId < nRofs; ++rofId) {
567 for (auto& vertex : rofVertices[rofId]) {
568 mTimeFrame->addPrimaryVertex(vertex);
569 }
570 if (mTimeFrame->hasMCinformation()) {
571 for (auto& label : rofLabels[rofId]) {
572 mTimeFrame->addPrimaryVertexLabel(label);
573 }
574 }
575 }
576}
577
578template <int NLayers>
580{
581 LOGP(info, "Using truth seeds as vertices; will skip computations");
582 const auto dc = o2::steer::DigitizationContext::loadFromFile("collisioncontext.root");
583 const auto irs = dc->getEventRecords();
585 int64_t roFrameLengthInBC = o2::itsmft::DPLAlpideParam<o2::detectors::DetID::ITS>::Instance().getROFLengthInBC(1);
587 const int iSrc = 0; // take only events from collision generator
588 auto eveId2colId = dc->getCollisionIndicesForSource(iSrc);
589 for (int iEve{0}; iEve < mcReader.getNEvents(iSrc); ++iEve) {
590 const auto& ir = irs[eveId2colId[iEve]];
591 if (!ir.isDummy()) { // do we need this, is this for diffractive events?
592 const auto& eve = mcReader.getMCEventHeader(iSrc, iEve);
593 auto bc = (ir - raw::HBFUtils::Instance().getFirstSampledTFIR()).toLong() - roFrameBiasInBC;
594 if (bc < 0) { // event happened before TF
595 continue;
596 }
597 Vertex vert;
598 vert.getTimeStamp().setTimeStamp(bc);
599 vert.getTimeStamp().setTimeStampError(roFrameLengthInBC / 2);
600 // set minimum to 1 sometimes for diffractive events there is nothing acceptance
601 vert.setNContributors(std::max(1L, std::ranges::count_if(mcReader.getTracks(iSrc, iEve), [](const auto& trk) {
602 if (!trk.isPrimary() || trk.GetPt() < 0.05 || std::abs(trk.GetEta()) > 1.1) {
603 return false;
604 }
605 const auto* p = o2::O2DatabasePDG::Instance()->GetParticle(trk.GetPdgCode());
606 return (!p) ? false : p->Charge() != 0;
607 })));
608 vert.setXYZ((float)eve.GetX(), (float)eve.GetY(), (float)eve.GetZ());
609 vert.setChi2(1); // not used as constraint
610 constexpr float cov = 25e-4;
611 vert.setSigmaX(cov);
612 vert.setSigmaY(cov);
613 vert.setSigmaZ(cov);
614 mTimeFrame->addPrimaryVertex(vert);
615 o2::MCCompLabel mcLbl(o2::MCCompLabel::maxTrackID(), iEve, iSrc, false);
616 VertexLabel lbl(mcLbl, 1.0);
617 mTimeFrame->addPrimaryVertexLabel(lbl);
618 }
619 mcReader.releaseTracksForSourceAndEvent(iSrc, iEve);
620 }
621 LOGP(info, "Imposed {} pv collisions from mc-truth", mTimeFrame->getPrimaryVertices().size());
622}
623
624template <int NLayers>
625void VertexerTraits<NLayers>::setNThreads(int n, std::shared_ptr<tbb::task_arena>& arena)
626{
627 if (arena == nullptr) {
628 mTaskArena = std::make_shared<tbb::task_arena>(std::abs(n));
629 LOGP(info, "Setting seeding vertexer with {} threads.", n);
630 } else {
631 mTaskArena = arena;
632 }
633}
634
635template <int NLayers>
636bool VertexerTraits<NLayers>::skipROF(int iteration, int rof) const
637{
638 return mVrtParams[iteration].PassFlags[IterationStep::SkipROFsAboveThreshold] &&
639 (int)mTimeFrame->getROFVertexLookupTableView().getVertices(1, rof).getEntries() > mVrtParams[iteration].vertPerRofThreshold;
640}
641
642template class VertexerTraits<7>;
643} // namespace o2::its
std::vector< std::string > labels
size_t minSize
uint64_t vertex
Definition RawEventData.h:9
uint64_t bc
Definition RawEventData.h:5
int32_t i
Mode
Definition Utils.h:89
uint32_t j
Definition RawData.h:0
Class to compute the primary vertex in ITS from tracklets.
int nClusters
static constexpr int maxTrackID()
static TDatabasePDG * Instance()
HMPID cluster implementation.
Definition Cluster.h:27
virtual void initialise(const TrackingParameters &trackingParams)
virtual void updateVertexingParameters(const std::vector< VertexingParameters > &vrtPar)
virtual void computeTracklets(const int iteration)
static DigitizationContext * loadFromFile(std::string_view filename="")
size_t getNEvents(int source) const
Get number of events.
o2::dataformats::MCEventHeader const & getMCEventHeader(int source, int event) const
retrieves the MCEventHeader for a given eventID and sourceID
void releaseTracksForSourceAndEvent(int source, int event)
API to ask releasing tracks (freeing memory) for source + event.
std::vector< MCTrack > const & getTracks(int source, int event) const
variant returning all tracks for source and event at once
GLdouble n
Definition glcorearb.h:1982
GLuint index
Definition glcorearb.h:781
GLuint GLsizei const GLchar * label
Definition glcorearb.h:2519
constexpr float TwoPI
void deepVectorClear(std::vector< T > &vec)
std::pmr::vector< T > bounded_vector
std::pair< o2::MCCompLabel, float > VertexLabel
Definition Vertex.h:27
std::vector< Cluster > getClusters(int event)
Common utility functions.
void empty(int)
int32_t w
std::shared_ptr< BoundedMemoryResource > memoryPool
o2::InteractionRecord ir(0, 0)
std::vector< Cluster > clusters
std::vector< Tracklet64 > tracklets