HEALPix

class astropy_healpix.HEALPix(nside=None, order='ring', frame=None)[source]

Bases: object

A HEALPix pixellization.

Parameters
nsideint

Number of pixels along the side of each of the 12 top-level HEALPix tiles

order{ ‘nested’ | ‘ring’ }

Order of HEALPix pixels

frameBaseCoordinateFrame, optional

The celestial coordinate frame of the pixellization. This can be ommitted, in which case the pixellization will not be attached to any particular celestial frame, and the methods ending in _skycoord will not work (but the _lonlat methods will still work and continue to return generic longitudes/latitudes).

Attributes Summary

npix

The number of pixels in the pixellization of the sphere.

pixel_area

The area of a single HEALPix pixel.

pixel_resolution

The resolution of a single HEALPix pixel.

Methods Summary

bilinear_interpolation_weights(self, lon, lat)

Get the four neighbours for each (lon, lat) position and the weight associated with each one for bilinear interpolation.

boundaries_lonlat(self, healpix_index, step)

Return the longitude and latitude of the edges of HEALPix pixels

boundaries_skycoord(self, healpix_index, step)

Return the celestial coordinates of the edges of HEALPix pixels

cone_search_lonlat(self, lon, lat, radius)

Find all the HEALPix pixels within a given radius of a longitude/latitude.

cone_search_skycoord(self, skycoord, radius)

Find all the HEALPix pixels within a given radius of a celestial position.

from_header(input_data[, field, hdu_in, nested])

Parameters

healpix_to_lonlat(self, healpix_index[, dx, dy])

Convert HEALPix indices (optionally with offsets) to longitudes/latitudes

healpix_to_skycoord(self, healpix_index[, …])

Convert HEALPix indices (optionally with offsets) to celestial coordinates.

interpolate_bilinear_lonlat(self, lon, lat, …)

Interpolate values at specific longitudes/latitudes using bilinear interpolation

interpolate_bilinear_skycoord(self, …)

Interpolate values at specific celestial coordinates using bilinear interpolation.

lonlat_to_healpix(self, lon, lat[, …])

Convert longitudes/latitudes to HEALPix indices (optionally with offsets)

neighbours(self, healpix_index)

Find all the HEALPix pixels that are the neighbours of a HEALPix pixel

nested_to_ring(self, nested_index)

Convert a healpix ‘nested’ index to a healpix ‘ring’ index

ring_to_nested(self, ring_index)

Convert a healpix ‘ring’ index to a healpix ‘nested’ index

skycoord_to_healpix(self, skycoord[, …])

Convert celestial coordinates to HEALPix indices (optionally with offsets).

Attributes Documentation

npix

The number of pixels in the pixellization of the sphere.

pixel_area

The area of a single HEALPix pixel.

pixel_resolution

The resolution of a single HEALPix pixel.

Methods Documentation

bilinear_interpolation_weights(self, lon, lat)[source]

Get the four neighbours for each (lon, lat) position and the weight associated with each one for bilinear interpolation.

Parameters
lon, latQuantity

The longitude and latitude values as Quantity instances with angle units.

Returns
indicesndarray

2-D array with shape (4, N) giving the four indices to use for the interpolation

weightsndarray

2-D array with shape (4, N) giving the four weights to use for the interpolation

boundaries_lonlat(self, healpix_index, step)[source]

Return the longitude and latitude of the edges of HEALPix pixels

This returns the longitude and latitude of points along the edge of each HEALPIX pixel. The number of points returned for each pixel is 4 * step, so setting step to 1 returns just the corners.

Parameters
healpix_indexndarray

1-D array of HEALPix pixels

stepint

The number of steps to take along each edge.

Returns
lon, latQuantity

The longitude and latitude, as 2-D arrays where the first dimension is the same as the healpix_index input, and the second dimension has size 4 * step.

boundaries_skycoord(self, healpix_index, step)[source]

Return the celestial coordinates of the edges of HEALPix pixels

This returns the celestial coordinates of points along the edge of each HEALPIX pixel. The number of points returned for each pixel is 4 * step, so setting step to 1 returns just the corners.

This method requires that a celestial frame was specified when initializing HEALPix. If you don’t know or need the celestial frame, you can instead use boundaries_lonlat().

Parameters
healpix_indexndarray

1-D array of HEALPix pixels

stepint

The number of steps to take along each edge.

Returns
skycoordSkyCoord

The celestial coordinates of the HEALPix pixel boundaries

cone_search_lonlat(self, lon, lat, radius)[source]

Find all the HEALPix pixels within a given radius of a longitude/latitude.

Note that this returns all pixels that overlap, including partially, with the search cone. This function can only be used for a single lon/lat pair at a time, since different calls to the function may result in a different number of matches.

Parameters
lon, latQuantity

The longitude and latitude to search around

radiusQuantity

The search radius

Returns
healpix_indexndarray

1-D array with all the matching HEALPix pixel indices.

cone_search_skycoord(self, skycoord, radius)[source]

Find all the HEALPix pixels within a given radius of a celestial position.

Note that this returns all pixels that overlap, including partially, with the search cone. This function can only be used for a single celestial position at a time, since different calls to the function may result in a different number of matches.

This method requires that a celestial frame was specified when initializing HEALPix. If you don’t know or need the celestial frame, you can instead use cone_search_lonlat().

Parameters
skycoordSkyCoord

The celestial coordinates to use for the cone search

