Medusa  1.1
Coordinate Free Mehless Method implementation
mm::integrators::ExplicitMultistep Namespace Reference

Detailed Description

Namespace containing factory functions for explicit linear multistep integrators.

Usage example:

/*
* Solving:
* x' = -y
* y' = x
*
* with solution (x, y) = (cos(t), sin(t))
*/
std::function<Eigen::VectorXd(double, const Eigen::VectorXd&)> func =
[](double, const Eigen::VectorXd& y) {
Eigen::VectorXd r(2);
r(0) = -y(1);
r(1) = y(0);
return r;
};
double tmax = 2*PI;
double dt = 0.1;
Eigen::VectorXd y0(2); y0 << 1.0, 0.0;
auto solver = integrators::ExplicitMultistep::AB3().solve(func, 0.0, tmax, dt, y0);
std::cout << solver << std::endl;
for (auto step : solver) {
std::cout << step << std::endl;
if (step.is_last()) {
std::cout << "This is the last step." << std::endl;
}
// Read/write access to current time and value.
step.time();
step.value();
}

Functions

template<class scalar_t = double>
static AdamsBashforth< scalar_t, 1 > AB1 ()
 Standard Euler's method. More...
 
template<class scalar_t = double>
static AdamsBashforth< scalar_t, 2 > AB2 ()
 Two step AB method. More...
 
template<class scalar_t = double>
static AdamsBashforth< scalar_t, 3 > AB3 ()
 Three step AB method. More...
 
template<class scalar_t = double>
static AdamsBashforth< scalar_t, 4 > AB4 ()
 Four step AB method. More...
 
template<class scalar_t = double>
static AdamsBashforth< scalar_t, 5 > AB5 ()
 Five step AB method. More...
 
template<int order, class scalar_t = double>
static AdamsBashforth< scalar_t, order > of_order ()
 Returns Adams-Bashforth explicit method of requested order with given floating point type. More...
 

Function Documentation

◆ AB1()

template<class scalar_t = double>
static AdamsBashforth<scalar_t, 1> mm::integrators::ExplicitMultistep::AB1 ( )
static

Standard Euler's method.

Examples
test/end2end/diffusion_explicit.cpp, and test/integrators/AdamsBashforth_test.cpp.

Definition at line 248 of file AdamsBashforth_fwd.hpp.

◆ AB2()

template<class scalar_t = double>
static AdamsBashforth<scalar_t, 2> mm::integrators::ExplicitMultistep::AB2 ( )
static

Two step AB method.

Examples
test/integrators/AdamsBashforth_test.cpp.

Definition at line 255 of file AdamsBashforth_fwd.hpp.

◆ AB3()

template<class scalar_t = double>
static AdamsBashforth<scalar_t, 3> mm::integrators::ExplicitMultistep::AB3 ( )
static

Three step AB method.

Examples
test/integrators/AdamsBashforth_test.cpp.

Definition at line 262 of file AdamsBashforth_fwd.hpp.

◆ AB4()

template<class scalar_t = double>
static AdamsBashforth<scalar_t, 4> mm::integrators::ExplicitMultistep::AB4 ( )
static

Four step AB method.

Examples
test/integrators/AdamsBashforth_test.cpp.

Definition at line 269 of file AdamsBashforth_fwd.hpp.

◆ AB5()

template<class scalar_t = double>
static AdamsBashforth<scalar_t, 5> mm::integrators::ExplicitMultistep::AB5 ( )
static

Five step AB method.

Examples
test/integrators/AdamsBashforth_test.cpp.

Definition at line 276 of file AdamsBashforth_fwd.hpp.

◆ of_order()

template<int order, class scalar_t = double>
static AdamsBashforth<scalar_t, order> mm::integrators::ExplicitMultistep::of_order ( )
static

Returns Adams-Bashforth explicit method of requested order with given floating point type.

Definition at line 296 of file AdamsBashforth_fwd.hpp.

mm::integrators::ExplicitMultistep::AB3
static AdamsBashforth< scalar_t, 3 > AB3()
Three step AB method.
Definition: AdamsBashforth_fwd.hpp:262
mm::PI
static const double PI
Mathematical constant pi in double precision.
Definition: Config.hpp:44