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