Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-22 07:34:24

0001 // C++ headers
0002 #include <cassert>
0003 #include <cstdint>
0004 #include <type_traits>
0005 
0006 // Alpaka headers
0007 #include <alpaka/alpaka.hpp>
0008 
0009 // CMSSW headers
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 // local headers
0022 #include "PixelRecHitKernel.h"
0023 #include "PixelRecHits.h"
0024 
0025 //#define GPU_DEBUG
0026 
0027 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0028   using namespace cms::alpakatools;
0029   template <typename TrackerTraits>
0030   class setHitsLayerStart {
0031   public:
0032     ALPAKA_FN_ACC void operator()(Acc1D const& acc,
0033                                   uint32_t const* __restrict__ hitsModuleStart,
0034                                   pixelCPEforDevice::ParamsOnDeviceT<TrackerTraits> const* __restrict__ cpeParams,
0035                                   uint32_t* __restrict__ hitsLayerStart) const {
0036       ALPAKA_ASSERT_ACC(0 == hitsModuleStart[0]);
0037 
0038       for (int32_t i : cms::alpakatools::uniform_elements(acc, TrackerTraits::numberOfLayers + 1)) {
0039         hitsLayerStart[i] = hitsModuleStart[cpeParams->layerGeometry().layerStart[i]];
0040 #ifdef GPU_DEBUG
0041         int old = i == 0 ? 0 : hitsModuleStart[cpeParams->layerGeometry().layerStart[i - 1]];
0042         printf("LayerStart %d/%d at module %d: %d - %d\n",
0043                i,
0044                TrackerTraits::numberOfLayers,
0045                cpeParams->layerGeometry().layerStart[i],
0046                hitsLayerStart[i],
0047                hitsLayerStart[i] - old);
0048 #endif
0049       }
0050     }
0051   };
0052 
0053   namespace pixelgpudetails {
0054 
0055     template <typename TrackerTraits>
0056     TrackingRecHitsSoACollection<TrackerTraits> PixelRecHitKernel<TrackerTraits>::makeHitsAsync(
0057         SiPixelDigisSoACollection const& digis_d,
0058         SiPixelClustersSoACollection const& clusters_d,
0059         BeamSpotPOD const* bs_d,
0060         pixelCPEforDevice::ParamsOnDeviceT<TrackerTraits> const* cpeParams,
0061         Queue queue) const {
0062       using namespace pixelRecHits;
0063       auto nHits = clusters_d.nClusters();
0064       auto offsetBPIX2 = clusters_d.offsetBPIX2();
0065 
0066       TrackingRecHitsSoACollection<TrackerTraits> hits_d(queue, nHits, offsetBPIX2, clusters_d->clusModuleStart());
0067 
0068       int activeModulesWithDigis = digis_d.nModules();
0069 
0070       // protect from empty events
0071       if (activeModulesWithDigis) {
0072         int threadsPerBlock = 128;
0073         // note: the kernel should work with an arbitrary number of blocks
0074         int blocks = activeModulesWithDigis;
0075         const auto workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(blocks, threadsPerBlock);
0076 
0077 #ifdef GPU_DEBUG
0078         std::cout << "launching GetHits kernel on " << alpaka::core::demangled<Acc1D> << " with " << blocks << " blocks"
0079                   << std::endl;
0080 #endif
0081         alpaka::exec<Acc1D>(queue,
0082                             workDiv1D,
0083                             GetHits<TrackerTraits>{},
0084                             cpeParams,
0085                             bs_d,
0086                             digis_d.view(),
0087                             digis_d.nDigis(),
0088                             digis_d.nModules(),
0089                             clusters_d.view(),
0090                             hits_d.view());
0091 #ifdef GPU_DEBUG
0092         alpaka::wait(queue);
0093 #endif
0094 
0095         // assuming full warp of threads is better than a smaller number...
0096         if (nHits) {
0097           const auto workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(1, 32);
0098           alpaka::exec<Acc1D>(queue,
0099                               workDiv1D,
0100                               setHitsLayerStart<TrackerTraits>{},
0101                               clusters_d->clusModuleStart(),
0102                               cpeParams,
0103                               hits_d.view().hitsLayerStart().data());
0104           constexpr auto nLayers = TrackerTraits::numberOfLayers;
0105 
0106           // Use a view since it's runtime sized and can't use the implicit definition
0107           // see HeterogeneousCore/AlpakaInterface/interface/OneToManyAssoc.h:100
0108           typename TrackingRecHitSoA<TrackerTraits>::PhiBinnerView hrv_d;
0109           hrv_d.assoc = &(hits_d.view().phiBinner());
0110           hrv_d.offSize = -1;
0111           hrv_d.offStorage = nullptr;
0112           hrv_d.contentSize = nHits;
0113           hrv_d.contentStorage = hits_d.view().phiBinnerStorage();
0114 
0115           cms::alpakatools::fillManyFromVector<Acc1D>(&(hits_d.view().phiBinner()),
0116                                                       hrv_d,
0117                                                       nLayers,
0118                                                       hits_d.view().iphi(),
0119                                                       hits_d.view().hitsLayerStart().data(),
0120                                                       nHits,
0121                                                       (uint32_t)256,
0122                                                       queue);
0123 
0124 #ifdef GPU_DEBUG
0125           alpaka::wait(queue);
0126 #endif
0127         }
0128       }
0129 
0130 #ifdef GPU_DEBUG
0131       alpaka::wait(queue);
0132       std::cout << "PixelRecHitKernel -> DONE!" << std::endl;
0133 #endif
0134 
0135       return hits_d;
0136     }
0137 
0138     template class PixelRecHitKernel<pixelTopology::Phase1>;
0139     template class PixelRecHitKernel<pixelTopology::Phase2>;
0140     template class PixelRecHitKernel<pixelTopology::HIonPhase1>;
0141 
0142   }  // namespace pixelgpudetails
0143 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE