barycentric_correction

This task barycentric_correction calculates the barycentric correction between HST and the Solar System barycenter for times in HST/STIS observations. It updates the times in uncalibrated (tag, raw) and calibrated data. It stores the old times in a new column for tag data and in new header keywords for raw and calibrated data. The function puts the various delay terms and positions in the header as well.

The function works by first determining the position of the Earth with Astropy and HST with either JPL Horizons or the STScI-provided HST orbital files.

For *tag.fits files with many times to convert, performance can be slow due to the high number of coordinate transformations.

Stand-alone tasks calc_delay_jpl calc_delay_orbfile can be used to calculate barycentric corrections without file modifications.

Task combine_hst_orb_files combines two neighboring HST orbital files together for use on MAST data spanning their coverage.

Task odelay_file_compare shows differences in times between two STIS FITS files. This is useful in what barycentric corrections were calculated.

Examples

barycentric_correction with default values:

>>> import stistools
>>> stistools.barycentric_correction.barycentric_correction("oep502040_x1d.fits")
exception stistools.barycentric_correction.OrbFileError

Orbit file does not cover the requested time range.

stistools.barycentric_correction.barycentric_correction(table_names, verbose=True, distance=1000000000.0, hst_orb=None, time_script=False, outfiles=None, time_system='TDB')

Calculates time-delay barycentric corrections from HST's position to the Solar System barycenter. This correction includes the classic geometric Roemer delay, as well as the general relativistic Einstein delay.

This function uses modern Astropy tools to calculate the position of the barycenter, Earth's location, and deal with the time standards. The function has been tested against a Python implementation of the previous odelaytime IRAF task and found to be consistent. We have also tested against other barycentric correction tools, including barycorr, astroutils, and pintbary.

HST's changing location around the Earth can lead to time-delay differences of up to ~46 milliseconds. HST's location can be determined either through STScI-provided HST orbital files, or through a query to JPL Horizons.

These calculations are accurate to within 1 millisecond outside the Solar System, and to within 5 milliseconds inside the Solar System.

Parameters:
table_names: str or Iterable[str]

List of strings with the file names to be time-corrected.

verbose: bool

Prints completion messages during execution. Default=True

distance: float

Distance the object is from HST in AU. Default is a trillion AU. Most important for objects in our Solar System as it is repsonsible for second-order correction, up to minutes. At 1 parsec, the correction can be on the order of a few ms.

hst_orb: str or None

Name of HST orbital file (generally starts p, ends as .fit) that covers the time of the observations. If not provided, JPL Horizons will be used to get HST's orbital position.

time_script: bool

Set to True if you want to time how long the script takes. Useful for debugging, especially for tag files with large numbers of events.

outfiles: str or list[str] or None

If not None, then it is a list of output files for each table_names. Each table_name will be copied over to the corresponding outfile name. Default None.

time_system: str

Define as either "TDB" or "UTC" for final time standard conversion. They will be different by about 69.184 seconds plus or minus a few ms depending on where Earth is in its orbit. Final results will then be in either BJD_TDB or BJD_UTC. Default is "TDB".

Returns:
None

However, the file is written to output.

Notes

If an OrbFileError is raised, you may need to combine two neighboring HST orbital files using combine_hst_orb_files.

stistools.barycentric_correction.calc_delay_jpl(times, ra, dec, distance=1000000000.0, verbose=True)

Calculate the barycentric light-travel time correction for HST using JPL Horizons.

Compute the barycentric light-travel time delay for the Hubble Space Telescope (HST) at the provided observation times and a target sky position (RA, Dec) by querying JPL Horizons for HST's geocentric state-vector (or interpolating a regular-sampled vector set) and combining it with Earth's barycentric position. A finite-distance correction is applied to account for targets that are not at infinite distance.

Parameters:
times: array-like or `~astropy.time.Time`

Observation times in Modified Julian Date (MJD). Can be a scalar or an array.

ra: float or `~astropy.units.Quantity`

Right ascension of the target in degrees.

dec: float or `~astropy.units.Quantity`

Declination of the target in degrees.

distance: float, optional

Distance to the target (default 1e9). This value is used in the finite-distance correction term. (See Notes for how it enters the equation.)

verbose: bool, optional

If True (default) print diagnostics about interpolation, the finite-distance correction, and the computed light-travel times.

Returns:
lt_time: Quantity

Light-travel time correction(s) (Astropy Quantity) in days. This is the time that should be added to the observation times to obtain barycentric arrival times, and includes the finite-distance correction.

Raises:
ValueError

If the Horizons query or interpolation fails to produce vectors covering the requested time range, or if the vector transformation cannot be performed.

RemoteServiceError

If the Horizons query or interpolation fails to produce vectors covering the requested time range.

ConvertError

If the vector transformation cannot be performed.

Notes

  • Uses the JPL planetary ephemeris (jpl) via Astropy's solar_system_ephemeris for consistency with Horizons.

  • For multiple times the function queries Horizons for HST (-48, location 500) with a regular step (1 minute) over the time range (plus a 5-minute margin) and cubic-interpolates the returned vectors to the requested times. For a single scalar time it queries Horizons directly at that epoch.

  • The finite-distance correction implemented here follows the same algebraic form used elsewhere in this codebase:

    \mbox{correction\_term} = \frac{-0.5}{D} \left(|r|^2 - (n \cdot r)^2\right) / c

    where r is the HST barycentric position vector, n is the unit vector toward the target, D is the provided distance, and c is the speed of light. The expression is converted into days before being returned.

  • The HST vector returned from Horizons is treated as geocentric and then transformed to ITRS/ICRS and added to Earth's barycentric position from Astropy's get_body_barycentric('earth', ...) to obtain HST's barycentric position.

  • Requires astroquery (Horizons), Astropy with the jplephem ephemeris available, and a working internet connection for Horizons queries when not using a local orbit file.

