vhel_to_gal

gala.coordinates.vhel_to_gal(coordinate, pm, rv, vcirc=<Quantity 220.0 km / s>, vlsr=<Quantity [ 10., 5.25, 7.17] km / s>, galactocentric_frame=None)[source]

Convert a Heliocentric velocity in spherical coordinates (e.g., proper motion and radial velocity) in the ICRS or Galactic frame to a Galactocentric, cartesian velocity.

The frame of the input coordinate determines how to interpret the given proper motions. For example, if the input coordinate is in the ICRS frame, the proper motions are assumed to be \((\mu_\alpha\cos\delta,\mu_\delta)\). This function also handles array inputs (see examples below).

TODO: Roundtrip using galactic coordinates only maintains relative precision of ~1E-5. Why?

Parameters:

coordinate : SkyCoord, BaseCoordinateFrame

This is most commonly a SkyCoord object, but alternatively, it can be any coordinate frame object that is transformable to the Galactocentric frame.

pm : Quantity or iterable of Quantity objects

Proper motion in the same frame as the coordinate. For example, if your input coordinate is in ICRS, then the proper motion is assumed to be in this frame as well. The order of elements should always be proper motion in (longitude, latitude), and should have shape (2,N). The longitude component is assumed to have the cosine of the latitude already multiplied in, so that in ICRS, for example, this would be \(\mu_\alpha\cos\delta\).

rv : Quantity

Barycentric radial velocity. Should have shape (1,N) or (N,).

vcirc : Quantity (optional)

Circular velocity of the Sun.

vlsr : Quantity (optional)

Velocity of the Sun relative to the local standard of rest (LSR).

galactocentric_frame : Galactocentric (optional)

An instantiated Galactocentric frame object with custom parameters for the Galactocentric coordinates. For example, if you want to set your own position of the Galactic center, you can pass in a frame with custom galcen_ra and galcen_dec.

Returns:

vxyz : Quantity (optional)

Cartesian velocity components (U,V,W). A Quantity object with shape (3,N).

Examples

>>> import astropy.units as u
>>> import astropy.coordinates as coord
>>> c = coord.SkyCoord(ra=196.5*u.degree, dec=-10.33*u.deg, distance=16.2*u.kpc)
>>> pm = [-1.53, 3.5]*u.mas/u.yr
>>> rv = 161.4*u.km/u.s
>>> vhel_to_gal(c, pm=pm, rv=rv) 
<Quantity [-137.29984564, 262.64052249, 305.50786499] km / s>
>>> c = coord.SkyCoord(ra=[196.5,51.3]*u.degree, dec=[-10.33,2.1]*u.deg, distance=[16.2,11.]*u.kpc)
>>> pm = [[-1.53,4.5], [3.5,10.9]]*u.mas/u.yr
>>> rv = [161.4,-210.2]*u.km/u.s
>>> vhel_to_gal(c, pm=pm, rv=rv) 
<Quantity [[-137.29984564,-212.10415701],
           [ 262.64052249, 496.85687803],
           [ 305.50786499, 554.16562628]] km / s>