Page MenuHomec4science

test_grid.cpp
No OneTemporary

File Metadata

Created
Tue, May 28, 19:04

test_grid.cpp

/**
*
* @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 "grid_view.hh"
using namespace tamaas;
/* -------------------------------------------------------------------------- */
// Testing 1D creation
TEST(TestGridCreation, Creation1d) {
Grid<Real, 1> grid{{5}, 2};
std::array<UInt, 1> correct_size{{5}};
ASSERT_TRUE(grid.sizes() == correct_size) << "Wrong sizes";
std::array<UInt, 2> correct_strides{{2, 1}};
ASSERT_TRUE(grid.getStrides() == correct_strides) << "Wrong strides";
}
// Testing 2D creation
TEST(TestGridCreation, Creation2d) {
Grid<Real, 2> grid({5, 2}, 3);
std::array<UInt, 2> correct_size{{5, 2}};
ASSERT_TRUE(grid.sizes() == correct_size) << "Wrong sizes";
std::array<UInt, 3> correct_strides{{6, 3, 1}};
ASSERT_TRUE(grid.getStrides() == correct_strides) << "Wrong strides";
}
// Testing 3D creation
TEST(TestGridCreation, Creation3d) {
Grid<Real, 3> grid({8, 5, 2}, 3);
std::array<UInt, 3> correct_size{{8, 5, 2}};
ASSERT_TRUE(grid.sizes() == correct_size) << "Wrong sizes";
std::array<UInt, 4> correct_strides{{30, 15, 3, 1}};
ASSERT_TRUE(grid.getStrides() == correct_strides) << "Wrong strides";
}
/* -------------------------------------------------------------------------- */
// Testing iterators in STL function iota
TEST(TestGridIterators, Iota) {
constexpr UInt N = 20;
Grid<UInt, 1> grid({N}, 1);
std::iota(grid.begin(), grid.end(), 0);
const UInt * p = grid.getInternalData();
for (UInt i = 0 ; i < N ; ++i)
ASSERT_TRUE(p[i] == i) << "Iota fill failed";
}
// Testing filling grid with OpenMP loop
TEST(TestGridIterators, OpenMPLoops) {
Grid<UInt, 1> grid({20}, 1);
#ifndef USE_CUDA
#pragma omp parallel for
#endif
for (auto it = grid.begin() ; it < grid.end() ; ++it) {
UInt i = it - grid.begin();
*it = i;
}
std::vector<UInt> correct(20);
std::iota(correct.begin(), correct.end(), 0);
ASSERT_TRUE(compare(grid, correct)) << "Fill using OpenMP loop failed";
}
/* -------------------------------------------------------------------------- */
// Testing scalar simple loop-based operators
TEST(TestGridOperators, LoopOperators) {
Grid<UInt, 1> grid({20}, 1);
grid = 1;
EXPECT_TRUE(std::all_of(grid.begin(), grid.end(), [](UInt x) {
return x == 1;
})) << "Assignement operator failed";
grid += 1;
EXPECT_TRUE(std::all_of(grid.begin(), grid.end(), [](UInt x) {
return x == 2;
})) << "Plus-equal operator failed";
}
// Testing loop-based functions with reductions
TEST(TestGridOperators, Stats) {
constexpr UInt N = 20;
Grid<Real, 1> grid({N}, 1);
std::iota(grid.begin(), grid.end(), 0);
auto b = grid.begin(), e = grid.end();
Real min = *std::min_element(b, e);
Real max = *std::max_element(b, e);
Real acc = std::accumulate(b, e, 0);
Real mea = acc / N;
Real var = 35;
EXPECT_DOUBLE_EQ(grid.min(), min) << "Minimum function failed";
EXPECT_DOUBLE_EQ(grid.max(), max) << "Maximum function failed";
EXPECT_DOUBLE_EQ(grid.sum(), acc) <<"Sum function failed";
EXPECT_DOUBLE_EQ(grid.mean(), mea) << "Mean function failed";
EXPECT_DOUBLE_EQ(grid.var(), var) << "Var function failed";
}
/* -------------------------------------------------------------------------- */
TEST(TestGridView, View2D) {
Grid<UInt, 3> grid3d({10, 10, 10}, 3);
std::iota(grid3d.begin(), grid3d.end(), 0);
for (UInt i : Loop::range(10)) {
auto grid_view = make_view(grid3d, i);
std::vector<UInt> reference(300);
std::iota(reference.begin(), reference.end(), 300 * i);
EXPECT_TRUE(
std::equal(grid_view.begin(), grid_view.end(), reference.begin()));
}
}
/* -------------------------------------------------------------------------- */
TEST(TestGridOperators, BroadcastOperators) {
Grid<UInt, 2> grid({10, 10}, 3), solution({10, 10}, 3);
// Filling up solution
Loop::stridedLoop([](VectorProxy<UInt, 3>&& v) {
v(0) = 1;
v(1) = 2;
v(2) = 3;
}, solution);
Vector<UInt, 3> v;
v(0) = 1;
v(1) = 2;
v(2) = 3;
grid = v;
EXPECT_TRUE(std::equal(grid.begin(), grid.end(), solution.begin()))
<< "Broadcast assignement failure";
v(1) = 1;
v(2) = 1;
solution += 1;
grid += v;
EXPECT_TRUE(std::equal(grid.begin(), grid.end(), solution.begin()))
<< "Broadcast += failure";
solution -= 1;
grid -= v;
EXPECT_TRUE(std::equal(grid.begin(), grid.end(), solution.begin()))
<< "Broadcast -= failure";
v(1) = 2;
v(2) = 3;
solution *= solution;
grid *= v;
EXPECT_TRUE(std::equal(grid.begin(), grid.end(), solution.begin()))
<< "Broadcast *= failure";
Loop::stridedLoop([](VectorProxy<UInt, 3>&& v) {
v(0) = 1;
v(1) = 2;
v(2) = 3;
}, solution);
grid /= v;
EXPECT_TRUE(std::equal(grid.begin(), grid.end(), solution.begin()))
<< "Broadcast /= failure";
}

Event Timeline