pint.observatory.clock_file.ClockFile

class pint.observatory.clock_file.ClockFile(mjd, clock, *, filename=None, friendly_name=None, comments=None, header=None, leading_comment=None, valid_beyond_ends=False)[source]

Bases: object

A clock correction file.

The ClockFile class provides a way to read various formats of clock files. It will provide the clock information from the file as arrays of times and clock correction values via the ClockFile.time and ClockFile.clock properties. The file should be initially read using the ClockFile.read() method, for example:

>>> cf = ClockFile.read(os.getenv('TEMPO')+'/clock/time_gbt.dat')
>>> print cf.time
[ 51909.5  51910.5  51911.5 ...,  57475.5  57476.5  57477.5]
>>> print cf.clock
[-3.14  -3.139 -3.152 ...,  0.179  0.185  0.188] us

Or:

>>> cf = ClockFile.read(os.getenv('TEMPO2')+'/clock/gbt2gps.clk',
                            format='tempo2')
>>> print cf.time
[ 51909.5  51910.5  51911.5 ...,  57411.5  57412.5  57413.5]
>>> print cf.clock
[ -3.14000000e-06  -3.13900000e-06  -3.15200000e-06 ...,   1.80000000e-08
   2.10000000e-08   2.30000000e-08] s

Clock correction file objects preserve the comments in the original file, and can be written in TEMPO or TEMPO2 format.

If you want to provide a new clock correction file format, you should write a function that can read it and return an appropriate ClockFile object. It should take a filename, Path, or file-like object as its first argument, and probably at least a friendly_name keyword argument. You would then add it to ClockFile._formats as something like ClockFile._formats["haiku"] = my_read_haiku_function. Then ClockFile.read(filename, format="haiku", ...) will call my_read_haiku_function(filename, ...).

Parameters:
  • mjd (np.ndarray) – The MJDs at which clock corrections are measured

  • clock (astropy.units.Quantity) – The clock corrections at those MJDs (units of time)

  • comments (list of str or None) – The comments following each clock correction; should match mjd and clock in length. If not provided, a list of empty comments is used. The first line of these strings are normally on the same line as the corresponding clock correction, while any subsequent lines occur between it and the next and should probably start with ``# ``.

  • leading_comment (str) – A comment to put at the top of the file. Lines should probably start with ``# ``.

  • filename (str or None) – If present, a file that can be read to reproduce this data.

  • friendly_name (str or None) – A descriptive file name, in case the filename is long or uninformative. If not provided defaults to the filename.

  • header (str) – A header to include, if output in TEMPO2 format.

Methods

evaluate(t[, limits])

Evaluate the clock corrections at the times t.

export(filename)

Write this clock correction file to the specified location.

last_correction_mjd()

Last MJD for which corrections are available.

merge(clocks, *[, trim])

Compute clock correction information for a combination of clocks.

read(filename[, format])

Read file, selecting an appropriate subclass based on format.

write_tempo2_clock_file(filename[, hdrline, ...])

Write clock corrections as a TEMPO2-format clock correction file.

write_tempo_clock_file(filename, obscode[, ...])

Write clock corrections as a TEMPO-format clock correction file.

Attributes

clock

An astropy.units.Quantity recording the amounts of clock correction.

time

An astropy.time.Time recording the dates of clock corrections.

classmethod read(filename, format='tempo', **kwargs)[source]

Read file, selecting an appropriate subclass based on format.

Any additional keyword arguments are passed to the appropriate reader function.

property time

An astropy.time.Time recording the dates of clock corrections.

property clock

An astropy.units.Quantity recording the amounts of clock correction.

evaluate(t, limits='warn')[source]

Evaluate the clock corrections at the times t.

If values past the end are encountered, check for new clock corrections in the global repository. Delegates the actual computation to the included ClockFile object; anything still not covered is treated according to limits.

Parameters:
  • t (astropy.time.Time) – An array-valued Time object specifying the times at which to evaluate the clock correction.

  • limits ("warn" or "error") – If “error”, raise an exception if times outside the range in the clock file are presented (or if the clock file is empty); if “warn”, extrapolate by returning the value at the nearest endpoint but emit a warning.

Returns:

corrections – The corrections in units of microseconds.

Return type:

astropy.units.Quantity

last_correction_mjd()[source]

Last MJD for which corrections are available.

static merge(clocks, *, trim=True)[source]

Compute clock correction information for a combination of clocks.

This combines the clock correction files to produce a single file that produces clock corrections that are the sum of the original input files. For example, one could combine ao2gps.clk and gps2utc.clk to produce ao2utc.clk.

Note that clock correction files can specify discontinuities by repeating MJDs. This merged object should accurately propagate these discontinuities.

Parameters:

trim (bool) – Whether to trim the resulting clock file to the MJD range covered by both input clock files.

Returns:

A merged ClockFile object.

Return type:

ClockFile

write_tempo_clock_file(filename, obscode, extra_comment=None)[source]

Write clock corrections as a TEMPO-format clock correction file.

Parameters:
  • filename (str or pathlib.Path or file-like) – The destination

  • obscode (str) – The TEMPO observatory code. TEMPO effectively concatenates all its clock corrections and uses this field to determine which observatory the clock corrections are relevant to. TEMPO observatory codes are case-insensitive one-character values agreed upon by convention and occurring in tim files. PINT recognizes these as aliases of observatories for which they are agreed upon, and an Observatory object contains a field that can be used to retrieve this.

  • comments (str) – Additional comments to include. These will be included below the headings and each line should be prefaced with # to avoid conflicting with the clock corrections.

write_tempo2_clock_file(filename, hdrline=None, extra_comment=None)[source]

Write clock corrections as a TEMPO2-format clock correction file.

Parameters:
  • filename (str or pathlib.Path or file-like) – The destination

  • hdrline (str) – The first line of the file. Should start with # and consist of a pair of timescales, like UTC(AO) UTC(GPS) that this clock file transforms between. If no value is provided, the value of self.header is used.

  • extra_comment (str) – Additional comments to include. Lines should probably start with # so they will be interpreted as comments. This field frequently contains details of the origin of the file, or at least the commands used to convert it from its original format.

export(filename)[source]

Write this clock correction file to the specified location.