pint.event_toas
Generic functions to load TOAs from events files, along with specific implementations for different missions.
The versions that look like get_..._TOAs()
are preferred: the others are retained for backward compatibility.
Instrument-specific Functions
- pint.event_toas.get_NuSTAR_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a NuSTAR file as a
pint.toa.TOAs
objectCorrectly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
ephem (str, optional) – The name of the solar system ephemeris to use; defaults to “DE421”.
planets (bool, optional) – Whether to apply Shapiro delays based on planet positions. Note that a long-standing TEMPO2 bug in this feature went unnoticed for years. Defaults to False.
- Return type:
See also
- pint.event_toas.get_NICER_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a NICER file as a
pint.toa.TOAs
objectCorrectly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
ephem (str, optional) – The name of the solar system ephemeris to use; defaults to “DE421”.
planets (bool, optional) – Whether to apply Shapiro delays based on planet positions. Note that a long-standing TEMPO2 bug in this feature went unnoticed for years. Defaults to False.
- Return type:
See also
- pint.event_toas.get_RXTE_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a RXTE file as a
pint.toa.TOAs
objectCorrectly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
ephem (str, optional) – The name of the solar system ephemeris to use; defaults to “DE421”.
planets (bool, optional) – Whether to apply Shapiro delays based on planet positions. Note that a long-standing TEMPO2 bug in this feature went unnoticed for years. Defaults to False.
- Return type:
See also
- pint.event_toas.get_IXPE_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a IXPE file as a
pint.toa.TOAs
objectCorrectly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
ephem (str, optional) – The name of the solar system ephemeris to use; defaults to “DE421”.
planets (bool, optional) – Whether to apply Shapiro delays based on planet positions. Note that a long-standing TEMPO2 bug in this feature went unnoticed for years. Defaults to False.
- Return type:
See also
- pint.event_toas.get_Swift_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a Swift file as a
pint.toa.TOAs
objectCorrectly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
ephem (str, optional) – The name of the solar system ephemeris to use; defaults to “DE421”.
planets (bool, optional) – Whether to apply Shapiro delays based on planet positions. Note that a long-standing TEMPO2 bug in this feature went unnoticed for years. Defaults to False.
- Return type:
See also
- pint.event_toas.get_XMM_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a XMM file as a
pint.toa.TOAs
objectCorrectly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
ephem (str, optional) – The name of the solar system ephemeris to use; defaults to “DE421”.
planets (bool, optional) – Whether to apply Shapiro delays based on planet positions. Note that a long-standing TEMPO2 bug in this feature went unnoticed for years. Defaults to False.
- Return type:
See also
- pint.event_toas.load_NuSTAR_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a NuSTAR file as PINT
TOA
objects.Correctly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
- Returns:
toalist
- Return type:
list of
TOA
objects
Note
This list should be converted into a
TOAs
object withpint.toa.get_TOAs_list()
for most operationsSee also
- pint.event_toas.load_NICER_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a NICER file as PINT
TOA
objects.Correctly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
- Returns:
toalist
- Return type:
list of
TOA
objects
Note
This list should be converted into a
TOAs
object withpint.toa.get_TOAs_list()
for most operationsSee also
- pint.event_toas.load_RXTE_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a RXTE file as PINT
TOA
objects.Correctly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
- Returns:
toalist
- Return type:
list of
TOA
objects
Note
This list should be converted into a
TOAs
object withpint.toa.get_TOAs_list()
for most operationsSee also
- pint.event_toas.load_IXPE_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a IXPE file as PINT
TOA
objects.Correctly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
- Returns:
toalist
- Return type:
list of
TOA
objects
Note
This list should be converted into a
TOAs
object withpint.toa.get_TOAs_list()
for most operationsSee also
- pint.event_toas.load_Swift_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a Swift file as PINT
TOA
objects.Correctly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
- Returns:
toalist
- Return type:
list of
TOA
objects
Note
This list should be converted into a
TOAs
object withpint.toa.get_TOAs_list()
for most operationsSee also
- pint.event_toas.load_XMM_TOAs(eventname[, minmjd, maxmjd, errors, ephem, planets])
Read photon event times out of a XMM file as PINT
TOA
objects.Correctly handles raw event files, or ones processed with axBary to have barycentered TOAs. Different conditions may apply to different missions.
The minmjd/maxmjd parameters can be used to avoid instantiation of TOAs we don’t want, which can otherwise be very slow.
- Parameters:
eventname (str) – File name of the FITS event list
minmjd (float, default "-infinity") – minimum MJD timestamp to return
maxmjd (float, default "infinity") – maximum MJD timestamp to return
errors (astropy.units.Quantity or float, optional) – The uncertainty on the TOA; if it’s a float it is assumed to be in microseconds
- Returns:
toalist
- Return type:
list of
TOA
objects
Note
This list should be converted into a
TOAs
object withpint.toa.get_TOAs_list()
for most operationsSee also
Functions
|
|
|
|
|
Read photon event times out of a FITS file as a |
|
Read photon event times out of a FITS file as |
|
Read photon event times out of a FITS file as PINT |
|
Read photon event times out of a FITS file as a list of PINT |
Read all the relevant information about missions in xselect.mdb. |