CMS 3D CMS Logo

TrackingRecHitsSoACollection.h
Go to the documentation of this file.
1 #ifndef DataFormats_TrackingRecHitSoA_interface_alpaka_TrackingRecHitsSoACollection_h
2 #define DataFormats_TrackingRecHitSoA_interface_alpaka_TrackingRecHitsSoACollection_h
3 
4 #include <cstdint>
5 #include <type_traits>
6 
7 #include <alpaka/alpaka.hpp>
8 
16 
18 
19  template <typename TrackerTraits>
20  using TrackingRecHitsSoACollection = std::conditional_t<std::is_same_v<Device, alpaka::DevCpu>,
23 
24  // Classes definition for Phase1/Phase2, to make the classes_def lighter. Not actually used in the code.
28 
29 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
30 
31 namespace cms::alpakatools {
32  template <typename TrackerTraits, typename TDevice>
33  struct CopyToHost<TrackingRecHitDevice<TrackerTraits, TDevice>> {
34  template <typename TQueue>
35  static auto copyAsync(TQueue& queue, TrackingRecHitDevice<TrackerTraits, TDevice> const& deviceData) {
36  TrackingRecHitHost<TrackerTraits> hostData(queue, deviceData.view().metadata().size());
37 
38  // Don't bother if zero hits
39  if (deviceData.view().metadata().size() == 0) {
40  std::memset(hostData.buffer().data(),
41  0,
42  alpaka::getExtentProduct(hostData.buffer()) *
43  sizeof(alpaka::Elem<typename TrackingRecHitHost<TrackerTraits>::Buffer>));
44  return hostData;
45  }
46 
47  alpaka::memcpy(queue, hostData.buffer(), deviceData.buffer());
48 #ifdef GPU_DEBUG
49  printf("TrackingRecHitsSoACollection: I'm copying to host.\n");
51  assert(deviceData.nHits() == hostData.nHits());
52  assert(deviceData.offsetBPIX2() == hostData.offsetBPIX2());
53 #endif
54  // Update the contents address of the phiBinner histo container after the copy from device happened
57  pbv.assoc = &(hostData.view().phiBinner());
58  pbv.offSize = -1;
59  pbv.offStorage = nullptr;
60  pbv.contentSize = hostData.nHits();
61  pbv.contentStorage = hostData.view().phiBinnerStorage();
62  hostData.view().phiBinner().initStorage(pbv);
63 
64  return hostData;
65  }
66  };
67 } // namespace cms::alpakatools
68 
72 
73 #endif // DataFormats_TrackingRecHitSoA_interface_alpaka_TrackingRecHitsSoACollection_h
typename PhiBinner::View PhiBinnerView
static auto copyAsync(TQueue &queue, TrackingRecHitDevice< TrackerTraits, TDevice > const &deviceData)
#define ASSERT_DEVICE_MATCHES_HOST_COLLECTION(DEVICE_COLLECTION, HOST_COLLECTION)
assert(be >=bs)
TrackingRecHitsSoACollection< pixelTopology::Phase2 > TrackingRecHitSoAPhase2
TrackingRecHitsSoACollection< pixelTopology::Phase1 > TrackingRecHitSoAPhase1
TrackingRecHitsSoACollection< pixelTopology::HIonPhase1 > TrackingRecHitSoAHIonPhase1
std::conditional_t< std::is_same_v< Device, alpaka::DevCpu >, TrackingRecHitHost< TrackerTraits >, TrackingRecHitDevice< TrackerTraits, Device > > TrackingRecHitsSoACollection