The solve function

In this manual, we explain the solve function from OptimalControl.jl package.

Basic usage

Let us 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:

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

Note that we must import NLPModelsIpopt.jl before calling solve. This is because the default method uses a direct approach, which transforms the optimal control problem into a nonlinear program (NLP) of the form:

\[\text{minimize}\quad F(y), \quad\text{subject to the constraints}\quad g(y) \le 0, \quad h(y) = 0. \]

Warning

Calling solve without loading a NLP solver package first will notify the user:

julia> solve(ocp)
ERROR: ExtensionError. Please make: julia> using NLPModelsIpopt

Resolution methods and algorithms

OptimalControl.jl offers a list of methods. To get it, simply call methods.

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

Each line is a method, with priority going from top to bottom. This means that

solve(ocp)

is equivalent to

solve(ocp, :direct, :adnlp, :ipopt)
  1. The first symbol refers to the general class of method. The only possible value is:
    • :direct: currently, only the so-called direct approach is implemented. Direct methods discretise the original optimal control problem and solve the resulting NLP. In this case, the main solve method redirects to CTDirect.solve.
  2. The second symbol refers to the NLP modeler. The possible values are:
  3. The third symbol specifies the NLP solver. Possible values are:
    • :ipopt: calls NLPModelsIpopt.ipopt to solve the NLP problem.
    • :madnlp: creates a MadNLP.MadNLPSolver instance from the NLP problem and solve it. MadNLP.jl is an open-source solver in Julia implementing a filter line-search interior-point algorithm like Ipopt.
    • :knitro: uses the Knitro solver (license required).
Warning

When using :exa for more performance (in particular to solve on GPU), there are limitations on the syntax:

  • dynamics must be declared coordinate by coordinate (not globally as a vector valued expression)
  • nonlinear constraints (boundary, variable, control, state, mixed ones, see Constraints must also be scalar expressions (linear constraints aka. ranges, on the other hand, can be vectors)
  • all expressions must only involve algebraic operations that are known to ExaModels (check the documentation), although one can provide additional user defined functions through registration (check ExaModels API)

For instance, let us try MadNLP solver with ExaModel modeller.

using MadNLP

ocp = @def begin
    t ∈ [ t0, tf ], time
    x = (q, v) ∈ R², state
    u ∈ R, control
    x(t0) == x0
    x(tf) == [0, 0]
    ∂(q)(t) == v(t)
    ∂(v)(t) == u(t)
    0.5∫( u(t)^2 ) → min
end

solve(ocp, :exa, :madnlp; scheme=:trapeze)
nothing # hide
<< @example-block not executed in draft mode >>

Note that you can provide a partial description. If multiple full descriptions contain it, priority is given to the first one in the list. Hence, all of the following calls are equivalent:

solve(ocp)
solve(ocp, :direct                )
solve(ocp,          :adnlp        )
solve(ocp,                  :ipopt)
solve(ocp, :direct, :adnlp        )
solve(ocp, :direct,         :ipopt)
solve(ocp, :direct, :adnlp, :ipopt)

Direct method and options

The main options for the direct method, with their [default] values, are:

  • display ([true], false): setting display = false disables output.
  • init: information for the initial guess. It can be given as numerical values, functions, or an existing solution. See how to set an initial guess.
  • grid_size ([250]): number of time steps in the (uniform) time discretization grid. More precisely, if N = grid_size and the initial and final times are t0 and tf, then the step length Δt = (tf - t0) / N.
  • time_grid ([nothing]): explicit time grid (can be non-uniform). If time_grid = nothing, a uniform grid of length grid_size is used.
  • scheme (:trapeze, [:midpoint], :euler, :euler_implicit, :gauss_legendre_2, :gauss_legendre_3): the discretisation scheme to transform the dynamics into nonlinear equations. See the discretization method tutorial for more details.
  • adnlp_backend ([:optimized], :manual, :default): backend used for automatic differentiation to create the ADNLPModels.ADNLPModel.

For advanced usage, see:

Note

The main solve method from OptimalControl.jl simply redirects to CTDirect.solve in that case.

NLP solvers specific options

In addition to these options, any remaining keyword arguments passed to solve are forwarded to the NLP solver.

Warning

The option names and accepted values depend on the chosen solver. For example, in Ipopt, print_level expects an integer, whereas in MadNLP it must be a MadNLP.LogLevels value (valid options: MadNLP.{TRACE, DEBUG, INFO, NOTICE, WARN, ERROR}). Moreover, some options are solver-specific: for instance, mu_strategy exists in Ipopt but not in MadNLP.

Please refer to the Ipopt options list and the NLPModelsIpopt.jl documentation.

sol = solve(ocp; max_iter=0, tol=1e-6, print_level=0)
iterations(sol)
<< @example-block not executed in draft mode >>

Likewise, see the MadNLP.jl options and the MadNLP.jl documentation.

sol = solve(ocp, :madnlp; max_iter=0, tol=1e-6, print_level=MadNLP.ERROR)
iterations(sol)
<< @example-block not executed in draft mode >>