Skip to content

Commit ab01747

Browse files
committed
tests: Add AdaptiveErrorTest for AD example to MPI PFASST.
Signed-off-by: Matthew Emmett <[email protected]>
1 parent 3ebb641 commit ab01747

File tree

2 files changed

+29
-6
lines changed

2 files changed

+29
-6
lines changed

examples/advection_diffusion/mpi_pfasst.cpp

Lines changed: 3 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -31,12 +31,10 @@ namespace pfasst
3131
{
3232
namespace advection_diffusion
3333
{
34-
error_map run_mpi_pfasst()
34+
error_map run_mpi_pfasst(double abs_residual_tol, size_t niters=4)
3535
{
36-
// const size_t nsteps = 8;
3736
const size_t nsteps = 4;
3837
const double dt = 0.01;
39-
const size_t niters = 4;
4038

4139
vector<pair<size_t, quadrature::QuadratureType>> nodes = {
4240
{ 3, quadrature::QuadratureType::GaussLobatto },
@@ -68,7 +66,7 @@ namespace pfasst
6866
pf.set_comm(&comm);
6967
pf.set_duration(0.0, nsteps * dt, dt, niters);
7068
pf.set_nsweeps({2, 1});
71-
// pf.get_finest<AdvectionDiffusionSweeper<>>()->set_residual_tolerances(1e-8, 0.0);
69+
pf.get_finest<AdvectionDiffusionSweeper<>>()->set_residual_tolerances(abs_residual_tol, 0.0);
7270
pf.run();
7371

7472
auto fine = pf.get_finest<AdvectionDiffusionSweeper<>>();
@@ -83,7 +81,7 @@ namespace pfasst
8381
int main(int argc, char** argv)
8482
{
8583
MPI_Init(&argc, &argv);
86-
pfasst::examples::advection_diffusion::run_mpi_pfasst();
84+
pfasst::examples::advection_diffusion::run_mpi_pfasst(0.0);
8785
fftw_cleanup();
8886
MPI_Finalize();
8987
}

tests/examples/advection_diffusion/test_mpi_advection_diffusion.cpp

Lines changed: 26 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -24,7 +24,7 @@ TEST(ErrorTest, MPIPFASST)
2424
{
2525
typedef error_map::value_type vtype;
2626

27-
auto errors = run_mpi_pfasst();
27+
auto errors = run_mpi_pfasst(0.0);
2828
auto get_step = [](const vtype x) { return get<0>(get<0>(x)); };
2929
auto get_iter = [](const vtype x) { return get<1>(get<0>(x)); };
3030
auto get_error = [](const vtype x) { return get<1>(x); };
@@ -40,6 +40,31 @@ TEST(ErrorTest, MPIPFASST)
4040
}
4141
}
4242

43+
TEST(AdaptiveErrorTest, MPIPFASST)
44+
{
45+
typedef error_map::value_type vtype;
46+
47+
auto errors = run_mpi_pfasst(1.e-8, 12);
48+
auto get_step = [](const vtype x) { return get<0>(get<0>(x)); };
49+
auto get_iter = [](const vtype x) { return get<1>(get<0>(x)); };
50+
auto get_error = [](const vtype x) { return get<1>(x); };
51+
52+
auto max_iter = get_iter(*std::max_element(errors.begin(), errors.end(),
53+
[get_iter](const vtype p1, const vtype p2) { return get_iter(p1) < get_iter(p2); }));
54+
55+
vector<double> ub = { 5e-8, 5e-8, 5e-8, 5e-8 };
56+
for (auto& x: errors) {
57+
if (get_iter(x) == max_iter) {
58+
EXPECT_LE(get_error(x), ub[get_step(x)]);
59+
}
60+
}
61+
62+
int rank;
63+
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
64+
vector<size_t> ic = { 1, 1, 2, 2 };
65+
ASSERT_EQ(max_iter, (size_t) ic[rank]);
66+
}
67+
4368
int main(int argc, char** argv)
4469
{
4570
testing::InitGoogleTest(&argc, argv);

0 commit comments

Comments
 (0)