Discretization
This page gives details on the methods for evaluating discretized dynamics, as well as instructions on how to define a custom integration method.
Model Discretization
With a model defined, we can compute the discrete dynamics and discrete dynamics Jacobians for an explicit integration rule with the following methods
x′ = discrete_dynamics(Q, model, x, u, t, dt)
x′ = discrete_dynamics(Q, model, z::KnotPoint)
discrete_jacobian!(Q, ∇f, model, z::AbstractKnotPoint)
For more information, see the docstrings for discrete_dynamics
and discrete_jacobian!
.
Integration Schemes
RobotDynamics.jl has already defined a handful of integration schemes for computing discrete dynamics. The integration schemes are specified as abstract types, so that methods can efficiently dispatch based on the integration scheme selected. Here is the current set of implemented types:
RobotDynamics.QuadratureRule
— TypeIntegration rule for approximating the continuous integrals for the equations of motion
RobotDynamics.Explicit
— TypeIntegration rules of the form x′ = f(x,u), where x′ is the next state
RobotDynamics.RK2
— TypeSecond-order Runge-Kutta method with zero-order-old on the controls (i.e. midpoint)
RobotDynamics.RK3
— TypeSecond-order Runge-Kutta method with zero-order-old on the controls
RobotDynamics.RK4
— TypeFourth-order Runge-Kutta method with zero-order-old on the controls
RobotDynamics.Exponential
— TypeExponential integration for linear systems with ZOH on controls.
RobotDynamics.PassThrough
— TypeIntegration type for systems with user defined discrete dynamics.
RobotDynamics.Implicit
— TypeIntegration rules of the form x′ = f(x,u,x′,u′), where x′,u′ are the states and controls at the next time step.
RobotDynamics.HermiteSimpson
— TypeThird-order Runge-Kutta method with first-order-hold on the controls
Defining a New Integration Scheme
Explicit Methods
Explicit integration schemes are understandably simpler, since the output is not a function of itself, as is the case with implict schemes. As such, as a minimum, the user only needs to define the following method for a new rule MyQ
:
abstract type MyQ <: RobotDynamics.Explicit end
x′ = discrete_dynamics(::Type{MyQ}, model::AbstractModel, x, u, t, dt)
which will make calls to the continuous-time dynamics function dynamics(model, x, u, t)
.
Below is an example of the default integration method RK3
, a third-order Runge-Kutta method:
function discrete_dynamics(::Type{RK3}, model::AbstractModel,
x::StaticVector, u::StaticVector, t, dt)
k1 = dynamics(model, x, u, t )*dt;
k2 = dynamics(model, x + k1/2, u, t + dt/2)*dt;
k3 = dynamics(model, x - k1 + 2*k2, u, t + dt )*dt;
x + (k1 + 4*k2 + k3)/6
end
Implicit Methods (experimental)
Incorporating implicit integration methods is still under development (great option for someone looking to contribute!).