EventList#

class gammapy.data.EventList(table, meta=None)[source]#

Bases: object

Event list.

Event list data is stored as 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)

Note that TIME is usually sorted, but sometimes it is not. E.g. when simulating data, or processing it in certain ways. So generally any analysis code should assume TIME is not sorted.

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:
tableTable

Event list table.

metaEventListMetaData

The metadata. Default is None.

Examples

>>> from gammapy.data import EventList
>>> events = EventList.read("$GAMMAPY_DATA/cta-1dc/data/baseline/gps/gps_baseline_110380.fits")
>>> print(events)
EventList
---------

  Instrument       : None
  Telescope        : CTA
  Obs. ID          : 110380

  Number of events : 106217
  Event rate       : 59.273 1 / s

  Time start       : 59235.5
  Time stop        : 59235.52074074074

  Min. energy      : 3.00e-02 TeV
  Max. energy      : 1.46e+02 TeV
  Median energy    : 1.02e-01 TeV

  Max. offset      : 5.0 deg

Attributes Summary

altaz

ALT / AZ position computed from RA / DEC as a SkyCoord object.

altaz_frame

ALT / AZ frame as an AltAz object.

altaz_from_table

ALT / AZ position from table as a SkyCoord object.

energy

Event energies as a Quantity.

galactic

Event Galactic sky coordinates as a SkyCoord object.

galactic_median

Median position as a SkyCoord object.

is_pointed_observation

Whether observation is pointed.

observation_dead_time_fraction

Dead-time fraction as a float.

observation_live_time_duration

Live-time duration in seconds as a Quantity.

observation_time_duration

Observation time duration in seconds as a Quantity.

observation_time_start

Observation start time as a Time object.

observation_time_stop

Observation stop time as a Time object.

observatory_earth_location

Observatory location as an EarthLocation object.

offset

Event offset from the array pointing position as an Angle.

offset_from_median

Event offset from the median position as an Angle.

pointing_radec

Pointing RA / DEC sky coordinates as a SkyCoord object.

radec

Event RA / DEC sky coordinates as a SkyCoord object.

time

Event times as a Time object.

time_ref

Time reference as a Time object.

Methods Summary

check([checks])

Run checks.

copy()

Copy event list (EventList).

from_stack(event_lists, **kwargs)

Stack (concatenate) list of event lists.

map_coord(geom)

Event map coordinates for a given geometry.

peek([allsky])

Quick look plots.

plot_energy([ax])

Plot counts as a function of energy.

plot_energy_offset([ax, center])

Plot counts histogram with energy and offset axes.

plot_image([ax, allsky])

Quick look counts map sky plot.

plot_offset2_distribution([ax, center, ...])

Plot offset^2 distribution of the events.

plot_time([ax])

Plot an event rate time curve.

read(filename[, hdu, checksum])

Read from FITS file.

select_energy(energy_range)

Select events in energy band.

select_mask(mask)

Select events inside a mask (EventList).

select_offset(offset_band)

Select events in offset band.

select_parameter(parameter, band)

Select events with respect to a specified parameter.

select_rad_max(rad_max[, position])

Select energy dependent offset.

select_region(regions[, wcs])

Select events in given region.

select_row_subset(row_specifier)

Select table row subset.

select_time(time_interval)

Select events in time interval.

stack(other)

Stack with another EventList in place.

to_table_hdu([format])

Convert event list to a BinTableHDU.

write(filename[, gti, overwrite, format, ...])

Deprecated since version v1.2.

Attributes Documentation

altaz#

ALT / AZ position computed from RA / DEC as a SkyCoord object.

altaz_frame#

ALT / AZ frame as an AltAz object.

altaz_from_table#

ALT / AZ position from table as a SkyCoord object.

energy#

Event energies as a Quantity.

galactic#

Event Galactic sky coordinates as a SkyCoord object.

Always computed from RA / DEC using Astropy.

galactic_median#

Median position as a SkyCoord object.

is_pointed_observation#

Whether observation is pointed.

observation_dead_time_fraction#

Dead-time fraction as a float.

This is a keyword related to IACTs. Defined as dead-time over observation time.

Dead-time is defined as the time during the observation where the detector didn’t record events: http://en.wikipedia.org/wiki/Dead_time https://ui.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 as a 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_duration#

Observation time duration in seconds as a Quantity.

This is a keyword related to IACTs. The wall time, including dead-time.

