Expand description
Symplectic integrators for Hamiltonian systems
This module provides specialized integrators for Hamiltonian systems that preserve important geometric properties such as energy conservation (to within bounded error), phase space volume, and symplectic structure.
Symplectic integrators are particularly useful for:
- Long-time integration of mechanical systems
- Molecular dynamics simulations
- Orbital mechanics
- Any conservative system where energy conservation is critical
§Methods
Several symplectic integration methods are provided:
- Symplectic Euler: First-order, simple method (forward and backward variants)
- Leapfrog/Störmer-Verlet: Second-order, widely used for molecular dynamics
- Symplectic Runge-Kutta: Higher-order methods with improved accuracy
- Composition methods: Methods constructed by composing lower-order integrators
§Basic Usage
use scirs2_integrate::symplectic::potential::HamiltonianSystem;
use scirs2_integrate::symplectic::leapfrog::StormerVerlet;
use scirs2_integrate::symplectic::SymplecticIntegrator;
use scirs2_core::ndarray::array;
// Define a simple harmonic oscillator: H = p²/2 + q²/2
let system = HamiltonianSystem::new(
|_t, _q, p| p.clone(), // dq/dt = p
|_t, q, _p| -q.clone(), // dp/dt = -q
);
// Initial conditions: (q0, p0) = (1.0, 0.0)
let q0 = array![1.0];
let p0 = array![0.0];
let t = 0.0;
let dt = 0.01;
// Create integrator
let integrator = StormerVerlet::<f64>::new();
// Take one step
let (q1, p1) = integrator.step(&system, t, &q0, &p0, dt).unwrap();
// Energy should be conserved (approximately)
let initial_energy = 0.5_f64 * p0.dot(&p0) + 0.5_f64 * q0.dot(&q0);
let final_energy = 0.5_f64 * p1.dot(&p1) + 0.5_f64 * q1.dot(&q1);
// For a single step with dt=0.1, the energy error is acceptable
assert!((initial_energy - final_energy).abs() < 1e-6);Re-exports§
pub use composition::CompositionMethod;pub use euler::symplectic_euler;pub use euler::symplectic_euler_a;pub use euler::symplectic_euler_b;pub use leapfrog::position_verlet;pub use leapfrog::velocity_verlet;pub use leapfrog::StormerVerlet;pub use potential::HamiltonianSystem;pub use potential::SeparableHamiltonian;pub use runge_kutta::GaussLegendre4;pub use runge_kutta::GaussLegendre6;
Modules§
- composition
- Composition methods for symplectic integration
- euler
- Symplectic Euler methods
- leapfrog
- Leapfrog and Störmer-Verlet integrators
- potential
- Specialized symplectic integrators for separable Hamiltonian systems
- runge_
kutta - Symplectic Runge-Kutta methods
Structs§
- Symplectic
Result - Result of symplectic integration containing state history
Traits§
- Hamiltonian
Fn - Trait for Hamiltonian systems
- Symplectic
Integrator - Trait for symplectic integrators
Type Aliases§
- Hamiltonian
FnBox - Type alias for Hamiltonian function