SkyCube

class gammapy.cube.SkyCube(name=None, data=None, wcs=None, energy_axis=None, meta=None)[source]

Bases: gammapy.image.core.MapBase

Sky cube with dimensions lon, lat and energy.

Note

A new set of map and cube classes is being developed in gammapy.maps and long-term will replace the existing gammapy.image.SkyImage and gammapy.cube.SkyCube classes. Please consider trying out gammapy.maps and changing your scripts to use those new classes. See Data Structures for Images and Cubes (gammapy.maps).

The order of the sky cube axes is defined as following:

  • The data array axis order is (energy, lat, lon).
  • The wcs object is a two dimensional celestial WCS with axis order (lon, lat).

For further information, see Cube Style Analysis (gammapy.cube).

Parameters:

name : str

Name of the cube.

data : ndarray

Data array.

wcs : WCS

WCS transformation object.

energy_axis : LogEnergyAxis

Energy axis object, defining the energy transformation.

meta : OrderedDict

Dictionary to store meta data.

Attributes Summary

bin_size Sky cube element bin size (Quantity)
energy_width Energy bin width vector (Quantity)
is_mask Is this a mask (check values, not dtype).
sky_image_ref Empty reference SkyImage.

Methods Summary

assert_allclose(cube1, cube2)
convolve(kernels, **kwargs) Convolve cube with a given set of kernels.
cutout(position, size) Cut out rectangular piece of a cube.
empty([emin, emax, enumbins, eunit, mode]) Create empty sky cube with log equal energy binning from the scratch.
empty_like(reference[, energies, unit, fill]) Create an empty sky cube with a given WCS and energy specification.
energies([mode]) Energy coordinate vector.
fill_events(events[, weights]) Fill events (modifies data attribute).
from_wcs_nd_map(wcs_map_nd) Create from a gammapy.maps.WcsNDMap.
info() Print summary info about the cube.
lookup(position, energy[, interpolation]) Lookup value in the cube at given sky position and energy.
plot_rgb([ax, fig]) Plot sky cube as RGB image.
read(filename[, format]) Read sky cube from FITS file.
read_hdu(hdu_list) Read sky cube from HDU list.
region_mask(region) Create a boolean cube mask for a region.
reproject(reference[, mode]) Reproject spatial dimensions onto a reference image.
show([viewer, ds9options]) Show sky cube in image viewer.
sky_image(energy[, interpolation]) Slice a 2-dim SkyImage from the cube at a given energy.
sky_image_idx(idx) Slice a 2-dim SkyImage from the cube at a given index.
sky_image_integral(emin, emax[, nbins, …]) Integrate cube along the energy axes using the log-log trapezoidal rule.
sky_image_sum(emin, emax) Sum cube along the energy axis.
spectrum(region) Extract spectrum in a given sky region.
to_fits(format) Write to FITS HDU list.
to_images() Convert to SkyCubeImages.
to_sherpa_data3d([dstype]) Convert sky cube to sherpa Data3D or Data3DInt object.
to_wcs_nd_map([energy_axis_mode]) Convert to a gammapy.maps.WcsNDMap.
wcs_pixel_to_skycoord(x, y, z) Convert pixel to world coordinates.
wcs_skycoord_to_pixel(position, energy) Convert world to pixel coordinates.
write(filename[, format]) Write to FITS file.

Attributes Documentation

bin_size

Sky cube element bin size (Quantity)

This is a convenience method which computes this:

cube.energy_width * cube.sky_image_ref.solid_angle()

Units could be “TeV” (or whatever energy_width returns) times “sr”

energy_width

Energy bin width vector (Quantity)

is_mask

Is this a mask (check values, not dtype).

sky_image_ref

Empty reference SkyImage.

Examples

Can be used to access the spatial information of the cube:

>>> from gammapy.cube import SkyCube
>>> cube = SkyCube.empty()
>>> coords = cube.sky_image_ref.coordinates()
>>> solid_angle = cube.sky_image_ref.solid_angle()

Methods Documentation

static assert_allclose(cube1, cube2)[source]
convolve(kernels, **kwargs)[source]

Convolve cube with a given set of kernels.

Parameters:

kernels : list or ndarray

List of 2D convolution kernels or 3D array. The energy axis must correspond to array axis=0.

