Map

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

Bases: abc.ABC

Abstract map class.

This can represent WCS- or HEALPIX-based maps with 2 spatial dimensions and N non-spatial dimensions.

Parameters:
geom : MapGeom

Geometry

data : ndarray

Data array

meta : OrderedDict

Dictionary to store meta data.

unit : str or Unit

Data unit

Attributes Summary

data Data array (ndarray)
geom Map geometry (MapGeom)
meta Map meta (OrderedDict)
quantity Map data times unit (Quantity)
unit Map unit (Unit)

Methods Summary

coadd(self, map_in) Add the contents of map_in to this map.
copy(self, \*\*kwargs) Copy map instance and overwrite given attributes, except for geometry.
create(\*\*kwargs) Create an empty map object.
crop(self, crop_width) Crop the spatial dimensions of the map.
downsample(self, factor[, preserve_counts, axis]) 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.
from_geom(geom[, meta, data, map_type, unit]) Generate an empty map from a MapGeom instance.
from_hdulist(hdulist[, hdu, hdu_bands, map_type]) Create from astropy.io.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, fill_value]) Interpolate map values at the given pixel coordinates.
iter_by_image(self) Iterate over image planes of the map.
pad(self, pad_width[, mode, cval, order]) Pad the spatial dimensions of the map.
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.
reproject(self, geom[, order, mode]) Reproject this map to a different geometry.
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[, keepdims]) Reduce to a 2D image by summing over non-spatial dimensions.
upsample(self, factor[, order, …]) 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 (MapGeom)

meta

Map meta (OrderedDict)

quantity

Map data times unit (Quantity)

unit

Map unit (Unit)

Methods Documentation

coadd(self, map_in)[source]

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_in : Map

Input map.

copy(self, **kwargs)[source]

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

Parameters:
**kwargs : dict

Keyword arguments to overwrite in the map constructor.

Returns:
copy : Map

Copied Map.

static create(**kwargs)[source]

Create an empty map object.

This method accepts generic options listed below, as well as options for HpxMap and WcsMap objects. For WCS-specific options, see WcsMap.create and for HPX-specific options, see HpxMap.create.

Parameters:
coordsys : str

Coordinate system, either Galactic (‘GAL’) or Equatorial (‘CEL’).

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

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

binsz : float or ndarray

Pixel size in degrees.

skydir : SkyCoord

Coordinate of map center.

axes : list

List of MapAxis objects for each non-spatial dimension. If None then the map will be a 2D image.

dtype : str

Data type, default is ‘float32’

unit : str or Unit

Data unit.

meta : OrderedDict

Dictionary to store meta data.

Returns:
map : Map

Empty 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:
map : Map

Cropped map.

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

Downsample the spatial dimension by a given factor.

Parameters:
factor : int

Downsampling factor.

preserve_counts : bool

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

axis : str

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

Returns:
map : Map

Downsampled map.

fill_by_coord(self, coords, weights=None)[source]

Fill pixels at coords with given weights.

Parameters:
coords : tuple 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.

weights : ndarray

Weights vector. Default is weight of one.

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

Fill pixels at idx with given weights.

Parameters:
idx : tuple

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.

weights : ndarray

Weights vector. Default is weight of one.

fill_by_pix(self, pix, weights=None)[source]

Fill pixels at pix with given weights.

Parameters:
pix : tuple

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.

weights : ndarray

Weights vector. Default is weight of one.

static from_geom(geom, meta=None, data=None, map_type='auto', unit='')[source]

Generate an empty map from a MapGeom instance.

Parameters:
geom : MapGeom

Map geometry.

data : numpy.ndarray

data array

meta : OrderedDict

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.

unit : str or Unit

Data unit.

Returns:
map_out : Map

Map object

static from_hdulist(hdulist, hdu=None, hdu_bands=None, map_type='auto')[source]

Create from astropy.io.fits.HDUList.

get_by_coord(self, coords)[source]

Return map values at the given map coordinates.

Parameters:
coords : tuple 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:
vals : ndarray

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:
idx : tuple

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:
vals : ndarray

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

get_by_pix(self, pix)[source]

Return map values at the given pixel coordinates.

Parameters:
pix : tuple

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:
vals : ndarray

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

get_image_by_coord(self, coords)[source]

Return spatial map at the given axis coordinates.

Parameters:
coords : tuple 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_out : Map

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)[source]

Return spatial map at the given axis pixel indices.

Parameters:
idx : tuple

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

Returns:
map_out : Map

Map with spatial dimensions only.

get_image_by_pix(self, pix)[source]

Return spatial map at the given axis pixel coordinates

Parameters:
pix : tuple

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_out : Map

Map with spatial dimensions only.

interp_by_coord(self, coords, interp=None, fill_value=None)[source]

Interpolate map values at the given map coordinates.

Parameters:
coords : tuple 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_value : None or float value

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

Returns:
vals : ndarray

Interpolated pixel values.

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

Interpolate map values at the given pixel coordinates.

Parameters:
pix : tuple

Tuple of pixel coordinate arrays for each dimension of the map. Tuple should be ordered as (p_lon, p_lat, p_0, …, p_n) where p_i are pixel 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_value : None or float value

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

Returns:
vals : ndarray

Interpolated pixel values.

iter_by_image(self)[source]

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.

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.

cval : float

Padding value when mode=’consant’.

order : int

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

Returns:
map : Map

Padded map.

plot_interactive(self, rc_params=None, **kwargs)[source]

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

Parameters:
rc_params : dict

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

**kwargs : dict

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')[source]

Read a map from a FITS file.

Parameters:
filename : str or Path

Name of the FITS file.

hdu : str

Name or index of the HDU with the map data.

hdu_bands : str

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_out : Map

Map object

reproject(self, geom, order=1, mode='interp')[source]

Reproject this map to a different geometry.

Only spatial axes are reprojected, if you would like to reproject non-spatial axes consider using Map.interp_by_coord() instead.

Parameters:
geom : MapGeom

Geometry of projection.

mode : {‘interp’, ‘exact’}

Method for reprojection. ‘interp’ method interpolates at pixel centers. ‘exact’ method integrates over intersection of pixels.

order : int or str

Order of interpolating polynomial (0 = nearest-neighbor, 1 = linear, 2 = quadratic, 3 = cubic).

Returns:
map : Map

Reprojected map.

set_by_coord(self, coords, vals)[source]

Set pixels at coords with given vals.

Parameters:
coords : tuple 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.

vals : ndarray

Values vector.

set_by_idx(self, idx, vals)[source]

Set pixels at idx with given vals.

Parameters:
idx : tuple

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.

vals : ndarray

Values vector.

set_by_pix(self, pix, vals)[source]

Set pixels at pix with given vals.

Parameters:
pix : tuple

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.

vals : ndarray

Values vector.

slice_by_idx(self, slices)[source]

Slice sub map from map object.

For usage examples, see Indexing and Slicing.

Parameters:
slices : dict

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_out : Map

Sliced map object.

sum_over_axes(self, keepdims=False)[source]

Reduce to a 2D image by summing over non-spatial dimensions.

upsample(self, factor, order=0, preserve_counts=True, axis=None)[source]

Upsample the spatial dimension by a given factor.

Parameters:
factor : int

Upsampling factor.

order : int

Order of the interpolation used for upsampling.

preserve_counts : bool

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

axis : str

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

Returns:
map : Map

Upsampled map.

write(self, filename, overwrite=False, **kwargs)[source]

Write to a FITS file.

Parameters:
filename : str

Output file name.

overwrite : bool

Overwrite existing file?

hdu : str

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

hdu_bands : str

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

conv : str

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

sparse : bool

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