Medusa  1.1
Coordinate Free Mehless Method implementation
mm::AdamsBashforth< Scalar, num_steps > Class Template Reference

#include <AdamsBashforth_fwd.hpp>

Detailed Description

template<typename Scalar, int num_steps>
class mm::AdamsBashforth< Scalar, num_steps >

Class representing an AdamsBashforth method, an explicit linear multistep method.

Template Parameters
ScalarType of scalars to use for calculations, e.g. double.
num_stepsNumber of previous steps the method uses.

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();
}

Definition at line 29 of file AdamsBashforth_fwd.hpp.

Public Member Functions

 AdamsBashforth (const Eigen::Matrix< scalar_t, steps, 1 > &b)
 Construct an num_steps-step Adams Bashforth method with coefficients b, given as \(y_{n+s} = y_{n+s-1} + h \sum_{j=0}^{s-1} b_j f(t_j, y_j)\). More...
 
template<typename func_t >
Integrator< func_t > solve (const func_t &func, scalar_t t0, scalar_t t_max, scalar_t dt, const Eigen::VectorXd &y0) const
 Returns a solver using this method. More...
 
template<typename initial_method_t , typename func_t >
Integrator< func_t, initial_method_t > solve (const initial_method_t &method, const func_t &func, scalar_t t0, scalar_t t_max, scalar_t dt, const Eigen::VectorXd &y0) const
 Same as AdamsBashforth::solve, but offers the option of supplying your own initial method. More...
 

Public Types

enum  { steps = num_steps }
 
typedef Scalar scalar_t
 floating point type More...
 

Classes

class  Integrator
 Integrator class for AB methods. More...
 

Friends

template<typename scalar_t , int stages>
std::ostream & operator<< (std::ostream &os, const AdamsBashforth< scalar_t, stages > &)
 Output the method's tableau for debugging. More...
 

Private Member Functions

template<typename func_t >
Eigen::VectorXd step (const func_t &func, scalar_t t, const Eigen::Matrix< scalar_t, Eigen::Dynamic, steps > &last_ys, scalar_t dt) const
 Returns next value. More...
 

Private Attributes

Eigen::Matrix< scalar_t, steps, 1 > b_
 coefficients of the multistep method More...
 

Member Enumeration Documentation

◆ anonymous enum

template<typename Scalar , int num_steps>
anonymous enum
Enumerator
steps 

number of steps of the multistep method

Definition at line 32 of file AdamsBashforth_fwd.hpp.

Constructor & Destructor Documentation

◆ AdamsBashforth()

template<typename Scalar , int num_steps>
mm::AdamsBashforth< Scalar, num_steps >::AdamsBashforth ( const Eigen::Matrix< scalar_t, steps, 1 > &  b)
inline

Construct an num_steps-step Adams Bashforth method with coefficients b, given as \(y_{n+s} = y_{n+s-1} + h \sum_{j=0}^{s-1} b_j f(t_j, y_j)\).

Parameters
bVector of size num_steps representing the coefficients of the method.

Definition at line 44 of file AdamsBashforth_fwd.hpp.

Member Function Documentation

◆ solve() [1/2]

template<typename Scalar , int num_steps>
template<typename func_t >
Integrator<func_t> mm::AdamsBashforth< Scalar, num_steps >::solve ( const func_t &  func,
scalar_t  t0,
scalar_t  t_max,
scalar_t  dt,
const Eigen::VectorXd &  y0 
) const
inline

Returns a solver using this method.

Parameters
funcFunction to integrate.
t0Start time.
t_maxEnd time.
dtTime step.
y0Initial value.
Warning
Last step might not finish exactly on t_max. The number of steps is calculated as std::ceil((t_max - t0)/dt).
Returns
A solver object.

Definition at line 203 of file AdamsBashforth_fwd.hpp.

◆ solve() [2/2]

template<typename Scalar , int num_steps>
template<typename initial_method_t , typename func_t >
Integrator<func_t, initial_method_t> mm::AdamsBashforth< Scalar, num_steps >::solve ( const initial_method_t &  method,
const func_t &  func,
scalar_t  t0,
scalar_t  t_max,
scalar_t  dt,
const Eigen::VectorXd &  y0 
) const
inline

Same as AdamsBashforth::solve, but offers the option of supplying your own initial method.

Definition at line 211 of file AdamsBashforth_fwd.hpp.

◆ step()

template<typename Scalar , int num_steps>
template<typename func_t >
Eigen::VectorXd mm::AdamsBashforth< Scalar, num_steps >::step ( const func_t &  func,
scalar_t  t,
const Eigen::Matrix< scalar_t, Eigen::Dynamic, steps > &  last_ys,
scalar_t  dt 
) const
inlineprivate

Returns next value.

Definition at line 220 of file AdamsBashforth_fwd.hpp.

Friends And Related Function Documentation

◆ operator<<

template<typename Scalar , int num_steps>
template<typename scalar_t , int stages>
std::ostream& operator<< ( std::ostream &  os,
const AdamsBashforth< scalar_t, stages > &   
)
friend

Output the method's tableau for debugging.

Member Data Documentation

◆ b_

template<typename Scalar , int num_steps>
Eigen::Matrix<scalar_t, steps, 1> mm::AdamsBashforth< Scalar, num_steps >::b_
private

coefficients of the multistep method

Definition at line 36 of file AdamsBashforth_fwd.hpp.


The documentation for this class was generated from the following file:
mm::AdamsBashforth::step
Eigen::VectorXd step(const func_t &func, scalar_t t, const Eigen::Matrix< scalar_t, Eigen::Dynamic, steps > &last_ys, scalar_t dt) const
Returns next value.
Definition: AdamsBashforth_fwd.hpp:220
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