Data Module

The classes in gbm.data are interfaces to the data that GBM produces.

Primary GBM Data Classes

The primary data classes for GBM data:

Cspec

class gbm.data.Cspec[source]

Bases: gbm.data.data.DataFile, gbm.data.phaii.PHAII

Class for CSPEC data.

Attributes
  • data (TimeEnergyBins) – The PHAII data

  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • energy_range (float, float) – The energy range of the data

  • filename (str) – The filename

  • full_path (str) – The full path+filename

  • gti ([(float, float), …]) – A list of tuples defining the good time intervals

  • headers (dict) – The headers for each extension of the PHAII

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

  • numchans (int) – The number of energy channels

  • trigtime (float) – The trigger time of the data, if available.

  • time_range (float, float) – The time range of the data

classmethod from_data(data, detector=None, **kwargs)[source]

Create a PHAII object from TimeEnergyBins data object.

Parameters
  • data (TimeEnergyBins) – The PHAII data

  • gti ([(float, float), ...], optional) – The list of tuples representing the good time intervals (start, stop). If omitted, the GTI is assumed to be [(tstart, tstop)].

  • trigtime (float, optional) – The trigger time, if applicable. If provided, the data times will be shifted relative to the trigger time.

  • detector (str, optional) – The detector that produced the data

  • object (str, optional) – The object being observed

  • ra_obj (float, optional) – The RA of the object

  • dec_obj (float, optional) – The Dec of the object

  • err_rad (float, optional) – The localization error radius of the object

Returns

PHAII – The newly created PHAII object

get_exposure(time_ranges=None)

Calculate the total exposure of a time range or time ranges of data

Parameters

time_ranges ([(float, float), ...], optional) – The time range or time ranges over which to calculate the exposure. If omitted, calculates the total exposure of the data.

Returns

float – The exposure of the time selections

classmethod open(filename)

Open a PHAII FITS file and return the PHAII object

Parameters

filename (str) – The filename of the FITS file

Returns

PHAII – The PHAII object

rebin_energy(method, *args, energy_range=(None, None))

Rebin the PHAII in energy given a rebinning method. Produces a new PHAII object.

Parameters
  • method (<function>) – The rebinning function

  • *args – Arguments to be passed to the rebinning function

  • energy_range ((float, float), optional) – The starting and ending energy to rebin. If omitted, uses the full range of data. Setting start or end to None will use the data from the beginning or end of the data, respectively.

Returns

PHAII: A new PHAII object with the requested rebinned data

rebin_time(method, *args, time_range=(None, None))

Rebin the PHAII in time given a rebinning method. Produces a new PHAII object.

Parameters
  • method (<function>) – The rebinning function

  • *args – Arguments to be passed to the rebinning function

  • time_range ((float, float), optional) – The starting and ending time to rebin. If omitted, uses the full range of data. Setting start or end to None will use the data from the beginning or end of the data, respectively.

Returns

PHAII: A new PHAII object with the requested rebinned data

set_properties(detector=None, trigtime=None, tstart=None, datatype=None, extension=None, meta=None)

Set the properties of the data file. Useful for creating new data files. A standardized filename will be built from the parameters

Parameters
  • detector (str, optional) – The detector the data file belongs to

  • trigtime (float, optional) – The trigger time, if applicable

  • tstart (float) – The start time of the data file. Must be set if trigtime is not.

  • datatype (str, optional) – The type of data the file contains

  • extension (str, optional) – The extension of the data file

  • meta (str, optional) – A metadata atttribute to be added to the filename

slice_energy(energy_ranges)

Slice the PHAII by one or more energy range. Produces a new PHAII object.

Parameters

energy_ranges ([(float, float), ...]) – The energy ranges to slice the data to.

Returns

PHAII – A new PHAII object with the requested energy ranges

slice_time(time_ranges)

Slice the PHAII by one or more time range. Produces a new PHAII object.

Parameters

time_ranges ([(float, float), ...]) – The time ranges to slice the data to.

Returns

PHAII – A new PHAII object with the requested time ranges

to_lightcurve(time_range=None, energy_range=None, channel_range=None)

Integrate the PHAII data over energy to produce a lightcurve

Parameters
  • time_range ((float, float), optional) – The time range of the lightcurve. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the lightcurve. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the lightcurve. If omitted, uses the entire energy range of the data.

Returns

TimeBins – The lightcurve data

to_pha(time_ranges=None, energy_range=None, channel_range=None, **kwargs)

Integrate the PHAII data over time to produce a PHA object

Parameters
  • time_ranges ([(float, float), ...], optional) – The time range of the spectrum. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the spectrum. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the spectrum. If omitted, uses the entire energy range of the data.

  • **kwargs – Options passed to gbm.data.PHA.from_data()

Returns

PHA – The single spectrum PHA

to_spectrum(time_range=None, energy_range=None, channel_range=None)

Integrate the PHAII data over time to produce a count spectrum

Parameters
  • time_range ((float, float), optional) – The time range of the spectrum. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the spectrum. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the spectrum. If omitted, uses the entire energy range of the data.

Returns

EnergyBins – The count spectrum data

write(directory, filename=None)

Writes the data to a FITS file.

Parameters
  • directory (str) – The directory to write the file

  • filename (str, optional) – If set, will override the standardized name


Ctime

class gbm.data.Ctime[source]

Bases: gbm.data.data.DataFile, gbm.data.phaii.PHAII

Class for CTIME data.

Attributes
  • data (TimeEnergyBins) – The PHAII data

  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • energy_range (float, float) – The energy range of the data

  • filename (str) – The filename

  • full_path (str) – The full path+filename

  • gti ([(float, float), …]) – A list of tuples defining the good time intervals

  • headers (dict) – The headers for each extension of the PHAII

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

  • numchans (int) – The number of energy channels

  • trigtime (float) – The trigger time of the data, if available.

  • time_range (float, float) – The time range of the data

classmethod from_data(data, detector=None, **kwargs)[source]

Create a PHAII object from TimeEnergyBins data object.

Parameters
  • data (TimeEnergyBins) – The PHAII data

  • gti ([(float, float), ...], optional) – The list of tuples representing the good time intervals (start, stop). If omitted, the GTI is assumed to be [(tstart, tstop)].

  • trigtime (float, optional) – The trigger time, if applicable. If provided, the data times will be shifted relative to the trigger time.

  • detector (str, optional) – The detector that produced the data

  • object (str, optional) – The object being observed

  • ra_obj (float, optional) – The RA of the object

  • dec_obj (float, optional) – The Dec of the object

  • err_rad (float, optional) – The localization error radius of the object

Returns

PHAII – The newly created PHAII object

get_exposure(time_ranges=None)

Calculate the total exposure of a time range or time ranges of data

Parameters

time_ranges ([(float, float), ...], optional) – The time range or time ranges over which to calculate the exposure. If omitted, calculates the total exposure of the data.

Returns

float – The exposure of the time selections

classmethod open(filename)

Open a PHAII FITS file and return the PHAII object

Parameters

filename (str) – The filename of the FITS file

Returns

PHAII – The PHAII object

rebin_energy(method, *args, energy_range=(None, None))

Rebin the PHAII in energy given a rebinning method. Produces a new PHAII object.

Parameters
  • method (<function>) – The rebinning function

  • *args – Arguments to be passed to the rebinning function

  • energy_range ((float, float), optional) – The starting and ending energy to rebin. If omitted, uses the full range of data. Setting start or end to None will use the data from the beginning or end of the data, respectively.

Returns

PHAII: A new PHAII object with the requested rebinned data

rebin_time(method, *args, time_range=(None, None))

Rebin the PHAII in time given a rebinning method. Produces a new PHAII object.

Parameters
  • method (<function>) – The rebinning function

  • *args – Arguments to be passed to the rebinning function

  • time_range ((float, float), optional) – The starting and ending time to rebin. If omitted, uses the full range of data. Setting start or end to None will use the data from the beginning or end of the data, respectively.

Returns

PHAII: A new PHAII object with the requested rebinned data

set_properties(detector=None, trigtime=None, tstart=None, datatype=None, extension=None, meta=None)

Set the properties of the data file. Useful for creating new data files. A standardized filename will be built from the parameters

Parameters
  • detector (str, optional) – The detector the data file belongs to

  • trigtime (float, optional) – The trigger time, if applicable

  • tstart (float) – The start time of the data file. Must be set if trigtime is not.

  • datatype (str, optional) – The type of data the file contains

  • extension (str, optional) – The extension of the data file

  • meta (str, optional) – A metadata atttribute to be added to the filename

slice_energy(energy_ranges)

Slice the PHAII by one or more energy range. Produces a new PHAII object.

Parameters

energy_ranges ([(float, float), ...]) – The energy ranges to slice the data to.

Returns

PHAII – A new PHAII object with the requested energy ranges

slice_time(time_ranges)

Slice the PHAII by one or more time range. Produces a new PHAII object.

Parameters

time_ranges ([(float, float), ...]) – The time ranges to slice the data to.

Returns

PHAII – A new PHAII object with the requested time ranges

to_lightcurve(time_range=None, energy_range=None, channel_range=None)

Integrate the PHAII data over energy to produce a lightcurve

Parameters
  • time_range ((float, float), optional) – The time range of the lightcurve. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the lightcurve. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the lightcurve. If omitted, uses the entire energy range of the data.

Returns

TimeBins – The lightcurve data

to_pha(time_ranges=None, energy_range=None, channel_range=None, **kwargs)

Integrate the PHAII data over time to produce a PHA object

Parameters
  • time_ranges ([(float, float), ...], optional) – The time range of the spectrum. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the spectrum. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the spectrum. If omitted, uses the entire energy range of the data.

  • **kwargs – Options passed to gbm.data.PHA.from_data()

Returns

PHA – The single spectrum PHA

to_spectrum(time_range=None, energy_range=None, channel_range=None)

Integrate the PHAII data over time to produce a count spectrum

Parameters
  • time_range ((float, float), optional) – The time range of the spectrum. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the spectrum. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the spectrum. If omitted, uses the entire energy range of the data.

Returns

EnergyBins – The count spectrum data

write(directory, filename=None)

Writes the data to a FITS file.

Parameters
  • directory (str) – The directory to write the file

  • filename (str, optional) – If set, will override the standardized name


GbmHealPix

class gbm.data.GbmHealPix[source]

Bases: gbm.data.localization.HealPix

Class for GBM HEALPix localization files.

Attributes
  • <detector_name>_pointing (float, float) – The RA, Dec of the detector pointing (e.g. GbmHealPix.n0_pointing)

  • centroid (float, float) – The RA and Dec of the highest probability pixel

  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • filename (str) – The filename

  • full_path (str) – The full path+filename

  • geo_location (float, float) – The geocenter RA, Dec at trigtime

  • geo_probability (float) – The amount of localization probability on the Earth

  • geo_radius (float) – The apparent Earth radius as observed by Fermi

  • headers (dict) – The headers for each extension

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

  • npix (int) – Number of pixels in the HEALPix map

  • nside (int) – The HEALPix resolution

  • quaternion (np.array) – The spacecraft attitude quaternion

  • pixel_area (float) – The area of each pixel in square degrees

  • scpos (np.array) – The spacecraft position in Earth inertial coordinates

  • sun_location (float, float) – The Sun RA, Dec at trigtime

  • trigtime (float) – The time corresponding to the localization

area(clevel)

Calculate the sky area contained within a given confidence region

Parameters

clevel (float) – The localization confidence level (valid range 0-1)

Returns

float – The area contained in square degrees

confidence(ra, dec)

Calculate the localization confidence level for a given point. This function interpolates the map at the requested point rather than providing the value at the nearest pixel center.

Parameters
  • ra (float) – The RA

  • dec (float) – The Dec

Returns

float – The localization confidence level

confidence_region_path(clevel, numpts_ra=360, numpts_dec=180)

Return the bounding path for a given confidence region

Parameters
  • clevel (float) – The localization confidence level (valid range 0-1)

  • numpts_ra (int, optional) – The number of grid points along the RA axis. Default is 360.

  • numpts_dec (int, optional) – The number of grid points along the Dec axis. Default is 180.

Returns

[(np.array, np.array), …] – A list of RA, Dec points, where each item in the list is a continuous closed path.

convolve(model, *args)

Convolve the map with a model kernel. The model can be a Gaussian kernel or any mixture of Gaussian kernels. Uses healpy.smoothing.

An example of a model kernel with a 50%/50% mixture of two Gaussians, one with a 1-deg width, and the other with a 3-deg width:

def gauss_mix_example():
    sigma1 = np.deg2rad(1.0)
    sigma2 = np.deg2rad(3.0)
    frac1 = 0.50
    return ([sigma1, sigma2], [frac1])
Parameters
  • model (<function>) – The function representing the model kernel

  • *args – Arguments to be passed to the model kernel function

Returns

HealPix – A new HealPix object that is a result of the convolution with the model kernel

classmethod from_annulus(center_ra, center_dec, radius, sigma, nside=None, **kwargs)

Create a HealPix object of a Gaussian-width annulus