observation_time_start#

Observation start time as a Time object.

observation_time_stop#

Observation stop time as a Time object.

observatory_earth_location#

Observatory location as an EarthLocation object.

offset#

Event offset from the array pointing position as an Angle.

offset_from_median#

Event offset from the median position as an Angle.

pointing_radec#

Pointing RA / DEC sky coordinates as a SkyCoord object.

radec#

Event RA / DEC sky coordinates as a SkyCoord object.

time#

Event times as a Time object.

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 as a Time object.

Methods Documentation

check(checks='all')[source]#

Run checks.

This is a generator that yields a list of dicts.

copy()[source]#

Copy event list (EventList).

classmethod from_stack(event_lists, **kwargs)[source]#

Stack (concatenate) list of event lists.

Calls vstack.

Parameters:
event_listslist

List of EventList to stack.

**kwargsdict, optional

Keyword arguments passed to vstack.

map_coord(geom)[source]#

Event map coordinates for a given geometry.

Parameters:
geomGeom

Geometry.

Returns:
coordMapCoord

Coordinates.

peek(allsky=False)[source]#

Quick look plots.

Parameters:
allskybool, optional

Whether to look at the events all-sky. Default is False.

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

Plot counts as a function of energy.

Parameters:
axAxes, optional

Matplotlib axes. Default is None

**kwargsdict, optional

Keyword arguments passed to hist.

Returns:
axAxes

Matplotlib axes.

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

Plot counts histogram with energy and offset axes.

Parameters:
axAxis, optional

Plot axis. Default is None.

centerSkyCoord, optional

Sky coord from which offset is computed. Default is None.

**kwargsdict, optional

Keyword arguments forwarded to pcolormesh.

Returns:
axAxis

Plot axis.

plot_image(ax=None, allsky=False)[source]#

Quick look counts map sky plot.

Parameters:
axAxes, optional

Matplotlib axes.

allskybool, optional

Whether to plot on an all sky geom. Default is False.

plot_offset2_distribution(ax=None, center=None, max_percentile=98, **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:
axAxes, optional

Matplotlib axes. Default is None.

centerastropy.coordinates.SkyCoord, optional

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

max_percentilefloat, optional

Define the percentile of the offset^2 distribution used to define the maximum offset^2 value. Default is 98.

**kwargsdict, optional

Extra keyword arguments are passed to hist.

Returns:
axAxes

Matplotlib axes.

Examples

Load an example event list:

>>> from gammapy.data import EventList
>>> from astropy import units as u
>>> filename = "$GAMMAPY_DATA/hess-dl3-dr1/data/hess_dl3_dr1_obs_id_023523.fits.gz"
>>> events = EventList.read(filename)
>>> #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) * u.deg ** 2
>>> 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, **kwargs)[source]#

Plot an event rate time curve.

Parameters:
axAxes, optional

Matplotlib axes. Default is None.

**kwargsdict, optional

Keyword arguments passed to errorbar.

Returns:
axAxes

Matplotlib axes.

classmethod read(filename, hdu='EVENTS', checksum=False, **kwargs)[source]#

Read from FITS file.

Format specification: EVENTS

Parameters:
filenamepathlib.Path, str

Filename

hdustr

Name of events HDU. Default is “EVENTS”.

checksumbool

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

select_energy(energy_range)[source]#

Select events in energy band.

Parameters:
energy_rangeQuantity

Energy range [energy_min, energy_max).

Returns:
event_listEventList

Copy of event list with selection applied.

Examples

>>> from astropy import units as u
>>> from gammapy.data import EventList
>>> filename = "$GAMMAPY_DATA/fermi_3fhl/fermi_3fhl_events_selected.fits.gz"
>>> event_list = EventList.read(filename)
>>> energy_range =[1, 20] * u.TeV
>>> event_list = event_list.select_energy(energy_range=energy_range)
select_mask(mask)[source]#

Select events inside a mask (EventList).

Parameters:
maskMap

Mask.

Returns:
event_listEventList

Copy of event list with selection applied.

Examples

>>> from gammapy.data import EventList
>>> from gammapy.maps import WcsGeom, Map
>>> geom = WcsGeom.create(skydir=(0,0), width=(4, 4), frame="galactic")
>>> mask = geom.region_mask("galactic;circle(0, 0, 0.5)")
>>> filename = "$GAMMAPY_DATA/cta-1dc/data/baseline/gps/gps_baseline_110380.fits"
>>> events = EventList.read(filename)
>>> masked_event = events.select_mask(mask)
>>> len(masked_event.table)
5594
select_offset(offset_band)[source]#

