BatseTte

class gdt.missions.cgro.batse.tte.BatseTte[source]

Bases: PhotonList

Class for BATSE Time-Tagged Event data.

Note

The deadtime for BATSE TTE is assumed to be 3.3 microsec per event based on the analysis of Gjesteland et al. (2010).

References

`Gjesteland, T. et al. 2010, Journal of Geophysical Research: Space

Physics, 115, A00E21

<https://agupubs.onlinelibrary.wiley.com/doi/10.1029/2009JA014578>`_

Attributes Summary

data

The event data

detector

The BATSE detector(s)

ebounds

The energy-channel mapping

ecalib

Energy calibration data

energy_range

The energy range of the data

event_deadtime

Per event deadtime for non-overflow channels

filename

The filename

gti

The good time intervals

hdulist

The list of Header Data Units

headers

The headers

num_chans

The number of energy channels

num_hdus

The number of HDUs

overflow_deadtime

Per event deadtime for overflow channels

time_range

The time range of the data

trigtime

The trigger time of the data, if available.

Methods Summary

close()

Close the file

column(hdu_num, col_name)

Return a column from an HDU as an array.

columns_as_array(hdu_num, col_names[, dtype])

Return a list of columns from an HDU as an array.

from_data(data[, gti, trigger_time, ...])

Create a BATSE TTE object from data.

get_column_names(hdu_num)

Get the column names in a HDU.

get_exposure([time_ranges])

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

merge(ttes[, primary, force_unique])

Merge a list of PhotonList objects into a single new PhotonList o bject.

open(file_path, **kwargs)

Open a BATSE TTE FITS file and return either a BatseTteTrigger or BatseTteMulti object depending on the type of file.

rebin_energy(method, *args, **kwargs)

Rebin the PhotonList in energy given a rebinning method.

set_ebounds(ebounds)

Set the energy calibration (ebounds) of the data.

slice_energy(energy_ranges)

Slice the BatseTte by one or more energy range.

slice_time(time_ranges)

Slice the BatseTte by one or more time range.

to_pha([time_ranges, energy_range, ...])

Integrate the PhotonList data over one or more time ranges to produce a PHA object

to_phaii(bin_method, *args[, time_range, ...])

Convert the TTE data to PHAII data by binning the data in time.

to_spectrum([time_range, energy_range, ...])

Integrate the PhotonList data over time to produce a count spectrum

write(directory[, filename])

Write the file to disk.

Attributes Documentation

data

The event data

Type:

(EventList)

detector

The BATSE detector(s)

Type:

(BatseDetector) or list

ebounds

The energy-channel mapping

Type:

(Ebounds)

ecalib

Energy calibration data

Type:

(BatseEnergyCalib)

energy_range

The energy range of the data

Type:

(float, float)

event_deadtime = 3.3e-06

Per event deadtime for non-overflow channels

filename

The filename

Type:

(str)

gti

The good time intervals

Type:

(Gti)

hdulist

The list of Header Data Units

Type:

(astropy.io.fits.hdu.HDUList)

headers

The headers

Type:

(FileHeaders)

num_chans

The number of energy channels

Type:

(int)

num_hdus

The number of HDUs

Type:

(int)

overflow_deadtime = 3.3e-06

Per event deadtime for overflow channels

time_range

The time range of the data

Type:

(float, float)

trigtime

The trigger time of the data, if available.

Type:

(float)

Methods Documentation

close()

Close the file

column(hdu_num: int, col_name: str) array

Return a column from an HDU as an array.

Parameters:
  • hdu_num (int) – The HDU number

  • col_name (str) – The name of the column

Returns:

(np.array)

columns_as_array(hdu_num: int, col_names: List[str], dtype: dtype = None) array

Return a list of columns from an HDU as an array.

Parameters:
  • hdu_num (int) – The HDU number

  • col_names (list of str) – The names of the columns

  • dtype (np.dtype, optional) – The custom dtype of the output array

Returns:

(np.array)

classmethod from_data(data, gti=None, trigger_time=None, filename=None, headers=None, ecalib=None, detector=None, **kwargs)[source]

Create a BATSE TTE object from data.

Parameters:
  • data (TimeEnergyBins) – The PHAII data

  • gti (Gti, optional) – The Good Time Intervals object. If omitted, the GTI spans (tstart, tstop)

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

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

  • headers (FileHeaders) – The file headers

  • ecalib (BatseEnergyCalib) – The detector calibration

Returns:

(BatseTte)

get_column_names(hdu_num: int)

Get the column names in a HDU. Returns empty if there is no data extension in the HDU.

Parameters:

hdu_num (int) – The HDU number

Returns:

(tuple)

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)

classmethod merge(ttes, primary=0, force_unique=True)

Merge a list of PhotonList objects into a single new PhotonList o bject.

Warning

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

Parameters:
  • ttes (list) – A list of PhotonList objects to merge

  • primary (int) – The index into the list of PhotonList objects to designate the primary PhotonList object. The primary object will be the reference for the header information for the new merged PhotonList object. Default is the first PhotonList 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:

(PhotonList)

classmethod open(file_path, **kwargs)[source]

Open a BATSE TTE FITS file and return either a BatseTteTrigger or BatseTteMulti object depending on the type of file.

Parameters:

file_path (str) – The file path of the FITS file

Returns:

(BatseTteTrigger or BatseTteMulti)

rebin_energy(method, *args, **kwargs)

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

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

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

Returns

(PhotonList)

set_ebounds(ebounds)

Set the energy calibration (ebounds) of the data. If the data already has an energy calibration, this method will update the calibration to the new ebounds.

Parameters:

ebounds (Ebounds) – The ebounds

slice_energy(energy_ranges)[source]

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

Parameters:

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

Returns:

(BatseTte)

slice_time(time_ranges)[source]

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

Parameters:

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

Returns:

(BatseTte)

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

Integrate the PhotonList data over one or more time ranges to produce a PHA object

Note::

If the data does not have an energy calibration (ebounds), then a PHA object cannot be created and calling this method will raise an exception.

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.PHA.from_data()

Returns:

(PHA)

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

Convert the TTE data to PHAII data by binning the data in time.

Parameters:
  • bin_method (<function>) – A binning function for unbinned data

  • *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.

  • phaii_class (class) – The Phaii subclass that the data will be converted to. Default is the base Phaii class.

  • headers (FileHeaders, optional) – The PHAII headers

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

Returns:

(BatsePhaiiTrigger)

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

Integrate the PhotonList 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)

write(directory: Union[str, Path], filename: str = None, **kwargs)

Write the file to disk.

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

  • filename (str, optional) – The filename. If omitted, attempts to use the filename if set.