Parameters
  • center_ra (float) – The RA of the center of the annulus

  • center_dec (float) – The Dec of the center of the annulus

  • radius (float) – The radius of the annulus, in degrees, measured to the center of the of the annulus

  • sigma (float) – The Gaussian standard deviation width of the annulus, in degrees

  • nside (int, optional) – The nside of the HEALPix to make. By default, the nside is automatically determined by the sigma width. Set this argument to override the default.

  • **kwargs – Options to pass to from_data()

Returns

HealPix – The HEALPix annulus

classmethod from_chi2grid(chi2grid, nside=128, tcat=None)[source]

Create a GbmHealPix object from a chi2grid object

Parameters
  • (class (chi2grid) – Chi2Grid): The chi2grid object containing the chi-squared/log-likelihood info

  • nside (int, optional) – The nside resolution to use. Default is 128

  • tcat (Tcat, optional) – The associated Tcat to fill out the primary header info

Returns

GbmHealPix – The GBM HEALPix localization

classmethod from_data(prob_arr, sig_arr, tcat=None, trigtime=None, quaternion=None, scpos=None)[source]

Create a HealPix object from healpix arrays and optional metadata

Parameters
  • prob_arr (np.array) – The HEALPix array containing the probability/pixel

  • sig_arr (np.array) – The HEALPix array containing the signficance

  • tcat (Tcat, optional) – The associated Tcat to fill out the primary header info

  • trigtime (float, optional) – The time corresponding to the localization

  • quaternion (np.array, optional) – The associated spacecraft quaternion used to determine the detector pointings in equatorial coordinates

  • scpos (np.array, optional) – The associated spacecraft position in Earth inertial coordinates used to determine the geocenter location in equatorial coordinates

Returns

GbmHealPix – The HEALPix localization

classmethod from_gaussian(center_ra, center_dec, sigma, nside=None, **kwargs)

Create a HealPix object of a Gaussian

Parameters
  • center_ra (float) – The RA of the center of the Gaussian

  • center_dec (float) – The Dec of the center of the Gaussian

  • sigma (float) – The Gaussian standard deviation, in degrees

  • nside (int, optional) – The nside of the HEALPix to make. By default, the nside is automatically determined by the sigma of the Gaussian. Set this argument to override the default.

  • **kwargs – Options to pass to from_data()

Returns

HealPix – The HEALPix Gaussian

classmethod from_vertices(ra_pts, dec_pts, nside=64, **kwargs)

Create a HealPix object from a list of RA, Dec vertices. The probability within the vertices will be distributed uniformly and zero probability outside the vertices.

Parameters
  • ra_pts (np.array) – The array of RA coordinates

  • dec_pts (np.array) – The array of Dec coordinates

  • nside (int, optional) – The nside of the HEALPix to make. Default is 64.

  • **kwargs – Options to pass to from_data()

Returns

HealPix – The HEALPix object

classmethod multiply(healpix1, healpix2, primary=1, output_nside=128)[source]

Multiply two GbmHealPix maps and return a new GbmHealPix object

Note

Either healpix1 or healpix2 can be a non-GbmHealPix object, however at least one of them must be a GbmHealPix object and the primary argument must be set to the appropriate GbmHealPix object otherwise a TypeError will be raised.

Parameters
  • healpix1 (HealPix or GbmHealPix) – One of the HEALPix maps to multiply

  • healpix2 (HealPix or GbmHealPix) – The other HEALPix map to multiply

  • primary (int, optional) – If 1, use the first map header information, or if 2, use the second map header information. Default is 1.

  • output_nside (int, optional) – The nside of the multiplied map. Default is 128.

Returns

GbmHealPix: The multiplied map

observable_fraction(healpix)[source]

The observable fraction of a healpix probability region on the sky. Non-observable regions are ones that are behind the Earth.

Parameters

healpix (HealPix) – The healpix region for which to calculate the observable fraction.

Returns

float – The fraction of the map (based on probability) that is observable.

classmethod open(filename)[source]

Open a GBM HEALPix FITS file and return the GbmHealPix object

Parameters

filename (str) – The filename of the FITS file

Returns

GbmHealPix – The GBM HEALPix localization

prob_array(numpts_ra=360, numpts_dec=180, sqdegrees=True, sig=False)

Return the localization probability mapped to a grid on the sky

Parameters
  • numpts_ra (int, optional) – The number of grid points along the RA axis. Default is 360.

  • numpts_dec (int, optional) – The number of grid points along the Dec axis. Default is 180.

  • sqdegrees (bool, optional) – If True, the prob_array is in units of probability per square degrees, otherwise in units of probability per pixel. Default is True

  • sig (bool, optional) – Set True to retun the significance map on a grid instead of the probability. Default is False.

Returns

3-tuple containing

  • np.array: The probability (or significance) array with shape (numpts_dec, numpts_ra)

  • np.array: The RA grid points

  • np.array: The Dec grid points

probability(ra, dec, per_pixel=False)

Calculate the localization probability at a given point. This function interpolates the map at the requested point rather than providing the vale at the nearest pixel center.

Parameters
  • ra (float) – The RA

  • dec (float) – The Dec

  • per_pixel (bool, optional) – If True, return probability per pixel, otherwise return probability per square degree. Default is False.

Returns

float – The localization probability

region_probability(healpix, prior=0.5)[source]

The probability that the localization is associated with the localization region from another map. This is calculated against the null hypothesis that the two maps represent unassociated sources:

\(P(A | \mathcal{I}) = \frac{P(\mathcal{I} | A) \ P(A)} {P(\mathcal{I} | A) \ P(A) + P(\mathcal{I} | \neg A) \ P(\neg A)}\)

where

  • \(P(\mathcal{I} | A)\) is the integral over the overlap of the two maps once the Earth occultation has been removed for this map.

  • \(P(\mathcal{I} | \neg A)\) is the integral over the overlap of this map with a uniform distribution on the sky (i.e. the probability the localization is associated with a random point on the sky)

  • \(P(A)\) is the prior probability that this localization is associated with the other HEALPix map.

Note

The localization region of this map overlapping the Earth will be removed and the remaining unocculted region is used for the calculation. The other map is assumed to have no exclusionary region.

Parameters
  • healpix (HealPix) – The healpix map for which to calculate the spatial association

  • prior (float, optional) – The prior probability that the localization is associated with the source. Default is 0.5

Returns

float – The probability that the two HEALPix maps are associated.

classmethod remove_earth(healpix)[source]

Return a new GbmHealPix with the probability on the Earth masked out. The remaining probability on the sky is renormalized.

Note

The geo_location attribute must be available to use this function

Parameters

healpix (GbmHealPix) – The map for which the Earth will be removed

Returns

GbmHealPix – GBM HEALPix localization

set_properties(detector=None, trigtime=None, tstart=None, datatype=None, extension=None, meta=None)

Set the properties of the data file. Useful for creating new data files. A standardized filename will be built from the parameters

Parameters
  • detector (str, optional) – The detector the data file belongs to

  • trigtime (float, optional) – The trigger time, if applicable

  • tstart (float) – The start time of the data file. Must be set if trigtime is not.

  • datatype (str, optional) – The type of data the file contains

  • extension (str, optional) – The extension of the data file

  • meta (str, optional) – A metadata atttribute to be added to the filename

source_probability(ra, dec, prior=0.5)[source]

The probability that the GbmHealPix localization is associated with a known point location. This is calculated against the null hypothesis that the localization originates from an unassociated random source that has equal probability of origination anywhere in the sky:

\(P(A | \mathcal{I}) = \frac{P(\mathcal{I} | A) \ P(A)} {P(\mathcal{I} | A) \ P(A) + P(\mathcal{I} | \neg A) \ P(\neg A)}\)

where

  • \(P(\mathcal{I} | A)\) is the probability of the localization at the point source once the Earth occultation has been removed

  • \(P(\mathcal{I} | \neg A)\) is the probability per pixel assuming a uniform distribution on the sky (i.e. the probability the localization is associated with a random point on the sky)

  • \(P(A)\) is the prior probability that the localization is associated with the point source

Note

If the point source is behind the Earth, then it is assumed that GBM could not observe it, therefore the probability will be zero.

Parameters
  • ra (float) – The RA of the known source location

  • dec (float) – The Dec of the known source location

  • prior (float, optional) – The prior probability that the localization is associated with the source. Default is 0.5

Returns

float – The probability that the localization is spatially associated with the point source

write(directory, filename=None)[source]

Write the GbmHealPix object to a FITS file

Parameters
  • directory (str) – The directory to write to

  • filename (str, optional) – The filename of the FITS file


PosHist

class gbm.data.PosHist[source]

Bases: gbm.data.data.DataFile

Class for Fermi Position History

Attributes
  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • filename (str) – The filename

  • full_path (str) – The full path+filename

  • gti ([(float, float), …]) – A list of time intervals where Fermi is outside the SAA

  • headers (dict) – The headers for each extension of the poshist file

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

  • sun_occulted ([(float, float), …]) – A list of time intervals where the sun is occulted by the Earth

  • time_range (float, float) – The time range of the poshist

detector_angle(ra, dec, det, times)[source]

Determine the angle between a sky location and a GBM detector

Parameters
  • ra (float) – The RA of a sky position

  • dec (float) – The Dec of a sky position

  • det (str) – The GBM detector

  • times (float or np.array) – Time(s) in MET

Returns

np.array – The angle, in degrees, between the position and detector pointing

detector_pointing(det, times)[source]

Retrieve the pointing of a GBM detector in equatorial coordinates

Parameters
  • det (str) – The GBM detector

  • times (float or np.array) – Time(s) in MET

Returns

(np.array, np.array) – The RA, Dec of the detector pointing

get_altitude(times)[source]

Retrieve the altitude of Fermi in orbit

Parameters

times (float or np.array) – Time(s) in MET

Returns

np.array – An array of altitudes at the requested times

get_angular_velocity(times)[source]

Retrieve the angular velocity of Fermi

Parameters

times (float or np.array) – Time(s) in MET

Returns

np.array: A (3, n) array of angular velocity components at n times

get_earth_radius(times)[source]

Retrieve the angular radius of the Earth as observed by Fermi

Parameters

times (float or np.array) – Time(s) in MET

Returns

np.array – An array of the apparent angular radius of the Earth

get_eic(times)[source]

Retrieve the position of Fermi in Earth inertial coordinates

Parameters

times (float or np.array) – Time(s) in MET

Returns

np.array – A (3, n) array of coordinates at n times

get_geocenter_radec(times)[source]

Retrieve the Equatorial position of the Geocenter as observed by Fermi

Parameters

times (float or np.array) – Time(s) in MET

Returns

(np.array, np.array) – The RA and Dec of the geocenter

get_latitude(times)[source]

Retrieve the position of Fermi in Earth latitude

Parameters

times (float or np.array) – Time(s) in MET

Returns

np.array – An array of latitudes at the requested times

get_longitude(times)[source]

Retrieve the position of Fermi in Earth East longitude

Parameters
  • times (float or np.array) – Time(s) in MET

  • np.array – An array of East longitudes at the requested times

get_mcilwain_l(times)[source]

Retrieve the approximate McIlwain L value as determined by the orbital position of Fermi

Parameters

times (float or np.array) – Time(s) in MET

Returns

np.array – An array of McIlwain L values

get_quaternions(times)[source]

Retrieve the Fermi attitude quaternions

Parameters

times (float or np.array) – Time(s) in MET

Returns

np.array – A (4, n) array of quaternions at n times

get_saa_passage(times)[source]

Determine if Fermi is in an SAA passage

Parameters

times (float or np.array) – Time(s) in MET

Returns

np.array(dtype=bool) – A boolean array, where True indicates Fermi is in SAA

get_sun_visibility(times)[source]

Determine if the sun is visible (not occulted)

Parameters

times (float or np.array) – Time(s) in MET

Returns

np.array(dtype=bool) – A boolean array, where True indicates the sun is visible

get_velocity(times)[source]

Retrieve the velocity of Fermi in Earth inertial coordinates

Parameters

times (float or np.array) – Time(s) in MET

Returns

np.array – A (3, n) array of velocity components at n times

location_visible(ra, dec, times)[source]

Determine if a sky location is visible or occulted by the Earth

Parameters
  • ra (float) – The RA of a sky position

  • dec (float) – The Dec of a sky position

  • times (float or np.array) – Time(s) in MET

Returns

np.array(dtype=bool) – A boolean array where True indicates the position is visible.

classmethod merge(poshists)[source]

Merge multiple PosHists into a single PosHist object

Parameters

poshists (list of PosHist) – List of PosHist objects

Returns

PosHist – The merged PosHist object

classmethod open(filename)[source]

Open and read a position history file

Parameters

filename (str) – The filename of the FITS file

Returns

PosHist – The PosHist object

to_equatorial(fermi_az, fermi_zen, times)[source]

Convert a position in spacecraft coordinates to equatorial coordinates

Parameters
  • fermi_az (float) – The Fermi azimuth of a position

  • fermi_zen (float) – The Fermi zenith of a position

  • times (float or np.array) – Time(s) in MET

Returns

