llama.files.healpix.skymap module

Define class containing a HEALPix skymap and convenience methods for manipulating, creating, and reformating skymaps.

class llama.files.healpix.skymap.HEALPixPSF

Bases: llama.files.healpix.skymap.HEALPixRepresentableFileHandler

A class with utility functions for generating skymaps in-memory that follow the HEALPix parameters (NSIDE, ORDERING) of the file described by a template HEALPixSkyMapFileHandler. Used to dynamically generate skymaps that are compatible with the template skymap based on some PSF (Gaussian by default). Implementations of this class should contain point-like events with some sort of point-spread function that can be used to generate skymaps.

get_healpix(nest=None, template_skymap=None)

Generate a list of HEALPixSkyMap instances from the events contained in this FileHandler.

Parameters
  • nest (bool, optional) – If True, return skymaps with NESTED ordering. If False, return skymaps with RING ordering. Overrides the ordering (but not NSIDE) defined in template_skymap_filehandler.

  • template_skymap (HEALPixSkyMap, optional) – Provide a HEALPixSkyMap as a template to override the ones provided by template_skymap_filehandler. Also overrides nest (if provided).

Returns

skymaps – A list of HEALPixSkyMap instances of the triggers in this FileHandler that conform to the HEALPixSkyMap parameters in this class’s template_skymap_filehandler

Return type

list

abstract source_locations()

Returns a list of tuples of (RA, Dec, sigma) for all point-like sources, where (RA, Dec) is the central sky position and sigma is the standard deviation of the Gaussian PSF. Since this depends on the structure of data in the relevant file handler, it must be specified for each subclass.

abstract property template_skymap_filehandler

The HEALPixSkyMapFileHandler whose HEALPix parameters should be used as a template for the HEALPixSkyMap output by this FileHandler. It is assumed that this HEALPixSkyMapFileHandler only returns one skymap in its list; consequently, only the first skymap loaded by this file handler will be used. Note that the template skymap file handler is likely not a dependency of this file handler; it is only used as a default for formatting generated skymaps from this file handler’s data.

class llama.files.healpix.skymap.HEALPixRepresentableFileHandler

Bases: llama.filehandler.TriggerList

A FileHandler class containing triggers that can be represented in HEALPix. Provides get_healpix and get_healpix_descriptions methods for this purpose. Implementations need not store data in a HEALPix format natively, as long as they provide a get_healpix method for returning valid HEALPixSkyMap instances.

abstract get_healpix(nest=None, template_skymap=None)

Return a list of HEALPix skymaps representing the triggers in this FileHandler.

get_healpix_descriptions(skymap_indices=None)

Get a list of descriptive strings for each HEALPix skymap contained in this file.

Parameters

skymap_indices (list) – A list of tuples describing the index of each skymap used to produce this output skymap, with each entry in the tuple corresponding to a detector. In this case, the descriptive strings will be generated based on this FileHandler instance’s detector property and the indices of the skymap_indices in skymap_indices. The ordering of the returned descriptive strings will match the ordering of skymap_indices. Note that for single-detector skymaps, the elements of the skymap_indices list should just be single-element tuples corresponding to the index of each skymap.

Returns

descriptions – Descriptions of each HEALPix skymap in the order in which they appear in the file.

Return type

list

class llama.files.healpix.skymap.HEALPixSkyMap(pixels, nest=True, title=None)

Bases: llama.files.healpix.skymap.SkyMap, llama.files.healpix.skymap.MollviewMixin

A HEALPix skymap object residing in memory, including methods for working with the skymap and combining skymaps with each other.

ang2pix(ra, dec, degrees=True)

Return the pixel index on this skymap at right ascensions ra and declinations dec.

conform_skymaps(*args)

Make sure all skymaps contained in *args have the same ordering (NESTED or RING) as the current skymap (self), reordering the pixels if necessary, and returning a tuple of the (possibly reordered) input skymaps.

dangle(ra, dec, degrees=True)

Get a new skymap with the same pixelization where each pixel value is the distance from that pixel to the point with right ascension ra and declination dec. By default, ra and dec are in degrees, though this can be overridden by setting degrees=False, in which case they will be interpretted as radians. The measured angle differences in the returned skymap will always be in degrees. Dot product formula per-pixel is:

sin(DEC1)*sin(DEC2) + cos(DEC1)*cos(DEC2)*cos(RA2-RA1)

So we get the angle by taking the arccos of this function.

Examples

We should find that the distance from any pixel to the North pole is equal to 90 minus the declination (within some small error):

>>> skymap = HEALPixSkyMap([1]*12)
>>> skymap.dangle(32, 90).pixels + skymap.dec - 90 < 1e-12
array([ True,  True,  True,  True,  True,  True,  True,  True,  True,
        True,  True,  True])

Likewise, the distance from any pixel to the South pole should be equal to 90 plus the declination:

>>> skymap.dangle(359, -90).pixels - skymap.dec - 90. < 1e-12
array([ True,  True,  True,  True,  True,  True,  True,  True,  True,
        True,  True,  True])
property dec

Declinations in degrees.

mollview(**kwargs)

Return a Mollweide projection of this skymap as a matplotlib figure with the correct nesting. Remaining **kwargs are passed to healpy.mollview. Use the LIGO standard rotation where 180 degrees is at the center of the map.

property nside

The NSIDE HEALPix parameter indicating how many times each base pixel has been subdivided in each of the 2 dimensions on the sphere.

classmethod nside2ang(nside, nest=True, degrees=True)

Return (RA, Dec) in degrees for the pixels corresponding to a HEALPix skymap with NSIDE=nside where nest=True (DEFAULT) indicates NEST ordering and nest=False indicates RING ordering. If degrees=True (DEFAULT), the return value will be in degrees; otherwise, radians.

