Differential Geometry
Index
- CTFlows.:⅋
- CTFlows.:⅋
- CTFlows.:⋅
- CTFlows.:⋅
- CTFlows.:⋅
- CTFlows.Lie
- CTFlows.Lie
- CTFlows.Lie
- CTFlows.Lie
- CTFlows.Lift
- CTFlows.Lift
- CTFlows.Poisson
- CTFlows.Poisson
- CTFlows.Poisson
- CTFlows.Poisson
- CTFlows.Poisson
- CTFlows.Poisson
- CTFlows.∂ₜ
- CTFlows.@Lie
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.:⅋ — Method"Directional derivative" of a vector field in the autonomous case, used internally for computing the Lie bracket.
Example:
julia> X = VectorField(x -> [x[2], -x[1]])
julia> Y = VectorField(x -> [x[1], x[2]])
julia> (X ⅋ Y)([1, 2])
[2, -1]CTFlows.:⅋ — Method"Directional derivative" of a vector field in the nonautonomous case, used internally for computing the Lie bracket.
Example:
julia> X = VectorField((t, x, v) -> [t + v[1] + v[2] + x[2], -x[1]], NonFixed, NonAutonomous)
julia> Y = VectorField((t, x, v) ->  [v[1] + v[2] + x[1], x[2]], NonFixed, NonAutonomous)
julia> (X ⅋ Y)(1, [1, 2], [2, 3])
[8, -1]CTFlows.:⋅ — MethodLie derivative of a scalar function along a vector field in the autonomous case.
Example:
julia> φ = x -> [x[2], -x[1]]
julia> X = VectorField(φ)
julia> f = x -> x[1]^2 + x[2]^2
julia> (X⋅f)([1, 2])
0CTFlows.:⋅ — MethodLie derivative of a scalar function along a vector field in the nonautonomous case.
Example:
julia> φ = (t, x, v) -> [t + x[2] + v[1], -x[1] + v[2]]
julia> X = VectorField(φ, NonAutonomous, NonFixed)
julia> f = (t, x, v) -> t + x[1]^2 + x[2]^2
julia> (X⋅f)(1, [1, 2], [2, 1])
10CTFlows.:⋅ — MethodLie derivative of a scalar function along a function (considered autonomous and non-variable).
Example:
julia> φ = x -> [x[2], -x[1]]
julia> f = x -> x[1]^2 + x[2]^2
julia> (φ⋅f)([1, 2])
0
julia> φ = (t, x, v) -> [t + x[2] + v[1], -x[1] + v[2]]
julia> f = (t, x, v) -> t + x[1]^2 + x[2]^2
julia> (φ⋅f)(1, [1, 2], [2, 1])
MethodErrorCTFlows.Lie — MethodLie derivative of a scalar function along a vector field.
Example:
julia> φ = x -> [x[2], -x[1]]
julia> X = VectorField(φ)
julia> f = x -> x[1]^2 + x[2]^2
julia> Lie(X,f)([1, 2])
0
julia> φ = (t, x, v) -> [t + x[2] + v[1], -x[1] + v[2]]
julia> X = VectorField(φ, NonAutonomous, NonFixed)
julia> f = (t, x, v) -> t + x[1]^2 + x[2]^2
julia> Lie(X, f)(1, [1, 2], [2, 1])
10CTFlows.Lie — MethodLie derivative of a scalar function along a function with specified dependencies.
Example:
julia> φ = x -> [x[2], -x[1]]
julia> f = x -> x[1]^2 + x[2]^2
julia> Lie(φ,f)([1, 2])
0
julia> φ = (t, x, v) -> [t + x[2] + v[1], -x[1] + v[2]]
julia> f = (t, x, v) -> t + x[1]^2 + x[2]^2
julia> Lie(φ, f, autonomous=false, variable=true)(1, [1, 2], [2, 1])
10CTFlows.Lie — MethodLie bracket of two vector fields in the autonomous case.
Example:
julia> f = x -> [x[2], 2x[1]]
julia> g = x -> [3x[2], -x[1]]
julia> X = VectorField(f)
julia> Y = VectorField(g)
julia> Lie(X, Y)([1, 2])
[7, -14]CTFlows.Lie — MethodLie bracket of two vector fields in the nonautonomous case.
Example:
julia> f = (t, x, v) -> [t + x[2] + v, -2x[1] - v]
julia> g = (t, x, v) -> [t + 3x[2] + v, -x[1] - v]
julia> X = VectorField(f, NonAutonomous, NonFixed)
julia> Y = VectorField(g, NonAutonomous, NonFixed)
julia> Lie(X, Y)(1, [1, 2], 1)
[-7, 12]CTFlows.Lift — MethodLift(
    X::CTFlows.VectorField
) -> CTFlows.HamiltonianLift{CTFlows.VectorField{TF, TD, VD}} where {TF<:Function, TD<:CTFlows.TimeDependence, VD<:CTFlows.VariableDependence}
Construct the Hamiltonian lift of a VectorField.
Arguments
- X::VectorField: The vector field to lift. Its signature determines if it is autonomous and/or variable.
Returns
- A HamiltonianLiftcallable object representing the Hamiltonian lift ofX.
Examples
julia> HL = Lift(VectorField(x -> [x[1]^2, x[2]^2], autonomous=true, variable=false))
julia> HL([1, 0], [0, 1])  # returns 0
julia> HL2 = Lift(VectorField((t, x, v) -> [t + x[1]^2, x[2]^2 + v], autonomous=false, variable=true))
julia> HL2(1, [1, 0], [0, 1], 1)  # returns 1
julia> H = Lift(x -> 2x)
julia> H(1, 1)  # returns 2
julia> H2 = Lift((t, x, v) -> 2x + t - v, autonomous=false, variable=true)
julia> H2(1, 1, 1, 1)  # returns 2
# Alternative syntax using symbols for autonomy and variability
julia> H3 = Lift((t, x, v) -> 2x + t - v, NonAutonomous, NonFixed)
julia> H3(1, 1, 1, 1)  # returns 2CTFlows.Lift — MethodLift(
    X::Function;
    autonomous,
    variable
) -> CTFlows.var"#19#23"{<:Function}
Construct the Hamiltonian lift of a function.
Arguments
- X::Function: The function representing the vector field.
- autonomous::Bool=true: Whether the function is autonomous (time-independent).
- variable::Bool=false: Whether the function depends on an additional variable argument.
Returns
- A callable function computing the Hamiltonian lift,
(and variants depending on autonomous and variable).
Details
Depending on the autonomous and variable flags, the returned function has one of the following call signatures:
- (x, p)if- autonomous=trueand- variable=false
- (x, p, v)if- autonomous=trueand- variable=true
- (t, x, p)if- autonomous=falseand- variable=false
- (t, x, p, v)if- autonomous=falseand- variable=true
Examples
julia> H = Lift(x -> 2x)
julia> H(1, 1)  # returns 2
julia> H2 = Lift((t, x, v) -> 2x + t - v, autonomous=false, variable=true)
julia> H2(1, 1, 1, 1)  # returns 2CTFlows.Poisson — MethodPoisson(
    f::Function,
    g::Function;
    autonomous,
    variable
) -> CTFlows.Hamiltonian
Poisson bracket of two functions. The time and variable dependence are specified with keyword arguments.
Returns a Hamiltonian computed from the functions promoted as Hamiltonians.
Example
julia> f = (x, p) -> x[2]^2 + 2x[1]^2 + p[1]^2
julia> g = (x, p) -> 3x[2]^2 - x[1]^2 + p[2]^2 + p[1]
julia> Poisson(f, g)([1, 2], [2, 1])     # -20
julia> f = (t, x, p, v) -> t*v[1]*x[2]^2 + 2x[1]^2 + p[1]^2 + v[2]
julia> g = (t, x, p, v) -> 3x[2]^2 - x[1]^2 + p[2]^2 + p[1] + t - v[2]
julia> Poisson(f, g, autonomous=false, variable=true)(2, [1, 2], [2, 1], [4, 4])     # -76CTFlows.Poisson — MethodPoisson(
    f::CTFlows.AbstractHamiltonian{TD<:CTFlows.TimeDependence, VD<:CTFlows.VariableDependence},
    g::Function
) -> CTFlows.Hamiltonian
Poisson bracket of a Hamiltonian and a function.
Returns a Hamiltonian representing {f, g} where f is already a Hamiltonian.
Example
julia> f = (x, p) -> x[2]^2 + 2x[1]^2 + p[1]^2
julia> g = (x, p) -> 3x[2]^2 - x[1]^2 + p[2]^2 + p[1]
julia> F = Hamiltonian(f)
julia> Poisson(F, g)([1, 2], [2, 1])     # -20
julia> f = (t, x, p, v) -> t*v[1]*x[2]^2 + 2x[1]^2 + p[1]^2 + v[2]
julia> g = (t, x, p, v) -> 3x[2]^2 - x[1]^2 + p[2]^2 + p[1] + t - v[2]
julia> F = Hamiltonian(f, autonomous=false, variable=true)
julia> Poisson(F, g)(2, [1, 2], [2, 1], [4, 4])     # -76CTFlows.Poisson — MethodPoisson(
    f::Function,
    g::CTFlows.AbstractHamiltonian{TD<:CTFlows.TimeDependence, VD<:CTFlows.VariableDependence}
) -> CTFlows.Hamiltonian
Poisson bracket of a function and a Hamiltonian.
Returns a Hamiltonian representing {f, g} where g is already a Hamiltonian.
Example
julia> f = (x, p) -> x[2]^2 + 2x[1]^2 + p[1]^2
julia> g = (x, p) -> 3x[2]^2 - x[1]^2 + p[2]^2 + p[1]
julia> G = Hamiltonian(g)
julia> Poisson(f, G)([1, 2], [2, 1])     # -20
julia> f = (t, x, p, v) -> t*v[1]*x[2]^2 + 2x[1]^2 + p[1]^2 + v[2]
julia> g = (t, x, p, v) -> 3x[2]^2 - x[1]^2 + p[2]^2 + p[1] + t - v[2]
julia> G = Hamiltonian(g, autonomous=false, variable=true)
julia> Poisson(f, G)(2, [1, 2], [2, 1], [4, 4])     # -76CTFlows.Poisson — MethodPoisson(
    f::CTFlows.AbstractHamiltonian{CTFlows.Autonomous, V<:CTFlows.VariableDependence},
    g::CTFlows.AbstractHamiltonian{CTFlows.Autonomous, V<:CTFlows.VariableDependence}
) -> Any
Poisson bracket of two Hamiltonian functions (subtype of AbstractHamiltonian). Autonomous case.
Returns a Hamiltonian representing the Poisson bracket {f, g} of two autonomous Hamiltonian functions f and g.
Example
julia> f = (x, p) -> x[2]^2 + 2x[1]^2 + p[1]^2
julia> g = (x, p) -> 3x[2]^2 - x[1]^2 + p[2]^2 + p[1]
julia> F = Hamiltonian(f)
julia> G = Hamiltonian(g)
julia> Poisson(f, g)([1, 2], [2, 1])     # -20
julia> Poisson(f, G)([1, 2], [2, 1])     # -20
julia> Poisson(F, g)([1, 2], [2, 1])     # -20CTFlows.Poisson — MethodPoisson(
    f::CTFlows.AbstractHamiltonian{CTFlows.NonAutonomous, V<:CTFlows.VariableDependence},
    g::CTFlows.AbstractHamiltonian{CTFlows.NonAutonomous, V<:CTFlows.VariableDependence}
) -> Any
Poisson bracket of two Hamiltonian functions. Non-autonomous case.
Returns a Hamiltonian representing {f, g} where f and g are time-dependent.
Example
julia> f = (t, x, p, v) -> t*v[1]*x[2]^2 + 2x[1]^2 + p[1]^2 + v[2]
julia> g = (t, x, p, v) -> 3x[2]^2 - x[1]^2 + p[2]^2 + p[1] + t - v[2]
julia> F = Hamiltonian(f, autonomous=false, variable=true)
julia> G = Hamiltonian(g, autonomous=false, variable=true)
julia> Poisson(F, G)(2, [1, 2], [2, 1], [4, 4])     # -76
julia> Poisson(f, g, NonAutonomous, NonFixed)(2, [1, 2], [2, 1], [4, 4])     # -76CTFlows.Poisson — MethodPoisson(
    f::CTFlows.HamiltonianLift{T<:CTFlows.TimeDependence, V<:CTFlows.VariableDependence},
    g::CTFlows.HamiltonianLift{T<:CTFlows.TimeDependence, V<:CTFlows.VariableDependence}
)
Poisson bracket of two HamiltonianLift vector fields.
Returns the HamiltonianLift corresponding to the Lie bracket of vector fields f.X and g.X.
Example
julia> f = x -> [x[1]^2 + x[2]^2, 2x[1]^2]
julia> g = x -> [3x[2]^2, x[2] - x[1]^2]
julia> F = Lift(f)
julia> G = Lift(g)
julia> Poisson(F, G)([1, 2], [2, 1])     # -64
julia> f = (t, x, v) -> [t*v[1]*x[2]^2, 2x[1]^2 + v[2]]
julia> g = (t, x, v) -> [3x[2]^2 - x[1]^2, t - v[2]]
julia> F = Lift(f, NonAutonomous, NonFixed)
julia> G = Lift(g, NonAutonomous, NonFixed)
julia> Poisson(F, G)(2, [1, 2], [2, 1], [4, 4])     # 100CTFlows.∂ₜ — MethodPartial derivative with respect to time of a function.
Example:
julia> ∂ₜ((t,x) -> t*x)(0,8)
8CTFlows.@Lie — MacroCompute Lie or Poisson brackets.
This macro provides a unified notation to define recursively nested Lie brackets (for vector fields) or Poisson brackets (for Hamiltonians).
Syntax
- @Lie [F, G]: computes the Lie bracket- [F, G]of two vector fields.
- @Lie [[F, G], H]: supports arbitrarily nested Lie brackets.
- @Lie {H, K}: computes the Poisson bracket- {H, K}of two Hamiltonians.
- @Lie {{H, K}, L}: supports arbitrarily nested Poisson brackets.
- @Lie expr autonomous = false: specifies a non-autonomous system.
- @Lie expr variable = true: indicates presence of an auxiliary variable- v.
Keyword-like arguments can be provided to control the evaluation context for Poisson brackets with raw functions:
- autonomous = Bool: whether the system is time-independent (default:- true).
- variable = Bool: whether the system depends on an extra variable- v(default:- false).
Bracket type detection
- Square brackets [...]denote Lie brackets betweenVectorFieldobjects.
- Curly brackets {...}denote Poisson brackets betweenHamiltonianobjects or between raw functions.
- The macro automatically dispatches to LieorPoissondepending on the input pattern.
Return
A callable object representing the specified Lie or Poisson bracket expression. The returned function can be evaluated like any other vector field or Hamiltonian.
Examples
■ Lie brackets with VectorField (autonomous)
julia> F1 = VectorField(x -> [0, -x[3], x[2]])
julia> F2 = VectorField(x -> [x[3], 0, -x[1]])
julia> L = @Lie [F1, F2]
julia> L([1.0, 2.0, 3.0])
3-element Vector{Float64}:
  2.0
 -1.0
  0.0■ Lie brackets with VectorField (non-autonomous, with auxiliary variable)
