Drake
integrator_base.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include <algorithm>
4 #include <limits>
5 #include <memory>
6 #include <utility>
7 
14 
15 namespace drake {
16 namespace systems {
17 
18 /*
19 TODO(edrumwri): comments below in preparation of upcoming functionality
20  Note: we ensure algorithmically that no report time, scheduled time, or
21 final time t can occur *within* an event window, that is, we will never have
22  t_low < t < t_high for any interesting t. Further, t_report, t_scheduled and
23  t_final can coincide with t_high but only t_report can be at t_low. The
24  interior of t_low:t_high is a "no man's land" where we don't understand the
25 solution, so must be avoided.
26 */
27 
28 /**
29  * An abstract class for an integrator for ODEs and DAEs as represented by a
30  * Drake System. Integrators solve initial value problems of the form:<pre>
31  * ẋ(t) = f(t, x(t)) with f : ℝ × ℝⁿ → ℝⁿ
32  * </pre>
33  * (i.e., `f()` is an ordinary differential equation) given initial conditions
34  * (t₀, x₀). Thus, integrators advance the continuous state of a dynamical
35  * system forward in time.
36  *
37  * Apart from solving initial value problems, for which the integrator is a
38  * key component of a simulator, integrators can also be used to solve
39  * boundary value problems (via numerical methods like the Multiple Shooting
40  * Method) and trajectory optimization problems (via numerical methods like
41  * direct transcription). This class and its derivatives were developed
42  * primarily toward the former application (through IntegrateAtMost() and
43  * the Simulator class). However, the IntegratorBase architecture was developed
44  * to support these ancillary applications as well using the
45  * IntegrateWithMultipleSteps() and IntegrateWithSingleFixedStep() methods;
46  * the latter permits the caller to advance time using fixed steps in
47  * applications where variable stepping would be deleterious (e.g., direct
48  * transcription).
49  *
50  * A natural question for a user to ask of an integrator is: Which scheme
51  * (method) should be applied to a particular problem? The answer is whichever
52  * one most quickly computes the solution to the desired accuracy! Selecting
53  * an integration scheme for a particular problem is presently an artform. As
54  * examples of some selection criteria: multistep methods generally work poorly
55  * when events (that require state reinitializations) are common, symplectic
56  * methods generally work well at maintaining stability for large integration
57  * steps, and stiff integrators are often best for computationally stiff
58  * systems. If ignorant as to the characteristics of a particular problem, it
59  * is often best to start with an explicit, Runge-Kutta type method. Statistics
60  * collected by the integrator can help diagnose performance issues and possibly
61  * point to the use of a different integration scheme.
62  *
63  * Some systems are known to exhibit "computational stiffness", by which it is
64  * meant that (excessively) small integration steps are necessary for purposes
65  * of stability: in other words, steps must be taken smaller than that
66  * required to achieve a desired accuracy *over a particular interval*.
67  * Thus, the nature of computationally stiff problems is that the solution to
68  * the ODE is *smooth* in the interval of stiffness (in contrast, some problems
69  * possess such high frequency dynamics that very small steps are simply
70  * necessary to capture the solution accurately). Implicit integrators are the
71  * go-to approach for solving computationally stiff problems, but careful
72  * consideration is warranted. Implicit integrators typically require much more
73  * computation than non-implicit (explicit) integrators, stiffness might be an
74  * issue on only a very small time interval, and some problems might be only
75  * "moderately stiff". Put another way, applying an implicit integrator to a
76  * potentially stiff problem might not yield faster computation. The first
77  * chapter of [Hairer, 1996] illustrates the issues broached in this paragraph
78  * using various examples.
79  *
80  * Established methods for integrating ordinary differential equations
81  * invariably make provisions for estimating the "local error" (i.e., the
82  * error over a small time interval) of a solution. Although the relationship
83  * between local error and global error (i.e., the accumulated error over
84  * multiple time steps) can be tenuous, such error estimates can allow
85  * integrators to work adaptively, subdividing time intervals as necessary
86  * (if, e.g., the system is particularly dynamic or stationary in an interval).
87  * Even for applications that do not recommend such adaptive integration- like
88  * direct transcription methods for trajectory optimization- error estimation
89  * allows the user to assess the accuracy of the solution.
90  *
91  * IntegratorBase provides numerous settings and flags that can leverage
92  * problem-specific information to speed integration and/or improve integration
93  * accuracy. As an example, set_maximum_step_size() allows the user to prevent
94  * overly large integration steps (that integration error control alone might
95  * be insufficient to detect). As noted previously, IntegratorBase also collects
96  * a plethora of statistics that can be used to diagnose poor integration
97  * performance. For example, a large number of shrinkages due to error control
98  * could indicate that a system is computationally stiff.
99  *
100  * - [Hairer, 1996] E. Hairer and G. Wanner. Solving Ordinary Differential
101  * Equations II (Stiff and Differential-Algebraic Problems).
102  * Springer, 1996.
103 
104  * @tparam T The vector element type, which must be a valid Eigen scalar.
105  */
106 template <class T>
108  public:
110 
111  /**
112  * Status returned by StepOnceAtMost().
113  * When a step is successful, it will return an indication of what caused it
114  * to stop where it did. When unsuccessful it will throw an exception so you
115  * won't see any return value. When return of control is due ONLY to reaching
116  * a publish time, (status is kReachedPublishTime) the context may return an
117  * interpolated value at an earlier time.
118  *
119  * Note: the simulation step must always end at an update time but can end
120  * after a publish time.
121  */
122  // TODO(edrumwri): incorporate kReachedZeroCrossing into the simulator.
123  enum StepResult {
124  /**
125  * Indicates a publish time has been reached but not an update time.
126  */
128  /**
129  * Localized an event; this is the *before* state (interpolated).
130  */
132  /**
133  * Indicates that integration terminated at an update time.
134  */
136  /**
137  * User requested control whenever an internal step is successful.
138  */
140  /**
141  * Reached the desired integration time without reaching an update time.
142  */
144  /** Took maximum number of steps without finishing integrating over the
145  * interval. */
147  };
148 
149  /**
150  * Maintains references to the system being integrated and the context used
151  * to specify the initial conditions for that system (if any).
152  * @param system A reference to the system to be integrated; the integrator
153  * will maintain a reference to the system in perpetuity, so
154  * the integrator must not outlive the system.
155  * @param context A pointer to a writeable context (nullptr is ok, but a
156  * non-null pointer must be set before Initialize() is
157  * called). The integrator will advance the system state using
158  * the pointer to this context. The pointer to the context will
159  * be maintained internally. The integrator must not outlive
160  * the context.
161  */
162  explicit IntegratorBase(const System<T>& system,
163  Context<T>* context = nullptr)
164  : system_(system), context_(context) {
165  initialization_done_ = false;
166  }
167 
168  /** Destructor. */
169  virtual ~IntegratorBase() = default;
170 
171  /**
172  * Indicates whether an integrator supports error estimation.
173  * Without error estimation, target accuracy will be unused.
174  */
175  virtual bool supports_error_estimation() const = 0;
176 
177  /**
178  * Sets an integrator with error control to fixed step mode. If the integrator
179  * runs in fixed step mode, it will always take the maximum step size
180  * directed (which may be that determined by get_maximum_step_size(), or may
181  * be smaller, as directed by, e.g., Simulator for event handling purposes).
182  * @throws std::logic_error if integrator does not support error
183  * estimation and @p flag is set to `false`.
184  */
185  void set_fixed_step_mode(bool flag) {
186  if (!flag && !supports_error_estimation())
187  throw std::logic_error("Integrator does not support accuracy estimation");
188  fixed_step_mode_ = flag;
189  }
190 
191  /**
192  * Gets whether an integrator is running in fixed step mode. If the integrator
193  * does not support error estimation, this function will always return `true`.
194  * If the integrator runs in fixed step mode, it will always take the maximum
195  * step size directed (which may be that determined by get_maximum_step_size()
196  * or may be smaller, as directed by, e.g., Simulator for event handling
197  * purposes).
198  * @sa set_fixed_step_mode()
199  */
200  bool get_fixed_step_mode() const {
201  return (!supports_error_estimation() || fixed_step_mode_);
202  }
203 
204  /** Request that the integrator attempt to achieve a particular accuracy for
205  * the continuous portions of the simulation. Otherwise a default accuracy is
206  * chosen for you. This may be ignored for fixed-step integration since
207  * accuracy control requires variable step sizes. You should call
208  * supports_error_estimation() to ensure that the
209  * integrator supports this capability before calling this function; if
210  * the integrator does not support it, this method will throw an exception.
211  *
212  * Integrators vary in the range of accuracy (loosest to tightest) that they
213  * can support. If you request accuracy outside the supported range for the
214  * chosen integrator it will be quietly adjusted to be in range. You can find
215  * out the accuracy setting actually being used using `get_accuracy_in_use()`.
216  *
217  * The precise meaning of *accuracy* is a complicated discussion, but
218  * translates roughly to the number of significant digits you want in the
219  * results. By convention it is supplied as `10^-digits`, meaning that an
220  * accuracy of 1e-3 provides about three significant digits. For more
221  * information, see [Sherman 2011].
222  * - M. Sherman, et al. Procedia IUTAM 2:241-261 (2011), Section 3.3.
223  * http://dx.doi.org/10.1016/j.piutam.2011.04.023
224  * @throws std::logic_error if integrator does not support error
225  * estimation.
226  */
227  // TODO(edrumwri): complain if integrator with error estimation wants to drop
228  // below the minimum step size
229  void set_target_accuracy(double accuracy) {
231  throw std::logic_error(
232  "Integrator does not support accuracy estimation "
233  "and user has requested error control");
234  target_accuracy_ = accuracy;
235  accuracy_in_use_ = accuracy;
236  }
237 
238  /**
239  * Gets the target accuracy.
240  * @sa get_accuracy_in_use()
241  */
242  double get_target_accuracy() const { return target_accuracy_; }
243 
244  /**
245  * Gets the accuracy in use by the integrator. This number may differ from
246  * the target accuracy if, for example, the user has requested an accuracy
247  * not attainable or not recommended for the particular integrator.
248  */
249  double get_accuracy_in_use() const { return accuracy_in_use_; }
250 
251  /**
252  * Sets the maximum step size that may be taken by this integrator. The
253  * integrator may stretch the maximum step size by as much as 1% to reach a
254  * discrete event. For fixed step integrators, all steps will be taken at the
255  * maximum step size *unless* an event would be missed.
256  */
257  // TODO(edrumwri): Update this comment when stretch size is configurable.
258  void set_maximum_step_size(const T& max_step_size) {
259  DRAKE_ASSERT(max_step_size >= 0.0);
260  max_step_size_ = max_step_size;
261  }
262 
263  /**
264  * Gets the maximum step size that may be taken by this integrator. This is
265  * a soft maximum: the integrator may stretch it by as much as 1% to hit a
266  * discrete event.
267  * @sa set_requested_minimum_step_size()
268  */
269  // TODO(edrumwri): Update this comment when stretch size is configurable.
270  const T& get_maximum_step_size() const { return max_step_size_; }
271 
272  /**
273  * @anchor Minstep
274  * @name Methods for minimum integration step size selection and behavior
275  *
276  * Variable step integrators reduce their step sizes as needed to achieve
277  * requirements such as specified accuracy or step convergence. However, it is
278  * not possible to take an arbitrarily small step. Normally integrators choose
279  * an appropriate minimum step and throw an exception if the requirements
280  * can't be achieved without going below that. Methods in this section allow
281  * you to influence two aspects of this procedure:
282  * - you can increase the minimum step size, and
283  * - you can control whether an exception is thrown if a smaller step would
284  * have been needed to achieve the aforementioned integrator requirements.
285  *
286  * By default, integrators allow a very small minimum step which can
287  * result in long run times. Setting a larger minimum can be helpful as a
288  * diagnostic to figure out what aspect of your simulation is requiring small
289  * steps. You can set the minimum to what should be a "reasonable" minimum
290  * based on what you know about the physical system. You will then get an
291  * std::runtime_error exception thrown at any point in time where your model
292  * behaves unexpectedly (due to, e.g., a discontinuity in the derivative
293  * evaluation function).
294  *
295  * If you disable the exception (via
296  * `set_throw_on_minimum_step_size_violation(false)`), the integrator will
297  * simply proceed with a step of the minimum size: accuracy is guaranteed
298  * only when the minimum step size is not violated. Beware that there can be
299  * no guarantee about the magnitude of any errors introduced by violating the
300  * accuracy "requirements" in this manner, so disabling the exception should
301  * be done warily.
302  *
303  * #### Details
304  * Because time is maintained to finite precision, there is an absolute
305  * minimum step size `h_floor` required to avoid roundoff error. The
306  * integrator will never take a step smaller than `h_floor`. We calculate
307  * `h_floor=max(ε,ε⋅t)`, where t is the current time and ε is a small multiple
308  * of machine precision, typically a number like 1e-14. Note that `h_floor`
309  * necessarily grows with time; if that is a concern you should limit how
310  * long your simulations are allowed to run without resetting time.
311  *
312  * You may request a larger minimum step size `h_min`. Then at every time t,
313  * the integrator determines a "working" minimum `h_work=max(h_min,h_floor)`.
314  * If the step size selection algorithm determines that a step smaller than
315  * `h_work` is needed to meet accuracy or other needs, then a
316  * std::runtime_error exception will be thrown and the simulation halted. On
317  * the other hand, if you have suppressed the exception (again, via
318  * `set_throw_on_minimum_step_size_violation(false)`), the integration
319  * will continue, taking a step of size `h_work`.
320  *
321  * Under some circumstances the integrator may legitimately take a step of
322  * size `h` smaller than your specified `h_min`, although never smaller than
323  * `h_floor`. For example, occasionally the integrator may reach an event or
324  * time limit that occurs a very short time after the end of a previous step,
325  * necessitating that a tiny "sliver" of a step be taken to complete the
326  * interval. That does not indicate an error, and required accuracy and
327  * convergence goals are achieved. Larger steps can resume immediately
328  * afterwards. Another circumstance is when one of the integrator's stepping
329  * methods is called directly requesting a very small step, for example
330  * `IntegrateWithMultipleSteps(h)`. No exception will be thrown in either of
331  * these cases.
332  */
333 
334  //@{
335  /**
336  * Sets the requested minimum step size `h_min` that may be taken by this
337  * integrator. No step smaller than this will be taken except under
338  * circumstances as described @link Minstep above. @endlink This setting will
339  * be ignored if it is smaller than the absolute minimum `h_floor` also
340  * described above. Default value is zero.
341  * @param min_step_size a non-negative value. Setting this value to zero
342  * will cause the integrator to use a reasonable value
343  * instead (see get_working_minimum_step_size()).
344  * @sa get_requested_minimum_step_size()
345  * @sa get_working_minimum_step_size()
346  */
347  void set_requested_minimum_step_size(const T& min_step_size) {
348  DRAKE_ASSERT(min_step_size >= 0.0);
349  req_min_step_size_ = min_step_size;
350  }
351 
352  /**
353  * Gets the requested minimum step size `h_min` for this integrator.
354  * @sa set_requested_minimum_step_size()
355  * @sa get_working_minimum_step_size(T)
356  */
358  return req_min_step_size_; }
359 
360  /**
361  * Sets whether the integrator should throw a std::runtime_error exception
362  * when the integrator's step size selection algorithm determines that it
363  * must take a step smaller than the minimum step size (for, e.g., purposes
364  * of error control). Default is `true`. If `false`, the integrator will
365  * advance time and state using the minimum specified step size in such
366  * situations. See @link Minstep this section @endlink for more detail.
367  */
369  min_step_exceeded_throws_ = throws;
370  }
371 
372  /**
373  * Reports the current setting of the throw_on_minimum_step_size_violation
374  * flag.
375  * @sa set_throw_on_minimum_step_size_violation().
376  */
378  return min_step_exceeded_throws_;
379  }
380 
381  /**
382  * Gets the current value of the working minimum step size `h_work(t)` for
383  * this integrator, which may vary with the current time t as stored in the
384  * integrator's context.
385  * See @link Minstep this section @endlink for more detail.
386  */
388  using std::max;
389  // Tolerance is just a number close to machine epsilon.
390  const double tol = 1e-14;
391  const T smart_minimum = max(tol, get_context().get_time()*tol);
392  return max(smart_minimum, req_min_step_size_);
393  }
394  //@}
395 
396  /**
397  * Resets the integrator to initial values, i.e., default construction
398  * values.
399  */
400  void Reset() {
401  // Kill the error estimate and weighting matrices.
402  err_est_.reset();
403  qbar_weight_.setZero(0);
404  z_weight_.setZero(0);
405  pinvN_dq_change_.reset();
406  unweighted_substate_change_.setZero(0);
407  weighted_q_change_.reset();
408 
409  // Integrator no longer operates in fixed step mode.
410  fixed_step_mode_ = false;
411 
412  // Statistics no longer valid.
413  ResetStatistics();
414 
415  // Wipe out settings.
416  req_min_step_size_ = 0;
417  max_step_size_ = nan();
418  accuracy_in_use_ = nan();
419 
420  // Indicate values used for error controlled integration no longer valid.
421  prev_step_size_ = nan();
422  ideal_next_step_size_ = nan();
423 
424  // Call the derived integrator reset routine.
425  DoReset();
426 
427  // Indicate that initialization is necessary.
428  initialization_done_ = false;
429  }
430 
431  /**
432  * An integrator must be initialized before being used. The pointer to the
433  * context must be set before Initialize() is called (or an std::logic_error
434  * will be thrown). If Initialize() is not called, an exception will be
435  * thrown when attempting to call StepOnceAtMost(). To reinitialize the
436  * integrator, Reset() should be called followed by Initialize().
437  * @throws std::logic_error If the context has not been set or a user-set
438  * parameter has been set illogically (i.e., one of the
439  * weighting matrix coefficients is set to a negative value- this
440  * check is only performed for integrators that support error
441  * estimation; the maximum step size is smaller than the minimum
442  * step size; the requested initial step size is outside of the
443  * interval [minimum step size, maximum step size]).
444  * @sa Reset()
445  */
446  void Initialize() {
447  if (!context_) throw std::logic_error("Context has not been set.");
448 
449  // Verify that user settings are reasonable.
450  if (max_step_size_ < req_min_step_size_) {
451  throw std::logic_error("Integrator maximum step size is less than the "
452  "minimum step size");
453  }
454  if (req_initial_step_size_ > max_step_size_) {
455  throw std::logic_error("Requested integrator initial step size is larger "
456  "than the maximum step size.");
457  }
458  if (req_initial_step_size_ < req_min_step_size_) {
459  throw std::logic_error("Requested integrator initial step size is smaller"
460  " than the minimum step size.");
461  }
462 
463  // TODO(edrumwri): Compute qbar_weight_, z_weight_ automatically.
464  // Set error weighting vectors if not already done.
466  // Allocate space for the error estimate.
467  err_est_ = system_.AllocateTimeDerivatives();
468 
469  const auto& xc = context_->get_state().get_continuous_state();
470  const int gv_size = xc.get_generalized_velocity().size();
471  const int misc_size = xc.get_misc_continuous_state().size();
472  if (qbar_weight_.size() != gv_size) qbar_weight_.setOnes(gv_size);
473  if (z_weight_.size() != misc_size) z_weight_.setOnes(misc_size);
474 
475  // Verify that minimum values of the weighting matrices are non-negative.
476  if ((qbar_weight_.size() && qbar_weight_.minCoeff() < 0) ||
477  (z_weight_.size() && z_weight_.minCoeff() < 0))
478  throw std::logic_error("Scaling coefficient is less than zero.");
479  }
480 
481  // Statistics no longer valid.
482  ResetStatistics();
483 
484  // Call the derived integrator initialization routine (if any)
485  DoInitialize();
486 
487  initialization_done_ = true;
488  }
489 
490  /**
491  * Request that the first attempted integration step have a particular size.
492  * If no request is made, the integrator will estimate a suitable size
493  * for the initial step attempt. *If the integrator does not support error
494  * control*, this method will throw a std::logic_error (call
495  * supports_error_estimation() to verify before calling this method). For
496  * variable-step integration, the initial target will be treated as a maximum
497  * step size subject to accuracy requirements and event occurrences. You can
498  * find out what size *actually* worked with
499  * `get_actual_initial_step_size_taken()`.
500  * @throws std::logic_error If the integrator does not support error
501  * estimation.
502  */
503  void request_initial_step_size_target(const T& step_size) {
504  using std::isnan;
506  throw std::logic_error(
507  "Integrator does not support error estimation and "
508  "user has initial step size target");
509  req_initial_step_size_ = step_size;
510  }
511 
512  /**
513  * Gets the target size of the first integration step. You can find out what
514  * step size was *actually* used for the first integration step with
515  * `get_actual_initial_step_size_taken()`.
516  * @see request_initial_step_size_target()
517  */
518  const T& get_initial_step_size_target() const {
519  return req_initial_step_size_;
520  }
521 
522  /**
523  * Integrates the system forward in time by a single step with step size
524  * subject to integration error tolerances (assuming that the integrator
525  * supports error estimation). The integrator must already have
526  * been initialized or an exception will be thrown. The context will be
527  * integrated forward by an amount that will never exceed the minimum of
528  * `publish_dt`, `update_dt`, and `1.01 * get_maximum_step_size()`.
529  *
530  * @param publish_dt The step size, >= 0.0 (exception will be thrown
531  * if this is not the case) at which the next publish will occur.
532  * @param update_dt The step size, > 0.0 (exception will be thrown
533  * if this is not the case) at which the next update will occur.
534  * @param boundary_dt The step size, >= 0.0 (exception will be thrown
535  * if this is not the case) marking the end of the user-designated
536  * simulated interval.
537  * @throws std::logic_error If the integrator has not been initialized or one
538  * of publish_dt, update_dt, or boundary_dt is
539  * negative.
540  * @return The reason for the integration step ending.
541  * @warning Users should generally not call this function directly; within
542  * simulation circumstances, users will typically call
543  * `Simulator::StepTo()`. In other circumstances, users will
544  * typically call `IntegratorBase::IntegrateWithMultipleSteps()`.
545  *
546  * This method at a glance:
547  * - For integrating ODEs/DAEs via Simulator
548  * - Supports fixed step and variable step integration schemes
549  * - Takes only a single step forward.
550  */
551  // TODO(edrumwri): Make the stretch size configurable.
552  StepResult IntegrateAtMost(const T& publish_dt, const T& update_dt,
553  const T& boundary_dt);
554 
555  /// Gets the stretch factor (> 1), which is multiplied by the maximum
556  /// (typically user-designated) integration step size to obtain the amount
557  /// that the integrator is able to stretch the maximum time step toward
558  /// hitting an upcoming publish or update event in IntegrateAtMost().
559  /// @sa IntegrateAtMost()
560  double get_stretch_factor() const { return 1.01; }
561 
562  /// Stepping function for integrators operating outside of Simulator that
563  /// advances the continuous state exactly by @p dt. This method is designed
564  /// for integrator users that do not wish to consider publishing or
565  /// discontinuous, mid-interval updates. This method will step the integrator
566  /// multiple times, as necessary, to attain requested error tolerances and
567  /// to ensure the integrator converges.
568  /// @warning Users should simulate systems using `Simulator::StepTo()` in
569  /// place of this function (which was created for off-simulation
570  /// purposes), generally.
571  /// @param dt The non-negative integration step to take.
572  /// @throws std::logic_error If the integrator has not been initialized or
573  /// dt is negative.
574  /// @sa IntegrateAtMost(), which is designed to be operated by Simulator and
575  /// accounts for publishing and state reinitialization.
576  /// @sa IntegrateWithSingleStep(), which is also designed to be operated
577  /// *outside of* Simulator, but throws an exception if the integrator
578  /// cannot advance time by @p dt in a single step.
579  ///
580  /// This method at a glance:
581  /// - For integrating ODEs/DAEs not using Simulator
582  /// - Supports fixed step and variable step integration schemes
583  /// - Takes as many steps as necessary until time has advanced by @p dt
584  void IntegrateWithMultipleSteps(const T& dt) {
585  using std::max;
586  using std::min;
587 
588  const Context<T>& context = get_context();
589  const T inf = std::numeric_limits<double>::infinity();
590  T t_remaining = dt;
591 
592  // Note: A concern below is that the while loop while run forever because
593  // t_remaining could be small, but not quite zero, if dt is relatively
594  // small compared to the context time. In such a case, t_final will be
595  // equal to context.get_time() in the expression immediately below,
596  // context.get_time() will not change during the call to IntegrateAtMost(),
597  // and t_remaining will be equal to zero (meaning that the loop will
598  // indeed terminate, as desired).
599  const T t_final = context.get_time() + t_remaining;
600  do {
601  IntegrateAtMost(inf, inf, min(t_remaining, get_maximum_step_size()));
602  t_remaining = t_final - context.get_time();
603  } while (t_remaining > 0);
604  }
605 
606  /// Stepping function for integrators operating outside of Simulator that
607  /// advances the continuous state exactly by @p dt *and using a single fixed
608  /// step*. This method is designed for integrator users that do not wish to
609  /// consider publishing or discontinuous, mid-interval updates. One such
610  /// example application is that of direct transcription for trajectory
611  /// optimization, for which the integration process should be _consistent_: it
612  /// should execute the same sequence of arithmetic operations for all values
613  /// of the nonlinear programming variables. In keeping with the naming
614  /// semantics of this function, error controlled integration is not supported
615  /// (though error estimates will be computed for integrators that support that
616  /// feature), which is a minimal requirement for "consistency".
617  /// @warning Users should simulate systems using `Simulator::StepTo()` in
618  /// place of this function (which was created for off-simulation
619  /// purposes), generally.
620  /// @param dt The non-negative integration step to take.
621  /// @throws std::logic_error If the integrator has not been initialized or
622  /// dt is negative **or** if the integrator
623  /// is not operating in fixed step mode.
624  /// @throws std::runtime_error If the integrator was unable to take a step
625  /// of the requested size.
626  /// @sa IntegrateAtMost(), which is designed to be operated by Simulator and
627  /// accounts for publishing and state reinitialization.
628  /// @sa IntegrateWithMultipleSteps(), which is also designed to be operated
629  /// *outside of* Simulator, but will take as many integration steps as
630  /// necessary until time has been stepped forward by @p dt.
631  ///
632  /// This method at a glance:
633  /// - For integrating ODEs/DAEs not using Simulator
634  /// - Fixed step integration (no step size reductions for error control or
635  /// integrator convergence)
636  /// - Takes only a single step forward.
637  void IntegrateWithSingleFixedStep(const T& dt) {
638  if (dt < 0) {
639  throw std::logic_error("IntegrateWithSingleFixedStep() called with a "
640  "negative step size.");
641  }
642  if (!this->get_fixed_step_mode())
643  throw std::logic_error("IntegrateWithSingleFixedStep() requires fixed "
644  "stepping.");
645  if (!Step(dt)) {
646  throw std::runtime_error("Integrator was unable to take a single fixed "
647  "step of the requested size.");
648  }
649 
650  UpdateStepStatistics(dt);
651  }
652 
653  /**
654  * @name Integrator statistics methods.
655  *
656  * @{
657  * These methods allow the caller to manipulate and query integrator
658  * statistics. Generally speaking, the larger the integration step taken, the
659  * faster a simulation will run. These methods allow querying (and resetting)
660  * the integrator statistics as one means of determining how to make
661  * a simulation run faster.
662  */
663  /**
664  * Forget accumulated statistics. These are reset to the values they have
665  * post construction or immediately after `Initialize()`.
666  */
668  actual_initial_step_size_taken_ = nan();
669  smallest_adapted_step_size_taken_ = nan();
670  largest_step_size_taken_ = nan();
671  num_steps_taken_ = 0;
672  num_ode_evals_ = 0;
673  num_shrinkages_from_error_control_ = 0;
674  num_shrinkages_from_substep_failures_ = 0;
675  num_substep_failures_ = 0;
677  }
678 
679  /// Gets the number of failed sub-steps (implying one or more step size
680  /// reductions was required to permit solving the necessary nonlinear system
681  /// of equations).
682  int64_t get_num_substep_failures() const {
683  return num_substep_failures_;
684  }
685 
686  /// Gets the number of step size shrinkages due to sub-step failures (e.g.,
687  /// integrator convergence failures) since the last call to ResetStatistics()
688  /// or Initialize().
690  return num_shrinkages_from_substep_failures_;
691  }
692 
693  /// Gets the number of step size shrinkages due to failure to meet targeted
694  /// error tolerances, since the last call to ResetStatistics or Initialize().
696  return num_shrinkages_from_error_control_;
697  }
698 
699 /**
700  * Returns the number of ODE function evaluations (calls to
701  * CalcTimeDerivatives()) since the last call to ResetStatistics() or
702  * Initialize(). This count includes *all* such calls including (1)
703  * those necessary to compute Jacobian matrices; (2) those used in rejected
704  * integrated steps (for, e.g., purposes of error control); (3) those used
705  * strictly for integrator error estimation; and (4) calls that exhibit little
706  * cost (due to results being cached).
707  */
708  int64_t get_num_derivative_evaluations() const { return num_ode_evals_; }
709 
710  /**
711  * The actual size of the successful first step.
712  */
714  return actual_initial_step_size_taken_;
715  }
716 
717  /**
718  * The size of the smallest step taken *as the result of a controlled
719  * integration step adjustment* since the last Initialize() or
720  * ResetStatistics() call. This value will be NaN for integrators without
721  * error estimation.
722  */
724  return smallest_adapted_step_size_taken_;
725  }
726 
727  /**
728  * The size of the largest step taken since the last Initialize() or
729  * ResetStatistics() call.
730  */
731  const T& get_largest_step_size_taken() const {
732  return largest_step_size_taken_;
733  }
734 
735  /**
736  * The number of integration steps taken since the last Initialize()
737  * or ResetStatistics() call.
738  */
739  int64_t get_num_steps_taken() const { return num_steps_taken_; }
740 
741  /**
742  * @}
743  */
744 
745  /**
746  * Return the step size the integrator would like to take next, based
747  * primarily on the integrator's accuracy prediction. This value will not
748  * be computed for integrators that do not support error estimation and
749  * NaN will be returned.
750  */
751  const T& get_ideal_next_step_size() const { return ideal_next_step_size_; }
752 
753  /**
754  * Returns a const reference to the internally-maintained Context holding
755  * the most recent state in the trajectory. This is suitable for publishing or
756  * extracting information about this trajectory step.
757  */
758  const Context<T>& get_context() const { return *context_; }
759 
760  /**
761  * Returns a mutable pointer to the internally-maintained Context holding
762  * the most recent state in the trajectory.
763  */
764  Context<T>* get_mutable_context() { return context_; }
765 
766  /**
767  * Replace the pointer to the internally-maintained Context with a different
768  * one. This is useful for supplying a new set of initial conditions or
769  * wiping out the current context (by passing in a null pointer). You
770  * should invoke Initialize() after replacing the Context unless the
771  * context is null.
772  * @param context The pointer to the new context or nullptr to wipe out
773  * the current context without replacing it with another.
774  */
775  void reset_context(Context<T>* context) {
776  context_ = context;
777  initialization_done_ = false;
778  }
779 
780  /**
781  * Gets a constant reference to the system that is being integrated (and
782  * was provided to the constructor of the integrator).
783  */
784  const System<T>& get_system() const { return system_; }
785 
786  /// Indicates whether the integrator has been initialized.
787  bool is_initialized() const { return initialization_done_; }
788 
789  /**
790  * Derived classes must override this function to return the order of
791  * the integrator's error estimate. If the integrator does not provide an
792  * error estimate, the derived class implementation should return 0.
793  */
794  virtual int get_error_estimate_order() const = 0;
795 
796  /**
797  * Gets the size of the last (previous) integration step. If no integration
798  * steps have been taken, value will be NaN.
799  */
801  return prev_step_size_;
802  }
803 
804  /**
805  * Gets the error estimate (used only for integrators that support error
806  * estimation). If the integrator does not support error estimation, nullptr
807  * is returned.
808  */
810  return err_est_.get();
811  }
812 
813  /**
814  * @name Methods for weighting state variable errors
815  * @{
816  * This group of methods describes how errors for state variables with
817  * heterogeneous units are weighted in the context of error-controlled
818  * integration. This is an advanced topic and most users can simply specify
819  * desired accuracy and accept the default state variable weights.
820  *
821  * A collection of state variables is generally defined in heterogenous units
822  * (e.g. length, angles, velocities, energy). Some of the state
823  * variables cannot even be expressed in meaningful units, like
824  * quaternions. Certain integrators provide an estimate of the absolute error
825  * made in each state variable during an integration step. These errors must
826  * be properly weighted to obtain an "accuracy" _with respect to each
827  * particular variable_. These per-variable accuracy determinations can be
828  * compared against the user's requirements and used to select an appropriate
829  * size for the next step [Sherman 2011]. The weights are
830  * normally determined automatically using the system's characteristic
831  * dimensions, so *most users can stop reading now!* Custom weighting is
832  * primarily useful for performance improvement; an optimal weighting would
833  * allow an error-controlled integrator to provide the desired level of
834  * accuracy across all state variables without wasting computation
835  * achieving superfluous accuracy for some of those variables.
836  *
837  * Users interested in more precise control over state variable weighting may
838  * use the methods in this group to access and modify weighting factors for
839  * individual state variables. Changes to these weights can only be made prior
840  * to integrator initialization or as a result of an event being triggered
841  * and then followed by re-initialization.
842  *
843  * <h4>Relative versus absolute accuracy</h4>
844  *
845  * %State variable integration error, as estimated by an integrator, is an
846  * absolute quantity with the same
847  * units as the variable. At each time step we therefore need to determine
848  * an absolute error that would be deemed "good enough", i.e. satisfies
849  * the user's accuracy requirement. If a variable is maintained to a
850  * *relative* accuracy then that "good enough" value is defined to be the
851  * required accuracy `a` (a fraction like 0.001) times the current value of
852  * the variable, as long as that value
853  * is far from zero. For variables maintained to an *absolute* accuracy, or
854  * relative variables that are at or near zero (where relative accuracy would
855  * be undefined or too strict, respectively), we need a different way to
856  * determine the "good enough" absolute error. The methods in this section
857  * control how that absolute error value is calculated.
858  *
859  * <h4>How to choose weights</h4>
860  *
861  * The weight `wᵢ` for a state variable `xᵢ` should be
862  * chosen so that the product `wᵢ * dxᵢ` is unitless, and in particular is 1
863  * when `dxᵢ` represents a "unit effect" of state variable `xᵢ`; that is, the
864  * change in `xᵢ` that produces a unit change in some quantity of interest in
865  * the system being simulated. Why unity (1)? Aside from normalizing the
866  * values, unity "grounds" the weighted error to the user-specified accuracy.
867  * A weighting can be applied individually to each state variable, but
868  * typically it is done approximately by combining the known type of the
869  * variable (e.g. length, angle) with a "characteristic scale" for that
870  * quantity. For example, if a "characteristic length" for the system being
871  * simulated is 0.1 meters, and `x₀` is a length variable measured in meters,
872  * then `w₀` should be 10 so that `w₀*dx₀=1` when `dx₀=0.1`. For angles
873  * representing pointing accuracy (say a camera direction) we typically assume
874  * a "characteristic angle" is one radian (about 60 degrees), so if x₁ is a
875  * pointing direction then w₁=1 is an appropriate weight. We can now scale an
876  * error vector `e=[dx₀ dx₁]` to a unitless fractional error vector
877  * `f=[w₀*dx₀ w₁*dx₁]`. Now to achieve a given accuracy `a`, say `a=.0001`,
878  * we need only check that `|fᵢ|<=a` for each element `i` of `f`. Further,
879  * this gives us a quantitative measure of "worst accuracy" that we can use
880  * to increase or reduce size of the next attempted step, so that we will just
881  * achieve the required accuracy but not much more. We'll be more precise
882  * about this below.
883  *
884  * <h4>Some subtleties for second-order dynamic systems</h4>
885  *
886  * Systems governed by 2nd-order differential equations are typically split
887  * into second order (configuration) variables q, and rate (velocity)
888  * variables v, where the time derivatives qdot of q are linearly related to
889  * v by the kinematic differential equation `qdot = dq/dt = N(q)*v`.
890  * Velocity variables are
891  * chosen to be physically significant, but configuration variables
892  * may be chosen for convenience and do not necessarily have direct physical
893  * interpretation. For examples, quaternions are chosen as a numerically
894  * stable orientation representation. This is problematic for choosing weights
895  * which must be done by physical reasoning
896  * as sketched above. We resolve this by introducing
897  * the notion of "quasi-coordinates" ꝗ (pronounced "qbar") which are defined
898  * by the equation `ꝗdot = dꝗ/dt = v`. Other than time scaling,
899  * quasi-coordinates have the same units as their corresponding velocity
900  * variables. That is, for weighting we need to think
901  * of the configuration coordinates in the same physical space as the velocity
902  * variables; weight those by their physical significance; and then map back
903  * to an instantaneous weighting
904  * on the actual configuration variables q. This mapping is performed
905  * automatically; you need only to be concerned about physical weightings.
906  *
907  * Note that generalized quasi-coordinates `ꝗ` can only be defined locally for
908  * a particular configuration `q`. There is in general no meaningful set of
909  * `n` generalized
910  * coordinates which can be differentiated with respect to time to yield `v`.
911  * For example, the Hairy Ball Theorem implies that it is not possible for
912  * three orientation variables to represent all 3D rotations without
913  * singularities, yet three velocity variables can represent angular velocity
914  * in 3D without singularities.
915  *
916  * To summarize, separate weights can be provided for each of
917  * - `n` generalized quasi-coordinates `ꝗ` (configuration variables in the
918  * velocity variable space), and
919  * - `nz` miscellaneous continuous state variables `z`.
920  *
921  * Weights on the generalized velocity variables `v (= dꝗ/dt)` are derived
922  * directly from the weights on `ꝗ`, weighted by a characteristic time.
923  * Weights on the actual `nq` generalized coordinates can
924  * be calculated efficiently from weights on the quasi-coordinates (details
925  * below).
926  *
927  * <h4>How the weights are used</h4>
928  *
929  * The errors in the `ꝗ` and `z` variables are weighted by the diagonal
930  * elements
931  * of diagonal weighting matrices Wꝗ and Wz, respectively. (The block-diagonal
932  * weighting matrix `Wq` on the original generalized coordinates `q` is
933  * calculated from `N` and `Wꝗ`; see below.) In the absence of
934  * other information, the default for all weighting values is one, so `Wꝗ` and
935  * `Wz` are `n × n` and `nz × nz` identity matrices. The weighting matrix `Wv`
936  * for the velocity variables is just `Wv = τ*Wꝗ` where `τ` is a
937  * "characteristic time" for the system, that is, a quantity in time units
938  * that represents a significant evolution of the trajectory. This serves to
939  * control the accuracy with which velocity is determined relative to
940  * configuration. Note that larger values of `τ` are more conservative since
941  * they increase the velocity weights. Typically we use `τ=1.0` or `0.1`
942  * seconds for human-scale mechanical systems.
943  * <!-- TODO(sherm1): provide more guidance for velocity weighting. -->
944  *
945  * The weighting matrices `Wq`, `Wv`, and `Wz` are used to compute a weighted
946  * infinity norm as follows. Although `Wv` and `Wz` are constant, the actual
947  * weightings may be state dependent for relative-error calculations.
948  * Define block diagonal error weighting matrix `E=diag(Eq,Ev,Ez)` as follows:
949  * <pre>
950  * Eq = Wq
951  * Ev: Ev(i,i) = { min(Wv(i,i), 1/|vᵢ|) if vᵢ is relative
952  * { Wv(i,i) if vᵢ is absolute
953  * Ez: Ez(i,i) = { min(Wz(i,i), 1/|zᵢ|) if zᵢ is relative
954  * { Wz(i,i) if zᵢ is absolute
955  * </pre>
956  * (`Ev` and `Ez` are diagonal.) A `v` or `z` will be maintained to relative
957  * accuracy unless (a) it is "close" to zero (less than 1), or (b) the
958  * variable has been defined as requiring absolute accuracy. Position
959  * variables `q` are always maintained to absolute accuracy (see
960  * [Sherman 2011] for rationale).
961  *
962  * Now given an error estimate vector `e=[eq ev ez]`, the vector `f=E*e`
963  * can be considered to provide a unitless fractional error for each of the
964  * state variables. To achieve a given user-specified accuracy `a`, we require
965  * that norm_inf(`f`) <= `a`. That is, no element of `f` can have absolute
966  * value larger than `a`. We also use `f` to determine an ideal next step
967  * size using an appropriate integrator-specific computation.
968  *
969  * <h4>Determining weights for q</h4>
970  *
971  * The kinematic differential equations `qdot=N(q)*v` employ an `nq × n`
972  * matrix `N`. By construction, this relationship is invertible using `N`'s
973  * left pseudo-inverse `N⁺` so that `v=N⁺ qdot` and `N⁺ N = I` (the identity
974  * matrix); however, `N N⁺ != I`, as `N` has more rows than columns generally.
975  * [Nikravesh 1988] shows how such a matrix `N` can be determined and provides
976  * more information. Given this relationship between `N` and `N⁺`, we can
977  * relate weighted errors in configuration coordinates `q` to weighted errors
978  * in generalized quasi-coordinates `ꝗ`, as the following derivation shows:
979  * <pre>
980  * v = N⁺ qdot Inverse kinematic differential equation
981  * dꝗ/dt = N⁺ dq/dt Use synonyms for v and qdot
982  * dꝗ = N⁺ dq Change time derivatives to differentials
983  * Wꝗ dꝗ = Wꝗ N⁺ dq Pre-multiply both sides by Wꝗ
984  * N Wꝗ dꝗ = N Wꝗ N⁺ dq Pre-multiply both sides by N
985  * N Wꝗ dꝗ = Wq dq Define Wq := N Wꝗ N⁺
986  * N Wꝗ v = Wq qdot Back to time derivatives.
987  * </pre>
988  * The last two equations show that `Wq` as defined above provides the
989  * expected relationship between the weighted `ꝗ` or `v` variables in velocity
990  * space and the weighted `q` or `qdot` (resp.) variables in configuration
991  * space.
992  *
993  * Finally, note that a diagonal entry of one of the weighting matrices can
994  * be set to zero to disable error estimation for that state variable
995  * (i.e., auxiliary variable or configuration/velocity variable pair), but
996  * that setting an entry to a negative value will cause an exception to be
997  * thrown when the integrator is initialized.
998  *
999  * - [Nikravesh 1988] P. Nikravesh. Computer-Aided Analysis of Mechanical
1000  * Systems. Prentice Hall, 1988. Sec. 6.3.
1001  * - [Sherman 2011] M. Sherman, et al. Procedia IUTAM 2:241-261 (2011),
1002  * Section 3.3. http://dx.doi.org/10.1016/j.piutam.2011.04.023
1003  *
1004  * @sa CalcStateChangeNorm()
1005  */
1006  /**
1007  * Gets the weighting vector (equivalent to a diagonal matrix) applied to
1008  * weighting both generalized coordinate and velocity state variable errors,
1009  * as described in the group documentation. Only used for integrators that
1010  * support error estimation.
1011  */
1012  const Eigen::VectorXd& get_generalized_state_weight_vector() const {
1013  return qbar_weight_;
1014  }
1015 
1016  /**
1017  * Gets a mutable weighting vector (equivalent to a diagonal matrix) applied
1018  * to weighting both generalized coordinate and velocity state variable
1019  * errors, as described in the group documentation. Only used for
1020  * integrators that support error estimation. Returns a VectorBlock
1021  * to make the values mutable without permitting changing the size of
1022  * the vector. Requires re-initializing the integrator after calling
1023  * this method; if Initialize() is not called afterward, an exception will be
1024  * thrown when attempting to call StepOnceAtMost(). If the caller sets
1025  * one of the entries to a negative value, an exception will be thrown
1026  * when the integrator is initialized.
1027  */
1028  Eigen::VectorBlock<Eigen::VectorXd>
1030  initialization_done_ = false;
1031  return qbar_weight_.head(qbar_weight_.rows());
1032  }
1033 
1034  /**
1035  * Gets the weighting vector (equivalent to a diagonal matrix) for
1036  * weighting errors in miscellaneous continuous state variables `z`. Only used
1037  * for integrators that support error estimation.
1038  */
1039  const Eigen::VectorXd& get_misc_state_weight_vector() const {
1040  return z_weight_;
1041  }
1042 
1043  /**
1044  * Gets a mutable weighting vector (equivalent to a diagonal matrix) for
1045  * weighting errors in miscellaneous continuous state variables `z`. Only used
1046  * for integrators that support error estimation. Returns a VectorBlock
1047  * to make the values mutable without permitting changing the size of
1048  * the vector. Requires re-initializing the integrator after calling this
1049  * method. If Initialize() is not called afterward, an exception will be
1050  * thrown when attempting to call StepOnceAtMost(). If the caller sets
1051  * one of the entries to a negative value, an exception will be thrown
1052  * when the integrator is initialized.
1053  */
1054  Eigen::VectorBlock<Eigen::VectorXd> get_mutable_misc_state_weight_vector() {
1055  initialization_done_ = false;
1056  return z_weight_.head(z_weight_.rows());
1057  }
1058 
1059  /**
1060  * @}
1061  */
1062 
1063  protected:
1064  /// Resets any statistics particular to a specific integrator. The default
1065  /// implementation of this function does nothing. If your integrator
1066  /// collects its own statistics, you should re-implement this method and
1067  /// reset them there.
1068  virtual void DoResetStatistics() {}
1069 
1070  /// Evaluates the derivative function (and updates call statistics).
1071  /// Subclasses should call this function rather than calling
1072  /// system.CalcTimeDerivatives() directly.
1073  void CalcTimeDerivatives(const Context<T>& context,
1074  ContinuousState<T>* dxdt) {
1075  get_system().CalcTimeDerivatives(context, dxdt);
1076  ++num_ode_evals_;
1077  }
1078 
1079  /// Evaluates the derivative function (and updates call statistics).
1080  /// Subclasses should call this function rather than calling
1081  /// system.CalcTimeDerivatives() directly. This version of this function
1082  /// exists to allow integrators to count AutoDiff'd systems in derivative
1083  /// function evaluations.
1084  template <typename U>
1085  void CalcTimeDerivatives(const System<U>& system,
1086  const Context<U>& context,
1087  ContinuousState<U>* dxdt) {
1088  system.CalcTimeDerivatives(context, dxdt);
1089  ++num_ode_evals_;
1090  }
1091 
1092  /**
1093  * Sets the working ("in use") accuracy for this integrator. The working
1094  * accuracy may not be equivalent to the target accuracy when the latter is
1095  * too loose or tight for an integrator's capabilities.
1096  * @sa get_accuracy_in_use()
1097  * @sa get_target_accuracy()
1098  */
1099  void set_accuracy_in_use(double accuracy) { accuracy_in_use_ = accuracy; }
1100 
1101  /// Generic code for validating (and resetting, if need be) the integrator
1102  /// working accuracy for error controlled integrators. This method is
1103  /// intended to be called from an integrator's DoInitialize() method.
1104  /// @param default_accuracy a reasonable default accuracy setting for this
1105  /// integrator.
1106  /// @param loosest_accuracy the loosest accuracy that this integrator should
1107  /// support.
1108  /// @param max_step_fraction a fraction of the maximum step size to use when
1109  /// setting the integrator accuracy and the user has not specified
1110  /// accuracy directly.
1111  /// @throws std::logic_error if neither the initial step size target nor
1112  /// the maximum step size has been set.
1113  void InitializeAccuracy(double default_accuracy, double loosest_accuracy,
1114  double max_step_fraction) {
1115  using std::isnan;
1116 
1117  // Set an artificial step size target, if not set already.
1118  if (isnan(this->get_initial_step_size_target())) {
1119  // Verify that maximum step size has been set.
1120  if (isnan(this->get_maximum_step_size()))
1121  throw std::logic_error("Neither initial step size target nor maximum "
1122  "step size has been set!");
1123 
1125  this->get_maximum_step_size() * max_step_fraction);
1126  }
1127 
1128  // Sets the working accuracy to a good value.
1129  double working_accuracy = this->get_target_accuracy();
1130 
1131  // If the user asks for accuracy that is looser than the loosest this
1132  // integrator can provide, use the integrator's loosest accuracy setting
1133  // instead.
1134  if (isnan(working_accuracy)) {
1135  working_accuracy = default_accuracy;
1136  } else {
1137  if (working_accuracy > loosest_accuracy)
1138  working_accuracy = loosest_accuracy;
1139  }
1140  this->set_accuracy_in_use(working_accuracy);
1141  }
1142 
1143  /**
1144  * Default code for advancing the continuous state of the system by a single
1145  * step of @p dt_max (or smaller, depending on error control). This particular
1146  * function is designed to be called directly by an error estimating
1147  * integrator's DoStep() method to effect error-controlled integration.
1148  * The integrator can effect error controlled integration without calling this
1149  * method, if the implementer so chooses, but this default method is expected
1150  * to function well in most circumstances.
1151  * @param[in] dt_max The maximum step size to be taken. The integrator may
1152  * take a smaller step than specified to satisfy accuracy
1153  * requirements, to resolve integrator convergence problems, or
1154  * to respect the integrator's maximum step size.
1155  * @throws std::logic_error if integrator does not support error
1156  * estimation.
1157  * @note This function will shrink the integration step as necessary whenever
1158  * the integrator's DoStep() fails to take the requested step
1159  * e.g., due to integrator convergence failure.
1160  * @returns `true` if the full step of size @p dt_max is taken and `false`
1161  * otherwise (i.e., a smaller step than @p dt_max was taken).
1162  */
1163  bool StepOnceErrorControlledAtMost(const T& dt_max);
1164 
1165  /**
1166  * Computes the infinity norm of a change in continuous state. We use the
1167  * infinity norm to capture the idea that, by providing accuracy requirements,
1168  * the user can indirectly specify error tolerances that act to limit the
1169  * largest error in any state vector component.
1170  * @returns the norm (a non-negative value)
1171  */
1172  T CalcStateChangeNorm(const ContinuousState<T>& dx_state) const;
1173 
1174  /**
1175  * Calculates adjusted integrator step sizes toward keeping state variables
1176  * within error bounds on the next integration step. Note that it is not
1177  * guaranteed that the (possibly) reduced step size will keep state variables
1178  * within error bounds; however, the process of (1) taking a trial
1179  * integration step, (2) calculating the error, and (3) adjusting the step
1180  * size can be repeated until convergence.
1181  * @param err
1182  * The norm of the integrator error that was computed using
1183  * @p attempted_step_size.
1184  * @param attempted_step_size
1185  * The step size that was attempted.
1186  * @param[in,out] at_minimum_step_size
1187  * If `true` on entry, the error control mechanism is not allowed to
1188  * shrink the step because the integrator is stepping at the minimum
1189  * step size (note that this condition will only occur if
1190  * `get_throw_on_minimum_step_size_violation() == false`- an exception
1191  * would be thrown otherwise). If `true` on entry and `false` on exit,
1192  * the error control mechanism has managed to increase the step size
1193  * above the working minimum; if `true` on entry and `true` on exit,
1194  * error control would like to shrink the step size but cannot. If
1195  * `false` on entry and `true` on exit, error control shrank the step
1196  * to the working minimum step size.
1197  * @returns a pair of types bool and T; the bool will be set to `true` if
1198  * the integration step was to be considered successful and `false`
1199  * otherwise. The value of the T type will be set to the recommended next
1200  * step size.
1201  */
1202  std::pair<bool, T> CalcAdjustedStepSize(
1203  const T& err,
1204  const T& attempted_step_size,
1205  bool* at_minimum_step_size) const;
1206 
1207  /**
1208  * Derived classes can override this method to perform special
1209  * initialization. This method is called during the Initialize() method. This
1210  * default method does nothing.
1211  */
1212  virtual void DoInitialize() {}
1213 
1214  /**
1215  * Derived classes can override this method to perform routines when
1216  * Reset() is called. This default method does nothing.
1217  */
1218  virtual void DoReset() {}
1219 
1220  /**
1221  * Derived classes must implement this method to (1) integrate the continuous
1222  * portion of this system forward by a single step of size @p dt and
1223  * (2) set the error estimate (via get_mutable_error_estimate()). This
1224  * method is called during the default Step() method.
1225  * @param dt The integration step to take.
1226  * @returns `true` if successful, `false` if the integrator was unable to take
1227  * a single step of size @p dt (due to, e.g., an integrator
1228  * convergence failure).
1229  * @post If the time on entry is denoted `t`, the time and state will be
1230  * advanced to `t+dt` if the method returns `true`; otherwise, the
1231  * time and state should be reset to those at `t`.
1232  * @warning It is expected that DoStep() will return `true` for some, albeit
1233  * possibly very small, positive value of @p dt. The derived
1234  * integrator's stepping algorithm can make this guarantee, for
1235  * example, by switching to an algorithm not subject to convergence
1236  * failures (e.g., explicit Euler) for very small step sizes.
1237  */
1238  virtual bool DoStep(const T& dt) = 0;
1239 
1240  /**
1241  * Gets an error estimate of the state variables recorded by the last call
1242  * to StepOnceFixedSize(). If the integrator does not support error
1243  * estimation, this function will return nullptr.
1244  */
1245  ContinuousState<T>* get_mutable_error_estimate() { return err_est_.get(); }
1246 
1247  // Sets the actual initial step size taken.
1249  actual_initial_step_size_taken_ = dt;
1250  }
1251 
1252  /**
1253  * Sets the size of the smallest-step-taken statistic as the result of a
1254  * controlled integration step adjustment.
1255  */
1257  smallest_adapted_step_size_taken_ = dt;
1258  }
1259 
1260  // Sets the largest-step-size-taken statistic.
1261  void set_largest_step_size_taken(const T& dt) {
1262  largest_step_size_taken_ = dt;
1263  }
1264 
1265  // Sets the "ideal" next step size (typically done via error control).
1266  void set_ideal_next_step_size(const T& dt) { ideal_next_step_size_ = dt; }
1267 
1268  private:
1269  // Validates that a smaller step size does not fall below the working minimum
1270  // and throws an exception if desired.
1271  void ValidateSmallerStepSize(const T& current_step_size,
1272  const T& new_step_size) const {
1273  if (new_step_size < get_working_minimum_step_size() &&
1274  new_step_size < current_step_size && // Verify step adjusted downward.
1275  min_step_exceeded_throws_) {
1276  SPDLOG_DEBUG(drake::log(), "Integrator wants to select too small step "
1277  "size of {}; working minimum is ", new_step_size,
1279  std::ostringstream str;
1280  str << "Error control wants to select step smaller than minimum" <<
1281  " allowed (" << get_working_minimum_step_size() << ")";
1282  throw std::runtime_error(str.str());
1283  }
1284  }
1285 
1286  // Updates the integrator statistics, accounting for a step just taken of
1287  // size dt.
1288  void UpdateStepStatistics(const T& dt) {
1289  // Handle first step specially.
1290  if (++num_steps_taken_ == 1) {
1293  } else {
1295  }
1296 
1297  // Update the previous step size.
1298  prev_step_size_ = dt;
1299  }
1300 
1301  // Steps the system forward exactly by @p dt, if possible, by calling DoStep.
1302  // Does necessary pre-initialization and post-cleanup. This method does not
1303  // update general integrator statistics (which are updated in the calling
1304  // methods), because error control might decide that it does not like the
1305  // result of the step and might "rewind" and take a smaller one.
1306  // @returns `true` if successful, `false` otherwise (due to, e.g., integrator
1307  // convergence failure).
1308  // @sa DoStep()
1309  bool Step(const T& dt) {
1310  if (!DoStep(dt))
1311  return false;
1312  return true;
1313  }
1314 
1315  // Reference to the system being simulated.
1316  const System<T>& system_;
1317 
1318  // Pointer to the context.
1319  Context<T>* context_{nullptr}; // The trajectory Context.
1320 
1321  // Runtime variables.
1322  // For variable step integrators, this is set at the end of each step to guide
1323  // the next one.
1324  T ideal_next_step_size_{nan()}; // Indicates that the value is uninitialized.
1325 
1326  // The scaling factor to apply to an integration step size when an integrator
1327  // convergence failure occurs (to make convergence more likely on the next
1328  // attempt).
1329  // TODO(edrumwri): Allow subdivision factor to be user-tweakable.
1330  const double subdivision_factor_{0.5};
1331 
1332  // The accuracy being used.
1333  double accuracy_in_use_{nan()};
1334 
1335  // The maximum step size.
1336  T max_step_size_{nan()};
1337 
1338  // The minimum step size.
1339  T req_min_step_size_{0};
1340 
1341  // The last step taken by the integrator.
1342  T prev_step_size_{nan()};
1343 
1344  // Whether error-controlled integrator is running in fixed step mode. Value
1345  // is irrelevant for integrators without error estimation capabilities.
1346  bool fixed_step_mode_{false};
1347 
1348  // When the minimum step is exceeded, does the integrator throw an exception?
1349  bool min_step_exceeded_throws_{true};
1350 
1351  // Statistics.
1352  T actual_initial_step_size_taken_{nan()};
1353  T smallest_adapted_step_size_taken_{nan()};
1354  T largest_step_size_taken_{nan()};
1355  int64_t num_steps_taken_{0};
1356  int64_t num_ode_evals_{0};
1357  int64_t num_shrinkages_from_error_control_{0};
1358  int64_t num_shrinkages_from_substep_failures_{0};
1359  int64_t num_substep_failures_{0};
1360 
1361  // Applied as diagonal matrices to weight state change variables.
1362  Eigen::VectorXd qbar_weight_, z_weight_;
1363 
1364  // State copy for reversion during error-controlled integration.
1365  VectorX<T> xc0_save_;
1366 
1367  // The error estimate computed during integration with error control.
1368  std::unique_ptr<ContinuousState<T>> err_est_;
1369 
1370  // The pseudo-inverse of the matrix that converts time derivatives of
1371  // generalized coordinates to generalized velocities, multiplied by the
1372  // change in the generalized coordinates (used in state change norm
1373  // calculations).
1374  mutable std::unique_ptr<VectorBase<T>> pinvN_dq_change_;
1375 
1376  // Vectors used in state change norm calculations.
1377  mutable VectorX<T> unweighted_substate_change_;
1378  mutable std::unique_ptr<VectorBase<T>> weighted_q_change_;
1379 
1380  // Variable for indicating when an integrator has been initialized.
1381  bool initialization_done_{false};
1382 
1383  // This a workaround for an apparent bug in clang 3.8 in which
1384  // defining this as a static constexpr member kNaN failed to instantiate
1385  // properly for the AutoDiffXd instantiation (worked in gcc and MSVC).
1386  // Restore to sanity when some later clang is current.
1387  static constexpr double nan() {
1388  return std::numeric_limits<double>::quiet_NaN();
1389  }
1390 
1391  double target_accuracy_{nan()}; // means "unspecified, use default"
1392  T req_initial_step_size_{nan()}; // means "unspecified, use default"
1393 };
1394 
1395 template <class T>
1397  using std::isnan;
1398  using std::min;
1399 
1400  // Verify that the integrator supports error estimates.
1401  if (!supports_error_estimation()) {
1402  throw std::logic_error("StepOnceErrorControlledAtMost() requires error "
1403  "estimation.");
1404  }
1405 
1406  // Save time, continuous variables, and time derivative because we'll possibly
1407  // revert time and state.
1408  const Context<T>& context = get_context();
1409  const T current_time = context.get_time();
1410  VectorBase<T>& xc =
1411  get_mutable_context()->get_mutable_continuous_state_vector();
1412  xc0_save_ = xc.CopyToVector();
1413 
1414  // Set the step size to attempt.
1415  T step_size_to_attempt = get_ideal_next_step_size();
1416  if (isnan(step_size_to_attempt)) {
1417  // Integrator has not taken a step. Set the current step size to the
1418  // initial step size.
1419  step_size_to_attempt = get_initial_step_size_target();
1420  DRAKE_DEMAND(!isnan(step_size_to_attempt));
1421  }
1422 
1423  // This variable indicates when the integrator has been pushed to its minimum
1424  // step limit. It can only be "true" if minimum step exceptions have been
1425  // suppressed by the user via set_throw_on_minimum_step_size_violation(false),
1426  // and the error control mechanism determines that the step is as low as it
1427  // can go.
1428  bool at_minimum_step_size = false;
1429 
1430  bool step_succeeded = false;
1431  do {
1432  // Constants used to determine whether modifications to the step size are
1433  // close enough to the attempted step size to use the unadjusted originals,
1434  // or (1) whether the step size to be attempted is so small that we should
1435  // consider it to be artifically limited or (2) whether the step size to
1436  // be attempted is sufficiently close to that requested such that the step
1437  // size should be stretched slightly.
1438  const double near_enough_smaller = 0.95;
1439  const double near_enough_larger = 1.001;
1440 
1441  // If we lose more than a small fraction of the step size we wanted
1442  // to take due to a need to stop at dt_max, make a note of that so the
1443  // step size adjuster won't try to grow from the current step.
1444  bool dt_was_artificially_limited = false;
1445  if (dt_max < near_enough_smaller * step_size_to_attempt) {
1446  // dt_max much smaller than current step size.
1447  dt_was_artificially_limited = true;
1448  step_size_to_attempt = dt_max;
1449  } else {
1450  if (dt_max < near_enough_larger * step_size_to_attempt) {
1451  // dt_max is roughly current step. Make it the step size to prevent
1452  // creating a small sliver (the remaining step).
1453  step_size_to_attempt = dt_max;
1454  }
1455  }
1456 
1457  // Limit the current step size.
1458  step_size_to_attempt = min(step_size_to_attempt, get_maximum_step_size());
1459 
1460  // Keep adjusting the integration step size until any integrator
1461  // convergence failures disappear. Note: this loop's correctness is
1462  // predicated on the assumption that an integrator will always converge for
1463  // a sufficiently small, yet nonzero step size.
1464  T adjusted_step_size = step_size_to_attempt;
1465  while (!Step(adjusted_step_size)) {
1466  SPDLOG_DEBUG(drake::log(), "Sub-step failed at {}", adjusted_step_size);
1467  adjusted_step_size *= subdivision_factor_;
1468 
1469  // Note: we could give the user more rope to hang themselves by looking
1470  // for zero rather than machine epsilon, which might be advantageous if
1471  // the user were modeling systems over extremely small time scales.
1472  // However, that issue could be addressed instead by scaling units, and
1473  // using machine epsilon allows failure to be detected much more rapidly.
1474  if (adjusted_step_size < std::numeric_limits<double>::epsilon()) {
1475  throw std::runtime_error("Integrator has been directed to a near zero-"
1476  "length step in order to obtain convergence.");
1477  }
1478  ValidateSmallerStepSize(step_size_to_attempt, adjusted_step_size);
1479  ++num_shrinkages_from_substep_failures_;
1480  ++num_substep_failures_;
1481  }
1482  step_size_to_attempt = adjusted_step_size;
1483 
1484  //--------------------------------------------------------------------
1485  T err_norm = CalcStateChangeNorm(*get_error_estimate());
1486  T next_step_size;
1487  std::tie(step_succeeded, next_step_size) = CalcAdjustedStepSize(
1488  err_norm, step_size_to_attempt, &at_minimum_step_size);
1489  SPDLOG_DEBUG(drake::log(), "Next step size: {}", next_step_size);
1490 
1491  if (step_succeeded) {
1492  // Only update the next step size (retain the previous one) if the
1493  // step size was not artificially limited.
1494  if (!dt_was_artificially_limited)
1495  ideal_next_step_size_ = next_step_size;
1496 
1498  set_actual_initial_step_size_taken(step_size_to_attempt);
1499 
1500  // Record the adapted step size taken.
1502  (step_size_to_attempt < get_smallest_adapted_step_size_taken() &&
1503  step_size_to_attempt < dt_max))
1504  set_smallest_adapted_step_size_taken(step_size_to_attempt);
1505  } else {
1506  ++num_shrinkages_from_error_control_;
1507 
1508  // Set the next step size to attempt.
1509  step_size_to_attempt = next_step_size;
1510 
1511  // Reset the time, state, and time derivative at t0.
1512  get_mutable_context()->set_time(current_time);
1513  xc.SetFromVector(xc0_save_);
1514  }
1515  } while (!step_succeeded);
1516 
1517  return (step_size_to_attempt == dt_max);
1518 }
1519 
1520 template <class T>
1522  const ContinuousState<T>& dx_state) const {
1523  using std::max;
1524  const Context<T>& context = get_context();
1525  const auto& system = get_system();
1526 
1527  // Get weighting matrices.
1528  const auto& qbar_v_weight = this->get_generalized_state_weight_vector();
1529  const auto& z_weight = this->get_misc_state_weight_vector();
1530 
1531  // Get the differences in the generalized position, velocity, and
1532  // miscellaneous continuous state vectors.
1533  const VectorBase<T>& dgq = dx_state.get_generalized_position();
1534  const VectorBase<T>& dgv = dx_state.get_generalized_velocity();
1535  const VectorBase<T>& dgz = dx_state.get_misc_continuous_state();
1536 
1537  // (re-)Initialize pinvN_dq_change_ and weighted_q_change_, if necessary.
1538  // Reinitialization might be required if the system state variables can
1539  // change during the course of the simulation.
1540  if (pinvN_dq_change_ == nullptr) {
1541  pinvN_dq_change_ = std::make_unique<BasicVector<T>>(dgv.size());
1542  weighted_q_change_ = std::make_unique<BasicVector<T>>(dgq.size());
1543  }
1544  DRAKE_DEMAND(pinvN_dq_change_->size() == dgv.size());
1545  DRAKE_DEMAND(weighted_q_change_->size() == dgq.size());
1546 
1547  // TODO(edrumwri): Acquire characteristic time properly from the system
1548  // (i.e., modify the System to provide this value).
1549  const double characteristic_time = 1.0;
1550 
1551  // Computes the infinity norm of the weighted velocity variables.
1552  unweighted_substate_change_ = dgv.CopyToVector();
1553  T v_nrm = qbar_v_weight.cwiseProduct(unweighted_substate_change_).
1554  template lpNorm<Eigen::Infinity>() * characteristic_time;
1555 
1556  // Compute the infinity norm of the weighted auxiliary variables.
1557  unweighted_substate_change_ = dgz.CopyToVector();
1558  T z_nrm = (z_weight.cwiseProduct(unweighted_substate_change_))
1559  .template lpNorm<Eigen::Infinity>();
1560 
1561  // Compute N * Wq * dq = N * Wꝗ * N+ * dq.
1562  unweighted_substate_change_ = dgq.CopyToVector();
1563  system.MapQDotToVelocity(context, unweighted_substate_change_,
1564  pinvN_dq_change_.get());
1565  system.MapVelocityToQDot(
1566  context, qbar_v_weight.cwiseProduct(pinvN_dq_change_->CopyToVector()),
1567  weighted_q_change_.get());
1568  T q_nrm = weighted_q_change_->CopyToVector().
1569  template lpNorm<Eigen::Infinity>();
1570  SPDLOG_DEBUG(drake::log(), "dq norm: {}, dv norm: {}, dz norm: {}",
1571  q_nrm, v_nrm, z_nrm);
1572 
1573  // TODO(edrumwri): Record the worst offender (which of the norms resulted
1574  // in the largest value).
1575  // Infinity norm of the concatenation of multiple vectors is equal to the
1576  // maximum of the infinity norms of the individual vectors.
1577  return max(z_nrm, max(q_nrm, v_nrm));
1578 }
1579 
1580 template <class T>
1582  const T& err,
1583  const T& step_taken,
1584  bool* at_minimum_step_size) const {
1585  using std::pow;
1586  using std::min;
1587  using std::max;
1588  using std::isnan;
1589  using std::isinf;
1590 
1591  // Magic numbers come from Simbody.
1592  const double kSafety = 0.9;
1593  const double kMinShrink = 0.1;
1594  const double kMaxGrow = 5.0;
1595  const double kHysteresisLow = 0.9;
1596  const double kHysteresisHigh = 1.2;
1597 
1598  // Get the order for the integrator's error estimate.
1599  const int err_order = get_error_estimate_order();
1600 
1601  // Set value for new step size to invalid value initially.
1602  T new_step_size(-1);
1603 
1604  // First, make a guess at the next step size to use based on
1605  // the supplied error norm. Watch out for NaN. Further adjustments will be
1606  // made in blocks of code that follow.
1607  if (isnan(err) || isinf(err)) // e.g., integrand returned NaN.
1608  new_step_size = kMinShrink * step_taken;
1609  else if (err == 0) // A "perfect" step; can happen if no dofs for example.
1610  new_step_size = kMaxGrow * step_taken;
1611  else // Choose best step for skating just below the desired accuracy.
1612  new_step_size = kSafety * step_taken *
1613  pow(get_accuracy_in_use() / err, 1.0 / err_order);
1614 
1615  // Error indicates that the step size can be increased.
1616  if (new_step_size > step_taken) {
1617  // If the integrator has been directed down to the minimum step size, but
1618  // now error indicates that the step size can be increased, de-activate
1619  // at_minimum_step_size.
1620  *at_minimum_step_size = false;
1621 
1622  // If the new step is bigger than the old, don't make the change if the
1623  // old one was small for some unimportant reason (like reached a publishing
1624  // interval). Also, don't grow the step size if the change would be very
1625  // small; better to keep the step size stable in that case (maybe just
1626  // for aesthetic reasons).
1627  if (new_step_size < kHysteresisHigh * step_taken)
1628  new_step_size = step_taken;
1629  }
1630 
1631  // If error indicates that we should shrink the step size but are not allowed
1632  // to, quit and indicate that the step was successful.
1633  if (new_step_size < step_taken && *at_minimum_step_size) {
1634  return std::make_pair(true, step_taken);
1635  }
1636 
1637  // If we're supposed to shrink the step size but the one we have actually
1638  // achieved the desired accuracy last time, we won't change the step now.
1639  // Otherwise, if we are going to shrink the step, let's not be shy -- we'll
1640  // shrink it by at least a factor of kHysteresisLow.
1641  if (new_step_size < step_taken) {
1642  if (err <= get_accuracy_in_use()) {
1643  new_step_size = step_taken; // not this time
1644  } else {
1645  T test_value = kHysteresisLow * step_taken;
1646  new_step_size = min(new_step_size, test_value);
1647  }
1648  }
1649 
1650  // Keep the size change within the allowable bounds.
1651  T max_grow_step = kMaxGrow * step_taken;
1652  T min_shrink_step = kMinShrink * step_taken;
1653  new_step_size = min(new_step_size, max_grow_step);
1654  new_step_size = max(new_step_size, min_shrink_step);
1655 
1656  // Apply user-requested limits on min and max step size.
1657  // TODO(edrumwri): Introduce some feedback to the user when integrator wants
1658  // to take a smaller step than user has selected as the minimum. Options for
1659  // this feedback could include throwing a special exception, logging, setting
1660  // a flag in the integrator that allows throwing an exception, or returning
1661  // a special status from IntegrateAtMost().
1662  if (!isnan(get_maximum_step_size()))
1663  new_step_size = min(new_step_size, get_maximum_step_size());
1664  ValidateSmallerStepSize(step_taken, new_step_size);
1665 
1666  // Increase the next step size, as necessary.
1667  new_step_size = max(new_step_size, get_working_minimum_step_size());
1668  if (new_step_size == get_working_minimum_step_size()) {
1669  // Indicate that the step is integrator is now trying the minimum step
1670  // size.
1671  *at_minimum_step_size = true;
1672 
1673  // If the integrator wants to shrink the step size below the
1674  // minimum allowed and exceptions are suppressed, indicate that status.
1675  if (new_step_size < step_taken)
1676  return std::make_pair(false, new_step_size);
1677  }
1678 
1679  return std::make_pair(new_step_size >= step_taken, new_step_size);
1680 }
1681 
1682 template <class T>
1684  const T& publish_dt, const T& update_dt, const T& boundary_dt) {
1685 
1687  throw std::logic_error("Integrator not initialized.");
1688 
1689  // Now that integrator has been checked for initialization, get the current
1690  // time.
1691  const T t0 = context_->get_time();
1692 
1693  // Verify that update dt is positive.
1694  if (update_dt <= 0.0)
1695  throw std::logic_error("Update dt must be strictly positive.");
1696 
1697  // Verify that other dt's are non-negative.
1698  if (publish_dt < 0.0)
1699  throw std::logic_error("Publish dt is negative.");
1700  if (boundary_dt < 0.0)
1701  throw std::logic_error("Boundary dt is negative.");
1702 
1703  // The size of the integration step is the minimum of the time until the next
1704  // update event, the time until the next publish event, the boundary time
1705  // (i.e., the maximum time that the user wished to step to), and the maximum
1706  // step size (which may stretch slightly to hit a discrete event).
1707 
1708  // We report to the caller which event ultimately constrained the step size.
1709  // If multiple events constrained it equally, we prefer to report update
1710  // events over publish events, publish events over boundary step limits,
1711  // and boundary limits over maximum step size limits. The caller must
1712  // determine event simultaneity by inspecting the time.
1713 
1714  // The maintainer of this code is advised to consider that, while updates
1715  // and boundary times, may both conceptually be deemed events, the distinction
1716  // is made for a reason. If both an update and a boundary time occur
1717  // simultaneously, the following behavior should result:
1718  // (1) kReachedUpdateTime is returned, (2) Simulator::StepTo() performs the
1719  // necessary update, (3) IntegrateAtMost() is called with boundary_dt=0 and
1720  // returns kReachedBoundaryTime, and (4) the simulation terminates. This
1721  // sequence of operations will ensure that the simulation state is valid if
1722  // Simulator::StepTo() is called again to advance time further.
1723 
1724  // We now analyze the following simultaneous cases with respect to Simulator:
1725  //
1726  // { publish, update }
1727  // kReachedUpdateTime will be returned, an update will be followed by a
1728  // publish.
1729  //
1730  // { publish, update, max step }
1731  // kReachedUpdateTime will be returned, an update will be followed by a
1732  // publish.
1733  //
1734  // { publish, boundary time, max step }
1735  // kReachedPublishTime will be returned, a publish will be performed followed
1736  // by another call to this function, which should return kReachedBoundaryTime
1737  // (followed in rapid succession by StepTo(.) return).
1738  //
1739  // { publish, boundary time, max step }
1740  // kReachedPublishTime will be returned, a publish will be performed followed
1741  // by another call to this function, which should return kReachedBoundaryTime
1742  // (followed in rapid succession by StepTo(.) return).
1743  //
1744  // { publish, update, boundary time, maximum step size }
1745  // kUpdateTimeReached will be returned, an update followed by a publish
1746  // will then be performed followed by another call to this function, which
1747  // should return kReachedBoundaryTime (followed in rapid succession by
1748  // StepTo(.) return).
1749 
1750  // By default, the candidate dt is the next discrete update event.
1752  T dt = update_dt;
1753 
1754  // If the next discrete publish event is sooner than the next discrete update
1755  // event, the publish event becomes the candidate dt
1756  if (publish_dt < update_dt) {
1757  candidate_result = IntegratorBase<T>::kReachedPublishTime;
1758  dt = publish_dt;
1759  }
1760 
1761  // If the stop time (boundary time) is sooner than the candidate, use it
1762  // instead.
1763  if (boundary_dt < dt) {
1764  candidate_result = IntegratorBase<T>::kReachedBoundaryTime;
1765  dt = boundary_dt;
1766  }
1767 
1768  // If there is no continuous state, there will be no need to limit the
1769  // integration step size.
1770  if (get_context().get_continuous_state().size() == 0) {
1771  Context<T>* context = get_mutable_context();
1772  context->set_time(context->get_time() + dt);
1773  return candidate_result;
1774  }
1775 
1776  // If all events are farther into the future than the maximum step
1777  // size times a stretch factor of 1.01, the maximum step size becomes the
1778  // candidate dt. Put another way, if the maximum step occurs right before
1779  // an update or a publish, the update or publish is done instead. In contrast,
1780  // we never step past boundary_dt, even if doing so would allow hitting a
1781  // publish or an update.
1782  const bool reached_boundary =
1783  (candidate_result == IntegratorBase<T>::kReachedBoundaryTime);
1784  const T& max_dt = IntegratorBase<T>::get_maximum_step_size();
1785  if ((reached_boundary && max_dt < dt) ||
1786  (!reached_boundary && max_dt * get_stretch_factor() < dt)) {
1787  candidate_result = IntegratorBase<T>::kTimeHasAdvanced;
1788  dt = max_dt;
1789  }
1790 
1791  if (dt < 0.0) throw std::logic_error("Negative dt.");
1792 
1793  // If error control is disabled, call the generic stepper. Otherwise, use
1794  // the error controlled method.
1795  bool full_step = true;
1796  if (this->get_fixed_step_mode()) {
1797  T adjusted_dt = dt;
1798  while (!Step(adjusted_dt)) {
1799  ++num_shrinkages_from_substep_failures_;
1800  ++num_substep_failures_;
1801  adjusted_dt *= subdivision_factor_;
1802  ValidateSmallerStepSize(dt, adjusted_dt);
1803  full_step = false;
1804  }
1805  } else {
1806  full_step = StepOnceErrorControlledAtMost(dt);
1807  }
1808 
1809  // Update generic statistics.
1810  const T actual_dt = context_->get_time() - t0;
1811  UpdateStepStatistics(actual_dt);
1812 
1813  if (full_step) {
1814  // If the integrator took the entire maximum step size we allowed above,
1815  // we report to the caller that a step constraint was hit, which may
1816  // indicate a discrete event has arrived.
1817  return candidate_result;
1818  } else {
1819  // Otherwise, we report to the caller that time has advanced, but no
1820  // discrete event has arrived.
1822  }
1823 }
1824 
1825 } // namespace systems
1826 } // namespace drake
void CalcTimeDerivatives(const System< U > &system, const Context< U > &context, ContinuousState< U > *dxdt)
Evaluates the derivative function (and updates call statistics).
Definition: integrator_base.h:1085
T get_working_minimum_step_size() const
Gets the current value of the working minimum step size h_work(t) for this integrator, which may vary with the current time t as stored in the integrator&#39;s context.
Definition: integrator_base.h:387
virtual int size() const =0
Returns the number of elements in the vector.
double get_target_accuracy() const
Gets the target accuracy.
Definition: integrator_base.h:242
void set_largest_step_size_taken(const T &dt)
Definition: integrator_base.h:1261
void CalcTimeDerivatives(const Context< T > &context, ContinuousState< T > *dxdt)
Evaluates the derivative function (and updates call statistics).
Definition: integrator_base.h:1073
bool isnan(const Eigen::AutoDiffScalar< DerType > &x)
Overloads isnan to mimic std::isnan from <cmath>.
Definition: autodiff_overloads.h:52
Localized an event; this is the before state (interpolated).
Definition: integrator_base.h:131
void set_ideal_next_step_size(const T &dt)
Definition: integrator_base.h:1266
const Eigen::VectorXd & get_generalized_state_weight_vector() const
Gets the weighting vector (equivalent to a diagonal matrix) applied to weighting both generalized coo...
Definition: integrator_base.h:1012
const T & get_actual_initial_step_size_taken() const
The actual size of the successful first step.
Definition: integrator_base.h:713
virtual void DoReset()
Derived classes can override this method to perform routines when Reset() is called.
Definition: integrator_base.h:1218
#define SPDLOG_DEBUG(logger,...)
Definition: text_logging.h:93
Definition: automotive_demo.cc:88
Context< T > * get_mutable_context()
Returns a mutable pointer to the internally-maintained Context holding the most recent state in the t...
Definition: integrator_base.h:764
const T & get_previous_integration_step_size() const
Gets the size of the last (previous) integration step.
Definition: integrator_base.h:800
void Initialize()
An integrator must be initialized before being used.
Definition: integrator_base.h:446
double get_stretch_factor() const
Gets the stretch factor (> 1), which is multiplied by the maximum (typically user-designated) integra...
Definition: integrator_base.h:560
void InitializeAccuracy(double default_accuracy, double loosest_accuracy, double max_step_fraction)
Generic code for validating (and resetting, if need be) the integrator working accuracy for error con...
Definition: integrator_base.h:1113
Context is an abstract base class template that represents all the inputs to a System: time...
Definition: query_handle.h:10
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > VectorX
A column vector of any size, templated on scalar type.
Definition: eigen_types.h:46
VectorBase is an abstract base class that real-valued signals between Systems and real-valued System ...
Definition: vector_base.h:27
Eigen::VectorBlock< Eigen::VectorXd > get_mutable_generalized_state_weight_vector()
Gets a mutable weighting vector (equivalent to a diagonal matrix) applied to weighting both generaliz...
Definition: integrator_base.h:1029
const VectorBase< T > & get_misc_continuous_state() const
Returns a const reference to the subset of the continuous state vector that is other continuous state...
Definition: continuous_state.h:142
void set_accuracy_in_use(double accuracy)
Sets the working ("in use") accuracy for this integrator.
Definition: integrator_base.h:1099
void set_throw_on_minimum_step_size_violation(bool throws)
Sets whether the integrator should throw a std::runtime_error exception when the integrator&#39;s step si...
Definition: integrator_base.h:368
logging::logger * log()
Retrieve an instance of a logger to use for logging; for example: drake::log()->info("potato!") ...
Definition: text_logging.cc:38
void request_initial_step_size_target(const T &step_size)
Request that the first attempted integration step have a particular size.
Definition: integrator_base.h:503
void set_maximum_step_size(const T &max_step_size)
Sets the maximum step size that may be taken by this integrator.
Definition: integrator_base.h:258
const Context< T > & get_context() const
Returns a const reference to the internally-maintained Context holding the most recent state in the t...
Definition: integrator_base.h:758
virtual VectorX< T > CopyToVector() const
Copies the entire state to a vector with no semantics.
Definition: vector_base.h:98
void IntegrateWithSingleFixedStep(const T &dt)
Stepping function for integrators operating outside of Simulator that advances the continuous state e...
Definition: integrator_base.h:637
ContinuousState< T > * get_mutable_error_estimate()
Gets an error estimate of the state variables recorded by the last call to StepOnceFixedSize().
Definition: integrator_base.h:1245
const T & get_initial_step_size_target() const
Gets the target size of the first integration step.
Definition: integrator_base.h:518
#define DRAKE_ASSERT(condition)
DRAKE_ASSERT(condition) is similar to the built-in assert(condition) from the C++ system header <cas...
Definition: drake_assert.h:37
bool get_fixed_step_mode() const
Gets whether an integrator is running in fixed step mode.
Definition: integrator_base.h:200
StepResult
Status returned by StepOnceAtMost().
Definition: integrator_base.h:123
const ContinuousState< T > * get_error_estimate() const
Gets the error estimate (used only for integrators that support error estimation).
Definition: integrator_base.h:809
virtual void set_time(const T &time_sec)
Set the current time in seconds.
Definition: context.h:52
User requested control whenever an internal step is successful.
Definition: integrator_base.h:139
bool StepOnceErrorControlledAtMost(const T &dt_max)
Default code for advancing the continuous state of the system by a single step of dt_max (or smaller...
Definition: integrator_base.h:1396
Provides Drake&#39;s assertion implementation.
This is the entry point for all text logging within Drake.
Expression max(const Expression &e1, const Expression &e2)
Definition: symbolic_expression.cc:697
StepResult IntegrateAtMost(const T &publish_dt, const T &update_dt, const T &boundary_dt)
Integrates the system forward in time by a single step with step size subject to integration error to...
Definition: integrator_base.h:1683
Took maximum number of steps without finishing integrating over the interval.
Definition: integrator_base.h:146
#define DRAKE_DEMAND(condition)
Evaluates condition and iff the value is false will trigger an assertion failure with a message showi...
Definition: drake_assert.h:45
int64_t get_num_derivative_evaluations() const
Returns the number of ODE function evaluations (calls to CalcTimeDerivatives()) since the last call t...
Definition: integrator_base.h:708
const System< T > & get_system() const
Gets a constant reference to the system that is being integrated (and was provided to the constructor...
Definition: integrator_base.h:784
Expression min(const Expression &e1, const Expression &e2)
Definition: symbolic_expression.cc:685
void set_fixed_step_mode(bool flag)
Sets an integrator with error control to fixed step mode.
Definition: integrator_base.h:185
Reached the desired integration time without reaching an update time.
Definition: integrator_base.h:143
void Reset()
Resets the integrator to initial values, i.e., default construction values.
Definition: integrator_base.h:400
void set_actual_initial_step_size_taken(const T &dt)
Definition: integrator_base.h:1248
const T & get_ideal_next_step_size() const
Return the step size the integrator would like to take next, based primarily on the integrator&#39;s accu...
Definition: integrator_base.h:751
const T & get_time() const
Returns the current time in seconds.
Definition: context.h:49
double get_accuracy_in_use() const
Gets the accuracy in use by the integrator.
Definition: integrator_base.h:249
const T & get_maximum_step_size() const
Gets the maximum step size that may be taken by this integrator.
Definition: integrator_base.h:270
const VectorBase< T > & get_generalized_position() const
Returns a const reference to the subset of the state vector that is generalized position q...
Definition: continuous_state.h:118
A superclass template for systems that receive input, maintain state, and produce output of a given m...
Definition: input_port_descriptor.h:12
const Eigen::VectorXd & get_misc_state_weight_vector() const
Gets the weighting vector (equivalent to a diagonal matrix) for weighting errors in miscellaneous con...
Definition: integrator_base.h:1039
void set_requested_minimum_step_size(const T &min_step_size)
Sets the requested minimum step size h_min that may be taken by this integrator.
Definition: integrator_base.h:347
virtual void SetFromVector(const Eigen::Ref< const VectorX< T >> &value)
Replaces the entire vector with the contents of value.
Definition: vector_base.h:80
virtual ~IntegratorBase()=default
Destructor.
int64_t get_num_steps_taken() const
The number of integration steps taken since the last Initialize() or ResetStatistics() call...
Definition: integrator_base.h:739
const T & get_smallest_adapted_step_size_taken() const
The size of the smallest step taken as the result of a controlled integration step adjustment since t...
Definition: integrator_base.h:723
virtual void DoInitialize()
Derived classes can override this method to perform special initialization.
Definition: integrator_base.h:1212
void CalcTimeDerivatives(const Context< T > &context, ContinuousState< T > *derivatives) const
Calculates the time derivatives xcdot of the continuous state xc.
Definition: system.h:543
T CalcStateChangeNorm(const ContinuousState< T > &dx_state) const
Computes the infinity norm of a change in continuous state.
Definition: integrator_base.h:1521
int64_t get_num_step_shrinkages_from_substep_failures() const
Gets the number of step size shrinkages due to sub-step failures (e.g., integrator convergence failur...
Definition: integrator_base.h:689
void IntegrateWithMultipleSteps(const T &dt)
Stepping function for integrators operating outside of Simulator that advances the continuous state e...
Definition: integrator_base.h:584
Indicates that integration terminated at an update time.
Definition: integrator_base.h:135
void reset_context(Context< T > *context)
Replace the pointer to the internally-maintained Context with a different one.
Definition: integrator_base.h:775
virtual bool supports_error_estimation() const =0
Indicates whether an integrator supports error estimation.
IntegratorBase(const System< T > &system, Context< T > *context=nullptr)
Maintains references to the system being integrated and the context used to specify the initial condi...
Definition: integrator_base.h:162
int64_t get_num_step_shrinkages_from_error_control() const
Gets the number of step size shrinkages due to failure to meet targeted error tolerances, since the last call to ResetStatistics or Initialize().
Definition: integrator_base.h:695
bool get_throw_on_minimum_step_size_violation() const
Reports the current setting of the throw_on_minimum_step_size_violation flag.
Definition: integrator_base.h:377
virtual bool DoStep(const T &dt)=0
Derived classes must implement this method to (1) integrate the continuous portion of this system for...
virtual void DoResetStatistics()
Resets any statistics particular to a specific integrator.
Definition: integrator_base.h:1068
Eigen::VectorBlock< Eigen::VectorXd > get_mutable_misc_state_weight_vector()
Gets a mutable weighting vector (equivalent to a diagonal matrix) for weighting errors in miscellaneo...
Definition: integrator_base.h:1054
const VectorBase< T > & get_generalized_velocity() const
Returns a const reference to the subset of the continuous state vector that is generalized velocity v...
Definition: continuous_state.h:130
An abstract class for an integrator for ODEs and DAEs as represented by a Drake System.
Definition: integrator_base.h:107
bool is_initialized() const
Indicates whether the integrator has been initialized.
Definition: integrator_base.h:787
ContinuousState is a container for all the continuous state variables xc.
Definition: continuous_state.h:44
bool isinf(const Eigen::AutoDiffScalar< DerType > &x)
Overloads isinf to mimic std::isinf from <cmath>.
Definition: autodiff_overloads.h:45
Polynomial< CoefficientType > pow(const Polynomial< CoefficientType > &base, typename Polynomial< CoefficientType >::PowerType exponent)
Provides power function for Polynomial.
Definition: polynomial.h:449
const T & get_requested_minimum_step_size() const
Gets the requested minimum step size h_min for this integrator.
Definition: integrator_base.h:357
#define DRAKE_NO_COPY_NO_MOVE_NO_ASSIGN(Classname)
DRAKE_NO_COPY_NO_MOVE_NO_ASSIGN deletes the special member functions for copy-construction, copy-assignment, move-construction, and move-assignment.
Definition: drake_copyable.h:33
virtual int get_error_estimate_order() const =0
Derived classes must override this function to return the order of the integrator&#39;s error estimate...
void set_target_accuracy(double accuracy)
Request that the integrator attempt to achieve a particular accuracy for the continuous portions of t...
Definition: integrator_base.h:229
std::pair< bool, T > CalcAdjustedStepSize(const T &err, const T &attempted_step_size, bool *at_minimum_step_size) const
Calculates adjusted integrator step sizes toward keeping state variables within error bounds on the n...
Definition: integrator_base.h:1581
Indicates a publish time has been reached but not an update time.
Definition: integrator_base.h:127
int64_t get_num_substep_failures() const
Gets the number of failed sub-steps (implying one or more step size reductions was required to permit...
Definition: integrator_base.h:682
void ResetStatistics()
Forget accumulated statistics.
Definition: integrator_base.h:667
Provides careful macros to selectively enable or disable the special member functions for copy-constr...
void set_smallest_adapted_step_size_taken(const T &dt)
Sets the size of the smallest-step-taken statistic as the result of a controlled integration step adj...
Definition: integrator_base.h:1256
const T & get_largest_step_size_taken() const
The size of the largest step taken since the last Initialize() or ResetStatistics() call...
Definition: integrator_base.h:731