radiusQuantity

The search radius

Returns
healpix_indexndarray

1-D array with all the matching HEALPix pixel indices.

classmethod from_header(input_data, field=0, hdu_in=None, nested=None)[source]
Parameters
input_datastr or TableHDU or BinTableHDU or tuple

The input data to reproject. This can be:

hdu_inint or str, optional

If input_data is a FITS file, specifies the HDU to use. (the default HDU for HEALPIX data is 1, unlike with image files where it is generally 0)

nestedbool, optional

The order of the healpix_data, either nested (True) or ring (False). If a FITS file is passed in, this is determined from the header.

Returns
healpixHEALPix

A HEALPix pixellization corresponding to the input data.

healpix_to_lonlat(self, healpix_index, dx=None, dy=None)[source]

Convert HEALPix indices (optionally with offsets) to longitudes/latitudes

Parameters
healpix_indexndarray

1-D array of HEALPix indices

dx, dyndarray, optional

1-D arrays of offsets inside the HEALPix pixel, which must be in the range [0:1] (0.5 is the center of the HEALPix pixels). If not specified, the position at the center of the pixel is used.

Returns
lonLongitude

The longitude values

latLatitude

The latitude values

healpix_to_skycoord(self, healpix_index, dx=None, dy=None)[source]

Convert HEALPix indices (optionally with offsets) to celestial coordinates.

Note that this method requires that a celestial frame was specified when initializing HEALPix. If you don’t know or need the celestial frame, you can instead use healpix_to_lonlat().

Parameters
healpix_indexndarray

1-D array of HEALPix indices

dx, dyndarray, optional

1-D arrays of offsets inside the HEALPix pixel, which must be in the range [0:1] (0.5 is the center of the HEALPix pixels). If not specified, the position at the center of the pixel is used.

Returns
coordSkyCoord

The resulting celestial coordinates

interpolate_bilinear_lonlat(self, lon, lat, values)[source]

Interpolate values at specific longitudes/latitudes using bilinear interpolation

If a position does not have four neighbours, this currently returns NaN.

Parameters
lon, latQuantity

The longitude and latitude values as Quantity instances with angle units.

valuesndarray

1-D array with the values in each HEALPix pixel. This must have a length of the form 12 * nside ** 2 (and nside is determined automatically from this).

Returns
resultndarray

1-D array of interpolated values

interpolate_bilinear_skycoord(self, skycoord, values)[source]

Interpolate values at specific celestial coordinates using bilinear interpolation.

If a position does not have four neighbours, this currently returns NaN.

Note that this method requires that a celestial frame was specified when initializing HEALPix. If you don’t know or need the celestial frame, you can instead use interpolate_bilinear_lonlat().

Parameters
skycoordSkyCoord

The celestial coordinates at which to interpolate

valuesndarray

1-D array with the values in each HEALPix pixel. This must have a length of the form 12 * nside ** 2 (and nside is determined automatically from this).

Returns
resultndarray

1-D array of interpolated values

lonlat_to_healpix(self, lon, lat, return_offsets=False)[source]

Convert longitudes/latitudes to HEALPix indices (optionally with offsets)

Parameters
lon, latQuantity

The longitude and latitude values as Quantity instances with angle units.

return_offsetsbool

If True, the returned values are the HEALPix pixel as well as dx and dy, the fractional positions inside the pixel. If False (the default), only the HEALPix pixel is returned.

Returns
healpix_indexndarray

1-D array of HEALPix indices

dx, dyndarray

1-D arrays of offsets inside the HEALPix pixel in the range [0:1] (0.5 is the center of the HEALPix pixels). This is returned if return_offsets is True.

neighbours(self, healpix_index)[source]

Find all the HEALPix pixels that are the neighbours of a HEALPix pixel

Parameters
healpix_indexndarray

Array of HEALPix pixels

Returns
neighndarray

Array giving the neighbours starting SW and rotating clockwise. This has one extra dimension compared to healpix_index - the first dimension - which is set to 8. For example if healpix_index has shape (2, 3), neigh has shape (8, 2, 3).

nested_to_ring(self, nested_index)[source]

Convert a healpix ‘nested’ index to a healpix ‘ring’ index

Parameters
nested_indexndarray

Healpix index using the ‘nested’ ordering

Returns
ring_indexndarray

Healpix index using the ‘ring’ ordering

ring_to_nested(self, ring_index)[source]

Convert a healpix ‘ring’ index to a healpix ‘nested’ index

Parameters
ring_indexndarray

Healpix index using the ‘ring’ ordering

Returns
nested_indexndarray

Healpix index using the ‘nested’ ordering

skycoord_to_healpix(self, skycoord, return_offsets=False)[source]

Convert celestial coordinates to HEALPix indices (optionally with offsets).

Note that this method requires that a celestial frame was specified when initializing HEALPix. If you don’t know or need the celestial frame, you can instead use lonlat_to_healpix().

Parameters
skycoordSkyCoord

The celestial coordinates to convert

return_offsetsbool

If True, the returned values are the HEALPix pixel as well as dx and dy, the fractional positions inside the pixel. If False (the default), only the HEALPix pixel is returned.

Returns
healpix_indexndarray

1-D array of HEALPix indices

dx, dyndarray

1-D arrays of offsets inside the HEALPix pixel in the range [0:1] (0.5 is the center of the HEALPix pixels). This is returned if return_offsets is True.