CMS 3D CMS Logo

getDeviceCachingAllocator.h
Go to the documentation of this file.
1 #ifndef HeterogeneousCore_AlpakaInterface_interface_getDeviceCachingAllocator_h
2 #define HeterogeneousCore_AlpakaInterface_interface_getDeviceCachingAllocator_h
3 
4 #include <cassert>
5 #include <memory>
6 
7 #include <alpaka/alpaka.hpp>
8 
14 
15 namespace cms::alpakatools {
16 
17  namespace detail {
18 
19  template <typename TDev,
20  typename TQueue,
21  typename = std::enable_if_t<alpaka::isDevice<TDev> and alpaka::isQueue<TQueue>>>
23  using Allocator = CachingAllocator<TDev, TQueue>;
24  auto const& devices = cms::alpakatools::devices<alpaka::Platform<TDev>>();
25  ssize_t const size = devices.size();
26 
27  // allocate the storage for the objects
28  auto ptr = std::allocator<Allocator>().allocate(size);
29 
30  // construct the objects in the storage
31  ptrdiff_t index = 0;
32  try {
33  for (; index < size; ++index) {
34 #if __cplusplus >= 202002L
35  std::construct_at(
36 #else
37  std::allocator<Allocator>().construct(
38 #endif
39  ptr + index,
40  devices[index],
41  config,
42  true, // reuseSameQueueAllocations
43  debug);
44  }
45  } catch (...) {
46  --index;
47  // destroy any object that had been succesfully constructed
48  while (index >= 0) {
49  std::destroy_at(ptr + index);
50  --index;
51  }
52  // deallocate the storage
53  std::allocator<Allocator>().deallocate(ptr, size);
54  // rethrow the exception
55  throw;
56  }
57 
58  // use a custom deleter to destroy all objects and deallocate the memory
59  auto deleter = [size](Allocator* allocators) {
60  for (size_t i = size; i > 0; --i) {
61  std::destroy_at(allocators + i - 1);
62  }
63  std::allocator<Allocator>().deallocate(allocators, size);
64  };
65 
66  return std::unique_ptr<Allocator[], decltype(deleter)>(ptr, deleter);
67  }
68 
69  } // namespace detail
70 
71  template <typename TDev,
72  typename TQueue,
73  typename = std::enable_if_t<alpaka::isDevice<TDev> and alpaka::isQueue<TQueue>>>
76  bool debug = false) {
77  // initialise all allocators, one per device
78  CMS_THREAD_SAFE static auto allocators = detail::allocate_device_allocators<TDev, TQueue>(config, debug);
79 
80  size_t const index = alpaka::getNativeHandle(device);
81  assert(index < cms::alpakatools::devices<alpaka::Platform<TDev>>().size());
82 
83  // the public interface is thread safe
84  return allocators[index];
85  }
86 
87 } // namespace cms::alpakatools
88 
89 #endif // HeterogeneousCore_AlpakaInterface_interface_getDeviceCachingAllocator_h
def config(tmpl, pkg_help)
Definition: cms.py:19
Definition: config.py:1
assert(be >=bs)
#define CMS_THREAD_SAFE
auto allocate_device_allocators(AllocatorConfig const &config, bool debug)
#define debug
Definition: HDRShower.cc:19
std::vector< alpaka::Dev< TPlatform > > const & devices()
Definition: devices.h:22
CachingAllocator< TDev, TQueue > & getDeviceCachingAllocator(TDev const &device, AllocatorConfig const &config=AllocatorConfig{}, bool debug=false)
std::unique_ptr< GeometricDet > construct(DDCompactView const &cpv, std::vector< int > const &detidShifts)