FOVCube

class gammapy.background.FOVCube(coordx_edges=None, coordy_edges=None, energy_edges=None, data=None, scheme=None)[source]

Bases: object

Field of view cube.

Container class for cubes (X, Y, energy).

The class has methods for reading a cube from a FITS file, write a cube to a FITS file and plot the cubes among others.

The order of the axes in the cube is (E, y, x), so in order to access the data correctly, the call is cube.data[energy_bin, coordy_bin, coordx_bin].

This class is very generic and can be used to contain cubes of different kinds of data. However, for the FITS reading/writing methods, special parameter names have to be defined, following the corresponding specifications.

This is taken care of by the define_scheme method. The user only has to specify the correct scheme parameter.

Currently accepted schemes are:

  • bg_cube: scheme for background cubes; spatial coordinates (X, Y) are in detector coordinates (a.k.a. nominal system coordinates).
  • bg_counts_cube: scheme for count cubes specific for background cube determination
  • bg_livetime_cube: scheme for livetime cubes specific for background cube determination

If no scheme is specified, a generic one is applied. New ones can be defined in define_scheme. The method also defines useful parameter names for the plots axis/title labels specific to each scheme.

Parameters:

coordx_edges : Angle, optional

Spatial bin edges vector (low and high). X coordinate.

coordy_edges : Angle, optional

Spatial bin edges vector (low and high). Y coordinate.

energy_edges : EnergyBounds, optional

Energy bin edges vector (low and high).

data : Quantity, optional

Data cube matrix in (energy, X, Y) format.

scheme : str, optional

String identifying parameter naming scheme for FITS files and plots.

Examples

Access cube data:

energy_bin = cube.energy_edges.find_energy_bin('2 TeV')
coord_bin = cube.find_coord_bin(coord=Angle([0., 0.], 'deg'))
cube.data[energy_bin, coord_bin[1], coord_bin[0]]

Attributes Summary

bin_volume Per-pixel bin volume.
coord_wcs WCS object describing the coordinates of the coord (X, Y) bins (WCS).
energy_edges Energy binning (EnergyBounds).
image_bin_centers Image bin centers (x, y) (2x Angle).
image_extent Image extent (Angle).
integral Integral of the cube (Quantity).
integral_images Integral of the cube images (Quantity).
scheme
scheme_dict Naming scheme, depending on the kind of cube (dict).
spectrum_extent Spectrum extent (Quantity).

Methods Summary

define_scheme([scheme]) Define naming scheme, depending on the kind of cube.
fill_events(event_lists) Fill events histogram.
find_coord_bin(coord) Find the bins that contain the specified coord (X, Y) pairs.
find_coord_bin_edges(coord) Find the bin edges of the specified coord (X, Y) pairs.
from_fits_image(image_hdu, energy_hdu[, scheme]) Read cube from a FITS image.
from_fits_table(hdu[, scheme]) Read cube from a FITS binary table.
make_spectrum(coord[, ebounds]) Generate energy spectrum at a certain position in the FOV.
plot_image(energy[, ax, style_kwargs]) Plot image for the energy bin containing the specified energy.
plot_spectrum(coord[, ebounds, ax, style_kwargs]) Plot spectra for the coord bin containing the specified coord (X, Y) pair.
read(filename[, format, scheme, hdu]) Read cube from FITS file.
to_fits_image() Convert cube to image FITS format.
to_fits_table() Convert cube to binary table FITS format.
to_table() Convert cube to astropy table format.
write(outfile[, format]) Write cube to fits file.

Attributes Documentation

bin_volume

Per-pixel bin volume.

TODO: explain with formula and units

coord_wcs

WCS object describing the coordinates of the coord (X, Y) bins (WCS).

This method gives the correct answer only for linear X, Y binning.

energy_edges

Energy binning (EnergyBounds).

image_bin_centers

Image bin centers (x, y) (2x Angle).

Returning two separate elements for the X and Y bin centers.

image_extent

Image extent (Angle).

The output array format is (x_lo, x_hi, y_lo, y_hi).

integral

Integral of the cube (Quantity).

The returned quantity has dimension of the data in the cube times solid angle times energy.

integral_images

Integral of the cube images (Quantity).

Calculate the integral of each energy bin (slice) in the cube. Returns an array of integrals.

The returned quantities have dimensions of the data in the cube times solid angle.

