CMS 3D CMS Logo

CAHitNtupletGeneratorKernelsImpl.h
Go to the documentation of this file.
1 //
2 // Original Author: Felice Pantaleo, CERN
3 //
4 
5 // #define NTUPLE_DEBUG
6 
7 #include <cmath>
8 #include <cstdint>
9 
10 #include <cuda_runtime.h>
11 
15 
16 #include "CAConstants.h"
18 #include "GPUCACell.h"
19 #include "gpuFishbone.h"
20 #include "gpuPixelDoublets.h"
21 
24 
27 
31 
32 __global__ void kernel_checkOverflows(HitContainer const *foundNtuplets,
36  GPUCACell const *__restrict__ cells,
37  uint32_t const *__restrict__ nCells,
40  GPUCACell::OuterHitOfCell const *__restrict__ isOuterHitOfCell,
41  uint32_t nHits,
42  uint32_t maxNumberOfDoublets,
44  auto first = threadIdx.x + blockIdx.x * blockDim.x;
45 
46  auto &c = *counters;
47  // counters once per event
48  if (0 == first) {
49  atomicAdd(&c.nEvents, 1);
50  atomicAdd(&c.nHits, nHits);
51  atomicAdd(&c.nCells, *nCells);
52  atomicAdd(&c.nTuples, apc->get().m);
53  atomicAdd(&c.nFitTracks, tupleMultiplicity->size());
54  }
55 
56 #ifdef NTUPLE_DEBUG
57  if (0 == first) {
58  printf("number of found cells %d, found tuples %d with total hits %d out of %d\n",
59  *nCells,
60  apc->get().m,
61  apc->get().n,
62  nHits);
63  if (apc->get().m < caConstants::maxNumberOfQuadruplets()) {
64  assert(foundNtuplets->size(apc->get().m) == 0);
65  assert(foundNtuplets->size() == apc->get().n);
66  }
67  }
68 
69  for (int idx = first, nt = foundNtuplets->nbins(); idx < nt; idx += gridDim.x * blockDim.x) {
70  if (foundNtuplets->size(idx) > 5)
71  printf("ERROR %d, %d\n", idx, foundNtuplets->size(idx));
72  assert(foundNtuplets->size(idx) < 6);
73  for (auto ih = foundNtuplets->begin(idx); ih != foundNtuplets->end(idx); ++ih)
74  assert(*ih < nHits);
75  }
76 #endif
77 
78  if (0 == first) {
80  printf("Tuples overflow\n");
82  printf("Cells overflow\n");
83  if (cellNeighbors && cellNeighbors->full())
84  printf("cellNeighbors overflow\n");
85  if (cellTracks && cellTracks->full())
86  printf("cellTracks overflow\n");
87  }
88 
89  for (int idx = first, nt = (*nCells); idx < nt; idx += gridDim.x * blockDim.x) {
90  auto const &thisCell = cells[idx];
91  if (thisCell.outerNeighbors().full()) //++tooManyNeighbors[thisCell.theLayerPairId];
92  printf("OuterNeighbors overflow %d in %d\n", idx, thisCell.layerPairId());
93  if (thisCell.tracks().full()) //++tooManyTracks[thisCell.theLayerPairId];
94  printf("Tracks overflow %d in %d\n", idx, thisCell.layerPairId());
95  if (thisCell.isKilled())
96  atomicAdd(&c.nKilledCells, 1);
97  if (thisCell.unused())
98  atomicAdd(&c.nEmptyCells, 1);
99  if (0 == hitToTuple->size(thisCell.inner_hit_id()) && 0 == hitToTuple->size(thisCell.outer_hit_id()))
100  atomicAdd(&c.nZeroTrackCells, 1);
101  }
102 
103  for (int idx = first, nt = nHits; idx < nt; idx += gridDim.x * blockDim.x) {
104  if (isOuterHitOfCell[idx].full()) // ++tooManyOuterHitOfCell;
105  printf("OuterHitOfCell overflow %d\n", idx);
106  }
107 }
108 
109 __global__ void kernel_fishboneCleaner(GPUCACell const *cells, uint32_t const *__restrict__ nCells, Quality *quality) {
110  constexpr auto bad = pixelTrack::Quality::bad;
111 
112  auto first = threadIdx.x + blockIdx.x * blockDim.x;
113  for (int idx = first, nt = (*nCells); idx < nt; idx += gridDim.x * blockDim.x) {
114  auto const &thisCell = cells[idx];
115  if (!thisCell.isKilled())
116  continue;
117 
118  for (auto it : thisCell.tracks())
119  quality[it] = bad;
120  }
121 }
122 
123 __global__ void kernel_earlyDuplicateRemover(GPUCACell const *cells,
124  uint32_t const *__restrict__ nCells,
126  Quality *quality) {
127  // constexpr auto bad = trackQuality::bad;
128  constexpr auto dup = pixelTrack::Quality::dup;
129  // constexpr auto loose = trackQuality::loose;
130 
131  assert(nCells);
132  auto first = threadIdx.x + blockIdx.x * blockDim.x;
133  for (int idx = first, nt = (*nCells); idx < nt; idx += gridDim.x * blockDim.x) {
134  auto const &thisCell = cells[idx];
135 
136  if (thisCell.tracks().size() < 2)
137  continue;
138  //if (0==thisCell.theUsed) continue;
139  // if (thisCell.theDoubletId < 0) continue;
140 
141  uint32_t maxNh = 0;
142 
143  // find maxNh
144  for (auto it : thisCell.tracks()) {
145  auto nh = foundNtuplets->size(it);
146  maxNh = std::max(nh, maxNh);
147  }
148 
149  for (auto it : thisCell.tracks()) {
150  if (foundNtuplets->size(it) != maxNh)
151  quality[it] = dup; //no race: simple assignment of the same constant
152  }
153  }
154 }
155 
156 __global__ void kernel_fastDuplicateRemover(GPUCACell const *__restrict__ cells,
157  uint32_t const *__restrict__ nCells,
158  HitContainer const *__restrict__ foundNtuplets,
159  TkSoA *__restrict__ tracks) {
160  constexpr auto bad = pixelTrack::Quality::bad;
161  constexpr auto dup = pixelTrack::Quality::dup;
163 
164  assert(nCells);
165 
166  auto first = threadIdx.x + blockIdx.x * blockDim.x;
167  for (int idx = first, nt = (*nCells); idx < nt; idx += gridDim.x * blockDim.x) {
168  auto const &thisCell = cells[idx];
169  if (thisCell.tracks().size() < 2)
170  continue;
171  // if (thisCell.theDoubletId < 0) continue;
172 
173  float mc = 10000.f;
174  uint16_t im = 60000;
175 
176  auto score = [&](auto it) {
177  return std::abs(tracks->tip(it)); // tip
178  // return tracks->chi2(it); //chi2
179  };
180 
181  // find min score
182  for (auto it : thisCell.tracks()) {
183  if (tracks->quality(it) == loose && score(it) < mc) {
184  mc = score(it);
185  im = it;
186  }
187  }
188  // mark all other duplicates
189  for (auto it : thisCell.tracks()) {
190  if (tracks->quality(it) != bad && it != im)
191  tracks->quality(it) = dup; //no race: simple assignment of the same constant
192  }
193  }
194 }
195 
196 __global__ void kernel_connect(cms::cuda::AtomicPairCounter *apc1,
197  cms::cuda::AtomicPairCounter *apc2, // just to zero them,
198  GPUCACell::Hits const *__restrict__ hhp,
199  GPUCACell *cells,
200  uint32_t const *__restrict__ nCells,
202  GPUCACell::OuterHitOfCell const *__restrict__ isOuterHitOfCell,
203  float hardCurvCut,
204  float ptmin,
205  float CAThetaCutBarrel,
206  float CAThetaCutForward,
207  float dcaCutInnerTriplet,
208  float dcaCutOuterTriplet) {
209  auto const &hh = *hhp;
210 
211  auto firstCellIndex = threadIdx.y + blockIdx.y * blockDim.y;
212  auto first = threadIdx.x;
213  auto stride = blockDim.x;
214 
215  if (0 == (firstCellIndex + first)) {
216  (*apc1) = 0;
217  (*apc2) = 0;
218  } // ready for next kernel
219 
220  for (int idx = firstCellIndex, nt = (*nCells); idx < nt; idx += gridDim.y * blockDim.y) {
221  auto cellIndex = idx;
222  auto &thisCell = cells[idx];
223  auto innerHitId = thisCell.inner_hit_id();
224  int numberOfPossibleNeighbors = isOuterHitOfCell[innerHitId].size();
225  auto vi = isOuterHitOfCell[innerHitId].data();
226 
227  auto ri = thisCell.inner_r(hh);
228  auto zi = thisCell.inner_z(hh);
229 
230  auto ro = thisCell.outer_r(hh);
231  auto zo = thisCell.outer_z(hh);
232  auto isBarrel = thisCell.inner_detIndex(hh) < caConstants::last_barrel_detIndex;
233 
234  for (int j = first; j < numberOfPossibleNeighbors; j += stride) {
235  auto otherCell = __ldg(vi + j);
236  auto &oc = cells[otherCell];
237  auto r1 = oc.inner_r(hh);
238  auto z1 = oc.inner_z(hh);
239  bool aligned = GPUCACell::areAlignedRZ(
240  r1,
241  z1,
242  ri,
243  zi,
244  ro,
245  zo,
246  ptmin,
247  isBarrel ? CAThetaCutBarrel : CAThetaCutForward); // 2.f*thetaCut); // FIXME tune cuts
248  if (aligned && thisCell.dcaCut(hh,
249  oc,
250  oc.inner_detIndex(hh) < caConstants::last_bpix1_detIndex ? dcaCutInnerTriplet
251  : dcaCutOuterTriplet,
252  hardCurvCut)) { // FIXME tune cuts
253  oc.addOuterNeighbor(cellIndex, *cellNeighbors);
254  thisCell.setUsedBit(1);
255  oc.setUsedBit(1);
256  }
257  } // loop on inner cells
258  } // loop on outer cells
259 }
260 
261 __global__ void kernel_find_ntuplets(GPUCACell::Hits const *__restrict__ hhp,
262  GPUCACell *__restrict__ cells,
263  uint32_t const *nCells,
267  Quality *__restrict__ quality,
268  unsigned int minHitsPerNtuplet) {
269  // recursive: not obvious to widen
270  auto const &hh = *hhp;
271 
272  auto first = threadIdx.x + blockIdx.x * blockDim.x;
273  for (int idx = first, nt = (*nCells); idx < nt; idx += gridDim.x * blockDim.x) {
274  auto const &thisCell = cells[idx];
275  if (thisCell.isKilled())
276  continue; // cut by earlyFishbone
277 
278  auto pid = thisCell.layerPairId();
279  auto doit = minHitsPerNtuplet > 3 ? pid < 3 : pid < 8 || pid > 12;
280  if (doit) {
282  stack.reset();
283  thisCell.find_ntuplets(hh, cells, *cellTracks, *foundNtuplets, *apc, quality, stack, minHitsPerNtuplet, pid < 3);
284  assert(stack.empty());
285  // printf("in %d found quadruplets: %d\n", cellIndex, apc->get());
286  }
287  }
288 }
289 
290 __global__ void kernel_mark_used(GPUCACell::Hits const *__restrict__ hhp,
291  GPUCACell *__restrict__ cells,
292  uint32_t const *nCells) {
293  auto first = threadIdx.x + blockIdx.x * blockDim.x;
294  for (int idx = first, nt = (*nCells); idx < nt; idx += gridDim.x * blockDim.x) {
295  auto &thisCell = cells[idx];
296  if (!thisCell.tracks().empty())
297  thisCell.setUsedBit(2);
298  }
299 }
300 
301 __global__ void kernel_countMultiplicity(HitContainer const *__restrict__ foundNtuplets,
302  Quality const *__restrict__ quality,
304  auto first = blockIdx.x * blockDim.x + threadIdx.x;
305  for (int it = first, nt = foundNtuplets->nbins(); it < nt; it += gridDim.x * blockDim.x) {
306  auto nhits = foundNtuplets->size(it);
307  if (nhits < 3)
308  continue;
310  continue;
312  if (nhits > 5)
313  printf("wrong mult %d %d\n", it, nhits);
314  assert(nhits < 8);
315  tupleMultiplicity->countDirect(nhits);
316  }
317 }
318 
319 __global__ void kernel_fillMultiplicity(HitContainer const *__restrict__ foundNtuplets,
320  Quality const *__restrict__ quality,
322  auto first = blockIdx.x * blockDim.x + threadIdx.x;
323  for (int it = first, nt = foundNtuplets->nbins(); it < nt; it += gridDim.x * blockDim.x) {
324  auto nhits = foundNtuplets->size(it);
325  if (nhits < 3)
326  continue;
328  continue;
330  if (nhits > 5)
331  printf("wrong mult %d %d\n", it, nhits);
332  assert(nhits < 8);
333  tupleMultiplicity->fillDirect(nhits, it);
334  }
335 }
336 
337 __global__ void kernel_classifyTracks(HitContainer const *__restrict__ tuples,
338  TkSoA const *__restrict__ tracks,
340  Quality *__restrict__ quality) {
341  int first = blockDim.x * blockIdx.x + threadIdx.x;
342  for (int it = first, nt = tuples->nbins(); it < nt; it += gridDim.x * blockDim.x) {
343  auto nhits = tuples->size(it);
344  if (nhits == 0)
345  break; // guard
346 
347  // if duplicate: not even fit
349  continue;
350 
352 
353  // mark doublets as bad
354  if (nhits < 3)
355  continue;
356 
357  // if the fit has any invalid parameters, mark it as bad
358  bool isNaN = false;
359  for (int i = 0; i < 5; ++i) {
360  isNaN |= std::isnan(tracks->stateAtBS.state(it)(i));
361  }
362  if (isNaN) {
363 #ifdef NTUPLE_DEBUG
364  printf("NaN in fit %d size %d chi2 %f\n", it, tuples->size(it), tracks->chi2(it));
365 #endif
366  continue;
367  }
368 
369  // compute a pT-dependent chi2 cut
370  // default parameters:
371  // - chi2MaxPt = 10 GeV
372  // - chi2Coeff = { 0.68177776, 0.74609577, -0.08035491, 0.00315399 }
373  // - chi2Scale = 30 for broken line fit, 45 for Riemann fit
374  // (see CAHitNtupletGeneratorGPU.cc)
375  float pt = std::min<float>(tracks->pt(it), cuts.chi2MaxPt);
376  float chi2Cut = cuts.chi2Scale *
377  (cuts.chi2Coeff[0] + pt * (cuts.chi2Coeff[1] + pt * (cuts.chi2Coeff[2] + pt * cuts.chi2Coeff[3])));
378  // above number were for Quads not normalized so for the time being just multiple by ndof for Quads (triplets to be understood)
379  if (3.f * tracks->chi2(it) >= chi2Cut) {
380 #ifdef NTUPLE_DEBUG
381  printf("Bad fit %d size %d pt %f eta %f chi2 %f\n",
382  it,
383  tuples->size(it),
384  tracks->pt(it),
385  tracks->eta(it),
386  3.f * tracks->chi2(it));
387 #endif
388  continue;
389  }
390 
391  // impose "region cuts" based on the fit results (phi, Tip, pt, cotan(theta)), Zip)
392  // default cuts:
393  // - for triplets: |Tip| < 0.3 cm, pT > 0.5 GeV, |Zip| < 12.0 cm
394  // - for quadruplets: |Tip| < 0.5 cm, pT > 0.3 GeV, |Zip| < 12.0 cm
395  // (see CAHitNtupletGeneratorGPU.cc)
396  auto const &region = (nhits > 3) ? cuts.quadruplet : cuts.triplet;
397  bool isOk = (std::abs(tracks->tip(it)) < region.maxTip) and (tracks->pt(it) > region.minPt) and
398  (std::abs(tracks->zip(it)) < region.maxZip);
399 
400  if (isOk)
402  }
403 }
404 
405 __global__ void kernel_doStatsForTracks(HitContainer const *__restrict__ tuples,
406  Quality const *__restrict__ quality,
408  int first = blockDim.x * blockIdx.x + threadIdx.x;
409  for (int idx = first, ntot = tuples->nbins(); idx < ntot; idx += gridDim.x * blockDim.x) {
410  if (tuples->size(idx) == 0)
411  break; //guard
413  continue;
414  atomicAdd(&(counters->nGoodTracks), 1);
415  }
416 }
417 
418 __global__ void kernel_countHitInTracks(HitContainer const *__restrict__ tuples,
419  Quality const *__restrict__ quality,
421  int first = blockDim.x * blockIdx.x + threadIdx.x;
422  for (int idx = first, ntot = tuples->nbins(); idx < ntot; idx += gridDim.x * blockDim.x) {
423  if (tuples->size(idx) == 0)
424  break; // guard
426  continue;
427  for (auto h = tuples->begin(idx); h != tuples->end(idx); ++h)
428  hitToTuple->countDirect(*h);
429  }
430 }
431 
432 __global__ void kernel_fillHitInTracks(HitContainer const *__restrict__ tuples,
433  Quality const *__restrict__ quality,
435  int first = blockDim.x * blockIdx.x + threadIdx.x;
436  for (int idx = first, ntot = tuples->nbins(); idx < ntot; idx += gridDim.x * blockDim.x) {
437  if (tuples->size(idx) == 0)
438  break; // guard
440  continue;
441  for (auto h = tuples->begin(idx); h != tuples->end(idx); ++h)
442  hitToTuple->fillDirect(*h, idx);
443  }
444 }
445 
446 __global__ void kernel_fillHitDetIndices(HitContainer const *__restrict__ tuples,
447  TrackingRecHit2DSOAView const *__restrict__ hhp,
448  HitContainer *__restrict__ hitDetIndices) {
449  int first = blockDim.x * blockIdx.x + threadIdx.x;
450  // copy offsets
451  for (int idx = first, ntot = tuples->totbins(); idx < ntot; idx += gridDim.x * blockDim.x) {
452  hitDetIndices->off[idx] = tuples->off[idx];
453  }
454  // fill hit indices
455  auto const &hh = *hhp;
456  auto nhits = hh.nHits();
457  for (int idx = first, ntot = tuples->size(); idx < ntot; idx += gridDim.x * blockDim.x) {
458  assert(tuples->bins[idx] < nhits);
459  hitDetIndices->bins[idx] = hh.detectorIndex(tuples->bins[idx]);
460  }
461 }
462 
463 __global__ void kernel_doStatsForHitInTracks(CAHitNtupletGeneratorKernelsGPU::HitToTuple const *__restrict__ hitToTuple,
465  auto &c = *counters;
466  int first = blockDim.x * blockIdx.x + threadIdx.x;
467  for (int idx = first, ntot = hitToTuple->nbins(); idx < ntot; idx += gridDim.x * blockDim.x) {
468  if (hitToTuple->size(idx) == 0)
469  continue; // SHALL NOT BE break
470  atomicAdd(&c.nUsedHits, 1);
471  if (hitToTuple->size(idx) > 1)
472  atomicAdd(&c.nDupHits, 1);
473  }
474 }
475 
476 __global__ void kernel_sharedHitCleaner(TrackingRecHit2DSOAView const *__restrict__ hhp,
477  HitContainer const *__restrict__ ptuples,
478  TkSoA const *__restrict__ ptracks,
479  Quality *__restrict__ quality,
480  uint16_t nmin,
482  constexpr auto bad = pixelTrack::Quality::bad;
483  constexpr auto dup = pixelTrack::Quality::dup;
484  // constexpr auto loose = trackQuality::loose;
485 
486  auto &hitToTuple = *phitToTuple;
487  auto const &foundNtuplets = *ptuples;
488  auto const &tracks = *ptracks;
489 
490  auto const &hh = *hhp;
491  int l1end = hh.hitsLayerStart()[1];
492 
493  int first = blockDim.x * blockIdx.x + threadIdx.x;
494  for (int idx = first, ntot = hitToTuple.nbins(); idx < ntot; idx += gridDim.x * blockDim.x) {
495  if (hitToTuple.size(idx) < 2)
496  continue;
497 
498  float mc = 10000.f;
499  uint16_t im = 60000;
500  uint32_t maxNh = 0;
501 
502  // find maxNh
503  for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
504  uint32_t nh = foundNtuplets.size(*it);
505  maxNh = std::max(nh, maxNh);
506  }
507  // kill all tracks shorter than maxHn (only triplets???)
508  for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
509  uint32_t nh = foundNtuplets.size(*it);
510 
511  //checking if shared hit is on bpix1 and if the tuple is short enough
512  if (idx < l1end and nh > nmin)
513  continue;
514 
515  if (maxNh != nh)
516  quality[*it] = dup;
517  }
518 
519  if (maxNh > 3)
520  continue;
521  // for triplets choose best tip!
522  for (auto ip = hitToTuple.begin(idx); ip != hitToTuple.end(idx); ++ip) {
523  auto const it = *ip;
524  if (quality[it] != bad && std::abs(tracks.tip(it)) < mc) {
525  mc = std::abs(tracks.tip(it));
526  im = it;
527  }
528  }
529  // mark duplicates
530  for (auto ip = hitToTuple.begin(idx); ip != hitToTuple.end(idx); ++ip) {
531  auto const it = *ip;
532  if (quality[it] != bad && it != im)
533  quality[it] = dup; //no race: simple assignment of the same constant
534  }
535  } // loop over hits
536 }
537 
538 __global__ void kernel_print_found_ntuplets(TrackingRecHit2DSOAView const *__restrict__ hhp,
539  HitContainer const *__restrict__ ptuples,
540  TkSoA const *__restrict__ ptracks,
541  Quality const *__restrict__ quality,
543  uint32_t maxPrint,
544  int iev) {
545  auto const &foundNtuplets = *ptuples;
546  auto const &tracks = *ptracks;
547  int first = blockDim.x * blockIdx.x + threadIdx.x;
548  for (int i = first, np = std::min(maxPrint, foundNtuplets.nbins()); i < np; i += blockDim.x * gridDim.x) {
549  auto nh = foundNtuplets.size(i);
550  if (nh < 3)
551  continue;
552  printf("TK: %d %d %d %f %f %f %f %f %f %f %d %d %d %d %d\n",
553  10000 * iev + i,
554  int(quality[i]),
555  nh,
556  tracks.charge(i),
557  tracks.pt(i),
558  tracks.eta(i),
559  tracks.phi(i),
560  tracks.tip(i),
561  tracks.zip(i),
562  // asinhf(fit_results[i].par(3)),
563  tracks.chi2(i),
564  *foundNtuplets.begin(i),
565  *(foundNtuplets.begin(i) + 1),
566  *(foundNtuplets.begin(i) + 2),
567  nh > 3 ? int(*(foundNtuplets.begin(i) + 3)) : -1,
568  nh > 4 ? int(*(foundNtuplets.begin(i) + 4)) : -1);
569  }
570 }
571 
572 __global__ void kernel_printCounters(cAHitNtupletGenerator::Counters const *counters) {
573  auto const &c = *counters;
574  printf(
575  "||Counters | nEvents | nHits | nCells | nTuples | nFitTacks | nGoodTracks | nUsedHits | nDupHits | "
576  "nKilledCells | "
577  "nEmptyCells | nZeroTrackCells ||\n");
578  printf("Counters Raw %lld %lld %lld %lld %lld %lld %lld %lld %lld %lld %lld\n",
579  c.nEvents,
580  c.nHits,
581  c.nCells,
582  c.nTuples,
583  c.nGoodTracks,
584  c.nFitTracks,
585  c.nUsedHits,
586  c.nDupHits,
587  c.nKilledCells,
588  c.nEmptyCells,
589  c.nZeroTrackCells);
590  printf("Counters Norm %lld || %.1f| %.1f| %.1f| %.1f| %.1f| %.1f| %.1f| %.1f| %.3f| %.3f||\n",
591  c.nEvents,
592  c.nHits / double(c.nEvents),
593  c.nCells / double(c.nEvents),
594  c.nTuples / double(c.nEvents),
595  c.nFitTracks / double(c.nEvents),
596  c.nGoodTracks / double(c.nEvents),
597  c.nUsedHits / double(c.nEvents),
598  c.nDupHits / double(c.nEvents),
599  c.nKilledCells / double(c.nEvents),
600  c.nEmptyCells / double(c.nCells),
601  c.nZeroTrackCells / double(c.nCells));
602 }
iev
const HitContainer *__restrict__ const TkSoA *__restrict__ const Quality *__restrict__ const CAHitNtupletGeneratorKernelsGPU::HitToTuple *__restrict__ uint32_t int iev
Definition: CAHitNtupletGeneratorKernelsImpl.h:544
pixelCPEforGPU.h
cellNeighbors
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple cms::cuda::AtomicPairCounter const GPUCACell *__restrict__ const uint32_t *__restrict__ const gpuPixelDoublets::CellNeighborsVector * cellNeighbors
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
mps_fire.i
i
Definition: mps_fire.py:428
cuts
const TkSoA *__restrict__ CAHitNtupletGeneratorKernelsGPU::QualityCuts cuts
Definition: CAHitNtupletGeneratorKernelsImpl.h:338
loose
constexpr auto loose
Definition: CAHitNtupletGeneratorKernelsImpl.h:162
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
CaloTowersParam_cfi.mc
mc
Definition: CaloTowersParam_cfi.py:8
nt
int nt
Definition: AMPTWrapper.h:42
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
min
T min(T a, T b)
Definition: MathUtil.h:58
pixelTrack::Quality::bad
cAHitNtupletGenerator::QualityCuts
Definition: CAHitNtupletGeneratorKernels.h:36
isOuterHitOfCell
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple cms::cuda::AtomicPairCounter const GPUCACell *__restrict__ const uint32_t *__restrict__ const gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell *__restrict__ isOuterHitOfCell
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
caConstants::last_bpix1_detIndex
constexpr uint32_t last_bpix1_detIndex
Definition: CAConstants.h:62
h
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
Definition: L1TUtmAlgorithmRcd.h:4
np
int np
Definition: AMPTWrapper.h:43
CAHitNtupletGeneratorKernels.h
TrackingRecHit2DHeterogeneous
Definition: TrackingRecHit2DHeterogeneous.h:8
cms::cuda::HistoContainer::m
return c m
Definition: HistoContainer.h:242
cells
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple cms::cuda::AtomicPairCounter const GPUCACell *__restrict__ cells
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
gpuPixelDoublets.h
gpuPixelDoublets::ntot
__shared__ uint32_t ntot
Definition: gpuPixelDoubletsAlgos.h:67
dup
constexpr auto dup
Definition: CAHitNtupletGeneratorKernelsImpl.h:161
TrackingRecHit2DSOAView
Definition: TrackingRecHit2DSOAView.h:15
full
Definition: GenABIO.cc:168
ptracks
const HitContainer *__restrict__ const TkSoA *__restrict__ ptracks
Definition: CAHitNtupletGeneratorKernelsImpl.h:477
caConstants::maxNumberOfQuadruplets
constexpr uint32_t maxNumberOfQuadruplets
Definition: CAConstants.h:41
nCells
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple cms::cuda::AtomicPairCounter const GPUCACell *__restrict__ const uint32_t *__restrict__ nCells
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
TrackingRecHit2DCUDA
TrackingRecHit2DHeterogeneous< cms::cudacompat::GPUTraits > TrackingRecHit2DCUDA
Definition: TrackingRecHit2DHeterogeneous.h:137
gpuFishbone.h
hhp
const TrackingRecHit2DSOAView *__restrict__ hhp
Definition: CAHitNtupletGeneratorKernelsImpl.h:447
pixelTrack::Quality::dup
caConstants::last_barrel_detIndex
constexpr uint32_t last_barrel_detIndex
Definition: CAConstants.h:63
__global__
#define __global__
Definition: cudaCompat.h:19
cms::cuda::SimpleVector
Definition: SimpleVector.h:15
heavyIonCSV_trainingSettings.idx
idx
Definition: heavyIonCSV_trainingSettings.py:5
quality
const uint32_t *__restrict__ Quality * quality
Definition: CAHitNtupletGeneratorKernelsImpl.h:109
pixelTrack::Quality
Quality
Definition: TrackSoAHeterogeneousT.h:10
h
counters
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple cms::cuda::AtomicPairCounter const GPUCACell *__restrict__ const uint32_t *__restrict__ const gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell *__restrict__ uint32_t uint32_t CAHitNtupletGeneratorKernelsGPU::Counters * counters
Definition: CAHitNtupletGeneratorKernelsImpl.h:43
caConstants::TupleMultiplicity
cms::cuda::OneToManyAssoc< tindex_type, 8, maxTuples > TupleMultiplicity
Definition: CAConstants.h:79
cms::cuda::nh
uint32_t nh
Definition: HistoContainer.h:23
CommonMethods.isnan
def isnan(num)
Definition: CommonMethods.py:98
cms::cudacompat::atomicAdd
T1 atomicAdd(T1 *a, T2 b)
Definition: cudaCompat.h:51
first
auto first
Definition: CAHitNtupletGeneratorKernelsImpl.h:112
cms::cudacompat::gridDim
const dim3 gridDim
Definition: cudaCompat.h:33
TrackSoAHeterogeneousT
Definition: TrackSoAHeterogeneousT.h:14
PixelPluginsPhase0_cfi.isBarrel
isBarrel
Definition: PixelPluginsPhase0_cfi.py:17
svgfig.stack
stack
Definition: svgfig.py:559
hh
const auto & hh
Definition: CAHitNtupletGeneratorKernelsImpl.h:455
nmin
const HitContainer *__restrict__ const TkSoA *__restrict__ Quality *__restrict__ uint16_t nmin
Definition: CAHitNtupletGeneratorKernelsImpl.h:477
cms::cuda::AtomicPairCounter
Definition: AtomicPairCounter.h:11
l1end
int l1end
Definition: CAHitNtupletGeneratorKernelsImpl.h:491
cms::cudacompat::blockDim
const dim3 blockDim
Definition: cudaCompat.h:30
tracks
const uint32_t *__restrict__ const HitContainer *__restrict__ TkSoA *__restrict__ tracks
Definition: CAHitNtupletGeneratorKernelsImpl.h:159
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
tupleMultiplicity
const caConstants::TupleMultiplicity * tupleMultiplicity
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
HLT_FULL_cff.region
region
Definition: HLT_FULL_cff.py:88271
createfilelist.int
int
Definition: createfilelist.py:10
cms::cuda::VecArray
Definition: VecArray.h:14
foundNtuplets
const uint32_t *__restrict__ HitContainer * foundNtuplets
Definition: CAHitNtupletGeneratorKernelsImpl.h:124
cudaCheck.h
nhits
auto nhits
Definition: CAHitNtupletGeneratorKernelsImpl.h:456
cms::cudacompat::threadIdx
const dim3 threadIdx
Definition: cudaCompat.h:29
maxPrint
const HitContainer *__restrict__ const TkSoA *__restrict__ const Quality *__restrict__ const CAHitNtupletGeneratorKernelsGPU::HitToTuple *__restrict__ uint32_t maxPrint
Definition: CAHitNtupletGeneratorKernelsImpl.h:539
cAHitNtupletGenerator::Counters
Definition: CAHitNtupletGeneratorKernels.h:12
cms::cuda::HistoContainer::nbins
static constexpr uint32_t nbins()
Definition: HistoContainer.h:175
diffTwoXMLs.r1
r1
Definition: diffTwoXMLs.py:53
apc
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple cms::cuda::AtomicPairCounter * apc
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
GPUCACell.h
phitToTuple
const HitContainer *__restrict__ const TkSoA *__restrict__ Quality *__restrict__ uint16_t const CAHitNtupletGeneratorKernelsGPU::HitToTuple *__restrict__ phitToTuple
Definition: CAHitNtupletGeneratorKernelsImpl.h:481
ptmin
double ptmin
Definition: HydjetWrapper.h:84
maxNumberOfDoublets
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple cms::cuda::AtomicPairCounter const GPUCACell *__restrict__ const uint32_t *__restrict__ const gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell *__restrict__ uint32_t uint32_t maxNumberOfDoublets
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
ptuples
const HitContainer *__restrict__ ptuples
Definition: CAHitNtupletGeneratorKernelsImpl.h:477
hitToTuple
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple * hitToTuple
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
hitDetIndices
const TrackingRecHit2DSOAView *__restrict__ HitContainer *__restrict__ hitDetIndices
Definition: CAHitNtupletGeneratorKernelsImpl.h:448
cellTracks
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple cms::cuda::AtomicPairCounter const GPUCACell *__restrict__ const uint32_t *__restrict__ const gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector * cellTracks
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
PixelTripletNoTipGenerator_cfi.chi2Cut
chi2Cut
Definition: PixelTripletNoTipGenerator_cfi.py:10
pixelTrack::HitContainer
TrackSoA::HitContainer HitContainer
Definition: TrackSoAHeterogeneousT.h:69
GPUCACell
Definition: GPUCACell.h:20
caConstants::HitToTuple
cms::cuda::OneToManyAssoc< tindex_type, pixelGPUConstants::maxNumberOfHits, 4 *maxTuples > HitToTuple
Definition: CAConstants.h:78
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
pixelTrack::TrackSoA
TrackSoAHeterogeneousT< maxNumber()> TrackSoA
Definition: TrackSoAHeterogeneousT.h:67
cms::cuda::HistoContainer
Definition: HistoContainer.h:152
offlineSlimmedPrimaryVertices_cfi.score
score
Definition: offlineSlimmedPrimaryVertices_cfi.py:6
cuda_assert.h
c
auto & c
Definition: CAHitNtupletGeneratorKernelsImpl.h:46
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
cms::cudacompat::__ldg
T __ldg(T const *x)
Definition: cudaCompat.h:82
pixelTrack::Quality::loose
CAConstants.h
assert
assert(nCells)
cms::cuda::HistoContainer::n
__host__ __device__ const index_type uint32_t n
Definition: HistoContainer.h:235
gpuPixelDoublets::stride
auto stride
Definition: gpuPixelDoubletsAlgos.h:80
nHits
const caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple cms::cuda::AtomicPairCounter const GPUCACell *__restrict__ const uint32_t *__restrict__ const gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell *__restrict__ uint32_t nHits
Definition: CAHitNtupletGeneratorKernelsImpl.h:33
cms::cudacompat::blockIdx
const dim3 blockIdx
Definition: cudaCompat.h:32