WcsNDMap¶
-
class
gammapy.maps.WcsNDMap(geom, data=None, dtype='float32', meta=None, unit='')[source]¶ Bases:
gammapy.maps.WcsMapRepresentation of a N+2D map using WCS with two spatial dimensions and N 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: geom :
WcsGeomWCS geometry object.
data :
ndarrayData array. If none then an empty array will be allocated.
dtype : str, optional
Data type, default is float32
meta :
OrderedDictDictionary to store meta data.
unit : str or
UnitThe map unit
Attributes Summary
dataData array ( ndarray)geomMap geometry ( MapGeom)metaMap meta ( OrderedDict)quantityMap data times unit ( Quantity)unitMap unit ( Unit)Methods Summary
coadd(map_in)Add the contents of map_into this map.convolve(kernel[, use_fft])Convolve map with a kernel. copy(**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(crop_width)Crop the spatial dimension of the map by removing a number of pixels from the edge of the map. cutout(position, width[, mode])Create a cutout around a given position. downsample(factor[, preserve_counts])Downsample the spatial dimension by a given factor. fill_by_coord(coords[, weights])Fill pixels at coordswith givenweights.fill_by_idx(idx[, weights])Fill pixels at idxwith givenweights.fill_by_pix(pix[, weights])Fill pixels at pixwith givenweights.from_geom(geom[, meta, data, map_type, unit])Generate an empty map from a MapGeominstance.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(coords)Return map values at the given map coordinates. get_by_idx(idx)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 interp_by_coord(coords[, interp, fill_value])Interpolate map values at the given map coordinates. interp_by_pix(pix[, interp, fill_value])Interpolate map values at the given pixel coordinates. iter_by_coord([buffersize])Iterate over elements of the map returning a tuple with values and map coordinates. iter_by_image()Iterate over image planes of the map returning a tuple with the image array and image plane index. iter_by_pix([buffersize])Iterate over elements of the map returning a tuple with values and pixel coordinates. make_hdu([hdu, hdu_bands, sparse, conv])Make a FITS HDU from this map. pad(pad_width[, mode, cval, order])Pad the spatial dimension of the map by extending the edge of the map by the given number of pixels. plot([ax, fig, add_cbar, stretch])Plot image on matplotlib WCS axes. plot_interactive([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(geom[, order, mode])Reproject this map to a different geometry. set_by_coord(coords, vals)Set pixels at coordswith givenvals.set_by_idx(idx, vals)Set pixels at idxwith givenvals.set_by_pix(pix, vals)Set pixels at pixwith givenvals.slice_by_idx(slices)Slice sub map from map object. smooth(width[, kernel])Smooth the image (works on a 2D image and returns a copy). sum_over_axes()Reduce to a 2D image by summing over non-spatial dimensions. to_hdulist([hdu, hdu_bands, sparse, conv])Convert to HDUList.upsample(factor[, order, preserve_counts])Upsample the spatial dimension by a given factor. write(filename[, overwrite])Write to a FITS file. Attributes Documentation
-
meta¶ Map meta (
OrderedDict)
Methods Documentation
-
coadd(map_in)¶ Add the contents of
map_into 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 :
MapInput map.
-
convolve(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 :
PSFKernelornumpy.ndarrayConvolution kernel.
use_fft : bool
kwargs : dict
Keyword arguments passed to
scipy.signal.fftconvolveorscipy.ndimage.convolve.Returns: map :
WcsNDMapConvolved map.
-
copy(**kwargs)¶ Copy map instance and overwrite given attributes, except for geometry.
Parameters: **kwargs : dict
Keyword arguments to overwrite in the map constructor.
Returns: copy :
MapCopied Map.
-
classmethod
create(map_type='wcs', npix=None, binsz=0.1, width=None, proj='CAR', coordsys='CEL', refpix=None, axes=None, skydir=None, dtype='float32', conv='gadf', 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.
npix : int 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.
width : float 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.
binsz : float 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.
skydir : tuple or
SkyCoordSky 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.
coordsys : {‘CEL’, ‘GAL’}, optional
Coordinate system, either Galactic (‘GAL’) or Equatorial (‘CEL’).
axes : list
List of non-spatial axes.
proj : string, optional
Any valid WCS projection type. Default is ‘CAR’ (cartesian).
refpix : tuple
Reference pixel of the projection. If None then this will be chosen to be center of the map.
dtype : str, optional
Data type, default is float32
conv : {‘fgst-ccube’,’fgst-template’,’gadf’}, optional
FITS format convention. Default is ‘gadf’.
meta :
OrderedDictDictionary to store meta data.
unit : str or
UnitThe unit of the map
Returns: map :
WcsMapA WCS map object.
-
crop(crop_width)[source]¶ Crop the spatial dimension of the map by removing a number of pixels from the edge of the map.
Parameters: crop_width : {sequence, array_like, int}
Number of pixels cropped from the edges of each axis. Defined analogously to
pad_withfromnumpy.pad.Returns: map :
MapCropped map.
-
cutout(position, width, mode='trim')[source]¶ Create a cutout around a given position.
Parameters: position :
SkyCoordCenter position of the cutout region.
width : tuple of
AngleAngular sizes of the region in (lon, lat) in that specific order. If only one value is passed, a square region is extracted.
mode : {‘trim’, ‘partial’, ‘strict’}
Mode option for Cutout2D, for details see
Cutout2D.Returns: cutout :
WcsNDMapCutout map
-
downsample(factor, preserve_counts=True)[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).
Returns: map :
MapDownsampled map.
-
fill_by_coord(coords, weights=None)¶ Fill pixels at
coordswith givenweights.Parameters: coords : tuple or
MapCoordCoordinate 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 :
ndarrayWeights vector. Default is weight of one.
-
fill_by_idx(idx, weights=None)[source]¶ Fill pixels at
idxwith givenweights.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 :
ndarrayWeights vector. Default is weight of one.
-
fill_by_pix(pix, weights=None)¶ Fill pixels at
pixwith givenweights.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 :
ndarrayWeights vector. Default is weight of one.
-
static
from_geom(geom, meta=None, data=None, map_type='auto', unit='')¶ Generate an empty map from a
MapGeominstance.Parameters: geom :
MapGeomMap geometry.
data :
numpy.ndarraydata array
meta :
OrderedDictDictionary 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
UnitData unit.
Returns: map_out :
MapMap object
-
classmethod
from_hdu(hdu, hdu_bands=None)[source]¶ Make a WcsNDMap object from a FITS HDU.
Parameters: hdu :
BinTableHDUorImageHDUThe map FITS HDU.
hdu_bands :
BinTableHDUThe BANDS table HDU.
-
classmethod
from_hdulist(hdu_list, hdu=None, hdu_bands=None)¶ Make a WcsMap object from a FITS HDUList.
Parameters: hdu_list :
HDUListHDU list containing HDUs for map data and bands.
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.
Returns: wcs_map :
WcsMapMap object
-
get_by_coord(coords)¶ Return map values at the given map coordinates.
Parameters: coords : tuple or
MapCoordCoordinate 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 :
ndarrayValues of pixels in the map. np.nan used to flag coords outside of map.
-
get_by_idx(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 :
ndarrayArray of pixel values. np.nan used to flag coordinate outside of map
-
get_by_pix(pix)¶ 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 :
ndarrayArray of pixel values. np.nan used to flag coordinates outside of map
-
get_image_by_coord(coords)¶ 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 :
MapMap with spatial dimensions only.
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: 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 :
MapMap with spatial dimensions only.
See also
-
get_image_by_pix(pix)¶ 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 :
MapMap with spatial dimensions only.
See also
-
interp_by_coord(coords, interp=None, fill_value=None)[source]¶ Interpolate map values at the given map coordinates.
Parameters: coords : tuple or
MapCoordCoordinate 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 :
ndarrayInterpolated pixel values.
-
interp_by_pix(pix, interp=None, fill_value=None)[source]¶ Interpolate map values at the given pixel coordinates.
-
iter_by_coord(buffersize=1)¶ Iterate over elements of the map returning a tuple with values and map coordinates.
Parameters: buffersize : int
Set the size of the buffer. The map will be returned in chunks of the given size.
Returns: val :
ndarrayMap values.
coords : tuple
Tuple of map coordinates.
-
iter_by_image()¶ Iterate over image planes of the map returning a tuple with the image array and image plane index.
The image plane index is in data order, so that the data array can be indexed directly. See Iterating on a Map for further information.
Returns: val :
ndarrayArray of image plane values.
idx : tuple
Index of image plane.
-
iter_by_pix(buffersize=1)¶ Iterate over elements of the map returning a tuple with values and pixel coordinates.
Parameters: buffersize : int
Set the size of the buffer. The map will be returned in chunks of the given size.
Returns: val :
ndarrayMap values.
pix : tuple
Tuple of pixel coordinates.
-
make_hdu(hdu='SKYMAP', hdu_bands=None, sparse=False, conv=None)¶ Make a FITS HDU from this map.
Parameters: hdu : str
The HDU extension name.
hdu_bands : str
The HDU extension name for BANDS table.
sparse : bool
Set INDXSCHM to SPARSE and sparsify the map by only writing pixels with non-zero amplitude.
Returns: hdu :
BinTableHDUorImageHDUHDU containing the map data.
-
pad(pad_width, mode='constant', cval=0, order=1)[source]¶ Pad the spatial dimension of the map by extending the edge of the map by the given number of pixels.
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 :
MapPadded map.
-
plot(ax=None, fig=None, add_cbar=False, stretch='linear', **kwargs)[source]¶ Plot image on matplotlib WCS axes.
Parameters: ax :
WCSAxes, optionalWCS axis object to plot on.
fig :
FigureFigure object.
add_cbar : bool
Add color bar?
stretch : str
Passed to
astropy.visualization.simple_norm.**kwargs : dict
Keyword arguments passed to
imshow.Returns: fig :
FigureFigure object.
ax :
WCSAxesWCS axis object
cbar :
Colorbaror NoneColorbar object.
-
plot_interactive(rc_params=None, **kwargs)¶ 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_EXTRA/datasets/vela_region/gll_iem_v05_rev1_cutout.fits") m.plot_interactive(cmap='gnuplot2')
If you would like to adjust the figure size you can use the
rc_paramsargument: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: filename : str or
PathName 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 :
MapMap object
-
reproject(geom, order=1, mode='interp')¶ 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 :
MapGeomGeometry 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 :
MapReprojected map.
-
set_by_coord(coords, vals)¶ Set pixels at
coordswith givenvals.Parameters: coords : tuple or
MapCoordCoordinate 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 :
ndarrayValues vector.
-
set_by_idx(idx, vals)[source]¶ Set pixels at
idxwith givenvals.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 :
ndarrayValues vector.
-
set_by_pix(pix, vals)¶ Set pixels at
pixwith givenvals.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 :
ndarrayValues vector.
-
slice_by_idx(slices)¶ Slice sub map from map object.
For usage examples, see Indexing and Slicing.
Parameters: slices : dict
Dict of axes names and integers or
sliceobject 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 :
MapSliced map object.
-
smooth(width, kernel='gauss', **kwargs)[source]¶ Smooth the image (works on a 2D image and returns a copy).
Parameters: width :
Quantityor floatSmoothing 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
kwargs : dict
Keyword arguments passed to
uniform_filter(‘box’),gaussian_filter(‘gauss’) orconvolve(‘disk’).Returns: image :
WcsNDMapSmoothed image (a copy, the original object is unchanged).
-
to_hdulist(hdu=None, hdu_bands=None, sparse=False, conv=None)¶ Convert to
HDUList.Parameters: 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.
sparse : bool
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_list :
HDUList
-
upsample(factor, order=0, preserve_counts=True)[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).
Returns: map :
MapUpsampled map.
-
write(filename, overwrite=False, **kwargs)¶ 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.
-