HpxNDMap

class gammapy.maps.HpxNDMap(geom, data=None, dtype='float32', meta=None, unit='')[source]

Bases: gammapy.maps.HpxMap

HEALPix map with any number of non-spatial dimensions.

This class uses a N+1D numpy array to represent the sequence of HEALPix image planes. Following the convention of WCS-based maps this class uses a column-wise ordering for the data array with the spatial dimension being tied to the last index of the array.

Parameters
geomHpxGeom

HEALPIX geometry object.

datandarray

HEALPIX data array. If none then an empty array will be allocated.

metadict

Dictionary to store meta data.

unitstr or Unit

The map unit

Attributes Summary

data

Data array (ndarray)

geom

Map geometry (Geom)

meta

Map meta (dict)

quantity

Map data times unit (Quantity)

unit

Map unit (Unit)

Methods Summary

apply_edisp(self, edisp)

Apply energy dispersion to map.

coadd(self, map_in[, weights])

Add the contents of map_in to this map.

copy(self, \*\*kwargs)

Copy map instance and overwrite given attributes, except for geometry.

create([nside, binsz, nest, map_type, …])

Factory method to create an empty HEALPix map.

crop(self, crop_width)

Crop the spatial dimensions of the map.

downsample(self, factor[, preserve_counts])

Downsample the spatial dimension by a given factor.

fill_by_coord(self, coords[, weights])

Fill pixels at coords with given weights.

fill_by_idx(self, idx[, weights])

Fill pixels at idx with given weights.

fill_by_pix(self, pix[, weights])

Fill pixels at pix with given weights.

fill_events(self, events)

Fill event coordinates (EventList).

from_geom(geom[, meta, data, map_type, …])

Generate an empty map from a Geom instance.

from_hdu(hdu[, hdu_bands])

Make a HpxNDMap object from a FITS HDU.

from_hdulist(hdu_list[, hdu, hdu_bands])

Make a HpxMap object from a FITS HDUList.

get_by_coord(self, coords)

Return map values at the given map coordinates.

get_by_idx(self, idx)

Return map values at the given pixel indices.

get_by_pix(self, pix)

Return map values at the given pixel coordinates.

get_image_by_coord(self, coords)

Return spatial map at the given axis coordinates.

get_image_by_idx(self, idx)

Return spatial map at the given axis pixel indices.

get_image_by_pix(self, pix)

Return spatial map at the given axis pixel coordinates

interp_by_coord(self, coords[, interp])

Interpolate map values at the given map coordinates.

interp_by_pix(self, pix[, interp])

Interpolate map values at the given pixel coordinates.

iter_by_image(self)

Iterate over image planes of the map.

make_hdu(self[, hdu, hdu_bands, sparse, conv])

Make a FITS HDU with input data.

make_wcs_mapping(self[, sum_bands, proj, …])

Make a HEALPix to WCS mapping object.

pad(self, pad_width[, mode, cval, order])

Pad the spatial dimensions of the map.

plot(self[, method, ax, normalize, proj, …])

Quickplot method.

plot_interactive(self[, rc_params])

Plot map with interactive widgets to explore the non spatial axes.

read(filename[, hdu, hdu_bands, map_type])

Read a map from a FITS file.

set_by_coord(self, coords, vals)

Set pixels at coords with given vals.

set_by_idx(self, idx, vals)

Set pixels at idx with given vals.

set_by_pix(self, pix, vals)

Set pixels at pix with given vals.

slice_by_idx(self, slices)

Slice sub map from map object.

sum_over_axes(self)

Sum over all non-spatial dimensions.

to_hdulist(self[, hdu, hdu_bands, sparse, conv])

Convert to HDUList.

to_swapped(self)

Return a new map with the opposite scheme (ring or nested).

to_ud_graded(self, nside[, preserve_counts])

Upgrade or downgrade the resolution of the map to the chosen nside.

to_wcs(self[, sum_bands, normalize, proj, …])

Make a WCS object and convert HEALPIX data into WCS projection.

upsample(self, factor[, preserve_counts])

Upsample the spatial dimension by a given factor.

write(self, filename[, overwrite])

Write to a FITS file.

Attributes Documentation

data

Data array (ndarray)

geom

Map geometry (Geom)

meta

Map meta (dict)

quantity