(np.array, np.array) – The RA, Dec of the position

to_fermi_frame(ra, dec, times)[source]

Convert an equatorial position to a position in spacecraft coordinates

Parameters
  • ra (float) – The RA of a sky position

  • dec (float) – The Dec of a sky position

  • times (float or np.array) – Time(s) in MET

Returns

(np.array, np.array) – The Fermi azimuth, zenith


RSP

class gbm.data.RSP[source]

Bases: gbm.data.data.DataFile

Class for single- and multi-DRM response files

Attributes
  • channel_centroids (np.array) – Geometric centroids of the energy channel bins

  • channel_widths (np.array) – Widths of the energy channel bins

  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • ebounds (np.recarray) – The energy channel bounds

  • filename (str) – The filename

  • full_path (str) – The full path+filename

  • headers (dict) – The headers of the RSP file

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

  • numchans (int) – Number of energy channels

  • numdrms (int) – Number of DRMs

  • numebins (int) – Number of photon bins

  • photon_bins (np.array, np.array) – The photon bin edges

  • photon_bin_centroids (np.array) – Geometric centroids of the photon bins

  • photon_bin_widths (np.array) – Widths of the photon bins

  • trigtime (float) – The trigger time, if applicable

  • tcent (np.array) – The center times of the intervals over which the DRMs are valid

  • tstart (np.array) – The start times of the intervals over which the DRMs are valid

  • tstop (np.array) – The end times of the intervals over which the DRMs are valid

channel_effective_area(index=None, atime=None)[source]

Returns the effective area as a function of recorded channel energy (integrated over incident photon bins).

Parameters
  • index (int, optional) – The DRM index

  • atime (float, optional) – The time corresponding to a DRM

Returns

Bins – A Bins object containing the effective area histogram

drm(index)[source]

Return the DRM based on the index into the DRM list

Parameters

index (int) – The DRM index

Returns

np.array – The requested DRM

drm_index(time_range)[source]

Return the indices of the DRMs that cover the requested time range

If the time range is before (after) the times covered by the DRMs in the RSP object, then the first (last) index will be returned.

Parameters

time_range (float, float) – The time range

Returns

np.array – An array of the indices of the DRMs covering the time range

effective_area(energy, index=None, atime=None, interp_kind='linear')[source]

Calculate the effective area at a given energy or an array of energies. This function interpolates the DRM (in log10(cm^2)) to provide the effective area for any energy within the photon energy bounds of the DRM.

Parameters
  • energy (float or np.array) – The photon energy or energies at which to calculate the effective area.

  • index (int, optional) – The DRM index. Default is 0

  • atime (float, optional) – The time corresponding to a DRM. The nearest DRM to the time will be chosen

  • interp_kind (str, optional) – The kind of interpolation to be passed to scipy.interp1d. Default is ‘linear’.

Returns

(np.array) – The effective area

extract_drm(index=None, atime=None)[source]

Extract a single DRM from a multi-DRM RSP2 object and return a new single-DRM RSP object.

Parameters
  • index (int, optional) – The DRM index to retrieve

  • atime (float, optional) – The time corresponding to a DRM. The nearest DRM to the time will be chosen

Returns

RSP – The single-DRM RSP object

fold_spectrum(function, params, index=0, atime=None, channel_mask=None)[source]

Fold a photon spectrum through a DRM to get a count spectrum

Parameters
  • function (<function>) – A photon spectrum function. The function must accept a list of function parameters as its first argument and an array of photon energies as its second argument. The function must return the evaluation of the photon model in units of ph/s-cm^2-keV.

  • params (list of float) – A list of parameter values to be passed to the photon spectrum function

  • index (int, optional) – The index of the DRM to be used from the list of DRMs in the response. Default is 0.

  • atime (float, optional) – A time associated with a DRM in the list of available DRMs. The DRM covering the time closest to atime will be used. If set, this overrides the index argument.

  • channel_mask (np.array, optional) – A boolean mask where True indicates the channel is to be used for folding and False indicates the channel is to not be used for folding. If omitted, all channels are used.

Returns

np.array – The count spectrum

classmethod from_arrays(chan_lo, chan_hi, energies_lo, energies_hi, eff_area, trigtime=None, tstart=None, tstop=None, detnam=None, object=None, ra_obj=None, dec_obj=None, err_rad=None, infile=None, ch2e_ver=None, gain=None, mat_type=None, src_az=None, src_el=None, geo_az=None, geo_el=None, det_ang=None, geo_ang=None, atscat=None, datatype=None, filetype='DRM', filename=None)[source]

Create a single-DRM RSP object from arrays of channel energies, incident photon energies and effective area.

Parameters
  • chan_lo (np.array) – The low edges of the energy channels

  • chan_hi (np.array) – The high edges of the energy channels

  • energies_lo (np.array) – The low edges of the photon bins

  • energies_hi (np.array) – The high edges of the photon bins

  • eff_area (np.array) – The effective area of each photon bin mapped to each energy channel

  • trigtime (float, optional) – The trigger time, if applicable.

  • tstart (float, optional) – The start time of the observation for which the response is valid

  • tstop (float, optional) – The stop time of the observation for which the response is valid

  • detnam (str, optional) – The detector name the corresponds to the response

  • object (str, optional) – The object being observed

  • ra_obj (float, optional) – The RA of the object

  • dec_obj (float, optional) – The Dec of the object

  • err_rad (float, optional) – The localization error radius of the object

  • infile (str, optional) – The file from which the response was made

  • ch2e_ver (str, optional) – The version of the channel to energy calibration

  • gain (float, optional) – The gain, if applicable

  • mat_type (int, optional) – If set: 0=instrument scattering only, 1=atmospheric scattering only, 2=instrument+atmospheric scattering

  • src_az (float, optional) – The azimuth in spacecraft coordinates for which the response is valid

  • src_el (float, optional) – The elevation in spacecraft coordinates for which the response is valid

  • geo_az (float, optional) – The location of the geocenter in spacecraft azimuth

  • geo_el (float, optional) – The location of the geocenter in spacecraft elevation

  • det_ang (float, optional) – The angle between detector normal and position for which the response is valid

  • geo_ang (float, optional) – The angle between the geocenter and the position for which the response is valid

  • atscat (str, optional) – The atmospheric scattering database used

  • datatype (str, optional) – The datatype to which the response corresponds

  • filetype (str, optional) – If not set, default is ‘DRM’

  • filename (str, optional) – If not set, a default name will be constructed from other inputs

Returns

RSP – The newly created RSP object

interpolate(atime, **kwargs)[source]

Interpolate over a multi-DRM object for a given time and return a new single-DRM RSP object. This function uses scipy.interpolate.interp1d.

If the given time is before (after) the times covered by the DRMs within the RSP object, then the first (last) DRM will be returned.

Parameters
  • atime (float) – The time corresponding to a DRM. The nearest DRM to the time will be chosen

  • **kwargs – Keyword arguments to be passed to scipy.interpolate.interp1d

Returns

RSP – The interpolated single-DRM RSP object

nearest_drm(atime)[source]

Return the nearest DRM to the requested time

Parameters

atime (float) – The requested time

Returns

np.array – The nearest DRM

classmethod open(filename)[source]

Read a RSP or RSP2 file from disk

Parameters

filename (str) – The filename

Returns

RSP The RSP object

photon_effective_area(index=None, atime=None)[source]

Returns the effective area as a function of incident photon energy (integrated over recorded energy channels).

Parameters
  • index (int, optional) – The DRM index

  • atime (float, optional) – The time corresponding to a DRM

Returns

Bins – A Bins object containing the effective area histogram

rebin(index=0, factor=None, edge_indices=None)[source]

Rebins the channel energy axis of a DRM and returns a new response object. Rebinning can only be used to downgrade the channel resolution and is constrained to the channel edges of the current DRM.

Rebinning can be performed by either defining an integral factor of the number of current energy channels to combine (e.g. factor=4 for 128 energy channels would result in 32 energy channels), or by defining an index array into the current channel edges that will define the new set of edges.

Parameters
  • index (int, optional) – The DRM index. Default is 0

  • factor (int, optional) – The rebinning factor. Must set either this or edge_indices

  • edge_indices (np.array, optional) – The index array that represents which energy edges should remain in the rebinned DRM.

Returns

(RSP) – A new rebinned single-DRM response object

resample(index=0, num_photon_bins=None, photon_bin_edges=None, num_interp_points=20, interp_kind='linear')[source]

Resamples the incident photon axis of a DRM and returns a new response object. Resampling can be used to downgrade the photon energy resolution, upgrade the resolution, and/or redefine the edges of the incident photon bins. By definition, the resampling can only be performed within the photon energy range of the current object.

The process for resampling is to create a high-resolution grid for each new photon bin, interpolate the differential effective area onto that grid (interpolation is performed on log10(effective area)), and then integrate the differential effective area over the grid points in each bin. The final effective area is then scaled by the ratio of the initial number of photon bins to the requested number of photon bins.

Parameters
  • index (int, optional) – The DRM index. Default is 0

  • num_photon_bins (int, optional) – The number of photon bins in the new DRM. The bin edges will be generated logarithmically. Only set this or photon_bin_edges.

  • photon_bin_edges (np.array, optional) – The array of photon bin edges. Only set this or num_photon_bins

  • num_interp_points (int, optional) – The number of interpolation points used to integrate over for each new photon bin. Default is 20.

  • interp_kind (str, optional) – The kind of interpolation to be passed to scipy.interp1d. Default is ‘linear’.

Returns

(RSP) – A new resampled single-DRM response object

set_properties(detector=None, trigtime=None, tstart=None, datatype=None, extension=None, meta=None)

Set the properties of the data file. Useful for creating new data files. A standardized filename will be built from the parameters

Parameters
  • detector (str, optional) – The detector the data file belongs to

  • trigtime (float, optional) – The trigger time, if applicable

  • tstart (float) – The start time of the data file. Must be set if trigtime is not.

  • datatype (str, optional) – The type of data the file contains

  • extension (str, optional) – The extension of the data file

  • meta (str, optional) – A metadata atttribute to be added to the filename

weighted(weights, bin_centers)[source]

Return a counts-weighted DRM. For long spectral accumulations, the time series of responses needs to be weighted by the counts history

Parameters
  • weights (np.array) – The weights of the counts history. Must sum to 1.0

  • bin_centers (np.array) – The times at the center of the time bins

Returns

np.array – The weighted DRM

write(directory, filename=None, **kwargs)[source]

Writes a single-DRM RSP object to a FITS file.

Parameters
  • directory (str) – Where to write the file

  • filename (str, optional) – The filename to write to

write_weighted(tstart, tstop, weights, directory, filename=None)[source]

Writes a counts-weighted DRM to a FITS file.

Parameters
  • tstart (np.array) – The start times of the count history

  • tstop (np.array) – The end times of the count history

  • weights (np.array) – The normalized weights

  • directory (str) – The directory to write the file

  • filename (str, optional) – If set, will override the standardized name


Tcat

class gbm.data.Tcat[source]

Bases: gbm.data.data.DataFile

Class for Trigger Catalog (TCAT) files

Attributes
  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • fermi_location (float, float) – Fermi’s orbital longitude and latitude

  • filename (str) – The filename

  • full_path (str) – The full path+filename

  • headers (dict) – The headers for each extension of the file

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

  • localizing_instrument (str) – The localizing instrument

  • location (float, float, float) – RA, Dec, and localization uncertainty

  • location_fermi_frame (float, float) – Location in Fermi azimuth and zenith

  • name (str) – Name of the trigger

  • time_range (float, float) – The time range

  • trigtime (float) – The trigger time

classmethod open(filename)[source]

Open a TCAT file and return the Tcat object

Parameters

filename (str) – The filename of the TCAT file

Returns

Tcat – The Tcat object


Trigdat

class gbm.data.Trigdat[source]

Bases: gbm.data.poshist.PosHist

Class for the GBM Trigger Data.

Attributes
  • backrates (BackRates) – A BackRates object containing the info from the on-board background

  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • filename (str) – The filename

  • fsw_locations – (FswLocation): A list of flight-software-determined locations for the event

  • full_path (str) – The full path+filename

  • headers (dict) – The headers for each extension of the file

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

  • maxrates (list of MaxRates) – A list of MaxRates objects, each containing maxrates info

  • num_maxrates (int) – The number of MaxRates issued by the flight software

  • time_range (float, float) – The time range of the data

  • triggered_detectors – (list of str): The detectors that were triggered

  • trigrates (MaxRates) – A MaxRates object containing the trigger information and rates

  • trigtime (float) – The trigger time

get_fsw_locations(index)[source]

Retrieve a flight software localization

Parameters

index (int) – The index of the localization to retrieve. Not to exceed num_maxrates-1

Returns

FswLocation – The flight software localization info

get_maxrates(index)[source]

Retrieve a MaxRates

Parameters

index (int) – The index of the MaxRates to retrieve. Not to exceed num_maxrates-1

Returns

MaxRates – The MaxRates object

get_saa_passage(times)[source]

Determine if Fermi is in an SAA passage

Parameters

times (float or np.array) – Time(s) in MET

