CMS 3D CMS Logo

PixelQuintuplet.h
Go to the documentation of this file.
1 #ifndef RecoTracker_LSTCore_src_alpaka_PixelQuintuplet_h
2 #define RecoTracker_LSTCore_src_alpaka_PixelQuintuplet_h
3 
12 
13 #include "Hit.h"
14 #include "PixelTriplet.h"
15 
17  ALPAKA_FN_ACC ALPAKA_FN_INLINE void addPixelQuintupletToMemory(ModulesConst modules,
19  SegmentsConst segments,
20  QuintupletsConst quintuplets,
21  PixelQuintuplets pixelQuintuplets,
22  unsigned int pixelIndex,
23  unsigned int t5Index,
24  unsigned int pixelQuintupletIndex,
25  float rzChiSquared,
26  float rPhiChiSquared,
27  float rPhiChiSquaredInwards,
28  float score,
29  float eta,
30  float phi,
31  float pixelRadius,
32  float quintupletRadius,
33  float centerX,
34  float centerY) {
35  pixelQuintuplets.pixelSegmentIndices()[pixelQuintupletIndex] = pixelIndex;
36  pixelQuintuplets.quintupletIndices()[pixelQuintupletIndex] = t5Index;
37  pixelQuintuplets.isDup()[pixelQuintupletIndex] = false;
38  pixelQuintuplets.score()[pixelQuintupletIndex] = __F2H(score);
39  pixelQuintuplets.eta()[pixelQuintupletIndex] = __F2H(eta);
40  pixelQuintuplets.phi()[pixelQuintupletIndex] = __F2H(phi);
41 
42  pixelQuintuplets.pixelRadius()[pixelQuintupletIndex] = __F2H(pixelRadius);
43  pixelQuintuplets.quintupletRadius()[pixelQuintupletIndex] = __F2H(quintupletRadius);
44  pixelQuintuplets.centerX()[pixelQuintupletIndex] = __F2H(centerX);
45  pixelQuintuplets.centerY()[pixelQuintupletIndex] = __F2H(centerY);
46 
47  pixelQuintuplets.logicalLayers()[pixelQuintupletIndex][0] = 0;
48  pixelQuintuplets.logicalLayers()[pixelQuintupletIndex][1] = 0;
49  pixelQuintuplets.logicalLayers()[pixelQuintupletIndex][2] = quintuplets.logicalLayers()[t5Index][0];
50  pixelQuintuplets.logicalLayers()[pixelQuintupletIndex][3] = quintuplets.logicalLayers()[t5Index][1];
51  pixelQuintuplets.logicalLayers()[pixelQuintupletIndex][4] = quintuplets.logicalLayers()[t5Index][2];
52  pixelQuintuplets.logicalLayers()[pixelQuintupletIndex][5] = quintuplets.logicalLayers()[t5Index][3];
53  pixelQuintuplets.logicalLayers()[pixelQuintupletIndex][6] = quintuplets.logicalLayers()[t5Index][4];
54 
55  pixelQuintuplets.lowerModuleIndices()[pixelQuintupletIndex][0] = segments.innerLowerModuleIndices()[pixelIndex];
56  pixelQuintuplets.lowerModuleIndices()[pixelQuintupletIndex][1] = segments.outerLowerModuleIndices()[pixelIndex];
57  pixelQuintuplets.lowerModuleIndices()[pixelQuintupletIndex][2] = quintuplets.lowerModuleIndices()[t5Index][0];
58  pixelQuintuplets.lowerModuleIndices()[pixelQuintupletIndex][3] = quintuplets.lowerModuleIndices()[t5Index][1];
59  pixelQuintuplets.lowerModuleIndices()[pixelQuintupletIndex][4] = quintuplets.lowerModuleIndices()[t5Index][2];
60  pixelQuintuplets.lowerModuleIndices()[pixelQuintupletIndex][5] = quintuplets.lowerModuleIndices()[t5Index][3];
61  pixelQuintuplets.lowerModuleIndices()[pixelQuintupletIndex][6] = quintuplets.lowerModuleIndices()[t5Index][4];
62 
63  unsigned int pixelInnerMD = segments.mdIndices()[pixelIndex][0];
64  unsigned int pixelOuterMD = segments.mdIndices()[pixelIndex][1];
65 
66  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][0] = mds.anchorHitIndices()[pixelInnerMD];
67  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][1] = mds.outerHitIndices()[pixelInnerMD];
68  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][2] = mds.anchorHitIndices()[pixelOuterMD];
69  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][3] = mds.outerHitIndices()[pixelOuterMD];
70 
71  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][4] = quintuplets.hitIndices()[t5Index][0];
72  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][5] = quintuplets.hitIndices()[t5Index][1];
73  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][6] = quintuplets.hitIndices()[t5Index][2];
74  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][7] = quintuplets.hitIndices()[t5Index][3];
75  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][8] = quintuplets.hitIndices()[t5Index][4];
76  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][9] = quintuplets.hitIndices()[t5Index][5];
77  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][10] = quintuplets.hitIndices()[t5Index][6];
78  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][11] = quintuplets.hitIndices()[t5Index][7];
79  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][12] = quintuplets.hitIndices()[t5Index][8];
80  pixelQuintuplets.hitIndices()[pixelQuintupletIndex][13] = quintuplets.hitIndices()[t5Index][9];
81 
82  pixelQuintuplets.rzChiSquared()[pixelQuintupletIndex] = rzChiSquared;
83  pixelQuintuplets.rPhiChiSquared()[pixelQuintupletIndex] = rPhiChiSquared;
84  pixelQuintuplets.rPhiChiSquaredInwards()[pixelQuintupletIndex] = rPhiChiSquaredInwards;
85  }
86 
87  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passPT5RZChiSquaredCuts(ModulesConst modules,
88  uint16_t lowerModuleIndex1,
89  uint16_t lowerModuleIndex2,
90  uint16_t lowerModuleIndex3,
91  uint16_t lowerModuleIndex4,
92  uint16_t lowerModuleIndex5,
93  float rzChiSquared) {
94  const int layer1 =
95  modules.layers()[lowerModuleIndex1] + 6 * (modules.subdets()[lowerModuleIndex1] == Endcap) +
96  5 * (modules.subdets()[lowerModuleIndex1] == Endcap and modules.moduleType()[lowerModuleIndex1] == TwoS);
97  const int layer2 =
98  modules.layers()[lowerModuleIndex2] + 6 * (modules.subdets()[lowerModuleIndex2] == Endcap) +
99  5 * (modules.subdets()[lowerModuleIndex2] == Endcap and modules.moduleType()[lowerModuleIndex2] == TwoS);
100  const int layer3 =
101  modules.layers()[lowerModuleIndex3] + 6 * (modules.subdets()[lowerModuleIndex3] == Endcap) +
102  5 * (modules.subdets()[lowerModuleIndex3] == Endcap and modules.moduleType()[lowerModuleIndex3] == TwoS);
103  const int layer4 =
104  modules.layers()[lowerModuleIndex4] + 6 * (modules.subdets()[lowerModuleIndex4] == Endcap) +
105  5 * (modules.subdets()[lowerModuleIndex4] == Endcap and modules.moduleType()[lowerModuleIndex4] == TwoS);
106  const int layer5 =
107  modules.layers()[lowerModuleIndex5] + 6 * (modules.subdets()[lowerModuleIndex5] == Endcap) +
108  5 * (modules.subdets()[lowerModuleIndex5] == Endcap and modules.moduleType()[lowerModuleIndex5] == TwoS);
109 
110  if (layer1 == 1 and layer2 == 2 and layer3 == 3) {
111  if (layer4 == 12 and layer5 == 13) {
112  return rzChiSquared < 451.141f;
113  } else if (layer4 == 4 and layer5 == 12) {
114  return rzChiSquared < 392.654f;
115  } else if (layer4 == 4 and layer5 == 5) {
116  return rzChiSquared < 225.322f;
117  } else if (layer4 == 7 and layer5 == 13) {
118  return rzChiSquared < 595.546f;
119  } else if (layer4 == 7 and layer5 == 8) {
120  return rzChiSquared < 196.111f;
121  }
122  } else if (layer1 == 1 and layer2 == 2 and layer3 == 7) {
123  if (layer4 == 13 and layer5 == 14) {
124  return rzChiSquared < 297.446f;
125  } else if (layer4 == 8 and layer5 == 14) {
126  return rzChiSquared < 451.141f;
127  } else if (layer4 == 8 and layer5 == 9) {
128  return rzChiSquared < 518.339f;
129  }
130  } else if (layer1 == 1 and layer2 == 7 and layer3 == 8) {
131  if (layer4 == 9 and layer5 == 10) {
132  return rzChiSquared < 341.75f;
133  } else if (layer4 == 9 and layer5 == 15) {
134  return rzChiSquared < 341.75f;
135  }
136  } else if (layer1 == 2 and layer2 == 3 and layer3 == 4) {
137  if (layer4 == 12 and layer5 == 13) {
138  return rzChiSquared < 392.655f;
139  } else if (layer4 == 5 and layer5 == 12) {
140  return rzChiSquared < 341.75f;
141  } else if (layer4 == 5 and layer5 == 6) {
142  return rzChiSquared < 112.537f;
143  }
144  } else if (layer1 == 2 and layer2 == 3 and layer4 == 7) {
145  if (layer4 == 13 and layer5 == 14) {
146  return rzChiSquared < 595.545f;
147  } else if (layer4 == 8 and layer5 == 14) {
148  return rzChiSquared < 74.198f;
149  }
150  } else if (layer1 == 2 and layer2 == 7 and layer3 == 8) {
151  if (layer4 == 14 and layer5 == 15) {
152  return rzChiSquared < 518.339f;
153  } else if (layer4 == 9 and layer5 == 10) {
154  return rzChiSquared < 8.046f;
155  } else if (layer4 == 9 and layer5 == 15) {
156  return rzChiSquared < 451.141f;
157  }
158  } else if (layer1 == 3 and layer2 == 7 and layer3 == 8 and layer4 == 14 and layer5 == 15) {
159  return rzChiSquared < 56.207f;
160  } else if (layer1 == 7 and layer2 == 8 and layer3 == 9) {
161  if (layer4 == 10 and layer5 == 11) {
162  return rzChiSquared < 64.578f;
163  } else if (layer4 == 10 and layer5 == 16) {
164  return rzChiSquared < 85.250f;
165  } else if (layer4 == 15 and layer5 == 16) {
166  return rzChiSquared < 85.250f;
167  }
168  }
169  return true;
170  }
171 
172  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passPT5RPhiChiSquaredCuts(ModulesConst modules,
173  uint16_t lowerModuleIndex1,
174  uint16_t lowerModuleIndex2,
175  uint16_t lowerModuleIndex3,
176  uint16_t lowerModuleIndex4,
177  uint16_t lowerModuleIndex5,
178  float rPhiChiSquared) {
179  const int layer1 =
180  modules.layers()[lowerModuleIndex1] + 6 * (modules.subdets()[lowerModuleIndex1] == Endcap) +
181  5 * (modules.subdets()[lowerModuleIndex1] == Endcap and modules.moduleType()[lowerModuleIndex1] == TwoS);
182  const int layer2 =
183  modules.layers()[lowerModuleIndex2] + 6 * (modules.subdets()[lowerModuleIndex2] == Endcap) +
184  5 * (modules.subdets()[lowerModuleIndex2] == Endcap and modules.moduleType()[lowerModuleIndex2] == TwoS);
185  const int layer3 =
186  modules.layers()[lowerModuleIndex3] + 6 * (modules.subdets()[lowerModuleIndex3] == Endcap) +
187  5 * (modules.subdets()[lowerModuleIndex3] == Endcap and modules.moduleType()[lowerModuleIndex3] == TwoS);
188  const int layer4 =
189  modules.layers()[lowerModuleIndex4] + 6 * (modules.subdets()[lowerModuleIndex4] == Endcap) +
190  5 * (modules.subdets()[lowerModuleIndex4] == Endcap and modules.moduleType()[lowerModuleIndex4] == TwoS);
191  const int layer5 =
192  modules.layers()[lowerModuleIndex5] + 6 * (modules.subdets()[lowerModuleIndex5] == Endcap) +
193  5 * (modules.subdets()[lowerModuleIndex5] == Endcap and modules.moduleType()[lowerModuleIndex5] == TwoS);
194 
195  if (layer1 == 1 and layer2 == 2 and layer3 == 3) {
196  if (layer4 == 12 and layer5 == 13) {
197  return rPhiChiSquared < 48.921f;
198  } else if (layer4 == 4 and layer5 == 12) {
199  return rPhiChiSquared < 97.948f;
200  } else if (layer4 == 4 and layer5 == 5) {
201  return rPhiChiSquared < 129.3f;
202  } else if (layer4 == 7 and layer5 == 13) {
203  return rPhiChiSquared < 56.21f;
204  } else if (layer4 == 7 and layer5 == 8) {
205  return rPhiChiSquared < 74.198f;
206  }
207  } else if (layer1 == 1 and layer2 == 2 and layer3 == 7) {
208  if (layer4 == 13 and layer5 == 14) {
209  return rPhiChiSquared < 21.265f;
210  } else if (layer4 == 8 and layer5 == 14) {
211  return rPhiChiSquared < 37.058f;
212  } else if (layer4 == 8 and layer5 == 9) {
213  return rPhiChiSquared < 42.578f;
214  }
215  } else if (layer1 == 1 and layer2 == 7 and layer3 == 8) {
216  if (layer4 == 9 and layer5 == 10) {
217  return rPhiChiSquared < 32.253f;
218  } else if (layer4 == 9 and layer5 == 15) {
219  return rPhiChiSquared < 37.058f;
220  }
221  } else if (layer1 == 2 and layer2 == 3 and layer3 == 4) {
222  if (layer4 == 12 and layer5 == 13) {
223  return rPhiChiSquared < 97.947f;
224  } else if (layer4 == 5 and layer5 == 12) {
225  return rPhiChiSquared < 129.3f;
226  } else if (layer4 == 5 and layer5 == 6) {
227  return rPhiChiSquared < 170.68f;
228  }
229  } else if (layer1 == 2 and layer2 == 3 and layer3 == 7) {
230  if (layer4 == 13 and layer5 == 14) {
231  return rPhiChiSquared < 48.92f;
232  } else if (layer4 == 8 and layer5 == 14) {
233  return rPhiChiSquared < 74.2f;
234  }
235  } else if (layer1 == 2 and layer2 == 7 and layer3 == 8) {
236  if (layer4 == 14 and layer5 == 15) {
237  return rPhiChiSquared < 42.58f;
238  } else if (layer4 == 9 and layer5 == 10) {
239  return rPhiChiSquared < 37.06f;
240  } else if (layer4 == 9 and layer5 == 15) {
241  return rPhiChiSquared < 48.92f;
242  }
243  } else if (layer1 == 3 and layer2 == 7 and layer3 == 8 and layer4 == 14 and layer5 == 15) {
244  return rPhiChiSquared < 85.25f;
245  } else if (layer1 == 7 and layer2 == 8 and layer3 == 9) {
246  if (layer4 == 10 and layer5 == 11) {
247  return rPhiChiSquared < 42.58f;
248  } else if (layer4 == 10 and layer5 == 16) {
249  return rPhiChiSquared < 37.06f;
250  } else if (layer4 == 15 and layer5 == 16) {
251  return rPhiChiSquared < 37.06f;
252  }
253  }
254  return true;
255  }
256 
257  template <typename TAcc>
258  ALPAKA_FN_ACC ALPAKA_FN_INLINE float computeChiSquaredpT5(TAcc const& acc,
259  unsigned int nPoints,
260  float* xs,
261  float* ys,
262  float* delta1,
263  float* delta2,
264  float* slopes,
265  bool* isFlat,
266  float g,
267  float f,
268  float radius) {
269  /*
270  Given values of (g, f, radius) and a set of points (and its uncertainties) compute chi squared
271  */
272  float c = g * g + f * f - radius * radius;
273  float chiSquared = 0.f;
274  float absArctanSlope, angleM, xPrime, yPrime, sigma2;
275  for (size_t i = 0; i < nPoints; i++) {
276  absArctanSlope = ((slopes[i] != kVerticalModuleSlope) ? alpaka::math::abs(acc, alpaka::math::atan(acc, slopes[i]))
277  : kPi / 2.f);
278  if (xs[i] > 0 and ys[i] > 0) {
279  angleM = kPi / 2.f - absArctanSlope;
280  } else if (xs[i] < 0 and ys[i] > 0) {
281  angleM = absArctanSlope + kPi / 2.f;
282  } else if (xs[i] < 0 and ys[i] < 0) {
283  angleM = -(absArctanSlope + kPi / 2.f);
284  } else if (xs[i] > 0 and ys[i] < 0) {
285  angleM = -(kPi / 2.f - absArctanSlope);
286  } else {
287  angleM = 0;
288  }
289  if (not isFlat[i]) {
290  xPrime = xs[i] * alpaka::math::cos(acc, angleM) + ys[i] * alpaka::math::sin(acc, angleM);
291  yPrime = ys[i] * alpaka::math::cos(acc, angleM) - xs[i] * alpaka::math::sin(acc, angleM);
292  } else {
293  xPrime = xs[i];
294  yPrime = ys[i];
295  }
296  sigma2 = 4 * ((xPrime * delta1[i]) * (xPrime * delta1[i]) + (yPrime * delta2[i]) * (yPrime * delta2[i]));
297  chiSquared += (xs[i] * xs[i] + ys[i] * ys[i] - 2 * g * xs[i] - 2 * f * ys[i] + c) *
298  (xs[i] * xs[i] + ys[i] * ys[i] - 2 * g * xs[i] - 2 * f * ys[i] + c) / (sigma2);
299  }
300  return chiSquared;
301  }
302 
303  template <typename TAcc>
304  ALPAKA_FN_ACC ALPAKA_FN_INLINE void computeSigmasForRegression_pT5(TAcc const& acc,
306  const uint16_t* lowerModuleIndices,
307  float* delta1,
308  float* delta2,
309  float* slopes,
310  bool* isFlat,
311  unsigned int nPoints = 5,
312  bool anchorHits = true) {
313  /*
314  bool anchorHits required to deal with a weird edge case wherein
315  the hits ultimately used in the regression are anchor hits, but the
316  lower modules need not all be Pixel Modules (in case of PS). Similarly,
317  when we compute the chi squared for the non-anchor hits, the "partner module"
318  need not always be a PS strip module, but all non-anchor hits sit on strip
319  modules.
320  */
321  ModuleType moduleType;
322  short moduleSubdet, moduleSide;
323  float inv1 = kWidthPS / kWidth2S;
324  float inv2 = kPixelPSZpitch / kWidth2S;
325  float inv3 = kStripPSZpitch / kWidth2S;
326  for (size_t i = 0; i < nPoints; i++) {
327  moduleType = modules.moduleType()[lowerModuleIndices[i]];
328  moduleSubdet = modules.subdets()[lowerModuleIndices[i]];
329  moduleSide = modules.sides()[lowerModuleIndices[i]];
330  const float& drdz = modules.drdzs()[lowerModuleIndices[i]];
331  slopes[i] = modules.dxdys()[lowerModuleIndices[i]];
332  //category 1 - barrel PS flat
333  if (moduleSubdet == Barrel and moduleType == PS and moduleSide == Center) {
334  delta1[i] = inv1;
335  delta2[i] = inv1;
336  slopes[i] = -999.f;
337  isFlat[i] = true;
338  }
339  //category 2 - barrel 2S
340  else if (moduleSubdet == Barrel and moduleType == TwoS) {
341  delta1[i] = 1.f;
342  delta2[i] = 1.f;
343  slopes[i] = -999.f;
344  isFlat[i] = true;
345  }
346  //category 3 - barrel PS tilted
347  else if (moduleSubdet == Barrel and moduleType == PS and moduleSide != Center) {
348  delta1[i] = inv1;
349  isFlat[i] = false;
350 
351  if (anchorHits) {
352  delta2[i] = (inv2 * drdz / alpaka::math::sqrt(acc, 1 + drdz * drdz));
353  } else {
354  delta2[i] = (inv3 * drdz / alpaka::math::sqrt(acc, 1 + drdz * drdz));
355  }
356  }
357  //category 4 - endcap PS
358  else if (moduleSubdet == Endcap and moduleType == PS) {
359  delta1[i] = inv1;
360  isFlat[i] = false;
361  /*
362  despite the type of the module layer of the lower module index,
363  all anchor hits are on the pixel side and all non-anchor hits are
364  on the strip side!
365  */
366  if (anchorHits) {
367  delta2[i] = inv2;
368  } else {
369  delta2[i] = inv3;
370  }
371  }
372  //category 5 - endcap 2S
373  else if (moduleSubdet == Endcap and moduleType == TwoS) {
374  delta1[i] = 1.f;
375  delta2[i] = 500.f * inv1;
376  isFlat[i] = false;
377  }
378 #ifdef WARNINGS
379  else {
380  printf("ERROR!!!!! I SHOULDN'T BE HERE!!!! subdet = %d, type = %d, side = %d\n",
381  moduleSubdet,
382  moduleType,
383  moduleSide);
384  }
385 #endif
386  }
387  }
388 
389  template <typename TAcc>
390  ALPAKA_FN_ACC ALPAKA_FN_INLINE float computePT5RPhiChiSquared(TAcc const& acc,
392  uint16_t* lowerModuleIndices,
393  float g,
394  float f,
395  float radius,
396  float* xs,
397  float* ys) {
398  /*
399  Compute circle parameters from 3 pixel hits, and then use them to compute the chi squared for the outer hits
400  */
401 
402  float delta1[5], delta2[5], slopes[5];
403  bool isFlat[5];
404  float chiSquared = 0;
405 
406  computeSigmasForRegression_pT5(acc, modules, lowerModuleIndices, delta1, delta2, slopes, isFlat);
407  chiSquared = computeChiSquaredpT5(acc, 5, xs, ys, delta1, delta2, slopes, isFlat, g, f, radius);
408 
409  return chiSquared;
410  }
411 
412  ALPAKA_FN_ACC ALPAKA_FN_INLINE float computePT5RPhiChiSquaredInwards(
413  float g, float f, float r, float* xPix, float* yPix) {
414  /*
415  Using the computed regression center and radius, compute the chi squared for the pixels
416  */
417 
418  float chiSquared = 0;
419  for (size_t i = 0; i < 2; i++) {
420  float residual = (xPix[i] - g) * (xPix[i] - g) + (yPix[i] - f) * (yPix[i] - f) - r * r;
421  chiSquared += residual * residual;
422  }
423  chiSquared *= 0.5f;
424  return chiSquared;
425  }
426 
427  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passPT5RPhiChiSquaredInwardsCuts(ModulesConst modules,
428  uint16_t lowerModuleIndex1,
429  uint16_t lowerModuleIndex2,
430  uint16_t lowerModuleIndex3,
431  uint16_t lowerModuleIndex4,
432  uint16_t lowerModuleIndex5,
433  float rPhiChiSquared) {
434  const int layer1 =
435  modules.layers()[lowerModuleIndex1] + 6 * (modules.subdets()[lowerModuleIndex1] == Endcap) +
436  5 * (modules.subdets()[lowerModuleIndex1] == Endcap and modules.moduleType()[lowerModuleIndex1] == TwoS);
437  const int layer2 =
438  modules.layers()[lowerModuleIndex2] + 6 * (modules.subdets()[lowerModuleIndex2] == Endcap) +
439  5 * (modules.subdets()[lowerModuleIndex2] == Endcap and modules.moduleType()[lowerModuleIndex2] == TwoS);
440  const int layer3 =
441  modules.layers()[lowerModuleIndex3] + 6 * (modules.subdets()[lowerModuleIndex3] == Endcap) +
442  5 * (modules.subdets()[lowerModuleIndex3] == Endcap and modules.moduleType()[lowerModuleIndex3] == TwoS);
443  const int layer4 =
444  modules.layers()[lowerModuleIndex4] + 6 * (modules.subdets()[lowerModuleIndex4] == Endcap) +
445  5 * (modules.subdets()[lowerModuleIndex4] == Endcap and modules.moduleType()[lowerModuleIndex4] == TwoS);
446  const int layer5 =
447  modules.layers()[lowerModuleIndex5] + 6 * (modules.subdets()[lowerModuleIndex5] == Endcap) +
448  5 * (modules.subdets()[lowerModuleIndex5] == Endcap and modules.moduleType()[lowerModuleIndex5] == TwoS);
449 
450  if (layer1 == 1 and layer2 == 2 and layer3 == 3) {
451  if (layer4 == 12 and layer5 == 13) {
452  return rPhiChiSquared < 451.141f;
453  } else if (layer4 == 4 and layer5 == 12) {
454  return rPhiChiSquared < 786.173f;
455  } else if (layer4 == 4 and layer5 == 5) {
456  return rPhiChiSquared < 595.545f;
457  } else if (layer4 == 7 and layer5 == 13) {
458  return rPhiChiSquared < 581.339f;
459  } else if (layer4 == 7 and layer5 == 8) {
460  return rPhiChiSquared < 112.537f;
461  }
462  } else if (layer1 == 1 and layer2 == 2 and layer3 == 7) {
463  if (layer4 == 13 and layer5 == 14) {
464  return rPhiChiSquared < 225.322f;
465  } else if (layer4 == 8 and layer5 == 14) {
466  return rPhiChiSquared < 1192.402f;
467  } else if (layer4 == 8 and layer5 == 9) {
468  return rPhiChiSquared < 786.173f;
469  }
470  } else if (layer1 == 1 and layer2 == 7 and layer3 == 8) {
471  if (layer4 == 9 and layer5 == 10) {
472  return rPhiChiSquared < 1037.817f;
473  } else if (layer4 == 9 and layer5 == 15) {
474  return rPhiChiSquared < 1808.536f;
475  }
476  } else if (layer1 == 2 and layer2 == 3 and layer3 == 4) {
477  if (layer4 == 12 and layer5 == 13) {
478  return rPhiChiSquared < 684.253f;
479  } else if (layer4 == 5 and layer5 == 12) {
480  return rPhiChiSquared < 684.253f;
481  } else if (layer4 == 5 and layer5 == 6) {
482  return rPhiChiSquared < 684.253f;
483  }
484  } else if (layer1 == 2 and layer2 == 3 and layer3 == 7) {
485  if (layer4 == 13 and layer5 == 14) {
486  return rPhiChiSquared < 451.141f;
487  } else if (layer4 == 8 and layer5 == 14) {
488  return rPhiChiSquared < 518.34f;
489  }
490  } else if (layer1 == 2 and layer2 == 7 and layer3 == 8) {
491  if (layer4 == 14 and layer5 == 15) {
492  return rPhiChiSquared < 2077.92f;
493  } else if (layer4 == 9 and layer5 == 10) {
494  return rPhiChiSquared < 74.20f;
495  } else if (layer4 == 9 and layer5 == 15) {
496  return rPhiChiSquared < 1808.536f;
497  }
498  } else if (layer1 == 3 and layer2 == 7 and layer3 == 8 and layer4 == 14 and layer5 == 15) {
499  return rPhiChiSquared < 786.173f;
500  } else if (layer1 == 7 and layer2 == 8 and layer3 == 9) {
501  if (layer4 == 10 and layer5 == 11) {
502  return rPhiChiSquared < 1574.076f;
503  } else if (layer4 == 10 and layer5 == 16) {
504  return rPhiChiSquared < 5492.11f;
505  } else if (layer4 == 15 and layer5 == 16) {
506  return rPhiChiSquared < 2743.037f;
507  }
508  }
509  return true;
510  }
511 
512  template <typename TAcc>
513  ALPAKA_FN_ACC ALPAKA_FN_INLINE float computePT5RZChiSquared(TAcc const& acc,
515  uint16_t* lowerModuleIndices,
516  float* rtPix,
517  float* zPix,
518  float* rts,
519  float* zs) {
520  //use the two anchor hits of the pixel segment to compute the slope
521  //then compute the pseudo chi squared of the five outer hits
522 
523  float slope = (zPix[1] - zPix[0]) / (rtPix[1] - rtPix[0]);
524  float residual = 0;
525  float error2 = 0;
526  //hardcoded array indices!!!
527  float RMSE = 0;
528  for (size_t i = 0; i < Params_T5::kLayers; i++) {
529  uint16_t& lowerModuleIndex = lowerModuleIndices[i];
530  const int moduleType = modules.moduleType()[lowerModuleIndex];
531  const int moduleSide = modules.sides()[lowerModuleIndex];
532  const int moduleSubdet = modules.subdets()[lowerModuleIndex];
533 
534  residual = (moduleSubdet == Barrel) ? (zs[i] - zPix[0]) - slope * (rts[i] - rtPix[0])
535  : (rts[i] - rtPix[0]) - (zs[i] - zPix[0]) / slope;
536  const float& drdz = modules.drdzs()[lowerModuleIndex];
537  //PS Modules
538  if (moduleType == 0) {
539  error2 = kPixelPSZpitch * kPixelPSZpitch;
540  } else //2S modules
541  {
542  error2 = kStrip2SZpitch * kStrip2SZpitch;
543  }
544 
545  //special dispensation to tilted PS modules!
546  if (moduleType == 0 and moduleSubdet == Barrel and moduleSide != Center) {
547  error2 /= (1.f + drdz * drdz);
548  }
549  RMSE += (residual * residual) / error2;
550  }
551 
552  RMSE = alpaka::math::sqrt(acc, 0.2f * RMSE); // Divided by the degree of freedom 5.
553  return RMSE;
554  }
555 
556  template <typename TAcc>
557  ALPAKA_FN_ACC ALPAKA_FN_INLINE bool runPixelQuintupletDefaultAlgo(TAcc const& acc,
560  MiniDoubletsConst mds,
561  SegmentsConst segments,
562  SegmentsPixelConst segmentsPixel,
563  TripletsConst triplets,
564  QuintupletsConst quintuplets,
565  unsigned int pixelSegmentIndex,
566  unsigned int quintupletIndex,
567  float& rzChiSquared,
568  float& rPhiChiSquared,
569  float& rPhiChiSquaredInwards,
570  float& pixelRadius,
571  float& quintupletRadius,
572  float& centerX,
573  float& centerY,
574  unsigned int pixelSegmentArrayIndex) {
575  unsigned int t5InnerT3Index = quintuplets.tripletIndices()[quintupletIndex][0];
576  unsigned int t5OuterT3Index = quintuplets.tripletIndices()[quintupletIndex][1];
577 
578  float pixelRadiusTemp, tripletRadius, rPhiChiSquaredTemp, rzChiSquaredTemp, rPhiChiSquaredInwardsTemp, centerXTemp,
579  centerYTemp;
580 
581  if (not runPixelTripletDefaultAlgo(acc,
582  modules,
583  ranges,
584  mds,
585  segments,
586  segmentsPixel,
587  triplets,
588  pixelSegmentIndex,
589  t5InnerT3Index,
590  pixelRadiusTemp,
591  tripletRadius,
592  centerXTemp,
593  centerYTemp,
594  rzChiSquaredTemp,
595  rPhiChiSquaredTemp,
596  rPhiChiSquaredInwardsTemp,
597  false))
598  return false;
599 
600  unsigned int firstSegmentIndex = triplets.segmentIndices()[t5InnerT3Index][0];
601  unsigned int secondSegmentIndex = triplets.segmentIndices()[t5InnerT3Index][1];
602  unsigned int thirdSegmentIndex = triplets.segmentIndices()[t5OuterT3Index][0];
603  unsigned int fourthSegmentIndex = triplets.segmentIndices()[t5OuterT3Index][1];
604 
605  unsigned int pixelInnerMDIndex = segments.mdIndices()[pixelSegmentIndex][0];
606  unsigned int pixelOuterMDIndex = segments.mdIndices()[pixelSegmentIndex][1];
607  unsigned int firstMDIndex = segments.mdIndices()[firstSegmentIndex][0];
608  unsigned int secondMDIndex = segments.mdIndices()[secondSegmentIndex][0];
609  unsigned int thirdMDIndex = segments.mdIndices()[secondSegmentIndex][1];
610  unsigned int fourthMDIndex = segments.mdIndices()[thirdSegmentIndex][1];
611  unsigned int fifthMDIndex = segments.mdIndices()[fourthSegmentIndex][1];
612 
613  uint16_t lowerModuleIndex1 = quintuplets.lowerModuleIndices()[quintupletIndex][0];
614  uint16_t lowerModuleIndex2 = quintuplets.lowerModuleIndices()[quintupletIndex][1];
615  uint16_t lowerModuleIndex3 = quintuplets.lowerModuleIndices()[quintupletIndex][2];
616  uint16_t lowerModuleIndex4 = quintuplets.lowerModuleIndices()[quintupletIndex][3];
617  uint16_t lowerModuleIndex5 = quintuplets.lowerModuleIndices()[quintupletIndex][4];
618 
619  uint16_t lowerModuleIndices[Params_T5::kLayers] = {
620  lowerModuleIndex1, lowerModuleIndex2, lowerModuleIndex3, lowerModuleIndex4, lowerModuleIndex5};
621 
622  float zPix[Params_pLS::kLayers] = {mds.anchorZ()[pixelInnerMDIndex], mds.anchorZ()[pixelOuterMDIndex]};
623  float rtPix[Params_pLS::kLayers] = {mds.anchorRt()[pixelInnerMDIndex], mds.anchorRt()[pixelOuterMDIndex]};
624  float zs[Params_T5::kLayers] = {mds.anchorZ()[firstMDIndex],
625  mds.anchorZ()[secondMDIndex],
626  mds.anchorZ()[thirdMDIndex],
627  mds.anchorZ()[fourthMDIndex],
628  mds.anchorZ()[fifthMDIndex]};
629  float rts[Params_T5::kLayers] = {mds.anchorRt()[firstMDIndex],
630  mds.anchorRt()[secondMDIndex],
631  mds.anchorRt()[thirdMDIndex],
632  mds.anchorRt()[fourthMDIndex],
633  mds.anchorRt()[fifthMDIndex]};
634 
635  rzChiSquared = computePT5RZChiSquared(acc, modules, lowerModuleIndices, rtPix, zPix, rts, zs);
636 
637  if (/*pixelRadius*/ 0 < 5.0f * kR1GeVf) { // FIXME: pixelRadius is not defined yet
639  lowerModuleIndex1,
640  lowerModuleIndex2,
641  lowerModuleIndex3,
642  lowerModuleIndex4,
643  lowerModuleIndex5,
644  rzChiSquared))
645  return false;
646  }
647 
648  //outer T5
649  float xs[Params_T5::kLayers] = {mds.anchorX()[firstMDIndex],
650  mds.anchorX()[secondMDIndex],
651  mds.anchorX()[thirdMDIndex],
652  mds.anchorX()[fourthMDIndex],
653  mds.anchorX()[fifthMDIndex]};
654  float ys[Params_T5::kLayers] = {mds.anchorY()[firstMDIndex],
655  mds.anchorY()[secondMDIndex],
656  mds.anchorY()[thirdMDIndex],
657  mds.anchorY()[fourthMDIndex],
658  mds.anchorY()[fifthMDIndex]};
659 
660  //get the appropriate radii and centers
661  centerX = segmentsPixel.circleCenterX()[pixelSegmentArrayIndex];
662  centerY = segmentsPixel.circleCenterY()[pixelSegmentArrayIndex];
663  pixelRadius = segmentsPixel.circleRadius()[pixelSegmentArrayIndex];
664 
665  float T5CenterX = quintuplets.regressionG()[quintupletIndex];
666  float T5CenterY = quintuplets.regressionF()[quintupletIndex];
667  quintupletRadius = quintuplets.regressionRadius()[quintupletIndex];
668 
669  rPhiChiSquared = computePT5RPhiChiSquared(acc, modules, lowerModuleIndices, centerX, centerY, pixelRadius, xs, ys);
670 
671  if (pixelRadius < 5.0f * kR1GeVf) {
673  lowerModuleIndex1,
674  lowerModuleIndex2,
675  lowerModuleIndex3,
676  lowerModuleIndex4,
677  lowerModuleIndex5,
678  rPhiChiSquared))
679  return false;
680  }
681 
682  float xPix[] = {mds.anchorX()[pixelInnerMDIndex], mds.anchorX()[pixelOuterMDIndex]};
683  float yPix[] = {mds.anchorY()[pixelInnerMDIndex], mds.anchorY()[pixelOuterMDIndex]};
684  rPhiChiSquaredInwards = computePT5RPhiChiSquaredInwards(T5CenterX, T5CenterY, quintupletRadius, xPix, yPix);
685 
686  if (quintuplets.regressionRadius()[quintupletIndex] < 5.0f * kR1GeVf) {
688  lowerModuleIndex1,
689  lowerModuleIndex2,
690  lowerModuleIndex3,
691  lowerModuleIndex4,
692  lowerModuleIndex5,
693  rPhiChiSquaredInwards))
694  return false;
695  }
696  //trusting the T5 regression center to also be a good estimate..
697  centerX = (centerX + T5CenterX) / 2;
698  centerY = (centerY + T5CenterY) / 2;
699 
700  return true;
701  }
702 
704  template <typename TAcc>
705  ALPAKA_FN_ACC void operator()(TAcc const& acc,
707  ModulesPixelConst modulesPixel,
708  MiniDoubletsConst mds,
709  SegmentsConst segments,
710  SegmentsPixel segmentsPixel,
711  Triplets triplets,
712  Quintuplets quintuplets,
713  QuintupletsOccupancyConst quintupletsOccupancy,
714  PixelQuintuplets pixelQuintuplets,
715  unsigned int* connectedPixelSize,
716  unsigned int* connectedPixelIndex,
717  unsigned int nPixelSegments,
718  ObjectRangesConst ranges) const {
719  auto const globalBlockIdx = alpaka::getIdx<alpaka::Grid, alpaka::Blocks>(acc);
720  auto const globalThreadIdx = alpaka::getIdx<alpaka::Grid, alpaka::Threads>(acc);
721  auto const gridBlockExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Blocks>(acc);
722  auto const gridThreadExtent = alpaka::getWorkDiv<alpaka::Grid, alpaka::Threads>(acc);
723 
724  for (unsigned int i_pLS = globalThreadIdx[1]; i_pLS < nPixelSegments; i_pLS += gridThreadExtent[1]) {
725  auto iLSModule_max = connectedPixelIndex[i_pLS] + connectedPixelSize[i_pLS];
726  for (unsigned int iLSModule = connectedPixelIndex[i_pLS] + globalBlockIdx[0]; iLSModule < iLSModule_max;
727  iLSModule += gridBlockExtent[0]) {
728  //these are actual module indices
729  uint16_t quintupletLowerModuleIndex = modulesPixel.connectedPixels()[iLSModule];
730  if (quintupletLowerModuleIndex >= modules.nLowerModules())
731  continue;
732  if (modules.moduleType()[quintupletLowerModuleIndex] == TwoS)
733  continue;
734  uint16_t pixelModuleIndex = modules.nLowerModules();
735  if (segmentsPixel.isDup()[i_pLS])
736  continue;
737  unsigned int nOuterQuintuplets = quintupletsOccupancy.nQuintuplets()[quintupletLowerModuleIndex];
738 
739  if (nOuterQuintuplets == 0)
740  continue;
741 
742  unsigned int pixelSegmentIndex = ranges.segmentModuleIndices()[pixelModuleIndex] + i_pLS;
743 
744  //fetch the quintuplet
745  for (unsigned int outerQuintupletArrayIndex = globalThreadIdx[2];
746  outerQuintupletArrayIndex < nOuterQuintuplets;
747  outerQuintupletArrayIndex += gridThreadExtent[2]) {
748  unsigned int quintupletIndex =
749  ranges.quintupletModuleIndices()[quintupletLowerModuleIndex] + outerQuintupletArrayIndex;
750 
751  if (quintuplets.isDup()[quintupletIndex])
752  continue;
753 
754  float rzChiSquared, rPhiChiSquared, rPhiChiSquaredInwards, pixelRadius, quintupletRadius, centerX, centerY;
755 
757  modules,
758  ranges,
759  mds,
760  segments,
761  segmentsPixel,
762  triplets,
763  quintuplets,
764  pixelSegmentIndex,
765  quintupletIndex,
766  rzChiSquared,
767  rPhiChiSquared,
768  rPhiChiSquaredInwards,
769  pixelRadius,
770  quintupletRadius,
771  centerX,
772  centerY,
773  static_cast<unsigned int>(i_pLS));
774  if (success) {
775  unsigned int totOccupancyPixelQuintuplets = alpaka::atomicAdd(
776  acc, &pixelQuintuplets.totOccupancyPixelQuintuplets(), 1u, alpaka::hierarchy::Threads{});
777  if (totOccupancyPixelQuintuplets >= n_max_pixel_quintuplets) {
778 #ifdef WARNINGS
779  printf("Pixel Quintuplet excess alert!\n");
780 #endif
781  } else {
782  unsigned int pixelQuintupletIndex =
783  alpaka::atomicAdd(acc, &pixelQuintuplets.nPixelQuintuplets(), 1u, alpaka::hierarchy::Threads{});
784  float eta = __H2F(quintuplets.eta()[quintupletIndex]);
785  float phi = __H2F(quintuplets.phi()[quintupletIndex]);
786 
788  mds,
789  segments,
790  quintuplets,
791  pixelQuintuplets,
792  pixelSegmentIndex,
793  quintupletIndex,
794  pixelQuintupletIndex,
795  rzChiSquared,
796  rPhiChiSquared,
797  rPhiChiSquaredInwards,
798  rPhiChiSquared,
799  eta,
800  phi,
801  pixelRadius,
802  quintupletRadius,
803  centerX,
804  centerY);
805 
806  triplets.partOfPT5()[quintuplets.tripletIndices()[quintupletIndex][0]] = true;
807  triplets.partOfPT5()[quintuplets.tripletIndices()[quintupletIndex][1]] = true;
808  segmentsPixel.partOfPT5()[i_pLS] = true;
809  quintuplets.partOfPT5()[quintupletIndex] = true;
810  } // tot occupancy
811  } // end success
812  } // end T5
813  } // end iLS
814  } // end i_pLS
815  }
816  };
817 } // namespace ALPAKA_ACCELERATOR_NAMESPACE::lst
818 #endif
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kR1GeVf
Definition: Common.h:50
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passPT5RPhiChiSquaredCuts(ModulesConst modules, uint16_t lowerModuleIndex1, uint16_t lowerModuleIndex2, uint16_t lowerModuleIndex3, uint16_t lowerModuleIndex4, uint16_t lowerModuleIndex5, float rPhiChiSquared)
static const double slope[3]
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
ALPAKA_FN_ACC ALPAKA_FN_INLINE float computePT5RZChiSquared(TAcc const &acc, ModulesConst modules, uint16_t *lowerModuleIndices, float *rtPix, float *zPix, float *rts, float *zs)
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kWidthPS
Definition: Common.h:58
ALPAKA_FN_ACC ALPAKA_FN_INLINE void computeSigmasForRegression_pT5(TAcc const &acc, ModulesConst modules, const uint16_t *lowerModuleIndices, float *delta1, float *delta2, float *slopes, bool *isFlat, unsigned int nPoints=5, bool anchorHits=true)
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kWidth2S
Definition: Common.h:57
ALPAKA_FN_ACC ALPAKA_FN_INLINE float computePT5RPhiChiSquaredInwards(float g, float f, float r, float *xPix, float *yPix)
ALPAKA_FN_ACC ALPAKA_FN_INLINE float computePT5RPhiChiSquared(TAcc const &acc, ModulesConst modules, uint16_t *lowerModuleIndices, float g, float f, float radius, float *xs, float *ys)
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kVerticalModuleSlope
Definition: Common.h:61
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passPT5RZChiSquaredCuts(ModulesConst modules, uint16_t lowerModuleIndex1, uint16_t lowerModuleIndex2, uint16_t lowerModuleIndex3, uint16_t lowerModuleIndex4, uint16_t lowerModuleIndex5, float rzChiSquared)
#define __H2F
Definition: Common.h:50
#define __F2H
Definition: Common.h:49
PixelQuintupletsSoA::View PixelQuintuplets
QuintupletsSoA::View Quintuplets
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool runPixelTripletDefaultAlgo(TAcc const &acc, ModulesConst modules, ObjectRangesConst ranges, MiniDoubletsConst mds, SegmentsConst segments, SegmentsPixelConst segmentsPixel, TripletsConst triplets, unsigned int pixelSegmentIndex, unsigned int tripletIndex, float &pixelRadius, float &tripletRadius, float &centerX, float &centerY, float &rzChiSquared, float &rPhiChiSquared, float &rPhiChiSquaredInwards, bool runChiSquaredCuts=true)
Definition: PixelTriplet.h:670
TripletsSoA::ConstView TripletsConst
Definition: TripletsSoA.h:31
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
ModulesPixelSoA::ConstView ModulesPixelConst
Definition: ModulesSoA.h:49
SegmentsPixelSoA::View SegmentsPixel
Definition: SegmentsSoA.h:52
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ModuleType
Definition: DQModule.h:24
double f[11][100]
MiniDoubletsSoA::ConstView MiniDoubletsConst
ModulesSoA::ConstView ModulesConst
Definition: ModulesSoA.h:47
string ranges
Definition: diffTwoXMLs.py:79
SegmentsPixelSoA::ConstView SegmentsPixelConst
Definition: SegmentsSoA.h:53
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kPixelPSZpitch
Definition: Common.h:54
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kStripPSZpitch
Definition: Common.h:55
ALPAKA_FN_ACC void operator()(TAcc const &acc, ModulesConst modules, ModulesPixelConst modulesPixel, MiniDoubletsConst mds, SegmentsConst segments, SegmentsPixel segmentsPixel, Triplets triplets, Quintuplets quintuplets, QuintupletsOccupancyConst quintupletsOccupancy, PixelQuintuplets pixelQuintuplets, unsigned int *connectedPixelSize, unsigned int *connectedPixelIndex, unsigned int nPixelSegments, ObjectRangesConst ranges) const
ALPAKA_FN_ACC ALPAKA_FN_INLINE float computeChiSquaredpT5(TAcc const &acc, unsigned int nPoints, float *xs, float *ys, float *delta1, float *delta2, float *slopes, bool *isFlat, float g, float f, float radius)
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kPi
Definition: Common.h:39
ALPAKA_STATIC_ACC_MEM_GLOBAL constexpr float kStrip2SZpitch
Definition: Common.h:56
TripletsSoA::View Triplets
Definition: TripletsSoA.h:30
QuintupletsSoA::ConstView QuintupletsConst
ALPAKA_FN_ACC ALPAKA_FN_INLINE void addPixelQuintupletToMemory(ModulesConst modules, MiniDoubletsConst mds, SegmentsConst segments, QuintupletsConst quintuplets, PixelQuintuplets pixelQuintuplets, unsigned int pixelIndex, unsigned int t5Index, unsigned int pixelQuintupletIndex, float rzChiSquared, float rPhiChiSquared, float rPhiChiSquaredInwards, float score, float eta, float phi, float pixelRadius, float quintupletRadius, float centerX, float centerY)
ObjectRangesSoA::ConstView ObjectRangesConst
SegmentsSoA::ConstView SegmentsConst
Definition: SegmentsSoA.h:49
QuintupletsOccupancySoA::ConstView QuintupletsOccupancyConst
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool runPixelQuintupletDefaultAlgo(TAcc const &acc, ModulesConst modules, ObjectRangesConst ranges, MiniDoubletsConst mds, SegmentsConst segments, SegmentsPixelConst segmentsPixel, TripletsConst triplets, QuintupletsConst quintuplets, unsigned int pixelSegmentIndex, unsigned int quintupletIndex, float &rzChiSquared, float &rPhiChiSquared, float &rPhiChiSquaredInwards, float &pixelRadius, float &quintupletRadius, float &centerX, float &centerY, unsigned int pixelSegmentArrayIndex)
constexpr unsigned int n_max_pixel_quintuplets
Definition: Common.h:36
ALPAKA_FN_ACC ALPAKA_FN_INLINE bool passPT5RPhiChiSquaredInwardsCuts(ModulesConst modules, uint16_t lowerModuleIndex1, uint16_t lowerModuleIndex2, uint16_t lowerModuleIndex3, uint16_t lowerModuleIndex4, uint16_t lowerModuleIndex5, float rPhiChiSquared)
T1 atomicAdd(T1 *a, T2 b)
Definition: cudaCompat.h:61
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE float eta(TAcc const &acc, float x, float y, float z)
Definition: Hit.h:11