HpxMap

class gammapy.maps.HpxMap(geom, data, meta=None, unit='')[source]

Bases: gammapy.maps.Map

Base class for HEALPIX map classes.

Parameters
geomHpxGeom

HEALPix geometry object.

datandarray

Data array.

metadict

Dictionary to store meta data.

unitUnit

The map unit

Attributes Summary

data

Data array (ndarray)

geom

Map geometry (Geom)

is_mask

Whether map is mask with bool dtype

meta

Map meta (dict)

quantity

Map data times unit (Quantity)

tag

unit

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

downsample(factor[, preserve_counts, axis_name])

Downsample the spatial dimension by a given factor.

fill_by_coord(coords[, weights])

Fill pixels at coords with given weights.

fill_by_idx(idx[, weights])

Fill pixels at idx with given weights.

fill_by_pix(pix[, weights])

Fill pixels at pix with given weights.

fill_events(events)

Fill event coordinates (EventList).

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

Generate an empty map from a Geom instance.

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.

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.

iter_by_image()

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

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

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 given vals.

set_by_idx(idx, vals)

Set pixels at idx with given vals.

set_by_pix(pix, vals)

Set pixels at pix with given vals.

slice_by_idx(slices)

Slice sub map from map object.

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_swapped()

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.

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

Upsample the spatial dimension by a given factor.

write(filename[, overwrite])

Write to a FITS file.

Attributes Documentation

data

Data array (ndarray)

geom

Map geometry (Geom)

is_mask

Whether map is mask with bool dtype

meta

Map meta (dict)

quantity

Map data times unit (Quantity)

tag = 'map'
unit

Map unit (Unit)

Methods Documentation

apply_edisp(edisp)

Apply energy dispersion to map. Requires energy axis.

Parameters
edispgammapy.irf.EDispKernel

Energy dispersion matrix

Returns
mapWcsNDMap

Map with energy dispersion applied.

coadd(map_in, weights=None)

Add the contents of map_in to this map.

This method can be used to combine maps containing integral quantities (e.g. counts) or differential quantities if the maps have the same binning.

Parameters
map_inMap

Input map.

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

The weight factors while adding

copy(**kwargs)

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

Parameters
**kwargsdict

Keyword arguments to overwrite in the map constructor.

Returns
copyMap

Copied Map.

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='')[source]

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.

metadict

Dictionary to store meta data.

unitstr or Unit

The map unit

Returns
mapHpxMap

A HPX map object.

abstract crop(crop_width)

Crop the spatial dimensions of the map.

Parameters
crop_width{sequence, array_like, int}

Number of pixels cropped from the edges of each axis. Defined analogously to pad_with from numpy.pad.

Returns
mapMap

Cropped map.

cumsum(axis_name)

Compute cumulative sum along a given axis

Parameters
axis_namestr

Along which axis to integrate.

Returns
cumsumMap

Map with cumulative sum

abstract downsample(factor, preserve_counts=True, axis_name=None)

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
mapMap

Downsampled map.

fill_by_coord(coords, weights=None)

Fill pixels at coords with given weights.

Parameters
coordstuple or MapCoord

Coordinate arrays for each dimension of the map. Tuple should be ordered as (lon, lat, x_0, …, x_n) where x_i are coordinates for non-spatial dimensions of the map.

weightsndarray

Weights vector. Default is weight of one.

abstract fill_by_idx(idx, weights=None)

Fill pixels at idx with given weights.

Parameters
idxtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps.

weightsndarray

Weights vector. Default is weight of one.

fill_by_pix(pix, weights=None)

Fill pixels at pix with given weights.

Parameters
pixtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps. Pixel indices can be either float or integer type. Float indices will be rounded to the nearest integer.

weightsndarray

Weights vector. Default is weight of one.

fill_events(events)

Fill event coordinates (EventList).

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

Generate an empty map from a Geom instance.