Returns

np.array(dtype=bool) – A boolean array, where True indicates Fermi is in SAA

classmethod open(filename)[source]

Open and read a trigdat file

Parameters

filename (str) – The filename of the trigdat file

Returns

Trigdat – The Trigdat object

sum_detectors(detectors, timescale=1024)[source]

Sum the data from a list of detectors and convert to a CTIME-like PHAII object

Parameters
  • detectors (list of str) – The detectors to sum

  • timescale (int, optional) – The minimum timescale in ms of the data to return. Available options are 1024, 256, and 64.

Returns

Ctime – The CTIME-like PHAII object with the detector-summed data

to_ctime(detector, timescale=1024)[source]

Convert the data for a detector to a CTIME-like PHAII object

Parameters
  • detector (str) – The detector to convert

  • timescale (int, optional) – The minimum timescale in ms of the data to return. Available options are 1024, 256, and 64.

Returns

Ctime – The CTIME-like PHAII object with the trigdat data


TTE

class gbm.data.TTE[source]

Bases: gbm.data.data.DataFile

Class for Time-Tagged Event data

Attributes
  • data (EventList) – The event data

  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • energy_range (float, float) – The energy range of the data

  • filename (str) – The filename

  • full_path (str) – The full path+filename

  • gti ([(float, float), …]) – A list of tuples defining the good time intervals

  • headers (dict) – The headers for each extension of the TTE

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

  • numchans (int) – The number of energy channels

  • trigtime (float) – The trigger time of the data, if available.

  • time_range (float, float) – The time range of the data

classmethod from_data(data, gti=None, trigtime=0.0, detector=None, object=None, ra_obj=None, dec_obj=None, err_rad=None)[source]

Create a TTE object from an EventList data object.

Parameters
  • data (EventList) – The event data

  • gti ([(float, float), ...], optional) – The list of tuples representing the good time intervals (start, stop). If omitted, the GTI is assumed to be [(tstart, tstop)].

  • trigtime (float, optional) – The trigger time, if applicable. If provided, the data times will be shifted relative to the trigger time.

  • detector (str, optional) – The detector that produced the data

  • object (str, optional) – The object being observed

  • ra_obj (float, optional) – The RA of the object

  • dec_obj (float, optional) – The Dec of the object

  • err_rad (float, optional) – The localization error radius of the object

Returns

TTE – The newly created TTE object

get_exposure(time_ranges=None)[source]

Calculate the total exposure of a time range or time ranges of data

Parameters

time_ranges ([(float, float), ...], optional) – The time range or time ranges over which to calculate the exposure. If omitted, calculates the total exposure of the data.

Returns

float – The exposure of the time selections

classmethod merge(ttes, primary=0, force_unique=True)[source]

Merge multiple TTE objects into a single new TTE object.

Note

The amount of data in a single TTE object can be quite large. It is up to you to determine if you have enough memory to support the merge.

Parameters
  • ttes (lists) – A list of TTE objects to merge

  • primary (int) – The index into the list of TTE objects to designate the primary TTE object. The primary object will be the reference for the header information for the new merged TTE object. Default is the first TTE object in the list (primary=0).

  • force_unique (bool, optional) – If True, force all events to be unique via brute force sorting. If False, the EventLists will only be checked and masked for overlapping time ranges. Events can potentially be lost if the merged EventLists contain overlapping times (but not necessarily duplicate events), however this method is much faster. Default is True.

Returns

TTE – The new merged TTE object

classmethod open(filename)[source]

Open a TTE FITS file and return the TTE object

Parameters

filename (str) – The filename of the FITS file

Returns

TTE – The TTE object

rebin_energy(bin_method, *args)[source]

Produce a TTE object with rebinned PHA channels

Parameters
  • bin_method (<function>) – A rebinning function

  • *args – Arguments to be passed to the bin method

Returns

TTE: A new TTE object with the rebinned PHA channels

set_properties(detector=None, trigtime=None, tstart=None, datatype=None, extension=None, meta=None)

Set the properties of the data file. Useful for creating new data files. A standardized filename will be built from the parameters

Parameters
  • detector (str, optional) – The detector the data file belongs to

  • trigtime (float, optional) – The trigger time, if applicable

  • tstart (float) – The start time of the data file. Must be set if trigtime is not.

  • datatype (str, optional) – The type of data the file contains

  • extension (str, optional) – The extension of the data file

  • meta (str, optional) – A metadata atttribute to be added to the filename

slice_energy(energy_ranges)[source]

Slice the TTE by one or more energy range. Produces a new TTE object.

Parameters

energy_ranges ([(float, float), ...]) – The energy ranges to slice the data to.

Returns

TTE – A new TTE object with the requested energy ranges

slice_time(time_ranges)[source]

Slice the TTE by one or more time range. Produces a new TTE object.

Parameters

time_ranges ([(float, float), ...]) – The time ranges to slice the data to.

Returns

TTE – A new TTE object with the requested time ranges

to_pha(time_ranges=None, energy_range=None, channel_range=None, **kwargs)[source]

Integrate the TTE data over time to produce a PHA object

Parameters
  • time_ranges ([(float, float), ...], optional) – The time range of the spectrum. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the spectrum. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the spectrum. If omitted, uses the entire energy range of the data.

  • **kwargs – Options passed to PHA.from_data()

Returns

PHA – The PHA object

to_phaii(bin_method, *args, time_range=None, energy_range=None, channel_range=None, **kwargs)[source]

Utilizing a binning function, convert the data to a PHAII object

Parameters
  • bin_method (<function>) – A binning function

  • *args – Arguments to pass to the binning function

  • time_range ([(float, float), ...], optional) – The time range of the spectrum. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the spectrum. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the spectrum. If omitted, uses the entire energy range of the data.

  • **kwargs – Options to pass to the binning function

Returns

Cspec – The PHAII object

to_spectrum(time_range=None, energy_range=None, channel_range=None)[source]

Integrate the TTE data over time to produce a count spectrum

Parameters
  • time_range ([(float, float), ...], optional) – The time range of the spectrum. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the spectrum. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the spectrum. If omitted, uses the entire energy range of the data.

Returns

EnergyBins – The spectrum data

write(directory, filename=None)[source]

Writes the data to a FITS file.

Parameters
  • directory (str) – The directory to write the file

  • filename (str, optional) – If set, will override the standardized name


Scat

class gbm.data.Scat[source]

Bases: gbm.data.data.DataFile

A container class for the spectral fit data in an SCAT file

Attributes
  • detectors (list) – The DetectorData objects used in the analysis

  • headers (dict) – The SCAT file headers

  • model_fits (list) – The GbmModelFit objects, one for each model fit

  • num_detectors (int) – The number of detectors in the SCAT file

  • num_fits (int) – The number of model fits

add_detector_data(detector_data)[source]

Add a new detector to the Scat

Parameters

detector_data (DetectorData) – The detector data

add_model_fit(model_fit)[source]

Add a new model fit to the Scat

Parameters

model_fit (GbmModelFit) – The model fit data

classmethod open(filename)[source]

Open a SCAT FITS file and create a Scat object

Parameters

filename (str) – The file to open

Returns

(Scat)

set_properties(detector=None, trigtime=None, tstart=None, datatype=None, extension=None, meta=None)

Set the properties of the data file. Useful for creating new data files. A standardized filename will be built from the parameters

Parameters
  • detector (str, optional) – The detector the data file belongs to

  • trigtime (float, optional) – The trigger time, if applicable

  • tstart (float) – The start time of the data file. Must be set if trigtime is not.

  • datatype (str, optional) – The type of data the file contains

  • extension (str, optional) – The extension of the data file

  • meta (str, optional) – A metadata atttribute to be added to the filename

write(directory, filename=None)[source]

Write a Scat object to a FITS file

Parameters
  • directory (str) – The directory where the file is to be written

  • filename (str, optional) – The filename. If not set, a default filename will be used.


Base Data Classes

Some base data classes. These are not intended to be instantiated, but can be useful to inherit if designing new data classes.

DataFile

class gbm.data.data.DataFile[source]

Base class for an interface to data files

Note

This class should not be used directly, instead use one that inherits from it and is specific to your data type e.g. TTE

Attributes
  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • filename (str) – The filename

  • full_path (str) – The full path+filename

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

set_properties(detector=None, trigtime=None, tstart=None, datatype=None, extension=None, meta=None)[source]

Set the properties of the data file. Useful for creating new data files. A standardized filename will be built from the parameters

Parameters
  • detector (str, optional) – The detector the data file belongs to

  • trigtime (float, optional) – The trigger time, if applicable

  • tstart (float) – The start time of the data file. Must be set if trigtime is not.

  • datatype (str, optional) – The type of data the file contains

  • extension (str, optional) – The extension of the data file

  • meta (str, optional) – A metadata atttribute to be added to the filename


HealPix

class gbm.data.HealPix[source]

Bases: gbm.data.data.DataFile

Base class for HEALPix localization files.

Attributes
  • centroid (float, float) – The RA and Dec of the highest probability pixel

  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • filename (str) – The filename

  • full_path (str) – The full path+filename

  • headers (dict) – The headers for each extension

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

  • npix (int) – Number of pixels in the HEALPix map

  • nside (int) – The HEALPix resolution

  • pixel_area (float) – The area of each pixel in square degrees

  • trigtime (float) – The time corresponding to the localization

area(clevel)[source]

Calculate the sky area contained within a given confidence region

Parameters

clevel (float) – The localization confidence level (valid range 0-1)

Returns

float – The area contained in square degrees

confidence(ra, dec)[source]

Calculate the localization confidence level for a given point. This function interpolates the map at the requested point rather than providing the value at the nearest pixel center.

Parameters
  • ra (float) – The RA

  • dec (float) – The Dec

Returns

float – The localization confidence level

confidence_region_path(clevel, numpts_ra=360, numpts_dec=180)[source]

Return the bounding path for a given confidence region

Parameters
  • clevel (float) – The localization confidence level (valid range 0-1)

  • numpts_ra (int, optional) – The number of grid points along the RA axis. Default is 360.

  • numpts_dec (int, optional) – The number of grid points along the Dec axis. Default is 180.

Returns

[(np.array, np.array), …] – A list of RA, Dec points, where each item in the list is a continuous closed path.

convolve(model, *args)[source]

Convolve the map with a model kernel. The model can be a Gaussian kernel or any mixture of Gaussian kernels. Uses healpy.smoothing.

An example of a model kernel with a 50%/50% mixture of two Gaussians, one with a 1-deg width, and the other with a 3-deg width:

def gauss_mix_example():
    sigma1 = np.deg2rad(1.0)
    sigma2 = np.deg2rad(3.0)
    frac1 = 0.50
    return ([sigma1, sigma2], [frac1])
Parameters
  • model (<function>) – The function representing the model kernel

  • *args – Arguments to be passed to the model kernel function

Returns

HealPix – A new HealPix object that is a result of the convolution with the model kernel

classmethod from_annulus(center_ra, center_dec, radius, sigma, nside=None, **kwargs)[source]

Create a HealPix object of a Gaussian-width annulus

Parameters
  • center_ra (float) – The RA of the center of the annulus

  • center_dec (float) – The Dec of the center of the annulus

  • radius (float) – The radius of the annulus, in degrees, measured to the center of the of the annulus

  • sigma (float) – The Gaussian standard deviation width of the annulus, in degrees

  • nside (int, optional) – The nside of the HEALPix to make. By default, the nside is automatically determined by the sigma width. Set this argument to override the default.

  • **kwargs – Options to pass to from_data()

Returns

HealPix – The HEALPix annulus

classmethod from_data(prob_arr, sig_arr, trigtime=None)[source]

Create a HealPix object from healpix arrays

Parameters
  • prob_arr (np.array) – The HEALPix array containing the probability/pixel

  • sig_arr (np.array) – The HEALPix array containing the signficance

  • trigtime (float, optional) – The time corresponding to the localization

Returns

HealPix – The HEALPix localization

classmethod from_gaussian(center_ra, center_dec, sigma, nside=None, **kwargs)[source]

Create a HealPix object of a Gaussian

Parameters
  • center_ra (float) – The RA of the center of the Gaussian

  • center_dec (float) – The Dec of the center of the Gaussian

  • sigma (float) – The Gaussian standard deviation, in degrees

  • nside (int, optional) – The nside of the HEALPix to make. By default, the nside is automatically determined by the sigma of the Gaussian. Set this argument to override the default.

  • **kwargs – Options to pass to from_data()

Returns

HealPix – The HEALPix Gaussian

classmethod from_vertices(ra_pts, dec_pts, nside=64, **kwargs)[source]

Create a HealPix object from a list of RA, Dec vertices. The probability within the vertices will be distributed uniformly and zero probability outside the vertices.

Parameters
  • ra_pts (np.array) – The array of RA coordinates

  • dec_pts (np.array) – The array of Dec coordinates

  • nside (int, optional) – The nside of the HEALPix to make. Default is 64.

  • **kwargs – Options to pass to from_data()

Returns

HealPix – The HEALPix object

classmethod multiply(healpix1, healpix2, primary=1, output_nside=128)[source]

