astrobase.coordutils module

Contains various useful tools for coordinate conversion, etc.

astrobase.coordutils.angle_wrap(angle, radians=False)[source]

Wraps the input angle to 360.0 degrees.

Parameters:
  • angle (float) – The angle to wrap around 360.0 deg.
  • radians (bool) – If True, will assume that the input is in radians. The output will then also be in radians.
Returns:

Wrapped angle. If radians is True: input is assumed to be in radians, output is also in radians.

Return type:

float

astrobase.coordutils.decimal_to_dms(decimal_value)[source]

Converts from decimal degrees (for declination coords) to DD:MM:SS.

Parameters:decimal_value (float) – A decimal value to convert to degrees, minutes, seconds sexagesimal format.
Returns:A four element tuple is returned: (sign, HH, MM, SS.ssss…)
Return type:tuple
astrobase.coordutils.decimal_to_hms(decimal_value)[source]

Converts from decimal degrees (for RA coords) to HH:MM:SS.

Parameters:decimal_value (float) – A decimal value to convert to hours, minutes, seconds. Negative values will be wrapped around 360.0.
Returns:A three element tuple is returned: (HH, MM, SS.ssss…)
Return type:tuple
astrobase.coordutils.hms_str_to_tuple(hms_string)[source]

Converts a string of the form HH:MM:SS or HH MM SS to a tuple of the form (HH, MM, SS).

Parameters:hms_string (str) – A RA coordinate string of the form ‘HH:MM:SS.sss’ or ‘HH MM SS.sss’.
Returns:A three element tuple is returned (HH, MM, SS.ssss…)
Return type:tuple
astrobase.coordutils.dms_str_to_tuple(dms_string)[source]

Converts a string of the form [+-]DD:MM:SS or [+-]DD MM SS to a tuple of the form (sign, DD, MM, SS).

Parameters:dms_string (str) – A declination coordinate string of the form ‘[+-]DD:MM:SS.sss’ or ‘[+-]DD MM SS.sss’. The sign in front of DD is optional. If it’s not there, this function will assume that the coordinate string is a positive value.
Returns:A four element tuple of the form: (sign, DD, MM, SS.ssss…).
Return type:tuple
astrobase.coordutils.hms_str_to_decimal(hms_string)[source]

Converts a HH:MM:SS string to decimal degrees.

Parameters:hms_string (str) – A right ascension coordinate string of the form: ‘HH:MM:SS.sss’ or ‘HH MM SS.sss’.
Returns:The RA value in decimal degrees (wrapped around 360.0 deg if necessary.)
Return type:float
astrobase.coordutils.dms_str_to_decimal(dms_string)[source]

Converts a DD:MM:SS string to decimal degrees.

Parameters:dms_string (str) – A declination coordinate string of the form: ‘[+-]DD:MM:SS.sss’ or ‘[+-]DD MM SS.sss’.
Returns:The declination value in decimal degrees.
Return type:float
astrobase.coordutils.hms_to_decimal(hours, minutes, seconds, returndeg=True)[source]

Converts from HH, MM, SS to a decimal value.

Parameters:
  • hours (int) – The HH part of a RA coordinate.
  • minutes (int) – The MM part of a RA coordinate.
  • seconds (float) – The SS.sss part of a RA coordinate.
  • returndeg (bool) – If this is True, then will return decimal degrees as the output. If this is False, then will return decimal HOURS as the output. Decimal hours are sometimes used in FITS headers.
Returns:

The right ascension value in either decimal degrees or decimal hours depending on returndeg.

Return type:

float

astrobase.coordutils.dms_to_decimal(sign, degrees, minutes, seconds)[source]

Converts from DD:MM:SS to a decimal value.

Parameters:
  • sign ({'+', '-', ''}) – The sign part of a Dec coordinate.
  • degrees (int) – The DD part of a Dec coordinate.
  • minutes (int) – The MM part of a Dec coordinate.
  • seconds (float) – The SS.sss part of a Dec coordinate.
Returns:

The declination value in decimal degrees.

Return type:

float

astrobase.coordutils.great_circle_dist(ra1, dec1, ra2, dec2)[source]

Calculates the great circle angular distance between two coords.

This calculates the great circle angular distance in arcseconds between two coordinates (ra1,dec1) and (ra2,dec2). This is basically a clone of GCIRC from the IDL Astrolib.

