Project
Loading...
Searching...
No Matches
LineVertexerHelpers.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.
11
12#include <algorithm>
13#include <array>
14#include <cmath>
15#include <cstdint>
16#include <limits>
17#include <numeric>
18#include <queue>
19#include <unordered_map>
20#include <utility>
21#include <vector>
22
23#include <Math/SMatrix.h>
24#include <Math/SVector.h>
25
29
31{
32namespace
33{
35using SVector3 = ROOT::Math::SVector<float, 3>;
36
37constexpr float TukeyC = 4.685f;
38constexpr float TukeyC2 = TukeyC * TukeyC;
39constexpr float InitialScale2 = 5.f;
40constexpr float MinScale2 = 1.f;
41constexpr float MedianToSigma = 1.4826f;
42constexpr float VertexShiftZTol = 0.01f;
43constexpr float VertexShiftR2Tol = 1.e-4f;
44constexpr int MaxFitIterations = 10;
45constexpr int MaxSeedsPerCluster = 32;
46constexpr float MinRelativePeakSupport = 0.1f;
47constexpr int MaxHistogramBins = 0x7fff;
48constexpr float TieTolerance = 1e-5f;
49
50struct LineRef {
51 LineRef(const Line& line, const int index, const float beamX, const float beamY, const float maxZ) : lineIndex(index)
52 {
53 const auto symTime = line.mTime.makeSymmetrical();
54 tCenter = symTime.getTimeStamp();
55 tHalfWidth = symTime.getTimeStampError();
56 const auto dx = line.originPoint(0) - beamX;
57 const auto dy = line.originPoint(1) - beamY;
58 const auto ux = line.cosinesDirector(0);
59 const auto uy = line.cosinesDirector(1);
60 const auto uz = line.cosinesDirector(2);
61 const auto den = math_utils::SqSum(ux, uy);
62 if (den <= constants::Tolerance) {
63 lineIndex = constants::UnusedIndex;
64 return;
65 }
66 const auto s0 = -((dx * ux) + (dy * uy)) / den;
67 const auto xb = dx + (s0 * ux);
68 const auto yb = dy + (s0 * uy);
69 zBeam = line.originPoint(2) + s0 * uz;
70 if (!std::isfinite(zBeam) || o2::gpu::CAMath::Abs(zBeam) > maxZ) {
71 lineIndex = constants::UnusedIndex;
72 }
73 }
74 bool isDead() const noexcept { return lineIndex == constants::UnusedIndex; }
75
76 int lineIndex = constants::UnusedIndex;
77 float zBeam = 0.f;
78 float tCenter = 0.f;
79 float tHalfWidth = 0.f;
80};
81
82struct VertexSeed {
83 explicit VertexSeed(const std::shared_ptr<BoundedMemoryResource>& mr) : contributors(mr.get()), assigned(mr.get()) {}
84
85 std::array<float, 3> vertex = {};
86 TimeEstBC time;
87 float scale2 = InitialScale2;
88 bounded_vector<int> contributors;
89 bounded_vector<int> assigned;
90 bool valid = false;
91 bool isUsableSeed() const noexcept
92 {
93 return valid && contributors.size() >= 2;
94 }
95};
96
97void compactSeeds(bounded_vector<VertexSeed>& seeds)
98{
99 seeds.erase(std::remove_if(seeds.begin(), seeds.end(), [](const VertexSeed& seed) {
100 return !seed.isUsableSeed();
101 }),
102 seeds.end());
103}
104
105struct Histogram2D {
106 explicit Histogram2D(const std::shared_ptr<BoundedMemoryResource>& mr) : bins(mr.get()) {}
107
108 int nTimeBins = 0;
109 int nZBins = 0;
110 float timeMin = 0.f;
111 float zMin = 0.f;
112 float timeBinSize = 1.f;
113 float zBinSize = 1.f;
114 bounded_vector<float> bins;
115
116 int getIndex(const int tBin, const int zBin) const noexcept
117 {
118 return (tBin * nZBins) + zBin;
119 }
120
121 std::pair<int, int> decodeIndex(const int index) const noexcept
122 {
123 return {index / nZBins, index % nZBins};
124 }
125
126 int getTimeBin(const float time) const noexcept
127 {
128 if (time < timeMin) {
129 return -1;
130 }
131 const auto bin = static_cast<int>((time - timeMin) / timeBinSize);
132 return (bin >= 0 && bin < nTimeBins) ? bin : -1;
133 }
134
135 int getZBin(const float z) const noexcept
136 {
137 if (z < zMin) {
138 return -1;
139 }
140 const auto bin = static_cast<int>((z - zMin) / zBinSize);
141 return (bin >= 0 && bin < nZBins) ? bin : -1;
142 }
143
144 void fill(const float time, const float z, const float weight) noexcept
145 {
146 const auto tBin = getTimeBin(time);
147 const auto zBin = getZBin(z);
148 if (tBin < 0 || zBin < 0) {
149 return;
150 }
151 bins[getIndex(tBin, zBin)] += weight;
152 }
153
154 int findPeakBin() const noexcept
155 {
156 float bestWeight = 0.f;
157 int bestIndex = -1;
158 for (int index = 0; index < static_cast<int>(bins.size()); ++index) {
159 if (bins[index] > bestWeight) {
160 bestWeight = bins[index];
161 bestIndex = index;
162 }
163 }
164 return bestIndex;
165 }
166
167 void suppressBin(const int index) noexcept
168 {
169 if (index >= 0 && index < static_cast<int>(bins.size())) {
170 bins[index] = -1.f;
171 }
172 }
173
174 void suppressNeighborhood(const int index, const int radiusTime, const int radiusZ) noexcept
175 {
176 if (index < 0) {
177 return;
178 }
179 const auto [tBin, zBin] = decodeIndex(index);
180 for (int dt = -radiusTime; dt <= radiusTime; ++dt) {
181 const auto tt = tBin + dt;
182 if (tt < 0 || tt >= nTimeBins) {
183 continue;
184 }
185 for (int dz = -radiusZ; dz <= radiusZ; ++dz) {
186 const auto zz = zBin + dz;
187 if (zz < 0 || zz >= nZBins) {
188 continue;
189 }
190 bins[getIndex(tt, zz)] = -1.f;
191 }
192 }
193 }
194
195 float getNeighborhoodSum(const int index, const int radiusTime, const int radiusZ) const noexcept
196 {
197 if (index < 0) {
198 return 0.f;
199 }
200 const auto [tBin, zBin] = decodeIndex(index);
201 float sum = 0.f;
202 for (int dt = -radiusTime; dt <= radiusTime; ++dt) {
203 const auto tt = tBin + dt;
204 if (tt < 0 || tt >= nTimeBins) {
205 continue;
206 }
207 for (int dz = -radiusZ; dz <= radiusZ; ++dz) {
208 const auto zz = zBin + dz;
209 if (zz < 0 || zz >= nZBins) {
210 continue;
211 }
212 const auto value = bins[getIndex(tt, zz)];
213 if (value > 0.f) {
214 sum += value;
215 }
216 }
217 }
218 return sum;
219 }
220
221 float getTimeBinCenter(const int tBin) const noexcept
222 {
223 return timeMin + ((static_cast<float>(tBin) + 0.5f) * timeBinSize);
224 }
225
226 float getZBinCenter(const int zBin) const noexcept
227 {
228 return zMin + ((static_cast<float>(zBin) + 0.5f) * zBinSize);
229 }
230
231 TimeEstBC getTimeInterval(const int tBin) const noexcept
232 {
233 const auto lowFloat = timeMin + (static_cast<float>(tBin) * timeBinSize);
234 const auto highFloat = lowFloat + timeBinSize;
235 const auto low = std::max<double>(0., std::floor(lowFloat));
236 const auto high = std::max(low + 1., (double)std::ceil(highFloat));
237 constexpr auto maxTS = std::numeric_limits<TimeStampType>::max();
238 const auto clampedLow = std::min<double>(low, maxTS - 1.);
239 const auto width = std::min<double>(high - clampedLow, std::numeric_limits<TimeStampErrorType>::max());
240 return {static_cast<TimeStampType>(clampedLow), static_cast<TimeStampErrorType>(std::max<double>(1., width))};
241 }
242
243 TimeEstBC getTimeNeighborhoodInterval(const int tBin, const int radius) const noexcept
244 {
245 const auto lowBin = std::max(0, tBin - radius);
246 const auto highBin = std::min(nTimeBins - 1, tBin + radius);
247 const auto lowFloat = timeMin + (static_cast<float>(lowBin) * timeBinSize);
248 const auto highFloat = timeMin + (static_cast<float>(highBin + 1) * timeBinSize);
249 const auto low = std::max<double>(0., std::floor(lowFloat));
250 const auto high = std::max(low + 1., (double)std::ceil(highFloat));
251 constexpr auto maxTS = std::numeric_limits<TimeStampType>::max();
252 const auto clampedLow = std::min<double>(low, maxTS - 1.);
253 const auto width = std::min<double>(high - clampedLow, std::numeric_limits<TimeStampErrorType>::max());
254 return {static_cast<TimeStampType>(clampedLow), static_cast<TimeStampErrorType>(std::max<double>(1., width))};
255 }
256};
257
258class SeedHistogram
259{
260 public:
261 SeedHistogram(std::span<const int> members,
262 std::span<const LineRef> lineRefs,
263 std::span<const Line> lines,
264 const Settings& settings)
265 : mMembers(members), mLineRefs(lineRefs), mSeedMemberRadiusTime(settings.seedMemberRadiusTime), mSeedMemberRadiusZ(settings.seedMemberRadiusZ), mMemoryPool(settings.memoryPool), mHistogram(mMemoryPool)
266 {
267 const auto zBinSize = 0.25f * settings.clusterCut;
268 const auto timeBinSize = medianTimeError(lines);
269
270 float minZ = std::numeric_limits<float>::max();
271 float maxZ = std::numeric_limits<float>::lowest();
272 float minTime = std::numeric_limits<float>::max();
273 float maxTime = std::numeric_limits<float>::lowest();
274 for (const auto lineRefIdx : mMembers) {
275 minZ = std::min(minZ, mLineRefs[lineRefIdx].zBeam);
276 maxZ = std::max(maxZ, mLineRefs[lineRefIdx].zBeam);
277 minTime = std::min(minTime, mLineRefs[lineRefIdx].tCenter);
278 maxTime = std::max(maxTime, mLineRefs[lineRefIdx].tCenter);
279 }
280
281 const auto dz = std::max(0.f, maxZ - minZ);
282 const auto dt = std::max(0.f, maxTime - minTime);
283 mHistogram.nZBins = 1 + static_cast<int>(dz / zBinSize);
284 mHistogram.nTimeBins = 1 + static_cast<int>(dt / timeBinSize);
285 if (mHistogram.nTimeBins * mHistogram.nZBins > MaxHistogramBins) {
286 if (mHistogram.nTimeBins > mHistogram.nZBins) {
287 mHistogram.nTimeBins = std::max(1, (MaxHistogramBins - 1) / std::max(1, mHistogram.nZBins));
288 } else {
289 mHistogram.nZBins = std::max(1, (MaxHistogramBins - 1) / std::max(1, mHistogram.nTimeBins));
290 }
291 }
292
293 mHistogram.timeBinSize = std::max(timeBinSize, dt / (float)std::max(1, mHistogram.nTimeBins));
294 mHistogram.zBinSize = std::max(zBinSize, dz / (float)std::max(1, mHistogram.nZBins));
295 const auto paddedTime = 0.5f * ((float)mHistogram.nTimeBins * mHistogram.timeBinSize - dt);
296 const auto paddedZ = 0.5f * ((float)mHistogram.nZBins * mHistogram.zBinSize - dz);
297 mHistogram.timeMin = minTime - paddedTime;
298 mHistogram.zMin = minZ - paddedZ;
299 mHistogram.bins.assign((size_t)mHistogram.nTimeBins * (size_t)mHistogram.nZBins, 0.f);
300
301 for (const auto lineRefIdx : mMembers) {
302 mHistogram.fill(mLineRefs[lineRefIdx].tCenter, mLineRefs[lineRefIdx].zBeam, 1.f);
303 }
304 }
305
306 int findPeakBin() const noexcept
307 {
308 return mHistogram.findPeakBin();
309 }
310
311 float getPeakSupport(const int peakIndex) const noexcept
312 {
313 return mHistogram.getNeighborhoodSum(peakIndex, mSeedMemberRadiusTime, mSeedMemberRadiusZ);
314 }
315
316 bounded_vector<int> collectLocalMembers(const int peakIndex, const int radiusTime, const int radiusZ) const
317 {
318 bounded_vector<int> localMembers(mMemoryPool.get());
319 localMembers.reserve(mMembers.size());
320 const auto [timeBin, zBin] = mHistogram.decodeIndex(peakIndex);
321 for (const auto lineRefIdx : mMembers) {
322 const auto memberTimeBin = mHistogram.getTimeBin(mLineRefs[lineRefIdx].tCenter);
323 const auto memberZBin = mHistogram.getZBin(mLineRefs[lineRefIdx].zBeam);
324 if (memberTimeBin < 0 || memberZBin < 0) {
325 continue;
326 }
327 if (o2::gpu::GPUCommonMath::Abs(memberTimeBin - timeBin) > radiusTime) {
328 continue;
329 }
330 if (o2::gpu::GPUCommonMath::Abs(memberZBin - zBin) > radiusZ) {
331 continue;
332 }
333 localMembers.push_back(lineRefIdx);
334 }
335 return localMembers;
336 }
337
338 TimeEstBC getPeakTimeInterval(const int peakIndex, const int radius = 0) const noexcept
339 {
340 return mHistogram.getTimeNeighborhoodInterval(mHistogram.decodeIndex(peakIndex).first, radius);
341 }
342
343 float getPeakZCenter(const int peakIndex) const noexcept
344 {
345 return mHistogram.getZBinCenter(mHistogram.decodeIndex(peakIndex).second);
346 }
347
348 void suppressPeak(const int peakIndex) noexcept
349 {
350 mHistogram.suppressBin(peakIndex);
351 }
352
353 void suppressPeakNeighborhood(const int peakIndex) noexcept
354 {
355 mHistogram.suppressNeighborhood(peakIndex, mSeedMemberRadiusTime, mSeedMemberRadiusZ);
356 }
357
358 private:
359 float medianTimeError(std::span<const Line> lines) const
360 {
361 bounded_vector<float> errors(mMemoryPool.get());
362 errors.reserve(mMembers.size());
363 for (const auto lineRefIdx : mMembers) {
364 errors.push_back(static_cast<float>(lines[mLineRefs[lineRefIdx].lineIndex].mTime.getTimeStampError()));
365 }
366 std::sort(errors.begin(), errors.end());
367 return errors.empty() ? 1.f : std::max(1.f, errors[errors.size() / 2]);
368 }
369
370 std::span<const int> mMembers;
371 std::span<const LineRef> mLineRefs;
372 int mSeedMemberRadiusTime = 1;
373 int mSeedMemberRadiusZ = 2;
374 std::shared_ptr<BoundedMemoryResource> mMemoryPool;
375 Histogram2D mHistogram;
376};
377
378float updateScale2(const std::span<const float> chi2s, const std::shared_ptr<BoundedMemoryResource>& mr) noexcept
379{
380 if (chi2s.empty()) {
381 return MinScale2;
382 }
383
384 bounded_vector<float> sorted(chi2s.begin(), chi2s.end(), mr.get());
385 std::sort(sorted.begin(), sorted.end());
386 const auto median = sorted[sorted.size() / 2];
387
388 for (auto& value : sorted) {
389 value = o2::gpu::GPUCommonMath::Abs(value - median);
390 }
391 std::sort(sorted.begin(), sorted.end());
392 const auto mad = sorted[sorted.size() / 2];
393 if (!std::isfinite(mad) || mad <= constants::Tolerance) {
394 return MinScale2;
395 }
396 return std::max(MinScale2, MedianToSigma * mad);
397}
398
399class VertexFit
400{
401 public:
402 void add(const Line& line, const float weight) noexcept
403 {
404 const auto& direction = line.cosinesDirector;
405 const auto& origin = line.originPoint;
406 const auto det = ROOT::Math::Dot(direction, direction);
407 if (det <= constants::Tolerance) {
408 return;
409 }
410
411 for (int i = 0; i < 3; ++i) {
412 for (int j = i; j < 3; ++j) {
413 mMatrix(i, j) += weight * (((i == j ? det : 0.f) - direction(i) * direction(j)) / det);
414 }
415 }
416
417 const auto dDotO = ROOT::Math::Dot(direction, origin);
418 for (int i = 0; i < 3; ++i) {
419 mRhs(i) += weight * ((direction(i) * dDotO - det * origin(i)) / det);
420 }
421 }
422
423 bool solve(std::array<float, 3>& vertexOut) const noexcept
424 {
425 SymMatrix3 inv{mMatrix};
426 if (!inv.InvertFast()) {
427 return false;
428 }
429 const auto solution = inv * mRhs;
430 vertexOut[0] = static_cast<float>(-solution(0));
431 vertexOut[1] = static_cast<float>(-solution(1));
432 vertexOut[2] = static_cast<float>(-solution(2));
433 return std::isfinite(vertexOut[0]) && std::isfinite(vertexOut[1]) && std::isfinite(vertexOut[2]);
434 }
435
436 private:
437 SymMatrix3 mMatrix;
438 SVector3 mRhs;
439};
440
441VertexSeed fitSeed(const VertexSeed& initialSeed,
442 std::span<const int> members,
443 std::span<const LineRef> lineRefs,
444 std::span<const Line> lines,
445 const std::shared_ptr<BoundedMemoryResource>& mr,
446 const float pairCut2)
447{
448 VertexSeed seed{mr};
449 seed.vertex = initialSeed.vertex;
450 seed.time = initialSeed.time;
451 seed.scale2 = initialSeed.scale2;
452 seed.valid = false;
453 seed.contributors.clear();
454 seed.assigned.clear();
455 if (members.size() < 2) {
456 return seed;
457 }
458
459 for (int iteration = 0; iteration < MaxFitIterations; ++iteration) {
460 VertexFit vertexFit;
461 TimeEstBC commonTime{};
462 bool hasCommonTime = false;
463 bounded_vector<int> contributors{mr.get()};
464 const auto scale2 = std::max(seed.scale2, MinScale2);
465 const auto tukeyFactor = 1.f / (scale2 * TukeyC2);
466
467 for (const auto lineRefIdx : members) {
468 const auto lineIdx = lineRefs[lineRefIdx].lineIndex;
469 const auto& line = lines[lineIdx];
470 if (!line.mTime.isCompatible(seed.time)) {
471 continue;
472 }
473 if (hasCommonTime && !line.mTime.isCompatible(commonTime)) {
474 continue;
475 }
476
477 const auto chi2 = Line::getDistance2FromPoint(line, seed.vertex) / pairCut2;
478 auto weight = 1.f - (chi2 * tukeyFactor);
479 if (weight <= 0.f) {
480 continue;
481 }
482 weight *= weight;
483
484 if (!hasCommonTime) {
485 commonTime = line.mTime;
486 hasCommonTime = true;
487 } else {
488 commonTime += line.mTime;
489 }
490
491 contributors.push_back(lineRefIdx);
492 vertexFit.add(line, weight);
493 }
494
495 if (!hasCommonTime || contributors.size() < 2) {
496 return seed;
497 }
498
499 std::sort(contributors.begin(), contributors.end());
500
501 std::array<float, 3> updatedVertex{};
502 if (!vertexFit.solve(updatedVertex)) {
503 return seed;
504 }
505
506 const auto sameContributors = contributors == seed.contributors;
507 const auto dz = o2::gpu::GPUCommonMath::Abs(updatedVertex[2] - seed.vertex[2]);
508 const auto oldR2 = (seed.vertex[0] * seed.vertex[0]) + (seed.vertex[1] * seed.vertex[1]);
509 const auto newR2 = (updatedVertex[0] * updatedVertex[0]) + (updatedVertex[1] * updatedVertex[1]);
510 const auto dr2 = o2::gpu::GPUCommonMath::Abs(newR2 - oldR2);
511
512 seed.vertex = updatedVertex;
513 seed.time = commonTime;
514 bounded_vector<float> updatedChi2s{mr.get()};
515 updatedChi2s.reserve(contributors.size());
516 for (const auto lineRefIx : contributors) {
517 updatedChi2s.push_back(Line::getDistance2FromPoint(lines[lineRefs[lineRefIx].lineIndex], seed.vertex) / pairCut2);
518 }
519 seed.scale2 = updateScale2(updatedChi2s, mr);
520 seed.contributors = std::move(contributors);
521 seed.valid = true;
522
523 if (sameContributors && dz < VertexShiftZTol && dr2 < VertexShiftR2Tol) {
524 break;
525 }
526 }
527
528 return seed;
529}
530
531size_t countSharedContributors(std::span<const int> lhs, std::span<const int> rhs) noexcept
532{
533 size_t shared = 0;
534 auto lhsIt = lhs.begin();
535 auto rhsIt = rhs.begin();
536 while (lhsIt != lhs.end() && rhsIt != rhs.end()) {
537 if (*lhsIt == *rhsIt) {
538 ++shared;
539 ++lhsIt;
540 ++rhsIt;
541 } else if (*lhsIt < *rhsIt) {
542 ++lhsIt;
543 } else {
544 ++rhsIt;
545 }
546 }
547 return shared;
548}
549
550bounded_vector<int> collectCompatibleContributors(const VertexSeed& seed,
551 std::span<const int> members,
552 std::span<const LineRef> lineRefs,
553 std::span<const Line> lines,
554 const std::shared_ptr<BoundedMemoryResource>& mr,
555 const float pairCut2)
556{
557 bounded_vector<int> contributors{mr.get()};
558 contributors.reserve(members.size());
559 for (const auto lineRefIdx : members) {
560 const auto lineIdx = lineRefs[lineRefIdx].lineIndex;
561 const auto& line = lines[lineIdx];
562 if (!line.mTime.isCompatible(seed.time)) {
563 continue;
564 }
565 if (Line::getDistance2FromPoint(line, seed.vertex) >= pairCut2) {
566 continue;
567 }
568 contributors.push_back(lineRefIdx);
569 }
570 std::sort(contributors.begin(), contributors.end());
571 return contributors;
572}
573
574void deduplicateSeeds(bounded_vector<VertexSeed>& seeds, const Settings& settings)
575{
576 if (seeds.size() < 2) {
577 return;
578 }
579
580 std::sort(seeds.begin(), seeds.end(), [](const VertexSeed& lhs, const VertexSeed& rhs) {
581 if (lhs.contributors.size() != rhs.contributors.size()) {
582 return lhs.contributors.size() > rhs.contributors.size();
583 }
584 if (o2::gpu::GPUCommonMath::Abs(lhs.scale2 - rhs.scale2) > constants::Tolerance) {
585 return lhs.scale2 < rhs.scale2;
586 }
587 return lhs.vertex[2] < rhs.vertex[2];
588 });
589
590 const auto dedupZCut = settings.seedDedupZCut > 0.f ? settings.seedDedupZCut : 0.25f * settings.clusterCut;
591 for (size_t i = 0; i < seeds.size(); ++i) {
592 auto& candidate = seeds[i];
593 if (!candidate.isUsableSeed()) {
594 candidate.valid = false;
595 continue;
596 }
597 bool duplicate = false;
598 for (size_t j = 0; j < i; ++j) {
599 const auto& kept = seeds[j];
600 if (!kept.isUsableSeed()) {
601 continue;
602 }
603 if (!candidate.time.isCompatible(kept.time)) {
604 continue;
605 }
606 const auto shared = countSharedContributors(candidate.contributors, kept.contributors);
607 const auto minSize = std::min(candidate.contributors.size(), kept.contributors.size());
608 const auto zDelta = o2::gpu::GPUCommonMath::Abs(candidate.vertex[2] - kept.vertex[2]);
609 const bool clearlyWorse = kept.contributors.size() > candidate.contributors.size() ||
610 kept.scale2 + constants::Tolerance < 0.9f * candidate.scale2;
611 const bool overlapDuplicate = shared > 0 && shared * 2 >= minSize;
612 const bool nearbyDuplicate = zDelta < dedupZCut && (shared > 0 || clearlyWorse);
613 if (overlapDuplicate || nearbyDuplicate) {
614 duplicate = true;
615 break;
616 }
617 }
618 if (duplicate) {
619 candidate.valid = false;
620 }
621 }
622 compactSeeds(seeds);
623}
624
625void deduplicateRefittedSeeds(bounded_vector<VertexSeed>& seeds, const Settings& settings)
626{
627 if (seeds.size() < 2) {
628 return;
629 }
630
631 std::sort(seeds.begin(), seeds.end(), [](const VertexSeed& lhs, const VertexSeed& rhs) {
632 if (lhs.contributors.size() != rhs.contributors.size()) {
633 return lhs.contributors.size() > rhs.contributors.size();
634 }
635 if (o2::gpu::GPUCommonMath::Abs(lhs.scale2 - rhs.scale2) > constants::Tolerance) {
636 return lhs.scale2 < rhs.scale2;
637 }
638 return lhs.vertex[2] < rhs.vertex[2];
639 });
640
641 const auto zCut = settings.refitDedupZCut > 0.f ? settings.refitDedupZCut : 0.25f * settings.clusterCut;
642 for (size_t i = 0; i < seeds.size(); ++i) {
643 auto& candidate = seeds[i];
644 if (!candidate.isUsableSeed()) {
645 candidate.valid = false;
646 continue;
647 }
648 bool duplicate = false;
649 for (size_t j = 0; j < i; ++j) {
650 const auto& kept = seeds[j];
651 if (!kept.isUsableSeed()) {
652 continue;
653 }
654 if (!candidate.time.isCompatible(kept.time)) {
655 continue;
656 }
657 const auto shared = countSharedContributors(candidate.contributors, kept.contributors);
658 const auto minSize = std::min(candidate.contributors.size(), kept.contributors.size());
659 const auto zDelta = o2::gpu::GPUCommonMath::Abs(candidate.vertex[2] - kept.vertex[2]);
660 const bool overlapDuplicate = shared > 0 && shared * 2 >= minSize;
661 const bool lowSupportPair = std::min(candidate.contributors.size(), kept.contributors.size()) < 4;
662 const bool clearlyWorse = kept.contributors.size() > candidate.contributors.size() ||
663 kept.scale2 + constants::Tolerance < 0.9f * candidate.scale2;
664 const bool geometricDuplicate = zDelta < zCut && (lowSupportPair || clearlyWorse);
665 if (overlapDuplicate || geometricDuplicate) {
666 duplicate = true;
667 break;
668 }
669 }
670 if (duplicate) {
671 candidate.valid = false;
672 }
673 }
674 compactSeeds(seeds);
675}
676
677struct OrderedComponent {
678 explicit OrderedComponent(const std::shared_ptr<BoundedMemoryResource>& mr) : members(mr.get()) {}
679 float center = 0.f;
680 bounded_vector<int> members;
681};
682
683bounded_vector<bounded_vector<int>> buildCoarseClusters(std::span<const LineRef> lineRefs,
684 std::span<const Line> lines,
685 const Settings& settings)
686{
687 bounded_vector<bounded_vector<int>> clusters(settings.memoryPool.get());
688 if (lineRefs.size() < 2) {
689 return clusters;
690 }
691
692 bounded_vector<int> sortedByLower(lineRefs.size(), settings.memoryPool.get());
693 std::iota(sortedByLower.begin(), sortedByLower.end(), 0);
694 std::sort(sortedByLower.begin(), sortedByLower.end(), [&](const int lhs, const int rhs) {
695 const auto lhsLower = lines[lineRefs[lhs].lineIndex].mTime.lower();
696 const auto rhsLower = lines[lineRefs[rhs].lineIndex].mTime.lower();
697 if (lhsLower != rhsLower) {
698 return lhsLower < rhsLower;
699 }
700 return lineRefs[lhs].lineIndex < lineRefs[rhs].lineIndex;
701 });
702
703 const auto coarseZWindow = settings.coarseZWindow > 0.f ? settings.coarseZWindow : settings.clusterCut;
704 bounded_vector<int> parent(lineRefs.size(), settings.memoryPool.get());
705 bounded_vector<int> componentSize(lineRefs.size(), 1, settings.memoryPool.get());
706 std::iota(parent.begin(), parent.end(), 0);
707 float minZ = std::numeric_limits<float>::max();
708 float maxZ = std::numeric_limits<float>::lowest();
709 for (const auto& lineRef : lineRefs) {
710 minZ = std::min(minZ, lineRef.zBeam);
711 maxZ = std::max(maxZ, lineRef.zBeam);
712 }
713 const auto nZBins = std::max(1, 1 + static_cast<int>((maxZ - minZ) / coarseZWindow));
714 auto getZBin = [&](const float z) {
715 return std::clamp(static_cast<int>((z - minZ) / coarseZWindow), 0, nZBins - 1);
716 };
717
718 auto findRoot = [&](int idx) {
719 int root = idx;
720 while (parent[root] != root) {
721 root = parent[root];
722 }
723 while (parent[idx] != idx) {
724 const auto next = parent[idx];
725 parent[idx] = root;
726 idx = next;
727 }
728 return root;
729 };
730
731 auto unite = [&](const int lhs, const int rhs) {
732 auto lhsRoot = findRoot(lhs);
733 auto rhsRoot = findRoot(rhs);
734 if (lhsRoot == rhsRoot) {
735 return;
736 }
737 if (componentSize[lhsRoot] < componentSize[rhsRoot]) {
738 std::swap(lhsRoot, rhsRoot);
739 }
740 parent[rhsRoot] = lhsRoot;
741 componentSize[lhsRoot] += componentSize[rhsRoot];
742 };
743
744 using ActiveEntry = std::pair<TimeStampType, int>;
745 bounded_vector<ActiveEntry> activeEntries(settings.memoryPool.get());
746 std::priority_queue<ActiveEntry, bounded_vector<ActiveEntry>, std::greater<ActiveEntry>> activeByUpper(std::greater<ActiveEntry>{}, std::move(activeEntries));
747 bounded_vector<uint8_t> activeMask(lineRefs.size(), 0, settings.memoryPool.get());
748 bounded_vector<bounded_vector<int>> activeByZBin(settings.memoryPool.get());
749 activeByZBin.reserve(nZBins);
750 for (int iBin = 0; iBin < nZBins; ++iBin) {
751 activeByZBin.emplace_back();
752 }
753 for (const auto lineRefIdx : sortedByLower) {
754 const auto& lineRef = lineRefs[lineRefIdx];
755 const auto& line = lines[lineRef.lineIndex];
756 const auto currentLower = line.mTime.lower();
757
758 while (!activeByUpper.empty() && activeByUpper.top().first < currentLower) {
759 activeMask[activeByUpper.top().second] = 0;
760 activeByUpper.pop();
761 }
762
763 const auto zBin = getZBin(lineRef.zBeam);
764 for (int neighborBin = std::max(0, zBin - 1); neighborBin <= std::min(nZBins - 1, zBin + 1); ++neighborBin) {
765 auto& bucket = activeByZBin[neighborBin];
766 size_t writePos = 0;
767 for (size_t readPos = 0; readPos < bucket.size(); ++readPos) {
768 const auto oLineRefIdx = bucket[readPos];
769 if (!activeMask[oLineRefIdx]) {
770 continue;
771 }
772 bucket[writePos++] = oLineRefIdx;
773 const auto& oLineRef = lineRefs[oLineRefIdx];
774 if (o2::gpu::GPUCommonMath::Abs(lineRef.zBeam - oLineRef.zBeam) >= coarseZWindow) {
775 continue;
776 }
777 const auto& otherLine = lines[oLineRef.lineIndex];
778 if (line.mTime.isCompatible(otherLine.mTime)) {
779 unite(lineRefIdx, oLineRefIdx);
780 }
781 }
782 bucket.resize(writePos);
783 }
784
785 activeMask[lineRefIdx] = 1;
786 activeByUpper.emplace(line.mTime.upper(), lineRefIdx);
787 activeByZBin[zBin].push_back(lineRefIdx);
788 }
789
790 std::unordered_map<int, bounded_vector<int>> components;
791 components.reserve(lineRefs.size());
792 for (int lineRefIdx = 0; lineRefIdx < static_cast<int>(lineRefs.size()); ++lineRefIdx) {
793 const auto root = findRoot(lineRefIdx);
794 auto [it, inserted] = components.try_emplace(root, std::pmr::polymorphic_allocator<int>{settings.memoryPool.get()});
795 (void)inserted;
796 it->second.push_back(lineRefIdx);
797 }
798
799 bounded_vector<OrderedComponent> orderedComponents(settings.memoryPool.get());
800 orderedComponents.reserve(components.size());
801 for (auto& [root, members] : components) {
802 (void)root;
803 if (members.size() < 2) {
804 continue;
805 }
806 std::sort(members.begin(), members.end(), [&](const int lhs, const int rhs) {
807 const auto lhsLower = lines[lineRefs[lhs].lineIndex].mTime.lower();
808 const auto rhsLower = lines[lineRefs[rhs].lineIndex].mTime.lower();
809 if (lhsLower != rhsLower) {
810 return lhsLower < rhsLower;
811 }
812 return lineRefs[lhs].lineIndex < lineRefs[rhs].lineIndex;
813 });
814 orderedComponents.emplace_back(settings.memoryPool);
815 orderedComponents.back().center = lineRefs[members.front()].tCenter;
816 orderedComponents.back().members = std::move(members);
817 }
818
819 std::sort(orderedComponents.begin(), orderedComponents.end(), [](const auto& lhs, const auto& rhs) {
820 if (o2::gpu::GPUCommonMath::Abs(lhs.center - rhs.center) > TieTolerance) {
821 return lhs.center < rhs.center;
822 }
823 return lhs.members.front() < rhs.members.front();
824 });
825 clusters.reserve(orderedComponents.size());
826 for (auto& component : orderedComponents) {
827 clusters.push_back(std::move(component.members));
828 }
829 return clusters;
830}
831
832bounded_vector<VertexSeed> buildSeeds(std::span<const int> members,
833 std::span<const LineRef> lineRefs,
834 std::span<const Line> lines,
835 const Settings& settings)
836{
837 SeedHistogram histogram(members, lineRefs, lines, settings);
838 bounded_vector<VertexSeed> seeds(settings.memoryPool.get());
839 seeds.reserve(MaxSeedsPerCluster);
840 float leadingPeakSupport = 0.f;
841
842 while (static_cast<int>(seeds.size()) < MaxSeedsPerCluster) {
843 const auto peak = histogram.findPeakBin();
844 if (peak < 0) {
845 break;
846 }
847 const auto peakSupport = histogram.getPeakSupport(peak);
848 if (peakSupport < 2.f) {
849 break;
850 }
851 if (leadingPeakSupport <= 0.f) {
852 leadingPeakSupport = peakSupport;
853 } else if (peakSupport < std::max(2.f, MinRelativePeakSupport * leadingPeakSupport)) {
854 break;
855 }
856 auto localMembers = histogram.collectLocalMembers(peak, 0, 0);
857 if (localMembers.size() < 2) {
858 localMembers = histogram.collectLocalMembers(peak, settings.seedMemberRadiusTime, settings.seedMemberRadiusZ);
859 }
860 if (localMembers.size() < 2) {
861 histogram.suppressPeak(peak);
862 continue;
863 }
864
865 VertexSeed seed(settings.memoryPool);
866 seed.vertex = {settings.beamX, settings.beamY, histogram.getPeakZCenter(peak)};
867 seed.time = histogram.getPeakTimeInterval(peak);
868 seed.scale2 = InitialScale2;
869
870 auto fitted = fitSeed(seed, localMembers, lineRefs, lines, settings.memoryPool, settings.pairCut2);
871 if (fitted.valid && fitted.contributors.size() >= 2) {
872 seeds.push_back(std::move(fitted));
873 histogram.suppressPeakNeighborhood(peak);
874 } else {
875 histogram.suppressPeak(peak);
876 }
877 }
878
879 return seeds;
880}
881
882void assignLinesToSeeds(bounded_vector<VertexSeed>& seeds,
883 std::span<const int> members,
884 std::span<const LineRef> lineRefs,
885 std::span<const Line> lines,
886 const float pairCut2)
887{
888 for (auto& seed : seeds) {
889 seed.assigned.clear();
890 }
891
892 for (const auto lineRefIdx : members) {
893 const auto lineIdx = lineRefs[lineRefIdx].lineIndex;
894 const auto& line = lines[lineIdx];
895
896 int bestSeed = -1;
897 float bestScore = std::numeric_limits<float>::max();
898 size_t bestMult = 0;
899 float bestZResidual = std::numeric_limits<float>::max();
900
901 for (int seedIdx = 0; seedIdx < static_cast<int>(seeds.size()); ++seedIdx) {
902 const auto& seed = seeds[seedIdx];
903 if (!seed.valid || seed.contributors.size() < 2) {
904 continue;
905 }
906 if (!line.mTime.isCompatible(seed.time)) {
907 continue;
908 }
909
910 const auto distance2 = Line::getDistance2FromPoint(line, seed.vertex);
911 if (distance2 >= pairCut2) {
912 continue;
913 }
914
915 const auto score = distance2 / std::max(seed.scale2, MinScale2);
916 const auto zResidual = o2::gpu::GPUCommonMath::Abs(lineRefs[lineRefIdx].zBeam - seed.vertex[2]);
917 const auto multiplicity = seed.contributors.size();
918
919 const auto betterScore = score + TieTolerance < bestScore;
920 const auto betterMultiplicity = o2::gpu::GPUCommonMath::Abs(score - bestScore) <= TieTolerance && multiplicity > bestMult;
921 const auto betterZ = o2::gpu::GPUCommonMath::Abs(score - bestScore) <= TieTolerance &&
922 multiplicity == bestMult && zResidual + constants::Tolerance < bestZResidual;
923 if (betterScore || betterMultiplicity || betterZ) {
924 bestSeed = seedIdx;
925 bestScore = score;
926 bestMult = multiplicity;
927 bestZResidual = zResidual;
928 }
929 }
930
931 if (bestSeed >= 0) {
932 seeds[bestSeed].assigned.push_back(lineRefIdx);
933 }
934 }
935}
936
937ClusterLines materializeCluster(const VertexSeed& seed,
938 std::span<const LineRef> lineRefs,
939 std::span<const Line> lines,
940 const std::shared_ptr<BoundedMemoryResource>& mr)
941{
942 bounded_vector<int> lineIndices{mr.get()};
943 lineIndices.reserve(seed.contributors.size());
944 for (const auto lineRefIdx : seed.contributors) {
945 lineIndices.push_back(lineRefs[lineRefIdx].lineIndex);
946 }
947 std::sort(lineIndices.begin(), lineIndices.end());
948 lineIndices.erase(std::unique(lineIndices.begin(), lineIndices.end()), lineIndices.end());
949
950 if (lineIndices.size() < 2) {
951 return {};
952 }
953
954 return {std::span<const int>{lineIndices.data(), lineIndices.size()}, lines};
955}
956
957} // namespace
958
959bounded_vector<ClusterLines> buildClusters(std::span<const Line> lines, const Settings& settings)
960{
962 if (lines.size() < 2) {
963 return clusters;
964 }
965
966 bounded_vector<LineRef> refs(settings.memoryPool.get());
967 refs.reserve(lines.size());
968 for (int lineIdx = 0; lineIdx < static_cast<int>(lines.size()); ++lineIdx) {
969 LineRef ref(lines[lineIdx], lineIdx, settings.beamX, settings.beamY, settings.maxZ);
970 if (!ref.isDead()) {
971 refs.push_back(ref);
972 }
973 }
974
975 if (refs.size() < 2) {
976 return clusters;
977 }
978
979 const auto coarseClusters = buildCoarseClusters(refs, lines, settings);
980
981 for (const auto& members : coarseClusters) {
982 auto seeds = buildSeeds(members, refs, lines, settings);
983 if (seeds.empty()) {
984 continue;
985 }
986
987 for (auto& seed : seeds) {
988 if (!seed.isUsableSeed()) {
989 seed.valid = false;
990 continue;
991 }
992 auto contributors = collectCompatibleContributors(seed, members, refs, lines, settings.memoryPool, settings.pairCut2);
993 if (contributors.size() < 2) {
994 seed.valid = false;
995 continue;
996 }
997 seed.contributors = std::move(contributors);
998 }
999 compactSeeds(seeds);
1000 if (seeds.empty()) {
1001 continue;
1002 }
1003 deduplicateSeeds(seeds, settings);
1004 if (seeds.empty()) {
1005 continue;
1006 }
1007 assignLinesToSeeds(seeds, members, refs, lines, settings.pairCut2);
1008 for (auto& seed : seeds) {
1009 if (seed.assigned.size() < 2) {
1010 seed.valid = false;
1011 continue;
1012 }
1013 seed = fitSeed(seed, seed.assigned, refs, lines, settings.memoryPool, settings.pairCut2);
1014 if (!seed.isUsableSeed()) {
1015 seed.valid = false;
1016 continue;
1017 }
1018 }
1019 compactSeeds(seeds);
1020 deduplicateRefittedSeeds(seeds, settings);
1021 for (auto& refit : seeds) {
1022 auto cluster = materializeCluster(refit, refs, lines, settings.memoryPool);
1023 if (cluster.getSize() < 2) {
1024 continue;
1025 }
1026 if (!cluster.isValid()) {
1027 continue;
1028 }
1029 clusters.push_back(std::move(cluster));
1030 }
1031 }
1032
1033 return clusters;
1034}
1035
1036} // namespace o2::its::line_vertexer
header::DataOrigin origin
size_t minSize
uint64_t vertex
Definition RawEventData.h:9
int16_t time
Definition RawEventData.h:4
int32_t i
bounded_vector< int > members
bool valid
float tCenter
float tHalfWidth
bounded_vector< int > assigned
float timeMin
int nTimeBins
float scale2
int lineIndex
float zBeam
float center
bounded_vector< int > contributors
float zMin
bounded_vector< float > bins
float timeBinSize
float zBinSize
uint32_t j
Definition RawData.h:0
float sum(float s, o2::dcs::DataPointValue v)
Definition dcs-ccdb.cxx:39
GLuint index
Definition glcorearb.h:781
GLint GLsizei width
Definition glcorearb.h:270
GLuint GLuint GLfloat weight
Definition glcorearb.h:5477
GLint GLenum GLint components
Definition glcorearb.h:5520
GLsizei const GLfloat * value
Definition glcorearb.h:819
typedef void(APIENTRYP PFNGLCULLFACEPROC)(GLenum mode)
GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat maxZ
Definition glcorearb.h:2910
GLint ref
Definition glcorearb.h:291
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
GLfloat GLfloat minZ
Definition glcorearb.h:2910
GLuint GLfloat GLfloat GLfloat GLfloat GLfloat GLfloat s0
Definition glcorearb.h:5034
int32_t const char int32_t line
auto get(const std::byte *buffer, size_t=0)
Definition DataHeader.h:454
bounded_vector< ClusterLines > buildClusters(std::span< const Line > lines, const Settings &settings)
std::pmr::vector< T > bounded_vector
uint32_t TimeStampType
Definition TimeEstBC.h:26
D const SVectorGPU< T, D > & rhs
Definition SMatrixGPU.h:193
if(!okForPhiMin(phi0, phi1))
R median(std::vector< T > v)
Definition fit.h:520
constexpr auto getIndex(const container_T &container, typename container_T::const_iterator iter) -> typename container_T::source_type
Definition algorithm.h:65
static float getDistance2FromPoint(const Line &line, const std::array< float, 3 > &point)
std::shared_ptr< BoundedMemoryResource > memoryPool
std::vector< Cluster > clusters