Multiply two HealPix maps and return a new HealPix object

Parameters
  • healpix1 (HealPix) – One of the HEALPix maps to multiply

  • healpix2 (HealPix) – The other HEALPix map to multiply

  • primary (int, optional) – If 1, use the first map header information, or if 2, use the second map header information. Default is 1.

  • output_nside (int, optional) – The nside of the multiplied map. Default is 128.

Returns

HealPix: The multiplied map

prob_array(numpts_ra=360, numpts_dec=180, sqdegrees=True, sig=False)[source]

Return the localization probability mapped to a grid on the sky

Parameters
  • numpts_ra (int, optional) – The number of grid points along the RA axis. Default is 360.

  • numpts_dec (int, optional) – The number of grid points along the Dec axis. Default is 180.

  • sqdegrees (bool, optional) – If True, the prob_array is in units of probability per square degrees, otherwise in units of probability per pixel. Default is True

  • sig (bool, optional) – Set True to retun the significance map on a grid instead of the probability. Default is False.

Returns

3-tuple containing

  • np.array: The probability (or significance) array with shape (numpts_dec, numpts_ra)

  • np.array: The RA grid points

  • np.array: The Dec grid points

probability(ra, dec, per_pixel=False)[source]

Calculate the localization probability at a given point. This function interpolates the map at the requested point rather than providing the vale at the nearest pixel center.

Parameters
  • ra (float) – The RA

  • dec (float) – The Dec

  • per_pixel (bool, optional) – If True, return probability per pixel, otherwise return probability per square degree. Default is False.

Returns

float – The localization probability

region_probability(healpix, prior=0.5)[source]

The probability that the HealPix localization is associated with another HealPix map. This is calculated against the null hypothesis that the two HealPix maps are unassociated:

\(P(A | \mathcal{I}) = \frac{P(\mathcal{I} | A) \ P(A)} {P(\mathcal{I} | A) \ P(A) + P(\mathcal{I} | \neg A) \ P(\neg A)}\)

where

  • \(P(\mathcal{I} | A)\) is the integral over the overlap of the two maps once the Earth occultation has been removed for this map.

  • \(P(\mathcal{I} | \neg A)\) is the integral over the overlap of this map with a uniform distribution on the sky (i.e. the probability the localization is associated with a random point on the sky)

  • \(P(A)\) is the prior probability that this localization is associated with the other HEALPix map.

Parameters
  • healpix (HealPix) – The healpix map for which to calculate the spatial association

  • prior (float, optional) – The prior probability that the localization is associated with the source. Default is 0.5

Returns

float – The probability that this HealPix localization is associated with the input HealPix map

set_properties(detector=None, trigtime=None, tstart=None, datatype=None, extension=None, meta=None)

Set the properties of the data file. Useful for creating new data files. A standardized filename will be built from the parameters

Parameters
  • detector (str, optional) – The detector the data file belongs to

  • trigtime (float, optional) – The trigger time, if applicable

  • tstart (float) – The start time of the data file. Must be set if trigtime is not.

  • datatype (str, optional) – The type of data the file contains

  • extension (str, optional) – The extension of the data file

  • meta (str, optional) – A metadata atttribute to be added to the filename

source_probability(ra, dec, prior=0.5)[source]

The probability that the HealPix localization is associated with a known point location. This is calculated against the null hypothesis that the HealPix localization originates from an unassociated random source that has equal probability of origination anywhere in the sky:

\(P(A | \mathcal{I}) = \frac{P(\mathcal{I} | A) \ P(A)} {P(\mathcal{I} | A) \ P(A) + P(\mathcal{I} | \neg A) \ P(\neg A)}\)

where

  • \(P(\mathcal{I} | A)\) is the probability of the localization at the point source once

  • \(P(\mathcal{I} | \neg A)\) is the probability per pixel assuming a uniform distribution on the sky (i.e. the probability the localization is associated with a random point on the sky)

  • \(P(A)\) is the prior probability that the localization is associated with the point source

Parameters
  • ra (float) – The RA of the known source location

  • dec (float) – The Dec of the known source location

  • prior (float, optional) – The prior probability that the localization is associated with the source. Default is 0.5

Returns

float – The probability that the HealPix localization is spatially associated with the point source


PHAII

class gbm.data.phaii.PHAII[source]

Bases: object

PHAII class for time history and spectra.

Attributes
  • data (TimeEnergyBins) – The PHAII data

  • gti ([(float, float), …]) – A list of tuples defining the good time intervals

  • headers (dict) – The headers for each extension of the PHAII

  • trigtime (float) – The trigger time of the data, if available.

  • time_range (float, float) – The time range of the data

  • energy_range (float, float) – The energy range of the data

  • numchans (int) – The number of energy channels

classmethod from_data(data, gti=None, trigtime=0.0, detector=None, object=None, ra_obj=None, dec_obj=None, err_rad=None)[source]

Create a PHAII object from TimeEnergyBins data object.

Parameters
  • data (TimeEnergyBins) – The PHAII data

  • gti ([(float, float), ...], optional) – The list of tuples representing the good time intervals (start, stop). If omitted, the GTI is assumed to be [(tstart, tstop)].

  • trigtime (float, optional) – The trigger time, if applicable. If provided, the data times will be shifted relative to the trigger time.

  • detector (str, optional) – The detector that produced the data

  • object (str, optional) – The object being observed

  • ra_obj (float, optional) – The RA of the object

  • dec_obj (float, optional) – The Dec of the object

  • err_rad (float, optional) – The localization error radius of the object

Returns

PHAII – The newly created PHAII object

get_exposure(time_ranges=None)[source]

Calculate the total exposure of a time range or time ranges of data

Parameters

time_ranges ([(float, float), ...], optional) – The time range or time ranges over which to calculate the exposure. If omitted, calculates the total exposure of the data.

Returns

float – The exposure of the time selections

classmethod open(filename)[source]

Open a PHAII FITS file and return the PHAII object

Parameters

filename (str) – The filename of the FITS file

Returns

PHAII – The PHAII object

rebin_energy(method, *args, energy_range=(None, None))[source]

Rebin the PHAII in energy given a rebinning method. Produces a new PHAII object.

Parameters
  • method (<function>) – The rebinning function

  • *args – Arguments to be passed to the rebinning function

  • energy_range ((float, float), optional) – The starting and ending energy to rebin. If omitted, uses the full range of data. Setting start or end to None will use the data from the beginning or end of the data, respectively.

Returns

PHAII: A new PHAII object with the requested rebinned data

rebin_time(method, *args, time_range=(None, None))[source]

Rebin the PHAII in time given a rebinning method. Produces a new PHAII object.

Parameters
  • method (<function>) – The rebinning function

  • *args – Arguments to be passed to the rebinning function

  • time_range ((float, float), optional) – The starting and ending time to rebin. If omitted, uses the full range of data. Setting start or end to None will use the data from the beginning or end of the data, respectively.

Returns

PHAII: A new PHAII object with the requested rebinned data

slice_energy(energy_ranges)[source]

Slice the PHAII by one or more energy range. Produces a new PHAII object.

Parameters

energy_ranges ([(float, float), ...]) – The energy ranges to slice the data to.

Returns

PHAII – A new PHAII object with the requested energy ranges

slice_time(time_ranges)[source]

Slice the PHAII by one or more time range. Produces a new PHAII object.

Parameters

time_ranges ([(float, float), ...]) – The time ranges to slice the data to.

Returns

PHAII – A new PHAII object with the requested time ranges

to_lightcurve(time_range=None, energy_range=None, channel_range=None)[source]

Integrate the PHAII data over energy to produce a lightcurve

Parameters
  • time_range ((float, float), optional) – The time range of the lightcurve. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the lightcurve. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the lightcurve. If omitted, uses the entire energy range of the data.

Returns

TimeBins – The lightcurve data

to_pha(time_ranges=None, energy_range=None, channel_range=None, **kwargs)[source]

Integrate the PHAII data over time to produce a PHA object

Parameters
  • time_ranges ([(float, float), ...], optional) – The time range of the spectrum. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the spectrum. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the spectrum. If omitted, uses the entire energy range of the data.

  • **kwargs – Options passed to gbm.data.PHA.from_data()

Returns

PHA – The single spectrum PHA

to_spectrum(time_range=None, energy_range=None, channel_range=None)[source]

Integrate the PHAII data over time to produce a count spectrum

Parameters
  • time_range ((float, float), optional) – The time range of the spectrum. If omitted, uses the entire time range of the data.

  • energy_range ((float, float), optional) – The energy range of the spectrum. If omitted, uses the entire energy range of the data.

  • channel_range ((int, int), optional) – The channel range of the spectrum. If omitted, uses the entire energy range of the data.

Returns

EnergyBins – The count spectrum data

write(directory, filename=None)[source]

Writes the data to a FITS file.

Parameters
  • directory (str) – The directory to write the file

  • filename (str, optional) – If set, will override the standardized name


Data Class Miscellania

Miscellaneous data classes: Collections of data, single-spectrum PHA and BAK classes, and Trigdat info classes

DataCollection

class gbm.data.DataCollection[source]

Bases: object

A container for a collection of like data objects, such as a collection of Ctime objects. This class exposes the individual objects’ attributes, and it exposes the methods of the object class so that methods can be called on the collection as a whole. For that reason, each object in the DataCollection must be of the same type, otherwise an error is raised. The type of the collection is set by the first object inserted into the collection and the collection type is immutable thereafter.

Objects are stored in the collection in the order they are added.

The number of items in the collection can be retrieved by len() and one can iterate over the items:

[data_item for data_item in DataCollection]

In addition to the DataCollection methods, all of the individual object attributes and methods are exposed, and they become methods of the DataCollection. Note that individual object attributes become methods i.e. if you have an item attribute called item.name, then the corresponding DataCollection method would be item.name().

Attributes
  • items (list) – The names of the items in the DataCollection

  • types (str) – The type of the objects in the DataCollection

classmethod from_list(data_list, names=None)[source]

Given a list of objects and optionally a list of corresponding names, create a new DataCollection.

Parameters
  • data_list (list of objects) – The list of objects to be in the collection

  • names (list of str, optional) – The list of corresponding names to the objects. If not set, will try to retrieve a name from object.filename (assuming it’s a data object). If that fails, each item will be named ambiguously ‘item1’, ‘item2’, etc.

Returns

DataCollection: The newly created collection

get_item(item_name)[source]

Retrieve an object from the DataCollection by name

Parameters

item_name (str) – The name of the item to retrieve

Returns

object – The retrieved data item

include(data_item, name=None)[source]

Insert an object into the collection. The first item inserted will set the immutable type.

Parameters
  • data_item (object) – A data object to include

  • name (str, optional) – An optional corresponding name. If not set, will try to retrieve a name from object.filename (assuming it’s a data object). If that fails, each item will be named ambiguously ‘item1’, ‘item2’, etc.

remove(item_name)[source]

Remove an object from the collection given the name

Parameters

item_name (str) – The name of the item to remove

to_list()[source]

Return the objects contained in the DataCollection as a list.

Returns

(list of objects) – The list of objects, in the order that they were inserted


GbmDetectorCollection

class gbm.data.GbmDetectorCollection[source]

Bases: gbm.data.collection.DataCollection

A container for a collection of GBM-specific data objects, such as a collection of Ctime objects from different detectors.

The special behavior of this class is to provide a way to interact with a collection of detector data that may contain a mix of different types of detectors. For example, many times we want a collection of GBM NaI and GBM BGO detectors. These detectors have very different energy ranges, and so may require different inputs for a variety of functions. This collection allows one to specify the different arguments for NaI and BGO data without having to implement many ugly and space-wasting loops and if...else decisions.

In addition to the GbmDetectorCollection methods, all of the individual object attributes and methods are exposed, and they become methods of the DataCollection. Note that individual object attributes become methods i.e. if you have an item attribute called item.name, then the corresponding DataCollection method would be item.name().

Attributes
  • items (list) – The names of the items in the DataCollection

  • types (str) – The type of the objects in the DataCollection

classmethod from_list(data_list, names=None, dets=None)[source]

Given a list of objects and optionally a list of corresponding names and corresponding detector names, create a new GbmDetectorCollection.

Parameters
  • data_list (list of objects) – The list of objects to be in the collection

  • names (list of str, optional) – The list of corresponding names to the objects. If not set, will try to retrieve a name from object.filename (assuming it’s a data object). If that fails, each item will be named ambiguously ‘item1’, ‘item2’, etc.

  • dets (list of str, optional) – The detector names for each object. If not set, will try to retrieve from the object.detector attribute. If that attribute doesn’t exist, an error will be raised, and the user will need to specify this list.

Returns

GbmDetectorCollection: The newly created collection

get_item(item_name)

Retrieve an object from the DataCollection by name

Parameters

item_name (str) – The name of the item to retrieve

Returns

object – The retrieved data item

include(data_item, det, name=None)[source]

Insert an object into the GbmDetectorCollection. The first item inserted will set the immutable type.

