stepper_wrapper.h Source File

CPP API: stepper_wrapper.h Source File
stepper_wrapper.h
Go to the documentation of this file.
1 /*
2 * Copyright (C) 2020-2026 MEmilio
3 *
4 * Authors: Rene Schmieding
5 *
6 * Contact: Martin J. Kuehn <Martin.Kuehn@DLR.de>
7 *
8 * Licensed under the Apache License, Version 2.0 (the "License");
9 * you may not use this file except in compliance with the License.
10 * You may obtain a copy of the License at
11 *
12 * http://www.apache.org/licenses/LICENSE-2.0
13 *
14 * Unless required by applicable law or agreed to in writing, software
15 * distributed under the License is distributed on an "AS IS" BASIS,
16 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17 * See the License for the specific language governing permissions and
18 * limitations under the License.
19 */
20 #ifndef MIO_MATH_STEPPER_WRAPPER_H
21 #define MIO_MATH_STEPPER_WRAPPER_H
22 
24 #include "memilio/utils/logging.h"
25 
26 #include "boost/numeric/odeint/external/eigen/eigen_algebra.hpp"
27 #include "boost/numeric/odeint/external/eigen/eigen_resize.hpp"
28 #include "boost/numeric/odeint/stepper/controlled_runge_kutta.hpp"
29 #include "boost/numeric/odeint/stepper/runge_kutta_fehlberg78.hpp"
30 #include "boost/numeric/odeint/stepper/runge_kutta_cash_karp54.hpp"
31 #include <boost/core/ref.hpp>
32 // #include "boost/numeric/odeint/stepper/runge_kutta_dopri5.hpp" // TODO: reenable once boost bug is fixed
33 
34 namespace mio
35 {
36 
37 namespace details
38 {
39 
41 template <class Value, class Time>
42 struct step_adjuster : public boost::numeric::odeint::default_step_adjuster<Value, Time> {
43  using boost::numeric::odeint::default_step_adjuster<Value, Time>::default_step_adjuster;
44  void set_dt_max(const Time& dt_max)
45  {
46  this->m_max_dt = dt_max;
47  }
48 };
49 
50 } // namespace details
51 
56 template <typename FP, template <class State, class Value, class Deriv, class Time, class Algebra, class Operations,
57  class Resizer> class ControlledStepper>
59 {
60  using Algebra = boost::numeric::odeint::vector_space_algebra;
61  using Operations = typename boost::numeric::odeint::operations_dispatcher<Eigen::VectorX<FP>>::operations_type;
62  using Resizer = boost::numeric::odeint::initially_resizer;
63  using ErrorChecker = boost::numeric::odeint::default_error_checker<FP, Algebra, Operations>;
65  // Note: use a reference_wrapper so we can both update dt_max, and replace the stepper to change tolerances
66  using Stepper = boost::numeric::odeint::controlled_runge_kutta<
67  ControlledStepper<Eigen::VectorX<FP>, FP, Eigen::VectorX<FP>, FP, Algebra, Operations, Resizer>, ErrorChecker,
68  boost::reference_wrapper<StepAdjuster>>;
69  static constexpr bool is_fsal_stepper = std::is_same_v<typename Stepper::stepper_type::stepper_category,
70  boost::numeric::odeint::explicit_error_stepper_fsal_tag>;
71  static_assert(!is_fsal_stepper,
72  "FSAL steppers cannot be used until https://github.com/boostorg/odeint/issues/72 is resolved.");
73 
74 public:
82  ControlledStepperWrapper(FP abs_tol = 1e-10, FP rel_tol = 1e-5, FP dt_min = std::numeric_limits<FP>::min(),
83  FP dt_max = std::numeric_limits<FP>::max())
84  : OdeIntegratorCore<FP>(dt_min, dt_max)
85  , m_abs_tol(abs_tol)
86  , m_rel_tol(rel_tol)
87  , m_step_adjuster(StepAdjuster{this->get_dt_max()})
88  , m_stepper(create_stepper())
89  {
90  }
91 
92  // if needed, make sure to create a new m_stepper
97 
98  std::unique_ptr<OdeIntegratorCore<FP>> clone() const override
99  {
100  return std::make_unique<ControlledStepperWrapper>(m_abs_tol, m_rel_tol, this->get_dt_min(), this->get_dt_max());
101  }
102 
111  bool step(const mio::DerivFunction<FP>& f, Eigen::Ref<const Eigen::VectorX<FP>> yt, FP& t, FP& dt,
112  Eigen::Ref<Eigen::VectorX<FP>> ytp1) const override
113  {
114  using boost::numeric::odeint::fail;
115  using std::max;
116  assert(0 <= this->get_dt_min());
117  assert(this->get_dt_min() <= this->get_dt_max());
118  // synchronise dt_max of the base class with the stepper
119  m_step_adjuster.set_dt_max(this->get_dt_max());
120  // warn about (upcoming) restrictions to dt
121  if (dt < this->get_dt_min() || dt > this->get_dt_max()) {
122  mio::log_warning("IntegratorCore: Restricting given step size dt = {} to [{}, {}].", dt, this->get_dt_min(),
123  this->get_dt_max());
124  }
125  // set initial values for exit conditions
126  auto step_result = fail;
127  bool is_dt_valid = true;
128  // make a integration step, adapting dt to a possibly larger value on success,
129  // or a strictly smaller value on fail.
130  // stop only on a successful step or a failed step size adaption (w.r.t. the minimal step size dt_min)
131  while (step_result == fail && is_dt_valid) {
132  if (dt < this->get_dt_min()) {
133  is_dt_valid = false;
134  dt = this->get_dt_min();
135  }
136  // we use the scheme try_step(sys, in, t, out, dt) with sys=f, in=y(t), out=y(t').
137  // this is similiar to do_step, but it can adapt the step size dt. If successful, it also updates t.
138  // Note: the resizer used by m_stepper restricts dt to dt_max (via making a failed step)
139  if constexpr (!is_fsal_stepper) { // prevent compile time errors with fsal steppers
140  step_result = m_stepper.try_step(
141  // reorder arguments of the DerivFunction f for the stepper
142  [&](const Eigen::VectorX<FP>& x, Eigen::VectorX<FP>& dxds, FP s) {
143  f(x, s, dxds);
144  },
145  yt, t, ytp1, dt);
146  }
147  }
148  // bound dt from below
149  // the last adaptive step (successful or not) may have calculated a new step size smaller than m_dt_min
150 
151  dt = max<FP>(dt, this->get_dt_min());
152  // check whether the last step failed (which means that m_dt_min was still too large to suffice tolerances)
153  if (step_result == fail) {
154  // adaptive stepping failed, but we still return the result of the last attempt
155  t += this->get_dt_min();
156  return false;
157  }
158  else {
159  // successfully made an integration step
160  return true;
161  }
162  }
163 
165  void set_abs_tolerance(FP abs_tol)
166  {
167  m_abs_tol = abs_tol;
168  m_stepper = create_stepper();
169  }
170 
172  void set_rel_tolerance(FP rel_tol)
173  {
174  m_rel_tol = rel_tol;
175  m_stepper = create_stepper();
176  }
177 
179  void set_dt_min(FP dt_min)
180  {
181  this->get_dt_min() = dt_min;
182  }
183 
185  void set_dt_max(FP dt_max)
186  {
187  this->get_dt_max() = dt_max;
188  }
189 
190 private:
193  {
194  // for more options see: boost/boost/numeric/odeint/stepper/controlled_runge_kutta.hpp
195  return Stepper(ErrorChecker(m_abs_tol, m_rel_tol), boost::ref(m_step_adjuster));
196  }
197 
198  FP m_abs_tol, m_rel_tol;
200  mutable Stepper m_stepper;
201 };
202 
207 template <typename FP, template <class State, class Value, class Deriv, class Time, class Algebra, class Operations,
208  class Resizer> class ExplicitStepper>
210 {
211 public:
212  using Stepper =
213  ExplicitStepper<Eigen::VectorX<FP>, FP, Eigen::VectorX<FP>, FP, boost::numeric::odeint::vector_space_algebra,
214  typename boost::numeric::odeint::operations_dispatcher<Eigen::VectorX<FP>>::operations_type,
215  boost::numeric::odeint::initially_resizer>;
216 
221  : mio::OdeIntegratorCore<FP>(FP{}, FP{})
222  {
223  }
224 
225  std::unique_ptr<OdeIntegratorCore<FP>> clone() const override
226  {
227  return std::make_unique<ExplicitStepperWrapper>(*this);
228  }
229 
238  bool step(const mio::DerivFunction<FP>& f, Eigen::Ref<const Eigen::VectorX<FP>> yt, FP& t, FP& dt,
239  Eigen::Ref<Eigen::VectorX<FP>> ytp1) const override
240  {
241  // copy the values from y(t) to ytp1, since we use the scheme do_step(sys, inout, t, dt) with
242  // sys=f, inout=y(t) for in-place computation - also, this form is shared by several steppers in boost
243  ytp1 = yt;
244  m_stepper.do_step(
245  // reorder arguments of the DerivFunction f for the stepper
246  [&](const Eigen::VectorX<FP>& x, Eigen::VectorX<FP>& dxds, FP s) {
247  f(x, s, dxds);
248  },
249  ytp1, t, dt);
250  // update time (it is not modified by do_step)
251  t += dt;
252  return true; // no step size adaption
253  }
254 
255 private:
256  mutable Stepper m_stepper;
257 };
258 
259 } // namespace mio
260 
261 #endif // MIO_MATH_STEPPER_WRAPPER_H
This is an adaptive IntegratorCore.
Definition: stepper_wrapper.h:59
void set_dt_max(FP dt_max)
Definition: stepper_wrapper.h:185
Stepper m_stepper
A stepper instance used for integration.
Definition: stepper_wrapper.h:200
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 and adapt the step size dt.
Definition: stepper_wrapper.h:111
boost::numeric::odeint::default_error_checker< FP, Algebra, Operations > ErrorChecker
Definition: stepper_wrapper.h:63
FP m_abs_tol
Definition: stepper_wrapper.h:198
typename boost::numeric::odeint::operations_dispatcher< Eigen::VectorX< FP > >::operations_type Operations
Definition: stepper_wrapper.h:61
boost::numeric::odeint::controlled_runge_kutta< ControlledStepper< Eigen::VectorX< FP >, FP, Eigen::VectorX< FP >, FP, Algebra, Operations, Resizer >, ErrorChecker, boost::reference_wrapper< StepAdjuster > > Stepper
Definition: stepper_wrapper.h:68
ControlledStepperWrapper & operator=(ControlledStepperWrapper &&other)=delete
boost::numeric::odeint::initially_resizer Resizer
Definition: stepper_wrapper.h:62
void set_abs_tolerance(FP abs_tol)
Definition: stepper_wrapper.h:165
ControlledStepperWrapper & operator=(const ControlledStepperWrapper &other)=delete
Stepper create_stepper()
(Re)initialize the internal stepper.
Definition: stepper_wrapper.h:192
ControlledStepperWrapper(ControlledStepperWrapper &&other)=delete
void set_dt_min(FP dt_min)
Definition: stepper_wrapper.h:179
void set_rel_tolerance(FP rel_tol)
Definition: stepper_wrapper.h:172
std::unique_ptr< OdeIntegratorCore< FP > > clone() const override
Definition: stepper_wrapper.h:98
ControlledStepperWrapper(const ControlledStepperWrapper &other)=delete
boost::numeric::odeint::vector_space_algebra Algebra
Definition: stepper_wrapper.h:60
StepAdjuster m_step_adjuster
Defines step sizing. Holds a copy of dt_max that has to be updated.
Definition: stepper_wrapper.h:199
ControlledStepperWrapper(FP abs_tol=1e-10, FP rel_tol=1e-5, FP dt_min=std::numeric_limits< FP >::min(), FP dt_max=std::numeric_limits< FP >::max())
Set up the integrator.
Definition: stepper_wrapper.h:82
This is a non-adaptive IntegratorCore.
Definition: stepper_wrapper.h:210
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.
Definition: stepper_wrapper.h:238
Stepper m_stepper
A stepper instance used for integration.
Definition: stepper_wrapper.h:256
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 > Stepper
Definition: stepper_wrapper.h:215
std::unique_ptr< OdeIntegratorCore< FP > > clone() const override
Definition: stepper_wrapper.h:225
ExplicitStepperWrapper()
Set up the integrator.
Definition: stepper_wrapper.h:220
Interface class defining the integration step used in a SystemIntegrator.
Definition: integrator.h:48
static min_max_return_type< ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 >, ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > >::type min(const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &a, const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &b)
Definition: ad.hpp:2599
static min_max_return_type< ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 >, ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > >::type max(const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &a, const ad::internal::active_type< AD_TAPE_REAL, DATA_HANDLER_1 > &b)
Definition: ad.hpp:2596
A collection of classes to simplify handling of matrix shapes in meta programming.
Definition: models/abm/analyze_result.h:30
void log_warning(spdlog::string_view_t fmt, const Args &... args)
Definition: logging.h:112
std::function< void(Eigen::Ref< const Eigen::VectorX< FP > > y, FP t, Eigen::Ref< Eigen::VectorX< FP > > dydt)> DerivFunction
Function template to be integrated.
Definition: integrator.h:39
Extends the default_step_adjuster with a setter for dt_max.
Definition: stepper_wrapper.h:42
void set_dt_max(const Time &dt_max)
Definition: stepper_wrapper.h:44