Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F71014414
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
Mon, Jul 8, 19:05
Size
6 KB
Mime Type
text/x-c++
Expires
Wed, Jul 10, 19:05 (2 d)
Engine
blob
Format
Raw Data
Handle
18896457
Attached To
rTAMAAS tamaas
test_loop.cpp
View Options
/**
*
* @author Lucas Frérot <lucas.frerot@epfl.ch>
*
* @section LICENSE
*
* Copyright (©) 2017 EPFL (Ecole Polytechnique Fédérale de
* Lausanne) Laboratory (LSMS - Laboratoire de Simulation en Mécanique des
* Solides)
*
* Tamaas is free software: you can redistribute it and/or modify it under the
* terms of the GNU Lesser General Public License as published by the Free
* Software Foundation, either version 3 of the License, or (at your option) any
* later version.
*
* Tamaas 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 Lesser General Public License for more
* details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with Tamaas. If not, see <http://www.gnu.org/licenses/>.
*
*/
/* -------------------------------------------------------------------------- */
#include "test.hh"
#include "grid.hh"
#include "static_types.hh"
#include "loop.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>());
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<UInt>());
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";
}
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(TestStridedLoops, OneArgument) {
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 = AddOneInplace<UInt>();
Loop::stridedLoop(add_one_inplace, grid);
ASSERT_TRUE(compare(solution, grid)) << "One argument strided loop failed";
}
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::stridedLoop(add_one_inplace, 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::stridedLoop(set_one, 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::stridedLoop(BroadcastSet123(), grid);
auto res = Loop::stridedReduce<operation::plus>(VectorReduction(), grid);
ASSERT_EQ(res(0), 100);
ASSERT_EQ(res(1), 200);
ASSERT_EQ(res(2), 300);
}
Event Timeline
Log In to Comment