EventList

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

Bases: object

Event list.

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)

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

Attributes Summary

altaz

ALT / AZ position computed from RA / DEC (SkyCoord).

altaz_frame

ALT / AZ frame (AltAz).

altaz_from_table

ALT / AZ position from table (SkyCoord).

energy

Event energies (Quantity).

galactic

Event Galactic sky coordinates (SkyCoord).

galactic_median

Median position in radec

is_pointed_observation

Whether observation is pointed

observation_dead_time_fraction

Dead-time fraction (float).

observation_live_time_duration

Live-time duration in seconds (Quantity).

observation_time_duration

Observation time duration in seconds (Quantity).

observation_time_start

Observation start time (Time).

observation_time_stop

Observation stop time (Time).

observatory_earth_location

Observatory location (EarthLocation).

offset

Event offset from the array pointing position (Angle).

offset_from_median

Event offset from the median 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.

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([allsky])

Quick look counts map sky plot.

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

Attributes Documentation

altaz

ALT / AZ position computed from RA / DEC (SkyCoord).

altaz_frame

ALT / AZ frame (AltAz).

altaz_from_table

ALT / AZ position from table (SkyCoord).

energy

Event energies (Quantity).

galactic

Event Galactic sky coordinates (SkyCoord).

Always computed from RA / DEC using Astropy.

galactic_median

Median position in radec

is_pointed_observation

Whether observation is pointed

observation_dead_time_fraction

Dead-time fraction (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 (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 (Quantity).

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

observation_time_start

Observation start time (Time).

observation_time_stop

Observation stop time (Time).

observatory_earth_location

Observatory location (EarthLocation).

offset

Event offset from the array pointing position (Angle).

offset_from_median

Event offset from the median 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.

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

Stack (concatenate) list of event lists.

Calls vstack.

Parameters
event_listslist

list of EventList to stack

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

Whether to look at the events allsky

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

Plot counts as a function of energy.

Parameters
axAxes or None

Axes

**kwargsdict

Keyword arguments passed to hist

Returns
axAxes or None

Axes

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

Plot counts histogram with energy and offset axes

Parameters
axAxis

Plot axis

centerSkyCoord

Sky coord from which offset is computed

**kwargsdict

Keyword arguments forwarded to pcolormesh

Returns
axAxis

Plot axis

plot_image(allsky=False)[source]

Quick look counts map sky plot.

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
axAxes (optional)

Axes

centerastropy.coordinates.SkyCoord

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

**kwargs :

Extra keyword arguments are passed to hist.

Returns
axAxes

Axes

Examples

Load an example event list:

>>> from gammapy.data import EventList
>>> from astropy import units as u
>>> import matplotlib.pyplot as plt
>>> 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):

>>> plt.cla()
>>> events.plot_offset2_distribution()
<AxesSubplot:xlabel='Offset^2 (deg2)', ylabel='Counts'>

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
>>> plt.cla()
>>> events.plot_offset2_distribution(center=center, bins=bins)
<AxesSubplot:xlabel='Offset^2 (deg2)', ylabel='Counts'>

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]

Plots an event rate time curve.

Parameters
axAxes or None

Axes

**kwargsdict

Keyword arguments passed to errorbar

Returns
axAxes

Axes

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

Read from FITS file.

Format specification: EVENTS

Parameters
filenamepathlib.Path, str

Filename

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
>>> event_list = EventList.read('$GAMMAPY_DATA/fermi_3fhl/fermi_3fhl_events_selected.fits.gz')
>>> 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

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.

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

Min and max value for the parameter to be selected (min <= parameter < max). If parameter is not dimensionless you have to provide a Quantity.

Returns
event_listEventList

Copy of event list with selection applied.

Examples

>>> from astropy import units as u
>>> from gammapy.data import EventList
>>> event_list = EventList.read('$GAMMAPY_DATA/fermi_3fhl/fermi_3fhl_events_selected.fits.gz')
>>> zd = (0, 30) * u.deg
>>> event_list = event_list.select_parameter(parameter='ZENITH_ANGLE', band=zd)
select_region(regions, wcs=None)[source]

Select events in given region.

Parameters
regionsstr, Region or list of Region

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

wcsWCS

World coordinate system transformation

Returns
event_listEventList

Copy of event list with selection applied.

select_row_subset(row_specifier)[source]

Select table row subset.

Parameters
row_specifierslice, int, or array of ints

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

Returns
event_listEventList

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_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