CMS 3D CMS Logo

AlpakaTestWrapperAdditionModule.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 
19 
21 
23  public:
25  ~AlpakaTestWrapperAdditionModule() override = default;
26 
27  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
28 
29  void analyze(edm::StreamID, edm::Event const& event, edm::EventSetup const& setup) const override;
30 
31  private:
32  const uint32_t size_;
33  };
34 
36  : size_(config.getParameter<uint32_t>("size")) {}
37 
40  desc.add<uint32_t>("size", 1024 * 1024);
41 
42  // ignore the alpaka = cms.untracked.PSet(...) injected by the framework
44  alpaka.setAllowAnything();
45  desc.addUntracked<edm::ParameterSetDescription>("alpaka", alpaka);
46 
47  descriptions.addWithDefaultLabel(desc);
48  }
49 
51  edm::Event const& event,
52  edm::EventSetup const& setup) const {
53  // require a valid Alpaka backend for running
55  if (not service or not service->enabled()) {
56  std::cout << "The " << ALPAKA_TYPE_ALIAS_NAME(AlpakaService)
57  << " is not available or disabled, the test will be skipped.\n";
58  return;
59  }
60 
61  // random number generator with a gaussian distribution
62  std::random_device rd{};
63  std::default_random_engine rand{rd()};
64  std::normal_distribution<float> dist{0., 1.};
65 
66  // tolerance
67  constexpr float epsilon = 0.000001;
68 
69  // allocate input and output host buffers
70  std::vector<float> in1_h(size_);
71  std::vector<float> in2_h(size_);
72  std::vector<float> out_h(size_);
73 
74  // fill the input buffers with random data, and the output buffer with zeros
75  for (uint32_t i = 0; i < size_; ++i) {
76  in1_h[i] = dist(rand);
77  in2_h[i] = dist(rand);
78  out_h[i] = 0.;
79  }
80 
81  // run the test on all available devices
82  for (auto const& device : cms::alpakatools::devices<Platform>()) {
83  Queue queue{device};
84 
85  // allocate input and output buffers on the device
86  auto in1_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
87  auto in2_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
88  auto out_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
89 
90  // copy the input data to the device
91  // FIXME: pass the explicit size of type uint32_t to avoid compilation error
92  // The destination view and the extent are required to have compatible index types!
93  alpaka::memcpy(queue, in1_d, in1_h, size_);
94  alpaka::memcpy(queue, in2_d, in2_h, size_);
95 
96  // fill the output buffer with zeros
97  alpaka::memset(queue, out_d, 0);
98 
99  // launch the 1-dimensional kernel for vector addition
100  test::wrapper_add_vectors_f(queue, in1_d.data(), in2_d.data(), out_d.data(), size_);
101 
102  // copy the results from the device to the host
103  alpaka::memcpy(queue, out_h, out_d);
104 
105  // wait for all the operations to complete
107 
108  // check the results
109  for (uint32_t i = 0; i < size_; ++i) {
110  float sum = in1_h[i] + in2_h[i];
111  assert(out_h[i] < sum + epsilon);
112  assert(out_h[i] > sum - epsilon);
113  }
114  }
115 
116  std::cout << "All tests passed.\n";
117  }
118 
119 } // namespace ALPAKA_ACCELERATOR_NAMESPACE
120 
122 DEFINE_FWK_ALPAKA_MODULE(AlpakaTestWrapperAdditionModule);
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void wrapper_add_vectors_f(Queue &queue, const float *__restrict__ in1, const float *__restrict__ in2, float *__restrict__ out, uint32_t size)
Definition: config.py:1
assert(be >=bs)
void analyze(edm::StreamID, edm::Event const &event, edm::EventSetup const &setup) const override
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
#define DEFINE_FWK_ALPAKA_MODULE(name)
Definition: MakerMacros.h:16
Definition: event.py:1