22 template <
typename MATRIX,
typename SCALAR =
double>
26 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
35 for (
size_t i = 0;
i < numSteps; ++
i)
43 const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
51 for (
size_t i = 0;
i < numSteps; ++
i)
67 virtual void do_step(
const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
80 template <
typename MATRIX,
typename SCALAR =
double>
84 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
88 virtual void do_step(
const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
93 rhs(stateInOut, derivative_, time);
94 stateInOut += dt * derivative_;
106 template <
typename MATRIX,
typename SCALAR =
double>
110 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
114 virtual void do_step(
const std::function<
void(
const MATRIX&, MATRIX&,
SCALAR)>& rhs,
120 SCALAR timePlusHalfStep = time + halfStep;
121 rhs(stateInOut, k1_, time);
122 rhs(stateInOut + halfStep * k1_, k2_, timePlusHalfStep);
123 rhs(stateInOut + halfStep * k2_, k3_, timePlusHalfStep);
124 rhs(stateInOut + dt * k3_, k4_, time + dt);
125 stateInOut += oneSixth_ * dt * (k1_ +
SCALAR(2.0) * k2_ +
SCALAR(2.0) * k3_ + k4_);
virtual void integrate_n_steps(std::function< void(const MATRIX &x, const SCALAR &t)> observe, 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: SteppersCT.h:42
Custom implementation of the euler stepper.
Definition: SteppersCT.h:81
This class serves as a common interface between the ODEInt and our custom integrators.
Definition: StepperBase.h:20
EIGEN_MAKE_ALIGNED_OPERATOR_NEW StepperRK4CT()
Definition: SteppersCT.h:112
CppAD::AD< CppAD::cg::CG< double > > SCALAR
virtual EIGEN_MAKE_ALIGNED_OPERATOR_NEW 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: SteppersCT.h:28
Custom implementation of the rk4 integration scheme.
Definition: SteppersCT.h:107
EIGEN_MAKE_ALIGNED_OPERATOR_NEW StepperEulerCT()
Definition: SteppersCT.h:86
virtual void do_step(const std::function< void(const MATRIX &, MATRIX &, SCALAR)> &rhs, MATRIX &stateInOut, const SCALAR time, const SCALAR dt)=0
Implements a single step of the integration scheme.
The stepper interface for custom steppers.
Definition: SteppersCT.h:23