galpy.orbit.Orbit.integrate_dxdv¶
Supported for planar (phase-space dimension 4) and full 3D (phase-space
dimension 6) Orbit instances.
- Orbit.integrate_dxdv(dxdv, t, pot, method='dopr54_c', progressbar=True, dt=None, numcores=2, force_map=False, rectIn=False, rectOut=False, rtol=None, atol=None)[source]¶
Integrate the orbit and a small area of phase space.
- Parameters:
dxdv (numpy.ndarray) – Initial conditions for the phase-space deviation in cylindrical or rectangular coordinates. The shape of the array should be (*input_shape, 4) for planar (4D) orbits and (*input_shape, 6) for 3D (6D) orbits.
t (list, numpy.ndarray or Quantity) – List of equispaced times at which to compute the orbit. The initial condition is t[0]. (note that for method=’odeint’, method=’dop853’, and method=’dop853_c’, the time array can be non-equispaced).
pot (Potential, DissipativeForce, or a combined force/potential formed using addition (pot1+pot2+force3+…)) – Gravitational field to integrate the orbit in.
method (str, optional) – Integration method. Default is ‘dopr54_c’. See Notes for more information.
progressbar (bool, optional) – If True, display a tqdm progress bar when integrating multiple orbits (requires tqdm to be installed!). Default is True.
dt (float, optional) – If set, force the integrator to use this basic stepsize; must be an integer divisor of output stepsize (only works for the C integrators that use a fixed stepsize) (can be Quantity).
numcores (int, optional) – Number of cores to use for Python-based multiprocessing (pure Python or using force_map=True); default = OMP_NUM_THREADS.
force_map (bool, optional) – If True, force use of Python-based multiprocessing (not recommended). Default is False.
rectIn (bool, optional) – If True, input dxdv is in rectangular coordinates. Default is False.
rectOut (bool, optional) – If True, output dxdv (that in orbit_dxdv) is in rectangular coordinates. Default is False.
rtol (float, optional) – Relative tolerance. Default is None.
atol (float, optional) – Absolute tolerance. Default is None.
- Returns:
Get the actual orbit using getOrbit_dxdv(), the orbit that is integrated alongside with dxdv is stored as usual, any previous regular orbit integration will be erased!
- Return type:
None
Notes
Possible integration methods are the non-symplectic ones in galpy:
‘odeint’ for scipy’s odeint
‘rk4_c’ for a 4th-order Runge-Kutta integrator in C
‘rk6_c’ for a 6-th order Runge-Kutta integrator in C
‘dopr54_c’ for a 5-4 Dormand-Prince integrator in C
‘dop853’ for a 8-5-3 Dormand-Prince integrator in Python
‘dop853_c’ for a 8-5-3 Dormand-Prince integrator in C
For 3D (6D) orbits, dissipative (velocity-dependent) forces are supported by the C-based methods for forces with a C implementation of the velocity-dependent force Jacobian (dF/dx, dF/dv), advertised by
hasC_dxdv3d=True. The phase-space-volume evolution follows det M(t) = exp(int tr(dF/dv) dt’): < 1 for friction, but exactly 1 for non-inertial frames (NonInertialFrameForce’s dF/dv = -2 [Omega]_x is antisymmetric, so a rotating frame preserves phase-space volume). The pure-Python methods (‘odeint’, ‘dop853’) raise aNotImplementedErrorfor dissipative forces.Rotated
EllipsoidalPotentialinstances (non-trivialzvec/pa) do not support C-based variational integration directly (hasC_dxdv3d=False); wrap the aligned potential in aRotateAndTiltWrapperPotentialinstead (identical physics, full 3D dxdv support).2011-10-17 - Written - Bovy (IAS)
2014-06-29 - Added rectIn and rectOut - Bovy (IAS)
2019-05-21 - Parallelized and incorporated into new Orbits class - Bovy (UofT)
2026-06-09 - Dissipative-force support in C (3D) - Bovy (UofT)