File indexing completed on 2024-12-12 03:12:13
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
0013 template <typename TrackerTraits, typename TDev>
0014 class TrackingRecHitDevice : public PortableDeviceCollection<TrackingRecHitLayout<TrackerTraits>, TDev> {
0015 public:
0016 using hitSoA = TrackingRecHitSoA<TrackerTraits>;
0017
0018
0019 using PortableDeviceCollection<TrackingRecHitLayout<TrackerTraits>, TDev>::view;
0020 using PortableDeviceCollection<TrackingRecHitLayout<TrackerTraits>, TDev>::const_view;
0021 using PortableDeviceCollection<TrackingRecHitLayout<TrackerTraits>, TDev>::buffer;
0022
0023 TrackingRecHitDevice(edm::Uninitialized)
0024 : PortableDeviceCollection<TrackingRecHitLayout<TrackerTraits>, TDev>{edm::kUninitialized} {}
0025
0026
0027 template <typename TQueue>
0028 explicit TrackingRecHitDevice(TQueue queue, uint32_t nHits, int32_t offsetBPIX2, uint32_t const* hitsModuleStart)
0029 : PortableDeviceCollection<TrackingRecHitLayout<TrackerTraits>, TDev>(nHits, queue), offsetBPIX2_{offsetBPIX2} {
0030 auto start_h = cms::alpakatools::make_device_view(queue, hitsModuleStart, TrackerTraits::numberOfModules + 1);
0031 auto start_d =
0032 cms::alpakatools::make_device_view(queue, view().hitsModuleStart().data(), TrackerTraits::numberOfModules + 1);
0033 alpaka::memcpy(queue, start_d, start_h);
0034
0035 auto off_h = cms::alpakatools::make_host_view(offsetBPIX2_);
0036 auto off_d = cms::alpakatools::make_device_view(queue, view().offsetBPIX2());
0037 alpaka::memcpy(queue, off_d, off_h);
0038 }
0039
0040 uint32_t nHits() const { return view().metadata().size(); }
0041
0042 int32_t offsetBPIX2() const { return offsetBPIX2_; }
0043
0044 uint32_t const* hitsModuleStart() const { return view().hitsModuleStart().data(); }
0045
0046
0047 template <typename TQueue>
0048 void updateFromDevice(TQueue queue) {
0049 auto off_h = cms::alpakatools::make_host_view(offsetBPIX2_);
0050 auto off_d = cms::alpakatools::make_device_view(queue, view().offsetBPIX2());
0051 alpaka::memcpy(queue, off_h, off_d);
0052 }
0053
0054 private:
0055
0056 int32_t offsetBPIX2_ = 0;
0057 };
0058
0059 #endif