Problem Types

The following data structures are all implemented in GeometricEquations.jl. Each problem type is derived from EquationProblem.

Geometric Equation Problems

GeometricEquations.EquationProblemType

EquationProblem: stores a GeometricEquation together with initial conditions, parameters, time span and time step size.

Parameters

  • ST <: GeometricEquation: super type, used for dispatch
  • DT <: Number: data type
  • TT <: Real: time step type
  • AT <: AbstractArray{DT}: array type of state variable
  • equType <: GeometricEquation: equation type
  • functionsType <: NamedTuple: types of all function methods
  • solutionsType <: NamedTuple: types of all solution methods
  • icsType <: NamedTuple: types of all initial conditions
  • parType <: OptionalParameters: parameters type

Fields

  • equation: reference to the parent equation object holding the vector fields, etc.
  • functions: methods for all vector fields, etc., that define the problem
  • solutions: methods for all solutions, etc., if defined
  • tspan: time span for problem (t₀,t₁)
  • tstep: time step to be used in simulation
  • ics: NamedTuple containing the initial conditions, must contain one field for each state variable
  • parameters: either a NamedTuple containing the equation's parameters or NullParameters indicating that the equation does not have any parameters

Subtypes

The EquationProblem type has various subtypes for the different equations types, that are defined e.g. via

const ODEProblem = EquationProblem{ODE}

and provide convenience constructors to construct an equation and the corresponding problem in one step, e.g.,

ODEProblem(v, tspan, tstep, ics::NamedTuple; kwargs...)
ODEProblem(v, tspan, tstep, q₀::StateVariable; kwargs...)

All problem subtypes take the following keyword arguments:

  • invariants = NullInvariants()
  • parameters = NullParameters()
  • periodicity = NullPeriodicity()

If not set to their corresponding Null types, the user needs to pass a NamedTuple whose values are

  • functions for invariants,
  • arbitrary data structures for parameters,
  • the same data structure as the solution for periodicity.

The latter should be zero everywhere, except for those components, that are periodic, i.e., whose value are supposed to stay within a range (0, max). Support for ranges starting with other values than zero is currently missing but can be added if demand arises.

source
GeometricEquations.EnsembleProblemType

EnsembleProblem: stores a GeometricEquation together with multiple sets of initial conditions and/or parameters, a time span for integration and a time step size.

An EnsembleProblem is initialized by providing a GeometricEquation, an integration time span (typically a tuple with two values, start and end time, respectively), a timestep, and one of the three following options:

  • a vector of initial conditions and a vector of parameter sets, with both vectors having the same length,
  • a vector of initial conditions and a single set of parameters,
  • a single initial condition and a vector of parameter sets.

Each initial condition is a NamedTuple that contains one field for each state variable. Each parameter set is either a NamedTuple containing the equation's parameters or NullParameters indicating that the equation does not have any parameters.

The different constructors then generate two vectors, one for the initial conditions and one for the parameters, where each pair of the corresponding entries defines one problem. In the first case, the respective constructor just checks if both vectors are of the same size. In the second case, the respective constructor creates a parameter vector, where each entry holds the same parameter set. In the third case, the respective constructor creates an initial condition vector, where each entry holds the same initial conditions. One may be inclined to think that the first and second constructor lead to a waste of memory, but in reality the respective vectors only hold references to the same initial conditons or parameters. Thus the data is not actually duplicated.

Each pair of initial conditions and parameters is referred to as sample. The methods

length(::EnsembleProblem)
nsamples(::EnsembleProblem)

return the number of samples in a EnsembleProblem. A single initial condition or parameter set can be retrieved by the methods

initial_condition(::EnsembleProblem, i)
parameter(::EnsembleProblem, i)

where i is the index of the sample. Typically, however, e.g. when integrating all samples of an EnsembleProblem with GeometricIntegrators, it is more convenient to retrieve the corresponding GeometricProblem via the method

problem(::EnsembleProblem, i)

The EnsembleProblem also allows to iterate over all samples, e.g.

for problem in ensemble
    # ...
    # integrate problem
    # ...
end

where ensemble is an EnsembleProblem and problem is the corresponding GeometricProblem.

Parameters

  • ST <: GeometricEquation: super type, used for dispatch
  • DT <: Number: data type
  • TT <: Real: time step type
  • AT <: AbstractArray{DT}: array type of state variable
  • equType <: GeometricEquation: equation type
  • functionsType <: NamedTuple: types of all function methods
  • solutionsType <: NamedTuple: types of all solution methods
  • icsType <: AbstractVector{<:NamedTuple}: types of all initial conditions
  • parType <: AbstractVector{<:OptionalParameters}: parameters type

