pint.toa.TOAs

class pint.toa.TOAs(toafile=None, toalist=None, toatable=None, tzr=False)[source]

Bases: object

A class of multiple TOAs, loaded from zero or more files.

Normally these objects should be read from a file with pint.toa.get_TOAs(). Constructing them with the constructor here does not apply the clock corrections and the resulting TOAs object may not be set up the way one would normally expect.

The contents are stored in an astropy.table.Table. Not all columns of the table are computed automatically, as their computation can be expensive. Methods of this class are available to populate these additional columns. Methods that return data from the columns do so in the internal order.

TOAs objects can accept indices that are boolean, list-of-integer, or slice, to produce a new TOAs object containing a subset of the TOAs in the original. For example, to obtain a new TOAs object containing the entries above 1 GHz:

>>> t[t.table['freq'] > 1*u.GHz]

TOAs objects also accept indexing to select columns or flags, so that for example t['mjd'] returns the Column contained in the object, and t['fish'] returns an array of strings that has the value associated with the flag -fish for each TOA that has that flag or the empty string for each TOA that doesn’t. TOAs objects also support assignment through these methods:

>>> t['high', t['freq'] > 1*u.GHz] = "1"

TOAs used to be grouped by observatory, but currently are not. To iterate over all observatory groups:

>>> for obs,idx in t.get_obs_groups():

Or to actually reorder the TOAs by observatory:

>>> t.table.group_by("obs")

Other sorting is fine too:

>>> t.table.sort("freq")
Columns in .table

Name

Contents

index

location of the TOA in the original input

mjd

the exact time of arrival (an astropy.time.Time object)

mjd_float

the time of arrival in floating-point (may be microseconds off)

error

the standard error (an astropy.units.Quantity describing the claimed uncertainty on the pulse arrival time)

freq

the observing frequency (an astropy.units.Quantity)

obs

the observatory at which the TOA was acquired (a pint.observatory.Observatory object)

flags

free-form flags associated with the TOA (a dictionary mapping flag to value)

tdb

the pulse arrival time converted to TDB (but not barycentered, that is, not corrected for light travel time; an astropy.time.Time object); computed by pint.toa.TOAs.compute_TDBs()

tdbld

a longdouble version of tdb for computational convenience

ssb_obs_pos, ssb_obs_vel

position and velocity of the observatory at the time of the TOA; computed by pint.toa.TOAs.compute_posvels()

ssb_obs_vel_ecl

velocity of the observatory in ecliptic coordinates at the time of the TOA; computed by pint.toa.TOAs.add_vel_ecl()

obs_sun_pos, obs_jupiter_pos, obs_saturn_pos, obs_venus_pos, obs_uranus_pos, obs_neptune_pos, obs_earth_pos

position of various celestial objects at the time of the TOA; computed by pint.toa.TOAs.compute_posvels()

pulse_number

integer number of turns since a fiducial moment; optional; can be computed from a model with pint.toa.TOAs.compute_pulse_numbers() or extracted from the pn entry in flags with pint.toa.TOAs.phase_columns_from_flags(). If it is present for some TOAs but not all, missing values are filled with NaN.

delta_pulse_number

number of turns to adjust pulse number by, compared to the model; PHASE statements in the .tim file or the padd entry in flags carry this information, and pint.toa.TOAs.phase_columns_from_flags() creates the column.

Parameters:
  • toafile (str, optional) – Filename to load TOAs from.

  • toalist (list of TOA objects, optional) – The TOA objects this TOAs should contain.

  • toatable (astropy.table.Table, optional) – An existing TOA table

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

  • provided. (Exactly one of these three parameters must be) –

Variables:
  • table (astropy.table.Table) – The data for all the TOAs.

  • commands (list of str) – “Commands” that were written in the file; these will have affected how some or all TOAs were interpreted during loading.

  • filename (str, optional) – The filename (if any) that the TOAs were loaded from.

  • planets (bool) – Whether planetary Shapiro delay should be considered.

  • ephem (object) – The Solar System ephemeris in use.

  • clock_corr_info (dict) – Information about the clock correction chains in use.

  • merged (bool) – If this object was merged from several files (and thus the filename of the first is not useful for referring to the whole object).

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

  • was_pickled (bool) – Whether this file was loaded from a pickle.

  • alias_translation (dict or None) – Translate observatory names by looking them up in this dictionary; this may be necessary to convert observatory names into something TEMPO can understand, or to cope with different setups using different names for the same observatory or the same name for different observatories. There is a dictionary tempo_aliases available to use names as compatible with TEMPO as possible.

  • wideband (bool) – Whether the TOAs also have wideband DM information

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

Methods

add_vel_ecl(obliquity)

Compute and add a column to self.table with velocities in ecliptic coordinates.

