EDispKernel#

class gammapy.irf.EDispKernel(axes, data=0, unit='', is_pointlike=False, fov_alignment=FoVAlignment.RADEC, meta=None, interp_kwargs=None)[source]#

Bases: IRF

Energy dispersion matrix.

Data format specification: RMF.

Parameters:
energy_axis_trueMapAxis

True energy axis. Its name must be “energy_true”.

energy_axisMapAxis

Reconstructed energy axis. Its name must be “energy”.

dataarray_like

2D energy dispersion matrix.

Examples

Create a Gaussian energy dispersion matrix:

>>> from gammapy.maps import MapAxis
>>> from gammapy.irf import EDispKernel
>>> energy = MapAxis.from_energy_bounds(0.1, 10, 10, unit='TeV')
>>> energy_true = MapAxis.from_energy_bounds(0.1, 10, 10, unit='TeV', name='energy_true')
>>> edisp = EDispKernel.from_gauss(energy_axis_true=energy_true, energy_axis=energy, sigma=0.1, bias=0)

Have a quick look:

>>> print(edisp)
EDispKernel
-----------

  axes  : ['energy_true', 'energy']
  shape : (10, 10)
  ndim  : 2
  unit  :
  dtype : float64

>>> edisp.peek()

Attributes Summary

axes

MapAxes.

data

default_interp_kwargs

Default Interpolation kwargs for IRF.

fov_alignment

Alignment of the field of view coordinate axes, see FoVAlignment.

has_offset_axis

Whether the IRF explicitly depends on offset.

is_pointlike

Whether the IRF is pointlike of full containment.

pdf_matrix

Energy dispersion PDF matrix as a ndarray.

quantity

Quantity as a Quantity object.

required_axes

tag

unit

Map unit as a Unit object.

Methods Summary

cumsum(axis_name)

Compute cumsum along a given axis.

evaluate([method])

Evaluate IRF.

from_diagonal_response(energy_axis_true[, ...])

Create energy dispersion from a diagonal response, i.e. perfect energy resolution.

from_gauss(energy_axis_true, energy_axis, ...)

Create Gaussian energy dispersion matrix (EnergyDispersion).

from_hdulist(hdulist[, hdu1, hdu2])

Create EnergyDispersion object from HDUList.

from_table(table[, format])

Read from Table.

get_bias(energy_true)

Get reconstruction bias for a given true energy.

get_bias_energy(bias[, energy_min, energy_max])

Find energy corresponding to a given bias.

get_mean(energy_true)

Get mean reconstructed energy for a given true energy.

get_resolution(energy_true)

Get energy resolution for a given true energy.

integral(axis_name, **kwargs)

Compute integral along a given axis.

integrate_log_log(axis_name, **kwargs)

Integrate along a given axis.

interp_missing_data(axis_name)

Interpolate missing data along a given axis.

is_allclose(other[, rtol_axes, atol_axes])

Compare two data IRFs for equivalency.

normalize(axis_name)

Normalise data in place along a given axis.

pad(pad_width, axis_name, **kwargs)

Pad IRF along a given axis.

pdf_in_safe_range(lo_threshold, hi_threshold)

PDF matrix with bins outside threshold set to 0.

peek([figsize])

Quick-look summary plots.

plot_bias([ax])

Plot reconstruction bias.

plot_matrix([ax, add_cbar, axes_loc, ...])

Plot PDF matrix.

read(filename[, hdu1, hdu2, checksum])

Read from file.

slice_by_idx(slices)

Slice sub IRF from IRF object.

to_hdulist([format])

Convert RMF to FITS HDU list format.

to_image([lo_threshold, hi_threshold])

Return a 2D edisp by summing the pdf matrix over the ereco axis.

to_table([format])

Convert to Table.

to_table_hdu([format])

Convert to BinTableHDU.

to_unit(unit)

Convert IRF to different unit.

write(filename[, format, checksum])

Write to file.

Attributes Documentation

axes#

MapAxes.

data#
default_interp_kwargs = {'bounds_error': False, 'fill_value': 0, 'method': 'nearest'}#

Default Interpolation kwargs for IRF. Fill zeros and do not interpolate

fov_alignment#

Alignment of the field of view coordinate axes, see FoVAlignment.

has_offset_axis#

Whether the IRF explicitly depends on offset.

is_pointlike#

Whether the IRF is pointlike of full containment.

pdf_matrix#

Energy dispersion PDF matrix as a ndarray.

Rows (first index): True Energy Columns (second index): Reco Energy

quantity#

Quantity as a Quantity object.

required_axes = ['energy_true', 'energy']#
tag = 'edisp_kernel'#
unit#

Map unit as a Unit object.

Methods Documentation

cumsum(axis_name)#

Compute cumsum along a given axis.

Parameters:
axis_namestr

Along which axis to integrate.

Returns:
irfIRF

Cumsum IRF.

