Skip to content

Commit

Permalink
Add HEALPix.from_header
Browse files Browse the repository at this point in the history
Fixes #126.
  • Loading branch information
lpsinger committed May 15, 2019
1 parent 11c89b4 commit 39e166d
Show file tree
Hide file tree
Showing 2 changed files with 59 additions and 4 deletions.
43 changes: 39 additions & 4 deletions astropy_healpix/high_level.py
Original file line number Diff line number Diff line change
Expand Up @@ -7,10 +7,11 @@
from astropy.coordinates.representation import UnitSphericalRepresentation

from .core import (nside_to_pixel_area, nside_to_pixel_resolution,
nside_to_npix, healpix_to_lonlat, lonlat_to_healpix,
bilinear_interpolation_weights, interpolate_bilinear_lonlat,
ring_to_nested, nested_to_ring, healpix_cone_search,
boundaries_lonlat, neighbours)
nside_to_npix, npix_to_nside, healpix_to_lonlat,
lonlat_to_healpix, bilinear_interpolation_weights,
interpolate_bilinear_lonlat, ring_to_nested, nested_to_ring,
healpix_cone_search, boundaries_lonlat, neighbours)
from .utils import parse_input_healpix_data

__all__ = ['HEALPix']

Expand Down Expand Up @@ -52,6 +53,40 @@ def __init__(self, nside=None, order='ring', frame=None):
self.order = order
self.frame = frame

@classmethod
def from_header(cls, input_data, field=0, hdu_in=None, nested=None):
"""
Parameters
----------
input_data : str or `~astropy.io.fits.TableHDU` or `~astropy.io.fits.BinTableHDU` or tuple
The input data to reproject. This can be:
* The name of a HEALPIX FITS file
* A `~astropy.io.fits.TableHDU` or `~astropy.io.fits.BinTableHDU`
instance
* A tuple where the first element is a `~numpy.ndarray` and the
second element is a `~astropy.coordinates.BaseCoordinateFrame`
instance or a string alias for a coordinate frame.
hdu_in : int 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)
nested : bool, 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
-------
healpix : `~astropy_healpix.HEALPix`
A HEALPix pixellization corresponding to the input data.
"""
array_in, frame, nested = parse_input_healpix_data(
input_data, field=field, hdu_in=hdu_in, nested=nested)
nside = npix_to_nside(len(array_in))
order = 'nested' if nested else 'ring'
return cls(nside=nside, order=order, frame=frame)

@property
def pixel_area(self):
"""
Expand Down
20 changes: 20 additions & 0 deletions astropy_healpix/tests/test_high_level.py
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,7 @@

from astropy import units as u
from astropy.coordinates import Longitude, Latitude, Galactic, SkyCoord
from astropy.io.fits import Header

from ..high_level import HEALPix

Expand Down Expand Up @@ -108,6 +109,25 @@ class TestCelestialHEALPix:
def setup_class(self):
self.pix = HEALPix(nside=256, order='nested', frame=Galactic())

def test_healpix_from_header(self):
"""Test instantiation from a FITS header.
Notes
-----
We don't need to test all possible options, because
:meth:`~astropy_healpix.HEALPix.from_header` is just a wrapper around
:meth:`~astropy_healpix.utils.parse_input_healpix_data`, which is
tested exhaustively in :mod:`~astropy_healpix.tests.test_utils`.
"""

pix = HEALPix.from_header(
(np.empty(self.pix.npix), 'G'),
nested=self.pix.order == 'nested')

assert pix.nside == self.pix.nside
assert type(pix.frame) == type(self.pix.frame)
assert pix.order == self.pix.order

def test_healpix_to_skycoord(self):
coord = self.pix.healpix_to_skycoord([1, 2, 3])

Expand Down

0 comments on commit 39e166d

Please sign in to comment.