Returns:

convolved : SkyCube

Convolved cube.

cutout(position, size)[source]

Cut out rectangular piece of a cube.

See cutout() for details.

classmethod empty(emin=0.5, emax=100, enumbins=10, eunit='TeV', mode='edges', **kwargs)[source]

Create empty sky cube with log equal energy binning from the scratch.

Parameters:

emin : float

Minimum energy.

emax : float

Maximum energy.

enumbins : int

Number of energy bins.

eunit : str

Energy unit.

mode : {‘edges’, ‘center’}

Whether emin and emax denote the bin centers or edges.

kwargs : dict

Keyword arguments passed to empty to create the spatial part of the cube.

Returns:

empty_cube : SkyCube

Empty sky cube object.

Examples

Create an empty sky cube:

from gammapy.cube import SkyCube
cube = SkyCube.empty(nxpix=11, nypix=7, enumbins=3, mode='center',
                emin=1, emax=100, eunit='TeV')
classmethod empty_like(reference, energies=None, unit='', fill=0)[source]

Create an empty sky cube with a given WCS and energy specification.

Parameters:

reference : SkyCube or SkyImage

Reference sky cube or image.

energies : Energy or EnergyBounds (optional)

Reference energies, mandatory when a SkyImage is passed.

unit : str

String specifying the data units.

fill : float

Value to fill the data array with.

Returns:

empty_cube : SkyCube

Empty sky cube object.

Examples

Create an empty sky cube from an image and energy center specification:

from astropy import units as u
from gammapy.image import SkyImage
from gammapy.cube import SkyCube
from gammapy.utils.energy import Energy, EnergyBounds

# define reference image
image = SkyImage.empty(nxpix=11, nypix=7)

# define energy binning centers
energies = Energy.equal_log_spacing(1 * u.TeV, 100 * u.TeV, 3)
cube = SkyCube.empty_like(reference=image, energies=energies)

# define energy binning bounds
energies = EnergyBounds.equal_log_spacing(1 * u.TeV, 100 * u.TeV, 3)
cube = SkyCube.empty_like(reference=image, energies=energies)
energies(mode='center')[source]

Energy coordinate vector.

Parameters:

mode : {‘center’, ‘edges’}

Return coordinate values at the pixels edges or pixel centers.

Returns:

coordinates : Quantity

Energy

fill_events(events, weights=None)[source]

Fill events (modifies data attribute).

Parameters:

events : EventList

Event list

weights : str, optional

Column to use as weights (none by default)

classmethod from_wcs_nd_map(wcs_map_nd)[source]

Create from a gammapy.maps.WcsNDMap.

There is no copy of the data or wcs object, this conversion is cheap.

This is meant to help migrate code using SkyCube over to the new maps classes.

info()[source]

Print summary info about the cube.

lookup(position, energy, interpolation=None)[source]

Lookup value in the cube at given sky position and energy.

Parameters:

position : SkyCoord

Position on the sky.

energy : Quantity

Energy

interpolation : {None, ‘linear’, ‘nearest’}

Interpolate data values between energies.

Returns:

value : Quantity

Value at the given sky position and energy.

plot_rgb(ax=None, fig=None, **kwargs)[source]

Plot sky cube as RGB image.

Parameters:

ax : WCSAxes, optional

WCS axis object to plot on.

fig : Figure, optional

Figure

Returns:

ax : WCSAxes

WCS axis object

classmethod read(filename, format='fermi-counts')[source]

Read sky cube from FITS file.

Parameters:

filename : str

File name

format : {‘fermi-counts’, ‘fermi-background’, ‘fermi-exposure’}

Fits file format.

Returns:

sky_cube : SkyCube

Sky cube

classmethod read_hdu(hdu_list)[source]

Read sky cube from HDU list.

Parameters:

object_hdu : ImageHDU

Image HDU object to be read

energy_table_hdu : TableHDU

Table HDU object giving energies of each slice of the Image HDU object_hdu

Returns:

sky_cube : SkyCube

Sky cube

region_mask(region)[source]

Create a boolean cube mask for a region.

The mask is:

  • True for pixels inside the region
  • False for pixels outside the region
Parameters:

region : PixelRegion or SkyRegion

Region in pixel or sky coordinates.

