Extension Types
Index
CTFlowsODE.HamiltonianFlow
CTFlowsODE.ODEFlow
CTFlowsODE.OptimalControlFlow
CTFlowsODE.OptimalControlFlowSolution
CTFlowsODE.VectorFieldFlow
CTModels.Solution
In the examples in the documentation below, the methods are not prefixed by the module name even if they are private.
julia> using CTFlows
julia> x = 1
julia> private_fun(x) # throw an error
must be replaced by
julia> using CTFlows
julia> x = 1
julia> CTFlows.private_fun(x)
However, if the method is reexported by another package, then, there is no need of prefixing.
julia> module OptimalControl
import CTFlows: private_fun
export private_fun
end
julia> using OptimalControl
julia> x = 1
julia> private_fun(x)
Documentation
CTFlowsODE.HamiltonianFlow
— Typestruct HamiltonianFlow <: CTFlowsODE.AbstractFlow{Union{Real, AbstractVector{<:Real}}, Union{Real, AbstractVector{<:Real}}}
A flow object for integrating Hamiltonian dynamics in optimal control.
Represents the time evolution of a Hamiltonian system using the canonical form of Hamilton's equations. The struct holds the numerical integration setup and metadata for event handling.
Fields
f::Function
: Flow integrator function, called likef(t0, x0, p0, tf; ...)
.rhs!::Function
: Right-hand side of the ODE system, used by solvers.tstops::Times
: List of times at which integration should pause or apply discrete effects.jumps::Vector{Tuple{Time,Costate}}
: List of jump discontinuities for the costate at given times.
Usage
Instances of HamiltonianFlow
are callable and forward arguments to the underlying flow function f
.
Example
julia> flow = HamiltonianFlow(f, rhs!)
julia> xf, pf = flow(0.0, x0, p0, 1.0)
CTFlowsODE.ODEFlow
— Typestruct ODEFlow <: CTFlowsODE.AbstractFlow{Any, Any}
Generic flow object for arbitrary ODE systems with jumps and events.
A catch-all flow for general-purpose ODE integration. Supports dynamic typing and arbitrary state structures.
Fields
f::Function
: Integrator function called with time span and initial conditions.rhs::Function
: Right-hand side for the differential equation.tstops::Times
: Times at which the integrator is forced to stop.jumps::Vector{Tuple{Time,Any}}
: User-defined jumps applied to the state during integration.
Example
julia> flow = ODEFlow(f, rhs)
julia> result = flow(0.0, u0, 1.0)
CTFlowsODE.OptimalControlFlow
— Typestruct OptimalControlFlow{VD} <: CTFlowsODE.AbstractFlow{Union{Real, AbstractVector{<:Real}}, Union{Real, AbstractVector{<:Real}}}
A flow object representing the solution of an optimal control problem.
Supports Hamiltonian-based and classical formulations. Provides call overloads for different control settings:
- Fixed external variables
- Parametric (non-fixed) control problems
Fields
f::Function
: Main integrator that receives the RHS and other arguments.rhs!::Function
: ODE right-hand side.tstops::Times
: Times where the solver should stop (e.g., nonsmooth dynamics).jumps::Vector{Tuple{Time,Costate}}
: Costate jump conditions.feedback_control::ControlLaw
: Feedback lawu(t, x, p, v)
.ocp::Model
: The optimal control problem definition.kwargs_Flow::Any
: Extra solver arguments.
Call Signatures
F(t0, x0, p0, tf; kwargs...)
: Solves with fixed variable dimension.F(t0, x0, p0, tf, v; kwargs...)
: Solves with parameterv
.F(tspan, x0, p0; ...)
: Solves and returns a fullOptimalControlSolution
.
Example
julia> flow = OptimalControlFlow(...)
julia> sol = flow(0.0, x0, p0, 1.0)
julia> opt_sol = flow((0.0, 1.0), x0, p0)
CTFlowsODE.OptimalControlFlowSolution
— Typestruct OptimalControlFlowSolution
Wraps the low-level ODE solution, control feedback law, model structure, and problem parameters.
Fields
ode_sol::Any
: The ODE solution (from DifferentialEquations.jl).feedback_control::ControlLaw
: Feedback control lawu(t, x, p, v)
.ocp::Model
: The optimal control model used.variable::Variable
: External or design parameters of the control problem.
Usage
You can evaluate the flow solution like a callable ODE solution.
Example
julia> sol = OptimalControlFlowSolution(ode_sol, u, model, v)
julia> x = sol(t)
CTFlowsODE.VectorFieldFlow
— Typestruct VectorFieldFlow <: CTFlowsODE.AbstractFlow{Union{Real, AbstractVector{<:Real}}, Union{Real, AbstractVector{<:Real}}}
A flow object for integrating general vector field dynamics.
Used for systems where the vector field is given explicitly, rather than derived from a Hamiltonian. Useful in settings like controlled systems or classical mechanics outside the Hamiltonian framework.
Fields
f::Function
: Flow integrator function.rhs::Function
: ODE right-hand side function.tstops::Times
: Event times (e.g., to trigger callbacks).jumps::Vector{Tuple{Time,State}}
: Discrete jump events on the state trajectory.
Example
julia> flow = VectorFieldFlow(f, rhs)
julia> xf = flow(0.0, x0, 1.0)
CTModels.Solution
— MethodSolution(
ocfs::CTFlowsODE.OptimalControlFlowSolution;
kwargs...
) -> 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{CTModels.var"#114#136"{CTFlowsODE.var"#8#15"{CTFlowsODE.var"#x#11"{CTFlowsODE.OptimalControlFlowSolution, Int64}}}}, CTModels.StateModelSolution{CTModels.var"#115#137"{CTFlowsODE.var"#8#15"{CTFlowsODE.var"#x#11"{CTFlowsODE.OptimalControlFlowSolution, Int64}}}}}, ControlModelType<:Union{CTModels.ControlModelSolution{CTModels.var"#116#138"{CTFlowsODE.var"#9#16"{CTFlowsODE.var"#u#13"{CTFlowsODE.OptimalControlFlowSolution, CTFlowsODE.var"#p#12"{CTFlowsODE.OptimalControlFlowSolution, Int64}, CTFlowsODE.var"#x#11"{CTFlowsODE.OptimalControlFlowSolution, Int64}}}}}, CTModels.ControlModelSolution{CTModels.var"#117#139"{CTFlowsODE.var"#9#16"{CTFlowsODE.var"#u#13"{CTFlowsODE.OptimalControlFlowSolution, CTFlowsODE.var"#p#12"{CTFlowsODE.OptimalControlFlowSolution, Int64}, CTFlowsODE.var"#x#11"{CTFlowsODE.OptimalControlFlowSolution, Int64}}}}}}, VariableModelType<:Union{CTModels.VariableModelSolution{Vector{Float64}}, CTModels.VariableModelSolution{Float64}}, CostateModelType<:Union{CTModels.var"#118#140"{CTFlowsODE.var"#10#17"{CTFlowsODE.var"#p#12"{CTFlowsODE.OptimalControlFlowSolution, Int64}}}, CTModels.var"#119#141"{CTFlowsODE.var"#10#17"{CTFlowsODE.var"#p#12"{CTFlowsODE.OptimalControlFlowSolution, Int64}}}}, 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"#121#143", CTModels.var"#122#144"}, BC_Dual<:Union{Nothing, Vector{Float64}}, SC_LB_Dual<:Union{Nothing, CTModels.var"#124#146", CTModels.var"#125#147"}, SC_UB_Dual<:Union{Nothing, CTModels.var"#127#149", CTModels.var"#128#150"}, CC_LB_Dual<:Union{Nothing, CTModels.var"#130#152", CTModels.var"#131#153"}, CC_UB_Dual<:Union{Nothing, CTModels.var"#133#155", CTModels.var"#134#156"}, VC_LB_Dual<:Union{Nothing, Vector{Float64}}, VC_UB_Dual<:Union{Nothing, Vector{Float64}}}), ModelType<:CTModels.Model}
Constructs an OptimalControlSolution
from an OptimalControlFlowSolution
.
This evaluates the objective (Mayer and/or Lagrange costs), extracts the time-dependent state, costate, and control trajectories, and builds a full CTModels.Solution
.
Returns a CTModels.Solution
ready for evaluation, reporting, or analysis.
Keyword Arguments
alg
: Optional solver for computing Lagrange integral, if needed.- Additional kwargs passed to the internal solver.
Example
julia> sol = Solution(optflow_solution)