API of the Flow function

CTFlows.FlowFunction
Flow(
    vf::VectorField;
    alg,
    abstol,
    reltol,
    saveat,
    internalnorm,
    kwargs_Flow...
) -> CTFlowsODE.VectorFieldFlow

Constructs a flow object for a classical (non-Hamiltonian) vector field.

This creates a VectorFieldFlow that integrates the ODE system dx/dt = vf(t, x, v) using DifferentialEquations.jl. It handles both fixed and parametric dynamics, as well as jump discontinuities and event stopping.

Keyword Arguments

  • alg, abstol, reltol, saveat, internalnorm: Solver options.
  • kwargs_Flow...: Additional arguments passed to the solver configuration.

Example

julia> vf(t, x, v) = -v * x
julia> flow = CTFlows.Flow(CTFlows.VectorField(vf))
julia> x1 = flow(0.0, 1.0, 1.0)
source
Flow(
    h::CTFlows.AbstractHamiltonian;
    alg,
    abstol,
    reltol,
    saveat,
    internalnorm,
    kwargs_Flow...
) -> CTFlowsODE.HamiltonianFlow

Constructs a Hamiltonian flow from a scalar Hamiltonian.

This method builds a numerical integrator that simulates the evolution of a Hamiltonian system given a Hamiltonian function h(t, x, p, l) or h(x, p).

Internally, it computes the right-hand side of Hamilton’s equations via automatic differentiation and returns a HamiltonianFlow object.

Keyword Arguments

  • alg, abstol, reltol, saveat, internalnorm: solver options.
  • kwargs_Flow...: forwarded to the solver.

Example

julia> H(x, p) = dot(p, p) + dot(x, x)
julia> flow = CTFlows.Flow(CTFlows.Hamiltonian(H))
julia> xf, pf = flow(0.0, x0, p0, 1.0)
source
Flow(
    hv::HamiltonianVectorField;
    alg,
    abstol,
    reltol,
    saveat,
    internalnorm,
    kwargs_Flow...
) -> CTFlowsODE.HamiltonianFlow

Constructs a Hamiltonian flow from a precomputed Hamiltonian vector field.

This method assumes you already provide the Hamiltonian vector field (dx/dt, dp/dt) instead of deriving it from a scalar Hamiltonian.

Returns a HamiltonianFlow object that integrates the given system.

Keyword Arguments

  • alg, abstol, reltol, saveat, internalnorm: solver options.
  • kwargs_Flow...: forwarded to the solver.

Example

julia> hv(t, x, p, l) = (∇ₚH, -∇ₓH)
julia> flow = CTFlows.Flow(CTFlows.HamiltonianVectorField(hv))
julia> xf, pf = flow(0.0, x0, p0, 1.0, l)
source
Flow(
    ocp::Model,
    u::CTFlows.ControlLaw;
    alg,
    abstol,
    reltol,
    saveat,
    internalnorm,
    kwargs_Flow...
) -> CTFlowsODE.OptimalControlFlow

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::CTModels.Model: An optimal control problem defined using CTModels.
  • u::CTFlows.ControlLaw: A feedback control law generated by ControlLaw(...) 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 from t0 to tf.
  • f((t0, tf), x0, p0) returns the full trajectory over the interval.

Example

julia> u = (x, p) -> p
julia> f = Flow(ocp, ControlLaw(u))
source
Flow(
    ocp::Model,
    u::Function;
    autonomous,
    variable,
    alg,
    abstol,
    reltol,
    saveat,
    internalnorm,
    kwargs_Flow...
) -> CTFlowsODE.OptimalControlFlow

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::CTModels.Model: The optimal control problem.
  • u::Function: A feedback control function:
    • If ocp is autonomous: u(x, p)
    • If non-autonomous: u(t, x, p)
  • 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)
source
Flow(
    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...
) -> CTFlowsODE.OptimalControlFlow

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)
Warning

All input functions must match the autonomous/non-autonomous nature of the problem.

source
Flow(
    ocp::Model,
    u::Function,
    g::Function,
    μ::Function;
    autonomous,
    variable,
    alg,
    abstol,
    reltol,
    saveat,
    internalnorm,
    kwargs_Flow...
) -> CTFlowsODE.OptimalControlFlow

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.

source
Flow(
    dyn::Function;
    autonomous,
    variable,
    alg,
    abstol,
    reltol,
    saveat,
    internalnorm,
    kwargs_Flow...
) -> CTFlowsODE.ODEFlow

Constructs a Flow from a user-defined dynamical system given as a Julia function.

This high-level interface handles:

  • autonomous and non-autonomous systems,
  • presence or absence of additional variables (v),
  • selection of ODE solvers and tolerances,
  • and integrates with the CTFlows event system (e.g., jumps, callbacks).

Arguments

  • dyn: A function defining the vector field. Its signature must match the values of autonomous and variable.
  • autonomous: Whether the dynamics are time-independent (false by default).
  • variable: Whether the dynamics depend on a control or parameter v.
  • alg, abstol, reltol, saveat, internalnorm: Solver settings passed to OrdinaryDiffEq.solve.
  • kwargs_Flow: Additional keyword arguments passed to the solver.

Returns

An ODEFlow object, wrapping both the full solver and its right-hand side (RHS).

Supported Function Signatures for dyn

Depending on the (autonomous, variable) flags:

  • (false, false): dyn(x)
  • (false, true): dyn(x, v)
  • (true, false): dyn(t, x)
  • (true, true): dyn(t, x, v)

Example

julia> dyn(t, x, v) = [-x[1] + v[1] * sin(t)]
julia> flow = CTFlows.Flow(dyn; autonomous=true, variable=true)
julia> xT = flow((0.0, 1.0), [1.0], [0.1])
source