ExplicitStepperWrapper< FP, ExplicitStepper > Class Template Reference

CPP API: mio::ExplicitStepperWrapper< FP, ExplicitStepper > Class Template Reference
mio::ExplicitStepperWrapper< FP, ExplicitStepper > Class Template Reference

This is a non-adaptive IntegratorCore. More...

#include <stepper_wrapper.h>

Inheritance diagram for mio::ExplicitStepperWrapper< FP, ExplicitStepper >:
Collaboration diagram for mio::ExplicitStepperWrapper< FP, ExplicitStepper >:

Public Types

using Stepper = ExplicitStepper< Eigen::VectorX< FP >, FP, Eigen::VectorX< FP >, FP, boost::numeric::odeint::vector_space_algebra, typename boost::numeric::odeint::operations_dispatcher< Eigen::VectorX< FP > >::operations_type, boost::numeric::odeint::initially_resizer >
 

Public Member Functions

std::unique_ptr< OdeIntegratorCore< FP > > clone () const override
 
 ExplicitStepperWrapper ()
 Set up the integrator. More...
 
bool step (const mio::DerivFunction< FP > &f, Eigen::Ref< const Eigen::VectorX< FP >> yt, FP &t, FP &dt, Eigen::Ref< Eigen::VectorX< FP >> ytp1) const override
 Make a single integration step on a system of ODEs with fixed step size dt. More...
 
- Public Member Functions inherited from mio::IntegratorCore< FP, Integrands >
 IntegratorCore (const FP &dt_min, const FP &dt_max)
 Initialize an IntegratorCore. More...
 
virtual bool step (const Integrands &... fs, Eigen::Ref< const Eigen::VectorX< FP >> yt, FP &t, FP &dt, Eigen::Ref< Eigen::VectorX< FP >> ytp1) const =0
 Make a single integration step. More...
 
virtual ~IntegratorCore ()
 
FP & get_dt_min ()
 Access lower bound to the step size dt. More...
 
const FP & get_dt_min () const
 Access lower bound to the step size dt. More...
 
FP & get_dt_max ()
 Access upper bound to the step size dt. More...
 
const FP & get_dt_max () const
 Access upper bound to the step size dt. More...
 

Private Attributes

Stepper m_stepper
 A stepper instance used for integration. More...
 

Detailed Description

template<typename FP, template< class State, class Value, class Deriv, class Time, class Algebra, class Operations, class Resizer > class ExplicitStepper>
class mio::ExplicitStepperWrapper< FP, ExplicitStepper >

This is a non-adaptive IntegratorCore.

It creates and manages an instance of an explicit stepper from boost::numeric::odeint, wrapped as mio::IntegratorCore.

Member Typedef Documentation

◆ Stepper

template<typename FP , template< class State, class Value, class Deriv, class Time, class Algebra, class Operations, class Resizer > class ExplicitStepper>
using mio::ExplicitStepperWrapper< FP, ExplicitStepper >::Stepper = ExplicitStepper<Eigen::VectorX<FP>, FP, Eigen::VectorX<FP>, FP, boost::numeric::odeint::vector_space_algebra, typename boost::numeric::odeint::operations_dispatcher<Eigen::VectorX<FP> >::operations_type, boost::numeric::odeint::initially_resizer>

Constructor & Destructor Documentation

◆ ExplicitStepperWrapper()

template<typename FP , template< class State, class Value, class Deriv, class Time, class Algebra, class Operations, class Resizer > class ExplicitStepper>
mio::ExplicitStepperWrapper< FP, ExplicitStepper >::ExplicitStepperWrapper ( )
inline

Set up the integrator.

Member Function Documentation

◆ clone()

template<typename FP , template< class State, class Value, class Deriv, class Time, class Algebra, class Operations, class Resizer > class ExplicitStepper>
std::unique_ptr<OdeIntegratorCore<FP> > mio::ExplicitStepperWrapper< FP, ExplicitStepper >::clone ( ) const
inlineoverridevirtual

◆ step()

template<typename FP , template< class State, class Value, class Deriv, class Time, class Algebra, class Operations, class Resizer > class ExplicitStepper>
bool mio::ExplicitStepperWrapper< FP, ExplicitStepper >::step ( const mio::DerivFunction< FP > &  f,
Eigen::Ref< const Eigen::VectorX< FP >>  yt,
FP &  t,
FP &  dt,
Eigen::Ref< Eigen::VectorX< FP >>  ytp1 
) const
inlineoverride

Make a single integration step on a system of ODEs with fixed step size dt.

Parameters
[in]ytValue of y at t, y(t).
[in,out]tCurrent time. Overwritten with t+dt.
[in]dtCurrent time step size h=dt.
[out]ytp1The approximated value of y(t+dt).

Member Data Documentation

◆ m_stepper

template<typename FP , template< class State, class Value, class Deriv, class Time, class Algebra, class Operations, class Resizer > class ExplicitStepper>
Stepper mio::ExplicitStepperWrapper< FP, ExplicitStepper >::m_stepper
mutableprivate

A stepper instance used for integration.