RegionNDMap#

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

Bases: Map

N-dimensional region map.

A RegionNDMap owns a RegionGeom instance as well as a data array containing the values associated to that region in the sky along the non-spatial axis, usually an energy axis. The spatial dimensions of a RegionNDMap are reduced to a single spatial bin with an arbitrary shape, and any extra dimensions are described by an arbitrary number of non-spatial axes.

Parameters:
geomRegionGeom

Region geometry object.

datandarray

Data array. If None then an empty array will be allocated.

dtypestr, optional

Data type. Default is “float32”.

metadict, optional

Dictionary to store metadata. Default is None.

unitstr or Unit, optional

The map unit. Default is “”.

Attributes Summary

data

Data array as a ndarray object.

geom

Map geometry as a Geom object.

is_mask

Whether map is a mask with boolean data type.

meta

Map metadata as a dict.

quantity

Map data as a Quantity object.

tag

unit

Map unit as an Unit object.

Methods Summary

coadd(map_in[, weights])

Add the contents of map_in to this map.

copy(**kwargs)

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

create(region[, axes, dtype, meta, unit, ...])

Create an empty region map object.

crop()

Crop the spatial dimensions of the map.

cumsum(axis_name)

Compute cumulative sum along a given axis.

cutout(*args, **kwargs)

Return self.

dot(other)

Apply dot product with the input map.

downsample(factor[, preserve_counts, ...])

Downsample the non-spatial dimension by a given factor.

fill_by_coord(coords[, weights])

Fill pixels at coords with given weights.

fill_by_idx(idx[, weights])

Fill pixels at idx with given weights.

fill_by_pix(pix[, weights])

Fill pixels at pix with given weights.

fill_events(events[, weights])

Fill the map from an EventList object.

from_geom(geom[, meta, data, unit, dtype])

Generate an empty map from a Geom instance.

from_hdulist(hdulist[, format, ogip_column, hdu])

Create from HDUList.

from_stack(maps[, axis, axis_name])

Create Map from a list of images and a non-spatial axis.

from_table(table[, format, colname])

Create region map from table.

get_by_coord(coords[, fill_value])

Return map values at the given map coordinates.

get_by_idx(idxs)

Return map values at the given pixel indices.

get_by_pix(pix[, fill_value])

Return map values at the given pixel coordinates.

get_image_by_coord(coords)

Return spatial map at the given axis coordinates.

get_image_by_idx(idx)

Return spatial map at the given axis pixel indices.

get_image_by_pix(pix)

Return spatial map at the given axis pixel coordinates

get_spectrum(*args, **kwargs)

Return self.

integral(axis_name, coords, **kwargs)

Compute integral along a given axis.

interp_by_coord(coords, **kwargs)

Interpolate map values at the given map coordinates.

interp_by_pix(pix, **kwargs)

Interpolate map values at the given pixel coordinates.

interp_to_geom(geom[, preserve_counts, ...])

Interpolate map to input geometry.

is_allclose(other[, rtol_axes, atol_axes])

Compare two Maps for close equivalency.

iter_by_axis(axis_name[, keepdims])

Iterate over a given axis.

iter_by_axis_data(axis_name)

Iterate data by axis.

iter_by_image([keepdims])

Iterate over image planes of a map.

iter_by_image_data()

Iterate over image planes of the map.

iter_by_image_index()

Iterate over image planes of the map.

mask_nearest_position(position)

Given a sky coordinate return nearest valid position in the mask.

normalize([axis_name])

Normalise data in place along a given axis.

pad(pad_width[, axis_name, mode, cval, method])

Pad the spatial dimensions of the map.

plot([ax, axis_name])

Plot the data contained in region map along the non-spatial axis.

plot_grid([figsize, ncols])

Plot map as a grid of subplots for non-spatial axes.

plot_hist([ax])

Plot as histogram.

plot_interactive()

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

plot_mask([ax])

Plot the mask as a shaded area in a xmin-xmax range.

plot_region([ax])

Plot region.

read(filename[, format, ogip_column, hdu, ...])

Read from file.

reduce(axis_name[, func, keepdims, weights])

Reduce map over a single non-spatial axis.

reduce_over_axes([func, keepdims, ...])

Reduce map over non-spatial axes.

rename_axes(names, new_names)