Fields

  • equation: reference to the parent equation object holding the vector fields, etc.
  • functions: methods for all vector fields, etc., that define the problem
  • solutions: methods for all solutions, etc., if defined
  • tspan: time span for problem (t₀,t₁)
  • tstep: time step to be used in simulation
  • ics: vector of NamedTuple containing the initial conditions, each NamedTuple must contain one field for each state variable
  • parameters: vector of either NamedTuple containing the equation's parameters or NullParameters indicating that the equation does not have any parameters

Constructors

The EnsembleProblem provides the following constructors:

EnsembleProblem(equ, tspan, tstep, ics::AbstractVector{<:NamedTuple}, parameters::AbstractVector{<:OptionalParameters})
EnsembleProblem(equ, tspan, tstep, ics::AbstractVector{<:NamedTuple}, parameters::OptionalParameters=NullParameters())
EnsembleProblem(equ, tspan, tstep, ics::NamedTuple, parameters::AbstractVector{<:OptionalParameters})
EnsembleProblem(equ, tspan, tstep, ics, ::Nothing) = 
    EnsembleProblem(equ, tspan, tstep, ics, NullParameters())
EnsembleProblem(equ, tspan, tstep, ics; parameters = NullParameters()) = 
    EnsembleProblem(equ, tspan, tstep, ics, parameters)
  • equ is a subtype of GeometricEquation
  • tspan is a tuple (t₀,t₁) of the integration time span with t₀ the start time and t₁ the end time
  • tstep is the time step, typically a value of some AbstractFloat subtype
  • ics are the initial conditions, either a single set or a vector of multiple sets
  • parameters are the static parameters of the problem, either a single set or a vector of multiple sets
source

Ordinary Differential Equations

GeometricEquations.ODEProblemType

ODEProblem: Ordinary Differential Equation Problem

Ordinary differential equations define an initial value problem of the form

\[\dot{q} (t) = v(t, q(t)) ,\]

with vector field $v$.

The dynamical variables $q$ with initial condition $q_{0}$ take values in $\mathbb{R}^{d}$.

Constructors

ODEProblem(v, tspan, tstep, ics::NamedTuple; kwargs...)
ODEProblem(v, tspan, tstep, q₀::StateVariable; kwargs...)
ODEProblem(v, tspan, tstep, q₀::AbstractArray; kwargs...)

where v is the function computing the vector field, tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entry q of type StateVariable. The initial condition q₀ can also be prescribed directly, as a StateVariable or an AbstractArray{<:Number}.

For possible keyword arguments see the documentation on EquationProblem subtypes.

Function Definitions

The function v providing the vector field must have the interface

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

where t is the current time, q is the current solution vector, v is the vector which holds the result of evaluating the vector field $v$ on t and q, and params is a NamedTuple of additional parameters on which the vector field may depend.

source
GeometricEquations.PODEProblemType

PODEProblem: Partitioned Ordinary Differential Equation Problem

A partitioned ordinary differential equation is an initial value problem of the form

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

with vector fields $v$ and $f$.

The dynamical variables $(q,p)$ with initial conditions $(q(t_{0}) = q_{0}, p(t_{0}) = p_{0})$ take values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$.

Constructors

PODEProblem(v, f, tspan, tstep, ics; kwargs...)
PODEProblem(v, f, tspan, tstep, q₀::StateVariable, p₀::StateVariable; kwargs...)
PODEProblem(v, f, tspan, tstep, q₀::AbstractArray, p₀::AbstractArray; kwargs...)

where v and f are the function computing the vector fields, tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entries q and p. The initial conditions q₀ and p₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}. For the interfaces of the functions v and f see PODE.

For possible keyword arguments see the documentation on EquationProblem subtypes.

Function Definitions

The functions v and f must have the interface

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

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

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

source
GeometricEquations.HODEProblemType

HODEProblem: Hamiltonian Ordinary Differential Equation Problem

A canonical Hamiltonian system of equations is special case of a partitioned ordinary differential equation,

\[\begin{aligned} \dot{q} (t) &= v(t, q(t), p(t)) , \\ \dot{p} (t) &= f(t, q(t), p(t)) , \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}\]

The dynamical variables $(q,p)$ with initial conditions $(q(t_{0}) = q_{0}, p(t_{0}) = p_{0})$ take values in $T^{*} Q \simeq \mathbb{R}^{d} \times \mathbb{R}^{d}$.

Constructors

HODEProblem(v, f, hamiltonian, tspan, tstep, ics; kwargs...)
HODEProblem(v, f, hamiltonian, tspan, tstep, q₀::StateVariable, p₀::StateVariable; kwargs...)
HODEProblem(v, f, hamiltonian, tspan, tstep, q₀::AbstractArray, p₀::AbstractArray; kwargs...)

