Project
Loading...
Searching...
No Matches
GPUTPCGMMerger.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
14
15#define GPUCA_CADEBUG 0
16#define GPUCA_MERGE_LOOPER_MC 0
17
18#include "GPUCommonDef.h"
19
20#if !defined(GPUCA_GPUCODE) && (defined(GPUCA_MERGER_BY_MC_LABEL) || defined(GPUCA_CADEBUG_ENABLED) || GPUCA_MERGE_LOOPER_MC)
22#include "GPUROOTDump.h"
23#endif
24
25#ifndef GPUCA_GPUCODE_DEVICE
26#include <cstdio>
27#include <cstring>
28#include <cmath>
29#include "GPUReconstruction.h"
30#endif
31
32#include "GPUTPCTracker.h"
33#include "GPUTPCClusterData.h"
34#include "GPUTPCTrackParam.h"
35#include "GPUTPCGMMerger.h"
36#include "GPUO2DataTypes.h"
37#include "TPCFastTransform.h"
38#include "GPUTPCConvertImpl.h"
39#include "GPUTPCGeometry.h"
40#include "GPUDefParametersRuntime.h"
41
42#include "GPUCommonMath.h"
43#include "GPUCommonAlgorithm.h"
44#include "GPUCommonConstants.h"
45
46#include "GPUTPCTrackParam.h"
47#include "GPUTPCGMMergedTrack.h"
48#include "GPUParam.h"
50
51#include "GPUTPCGMTrackParam.h"
52#include "GPUTPCGMSectorTrack.h"
53#include "GPUTPCGMBorderTrack.h"
54
57#ifndef GPUCA_GPUCODE
60#endif
61
62static constexpr int32_t kMaxParts = 400;
63static constexpr int32_t kMaxClusters = GPUCA_MERGER_MAX_TRACK_CLUSTERS;
64
65using namespace o2::gpu;
66using namespace o2::tpc;
67using namespace gputpcgmmergertypes;
68
69namespace o2::gpu::internal
70{
72 float refz;
73 float x;
74 float y;
75 uint32_t id;
76};
77
79 GPUd() bool operator()(const GPUTPCGMBorderRange& a, const GPUTPCGMBorderRange& b)
80 {
81 return GPUCA_DETERMINISTIC_CODE((a.fMax != b.fMax) ? (a.fMax < b.fMax) : (a.fId < b.fId), a.fMax < b.fMax);
82 }
83};
85 GPUd() bool operator()(const GPUTPCGMBorderRange& a, const GPUTPCGMBorderRange& b)
86 {
87 return GPUCA_DETERMINISTIC_CODE((a.fMin != b.fMin) ? (a.fMin < b.fMin) : (a.fId < b.fId), a.fMin < b.fMin);
88 }
89};
90
94 GPUd() bool operator()(const int32_t aa, const int32_t bb)
95 {
98 if (a.CCE() != b.CCE()) {
99 return a.CCE() > b.CCE();
100 }
101 if (a.Legs() != b.Legs()) {
102 return a.Legs() > b.Legs();
103 }
104 GPUCA_DETERMINISTIC_CODE( // clang-format off
105 if (a.NClusters() != b.NClusters()) {
106 return a.NClusters() > b.NClusters();
107 } if (CAMath::Abs(a.GetParam().GetQPt()) != CAMath::Abs(b.GetParam().GetQPt())) {
108 return CAMath::Abs(a.GetParam().GetQPt()) > CAMath::Abs(b.GetParam().GetQPt());
109 } if (a.GetParam().GetY() != b.GetParam().GetY()) {
110 return a.GetParam().GetY() > b.GetParam().GetY();
111 }
112 return aa > bb;
113 , // !GPUCA_DETERMINISTIC_CODE
114 return a.NClusters() > b.NClusters();
115 ) // clang-format on
116 }
117};
118
122 GPUd() bool operator()(const int32_t aa, const int32_t bb)
123 {
126 GPUCA_DETERMINISTIC_CODE( // clang-format off
127 if (CAMath::Abs(a.GetParam().GetQPt()) != CAMath::Abs(b.GetParam().GetQPt())) {
128 return CAMath::Abs(a.GetParam().GetQPt()) > CAMath::Abs(b.GetParam().GetQPt());
129 } if (a.GetParam().GetY() != b.GetParam().GetY()) {
130 return a.GetParam().GetY() > b.GetParam().GetY();
131 }
132 return a.GetParam().GetZ() > b.GetParam().GetZ();
133 , // !GPUCA_DETERMINISTIC_CODE
134 return CAMath::Abs(a.GetParam().GetQPt()) > CAMath::Abs(b.GetParam().GetQPt());
135 ) // clang-format on
136 }
137};
138
140 GPUd() bool operator()(const MergeLooperParam& a, const MergeLooperParam& b)
141 {
142 return GPUCA_DETERMINISTIC_CODE(CAMath::Abs(a.refz) != CAMath::Abs(b.refz) ? CAMath::Abs(a.refz) < CAMath::Abs(b.refz) : a.id < b.id, CAMath::Abs(a.refz) < CAMath::Abs(b.refz));
143 }
144};
145
146} // namespace o2::gpu::internal
147
148using namespace o2::gpu::internal;
149
150#ifndef GPUCA_GPUCODE
151
152#include "GPUQA.h"
153#include "GPUMemorySizeScalers.h"
154
156{
157 for (int32_t iSector = 0; iSector < NSECTORS; iSector++) {
158 mNextSectorInd[iSector] = iSector + 1;
159 mPrevSectorInd[iSector] = iSector - 1;
160 }
161 int32_t mid = NSECTORS / 2 - 1;
162 int32_t last = NSECTORS - 1;
163 mNextSectorInd[mid] = 0;
164 mPrevSectorInd[0] = mid;
165 mNextSectorInd[last] = NSECTORS / 2;
166 mPrevSectorInd[NSECTORS / 2] = last;
167}
168
169// DEBUG CODE
170#if !defined(GPUCA_GPUCODE) && (defined(GPUCA_MERGER_BY_MC_LABEL) || defined(GPUCA_CADEBUG_ENABLED) || GPUCA_MERGE_LOOPER_MC)
171#include "GPUQAHelper.h"
172
173void GPUTPCGMMerger::CheckMergedTracks()
174{
175 std::vector<bool> trkUsed(SectorTrackInfoLocalTotal());
176 for (int32_t i = 0; i < SectorTrackInfoLocalTotal(); i++) {
177 trkUsed[i] = false;
178 }
179
180 for (int32_t itr = 0; itr < SectorTrackInfoLocalTotal(); itr++) {
181 GPUTPCGMSectorTrack& track = mSectorTrackInfos[itr];
182 if (track.PrevSegmentNeighbour() >= 0) {
183 continue;
184 }
185 if (track.PrevNeighbour() >= 0) {
186 continue;
187 }
188 int32_t leg = 0;
189 GPUTPCGMSectorTrack *trbase = &track, *tr = &track;
190 while (true) {
191 int32_t iTrk = tr - mSectorTrackInfos;
192 if (trkUsed[iTrk]) {
193 GPUError("FAILURE: double use");
194 }
195 trkUsed[iTrk] = true;
196
197 int32_t jtr = tr->NextSegmentNeighbour();
198 if (jtr >= 0) {
199 tr = &(mSectorTrackInfos[jtr]);
200 continue;
201 }
202 jtr = trbase->NextNeighbour();
203 if (jtr >= 0) {
204 trbase = &(mSectorTrackInfos[jtr]);
205 tr = trbase;
206 if (tr->PrevSegmentNeighbour() >= 0) {
207 break;
208 }
209 leg++;
210 continue;
211 }
212 break;
213 }
214 }
215 for (int32_t i = 0; i < SectorTrackInfoLocalTotal(); i++) {
216 if (trkUsed[i] == false) {
217 GPUError("FAILURE: trk missed");
218 }
219 }
220}
221
222template <class T>
223inline const auto* resolveMCLabels(const o2::dataformats::ConstMCTruthContainerView<o2::MCCompLabel>* a, const AliHLTTPCClusterMCLabel* b)
224{
225 return a;
226}
227template <>
228inline const auto* resolveMCLabels<AliHLTTPCClusterMCLabel>(const o2::dataformats::ConstMCTruthContainerView<o2::MCCompLabel>* a, const AliHLTTPCClusterMCLabel* b)
229{
230 return b;
231}
232
233template <class T, class S>
234int64_t GPUTPCGMMerger::GetTrackLabelA(const S& trk) const
235{
236 GPUTPCGMSectorTrack* sectorTrack = nullptr;
237 int32_t nClusters = 0;
238 if constexpr (std::is_same_v<S, GPUTPCGMBorderTrack&>) {
239 sectorTrack = &mSectorTrackInfos[trk.TrackID()];
240 nClusters = sectorTrack->OrigTrack()->NHits();
241 } else {
242 nClusters = trk.NClusters();
243 }
244 auto acc = GPUTPCTrkLbl<false, GPUTPCTrkLbl_ret>(resolveMCLabels<T>(GetConstantMem()->ioPtrs.clustersNative ? GetConstantMem()->ioPtrs.clustersNative->clustersMCTruth : nullptr, GetConstantMem()->ioPtrs.mcLabelsTPC), 0.5f);
245 for (int32_t i = 0; i < nClusters; i++) {
246 int32_t id;
247 if constexpr (std::is_same_v<S, GPUTPCGMBorderTrack&>) {
248 const GPUTPCTracker& tracker = GetConstantMem()->tpcTrackers[sectorTrack->Sector()];
249 const GPUTPCHitId& ic = tracker.TrackHits()[sectorTrack->OrigTrack()->FirstHitID() + i];
250 id = tracker.Data().ClusterDataIndex(tracker.Data().Row(ic.RowIndex()), ic.HitIndex()) + GetConstantMem()->ioPtrs.clustersNative->clusterOffset[sectorTrack->Sector()][0];
251 } else {
252 id = mClusters[trk.FirstClusterRef() + i].num;
253 }
254 acc.addLabel(id);
255 }
256 return acc.computeLabel().id;
257}
258
259template <class S>
260int64_t GPUTPCGMMerger::GetTrackLabel(const S& trk) const
261{
262#ifdef GPUCA_TPC_GEOMETRY_O2
263 if (GetConstantMem()->ioPtrs.clustersNative->clustersMCTruth) {
264 return GetTrackLabelA<o2::dataformats::ConstMCTruthContainerView<o2::MCCompLabel>, S>(trk);
265 } else
266#endif
267 {
268 return GetTrackLabelA<AliHLTTPCClusterMCLabel, S>(trk);
269 }
270}
271
272#endif
273// END DEBUG CODE
274
275void GPUTPCGMMerger::PrintMergeGraph(const GPUTPCGMSectorTrack* trk, std::ostream& out) const
276{
277 const GPUTPCGMSectorTrack* orgTrack = trk;
278 while (trk->PrevSegmentNeighbour() >= 0) {
279 trk = &mSectorTrackInfos[trk->PrevSegmentNeighbour()];
280 }
281 const GPUTPCGMSectorTrack* orgTower = trk;
282 while (trk->PrevNeighbour() >= 0) {
283 trk = &mSectorTrackInfos[trk->PrevNeighbour()];
284 }
285
286 int32_t nextId = trk - mSectorTrackInfos;
287 out << "Graph of track " << (orgTrack - mSectorTrackInfos) << "\n";
288 while (nextId >= 0) {
289 trk = &mSectorTrackInfos[nextId];
290 if (trk->PrevSegmentNeighbour() >= 0) {
291 out << "TRACK TREE INVALID!!! " << trk->PrevSegmentNeighbour() << " --> " << nextId << "\n";
292 }
293 out << (trk == orgTower ? "--" : " ");
294 while (nextId >= 0) {
295 GPUTPCGMSectorTrack* trk2 = &mSectorTrackInfos[nextId];
296 if (trk != trk2 && (trk2->PrevNeighbour() >= 0 || trk2->NextNeighbour() >= 0)) {
297 out << " (TRACK TREE INVALID!!! " << trk2->PrevNeighbour() << " <-- " << nextId << " --> " << trk2->NextNeighbour() << ") ";
298 }
299 char tmp[128];
300 snprintf(tmp, 128, " %s%5d(%5.2f)", trk2 == orgTrack ? "!" : " ", nextId, trk2->QPt());
301 out << tmp;
302 nextId = trk2->NextSegmentNeighbour();
303 }
304 out << "\n";
305 nextId = trk->NextNeighbour();
306 }
307}
308
310
312{
313 computePointerWithAlignment(mem, mSectorTrackInfos, mNTotalSectorTracks);
314 computePointerWithAlignment(mem, mSectorTrackInfoIndex, NSECTORS * 2 + 1);
315 if (mRec->GetProcessingSettings().deterministicGPUReconstruction) {
316 computePointerWithAlignment(mem, mTmpSortMemory, std::max(mNTotalSectorTracks, mNMaxTracks * 2));
317 }
318
319 void* memBase = mem;
320 computePointerWithAlignment(mem, mBorderMemory, 2 * mNTotalSectorTracks); // MergeBorders & Resolve
321 computePointerWithAlignment(mem, mBorderRangeMemory, 2 * mNTotalSectorTracks);
322 int32_t nTracks = 0;
323 for (int32_t iSector = 0; iSector < NSECTORS; iSector++) {
324 const int32_t n = *mRec->GetConstantMem().tpcTrackers[iSector].NTracks();
325 mBorder[iSector] = mBorderMemory + 2 * nTracks;
326 mBorder[NSECTORS + iSector] = mBorderMemory + 2 * nTracks + n;
327 mBorderRange[iSector] = mBorderRangeMemory + 2 * nTracks;
328 nTracks += n;
329 }
330 computePointerWithAlignment(mem, mTrackLinks, mNTotalSectorTracks);
331 computePointerWithAlignment(mem, mTrackCCRoots, mNTotalSectorTracks);
332 void* memMax = mem;
333 mem = memBase;
334 computePointerWithAlignment(mem, mTrackIDs, GPUCA_NSECTORS * mNMaxSingleSectorTracks); // UnpackResetIds - RefitSectorTracks - UnpackSectorGlobal
335 memMax = (void*)std::max((size_t)mem, (size_t)memMax);
336 mem = memBase;
337 computePointerWithAlignment(mem, mTrackSort, mNMaxTracks); // PrepareClustersForFit0 - SortTracksQPt - PrepareClustersForFit1 - PrepareClustersForFit1 / Finalize0 - Finalize2
338 computePointerWithAlignment(mem, mSharedCount, mNMaxClusters);
339 memMax = (void*)std::max((size_t)mem, (size_t)memMax);
340 mem = memBase;
341 computePointerWithAlignment(mem, mLoopData, mNMaxTracks); // GPUTPCGMMergerTrackFit - GPUTPCGMMergerFollowLoopers
342 computePointerWithAlignment(mem, mRetryRefitIds, mNMaxTracks); // Reducing mNMaxTracks for mLoopData / mRetryRefitIds does not save memory, since the other parts are larger anyway
343 memMax = (void*)std::max((size_t)mem, (size_t)memMax);
344 mem = memBase;
345 computePointerWithAlignment(mem, mLooperCandidates, mNMaxLooperMatches); // MergeLoopers 1-3
346 memMax = (void*)std::max((size_t)mem, (size_t)memMax);
347 return memMax;
348}
349
351{
352 computePointerWithAlignment(mem, mMemory);
353 return mem;
354}
355
357{
358 computePointerWithAlignment(mem, mTrackOrderAttach, mNMaxTracks);
359 const bool mergerSortTracks = mRec->GetProcessingSettings().mergerSortTracks == -1 ? mRec->getGPUParameters(mRec->GetRecoStepsGPU() & GPUDataTypes::RecoStep::TPCMerging).par_SORT_BEFORE_FIT : mRec->GetProcessingSettings().mergerSortTracks;
360 if (mergerSortTracks) {
361 computePointerWithAlignment(mem, mTrackOrderProcess, mNMaxTracks);
362 }
363 return mem;
364}
365
367{
368 computePointerWithAlignment(mem, mMergedTracks, mNMaxTracks);
369 if (mRec->GetParam().dodEdxEnabled) {
370 computePointerWithAlignment(mem, mMergedTracksdEdx, mNMaxTracks);
371 if (mRec->GetParam().rec.tpc.dEdxClusterRejectionFlagMask != mRec->GetParam().rec.tpc.dEdxClusterRejectionFlagMaskAlt) {
372 computePointerWithAlignment(mem, mMergedTracksdEdxAlt, mNMaxTracks);
373 }
374 }
375 computePointerWithAlignment(mem, mClusters, mNMaxOutputTrackClusters);
376 if (mRec->GetParam().par.earlyTpcTransform) {
377 computePointerWithAlignment(mem, mClustersXYZ, mNMaxOutputTrackClusters);
378 }
379 computePointerWithAlignment(mem, mClusterAttachment, mNMaxClusters);
380 return mem;
381}
382
384{
385 if ((mRec->GetRecoSteps() & GPUDataTypes::RecoStep::Refit) || mRec->GetProcessingSettings().outputSharedClusterMap) {
386 computePointerWithAlignment(mem, mClusterStateExt, mNMaxClusters);
387 } else {
388 mClusterStateExt = nullptr;
389 }
390 return mem;
391}
392
394{
395 computePointerWithAlignment(mem, mOutputTracksTPCO2, HostProcessor(this).NOutputTracksTPCO2());
396 return mem;
397}
398
400{
401 computePointerWithAlignment(mem, mOutputClusRefsTPCO2, HostProcessor(this).NOutputClusRefsTPCO2());
402 return mem;
403}
404
406{
407 computePointerWithAlignment(mem, mOutputTracksTPCO2MC, HostProcessor(this).NOutputTracksTPCO2());
408 return mem;
409}
410
412{
413 computePointerWithAlignment(mem, mTrackSortO2, mNMaxTracks);
414 computePointerWithAlignment(mem, mClusRefTmp, mNMaxTracks);
415 return mem;
416}
417
419{
425 if (mRec->GetProcessingSettings().createO2Output) {
429 if (mRec->GetProcessingSettings().runMC) {
431 }
432 }
434}
435
437{
438 mNTotalSectorTracks = 0;
439 mNClusters = 0;
440 mNMaxSingleSectorTracks = 0;
441 for (int32_t iSector = 0; iSector < NSECTORS; iSector++) {
442 uint32_t ntrk = *mRec->GetConstantMem().tpcTrackers[iSector].NTracks();
443 mNTotalSectorTracks += ntrk;
444 mNClusters += *mRec->GetConstantMem().tpcTrackers[iSector].NTrackHits();
445 if (mNMaxSingleSectorTracks < ntrk) {
446 mNMaxSingleSectorTracks = ntrk;
447 }
448 }
449 mNMaxOutputTrackClusters = mRec->MemoryScalers()->NTPCMergedTrackHits(mNClusters);
450 if (CAMath::Abs(Param().polynomialField.GetNominalBz()) < (gpu_common_constants::kZeroFieldCut * gpu_common_constants::kCLight)) {
451 mNMaxTracks = mRec->MemoryScalers()->getValue(mNTotalSectorTracks, mNTotalSectorTracks); // 0 magnetic field
452 } else {
453 mNMaxTracks = mRec->MemoryScalers()->NTPCMergedTracks(mNTotalSectorTracks);
454 }
455 if (io.clustersNative) {
456 mNMaxClusters = io.clustersNative->nClustersTotal;
458 mNMaxClusters = 0;
459 for (int32_t i = 0; i < NSECTORS; i++) {
460 mNMaxClusters += mRec->GetConstantMem().tpcTrackers[i].NHitsTotal();
461 }
462 } else {
463 mNMaxClusters = mNClusters;
464 }
465 mNMaxLooperMatches = mNMaxClusters / 4; // We have that much scratch memory anyway
466}
467
469{
470 for (int32_t i = 0; i < NSECTORS; i++) {
471 if (mRec->GetConstantMem().tpcTrackers[i].CommonMemory()->nLocalTracks > (int32_t)mNMaxSingleSectorTracks) {
472 throw std::runtime_error("mNMaxSingleSectorTracks too small");
473 }
474 }
476 throw std::runtime_error("Must run also sector tracking");
477 }
478 return 0;
479}
480
481#endif // GPUCA_GPUCODE
482
483GPUd() void GPUTPCGMMerger::ClearTrackLinks(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, bool output)
484{
485 const int32_t n = output ? mMemory->nMergedTracks : SectorTrackInfoLocalTotal();
486 for (int32_t i = iBlock * nThreads + iThread; i < n; i += nThreads * nBlocks) {
487 mTrackLinks[i] = -1;
488 }
489}
490
491GPUd() int32_t GPUTPCGMMerger::RefitSectorTrack(GPUTPCGMSectorTrack& sectorTrack, const GPUTPCTrack* inTrack, float alpha, int32_t sector)
492{
494 prop.SetMaterialTPC();
495 prop.SetMaxSinPhi(GPUCA_MAX_SIN_PHI);
496 prop.SetToyMCEventsFlag(false);
497 prop.SetSeedingErrors(true); // Larger errors for seeds, better since we don't start with good hypothesis
498 prop.SetFitInProjections(false);
499 prop.SetPolynomialField(&Param().polynomialField);
501 trk.X() = inTrack->Param().GetX();
502 trk.Y() = inTrack->Param().GetY();
503 trk.Z() = inTrack->Param().GetZ();
504 trk.SinPhi() = inTrack->Param().GetSinPhi();
505 trk.DzDs() = inTrack->Param().GetDzDs();
506 trk.QPt() = inTrack->Param().GetQPt();
507 trk.TZOffset() = Param().par.earlyTpcTransform ? inTrack->Param().GetZOffset() : GetConstantMem()->calibObjects.fastTransformHelper->getCorrMap()->convZOffsetToVertexTime(sector, inTrack->Param().GetZOffset(), Param().continuousMaxTimeBin);
508 trk.ShiftZ(this, sector, sectorTrack.ClusterZT0(), sectorTrack.ClusterZTN(), inTrack->Param().GetX(), inTrack->Param().GetX()); // We do not store the inner / outer cluster X, so we just use the track X instead
509 sectorTrack.SetX2(0.f);
510 for (int32_t way = 0; way < 2; way++) {
511 if (way) {
512 prop.SetFitInProjections(true);
513 prop.SetPropagateBzOnly(true);
514 }
515 trk.ResetCovariance();
516 prop.SetTrack(&trk, alpha);
517 int32_t start = way ? inTrack->NHits() - 1 : 0;
518 int32_t end = way ? 0 : (inTrack->NHits() - 1);
519 int32_t incr = way ? -1 : 1;
520 for (int32_t i = start; i != end; i += incr) {
521 float x, y, z;
522 int32_t row, flags;
523 const GPUTPCTracker& tracker = GetConstantMem()->tpcTrackers[sector];
524 const GPUTPCHitId& ic = tracker.TrackHits()[inTrack->FirstHitID() + i];
525 int32_t clusterIndex = tracker.Data().ClusterDataIndex(tracker.Data().Row(ic.RowIndex()), ic.HitIndex());
526 row = ic.RowIndex();
527 const ClusterNative& cl = GetConstantMem()->ioPtrs.clustersNative->clustersLinear[GetConstantMem()->ioPtrs.clustersNative->clusterOffset[sector][0] + clusterIndex];
528 flags = cl.getFlags();
529 if (Param().par.earlyTpcTransform) {
530 x = tracker.Data().ClusterData()[clusterIndex].x;
531 y = tracker.Data().ClusterData()[clusterIndex].y;
532 z = tracker.Data().ClusterData()[clusterIndex].z - trk.TZOffset();
533 } else {
534 GetConstantMem()->calibObjects.fastTransformHelper->Transform(sector, row, cl.getPad(), cl.getTime(), x, y, z, trk.TZOffset());
535 }
536 if (prop.PropagateToXAlpha(x, alpha, true)) {
537 return way == 0;
538 }
539 trk.ConstrainSinPhi();
540 if (prop.Update(y, z, row, Param(), flags & GPUTPCGMMergedTrackHit::clustererAndSharedFlags, 0, nullptr, false, sector, -1.f, 0.f, 0.f)) { // TODO: Use correct time / avgCharge
541 return way == 0;
542 }
543 trk.ConstrainSinPhi();
544 }
545 if (way) {
546 sectorTrack.SetParam2(trk);
547 } else {
548 sectorTrack.Set(trk, inTrack, alpha, sector);
549 }
550 }
551 return 0;
552}
553
554GPUd() void GPUTPCGMMerger::SetTrackClusterZT(GPUTPCGMSectorTrack& track, int32_t iSector, const GPUTPCTrack* sectorTr)
555{
556 const GPUTPCTracker& trk = GetConstantMem()->tpcTrackers[iSector];
557 const GPUTPCHitId& ic1 = trk.TrackHits()[sectorTr->FirstHitID()];
558 const GPUTPCHitId& ic2 = trk.TrackHits()[sectorTr->FirstHitID() + sectorTr->NHits() - 1];
559 int32_t clusterIndex1 = trk.Data().ClusterDataIndex(trk.Data().Row(ic1.RowIndex()), ic1.HitIndex());
560 int32_t clusterIndex2 = trk.Data().ClusterDataIndex(trk.Data().Row(ic2.RowIndex()), ic2.HitIndex());
561 if (Param().par.earlyTpcTransform) {
562 track.SetClusterZT(trk.Data().ClusterData()[clusterIndex1].z, trk.Data().ClusterData()[clusterIndex2].z);
563 } else {
564 const ClusterNative* cl = GetConstantMem()->ioPtrs.clustersNative->clustersLinear + GetConstantMem()->ioPtrs.clustersNative->clusterOffset[iSector][0];
565 track.SetClusterZT(cl[clusterIndex1].getTime(), cl[clusterIndex2].getTime());
566 }
567}
568
569GPUd() void GPUTPCGMMerger::UnpackSaveNumber(int32_t id)
570{
571 mSectorTrackInfoIndex[id] = mMemory->nUnpackedTracks;
572}
573
574GPUd() void GPUTPCGMMerger::UnpackSectorGlobal(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t iSector)
575{
576 const GPUTPCTracker& trk = GetConstantMem()->tpcTrackers[iSector];
577 float alpha = Param().Alpha(iSector);
578 const GPUTPCTrack* sectorTr = mMemory->firstExtrapolatedTracks[iSector];
579 uint32_t nLocalTracks = trk.CommonMemory()->nLocalTracks;
580 uint32_t nTracks = *trk.NTracks();
581 for (uint32_t itr = nLocalTracks + iBlock * nThreads + iThread; itr < nTracks; itr += nBlocks * nThreads) {
582 sectorTr = &trk.Tracks()[itr];
583 int32_t localId = mTrackIDs[(sectorTr->LocalTrackId() >> 24) * mNMaxSingleSectorTracks + (sectorTr->LocalTrackId() & 0xFFFFFF)];
584 if (localId == -1) {
585 continue;
586 }
587 uint32_t myTrack = CAMath::AtomicAdd(&mMemory->nUnpackedTracks, 1u);
588 GPUTPCGMSectorTrack& track = mSectorTrackInfos[myTrack];
589 SetTrackClusterZT(track, iSector, sectorTr);
590 track.Set(this, sectorTr, alpha, iSector);
591 track.SetGlobalSectorTrackCov();
592 track.SetPrevNeighbour(-1);
593 track.SetNextNeighbour(-1);
594 track.SetNextSegmentNeighbour(-1);
595 track.SetPrevSegmentNeighbour(-1);
596 track.SetLocalTrackId(localId);
597 }
598}
599
600GPUd() void GPUTPCGMMerger::UnpackResetIds(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t iSector)
601{
602 const GPUTPCTracker& trk = GetConstantMem()->tpcTrackers[iSector];
603 uint32_t nLocalTracks = trk.CommonMemory()->nLocalTracks;
604 for (uint32_t i = iBlock * nThreads + iThread; i < nLocalTracks; i += nBlocks * nThreads) {
605 mTrackIDs[iSector * mNMaxSingleSectorTracks + i] = -1;
606 }
607}
608
609GPUd() void GPUTPCGMMerger::RefitSectorTracks(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t iSector)
610{
611 const GPUTPCTracker& trk = GetConstantMem()->tpcTrackers[iSector];
612 uint32_t nLocalTracks = trk.CommonMemory()->nLocalTracks;
613
614 float alpha = Param().Alpha(iSector);
615 const GPUTPCTrack* sectorTr = nullptr;
616
617 for (uint32_t itr = iBlock * nThreads + iThread; itr < nLocalTracks; itr += nBlocks * nThreads) {
618 sectorTr = &trk.Tracks()[itr];
620 SetTrackClusterZT(track, iSector, sectorTr);
621 if (Param().rec.tpc.mergerCovSource == 0) {
622 track.Set(this, sectorTr, alpha, iSector);
623 if (!track.FilterErrors(this, iSector, GPUCA_MAX_SIN_PHI, 0.1f)) {
624 continue;
625 }
626 } else if (Param().rec.tpc.mergerCovSource == 1) {
627 track.Set(this, sectorTr, alpha, iSector);
628 track.CopyBaseTrackCov();
629 } else if (Param().rec.tpc.mergerCovSource == 2) {
630 if (RefitSectorTrack(track, sectorTr, alpha, iSector)) {
631 track.Set(this, sectorTr, alpha, iSector); // TODO: Why does the refit fail, it shouldn't, this workaround should be removed
632 if (!track.FilterErrors(this, iSector, GPUCA_MAX_SIN_PHI, 0.1f)) {
633 continue;
634 }
635 }
636 }
637
638 CADEBUG(GPUInfo("INPUT Sector %d, Track %u, QPt %f DzDs %f", iSector, itr, track.QPt(), track.DzDs()));
639 track.SetPrevNeighbour(-1);
640 track.SetNextNeighbour(-1);
641 track.SetNextSegmentNeighbour(-1);
642 track.SetPrevSegmentNeighbour(-1);
643 track.SetExtrapolatedTrackId(0, -1);
644 track.SetExtrapolatedTrackId(1, -1);
645 uint32_t myTrack = CAMath::AtomicAdd(&mMemory->nUnpackedTracks, 1u);
646 mTrackIDs[iSector * mNMaxSingleSectorTracks + sectorTr->LocalTrackId()] = myTrack;
647 mSectorTrackInfos[myTrack] = track;
648 }
649}
650
651GPUd() void GPUTPCGMMerger::LinkExtrapolatedTracks(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
652{
653 for (int32_t itr = SectorTrackInfoGlobalFirst(0) + iBlock * nThreads + iThread; itr < SectorTrackInfoGlobalLast(NSECTORS - 1); itr += nThreads * nBlocks) {
654 GPUTPCGMSectorTrack& extrapolatedTrack = mSectorTrackInfos[itr];
655 GPUTPCGMSectorTrack& localTrack = mSectorTrackInfos[extrapolatedTrack.LocalTrackId()];
656 if (localTrack.ExtrapolatedTrackId(0) != -1 || !CAMath::AtomicCAS(&localTrack.ExtrapolatedTrackIds()[0], -1, itr)) {
657 localTrack.SetExtrapolatedTrackId(1, itr);
658 }
659 }
660}
661
662GPUd() void GPUTPCGMMerger::MergeSectorsPrepareStep2(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t iBorder, GPUTPCGMBorderTrack** B, GPUAtomic(uint32_t) * nB, bool useOrigTrackParam)
663{
664 //* prepare sector tracks for merging with next/previous/same sector
665 //* each track transported to the border line
666
667 float fieldBz = Param().bzCLight;
668
669 float dAlpha = Param().par.dAlpha / 2;
670 float x0 = 0;
671
672 if (iBorder == 0) { // transport to the left edge of the sector and rotate horizontally
673 dAlpha = dAlpha - CAMath::Pi() / 2;
674 } else if (iBorder == 1) { // transport to the right edge of the sector and rotate horizontally
675 dAlpha = -dAlpha - CAMath::Pi() / 2;
676 } else if (iBorder == 2) { // transport to the middle of the sector and rotate vertically to the border on the left
677 x0 = GPUTPCGeometry::Row2X(63);
678 } else if (iBorder == 3) { // transport to the middle of the sector and rotate vertically to the border on the right
679 dAlpha = -dAlpha;
680 x0 = GPUTPCGeometry::Row2X(63);
681 } else if (iBorder == 4) { // transport to the middle of the sßector, w/o rotation
682 dAlpha = 0;
683 x0 = GPUTPCGeometry::Row2X(63);
684 }
685
686 const float maxSin = CAMath::Sin(60.f / 180.f * CAMath::Pi());
687 float cosAlpha = CAMath::Cos(dAlpha);
688 float sinAlpha = CAMath::Sin(dAlpha);
689
690 GPUTPCGMSectorTrack trackTmp;
691 for (int32_t itr = iBlock * nThreads + iThread; itr < SectorTrackInfoLocalTotal(); itr += nThreads * nBlocks) {
692 const GPUTPCGMSectorTrack* track = &mSectorTrackInfos[itr];
693 int32_t iSector = track->Sector();
694
695 if (track->PrevSegmentNeighbour() >= 0 && track->Sector() == mSectorTrackInfos[track->PrevSegmentNeighbour()].Sector()) {
696 continue;
697 }
698 if (useOrigTrackParam) { // TODO: Check how far this makes sense with sector track refit
699 if (CAMath::Abs(track->QPt()) * Param().qptB5Scaler < Param().rec.tpc.mergerLooperQPtB5Limit) {
700 continue;
701 }
702 const GPUTPCGMSectorTrack* trackMin = track;
703 while (track->NextSegmentNeighbour() >= 0 && track->Sector() == mSectorTrackInfos[track->NextSegmentNeighbour()].Sector()) {
704 track = &mSectorTrackInfos[track->NextSegmentNeighbour()];
705 if (track->OrigTrack()->Param().X() < trackMin->OrigTrack()->Param().X()) {
706 trackMin = track;
707 }
708 }
709 trackTmp = *trackMin;
710 track = &trackTmp;
711 if (Param().rec.tpc.mergerCovSource == 2 && trackTmp.X2() != 0.f) {
712 trackTmp.UseParam2();
713 } else {
714 trackTmp.Set(this, trackMin->OrigTrack(), trackMin->Alpha(), trackMin->Sector());
715 }
716 } else {
717 if (CAMath::Abs(track->QPt()) * Param().qptB5Scaler < Param().rec.tpc.mergerLooperSecondHorizontalQPtB5Limit) {
718 if (iBorder == 0 && track->NextNeighbour() >= 0) {
719 continue;
720 }
721 if (iBorder == 1 && track->PrevNeighbour() >= 0) {
722 continue;
723 }
724 }
725 }
727
728 if (track->TransportToXAlpha(this, x0, sinAlpha, cosAlpha, fieldBz, b, maxSin)) {
729 b.SetTrackID(itr);
730 b.SetNClusters(track->NClusters());
731 for (int32_t i = 0; i < 4; i++) {
732 if (CAMath::Abs(b.Cov()[i]) >= 5.0f) {
733 b.SetCov(i, 5.0f);
734 }
735 }
736 if (CAMath::Abs(b.Cov()[4]) >= 0.5f) {
737 b.SetCov(4, 0.5f);
738 }
739 uint32_t myTrack = CAMath::AtomicAdd(&nB[iSector], 1u);
740 B[iSector][myTrack] = b;
741 }
742 }
743}
744
745template <>
746GPUd() void GPUTPCGMMerger::MergeBorderTracks<0>(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t iSector1, GPUTPCGMBorderTrack* B1, int32_t N1, int32_t iSector2, GPUTPCGMBorderTrack* B2, int32_t N2, int32_t mergeMode)
747{
748 CADEBUG(GPUInfo("\nMERGING Sectors %d %d NTracks %d %d CROSS %d", iSector1, iSector2, N1, N2, mergeMode));
749 GPUTPCGMBorderRange* range1 = mBorderRange[iSector1];
750 GPUTPCGMBorderRange* range2 = mBorderRange[iSector2] + *GetConstantMem()->tpcTrackers[iSector2].NTracks();
751 bool sameSector = (iSector1 == iSector2);
752 for (int32_t itr = iBlock * nThreads + iThread; itr < N1; itr += nThreads * nBlocks) {
753 GPUTPCGMBorderTrack& b = B1[itr];
754 float d = CAMath::Max(0.5f, 3.5f * CAMath::Sqrt(b.Cov()[1]));
755 if (CAMath::Abs(b.Par()[4]) * Param().qptB5Scaler >= 20) {
756 d *= 2;
757 } else if (d > 3) {
758 d = 3;
759 }
760 CADEBUG(printf(" Input Sector 1 %d Track %d: ", iSector1, itr); for (int32_t i = 0; i < 5; i++) { printf("%8.3f ", b.Par()[i]); } printf(" - "); for (int32_t i = 0; i < 5; i++) { printf("%8.3f ", b.Cov()[i]); } printf(" - D %8.3f\n", d));
762 range.fId = itr;
763 range.fMin = b.Par()[1] + b.ZOffsetLinear() - d;
764 range.fMax = b.Par()[1] + b.ZOffsetLinear() + d;
765 range1[itr] = range;
766 if (sameSector) {
767 range2[itr] = range;
768 }
769 }
770 if (!sameSector) {
771 for (int32_t itr = iBlock * nThreads + iThread; itr < N2; itr += nThreads * nBlocks) {
772 GPUTPCGMBorderTrack& b = B2[itr];
773 float d = CAMath::Max(0.5f, 3.5f * CAMath::Sqrt(b.Cov()[1]));
774 if (CAMath::Abs(b.Par()[4]) * Param().qptB5Scaler >= 20) {
775 d *= 2;
776 } else if (d > 3) {
777 d = 3;
778 }
779 CADEBUG(printf(" Input Sector 2 %d Track %d: ", iSector2, itr); for (int32_t i = 0; i < 5; i++) { printf("%8.3f ", b.Par()[i]); } printf(" - "); for (int32_t i = 0; i < 5; i++) { printf("%8.3f ", b.Cov()[i]); } printf(" - D %8.3f\n", d));
781 range.fId = itr;
782 range.fMin = b.Par()[1] + b.ZOffsetLinear() - d;
783 range.fMax = b.Par()[1] + b.ZOffsetLinear() + d;
784 range2[itr] = range;
785 }
786 }
787}
788
789template <>
790GPUd() void GPUTPCGMMerger::MergeBorderTracks<1>(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t iSector1, GPUTPCGMBorderTrack* B1, int32_t N1, int32_t iSector2, GPUTPCGMBorderTrack* B2, int32_t N2, int32_t mergeMode)
791{
792#if !defined(GPUCA_GPUCODE_COMPILEKERNELS)
793 GPUTPCGMBorderRange* range1 = mBorderRange[iSector1];
794 GPUTPCGMBorderRange* range2 = mBorderRange[iSector2] + *GetConstantMem()->tpcTrackers[iSector2].NTracks();
795
796 if (iThread == 0) {
797 if (iBlock == 0) {
798 GPUCommonAlgorithm::sortDeviceDynamic(range1, range1 + N1, [](const GPUTPCGMBorderRange& a, const GPUTPCGMBorderRange& b) { return GPUCA_DETERMINISTIC_CODE((a.fMin != b.fMin) ? (a.fMin < b.fMin) : (a.fId < b.fId), a.fMin < b.fMin); });
799 } else if (iBlock == 1) {
800 GPUCommonAlgorithm::sortDeviceDynamic(range2, range2 + N2, [](const GPUTPCGMBorderRange& a, const GPUTPCGMBorderRange& b) { return GPUCA_DETERMINISTIC_CODE((a.fMax != b.fMax) ? (a.fMax < b.fMax) : (a.fId < b.fId), a.fMax < b.fMax); });
801 }
802 }
803#else
804 printf("This sorting variant is disabled for RTC");
805#endif
806}
807
808template <>
809GPUd() void GPUTPCGMMerger::MergeBorderTracks<3>(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, GPUTPCGMBorderRange* range, int32_t N, int32_t cmpMax)
810{
811#ifndef GPUCA_SPECIALIZE_THRUST_SORTS
812 if (iThread == 0 && iBlock == 0) {
813 if (cmpMax) {
814 GPUCommonAlgorithm::sortDeviceDynamic(range, range + N, MergeBorderTracks_compMax());
815 } else {
816 GPUCommonAlgorithm::sortDeviceDynamic(range, range + N, MergeBorderTracks_compMin());
817 }
818 }
819#endif
820}
821
822template <>
823GPUd() void GPUTPCGMMerger::MergeBorderTracks<2>(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t iSector1, GPUTPCGMBorderTrack* B1, int32_t N1, int32_t iSector2, GPUTPCGMBorderTrack* B2, int32_t N2, int32_t mergeMode)
824{
825 // int32_t statAll = 0, statMerged = 0;
826 float factor2ys = Param().rec.tpc.trackMergerFactor2YS;
827 float factor2zt = Param().rec.tpc.trackMergerFactor2ZT;
828 float factor2k = Param().rec.tpc.trackMergerFactor2K;
829 float factor2General = Param().rec.tpc.trackMergerFactor2General;
830
831 factor2k = factor2General * factor2k;
832 factor2ys = factor2General * factor2ys;
833 factor2zt = factor2General * factor2zt;
834
835 int32_t minNPartHits = Param().rec.tpc.trackMergerMinPartHits;
836 int32_t minNTotalHits = Param().rec.tpc.trackMergerMinTotalHits;
837
838 bool sameSector = (iSector1 == iSector2);
839
840 GPUTPCGMBorderRange* range1 = mBorderRange[iSector1];
841 GPUTPCGMBorderRange* range2 = mBorderRange[iSector2] + *GetConstantMem()->tpcTrackers[iSector2].NTracks();
842
843 int32_t i2 = 0;
844 for (int32_t i1 = iBlock * nThreads + iThread; i1 < N1; i1 += nThreads * nBlocks) {
845 GPUTPCGMBorderRange r1 = range1[i1];
846 while (i2 < N2 && range2[i2].fMax < r1.fMin) {
847 i2++;
848 }
849
851 if (b1.NClusters() < minNPartHits) {
852 continue;
853 }
854 int32_t iBest2 = -1;
855 int32_t lBest2 = 0;
856 // statAll++;
857 for (int32_t k2 = i2; k2 < N2; k2++) {
858 GPUTPCGMBorderRange r2 = range2[k2];
859 if (r2.fMin > r1.fMax) {
860 break;
861 }
862 if (sameSector && (r1.fId >= r2.fId)) {
863 continue;
864 }
865 // do check
866
867 GPUTPCGMBorderTrack& b2 = B2[r2.fId];
868#if defined(GPUCA_MERGER_BY_MC_LABEL) && !defined(GPUCA_GPUCODE)
869 int64_t label1 = GetTrackLabel(b1);
870 int64_t label2 = GetTrackLabel(b2);
871 if (label1 != label2 && label1 != -1) // DEBUG CODE, match by MC label
872#endif
873 {
874 CADEBUG(if (GetConstantMem()->ioPtrs.mcLabelsTPC) {printf("Comparing track %3d to %3d: ", r1.fId, r2.fId); for (int32_t i = 0; i < 5; i++) { printf("%8.3f ", b1.Par()[i]); } printf(" - "); for (int32_t i = 0; i < 5; i++) { printf("%8.3f ", b1.Cov()[i]); } printf("\n%28s", ""); });
875 CADEBUG(if (GetConstantMem()->ioPtrs.mcLabelsTPC) {for (int32_t i = 0; i < 5; i++) { printf("%8.3f ", b2.Par()[i]); } printf(" - "); for (int32_t i = 0; i < 5; i++) { printf("%8.3f ", b2.Cov()[i]); } printf(" - %5s - ", GetTrackLabel(b1) == GetTrackLabel(b2) ? "CLONE" : "FAKE"); });
876 if (b2.NClusters() < lBest2) {
877 CADEBUG2(continue, printf("!NCl1\n"));
878 }
879 if (mergeMode > 0) {
880 // Merging CE tracks
881 int32_t maxRowDiff = mergeMode == 2 ? 1 : 3; // TODO: check cut
882 if (CAMath::Abs(b1.Row() - b2.Row()) > maxRowDiff) {
883 CADEBUG2(continue, printf("!ROW\n"));
884 }
885 if (CAMath::Abs(b1.Par()[2] - b2.Par()[2]) > 0.5f || CAMath::Abs(b1.Par()[3] - b2.Par()[3]) > 0.5f) {
886 CADEBUG2(continue, printf("!CE SinPhi/Tgl\n")); // Crude cut to avoid totally wrong matches, TODO: check cut
887 }
888 }
889
890 GPUCA_DEBUG_STREAMER_CHECK(float weight = b1.Par()[4] * b1.Par()[4]; if (o2::utils::DebugStreamer::checkStream(o2::utils::StreamFlags::streamMergeBorderTracksAll, b1.TrackID(), weight)) { MergedTrackStreamer(b1, b2, "merge_all_tracks", iSector1, iSector2, mergeMode, weight, o2::utils::DebugStreamer::getSamplingFrequency(o2::utils::StreamFlags::streamMergeBorderTracksAll)); });
891
892 if (!b1.CheckChi2Y(b2, factor2ys)) {
893 CADEBUG2(continue, printf("!Y\n"));
894 }
895 // if( !b1.CheckChi2Z(b2, factor2zt ) ) CADEBUG2(continue, printf("!NCl1\n"));
896 if (!b1.CheckChi2QPt(b2, factor2k)) {
897 CADEBUG2(continue, printf("!QPt\n"));
898 }
899 float fys = CAMath::Abs(b1.Par()[4]) * Param().qptB5Scaler < 20 ? factor2ys : (2.f * factor2ys);
900 float fzt = CAMath::Abs(b1.Par()[4]) * Param().qptB5Scaler < 20 ? factor2zt : (2.f * factor2zt);
901 if (!b1.CheckChi2YS(b2, fys)) {
902 CADEBUG2(continue, printf("!YS\n"));
903 }
904 if (!b1.CheckChi2ZT(b2, fzt)) {
905 CADEBUG2(continue, printf("!ZT\n"));
906 }
907 if (CAMath::Abs(b1.Par()[4]) * Param().qptB5Scaler < 20) {
908 if (b2.NClusters() < minNPartHits) {
909 CADEBUG2(continue, printf("!NCl2\n"));
910 }
911 if (b1.NClusters() + b2.NClusters() < minNTotalHits) {
912 CADEBUG2(continue, printf("!NCl3\n"));
913 }
914 }
915 CADEBUG(printf("OK: dZ %8.3f D1 %8.3f D2 %8.3f\n", CAMath::Abs(b1.Par()[1] - b2.Par()[1]), 3.5f * sqrt(b1.Cov()[1]), 3.5f * sqrt(b2.Cov()[1])));
916 } // DEBUG CODE, match by MC label
917 lBest2 = b2.NClusters();
918 iBest2 = b2.TrackID();
919 }
920
921 if (iBest2 < 0) {
922 continue;
923 }
924 GPUCA_DEBUG_STREAMER_CHECK(float weight = b1.Par()[4] * b1.Par()[4]; if (o2::utils::DebugStreamer::checkStream(o2::utils::StreamFlags::streamMergeBorderTracksBest, b1.TrackID(), weight)) { MergedTrackStreamer(b1, MergedTrackStreamerFindBorderTrack(B2, N2, iBest2), "merge_best_track", iSector1, iSector2, mergeMode, weight, o2::utils::DebugStreamer::getSamplingFrequency(o2::utils::StreamFlags::streamMergeBorderTracksBest)); });
925
926 // statMerged++;
927
928 CADEBUG(GPUInfo("Found match %d %d", b1.TrackID(), iBest2));
929
930 mTrackLinks[b1.TrackID()] = iBest2;
931 if (mergeMode > 0) {
932 GPUCA_DETERMINISTIC_CODE(CAMath::AtomicMax(&mTrackLinks[iBest2], b1.TrackID()), mTrackLinks[iBest2] = b1.TrackID());
933 }
934 }
935 // GPUInfo("STAT: sectors %d, %d: all %d merged %d", iSector1, iSector2, statAll, statMerged);
936}
937
938GPUdii() void GPUTPCGMMerger::MergeBorderTracksSetup(int32_t& n1, int32_t& n2, GPUTPCGMBorderTrack*& b1, GPUTPCGMBorderTrack*& b2, int32_t& jSector, int32_t iSector, int8_t withinSector, int8_t mergeMode) const
939{
940 if (withinSector == 1) { // Merge tracks within the same sector
941 jSector = iSector;
942 n1 = n2 = mMemory->tmpCounter[iSector];
943 b1 = b2 = mBorder[iSector];
944 } else if (withinSector == -1) { // Merge tracks accross the central electrode
945 jSector = (iSector + NSECTORS / 2);
946 const int32_t offset = mergeMode == 2 ? NSECTORS : 0;
947 n1 = mMemory->tmpCounter[iSector + offset];
948 n2 = mMemory->tmpCounter[jSector + offset];
949 b1 = mBorder[iSector + offset];
950 b2 = mBorder[jSector + offset];
951 } else { // Merge tracks of adjacent sectors
952 jSector = mNextSectorInd[iSector];
953 n1 = mMemory->tmpCounter[iSector];
954 n2 = mMemory->tmpCounter[NSECTORS + jSector];
955 b1 = mBorder[iSector];
956 b2 = mBorder[NSECTORS + jSector];
957 }
958}
959
960template <int32_t I>
961GPUd() void GPUTPCGMMerger::MergeBorderTracks(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t iSector, int8_t withinSector, int8_t mergeMode)
962{
963 int32_t n1, n2;
965 int32_t jSector;
966 MergeBorderTracksSetup(n1, n2, b1, b2, jSector, iSector, withinSector, mergeMode);
967 MergeBorderTracks<I>(nBlocks, nThreads, iBlock, iThread, iSector, b1, n1, jSector, b2, n2, mergeMode);
968}
969
970#if !defined(GPUCA_GPUCODE) || defined(GPUCA_GPUCODE_DEVICE) // FIXME: DR: WORKAROUND to avoid CUDA bug creating host symbols for device code.
971template GPUdni() void GPUTPCGMMerger::MergeBorderTracks<0>(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t iSector, int8_t withinSector, int8_t mergeMode);
972template GPUdni() void GPUTPCGMMerger::MergeBorderTracks<1>(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t iSector, int8_t withinSector, int8_t mergeMode);
973template GPUdni() void GPUTPCGMMerger::MergeBorderTracks<2>(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t iSector, int8_t withinSector, int8_t mergeMode);
974#endif
975
976GPUd() void GPUTPCGMMerger::MergeWithinSectorsPrepare(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
977{
978 float x0 = GPUTPCGeometry::Row2X(63);
979 const float maxSin = CAMath::Sin(60.f / 180.f * CAMath::Pi());
980
981 for (int32_t itr = iBlock * nThreads + iThread; itr < SectorTrackInfoLocalTotal(); itr += nThreads * nBlocks) {
982 GPUTPCGMSectorTrack& track = mSectorTrackInfos[itr];
983 int32_t iSector = track.Sector();
985 if (track.TransportToX(this, x0, Param().bzCLight, b, maxSin)) {
986 b.SetTrackID(itr);
987 CADEBUG(printf("WITHIN SECTOR %d Track %d - ", iSector, itr); for (int32_t i = 0; i < 5; i++) { printf("%8.3f ", b.Par()[i]); } printf(" - "); for (int32_t i = 0; i < 5; i++) { printf("%8.3f ", b.Cov()[i]); } printf("\n"));
988 b.SetNClusters(track.NClusters());
989 uint32_t myTrack = CAMath::AtomicAdd(&mMemory->tmpCounter[iSector], 1u);
990 mBorder[iSector][myTrack] = b;
991 }
992 }
993}
994
995GPUd() void GPUTPCGMMerger::MergeSectorsPrepare(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int32_t border0, int32_t border1, int8_t useOrigTrackParam)
996{
997 bool part2 = iBlock & 1;
998 int32_t border = part2 ? border1 : border0;
999 GPUAtomic(uint32_t)* n = mMemory->tmpCounter;
1000 GPUTPCGMBorderTrack** b = mBorder;
1001 if (part2) {
1002 n += NSECTORS;
1003 b += NSECTORS;
1004 }
1005 MergeSectorsPrepareStep2((nBlocks + !part2) >> 1, nThreads, iBlock >> 1, iThread, border, b, n, useOrigTrackParam);
1006}
1007
1008GPUdi() void GPUTPCGMMerger::setBlockRange(int32_t elems, int32_t nBlocks, int32_t iBlock, int32_t& start, int32_t& end)
1009{
1010 start = (elems + nBlocks - 1) / nBlocks * iBlock;
1011 end = (elems + nBlocks - 1) / nBlocks * (iBlock + 1);
1012 end = CAMath::Min(elems, end);
1013}
1014
1015GPUd() void GPUTPCGMMerger::hookEdge(int32_t u, int32_t v)
1016{
1017 if (v < 0) {
1018 return;
1019 }
1020 while (true) {
1021 u = mTrackCCRoots[u];
1022 v = mTrackCCRoots[v];
1023 if (u == v) {
1024 break;
1025 }
1026 int32_t h = CAMath::Max(u, v);
1027 int32_t l = CAMath::Min(u, v);
1028
1029 int32_t old = CAMath::AtomicCAS(&mTrackCCRoots[h], h, l);
1030 if (old == h) {
1031 break;
1032 }
1033
1034 u = mTrackCCRoots[h];
1035 v = l;
1036 }
1037}
1038
1039GPUd() void GPUTPCGMMerger::ResolveFindConnectedComponentsSetup(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1040{
1041 int32_t start, end;
1042 setBlockRange(SectorTrackInfoLocalTotal(), nBlocks, iBlock, start, end);
1043 for (int32_t i = start + iThread; i < end; i += nThreads) {
1044 mTrackCCRoots[i] = i;
1045 }
1046}
1047
1048GPUd() void GPUTPCGMMerger::ResolveFindConnectedComponentsHookLinks(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1049{
1050 // Compute connected components in parallel, step 1.
1051 // Source: Adaptive Work-Efficient Connected Components on the GPU, Sutton et al, 2016 (https://arxiv.org/pdf/1612.01178.pdf)
1052 int32_t start, end;
1053 setBlockRange(SectorTrackInfoLocalTotal(), nBlocks, iBlock, start, end);
1054 for (int32_t itr = start + iThread; itr < end; itr += nThreads) {
1055 hookEdge(itr, mTrackLinks[itr]);
1056 }
1057}
1058
1059GPUd() void GPUTPCGMMerger::ResolveFindConnectedComponentsHookNeighbors(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1060{
1061 // Compute connected components in parallel, step 1 - Part 2.
1062 nBlocks = nBlocks / 4 * 4;
1063 if (iBlock >= nBlocks) {
1064 return;
1065 }
1066
1067 int32_t start, end;
1068 setBlockRange(SectorTrackInfoLocalTotal(), nBlocks / 4, iBlock / 4, start, end);
1069
1070 int32_t myNeighbor = iBlock % 4;
1071
1072 for (int32_t itr = start + iThread; itr < end; itr += nThreads) {
1073 int32_t v = mSectorTrackInfos[itr].AnyNeighbour(myNeighbor);
1074 hookEdge(itr, v);
1075 }
1076}
1077
1078GPUd() void GPUTPCGMMerger::ResolveFindConnectedComponentsMultiJump(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1079{
1080 // Compute connected components in parallel, step 2.
1081 int32_t start, end;
1082 setBlockRange(SectorTrackInfoLocalTotal(), nBlocks, iBlock, start, end);
1083 for (int32_t itr = start + iThread; itr < end; itr += nThreads) {
1084 int32_t root = itr;
1085 int32_t next = mTrackCCRoots[root];
1086 if (root == next) {
1087 continue;
1088 }
1089 do {
1090 root = next;
1091 next = mTrackCCRoots[next];
1092 } while (root != next);
1093 mTrackCCRoots[itr] = root;
1094 }
1095}
1096
1097GPUd() void GPUTPCGMMerger::ResolveMergeSectors(GPUResolveSharedMemory& smem, int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread, int8_t useOrigTrackParam, int8_t mergeAll)
1098{
1099 if (!mergeAll) {
1100 /*int32_t neighborType = useOrigTrackParam ? 1 : 0;
1101 int32_t old1 = newTrack2.PrevNeighbour(0);
1102 int32_t old2 = newTrack1.NextNeighbour(0);
1103 if (old1 < 0 && old2 < 0) neighborType = 0;
1104 if (old1 == itr) continue;
1105 if (neighborType) old1 = newTrack2.PrevNeighbour(1);
1106 if ( old1 >= 0 )
1107 {
1108 GPUTPCGMSectorTrack &oldTrack1 = mSectorTrackInfos[old1];
1109 if ( oldTrack1.NClusters() < newTrack1.NClusters() ) {
1110 newTrack2.SetPrevNeighbour( -1, neighborType );
1111 oldTrack1.SetNextNeighbour( -1, neighborType );
1112 } else continue;
1113 }
1114
1115 if (old2 == itr2) continue;
1116 if (neighborType) old2 = newTrack1.NextNeighbour(1);
1117 if ( old2 >= 0 )
1118 {
1119 GPUTPCGMSectorTrack &oldTrack2 = mSectorTrackInfos[old2];
1120 if ( oldTrack2.NClusters() < newTrack2.NClusters() )
1121 {
1122 oldTrack2.SetPrevNeighbour( -1, neighborType );
1123 } else continue;
1124 }
1125 newTrack1.SetNextNeighbour( itr2, neighborType );
1126 newTrack2.SetPrevNeighbour( itr, neighborType );*/
1127 }
1128
1129 int32_t start, end;
1130 setBlockRange(SectorTrackInfoLocalTotal(), nBlocks, iBlock, start, end);
1131
1132 for (int32_t baseIdx = 0; baseIdx < SectorTrackInfoLocalTotal(); baseIdx += nThreads) {
1133 int32_t itr = baseIdx + iThread;
1134 bool inRange = itr < SectorTrackInfoLocalTotal();
1135
1136 int32_t itr2 = -1;
1137 if (inRange) {
1138 itr2 = mTrackLinks[itr];
1139 }
1140
1141 bool resolveSector = (itr2 > -1);
1142 if (resolveSector) {
1143 int32_t root = mTrackCCRoots[itr];
1144 resolveSector &= (start <= root) && (root < end);
1145 }
1146
1147 int16_t smemIdx = work_group_scan_inclusive_add(int16_t(resolveSector));
1148
1149 if (resolveSector) {
1150 smem.iTrack1[smemIdx - 1] = itr;
1151 smem.iTrack2[smemIdx - 1] = itr2;
1152 }
1153 GPUbarrier();
1154
1155 if (iThread < nThreads - 1) {
1156 continue;
1157 }
1158
1159 const int32_t nSectors = smemIdx;
1160
1161 for (int32_t i = 0; i < nSectors; i++) {
1162 itr = smem.iTrack1[i];
1163 itr2 = smem.iTrack2[i];
1164
1165 GPUTPCGMSectorTrack* track1 = &mSectorTrackInfos[itr];
1166 GPUTPCGMSectorTrack* track2 = &mSectorTrackInfos[itr2];
1167 GPUTPCGMSectorTrack* track1Base = track1;
1168 GPUTPCGMSectorTrack* track2Base = track2;
1169
1170 bool sameSegment = CAMath::Abs(track1->NClusters() > track2->NClusters() ? track1->QPt() : track2->QPt()) * Param().qptB5Scaler < 2 || track1->QPt() * track2->QPt() > 0;
1171 // GPUInfo("\nMerge %d with %d - same segment %d", itr, itr2, (int32_t) sameSegment);
1172 // PrintMergeGraph(track1, std::cout);
1173 // PrintMergeGraph(track2, std::cout);
1174
1175 while (track2->PrevSegmentNeighbour() >= 0) {
1176 track2 = &mSectorTrackInfos[track2->PrevSegmentNeighbour()];
1177 }
1178 if (sameSegment) {
1179 if (track1 == track2) {
1180 continue;
1181 }
1182 while (track1->PrevSegmentNeighbour() >= 0) {
1183 track1 = &mSectorTrackInfos[track1->PrevSegmentNeighbour()];
1184 if (track1 == track2) {
1185 goto NextTrack;
1186 }
1187 }
1188 GPUCommonAlgorithm::swap(track1, track1Base);
1189 for (int32_t k = 0; k < 2; k++) {
1190 GPUTPCGMSectorTrack* tmp = track1Base;
1191 while (tmp->Neighbour(k) >= 0) {
1192 tmp = &mSectorTrackInfos[tmp->Neighbour(k)];
1193 if (tmp == track2) {
1194 goto NextTrack;
1195 }
1196 }
1197 }
1198
1199 while (track1->NextSegmentNeighbour() >= 0) {
1200 track1 = &mSectorTrackInfos[track1->NextSegmentNeighbour()];
1201 if (track1 == track2) {
1202 goto NextTrack;
1203 }
1204 }
1205 } else {
1206 while (track1->PrevSegmentNeighbour() >= 0) {
1207 track1 = &mSectorTrackInfos[track1->PrevSegmentNeighbour()];
1208 }
1209
1210 if (track1 == track2) {
1211 continue;
1212 }
1213 for (int32_t k = 0; k < 2; k++) {
1214 GPUTPCGMSectorTrack* tmp = track1;
1215 while (tmp->Neighbour(k) >= 0) {
1216 tmp = &mSectorTrackInfos[tmp->Neighbour(k)];
1217 if (tmp == track2) {
1218 goto NextTrack;
1219 }
1220 }
1221 }
1222
1223 float z1min, z1max, z2min, z2max;
1224 z1min = track1->MinClusterZT();
1225 z1max = track1->MaxClusterZT();
1226 z2min = track2->MinClusterZT();
1227 z2max = track2->MaxClusterZT();
1228 if (track1 != track1Base) {
1229 z1min = CAMath::Min(z1min, track1Base->MinClusterZT());
1230 z1max = CAMath::Max(z1max, track1Base->MaxClusterZT());
1231 }
1232 if (track2 != track2Base) {
1233 z2min = CAMath::Min(z2min, track2Base->MinClusterZT());
1234 z2max = CAMath::Max(z2max, track2Base->MaxClusterZT());
1235 }
1236 bool goUp = z2max - z1min > z1max - z2min;
1237
1238 if (track1->Neighbour(goUp) < 0 && track2->Neighbour(!goUp) < 0) {
1239 track1->SetNeighbor(track2 - mSectorTrackInfos, goUp);
1240 track2->SetNeighbor(track1 - mSectorTrackInfos, !goUp);
1241 // GPUInfo("Result (simple neighbor)");
1242 // PrintMergeGraph(track1, std::cout);
1243 continue;
1244 } else if (track1->Neighbour(goUp) < 0) {
1245 track2 = &mSectorTrackInfos[track2->Neighbour(!goUp)];
1246 GPUCommonAlgorithm::swap(track1, track2);
1247 } else if (track2->Neighbour(!goUp) < 0) {
1248 track1 = &mSectorTrackInfos[track1->Neighbour(goUp)];
1249 } else { // Both would work, but we use the simpler one
1250 track1 = &mSectorTrackInfos[track1->Neighbour(goUp)];
1251 }
1252 track1Base = track1;
1253 }
1254
1255 track2Base = track2;
1256 if (!sameSegment) {
1257 while (track1->NextSegmentNeighbour() >= 0) {
1258 track1 = &mSectorTrackInfos[track1->NextSegmentNeighbour()];
1259 }
1260 }
1261 track1->SetNextSegmentNeighbour(track2 - mSectorTrackInfos);
1262 track2->SetPrevSegmentNeighbour(track1 - mSectorTrackInfos);
1263 // k = 0: Merge right side
1264 // k = 1: Merge left side
1265 for (int32_t k = 0; k < 2; k++) {
1266 track1 = track1Base;
1267 track2 = track2Base;
1268 while (track2->Neighbour(k) >= 0) {
1269 if (track1->Neighbour(k) >= 0) {
1270 GPUTPCGMSectorTrack* track1new = &mSectorTrackInfos[track1->Neighbour(k)];
1271 GPUTPCGMSectorTrack* track2new = &mSectorTrackInfos[track2->Neighbour(k)];
1272 track2->SetNeighbor(-1, k);
1273 track2new->SetNeighbor(-1, k ^ 1);
1274 track1 = track1new;
1275 while (track1->NextSegmentNeighbour() >= 0) {
1276 track1 = &mSectorTrackInfos[track1->NextSegmentNeighbour()];
1277 }
1278 track1->SetNextSegmentNeighbour(track2new - mSectorTrackInfos);
1279 track2new->SetPrevSegmentNeighbour(track1 - mSectorTrackInfos);
1280 track1 = track1new;
1281 track2 = track2new;
1282 } else {
1283 GPUTPCGMSectorTrack* track2new = &mSectorTrackInfos[track2->Neighbour(k)];
1284 track1->SetNeighbor(track2->Neighbour(k), k);
1285 track2->SetNeighbor(-1, k);
1286 track2new->SetNeighbor(track1 - mSectorTrackInfos, k ^ 1);
1287 }
1288 }
1289 }
1290 // GPUInfo("Result");
1291 // PrintMergeGraph(track1, std::cout);
1292 NextTrack:;
1293 }
1294 }
1295}
1296
1297GPUd() void GPUTPCGMMerger::MergeCEFill(const GPUTPCGMSectorTrack* track, const GPUTPCGMMergedTrackHit& cls, const GPUTPCGMMergedTrackHitXYZ* clsXYZ, int32_t itr)
1298{
1299 if (Param().rec.tpc.mergerCERowLimit > 0 && CAMath::Abs(track->QPt()) * Param().qptB5Scaler < 0.3f && (cls.row < Param().rec.tpc.mergerCERowLimit || cls.row >= GPUCA_ROW_COUNT - Param().rec.tpc.mergerCERowLimit)) {
1300 return;
1301 }
1302
1303 float z = 0;
1304 if (Param().par.earlyTpcTransform) {
1305 z = clsXYZ->z;
1306 } else {
1307 float x, y;
1308 auto& cln = mConstantMem->ioPtrs.clustersNative->clustersLinear[cls.num];
1309 GPUTPCConvertImpl::convert(*mConstantMem, cls.sector, cls.row, cln.getPad(), cln.getTime(), x, y, z);
1310 }
1311
1312 if (!Param().par.continuousTracking && CAMath::Abs(z) > 10) {
1313 return;
1314 }
1315 int32_t sector = track->Sector();
1316 for (int32_t attempt = 0; attempt < 2; attempt++) {
1318 const float x0 = GPUTPCGeometry::Row2X(attempt == 0 ? 63 : cls.row);
1319 if (track->TransportToX(this, x0, Param().bzCLight, b, GPUCA_MAX_SIN_PHI_LOW)) {
1320 b.SetTrackID(itr);
1321 b.SetNClusters(mMergedTracks[itr].NClusters());
1322 if (CAMath::Abs(b.Cov()[4]) >= 0.5f) {
1323 b.SetCov(4, 0.5f); // TODO: Is this needed and better than the cut in BorderTrack?
1324 }
1325 if (track->CSide()) {
1326 b.SetPar(1, b.Par()[1] - 2 * (z - b.ZOffsetLinear()));
1327 b.SetZOffsetLinear(-b.ZOffsetLinear());
1328 }
1329 b.SetRow(cls.row);
1330 uint32_t id = sector + attempt * NSECTORS;
1331 uint32_t myTrack = CAMath::AtomicAdd(&mMemory->tmpCounter[id], 1u);
1332 mBorder[id][myTrack] = b;
1333 break;
1334 }
1335 }
1336}
1337
1338GPUd() void GPUTPCGMMerger::MergeCE(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1339{
1340 const ClusterNative* cls = Param().par.earlyTpcTransform ? nullptr : mConstantMem->ioPtrs.clustersNative->clustersLinear;
1341 for (uint32_t i = iBlock * nThreads + iThread; i < mMemory->nMergedTracks; i += nThreads * nBlocks) {
1342 if (mMergedTracks[i].CSide() == 0 && mTrackLinks[i] >= 0) {
1343 if (mTrackLinks[mTrackLinks[i]] != (int32_t)i) {
1344 continue;
1345 }
1346 GPUTPCGMMergedTrack* trk[2] = {&mMergedTracks[i], &mMergedTracks[mTrackLinks[i]]};
1347
1348 if (!trk[1]->OK() || trk[1]->CCE()) {
1349 continue;
1350 }
1351 bool celooper = (trk[0]->GetParam().GetQPt() * Param().qptB5Scaler > 1 && trk[0]->GetParam().GetQPt() * trk[1]->GetParam().GetQPt() < 0);
1352 bool looper = trk[0]->Looper() || trk[1]->Looper() || celooper;
1353 if (!looper && trk[0]->GetParam().GetPar(3) * trk[1]->GetParam().GetPar(3) < 0) {
1354 continue;
1355 }
1356
1357 uint32_t newRef = CAMath::AtomicAdd(&mMemory->nOutputTrackClusters, trk[0]->NClusters() + trk[1]->NClusters());
1358 if (newRef + trk[0]->NClusters() + trk[1]->NClusters() >= mNMaxOutputTrackClusters) {
1359 raiseError(GPUErrors::ERROR_MERGER_CE_HIT_OVERFLOW, newRef + trk[0]->NClusters() + trk[1]->NClusters(), mNMaxOutputTrackClusters);
1360 for (uint32_t k = newRef; k < mNMaxOutputTrackClusters; k++) {
1361 mClusters[k].num = 0;
1362 mClusters[k].state = 0;
1363 }
1364 CAMath::AtomicExch(&mMemory->nOutputTrackClusters, mNMaxOutputTrackClusters);
1365 return;
1366 }
1367
1368 bool needswap = false;
1369 if (looper) {
1370 float z0max, z1max;
1371 if (Param().par.earlyTpcTransform) {
1372 z0max = CAMath::Max(CAMath::Abs(mClustersXYZ[trk[0]->FirstClusterRef()].z), CAMath::Abs(mClustersXYZ[trk[0]->FirstClusterRef() + trk[0]->NClusters() - 1].z));
1373 z1max = CAMath::Max(CAMath::Abs(mClustersXYZ[trk[1]->FirstClusterRef()].z), CAMath::Abs(mClustersXYZ[trk[1]->FirstClusterRef() + trk[1]->NClusters() - 1].z));
1374 } else {
1375 z0max = -CAMath::Min(cls[mClusters[trk[0]->FirstClusterRef()].num].getTime(), cls[mClusters[trk[0]->FirstClusterRef() + trk[0]->NClusters() - 1].num].getTime());
1376 z1max = -CAMath::Min(cls[mClusters[trk[1]->FirstClusterRef()].num].getTime(), cls[mClusters[trk[1]->FirstClusterRef() + trk[1]->NClusters() - 1].num].getTime());
1377 }
1378 if (z1max < z0max) {
1379 needswap = true;
1380 }
1381 } else {
1382 if (mClusters[trk[0]->FirstClusterRef()].row > mClusters[trk[1]->FirstClusterRef()].row) {
1383 needswap = true;
1384 }
1385 }
1386 if (needswap) {
1387 GPUCommonAlgorithm::swap(trk[0], trk[1]);
1388 }
1389
1390 bool reverse[2] = {false, false};
1391 if (looper) {
1392 if (Param().par.earlyTpcTransform) {
1393 reverse[0] = (mClustersXYZ[trk[0]->FirstClusterRef()].z > mClustersXYZ[trk[0]->FirstClusterRef() + trk[0]->NClusters() - 1].z) ^ (trk[0]->CSide() > 0);
1394 reverse[1] = (mClustersXYZ[trk[1]->FirstClusterRef()].z < mClustersXYZ[trk[1]->FirstClusterRef() + trk[1]->NClusters() - 1].z) ^ (trk[1]->CSide() > 0);
1395 } else {
1396 reverse[0] = cls[mClusters[trk[0]->FirstClusterRef()].num].getTime() < cls[mClusters[trk[0]->FirstClusterRef() + trk[0]->NClusters() - 1].num].getTime();
1397 reverse[1] = cls[mClusters[trk[1]->FirstClusterRef()].num].getTime() > cls[mClusters[trk[1]->FirstClusterRef() + trk[1]->NClusters() - 1].num].getTime();
1398 }
1399 }
1400
1401 if (Param().par.continuousTracking) {
1402 if (Param().par.earlyTpcTransform) {
1403 const float z0 = trk[0]->CSide() ? CAMath::Max(mClustersXYZ[trk[0]->FirstClusterRef()].z, mClustersXYZ[trk[0]->FirstClusterRef() + trk[0]->NClusters() - 1].z) : CAMath::Min(mClustersXYZ[trk[0]->FirstClusterRef()].z, mClustersXYZ[trk[0]->FirstClusterRef() + trk[0]->NClusters() - 1].z);
1404 const float z1 = trk[1]->CSide() ? CAMath::Max(mClustersXYZ[trk[1]->FirstClusterRef()].z, mClustersXYZ[trk[1]->FirstClusterRef() + trk[1]->NClusters() - 1].z) : CAMath::Min(mClustersXYZ[trk[1]->FirstClusterRef()].z, mClustersXYZ[trk[1]->FirstClusterRef() + trk[1]->NClusters() - 1].z);
1405 const float offset = CAMath::Abs(z1) > CAMath::Abs(z0) ? -z0 : z1;
1406 trk[1]->Param().Z() += trk[1]->Param().TZOffset() - offset;
1407 trk[1]->Param().TZOffset() = offset;
1408 } else {
1409 GPUTPCGMMergedTrackHit* clsmax;
1410 const float tmax = CAMath::MaxWithRef(cls[mClusters[trk[0]->FirstClusterRef()].num].getTime(), cls[mClusters[trk[0]->FirstClusterRef() + trk[0]->NClusters() - 1].num].getTime(),
1411 cls[mClusters[trk[1]->FirstClusterRef()].num].getTime(), cls[mClusters[trk[1]->FirstClusterRef() + trk[1]->NClusters() - 1].num].getTime(),
1412 &mClusters[trk[0]->FirstClusterRef()], &mClusters[trk[0]->FirstClusterRef() + trk[0]->NClusters() - 1],
1413 &mClusters[trk[1]->FirstClusterRef()], &mClusters[trk[1]->FirstClusterRef() + trk[1]->NClusters() - 1], clsmax);
1414 const float offset = CAMath::Max(tmax - mConstantMem->calibObjects.fastTransformHelper->getCorrMap()->getMaxDriftTime(clsmax->sector, clsmax->row, cls[clsmax->num].getPad()), 0.f);
1415 trk[1]->Param().Z() += mConstantMem->calibObjects.fastTransformHelper->getCorrMap()->convDeltaTimeToDeltaZinTimeFrame(trk[1]->CSide() * NSECTORS / 2, trk[1]->Param().TZOffset() - offset);
1416 trk[1]->Param().TZOffset() = offset;
1417 }
1418 }
1419
1420 int32_t pos = newRef;
1421 int32_t leg = -1;
1422 int32_t lastLeg = -1;
1423#pragma unroll
1424 for (int32_t k = 1; k >= 0; k--) {
1425 int32_t loopstart = reverse[k] ? (trk[k]->NClusters() - 1) : 0;
1426 int32_t loopend = reverse[k] ? -1 : (int32_t)trk[k]->NClusters();
1427 int32_t loopinc = reverse[k] ? -1 : 1;
1428 for (int32_t j = loopstart; j != loopend; j += loopinc) {
1429 if (Param().par.earlyTpcTransform) {
1430 mClustersXYZ[pos] = mClustersXYZ[trk[k]->FirstClusterRef() + j];
1431 }
1432 mClusters[pos] = mClusters[trk[k]->FirstClusterRef() + j];
1433 if (looper) {
1434 if (mClusters[trk[k]->FirstClusterRef() + j].leg != lastLeg) {
1435 leg++;
1436 lastLeg = mClusters[trk[k]->FirstClusterRef() + j].leg;
1437 }
1438 mClusters[pos].leg = leg;
1439 }
1440 pos++;
1441 }
1442 if (celooper) {
1443 lastLeg = -1;
1444 }
1445 }
1446 trk[1]->SetFirstClusterRef(newRef);
1447 trk[1]->SetNClusters(trk[0]->NClusters() + trk[1]->NClusters());
1448 if (trk[1]->NClusters() > GPUCA_MERGER_MAX_TRACK_CLUSTERS) {
1449 trk[1]->SetFirstClusterRef(trk[1]->FirstClusterRef() + trk[1]->NClusters() - GPUCA_MERGER_MAX_TRACK_CLUSTERS);
1450 trk[1]->SetNClusters(GPUCA_MERGER_MAX_TRACK_CLUSTERS);
1451 }
1452 trk[1]->SetCCE(true);
1453 if (looper) {
1454 trk[1]->SetLooper(true);
1455 trk[1]->SetLegs(leg + 1);
1456 }
1457 trk[0]->SetNClusters(0);
1458 trk[0]->SetOK(false);
1459 }
1460 }
1461
1462 // for (int32_t i = 0;i < mMemory->nMergedTracks;i++) {if (mMergedTracks[i].CCE() == false) {mMergedTracks[i].SetNClusters(0);mMergedTracks[i].SetOK(false);}} //Remove all non-CE tracks
1463}
1464
1465namespace o2::gpu::internal
1466{
1467namespace // anonymous
1468{
1469struct GPUTPCGMMerger_CompareClusterIdsLooper {
1470 struct clcomparestruct {
1471 uint8_t leg;
1472 };
1473
1474 const uint8_t leg;
1475 const bool outwards;
1477 const clcomparestruct* const cmp2;
1478 GPUd() GPUTPCGMMerger_CompareClusterIdsLooper(uint8_t l, bool o, const GPUTPCGMMerger::trackCluster* c1, const clcomparestruct* c2) : leg(l), outwards(o), cmp1(c1), cmp2(c2) {}
1479 GPUd() bool operator()(const int16_t aa, const int16_t bb)
1480 {
1481 const clcomparestruct& a = cmp2[aa];
1482 const clcomparestruct& b = cmp2[bb];
1485 if (a.leg != b.leg) {
1486 return ((leg > 0) ^ (a.leg > b.leg));
1487 }
1488 if (a1.row != b1.row) {
1489 return ((a1.row > b1.row) ^ ((a.leg - leg) & 1) ^ outwards);
1490 }
1491 return GPUCA_DETERMINISTIC_CODE((a1.id != b1.id) ? (a1.id > b1.id) : (aa > bb), a1.id > b1.id);
1492 }
1493};
1494
1495struct GPUTPCGMMerger_CompareClusterIds {
1497 GPUd() GPUTPCGMMerger_CompareClusterIds(const GPUTPCGMMerger::trackCluster* cmp) : mCmp(cmp) {}
1498 GPUd() bool operator()(const int16_t aa, const int16_t bb)
1499 {
1502 if (a.row != b.row) {
1503 return (a.row > b.row);
1504 }
1505 return GPUCA_DETERMINISTIC_CODE((a.id != b.id) ? (a.id > b.id) : (aa > bb), a.id > b.id);
1506 }
1507};
1508} // anonymous namespace
1509} // namespace o2::gpu::internal
1510
1511GPUd() void GPUTPCGMMerger::CollectMergedTracks(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1512{
1513 GPUTPCGMSectorTrack* trackParts[kMaxParts];
1514
1515 for (int32_t itr = iBlock * nThreads + iThread; itr < SectorTrackInfoLocalTotal(); itr += nThreads * nBlocks) {
1516
1517 GPUTPCGMSectorTrack& track = mSectorTrackInfos[itr];
1518
1519 if (track.PrevSegmentNeighbour() >= 0) {
1520 continue;
1521 }
1522 if (track.PrevNeighbour() >= 0) {
1523 continue;
1524 }
1525 int32_t nParts = 0;
1526 int32_t nHits = 0;
1527 int32_t leg = 0;
1528 GPUTPCGMSectorTrack *trbase = &track, *tr = &track;
1529 tr->SetPrevSegmentNeighbour(1000000000);
1530 while (true) {
1531 if (nParts >= kMaxParts) {
1532 break;
1533 }
1534 if (nHits + tr->NClusters() > kMaxClusters) {
1535 break;
1536 }
1537 nHits += tr->NClusters();
1538
1539 tr->SetLeg(leg);
1540 trackParts[nParts++] = tr;
1541 for (int32_t i = 0; i < 2; i++) {
1542 if (tr->ExtrapolatedTrackId(i) != -1) {
1543 if (nParts >= kMaxParts) {
1544 break;
1545 }
1546 if (nHits + mSectorTrackInfos[tr->ExtrapolatedTrackId(i)].NClusters() > kMaxClusters) {
1547 break;
1548 }
1549 trackParts[nParts] = &mSectorTrackInfos[tr->ExtrapolatedTrackId(i)];
1550 trackParts[nParts++]->SetLeg(leg);
1551 nHits += mSectorTrackInfos[tr->ExtrapolatedTrackId(i)].NClusters();
1552 }
1553 }
1554 int32_t jtr = tr->NextSegmentNeighbour();
1555 if (jtr >= 0) {
1556 tr = &(mSectorTrackInfos[jtr]);
1557 tr->SetPrevSegmentNeighbour(1000000002);
1558 continue;
1559 }
1560 jtr = trbase->NextNeighbour();
1561 if (jtr >= 0) {
1562 trbase = &(mSectorTrackInfos[jtr]);
1563 tr = trbase;
1564 if (tr->PrevSegmentNeighbour() >= 0) {
1565 break;
1566 }
1567 tr->SetPrevSegmentNeighbour(1000000001);
1568 leg++;
1569 continue;
1570 }
1571 break;
1572 }
1573
1574 // unpack and sort clusters
1575 if (nParts > 1 && leg == 0) {
1576 GPUCommonAlgorithm::sort(trackParts, trackParts + nParts, [](const GPUTPCGMSectorTrack* a, const GPUTPCGMSectorTrack* b) {
1577 GPUCA_DETERMINISTIC_CODE( // clang-format off
1578 if (a->X() != b->X()) {
1579 return (a->X() > b->X());
1580 }
1581 if (a->Y() != b->Y()) {
1582 return (a->Y() > b->Y());
1583 }
1584 if (a->Z() != b->Z()) {
1585 return (a->Z() > b->Z());
1586 }
1587 return a->QPt() > b->QPt();
1588 , // !GPUCA_DETERMINISTIC_CODE
1589 return (a->X() > b->X());
1590 ) // clang-format on
1591 });
1592 }
1593
1594 if (Param().rec.tpc.dropLoopers && leg > 0) {
1595 nParts = 1;
1596 leg = 0;
1597 }
1598
1599 trackCluster trackClusters[kMaxClusters];
1600 nHits = 0;
1601 for (int32_t ipart = 0; ipart < nParts; ipart++) {
1602 const GPUTPCGMSectorTrack* t = trackParts[ipart];
1603 CADEBUG(printf("Collect Track %d Part %d QPt %f DzDs %f\n", mMemory->nMergedTracks, ipart, t->QPt(), t->DzDs()));
1604 int32_t nTrackHits = t->NClusters();
1605 trackCluster* c2 = trackClusters + nHits + nTrackHits - 1;
1606 for (int32_t i = 0; i < nTrackHits; i++, c2--) {
1607 const GPUTPCTracker& trk = GetConstantMem()->tpcTrackers[t->Sector()];
1608 const GPUTPCHitId& ic = trk.TrackHits()[t->OrigTrack()->FirstHitID() + i];
1609 uint32_t id = trk.Data().ClusterDataIndex(trk.Data().Row(ic.RowIndex()), ic.HitIndex()) + GetConstantMem()->ioPtrs.clustersNative->clusterOffset[t->Sector()][0];
1610 *c2 = trackCluster{id, (uint8_t)ic.RowIndex(), t->Sector(), t->Leg()};
1611 }
1612 nHits += nTrackHits;
1613 }
1614 if (nHits < GPUCA_TRACKLET_SELECTOR_MIN_HITS_B5(track.QPt() * Param().qptB5Scaler)) {
1615 continue;
1616 }
1617
1618 int32_t ordered = leg == 0;
1619 if (ordered) {
1620 for (int32_t i = 1; i < nHits; i++) {
1621 if (trackClusters[i].row > trackClusters[i - 1].row || trackClusters[i].id == trackClusters[i - 1].id) {
1622 ordered = 0;
1623 break;
1624 }
1625 }
1626 }
1627 int32_t firstTrackIndex = 0;
1628 int32_t lastTrackIndex = nParts - 1;
1629 if (ordered == 0) {
1630 int32_t nTmpHits = 0;
1631 trackCluster trackClustersUnsorted[kMaxClusters];
1632 int16_t clusterIndices[kMaxClusters];
1633 for (int32_t i = 0; i < nHits; i++) {
1634 trackClustersUnsorted[i] = trackClusters[i];
1635 clusterIndices[i] = i;
1636 }
1637
1638 if (leg > 0) {
1639 // Find QPt and DzDs for the segment closest to the vertex, if low/mid Pt
1640 float baseZT = 1e9;
1641 uint8_t baseLeg = 0;
1642 for (int32_t i = 0; i < nParts; i++) {
1643 if (trackParts[i]->Leg() == 0 || trackParts[i]->Leg() == leg) {
1644 float zt;
1645 if (Param().par.earlyTpcTransform) {
1646 zt = CAMath::Min(CAMath::Abs(trackParts[i]->ClusterZT0()), CAMath::Abs(trackParts[i]->ClusterZTN()));
1647 } else {
1648 zt = -trackParts[i]->MinClusterZT(); // Negative time ~ smallest z, to behave the same way // TODO: Check all these min / max ZT
1649 }
1650 if (zt < baseZT) {
1651 baseZT = zt;
1652 baseLeg = trackParts[i]->Leg();
1653 }
1654 }
1655 }
1656 int32_t iLongest = 1e9;
1657 int32_t length = 0;
1658 for (int32_t i = (baseLeg ? (nParts - 1) : 0); baseLeg ? (i >= 0) : (i < nParts); baseLeg ? i-- : i++) {
1659 if (trackParts[i]->Leg() != baseLeg) {
1660 break;
1661 }
1662 if (trackParts[i]->OrigTrack()->NHits() > length) {
1663 iLongest = i;
1664 length = trackParts[i]->OrigTrack()->NHits();
1665 }
1666 }
1667 bool outwards;
1668 if (Param().par.earlyTpcTransform) {
1669 outwards = (trackParts[iLongest]->ClusterZT0() > trackParts[iLongest]->ClusterZTN()) ^ trackParts[iLongest]->CSide();
1670 } else {
1671 outwards = trackParts[iLongest]->ClusterZT0() < trackParts[iLongest]->ClusterZTN();
1672 }
1673 GPUTPCGMMerger_CompareClusterIdsLooper::clcomparestruct clusterSort[kMaxClusters];
1674 for (int32_t iPart = 0; iPart < nParts; iPart++) {
1675 const GPUTPCGMSectorTrack* t = trackParts[iPart];
1676 int32_t nTrackHits = t->NClusters();
1677 for (int32_t j = 0; j < nTrackHits; j++) {
1678 int32_t i = nTmpHits + j;
1679 clusterSort[i].leg = t->Leg();
1680 }
1681 nTmpHits += nTrackHits;
1682 }
1683
1684 GPUCommonAlgorithm::sort(clusterIndices, clusterIndices + nHits, GPUTPCGMMerger_CompareClusterIdsLooper(baseLeg, outwards, trackClusters, clusterSort));
1685 } else {
1686 GPUCommonAlgorithm::sort(clusterIndices, clusterIndices + nHits, GPUTPCGMMerger_CompareClusterIds(trackClusters));
1687 }
1688 nTmpHits = 0;
1689 firstTrackIndex = lastTrackIndex = -1;
1690 for (int32_t i = 0; i < nParts; i++) {
1691 nTmpHits += trackParts[i]->NClusters();
1692 if (nTmpHits > clusterIndices[0] && firstTrackIndex == -1) {
1693 firstTrackIndex = i;
1694 }
1695 if (nTmpHits > clusterIndices[nHits - 1] && lastTrackIndex == -1) {
1696 lastTrackIndex = i;
1697 }
1698 }
1699
1700 int32_t nFilteredHits = 0;
1701 int32_t indPrev = -1;
1702 for (int32_t i = 0; i < nHits; i++) {
1703 int32_t ind = clusterIndices[i];
1704 if (indPrev >= 0 && trackClustersUnsorted[ind].id == trackClustersUnsorted[indPrev].id) {
1705 continue;
1706 }
1707 indPrev = ind;
1708 trackClusters[nFilteredHits] = trackClustersUnsorted[ind];
1709 nFilteredHits++;
1710 }
1711 nHits = nFilteredHits;
1712 }
1713
1714 const uint32_t iOutTrackFirstCluster = CAMath::AtomicAdd(&mMemory->nOutputTrackClusters, (uint32_t)nHits);
1715 if (iOutTrackFirstCluster >= mNMaxOutputTrackClusters) {
1716 raiseError(GPUErrors::ERROR_MERGER_HIT_OVERFLOW, iOutTrackFirstCluster, mNMaxOutputTrackClusters);
1717 CAMath::AtomicExch(&mMemory->nOutputTrackClusters, mNMaxOutputTrackClusters);
1718 continue;
1719 }
1720
1721 GPUTPCGMMergedTrackHit* const cl = mClusters + iOutTrackFirstCluster;
1722
1723 for (int32_t i = 0; i < nHits; i++) {
1724 uint8_t state;
1725 if (Param().par.earlyTpcTransform) {
1726 const GPUTPCClusterData& c = GetConstantMem()->tpcTrackers[trackClusters[i].sector].ClusterData()[trackClusters[i].id - GetConstantMem()->tpcTrackers[trackClusters[i].sector].Data().ClusterIdOffset()];
1727 GPUTPCGMMergedTrackHitXYZ* const clXYZ = mClustersXYZ + iOutTrackFirstCluster;
1728 clXYZ[i].x = c.x;
1729 clXYZ[i].y = c.y;
1730 clXYZ[i].z = c.z;
1731 clXYZ[i].amp = c.amp;
1732#ifdef GPUCA_TPC_RAW_PROPAGATE_PAD_ROW_TIME
1733 clXYZ[i].pad = c.mPad;
1734 clXYZ[i].time = c.mTime;
1735#endif
1736 state = c.flags;
1737 } else {
1738 const ClusterNative& c = GetConstantMem()->ioPtrs.clustersNative->clustersLinear[trackClusters[i].id];
1739 state = c.getFlags();
1740 }
1741 cl[i].state = state & GPUTPCGMMergedTrackHit::clustererAndSharedFlags; // Only allow edge, deconvoluted, and shared flags
1742 cl[i].row = trackClusters[i].row;
1743 cl[i].num = trackClusters[i].id;
1744 cl[i].sector = trackClusters[i].sector;
1745 cl[i].leg = trackClusters[i].leg;
1746 }
1747
1748 uint32_t iOutputTrack = CAMath::AtomicAdd(&mMemory->nMergedTracks, 1u);
1749 if (iOutputTrack >= mNMaxTracks) {
1750 raiseError(GPUErrors::ERROR_MERGER_TRACK_OVERFLOW, iOutputTrack, mNMaxTracks);
1751 CAMath::AtomicExch(&mMemory->nMergedTracks, mNMaxTracks);
1752 continue;
1753 }
1754
1755 GPUTPCGMMergedTrack& mergedTrack = mMergedTracks[iOutputTrack];
1756
1757 mergedTrack.SetFlags(0);
1758 mergedTrack.SetOK(1);
1759 mergedTrack.SetLooper(leg > 0);
1760 mergedTrack.SetLegs(leg);
1761 mergedTrack.SetNClusters(nHits);
1762 mergedTrack.SetFirstClusterRef(iOutTrackFirstCluster);
1763 GPUTPCGMTrackParam& p1 = mergedTrack.Param();
1764 const GPUTPCGMSectorTrack& p2 = *trackParts[firstTrackIndex];
1765 mergedTrack.SetCSide(p2.CSide());
1766
1768 const float toX = Param().par.earlyTpcTransform ? mClustersXYZ[iOutTrackFirstCluster].x : GPUTPCGeometry::Row2X(cl[0].row);
1769 if (p2.TransportToX(this, toX, Param().bzCLight, b, GPUCA_MAX_SIN_PHI, false)) {
1770 p1.X() = toX;
1771 p1.Y() = b.Par()[0];
1772 p1.Z() = b.Par()[1];
1773 p1.SinPhi() = b.Par()[2];
1774 } else {
1775 p1.X() = p2.X();
1776 p1.Y() = p2.Y();
1777 p1.Z() = p2.Z();
1778 p1.SinPhi() = p2.SinPhi();
1779 }
1780 p1.TZOffset() = p2.TZOffset();
1781 p1.DzDs() = p2.DzDs();
1782 p1.QPt() = p2.QPt();
1783 mergedTrack.SetAlpha(p2.Alpha());
1784 if (CAMath::Abs(Param().polynomialField.GetNominalBz()) < (gpu_common_constants::kZeroFieldCut * gpu_common_constants::kCLight)) {
1785 p1.QPt() = 100.f / Param().rec.bz0Pt10MeV;
1786 }
1787
1788 // if (nParts > 1) printf("Merged %d: QPt %f %d parts %d hits\n", mMemory->nMergedTracks, p1.QPt(), nParts, nHits);
1789
1790 /*if (GPUQA::QAAvailable() && mRec->GetQA() && mRec->GetQA()->SuppressTrack(mMemory->nMergedTracks))
1791 {
1792 mergedTrack.SetOK(0);
1793 mergedTrack.SetNClusters(0);
1794 }
1795 if (mergedTrack.NClusters() && mergedTrack.OK()) */
1796 if (Param().rec.tpc.mergeCE) {
1797 bool CEside;
1798 if (Param().par.earlyTpcTransform) {
1799 const GPUTPCGMMergedTrackHitXYZ* const clXYZ = mClustersXYZ + iOutTrackFirstCluster;
1800 CEside = (mergedTrack.CSide() != 0) ^ (clXYZ[0].z > clXYZ[nHits - 1].z);
1801 } else {
1802 auto& cls = mConstantMem->ioPtrs.clustersNative->clustersLinear;
1803 CEside = cls[cl[0].num].getTime() < cls[cl[nHits - 1].num].getTime();
1804 }
1805 MergeCEFill(trackParts[CEside ? lastTrackIndex : firstTrackIndex], cl[CEside ? (nHits - 1) : 0], Param().par.earlyTpcTransform ? &(mClustersXYZ + iOutTrackFirstCluster)[CEside ? (nHits - 1) : 0] : nullptr, iOutputTrack);
1806 }
1807 } // itr
1808}
1809
1810GPUd() void GPUTPCGMMerger::SortTracksPrepare(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1811{
1812 for (uint32_t i = iBlock * nThreads + iThread; i < mMemory->nMergedTracks; i += nThreads * nBlocks) {
1813 mTrackOrderProcess[i] = i;
1814 }
1815}
1816
1817GPUd() void GPUTPCGMMerger::PrepareClustersForFit0(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1818{
1819 for (uint32_t i = iBlock * nThreads + iThread; i < mMemory->nMergedTracks; i += nBlocks * nThreads) {
1820 mTrackSort[i] = i;
1821 }
1822}
1823
1824GPUd() void GPUTPCGMMerger::SortTracks(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1825{
1826#ifndef GPUCA_SPECIALIZE_THRUST_SORTS
1827 if (iThread == 0 && iBlock == 0) {
1828 GPUCommonAlgorithm::sortDeviceDynamic(mTrackOrderProcess, mTrackOrderProcess + mMemory->nMergedTracks, GPUTPCGMMergerSortTracks_comp(mMergedTracks));
1829 }
1830#endif
1831}
1832
1833GPUd() void GPUTPCGMMerger::SortTracksQPt(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1834{
1835#ifndef GPUCA_SPECIALIZE_THRUST_SORTS
1836 if (iThread == 0 && iBlock == 0) {
1837 GPUCommonAlgorithm::sortDeviceDynamic(mTrackSort, mTrackSort + mMemory->nMergedTracks, GPUTPCGMMergerSortTracksQPt_comp(mMergedTracks));
1838 }
1839#endif
1840}
1841
1842GPUd() void GPUTPCGMMerger::PrepareClustersForFit1(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1843{
1844 for (uint32_t i = iBlock * nThreads + iThread; i < mMemory->nMergedTracks; i += nBlocks * nThreads) {
1845 mTrackOrderAttach[mTrackSort[i]] = i;
1846 const GPUTPCGMMergedTrack& trk = mMergedTracks[i];
1847 if (trk.OK()) {
1848 for (uint32_t j = 0; j < trk.NClusters(); j++) {
1849 mClusterAttachment[mClusters[trk.FirstClusterRef() + j].num] = attachAttached | attachGood;
1850 CAMath::AtomicAdd(&mSharedCount[mClusters[trk.FirstClusterRef() + j].num], 1u);
1851 }
1852 }
1853 }
1854}
1855
1856GPUd() void GPUTPCGMMerger::PrepareClustersForFit2(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1857{
1858 for (uint32_t i = iBlock * nThreads + iThread; i < mMemory->nOutputTrackClusters; i += nBlocks * nThreads) {
1859 if (mSharedCount[mClusters[i].num] > 1) {
1860 mClusters[i].state |= GPUTPCGMMergedTrackHit::flagShared;
1861 }
1862 }
1863 if (mClusterStateExt) {
1864 for (uint32_t i = iBlock * nThreads + iThread; i < mNMaxClusters; i += nBlocks * nThreads) {
1865 uint8_t state = GetConstantMem()->ioPtrs.clustersNative->clustersLinear[i].getFlags();
1866 if (mSharedCount[i] > 1) {
1868 }
1869 mClusterStateExt[i] = state;
1870 }
1871 }
1872}
1873
1874GPUd() void GPUTPCGMMerger::Finalize0(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1875{
1876 for (uint32_t i = iBlock * nThreads + iThread; i < mMemory->nMergedTracks; i += nThreads * nBlocks) {
1877 mTrackSort[mTrackOrderAttach[i]] = i;
1878 }
1879 for (uint32_t i = iBlock * nThreads + iThread; i < mMemory->nOutputTrackClusters; i += nThreads * nBlocks) {
1880 mClusterAttachment[mClusters[i].num] = 0; // Reset adjacent attachment for attached clusters, set correctly below
1881 }
1882}
1883
1884GPUd() void GPUTPCGMMerger::Finalize1(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1885{
1886 for (uint32_t i = iBlock * nThreads + iThread; i < mMemory->nMergedTracks; i += nThreads * nBlocks) {
1887 const GPUTPCGMMergedTrack& trk = mMergedTracks[i];
1888 if (!trk.OK() || trk.NClusters() == 0) {
1889 continue;
1890 }
1891 uint8_t goodLeg = mClusters[trk.FirstClusterRef() + trk.NClusters() - 1].leg;
1892 for (uint32_t j = 0; j < trk.NClusters(); j++) {
1893 int32_t id = mClusters[trk.FirstClusterRef() + j].num;
1894 uint32_t weight = mTrackOrderAttach[i] | attachAttached;
1895 uint8_t clusterState = mClusters[trk.FirstClusterRef() + j].state;
1896 if (!(clusterState & GPUTPCGMMergedTrackHit::flagReject)) {
1897 weight |= attachGood;
1898 } else if (clusterState & GPUTPCGMMergedTrackHit::flagHighIncl) {
1900 }
1901 if (mClusters[trk.FirstClusterRef() + j].leg == goodLeg) {
1903 }
1904 CAMath::AtomicMax(&mClusterAttachment[id], weight);
1905 }
1906 }
1907}
1908
1909GPUd() void GPUTPCGMMerger::Finalize2(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1910{
1911 for (uint32_t i = iBlock * nThreads + iThread; i < mNMaxClusters; i += nThreads * nBlocks) {
1912 if (mClusterAttachment[i] != 0) {
1913 mClusterAttachment[i] = (mClusterAttachment[i] & attachFlagMask) | mTrackSort[mClusterAttachment[i] & attachTrackMask];
1914 }
1915 }
1916}
1917
1918GPUd() void GPUTPCGMMerger::MergeLoopersInit(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1919{
1920 const float lowPtThresh = Param().rec.tpc.rejectQPtB5 * 1.1f; // Might need to merge tracks above the threshold with parts below the threshold
1921 for (uint32_t i = get_global_id(0); i < mMemory->nMergedTracks; i += get_global_size(0)) {
1922 const auto& trk = mMergedTracks[i];
1923 const auto& p = trk.GetParam();
1924 const float qptabs = CAMath::Abs(p.GetQPt());
1925 if (trk.NClusters() && qptabs * Param().qptB5Scaler > 5.f && qptabs * Param().qptB5Scaler <= lowPtThresh) {
1926 const int32_t sector = mClusters[trk.FirstClusterRef() + trk.NClusters() - 1].sector;
1927 const float refz = p.GetZ() + (Param().par.earlyTpcTransform ? p.GetTZOffset() : GetConstantMem()->calibObjects.fastTransformHelper->getCorrMap()->convVertexTimeToZOffset(sector, p.GetTZOffset(), Param().continuousMaxTimeBin)) + (trk.CSide() ? -100 : 100);
1928 float sinA, cosA;
1929 CAMath::SinCos(trk.GetAlpha(), sinA, cosA);
1930 float gx = cosA * p.GetX() - sinA * p.GetY();
1931 float gy = cosA * p.GetY() + sinA * p.GetX();
1932 float bz = Param().polynomialField.GetFieldBz(gx, gy, p.GetZ());
1933 const float r1 = p.GetQPt() * bz;
1934 const float r = CAMath::Abs(r1) > 0.0001f ? (1.f / r1) : 10000;
1935 const float mx = p.GetX() + r * p.GetSinPhi();
1936 const float my = p.GetY() - r * CAMath::Sqrt(1 - p.GetSinPhi() * p.GetSinPhi());
1937 const float gmx = cosA * mx - sinA * my;
1938 const float gmy = cosA * my + sinA * mx;
1939 uint32_t myId = CAMath::AtomicAdd(&mMemory->nLooperMatchCandidates, 1u);
1940 if (myId >= mNMaxLooperMatches) {
1941 raiseError(GPUErrors::ERROR_LOOPER_MATCH_OVERFLOW, myId, mNMaxLooperMatches);
1942 CAMath::AtomicExch(&mMemory->nLooperMatchCandidates, mNMaxLooperMatches);
1943 return;
1944 }
1945 mLooperCandidates[myId] = MergeLooperParam{refz, gmx, gmy, i};
1946
1947 /*printf("Track %u Sanity qpt %f snp %f bz %f\n", mMemory->nLooperMatchCandidates, p.GetQPt(), p.GetSinPhi(), bz);
1948 for (uint32_t k = 0;k < trk.NClusters();k++) {
1949 float xx, yy, zz;
1950 if (Param().par.earlyTpcTransform) {
1951 const float zOffset = (mClusters[trk.FirstClusterRef() + k].sector < 18) == (mClusters[trk.FirstClusterRef() + 0].sector < 18) ? p.GetTZOffset() : -p.GetTZOffset();
1952 xx = mClustersXYZ[trk.FirstClusterRef() + k].x;
1953 yy = mClustersXYZ[trk.FirstClusterRef() + k].y;
1954 zz = mClustersXYZ[trk.FirstClusterRef() + k].z - zOffset;
1955 } else {
1956 const ClusterNative& GPUrestrict() cl = GetConstantMem()->ioPtrs.clustersNative->clustersLinear[mClusters[trk.FirstClusterRef() + k].num];
1957 GetConstantMem()->calibObjects.fastTransformHelper->Transform(mClusters[trk.FirstClusterRef() + k].sector, mClusters[trk.FirstClusterRef() + k].row, cl.getPad(), cl.getTime(), xx, yy, zz, p.GetTZOffset());
1958 }
1959 float sa2, ca2;
1960 CAMath::SinCos(Param().Alpha(mClusters[trk.FirstClusterRef() + k].sector), sa2, ca2);
1961 float cx = ca2 * xx - sa2 * yy;
1962 float cy = ca2 * yy + sa2 * xx;
1963 float dist = CAMath::Sqrt((cx - gmx) * (cx - gmx) + (cy - gmy) * (cy - gmy));
1964 printf("Hit %3d/%3d sector %d xy %f %f R %f\n", k, trk.NClusters(), (int32_t)mClusters[trk.FirstClusterRef() + k].sector, cx, cy, dist);
1965 }*/
1966 }
1967 }
1968}
1969
1970GPUd() void GPUTPCGMMerger::MergeLoopersSort(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1971{
1972#ifndef GPUCA_SPECIALIZE_THRUST_SORTS
1973 if (iThread == 0 && iBlock == 0) {
1974 GPUCommonAlgorithm::sortDeviceDynamic(mLooperCandidates, mLooperCandidates + mMemory->nLooperMatchCandidates, GPUTPCGMMergerMergeLoopers_comp());
1975 }
1976#endif
1977}
1978
1979GPUd() void GPUTPCGMMerger::MergeLoopersMain(int32_t nBlocks, int32_t nThreads, int32_t iBlock, int32_t iThread)
1980{
1981 const MergeLooperParam* params = mLooperCandidates;
1982
1983#if GPUCA_MERGE_LOOPER_MC && !defined(GPUCA_GPUCODE)
1984 std::vector<int64_t> paramLabels(mMemory->nLooperMatchCandidates);
1985 for (uint32_t i = 0; i < mMemory->nLooperMatchCandidates; i++) {
1986 paramLabels[i] = GetTrackLabel(mMergedTracks[params[i].id]);
1987 }
1988 /*std::vector<bool> dropped(mMemory->nLooperMatchCandidates);
1989 std::vector<bool> droppedMC(mMemory->nLooperMatchCandidates);
1990 std::vector<int32_t> histMatch(101);
1991 std::vector<int32_t> histFail(101);*/
1992 if (!mRec->GetProcessingSettings().runQA) {
1993 throw std::runtime_error("Need QA enabled for the Merge Loopers MC QA");
1994 }
1995#endif
1996
1997 for (uint32_t i = get_global_id(0); i < mMemory->nLooperMatchCandidates; i += get_global_size(0)) {
1998 for (uint32_t j = i + 1; j < mMemory->nLooperMatchCandidates; j++) {
1999 // int32_t bs = 0;
2000 if (CAMath::Abs(params[j].refz) > CAMath::Abs(params[i].refz) + 100.f) {
2001 break;
2002 }
2003 const float d2xy = CAMath::Sum2(params[i].x - params[j].x, params[i].y - params[j].y);
2004 if (d2xy > 15.f) {
2005 // bs |= 1;
2006 continue;
2007 }
2008 const auto& trk1 = mMergedTracks[params[i].id];
2009 const auto& trk2 = mMergedTracks[params[j].id];
2010 const auto& param1 = trk1.GetParam();
2011 const auto& param2 = trk2.GetParam();
2012 if (CAMath::Abs(param1.GetDzDs()) > 0.03f && CAMath::Abs(param2.GetDzDs()) > 0.03f && param1.GetDzDs() * param2.GetDzDs() * param1.GetQPt() * param2.GetQPt() < 0) {
2013 // bs |= 2;
2014 continue;
2015 }
2016
2017 const float dznormalized = (CAMath::Abs(params[j].refz) - CAMath::Abs(params[i].refz)) / (CAMath::TwoPi() * 0.5f * (CAMath::Abs(param1.GetDzDs()) + CAMath::Abs(param2.GetDzDs())) * 1.f / (0.5f * (CAMath::Abs(param1.GetQPt()) + CAMath::Abs(param2.GetQPt())) * CAMath::Abs(Param().polynomialField.GetNominalBz())));
2018 const float phasecorr = CAMath::Modf((CAMath::ASin(param1.GetSinPhi()) + trk1.GetAlpha() - CAMath::ASin(param2.GetSinPhi()) - trk2.GetAlpha()) / CAMath::TwoPi() + 5.5f, 1.f) - 0.5f;
2019 const float phasecorrdirection = (params[j].refz * param1.GetQPt() * param1.GetDzDs()) > 0 ? 1 : -1;
2020 const float dzcorr = dznormalized + phasecorr * phasecorrdirection;
2021 const bool sameside = !(trk1.CSide() ^ trk2.CSide());
2022 const float dzcorrlimit[4] = {sameside ? 0.018f : 0.012f, sameside ? 0.12f : 0.025f, 0.14f, 0.15f};
2023 const int32_t dzcorrcount = sameside ? 4 : 2;
2024 bool dzcorrok = false;
2025 float dznorm = 0.f;
2026 for (int32_t k = 0; k < dzcorrcount; k++) {
2027 const float d = CAMath::Abs(dzcorr - 0.5f * k);
2028 if (d <= dzcorrlimit[k]) {
2029 dzcorrok = true;
2030 dznorm = d / dzcorrlimit[k];
2031 break;
2032 }
2033 }
2034 if (!dzcorrok) {
2035 // bs |= 4;
2036 continue;
2037 }
2038
2039 const float dtgl = param1.GetDzDs() - (param1.GetQPt() * param2.GetQPt() > 0 ? param2.GetDzDs() : -param2.GetDzDs());
2040 const float dqpt = (CAMath::Abs(param1.GetQPt()) - CAMath::Abs(param2.GetQPt())) / CAMath::Min(param1.GetQPt(), param2.GetQPt());
2041 float d = CAMath::Sum2(dtgl * (1.f / 0.03f), dqpt * (1.f / 0.04f)) + d2xy * (1.f / 4.f) + dznorm * (1.f / 0.3f);
2042 bool EQ = d < 6.f;
2043#if GPUCA_MERGE_LOOPER_MC && !defined(GPUCA_GPUCODE)
2044 const int64_t label1 = paramLabels[i];
2045 const int64_t label2 = paramLabels[j];
2046 bool labelEQ = label1 != -1 && label1 == label2;
2047 if (1 || EQ || labelEQ) {
2048 // printf("Matching track %d/%d %u-%u (%ld/%ld): dist %f side %d %d, tgl %f %f, qpt %f %f, x %f %f, y %f %f\n", (int32_t)EQ, (int32_t)labelEQ, i, j, label1, label2, d, (int32_t)mMergedTracks[params[i].id].CSide(), (int32_t)mMergedTracks[params[j].id].CSide(), params[i].tgl, params[j].tgl, params[i].qpt, params[j].qpt, params[i].x, params[j].x, params[i].y, params[j].y);
2049 static auto& tup = GPUROOTDump<TNtuple>::get("mergeloopers", "labeleq:sides:d2xy:tgl1:tgl2:qpt1:qpt2:dz:dzcorr:dtgl:dqpt:dznorm:bs");
2050 tup.Fill((float)labelEQ, (trk1.CSide() ? 1 : 0) | (trk2.CSide() ? 2 : 0), d2xy, param1.GetDzDs(), param2.GetDzDs(), param1.GetQPt(), param2.GetQPt(), CAMath::Abs(params[j].refz) - CAMath::Abs(params[i].refz), dzcorr, dtgl, dqpt, dznorm, bs);
2051 static auto tup2 = GPUROOTDump<TNtuple>::getNew("mergeloopers2", "labeleq:refz1:refz2:tgl1:tgl2:qpt1:qpt2:snp1:snp2:a1:a2:dzn:phasecor:phasedir:dzcorr");
2052 tup2.Fill((float)labelEQ, params[i].refz, params[j].refz, param1.GetDzDs(), param2.GetDzDs(), param1.GetQPt(), param2.GetQPt(), param1.GetSinPhi(), param2.GetSinPhi(), trk1.GetAlpha(), trk2.GetAlpha(), dznormalized, phasecorr, phasecorrdirection, dzcorr);
2053 }
2054 /*if (EQ) {
2055 dropped[j] = true;
2056 }
2057 if (labelEQ) {
2058 droppedMC[j] = true;
2059 histMatch[CAMath::Min<int32_t>(100, d * 10.f)]++;
2060 }
2061 if (d < 10.f && !labelEQ) {
2062 histFail[CAMath::Min<int32_t>(100, d * 10.f)]++;
2063 }*/
2064#endif
2065 if (EQ) {
2066 mMergedTracks[params[j].id].SetMergedLooper(true);
2067 if (CAMath::Abs(param2.GetQPt() * Param().qptB5Scaler) >= Param().rec.tpc.rejectQPtB5) {
2068 mMergedTracks[params[i].id].SetMergedLooper(true);
2069 }
2070 }
2071 }
2072 }
2073 /*#if GPUCA_MERGE_LOOPER_MC && !defined(GPUCA_GPUCODE)
2074 int32_t total = 0, totalmc = 0, good = 0, missed = 0, fake = 0;
2075 for (uint32_t i = 0; i < mMemory->nLooperMatchCandidates; i++) {
2076 total += dropped[i];
2077 totalmc += droppedMC[i];
2078 good += dropped[i] && droppedMC[i];
2079 missed += droppedMC[i] && !dropped[i];
2080 fake += dropped[i] && !droppedMC[i];
2081 }
2082 if (good) {
2083 printf("%20s: %8d\n", "Total", total);
2084 printf("%20s: %8d\n", "TotalMC", totalmc);
2085 printf("%20s: %8d (%8.3f%% %8.3f%%)\n", "Good", good, 100.f * good / total, 100.f * good / totalmc);
2086 printf("%20s: %8d (%8.3f%%)\n", "Missed", missed, 100.f * missed / totalmc);
2087 printf("%20s: %8d (%8.3f%%)\n", "Fake", fake, 100.f * fake / total);
2088 }
2089 printf("Match histo\n");
2090 for (uint32_t i = 0; i < histMatch.size(); i++) {
2091 printf("%8.3f: %3d\n", i / 10.f + 0.05f, histMatch[i]);
2092 }
2093 printf("Fake histo\n");
2094 for (uint32_t i = 0; i < histFail.size(); i++) {
2095 printf("%8.3f: %3d\n", i / 10.f + 0.05f, histFail[i]);
2096 }
2097#endif*/
2098}
benchmark::State & state
Class of a TPC cluster in TPC-native coordinates (row, time)
A const (ready only) version of MCTruthContainer.
int32_t i
#define GPUdii()
#define get_global_size(dim)
#define GPUAtomic(type)
#define GPUbarrier()
#define GPUdni()
#define get_global_id(dim)
#define GPUCA_DEBUG_STREAMER_CHECK(...)
#define GPUCA_DETERMINISTIC_CODE(det, indet)
#define GPUCA_MERGER_MAX_TRACK_CLUSTERS
#define GPUCA_MAX_SIN_PHI_LOW
#define GPUCA_TRACKLET_SELECTOR_MIN_HITS_B5(QPTB5)
#define GPUCA_MAX_SIN_PHI
#define CADEBUG(...)
Definition GPUDef.h:79
#define CADEBUG2(cmd,...)
Definition GPUDef.h:80
bool o
uint8_t leg
const GPUTPCGMMerger::trackCluster & b1
const GPUTPCGMMerger::trackCluster & a1
const GPUTPCGMMerger::trackCluster *const mCmp
const bool outwards
const clcomparestruct *const cmp2
const GPUTPCGMMerger::trackCluster *const cmp1
bool const GPUTPCGMMerger::trackCluster * c1
bool const GPUTPCGMMerger::trackCluster const clcomparestruct * c2
#define GPUCA_NSECTORS
#define GPUCA_ROW_COUNT
constexpr int p2()
constexpr int p1()
constexpr to accelerate the coordinates changing
void output(const std::map< std::string, ChannelStat > &channels)
Definition rawdump.cxx:197
uint16_t pos
Definition RawData.h:3
uint32_t j
Definition RawData.h:0
uint32_t c
Definition RawData.h:2
Definition of TPCFastTransform class.
std::array< int, 64 > reverse(std::array< int, 64 > a)
int nClusters
double num
Definition B.h:16
Class for time synchronization of RawReader instances.
GPUReconstruction * mRec
static void computePointerWithAlignment(T *&basePtr, S *&objPtr, size_t nEntries=1)
static GPUROOTDump< T, Args... > getNew(const char *name1, Names... names)
Definition GPUROOTDump.h:64
static GPUROOTDump< T, Args... > & get(const char *name1, Names... names)
Definition GPUROOTDump.h:58
RecoStepField GetRecoStepsGPU() const
int16_t RegisterMemoryAllocation(T *proc, void *(T::*setPtr)(void *), int32_t type, const char *name="", const GPUMemoryReuse &re=GPUMemoryReuse())
RecoStepField GetRecoSteps() const
const GPUParam & GetParam() const
const GPUConstantMem & GetConstantMem() const
GPUMemorySizeScalers * MemoryScalers()
const GPUSettingsProcessing & GetProcessingSettings() const
virtual const GPUDefParameters & getGPUParameters(bool doGPU) const =0
void * SetPointersRefitScratch(void *mem)
void * SetPointersMemory(void *mem)
static constexpr const int32_t NSECTORS
void * SetPointersOutput(void *mem)
void * SetPointersOutputO2MC(void *mem)
void * SetPointersOutputO2Scratch(void *mem)
void * SetPointersMerger(void *mem)
void * SetPointersOutputO2Clus(void *mem)
void * SetPointersOutputState(void *mem)
void SetMaxData(const GPUTrackingInOutPointers &io)
void * SetPointersOutputO2(void *mem)
int32_t int32_t int32_t bool float float Z
int32_t int32_t int32_t bool float Y
GLint GLint GLsizei GLint border
Definition glcorearb.h:275
GLdouble n
Definition glcorearb.h:1982
GLfloat GLfloat GLfloat alpha
Definition glcorearb.h:279
GLint GLenum GLint x
Definition glcorearb.h:403
GLuint GLuint end
Definition glcorearb.h:469
const GLdouble * v
Definition glcorearb.h:832
GLuint GLuint GLfloat weight
Definition glcorearb.h:5477
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLenum GLint * range
Definition glcorearb.h:1899
GLint y
Definition glcorearb.h:270
GLenum const GLfloat * params
Definition glcorearb.h:272
GLintptr offset
Definition glcorearb.h:660
GLuint GLsizei GLsizei * length
Definition glcorearb.h:790
typedef void(APIENTRYP PFNGLCULLFACEPROC)(GLenum mode)
GLuint GLfloat x0
Definition glcorearb.h:5034
GLbitfield flags
Definition glcorearb.h:1570
GLboolean r
Definition glcorearb.h:1233
GLuint start
Definition glcorearb.h:469
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLuint id
Definition glcorearb.h:650
GLdouble GLdouble GLdouble z
Definition glcorearb.h:843
uint8_t itsSharedClusterMap uint8_t
Node par(int index)
Parameters.
GPUdi() o2
Definition TrackTRD.h:38
GPUd() const expr uint32_t MultivariatePolynomialHelper< Dim
float float float float z1
Definition MathUtils.h:44
if(!okForPhiMin(phi0, phi1))
Global TPC definitions and constants.
Definition SimTraits.h:167
@ streamMergeBorderTracksBest
stream MergeBorderTracks best track
@ streamMergeBorderTracksAll
stream MergeBorderTracks all tracks
std::string getTime(uint64_t ts)
GPUReconstruction * rec
GPUTPCTracker tpcTrackers[GPUCA_NSECTORS]
size_t getValue(size_t maxVal, size_t val)
size_t NTPCMergedTrackHits(size_t tpcSectorTrackHitss)
size_t NTPCMergedTracks(size_t tpcSectorTracks)
const o2::tpc::ClusterNativeAccess * clustersNative
const GPUTPCGMMergedTrack & GPUrestrict() b
GPUd() bool operator()(const int32_t aa
GPUhd() GPUTPCGMMergerSortTracksQPt_comp(GPUTPCGMMergedTrack *cmp)
const GPUTPCGMMergedTrack & GPUrestrict() b
GPUhd() GPUTPCGMMergerSortTracks_comp(GPUTPCGMMergedTrack *cmp)
GPUd() bool operator()(const int32_t aa
std::vector< int > row
char const *restrict const cmp
Definition x9.h:96