Coordinate Systems (gala.coordinates)#

Introduction#

The coordinates subpackage primarily provides specialty astropy.coordinates frame classes for coordinate systems defined by the stellar streams, and for other common Galactic dynamics tasks like removing solar reflex motion from proper motions or radial velocities, and transforming a proper motion covariance matrix from one frame to another.

For the examples below the following imports have already been executed:

>>> import numpy as np
>>> import astropy.coordinates as coord
>>> import astropy.units as u
>>> import gala.coordinates as gc

We will also set the default Astropy Galactocentric frame parameters to the values adopted in Astropy v4.0:

>>> _ = coord.galactocentric_frame_defaults.set('v4.0')

Stellar stream coordinate frames#

gala provides Astropy coordinate frame classes for transforming to several built-in stellar stream stream coordinate frames (as defined in the references below), and for transforming positions and velocities to and from coordinate systems defined by great circles or poles. These classes behave like the built-in astropy coordinates frames (e.g., ICRS or Galactic) and can be transformed to and from other astropy coordinate frames. For example, to convert a set of ICRS (RA, Dec) coordinates to a coordinate system aligned with the Sagittarius stream with the SagittariusLaw10 frame:

>>> c = coord.ICRS(ra=100.68458*u.degree, dec=41.26917*u.degree)
>>> sgr = c.transform_to(gc.SagittariusLaw10())
>>> (sgr.Lambda, sgr.Beta) 
(<Longitude 179.58511053544734 deg>, <Latitude -12.558450192162654 deg>)

Or, to transform from SagittariusLaw10 coordinates to the Galactic frame:

>>> sgr = gc.SagittariusLaw10(Lambda=156.342*u.degree, Beta=1.1*u.degree)
>>> c = sgr.transform_to(coord.Galactic())
>>> (c.l, c.b) 
(<Longitude 182.5922090437946 deg>, <Latitude -9.539692094685893 deg>)

These transformations also handle velocities so that proper motion components can be transformed between the systems. For example, to transform from GD1Koposov10 proper motions to Galactic proper motions:

>>> gd1 = gc.GD1Koposov10(phi1=-35*u.degree, phi2=0*u.degree,
...                       pm_phi1_cosphi2=-12.20*u.mas/u.yr,
...                       pm_phi2=-3.10*u.mas/u.yr)
>>> gd1.transform_to(coord.Galactic()) 
<Galactic Coordinate: (l, b) in deg
    (181.28968151, 54.84972806)
 (pm_l_cosb, pm_b) in mas / yr
    (12.03209393, -3.69847479)>

As with the other Astropy coordinate frames, with a full specification of the 3D position and velocity, we can transform to a Galactocentric frame:

>>> gd1 = gc.GD1Koposov10(phi1=-35.00*u.degree, phi2=0.04*u.degree,
...                       distance=7.83*u.kpc,
...                       pm_phi1_cosphi2=-12.20*u.mas/u.yr,
...                       pm_phi2=-3.10*u.mas/u.yr,
...                       radial_velocity=-32*u.km/u.s)
>>> gd1.transform_to(coord.Galactocentric()) 
<Galactocentric Coordinate (galcen_coord=<ICRS Coordinate: (ra, dec) in deg
    (266.4051, -28.936175)>, galcen_distance=8.122 kpc, galcen_v_sun=(12.9, 245.6, 7.78) km / s, z_sun=20.8 pc, roll=0.0 deg): (x, y, z) in kpc
    (-12.61622659, -0.09870921, 6.43179403)
(v_x, v_y, v_z) in km / s
    (-71.14675268, -203.01648654, -97.12884319)>

For custom great circle coordinate systems, and for more information about the stellar stream frames, see Great circle and stellar stream coordinate frames.

Correcting velocities for solar reflex motion#

The reflex_correct function accepts an Astropy SkyCoord object with position and velocity information, and returns a coordinate object with the solar motion added back in to the velocity components. This is useful for computing velocities in a Galactocentric reference frame, rather than a solar system barycentric frame.

The reflex_correct function accepts a coordinate object with scalar or array values:

>>> c = coord.SkyCoord(ra=[180.323, 1.523]*u.deg,
...                    dec=[-17, 29]*u.deg,
...                    distance=[172, 412]*u.pc,
...                    pm_ra_cosdec=[-11, 3]*u.mas/u.yr,
...                    pm_dec=[4, 8]*u.mas/u.yr,
...                    radial_velocity=[114, -21]*u.km/u.s)
>>> gc.reflex_correct(c) 
<SkyCoord (ICRS): (ra, dec, distance) in (deg, deg, pc)
    [(180.323, -17., 172.), (  1.523,  29., 412.)]
(pm_ra_cosdec, pm_dec, radial_velocity) in (mas / yr, mas / yr, km / s)
    [(139.47001884, 175.45769809, -47.09032586),
    (-61.01738781,  61.51055793, 163.36721898)]>

By default, this uses the solar location and velocity from the astropy.coordinates.Galactocentric frame class. To modify these parameters, for example, to change the solar velocity, or the sun’s height above the Galactic midplane, use the arguments of the astropy.coordinates.Galactocentric class and pass in an instance of the astropy.coordinates.Galactocentric frame:

>>> vsun = coord.CartesianDifferential([11., 245., 7.]*u.km/u.s)
>>> gc_frame = coord.Galactocentric(galcen_v_sun=vsun, z_sun=0*u.pc)
>>> gc.reflex_correct(c, gc_frame) 
<SkyCoord (ICRS): (ra, dec, distance) in (deg, deg, pc)
    [(180.323, -17., 172.), (  1.523,  29., 412.)]
(pm_ra_cosdec, pm_dec, radial_velocity) in (mas / yr, mas / yr, km / s)
    [(136.93481249, 175.37627916, -47.6177433 ),
    (-59.96484921,  61.41044742, 163.90707073)]>

If you don’t have radial velocity information and want to correct the proper motions, pass in zeros for the radial velocity (and ignore the output value of the radial velocity):

>>> c = coord.SkyCoord(ra=162*u.deg,
...                    dec=-17*u.deg,
...                    distance=172*u.pc,
...                    pm_ra_cosdec=-11*u.mas/u.yr,
...                    pm_dec=4*u.mas/u.yr,
...                    radial_velocity=0*u.km/u.s)
>>> gc.reflex_correct(c) 
<SkyCoord (ICRS): (ra, dec, distance) in (deg, deg, pc)
    (162., -17., 172.)
(pm_ra_cosdec, pm_dec, radial_velocity) in (mas / yr, mas / yr, km / s)
    (88.20380175, 163.68500525, -192.48721942)>

Similarly, if you don’t have proper motion information and want to correct the proper motions, pass in zeros for the proper motions (and ignore the output values of the proper motions) – this is sometimes called “v_GSR”:

>>> c = coord.SkyCoord(ra=162*u.deg,
...                    dec=-17*u.deg,
...                    distance=172*u.pc,
...                    pm_ra_cosdec=0*u.mas/u.yr,
...                    pm_dec=0*u.mas/u.yr,
...                    radial_velocity=127*u.km/u.s)
>>> gc.reflex_correct(c) 
<SkyCoord (ICRS): (ra, dec, distance) in (deg, deg, pc)
    (162., -17., 172.)
(pm_ra_cosdec, pm_dec, radial_velocity) in (mas / yr, mas / yr, km / s)
    (99.20380175, 159.68500525, -65.48721942)>

Transforming a proper motion covariance matrix to a new coordinate frame#

When working with Gaia or other astrometric data sets, you may need to transform the reported covariance matrix between proper motion components into a new coordinate system. For example, Gaia data are provided in the ICRS (equatorial) coordinate frame, but for Galactic science, we often want to instead work in the Galactic coordinate system. For this and other transformations that only require a rotation (i.e. the origin doesn’t change), the astrometric covariance matrix can be transformed exactly through a projection of the rotation onto the tangent plane at a given location. The details of this procedure are explained in this document from the Gaia data processing team, and this functionality is implemented in gala. Let’s first create a coordinate object to transform:

>>> c = coord.SkyCoord(ra=62*u.deg,
...                    dec=17*u.deg,
...                    pm_ra_cosdec=1*u.mas/u.yr,
...                    pm_dec=3*u.mas/u.yr)