Parameters
  • data_item (object) – A data object to include

  • det (str) – The corresponding detector for the item

  • name (str, optional) – An optional corresponding name. If not set, will try to retrieve a name from object.filename (assuming it’s a data object). If that fails, each item will be named ambiguously ‘item1’, ‘item2’, etc.

remove(item_name)[source]

Remove an object from the collection given the name

Parameters

item_name (str) – The name of the item to remove

to_list()

Return the objects contained in the DataCollection as a list.

Returns

(list of objects) – The list of objects, in the order that they were inserted


BAK

class gbm.data.BAK[source]

Bases: gbm.data.pha.PHA

Class for a PHA background spectrum.

Attributes
  • data (EnergyBins) – The PHA data

  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • energy_range (float, float) – The energy range of the spectrum

  • exposure (float) – The exposure of the PHA data

  • filename (str) – The filename

  • full_path (str) – The full path+filename

  • gti ([(float, float), …]) – The good time intervals

  • headers (dict) – The headers for each extension of the PHA

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

  • numchans (int) – The number of energy channels

  • tcent (float) – The center time of the data

  • time_range (float, float) – The time range of the spectrum

  • trigtime (float) – The trigger time of the data, if available

classmethod from_data(data, tstart, tstop, gti=None, trigtime=0.0, detector=None, object=None, ra_obj=None, dec_obj=None, err_rad=None, datatype=None, meta=None)[source]

Create a BAK object from a BackgroundSpectrum data object.

Parameters
  • data (BackgroundSpectrum) – The background count rate spectrum

  • tstart (float) – The start time of the data

  • tstop (float) – The end time of the data

  • gti ([(float, float), ...], optional) – The list of tuples representing the good time intervals (start, stop). If omitted, the GTI is assumed to be [(tstart, tstop)].

  • trigtime (float, optional) – The trigger time, if applicable. If provided, the data times will be shifted relative to the trigger time. Default is zero.

  • detector (str, optional) – The detector that produced the data

  • object (str, optional) – The object being observed

  • ra_obj (float, optional) – The RA of the object

  • dec_obj (float, optional) – The Dec of the object

  • err_rad (float, optional) – The localization error radius of the object

  • datatype (str, optional) – The datatype from which the PHA is created

  • meta (str, optional) – Additional metadata string to be added to standard filename

Returns

BAK – The BAK object

classmethod open(filename)[source]

Open a BAK FITS file and return the BAK object

Parameters

filename (str) – The filename of the FITS file

Returns

BAK – The BAK object

rebin_energy(method, *args, emin=None, emax=None)[source]

Not Implemented

set_properties(detector=None, trigtime=None, tstart=None, datatype=None, extension=None, meta=None)

Set the properties of the data file. Useful for creating new data files. A standardized filename will be built from the parameters

Parameters
  • detector (str, optional) – The detector the data file belongs to

  • trigtime (float, optional) – The trigger time, if applicable

  • tstart (float) – The start time of the data file. Must be set if trigtime is not.

  • datatype (str, optional) – The type of data the file contains

  • extension (str, optional) – The extension of the data file

  • meta (str, optional) – A metadata atttribute to be added to the filename

slice_energy(method, *args, emin=None, emax=None)[source]

Not Implemented

write(directory, filename=None)[source]

Writes the data to a FITS file.

Parameters
  • directory (str) – The directory to write the file

  • filename (str, optional) – If set, will override the standardized name


PHA

class gbm.data.PHA[source]

Bases: gbm.data.data.DataFile

PHA class for count spectra.

Attributes
  • data (EnergyBins) – The PHA data

  • datatype (str) – The datatype of the file

  • detector (str) – The GBM detector the file is associated with

  • directory (str) – The directory the file is located in

  • energy_range (float, float) – The energy range of the spectrum

  • exposure (float) – The exposure of the PHA data

  • filename (str) – The filename

  • full_path (str) – The full path+filename

  • gti ([(float, float), …]) – The good time intervals

  • headers (dict) – The headers for each extension of the PHA

  • id (str) – The GBM file ID

  • is_gbm_file (bool) – True if the file is a valid GBM standard file, False if it is not.

  • is_trigger (bool) – True if the file is a GBM trigger file, False if not

  • numchans (int) – The number of energy channels

  • tcent (float) – The center time of the data

  • time_range (float, float) – The time range of the spectrum

  • trigtime (float) – The trigger time of the data, if available

  • valid_channels (np.array(dtype=bool)) – The channels that are marked valid and available for fitting

classmethod from_data(data, tstart, tstop, gti=None, trigtime=0.0, detector=None, object=None, ra_obj=None, dec_obj=None, err_rad=None, datatype=None, backfile=None, rspfile=None, meta=None, channel_mask=None)[source]

Create a PHA object from an EnergyBins data object.

Parameters
  • data (EnergyBins) – The PHA count spectrum data

  • tstart (float) – The start time of the data

  • tstop (float) – The end time of the data

  • gti ([(float, float), ...], optional) – The list of tuples representing the good time intervals (start, stop). If omitted, the GTI is assumed to be [(tstart, tstop)].

  • trigtime (float, optional) – The trigger time, if applicable. If provided, the data times will be shifted relative to the trigger time. Default is zero.

  • detector (str, optional) – The detector that produced the data

  • object (str, optional) – The object being observed

  • ra_obj (float, optional) – The RA of the object

  • dec_obj (float, optional) – The Dec of the object

  • err_rad (float, optional) – The localization error radius of the object

  • datatype (str, optional) – The datatype from which the PHA is created

  • backfile (str, optional) – The associated background file

  • rspfile (str, optional) – The associated response file

  • meta (str, optional) – Additional metadata string to be added to standard filename

  • channel_mask (np.array(dtype=bool)) – A boolean array representing the valid channels to be used for fitting. If omitted, assumes all non-zero count channels are valid.

Returns

PHA – The PHA object

classmethod open(filename)[source]

Open a PHA FITS file and return the PHA object

Parameters

filename (str) – The filename of the FITS file

Returns

PHA – The PHA object

rebin_energy(method, *args, emin=None, emax=None)[source]

Rebin the PHA in energy given a rebinning method. Produces a new PHA object.

Parameters
  • method (<function>) – The rebinning function

  • *args – Arguments to be passed to the rebinning function

  • emin (float, optional) – The minimum energy to rebin. If omitted, uses the lowest energy edge of the data

  • emax (float, optional) – The maximum energy to rebin. If omitted, uses the highest energy edge of the data

Returns

PHA – A new PHA object with the requested rebinned data

set_properties(detector=None, trigtime=None, tstart=None, datatype=None, extension=None, meta=None)

Set the properties of the data file. Useful for creating new data files. A standardized filename will be built from the parameters

Parameters
  • detector (str, optional) – The detector the data file belongs to

  • trigtime (float, optional) – The trigger time, if applicable

  • tstart (float) – The start time of the data file. Must be set if trigtime is not.

  • datatype (str, optional) – The type of data the file contains

  • extension (str, optional) – The extension of the data file

  • meta (str, optional) – A metadata atttribute to be added to the filename

slice_energy(energy_ranges)[source]

Slice the PHA by one or more energy range. Produces a new PHA object.

Parameters

energy_ranges ([(float, float), ...]) – The energy ranges to slice over

Returns

PHA – A new PHA object with the requested energy range

write(directory, filename=None, backfile=None)[source]

Writes the data to a FITS file.

Parameters
  • directory (str) – The directory to write the file

  • filename (str, optional) – If set, will override the standardized name

  • backfile (str, optional) – If set, will add the associated BAK file name to the PHA SPECTRUM header


scat.Parameter

class gbm.data.scat.Parameter(value, uncert, name='', units=None, support=(- inf, inf))[source]

Bases: object

A fit parameter class

Parameters
  • value (float) – The central fit value

  • uncert (float or 2-tuple) – The 1-sigma uncertainty. If a 2-tuple, then is of the form (low, high)

  • name (str, optional) – The name of the parameter

  • units (str, optional) – The units of the parameter

  • support (2-tuple, optional) – The valid support of the parameter

Attributes
  • name (str) – The name of the parameter

  • support (2-tuple) – The valid support of the parameter

  • uncertainty (2-tuple) – The 1-sigma uncertainty

  • units (str) – The units of the parameter

  • value (float) – The central fit value

one_sigma_range()[source]

Return the 1 sigma range of the parameter fit

to_fits_value()[source]

Return as a tuple to be used for a FITS file

Returns

(tuple)

2-value tuple (value, uncertainty) or 3-value tuple

(value, +uncertainty, -uncertainty)

valid_value()[source]

Check if the parameter value is within the allowed parameter range


scat.GbmModelFit

class gbm.data.scat.GbmModelFit(name, time_range, **kwargs)[source]

Bases: gbm.data.scat.ModelFit

A container for the info from a model fit, with values used in the GBM SCAT files. Inherits from ModelFit.

Attributes
  • covariance (np.array) – The covariance matrix of the fit

  • dof (int) – The degrees-of-freedom of the fit

  • duration_fluence (EnergyFluence) – The energy fluence over the duration energy range, nominally 50-300 keV

  • energy_fluence (EnergyFluence) – The energy fluence, nominally over 10-1000 keV

  • energy_fluence_50_300 (EnergyFluence) – The energy fluence over 50-300 keV

  • energy_flux (EnergyFlux) – The energy flux, nominally over 10-1000 keV

  • flux_energy_range (tuple) – The energy range of the flux and fluence, (low, high)

  • name (str) – The name of the model

  • parameters (list) – A list of model parameters

  • photon_fluence (PhotonFluence) – The photon fluence, nominally over 10-1000 keV

  • photon_flux (PhotonFlux) – The photon flux, nominally over 10-1000 keV

  • photon_flux_50_300 (PhotonFlux) – The photon flux over 50-300 keV

  • stat_name (str) – The name of the fit statistic

  • stat_value (float) – The fit statistic value

  • time_range (float, float) – The time range of the model fit, (low, high)

classmethod from_fits_row(fits_row, model_name, param_names=None, flux_range=(10.0, 1000.0), dur_range=(50.0, 300.0))[source]

Read a FITS row and return a GbmModelFit object

Returns

(GbmModelFit)

parameter_list()

Return the list of parameter names

Returns

(list) – The parameter names

to_fits_row()[source]

Return the contained data as a FITS table row

Returns

(astropy.io.fits.BinTableHDU) – The FITS table


scat.DetectorData

class gbm.data.scat.DetectorData(instrument, detector, datatype, filename, numchans, **kwargs)[source]

Bases: object

A container for detector info used in a fit

Parameters
  • instrument (str) – The name of the instrument

  • detector (str) – The name of the detector

  • datatype (str) – The name of the datatype

  • filename (str) – The filename of the data file

  • numchans (int) – Number of energy channels used

  • active (bool, optional) – True if the detector is used in the fit

  • response (str, optional) – The filename of the detector response

  • time_range (tuple, optional) – The time range of the data used

  • energy_range (tuple, optional) – The energy range of the data used

  • channel_range (tuple, optional) – The energy channel range of the data

  • energy_edges (np.array, optional) – The edges of the energy channels

  • photon_counts (np.array, optional) – The deconvolved photon counts for the detector

  • photon_model (np.array, optional) – The photon model for the detector

  • photon_errors (np.array, optional) – The deconvolved photon count errors for the detector

Attributes
  • active (bool, optional) – True if the detector is used in the fit

  • channel_range = (int, int) – The energy channel range of the data

  • datatype (str) – The name of the datatype

  • detector (str) – The name of the detector

  • energy_edges (np.array) – The edges of the energy channels

  • energy_range (float, float) – The energy range of the data used

  • filename (str) – The filename of the data file

  • instrument (str) – The name of the instrument

  • numchans (int) – Number of energy channels used

  • photon_counts (np.array) – The deconvolved photon counts for the detector

  • photon_errors (np.array) – The deconvolved photon count errors for the detector

  • photon_model (np.array) – The photon model for the detector

  • response (str) – The filename of the detector response

  • time_range (float, float) – The time range of the data used

classmethod from_fits_row(fits_row)[source]

Read a FITS row and return a DetectorData object

Returns

(DetectorData)

to_fits_row()[source]

Return the contained data as a FITS table row

Returns

(astropy.io.fits.BinTableHDU) – The FITS row


trigdat.MaxRates

class gbm.data.trigdat.MaxRates(rec_array)[source]

Bases: object

Class for the MaxRates data in Trigdat.

Parameters

rec_array (np.recarray) – The FITS TRIGRATE or MAXRATES record array from the trigdat file

Attributes
  • all_rates (np.array) – An array (numchans, numdets) of the maxrates

  • eic (np.array) – The position of Fermi in Earth inertial coordinates

  • numchans (int) – The number of energy channels

  • numdets (int) – The number of detectors

  • quaternions (np.array) – The quaternions at the maxrates time

  • timescale (float) – The timescale of the maxrates

  • time_range (float, float) – The time range of the maxrates

get_detector(det)[source]

Retrieve the maxrates for a detector

Parameters

det (str) – The detector

Returns

np.array – An array of size (numchans,) of rates for the detector


trigdat.BackRates