adjust_TOAs(delta)

Apply a time delta to TOAs.

apply_clock_corrections([include_bipm, ...])

Apply observatory clock corrections and TIME statments.

check_hashes([timfile])

Determine whether the input files are the same as when loaded.

compute_TDBs([method, ephem])

Compute and add TDB and TDB long double columns to the TOA table.

compute_posvels([ephem, planets])

Compute positions and velocities of the observatories and Earth.

compute_pulse_numbers(model)

Set pulse numbers (in TOA table column pulse_numbers) based on model.

get_all_flags()

Return a list of all the flags used by any TOA.

get_clusters([gap_limit, add_column, add_flag])

Identify toas within gap limit (default 2h = 0.0833d) of each other as the same cluster.

get_dm_errors()

Get the Wideband DM data error.

get_dms()

Get the Wideband DM data.

get_errors()

Return a numpy array of the TOA errors in us.

get_flag_value(flag[, fill_value, as_type])

Get the requested TOA flag values.

get_flags()

Return a numpy array of the TOA flags.

get_freqs()

Return a Quantity of the observing frequencies for the TOAs.

get_highest_density_range([ndays])

Print the range of mjds (default 7 days) with the most toas

get_mjds([high_precision])

Array of MJDs in the TOAs object.

get_obs_groups()

Return an iterator over the different observatories

get_obss()

Return a numpy array of the observatories for each TOA.

get_pulse_numbers()

Return a numpy array of the pulse numbers for each TOA if they exist.

get_summary()

Return a short ASCII summary of the TOAs.

is_wideband()

Whether or not the data have wideband TOA values

merge(t, *args[, strict])

Merge TOAs instances into the existing object

phase_columns_from_flags()

Create and/or modify pulse_number and delta_pulse_number columns.

print_summary()

Prints self.get_summary().

remove_pulse_numbers()

Delete pulse numbers from TOA data

renumber([index_order])

Recreate the index column so the values go from 0 to len(self)-1.

select(selectarray)

Apply a boolean selection or mask array to the TOA table.

to_TOA_list([clkcorr])

Turn a pint.toa.TOAs object into a list of pint.toa.TOA objects

unselect()

Return to previous selected version of the TOA table (stored in stack).

update_all_times([tdb_method])

Update the various derived time columns

update_mjd_float()

Update the mjd_float column from the mjd column

write_TOA_file(filename[, name, format, ...])

Write this object to a .tim file.

Attributes

first_MJD

The first MJD, in Time format.

last_MJD

The last MJD, in Time format.

ntoas

The number of TOAs.

observatories

The set of observatories in use by these TOAs.

wideband

Whether or not the data have wideband TOA values

__add__(other)[source]

Addition operator, allowing merging/concatenation of TOAs

property ntoas

The number of TOAs. Also available as len(toas).

property observatories

The set of observatories in use by these TOAs.

property first_MJD

The first MJD, in Time format.

property last_MJD

The last MJD, in Time format.

property wideband

Whether or not the data have wideband TOA values

to_TOA_list(clkcorr=False)[source]

Turn a pint.toa.TOAs object into a list of pint.toa.TOA objects

This effectively undoes pint.toa.get_TOAs_list(), optionally undoing clock corrections too

Parameters:

clkcorr (bool, optional) – Whether or not to undo any clock corrections

Returns:

Of pint.toa.TOA objects

Return type:

list

is_wideband()[source]

Whether or not the data have wideband TOA values

get_all_flags()[source]

Return a list of all the flags used by any TOA.

get_freqs()[source]

Return a Quantity of the observing frequencies for the TOAs.

get_mjds(high_precision=False)[source]

Array of MJDs in the TOAs object.

With high_precision is True Return an array of the astropy.times (UTC) of the TOAs

With high_precision is False Return an array of toas in mjd as double precision floats

WARNING: Depending on the situation, you may get MJDs in a different scales (e.g. UTC, TT, or TDB) or even a mixture of scales if some TOAs are barycentred and some are not (a perfectly valid situation when fitting both Fermi and radio TOAs)

get_errors()[source]

Return a numpy array of the TOA errors in us.

get_obss()[source]

Return a numpy array of the observatories for each TOA.

get_obs_groups()[source]

Return an iterator over the different observatories

get_pulse_numbers()[source]

Return a numpy array of the pulse numbers for each TOA if they exist.

get_flags()[source]

Return a numpy array of the TOA flags.

get_flag_value(flag, fill_value=None, as_type=None)[source]

Get the requested TOA flag values.

Parameters:
  • flag_name (str) – The request flag name.

  • fill_value – The value to include for missing flags.

  • as_type (callable or None) – If provided, this is called on each value to convert them from to the desired type. All flag values are stored as strings internally.

