HpxGeom#
- class gammapy.maps.HpxGeom[source]#
Bases:
Geom
Geometry class for HEALPix maps.
This class performs mapping between partial-sky indices (pixel number within a HEALPix region) and all-sky indices (pixel number within an all-sky HEALPix map). Multi-band HEALPix geometries use a global indexing scheme that assigns a unique pixel number based on the all-sky index and band index. In the single-band case the global index is the same as the HEALPix index.
By default, the constructor will return an all-sky map. Partial-sky maps can be defined with the
region
argument.- Parameters:
- nside
ndarray
HEALPix NSIDE parameter, the total number of pixels is 12*nside*nside. For multi-dimensional maps one can pass either a single
nside
value or a vector ofnside
values defining the pixel size for each image plane. Ifnside
is not a scalar then its dimensionality should match that of the non-spatial axes. If nest is True,nside
must be a power of 2, less than 2**30.- nestbool, optional
Indexing scheme. If True, “NESTED” scheme. If False, “RING” scheme. Default is True.
- frame{“icrs”, “galactic”}
Coordinate system. Default is “icrs”.
- regionstr or tuple, optional
Spatial geometry for partial-sky maps. If None, the map will encompass the whole sky. String input will be parsed according to HPX_REG header keyword conventions. Tuple input can be used to define an explicit list of pixels encompassed by the geometry. Default is None.
- axeslist
Axes for non-spatial dimensions.
- nside
Notes
For examples, see HEALPix-based maps
Attributes Summary
List of non-spatial axes.
All axes names.
Map coordinates of the center of the geometry as a tuple.
Pixel coordinates of the center of the geometry as a tuple.
Sky coordinate of the center of the geometry.
Shape of the
ndarray
matching this geometry.Shape of data of the non-spatial axes and unit spatial axes.
HEALPix pixel and band indices for every pixel in the map.
Flag for all-sky maps.
Flag identifying whether this geometry is regular in non-spatial dimensions.
Number of dimensions as an integer.
Whether HEALPix order is nested as a boolean.
Number of pixels in each band.
Maximum number of pixels.
NSIDE in each band.
The order in each band (
NSIDE = 2 ** ORDER
).HEALPix ordering ('NESTED' or 'RING').
Map projection.
Region string.
Shape of non-spatial axes.
Width of the map.
Methods Summary
contains
(coords)Check if a given map coordinate is contained in the geometry.
coord_to_pix
(coords)Convert map coordinates to pixel coordinates.
create
([nside, binsz, nest, frame, region, ...])Create an HpxGeom object.
crop
(crop_width)Crop the geometry at the edges.
cutout
(position, width, **kwargs)Create a cutout around a given position.
downsample
(factor[, axis_name])Downsample the spatial dimension of the geometry by a given factor.
from_hdu
(hdu[, hdu_bands])Create an HPX object from a BinTable HDU.
from_header
(header[, hdu_bands, format])Create an HPX object from a FITS header.
get_coord
([idx, flat, sparse, mode, axis_name])Get the coordinate array for this geometry.
get_idx
([idx, local, flat, sparse, mode, ...])Get tuple of pixel indices for this geometry.
get_index_list
(nside, nest, region)Get list of pixels indices for all the pixels in a region.
global_to_local
(idx_global[, ravel])Compute global (all-sky) index from a local (partial-sky) index.
interp_weights
(coords[, idxs])Get interpolation weights for given coordinates.
is_aligned
(other)Check if HEALPix geoms and extra axes are aligned.
is_allclose
(other[, rtol_axes, atol_axes])Compare two data IRFs for equivalency.
local_to_global
(idx_local)Compute a global index (all-sky) from a local (partial-sky) index.
pix_to_coord
(pix)Convert pixel coordinates to map coordinates.
pix_to_idx
(pix[, clip])Convert pixel coordinates to pixel indices.
region_mask
(regions)Create a mask from a given list of regions.
separation
(center)Compute sky separation with respect to a given center.
Solid angle array as a
Quantity
insr
.to_binsz
(binsz)Change pixel size of the geometry.
to_cube
(axes)Append non-spatial axes to create a higher-dimensional geometry.
to_header
([format])Build and return FITS header for this HEALPix map.
to_image
()Create a 2D image geometry (drop non-spatial dimensions).
to_nside
(nside)Upgrade or downgrade the resolution to a given NSIDE.
Geometry copy with swapped ORDERING (NEST->RING or vice versa).
to_wcs_geom
([proj, oversample, width_pix])Make a WCS projection appropriate for this HEALPix pixelization.
to_wcs_tiles
([nside_tiles, margin])Create WCS tiles geometries from HPX geometry with given nside.
upsample
(factor)Upsample the spatial dimension of the geometry by a given factor.
Attributes Documentation
- axes#
List of non-spatial axes.
- axes_names#
All axes names.
- center_coord#
Map coordinates of the center of the geometry as a tuple.
- center_pix#
Pixel coordinates of the center of the geometry as a tuple.
- center_skydir#
Sky coordinate of the center of the geometry.
- Returns:
- center
SkyCoord
Center position.
- center
- data_shape_axes#
Shape of data of the non-spatial axes and unit spatial axes.
- frame#
- ipix#
HEALPix pixel and band indices for every pixel in the map.
- is_allsky#
Flag for all-sky maps.
- is_hpx = True#
- is_region = False#
- is_regular#
Flag identifying whether this geometry is regular in non-spatial dimensions.
False for multi-resolution or irregular geometries. If True, all image planes have the same pixel geometry.
- ndim#
Number of dimensions as an integer.
- nest#
Whether HEALPix order is nested as a boolean.
- npix#
Number of pixels in each band.
For partial-sky geometries this can be less than the number of pixels for the band NSIDE.
- npix_max#
Maximum number of pixels.
- nside#
NSIDE in each band.
- order#
The order in each band (
NSIDE = 2 ** ORDER
).Set to -1 for bands with NSIDE that is not a power of 2.
- ordering#
HEALPix ordering (‘NESTED’ or ‘RING’).
- pixel_scales#
- projection#
Map projection.
- region#
Region string.
- shape_axes#
Shape of non-spatial axes.
- width#
Width of the map.
Methods Documentation
- coord_to_pix(coords)[source]#
Convert map coordinates to pixel coordinates.
- Parameters:
- coordstuple
Coordinate values in each dimension of the map. This can either be a tuple of numpy arrays or a MapCoord object. If passed as a tuple then the ordering should be (longitude, latitude, c_0, …, c_N) where c_i is the coordinate vector for axis i.
- Returns:
- pixtuple
Tuple of pixel coordinates in image and band dimensions.
- classmethod create(nside=None, binsz=None, nest=True, frame='icrs', region=None, axes=None, skydir=None, width=None)[source]#
Create an HpxGeom object.
- Parameters:
- nsideint or
ndarray
, optional HEALPix NSIDE parameter. This parameter sets the size of the spatial pixels in the map. If nest is True,
nside
must be a power of 2, less than 2**30. Default is None.- binszfloat or
ndarray
, optional Approximate pixel size in degrees. An
nside
will be chosen that corresponds to a pixel size closest to this value. This option is superseded bynside
. Default is None.- nestbool, optional
Indexing scheme. If True, “NESTED” scheme. If False, “RING” scheme. Default is True.
- frame{“icrs”, “galactic”}
Coordinate system, either Galactic (“galactic”) or Equatorial (“icrs”). Default is “icrs”.
- regionstr, optional
HEALPix region string. Allows for partial-sky maps. Default is None.
- axeslist, optional
List of axes for non-spatial dimensions. Default is None.
- skydirtuple or
SkyCoord
, optional Sky position of map center. Can be either a SkyCoord object or a tuple of longitude and latitude in deg in the coordinate system of the map. Default is None.
- widthfloat, optional
Diameter of the map in degrees. If set the map will encompass all pixels within a circular region centered on
skydir
. Default is None.
- nsideint or
- Returns:
- geom
HpxGeom
A HEALPix geometry object.
- geom
Examples
>>> from gammapy.maps import HpxGeom, MapAxis >>> axis = MapAxis.from_bounds(0,1,2) >>> geom = HpxGeom.create(nside=16) >>> geom = HpxGeom.create(binsz=0.1, width=10.0) >>> geom = HpxGeom.create(nside=64, width=10.0, axes=[axis]) >>> geom = HpxGeom.create(nside=[32,64], width=10.0, axes=[axis])
- crop(crop_width)[source]#
Crop the geometry at the edges.
- Parameters:
- crop_width{sequence, array_like, int}
Number of values cropped from the edges of each axis.
- Returns:
- geom
Geom
Cropped geometry.
- geom
- downsample(factor, axis_name=None)[source]#
Downsample the spatial dimension of the geometry by a given factor.
- Parameters:
- factorint
Downsampling factor.
- axis_namestr
Axis to downsample.
- Returns:
- geom
Geom
Downsampled geometry.
- geom
- classmethod from_hdu(hdu, hdu_bands=None)[source]#
Create an HPX object from a BinTable HDU.
- Parameters:
- hdu
BinTableHDU
The FITS HDU.
- hdu_bands
BinTableHDU
, optional The BANDS table HDU. Default is None.
- hdu
- Returns:
- hpx
HpxGeom
HEALPix geometry.
- hpx
- classmethod from_header(header, hdu_bands=None, format=None)[source]#
Create an HPX object from a FITS header.
- Parameters:
- header
Header
The FITS header.
- hdu_bands
BinTableHDU
, optional The BANDS table HDU. Default is None.
- formatstr, optional
FITS convention. Default is None. If None the format is guessed. The following formats are supported:
“gadf”
“fgst-ccube”
“fgst-ltcube”
“fgst-bexpcube”
“fgst-srcmap”
“fgst-template”
“fgst-srcmap-sparse”
“galprop”
“galprop2”
“healpy”
- header
- Returns:
- hpx
HpxGeom
HEALPix geometry.
- hpx
- get_coord(idx=None, flat=False, sparse=False, mode='center', axis_name=None)[source]#
Get the coordinate array for this geometry.
Returns a coordinate array with the same shape as the data array. Pixels outside the geometry are set to NaN. Coordinates for a single image plane can be accessed by setting
idx
to the index tuple of a plane.- Parameters:
- idxtuple, optional
A tuple of indices with one index for each non-spatial dimension. If defined only coordinates for the image plane with this index will be returned. If none then coordinates for all pixels will be returned. Default is None.
- flatbool, optional
Return a flattened array containing only coordinates for pixels contained in the geometry. Default is False.
- Returns:
- coordstuple
Tuple of coordinate vectors with one vector for each dimension.
- get_idx(idx=None, local=False, flat=False, sparse=False, mode='center', axis_name=None)[source]#
Get tuple of pixel indices for this geometry.
Returns all pixels in the geometry by default. Pixel indices for a single image plane can be accessed by setting
idx
to the index tuple of a plane.- Parameters:
- idxtuple, optional
A tuple of indices with one index for each non-spatial dimension. If defined only pixels for the image plane with this index will be returned. If none then all pixels will be returned. Default is None.
- localbool, optional
Flag to return local or global pixel indices. Local indices run from 0 to the number of pixels in a given image plane. Default is False.
- flatbool, optional
Return a flattened array containing only indices for pixels contained in the geometry. Default is False.
- Returns:
- idxtuple
Tuple of pixel index vectors with one vector for each dimension.
- static get_index_list(nside, nest, region)[source]#
Get list of pixels indices for all the pixels in a region.
- Parameters:
- nsideint
HEALPix NSIDE parameter.
- nestbool
Indexing scheme. If True, “NESTED” scheme. If False, “RING” scheme.
- regionstr
HEALPix region string.
- Returns:
- ilist
ndarray
List of pixel indices.
- ilist
- global_to_local(idx_global, ravel=False)[source]#
Compute global (all-sky) index from a local (partial-sky) index.
- Parameters:
- idx_globaltuple
A tuple of pixel indices with global HEALPix pixel indices.
- ravelbool, optional
Return a raveled index. Default is False.
- Returns:
- idx_localtuple
A tuple of pixel indices with local HEALPix pixel indices.
- is_aligned(other)[source]#
Check if HEALPix geoms and extra axes are aligned.
- Parameters:
- other
HpxGeom
Other geometry.
- other
- Returns:
- alignedbool
Whether geometries are aligned.
- is_allclose(other, rtol_axes=1e-06, atol_axes=1e-06)[source]#
Compare two data IRFs for equivalency.
- Parameters:
- other
HpxGeom
Geometry to compare against.
- rtol_axesfloat, optional
Relative tolerance for axes comparison. Default is 1e-6.
- atol_axesfloat, optional
Relative tolerance for axes comparison. Default is 1e-6.
- other
- Returns:
- is_allclosebool
Whether the geometry is all close.
- local_to_global(idx_local)[source]#
Compute a global index (all-sky) from a local (partial-sky) index.
- Parameters:
- idx_localtuple
A tuple of pixel indices with local HEALPix pixel indices.
- Returns:
- idx_globaltuple
A tuple of pixel index vectors with global HEALPix pixel indices.
- pix_to_coord(pix)[source]#
Convert pixel coordinates to map coordinates.
- Parameters:
- pixtuple
Tuple of pixel coordinates.
- Returns:
- coordstuple
Tuple of map coordinates.
- pix_to_idx(pix, clip=False)[source]#
Convert pixel coordinates to pixel indices.
Returns -1 for pixel coordinates that lie outside the map.
- Parameters:
- pixtuple
Tuple of pixel coordinates.
- clipbool
Choose whether to clip indices to the valid range of the geometry. If False then indices for coordinates outside the geometry range will be set -1. Default is False.
- Returns:
- idxtuple
Tuple of pixel indices.
- region_mask(regions)[source]#
Create a mask from a given list of regions.
The mask is filled such that a pixel inside the region is filled with “True”. To invert the mask, e.g. to create a mask with exclusion regions the tilde (~) operator can be used (see example below).
- Parameters:
- regionsstr,
Region
or list ofRegion
Region or list of regions (pixel or sky regions accepted). A region can be defined as a string ind DS9 format as well. See http://ds9.si.edu/doc/ref/region.html for details.
- regionsstr,
- Returns:
- mask_map
WcsNDMap
of boolean type Boolean region mask.
- mask_map
- solid_angle()[source]#
Solid angle array as a
Quantity
insr
.The array has the same dimensionality as
map.nside
as all pixels have the same solid angle.
- to_cube(axes)[source]#
Append non-spatial axes to create a higher-dimensional geometry.
This will result in a new geometry with N+M dimensions where N is the number of current dimensions and M is the number of axes in the list.
- Parameters:
- axeslist
Axes that will be appended to this geometry.
- Returns:
- geom
Geom
Map geometry.
- geom
- to_image()[source]#
Create a 2D image geometry (drop non-spatial dimensions).
- Returns:
- geom
Geom
Image geometry.
- geom
- to_nside(nside)[source]#
Upgrade or downgrade the resolution to a given NSIDE.
- Parameters:
- nsideint
HEALPix NSIDE parameter.
- Returns:
- geom
HpxGeom
A HEALPix geometry object.
- geom
- to_swapped()[source]#
Geometry copy with swapped ORDERING (NEST->RING or vice versa).
- Returns:
- geom
HpxGeom
A HEALPix geometry object.
- geom
- to_wcs_geom(proj='AIT', oversample=2, width_pix=None)[source]#
Make a WCS projection appropriate for this HEALPix pixelization.
- Parameters:
- projstr, optional
Projection type of WCS geometry. Default is “AIT”.
- oversamplefloat, optional
Oversampling factor for WCS map. This will be the approximate ratio of the width of a HEALPix pixel to a WCS pixel. If this parameter is None then the width will be set from
width_pix
. Default is 2.- width_pixint, optional
Width of the WCS geometry in pixels. The pixel size will be set to the number of pixels satisfying
oversample
orwidth_pix
whichever is smaller. If this parameter is None then the width will be set fromoversample
. Default is None.
- Returns:
- wcs
WcsGeom
WCS geometry.
- wcs
- to_wcs_tiles(nside_tiles=4, margin='0 deg')[source]#
Create WCS tiles geometries from HPX geometry with given nside.
The HEALPix geom is divide into superpixels defined by
nside_tiles
, which are then represented by a WCS geometry using a tangential projection. The number of WCS tiles is given by the number of pixels for the givennside_tiles
.- Parameters:
- nside_tilesint, optional
HEALPix NSIDE parameter for super pixel tiles. Default is 4.
- margin
Quantity
, optional Width margin of the WCS tile. Default is “0 deg”.
- Returns:
- wcs_tileslist
List of WCS tile geometries.
- upsample(factor)[source]#
Upsample the spatial dimension of the geometry by a given factor.
- Parameters:
- factorint
Upsampling factor.
- axis_namestr
Axis to upsample.
- Returns:
- geom
Geom
Upsampled geometry.
- geom
- classmethod __new__(*args, **kwargs)#