File indexing completed on 2025-07-03 00:42:11
0001 #ifndef DataFormats_TrackingRecHitSoA_interface_TrackingRecHitSoADevice_h
0002 #define DataFormats_TrackingRecHitSoA_interface_TrackingRecHitSoADevice_h
0003
0004 #include <cstdint>
0005
0006 #include <alpaka/alpaka.hpp>
0007
0008 #include "DataFormats/Common/interface/Uninitialized.h"
0009 #include "DataFormats/Portable/interface/PortableDeviceCollection.h"
0010 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsHost.h"
0011 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsSoA.h"
0012 #include "DataFormats/SiPixelClusterSoA/interface/SiPixelClustersDevice.h"
0013
0014
0015
0016
0017
0018 namespace reco {
0019
0020 template <typename TDev>
0021 using HitPortableCollectionDevice = PortableDeviceMultiCollection<TDev, reco::TrackingRecHitSoA, reco::HitModuleSoA>;
0022
0023 template <typename TDev>
0024 class TrackingRecHitDevice : public HitPortableCollectionDevice<TDev> {
0025 public:
0026 TrackingRecHitDevice() = default;
0027
0028 TrackingRecHitDevice(edm::Uninitialized) : HitPortableCollectionDevice<TDev>{edm::kUninitialized} {}
0029
0030
0031 template <typename TQueue>
0032 explicit TrackingRecHitDevice(TQueue queue, uint32_t nHits, uint32_t nModules)
0033 : HitPortableCollectionDevice<TDev>({{int(nHits), int(nModules + 1)}}, queue) {}
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043 template <typename TQueue>
0044 explicit TrackingRecHitDevice(TQueue queue, SiPixelClustersDevice<TDev> const &clusters)
0045 : HitPortableCollectionDevice<TDev>({{int(clusters.nClusters()), clusters.view().metadata().size()}}, queue),
0046 offsetBPIX2_{clusters.offsetBPIX2()} {
0047 auto hitsView = this->template view<TrackingRecHitSoA>();
0048 auto modsView = this->template view<HitModuleSoA>();
0049
0050 auto nModules = clusters.view().metadata().size();
0051
0052 auto clusters_m = cms::alpakatools::make_device_view(queue, clusters.view().clusModuleStart(), nModules);
0053 auto hits_m = cms::alpakatools::make_device_view(queue, modsView.moduleStart(), nModules);
0054
0055 alpaka::memcpy(queue, hits_m, clusters_m);
0056
0057 auto off_h = cms::alpakatools::make_host_view(offsetBPIX2_);
0058 auto off_d = cms::alpakatools::make_device_view(queue, hitsView.offsetBPIX2());
0059 alpaka::memcpy(queue, off_d, off_h);
0060 }
0061
0062 uint32_t nHits() const { return this->template view<TrackingRecHitSoA>().metadata().size(); }
0063 uint32_t nModules() const { return this->template view<HitModuleSoA>().metadata().size() - 1; }
0064
0065 int32_t offsetBPIX2() const { return offsetBPIX2_; }
0066
0067
0068 template <typename TQueue>
0069 void updateFromDevice(TQueue queue) {
0070 auto off_h = cms::alpakatools::make_host_view(offsetBPIX2_);
0071 auto off_d = cms::alpakatools::make_device_view(queue, this->template view<TrackingRecHitSoA>().offsetBPIX2());
0072 alpaka::memcpy(queue, off_h, off_d);
0073 }
0074
0075 private:
0076
0077 int32_t offsetBPIX2_ = 0;
0078 };
0079 }
0080
0081 #endif