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