SkyImage

class gammapy.image.SkyImage(name=None, data=None, wcs=None, unit='', meta=None)[source]

Bases: gammapy.image.core.MapBase

Sky image.

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

For further information, see Image processing and analysis tools (gammapy.image).

Parameters:

name : str

Name of the image.

data : ndarray

Data array.

wcs : WCS

WCS transformation object.

unit : str

String specifying the data units.

meta : OrderedDict

Dictionary to store meta data.

Attributes Summary

center Center sky coordinate of the image (SkyCoord).
center_pix Center pixel coordinate of the image (PixCoord).
distance_image Distance to nearest exclusion region.
height Maximum angular height of the image.
is_mask Is this a mask (check values, not dtype).
width Maximum angular width of the image.

Methods Summary

assert_allclose(image1, image2[, check_wcs, …]) Assert all-close for SkyImage.
contains(position) Check if given position on the sky is contained in the image.
convolve(kernel[, use_fft]) Convolve sky image with kernel.
coordinates([mode]) Sky coordinate images.
coordinates_pix([mode]) Pixel sky coordinate images.
copy() Copy image.
crop(crop_width) Crop sky image at the edges with given crop width.
cutout(position, size[, copy]) Cut out rectangular piece of a image.
downsample(factor[, method]) Down sample image by a given factor.
empty([name, nxpix, nypix, binsz, xref, …]) Create an empty image from scratch.
empty_like(image[, name, unit, fill, meta]) Create an empty image like the given image.
fill_events(events[, weights]) Fill events (modifies data attribute).
footprint([mode]) Footprint of the image on the sky.
from_image_hdu(image_hdu) Create image from ImageHDU.
from_wcs_nd_map(wcs_map_nd) Create from a gammapy.maps.WcsNDMap.
info() Print summary info about the image.
lookup(position[, interpolation]) Lookup value at given sky position.
lookup_max([region]) Find position of maximum in a image.
lookup_pix(position[, interpolation]) Lookup value at given pixel position.
pad(pad_width[, mode]) Pad sky image at the edges.
paste(image[, method, wcs_check]) Paste smaller image into image.
plot([ax, fig, add_cbar, stretch]) Plot image on matplotlib WCS axes.
read(filename[, hdu]) Read image from FITS file (SkyImage).
region_mask(region) Create a boolean mask for a region.
reproject(reference[, mode]) Reproject image to given reference.
show([viewer, ds9options]) Show image in image viewer.
smooth([kernel, radius]) Smooth the image (works on and returns a copy).
solid_angle() Solid angle image (2-dim Quantity in sr).
threshold(threshold) Threshold this image, creating a mask.
to_image_hdu() Convert image to a PrimaryHDU.
to_quantity() Convert image to Quantity.
to_sherpa_data2d([dstype]) Convert image to Data2D or Data2DInt class.
to_wcs_nd_map() Convert to a gammapy.maps.WcsNDMap.
upsample(factor, **kwargs) Up sample image by a given factor.
wcs_pixel_scale([method]) Pixel scale.
wcs_pixel_to_skycoord(xp, yp) Convert a set of pixel coordinates into a SkyCoord coordinate.
wcs_skycoord_to_pixel(coords) Convert a set of SkyCoord coordinates into pixels.
write(filename, *args, **kwargs) Write image to FITS file.

Attributes Documentation

center

Center sky coordinate of the image (SkyCoord).

center_pix

Center pixel coordinate of the image (PixCoord).

distance_image

Distance to nearest exclusion region.

Compute distance image, i.e. the Euclidean (=Cartesian 2D) distance (in pixels) to the nearest exclusion region.

We need to call distance_transform_edt twice because it only computes dist for pixels outside exclusion regions, so to get the distances for pixels inside we call it on the inverted mask and then combine both distance images into one, using negative distances (note the minus sign) for pixels inside exclusion regions.

If data consist only of ones, it’ll be supposed to be far away from zero pixels, so in capacity of answer it should be return the matrix with the shape as like as data but packed by constant value Max_Value (MAX_VALUE = 1e10).

If data consist only of zeros, it’ll be supposed to be deep inside an exclusion region, so in capacity of answer it should be return the matrix with the shape as like as data but packed by constant value -Max_Value (MAX_VALUE = 1e10).

Returns:

distance : SkyImage

Sky image of distance to nearest exclusion region.

Examples

>>> from gammapy.image import SkyImage
>>> data = np.array([[0., 0., 1.], [1., 1., 1.]])
>>> mask = SkyImage(data=data)
>>> print(mask.distance_image.data)
[[-1, -1, 1], [1, 1, 1.41421356]]
height

Maximum angular height of the image.

is_mask

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

