Dynamics (gala.dynamics)

Introduction

This subpackage contains functions and classes useful for gravitational dynamics. The fundamental objects used by many of the functions and utilities in this and other subpackages are the PhaseSpacePosition and Orbit subclasses.

There are utilities for transforming orbits in phase-space to action-angle coordinates, tools for visualizing and computing dynamical quantities from orbits, tools to generate mock stellar streams, and tools useful for nonlinear dynamics such as Lyapunov exponent estimation.

For code blocks below and any pages linked below, I assume the following imports have already been excuted:

>>> import astropy.units as u
>>> import numpy as np
>>> import gala.potential as gp
>>> import gala.dynamics as gd
>>> from gala.units import galactic

Getting started: Working with orbits

Some simple tools are provided for inspecting and plotting orbits. For example, we’ll start by integrating an orbit in Cartesian coordinates using the gala.potential and gala.integrate subpackages:

>>> pot = gp.MiyamotoNagaiPotential(m=2.5E11*u.Msun, a=6.5*u.kpc,
...                                 b=0.26*u.kpc, units=galactic)
>>> w0 = gd.CartesianPhaseSpacePosition(pos=[11., 0., 0.2]*u.kpc,
...                                     vel=[0., 200, 100]*u.km/u.s)
>>> orbit = pot.integrate_orbit(w0, dt=1., n_steps=1000)

This will integrate an orbit from the specified initial conditions (w0) and return an orbit object. There are many useful methods of the Orbit subclasses and many functions that accept Orbit objects. For example, we can easily visualize the orbit by plotting the time series in all projections using the plot() method:

>>> fig = orbit.plot()

(Source code, png, hires.png, pdf)

../_images/index-1.png

From this object, we can easily compute dynamical quantities such as the energy or angular momentum (I take the 0th element because these functions return the quantities computed at every timestep):

>>> orbit.energy()[0] 
<Quantity -0.06074019848886105 kpc2 / Myr2>

Let’s see how well the integrator conserves energy and the z component of angular momentum:

>>> E = orbit.energy()
>>> Lz = orbit.angular_momentum()[2]
>>> np.std(E), np.std(Lz)
(<Quantity 4.654233175716351e-06 kpc2 / Myr2>,
 <Quantity 9.675900603446092e-16 kpc2 / Myr>)

Using gala.dynamics

More details are provided in the linked pages below:

gala.dynamics Package

Functions

check_angle_sampling(nvecs, angles) Returns a list of the index of elements of n which do not have adequate toy angle coverage.
estimate_dt_n_steps(w0, potential, ...[, ...]) Estimate the timestep and number of steps to integrate an orbit for given its initial conditions and a potential object.
fast_lyapunov_max(w0, potential, dt, n_steps) Compute the maximum Lyapunov exponent using a C-implemented estimator that uses the DOPRI853 integrator.
find_actions(orbit, N_max[, ...]) Find approximate actions and angles for samples of a phase-space orbit.
fit_harmonic_oscillator(orbit[, omega0]) Fit the toy harmonic oscillator potential to the sum of the energy
fit_isochrone(orbit[, m0, b0]) Fit the toy Isochrone potential to the sum of the energy residuals relative
fit_toy_potential(orbit[, ...]) Fit a best fitting toy potential to the orbit provided.
generate_n_vectors(N_max[, dx, dy, dz, ...]) Generate integer vectors, \(\boldsymbol{n}\), with \(|\boldsymbol{n}| < N_{\rm max}\).
harmonic_oscillator_to_aa(w, potential) Transform the input cartesian position and velocity to action-angle coordinates for the Harmonic Oscillator potential.
isochrone_to_aa(w, potential) Transform the input cartesian position and velocity to action-angle coordinates in the Isochrone potential.
lyapunov_max(w0, integrator, dt, n_steps[, ...]) Compute the maximum Lyapunov exponent of an orbit by integrating many nearby orbits (noffset) separated with isotropically distributed directions but the same initial deviation length, d0.
peak_to_peak_period(t, f[, amplitude_threshold]) Estimate the period of the input time series by measuring the average peak-to-peak time.
plot_orbits(x[, t, ix, axes, triangle, ...]) Given time series of positions, x, make nice plots of the orbit in cartesian projections.
surface_of_section(orbit, plane_ix[, ...]) Generate and return a surface of section from the given orbit.
three_panel(q[, relative_to, autolim, axes, ...]) Given 3D quantities, q, make a nice three-panel or triangle plot of projections of the values.

Classes

CartesianOrbit(pos, vel[, t, potential]) Represents an orbit in Cartesian coordinates – positions and velocities (conjugate momenta) at different times.
CartesianPhaseSpacePosition(pos, vel) Represents phase-space positions in Cartesian coordinates, e.g., positions and conjugate momenta (velocities).
Orbit
PhaseSpacePosition

Class Inheritance Diagram

Inheritance diagram of gala.dynamics.orbit.CartesianOrbit, gala.dynamics.core.CartesianPhaseSpacePosition, gala.dynamics.orbit.Orbit, gala.dynamics.core.PhaseSpacePosition