Differential geometry
Index
CTBase.:⋅
CTBase.:⋅
CTBase.:⋅
CTBase.Lie
CTBase.Lie
CTBase.Lie
CTBase.Lie
CTBase.Lie
CTBase.Lift
CTBase.Lift
CTBase.Lift
CTBase.Poisson
CTBase.Poisson
CTBase.Poisson
CTBase.Poisson
CTBase.Poisson
CTBase.Poisson
CTBase.Poisson
CTBase.∂ₜ
CTBase.@Lie
CTBase.@Lie
CTBase.@Lie
Documentation
CTBase.:⋅
— Method⋅(
X::Function,
f::Function
) -> CTBase.var"#90#92"{VectorField{Autonomous, Fixed}, <:Function}
Lie derivative of a scalar function along a function. In this case both functions will be 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])
MethodError
CTBase.:⋅
— Method⋅(
X::VectorField{Autonomous, <:VariableDependence},
f::Function
) -> CTBase.var"#90#92"{VectorField{Autonomous, var"#s205"}, <:Function} where var"#s205"<:VariableDependence
Lie derivative of a scalar function along a vector field : L_X(f) = X⋅f, in 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])
0
CTBase.:⋅
— Method⋅(
X::VectorField{NonAutonomous, <:VariableDependence},
f::Function
) -> CTBase.var"#94#96"{VectorField{NonAutonomous, var"#s205"}, <:Function} where var"#s205"<:VariableDependence
Lie derivative of a scalar function along a vector field : L_X(f) = X⋅f, in 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])
10
CTBase.Lie
— MethodLie(
X::Function,
f::Function,
dependences::DataType...
) -> Union{CTBase.var"#90#92"{VectorField{Autonomous, var"#s205"}, <:Function} where var"#s205"<:VariableDependence, CTBase.var"#94#96"{VectorField{NonAutonomous, var"#s205"}, <:Function} where var"#s205"<:VariableDependence}
Lie derivative of a scalar function along a vector field or a function. Dependencies are specified with DataType : Autonomous, NonAutonomous and Fixed, NonFixed.
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, NonAutonomous, NonFixed)(1, [1, 2], [2, 1])
10
CTBase.Lie
— MethodLie(
X::Function,
f::Function;
autonomous,
variable
) -> Union{CTBase.var"#90#92"{VectorField{Autonomous, var"#s205"}, <:Function} where var"#s205"<:VariableDependence, CTBase.var"#94#96"{VectorField{NonAutonomous, var"#s205"}, <:Function} where var"#s205"<:VariableDependence}
Lie derivative of a scalar function along a function. Dependencies are specified with boolean : autonomous and variable.
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])
10
CTBase.Lie
— MethodLie(
X::VectorField,
f::Function
) -> Union{CTBase.var"#90#92"{VectorField{Autonomous, var"#s205"}, <:Function} where var"#s205"<:VariableDependence, CTBase.var"#94#96"{VectorField{NonAutonomous, var"#s205"}, <:Function} where var"#s205"<:VariableDependence}
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])
10
CTBase.Lie
— MethodLie(
X::VectorField{Autonomous, V<:VariableDependence},
Y::VectorField{Autonomous, V<:VariableDependence}
) -> VectorField
Lie bracket of two vector fields: [X, Y] = Lie(X, Y), 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]
CTBase.Lie
— MethodLie(
X::VectorField{NonAutonomous, V<:VariableDependence},
Y::VectorField{NonAutonomous, V<:VariableDependence}
) -> VectorField{NonAutonomous}
Lie bracket of two vector fields: [X, Y] = Lie(X, Y), 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]
CTBase.Lift
— MethodLift(
X::Function,
dependences::DataType...
) -> CTBase.var"#82#86"{<:Function}
Return the Lift of a function. Dependencies are specified with DataType : Autonomous, NonAutonomous and Fixed, NonFixed.
Example
julia> H = Lift(x -> 2x)
julia> H(1, 1)
2
julia> H = Lift((t, x, v) -> 2x + t - v, NonAutonomous, NonFixed)
julia> H(1, 1, 1, 1)
2
CTBase.Lift
— MethodLift(
X::Function;
autonomous,
variable
) -> CTBase.var"#82#86"{<:Function}
Return the Lift of a function. Dependencies are specified with boolean : autonomous and variable.
Example
julia> H = Lift(x -> 2x)
julia> H(1, 1)
2
julia> H = Lift((t, x, v) -> 2x + t - v, autonomous=false, variable=true)
julia> H(1, 1, 1, 1)
2
CTBase.Lift
— MethodLift(X::VectorField) -> HamiltonianLift
Return the HamiltonianLift of a VectorField.
Example
julia> HL = Lift(VectorField(x -> [x[1]^2,x[2]^2], autonomous=true, variable=false))
julia> HL([1, 0], [0, 1])
0
julia> HL = Lift(VectorField((t, x, v) -> [t+x[1]^2,x[2]^2+v], autonomous=false, variable=true))
julia> HL(1, [1, 0], [0, 1], 1)
1
julia> H = Lift(x -> 2x)
julia> H(1, 1)
2
julia> H = Lift((t, x, v) -> 2x + t - v, autonomous=false, variable=true)
julia> H(1, 1, 1, 1)
2
julia> H = Lift((t, x, v) -> 2x + t - v, NonAutonomous, NonFixed)
julia> H(1, 1, 1, 1)
2
CTBase.Poisson
— MethodPoisson(
f::Function,
g::Function,
dependences::DataType...
) -> Hamiltonian
Poisson bracket of two functions : {f, g} = Poisson(f, g) Dependencies are specified with DataType : Autonomous, NonAutonomous and Fixed, NonFixed.
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, NonAutonomous, NonFixed)(2, [1, 2], [2, 1], [4, 4])
-76
CTBase.Poisson
— MethodPoisson(
f::Function,
g::Function;
autonomous,
variable
) -> Hamiltonian
Poisson bracket of two functions : {f, g} = Poisson(f, g) Dependencies are specified with boolean : autonomous and variable.
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])
-76
CTBase.Poisson
— MethodPoisson(
f::AbstractHamiltonian{Autonomous, V<:VariableDependence},
g::AbstractHamiltonian{Autonomous, V<:VariableDependence}
) -> HamiltonianLift
Poisson bracket of two Hamiltonian functions (subtype of AbstractHamiltonian) : {f, g} = Poisson(f, g), autonomous case
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])
-20
CTBase.Poisson
— MethodPoisson(
f::AbstractHamiltonian{NonAutonomous, V<:VariableDependence},
g::AbstractHamiltonian{NonAutonomous, V<:VariableDependence}
) -> HamiltonianLift
Poisson bracket of two Hamiltonian functions (subtype of AbstractHamiltonian) : {f, g} = Poisson(f, g), non autonomous case
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])
-76
CTBase.Poisson
— MethodPoisson(
f::AbstractHamiltonian{T<:TimeDependence, V<:VariableDependence},
g::Function
) -> Hamiltonian
Poisson bracket of an Hamiltonian function (subtype of AbstractHamiltonian) and a function : {f, g} = Poisson(f, g), autonomous case
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])
-76
CTBase.Poisson
— MethodPoisson(
f::Function,
g::AbstractHamiltonian{T<:TimeDependence, V<:VariableDependence}
) -> Hamiltonian
Poisson bracket of a function and an Hamiltonian function (subtype of AbstractHamiltonian) : {f, g} = Poisson(f, 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> 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])
-76
CTBase.Poisson
— MethodPoisson(
f::HamiltonianLift{T<:TimeDependence, V<:VariableDependence},
g::HamiltonianLift{T<:TimeDependence, V<:VariableDependence}
) -> HamiltonianLift
Poisson bracket of two HamiltonianLift functions : {f, g} = Poisson(f, g)
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])
100
CTBase.∂ₜ
— Method∂ₜ(f) -> CTBase.var"#99#101"
Partial derivative wrt time of a function.
Example
julia> ∂ₜ((t,x) -> t*x)(0,8)
8
CTBase.@Lie
— MacroMacros for Poisson brackets
Example
julia> H0 = (x, p) -> 0.5*(x[1]^2+x[2]^2+p[1]^2)
julia> H1 = (x, p) -> 0.5*(x[1]^2+x[2]^2+p[2]^2)
julia> @Lie {H0, H1}([1, 2, 3], [1, 0, 7]) autonomous=true variable=false
#
julia> H0 = (t, x, p) -> 0.5*(x[1]^2+x[2]^2+p[1]^2)
julia> H1 = (t, x, p) -> 0.5*(x[1]^2+x[2]^2+p[2]^2)
julia> @Lie {H0, H1}(1, [1, 2, 3], [1, 0, 7]) autonomous=false variable=false
#
julia> H0 = (x, p, v) -> 0.5*(x[1]^2+x[2]^2+p[1]^2+v)
julia> H1 = (x, p, v) -> 0.5*(x[1]^2+x[2]^2+p[2]^2+v)
julia> @Lie {H0, H1}([1, 2, 3], [1, 0, 7], 2) autonomous=true variable=true
#
julia> H0 = (t, x, p, v) -> 0.5*(x[1]^2+x[2]^2+p[1]^2+v)
julia> H1 = (t, x, p, v) -> 0.5*(x[1]^2+x[2]^2+p[2]^2+v)
julia> @Lie {H0, H1}(1, [1, 2, 3], [1, 0, 7], 2) autonomous=false variable=true
CTBase.@Lie
— MacroMacros for Lie and Poisson brackets
Example
julia> H0 = (t, x, p) -> 0.5*(x[1]^2+x[2]^2+p[1]^2)
julia> H1 = (t, x, p) -> 0.5*(x[1]^2+x[2]^2+p[2]^2)
julia> @Lie {H0, H1}(1, [1, 2, 3], [1, 0, 7]) autonomous=false
#
julia> H0 = (x, p, v) -> 0.5*(x[1]^2+x[2]^2+p[1]^2+v)
julia> H1 = (x, p, v) -> 0.5*(x[1]^2+x[2]^2+p[2]^2+v)
julia> @Lie {H0, H1}([1, 2, 3], [1, 0, 7], 2) variable=true
#
CTBase.@Lie
— MacroMacros for Lie and Poisson brackets
Example
julia> F0 = VectorField(x -> [x[1], x[2], (1-x[3])])
julia> F1 = VectorField(x -> [0, -x[3], x[2]])
julia> @Lie [F0, F1]([1, 2, 3])
[0, 5, 4]
#
julia> F0 = VectorField((t, x) -> [t+x[1], x[2], (1-x[3])], autonomous=false)
julia> F1 = VectorField((t, x) -> [t, -x[3], x[2]], autonomous=false)
julia> @Lie [F0, F1](1, [1, 2, 3])
#
julia> F0 = VectorField((x, v) -> [x[1]+v, x[2], (1-x[3])], variable=true)
julia> F1 = VectorField((x, v) -> [0, -x[3]-v, x[2]], variable=true)
julia> @Lie [F0, F1]([1, 2, 3], 2)
#
julia> F0 = VectorField((t, x, v) -> [t+x[1]+v, x[2], (1-x[3])], autonomous=false, variable=true)
julia> F1 = VectorField((t, x, v) -> [t, -x[3]-v, x[2]], autonomous=false, variable=true)
julia> @Lie [F0, F1](1, [1, 2, 3], 2)
#
julia> H0 = Hamiltonian((x, p) -> 0.5*(2x[1]^2+x[2]^2+p[1]^2))
julia> H1 = Hamiltonian((x, p) -> 0.5*(3x[1]^2+x[2]^2+p[2]^2))
julia> @Lie {H0, H1}([1, 2, 3], [1, 0, 7])
3.0
#
julia> H0 = Hamiltonian((t, x, p) -> 0.5*(x[1]^2+x[2]^2+p[1]^2), autonomous=false)
julia> H1 = Hamiltonian((t, x, p) -> 0.5*(x[1]^2+x[2]^2+p[2]^2), autonomous=false)
julia> @Lie {H0, H1}(1, [1, 2, 3], [1, 0, 7])
#
julia> H0 = Hamiltonian((x, p, v) -> 0.5*(x[1]^2+x[2]^2+p[1]^2+v), variable=true)
julia> H1 = Hamiltonian((x, p, v) -> 0.5*(x[1]^2+x[2]^2+p[2]^2+v), variable=true)
julia> @Lie {H0, H1}([1, 2, 3], [1, 0, 7], 2)
#
julia> H0 = Hamiltonian((t, x, p, v) -> 0.5*(x[1]^2+x[2]^2+p[1]^2+v), autonomous=false, variable=true)
julia> H1 = Hamiltonian((t, x, p, v) -> 0.5*(x[1]^2+x[2]^2+p[2]^2+v), autonomous=false, variable=true)
julia> @Lie {H0, H1}(1, [1, 2, 3], [1, 0, 7], 2)
#