width

Maximum angular width of the image.

Methods Documentation

static assert_allclose(image1, image2, check_wcs=True, check_name=True, check_unit=True)[source]

Assert all-close for SkyImage.

A useful helper function to implement tests.

contains(position)[source]

Check if given position on the sky is contained in the image.

Parameters:

position : SkyCoord

Position on the sky.

Returns:

containment : array

Bool array

convolve(kernel, use_fft=False, **kwargs)[source]

Convolve sky image with kernel.

Parameters:

kernel : ndarray

2D array representing the convolution kernel.

**kwargs : dict

Further keyword arguments passed to convolve.

coordinates(mode='center')[source]

Sky coordinate images.

Parameters:

mode : {‘center’, ‘edges’}

Return coordinate values at the pixels edges or pixel centers.

Returns:

coordinates : SkyCoord

Position on the sky.

coordinates_pix(mode='center')[source]

Pixel sky coordinate images.

Parameters:

mode : {‘center’, ‘edges’}

Return coordinate values at the pixels edges or pixel centers.

Returns:

x, y : tuple

Return arrays representing the coordinates of a sky grid.

copy()[source]

Copy image.

crop(crop_width)[source]

Crop sky image at the edges with given crop width.

Analogous method to SkyImage.pad() to crop the sky image at the edges. Adapts the WCS specification accordingly.

Parameters:

crop_width : {sequence, array_like, int}

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

Returns:

image : SkyImage

Cropped image

cutout(position, size, copy=True)[source]

Cut out rectangular piece of a image.

See Cutout and paste for more information how to cut and paste sky images.

Parameters:

position : SkyCoord

Position of the center of the image to cut out.

size : int, array-like, Quantity

The size of the cutout array along each axis. If size is a scalar number or a scalar Quantity, then a square cutout of size will be created. If size has two elements, they should be in (ny, nx) order. Scalar numbers in size are assumed to be in units of pixels. size can also be a Quantity object or contain Quantity objects. Such Quantity objects must be in pixel or angular units. For all cases, size will be converted to an integer number of pixels, rounding the the nearest integer. See the mode keyword for additional details on the final cutout size.

Note

If size is in angular units, the cutout size is converted to pixels using the pixel scales along each axis of the image at the CRPIX location. Projection and other non-linear distortions are not taken into account.

Returns:

cutout : SkyImage

Cut out image.

downsample(factor, method=<function nansum>)[source]

Down sample image by a given factor.

The image is down sampled using skimage.measure.block_reduce. If the shape of the data is not divisible by the down sampling factor, the image must be padded beforehand to the correct shape.

Parameters:

factor : int

Down sampling factor.

method : np.ufunc (np.nansum), optional

Method how to combine the image blocks.

Returns:

image : SkyImage

Down sampled image.

classmethod empty(name=None, nxpix=200, nypix=200, binsz=0.02, xref=0, yref=0, fill=0, proj='CAR', coordsys='GAL', xrefpix=None, yrefpix=None, dtype='float64', unit='', meta=None)[source]

Create an empty image from scratch.

