Hamiltonian

Index

Warning

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

CTFlows.FlowMethod
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
CTFlows.FlowMethod
Flow(
    hv::CTFlows.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
CTFlowsODE.hamiltonian_usageMethod
hamiltonian_usage(
    alg,
    abstol,
    reltol,
    saveat,
    internalnorm;
    kwargs_Flow...
) -> Any

Constructs a solver function for Hamiltonian systems, with configurable solver options.

Returns a callable object that integrates Hamilton's equations from an initial state (x0, p0) over a time span tspan = (t0, tf), with optional external parameters.

The returned function has two methods:

  • f(tspan, x0, p0, v=default_variable; kwargs...) → returns the full trajectory (solution object).
  • f(t0, x0, p0, tf, v=default_variable; kwargs...) → returns only the final (x, p) state.

Internally, it uses OrdinaryDiffEq.solve and supports events and callbacks.

Arguments

  • alg: integration algorithm (e.g. Tsit5()).
  • abstol, reltol: absolute and relative tolerances.
  • saveat: time points for saving.
  • internalnorm: norm used for adaptive integration.
  • kwargs_Flow...: additional keyword arguments for the solver.

Example

julia> flowfun = hamiltonian_usage(Tsit5(), 1e-8, 1e-8, 0.1, norm)
julia> xf, pf = flowfun(0.0, x0, p0, 1.0)
source
CTFlowsODE.rhsMethod
rhs(
    h::CTFlows.AbstractHamiltonian
) -> CTFlowsODE.var"#rhs!#43"{<:CTFlows.AbstractHamiltonian{TD, VD}} where {TD<:CTFlows.TimeDependence, VD<:CTFlows.VariableDependence}

Constructs the right-hand side of Hamilton's equations from a scalar Hamiltonian function.

Given a Hamiltonian h(t, x, p, l) (or h(x, p) in the autonomous case), returns an in-place function rhs!(dz, z, v, t) suitable for numerical integration.

This function computes the canonical Hamiltonian vector field using automatic differentiation:

dz[1:n]     =  ∂H/∂p
dz[n+1:2n]  = -∂H/∂x

Arguments

  • h: a subtype of CTFlows.AbstractHamiltonian defining the scalar Hamiltonian.

Returns

  • rhs!: a function for use in an ODE solver.
source