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