Parameters
geomGeom

Map geometry.

datanumpy.ndarray

data array

metadict

Dictionary to store meta data.

unitstr or Unit

Data unit.

Returns
map_outMap

Map object

classmethod from_hdulist(hdu_list, hdu=None, hdu_bands=None, format=None, colname=None)[source]

Make a HpxMap object from a FITS HDUList.

Parameters
hdu_listHDUList

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”

Returns
hpx_mapHpxMap

Map object

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
mapslist of Map objects

List of maps

axisMapAxis

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

axis_namestr

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

Returns
mapMap

Map with additional non-spatial axis.

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

Returns
valsndarray

Values of pixels in the map. np.nan used to flag coords outside of map.

abstract get_by_idx(idx)

Return map values at the given pixel indices.

Parameters
idxtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps.

Returns
valsndarray

Array of pixel values. np.nan used to flag coordinate outside of map

get_by_pix(pix, fill_value=nan)

Return map values at the given pixel coordinates.

Parameters
pixtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps. Pixel indices can be either float or integer type.

fill_valuefloat

Value which is returned if the position is outside of the projection footprint

Returns
valsndarray

Array 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
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_outMap

Map with spatial dimensions only.

Examples

import numpy as np
from gammapy.maps import Map, MapAxis
from astropy.coordinates import SkyCoord
from astropy import units as u

# Define map axes
energy_axis = MapAxis.from_edges(
    np.logspace(-1., 1., 4), unit='TeV', name='energy',
)

time_axis = MapAxis.from_edges(
    np.linspace(0., 10, 20), unit='h', name='time',
)

# Define map center
skydir = SkyCoord(0, 0, frame='galactic', unit='deg')

# Create map
m_wcs = Map.create(
    map_type='wcs',
    binsz=0.02,
    skydir=skydir,
    width=10.0,
    axes=[energy_axis, time_axis],
)

# Get image by coord tuple
image = m_wcs.get_image_by_coord(('500 GeV', '1 h'))

# Get image by coord dict with strings
image = m_wcs.get_image_by_coord({'energy': '500 GeV', 'time': '1 h'})

# Get image by coord dict with quantities
image = m_wcs.get_image_by_coord({'energy': 0.5 * u.TeV, 'time': 1 * u.h})
get_image_by_idx(idx)

Return spatial map at the given axis pixel indices.

Parameters
idxtuple

Tuple of scalar indices for each non spatial dimension of the map. Tuple should be ordered as (I_0, …, I_n).

Returns
map_outMap

Map with spatial dimensions only.

get_image_by_pix(pix)

Return spatial map at the given axis pixel coordinates

Parameters
pixtuple

Tuple of scalar pixel coordinates for each non-spatial dimension of the map. Tuple should be ordered as (I_0, …, I_n). Pixel coordinates can be either float or integer type.

Returns
map_outMap

Map with spatial dimensions only.

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 the region, 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.

weightsWcsNDMap

Array to be used as weights. The geometry must be equivalent.

Returns
spectrumRegionNDMap

Spectrum in the given region.

integral(axis_name, coords, **kwargs)

Compute integral along a given axis

This method uses interpolation of the cumulative sum.

Parameters
axis_namestr

Along which axis to integrate.

coordsdict or MapCoord

Map coordinates

**kwargsdict

Coordinates at which to evaluate the IRF

Returns
arrayQuantity

Returns 2D array with axes offset

abstract interp_by_coord(coords, method='nearest', fill_value=None)

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{“nearest”, “linear”}

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.

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
valsndarray

Interpolated pixel values.

abstract interp_by_pix(pix, method='nearest', fill_value=None)

Interpolate map values at the given pixel coordinates.

Parameters
pixtuple

Tuple of pixel coordinate arrays for each dimension of the map. Tuple should be ordered as (p_lon, p_lat, p_0, …, p_n) where p_i are pixel coordinates for non-spatial dimensions of the map.

