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 "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/Frameworkfwd.h"
0008 #include "FWCore/Framework/interface/global/EDAnalyzer.h"
0009 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0012 #include "FWCore/ServiceRegistry/interface/Service.h"
0013 #include "HeterogeneousCore/AlpakaInterface/interface/config.h"
0014 #include "HeterogeneousTest/AlpakaOpaque/interface/alpaka/DeviceAdditionOpaque.h"
0015 #include "HeterogeneousCore/AlpakaServices/interface/alpaka/AlpakaService.h"
0016
0017 namespace ALPAKA_ACCELERATOR_NAMESPACE {
0018
0019 class AlpakaTestOpaqueAdditionModule : public edm::global::EDAnalyzer<> {
0020 public:
0021 explicit AlpakaTestOpaqueAdditionModule(edm::ParameterSet const& config);
0022 ~AlpakaTestOpaqueAdditionModule() override = default;
0023
0024 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0025
0026 void analyze(edm::StreamID, edm::Event const& event, edm::EventSetup const& setup) const override;
0027
0028 private:
0029 const uint32_t size_;
0030 };
0031
0032 AlpakaTestOpaqueAdditionModule::AlpakaTestOpaqueAdditionModule(edm::ParameterSet const& config)
0033 : size_(config.getParameter<uint32_t>("size")) {}
0034
0035 void AlpakaTestOpaqueAdditionModule::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0036 edm::ParameterSetDescription desc;
0037 desc.add<uint32_t>("size", 1024 * 1024);
0038
0039
0040 edm::ParameterSetDescription alpaka;
0041 alpaka.setAllowAnything();
0042 desc.addUntracked<edm::ParameterSetDescription>("alpaka", alpaka);
0043
0044 descriptions.addWithDefaultLabel(desc);
0045 }
0046
0047 void AlpakaTestOpaqueAdditionModule::analyze(edm::StreamID,
0048 edm::Event const& event,
0049 edm::EventSetup const& setup) const {
0050
0051 edm::Service<ALPAKA_TYPE_ALIAS(AlpakaService)> service;
0052 if (not service or not service->enabled()) {
0053 std::cout << "The " << ALPAKA_TYPE_ALIAS_NAME(AlpakaService)
0054 << " is not available or disabled, the test will be skipped.\n";
0055 return;
0056 }
0057
0058
0059 std::random_device rd{};
0060 std::default_random_engine rand{rd()};
0061 std::normal_distribution<float> dist{0., 1.};
0062
0063
0064 constexpr float epsilon = 0.000001;
0065
0066
0067 std::vector<float> in1(size_);
0068 std::vector<float> in2(size_);
0069 std::vector<float> out(size_);
0070
0071
0072 for (uint32_t i = 0; i < size_; ++i) {
0073 in1[i] = dist(rand);
0074 in2[i] = dist(rand);
0075 out[i] = 0.;
0076 }
0077
0078
0079 test::opaque_add_vectors_f(in1.data(), in2.data(), out.data(), size_);
0080
0081
0082 for (uint32_t i = 0; i < size_; ++i) {
0083 float sum = in1[i] + in2[i];
0084 assert(out[i] < sum + epsilon);
0085 assert(out[i] > sum - epsilon);
0086 }
0087
0088 std::cout << "All tests passed.\n";
0089 }
0090
0091 }
0092
0093 #include "HeterogeneousCore/AlpakaCore/interface/alpaka/MakerMacros.h"
0094 DEFINE_FWK_ALPAKA_MODULE(AlpakaTestOpaqueAdditionModule);