Contouring (ligo.skymap.postprocess.contour)

ligo.skymap.postprocess.contour.contour(m, levels, nest=False, degrees=False, simplify=True)[source] [edit on github]

Calculate contours from a HEALPix dataset.

Parameters:
mnumpy.ndarray

The HEALPix dataset.

levelslist

The list of contour values.

nestbool, default=False

Indicates whether the input sky map is in nested rather than ring-indexed HEALPix coordinates (default: ring).

degreesbool, default=False

Whether the contours are in degrees instead of radians.

simplifybool, default=True

Whether to simplify the paths.

Returns:
list

A list with the same length as levels. Each item is a list of disjoint polygons, of which each item is a list of points, of which each is a list consisting of the right ascension and declination.

Examples

A very simply example sky map…

>>> nside = 32
>>> npix = ah.nside_to_npix(nside)
>>> ra, dec = hp.pix2ang(nside, np.arange(npix), lonlat=True)
>>> m = dec
>>> contour(m, [10, 20, 30], degrees=True)
[[[[..., ...], ...], ...], ...]
ligo.skymap.postprocess.contour.simplify(vertices, min_area)[source] [edit on github]

Simplify a polygon on the unit sphere.

This is a naive, slow implementation of Visvalingam’s algorithm (see http://bost.ocks.org/mike/simplify/), adapted for for linear rings on a sphere.

Parameters:
verticesnp.ndarray

An Nx3 array of Cartesian vertex coordinates. Each vertex should be a unit vector.

min_areafloat

The minimum area of triangles formed by adjacent triplets of vertices.

Returns:
verticesnp.ndarray