Map data times unit (Quantity)

unit

Map unit (Unit)

Methods Documentation

apply_edisp(self, edisp)

Apply energy dispersion to map. Requires energy axis.

Parameters
edispgammapy.irf.EDispKernel

Energy dispersion matrix

Returns
mapWcsNDMap

Map with energy dispersion applied.

coadd(self, map_in, weights=None)

Add the contents of map_in to this map.

This method can be used to combine maps containing integral quantities (e.g. counts) or differential quantities if the maps have the same binning.

Parameters
map_inMap

Input map.

weights: `Map` or `~numpy.ndarray`

The weight factors while adding

copy(self, **kwargs)

Copy map instance and overwrite given attributes, except for geometry.

Parameters
**kwargsdict

Keyword arguments to overwrite in the map constructor.

Returns
copyMap

Copied Map.

classmethod create(nside=None, binsz=None, nest=True, map_type='hpx', frame='icrs', data=None, skydir=None, width=None, dtype='float32', region=None, axes=None, meta=None, unit='')

Factory method to create an empty HEALPix map.

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 correponds 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.

map_type{‘hpx’, ‘hpx-sparse’}

Map type. Selects the class that will be used to instantiate the map.

widthfloat

Diameter of the map in degrees. If None then an all-sky geometry will be created.

axeslist

List of MapAxis objects for each non-spatial dimension.

metadict

Dictionary to store meta data.

unitstr or Unit

The map unit

Returns
mapHpxMap

A HPX map object.

crop(self, crop_width)[source]

Crop the spatial dimensions of the map.

Parameters
crop_width{sequence, array_like, int}

Number of pixels cropped from the edges of each axis. Defined analogously to pad_with from numpy.pad.

Returns
mapMap

Cropped map.

downsample(self, factor, preserve_counts=True)[source]

Downsample the spatial dimension by a given factor.

Parameters
factorint

Downsampling factor.

preserve_countsbool

Preserve the integral over each bin. This should be true if the map is an integral quantity (e.g. counts) and false if the map is a differential quantity (e.g. intensity).

axisstr

Which axis to downsample. By default spatial axes are downsampled.

Returns
mapMap

Downsampled map.

fill_by_coord(self, coords, weights=None)

Fill pixels at coords with given weights.

Parameters
coordstuple or MapCoord

Coordinate arrays for each dimension of the map. Tuple should be ordered as (lon, lat, x_0, …, x_n) where x_i are coordinates for non-spatial dimensions of the map.

weightsndarray

Weights vector. Default is weight of one.

fill_by_idx(self, idx, weights=None)[source]

Fill pixels at idx with given weights.

Parameters
idxtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps.

weightsndarray

Weights vector. Default is weight of one.

fill_by_pix(self, pix, weights=None)

Fill pixels at pix with given weights.

Parameters
pixtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps. Pixel indices can be either float or integer type. Float indices will be rounded to the nearest integer.

weightsndarray

Weights vector. Default is weight of one.

fill_events(self, events)

Fill event coordinates (EventList).

static from_geom(geom, meta=None, data=None, map_type='auto', unit='', dtype='float32')

Generate an empty map from a Geom instance.

Parameters
geomGeom

Map geometry.

datanumpy.ndarray

data array

metadict

Dictionary to store meta data.

map_type{‘wcs’, ‘wcs-sparse’, ‘hpx’, ‘hpx-sparse’, ‘auto’}

Map type. Selects the class that will be used to instantiate the map. The map type should be consistent with the geometry. If map_type is ‘auto’ then an appropriate map type will be inferred from type of geom.

unitstr or Unit

Data unit.

Returns
map_outMap

Map object

classmethod from_hdu(hdu, hdu_bands=None)[source]

Make a HpxNDMap object from a FITS HDU.

Parameters
hduBinTableHDU

The FITS HDU

hdu_bandsBinTableHDU

The BANDS table HDU

classmethod from_hdulist(hdu_list, hdu=None, hdu_bands=None)

Make a HpxMap object from a FITS HDUList.

Parameters
hdu_listHDUList

HDU list containing HDUs for map data and bands.

hdustr

Name or index of the HDU with the map data. If None then the method will try to load map data from the first BinTableHDU in the file.

hdu_bandsstr

Name or index of the HDU with the BANDS table.

