1 #ifndef RecoTracker_LSTCore_src_alpaka_Segment_h 2 #define RecoTracker_LSTCore_src_alpaka_Segment_h 22 short subdet =
modules.subdets()[moduleIndex];
23 short layer =
modules.layers()[moduleIndex];
25 short rod =
modules.rods()[moduleIndex];
27 return (subdet == Barrel) && (((
side !=
Center) && (layer == 3)) ||
28 ((
side ==
NegZ) && (((layer == 2) && (rod > 5)) || ((layer == 1) && (rod > 9)))) ||
29 ((
side ==
PosZ) && (((layer == 2) && (rod < 8)) || ((layer == 1) && (rod < 4)))));
36 return (subdet == Barrel) && (((
side !=
Center) && (layer == 3)) ||
37 ((
side ==
NegZ) && (((layer == 2) && (rod > 5)) || ((layer == 1) && (rod > 9)))) ||
38 ((
side ==
PosZ) && (((layer == 2) && (rod < 8)) || ((layer == 1) && (rod < 4)))));
42 static constexpr float miniDeltaTilted[3] = {0.26f, 0.26f, 0.26f};
43 static constexpr float miniDeltaFlat[6] = {0.26f, 0.16f, 0.16f, 0.18f, 0.18f, 0.18f};
44 static constexpr float miniDeltaLooseTilted[3] = {0.4f, 0.4f, 0.4f};
45 static constexpr float miniDeltaEndcap[5][15] = {
46 {0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f},
47 {0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f},
48 {0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f},
49 {0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f},
50 {0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f}};
52 unsigned int iL = layer - 1;
53 unsigned int iR =
ring - 1;
55 float moduleSeparation = 0;
58 moduleSeparation = miniDeltaFlat[iL];
60 moduleSeparation = miniDeltaTilted[iL];
61 }
else if (subdet ==
Endcap) {
62 moduleSeparation = miniDeltaEndcap[iL][iR];
65 moduleSeparation = miniDeltaLooseTilted[iL];
68 return moduleSeparation;
72 static constexpr float miniDeltaTilted[3] = {0.26f, 0.26f, 0.26f};
73 static constexpr float miniDeltaFlat[6] = {0.26f, 0.16f, 0.16f, 0.18f, 0.18f, 0.18f};
74 static constexpr float miniDeltaLooseTilted[3] = {0.4f, 0.4f, 0.4f};
75 static constexpr float miniDeltaEndcap[5][15] = {
76 {0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f},
77 {0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f},
78 {0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f},
79 {0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f},
80 {0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.4f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f, 0.18f}};
82 unsigned int iL =
modules.layers()[moduleIndex] - 1;
83 unsigned int iR =
modules.rings()[moduleIndex] - 1;
84 short subdet =
modules.subdets()[moduleIndex];
87 float moduleSeparation = 0;
90 moduleSeparation = miniDeltaFlat[iL];
92 moduleSeparation = miniDeltaTilted[iL];
93 }
else if (subdet ==
Endcap) {
94 moduleSeparation = miniDeltaEndcap[iL][iR];
97 moduleSeparation = miniDeltaLooseTilted[iL];
100 return moduleSeparation;
103 template <
typename TAcc>
105 float* dAlphaThresholdValues,
116 uint16_t innerLowerModuleIndex,
117 uint16_t outerLowerModuleIndex,
118 unsigned int innerMDIndex,
119 unsigned int outerMDIndex) {
120 float sdMuls = (
modules.subdets()[innerLowerModuleIndex] == Barrel)
125 float segmentDr =
alpaka::math::sqrt(acc, (yOut - yIn) * (yOut - yIn) + (xOut - xIn) * (xOut - xIn));
127 const float dAlpha_Bfield =
131 modules.subdets()[innerLowerModuleIndex] == Barrel and
modules.sides()[innerLowerModuleIndex] !=
Center;
133 modules.subdets()[outerLowerModuleIndex] == Barrel and
modules.sides()[outerLowerModuleIndex] !=
Center;
135 float drdzInner =
modules.drdzs()[innerLowerModuleIndex];
136 float drdzOuter =
modules.drdzs()[outerLowerModuleIndex];
139 const float innerminiTilt2 = isInnerTilted
141 (1.
f + drdzInner * drdzInner) / (innerModuleGapSize * innerModuleGapSize))
144 const float outerminiTilt2 = isOuterTilted
146 (1.
f + drdzOuter * drdzOuter) / (outerModuleGapSize * outerModuleGapSize))
149 float miniDelta = innerModuleGapSize;
151 float sdLumForInnerMini2;
152 float sdLumForOuterMini2;
154 if (
modules.subdets()[innerLowerModuleIndex] == Barrel) {
155 sdLumForInnerMini2 = innerminiTilt2 * (dAlpha_Bfield * dAlpha_Bfield);
157 sdLumForInnerMini2 = (mds.dphis()[innerMDIndex] * mds.dphis()[innerMDIndex]) * (
kDeltaZLum *
kDeltaZLum) /
158 (mds.dzs()[innerMDIndex] * mds.dzs()[innerMDIndex]);
161 if (
modules.subdets()[outerLowerModuleIndex] == Barrel) {
162 sdLumForOuterMini2 = outerminiTilt2 * (dAlpha_Bfield * dAlpha_Bfield);
164 sdLumForOuterMini2 = (mds.dphis()[outerMDIndex] * mds.dphis()[outerMDIndex]) * (
kDeltaZLum *
kDeltaZLum) /
165 (mds.dzs()[outerMDIndex] * mds.dzs()[outerMDIndex]);
169 float dAlpha_res_inner =
172 float dAlpha_res_outer =
176 float dAlpha_res = dAlpha_res_inner + dAlpha_res_outer;
178 if (
modules.subdets()[innerLowerModuleIndex] == Barrel and
modules.sides()[innerLowerModuleIndex] ==
Center) {
179 dAlphaThresholdValues[0] = dAlpha_Bfield +
alpaka::math::sqrt(acc, dAlpha_res * dAlpha_res + sdMuls * sdMuls);
181 dAlphaThresholdValues[0] =
182 dAlpha_Bfield +
alpaka::math::sqrt(acc, dAlpha_res * dAlpha_res + sdMuls * sdMuls + sdLumForInnerMini2);
185 if (
modules.subdets()[outerLowerModuleIndex] == Barrel and
modules.sides()[outerLowerModuleIndex] ==
Center) {
186 dAlphaThresholdValues[1] = dAlpha_Bfield +
alpaka::math::sqrt(acc, dAlpha_res * dAlpha_res + sdMuls * sdMuls);
188 dAlphaThresholdValues[1] =
189 dAlpha_Bfield +
alpaka::math::sqrt(acc, dAlpha_res * dAlpha_res + sdMuls * sdMuls + sdLumForOuterMini2);
193 dAlphaThresholdValues[2] = dAlpha_Bfield +
alpaka::math::sqrt(acc, dAlpha_res * dAlpha_res + sdMuls * sdMuls);
197 unsigned int lowerMDIndex,
198 unsigned int upperMDIndex,
199 uint16_t innerLowerModuleIndex,
200 uint16_t outerLowerModuleIndex,
201 unsigned int innerMDAnchorHitIndex,
202 unsigned int outerMDAnchorHitIndex,
210 segments.mdIndices()[
idx][0] = lowerMDIndex;
211 segments.mdIndices()[
idx][1] = upperMDIndex;
212 segments.innerLowerModuleIndices()[
idx] = innerLowerModuleIndex;
213 segments.outerLowerModuleIndices()[
idx] = outerLowerModuleIndex;
214 segments.innerMiniDoubletAnchorHitIndices()[
idx] = innerMDAnchorHitIndex;
215 segments.outerMiniDoubletAnchorHitIndices()[
idx] = outerMDAnchorHitIndex;
217 segments.dPhis()[
idx] =
__F2H(dPhi);
220 segments.dPhiChanges()[
idx] =
__F2H(dPhiChange);
221 segments.dPhiChangeMins()[
idx] =
__F2H(dPhiChangeMin);
222 segments.dPhiChangeMaxs()[
idx] =
__F2H(dPhiChangeMax);
225 template <
typename TAcc>
230 unsigned int innerMDIndex,
231 unsigned int outerMDIndex,
232 uint16_t pixelModuleIndex,
233 unsigned int hitIdxs[4],
234 unsigned int innerAnchorHitIndex,
235 unsigned int outerAnchorHitIndex,
238 unsigned int pixelSegmentArrayIndex,
240 segments.mdIndices()[
idx][0] = innerMDIndex;
241 segments.mdIndices()[
idx][1] = outerMDIndex;
242 segments.innerLowerModuleIndices()[
idx] = pixelModuleIndex;
243 segments.outerLowerModuleIndices()[
idx] = pixelModuleIndex;
244 segments.innerMiniDoubletAnchorHitIndices()[
idx] = innerAnchorHitIndex;
245 segments.outerMiniDoubletAnchorHitIndices()[
idx] = outerAnchorHitIndex;
246 segments.dPhiChanges()[
idx] =
__F2H(dPhiChange);
248 segmentsPixel.isDup()[pixelSegmentArrayIndex] =
false;
249 segmentsPixel.partOfPT5()[pixelSegmentArrayIndex] =
false;
250 segmentsPixel.score()[pixelSegmentArrayIndex] =
score;
251 segmentsPixel.pLSHitsIdxs()[pixelSegmentArrayIndex].x = hitIdxs[0];
252 segmentsPixel.pLSHitsIdxs()[pixelSegmentArrayIndex].y = hitIdxs[1];
253 segmentsPixel.pLSHitsIdxs()[pixelSegmentArrayIndex].z = hitIdxs[2];
254 segmentsPixel.pLSHitsIdxs()[pixelSegmentArrayIndex].w = hitIdxs[3];
260 float circleRadius = mds.outerX()[innerMDIndex] / (2 *
k2Rinv1GeVf);
261 float circlePhi = mds.outerZ()[innerMDIndex];
262 float candidateCenterXs[] = {mds.anchorX()[innerMDIndex] + circleRadius *
alpaka::math::sin(acc, circlePhi),
264 float candidateCenterYs[] = {mds.anchorY()[innerMDIndex] - circleRadius *
alpaka::math::cos(acc, circlePhi),
271 for (
size_t i = 0;
i < 2;
i++) {
274 (mds.anchorX()[outerMDIndex] - candidateCenterXs[
i]) *
275 (mds.anchorX()[outerMDIndex] - candidateCenterXs[
i]) +
276 (mds.anchorY()[outerMDIndex] - candidateCenterYs[
i]) *
277 (mds.anchorY()[outerMDIndex] - candidateCenterYs[
i])) -
279 if (chiSquared < bestChiSquared) {
280 bestChiSquared = chiSquared;
284 segmentsPixel.circleCenterX()[pixelSegmentArrayIndex] = candidateCenterXs[bestIndex];
285 segmentsPixel.circleCenterY()[pixelSegmentArrayIndex] = candidateCenterYs[bestIndex];
286 segmentsPixel.circleRadius()[pixelSegmentArrayIndex] = circleRadius;
289 template <
typename TAcc>
293 uint16_t innerLowerModuleIndex,
294 uint16_t outerLowerModuleIndex,
295 unsigned int innerMDIndex,
296 unsigned int outerMDIndex,
301 float& dPhiChangeMin,
302 float& dPhiChangeMax) {
303 float sdMuls = (
modules.subdets()[innerLowerModuleIndex] == Barrel)
307 float xIn, yIn, zIn, rtIn, xOut, yOut, zOut, rtOut;
309 xIn = mds.anchorX()[innerMDIndex];
310 yIn = mds.anchorY()[innerMDIndex];
311 zIn = mds.anchorZ()[innerMDIndex];
312 rtIn = mds.anchorRt()[innerMDIndex];
314 xOut = mds.anchorX()[outerMDIndex];
315 yOut = mds.anchorY()[outerMDIndex];
316 zOut = mds.anchorZ()[outerMDIndex];
317 rtOut = mds.anchorRt()[outerMDIndex];
320 float sdPVoff = 0.1f / rtOut;
325 float zLo = zIn + (zIn -
kDeltaZLum) * (rtOut / rtIn - 1.
f) * (zIn > 0.f ? 1.f : dzDrtScale) -
327 float zHi = zIn + (zIn +
kDeltaZLum) * (rtOut / rtIn - 1.
f) * (zIn < 0.f ? 1.f : dzDrtScale) + zGeom;
329 if ((zOut < zLo) || (zOut > zHi))
334 dPhi =
phi_mpi_pi(acc, mds.anchorPhi()[outerMDIndex] - mds.anchorPhi()[innerMDIndex]);
339 dPhiChange =
phi_mpi_pi(acc,
phi(acc, xOut - xIn, yOut - yIn) - mds.anchorPhi()[innerMDIndex]);
344 float dAlphaThresholdValues[3];
346 dAlphaThresholdValues,
357 innerLowerModuleIndex,
358 outerLowerModuleIndex,
362 float innerMDAlpha = mds.dphichanges()[innerMDIndex];
363 float outerMDAlpha = mds.dphichanges()[outerMDIndex];
364 float dAlphaInnerMDSegment = innerMDAlpha - dPhiChange;
365 float dAlphaOuterMDSegment = outerMDAlpha - dPhiChange;
366 float dAlphaInnerMDOuterMD = innerMDAlpha - outerMDAlpha;
368 float dAlphaInnerMDSegmentThreshold = dAlphaThresholdValues[0];
369 float dAlphaOuterMDSegmentThreshold = dAlphaThresholdValues[1];
370 float dAlphaInnerMDOuterMDThreshold = dAlphaThresholdValues[2];
372 if (
alpaka::math::abs(acc, dAlphaInnerMDSegment) >= dAlphaInnerMDSegmentThreshold)
374 if (
alpaka::math::abs(acc, dAlphaOuterMDSegment) >= dAlphaOuterMDSegmentThreshold)
376 return alpaka::math::abs(acc, dAlphaInnerMDOuterMD) < dAlphaInnerMDOuterMDThreshold;
379 template <
typename TAcc>
383 uint16_t innerLowerModuleIndex,
384 uint16_t outerLowerModuleIndex,
385 unsigned int innerMDIndex,
386 unsigned int outerMDIndex,
391 float& dPhiChangeMin,
392 float& dPhiChangeMax) {
393 float xIn, yIn, zIn, rtIn, xOut, yOut, zOut, rtOut;
395 xIn = mds.anchorX()[innerMDIndex];
396 yIn = mds.anchorY()[innerMDIndex];
397 zIn = mds.anchorZ()[innerMDIndex];
398 rtIn = mds.anchorRt()[innerMDIndex];
400 xOut = mds.anchorX()[outerMDIndex];
401 yOut = mds.anchorY()[outerMDIndex];
402 zOut = mds.anchorZ()[outerMDIndex];
403 rtOut = mds.anchorRt()[outerMDIndex];
405 bool outerLayerEndcapTwoS =
409 float disks2SMinRadius = 60.f;
411 float rtGeom = ((rtIn < disks2SMinRadius && rtOut < disks2SMinRadius)
420 float dz = zOut - zIn;
421 float dLum = alpaka::math::copysign(acc,
kDeltaZLum, zIn);
425 acc, rtIn * (1.
f +
dz / (zIn + dLum) * drtDzScale) - rtGeom, rtIn - 0.5
f * rtGeom);
426 float rtHi = rtIn * (zOut - dLum) / (zIn - dLum) +
430 if ((rtOut < rtLo) || (rtOut > rtHi))
433 dPhi =
phi_mpi_pi(acc, mds.anchorPhi()[outerMDIndex] - mds.anchorPhi()[innerMDIndex]);
435 float sdCut = sdSlope;
436 if (outerLayerEndcapTwoS) {
437 float dPhiPos_high =
phi_mpi_pi(acc, mds.anchorHighEdgePhi()[outerMDIndex] - mds.anchorPhi()[innerMDIndex]);
438 float dPhiPos_low =
phi_mpi_pi(acc, mds.anchorLowEdgePhi()[outerMDIndex] - mds.anchorPhi()[innerMDIndex]);
449 float dzFrac =
dz / zIn;
450 dPhiChange = dPhi / dzFrac * (1.f + dzFrac);
451 dPhiChangeMin =
dPhiMin / dzFrac * (1.f + dzFrac);
452 dPhiChangeMax =
dPhiMax / dzFrac * (1.f + dzFrac);
457 float dAlphaThresholdValues[3];
459 dAlphaThresholdValues,
470 innerLowerModuleIndex,
471 outerLowerModuleIndex,
475 float innerMDAlpha = mds.dphichanges()[innerMDIndex];
476 float outerMDAlpha = mds.dphichanges()[outerMDIndex];
477 float dAlphaInnerMDSegment = innerMDAlpha - dPhiChange;
478 float dAlphaOuterMDSegment = outerMDAlpha - dPhiChange;
479 float dAlphaInnerMDOuterMD = innerMDAlpha - outerMDAlpha;
481 float dAlphaInnerMDSegmentThreshold = dAlphaThresholdValues[0];
482 float dAlphaOuterMDSegmentThreshold = dAlphaThresholdValues[1];
483 float dAlphaInnerMDOuterMDThreshold = dAlphaThresholdValues[2];
485 if (
alpaka::math::abs(acc, dAlphaInnerMDSegment) >= dAlphaInnerMDSegmentThreshold)
487 if (
alpaka::math::abs(acc, dAlphaOuterMDSegment) >= dAlphaOuterMDSegmentThreshold)
489 return alpaka::math::abs(acc, dAlphaInnerMDOuterMD) < dAlphaInnerMDOuterMDThreshold;
492 template <
typename TAcc>
496 uint16_t innerLowerModuleIndex,
497 uint16_t outerLowerModuleIndex,
498 unsigned int innerMDIndex,
499 unsigned int outerMDIndex,
504 float& dPhiChangeMin,
505 float& dPhiChangeMax) {
506 if (
modules.subdets()[innerLowerModuleIndex] == Barrel and
modules.subdets()[outerLowerModuleIndex] == Barrel) {
510 innerLowerModuleIndex,
511 outerLowerModuleIndex,
524 innerLowerModuleIndex,
525 outerLowerModuleIndex,
538 template <
typename TAcc>
546 auto const globalBlockIdx = alpaka::getIdx<alpaka::Grid, alpaka::Blocks>(acc);
547 auto const blockThreadIdx = alpaka::getIdx<alpaka::Block, alpaka::Threads>(acc);
548 auto const gridBlockExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Blocks>(acc);
549 auto const blockThreadExtent = alpaka::getWorkDiv<alpaka::Block, alpaka::Threads>(acc);
551 for (uint16_t innerLowerModuleIndex = globalBlockIdx[2]; innerLowerModuleIndex <
modules.nLowerModules();
552 innerLowerModuleIndex += gridBlockExtent[2]) {
553 unsigned int nInnerMDs = mdsOccupancy.nMDs()[innerLowerModuleIndex];
557 unsigned int nConnectedModules =
modules.nConnectedModules()[innerLowerModuleIndex];
559 for (uint16_t outerLowerModuleArrayIdx = blockThreadIdx[1]; outerLowerModuleArrayIdx < nConnectedModules;
560 outerLowerModuleArrayIdx += blockThreadExtent[1]) {
561 uint16_t outerLowerModuleIndex =
modules.moduleMap()[innerLowerModuleIndex][outerLowerModuleArrayIdx];
563 unsigned int nOuterMDs = mdsOccupancy.nMDs()[outerLowerModuleIndex];
565 unsigned int limit = nInnerMDs * nOuterMDs;
569 for (
unsigned int hitIndex = blockThreadIdx[2]; hitIndex <
limit; hitIndex += blockThreadExtent[2]) {
570 unsigned int innerMDArrayIdx = hitIndex / nOuterMDs;
571 unsigned int outerMDArrayIdx = hitIndex % nOuterMDs;
572 if (outerMDArrayIdx >= nOuterMDs)
575 unsigned int innerMDIndex =
ranges.mdRanges()[innerLowerModuleIndex][0] + innerMDArrayIdx;
576 unsigned int outerMDIndex =
ranges.mdRanges()[outerLowerModuleIndex][0] + outerMDArrayIdx;
578 float dPhi,
dPhiMin,
dPhiMax, dPhiChange, dPhiChangeMin, dPhiChangeMax;
580 unsigned int innerMiniDoubletAnchorHitIndex = mds.anchorHitIndices()[innerMDIndex];
581 unsigned int outerMiniDoubletAnchorHitIndex = mds.anchorHitIndices()[outerMDIndex];
589 innerLowerModuleIndex,
590 outerLowerModuleIndex,
599 unsigned int totOccupancySegments =
601 &segmentsOccupancy.totOccupancySegments()[innerLowerModuleIndex],
603 alpaka::hierarchy::Threads{});
604 if (static_cast<int>(totOccupancySegments) >=
ranges.segmentModuleOccupancy()[innerLowerModuleIndex]) {
606 printf(
"Segment excess alert! Module index = %d\n", innerLowerModuleIndex);
610 acc, &segmentsOccupancy.nSegments()[innerLowerModuleIndex], 1u, alpaka::hierarchy::Threads{});
611 unsigned int segmentIdx =
ranges.segmentModuleIndices()[innerLowerModuleIndex] + segmentModuleIdx;
616 innerLowerModuleIndex,
617 outerLowerModuleIndex,
618 innerMiniDoubletAnchorHitIndex,
619 outerMiniDoubletAnchorHitIndex,
636 template <
typename TAcc>
642 static_assert(std::is_same_v<TAcc, ALPAKA_ACCELERATOR_NAMESPACE::Acc1D>,
"Should be Acc1D");
643 ALPAKA_ASSERT_ACC((alpaka::getWorkDiv<alpaka::Grid, alpaka::Blocks>(acc)[0] == 1));
645 auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
646 auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
649 int& nTotalSegments = alpaka::declareSharedVar<int, __COUNTER__>(acc);
653 alpaka::syncBlockThreads(acc);
655 for (uint16_t
i = globalThreadIdx[0];
i <
modules.nLowerModules();
i += gridThreadExtent[0]) {
656 if (
modules.nConnectedModules()[
i] == 0) {
657 ranges.segmentModuleIndices()[
i] = nTotalSegments;
658 ranges.segmentModuleOccupancy()[
i] = 0;
672 else if (module_layers <= 2 && module_subdets == 4 && module_rings >= 11)
681 category_number = -1;
684 if (module_eta < 0.75
f)
686 else if (module_eta < 1.5
f)
688 else if (module_eta < 2.25
f)
690 else if (module_eta < 3.0
f)
696 if (category_number == 0 && eta_number == 0)
698 else if (category_number == 0 && eta_number == 1)
700 else if (category_number == 0 && eta_number == 2)
702 else if (category_number == 0 && eta_number == 3)
704 else if (category_number == 1 && eta_number == 0)
706 else if (category_number == 1 && eta_number == 1)
708 else if (category_number == 2 && eta_number == 1)
710 else if (category_number == 2 && eta_number == 2)
712 else if (category_number == 3 && eta_number == 1)
714 else if (category_number == 3 && eta_number == 2)
716 else if (category_number == 3 && eta_number == 3)
721 printf(
"Unhandled case in createSegmentArrayRanges! Module index = %i\n",
i);
725 int nTotSegs =
alpaka::atomicAdd(acc, &nTotalSegments, occupancy, alpaka::hierarchy::Threads{});
726 ranges.segmentModuleIndices()[
i] = nTotSegs;
727 ranges.segmentModuleOccupancy()[
i] = occupancy;
731 alpaka::syncBlockThreads(acc);
733 ranges.segmentModuleIndices()[
modules.nLowerModules()] = nTotalSegments;
734 ranges.nTotalSegs() = nTotalSegments;
740 template <
typename TAcc>
746 static_assert(std::is_same_v<TAcc, ALPAKA_ACCELERATOR_NAMESPACE::Acc1D>,
"Should be Acc1D");
747 ALPAKA_ASSERT_ACC((alpaka::getWorkDiv<alpaka::Grid, alpaka::Blocks>(acc)[0] == 1));
749 auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
750 auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
752 for (uint16_t
i = globalThreadIdx[0];
i <
modules.nLowerModules();
i += gridThreadExtent[0]) {
753 if (segmentsOccupancy.nSegments()[
i] == 0) {
754 ranges.segmentRanges()[
i][0] = -1;
755 ranges.segmentRanges()[
i][1] = -1;
758 ranges.segmentRanges()[
i][1] =
ranges.segmentModuleIndices()[
i] + segmentsOccupancy.nSegments()[
i] - 1;
765 template <
typename TAcc>
773 unsigned int* hitIndices0,
774 unsigned int* hitIndices1,
775 unsigned int* hitIndices2,
776 unsigned int* hitIndices3,
778 uint16_t pixelModuleIndex,
780 auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
781 auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
783 for (
int tid = globalThreadIdx[2]; tid < size; tid += gridThreadExtent[2]) {
784 unsigned int innerMDIndex =
ranges.miniDoubletModuleIndices()[pixelModuleIndex] + 2 * (tid);
785 unsigned int outerMDIndex =
ranges.miniDoubletModuleIndices()[pixelModuleIndex] + 2 * (tid) + 1;
786 unsigned int pixelSegmentIndex =
ranges.segmentModuleIndices()[pixelModuleIndex] + tid;
822 float slope = alpaka::math::sinh(acc,
hits.ys()[mds.outerHitIndices()[innerMDIndex]]);
824 hits.zs()[mds.anchorHitIndices()[innerMDIndex]] -
slope *
hits.rts()[mds.anchorHitIndices()[innerMDIndex]];
825 float score_lsq = (
hits.rts()[mds.anchorHitIndices()[outerMDIndex]] *
slope + intercept) -
826 (
hits.zs()[mds.anchorHitIndices()[outerMDIndex]]);
827 score_lsq = score_lsq * score_lsq;
829 unsigned int hits1[Params_pLS::kHits];
830 hits1[0] =
hits.idxs()[mds.anchorHitIndices()[innerMDIndex]];
831 hits1[1] =
hits.idxs()[mds.anchorHitIndices()[outerMDIndex]];
832 hits1[2] =
hits.idxs()[mds.outerHitIndices()[innerMDIndex]];
833 hits1[3] =
hits.idxs()[mds.outerHitIndices()[outerMDIndex]];
ALPAKA_FN_ACC void operator()(TAcc const &acc, ModulesConst modules, ObjectRangesConst ranges, HitsConst hits, MiniDoublets mds, Segments segments, SegmentsPixel segmentsPixel, unsigned int *hitIndices0, unsigned int *hitIndices1, unsigned int *hitIndices2, unsigned int *hitIndices3, float *dPhiChange, uint16_t pixelModuleIndex, int size) const
MiniDoubletsSoA::View MiniDoublets
ALPAKA_FN_ACC ALPAKA_FN_INLINE void dAlphaThreshold(TAcc const &acc, float *dAlphaThresholdValues, ModulesConst modules, MiniDoubletsConst mds, float xIn, float yIn, float zIn, float rtIn, float xOut, float yOut, float zOut, float rtOut, uint16_t innerLowerModuleIndex, uint16_t outerLowerModuleIndex, unsigned int innerMDIndex, unsigned int outerMDIndex)
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kMiniMulsPtScaleBarrel[6]
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float k2Rinv1GeVf
ObjectRangesSoA::View ObjectRanges
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kMiniMulsPtScaleEndcap[5]
static const double slope[3]
const std::vector< int > & module_rings()
Sin< T >::type sin(const T &t)
ALPAKA_FN_ACC void operator()(TAcc const &acc, ModulesConst modules, SegmentsOccupancyConst segmentsOccupancy, ObjectRanges ranges) const
ALPAKA_FN_ACC int side(int ieta, int iphi)
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool isTighterTiltedModules_seg(ModulesConst modules, unsigned int moduleIndex)
ALPAKA_FN_ACC void operator()(TAcc const &acc, ModulesConst modules, MiniDoubletsConst mds, MiniDoubletsOccupancyConst mdsOccupancy, Segments segments, SegmentsOccupancy segmentsOccupancy, ObjectRangesConst ranges) const
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool runSegmentDefaultAlgoBarrel(TAcc const &acc, ModulesConst modules, MiniDoubletsConst mds, uint16_t innerLowerModuleIndex, uint16_t outerLowerModuleIndex, unsigned int innerMDIndex, unsigned int outerMDIndex, float &dPhi, float &dPhiMin, float &dPhiMax, float &dPhiChange, float &dPhiChangeMin, float &dPhiChangeMax)
ALPAKA_FN_ACC ALPAKA_FN_INLINE void addPixelSegmentToMemory(TAcc const &acc, Segments segments, SegmentsPixel segmentsPixel, MiniDoubletsConst mds, unsigned int innerMDIndex, unsigned int outerMDIndex, uint16_t pixelModuleIndex, unsigned int hitIdxs[4], unsigned int innerAnchorHitIndex, unsigned int outerAnchorHitIndex, float dPhiChange, unsigned int idx, unsigned int pixelSegmentArrayIndex, float score)
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kVerticalModuleSlope
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float ptCut
SegmentsOccupancySoA::View SegmentsOccupancy
ALPAKA_FN_ACC ALPAKA_FN_INLINE void addSegmentToMemory(Segments segments, unsigned int lowerMDIndex, unsigned int upperMDIndex, uint16_t innerLowerModuleIndex, uint16_t outerLowerModuleIndex, unsigned int innerMDAnchorHitIndex, unsigned int outerMDAnchorHitIndex, float dPhi, float dPhiMin, float dPhiMax, float dPhiChange, float dPhiChangeMin, float dPhiChangeMax, unsigned int idx)
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE float phi_mpi_pi(TAcc const &acc, float x)
ALPAKA_FN_ACC ALPAKA_FN_INLINE void addMDToMemory(TAcc const &acc, MiniDoublets mds, HitsConst hits, ModulesConst modules, unsigned int lowerHitIdx, unsigned int upperHitIdx, uint16_t lowerModuleIdx, float dz, float dPhi, float dPhiChange, float shiftedX, float shiftedY, float shiftedZ, float noShiftedDphi, float noShiftedDPhiChange, unsigned int idx)
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE float phi(TAcc const &acc, float x, float y)
SegmentsPixelSoA::View SegmentsPixel
Cos< T >::type cos(const T &t)
Tan< T >::type tan(const T &t)
Abs< T >::type abs(const T &t)
MiniDoubletsSoA::ConstView MiniDoubletsConst
HitsSoA::ConstView HitsConst
ModulesSoA::ConstView ModulesConst
const std::vector< int > & module_layers()
SegmentsSoA::View Segments
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kPixelPSZpitch
const std::vector< int > & module_subdets()
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool runSegmentDefaultAlgo(TAcc const &acc, ModulesConst modules, MiniDoubletsConst mds, uint16_t innerLowerModuleIndex, uint16_t outerLowerModuleIndex, unsigned int innerMDIndex, unsigned int outerMDIndex, float &dPhi, float &dPhiMin, float &dPhiMax, float &dPhiChange, float &dPhiChangeMin, float &dPhiChangeMax)
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kStrip2SZpitch
SegmentsOccupancySoA::ConstView SegmentsOccupancyConst
ObjectRangesSoA::ConstView ObjectRangesConst
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool runSegmentDefaultAlgoEndcap(TAcc const &acc, ModulesConst modules, MiniDoubletsConst mds, uint16_t innerLowerModuleIndex, uint16_t outerLowerModuleIndex, unsigned int innerMDIndex, unsigned int outerMDIndex, float &dPhi, float &dPhiMin, float &dPhiMax, float &dPhiChange, float &dPhiChangeMin, float &dPhiChangeMax)
ALPAKA_FN_ACC void operator()(TAcc const &acc, ModulesConst modules, ObjectRanges ranges, MiniDoubletsConst mds) const
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kSinAlphaMax
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kDeltaZLum
ALPAKA_FN_ACC ALPAKA_FN_INLINE float moduleGapSize_seg(short layer, short ring, short subdet, short side, short rod)
ALPAKA_ASSERT_ACC(offsets)
MiniDoubletsOccupancySoA::ConstView MiniDoubletsOccupancyConst
T1 atomicAdd(T1 *a, T2 b)