HpxGeom¶
- 
class gammapy.maps.HpxGeom(nside, nest=True, frame='icrs', region=None, axes=None)[source]¶
- Bases: - gammapy.maps.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 - regionargument.- Parameters
- nsidendarray
- 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 of nside values defining the pixel size for each image plane. If nside is not a scalar then its dimensionality should match that of the non-spatial axes. 
- nestbool
- True -> ‘NESTED’, False -> ‘RING’ indexing scheme 
- framestr
- Coordinate system, “icrs” | “galactic” 
- regionstr or tuple
- 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. 
- axeslist
- Axes for non-spatial dimensions. 
 
- nside
 - Attributes Summary - If the geom contains an energy axis rename it to energy true - List of non-spatial axes. - All axes names - Map coordinates of the center of the geometry (tuple). - Pixel coordinates of the center of the geometry (tuple). - Sky coordinate of the center of the geometry. - Shape of the Numpy data array matching this geometry. - Shape of data of the non-spatial axes and unit spatial axes. - Whether geom has an energy axis - HEALPIX pixel and band indices for every pixel in the map. - Flag for all-sky maps. - Whether the geom non spatial axes have length 1, i.e. - Whether the geom is an image without extra dimensions. - Flag identifying whether this geometry is regular in non-spatial dimensions. - Number of dimensions (int). - Is HEALPix order nested? (bool). - Number of pixels in each band. - Max. - NSIDE in each band. - 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. - contains_pix(pix)- Check if a given pixel coordinate is contained in the geometry. - coord_to_idx(coords[, clip])- Convert map coordinates to pixel indices. - coord_to_pix(coords)- Convert map coordinates to pixel coordinates. - copy(**kwargs)- Copy and overwrite given attributes. - 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. - data_nbytes([dtype])- Estimate memory usage in megabytes of the Numpy data array matching this geometry depending on the given type. - downsample(factor[, axis_name])- Downsample the spatial dimension of the geometry by a given factor. - drop(axis_name)- Drop an axis from the geom. - energy_mask([energy_min, energy_max, …])- Create a mask for a given energy range. - from_hdu(hdu[, hdu_bands])- Create an HPX object from a BinTable HDU. - from_hdulist(hdulist[, hdu, hdu_bands])- Load a geometry object from a FITS HDUList. - 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 coords - is_aligned(other)- Check if HEALPIx geoms and extra axes are aligned. - local_to_global(idx_local)- Compute a local index (partial-sky) from a global (all-sky) index. - pad(pad_width, axis_name)- Pad the geometry at the edges. - 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 - replace_axis(axis)- Replace axis with a new one. - resample_axis(axis)- Resample geom to a new axis binning. - separation(center)- Compute sky separation wrt a given center. - slice_by_idx(slices)- Create a new geometry by slicing the non-spatial axes. - Solid angle array ( - Quantityin- sr).- squash(axis_name)- Squash geom axis. - to_bands_hdu([hdu_bands, format])- 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 2D image geometry (drop non-spatial dimensions). - to_nside(nside)- Upgrade or downgrade the reoslution 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 HPX 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 - 
as_energy_true¶
- If the geom contains an energy axis rename it to energy true 
 - 
axes¶
- List of non-spatial axes. 
 - 
axes_names¶
- All axes names 
 - 
center_coord¶
- Map coordinates of the center of the geometry (tuple). 
 - 
center_pix¶
- Pixel coordinates of the center of the geometry (tuple). 
 - 
data_shape¶
- Shape of the Numpy data array matching this geometry. 
 - 
data_shape_axes¶
- Shape of data of the non-spatial axes and unit spatial axes. 
 - 
frame¶
 - 
has_energy_axis¶
- Whether geom has an energy axis 
 - 
ipix¶
- HEALPIX pixel and band indices for every pixel in the map. 
 - 
is_allsky¶
- Flag for all-sky maps. 
 - 
is_flat¶
- Whether the geom non spatial axes have length 1, i.e. if the geom is equivalent to an image. 
 - 
is_hpx= True¶
 - 
is_image¶
- Whether the geom is an image without extra dimensions. 
 - 
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 (int). 
 - 
nest¶
- Is HEALPix order nested? (bool). 
 - 
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¶
- Max. number of pixels 
 - 
nside¶
- NSIDE in each band. 
 - 
order¶
- 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 - 
contains_pix(pix)¶
- Check if a given pixel coordinate is contained in the geometry. - Parameters
- pixtuple
- Tuple of pixel coordinates. 
 
- Returns
- containmentndarray
- Bool array. 
 
- containment
 
 - 
coord_to_idx(coords, clip=False)¶
- Convert map coordinates to pixel indices. - Parameters
- coordstuple or MapCoord
- 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. 
- 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. 
 
- coordstuple or 
- Returns
- pixtuple
- Tuple of pixel indices in image and band dimensions. Elements set to -1 correspond to coordinates outside the map. 
 
 
 - 
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. 
 
 
 - 
copy(**kwargs)¶
- Copy and overwrite given attributes. - Parameters
- **kwargsdict
- Keyword arguments to overwrite in the map geometry constructor. 
 
- Returns
- copyGeom
- Copied map geometry. 
 
- copy
 
 - 
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
- HEALPix NSIDE parameter. This parameter sets the size of the spatial pixels in the map. 
- binszfloat or ndarray
- Approximate pixel size in degrees. An NSIDE will be chosen that corresponds to a pixel size closest to this value. This option is superseded by nside. 
- nestbool
- True for HEALPIX “NESTED” indexing scheme, False for “RING” scheme 
- frame{“icrs”, “galactic”}, optional
- Coordinate system, either Galactic (“galactic”) or Equatorial (“icrs”). 
- skydirtuple or SkyCoord
- 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. 
- regionstr
- HPX region string. Allows for partial-sky maps. 
- widthfloat
- Diameter of the map in degrees. If set the map will encompass all pixels within a circular region centered on - skydir.
- axeslist
- List of axes for non-spatial dimensions. 
 