Parameters:
  • ra1,dec1 (float or array-like) – The first coordinate’s right ascension and declination value(s) in decimal degrees.
  • ra2,dec2 (float or array-like) – The second coordinate’s right ascension and declination value(s) in decimal degrees.
Returns:

Great circle distance between the two coordinates in arseconds.

Return type:

float or array-like

Notes

If (ra1, dec1) is scalar and (ra2, dec2) is scalar: the result is a float distance in arcseconds.

If (ra1, dec1) is scalar and (ra2, dec2) is array-like: the result is an np.array with distance in arcseconds between (ra1, dec1) and each element of (ra2, dec2).

If (ra1, dec1) is array-like and (ra2, dec2) is scalar: the result is an np.array with distance in arcseconds between (ra2, dec2) and each element of (ra1, dec1).

If (ra1, dec1) and (ra2, dec2) are both array-like: the result is an np.array with the pair-wise distance in arcseconds between each element of the two coordinate lists. In this case, if the input array-likes are not the same length, then excess elements of the longer one will be ignored.

astrobase.coordutils.xmatch_basic(ra1, dec1, ra2, dec2, match_radius=5.0)[source]

Finds the closest object in (ra2, dec2) to scalar coordinate pair (ra1, dec1) and returns the distance in arcseconds.

This is a quick matcher that uses the great_circle_dist function to find the closest object in (ra2, dec2) within match_radius arcseconds to (ra1, dec1). (ra1, dec1) must be a scalar pair, while (ra2, dec2) must be array-likes of the same lengths.

Parameters:
  • ra1,dec1 (float) – Coordinate of the object to find matches to. In decimal degrees.
  • ra2,dec2 (array-like) – The coordinates that will be searched for matches. In decimal degrees.
  • match_radius (float) – The match radius in arcseconds to use for the match.
Returns:

A two element tuple like the following:

(True -> no match found or False -> found a match,
 minimum distance between target and list in arcseconds)

Return type:

tuple

astrobase.coordutils.xmatch_neighbors(ra1, dec1, ra2, dec2, match_radius=60.0, includeself=False, sortresults=True)[source]

Finds the closest objects in (ra2, dec2) to scalar coordinate pair (ra1, dec1) and returns the indices of the objects that match.

This is a quick matcher that uses the great_circle_dist function to find the closest object in (ra2, dec2) within match_radius arcseconds to (ra1, dec1). (ra1, dec1) must be a scalar pair, while (ra2, dec2) must be array-likes of the same lengths.

Parameters:
  • ra1,dec1 (float) – Coordinate of the object to find matches to. In decimal degrees.
  • ra2,dec2 (array-like) – The coordinates that will be searched for matches. In decimal degrees.
  • match_radius (float) – The match radius in arcseconds to use for the match.
  • includeself (bool) – If this is True, the object itself will be included in the match results.
  • sortresults (bool) – If this is True, the match indices will be sorted by distance.
Returns:

A tuple like the following is returned:

(True -> matches found or False -> no matches found,
 minimum distance between target and list,
 np.array of indices where list of coordinates is
 closer than `match_radius` arcseconds from the target,
 np.array of distances in arcseconds)

Return type:

tuple

astrobase.coordutils.make_kdtree(ra, decl)[source]

This makes a scipy.spatial.CKDTree on (ra, decl).

Parameters:ra,decl (array-like) – The right ascension and declination coordinate pairs in decimal degrees.
Returns:The cKDTRee object generated by this function is returned and can be used to run various spatial queries.
Return type:scipy.spatial.CKDTree
astrobase.coordutils.conesearch_kdtree(kdtree, racenter, declcenter, searchradiusdeg, conesearchworkers=1)[source]

This does a cone-search around (racenter, declcenter) in kdtree.

Parameters:
  • kdtree (scipy.spatial.CKDTree) – This is a kdtree object generated by the make_kdtree function.
  • racenter,declcenter (float or array-like) – This is the center coordinate to run the cone-search around in decimal degrees. If this is an np.array, will search for all coordinate pairs in the array.
  • searchradiusdeg (float) – The search radius to use for the cone-search in decimal degrees.
  • conesearchworkers (int) – The number of parallel workers to launch for the cone-search.
Returns:

