join_skycoord¶
- astropy.table.join_skycoord(distance, distance_func='search_around_sky')[source]¶
Helper function to join on SkyCoord columns using distance matching.
This function is intended for use in
table.join()
to allow performing a table join where the key columns are bothSkyCoord
objects, matched by computing the distance between points and accepting values belowdistance
.The distance cross-matching is done using either
search_around_sky
orsearch_around_3d
, depending on the value ofdistance_func
. The default is'search_around_sky'
.One can also provide a function object for
distance_func
, in which case it must be a function that follows the same input and output API assearch_around_sky
. In this case the function will be called with(skycoord1, skycoord2, distance)
as arguments.- Parameters:
- distance
Quantity
[:ref: ‘angle’, :ref: ‘length’] Maximum distance between points to be considered a join match. Must have angular or distance units.
- distance_func
python:str
or python:function Specifies the function for performing the cross-match based on
distance
. If supplied as a string this specifies the name of a function inastropy.coordinates
. If supplied as a function then that function is called directly.
- distance
- Returns:
- join_funcpython:function
Function that accepts two
SkyCoord
columns (col1, col2) and returns the tuple (ids1, ids2) of pair-matched unique identifiers.
Examples
This example shows an inner join of two
SkyCoord
columns, taking any sources within 0.2 deg to be a match. Note the newsc_id
column which is added and provides a unique source identifier for the matches.>>> from astropy.coordinates import SkyCoord >>> import astropy.units as u >>> from astropy.table import Table, join_skycoord >>> from astropy import table
>>> sc1 = SkyCoord([0, 1, 1.1, 2], [0, 0, 0, 0], unit='deg') >>> sc2 = SkyCoord([0.5, 1.05, 2.1], [0, 0, 0], unit='deg')
>>> join_func = join_skycoord(0.2 * u.deg) >>> join_func(sc1, sc2) # Associate each coordinate with unique source ID (array([3, 1, 1, 2]), array([4, 1, 2]))
>>> t1 = Table([sc1], names=['sc']) >>> t2 = Table([sc2], names=['sc']) >>> t12 = table.join(t1, t2, join_funcs={'sc': join_skycoord(0.2 * u.deg)}) >>> print(t12) # Note new `sc_id` column with the IDs from join_func() sc_id sc_1 sc_2 deg,deg deg,deg ----- ------- -------- 1 1.0,0.0 1.05,0.0 1 1.1,0.0 1.05,0.0 2 2.0,0.0 2.1,0.0