Returns:

  • values (list) – A list of flag values from each TOA. If the TOA does not have the flag, it will fill up with the fill_value.

  • valid_index (list) – The indices, in self.table, of the places where the flag values occur.

get_dms()[source]

Get the Wideband DM data.

Note

This does not handle situations where some but not all TOAs have DM information.

get_dm_errors()[source]

Get the Wideband DM data error.

Note

This does not handle situations where some but not all TOAs have DM information.

get_clusters(gap_limit=<Quantity 2. h>, add_column=False, add_flag=None)[source]

Identify toas within gap limit (default 2h = 0.0833d) of each other as the same cluster.

Clusters can be larger than the gap limit - if toas are separated by a gap larger than the gap limit, a new cluster starts and continues until another such gap is found.

Cluster info can be added as a clusters column to the pint.toa.TOAs.table object if add_column is True. Cluster info can also be added as a flag with name specified. In those cases self.table.meta['cluster_gap'] will be set to the gap_limit. If the desired clustering corresponds to that in pint.toa.TOAs.table then that column is returned.

Parameters:
  • gap_limit (astropy.units.Quantity, optional) – The minimum size of gap to create a new cluster. Defaults to two hours.

  • add_column (bool, optional) – Whether or not to add a clusters column to the TOA table (default: False)

  • add_flag (str, optional) – If not None, will add a flag with that name to the TOA table whose value is the cluster number (as a string, starting at 0) (default: None)

Returns:

clusters – The cluster number associated to each TOA. Clusters are numbered chronologically from zero.

Return type:

numpy.ndarray

get_highest_density_range(ndays=<Quantity 7. d>)[source]

Print the range of mjds (default 7 days) with the most toas

check_hashes(timfile=None)[source]

Determine whether the input files are the same as when loaded.

Parameters:

timfile (str or list of str or file-like or None) – If provided this should match the list of files the TOAs object was loaded from. If this is a string or list of strings, and the number matches the number of files this TOAs object was loaded from, it is assumed that these are supposed to be the same files, re-named or moved; their contents are then checked. If the contents or the number doesn’t match, this function returns False.

Returns:

True if the contents of the TOAs object matches the content of the files.

Return type:

bool

select(selectarray)[source]

Apply a boolean selection or mask array to the TOA table.

Deprecated. Use toas[selectarray] to get a new object instead.

This operation modifies the TOAs object in place, shrinking its table down to just those TOAs where selectarray is True. This function also stores the old table in a stack.

unselect()[source]

Return to previous selected version of the TOA table (stored in stack).

Deprecated. Use toas[selectarray] to get a new object instead.

get_summary()[source]

Return a short ASCII summary of the TOAs.

This includes summary information about the errors and frequencies but is never more than a few lines regardless of how many TOAs there are.

print_summary()[source]

Prints self.get_summary().

phase_columns_from_flags()[source]

Create and/or modify pulse_number and delta_pulse_number columns.

Scans pulse numbers from the table flags and creates a new table column. Modifes the delta_pulse_number column, if required. Removes the pulse numbers from the flags.

compute_pulse_numbers(model)[source]

Set pulse numbers (in TOA table column pulse_numbers) based on model.

Replace any existing pulse numbers by computing phases according to model and then setting the pulse number of each to their integer part, which the nearest integer since Phase objects ensure that.

Parameters:

model (pint.models.timing_model.TimingModel) – The model defining times of arrival; the pulse numbers assigned will be the nearest integer number of turns to that predicted by the model.

remove_pulse_numbers()[source]

Delete pulse numbers from TOA data

adjust_TOAs(delta)[source]

Apply a time delta to TOAs.

Adjusts the time (MJD) of the TOAs by applying delta, which should be a scalar or have the same shape as self.table['mjd']. This function does not change the pulse numbers column, if present, but does recompute mjd_float, the TDB times, and the observatory positions and velocities.

Parameters:

delta (astropy.time.TimeDelta or astropy.units.Quantity) – The time difference to add to the MJD of each TOA

renumber(index_order=True)[source]

Recreate the index column so the values go from 0 to len(self)-1.

This modifies the TOAs object and also returns it, for calling convenience.

Parameters:

index_order (bool) – If True, preserve the order of the index column, but renumber so there are no gaps. If False, number according to the order TOAs occur in the object.

Return type:

self

write_TOA_file(filename, name='unk', format='tempo2', commentflag=None, order_by_index=True, *, include_pn=True, include_info=True, comment=None)[source]

Write this object to a .tim file.

This function writes the contents of this object to a (single) .tim file. If TEMPO2 format is used, this file is able to represent the contents of this object to nanosecond level. No TIME or EFAC commands are emitted.

