OrdinaryDiffEqSymplecticRK

Symplectic integrators are specialized methods for solving Hamiltonian systems and second-order differential equations that preserve important geometric properties of the phase space. These methods are essential for long-time integration of conservative mechanical systems.

Key Properties

Symplectic integrators provide:

  • Exact conservation of symplectic structure in phase space
  • Bounded energy error over long time periods
  • Excellent long-time stability without secular drift
  • Preservation of periodic orbits and other geometric structures
  • Linear energy drift instead of quadratic (much better than standard methods)

When to Use Symplectic Methods

Symplectic integrators are essential for:

  • Hamiltonian systems and conservative mechanical problems
  • Molecular dynamics and N-body simulations
  • Celestial mechanics and orbital computations
  • Plasma physics and charged particle dynamics
  • Long-time integration where energy conservation is critical
  • Oscillatory problems requiring preservation of periodic structure
  • Classical mechanics problems with known analytical properties

Mathematical Background

For a Hamiltonian system with energy H(p,q), symplectic integrators preserve the symplectic structure dp ∧ dq. While standard integrators have energy error growing quadratically over time, symplectic methods maintain bounded energy with only linear drift, making them superior for long-time integration.

Solver Selection Guide

First-order methods

  • SymplecticEuler: First-order, simplest symplectic method. Only recommended when the dynamics function f is not differentiable.

Second-order methods

  • McAte2: Optimized second-order McLachlan-Atela method, recommended for most applications
  • VelocityVerlet: Second-order, common choice for molecular dynamics but less efficient in terms of accuracy than McAte2
  • VerletLeapfrog: Second-order, kick-drift-kick formulation
  • LeapfrogDriftKickDrift: Alternative second-order leapfrog
  • PseudoVerletLeapfrog: Modified Verlet scheme

Third-order methods

  • Ruth3: Third-order method
  • McAte3: Optimized third-order McLachlan-Atela method

Fourth-order methods

  • CandyRoz4: Fourth-order method
  • McAte4: Fourth-order McLachlan-Atela (requires quadratic kinetic energy)
  • CalvoSanz4: Optimized fourth-order method
  • McAte42: Alternative fourth-order method (BROKEN)

Higher-order methods

  • McAte5: Fifth-order McLachlan-Atela method
  • Yoshida6: Sixth-order method
  • KahanLi6: Optimized sixth-order method
  • McAte8: Eighth-order McLachlan-Atela method
  • KahanLi8: Optimized eighth-order method
  • SofSpa10: Tenth-order method for highest precision

Method Selection Guidelines

  • For most applications: McAte2 (second-order, optimal efficiency)
  • For molecular dynamics (common choice): VelocityVerlet (less efficient than McAte2 but widely used)
  • For non-differentiable dynamics: SymplecticEuler (first-order, only when necessary)
  • For computational efficiency: McAte2 or McAte3

Important Note on Chaotic Systems

Most N-body problems (molecular dynamics, astrophysics) are chaotic systems where solutions diverge onto shadow trajectories. In such cases, higher-order methods provide no practical advantage because the true error remains O(1) for sufficiently long integrations - exactly the scenarios where symplectic methods are most needed. The geometric properties preserved by symplectic integrators are more important than high-order accuracy for chaotic systems.

For more information on chaos and accuracy in numerical integration, see: How Chaotic is Chaos? How Some AI for Science (SciML) Papers are Overstating Accuracy Claims

Installation

To be able to access the solvers in OrdinaryDiffEqSymplecticRK, you must first install them use the Julia package manager:

using Pkg
Pkg.add("OrdinaryDiffEqSymplecticRK")

This will only install the solvers listed at the bottom of this page. If you want to explore other solvers for your problem, you will need to install some of the other libraries listed in the navigation bar on the left.

Example usage

using OrdinaryDiffEqSymplecticRK
function HH_acceleration!(dv, v, u, p, t)
    x, y = u
    dx, dy = dv
    dv[1] = -x - 2x * y
    dv[2] = y^2 - y - x^2
end
initial_positions = [0.0, 0.1]
initial_velocities = [0.5, 0.0]
tspan = (0.0, 1.0)
prob = SecondOrderODEProblem(HH_acceleration!, initial_velocities, initial_positions, tspan)
sol = solve(prob, KahanLi8(), dt = 1 / 10)