where v and f are the function computing the vector fields, hamiltonian returns the value of the Hamiltonian (i.e. the total energy), tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entries q and p. The initial conditions q₀ and p₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}. For the interfaces of the functions v, f and hamiltonian see HODE.

For possible keyword arguments see the documentation on EquationProblem subtypes.

Function Definitions

The functions v, f and hamiltonian must have the interface

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

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

function hamiltonian(t, q, p, params)
    return ...
end

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

source
GeometricEquations.IODEProblemType

IODEProblem: Implicit Ordinary Differential Equation Problem

An implicit ordinary differential equations is an initial value problem of the form

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

with force field $f$, the momentum defined by $p$. This is a special case of a differential algebraic equation with dynamical variables $(q,p)$ and algebraic variable $v$, that is determined such that the constraint $p(t) = ϑ(t, q(t), v(t))$ is satisfied.

Most integrators perform a projection step in order to enforce this constraint. To this end, the system is extended to

\[\begin{aligned} \dot{q} (t) &= v(t) + λ(t) , \\ \dot{p} (t) &= f(t, q(t), v(t)) + g(t, q(t), v(t), λ(t)) , \\ p(t) &= ϑ(t, q(t), v(t)) , \end{aligned}\]

where the vector field defining the projection step is usually given as

\[\begin{aligned} g(t, q(t), v(t), λ(t)) &= λ(t) \cdot \nabla ϑ(t, q(t), v(t)) . \end{aligned}\]

The dynamical variables $(q,p)$ with initial conditions $(q(t_{0}) = q_{0}, p(t_{0}) = p_{0})$ take values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$. The algebraic variable $λ$ with initial condition $λ(t_{0}) = λ_{0}$ takes values in $\mathbb{R}^{m}$.

Constructors

IODEProblem(ϑ, f, tspan, tstep, ics; kwargs...)
IODEProblem(ϑ, f, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::AlgebraicVariable; kwargs...)
IODEProblem(ϑ, f, tspan, tstep, q₀::AbstractArray, p₀::AbstractArray, λ₀::AbstractArray = zero(q₀); kwargs...)
IODEProblem(ϑ, f, g, tspan, tstep, ics; kwargs...)
IODEProblem(ϑ, f, g, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::AlgebraicVariable; kwargs...)
IODEProblem(ϑ, f, g, tspan, tstep, q₀::AbstractArray, p₀::AbstractArray, λ₀::AbstractArray = zero(q₀); kwargs...)

The functions ϑ, f and g compute the momentum and the vector fields, respectively.

tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entries q, p and λ. The initial conditions q₀ and p₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}. For the interfaces of the functions ϑ, f and g see IODE.

In addition to the standard keyword arguments for EquationProblem subtypes, an IODEProblem accepts functions and for the computation of initial guesses for the vector fields with default values v̄ = _iode_default_v̄ and f̄ = f.

source
GeometricEquations.LODEProblemType

LODEProblem: Lagrangian Ordinary Differential Equation Problem

A Lagrangian system of equations is a special case of an implicit ordinary differential equations, that is an implicit initial value problem of the form

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

with momentum $p$ and force field $f$, given by

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

This is a special case of an implicit ordinary differential equation, that is defined by a Lagrangian, as well as a special case of a differential algebraic equation with dynamical variables $(q,p)$ and algebraic variable $v$, that is determined such that the constraint $p(t) = ϑ(t, q(t), v(t))$ is satisfied.

Many integrators perform a projection step in order to enforce this constraint. To this end, the system is extended to

\[\begin{aligned} \dot{q} (t) &= v(t) + \lambda(t) , \\ \dot{p} (t) &= f(t, q(t), v(t)) + g(t, q(t), v(t), \lambda(t)) , \\ p(t) &= ϑ(t, q(t), v(t)) , \end{aligned}\]

where the vector field defining the projection step is usually given as

\[\begin{aligned} g(t, q(t), v(t), λ(t)) &= λ(t) \cdot \nabla ϑ(t, q(t), v(t)) . \end{aligned}\]

The dynamical variables $(q,p)$ with initial conditions $(q(t_{0}) = q_{0}, p(t_{0}) = p_{0})$ take values in $T^{*} Q \simeq \mathbb{R}^{d} \times \mathbb{R}^{d}$. The algebraic variable $λ$ with initial condition $λ(t_{0}) = λ_{0}$ takes values in $\mathbb{R}^{m}$.

Constructors