Returns:

mask : SkyCube

A boolean sky cube mask.

reproject(reference, mode='interp', *args, **kwargs)[source]

Reproject spatial dimensions onto a reference image.

Parameters:

reference : Header, SkyImage or SkyCube

Reference wcs specification to reproject the data on.

mode : {‘interp’, ‘exact’}

Interpolation mode.

*args : list

Arguments passed to reproject_interp or reproject_exact.

**kwargs : dict

Keyword arguments passed to reproject_interp or reproject_exact.

Returns:

reprojected_cube : SkyCube

Cube spatially reprojected to the reference.

show(viewer='mpl', ds9options=None, **kwargs)[source]

Show sky cube in image viewer.

Parameters:

viewer : {‘mpl’, ‘ds9’}

Which image viewer to use. Option ‘ds9’ requires ds9 to be installed.

ds9options : list, optional

List of options passed to ds9. E.g. [‘-cmap’, ‘heat’, ‘-scale’, ‘log’]. Any valid ds9 command line option can be passed. See http://ds9.si.edu/doc/ref/command.html for details.

**kwargs : dict

Keyword arguments passed to matplotlib.pyplot.imshow.

sky_image(energy, interpolation=None)[source]

Slice a 2-dim SkyImage from the cube at a given energy.

Parameters:

energy : Quantity

Energy value

interpolation : {None, ‘linear’, ‘nearest’}

Interpolate data values between energies. None corresponds to ‘nearest’, but might have advantages in performance, because no interpolator is set up.

Returns:

image : SkyImage

2-dim sky image

sky_image_idx(idx)[source]

Slice a 2-dim SkyImage from the cube at a given index.

Parameters:

idx : int

Index of the sky image.

Returns:

image : SkyImage

2-dim sky image

sky_image_integral(emin, emax, nbins=10, per_decade=False, interpolation='linear')[source]

Integrate cube along the energy axes using the log-log trapezoidal rule.

Parameters:

emin : Quantity

Integration range minimum.

emax : Quantity

Integration range maximum.

nbins : int, optional

Number of grid points used for the integration.

per_decade : bool

Whether nbins is per decade.

interpolation : {None, ‘linear’, ‘nearest’}

Interpolate data values between energies.

Returns:

image : SkyImage

Integral image.

sky_image_sum(emin, emax)[source]

Sum cube along the energy axis.

Similar to the sky_image_integral method, but not doing interpolation / integration.

Just selects a subset of energy bins and sums those. This is useful for counts.

spectrum(region)[source]

Extract spectrum in a given sky region.

Parameters:

region : SkyRegion

Sky region to extract the spectrum from.

Returns:

spectrum : Table

Summed spectrum of pixels in the mask.

to_fits(format)[source]

Write to FITS HDU list.

Parameters:

format : {‘fermi-counts’, ‘fermi-background’, ‘fermi-exposure’}

Fits file format.

Returns:

hdu_list : HDUList

to_images()[source]

Convert to SkyCubeImages.

to_sherpa_data3d(dstype='Data3D')[source]

Convert sky cube to sherpa Data3D or Data3DInt object.

Parameters:

dstype : {‘Data3D’, ‘Data3DInt’}

Sherpa data type.

to_wcs_nd_map(energy_axis_mode='center')[source]

Convert to a gammapy.maps.WcsNDMap.

There is no copy of the data or wcs object, this conversion is cheap.

This is meant to help migrate code using SkyCube over to the new maps classes.

wcs_pixel_to_skycoord(x, y, z)[source]

Convert pixel to world coordinates.

Parameters:

x : ndarry

x coordinate array

y : ndarry

y coordinate array

z : ndarry

z coordinate array

Returns:

(position, energy) : tuple

Tuple of (SkyCoord, Quantity).

wcs_skycoord_to_pixel(position, energy)[source]

Convert world to pixel coordinates.

Parameters:

position : SkyCoord

Position on the sky.

energy : Quantity

Energy

Returns:

(x, y, z) : tuple

Tuple of x, y, z coordinates.

write(filename, format='fermi-counts', **kwargs)[source]

Write to FITS file.

Parameters:

filename : str

Filename

format : {‘fermi-counts’, ‘fermi-background’, ‘fermi-exposure’}

Fits file format.