Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef DataFormats_TrackingRecHitSoA_interface_TrackingRecHitsHost_h
0002 #define DataFormats_TrackingRecHitSoA_interface_TrackingRecHitsHost_h
0003 
0004 #include <cstdint>
0005 
0006 #include <alpaka/alpaka.hpp>
0007 
0008 #include "DataFormats/Common/interface/Uninitialized.h"
0009 #include "DataFormats/Portable/interface/PortableHostCollection.h"
0010 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsSoA.h"
0011 #include "DataFormats/SiPixelClusterSoA/interface/SiPixelClustersHost.h"
0012 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0013 
0014 // TODO: The class is created via inheritance of the PortableCollection.
0015 // This is generally discouraged, and should be done via composition.
0016 // See: https://github.com/cms-sw/cmssw/pull/40465#discussion_r1067364306
0017 
0018 namespace reco {
0019 
0020   using HitPortableCollectionHost = PortableHostMultiCollection<reco::TrackingRecHitSoA, reco::HitModuleSoA>;
0021 
0022   class TrackingRecHitHost : public HitPortableCollectionHost {
0023   public:
0024     TrackingRecHitHost(edm::Uninitialized)
0025         : PortableHostMultiCollection<reco::TrackingRecHitSoA, reco::HitModuleSoA>{edm::kUninitialized} {}
0026 
0027     // Constructor which specifies only the SoA size, to be used when copying the results from the device to the host
0028     template <typename TQueue>
0029     explicit TrackingRecHitHost(TQueue queue, uint32_t nHits, uint32_t nModules)
0030         : HitPortableCollectionHost({{int(nHits), int(nModules + 1)}}, queue) {}
0031     // Why this +1? See TrackingRecHitDevice.h constructor for an explanation
0032 
0033     // Constructor from clusters
0034     template <typename TQueue>
0035     explicit TrackingRecHitHost(TQueue queue, SiPixelClustersHost const &clusters)
0036         : HitPortableCollectionHost({{int(clusters.nClusters()), clusters.view().metadata().size()}}, queue) {
0037       auto hitsView = this->template view<TrackingRecHitSoA>();
0038       auto modsView = this->template view<HitModuleSoA>();
0039 
0040       auto nModules = clusters.view().metadata().size();
0041 
0042       auto clusters_m = cms::alpakatools::make_host_view(clusters.view().clusModuleStart(), nModules);
0043       auto hits_m = cms::alpakatools::make_host_view(modsView.moduleStart(), nModules);
0044 
0045       alpaka::memcpy(queue, hits_m, clusters_m);
0046 
0047       hitsView.offsetBPIX2() = clusters.offsetBPIX2();
0048     }
0049 
0050     uint32_t nHits() const { return this->template view<TrackingRecHitSoA>().metadata().size(); }
0051     uint32_t nModules() const { return this->template view<HitModuleSoA>().metadata().size() - 1; }
0052 
0053     int32_t offsetBPIX2() const { return this->template view<TrackingRecHitSoA>().offsetBPIX2(); }
0054 
0055     // do nothing for a host collection
0056     template <typename TQueue>
0057     void updateFromDevice(TQueue) {}
0058   };
0059 
0060 }  // namespace reco
0061 
0062 #endif  // DataFormats_TrackingRecHitSoA_interface_TrackingRecHitsHost_h