File indexing completed on 2024-04-23 22:56:27
0001
0002 #include <cassert>
0003 #include <cstdint>
0004 #include <type_traits>
0005
0006
0007 #include <alpaka/alpaka.hpp>
0008
0009
0010 #include "DataFormats/BeamSpot/interface/BeamSpotPOD.h"
0011 #include "DataFormats/SiPixelClusterSoA/interface/alpaka/SiPixelClustersSoACollection.h"
0012 #include "DataFormats/SiPixelDigiSoA/interface/alpaka/SiPixelDigisSoACollection.h"
0013 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsSoA.h"
0014 #include "DataFormats/TrackingRecHitSoA/interface/alpaka/TrackingRecHitsSoACollection.h"
0015 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0016 #include "HeterogeneousCore/AlpakaInterface/interface/HistoContainer.h"
0017 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0018 #include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h"
0019 #include "RecoLocalTracker/SiPixelRecHits/interface/pixelCPEforDevice.h"
0020
0021
0022 #include "PixelRecHitKernel.h"
0023 #include "PixelRecHits.h"
0024
0025
0026
0027 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0028 using namespace cms::alpakatools;
0029 template <typename TrackerTraits>
0030 class setHitsLayerStart {
0031 public:
0032 template <typename TAcc, typename = std::enable_if_t<alpaka::isAccelerator<TAcc>>>
0033 ALPAKA_FN_ACC void operator()(TAcc const& acc,
0034 uint32_t const* __restrict__ hitsModuleStart,
0035 pixelCPEforDevice::ParamsOnDeviceT<TrackerTraits> const* __restrict__ cpeParams,
0036 uint32_t* __restrict__ hitsLayerStart) const {
0037 ALPAKA_ASSERT_ACC(0 == hitsModuleStart[0]);
0038
0039 for (int32_t i : cms::alpakatools::uniform_elements(acc, TrackerTraits::numberOfLayers + 1)) {
0040 hitsLayerStart[i] = hitsModuleStart[cpeParams->layerGeometry().layerStart[i]];
0041 #ifdef GPU_DEBUG
0042 int old = i == 0 ? 0 : hitsModuleStart[cpeParams->layerGeometry().layerStart[i - 1]];
0043 printf("LayerStart %d/%d at module %d: %d - %d\n",
0044 i,
0045 TrackerTraits::numberOfLayers,
0046 cpeParams->layerGeometry().layerStart[i],
0047 hitsLayerStart[i],
0048 hitsLayerStart[i] - old);
0049 #endif
0050 }
0051 }
0052 };
0053
0054 namespace pixelgpudetails {
0055
0056 template <typename TrackerTraits>
0057 TrackingRecHitsSoACollection<TrackerTraits> PixelRecHitKernel<TrackerTraits>::makeHitsAsync(
0058 SiPixelDigisSoACollection const& digis_d,
0059 SiPixelClustersSoACollection const& clusters_d,
0060 BeamSpotPOD const* bs_d,
0061 pixelCPEforDevice::ParamsOnDeviceT<TrackerTraits> const* cpeParams,
0062 Queue queue) const {
0063 using namespace pixelRecHits;
0064 auto nHits = clusters_d.nClusters();
0065 auto offsetBPIX2 = clusters_d.offsetBPIX2();
0066
0067 TrackingRecHitsSoACollection<TrackerTraits> hits_d(queue, nHits, offsetBPIX2, clusters_d->clusModuleStart());
0068
0069 int activeModulesWithDigis = digis_d.nModules();
0070
0071
0072 if (activeModulesWithDigis) {
0073 int threadsPerBlock = 128;
0074
0075 int blocks = activeModulesWithDigis;
0076 const auto workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(blocks, threadsPerBlock);
0077
0078 #ifdef GPU_DEBUG
0079 std::cout << "launching GetHits kernel on " << alpaka::core::demangled<Acc1D> << " with " << blocks << " blocks"
0080 << std::endl;
0081 #endif
0082 alpaka::exec<Acc1D>(queue,
0083 workDiv1D,
0084 GetHits<TrackerTraits>{},
0085 cpeParams,
0086 bs_d,
0087 digis_d.view(),
0088 digis_d.nDigis(),
0089 digis_d.nModules(),
0090 clusters_d.view(),
0091 hits_d.view());
0092 #ifdef GPU_DEBUG
0093 alpaka::wait(queue);
0094 #endif
0095
0096
0097 if (nHits) {
0098 const auto workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(1, 32);
0099 alpaka::exec<Acc1D>(queue,
0100 workDiv1D,
0101 setHitsLayerStart<TrackerTraits>{},
0102 clusters_d->clusModuleStart(),
0103 cpeParams,
0104 hits_d.view().hitsLayerStart().data());
0105 constexpr auto nLayers = TrackerTraits::numberOfLayers;
0106
0107
0108
0109 typename TrackingRecHitSoA<TrackerTraits>::PhiBinnerView hrv_d;
0110 hrv_d.assoc = &(hits_d.view().phiBinner());
0111 hrv_d.offSize = -1;
0112 hrv_d.offStorage = nullptr;
0113 hrv_d.contentSize = nHits;
0114 hrv_d.contentStorage = hits_d.view().phiBinnerStorage();
0115
0116 cms::alpakatools::fillManyFromVector<Acc1D>(&(hits_d.view().phiBinner()),
0117 hrv_d,
0118 nLayers,
0119 hits_d.view().iphi(),
0120 hits_d.view().hitsLayerStart().data(),
0121 nHits,
0122 (uint32_t)256,
0123 queue);
0124
0125 #ifdef GPU_DEBUG
0126 alpaka::wait(queue);
0127 #endif
0128 }
0129 }
0130
0131 #ifdef GPU_DEBUG
0132 alpaka::wait(queue);
0133 std::cout << "PixelRecHitKernel -> DONE!" << std::endl;
0134 #endif
0135
0136 return hits_d;
0137 }
0138
0139 template class PixelRecHitKernel<pixelTopology::Phase1>;
0140 template class PixelRecHitKernel<pixelTopology::Phase2>;
0141 template class PixelRecHitKernel<pixelTopology::HIonPhase1>;
0142
0143 }
0144 }