SPARK
GeometricIntegrators.SPARK.AbstractTableauSPARK
— TypeHolds the tableau of an Specialised Partitioned Additive Runge-Kutta method for Variational systems.
GeometricIntegrators.SPARK.CoefficientsARK
— TypeHolds the coefficients of an additive Runge-Kutta method.
GeometricIntegrators.SPARK.CoefficientsIRK
— TypeHolds the coefficients of an additive Runge-Kutta method.
GeometricIntegrators.SPARK.CoefficientsMRK
— TypeHolds the multiplier Runge-Kutta coefficients.
GeometricIntegrators.SPARK.CoefficientsPRK
— TypeHolds the coefficients of a projective Runge-Kutta method.
GeometricIntegrators.SPARK.CoefficientsSPARK
— TypeHolds the coefficients of a SPARK method.
GeometricIntegrators.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.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.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.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.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.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.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.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.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.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.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.SPARK.SLRK
— TypeHolds all parameters of an Specialised Partitioned Additive Runge-Kutta method for variational systems subject to constraints.
GeometricIntegrators.SPARK.TableauHPARK
— TypeHolds the tableau of a Hamiltonian Partitioned Additive Runge-Kutta methods.
GeometricIntegrators.SPARK.TableauSPARK
— TypeHolds the tableau of an Specialised Partitioned Additive Runge-Kutta method.
GeometricIntegrators.SPARK.TableauVPARK
— TypeHolds the tableau of an Variational Partitioned Additive Runge-Kutta method.
GeometricIntegrators.SPARK.VSPARKsecondary
— 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.SPARK.SPARKGLRK
— MethodSPARK tableau for Gauss-Legendre Runge-Kutta method with s stages.
GeometricIntegrators.SPARK.SPARKGLRKLobattoIIIAIIIB
— FunctionSPARK tableau for Gauss-Legendre/Gauss-Lobatto-IIIA-IIIB methods.
GeometricIntegrators.SPARK.SPARKGLRKLobattoIIIBIIIA
— FunctionSPARK tableau for Gauss-Legendre/Gauss-Lobatto-IIIB-IIIA methods.
GeometricIntegrators.SPARK.SPARKGLVPRK
— MethodTableau for Variational Gauss-Legendre method with s stages.
GeometricIntegrators.SPARK.SPARKLobABC
— MethodTableau for Gauss-Lobatto IIIA-IIIB-IIIC method with s stages.
GeometricIntegrators.SPARK.SPARKLobABD
— MethodTableau for Gauss-Lobatto IIIA-IIIB-IIID method with s stages.
GeometricIntegrators.SPARK.SPARKLobatto
— MethodSPARK tableau for Gauss-Lobatto methods.
GeometricIntegrators.SPARK.SPARKLobattoIIIAIIIB
— MethodSPARK tableau for Gauss-Lobatto IIIA-IIIB method with s stages.
GeometricIntegrators.SPARK.SPARKLobattoIIIBIIIA
— MethodSPARK tableau for Gauss-Lobatto IIIB-IIIA method with s stages.
GeometricIntegrators.SPARK.SPARKVPRK
— MethodSPARK Tableau for Variational Partitioned Runge-Kutta Methods.
GeometricIntegrators.SPARK.TableauGausspSymplectic
— MethodTableau for Gauss-Legendre method with s stages and symplectic projection.
GeometricIntegrators.SPARK.TableauHPARKGLRK
— MethodTableau for Gauss-Legendre HPARK method with s stages.
GeometricIntegrators.SPARK.TableauHPARKLobattoIIIAIIIB
— MethodSPARK tableau for Gauss-Lobatto IIIA-IIIB HPARK method with s stages.
GeometricIntegrators.SPARK.TableauHPARKLobattoIIIBIIIA
— MethodSPARK tableau for Gauss-Lobatto IIIB-IIIA method with s stages.
GeometricIntegrators.SPARK.TableauHSPARKGLRKpSymmetric
— MethodTableau for Gauss-Legendre method with s stages and symplectic projection.
GeometricIntegrators.SPARK.TableauHSPARKLobattoIIIAIIIBpSymmetric
— MethodTableau for Gauss-Lobatto IIIA-IIIB method with s stages and symmetric projection.
GeometricIntegrators.SPARK.TableauHSPARKLobattoIIIBIIIApSymmetric
— MethodTableau for Gauss-Lobatto IIIB-IIIA method with s stages and symmetric projection.
GeometricIntegrators.SPARK.TableauLobattoIIIAIIIBpSymplectic
— MethodTableau for Gauss-Lobatto IIIA-IIIB method with s stages and symplectic projection.
GeometricIntegrators.SPARK.TableauLobattoIIIBIIIApSymplectic
— MethodTableau for Gauss-Lobatto IIIB-IIIA method with s stages and symplectic projection.
GeometricIntegrators.SPARK.TableauVSPARKGLRKpInternal
— MethodTableau for Gauss-Legendre method with s stages and symplectic projection.
GeometricIntegrators.SPARK.TableauVSPARKGLRKpLobattoIIIAIIIB
— MethodTableau for Gauss-Legendre method with s stages and Lobatto-IIIA-IIIB projection.
GeometricIntegrators.SPARK.TableauVSPARKGLRKpLobattoIIIBIIIA
— MethodTableau for Gauss-Legendre method with s stages and Lobatto-IIIB-IIIA projection.
GeometricIntegrators.SPARK.TableauVSPARKGLRKpMidpoint
— MethodTableau for Gauss-Legendre method with s stages and midpoint projection.
GeometricIntegrators.SPARK.TableauVSPARKGLRKpModifiedInternal
— MethodTableau for Gauss-Legendre method with s stages and symplectic projection.
GeometricIntegrators.SPARK.TableauVSPARKGLRKpModifiedLobattoIIIAIIIB
— MethodTableau for Gauss-Legendre method with s stages and symplectic projection.
GeometricIntegrators.SPARK.TableauVSPARKGLRKpModifiedLobattoIIIBIIIA
— MethodTableau for Gauss-Legendre method with s stages and symplectic projection.
GeometricIntegrators.SPARK.TableauVSPARKGLRKpModifiedMidpoint
— MethodTableau for Gauss-Legendre method with s stages and midpoint projection.
GeometricIntegrators.SPARK.TableauVSPARKGLRKpSymmetric
— MethodTableau for Gauss-Legendre method with s stages and symplectic projection.
GeometricIntegrators.SPARK.TableauVSPARKInternalProjection
— MethodConsider a symplectic pair of tableaus $(a^{1}, b^{1}, c^{1})$ and $(a^{3}, b^{3}, c^{3})$, i.e., satsifying $b^{1}_{i} b^{3}_{j} = b^{1}_{i} a^{3}_{ij} + b^{3}_{j} a^{1}_{ji}$, with an arbitrary number of stages $s$. Use the same tableaus for $\tilde{a}^{1}$ and $\tilde{a}^{3}$, so that $\tilde{s} = s$, as well as
\[\begin{aligned} \begin{array}{c|cc} & \tfrac{1}{2} b^{1} \\ & \vdots \\ & \tfrac{1}{2} b^{1} \\ \hline a^{2} & \\ \end{array} && \begin{array}{c|cc} & \tfrac{1}{2} b^{3} \\ & \vdots \\ & \tfrac{1}{2} b^{3} \\ \hline a^{4} & \\ \end{array} && \begin{array}{c|cc} & \tfrac{1}{2} b^{1} \\ c^{1} & \vdots \\ & \tfrac{1}{2} b^{1} \\ \hline \tilde{a}^{2} & \tfrac{1}{2} (1 + R(\infty)) \, b^{1} \\ \end{array} && \begin{array}{c|cc} & \tfrac{1}{2} b^{3} \\ c^{3} & \vdots \\ & \tfrac{1}{2} b^{3} \\ \hline \tilde{a}^{4} & \tfrac{1}{2} (1 + R(\infty)) \, b^{3} \\ \end{array} \end{aligned}\]
Set $\omega = [0, ..., 0, 1]$ and
\[\delta_{ij} = \begin{cases} +1 & j = i , \\ -1 & j = \tilde{s} , \\ 0 & \text{else} , \end{cases}\]
so that $\Lambda_{1} = \Lambda_{2} = ... = \Lambda_{\tilde{s}}$.
This methods is constructed to satisfy the constraint on the projective stages, $\phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) = 0$ for $i = 1, \, ..., \, \tilde{s}$. Note, however, that it violates the symplecticity conditions $b^{1}_{i} b^{4}_{j} = b^{1}_{i} a^{4}_{ij} + b^{4}_{j} \tilde{a}^{1}_{ji}$ and $b^{2}_{i} b^{3}_{j} = b^{2}_{i} \tilde{a}^{3}_{ij} + b^{3}_{j} a^{2}_{ji}$.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIAIIIBProjection
— MethodConsider a symplectic pair of tableaus $(a^{1}, b^{1}, c^{1})$ and $(a^{3}, b^{3}, c^{3})$, i.e., satsifying $b^{1}_{i} b^{3}_{j} = b^{1}_{i} a^{3}_{ij} + b^{3}_{j} a^{1}_{ji}$, with an arbitrary number of stages $s$. For the projection, choose the Lobatto-IIIA and IIIB tableaus with $\tilde{s} = 2$ stages for $(\tilde{a}^{4}, b^{4})$ and $(\tilde{a}^{2}, b^{2})$, respectively, and choose $\tilde{a}^{1}$ and $\tilde{a}^{3}$ such that the projective stages correspond to the initial condition and the solution, i.e.,
\[\begin{aligned} \begin{array}{c|cc} 0 & 0 \\ 1 & b^{1} \\ \hline \tilde{a}^{1} & \\ \end{array} && \begin{array}{c|cc} 0 & \tfrac{1}{2} & 0 \\ 1 & \tfrac{1}{2} & 0 \\ \hline \tilde{a}^{2} & \\ \end{array} && \begin{array}{c|cc} 0 & 0 \\ 1 & b^{3} \\ \hline \tilde{a}^{3} & \\ \end{array} && \begin{array}{c|cc} 0 & 0 & 0 \\ 1 & \tfrac{1}{2} & \tfrac{1}{2} \\ \hline \tilde{a}^{4} & \\ \end{array} \end{aligned}\]
and compute $a^{2}$ and $a^{4}$ by the symplecticity conditions, that is $a^{2}_{ij} = b^{2}_{j} ( b^{3}_{i} - \tilde{a}^{3}_{ji} ) / b^{3}_{i}$ and $a^{4}_{ij} = b^{4}_{j} ( b^{1}_{i} - \tilde{a}^{1}_{ji}) / b^{1}_{i}$. Finally choose $\omega = [0, 0, 1]$ and $\delta = [-1, R_{\infty}]$, implying that $\rho = 1$. By construction, this method satisfies all symplecticity conditions, but the constraint on the projection stages, $\phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) = 0$ for $i = 1, \, ..., \, \tilde{s}$, is not satisfied exactly, but only approximately, although with bounded error.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIAIIIBpLobattoIIIAIIIB
— MethodTableau for Gauss-Lobatto IIIA-IIIB method with s stages and Lobatto-IIIA-IIIB projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIAIIIBpLobattoIIIBIIIA
— MethodTableau for Gauss-Lobatto IIIA-IIIB method with s stages and Lobatto-IIIB-IIIA projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIAIIIBpMidpoint
— MethodTableau for Gauss-Lobatto IIIA-IIIB method with s stages and midpoint projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIAIIIBpModifiedLobattoIIIAIIIB
— MethodTableau for Gauss-Lobatto IIIA-IIIB method with s stages and Lobatto-IIIA-IIIB projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIAIIIBpModifiedLobattoIIIBIIIA
— MethodTableau for Gauss-Lobatto IIIA-IIIB method with s stages and Lobatto-IIIB-IIIA projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIAIIIBpModifiedMidpoint
— MethodTableau for Gauss-Lobatto IIIA-IIIB method with s stages and midpoint projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIAIIIBpSymmetric
— MethodTableau for Gauss-Lobatto IIIA-IIIB method with s stages and symmetric projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIBIIIAProjection
— MethodThis methods is the same as TableauVSPARKLobattoIIIAIIIBProjection
, except for using Lobatto-IIIA and IIIB tableaus with $\tilde{s} = 2$ stages for $(\tilde{a}^{2}, b^{2})$, and $(\tilde{a}^{4}, b^{4})$ respectively, instead of the other way around.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIBIIIApLobattoIIIAIIIB
— MethodTableau for Gauss-Lobatto IIIB-IIIA method with s stages and Lobatto-IIIA-IIIB projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIBIIIApLobattoIIIBIIIA
— MethodTableau for Gauss-Lobatto IIIB-IIIA method with s stages and Lobatto-IIIB-IIIA projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIBIIIApMidpoint
— MethodTableau for Gauss-Lobatto IIIB-IIIA method with s stages and midpoint projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIBIIIApModifiedLobattoIIIAIIIB
— MethodTableau for Gauss-Lobatto IIIB-IIIA method with s stages and Lobatto-IIIA-IIIB projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIBIIIApModifiedLobattoIIIBIIIA
— MethodTableau for Gauss-Lobatto IIIB-IIIA method with s stages and Lobatto-IIIB-IIIA projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIBIIIApModifiedMidpoint
— MethodTableau for Gauss-Lobatto IIIB-IIIA method with s stages and midpoint projection.
GeometricIntegrators.SPARK.TableauVSPARKLobattoIIIBIIIApSymmetric
— MethodTableau for Gauss-Lobatto IIIB-IIIA method with s stages and symmetric projection.
GeometricIntegrators.SPARK.TableauVSPARKMidpointProjection
— MethodConsider a symplectic pair of tableaus $(a^{1}, b^{1}, c^{1})$ and $(a^{3}, b^{3}, c^{3})$, i.e., satsifying $b^{1}_{i} b^{3}_{j} = b^{1}_{i} a^{3}_{ij} + b^{3}_{j} a^{1}_{ji}$, with an arbitrary number of stages $s$. For the projection, choose the tableau with $\tilde{s} = 1$ and $\rho = 0$, such that $\tilde{Q}_{n,1} = \tfrac{1}{2} ( q_{n} + q_{n+1})$, $\tilde{P}_{n,1} = \tfrac{1}{2} ( p_{n} + p_{n+1})$, i.e.,
\[\begin{aligned} \begin{array}{c|cc} & \tfrac{1}{2} \\ & \vdots \\ & \tfrac{1}{2} \\ \hline a^{2} & \\ \end{array} && \begin{array}{c|cc} & \tfrac{1}{2} \\ & \vdots \\ & \tfrac{1}{2} \\ \hline a^{4} & \\ \end{array} && \begin{array}{c|cc} \tfrac{1}{2} & \tfrac{1}{2} \\ \hline \tilde{a}^{2} & \tfrac{1}{2} (1 + R(\infty))\\ \end{array} && \begin{array}{c|cc} \tfrac{1}{2} & \tfrac{1}{2} \\ \hline \tilde{a}^{4} & \tfrac{1}{2} (1 + R(\infty))\\ \end{array} \end{aligned}\]
The coefficients $\tilde{a}^{1}$ and $\tilde{a}^{3}$ are determined by the symplecticity conditions, specifically $a^{4}_{ij} = b^{4}_{j} ( b^{1}_{i} - \tilde{a}^{1}_{ji}) / b^{1}_{i}$ and $a^{2}_{ij} = b^{2}_{j} ( b^{3}_{i} - \tilde{a}^{3}_{ji} ) / b^{3}_{i}$, and $\omega = [0, 1]$.
GeometricIntegrators.SPARK.TableauVSPARKModifiedInternalProjection
— FunctionConsider a symplectic pair of tableaus $(a^{1}, b^{1}, c^{1})$ and $(a^{3}, b^{3}, c^{3})$, i.e., satsifying $b^{1}_{i} b^{3}_{j} = b^{1}_{i} a^{3}_{ij} + b^{3}_{j} a^{1}_{ji}$, with an arbitrary number of stages $s$, and set
\[\begin{aligned} \begin{array}{c|cc} & \tfrac{1}{2} b^{1} \\ & \vdots \\ & \tfrac{1}{2} b^{1} \\ \hline a^{2} & \\ \end{array} && \begin{array}{c|cc} & \tfrac{1}{2} b^{4} \\ & \vdots \\ & \tfrac{1}{2} b^{4} \\ \hline a^{4} & \\ \end{array} && \begin{array}{c|cc} & \tfrac{1}{2} b^{1} \\ c^{1} & \vdots \\ & \tfrac{1}{2} b^{1} \\ \hline \tilde{a}^{2} & \tfrac{1}{2} (1 + R(\infty)) \, b^{1} \\ \end{array} && \begin{array}{c|cc} & \tfrac{1}{2} b^{3} \\ c^{3} & \vdots \\ & \tfrac{1}{2} b^{3} \\ \hline \tilde{a}^{4} & \tfrac{1}{2} (1 + R(\infty)) \, b^{3} \\ \end{array} \end{aligned}\]
Note that by this definition $\tilde{s} = s$. The coefficients $\tilde{a}^{1}$ and $\tilde{a}^{3}$ are determined by the (modified) symplecticity conditions, specifically $a^{4}_{ij} = b^{3}_{j} ( b^{1}_{i} - \tilde{a}^{1}_{ji}) / b^{1}_{i}$ and $a^{2}_{ij} = b^{1}_{j} ( b^{3}_{i} - \tilde{a}^{3}_{ji} ) / b^{3}_{i}$, where $b^{2}$ has been replaced with $b^{1}$ and $b^{4}$ with $b^{3}$, respectively. Set $\omega = [0, ..., 0, 1]$ and
\[\delta_{ij} = \begin{cases} +1 & j = i , \\ -1 & j = \tilde{s} , \\ 0 & \text{else} , \end{cases}\]
so that $\Lambda_{1} = \Lambda_{2} = ... = \Lambda_{\tilde{s}}$.
Note that this method satisfies the symplecticity conditions $b^{1}_{i} b^{4}_{j} = b^{1}_{i} a^{4}_{ij} + b^{4}_{j} \tilde{a}^{1}_{ji}$ and $b^{2}_{i} b^{3}_{j} = b^{2}_{i} \tilde{a}^{3}_{ij} + b^{3}_{j} a^{2}_{ji}$ only if $R(\infty) = 1$ due to the definitions of $b^{2}$ and $b^{4}$. Moreover, it does usually not satisfy the constraint on the projective stages, $\phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) = 0$ for $i = 1, \, ..., \, \tilde{s}$, exactly, but only approximately with bounded error, thus implying a residual in the symplecticity equation even if $R(\infty) = 1$.
GeometricIntegrators.SPARK.TableauVSPARKModifiedLobattoIIIAIIIBProjection
— MethodConsider a symplectic pair of tableaus $(a^{1}, b^{1}, c^{1})$ and $(a^{3}, b^{3}, c^{3})$, i.e., satsifying $b^{1}_{i} b^{3}_{j} = b^{1}_{i} a^{3}_{ij} + b^{3}_{j} a^{1}_{ji}$, with an arbitrary number of stages $s$. For the projection, choose the Lobatto-IIIA and IIIB tableaus with $\tilde{s} = 2$ stages for $(\tilde{a}^{4}, b^{4})$ and $(\tilde{a}^{2}, b^{2})$, respectively.
The coefficients $\tilde{a}^{1}$ and $\tilde{a}^{3}$ are determined by the relations
\[\begin{aligned} \sum \limits_{j=1}^{s} \tilde{a}^{1}_{ij} (c_{j}^{1})^{k-1} &= \frac{(c_{i}^{2})^k}{k} , \qquad & \sum \limits_{j=1}^{s} \tilde{a}^{3}_{ij} (c_{j}^{3})^{k-1} &= \frac{(c_{i}^{4})^k}{k} , \qquad & i &= 1 , \, ... , \, \tilde{s} , \qquad & k &= 1 , \, ... , \, s . \end{aligned}\]
The coefficients $a^{2}$ and $a^{4}$ by the symplecticity conditions, that is $a^{2}_{ij} = b^{2}_{j} ( b^{3}_{i} - \tilde{a}^{3}_{ji} ) / b^{3}_{i}$ and $a^{4}_{ij} = b^{4}_{j} ( b^{1}_{i} - \tilde{a}^{1}_{ji}) / b^{1}_{i}$. Finally choose $\omega = [0, 0, 1]$ and $\delta = [-1, R_{\infty}]$, implying that $\rho = 1$. By construction, this method satisfies all symplecticity conditions, but the constraint on the projection stages, $\phi(\tilde{Q}_{n,i}, \tilde{P}_{n,i}) = 0$ for $i = 1, \, ..., \, \tilde{s}$, is not satisfied exactly, but only approximately, although with bounded error.
GeometricIntegrators.SPARK.TableauVSPARKModifiedLobattoIIIBIIIAProjection
— MethodThis methods is the same as TableauVSPARKModifiedLobattoIIIAIIIBProjection
, except for using Lobatto-IIIA and IIIB tableaus with $\tilde{s} = 2$ stages for $(\tilde{a}^{2}, b^{2})$, and $(\tilde{a}^{4}, b^{4})$ respectively, instead of the other way around.
GeometricIntegrators.SPARK.TableauVSPARKModifiedMidpointProjection
— MethodConsider a symplectic pair of tableaus $(a^{1}, b^{1}, c^{1})$ and $(a^{3}, b^{3}, c^{3})$, i.e., satsifying $b^{1}_{i} b^{3}_{j} = b^{1}_{i} a^{3}_{ij} + b^{3}_{j} a^{1}_{ji}$, with an arbitrary number of stages $s$. For the projection, choose the tableau with $\tilde{s} = 1$ and $\rho = 0$, such that $\tilde{Q}_{n,1} = \tfrac{1}{2} ( q_{n} + q_{n+1})$, $\tilde{P}_{n,1} = \tfrac{1}{2} ( p_{n} + p_{n+1})$, i.e.,
\[\begin{aligned} \begin{array}{c|cc} \tfrac{1}{2} & \tfrac{1}{2} b^{1} \\ \hline \tilde{a}^{1} & \\ \end{array} && \begin{array}{c|cc} \tfrac{1}{2} & \tfrac{1}{2} \\ \hline \tilde{a}^{2} & \tfrac{1}{2} ( 1 + R (\infty) ) \\ \end{array} && \begin{array}{c|cc} \tfrac{1}{2} & \tfrac{1}{2} b^{3} \\ \hline \tilde{a}^{3} & \\ \end{array} && \begin{array}{c|cc} \tfrac{1}{2} & \tfrac{1}{2} \\ \hline \tilde{a}^{4} & \tfrac{1}{2} ( 1 + R (\infty) ) \\ \end{array} \end{aligned}\]
The coefficients $a^{2}$ and $a^{4}$ are determined by the symplecticity conditions, specifically $a^{4}_{ij} = b^{4}_{j} ( b^{1}_{i} - \tilde{a}^{1}_{ji}) / b^{1}_{i}$ and $a^{2}_{ij} = b^{2}_{j} ( b^{3}_{i} - \tilde{a}^{3}_{ji} ) / b^{3}_{i}$, and $\omega = [0, 1]$.
GeometricIntegrators.SPARK.TableauVSPARKSymmetricProjection
— MethodConsider a symplectic pair of tableaus $(a^{1}, b^{1}, c^{1})$ and $(a^{3}, b^{3}, c^{3})$, i.e., satsifying $b^{1}_{i} b^{3}_{j} = b^{1}_{i} a^{3}_{ij} + b^{3}_{j} a^{1}_{ji}$, with an arbitrary number of stages $s$. For the projection, choose the tableau with $\tilde{s} = 2$ and $\rho = 1$, such that $\tilde{Q}_{n,1} = q_{n}$, $\tilde{Q}_{n,2} = q_{n+1}$, $\tilde{P}_{n,1} = p_{n}$, $\tilde{P}_{n,2} = p_{n+1}$, i.e.,
\[\begin{aligned} \begin{array}{c|cc} 0 & 0 \\ 1 & b^{1} \\ \hline \tilde{a}^{1} & \\ \end{array} && \begin{array}{c|cc} & 0 & 0 \\ & \tfrac{1}{2} & \tfrac{1}{2} \\ \hline \tilde{a}^{2} & \tfrac{1}{2} & \tfrac{1}{2} \\ \end{array} && \begin{array}{c|cc} 0 & 0 \\ 1 & b^{3} \\ \hline \tilde{a}^{3} & \\ \end{array} && \begin{array}{c|cc} & 0 & 0 \\ & \tfrac{1}{2} & \tfrac{1}{2} \\ \hline \tilde{a}^{4} & \tfrac{1}{2} & \tfrac{1}{2} \\ \end{array} \end{aligned}\]
The coefficients $a^{2}$ and $a^{4}$ are determined by the symplecticity conditions, specifically $a^{4}_{ij} = b^{4}_{j} ( b^{1}_{i} - \tilde{a}^{1}_{ji}) / b^{1}_{i}$ and $a^{2}_{ij} = b^{2}_{j} ( b^{3}_{i} - \tilde{a}^{3}_{ji} ) / b^{3}_{i}$. Further choose $\omega = [1, 1, 0]$ and $\delta = [-1, R_{\infty}]$, so that $\tilde{\Lambda}_{n,1} = R_{\infty} \tilde{\Lambda}_{n,2}$ and
\[\tilde{P}_{n,1} - \vartheta (\tilde{Q}_{n,1}) + R_{\infty} ( \tilde{P}_{n,2} - \vartheta (\tilde{Q}_{n,2}) ) = 0 .\]
Due to the particular choice of projective stages, this is equivalent to
\[p_{n} - \vartheta (q_{n}) + R_{\infty} ( p_{n+1} - \vartheta (q_{n+1}) ) = 0 ,\]
so that the constraint $\phi(q_{n+1}, p_{n+1}) = 0$ is satisfied if $\phi(q_{n}, p_{n}) = 0$. Note that the choice of $\tilde{a}^{2}$ and $\tilde{a}^{4}$ violates the symplecticity condition $b^{2}_{i} b^{4}_{j} = b^{2}_{i} \tilde{a}^{4}_{ij} + b^{4}_{j} \tilde{a}^{2}_{ji}$.
GeometricIntegrators.SPARK.lobatto_gauss_coefficients
— FunctionThe projective Lobatto-GLRK coefficients are implicitly given by
\[\sum \limits_{j=1}^{s} a_{ij} c_{j}^{k-1} = \frac{\bar{c}_i^k}{k} \qquad i = 1 , \, ... , \, \sigma , \; k = 1 , \, ... , \, s ,\]
where $c$ are Gauß-Legendre nodes with $s$ stages and $\bar{c}$ are Gauß-Lobatto nodes with $\sigma$ stages.