CMS 3D CMS Logo

CAHitNtupletGeneratorKernelsImpl.h
Go to the documentation of this file.
1 #ifndef RecoTracker_PixelSeeding_plugins_alpaka_CAHitNtupletGeneratorKernelsImpl_h
2 #define RecoTracker_PixelSeeding_plugins_alpaka_CAHitNtupletGeneratorKernelsImpl_h
3 
4 //#define GPU_DEBUG
5 //#define NTUPLE_DEBUG
6 
7 // C++ includes
8 #include <cmath>
9 #include <cstdint>
10 #include <cstdio>
11 #include <limits>
12 #include <type_traits>
13 
14 // Alpaka includes
15 #include <alpaka/alpaka.hpp>
16 
17 // CMSSW includes
25 
26 // local includes
27 #include "CACell.h"
29 #include "CAStructures.h"
30 
32 
35  constexpr float nSigma2 = 25.f;
36 
37  // all of these below are mostly to avoid brining around the relative namespace
38 
39  template <typename TrackerTraits>
41 
42  template <typename TrackerTraits>
44 
45  template <typename TrackerTraits>
47 
48  template <typename TrackerTraits>
50 
51  template <typename TrackerTraits>
53 
55 
56  template <typename TrackerTraits>
58 
59  template <typename TrackerTraits>
61 
62  template <typename TrackerTraits>
64 
65  template <typename TrackerTraits>
67 
68  template <typename TrackerTraits>
70 
72 
73  template <typename TrackerTraits>
75  public:
76  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
77  ALPAKA_FN_ACC void operator()(TAcc const &acc,
82  CACellT<TrackerTraits> const *__restrict__ cells,
83  uint32_t const *__restrict__ nCells,
87  int32_t nHits,
88  uint32_t maxNumberOfDoublets,
89  Counters *counters) const {
90  auto &c = *counters;
91  // counters once per event
93  alpaka::atomicAdd(acc, &c.nEvents, 1ull, alpaka::hierarchy::Blocks{});
94  alpaka::atomicAdd(acc, &c.nHits, static_cast<unsigned long long>(nHits), alpaka::hierarchy::Blocks{});
95  alpaka::atomicAdd(acc, &c.nCells, static_cast<unsigned long long>(*nCells), alpaka::hierarchy::Blocks{});
97  acc, &c.nTuples, static_cast<unsigned long long>(apc->get().first), alpaka::hierarchy::Blocks{});
99  &c.nFitTracks,
100  static_cast<unsigned long long>(tupleMultiplicity->size()),
102  }
103 
104 #ifdef NTUPLE_DEBUGS
106  printf("number of found cells %d \n found tuples %d with total hits %d out of %d\n",
107  *nCells,
108  apc->get().first,
109  apc->get().second,
110  nHits);
111  if (apc->get().first < TrackerTraits::maxNumberOfQuadruplets) {
112  ALPAKA_ASSERT_ACC(tracks_view.hitIndices().size(apc->get().first) == 0);
113  ALPAKA_ASSERT_ACC(tracks_view.hitIndices().size() == apc->get().second);
114  }
115  }
116 
117  for (auto idx : cms::alpakatools::uniform_elements(acc, tracks_view.hitIndices().nOnes())) {
118  if (tracks_view.hitIndices().size(idx) > TrackerTraits::maxHitsOnTrack) // current real limit
119  printf("ERROR %d, %d\n", idx, tracks_view.hitIndices().size(idx));
120  ALPAKA_ASSERT_ACC(ftracks_view.hitIndices().size(idx) <= TrackerTraits::maxHitsOnTrack);
121  for (auto ih = tracks_view.hitIndices().begin(idx); ih != tracks_view.hitIndices().end(idx); ++ih)
122  ALPAKA_ASSERT_ACC(int(*ih) < nHits);
123  }
124 #endif
125 
127  if (apc->get().first >= TrackerTraits::maxNumberOfQuadruplets)
128  printf("Tuples overflow\n");
129  if (*nCells >= maxNumberOfDoublets)
130  printf("Cells overflow\n");
131  if (cellNeighbors && cellNeighbors->full())
132  printf("cellNeighbors overflow %d %d \n", cellNeighbors->capacity(), cellNeighbors->size());
133  if (cellTracks && cellTracks->full())
134  printf("cellTracks overflow\n");
135  if (int(hitToTuple->nOnes()) < nHits)
136  printf("ERROR hitToTuple overflow %d %d\n", hitToTuple->nOnes(), nHits);
137 #ifdef GPU_DEBUG
138  printf("size of cellNeighbors %d \n cellTracks %d \n hitToTuple %d \n",
139  cellNeighbors->size(),
140  cellTracks->size(),
141  hitToTuple->size());
142 #endif
143  }
144 
145  for (auto idx : cms::alpakatools::uniform_elements(acc, *nCells)) {
146  auto const &thisCell = cells[idx];
147  if (thisCell.hasFishbone() && !thisCell.isKilled())
148  alpaka::atomicAdd(acc, &c.nFishCells, 1ull, alpaka::hierarchy::Blocks{});
149  if (thisCell.outerNeighbors().full()) //++tooManyNeighbors[thisCell.theLayerPairId];
150  printf("OuterNeighbors overflow %d in %d\n", idx, thisCell.layerPairId());
151  if (thisCell.tracks().full()) //++tooManyTracks[thisCell.theLayerPairId];
152  printf("Tracks overflow %d in %d\n", idx, thisCell.layerPairId());
153  if (thisCell.isKilled())
154  alpaka::atomicAdd(acc, &c.nKilledCells, 1ull, alpaka::hierarchy::Blocks{});
155  if (!thisCell.unused())
156  alpaka::atomicAdd(acc, &c.nEmptyCells, 1ull, alpaka::hierarchy::Blocks{});
157  if ((0 == hitToTuple->size(thisCell.inner_hit_id())) && (0 == hitToTuple->size(thisCell.outer_hit_id())))
158  alpaka::atomicAdd(acc, &c.nZeroTrackCells, 1ull, alpaka::hierarchy::Blocks{});
159  }
160 
161  // FIXME this loop was up to nHits - isOuterHitOfCell.offset in the CUDA version
162  for (auto idx : cms::alpakatools::uniform_elements(acc, nHits))
163  if ((*isOuterHitOfCell).container[idx].full()) // ++tooManyOuterHitOfCell;
164  printf("OuterHitOfCell overflow %d\n", idx);
165  }
166  };
167 
168  template <typename TrackerTraits>
170  public:
171  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
172  ALPAKA_FN_ACC void operator()(TAcc const &acc,
174  uint32_t const *__restrict__ nCells,
177 
178  for (auto idx : cms::alpakatools::uniform_elements(acc, *nCells)) {
179  auto const &thisCell = cells[idx];
180  if (!thisCell.isKilled())
181  continue;
182 
183  for (auto it : thisCell.tracks())
185  }
186  }
187  };
188 
189  // remove shorter tracks if sharing a cell
190  // It does not seem to affect efficiency in any way!
191  template <typename TrackerTraits>
193  public:
194  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
195  ALPAKA_FN_ACC void operator()(TAcc const &acc,
197  uint32_t const *__restrict__ nCells,
199  bool dupPassThrough) const {
200  // quality to mark rejected
201  constexpr auto reject = Quality::edup;
203  for (auto idx : cms::alpakatools::uniform_elements(acc, *nCells)) {
204  auto const &thisCell = cells[idx];
205 
206  if (thisCell.tracks().size() < 2)
207  continue;
208 
209  int8_t maxNl = 0;
210 
211  // find maxNl
212  for (auto it : thisCell.tracks()) {
213  auto nl = tracks_view[it].nLayers();
214  maxNl = std::max(nl, maxNl);
215  }
216 
217  // if (maxNl<4) continue;
218  // quad pass through (leave it here for tests)
219  // maxNl = std::min(4, maxNl);
220 
221  for (auto it : thisCell.tracks()) {
222  if (tracks_view[it].nLayers() < maxNl)
223  tracks_view[it].quality() = reject; // no race: simple assignment of the same constant
224  }
225  }
226  }
227  };
228 
229  // assume the above (so, short tracks already removed)
230  template <typename TrackerTraits>
232  public:
233  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
234  ALPAKA_FN_ACC void operator()(TAcc const &acc,
235  CACellT<TrackerTraits> const *__restrict__ cells,
236  uint32_t const *__restrict__ nCells,
238  bool dupPassThrough) const {
239  // quality to mark rejected
242 
244  const auto ntNCells = (*nCells);
245 
246  for (auto idx : cms::alpakatools::uniform_elements(acc, ntNCells)) {
247  auto const &thisCell = cells[idx];
248  if (thisCell.tracks().size() < 2)
249  continue;
250 
251  float mc = maxScore;
252  uint16_t im = tkNotFound;
253 
254  auto score = [&](auto it) { return std::abs(reco::tip(tracks_view, it)); };
255 
256  // full crazy combinatorics
257  int ntr = thisCell.tracks().size();
258  for (int i = 0; i < ntr - 1; ++i) {
259  auto it = thisCell.tracks()[i];
260  auto qi = tracks_view[it].quality();
261  if (qi <= reject)
262  continue;
263  auto opi = tracks_view[it].state()(2);
264  auto e2opi = tracks_view[it].covariance()(9);
265  auto cti = tracks_view[it].state()(3);
266  auto e2cti = tracks_view[it].covariance()(12);
267  for (auto j = i + 1; j < ntr; ++j) {
268  auto jt = thisCell.tracks()[j];
269  auto qj = tracks_view[jt].quality();
270  if (qj <= reject)
271  continue;
272  auto opj = tracks_view[jt].state()(2);
273  auto ctj = tracks_view[jt].state()(3);
274  auto dct = nSigma2 * (tracks_view[jt].covariance()(12) + e2cti);
275  if ((cti - ctj) * (cti - ctj) > dct)
276  continue;
277  auto dop = nSigma2 * (tracks_view[jt].covariance()(9) + e2opi);
278  if ((opi - opj) * (opi - opj) > dop)
279  continue;
280  if ((qj < qi) || (qj == qi && score(it) < score(jt)))
281  tracks_view[jt].quality() = reject;
282  else {
283  tracks_view[it].quality() = reject;
284  break;
285  }
286  }
287  }
288 
289  // find maxQual
290  auto maxQual = reject; // no duplicate!
291  for (auto it : thisCell.tracks()) {
292  if (tracks_view[it].quality() > maxQual)
293  maxQual = tracks_view[it].quality();
294  }
295 
296  if (maxQual <= loose)
297  continue;
298 
299  // find min score
300  for (auto it : thisCell.tracks()) {
301  if (tracks_view[it].quality() == maxQual && score(it) < mc) {
302  mc = score(it);
303  im = it;
304  }
305  }
306 
307  if (tkNotFound == im)
308  continue;
309 
310  // mark all other duplicates (not yet, keep it loose)
311  for (auto it : thisCell.tracks()) {
312  if (tracks_view[it].quality() > loose && it != im)
313  tracks_view[it].quality() = loose; //no race: simple assignment of the same constant
314  }
315  }
316  }
317  };
318 
319  template <typename TrackerTraits>
321  public:
322  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
323  ALPAKA_FN_ACC void operator()(TAcc const &acc,
325  cms::alpakatools::AtomicPairCounter *apc2, // just to zero them
328  uint32_t *nCells,
333 
335  *apc1 = 0;
336  *apc2 = 0;
337  } // ready for next kernel
338 
339  // loop on outer cells
340  for (uint32_t idx : cms::alpakatools::uniform_elements_y(acc, *nCells)) {
341  auto cellIndex = idx;
342  auto &thisCell = cells[idx];
343  auto innerHitId = thisCell.inner_hit_id();
344  if (int(innerHitId) < isOuterHitOfCell->offset)
345  continue;
346 
347  uint32_t numberOfPossibleNeighbors = (*isOuterHitOfCell)[innerHitId].size();
348  auto vi = (*isOuterHitOfCell)[innerHitId].data();
349  auto ri = thisCell.inner_r(hh);
350  auto zi = thisCell.inner_z(hh);
351  auto ro = thisCell.outer_r(hh);
352  auto zo = thisCell.outer_z(hh);
353  auto isBarrel = thisCell.inner_detIndex(hh) < TrackerTraits::last_barrel_detIndex;
354 
355  // loop on inner cells
356  for (uint32_t j : cms::alpakatools::independent_group_elements_x(acc, numberOfPossibleNeighbors)) {
357  auto otherCell = (vi[j]);
358  auto &oc = cells[otherCell];
359  auto r1 = oc.inner_r(hh);
360  auto z1 = oc.inner_z(hh);
361  bool aligned = Cell::areAlignedRZ(
362  r1,
363  z1,
364  ri,
365  zi,
366  ro,
367  zo,
368  params.ptmin_,
369  isBarrel ? params.CAThetaCutBarrel_ : params.CAThetaCutForward_); // 2.f*thetaCut); // FIXME tune cuts
370  if (aligned &&
371  thisCell.dcaCut(hh,
372  oc,
373  oc.inner_detIndex(hh) < TrackerTraits::last_bpix1_detIndex ? params.dcaCutInnerTriplet_
374  : params.dcaCutOuterTriplet_,
375  params.hardCurvCut_)) { // FIXME tune cuts
376  oc.addOuterNeighbor(acc, cellIndex, *cellNeighbors);
377  thisCell.setStatusBits(Cell::StatusBit::kUsed);
378  oc.setStatusBits(Cell::StatusBit::kUsed);
379  }
380  } // loop on inner cells
381  } // loop on outer cells
382  }
383  };
384  template <typename TrackerTraits>
386  public:
387  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
388  ALPAKA_FN_ACC void operator()(TAcc const &acc,
391  CACellT<TrackerTraits> *__restrict__ cells,
392  uint32_t const *nCells,
396  // recursive: not obvious to widen
397 
399 
400 #ifdef GPU_DEBUG
402  printf("starting producing ntuplets from %d cells \n", *nCells);
403 #endif
404 
405  for (auto idx : cms::alpakatools::uniform_elements(acc, (*nCells))) {
406  auto const &thisCell = cells[idx];
407 
408  // cut by earlyFishbone
409  if (thisCell.isKilled())
410  continue;
411 
412  // we require at least three hits
413  if (thisCell.outerNeighbors().empty())
414  continue;
415 
416  auto pid = thisCell.layerPairId();
417  bool doit = params.startingLayerPair(pid);
418 
420 
421  if (doit) {
422  typename Cell::TmpTuple stack;
423  stack.reset();
424  bool bpix1Start = params.startAt0(pid);
425  thisCell.template find_ntuplets<maxDepth, TAcc>(acc,
426  hh,
427  cells,
428  *cellTracks,
429  tracks_view.hitIndices(),
430  *apc,
431  tracks_view.quality(),
432  stack,
433  params.minHitsPerNtuplet_,
434  bpix1Start);
435  ALPAKA_ASSERT_ACC(stack.empty());
436  }
437  }
438  }
439  };
440 
441  template <typename TrackerTraits>
443  public:
444  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
445  ALPAKA_FN_ACC void operator()(TAcc const &acc,
446  CACellT<TrackerTraits> *__restrict__ cells,
447  uint32_t const *nCells) const {
449  for (auto idx : cms::alpakatools::uniform_elements(acc, (*nCells))) {
450  auto &thisCell = cells[idx];
451  if (!thisCell.tracks().empty())
452  thisCell.setStatusBits(Cell::StatusBit::kInTrack);
453  }
454  }
455  };
456 
457  template <typename TrackerTraits>
459  public:
460  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
461  ALPAKA_FN_ACC void operator()(TAcc const &acc,
464  for (auto it : cms::alpakatools::uniform_elements(acc, tracks_view.hitIndices().nOnes())) {
465  auto nhits = tracks_view.hitIndices().size(it);
466  if (nhits < 3)
467  continue;
468  if (tracks_view[it].quality() == Quality::edup)
469  continue;
470  ALPAKA_ASSERT_ACC(tracks_view[it].quality() == Quality::bad);
471  if (nhits > TrackerTraits::maxHitsOnTrack) // current limit
472  printf("wrong mult %d %d\n", it, nhits);
473  ALPAKA_ASSERT_ACC(nhits <= TrackerTraits::maxHitsOnTrack);
474  tupleMultiplicity->count(acc, nhits);
475  }
476  }
477  };
478 
479  template <typename TrackerTraits>
481  public:
482  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
483  ALPAKA_FN_ACC void operator()(TAcc const &acc,
486  for (auto it : cms::alpakatools::uniform_elements(acc, tracks_view.hitIndices().nOnes())) {
487  auto nhits = tracks_view.hitIndices().size(it);
488  if (nhits < 3)
489  continue;
490  if (tracks_view[it].quality() == Quality::edup)
491  continue;
492  ALPAKA_ASSERT_ACC(tracks_view[it].quality() == Quality::bad);
493  if (nhits > TrackerTraits::maxHitsOnTrack)
494  printf("wrong mult %d %d\n", it, nhits);
495  ALPAKA_ASSERT_ACC(nhits <= TrackerTraits::maxHitsOnTrack);
496  tupleMultiplicity->fill(acc, nhits, it);
497  }
498  }
499  };
500 
501  template <typename TrackerTraits>
503  public:
504  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
505  ALPAKA_FN_ACC void operator()(TAcc const &acc,
508  for (auto it : cms::alpakatools::uniform_elements(acc, tracks_view.hitIndices().nOnes())) {
509  auto nhits = tracks_view.hitIndices().size(it);
510  if (nhits == 0)
511  break; // guard
512 
513  // if duplicate: not even fit
514  if (tracks_view[it].quality() == Quality::edup)
515  continue;
516 
517  ALPAKA_ASSERT_ACC(tracks_view[it].quality() == Quality::bad);
518 
519  // mark doublets as bad
520  if (nhits < 3)
521  continue;
522 
523  // if the fit has any invalid parameters, mark it as bad
524  bool isNaN = false;
525  for (int i = 0; i < 5; ++i) {
526  isNaN |= edm::isNotFinite(tracks_view[it].state()(i));
527  }
528  if (isNaN) {
529 #ifdef NTUPLE_DEBUG
530  printf("NaN in fit %d size %d chi2 %f\n", it, tracks_view.hitIndices().size(it), tracks_view[it].chi2());
531 #endif
532  continue;
533  }
534 
535  tracks_view[it].quality() = Quality::strict;
536 
537  if (cuts.strictCut(tracks_view, it))
538  continue;
539 
540  tracks_view[it].quality() = Quality::tight;
541 
542  if (cuts.isHP(tracks_view, nhits, it))
544  }
545  }
546  };
547 
548  template <typename TrackerTraits>
550  public:
551  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
552  ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView<TrackerTraits> tracks_view, Counters *counters) const {
553  for (auto idx : cms::alpakatools::uniform_elements(acc, tracks_view.hitIndices().nOnes())) {
554  if (tracks_view.hitIndices().size(idx) == 0)
555  break; //guard
557  continue;
558  alpaka::atomicAdd(acc, &(counters->nLooseTracks), 1ull, alpaka::hierarchy::Blocks{});
559  if (tracks_view[idx].quality() < Quality::strict)
560  continue;
561  alpaka::atomicAdd(acc, &(counters->nGoodTracks), 1ull, alpaka::hierarchy::Blocks{});
562  }
563  }
564  };
565 
566  template <typename TrackerTraits>
568  public:
569  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
570  ALPAKA_FN_ACC void operator()(TAcc const &acc,
573  for (auto idx : cms::alpakatools::uniform_elements(acc, tracks_view.hitIndices().nOnes())) {
574  if (tracks_view.hitIndices().size(idx) == 0)
575  break; // guard
576  for (auto h = tracks_view.hitIndices().begin(idx); h != tracks_view.hitIndices().end(idx); ++h)
577  hitToTuple->count(acc, *h);
578  }
579  }
580  };
581 
582  template <typename TrackerTraits>
584  public:
585  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
586  ALPAKA_FN_ACC void operator()(TAcc const &acc,
589  for (auto idx : cms::alpakatools::uniform_elements(acc, tracks_view.hitIndices().nOnes())) {
590  if (tracks_view.hitIndices().size(idx) == 0)
591  break; // guard
592  for (auto h = tracks_view.hitIndices().begin(idx); h != tracks_view.hitIndices().end(idx); ++h)
593  hitToTuple->fill(acc, *h, idx);
594  }
595  }
596  };
597 
598  template <typename TrackerTraits>
600  public:
601  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
602  ALPAKA_FN_ACC void operator()(TAcc const &acc,
605  // copy offsets
606  for (auto idx : cms::alpakatools::uniform_elements(acc, tracks_view.hitIndices().nOnes())) {
607  tracks_view.detIndices().off[idx] = tracks_view.hitIndices().off[idx];
608  }
609  // fill hit indices
610  for (auto idx : cms::alpakatools::uniform_elements(acc, tracks_view.hitIndices().size())) {
611  ALPAKA_ASSERT_ACC(tracks_view.hitIndices().content[idx] < (uint32_t)hh.metadata().size());
612  tracks_view.detIndices().content[idx] = hh[tracks_view.hitIndices().content[idx]].detectorIndex();
613  }
614  }
615  };
616 
617  template <typename TrackerTraits>
619  public:
620  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
621  ALPAKA_FN_ACC void operator()(TAcc const &acc,
624  // clamp the number of tracks to the capacity of the SoA
625  auto ntracks = std::min<int>(apc->get().first, tracks_view.metadata().size() - 1);
626 
628  tracks_view.nTracks() = ntracks;
629  for (auto idx : cms::alpakatools::uniform_elements(acc, ntracks)) {
632  }
633  }
634  };
635 
636  template <typename TrackerTraits>
638  public:
639  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
640  ALPAKA_FN_ACC void operator()(TAcc const &acc,
641  HitToTuple<TrackerTraits> const *__restrict__ hitToTuple,
642  Counters *counters) const {
643  auto &c = *counters;
644  for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple->nOnes())) {
645  if (hitToTuple->size(idx) == 0)
646  continue; // SHALL NOT BE break
647  alpaka::atomicAdd(acc, &c.nUsedHits, 1ull, alpaka::hierarchy::Blocks{});
648  if (hitToTuple->size(idx) > 1)
649  alpaka::atomicAdd(acc, &c.nDupHits, 1ull, alpaka::hierarchy::Blocks{});
650  }
651  }
652  };
653 
654  template <typename TrackerTraits>
656  public:
657  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
658  ALPAKA_FN_ACC void operator()(TAcc const &acc,
659  int *__restrict__ nshared,
660  HitContainer<TrackerTraits> const *__restrict__ ptuples,
661  Quality const *__restrict__ quality,
662  HitToTuple<TrackerTraits> const *__restrict__ phitToTuple) const {
664 
665  auto &hitToTuple = *phitToTuple;
666  auto const &foundNtuplets = *ptuples;
667  for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple->nbins())) {
668  if (hitToTuple.size(idx) < 2)
669  continue;
670 
671  int nt = 0;
672 
673  // count "good" tracks
674  for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
675  if (quality[*it] < loose)
676  continue;
677  ++nt;
678  }
679 
680  if (nt < 2)
681  continue;
682 
683  // now mark each track triplet as sharing a hit
684  for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
685  if (foundNtuplets.size(*it) > 3)
686  continue;
687  alpaka::atomicAdd(acc, &nshared[*it], 1ull, alpaka::hierarchy::Blocks{});
688  }
689 
690  } // hit loop
691  }
692  };
693 
694  template <typename TrackerTraits>
696  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
697  ALPAKA_FN_ACC void operator()(TAcc const &acc,
698  int const *__restrict__ nshared,
699  HitContainer<TrackerTraits> const *__restrict__ tuples,
700  Quality *__restrict__ quality,
701  bool dupPassThrough) const {
702  // constexpr auto bad = Quality::bad;
703  constexpr auto dup = Quality::dup;
705  // constexpr auto strict = Quality::strict;
706 
707  // quality to mark rejected
708  auto const reject = dupPassThrough ? loose : dup;
709  for (auto idx : cms::alpakatools::uniform_elements(acc, tuples->nbins())) {
710  if (tuples->size(idx) == 0)
711  break; //guard
712  if (quality[idx] <= reject)
713  continue;
714  if (nshared[idx] > 2)
715  quality[idx] = reject;
716  }
717  }
718  };
719 
720  // mostly for very forward triplets.....
721  template <typename TrackerTraits>
723  public:
724  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
725  ALPAKA_FN_ACC void operator()(TAcc const &acc,
727  uint16_t nmin,
728  bool dupPassThrough,
729  HitToTuple<TrackerTraits> const *__restrict__ phitToTuple) const {
730  // quality to mark rejected
732 
733  auto &hitToTuple = *phitToTuple;
734 
735  for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple.nOnes())) {
736  if (hitToTuple.size(idx) < 2)
737  continue;
738 
739  auto score = [&](auto it, auto nl) { return std::abs(reco::tip(tracks_view, it)); };
740 
741  // full combinatorics
742  for (auto ip = hitToTuple.begin(idx); ip < hitToTuple.end(idx) - 1; ++ip) {
743  auto const it = *ip;
744  auto qi = tracks_view[it].quality();
745  if (qi <= reject)
746  continue;
747  auto opi = tracks_view[it].state()(2);
748  auto e2opi = tracks_view[it].covariance()(9);
749  auto cti = tracks_view[it].state()(3);
750  auto e2cti = tracks_view[it].covariance()(12);
751  auto nli = tracks_view[it].nLayers();
752  for (auto jp = ip + 1; jp < hitToTuple.end(idx); ++jp) {
753  auto const jt = *jp;
754  auto qj = tracks_view[jt].quality();
755  if (qj <= reject)
756  continue;
757  auto opj = tracks_view[jt].state()(2);
758  auto ctj = tracks_view[jt].state()(3);
759  auto dct = nSigma2 * (tracks_view[jt].covariance()(12) + e2cti);
760  if ((cti - ctj) * (cti - ctj) > dct)
761  continue;
762  auto dop = nSigma2 * (tracks_view[jt].covariance()(9) + e2opi);
763  if ((opi - opj) * (opi - opj) > dop)
764  continue;
765  auto nlj = tracks_view[jt].nLayers();
766  if (nlj < nli || (nlj == nli && (qj < qi || (qj == qi && score(it, nli) < score(jt, nlj)))))
767  tracks_view[jt].quality() = reject;
768  else {
769  tracks_view[it].quality() = reject;
770  break;
771  }
772  }
773  }
774  }
775  }
776  };
777 
778  template <typename TrackerTraits>
780  public:
781  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
782  ALPAKA_FN_ACC void operator()(TAcc const &acc,
785  int nmin,
786  bool dupPassThrough,
787  HitToTuple<TrackerTraits> const *__restrict__ phitToTuple) const {
788  // quality to mark rejected
790  // quality of longest track
791  auto const longTqual = Quality::highPurity;
792 
793  auto &hitToTuple = *phitToTuple;
794 
795  uint32_t l1end = hh.hitsLayerStart()[1];
796 
797  for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple.nOnes())) {
798  if (hitToTuple.size(idx) < 2)
799  continue;
800 
801  int8_t maxNl = 0;
802 
803  // find maxNl
804  for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
805  if (tracks_view[*it].quality() < longTqual)
806  continue;
807  // if (tracks_view[*it].nHits()==3) continue;
808  auto nl = tracks_view[*it].nLayers();
809  maxNl = std::max(nl, maxNl);
810  }
811 
812  if (maxNl < 4)
813  continue;
814 
815  // quad pass through (leave for tests)
816  // maxNl = std::min(4, maxNl);
817 
818  // kill all tracks shorter than maxHl (only triplets???
819  for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
820  auto nl = tracks_view[*it].nLayers();
821 
822  //checking if shared hit is on bpix1 and if the tuple is short enough
823  if (idx < l1end and nl > nmin)
824  continue;
825 
826  if (nl < maxNl && tracks_view[*it].quality() > reject)
827  tracks_view[*it].quality() = reject;
828  }
829  }
830  }
831  };
832  template <typename TrackerTraits>
834  public:
835  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
836  ALPAKA_FN_ACC void operator()(TAcc const &acc,
838  uint16_t nmin,
839  bool dupPassThrough,
840  HitToTuple<TrackerTraits> const *__restrict__ phitToTuple) const {
841  // quality to mark rejected
842  auto const reject = Quality::loose;
844  auto const good = Quality::strict;
845 
846  auto &hitToTuple = *phitToTuple;
847 
848  for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple.nOnes())) {
849  if (hitToTuple.size(idx) < 2)
850  continue;
851 
852  float mc = maxScore;
853  uint16_t im = tkNotFound;
854  bool onlyTriplets = true;
855 
856  // check if only triplets
857  for (auto it = hitToTuple.begin(idx); it != hitToTuple.end(idx); ++it) {
858  if (tracks_view[*it].quality() <= good)
859  continue;
860  onlyTriplets &= reco::isTriplet(tracks_view, *it);
861  if (!onlyTriplets)
862  break;
863  }
864 
865  // only triplets
866  if (!onlyTriplets)
867  continue;
868 
869  // for triplets choose best tip! (should we first find best quality???)
870  for (auto ip = hitToTuple.begin(idx); ip != hitToTuple.end(idx); ++ip) {
871  auto const it = *ip;
874  im = it;
875  }
876  }
877 
878  if (tkNotFound == im)
879  continue;
880 
881  // mark worse ambiguities
882  for (auto ip = hitToTuple.begin(idx); ip != hitToTuple.end(idx); ++ip) {
883  auto const it = *ip;
884  if (tracks_view[it].quality() > reject && it != im)
885  tracks_view[it].quality() = reject; //no race: simple assignment of the same constant
886  }
887 
888  } // loop over hits
889  }
890  };
891 
892  template <typename TrackerTraits>
894  public:
895  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
896  ALPAKA_FN_ACC void operator()(TAcc const &acc,
898  uint16_t nmin,
899  bool dupPassThrough,
900  HitToTuple<TrackerTraits> const *__restrict__ phitToTuple) const {
901  // quality to mark rejected
902  auto const reject = Quality::loose;
904  auto const good = Quality::loose;
905 
906  auto &hitToTuple = *phitToTuple;
907 
908  for (auto idx : cms::alpakatools::uniform_elements(acc, hitToTuple.nOnes())) {
909  if (hitToTuple.size(idx) < 2)
910  continue;
911 
912  float mc = maxScore;
913  uint16_t im = tkNotFound;
914 
915  // choose best tip! (should we first find best quality???)
916  for (auto ip = hitToTuple.begin(idx); ip != hitToTuple.end(idx); ++ip) {
917  auto const it = *ip;
920  im = it;
921  }
922  }
923 
924  if (tkNotFound == im)
925  continue;
926 
927  // mark worse ambiguities
928  for (auto ip = hitToTuple.begin(idx); ip != hitToTuple.end(idx); ++ip) {
929  auto const it = *ip;
930  if (tracks_view[it].quality() > reject && reco::isTriplet(tracks_view, it) && it != im)
931  tracks_view[it].quality() = reject; //no race: simple assignment of the same constant
932  }
933 
934  } // loop over hits
935  }
936  };
937 
938  template <typename TrackerTraits>
940  public:
941  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
942  ALPAKA_FN_ACC void operator()(TAcc const &acc,
945  HitToTuple<TrackerTraits> const *__restrict__ phitToTuple,
946  int32_t firstPrint,
947  int32_t lastPrint,
948  int iev) const {
950 
951  for (auto i :
953  auto nh = tracks_view.hitIndices().size(i);
954  if (nh < 3)
955  continue;
956  if (tracks_view[i].quality() < loose)
957  continue;
958  printf("TK: %d %d %d %d %f %f %f %f %f %f %f %.3f %.3f %.3f %.3f %.3f %.3f %.3f\n",
959  10000 * iev + i,
960  int(tracks_view[i].quality()),
961  nh,
962  tracks_view[i].nLayers(),
964  tracks_view[i].pt(),
965  tracks_view[i].eta(),
969  tracks_view[i].chi2(),
970  hh[*tracks_view.hitIndices().begin(i)].zGlobal(),
971  hh[*(tracks_view.hitIndices().begin(i) + 1)].zGlobal(),
972  hh[*(tracks_view.hitIndices().begin(i) + 2)].zGlobal(),
973  nh > 3 ? hh[int(*(tracks_view.hitIndices().begin(i) + 3))].zGlobal() : 0,
974  nh > 4 ? hh[int(*(tracks_view.hitIndices().begin(i) + 4))].zGlobal() : 0,
975  nh > 5 ? hh[int(*(tracks_view.hitIndices().begin(i) + 5))].zGlobal() : 0,
976  nh > 6 ? hh[int(*(tracks_view.hitIndices().begin(i) + nh - 1))].zGlobal() : 0);
977  }
978  }
979  };
980 
982  public:
983  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
984  ALPAKA_FN_ACC void operator()(TAcc const &acc, Counters const *counters) const {
985  auto const &c = *counters;
986  printf(
987  "||Counters | nEvents | nHits | nCells | nTuples | nFitTacks | nLooseTracks | nGoodTracks | nUsedHits | "
988  "nDupHits | nFishCells | nKilledCells | nUsedCells | nZeroTrackCells ||\n");
989  printf("Counters Raw %lld %lld %lld %lld %lld %lld %lld %lld %lld %lld %lld %lld %lld\n",
990  c.nEvents,
991  c.nHits,
992  c.nCells,
993  c.nTuples,
994  c.nFitTracks,
995  c.nLooseTracks,
996  c.nGoodTracks,
997  c.nUsedHits,
998  c.nDupHits,
999  c.nFishCells,
1000  c.nKilledCells,
1001  c.nEmptyCells,
1002  c.nZeroTrackCells);
1003  printf(
1004  "Counters Norm %lld || %.1f| %.1f| %.1f| %.1f| %.1f| %.1f| %.1f| %.1f| %.3f| %.3f| %.3f| "
1005  "%.3f||\n",
1006  c.nEvents,
1007  c.nHits / double(c.nEvents),
1008  c.nCells / double(c.nEvents),
1009  c.nTuples / double(c.nEvents),
1010  c.nFitTracks / double(c.nEvents),
1011  c.nLooseTracks / double(c.nEvents),
1012  c.nGoodTracks / double(c.nEvents),
1013  c.nUsedHits / double(c.nEvents),
1014  c.nDupHits / double(c.nEvents),
1015  c.nFishCells / double(c.nCells),
1016  c.nKilledCells / double(c.nCells),
1017  c.nEmptyCells / double(c.nCells),
1018  c.nZeroTrackCells / double(c.nCells));
1019  }
1020  };
1021 
1022 } // namespace ALPAKA_ACCELERATOR_NAMESPACE::caHitNtupletGeneratorKernels
1023 
1024 #endif // RecoTracker_PixelSeeding_plugins_alpaka_CAHitNtupletGeneratorKernelsImpl_h
HitContainer< TrackerTraits > const *__restrict__ tuples
TupleMultiplicity< TrackerTraits > const HitToTuple< TrackerTraits > const cms::cuda::AtomicPairCounter GPUCACellT< TrackerTraits > const *__restrict__ uint32_t const *__restrict__ CellNeighborsVector< TrackerTraits > const CellTracksVector< TrackerTraits > const OuterHitOfCell< TrackerTraits > const isOuterHitOfCell
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, cms::alpakatools::AtomicPairCounter *apc) const
ALPAKA_FN_ACC auto uniform_elements(TAcc const &acc, TArgs... args)
Definition: workdivision.h:311
TupleMultiplicity< TrackerTraits > const HitToTuple< TrackerTraits > const cms::cuda::AtomicPairCounter GPUCACellT< TrackerTraits > const *__restrict__ uint32_t const *__restrict__ CellNeighborsVector< TrackerTraits > const * cellNeighbors
TupleMultiplicity< TrackerTraits > const HitToTuple< TrackerTraits > const cms::cuda::AtomicPairCounter GPUCACellT< TrackerTraits > const *__restrict__ uint32_t const *__restrict__ CellNeighborsVector< TrackerTraits > const CellTracksVector< TrackerTraits > const * cellTracks
uint32_t const *__restrict__ TkSoAView< TrackerTraits > tracks_view
ALPAKA_FN_ACC void operator()(TAcc const &acc, cms::alpakatools::AtomicPairCounter *apc1, cms::alpakatools::AtomicPairCounter *apc2, HitsConstView< TrackerTraits > hh, CACellT< TrackerTraits > *cells, uint32_t *nCells, CellNeighborsVector< TrackerTraits > *cellNeighbors, OuterHitOfCell< TrackerTraits > const *isOuterHitOfCell, CAParams< TrackerTraits > params) const
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:90
caHitNtupletGenerator::Counters Counters
typename reco::TrackSoA< TrackerTraits >::HitContainer HitContainer
ALPAKA_FN_ACC void operator()(TAcc const &acc, HitToTuple< TrackerTraits > const *__restrict__ hitToTuple, Counters *counters) const
TkSoAView< TrackerTraits > HitToTuple< TrackerTraits > const *__restrict__ int32_t int32_t int iev
TupleMultiplicity< TrackerTraits > const HitToTuple< TrackerTraits > const cms::cuda::AtomicPairCounter GPUCACellT< TrackerTraits > const *__restrict__ cells
TupleMultiplicity< TrackerTraits > const * tupleMultiplicity
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, TupleMultiplicity< TrackerTraits > const *tupleMultiplicity, HitToTuple< TrackerTraits > const *hitToTuple, cms::alpakatools::AtomicPairCounter *apc, CACellT< TrackerTraits > const *__restrict__ cells, uint32_t const *__restrict__ nCells, CellNeighborsVector< TrackerTraits > const *cellNeighbors, CellTracksVector< TrackerTraits > const *cellTracks, OuterHitOfCell< TrackerTraits > const *isOuterHitOfCell, int32_t nHits, uint32_t maxNumberOfDoublets, Counters *counters) const
uint32_t const *__restrict__ TkSoAView< TrackerTraits > bool dupPassThrough
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float charge(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:73
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, HitToTuple< TrackerTraits > *hitToTuple) const
TupleMultiplicity< TrackerTraits > const HitToTuple< TrackerTraits > const cms::cuda::AtomicPairCounter GPUCACellT< TrackerTraits > const *__restrict__ uint32_t const *__restrict__ CellNeighborsVector< TrackerTraits > const CellTracksVector< TrackerTraits > const OuterHitOfCell< TrackerTraits > const int32_t uint32_t maxNumberOfDoublets
TkSoAView< TrackerTraits > GPUCACellT< TrackerTraits > *__restrict__ uint32_t const CellTracksVector< TrackerTraits > cms::cuda::AtomicPairCounter CAParams< TrackerTraits > params
string quality
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float tip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:85
ALPAKA_FN_ACC void operator()(TAcc const &acc, CACellT< TrackerTraits > const *__restrict__ cells, uint32_t const *__restrict__ nCells, TkSoAView< TrackerTraits > tracks_view, bool dupPassThrough) const
TupleMultiplicity< TrackerTraits > const HitToTuple< TrackerTraits > const cms::cuda::AtomicPairCounter GPUCACellT< TrackerTraits > const *__restrict__ uint32_t const *__restrict__ CellNeighborsVector< TrackerTraits > const CellTracksVector< TrackerTraits > const OuterHitOfCell< TrackerTraits > const int32_t uint32_t Counters * counters
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, uint16_t nmin, bool dupPassThrough, HitToTuple< TrackerTraits > const *__restrict__ phitToTuple) const
ALPAKA_FN_ACC void operator()(TAcc const &acc, HitsConstView< TrackerTraits > hh, TkSoAView< TrackerTraits > tracks_view, CACellT< TrackerTraits > *__restrict__ cells, uint32_t const *nCells, CellTracksVector< TrackerTraits > *cellTracks, cms::alpakatools::AtomicPairCounter *apc, CAParams< TrackerTraits > params) const
static constexpr __host__ __device__ int computeNumberOfLayers(const TrackSoAConstView &tracks, int32_t i)
ALPAKA_FN_ACC auto uniform_elements_y(TAcc const &acc, TArgs... args)
Definition: workdivision.h:344
ALPAKA_FN_ACC void operator()(TAcc const &acc, HitsConstView< TrackerTraits > hh, TkSoAView< TrackerTraits > tracks_view, int nmin, bool dupPassThrough, HitToTuple< TrackerTraits > const *__restrict__ phitToTuple) const
TkSoAView< TrackerTraits > HitToTuple< TrackerTraits > const *__restrict__ int32_t firstPrint
__device__ __host__ Counters get() const
stack
Definition: svgfig.py:559
std::vector< Block > Blocks
Definition: Block.h:99
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float phi(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:80
typename CACellT< TrackerTraits >::HitsConstView HitsConstView
ALPAKA_FN_ACC void operator()(TAcc const &acc, CACellT< TrackerTraits > *__restrict__ cells, uint32_t const *nCells) const
ALPAKA_FN_ACC auto independent_group_elements_x(TAcc const &acc, TArgs... args)
TupleMultiplicity< TrackerTraits > const HitToTuple< TrackerTraits > const * hitToTuple
ALPAKA_FN_ACC void operator()(TAcc const &acc, int *__restrict__ nshared, HitContainer< TrackerTraits > const *__restrict__ ptuples, Quality const *__restrict__ quality, HitToTuple< TrackerTraits > const *__restrict__ phitToTuple) const
uint32_t nh
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr bool isTriplet(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:95
ALPAKA_FN_ACC void operator()(TAcc const &acc, int const *__restrict__ nshared, HitContainer< TrackerTraits > const *__restrict__ tuples, Quality *__restrict__ quality, bool dupPassThrough) const
int nt
Definition: AMPTWrapper.h:42
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, uint16_t nmin, bool dupPassThrough, HitToTuple< TrackerTraits > const *__restrict__ phitToTuple) const
TupleMultiplicity< TrackerTraits > const HitToTuple< TrackerTraits > const cms::cuda::AtomicPairCounter GPUCACellT< TrackerTraits > const *__restrict__ uint32_t const *__restrict__ nCells
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, TupleMultiplicity< TrackerTraits > *tupleMultiplicity) const
ALPAKA_FN_ACC constexpr bool once_per_grid(TAcc const &acc)
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, Counters *counters) const
ALPAKA_FN_ACC void operator()(TAcc const &acc, CACellT< TrackerTraits > const *cells, uint32_t const *__restrict__ nCells, TkSoAView< TrackerTraits > tracks_view, bool dupPassThrough) const
ALPAKA_FN_ACC void operator()(TAcc const &acc, CACellT< TrackerTraits > const *cells, uint32_t const *__restrict__ nCells, TkSoAView< TrackerTraits > tracks_view) const
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, HitToTuple< TrackerTraits > *hitToTuple) const
TupleMultiplicity< TrackerTraits > const HitToTuple< TrackerTraits > const cms::cuda::AtomicPairCounter * apc
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, QualityCuts< TrackerTraits > cuts) const
TkSoAView< TrackerTraits > HitToTuple< TrackerTraits > const *__restrict__ int32_t int32_t lastPrint
ALPAKA_FN_ACC void operator()(TAcc const &acc, Counters const *counters) const
ALPAKA_FN_ACC void operator()(TAcc const &acc, HitsConstView< TrackerTraits > hh, TkSoAView< TrackerTraits > tracks_view, HitToTuple< TrackerTraits > const *__restrict__ phitToTuple, int32_t firstPrint, int32_t lastPrint, int iev) const
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, HitsConstView< TrackerTraits > hh) const
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, uint16_t nmin, bool dupPassThrough, HitToTuple< TrackerTraits > const *__restrict__ phitToTuple) const
TrackingRecHitSoAConstView< TrackerTraits > HitsConstView
Definition: CACell.h:36
ALPAKA_FN_ACC void operator()(TAcc const &acc, TkSoAView< TrackerTraits > tracks_view, TupleMultiplicity< TrackerTraits > *tupleMultiplicity) const
HitContainer< TrackerTraits > const *__restrict__ Quality const *__restrict__ HitToTuple< TrackerTraits > const *__restrict__ phitToTuple
T1 atomicAdd(T1 *a, T2 b)
Definition: cudaCompat.h:61
HitContainer< TrackerTraits > const *__restrict__ ptuples
typename reco::TrackSoA< TrackerTraits >::template Layout<>::View TrackSoAView
Definition: TracksSoA.h:45