Multi-Order HEALPix Maps (ligo.skymap.moc
)¶
Support for HEALPix UNIQ pixel indexing [1] and multi-order coverage (MOC) maps [2].
References¶
Reinecke & Hivon, 2015. “Efficient data structures for masks on 2D grids.” AA 580, A132. doi:10.1051/0004-6361/201526549
Boch et al., 2014. “MOC - HEALPix Multi-Order Coverage map.” IVOA Recommendation <http://ivoa.net/documents/MOC/>.
- ligo.skymap.moc.nest2uniq(x1, x2, /, out=None, *, where=True, casting='same_kind', order='K', dtype=None, subok=True[, signature, extobj])¶
Convert a pixel index from NESTED to NUNIQ ordering.
- Parameters:
- order
numpy.ndarray
HEALPix resolution order, the logarithm base 2 of
nside
- ipix
numpy.ndarray
NESTED pixel index
- order
- Returns:
- uniq
numpy.ndarray
NUNIQ pixel index
- uniq
- ligo.skymap.moc.uniq2nest(x, [out1, out2, ]/, [out=(None, None), ]*, where=True, casting='same_kind', order='K', dtype=None, subok=True[, signature, extobj])¶
Convert a pixel index from NUNIQ to NESTED ordering.
- Parameters:
- uniq
numpy.ndarray
NUNIQ pixel index
- uniq
- Returns:
- order
numpy.ndarray
HEALPix resolution order (logarithm base 2 of
nside
)- ipix
numpy.ndarray
NESTED pixel index
- order
- ligo.skymap.moc.uniq2order(x, /, out=None, *, where=True, casting='same_kind', order='K', dtype=None, subok=True[, signature, extobj])¶
Determine the HEALPix resolution order of a HEALPix NESTED index.
- Parameters:
- uniq
numpy.ndarray
NUNIQ pixel index
- uniq
- Returns:
- order
numpy.ndarray
HEALPix resolution order, the logarithm base 2 of
nside
- order
- ligo.skymap.moc.uniq2pixarea(x, /, out=None, *, where=True, casting='same_kind', order='K', dtype=None, subok=True[, signature, extobj])¶
Determine the area of a HEALPix NESTED index.
- Parameters:
- uniq
numpy.ndarray
NUNIQ pixel index
- uniq
- Returns:
- area
numpy.ndarray
The pixel’s area in steradians
- area
- ligo.skymap.moc.rasterize(moc_data, order=None)[source]¶
Convert a multi-order HEALPix dataset to fixed-order NESTED ordering.
- Parameters:
- moc_data
numpy.ndarray
A multi-order HEALPix dataset stored as a Numpy record array whose first column is called UNIQ and contains the NUNIQ pixel index. Every point on the unit sphere must be contained in exactly one pixel in the dataset.
- orderint, optional
The desired output resolution order, or
None
for the maximum resolution present in the dataset.
- moc_data
- Returns:
- nested_data
numpy.ndarray
A fixed-order, NESTED-ordering HEALPix dataset with all of the columns that were in moc_data, with the exception of the UNIQ column.
- nested_data
- ligo.skymap.moc.bayestar_adaptive_grid(probdensity, *args, top_nside=16, rounds=8, **kwargs)[source]¶
Create a sky map by evaluating a function on an adaptive grid.
Perform the BAYESTAR adaptive mesh refinement scheme as described in Section VI of Singer & Price 2016, PRD, 93, 024013 doi:10.1103/PhysRevD.93.024013. This computes the sky map using a provided analytic function and refines the grid, dividing the highest 25% into subpixels and then recalculating their values. The extra given args and kwargs will be passed to the given probdensity function.
- Parameters:
- probdensitycallable
Probability density function. The first argument consists of column-stacked array of right ascension and declination in radians. The return value must be a 1D array of the probability density in inverse steradians with the same length as the argument.
- top_nsideint
HEALPix NSIDE resolution of initial evaluation of the sky map
- roundsint
Number of refinement rounds, including the initial sky map evaluation
- Returns:
- skymapastropy.table.Table
An astropy Table with UNIQ and PROBDENSITY columns, representing a multi-ordered sky map