CMS 3D CMS Logo

CAHitNtupletGeneratorKernels.dev.cc
Go to the documentation of this file.
1 // C++ headers
2 #ifdef DUMP_GPU_TK_TUPLES
3 #include <mutex>
4 #endif
5 
6 // Alpaka headers
7 #include <alpaka/alpaka.hpp>
8 
9 // CMSSW headers
16 
17 // local headers
18 #include "CAFishbone.h"
21 
22 //#define GPU_DEBUG
23 //#define NTUPLE_DEBUG
24 
26 
27  template <typename TrackerTraits>
29  uint32_t nhits,
30  uint32_t offsetBPIX2,
31  Queue &queue)
32  : m_params(params),
34  // ALLOCATIONS FOR THE INTERMEDIATE RESULTS (STAYS ON WORKER)
36  counters_{cms::alpakatools::make_device_buffer<Counters>(queue)},
37 
38  // workspace
39  device_hitToTuple_{cms::alpakatools::make_device_buffer<HitToTuple>(queue)},
40  device_hitToTupleStorage_{
41  cms::alpakatools::make_device_buffer<typename HitToTuple::Counter[]>(queue, nhits + 1)},
42  device_tupleMultiplicity_{cms::alpakatools::make_device_buffer<TupleMultiplicity>(queue)},
43 
44  // NB: In legacy, device_theCells_ and device_isOuterHitOfCell_ were allocated inside buildDoublets
45  device_theCells_{
46  cms::alpakatools::make_device_buffer<CACell[]>(queue, m_params.caParams_.maxNumberOfDoublets_)},
47  // in principle we can use "nhits" to heuristically dimension the workspace...
48  device_isOuterHitOfCell_{
49  cms::alpakatools::make_device_buffer<OuterHitOfCellContainer[]>(queue, std::max(1u, nhits - offsetBPIX2))},
50  isOuterHitOfCell_{cms::alpakatools::make_device_buffer<OuterHitOfCell>(queue)},
51 
52  device_theCellNeighbors_{cms::alpakatools::make_device_buffer<CellNeighborsVector>(queue)},
53  device_theCellTracks_{cms::alpakatools::make_device_buffer<CellTracksVector>(queue)},
54  // NB: In legacy, cellStorage_ was allocated inside buildDoublets
55  cellStorage_{cms::alpakatools::make_device_buffer<unsigned char[]>(
56  queue,
57  TrackerTraits::maxNumOfActiveDoublets * sizeof(CellNeighbors) +
58  TrackerTraits::maxNumOfActiveDoublets * sizeof(CellTracks))},
59  device_cellCuts_{cms::alpakatools::make_device_buffer<CellCuts>(queue)},
60  device_theCellNeighborsContainer_{reinterpret_cast<CellNeighbors *>(cellStorage_.data())},
61  device_theCellTracksContainer_{reinterpret_cast<CellTracks *>(
62  cellStorage_.data() + TrackerTraits::maxNumOfActiveDoublets * sizeof(CellNeighbors))},
63 
64  // NB: In legacy, device_storage_ was allocated inside allocateOnGPU
65  device_storage_{
66  cms::alpakatools::make_device_buffer<cms::alpakatools::AtomicPairCounter::DoubleWord[]>(queue, 3u)},
67  device_hitTuple_apc_{reinterpret_cast<cms::alpakatools::AtomicPairCounter *>(device_storage_.data())},
68  device_hitToTuple_apc_{reinterpret_cast<cms::alpakatools::AtomicPairCounter *>(device_storage_.data() + 1)},
69  device_nCells_{cms::alpakatools::make_device_view(alpaka::getDev(queue),
70  *reinterpret_cast<uint32_t *>(device_storage_.data() + 2))} {
71 #ifdef GPU_DEBUG
72  std::cout << "Allocation for tuple building. N hits " << nhits << std::endl;
73 #endif
74 
75  alpaka::memset(queue, counters_, 0);
76  alpaka::memset(queue, device_nCells_, 0);
77  alpaka::memset(queue, cellStorage_, 0);
78 
79  auto cellCuts_h = cms::alpakatools::make_host_view(m_params.cellCuts_);
80  alpaka::memcpy(queue, device_cellCuts_, cellCuts_h);
81 
82  [[maybe_unused]] TupleMultiplicity *tupleMultiplicityDeviceData = device_tupleMultiplicity_.data();
83  using TM = cms::alpakatools::OneToManyAssocRandomAccess<typename TrackerTraits::tindex_type,
84  TrackerTraits::maxHitsOnTrack + 1,
85  TrackerTraits::maxNumberOfTuples>;
86  TM *tm = device_tupleMultiplicity_.data();
87  TM::template launchZero<Acc1D>(tm, queue);
88  TupleMultiplicity::template launchZero<Acc1D>(tupleMultiplicityDeviceData, queue);
89 
90  device_hitToTupleView_.assoc = device_hitToTuple_.data();
91  device_hitToTupleView_.offStorage = device_hitToTupleStorage_.data();
92  device_hitToTupleView_.offSize = nhits + 1;
93 
94  HitToTuple::template launchZero<Acc1D>(device_hitToTupleView_, queue);
95 #ifdef GPU_DEBUG
96  std::cout << "Allocations for CAHitNtupletGeneratorKernels: done!" << std::endl;
97 #endif
98  }
99 
100  template <typename TrackerTraits>
102  uint32_t offsetBPIX2,
104  Queue &queue) {
105  using namespace caPixelDoublets;
106  using namespace caHitNtupletGeneratorKernels;
107 
108  // zero tuples
109  HitContainer::template launchZero<Acc1D>(&(tracks_view.hitIndices()), queue);
110 
111  uint32_t nhits = hh.metadata().size();
112 
113 #ifdef NTUPLE_DEBUG
114  std::cout << "start tuple building. N hits " << nhits << std::endl;
115  if (nhits < 2)
116  std::cout << "too few hits " << nhits << std::endl;
117 #endif
118 
119  //
120  // applying combinatoric cleaning such as fishbone at this stage is too expensive
121  //
122 
123  const auto nthTot = 64;
124  const auto stride = 4;
125  auto blockSize = nthTot / stride;
126  auto numberOfBlocks = cms::alpakatools::divide_up_by(3 * m_params.caParams_.maxNumberOfDoublets_ / 4, blockSize);
127  const auto rescale = numberOfBlocks / 65536;
128  blockSize *= (rescale + 1);
129  numberOfBlocks = cms::alpakatools::divide_up_by(3 * m_params.caParams_.maxNumberOfDoublets_ / 4, blockSize);
130  assert(numberOfBlocks < 65536);
131  assert(blockSize > 0 && 0 == blockSize % 16);
132  const Vec2D blks{numberOfBlocks, 1u};
133  const Vec2D thrs{blockSize, stride};
134  const auto kernelConnectWorkDiv = cms::alpakatools::make_workdiv<Acc2D>(blks, thrs);
135 
136  alpaka::exec<Acc2D>(queue,
137  kernelConnectWorkDiv,
138  Kernel_connect<TrackerTraits>{},
139  this->device_hitTuple_apc_,
140  this->device_hitToTuple_apc_, // needed only to be reset, ready for next kernel
141  hh,
142  this->device_theCells_.data(),
143  this->device_nCells_.data(),
144  this->device_theCellNeighbors_.data(),
145  this->isOuterHitOfCell_.data(),
146  this->m_params.caParams_);
147 
148  // do not run the fishbone if there are hits only in BPIX1
149  if (this->m_params.earlyFishbone_ and nhits > offsetBPIX2) {
150  const auto nthTot = 128;
151  const auto stride = 16;
152  const auto blockSize = nthTot / stride;
153  const auto numberOfBlocks = cms::alpakatools::divide_up_by(nhits - offsetBPIX2, blockSize);
154  const Vec2D blks{numberOfBlocks, 1u};
155  const Vec2D thrs{blockSize, stride};
156  const auto fishboneWorkDiv = cms::alpakatools::make_workdiv<Acc2D>(blks, thrs);
157  alpaka::exec<Acc2D>(queue,
158  fishboneWorkDiv,
160  hh,
161  this->device_theCells_.data(),
162  this->device_nCells_.data(),
163  this->isOuterHitOfCell_.data(),
164  nhits,
165  false);
166  }
167  blockSize = 64;
168  numberOfBlocks = cms::alpakatools::divide_up_by(3 * m_params.caParams_.maxNumberOfDoublets_ / 4, blockSize);
169  auto workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
170  alpaka::exec<Acc1D>(queue,
171  workDiv1D,
172  Kernel_find_ntuplets<TrackerTraits>{},
173  hh,
174  tracks_view,
175  this->device_theCells_.data(),
176  this->device_nCells_.data(),
177  this->device_theCellTracks_.data(),
178  this->device_hitTuple_apc_,
179  this->m_params.caParams_);
180 #ifdef GPU_DEBUG
182 #endif
183 
184  if (this->m_params.doStats_)
185  alpaka::exec<Acc1D>(queue,
186  workDiv1D,
187  Kernel_mark_used<TrackerTraits>{},
188  this->device_theCells_.data(),
189  this->device_nCells_.data());
190 
191 #ifdef GPU_DEBUG
193 #endif
194 
195  blockSize = 128;
196  numberOfBlocks = cms::alpakatools::divide_up_by(HitContainer{}.totOnes(), blockSize);
197  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
198 
199  alpaka::exec<Acc1D>(
200  queue, workDiv1D, typename HitContainer::finalizeBulk{}, this->device_hitTuple_apc_, &tracks_view.hitIndices());
201 
202 #ifdef GPU_DEBUG
204 #endif
205 
206  alpaka::exec<Acc1D>(queue, workDiv1D, Kernel_fillHitDetIndices<TrackerTraits>{}, tracks_view, hh);
207 
208 #ifdef GPU_DEBUG
210 #endif
211  alpaka::exec<Acc1D>(queue, workDiv1D, Kernel_fillNLayers<TrackerTraits>{}, tracks_view, this->device_hitTuple_apc_);
212 
213 #ifdef GPU_DEBUG
215 #endif
216 
217  // remove duplicates (tracks that share a doublet)
218  numberOfBlocks = cms::alpakatools::divide_up_by(3 * m_params.caParams_.maxNumberOfDoublets_ / 4, blockSize);
219  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
220 
221  alpaka::exec<Acc1D>(queue,
222  workDiv1D,
223  Kernel_earlyDuplicateRemover<TrackerTraits>{},
224  this->device_theCells_.data(),
225  this->device_nCells_.data(),
226  tracks_view,
227  this->m_params.dupPassThrough_);
228 #ifdef GPU_DEBUG
230 #endif
231 
232  blockSize = 128;
233  numberOfBlocks = cms::alpakatools::divide_up_by(3 * TrackerTraits::maxNumberOfTuples / 4, blockSize);
234  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
235 
236  alpaka::exec<Acc1D>(queue,
237  workDiv1D,
238  Kernel_countMultiplicity<TrackerTraits>{},
239  tracks_view,
240  this->device_tupleMultiplicity_.data());
241  TupleMultiplicity::template launchFinalize<Acc1D>(this->device_tupleMultiplicity_.data(), queue);
242 
243  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
244  alpaka::exec<Acc1D>(
245  queue, workDiv1D, Kernel_fillMultiplicity<TrackerTraits>{}, tracks_view, this->device_tupleMultiplicity_.data());
246 #ifdef GPU_DEBUG
248 #endif
249  // do not run the fishbone if there are hits only in BPIX1
250  if (this->m_params.lateFishbone_ and nhits > offsetBPIX2) {
251  const auto nthTot = 128;
252  const auto stride = 16;
253  const auto blockSize = nthTot / stride;
254  const auto numberOfBlocks = cms::alpakatools::divide_up_by(nhits - offsetBPIX2, blockSize);
255  const Vec2D blks{numberOfBlocks, 1u};
256  const Vec2D thrs{blockSize, stride};
257  const auto workDiv2D = cms::alpakatools::make_workdiv<Acc2D>(blks, thrs);
258 
259  alpaka::exec<Acc2D>(queue,
260  workDiv2D,
262  hh,
263  this->device_theCells_.data(),
264  this->device_nCells_.data(),
265  this->isOuterHitOfCell_.data(),
266  nhits,
267  true);
268  }
269 
270 #ifdef GPU_DEBUG
272 #endif
273  }
274 
275  template <typename TrackerTraits>
277  uint32_t offsetBPIX2,
278  Queue &queue) {
279  using namespace caPixelDoublets;
281  using OuterHitOfCell = typename CACell::OuterHitOfCell;
282  using CellNeighbors = typename CACell::CellNeighbors;
283  using CellTracks = typename CACell::CellTracks;
285 
286  auto nhits = hh.metadata().size();
287 #ifdef NTUPLE_DEBUG
288  std::cout << "building Doublets out of " << nhits << " Hits" << std::endl;
289 #endif
290 
291 #ifdef GPU_DEBUG
293 #endif
294 
295  // in principle we can use "nhits" to heuristically dimension the workspace...
296  ALPAKA_ASSERT_ACC(this->device_isOuterHitOfCell_.data());
297 
298  alpaka::exec<Acc1D>(
299  queue,
300  cms::alpakatools::make_workdiv<Acc1D>(1, 1),
301  [] ALPAKA_FN_ACC(Acc1D const &acc,
303  OuterHitOfCellContainer *container,
304  int32_t const *offset) {
305  // this code runs on the device
306  isOuterHitOfCell->container = container;
307  isOuterHitOfCell->offset = *offset;
308  },
309  this->isOuterHitOfCell_.data(),
310  this->device_isOuterHitOfCell_.data(),
311  &hh.offsetBPIX2());
312 
313  {
314  int threadsPerBlock = 128;
315  // at least one block!
316  int blocks = std::max(1u, cms::alpakatools::divide_up_by(nhits - offsetBPIX2, threadsPerBlock));
317  const auto workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(blocks, threadsPerBlock);
318 
319  alpaka::exec<Acc1D>(queue,
320  workDiv1D,
322  this->isOuterHitOfCell_.data(),
323  nhits,
324  this->device_theCellNeighbors_.data(),
325  this->device_theCellNeighborsContainer_,
326  this->device_theCellTracks_.data(),
327  this->device_theCellTracksContainer_);
328  }
329 
330 #ifdef GPU_DEBUG
332 #endif
333 
334  if (0 == nhits)
335  return; // protect against empty events
336 
337  // take all layer pairs into account
338  auto nActualPairs = this->m_params.nPairs();
339 
340  const int stride = 4;
341  const int threadsPerBlock = TrackerTraits::getDoubletsFromHistoMaxBlockSize / stride;
342  int blocks = (4 * nhits + threadsPerBlock - 1) / threadsPerBlock;
343  const Vec2D blks{blocks, 1u};
344  const Vec2D thrs{threadsPerBlock, stride};
345  const auto workDiv2D = cms::alpakatools::make_workdiv<Acc2D>(blks, thrs);
346 
347  alpaka::exec<Acc2D>(queue,
348  workDiv2D,
350  this->device_theCells_.data(),
351  this->device_nCells_.data(),
352  this->device_theCellNeighbors_.data(),
353  this->device_theCellTracks_.data(),
354  hh,
355  this->isOuterHitOfCell_.data(),
356  nActualPairs,
357  this->m_params.caParams_.maxNumberOfDoublets_,
358  this->m_params.cellCuts_);
359 
360 #ifdef GPU_DEBUG
362 #endif
363  }
364 
365  template <typename TrackerTraits>
368  Queue &queue) {
369  using namespace caHitNtupletGeneratorKernels;
370 
371  uint32_t nhits = hh.metadata().size();
372 
373  auto blockSize = 64;
374 
375  // classify tracks based on kinematics
376  auto numberOfBlocks = cms::alpakatools::divide_up_by(3 * TrackerTraits::maxNumberOfQuadruplets / 4, blockSize);
377  auto workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
378  alpaka::exec<Acc1D>(
379  queue, workDiv1D, Kernel_classifyTracks<TrackerTraits>{}, tracks_view, this->m_params.qualityCuts_);
380 
381  if (this->m_params.lateFishbone_) {
382  // apply fishbone cleaning to good tracks
383  numberOfBlocks = cms::alpakatools::divide_up_by(3 * m_params.caParams_.maxNumberOfDoublets_ / 4, blockSize);
384  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
385  alpaka::exec<Acc1D>(queue,
386  workDiv1D,
387  Kernel_fishboneCleaner<TrackerTraits>{},
388  this->device_theCells_.data(),
389  this->device_nCells_.data(),
390  tracks_view);
391  }
392 
393  // mark duplicates (tracks that share a doublet)
394  numberOfBlocks = cms::alpakatools::divide_up_by(3 * m_params.caParams_.maxNumberOfDoublets_ / 4, blockSize);
395  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
396  alpaka::exec<Acc1D>(queue,
397  workDiv1D,
398  Kernel_fastDuplicateRemover<TrackerTraits>{},
399  this->device_theCells_.data(),
400  this->device_nCells_.data(),
401  tracks_view,
402  this->m_params.dupPassThrough_);
403 #ifdef GPU_DEBUG
405 #endif
406 
407  if (this->m_params.doSharedHitCut_ || this->m_params.doStats_) {
408  // fill hit->track "map"
409  numberOfBlocks = cms::alpakatools::divide_up_by(3 * TrackerTraits::maxNumberOfQuadruplets / 4, blockSize);
410  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
411  alpaka::exec<Acc1D>(queue,
412  workDiv1D,
413  Kernel_countHitInTracks<TrackerTraits>{},
414  tracks_view,
415  this->device_hitToTuple_.data()); //CHECK
416 
417  HitToTuple::template launchFinalize<Acc1D>(this->device_hitToTupleView_, queue);
418  alpaka::exec<Acc1D>(
419  queue, workDiv1D, Kernel_fillHitInTracks<TrackerTraits>{}, tracks_view, this->device_hitToTuple_.data());
420 #ifdef GPU_DEBUG
422 #endif
423  }
424 
425  if (this->m_params.doSharedHitCut_) {
426  // mark duplicates (tracks that share at least one hit)
427  numberOfBlocks = cms::alpakatools::divide_up_by(3 * TrackerTraits::maxNumberOfQuadruplets / 4,
428  blockSize); // TODO: Check if correct
429  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
430  alpaka::exec<Acc1D>(queue,
431  workDiv1D,
432  Kernel_rejectDuplicate<TrackerTraits>{},
433  tracks_view,
434  this->m_params.minHitsForSharingCut_,
435  this->m_params.dupPassThrough_,
436  this->device_hitToTuple_.data());
437 
438  alpaka::exec<Acc1D>(queue,
439  workDiv1D,
440  Kernel_sharedHitCleaner<TrackerTraits>{},
441  hh,
442  tracks_view,
443  this->m_params.minHitsForSharingCut_,
444  this->m_params.dupPassThrough_,
445  this->device_hitToTuple_.data());
446 
447  if (this->m_params.useSimpleTripletCleaner_) {
448  // (typename HitToTuple{}::capacity(),
449  numberOfBlocks = cms::alpakatools::divide_up_by(HitToTuple{}.capacity(), blockSize);
450  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
451  alpaka::exec<Acc1D>(queue,
452  workDiv1D,
453  Kernel_simpleTripletCleaner<TrackerTraits>{},
454  tracks_view,
455  this->m_params.minHitsForSharingCut_,
456  this->m_params.dupPassThrough_,
457  this->device_hitToTuple_.data());
458  } else {
459  numberOfBlocks = cms::alpakatools::divide_up_by(HitToTuple{}.capacity(), blockSize);
460  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
461  alpaka::exec<Acc1D>(queue,
462  workDiv1D,
463  Kernel_tripletCleaner<TrackerTraits>{},
464  tracks_view,
465  this->m_params.minHitsForSharingCut_,
466  this->m_params.dupPassThrough_,
467  this->device_hitToTuple_.data());
468  }
469 #ifdef GPU_DEBUG
471 #endif
472  }
473 
474  if (this->m_params.doStats_) {
475  numberOfBlocks =
476  cms::alpakatools::divide_up_by(std::max(nhits, m_params.caParams_.maxNumberOfDoublets_), blockSize);
477  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
478 
479  alpaka::exec<Acc1D>(queue,
480  workDiv1D,
481  Kernel_checkOverflows<TrackerTraits>{},
482  tracks_view,
483  this->device_tupleMultiplicity_.data(),
484  this->device_hitToTuple_.data(),
485  this->device_hitTuple_apc_,
486  this->device_theCells_.data(),
487  this->device_nCells_.data(),
488  this->device_theCellNeighbors_.data(),
489  this->device_theCellTracks_.data(),
490  this->isOuterHitOfCell_.data(),
491  nhits,
492  this->m_params.caParams_.maxNumberOfDoublets_,
493  this->counters_.data());
494  }
495 
496  if (this->m_params.doStats_) {
497  // counters (add flag???)
498 
499  numberOfBlocks = cms::alpakatools::divide_up_by(HitToTuple{}.capacity(), blockSize);
500  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
501  alpaka::exec<Acc1D>(queue,
502  workDiv1D,
503  Kernel_doStatsForHitInTracks<TrackerTraits>{},
504  this->device_hitToTuple_.data(),
505  this->counters_.data());
506 
507  numberOfBlocks = cms::alpakatools::divide_up_by(3 * TrackerTraits::maxNumberOfQuadruplets / 4, blockSize);
508  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(numberOfBlocks, blockSize);
509  alpaka::exec<Acc1D>(
510  queue, workDiv1D, Kernel_doStatsForTracks<TrackerTraits>{}, tracks_view, this->counters_.data());
511  }
512 #ifdef GPU_DEBUG
514 #endif
515 
516 #ifdef DUMP_GPU_TK_TUPLES
517  static std::atomic<int> iev(0);
518  static std::mutex lock;
519  workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(1u, 32u);
520  {
521  std::lock_guard<std::mutex> guard(lock);
522  ++iev;
523  for (int k = 0; k < 20000; k += 500) {
524  alpaka::exec<Acc1D>(queue,
525  workDiv1D,
526  Kernel_print_found_ntuplets<TrackerTraits>{},
527  hh,
528  tracks_view,
529  this->device_hitToTuple_.data(),
530  k,
531  k + 500,
532  iev);
534  }
535  alpaka::exec<Acc1D>(queue,
536  workDiv1D,
537  Kernel_print_found_ntuplets<TrackerTraits>{},
538  hh,
539  tracks_view,
540  this->device_hitToTuple_.data(),
541  20000,
542  1000000,
543  iev);
544 
546  }
547 #endif
548  }
549 
550  /* This will make sense when we will be able to run this once per job in Alpaka
551 
552  template <typename TrackerTraits>
553  void CAHitNtupletGeneratorKernels<TrackerTraits>::printCounters() {
554  auto workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(1,1);
555  alpaka::exec<Acc1D>(queue_, workDiv1D, Kernel_printCounters{}, this->counters_.data());
556  }
557  */
558 
562 
563 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
constexpr Idx divide_up_by(Idx value, Idx divisor)
Definition: workdivision.h:20
uint32_t const *__restrict__ TkSoAView< TrackerTraits > tracks_view
Vec< Dim2D > Vec2D
Definition: config.h:26
caStructures::CellNeighborsT< TrackerTraits > CellNeighbors
Definition: CAFishbone.h:23
static std::mutex mutex
Definition: Proxy.cc:8
CAHitNtupletGeneratorKernels(Params const &params, uint32_t nhits, uint32_t offsetBPIX2, Queue &queue)
TkSoAView< TrackerTraits > HitToTuple< TrackerTraits > const *__restrict__ int32_t int32_t int iev
caStructures::template HitToTupleT< TrackerTraits > HitToTuple
std::enable_if_t< not std::is_array_v< T >, device_view< TDev, T > > make_device_view(TDev const &device, T &data)
Definition: memory.h:260
assert(be >=bs)
typename reco::TrackSoA< TrackerTraits >::HitContainer HitContainer
constexpr uint32_t stride
Definition: HelixFit.h:22
def template(fileName, svg, replaceme="REPLACEME")
Definition: svgfig.py:521
uint32_t CellNeighborsVector< TrackerTraits > CellTracksVector< TrackerTraits > HitsConstView< TrackerTraits > OuterHitOfCell< TrackerTraits > int nActualPairs
ALPAKA_FN_ACC ALPAKA_FN_INLINE void uint32_t const uint32_t CACellT< TrackerTraits > uint32_t CellNeighborsVector< TrackerTraits > CellTracksVector< TrackerTraits > HitsConstView< TrackerTraits > hh
void buildDoublets(const HitsConstView &hh, uint32_t offsetBPIX2, Queue &queue)
caStructures::OuterHitOfCellT< TrackerTraits > OuterHitOfCell
Definition: CAFishbone.h:31
constexpr auto getDoubletsFromHistoMaxBlockSize
caStructures::CellTracksT< TrackerTraits > CellTracks
Definition: CAFishbone.h:25
void launchKernels(const HitsConstView &hh, uint32_t offsetBPIX2, TkSoAView &track_view, Queue &queue)
std::enable_if_t< not std::is_array_v< T >, host_view< T > > make_host_view(T &data)
Definition: memory.h:153
void classifyTuples(const HitsConstView &hh, TkSoAView &track_view, Queue &queue)
ALPAKA_FN_ACC ALPAKA_FN_INLINE void uint32_t const uint32_t CACellT< TrackerTraits > uint32_t CellNeighborsVector< TrackerTraits > CellTracksVector< TrackerTraits > HitsConstView< TrackerTraits > OuterHitOfCell< TrackerTraits > isOuterHitOfCell