TimeDelta#

class astropy.time.TimeDelta(val, val2=None, format=None, scale=None, precision=None, in_subfmt=None, out_subfmt=None, location=None, copy=False)[source]#

Bases: TimeBase

Represent the time difference between two times.

A TimeDelta object is initialized with one or more times in the val argument. The input times in val must conform to the specified format. 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:

>>> list(TimeDelta.FORMATS)
['sec', 'jd', 'datetime', 'quantity_str']

Note that for time differences, the scale can be among three groups: geocentric (‘tai’, ‘tt’, ‘tcg’), barycentric (‘tcb’, ‘tdb’), and rotational (‘ut1’). Within each of these, the scales for time differences are the same. Conversion between geocentric and barycentric is possible, as there is only a scale factor change, but one cannot convert to or from ‘ut1’, as this requires knowledge of the actual times, not just their difference. For a similar reason, ‘utc’ is not a valid scale for a time difference: a UTC day is not always 86400 seconds.

For more information see:

Parameters:
valsequence, ndarray, number, Quantity or TimeDelta object

Value(s) to initialize the time difference(s). Any quantities will be converted appropriately (with care taken to avoid rounding errors for regular time units).

val2sequence, ndarray, number, or Quantity; optional

Additional values, as needed to preserve precision.

formatstr, optional

Format of input value(s). For numerical inputs without units, “jd” is assumed and values are interpreted as days. A deprecation warning is raised in this case. To avoid the warning, either specify the format or add units to the input values.

scalestr, optional

Time scale of input value(s), must be one of the following values: (‘tdb’, ‘tt’, ‘ut1’, ‘tcg’, ‘tcb’, ‘tai’). If not given (or None), the scale is arbitrary; when added or subtracted from a Time instance, it will be used without conversion.

precisionint, optional

Digits of precision in string representation of time

in_subfmtstr, optional

Unix glob to select subformats for parsing input times

out_subfmtstr, optional

Unix glob to select subformat for outputting times

copybool, optional

Make a copy of the input values

Attributes Summary

FORMATS

Dict of time delta formats.

SCALES

List of time delta scales.

T

Return an instance with the data transposed.

cache

Return the cache associated with this instance.

format

Get or set time format.

in_subfmt

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

info

Container for meta information like name, description, format.

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.

mask

masked

ndim

The number of dimensions of the instance and underlying arrays.

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

The shape of the time instances.

size

The size of the object, as calculated from its shape.

unmasked

Get an instance without the mask.

value

Time value(s) in current format.

writeable

Methods Summary

argmax([axis, out])

Return indices of the maximum values along the given axis.

argmin([axis, out])

Return indices of the minimum values along the given axis.

argsort([axis, kind])

Returns the indices that would sort the time array.

copy([format])

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

diagonal(*args, **kwargs)

Return an instance with the specified diagonals.

filled(fill_value)

Get a copy of the underlying data, with masked values filled in.

flatten(*args, **kwargs)

Return a copy with the array collapsed into one dimension.

insert(obj, values[, axis])

Insert values before the given indices in the column and return a new Time or TimeDelta object.

isclose(other[, atol, rtol])

Returns a boolean or boolean array where two TimeDelta objects are element-wise equal within a time tolerance.

max([axis, out, keepdims])

Maximum along a given axis.

mean([axis, dtype, out, keepdims, where])

Mean along a given axis.

min([axis, out, keepdims])

Minimum along a given axis.

ptp([axis, out, keepdims])

Peak to peak (maximum - minimum) along a given axis.

ravel(*args, **kwargs)

Return an instance with the array collapsed into one dimension.

replicate(*args, **kwargs)

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

reshape(*args, **kwargs)

Returns an instance containing the same data with a new shape.

sort([axis])

Return a copy sorted along the specified axis.

squeeze(*args, **kwargs)

Return an instance with single-dimensional shape entries removed.

swapaxes(*args, **kwargs)

Return an instance with the given axes interchanged.

take(indices[, axis, out, mode])

Return a new instance formed from the elements at the given indices.

to(unit[, equivalencies])

Convert to a quantity in the specified unit.

