Internal functions

CTDirect.DOCPType

Struct for discretized optimal control problem DOCP

Contains:

  • a copy of the original OCP
  • data required to link the OCP with the discretized DOCP
source
CTDirect.Gauss_Legendre_1Type

Implicit Midpoint discretization, formulated as a generic IRK (ie Gauss Legendre 1) For testing purpose only, use :midpoint instead (cf midpoint.jl) !

source
CTDirect.DOCP_Hessian_patternMethod
DOCP_Hessian_pattern(
    docp::CTDirect.DOCP{<:CTDirect.GenericIRK}
) -> SparseArrays.SparseMatrixCSC{Bool, Int64}

Build sparsity pattern for Hessian of Lagrangian

source
CTDirect.DOCP_Hessian_patternMethod
DOCP_Hessian_pattern(
    docp::CTDirect.DOCP{CTDirect.Midpoint}
) -> SparseArrays.SparseMatrixCSC{Bool, Int64}

Build sparsity pattern for Hessian of Lagrangian

source
CTDirect.DOCP_Hessian_patternMethod
DOCP_Hessian_pattern(
    docp::CTDirect.DOCP{CTDirect.Trapeze}
) -> SparseArrays.SparseMatrixCSC{Bool, Int64}

Build sparsity pattern for Hessian of Lagrangian

source
CTDirect.DOCP_Hessian_patternMethod
DOCP_Hessian_pattern(
    docp::CTDirect.DOCP{D<:CTDirect.Discretization}
) -> SparseArrays.SparseMatrixCSC{Bool, Int64}

Build sparsity pattern for Hessian of Lagrangian

source
CTDirect.DOCP_Jacobian_patternMethod
DOCP_Jacobian_pattern(
    docp::CTDirect.DOCP{<:CTDirect.GenericIRK}
) -> SparseArrays.SparseMatrixCSC{Bool, Int64}

Build sparsity pattern for Jacobian of constraints

source
CTDirect.DOCP_Jacobian_patternMethod
DOCP_Jacobian_pattern(
    docp::CTDirect.DOCP{CTDirect.Midpoint}
) -> SparseArrays.SparseMatrixCSC{Bool, Int64}

Build sparsity pattern for Jacobian of constraints

source
CTDirect.DOCP_Jacobian_patternMethod
DOCP_Jacobian_pattern(
    docp::CTDirect.DOCP{CTDirect.Trapeze}
) -> SparseArrays.SparseMatrixCSC{Bool, Int64}

Build sparsity pattern for Jacobian of constraints

source
CTDirect.DOCP_Jacobian_patternMethod
DOCP_Jacobian_pattern(
    docp::CTDirect.DOCP{D<:CTDirect.Discretization}
) -> SparseArrays.SparseMatrixCSC{Bool, Int64}

Build sparsity pattern for Jacobian of constraints

source
CTDirect.DOCP_constraints!Method
DOCP_constraints!(c, xu, docp::CTDirect.DOCP) -> Any

Compute the constraints C for the DOCP problem (modeled as LB <= C(X) <= UB).

source
CTDirect.DOCP_initial_guessFunction
DOCP_initial_guess(docp::CTDirect.DOCP) -> Vector{Float64}
DOCP_initial_guess(
    docp::CTDirect.DOCP,
    init::CTModels.Init
) -> Vector{Float64}

Build initial guess for discretized problem

source
CTDirect.IRK_dimsMethod
IRK_dims(
    dim_NLP_steps,
    dim_NLP_x,
    dim_NLP_u,
    dim_NLP_v,
    dim_path_cons,
    dim_boundary_cons,
    stage
) -> NTuple{5, Any}

Return the dimension of the NLP variables and constraints for a generic IRK discretizion, with the control taken constant per step (ie not distinct controls at time stages)

source
CTDirect.SolverInfosMethod
SolverInfos() -> Tuple{Int64, Float64, String, Symbol, Bool}

Retrieve convergence information (Ipopt version)

source
CTDirect.__adnlp_backendMethod
__adnlp_backend() -> Symbol

Used to set the default backend for AD in ADNLPModels. The default value is :optimized.

source
CTDirect.__disc_methodMethod
__disc_method() -> Symbol

Used to set the default discretization method. The default value is trapeze.

source
CTDirect.__displayMethod
__display() -> Bool

Used to set the default display toggle. The default value is true.

source
CTDirect.__ipopt_linear_solverMethod
__ipopt_linear_solver() -> String

Used to set the default value of the linear solver of Ipopt for the direct method. The default value is mumps.