If (racenter, declcenter) is a single coordinate, this will return a list of the indices of the matching objects in the kdtree. If (racenter, declcenter) are array-likes, this will return an object array containing lists of matching object indices for each coordinate searched.

Return type:

list or np.array of lists

astrobase.coordutils.xmatch_kdtree(kdtree, extra, extdecl, xmatchdistdeg, closestonly=True)[source]

This cross-matches between kdtree and (extra, extdecl) arrays.

Returns the indices of the kdtree and the indices of extra, extdecl that xmatch successfully.

Parameters:
  • kdtree (scipy.spatial.CKDTree) – This is a kdtree object generated by the make_kdtree function.
  • extra,extdecl (array-like) – These are np.arrays of ‘external’ coordinates in decimal degrees that will be cross-matched against the objects in kdtree.
  • xmatchdistdeg (float) – The match radius to use for the cross-match in decimal degrees.
  • closestonly (bool) – If closestonly is True, then this function returns only the closest matching indices in (extra, extdecl) for each object in kdtree if there are any matches. Otherwise, it returns a list of indices in (extra, extdecl) for all matches within xmatchdistdeg between kdtree and (extra, extdecl).
Returns:

Returns a tuple of the form:

(list of `kdtree` indices matching to external objects,
 list of all `extra`/`extdecl` indices that match to each
 element in `kdtree` within the specified cross-match distance)

Return type:

tuple of lists

astrobase.coordutils.total_proper_motion(pmra, pmdecl, decl)[source]

This calculates the total proper motion of an object.

Parameters:
  • pmra (float or array-like) – The proper motion(s) in right ascension, measured in mas/yr.
  • pmdecl (float or array-like) – The proper motion(s) in declination, measured in mas/yr.
  • decl (float or array-like) – The declination of the object(s) in decimal degrees.
Returns:

The total proper motion(s) of the object(s) in mas/yr.

Return type:

float or array-like

astrobase.coordutils.reduced_proper_motion(mag, propermotion)[source]

This calculates the reduced proper motion using the mag measurement provided.

Parameters:
  • mag (float or array-like) – The magnitude(s) to use to calculate the reduced proper motion(s).
  • propermotion (float or array-like) – The total proper motion of the object(s). Use the total_proper_motion function to calculate this if you have pmra, pmdecl, and decl values. propermotion should be in mas/yr.
Returns:

The reduced proper motion for the object(s). This is effectively a measure of the absolute magnitude in the band provided.

Return type:

float or array-like

astrobase.coordutils.equatorial_to_galactic(ra, decl, equinox='J2000')[source]

This converts from equatorial coords to galactic coords.

Parameters:
  • ra (float or array-like) – Right ascension values(s) in decimal degrees.
  • decl (float or array-like) – Declination value(s) in decimal degrees.
  • equinox (str) – The equinox that the coordinates are measured at. This must be recognizable by Astropy’s SkyCoord class.
Returns:

The galactic coordinates (l, b) for each element of the input (ra, decl).

Return type:

tuple of (float, float) or tuple of (np.array, np.array)

astrobase.coordutils.galactic_to_equatorial(gl, gb)[source]

This converts from galactic coords to equatorial coordinates.

Parameters:
  • gl (float or array-like) – Galactic longitude values(s) in decimal degrees.
  • gb (float or array-like) – Galactic latitude value(s) in decimal degrees.
Returns:

The equatorial coordinates (RA, DEC) for each element of the input (gl, gb) in decimal degrees. These are reported in the ICRS frame.

Return type:

tuple of (float, float) or tuple of (np.array, np.array)

astrobase.coordutils.xieta_from_radecl(inra, indecl, incenterra, incenterdecl, deg=True)[source]

This returns the image-plane projected xi-eta coords for inra, indecl.

Parameters:
  • inra,indecl (array-like) – The equatorial coordinates to get the xi, eta coordinates for in decimal degrees or radians.
  • incenterra,incenterdecl (float) – The center coordinate values to use to calculate the plane-projected coordinates around.
  • deg (bool) – If this is True, the input angles are assumed to be in degrees and the output is in degrees as well.
Returns:

This is the (xi, eta) coordinate pairs corresponding to the image-plane projected coordinates for each pair of input equatorial coordinates in (inra, indecl).

Return type:

tuple of np.arrays