CartesianPhaseSpacePosition

class gala.dynamics.core.CartesianPhaseSpacePosition(pos, vel)[source]

Bases: gala.dynamics.PhaseSpacePosition

Represents phase-space positions in Cartesian coordinates, e.g., positions and conjugate momenta (velocities).

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 axis (0, 1) has special meaning:

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

So if the input position array, pos, has shape pos.shape = (3, 100), this would represent 100 3D positions (pos[0] is x, pos[1] is y, etc.). The same is true for velocity. The position and velocity arrays must have the same shape.

If the input position and velocity are arrays rather than Quantity objects, they are internally stored with dimensionles units.

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.

Methods Summary

angular_momentum()
energy(potential) The total energy per unit mass (e.g., kinetic + potential):
from_w(w[, units])
kinetic_energy() The kinetic energy per unit mass:
plot(**kwargs) Plot the positions in all projections.
potential_energy(potential) The potential energy per unit mass:
represent_as(Representation) Represent the position and velocity of the orbit in an alternate coordinate system.
to_frame(frame[, galactocentric_frame, ...]) Transform the orbit from Galactocentric, cartesian coordinates to Heliocentric coordinates in the specified Astropy coordinate frame.
w([units]) This returns a single array containing the phase-space positions.

Methods Documentation

angular_momentum()[source]
energy(potential)[source]

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

Parameters:

potential : gala.potential.PotentialBase

The potential object to compute the energy from.

Returns:

E : Quantity

The total energy.

classmethod from_w(w, units=None, **kwargs)[source]
kinetic_energy()[source]

The kinetic energy per unit mass:

\[E_K = \frac{1}{2} \, |\boldsymbol{v}|^2\]
Returns:

E : Quantity

The kinetic energy.

plot(**kwargs)[source]

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

Warning

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

Parameters:

relative_to : bool (optional)

Plot the values relative to this value or values.

autolim : bool (optional)

Automatically set the plot limits to be something sensible.

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 scatter(). You can pass in any of the usual style kwargs like color=..., marker=..., etc.

Returns:

fig : Figure

potential_energy(potential)[source]

The potential energy per unit mass:

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

potential : gala.potential.PotentialBase

The potential object to compute the energy from.

Returns:

E : Quantity

The potential energy.

represent_as(Representation)[source]

Represent the position and velocity of the orbit in an alternate coordinate system. Supports any of the Astropy coordinates representation classes.

Parameters:

Representation : BaseRepresentation

The class for the desired representation.

Returns:

pos : BaseRepresentation

vel : Quantity

The velocity in the new representation. All components have units of velocity – e.g., to get an angular velocity in spherical representations, you’ll need to divide by the radius.

to_frame(frame, galactocentric_frame=<Galactocentric Frame (galcen_distance=8.3 kpc, galcen_ra=266d24m18.36s, galcen_dec=-28d56m10.23s, z_sun=27.0 pc, roll=0.0 deg)>, vcirc=None, vlsr=None)[source]

Transform the orbit from Galactocentric, cartesian coordinates to Heliocentric coordinates in the specified Astropy coordinate frame.

Parameters:

frame : BaseCoordinateFrame

galactocentric_frame : Galactocentric

vcirc : Quantity

Circular velocity of the Sun. Passed to velocity transformation.

vlsr : Quantity

Velocity of the Sun relative to the LSR. Passed to velocity transformation.

Returns:

c : BaseCoordinateFrame

An instantiated coordinate frame.

v : tuple

A tuple of velocities represented as Quantity objects.

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,...).