source
CTDirect.__ipopt_mu_strategyMethod
__ipopt_mu_strategy() -> String

Used to set the default value of the μ strategy of Ipopt for the direct method. The default value is adaptive.

source
CTDirect.__ipopt_print_levelMethod
__ipopt_print_level() -> Int64

Used to set the default value of the print level of Ipopt for the direct method. The default value is 5.

source
CTDirect.__knitro_print_levelMethod
__knitro_print_level() -> Int64

Used to set the default value of the print level of Knitro for the direct method. The default value is 3.

source
CTDirect.__madnlp_linear_solverMethod
__madnlp_linear_solver() -> String

Used to set the default value of the linear solver of MadNLP for the direct method. The default value is umfpack.

source
CTDirect.__ocp_initMethod
__ocp_init()

Used to set the default initial guess. The default value is nothing and will correspond to 0.1 for all variables.

source
CTDirect.add_nonzero_block!Method
add_nonzero_block!(M, i_start, i_end, j_start, j_end; sym)

Add block of nonzeros elements to a sparsity pattern Format: boolean matrix (M) or index vectors (Is, Js) Includes a more compact method for single element case Option to add the symmetric block also (eg for Hessian) Note: independent from discretization scheme

source
CTDirect.available_methodsMethod
available_methods(

) -> Tuple{Tuple{Symbol, Symbol}, Tuple{Symbol, Symbol}, Tuple{Symbol, Symbol}}

Return the list of available methods to solve the optimal control problem.

source
CTDirect.build_OCP_solutionMethod
build_OCP_solution(
    docp,
    docp_solution
) -> CTModels.Solution{TimeGridModelType, TimesModelType, StateModelType, ControlModelType, VariableModelType, CostateModelType, Float64, CTModels.DualModel{CTModels.var"#69#81"{Matrix{Float64}}, CTModels.var"#70#82"{Matrix{Float64}}, Vector{Float64}, Vector{Float64}, CTModels.var"#71#83"{Matrix{Float64}, Int64}, CTModels.var"#72#84"{Matrix{Float64}, Int64}, CTModels.var"#73#85"{Matrix{Float64}, Int64}, CTModels.var"#74#86"{Matrix{Float64}, Int64}, Vector{Float64}, Vector{Float64}}, CTModels.SolverInfos{Dict{Symbol, Any}}} where {TimeGridModelType<:CTModels.TimeGridModel, TimesModelType<:CTModels.TimesModel, StateModelType<:Union{CTModels.StateModelSolution{TS} where TS<:CTModels.var"#63#75", CTModels.StateModelSolution{TS} where TS<:CTModels.var"#64#76"}, ControlModelType<:Union{CTModels.ControlModelSolution{TS} where TS<:CTModels.var"#65#77", CTModels.ControlModelSolution{TS} where TS<:CTModels.var"#66#78"}, VariableModelType<:Union{CTModels.VariableModelSolution{Vector{Float64}}, CTModels.VariableModelSolution{Float64}}, CostateModelType<:Union{CTModels.var"#67#79", CTModels.var"#68#80"}}

Build OCP functional solution from DOCP discrete solution (given as a SolverCore.GenericExecutionStats)

source
CTDirect.build_OCP_solutionMethod
build_OCP_solution(docp; primal, dual, mult_LB, mult_UB)

Build OCP functional solution from DOCP discrete solution (given as a SolverCore.GenericExecutionStats)

source
CTDirect.build_boundsMethod
build_bounds(
    dim_var,
    dim_box,
    box_triplet
) -> Tuple{Any, Any}

Build full, ordered sets of bounds for state, control or optimization variables

source
CTDirect.constraints_bounds!Method
constraints_bounds!(
    docp::CTDirect.DOCP
) -> Tuple{Vector{Float64}, Vector{Float64}}

Build upper and lower bounds vectors for the DOCP nonlinear constraints.

source
CTDirect.direct_transcriptionMethod
direct_transcription(
    ocp::CTModels.Model,
    description...;
    grid_size,
    disc_method,
    time_grid,
    init,
    adnlp_backend,
    solver_backend,
    show_time,
    matrix_free
) -> Tuple{CTDirect.DOCP{_A, CTModels.Model{TimesModelType, StateModelType, ControlModelType, VariableModelType, DynamicsModelType, ObjectiveModelType, ConstraintsModelType}} where {_A<:CTDirect.Discretization, TimesModelType<:CTModels.AbstractTimesModel, StateModelType<:CTModels.AbstractStateModel, ControlModelType<:CTModels.AbstractControlModel, VariableModelType<:CTModels.AbstractVariableModel, DynamicsModelType<:Function, ObjectiveModelType<:CTModels.AbstractObjectiveModel, ConstraintsModelType<:CTModels.AbstractConstraintsModel}, ADNLPModels.ADNLPModel{Float64, Vector{Float64}, Vector{Int64}}}