evaluate(method=None, **kwargs)#

Evaluate IRF.

Parameters:
**kwargsdict

Coordinates at which to evaluate the IRF.

methodstr {‘linear’, ‘nearest’}, optional

Interpolation method.

Returns:
arrayQuantity

Interpolated values.

classmethod from_diagonal_response(energy_axis_true, energy_axis=None)[source]#

Create energy dispersion from a diagonal response, i.e. perfect energy resolution.

This creates the matrix corresponding to a perfect energy response. It contains ones where the energy_true center is inside the e_reco bin. It is a square diagonal matrix if energy_true = e_reco.

This is useful in cases where code always applies an edisp, but you don’t want it to do anything.

Parameters:
energy_axis_trueMapAxis

True energy axis.

energy_axisMapAxis, optional

Reconstructed energy axis. Default is None.

Examples

If energy_true equals energy, you get a diagonal matrix:

>>> from gammapy.irf import EDispKernel
>>> from gammapy.maps import MapAxis
>>> import astropy.units as u
>>> energy_true_axis = MapAxis.from_energy_edges(
...            [0.5, 1, 2, 4, 6] * u.TeV, name="energy_true"
...        )
>>> edisp = EDispKernel.from_diagonal_response(energy_true_axis)
>>> edisp.plot_matrix() 

Example with different energy binnings:

>>> energy_true_axis = MapAxis.from_energy_edges(
...     [0.5, 1, 2, 4, 6] * u.TeV, name="energy_true"
... )
>>> energy_axis = MapAxis.from_energy_edges([2, 4, 6] * u.TeV)
>>> edisp = EDispKernel.from_diagonal_response(energy_true_axis, energy_axis)
>>> edisp.plot_matrix() 
classmethod from_gauss(energy_axis_true, energy_axis, sigma, bias, pdf_threshold=1e-06)[source]#

Create Gaussian energy dispersion matrix (EnergyDispersion).

Calls gammapy.irf.EnergyDispersion2D.from_gauss().

Parameters:
energy_axis_trueQuantity

Bin edges of true energy axis.

energy_axisQuantity

Bin edges of reconstructed energy axis.

biasfloat or ndarray

Center of Gaussian energy dispersion, bias.

sigmafloat or ndarray

RMS width of Gaussian energy dispersion, resolution.

pdf_thresholdfloat, optional

Zero suppression threshold. Default is 1e-6.

Returns:
edispEDispKernel

Energy dispersion kernel.

classmethod from_hdulist(hdulist, hdu1='MATRIX', hdu2='EBOUNDS')[source]#

Create EnergyDispersion object from HDUList.

Parameters:
hdulistHDUList

HDU list with MATRIX and EBOUNDS extensions.

hdu1str, optional

HDU containing the energy dispersion matrix. Default is “MATRIX”.

hdu2str, optional

HDU containing the energy axis information. Default is “EBOUNDS”.

classmethod from_table(table, format='gadf-dl3')#

Read from Table.

Parameters:
tableTable

Table with IRF data.

format{“gadf-dl3”}, optional

Format specification. Default is “gadf-dl3”.

Returns:
irfIRF

IRF class.

get_bias(energy_true)[source]#

Get reconstruction bias for a given true energy.

Bias is defined as

\[\frac{E_{reco}-E_{true}}{E_{true}}\]
Parameters:
energy_trueQuantity

True energy.

get_bias_energy(bias, energy_min=None, energy_max=None)[source]#

Find energy corresponding to a given bias.

In case the solution is not unique, provide the energy_min or energy_max arguments to limit the solution to the given range. By default, the peak energy of the bias is chosen as energy_min.

Parameters:
biasfloat

Bias value.

energy_minQuantity

Lower bracket value in case solution is not unique.

energy_maxQuantity

Upper bracket value in case solution is not unique.

Returns:
bias_energyQuantity

Reconstructed energy corresponding to the given bias.

get_mean(energy_true)[source]#

Get mean reconstructed energy for a given true energy.

get_resolution(energy_true)[source]#

Get energy resolution for a given true energy.

The resolution is given as a percentage of the true energy.

Parameters:
energy_trueQuantity

True energy.

integral(axis_name, **kwargs)#

Compute integral along a given axis.

This method uses interpolation of the cumulative sum.

Parameters:
axis_namestr

Along which axis to integrate.

**kwargsdict

Coordinates at which to evaluate the IRF.

Returns:
arrayQuantity

Returns 2D array with axes offset.

integrate_log_log(axis_name, **kwargs)#

Integrate along a given axis.

This method uses log-log trapezoidal integration.

Parameters:
axis_namestr

Along which axis to integrate.

**kwargsdict

Coordinates at which to evaluate the IRF.

Returns:
arrayQuantity

Returns 2D array with axes offset.

interp_missing_data(axis_name)#

Interpolate missing data along a given axis.