Uses the same parameter names as the Fermi tool gtbin (see https://fermi.gsfc.nasa.gov/ssc/data/analysis/scitools/help/gtbin.txt).

If no reference pixel position is given it is assumed to be at the center of the image.

Parameters:

name : str

Name of the image.

nxpix : int, optional

Number of pixels in x axis. Default is 200.

nypix : int, optional

Number of pixels in y axis. Default is 200.

binsz : float, optional

Bin size for x and y axes in units of degrees. Default is 0.02.

xref : float, optional

Coordinate system value at reference pixel for x axis. Default is 0.

yref : float, optional

Coordinate system value at reference pixel for y axis. Default is 0.

fill : float, optional

Fill image with constant value. Default is 0.

proj : string, optional

Any valid WCS projection type. Default is ‘CAR’ (cartesian).

coordsys : {‘CEL’, ‘GAL’}, optional

Coordinate system, either Galactic (‘GAL’) or Equatorial (‘CEL’). Default is ‘GAL’ (Galactic).

xrefpix : float, optional

Coordinate system reference pixel for x axis. Default is None.

yrefpix: float, optional

Coordinate system reference pixel for y axis. Default is None.

dtype : str, optional

Data type, default is float32

unit : str or Unit

Data unit.

meta : OrderedDict

Meta data attached to the image.

Returns:

image : SkyImage

Empty image.

classmethod empty_like(image, name=None, unit='', fill=0, meta=None)[source]

Create an empty image like the given image.

The WCS is copied over, the data array is filled with the fill value.

Parameters:

image : SkyImage or ImageHDU

Instance of SkyImage.

fill : float, optional

Fill image with constant value. Default is 0.

name : str

Name of the image.

unit : str

String specifying the data units.

meta : OrderedDict

Dictionary to store meta data.

fill_events(events, weights=None)[source]

Fill events (modifies data attribute).

Calls numpy.histogramdd

Parameters:

events : EventList

Event list

weights : str, optional

Column to use as weights (none by default)

Examples

Show example how to make an empty image and fill it.

footprint(mode='edges')[source]

Footprint of the image on the sky.

Parameters:

mode : {‘center’, ‘edges’}

Use corner pixel centers or corners?

Returns:

coordinates : OrderedDict

Dictionary of the positions of the corners of the image with keys {‘lower left’, ‘upper left’, ‘upper right’, ‘lower right’} and SkyCoord objects as values.

Examples

>>> from gammapy.image import SkyImage
>>> image = SkyImage.empty(nxpix=3, nypix=2)
>>> coord = image.footprint(mode='corner')
>>> coord['lower left']
<SkyCoord (Galactic): (l, b) in deg
    (0.03, -0.02)>
classmethod from_image_hdu(image_hdu)[source]

Create image from ImageHDU.

Parameters:

image_hdu : astropy.io.fits.ImageHDU

Source image HDU.

Examples

>>> from astropy.io import fits
>>> from gammapy.image import SkyImage
>>> hdu_list = fits.open('data.fits')
>>> image = SkyImage.from_image_hdu(hdu_list['myimage'])
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 SkyImage over to the new maps classes.

info()[source]

Print summary info about the image.

lookup(position, interpolation=None)[source]

Lookup value at given sky position.

Parameters:

position : SkyCoord

Position on the sky.

interpolation : {‘None’}

Interpolation mode.

lookup_max(region=None)[source]

Find position of maximum in a image.

Parameters:

region : SkyRegion (optional)

Limit lookup of maximum to that given sky region.

Returns:

(position, value): SkyCoord, float

Position and value of the maximum.

lookup_pix(position, interpolation=None)[source]

Lookup value at given pixel position.

Parameters:

position : PixCoord

Pixel coordinate position

interpolation : {‘None’}

Interpolation mode.

pad(pad_width, mode='reflect', **kwargs)[source]

Pad sky image at the edges.

Calls numpy.pad, passing mode and kwargs to it and adapts the wcs specifcation.

Parameters:

pad_width : {sequence, array_like, int}

Number of values padded to the edges of each axis, passed to numpy.pad

mode : str (‘reflect’)

Padding mode, passed to numpy.pad.

Returns:

image : SkyImage

Padded image

Examples

>>> from gammapy.image import SkyImage
>>> image = SkyImage.empty(nxpix=10, nypix=13)
>>> print(image.data.shape)
(13, 10)
>>> image2 = image.pad(pad_width=4, mode='reflect')
>>> image2.data.shape
(18, 21)
paste(image, method='sum', wcs_check=True)[source]

Paste smaller image into image.

WCS specifications of both images must be aligned. If not call SkyImage.reproject() on one of the images first. See Cutout and paste more for information how to cut and paste sky images.

Parameters:

image : SkyImage

Smaller image to paste.

method : {‘sum’, ‘replace’}, optional

Sum or replace total values with cutout values.

wcs_check : bool

Check if both WCS are aligned. Raises WcsError if not. Disable for performance critical computations.

plot(ax=None, fig=None, add_cbar=False, stretch='linear', **kwargs)[source]

Plot image on matplotlib WCS axes.

Parameters:

ax : WCSAxes, optional

WCS axis object to plot on.

fig : Figure, optional

Figure

stretch : str, optional

Scaling for image (‘linear’, ‘sqrt’, ‘log’). Similar to normalize and stretch functions in ds9. See http://docs.astropy.org/en/stable/visualization/normalization.html

Returns

——-

fig : Figure, optional

Figure

ax : WCSAxes, optional

WCS axis object

cbar : ?

Colorbar object (if add_cbar=True was set)

Examples

——–

>>> from astropy.visualization import simple_norm

>>> from gammapy.image import SkyImage

>>> filename = ‘$GAMMAPY_EXTRA/datasets/fermi_2fhl/fermi_2fhl_vela.fits.gz’

>>> image = SkyImage.read(filename, hdu=2)

>>> norm = simple_norm(image, ‘sqrt’)

>>> plt.imshow(image, norm = norm)

>>> plt.show()

>>> #Equivalent to :

>>> image.plot(stretch=’sqrt’)

classmethod read(filename, hdu=None, **kwargs)[source]

Read image from FITS file (SkyImage).

Parameters are passed to SmartHDUList.

region_mask(region)[source]

Create a boolean mask for a region.

The data of this image is unchanged, a new mask is returned.

The mask is:

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

region : PixelRegion or SkyRegion object

A region on the sky could be defined in pixel or sky coordinates.

Returns:

mask : SkyImage

A boolean sky mask.

Examples

>>> from gammapy.image import SkyImage
>>> from regions import CirclePixelRegion, PixCoord
>>> region = CirclePixelRegion(center=PixCoord(x=2, y=1), radius=1.1)
>>> image = SkyImage.empty(nxpix=5, nypix=4)
>>> mask = image.region_mask(region)
>>> print(mask.data.astype(int))
[[0 0 1 0 0]
 [0 1 1 1 0]
 [0 0 1 0 0]
 [0 0 0 0 0]]
reproject(reference, mode='interp', *args, **kwargs)[source]

Reproject image to given reference.

Parameters:

reference : Header, or SkyImage

Reference image 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:

image : SkyImage

Image reprojected onto reference.

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

Show image 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 imshow.

smooth(kernel='gauss', radius=<Quantity 0.1 deg>, **kwargs)[source]

Smooth the image (works on and returns a copy).

The definition of the smoothing parameter radius is equivalent to the one that is used in ds9 (see ds9 smoothing).

Parameters:

kernel : {‘gauss’, ‘disk’, ‘box’}

Kernel shape

radius : Quantity 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 SkyImage.wcs_pixel_scale().

kwargs : dict

Keyword arguments passed to uniform_filter (‘box’), gaussian_filter (‘gauss’) or convolve (‘disk’).

Returns:

image : SkyImage

Smoothed image (a copy, the original object is unchanged).

solid_angle()[source]

Solid angle image (2-dim Quantity in sr).

threshold(threshold)[source]

Threshold this image, creating a mask.

Parameters:

threshold : float

Threshold value.

Returns:

mask : SkyImage

Mask with 0 where data > threshold and 1 otherwise

Examples

TODO: some more docs and example

to_image_hdu()[source]

Convert image to a PrimaryHDU.

Returns:

hdu : PrimaryHDU

Primary image hdu object.

to_quantity()[source]

Convert image to Quantity.

to_sherpa_data2d(dstype='Data2D')[source]

Convert image to Data2D or Data2DInt class.

Parameters:

dstype : {‘Data2D’, ‘Data2DInt’}

Sherpa data type.

to_wcs_nd_map()[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 SkyImage over to the new maps classes.

upsample(factor, **kwargs)[source]

Up sample image by a given factor.

The image is up sampled using scipy.ndimage.zoom.

Parameters:

factor : int

Up sampling factor.

order : int

Order of the interpolation used for upsampling.

Returns:

image : SkyImage

Up sampled image.

wcs_pixel_scale(method='cdelt')[source]

Pixel scale.

Returns angles along each axis of the image pixel at the CRPIX location once it is projected onto the plane of intermediate world coordinates.

Calls proj_plane_pixel_scales.

Parameters:

method : {‘cdelt’, ‘proj_plane’} (default ‘cdelt’)

Result is calculated according to the ‘cdelt’ or ‘proj_plane’ methods.

Returns:

angle : Angle

An angle of projection plane increments corresponding to each pixel side (axis).

Examples

>>> from gammapy.image import SkyImage
>>> image = SkyImage.empty(nxpix=3, nypix=2)
>>> image.wcs_pixel_scale()
<Angle [ 0.02, 0.02] deg>
wcs_pixel_to_skycoord(xp, yp)[source]

Convert a set of pixel coordinates into a SkyCoord coordinate.

Calls pixel_to_skycoord, passing xp, yp to it.

Parameters:

xp, yp : float or ndarray

The coordinates to convert.

Returns:

coordinates : SkyCoord

The celestial coordinates.

Examples

>>> from gammapy.image import SkyImage
>>> image = SkyImage.empty(nxpix=10, nypix=15)
>>> x, y = [5, 3.4], [8, 11.2]
>>> image.wcs_pixel_to_skycoord(xp=x, yp=y)
<SkyCoord (Galactic): (l, b) in deg
    [(359.99, 0.02), (0.022, 0.084)]>
wcs_skycoord_to_pixel(coords)[source]

Convert a set of SkyCoord coordinates into pixels.

Calls skycoord_to_pixel, passing coords to it.

Parameters:

coords : SkyCoord

The coordinates to convert.

Returns:

xp, yp : ndarray

The pixel coordinates.

write(filename, *args, **kwargs)[source]

Write image to FITS file.

Parameters:

filename : str

Name of the FITS file.

*args : list

Arguments passed to writeto.

**kwargs : dict

Keyword arguments passed to writeto.