Discretize an optimal control problem into a nonlinear optimization problem (ie direct transcription)

Arguments

  • ocp: optimal control problem as defined in CTModels
  • [description]: can specifiy for instance the NLP model and / or solver (:ipopt, :madnlp or :knitro)

Keyword arguments (optional)

  • grid_size: number of time steps for the discretized problem ([250])
  • disc_method: discretization method ([:trapeze], :euler, :euler_implicit, :midpoint, gauss_legendre_2, gauss_legendre_3)
  • time_grid: explicit time grid (can be non uniform)
  • init: info for the starting guess (values as named tuple or existing solution)
  • adnlp_backend: backend for automatic differentiation in ADNLPModels ([:optimized], :manual, :default)
  • show_time: (:true, [:false]) show timing details from ADNLPModels
source
CTDirect.get_OCP_control_at_time_stepMethod
get_OCP_control_at_time_step(
    xu,
    docp::CTDirect.DOCP,
    i
) -> Any

Retrieve control variables at given time step from the NLP variables. Convention: 1 <= i <= dimNLPsteps(+1), with convention u(tf) = U_N Vector output

source
CTDirect.get_OCP_control_at_time_stepMethod
get_OCP_control_at_time_step(
    xu,
    docp::CTDirect.DOCP{CTDirect.Euler},
    i
) -> Any

Retrieve control variables at given time step from the NLP variables. Convention: see above for acplicit / implicit versions Vector output

source
CTDirect.get_OCP_state_at_time_stepMethod
get_OCP_state_at_time_step(
    xu,
    docp::CTDirect.DOCP,
    i
) -> Any

Retrieve state variables at given time step from the NLP variables. Convention: 1 <= i <= dimNLPsteps+1 Vector output

source
CTDirect.get_OCP_variableMethod
get_OCP_variable(xu, docp::CTDirect.DOCP) -> Any

Retrieve optimization variables from the NLP variables. Convention: stored at the end, hence not dependent on the discretization method Vector output

source
CTDirect.get_lagrange_state_at_time_stepMethod
get_lagrange_state_at_time_step(
    xu,
    docp::CTDirect.DOCP,
    i
) -> Any

Retrieve state variable for lagrange cost at given time step from the NLP variables. Convention: 1 <= i <= dimNLPsteps+1 (no check for actual lagrange cost presence !)

source
CTDirect.get_stagevars_at_time_stepMethod
get_stagevars_at_time_step(
    xu,
    docp::CTDirect.DOCP,
    i,
    j
) -> Any

Retrieve stage variables at given time step/stage from the NLP variables. Convention: 1 <= i <= dimNLPsteps(+1), 1 <= j <= s Vector output Note that passing correct indices is up to the caller, no checks are made here.

source
CTDirect.get_time_gridMethod
get_time_grid(xu, docp::CTDirect.DOCP) -> Any

Return time grid for variable time problems (times are then dependent on NLP variables)

source
CTDirect.parse_DOCP_solution_dualMethod
parse_DOCP_solution_dual(
    docp,
    multipliers,
    constraints
) -> NTuple{5, Any}

Recover OCP costate and constraints multipliers from DOCP multipliers

source
CTDirect.setStepConstraints!Method
setStepConstraints!(
    docp::CTDirect.DOCP{<:CTDirect.GenericIRK},
    c,
    xu,
    v,
    time_grid,
    i,
    work
) -> Any

Set the constraints corresponding to the state equation Convention: 1 <= i <= dimNLPsteps (+1)

source
CTDirect.setStepConstraints!Method
setStepConstraints!(
    docp::CTDirect.DOCP{CTDirect.Euler},
    c,
    xu,
    v,
    time_grid,
    i,
    work
) -> Any

Set the constraints corresponding to the state equation Convention: 1 <= i <= dimNLPsteps+1

source
CTDirect.setStepConstraints!Method
setStepConstraints!(
    docp::CTDirect.DOCP{CTDirect.Midpoint},
    c,
    xu,
    v,
    time_grid,
    i,
    work
) -> Any

Set the constraints corresponding to the state equation Convention: 1 <= i <= dimNLPsteps+1

source
CTDirect.setStepConstraints!Method
setStepConstraints!(
    docp::CTDirect.DOCP{CTDirect.Trapeze},
    c,
    xu,
    v,
    time_grid,
    i,
    work
) -> Any

