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

#include <RKExplicit_fwd.hpp>

Detailed Description

template<typename Scalar, int num_stages>
class mm::RKExplicit< Scalar, num_stages >

Class representing an explicit Runge-Kutta method.

Template Parameters
Scalartype of scalars to use for calculations, e.g. double
num_stagesnumber of stages of the method

Usage example:

std::function<Eigen::VectorXd(double, const Eigen::VectorXd&)> func =
[](double, const Eigen::VectorXd& y) {
return -y;
};
Eigen::VectorXd y0(1);
y0 << 1.0;
double tmax = 10.0;
double dt = 0.1;
auto integrator = integrators::Explicit::RK4().solve(func, 0.0, tmax, dt, y0);
std::cout << integrator << std::endl;
for (auto& step : integrator) {
// You can use read write access to:
step.value();
step.time();
// and check if this is the last step:
step.is_last();
}
// Aditionally, one can iterate manually
auto step = integrator.begin();
while (step) {
// do something with step
++step;
}
Eigen::VectorXd value = step.value(); // do something with value

Definition at line 32 of file RKExplicit_fwd.hpp.

Public Member Functions

 RKExplicit (const Eigen::Matrix< scalar_t, stages, 1 > &alpha, const Eigen::Matrix< scalar_t, stages, stages > &beta, const Eigen::Matrix< scalar_t, stages, 1 > &gamma)
 Construct a Runge-Kutta method with tableau \( \begin{array}{l|l} \alpha & \beta \\ \hline & \gamma^\mathsf{T} \end{array} \). 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...
 

Public Types

enum  { stages = num_stages }
 
typedef Scalar scalar_t
 Floating point type used in the computations. More...
 

Classes

class  Integrator
 Integrator class for RK methods. More...
 

Friends

template<typename Scalar2 , int num_steps2>
class AdamsBashforth
 Implements Adams Bashforth multistep methods. More...
 
template<typename scalar_t , int stages>
std::ostream & operator<< (std::ostream &os, const RKExplicit< 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::VectorXd &y, scalar_t dt) const
 Returns next value. More...
 

Private Attributes

Eigen::Matrix< scalar_t, stages, 1 > alpha_
 Left row of the tableau. More...
 
Eigen::Matrix< scalar_t, stages, stagesbeta_
 Central matrix of the tableau. More...
 
Eigen::Matrix< scalar_t, stages, 1 > gamma_
 Bottom row of the tableau. More...
 

Member Enumeration Documentation

◆ anonymous enum

template<typename Scalar , int num_stages>
anonymous enum
Enumerator
stages 

number of stages

Definition at line 35 of file RKExplicit_fwd.hpp.

Constructor & Destructor Documentation

◆ RKExplicit()

template<typename Scalar , int num_stages>
mm::RKExplicit< Scalar, num_stages >::RKExplicit ( const Eigen::Matrix< scalar_t, stages, 1 > &  alpha,
const Eigen::Matrix< scalar_t, stages, stages > &  beta,
const Eigen::Matrix< scalar_t, stages, 1 > &  gamma 
)
inline

Construct a Runge-Kutta method with tableau \( \begin{array}{l|l} \alpha & \beta \\ \hline & \gamma^\mathsf{T} \end{array} \).

Parameters
alphaVector of size stages.
betaMatrix of size stages times stages. Only strictly lower triangular part is used.
gammaVector of size stages.

Definition at line 56 of file RKExplicit_fwd.hpp.

Member Function Documentation

◆ solve()

template<typename Scalar , int num_stages>
template<typename func_t >
Integrator<func_t> mm::RKExplicit< Scalar, num_stages >::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 204 of file RKExplicit_fwd.hpp.

◆ step()

template<typename Scalar , int num_stages>
template<typename func_t >
Eigen::VectorXd mm::RKExplicit< Scalar, num_stages >::step ( const func_t &  func,
scalar_t  t,
const Eigen::VectorXd &  y,
scalar_t  dt 
) const
inlineprivate

Returns next value.

Definition at line 212 of file RKExplicit_fwd.hpp.

Friends And Related Function Documentation

◆ AdamsBashforth

template<typename Scalar , int num_stages>
template<typename Scalar2 , int num_steps2>
friend class AdamsBashforth
friend

Implements Adams Bashforth multistep methods.

Definition at line 230 of file RKExplicit_fwd.hpp.

◆ operator<<

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

Output the method's tableau for debugging.

Member Data Documentation

◆ alpha_

template<typename Scalar , int num_stages>
Eigen::Matrix<scalar_t, stages, 1> mm::RKExplicit< Scalar, num_stages >::alpha_
private

Left row of the tableau.

Definition at line 38 of file RKExplicit_fwd.hpp.

◆ beta_

template<typename Scalar , int num_stages>
Eigen::Matrix<scalar_t, stages, stages> mm::RKExplicit< Scalar, num_stages >::beta_
private

Central matrix of the tableau.

Definition at line 39 of file RKExplicit_fwd.hpp.

◆ gamma_

template<typename Scalar , int num_stages>
Eigen::Matrix<scalar_t, stages, 1> mm::RKExplicit< Scalar, num_stages >::gamma_
private

Bottom row of the tableau.

Definition at line 40 of file RKExplicit_fwd.hpp.


The documentation for this class was generated from the following file:
mm::RKExplicit::step
Eigen::VectorXd step(const func_t &func, scalar_t t, const Eigen::VectorXd &y, scalar_t dt) const
Returns next value.
Definition: RKExplicit_fwd.hpp:212
mm::integrators::Explicit::RK4
RKExplicit< scalar_t, 4 > RK4()
Standard RK4 4th order method.
Definition: RKExplicit.hpp:16