RegionNDMap¶
-
class
gammapy.maps.
RegionNDMap
(geom, data=None, dtype='float32', meta=None, unit='')[source]¶ Bases:
gammapy.maps.Map
Region ND map
- Parameters
- geom
RegionGeom
Region geometry object.
- data
ndarray
Data array. If none then an empty array will be allocated.
- dtypestr, optional
Data type, default is float32
- meta
dict
Dictionary to store meta data.
- unitstr or
Unit
The map unit
- geom
Attributes Summary
Data array (
ndarray
)Map geometry (
Geom
)Map meta (
dict
)Map data times unit (
Quantity
)Map unit (
Unit
)Methods Summary
apply_edisp
(edisp)Apply energy dispersion to map.
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, wcs])- Parameters
crop
()Crop the spatial dimensions of the map.
cutout
(*args, **kwargs)Return self
downsample
(factor[, preserve_counts, …])Downsample the spatial dimension by a given factor.
fill_by_coord
(coords[, weights])Fill pixels at
coords
with givenweights
.fill_by_idx
(idx[, weights])Fill pixels at
idx
with givenweights
.fill_by_pix
(pix[, weights])Fill pixels at
pix
with givenweights
.fill_events
(events)Fill event coordinates (
EventList
).from_geom
(geom[, meta, data, unit, dtype])Generate an empty map from a
Geom
instance.from_hdulist
(hdulist[, format, ogip_column])Create from
HDUList
.from_images
(images[, axis])Create Map from list of images and a non-spatial axis.
get_by_coord
(coords)Return map values at the given map coordinates.
get_by_idx
(idxs)Return map values at the given pixel indices.
get_by_pix
(pix)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
interp_by_coord
(coords[, interp])Interpolate map values at the given map coordinates.
interp_by_pix
(pix[, method, fill_value])Interpolate map values at the given pixel coordinates.
interp_to_geom
(geom[, preserve_counts])Interpolate map to input geometry.
Iterate over image planes of the map.
pad
()Pad the spatial dimensions of the map.
plot
([ax])Plot region map.
plot_grid
([figsize, ncols])Plot map as a grid of subplots for non-spatial axes
plot_hist
([ax])Plot as histogram.
Plot map with interactive widgets to explore the non spatial axes.
plot_region
([ax])Plot region
read
(filename[, format, ogip_column])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
resample_axis
(axis[, weights, ufunc])Resample map to a new axis binning by grouping over smaller bins and apply ufunc to the bin contents.
set_by_coord
(coords, vals)Set pixels at
coords
with givenvals
.set_by_idx
(idx, value)Set pixels at
idx
with givenvals
.set_by_pix
(pix, vals)Set pixels at
pix
with givenvals
.slice_by_idx
(slices)Slice sub map from map object.
stack
(other[, weights])Stack other region map into map.
sum_over_axes
([axes_names, keepdims, weights])To sum map values over all non-spatial axes.
to_cube
(axes)Append non-spatial axes to create a higher-dimensional Map.
to_hdulist
([format, ogip_column])Convert to
HDUList
.to_table
([format, ogip_column])Convert to
Table
.upsample
(factor[, preserve_counts, axis_name])Upsample the spatial dimension by a given factor.
write
(filename[, overwrite, format, ogip_column])Attributes Documentation
-
tag
= 'map'¶
Methods Documentation
-
apply_edisp
(edisp)¶ Apply energy dispersion to map. Requires energy axis.
- Parameters
- edisp
gammapy.irf.EDispKernel
Energy dispersion matrix
- edisp
- Returns
- map
WcsNDMap
Map with energy dispersion applied.
- map
-
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_in
Map
Input map.
- weights: `Map` or `~numpy.ndarray`
The weight factors while adding
- map_in
-
copy
(**kwargs)¶ Copy map instance and overwrite given attributes, except for geometry.
- Parameters
- **kwargsdict
Keyword arguments to overwrite in the map constructor.
- Returns
- copy
Map
Copied Map.
- copy
-
downsample
(factor, preserve_counts=True, axis_name='energy', weights=None)[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
- map
Map
Downsampled map.
- map
-
fill_by_coord
(coords, weights=None)¶ Fill pixels at
coords
with givenweights
.
-
fill_by_idx
(idx, weights=None)[source]¶ Fill pixels at
idx
with givenweights
.- 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.
- weights
ndarray
Weights vector. Default is weight of one.
-
fill_by_pix
(pix, weights=None)¶ Fill pixels at
pix
with givenweights
.- 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.
- weights
ndarray
Weights vector. Default is weight of one.
-
static
from_geom
(geom, meta=None, data=None, unit='', dtype='float32')¶ Generate an empty map from a
Geom
instance.- Parameters
- geom
Geom
Map geometry.
- data
numpy.ndarray
data array
- meta
dict
Dictionary to store meta data.
- unitstr or
Unit
Data unit.
- geom
- Returns
- map_out
Map
Map object
- map_out
-
classmethod
from_hdulist
(hdulist, format='ogip', ogip_column='COUNTS')[source]¶ Create from
HDUList
.- Parameters
- hdulist
HDUList
HDU list.
- format{“ogip”, “ogip-arf”}
Format specification
- ogip_column{“COUNTS”}
OGIP data format column
- hdulist
- Returns
- region_nd_map
RegionNDMap
Region map.
- region_nd_map
-
classmethod
from_images
(images, axis=None)¶ Create Map from list of images and a non-spatial axis.
If the images have a non-spatial axis of length 1 a new axes is generated from by merging the individual axes. The image geometries must be aligned.
-
get_by_coord
(coords)¶ Return map values at the given map coordinates.
-
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
- vals
ndarray
Array of pixel values. np.nan used to flag coordinate outside of map
- vals
-
get_by_pix
(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
- vals
ndarray
Array of pixel values. np.nan used to flag coordinates outside of map
- vals
-
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. 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.
- map_out
See also
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_out
Map
Map with spatial dimensions only.
- map_out
See also
-
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_out
Map
Map with spatial dimensions only.
- map_out
See also
-
interp_by_coord
(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.
- coordstuple or
- Returns
- vals
ndarray
Interpolated pixel values.
- vals
-
interp_by_pix
(pix, method='linear', fill_value=None)[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.
- 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
- vals
ndarray
Interpolated pixel values.
- vals
-
interp_to_geom
(geom, preserve_counts=False, **kwargs)¶ Interpolate map to input geometry.
- Parameters
- geom
Geom
Target Map geometry
- 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)
- **kwargsdict
Keyword arguments passed to
Map.interp_by_coord
- geom
- Returns
- interp_map
Map
Interpolated Map
- interp_map
-
iter_by_image
()¶ Iterate over image planes of the map.
This is a generator yielding
(data, idx)
tuples, wheredata
is anumpy.ndarray
view of the image plane data, andidx
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
()[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
- map
Map
Padded map.
- map
-
plot
(ax=None, **kwargs)[source]¶ Plot region map.
- Parameters
- ax
Axis
Axis used for plotting
- **kwargsdict
Keyword arguments passed to
errorbar
- ax
- Returns
- ax
Axis
Axis used for plotting
- ax
-
plot_grid
(figsize=None, ncols=3, **kwargs)¶ Plot map as a grid of subplots for non-spatial axes
- Parameters
- figsizetuple of int
Figsize to plot on
- ncolsint
Number of columns to plot
- **kwargsdict
Keyword arguments passed to
Map.plot
.
- Returns
- axes
ndarray
ofAxes
Axes grid
- axes
-
plot_interactive
()[source]¶ 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)
-
plot_region
(ax=None, **kwargs)[source]¶ Plot region
- Parameters
- ax
WCSAxes
Axes to plot on.
- **kwargsdict
Keyword arguments forwarded to
as_artist
- ax
-
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
- func
ufunc
Function to use for reducing the data.
- keepdimsbool, optional
If this is set to true, the axes which are summed over are left in the map with a single bin
- weights
Map
Weights to be applied.
- Returns
- map_out
Map
Map with the given non-spatial axes reduced
- map_out
-
reduce_over_axes
(func=<ufunc 'add'>, keepdims=False, axes_names=None, weights=None)¶ Reduce map over non-spatial axes
- Parameters
- func
ufunc
Function to use for reducing the data.
- keepdimsbool, optional
If this is set to true, the axes which are summed over are left in the map with a single bin
- axes_names: list
Names of MapAxis to reduce over If None, all will reduced
- weights
Map
Weights to be applied.
- func
- Returns
- map_out
Map
Map with non-spatial axes reduced
- map_out
-
resample_axis
(axis, weights=None, ufunc=<ufunc 'add'>)¶ Resample map to a new axis binning by grouping over smaller bins and apply ufunc to the bin contents.
By default, the map content are summed over the smaller bins. Other numpy ufunc can be used, e.g. np.logical_and, np.logical_or
-
set_by_coord
(coords, vals)¶ Set pixels at
coords
with givenvals
.
-
set_by_idx
(idx, value)[source]¶ Set pixels at
idx
with givenvals
.- 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.
- vals
ndarray
Values vector.
-
set_by_pix
(pix, vals)¶ Set pixels at
pix
with givenvals
.- 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.
- vals
ndarray
Values vector.
-
slice_by_idx
(slices)¶ Slice sub map from map object.
For usage examples, see Indexing and Slicing.
-
stack
(other, weights=None)[source]¶ Stack other region map into map.
- Parameters
- other
RegionNDMap
Other map to stack
- weights
RegionNDMap
Array to be used as weights. The spatial geometry must be equivalent to
other
and additional axes must be broadcastable.
- other
-
sum_over_axes
(axes_names=None, keepdims=True, weights=None)¶ To sum map values over all non-spatial axes.
- Parameters
- keepdimsbool, optional
If this is set to true, the axes which are summed over are left in the map with a single bin
- axes_names: list of str
Names of MapAxis to reduce over. If None, all will summed over
- weights
Map
Weights to be applied. The Map should have the same geometry.
- Returns
- map_out
Map
Map with non-spatial axes summed over
- map_out
-
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 the new geometry
- Parameters
- axeslist
Axes that will be appended to this Map. The axes should have only one bin
- Returns
- map
WcsNDMap
new map
- map
-
to_hdulist
(format='ogip', ogip_column='COUNTS')[source]¶ Convert to
HDUList
.- Parameters
- format{“ogip”, “ogip-sherpa”}
Format specification
- ogip_column{“COUNTS”, “SPECRESP”}
Ogip column format
- Returns
- hdulist
HDUList
HDU list
- hdulist
-
to_table
(format='ogip', ogip_column='COUNTS')[source]¶ Convert to
Table
.Data format specification: PHA
- Parameters
- format{“ogip”, “ogip-sherpa”}
Format specification
- ogip_column{“COUNTS”, “SPECRESP”}
Ogip column format
- Returns
- table
Table
Table
- table
-
upsample
(factor, preserve_counts=True, axis_name='energy')[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
- map
Map
Upsampled map.
- map