Weather and climate data
xarray
can leverage metadata that follows the Climate and Forecast (CF) conventions if present. Examples include automatic labelling of plots with descriptive names and units if proper metadata is present (see Plotting) and support for non-standard calendars used in climate science through the cftime
module (see Non-standard calendars and dates outside the Timestamp-valid range). There are also a number of geosciences-focused projects that build on xarray (see Xarray related projects).
CF-compliant coordinate variables
MetPy adds a metpy
accessor that allows accessing coordinates with appropriate CF metadata using generic names x
, y
, vertical
and time
. There is also a cartopy_crs attribute that provides projection information, parsed from the appropriate CF metadata, as a Cartopy projection object. See their documentation for more information.
Non-standard calendars and dates outside the Timestamp-valid range
Through the standalone cftime
library and a custom subclass ofpandas.Index
, xarray supports a subset of the indexingfunctionality enabled through the standard pandas.DatetimeIndex
fordates from non-standard calendars commonly used in climate science or datesusing a standard calendar, but outside the Timestamp-valid range(approximately between years 1678 and 2262).
Note
As of xarray version 0.11, by default, cftime.datetime
objectswill be used to represent times (either in indexes, as aCFTimeIndex
, or in data arrays with dtype object) ifany of the following are true:
The dates are from a non-standard calendar
Any dates are outside the Timestamp-valid range.
Otherwise pandas-compatible dates from a standard calendar will berepresented with the np.datetime64[ns]
data type, enabling the use of apandas.DatetimeIndex
or arrays with dtype np.datetime64[ns]
and their full set of associated features.
For example, you can create a DataArray indexed by a timecoordinate with dates from a no-leap calendar and aCFTimeIndex
will automatically be used:
- In [1]: from itertools import product
- In [2]: from cftime import DatetimeNoLeap
- In [3]: dates = [DatetimeNoLeap(year, month, 1) for year, month in
- ...: product(range(1, 3), range(1, 13))]
- ...:
- In [4]: da = xr.DataArray(np.arange(24), coords=[dates], dims=['time'], name='foo')
xarray also includes a cftime_range()
function, which enablescreating a CFTimeIndex
with regularly-spaced dates. Forinstance, we can create the same dates and DataArray we created above using:
- In [5]: dates = xr.cftime_range(start='0001', periods=24, freq='MS', calendar='noleap')
- In [6]: da = xr.DataArray(np.arange(24), coords=[dates], dims=['time'], name='foo')
With strftime()
we can also easily generate formatted strings fromthe datetime values of a CFTimeIndex
directly or through thedt()
accessor for a DataArray
using the same formatting as the standard datetime.strftime convention .
- In [7]: dates.strftime('%c')
- Out[7]:
- Index(['Tue Jan 1 00:00:00 1', 'Fri Feb 1 00:00:00 1',
- 'Fri Mar 1 00:00:00 1', 'Mon Apr 1 00:00:00 1',
- 'Wed May 1 00:00:00 1', 'Sat Jun 1 00:00:00 1',
- 'Mon Jul 1 00:00:00 1', 'Thu Aug 1 00:00:00 1',
- 'Sun Sep 1 00:00:00 1', 'Tue Oct 1 00:00:00 1',
- 'Fri Nov 1 00:00:00 1', 'Sun Dec 1 00:00:00 1',
- 'Wed Jan 1 00:00:00 2', 'Sat Feb 1 00:00:00 2',
- 'Sat Mar 1 00:00:00 2', 'Tue Apr 1 00:00:00 2',
- 'Thu May 1 00:00:00 2', 'Sun Jun 1 00:00:00 2',
- 'Tue Jul 1 00:00:00 2', 'Fri Aug 1 00:00:00 2',
- 'Mon Sep 1 00:00:00 2', 'Wed Oct 1 00:00:00 2',
- 'Sat Nov 1 00:00:00 2', 'Mon Dec 1 00:00:00 2'],
- dtype='object')
- In [8]: da['time'].dt.strftime('%Y%m%d')
- Out[8]:
- <xarray.DataArray 'strftime' (time: 24)>
- array([' 10101', ' 10201', ' 10301', ' 10401', ' 10501', ' 10601',
- ' 10701', ' 10801', ' 10901', ' 11001', ' 11101', ' 11201',
- ' 20101', ' 20201', ' 20301', ' 20401', ' 20501', ' 20601',
- ' 20701', ' 20801', ' 20901', ' 21001', ' 21101', ' 21201'],
- dtype=object)
- Coordinates:
- * time (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00
For data indexed by a CFTimeIndex
xarray currently supports:
- Partial datetime string indexing using strictly ISO 8601-format partialdatetime strings:
- In [9]: da.sel(time='0001')
- Out[9]:
- <xarray.DataArray 'foo' (time: 12)>
- array([ 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11])
- Coordinates:
- * time (time) object 0001-01-01 00:00:00 ... 0001-12-01 00:00:00
- In [10]: da.sel(time=slice('0001-05', '0002-02'))
- Out[10]:
- <xarray.DataArray 'foo' (time: 10)>
- array([ 4, 5, 6, 7, 8, 9, 10, 11, 12, 13])
- Coordinates:
- * time (time) object 0001-05-01 00:00:00 ... 0002-02-01 00:00:00
- Access of basic datetime components via the
dt
accessor (in this casejust “year”, “month”, “day”, “hour”, “minute”, “second”, “microsecond”,“season”, “dayofyear”, and “dayofweek”):
- In [11]: da.time.dt.year
- Out[11]:
- <xarray.DataArray 'year' (time: 24)>
- array([1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2])
- Coordinates:
- * time (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00
- In [12]: da.time.dt.month
- Out[12]:
- <xarray.DataArray 'month' (time: 24)>
- array([ 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 1, 2, 3, 4, 5, 6,
- 7, 8, 9, 10, 11, 12])
- Coordinates:
- * time (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00
- In [13]: da.time.dt.season
- Out[13]:
- <xarray.DataArray 'season' (time: 24)>
- array(['DJF', 'DJF', 'MAM', 'MAM', 'MAM', 'JJA', 'JJA', 'JJA', 'SON', 'SON',
- 'SON', 'DJF', 'DJF', 'DJF', 'MAM', 'MAM', 'MAM', 'JJA', 'JJA', 'JJA',
- 'SON', 'SON', 'SON', 'DJF'], dtype='<U3')
- Coordinates:
- * time (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00
- In [14]: da.time.dt.dayofyear
- Out[14]:
- <xarray.DataArray 'dayofyear' (time: 24)>
- array([ 1, 32, 60, 91, 121, 152, 182, 213, 244, 274, 305, 335, 1, 32,
- 60, 91, 121, 152, 182, 213, 244, 274, 305, 335])
- Coordinates:
- * time (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00
- In [15]: da.time.dt.dayofweek
- Out[15]:
- <xarray.DataArray 'dayofweek' (time: 24)>
- array([1, 4, 4, 0, 2, 5, 0, 3, 6, 1, 4, 6, 2, 5, 5, 1, 3, 6, 1, 4, 0, 2, 5, 0])
- Coordinates:
- * time (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00
- Group-by operations based on datetime accessor attributes (e.g. by month ofthe year):
- In [16]: da.groupby('time.month').sum()
- Out[16]:
- <xarray.DataArray 'foo' (month: 12)>
- array([12, 14, 16, 18, 20, 22, 24, 26, 28, 30, 32, 34])
- Coordinates:
- * month (month) int64 1 2 3 4 5 6 7 8 9 10 11 12
- Interpolation using
cftime.datetime
objects:
- In [17]: da.interp(time=[DatetimeNoLeap(1, 1, 15), DatetimeNoLeap(1, 2, 15)])
- Out[17]:
- <xarray.DataArray 'foo' (time: 2)>
- array([0.451613, 1.5 ])
- Coordinates:
- * time (time) object 0001-01-15 00:00:00 0001-02-15 00:00:00
- Interpolation using datetime strings:
- In [18]: da.interp(time=['0001-01-15', '0001-02-15'])
- Out[18]:
- <xarray.DataArray 'foo' (time: 2)>
- array([0.451613, 1.5 ])
- Coordinates:
- * time (time) object 0001-01-15 00:00:00 0001-02-15 00:00:00
- Differentiation:
- In [19]: da.differentiate('time')
- Out[19]:
- <xarray.DataArray 'foo' (time: 24)>
- array([3.733572e-07, 3.943755e-07, 3.943755e-07, 3.796819e-07, 3.796819e-07,
- 3.796819e-07, 3.796819e-07, 3.733572e-07, 3.796819e-07, 3.796819e-07,
- 3.796819e-07, 3.796819e-07, 3.733572e-07, 3.943755e-07, 3.943755e-07,
- 3.796819e-07, 3.796819e-07, 3.796819e-07, 3.796819e-07, 3.733572e-07,
- 3.796819e-07, 3.796819e-07, 3.796819e-07, 3.858025e-07])
- Coordinates:
- * time (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00
- Serialization:
- In [20]: da.to_netcdf('example-no-leap.nc')
- In [21]: xr.open_dataset('example-no-leap.nc')
- Out[21]:
- <xarray.Dataset>
- Dimensions: (time: 24)
- Coordinates:
- * time (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00
- Data variables:
- foo (time) int64 ...
- And resampling along the time dimension for data indexed by a
CFTimeIndex
:
- In [22]: da.resample(time='81T', closed='right', label='right', base=3).mean()
- Out[22]:
- <xarray.DataArray 'foo' (time: 12428)>
- array([ 0., nan, nan, ..., nan, nan, 23.])
- Coordinates:
- * time (time) object 0001-01-01 00:03:00 ... 0002-12-01 00:30:00
Note
For some use-cases it may still be useful to convert froma CFTimeIndex
to a pandas.DatetimeIndex
,despite the difference in calendar types. The recommended way of doing thisis to use the built-in to_datetimeindex()
method:
- In [23]: modern_times = xr.cftime_range('2000', periods=24, freq='MS', calendar='noleap')
- In [24]: da = xr.DataArray(range(24), [('time', modern_times)])
- In [25]: da
- Out[25]:
- <xarray.DataArray (time: 24)>
- array([ 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17,
- 18, 19, 20, 21, 22, 23])
- Coordinates:
- * time (time) object 2000-01-01 00:00:00 ... 2001-12-01 00:00:00
- In [26]: datetimeindex = da.indexes['time'].to_datetimeindex()
- In [27]: da['time'] = datetimeindex
However in this case one should use caution to only perform operations whichdo not depend on differences between dates (e.g. differentiation,interpolation, or upsampling with resample), as these could introduce subtleand silent errors due to the difference in calendar types between the datesencoded in your data and the dates stored in memory.