Examples

>>> import numpy as np
>>> ra, dec = HEALPixSkyMap.nside2ang(1, nest=True, degrees=True)
>>> all(1e-13 > ra - np.array([
...     45., 135., 225., 315.,   0.,  90., 180., 270.,  45., 135.,
...     225., 315.
... ]))
True
psf_gaussian(ra, dec, sigma, degrees=True, normalize=True)

Calculate a Gaussian PSF at all pixels in this skymap.

Parameters
  • ra (float) – Right Ascension of the Gaussian mean.

  • dec (float) – Declination of the Gaussian mean.

  • sigma (float) – Standard deviation of the Gaussian PSF.

  • degrees (bool, optional) – If true, input arguments (ra, dec, sigma) are assumed to be in degrees. Otherwise, they are assumed to be in radians.

  • normalize (bool, optional) – If true, normalize the returned skymap so that the integral over area is 1 (useful for when the Gaussian is meant to represent a probability density function).

Returns

skymap – A new skymap with the same pixelization where a 2-D Gaussian PSF (point spread function) centered at (ra, dec) with standard deviation sigma has been calculated for each pixel in this skymap.

Return type

HEALPixSkyMap

property ra

Right Ascensions in degrees.

to_nest()

Return a skymap with the same pixel values and NESTED ordering.

to_ring()

Return a skymap with the same pixel values and RING ordering.

class llama.files.healpix.skymap.MollviewMixin

Bases: abc.ABC

Defines an interface for plotting Mollweide projections.

abstract mollview(**kwargs)

Return a Mollweide projection of this skymap as a matplotlib figure.

class llama.files.healpix.skymap.SkyMap(pixels, ra, dec, degrees=True, area_per_pixel=None)

Bases: abc.ABC

An array representing a skymap.

property area_per_pixel

Get the solid angle per pixel. Can be set to a scalar value or a vector with the same length as the list of pixels; in either case, the pixel areas are returned as a vector for each pixel. If not set, assumes that all pixel areas are identical and that the pixel areas sum to a full sphere.

property as_pdf

Return a skymap whose pixels represent proability/pixel (rather than probability per solid angle, i.e. a proper PDF) and putting them into units of probability per solid angle. Returns a new skymap with units of inverse solid angle. If this skymap is already a PDF in units of inverse solid angle, this operation returns the original skymap.

confidence_region_size(percent=90)

Calculate the size of the confidence region in square degrees. By default, finds the 90% confidence region, though the percentage can be manually specified. Does NOT normalize pixels before summing them (in case the total probability does not add to 1, e.g. if we are already dealing with a partial skymap).

abstract dangle(ra, dec, degrees=True)

Get a new skymap with the same pixelization where each pixel value is the distance from that pixel to the point with right ascension ra and declination dec. By default, ra and dec are in degrees, though this can be overridden by setting degrees=False, in which case they will be interpretted as radians. The measured angle differences in the returned skymap will always be in degrees. Dot product formula per-pixel is:

sin(DEC1)*sin(DEC2) + cos(DEC1)*cos(DEC2)*cos(RA2-RA1)

property dec

Declinations in degrees.

integrate()

If pixels is given in units of inverse solid angle, return the integral of the pixel values over the whole sky. If pixels is dimensionless, assume the pixel values are probability per-pixel and simply return the sum of pixel values. Either way, the integral is dimensionless and equals the total probability of the SkyMap (given the stated assumptions).

map(function)

Map a function to the pixel locations of this skymap and return an identically sized and ordered skymap with the output values of that function as its pixel values. function must take vector arguments ra and dec and return a list of pixel values of the same length. No normalization takes place.

normalize(norm_to=1)

Scale the pixel values of this skymap so that they integrate to the value specified by norm_to (DEFAULT: 1). Returns a new skymap with the same units as self.pixels (assuming self.integrate is defined for self.units used in self.pixels.

property npix
property ra

Right Ascensions in degrees.

property unit

Get the units of the pixels of this skymap. If pixels is already an astropy.Quantity, simply return its unit. Otherwise, return astropy.units.Unit("") (i.e. dimensionless).

llama.files.healpix.skymap.psf_gaussian(delta_angle, sigma, degrees=True)

Return a 2-D Gaussian PDF function to be used as a PSF (point spread function) on the sphere. Note that this Gaussian PDF does not take account of curvature and is accurate only for small angular spreads.

Parameters
  • delta_angle (array_like) – A scalar or array of angular distances from the center of the Gaussian. Must be in same units (degrees or radians) as sigma.

  • sigma (float) – The standard deviation of the Gaussian. Must be in same units (degrees or radians) as delta_angle.

  • degrees (bool, optional) – Whether angles are given in degrees. If False, units are assumed to be in radians.

Returns

  • densities (array_like) – The values of the Gaussian at the distances given by delta_angle in units of [1/steradian] (even if inputs are given in degrees!).

  • >>> import numpy as np

  • >>> psf_gaussian(1, 1, degrees=False) == 1/(2*np.pi) * np.exp(-1/2)

  • True

llama.files.healpix.skymap.psf_gaussian_scales(sigma)

Get the characteristic maximum and minimum scales of a Gaussian PSF for a given angle sigma.

Parameters

sigma (int or float) – The angular standard deviation of the Gaussian PSF in degrees or radians.

Returns

  • min_scale (float) – The scale of inter-pixel distances required to capture the features of a Gaussian with the given standard deviation. In same units as sigma.

  • max_scale (float) – The diameter of a pixel mesh centered on the Gaussian mean that is large enough to capture almost all of the probability in the distribution. In same units as sigma.