Set the constraints corresponding to the state equation Convention: 1 <= i <= dimNLPsteps+1

source
CTDirect.setWorkArrayMethod
setWorkArray(docp::CTDirect.DOCP, xu, time_grid, v) -> Any

Set work array for all dynamics and lagrange cost evaluations

source
CTDirect.setWorkArrayMethod
setWorkArray(
    docp::CTDirect.DOCP{<:CTDirect.GenericIRK},
    xu,
    time_grid,
    v
) -> Any

Set work array for all dynamics and lagrange cost evaluations

source
CTDirect.setWorkArrayMethod
setWorkArray(
    docp::CTDirect.DOCP{CTDirect.Euler},
    xu,
    time_grid,
    v
) -> Any

Set work array for all dynamics and lagrange cost evaluations

source
CTDirect.setWorkArrayMethod
setWorkArray(
    docp::CTDirect.DOCP{CTDirect.Midpoint},
    xu,
    time_grid,
    v
) -> Any

Set work array for all dynamics and lagrange cost evaluations

source
CTDirect.setWorkArrayMethod
setWorkArray(
    docp::CTDirect.DOCP{CTDirect.Trapeze},
    xu,
    time_grid,
    v
) -> Any

Set work array for all dynamics and lagrange cost evaluations

source
CTDirect.set_control_at_time_step!Method
set_control_at_time_step!(
    xu,
    u_init,
    docp::CTDirect.DOCP,
    i
) -> Any

Set initial guess for control variables at given time step Convention: 1 <= i <= dimNLPsteps(+1)

source
CTDirect.set_state_at_time_step!Method
set_state_at_time_step!(
    xu,
    x_init,
    docp::CTDirect.DOCP,
    i
) -> Any

Set initial guess for state variables at given time step Convention: 1 <= i <= dimNLPsteps+1

source
CTDirect.solveMethod
solve(
    ocp::CTModels.Model,
    description::Symbol...;
    display,
    grid_size,
    disc_method,
    time_grid,
    init,
    adnlp_backend,
    kwargs...
) -> CTModels.Solution{TimeGridModelType, TimesModelType, StateModelType, ControlModelType, VariableModelType, CostateModelType, Float64, CTModels.DualModel{CTModels.var"#69#81"{Matrix{Float64}}, CTModels.var"#70#82"{Matrix{Float64}}, Vector{Float64}, Vector{Float64}, CTModels.var"#71#83"{Matrix{Float64}, Int64}, CTModels.var"#72#84"{Matrix{Float64}, Int64}, CTModels.var"#73#85"{Matrix{Float64}, Int64}, CTModels.var"#74#86"{Matrix{Float64}, Int64}, Vector{Float64}, Vector{Float64}}, CTModels.SolverInfos{Dict{Symbol, Any}}} where {TimeGridModelType<:CTModels.TimeGridModel, TimesModelType<:CTModels.TimesModel, StateModelType<:Union{CTModels.StateModelSolution{TS} where TS<:CTModels.var"#63#75", CTModels.StateModelSolution{TS} where TS<:CTModels.var"#64#76"}, ControlModelType<:Union{CTModels.ControlModelSolution{TS} where TS<:CTModels.var"#65#77", CTModels.ControlModelSolution{TS} where TS<:CTModels.var"#66#78"}, VariableModelType<:Union{CTModels.VariableModelSolution{Vector{Float64}}, CTModels.VariableModelSolution{Float64}}, CostateModelType<:Union{CTModels.var"#67#79", CTModels.var"#68#80"}}

Solve an OCP with a direct method

Arguments

  • ocp: optimal control problem as defined in CTBase
  • [description]: can specifiy for instance the NLP model and / or solver (:ipopt, :madnlp or :knitro)

Keyword arguments (optional)

  • display: ([true], false) will disable output if set to false
  • grid_size: number of time steps for the discretized problem ([250])
  • disc_method: discretization method ([:trapeze], :midpoint, gauss_legendre_2)
  • time_grid: explicit time grid (can be non uniform)
  • init: info for the starting guess (values or existing solution)
  • adnlp_backend: backend for automatic differentiation in ADNLPModels ([:optimized], :manual, :default)

All further keywords are passed to the inner call of solve_docp

source
CTDirect.variables_bounds!Method
variables_bounds!(
    docp::CTDirect.DOCP
) -> Tuple{Vector{Float64}, Vector{Float64}}

Build upper and lower bounds vectors for the DOCP variable box constraints.

source