An abstract class for an integrator for ODEs and DAEs as represented by a Drake System.
Integrators solve initial value problems of the form:
ẋ(t) = f(t, x(t)) with f : ℝ × ℝⁿ → ℝⁿ
(i.e., f()
is an ordinary differential equation) given initial conditions (t₀, x₀). Thus, integrators advance the continuous state of a dynamical system forward in time.
Drake's subclasses of IntegratorBase<T> should follow the naming pattern FooIntegrator<T>
by convention.
T | The scalar type, which must be one of the default scalars. |
#include <drake/systems/analysis/integrator_base.h>
Public Types | |
enum | StepResult { kReachedPublishTime = 1, kReachedZeroCrossing = 2, kReachedUpdateTime = 3, kTimeHasAdvanced = 4, kReachedBoundaryTime = 5, kReachedStepLimit = 6 } |
Status returned by IntegrateNoFurtherThanTime(). More... | |
Public Member Functions | |
IntegratorBase (const System< T > &system, Context< T > *context=nullptr) | |
Maintains references to the system being integrated and the context used to specify the initial conditions for that system (if any). More... | |
virtual | ~IntegratorBase ()=default |
void | Reset () |
Resets the integrator to initial values, i.e., default construction values. More... | |
void | Initialize () |
An integrator must be initialized before being used. More... | |
StepResult | IntegrateNoFurtherThanTime (const T &publish_time, const T &update_time, const T &boundary_time) |
(Internal use only) Integrates the system forward in time by a single step with step size subject to integration error tolerances (assuming that the integrator supports error estimation). More... | |
void | IntegrateWithMultipleStepsToTime (const T &t_final) |
Stepping function for integrators operating outside of Simulator that advances the continuous state exactly to t_final . More... | |
bool | IntegrateWithSingleFixedStepToTime (const T &t_target) |
Stepping function for integrators operating outside of Simulator that advances the continuous state using a single step to t_target . More... | |
const Context< T > & | get_context () const |
Returns a const reference to the internally-maintained Context holding the most recent state in the trajectory. More... | |
Context< T > * | get_mutable_context () |
Returns a mutable pointer to the internally-maintained Context holding the most recent state in the trajectory. More... | |
void | reset_context (Context< T > *context) |
Replace the pointer to the internally-maintained Context with a different one. More... | |
const System< T > & | get_system () const |
Gets a constant reference to the system that is being integrated (and was provided to the constructor of the integrator). More... | |
bool | is_initialized () const |
Indicates whether the integrator has been initialized. More... | |
const T & | get_previous_integration_step_size () const |
Gets the size of the last (previous) integration step. More... | |
Does not allow copy, move, or assignment | |
IntegratorBase (const IntegratorBase &)=delete | |
IntegratorBase & | operator= (const IntegratorBase &)=delete |
IntegratorBase (IntegratorBase &&)=delete | |
IntegratorBase & | operator= (IntegratorBase &&)=delete |
Methods for getting and setting integrator accuracy | |
The precise meaning of accuracy is a complicated discussion, but it translates roughly to the number of significant digits you want in the results. By convention it is supplied as Integrators vary in the range of accuracy (loosest to tightest) that they can support, and each integrator will choose a default accuracy to be used that lies somewhere within this range and attempts to balance computation and accuracy. If you request accuracy outside the supported range for the chosen integrator it will be quietly adjusted to be in range. You can find out the accuracy setting actually being used using Implicit integrators additionally use the accuracy setting for determining when the underlying Newton-Raphson root finding process has converged. For those integrators, the accuracy setting also limits the allowable iteration error in the Newton-Raphson process. Looser accuracy in that process certainly implies greater error in the ODE solution and might impact the stability of the solution negatively as well.
| |
void | set_target_accuracy (double accuracy) |
Request that the integrator attempt to achieve a particular accuracy for the continuous portions of the simulation. More... | |
double | get_target_accuracy () const |
Gets the target accuracy. More... | |
double | get_accuracy_in_use () const |
Gets the accuracy in use by the integrator. More... | |
Methods related to error estimation and control | |
Established methods for integrating ordinary differential equations invariably make provisions for estimating the "local error" (i.e., the error over a small time interval) of a solution. Although the relationship between local error and global error (i.e., the accumulated error over multiple time steps) can be tenuous, such error estimates can allow integrators to work adaptively, subdividing time intervals as necessary (if, e.g., the system is particularly dynamic or stationary in an interval). Even for applications that do not recommend such adaptive integration- like direct transcription methods for trajectory optimization- error estimation allows the user to assess the accuracy of the solution. | |
virtual bool | supports_error_estimation () const =0 |
Derived classes must override this function to indicate whether the integrator supports error estimation. More... | |
virtual int | get_error_estimate_order () const =0 |
Derived classes must override this function to return the order of the asymptotic term in the integrator's error estimate. More... | |
const ContinuousState< T > * | get_error_estimate () const |
Gets the error estimate (used only for integrators that support error estimation). More... | |
const T & | get_ideal_next_step_size () const |
Return the step size the integrator would like to take next, based primarily on the integrator's accuracy prediction. More... | |
void | set_fixed_step_mode (bool flag) |
Sets an integrator with error control to fixed step mode. More... | |
bool | get_fixed_step_mode () const |
Gets whether an integrator is running in fixed step mode. More... | |
Methods for weighting state variable errors (in the context of error control) | |
This group of methods describes how errors for state variables with heterogeneous units are weighted in the context of error-controlled integration. This is an advanced topic and most users can simply specify desired accuracy and accept the default state variable weights. A collection of state variables is generally defined in heterogeneous units (e.g. length, angles, velocities, energy). Some of the state variables cannot even be expressed in meaningful units, like quaternions. Certain integrators provide an estimate of the absolute error made in each state variable during an integration step. These errors must be properly weighted to obtain an "accuracy" with respect to each particular variable. These per-variable accuracy determinations can be compared against the user's requirements and used to select an appropriate size for the next step [Sherman 2011]. The weights are normally determined automatically using the system's characteristic dimensions, so most users can stop reading now! Custom weighting is primarily useful for performance improvement; an optimal weighting would allow an error-controlled integrator to provide the desired level of accuracy across all state variables without wasting computation achieving superfluous accuracy for some of those variables. Users interested in more precise control over state variable weighting may use the methods in this group to access and modify weighting factors for individual state variables. Changes to these weights can only be made prior to integrator initialization or as a result of an event being triggered and then followed by re-initialization. Relative versus absolute accuracy: State variable integration error, as estimated by an integrator, is an absolute quantity with the same units as the variable. At each time step we therefore need to determine an absolute error that would be deemed "good enough", i.e. satisfies the user's accuracy requirement. If a variable is maintained to a relative accuracy then that "good enough" value is defined to be the required accuracy How to choose weights: The weight Some subtleties for second-order dynamic systems: Systems governed by 2nd-order differential equations are typically split into second order (configuration) variables q, and rate (velocity) variables v, where the time derivatives qdot of q are linearly related to v by the kinematic differential equation Note that generalized quasi-coordinates To summarize, separate weights can be provided for each of
Weights on the generalized velocity variables How the weights are used: The errors in the The weighting matrices Eq = Wq Ev: Ev(i,i) = { min(Wv(i,i), 1/|vᵢ|) if vᵢ is relative { Wv(i,i) if vᵢ is absolute Ez: Ez(i,i) = { min(Wz(i,i), 1/|zᵢ|) if zᵢ is relative { Wz(i,i) if zᵢ is absolute ( Now given an error estimate vector Determining weights for q: The kinematic differential equations v = N⁺ qdot Inverse kinematic differential equation dꝗ/dt = N⁺ dq/dt Use synonyms for v and qdot dꝗ = N⁺ dq Change time derivatives to differentials Wꝗ dꝗ = Wꝗ N⁺ dq Pre-multiply both sides by Wꝗ N Wꝗ dꝗ = N Wꝗ N⁺ dq Pre-multiply both sides by N N Wꝗ dꝗ = Wq dq Define Wq := N Wꝗ N⁺ N Wꝗ v = Wq qdot Back to time derivatives. The last two equations show that Finally, note that a diagonal entry of one of the weighting matrices can be set to zero to disable error estimation for that state variable (i.e., auxiliary variable or configuration/velocity variable pair), but that setting an entry to a negative value will cause an exception to be thrown when the integrator is initialized.
| |
const Eigen::VectorXd & | get_generalized_state_weight_vector () const |
Gets the weighting vector (equivalent to a diagonal matrix) applied to weighting both generalized coordinate and velocity state variable errors, as described in the group documentation. More... | |
Eigen::VectorBlock< Eigen::VectorXd > | get_mutable_generalized_state_weight_vector () |
Gets a mutable weighting vector (equivalent to a diagonal matrix) applied to weighting both generalized coordinate and velocity state variable errors, as described in the group documentation. More... | |
const Eigen::VectorXd & | get_misc_state_weight_vector () const |
Gets the weighting vector (equivalent to a diagonal matrix) for weighting errors in miscellaneous continuous state variables z . More... | |
Eigen::VectorBlock< Eigen::VectorXd > | get_mutable_misc_state_weight_vector () |
Gets a mutable weighting vector (equivalent to a diagonal matrix) for weighting errors in miscellaneous continuous state variables z . More... | |
Methods related to initial step size | |
From [Watts 1983], "One of the more critical issues in solving ordinary differential equations by a step-by-step process occurs in the starting phase. Somehow the procedure must be supplied with an initial step size that is on scale for the problem at hand. It must be small enough to yield a reliable solution by the process, but it should not be so small as to significantly affect the efficiency of solution. The more important of these two possibilities is obviously the reliability question. The first step taken by the code must reflect how fast the solution changes near the initial point. For general purpose computing, an automatic step size adjustment procedure for choosing subsequent steps is essential to produce an accurate solution efficiently. This step size control is usually based on estimates of the local errors incurred by the numerical method. Because most codes also employ algorithmic devices which restrict the step size control to be moderately varying (for reliability), subsequent steps usually tend to stay on scale of the problem. This is not always so, as sometimes happens when working with crude tolerances on problems having rapidly varying components. Nevertheless, most step size adjustment procedures deal reasonably well with all but the most abrupt changes, leaving the most serious danger confined to the starting step size." Users may not have a good idea of an initial step size to take, so integration codes usually attempt to automatically select an initial step size. Sophisticated algorithms for initial step size selection are described in [Hindmarsh 1980], [Watts 1983], [Gladwell 1987], and [Hairer 2008]. These algorithms can fail to produce a good initial step size as well (see discussion in [Watts 1983]). Drake's integrators use a fraction (generally 1/10th of the maximum step size) to set the initial step size. If you have a problem that operates at wildly varying time scales, e.g., Robertson's canonical stiff system problem (that requires a large maximum step size to be efficient), consider setting both the initial and maximum step sizes (i.e., not using the defaults) to keep from missing phenomena that occur over small time scales near the beginning of the time interval being integrated.
| |
void | request_initial_step_size_target (const T &step_size) |
Request that the first attempted integration step have a particular size. More... | |
const T & | get_initial_step_size_target () const |
Gets the target size of the first integration step. More... | |
Methods related to maximum integration step size | |
Sets the nominal maximum step size- the actual maximum step size taken may be slightly larger (see set_maximum_step_size() and get_stretch_factor())- that an integrator will take. Each integrator has a default maximum step size, which might be infinite. | |
void | set_maximum_step_size (const T &max_step_size) |
Sets the maximum step size that may be taken by this integrator. More... | |
const T & | get_maximum_step_size () const |
Gets the maximum step size that may be taken by this integrator. More... | |
double | get_stretch_factor () const |
Gets the stretch factor (> 1), which is multiplied by the maximum (typically user-designated) integration step size to obtain the amount that the integrator is able to stretch the maximum time step toward hitting an upcoming publish or update event in IntegrateNoFurtherThanTime(). More... | |
Methods related to minimum integration step size selection and behavior | |
Variable step integrators reduce their step sizes as needed to achieve requirements such as specified accuracy or step convergence. However, it is not possible to take an arbitrarily small step. Normally integrators choose an appropriate minimum step and throw an exception if the requirements can't be achieved without going below that. Methods in this section allow you to influence two aspects of this procedure:
By default, integrators allow a very small minimum step which can result in long run times. Setting a larger minimum can be helpful as a diagnostic to figure out what aspect of your simulation is requiring small steps. You can set the minimum to what should be a "reasonable" minimum based on what you know about the physical system. You will then get an std::runtime_error exception thrown at any point in time where your model behaves unexpectedly (due to, e.g., a discontinuity in the derivative evaluation function). If you disable the exception (via DetailsBecause time is maintained to finite precision, the integrator uses a scalar You may request a larger minimum step size Under some circumstances the integrator may legitimately take a step of size | |
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. More... | |
const T & | get_requested_minimum_step_size () const |
Gets the requested minimum step size h_min for this integrator. More... | |
void | set_throw_on_minimum_step_size_violation (bool throws) |
Sets whether the integrator should throw a std::exception when the integrator's step size selection algorithm determines that it must take a step smaller than the minimum step size (for, e.g., purposes of error control). More... | |
bool | get_throw_on_minimum_step_size_violation () const |
Reports the current setting of the throw_on_minimum_step_size_violation flag. More... | |
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's context. More... | |
Integrator statistics methods | |
These methods allow the caller to manipulate and query integrator statistics. Generally speaking, the larger the integration step taken, the faster a simulation will run. These methods allow querying (and resetting) the integrator statistics as one means of determining how to make a simulation run faster. | |
void | ResetStatistics () |
Forget accumulated statistics. More... | |
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 solving the necessary nonlinear system of equations). More... | |
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 failures) since the last call to ResetStatistics() or Initialize(). More... | |
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(). More... | |
int64_t | get_num_derivative_evaluations () const |
Returns the number of ODE function evaluations (calls to CalcTimeDerivatives()) since the last call to ResetStatistics() or Initialize(). More... | |
const T & | get_actual_initial_step_size_taken () const |
The actual size of the successful first step. More... | |
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 the last Initialize() or ResetStatistics() call. More... | |
const T & | get_largest_step_size_taken () const |
The size of the largest step taken since the last Initialize() or ResetStatistics() call. More... | |
int64_t | get_num_steps_taken () const |
The number of integration steps taken since the last Initialize() or ResetStatistics() call. More... | |
void | add_derivative_evaluations (double evals) |
Manually increments the statistic for the number of ODE evaluations. More... | |
Methods for dense output computation | |
In general, dense output computations entail both CPU load and memory footprint increases during numerical integration. For some applications, the performance penalty may be prohibitive. As such, these computations are only carried out by explicit user request. The API to start and stop a dense integration process (i.e. a numerical integration process that also computes dense output) is consistent with this design choice. Once dense integration is started, and until it is stopped, all subsequent integration steps taken will update the allocated dense output. | |
void | StartDenseIntegration () |
Starts dense integration, allocating a new dense output for this integrator to use. More... | |
const trajectories::PiecewisePolynomial< T > * | get_dense_output () const |
Returns a const pointer to the integrator's current PiecewisePolynomial instance, holding a representation of the continuous state trajectory since the last StartDenseIntegration() call. More... | |
std::unique_ptr< trajectories::PiecewisePolynomial< T > > | StopDenseIntegration () |
Stops dense integration, yielding ownership of the current dense output to the caller. More... | |
Protected Member Functions | |
virtual void | DoResetStatistics () |
Resets any statistics particular to a specific integrator. More... | |
const ContinuousState< T > & | EvalTimeDerivatives (const Context< T > &context) |
Evaluates the derivative function and updates call statistics. More... | |
template<typename U > | |
const ContinuousState< U > & | EvalTimeDerivatives (const System< U > &system, const Context< U > &context) |
Evaluates the derivative function (and updates call statistics). More... | |
void | set_accuracy_in_use (double accuracy) |
Sets the working ("in use") accuracy for this integrator. More... | |
bool | StepOnceErrorControlledAtMost (const T &h_max) |
Default code for advancing the continuous state of the system by a single step of h_max (or smaller, depending on error control). More... | |
T | CalcStateChangeNorm (const ContinuousState< T > &dx_state) const |
Computes the infinity norm of a change in continuous state. More... | |
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 next integration step. More... | |
virtual void | DoInitialize () |
Derived classes can override this method to perform special initialization. More... | |
virtual void | DoReset () |
Derived classes can override this method to perform routines when Reset() is called. More... | |
trajectories::PiecewisePolynomial< T > * | get_mutable_dense_output () |
Returns a mutable pointer to the internally-maintained PiecewisePolynomial instance, holding a representation of the continuous state trajectory since the last time StartDenseIntegration() was called. More... | |
virtual bool | DoStep (const T &h)=0 |
Derived classes must implement this method to (1) integrate the continuous portion of this system forward by a single step of size h and (2) set the error estimate (via get_mutable_error_estimate()). More... | |
bool | DoDenseStep (const T &h) |
Calls DoStep(h) while recording the resulting step in the dense output. More... | |
ContinuousState< T > * | get_mutable_error_estimate () |
Gets an error estimate of the state variables recorded by the last call to StepOnceFixedSize(). More... | |
void | set_actual_initial_step_size_taken (const T &h) |
void | set_smallest_adapted_step_size_taken (const T &h) |
Sets the size of the smallest-step-taken statistic as the result of a controlled integration step adjustment. More... | |
void | set_largest_step_size_taken (const T &h) |
void | set_ideal_next_step_size (const T &h) |
enum StepResult |
Status returned by IntegrateNoFurtherThanTime().
When a step is successful, it will return an indication of what caused it to stop where it did. When unsuccessful it will throw an exception so you won't see any return value. When return of control is due ONLY to reaching a publish time, (status is kReachedPublishTime) the context may return an interpolated value at an earlier time.
|
delete |
|
delete |
|
explicit |
Maintains references to the system being integrated and the context used to specify the initial conditions for that system (if any).
system | A reference to the system to be integrated; the integrator will maintain a reference to the system in perpetuity, so the integrator must not outlive the system. |
context | A pointer to a writeable context (nullptr is ok, but a non-null pointer must be set before Initialize() is called). The integrator will advance the system state using the pointer to this context. The pointer to the context will be maintained internally. The integrator must not outlive the context. |
|
virtualdefault |
void add_derivative_evaluations | ( | double | evals | ) |
Manually increments the statistic for the number of ODE evaluations.
|
protected |
Calculates adjusted integrator step sizes toward keeping state variables within error bounds on the next integration step.
Note that it is not guaranteed that the (possibly) reduced step size will keep state variables within error bounds; however, the process of (1) taking a trial integration step, (2) calculating the error, and (3) adjusting the step size can be repeated until convergence.
err | The norm of the integrator error that was computed using attempted_step_size . | |
attempted_step_size | The step size that was attempted. | |
[in,out] | at_minimum_step_size | If true on entry, the error control mechanism is not allowed to shrink the step because the integrator is stepping at the minimum step size (note that this condition will only occur if get_throw_on_minimum_step_size_violation() == false - an exception would be thrown otherwise). If true on entry and false on exit, the error control mechanism has managed to increase the step size above the working minimum; if true on entry and true on exit, error control would like to shrink the step size but cannot. If false on entry and true on exit, error control shrank the step to the working minimum step size. |
true
if the integration step was to be considered successful and false
otherwise. The value of the T type will be set to the recommended next step size.
|
protected |
Computes the infinity norm of a change in continuous state.
We use the infinity norm to capture the idea that, by providing accuracy requirements, the user can indirectly specify error tolerances that act to limit the largest error in any state vector component.
|
protected |
Calls DoStep(h) while recording the resulting step in the dense output.
If the current dense output is already non-empty, then the time in the current context must match either the final segment time of the dense output, or the penultimate segment time (to support the case where the same integration step is attempted multiple times, which occurs e.g. in witness function isolation).
h | The integration step to take. |
true
if successful, false
if either the integrator was unable to take a single step of size h
or to advance its dense output an equal step.
|
protectedvirtual |
Derived classes can override this method to perform special initialization.
This method is called during the Initialize() method. This default method does nothing.
|
protectedvirtual |
Derived classes can override this method to perform routines when Reset() is called.
This default method does nothing.
Reimplemented in ImplicitIntegrator< T >.
|
protectedvirtual |
Resets any statistics particular to a specific integrator.
The default implementation of this function does nothing. If your integrator collects its own statistics, you should re-implement this method and reset them there.
Reimplemented in ImplicitIntegrator< T >.
|
protectedpure virtual |
Derived classes must implement this method to (1) integrate the continuous portion of this system forward by a single step of size h
and (2) set the error estimate (via get_mutable_error_estimate()).
This method is called during the integration process (via StepOnceErrorControlledAtMost(), IntegrateNoFurtherThanTime(), and IntegrateWithSingleFixedStepToTime()).
h | The integration step to take. |
true
if successful, false
if the integrator was unable to take a single step of size h
(due to, e.g., an integrator convergence failure). t
, the time and state will be advanced to t+h
if the method returns true
; otherwise, the time and state should be reset to those at t
. true
for some, albeit possibly very small, positive value of h
. The derived integrator's stepping algorithm can make this guarantee, for example, by switching to an algorithm not subject to convergence failures (e.g., explicit Euler) for very small step sizes.
|
protected |
Evaluates the derivative function and updates call statistics.
Subclasses should call this function rather than calling system.EvalTimeDerivatives() directly.
|
protected |
Evaluates the derivative function (and updates call statistics).
Subclasses should call this function rather than calling system.EvalTimeDerivatives() directly. This version of this function exists to allow integrators to include AutoDiff'd systems in derivative function evaluations.
double get_accuracy_in_use | ( | ) | const |
Gets the accuracy in use by the integrator.
This number may differ from the target accuracy if, for example, the user has requested an accuracy not attainable or not recommended for the particular integrator.
const T& get_actual_initial_step_size_taken | ( | ) | const |
The actual size of the successful first step.
const Context<T>& get_context | ( | ) | const |
Returns a const reference to the internally-maintained Context holding the most recent state in the trajectory.
This is suitable for publishing or extracting information about this trajectory step.
const trajectories::PiecewisePolynomial<T>* get_dense_output | ( | ) | const |
Returns a const pointer to the integrator's current PiecewisePolynomial instance, holding a representation of the continuous state trajectory since the last StartDenseIntegration() call.
This is suitable to query the integrator's current dense output, if any (may be nullptr).
const ContinuousState<T>* get_error_estimate | ( | ) | const |
Gets the error estimate (used only for integrators that support error estimation).
If the integrator does not support error estimation, nullptr is returned.
|
pure virtual |
Derived classes must override this function to return the order of the asymptotic term in the integrator's error estimate.
An error estimator approximates the truncation error in an integrator's solution. That truncation error e(.) is approximated by a Taylor Series expansion in the neighborhood around t:
e(t+h) ≈ e(t) + he(t) + he'(t) + ½h²e''(t) + ... ≈ e(t) + he(t) + he'(t) + ½h²e''(t) + O(h³)
where we have replaced the "..." with the asymptotic error of all terms truncated from the series.
Implementions should return the order of the asymptotic term in the Taylor Series expansion around the expression for the error. For an integrator that propagates a second-order solution and provides an estimate of the error using an embedded first-order method, this method should return "2", as can be seen in the derivation below, using y* as the true solution:
y̅ = y* + O(h³) [second order solution] ŷ = y* + O(h²) [embedded first-order method] e = (y̅ - ŷ) = O(h²)
If the integrator does not provide an error estimate, the derived class implementation should return 0.
Implemented in VelocityImplicitEulerIntegrator< T >, ImplicitEulerIntegrator< T >, SemiExplicitEulerIntegrator< T >, RadauIntegrator< T, num_stages >, RungeKutta3Integrator< T >, RungeKutta5Integrator< T >, BogackiShampine3Integrator< T >, ExplicitEulerIntegrator< T >, and RungeKutta2Integrator< T >.
bool get_fixed_step_mode | ( | ) | const |
Gets whether an integrator is running in fixed step mode.
If the integrator does not support error estimation, this function will always return true
.
const Eigen::VectorXd& get_generalized_state_weight_vector | ( | ) | const |
Gets the weighting vector (equivalent to a diagonal matrix) applied to weighting both generalized coordinate and velocity state variable errors, as described in the group documentation.
Only used for integrators that support error estimation.
const T& get_ideal_next_step_size | ( | ) | const |
Return the step size the integrator would like to take next, based primarily on the integrator's accuracy prediction.
This value will not be computed for integrators that do not support error estimation and NaN will be returned.
const T& get_initial_step_size_target | ( | ) | const |
Gets the target size of the first integration step.
You can find out what step size was actually used for the first integration step with get_actual_initial_step_size_taken()
.
const T& get_largest_step_size_taken | ( | ) | const |
The size of the largest step taken since the last Initialize() or ResetStatistics() call.
const T& get_maximum_step_size | ( | ) | const |
Gets the maximum step size that may be taken by this integrator.
This is a soft maximum: the integrator may stretch it by as much as 1% to hit a discrete event.
const Eigen::VectorXd& get_misc_state_weight_vector | ( | ) | const |
Gets the weighting vector (equivalent to a diagonal matrix) for weighting errors in miscellaneous continuous state variables z
.
Only used for integrators that support error estimation.
Context<T>* get_mutable_context | ( | ) |
Returns a mutable pointer to the internally-maintained Context holding the most recent state in the trajectory.
|
protected |
Returns a mutable pointer to the internally-maintained PiecewisePolynomial instance, holding a representation of the continuous state trajectory since the last time StartDenseIntegration() was called.
This is useful for derived classes to update the integrator's current dense output, if any (may be nullptr).
|
protected |
Gets an error estimate of the state variables recorded by the last call to StepOnceFixedSize().
If the integrator does not support error estimation, this function will return nullptr.
Eigen::VectorBlock<Eigen::VectorXd> get_mutable_generalized_state_weight_vector | ( | ) |
Gets a mutable weighting vector (equivalent to a diagonal matrix) applied to weighting both generalized coordinate and velocity state variable errors, as described in the group documentation.
Only used for integrators that support error estimation. Returns a VectorBlock to make the values mutable without permitting changing the size of the vector. Requires re-initializing the integrator after calling this method; if Initialize() is not called afterward, an exception will be thrown when attempting to call IntegrateNoFurtherThanTime(). If the caller sets one of the entries to a negative value, an exception will be thrown when the integrator is initialized.
Eigen::VectorBlock<Eigen::VectorXd> get_mutable_misc_state_weight_vector | ( | ) |
Gets a mutable weighting vector (equivalent to a diagonal matrix) for weighting errors in miscellaneous continuous state variables z
.
Only used for integrators that support error estimation. Returns a VectorBlock to make the values mutable without permitting changing the size of the vector. Requires re-initializing the integrator after calling this method. If Initialize() is not called afterward, an exception will be thrown when attempting to call IntegrateNoFurtherThanTime(). If the caller sets one of the entries to a negative value, an exception will be thrown when the integrator is initialized.
int64_t get_num_derivative_evaluations | ( | ) | const |
Returns the number of ODE function evaluations (calls to CalcTimeDerivatives()) since the last call to ResetStatistics() or Initialize().
This count includes all such calls including (1) those necessary to compute Jacobian matrices; (2) those used in rejected integrated steps (for, e.g., purposes of error control); (3) those used strictly for integrator error estimation; and (4) calls that exhibit little cost (due to results being cached).
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().
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 failures) since the last call to ResetStatistics() or Initialize().
int64_t get_num_steps_taken | ( | ) | const |
The number of integration steps taken since the last Initialize() or ResetStatistics() call.
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 solving the necessary nonlinear system of equations).
const T& get_previous_integration_step_size | ( | ) | const |
Gets the size of the last (previous) integration step.
If no integration steps have been taken, value will be NaN.
const T& get_requested_minimum_step_size | ( | ) | const |
Gets the requested minimum step size h_min
for this integrator.
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 the last Initialize() or ResetStatistics() call.
This value will be NaN for integrators without error estimation.
double get_stretch_factor | ( | ) | const |
Gets the stretch factor (> 1), which is multiplied by the maximum (typically user-designated) integration step size to obtain the amount that the integrator is able to stretch the maximum time step toward hitting an upcoming publish or update event in IntegrateNoFurtherThanTime().
const System<T>& get_system | ( | ) | const |
Gets a constant reference to the system that is being integrated (and was provided to the constructor of the integrator).
double get_target_accuracy | ( | ) | const |
Gets the target accuracy.
bool get_throw_on_minimum_step_size_violation | ( | ) | const |
Reports the current setting of the throw_on_minimum_step_size_violation flag.
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's context.
See this section for more detail.
void Initialize | ( | ) |
An integrator must be initialized before being used.
The pointer to the context must be set before Initialize() is called (or an std::exception will be thrown). If Initialize() is not called, an exception will be thrown when attempting to call IntegrateNoFurtherThanTime(). To reinitialize the integrator, Reset() should be called followed by Initialize().
std::exception | If the context has not been set or a user-set parameter has been set illogically (i.e., one of the weighting matrix coefficients is set to a negative value- this check is only performed for integrators that support error estimation; the maximum step size is smaller than the minimum step size; the requested initial step size is outside of the interval [minimum step size, maximum step size]). |
StepResult IntegrateNoFurtherThanTime | ( | const T & | publish_time, |
const T & | update_time, | ||
const T & | boundary_time | ||
) |
(Internal use only) Integrates the system forward in time by a single step with step size subject to integration error tolerances (assuming that the integrator supports error estimation).
The integrator must already have been initialized or an exception will be thrown. The context will be integrated to a time that will never exceed the minimum of publish_time
, update_time
, and the current time plus 1.01 * get_maximum_step_size()
.
publish_time | The present or future time (exception will be thrown if this is not the case) at which the next publish will occur. |
update_time | The present or future time (exception will be thrown if this is not the case) at which the next update will occur. |
boundary_time | The present or future time (exception will be thrown if this is not the case) marking the end of the user-designated simulated interval. |
std::exception | If the integrator has not been initialized or one of publish_time, update_time, or boundary_time is in the past. |
min(publish_time, update_time, boundary_time)
. Simulator::AdvanceTo()
. In other circumstances, users will typically call IntegratorBase::IntegrateWithMultipleStepsToTime()
.This method at a glance:
void IntegrateWithMultipleStepsToTime | ( | const T & | t_final | ) |
Stepping function for integrators operating outside of Simulator that advances the continuous state exactly to t_final
.
This method is designed for integrator users that do not wish to consider publishing or discontinuous, mid-interval updates. This method will step the integrator multiple times, as necessary, to attain requested error tolerances and to ensure the integrator converges.
Simulator::AdvanceTo()
in place of this function (which was created for off-simulation purposes), generally. t_final | The current or future time to integrate to. |
std::exception | If the integrator has not been initialized or t_final is in the past. |
t_final
in a single step.This method at a glance:
t_final
bool IntegrateWithSingleFixedStepToTime | ( | const T & | t_target | ) |
Stepping function for integrators operating outside of Simulator that advances the continuous state using a single step to t_target
.
This method is designed for integrator users that do not wish to consider publishing or discontinuous, mid-interval updates. One such example application is that of direct transcription for trajectory optimization, for which the integration process should be consistent: it should execute the same sequence of arithmetic operations for all values of the nonlinear programming variables. In keeping with the naming semantics of this function, error controlled integration is not supported (though error estimates will be computed for integrators that support that feature), which is a minimal requirement for "consistency".
Simulator::AdvanceTo()
in place of this function (which was created for off-simulation purposes), generally. t_target | The current or future time to integrate to. |
std::exception | If the integrator has not been initialized or t_target is in the past or the integrator is not operating in fixed step mode. |
t_target
. true
if the integrator was able to take a single fixed step to t_target
.This method at a glance:
bool is_initialized | ( | ) | const |
Indicates whether the integrator has been initialized.
|
delete |
|
delete |
void request_initial_step_size_target | ( | const T & | step_size | ) |
Request that the first attempted integration step have a particular size.
If no request is made, the integrator will estimate a suitable size for the initial step attempt. If the integrator does not support error control, this method will throw a std::exception (call supports_error_estimation() to verify before calling this method). For variable-step integration, the initial target will be treated as a maximum step size subject to accuracy requirements and event occurrences. You can find out what size actually worked with get_actual_initial_step_size_taken()
.
std::exception | If the integrator does not support error estimation. |
void Reset | ( | ) |
Resets the integrator to initial values, i.e., default construction values.
void reset_context | ( | Context< T > * | context | ) |
Replace the pointer to the internally-maintained Context with a different one.
This is useful for supplying a new set of initial conditions or wiping out the current context (by passing in a null pointer). You should invoke Initialize() after replacing the Context unless the context is null.
context | The pointer to the new context or nullptr to wipe out the current context without replacing it with another. |
void ResetStatistics | ( | ) |
Forget accumulated statistics.
These are reset to the values they have post construction or immediately after Initialize()
.
|
protected |
Sets the working ("in use") accuracy for this integrator.
The working accuracy may not be equivalent to the target accuracy when the latter is too loose or tight for an integrator's capabilities.
|
protected |
void set_fixed_step_mode | ( | bool | flag | ) |
Sets an integrator with error control to fixed step mode.
If the integrator runs in fixed step mode, it will always take the maximum step size directed (which may be that determined by get_maximum_step_size(), or may be smaller, as directed by, e.g., Simulator for event handling purposes).
std::exception | if integrator does not support error estimation and flag is set to false . |
|
protected |
|
protected |
void set_maximum_step_size | ( | const T & | max_step_size | ) |
Sets the maximum step size that may be taken by this integrator.
This setting should be used if you know the maximum time scale of your problem. The integrator may stretch the maximum step size by as much as 1% to reach a discrete event. For fixed step integrators, all steps will be taken at the maximum step size unless an event would be missed.
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.
No step smaller than this will be taken except under circumstances as described above. This setting will be ignored if it is smaller than the absolute minimum h_floor
also described above. Default value is zero.
min_step_size | a non-negative value. Setting this value to zero will cause the integrator to use a reasonable value instead (see get_working_minimum_step_size()). |
|
protected |
Sets the size of the smallest-step-taken statistic as the result of a controlled integration step adjustment.
void set_target_accuracy | ( | double | accuracy | ) |
Request that the integrator attempt to achieve a particular accuracy for the continuous portions of the simulation.
Otherwise a default accuracy is chosen for you. This may be ignored for fixed-step integration since accuracy control requires variable step sizes. You should call supports_error_estimation() to ensure that the integrator supports this capability before calling this function; if the integrator does not support it, this method will throw an exception.
std::exception | if integrator does not support error estimation. |
void set_throw_on_minimum_step_size_violation | ( | bool | throws | ) |
Sets whether the integrator should throw a std::exception when the integrator's step size selection algorithm determines that it must take a step smaller than the minimum step size (for, e.g., purposes of error control).
Default is true
. If false
, the integrator will advance time and state using the minimum specified step size in such situations. See this section for more detail.
void StartDenseIntegration | ( | ) |
Starts dense integration, allocating a new dense output for this integrator to use.
std::exception | if any of the preconditions is not met. |
|
protected |
Default code for advancing the continuous state of the system by a single step of h_max
(or smaller, depending on error control).
This particular function is designed to be called directly by an error estimating integrator's DoStep() method to effect error-controlled integration. The integrator can effect error controlled integration without calling this method, if the implementer so chooses, but this default method is expected to function well in most circumstances.
[in] | h_max | The maximum step size to be taken. The integrator may take a smaller step than specified to satisfy accuracy requirements, to resolve integrator convergence problems, or to respect the integrator's maximum step size. |
std::exception | if integrator does not support error estimation. |
true
if the full step of size h_max
is taken and false
otherwise (i.e., a smaller step than h_max
was taken). std::unique_ptr<trajectories::PiecewisePolynomial<T> > StopDenseIntegration | ( | ) |
Stops dense integration, yielding ownership of the current dense output to the caller.
std::exception | if any of the preconditions is not met. |
|
pure virtual |
Derived classes must override this function to indicate whether the integrator supports error estimation.
Without error estimation, the target accuracy setting (see accuracy settings) will be unused.
Implemented in VelocityImplicitEulerIntegrator< T >, ImplicitEulerIntegrator< T >, SemiExplicitEulerIntegrator< T >, RadauIntegrator< T, num_stages >, RungeKutta3Integrator< T >, RungeKutta5Integrator< T >, BogackiShampine3Integrator< T >, ExplicitEulerIntegrator< T >, and RungeKutta2Integrator< T >.