to_datetime()

Convert to datetime.timedelta object.

to_string()

Output a string representation of the Time or TimeDelta object.

to_value(*args, **kwargs)

Get time delta values expressed in specified output format or unit.

transpose(*args, **kwargs)

Return an instance with the data transposed.

Attributes Documentation

FORMATS = {'datetime': <class 'astropy.time.formats.TimeDeltaDatetime'>, 'jd': <class 'astropy.time.formats.TimeDeltaJD'>, 'quantity_str': <class 'astropy.time.formats.TimeDeltaQuantityString'>, 'sec': <class 'astropy.time.formats.TimeDeltaSec'>}#

Dict of time delta formats.

SCALES = ('tai', 'tt', 'tcg', 'tcb', 'tdb', 'ut1', 'local')#

List of time delta scales.

T#

Return an instance with the data transposed.

Parameters are as for T. All internal data are views of the data of the original.

cache#

Return the cache associated with this instance.

format#

Get or set time format.

The format defines the way times are represented when accessed via the .value attribute. By default it is the same as the format used for initializing the Time instance, but it can be set to any other value that could be used for initialization. These can be listed with:

>>> list(Time.FORMATS)
['jd', 'mjd', 'decimalyear', 'unix', 'unix_tai', 'cxcsec', 'gps', 'plot_date',
 'stardate', 'datetime', 'ymdhms', 'iso', 'isot', 'yday', 'datetime64',
 'fits', 'byear', 'jyear', 'byear_str', 'jyear_str']
in_subfmt#

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

info#

Container for meta information like name, description, format. This is required when the object is used as a mixin column within a table, but can be used as a general way to store meta information.

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.

mask#
masked#
ndim#

The number of dimensions of the instance and underlying arrays.

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#

The shape of the time instances.

Like shape, can be set to a new shape by assigning a tuple. Note that if different instances share some but not all underlying data, setting the shape of one instance can make the other instance unusable. Hence, it is strongly recommended to get new, reshaped instances with the reshape method.

Raises:
ValueError

If the new shape has the wrong total number of elements.

AttributeError

If the shape of the jd1, jd2, location, delta_ut1_utc, or delta_tdb_tt attributes cannot be changed without the arrays being copied. For these cases, use the Time.reshape method (which copies any arrays that cannot be reshaped in-place).

size#

The size of the object, as calculated from its shape.

unmasked#

Get an instance without the mask.

Note that while one gets a new instance, the underlying data will be shared.

value#

Time value(s) in current format.

writeable#

Methods Documentation

argmax(axis=None, out=None)#

Return indices of the maximum values along the given axis.

This is similar to argmax(), but adapted to ensure that the full precision given by the two doubles jd1 and jd2 is used. See argmax() for detailed documentation.

argmin(axis=None, out=None)#

Return indices of the minimum values along the given axis.

This is similar to argmin(), but adapted to ensure that the full precision given by the two doubles jd1 and jd2 is used. See argmin() for detailed documentation.

argsort(axis=-1, kind='stable')#

Returns the indices that would sort the time array.

This is similar to argsort(), but adapted to ensure that the full precision given by the two doubles jd1 and jd2 is used, and that corresponding attributes are copied. Internally, it uses lexsort(), and hence no sort method can be chosen.

Parameters:
axisint, optional

Axis along which to sort. Default is -1, which means sort along the last axis.

kind‘stable’, optional

Sorting is done with lexsort() so this argument is ignored, but kept for compatibility with argsort(). The sorting is stable, meaning that the order of equal elements is preserved.

Returns:
indicesndarray

An array of indices that sort the time array.

copy(format=None)#

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:
formatstr, optional

Time format of the copy.

Returns:
tmTime object

Copy of this object

diagonal(*args, **kwargs)#

Return an instance with the specified diagonals.

Parameters are as for diagonal(). All internal data are views of the data of the original.

filled(fill_value)#

Get a copy of the underlying data, with masked values filled in.

Parameters:
fill_valueobject

Value to replace masked values with. Note that if this value is masked

flatten(*args, **kwargs)#