and a covariance matrix for the proper motion components, for example, as would be constructed from a single row from a Gaia data release source catalog:

>>> cov = np.array([[0.53510132, 0.16637034],
...                 [0.16637034, 1.1235292 ]])

This matrix specifies the 2D error distribution for the proper motion measurement in the ICRS frame. To transform this matrix to, e.g., the Galactic coordinate system, use the function transform_pm_cov:

>>> gc.transform_pm_cov(c, cov, coord.Galactic()) 
array([[ 0.69450047, -0.309945  ],
       [-0.309945  ,  0.96413005]])

Note that this also works for all of the great circle or stellar stream coordinate frames implemented in gala:

>>> gc.transform_pm_cov(c, cov, gc.GD1Koposov10()) 
array([[1.10838914, 0.19067958],
       [0.19067958, 0.55024138]])

This works for array-valued coordinates as well, so try to avoid looping over this function and instead apply it to array-valued coordinate objects.

References#

Using gala.coordinates#

More details are provided in the linked pages below:

API#

gala.coordinates Package#

Functions#

cartesian_to_poincare_polar(w)

Convert an array of 6D Cartesian positions to Poincaré symplectic polar coordinates.

make_greatcircle_cls(cls_name[, ...])

pole_from_endpoints(coord1, coord2)

Compute the pole from a great circle that connects the two specified coordinates.

reflex_correct(coords[, galactocentric_frame])

Correct the input Astropy coordinate object for solar reflex motion.

transform_pm_cov(c, cov, to_frame)

Transform a proper motion covariance matrix to a new frame.

vgsr_to_vhel(coordinate, vgsr[, vsun])

Convert a radial velocity in the Galactic standard of rest (GSR) to a barycentric radial velocity.

vhel_to_vgsr(coordinate, vhel, vsun)

Convert a velocity from a heliocentric radial velocity to the Galactic standard of rest (GSR).

Classes#

GD1Koposov10(*args, **kwargs)

A Heliocentric spherical coordinate system defined by the orbit of the GD1 stream, as described in Koposov et al. 2010 (see: http://arxiv.org/abs/0907.1085).

GreatCircleICRSFrame(*args, **kwargs)

A coordinate frame defined by a pole and origin.

JhelumBonaca19(*args, **kwargs)

A Heliocentric spherical coordinate system defined by the orbit of the Jhelum stream, as described in Bonaca et al. 2019.

MagellanicStreamNidever08(*args, **kwargs)

A coordinate or frame aligned with the Magellanic Stream, as defined by Nidever et al. (2008, see: http://adsabs.harvard.edu/abs/2008ApJ...679..432N).

OphiuchusPriceWhelan16(*args, **kwargs)

A Heliocentric spherical coordinate system defined by the orbit of the Ophiuchus stream, as described in Price-Whelan et al. 2016 (see: https://arxiv.org/abs/1601.06790).

OrphanKoposov19(*args, **kwargs)

A coordinate frame for the Orphan stream defined by Sergey Koposov.

OrphanNewberg10(*args, **kwargs)

A Heliocentric spherical coordinate system defined by the orbit of the Orphan stream, as described in Newberg et al. 2010 (see: http://arxiv.org/abs/1001.0576).

Pal13Shipp20(*args, **kwargs)

A Heliocentric spherical coordinate system defined by the orbit of the Pal 13 stream by Shipp et al. (2020).

Pal5PriceWhelan18(*args, **kwargs)

A Heliocentric spherical coordinate system defined by the orbit of the Pal 5 stream by A.

SagittariusLaw10(*args, **kwargs)

A Heliocentric spherical coordinate system defined by the orbit of the Sagittarius dwarf galaxy, as described in http://adsabs.harvard.edu/abs/2003ApJ...599.1082M and further explained in http://www.stsci.edu/~dlaw/Sgr/.

SagittariusVasiliev21(*args, **kwargs)

A Heliocentric, right-handed spherical coordinate system defined by the orbit of the Sagittarius dwarf galaxy, as described in Vasiliev et al. 2021, https://ui.adsabs.harvard.edu/abs/2021MNRAS.501.2279V/abstract.