Rename the Map axes.

reorder_axes(axes_names)

Return a new map re-ordering the non-spatial axes.

reproject_by_image(geom[, preserve_counts, ...])

Reproject each image of a ND map to input 2d geometry.

reproject_to_geom(geom[, preserve_counts, ...])

Reproject map to input geometry.

resample(geom[, weights, preserve_counts])

Resample pixels to geom with given weights.

resample_axis(axis[, weights, ufunc])

Resample map to a new axis by grouping and reducing smaller bins by a given function ufunc.

sample_coord(n_events[, random_state])

Sample position and energy of events.

set_by_coord(coords, vals)

Set pixels at coords with given vals.

set_by_idx(idx, value)

Set pixels at idx with given vals.

set_by_pix(pix, vals)

Set pixels at pix with given vals.

slice_by_idx(slices)

Slice sub map from map object.

split_by_axis(axis_name)

Split a Map along an axis into multiple maps.

stack(other[, weights, nan_to_num])

Stack other region map into map.

sum_over_axes([axes_names, keepdims, weights])

Sum map values over all non-spatial axes.

to_cube(axes)

Append non-spatial axes to create a higher-dimensional Map.

to_hdulist([format, hdu, hdu_bands, hdu_region])

Convert to HDUList.

to_region_nd_map(*args, **kwargs)

Return self.

to_table([format])

Convert to Table.

to_unit(unit)

Convert map to a different unit.

upsample(factor[, order, preserve_counts, ...])

Upsample the non-spatial dimension by a given factor.

write(filename[, overwrite, format, hdu, ...])

Write map to file.

Attributes Documentation

data#

Data array as a ndarray object.

geom#

Map geometry as a Geom object.

is_mask#

Whether map is a mask with boolean data type.

meta#

Map metadata as a dict.

quantity#

Map data as a Quantity object.

tag = 'map'#
unit#

Map unit as an Unit object.

Methods Documentation

