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. This includes essential operations like:
- Plotting a solution: How to plot the optimal solution for your defined problem.
- Printing a solution: How to display a summary of your solution.
After covering these core functionalities, we'll delve into the structure of a solution. Since a solution is structured as a OptimalControl.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 solution.
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 >>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 >>For more details about plotting a solution, visit the plot manual.
Solution struct
The solution sol is a OptimalControl.Solution struct.
CTModels.OCP.Solution — Type
struct Solution{TimeGridModelType<:CTModels.OCP.AbstractTimeGridModel, TimesModelType<:CTModels.OCP.AbstractTimesModel, StateModelType<:CTModels.OCP.AbstractStateModel, ControlModelType<:CTModels.OCP.AbstractControlModel, VariableModelType<:CTModels.OCP.AbstractVariableModel, ModelType<:CTModels.OCP.AbstractModel, CostateModelType<:Function, ObjectiveValueType<:Real, DualModelType<:CTModels.OCP.AbstractDualModel, SolverInfosType<:CTModels.OCP.AbstractSolverInfos} <: CTModels.OCP.AbstractSolutionComplete solution of an optimal control problem.
Stores the optimal state, control, and costate trajectories, the optimisation variable value, objective value, dual variables, and solver information.
Fields
time_grid::TimeGridModelType: Discretised time points.times::TimesModelType: Initial and final time specification.state::StateModelType: State trajectoryt -> x(t)with metadata.control::ControlModelType: Control trajectoryt -> u(t)with metadata.variable::VariableModelType: Optimisation variable value with metadata.model::ModelType: Reference to the optimal control problem model.costate::CostateModelType: Costate (adjoint) trajectoryt -> p(t).objective::ObjectiveValueType: Optimal objective value.dual::DualModelType: Dual variables for all constraints.solver_infos::SolverInfosType: Solver statistics and status.
Example
julia> using CTModels
julia> # Solutions are typically returned by solvers
julia> sol = solve(ocp, ...) # Returns a Solution
julia> CTModels.objective(sol)Each field can be accessed directly (sol.state, 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 >>You can also retrieve the original optimal control problem from the solution:
model(sol) # returns the original OCP model<< @example-block not executed in draft mode >>Trajectories
The trajectory component provides access to the state, control, variable, and costate trajectories.
State trajectory
Get the state trajectory as a function of time:
x = state(sol) # returns a function of time<< @example-block not executed in draft mode >>Evaluate the state at any time (not just grid points):
t = 0.25
x(t) # returns state vector at t=0.25<< @example-block not executed in draft mode >>The state function can be evaluated at any time within the problem horizon, even if it's not a discretization grid point:
0.25 ∈ time_grid(sol) # false: not a grid point<< @example-block not executed in draft mode >>x(0.25) # still works: interpolated value<< @example-block not executed in draft mode >>Control trajectory
Get the control trajectory as a function of time:
u = control(sol) # returns a function of time<< @example-block not executed in draft mode >>u(t) # returns control value at t<< @example-block not executed in draft mode >>Variable values
Get the optimization variable values:
v = variable(sol) # returns an empty vector if no variable<< @example-block not executed in draft mode >>Costate trajectory
Get the costate (adjoint) trajectory as a function of time:
p = costate(sol) # returns a function of time<< @example-block not executed in draft mode >>p(t) # returns costate vector at t<< @example-block not executed in draft mode >>Time information
Get time-related information from the solution:
time_grid(sol) # returns the discretization time grid<< @example-block not executed in draft mode >>times(sol) # returns the TimesModel struct containing time information<< @example-block not executed in draft mode >>Unified vs. multiple grids:
With a standard collocation method, there is a single time grid that can be retrieved via time_grid(sol). The solution internally uses a UnifiedTimeGridModel for memory efficiency.
For discretization methods that use multiple grids (one per component), the solution uses a MultipleTimeGridModel. In this case, you must specify which component's grid you want:
time_grid(sol, :state)— state trajectory and state box constraint dualstime_grid(sol, :control)— control trajectory and control box constraint dualstime_grid(sol, :costate)— costate trajectory (maps to:stategrid)time_grid(sol, :path)— path constraint duals
Aliases are accepted: :costate/:costates map to :state, :dual/:duals map to :path, and plural forms (:states, :controls) are also valid.
All grids must be strictly increasing, finite, and non-empty.
Trajectories (state, control, costate, path_constraints_dual) can be provided either as matrices (rows = time points, columns = components) or as functions t -> vector for interpolated or analytical data.
Summary table
| Method | Returns | Description |
|---|---|---|
state(sol) | Function | State trajectory x(t) |
control(sol) | Function | Control trajectory u(t) |
variable(sol) | Vector | Variable values |
costate(sol) | Function | Costate trajectory p(t) |
time_grid(sol) | Vector{Float64} | Discretization time grid |
times(sol) | TimesModel | TimesModel struct containing time information |
Objective
The objective component provides access to the objective value.
Objective value
Get the optimal objective value:
objective(sol) # returns the objective value<< @example-block not executed in draft mode >>Summary table
| Method | Returns | Description |
|---|---|---|
objective(sol) | Float64 | Objective value |
Dual variables
The dual variables (Lagrange multipliers) provide sensitivity information about constraints.
To illustrate dual variables, we define a problem with various 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 of labeled constraints
Get the dual variable for a specific labeled constraint:
dual(sol, ocp, :eq_tf) # dual for variable constraint<< @example-block not executed in draft mode >>dual(sol, ocp, :eq_x0) # dual for boundary constraint<< @example-block not executed in draft mode >>For path constraints, the dual is a function 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 >>In all cases, dual(sol, ocp, :label) returns a signed multiplier μ (scalar, vector, or function of time, depending on the constraint type). The sign convention is, component-wise:
μ > 0⇒ the lower-side constraint is active (e.g.lb ≤ ...),μ < 0⇒ the upper-side constraint is active (e.g.... ≤ ub),μ = 0⇒ the constraint is inactive (or the component is never constrained).
For nonlinear path and boundary constraints, the solver already returns a signed multiplier natively; CTModels simply forwards it via path_constraints_dual(sol) / boundary_constraints_dual(sol).
For box constraints (state/control/variable components), the solver stores lower- and upper-bound multipliers separately as non-negative quantities. CTModels combines them into the signed multiplier explicitly:
μ = μ_lb − μ_ubcomputed per targeted primal component. This is the value returned by dual(sol, ocp, :label) for a box-constraint label.
Rationale for box constraints: after intersection of duplicate box declarations, the solver only sees a single effective bound per component, hence a single signed multiplier per component. If several labels target the same component, each label returns the same per-component multiplier (via the aliases mechanism; see the OCP manual and Duplicate box constraints).
The raw non-negative *_lb_dual(sol) / *_ub_dual(sol) accessors (see below) remain available if the unsigned components are needed separately.
Box constraint duals
Get dual variables for box constraints on state, control, and variable:
state_constraints_lb_dual(sol) # lower bound duals for state<< @example-block not executed in draft mode >>state_constraints_ub_dual(sol) # upper bound duals for state<< @example-block not executed in draft mode >>control_constraints_lb_dual(sol) # lower bound duals for control<< @example-block not executed in draft mode >>control_constraints_ub_dual(sol) # upper bound duals for control<< @example-block not executed in draft mode >>variable_constraints_lb_dual(sol) # lower bound duals for variable<< @example-block not executed in draft mode >>variable_constraints_ub_dual(sol) # upper bound duals for variable<< @example-block not executed in draft mode >>Box constraint dual dimensions
Get the dimensions of the box-constraint dual vectors (one entry per primal component that is box-constrained):
dim_dual_state_constraints_box(sol) # dimension of state box constraint duals<< @example-block not executed in draft mode >>dim_dual_control_constraints_box(sol) # dimension of control box constraint duals<< @example-block not executed in draft mode >>dim_dual_variable_constraints_box(sol) # dimension of variable box constraint duals<< @example-block not executed in draft mode >>Nonlinear constraint duals
Get dual variables for nonlinear path and boundary constraints:
path_constraints_dual(sol) # duals for nonlinear path constraints<< @example-block not executed in draft mode >>boundary_constraints_dual(sol) # duals for nonlinear boundary constraints<< @example-block not executed in draft mode >>Get the dimensions of the nonlinear constraints:
dim_path_constraints_nl(sol) # number of nonlinear path constraints<< @example-block not executed in draft mode >>dim_boundary_constraints_nl(sol) # number of nonlinear boundary constraints<< @example-block not executed in draft mode >>Summary table
| Method | Returns | Description |
|---|---|---|
dual(sol, ocp, label) | Real or Function | Signed dual for labeled constraint |
state_constraints_lb_dual(sol) | Dual values | State lower bound duals (non-negative) |
state_constraints_ub_dual(sol) | Dual values | State upper bound duals (non-negative) |
control_constraints_lb_dual(sol) | Dual values | Control lower bound duals (non-negative) |
control_constraints_ub_dual(sol) | Dual values | Control upper bound duals (non-negative) |
variable_constraints_lb_dual(sol) | Dual values | Variable lower bound duals (non-negative) |
variable_constraints_ub_dual(sol) | Dual values | Variable upper bound duals (non-negative) |
dim_dual_state_constraints_box(sol) | Int | Dimension of state box constraint duals |
dim_dual_control_constraints_box(sol) | Int | Dimension of control box constraint duals |
dim_dual_variable_constraints_box(sol) | Int | Dimension of variable box constraint duals |
path_constraints_dual(sol) | Dual values | Nonlinear path constraint duals (signed) |
boundary_constraints_dual(sol) | Dual values | Nonlinear boundary constraint duals (signed) |
dim_path_constraints_nl(sol) | Int | Number of nonlinear path constraints |
dim_boundary_constraints_nl(sol) | Int | Number of nonlinear boundary constraints |
Solution metadata
The solution metadata provides information about the solver performance and status.
Solver status
Check if the solution was successful:
successful(sol) # returns true if solver succeeded<< @example-block not executed in draft mode >>Get the solver status symbol:
status(sol) # returns solver status (e.g., :first_order)<< @example-block not executed in draft mode >>Get the solver message:
message(sol) # returns solver message string<< @example-block not executed in draft mode >>Iteration count
Get the number of solver iterations:
iterations(sol) # returns iteration count<< @example-block not executed in draft mode >>Constraints violation
Get the maximum constraint violation:
constraints_violation(sol) # returns max violation<< @example-block not executed in draft mode >>Additional solver information
Get additional solver-specific information:
infos(sol) # returns dictionary of solver info<< @example-block not executed in draft mode >>Summary table
| Method | Returns | Description |
|---|---|---|
successful(sol) | Bool | True if solver succeeded |
status(sol) | Symbol | Solver status |
message(sol) | String | Solver message |
iterations(sol) | Int | Number of iterations |
constraints_violation(sol) | Float64 | Maximum constraint violation |
infos(sol) | Dict | Additional solver information |