CartesianOrbit

class gala.dynamics.CartesianOrbit(pos, vel, t=None, potential=None)[source]

Bases: gala.dynamics.core.CartesianPhaseSpacePosition, gala.dynamics.Orbit

Represents an orbit in Cartesian coordinates – positions and velocities (conjugate momenta) at different times.

Warning

This is an experimental class. The API may change in a future release!

The position and velocity quantities (arrays) can have an arbitrary number of dimensions, but the first two axes (0, 1) have special meaning:

- `axis=0` is the coordinate dimension (e.g., x, y, z)
- `axis=1` is the time dimension

So if the input position array, pos, has shape pos.shape = (3, 100), this would be a 3D orbit (pos[0] is x, pos[1] is y, etc.). For representing multiple orbits, the position array could have 3 axes, e.g., it might have shape pos.shape = (3, 100, 8), where this is interpreted as a 3D position at 100 times for 8 different orbits. The same is true for velocity. The position and velocity arrays must have the same shape.

If a time argument is specified, the position and velocity arrays must have the same number of timesteps as the length of the time object:

len(t) == pos.shape[1]
Parameters:

pos : Quantity, array_like

Positions. If a numpy array (e.g., has no units), this will be stored as a dimensionless Quantity. See the note above about the assumed meaning of the axes of this object.

vel : Quantity, array_like

Velocities. If a numpy array (e.g., has no units), this will be stored as a dimensionless Quantity. See the note above about the assumed meaning of the axes of this object.

t : array_like, Quantity (optional)

Array of times. If a numpy array (e.g., has no units), this will be stored as a dimensionless Quantity.

potential : PotentialBase (optional)

The potential that the orbit was integrated in.

Attributes Summary

r The orbital radius.

Methods Summary

align_circulation_with_z([circulation]) If the input orbit is a tube orbit, this function aligns the circulation axis with the z axis and returns a copy.
circulation() Determine which axes the Orbit circulates around by checking whether there is a change of sign of the angular momentum about an axis.
energy([potential]) The total energy per unit mass (e.g., kinetic + potential):
plot(**kwargs) Plot the orbit in all projections.
potential_energy([potential]) The potential energy per unit mass:
w([units]) This returns a single array containing the phase-space positions.

Attributes Documentation

r

The orbital radius.

Methods Documentation

align_circulation_with_z(circulation=None)[source]

If the input orbit is a tube orbit, this function aligns the circulation axis with the z axis and returns a copy.

Parameters:

circulation : array_like (optional)

Array of bits that specify the axis about which the orbit circulates. If not provided, will compute this using circulation(). See that method for more information.

Returns:

orb : CartesianOrbit

A copy of the original orbit object with circulation aligned with the z axis.

circulation()[source]

Determine which axes the Orbit circulates around by checking whether there is a change of sign of the angular momentum about an axis. Returns a 2D array with ndim integers per orbit point. If a box orbit, all integers will be 0. A 1 indicates circulation about the corresponding axis.

TODO: clockwise / counterclockwise?

For example, for a single 3D orbit:

  • Box and boxlet = [0,0,0]
  • z-axis (short-axis) tube = [0,0,1]
  • x-axis (long-axis) tube = [1,0,0]
Returns:

circulation : numpy.ndarray

An array that specifies whether there is circulation about any of the axes of the input orbit. For a single orbit, will return a 1D array, but for multiple orbits, the shape will be (3, norbits).

energy(potential=None)[source]

The total energy per unit mass (e.g., kinetic + potential):

Returns:

E : Quantity

The total energy.

plot(**kwargs)[source]

Plot the orbit in all projections. This is a thin wrapper around plot_orbits – the docstring for this function is included here.

Warning

This will currently fail for orbits with fewer than 3 dimensions.

Parameters:

ix : int, array_like (optional)

Index or array of indices of orbits to plot. For example, if x is an array of shape (3,1024,32) - 1024 timesteps for 32 orbits in 3D positions – ix would specify which of the 32 orbits to plot.

axes : array_like (optional)

Array of matplotlib Axes objects.

triangle : bool (optional)

Make a triangle plot instead of plotting all projections in a single row.

subplots_kwargs : dict (optional)

Dictionary of kwargs passed to subplots().

labels : iterable (optional)

List or iterable of axis labels as strings. They should correspond to the dimensions of the input orbit.

**kwargs

All other keyword arguments are passed to plot(). You can pass in any of the usual style kwargs like color=..., marker=..., etc.

Returns:

fig : Figure

potential_energy(potential=None)[source]

The potential energy per unit mass:

\[E_\Phi = \Phi(\boldsymbol{q})\]
Returns:

E : Quantity

The potential energy.

w(units=None)[source]

This returns a single array containing the phase-space positions.

Parameters:

units : UnitSystem (optional)

The unit system to represent the position and velocity in before combining into the full array.

Returns:

w : ndarray

A numpy array of all positions and velocities, without units. Will have shape (2*ndim,...).