LODEProblem(ϑ, f, ω, l, tspan, tstep, ics; kwargs...)
LODEProblem(ϑ, f, ω, l, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::AlgebraicVariable; kwargs...)
LODEProblem(ϑ, f, ω, l, tspan, tstep, q₀::AbstractArray, p₀::AbstractArray, λ₀::AbstractArray = zero(q₀); kwargs...)
LODEProblem(ϑ, f, g, ω, l, tspan, tstep, ics; kwargs...)
LODEProblem(ϑ, f, g, ω, l, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::AlgebraicVariable; kwargs...)
LODEProblem(ϑ, f, g, ω, l, tspan, tstep, q₀::AbstractArray, p₀::AbstractArray, λ₀::AbstractArray = zero(q₀); kwargs...)

where ϑ, f and g are the functions computing the momentum and the vector fields, respectively, ω determines the symplectic matrix, and l returns the Lagrangian, tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entries q, p and λ. The initial conditions q₀, p₀ and λ₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}, where λ₀ can also be omitted. For the interfaces of the functions ϑ, f, g, ω and l see LODE.

In addition to the standard keyword arguments for EquationProblem subtypes, a LODEProblem accepts functions and for the computation of initial guesses for the vector fields with default values v̄ = _lode_default_v̄ and f̄ = f.

source
GeometricEquations.SODEProblemType

SODEProblem: Split Ordinary Differential Equation Problem

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) .\]

The dynamical variables $q$ with initial condition $q_{0}$ take values in $\mathbb{R}^{d}$.

Constructors

SODEProblem(v, q, tspan, tstep, ics::NamedTuple; kwargs...)
SODEProblem(v, q, tspan, tstep, q₀::StateVariable; kwargs...)
SODEProblem(v, q, tspan, tstep, q₀::AbstractArray; kwargs...)
SODEProblem(v, tspan, tstep, ics::NamedTuple; kwargs...)
SODEProblem(v, tspan, tstep, q₀::StateVariable; kwargs...)
SODEProblem(v, tspan, tstep, q₀::AbstractArray; kwargs...)

where v is a tuple of functions computing the vector fields for each substep, q is an optional tuple of functions computing the solution for each substep, tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entry q. The initial condition q₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}.

For possible keyword arguments see the documentation on EquationProblem subtypes.

Function Definitions

The functions v_i providing the vector field must have the interface

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

and the functions q_i providing the solutions must have the interface

function q_i(q₁, t₁, q₀, t₀, params)
    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 at time t₁, holding the result of computing one substep

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. with the vector field $v_i$.

source

Differential Algebraic Equations

GeometricEquations.DAEProblemType

DAEProblem: Differential Algebraic Equation Problem

Defines a differential algebraic initial value problem

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

with vector field $v$, projection $u$, algebraic constraint $\phi=0$.

Some integrators also enforce the secondary constraint $\psi$, that is the time derivative of the algebraic constraint $\phi$. In this case, the system of equations is modified as follows

\[\begin{aligned} \dot{q} (t) &= v(t, q(t)) + u(t, q(t), \lambda(t)) + \bar{u} (t, q(t), \dot{q} (t), \dot{p} (t), \mu(t)) , \\ 0 &= \phi (t, q(t)) , \\ 0 &= \psi (t, q(t), \dot{q} (t)) . \end{aligned}\]

The dynamical variable $q$ with initial conditions $q(t_{0}) = q_{0}$ takes values in $\mathbb{R}^{d}$. The algebraic variables $(λ,μ)$ with initial condition $(λ(t_{0}) = λ_{0}, μ(t_{0}) = μ_{0})$ take values in $\mathbb{R}^{m} \times \mathbb{R}^{m}$.

Constructors

DAEProblem(v, u, ϕ, ū, ψ, tspan, tstep, ics::NamedTuple; kwargs...)
DAEProblem(v, u, ϕ, ū, ψ, tspan, tstep, q₀::StateVariable, λ₀::StateVariable, μ₀::StateVariable = zero(λ₀); kwargs...)
DAEProblem(v, u, ϕ, tspan, tstep, ics::NamedTuple; kwargs...)
DAEProblem(v, u, ϕ, tspan, tstep, q₀::StateVariable, λ₀::StateVariable; kwargs...)

The functions v and u compute the vector field and the projection, respectively, ϕ provides the algebraic constraint. The functions ψ and are optional and provide the secondary constraint, that is the time derivative of the algebraic constraint, and the corresponding projection.

tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entries q, λ and μ. The initial conditions q₀, λ₀ and μ₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}. For the interfaces of the functions v, u, ϕ, , ψ see DAE.

In addition to the standard keyword arguments for EquationProblem subtypes, a DAEProblem accepts a function for the computation of an initial guess for the vector field with default value v̄ = v.

Function Definitions

