24 template<
typename scalar>
26 const size_t nparticles,
const size_t niters,
27 const double abs_res_tol,
const double rel_res_tol)
33 const double mass = 1.0;
34 const double charge = 1.0;
36 shared_ptr<bindings::WrapperInterface<double, double>> impl_solver = \
37 make_shared<bindings::WrapperSimplePhysicsSolver<double, double>>();
41 auto quad1 = quadrature::quadrature_factory<double>(nnodes,
43 auto factory1 = make_shared<ParticleCloudFactory<double>>(nparticles, 3, mass, charge);
44 string data_file1 =
"s" + to_string(nsteps) +
"_i" + to_string(niters)
45 +
"_dt" + to_string(dt) +
"_m" + to_string(nnodes)
46 +
"_p" + to_string(nparticles) +
"_level1.csv";
47 auto sweeper1 = make_shared<BorisSweeper<double, double>>(impl_solver, data_file1);
48 auto transfer1 = make_shared<InjectiveTransfer<double, double>>();
49 sweeper1->set_quadrature(quad1);
50 sweeper1->set_factory(factory1);
51 sweeper1->set_residual_tolerances(abs_res_tol, rel_res_tol);
52 controller.
add_level(sweeper1, transfer1);
55 auto quad2 = quadrature::quadrature_factory<double>(nnodes,
57 auto factory2 = make_shared<ParticleCloudFactory<double>>(nparticles, 3, mass, charge);
58 string data_file2 =
"s" + to_string(nsteps) +
"_i" + to_string(niters)
59 +
"_dt" + to_string(dt) +
"_m" + to_string(nnodes)
60 +
"_p" + to_string(nparticles) +
"_level2.csv";
61 auto sweeper2 = make_shared<BorisSweeper<double, double>>(impl_solver, data_file2);
62 auto transfer2 = make_shared<InjectiveTransfer<double, double>>();
63 sweeper2->set_quadrature(quad2);
64 sweeper2->set_factory(factory2);
65 controller.
add_level(sweeper2, transfer2);
71 shared_ptr<Particle<double>> center = make_shared<Particle<double>>();
72 center->pos()[0] = 10;
73 center->vel()[0] = 100;
74 center->vel()[2] = 100;
77 shared_ptr<ParticleCloud<double>> q0 = \
79 q0->distribute_around_center(center);
80 ML_CLOG(INFO,
"Boris",
"Initial Particle (fine) : "
82 fine_sweeper->set_initial_energy();
86 return fine_sweeper->get_errors();
92 #ifndef PFASST_UNIT_TESTING
93 int main(
int argc,
char** argv)
95 MPI_Init(&argc, &argv);
97 pfasst::examples::boris::init_opts<>,
98 pfasst::examples::boris::init_logs<>);
100 const size_t nsteps = pfasst::config::get_value<size_t>(
"num_steps", 1);
101 const double dt = pfasst::config::get_value<double>(
"delta_step", 0.015625);
102 const size_t nnodes = pfasst::config::get_value<size_t>(
"num_nodes", 5);
103 const size_t nparticles = pfasst::config::get_value<size_t>(
"num_particles", 1);
104 const size_t niters = pfasst::config::get_value<size_t>(
"num_iter", 2);
105 const double abs_res_tol = pfasst::config::get_value<double>(
"abs_res_tol", 0.0);
106 const double rel_res_tol = pfasst::config::get_value<double>(
"rel_res_tol", 0.0);
108 ML_CLOG(INFO,
"Boris",
"nsteps=" << nsteps <<
", "
109 <<
"dt=" << dt <<
", "
110 <<
"nnodes=" << nnodes <<
", "
111 <<
"nparticles=" << nparticles <<
", "
112 <<
"niter=" << niters <<
", "
113 <<
"abs res=" << abs_res_tol <<
", "
114 <<
"rel res=" << rel_res_tol);
116 pfasst::examples::boris::run_boris_pfasst<double>(nsteps,
dt, nnodes, nparticles, niters, abs_res_tol, rel_res_tol);
int main(int argc, char **argv)
map< error_index, ErrorTuple< scalar >> error_map
void setup(shared_ptr< WrapperInterface< scalar, time >> wrapper)
virtual void set_options(bool all_sweepers=true)
Set options from command line etc.
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++
Implementation of the PFASST algorithm as described in .
static void init(int argc, char **argv, std::function< void()> opts=nullptr, std::function< void()> logs=nullptr)
virtual void set_duration(time t0, time tend, time dt, size_t niters)
Set basic time scope of the Controller.
virtual void setup() override
Basic setup routine for this controller.
shared_ptr< R > get_finest()
Get coarsest level.
virtual void set_comm(ICommunicator *comm)
Set communicator.
error_map< scalar > run_boris_pfasst(const size_t nsteps, const scalar dt, const size_t nnodes, const size_t nparticles, const size_t niters, const double abs_res_tol, const double rel_res_tol)
virtual void run() override
Solve ODE using PFASST.