Return a copy with the array collapsed into one dimension.

Parameters are as for flatten().

insert(obj, values, axis=0)#

Insert values before the given indices in the column and return a new Time or TimeDelta object.

The values to be inserted must conform to the rules for in-place setting of Time objects (see Get and set values in the Time documentation).

The API signature matches the np.insert API, but is more limited. The specification of insert index obj must be a single integer, and the axis must be 0 for simple row insertion before the index.

Parameters:
objint

Integer index before which values is inserted.

valuesarray_like

Value(s) to insert. If the type of values is different from that of quantity, values is converted to the matching type.

axisint, optional

Axis along which to insert values. Default is 0, which is the only allowed value and will insert a row.

Returns:
outTime subclass

New time object with inserted value(s)

isclose(other, atol=None, rtol=0.0)[source]#

Returns a boolean or boolean array where two TimeDelta objects are element-wise equal within a time tolerance.

This effectively evaluates the expression below:

abs(self - other) <= atol + rtol * abs(other)
Parameters:
otherQuantity or TimeDelta

Quantity or TimeDelta object for comparison.

atolQuantity or TimeDelta

Absolute tolerance for equality with units of time (e.g. u.s or u.day). Default is one bit in the 128-bit JD time representation, equivalent to about 20 picosecs.

rtolfloat

Relative tolerance for equality

max(axis=None, out=None, keepdims=False)#

Maximum along a given axis.

This is similar to max(), but adapted to ensure that the full precision given by the two doubles jd1 and jd2 is used, and that corresponding attributes are copied.

Note that the out argument is present only for compatibility with np.max; since Time instances are immutable, it is not possible to have an actual out to store the result in.

mean(axis=None, dtype=None, out=None, keepdims=False, *, where=True)#

Mean along a given axis.

This is similar to mean(), but adapted to ensure that the full precision given by the two doubles jd1 and jd2 is used, and that corresponding attributes are copied.

Note that the out argument is present only for compatibility with np.mean; since Time instances are immutable, it is not possible to have an actual out to store the result in.

Similarly, the dtype argument is also present for compatibility only; it has no meaning for Time.

Parameters:
axisNone or int or tuple of int, optional

Axis or axes along which the means are computed. The default is to compute the mean of the flattened array.

dtypeNone

Only present for compatibility with mean(), must be None.

outNone

Only present for compatibility with mean(), must be None.

keepdimsbool, optional

If this is set to True, the axes which are reduced are left in the result as dimensions with size one. With this option, the result will broadcast correctly against the input array.

wherearray_like of bool, optional

Elements to include in the mean. See reduce for details.

Returns:
mTime

A new Time instance containing the mean values

min(axis=None, out=None, keepdims=False)#

Minimum along a given axis.

This is similar to min(), but adapted to ensure that the full precision given by the two doubles jd1 and jd2 is used, and that corresponding attributes are copied.

Note that the out argument is present only for compatibility with np.min; since Time instances are immutable, it is not possible to have an actual out to store the result in.

ptp(axis=None, out=None, keepdims=False)#

Peak to peak (maximum - minimum) along a given axis.

This is similar to ptp(), but adapted to ensure that the full precision given by the two doubles jd1 and jd2 is used.

Note that the out argument is present only for compatibility with ptp; since Time instances are immutable, it is not possible to have an actual out to store the result in.

ravel(*args, **kwargs)#

Return an instance with the array collapsed into one dimension.

Parameters are as for ravel(). Note that it is not always possible to unravel an array without copying the data. If you want an error to be raise if the data is copied, you should should assign shape (-1,) to the shape attribute.

replicate(*args, **kwargs)[source]#

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:
formatstr, optional

Time format of the replica.

copybool, optional

Return a true copy instead of using references where possible.

Returns:
tmTime object

Replica of this object

reshape(*args, **kwargs)#

Returns an instance containing the same data with a new shape.

Parameters are as for reshape(). Note that it is not always possible to change the shape of an array without copying the data (see reshape() documentation). If you want an error to be raise if the data is copied, you should assign the new shape to the shape attribute (note: this may not be implemented for all classes using NDArrayShapeMethods).

