Cross Match Catalogs with HEALPix Sky Maps (ligo.skymap.postprocess.crossmatch)

Catalog cross matching for HEALPix sky maps.

class ligo.skymap.postprocess.crossmatch.CrossmatchResult(searched_area, searched_prob, offset, searched_modes, contour_areas, area_probs, contour_modes, searched_prob_dist, contour_dists, searched_vol, searched_prob_vol, contour_vols, probdensity, probdensity_vol)

Bases: tuple

Cross match result as returned by crossmatch().

Notes

  • All probabilities returned are between 0 and 1.

  • All angles returned are in degrees.

  • All areas returned are in square degrees.

  • All distances are luminosity distances in units of Mpc.

  • All volumes are in units of Mpc³. If crossmatch() was run with cosmology=False, then all volumes are Euclidean volumes in luminosity distance. If crossmatch() was run with cosmology=True, then all volumes are comoving volumes.

area_probs

Probability within the 2D credible regions of the areas specified by the areas argument passed to crossmatch().

contour_areas

Area within the 2D credible regions of the probabilities specified by the contour argument passed to crossmatch().

contour_dists

Distance credible interval, marginalized over sky position, Same length as the coordinates argument passed to crossmatch().

contour_modes

Number of disconnected regions within the 2D credible regions of the probabilities specified by the contour argument passed to crossmatch().

contour_vols

Volume within the 3D credible regions of the probabilities specified by the contour argument passed to crossmatch().

offset

Angles on the sky between the target positions and the maximum a posteriori position. Same length as the coordinates argument passed to crossmatch().

probdensity

2D probability density per steradian at the positions of each of the targets. Same length as the coordinates argument passed to crossmatch().

probdensity_vol

3D probability density per cubic megaparsec at the positions of each of the targets. Same length as the coordinates argument passed to crossmatch().

searched_area

Area within the 2D credible region containing each target position. Same length as the coordinates argument passed to crossmatch().

searched_modes

Number of disconnected regions within the 2D credible regions containing each target position. Same length as the coordinates argument passed to crossmatch().

searched_prob

Probability within the 2D credible region containing each target position. Same length as the coordinates argument passed to crossmatch().

searched_prob_dist

Cumulative CDF of distance, marginalized over sky position, at the distance of each of the targets. Same length as the coordinates argument passed to crossmatch().

searched_prob_vol

Probability within the 3D credible region containing each target position. Same length as the coordinates argument passed to crossmatch().

searched_vol

Volume within the 3D credible region containing each target position. Same length as the coordinates argument passed to crossmatch().

ligo.skymap.postprocess.crossmatch.crossmatch(sky_map, coordinates=None, contours=(), areas=(), modes=False, cosmology=False)[source] [edit on github]

Cross match a sky map with a catalog of points.

Given a sky map and the true right ascension and declination (in radians), find the smallest area in deg^2 that would have to be searched to find the source, the smallest posterior mass, and the angular offset in degrees from the true location to the maximum (mode) of the posterior. Optionally, also compute the areas of and numbers of modes within the smallest contours containing a given total probability.

Parameters:
sky_mapastropy.table.Table

A multiresolution sky map, as returned by ligo.skymap.io.fits.read_sky_map() called with the keyword argument moc=True.

coordinatesastropy.coordinates.SkyCoord, optional

The catalog of target positions to match against.

contourstuple, optional

Credible levels between 0 and 1. If this argument is present, then calculate the areas and volumes of the 2D and 3D credible regions that contain these probabilities. For example, for contours=(0.5, 0.9), then areas and volumes of the 50% and 90% credible regions.

areastuple, optional

Credible areas in square degrees. If this argument is present, then calculate the probability contained in the 2D credible levels that have these areas. For example, for areas=(20, 100), then compute the probability within the smallest credible levels of 20 deg² and 100 deg², respectively.

modesbool, optional

If True, then enable calculation of the number of distinct modes or islands of probability. Note that this option may be computationally expensive.

cosmologybool, optional

If True, then search space by descending probability density per unit comoving volume. If False, then search space by descending probability per luminosity distance cubed.

Returns:
resultCrossmatchResult

Notes

This function is also be used for injection finding; see Gather Summary Statistics (ligo-skymap-stats).

Examples

First, some imports:

>>> from astroquery.vizier import VizierClass
>>> from astropy.coordinates import SkyCoord
>>> from ligo.skymap.io import read_sky_map
>>> from ligo.skymap.postprocess import crossmatch

Next, retrieve the GLADE catalog using Astroquery and get the coordinates of all its entries:

>>> vizier = VizierClass(
...     row_limit=-1,
...     columns=['recno', 'GWGC', '_RAJ2000', '_DEJ2000', 'Dist'])
>>> cat, = vizier.get_catalogs('VII/281/glade2')
>>> cat.sort('recno')  # sort catalog so that doctest output is stable
>>> del cat['recno']
>>> coordinates = SkyCoord(cat['_RAJ2000'], cat['_DEJ2000'], cat['Dist'])

Load the multiresolution sky map for S190814bv:

>>> url = 'https://gracedb.ligo.org/api/superevents/S190814bv/files/bayestar.multiorder.fits'
>>> skymap = read_sky_map(url, moc=True)

Perform the cross match:

>>> result = crossmatch(skymap, coordinates)

Using the cross match results, we can list the galaxies within the 90% credible volume:

>>> print(cat[result.searched_prob_vol < 0.9])
     _RAJ2000          _DEJ2000        GWGC            Dist
       deg               deg                           Mpc
----------------- ----------------- ---------- --------------------
  9.3396700000000 -19.9342460000000    NGC0171    57.56212553960000
 20.2009090000000 -31.1146050000000        ---   137.16022925600001
  8.9144680000000 -20.1252980000000 ESO540-003    49.07809291930000
 10.6762720000000 -21.7740820000000        ---   276.46938505499998
 13.5855170000000 -23.5523850000000        ---   138.44550704800000
 20.6362970000000 -29.9825150000000        ---   160.23313164900000
 13.1923880000000 -22.9750180000000        ---   236.96795954500001
 11.7813630000000 -24.3706470000000        ---   244.25031189699999
 19.1711120000000 -31.4339490000000        ---   152.13614001400001
 13.6367060000000 -23.4948790000000        ---   141.25162979500001
              ...               ...        ...                  ...
 11.3517000000000 -25.8597000000000        ---   335.73800000000000
 11.2074000000000 -25.7149000000000        ---   309.02999999999997
 11.1875000000000 -25.7504000000000        ---   295.12099999999998
 10.8609000000000 -25.6904000000000        ---   291.07200000000000
 10.6939000000000 -25.6778300000000        ---   323.59399999999999
 15.4935000000000 -26.0305000000000        ---   304.78899999999999
 15.2794000000000 -27.0411000000000        ---   320.62700000000001
 14.8324000000000 -27.0460000000000        ---   320.62700000000001
 14.5341000000000 -26.0949000000000        ---   307.61000000000001
 23.1281000000000 -31.1109200000000        ---   320.62700000000001
Length = 1479 rows