CMS 3D CMS Logo

Segment.h
Go to the documentation of this file.
1 #ifndef RecoTracker_LSTCore_src_alpaka_Segment_h
2 #define RecoTracker_LSTCore_src_alpaka_Segment_h
3 
5 
12 
13 #include "MiniDoublet.h"
14 #include "Hit.h"
15 
17 
18  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool isTighterTiltedModules_seg(ModulesConst modules, unsigned int moduleIndex) {
19  // The "tighter" tilted modules are the subset of tilted modules that have smaller spacing
20  // This is the same as what was previously considered as"isNormalTiltedModules"
21  // See Figure 9.1 of https://cds.cern.ch/record/2272264/files/CMS-TDR-014.pdf
22  short subdet = modules.subdets()[moduleIndex];
23  short layer = modules.layers()[moduleIndex];
24  short side = modules.sides()[moduleIndex];
25  short rod = modules.rods()[moduleIndex];
26 
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)))));
30  }
31 
32  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool isTighterTiltedModules_seg(short subdet, short layer, short side, short rod) {
33  // The "tighter" tilted modules are the subset of tilted modules that have smaller spacing
34  // This is the same as what was previously considered as"isNormalTiltedModules"
35  // See Figure 9.1 of https://cds.cern.ch/record/2272264/files/CMS-TDR-014.pdf
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)))));
39  }
40 
41  ALPAKA_FN_ACC ALPAKA_FN_INLINE float moduleGapSize_seg(short layer, short ring, short subdet, short side, short rod) {
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, /*10*/ 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, /*10*/ 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, /*10*/ 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, /*10*/ 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, /*10*/ 0.18f, 0.18f, 0.18f, 0.18f, 0.18f}};
51 
52  unsigned int iL = layer - 1;
53  unsigned int iR = ring - 1;
54 
55  float moduleSeparation = 0;
56 
57  if (subdet == Barrel and side == Center) {
58  moduleSeparation = miniDeltaFlat[iL];
59  } else if (isTighterTiltedModules_seg(subdet, layer, side, rod)) {
60  moduleSeparation = miniDeltaTilted[iL];
61  } else if (subdet == Endcap) {
62  moduleSeparation = miniDeltaEndcap[iL][iR];
63  } else //Loose tilted modules
64  {
65  moduleSeparation = miniDeltaLooseTilted[iL];
66  }
67 
68  return moduleSeparation;
69  }
70 
71  ALPAKA_FN_ACC ALPAKA_FN_INLINE float moduleGapSize_seg(ModulesConst modules, unsigned int moduleIndex) {
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, /*10*/ 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, /*10*/ 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, /*10*/ 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, /*10*/ 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, /*10*/ 0.18f, 0.18f, 0.18f, 0.18f, 0.18f}};
81 
82  unsigned int iL = modules.layers()[moduleIndex] - 1;
83  unsigned int iR = modules.rings()[moduleIndex] - 1;
84  short subdet = modules.subdets()[moduleIndex];
85  short side = modules.sides()[moduleIndex];
86 
87  float moduleSeparation = 0;
88 
89  if (subdet == Barrel and side == Center) {
90  moduleSeparation = miniDeltaFlat[iL];
91  } else if (isTighterTiltedModules_seg(modules, moduleIndex)) {
92  moduleSeparation = miniDeltaTilted[iL];
93  } else if (subdet == Endcap) {
94  moduleSeparation = miniDeltaEndcap[iL][iR];
95  } else //Loose tilted modules
96  {
97  moduleSeparation = miniDeltaLooseTilted[iL];
98  }
99 
100  return moduleSeparation;
101  }
102 
103  template <typename TAcc>
104  ALPAKA_FN_ACC ALPAKA_FN_INLINE void dAlphaThreshold(TAcc const& acc,
105  float* dAlphaThresholdValues,
107  MiniDoubletsConst mds,
108  float xIn,
109  float yIn,
110  float zIn,
111  float rtIn,
112  float xOut,
113  float yOut,
114  float zOut,
115  float rtOut,
116  uint16_t innerLowerModuleIndex,
117  uint16_t outerLowerModuleIndex,
118  unsigned int innerMDIndex,
119  unsigned int outerMDIndex) {
120  float sdMuls = (modules.subdets()[innerLowerModuleIndex] == Barrel)
121  ? kMiniMulsPtScaleBarrel[modules.layers()[innerLowerModuleIndex] - 1] * 3.f / ptCut
122  : kMiniMulsPtScaleEndcap[modules.layers()[innerLowerModuleIndex] - 1] * 3.f / ptCut;
123 
124  //more accurate then outer rt - inner rt
125  float segmentDr = alpaka::math::sqrt(acc, (yOut - yIn) * (yOut - yIn) + (xOut - xIn) * (xOut - xIn));
126 
127  const float dAlpha_Bfield =
128  alpaka::math::asin(acc, alpaka::math::min(acc, segmentDr * k2Rinv1GeVf / ptCut, kSinAlphaMax));
129 
130  bool isInnerTilted =
131  modules.subdets()[innerLowerModuleIndex] == Barrel and modules.sides()[innerLowerModuleIndex] != Center;
132  bool isOuterTilted =
133  modules.subdets()[outerLowerModuleIndex] == Barrel and modules.sides()[outerLowerModuleIndex] != Center;
134 
135  float drdzInner = modules.drdzs()[innerLowerModuleIndex];
136  float drdzOuter = modules.drdzs()[outerLowerModuleIndex];
137  float innerModuleGapSize = moduleGapSize_seg(modules, innerLowerModuleIndex);
138  float outerModuleGapSize = moduleGapSize_seg(modules, outerLowerModuleIndex);
139  const float innerminiTilt2 = isInnerTilted
140  ? ((0.5f * 0.5f) * (kPixelPSZpitch * kPixelPSZpitch) * (drdzInner * drdzInner) /
141  (1.f + drdzInner * drdzInner) / (innerModuleGapSize * innerModuleGapSize))
142  : 0;
143 
144  const float outerminiTilt2 = isOuterTilted
145  ? ((0.5f * 0.5f) * (kPixelPSZpitch * kPixelPSZpitch) * (drdzOuter * drdzOuter) /
146  (1.f + drdzOuter * drdzOuter) / (outerModuleGapSize * outerModuleGapSize))
147  : 0;
148 
149  float miniDelta = innerModuleGapSize;
150 
151  float sdLumForInnerMini2;
152  float sdLumForOuterMini2;
153 
154  if (modules.subdets()[innerLowerModuleIndex] == Barrel) {
155  sdLumForInnerMini2 = innerminiTilt2 * (dAlpha_Bfield * dAlpha_Bfield);
156  } else {
157  sdLumForInnerMini2 = (mds.dphis()[innerMDIndex] * mds.dphis()[innerMDIndex]) * (kDeltaZLum * kDeltaZLum) /
158  (mds.dzs()[innerMDIndex] * mds.dzs()[innerMDIndex]);
159  }
160 
161  if (modules.subdets()[outerLowerModuleIndex] == Barrel) {
162  sdLumForOuterMini2 = outerminiTilt2 * (dAlpha_Bfield * dAlpha_Bfield);
163  } else {
164  sdLumForOuterMini2 = (mds.dphis()[outerMDIndex] * mds.dphis()[outerMDIndex]) * (kDeltaZLum * kDeltaZLum) /
165  (mds.dzs()[outerMDIndex] * mds.dzs()[outerMDIndex]);
166  }
167 
168  // Unique stuff for the segment dudes alone
169  float dAlpha_res_inner =
170  0.02f / miniDelta *
171  (modules.subdets()[innerLowerModuleIndex] == Barrel ? 1.0f : alpaka::math::abs(acc, zIn) / rtIn);
172  float dAlpha_res_outer =
173  0.02f / miniDelta *
174  (modules.subdets()[outerLowerModuleIndex] == Barrel ? 1.0f : alpaka::math::abs(acc, zOut) / rtOut);
175 
176  float dAlpha_res = dAlpha_res_inner + dAlpha_res_outer;
177 
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);
180  } else {
181  dAlphaThresholdValues[0] =
182  dAlpha_Bfield + alpaka::math::sqrt(acc, dAlpha_res * dAlpha_res + sdMuls * sdMuls + sdLumForInnerMini2);
183  }
184 
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);
187  } else {
188  dAlphaThresholdValues[1] =
189  dAlpha_Bfield + alpaka::math::sqrt(acc, dAlpha_res * dAlpha_res + sdMuls * sdMuls + sdLumForOuterMini2);
190  }
191 
192  //Inner to outer
193  dAlphaThresholdValues[2] = dAlpha_Bfield + alpaka::math::sqrt(acc, dAlpha_res * dAlpha_res + sdMuls * sdMuls);
194  }
195 
196  ALPAKA_FN_ACC ALPAKA_FN_INLINE void addSegmentToMemory(Segments segments,
197  unsigned int lowerMDIndex,
198  unsigned int upperMDIndex,
199  uint16_t innerLowerModuleIndex,
200  uint16_t outerLowerModuleIndex,
201  unsigned int innerMDAnchorHitIndex,
202  unsigned int outerMDAnchorHitIndex,
203  float dPhi,
204  float dPhiMin,
205  float dPhiMax,
206  float dPhiChange,
207  float dPhiChangeMin,
208  float dPhiChangeMax,
209  unsigned int idx) {
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;
216 
217  segments.dPhis()[idx] = __F2H(dPhi);
218  segments.dPhiMins()[idx] = __F2H(dPhiMin);
219  segments.dPhiMaxs()[idx] = __F2H(dPhiMax);
220  segments.dPhiChanges()[idx] = __F2H(dPhiChange);
221  segments.dPhiChangeMins()[idx] = __F2H(dPhiChangeMin);
222  segments.dPhiChangeMaxs()[idx] = __F2H(dPhiChangeMax);
223  }
224 
225  template <typename TAcc>
226  ALPAKA_FN_ACC ALPAKA_FN_INLINE void addPixelSegmentToMemory(TAcc const& acc,
227  Segments segments,
228  SegmentsPixel segmentsPixel,
229  MiniDoubletsConst mds,
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,
236  float dPhiChange,
237  unsigned int idx,
238  unsigned int pixelSegmentArrayIndex,
239  float score) {
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);
247 
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];
255 
256  //computing circle parameters
257  /*
258  The two anchor hits are r3PCA and r3LH. p3PCA pt, eta, phi is hitIndex1 x, y, z
259  */
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),
263  mds.anchorX()[innerMDIndex] - circleRadius * alpaka::math::sin(acc, circlePhi)};
264  float candidateCenterYs[] = {mds.anchorY()[innerMDIndex] - circleRadius * alpaka::math::cos(acc, circlePhi),
265  mds.anchorY()[innerMDIndex] + circleRadius * alpaka::math::cos(acc, circlePhi)};
266 
267  //check which of the circles can accommodate r3LH better (we won't get perfect agreement)
268  float bestChiSquared = kVerticalModuleSlope;
269  float chiSquared;
270  size_t bestIndex;
271  for (size_t i = 0; i < 2; i++) {
272  chiSquared = alpaka::math::abs(acc,
273  alpaka::math::sqrt(acc,
274  (mds.anchorX()[outerMDIndex] - candidateCenterXs[i]) *
275  (mds.anchorX()[outerMDIndex] - candidateCenterXs[i]) +
276  (mds.anchorY()[outerMDIndex] - candidateCenterYs[i]) *
277  (mds.anchorY()[outerMDIndex] - candidateCenterYs[i])) -
278  circleRadius);
279  if (chiSquared < bestChiSquared) {
280  bestChiSquared = chiSquared;
281  bestIndex = i;
282  }
283  }
284  segmentsPixel.circleCenterX()[pixelSegmentArrayIndex] = candidateCenterXs[bestIndex];
285  segmentsPixel.circleCenterY()[pixelSegmentArrayIndex] = candidateCenterYs[bestIndex];
286  segmentsPixel.circleRadius()[pixelSegmentArrayIndex] = circleRadius;
287  }
288 
289  template <typename TAcc>
290  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool runSegmentDefaultAlgoBarrel(TAcc const& acc,
292  MiniDoubletsConst mds,
293  uint16_t innerLowerModuleIndex,
294  uint16_t outerLowerModuleIndex,
295  unsigned int innerMDIndex,
296  unsigned int outerMDIndex,
297  float& dPhi,
298  float& dPhiMin,
299  float& dPhiMax,
300  float& dPhiChange,
301  float& dPhiChangeMin,
302  float& dPhiChangeMax) {
303  float sdMuls = (modules.subdets()[innerLowerModuleIndex] == Barrel)
304  ? kMiniMulsPtScaleBarrel[modules.layers()[innerLowerModuleIndex] - 1] * 3.f / ptCut
305  : kMiniMulsPtScaleEndcap[modules.layers()[innerLowerModuleIndex] - 1] * 3.f / ptCut;
306 
307  float xIn, yIn, zIn, rtIn, xOut, yOut, zOut, rtOut;
308 
309  xIn = mds.anchorX()[innerMDIndex];
310  yIn = mds.anchorY()[innerMDIndex];
311  zIn = mds.anchorZ()[innerMDIndex];
312  rtIn = mds.anchorRt()[innerMDIndex];
313 
314  xOut = mds.anchorX()[outerMDIndex];
315  yOut = mds.anchorY()[outerMDIndex];
316  zOut = mds.anchorZ()[outerMDIndex];
317  rtOut = mds.anchorRt()[outerMDIndex];
318 
319  float sdSlope = alpaka::math::asin(acc, alpaka::math::min(acc, rtOut * k2Rinv1GeVf / ptCut, kSinAlphaMax));
320  float sdPVoff = 0.1f / rtOut;
321  float dzDrtScale = alpaka::math::tan(acc, sdSlope) / sdSlope; //FIXME: need appropriate value
322 
323  const float zGeom = modules.layers()[innerLowerModuleIndex] <= 2 ? 2.f * kPixelPSZpitch : 2.f * kStrip2SZpitch;
324 
325  float zLo = zIn + (zIn - kDeltaZLum) * (rtOut / rtIn - 1.f) * (zIn > 0.f ? 1.f : dzDrtScale) -
326  zGeom; //slope-correction only on outer end
327  float zHi = zIn + (zIn + kDeltaZLum) * (rtOut / rtIn - 1.f) * (zIn < 0.f ? 1.f : dzDrtScale) + zGeom;
328 
329  if ((zOut < zLo) || (zOut > zHi))
330  return false;
331 
332  float sdCut = sdSlope + alpaka::math::sqrt(acc, sdMuls * sdMuls + sdPVoff * sdPVoff);
333 
334  dPhi = phi_mpi_pi(acc, mds.anchorPhi()[outerMDIndex] - mds.anchorPhi()[innerMDIndex]);
335 
336  if (alpaka::math::abs(acc, dPhi) > sdCut)
337  return false;
338 
339  dPhiChange = phi_mpi_pi(acc, phi(acc, xOut - xIn, yOut - yIn) - mds.anchorPhi()[innerMDIndex]);
340 
341  if (alpaka::math::abs(acc, dPhiChange) > sdCut)
342  return false;
343 
344  float dAlphaThresholdValues[3];
345  dAlphaThreshold(acc,
346  dAlphaThresholdValues,
347  modules,
348  mds,
349  xIn,
350  yIn,
351  zIn,
352  rtIn,
353  xOut,
354  yOut,
355  zOut,
356  rtOut,
357  innerLowerModuleIndex,
358  outerLowerModuleIndex,
359  innerMDIndex,
360  outerMDIndex);
361 
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;
367 
368  float dAlphaInnerMDSegmentThreshold = dAlphaThresholdValues[0];
369  float dAlphaOuterMDSegmentThreshold = dAlphaThresholdValues[1];
370  float dAlphaInnerMDOuterMDThreshold = dAlphaThresholdValues[2];
371 
372  if (alpaka::math::abs(acc, dAlphaInnerMDSegment) >= dAlphaInnerMDSegmentThreshold)
373  return false;
374  if (alpaka::math::abs(acc, dAlphaOuterMDSegment) >= dAlphaOuterMDSegmentThreshold)
375  return false;
376  return alpaka::math::abs(acc, dAlphaInnerMDOuterMD) < dAlphaInnerMDOuterMDThreshold;
377  }
378 
379  template <typename TAcc>
380  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool runSegmentDefaultAlgoEndcap(TAcc const& acc,
382  MiniDoubletsConst mds,
383  uint16_t innerLowerModuleIndex,
384  uint16_t outerLowerModuleIndex,
385  unsigned int innerMDIndex,
386  unsigned int outerMDIndex,
387  float& dPhi,
388  float& dPhiMin,
389  float& dPhiMax,
390  float& dPhiChange,
391  float& dPhiChangeMin,
392  float& dPhiChangeMax) {
393  float xIn, yIn, zIn, rtIn, xOut, yOut, zOut, rtOut;
394 
395  xIn = mds.anchorX()[innerMDIndex];
396  yIn = mds.anchorY()[innerMDIndex];
397  zIn = mds.anchorZ()[innerMDIndex];
398  rtIn = mds.anchorRt()[innerMDIndex];
399 
400  xOut = mds.anchorX()[outerMDIndex];
401  yOut = mds.anchorY()[outerMDIndex];
402  zOut = mds.anchorZ()[outerMDIndex];
403  rtOut = mds.anchorRt()[outerMDIndex];
404 
405  bool outerLayerEndcapTwoS =
406  (modules.subdets()[outerLowerModuleIndex] == Endcap) && (modules.moduleType()[outerLowerModuleIndex] == TwoS);
407 
408  float sdSlope = alpaka::math::asin(acc, alpaka::math::min(acc, rtOut * k2Rinv1GeVf / ptCut, kSinAlphaMax));
409  float disks2SMinRadius = 60.f;
410 
411  float rtGeom = ((rtIn < disks2SMinRadius && rtOut < disks2SMinRadius)
412  ? (2.f * kPixelPSZpitch)
413  : ((rtIn < disks2SMinRadius || rtOut < disks2SMinRadius) ? (kPixelPSZpitch + kStrip2SZpitch)
414  : (2.f * kStrip2SZpitch)));
415 
416  //cut 0 - z compatibility
417  if (zIn * zOut < 0)
418  return false;
419 
420  float dz = zOut - zIn;
421  float dLum = alpaka::math::copysign(acc, kDeltaZLum, zIn);
422  float drtDzScale = sdSlope / alpaka::math::tan(acc, sdSlope);
423 
424  float rtLo = alpaka::math::max(
425  acc, rtIn * (1.f + dz / (zIn + dLum) * drtDzScale) - rtGeom, rtIn - 0.5f * rtGeom); //rt should increase
426  float rtHi = rtIn * (zOut - dLum) / (zIn - dLum) +
427  rtGeom; //dLum for luminous; rGeom for measurement size; no tanTheta_loc(pt) correction
428 
429  // Completeness
430  if ((rtOut < rtLo) || (rtOut > rtHi))
431  return false;
432 
433  dPhi = phi_mpi_pi(acc, mds.anchorPhi()[outerMDIndex] - mds.anchorPhi()[innerMDIndex]);
434 
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]);
439 
440  dPhiMax = alpaka::math::abs(acc, dPhiPos_high) > alpaka::math::abs(acc, dPhiPos_low) ? dPhiPos_high : dPhiPos_low;
441  dPhiMin = alpaka::math::abs(acc, dPhiPos_high) > alpaka::math::abs(acc, dPhiPos_low) ? dPhiPos_low : dPhiPos_high;
442  } else {
443  dPhiMax = dPhi;
444  dPhiMin = dPhi;
445  }
446  if (alpaka::math::abs(acc, dPhi) > sdCut)
447  return false;
448 
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);
453 
454  if (alpaka::math::abs(acc, dPhiChange) > sdCut)
455  return false;
456 
457  float dAlphaThresholdValues[3];
458  dAlphaThreshold(acc,
459  dAlphaThresholdValues,
460  modules,
461  mds,
462  xIn,
463  yIn,
464  zIn,
465  rtIn,
466  xOut,
467  yOut,
468  zOut,
469  rtOut,
470  innerLowerModuleIndex,
471  outerLowerModuleIndex,
472  innerMDIndex,
473  outerMDIndex);
474 
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;
480 
481  float dAlphaInnerMDSegmentThreshold = dAlphaThresholdValues[0];
482  float dAlphaOuterMDSegmentThreshold = dAlphaThresholdValues[1];
483  float dAlphaInnerMDOuterMDThreshold = dAlphaThresholdValues[2];
484 
485  if (alpaka::math::abs(acc, dAlphaInnerMDSegment) >= dAlphaInnerMDSegmentThreshold)
486  return false;
487  if (alpaka::math::abs(acc, dAlphaOuterMDSegment) >= dAlphaOuterMDSegmentThreshold)
488  return false;
489  return alpaka::math::abs(acc, dAlphaInnerMDOuterMD) < dAlphaInnerMDOuterMDThreshold;
490  }
491 
492  template <typename TAcc>
493  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool runSegmentDefaultAlgo(TAcc const& acc,
495  MiniDoubletsConst mds,
496  uint16_t innerLowerModuleIndex,
497  uint16_t outerLowerModuleIndex,
498  unsigned int innerMDIndex,
499  unsigned int outerMDIndex,
500  float& dPhi,
501  float& dPhiMin,
502  float& dPhiMax,
503  float& dPhiChange,
504  float& dPhiChangeMin,
505  float& dPhiChangeMax) {
506  if (modules.subdets()[innerLowerModuleIndex] == Barrel and modules.subdets()[outerLowerModuleIndex] == Barrel) {
507  return runSegmentDefaultAlgoBarrel(acc,
508  modules,
509  mds,
510  innerLowerModuleIndex,
511  outerLowerModuleIndex,
512  innerMDIndex,
513  outerMDIndex,
514  dPhi,
515  dPhiMin,
516  dPhiMax,
517  dPhiChange,
518  dPhiChangeMin,
519  dPhiChangeMax);
520  } else {
521  return runSegmentDefaultAlgoEndcap(acc,
522  modules,
523  mds,
524  innerLowerModuleIndex,
525  outerLowerModuleIndex,
526  innerMDIndex,
527  outerMDIndex,
528  dPhi,
529  dPhiMin,
530  dPhiMax,
531  dPhiChange,
532  dPhiChangeMin,
533  dPhiChangeMax);
534  }
535  }
536 
537  struct CreateSegments {
538  template <typename TAcc>
539  ALPAKA_FN_ACC void operator()(TAcc const& acc,
541  MiniDoubletsConst mds,
542  MiniDoubletsOccupancyConst mdsOccupancy,
543  Segments segments,
544  SegmentsOccupancy segmentsOccupancy,
545  ObjectRangesConst ranges) const {
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);
550 
551  for (uint16_t innerLowerModuleIndex = globalBlockIdx[2]; innerLowerModuleIndex < modules.nLowerModules();
552  innerLowerModuleIndex += gridBlockExtent[2]) {
553  unsigned int nInnerMDs = mdsOccupancy.nMDs()[innerLowerModuleIndex];
554  if (nInnerMDs == 0)
555  continue;
556 
557  unsigned int nConnectedModules = modules.nConnectedModules()[innerLowerModuleIndex];
558 
559  for (uint16_t outerLowerModuleArrayIdx = blockThreadIdx[1]; outerLowerModuleArrayIdx < nConnectedModules;
560  outerLowerModuleArrayIdx += blockThreadExtent[1]) {
561  uint16_t outerLowerModuleIndex = modules.moduleMap()[innerLowerModuleIndex][outerLowerModuleArrayIdx];
562 
563  unsigned int nOuterMDs = mdsOccupancy.nMDs()[outerLowerModuleIndex];
564 
565  unsigned int limit = nInnerMDs * nOuterMDs;
566 
567  if (limit == 0)
568  continue;
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)
573  continue;
574 
575  unsigned int innerMDIndex = ranges.mdRanges()[innerLowerModuleIndex][0] + innerMDArrayIdx;
576  unsigned int outerMDIndex = ranges.mdRanges()[outerLowerModuleIndex][0] + outerMDArrayIdx;
577 
578  float dPhi, dPhiMin, dPhiMax, dPhiChange, dPhiChangeMin, dPhiChangeMax;
579 
580  unsigned int innerMiniDoubletAnchorHitIndex = mds.anchorHitIndices()[innerMDIndex];
581  unsigned int outerMiniDoubletAnchorHitIndex = mds.anchorHitIndices()[outerMDIndex];
582  dPhiMin = 0;
583  dPhiMax = 0;
584  dPhiChangeMin = 0;
585  dPhiChangeMax = 0;
586  if (runSegmentDefaultAlgo(acc,
587  modules,
588  mds,
589  innerLowerModuleIndex,
590  outerLowerModuleIndex,
591  innerMDIndex,
592  outerMDIndex,
593  dPhi,
594  dPhiMin,
595  dPhiMax,
596  dPhiChange,
597  dPhiChangeMin,
598  dPhiChangeMax)) {
599  unsigned int totOccupancySegments =
600  alpaka::atomicAdd(acc,
601  &segmentsOccupancy.totOccupancySegments()[innerLowerModuleIndex],
602  1u,
603  alpaka::hierarchy::Threads{});
604  if (static_cast<int>(totOccupancySegments) >= ranges.segmentModuleOccupancy()[innerLowerModuleIndex]) {
605 #ifdef WARNINGS
606  printf("Segment excess alert! Module index = %d\n", innerLowerModuleIndex);
607 #endif
608  } else {
609  unsigned int segmentModuleIdx = alpaka::atomicAdd(
610  acc, &segmentsOccupancy.nSegments()[innerLowerModuleIndex], 1u, alpaka::hierarchy::Threads{});
611  unsigned int segmentIdx = ranges.segmentModuleIndices()[innerLowerModuleIndex] + segmentModuleIdx;
612 
613  addSegmentToMemory(segments,
614  innerMDIndex,
615  outerMDIndex,
616  innerLowerModuleIndex,
617  outerLowerModuleIndex,
618  innerMiniDoubletAnchorHitIndex,
619  outerMiniDoubletAnchorHitIndex,
620  dPhi,
621  dPhiMin,
622  dPhiMax,
623  dPhiChange,
624  dPhiChangeMin,
625  dPhiChangeMax,
626  segmentIdx);
627  }
628  }
629  }
630  }
631  }
632  }
633  };
634 
636  template <typename TAcc>
637  ALPAKA_FN_ACC void operator()(TAcc const& acc,
640  MiniDoubletsConst mds) const {
641  // implementation is 1D with a single block
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));
644 
645  auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
646  auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
647 
648  // Initialize variables in shared memory and set to 0
649  int& nTotalSegments = alpaka::declareSharedVar<int, __COUNTER__>(acc);
651  nTotalSegments = 0;
652  }
653  alpaka::syncBlockThreads(acc);
654 
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;
659  continue;
660  }
661 
662  short module_rings = modules.rings()[i];
663  short module_layers = modules.layers()[i];
664  short module_subdets = modules.subdets()[i];
665  float module_eta = alpaka::math::abs(acc, modules.eta()[i]);
666 
667  int category_number;
668  if (module_layers <= 3 && module_subdets == 5)
669  category_number = 0;
670  else if (module_layers >= 4 && module_subdets == 5)
671  category_number = 1;
672  else if (module_layers <= 2 && module_subdets == 4 && module_rings >= 11)
673  category_number = 2;
674  else if (module_layers >= 3 && module_subdets == 4 && module_rings >= 8)
675  category_number = 2;
676  else if (module_layers <= 2 && module_subdets == 4 && module_rings <= 10)
677  category_number = 3;
678  else if (module_layers >= 3 && module_subdets == 4 && module_rings <= 7)
679  category_number = 3;
680  else
681  category_number = -1;
682 
683  int eta_number;
684  if (module_eta < 0.75f)
685  eta_number = 0;
686  else if (module_eta < 1.5f)
687  eta_number = 1;
688  else if (module_eta < 2.25f)
689  eta_number = 2;
690  else if (module_eta < 3.0f)
691  eta_number = 3;
692  else
693  eta_number = -1;
694 
695  int occupancy;
696  if (category_number == 0 && eta_number == 0)
697  occupancy = 572;
698  else if (category_number == 0 && eta_number == 1)
699  occupancy = 300;
700  else if (category_number == 0 && eta_number == 2)
701  occupancy = 183;
702  else if (category_number == 0 && eta_number == 3)
703  occupancy = 62;
704  else if (category_number == 1 && eta_number == 0)
705  occupancy = 191;
706  else if (category_number == 1 && eta_number == 1)
707  occupancy = 128;
708  else if (category_number == 2 && eta_number == 1)
709  occupancy = 107;
710  else if (category_number == 2 && eta_number == 2)
711  occupancy = 102;
712  else if (category_number == 3 && eta_number == 1)
713  occupancy = 64;
714  else if (category_number == 3 && eta_number == 2)
715  occupancy = 79;
716  else if (category_number == 3 && eta_number == 3)
717  occupancy = 85;
718  else {
719  occupancy = 0;
720 #ifdef WARNINGS
721  printf("Unhandled case in createSegmentArrayRanges! Module index = %i\n", i);
722 #endif
723  }
724 
725  int nTotSegs = alpaka::atomicAdd(acc, &nTotalSegments, occupancy, alpaka::hierarchy::Threads{});
726  ranges.segmentModuleIndices()[i] = nTotSegs;
727  ranges.segmentModuleOccupancy()[i] = occupancy;
728  }
729 
730  // Wait for all threads to finish before reporting final values
731  alpaka::syncBlockThreads(acc);
733  ranges.segmentModuleIndices()[modules.nLowerModules()] = nTotalSegments;
734  ranges.nTotalSegs() = nTotalSegments;
735  }
736  }
737  };
738 
740  template <typename TAcc>
741  ALPAKA_FN_ACC void operator()(TAcc const& acc,
743  SegmentsOccupancyConst segmentsOccupancy,
744  ObjectRanges ranges) const {
745  // implementation is 1D with a single block
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));
748 
749  auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
750  auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
751 
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;
756  } else {
757  ranges.segmentRanges()[i][0] = ranges.segmentModuleIndices()[i];
758  ranges.segmentRanges()[i][1] = ranges.segmentModuleIndices()[i] + segmentsOccupancy.nSegments()[i] - 1;
759  }
760  }
761  }
762  };
763 
765  template <typename TAcc>
766  ALPAKA_FN_ACC void operator()(TAcc const& acc,
769  HitsConst hits,
770  MiniDoublets mds,
771  Segments segments,
772  SegmentsPixel segmentsPixel,
773  unsigned int* hitIndices0,
774  unsigned int* hitIndices1,
775  unsigned int* hitIndices2,
776  unsigned int* hitIndices3,
777  float* dPhiChange,
778  uint16_t pixelModuleIndex,
779  int size) const {
780  auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
781  auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
782 
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;
787 
788  addMDToMemory(acc,
789  mds,
790  hits,
791  modules,
792  hitIndices0[tid],
793  hitIndices1[tid],
794  pixelModuleIndex,
795  0,
796  0,
797  0,
798  0,
799  0,
800  0,
801  0,
802  0,
803  innerMDIndex);
804  addMDToMemory(acc,
805  mds,
806  hits,
807  modules,
808  hitIndices2[tid],
809  hitIndices3[tid],
810  pixelModuleIndex,
811  0,
812  0,
813  0,
814  0,
815  0,
816  0,
817  0,
818  0,
819  outerMDIndex);
820 
821  //in outer hits - pt, eta, phi
822  float slope = alpaka::math::sinh(acc, hits.ys()[mds.outerHitIndices()[innerMDIndex]]);
823  float intercept =
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;
828 
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]];
835  segments,
836  segmentsPixel,
837  mds,
838  innerMDIndex,
839  outerMDIndex,
840  pixelModuleIndex,
841  hits1,
842  hitIndices0[tid],
843  hitIndices2[tid],
844  dPhiChange[tid],
845  pixelSegmentIndex,
846  tid,
847  score_lsq);
848  }
849  }
850  };
851 } // namespace ALPAKA_ACCELERATOR_NAMESPACE::lst
852 
853 #endif
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
Definition: Segment.h:766
MiniDoubletsSoA::View MiniDoublets
ALPAKA_FN_ACC constexpr bool once_per_block(TAcc const &acc)
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)
Definition: Segment.h:104
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kMiniMulsPtScaleBarrel[6]
Definition: Common.h:42
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float k2Rinv1GeVf
Definition: Common.h:49
ObjectRangesSoA::View ObjectRanges
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kMiniMulsPtScaleEndcap[5]
Definition: Common.h:44
static const double slope[3]
const std::vector< int > & module_rings()
Definition: LSTEff.cc:3384
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
ALPAKA_FN_ACC void operator()(TAcc const &acc, ModulesConst modules, SegmentsOccupancyConst segmentsOccupancy, ObjectRanges ranges) const
Definition: Segment.h:741
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool isTighterTiltedModules_seg(ModulesConst modules, unsigned int moduleIndex)
Definition: Segment.h:18
ALPAKA_FN_ACC void operator()(TAcc const &acc, ModulesConst modules, MiniDoubletsConst mds, MiniDoubletsOccupancyConst mdsOccupancy, Segments segments, SegmentsOccupancy segmentsOccupancy, ObjectRangesConst ranges) const
Definition: Segment.h:539
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)
Definition: Segment.h:290
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)
Definition: Segment.h:226
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kVerticalModuleSlope
Definition: Common.h:61
#define __F2H
Definition: Common.h:49
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float ptCut
Definition: Common.h:52
SegmentsOccupancySoA::View SegmentsOccupancy
Definition: SegmentsSoA.h:50
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)
Definition: Segment.h:196
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE float phi_mpi_pi(TAcc const &acc, float x)
Definition: Hit.h:19
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)
Definition: MiniDoublet.h:17
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE float phi(TAcc const &acc, float x, float y)
Definition: Hit.h:29
T sqrt(T t)
Definition: SSEVec.h:23
SegmentsPixelSoA::View SegmentsPixel
Definition: SegmentsSoA.h:52
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
MiniDoubletsSoA::ConstView MiniDoubletsConst
HitsSoA::ConstView HitsConst
Definition: HitsSoA.h:33
ModulesSoA::ConstView ModulesConst
Definition: ModulesSoA.h:47
const std::vector< int > & module_layers()
Definition: LSTEff.cc:3388
SegmentsSoA::View Segments
Definition: SegmentsSoA.h:48
string ranges
Definition: diffTwoXMLs.py:79
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kPixelPSZpitch
Definition: Common.h:54
const std::vector< int > & module_subdets()
Definition: LSTEff.cc:3288
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)
Definition: Segment.h:493
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kStrip2SZpitch
Definition: Common.h:56
SegmentsOccupancySoA::ConstView SegmentsOccupancyConst
Definition: SegmentsSoA.h:51
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)
Definition: Segment.h:380
ALPAKA_FN_ACC void operator()(TAcc const &acc, ModulesConst modules, ObjectRanges ranges, MiniDoubletsConst mds) const
Definition: Segment.h:637
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kSinAlphaMax
Definition: Common.h:51
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kDeltaZLum
Definition: Common.h:53
ALPAKA_FN_ACC ALPAKA_FN_INLINE float moduleGapSize_seg(short layer, short ring, short subdet, short side, short rod)
Definition: Segment.h:41
MiniDoubletsOccupancySoA::ConstView MiniDoubletsOccupancyConst
T1 atomicAdd(T1 *a, T2 b)
Definition: cudaCompat.h:61