Full list of solvers

OrdinaryDiffEqSymplecticRK.VelocityVerletType
VelocityVerlet()

Symplectic Runge-Kutta Methods 2nd order explicit symplectic integrator. Requires f_2(t,u) = v, i.e. a second order ODE.

Keyword Arguments

References

@article{verlet1967computer, title={Computer" experiments" on classical fluids. I. Thermodynamical properties of Lennard-Jones molecules}, author={Verlet, Loup}, journal={Physical review}, volume={159}, number={1}, pages={98}, year={1967}, publisher={APS} }

source
OrdinaryDiffEqSymplecticRK.VerletLeapfrogType
VerletLeapfrog()

Symplectic Runge-Kutta Methods 2nd order explicit symplectic integrator. Kick-drift-kick form. Requires only one evaluation of f1 per step.

Keyword Arguments

References

@article{monaghan2005, title = {Smoothed particle hydrodynamics}, author = {Monaghan, Joseph J.}, year = {2005}, journal = {Reports on Progress in Physics}, volume = {68}, number = {8}, pages = {1703–1759}, doi = {10.1088/0034-4885/68/8/R01}, }

source
OrdinaryDiffEqSymplecticRK.LeapfrogDriftKickDriftType
LeapfrogDriftKickDrift()

Symplectic Runge-Kutta Methods 2nd order explicit symplectic integrator. Drift-kick-drift form of VerletLeapfrog designed to work when f1 depends on v. Requires two evaluation of f1 per step.

Keyword Arguments

References

@article{monaghan2005, title = {Smoothed particle hydrodynamics}, author = {Monaghan, Joseph J.}, year = {2005}, journal = {Reports on Progress in Physics}, volume = {68}, number = {8}, pages = {1703–1759}, doi = {10.1088/0034-4885/68/8/R01}, }

source
OrdinaryDiffEqSymplecticRK.PseudoVerletLeapfrogType
PseudoVerletLeapfrog()

Symplectic Runge-Kutta Methods 2nd order explicit symplectic integrator.

Keyword Arguments

References

@article{verlet1967computer, title={Computer" experiments" on classical fluids. I. Thermodynamical properties of Lennard-Jones molecules}, author={Verlet, Loup}, journal={Physical review}, volume={159}, number={1}, pages={98}, year={1967}, publisher={APS} }

source
OrdinaryDiffEqSymplecticRK.McAte2Type
McAte2()

Symplectic Runge-Kutta Methods Optimized efficiency 2nd order explicit symplectic integrator.

Keyword Arguments

References

@article{mclachlan1992accuracy, title={The accuracy of symplectic integrators}, author={McLachlan, Robert I and Atela, Pau}, journal={Nonlinearity}, volume={5}, number={2}, pages={541}, year={1992}, publisher={IOP Publishing} }

source
OrdinaryDiffEqSymplecticRK.Ruth3Type
Ruth3()

Symplectic Runge-Kutta Methods 3rd order explicit symplectic integrator.

Keyword Arguments

References

@article{ruth1983canonical, title={A canonical integration technique}, author={Ruth, Ronald D}, journal={IEEE Trans. Nucl. Sci.}, volume={30}, number={CERN-LEP-TH-83-14}, pages={2669–2671}, year={1983}}

source
OrdinaryDiffEqSymplecticRK.McAte3Type
McAte3()

Symplectic Runge-Kutta Methods Optimized efficiency 3rd order explicit symplectic integrator.

Keyword Arguments

References

@article{mclachlan1992accuracy, title={The accuracy of symplectic integrators}, author={McLachlan, Robert I and Atela, Pau}, journal={Nonlinearity}, volume={5}, number={2}, pages={541}, year={1992}, publisher={IOP Publishing} }

source
OrdinaryDiffEqSymplecticRK.CandyRoz4Type
CandyRoz4()

Symplectic Runge-Kutta Methods 4th order explicit symplectic integrator.

Keyword Arguments

References

@article{candy1991symplectic, itle={A symplectic integration algorithm for separable Hamiltonian functions}, uthor={Candy, J and Rozmus, W}, ournal={Journal of Computational Physics}, olume={92}, umber={1}, ages={230–256}, ear={1991}, ublisher={Elsevier}}

source
OrdinaryDiffEqSymplecticRK.McAte4Type
McAte4()