coadd(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

Map to add.

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

The weight factors while adding. Default is None.

copy(**kwargs)#

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

Parameters:
**kwargsdict, optional

Keyword arguments to overwrite in the map constructor.

Returns:
copyMap

Copied Map.

classmethod create(region, axes=None, dtype='float32', meta=None, unit='', wcs=None, binsz_wcs='0.1 deg', data=None)[source]#

Create an empty region map object.

Parameters:
regionstr or SkyRegion

Region specification.

axeslist of MapAxis, optional

Non-spatial axes. Default is None.

dtypestr, optional

Data type. Default is ‘float32’.

unitstr or Unit, optional

Data unit. Default is “”.

metadict, optional

Dictionary to store metadata. Default is None.

wcsWCS, optional

WCS projection to use for local projections of the region. Default is None.

binsz_wcs: `~astropy.units.Quantity` or str, optional

Bin size used for the default WCS, if wcs=None. Default is “0.1 deg”.

datandarray, optional

Data array. Default is None.

Returns:
mapRegionNDMap

Region map.

crop()[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.

cumsum(axis_name)#

Compute cumulative sum along a given axis.

Parameters:
axis_namestr

Along which axis to sum.

Returns:
cumsumMap

Map with cumulative sum.

cutout(*args, **kwargs)[source]#

Return self.

dot(other)#

Apply dot product with the input map.

The input Map has to share a single MapAxis with the current Map. Because it has no spatial dimension, it must be a RegionNDMap.

Parameters:
otherRegionNDMap

Map to apply the dot product to. It must share a unique non-spatial MapAxis with the current Map.

Returns:
mapMap

Map with dot product applied.

downsample(factor, preserve_counts=True, axis_name=None, weights=None)[source]#

Downsample the non-spatial dimension by a given factor.

By default, the first axes is downsampled.

Parameters:
factorint

Downsampling factor.

preserve_countsbool, optional

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). Default is True.

axis_namestr, optional

Which axis to downsample. Default is “energy”.

weightsRegionNDMap, optional

Contains the weights to apply to the axis to reduce. Default weights of one.

Returns:
mapRegionNDMap

Downsampled region map.

fill_by_coord(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, optional

Weights vector. If None, weights are set to 1. Default is None.

fill_by_idx(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, optional

Weights vector. If None, weights are set to 1. Default is None.

fill_by_pix(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, optional

Weights vector. If None, weights are set to 1. Default is None.

fill_events(events, weights=None)#

Fill the map from an EventList object.

Parameters:
eventsEventList

Events to fill in the map with.

weightsndarray, optional

Weights vector. The weights vector must be of the same length as the events column length. If None, weights are set to 1. Default is None.

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

Generate an empty map from a Geom instance.

Parameters:
geomGeom

Map geometry.

metadict, optional

Dictionary to store metadata. Default is None.

datanumpy.ndarray, optional

Data array. Default is None.

unitstr or Unit

Data unit.

dtypestr, optional

Data type. Default is ‘float32’.

Returns:
map_outMap

Map object.

classmethod from_hdulist(hdulist, format='gadf', ogip_column=None, hdu=None, **kwargs)[source]#

Create from HDUList.

Parameters:
hdulistHDUList

HDU list.

format{“gadf”, “ogip”, “ogip-arf”}

Format specification. Default is “gadf”.

ogip_column{“COUNTS”, “QUALITY”, “BACKSCAL”}, optional

If format ‘ogip’ is chosen which table HDU column to read. Default is None.

hdustr, optional

Name or index of the HDU with the map data. Default is None.

Returns:
region_nd_mapRegionNDMap

Region map.

classmethod from_stack(maps, axis=None, axis_name=None)#

Create Map from a list of images and a non-spatial axis.

The image geometries must be aligned, except for the axis that is stacked.

Parameters:
mapslist of Map objects

List of maps.

axisMapAxis, optional

If a MapAxis is provided the maps are stacked along the last data axis and the new axis is introduced. Default is None.

axis_namestr, optional

If an axis name is as string the given the maps are stacked along the given axis name.

Returns:
mapMap

Map with additional non-spatial axis.

classmethod from_table(table, format='', colname=None)[source]#

Create region map from table.

Parameters:
tableTable

Table with input data.

format{“gadf-sed”, “lightcurve”, “profile”}

Format to use.

colnamestr

Column name to take the data from.

Returns:
region_mapRegionNDMap

Region map.

get_by_coord(coords, fill_value=nan)#

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.

fill_valuefloat

Value which is returned if the position is outside the projection footprint. Default is numpy.nan.

Returns:
valsndarray

Values of pixels in the map. numpy.nan is used to flag coordinates outside the map.

get_by_idx(idxs)[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. numpy.nan is used to flag coordinates outside the map.

get_by_pix(pix, fill_value=nan)#

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.

fill_valuefloat

Value which is returned if the position is outside the projection footprint. Default is numpy.nan.

Returns:
valsndarray

Array of pixel values. numpy.nan is used to flag coordinates outside the map.

get_image_by_coord(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. Dictionary 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.

See also

get_image_by_idx, get_image_by_pix.

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

See also

get_image_by_coord, get_image_by_pix.
get_image_by_pix(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.

See also

get_image_by_coord, get_image_by_idx.
get_spectrum(*args, **kwargs)[source]#

Return self.

integral(axis_name, coords, **kwargs)#

Compute integral along a given axis.

This method uses interpolation of the cumulative sum.

Parameters:
axis_namestr

Along which axis to integrate.

coordsdict or MapCoord

Map coordinates.

**kwargsdict, optional

Keyword arguments passed to Map.interp_by_coord.

Returns:
arrayQuantity

2D array with axes offset.

interp_by_coord(coords, **kwargs)[source]#

Interpolate map values at the given map coordinates.

Parameters:
coordstuple, dict 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. “lon” and “lat” are optional and will be taken at the center of the region by default.

method{“linear”, “nearest”}

Method to interpolate data values. Default is “linear”.

fill_valueNone or float value

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

values_scale{“lin”, “log”, “sqrt”}

Optional value scaling.

Returns:
valsndarray

Interpolated pixel values.

interp_by_pix(pix, **kwargs)[source]#

Interpolate map values at the given pixel coordinates.

Parameters:
pixtuple

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.

method{“linear”, “nearest”}

Method to interpolate data values. Default is “linear”.

fill_valuefloat, optional

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

Returns:
valsndarray

Interpolated pixel values.

interp_to_geom(geom, preserve_counts=False, fill_value=0, **kwargs)#

Interpolate map to input geometry.

Parameters:
geomGeom

Target Map geometry.

preserve_countsbool, optional

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). Default is False.

fill_valuefloat, optional

The value to use for points outside the interpolation domain. Default is 0.

**kwargsdict, optional

Keyword arguments passed to Map.interp_by_coord.

Returns:
interp_mapMap

Interpolated Map.

is_allclose(other, rtol_axes=0.001, atol_axes=1e-06, **kwargs)#

Compare two Maps for close equivalency.

Parameters:
othergammapy.maps.Map

The Map to compare against.

rtol_axesfloat, optional

Relative tolerance for the axes’ comparison. Default is 1e-3.

atol_axesfloat, optional

Absolute tolerance for the axes’ comparison. Default is 1e-6.

**kwargsdict, optional

Keywords passed to allclose.

Returns:
is_allclosebool

Whether the Map is all close.

iter_by_axis(axis_name, keepdims=False)#

Iterate over a given axis.

Yields:
mapMap

Map iteration.

See also

iter_by_image

iterate by image returning a map.

iter_by_axis_data(axis_name)[source]#

Iterate data by axis.

Parameters:
axis_namestr

Axis name.

Returns:
idx, datatuple, Quantity

Data and index.

iter_by_image(keepdims=False)#

Iterate over image planes of a map.

Parameters:
keepdimsbool, optional

Keep dimensions. Default is False.

Yields:
mapMap

Map iteration.

See also

iter_by_image_data

iterate by image returning data and index.

iter_by_image_data()#

Iterate over image planes of the map.

The image plane index is in data order, so that the data array can be indexed directly.

Yields:
(data, idx)tuple

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.

See also

iter_by_image

iterate by image returning a map.

iter_by_image_index()#

Iterate over image planes of the map.

The image plane index is in data order, so that the data array can be indexed directly.

Yields:
idxtuple

idx is a tuple of int, the index of the image plane.

See also

iter_by_image

iterate by image returning a map.

mask_nearest_position(position)#

Given a sky coordinate return nearest valid position in the mask.

If the mask contains additional axes, the mask is reduced over those axes.

Parameters:
positionSkyCoord

Test position.

Returns:
positionSkyCoord

The nearest position in the mask.

normalize(axis_name=None)#

Normalise data in place along a given axis.

Parameters:
axis_namestr, optional

Along which axis to normalise.

pad(pad_width, axis_name=None, mode='constant', cval=0, method='linear')#

Pad the spatial dimensions of the map.

Parameters:
pad_width{sequence, array_like, int}

Number of pixels padded to the edges of each axis.

axis_namestr, optional

Which axis to downsample. By default, spatial axes are padded. Default is None.

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

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

cvalfloat, optional

Padding value when mode='consant'. Default is 0.

Returns:
mapMap

Padded map.

plot(ax=None, axis_name=None, **kwargs)[source]#

Plot the data contained in region map along the non-spatial axis.

Parameters:
axAxis, optional

Axis used for plotting. Default is None.

axis_namestr, optional

Which axis to plot on the x-axis. Extra axes will be plotted as additional lines. Default is None.

**kwargsdict

Keyword arguments passed to errorbar.

Returns:
axAxis

Axis used for plotting.

plot_grid(figsize=None, ncols=3, **kwargs)#

Plot map as a grid of subplots for non-spatial axes.

Parameters:
figsizetuple of int, optional

Figsize to plot on. Default is None.

ncolsint, optional

Number of columns to plot. Default is 3.

**kwargsdict, optional

Keyword arguments passed to WcsNDMap.plot.

Returns:
axesndarray of Axes

Axes grid.

plot_hist(ax=None, **kwargs)[source]#

Plot as histogram.

Parameters:
axaxis, optional

Axis instance to be used for the plot. Default is None.

**kwargsdict

Keyword arguments passed to hist.

Returns:
axAxis

Axis used for plotting.

plot_interactive()[source]#

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

Parameters:
rc_paramsdict, optional

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

**kwargsdict, optional

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)
plot_mask(ax=None, **kwargs)[source]#

Plot the mask as a shaded area in a xmin-xmax range.

Parameters:
axaxis, optional

Axis instance to be used for the plot. Default is None.

**kwargsdict

Keyword arguments passed to axvspan.

Returns:
axAxis

Axis used for plotting.

plot_region(ax=None, **kwargs)[source]#

Plot region.

Parameters:
axWCSAxes, optional

Axes to plot on. If no axes are given, the region is shown using the minimal equivalent WCS geometry. Default is None.

**kwargsdict

Keyword arguments forwarded to as_artist.

classmethod read(filename, format='gadf', ogip_column=None, hdu=None, checksum=False)[source]#

Read from file.

Parameters:
filenamepathlib.Path or str

Filename.

format{“gadf”, “ogip”, “ogip-arf”}

Which format to use. Default is “gadf”.

ogip_column{None, “COUNTS”, “QUALITY”, “BACKSCAL”}, optional

If format ‘ogip’ is chosen which table hdu column to read. Default is None.

hdustr, optional

Name or index of the HDU with the map data. Default is None.

checksumbool

If True checks both DATASUM and CHECKSUM cards in the file headers. Default is False.

Returns:
region_mapRegionNDMap

Region map.

reduce(axis_name, func=<ufunc 'add'>, keepdims=False, weights=None)#

Reduce map over a single non-spatial axis.

Parameters:
axis_name: str

The name of the axis to reduce over.

funcufunc, optional

Function to use for reducing the data. Default is numpy.add.

keepdimsbool, optional

If this is set to true, the axes which are summed over are left in the map with a single bin. Default is False.

weightsMap

Weights to be applied. The map should have the same geometry as this one. Default is None.

Returns:
map_outMap

Map with the given non-spatial axes reduced.

reduce_over_axes(func=<ufunc 'add'>, keepdims=False, axes_names=None, weights=None)#

Reduce map over non-spatial axes.

Parameters:
funcufunc, optional

Function to use for reducing the data. Default is numpy.add.

keepdimsbool, optional

If this is set to true, the axes which are summed over are left in the map with a single bin. Default is False.

axes_names: list, optional

Names of axis to reduce over. If None, all non-spatial axis will be reduced.

weightsMap, optional

Weights to be applied. The map should have the same geometry as this one. Default is None.

Returns:
map_outMap

Map with non-spatial axes reduced.

rename_axes(names, new_names)#

Rename the Map axes.

Parameters:
namesstr or list of str

Names of the axes.

new_namesstr or list of str

New names of the axes. The list must be of the same length as names).

Returns:
geomMap

Map with renamed axes.

reorder_axes(axes_names)#

Return a new map re-ordering the non-spatial axes.

Parameters:
axes_nameslist of str

The list of axes names in the required order.

Returns:
mapMap

The map with axes re-ordered.

reproject_by_image(geom, preserve_counts=False, precision_factor=10)#

Reproject each image of a ND map to input 2d geometry.

For large maps this method is faster than reproject_to_geom.

Parameters:
geomGeom

Target slice geometry in 2D.

preserve_countsbool, optional

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). Default is False.

precision_factorint, optional

Minimal factor between the bin size of the output map and the oversampled base map. Used only for the oversampling method. Default is 10.

Returns:
output_mapMap

Reprojected Map.

reproject_to_geom(geom, preserve_counts=False, precision_factor=10)#

Reproject map to input geometry.

Parameters:
geomGeom

Target Map geometry.

preserve_countsbool, optional

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). Default is False.

precision_factorint, optional

Minimal factor between the bin size of the output map and the oversampled base map. Used only for the oversampling method. Default is 10.

Returns:
output_mapMap

Reprojected Map.

resample(geom, weights=None, preserve_counts=True)#

Resample pixels to geom with given weights.

Parameters:
geomGeom

Target Map geometry.

weightsndarray, optional

Weights vector. It must have same shape as the data of the map. If set to None, weights will be set to 1. Default is None.

preserve_countsbool, optional

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). Default is True.

Returns:
resampled_mapMap

Resampled map.

resample_axis(axis, weights=None, ufunc=<ufunc 'add'>)#

Resample map to a new axis by grouping and reducing smaller bins by a given function ufunc.

By default, the map content are summed over the smaller bins. Other numpy.ufunc can be used, e.g. numpy.logical_and or numpy.logical_or.

Parameters:
axisMapAxis

New map axis.

weightsMap, optional

Array to be used as weights. The spatial geometry must be equivalent to other and additional axes must be broadcastable. If set to None, weights will be set to 1. Default is None.

ufuncufunc

Universal function to use to resample the axis. Default is numpy.add.

Returns:
mapMap

Map with resampled axis.

sample_coord(n_events, random_state=0)#

Sample position and energy of events.

Parameters:
n_eventsint

Number of events to sample.

random_state{int, ‘random-seed’, ‘global-rng’, RandomState}

Defines random number generator initialisation. Passed to get_random_state. Default is 0.

Returns:
coordsMapCoord object.

Sequence of coordinates and energies of the sampled events.

set_by_coord(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(idx, value)[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(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(slices)#

Slice sub map from map object.

Parameters:
slicesdict

Dictionary 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 dictionary are kept unchanged.

Returns:
map_outMap

Sliced map object.

split_by_axis(axis_name)#

Split a Map along an axis into multiple maps.

Parameters:
axis_namestr

Name of the axis to split.

Returns:
mapslist

A list of Map.

stack(other, weights=None, nan_to_num=True)[source]#

Stack other region map into map.

Parameters:
otherRegionNDMap

Other map to stack.

weightsRegionNDMap, optional

Array to be used as weights. The spatial geometry must be equivalent to other and additional axes must be broadcastable. Default is None.

nan_to_num: bool, optional

Non-finite values are replaced by zero if True. Default is True.

sum_over_axes(axes_names=None, keepdims=True, weights=None)#

Sum map values over all non-spatial axes.

Parameters:
axes_names: list of str

Names of the axis to reduce over. If None, all non-spatial axis will be summed over. Default is None.

keepdimsbool, optional

If this is set to true, the axes which are summed over are left in the map with a single bin. Default is True.

weightsMap, optional

Weights to be applied. The map should have the same geometry as this one. Default is None.

Returns:
map_outMap

Map with non-spatial axes summed over.

to_cube(axes)#

Append non-spatial axes to create a higher-dimensional Map.

This will result in a Map with 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. The data is reshaped onto this new geometry.

Parameters:
axeslist

Axes that will be appended to this Map. The axes should have only one bin.

Returns:
mapWcsNDMap

New map.

to_hdulist(format='gadf', hdu='SKYMAP', hdu_bands=None, hdu_region=None)[source]#

Convert to HDUList.

Parameters:
format{“gadf”, “ogip”, “ogip-sherpa”, “ogip-arf”, “ogip-arf-sherpa”}

Format specification. Default is “gadf”.

hdustr, optional

Name of the HDU with the map data, used for “gadf” format. Default is “SKYMAP”.

hdu_bandsstr, optional

Name or index of the HDU with the BANDS table, used for “gadf” format. Default is None.

hdu_regionstr, optional

Name or index of the HDU with the region table. Default is None.

Returns:
hdulistHDUList

HDU list.

to_region_nd_map(*args, **kwargs)[source]#

Return self.

to_table(format='gadf')[source]#

Convert to Table.

Data format specification: PHA.

Parameters:
format{“gadf”, “ogip”, “ogip-arf”, “ogip-arf-sherpa”}

Format specification. Default is “gadf”.

Returns:
tableTable

Table.

to_unit(unit)#

Convert map to a different unit.

Parameters:
unitstr or Unit

New unit.

Returns:
mapMap

Map with new unit and converted data.

upsample(factor, order=0, preserve_counts=True, axis_name=None)[source]#

Upsample the non-spatial dimension by a given factor.

By default, the first axes is upsampled.

Parameters:
factorint

Upsampling factor.

orderint, optional

Order of the interpolation used for upsampling. Default is 0.

preserve_countsbool, optional

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

axis_namestr, optional

Which axis to upsample. Default is “energy”.

Returns:
mapRegionNDMap

Upsampled region map.

write(filename, overwrite=False, format='gadf', hdu='SKYMAP', checksum=False)[source]#

Write map to file.

Parameters:
filenamepathlib.Path or str

Filename.

overwritebool, optional

Overwrite existing file. Default is False.

format{“gadf”, “ogip”, “ogip-sherpa”, “ogip-arf”, “ogip-arf-sherpa”}

Which format to use. Default is “gadf”.

hdustr, optional

Name of the HDU. Default is “SKYMAP”.

checksumbool, optional

When True adds both DATASUM and CHECKSUM cards to the headers written to the file. Default is False.