Lagrangian Systems
The Euler-Lagrange equations, that is the dynamical equations of a Lagrangian system, are given in terms of the Lagrangian $L(x,v)$ by
\[\frac{d}{dt} \frac{\partial L}{\partial v} - \frac{\partial L}{\partial x} = 0 .\]
For regular (i.e. non-degenerate) Lagrangians, this is a set of second-order ordinary differential equations. In many numerical applications, it is advantageous to solve the implicit form of these equations, given by
\[\begin{align*} \frac{d \vartheta}{dt} &= f , & \vartheta &= \frac{\partial L}{\partial v} , & f = \frac{\partial L}{\partial x} . \end{align*}\]
In the following, we show how these equations can be obtained for the example of a particle in a square potential.
Particle in a potential
Before any use, we need to load EulerLagrange
:
using EulerLagrange
Next, we generate symbolic variables for a two-dimensional system:
t, x, v = lagrangian_variables(2)
(t, (x(t))[1:2], (v(t))[1:2])
With those variables, we can construct a Lagrangian
using LinearAlgebra
L = v ⋅ v / 2 - x ⋅ x / 2
\[ \begin{equation} \frac{1}{2} \left( v\left( t \right)_{1}^{2} + v\left( t \right)_{2}^{2} \right) - \frac{1}{2} \left( x\left( t \right)_{1}^{2} + x\left( t \right)_{2}^{2} \right) \end{equation} \]
This Lagrangian together with the symbolic variables is then used to construct a LagrangianSystem
:
lag_sys = LagrangianSystem(L, t, x, v)
Lagrangian system with
L = (1//2)*((v(t))[1]^2 + (v(t))[2]^2) - (1//2)*((x(t))[1]^2 + (x(t))[2]^2)
The constructor computes the Euler-Lagrange equations and generates the corresponding Julia code. In the last step, we can now construct a LODEProblem
from the LagrangianSystem
and some appropriate initial conditions, a time span to integrate over and a time step:
tspan = (0.0, 10.0)
tstep = 0.01
q₀ = [1.0, 1.0]
p₀ = [0.5, 2.0]
lprob = LODEProblem(lag_sys, tspan, tstep, q₀, p₀)
Geometric Equation Problem for Lagrangian Ordinary Differential Equation (LODE)
with vector fields
ϑ = RuntimeGeneratedFunctions.RuntimeGeneratedFunction{(:ˍ₋out, :t, :X, :V, :params), EulerLagrange.var"#_RGF_ModTag", EulerLagrange.var"#_RGF_ModTag", (0x9074dbae, 0x9dbd5857, 0x29f62319, 0xac7fc940, 0x64f9e588), Expr}(quote
#= none:1 =#
#= none:5 =#
begin
#= none:7 =#
#= none:7 =# @inbounds begin
#= none:9 =#
ˍ₋out[1] = getindex(V, 1)
#= none:10 =#
ˍ₋out[2] = getindex(V, 2)
#= none:12 =#
nothing
end
end
end)
f = RuntimeGeneratedFunctions.RuntimeGeneratedFunction{(:ˍ₋out, :t, :X, :V, :params), EulerLagrange.var"#_RGF_ModTag", EulerLagrange.var"#_RGF_ModTag", (0xcb91fc44, 0xcf8fc5e3, 0xe0155e43, 0x302a3ab5, 0xe40c74b9), Expr}(quote
#= none:1 =#
#= none:5 =#
begin
#= none:7 =#
#= none:7 =# @inbounds begin
#= none:9 =#
ˍ₋out[1] = -1 // 1 * getindex(X, 1)
#= none:10 =#
ˍ₋out[2] = -1 // 1 * getindex(X, 2)
#= none:12 =#
nothing
end
end
end)
g = RuntimeGeneratedFunctions.RuntimeGeneratedFunction{(:ˍ₋out, :t, :X, :Λ, :V, :params), EulerLagrange.var"#_RGF_ModTag", EulerLagrange.var"#_RGF_ModTag", (0x6f591718, 0x213c7aee, 0xa402e277, 0xe90cbccb, 0xda71f8c3), Expr}(quote
#= none:1 =#
#= none:5 =#
begin
#= none:7 =#
#= none:7 =# @inbounds begin
#= none:9 =#
ˍ₋out[1] = (Differential(t))(getindex(V, 1))
#= none:10 =#
ˍ₋out[2] = (Differential(t))(getindex(V, 2))
#= none:12 =#
nothing
end
end
end)
Lagrangian: L = RuntimeGeneratedFunctions.RuntimeGeneratedFunction{(:t, :X, :V, :params), EulerLagrange.var"#_RGF_ModTag", EulerLagrange.var"#_RGF_ModTag", (0x79cf434b, 0x0601728f, 0x22d644f8, 0x940a3832, 0x772811c7), Expr}(quote
#= none:1 =#
#= none:5 =#
1 // 2 * (getindex(V, 1) ^ 2 + getindex(V, 2) ^ 2) + -1 // 2 * (getindex(X, 1) ^ 2 + getindex(X, 2) ^ 2)
end)
Invariants:
GeometricBase.NullInvariants()
Timespan: (0.0, 10.0)
Timestep: 0.01
Initial conditions:
(t = 0.0, q = [1.0, 1.0], p = [0.5, 2.0], v = [0.0, 0.0])
Parameters:
GeometricBase.NullParameters()
We can integrate this system using GeometricIntegrators:
using GeometricIntegrators
sol = integrate(lprob, Gauss(1))
using CairoMakie
fig = lines(parent(sol.q[:,1]), parent(sol.q[:,2]);
axis = (; xlabel = "x₁", ylabel = "x₂", title = "Particle moving in a square potential"),
figure = (; size = (800,600), fontsize = 22))
┌ Warning: Hermite Extrapolation: q's history[1] and history[2] are identical!
└ @ GeometricIntegrators.Extrapolators ~/.julia/packages/GeometricIntegrators/xUuFd/src/extrapolation/hermite.jl:211
Parameters
We can also include parametric dependencies in the Lagrangian. Consider, for example, a parameter α
that determines the strength of the potential.
The easiest way, to account for parameters, is to create a named tuple with typical values for each parameter, e.g.,
params = (α = 5.0,)
(α = 5.0,)
In the next step, we use the function symbolize
to generate a symbolic version of the parameters:
sparams = symbolize(params)
(α = αₚ,)
Now we modify the Lagrangian to account for the parameter:
L = v ⋅ v / 2 - sparams.α * (x ⋅ x) / 2
\[ \begin{equation} \frac{1}{2} \left( v\left( t \right)_{1}^{2} + v\left( t \right)_{2}^{2} \right) - \frac{1}{2} \left( x\left( t \right)_{1}^{2} + x\left( t \right)_{2}^{2} \right) \mathtt{\alpha{_p}} \end{equation} \]
From here on, everything follows along the same lines as before, the only difference being that we also need to pass the symbolic parameters sparams
to the LagrangianSystem
constructor:
lag_sys = LagrangianSystem(L, t, x, v, sparams)
Lagrangian system with
L = (1//2)*((v(t))[1]^2 + (v(t))[2]^2) - (1//2)*((x(t))[1]^2 + (x(t))[2]^2)*αₚ
Analogously, we need to pass actual parameter values params
to the LODEProblem
constructor via the parameters
keyword argument:
lprob = LODEProblem(lag_sys, tspan, tstep, q₀, p₀; parameters = params)
Geometric Equation Problem for Lagrangian Ordinary Differential Equation (LODE)
with vector fields
ϑ = RuntimeGeneratedFunctions.RuntimeGeneratedFunction{(:ˍ₋out, :t, :X, :V, :params), EulerLagrange.var"#_RGF_ModTag", EulerLagrange.var"#_RGF_ModTag", (0x9074dbae, 0x9dbd5857, 0x29f62319, 0xac7fc940, 0x64f9e588), Expr}(quote
#= none:1 =#
#= none:5 =#
begin
#= none:7 =#
#= none:7 =# @inbounds begin
#= none:9 =#
ˍ₋out[1] = getindex(V, 1)
#= none:10 =#
ˍ₋out[2] = getindex(V, 2)
#= none:12 =#
nothing
end
end
end)
f = RuntimeGeneratedFunctions.RuntimeGeneratedFunction{(:ˍ₋out, :t, :X, :V, :params), EulerLagrange.var"#_RGF_ModTag", EulerLagrange.var"#_RGF_ModTag", (0x7e9a6a94, 0xb6bea226, 0x391e4a0d, 0x03fd8e57, 0xc319b875), Expr}(quote
#= none:1 =#
#= none:5 =#
begin
#= none:7 =#
#= none:7 =# @inbounds begin
#= none:9 =#
ˍ₋out[1] = (-1 // 1 * getindex(X, 1)) * params.α
#= none:10 =#
ˍ₋out[2] = (-1 // 1 * getindex(X, 2)) * params.α
#= none:12 =#
nothing
end
end
end)
g = RuntimeGeneratedFunctions.RuntimeGeneratedFunction{(:ˍ₋out, :t, :X, :Λ, :V, :params), EulerLagrange.var"#_RGF_ModTag", EulerLagrange.var"#_RGF_ModTag", (0x6f591718, 0x213c7aee, 0xa402e277, 0xe90cbccb, 0xda71f8c3), Expr}(quote
#= none:1 =#
#= none:5 =#
begin
#= none:7 =#
#= none:7 =# @inbounds begin
#= none:9 =#
ˍ₋out[1] = (Differential(t))(getindex(V, 1))
#= none:10 =#
ˍ₋out[2] = (Differential(t))(getindex(V, 2))
#= none:12 =#
nothing
end
end
end)
Lagrangian: L = RuntimeGeneratedFunctions.RuntimeGeneratedFunction{(:t, :X, :V, :params), EulerLagrange.var"#_RGF_ModTag", EulerLagrange.var"#_RGF_ModTag", (0xd47b589b, 0x156de4f8, 0xdb8bbeeb, 0x8817e5e0, 0x6c7ea09c), Expr}(quote
#= none:1 =#
#= none:5 =#
1 // 2 * (getindex(V, 1) ^ 2 + getindex(V, 2) ^ 2) + (-1 // 2 * (getindex(X, 1) ^ 2 + getindex(X, 2) ^ 2)) * params.α
end)
Invariants:
GeometricBase.NullInvariants()
Timespan: (0.0, 10.0)
Timestep: 0.01
Initial conditions:
(t = 0.0, q = [1.0, 1.0], p = [0.5, 2.0], v = [0.0, 0.0])
Parameters:
(α = 5.0,)
This problem can again be integrated using GeometricIntegrators:
sol = integrate(lprob, Gauss(1))
fig = lines(parent(sol.q[:,1]), parent(sol.q[:,2]);
axis = (; xlabel = "x₁", ylabel = "x₂", title = "Particle moving in a square potential"),
figure = (; size = (800,600), fontsize = 22))
┌ Warning: Hermite Extrapolation: q's history[1] and history[2] are identical!
└ @ GeometricIntegrators.Extrapolators ~/.julia/packages/GeometricIntegrators/xUuFd/src/extrapolation/hermite.jl:211