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 of pandas.Index, xarray supports a subset of the indexing functionality enabled through the standard pandas.DatetimeIndex for dates from non-standard calendars commonly used in climate science or dates using 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 objects will be used to represent times (either in indexes, as a CFTimeIndex, or in data arrays with dtype object) if any 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 be represented with the np.datetime64[ns] data type, enabling the use of a pandas.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 time coordinate with dates from a no-leap calendar and a CFTimeIndex 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 enables creating a CFTimeIndex with regularly-spaced dates. For instance, 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")

Mirroring pandas’ method with the same name, infer_freq() allows one to infer the sampling frequency of a CFTimeIndex or a 1-D DataArray containing cftime objects. It also works transparently with np.datetime64[ns] and np.timedelta64[ns] data.

In [7]: xr.infer_freq(dates)
Out[7]: 'MS'

With strftime() we can also easily generate formatted strings from the datetime values of a CFTimeIndex directly or through the dt() accessor for a DataArray using the same formatting as the standard datetime.strftime convention .

In [8]: dates.strftime("%c")
Out[8]: 
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 [9]: da["time"].dt.strftime("%Y%m%d")
Out[9]: 
<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:

In [10]: da.sel(time="0001")
Out[10]: 
<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 [11]: da.sel(time=slice("0001-05", "0002-02"))
Out[11]: 
<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 case just “year”, “month”, “day”, “hour”, “minute”, “second”, “microsecond”, “season”, “dayofyear”, “dayofweek”, and “days_in_month”):

In [12]: da.time.dt.year
Out[12]: 
<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 [13]: da.time.dt.month
Out[13]: 
<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 [14]: da.time.dt.season
Out[14]: 
<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 [15]: da.time.dt.dayofyear
Out[15]: 
<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 [16]: da.time.dt.dayofweek
Out[16]: 
<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

In [17]: da.time.dt.days_in_month
Out[17]: 
<xarray.DataArray 'days_in_month' (time: 24)>
array([31, 28, 31, 30, 31, 30, 31, 31, 30, 31, 30, 31, 31, 28, 31, 30, 31,
       30, 31, 31, 30, 31, 30, 31])
Coordinates:
  * time     (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00
  • Rounding of datetimes to fixed frequencies via the dt accessor:

In [18]: da.time.dt.ceil("3D")
Out[18]: 
<xarray.DataArray 'ceil' (time: 24)>
array([cftime.DatetimeNoLeap(1, 1, 3, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 2, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 3, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 4, 3, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 5, 3, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 6, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 7, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 8, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 9, 3, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 10, 3, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 11, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 12, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 1, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 2, 3, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 3, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 4, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 5, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 6, 3, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 7, 3, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 8, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 9, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 10, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 11, 3, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 12, 3, 0, 0, 0, 0)], dtype=object)
Coordinates:
  * time     (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00

In [19]: da.time.dt.floor("5D")
Out[19]: 
<xarray.DataArray 'floor' (time: 24)>
array([cftime.DatetimeNoLeap(1, 1, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 1, 31, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 2, 25, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 4, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 5, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 5, 31, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 6, 30, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 7, 30, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 8, 29, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 9, 28, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 10, 28, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 11, 27, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 1, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 1, 31, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 2, 25, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 4, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 5, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 5, 31, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 6, 30, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 7, 30, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 8, 29, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 9, 28, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 10, 28, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 11, 27, 0, 0, 0, 0)], dtype=object)
Coordinates:
  * time     (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00

In [20]: da.time.dt.round("2D")
Out[20]: 
<xarray.DataArray 'round' (time: 24)>
array([cftime.DatetimeNoLeap(1, 1, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 2, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 3, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 3, 31, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 5, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 6, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 7, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 8, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 9, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 10, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 11, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(1, 11, 30, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 1, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 2, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 3, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 4, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 5, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 6, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 6, 30, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 8, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 9, 2, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 9, 30, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 11, 1, 0, 0, 0, 0),
       cftime.DatetimeNoLeap(2, 12, 1, 0, 0, 0, 0)], dtype=object)
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 of the year):

In [21]: da.groupby("time.month").sum()
Out[21]: 
<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
In [22]: da.interp(time=[DatetimeNoLeap(1, 1, 15), DatetimeNoLeap(1, 2, 15)])
Out[22]: 
<xarray.DataArray 'foo' (time: 2)>
array([0.4516129, 1.5      ])
Coordinates:
  * time     (time) object 0001-01-15 00:00:00 0001-02-15 00:00:00
  • Interpolation using datetime strings:

In [23]: da.interp(time=["0001-01-15", "0001-02-15"])
Out[23]: 
<xarray.DataArray 'foo' (time: 2)>
array([0.4516129, 1.5      ])
Coordinates:
  * time     (time) object 0001-01-15 00:00:00 0001-02-15 00:00:00
  • Differentiation:

In [24]: da.differentiate("time")
Out[24]: 
<xarray.DataArray 'foo' (time: 24)>
array([3.73357228e-07, 3.94375523e-07, 3.94375523e-07, 3.79681859e-07,
       3.79681859e-07, 3.79681859e-07, 3.79681859e-07, 3.73357228e-07,
       3.79681859e-07, 3.79681859e-07, 3.79681859e-07, 3.79681859e-07,
       3.73357228e-07, 3.94375523e-07, 3.94375523e-07, 3.79681859e-07,
       3.79681859e-07, 3.79681859e-07, 3.79681859e-07, 3.73357228e-07,
       3.79681859e-07, 3.79681859e-07, 3.79681859e-07, 3.85802469e-07])
Coordinates:
  * time     (time) object 0001-01-01 00:00:00 ... 0002-12-01 00:00:00
  • Serialization:

In [25]: da.to_netcdf("example-no-leap.nc")

In [26]: xr.open_dataset("example-no-leap.nc")
Out[26]: 
<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 [27]: da.resample(time="81T", closed="right", label="right", base=3).mean()
Out[27]: 
<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 from a CFTimeIndex to a pandas.DatetimeIndex, despite the difference in calendar types. The recommended way of doing this is to use the built-in to_datetimeindex() method:

In [28]: modern_times = xr.cftime_range("2000", periods=24, freq="MS", calendar="noleap")

In [29]: da = xr.DataArray(range(24), [("time", modern_times)])

In [30]: da
Out[30]: 
<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 [31]: datetimeindex = da.indexes["time"].to_datetimeindex()

In [32]: da["time"] = datetimeindex

However in this case one should use caution to only perform operations which do not depend on differences between dates (e.g. differentiation, interpolation, or upsampling with resample), as these could introduce subtle and silent errors due to the difference in calendar types between the dates encoded in your data and the dates stored in memory.