WcsNDMap¶
-
class
gammapy.maps.
WcsNDMap
(geom, data=None, dtype='float32', meta=None, unit='')[source]¶ Bases:
gammapy.maps.WcsMap
HEALPix map with any number of non-spatial dimensions.
This class uses an ND numpy array to store map values. For maps with non-spatial dimensions and variable pixel size it will allocate an array with dimensions commensurate with the largest image plane.
- Parameters
Attributes Summary
Data array (
ndarray
)Map geometry (
Geom
)Map meta (
dict
)Map data times unit (
Quantity
)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.convolve
(self, kernel[, use_fft])Convolve map with a kernel.
copy
(self, \*\*kwargs)Copy map instance and overwrite given attributes, except for geometry.
create
([map_type, npix, binsz, width, proj, …])Factory method to create an empty WCS map.
crop
(self, crop_width)Crop the spatial dimensions of the map.
cutout
(self, position, width[, mode])Create a cutout around a given position.
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 givenweights
.fill_by_idx
(self, idx[, weights])Fill pixels at
idx
with givenweights
.fill_by_pix
(self, pix[, weights])Fill pixels at
pix
with givenweights
.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 WcsNDMap object from a FITS HDU.
from_hdulist
(hdu_list[, hdu, hdu_bands])Make a WcsMap 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
get_spectrum
(self[, region, func])Extract spectrum in a given region.
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.
make_hdu
(self[, hdu, hdu_bands, sparse, conv])Make a FITS HDU from this map.
pad
(self, pad_width[, mode, cval, order])Pad the spatial dimensions of the map.
plot
(self[, ax, fig, add_cbar, stretch])Plot image on matplotlib WCS axes.
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.
reduce
(self, axis[, func, keepdims])Reduce map over a single non-spatial axis
reduce_over_axes
(self[, func, keepdims, axes])Reduce map over non-spatial axes
sample_coord
(self, n_events[, random_state])Sample position and energy of events.
set_by_coord
(self, coords, vals)Set pixels at
coords
with givenvals
.set_by_idx
(self, idx, vals)Set pixels at
idx
with givenvals
.set_by_pix
(self, pix, vals)Set pixels at
pix
with givenvals
.slice_by_idx
(self, slices)Slice sub map from map object.
smooth
(self, width[, kernel])Smooth the map.
stack
(self, other[, weights])Stack cutout into map.
sum_over_axes
(self[, axes, keepdims])To sum map values over all non-spatial axes.
to_hdulist
(self[, hdu, hdu_bands, sparse, conv])Convert to
HDUList
.upsample
(self, factor[, order, …])Upsample the spatial dimension by a given factor.
write
(self, filename[, overwrite])Write to a FITS file.
Attributes Documentation
Methods Documentation
-
apply_edisp
(self, 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
(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_in
Map
Input map.
- weights: `Map` or `~numpy.ndarray`
The weight factors while adding
- map_in
-
convolve
(self, kernel, use_fft=True, **kwargs)[source]¶ Convolve map with a kernel.
If the kernel is two dimensional, it is applied to all image planes likewise. If the kernel is higher dimensional it must match the map in the number of dimensions and the corresponding kernel is selected for every image plane.
- Parameters
- kernel
PSFKernel
ornumpy.ndarray
Convolution kernel.
- use_fftbool
- kwargsdict
Keyword arguments passed to
scipy.signal.fftconvolve
orscipy.ndimage.convolve
.
- kernel
- Returns
- map
WcsNDMap
Convolved map.
- map
-
copy
(self, **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
-
classmethod
create
(map_type='wcs', npix=None, binsz=0.1, width=None, proj='CAR', frame='icrs', refpix=None, axes=None, skydir=None, dtype='float32', meta=None, unit='')¶ Factory method to create an empty WCS map.
- Parameters
- map_type{‘wcs’, ‘wcs-sparse’}
Map type. Selects the class that will be used to instantiate the map.
- npixint or tuple or list
Width of the map in pixels. A tuple will be interpreted as parameters for longitude and latitude axes. For maps with non-spatial dimensions, list input can be used to define a different map width in each image plane. This option supersedes width.
- widthfloat or tuple or list
Width of the map in degrees. A tuple will be interpreted as parameters for longitude and latitude axes. For maps with non-spatial dimensions, list input can be used to define a different map width in each image plane.
- binszfloat or tuple or list
Map pixel size in degrees. A tuple will be interpreted as parameters for longitude and latitude axes. For maps with non-spatial dimensions, list input can be used to define a different bin size in each image plane.
- 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.
- frame{“icrs”, “galactic”}, optional
Coordinate system, either Galactic (“galactic”) or Equatorial (“icrs”).
- axeslist
List of non-spatial axes.
- projstring, optional
Any valid WCS projection type. Default is ‘CAR’ (cartesian).
- refpixtuple
Reference pixel of the projection. If None then this will be chosen to be center of the map.
- dtypestr, optional
Data type, default is float32
- conv{‘fgst-ccube’,’fgst-template’,’gadf’}, optional
FITS format convention. Default is ‘gadf’.
- meta
dict
Dictionary to store meta data.
- unitstr or
Unit
The unit of the map
- Returns
- map
WcsMap
A WCS map object.
- map
-
cutout
(self, position, width, mode='trim')[source]¶ Create a cutout around a given position.
- Parameters
- Returns
- cutout
WcsNDMap
Cutout map
- cutout
-
downsample
(self, factor, preserve_counts=True, axis=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
(self, coords, weights=None)¶ Fill pixels at
coords
with givenweights
.
-
fill_by_idx
(self, 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
(self, 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, map_type='auto', 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.
- 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.
- geom
- Returns
- map_out
Map
Map object
- map_out
-
classmethod
from_hdu
(hdu, hdu_bands=None)[source]¶ Make a WcsNDMap object from a FITS HDU.
- Parameters
- hdu
BinTableHDU
orImageHDU
The map FITS HDU.
- hdu_bands
BinTableHDU
The BANDS table HDU.
- hdu
-
classmethod
from_hdulist
(hdu_list, hdu=None, hdu_bands=None)¶ Make a WcsMap object from a FITS HDUList.
-
get_by_coord
(self, coords)¶ Return map values at the given map coordinates.
-
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
- vals
ndarray
Array of pixel values. np.nan used to flag coordinate outside of map
- vals
-
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
- vals
ndarray
Array of pixel values. np.nan used to flag coordinates outside of map
- vals
-
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_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
(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_out
Map
Map with spatial dimensions only.
- map_out
See also
-
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_out
Map
Map with spatial dimensions only.
- map_out
See also
-
get_spectrum
(self, region=None, func=<function nansum at 0x10a7d0ae8>)[source]¶ Extract spectrum in a given region.
The spectrum can be computed by summing (or, more generally, applying
func
) along the spatial axes in each energy bin. This occurs only inside theregion
, which by default is assumed to be the whole spatial extension of the map.- Parameters
- region: `~regions.Region`
Region (pixel or sky regions accepted).
- funcnumpy.ufunc
Function to reduce the data.
- Returns
- spectrum
RegionNDMap
Spectrum in the given region.
- spectrum
-
interp_by_coord
(self, coords, interp=None, fill_value=None)[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
(self, pix, interp=None, fill_value=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, 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.
-
make_hdu
(self, hdu='SKYMAP', hdu_bands=None, sparse=False, conv=None)¶ Make a FITS HDU from this map.
- 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.
- Returns
- hdu
BinTableHDU
orImageHDU
HDU containing the map data.
- hdu
-
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
- map
Map
Padded map.
- map
-
plot
(self, ax=None, fig=None, add_cbar=False, stretch='linear', **kwargs)[source]¶ Plot image on matplotlib WCS axes.
- Parameters
- ax
WCSAxes
, optional WCS axis object to plot on.
- fig
Figure
Figure object.
- add_cbarbool
Add color bar?
- stretchstr
Passed to
astropy.visualization.simple_norm
.- **kwargsdict
Keyword arguments passed to
imshow
.
- ax
- Returns
-
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.
- filenamestr or
- Returns
- map_out
Map
Map object
- map_out
-
reduce
(self, axis, func=<ufunc 'add'>, keepdims=False)[source]¶ Reduce map over a single non-spatial axis
- Parameters
- axis: 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
- Returns
- map_out
WcsNDMap
Map with the given non-spatial axes reduced
- map_out
-
reduce_over_axes
(self, func=<ufunc 'add'>, keepdims=False, axes=None)[source]¶ 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: list
Names of MapAxis to reduce over If None, all will reduced
- func
- Returns
- map_out
WcsNDMap
Map with non-spatial axes reduced
- map_out
-
sample_coord
(self, n_events, random_state=0)[source]¶ 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
.
- Returns
- coords
MapCoord
object. Sequence of coordinates and energies of the sampled events.
- coords
-
set_by_coord
(self, coords, vals)¶ Set pixels at
coords
with givenvals
.
-
set_by_idx
(self, idx, vals)[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
(self, 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
(self, slices)¶ Slice sub map from map object.
For usage examples, see Indexing and Slicing.
-
smooth
(self, width, kernel='gauss', **kwargs)[source]¶ Smooth the map.
Iterates over 2D image planes, processing one at a time.
- Parameters
- width
Quantity
, str or float Smoothing width given as quantity or float. If a float is given it interpreted as smoothing width in pixels. If an (angular) quantity is given it converted to pixels using
geom.wcs.wcs.cdelt
. It corresponds to the standard deviation in case of a Gaussian kernel, the radius in case of a disk kernel, and the side length in case of a box kernel.- kernel{‘gauss’, ‘disk’, ‘box’}
Kernel shape
- kwargsdict
Keyword arguments passed to
uniform_filter
(‘box’),gaussian_filter
(‘gauss’) orconvolve
(‘disk’).
- width
- Returns
- image
WcsNDMap
Smoothed image (a copy, the original object is unchanged).
- image
-
sum_over_axes
(self, axes=None, keepdims=False)[source]¶ 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: list
Names of MapAxis to reduce over If None, all will summed over
- Returns
- map_out
WcsNDMap
Map with non-spatial axes summed over
- map_out
-
to_hdulist
(self, hdu=None, hdu_bands=None, sparse=False, conv='gadf')¶ Convert to
HDUList
.- Parameters
- hdustr
Name or index of the HDU with the map data.
- hdu_bandsstr
Name or index of the HDU with the BANDS table.
- sparsebool
Sparsify the map by only writing pixels with non-zero amplitude.
- conv{‘gadf’, ‘fgst-ccube’,’fgst-template’}
FITS format convention.
- Returns
- hdu_list
HDUList
- hdu_list
-
upsample
(self, factor, order=0, preserve_counts=True, axis=None)[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
-
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.