Performance / Accuracy Considerations

  • Interpolation at 1-minute sampling is a compromise between speed and accuracy; a timing error of ~1 minute corresponds to a maximum light-travel time error of a few milliseconds. Increase sampling density or query exact epochs if sub-ms accuracy is required.

  • Converting large arrays of epochs and coordinate transforms may be time-consuming.

Examples

>>> # single time (MJD), RA/Dec in degrees
>>> barycentric_correction.calc_delay_jpl(55521.123, 210.8023, -47.393)
Finite distance correction: -7.804152203899432e-08 s
Light travel times: -405.56807173759836 s
>>> # multiple times
>>> times = [60200.123, 60200.124, 60200.125]
>>> lt_array = calc_delay_jpl(times, 210.8023, -47.393, distance=1e9, verbose=False)
stistools.barycentric_correction.calc_delay_orbfile(times, ra, dec, hst_orb, distance=1000000000.0, verbose=True)

Calculate the light-travel time correction for HST observations using an orbit file.

This function computes the barycentric light-travel time delay for the Hubble Space Telescope (HST) given a set of observation times and a target sky position (RA, Dec). It interpolates HST's position from a provided orbit file, combines it with Earth's barycentric position, and calculates the light-travel time correction to the Solar System barycenter, including a finite-distance correction term.

Parameters:
times: array-like or float

Observation times in Modified Julian Date (MJD), corresponding to HST exposures.

ra: float

Right ascension of the target in degrees.

dec: float

Declination of the target in degrees.

hst_orb: str

Path to the HST orbit FITS file. This file must contain columns TIME, X, Y, and Z giving HST's position (in km) relative to the Earth's center.

distance: float, optional

Distance to the target in kilometers (default is 1e9, effectively infinite distance). Used to apply the finite-distance light-travel time correction.

verbose: bool, optional

If True (default), print information about the finite-distance correction and the calculated light-travel times.

Returns:
lt_time: Quantity

The barycentric light-travel time correction(s) in days, including the finite-distance correction term.

Raises:
OrbFileError

If the orbit file does not cover the requested time range.

Notes

  • Uses the JPL planetary ephemeris (jplephem) for consistency with NASA Horizons.

  • Interpolates HST's orbital position at each observation time using cubic interpolation to avoid excessive Horizons queries.

  • The resulting correction term accounts for HST's motion relative to the Solar System barycenter.

  • The finite-distance correction term scales as ~(r^2 - (r·n)^2) / (2cD), where r is HST's barycentric position vector, n is the unit vector toward the target, c is the speed of light, and D is the target distance.

  • HST orbit files may be downloaded using information provided at https://www.stsci.edu/~STIS/monitors/ephemeris_files.html

Examples

>>> barycentric_correction.calc_delay_orbfile([55521.123], 210.8023, -47.393, hst_orb='pubj0000r.fit')
Finite distance correction: [-7.80563708e-08] s
Light travel times: [-405.56687377] s
stistools.barycentric_correction.combine_hst_orb_files(file1, file2, outname, overwrite=False)

Combine two neighboring hst_orb ephemeris files (ORBs) together. Useful for correcting MAST datasets that are observed over an ORB file boundary.

Parameters:
file1: str

Path to the first hst_orb FITS file.

file2: str

Path to the second hst_orb FITS file.

outname: str

Name/path of combined hst_orb FITS file to output.

overwrite: bool

Overwrite a pre-existing output file? Default=False

Returns:
None
Raises:
OrbFileError

The time gap between orbit files is too large (>120 s)

FileExistsError

Output file already exists and overwrite=False

ValueError

Unexpected duplicate times detected after truncation of later file

stistools.barycentric_correction.odelay_file_compare(file1, file2)

Compare timing information between two FITS files.

Computes and prints the differences in exposure start times and data timestamps between two FITS files, typically used for verifying time coordinate consistency in astronomical observations.

Parameters:
file1: str

Path to the first FITS file.

file2: str

Path to the second FITS file to compare against file1.

Returns:
None

This function only prints comparison results to stdout.

Notes

The function prints the following comparisons (file2 - file1):

  • TEXPSTRT header keyword difference in days and seconds

  • EXPSTART header keyword difference in days and seconds (from extension 1)

  • First and last TIME values difference in seconds (only if 'tag' is in file1 name)

  • TT to TDB time scale conversion difference for reference

The function assumes:

  • Both files have TEXPSTRT in the primary header (extension 0)

  • Both files have EXPSTART in extension 1 header

  • If tag appears in file1 name, extension 1 contains a table. The script will print the difference between the first and last times in each file.

Examples

>>> odelay_file_compare('observation1_tag.fits', 'observation2_tag.fits')
TEXPSTRT f2-f1: -1.5854311641305685e-08 days
TEXPSTRT f2-f1: -0.0013698125258088112 seconds
EXPSTART f2-f1: -1.5854311641305685e-08 days
EXPSTART f2-f1: -0.0013698125258088112 seconds
First TIME f2-f1: 0.0 seconds
Last TIME f2-f1: -0.23225000000002183 seconds
In case it is helpful, the difference between TT and TDB_BJD is -0.001497 s