Symplectic Runge-Kutta Methods 4th order explicit symplectic integrator. Requires quadratic kinetic energy.

Keyword Arguments

References

@article{mclachlan1992accuracy, title={The accuracy of symplectic integrators}, author={McLachlan, Robert I and Atela, Pau}, journal={Nonlinearity}, volume={5}, number={2}, pages={541}, year={1992}, publisher={IOP Publishing} }

source
OrdinaryDiffEqSymplecticRK.CalvoSanz4Type
CalvoSanz4()

Symplectic Runge-Kutta Methods Optimized efficiency 4th order explicit symplectic integrator.

Keyword Arguments

References

@article{sanz1993symplectic, title={Symplectic numerical methods for Hamiltonian problems}, author={Sanz-Serna, Jes{'u}s Maria and Calvo, Mari-Paz}, journal={International Journal of Modern Physics C}, volume={4}, number={02}, pages={385–392}, year={1993}, publisher={World Scientific} }

source
OrdinaryDiffEqSymplecticRK.McAte42Type
McAte42()

Symplectic Runge-Kutta Methods 4th order explicit symplectic integrator. BROKEN

Keyword Arguments

References

@article{mclachlan1992accuracy, title={The accuracy of symplectic integrators}, author={McLachlan, Robert I and Atela, Pau}, journal={Nonlinearity}, volume={5}, number={2}, pages={541}, year={1992}, publisher={IOP Publishing} }

source
OrdinaryDiffEqSymplecticRK.McAte5Type
McAte5()

Symplectic Runge-Kutta Methods Optimized efficiency 5th order explicit symplectic integrator. Requires quadratic kinetic energy.

Keyword Arguments

References

@article{mclachlan1992accuracy, title={The accuracy of symplectic integrators}, author={McLachlan, Robert I and Atela, Pau}, journal={Nonlinearity}, volume={5}, number={2}, pages={541}, year={1992}, publisher={IOP Publishing} }

source
OrdinaryDiffEqSymplecticRK.Yoshida6Type
Yoshida6()

Symplectic Runge-Kutta Methods 6th order explicit symplectic integrator.

Keyword Arguments

References

@article{yoshida1990construction, title={Construction of higher order symplectic integrators}, author={Yoshida, Haruo}, journal={Physics letters A}, volume={150}, number={5-7}, pages={262–268}, year={1990}, publisher={Elsevier}}

source
OrdinaryDiffEqSymplecticRK.KahanLi6Type
KahanLi6()

Symplectic Runge-Kutta Methods Optimized efficiency 6th order explicit symplectic integrator.

Keyword Arguments

References

@article{yoshida1990construction, title={Construction of higher order symplectic integrators}, author={Yoshida, Haruo}, journal={Physics letters A}, volume={150}, number={5-7}, pages={262–268}, year={1990}, publisher={Elsevier}}

source
OrdinaryDiffEqSymplecticRK.McAte8Type
McAte8()

Symplectic Runge-Kutta Methods 8th order explicit symplectic integrator.

Keyword Arguments

References

@article{mclachlan1995numerical, title={On the numerical integration of ordinary differential equations by symmetric composition methods}, author={McLachlan, Robert I}, journal={SIAM Journal on Scientific Computing}, volume={16}, number={1}, pages={151–168}, year={1995}, publisher={SIAM} }

source
OrdinaryDiffEqSymplecticRK.KahanLi8Type
KahanLi8()

Symplectic Runge-Kutta Methods Optimized efficiency 8th order explicit symplectic integrator.

Keyword Arguments

References

@article{kahan1997composition, title={Composition constants for raising the orders of unconventional schemes for ordinary differential equations}, author={Kahan, William and Li, Ren-Cang}, journal={Mathematics of computation}, volume={66}, number={219}, pages={1089–1099}, year={1997}}

source
OrdinaryDiffEqSymplecticRK.SofSpa10Type
SofSpa10()

Symplectic Runge-Kutta Methods 10th order explicit symplectic integrator.

Keyword Arguments

References

@article{sofroniou2005derivation, title={Derivation of symmetric composition constants for symmetric integrators}, author={Sofroniou, Mark and Spaletta, Giulia}, journal={Optimization Methods and Software}, volume={20}, number={4-5}, pages={597–613}, year={2005}, publisher={Taylor \& Francis}}

source