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
0019 #include "AlpakaTestDeviceAdditionAlgo.h"
0020
0021 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0022
0023 class AlpakaTestDeviceAdditionModule : public edm::global::EDAnalyzer<> {
0024 public:
0025 explicit AlpakaTestDeviceAdditionModule(edm::ParameterSet const& config);
0026 ~AlpakaTestDeviceAdditionModule() override = default;
0027
0028 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0029
0030 void analyze(edm::StreamID, edm::Event const& event, edm::EventSetup const& setup) const override;
0031
0032 private:
0033 const uint32_t size_;
0034 };
0035
0036 AlpakaTestDeviceAdditionModule::AlpakaTestDeviceAdditionModule(edm::ParameterSet const& config)
0037 : size_(config.getParameter<uint32_t>("size")) {}
0038
0039 void AlpakaTestDeviceAdditionModule::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0040 edm::ParameterSetDescription desc;
0041 desc.add<uint32_t>("size", 1024 * 1024);
0042
0043
0044 edm::ParameterSetDescription alpaka;
0045 alpaka.setAllowAnything();
0046 desc.addUntracked<edm::ParameterSetDescription>("alpaka", alpaka);
0047
0048 descriptions.addWithDefaultLabel(desc);
0049 }
0050
0051 void AlpakaTestDeviceAdditionModule::analyze(edm::StreamID,
0052 edm::Event const& event,
0053 edm::EventSetup const& setup) const {
0054
0055 edm::Service<ALPAKA_TYPE_ALIAS(AlpakaService)> service;
0056 if (not service or not service->enabled()) {
0057 std::cout << "The " << ALPAKA_TYPE_ALIAS_NAME(AlpakaService)
0058 << " is not available or disabled, the test will be skipped.\n";
0059 return;
0060 }
0061
0062
0063 std::random_device rd{};
0064 std::default_random_engine rand{rd()};
0065 std::normal_distribution<float> dist{0., 1.};
0066
0067
0068 constexpr float epsilon = 0.000001;
0069
0070
0071 std::vector<float> in1_h(size_);
0072 std::vector<float> in2_h(size_);
0073 std::vector<float> out_h(size_);
0074
0075
0076 for (uint32_t i = 0; i < size_; ++i) {
0077 in1_h[i] = dist(rand);
0078 in2_h[i] = dist(rand);
0079 out_h[i] = 0.;
0080 }
0081
0082
0083 for (auto const& device : cms::alpakatools::devices<Platform>()) {
0084 Queue queue{device};
0085
0086
0087 auto in1_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
0088 auto in2_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
0089 auto out_d = cms::alpakatools::make_device_buffer<float[]>(queue, size_);
0090
0091
0092
0093
0094 alpaka::memcpy(queue, in1_d, in1_h, size_);
0095 alpaka::memcpy(queue, in2_d, in2_h, size_);
0096
0097
0098 alpaka::memset(queue, out_d, 0);
0099
0100
0101 HeterogeneousTestAlpakaDevicePlugins::wrapper_add_vectors_f(
0102 queue, in1_d.data(), in2_d.data(), out_d.data(), size_);
0103
0104
0105 alpaka::memcpy(queue, out_h, out_d);
0106
0107
0108 alpaka::wait(queue);
0109
0110
0111 for (uint32_t i = 0; i < size_; ++i) {
0112 float sum = in1_h[i] + in2_h[i];
0113 assert(out_h[i] < sum + epsilon);
0114 assert(out_h[i] > sum - epsilon);
0115 }
0116 }
0117
0118 std::cout << "All tests passed.\n";
0119 }
0120
0121 }
0122
0123 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/MakerMacros.h"
0124 DEFINE_FWK_ALPAKA_MODULE(AlpakaTestDeviceAdditionModule);