Skip to content

Commit c955174

Browse files
authored
Merge pull request #65 from acts-project/tests/expand_cuda_strided_tests
Expand CUDA strided array tests
2 parents 8b67233 + 91f9d01 commit c955174

File tree

3 files changed

+172
-92
lines changed

3 files changed

+172
-92
lines changed

tests/cuda/CMakeLists.txt

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -11,7 +11,7 @@ find_package(CUDAToolkit REQUIRED)
1111
include(covfie-compiler-options-cuda)
1212

1313
# Create the test executable from the individual test groups.
14-
add_executable(test_cuda test_cuda_array.cu)
14+
add_executable(test_cuda test_cuda_strided_array.cu)
1515

1616
# Ensure that the tests are linked against the required libraries.
1717
target_link_libraries(

tests/cuda/test_cuda_array.cu

Lines changed: 0 additions & 91 deletions
This file was deleted.
Lines changed: 171 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,171 @@
1+
/*
2+
* SPDX-PackageName: "covfie, a part of the ACTS project"
3+
* SPDX-FileCopyrightText: 2022 CERN
4+
* SPDX-License-Identifier: MPL-2.0
5+
*/
6+
7+
#include <optional>
8+
#include <utility>
9+
10+
#include <gtest/gtest.h>
11+
12+
#include <covfie/core/backend/primitive/array.hpp>
13+
#include <covfie/core/backend/transformer/strided.hpp>
14+
#include <covfie/core/field.hpp>
15+
#include <covfie/core/field_view.hpp>
16+
#include <covfie/core/utility/nd_map.hpp>
17+
#include <covfie/core/vector.hpp>
18+
#include <covfie/cuda/backend/primitive/cuda_device_array.hpp>
19+
20+
#include "retrieve_vector.hpp"
21+
22+
namespace {
23+
class NameGenerator
24+
{
25+
public:
26+
template <typename B>
27+
static std::string GetName(int)
28+
{
29+
using scalar_t = typename B::covariant_output_t::scalar_t;
30+
constexpr std::size_t dims = B::covariant_output_t::dimensions;
31+
32+
if constexpr (std::is_same_v<scalar_t, float>) {
33+
return "Float[" + std::to_string(dims) + "]";
34+
} else if constexpr (std::is_same_v<scalar_t, std::size_t>) {
35+
return "UInt64[" + std::to_string(dims) + "]";
36+
}
37+
}
38+
};
39+
40+
template <std::size_t N, std::size_t... Is>
41+
covfie::array::array<std::size_t, N>
42+
make_size_array_helper(std::index_sequence<Is...> &&)
43+
{
44+
covfie::array::array<std::size_t, N> rv;
45+
((rv[Is] = 10UL), ...);
46+
47+
return rv;
48+
}
49+
50+
template <std::size_t N>
51+
covfie::array::array<std::size_t, N> make_size_array()
52+
{
53+
return make_size_array_helper<N>(std::make_index_sequence<N>());
54+
}
55+
56+
template <std::size_t N, typename T, typename U, std::size_t... Is>
57+
covfie::array::array<U, N>
58+
convert_array_helper(const covfie::array::array<T, N> & a, std::index_sequence<Is...> &&)
59+
{
60+
covfie::array::array<U, N> rv;
61+
((rv[Is] = static_cast<U>(a[Is])), ...);
62+
return rv;
63+
}
64+
65+
template <std::size_t N, typename T, typename U>
66+
covfie::array::array<U, N> convert_array(const covfie::array::array<T, N> & a)
67+
{
68+
return convert_array_helper<N, T, U>(a, std::make_index_sequence<N>());
69+
}
70+
}
71+
72+
template <std::size_t N, typename B>
73+
class TestLookupGeneric : public ::testing::Test
74+
{
75+
protected:
76+
void SetUp() override
77+
{
78+
using scalar_t = typename B::covariant_output_t::scalar_t;
79+
80+
using canonical_backend_t = covfie::backend::strided<
81+
covfie::vector::vector_d<std::size_t, N>,
82+
covfie::backend::array<covfie::vector::vector_d<scalar_t, N>>>;
83+
84+
m_sizes = make_size_array<N>();
85+
86+
covfie::field<canonical_backend_t> f(covfie::make_parameter_pack(
87+
typename canonical_backend_t::configuration_t(m_sizes)
88+
));
89+
covfie::field_view<canonical_backend_t> fv(f);
90+
91+
covfie::utility::nd_map<decltype(m_sizes)>(
92+
[&fv](decltype(m_sizes) t) {
93+
fv.at(t) = convert_array<N, std::size_t, scalar_t>(t);
94+
},
95+
m_sizes
96+
);
97+
98+
m_field = covfie::field<B>(f);
99+
}
100+
101+
covfie::array::array<std::size_t, N> m_sizes;
102+
std::optional<covfie::field<B>> m_field;
103+
};
104+
105+
template <std::size_t N, typename T>
106+
using CudaStridedArrayBackend = covfie::backend::strided<
107+
covfie::vector::vector_d<std::size_t, N>,
108+
covfie::backend::cuda_device_array<covfie::vector::vector_d<T, N>>>;
109+
110+
// 1D tests
111+
template <typename B>
112+
using TestCudaStrided1D = TestLookupGeneric<1, B>;
113+
using Backends1D = ::testing::Types<
114+
CudaStridedArrayBackend<1, std::size_t>,
115+
CudaStridedArrayBackend<1, float>>;
116+
TYPED_TEST_SUITE(TestCudaStrided1D, Backends1D, NameGenerator);
117+
TYPED_TEST(TestCudaStrided1D, LookUp)
118+
{
119+
for (std::size_t x = 0; x < this->m_sizes[0]; ++x) {
120+
std::decay_t<typename covfie::field<TypeParam>::output_t> o =
121+
retrieve_vector<covfie::field<TypeParam>>(*(this->m_field), {x});
122+
123+
EXPECT_EQ(o[0], x);
124+
}
125+
}
126+
127+
// 2D tests
128+
template <typename B>
129+
using TestCudaStrided2D = TestLookupGeneric<2, B>;
130+
using Backends2D = ::testing::Types<
131+
CudaStridedArrayBackend<2, std::size_t>,
132+
CudaStridedArrayBackend<2, float>>;
133+
TYPED_TEST_SUITE(TestCudaStrided2D, Backends2D, NameGenerator);
134+
TYPED_TEST(TestCudaStrided2D, LookUp)
135+
{
136+
for (std::size_t x = 0; x < this->m_sizes[0]; ++x) {
137+
for (std::size_t y = 0; y < this->m_sizes[1]; ++y) {
138+
std::decay_t<typename covfie::field<TypeParam>::output_t> o =
139+
retrieve_vector<covfie::field<TypeParam>>(
140+
*(this->m_field), {x, y}
141+
);
142+
143+
EXPECT_EQ(o[0], x);
144+
EXPECT_EQ(o[1], y);
145+
}
146+
}
147+
}
148+
149+
template <typename B>
150+
using TestCudaStrided3D = TestLookupGeneric<3, B>;
151+
using Backends3D = ::testing::Types<
152+
CudaStridedArrayBackend<3, std::size_t>,
153+
CudaStridedArrayBackend<3, float>>;
154+
TYPED_TEST_SUITE(TestCudaStrided3D, Backends3D, NameGenerator);
155+
TYPED_TEST(TestCudaStrided3D, LookUp)
156+
{
157+
for (std::size_t x = 0; x < this->m_sizes[0]; ++x) {
158+
for (std::size_t y = 0; y < this->m_sizes[1]; ++y) {
159+
for (std::size_t z = 0; z < this->m_sizes[2]; ++z) {
160+
std::decay_t<typename covfie::field<TypeParam>::output_t> o =
161+
retrieve_vector<covfie::field<TypeParam>>(
162+
*(this->m_field), {x, y, z}
163+
);
164+
165+
EXPECT_EQ(o[0], x);
166+
EXPECT_EQ(o[1], y);
167+
EXPECT_EQ(o[2], z);
168+
}
169+
}
170+
}
171+
}

0 commit comments

Comments
 (0)