Equations

GeometricIntegrators.Equations.DAEType

DAE: Differential Algebraic Equation

Defines a differential algebraic initial value problem

\[\begin{aligned} \dot{q} (t) &= v(t, q(t)) + u(t, q(t), \lambda(t)) , & q(t_{0}) &= q_{0} , \\ 0 &= \phi (t, q(t), \lambda(t)) , & \lambda(t_{0}) &= \lambda_{0} , \end{aligned}\]

with vector field $v$, projection $u$, algebraic constraint $\phi=0$, initial conditions $q_{0}$ and $\lambda_{0}$, the dynamical variable $q$ taking values in $\mathbb{R}^{m}$ and the algebraic variable $\lambda$ taking values in $\mathbb{R}^{n}$.

Fields

  • d: dimension of dynamical variable $q$ and the vector field $v$
  • m: dimension of algebraic variable $\lambda$ and the constraint $\phi$
  • n: number of initial conditions
  • v: function computing the vector field
  • u: function computing the projection
  • ϕ: algebraic constraint
  • : function computing an initial guess for the velocity field $v$ (optional)
  • h: function computing the Hamiltonian (optional)
  • t₀: initial time
  • q₀: initial condition for dynamical variable $q$
  • λ₀: initial condition for algebraic variable $\lambda$

The function v, providing the vector field, takes three arguments, v(t, q, v), the functions u and ϕ, providing the projection and the algebraic constraint take four arguments, u(t, q, λ, u) and ϕ(t, q, λ, ϕ), where t is the current time, q and λ are the current solution vectors, and v, u and ϕ are the vectors which hold the result of evaluating the vector field $v$, the projection $u$ and the algebraic constraint $\phi$ on t, q and λ.

Example

    function v(t, q, v)
        v[1] = q[1]
        v[2] = q[2]
    end

    function u(t, q, λ, u)
        u[1] = +λ[1]
        u[2] = -λ[1]
    end

    function ϕ(t, q, λ, ϕ)
        ϕ[1] = q[2] - q[1]
    end

    t₀ = 0.
    q₀ = [1., 1.]
    λ₀ = [0.]

    dae = DAE(v, u, ϕ, t₀, q₀, λ₀)
source
GeometricIntegrators.Equations.HDAEType

HDAE: Hamiltonian Differential Algebraic Equation EXPERIMENTAL

Defines a Hamiltonian differential algebraic initial value problem, that is a canonical Hamiltonian system of equations subject to Dirac constraints,

\[\begin{aligned} \dot{q} (t) &= v(t, q(t), p(t)) + u(t, q(t), p(t), \lambda(t)) + \bar{g}(t, q(t), p(t), \lambda(t), \gamma(t)) , & q(t_{0}) &= q_{0} , \\ \dot{p} (t) &= f(t, q(t), p(t)) + g(t, q(t), p(t), \lambda(t)) + \bar{f}(t, q(t), p(t), \lambda(t), \gamma(t)) , & p(t_{0}) &= p_{0} , \\ 0 &= \phi (t, q(t), p(t)) , \\ 0 &= \psi (t, q(t), p(t), \dot{q}(t), \dot{p}(t)) , \end{aligned}\]

with vector fields $v$, $u$, $\bar{u}$ and $f$, $g$, $\bar{g}$, primary constraint $\phi(q,p)=0$ and secondary constraint $\psi(q,p,\lambda)=0$, initial conditions $(q_{0}, p_{0})$, the dynamical variables $(q,p)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$ and the algebraic variables $(\lambda, \gamma)$ taking values in $\mathbb{R}^{n} \times \mathbb{R}^{d}$.

