EventListLAT

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

Bases: gammapy.data.EventListBase

Event list for Fermi-LAT dataset.

Fermi-LAT data products https://fermi.gsfc.nasa.gov/ssc/data/analysis/documentation/Cicerone/Cicerone_Data/LAT_DP.html Data format specification (columns) https://fermi.gsfc.nasa.gov/ssc/data/analysis/documentation/Cicerone/Cicerone_Data/LAT_Data_Columns.html

For further information, see the base class: EventListBase.

Parameters
tableTable

Event list table

Examples

To load an example Fermi-LAT event list:

>>> from gammapy.data import EventListLAT
>>> filename = "$GAMMAPY_DATA/fermi-3fhl-gc/fermi-3fhl-gc-events.fits.gz"
>>> events = EventListLAT.read(filename)

Attributes Summary

energy

Event energies (Quantity).

galactic

Event Galactic sky coordinates (SkyCoord).

observation_time_end

Observation stop time (Time).

observation_time_start

Observation start time (Time).

radec

Event RA / DEC sky coordinates (SkyCoord).

time

Event times (Time).

time_ref

Time reference (Time).

Methods Summary

check(self[, checks])

Run checks.

map_coord(self, geom)

Event map coordinates for a given geometry.

plot_energy(self[, ax, ebounds])

Plot counts as a function of energy.

plot_energy_offset(self[, ax])

Plot counts histogram with energy and offset axes.

plot_image(self)

Quick look counts map sky plot.

plot_offset2_distribution(self[, ax, center])

Plot offset^2 distribution of the events.

plot_time(self[, ax])

Plots an event rate time curve.

read(filename, \*\*kwargs)

Read from FITS file.

select_energy(self, energy_band)

Select events in energy band.

select_map_mask(self, mask)

Select events inside a mask (EventList).

select_parameter(self, parameter, band)

Select events with respect to a specified parameter.

select_region(self, region[, wcs])

Select events in given region.

select_row_subset(self, row_specifier)

Select table row subset.

select_time(self, 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_time_end

Observation stop time (Time).

observation_time_start

Observation start time (Time).

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(self, checks='all')

Run checks.

This is a generator that yields a list of dicts.

map_coord(self, geom)

Event map coordinates for a given geometry.

Parameters
geomGeom

Geometry

Returns
coordMapCoord

Coordinates

plot_energy(self, ax=None, ebounds=None, **kwargs)

Plot counts as a function of energy.

plot_energy_offset(self, ax=None)

Plot counts histogram with energy and offset axes.

plot_image(self)[source]

Quick look counts map sky plot.

plot_offset2_distribution(self, ax=None, center=None, **kwargs)

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 matplotlib.pyplot.hist.

Returns
axAxes

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(self, ax=None)

Plots an event rate time curve.

Parameters
axAxes or None

Axes

Returns
axAxes

Axes

classmethod read(filename, **kwargs)

Read from FITS file.

Format specification: EVENTS

Parameters
filenamepathlib.Path, str

Filename

select_energy(self, energy_band)

Select events in energy band.

Parameters
energy_bandQuantity

Energy band [energy_min, energy_max)

Returns
event_listEventList

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_map_mask(self, mask)

Select events inside a mask (EventList).

Parameters
maskMap

Mask

select_parameter(self, parameter, band)

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 gammapy.data import EventList
>>> event_list = EventList.read('events.fits')
>>> phase_region = (0.3, 0.5)
>>> event_list = event_list.select_parameter(parameter='PHASE', band=phase_region)
select_region(self, region, wcs=None)

Select events in given region.

Parameters
regionSkyRegion or str

Sky region or string defining a sky region

wcsWCS

World coordinate system transformation

Returns
event_listEventList

Copy of event list with selection applied.

select_row_subset(self, row_specifier)

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(self, time_interval)

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.

classmethod stack(event_lists, **kwargs)

Stack (concatenate) list of event lists.

Calls vstack.

Parameters
event_listslist

list of EventList to stack