Select events in offset band.

Parameters:
offset_bandAngle

offset band [offset_min, offset_max).

Returns:
event_listEventList

Copy of event list with selection applied.

Examples

>>> from gammapy.data import EventList
>>> import astropy.units as u
>>> filename = "$GAMMAPY_DATA/cta-1dc/data/baseline/gps/gps_baseline_110380.fits"
>>> events = EventList.read(filename)
>>> selected_events = events.select_offset([0.3, 0.9]*u.deg)
>>> len(selected_events.table)
12688
select_parameter(parameter, band)[source]#

Select events with respect to a specified parameter.

Parameters:
parameterstr

Parameter used for the selection. Must be present in self.table.

bandtuple or astropy.units.Quantity

Minimum and maximum value for the parameter to be selected (minimum <= parameter < maximum). If parameter is not dimensionless, a Quantity is required.

Returns:
event_listEventList

Copy of event list with selection applied.

Examples

>>> from astropy import units as u
>>> from gammapy.data import EventList
>>> filename = "$GAMMAPY_DATA/fermi_3fhl/fermi_3fhl_events_selected.fits.gz"
>>> event_list = EventList.read(filename)
>>> zd = (0, 30) * u.deg
>>> event_list = event_list.select_parameter(parameter='ZENITH_ANGLE', band=zd)
>>> print(len(event_list.table))
123944
select_rad_max(rad_max, position=None)[source]#

Select energy dependent offset.

Parameters:
rad_maxRadMax2D

Rad max definition.

positionSkyCoord, optional

Center position. Default is the pointing position.

Returns:
event_listEventList

Copy of event list with selection applied.

select_region(regions, wcs=None)[source]#

Select events in given region.

Parameters:
regionsstr or Region or list of Region

Region or list of regions (pixel or sky regions accepted). A region can be defined as a string in the DS9 format as well. See http://ds9.si.edu/doc/ref/region.html for details.

wcsWCS, optional

World coordinate system transformation. Default is None.

Returns:
event_listEventList

Copy of event list with selection applied.

select_row_subset(row_specifier)[source]#

Select table row subset.

Parameters:
row_specifierslice or int or array of int

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

Returns:
event_listEventList

New event list with table row subset selected.

Examples

>>> from gammapy.data import EventList
>>> import numpy as np
>>> filename = "$GAMMAPY_DATA/cta-1dc/data/baseline/gps/gps_baseline_110380.fits"
>>> events = EventList.read(filename)
>>> #Use a boolean mask as ``row_specifier``:
>>> mask = events.table['MC_ID'] == 1
>>> events2 = events.select_row_subset(mask)
>>> print(len(events2.table))
97978
>>> #Use row index array as ``row_specifier``:
>>> idx = np.where(events.table['MC_ID'] == 1)[0]
>>> events2 = events.select_row_subset(idx)
>>> print(len(events2.table))
97978
select_time(time_interval)[source]#

Select events in time interval.

Parameters:
time_intervalastropy.time.Time

Start time (inclusive) and stop time (exclusive) for the selection.

Returns:
eventsEventList

Copy of event list with selection applied.

stack(other)[source]#

Stack with another EventList in place.

Calls vstack.

Parameters:
otherEventList

Event list to stack to self.

to_table_hdu(format='gadf')[source]#

Convert event list to a BinTableHDU.

Parameters:
formatstr, optional

Output format, currently only “gadf” is supported. Default is “gadf”.

Returns:
hduastropy.io.fits.BinTableHDU

EventList converted to FITS representation.

write(filename, gti=None, overwrite=False, format='gadf', checksum=False)[source]#

Deprecated since version v1.2: To write an EventList utilise Observation.write() with include_irfs=False

Write the event list to a FITS file.

If a GTI object is provided, it is saved into a second extension in the file.

Parameters:
filenamepathlib.Path or str

Filename.

gtiGTI, optional

Good Time Intervals object to save to the same file. Default is None.

overwritebool, optional

Overwrite existing file. Default is False.

formatstr, optional

FITS format convention. Default is “gadf”.

checksumbool, optional

When True adds both DATASUM and CHECKSUM cards to the headers written to the file. Default is False.