Integrators
Common
GeometricIntegrators.Solutions.AtomicSolution
— MethodCreate AtomicSolution for DAE.
GeometricIntegrators.Solutions.AtomicSolution
— MethodCreate AtomicSolution for ODE.
GeometricIntegrators.Solutions.AtomicSolution
— MethodCreate AtomicSolution for partitioned DAE.
GeometricIntegrators.Solutions.AtomicSolution
— MethodCreate AtomicSolution for partitioned ODE.
GeometricIntegrators.Solutions.AtomicSolution
— MethodCreate AtomicSolution for partitioned SDE.
GeometricIntegrators.Solutions.AtomicSolution
— MethodCreate AtomicSolution for SDE.
GeometricIntegrators.Integrators.AbstractTableau
— TypeHolds the information for the various methods' tableaus.
GeometricIntegrators.Integrators.create_internal_stage_matrix
— MethodCreate a vector of S solution matrices of type DT to store the solution of S internal stages for a problem with DxM
dimensions.
GeometricIntegrators.Integrators.create_internal_stage_vector
— MethodCreate a vector of S solution vectors of type DT to store the solution of S internal stages for a problem with D
dimensions.
GeometricIntegrators.Integrators.create_internal_stage_vector_with_zero
— MethodCreate a vector of (S,M+1) solution vectors of type DT to store the solution of S internal stages and M random processes for a problem with D
dimensions.
GeometricIntegrators.Integrators.create_internal_stage_vector_with_zero
— MethodCreate a vector of S+1 solution vectors of type DT to store the solution of S internal stages and the solution of the previous timestep for a problem with D
dimensions.
GeometricIntegrators.Integrators.create_nonlinear_solver
— FunctionCreate nonlinear solver object for a system of N
equations with data type DT
. The function $f(x)=0$ to be solved for is determined by a julia function function_stages!(x, b, params)
, where x
is the current solution and b
is the output vector, s.th. $b = f(x)$. params
are a set of parameters depending on the equation and integrator that is used. The solver type is obtained from the config dictionary (:nls_solver
).
GeometricIntegrators.Integrators.Integrator
— MethodPrint error for integrators not implemented, yet.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for special partitioned additive Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for variational partitioned additive Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for variational special partitioned additive Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for variational special partitioned additive Runge-Kutta tableau with projection on primary constraint.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for Projected Gauss-Legendre Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for variational partitioned Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for diagonally implicit Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for explicit Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for fully implicit Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for Hamiltonian partitioned additive Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for Hamiltonian special partitioned additive Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for Hamiltonian special partitioned additive Runge-Kutta tableau with projection on primary constraint.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for explicit partitioned Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for implicit partitioned Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for stochastic fully implicit partitioned Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for stochastic explicit Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for stochastic fully implicit Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for weak explicit Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for weak fully implicit Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for splitting tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for stochastic fully implicit split partitioned Runge-Kutta tableau.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for variational special partitioned additive Runge-Kutta tableau with projection on secondary constraint.
GeometricIntegrators.Integrators.Integrator
— MethodCreate integrator for formal Lagrangian Runge-Kutta tableau.
GeometricIntegrators.Integrators.IntegratorConstructor
— MethodCreate integrator constructor for diagonally implicit Runge-Kutta tableau.
GeometricIntegrators.Integrators.IntegratorConstructor
— MethodCreate integrator constructor for explicit Runge-Kutta tableau.
GeometricIntegrators.Integrators.IntegratorConstructor
— MethodCreate integrator constructor for fully implicit Runge-Kutta tableau.
GeometricIntegrators.Integrators.IntegratorConstructor
— MethodCreate integrator constructor for exact solution.
GeometricIntegrators.Integrators.integrate!
— MethodIntegrate ODE for initial conditions m with m₁ ≤ m ≤ m₂.
GeometricIntegrators.Integrators.integrate!
— MethodIntegrate equation for all initial conditions.
GeometricIntegrators.Integrators.integrate!
— MethodIntegrate ODE for initial conditions m with m₁ ≤ m ≤ m₂ for time steps n with n₁ ≤ n ≤ n₂.
GeometricIntegrators.Integrators.integrate
— MethodIntegrate given equation with given tableau for ntime time steps and return solution.
GeometricIntegrators.Integrators.integrate
— MethodApply integrator for ntime time steps and return solution.
GeometricIntegrators.Integrators.integrate
— MethodIntegrate ODE specified by vector field and initial condition with given tableau for ntime time steps and return solution.
GeometricIntegrators.Integrators.integrate
— MethodIntegrate PODE specified by two vector fields and initial conditions with given tableau for ntime time steps and return solution.
Initial Guesses
GeometricIntegrators.Integrators.aitken_neville
— MethodCompute p(x) where p is the unique polynomial of degree length(xi), such that p(x[i]) = y[i]) for all i.
ti: interpolation nodes
xi: interpolation values
t: evaluation point
x: evaluation value
GeometricIntegrators.Integrators.euler_extrapolation
— MethodEuler extrapolation method with arbitrary order p.
v: function to compute vector field
t₀: initial time
t₁: final time
x₀: initial value
x₁: final value
s: number of interpolations (order p=s+1)
TODO This is probably broken!
GeometricIntegrators.Integrators.midpoint_extrapolation
— MethodMidpoint extrapolation method with arbitrary order p.
v: function to compute vector field
f: function to compute force field
t₀: initial time
t₁: final time
q₀: initial positions
p₀: initial momenta
q₁: final positions
p₁: final momenta
s: number of interpolations (order p=2s+2)
GeometricIntegrators.Integrators.midpoint_extrapolation
— MethodMidpoint extrapolation method with arbitrary order p.
v: function to compute vector field
t₀: initial time
t₁: final time
x₀: initial value
x₁: final value
s: number of interpolations (order p=2s+2)
GeometricIntegrators.Integrators.InitialGuessODE
— TypeInitialGuessODE
: Initial guess for ordinary differential equations
Fields
int
: interpolation structurev
: vector fieldΔt
: time steps
: number of extrapolation stages (for initialisation)
GeometricIntegrators.Integrators.initialize!
— MethodInitialise initial guess, i.e., given t₀, t₁, q₁ compute q₀, v₀, v₁.
GeometricIntegrators.Integrators.update_vector_fields!
— Methodcompute vector field of new solution
GeometricIntegrators.Integrators.InitialGuessPODE
— TypeInitialGuessPODE
: Initial guess for partitioned ordinary differential equations
Fields
int
: interpolation structurev
: vector field for qf
: vector field for pΔt
: time steps
: number of extrapolation stages (for initialisation)
GeometricIntegrators.Integrators.initialize!
— MethodInitialise initial guess, i.e., given t₀, t₁, q₁ compute q₀, v₀, v₁.
Splitting Methods
GeometricIntegrators.Integrators.IntegratorSplitting
— TypeSplitting integrator.
GeometricIntegrators.Integrators.IntegratorSplitting
— MethodConstruct splitting integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with splitting integrator.
Runge-Kutta Methods
GeometricIntegrators.Integrators.CoefficientsRK
— TypeHolds the coefficients of a Runge-Kutta method.
Base.show
— MethodPrint Runge-Kutta coefficients.
GeometricIntegrators.Integrators.AbstractTableauERK
— TypeHolds the tableau of an explicit Runge-Kutta method.
GeometricIntegrators.Integrators.AbstractTableauIRK
— TypeHolds the tableau of an implicit Runge-Kutta method.
GeometricIntegrators.Integrators.AbstractTableauPRK
— TypeHolds the tableau of a partitioned Runge-Kutta method.
GeometricIntegrators.Integrators.AbstractTableauRK
— TypeHolds the tableau of a Runge-Kutta method.
GeometricIntegrators.Integrators.readTableauRKHeaderFromFile
— MethodReads and parses Tableau metadata from file.
GeometricIntegrators.Integrators.writeTableauToFile
— MethodWrite Runge-Kutta tableau to file.
GeometricIntegrators.Integrators.IntegratorCacheERK
— TypeExplicit Runge-Kutta integrator cache.
GeometricIntegrators.Integrators.IntegratorERK
— TypeExplicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.ParametersERK
— TypeParameters for right-hand side function of explicit Runge-Kutta methods.
GeometricIntegrators.Integrators.TableauERK
— TypeHolds the tableau of an explicit Runge-Kutta method.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with explicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.readTableauERKFromFile
— MethodRead explicit Runge-Kutta tableau from file.
GeometricIntegrators.Integrators.IntegratorCacheDIRK
— TypeDiagonally implicit Runge-Kutta integrator cache.
GeometricIntegrators.Integrators.IntegratorDIRK
— TypeDiagonally implicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.ParametersDIRK
— TypeParameters for right-hand side function of diagonally implicit Runge-Kutta methods.
GeometricIntegrators.Integrators.TableauDIRK
— TypeHolds the tableau of a diagonally implicit Runge-Kutta method.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of fully implicit Runge-Kutta methods.
GeometricIntegrators.Integrators.initial_guess!
— MethodCompute initial guess for internal stages.
GeometricIntegrators.Integrators.initialize!
— MethodInitialise initial guess
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with diagonally implicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.IntegratorCacheFIRK
— TypeFully implicit Runge-Kutta integrator cache.
Fields
q̃
: initial guess of solutionṽ
: initial guess of vector fields̃
: holds shift due to periodicity of solutionQ
: internal stages of solutionV
: internal stages of vector fieldY
: vector field of internal stages
GeometricIntegrators.Integrators.IntegratorFIRK
— TypeFully implicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.ParametersFIRK
— TypeParameters for right-hand side function of fully implicit Runge-Kutta methods.
GeometricIntegrators.Integrators.TableauFIRK
— TypeHolds the tableau of a fully implicit Runge-Kutta method.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of fully implicit Runge-Kutta methods.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with fully implicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.jacobian!
— MethodCompute stages of fully implicit Runge-Kutta methods.
GeometricIntegrators.Integrators.IntegratorCacheEPRK
— TypeExplicit Runge-Kutta integrator cache.
GeometricIntegrators.Integrators.IntegratorEPRK
— TypeExplicit partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.ParametersEPRK
— TypeParameters for right-hand side function of explicit partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.TableauEPRK
— TypeTableauEPRK
: Tableau of an Explicit Partitioned Runge-Kutta method
\[\begin{aligned} V_{n,i} &= \hphantom{-} \dfrac{\partial H}{\partial p} (Q_{n,i}, P_{n,i}) , & Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} \, V_{n,j} , & q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} \, V_{n,i} , \\ F_{n,i} &= - \dfrac{\partial H}{\partial q} (Q_{n,i}, P_{n,i}) , & P_{n,i} &= p_{n} + h \sum \limits_{i=1}^{s} \bar{a}_{ij} \, F_{n,j} , & p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} \bar{b}_{i} \, F_{n,i} , \end{aligned}\]
usually satisfying the symplecticity conditions
\[\begin{aligned} b_{i} \bar{a}_{ij} + b_{j} a_{ji} &= b_{i} b_{j} , & \bar{b}_i &= b_i . \end{aligned}\]
GeometricIntegrators.Integrators.computeStageP!
— MethodCompute P stages of explicit partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.computeStageQ!
— MethodCompute Q stages of explicit partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate partitioned ODE with explicit partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.IntegratorCacheIPRK
— TypeImplicit partitioned Runge-Kutta integrator cache.
Fields
q̃
: initial guess of qp̃
: initial guess of pṽ
: initial guess of vf̃
: initial guess of fs̃
: holds shift due to periodicity of solutionQ
: internal stages of qP
: internal stages of pV
: internal stages of vF
: internal stages of fY
: vector field of internal stages of qZ
: vector field of internal stages of p
GeometricIntegrators.Integrators.IntegratorIPRK
— TypeImplicit partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.ParametersIPRK
— TypeParameters for right-hand side function of implicit partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.TableauIPRK
— TypeTableauIPRK
: Tableau of an Implicit Partitioned Runge-Kutta method
\[\begin{aligned} P_{n,i} &= \dfrac{\partial L}{\partial v} (Q_{n,i}, V_{n,i}) , & Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} \, V_{n,j} , & q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} \, V_{n,i} , \\ F_{n,i} &= \dfrac{\partial L}{\partial q} (Q_{n,i}, V_{n,i}) , & P_{n,i} &= p_{n} + h \sum \limits_{i=1}^{s} \bar{a}_{ij} \, F_{n,j} , & p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} \bar{b}_{i} \, F_{n,i} , \end{aligned}\]
usually satisfying the symplecticity conditions
\[\begin{aligned} b_{i} \bar{a}_{ij} + b_{j} a_{ji} &= b_{i} b_{j} , & \bar{b}_i &= b_i . \end{aligned}\]
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of implicit partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with implicit partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.IntegratorFLRK
— TypeFormal Lagrangian Runge-Kutta integrator.
GeometricIntegrators.Integrators.ParametersFLRK
— TypeParameters for right-hand side function of formal Lagrangian Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of formal Lagrangian Runge-Kutta methods.
GeometricIntegrators.Integrators.CoefficientsPGLRK
— TypeHolds the coefficients of a projected Gauss-Legendre Runge-Kutta method.
GeometricIntegrators.Integrators.IntegratorPGLRK
— TypeProjected Gauss-Legendre Runge-Kutta integrator.
Reference: LUIGI BRUGNANO, FELICE IAVERNARO, AND DONATO TRIGIANTE.
ENERGY- AND QUADRATIC INVARIANTS–PRESERVING INTEGRATORS BASED
UPON GAUSS COLLOCATION FORMULAE.
SIAM J. NUMER. ANAL. Vol. 50, No. 6, pp. 2897–2916, 2012.
GeometricIntegrators.Integrators.ParametersPGLRK
— TypeParameters for right-hand side function of projected Gauss-Legendre Runge-Kutta methods.
Base.show
— MethodPrint Runge-Kutta coefficients.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of projected Gauss-Legendre Runge-Kutta methods.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with projected Gauss-Legendre Runge-Kutta integrator.
SPARK Methods
GeometricIntegrators.Integrators.SPARK.AbstractParametersSPARK
— TypeParameters for right-hand side function of Specialised Partitioned Additive Runge-Kutta methods.
GeometricIntegrators.Integrators.SPARK.AbstractTableauSPARK
— TypeHolds the tableau of an Specialised Partitioned Additive Runge-Kutta method for Variational systems.
GeometricIntegrators.Integrators.SPARK.CoefficientsARK
— TypeHolds the coefficients of an additive Runge-Kutta method.
GeometricIntegrators.Integrators.SPARK.CoefficientsIRK
— TypeHolds the coefficients of an additive Runge-Kutta method.
GeometricIntegrators.Integrators.SPARK.CoefficientsMRK
— TypeHolds the multiplier Runge-Kutta coefficients.
GeometricIntegrators.Integrators.SPARK.CoefficientsPRK
— TypeHolds the coefficients of a projective Runge-Kutta method.
GeometricIntegrators.Integrators.SPARK.CoefficientsSPARK
— TypeHolds the coefficients of a SPARK method.
GeometricIntegrators.Integrators.SPARK.IntegratorCacheSPARK
— TypeCache of a Specialised Partitioned Additive Runge-Kutta integrator.
Fields
n
: time step numbert
: time of current time stept̅
: time of previous time stepq
: current solution of qq̅
: previous solution of qp
: current solution of pp̅
: previous solution of pv
: vector field of qv̅
: vector field of q̅f
: vector field of pf̅
: vector field of p̅q̃
: initial guess of qp̃
: initial guess of pṽ
: initial guess of vf̃
: initial guess of fs̃
: holds shift due to periodicity of solutionQ
: internal stages of qP
: internal stages of pV
: internal stages of vF
: internal stages of fY
: vector field of internal stages of qZ
: vector field of internal stages of p
GeometricIntegrators.Integrators.SPARK.IntegratorHPARK
— TypePartitioned Additive Runge-Kutta integrator for Hamiltonian systems subject to Dirac constraints EXPERIMENTAL.
This integrator solves the following system of equations for the internal stages,
\[\begin{aligned} Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} U_{n,j} , & i &= 1, ..., s , \\ P_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} a_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} G_{n,j} , & i &= 1, ..., s , \\ \tilde{Q}_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} U_{n,j} , & i &= 1, ..., r , \\ \tilde{P}_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} G_{n,j} , & i &= 1, ..., r , \\ \tilde{\Phi}_{n,i} &= 0 , & i &= 1, ..., r , \end{aligned}\]
with definitions
\[\begin{aligned} V_{n,i} &= \hphantom{-} \frac{\partial H}{\partial p} (Q_{n,i}, P_{n,i}) , & i &= 1, ..., s , \\ F_{n,i} &= - \frac{\partial H}{\partial q} (Q_{n,i}, P_{n,i}) , & i &= 1, ..., s , \\ U_{n,i} &= \hphantom{-} \frac{\partial \phi}{\partial p} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ G_{n,i} &= - \frac{\partial \phi}{\partial q} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ \tilde{\Phi}_{n,i} &= \phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) , & i &= 1, ..., r , \end{aligned}\]
and update rule
\[\begin{aligned} q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} V_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} U_{n,i} , \\ p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} b_{i} F_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} G_{n,i} . \end{aligned}\]
GeometricIntegrators.Integrators.SPARK.IntegratorHSPARK
— TypeSpecialised Partitioned Additive Runge-Kutta integrator for Hamiltonian systems EXPERIMENTAL.
This integrator solves the following system of equations for the internal stages,
\[\begin{aligned} Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} U_{n,j} , & i &= 1, ..., s , \\ P_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} a_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} G_{n,j} , & i &= 1, ..., s , \\ \tilde{Q}_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} U_{n,j} , & i &= 1, ..., r , \\ \tilde{P}_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} G_{n,j} , & i &= 1, ..., r , \\ 0 &= \sum \limits_{j=1}^{r} \omega_{ij} \tilde{\Phi}_{n,j} , & i &= 1, ..., r-1 , \end{aligned}\]
with definitions
\[\begin{aligned} V_{n,i} &= \hphantom{-} \frac{\partial H}{\partial p} (Q_{n,i}, P_{n,i}) , & i &= 1, ..., s , \\ F_{n,i} &= - \frac{\partial H}{\partial q} (Q_{n,i}, P_{n,i}) , & i &= 1, ..., s , \\ U_{n,i} &= \hphantom{-} \frac{\partial \phi}{\partial p} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ G_{n,i} &= - \frac{\partial \phi}{\partial q} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ \tilde{\Phi}_{n,i} &= \phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) , & i &= 1, ..., r , \end{aligned}\]
and update rule
\[\begin{aligned} q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} V_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} U_{n,i} , \\ p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} b_{i} F_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} G_{n,i} , \\ 0 &= \phi (q_{n+1}, p_{n+1}) . \end{aligned}\]
GeometricIntegrators.Integrators.SPARK.IntegratorHSPARKprimary
— TypeSpecialised Partitioned Additive Runge-Kutta integrator for Hamiltonian systems EXPERIMENTAL.
This integrator solves the following system of equations for the internal stages,
\[\begin{aligned} Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} U_{n,j} , & i &= 1, ..., s , \\ P_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} a_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} G_{n,j} , & i &= 1, ..., s , \\ \tilde{Q}_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} U_{n,j} , & i &= 1, ..., r , \\ \tilde{P}_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} G_{n,j} , & i &= 1, ..., r , \\ 0 &= \sum \limits_{j=1}^{r} \omega_{ij} \tilde{\Phi}_{n,j} , & i &= 1, ..., r-1 , \end{aligned}\]
with definitions
\[\begin{aligned} V_{n,i} &= \hphantom{-} \frac{\partial H}{\partial p} (Q_{n,i}, P_{n,i}) , & i &= 1, ..., s , \\ F_{n,i} &= - \frac{\partial H}{\partial q} (Q_{n,i}, P_{n,i}) , & i &= 1, ..., s , \\ U_{n,i} &= \hphantom{-} \frac{\partial \phi}{\partial p} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ G_{n,i} &= - \frac{\partial \phi}{\partial q} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ \tilde{\Phi}_{n,i} &= \phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) , & i &= 1, ..., r , \end{aligned}\]
and update rule
\[\begin{aligned} q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} V_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} U_{n,i} , \\ p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} b_{i} F_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} G_{n,i} , \\ 0 &= \phi (q_{n+1}, p_{n+1}) . \end{aligned}\]
GeometricIntegrators.Integrators.SPARK.IntegratorHSPARKsecondary
— TypeSpecialised Partitioned Additive Runge-Kutta integrator for Hamiltonian systems subject to Dirac constraints with projection on secondary constraint EXPERIMENTAL.
This integrator solves the following system of equations for the internal stages,
\[\begin{aligned} Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} U_{n,j} , & i &= 1, ..., s , \\ P_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} a_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} G_{n,j} , & i &= 1, ..., s , \\ \tilde{Q}_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} U_{n,j} , & i &= 1, ..., r , \\ \tilde{P}_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} G_{n,j} , & i &= 1, ..., r , \\ 0 &= \sum \limits_{j=1}^{r} \omega_{ij} \tilde{\Phi}_{n,j} , & i &= 1, ..., r-1 , \end{aligned}\]
with definitions
\[\begin{aligned} V_{n,i} &= \hphantom{-} \frac{\partial H}{\partial p} (Q_{n,i}, P_{n,i}) , & i &= 1, ..., s , \\ F_{n,i} &= - \frac{\partial H}{\partial q} (Q_{n,i}, P_{n,i}) , & i &= 1, ..., s , \\ U_{n,i} &= \hphantom{-} \frac{\partial \phi}{\partial p} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ G_{n,i} &= - \frac{\partial \phi}{\partial q} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ \tilde{\Phi}_{n,i} &= \phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) , & i &= 1, ..., r , \end{aligned}\]
and update rule
\[\begin{aligned} q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} V_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} U_{n,i} , \\ p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} b_{i} F_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} G_{n,i} , \\ 0 &= \phi (q_{n+1}, p_{n+1}) . \end{aligned}\]
GeometricIntegrators.Integrators.SPARK.IntegratorSLRK
— TypeSpecialised Partitioned Additive Runge-Kutta integrator for degenerate variational systems with projection on secondary constraint.
This integrator solves the following system of equations for the internal stages,
\[\begin{aligned} Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a^1_{ij} V_{n,j} + h \sum \limits_{j=1}^{\sigma} a^2_{ij} \tilde{\Lambda}_{n,j} , & i &= 1, ..., s , \\ P_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \bar{a}^1_{ij} F^1_{n,j} + h \sum \limits_{j=1}^{s} \bar{a}^2_{ij} F^2_{n,j} + h \sum \limits_{j=1}^{\sigma} \bar{a}^3_{ij} \tilde{F}^3_{n,j} , & i &= 1, ..., s , \\ 0 &= \Phi_{n,i} , & i &= 1, ..., s , \end{aligned}\]
the projective stages
\[\begin{aligned} \tilde{Q}_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} \alpha^1_{ij} \tilde{V}_{n,j} + h \sum \limits_{j=1}^{\sigma} \alpha^2_{ij} \tilde{\Lambda}_{n,j} , & i &= 1, ..., \sigma , \\ \tilde{P}_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \bar{\alpha}^1_{ij} F^1_{n,j} + h \sum \limits_{j=1}^{s} \bar{\alpha}^2_{ij} \tilde{F}^2_{n,j} + h \sum \limits_{j=1}^{\sigma} \bar{\alpha}^3_{ij} \tilde{F}^3_{n,j} , & i &= 1, ..., \sigma , \\ 0 &= \tilde{\Phi}_{n,i} , & i &= 1, ..., \sigma , \\ 0 &= \sum \limits_{j=1}^{\sigma} \omega_{ij} \tilde{\Psi}_{n,i} , & i &= 1, ..., \sigma-1 , \end{aligned}\]
and update rule
\[\begin{aligned} q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b^1_{i} V_{n,i} + h \sum \limits_{i=1}^{\sigma} b^2_{i} \tilde{\Lambda}_{n,i} , \\ p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} b^1_{i} F^1_{n,i} + h \sum \limits_{i=1}^{s} b^2_{i} F^2_{n,i} + h \sum \limits_{i=1}^{\sigma} b^3_{i} \tilde{F}^3_{n,i} , \\ 0 &= \phi (q_{n+1}, p_{n+1}) , \end{aligned}\]
with definitions
\[\begin{aligned} F^1_{n,i} &= \nabla H (Q_{n,i}) , & i &= 1, ..., s , \\ F^2_{n,i} &= \nabla \vartheta (Q_{n,i}) \cdot V_{n,i} , & i &= 1, ..., s , \\ \tilde{F}^3_{n,i} &= \nabla \phi (\tilde{Q}_{n,i}, \tilde{P}_{n,i}) \cdot \tilde{\Lambda}_{n,i} , & i &= 1, ..., \sigma , \\ \Phi_{n,i} &= \phi (Q_{n,i}, P_{n,i}) = P_{n,i} - \vartheta (Q_{n,i}) , & i &= 1, ..., s , \\ \tilde{\Phi}_{n,i} &= \phi (\tilde{Q}_{n,i}, \tilde{P}_{n,i}) = \tilde{P}_{n,i} - \vartheta (\tilde{Q}_{n,i}) , & i &= 1, ..., \sigma , \\ \tilde{\Psi}_{n,i} &= \psi (\tilde{Q}_{n,i}, \tilde{V}_{n,i}, \tilde{P}_{n,i}, \tilde{F}_{n,i}) = \tilde{F}_{n,i} - \tilde{V}_{n,i} \cdot \nabla \vartheta (\tilde{Q}_{n,i}) , & i &= 1, ..., \sigma , \end{aligned}\]
so that
\[\begin{aligned} F^1_{n,i} + F^2_{n,i} &= \frac{\partial L}{\partial q} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ \phi (Q_{n,i}, P_{n,i}) &= P_{n,i} - \frac{\partial L}{\partial v} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ \psi (\tilde{Q}_{n,i}, \tilde{V}_{n,i}, \tilde{P}_{n,i}, \tilde{F}_{n,i}) &= \tilde{F}_{n,i} - \tilde{V}_{n,i} \cdot \nabla \frac{\partial L}{\partial v} (\tilde{Q}_{n,i}, \tilde{V}_{n,i}) , & i &= 1, ..., \sigma . \end{aligned}\]
GeometricIntegrators.Integrators.SPARK.IntegratorSPARK
— TypeSpecialised Partitioned Additive Runge-Kutta integrator for index-two DAE systems.
This integrator solves the following system of equations for the internal stages,
\[\begin{aligned} Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} U_{n,j} , & i &= 1, ..., s , \\ P_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} a_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} G_{n,j} , & i &= 1, ..., s , \\ \tilde{Q}_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} U_{n,j} , & i &= 1, ..., r , \\ \tilde{P}_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} G_{n,j} , & i &= 1, ..., r , \\ 0 &= \sum \limits_{j=1}^{r} \omega_{ij} \tilde{\Phi}_{n,j} , & i &= 1, ..., r-1 , \end{aligned}\]
with definitions
\[\begin{aligned} P_{n,i} &= \frac{\partial L}{\partial v} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ F_{n,i} &= \frac{\partial L}{\partial q} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ U_{n,i} &= \hphantom{-} \frac{\partial \phi}{\partial p} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ G_{n,i} &= - \frac{\partial \phi}{\partial q} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ \tilde{\Phi}_{n,i} &= \phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) , & i &= 1, ..., r , \end{aligned}\]
and update rule
\[\begin{aligned} q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} V_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} U_{n,i} , \\ p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} b_{i} F_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} G_{n,i} , \\ 0 &= \phi (q_{n+1}, p_{n+1}) . \end{aligned}\]
GeometricIntegrators.Integrators.SPARK.IntegratorVPARK
— TypeVariational partitioned additive Runge-Kutta integrator.
This integrator solves the following system of equations for the internal stages,
\[\begin{aligned} Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} U_{n,j} , & i &= 1, ..., s , \\ P_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} a_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} G_{n,j} , & i &= 1, ..., s , \\ \tilde{Q}_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} U_{n,j} , & i &= 1, ..., r , \\ \tilde{P}_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} G_{n,j} , & i &= 1, ..., r , \\ \tilde{\Phi}_{n,i} &= 0 , & i &= 1, ..., r , \end{aligned}\]
with definitions
\[\begin{aligned} P_{n,i} &= \frac{\partial L}{\partial v} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ F_{n,i} &= \frac{\partial L}{\partial q} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ U_{n,i} &= \hphantom{-} \frac{\partial \phi}{\partial p} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ G_{n,i} &= - \frac{\partial \phi}{\partial q} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ \tilde{\Phi}_{n,i} &= \phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) , & i &= 1, ..., r , \end{aligned}\]
and update rule
\[\begin{aligned} q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} V_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} U_{n,i} , \\ p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} b_{i} F_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} G_{n,i} . \end{aligned}\]
GeometricIntegrators.Integrators.SPARK.IntegratorVSPARK
— TypeSpecialised Partitioned Additive Runge-Kutta integrator for Variational systems.
This integrator solves the following system of equations for the internal stages,
\[\begin{aligned} Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} U_{n,j} , & i &= 1, ..., s , \\ P_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} a_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} G_{n,j} , & i &= 1, ..., s , \\ \tilde{Q}_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} U_{n,j} , & i &= 1, ..., r , \\ \tilde{P}_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} G_{n,j} , & i &= 1, ..., r , \\ 0 &= \sum \limits_{j=1}^{r} \omega_{ij} \tilde{\Phi}_{n,j} , & i &= 1, ..., r-1 , \end{aligned}\]
with definitions
\[\begin{aligned} P_{n,i} &= \frac{\partial L}{\partial v} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ F_{n,i} &= \frac{\partial L}{\partial q} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ U_{n,i} &= \hphantom{-} \frac{\partial \phi}{\partial p} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ G_{n,i} &= - \frac{\partial \phi}{\partial q} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ \tilde{\Phi}_{n,i} &= \phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) , & i &= 1, ..., r , \end{aligned}\]
and update rule
\[\begin{aligned} q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} V_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} U_{n,i} , \\ p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} b_{i} F_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} G_{n,i} , \\ 0 &= \phi (q_{n+1}, p_{n+1}) . \end{aligned}\]
GeometricIntegrators.Integrators.SPARK.IntegratorVSPARKprimary
— TypeSpecialised Partitioned Additive Runge-Kutta integrator for Variational systems.
This integrator solves the following system of equations for the internal stages,
\[\begin{aligned} Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} U_{n,j} , & i &= 1, ..., s , \\ P_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} a_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \alpha_{ij} G_{n,j} , & i &= 1, ..., s , \\ \tilde{Q}_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} V_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} U_{n,j} , & i &= 1, ..., r , \\ \tilde{P}_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \tilde{a}_{ij} F_{n,j} + h \sum \limits_{j=1}^{r} \tilde{\alpha}_{ij} G_{n,j} , & i &= 1, ..., r , \\ 0 &= \sum \limits_{j=1}^{r} \omega_{ij} \tilde{\Phi}_{n,j} , & i &= 1, ..., r-1 , \\ 0 &= \sum \limits_{i=1}^{r} \tilde{d}_i \, \Lambda_{n,i} , \end{aligned}\]
with definitions
\[\begin{aligned} P_{n,i} &= \frac{\partial L}{\partial v} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ F_{n,i} &= \frac{\partial L}{\partial q} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ U_{n,i} &= \hphantom{-} \frac{\partial \phi}{\partial p} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ G_{n,i} &= - \frac{\partial \phi}{\partial q} (\tilde{Q}_{n,i}, \tilde{P}_{n,i})^{T} \Lambda_{n,i} , & i &= 1, ..., r , \\ \tilde{\Phi}_{n,i} &= \phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) , & i &= 1, ..., r , \end{aligned}\]
and update rule
\[\begin{aligned} q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} V_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} U_{n,i} , \\ p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} b_{i} F_{n,i} + h \sum \limits_{i=1}^{r} \beta_{i} G_{n,i} . \end{aligned}\]
GeometricIntegrators.Integrators.SPARK.IntegratorVSPARKsecondary
— TypeSpecialised Partitioned Additive Runge-Kutta integrator for degenerate variational systems with projection on secondary constraint.
This integrator solves the following system of equations for the internal stages,
\[\begin{aligned} Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a^1_{ij} V_{n,j} + h \sum \limits_{j=1}^{\sigma} a^2_{ij} \tilde{\Lambda}_{n,j} , & i &= 1, ..., s , \\ P_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \bar{a}^1_{ij} F^1_{n,j} + h \sum \limits_{j=1}^{s} \bar{a}^2_{ij} F^2_{n,j} + h \sum \limits_{j=1}^{\sigma} \bar{a}^3_{ij} \tilde{F}^3_{n,j} , & i &= 1, ..., s , \\ 0 &= \Phi_{n,i} , & i &= 1, ..., s , \end{aligned}\]
the projective stages
\[\begin{aligned} \tilde{Q}_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} \alpha^1_{ij} \tilde{V}_{n,j} + h \sum \limits_{j=1}^{\sigma} \alpha^2_{ij} \tilde{\Lambda}_{n,j} , & i &= 1, ..., \sigma , \\ \tilde{P}_{n,i} &= p_{n} + h \sum \limits_{j=1}^{s} \bar{\alpha}^1_{ij} F^1_{n,j} + h \sum \limits_{j=1}^{s} \bar{\alpha}^2_{ij} \tilde{F}^2_{n,j} + h \sum \limits_{j=1}^{\sigma} \bar{\alpha}^3_{ij} \tilde{F}^3_{n,j} , & i &= 1, ..., \sigma , \\ 0 &= \tilde{\Phi}_{n,i} , & i &= 1, ..., \sigma , \\ 0 &= \sum \limits_{j=1}^{\sigma} \omega_{ij} \tilde{\Psi}_{n,i} , & i &= 1, ..., \sigma-1 , \end{aligned}\]
and update rule
\[\begin{aligned} q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b^1_{i} V_{n,i} + h \sum \limits_{i=1}^{\sigma} b^2_{i} \tilde{\Lambda}_{n,i} , \\ p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} b^1_{i} F^1_{n,i} + h \sum \limits_{i=1}^{s} b^2_{i} F^2_{n,i} + h \sum \limits_{i=1}^{\sigma} b^3_{i} \tilde{F}^3_{n,i} , \\ 0 &= \phi (q_{n+1}, p_{n+1}) , \end{aligned}\]
with definitions
\[\begin{aligned} F^1_{n,i} &= \nabla H (Q_{n,i}) , & i &= 1, ..., s , \\ F^2_{n,i} &= \nabla \vartheta (Q_{n,i}) \cdot V_{n,i} , & i &= 1, ..., s , \\ \tilde{F}^3_{n,i} &= \nabla \phi (\tilde{Q}_{n,i}, \tilde{P}_{n,i}) \cdot \tilde{\Lambda}_{n,i} , & i &= 1, ..., \sigma , \\ \Phi_{n,i} &= \phi (Q_{n,i}, P_{n,i}) = P_{n,i} - \vartheta (Q_{n,i}) , & i &= 1, ..., s , \\ \tilde{\Phi}_{n,i} &= \phi (\tilde{Q}_{n,i}, \tilde{P}_{n,i}) = \tilde{P}_{n,i} - \vartheta (\tilde{Q}_{n,i}) , & i &= 1, ..., \sigma , \\ \tilde{\Psi}_{n,i} &= \psi (\tilde{Q}_{n,i}, \tilde{V}_{n,i}, \tilde{P}_{n,i}, \tilde{F}_{n,i}) = \tilde{F}_{n,i} - \tilde{V}_{n,i} \cdot \nabla \vartheta (\tilde{Q}_{n,i}) , & i &= 1, ..., \sigma , \end{aligned}\]
so that
\[\begin{aligned} F^1_{n,i} + F^2_{n,i} &= \frac{\partial L}{\partial q} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ \phi (Q_{n,i}, P_{n,i}) &= P_{n,i} - \frac{\partial L}{\partial v} (Q_{n,i}, V_{n,i}) , & i &= 1, ..., s , \\ \psi (\tilde{Q}_{n,i}, \tilde{V}_{n,i}, \tilde{P}_{n,i}, \tilde{F}_{n,i}) &= \tilde{F}_{n,i} - \tilde{V}_{n,i} \cdot \nabla \frac{\partial L}{\partial v} (\tilde{Q}_{n,i}, \tilde{V}_{n,i}) , & i &= 1, ..., \sigma . \end{aligned}\]
GeometricIntegrators.Integrators.SPARK.ParametersHPARK
— TypeParameters for right-hand side function of Hamiltonian Partitioned Additive Runge-Kutta methods.
GeometricIntegrators.Integrators.SPARK.ParametersSLRK
— TypeParameters for right-hand side function of Specialised Partitioned Additive Runge-Kutta methods for Variational systems.
GeometricIntegrators.Integrators.SPARK.ParametersVPARK
— TypeParameters for right-hand side function of Variational Partitioned Additive Runge-Kutta methods.
GeometricIntegrators.Integrators.SPARK.ParametersVSPARKsecondary
— TypeParameters for right-hand side function of Specialised Partitioned Additive Runge-Kutta methods for Variational systems.
GeometricIntegrators.Integrators.SPARK.TableauHPARK
— TypeHolds the tableau of a Hamiltonian Partitioned Additive Runge-Kutta methods.
GeometricIntegrators.Integrators.SPARK.TableauSLRK
— TypeHolds the tableau of an Specialised Partitioned Additive Runge-Kutta method for Variational systems.
GeometricIntegrators.Integrators.SPARK.TableauSPARK
— TypeHolds the tableau of an Specialised Partitioned Additive Runge-Kutta method.
GeometricIntegrators.Integrators.SPARK.TableauVPARK
— TypeHolds the tableau of an Variational Partitioned Additive Runge-Kutta method.
GeometricIntegrators.Integrators.SPARK.TableauVSPARKsecondary
— TypeHolds the tableau of an Specialised Partitioned Additive Runge-Kutta method for Variational systems.
Base.show
— MethodPrint additive Runge-Kutta coefficients.
Base.show
— MethodPrint additive Runge-Kutta coefficients.
Base.show
— MethodPrint multiplier Runge-Kutta coefficients.
Base.show
— MethodPrint projective Runge-Kutta coefficients.
Base.show
— MethodPrint SPARK coefficients.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of specialised partitioned additive Runge-Kutta methods for variational systems.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of specialised partitioned additive Runge-Kutta methods for variational systems.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of specialised partitioned additive Runge-Kutta methods for variational systems.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of specialised partitioned additive Runge-Kutta methods for variational systems.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of specialised partitioned additive Runge-Kutta methods for variational systems.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of variational partitioned additive Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of specialised partitioned additive Runge-Kutta methods for variational systems.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of variational partitioned additive Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of specialised partitioned additive Runge-Kutta methods for variational systems.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of specialised partitioned additive Runge-Kutta methods for variational systems.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate an implicit DAE with a specialised partitioned additive Runge-Kutta integrator.
VPRK Methods
GeometricIntegrators.Integrators.VPRK.AbstractParametersVPRK
— TypeParameters for right-hand side function of variational partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.VPRK.IntegratorCacheVPRK
— TypeVariational partitioned Runge-Kutta integrator cache.
Fields
q̃
: initial guess of qp̃
: initial guess of pṽ
: initial guess of vf̃
: initial guess of fs̃
: holds shift due to periodicity of solutionQ
: internal stages of qP
: internal stages of pV
: internal stages of vF
: internal stages of fY
: integral of vector field of internal stages of qZ
: integral of vector field of internal stages of p
GeometricIntegrators.Integrators.VPRK.IntegratorVPRK
— TypeVariational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKdegenerate
— TypeVariational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpInternal
— TypeVariational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpLegendre
— TypeVariational special partitioned additive Runge-Kutta integrator.
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpMidpoint
— TypeVariational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpSecondary
— TypeVariational partitioned Runge-Kutta integrator with projection on secondary constraint.
The VPRK integrator solves the following system of equations for the internal stages,
\[\begin{aligned} Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} \, \big( V_{n,j} + \Lambda_{n,j} \big) , \\ P_{n,i} &= p_{n} + h \sum \limits_{i=1}^{s} \bar{a}_{ij} \, \big( F_{n,j} + \nabla \vartheta (Q_{n,j})^{T} \cdot \Lambda_{n,j} \big) - d_i \lambda , \\ 0 &= \sum \limits_{i=1}^{s} d_i V_i , \\ 0 &= \sum \limits_{j=1}^{s} \omega_{ij} \Psi_{n,j} , \end{aligned}\]
with definitions
\[\begin{aligned} P_{n,i} &= \dfrac{\partial L}{\partial v} (Q_{n,i}, V_{n,i}) , \\ F_{n,i} &= \dfrac{\partial L}{\partial q} (Q_{n,i}, V_{n,i}) , \\ \Psi_{n,i} &= \psi(Q_{n,i}, V_{n,i}, P_{n,i}, F_{n,i}) , \end{aligned}\]
and update rule
\[\begin{aligned} q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} \, \big( V_{n,i} + \Lambda_{n,i} \big) , \\ p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} \bar{b}_{i} \, \big( F_{n,i} + \nabla \vartheta (Q_{n,j})^{T} \cdot \Lambda_{n,j} \big) , \\ 0 &= \phi (q_{n+1}, p_{n+1}) , \end{aligned}\]
satisfying the symplecticity conditions
\[\begin{aligned} b_{i} \bar{a}_{ij} + \bar{b}_{j} a_{ji} &= b_{i} \bar{b}_{j} , & \bar{b}_i &= b_i , \end{aligned}\]
the primary constraint,
\[\begin{aligned} \phi(q,p) = p - \vartheta (q) = 0 , \end{aligned}\]
at the final solution $(q_{n+1}, p_{n+1})$, and super positions of the secondary constraints,
\[\begin{aligned} \psi(q,\dot{q},p,\dot{p}) = \dot{p} - \dot{q} \cdot \nabla \vartheta (q) = \big( \nabla \vartheta (q) - \nabla \vartheta^{T} (q) \big) \cdot \dot{q} - \nabla H (q) = 0, \end{aligned}\]
which, evaluated at the internal stages, read
\[\begin{aligned} \Psi_{n,j} = \big( \nabla \vartheta (Q_{n,j}) - \nabla \vartheta^{T} (Q_{n,j}) \big) \cdot V_{n,j} - \nabla H (Q_{n,j}) . \end{aligned}\]
Here, $\omega$ is a $(s-1) \times s$ matrix, chosen such that the resulting method has optimal order. The vector $d$ is zero for Gauss-Legendre methods and needs to be chosen appropriately for Gauss-Lobatto methods (for details see documentation of VPRK methods).
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpStandard
— TypeVariational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpStandard
— MethodVariational partitioned Runge-Kutta integrator with standard projection.
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpSymmetric
— TypeVariational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpTableau
— TypeProjected Variational Gauss-Legendre Runge-Kutta integrator.
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpVariational
— TypeVariational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.VPRK.ParametersVPRK
— TypeParameters for right-hand side function of Variational Partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.VPRK.ParametersVPRKdegenerate
— TypeParameters for right-hand side function of Variational Partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.VPRK.ParametersVPRKpInternal
— TypeParameters for right-hand side function of Variational Partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.VPRK.ParametersVPRKpLegendre
— TypeParameters for right-hand side function of Variational Partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.VPRK.ParametersVPRKpMidpoint
— TypeParameters for right-hand side function of Variational Partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.VPRK.ParametersVPRKpSecondary
— TypeParameters for right-hand side function of Variational Partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.VPRK.ParametersVPRKpStandard
— TypeParameters for right-hand side function of Variational Partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.VPRK.ParametersVPRKpSymmetric
— TypeParameters for right-hand side function of Variational Partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.VPRK.ParametersVPRKpTableau
— TypeParameters for right-hand side function of projected Gauss-Legendre Runge-Kutta methods.
GeometricIntegrators.Integrators.VPRK.ParametersVPRKpVariational
— TypeParameters for right-hand side function of Variational Partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.VPRK.TableauVPRK
— TypeTableauVPRK
: Tableau of a Variational Partitioned Runge-Kutta method
\[\begin{aligned} P_{n,i} &= \dfrac{\partial L}{\partial v} (Q_{n,i}, V_{n,i}) , & Q_{n,i} &= q_{n} + h \sum \limits_{j=1}^{s} a_{ij} \, V_{n,j} , & q_{n+1} &= q_{n} + h \sum \limits_{i=1}^{s} b_{i} \, V_{n,i} , \\ F_{n,i} &= \dfrac{\partial L}{\partial q} (Q_{n,i}, V_{n,i}) , & P_{n,i} &= p_{n} + h \sum \limits_{i=1}^{s} \bar{a}_{ij} \, F_{n,j} - d_i \lambda , & p_{n+1} &= p_{n} + h \sum \limits_{i=1}^{s} \bar{b}_{i} \, F_{n,i} , \\ && 0 &= \sum \limits_{i=1}^{s} d_i V_i , && \end{aligned}\]
satisfying the symplecticity conditions
\[\begin{aligned} b_{i} \bar{a}_{ij} + b_{j} a_{ji} &= b_{i} b_{j} , & \bar{b}_i &= b_i . \end{aligned}\]
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpSymplectic
— MethodVariational partitioned Runge-Kutta integrator with symplectic projection.
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpVariationalP
— MethodVariational partitioned Runge-Kutta integrator with variational projection on $(p_{n}, q_{n+1})$.
GeometricIntegrators.Integrators.VPRK.IntegratorVPRKpVariationalQ
— MethodVariational partitioned Runge-Kutta integrator with variational projection on $(q_{n}, p_{n+1})$.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of variational partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute solution of degenerate symplectic partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of variational partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of variational special partitioned additive Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of variational partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of variational partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of projected variational partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of variational partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of projected variational partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of projected Gauss-Legendre Runge-Kutta methods.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with variational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with variational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate DAE with variational special partitioned additive Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with variational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with variational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with variational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with variational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with projected Gauss-Legendre Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with variational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with variational partitioned Runge-Kutta integrator.
Galerkin Variational Integrators
GeometricIntegrators.Integrators.IntegratorCGVI
— TypeContinuous Galerkin Variational Integrator.
GeometricIntegrators.Integrators.ParametersCGVI
— TypeParametersCGVI
: Parameters for right-hand side function of continuous Galerkin variational Integrator.
Parameters
Θ
: function of the noncanonical one-form (∂L/∂v)f
: function of the force (∂L/∂q)Δt
: time stepb
: weights of the quadrature rulec
: nodes of the quadrature rulex
: nodes of the basism
: mass matrixa
: derivative matrixr₀
: reconstruction coefficients at the beginning of the intervalr₁
: reconstruction coefficients at the end of the intervalt
: initial timeq
: solution of q at time tp
: solution of p at time t
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of variational partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate ODE with variational partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.IntegratorCacheDGVI
— TypeNonlinear function cache for Discontinuous Galerkin Variational Integrator.
Parameters
ST
: data typeD
: number of dimensionsS
: number of degrees of freedomR
: number of nodes of quadrature formula
Fields
X
: degrees of freedomQ
: solution at quadrature nodesV
: velocity at quadrature nodesP
: one-form at quadrature nodesF
: forces at quadrature nodesq
: current solution of $q_{n}$q⁻
: current solution of $q_{n}^{-}$q⁺
: current solution of $q_{n}^{+}$q̅
: current solution of $q_{n+1}$q̅⁻
: current solution of $q_{n+1}^{-}$q̅⁺
: current solution of $q_{n+1}^{+}$ϕ
: average of the solution at $t_{n}$ϕ̅
: average of the solution at $t_{n+1}$λ
: jump of the solution at $t_{n}$λ̅
: jump of the solution at $t_{n+1}$θ
: one-form evaluated across at $t_{n}$Θ̅
: one-form evaluated across at $t_{n+1}$g
: projection evaluated across at $t_{n}$g̅
: projection evaluated across at $t_{n+1}$
GeometricIntegrators.Integrators.IntegratorDGVI
— TypeIntegratorDGVI
: Discontinuous Galerkin Variational Integrator.
The DGVI integrators arise from the discretization of the action integral
\[\mathcal{A} [q] = \int \limits_{0}^{T} L(q(t), \dot{q}(t)) \, dt ,\]
with $L$ a fully degenerate Lagrangian of the form
\[L(q, \dot{q}) = \vartheta (q) \cdot \dot{q} - H(q) ,\]
where $\vartheta (q)$ denotes the Cartan one-form and $H(q)$ the Hamiltonian, which is usually given by the total energy of the system.
Discretization
Within each interval $(t_{n}, t_{n+1})$ a piecewise-polynomial approximation $q_h$ of the trajectory $q$ is constructed using $S$ basis functions $\varphi_{i}$,
\[q_h(t) \vert_{(t_{n}, t_{n+1})} = \sum \limits_{i=1}^{S} x_{n,i} \, \bar{\varphi}_{n,i} (t) ,\]
where $\bar{\varphi}_{n,i} (t)$ is a rescaled basis function, defined by
\[\bar{\varphi}_{n,i} (t) = \varphi_{i} \bigg( \frac{t - t_{n}}{t_{n+1} - t_{n}} \bigg) ,\]
and it is assumed that $\varphi_{i}$ is compactly supported on $[0,1]$. These approximations $q_h(t)$ are not assumed to be continuous across interval boundaries $t_{n}$ but usuaslly have jumps.
The integral over $(t_{n}, t_{n+1})$ is approximated by a quadrature rule with $R$ nodes $c_i$ and weights $b_i$. Denote by $m$ and $a$ mass and derivative matrices, respectively, whose elements are given by
\[m_{ij} = \varphi_j (c_i) , \qquad a_{ij} = \varphi_j' (c_i) , \qquad i = 1, ..., R , \; j = 1, ..., S .\]
With that, the solution and its time derivative at the quadrature points can be written as
\[Q_{n,i} \equiv q_h(t_n + c_i h) = m_{ij} x_{n,j} , \qquad V_{n,i} \equiv \dot{q}_h (t_n + c_i h) = a_{ij} x_{n,j} ,\]
where
\[x_{n} = ( x_{n,1}, ..., x_{n,S} )^T\]
is the vector containing the degrees of freedom of $q_h \vert_{[t_{n}, t_{n+1}]}$. The limits of $q_h(t)$ at $t_{n}$ and $t_{n+1}$ are given by
\[q_{n}^{+} = \lim \limits_{t \downarrow t_{n}} q_h(t) = \sum \limits_{j=1}^{S} r^{+}_{j} \delta x_{n,j} , \hspace{3em} q_{n+1}^{-} = \lim \limits_{t \uparrow t_{n+1}} q_h(t) = \sum \limits_{j=1}^{S} r^{-}_{j} \delta x_{n,j} .\]
The discrete action reads
\[\mathcal{A}_d [x_d] = h \sum \limits_{n=0}^{N-1} \bigg[ \sum \limits_{i=1}^{R} b_i \big[ \vartheta (Q_{n,i}) \cdot V_{n,i} - H(Q_{n,i}) \big] \\ + \frac{\vartheta (q_n) + \vartheta (q_n^+)}{2} \cdot (q_n^+ - q_n) + \frac{\vartheta (q_{n+1}^-) + \vartheta (q_{n+1})}{2} \cdot (q_{n+1} - q_{n+1}^-) \bigg] ,\]
so that using the relations
\[\delta Q_{n,i} = m_{ij} \delta x_{n,j} , \qquad \delta V_{n,i} = \frac{a_{ij}}{h} \delta x_{n,j} , \qquad \delta q_{n}^- = r^{-}_{j} \delta x_{n-1,j} , \qquad \delta q_{n}^+ = r^{+}_{j} \delta x_{n,j} ,\]
the discrete action principle leads to the discrete equations of motion,
\[0 = \sum \limits_{i=1}^{R} b_i \big[ h m_{ij} \nabla \vartheta (Q_{n,i}) \cdot V_{n,i} + a_{ij} \vartheta (Q_{n,i}) - h m_{ij} \nabla H(Q_{n,i}) \big] \\ + r^{+}_{j} \, \frac{ \vartheta ( q_{n } ) + \vartheta( q_{n }^+ ) }{2} - r^{-}_{j} \, \frac{ \vartheta ( q_{n+1}^- ) + \vartheta( q_{n+1} ) }{2} \big] \\ + h r^{+}_{j} \, \nabla \vartheta (q_{n }^+) \cdot (q_{n }^+ - q_{n } ) + h r^{-}_{j} \, \nabla \vartheta (q_{n+1}^-) \cdot (q_{n+1} - q_{n+1}^-) ,\]
and
\[\vartheta(q_{n}^+) = \vartheta (q_{n}^-) + \nabla \vartheta (q_{n}) \cdot (q_{n}^+ - q_{n}^-) ,\]
for all $n$ and all $j$. Let us introduce the variable $p_n$ as
\[p_{n} = \vartheta (q_{n}^-) + \nabla \vartheta (q_{n}) \cdot (q_{n} - q_{n}^-) ,\]
so that
\[\vartheta(q_{n}^+) = p_{n} + \nabla \vartheta (q_{n}) \cdot (q_{n}^+ - q_{n}) .\]
Then the above equations provide a map $(q_{n}, p_{n}) \mapsto (q_{n+1}, p_{n+1})$. In order to solve these equations, initial conditions $q_{0}$ and $p_{0} = \vartheta(q_{0})$ have to be prescribed.
Fields
equation
: Implicit Ordinary Differential Equationbasis
: piecewise polynomial basisquadrature
: numerical quadrature ruleΔt
: time stepparams
: ParametersDGVIsolver
: nonlinear solveriguess
: initial guessq
: current solution vector for trajectoryp
: current solution vector for one-formcache
: temporary variables for nonlinear solver
GeometricIntegrators.Integrators.ParametersDGVI
— TypeParametersDGVI
: Parameters for right-hand side function of Discontinuous Galerkin Variational Integrator.
Parameters
DT
: data typeTT
: parameter typeD
: dimension of the systemS
: number of basis nodesR
: number of quadrature nodes
Fields
Θ
: function of the noncanonical one-form (∂L/∂v)f
: function of the force (∂L/∂q)g
: function of the projection ∇ϑ(q)⋅vΔt
: time stepb
: quadrature weightsc
: quadrature nodesm
: mass matrixa
: derivative matrixr⁻
: reconstruction coefficients, jump lhs valuer⁺
: reconstruction coefficients, jump rhs valuet
: initial timeq
: solution of q at time tq⁻
: solution of q⁻ at time tq⁺
: solution of q⁺ at time t
GeometricIntegrators.Integrators.compute_stages_p!
— MethodCompute one-form and forces at quadrature nodes.
GeometricIntegrators.Integrators.compute_stages_q!
— MethodCompute solution at quadrature nodes and across jump.
GeometricIntegrators.Integrators.compute_stages_v!
— MethodCompute velocities at quadrature nodes.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of variational partitioned Runge-Kutta methods.
Stochastic Integrators
GeometricIntegrators.Integrators.Stochastic.IntegratorCacheSERK
— TypeStochastic Explicit Runge-Kutta integrator cache.
GeometricIntegrators.Integrators.Stochastic.IntegratorCacheSIPRK
— TypeStructure for holding the internal stages Q, the values of the drift vector and the diffusion matrix evaluated at the internal stages V=v(Q)
, B=B(Q)
, and the increments Y = Δt*a_drift*v(Q) + a_diff*B(Q)*ΔW
.
GeometricIntegrators.Integrators.Stochastic.IntegratorCacheSIRK
— TypeStructure for holding the internal stages Q
, the values of the drift vector and the diffusion matrix evaluated at the internal stages V=v(Q)
, B=B(Q)
, and the increments Y = a_drift*v(Q)*Δt + a_diff*B(Q)*ΔW
.
GeometricIntegrators.Integrators.Stochastic.IntegratorCacheSISPRK
— TypeStructure for holding the internal stages Q
, the values of the drift vector and the diffusion matrix evaluated at the internal stages V=v(Q)
, B=B(Q)
, and the increments Y = Δt*a_drift*v(Q) + a_diff*B(Q)*ΔW
.
GeometricIntegrators.Integrators.Stochastic.IntegratorCacheWERK
— TypeWeak Stochastic Explicit Runge-Kutta integrator cache.
GeometricIntegrators.Integrators.Stochastic.IntegratorCacheWIRK
— TypeStructure for holding the internal stages Q0
, and Q1
the values of the drift vector and the diffusion matrix evaluated at the internal stages V=v(Q0)
, B=B(Q1)
, and the increments Y = Δt*a_drift*v(Q) + a_diff*B(Q)*ΔW
.
GeometricIntegrators.Integrators.Stochastic.IntegratorSERK
— TypeStochastic Explicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.Stochastic.IntegratorSIPRK
— TypeStochastic implicit partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.Stochastic.IntegratorSIRK
— TypeStochastic implicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.Stochastic.IntegratorSISPRK
— TypeStochastic implicit partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.Stochastic.IntegratorWERK
— TypeWeak Stochastic Explicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.Stochastic.IntegratorWIRK
— TypeStochastic implicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.Stochastic.ParametersSERK
— TypeParameters for stochastic explicit Runge-Kutta methods.
GeometricIntegrators.Integrators.Stochastic.ParametersSIPRK
— TypeParameters for right-hand side function of implicit Runge-Kutta methods.
A
- if positive, the upper bound of the Wiener process increments; if A=0.0, no truncation
GeometricIntegrators.Integrators.Stochastic.ParametersSIRK
— TypeParameters for right-hand side function of implicit Runge-Kutta methods.
A
- if positive, the upper bound of the Wiener process increments; if A=0.0
, no truncation
GeometricIntegrators.Integrators.Stochastic.ParametersSISPRK
— TypeParameters for right-hand side function of implicit Runge-Kutta methods.
A
- if positive, the upper bound of the Wiener process increments; if A=0.0
, no truncation.
GeometricIntegrators.Integrators.Stochastic.ParametersWERK
— TypeParameters for weak stochastic explicit Runge-Kutta methods.
GeometricIntegrators.Integrators.Stochastic.ParametersWIRK
— TypeParameters for right-hand side function of weak implicit Runge-Kutta methods.
GeometricIntegrators.Integrators.Stochastic.TableauSERK
— TypeHolds the tableau of a stochastic explicit Runge-Kutta method.
Order of the tableau is not included, because unlike in the deterministic setting, it depends on the properties of the noise (e.g., the dimension of the Wiener process and the commutativity properties of the diffusion matrix).
Orders stored in qdrift
, qdiff
and qdiff2
are understood as the classical orders of these methods.
GeometricIntegrators.Integrators.Stochastic.TableauSIPRK
— TypeHolds the tableau of a stochastic implicit partitioned Runge-Kutta method. qdrift, pdrift hold the RK coefficients for the drift part, and qdiff, pdiff hold the RK coefficients for the diffusion part of the SDE.
Order of the tableau is not included, because unlike in the deterministic setting, it depends on the properties of the noise (e.g., the dimension of the Wiener process and the commutativity properties of the diffusion matrix).
Orders stored in qdrift and qdiff are understood as the classical orders of these methods.
GeometricIntegrators.Integrators.Stochastic.TableauSIRK
— TypeHolds the tableau of a stochastic implicit Runge-Kutta method.
qdrift
holds the RK coefficients for the drift part, qdiff
holds the RK coefficients for the diffusion part of the SDE.
Order of the tableau is not included, because unlike in the deterministic setting, it depends on the properties of the noise (e.g., the dimension of the Wiener process and the commutativity properties of the diffusion matrix).
Orders stored in qdrift
and qdiff
are understood as the classical orders of these methods.
GeometricIntegrators.Integrators.Stochastic.TableauSISPRK
— TypeHolds the tableau of a stochastic implicit split partitioned Runge-Kutta method.
qdrift
, pdrift1
, pdrift2
hold the RK coefficients for the drift parts, and qdiff
, pdiff1
, pdiff2
hold the RK coefficients for the diffusion part of the SDE.
Order of the tableau is not included, because unlike in the deterministic setting, it depends on the properties of the noise (e.g., the dimension of the Wiener process and the commutativity properties of the diffusion matrix)
Orders stored in qdrift
and qdiff
are understood as the classical orders of these methods.
GeometricIntegrators.Integrators.Stochastic.TableauWERK
— TypeHolds the tableau of a weak explicit Runge-Kutta method.
Reference: Andreas Rossler, "Second order Runge-Kutta methods for Stratonovich stochastic differential equations",
BIT Numerical Mathematics (2007) 47, equation (5.1).
Order of the tableau is not included, because unlike in the deterministic setting, it depends on the properties of the noise (e.g., the dimension of the Wiener process and the commutativity properties of the diffusion matrix)
Orders stored in qdrift
, qdiff
are irrelevant and set to 0.
qdrift0
, qdrift1
, qdrift2
correspond to A0, A1, A2 in the paper qdiff0
, qdiff1
, qdiff2
, qdiff3
correspond to B0, B1, B2, B3 qdrift0.b = alpha
qdiff0.b = beta1
qdiff3.b = beta2
qdrift0.c = c0
qdrift1.c = c1
qdrift2.c = c2
GeometricIntegrators.Integrators.Stochastic.TableauWIRK
— TypeHolds the tableau of a weak implicit Runge-Kutta method.
Reference: Wang, Hong, Xu, "Construction of Symplectic Runge-Kutta Methods for Stochastic Hamiltonian Systems",
Commun. Comput. Phys. 21(1), 2017.
Order of the tableau is not included, because unlike in the deterministic setting, it depends on the properties of the noise (e.g., the dimension of the Wiener process and the commutativity properties of the diffusion matrix)
Orders stored in qdrift and qdiff are understood as the classical orders of these methods.
qdrift0
, qdrift1
correspond to A0, A1 in the paper qdiff0
, qdiff1
, qdiff3
correspond to B0, B1, B3 qdrift0.b = alpha
qdiff0.b = beta
qdrift0.c = c0
qdrift1.c = c1
GeometricIntegrators.Integrators.Stochastic.compute_stages!
— MethodUnpacks the data stored in x = (Y[1][1], Y[1][2], ... Y[1][D], Y[2][1], ...)
into Y
, calculates the internal stages Q
, the values of the RHS of the SDE ( v(Q) and B(Q) ), and assigns them to V
and B
. Unlike for FIRK, here Y = a v(Q) Δt + â B(Q) ΔW
GeometricIntegrators.Integrators.Stochastic.compute_stages!
— MethodUnpacks the data stored in x = (Y0[1][1], Y0[1][2]], ... Y0[1][D], ... Y0[S][D], Y1[1][1,1], Y1[1][2,1], ... Y1[1][D,1], Y1[1][1,2], Y1[1][2,2], ... Y1[1][D,2], ... Y1[S][D,M] )
into Y0
and Y1
, calculates the internal stages Q0
and Q1
, the values of the RHS of the SDE ( v(Q0)
and B(Q1)
), and assigns them to V
and B
. Unlike for FIRK, here Y = Δt a v(Q) + â B(Q) ΔW
.
GeometricIntegrators.Integrators.Stochastic.compute_stages!
— MethodUnpacks the data stored in x = (Y[1][1], Y[1][2], ... Y[1][D], Y[2][1], ..., Z[1][1], Z[1][2], ... Z[1][D], Z[2][1], ...)
into Y
, Z
, calculates the internal stages Q
, P
, the values of the RHS of the SDE ( v(Q,P)
, f(Q,P)
, B(Q,P)
and G(Q,P)
), and assigns them to V
, F
, B
and G
. Unlike for FIRK, here Y = Δt a_drift v(Q,P) + a_diff B(Q,P) ΔW
, Z = Δt â_drift v(Q,P) + â_diff B(Q,P) ΔW
.
GeometricIntegrators.Integrators.Stochastic.compute_stages!
— MethodUnpacks the data stored in x = (Y[1][1], Y[1][2], ... Y[1][D], Y[2][1], ..., Z[1][1], Z[1][2], ... Z[1][D], Z[2][1], ...)
into Y
, Z
, calculates the internal stages Q
, P
, the values of the RHS of the SDE ( vi(Q,P)
, fi(Q,P)
, Bi(Q,P)
and Gi(Q,P)
), and assigns them to V[i]
, F[i]
, B[i]
and G[i]
. Unlike for FIRK, here Y = Δt a_drift v(Q,P) + a_diff B(Q,P) ΔW
, Z = Δt â1_drift f1(Q,P) + Δt â2_drift f2(Q,P) + â1_diff G1(Q,P) ΔW + â2_diff G2(Q,P) ΔW
.
GeometricIntegrators.Integrators.Stochastic.initial_guess!
— MethodThis function computes initial guesses for Y
, Z
and assigns them to int.solver.x The prediction is calculated using an explicit integrator.
SIMPLE SOLUTION The simplest initial guess for Y
, Z
is 0: int.solver.x .= zeros(eltype(int), 2*tableau(int).s*ndims(int))
USING AN EXPLICIT INTEGRATOR TO COMPUTE AN INITIAL GUESS Below we use the R2 method of Burrage & Burrage to calculate the internal stages at the times c[1]...c[s]. This approach seems to give very good approximations if the time step and magnitude of noise are not too large. If the noise intensity is too big, one may have to perform a few iterations of the explicit method with a smaller time step, use a higher-order explicit method (e.g. CL or G5), or use the simple solution above.
When calling this function, int.params should contain the data: int.params.q
- the q solution at the previous time step int.params.p
- the p solution at the previous time step int.params.t
- the time of the previous step int.params.ΔW
- the increment of the Brownian motion for the current step
GeometricIntegrators.Integrators.Stochastic.initial_guess!
— MethodThis function computes initial guesses for Y and assigns them to int.solver.x
. The prediction is calculated using an explicit integrator.
SIMPLE SOLUTION The simplest initial guess for Y
is 0: int.solver.x .= zeros(eltype(int), tableau(int).s*ndims(int))
USING AN EXPLICIT INTEGRATOR TO COMPUTE AN INITIAL GUESS Below we use the R2 method of Burrage & Burrage to calculate the internal stages at the times c[1]...c[s]
. This approach seems to give very good approximations if the time step and magnitude of noise are not too large. If the noise intensity is too big, one may have to perform a few iterations of the explicit method with a smaller time step, use a higher-order explicit method (e.g. CL or G5), or use the simple solution above.
When calling this function, int.params should contain the data: int.params.q
- the solution at the previous time step int.params.t
- the time of the previous step int.params.ΔW
- the increment of the Brownian motion for the current step
GeometricIntegrators.Integrators.Stochastic.initial_guess!
— MethodThis function computes initial guesses for Y
, Z
and assigns them to int.solver.x
.
For SISPRK we are NOT IMPLEMENTING an InitialGuess.
SIMPLE SOLUTION The simplest initial guess for Y
, Z
is 0.
GeometricIntegrators.Integrators.Stochastic.initial_guess!
— MethodThis function computes initial guesses for Y
and assigns them to int.solver.x For WIRK we are NOT IMPLEMENTING an InitialGuess.
Using an explicit integrator to predict the next step's value (like in SIRK) does not seem to be a good idea here, because the integrators are convergent in the weak sense only, and there is no guarantee that the explicit integrator will produce anything close to the desired solution...
The simplest initial guess for Y
is 0.
GeometricIntegrators.Integrators.Stochastic.update_solution!
— MethodUpdate solution for stochastic split partitioned Runge-Kutta methods
q
,p
: the solution vector to be updatedV
,F1
,F2
: the matrix containing the drift vectors evaluated at the internal stages v(Qi), fi(Qi)B
,G1
,G2
: the array containing the diffusion matrices evaluated at the internal stages B(Qi), Gi(Qi)bqdrift
,bpdrift1
,bpdrift2
: the Runge-Kutta coefficients for the drift parts of the q and p equationsbqdiff
,bpdiff1
,bpdiff2
: the Runge-Kutta coefficients for the diffusion parts of the q and p equationsΔt
: the time stepΔW
: the increments of the Brownian motion
GeometricIntegrators.Integrators.Stochastic.update_solution!
— MethodUpdate solution for stochastic partitioned Runge-Kutta methods
q
,p
: the solution vector to be updatedV
,F
: the matrix containing the drift vectors evaluated at the internal stages v(Qi), f(Qi)B
,G
: the array containing the diffusion matrices evaluated at the internal stages B(Qi), G(Qi)bqdrift
,bpdrift
: the Runge-Kutta coefficients for the drift parts of the q and p equationsbqdiff
,bpdiff
: the Runge-Kutta coefficients for the diffusion parts of the q and p equationsΔt
: the time stepΔW
: the increments of the Brownian motion
GeometricIntegrators.Integrators.Stochastic.update_solution!
— MethodUpdate solution for weak Runge-Kutta methods WERK
x
: the solution vector to be updatedV
: the matrix containing the drift vector evaluated at the internal stages v(Q_i)B1
: the array containing the diffusion matrix evaluated at the internal stages H^(l)i, such that B1[:,l,i] is evaluated at H^(l)iB2
: the array containing the diffusion matrix evaluated at the internal stages Ĥ^(l)i, such that B2[:,l,i] is evaluated at Ĥ^(l)iα
: the Runge-Kutta coefficients for the drift partβ1
: the Runge-Kutta coefficients for the diffusion term with the random incrementsβ2
: the Runge-Kutta coefficients for the second diffusion termΔt
: the time stepΔW
: the increments of the Brownian motion represented by the random variables Î^(k)
GeometricIntegrators.Integrators.Stochastic.update_solution!
— MethodUpdate solution for stochastic Runge-Kutta methods (SERK)
x
: the solution vector to be updatedV
: the matrix containing the drift vector evaluated at the internal stages v(Q_i)B
: the array containing the diffusion matrix evaluated at the internal stages B(Q_i)bdrift
: the Runge-Kutta coefficients for the drift partbdiff
: the Runge-Kutta coefficients for the ΔW terms of the diffusion partbdiff2
: the Runge-Kutta coefficients for the ΔZ terms of the diffusion partΔt
: the time stepΔW
: the increments of the Brownian motionΔZ
: the integrals of the increments of the Brownian motion
GeometricIntegrators.Integrators.Stochastic.update_solution!
— MethodUpdate solution for stochastic Runge-Kutta methods (SIRK and WIRK)
x
: the solution vector to be updatedV
: the matrix containing the drift vector evaluated at the internal stages v(Qi) (SIRK) or v(Q0i) (WIRK)B
: the array containing the diffusion matrix evaluated at the internal stages B(Qi) (SIRK) or B(Q1^(l)i) (WIRK)bdrift
: the Runge-Kutta coefficients for the drift partbdiff
: the Runge-Kutta coefficients for the diffusion partΔt
: the time stepΔW
: the increments of the Brownian motion (SFIRK) or the increments represented by the random variables Î^(k) (WFIRK)
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of implicit Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of stochastic implicit Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of stochastic implicit split partitioned Runge-Kutta methods.
GeometricIntegrators.Integrators.function_stages!
— MethodCompute stages of weak implicit Runge-Kutta methods.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate SDE with explicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate PSDE with a stochastic implicit partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate SDE with a stochastic implicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate PSDE with a stochastic implicit partitioned Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate SDE with explicit Runge-Kutta integrator.
GeometricIntegrators.Integrators.integrate_step!
— MethodIntegrate SDE with a stochastic implicit Runge-Kutta integrator.