CMS 3D CMS Logo

AlpakaTestDeviceAdditionModule.cc
Go to the documentation of this file.
1 #include <cstdint>
2 #include <iostream>
3 #include <random>
4 #include <vector>
5 
6 #include <alpaka/alpaka.hpp>
7 
18 
20 
22 
24  public:
26  ~AlpakaTestDeviceAdditionModule() override = default;
27 
28  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
29 
30  void analyze(edm::StreamID, edm::Event const& event, edm::EventSetup const& setup) const override;
31 
32  private:
33  const uint32_t size_;
34  };
35 
37  : size_(config.getParameter<uint32_t>("size")) {}
38 
41  desc.add<uint32_t>("size", 1024 * 1024);
42 
43  // ignore the alpaka = cms.untracked.PSet(...) injected by the framework
45  alpaka.setAllowAnything();
46  desc.addUntracked<edm::ParameterSetDescription>("alpaka", alpaka);
47 
48  descriptions.addWithDefaultLabel(desc);
49  }
50 
52  edm::Event const& event,
53  edm::EventSetup const& setup) const {
54  // require a valid Alpaka backend for running
56  if (not service or not service->enabled()) {
57  std::cout << "The " << ALPAKA_TYPE_ALIAS_NAME(AlpakaService)
58  << " is not available or disabled, the test will be skipped.\n";
59  return;
60  }
61 
62  // random number generator with a gaussian distribution
63  std::random_device rd{};
64  std::default_random_engine rand{rd()};
65  std::normal_distribution<float> dist{0., 1.};
66 
67  // tolerance
68  constexpr float epsilon = 0.000001;
69 
70  // allocate input and output host buffers
71  std::vector<float> in1_h(size_);
72  std::vector<float> in2_h(size_);
73  std::vector<float> out_h(size_);
74 
75  // fill the input buffers with random data, and the output buffer with zeros
76  for (uint32_t i = 0; i < size_; ++i) {
77  in1_h[i] = dist(rand);
78  in2_h[i] = dist(rand);
79  out_h[i] = 0.;
80  }
81 
82  // run the test on all available devices
83  for (auto const& device : cms::alpakatools::devices<Platform>()) {
84  Queue queue{device};
85 
86  // allocate input and output buffers on the device
87  auto in1_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
88  auto in2_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
89  auto out_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
90 
91  // copy the input data to the device
92  // FIXME: pass the explicit size of type uint32_t to avoid compilation error
93  // The destination view and the extent are required to have compatible index types!
94  alpaka::memcpy(queue, in1_d, in1_h, size_);
95  alpaka::memcpy(queue, in2_d, in2_h, size_);
96 
97  // fill the output buffer with zeros
98  alpaka::memset(queue, out_d, 0);
99 
100  // launch the 1-dimensional kernel for vector addition
102  queue, in1_d.data(), in2_d.data(), out_d.data(), size_);
103 
104  // copy the results from the device to the host
105  alpaka::memcpy(queue, out_h, out_d);
106 
107  // wait for all the operations to complete
109 
110  // check the results
111  for (uint32_t i = 0; i < size_; ++i) {
112  float sum = in1_h[i] + in2_h[i];
113  assert(out_h[i] < sum + epsilon);
114  assert(out_h[i] > sum - epsilon);
115  }
116  }
117 
118  std::cout << "All tests passed.\n";
119  }
120 
121 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
122 
124 DEFINE_FWK_ALPAKA_MODULE(AlpakaTestDeviceAdditionModule);
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: config.py:1
assert(be >=bs)
ALPAKA_ACCELERATOR_NAMESPACE::Queue Queue
Definition: LSTEvent.dev.cc:14
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
void analyze(edm::StreamID, edm::Event const &event, edm::EventSetup const &setup) const override
#define DEFINE_FWK_ALPAKA_MODULE(name)
Definition: MakerMacros.h:16
void wrapper_add_vectors_f(Queue &queue, const float *__restrict__ in1, const float *__restrict__ in2, float *__restrict__ out, uint32_t size)
Definition: event.py:1