Returns
hpx_mapHpxMap

Map object

get_by_coord(self, coords)

Return map values at the given map coordinates.

Parameters
coordstuple or MapCoord

Coordinate arrays for each dimension of the map. Tuple should be ordered as (lon, lat, x_0, …, x_n) where x_i are coordinates for non-spatial dimensions of the map.

Returns
valsndarray

Values of pixels in the map. np.nan used to flag coords outside of map.

get_by_idx(self, idx)[source]

Return map values at the given pixel indices.

Parameters
idxtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps.

Returns
valsndarray

Array of pixel values. np.nan used to flag coordinate outside of map

get_by_pix(self, pix)

Return map values at the given pixel coordinates.

Parameters
pixtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps. Pixel indices can be either float or integer type.

Returns
valsndarray

Array of pixel values. np.nan used to flag coordinates outside of map

get_image_by_coord(self, coords)

Return spatial map at the given axis coordinates.

Parameters
coordstuple or dict

Tuple should be ordered as (x_0, …, x_n) where x_i are coordinates for non-spatial dimensions of the map. Dict should specify the axis names of the non-spatial axes such as {‘axes0’: x_0, …, ‘axesn’: x_n}.

Returns
map_outMap

Map with spatial dimensions only.

Examples

import numpy as np
from gammapy.maps import Map, MapAxis
from astropy.coordinates import SkyCoord
from astropy import units as u

# Define map axes
energy_axis = MapAxis.from_edges(
    np.logspace(-1., 1., 4), unit='TeV', name='energy',
)

time_axis = MapAxis.from_edges(
    np.linspace(0., 10, 20), unit='h', name='time',
)

# Define map center
skydir = SkyCoord(0, 0, frame='galactic', unit='deg')

# Create map
m_wcs = Map.create(
    map_type='wcs',
    binsz=0.02,
    skydir=skydir,
    width=10.0,
    axes=[energy_axis, time_axis],
)

# Get image by coord tuple
image = m_wcs.get_image_by_coord(('500 GeV', '1 h'))

# Get image by coord dict with strings
image = m_wcs.get_image_by_coord({'energy': '500 GeV', 'time': '1 h'})

# Get image by coord dict with quantities
image = m_wcs.get_image_by_coord({'energy': 0.5 * u.TeV, 'time': 1 * u.h})
get_image_by_idx(self, idx)

Return spatial map at the given axis pixel indices.

Parameters
idxtuple

Tuple of scalar indices for each non spatial dimension of the map. Tuple should be ordered as (I_0, …, I_n).

Returns
map_outMap

Map with spatial dimensions only.

get_image_by_pix(self, pix)

Return spatial map at the given axis pixel coordinates

Parameters
pixtuple

Tuple of scalar pixel coordinates for each non-spatial dimension of the map. Tuple should be ordered as (I_0, …, I_n). Pixel coordinates can be either float or integer type.

Returns
map_outMap

Map with spatial dimensions only.

interp_by_coord(self, coords, interp=1)[source]

Interpolate map values at the given map coordinates.

Parameters
coordstuple or MapCoord

Coordinate arrays for each dimension of the map. Tuple should be ordered as (lon, lat, x_0, …, x_n) where x_i are coordinates for non-spatial dimensions of the map.

interp{None, ‘nearest’, ‘linear’, ‘cubic’, 0, 1, 2, 3}

Method to interpolate data values. By default no interpolation is performed and the return value will be the amplitude of the pixel encompassing the given coordinate. Integer values can be used in lieu of strings to choose the interpolation method of the given order (0=’nearest’, 1=’linear’, 2=’quadratic’, 3=’cubic’). Note that only ‘nearest’ and ‘linear’ methods are supported for all map types.

fill_valueNone or float value

The value to use for points outside of the interpolation domain. If None, values outside the domain are extrapolated.

Returns
valsndarray

Interpolated pixel values.

interp_by_pix(self, pix, interp=None)[source]

Interpolate map values at the given pixel coordinates.

iter_by_image(self)

Iterate over image planes of the map.

This is a generator yielding (data, idx) tuples, where data is a numpy.ndarray view of the image plane data, and idx is a tuple of int, the index of the image plane.

The image plane index is in data order, so that the data array can be indexed directly. See Iterating by image for further information.

