HpxNDMap#
- class gammapy.maps.HpxNDMap(geom, data=None, dtype='float32', meta=None, unit='')[source]#
Bases:
gammapy.maps.hpx.core.HpxMap
HEALPix map with any number of non-spatial dimensions.
This class uses a N+1D numpy array to represent the sequence of HEALPix image planes. Following the convention of WCS-based maps this class uses a column-wise ordering for the data array with the spatial dimension being tied to the last index of the array.
- Parameters
Attributes Summary
Data array (
ndarray
)Map geometry (
Geom
)Whether map is mask with bool dtype
Map meta (
dict
)Map data times unit (
Quantity
)Map unit (
Unit
)Methods Summary
apply_edisp
(edisp)Deprecated since version v1.1.
coadd
(map_in[, weights])Add the contents of
map_in
to this map.convolve
(kernel[, convolution_method])Convolve map with a WCS kernel.
convolve_full
(kernel)Convolve map with a symmetrical WCS kernel.
convolve_wcs
(kernel, **kwargs)Convolve map with a WCS kernel.
copy
(**kwargs)Copy map instance and overwrite given attributes, except for geometry.
create
([nside, binsz, nest, map_type, ...])Factory method to create an empty HEALPix map.
crop
(crop_width)Crop the spatial dimensions of the map.
cumsum
(axis_name)Compute cumulative sum along a given axis
cutout
(position, width, *args, **kwargs)Create a cutout around a given position.
dot
(other)Apply dot product with the input map.
downsample
(factor[, preserve_counts, axis_name])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[, weights])Fill event coordinates (
EventList
).from_geom
(geom[, meta, data, unit, dtype])Generate an empty map from a
Geom
instance.from_hdu
(hdu[, hdu_bands, format, colname])Make a HpxNDMap object from a FITS HDU.
from_hdulist
(hdu_list[, hdu, hdu_bands, ...])Make a HpxMap object from a FITS HDUList.
from_stack
(maps[, axis, axis_name])Create Map from list of images and a non-spatial axis.
from_wcs_tiles
(wcs_tiles[, nest])Create HEALPix map from WCS tiles.
get_by_coord
(coords[, fill_value])Return map values at the given map coordinates.
get_by_idx
(idx)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
([region, func, weights])Extract spectrum in a given region.
integral
(axis_name, coords, **kwargs)Compute integral along a given axis
interp_by_coord
(coords[, method, fill_value])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.
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_image
([keepdims])Iterate over image planes of a map.
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
([method, ax, normalize, proj, ...])Quickplot method.
plot_grid
([figsize, ncols])Plot map as a grid of subplots for non-spatial axes
plot_interactive
([rc_params])Plot map with interactive widgets to explore the non spatial axes.
plot_mask
([method, ax, proj, oversample, ...])Plot the mask as a shaded area
read
(filename[, hdu, hdu_bands, map_type, ...])Read a map from a FITS 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 order.
reproject_to_geom
(geom[, preserve_counts, ...])Reproject map to input geometry.
resample
(geom[, weights, preserve_counts])Resample pixels to
geom
with givenweights
.resample_axis
(axis[, weights, ufunc])Resample map to a new axis by grouping and reducing smaller bins by a given ufunc
sample_coord
(n_events[, random_state])Sample position and energy of events.
set_by_coord
(coords, vals)Set pixels at
coords
with givenvals
.set_by_idx
(idx, vals)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.
smooth
(width[, kernel])Smooth the map.
split_by_axis
(axis_name)Split a Map along an axis into multiple maps.
stack
(other[, weights, nan_to_num])Stack cutout 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_hdu
([hdu, hdu_bands, sparse, format])Make a FITS HDU with input data.
to_hdulist
([hdu, hdu_bands, sparse, format])Convert to
HDUList
.to_nside
(nside[, preserve_counts])Upsample or downsample the map to a given nside
to_region_nd_map
(region[, func, weights, method])Get region ND map in a given region.
Return a new map with the opposite scheme (ring or nested).
to_unit
(unit)Convert map to different unit
to_wcs
([sum_bands, normalize, proj, ...])Make a WCS object and convert HEALPIX data into WCS projection.
to_wcs_tiles
([nside_tiles, margin, method, ...])Convert HpxNDMap to a list of WCS tiles
upsample
(factor[, order, preserve_counts, ...])Upsample the spatial dimension by a given factor.
write
(filename[, overwrite])Write to a FITS file.
Attributes Documentation
- is_mask#
Whether map is mask with bool dtype
- tag = 'map'#
Methods Documentation
- apply_edisp(edisp)#
Deprecated since version v1.1: The apply_edisp function is deprecated and may be removed in a future version. Use gammapy.datasets.apply_edisp instead.
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
- convolve(kernel, convolution_method='wcs-tan', **kwargs)[source]#
Convolve map with a WCS kernel.
It projects the map into a WCS geometry, convolves with a WCS kernel and projects back into the initial Healpix geometry.
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
Convolution kernel. The pixel size must be upsampled by a factor 2 or bigger with respect to the input map to prevent artifacts in the projection.
- convolution_methodstr
Supported methods are : ‘wcs-tan’: project on WCS geometry and convolve with WCS kernel. See
convolve_wcs
.- **kwargsdict
Keyword arguments passed to
convolve
.
- kernel
- Returns
- map
HpxNDMap
Convolved map.
- map
- convolve_full(kernel)[source]#
Convolve map with a symmetrical WCS kernel.
It extracts the radial profile of the kernel (assuming radial symmetry) and convolves via
hp.sphtfunc.smoothing
. Since no projection is applied, this is suited for full-sky and large maps.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.
- convolve_wcs(kernel, **kwargs)[source]#
Convolve map with a WCS kernel.
It projects the map into a WCS geometry, convolves with a WCS kernel and projects back into the initial Healpix geometry.
If the kernel is two dimensional, it is applied to all image planes likewise. If the kernel is higher dimensional should either match the map in the number of dimensions or the map must be an image (no non-spatial axes). In that case, the corresponding kernel is selected and applied to every image plane or to the single input image respectively.
- 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
- classmethod create(nside=None, binsz=None, nest=True, map_type='hpx', frame='icrs', data=None, skydir=None, width=None, dtype='float32', region=None, axes=None, meta=None, unit='')#
Factory method to create an empty HEALPix map.
- Parameters
- nsideint or
ndarray
HEALPix NSIDE parameter. This parameter sets the size of the spatial pixels in the map.
- binszfloat or
ndarray
Approximate pixel size in degrees. An NSIDE will be chosen that corresponds to a pixel size closest to this value. This option is superseded by nside.
- nestbool
True for HEALPix “NESTED” indexing scheme, False for “RING” scheme.
- frame{“icrs”, “galactic”}, optional
Coordinate system, either Galactic (“galactic”) or Equatorial (“icrs”).
- 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.
- map_type{‘hpx’, ‘hpx-sparse’}
Map type. Selects the class that will be used to instantiate the map.
- widthfloat
Diameter of the map in degrees. If None then an all-sky geometry will be created.
- axeslist
List of
MapAxis
objects for each non-spatial dimension.- meta
dict
Dictionary to store meta data.
- unitstr or
Unit
The map unit
- nsideint or
- Returns
- map
HpxMap
A HPX map object.
- map
- cumsum(axis_name)#
Compute cumulative sum along a given axis
- Parameters
- axis_namestr
Along which axis to integrate.
- Returns
- cumsum
Map
Map with cumulative sum
- cumsum
- 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
- other
RegionNDMap
Map to apply the dot product to. It must share a unique non-spatial MapAxis with the current Map.
- other
- Returns
- map
Map
Map with dot product applied.
- map
- downsample(factor, preserve_counts=True, axis_name=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).
- axis_namestr
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_hdu(hdu, hdu_bands=None, format=None, colname=None)[source]#
Make a HpxNDMap object from a FITS HDU.
- Parameters
- hdu
BinTableHDU
The FITS HDU
- hdu_bands
BinTableHDU
The BANDS table HDU
- formatstr, optional
FITS convention. If None the format is guessed. The following formats are supported:
“gadf”
“fgst-ccube”
“fgst-ltcube”
“fgst-bexpcube”
“fgst-srcmap”
“fgst-template”
“fgst-srcmap-sparse”
“galprop”
“galprop2”
- colnamestr, optional
Data column name to be used for the HEALPix map.
- hdu
- Returns
- map
HpxMap
HEALPix map
- map
- classmethod from_hdulist(hdu_list, hdu=None, hdu_bands=None, format=None, colname=None)#
Make a HpxMap object from a FITS HDUList.
- Parameters
- hdu_list
HDUList
HDU list containing HDUs for map data and bands.
- hdustr
Name or index of the HDU with the map data. If None then the method will try to load map data from the first BinTableHDU in the file.
- hdu_bandsstr
Name or index of the HDU with the BANDS table.
- formatstr, optional
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). The following formats are supported:
“gadf” (default)
“fgst-ccube”
“fgst-ltcube”
“fgst-bexpcube”
“fgst-srcmap”
“fgst-template”
“fgst-srcmap-sparse”
“galprop”
“galprop2”
- hdu_list
- Returns
- hpx_map
HpxMap
Map object
- hpx_map
- classmethod from_stack(maps, axis=None, axis_name=None)#
Create Map from list of images and a non-spatial axis.
The image geometries must be aligned, except for the axis that is stacked.
- Parameters
- Returns
- map
Map
Map with additional non-spatial axis.
- 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 of the projection footprint
- coordstuple or
- Returns
- vals
ndarray
Values of pixels in the map. np.nan used to flag coords outside of map.
- vals
- get_by_idx(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(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 of the projection footprint
- 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
- get_spectrum(region=None, func=<function nansum>, weights=None)#
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.func
Function to reduce the data. Default is np.nansum. For a boolean Map, use np.any or np.all.
- weights
WcsNDMap
Array to be used as weights. The geometry must be equivalent.
- Returns
- spectrum
RegionNDMap
Spectrum in the given region.
- spectrum
- integral(axis_name, coords, **kwargs)#
Compute integral along a given axis
This method uses interpolation of the cumulative sum.
- interp_by_coord(coords, method='linear', 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.
- method{“linear”, “nearest”}
Method to interpolate data values. By default linear interpolation is performed.
- 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=None, fill_value=None)[source]#
Interpolate map values at the given pixel coordinates.
- interp_to_geom(geom, preserve_counts=False, fill_value=0, **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
- is_allclose(other, rtol_axes=0.001, atol_axes=1e-06, **kwargs)#
Compare two Maps for close equivalency
- Parameters
- other
gammapy.maps.Map
The Map to compare against
- rtol_axesfloat
Relative tolerance for the axes comparison.
- atol_axesfloat
Relative tolerance for the axes comparison.
- **kwargsdict
keywords passed to
numpy.allclose
- other
- Returns
- is_allclosebool
Whether the Map is all close.
- iter_by_axis(axis_name, keepdims=False)#
“Iterate over a given axis
- Yields
- map
Map
Map iteration.
- map
See also
iter_by_image
iterate by image returning a map
- iter_by_image(keepdims=False)#
Iterate over image planes of a map.
- Parameters
- keepdimsbool
Keep dimensions.
- Yields
- map
Map
Map iteration.
- map
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 anumpy.ndarray
view of the image plane data, andidx
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.
- normalize(axis_name=None)#
Normalise data in place along a given axis.
- Parameters
- axis_namestr
Along which axis to normalize.
- 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
Which axis to downsample. By default spatial axes are padded.
- 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’.
- Returns
- map
Map
Padded map.
- map
- plot(method='raster', ax=None, normalize=False, proj='AIT', oversample=2, width_pix=1000, **kwargs)[source]#
Quickplot method.
This will generate a visualization of the map by converting to a rasterized WCS image (method=’raster’) or drawing polygons for each pixel (method=’poly’).
- Parameters
- method{‘raster’,’poly’}
Method for mapping HEALPix pixels to a two-dimensional image. Can be set to ‘raster’ (rasterization to cartesian image plane) or ‘poly’ (explicit polygons for each pixel). WARNING: The ‘poly’ method is much slower than ‘raster’ and only suitable for maps with less than ~10k pixels.
- projstring, optional
Any valid WCS projection type.
- oversamplefloat
Oversampling factor for WCS map. This will be the approximate ratio of the width of a HPX pixel to a WCS pixel. If this parameter is None then the width will be set from
width_pix
.- width_pixint
Width of the WCS geometry in pixels. The pixel size will be set to the number of pixels satisfying
oversample
orwidth_pix
whichever is smaller. If this parameter is None then the width will be set fromoversample
.- **kwargsdict
Keyword arguments passed to
imshow
.
- Returns
- ax
WCSAxes
WCS axis object
- 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
WcsNDMap.plot
.
- Returns
- axes
ndarray
ofAxes
Axes grid
- axes
- plot_interactive(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)
- plot_mask(method='raster', ax=None, proj='AIT', oversample=2, width_pix=1000, **kwargs)[source]#
Plot the mask as a shaded area
- Parameters
- method{‘raster’,’poly’}
Method for mapping HEALPix pixels to a two-dimensional image. Can be set to ‘raster’ (rasterization to cartesian image plane) or ‘poly’ (explicit polygons for each pixel). WARNING: The ‘poly’ method is much slower than ‘raster’ and only suitable for maps with less than ~10k pixels.
- projstring, optional
Any valid WCS projection type.
- oversamplefloat
Oversampling factor for WCS map. This will be the approximate ratio of the width of a HPX pixel to a WCS pixel. If this parameter is None then the width will be set from
width_pix
.- width_pixint
Width of the WCS geometry in pixels. The pixel size will be set to the number of pixels satisfying
oversample
orwidth_pix
whichever is smaller. If this parameter is None then the width will be set fromoversample
.- **kwargsdict
Keyword arguments passed to
imshow
.
- Returns
- ax
WCSAxes
WCS axis object
- ax
- static read(filename, hdu=None, hdu_bands=None, map_type='auto', format=None, colname=None)#
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’, ‘region’}
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.
- colnamestr, optional
data column name to be used of healix map.
- filenamestr or
- Returns
- map_out
Map
Map object
- map_out
- reduce(axis_name, func=<ufunc 'add'>, keepdims=False, weights=None)#
Reduce map over a single non-spatial axis
- Parameters
- 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
- Returns
- map_out
Map
Map with non-spatial axes reduced
- map_out
- rename_axes(names, new_names)#
Rename the Map axes.
- Parameters
- nameslist or str
Names of the axes.
- new_nameslist or str
New names of the axes (list must be of same length than
names
).
- Returns
- geom
Map
Renamed Map.
- geom
- reorder_axes(axes_names)#
Return a new map re-ordering the non-spatial axes order.
- Parameters
- axes_nameslist of str
the list of axes names in the required order
- Returns
- map
Map
the map with axes re-ordered
- map
- reproject_to_geom(geom, preserve_counts=False, precision_factor=10)#
Reproject 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)
- precision_factorint
Minimal factor between the bin size of the output map and the oversampled base map. Used only for the oversampling method.
- geom
- Returns
- output_map
Map
Reprojected Map
- output_map
- resample(geom, weights=None, preserve_counts=True)#
Resample pixels to
geom
with givenweights
.- Parameters
- geom
Geom
Target Map geometry
- weights
ndarray
Weights vector. Default is weight of one. Must have same shape as the data of the map.
- 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)
- geom
- Returns
- resampled_map
Map
Resampled map
- 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 ufunc
By default, the map content are summed over the smaller bins. Other numpy ufunc can be used, e.g.
numpy.logical_and
ornumpy.logical_or
.
- sample_coord(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(coords, vals)#
Set pixels at
coords
with givenvals
.
- set_by_idx(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(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.
- smooth(width, kernel='gauss')[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
healpy.nside2resol
. It corresponds to the standard deviation in case of a Gaussian kernel, and the radius in case of a disk kernel.- kernel{‘gauss’, ‘disk’}
Kernel shape
- width
- Returns
- image
HpxNDMap
Smoothed image (a copy, the original object is unchanged).
- image
- 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
- 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_hdu(hdu=None, hdu_bands=None, sparse=False, format=None)#
Make a FITS HDU with input data.
- 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.
- format{‘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_out
BinTableHDU
orImageHDU
Output HDU containing map data.
- hdu_out
- to_hdulist(hdu='SKYMAP', hdu_bands=None, sparse=False, format='gadf')#
Convert to
HDUList
.- 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.
- formatstr, optional
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). The following formats are supported:
“gadf” (default)
“fgst-ccube”
“fgst-ltcube”
“fgst-bexpcube”
“fgst-srcmap”
“fgst-template”
“fgst-srcmap-sparse”
“galprop”
“galprop2”
- Returns
- hdu_list
HDUList
- hdu_list
- to_nside(nside, preserve_counts=True)[source]#
Upsample or downsample the map to a given nside
- Parameters
- nsideint
Nside
- 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).
- Returns
- geom
HpxNDMap
Healpix map with new nside.
- geom
- to_region_nd_map(region, func=<function nansum>, weights=None, method='nearest')[source]#
Get region ND map in a given region.
By default the whole map region is considered.
- Parameters
- region: `~regions.Region` or `~astropy.coordinates.SkyCoord`
Region.
- funcnumpy.func
Function to reduce the data. Default is np.nansum. For boolean Map, use np.any or np.all.
- weights
WcsNDMap
Array to be used as weights. The geometry must be equivalent.
- method{“nearest”, “linear”}
How to interpolate if a position is given.
- Returns
- spectrum
RegionNDMap
Spectrum in the given region.
- spectrum
- to_swapped()[source]#
Return a new map with the opposite scheme (ring or nested).
- Returns
- map
HpxMap
Map object.
- map
- to_unit(unit)#
Convert map to different unit
- Parameters
- unit
Unit
or str New unit
- unit
- Returns
- map
Map
Map with new unit and converted data
- map
- to_wcs(sum_bands=False, normalize=True, proj='AIT', oversample=2, width_pix=None, hpx2wcs=None, fill_nan=True)[source]#
Make a WCS object and convert HEALPIX data into WCS projection.
- Parameters
- sum_bandsbool
Sum over non-spatial axes before reprojecting. If False then the WCS map will have the same dimensionality as the HEALPix one.
- normalizebool
Preserve integral by splitting HEALPIX values between bins?
- projstr
WCS-projection
- oversamplefloat
Oversampling factor for WCS map. This will be the approximate ratio of the width of a HPX pixel to a WCS pixel. If this parameter is None then the width will be set from
width_pix
.- width_pixint
Width of the WCS geometry in pixels. The pixel size will be set to the number of pixels satisfying
oversample
orwidth_pix
whichever is smaller. If this parameter is None then the width will be set fromoversample
.- hpx2wcs
HpxToWcsMapping
Set the HPX to WCS mapping object that will be used to generate the WCS map. If none then a new mapping will be generated based on
proj
andoversample
arguments.
- Returns
- map_out
WcsMap
WCS map object.
- map_out
- to_wcs_tiles(nside_tiles=4, margin='0 deg', method='nearest', oversampling_factor=1)[source]#
Convert HpxNDMap to a list of WCS tiles
- Parameters
- nside_tilesint
Nside for super pixel tiles. Usually nsi
- marginAngle
Width margin of the wcs tile
- method{‘nearest’, ‘linear’}
Interpolation method
- oversampling_factorint
Oversampling factor.
- Returns
- wcs_tileslist of
WcsNDMap
WCS tiles.
- wcs_tileslist of
- upsample(factor, order=0, preserve_counts=True, axis_name=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).
- axis_namestr
Which axis to upsample. By default spatial axes are upsampled.
- Returns
- map
Map
Upsampled map.
- map
- write(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.
- formatstr, optional
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). The following formats are supported:
“gadf” (default)
“fgst-ccube”
“fgst-ltcube”
“fgst-bexpcube”
“fgst-srcmap”
“fgst-template”
“fgst-srcmap-sparse”
“galprop”
“galprop2”
- sparsebool
Sparsify the map by dropping pixels with zero amplitude. This option is only compatible with the ‘gadf’ format.