sort(axis=-1)#

Return a copy sorted along the specified axis.

This is similar to sort(), but internally uses indexing with lexsort() to ensure that the full precision given by the two doubles jd1 and jd2 is kept, and that corresponding attributes are properly sorted and copied as well.

Parameters:
axisint or None

Axis to be sorted. If None, the flattened array is sorted. By default, sort over the last axis.

squeeze(*args, **kwargs)#

Return an instance with single-dimensional shape entries removed.

Parameters are as for squeeze(). All internal data are views of the data of the original.

swapaxes(*args, **kwargs)#

Return an instance with the given axes interchanged.

Parameters are as for swapaxes(): axis1, axis2. All internal data are views of the data of the original.

take(indices, axis=None, out=None, mode='raise')#

Return a new instance formed from the elements at the given indices.

Parameters are as for take(), except that, obviously, no output array can be given.

to(unit, equivalencies=[])[source]#

Convert to a quantity in the specified unit.

Parameters:
unitastropy:unit-like

The unit to convert to.

equivalencieslist of tuple

A list of equivalence pairs to try if the units are not directly convertible (see Equivalencies). If None, no equivalencies will be applied at all, not even any set globallyq or within a context.

Returns:
quantityQuantity

The quantity in the units specified.

See also

to_value

get the numerical value in a given unit.

to_datetime()[source]#

Convert to datetime.timedelta object.

to_string()#

Output a string representation of the Time or TimeDelta object.

Similar to str(self.value) (which uses numpy array formatting) but array values are evaluated only for the items that actually are output. For large arrays this can be a substantial performance improvement.

Returns:
outstr

String representation of the time values.

to_value(*args, **kwargs)[source]#

Get time delta values expressed in specified output format or unit.

This method is flexible and handles both conversion to a specified TimeDelta format / sub-format AND conversion to a specified unit. If positional argument(s) are provided then the first one is checked to see if it is a valid TimeDelta format, and next it is checked to see if it is a valid unit or unit string.

To convert to a TimeDelta format and optional sub-format the options are:

tm = TimeDelta(1.0 * u.s)
tm.to_value('jd')  # equivalent of tm.jd
tm.to_value('jd', 'decimal')  # convert to 'jd' as a Decimal object
tm.to_value('jd', subfmt='decimal')
tm.to_value(format='jd', subfmt='decimal')

To convert to a unit with optional equivalencies, the options are:

tm.to_value('hr')  # convert to u.hr (hours)
tm.to_value('hr', equivalencies=[])
tm.to_value(unit='hr', equivalencies=[])

The built-in TimeDelta options for format are shown below:

>>> list(TimeDelta.FORMATS)
['sec', 'jd', 'datetime', 'quantity_str']

For the two numerical formats ‘jd’ and ‘sec’, the available subfmt options are: {‘float’, ‘long’, ‘decimal’, ‘str’, ‘bytes’}. Here, ‘long’ uses numpy.longdouble for somewhat enhanced precision (with the enhancement depending on platform), and ‘decimal’ instances of decimal.Decimal for full precision. For the ‘str’ and ‘bytes’ sub-formats, the number of digits is also chosen such that time values are represented accurately. Default: as set by out_subfmt (which by default picks the first available for a given format, i.e., ‘float’).

Parameters:
formatstr, optional

The format in which one wants the TimeDelta values. Default: the current format.

subfmtstr, optional

Possible sub-format in which the values should be given. Default: as set by out_subfmt (which by default picks the first available for a given format, i.e., ‘float’ or ‘date_hms’).

unitUnitBase instance or str, optional

The unit in which the value should be given.

equivalencieslist of tuple

A list of equivalence pairs to try if the units are not directly convertible (see Equivalencies). If None, no equivalencies will be applied at all, not even any set globally or within a context.

Returns:
valuendarray or scalar

The value in the format or units specified.

See also

to

Convert to a Quantity instance in a given unit.

value

The time value in the current format.

transpose(*args, **kwargs)#

Return an instance with the data transposed.

Parameters are as for transpose(). All internal data are views of the data of the original.