make_hdu(self, hdu=None, hdu_bands=None, sparse=False, conv=None)

Make a FITS HDU with input data.

Parameters
hdustr

The HDU extension name.

hdu_bandsstr

The HDU extension name for BANDS table.

sparsebool

Set INDXSCHM to SPARSE and sparsify the map by only writing pixels with non-zero amplitude.

conv{‘fgst-ccube’, ‘fgst-template’, ‘gadf’, None}, optional

FITS format convention. If None this will be set to the default convention of the map.

Returns
hdu_outBinTableHDU or ImageHDU

Output HDU containing map data.

make_wcs_mapping(self, sum_bands=False, proj='AIT', oversample=2, width_pix=None)[source]

Make a HEALPix to WCS mapping object.

Parameters
sum_bandsbool

sum over non-spatial dimensions before reprojecting

projstr

WCS-projection

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 oversample or width_pix whichever is smaller. If this parameter is None then the width will be set from oversample.

Returns
hpx2wcsHpxToWcsMapping
pad(self, pad_width, mode='constant', cval=0, order=1)[source]

Pad the spatial dimensions of the map.

Parameters
pad_width{sequence, array_like, int}

Number of pixels padded to the edges of each axis.

mode{‘edge’, ‘constant’, ‘interp’}

Padding mode. ‘edge’ pads with the closest edge value. ‘constant’ pads with a constant value. ‘interp’ pads with an extrapolated value.

cvalfloat

Padding value when mode=’consant’.

orderint

Order of interpolation when mode=’constant’ (0 = nearest-neighbor, 1 = linear, 2 = quadratic, 3 = cubic).

Returns
mapMap

Padded map.

plot(self, method='raster', ax=None, normalize=False, proj='AIT', oversample=2, width_pix=1000, **kwargs)[source]

Quickplot method.

This will generate a visualization of the map by converting to a rasterized WCS image (method=’raster’) or drawing polygons for each pixel (method=’poly’).

Parameters
method{‘raster’,’poly’}

Method for mapping HEALPix pixels to a two-dimensional image. Can be set to ‘raster’ (rasterization to cartesian image plane) or ‘poly’ (explicit polygons for each pixel). WARNING: The ‘poly’ method is much slower than ‘raster’ and only suitable for maps with less than ~10k pixels.

projstring, optional

Any valid WCS projection type.

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 oversample or width_pix whichever is smaller. If this parameter is None then the width will be set from oversample.

**kwargsdict

Keyword arguments passed to imshow.

Returns
figFigure

Figure object.

axWCSAxes

WCS axis object

imAxesImage or PatchCollection

Image object.

plot_interactive(self, rc_params=None, **kwargs)

Plot map with interactive widgets to explore the non spatial axes.

Parameters
rc_paramsdict

Passed to matplotlib.rc_context(rc=rc_params) to style the plot.

**kwargsdict

Keyword arguments passed to WcsNDMap.plot.

Examples

You can try this out e.g. using a Fermi-LAT diffuse model cube with an energy axis:

from gammapy.maps import Map

m = Map.read("$GAMMAPY_DATA/fermi_3fhl/gll_iem_v06_cutout.fits")
m.plot_interactive(add_cbar=True, stretch="sqrt")

If you would like to adjust the figure size you can use the rc_params argument:

rc_params = {'figure.figsize': (12, 6), 'font.size': 12}
m.plot_interactive(rc_params=rc_params)
static read(filename, hdu=None, hdu_bands=None, map_type='auto')

Read a map from a FITS file.

Parameters
filenamestr or Path

Name of the FITS file.

hdustr

Name or index of the HDU with the map data.

hdu_bandsstr

Name or index of the HDU with the BANDS table. If not defined this will be inferred from the FITS header of the map HDU.

map_type{‘wcs’, ‘wcs-sparse’, ‘hpx’, ‘hpx-sparse’, ‘auto’}

Map type. Selects the class that will be used to instantiate the map. The map type should be consistent with the format of the input file. If map_type is ‘auto’ then an appropriate map type will be inferred from the input file.

Returns
map_outMap

Map object

set_by_coord(self, coords, vals)

Set pixels at coords with given vals.

Parameters
coordstuple or MapCoord

