PFASST++
mpi_pfasst.cpp
Go to the documentation of this file.
1 
9 #include <cassert>
10 #include <memory>
11 #include <vector>
12 #include <utility>
13 using namespace std;
14 
15 #include <mpi.h>
16 #include <fftw3.h>
17 
18 #include <pfasst.hpp>
19 #include <pfasst/logging.hpp>
22 #include <pfasst/encap/vector.hpp>
23 
25 #include "spectral_transfer_1d.hpp"
26 
27 using namespace pfasst::encap;
28 using namespace pfasst::mpi;
29 
30 namespace pfasst
31 {
32  namespace examples
33  {
34  namespace advection_diffusion
35  {
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)
47  {
48  ML_CLOG(INFO, "Advec", "abs_res_tol: " << abs_res_tol << ", "
49  << "rel_res_tol: " << rel_res_tol << ", "
50  << "niter: " << niters << ", "
51  << "nsteps: " << nsteps << ", "
52  << "dt: " << dt << ", "
53  << "ndofs (f-c): " << ndofs_f << "-" << ndofs_c << ", "
54  << "nnodes (f-c): " << nnodes_f << "-" << nnodes_c);
55 
56  MPICommunicator comm(MPI_COMM_WORLD);
57  PFASST<> pf;
58 
59  auto quad_c = quadrature::quadrature_factory(nnodes_c, quadrature::QuadratureType::GaussLobatto);
60  auto factory_c = make_shared<VectorFactory<double>>(ndofs_c);
61  auto sweeper_c = make_shared<AdvectionDiffusionSweeper<>>(ndofs_c);
62  auto transfer_c = make_shared<SpectralTransfer1D<>>();
63 
64  sweeper_c->set_quadrature(quad_c);
65  sweeper_c->set_factory(factory_c);
66  sweeper_c->set_residual_tolerances(abs_res_tol, rel_res_tol);
67 
68  auto quad_f = quadrature::quadrature_factory(nnodes_f, quadrature::QuadratureType::GaussLobatto);
69  auto factory_f = make_shared<VectorFactory<double>>(ndofs_f);
70  auto sweeper_f = make_shared<AdvectionDiffusionSweeper<>>(ndofs_f);
71  auto transfer_f = make_shared<SpectralTransfer1D<>>();
72 
73  sweeper_f->set_quadrature(quad_f);
74  sweeper_f->set_factory(factory_f);
75  sweeper_f->set_residual_tolerances(abs_res_tol, rel_res_tol);
76 
77  ML_LOG(INFO, "expected quadrature error: " << quad_c->expected_error() << " (" << nnodes_c << ")");
78  ML_LOG(INFO, "expected quadrature error: " << quad_f->expected_error() << " (" << nnodes_f << ")");
79 
80  pf.add_level(sweeper_f, transfer_f);
81  pf.add_level(sweeper_c, transfer_c);
82  pf.setup();
83 
84  auto q0 = sweeper_f->get_start_state();
85  sweeper_f->exact(q0, 0.0);
86 
87  pf.set_comm(&comm);
88  pf.set_duration(0.0, nsteps * dt, dt, niters);
89  pf.set_nsweeps({2, 1});
90  pf.get_finest<AdvectionDiffusionSweeper<>>()->set_residual_tolerances(abs_res_tol, rel_res_tol);
91  pf.set_options();
92  pf.run();
93 
94  auto fine = pf.get_finest<AdvectionDiffusionSweeper<>>();
95  return fine->get_errors();
96  }
97  } // ::pfasst::examples::advection_diffusion
98  } // ::pfasst::examples
99 } // ::pfasst
100 
101 
102 #ifndef PFASST_UNIT_TESTING
103 int main(int argc, char** argv)
104 {
105  MPI_Init(&argc, &argv);
106  pfasst::init(argc, argv,
109 
110  const double tend = pfasst::config::get_value<double>("tend", 0.04);
111  const double dt = pfasst::config::get_value<double>("dt", 0.01);
112  const size_t nnodes_f = pfasst::config::get_value<size_t>("num_nodes", 5);
113  const size_t ndofs_f = pfasst::config::get_value<size_t>("spatial_dofs", 128);
114  const size_t niters = pfasst::config::get_value<size_t>("num_iter", 4);
115  const double abs_res_tol = pfasst::config::get_value<double>("abs_res_tol", 0.0);
116  const double rel_res_tol = pfasst::config::get_value<double>("rel_res_tol", 0.0);
117 
118  const size_t nsteps = tend / dt;
119  const size_t nnodes_c = (nnodes_f + 1) / 2;
120  const size_t ndofs_c = ndofs_f / 2;
121 
123  niters, nsteps, dt,
124  ndofs_f, ndofs_c, nnodes_f, nnodes_c);
125  fftw_cleanup();
126  MPI_Finalize();
127 }
128 #endif
virtual void set_options(bool all_sweepers=true)
Set options from command line etc.
STL namespace.
map< tuple< size_t, size_t >, double > error_map
virtual void add_level(shared_ptr< ISweeper< time >> sweeper, shared_ptr< ITransfer< time >> transfer=shared_ptr< ITransfer< time >>(nullptr), bool coarse=true)
Adding a level to the controller.
#define ML_CLOG(level, logger_id, x)
same as CLOG(level, logger, x) from easylogging++
Definition: logging.hpp:117
Implementation of the PFASST algorithm as described in .
Definition: pfasst.hpp:19
static void init(int argc, char **argv, std::function< void()> opts=nullptr, std::function< void()> logs=nullptr)
Definition: pfasst.hpp:13
virtual void set_nsweeps(vector< size_t > nsweeps)
Set desired number of sweeps for each level independently.
Definition: mlsdc_impl.hpp:56
#define ML_LOG(level, x)
same as LOG(level, x) from easylogging++
Definition: logging.hpp:110
virtual void set_duration(time t0, time tend, time dt, size_t niters)
Set basic time scope of the Controller.
error_map run_mpi_pfasst(const double abs_res_tol, const double rel_res_tol, const size_t niters, const size_t nsteps, const double dt, const size_t ndofs_f, const size_t ndofs_c, const size_t nnodes_f, const size_t nnodes_c)
Advection/diffusion example using an encapsulated IMEX sweeper.
Definition: mpi_pfasst.cpp:43
virtual void setup() override
Basic setup routine for this controller.
Definition: mlsdc_impl.hpp:45
shared_ptr< IQuadrature< precision > > quadrature_factory(const size_t nnodes, const QuadratureType qtype)
Instantiates quadrature handler for given number of nodes and type descriptor.
Definition: quadrature.hpp:53
shared_ptr< R > get_finest()
Get coarsest level.
Definition: interface.hpp:163
virtual void set_comm(ICommunicator *comm)
Set communicator.
Definition: pfasst_impl.hpp:30
Encapsulations (short encaps) are the central data type for all PFASST++ algorithms.
int main(int argc, char **argv)
Definition: mpi_pfasst.cpp:103
advection-diffusion sweeper with semi-implicit time-integration.
virtual void run() override
Solve ODE using PFASST.
Definition: pfasst_impl.hpp:40
float dt
Definition: plot.py:10