CartesianOrbit¶
-
class
gala.dynamics.orbit.CartesianOrbit(pos, vel, t=None, potential=None)[source]¶ Bases:
gala.dynamics.core.CartesianPhaseSpacePosition,gala.dynamics.OrbitRepresents 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 shapepos.shape = (3, 100), this would be a 3D orbit (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 :
Quantity, array_likePositions. 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_likeVelocities. 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
ndimnorbitsntimesrThe orbital radius. shapeMethods 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. angular_momentum()apocenter([type])Estimate the apocenter(s) of the orbit. circulation()Determine which axes the Orbit circulates around by checking whether there is a change of sign of the angular momentum about an axis. eccentricity()Returns the eccentricity computed from the mean apocenter and mean pericenter. energy([potential])The total energy per unit mass (e.g., kinetic + potential): estimate_period([radial])Estimate the period of the orbit. from_w(w[, units])kinetic_energy()The kinetic energy per unit mass: pericenter([type])Estimate the pericenter(s) of the orbit. plot(**kwargs)Plot the orbit 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. Attributes Documentation
-
ndim¶
-
norbits¶
-
ntimes¶
-
r¶ The orbital radius.
-
shape¶
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 :
CartesianOrbitA copy of the original orbit object with circulation aligned with the z axis.
-
angular_momentum()¶
-
apocenter(type=<function mean>)¶ Estimate the apocenter(s) of the orbit. By default, this returns the mean apocenter. To get, e.g., the minimum apocenter, pass in
type=np.min. To get all apocenters, pass intype=None.Parameters: type : func (optional)
By default, this returns the mean apocenter. To return all apocenters, pass in
None. To get, e.g., the minimum or maximum apocenter, pass innp.minornp.max.Returns: apo : float,
ndarrayEither a single number or an array of apocenters.
-
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
ndimintegers 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.ndarrayAn 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).
-
eccentricity()¶ 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}}\]Returns: ecc : float
The orbital eccentricity.
-
energy(potential=None)[source]¶ The total energy per unit mass (e.g., kinetic + potential):
Returns: E :
QuantityThe total energy.
-
estimate_period(radial=True)¶ 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.Parameters: radial : bool (optional)
What period to estimate. If
True, estimates the radial period. IfFalse, estimates period in each dimension, e.g., if the orbit is 3D, along x, y, and z.Returns: T :
QuantityThe period or periods.
-
from_w(w, units=None, **kwargs)¶
-
kinetic_energy()¶ The kinetic energy per unit mass:
\[E_K = \frac{1}{2} \, |\boldsymbol{v}|^2\]Returns: E :
QuantityThe kinetic energy.
-
pericenter(type=<function mean>)¶ Estimate the pericenter(s) of the orbit. By default, this returns the mean pericenter. To get, e.g., the minimum pericenter, pass in
type=np.min. To get all pericenters, pass intype=None.Parameters: type : func (optional)
By default, this returns the mean pericenter. To return all pericenters, pass in
None. To get, e.g., the minimum or maximum pericenter, pass innp.minornp.max.Returns: peri : float,
ndarrayEither a single number or an array of pericenters.
-
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
xis an array of shape(3,1024,32)- 1024 timesteps for 32 orbits in 3D positions –ixwould 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 likecolor=...,marker=..., etc.Returns: fig :
Figure
-
potential_energy(potential=None)[source]¶ The potential energy per unit mass:
\[E_\Phi = \Phi(\boldsymbol{q})\]Returns: E :
QuantityThe potential energy.
-
represent_as(Representation)¶ Represent the position and velocity of the orbit in an alternate coordinate system. Supports any of the Astropy coordinates representation classes.
Parameters: Representation :
BaseRepresentationThe class for the desired representation.
Returns: pos :
BaseRepresentationvel :
QuantityThe 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)¶ Transform the orbit from Galactocentric, cartesian coordinates to Heliocentric coordinates in the specified Astropy coordinate frame.
Parameters: frame :
BaseCoordinateFramegalactocentric_frame :
Galactocentricvcirc :
QuantityCircular velocity of the Sun. Passed to velocity transformation.
vlsr :
QuantityVelocity of the Sun relative to the LSR. Passed to velocity transformation.
Returns: An instantiated coordinate frame.
v : tuple
A tuple of velocities represented as
Quantityobjects.
-
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 :
ndarrayA numpy array of all positions and velocities, without units. Will have shape
(2*ndim,...).
-