Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-09 02:22:22

0001 #include <cstdint>
0002 #include <iostream>
0003 #include <random>
0004 #include <vector>
0005 
0006 #include <alpaka/alpaka.hpp>
0007 
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/Framework/interface/Frameworkfwd.h"
0010 #include "FWCore/Framework/interface/global/EDAnalyzer.h"
0011 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 #include "FWCore/ServiceRegistry/interface/Service.h"
0015 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0016 #include "HeterogeneousCore/AlpakaInterface/interface/memory.h"
0017 #include "HeterogeneousCore/AlpakaServices/interface/alpaka/AlpakaService.h"
0018 #include "HeterogeneousTest/AlpakaWrapper/interface/alpaka/DeviceAdditionWrapper.h"
0019 
0020 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0021 
0022   class AlpakaTestWrapperAdditionModule : public edm::global::EDAnalyzer<> {
0023   public:
0024     explicit AlpakaTestWrapperAdditionModule(edm::ParameterSet const& config);
0025     ~AlpakaTestWrapperAdditionModule() override = default;
0026 
0027     static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0028 
0029     void analyze(edm::StreamID, edm::Event const& event, edm::EventSetup const& setup) const override;
0030 
0031   private:
0032     const uint32_t size_;
0033   };
0034 
0035   AlpakaTestWrapperAdditionModule::AlpakaTestWrapperAdditionModule(edm::ParameterSet const& config)
0036       : size_(config.getParameter<uint32_t>("size")) {}
0037 
0038   void AlpakaTestWrapperAdditionModule::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0039     edm::ParameterSetDescription desc;
0040     desc.add<uint32_t>("size", 1024 * 1024);
0041 
0042     // ignore the alpaka = cms.untracked.PSet(...) injected by the framework
0043     edm::ParameterSetDescription alpaka;
0044     alpaka.setAllowAnything();
0045     desc.addUntracked<edm::ParameterSetDescription>("alpaka", alpaka);
0046 
0047     descriptions.addWithDefaultLabel(desc);
0048   }
0049 
0050   void AlpakaTestWrapperAdditionModule::analyze(edm::StreamID,
0051                                                 edm::Event const& event,
0052                                                 edm::EventSetup const& setup) const {
0053     // require a valid Alpaka backend for running
0054     edm::Service<ALPAKA_TYPE_ALIAS(AlpakaService)> service;
0055     if (not service or not service->enabled()) {
0056       std::cout << "The " << ALPAKA_TYPE_ALIAS_NAME(AlpakaService)
0057                 << " is not available or disabled, the test will be skipped.\n";
0058       return;
0059     }
0060 
0061     // random number generator with a gaussian distribution
0062     std::random_device rd{};
0063     std::default_random_engine rand{rd()};
0064     std::normal_distribution<float> dist{0., 1.};
0065 
0066     // tolerance
0067     constexpr float epsilon = 0.000001;
0068 
0069     // allocate input and output host buffers
0070     std::vector<float> in1_h(size_);
0071     std::vector<float> in2_h(size_);
0072     std::vector<float> out_h(size_);
0073 
0074     // fill the input buffers with random data, and the output buffer with zeros
0075     for (uint32_t i = 0; i < size_; ++i) {
0076       in1_h[i] = dist(rand);
0077       in2_h[i] = dist(rand);
0078       out_h[i] = 0.;
0079     }
0080 
0081     // run the test on all available devices
0082     for (auto const& device : cms::alpakatools::devices<Platform>()) {
0083       Queue queue{device};
0084 
0085       // allocate input and output buffers on the device
0086       auto in1_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
0087       auto in2_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
0088       auto out_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
0089 
0090       // copy the input data to the device
0091       // FIXME: pass the explicit size of type uint32_t to avoid compilation error
0092       // The destination view and the extent are required to have compatible index types!
0093       alpaka::memcpy(queue, in1_d, in1_h, size_);
0094       alpaka::memcpy(queue, in2_d, in2_h, size_);
0095 
0096       // fill the output buffer with zeros
0097       alpaka::memset(queue, out_d, 0);
0098 
0099       // launch the 1-dimensional kernel for vector addition
0100       test::wrapper_add_vectors_f(queue, in1_d.data(), in2_d.data(), out_d.data(), size_);
0101 
0102       // copy the results from the device to the host
0103       alpaka::memcpy(queue, out_h, out_d);
0104 
0105       // wait for all the operations to complete
0106       alpaka::wait(queue);
0107 
0108       // check the results
0109       for (uint32_t i = 0; i < size_; ++i) {
0110         float sum = in1_h[i] + in2_h[i];
0111         assert(out_h[i] < sum + epsilon);
0112         assert(out_h[i] > sum - epsilon);
0113       }
0114     }
0115 
0116     std::cout << "All tests passed.\n";
0117   }
0118 
0119 }  // namespace ALPAKA_ACCELERATOR_NAMESPACE
0120 
0121 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/MakerMacros.h"
0122 DEFINE_FWK_ALPAKA_MODULE(AlpakaTestWrapperAdditionModule);