astropy:docs

Time

class astropy.time.Time(val, val2=None, format=None, scale=None, precision=None, in_subfmt=None, out_subfmt=None, location=None, copy=False)[source] [edit on github]

Bases: object

Represent and manipulate times and dates for astronomy.

A Time object is initialized with one or more times in the val argument. The input times in val must conform to the specified format and must correspond to the specified time scale. The optional val2 time input should be supplied only for numeric input formats (e.g. JD) where very high precision (better than 64-bit precision) is required.

The allowed values for format can be listed with:

>>> sorted(Time.FORMATS)
['astropy_time', 'byear', 'byear_str', 'cxcsec', 'datetime',
'decimalyear', 'gps', 'iso', 'isot', 'jd', 'jyear', 'jyear_str', 'mjd',
'plot_date', 'unix', 'yday']
Parameters:

val : sequence, str, number, or Time object

Value(s) to initialize the time or times.

val2 : sequence, str, or number; optional

Value(s) to initialize the time or times.

format : str, optional

Format of input value(s)

scale : str, optional

Time scale of input value(s), must be one of the following: (‘tai’, ‘tcb’, ‘tcg’, ‘tdb’, ‘tt’, ‘ut1’, ‘utc’)

precision : int, optional

Digits of precision in string representation of time

in_subfmt : str, optional

Subformat for inputting string times

out_subfmt : str, optional

Subformat for outputting string times

location : EarthLocation or tuple, optional

If given as an tuple, it should be able to initialize an an EarthLocation instance, i.e., either contain 3 items with units of length for geocentric coordinates, or contain a longitude, latitude, and an optional height for geodetic coordinates. Can be a single location, or one for each input time.

copy : bool, optional

Make a copy of the input values

Attributes Summary

FORMATS Dict of time formats
SCALES List of time scales
delta_tdb_tt TDB - TT time scale offset
delta_ut1_utc UT1 - UTC time scale offset
format Time format
in_subfmt Unix wildcard pattern to select subformats for parsing string input times.
isscalar
jd1 First of the two doubles that internally store time value(s) in JD.
jd2 Second of the two doubles that internally store time value(s) in JD.
out_subfmt Unix wildcard pattern to select subformats for outputting times.
precision Decimal precision when outputting seconds as floating point (int value between 0 and 9 inclusive).
scale Time scale
shape
size
value Time value(s) in current format

Methods Summary

copy([format]) Return a fully independent copy the Time object, optionally changing the format.
get_delta_ut1_utc([iers_table, return_status]) Find UT1 - UTC differences by interpolating in IERS Table.
now() Creates a new object corresponding to the instant in time this method is called.
replicate([format, copy]) Return a replica of the Time object, optionally changing the format.
sidereal_time(kind[, longitude, model]) Calculate sidereal time.

Attributes Documentation

FORMATS = {u'mjd': <class 'astropy.time.core.TimeMJD'>, u'cxcsec': <class 'astropy.time.core.TimeCxcSec'>, u'jyear': <class 'astropy.time.core.TimeJulianEpoch'>, u'plot_date': <class 'astropy.time.core.TimePlotDate'>, u'byear': <class 'astropy.time.core.TimeBesselianEpoch'>, u'astropy_time': <class 'astropy.time.core.TimeAstropyTime'>, u'decimalyear': <class 'astropy.time.core.TimeDecimalYear'>, u'datetime': <class 'astropy.time.core.TimeDatetime'>, u'jyear_str': <class 'astropy.time.core.TimeJulianEpochString'>, u'unix': <class 'astropy.time.core.TimeUnix'>, u'iso': <class 'astropy.time.core.TimeISO'>, u'isot': <class 'astropy.time.core.TimeISOT'>, u'jd': <class 'astropy.time.core.TimeJD'>, u'yday': <class 'astropy.time.core.TimeYearDayTime'>, u'byear_str': <class 'astropy.time.core.TimeBesselianEpochString'>, u'gps': <class 'astropy.time.core.TimeGPS'>}

Dict of time formats

SCALES = (u'tai', u'tcb', u'tcg', u'tdb', u'tt', u'ut1', u'utc')

List of time scales

delta_tdb_tt

TDB - TT time scale offset

delta_ut1_utc

UT1 - UTC time scale offset

format

Time format

in_subfmt

Unix wildcard pattern to select subformats for parsing string input times.

isscalar
jd1

First of the two doubles that internally store time value(s) in JD.

jd2

Second of the two doubles that internally store time value(s) in JD.

out_subfmt

Unix wildcard pattern to select subformats for outputting times.

precision

