Solution

Index

Documentation

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

Retrieve convergence information from an NLP solution.

Arguments

  • nlp_solution: A solver execution statistics object.

Returns

  • (objective, iterations, constraints_violation, message, status, successful): A tuple containing the final objective value, iteration count, primal feasibility, solver message, solver status, and success flag.

Example

julia> SolverInfos(nlp_solution)
(1.23, 15, 1.0e-6, "Ipopt/generic", :first_order, true)
source
CTDirect.SolverInfosMethod
SolverInfos(

) -> Tuple{Float64, Int64, Float64, String, Symbol, Bool}

Return default convergence information for an NLP solution.

Returns

  • (objective, iterations, constraints_violation, message, status, successful): Default values representing an undefined solver state.

Example

julia> SolverInfos()
(0.0, 0, 0.0, "undefined", :undefined, true)
source
CTDirect.build_OCP_solutionMethod
build_OCP_solution(
    docp,
    nlp_solution;
    nlp_model,
    nlp_solver
) -> CTModels.Solution{TimeGridModelType, TimesModelType, StateModelType, ControlModelType, VariableModelType, CostateModelType, Float64, DualModelType, CTModels.SolverInfos{Dict{Symbol, Any}}, ModelType} where {TimeGridModelType<:CTModels.TimeGridModel, TimesModelType<:CTModels.TimesModel, StateModelType<:Union{CTModels.StateModelSolution{TS} where TS<:CTModels.var"#98#120", CTModels.StateModelSolution{TS} where TS<:CTModels.var"#99#121"}, ControlModelType<:Union{CTModels.ControlModelSolution{TS} where TS<:CTModels.var"#100#122", CTModels.ControlModelSolution{TS} where TS<:CTModels.var"#101#123"}, VariableModelType<:Union{CTModels.VariableModelSolution{Vector{Float64}}, CTModels.VariableModelSolution{Float64}}, CostateModelType<:Union{CTModels.var"#102#124", CTModels.var"#103#125"}, DualModelType<:(CTModels.DualModel{PC_Dual, BC_Dual, SC_LB_Dual, SC_UB_Dual, CC_LB_Dual, CC_UB_Dual, VC_LB_Dual, VC_UB_Dual} where {PC_Dual<:Union{Nothing, CTModels.var"#105#127", CTModels.var"#106#128"}, BC_Dual<:Union{Nothing, Vector{Float64}}, SC_LB_Dual<:Union{Nothing, CTModels.var"#108#130", CTModels.var"#109#131"}, SC_UB_Dual<:Union{Nothing, CTModels.var"#111#133", CTModels.var"#112#134"}, CC_LB_Dual<:Union{Nothing, CTModels.var"#114#136", CTModels.var"#115#137"}, CC_UB_Dual<:Union{Nothing, CTModels.var"#117#139", CTModels.var"#118#140"}, VC_LB_Dual<:Union{Nothing, Vector{Float64}}, VC_UB_Dual<:Union{Nothing, Vector{Float64}}}), ModelType<:CTModels.Model}

Build an OCP functional solution from a DOCP discrete solution given as a SolverCore.GenericExecutionStats object.

Arguments

  • docp: The discretized optimal control problem (DOCP).
  • nlp_solution: A solver execution statistics object.
  • nlp_model: The NLP model backend (default: ADNLPBackend()).
  • nlp_solver: The NLP solver backend (default: IpoptBackend()).

Returns

  • solution::CTModels.Solution: A functional OCP solution containing trajectories, multipliers, and solver information.

Example

julia> build_OCP_solution(docp, nlp_solution)
CTModels.Solution(...)
source
CTDirect.build_OCP_solutionMethod
build_OCP_solution(
    docp;
    primal,
    dual,
    multipliers_L,
    multipliers_U,
    nlp_model,
    nlp_solution
)

Build an OCP functional solution from a DOCP discrete solution, given explicit primal variables, and optionally dual variables and bound multipliers.

Arguments

  • docp: The discretized optimal control problem (DOCP).
  • primal: Array of primal decision variables.
  • dual: Array of dual variables (default: nothing).
  • multipliers_L: Lower bound multipliers (default: nothing).
  • multipliers_U: Upper bound multipliers (default: nothing).
  • nlp_model: The NLP model backend (default: ADNLPBackend()).
  • nlp_solution: A solver execution statistics object.

Returns

  • solution::CTModels.Solution: A functional OCP solution with trajectories, multipliers, and solver information.

Example

julia> build_OCP_solution(docp; primal=primal_vars, nlp_solution=nlp_solution)
CTModels.Solution(...)
source
CTDirect.is_emptyMethod
is_empty(t) -> Any

Check whether a collection t is empty or not defined.

Arguments

  • t: Any object that may be nothing or support length.

Returns

  • ::Bool: true if t is nothing or has length zero, otherwise false.

Example

julia> is_empty([])
true

julia> is_empty([1, 2, 3])
false

julia> is_empty(nothing)
true
source
CTDirect.parse_DOCP_solution_dualMethod
parse_DOCP_solution_dual(
    docp,
    multipliers;
    nlp_model,
    nlp_solution
)

Recover OCP costates and constraint multipliers from DOCP dual variables.

Arguments

  • docp: The discretized optimal control problem (DOCP).
  • multipliers: Array of dual variables (may be nothing).
  • nlp_model: The NLP model backend (default: ADNLPBackend()).
  • nlp_solution: A solver execution statistics object.

Returns

  • (P, path_constraints_dual, boundary_constraints_dual):
    • P: Costate trajectory.
    • path_constraints_dual: Path constraint multipliers.
    • boundary_constraints_dual: Boundary constraint multipliers.

Example

julia> P, path_dual, bound_dual = parse_DOCP_solution_dual(docp, duals; nlp_model=nlp_model, nlp_solution=nlp_solution)
([...] , [...], [...])
source
CTDirect.parse_DOCP_solution_primalMethod
parse_DOCP_solution_primal(
    docp,
    solution;
    multipliers_L,
    multipliers_U,
    nlp_model,
    nlp_solution
)

Recover OCP state, control, and optimization variables from DOCP primal variables. Bound multipliers are also parsed if available.

Arguments

  • docp: The discretized optimal control problem (DOCP).
  • solution: Array of primal decision variables.
  • multipliers_L: Lower bound multipliers.
  • multipliers_U: Upper bound multipliers.
  • nlp_model: The NLP model backend.
  • nlp_solution: A solver execution statistics object.

Returns

  • (X, U, v, box_multipliers):
    • X: State trajectory.
    • U: Control trajectory.
    • v: Optimization variables.
    • box_multipliers: Tuple of bound multipliers for states, controls, and variables.

Example

julia> X, U, v, box_mults = parse_DOCP_solution_primal(docp, primal;
       multipliers_L=mL, multipliers_U=mU, nlp_model=nlp_model, nlp_solution=nlp_solution)
([...] , [...], [...], (...))
source