xcdat.decode_non_cf_time
xcdat.decode_non_cf_time#
- xcdat.decode_non_cf_time(dataset)[source]#
Decodes time coordinates and time bounds with non-CF compliant units.
By default,
xarrayonly supports decoding time with CF compliant units 3. This function enables decoding time with non-CF compliant units.The time coordinates must have a “calendar” attribute set to a CF calendar type supported by
cftime(“noleap”, “360_day”, “365_day”, “366_day”, “gregorian”, “proleptic_gregorian”, “julian”, “all_leap”, or “standard”) and a “units” attribute set to a supported format (“months since …” or “years since …”).The logic for this function:
Extract units and reference date strings from the “units” attribute.
For example with “months since 1800-01-01”, the units are “months” and reference date is “1800-01-01”.
Using the reference date, create a reference
datetimeobject.Starting from the reference
datetimeobject, use the numerically encoded time coordinate values (each representing an offset) to create an array ofcftimeobjects based on the calendar type.Using the array of
cftimeobjects, create a new xr.DataArray of time coordinates to replace the numerically encoded ones.If it exists, create a time bounds DataArray using steps 3 and 4.
- Parameters
dataset (xr.Dataset) – Dataset with numerically encoded time coordinates and time bounds (if they exist). If the time coordinates cannot be decoded then the original dataset is returned.
- Returns
xr.Dataset – Dataset with decoded time coordinates and time bounds (if they exist) as
cftimeobjects.
Notes
Time coordinates are represented by
cftime.datetimeobjects because it is not restricted by thepandas.Timestamprange (years 1678 through 2262). Refer to 4 and 5 for more information on this limitation.References
- 3
https://cfconventions.org/cf-conventions/cf-conventions.html#time-coordinate
- 4
- 5
https://pandas.pydata.org/pandas-docs/stable/user_guide/timeseries.html#timestamp-limitations
Examples
Decode the time coordinates with non-CF units in a Dataset:
>>> from xcdat.dataset import decode_non_cf_time >>> >>> ds.time <xarray.DataArray 'time' (time: 3)> array([0, 1, 2]) Coordinates: * time (time) int64 0 1 2 Attributes: units: years since 2000-01-01 bounds: time_bnds axis: T long_name: time standard_name: time calendar: noleap >>> >>> ds_decoded = decode_non_cf_time(ds) >>> ds_decoded.time <xarray.DataArray 'time' (time: 3)> array([cftime.DatetimeNoLeap(1850, 1, 1, 0, 0, 0, 0, has_year_zero=True), cftime.DatetimeNoLeap(1850, 1, 1, 0, 0, 0, 0, has_year_zero=True), cftime.DatetimeNoLeap(1850, 1, 1, 0, 0, 0, 0, has_year_zero=True)], dtype='object') Coordinates: * time (time) datetime64[ns] 2000-01-01 2001-01-01 2002-01-01 Attributes: units: years since 2000-01-01 bounds: time_bnds axis: T long_name: time standard_name: time calendar: noleap
View time encoding information:
>>> ds_decoded.time.encoding {'source': None, 'dtype': dtype('int64'), 'original_shape': (3,), 'units': 'years since 2000-01-01', 'calendar': 'noleap'}