Some tools for Hamiltonian systems
Project description
pyHamSys
pyHamSys is a Python package for scientific computing involving Hamiltonian systems
Installation within a Python virtual environment:
python3 -m pip install pyhamsys
For more information on creating a Python virtual environment, click here.
Symplectic Integrators
pyHamSys includes a class SymplecticIntegrator containing the following symplectic splitting integrators:
Verlet
(order 2, all purpose), also referred to as Strang or Störmer-Verlet splitting- From Forest, Ruth, Physica D 43, 105 (1990):
FR
(order 4, all purpose)
- From Yoshida, Phys. Lett. A 150, 262 (1990):
Yo#
: # should be replaced by an even integer, e.g.,Yo6
for 6th order symplectic integrator (all purpose)Yos6
: (order 6, all purpose) optimized symplectic integrator (solution A from Table 1)
- From McLachlan, SIAM J. Sci. Comp. 16, 151 (1995):
M2
(order 2, all purpose)M4
(order 4, all purpose)
- From Omelyan, Mryglod, Folk, Comput. Phys. Commun. 146, 188 (2002):
EFRL
(order 4) optimized for H = A + BPEFRL
andVEFRL
(order 4) optimized for H = A(p) + B(q). ForPEFRL
, chi should be exp(h XA)exp(h XB). ForVEFRL
, chi should be exp(h XB)exp(h XA).
- From Blanes, Moan, J. Comput. Appl. Math. 142, 313 (2002):
BM4
(order 4, all purpose) refers to S6BM6
(order 6, all purpose) refers to S10RKN4b
(order 4) refers to SRKN6b optimized for H = A(p) + B(q). Here chi should be exp(h XB)exp(h XA).RKN6b
(order 6) refers to SRKN11b optimized for H = A(p) + B(q). Here chi should be exp(h XB)exp(h XA).RKN6a
(order 6) refers to SRKN14a optimized for H = A(p) + B(q). Here chi should be exp(h XA)exp(h XB).
- From Blanes, Casas, Farrés, Laskar, Makazaga, Murua, Appl. Numer. Math. 68, 58 (2013):
ABA104
(order (10,4)) optimized for H = A + ε B. Here chi should be exp(h XA)exp(h XB).ABA864
(order (8,6,4)) optimized for H = A + ε B. Here chi should be exp(h XA)exp(h XB).ABA1064
(order (10,6,4)) optimized for H = A + ε B. Here chi should be exp(h XA)exp(h XB).
All purpose integrators are for any splitting of the Hamiltonian H=∑k Ak in any order of the functions Ak. Otherwise, the order of the operators is specified for each integrator. These integrators are used in the functions solve_ivp_symp
and solve_ivp_sympext
by specifying the entry method
(default is BM4
).
HamSys class
Parameters
ndof
: number of degrees of freedom of the Hamiltonian system 'ndof' should be an integer or half an integer. Half integers denote an explicit time dependence.complex
: bolean
If False, the dynamical variables (q, p) are real and canonically conjugate. If True, the dynamical variables are (ψ, ψ*) where $\psi=(q + i p)/\sqrt{2}$. Default is False.
Attributes
hamiltonian
: callable
A function of (t, y) which returns the Hamiltonian H(t,y) where y is the state vector.y_dot
: callable
A function of (t, y) which returns {y,H(t,y)} where y is the state vector and H is the Hamiltonian. In canonical coordinates (used, e.g., insolve_ivp_sympext
) where y = (q, p), this function returns (∂H/∂p, -∂H/∂q).k_dot
: callable
A function of (t, y) which returns {k,H(t,y)} = -∂H/∂t where k is canonically conjugate to t and H is the Hamiltonian.
Functions
-
compute_vector_field
: from a callable function (Hamiltonian in canonical coordinates) written with symbolic variables (SymPy), computes the vector fields,y_dot
andk_dot
.Determine Hamilton's equations of motion from a given scalar function –the Hamiltonian– H(q, p, t) where q and p are respectively positions and momenta.
Parameters
hamiltonian
: callable Function H(q, p, t) –the Hamiltonian expressed in symbolic variables–, expressed using SymPy functions.output
: bool, optional If True, displays the equations of motion. Default is False.
The function
compute_vector_field
determines the HamSys function attributesy_dot
andk_dot
to be used insolve_ivp_sympext
. The derivatives are computed symbolically using SymPy. -
compute_energy
: callable A function ofsol
–a solution provided bysolve_ivp_sympext
– andmaxerror
, a boolean indicating whether the maximum error in total energy is given (if True) or all the values of the total energy (if False).Parameters
sol
: OdeSolution
Solution provided bysolve_ivp_sympext
.maxerror
: bool, optional
Default is True.
solve_ivp_symp and solve_ivp_sympext
The functions solve_ivp_symp
and solve_ivp_sympext
solve an initial value problem for a Hamiltonian system using an element of the class SymplecticIntegrator, an explicit symplectic splitting scheme (see [1]). These functions numerically integrate a system of ordinary differential equations given an initial value:
dy / dt = {y, H(t, y)}
y(t0) = y0
Here t is a 1-D independent variable (time), y(t) is an N-D vector-valued function (state). A Hamiltonian H(t, y) and a Poisson bracket {. , .} determine the differential equations. The goal is to find y(t) approximately satisfying the differential equations, given an initial value y(t0) = y0.
The function solve_ivp_symp
solves an initial value problem using an explicit symplectic integration. The Hamiltonian flow is defined by two functions chi
and chi_star
of (h, t, y) (see [2]). This function works for any set of coordinates, canonical or non-canonical, provided that the splitting H=∑k Ak leads to facilitated expressions for the operators exp(h Xk) where Xk = {Ak , ·}.
The function solve_ivp_sympext
solves an initial value problem using an explicit symplectic approximation obtained by an extension in phase space (see [3]). This symplectic approximation works for canonical Poisson brackets, and the state vector should be of the form y = (q, p).
Parameters:
chi
(forsolve_ivp_symp
) : callable
Function of (h, t, y) returning exp(h Xn)...exp(h X1) y at time t. If the selected integrator is not all purpose, refer to the list above for the specific ordering of the operators. The operator Xk is the Liouville operator associated with the function Ak, i.e., for Hamiltonian flows Xk = {Ak , ·} where {· , ·} is the Poisson bracket.chi
must return an array of the same shape asy
.chi_star
(forsolve_ivp_symp
) : callable
Function of (h, t, y) returning exp(h X1)...exp(h Xn) y at time t.chi_star
must return an array of the same shape asy
.hs
(forsolve_ivp_sympext
) : element of class HamSys
The attributesy_dot
ofhs
should be defined. Ifcheck_energy
is True. It the Hamiltonian system has an explicit time dependence (i.e., the parameterndof
ofhs
is half an integer), the attributek_dot
ofhs
should be specified.t_span
: 2-member sequence
Interval of integration (t0, tf). The solver starts with t=t0 and integrates until it reaches t=tf. Both t0 and tf must be floats or values interpretable by the float conversion function.y0
: array_like, shape (n,)
Initial state.step
: float
Step size.t_eval
: array_like or None, optional
Times at which to store the computed solution, must be sorted and equally spaced, and lie withint_span
. If None (default), use points selected by the solver.method
: string, optional
Integration methods are listed on pyhamsys.
'BM4' is the default.omega
(forsolve_ivp_sympext
) : float, optional
Coupling parameter in the extended phase space (see [3]). Default = 10.command
: function of (t, y)
Function to be run at each step size (e.g., plotting an observable associated with the state vector y, or register specific events).check_energy
(forsolve_ivp_sympext
) : bool, optional
If True, the attributehamiltonian
ofhs
should be defined. Default is False.
Returns:
Bunch object with the following fields defined:
t
: ndarray, shape (n_points,)
Time points.y
: ndarray, shape (n, n_points)
Values of the solutiony
att
.k
(forsolve_ivp_sympext
) : ndarray, shape (n_points,) Values ofk
att
. Only forsolve_ivp_sympext
and ifcheck_energy
is True for a Hamiltonian system with an explicit time dependence (i.e., the parameterndof
ofhs
is half an integer).err
(forsolve_ivp_sympext
) : float Error in the computation of the total energy. Only forsolve_ivp_sympext
and ifcheck_energy
is True.step
: step size used in the computation.
Remarks:
- Use
solve_ivp_symp
is the Hamiltonian can be split and if each partial operator exp(h Xk) can be easily expressed/computed. Otherwise usesolve_ivp_sympext
if your coordinates are canonical. - If
t_eval
is a linearly spaced list or array, or ift_eval
is None (default), the step size is slightly readjusted so that the output times contain the values int_eval
, or the final time tf corresponds to an integer number of step sizes. The step size used in the computation is recorded in the solution assol.step
. - For integrating multiple trajectories at the same time, extend phase space and define a state vector y = (y1, y2,...yN) where N is the number of trajectories. The Hamiltonian is given by $H(t,\mathbf{y})=\sum_{i=1}^N h(t, y_i)$.
References:
- [1] Hairer, Lubich, Wanner, 2003, Geometric Numerical Integration: Structure-Preserving Algorithms for Ordinary Differential Equations (Springer)
- [2] McLachlan, Tuning symplectic integrators is easy and worthwhile, Commun. Comput. Phys. 31, 987 (2022); arxiv:2104.10269
- [3] Tao, M., Explicit symplectic approximation of nonseparable Hamiltonians: Algorithm and long time performance, Phys. Rev. E 94, 043303 (2016)
Example
import numpy as xp
import sympy as sp
import matplotlib.pyplot as plt
from pyhamsys import HamSys, solve_ivp_sympext
hs = HamSys()
hamiltonian = lambda q, p, t: p**2 / 2 - sp.cos(q)
hs.compute_vector_field(hamiltonian, output=True)
sol = solve_ivp_sympext(hs, (0, 20), xp.asarray([3, 0]), step=1e-1, check_energy=True)
print(f"Error in energy : {sol.err}")
plt.plot(sol.y[0], sol.y[1])
plt.show()
For more information: cristel.chandre@cnrs.fr
Project details
Release history Release notifications | RSS feed
Download files
Download the file for your platform. If you're not sure which to choose, learn more about installing packages.
Source Distribution
Built Distribution
File details
Details for the file pyhamsys-0.43.tar.gz
.
File metadata
- Download URL: pyhamsys-0.43.tar.gz
- Upload date:
- Size: 14.7 kB
- Tags: Source
- Uploaded using Trusted Publishing? No
- Uploaded via: twine/5.1.1 CPython/3.9.15
File hashes
Algorithm | Hash digest | |
---|---|---|
SHA256 | e8f6e018a5312772e27031da1497dd01f61876743d82d64b6e4418e6859168af |
|
MD5 | f5c090f3d85ae6fc45dff955236ce656 |
|
BLAKE2b-256 | aef38e6ef9549d029986d682211430e17cdbb08951182c96571859516cd0cc27 |
File details
Details for the file pyhamsys-0.43-py3-none-any.whl
.
File metadata
- Download URL: pyhamsys-0.43-py3-none-any.whl
- Upload date:
- Size: 14.9 kB
- Tags: Python 3
- Uploaded using Trusted Publishing? No
- Uploaded via: twine/5.1.1 CPython/3.9.15
File hashes
Algorithm | Hash digest | |
---|---|---|
SHA256 | 9ecabb0e1b34404f5fb16f96bf65484735b61da61810496b4c568ef7d61fa94b |
|
MD5 | 96c7a600b8b5a6ee57dac3b6bc3e2c75 |
|
BLAKE2b-256 | 3e110b6391346ba918b5e80056e7703933d8097af51bc60ca17a065d2f668d8e |