Skip to content

Commit 9232697

Browse files
committed
examples: advec: make mpi example take cmd options
1 parent 690b095 commit 9232697

File tree

2 files changed

+28
-12
lines changed

2 files changed

+28
-12
lines changed

examples/advection_diffusion/mpi_pfasst.cpp

Lines changed: 22 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -40,17 +40,17 @@ namespace pfasst
4040
*
4141
* @ingroup AdvectionDiffusion
4242
*/
43-
error_map run_mpi_pfasst(double abs_residual_tol, size_t niters=4)
43+
error_map run_mpi_pfasst(const double abs_res_tol, const double rel_res_tol,
44+
const size_t niters, const size_t nsteps, const double dt,
45+
const size_t ndofs_f, const size_t ndofs_c,
46+
const size_t nnodes_f, const size_t nnodes_c)
4447
{
45-
const size_t nsteps = 4;
46-
const double dt = 0.01;
47-
4848
vector<pair<size_t, quadrature::QuadratureType>> nodes = {
49-
{ 3, quadrature::QuadratureType::GaussLobatto },
50-
{ 5, quadrature::QuadratureType::GaussLobatto }
49+
{ nnodes_c, quadrature::QuadratureType::GaussLobatto },
50+
{ nnodes_f, quadrature::QuadratureType::GaussLobatto }
5151
};
5252

53-
vector<size_t> ndofs = { 64, 128 };
53+
vector<size_t> ndofs = { ndofs_c, ndofs_f };
5454

5555
auto build_level = [ndofs](size_t level) {
5656
auto factory = make_shared<MPIVectorFactory<double>>(ndofs[level]);
@@ -75,7 +75,7 @@ namespace pfasst
7575
pf.set_comm(&comm);
7676
pf.set_duration(0.0, nsteps * dt, dt, niters);
7777
pf.set_nsweeps({2, 1});
78-
pf.get_finest<AdvectionDiffusionSweeper<>>()->set_residual_tolerances(abs_residual_tol, 0.0);
78+
pf.get_finest<AdvectionDiffusionSweeper<>>()->set_residual_tolerances(abs_res_tol, rel_res_tol);
7979
pf.set_options();
8080
pf.run();
8181

@@ -94,7 +94,20 @@ int main(int argc, char** argv)
9494
pfasst::init(argc, argv,
9595
pfasst::examples::advection_diffusion::AdvectionDiffusionSweeper<>::init_opts,
9696
pfasst::examples::advection_diffusion::AdvectionDiffusionSweeper<>::init_logs);
97-
pfasst::examples::advection_diffusion::run_mpi_pfasst(0.0);
97+
98+
const double tend = pfasst::config::get_value<double>("tend", 0.04);
99+
const double dt = pfasst::config::get_value<double>("dt", 0.01);
100+
const size_t nnodes_f = pfasst::config::get_value<size_t>("num_nodes", 5);
101+
const size_t ndofs_f = pfasst::config::get_value<size_t>("spatial_dofs", 128);
102+
const size_t niters = pfasst::config::get_value<size_t>("num_iter", 4);
103+
const double abs_res_tol = pfasst::config::get_value<double>("abs_res_tol", 0.0);
104+
const double rel_res_tol = pfasst::config::get_value<double>("rel_res_tol", 0.0);
105+
106+
const size_t nsteps = tend / dt;
107+
const size_t nnodes_c = (nnodes_f + 1) / 2;
108+
const size_t ndofs_c = ndofs_f / 2;
109+
110+
pfasst::examples::advection_diffusion::run_mpi_pfasst(abs_res_tol, rel_res_tol, niters, nsteps, dt, ndofs_f, ndofs_c, nnodes_f, nnodes_c);
98111
fftw_cleanup();
99112
MPI_Finalize();
100113
}

tests/examples/advection_diffusion/test_mpi_advection_diffusion.cpp

Lines changed: 6 additions & 3 deletions
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(0.0);
27+
auto errors = run_mpi_pfasst(0.0, 0.0, 4, 4, 0.01, 128, 64, 5, 3);
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); };
@@ -44,7 +44,7 @@ TEST(AdaptiveErrorTest, MPIPFASST)
4444
{
4545
typedef error_map::value_type vtype;
4646

47-
auto errors = run_mpi_pfasst(1.e-8, 12);
47+
auto errors = run_mpi_pfasst(1.e-8, 0.0, 12, 4, 0.01, 128, 64, 5, 3);
4848
auto get_step = [](const vtype x) { return get<0>(get<0>(x)); };
4949
auto get_iter = [](const vtype x) { return get<1>(get<0>(x)); };
5050
auto get_error = [](const vtype x) { return get<1>(x); };
@@ -68,7 +68,10 @@ TEST(AdaptiveErrorTest, MPIPFASST)
6868
int main(int argc, char** argv)
6969
{
7070
testing::InitGoogleTest(&argc, argv);
71-
MPI_Init(&argc, &argv);
71+
pfasst::init(argc, argv,
72+
pfasst::examples::advection_diffusion::AdvectionDiffusionSweeper<>::init_opts,
73+
pfasst::examples::advection_diffusion::AdvectionDiffusionSweeper<>::init_logs,
74+
true);
7275
int result = 1, max_result; // GTest return value 1 (failure), 0 (success)
7376
result = RUN_ALL_TESTS();
7477
MPI_Allreduce(&result, &max_result, 1, MPI_INT, MPI_MAX, MPI_COMM_WORLD);

0 commit comments

Comments
 (0)