footprint_healpix#
- m4opt.fov.footprint_healpix(hpx, region, target_coord=None, rotation=None)[source] [edit on github]#
Calculate the HEALPix pixels inside an observing footprint.
The region is expected to represent the field of view of an instrument at a standard orientation of R.A.=0, Dec.=0, P.A.=0. This function rotates the region as if the instrument is pointed at the given target coordinate and optional rotation.
The target coordinate and rotation may be arrays; in that case the return value is a Numpy array of arrays of uneven length.
- Parameters:
hpx (HEALPix) – The HEALPix object specifying the ordering and resolution.
region (Region | Regions) – The shape of the field of view in the standard orientation.
target_coord (SkyCoord | None) – The position for the center of the field of view.
rotation (Annotated[Quantity, PhysicalType('angle')] | None) – The rotation of the field of view about its center.
Examples
First, some imports:
>>> from astropy.coordinates import ICRS, SkyCoord >>> from astropy import units as u >>> from astropy_healpix import HEALPix >>> from m4opt.fov import footprint_healpix >>> import numpy as np >>> from regions import CircleSkyRegion, EllipseSkyRegion, PointSkyRegion, PolygonSkyRegion, RectangleSkyRegion, Regions
We support pointlike FOVs:
>>> hpx = HEALPix(nside=32, frame=ICRS()) >>> region = PointSkyRegion(SkyCoord(0 * u.deg, 0 * u.deg)) >>> target_coord = SkyCoord(5 * u.deg, -5 * u.deg) >>> footprint_healpix(hpx, region, target_coord) array([6593])
We support circular FOVs:
>>> hpx = HEALPix(nside=32, frame=ICRS()) >>> region = CircleSkyRegion(SkyCoord(0 * u.deg, 0 * u.deg), 3 * u.deg) >>> target_coord = SkyCoord(5 * u.deg, -5 * u.deg) >>> footprint_healpix(hpx, region, target_coord) array([6337, 6465, 6466, 6593, 6594, 6721, 6722, 6849, 6850])
and rectangular FOVs:
>>> region2 = RectangleSkyRegion(SkyCoord(0 * u.deg, 0 * u.deg), 6 * u.deg, 8 * u.deg) >>> target_coord = SkyCoord(5 * u.deg, -5 * u.deg) >>> footprint_healpix(hpx, region2, target_coord) array([6209, 6210, 6337, 6338, 6465, 6466, 6593, 6594, 6721, 6722, 6849, 6850, 6977, 6978])
We can compute the footprints for an array of target coordinates:
>>> ras, decs = np.meshgrid([0, 1], [2, 3]) * u.deg >>> target_coords = SkyCoord(ras, decs) >>> footprint_healpix(hpx, region, target_coords) array([[array([5696, 5824, 5951, 5952, 5953, 6079, 6080, 6207]), array([5568, 5696, 5697, 5824, 5951, 5952, 5953, 6080])], [array([5440, 5568, 5695, 5696, 5697, 5823, 5824, 5951, 5952]), array([5568, 5695, 5696, 5697, 5824, 5951, 5952, 5953])]], dtype=object)
We support polygon regions:
>>> region = PolygonSkyRegion(SkyCoord([-2, 2, 0] * u.deg, [0, 0, 2] * u.deg)) >>> footprint_healpix(hpx, region, target_coord) array([6593])
And arrays of target coordinates:
>>> footprint_healpix(hpx, region, target_coords) array([[array([5696, 5824, 5951]), array([5824])], [array([5696]), array([5696])]], dtype=object)
And even non-convex polygon regions:
>>> region = PolygonSkyRegion(SkyCoord([-2, 0, 2, 2, -2] * u.deg, [2, 0, 2, -2, -2] * u.deg)) >>> footprint_healpix(hpx, region, target_coord) array([6593, 6722])
Compound regions are also fine:
>>> regions = Regions([ ... CircleSkyRegion(SkyCoord(0 * u.deg, 0 * u.deg), 3 * u.deg), ... PolygonSkyRegion(SkyCoord([-2, 2, 0] * u.deg, [0, 0, 2] * u.deg))]) >>> footprint_healpix(hpx, regions, target_coord) array([6337, 6465, 6466, 6593, 6594, 6721, 6722, 6849, 6850])
Not all region types are supported:
>>> region = EllipseSkyRegion(SkyCoord(0 * u.deg, 0 * u.deg), 5 * u.deg, 2 * u.deg) >>> footprint_healpix(hpx, region, target_coord) Traceback (most recent call last): ... NotImplementedError: Footprint transformations are not implemented for EllipseSkyRegion