The functions v, u and ϕ must have the interface

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

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

function ϕ(ϕ, t, q, params)
    ϕ[1] = ...
end

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 λ.

Some integrators also enforce the secondary constraint $\psi$ and require the following additional functions

function ū(u, t, q, μ, params)
    u[1] = ...
    u[2] = ...
    ...
end

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

With the above function definitions the DAEProblem can be created by

tspan = (0.0, 1.0)
tstep = 0.1
q₀ = [1., 1.]
λ₀ = [0.]
μ₀ = [0.]

prob = DAEProblem(v, u, ϕ, tspan, tstep, q₀, λ₀)

or

prob = DAEProblem(v, u, ϕ, ū, ψ, tspan, tstep, q₀, λ₀, μ₀)
source
GeometricEquations.PDAEProblemType

PDAEProblem: Partitioned Differential Algebraic Equation Problem

A partitioned differential algebraic equation has the form

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

with vector fields $v$ and $f$, projection $u$ and $g$, algebraic constraint $\phi=0$.

Some integrators also enforce the secondary constraint $\psi$, that is the time derivative of the algebraic constraint $\phi$. In this case, the system of equations is modified as follows

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

The dynamical variables $(q,p)$ with initial conditions $(q(t_{0}) = q_{0}, p(t_{0}) = p_{0})$ take values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$. The algebraic variables $(λ,μ)$ with initial condition $(λ(t_{0}) = λ_{0}, μ(t_{0}) = μ_{0})$ take values in $\mathbb{R}^{m} \times \mathbb{R}^{m}$.

Constructors

PDAEProblem(v, f, u, g, ϕ, ū, ḡ, ψ, tspan, tstep, ics::NamedTuple; kwargs...)
PDAEProblem(v, f, u, g, ϕ, ū, ḡ, ψ, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::StateVariable, μ₀::StateVariable = zero(λ₀); kwargs...)
PDAEProblem(v, f, u, g, ϕ, tspan, tstep, ics::NamedTuple; kwargs...)
PDAEProblem(v, f, u, g, ϕ, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::StateVariable; kwargs...)

The functions v and f compute the vector field, u and g compute the projections, and ϕ provides the algebraic constraint. The functions ψ, and are optional and provide the secondary constraint and the corresponding projection.

tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entries q, p, λ and μ. The initial conditions q₀, p₀, λ₀ and μ₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}. For the interfaces of the functions v, f, u, g, ϕ, , , ψ see PDAE.

In addition to the standard keyword arguments for EquationProblem subtypes, a PDAEProblem accepts functions and for the computation of initial guesses for the vector fields with default values v̄ = v and f̄ = f.

Function Definitions

The functions v, f, u, g and ϕ must have the interface

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

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

function u(u, t, q, p, λ, params)
    u[1] = ...
    u[2] = ...
    ...
end

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

function ϕ(ϕ, t, q, p, params)
    ϕ[1] = ...
end

where t is the current time, q, p and λ are the current solution vectors, v, f, u and g are the vectors which hold the result of evaluating the vector fields $v$ and $f$, the projections $u$ and $g$, and ϕ holds the algebraic constraint $\phi$, evaluated on t, q, p and λ.

Some integrators also enforce the secondary constraint $\psi$ and require the following additional functions

function ū(u, t, q, p, μ, params)
    u[1] = ...
    u[2] = ...
    ...
end

function ḡ(g, t, q, p, μ, params)
    g[1] = ...
    g[2] = ...
    ...
end

function ψ(ψ, t, q, p, v, f, params)
    ψ[1] = ...
end

With the above function definitions the PDAEProblem can be created by

tspan = (0.0, 1.0)
tstep = 0.1
q₀ = [1., 1.]
p₀ = [1., 0.]
λ₀ = [0.]
μ₀ = [0.]

prob = PDAEProblem(v, f, u, g, ϕ, tspan, tstep, q₀, p₀, λ₀)

or

prob = PDAEProblem(v, f, u, g, ϕ, ū, ḡ, ψ, tspan, tstep, q₀, p₀, λ₀, μ₀)
source
GeometricEquations.HDAEProblemType

HDAEProblem: Hamiltonian Differential Algebraic Equation

