PhaseSpacePosition#

class gala.dynamics.PhaseSpacePosition(pos, vel=None, frame=None)[source]#

Bases: object

Represents phase-space positions, i.e. positions and conjugate momenta (velocities).

The class can be instantiated with Astropy representation objects (e.g., CartesianRepresentation), Astropy Quantity objects, or plain Numpy arrays.

If passing in representation objects, the default representation is taken to be the class that is passed in.

If passing in Quantity or Numpy array instances for both position and velocity, they are assumed to be Cartesian. Array inputs are interpreted as dimensionless quantities. The input position and velocity objects can have an arbitrary number of (broadcastable) dimensions. For Quantity or array inputs, the first axis (0) has special meaning:

  • axis=0 is the coordinate dimension (e.g., x, y, z for Cartesian)

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.

Parameters:
posrepresentation, quantity_like, or 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.

veldifferential, quantity_like, or 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.

frameFrameBase (optional)

The reference frame of the input phase-space positions.

Attributes Summary

data

pos_components

representation_mappings

shape

This is not the shape of the position or velocity arrays.

vel_components

Methods Summary

angular_momentum()

Compute the angular momentum for the phase-space positions contained in this object.

energy(hamiltonian)

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

from_hdf5(f)

Load an object from an HDF5 file.

from_w(w[, units])

get_components(which)

Get the component name dictionary for the desired object.

guiding_radius(potential[, t])

Compute the guiding-center radius

kinetic_energy()

The kinetic energy per unit mass:

plot([components, units, auto_aspect])

Plot the positions in all projections.

potential_energy(potential)

The potential energy per unit mass:

represent_as(new_pos[, new_vel])

Represent the position and velocity of the orbit in an alternate coordinate system.

reshape(new_shape)

Reshape the underlying position and velocity arrays.

to_coord_frame(frame[, galactocentric_frame])

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

to_frame(frame[, current_frame])

Transform to a new reference frame.

to_hdf5(f)

Serialize this object to an HDF5 file.

w([units])

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

Attributes Documentation

data#
pos_components#
representation_mappings = {<class 'astropy.coordinates.representation.cartesian.CartesianDifferential'>: [('d_x', 'v_x', 'recommended'), ('d_y', 'v_y', 'recommended'), ('d_z', 'v_z', 'recommended'), ('d_xyz', 'v_xyz', 'recommended')], <class 'astropy.coordinates.representation.cartesian.CartesianRepresentation'>: [('xyz', 'xyz', 'recommended')], <class 'astropy.coordinates.representation.cylindrical.CylindricalDifferential'>: [('d_rho', 'v_rho', 'recommended'), ('d_phi', 'pm_phi', 'recommended'), ('d_z', 'v_z', 'recommended')], <class 'astropy.coordinates.representation.spherical.PhysicsSphericalDifferential'>: [('d_phi', 'pm_phi', Unit("mas / yr")), ('d_theta', 'pm_theta', Unit("mas / yr")), ('d_r', 'radial_velocity', 'recommended')], <class 'astropy.coordinates.representation.spherical.SphericalCosLatDifferential'>: [('d_lon_coslat', 'pm_lon_coslat', Unit("mas / yr")), ('d_lat', 'pm_lat', Unit("mas / yr")), ('d_distance', 'radial_velocity', 'recommended')], <class 'astropy.coordinates.representation.spherical.SphericalDifferential'>: [('d_lon', 'pm_lon', Unit("mas / yr")), ('d_lat', 'pm_lat', Unit("mas / yr")), ('d_distance', 'radial_velocity', 'recommended')], <class 'astropy.coordinates.representation.spherical.UnitSphericalCosLatDifferential'>: [('d_lon_coslat', 'pm_lon_coslat', Unit("mas / yr")), ('d_lat', 'pm_lat', Unit("mas / yr")), ('d_distance', 'radial_velocity', 'recommended')], <class 'astropy.coordinates.representation.spherical.UnitSphericalDifferential'>: [('d_lon', 'pm_lon', Unit("mas / yr")), ('d_lat', 'pm_lat', Unit("mas / yr")), ('d_distance', 'radial_velocity', 'recommended')], <class 'gala.dynamics.representation_nd.NDCartesianDifferential'>: [('d_xyz', 'v_xyz', 'recommended'), ('d_x([0-9])', 'v_x{0}', 'recommended')], <class 'gala.dynamics.representation_nd.NDCartesianRepresentation'>: [('xyz', 'xyz', 'recommended')]}#
shape#

This is not the shape of the position or velocity arrays. That is accessed by doing, e.g., obj.x.shape.

vel_components#

Methods Documentation

angular_momentum()[source]#

Compute the angular momentum for the phase-space positions contained in this object:

.. math::

boldsymbol{{L}} = boldsymbol{{q}} times boldsymbol{{p}}