scheme = ''
scheme_dict

Naming scheme, depending on the kind of cube (dict).

spectrum_extent

Spectrum extent (Quantity).

The output array format is (e_lo, e_hi).

Methods Documentation

static define_scheme(scheme=None)[source]

Define naming scheme, depending on the kind of cube.

Parameters:

scheme : str, optional

String identifying parameter naming scheme for FITS files and plots.

Returns:

scheme_dict : dict

Dictionary containing parameter naming scheme for FITS files and plots.

fill_events(event_lists)[source]

Fill events histogram.

This add the counts to the existing value array.

Parameters:

event_lists : list of EventList

Python list of event list objects.

find_coord_bin(coord)[source]

Find the bins that contain the specified coord (X, Y) pairs.

Parameters:

coord : Angle

Array of coord (X, Y) pairs to search for.

Returns:

bin_index : ndarray

Array of integers with the indices (x, y) of the coord bin containing the specified coord (X, Y) pair.

find_coord_bin_edges(coord)[source]

Find the bin edges of the specified coord (X, Y) pairs.

Parameters:

coord : Angle

Array of coord (X, Y) pairs to search for.

Returns:

bin_edges : Angle

Coord bin edges (x_lo, x_hi, y_lo, y_hi).

classmethod from_fits_image(image_hdu, energy_hdu, scheme=None)[source]

Read cube from a FITS image.

Parameters:

image_hdu : PrimaryHDU

FOVCube image HDU.

energy_hdu : BinTableHDU

Energy binning table.

scheme : str, optional

String identifying parameter naming scheme for FITS files and plots.

Returns:

cube : FOVCube

FOVCube object.

classmethod from_fits_table(hdu, scheme=None)[source]

Read cube from a FITS binary table.

Parameters:

hdu : BinTableHDU

HDU binary table for the cube.

scheme : str, optional

String identifying parameter naming scheme for FITS files and plots.

Returns:

cube : FOVCube

FOVCube object.

make_spectrum(coord, ebounds=None)[source]

Generate energy spectrum at a certain position in the FOV.

Parameters:

coord : Quantity

Coord (X,Y) pair of cube bin to plot.

ebounds : EnergyBounds, optional

Energy binning for the spectrum

Returns:

spectrum : Quantity

Energy spectrum

plot_image(energy, ax=None, style_kwargs=None)[source]

Plot image for the energy bin containing the specified energy.

Parameters:

energy : Energy

Energy of cube bin to plot.

ax : Axes, optional

Axes of the figure for the plot.

style_kwargs : dict, optional

Style options for the plot.

Returns:

ax : Axes

Axes of the figure containing the plot.

plot_spectrum(coord, ebounds=None, ax=None, style_kwargs=None)[source]

Plot spectra for the coord bin containing the specified coord (X, Y) pair.

Parameters:

coord : Quantity

Coord (X,Y) pair of cube bin to plot.

ax : Axes, optional

Axes of the figure for the plot.

ebounds : EnergyBounds, optional

Energy binning for the spectrum

style_kwargs : dict, optional

Style options for the plot.

Returns:

ax : Axes

Axes of the figure containing the plot.

classmethod read(filename, format='table', scheme=None, hdu='bkg_3d')[source]

Read cube from FITS file.

Several input formats are accepted, depending on the value of the format parameter:

Parameters:

filename : str

Name of file with the cube.

format : str, optional

Format of the cube to read.

scheme : str, optional

String identifying parameter naming scheme for FITS files and plots.

Returns:

cube : FOVCube

FOVCube object.

to_fits_image()[source]

Convert cube to image FITS format.

Returns:

hdu_list : HDUList

HDU list with:

to_fits_table()[source]

Convert cube to binary table FITS format.

Returns:

tbhdu : BinTableHDU

Table containing the cube.

to_table()[source]

Convert cube to astropy table format.

The name of the table is stored in the table meta information under the keyword ‘name’.

Returns:

table : Table

Table containing the cube.

write(outfile, format='table', **kwargs)[source]

Write cube to fits file.

Several output formats are accepted, depending on the value of the format parameter:

Depending on the value of the format parameter, this method calls either writeto or writeto, forwarding the kwargs arguments.

Parameters:

outfile : str

Name of file to write.

format : str, optional

Format of the cube to write.

kwargs

Extra arguments for the corresponding astropy.io.fits writeto method.