A Hamiltonian differential algebraic is an 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{u} (t, q(t), p(t), \mu(t)) , \\ \dot{p} (t) &= f(t, q(t), p(t)) + g(t, q(t), p(t), \lambda(t)) + \bar{g} (t, q(t), p(t), \mu(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 fields $v$, $u$, $\bar{u}$ and $f$, $g$, $\bar{g}$, primary constraint $\phi(q,p)=0$ and secondary constraint $\psi(q,p,\dot{q},\dot{p})=0$.

The dynamical variables $(q,p)$ with initial conditions $(q(t_{0}) = q_{0}, p(t_{0}) = p_{0})$ take values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$. The algebraic variables $(λ,μ)$ with initial condition $(λ(t_{0}) = λ_{0}, μ(t_{0}) = μ_{0})$ take values in $\mathbb{R}^{m} \times \mathbb{R}^{m}$.

Constructors

HDAEProblem(v, f, u, g, ϕ, ū, ḡ, ψ, h, tspan, tstep, ics::NamedTuple; kwargs...)
HDAEProblem(v, f, u, g, ϕ, ū, ḡ, ψ, h, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::StateVariable, μ₀::StateVariable = zero(λ₀); kwargs...)
HDAEProblem(v, f, u, g, ϕ, h, tspan, tstep, ics::NamedTuple; kwargs...)
HDAEProblem(v, f, u, g, ϕ, h, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::StateVariable; kwargs...)

The functions v and f compute the vector field, u and g compute the projections, ϕ provides the algebraic constraint and h the Hamiltonian. The functions ψ, and are optional and provide the secondary constraint, that is the time derivative of the algebraic constraint, and the corresponding projection.

tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entries q, p, λ and μ. The initial conditions q₀, p₀, λ₀ and μ₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}. For the interfaces of the functions v, f, u, g, ϕ, , , ψ, and h see HDAE.

In addition to the standard keyword arguments for EquationProblem subtypes, a HDAEProblem accepts functions and for the computation of initial guesses for the vector fields with default values v̄ = v and f̄ = f.

Function Definitions

The functions v, f, u, g, ϕ and h must have the interface

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

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

function u(u, t, q, p, λ, params)
    u[1] = ...
    u[2] = ...
    ...
end

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

function ϕ(ϕ, t, q, p, params)
    ϕ[1] = ...
end

function h(t, q, p, params)
    ...
end

where t is the current time, q, p, λ and μ are the current solution vectors, v, f, u and g are the vectors which hold the result of evaluating the vector fields $v$ and $f$, the projections on the primary constraint $u$ and $g$, ϕ holds the algebraic constraint $\phi$, and h returns the Hamiltonian of the system, all evaluated on t, q, p and λ.

Some integrators also enforce the secondary constraint $\psi$ and require the following additional functions

function ū(u, t, q, p, μ, params)
    u[1] = ...
    u[2] = ...
    ...
end

function ḡ(g, t, q, p, μ, params)
    g[1] = ...
    g[2] = ...
    ...
end

function ψ(ψ, t, q, p, v, f, params)
    ψ[1] = ...
end

With the above function definitions the HDAEProblem can be created by

tspan = (0.0, 1.0)
tstep = 0.1
q₀ = [1., 1.]
p₀ = [1., 0.]
λ₀ = [0.]
μ₀ = [0.]

prob = HDAEProblem(v, f, u, g, ϕ, h, tspan, tstep, q₀, p₀, λ₀)

or

prob = HDAEProblem(v, f, u, g, ϕ, ū, ḡ, ψ, h, tspan, tstep, q₀, p₀, λ₀, μ₀)
source
GeometricEquations.IDAEProblemType

IDAEProblem: Implicit Differential Algebraic Equation Problem

An implicit differential algebraic initial value problem takes the form

\[\begin{aligned} \dot{q} (t) &= v(t) + u(t, q(t), v(t), p(t), \lambda(t)) , \\ \dot{p} (t) &= f(t, q(t), v(t)) + g(t, q(t), v(t), p(t), \lambda(t)) , \\ p(t) &= ϑ(t, q(t), v(t)) , \\ 0 &= \phi (t, q(t), v(t), p(t)) , \end{aligned}\]

with force field $f$, the momentum defined by $ϑ$, projections $u$ and $g$, algebraic constraint $\phi(t,q,v,p)=0$.

Some integrators also enforce the secondary constraint $\psi$, that is the time derivative of the algebraic constraint $\phi$. In this case, the system of equations is modified as follows

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

The dynamical variables $(q,p)$ with initial conditions $(q(t_{0}) = q_{0}, p(t_{0}) = p_{0})$ take values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$. The algebraic variables $(λ,μ)$ with initial condition $(λ(t_{0}) = λ_{0}, μ(t_{0}) = μ_{0})$ take values in $\mathbb{R}^{m} \times \mathbb{R}^{m}$.

Constructors

IDAEProblem(ϑ, f, u, g, ϕ, ū, ḡ, ψ, tspan, tstep, ics; kwargs...)
IDAEProblem(ϑ, f, u, g, ϕ, ū, ḡ, ψ, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::StateVariable = zero(q₀), μ₀::StateVariable = zero(λ₀); kwargs...)
IDAEProblem(ϑ, f, u, g, ϕ, tspan, tstep, ics; kwargs...)
IDAEProblem(ϑ, f, u, g, ϕ, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::StateVariable = zero(q₀); kwargs...)

The function ϑ computes the momentum, f computes the force field, u and g compute the projections, and ϕ provides the algebraic constraint. The functions ψ, and are optional and provide the secondary constraint, that is the time derivative of the algebraic constraint, and the corresponding projection.

tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entries q and p. The initial conditions q₀, p₀, λ₀ and μ₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}. For the interfaces of the functions ϑ, f, u, g, ϕ, , , ψ see IDAE.

