-
Notifications
You must be signed in to change notification settings - Fork 4.4k
Commit
This commit does not belong to any branch on this repository, and may belong to a fork outside of the repository.
Implement a simple vector add test with alpaka
Implement a simple vector add test with alpaka and cms::alpakatools utilities.
- Loading branch information
Showing
2 changed files
with
100 additions
and
0 deletions.
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
93 changes: 93 additions & 0 deletions
93
HeterogeneousCore/AlpakaInterface/test/alpaka/testKernel.dev.cc
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,93 @@ | ||
#include <random> | ||
#include <cmath> | ||
|
||
#include <alpaka/alpaka.hpp> | ||
|
||
#define CATCH_CONFIG_MAIN | ||
#include <catch.hpp> | ||
|
||
#include "HeterogeneousCore/AlpakaInterface/interface/config.h" | ||
#include "HeterogeneousCore/AlpakaInterface/interface/memory.h" | ||
#include "HeterogeneousCore/AlpakaInterface/interface/vec.h" | ||
#include "HeterogeneousCore/AlpakaInterface/interface/workdivision.h" | ||
|
||
// each test binary is built for a single Alpaka backend | ||
using namespace ALPAKA_ACCELERATOR_NAMESPACE; | ||
|
||
static constexpr auto s_tag = "[" ALPAKA_TYPE_ALIAS_NAME(alpakaTestKernel) "]"; | ||
|
||
struct VectorAddKernel { | ||
template <typename TAcc, typename T> | ||
ALPAKA_FN_ACC void operator()( | ||
TAcc const& acc, T const* __restrict__ in1, T const* __restrict__ in2, T* __restrict__ out, size_t size) const { | ||
for (auto index : cms::alpakatools::elements_with_stride(acc, size)) { | ||
out[index] = in1[index] + in2[index]; | ||
} | ||
} | ||
}; | ||
|
||
TEST_CASE("Standard checks of " ALPAKA_TYPE_ALIAS_NAME(alpakaTestKernel), s_tag) { | ||
SECTION("Kernel1D") { | ||
// get the list of devices on the current platform | ||
auto const& devices = cms::alpakatools::devices<Platform>(); | ||
if (devices.empty()) { | ||
std::cout << "No devices on platform\n"; | ||
return; | ||
} | ||
|
||
// random number generator with a gaussian distribution | ||
std::random_device rd{}; | ||
std::mt19937 gen{rd()}; | ||
std::normal_distribution<> dist{0., 1.}; | ||
|
||
// tolerance | ||
constexpr float epsilon = 0.0001; | ||
|
||
// buffer size | ||
constexpr size_t size = 1024 * 1024; | ||
|
||
// create a default queue on the first device | ||
auto const& device = devices[0]; | ||
auto queue = Queue(device); | ||
|
||
// allocate buffers for the input data on the host and fill them with random data | ||
auto in1_h = cms::alpakatools::make_host_buffer<float[]>(queue, size); | ||
auto in2_h = cms::alpakatools::make_host_buffer<float[]>(queue, size); | ||
alpaka::wait(queue); | ||
for (size_t i = 0; i < size; ++i) { | ||
in1_h[i] = dist(gen); | ||
in2_h[i] = dist(gen); | ||
} | ||
|
||
// allocate buffers for the input data on the device | ||
auto in1_d = cms::alpakatools::make_device_buffer<float[]>(queue, size); | ||
auto in2_d = cms::alpakatools::make_device_buffer<float[]>(queue, size); | ||
|
||
// copy the input data to the device; the size is known from the buffer objects | ||
alpaka::memcpy(queue, in1_d, in1_h); | ||
alpaka::memcpy(queue, in2_d, in2_h); | ||
|
||
// allocate a buffer for the results on the device, and fill it with zeros | ||
auto out_d = cms::alpakatools::make_device_buffer<float[]>(queue, size); | ||
alpaka::memset(queue, out_d, 0.); | ||
|
||
// launch the 1-dimensional kernel | ||
auto div = cms::alpakatools::make_workdiv<Acc1D>(16, 1024); | ||
alpaka::exec<Acc1D>(queue, div, VectorAddKernel{}, in1_d.data(), in2_d.data(), out_d.data(), size); | ||
|
||
// copy the results from the device to the host | ||
auto out_h = cms::alpakatools::make_host_buffer<float[]>(queue, size); | ||
alpaka::memcpy(queue, out_h, out_d); | ||
|
||
// wait for all the operations to complete | ||
alpaka::wait(queue); | ||
|
||
// check the results | ||
for (size_t i = 0; i < size; ++i) { | ||
//REQUIRE(std::abs(out_h[i] - (in1_h[i] + in2_h[i])) < epsilon); | ||
float sum = in1_h[i] + in2_h[i]; | ||
REQUIRE(out_h[i] < sum + epsilon); | ||
REQUIRE(out_h[i] > sum - epsilon); | ||
} | ||
} | ||
} |