From d4e6098fab576a5768b3f19222bebf92aa7da3c9 Mon Sep 17 00:00:00 2001 From: Stephen Nicholas Swatman Date: Wed, 20 Nov 2024 14:31:58 +0100 Subject: [PATCH] Replace uses of STL types in device code This enables us to use the code in CUDA without the use of the experimental constexpr relaxation. --- cmake/covfie-compiler-options-cuda.cmake | 3 - examples/core/generate_test_field.cpp | 2 +- examples/core/scaleup_bfield.cpp | 2 +- examples/core/slice3dto2d.cpp | 2 +- examples/cuda/render_slice.cu | 12 +- examples/cuda/render_slice_texture.cu | 10 +- lib/core/covfie/core/algebra/affine.hpp | 5 +- lib/core/covfie/core/algebra/matrix.hpp | 4 +- lib/core/covfie/core/algebra/vector.hpp | 5 +- lib/core/covfie/core/array.hpp | 125 ++++++++++++++++ .../covfie/core/backend/primitive/array.hpp | 1 - .../core/backend/primitive/constant.hpp | 1 - .../core/backend/transformer/hilbert.hpp | 17 +-- .../core/backend/transformer/morton.hpp | 12 +- .../core/backend/transformer/shuffle.hpp | 2 +- .../core/backend/transformer/strided.hpp | 17 +-- lib/core/covfie/core/utility/nd_map.hpp | 59 ++++++-- lib/core/covfie/core/utility/nd_size.hpp | 4 +- lib/core/covfie/core/utility/tuple.hpp | 28 ---- lib/core/covfie/core/vector.hpp | 7 +- .../cuda/backend/primitive/cuda_texture.hpp | 17 +-- tests/core/test_algebra.cpp | 137 +++++++++++------- tests/core/test_atlas_like_io.cpp | 4 +- tests/core/test_canonical.cpp | 6 +- tests/core/test_transformer_affine.cpp | 2 +- tests/core/test_utility.cpp | 16 +- tests/cuda/test_cuda_array.cu | 12 +- 27 files changed, 324 insertions(+), 188 deletions(-) create mode 100644 lib/core/covfie/core/array.hpp delete mode 100644 lib/core/covfie/core/utility/tuple.hpp diff --git a/cmake/covfie-compiler-options-cuda.cmake b/cmake/covfie-compiler-options-cuda.cmake index f24f4b2..43bf553 100644 --- a/cmake/covfie-compiler-options-cuda.cmake +++ b/cmake/covfie-compiler-options-cuda.cmake @@ -31,9 +31,6 @@ if("${CMAKE_CUDA_COMPILER_ID}" MATCHES "NVIDIA") # Make CUDA generate debug symbols for the device code as well in a debug # build. covfie_add_flag( CMAKE_CUDA_FLAGS_DEBUG "-G" ) - # Allow to use functions in device code that are constexpr, even if they are - # not marked with __device__. - covfie_add_flag( CMAKE_CUDA_FLAGS "--expt-relaxed-constexpr" ) endif() covfie_add_flag( CMAKE_CUDA_FLAGS "-Wfloat-conversion" ) diff --git a/examples/core/generate_test_field.cpp b/examples/core/generate_test_field.cpp index 1fc1b1f..f75634f 100644 --- a/examples/core/generate_test_field.cpp +++ b/examples/core/generate_test_field.cpp @@ -74,7 +74,7 @@ int main(int argc, char ** argv) BOOST_LOG_TRIVIAL(info) << "Allocating space for new vector field..."; core_field_t cf(covfie::make_parameter_pack(core_backend_t::configuration_t{ - 201, 201, 301})); + 201u, 201u, 301u})); BOOST_LOG_TRIVIAL(info) << "Filling new vector field..."; diff --git a/examples/core/scaleup_bfield.cpp b/examples/core/scaleup_bfield.cpp index 0e95689..bcbb0ff 100644 --- a/examples/core/scaleup_bfield.cpp +++ b/examples/core/scaleup_bfield.cpp @@ -99,7 +99,7 @@ int main(int argc, char ** argv) output_field_t new_field(covfie::make_parameter_pack( scaling * translation, std::monostate{}, - covfie::utility::nd_size<3>{601, 601, 901} + covfie::utility::nd_size<3>{601u, 601u, 901u} )); BOOST_LOG_TRIVIAL(info) << "Copying data from old field to new field..."; diff --git a/examples/core/slice3dto2d.cpp b/examples/core/slice3dto2d.cpp index a638c85..525588b 100644 --- a/examples/core/slice3dto2d.cpp +++ b/examples/core/slice3dto2d.cpp @@ -122,7 +122,7 @@ int main(int argc, char ** argv) covfie::utility::nd_size<3> in_size = f.backend().get_backend().get_backend().get_configuration(); - covfie::utility::nd_size<2> out_size{0, 0}; + covfie::utility::nd_size<2> out_size{0u, 0u}; if (vm["axis"].as() == "x") { out_size = {in_size[1], in_size[2]}; diff --git a/examples/cuda/render_slice.cu b/examples/cuda/render_slice.cu index 89b593d..0a5e393 100644 --- a/examples/cuda/render_slice.cu +++ b/examples/cuda/render_slice.cu @@ -95,12 +95,12 @@ __global__ void render( typename field_t::output_t p = vf.at(fx * 20000.f - 10000.f, fy * 20000.f - 10000.f, z); - out[height * x + y] = static_cast(std::lround( - 255.f * std::min( - std::sqrt( - std::pow(p[0] / 0.000299792458f, 2.f) + - std::pow(p[1] / 0.000299792458f, 2.f) + - std::pow(p[2] / 0.000299792458f, 2.f) + out[height * x + y] = static_cast(::lroundf( + 255.f * ::fmin( + ::sqrtf( + ::powf(p[0] / 0.000299792458f, 2.f) + + ::powf(p[1] / 0.000299792458f, 2.f) + + ::powf(p[2] / 0.000299792458f, 2.f) ), 1.0f ) diff --git a/examples/cuda/render_slice_texture.cu b/examples/cuda/render_slice_texture.cu index 3f00d12..f728895 100644 --- a/examples/cuda/render_slice_texture.cu +++ b/examples/cuda/render_slice_texture.cu @@ -97,11 +97,11 @@ __global__ void render( typename field_t::output_t p = vf.at(fx * 20000.f - 10000.f, fy * 20000.f - 10000.f, z); out[height * x + y] = static_cast(std::lround( - 255.f * std::min( - std::sqrt( - std::pow(p[0] / 0.000299792458f, 2.f) + - std::pow(p[1] / 0.000299792458f, 2.f) + - std::pow(p[2] / 0.000299792458f, 2.f) + 255.f * fmin( + ::sqrtf( + ::powf(p[0] / 0.000299792458f, 2.f) + + ::powf(p[1] / 0.000299792458f, 2.f) + + ::powf(p[2] / 0.000299792458f, 2.f) ), 1.0f ) diff --git a/lib/core/covfie/core/algebra/affine.hpp b/lib/core/covfie/core/algebra/affine.hpp index 04ea2ad..c8f631d 100644 --- a/lib/core/covfie/core/algebra/affine.hpp +++ b/lib/core/covfie/core/algebra/affine.hpp @@ -14,6 +14,7 @@ #include #include +#include #include namespace covfie::algebra { @@ -90,7 +91,7 @@ struct affine : public matrix { "of the matrix." ); - std::array arr{args...}; + array::array arr{args...}; matrix result = matrix::identity(); @@ -115,7 +116,7 @@ struct affine : public matrix { "the matrix." ); - std::array arr{args...}; + array::array arr{args...}; matrix result = matrix::identity(); diff --git a/lib/core/covfie/core/algebra/matrix.hpp b/lib/core/covfie/core/algebra/matrix.hpp index bf6cbc8..0089c89 100644 --- a/lib/core/covfie/core/algebra/matrix.hpp +++ b/lib/core/covfie/core/algebra/matrix.hpp @@ -10,9 +10,9 @@ #pragma once -#include #include +#include #include namespace covfie::algebra { @@ -26,7 +26,7 @@ struct matrix { { } - COVFIE_DEVICE matrix(std::array, N> l) + COVFIE_DEVICE matrix(array::array, N> l) { for (I i = 0; i < l.size(); ++i) { for (I j = 0; j < l[i].size(); ++j) { diff --git a/lib/core/covfie/core/algebra/vector.hpp b/lib/core/covfie/core/algebra/vector.hpp index e663460..88a66d3 100644 --- a/lib/core/covfie/core/algebra/vector.hpp +++ b/lib/core/covfie/core/algebra/vector.hpp @@ -16,6 +16,7 @@ #include #include +#include #include namespace covfie::algebra { @@ -26,7 +27,7 @@ struct vector : public matrix { { } - COVFIE_DEVICE vector(std::array l) + COVFIE_DEVICE vector(array::array l) : matrix() { for (I i = 0; i < N; ++i) { @@ -47,7 +48,7 @@ struct vector : public matrix { sizeof...(Args) == N, bool> = true> COVFIE_DEVICE vector(Args... args) - : vector(std::array{std::forward(args)...}) + : vector(array::array{std::forward(args)...}) { } diff --git a/lib/core/covfie/core/array.hpp b/lib/core/covfie/core/array.hpp new file mode 100644 index 0000000..2c682bd --- /dev/null +++ b/lib/core/covfie/core/array.hpp @@ -0,0 +1,125 @@ +/* + * This file is part of covfie, a part of the ACTS project + * + * Copyright (c) 2022 CERN + * + * This Source Code Form is subject to the terms of the Mozilla Public License, + * v. 2.0. If a copy of the MPL was not distributed with this file, You can + * obtain one at http://mozilla.org/MPL/2.0/. + */ + +#pragma once + +#include +#include +#include + +#include + +namespace covfie::array { +template +requires(_size > 0) struct array { + using scalar_t = _scalar_t; + using value_type = _scalar_t; + static constexpr std::size_t dimensions = _size; + + COVFIE_DEVICE array() = default; + + COVFIE_DEVICE array(const scalar_t (&arr)[dimensions]) + requires(dimensions > 1) + : array(arr, std::make_index_sequence()) + { + } + + COVFIE_DEVICE array(const scalar_t & val) + : array(val, std::make_index_sequence()) + { + } + + template + requires(sizeof...(Ts) == dimensions) COVFIE_DEVICE array(Ts... args) + : m_data{std::forward(args)...} + { + } + + COVFIE_DEVICE constexpr scalar_t & at(const std::size_t & n) + { + assert(n < dimensions); + + return m_data[n]; + } + + COVFIE_DEVICE constexpr const scalar_t & at(const std::size_t & n) const + { + assert(n < dimensions); + + return m_data[n]; + } + + COVFIE_DEVICE constexpr scalar_t & operator[](const std::size_t & n) + { + assert(n < dimensions); + + return m_data[n]; + } + + COVFIE_DEVICE constexpr const scalar_t & operator[](const std::size_t & n + ) const + { + assert(n < dimensions); + + return m_data[n]; + } + + COVFIE_DEVICE constexpr std::size_t size() const + { + return dimensions; + } + + COVFIE_DEVICE constexpr scalar_t * begin() + { + return m_data + 0; + } + + COVFIE_DEVICE constexpr const scalar_t * begin() const + { + return m_data + 0; + } + + COVFIE_DEVICE constexpr const scalar_t * cbegin() const + { + return m_data + 0; + } + + COVFIE_DEVICE constexpr scalar_t * end() + { + return m_data + dimensions; + } + + COVFIE_DEVICE constexpr const scalar_t * end() const + { + return m_data + dimensions; + } + + COVFIE_DEVICE constexpr const scalar_t * cend() const + { + return m_data + dimensions; + } + +private: + template + COVFIE_DEVICE + array(const scalar_t (&arr)[dimensions], std::index_sequence) + : m_data{arr[Is]...} + { + } + + template + COVFIE_DEVICE array(const scalar_t & val, std::index_sequence) + : m_data{((void)Is, val)...} + { + } + + scalar_t m_data[dimensions]; +}; +} diff --git a/lib/core/covfie/core/backend/primitive/array.hpp b/lib/core/covfie/core/backend/primitive/array.hpp index 7138d41..039d627 100644 --- a/lib/core/covfie/core/backend/primitive/array.hpp +++ b/lib/core/covfie/core/backend/primitive/array.hpp @@ -13,7 +13,6 @@ #include #include #include -#include #include #include diff --git a/lib/core/covfie/core/backend/primitive/constant.hpp b/lib/core/covfie/core/backend/primitive/constant.hpp index 26384e7..2b00359 100644 --- a/lib/core/covfie/core/backend/primitive/constant.hpp +++ b/lib/core/covfie/core/backend/primitive/constant.hpp @@ -10,7 +10,6 @@ #pragma once -#include #include #include diff --git a/lib/core/covfie/core/backend/transformer/hilbert.hpp b/lib/core/covfie/core/backend/transformer/hilbert.hpp index 4931f0e..9d47d97 100644 --- a/lib/core/covfie/core/backend/transformer/hilbert.hpp +++ b/lib/core/covfie/core/backend/transformer/hilbert.hpp @@ -17,6 +17,7 @@ #include #include +#include #include #include #include @@ -25,7 +26,6 @@ #include #include #include -#include #include namespace covfie::backend { @@ -118,20 +118,13 @@ struct hilbert { ); typename T::parent_t::non_owning_data_t nother(other); - using tuple_t = decltype(std::tuple_cat( - std::declval>() - )); - - utility::nd_map( - [¬her, &res](tuple_t t) { - auto old_c = utility::to_array(t); + utility::nd_map( + [¬her, &res](decltype(sizes) t) { coordinate_t c; for (std::size_t i = 0; i < contravariant_input_t::dimensions; ++i) { - c[i] = old_c[i]; + c[i] = t[i]; } std::size_t idx = calculate_index(c); @@ -141,7 +134,7 @@ struct hilbert { res[idx][i] = nother.at(c)[i]; } }, - std::tuple_cat(sizes) + sizes ); return res; diff --git a/lib/core/covfie/core/backend/transformer/morton.hpp b/lib/core/covfie/core/backend/transformer/morton.hpp index a20f0da..c58ad10 100644 --- a/lib/core/covfie/core/backend/transformer/morton.hpp +++ b/lib/core/covfie/core/backend/transformer/morton.hpp @@ -32,7 +32,6 @@ #include #include #include -#include #include namespace covfie::backend { @@ -165,18 +164,15 @@ struct morton { ); typename T::parent_t::non_owning_data_t nother(other); - using tuple_t = decltype(std::tuple_cat(sizes)); - - utility::nd_map( - [¬her, &res](tuple_t t) { - auto old_c = utility::to_array(t); + utility::nd_map( + [¬her, &res](decltype(sizes) t) { typename contravariant_input_t::vector_t c; for (std::size_t i = 0; i < contravariant_input_t::dimensions; ++i) { c[i] = static_cast< typename contravariant_input_t::vector_t::value_type>( - old_c[i] + t[i] ); } @@ -187,7 +183,7 @@ struct morton { res[idx][i] = nother.at(c)[i]; } }, - std::tuple_cat(sizes) + sizes ); return res; diff --git a/lib/core/covfie/core/backend/transformer/shuffle.hpp b/lib/core/covfie/core/backend/transformer/shuffle.hpp index 815ee9f..88ada28 100644 --- a/lib/core/covfie/core/backend/transformer/shuffle.hpp +++ b/lib/core/covfie/core/backend/transformer/shuffle.hpp @@ -107,7 +107,7 @@ struct shuffle { shuffle(typename contravariant_input_t::vector_t c, std::index_sequence) const { - return {std::get(c)...}; + return {c.at(Is)...}; } COVFIE_DEVICE typename covariant_output_t::vector_t diff --git a/lib/core/covfie/core/backend/transformer/strided.hpp b/lib/core/covfie/core/backend/transformer/strided.hpp index f55c46e..82a396c 100644 --- a/lib/core/covfie/core/backend/transformer/strided.hpp +++ b/lib/core/covfie/core/backend/transformer/strided.hpp @@ -22,7 +22,6 @@ #include #include #include -#include #include namespace covfie::backend { @@ -67,19 +66,13 @@ struct strided { ); typename T::parent_t::non_owning_data_t nother(other); - using tuple_t = decltype(std::tuple_cat( - std::declval< - std::array>() - )); - - utility::nd_map( - [&sizes, ¬her, &res](tuple_t t) { - coordinate_t c = utility::to_array(t); + utility::nd_map( + [&sizes, ¬her, &res](decltype(sizes) t) { typename contravariant_input_t::scalar_t idx = 0; for (std::size_t k = 0; k < contravariant_input_t::dimensions; ++k) { - typename contravariant_input_t::scalar_t tmp = c[k]; + typename contravariant_input_t::scalar_t tmp = t[k]; for (std::size_t l = k + 1; l < contravariant_input_t::dimensions; @@ -93,10 +86,10 @@ struct strided { for (std::size_t i = 0; i < covariant_output_t::dimensions; ++i) { - res[idx][i] = nother.at(c)[i]; + res[idx][i] = nother.at(t)[i]; } }, - std::tuple_cat(sizes) + sizes ); return res; diff --git a/lib/core/covfie/core/utility/nd_map.hpp b/lib/core/covfie/core/utility/nd_map.hpp index 48dc690..85d8fad 100644 --- a/lib/core/covfie/core/utility/nd_map.hpp +++ b/lib/core/covfie/core/utility/nd_map.hpp @@ -10,36 +10,71 @@ #pragma once -#include #include -#include +#include + +#include namespace covfie::utility { -template -auto tail_impl(std::index_sequence, [[maybe_unused]] std::tuple t) +template +auto tail_impl( + std::index_sequence, [[maybe_unused]] const array::array & t +) +{ + return array::array{t.at(Ns + 1u)...}; +} + +template +auto tail(const array::array & t) { - return std::make_tuple(std::get(t)...); + return tail_impl(std::make_index_sequence(), t); } -template -auto tail(std::tuple t) +template < + typename T, + std::size_t N1, + std::size_t N2, + std::size_t... Is1, + std::size_t... Is2> +array::array +cat_impl(const array::array & a1, const array::array & a2, std::index_sequence, std::index_sequence) { - return tail_impl(std::make_index_sequence(), t); + return {a1.at(Is1)..., a2.at(Is2)...}; +} + +template +array::array +cat(const array::array & a1, const array::array & a2) +{ + return cat_impl( + a1, a2, std::make_index_sequence(), std::make_index_sequence() + ); } template void nd_map(std::function f, Tuple s) { - if constexpr (std::tuple_size::value == 0u) { + if constexpr (Tuple::dimensions == 0u) { f({}); + } else if constexpr (Tuple::dimensions == 1u) { + for (typename Tuple::value_type i = + static_cast(0); + i < s.at(0); + ++i) + { + f(array::array{i}); + } } else { - using head_t = typename std::tuple_element<0u, Tuple>::type; using tail_t = decltype(tail(std::declval())); - for (head_t i = static_cast(0); i < std::get<0u>(s); ++i) { + for (typename Tuple::value_type i = + static_cast(0); + i < s.at(0); + ++i) + { nd_map( [f, i](tail_t r) { - f(std::tuple_cat(std::tuple{i}, r)); + f(cat(array::array{i}, r)); }, tail(s) ); diff --git a/lib/core/covfie/core/utility/nd_size.hpp b/lib/core/covfie/core/utility/nd_size.hpp index e0c2fc4..65b79de 100644 --- a/lib/core/covfie/core/utility/nd_size.hpp +++ b/lib/core/covfie/core/utility/nd_size.hpp @@ -10,9 +10,9 @@ #pragma once -#include +#include namespace covfie::utility { template -using nd_size = std::array; +using nd_size = array::array; } diff --git a/lib/core/covfie/core/utility/tuple.hpp b/lib/core/covfie/core/utility/tuple.hpp deleted file mode 100644 index 375b6b5..0000000 --- a/lib/core/covfie/core/utility/tuple.hpp +++ /dev/null @@ -1,28 +0,0 @@ -/* - * This file is part of covfie, a part of the ACTS project - * - * Copyright (c) 2022 CERN - * - * This Source Code Form is subject to the terms of the Mozilla Public License, - * v. 2.0. If a copy of the MPL was not distributed with this file, You can - * obtain one at http://mozilla.org/MPL/2.0/. - */ - -#pragma once - -#include -#include - -namespace covfie::utility { -template -std::array to_array(std::tuple i) -{ - return std::apply( - [](T t, Ts... ts) { - return std::array{ - static_cast(t), static_cast(ts)...}; - }, - i - ); -} -} diff --git a/lib/core/covfie/core/vector.hpp b/lib/core/covfie/core/vector.hpp index 721a3ac..cbdbc8a 100644 --- a/lib/core/covfie/core/vector.hpp +++ b/lib/core/covfie/core/vector.hpp @@ -10,9 +10,10 @@ #pragma once -#include #include +#include + namespace covfie::vector { template struct vector_d { @@ -25,7 +26,7 @@ struct array_vector_d { using vector_d = _vector_d; using scalar_t = typename vector_d::type; static constexpr std::size_t dimensions = vector_d::size; - using vector_t = std::array; + using vector_t = array::array; }; template @@ -34,7 +35,7 @@ struct array_reference_vector_d { using scalar_t = typename vector_d::type; static constexpr std::size_t dimensions = vector_d::size; using vector_t = std::add_lvalue_reference_t< - std::array>; + array::array>; }; template diff --git a/lib/cuda/covfie/cuda/backend/primitive/cuda_texture.hpp b/lib/cuda/covfie/cuda/backend/primitive/cuda_texture.hpp index 1909bc3..3d61f34 100644 --- a/lib/cuda/covfie/cuda/backend/primitive/cuda_texture.hpp +++ b/lib/cuda/covfie/cuda/backend/primitive/cuda_texture.hpp @@ -19,7 +19,6 @@ #include #include #include -#include #include #include #include @@ -96,17 +95,11 @@ struct cuda_texture { extent.width * extent.height * extent.depth ); - using tuple_t = decltype(std::tuple_cat( - std::declval< - std::array>( - ) - )); utility::nd_map( - std::function([&no, &stage, &srcSize](tuple_t i) -> void { - auto ia = utility::to_array(i); + std::function([&no, &stage, &srcSize](decltype(srcSize) i) -> void { typename T::parent_t::covariant_output_t::vector_t v = - no.at(ia); + no.at(i); std::size_t idx = 0; @@ -114,7 +107,7 @@ struct cuda_texture { k <= contravariant_input_t::dimensions; --k) { - std::size_t tmp = ia[k]; + std::size_t tmp = i[k]; for (std::size_t l = k - 1; l < k; --l) { tmp *= srcSize[l]; @@ -145,7 +138,7 @@ struct cuda_texture { stage[idx2].w = static_cast(v[3]); } }), - std::tuple_cat(srcSize) + srcSize ); if constexpr (_input_vector_t::size == 2) { diff --git a/tests/core/test_algebra.cpp b/tests/core/test_algebra.cpp index 543ec8b..0013fd9 100644 --- a/tests/core/test_algebra.cpp +++ b/tests/core/test_algebra.cpp @@ -16,7 +16,7 @@ TEST(TestAlgebra, VectorArrayInit1F) { - covfie::algebra::vector<1, float> v(std::array{17.2f}); + covfie::algebra::vector<1, float> v(covfie::array::array{17.2f}); EXPECT_FLOAT_EQ(v(0), 17.2f); } @@ -30,7 +30,7 @@ TEST(TestAlgebra, VectorVariadicInit1F) TEST(TestAlgebra, VectorArrayInit1D) { - covfie::algebra::vector<1, double> v(std::array{21.5}); + covfie::algebra::vector<1, double> v(covfie::array::array{21.5}); EXPECT_DOUBLE_EQ(v(0), 21.5); } @@ -44,7 +44,8 @@ TEST(TestAlgebra, VectorVariadicInit1D) TEST(TestAlgebra, VectorArrayInit2F) { - covfie::algebra::vector<2, float> v(std::array{17.2f, 92.4f}); + covfie::algebra::vector<2, float> v(covfie::array::array{ + 17.2f, 92.4f}); EXPECT_FLOAT_EQ(v(0), 17.2f); EXPECT_FLOAT_EQ(v(1), 92.4f); @@ -60,7 +61,8 @@ TEST(TestAlgebra, VectorVariadicInit2F) TEST(TestAlgebra, VectorArrayInit2D) { - covfie::algebra::vector<2, double> v(std::array{21.5, 11.8}); + covfie::algebra::vector<2, double> v(covfie::array::array{ + 21.5, 11.8}); EXPECT_DOUBLE_EQ(v(0), 21.5); EXPECT_DOUBLE_EQ(v(1), 11.8); @@ -76,7 +78,7 @@ TEST(TestAlgebra, VectorVariadicInit2D) TEST(TestAlgebra, VectorArrayInit3F) { - covfie::algebra::vector<3, float> v(std::array{ + covfie::algebra::vector<3, float> v(covfie::array::array{ 17.2f, 92.4f, 39.5f}); EXPECT_FLOAT_EQ(v(0), 17.2f); @@ -95,8 +97,8 @@ TEST(TestAlgebra, VectorVariadicInit3F) TEST(TestAlgebra, VectorArrayInit3D) { - covfie::algebra::vector<3, double> v(std::array{21.5, 11.8, 28.2} - ); + covfie::algebra::vector<3, double> v(covfie::array::array{ + 21.5, 11.8, 28.2}); EXPECT_DOUBLE_EQ(v(0), 21.5); EXPECT_DOUBLE_EQ(v(1), 11.8); @@ -180,16 +182,19 @@ TEST(TestAlgebra, VectorAssignment3D) TEST(TestAlgebra, MatrixInit1x1F) { - covfie::algebra::matrix<1, 1, float> m(std::array, 1>{ - {{5.2f}}}); + covfie::algebra::matrix<1, 1, float> m( + covfie::array::array, 1>{{{5.2f}}} + ); EXPECT_FLOAT_EQ(m(0, 0), 5.2f); } TEST(TestAlgebra, MatrixInit2x2F) { - covfie::algebra::matrix<2, 2, float> m(std::array, 2>{ - {{0.5f, 10.5f}, {1.5f, 11.5f}}}); + covfie::algebra::matrix<2, 2, float> m( + covfie::array::array, 2>{ + {{0.5f, 10.5f}, {1.5f, 11.5f}}} + ); EXPECT_FLOAT_EQ(m(0, 0), 0.5f); EXPECT_FLOAT_EQ(m(0, 1), 10.5f); @@ -199,8 +204,10 @@ TEST(TestAlgebra, MatrixInit2x2F) TEST(TestAlgebra, MatrixInit3x3F) { - covfie::algebra::matrix<3, 3, float> m(std::array, 3>{ - {{0.5f, 10.5f, 20.5f}, {1.5f, 11.5f, 21.5f}, {2.5f, 12.5f, 22.5f}}}); + covfie::algebra::matrix<3, 3, float> m( + covfie::array::array, 3>{ + {{0.5f, 10.5f, 20.5f}, {1.5f, 11.5f, 21.5f}, {2.5f, 12.5f, 22.5f}}} + ); EXPECT_FLOAT_EQ(m(0, 0), 0.5f); EXPECT_FLOAT_EQ(m(0, 1), 10.5f); @@ -215,8 +222,10 @@ TEST(TestAlgebra, MatrixInit3x3F) TEST(TestAlgebra, MatrixInit3x2F) { - covfie::algebra::matrix<3, 2, float> m(std::array, 3>{ - {{0.5f, 10.5f}, {1.5f, 11.5f}, {2.5f, 12.5f}}}); + covfie::algebra::matrix<3, 2, float> m( + covfie::array::array, 3>{ + {{0.5f, 10.5f}, {1.5f, 11.5f}, {2.5f, 12.5f}}} + ); EXPECT_FLOAT_EQ(m(0, 0), 0.5f); EXPECT_FLOAT_EQ(m(0, 1), 10.5f); @@ -228,8 +237,9 @@ TEST(TestAlgebra, MatrixInit3x2F) TEST(TestAlgebra, AffineInit1F) { - covfie::algebra::affine<1, float> m(std::array, 1>{ - {{0.5f, 10.5f}}}); + covfie::algebra::affine<1, float> m( + covfie::array::array, 1>{{{0.5f, 10.5f}}} + ); EXPECT_FLOAT_EQ(m(0, 0), 0.5f); EXPECT_FLOAT_EQ(m(0, 1), 10.5f); @@ -237,8 +247,10 @@ TEST(TestAlgebra, AffineInit1F) TEST(TestAlgebra, AffineInit2F) { - covfie::algebra::affine<2, float> m(std::array, 2>{ - {{0.5f, 10.5f, 20.5f}, {1.5f, 11.5f, 21.5f}}}); + covfie::algebra::affine<2, float> m( + covfie::array::array, 2>{ + {{0.5f, 10.5f, 20.5f}, {1.5f, 11.5f, 21.5f}}} + ); EXPECT_FLOAT_EQ(m(0, 0), 0.5f); EXPECT_FLOAT_EQ(m(0, 1), 10.5f); @@ -250,10 +262,12 @@ TEST(TestAlgebra, AffineInit2F) TEST(TestAlgebra, AffineInit3F) { - covfie::algebra::affine<3, float> m(std::array, 3>{ - {{0.5f, 10.5f, 20.5f, 30.5f}, - {1.5f, 11.5f, 21.5f, 31.5f}, - {2.5f, 12.5f, 22.5f, 32.5f}}}); + covfie::algebra::affine<3, float> m( + covfie::array::array, 3>{ + {{0.5f, 10.5f, 20.5f, 30.5f}, + {1.5f, 11.5f, 21.5f, 31.5f}, + {2.5f, 12.5f, 22.5f, 32.5f}}} + ); EXPECT_FLOAT_EQ(m(0, 0), 0.5f); EXPECT_FLOAT_EQ(m(0, 1), 10.5f); @@ -271,10 +285,12 @@ TEST(TestAlgebra, AffineInit3F) TEST(TestAlgebra, MatrixMatrixMultiplication1x1x1F) { - covfie::algebra::matrix<1, 1, float> m1(std::array, 1>{ - {{0.5f}}}); - covfie::algebra::matrix<1, 1, float> m2(std::array, 1>{ - {{1.5f}}}); + covfie::algebra::matrix<1, 1, float> m1( + covfie::array::array, 1>{{{0.5f}}} + ); + covfie::algebra::matrix<1, 1, float> m2( + covfie::array::array, 1>{{{1.5f}}} + ); covfie::algebra::matrix<1, 1, float> mr1 = m1 * m2; covfie::algebra::matrix<1, 1, float> mr2 = m2 * m1; @@ -285,10 +301,14 @@ TEST(TestAlgebra, MatrixMatrixMultiplication1x1x1F) TEST(TestAlgebra, MatrixMatrixMultiplication2x2x2F) { - covfie::algebra::matrix<2, 2, float> m1(std::array, 2>{ - {{0.5f, 10.5f}, {1.5f, 11.5f}}}); - covfie::algebra::matrix<2, 2, float> m2(std::array, 2>{ - {{10.5f, 110.5f}, {11.5f, 111.5f}}}); + covfie::algebra::matrix<2, 2, float> m1( + covfie::array::array, 2>{ + {{0.5f, 10.5f}, {1.5f, 11.5f}}} + ); + covfie::algebra::matrix<2, 2, float> m2( + covfie::array::array, 2>{ + {{10.5f, 110.5f}, {11.5f, 111.5f}}} + ); covfie::algebra::matrix<2, 2, float> mr1 = m1 * m2; covfie::algebra::matrix<2, 2, float> mr2 = m2 * m1; @@ -306,12 +326,16 @@ TEST(TestAlgebra, MatrixMatrixMultiplication2x2x2F) TEST(TestAlgebra, MatrixMatrixMultiplication3x3x3F) { - covfie::algebra::matrix<3, 3, float> m1(std::array, 3>{ - {{0.5f, 10.5f, 20.5f}, {1.5f, 11.5f, 21.5f}, {2.5f, 12.5f, 22.5f}}}); - covfie::algebra::matrix<3, 3, float> m2(std::array, 3>{ - {{10.5f, 110.5f, 120.5f}, - {11.5f, 111.5f, 121.5f}, - {12.5f, 112.5f, 122.5f}}}); + covfie::algebra::matrix<3, 3, float> m1( + covfie::array::array, 3>{ + {{0.5f, 10.5f, 20.5f}, {1.5f, 11.5f, 21.5f}, {2.5f, 12.5f, 22.5f}}} + ); + covfie::algebra::matrix<3, 3, float> m2( + covfie::array::array, 3>{ + {{10.5f, 110.5f, 120.5f}, + {11.5f, 111.5f, 121.5f}, + {12.5f, 112.5f, 122.5f}}} + ); covfie::algebra::matrix<3, 3, float> mr1 = m1 * m2; covfie::algebra::matrix<3, 3, float> mr2 = m2 * m1; @@ -339,10 +363,14 @@ TEST(TestAlgebra, MatrixMatrixMultiplication3x3x3F) TEST(TestAlgebra, MatrixMatrixMultiplication3x2x4F) { - covfie::algebra::matrix<3, 2, float> m1(std::array, 3>{ - {{0.5f, 10.5f}, {1.5f, 11.5f}, {2.5f, 12.5f}}}); - covfie::algebra::matrix<2, 4, float> m2(std::array, 2>{ - {{10.5f, 110.5f, 120.5f, 130.5f}, {11.5f, 111.5f, 121.5f, 131.5f}}}); + covfie::algebra::matrix<3, 2, float> m1( + covfie::array::array, 3>{ + {{0.5f, 10.5f}, {1.5f, 11.5f}, {2.5f, 12.5f}}} + ); + covfie::algebra::matrix<2, 4, float> m2( + covfie::array::array, 2>{ + {{10.5f, 110.5f, 120.5f, 130.5f}, {11.5f, 111.5f, 121.5f, 131.5f}}} + ); covfie::algebra::matrix<3, 4, float> mr = m1 * m2; @@ -362,8 +390,10 @@ TEST(TestAlgebra, MatrixMatrixMultiplication3x2x4F) TEST(TestAlgebra, MatrixVectorMultiplication3x3x1F) { - covfie::algebra::matrix<3, 3, float> m(std::array, 3>{ - {{0.5f, 10.5f, 20.5f}, {1.5f, 11.5f, 21.5f}, {2.5f, 12.5f, 22.5f}}}); + covfie::algebra::matrix<3, 3, float> m( + covfie::array::array, 3>{ + {{0.5f, 10.5f, 20.5f}, {1.5f, 11.5f, 21.5f}, {2.5f, 12.5f, 22.5f}}} + ); covfie::algebra::vector<3, float> v(4.2f, 9.1f, 5.5f); covfie::algebra::vector<3, float> r = m * v; @@ -375,8 +405,9 @@ TEST(TestAlgebra, MatrixVectorMultiplication3x3x1F) TEST(TestAlgebra, AffineVectorMultiplication1F) { - covfie::algebra::affine<1, float> m(std::array, 1>{ - {{0.5f, 10.5f}}}); + covfie::algebra::affine<1, float> m( + covfie::array::array, 1>{{{0.5f, 10.5f}}} + ); covfie::algebra::vector<1, float> v(4.2f); covfie::algebra::vector<1, float> r = m * v; @@ -386,8 +417,10 @@ TEST(TestAlgebra, AffineVectorMultiplication1F) TEST(TestAlgebra, AffineVectorMultiplication2F) { - covfie::algebra::affine<2, float> m(std::array, 2>{ - {{0.5f, 10.5f, 20.5f}, {1.5f, 11.5f, 21.5f}}}); + covfie::algebra::affine<2, float> m( + covfie::array::array, 2>{ + {{0.5f, 10.5f, 20.5f}, {1.5f, 11.5f, 21.5f}}} + ); covfie::algebra::vector<2, float> v(4.2f, 7.1f); covfie::algebra::vector<2, float> r = m * v; @@ -398,10 +431,12 @@ TEST(TestAlgebra, AffineVectorMultiplication2F) TEST(TestAlgebra, AffineVectorMultiplication3F) { - covfie::algebra::affine<3, float> m(std::array, 3>{ - {{0.5f, 10.5f, 20.5f, 30.5f}, - {1.5f, 11.5f, 21.5f, 31.5f}, - {2.5f, 12.5f, 22.5f, 32.5f}}}); + covfie::algebra::affine<3, float> m( + covfie::array::array, 3>{ + {{0.5f, 10.5f, 20.5f, 30.5f}, + {1.5f, 11.5f, 21.5f, 31.5f}, + {2.5f, 12.5f, 22.5f, 32.5f}}} + ); covfie::algebra::vector<3, float> v(4.2f, 7.1f, 5.9f); covfie::algebra::vector<3, float> r = m * v; diff --git a/tests/core/test_atlas_like_io.cpp b/tests/core/test_atlas_like_io.cpp index fd8f5d4..b890a59 100644 --- a/tests/core/test_atlas_like_io.cpp +++ b/tests/core/test_atlas_like_io.cpp @@ -35,7 +35,7 @@ TEST(TestAtlasLikeIO, WriteReadAtlasLike) { covfie::field bf( covfie::make_parameter_pack_for>( - {2, 2, 2}, {8} + {2u, 2u, 2u}, {8u} ) ); covfie::field::view_t bv(bf); @@ -111,7 +111,7 @@ TEST(TestAtlasLikeIO, WriteReadAtlasLikeChangeInterpolation) { covfie::field bf( covfie::make_parameter_pack_for>( - {3, 3, 3}, {27} + {3u, 3u, 3u}, {27u} ) ); covfie::field::view_t bv(bf); diff --git a/tests/core/test_canonical.cpp b/tests/core/test_canonical.cpp index 6ef6f08..94dc9a9 100644 --- a/tests/core/test_canonical.cpp +++ b/tests/core/test_canonical.cpp @@ -25,7 +25,7 @@ TEST(TestCanonicalLayout, Canonical1D) covfie::backend::identity>>; field_t f(covfie::make_parameter_pack( - field_t::backend_t::configuration_t{16}, + field_t::backend_t::configuration_t{16u}, field_t::backend_t::backend_t::configuration_t{} )); field_t::view_t fv(f); @@ -42,7 +42,7 @@ TEST(TestCanonicalLayout, Canonical2D) covfie::backend::identity>>; field_t f(covfie::make_parameter_pack( - field_t::backend_t::configuration_t{4, 4}, + field_t::backend_t::configuration_t{4u, 4u}, field_t::backend_t::backend_t::configuration_t{} )); @@ -66,7 +66,7 @@ TEST(TestCanonicalLayout, Canonical3D) covfie::backend::identity>>; field_t f(covfie::make_parameter_pack( - field_t::backend_t::configuration_t{4, 4, 4}, + field_t::backend_t::configuration_t{4u, 4u, 4u}, field_t::backend_t::backend_t::configuration_t{} )); diff --git a/tests/core/test_transformer_affine.cpp b/tests/core/test_transformer_affine.cpp index 07bd069..4206f8d 100644 --- a/tests/core/test_transformer_affine.cpp +++ b/tests/core/test_transformer_affine.cpp @@ -24,7 +24,7 @@ TEST(TestAffineTransformer, AffineConstant1Dto1D) field_t f(covfie::make_parameter_pack( field_t::backend_t::configuration_t(covfie::algebra::affine<1>( - std::array, 1>({{{0.f, 5.f}}}) + covfie::array::array, 1>{{0.f, 5.f}} )), field_t::backend_t::backend_t::configuration_t({5.f}) )); diff --git a/tests/core/test_utility.cpp b/tests/core/test_utility.cpp index 007b516..cd26b46 100644 --- a/tests/core/test_utility.cpp +++ b/tests/core/test_utility.cpp @@ -16,9 +16,11 @@ TEST(TestNDMap, Simple1DAdd) { std::size_t i = 0; - covfie::utility::nd_map>( - std::function([&i](std::tuple t) { i += std::get<0>(t); }), - {10} + covfie::utility::nd_map>( + std::function([&i](covfie::array::array t) { + i += t.at(0); + }), + {10u} ); EXPECT_EQ(i, 45); @@ -28,11 +30,11 @@ TEST(TestNDMap, Simple2DAdd) { std::size_t i = 0; - covfie::utility::nd_map>( - std::function([&i](std::tuple t) { - i += (std::get<0>(t) + std::get<1>(t)); + covfie::utility::nd_map>( + std::function([&i](covfie::array::array t) { + i += (t.at(0) + t.at(1)); }), - {5, 10} + {5u, 10u} ); EXPECT_EQ(i, 325); diff --git a/tests/cuda/test_cuda_array.cu b/tests/cuda/test_cuda_array.cu index 8ed8a70..b601150 100644 --- a/tests/cuda/test_cuda_array.cu +++ b/tests/cuda/test_cuda_array.cu @@ -17,7 +17,6 @@ #include #include #include -#include #include #include @@ -33,20 +32,15 @@ protected: covfie::vector::size3, covfie::backend::array>; - std::array sizes; - sizes.fill(10); + covfie::array::array sizes{10UL, 10UL, 10UL}; covfie::field f(covfie::make_parameter_pack( canonical_backend_t::configuration_t{sizes} )); covfie::field_view fv(f); - covfie::utility::nd_map( - [&fv](decltype(std::tuple_cat(sizes)) t) { - fv.at(covfie::utility::to_array(t)) = - covfie::utility::to_array(t); - }, - std::tuple_cat(sizes) + covfie::utility::nd_map( + [&fv](decltype(sizes) t) { fv.at(t) = t; }, sizes ); m_field = covfie::field(f);