The optimal control solution object: structure and usage

In this manual, we'll first recall the main functionalities you can use when working with a solution of an optimal control problem (SOL). This includes essential operations like:

  • Plotting a SOL: How to plot the optimal solution for your defined problem.
  • Printing a SOL: How to display a summary of your solution.

After covering these core functionalities, we'll delve into the structure of a SOL. Since a SOL is structured as a Solution struct, we'll first explain how to access its underlying attributes. Following this, we'll shift our focus to the simple properties inherent to a SOL.


Content


Main functionalities

Let's define a basic optimal control problem.

using OptimalControl

t0 = 0
tf = 1
x0 = [-1, 0]

ocp = @def begin
    t ∈ [ t0, tf ], time
    x = (q, v) ∈ R², state
    u ∈ R, control
    x(t0) == x0
    x(tf) == [0, 0]
    ẋ(t)  == [v(t), u(t)]
    0.5∫( u(t)^2 ) → min
end
nothing # hide
<< @example-block not executed in draft mode >>

We can now solve the problem (for more details, visit the solve manual):

using NLPModelsIpopt
sol = solve(ocp)
nothing # hide
<< @example-block not executed in draft mode >>
Note

You can export (or save) the solution in a Julia .jld2 data file and reload it later, and also export a discretised version of the solution in a more portable JSON format. Note that the optimal control problem is needed when loading a solution.

See the two functions:

To print sol, simply:

sol
<< @example-block not executed in draft mode >>

For complementary information, you can plot the solution:

using Plots
plot(sol)
<< @example-block not executed in draft mode >>
Note

For more details about plotting a solution, visit the plot manual.

Solution struct

The solution sol is a Solution struct.

CTModels.SolutionType
struct Solution{TimeGridModelType<:CTModels.AbstractTimeGridModel, TimesModelType<:CTModels.AbstractTimesModel, StateModelType<:CTModels.AbstractStateModel, ControlModelType<:CTModels.AbstractControlModel, VariableModelType<:CTModels.AbstractVariableModel, CostateModelType<:Function, ObjectiveValueType<:Real, DualModelType<:CTModels.AbstractDualModel, SolverInfosType<:CTModels.AbstractSolverInfos, ModelType<:CTModels.AbstractModel} <: CTModels.AbstractSolution

Fields

  • time_grid::CTModels.AbstractTimeGridModel

  • times::CTModels.AbstractTimesModel

  • state::CTModels.AbstractStateModel

  • control::CTModels.AbstractControlModel

  • variable::CTModels.AbstractVariableModel

  • costate::Function

  • objective::Real

  • dual::CTModels.AbstractDualModel

  • solver_infos::CTModels.AbstractSolverInfos

  • model::CTModels.AbstractModel

source

Each field can be accessed directly (ocp.times, etc) but we recommend to use the sophisticated getters we provide: the state(sol::Solution) method does not return sol.state but a function of time that can be called at any time, not only on the grid time_grid.

0.25 ∈ time_grid(sol)
<< @example-block not executed in draft mode >>
x = state(sol)
x(0.25)
<< @example-block not executed in draft mode >>

Attributes and properties

State, costate, control, variable and objective value

You can access the values of the state, costate, control and variable by eponymous functions. The returned values are functions of time for the state, costate and control and a scalar or a vector for the variable.

t = 0.25
x = state(sol)
p = costate(sol)
u = control(sol)
nothing # hide
<< @example-block not executed in draft mode >>

Since the state is of dimension 2, evaluating x(t) returns a vector:

x(t)
<< @example-block not executed in draft mode >>

It is the same for the costate:

p(t)
<< @example-block not executed in draft mode >>

But the control is one-dimensional:

u(t)
<< @example-block not executed in draft mode >>

There is no variable, hence, an empty vector is returned:

v = variable(sol)
<< @example-block not executed in draft mode >>

The objective value is accessed by:

objective(sol)
<< @example-block not executed in draft mode >>

Infos from the solver

The problem ocp is solved via a direct method (see solve manual for details). The solver stores data in sol, including the success of the optimization, the iteration count, the time grid used for discretisation, and other specific details within the solver_infos field.

time_grid(sol)
<< @example-block not executed in draft mode >>
constraints_violation(sol)
<< @example-block not executed in draft mode >>
infos(sol)
<< @example-block not executed in draft mode >>
iterations(sol)
<< @example-block not executed in draft mode >>
message(sol)
<< @example-block not executed in draft mode >>
status(sol)
<< @example-block not executed in draft mode >>
successful(sol)
<< @example-block not executed in draft mode >>

Dual variables

You can retrieved dual variables (or Lagrange multipliers) associated to labelled constraint. To illustrate this, we define a problem with constraints:

ocp = @def begin

    tf ∈ R,             variable
    t ∈ [0, tf],        time
    x = (q, v) ∈ R²,    state
    u ∈ R,              control

    tf ≥ 0,             (eq_tf)
    -1 ≤ u(t) ≤ 1,      (eq_u)
    v(t) ≤ 0.75,        (eq_v)

    x(0)  == [-1, 0],   (eq_x0)
    q(tf) == 0
    v(tf) == 0

    ẋ(t) == [v(t), u(t)]

    tf → min

end
sol = solve(ocp; display=false)
nothing # hide
<< @example-block not executed in draft mode >>

Dual variables corresponding to variable and boundary constraints are given as scalar or vectors.

dual(sol, ocp, :eq_tf)
<< @example-block not executed in draft mode >>
dual(sol, ocp, :eq_x0)
<< @example-block not executed in draft mode >>

The other type of constraints are associated to dual variables given as functions of time.

μ_u = dual(sol, ocp, :eq_u)
plot(time_grid(sol), μ_u)
<< @example-block not executed in draft mode >>
μ_v = dual(sol, ocp, :eq_v)
plot(time_grid(sol), μ_v)
<< @example-block not executed in draft mode >>