EDispKernel#

class gammapy.irf.EDispKernel[source]#

Bases: IRF

Energy dispersion matrix.

Data format specification: RMF.

Parameters:
axeslist of MapAxis or MapAxes
Required axes (in the given order) are:
  • energy_true (true energy axis)

  • energy (reconstructed energy axis)

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

default_interp_kwargs

Default Interpolation kwargs for IRF.

pdf_matrix

Energy dispersion PDF matrix as a ndarray.

required_axes

tag

Methods Summary

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.

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.

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, format])

Read from file.

to_hdulist([format, creation])

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.

write(filename[, format, checksum, creation])

Write to file.

Attributes Documentation

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

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

pdf_matrix#

Energy dispersion PDF matrix as a ndarray.

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

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

Methods Documentation

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

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.

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.

This method creates a figure with two subplots:

  • Bias plot : reconstruction bias as a function of true energy

  • Energy dispersion matrix plot : probability density function matrix to have energy as a function of energy_true

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, format='gadf')[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.

format{“gadf”, “gtdrm”}

FITS format convention. Default is “gadf”.

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

Convert RMF to FITS HDU list format.

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

Format to use. Default is “ogip”.

creationCreatorMetadata, optional

Creation metadata to add to the file. If None, default metadata is added. Default is None.

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.

write(filename, format='ogip', checksum=False, creation=None, **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.

creationCreatorMetadata, optional

Creation metadata to add to the file. If None, default metadata is added. Default is None.

__init__(axes, data=0, unit='', is_pointlike=False, fov_alignment=FoVAlignment.RADEC, meta=None, interp_kwargs=None)#
classmethod __new__(*args, **kwargs)#