Orbit#
- class gala.dynamics.Orbit(pos, vel, t=None, hamiltonian=None, potential=None, frame=None)[source]#
Bases:
PhaseSpacePosition
Represents an orbit: positions and velocities (conjugate momenta) as a function of time.
The class can be instantiated with Astropy representation objects (e.g.,
CartesianRepresentation
), AstropyQuantity
objects, or plain Numpy arrays.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 axes 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 shapepos.shape = (3, 100)
, this would be a 3D orbit at 100 times (pos[0]
isx
,pos[1]`
isy
, etc.). For representing multiple orbits, the position array could have 3 axes, e.g., it might have shapepos.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
representation
,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.- vel
differential
,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.- tarray_like,
Quantity
(optional) Array of times. If a numpy array (e.g., has no units), this will be stored as a dimensionless
Quantity
.- hamiltonian
Hamiltonian
(optional) The Hamiltonian that the orbit was integrated in.
- pos
Represents phase-space positions, i.e. positions and conjugate momenta (velocities).
The class can be instantiated with Astropy representation objects (e.g.,
CartesianRepresentation
), AstropyQuantity
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 shapepos.shape = (3, 100)
, this would represent 100 3D positions (pos[0]
isx
,pos[1]
isy
, etc.). The same is true for velocity.- Parameters:
- pos
representation
,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.- vel
differential
,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.- frame
FrameBase
(optional) The reference frame of the input phase-space positions.
- pos
Attributes Summary
This is not the shape of the position or velocity arrays.
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.
Compute the angular momentum for the phase-space positions contained in this object.
animate
([components, units, stride, ...])Animate an orbit or collection of orbits.
apocenter
([return_times, func, approximate])Estimate the apocenter(s) of the orbit by identifying local maxima in the spherical radius, fitting a parabola around these local maxima and then solving this parabola to find the apocenter(s).
Determine which axes the Orbit circulates around by checking whether there is a change of sign of the angular momentum about an axis.
eccentricity
(**kw)Returns the eccentricity computed from the mean apocenter and mean pericenter.
energy
([hamiltonian])The total energy per unit mass:
estimate_period
([radial])Estimate the period of the orbit.
from_galpy_orbit
(galpy_orbit)Create a Gala
PhaseSpacePosition
orOrbit
instance from agalpy.Orbit
instance.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
The kinetic energy per unit mass:
Generator for iterating over each orbit.
pericenter
([return_times, func, approximate])Estimate the pericenter(s) of the orbit by identifying local minima in the spherical radius, fitting a parabola around these local minima and then solving this parabola to find the pericenter(s).
plot
([components, units, auto_aspect])Plot the positions in all projections.
plot_3d
([components, units, auto_aspect, ...])Plot the specified 3D components.
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 different reference frame.
to_galpy_orbit
([ro, vo])Convert this object to a
galpy.Orbit
instance.to_hdf5
(f)Serialize this object to an HDF5 file.
w
([units])This returns a single array containing the phase-space positions.
zmax
([return_times, func, approximate])Estimate the maximum
z
height of the orbit by identifying local maxima in the absolute value of thez
position, fitting a parabola around these local maxima and then solving this parabola to find the maximumz
height.Attributes Documentation
- data#
- hamiltonian#
- norbits#
- ntimes#
- 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
- 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:
- circulationarray_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
Orbit
A copy of the original orbit object with circulation aligned with the z axis.
- orb
- angular_momentum()#
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:
- L
Quantity
Array of angular momentum vectors.
- L
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>
- animate(components=None, units=None, stride=1, segment_nsteps=10, underplot_full_orbit=True, show_time=True, marker_style=None, segment_style=None, FuncAnimation_kwargs=None, orbit_plot_kwargs=None, axes=None)[source]#
Animate an orbit or collection of orbits.
- 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['v_x', 'v_y', 'v_z']
. If the representation is different, the component names will be different. For example, for a Cylindrical representation, the components are['rho', 'phi', 'z']
and['v_rho', 'pm_phi', 'v_z']
.- units
UnitBase
, iterable (optional) A single unit or list of units to display the components in.
- stride
int
(optional) How often to draw a new frame, in terms of orbit timesteps.
- segment_nsteps
int
(optional) How many timesteps to draw in an orbit segment trailing the timestep marker. Set this to 0 or None to disable.
- underplot_full_orbitbool (optional)
Controls whether to under-plot the full orbit as a thin line.
- show_timebool (optional)
Controls whether to show a label of the current timestep
- marker_style
dict
orlist
ofdict
(optional) Matplotlib style arguments passed to
matplotlib.pyplot.plot
that control the plot style of the timestep marker. If a single dict is passed then the marker_style is applied to all orbits. If a list of dicts is passed then each dict will be applied to each orbit.- segment_style
dict
(optional) Matplotlib style arguments passed to
matplotlib.pyplot.plot
that control the plot style of the orbit segment. If a single dict is passed then the segment_style is applied to all orbits. If a list of dicts is passed then each dict will be applied to each orbit.- FuncAnimation_kwargs
dict
(optional) Keyword arguments passed through to
matplotlib.animation.FuncAnimation
.- orbit_plot_kwargs
dict
(optional) Keyword arguments passed through to
gala.dynamics.Orbit.plot
.- axes
matplotlib.axes.Axes
(optional) Where to draw the orbit.
- Returns:
- apocenter(return_times=False, func=<function mean>, approximate=False)[source]#
Estimate the apocenter(s) of the orbit by identifying local maxima in the spherical radius, fitting a parabola around these local maxima and then solving this parabola to find the apocenter(s).
By default, this returns the mean of all local maxima (apocenters). To get, e.g., the largest apocenter, pass in
func=np.max
. To get all apocenters, pass infunc=None
.- Parameters:
- Returns:
- 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)
.
- circulation
- eccentricity(**kw)[source]#
Returns the eccentricity computed from the mean apocenter and mean pericenter.
\[e = \frac{r_{\rm apo} - r_{\rm per}}{r_{\rm apo} + r_{\rm per}}\]- Parameters:
- **kw
Any keyword arguments passed to
apocenter()
andpericenter()
. For example,approximate=True
.
- Returns:
- ecc
float
The orbital eccentricity.
- ecc
- energy(hamiltonian=None)[source]#
The total energy per unit mass:
- Parameters:
- hamiltonian
gala.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.
- hamiltonian
- Returns:
- E
Quantity
The total energy.
- E
- estimate_period(radial=False)[source]#
Estimate the period of the orbit. By default, computes the radial period. If
radial==False
, this returns period estimates for each dimension of the orbit.
- classmethod from_galpy_orbit(galpy_orbit)[source]#
Create a Gala
PhaseSpacePosition
orOrbit
instance from agalpy.Orbit
instance.- Parameters:
- galpy_orbit
galpy.orbit.Orbit
- galpy_orbit
- Returns:
- orbit
Orbit
- orbit
- classmethod from_w(w, units=None, **kwargs)#
- get_components(which)#
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:
- which
str
Can either be
'pos'
or'vel'
to get the components for the position or velocity object.
- which
- kinetic_energy()#
The kinetic energy per unit mass:
\[E_K = \frac{1}{2} \, |\boldsymbol{v}|^2\]- Returns:
- E
Quantity
The kinetic energy.
- E
- pericenter(return_times=False, func=<function mean>, approximate=False)[source]#
Estimate the pericenter(s) of the orbit by identifying local minima in the spherical radius, fitting a parabola around these local minima and then solving this parabola to find the pericenter(s).
By default, this returns the mean of all local minima (pericenters). To get, e.g., the minimum pericenter, pass in
func=np.min
. To get all pericenters, pass infunc=None
.- Parameters:
- Returns:
- 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['v_x', 'v_y', 'v_z']
. If the representation is different, the component names will be different. For example, for a Cylindrical representation, the components are['rho', 'phi', 'z']
and['v_rho', 'pm_phi', 'v_z']
.- units
UnitBase
, 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_kwargs
dict
(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_function
callable()
(optional) The
matplotlib
plot function to use. By default, this isscatter()
, 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 likecolor=...
,marker=...
, etc.
- Returns:
- fig
Figure
- fig
- plot_3d(components=None, units=None, auto_aspect=True, subplots_kwargs=None, **kwargs)[source]#
Plot the specified 3D components.
- 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['v_x', 'v_y', 'v_z']
. If the representation is different, the component names will be different. For example, for a Cylindrical representation, the components are['rho', 'phi', 'z']
and['v_rho', 'pm_phi', 'v_z']
.- units
UnitBase
, iterable (optional) A single unit or list of units to display the components in.
- auto_aspectbool (optional)
Automatically enforce an equal aspect ratio.
- ax
matplotlib.axes.Axes
The matplotlib Axes object to draw on.
- subplots_kwargs
dict
(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_function
str
(optional) The
matplotlib
plot function to use. By default, this is ‘plot’ but can also be, e.g., ‘scatter’.- **kwargs
All other keyword arguments are passed to the
plot_function
. You can pass in any of the usual style kwargs likecolor=...
,marker=...
, etc.
- Returns:
- fig
Figure
- fig
- potential_energy(potential=None)[source]#
The potential energy per unit mass:
\[E_\Phi = \Phi(\boldsymbol{q})\]- Returns:
- E
Quantity
The potential energy.
- E
- 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_pos
BaseRepresentation
The type of representation to generate. Must be a class (not an instance), or the string name of the representation class.
- new_vel
BaseDifferential
(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.
- new_pos
- Returns:
- new_orbit
gala.dynamics.Orbit
- new_orbit
- to_coord_frame(frame, galactocentric_frame=None, **kwargs)#
Transform the orbit from Galactocentric, cartesian coordinates to Heliocentric coordinates in the specified Astropy coordinate frame.
- Parameters:
- frame
BaseCoordinateFrame
The frame instance specifying the desired output frame. For example,
ICRS
.- galactocentric_frame
Galactocentric
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.
- frame
- Returns:
- c
BaseCoordinateFrame
An instantiated coordinate frame containing the positions and velocities from this object transformed to the specified coordinate frame.
- c
- to_frame(frame, current_frame=None, **kwargs)[source]#
Transform to a different reference frame.
- Parameters:
- frame
gala.potential.CFrameBase
The frame to transform to.
- current_frame
gala.potential.CFrameBase
(optional) If the Orbit has no associated Hamiltonian, this specifies the current frame of the orbit.
- frame
- Returns:
- orbit
gala.dynamics.Orbit
The orbit in the new reference frame.
- orbit
- to_galpy_orbit(ro=None, vo=None)[source]#
Convert this object to a
galpy.Orbit
instance.- Parameters:
- ro
astropy.units.Quantity
orastropy.units.UnitBase
“Natural” length unit.
- vo
astropy.units.Quantity
orastropy.units.UnitBase
“Natural” velocity unit.
- ro
- Returns:
- galpy_orbit
galpy.orbit.Orbit
- galpy_orbit
- 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.
- units
- Returns:
- w
ndarray
A numpy array of all positions and velocities, without units. Will have shape
(2*ndim, ...)
.
- w
- zmax(return_times=False, func=<function mean>, approximate=False)[source]#
Estimate the maximum
z
height of the orbit by identifying local maxima in the absolute value of thez
position, fitting a parabola around these local maxima and then solving this parabola to find the maximumz
height.By default, this returns the mean of all local maxima. To get, e.g., the largest
z
excursion, pass infunc=np.max
. To get allz
maxima, pass infunc=None
.- Parameters:
- Returns: