19 template <
typename MATRIX,
typename SCALAR =
double>
23 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
42 throw std::runtime_error(
"Integrate_n_steps not implemented for the stepper type");
56 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
62 throw std::runtime_error(
"Integrate_n_steps not implemented for the stepper type");
81 throw std::runtime_error(
"integrate_const not implemented for the stepper type");
95 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
101 throw std::runtime_error(
"integrate_const not implemented for the stepper type");
123 throw std::runtime_error(
"integrate_adaptive not implemented for the stepper type");
141 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
147 throw std::runtime_error(
"integrate_adaptive not implemented for the stepper type");
160 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
165 throw std::runtime_error(
"integrate_times not implemented for the stepper type");
virtual void integrate_const(std::function< void(const MATRIX &x, const SCALAR &t)> observer, const std::function< void(const MATRIX &, MATRIX &, SCALAR)> &rhs, MATRIX &state, const SCALAR &startTime, const SCALAR &finalTime, SCALAR dt)
Equidistant integration based on initial and final time as well as step length.
Definition: StepperBase.h:94
virtual void integrate_n_steps(const std::function< void(const MATRIX &, MATRIX &, SCALAR)> &rhs, MATRIX &state, const SCALAR &startTime, size_t numSteps, SCALAR dt)
Performs numSteps integration steps.
Definition: StepperBase.h:36
virtual void integrate_n_steps(std::function< void(const MATRIX &x, const SCALAR &t)> observer, const std::function< void(const MATRIX &, MATRIX &, SCALAR)> &rhs, MATRIX &state, const SCALAR &startTime, size_t numSteps, SCALAR dt)
Performs numSteps integration steps.
Definition: StepperBase.h:55
This class serves as a common interface between the ODEInt and our custom integrators.
Definition: StepperBase.h:20
void setAdaptiveErrorTolerances(const SCALAR absErrTol, const SCALAR &relErrTol)
Sets the adaptive error tolerances.
Definition: StepperBase.h:174
CppAD::AD< CppAD::cg::CG< double > > SCALAR
virtual ~StepperBase()
Definition: StepperBase.h:26
An array in time.
Definition: TimeArray.h:22
virtual void integrate_adaptive(const std::function< void(const MATRIX &, MATRIX &, SCALAR)> &rhs, MATRIX &state, const SCALAR &startTime, const SCALAR &finalTime, SCALAR dtInitial=SCALAR(0.01))
Integrates forward in time from an initial to a final time. If an adaptive stepper is used...
Definition: StepperBase.h:117
virtual void integrate_const(const std::function< void(const MATRIX &, MATRIX &, SCALAR)> &rhs, MATRIX &state, const SCALAR &startTime, const SCALAR &finalTime, SCALAR dt)
Equidistant integration based on initial and final time as well as step length.
Definition: StepperBase.h:75
virtual void integrate_adaptive(std::function< void(const MATRIX &x, const SCALAR &t)> observer, const std::function< void(const MATRIX &, MATRIX &, SCALAR)> &rhs, MATRIX &state, const SCALAR &startTime, const SCALAR &finalTime, const SCALAR dtInitial=SCALAR(0.01))
Integrates forward in time from an initial to a final time. If an adaptive stepper is used...
Definition: StepperBase.h:140
virtual void integrate_times(std::function< void(const MATRIX &x, const SCALAR &t)> observer, const std::function< void(const MATRIX &, MATRIX &, SCALAR)> &rhs, MATRIX &state, const tpl::TimeArray< SCALAR > &timeTrajectory, SCALAR dtInitial=SCALAR(0.01))
Integrates a system using a given time sequence.
Definition: StepperBase.h:159
SCALAR absErrTol_
Definition: StepperBase.h:181
EIGEN_MAKE_ALIGNED_OPERATOR_NEW StepperBase()
Definition: StepperBase.h:25
SCALAR relErrTol_
Definition: StepperBase.h:182