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 Implicit integration rule with the following methods

RobotDynamics.discrete_dynamicsFunction

Compute the discretized dynamics of model using explicit integration scheme Q<:QuadratureRule.

Methods:

x′ = discrete_dynamics(model, model, z)  # uses RK3 as the default integration scheme
x′ = discrete_dynamics(Q, model, x, u, t, dt)
x′ = discrete_dynamics(Q, model, z::KnotPoint)

The default integration scheme is stored in TrajectoryOptimization.DEFAULT_Q

RobotDynamics.discrete_jacobian!Function
∇f = discrete_jacobian!(Q, ∇f, model, z::AbstractKnotPoint)

Compute the n × (n+m) discrete dynamics Jacobian ∇f of model using explicit integration scheme Q<:QuadratureRule.

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.RK2Type

Second-order Runge-Kutta method with zero-order-old on the controls (i.e. midpoint)

RobotDynamics.RK3Type

Second-order Runge-Kutta method with zero-order-old on the controls

RobotDynamics.RK4Type

Fourth-order Runge-Kutta method with zero-order-old on the controls

RobotDynamics.ImplicitType

Integration rules of the form x′ = f(x,u,x′,u′), where x′,u′ are the states and controls at the next time step.

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!).