EffectiveAreaTable2D#

class gammapy.irf.EffectiveAreaTable2D[source]#

Bases: IRF

2D effective area table.

Data format specification: AEFF_2D.

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

  • offset (field of view offset axis)

dataQuantity

Effective area.

metadict

Metadata dictionary.

Examples

Here’s an example you can use to learn about this class:

>>> from gammapy.irf import EffectiveAreaTable2D
>>> filename = '$GAMMAPY_DATA/cta-1dc/caldb/data/cta/1dc/bcf/South_z20_50h/irf_file.fits'
>>> aeff = EffectiveAreaTable2D.read(filename, hdu='EFFECTIVE AREA')
>>> print(aeff)
EffectiveAreaTable2D
--------------------

  axes  : ['energy_true', 'offset']
  shape : (42, 6)
  ndim  : 2
  unit  : m2
  dtype : >f4

Here’s another one, created from scratch, without reading a file:

>>> from gammapy.irf import EffectiveAreaTable2D
>>> from gammapy.maps import MapAxis
>>> energy_axis_true = MapAxis.from_energy_bounds(
...        "0.1 TeV", "100 TeV", nbin=30, name="energy_true"
...    )
>>> offset_axis = MapAxis.from_bounds(0, 5, nbin=4, name="offset")
>>> aeff = EffectiveAreaTable2D(axes=[energy_axis_true, offset_axis], data=1e10, unit="cm2")
>>> print(aeff)
EffectiveAreaTable2D
--------------------

  axes  : ['energy_true', 'offset']
  shape : (30, 4)
  ndim  : 2
  unit  : cm2
  dtype : float64

Attributes Summary

Methods Summary

from_parametrization([energy_axis_true, ...])

Create parametrized effective area.

peek([figsize])

Quick-look summary plots.

plot([ax, add_cbar, axes_loc, kwargs_colorbar])

Plot effective area image.

plot_energy_dependence([ax, offset])

Plot effective area versus energy for a given offset.

plot_offset_dependence([ax, energy])

Plot effective area versus offset for a given energy.

Attributes Documentation

default_unit = Unit("m2")#
required_axes = ['energy_true', 'offset']#
tag = 'aeff_2d'#

Methods Documentation

classmethod from_parametrization(energy_axis_true=None, instrument='HESS')[source]#

Create parametrized effective area.

Parametrizations of the effective areas of different Cherenkov telescopes taken from Appendix B of Abramowski et al. (2010), see https://ui.adsabs.harvard.edu/abs/2010MNRAS.402.1342A .

\[A_{eff}(E) = g_1 \left(\frac{E}{\mathrm{MeV}}\right)^{-g_2}\exp{\left(-\frac{g_3}{E}\right)}\]

This method does not model the offset dependence of the effective area, but just assumes that it is constant.

Parameters:
energy_axis_trueMapAxis, optional

Energy binning, analytic function is evaluated at log centers. Default is None.

instrument{‘HESS’, ‘HESS2’, ‘CTAO’}

Instrument name. Default is ‘HESS’.

Returns:
aeffEffectiveAreaTable2D

Effective area table.

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

Quick-look summary plots.

This method creates a figure with three subplots:

  • Energy dependence plot : effective area versus true energy for a given offset

  • Offset dependence plot : effective area versus true energy for a given offset

  • Effective area 2D map

Parameters:
figsizetuple, optional

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

plot(ax=None, add_cbar=True, axes_loc=None, kwargs_colorbar=None, **kwargs)[source]#

Plot effective area image.

Parameters:
axAxes, optional

Matplotlib axes. Default is None.

add_cbarbool, optional

Add a colorbar to the plot. Default is True.

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.

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

Plot effective area versus energy for a given offset.

Parameters:
axAxes, optional

Matplotlib axes. Default is None.

offsetlist of Angle, optional

Offset. Default is None.

kwargsdict

Forwarded to plt.plot().

Returns:
axAxes

Matplotlib axes.

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

Plot effective area versus offset for a given energy.

Parameters:
axAxes, optional

Matplotlib axes. Default is None.

energyQuantity

Energy.

**kwargsdict

Keyword argument passed to plot.

Returns:
axAxes

Matplotlib axes.

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