class gbm.data.trigdat.BackRates(rec_array)[source]

Bases: object

Class for the background rates data in Trigdat.

Parameters

rec_array (np.recarray) – The FITS BCKRATES record array from the trigdat file

Attributes
  • all_rates (np.array) – An array (numchans, numdets) of the maxrates

  • numchans (int) – The number of energy channels

  • numdets (int) – The number of detectors

  • quality (int, int) – The quality flags for the background

  • time_range (float, float) – The time range of the maxrates

get_detector(det)[source]

Retrieve the background rates for a detector

Parameters

det (str) – The detector

Returns

np.array – An array of size (numchans,) of background rates for the detector


trigdat.FswLocation

class gbm.data.trigdat.FswLocation(rec_array)[source]

Bases: object

Class for the flight software localization info

Parameters

rec_array (np.recarray) – The FITS OB_CALC record array from the trigdat file

Attributes
  • fluence (float) – The fluence of the localization interval

  • hardness_ratio (float) – The hardness ratio for the localization interval

  • intensity (float) – The brightness of the signal

  • location (float, float, float) – The RA, Dec, and statistical error of the onboard localization

  • location_sc (float) – The localization in spacecraft coordinates: Azimuth, Zenith

  • next_classification (str, float) – The next most likely classification of the trigger and the probability

  • significance (float) – The S/N ratio of the localization interval

  • spectrum (str) – The spectrum used in the localization

  • time (float) – Time at which the localization was calculated

  • timescale (float) – The localization interval timescale

  • top_classification (str, float) – The most likely classification of the trigger and the probability


Data Primitives

Primitive data classes. The GBM Primary classes act as wrappers around some of these classes with additional specialized metadata. These classes may be useful in designing new higher-level data classes.

Bins

class gbm.data.primitives.Bins(counts, lo_edges, hi_edges)[source]

Bases: object

A primitive class defining a set of histogram bins

Parameters
  • counts (np.array) – The array of counts in each bin

  • lo_edges (np.array) – The low-value edges of the bins

  • hi_edges (np.array) – The high-value edges of the bins

Attributes
  • centroids (np.array) – The centroids of the bins

  • counts (np.array) – The counts in each bin

  • count_uncertainty (np.array) – The count uncertainty in each bin

  • hi_edges (np.array) – The high-value edges of the bins

  • lo_edges (np.array) – The low-value edges of the bins

  • range (float, float) – The range of the bin edges

  • rates (np.array) – The count rate of each bin: counts/width

  • rate_uncertainty (np.array) – The count rate uncertainty of each bin

  • size (int) – Number of bins

  • widths (np.array) – The widths of the bins

closest_edge(val, which='either')[source]

Return the closest bin edge

Parameters
  • val (float) – Input value

  • which (str, optional) –

    Options are:

    • ’either’ - closest edge to val;

    • ’low’ - closest edge lower than val; or

    • ’high’ - closest edge higher than val. Default is ‘either’

Returns

float – The closest bin edge to the input value

contiguous_bins()[source]

Return a list of Bins, each one containing a continuous segment of data. This is done by comparing the edges of each bin, and if there is a gap between edges, the data is split into separate Bin objects, each containing a contiguous set of data.

Returns

list of Bins: A list of Bins

slice(lo_edge, hi_edge)[source]

Perform a slice over the range of the bins and return a new Bins object. Note that lo_edge and hi_edge values that fall inside a bin will result in that bin being included.

Parameters
  • lo_edge (float) – The start of the slice

  • hi_edge (float) – The end of the slice

Returns

Bins – A new Bins object containing the slice


EnergyBins

class gbm.data.primitives.EnergyBins(counts, lo_edges, hi_edges, exposure)[source]

Bases: gbm.data.primitives.TimeBins

A class defining a set of Energy (count spectra) bins.

Parameters
  • counts (np.array) – The array of counts in each bin

  • lo_edges (np.array) – The low-value edges of the bins

  • hi_edges (np.array) – The high-value edges of the bins

  • exposure (np.array) – The exposure of each bin

Attributes
  • centroids (np.array) – The geometric centroids of the bins

  • counts (np.array) – The counts in each bin

  • count_uncertainty (np.array) – The count uncertainty in each bin

  • exposure (np.array) – The exposure of each bin

  • hi_edges (np.array) – The high-value edges of the bins

  • lo_edges (np.array) – The low-value edges of the bins

  • range (float, float) – The range of the bin edges

  • rates (np.array) – The differential count rate of each bin: counts/(exposure*widths)

  • rate_uncertainty (np.array) – The count rate uncertainty of each bin

  • size (int) – Number of bins

  • widths (np.array) – The widths of the bins

closest_edge(val, which='either')

Return the closest bin edge

Parameters
  • val (float) – Input value

  • which (str, optional) –

    Options are:

    • ’either’ - closest edge to val;

    • ’low’ - closest edge lower than val; or

    • ’high’ - closest edge higher than val. Default is ‘either’

Returns

float – The closest bin edge to the input value

contiguous_bins()

Return a list of Bins, each one containing a continuous segment of data. This is done by comparing the edges of each bin, and if there is a gap between edges, the data is split into separate Bin objects, each containing a contiguous set of data.

Returns

list of Bins: A list of Bins

classmethod from_fits_array(arr, ebounds)[source]

Create an EnergyBins object from a FITS counts array

Parameters
  • arr (np.recarray or astropy.io.fits.fitsrec.FITS_rec) – The FITS counts array

  • ebounds (np.recarray or astropy.io.fits.fitsrec.FITS_rec) – The ebounds array, mapping channel numbers to energy bins

Returns

EnergyBins – The new EnergyBins object

classmethod merge(histos)

Merge multiple TimeBins together.

Parameters

histos (list of TimeBins) – A list containing the TimeBins to be merged

Returns

TimeBins – A new TimeBins object containing the merged TimeBins

rebin(method, *args, emin=None, emax=None)[source]

Rebin the EnergyBins object in given a binning function and return a a new EnergyBins object

The binning function should take as input an array of counts, array of exposures, and an array of bin edges. Additional arguments specific to the function are allowed. The function should return an array of the new counts, new exposure, and new edges.

Parameters
  • method (<function>) – A binning function

  • *args – Arguments to be passed to the binning function

  • emin (float, optional) – If set, defines the starting energy of the EnergyBins to be binned, otherwise binning will begin at the first bin edge.

  • emax (float, optional) – If set, defines the ending energy of the EnergyBins to be binned, otherwise binning will end at the last bin edge.

Returns

EnergyBins – The rebinned EnergyBins object

slice(tstart, tstop)

Perform a slice over a time range and return a new Bins object. Note that tstart and tstop values that fall inside a bin will result in that bin being included.

Parameters
  • lo_edge (float) – The start of the slice

  • hi_edge (float) – The end of the slice

Returns

TimeBins – A new TimeBins object containing the time slice

classmethod sum(histos)[source]

Sum multiple EnergyBins together if they have the same energy range (support). Example use would be summing two count spectra.

Parameters

histos (list of EnergyBins) – A list containing the EnergyBins to be summed

Returns

EnergyBins – A new EnergyBins object containing the summed EnergyBins


TimeBins

class gbm.data.primitives.TimeBins(counts, lo_edges, hi_edges, exposure)[source]

Bases: gbm.data.primitives.Bins

A class defining a set of Time History (lightcurve) bins.

Parameters
  • counts (np.array) – The array of counts in each bin

  • lo_edges (np.array) – The low-value edges of the bins

  • hi_edges (np.array) – The high-value edges of the bins

  • exposure (np.array) – The exposure of each bin

Attributes
  • centroids (np.array) – The centroids of the bins

  • counts (np.array) – The counts in each bin

  • count_uncertainty (np.array) – The count uncertainty in each bin

  • exposure (np.array) – The exposure of each bin

  • hi_edges (np.array) – The high-value edges of the bins

  • lo_edges (np.array) – The low-value edges of the bins

  • range (float, float) – The range of the bin edges

  • rates (np.array) – The count rate of each bin: counts/exposure

  • rate_uncertainty (np.array) – The count rate uncertainty of each bin

  • size (int) – Number of bins

  • widths (np.array) – The widths of the bins

closest_edge(val, which='either')

Return the closest bin edge

Parameters
  • val (float) – Input value

  • which (str, optional) –

    Options are:

    • ’either’ - closest edge to val;

    • ’low’ - closest edge lower than val; or

    • ’high’ - closest edge higher than val. Default is ‘either’

Returns

float – The closest bin edge to the input value

contiguous_bins()

Return a list of Bins, each one containing a continuous segment of data. This is done by comparing the edges of each bin, and if there is a gap between edges, the data is split into separate Bin objects, each containing a contiguous set of data.

Returns

list of Bins: A list of Bins

classmethod from_fits_array(arr)[source]

Create an TimeBins object from a FITS counts array

Parameters

arr (np.recarray or astropy.io.fits.fitsrec.FITS_rec) – The FITS counts array

Returns

TimeBins – The new TimeBins object

classmethod merge(histos)[source]

Merge multiple TimeBins together.

Parameters

histos (list of TimeBins) – A list containing the TimeBins to be merged

Returns

TimeBins – A new TimeBins object containing the merged TimeBins

rebin(method, *args, tstart=None, tstop=None)[source]

Rebin the TimeBins object in given a binning function and return a a new TimeBins object

The binning function should take as input an array of counts, array of exposures, and an array of bin edges. Additional arguments specific to the function are allowed. The function should return an array of the new counts, new exposure, and new edges.

Parameters
  • method (<function>) – A binning function

  • *args – Arguments to be passed to the binning function

  • tstart (float, optional) – If set, defines the start time of the TimeBins to be binned, otherwise binning will begin at the time of the first bin edge.

  • tstop (float, optional) – If set, defines the end time of the TimeBins to be binned, otherwise binning will end at the time of the last bin edge.

Returns

TimeBins – The rebinned TimeBins object

slice(tstart, tstop)[source]

Perform a slice over a time range and return a new Bins object. Note that tstart and tstop values that fall inside a bin will result in that bin being included.

Parameters
  • lo_edge (float) – The start of the slice

  • hi_edge (float) – The end of the slice

Returns

TimeBins – A new TimeBins object containing the time slice

classmethod sum(histos)[source]

Sum multiple TimeBins together if they have the same time range (support) and the same bin widths. Example use would be summing two histograms at different energies.

Parameters

histos (list of TimeBins) – A list containing the TimeBins to be summed

Returns

TimeBins – A new TimeBins object containing the summed TimeBins


TimeEnergyBins

class gbm.data.primitives.TimeEnergyBins(counts, tstart, tstop, exposure, emin, emax)[source]

Bases: object

A class defining a set of 2D Time-Energy bins.

Parameters
  • counts (np.array) – The array of counts in each bin

  • tstart (np.array) – The low-value edges of the time bins

  • tstop (np.array) – The high-value edges of the time bins

  • exposure (np.array) – The exposure of each bin

  • emin (np.array) – The low-value edges of the energy bins

  • emax (np.array) – The high-value edges of the energy bins

Attributes
  • chanwidths (np.array) – The bin widths along the energy axis

  • count_uncertainty (np.array) – The counts uncertainty in each bin

  • energy_centroids (np.array) – The bin centroids along the energy axis

  • energy_range (float, float) – The range of the data along the energy axis

  • numchans (int) – The number of energy channels along the energy axis

  • numtimes (int) – The number of bins along the time axis

  • rates (np.array) – The rates in each Time-Energy Bin

  • rates_per_kev (np.array) – The differential rates in units of counts/s/keV

  • rate_uncertainty (np.array) – The rate uncertainty in each bin

  • rate_uncertainty_per_kev (np.array) – The differential rate uncertainty in units of counts/s/keV

  • size (int, int) – The number of bins along both axes (numtimes, numchans)

  • time_centroids (np.array) – The bin centroids along the time axis

  • time_range (float, float) – The range of the data along the time axis

  • time_widths (np.array) – The bin widths along the time axis

closest_energy_edge(val, which='either')[source]

Return the closest energy bin edge

Parameters
  • val (float) – Input value

  • which (str, optional) –

    Options are:

    • ’either’ - closest edge to val;

    • ’low’ - closest edge lower than val;

    • ’high’ - closest edge higher than val. Default is ‘either’

Returns

float – The closest energy bin edge to the input value

closest_time_edge(val, which='either')[source]

Return the closest time bin edge

Parameters
  • val (float) – Input value

  • which (str, optional) –

    Options are:

    • ’either’ - closest edge to val;

    • ’low’ - closest edge lower than val;

    • ’high’ - closest edge higher than val. Default is ‘either’

Returns

float – The closest time bin edge to the input value

contiguous_energy_bins()[source]

Return a list of TimeEnergyBins, each one containing a contiguous energy segment of data. This is done by comparing the edges of each energy bin, and if thereis a gap between edges, the data is split into separate TimeEnergyBin objects, each containing an energy-contiguous set of data.

Returns

list of TimeEnergyBins – A list of TimeEnergyBins

contiguous_time_bins()[source]

