pint.toa.get_TOAs_array

pint.toa.get_TOAs_array(times, obs, scale=None, errors=<Quantity 1. us>, freqs=<Quantity inf MHz>, flags=None, model=None, ephem=None, include_bipm=True, bipm_version='BIPM2021', include_gps=True, planets=False, tdb_method='default', commands=None, hashes=None, limits='warn', tzr=False, **kwargs)[source]

Load and prepare TOAs for PINT use from an array of times.

Creates TOAs from a an array of times, applies clock corrections, computes key values (like TDB), computes the observatory position and velocity vectors.

Observatory clock corrections are also applied, and thus the exact result also depends on the precise values in observatory clock correction files; normally these do not change. See pint.toa.TOAs.apply_clock_corrections() f or further information on the meaning of the clock correction flags.

Note that this requires all TOAs to be measured at a single observatory (which can be a spacecraft). If you wish to combine multiple data-sets use pint.toa.merge_TOAs()

Parameters:
  • times (astropy.time.Time, float, numpy.ndarray, or tuple of floats/numpy.ndarray) – The times of the TOAs, which can be expressed as astropy Time values, floating point MJDs (64 or 80 bit precision), or a tuple of (MJD1,MJD2) values whose sum is the full precision MJD (usually the integer and fractional part of the MJD)

  • obs (str) – The observatory code for the TOA

  • errors (astropy.units.Quantity or numpy.ndarray or float, optional) – The uncertainty on the TOA; Assumed to be in microseconds if no units provided

  • freqs (astropy.units.Quantity or numpy.ndarray or float, optional) – Frequency corresponding to the TOAs; Assumed to be in MHz if no units provided

  • scale (str, optional) – Time scale for the TOAs. Defaults to the timescale appropriate to the site, but can be overridden

  • flags (dict or iterable) – Flags associated with the TOAs. If iterable, then must have same length as times. Any additional keyword arguments are interpreted as flags as well.

  • ephem (str) – The name of the solar system ephemeris to use; defaults to “DE421”.

  • include_bipm (bool or None) – Whether to apply the BIPM clock correction. Defaults to True.

  • bipm_version (str or None) – Which version of the BIPM tables to use for the clock correction. The format must be ‘BIPMXXXX’ where XXXX is a year.

  • include_gps (bool or None) – Whether to include the GPS clock correction. Defaults to True.

  • planets (bool or None) – 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.

  • model (pint.models.timing_model.TimingModel or None) – If a valid timing model is passed, model commands (such as BIPM version, planet shapiro delay, and solar system ephemeris) that affect TOA loading are applied.

  • tdb_method (str) – Which method to use for the clock correction to TDB. See pint.observatory.Observatory.get_TDBs() for details.

  • hashes (dict) – A dictionary of hashes of the files this data was read from (if any). This is used by check_hashes() to verify whether the data on disk has changed so that the file can be re-read if necessary.

  • limits ("warn" or "error") – What to do when encountering TOAs for which clock corrections are not available.

  • tzr (bool) – Whether the TOAs object corresponds to a TZR TOA

Returns:

TOAs – Completed TOAs object representing the data.

Return type:

pint.toa.TOAs

Examples

>>> from pint import pulsar_mjd
>>> from pint import toa
>>> from astropy import time
>>> t = pulsar_mjd.Time(np.linspace(55000, 56000, 100), scale="utc", format="pulsar_mjd")
>>> errors = np.random.normal(loc=1, scale=0.1, size=len(t)) << u.us
>>> freqs = np.ones(len(t)) * 1400 << u.MHz
>>> obs = "gbt"
>>> flags = {"quality": "good"}
>>> # additional kwargs are converted to flags
>>> # they are automatically made into strings
>>> toas = toa.get_TOAs_array(
        t,
        obs,
        errors=errors,
        freqs=freqs,
        ephem="DE440",
        flags=flags,
        int_time=np.ones(len(t)) * 30,
    )