Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-12-22 23:30:01

0001 #ifndef RecoTracker_LSTCore_src_alpaka_TrackCandidate_h
0002 #define RecoTracker_LSTCore_src_alpaka_TrackCandidate_h
0003 
0004 #include "RecoTracker/LSTCore/interface/alpaka/Common.h"
0005 #include "RecoTracker/LSTCore/interface/ModulesSoA.h"
0006 #include "RecoTracker/LSTCore/interface/MiniDoubletsSoA.h"
0007 #include "RecoTracker/LSTCore/interface/PixelQuintupletsSoA.h"
0008 #include "RecoTracker/LSTCore/interface/PixelTripletsSoA.h"
0009 #include "RecoTracker/LSTCore/interface/QuintupletsSoA.h"
0010 #include "RecoTracker/LSTCore/interface/SegmentsSoA.h"
0011 #include "RecoTracker/LSTCore/interface/TrackCandidatesSoA.h"
0012 #include "RecoTracker/LSTCore/interface/TripletsSoA.h"
0013 
0014 #include "Hit.h"
0015 
0016 namespace ALPAKA_ACCELERATOR_NAMESPACE::lst {
0017   ALPAKA_FN_ACC ALPAKA_FN_INLINE void addpLSTrackCandidateToMemory(TrackCandidates& cands,
0018                                                                    unsigned int trackletIndex,
0019                                                                    unsigned int trackCandidateIndex,
0020                                                                    uint4 hitIndices,
0021                                                                    int pixelSeedIndex) {
0022     cands.trackCandidateType()[trackCandidateIndex] = LSTObjType::pLS;
0023     cands.directObjectIndices()[trackCandidateIndex] = trackletIndex;
0024     cands.pixelSeedIndex()[trackCandidateIndex] = pixelSeedIndex;
0025 
0026     cands.objectIndices()[trackCandidateIndex][0] = trackletIndex;
0027     cands.objectIndices()[trackCandidateIndex][1] = trackletIndex;
0028 
0029     cands.hitIndices()[trackCandidateIndex][0] =
0030         hitIndices.x;  // Order explanation in https://github.com/SegmentLinking/TrackLooper/issues/267
0031     cands.hitIndices()[trackCandidateIndex][1] = hitIndices.z;
0032     cands.hitIndices()[trackCandidateIndex][2] = hitIndices.y;
0033     cands.hitIndices()[trackCandidateIndex][3] = hitIndices.w;
0034   }
0035 
0036   ALPAKA_FN_ACC ALPAKA_FN_INLINE void addTrackCandidateToMemory(TrackCandidates& cands,
0037                                                                 short trackCandidateType,
0038                                                                 unsigned int innerTrackletIndex,
0039                                                                 unsigned int outerTrackletIndex,
0040                                                                 const uint8_t* logicalLayerIndices,
0041                                                                 const uint16_t* lowerModuleIndices,
0042                                                                 const unsigned int* hitIndices,
0043                                                                 int pixelSeedIndex,
0044                                                                 float centerX,
0045                                                                 float centerY,
0046                                                                 float radius,
0047                                                                 unsigned int trackCandidateIndex,
0048                                                                 unsigned int directObjectIndex) {
0049     cands.trackCandidateType()[trackCandidateIndex] = trackCandidateType;
0050     cands.directObjectIndices()[trackCandidateIndex] = directObjectIndex;
0051     cands.pixelSeedIndex()[trackCandidateIndex] = pixelSeedIndex;
0052 
0053     cands.objectIndices()[trackCandidateIndex][0] = innerTrackletIndex;
0054     cands.objectIndices()[trackCandidateIndex][1] = outerTrackletIndex;
0055 
0056     size_t limits = trackCandidateType == LSTObjType::pT5 ? Params_pT5::kLayers : Params_pT3::kLayers;
0057 
0058     //send the starting pointer to the logicalLayer and hitIndices
0059     for (size_t i = 0; i < limits; i++) {
0060       cands.logicalLayers()[trackCandidateIndex][i] = logicalLayerIndices[i];
0061       cands.lowerModuleIndices()[trackCandidateIndex][i] = lowerModuleIndices[i];
0062     }
0063     for (size_t i = 0; i < 2 * limits; i++) {
0064       cands.hitIndices()[trackCandidateIndex][i] = hitIndices[i];
0065     }
0066     cands.centerX()[trackCandidateIndex] = __F2H(centerX);
0067     cands.centerY()[trackCandidateIndex] = __F2H(centerY);
0068     cands.radius()[trackCandidateIndex] = __F2H(radius);
0069   }
0070 
0071   ALPAKA_FN_ACC ALPAKA_FN_INLINE int checkPixelHits(
0072       unsigned int ix, unsigned int jx, MiniDoubletsConst mds, SegmentsConst segments, HitsConst hits) {
0073     int phits1[Params_pLS::kHits];
0074     int phits2[Params_pLS::kHits];
0075 
0076     phits1[0] = hits.idxs()[mds.anchorHitIndices()[segments.mdIndices()[ix][0]]];
0077     phits1[1] = hits.idxs()[mds.anchorHitIndices()[segments.mdIndices()[ix][1]]];
0078     phits1[2] = hits.idxs()[mds.outerHitIndices()[segments.mdIndices()[ix][0]]];
0079     phits1[3] = hits.idxs()[mds.outerHitIndices()[segments.mdIndices()[ix][1]]];
0080 
0081     phits2[0] = hits.idxs()[mds.anchorHitIndices()[segments.mdIndices()[jx][0]]];
0082     phits2[1] = hits.idxs()[mds.anchorHitIndices()[segments.mdIndices()[jx][1]]];
0083     phits2[2] = hits.idxs()[mds.outerHitIndices()[segments.mdIndices()[jx][0]]];
0084     phits2[3] = hits.idxs()[mds.outerHitIndices()[segments.mdIndices()[jx][1]]];
0085 
0086     int npMatched = 0;
0087 
0088     for (int i = 0; i < Params_pLS::kHits; i++) {
0089       bool pmatched = false;
0090       if (phits1[i] == -1)
0091         continue;
0092 
0093       for (int j = 0; j < Params_pLS::kHits; j++) {
0094         if (phits2[j] == -1)
0095           continue;
0096 
0097         if (phits1[i] == phits2[j]) {
0098           pmatched = true;
0099           break;
0100         }
0101       }
0102       if (pmatched)
0103         npMatched++;
0104     }
0105     return npMatched;
0106   }
0107 
0108   struct CrossCleanpT3 {
0109     template <typename TAcc>
0110     ALPAKA_FN_ACC void operator()(TAcc const& acc,
0111                                   ModulesConst modules,
0112                                   ObjectRangesConst ranges,
0113                                   PixelTriplets pixelTriplets,
0114                                   SegmentsPixelConst segmentsPixel,
0115                                   PixelQuintupletsConst pixelQuintuplets) const {
0116       auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
0117       auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
0118 
0119       unsigned int nPixelTriplets = pixelTriplets.nPixelTriplets();
0120       for (unsigned int pixelTripletIndex = globalThreadIdx[2]; pixelTripletIndex < nPixelTriplets;
0121            pixelTripletIndex += gridThreadExtent[2]) {
0122         if (pixelTriplets.isDup()[pixelTripletIndex])
0123           continue;
0124 
0125         // Cross cleaning step
0126         float eta1 = __H2F(pixelTriplets.eta_pix()[pixelTripletIndex]);
0127         float phi1 = __H2F(pixelTriplets.phi_pix()[pixelTripletIndex]);
0128 
0129         int pixelModuleIndex = modules.nLowerModules();
0130         unsigned int prefix = ranges.segmentModuleIndices()[pixelModuleIndex];
0131 
0132         unsigned int nPixelQuintuplets = pixelQuintuplets.nPixelQuintuplets();
0133         for (unsigned int pixelQuintupletIndex = globalThreadIdx[1]; pixelQuintupletIndex < nPixelQuintuplets;
0134              pixelQuintupletIndex += gridThreadExtent[1]) {
0135           unsigned int pLS_jx = pixelQuintuplets.pixelSegmentIndices()[pixelQuintupletIndex];
0136           float eta2 = segmentsPixel.eta()[pLS_jx - prefix];
0137           float phi2 = segmentsPixel.phi()[pLS_jx - prefix];
0138           float dEta = alpaka::math::abs(acc, (eta1 - eta2));
0139           float dPhi = calculate_dPhi(phi1, phi2);
0140 
0141           float dR2 = dEta * dEta + dPhi * dPhi;
0142           if (dR2 < 1e-5f)
0143             pixelTriplets.isDup()[pixelTripletIndex] = true;
0144         }
0145       }
0146     }
0147   };
0148 
0149   struct CrossCleanT5 {
0150     template <typename TAcc>
0151     ALPAKA_FN_ACC void operator()(TAcc const& acc,
0152                                   ModulesConst modules,
0153                                   Quintuplets quintuplets,
0154                                   QuintupletsOccupancyConst quintupletsOccupancy,
0155                                   PixelQuintupletsConst pixelQuintuplets,
0156                                   PixelTripletsConst pixelTriplets,
0157                                   ObjectRangesConst ranges) const {
0158       auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
0159       auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
0160 
0161       for (int innerInnerInnerLowerModuleArrayIndex = globalThreadIdx[0];
0162            innerInnerInnerLowerModuleArrayIndex < modules.nLowerModules();
0163            innerInnerInnerLowerModuleArrayIndex += gridThreadExtent[0]) {
0164         if (ranges.quintupletModuleIndices()[innerInnerInnerLowerModuleArrayIndex] == -1)
0165           continue;
0166 
0167         unsigned int nQuints = quintupletsOccupancy.nQuintuplets()[innerInnerInnerLowerModuleArrayIndex];
0168         for (unsigned int innerObjectArrayIndex = globalThreadIdx[1]; innerObjectArrayIndex < nQuints;
0169              innerObjectArrayIndex += gridThreadExtent[1]) {
0170           unsigned int quintupletIndex =
0171               ranges.quintupletModuleIndices()[innerInnerInnerLowerModuleArrayIndex] + innerObjectArrayIndex;
0172 
0173           // Don't add duplicate T5s or T5s that are accounted in pT5s
0174           if (quintuplets.isDup()[quintupletIndex] or quintuplets.partOfPT5()[quintupletIndex])
0175             continue;
0176           unsigned int loop_bound = pixelQuintuplets.nPixelQuintuplets() + pixelTriplets.nPixelTriplets();
0177           // Cross cleaning step
0178           float eta1 = __H2F(quintuplets.eta()[quintupletIndex]);
0179           float phi1 = __H2F(quintuplets.phi()[quintupletIndex]);
0180 
0181           for (unsigned int jx = globalThreadIdx[2]; jx < loop_bound; jx += gridThreadExtent[2]) {
0182             float eta2, phi2;
0183             if (jx < pixelQuintuplets.nPixelQuintuplets()) {
0184               eta2 = __H2F(pixelQuintuplets.eta()[jx]);
0185               phi2 = __H2F(pixelQuintuplets.phi()[jx]);
0186             } else {
0187               eta2 = __H2F(pixelTriplets.eta()[jx - pixelQuintuplets.nPixelQuintuplets()]);
0188               phi2 = __H2F(pixelTriplets.phi()[jx - pixelQuintuplets.nPixelQuintuplets()]);
0189             }
0190 
0191             float dEta = alpaka::math::abs(acc, eta1 - eta2);
0192             float dPhi = calculate_dPhi(phi1, phi2);
0193 
0194             float dR2 = dEta * dEta + dPhi * dPhi;
0195             if (dR2 < 1e-3f)
0196               quintuplets.isDup()[quintupletIndex] = true;
0197           }
0198         }
0199       }
0200     }
0201   };
0202 
0203   struct CrossCleanpLS {
0204     template <typename TAcc>
0205     ALPAKA_FN_ACC void operator()(TAcc const& acc,
0206                                   ModulesConst modules,
0207                                   ObjectRangesConst ranges,
0208                                   PixelTripletsConst pixelTriplets,
0209                                   TrackCandidates cands,
0210                                   SegmentsConst segments,
0211                                   SegmentsOccupancyConst segmentsOccupancy,
0212                                   SegmentsPixel segmentsPixel,
0213                                   MiniDoubletsConst mds,
0214                                   HitsConst hits,
0215                                   QuintupletsConst quintuplets) const {
0216       auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
0217       auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
0218 
0219       int pixelModuleIndex = modules.nLowerModules();
0220       unsigned int nPixels = segmentsOccupancy.nSegments()[pixelModuleIndex];
0221       for (unsigned int pixelArrayIndex = globalThreadIdx[2]; pixelArrayIndex < nPixels;
0222            pixelArrayIndex += gridThreadExtent[2]) {
0223         if (!segmentsPixel.isQuad()[pixelArrayIndex] || segmentsPixel.isDup()[pixelArrayIndex])
0224           continue;
0225 
0226         float eta1 = segmentsPixel.eta()[pixelArrayIndex];
0227         float phi1 = segmentsPixel.phi()[pixelArrayIndex];
0228         unsigned int prefix = ranges.segmentModuleIndices()[pixelModuleIndex];
0229 
0230         unsigned int nTrackCandidates = cands.nTrackCandidates();
0231         for (unsigned int trackCandidateIndex = globalThreadIdx[1]; trackCandidateIndex < nTrackCandidates;
0232              trackCandidateIndex += gridThreadExtent[1]) {
0233           short type = cands.trackCandidateType()[trackCandidateIndex];
0234           unsigned int innerTrackletIdx = cands.objectIndices()[trackCandidateIndex][0];
0235           if (type == LSTObjType::T5) {
0236             unsigned int quintupletIndex = innerTrackletIdx;  // T5 index
0237             float eta2 = __H2F(quintuplets.eta()[quintupletIndex]);
0238             float phi2 = __H2F(quintuplets.phi()[quintupletIndex]);
0239             float dEta = alpaka::math::abs(acc, eta1 - eta2);
0240             float dPhi = calculate_dPhi(phi1, phi2);
0241 
0242             float dR2 = dEta * dEta + dPhi * dPhi;
0243             if (dR2 < 1e-3f)
0244               segmentsPixel.isDup()[pixelArrayIndex] = true;
0245           }
0246           if (type == LSTObjType::pT3) {
0247             int pLSIndex = pixelTriplets.pixelSegmentIndices()[innerTrackletIdx];
0248             int npMatched = checkPixelHits(prefix + pixelArrayIndex, pLSIndex, mds, segments, hits);
0249             if (npMatched > 0)
0250               segmentsPixel.isDup()[pixelArrayIndex] = true;
0251 
0252             int pT3Index = innerTrackletIdx;
0253             float eta2 = __H2F(pixelTriplets.eta_pix()[pT3Index]);
0254             float phi2 = __H2F(pixelTriplets.phi_pix()[pT3Index]);
0255             float dEta = alpaka::math::abs(acc, eta1 - eta2);
0256             float dPhi = calculate_dPhi(phi1, phi2);
0257 
0258             float dR2 = dEta * dEta + dPhi * dPhi;
0259             if (dR2 < 0.000001f)
0260               segmentsPixel.isDup()[pixelArrayIndex] = true;
0261           }
0262           if (type == LSTObjType::pT5) {
0263             unsigned int pLSIndex = innerTrackletIdx;
0264             int npMatched = checkPixelHits(prefix + pixelArrayIndex, pLSIndex, mds, segments, hits);
0265             if (npMatched > 0) {
0266               segmentsPixel.isDup()[pixelArrayIndex] = true;
0267             }
0268 
0269             float eta2 = segmentsPixel.eta()[pLSIndex - prefix];
0270             float phi2 = segmentsPixel.phi()[pLSIndex - prefix];
0271             float dEta = alpaka::math::abs(acc, eta1 - eta2);
0272             float dPhi = calculate_dPhi(phi1, phi2);
0273 
0274             float dR2 = dEta * dEta + dPhi * dPhi;
0275             if (dR2 < 0.000001f)
0276               segmentsPixel.isDup()[pixelArrayIndex] = true;
0277           }
0278         }
0279       }
0280     }
0281   };
0282 
0283   struct AddpT3asTrackCandidates {
0284     template <typename TAcc>
0285     ALPAKA_FN_ACC void operator()(TAcc const& acc,
0286                                   uint16_t nLowerModules,
0287                                   PixelTripletsConst pixelTriplets,
0288                                   TrackCandidates cands,
0289                                   SegmentsPixelConst segmentsPixel,
0290                                   ObjectRangesConst ranges) const {
0291       // implementation is 1D with a single block
0292       static_assert(std::is_same_v<TAcc, ALPAKA_ACCELERATOR_NAMESPACE::Acc1D>, "Should be Acc1D");
0293       ALPAKA_ASSERT_ACC((alpaka::getWorkDiv<alpaka::Grid, alpaka::Blocks>(acc)[0] == 1));
0294 
0295       auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
0296       auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
0297 
0298       unsigned int nPixelTriplets = pixelTriplets.nPixelTriplets();
0299       unsigned int pLS_offset = ranges.segmentModuleIndices()[nLowerModules];
0300       for (unsigned int pixelTripletIndex = globalThreadIdx[0]; pixelTripletIndex < nPixelTriplets;
0301            pixelTripletIndex += gridThreadExtent[0]) {
0302         if ((pixelTriplets.isDup()[pixelTripletIndex]))
0303           continue;
0304 
0305         unsigned int trackCandidateIdx =
0306             alpaka::atomicAdd(acc, &cands.nTrackCandidates(), 1u, alpaka::hierarchy::Threads{});
0307         if (trackCandidateIdx >= n_max_pixel_track_candidates)  // This is done before any non-pixel TCs are added
0308         {
0309 #ifdef WARNINGS
0310           printf("Track Candidate excess alert! Type = pT3");
0311 #endif
0312           alpaka::atomicSub(acc, &cands.nTrackCandidates(), 1u, alpaka::hierarchy::Threads{});
0313           break;
0314 
0315         } else {
0316           alpaka::atomicAdd(acc, &cands.nTrackCandidatespT3(), 1u, alpaka::hierarchy::Threads{});
0317 
0318           float radius = 0.5f * (__H2F(pixelTriplets.pixelRadius()[pixelTripletIndex]) +
0319                                  __H2F(pixelTriplets.tripletRadius()[pixelTripletIndex]));
0320           unsigned int pT3PixelIndex = pixelTriplets.pixelSegmentIndices()[pixelTripletIndex];
0321           addTrackCandidateToMemory(cands,
0322                                     LSTObjType::pT3,
0323                                     pixelTripletIndex,
0324                                     pixelTripletIndex,
0325                                     pixelTriplets.logicalLayers()[pixelTripletIndex].data(),
0326                                     pixelTriplets.lowerModuleIndices()[pixelTripletIndex].data(),
0327                                     pixelTriplets.hitIndices()[pixelTripletIndex].data(),
0328                                     segmentsPixel.seedIdx()[pT3PixelIndex - pLS_offset],
0329                                     __H2F(pixelTriplets.centerX()[pixelTripletIndex]),
0330                                     __H2F(pixelTriplets.centerY()[pixelTripletIndex]),
0331                                     radius,
0332                                     trackCandidateIdx,
0333                                     pixelTripletIndex);
0334         }
0335       }
0336     }
0337   };
0338 
0339   struct AddT5asTrackCandidate {
0340     template <typename TAcc>
0341     ALPAKA_FN_ACC void operator()(TAcc const& acc,
0342                                   uint16_t nLowerModules,
0343                                   QuintupletsConst quintuplets,
0344                                   QuintupletsOccupancyConst quintupletsOccupancy,
0345                                   TrackCandidates cands,
0346                                   ObjectRangesConst ranges) const {
0347       auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
0348       auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
0349 
0350       for (int idx = globalThreadIdx[1]; idx < nLowerModules; idx += gridThreadExtent[1]) {
0351         if (ranges.quintupletModuleIndices()[idx] == -1)
0352           continue;
0353 
0354         unsigned int nQuints = quintupletsOccupancy.nQuintuplets()[idx];
0355         for (unsigned int jdx = globalThreadIdx[2]; jdx < nQuints; jdx += gridThreadExtent[2]) {
0356           unsigned int quintupletIndex = ranges.quintupletModuleIndices()[idx] + jdx;
0357           if (quintuplets.isDup()[quintupletIndex] or quintuplets.partOfPT5()[quintupletIndex])
0358             continue;
0359           if (!(quintuplets.tightCutFlag()[quintupletIndex]))
0360             continue;
0361 
0362           unsigned int trackCandidateIdx =
0363               alpaka::atomicAdd(acc, &cands.nTrackCandidates(), 1u, alpaka::hierarchy::Threads{});
0364           if (trackCandidateIdx - cands.nTrackCandidatespT5() - cands.nTrackCandidatespT3() >=
0365               n_max_nonpixel_track_candidates)  // pT5 and pT3 TCs have been added, but not pLS TCs
0366           {
0367 #ifdef WARNINGS
0368             printf("Track Candidate excess alert! Type = T5");
0369 #endif
0370             alpaka::atomicSub(acc, &cands.nTrackCandidates(), 1u, alpaka::hierarchy::Threads{});
0371             break;
0372           } else {
0373             alpaka::atomicAdd(acc, &cands.nTrackCandidatesT5(), 1u, alpaka::hierarchy::Threads{});
0374             addTrackCandidateToMemory(cands,
0375                                       LSTObjType::T5,
0376                                       quintupletIndex,
0377                                       quintupletIndex,
0378                                       quintuplets.logicalLayers()[quintupletIndex].data(),
0379                                       quintuplets.lowerModuleIndices()[quintupletIndex].data(),
0380                                       quintuplets.hitIndices()[quintupletIndex].data(),
0381                                       -1 /*no pixel seed index for T5s*/,
0382                                       quintuplets.regressionCenterX()[quintupletIndex],
0383                                       quintuplets.regressionCenterY()[quintupletIndex],
0384                                       quintuplets.regressionRadius()[quintupletIndex],
0385                                       trackCandidateIdx,
0386                                       quintupletIndex);
0387           }
0388         }
0389       }
0390     }
0391   };
0392 
0393   struct AddpLSasTrackCandidate {
0394     template <typename TAcc>
0395     ALPAKA_FN_ACC void operator()(TAcc const& acc,
0396                                   uint16_t nLowerModules,
0397                                   TrackCandidates cands,
0398                                   SegmentsOccupancyConst segmentsOccupancy,
0399                                   SegmentsPixelConst segmentsPixel,
0400                                   bool tc_pls_triplets) const {
0401       auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
0402       auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
0403 
0404       unsigned int nPixels = segmentsOccupancy.nSegments()[nLowerModules];
0405       for (unsigned int pixelArrayIndex = globalThreadIdx[2]; pixelArrayIndex < nPixels;
0406            pixelArrayIndex += gridThreadExtent[2]) {
0407         if ((tc_pls_triplets ? 0 : !segmentsPixel.isQuad()[pixelArrayIndex]) ||
0408             (segmentsPixel.isDup()[pixelArrayIndex]))
0409           continue;
0410 
0411         unsigned int trackCandidateIdx =
0412             alpaka::atomicAdd(acc, &cands.nTrackCandidates(), 1u, alpaka::hierarchy::Threads{});
0413         if (trackCandidateIdx - cands.nTrackCandidatesT5() >=
0414             n_max_pixel_track_candidates)  // T5 TCs have already been added
0415         {
0416 #ifdef WARNINGS
0417           printf("Track Candidate excess alert! Type = pLS");
0418 #endif
0419           alpaka::atomicSub(acc, &cands.nTrackCandidates(), 1u, alpaka::hierarchy::Threads{});
0420           break;
0421 
0422         } else {
0423           alpaka::atomicAdd(acc, &cands.nTrackCandidatespLS(), 1u, alpaka::hierarchy::Threads{});
0424           addpLSTrackCandidateToMemory(cands,
0425                                        pixelArrayIndex,
0426                                        trackCandidateIdx,
0427                                        segmentsPixel.pLSHitsIdxs()[pixelArrayIndex],
0428                                        segmentsPixel.seedIdx()[pixelArrayIndex]);
0429         }
0430       }
0431     }
0432   };
0433 
0434   struct AddpT5asTrackCandidate {
0435     template <typename TAcc>
0436     ALPAKA_FN_ACC void operator()(TAcc const& acc,
0437                                   uint16_t nLowerModules,
0438                                   PixelQuintupletsConst pixelQuintuplets,
0439                                   TrackCandidates cands,
0440                                   SegmentsPixelConst segmentsPixel,
0441                                   ObjectRangesConst ranges) const {
0442       // implementation is 1D with a single block
0443       static_assert(std::is_same_v<TAcc, ALPAKA_ACCELERATOR_NAMESPACE::Acc1D>, "Should be Acc1D");
0444       ALPAKA_ASSERT_ACC((alpaka::getWorkDiv<alpaka::Grid, alpaka::Blocks>(acc)[0] == 1));
0445 
0446       auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
0447       auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
0448 
0449       int nPixelQuintuplets = pixelQuintuplets.nPixelQuintuplets();
0450       unsigned int pLS_offset = ranges.segmentModuleIndices()[nLowerModules];
0451       for (int pixelQuintupletIndex = globalThreadIdx[0]; pixelQuintupletIndex < nPixelQuintuplets;
0452            pixelQuintupletIndex += gridThreadExtent[0]) {
0453         if (pixelQuintuplets.isDup()[pixelQuintupletIndex])
0454           continue;
0455 
0456         unsigned int trackCandidateIdx =
0457             alpaka::atomicAdd(acc, &cands.nTrackCandidates(), 1u, alpaka::hierarchy::Threads{});
0458         if (trackCandidateIdx >= n_max_pixel_track_candidates)  // No other TCs have been added yet
0459         {
0460 #ifdef WARNINGS
0461           printf("Track Candidate excess alert! Type = pT5");
0462 #endif
0463           alpaka::atomicSub(acc, &cands.nTrackCandidates(), 1u, alpaka::hierarchy::Threads{});
0464           break;
0465 
0466         } else {
0467           alpaka::atomicAdd(acc, &cands.nTrackCandidatespT5(), 1u, alpaka::hierarchy::Threads{});
0468 
0469           float radius = 0.5f * (__H2F(pixelQuintuplets.pixelRadius()[pixelQuintupletIndex]) +
0470                                  __H2F(pixelQuintuplets.quintupletRadius()[pixelQuintupletIndex]));
0471           unsigned int pT5PixelIndex = pixelQuintuplets.pixelSegmentIndices()[pixelQuintupletIndex];
0472           addTrackCandidateToMemory(cands,
0473                                     LSTObjType::pT5,
0474                                     pT5PixelIndex,
0475                                     pixelQuintuplets.quintupletIndices()[pixelQuintupletIndex],
0476                                     pixelQuintuplets.logicalLayers()[pixelQuintupletIndex].data(),
0477                                     pixelQuintuplets.lowerModuleIndices()[pixelQuintupletIndex].data(),
0478                                     pixelQuintuplets.hitIndices()[pixelQuintupletIndex].data(),
0479                                     segmentsPixel.seedIdx()[pT5PixelIndex - pLS_offset],
0480                                     __H2F(pixelQuintuplets.centerX()[pixelQuintupletIndex]),
0481                                     __H2F(pixelQuintuplets.centerY()[pixelQuintupletIndex]),
0482                                     radius,
0483                                     trackCandidateIdx,
0484                                     pixelQuintupletIndex);
0485         }
0486       }
0487     }
0488   };
0489 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE::lst
0490 
0491 ASSERT_DEVICE_MATCHES_HOST_COLLECTION(lst::TrackCandidatesDeviceCollection, lst::TrackCandidatesHostCollection);
0492 
0493 #endif