Coordinate arrays for each dimension of the map. Tuple should be ordered as (lon, lat, x_0, …, x_n) where x_i are coordinates for non-spatial dimensions of the map.

valsndarray

Values vector.

set_by_idx(self, idx, vals)[source]

Set pixels at idx with given vals.

Parameters
idxtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps.

valsndarray

Values vector.

set_by_pix(self, pix, vals)

Set pixels at pix with given vals.

Parameters
pixtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps. Pixel indices can be either float or integer type. Float indices will be rounded to the nearest integer.

valsndarray

Values vector.

slice_by_idx(self, slices)

Slice sub map from map object.

For usage examples, see Indexing and Slicing.

Parameters
slicesdict

Dict of axes names and integers or slice object pairs. Contains one element for each non-spatial dimension. For integer indexing the corresponding axes is dropped from the map. Axes not specified in the dict are kept unchanged.

Returns
map_outMap

Sliced map object.

sum_over_axes(self)[source]

Sum over all non-spatial dimensions.

Returns
map_outHpxNDMap

Summed map.

to_hdulist(self, hdu='SKYMAP', hdu_bands=None, sparse=False, conv='gadf')

Convert to HDUList.

Parameters
hdustr

The HDU extension name.

hdu_bandsstr

The HDU extension name for BANDS table.

sparsebool

Set INDXSCHM to SPARSE and sparsify the map by only writing pixels with non-zero amplitude.

conv{‘fgst-ccube’,’fgst-template’,’gadf’,None}, optional

FITS format convention. If None this will be set to the default convention of the map.

Returns
hdu_listHDUList
to_swapped(self)[source]

Return a new map with the opposite scheme (ring or nested).

Returns
mapHpxMap

Map object.

to_ud_graded(self, nside, preserve_counts=False)[source]

Upgrade or downgrade the resolution of the map to the chosen nside.

Parameters
nsideint

NSIDE parameter of the new map.

preserve_countsbool

Choose whether to preserve counts (total amplitude) or intensity (amplitude per unit solid angle).

Returns
mapHpxMap

Map object.

to_wcs(self, sum_bands=False, normalize=True, proj='AIT', oversample=2, width_pix=None, hpx2wcs=None)[source]

Make a WCS object and convert HEALPIX data into WCS projection.

Parameters
sum_bandsbool

Sum over non-spatial axes before reprojecting. If False then the WCS map will have the same dimensionality as the HEALPix one.

normalizebool

Preserve integral by splitting HEALPIX values between bins?

projstr

WCS-projection

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 oversample or width_pix whichever is smaller. If this parameter is None then the width will be set from oversample.

hpx2wcsHpxToWcsMapping

Set the HPX to WCS mapping object that will be used to generate the WCS map. If none then a new mapping will be generated based on proj and oversample arguments.

Returns
map_outWcsMap

WCS map object.

upsample(self, factor, preserve_counts=True)[source]

Upsample the spatial dimension by a given factor.

Parameters
factorint

Upsampling factor.

orderint

Order of the interpolation used for upsampling.

preserve_countsbool

Preserve the integral over each bin. This should be true if the map is an integral quantity (e.g. counts) and false if the map is a differential quantity (e.g. intensity).

axisstr

Which axis to upsample. By default spatial axes are upsampled.

Returns
mapMap

Upsampled map.

write(self, filename, overwrite=False, **kwargs)

Write to a FITS file.

Parameters
filenamestr

Output file name.

overwritebool

Overwrite existing file?

hdustr

Set the name of the image extension. By default this will be set to SKYMAP (for BINTABLE HDU) or PRIMARY (for IMAGE HDU).

hdu_bandsstr

Set the name of the bands table extension. By default this will be set to BANDS.

convstr

FITS format convention. By default files will be written to the gamma-astro-data-formats (GADF) format. This option can be used to write files that are compliant with format conventions required by specific software (e.g. the Fermi Science Tools). Supported conventions are ‘gadf’, ‘fgst-ccube’, ‘fgst-ltcube’, ‘fgst-bexpcube’, ‘fgst-template’, ‘fgst-srcmap’, ‘fgst-srcmap-sparse’, ‘galprop’, and ‘galprop2’.

sparsebool

Sparsify the map by dropping pixels with zero amplitude. This option is only compatible with the ‘gadf’ format.