EventListBase

class gammapy.data.EventListBase(table)[source]

Bases: object

Event list.

This class represents the base for two different event lists: - EventList: targeted for IACT event lists - EventListLAT: targeted for Fermi-LAT event lists

Event list data is stored in table (Table) data member.

The most important reconstructed event parameters are available as the following columns:

  • TIME - Mission elapsed time (sec)
  • RA, DEC - ICRS system position (deg)
  • ENERGY - Energy (usually MeV for Fermi and TeV for IACTs)

Other optional (columns) that are sometimes useful for high-level analysis:

  • GLON, GLAT - Galactic coordinates (deg)
  • DETX, DETY - Field of view coordinates (deg)

Note that when reading data for analysis you shouldn’t use those values directly, but access them via properties which create objects of the appropriate class:

Parameters:

table : Table

Event list table

Attributes Summary

energy Event energies (Quantity).
galactic Event Galactic sky coordinates (SkyCoord).
observation_dead_time_fraction Dead-time fraction (float).
observation_live_time_duration Live-time duration in seconds (Quantity).
observation_time_end Observation stop time (Time).
observation_time_start Observation start time (Time).
offset Event offset from the array pointing position (Angle).
pointing_radec Pointing RA / DEC sky coordinates (SkyCoord).
radec Event RA / DEC sky coordinates (SkyCoord).
time Event times (Time).
time_ref Time reference (Time)

Methods Summary

check([checks]) Run checks.
filter_circular_region(region) Create selection mask for event in given circular regions.
plot_energy([ax, ebounds]) Plot counts as a function of energy.
plot_energy_offset([ax]) Plot counts histogram with energy and offset axes.
plot_offset2_distribution([ax, center]) Plot offset^2 distribution of the events.
plot_time([ax]) Plots an event rate time curve.
read(filename, **kwargs) Read from FITS file.
select_circular_region(region) Select events in circular regions.
select_energy(energy_band) Select events in energy band.
select_row_subset(row_specifier) Select table row subset.
select_sky_box(lon_lim, lat_lim[, frame]) Select events in sky box.
select_sky_cone(center, radius) Select events in sky circle.
select_sky_ring(center, inner_radius, …) Select events in ring region on the sky.
select_time(time_interval) Select events in time interval.
stack(event_lists, **kwargs) Stack (concatenate) list of event lists.

Attributes Documentation

energy

Event energies (Quantity).

galactic

Event Galactic sky coordinates (SkyCoord).

Always computed from RA / DEC using Astropy.

observation_dead_time_fraction

Dead-time fraction (float).

Defined as dead-time over observation time.

Dead-time is defined as the time during the observation where the detector didn’t record events: https://en.wikipedia.org/wiki/Dead_time https://adsabs.harvard.edu/abs/2004APh….22..285F

The dead-time fraction is used in the live-time computation, which in turn is used in the exposure and flux computation.

observation_live_time_duration

Live-time duration in seconds (Quantity).

The dead-time-corrected observation time.

  • In Fermi-LAT it is automatically provided in the header of the event list.
  • In IACTs is computed as t_live = t_observation * (1 - f_dead)

where f_dead is the dead-time fraction.

observation_time_end

Observation stop time (Time).

observation_time_start

Observation start time (Time).

offset

Event offset from the array pointing position (Angle).

pointing_radec

Pointing RA / DEC sky coordinates (SkyCoord).

radec

Event RA / DEC sky coordinates (SkyCoord).

time

Event times (Time).

Notes

Times are automatically converted to 64-bit floats. With 32-bit floats times will be incorrect by a few seconds when e.g. adding them to the reference time.

time_ref

Time reference (Time)

Methods Documentation

check(checks='all')[source]

Run checks.

This is a generator that yields a list of dicts.

filter_circular_region(region)[source]

Create selection mask for event in given circular regions.

TODO: Extend to support generic regions

Parameters:

region : list of SkyRegion

List of sky regions

Returns:

index_array : numpy.ndarray

Index array of selected events

plot_energy(ax=None, ebounds=None, **kwargs)[source]

Plot counts as a function of energy.

plot_energy_offset(ax=None)[source]

Plot counts histogram with energy and offset axes.