In addition to the standard keyword arguments for EquationProblem subtypes, an IDAEProblem accepts functions and for the computation of initial guesses for the vector fields with default values v̄ = _idae_default_v̄ and f̄ = f.

Function Definitions

The functions ϑ and f must have the interface

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

function f(f, t, q, v, params)
    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 u(u, t, q, v, p, λ, params)
    u[1] = ...
    u[2] = ...
    ...
end

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

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

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

Some integrators also enforce the secondary constraint $\psi$ and require the following additional functions

function ū(u, t, q, v, p, μ, params)
    u[1] = ...
    u[2] = ...
    ...
end

function ḡ(g, t, q, v, p, μ, params)
    g[1] = ...
    g[2] = ...
    ...
end

function ψ(ψ, t, q, v, p, q̇, ṗ, params)
    ψ[1] = ...
end
source
GeometricEquations.LDAEProblemType

LDAEProblem: Lagrangian Differential Algebraic Equation Problem

A special case of an implicit initial value problem is a Lagrangian differential algebraic equation of the form

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

with momentum $p$ and force field $f$, given by

\[\begin{aligned} p &= \frac{\partial L}{\partial v} (q,v) , & f &= \frac{\partial L}{\partial q} (q,v) , \end{aligned}\]

projection fields $u$, $\bar{u}$ and $g$, $\bar{g}$. This is a special case of a differential algebraic equation with dynamical variables $(q,p)$ and algebraic variables $v$, $\lambda$ and $\mu$.

The dynamical variables $(q,p)$ with initial conditions $(q(t_{0}) = q_{0}, p(t_{0}) = p_{0})$ take values in $\mathbb{R}^{d} \times \mathbb{R}^{d}$. The algebraic variables $(λ,μ)$ with initial condition $(λ(t_{0}) = λ_{0}, μ(t_{0}) = μ_{0})$ take values in $\mathbb{R}^{m} \times \mathbb{R}^{m}$.

Constructors

LDAEProblem(ϑ, f, u, g, ϕ, ū, ḡ, ψ, ω, l, tspan, tstep, ics; kwargs...)
LDAEProblem(ϑ, f, u, g, ϕ, ū, ḡ, ψ, ω, l, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::StateVariable = zero(q₀), μ₀::StateVariable = zero(λ₀); kwargs...)
LDAEProblem(ϑ, f, u, g, ϕ, ω, l, tspan, tstep, ics; kwargs...)
LDAEProblem(ϑ, f, u, g, ϕ, ω, l, tspan, tstep, q₀::StateVariable, p₀::StateVariable, λ₀::StateVariable = zero(q₀); kwargs...)

The function ϑ computes the momentum, f computes the force field, u and g compute the projections, and ϕ provides the algebraic constraint. The functions ψ, and are optional and provide the secondary constraint, that is the time derivative of the algebraic constraint, and the corresponding projection.

tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entries q and p. The initial conditions q₀, p₀, λ₀ and μ₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}. For the interfaces of the functions ϑ, f, u, g, ϕ, , , ψ, ω and l see LDAE.

In addition to the standard keyword arguments for EquationProblem subtypes, a LDAEProblem accepts functions and for the computation of initial guesses for the vector fields with default values v̄ = _ldae_default_v̄ and f̄ = f.

Function Definitions

The functions ϑ and f must have the interface

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

function f(f, t, q, v, params)
    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 u(u, t, q, v, p, μ, params)
    u[1] = ...
    u[2] = ...
    ...
end

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

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

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

and the functions ω and l, computing the symplectic matrix and the Lagrangian, have the following signature

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

function l(t, q, v, params)
    return ...
end

Some integrators also enforce the secondary constraint $\psi$ and require the following additional functions

function ū(u, t, q, v, p, μ, params)
    u[1] = ...
    u[2] = ...
    ...
end

function ḡ(g, t, q, v, p, μ, params)
    g[1] = ...
    g[2] = ...
    ...
end