Return a list of TimeEnergyBins, each one containing a contiguous time segment of data. This is done by comparing the edges of each time bin, and if thereis a gap between edges, the data is split into separate TimeEnergyBin objects, each containing a time-contiguous set of data.

Returns

list of TimeEnergyBins – A list of TimeEnergyBins

get_exposure(time_ranges=None, scale=False)[source]

Calculate the total exposure of a time range or time ranges of data

Parameters
  • time_ranges ([(float, float), ...], optional) – The time range or time ranges over which to calculate the exposure. If omitted, calculates the total exposure of the data

  • scale (bool, optional) – If True and the time ranges don’t match up with the data binning, will scale the exposure based on the requested time range. Default is False.

Returns

float – The exposure of the time selections

integrate_energy(emin=None, emax=None)[source]

Integrate the histogram over the energy axis (producing a lightcurve). Limits on the integration smaller than the full range can be set.

Parameters
  • emin (float, optional) – The low end of the integration range. If not set, uses the lowest energy edge of the histogram

  • emax (float, optional) – The high end of the integration range. If not set, uses the highest energy edge of the histogram

Returns

TimeBins – A TimeBins object containing the lightcurve histogram

integrate_time(tstart=None, tstop=None)[source]

Integrate the histogram over the time axis (producing a count rate spectrum). Limits on the integration smaller than the full range can be set.

Parameters
  • tstart (float, optional) – The low end of the integration range. If not set, uses the lowest time edge of the histogram

  • tstop (float, optional) – The high end of the integration range. If not set, uses the highest time edge of the histogram

Returns

EnergyBins – A EnergyBins object containing the count rate spectrum histogram

classmethod merge_energy(histos)[source]

Merge multiple TimeEnergyBins together along the energy axis.

Parameters

histos (list of TimeEnergyBins) – A list containing the TimeEnergyBins to be merged

Returns

TimeEnergyBins – A new TimEnergyBins object containing the merged TimeEnergyBins

classmethod merge_time(histos)[source]

Merge multiple TimeEnergyBins together along the time axis.

Parameters

histos (list of TimeEnergyBins) – A list containing the TimeEnergyBins to be merged

Returns

TimeEnergyBins: A new TimEnergyBins object containing the merged TimeEnergyBins

rebin_energy(method, *args, emin=None, emax=None)[source]

Rebin the TimeEnergyBins object along the energy axis given a binning function and return a new TimeEnergyBins object

The binning function should take as input an array of counts, array of exposures, and an array of bin edges. Additional arguments specific to the function are allowed. The function should return an array of the new counts, new exposure, and new edges.

Parameters
  • method (<function>) – A binning function

  • *args – Arguments to be passed to the binning function

  • emin (float, optional) – If set, defines the starting edge of the TimeEnergyBins to be binned, otherwise binning will begin at the the first bin edge.

  • emax (float, optional) – If set, defines the ending edge of the TimeEnergyBins to be binned, otherwise binning will end at the last bin edge.

Returns

TimeEnergyBins – The rebinned TimeEnergyBins object

rebin_time(method, *args, tstart=None, tstop=None)[source]

Rebin the TimeEnergyBins object along the time axis given a binning function and return a new TimeEnergyBins object

The binning function should take as input an array of counts, array of exposures, and an array of bin edges. Additional arguments specific to the function are allowed. The function should return an array of the new counts, new exposure, and new edges.

Parameters
  • method (<function>) – A binning function

  • *args – Arguments to be passed to the binning function

  • tstart (float, optional) – If set, defines the start time of the TimeEnergyBins to be binned, otherwise binning will begin at the time of the first bin edge.

  • tstop (float, optional) – If set, defines the end time of the TimeEnergyBins to be binned, otherwise binning will end at the time of the last bin edge.

Returns

TimeEnergyBins – The rebinned TimeEnergyBins object

slice_energy(emin, emax)[source]

Perform a slice over an energy range and return a new TimeEnergyBins object. Note that emin and emax values that fall inside a bin will result in that bin being included.

Parameters
  • emin (float) – The start of the slice

  • emax (float) – The end of the slice

Returns

TimeEnergyBins – A new TimeEnergyBins object containing the energy slice

slice_time(tstart, tstop)[source]

Perform a slice over a time range and return a new TimeEnergyBins object. Note that tstart and tstop values that fall inside a bin will result in that bin being included.

Parameters
  • tstart (float) – The start of the slice

  • tstop (float) – The end of the slice

Returns

TimeEnergyBins – A new TimeEnergyBins object containing the time slice


EventList

class gbm.data.primitives.EventList[source]

Bases: object

A primitive class defining a TTE event list.

Attributes
  • channel_range (int, int) – The range of the channels in the list

  • count_spectrum (EnergyBins) – The count spectrum histogram

  • emax (np.array) – The PHA upper energy edges

  • emin (np.array) – The PHA lower energy edges

  • energy_range (float, float) – The range of the energy edges in the list

  • numchans (int) – The number of energy channels. Note that not all channels will necessarily have events, especially if a slice is made over energy.

  • pha (np.array) – The PHA channel array

  • size (int) – The number of events in the list

  • time (np.array) – The time array

  • time_range (float, float) – The range of the times in the list

bin(method, *args, tstart=None, tstop=None, event_deadtime=2.6e-06, overflow_deadtime=1e-05, **kwargs)[source]

Bin the EventList in time given a binning function and return a 2D time-energy channel histogram.

The binning function should take as input an array of times as well as a tstart and tstop keywords for partial list binning. Additional arguments and keyword arguments specific to the function are allowed. The function should return an array of time edges for the bins, such that, for n bins, there are n + 1 edges.

Parameters
  • method (<function>) – A binning function

  • *args – Arguments to be passed to the binning function

  • tstart (float, optional) – If set, defines the start time of the EventList to be binned, otherwise binning will begin at the time of the first event.

  • tstop (float, optional) – If set, defines the end time of the EventList to be binned, otherwise binning will end at the time of the last event.

  • event_deadtime (float, optional) – The deadtime per event in seconds. Default is 2.6e-6.

  • overflow_deadtime (float, optional) – The deadtime per event in the overflow channel in seconds. Default is 1e-5.

  • **kwargs – Options to be passed to the binning function

Returns

TimeEnergyBins – A Time-Energy Channel histogram

channel_slice(chanlo, chanhi)[source]

Perform a slice in energy channels of the EventList and return a new EventList

Parameters
  • chanlo (int) – The start of the channel slice

  • chanhi (int) – The end of the channel slice

Returns

EventList – A new EventList object containing the channel slice

energy_slice(emin, emax)[source]

Perform a slice in energy of the EventList and return a new EventList. Since energies are binned, an emin or emax falling inside of an energy channel bin will include that bin in the slice.

Parameters
  • emin (float) – The start of the energy slice

  • emax (float) – The end of the energy slice

Returns

EventList – A new EventList object containing the energy slice

classmethod from_fits_array(events_arr, ebounds)[source]

Create an EventList object from TTE FITS arrays

Parameters
  • events_arr (np.recarray or astropy.io.fits.fitsrec.FITS_rec) – The TTE events array

  • ebounds (np.recarray or astropy.io.fits.fitsrec.FITS_rec) – The TTE Ebounds array, mapping channel numbers to energy bins

Returns

EventList – The new EventList

classmethod from_lists(times_list, pha_list, chan_lo, chan_hi)[source]

Create an EventList object from lists of times, channels, and the channel bounds. The list of times and channels must be the same length, and the list of channel boundaries are used to map the PHA index number in pha_list to energy channels.

Parameters
  • times_list (list of float) – A list of event times

  • pha_list (list of int) – A list of PHA channels. Must be same length as times_list.

  • chan_lo (list of float) – A list of the lower edges for the energy channels.

  • chan_hi (list of float) – A list of the upper edges for the energy channels. Must be same length as chan_hi.

Returns

EventList – The new EventList

get_exposure(time_ranges=None, event_deadtime=2.6e-06, overflow_deadtime=1e-05)[source]

Calculate the total exposure of a time range or time ranges of data

Parameters
  • time_ranges ([(float, float), ...], optional) – The time range or time ranges over which to calculate the exposure. If omitted, calculates the total exposure of the data.

  • event_deadtime (float, optional) – The deadtime per event in seconds. Default is 2.6e-6.

  • overflow_deadtime (float, optional) – The deadtime per event in the overflow channel in seconds. Default is 1e-5.

Returns

float – The exposure of the time selections

classmethod merge(eventlists, sort_attrib=None, force_unique=False)[source]

Merge multiple EventLists together in time and optionally sort.

Parameters
  • eventlist (list of EventList) – A list containing the EventLists to be merged

  • sort_attrib (str, optional) – The name of the EventList attribute to sort, either ‘TIME’ or ‘PHA’

  • force_unique (bool, optional) – If True, force all events to be unique via brute force sorting. If False, the EventLists will only be checked and masked for overlapping time ranges. Events can potentially be lost if the merged EventLists contain overlapping times (but not necessarily duplicate events), however this method is much faster. Default is False.

Returns

EventList – A new EventList object containing the merged EventLists

rebin_energy(method, *args)[source]

Rebin the PHA channels using the specified binning algorithm. This does not change the number of events in the EventList, but changes their assignment to a PHA channel and bins the energy bounds mapping to those channels. A new EventList object is returned.

Parameters
  • method (<function>) – A binning function

  • *args – Arguments to be passed to the binning function

Returns

class:`EventList: A new EventList object with the rebinned PHA channels

sort(attrib)[source]

In-place sort by attribute. Either by ‘TIME’ or ‘PHA’

Parameters

attrib (str) – The name of the EventList attribute, either ‘TIME’ or ‘PHA’

time_slice(tstart, tstop)[source]

Perform a slice in time of the EventList and return a new EventList

Parameters
  • tstart (float) – The start of the time slice

  • tstop (float) – The end of the time slice

Returns

EventList – A new EventList object containing the time slice


GTI

class gbm.data.primitives.GTI(tstart=None, tstop=None)[source]

Bases: object

A primitive class defining a set of Good Time Intervals (GTIs)

Parameters
  • tstart (float) – The start time of the GTI

  • tstop (float) – The end time of the GTI

Attributes
  • num_intervals (int) – The number of intervals in the GTI

  • range (float, float) – The full range of the GTI

as_list()[source]

Return the GTI as a list of tuples.

Returns

[(float, float), …] – The list of GTI interval tuples

contains(a_time, inclusive=True)[source]

Determine if the GTI contains a time.

Parameters
  • a_time (float) – The input time to check

  • inclusive (bool, optional) – If True, then includes the edges of the range for the check, otherwise it is edge-exclusive. Default is True.

Returns

bool – True if the time is in the GTI, False otherwise

classmethod from_boolean_mask(times, mask)[source]

Create a new GTI object from a list of times and a Boolean mask Splits the boolean mask into segments of contiguous values and applies to array of times to create a GTI object.

Parameters
  • times (np.array) – An array of times

  • mask (np.array(dtype=bool)) – The boolean array. Must be the same size as times.

Returns

GTI – The new GTI object

classmethod from_list(gti_list)[source]

Create a new GTI object from a list of tuples.

Parameters

gti_list ([(float, float), ...]) – A list of interval tuples

Returns

GTI – The new GTI object

insert(tstart, tstop)[source]

Insert a new interval into the GTI

Parameters
  • tstart (float) – The start time of the new interval

  • tstop (float) – The end time of the new interval

classmethod merge(gti1, gti2)[source]

Return a new GTI object that is a merge of two existing GTI objects.

Parameters
  • gti1 (GTI) – A GTI to be merged

  • gti2 (GTI) – A GTI to be merged

Returns

GTI – The new merged GTI object

property range

The full range of the GTI

Type

(float, float)


TimeRange

class gbm.data.primitives.TimeRange(tstart, tstop)[source]

Bases: object

A primitive class defining a time range

Parameters
  • tstart (float) – The start time of the range

  • tstop (float) – The end time of the range

Attributes
  • center (float) – The center of the time range

  • duration (float) – The duration of the time range

  • tstart (float) – The start time of the range

  • tstop (float) – The end time of the range

as_tuple()[source]

Return the time range as a tuple.

Returns

(float, float) – The starting and ending time of the time range

contains(a_time, inclusive=True)[source]

Determine if the time range contains a time.

Parameters
  • a_time (float) – The input time to check

  • inclusive (bool, optional) – If True, then includes the edges of the range for the check, otherwise it is edge-exclusive. Default is True.

Returns

bool – True if the time is in the time range, False otherwise

classmethod intersection(range1, range2)[source]

Return a new TimeRange that is the intersection of two input TimeRanges. If the input TimeRanges do not intersect, then None is returned.

Parameters
  • range1 (TimeRange) – A time range used to calculate the intersection

  • range2 (TimeRange) – Another time range used to calculate the intersection

Returns

TimeRange – The intersected time range

classmethod union(range1, range2)[source]

Return a new TimeRange that is the union of two input TimeRanges

Parameters
  • range1 (TimeRange) – A time range used to calculate the union

  • range2 (TimeRange) – Another time range used to calculate the union

Returns

TimeRange – The unionized time range