- nsideint or 
- Returns
- geomHpxGeom
- 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
- geomGeom
- Cropped geometry. 
 
- geom
 
 - 
data_nbytes(dtype='float32')¶
- Estimate memory usage in megabytes of the Numpy data array matching this geometry depending on the given type. - Parameters
- dtypedata-type
- The desired data-type for the array. Default is “float32” 
 
- Returns
- memoryQuantity
- Estimated memory usage in megabytes (MB) 
 
- memory
 
 - 
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
- geomGeom
- Downsampled geometry. 
 
- geom
 
 - 
drop(axis_name)¶
- Drop an axis from the geom. - Parameters
- axis_namestr
- Name of the axis to remove. 
 
- Returns
- geomGeom
- New geom with the axis removed. 
 
- geom
 
 - 
energy_mask(energy_min=None, energy_max=None, round_to_edge=False)¶
- Create a mask for a given energy range. - The energy bin must be fully contained to be included in the mask. 
 - 
classmethod from_hdu(hdu, hdu_bands=None)[source]¶
- Create an HPX object from a BinTable HDU. - Parameters
- hduBinTableHDU
- The FITS HDU 
- hdu_bandsBinTableHDU
- The BANDS table HDU 
 
- hdu
- Returns
- hpxHpxGeom
- HEALPix geometry. 
 
- hpx
 
 - 
classmethod from_hdulist(hdulist, hdu=None, hdu_bands=None)¶
- Load a geometry object from a FITS HDUList. 
 - 
classmethod from_header(header, hdu_bands=None, format=None)[source]¶
- Create an HPX object from a FITS header. - Parameters
- headerHeader
- The FITS header 
- hdu_bandsBinTableHDU
- The BANDS table HDU. 
- formatstr, optional
- FITS convention. 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
- hpxHpxGeom
- 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 - idxto 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. 
- flatbool, optional
- Return a flattened array containing only coordinates for pixels contained in the geometry. 
 
- 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 - idxto 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. 
- localbool
- Flag to return local or global pixel indices. Local indices run from 0 to the number of pixels in a given image plane. 
- flatbool, optional
- Return a flattened array containing only indices for pixels contained in the geometry. 
 
- 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
- True for ‘NESTED’, False = ‘RING’ 
- regionstr
- HEALPIX region string 
 
- Returns
- ilistndarray
- 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
- Return a raveled index. 
 
- 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
- otherHpxGeom
- Other geom. 
 
- other
- Returns
- alignedbool
- Whether geometries are aligned 
 
 
 - 
local_to_global(idx_local)[source]¶
- Compute a local index (partial-sky) from a global (all-sky) index. - Returns
- idx_globaltuple
- A tuple of pixel index vectors with global HEALPIX pixel indices 
 
 
 - 
pad(pad_width, axis_name)¶
- Pad the geometry at the edges. - Parameters
- pad_width{sequence, array_like, int}
- Number of values padded to the edges of each axis. 
- axis_namestr
- Name of the axis to pad. 
 
- Returns
- geomGeom
- Padded geometry. 
 
- geom
 
 - 
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 of 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. 
 
- 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, Regionor 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_mapWcsNDMapof boolean type
- Boolean region mask 
 
- mask_map
 
 - 
replace_axis(axis)¶
- Replace axis with a new one. 
 - 
resample_axis(axis)¶
- Resample geom to a new axis binning. - This method groups the existing bins into a new binning. 
 - 
slice_by_idx(slices)¶
- Create a new geometry by slicing the non-spatial axes. 
 - 
solid_angle()[source]¶
- Solid angle array ( - Quantityin- sr).- The array has the same dimensionality as - map.nsidesince all pixels have the same solid angle.
 - 
squash(axis_name)¶
- Squash geom axis. - Parameters
- axis_namestr
- Axis to squash. 
 
- Returns
- geomGeom
- Geom with squashed axis. 
 
- geom
 
 - 
to_bands_hdu(hdu_bands=None, format='gadf')¶
 - 
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
- geomGeom
- Map geometry. 
 
- geom
 
 - 
to_image()[source]¶
- Create 2D image geometry (drop non-spatial dimensions). - Returns
- geomGeom
- Image geometry. 
 
- geom
 
 - 
to_nside(nside)[source]¶
- Upgrade or downgrade the reoslution to a given nside - Parameters
- nsideint
- Nside 
 
- Returns
- geomHpxGeom
- A HEALPix geometry object. 
 
- geom
 
 - 
to_swapped()[source]¶
- Geometry copy with swapped ORDERING (NEST->RING or vice versa). - Returns
- geomHpxGeom
- A HEALPix geometry object. 
 
- geom
 
 - 
to_wcs_geom(proj='AIT', oversample=2, width_pix=None)[source]¶
- Make a WCS projection appropriate for this HPX pixelization. - Parameters
- projstr
- Projection type of WCS geometry. 
- oversamplefloat
- Oversampling factor for WCS map. This will be the approximate ratio of the width of a HPX pixel to a WCS pixel. If this parameter is None then the width will be set from - width_pix.
- width_pixint
- Width of the WCS geometry in pixels. The pixel size will be set to the number of pixels satisfying - oversampleor- width_pixwhichever is smaller. If this parameter is None then the width will be set from- oversample.
 
- Returns
- wcsWcsGeom
- 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 given nside_tiles. - Parameters
- nside_tilesint
- Nside for super pixel tiles. Usually nsi 
- marginAngle
- Width margin of the wcs tile 
 
- Returns
- wcs_tileslist
- List of WCS tile geoms.