Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F61234014
test_grid.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, May 5, 09:27
Size
4 KB
Mime Type
text/x-c
Expires
Tue, May 7, 09:27 (1 d, 23 h)
Engine
blob
Format
Raw Data
Handle
17483655
Attached To
rTAMAAS tamaas
test_grid.cpp
View Options
/**
*
* @author Lucas Frérot <lucas.frerot@epfl.ch>
*
* @section LICENSE
*
* Copyright (©) 2016 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 <array>
#include <algorithm>
#include <numeric>
#include <criterion/criterion.h>
#include "grid.hh"
using namespace tamaas;
void tamaas_init() {
initialize();
}
void tamaas_teardown() {
finalize();
}
/* -------------------------------------------------------------------------- */
TestSuite(grid_creation, .init = tamaas_init, .fini = tamaas_teardown);
// Testing 1D creation
Test(grid_creation, creation_1d) {
Grid<Real, 1> grid({5}, 2);
std::array<UInt, 1> correct_size = {5};
cr_assert(grid.sizes() == correct_size, "Wrong sizes");
std::array<UInt, 2> correct_strides = {2, 1};
cr_assert(grid.getStrides() == correct_strides, "Wrong strides");
}
// Testing 2D creation
Test(grid_creation, creation_2d) {
Grid<Real, 2> grid({5, 2}, 3);
std::array<UInt, 2> correct_size = {5, 2};
cr_assert(grid.sizes() == correct_size, "Wrong sizes");
std::array<UInt, 3> correct_strides = {6, 3, 1};
cr_assert(grid.getStrides() == correct_strides, "Wrong strides");
}
// Testing 3D creation
Test(grid_creation, creation_3d) {
Grid<Real, 3> grid({8, 5, 2}, 3);
std::array<UInt, 3> correct_size = {8, 5, 2};
cr_assert(grid.sizes() == correct_size, "Wrong sizes");
std::array<UInt, 4> correct_strides = {30, 15, 3, 1};
cr_assert(grid.getStrides() == correct_strides, "Wrong strides");
}
/* -------------------------------------------------------------------------- */
TestSuite(grid_iterators, .init = tamaas_init, .fini = tamaas_teardown);
// Testing iterators in STL function iota
Test(grid_iterators, 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)
cr_assert(p[i] == i, "Iota fill failed");
}
// Testing filling grid with OpenMP loop
Test(grid_iterators, openmp_loops) {
Grid<UInt, 1> grid({20}, 1);
#pragma omp parallel for
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);
cr_assert(std::mismatch(correct.begin(), correct.end(), grid.begin())
== std::make_pair(correct.end(), grid.end()),
"Fill using OpenMP loop failed");
}
/* -------------------------------------------------------------------------- */
TestSuite(grid_operators, .init = tamaas_init, .fini = tamaas_teardown);
// Testing scalar simple loop-based operators
Test(grid_operators, loop_operators) {
Grid<UInt, 1> grid({20}, 1);
grid = 1;
cr_expect(std::all_of(grid.begin(), grid.end(), [](UInt x){return x == 1;}),
"Assignement operator failed");
grid += 1;
cr_expect(std::all_of(grid.begin(), grid.end(), [](UInt x){return x == 2;}),
"Plus-equal operator failed");
}
// Testing loop-based functions with reductions
Test(grid_operators, 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;
cr_expect(grid.min() == min,
"Minimum function failed");
cr_expect(grid.max() == max,
"Maximum function failed");
cr_expect(grid.sum() == acc,
"Sum function failed");
cr_expect(grid.mean() == mea,
"Mean function failed");
cr_expect(grid.var() == var,
"Var function failed");
}
Event Timeline
Log In to Comment