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