is_allclose(other, rtol_axes=0.001, atol_axes=1e-06, **kwargs)#

Compare two data IRFs for equivalency.

Parameters:
otherIRF

The IRF to compare against.

rtol_axesfloat, optional

Relative tolerance for the axis comparison. Default is 1e-3.

atol_axesfloat, optional

Absolute tolerance for the axis comparison. Default is 1e-6.

**kwargsdict

Keywords passed to numpy.allclose.

Returns:
is_allclosebool

Whether the IRF is all close.

normalize(axis_name)#

Normalise data in place along a given axis.

Parameters:
axis_namestr

Along which axis to normalize.

pad(pad_width, axis_name, **kwargs)#

Pad IRF along a given axis.

Parameters:
pad_width{sequence, array_like, int}

Number of pixels padded to the edges of each axis.

axis_namestr

Axis to downsample. By default, spatial axes are padded.

**kwargsdict

Keyword argument forwarded to pad.

Returns:
irfIRF

Padded IRF.

pdf_in_safe_range(lo_threshold, hi_threshold)[source]#

PDF matrix with bins outside threshold set to 0.

Parameters:
lo_thresholdQuantity

Low reconstructed energy threshold.

hi_thresholdQuantity

High reconstructed energy threshold.

peek(figsize=(15, 5))[source]#

Quick-look summary plots.

Parameters:
figsizetuple, optional

Size of the figure. Default is (15, 5).

plot_bias(ax=None, **kwargs)[source]#

Plot reconstruction bias.

See get_bias method.

Parameters:
axAxes, optional

Matplotlib axes. Default is None.

**kwargsdict

Keyword arguments.

Returns:
axAxes, optional

Matplotlib axes.

plot_matrix(ax=None, add_cbar=False, axes_loc=None, kwargs_colorbar=None, **kwargs)[source]#

Plot PDF matrix.

Parameters:
axAxes, optional

Matplotlib axes. Default is None.

add_cbarbool, optional

Add a colorbar to the plot. Default is False.

axes_locdict, optional

Keyword arguments passed to append_axes.

kwargs_colorbardict, optional

Keyword arguments passed to colorbar.

kwargsdict

Keyword arguments passed to pcolormesh.

Returns:
axAxes

Matplotlib axes.

classmethod read(filename, hdu1='MATRIX', hdu2='EBOUNDS', checksum=False)[source]#

Read from file.

Parameters:
filenamepathlib.Path or str

File to read.

hdu1str, optional

HDU containing the energy dispersion matrix. Default is “MATRIX”.

hdu2str, optional

HDU containing the energy axis information. Default is “EBOUNDS”.

checksumbool

If True checks both DATASUM and CHECKSUM cards in the file headers. Default is False.

slice_by_idx(slices)#

Slice sub IRF from IRF object.

Parameters:
slicesdict

Dictionary of axes names and slice object pairs. Contains one element for each non-spatial dimension. Axes not specified in the dictionary are kept unchanged.

Returns:
slicedIRF

Sliced IRF object.

to_hdulist(format='ogip', **kwargs)[source]#

Convert RMF to FITS HDU list format.

Parameters:
format{“ogip”, “ogip-sherpa”}

Format to use. Default is “ogip”.

Returns:
hdulistHDUList

RMF in HDU list format.

Notes

For more information on the RMF FITS file format see: https://heasarc.gsfc.nasa.gov/docs/heasarc/caldb/docs/summary/cal_gen_92_002_summary.html

to_image(lo_threshold=None, hi_threshold=None)[source]#

Return a 2D edisp by summing the pdf matrix over the ereco axis.

Parameters:
lo_thresholdQuantity, optional

Low reconstructed energy threshold. Default is None.

hi_thresholdQuantity, optional

High reconstructed energy threshold. Default is None.

to_table(format='ogip')[source]#

Convert to Table.

The output table is in the OGIP RMF format. https://heasarc.gsfc.nasa.gov/docs/heasarc/caldb/docs/memos/cal_gen_92_002/cal_gen_92_002.html#Tab:1 # noqa: E501

Parameters:
format{“ogip”, “ogip-sherpa”}

Format to use. Default is “ogip”.

Returns:
tableTable

Matrix table.

to_table_hdu(format='gadf-dl3')#

Convert to BinTableHDU.

Parameters:
format{“gadf-dl3”}, optional

Format specification. Default is “gadf-dl3”.

Returns:
hduBinTableHDU

IRF data table HDU.

to_unit(unit)#

Convert IRF to different unit.

Parameters:
unitUnit or str

New unit.

Returns:
irfIRF

IRF with new unit and converted data.

write(filename, format='ogip', checksum=False, **kwargs)[source]#

Write to file.

Parameters:
filenamestr

Filename.

format{“ogip”, “ogip-sherpa”}

Format to use. Default is “ogip”.

checksumbool

If True checks both DATASUM and CHECKSUM cards in the file headers. Default is False.