Optimal Control Problem
Index
CTFlows.FlowCTFlows.FlowCTFlows.FlowCTFlows.FlowCTFlows.FlowCTFlows.FlowCTFlows.FlowCTFlows.FlowCTFlows.FlowCTFlowsODE.__ocp_Flow
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 errormust 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
CTFlows.Flow — Method
Flow(
ocp::CTModels.OCP.Model,
u::Function,
g::Function,
μ::Function;
autonomous,
variable,
alg,
abstol,
reltol,
saveat,
internalnorm,
kwargs_Flow...
) -> Union{CTFlowsODE.OptimalControlFlow{CTFlows.Fixed}, CTFlowsODE.OptimalControlFlow{CTFlows.NonFixed}}
Construct a flow from a raw feedback control, constraint, and multiplier.
This version is for defining flows directly from user functions without wrapping them into ControlLaw, Constraint, or Multiplier types. Automatically wraps and adapts them based on time dependence.
Arguments
ocp::CTModels.Model: The optimal control problem.u::Function: Control law.g::Function: Constraint.μ::Function: Multiplier.autonomous::Bool: Whether the system is autonomous.variable::Bool: Whether time is a free variable.alg,abstol,reltol,saveat,internalnorm: Solver parameters.kwargs_Flow: Additional options.
Returns
A Flow object ready for trajectory integration.
CTFlows.Flow — Method
Flow(
ocp::CTModels.OCP.Model{<:CTModels.OCP.TimeDependence, <:CTModels.OCP.AbstractTimesModel, <:CTModels.OCP.AbstractStateModel, CTModels.OCP.ControlModel},
u::CTFlows.ControlLaw;
alg,
abstol,
reltol,
saveat,
internalnorm,
kwargs_Flow...
) -> Union{CTFlowsODE.OptimalControlFlow{CTFlows.Fixed}, CTFlowsODE.OptimalControlFlow{CTFlows.NonFixed}}
Construct a flow for an optimal control problem using a given control law.
This method builds the Hamiltonian system associated with the optimal control problem (ocp) and integrates the corresponding state–costate dynamics using the specified control law u.
Arguments
ocp::WithControlModel: An optimal control problem with control variables.u::CTFlows.ControlLaw: A feedback control law generated byControlLaw(...)or similar.alg: Integration algorithm (default inferred).abstol: Absolute tolerance for the ODE solver.reltol: Relative tolerance for the ODE solver.saveat: Time points at which to save the solution.internalnorm: Optional norm function used by the integrator.kwargs_Flow: Additional keyword arguments passed to the solver.
Returns
A flow object f such that:
f(t0, x0, p0, tf)integrates the state and costate fromt0totf.f((t0, tf), x0, p0)returns the full trajectory over the interval.
Throws
CTBase.Exceptions.PreconditionError: If called on a control-free problem (EmptyControlModel).
Example
julia> u = (x, p) -> p
julia> f = Flow(ocp, ControlLaw(u))Notes
For control-free problems (parameter estimation, optimal design), use Flow(ocp) instead.
CTFlows.Flow — Method
Flow(
ocp::CTModels.OCP.Model{<:CTModels.OCP.TimeDependence, <:CTModels.OCP.AbstractTimesModel, <:CTModels.OCP.AbstractStateModel, CTModels.OCP.ControlModel},
u::Function;
autonomous,
variable,
alg,
abstol,
reltol,
saveat,
internalnorm,
kwargs_Flow...
) -> Union{CTFlowsODE.OptimalControlFlow{CTFlows.Fixed}, CTFlowsODE.OptimalControlFlow{CTFlows.NonFixed}}
Construct a flow for an optimal control problem using a control function in feedback form.
This method constructs the Hamiltonian and integrates the associated state–costate dynamics using a raw function u. It automatically wraps u as a control law.
Arguments
ocp::WithControlModel: The optimal control problem with control variables.u::Function: A feedback control function:- If
ocpis autonomous:u(x, p) - If non-autonomous:
u(t, x, p)
- If
autonomous::Bool: Whether the control law depends on time.variable::Bool: Whether the OCP involves variable time (e.g., free final time).alg,abstol,reltol,saveat,internalnorm: ODE solver parameters.kwargs_Flow: Additional options.
Returns
A Flow object compatible with function call interfaces for state propagation.
Example
julia> u = (t, x, p) -> t + p
julia> f = Flow(ocp, u)CTFlows.Flow — Method
Flow(
ocp::CTModels.OCP.Model{<:CTModels.OCP.TimeDependence, <:CTModels.OCP.AbstractTimesModel, <:CTModels.OCP.AbstractStateModel, CTModels.OCP.EmptyControlModel},
u::CTFlows.ControlLaw;
kwargs...
)
Guard method that prevents providing a control law to a control-free problem.
Throws
CTBase.Exceptions.PreconditionError: Always throws with a clear error message.
CTFlows.Flow — Method
Flow(
ocp::CTModels.OCP.Model{<:CTModels.OCP.TimeDependence, <:CTModels.OCP.AbstractTimesModel, <:CTModels.OCP.AbstractStateModel, CTModels.OCP.EmptyControlModel},
g::Function,
μ::Function;
autonomous,
variable,
alg,
abstol,
reltol,
saveat,
internalnorm,
kwargs_Flow...
) -> Union{CTFlowsODE.OptimalControlFlow{CTFlows.Fixed}, CTFlowsODE.OptimalControlFlow{CTFlows.NonFixed}}
Construct a flow for a control-free problem with raw constraint and multiplier functions.
This version is for defining flows directly from user functions without wrapping them into StateConstraint or Multiplier types. Automatically wraps and adapts them based on time dependence.
Arguments
ocp::ControlFreeModel: The control-free optimal control problem.g::Function: State constraint function.μ::Function: Multiplier function.autonomous::Bool: Whether the system is autonomous.variable::Bool: Whether time is a free variable.alg,abstol,reltol,saveat,internalnorm: Solver parameters.kwargs_Flow: Additional options.
Returns
A Flow object ready for trajectory integration.
Example
julia> # Control-free with raw functions
julia> f = Flow(ocp, (x, v) -> x[1] - 1.0, (x, p, v) -> p[1]; autonomous=true, variable=true)CTFlows.Flow — Method
Flow(
ocp::CTModels.OCP.Model{<:CTModels.OCP.TimeDependence, <:CTModels.OCP.AbstractTimesModel, <:CTModels.OCP.AbstractStateModel, CTModels.OCP.EmptyControlModel},
u::Function;
kwargs...
)
Guard method that prevents providing a control function to a control-free problem.
Throws
CTBase.Exceptions.PreconditionError: Always throws with a clear error message.
CTFlows.Flow — Method
Flow(
ocp::CTModels.OCP.Model{<:CTModels.OCP.TimeDependence, <:CTModels.OCP.AbstractTimesModel, <:CTModels.OCP.AbstractStateModel, CTModels.OCP.EmptyControlModel};
alg,
abstol,
reltol,
saveat,
internalnorm,
kwargs_Flow...
) -> Union{CTFlowsODE.OptimalControlFlow{CTFlows.Fixed}, CTFlowsODE.OptimalControlFlow{CTFlows.NonFixed}}
Construct a flow for a control-free optimal control problem.
This method builds the Hamiltonian system for problems without control variables (parameter estimation, optimal design). The control is internally set to an empty array.
Arguments
ocp::ControlFreeModel: A control-free optimal control problem (withEmptyControlModel).alg: Integration algorithm (default inferred).abstol: Absolute tolerance for the ODE solver.reltol: Relative tolerance for the ODE solver.saveat: Time points at which to save the solution.internalnorm: Optional norm function used by the integrator.kwargs_Flow: Additional keyword arguments passed to the solver.
Returns
A flow object f such that:
f(t0, x0, p0, tf)integrates the state and costate fromt0totf.f((t0, tf), x0, p0)returns the full trajectory over the interval.
Example
julia> # Control-free problem (parameter estimation)
julia> ocp = @def begin
λ ∈ R, variable
t ∈ [0, 10], time
x ∈ R, state
x(0) == 2.0
ẋ(t) == λ * x(t)
∫((x(t) - data(t))^2) → min
end
julia> f = Flow(ocp)
julia> sol = f((0.0, 10.0), x0, p0, λ)Notes
This method is for control-free problems only. For problems with control variables, use Flow(ocp, u) instead.
CTFlows.Flow — Method
Flow(
ocp::CTModels.OCP.Model,
u::Union{CTFlows.ControlLaw{<:Function, T, V}, CTFlows.FeedbackControl{<:Function, T, V}},
g::Union{CTFlows.MixedConstraint{<:Function, T, V}, CTFlows.StateConstraint{<:Function, T, V}},
μ::CTFlows.Multiplier{<:Function, T, V};
alg,
abstol,
reltol,
saveat,
internalnorm,
kwargs_Flow...
) -> Union{CTFlowsODE.OptimalControlFlow{CTFlows.Fixed}, CTFlowsODE.OptimalControlFlow{CTFlows.NonFixed}}
Construct a flow for an optimal control problem with control and constraint multipliers in feedback form.
This variant constructs a Hamiltonian system incorporating both the control law and a multiplier law (e.g., for enforcing state or mixed constraints). All inputs must be consistent in time dependence.
Arguments
ocp::CTModels.Model: The optimal control problem.u::ControlLaw or FeedbackControl: Feedback control.g::StateConstraint or MixedConstraint: Constraint function.μ::Multiplier: Multiplier function.alg,abstol,reltol,saveat,internalnorm: Solver settings.kwargs_Flow: Additional options.
Returns
A Flow object that integrates the constrained Hamiltonian dynamics.
Example
julia> f = Flow(ocp, (x, p) -> p[1], (x, u) -> x[1] - 1, (x, p) -> x[1]+p[1])For non-autonomous cases:
julia> f = Flow(ocp, (t, x, p) -> t + p, (t, x, u) -> x - 1, (t, x, p) -> x+p)CTFlows.Flow — Method
Flow(
ocp::CTModels.OCP.Model{<:CTModels.OCP.TimeDependence, <:CTModels.OCP.AbstractTimesModel, <:CTModels.OCP.AbstractStateModel, CTModels.OCP.EmptyControlModel},
g::CTFlows.StateConstraint{<:Function, T, V},
μ::CTFlows.Multiplier{<:Function, T, V};
alg,
abstol,
reltol,
saveat,
internalnorm,
kwargs_Flow...
) -> Union{CTFlowsODE.OptimalControlFlow{CTFlows.Fixed}, CTFlowsODE.OptimalControlFlow{CTFlows.NonFixed}}
Construct a flow for a control-free optimal control problem with state constraint and multiplier.
This method builds the Hamiltonian system for control-free problems with state constraints. The control is internally set to an empty array.
Arguments
ocp::ControlFreeModel: A control-free optimal control problem (withEmptyControlModel).g::CTFlows.StateConstraint: State constraint function.μ::CTFlows.Multiplier: Multiplier function associated with the constraint.alg,abstol,reltol,saveat,internalnorm: Solver settings.kwargs_Flow: Additional options.
Returns
A Flow object that integrates the constrained Hamiltonian dynamics.
Example
julia> # Control-free problem with state constraint
julia> g = StateConstraint((x, v) -> x[1] - 1.0)
julia> μ = Multiplier((x, p, v) -> p[1])
julia> f = Flow(ocp, g, μ)CTFlowsODE.__ocp_Flow — Method
__ocp_Flow(
ocp::CTModels.OCP.Model,
h::CTFlows.Hamiltonian,
u::CTFlows.ControlLaw,
alg,
abstol,
reltol,
saveat,
internalnorm;
kwargs_Flow...
) -> Union{CTFlowsODE.OptimalControlFlow{CTFlows.Fixed}, CTFlowsODE.OptimalControlFlow{CTFlows.NonFixed}}
Internal helper: builds the OptimalControlFlow object from a Hamiltonian and control law.
This function assembles the RHS, constructs the integrator, and packages the flow object.
Arguments
ocp: The original optimal control problem.h: A Hamiltonian structure.u: A control law.alg,abstol,reltol,saveat,internalnorm: Integration parameters.kwargs_Flow: Additional parameters.
Returns
An OptimalControlFlow object, callable as a function for integration.