Decimal precision when outputting seconds as floating point (int value between 0 and 9 inclusive).

scale

Time scale

shape
size
value

Time value(s) in current format

Methods Documentation

copy(format=None)[source] [edit on github]

Return a fully independent copy the Time object, optionally changing the format.

If format is supplied then the time format of the returned Time object will be set accordingly, otherwise it will be unchanged from the original.

In this method a full copy of the internal time arrays will be made. The internal time arrays are normally not changeable by the user so in most cases the replicate() method should be used.

Parameters:

format : str, optional

Time format of the copy.

Returns:

tm : Time object

Copy of this object

get_delta_ut1_utc(iers_table=None, return_status=False)[source] [edit on github]

Find UT1 - UTC differences by interpolating in IERS Table.

Parameters:

iers_table : astropy.utils.iers.IERS table, optional

Table containing UT1-UTC differences from IERS Bulletins A and/or B. If None, use default version (see astropy.utils.iers)

return_status : bool

Whether to return status values. If False (default), iers raises IndexError if any time is out of the range covered by the IERS table.

Returns:

ut1_utc : float or float array

UT1-UTC, interpolated in IERS Table

status : int or int array

Status values (if return_status=`True`):: astropy.utils.iers.FROM_IERS_B astropy.utils.iers.FROM_IERS_A astropy.utils.iers.FROM_IERS_A_PREDICTION astropy.utils.iers.TIME_BEFORE_IERS_RANGE astropy.utils.iers.TIME_BEYOND_IERS_RANGE

Notes

In normal usage, UT1-UTC differences are calculated automatically on the first instance ut1 is needed.

Examples

To check in code whether any times are before the IERS table range:

>>> from astropy.utils.iers import TIME_BEFORE_IERS_RANGE
>>> t = Time(['1961-01-01', '2000-01-01'], scale='utc')
>>> delta, status = t.get_delta_ut1_utc(return_status=True)
>>> status == TIME_BEFORE_IERS_RANGE
array([ True, False], dtype=bool)

To use an updated IERS A bulletin to calculate UT1-UTC (see also astropy.utils.iers):

>>> from astropy.utils.iers import IERS_A, IERS_A_URL
>>> from astropy.utils.data import download_file
>>> t = Time(['1974-01-01', '2000-01-01'], scale='utc')
>>> iers_a_file = download_file(IERS_A_URL,
...                             cache=True)        
Downloading ... [Done]
>>> iers_a = IERS_A.open(iers_a_file)              
>>> t.delta_ut1_utc = t.get_delta_ut1_utc(iers_a)  

The delta_ut1_utc property will be used to calculate t.ut1; raises IndexError if any of the times is out of range.

classmethod now()[source] [edit on github]

Creates a new object corresponding to the instant in time this method is called.

Note

“Now” is determined using the utcnow function, so its accuracy and precision is determined by that function. Generally that means it is set by the accuracy of your system clock.

Returns:

nowtime

A new Time object (or a subclass of Time if this is called from such a subclass) at the current time.

replicate(format=None, copy=False)[source] [edit on github]

Return a replica of the Time object, optionally changing the format.

If format is supplied then the time format of the returned Time object will be set accordingly, otherwise it will be unchanged from the original.

If copy is set to True then a full copy of the internal time arrays will be made. By default the replica will use a reference to the original arrays when possible to save memory. The internal time arrays are normally not changeable by the user so in most cases it should not be necessary to set copy to True.

The convenience method copy() is available in which copy is True by default.

Parameters:

format : str, optional

Time format of the replica.

copy : bool, optional

Return a true copy instead of using references where possible.

Returns:

tm : Time object

Replica of this object

sidereal_time(kind, longitude=None, model=None)[source] [edit on github]

Calculate sidereal time.

Parameters:

kind : str

'mean' or 'apparent', i.e., accounting for precession only, or also for nutation.

longitude : Quantity, str, or None; optional

The longitude on the Earth at which to compute the sidereal time. Can be given as a Quantity with angular units (or an Angle or Longitude), or as a name of an observatory (currently, only 'greenwich' is supported, equivalent to 0 deg). If None (default), the lon attribute of the Time object is used.

model : str or None; optional

Precession (and nutation) model to use. The available ones are: - apparent: [u’IAU1994’, u’IAU2000A’, u’IAU2000B’, u’IAU2006A’] - mean: [u’IAU1982’, u’IAU2000’, u’IAU2006’] If None (default), the last (most recent) one from the appropriate list above is used.

Returns:

sidereal time : Longitude

Sidereal time as a quantity with units of hourangle

Page Contents