Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F84443533
test_loop.cpp
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Sun, Sep 22, 22:39
Size
10 KB
Mime Type
text/x-c++
Expires
Tue, Sep 24, 22:39 (2 d)
Engine
blob
Format
Raw Data
Handle
21020909
Attached To
rTAMAAS tamaas
test_loop.cpp
View Options
/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see <https://www.gnu.org/licenses/>.
*
*/
/* -------------------------------------------------------------------------- */
#include "grid.hh"
#include "grid_view.hh"
#include "loop.hh"
#include "static_types.hh"
#include "test.hh"
/* -------------------------------------------------------------------------- */
/* WARNING: here we cannot use lambdas for tests because GoogleTest declares */
/* test functions as private members of classes which is incompatible with */
/* cuda's extended lambdas. I know... it's f*cking stupid */
/* -------------------------------------------------------------------------- */
using namespace tamaas;
template <typename T>
struct AddOneInplace {
CUDA_LAMBDA void operator()(T& x) { x += 1; }
};
// Testing loops on one grid
TEST(TestLoops, OneArgument) {
Grid<Real, 1> grid({20}, 1);
Grid<Real, 1> solution({20}, 1);
auto add_one = [](Real& x) { return x + 1; };
std::iota(grid.begin(), grid.end(), 1);
// Makeing solution
std::transform(grid.begin(), grid.end(), solution.begin(), add_one);
auto add_one_inplace = AddOneInplace<Real>();
Loop::loop(add_one_inplace, grid);
ASSERT_TRUE(compare(grid, solution, AreFloatEqual()))
<< "One argument loop failed";
}
struct PrimalTest {
CUDA_LAMBDA void operator()(Int& primal, Int& val) {
val = (primal > 0) ? -1 : 1;
}
};
// Testing loops on two grids
TEST(TestLoops, TwoArguments) {
// Why no ints?
Grid<Int, 2> grid({20, 20}, 1);
Grid<Int, 2> primal({20, 20}, 1);
Grid<Int, 2> solution({20, 20}, 1);
primal(0, 0) = 1;
primal(0, 1) = 1;
primal(1, 0) = 1;
primal(1, 1) = 1;
std::transform(primal.begin(), primal.end(), solution.begin(),
[](Int& primal) { return (primal > 0) ? -1 : 1; });
auto primal_test = PrimalTest();
Loop::loop(primal_test, primal, grid);
ASSERT_TRUE(compare(solution, grid)) << "Two argument loop failed";
}
struct AssignUInt {
CUDA_LAMBDA void operator()(UInt& x, UInt i) { x = i; }
};
// Testing an enumeration
TEST(TestLoops, Enumeration) {
Grid<UInt, 1> grid({100}, 1);
Grid<UInt, 1> solution({100}, 1);
std::iota(solution.begin(), solution.end(), 0);
auto assign_uint = AssignUInt();
Loop::loop(assign_uint, grid, Loop::range(100));
ASSERT_TRUE(compare(solution, grid)) << "Enumeration loop failed";
}
/* -------------------------------------------------------------------------- */
struct Identity {
CUDA_LAMBDA UInt operator()(UInt& x) const { return x; }
};
// Testing one grid reductions
TEST(TestReductions, OneArgument) {
Grid<UInt, 1> grid({6}, 1);
std::iota(grid.begin(), grid.end(), 1);
const auto id = Identity();
// Sum reduction
UInt sol = std::accumulate(grid.begin(), grid.end(), 0, std::plus<>());
UInt red = Loop::reduce<operation::plus>(id, grid);
ASSERT_TRUE(sol == red) << "Addition reduction failed on one argument";
// Product reduction
sol = std::accumulate(grid.begin(), grid.end(), 1, std::multiplies<>());
red = Loop::reduce<operation::times>(id, grid);
ASSERT_TRUE(sol == red) << "Multiplication reduction failed on one argument";
// Min reduction
sol = *std::min_element(grid.begin(), grid.end());
red = Loop::reduce<operation::min>(id, grid);
ASSERT_TRUE(sol == red) << "Min reduction failed on one argument";
// Max reduction
sol = *std::max_element(grid.begin(), grid.end());
red = Loop::reduce<operation::max>(id, grid);
ASSERT_TRUE(sol == red) << "Max reduction failed on one argument";
}
TEST(TestReductions, ReduceAndTransform) {
UInt n = 20;
Grid<UInt, 1> grid({n}, 1), solution({n}, 1);
std::iota(solution.begin(), solution.end(), 0);
UInt sum_value = (n - 1) * n / 2;
UInt res = Loop::reduce<operation::plus>(
[](UInt& x, UInt i) {
x = i;
return x;
},
grid, Loop::range(n));
EXPECT_EQ(res, sum_value) << "Reduction failed";
EXPECT_TRUE(compare(grid, solution)) << "Assign failed";
}
struct PrimalReduce {
CUDA_LAMBDA UInt operator()(UInt& p, UInt& val) { return (p > 0) ? val : 0; }
};
TEST(TestReductions, TwoArguments) {
Grid<UInt, 1> grid({20}, 1);
Grid<UInt, 1> primal({20}, 1);
grid = 1;
primal(0) = 1;
primal(1) = 1;
auto primal_reduce = PrimalReduce();
// Reduce on values where primal > 0
UInt red = Loop::reduce<operation::plus>(primal_reduce, primal, grid);
ASSERT_TRUE(red == 2) << "Two args reduction failed";
}
/* -------------------------------------------------------------------------- */
TEST(TestRange, type_trait) {
Grid<UInt, 1> grid({1}, 1);
auto gridrange = range<VectorProxy<UInt, 1>>(grid);
static_assert(decltype(gridrange)::is_valid_container<Grid<UInt, 1>>::value,
"is_valid_container Type trait is wrong");
static_assert(
not decltype(gridrange)::is_valid_container<Grid<Real, 1>&>::value,
"is_valid_container Type trait is wrong");
static_assert(not Range<VectorProxy<Real, 1>, Real,
1>::is_valid_container<decltype(grid)>::value,
"is_valid_container Type trait is wrong");
}
TEST(TestRange, headless) {
Grid<UInt, 1> grid({10}, 1), solution({10}, 1);
std::fill(++solution.begin(), solution.end(), 1);
auto gridrange = range<VectorProxy<UInt, 1>>(grid).headless();
Loop::loop([](auto x) { x = 1; }, gridrange);
ASSERT_TRUE(compare(grid, solution)) << "Headless fail";
}
template <typename T>
using WrapVector = VectorProxy<T, 2>;
struct AddOneVector {
CUDA_LAMBDA void operator()(WrapVector<UInt> x) { x(0) += 1; }
};
TEST(TestStridedLoops, VectorStride) {
Grid<UInt, 2> grid({10, 10}, 2);
std::iota(grid.begin(), grid.end(), 1);
Grid<UInt, 2> solution({10, 10}, 2);
solution = grid;
std::for_each(solution.begin(), solution.end(), [](UInt& x) {
if (x % 2 == 1)
x += 1;
});
auto add_one_inplace = AddOneVector();
Loop::loop(add_one_inplace, range<WrapVector<UInt>>(grid));
ASSERT_TRUE(compare(solution, grid)) << "Static vector strided loop failed";
}
template <typename T>
using WrapMatrix = MatrixProxy<T, 2, 2>;
struct SetOneMatrix {
CUDA_LAMBDA void operator()(WrapMatrix<UInt> x) {
x(0, 0) = 1;
x(1, 1) = 1;
}
};
TEST(TestStridedLoops, MatrixStride) {
Grid<UInt, 2> grid({10, 10}, 4);
Grid<UInt, 2> solution({10, 10}, 4);
std::iota(solution.begin(), solution.end(), 0);
std::for_each(solution.begin(), solution.end(), [](UInt& x) {
if (x % 4 == 0 || x % 4 == 3)
x = 1;
else
x = 0;
});
auto set_one = SetOneMatrix();
Loop::loop(set_one, range<WrapMatrix<UInt>>(grid));
ASSERT_TRUE(compare(solution, grid)) << "Static matrix strided loop failed";
}
struct VectorReduction {
CUDA_LAMBDA Vector<UInt, 3> operator()(VectorProxy<UInt, 3> v) const {
return v;
}
};
struct BroadcastSet123 {
CUDA_LAMBDA inline void operator()(VectorProxy<UInt, 3> v) const {
v(0) = 1;
v(1) = 2;
v(2) = 3;
}
};
TEST(TestStridedReduction, VectorReduce) {
Grid<UInt, 2> grid({10, 10}, 3);
Loop::loop(BroadcastSet123(), range<VectorProxy<UInt, 3>>(grid));
auto res = Loop::reduce<operation::plus>(VectorReduction(),
range<VectorProxy<UInt, 3>>(grid));
ASSERT_EQ(res(0), 100);
ASSERT_EQ(res(1), 200);
ASSERT_EQ(res(2), 300);
}
TEST(TestViewReduction, ScalarReduce) {
Grid<UInt, 2> grid({10, 10}, 3);
Loop::loop(BroadcastSet123(), range<VectorProxy<UInt, 3>>(grid));
auto view = make_component_view(grid, 2);
UInt res = Loop::reduce<operation::plus>([](UInt& x) { return x; }, view);
EXPECT_EQ(res, 300) << "Reduce on component view fail";
}
TEST(TestViewReduction, VectorReduce) {
Grid<UInt, 2> grid({10, 10}, 3);
auto view2 = make_view(grid, 0);
Loop::loop(BroadcastSet123(), range<VectorProxy<UInt, 3>>(view2));
auto res2 = Loop::reduce<operation::plus>(VectorReduction(),
range<VectorProxy<UInt, 3>>(view2));
EXPECT_EQ(res2(0), 10);
EXPECT_EQ(res2(1), 20);
EXPECT_EQ(res2(2), 30);
}
TEST(TestViewLoop, ScalarLoop) {
Grid<UInt, 2> grid({10, 10}, 3), solution({10, 10}, 3);
auto view = make_component_view(grid, 2);
Loop::loop([](auto& x) { x = 1; }, view);
Loop::loop([](auto v) { v(2) = 1; }, range<VectorProxy<UInt, 3>>(solution));
ASSERT_TRUE(compare(grid, solution)) << "View loop fail";
}
TEST(TestLoopChecks, Components) {
Grid<UInt, 2> grid({10, 10}, 3);
EXPECT_THROW(
Loop::loop([](auto v) { v = 0; }, range<VectorProxy<UInt, 2>>(grid)),
Exception)
<< "Broken check on number of components";
}
TEST(TestLoopChecks, LoopSize) {
Grid<UInt, 1> grid({10}, 2), other({10}, 1);
EXPECT_THROW(Loop::loop([](auto& x, auto& y) { x = y; }, grid, other),
Exception)
<< "Check on loop size without ranges fail";
other.resize({11});
EXPECT_THROW(Loop::loop([](auto x, auto y) { x(0) = y(0); },
range<VectorProxy<UInt, 2>>(grid),
range<VectorProxy<UInt, 1>>(other)),
Exception)
<< "Check on loop size with ranges fail";
Grid<UInt, 2> twod({10, 11}, 2);
auto view = make_view(twod, 0);
EXPECT_THROW(Loop::loop([](auto& x, auto& y) { x = y; }, grid, view),
Exception)
<< "Check on loop size with view fail";
}
TEST(TestReductions, ReduceAndTransformVector) {
UInt n = 20;
Grid<UInt, 1> grid({n}, 2), solution({n}, 2);
std::iota(solution.begin(), solution.end(), 1);
std::iota(grid.begin(), grid.end(), 0);
UInt sum_value = (2 * n + 1) * 2 * n / 2;
UInt res = Loop::reduce<operation::plus>(
[](auto x, UInt) {
x += 1;
return x(0) + x(1);
},
range<VectorProxy<UInt, 2>>(grid), Loop::range(n));
EXPECT_EQ(res, sum_value) << "Reduction failed";
EXPECT_TRUE(compare(grid, solution)) << "Assign failed";
}
Event Timeline
Log In to Comment