Fields

  • d: dimension of dynamical variables $q$ and $p$ as well as the vector fields $v$ and $f$
  • m: dimension of algebraic variables $\lambda$ and $\gamma$ and the constraints $\phi$ and $\psi$
  • v: function computing the Hamiltonian vector field $v$
  • f: function computing the Hamiltonian vector field $f$
  • u: function computing the primary projection field $u$
  • g: function computing the primary projection field $g$
  • : function computing the secondary projection field $\bar{u}$
  • : function computing the secondary projection field $\bar{g}$
  • ϕ: primary constraints
  • ψ: secondary constraints
  • : function computing an initial guess for the velocity field $v$` (optional)
  • : function computing an initial guess for the force field $f$ (optional)
  • h: function computing the Hamiltonian $H$
  • t₀: initial time
  • q₀: initial condition for dynamical variable $q$
  • p₀: initial condition for dynamical variable $p$
  • λ₀: initial condition for algebraic variable $λ$
source
GeometricIntegrators.Equations.HODEType

HODE: Hamiltonian Ordinary Differential Equation EXPERIMENTAL

Defines a Hamiltonian ordinary differential initial value problem, that is a canonical Hamiltonian system of equations,

\[\begin{aligned} \dot{q} (t) &= v(t, q(t), p(t)) , & q(t_{0}) &= q_{0} , \\ \dot{p} (t) &= f(t, q(t), p(t)) , & p(t_{0}) &= p_{0} , \end{aligned}\]

with vector fields $v$ and $f$, given by

\[\begin{aligned} v &= \frac{\partial H}{\partial p} , & f &= - \frac{\partial H}{\partial q} , \end{aligned}\]

initial conditions $(q_{0}, p_{0})$ and the dynamical variables $(q,p)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$.

Fields

  • d: dimension of dynamical variables $q$ and $p$ as well as the vector fields $v$ and $f$
  • v: function computing the vector field $v$
  • f: function computing the vector field $f$
  • h: function computing the Hamiltonian $H$
  • t₀: initial time
  • q₀: initial condition for dynamical variable $q$
  • p₀: initial condition for dynamical variable $p$
source
GeometricIntegrators.Equations.IDAEType

IDAE: Implicit Differential Algebraic Equation

Defines a partitioned differential algebraic initial value problem

\[\begin{aligned} \dot{q} (t) &= v(t) + u(t, q(t), p(t), \lambda(t)) , & q(t_{0}) &= q_{0} , \\ \dot{p} (t) &= f(t, q(t), v(t)) + r(t, q(t), p(t), \lambda(t)) , & p(t_{0}) &= p_{0} , \\ p(t) &= p(t, q(t), v(t)) , && \\ 0 &= \phi (t, q(t), p(t), \lambda(t)) , & \lambda(t_{0}) &= \lambda_{0} , \end{aligned}\]

with vector field $f$, the momentum defined by $p$, projection $u$ and $r$, algebraic constraint $\phi=0$, conditions $(q_{0}, p_{0})$ and $\lambda_{0}$, the dynamical variables $(q,p)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$ and the algebraic variable $\lambda$ taking values in $\mathbb{R}^{n}$.

Fields

  • d: dimension of dynamical variables $q$ and $p$ as well as the vector fields $f$ and $p$
  • m: dimension of algebraic variable $\lambda$ and the constraint $\phi$
  • ϑ: function determining the momentum
  • f: function computing the vector field $f$
  • u: function computing the projection for $q$
  • g: function computing the projection for $p$
  • ϕ: algebraic constraint
  • : function computing an initial guess for the velocity field $v$ (optional)
  • : function computing an initial guess for the force field $f$ (optional)
  • h: function computing the Hamiltonian (optional)
  • t₀: initial time
  • q₀: initial condition for dynamical variable $q$
  • p₀: initial condition for dynamical variable $p$
  • λ₀: initial condition for algebraic variable $\lambda$
source
GeometricIntegrators.Equations.IODEType

IODE: Implicit Ordinary Differential Equation

Defines an implicit initial value problem

\[\begin{aligned} \dot{q} (t) &= v(t) , & q(t_{0}) &= q_{0} , \\ \dot{p} (t) &= f(t, q(t), v(t)) , & p(t_{0}) &= p_{0} , \\ p(t) &= ϑ(t, q(t), v(t)) \end{aligned}\]

with vector field $f$, the momentum defined by $p$, initial conditions $(q_{0}, p_{0})$ and the solution $(q,p)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$. This is a special case of a differential algebraic equation with dynamical variables $(q,p)$ and algebraic variable $v$.

Fields

  • d: dimension of dynamical variables $q$ and $p$ as well as the vector fields $f$ and $p$
  • ϑ: function determining the momentum
  • f: function computing the vector field
  • g: function determining the projection, given by $\nabla \vartheta (q) \cdot \lambda$
  • : function computing an initial guess for the velocity field $v$ (optional)
  • : function computing an initial guess for the force field $f$ (optional)
  • h: function computing the Hamiltonian (optional)
  • t₀: initial time (optional)
  • q₀: initial condition for q
  • p₀: initial condition for p

The functions ϑ and f must have the interface

    function ϑ(t, q, v, p)
        p[1] = ...
        p[2] = ...
        ...
    end

and

    function f(t, q, v, f)
        f[1] = ...
        f[2] = ...
        ...
    end

where t is the current time, q is the current solution vector, v is the current velocity and f and p are the vectors which hold the result of evaluating the functions $f$ and $ϑ$ on t, q and v. In addition, the functions g, and are specified by

    function g(t, q, λ, g)
        g[1] = ...
        g[2] = ...
        ...
    end

    function v̄(t, q, v)
        v[1] = ...
        v[2] = ...
        ...
    end

    function f̄(t, q, v, f)
        f[1] = ...
        f[2] = ...
        ...
    end

The function g is used in projection methods that enforce $p = ϑ(q)$. The functions and are used for initial guesses in nonlinear implicit solvers.

source
GeometricIntegrators.Equations.ODEType

ODE: Ordinary Differential Equation

Defines an initial value problem

\[\dot{q} (t) = v(t, q(t)) , \qquad q(t_{0}) = q_{0} ,\]

with vector field $v$, initial condition $q_{0}$ and the solution $q$ taking values in $\mathbb{R}^{d}$.

Fields

  • d: dimension of dynamical variable $q$ and the vector field $v$
  • v: function computing the vector field
  • h: function computing the Hamiltonian (optional)
  • t₀: initial time
  • q₀: initial condition

The function v providing the vector field must have the interface

    function v(t, q, v)
        v[1] = ...
        v[2] = ...
        ...
    end

where t is the current time, q is the current solution vector, and v is the vector which holds the result of evaluating the vector field $v$ on t and q.

source
GeometricIntegrators.Equations.PDAEType

PDAE: Partitioned Differential Algebraic Equation

Defines a partitioned differential algebraic initial value problem

\[\begin{aligned} \dot{q} (t) &= v(t, q(t), p(t)) + u(t, q(t), p(t), \lambda(t)) , & q(t_{0}) &= q_{0} , \\ \dot{p} (t) &= f(t, q(t), p(t)) + r(t, q(t), p(t), \lambda(t)) , & p(t_{0}) &= p_{0} , \\ 0 &= \phi (t, q(t), p(t), \lambda(t)) , & \lambda(t_{0}) &= \lambda_{0} , \end{aligned}\]

with vector fields $v$ and $f$, projection $u$ and $r$, algebraic constraint $\phi=0$, conditions $(q_{0}, p_{0})$ and $\lambda_{0}$, the dynamical variables $(q,p)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$ and the algebraic variable $\lambda$ taking values in $\mathbb{R}^{n}$.

Fields

  • d: dimension of dynamical variables $q$ and $p$ as well as the vector fields $f$ and $p$
  • m: dimension of algebraic variable $\lambda$ and the constraint $\phi$
  • v: function computing the vector field $v$
  • f: function computing the vector field $f$
  • u: function computing the projection for $q$
  • g: function computing the projection for $p$
  • ϕ: algebraic constraint
  • : function computing an initial guess for the velocity field $v$ (optional)
  • : function computing an initial guess for the force field $f$ (optional)
  • h: function computing the Hamiltonian (optional)
  • t₀: initial time
  • q₀: initial condition for dynamical variable $q$
  • p₀: initial condition for dynamical variable $p$
  • λ₀: initial condition for algebraic variable $\lambda$
source
GeometricIntegrators.Equations.PODEType

PODE: Partitioned Ordinary Differential Equation

Defines a partitioned initial value problem

\[\begin{aligned} \dot{q} (t) &= v(t, q(t), p(t)) , & q(t_{0}) &= q_{0} , \\ \dot{p} (t) &= f(t, q(t), p(t)) , & p(t_{0}) &= p_{0} , \end{aligned}\]

with vector fields $v$ and $f$, initial conditions $(q_{0}, p_{0})$ and the solution $(q,p)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$.

Fields

  • d: dimension of dynamical variables $q$ and $p$ as well as the vector fields $v$ and $f$
  • v: function computing the vector field $v$
  • f: function computing the vector field $f$
  • h: function computing the Hamiltonian (optional)
  • t₀: initial time
  • q₀: initial condition for q
  • p₀: initial condition for p

The functions v and f must have the interface

    function v(t, q, p, v)
        v[1] = ...
        v[2] = ...
        ...
    end

and

    function f(t, q, p, f)
        f[1] = ...
        f[2] = ...
        ...
    end

where t is the current time, q and p are the current solution vectors and v and f are the vectors which hold the result of evaluating the vector fields $v$ and $f$ on t, q and p.

source
GeometricIntegrators.Equations.PSDEType

PSDE: Stratonovich Partitioned Stochastic Differential Equation

Defines a partitioned stochastic differential initial value problem

\[\begin{aligned} dq (t) &= v(t, q(t)) \, dt + B(t, q(t)) \circ dW , & q(t_{0}) &= q_{0} , \\ dp (t) &= f(t, q(t)) \, dt + G(t, q(t)) \circ dW , & p(t_{0}) &= p_{0} \end{aligned}\]

with the drift vector fields $v$ and $f$, diffusion matrices $B$ and $G$, initial conditions $q_{0}$ and $p_{0}$, the dynamical variables $(q,p)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$, and the m-dimensional Wiener process W

Fields

  • d: dimension of dynamical variable $q$ and the vector field $v$
  • m: dimension of the Wiener process
  • ns: number of sample paths
  • v: function computing the drift vector field for the position variable $q$
  • f: function computing the drift vector field for the momentum variable $p$
  • B: function computing the d x m diffusion matrix for the position variable $q$
  • G: function computing the d x m diffusion matrix for the momentum variable $p$
  • t₀: initial time
  • q₀: initial condition for dynamical variable $q$ (may be a random variable itself)
  • p₀: initial condition for dynamical variable $p$ (may be a random variable itself)

The functions v, f, B and G, providing the drift vector fields and diffusion matrices, take four arguments, v(t, q, p, v), f(t, q, p, f), B(t, q, p, B) and G(t, q, p, G), where t is the current time, (q, p) is the current solution vector, and v, f, B and G are the variables which hold the result of evaluating the vector fields $v$, $f$ and the matrices $B$, $G$ on t and (q,p).

Example

    function v(λ, t, q, v)
        v[1] = λ*q[1]
        v[2] = λ*q[2]
    end

    function B(μ, t, q, B)
        B[1] = μ*q[1]
        B[2] = μ*q[2]
    end

    t₀ = 0.
    q₀ = [1., 1.]
    λ  = 2.
    μ  = 1.

    v_sde = (t, q, v) -> v(λ, t, q, v)
    B_sde = (t, q, B) -> B(μ, t, q, B)

    sde = SDE(v_sde, B_sde, t₀, q₀)
source
GeometricIntegrators.Equations.SDEType

SDE: Stratonovich Stochastic Differential Equation

Defines a stochastic differential initial value problem

\[\begin{aligned} dq (t) &= v(t, q(t)) \, dt + B(t, q(t)) \circ dW , & q(t_{0}) &= q_{0} , \end{aligned}\]

with drift vector field $v$, diffusion matrix $B$, initial conditions $q_{0}$, the dynamical variable $q$ taking values in $\mathbb{R}^{d}$, and the m-dimensional Wiener process W

Fields

  • d: dimension of dynamical variable $q$ and the vector field $v$
  • m: dimension of the Wiener process
  • ns: number of sample paths
  • v: function computing the deterministic vector field
  • B: function computing the d x m diffusion matrix
  • t₀: initial time
  • q₀: initial condition for dynamical variable $q$ (may be a random variable itself)

Parameters

  • N: dimension of nitial condition array: N=1 - single, N=2 - multiple

The functions v and B, providing the drift vector field and diffusion matrix, v(t, q, v) and B(t, q, B, col=0), where t is the current time, q is the current solution vector, and v and B are the variables which hold the result of evaluating the vector field $v$ and the matrix $B$ on t and q (if col==0), or the column col of the matrix B (if col>0).

Example

    function v(t, q, v, p)
        λ = p[:λ]
        v[1] = λ*q[1]
        v[2] = λ*q[2]
    end

    function B(t, q, B, p, col=0)
        μ = p[:μ]
        if col==0 #whole matrix
            B[1,1] = μ*q[1]
            B[2,1] = μ*q[2]
        elseif col==1
            #just first column
        end
    end

    t₀ = 0.
    q₀ = [1., 1.]
    λ  = 2.
    μ  = 1.
    p = (λ=λ, μ=μ)

    sde = SDE(v, B, t₀, q₀; parameters=p)
source
GeometricIntegrators.Equations.SODEType

SODE: Split Ordinary Differential Equation

Defines an initial value problem

\[\dot{q} (t) = v(t, q(t)) , \qquad q(t_{0}) = q_{0} ,\]

with vector field $v$, initial condition $q_{0}$ and the solution $q$ taking values in $\mathbb{R}^{d}$. Here, the vector field $v$ is given as a sum of vector fields

\[v (t) = v_1 (t) + ... + v_r (t) .\]

Fields

  • d: dimension of dynamical variable $q$ and the vector field $v$
  • v: tuple of functions computing the vector field
  • q: tuple of functions computing the solution
  • t₀: initial time
  • q₀: initial condition

The functions v_i providing the vector field must have the interface

    function v_i(t, q, v)
        v[1] = ...
        v[2] = ...
        ...
    end

and the functions q_i providing the solutions must have the interface

    function q_i(t, q₀, q₁, h)
        q₁[1] = q₀[1] + ...
        q₁[2] = q₀[2] + ...
        ...
    end

where t is the current time, q₀ is the current solution vector, q₁ is the new solution vector which holds the result of computing one substep with the vector field $v_i$ on t and q₀, and h is the (sub-)timestep to compute the update for.

The fact that the function v returns the solution and not just the vector field for each substep increases the flexibility for the use of splitting methods, e.g., it allows to use another integrator for solving substeps.

source
GeometricIntegrators.Equations.SPDAEType

SPDAE: Split Partitioned Differential Algebraic Equation EXPERIMENTAL

Defines a split differential algebraic initial value problem, that is a canonical Hamiltonian system of equations subject to Dirac constraints,

\[\begin{aligned} \dot{q} (t) &= v_1(t, q(t), p(t)) + v_2(t, q(t), p(t), \lambda(t)) + v_3(t, q(t), p(t), \lambda(t), \gamma(t)) , & q(t_{0}) &= q_{0} , \\ \dot{p} (t) &= f_1(t, q(t), p(t)) + f_2(t, q(t), p(t), \lambda(t)) + f_3(t, q(t), p(t), \lambda(t), \gamma(t)) , & p(t_{0}) &= p_{0} , \\ 0 &= \phi (t, q(t), p(t)) , \\ 0 &= \psi (t, q(t), p(t), \dot{q}(t), \dot{p}(t)) , \end{aligned}\]

with vector fields $v_i$ and $f_i$ for $i = 1 ... 3$, primary constraint $\phi(q,p)=0$ and secondary constraint $\psi(q,p,\lambda)=0$, initial conditions $(q_{0}, p_{0})$, the dynamical variables $(q,p)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$ and the algebraic variables $(\lambda, \gamma)$ taking values in $\mathbb{R}^{n} \times \mathbb{R}^{d}$.

Fields

  • d: dimension of dynamical variables $q$ and $p$ as well as the vector fields $v$ and $f$
  • m: dimension of algebraic variables $\lambda$ and $\gamma$ and the constraints $\phi$ and $\psi$
  • n: number of initial conditions
  • v: tuple of functions computing the vector fields $v_i$, $i = 1 ... 3$
  • f: tuple of functions computing the vector fields $f_i$, $i = 1 ... 3$
  • ϕ: primary constraints
  • ψ: secondary constraints
  • t₀: initial time
  • q₀: initial condition for dynamical variable $q$
  • p₀: initial condition for dynamical variable $p$
  • λ₀: initial condition for algebraic variable $λ$
source
GeometricIntegrators.Equations.SPSDEType

SPSDE: Stratonovich Split Partitioned Stochastic Differential Equation

Defines a partitioned stochastic differential initial value problem

\[\begin{aligned} dq (t) &= v(t, q(t)) \, dt + B(t, q(t)) \circ dW , & q(t_{0}) &= q_{0} , \\ dp (t) &= [ f_1(t, q(t)) + f_2(t, q(t)) ] \, dt + [ G_1(t, q(t)) + G_2(t, q(t)) ] \circ dW , & p(t_{0}) &= p_{0} , \end{aligned}\]

with the drift vector fields $v$ and $f_i$, diffusion matrices $B$ and $G_i$, initial conditions $q_{0}$ and $p_{0}$, the dynamical variables $(q,p)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$, and the m-dimensional Wiener process W

Fields

  • d: dimension of dynamical variable $q$ and the vector fields $vi$
  • m: dimension of the Wiener process
  • ni: number of initial conditions
  • ns: number of sample paths
  • v : function computing the drift vector field for the position variable $q$
  • f1: function computing the drift vector field for the momentum variable $p$
  • f2: function computing the drift vector field for the momentum variable $p$
  • B : function computing the d x m diffusion matrix for the position variable $q$
  • G1: function computing the d x m diffusion matrix for the momentum variable $p$
  • G2: function computing the d x m diffusion matrix for the momentum variable $p$
  • t₀: initial time
  • q₀: initial condition for dynamical variable $q$ (may be a random variable itself)
  • p₀: initial condition for dynamical variable $p$ (may be a random variable itself)

The functions v, f, B and G, providing the drift vector fields and diffusion matrices, take four arguments, v(t, q, p, v), f(t, q, p, f), B(t, q, p, B) and G(t, q, p, G), where t is the current time, (q, p) is the current solution vector, and v, f, B and G are the variables which hold the result of evaluating the vector fields $v$, $f$ and the matrices $B$, $G$ on t and (q,p).

Example

    function v(λ, t, q, v)
        v[1] = λ*q[1]
        v[2] = λ*q[2]
    end

    function B(μ, t, q, B)
        B[1] = μ*q[1]
        B[2] = μ*q[2]
    end

    t₀ = 0.
    q₀ = [1., 1.]
    λ  = 2.
    μ  = 1.

    v_sde = (t, q, v) -> v(λ, t, q, v)
    B_sde = (t, q, B) -> B(μ, t, q, B)

    sde = SDE(v_sde, B_sde, t₀, q₀)
source
GeometricIntegrators.Equations.VDAEType

VDAE: Variational Differential Algebraic Equation EXPERIMENTAL

Defines an implicit initial value problem

\[\begin{aligned} \dot{q} (t) &= v(t) + \lambda(t), & q(t_{0}) &= q_{0} , \\ \dot{p} (t) &= f(t, q(t), v(t)) + g(t, q(t), \lambda(t)) + \bar{g} (t, q(t), \mu(t)) , & p(t_{0}) &= p_{0} , \\ p(t) &= ϑ(t, q(t), v(t)) , \\ 0 &= \phi (t, q(t), p(t)) , \\ 0 &= \psi (t, q(t), p(t), \dot{q}(t), \dot{p}(t)) , \end{aligned}\]

with vector field $f$, the momentum defined by $p$, initial conditions $(q_{0}, p_{0})$ and the solution $(q,p)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$ and the algebraic variables $(v, \lambda, \mu)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d} \times \mathbb{R}^{m}$. This is a special case of a differential algebraic equation with dynamical variables $(q,p)$ and algebraic variables $v$, $\lambda$ and $\mu$.

Fields

  • d: dimension of dynamical variables $q$ and $p$ as well as the vector fields $f$ and $p$
  • ϑ: function determining the momentum
  • f: function computing the vector field
  • g: function determining the primary projection, usually given by $\nabla \vartheta (q) \cdot \lambda$
  • : function determining the secondary projection, usually given by $\lambda \cdot \nabla \vartheta (q)$
  • ϕ: primary constraints, usually given by $p - \vartheta (q)$
  • ψ: secondary constraints, usually given by $\dot{p} - \dot{q} \cdot \nabla \vartheta (q)$
  • : function computing an initial guess for the velocity field $v$ (optional)
  • : function computing an initial guess for the force field $f$ (optional)
  • h: function computing the Hamiltonian (optional)
  • Ω: symplectic matrix (optional)
  • ∇H: gradient of the Hamiltonian (optional)
  • t₀: initial time (optional)
  • q₀: initial condition for q
  • p₀: initial condition for p
  • λ₀: initial condition for λ (optional)
  • μ₀: initial condition for μ (optional)

The functions ϑ and f must have the interface

    function ϑ(t, q, v, p)
        p[1] = ...
        p[2] = ...
        ...
    end

and

    function f(t, q, v, f)
        f[1] = ...
        f[2] = ...
        ...
    end

where t is the current time, q is the current solution vector, v is the current velocity and f and p are the vectors which hold the result of evaluating the functions $f$ and $ϑ$ on t, q and v. The funtions g, and are specified by

    function g(t, q, λ, g)
        g[1] = ...
        g[2] = ...
        ...
    end

    function v̄(t, q, v)
        v[1] = ...
        v[2] = ...
        ...
    end

    function f̄(t, q, v, f)
        f[1] = ...
        f[2] = ...
        ...
    end
source
GeometricIntegrators.Equations.VODEType

VODE: Variational Ordinary Differential Equation EXPERIMENTAL

Defines an implicit initial value problem

\[\begin{aligned} \dot{q} (t) &= v(t) , & q(t_{0}) &= q_{0} , \\ \dot{p} (t) &= f(t, q(t), v(t)) , & p(t_{0}) &= p_{0} , \\ p(t) &= ϑ(t, q(t), v(t)) \end{aligned}\]

with vector field $f$, the momentum defined by $p$, initial conditions $(q_{0}, p_{0})$ and the solution $(q,p)$ taking values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$. This is a special case of a differential algebraic equation with dynamical variables $(q,p)$ and algebraic variable $v$.

Fields

  • d: dimension of dynamical variables $q$ and $p$ as well as the vector fields $f$ and $p$
  • ϑ: function determining the momentum
  • f: function computing the vector field
  • g: function determining the projection, given by $\nabla \vartheta (q) \cdot \lambda$
  • : function computing an initial guess for the velocity field $v$ (optional)
  • : function computing an initial guess for the force field $f$ (optional)
  • h: function computing the Hamiltonian (optional)
  • Ω: symplectic matrix (optional)
  • ∇H: gradient of the Hamiltonian (optional)
  • t₀: initial time (optional)
  • q₀: initial condition for q
  • p₀: initial condition for p
  • λ₀: initial condition for λ (optional)

The functions ϑ and f must have the interface

    function ϑ(t, q, v, p)
        p[1] = ...
        p[2] = ...
        ...
    end

and

    function f(t, q, v, f)
        f[1] = ...
        f[2] = ...
        ...
    end

where t is the current time, q is the current solution vector, v is the current velocity and f and p are the vectors which hold the result of evaluating the functions $f$ and $ϑ$ on t, q and v. The funtions g and v are specified by

    function g(t, q, λ, g)
        g[1] = ...
        g[2] = ...
        ...
    end

and

    function v(t, q, p, v)
        v[1] = ...
        v[2] = ...
        ...
    end
source