8 #include <boost/numeric/odeint.hpp> 22 template <
class STEPPER,
typename MATRIX,
typename SCALAR =
double>
26 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
35 boost::numeric::odeint::integrate_n_steps(stepper_, rhs, state, startTime, dt, numSteps);
39 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
45 boost::numeric::odeint::integrate_n_steps(stepper_, rhs, state, startTime, dt, numSteps, observer);
55 boost::numeric::odeint::integrate_const(stepper_, rhs, state, startTime, finalTime, dt);
59 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
65 boost::numeric::odeint::integrate_const(stepper_, rhs, state, startTime, finalTime, dt, observer);
74 boost::numeric::odeint::integrate_adaptive(stepper_, rhs, state, startTime, finalTime, dtInitial);
78 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
84 boost::numeric::odeint::integrate_adaptive(stepper_, rhs, state, startTime, finalTime, dtInitial, observer);
88 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
93 boost::numeric::odeint::integrate_times(
94 stepper_, rhs, state, &timeTrajectory.front(), &timeTrajectory.back() + 1, dtInitial, observer);
109 template <
class STEPPER,
typename MATRIX,
typename SCALAR =
double>
113 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
114 typedef typename boost::numeric::odeint::result_of::make_dense_output<STEPPER>::type
StepperDense;
118 stepperDense_ = boost::numeric::odeint::make_dense_output(this->
absErrTol_, this->
relErrTol_, stepper_);
127 stepperDense_.initialize(state, startTime, dtInitial);
128 boost::numeric::odeint::integrate_adaptive(stepperDense_, rhs, state, startTime, finalTime, dtInitial);
132 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
138 stepperDense_.initialize(state, startTime, dtInitial);
139 boost::numeric::odeint::integrate_adaptive(
140 stepperDense_, rhs, state, startTime, finalTime, dtInitial, observer);
144 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
149 stepperDense_.initialize(state, timeTrajectory.front(), dtInitial);
150 boost::numeric::odeint::integrate_times(
151 stepperDense_, rhs, state, &timeTrajectory.front(), &timeTrajectory.back() + 1, dtInitial, observer);
156 StepperDense stepperDense_;
166 template <
class STEPPER,
typename MATRIX,
typename SCALAR =
double>
170 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
171 typedef typename boost::numeric::odeint::result_of::make_controlled<STEPPER>::type
StepperControlled;
175 stepperControlled_ = boost::numeric::odeint::make_controlled(this->
absErrTol_, this->
relErrTol_, stepper_);
184 boost::numeric::odeint::integrate_adaptive(stepperControlled_, rhs, state, startTime, finalTime, dtInitial);
188 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
194 boost::numeric::odeint::integrate_adaptive(
195 stepperControlled_, rhs, state, startTime, finalTime, dtInitial, observer);
199 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
204 boost::numeric::odeint::integrate_times(
205 stepperControlled_, rhs, state, &timeTrajectory.front(), &timeTrajectory.back() + 1, dtInitial, observer);
210 StepperControlled stepperControlled_;
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)) override
Integrates forward in time from an initial to a final time. If an adaptive stepper is used...
Definition: SteppersODEInt.h:121
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)) override
Integrates forward in time from an initial to a final time. If an adaptive stepper is used...
Definition: SteppersODEInt.h:187
StepperODEIntControlled()
Definition: SteppersODEInt.h:173
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) override
Performs numSteps integration steps.
Definition: SteppersODEInt.h:38
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) override
Equidistant integration based on initial and final time as well as step length.
Definition: SteppersODEInt.h:58
EIGEN_MAKE_ALIGNED_OPERATOR_NEW typedef boost::numeric::odeint::result_of::make_dense_output< STEPPER >::type StepperDense
Definition: SteppersODEInt.h:114
The interface to call ODEInt Dense Output Integration routines.
Definition: SteppersODEInt.h:110
This class serves as a common interface between the ODEInt and our custom integrators.
Definition: StepperBase.h:20
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)) override
Integrates forward in time from an initial to a final time. If an adaptive stepper is used...
Definition: SteppersODEInt.h:68
CppAD::AD< CppAD::cg::CG< double > > SCALAR
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)) override
Integrates forward in time from an initial to a final time. If an adaptive stepper is used...
Definition: SteppersODEInt.h:131
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)) override
Integrates forward in time from an initial to a final time. If an adaptive stepper is used...
Definition: SteppersODEInt.h:77
virtual void integrate_const(const std::function< void(const MATRIX &, MATRIX &, SCALAR)> &rhs, MATRIX &state, const SCALAR &startTime, const SCALAR &finalTime, SCALAR dt) override
Equidistant integration based on initial and final time as well as step length.
Definition: SteppersODEInt.h:49
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)) override
Integrates a system using a given time sequence.
Definition: SteppersODEInt.h:143
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)) override
Integrates a system using a given time sequence.
Definition: SteppersODEInt.h:198
An array in time.
Definition: TimeArray.h:22
StepperODEIntDenseOutput()
Definition: SteppersODEInt.h:116
The interface to call ODEInt Controlled integration routines.
Definition: SteppersODEInt.h:167
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)) override
Integrates a system using a given time sequence.
Definition: SteppersODEInt.h:87
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)) override
Integrates forward in time from an initial to a final time. If an adaptive stepper is used...
Definition: SteppersODEInt.h:178
EIGEN_MAKE_ALIGNED_OPERATOR_NEW StepperODEInt()
Definition: SteppersODEInt.h:28
virtual void integrate_n_steps(const std::function< void(const MATRIX &, MATRIX &, SCALAR)> &rhs, MATRIX &state, const SCALAR &startTime, size_t numSteps, SCALAR dt) override
Performs numSteps integration steps.
Definition: SteppersODEInt.h:29
SCALAR absErrTol_
Definition: StepperBase.h:181
SCALAR relErrTol_
Definition: StepperBase.h:182
EIGEN_MAKE_ALIGNED_OPERATOR_NEW typedef boost::numeric::odeint::result_of::make_controlled< STEPPER >::type StepperControlled
Definition: SteppersODEInt.h:171
The interface to call the integration routines from ODEInt.
Definition: SteppersODEInt.h:23