CMS 3D CMS Logo

All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
PixelRecHitKernels.dev.cc
Go to the documentation of this file.
1 // C++ headers
2 #include <algorithm>
3 #include <numeric>
4 
5 // Alpaka headers
6 #include <alpaka/alpaka.hpp>
7 
8 // CMSSW headers
13 
14 #include "PixelRecHitKernel.h"
15 #include "PixelRecHits.h"
16 
17 //#define GPU_DEBUG
18 
20  using namespace cms::alpakatools;
21  template <typename TrackerTraits>
23  public:
24  template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
25  ALPAKA_FN_ACC void operator()(TAcc const& acc,
26  uint32_t const* __restrict__ hitsModuleStart,
27  pixelCPEforDevice::ParamsOnDeviceT<TrackerTraits> const* __restrict__ cpeParams,
28  uint32_t* __restrict__ hitsLayerStart) const {
29  assert(0 == hitsModuleStart[0]);
30 
32  hitsLayerStart[i] = hitsModuleStart[cpeParams->layerGeometry().layerStart[i]];
33 #ifdef GPU_DEBUG
34  int old = i == 0 ? 0 : hitsModuleStart[cpeParams->layerGeometry().layerStart[i - 1]];
35  printf("LayerStart %d/%d at module %d: %d - %d\n",
36  i,
38  cpeParams->layerGeometry().layerStart[i],
39  hitsLayerStart[i],
40  hitsLayerStart[i] - old);
41 #endif
42  }
43  }
44  };
45 
46  namespace pixelgpudetails {
47 
48  template <typename TrackerTraits>
50  SiPixelDigisSoACollection const& digis_d,
51  SiPixelClustersSoACollection const& clusters_d,
52  BeamSpotPOD const* bs_d,
54  Queue queue) const {
55  using namespace pixelRecHits;
56  auto nHits = clusters_d.nClusters();
57  auto offsetBPIX2 = clusters_d.offsetBPIX2();
58 
59  TrackingRecHitsSoACollection<TrackerTraits> hits_d(nHits, offsetBPIX2, clusters_d->clusModuleStart(), queue);
60 
61  int activeModulesWithDigis = digis_d.nModules();
62 
63  // protect from empty events
64  if (activeModulesWithDigis) {
65  int threadsPerBlock = 128;
66  int blocks = activeModulesWithDigis;
67  const auto workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(blocks, threadsPerBlock);
68 
69 #ifdef GPU_DEBUG
70  std::cout << "launching GetHits kernel on " << alpaka::core::demangled<Acc1D> << " with " << blocks << " blocks"
71  << std::endl;
72 #endif
73  alpaka::exec<Acc1D>(queue,
74  workDiv1D,
75  GetHits<TrackerTraits>{},
76  cpeParams,
77  bs_d,
78  digis_d.view(),
79  digis_d.nDigis(),
80  clusters_d.view(),
81  hits_d.view());
82 #ifdef GPU_DEBUG
84 #endif
85 
86  // assuming full warp of threads is better than a smaller number...
87  if (nHits) {
88  const auto workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(1, 32);
89  alpaka::exec<Acc1D>(queue,
90  workDiv1D,
92  clusters_d->clusModuleStart(),
93  cpeParams,
94  hits_d.view().hitsLayerStart().data());
96 
97  // Use a view since it's runtime sized and can't use the implicit definition
98  // see HeterogeneousCore/AlpakaInterface/interface/OneToManyAssoc.h:100
100  hrv_d.assoc = &(hits_d.view().phiBinner());
101  hrv_d.offSize = -1;
102  hrv_d.offStorage = nullptr;
103  hrv_d.contentSize = nHits;
104  hrv_d.contentStorage = hits_d.view().phiBinnerStorage();
105 
106  // fillManyFromVector<Acc1D>(h_d.data(), nParts, v_d.data(), offsets_d.data(), offsets[10], 256, queue);
107  /* cms::alpakatools::fillManyFromVector<Acc1D>(&(hits_d.view().phiBinner()),
108  nLayers,
109  hits_d.view().iphi(),
110  hits_d.view().hitsLayerStart().data(),
111  nHits,
112  (uint32_t)256,
113  queue);
114 */
115  cms::alpakatools::fillManyFromVector<Acc1D>(&(hits_d.view().phiBinner()),
116  hrv_d,
117  nLayers,
118  hits_d.view().iphi(),
119  hits_d.view().hitsLayerStart().data(),
120  nHits,
121  (uint32_t)256,
122  queue);
123 
124 #ifdef GPU_DEBUG
126 #endif
127  }
128  }
129 
130 #ifdef GPU_DEBUG
132  std::cout << "PixelRecHitKernel -> DONE!" << std::endl;
133 #endif
134 
135  return hits_d;
136  }
137 
141 
142  } // namespace pixelgpudetails
143 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
ALPAKA_FN_ACC void operator()(TAcc const &acc, uint32_t const *__restrict__ hitsModuleStart, pixelCPEforDevice::ParamsOnDeviceT< TrackerTraits > const *__restrict__ cpeParams, uint32_t *__restrict__ hitsLayerStart) const
typename PhiBinner::View PhiBinnerView
constexpr uint32_t numberOfLayers
assert(be >=bs)
std::conditional_t< std::is_same_v< Device, alpaka::DevCpu >, SiPixelDigisHost, SiPixelDigisDevice< Device > > SiPixelDigisSoACollection
constexpr LayerGeometry const &__restrict__ layerGeometry() const
std::conditional_t< std::is_same_v< Device, alpaka::DevCpu >, TrackingRecHitHost< TrackerTraits >, TrackingRecHitDevice< TrackerTraits, Device > > TrackingRecHitsSoACollection
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
std::conditional_t< std::is_same_v< Device, alpaka::DevCpu >, SiPixelClustersHost, SiPixelClustersDevice< Device > > SiPixelClustersSoACollection