Parameters:
  • filename (str or file-like) – File name to write to; can be an open file object

  • name (str) – Value to put in the “name” field of tempo2 files, if a “-name” flag is not available.

  • format (str) – Format specifier for file (‘TEMPO’ or ‘Princeton’) or (‘Tempo2’ or ‘1’); note that not all features may be supported in ‘TEMPO’ mode.

  • commentflag (str or None) – If a string, and that string is a TOA flag, that TOA will be commented in the output file. If None (or non-string), no TOAs will be commented.

  • order_by_index (bool) – If True, write the TOAs in the order specified in the “index” column (which is usually the same as the original file); if False, write them in the order they occur in the TOAs object.

  • include_pn (bool, optional) – Include pulse numbers (if available)

  • include_info (bool, optional) – Include information string if True

  • comment (str, optional) – Additional string to include in TOA file

apply_clock_corrections(include_bipm=True, bipm_version='BIPM2021', include_gps=True, limits='warn')[source]

Apply observatory clock corrections and TIME statments.

Apply clock corrections to all the TOAs where corrections are available. This routine actually changes the value of the TOA, although the correction is also listed as a new flag for the TOA called ‘clkcorr’ so that it can be reversed if necessary. This routine also applies all ‘TIME’ commands (-to flags) and treats them exactly as if they were a part of the observatory clock corrections.

If the clock corrections have already been applied they will not be re-applied.

Options to include GPS or BIPM clock corrections are set to True by default in order to give the most accurate clock corrections.

A description of how PINT handles clock corrections and timescales is here: https://github.com/nanograv/PINT/wiki/Clock-Corrections-and-Timescales-in-PINT

Parameters:
  • include_bipm (bool) – Whether or not to include BIPM correction

  • bipm_version (str) – BIPM version to use. The format must be ‘BIPMXXXX’ where XXXX is a year.

  • include_gps (bool) – Whether or not to include GPS corrections

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

compute_TDBs(method='default', ephem=None)[source]

Compute and add TDB and TDB long double columns to the TOA table.

This routine creates new columns ‘tdb’ and ‘tdbld’ in a TOA table for TDB times, using the Observatory locations and IERS Earth rotation corrections for UT1.

If these columns are already present, delete and replace them.

Parameters:
  • method (str) – Which method to use. See pint.observatory.Observatory.get_TDBs() for details.

  • ephem (str or None) – Solar System ephemeris to use for the computation. If not specified use the value in self.ephem if available, else use pint.toa.EPHEM_default If specified, replace self.ephem.

compute_posvels(ephem=None, planets=None)[source]

Compute positions and velocities of the observatories and Earth.

Compute the positions and velocities of the observatory (wrt the Geocenter) and the center of the Earth (referenced to the SSB) for each TOA. The JPL solar system ephemeris can be set using the ‘ephem’ parameter. The positions and velocities are set with PosVel class instances which have astropy units.

If the required columns already exist, they will be replaced.

Parameters:
  • ephem (str) – The Solar System ephemeris to use; if not specified, use the default ephemeris for the TOAs object. If specified, replace the TOAs object’s ephem attribute with this value and do the computation.

  • planets (bool) – Whether to compute positions for the Solar System planets. If not specified, use the value stored in self.planets; if specified, set self.planets to this value.

add_vel_ecl(obliquity)[source]

Compute and add a column to self.table with velocities in ecliptic coordinates.

Called in barycentric_radio_freq() in AstrometryEcliptic (astrometry.py) if ssb_obs_vel_ecl column does not already exist. If compute_posvels() called again for a TOAs object (aka TOAs modified), deletes this column so that this function will be called again and velocities will be calculated with updated TOAs.

update_mjd_float()[source]

Update the mjd_float column from the mjd column

update_all_times(tdb_method='default')[source]

Update the various derived time columns

Updates:

  • mjd_float

  • tdb

  • tdbld

  • ssb_obs_pos

  • ssb_obs_vel

  • obs_sun_pos

  • obs_*_pos if planets is True

  • ssb_obs_vel_ecl if obliquity is not None

Parameters:

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

merge(t, *args, strict=False)[source]

Merge TOAs instances into the existing object

In order for a merge to work, each TOAs instance needs to have been created using the same Solar System Ephemeris (EPHEM), the same reference timescale (i.e. CLOCK), the same value of .planets (i.e. whether planetary PosVel columns are in the tables or not), and the same value of the obliquity.

If pulse_number, [ssb_obs_pos, ssb_obs_vel, obs_sun_pos], [tdb, tdbld] columns are present in some but not all data try to put them in unless strict=True

Parameters:
  • t (pint.toa.TOAs) –

  • args (pint.toa.TOAs, optional) – Additional TOAs to merge

  • strict (bool, optional) – If strict, do not add in missing columns to facilitate merging