function ψ(ψ, t, q, v, p, q̇, ṗ, params)
    ψ[1] = ...
end

With the above function definitions the LDAEProblem can be created by

tspan = (0.0, 1.0)
tstep = 0.1
q₀ = [1., 1.]
p₀ = [1., 0.]
λ₀ = [0.]
μ₀ = [0.]

prob = LDAEProblem(ϑ, f, u, g, ϕ, ω, l, tspan, tstep, q₀, p₀, λ₀)

or

prob = LDAEProblem(ϑ, f, u, g, ϕ, ū, ḡ, ψ, ω, l, tspan, tstep, q₀, p₀, λ₀, μ₀)
source

Stochastic Differential Equations

GeometricEquations.SDEProblemType

SDEProblem: Stratonovich Stochastic Differential Equation Problem

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

Constructors

SDEProblem(v, B, tspan, tstep, ics::NamedTuple; kwargs...)
SDEProblem(v, B, tspan, tstep, q₀::StateVariable; kwargs...)

where v is the function computing the vector field and B computes the diffusion matrix tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entry q. The initial condition q₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}.

For possible keyword arguments see the documentation on EquationProblem subtypes.

Function Definitions

The functions v and B, providing the drift vector field and diffusion matrix. The function v must have the interface

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

where t is the current time, q is the current solution vector, v is the vector which holds the result of evaluating the vector field $v$ on t and q, and params are additional parameters. The function B should have a method with interface

function B(B, t, q, params)
    B[1,1] = ...
    ...
end
source
GeometricEquations.PSDEProblemType

PSDEProblem: Stratonovich Partitioned Stochastic Differential Equation Problem

A partitioned stochastic differential equations is an initial value problem of the form

\[\begin{aligned} dq (t) &= v(t, q(t), p(t)) \, dt + B(t, q(t), p(t)) \circ dW , & q(t_{0}) &= q_{0} , \\ dp (t) &= f(t, q(t), p(t)) \, dt + G(t, q(t), p(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

Constructors

PSDEProblem(v, f, B, G, tspan, tstep, ics::NamedTuple; kwargs...)
PSDEProblem(v, f, B, G, tspan, tstep, q₀::StateVariable; p₀::StateVariable; kwargs...)

where v and f are the functions computing the vector field and B and G compute the diffusion matrices, tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entry q. The initial condition q₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}.

For possible keyword arguments see the documentation on EquationProblem subtypes.

Function Definitions

The functions v, f, B and G, providing the drift vector fields and diffusion matrices, each take five arguments, v(v, t, q, p, params), f(f, t, q, p, params), B(B, t, q, p, params) and G(G, t, q, p, params), where t is the current time, (q, p) is the current solution, 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), and params are optional parameters.

The corresponding methods should have the following signatures:

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

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

function B(B, t, q, p, params)
    B[1,1] = ...
    ...
end

function G(G, t, q, p, params)
    G[1,1] = ...
    ...
end
source
GeometricEquations.SPSDEProblemType

SPSDEProblem: Stratonovich Split Partitioned Stochastic Differential Equation Problem

Defines a partitioned stochastic differential initial value problem

\[\begin{aligned} dq (t) &= v(t, q(t), p(t)) \, dt + B(t, q(t), p(t)) \circ dW , & q(t_{0}) &= q_{0} , \\ dp (t) &= [ f_1(t, q(t), p(t)) + f_2(t, q(t), p(t)) ] \, dt + [ G_1(t, q(t), p(t)) + G_2(t, q(t), p(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

Constructors

SPSDEProblem(v, f1, f2, B, G1, G2, tspan, tstep, ics::NamedTuple; kwargs...)
SPSDEProblem(v, f1, f2, B, G1, G2, tspan, tstep, q₀::StateVariable; p₀::StateVariable; kwargs...)

where v and f are the functions computing the vector field and Bᵢ and Gᵢ compute the diffusion matrices, tspan is the time interval (t₀,t₁) for the problem to be solved in, tstep is the time step to be used in the simulation, and ics is a NamedTuple with entry q. The initial condition q₀ can also be prescribed directly, with StateVariable an AbstractArray{<:Number}.

For possible keyword arguments see the documentation on EquationProblem subtypes.

Function Definitions

The functions v, f1, f2, B, G1 and G2, providing the drift vector fields and diffusion matrices, all take five arguments, (out, t, q, p, params).

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

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

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

function B(B, t, q, p, params)
    B[1,1] = ...
    ...
end

function G1(G, t, q, p, params)
    G[1,1] = ...
    ...
end

function G2(G, t, q, p, params)
    G[1,1] = ...
    ...
end

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_i$, $G_i$ on (t,q,p).

source