method{“nearest”, “linear”}

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.

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
valsndarray

Interpolated pixel values.

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

Interpolate map to input geometry.

Parameters
geomGeom

Target Map geometry

preserve_countsbool

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

Returns
interp_mapMap

Interpolated Map

iter_by_image()

Iterate over image planes of the map.

This is a generator yielding (data, idx) tuples, where data is a numpy.ndarray view of the image plane data, and idx is a tuple of int, the index of the image plane.

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

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.

Parameters
positionSkyCoord

Test position

Returns
positionSkyCoord

Nearest position in the mask

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
mapMap

Padded map.

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
axesndarray of Axes

Axes grid

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

Returns
map_outMap

Map object

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

Reduce map over a single non-spatial axis

Parameters
axis_name: str

The name of the axis to reduce over

funcufunc

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

weightsMap

Weights to be applied.

Returns
map_outMap

Map with the given non-spatial axes reduced

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

Reduce map over non-spatial axes

Parameters
funcufunc

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

weightsMap

Weights to be applied.

Returns
map_outMap

Map with non-spatial axes reduced

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

Parameters
axisMapAxis

New map axis.

weightsMap

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

ufuncufunc

ufunc to use to resample the axis. Default is numpy.add.

Returns
mapMap

Map with resampled axis.

set_by_coord(coords, vals)

Set pixels at coords with given vals.

Parameters
coordstuple or MapCoord

Coordinate arrays for each dimension of the map. Tuple should be ordered as (lon, lat, x_0, …, x_n) where x_i are coordinates for non-spatial dimensions of the map.

valsndarray

Values vector.

abstract set_by_idx(idx, vals)

Set pixels at idx with given vals.

Parameters
idxtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps.

valsndarray

Values vector.

set_by_pix(pix, vals)

Set pixels at pix with given vals.

Parameters
pixtuple

Tuple of pixel index arrays for each dimension of the map. Tuple should be ordered as (I_lon, I_lat, I_0, …, I_n) for WCS maps and (I_hpx, I_0, …, I_n) for HEALPix maps. Pixel indices can be either float or integer type. Float indices will be rounded to the nearest integer.

valsndarray

Values vector.

slice_by_idx(slices)

Slice sub map from map object.

Parameters
slicesdict

Dict of axes names and integers or slice object pairs. Contains one element for each non-spatial dimension. For integer indexing the corresponding axes is dropped from the map. Axes not specified in the dict are kept unchanged.

Returns
map_outMap

Sliced map object.

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

weightsMap

Weights to be applied. The Map should have the same geometry.

Returns
map_outMap

Map with non-spatial axes summed over

to_cube(axes)

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

This will result in a Map with a new geometry with N+M dimensions where N is the number of current dimensions and M is the number of axes in the list. The data is reshaped onto the new geometry

Parameters
axeslist

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

Returns
mapWcsNDMap

new map

to_hdu(hdu=None, hdu_bands=None, sparse=False, format=None)[source]

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_outBinTableHDU or ImageHDU

Output HDU containing map data.

to_hdulist(hdu='SKYMAP', hdu_bands=None, sparse=False, format='gadf')[source]

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_listHDUList
abstract to_swapped()[source]

Return a new map with the opposite scheme (ring or nested).

Returns
mapHpxMap

Map object.

to_unit(unit)

Convert map to different unit

Parameters
unitUnit or str

New unit

Returns
mapMap

Map with new unit and converted data

abstract to_wcs(sum_bands=False, normalize=True, proj='AIT', oversample=2, width_pix=None, hpx2wcs=None)[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 or width_pix whichever is smaller. If this parameter is None then the width will be set from oversample.

hpx2wcsHpxToWcsMapping

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 and oversample arguments.

Returns
map_outWcsMap

WCS map object.

abstract upsample(factor, order=0, preserve_counts=True, axis=None)

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
mapMap

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