PhotonList¶
- class gdt.core.tte.PhotonList[source]¶
Bases:
FitsFileContextManagerClass for Photon Lists or Time-Tagged Event data
Attributes Summary
The event data
The energy-channel mapping
The energy range of the data
The deadtime imposed per event
The filename
The good time intervals
The list of Header Data Units
The headers
The number of energy channels
The number of HDUs
The per-event deadtime imposed by events in the overflow channel
The time range of the data
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 PhotonList object from an
EventListobject.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.
hdu_index_from_name(ext_name)Returns the index into the HDU list based on the extension name.
merge(ttes[, primary, force_unique])Merge a list of PhotonList objects into a single new PhotonList o bject.
open(file_path[, mode, memmap])Open a FITS 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, **kwargs)Slice the PhotonList by one or more energy ranges.
slice_time(time_ranges, **kwargs)Slice the PhotonList 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 PhotonList 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
- energy_range¶
The energy range of the data
- Type:
(float, float)
- event_deadtime¶
The deadtime imposed per event
- Type:
(float)
- filename¶
The filename
- Type:
(str)
- hdulist¶
The list of Header Data Units
- Type:
(astropy.io.fits.hdu.HDUList)
- headers¶
The headers
- Type:
- num_chans¶
The number of energy channels
- Type:
(int)
- num_hdus¶
The number of HDUs
- Type:
(int)
- overflow_deadtime¶
The per-event deadtime imposed by events in the overflow channel
- Type:
(float)
- 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, event_deadtime=0.0, overflow_deadtime=0.0, **kwargs)[source]¶
Create a PhotonList object from an
EventListobject.- Parameters:
data (
EventList) – The event datagti (
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) – The filename
headers (
FileHeaders) – The file headersevent_deadtime (float, optional) – The deadtime imposed per event. Default is 0.0.
overflow_deadtime (float, optional) – The per-event deadtime imposed by events in the overflow channel. Default is 0.0.
- Returns:
- 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)[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)
- hdu_index_from_name(ext_name)¶
Returns the index into the HDU list based on the extension name. If there is no match, returns None.
- Parameters:
ext_name (str) – The extension name
- Returns:
(int)
- classmethod merge(ttes, primary=0, force_unique=True)[source]¶
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:
- classmethod open(file_path: Union[str, Path], mode: str = 'readonly', memmap: bool = None)¶
Open a FITS file.
- Parameters:
file_path (str) – The file path
mode (str) – The file access mode
memmap (bool) – If True, memory map when reading the file
- Returns:
(
FitsFileContextManager)
- rebin_energy(method, *args, **kwargs)[source]¶
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
- set_ebounds(ebounds)[source]¶
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, **kwargs)[source]¶
Slice the PhotonList by one or more energy ranges. Produces a new PhotonList object.
- Note::
If the data does not have an energy calibration (ebounds), then this function will slice by energy channels, and therefore the
energy_rangesargument should be a range(s) of energy channels instead of energies.
- Parameters:
energy_ranges ([(float, float), ...]) – The energy ranges to slice the data to.
- Returns:
- slice_time(time_ranges, **kwargs)[source]¶
Slice the PhotonList by one or more time range. Produces a new PhotonList object.
- Parameters:
time_ranges ([(float, float), ...]) – The time ranges to slice the data to.
- Returns:
- to_pha(time_ranges=None, energy_range=None, channel_range=None, **kwargs)[source]¶
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, phaii_class=<class 'gdt.core.phaii.Phaii'>, headers=None, **kwargs)[source]¶
Convert the PhotonList data to PHAII data by binning the data in time.
- Note::
If the data has no energy calibration, then
energy_rangeis ignored, and onlychannel_rangeis used.
- 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
Phaiiclass.headers (
FileHeaders, optional) – The PHAII headers**kwargs – Options to pass to the binning function
- Returns:
(
Phaii)
- to_spectrum(time_range=None, energy_range=None, channel_range=None)[source]¶
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:
- 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
filenameif set.