See Array shapes for more information about the shapes of input and output objects.

Returns:
LQuantity

Array of angular momentum vectors.

Examples

>>> import numpy as np
>>> import astropy.units as u
>>> pos = np.array([1., 0, 0]) * u.au
>>> vel = np.array([0, 2*np.pi, 0]) * u.au/u.yr
>>> w = PhaseSpacePosition(pos, vel)
>>> w.angular_momentum() 
<Quantity [0.        ,0.        ,6.28318531] AU2 / yr>
energy(hamiltonian)[source]#

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

Parameters:
hamiltoniangala.potential.Hamiltonian, gala.potential.PotentialBase instance

The Hamiltonian object to evaluate the energy. If a potential is passed in, this assumes a static reference frame.

Returns:
EQuantity

The total energy.

classmethod from_hdf5(f)[source]#

Load an object from an HDF5 file.

Requires h5py.

Parameters:
fstr, h5py.File

Either the filename or an open HDF5 file.

classmethod from_w(w, units=None, **kwargs)[source]#
get_components(which)[source]#

Get the component name dictionary for the desired object.

The returned dictionary maps component names on this class to component names on the desired object.

Parameters:
whichstr

Can either be 'pos' or 'vel' to get the components for the position or velocity object.

guiding_radius(potential, t=0.0, **root_kwargs)[source]#

Compute the guiding-center radius

Parameters:
potentialgala.potential.PotentialBase subclass instance

The potential to compute the guiding radius in.

tquantity-like (optional)

Time.

**root_kwargs

Any additional keyword arguments are passed to root.

Returns:
RgQuantity

Guiding-center radius.

kinetic_energy()[source]#

The kinetic energy per unit mass:

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

The kinetic energy.

plot(components=None, units=None, auto_aspect=True, **kwargs)[source]#

Plot the positions in all projections. This is a wrapper around plot_projections for fast access and quick visualization. All extra keyword arguments are passed to that function (the docstring for this function is included here for convenience).

Parameters:
componentsiterable (optional)

A list of component names (strings) to plot. By default, this is the Cartesian positions ['x', 'y', 'z']. To plot Cartesian velocities, pass in the velocity component names ['d_x', 'd_y', 'd_z'].

unitsUnitBase, iterable (optional)

A single unit or list of units to display the components in.

auto_aspectbool (optional)

Automatically enforce an equal aspect ratio.

relative_tobool (optional)

Plot the values relative to this value or values.

autolimbool (optional)

Automatically set the plot limits to be something sensible.

axesarray_like (optional)

Array of matplotlib Axes objects.

subplots_kwargsdict (optional)

Dictionary of kwargs passed to subplots().

labelsiterable (optional)

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

plot_functioncallable() (optional)

The matplotlib plot function to use. By default, this is scatter(), but can also be, e.g., plot().

**kwargs

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

Returns:
figFigure
potential_energy(potential)[source]#

The potential energy per unit mass:

\[E_\Phi = \Phi(\boldsymbol{q})\]
Parameters:
potentialgala.potential.PotentialBase

The potential object to compute the energy from.

Returns:
EQuantity

The potential energy.

represent_as(new_pos, new_vel=None)[source]#

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

Parameters:
new_posBaseRepresentation

The type of representation to generate. Must be a class (not an instance), or the string name of the representation class.

new_velBaseDifferential (optional)

Class in which any velocities should be represented. Must be a class (not an instance), or the string name of the differential class. If None, uses the default differential for the new position class.

Returns:
new_pspgala.dynamics.PhaseSpacePosition
reshape(new_shape)[source]#

Reshape the underlying position and velocity arrays.

to_coord_frame(frame, galactocentric_frame=None, **kwargs)[source]#

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

Parameters:
frameBaseCoordinateFrame

The frame instance specifying the desired output frame. For example, ICRS.

galactocentric_frameGalactocentric

This is the assumed frame that the position and velocity of this object are in. The Galactocentric instand should have parameters specifying the position and motion of the sun in the Galactocentric frame, but no data.

Returns:
cBaseCoordinateFrame

An instantiated coordinate frame containing the positions and velocities from this object transformed to the specified coordinate frame.

to_frame(frame, current_frame=None, **kwargs)[source]#

Transform to a new reference frame.

Parameters:
frameFrameBase

The frame to transform to.

current_framegala.potential.CFrameBase

The current frame the phase-space position is in.

**kwargs

Any additional arguments are passed through to the individual frame transformation functions (see: transformations).

Returns:
pspgala.dynamics.PhaseSpacePosition

The phase-space position in the new reference frame.

to_hdf5(f)[source]#

Serialize this object to an HDF5 file.

Requires h5py.

Parameters:
fstr, h5py.File

Either the filename or an open HDF5 file.

w(units=None)[source]#

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

Parameters:
unitsUnitSystem (optional)

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

Returns:
wndarray

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