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