- 3.0.2 core module.
SteppersODEInt.h
Go to the documentation of this file.
1 /**********************************************************************************************************************
2 This file is part of the Control Toolbox (https://github.com/ethz-adrl/control-toolbox), copyright by ETH Zurich.
3 Licensed under the BSD-2 license (see LICENSE file in main directory)
4 **********************************************************************************************************************/
5 
6 #pragma once
7 
8 #include <boost/numeric/odeint.hpp>
10 
11 namespace ct {
12 namespace core {
13 namespace internal {
14 
22 template <class STEPPER, typename MATRIX, typename SCALAR = double>
23 class StepperODEInt : public StepperBase<MATRIX, SCALAR>
24 {
25 public:
26  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
27 
29  virtual void integrate_n_steps(const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
30  MATRIX& state,
31  const SCALAR& startTime,
32  size_t numSteps,
33  SCALAR dt) override
34  {
35  boost::numeric::odeint::integrate_n_steps(stepper_, rhs, state, startTime, dt, numSteps);
36  }
37 
38  virtual void integrate_n_steps(std::function<void(const MATRIX& x, const SCALAR& t)> observer,
39  const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
40  MATRIX& state,
41  const SCALAR& startTime,
42  size_t numSteps,
43  SCALAR dt) override
44  {
45  boost::numeric::odeint::integrate_n_steps(stepper_, rhs, state, startTime, dt, numSteps, observer);
46  }
47 
48 
49  virtual void integrate_const(const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
50  MATRIX& state,
51  const SCALAR& startTime,
52  const SCALAR& finalTime,
53  SCALAR dt) override
54  {
55  boost::numeric::odeint::integrate_const(stepper_, rhs, state, startTime, finalTime, dt);
56  }
57 
58  virtual void integrate_const(std::function<void(const MATRIX& x, const SCALAR& t)> observer,
59  const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
60  MATRIX& state,
61  const SCALAR& startTime,
62  const SCALAR& finalTime,
63  SCALAR dt) override
64  {
65  boost::numeric::odeint::integrate_const(stepper_, rhs, state, startTime, finalTime, dt, observer);
66  }
67 
68  virtual void integrate_adaptive(const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
69  MATRIX& state,
70  const SCALAR& startTime,
71  const SCALAR& finalTime,
72  SCALAR dtInitial = SCALAR(0.01)) override
73  {
74  boost::numeric::odeint::integrate_adaptive(stepper_, rhs, state, startTime, finalTime, dtInitial);
75  }
76 
77  virtual void integrate_adaptive(std::function<void(const MATRIX& x, const SCALAR& t)> observer,
78  const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
79  MATRIX& state,
80  const SCALAR& startTime,
81  const SCALAR& finalTime,
82  const SCALAR dtInitial = SCALAR(0.01)) override
83  {
84  boost::numeric::odeint::integrate_adaptive(stepper_, rhs, state, startTime, finalTime, dtInitial, observer);
85  }
86 
87  virtual void integrate_times(std::function<void(const MATRIX& x, const SCALAR& t)> observer,
88  const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
89  MATRIX& state,
90  const tpl::TimeArray<SCALAR>& timeTrajectory,
91  SCALAR dtInitial = SCALAR(0.01)) override
92  {
93  boost::numeric::odeint::integrate_times(
94  stepper_, rhs, state, &timeTrajectory.front(), &timeTrajectory.back() + 1, dtInitial, observer);
95  }
96 
97 private:
98  STEPPER stepper_;
99 };
100 
101 
109 template <class STEPPER, typename MATRIX, typename SCALAR = double>
110 class StepperODEIntDenseOutput : public StepperODEInt<STEPPER, MATRIX, SCALAR>
111 {
112 public:
113  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
114  typedef typename boost::numeric::odeint::result_of::make_dense_output<STEPPER>::type StepperDense;
115 
117  {
118  stepperDense_ = boost::numeric::odeint::make_dense_output(this->absErrTol_, this->relErrTol_, stepper_);
119  }
120 
121  virtual void integrate_adaptive(const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
122  MATRIX& state,
123  const SCALAR& startTime,
124  const SCALAR& finalTime,
125  SCALAR dtInitial = SCALAR(0.01)) override
126  {
127  stepperDense_.initialize(state, startTime, dtInitial);
128  boost::numeric::odeint::integrate_adaptive(stepperDense_, rhs, state, startTime, finalTime, dtInitial);
129  }
130 
131  virtual void integrate_adaptive(std::function<void(const MATRIX& x, const SCALAR& t)> observer,
132  const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
133  MATRIX& state,
134  const SCALAR& startTime,
135  const SCALAR& finalTime,
136  const SCALAR dtInitial = SCALAR(0.01)) override
137  {
138  stepperDense_.initialize(state, startTime, dtInitial);
139  boost::numeric::odeint::integrate_adaptive(
140  stepperDense_, rhs, state, startTime, finalTime, dtInitial, observer);
141  }
142 
143  virtual void integrate_times(std::function<void(const MATRIX& x, const SCALAR& t)> observer,
144  const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
145  MATRIX& state,
146  const tpl::TimeArray<SCALAR>& timeTrajectory,
147  SCALAR dtInitial = SCALAR(0.01)) override
148  {
149  stepperDense_.initialize(state, timeTrajectory.front(), dtInitial);
150  boost::numeric::odeint::integrate_times(
151  stepperDense_, rhs, state, &timeTrajectory.front(), &timeTrajectory.back() + 1, dtInitial, observer);
152  }
153 
154 private:
155  STEPPER stepper_;
156  StepperDense stepperDense_;
157 };
158 
166 template <class STEPPER, typename MATRIX, typename SCALAR = double>
167 class StepperODEIntControlled : public StepperODEInt<STEPPER, MATRIX, SCALAR>
168 {
169 public:
170  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
171  typedef typename boost::numeric::odeint::result_of::make_controlled<STEPPER>::type StepperControlled;
172 
174  {
175  stepperControlled_ = boost::numeric::odeint::make_controlled(this->absErrTol_, this->relErrTol_, stepper_);
176  }
177 
178  virtual void integrate_adaptive(const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
179  MATRIX& state,
180  const SCALAR& startTime,
181  const SCALAR& finalTime,
182  SCALAR dtInitial = SCALAR(0.01)) override
183  {
184  boost::numeric::odeint::integrate_adaptive(stepperControlled_, rhs, state, startTime, finalTime, dtInitial);
185  }
186 
187  virtual void integrate_adaptive(std::function<void(const MATRIX& x, const SCALAR& t)> observer,
188  const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
189  MATRIX& state,
190  const SCALAR& startTime,
191  const SCALAR& finalTime,
192  const SCALAR dtInitial = SCALAR(0.01)) override
193  {
194  boost::numeric::odeint::integrate_adaptive(
195  stepperControlled_, rhs, state, startTime, finalTime, dtInitial, observer);
196  }
197 
198  virtual void integrate_times(std::function<void(const MATRIX& x, const SCALAR& t)> observer,
199  const std::function<void(const MATRIX&, MATRIX&, SCALAR)>& rhs,
200  MATRIX& state,
201  const tpl::TimeArray<SCALAR>& timeTrajectory,
202  SCALAR dtInitial = SCALAR(0.01)) override
203  {
204  boost::numeric::odeint::integrate_times(
205  stepperControlled_, rhs, state, &timeTrajectory.front(), &timeTrajectory.back() + 1, dtInitial, observer);
206  }
207 
208 private:
209  STEPPER stepper_;
210  StepperControlled stepperControlled_;
211 };
212 }
213 }
214 }
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