Page MenuHomec4science

test_mpi.cpp
No OneTemporary

File Metadata

Created
Sat, Jan 18, 06:50

test_mpi.cpp

/**
*
* @author Guillaume Anciaux <guillaume.anciaux@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 "fft_plan_manager.hh"
#include "grid.hh"
#include "grid_hermitian.hh"
#include "loop.hh"
#include "static_types.hh"
#include "tamaas.hh"
#include "test.hh"
#include <algorithm>
#include <fftw3-mpi.h>
/* -------------------------------------------------------------------------- */
using namespace tamaas;
namespace tamaas {
extern int mpi_rank;
}
#define TEST_INFO __PRETTY_FUNCTION__ << " nb_components = " << nb_components
#define TEST_SUCCESS \
MPI_Barrier(MPI_COMM_WORLD); \
std::cout << "success: " << TEST_INFO << std::endl;
#define TEST_START std::cout << "start: " << TEST_INFO << std::endl;
template <typename T, UInt dim>
void TEST1(UInt nb_components, std::initializer_list<UInt> dims) {
TEST_START;
Grid<T, dim> grid(dims, nb_components);
Grid<T, dim> solution(dims, nb_components);
std::iota(grid.begin(), grid.end(), 1);
// Makeing solution
auto add_one = [] CUDA_LAMBDA(T & x) { return x + 1; };
std::transform(grid.begin(), grid.end(), solution.begin(), add_one);
// apply tamaas loop
auto add_one_inplace = [] CUDA_LAMBDA(T & x) { x += 1; };
Loop::loop(add_one_inplace, grid);
assert(compare(grid, solution, AreFloatEqual()));
grid.gather();
solution.gather();
assert(compare(grid, solution, AreFloatEqual()));
TEST_SUCCESS;
}
template <typename T>
void TEST2(UInt nb_components) {
TEST_START;
Grid<T, 2> grid({10, 10}, nb_components);
std::iota(grid.begin(), grid.end(), 1);
Grid<T, 2> solution({10, 10}, nb_components);
solution = grid;
if (nb_components > 1) {
std::for_each(solution.begin(), solution.end(), [&] CUDA_LAMBDA(T & x) {
if (int(x) % nb_components == 1)
x += 1;
});
} else {
std::for_each(solution.begin(), solution.end(),
[&] CUDA_LAMBDA(T & x) { x += 1; });
}
Loop::stridedLoop([] CUDA_LAMBDA(T & x) { x += 1; }, grid);
assert(compare(solution, grid));
grid.gather();
solution.gather();
assert(compare(solution, grid));
TEST_SUCCESS;
}
template <typename T>
void TEST3(UInt nb_components) {
TEST_START;
Grid<T, 2> grid({2, 2}, nb_components);
std::iota(grid.begin(), grid.end(), 1);
const auto id = [] CUDA_LAMBDA(T & x) { return x; };
// Sum reduction
T red = Loop::reduce<operation::plus>(id, grid);
// gather the grid to test in sequential
grid.gather();
if (mpi_rank == 0) {
T sol = std::accumulate(grid.begin(), grid.end(), 0, std::plus<T>());
assert(sol == red);
}
// re-distribute the grid to test in sequential
grid.distribute();
// Product reduction
red = Loop::reduce<operation::times>(id, grid);
// gather the grid to test in sequential
grid.gather();
if (mpi_rank == 0) {
T sol =
std::accumulate(grid.begin(), grid.end(), T(1.), std::multiplies<T>());
assert(sol == red);
}
// re-distribute the grid to test in sequential
grid.distribute();
red = Loop::reduce<operation::min>(id, grid);
grid.gather();
if (mpi_rank == 0) {
// Min reduction
T sol = *std::min_element(grid.begin(), grid.end());
assert(sol == red);
}
// re-distribute the grid to test in sequential
grid.distribute();
red = Loop::reduce<operation::max>(id, grid);
if (mpi_rank == 0) {
// Max reduction
T sol = *std::max_element(grid.begin(), grid.end());
assert(sol == red);
}
TEST_SUCCESS;
}
void TEST4() {
constexpr UInt size = 2;
/* get local data size and allocate */
ptrdiff_t local_n0, local_0_start;
ptrdiff_t alloc_local = fftw_mpi_local_size_2d(size, size / 2 + 1, MPI_COMM_WORLD, &local_n0,
&local_0_start);
TAMAAS_LOG(alloc_local);
std::vector<double> data(local_n0 * size);
std::vector<fftw_complex> solution(local_n0 * (size / 2 + 1));
std::for_each(data.begin(), data.end(), [](double& x) { x = 0.; });
for (auto && x: solution){
x[0] = 0.;
x[1] = 0.;
};
std::iota(std::begin(data), std::end(data), 0);
fftw_plan solution_plan;
if (mpi_size > 1 || true) {
solution_plan = fftw_mpi_plan_dft_r2c_2d(size, size, &data[0], &solution[0],
MPI_COMM_WORLD, FFTW_ESTIMATE);
} else {
solution_plan =
fftw_plan_dft_r2c_2d(size, size, &data[0], &solution[0], FFTW_ESTIMATE);
}
fftw_execute(solution_plan);
Grid<Real, 2> grid({size, size}, 1);
std::iota(grid.begin(), grid.end(), 0);
TAMAAS_LOG(grid);
std::for_each(data.begin(), data.end(), [](double& x) { TAMAAS_LOG(x); });
GridHermitian<Real, 2> result({size, size / 2 + 1}, 1);
FFTPlanManager::get().createPlan(grid, result).forwardTransform();
#ifdef USE_CUDA
cudaDeviceSynchronize();
#endif
TAMAAS_LOG(result);
std::for_each(solution.begin(), solution.end(),
[](fftw_complex& x) { TAMAAS_LOG(x[0] << " , " << x[1]); });
assert(compare(result, solution, AreComplexEqual()));
fftw_destroy_plan(solution_plan);
FFTPlanManager::get().destroyPlan(grid, result);
}
void TEST5() {
constexpr UInt size = 10;
double data[size * size] = {0};
fftw_complex solution[size * (size / 2 + 1)] = {{0}};
std::iota(std::begin(data), std::end(data), 0);
fftw_plan solution_plan =
fftw_plan_dft_r2c_2d(size, size, data, solution, FFTW_ESTIMATE);
fftw_execute(solution_plan);
Grid<Real, 2> grid({size, size}, 1);
std::iota(grid.begin(), grid.end(), 0);
GridHermitian<Real, 2> result({size, size / 2 + 1}, 1);
FFTPlanManager::get().createPlan(grid, result).forwardTransform();
#ifdef USE_CUDA
cudaDeviceSynchronize();
#endif
assert(compare(result, solution, AreComplexEqual()));
fftw_destroy_plan(solution_plan);
FFTPlanManager::get().destroyPlan(grid, result);
}
////////////////////////////////
template <typename T>
void runGridTests() {
for (UInt nb_components = 1; nb_components < 10; ++nb_components) {
TEST1<T, 2>(nb_components, {10, 20});
TEST1<T, 3>(nb_components, {10, 20, 30});
TEST2<T>(nb_components);
TEST3<T>(nb_components);
}
}
void runFFTWTests() {
TEST4();
// TEST5();
}
int main() {
initialize();
runFFTWTests();
runGridTests<Real>();
runGridTests<UInt>();
runGridTests<int>();
finalize();
}

Event Timeline