Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-07-03 00:42:36

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 
0030   using namespace ALPAKA_ACCELERATOR_NAMESPACE::reco;
0031 
0032   namespace pixelgpudetails {
0033 
0034     template <typename TrackerTraits>
0035     TrackingRecHitsSoACollection PixelRecHitKernel<TrackerTraits>::makeHitsAsync(
0036         SiPixelDigisSoACollection const& digis_d,
0037         SiPixelClustersSoACollection const& clusters_d,
0038         BeamSpotPOD const* bs_d,
0039         pixelCPEforDevice::ParamsOnDeviceT<TrackerTraits> const* cpeParams,
0040         Queue queue) const {
0041       using namespace pixelRecHits;
0042 
0043       TrackingRecHitsSoACollection hits_d(queue, clusters_d);
0044 
0045       int activeModulesWithDigis = digis_d.nModules();
0046 
0047       // protect from empty events
0048       if (activeModulesWithDigis) {
0049         int threadsPerBlock = 128;
0050         // note: the kernel should work with an arbitrary number of blocks
0051         int blocks = activeModulesWithDigis;
0052         const auto workDiv1D = cms::alpakatools::make_workdiv<Acc1D>(blocks, threadsPerBlock);
0053 
0054 #ifdef GPU_DEBUG
0055         std::cout << "launching GetHits kernel on " << alpaka::core::demangled<Acc1D> << " with " << blocks << " blocks"
0056                   << std::endl;
0057 #endif
0058         alpaka::exec<Acc1D>(queue,
0059                             workDiv1D,
0060                             GetHits<TrackerTraits>{},
0061                             cpeParams,
0062                             bs_d,
0063                             digis_d.view(),
0064                             digis_d.nDigis(),
0065                             digis_d.nModules(),
0066                             clusters_d.view(),
0067                             hits_d.view());
0068 #ifdef GPU_DEBUG
0069         alpaka::wait(queue);
0070 #endif
0071       }
0072 
0073 #ifdef GPU_DEBUG
0074       alpaka::wait(queue);
0075       std::cout << "PixelRecHitKernel -> DONE!" << std::endl;
0076 #endif
0077 
0078       return hits_d;
0079     }
0080 
0081     template class PixelRecHitKernel<pixelTopology::Phase1>;
0082     template class PixelRecHitKernel<pixelTopology::Phase2>;
0083     template class PixelRecHitKernel<pixelTopology::HIonPhase1>;
0084 
0085   }  // namespace pixelgpudetails
0086 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE