Skip to content

Commit beec217

Browse files
author
Matthew Emmett
committed
tests: Add test_mpi_advection_diffusion.cpp.
This implements a basic test for the MPI-PFASST version of the advection/diffusion example.
1 parent e7a958a commit beec217

File tree

3 files changed

+106
-4
lines changed

3 files changed

+106
-4
lines changed

examples/advection_diffusion/mpi_pfasst.cpp

Lines changed: 11 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -24,10 +24,8 @@ using namespace pfasst;
2424
using namespace pfasst::encap;
2525
using namespace pfasst::mpi;
2626

27-
int main(int argc, char** argv)
27+
error_map run_mpi_pfasst()
2828
{
29-
MPI_Init(&argc, &argv);
30-
3129
const size_t nsteps = 4;
3230
const double dt = 0.01;
3331
const size_t niters = 4;
@@ -63,7 +61,16 @@ int main(int argc, char** argv)
6361
pf.set_duration(0.0, nsteps * dt, dt, niters);
6462
pf.run();
6563

66-
fftw_cleanup();
64+
auto fine = pf.get_finest<AdvectionDiffusionSweeper<>>();
65+
return fine->get_errors();
66+
}
6767

68+
#ifndef PFASST_UNIT_TESTING
69+
int main(int argc, char** argv)
70+
{
71+
MPI_Init(&argc, &argv);
72+
run_mpi_pfasst();
73+
fftw_cleanup();
6874
MPI_Finalize();
6975
}
76+
#endif

tests/examples/advection_diffusion/CMakeLists.txt

Lines changed: 36 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -9,6 +9,12 @@ set(TESTS
99
test_advection_diffusion
1010
)
1111

12+
if(${pfasst_WITH_MPI})
13+
set(MPI_TESTS
14+
test_mpi_advection_diffusion
15+
)
16+
endif()
17+
1218
foreach(test ${TESTS})
1319
message(STATUS " ${test}")
1420
add_executable(${test} ${test}.cpp)
@@ -28,3 +34,33 @@ foreach(test ${TESTS})
2834
COMMAND ${CMAKE_RUNTIME_OUTPUT_DIRECTORY}/tests/examples/advection_diffusion/${test} --gtest_output=xml:${test}_out.xml
2935
)
3036
endforeach(test)
37+
38+
if(${pfasst_WITH_MPI})
39+
include_directories(${MPI_CXX_INCLUDE_PATH})
40+
foreach(test ${MPI_TESTS})
41+
message(STATUS " ${test}")
42+
add_executable(${test} ${test}.cpp)
43+
add_dependencies(${test} googlemock)
44+
if(NOT FFTW_FOUND)
45+
add_dependencies(${test} googlemock fftw3)
46+
endif()
47+
if(MPI_COMPILE_FLAGS)
48+
set_target_properties(${test} PROPERTIES COMPILE_FLAGS "${MPI_COMPILE_FLAGS}")
49+
endif()
50+
if(MPI_LINK_FLAGS)
51+
set_target_properties(${test} PROPERTIES LINK_FLAGS "${MPI_LINK_FLAGS}")
52+
endif()
53+
target_link_libraries(${test}
54+
${3rdparty_DEPENDEND_LIBS}
55+
${FFTW_LIBRARIES}
56+
${pfasst_DEPENDED_LIBS}
57+
${MPI_CXX_LIBRARIES}
58+
)
59+
set_target_properties(${test}
60+
PROPERTIES RUNTIME_OUTPUT_DIRECTORY ${CMAKE_RUNTIME_OUTPUT_DIRECTORY}/tests/examples/advection_diffusion
61+
)
62+
add_test(NAME ${test}
63+
COMMAND ${MPIEXEC} ${MPIEXEC_NUMPROC_FLAG} 4 ${MPIEXEC_PREFLAGS} ${CMAKE_RUNTIME_OUTPUT_DIRECTORY}/tests/examples/advection_diffusion/${test} ${MPIEXEC_POSTFLAGS}
64+
)
65+
endforeach(test)
66+
endif()
Lines changed: 59 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,59 @@
1+
/*
2+
* Tests for the advection-diffusion examples.
3+
*/
4+
5+
#include <algorithm>
6+
#include <iostream>
7+
#include <tuple>
8+
#include <vector>
9+
10+
#include <gtest/gtest.h>
11+
#include <gmock/gmock.h>
12+
13+
#include <mpi.h>
14+
15+
#define PFASST_UNIT_TESTING
16+
#include "../examples/advection_diffusion/mpi_pfasst.cpp"
17+
#undef PFASST_UNIT_TESTING
18+
19+
20+
using namespace std;
21+
22+
// MATCHER(DoubleNear, "")
23+
// {
24+
// return abs(get<0>(arg) - get<1>(arg)) < 1e-15;
25+
// }
26+
27+
// MATCHER(DoubleLess, "")
28+
// {
29+
// return get<0>(arg) < get<1>(arg);
30+
// }
31+
32+
TEST(ErrorTest, MPIPFASST)
33+
{
34+
typedef error_map::value_type vtype;
35+
36+
auto errors = run_mpi_pfasst();
37+
auto get_step = [](const vtype x) { return get<0>(get<0>(x)); };
38+
auto get_iter = [](const vtype x) { return get<1>(get<0>(x)); };
39+
auto get_error = [](const vtype x) { return get<1>(x); };
40+
41+
auto max_iter = get_iter(*std::max_element(errors.begin(), errors.end(),
42+
[get_iter](const vtype p1, const vtype p2) { return get_iter(p1) < get_iter(p2); }));
43+
44+
vector<double> ex = { 1.224103e-10, 5.145808e-10, 3.389905e-9, 1.920198e-7 };
45+
for (auto& x: errors) {
46+
if (get_iter(x) == max_iter) {
47+
EXPECT_NEAR(get_error(x), ex[get_step(x)], 1e-12);
48+
}
49+
}
50+
}
51+
52+
int main(int argc, char** argv)
53+
{
54+
testing::InitGoogleTest(&argc, argv);
55+
MPI_Init(&argc, &argv);
56+
auto result = RUN_ALL_TESTS();
57+
MPI_Finalize();
58+
return result;
59+
}

0 commit comments

Comments
 (0)