julia> F1 = VectorField((t, x, v) -> [0, -x[3], x[2]]; autonomous=false, variable=true)
julia> F2 = VectorField((t, x, v) -> [x[3], 0, -x[1]]; autonomous=false, variable=true)
julia> L = @Lie [F1, F2]
julia> L(0.0, [1.0, 2.0, 3.0], 1.0)
3-element Vector{Float64}:
  2.0
 -1.0
  0.0■ Poisson brackets with Hamiltonian (autonomous)
julia> H1 = Hamiltonian((x, p) -> x[1]^2 + p[2]^2)
julia> H2 = Hamiltonian((x, p) -> x[2]^2 + p[1]^2)
julia> P = @Lie {H1, H2}
julia> P([1.0, 1.0], [3.0, 2.0])
-4.0■ Poisson brackets with Hamiltonian (non-autonomous, with variable)
julia> H1 = Hamiltonian((t, x, p, v) -> x[1]^2 + p[2]^2 + v; autonomous=false, variable=true)
julia> H2 = Hamiltonian((t, x, p, v) -> x[2]^2 + p[1]^2 + v; autonomous=false, variable=true)
julia> P = @Lie {H1, H2}
julia> P(1.0, [1.0, 3.0], [4.0, 2.0], 3.0)
8.0■ Poisson brackets from raw functions
julia> H1 = (x, p) -> x[1]^2 + p[2]^2
julia> H2 = (x, p) -> x[2]^2 + p[1]^2
julia> P = @Lie {H1, H2}
julia> P([1.0, 1.0], [3.0, 2.0])
-4.0■ Poisson bracket with non-autonomous raw functions
julia> H1 = (t, x, p) -> x[1]^2 + p[2]^2 + t
julia> H2 = (t, x, p) -> x[2]^2 + p[1]^2 + t
julia> P = @Lie {H1, H2} autonomous = false
julia> P(3.0, [1.0, 2.0], [4.0, 1.0])
-8.0■ Nested brackets
julia> F = VectorField(x -> [-x[1], x[2], x[3]])
julia> G = VectorField(x -> [x[3], -x[2], 0])
julia> H = VectorField(x -> [0, 0, -x[1]])
julia> nested = @Lie [[F, G], H]
julia> nested([1.0, 2.0, 3.0])
3-element Vector{Float64}:
  2.0
  0.0
 -6.0julia> H1 = (x, p) -> x[2]*x[1]^2 + p[1]^2
julia> H2 = (x, p) -> x[1]*p[2]^2
julia> H3 = (x, p) -> x[1]*p[2] + x[2]*p[1]
julia> nested_poisson = @Lie {{H1, H2}, H3}
julia> nested_poisson([1.0, 2.0], [0.5, 1.0])
14.0■ Mixed expressions with arithmetic
julia> F1 = VectorField(x -> [0, -x[3], x[2]])
julia> F2 = VectorField(x -> [x[3], 0, -x[1]])
julia> x = [1.0, 2.0, 3.0]
julia> @Lie [F1, F2](x) + 3 * [F1, F2](x)
3-element Vector{Float64}:
  8.0
 -4.0
  0.0julia> H1 = (x, p) -> x[1]^2
julia> H2 = (x, p) -> p[1]^2
julia> H3 = (x, p) -> x[1]*p[1]
julia> x = [1.0, 2.0, 3.0]
julia> p = [3.0, 2.0, 1.0]
julia> @Lie {H1, H2}(x, p) + 2 * {H2, H3}(x, p)
24.0