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