plot_offset2_distribution(ax=None, center=None, **kwargs)[source]

Plot offset^2 distribution of the events.

The distribution shown in this plot is for this quantity:

offset = center.separation(events.radec).deg
offset2 = offset ** 2

Note that this method is just for a quicklook plot.

If you want to do computations with the offset or offset^2 values, you can use the line above. As an example, here’s how to compute the 68% event containment radius using numpy.percentile:

import numpy as np
r68 = np.percentile(offset, q=68)
Parameters:

ax : Axes (optional)

Axes

center : astropy.coordinates.SkyCoord

Center position for the offset^2 distribution. Default is the observation pointing position.

**kwargs :

Extra keyword arguments are passed to matplotlib.pyplot.hist.

Returns:

ax : Axes

Axes

Examples

Load an example event list:

>>> from gammapy.data import EventList
>>> events = EventList.read('$GAMMAPY_DATA/hess-dl3-dr1/data/hess_dl3_dr1_obs_id_023523.fits.gz')

Plot the offset^2 distribution wrt. the observation pointing position (this is a commonly used plot to check the background spatial distribution):

>>> events.plot_offset2_distribution()

Plot the offset^2 distribution wrt. the Crab pulsar position (this is commonly used to check both the gamma-ray signal and the background spatial distribution):

>>> import numpy as np
>>> from astropy.coordinates import SkyCoord
>>> center = SkyCoord(83.63307, 22.01449, unit='deg')
>>> bins = np.linspace(start=0, stop=0.3 ** 2, num=30)
>>> events.plot_offset2_distribution(center=center, bins=bins)

Note how we passed the bins option of matplotlib.pyplot.hist to control the histogram binning, in this case 30 bins ranging from 0 to (0.3 deg)^2.

plot_time(ax=None)[source]

Plots an event rate time curve.

Parameters:

ax : Axes or None

Axes

Returns:

ax : Axes

Axes

classmethod read(filename, **kwargs)[source]

Read from FITS file.

Format specification: EVENTS

Parameters:

filename : Path, str

Filename

select_circular_region(region)[source]

Select events in circular regions.

TODO: Extend to support generic regions

Parameters:

region : CircleSkyRegion or list of CircleSkyRegion

(List of) sky region(s)

Returns:

event_list : EventList

Copy of event list with selection applied.

select_energy(energy_band)[source]

Select events in energy band.

Parameters:

energy_band : Quantity

Energy band [energy_min, energy_max)

Returns:

event_list : EventList

Copy of event list with selection applied.

Examples

>>> from astropy.units import Quantity
>>> from gammapy.data import EventList
>>> event_list = EventList.read('events.fits')
>>> energy_band = Quantity([1, 20], 'TeV')
>>> event_list = event_list.select_energy()
select_row_subset(row_specifier)[source]

Select table row subset.

Parameters:

row_specifier : slice, int, or array of ints

Specification for rows to select, passed on to self.table[row_specifier].

Returns:

event_list : EventList

New event list with table row subset selected

Examples

Use a boolean mask as row_specifier:

mask = events.table[‘FOO’] > 42 events2 = events.select_row_subset(mask)

Use row index array as row_specifier:

idx = np.where(events.table[‘FOO’] > 42)[0] events2 = events.select_row_subset(idx)
select_sky_box(lon_lim, lat_lim, frame='icrs')[source]

Select events in sky box.

TODO: move gammapy.catalog.select_sky_box to gammapy.utils.

select_sky_cone(center, radius)[source]

Select events in sky circle.

Parameters:

center : SkyCoord

Sky circle center

radius : Angle

Sky circle radius

Returns:

event_list : EventList

Copy of event list with selection applied.

select_sky_ring(center, inner_radius, outer_radius)[source]

Select events in ring region on the sky.

Parameters:

center : SkyCoord

Sky ring center

inner_radius, outer_radius : Angle

Sky ring inner and outer radius

Returns:

event_list : EventList

Copy of event list with selection applied.

select_time(time_interval)[source]

Select events in time interval.

classmethod stack(event_lists, **kwargs)[source]

Stack (concatenate) list of event lists.

Calls vstack.

Parameters:

event_lists : list

list of EventList to stack