File indexing completed on 2025-07-03 00:42:11
0001 #ifndef DataFormats_TrackingRecHitSoA_interface_alpaka_TrackingRecHitsSoACollection_h
0002 #define DataFormats_TrackingRecHitSoA_interface_alpaka_TrackingRecHitsSoACollection_h
0003
0004
0005
0006 #include <cstdint>
0007 #include <type_traits>
0008
0009 #include <alpaka/alpaka.hpp>
0010
0011 #include "DataFormats/Portable/interface/alpaka/PortableCollection.h"
0012 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsDevice.h"
0013 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsHost.h"
0014 #include "DataFormats/TrackingRecHitSoA/interface/TrackingRecHitsSoA.h"
0015 #include "Geometry/CommonTopologies/interface/SimplePixelTopology.h"
0016 #include "HeterogeneousCore/AlpakaInterface/interface/CopyToHost.h"
0017 #include "HeterogeneousCore/AlpakaInterface/interface/CopyToDevice.h"
0018 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0019
0020 namespace ALPAKA_ACCELERATOR_NAMESPACE::reco {
0021
0022 using TrackingRecHitsSoACollection = std::conditional_t<std::is_same_v<Device, alpaka::DevCpu>,
0023 ::reco::TrackingRecHitHost,
0024 ::reco::TrackingRecHitDevice<Device>>;
0025 }
0026
0027 namespace cms::alpakatools {
0028 template <typename TDevice>
0029 struct CopyToHost<::reco::TrackingRecHitDevice<TDevice>> {
0030 template <typename TQueue>
0031 static auto copyAsync(TQueue& queue, ::reco::TrackingRecHitDevice<TDevice> const& deviceData) {
0032 auto nHits = deviceData.nHits();
0033
0034 reco::TrackingRecHitHost hostData(queue, nHits, deviceData.nModules());
0035
0036
0037 if (nHits == 0) {
0038 std::memset(
0039 hostData.buffer().data(),
0040 0,
0041 alpaka::getExtentProduct(hostData.buffer()) * sizeof(alpaka::Elem<reco::TrackingRecHitHost::Buffer>));
0042 return hostData;
0043 }
0044
0045 alpaka::memcpy(queue, hostData.buffer(), deviceData.buffer());
0046 #ifdef GPU_DEBUG
0047 printf("TrackingRecHitsSoACollection: I'm copying to host.\n");
0048 alpaka::wait(queue);
0049 assert(deviceData.nHits() == hostData.nHits());
0050 assert(deviceData.nModules() == hostData.nModules());
0051 assert(deviceData.offsetBPIX2() == hostData.offsetBPIX2());
0052 #endif
0053
0054 return hostData;
0055 }
0056 };
0057
0058 template <>
0059 struct CopyToDevice<::reco::TrackingRecHitHost> {
0060 template <typename TQueue>
0061 static auto copyAsync(TQueue& queue, reco::TrackingRecHitHost const& hostData) {
0062 using TDevice = typename alpaka::trait::DevType<TQueue>::type;
0063
0064 auto nHits = hostData.nHits();
0065
0066 reco::TrackingRecHitDevice<TDevice> deviceData(queue, nHits, hostData.nModules());
0067
0068 if (nHits == 0) {
0069 std::memset(
0070 deviceData.buffer().data(),
0071 0,
0072 alpaka::getExtentProduct(deviceData.buffer()) * sizeof(alpaka::Elem<reco::TrackingRecHitHost::Buffer>));
0073 return deviceData;
0074 }
0075
0076 alpaka::memcpy(queue, deviceData.buffer(), hostData.buffer());
0077
0078 #ifdef GPU_DEBUG
0079 printf("TrackingRecHitsSoACollection: I'm copying to device.\n");
0080 alpaka::wait(queue);
0081 assert(deviceData.nHits() == hostData.nHits());
0082 assert(deviceData.nModules() == hostData.nModules());
0083 assert(deviceData.offsetBPIX2() == hostData.offsetBPIX2());
0084 #endif
0085 return deviceData;
0086 }
0087 };
0088
0089 }
0090
0091 ASSERT_DEVICE_MATCHES_HOST_COLLECTION(reco::TrackingRecHitsSoACollection, reco::TrackingRecHitHost);
0092
0093 #endif