Differential Geometry
Index
CTFlows.:⅋CTFlows.:⅋CTFlows.:⋅CTFlows.:⋅CTFlows.:⋅CTFlows.LieCTFlows.LieCTFlows.LieCTFlows.LieCTFlows.LiftCTFlows.LiftCTFlows.PoissonCTFlows.PoissonCTFlows.PoissonCTFlows.PoissonCTFlows.PoissonCTFlows.PoissonCTFlows._Lie_bracketCTFlows._Lie_bracketCTFlows._Lie_bracketCTFlows._Lie_bracketCTFlows._Lie_bracketCTFlows.__check_bracket_consistencyCTFlows.__is_mixed_usageCTFlows.__parse_lie_argsCTFlows.__transform_lie_poisson_expressionCTFlows._get_TDCTFlows._get_TDCTFlows._get_TDCTFlows._get_VDCTFlows._get_VDCTFlows._get_VDCTFlows.∂ₜ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.:⋅ — Method
Lie 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.:⋅ — Method
Lie 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.:⋅ — Method
Lie 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 — Method
Lie 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 — Method
Lie 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 — Method
Lie 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 — Method
Lie 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 — Method
Lift(
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 — Method
Lift(
X::Function;
autonomous,
variable
) -> CTFlows.var"#21#22"{<: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)ifautonomous=trueandvariable=false(x, p, v)ifautonomous=trueandvariable=true(t, x, p)ifautonomous=falseandvariable=false(t, x, p, v)ifautonomous=falseandvariable=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 — Method
Poisson(
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 — Method
Poisson(
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 — Method
Poisson(
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 — Method
Poisson(
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 — Method
Poisson(
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 — Method
Poisson(
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._Lie_bracket — Method
_Lie_bracket(
f::Function,
g::Function;
autonomous,
variable
) -> Any
Internal Lie bracket for two plain functions.
Automatically wraps both functions in VectorField objects with the specified time and variable dependence.
Arguments
f::Function: First function representing a vector fieldg::Function: Second function representing a vector fieldautonomous::Bool=true: Whether the functions are time-independentvariable::Bool=false: Whether the functions depend on an auxiliary variablev
Returns
- A
VectorFieldrepresenting the Lie bracket[f, g]
Example
julia> f = x -> [x[2], -x[1]]
julia> g = x -> [x[1], x[2]]
julia> bracket = CTFlows._Lie_bracket(f, g)
julia> bracket([1, 2])
2-element Vector{Float64}:
3.0
-3.0CTFlows._Lie_bracket — Method
_Lie_bracket(
f::CTFlows.AbstractVectorField{TD<:CTFlows.TimeDependence, VD<:CTFlows.VariableDependence},
g::Function
) -> Any
Internal Lie bracket for a vector field and a plain function.
Automatically wraps the function in a VectorField with time and variable dependence inferred from the first argument.
Arguments
f::AbstractVectorField{TD,VD}: Vector fieldg::Function: Function representing a vector field
Returns
- A
VectorFieldrepresenting the Lie bracket[f, g]
CTFlows._Lie_bracket — Method
_Lie_bracket(
f::Function,
g::CTFlows.AbstractVectorField{TD<:CTFlows.TimeDependence, VD<:CTFlows.VariableDependence}
) -> Any
Internal Lie bracket for a plain function and a vector field.
Automatically wraps the function in a VectorField with time and variable dependence inferred from the second argument.
Arguments
f::Function: Function representing a vector fieldg::AbstractVectorField{TD,VD}: Vector field
Returns
- A
VectorFieldrepresenting the Lie bracket[f, g]
CTFlows._Lie_bracket — Method
_Lie_bracket(
X::CTFlows.VectorField{<:Function, CTFlows.Autonomous, V<:CTFlows.VariableDependence},
Y::CTFlows.VectorField{<:Function, CTFlows.Autonomous, V<:CTFlows.VariableDependence}
) -> Any
Internal Lie bracket for two vector fields in the autonomous case.
This is the core implementation used by both Lie and the @Lie macro.
Arguments
X::VectorField{<:Function,Autonomous,V}: First vector fieldY::VectorField{<:Function,Autonomous,V}: Second vector field
Returns
VectorField{<:Function,Autonomous,V}: The Lie bracket[X, Y]
CTFlows._Lie_bracket — Method
_Lie_bracket(
X::CTFlows.VectorField{<:Function, CTFlows.NonAutonomous, V<:CTFlows.VariableDependence},
Y::CTFlows.VectorField{<:Function, CTFlows.NonAutonomous, V<:CTFlows.VariableDependence}
) -> Any
Internal Lie bracket for two vector fields in the non-autonomous case.
This is the core implementation used by both Lie and the @Lie macro.
Arguments
X::VectorField{<:Function,NonAutonomous,V}: First vector fieldY::VectorField{<:Function,NonAutonomous,V}: Second vector field
Returns
VectorField{<:Function,NonAutonomous,V}: The Lie bracket[X, Y]
CTFlows.__check_bracket_consistency — Method
__check_bracket_consistency(
a,
b,
has_autonomous::Bool,
autonomous::Bool,
has_variable::Bool,
variable::Bool
)
Check consistency between bracket operands and user-specified autonomous/variable arguments.
This function validates that:
- Both operands have matching
TimeDependence(if both are typed objects) - Both operands have matching
VariableDependence(if both are typed objects) - User-specified
autonomousmatches operands'TimeDependence(if provided) - User-specified
variablematches operands'VariableDependence(if provided)
Arguments
a: First bracket operand (Function, AbstractVectorField, or AbstractHamiltonian)b: Second bracket operand (Function, AbstractVectorField, or AbstractHamiltonian)has_autonomous::Bool: Whether user explicitly providedautonomousargumentautonomous::Bool: User-specified autonomous valuehas_variable::Bool: Whether user explicitly providedvariableargumentvariable::Bool: User-specified variable value
Throws
CTBase.Exceptions.IncorrectArgument: If any consistency check fails
Examples
julia> F1 = VectorField(x -> [x[2], -x[1]]) # Autonomous, Fixed
julia> F2 = VectorField((t, x) -> [x[2], -x[1]]; autonomous=false) # NonAutonomous, Fixed
julia> CTFlows.__check_bracket_consistency(F1, F2, false, true, false, false)
ERROR: IncorrectArgument: Mismatched time dependence: first operand is Autonomous, second is NonAutonomous.
julia> CTFlows.__check_bracket_consistency(F1, F1, true, false, false, false)
ERROR: IncorrectArgument: autonomous=false conflicts with first operand which is Autonomous.
julia> CTFlows.__check_bracket_consistency(F1, F1, false, true, false, false) # OK, no errorCTFlows.__is_mixed_usage — Method
__is_mixed_usage(x) -> Bool
Helper function to detect mixed usage of Lie and Poisson brackets in an expression.
This function walks through an expression tree and checks whether it contains both:
- Lie bracket notation
[_, _](square brackets) - Poisson bracket notation
{_, _}(curly braces)
Mixing these two types of brackets in the same expression is not allowed and will trigger an error in the @Lie macro.
Arguments
x::Expr: An expression to analyze
Returns
Bool:trueif the expression contains both Lie and Poisson brackets,falseotherwise
Examples
julia> CTFlows.__is_mixed_usage(:([F0, F1]))
false
julia> CTFlows.__is_mixed_usage(:({H0, H1}))
false
julia> CTFlows.__is_mixed_usage(:([F0, F1] + {H0, H1}))
trueCTFlows.__parse_lie_args — Method
__parse_lie_args(
args...
) -> Tuple{Any, Any, Bool, Bool, Union{Nothing, Expr}}
Helper function to parse keyword arguments for the @Lie macro.
This function parses the keyword arguments and returns the autonomous and variable flags, along with flags indicating whether each argument was explicitly provided by the user. If an invalid argument is encountered, it returns an error expression that can be thrown.
Arguments
args...: Variable arguments passed to the macro
Returns
Tuple{Bool, Bool, Bool, Bool, Union{Expr, Nothing}}:autonomous: Whether the system is time-independentvariable: Whether the system depends on an extra variablevhas_autonomous: Whether the user explicitly provided theautonomousargumenthas_variable: Whether the user explicitly provided thevariableargumenterror_expr: An expression to throw an error if invalid arguments are found, ornothingif parsing succeeded
Examples
julia> autonomous, variable, has_autonomous, has_variable, error = CTFlows.__parse_lie_args()
(true, false, false, false, nothing)
julia> autonomous, variable, has_autonomous, has_variable, error = CTFlows.__parse_lie_args(:(autonomous = false))
(false, false, true, false, nothing)
julia> autonomous, variable, has_autonomous, has_variable, error = CTFlows.__parse_lie_args(:(autonomous = false), :(variable = true))
(false, true, true, true, nothing)
julia> autonomous, variable, has_autonomous, has_variable, error = CTFlows.__parse_lie_args(:(invalid_arg = true))
(true, false, false, false, :(throw(CTBase.Exceptions.IncorrectArgument("Invalid argument: invalid_arg = true"))))CTFlows.__transform_lie_poisson_expression — Method
__transform_lie_poisson_expression(
x,
autonomous,
variable,
has_autonomous,
has_variable
) -> Any
Helper function to transform Lie and Poisson bracket expressions into their corresponding function calls.
This function walks through an expression and transforms:
- Lie bracket notation
[a, b]intoCTFlows._Lie_bracket(a, b)calls with runtime type checks - Poisson bracket notation
{c, d}intoCTFlows.Poisson(c, d)calls with runtime type checks
For both Lie and Poisson brackets with raw functions, it generates runtime type checks to determine whether to pass autonomous and variable keyword arguments.
Arguments
x::Expr: An expression to transformautonomous::Bool: Whether the system is time-independentvariable::Bool: Whether the system depends on an extra variablev
Returns
Expr: The transformed expression with_Lie_bracketorPoissonfunction calls
Examples
julia> expr = :([F0, F1])
julia> result = CTFlows.__transform_lie_poisson_expression(expr, true, false)
quote
if isa(F0, Function) && isa(F1, Function)
CTFlows._Lie_bracket(F0, F1; autonomous=true, variable=false)
else
CTFlows._Lie_bracket(F0, F1)
end
end
julia> expr = :({f, g})
julia> result = CTFlows.__transform_lie_poisson_expression(expr, true, false)
quote
if isa(f, Function) && isa(g, Function)
CTFlows.Poisson(f, g; autonomous=true, variable=false)
else
CTFlows.Poisson(f, g)
end
endCTFlows._get_TD — Method
_get_TD(_::Function)
Extract the TimeDependence type from a plain Function.
Returns nothing since plain functions don't have type parameters.
CTFlows._get_TD — Method
_get_TD(_::CTFlows.AbstractHamiltonian{TD, VD}) -> Any
Extract the TimeDependence type from an AbstractHamiltonian.
Returns the type parameter TD (either Autonomous or NonAutonomous).
CTFlows._get_TD — Method
_get_TD(_::CTFlows.AbstractVectorField{TD, VD}) -> Any
Extract the TimeDependence type from an AbstractVectorField.
Returns the type parameter TD (either Autonomous or NonAutonomous).
CTFlows._get_VD — Method
_get_VD(_::Function)
Extract the VariableDependence type from a plain Function.
Returns nothing since plain functions don't have type parameters.
CTFlows._get_VD — Method
_get_VD(_::CTFlows.AbstractHamiltonian{TD, VD}) -> Any
Extract the VariableDependence type from an AbstractHamiltonian.
Returns the type parameter VD (either Fixed or NonFixed).
CTFlows._get_VD — Method
_get_VD(_::CTFlows.AbstractVectorField{TD, VD}) -> Any
Extract the VariableDependence type from an AbstractVectorField.
Returns the type parameter VD (either Fixed or NonFixed).
CTFlows.∂ₜ — Method
Partial derivative with respect to time of a function.
Example:
julia> ∂ₜ((t,x) -> t*x)(0,8)
8CTFlows.@Lie — Macro
Compute 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 variablev.
Keyword-like arguments can be provided to control the evaluation context for both Lie and Poisson brackets with raw functions:
autonomous = Bool: whether the system is time-independent (default:true).variable = Bool: whether the system depends on an extra variablev(default:false).
Bracket type detection
- Square brackets
[...]denote Lie brackets betweenVectorFieldobjects or raw functions. - Curly brackets
{...}denote Poisson brackets betweenHamiltonianobjects or raw functions. - The macro automatically dispatches to
_Lie_bracketorPoissondepending on the input pattern. - Raw functions are automatically wrapped in
VectorFieldorHamiltonianobjects as needed.
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■ Lie brackets from raw functions (autonomous)
julia> f1 = x -> [0, -x[3], x[2]]
julia> f2 = 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 bracket with non-autonomous raw functions
julia> f1 = (t, x) -> [t, -x[3], x[2]]
julia> f2 = (t, x) -> [x[3], 0, -x[1]]
julia> L = @Lie [f1, f2] autonomous = false
julia> L(1.0, [1.0, 2.0, 3.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