Computation

The labels associated with DataArray andDataset objects enables some powerful shortcuts forcomputation, notably including aggregation and broadcasting by dimensionnames.

Basic array math

Arithmetic operations with a single DataArray automatically vectorize (likenumpy) over all array values:

  1. In [1]: arr = xr.DataArray(np.random.RandomState(0).randn(2, 3),
  2. ...: [('x', ['a', 'b']), ('y', [10, 20, 30])])
  3. ...:
  4.  
  5. In [2]: arr - 3
  6. Out[2]:
  7. <xarray.DataArray (x: 2, y: 3)>
  8. array([[-1.235948, -2.599843, -2.021262],
  9. [-0.759107, -1.132442, -3.977278]])
  10. Coordinates:
  11. * x (x) <U1 'a' 'b'
  12. * y (y) int64 10 20 30
  13.  
  14. In [3]: abs(arr)
  15. Out[3]:
  16. <xarray.DataArray (x: 2, y: 3)>
  17. array([[1.764052, 0.400157, 0.978738],
  18. [2.240893, 1.867558, 0.977278]])
  19. Coordinates:
  20. * x (x) <U1 'a' 'b'
  21. * y (y) int64 10 20 30

You can also use any of numpy’s or scipy’s many ufunc functions directly ona DataArray:

  1. In [4]: np.sin(arr)
  2. Out[4]:
  3. <xarray.DataArray (x: 2, y: 3)>
  4. array([[ 0.981384, 0.389563, 0.829794],
  5. [ 0.783762, 0.956288, -0.828978]])
  6. Coordinates:
  7. * x (x) <U1 'a' 'b'
  8. * y (y) int64 10 20 30

Use where() to conditionally switch between values:

  1. In [5]: xr.where(arr > 0, 'positive', 'negative')
  2. Out[5]:
  3. <xarray.DataArray (x: 2, y: 3)>
  4. array([['positive', 'positive', 'positive'],
  5. ['positive', 'positive', 'negative']], dtype='<U8')
  6. Coordinates:
  7. * x (x) <U1 'a' 'b'
  8. * y (y) int64 10 20 30

Use @ to perform matrix multiplication:

  1. In [6]: arr @ arr
  2. Out[6]:
  3. <xarray.DataArray ()>
  4. array(13.694382)

Data arrays also implement many numpy.ndarray methods:

  1. In [7]: arr.round(2)
  2. Out[7]:
  3. <xarray.DataArray (x: 2, y: 3)>
  4. array([[ 1.76, 0.4 , 0.98],
  5. [ 2.24, 1.87, -0.98]])
  6. Coordinates:
  7. * x (x) <U1 'a' 'b'
  8. * y (y) int64 10 20 30
  9.  
  10. In [8]: arr.T
  11. Out[8]:
  12. <xarray.DataArray (y: 3, x: 2)>
  13. array([[ 1.764052, 2.240893],
  14. [ 0.400157, 1.867558],
  15. [ 0.978738, -0.977278]])
  16. Coordinates:
  17. * x (x) <U1 'a' 'b'
  18. * y (y) int64 10 20 30

Missing values

xarray objects borrow the isnull(),notnull(), count(),dropna(), fillna(),ffill(), and bfill()methods for working with missing data from pandas:

  1. In [9]: x = xr.DataArray([0, 1, np.nan, np.nan, 2], dims=['x'])
  2.  
  3. In [10]: x.isnull()
  4. Out[10]:
  5. <xarray.DataArray (x: 5)>
  6. array([False, False, True, True, False])
  7. Dimensions without coordinates: x
  8.  
  9. In [11]: x.notnull()
  10. Out[11]:
  11. <xarray.DataArray (x: 5)>
  12. array([ True, True, False, False, True])
  13. Dimensions without coordinates: x
  14.  
  15. In [12]: x.count()
  16. Out[12]:
  17. <xarray.DataArray ()>
  18. array(3)
  19.  
  20. In [13]: x.dropna(dim='x')
  21. Out[13]:
  22. <xarray.DataArray (x: 3)>
  23. array([0., 1., 2.])
  24. Dimensions without coordinates: x
  25.  
  26. In [14]: x.fillna(-1)
  27. Out[14]:
  28. <xarray.DataArray (x: 5)>
  29. array([ 0., 1., -1., -1., 2.])
  30. Dimensions without coordinates: x
  31.  
  32. In [15]: x.ffill('x')
  33. Out[15]:
  34. <xarray.DataArray (x: 5)>
  35. array([0., 1., 1., 1., 2.])
  36. Dimensions without coordinates: x
  37.  
  38. In [16]: x.bfill('x')
  39. Out[16]:
  40. <xarray.DataArray (x: 5)>
  41. array([0., 1., 2., 2., 2.])
  42. Dimensions without coordinates: x

Like pandas, xarray uses the float value np.nan (not-a-number) to representmissing values.

xarray objects also have an interpolate_na() methodfor filling missing values via 1D interpolation.

  1. In [17]: x = xr.DataArray([0, 1, np.nan, np.nan, 2], dims=['x'],
  2. ....: coords={'xx': xr.Variable('x', [0, 1, 1.1, 1.9, 3])})
  3. ....:
  4.  
  5. In [18]: x.interpolate_na(dim='x', method='linear', use_coordinate='xx')
  6. Out[18]:
  7. <xarray.DataArray (x: 5)>
  8. array([0. , 1. , 1.05, 1.45, 2. ])
  9. Coordinates:
  10. xx (x) float64 0.0 1.0 1.1 1.9 3.0
  11. Dimensions without coordinates: x

Note that xarray slightly diverges from the pandas interpolate syntax byproviding the use_coordinate keyword which facilitates a clear specificationof which values to use as the index in the interpolation.

Aggregation

Aggregation methods have been updated to take a dim argument instead ofaxis. This allows for very intuitive syntax for aggregation methods that areapplied along particular dimension(s):

  1. In [19]: arr.sum(dim='x')
  2. Out[19]:
  3. <xarray.DataArray (y: 3)>
  4. array([4.004946e+00, 2.267715e+00, 1.460104e-03])
  5. Coordinates:
  6. * y (y) int64 10 20 30
  7.  
  8. In [20]: arr.std(['x', 'y'])
  9. Out[20]:
  10. <xarray.DataArray ()>
  11. array(1.090383)
  12.  
  13. In [21]: arr.min()
  14. Out[21]:
  15. <xarray.DataArray ()>
  16. array(-0.977278)

If you need to figure out the axis number for a dimension yourself (say,for wrapping code designed to work with numpy arrays), you can use theget_axis_num() method:

  1. In [22]: arr.get_axis_num('y')
  2. Out[22]: 1

These operations automatically skip missing values, like in pandas:

  1. In [23]: xr.DataArray([1, 2, np.nan, 3]).mean()
  2. Out[23]:
  3. <xarray.DataArray ()>
  4. array(2.)

If desired, you can disable this behavior by invoking the aggregation methodwith skipna=False.

Rolling window operations

DataArray objects include a rolling() method. Thismethod supports rolling window aggregation:

  1. In [24]: arr = xr.DataArray(np.arange(0, 7.5, 0.5).reshape(3, 5),
  2. ....: dims=('x', 'y'))
  3. ....:
  4.  
  5. In [25]: arr
  6. Out[25]:
  7. <xarray.DataArray (x: 3, y: 5)>
  8. array([[0. , 0.5, 1. , 1.5, 2. ],
  9. [2.5, 3. , 3.5, 4. , 4.5],
  10. [5. , 5.5, 6. , 6.5, 7. ]])
  11. Dimensions without coordinates: x, y

rolling() is applied along one dimension using thename of the dimension as a key (e.g. y) and the window size as the value(e.g. 3). We get back a Rolling object:

  1. In [26]: arr.rolling(y=3)
  2. Out[26]: DataArrayRolling [window->3,center->False,dim->y]

Aggregation and summary methods can be applied directly to the Rollingobject:

  1. In [27]: r = arr.rolling(y=3)
  2.  
  3. In [28]: r.reduce(np.std)
  4. Out[28]:
  5. <xarray.DataArray (x: 3, y: 5)>
  6. array([[ nan, nan, 0.408248, 0.408248, 0.408248],
  7. [ nan, nan, 0.408248, 0.408248, 0.408248],
  8. [ nan, nan, 0.408248, 0.408248, 0.408248]])
  9. Dimensions without coordinates: x, y
  10.  
  11. In [29]: r.mean()
  12. Out[29]:
  13. <xarray.DataArray (x: 3, y: 5)>
  14. array([[nan, nan, 0.5, 1. , 1.5],
  15. [nan, nan, 3. , 3.5, 4. ],
  16. [nan, nan, 5.5, 6. , 6.5]])
  17. Dimensions without coordinates: x, y

Aggregation results are assigned the coordinate at the end of each window bydefault, but can be centered by passing center=True when constructing theRolling object:

  1. In [30]: r = arr.rolling(y=3, center=True)
  2.  
  3. In [31]: r.mean()
  4. Out[31]:
  5. <xarray.DataArray (x: 3, y: 5)>
  6. array([[nan, 0.5, 1. , 1.5, nan],
  7. [nan, 3. , 3.5, 4. , nan],
  8. [nan, 5.5, 6. , 6.5, nan]])
  9. Dimensions without coordinates: x, y

As can be seen above, aggregations of windows which overlap the border of thearray produce nans. Settingmin_periods in the call to rollingchanges the minimum number of observations within the window required to havea value when aggregating:

  1. In [32]: r = arr.rolling(y=3, min_periods=2)
  2.  
  3. In [33]: r.mean()
  4. Out[33]:
  5. <xarray.DataArray (x: 3, y: 5)>
  6. array([[ nan, 0.25, 0.5 , 1. , 1.5 ],
  7. [ nan, 2.75, 3. , 3.5 , 4. ],
  8. [ nan, 5.25, 5.5 , 6. , 6.5 ]])
  9. Dimensions without coordinates: x, y
  10.  
  11. In [34]: r = arr.rolling(y=3, center=True, min_periods=2)
  12.  
  13. In [35]: r.mean()
  14. Out[35]:
  15. <xarray.DataArray (x: 3, y: 5)>
  16. array([[0.25, 0.5 , 1. , 1.5 , 1.75],
  17. [2.75, 3. , 3.5 , 4. , 4.25],
  18. [5.25, 5.5 , 6. , 6.5 , 6.75]])
  19. Dimensions without coordinates: x, y

Note that rolling window aggregations are faster when bottleneck is installed.

We can also manually iterate through Rolling objects:

  1. for label, arr_window in r:
  2. # arr_window is a view of x

While rolling provides a simple moving average, DataArray also supportsan exponential moving average with rolling_exp().This is similiar to pandas’ ewm method. numbagg is required.

  1. arr.rolling_exp(y=3).mean()

The rolling_exp method takes a window_type kwarg, which can be 'alpha','com' (for center-of-mass), 'span', and 'halflife'. The default isspan.

Finally, the rolling object has a construct method which returns aview of the original DataArray with the windowed dimension inthe last position.You can use this for more advanced rolling operations such as strided rolling,windowed rolling, convolution, short-time FFT etc.

  1. # rolling with 2-point stride
  2. In [36]: rolling_da = r.construct('window_dim', stride=2)
  3.  
  4. In [37]: rolling_da
  5. Out[37]:
  6. <xarray.DataArray (x: 3, y: 3, window_dim: 3)>
  7. array([[[nan, 0. , 0.5],
  8. [0.5, 1. , 1.5],
  9. [1.5, 2. , nan]],
  10.  
  11. [[nan, 2.5, 3. ],
  12. [3. , 3.5, 4. ],
  13. [4. , 4.5, nan]],
  14.  
  15. [[nan, 5. , 5.5],
  16. [5.5, 6. , 6.5],
  17. [6.5, 7. , nan]]])
  18. Dimensions without coordinates: x, y, window_dim
  19.  
  20. In [38]: rolling_da.mean('window_dim', skipna=False)
  21. Out[38]:
  22. <xarray.DataArray (x: 3, y: 3)>
  23. array([[nan, 1. , nan],
  24. [nan, 3.5, nan],
  25. [nan, 6. , nan]])
  26. Dimensions without coordinates: x, y

Because the DataArray given by r.construct('window_dim') is a viewof the original array, it is memory efficient.You can also use construct to compute a weighted rolling sum:

  1. In [39]: weight = xr.DataArray([0.25, 0.5, 0.25], dims=['window'])
  2.  
  3. In [40]: arr.rolling(y=3).construct('window').dot(weight)
  4. Out[40]:
  5. <xarray.DataArray (x: 3, y: 5)>
  6. array([[nan, nan, 0.5, 1. , 1.5],
  7. [nan, nan, 3. , 3.5, 4. ],
  8. [nan, nan, 5.5, 6. , 6.5]])
  9. Dimensions without coordinates: x, y

Note

numpy’s Nan-aggregation functions such as nansum copy the original array.In xarray, we internally use these functions in our aggregation methods(such as .sum()) if skipna argument is not specified or set to True.This means rolling_da.mean('window_dim') is memory inefficient.To avoid this, use skipna=False as the above example.

Coarsen large arrays

DataArray and Dataset objects include acoarsen() and coarsen()methods. This supports the block aggregation along multiple dimensions,

  1. In [41]: x = np.linspace(0, 10, 300)
  2.  
  3. In [42]: t = pd.date_range('15/12/1999', periods=364)
  4.  
  5. In [43]: da = xr.DataArray(np.sin(x) * np.cos(np.linspace(0, 1, 364)[:, np.newaxis]),
  6. ....: dims=['time', 'x'], coords={'time': t, 'x': x})
  7. ....:
  8.  
  9. In [44]: da
  10. Out[44]:
  11. <xarray.DataArray (time: 364, x: 300)>
  12. array([[ 0. , 0.033439, 0.06684 , ..., -0.486721, -0.51566 , -0.544021],
  13. [ 0. , 0.033438, 0.06684 , ..., -0.486719, -0.515658, -0.544019],
  14. [ 0. , 0.033438, 0.066839, ..., -0.486714, -0.515652, -0.544013],
  15. ...,
  16. [ 0. , 0.018222, 0.036423, ..., -0.265229, -0.280998, -0.296454],
  17. [ 0. , 0.018144, 0.036268, ..., -0.264104, -0.279806, -0.295196],
  18. [ 0. , 0.018067, 0.036114, ..., -0.262977, -0.278612, -0.293936]])
  19. Coordinates:
  20. * time (time) datetime64[ns] 1999-12-15 1999-12-16 ... 2000-12-12
  21. * x (x) float64 0.0 0.03344 0.06689 0.1003 ... 9.9 9.933 9.967 10.0

In order to take a block mean for every 7 days along time dimension andevery 2 points along x dimension,

  1. In [45]: da.coarsen(time=7, x=2).mean()
  2. Out[45]:
  3. <xarray.DataArray (time: 52, x: 150)>
  4. array([[ 0.016718, 0.083499, 0.149906, ..., -0.411988, -0.471957, -0.529814],
  5. [ 0.016713, 0.08347 , 0.149854, ..., -0.411846, -0.471794, -0.529631],
  6. [ 0.016701, 0.08341 , 0.149747, ..., -0.41155 , -0.471455, -0.529251],
  7. ...,
  8. [ 0.009682, 0.048356, 0.086814, ..., -0.238592, -0.273321, -0.306828],
  9. [ 0.009417, 0.047034, 0.084441, ..., -0.232071, -0.265851, -0.298441],
  10. [ 0.009149, 0.045695, 0.082037, ..., -0.225463, -0.258281, -0.289944]])
  11. Coordinates:
  12. * time (time) datetime64[ns] 1999-12-18 1999-12-25 ... 2000-12-09
  13. * x (x) float64 0.01672 0.08361 0.1505 0.2174 ... 9.849 9.916 9.983

coarsen() raises an ValueError if the datalength is not a multiple of the corresponding window size.You can choose boundary='trim' or boundary='pad' options for trimmingthe excess entries or padding nan to insufficient entries,

  1. In [46]: da.coarsen(time=30, x=2, boundary='trim').mean()
  2. Out[46]:
  3. <xarray.DataArray (time: 12, x: 150)>
  4. array([[ 0.016701, 0.083413, 0.149751, ..., -0.411563, -0.471469, -0.529267],
  5. [ 0.016589, 0.082853, 0.148746, ..., -0.4088 , -0.468305, -0.525715],
  6. [ 0.016364, 0.081727, 0.146725, ..., -0.403247, -0.461943, -0.518573],
  7. ...,
  8. [ 0.011838, 0.059126, 0.106149, ..., -0.291732, -0.334196, -0.375165],
  9. [ 0.010824, 0.05406 , 0.097053, ..., -0.266733, -0.305558, -0.343017],
  10. [ 0.009736, 0.048624, 0.087295, ..., -0.239913, -0.274835, -0.308527]])
  11. Coordinates:
  12. * time (time) datetime64[ns] 1999-12-29T12:00:00 ... 2000-11-23T12:00:00
  13. * x (x) float64 0.01672 0.08361 0.1505 0.2174 ... 9.849 9.916 9.983

If you want to apply a specific function to coordinate, you can pass thefunction or method name to coord_func option,

  1. In [47]: da.coarsen(time=7, x=2, coord_func={'time': 'min'}).mean()
  2. Out[47]:
  3. <xarray.DataArray (time: 52, x: 150)>
  4. array([[ 0.016718, 0.083499, 0.149906, ..., -0.411988, -0.471957, -0.529814],
  5. [ 0.016713, 0.08347 , 0.149854, ..., -0.411846, -0.471794, -0.529631],
  6. [ 0.016701, 0.08341 , 0.149747, ..., -0.41155 , -0.471455, -0.529251],
  7. ...,
  8. [ 0.009682, 0.048356, 0.086814, ..., -0.238592, -0.273321, -0.306828],
  9. [ 0.009417, 0.047034, 0.084441, ..., -0.232071, -0.265851, -0.298441],
  10. [ 0.009149, 0.045695, 0.082037, ..., -0.225463, -0.258281, -0.289944]])
  11. Coordinates:
  12. * time (time) datetime64[ns] 1999-12-15 1999-12-22 ... 2000-12-06
  13. * x (x) float64 0.01672 0.08361 0.1505 0.2174 ... 9.849 9.916 9.983

Computation using Coordinates

Xarray objects have some handy methods for the computation with theircoordinates. differentiate() computes derivatives bycentral finite differences using their coordinates,

  1. In [48]: a = xr.DataArray([0, 1, 2, 3], dims=['x'], coords=[[0.1, 0.11, 0.2, 0.3]])
  2.  
  3. In [49]: a
  4. Out[49]:
  5. <xarray.DataArray (x: 4)>
  6. array([0, 1, 2, 3])
  7. Coordinates:
  8. * x (x) float64 0.1 0.11 0.2 0.3
  9.  
  10. In [50]: a.differentiate('x')
  11. Out[50]:
  12. <xarray.DataArray (x: 4)>
  13. array([100. , 91.111111, 10.584795, 10. ])
  14. Coordinates:
  15. * x (x) float64 0.1 0.11 0.2 0.3

This method can be used also for multidimensional arrays,

  1. In [51]: a = xr.DataArray(np.arange(8).reshape(4, 2), dims=['x', 'y'],
  2. ....: coords={'x': [0.1, 0.11, 0.2, 0.3]})
  3. ....:
  4.  
  5. In [52]: a.differentiate('x')
  6. Out[52]:
  7. <xarray.DataArray (x: 4, y: 2)>
  8. array([[200. , 200. ],
  9. [182.222222, 182.222222],
  10. [ 21.169591, 21.169591],
  11. [ 20. , 20. ]])
  12. Coordinates:
  13. * x (x) float64 0.1 0.11 0.2 0.3
  14. Dimensions without coordinates: y

integrate() computes integration based ontrapezoidal rule using their coordinates,

  1. In [53]: a.integrate('x')
  2. Out[53]:
  3. <xarray.DataArray (y: 2)>
  4. array([0.78, 0.98])
  5. Dimensions without coordinates: y

Note

These methods are limited to simple cartesian geometry. Differentiationand integration along multidimensional coordinate are not supported.

Broadcasting by dimension name

DataArray objects are automatically align themselves (“broadcasting” inthe numpy parlance) by dimension name instead of axis order. With xarray, youdo not need to transpose arrays or insert dimensions of length 1 to get arrayoperations to work, as commonly done in numpy with np.reshape() ornp.newaxis.

This is best illustrated by a few examples. Consider two one-dimensionalarrays with different sizes aligned along different dimensions:

  1. In [54]: a = xr.DataArray([1, 2], [('x', ['a', 'b'])])
  2.  
  3. In [55]: a
  4. Out[55]:
  5. <xarray.DataArray (x: 2)>
  6. array([1, 2])
  7. Coordinates:
  8. * x (x) <U1 'a' 'b'
  9.  
  10. In [56]: b = xr.DataArray([-1, -2, -3], [('y', [10, 20, 30])])
  11.  
  12. In [57]: b
  13. Out[57]:
  14. <xarray.DataArray (y: 3)>
  15. array([-1, -2, -3])
  16. Coordinates:
  17. * y (y) int64 10 20 30

With xarray, we can apply binary mathematical operations to these arrays, andtheir dimensions are expanded automatically:

  1. In [58]: a * b
  2. Out[58]:
  3. <xarray.DataArray (x: 2, y: 3)>
  4. array([[-1, -2, -3],
  5. [-2, -4, -6]])
  6. Coordinates:
  7. * x (x) <U1 'a' 'b'
  8. * y (y) int64 10 20 30

Moreover, dimensions are always reordered to the order in which they firstappeared:

  1. In [59]: c = xr.DataArray(np.arange(6).reshape(3, 2), [b['y'], a['x']])
  2.  
  3. In [60]: c
  4. Out[60]:
  5. <xarray.DataArray (y: 3, x: 2)>
  6. array([[0, 1],
  7. [2, 3],
  8. [4, 5]])
  9. Coordinates:
  10. * y (y) int64 10 20 30
  11. * x (x) <U1 'a' 'b'
  12.  
  13. In [61]: a + c
  14. Out[61]:
  15. <xarray.DataArray (x: 2, y: 3)>
  16. array([[1, 3, 5],
  17. [3, 5, 7]])
  18. Coordinates:
  19. * x (x) <U1 'a' 'b'
  20. * y (y) int64 10 20 30

This means, for example, that you always subtract an array from its transpose:

  1. In [62]: c - c.T
  2. Out[62]:
  3. <xarray.DataArray (y: 3, x: 2)>
  4. array([[0, 0],
  5. [0, 0],
  6. [0, 0]])
  7. Coordinates:
  8. * y (y) int64 10 20 30
  9. * x (x) <U1 'a' 'b'

You can explicitly broadcast xarray data structures by using thebroadcast() function:

  1. In [63]: a2, b2 = xr.broadcast(a, b)
  2.  
  3. In [64]: a2
  4. Out[64]:
  5. <xarray.DataArray (x: 2, y: 3)>
  6. array([[1, 1, 1],
  7. [2, 2, 2]])
  8. Coordinates:
  9. * x (x) <U1 'a' 'b'
  10. * y (y) int64 10 20 30
  11.  
  12. In [65]: b2
  13. Out[65]:
  14. <xarray.DataArray (x: 2, y: 3)>
  15. array([[-1, -2, -3],
  16. [-1, -2, -3]])
  17. Coordinates:
  18. * y (y) int64 10 20 30
  19. * x (x) <U1 'a' 'b'

Automatic alignment

xarray enforces alignment between indexCoordinates (that is,coordinates with the same name as a dimension, marked by *) on objects usedin binary operations.

Similarly to pandas, this alignment is automatic for arithmetic on binaryoperations. The default result of a binary operation is by the intersection(not the union) of coordinate labels:

  1. In [66]: arr = xr.DataArray(np.arange(3), [('x', range(3))])
  2.  
  3. In [67]: arr + arr[:-1]
  4. Out[67]:
  5. <xarray.DataArray (x: 2)>
  6. array([0, 2])
  7. Coordinates:
  8. * x (x) int64 0 1

If coordinate values for a dimension are missing on either argument, allmatching dimensions must have the same size:

  1. In [68]: arr + xr.DataArray([1, 2], dims='x')
  2. ValueError: arguments without labels along dimension 'x' cannot be aligned because they have different dimension size(s) {2} than the size of the aligned dimension labels: 3

However, one can explicitly change this default automatic alignment type (“inner”)via set_options() in context manager:

  1. In [69]: with xr.set_options(arithmetic_join="outer"):
  2. ....: arr + arr[:1]
  3. ....:
  4.  
  5. In [70]: arr + arr[:1]
  6. Out[70]:
  7. <xarray.DataArray (x: 1)>
  8. array([0])
  9. Coordinates:
  10. * x (x) int64 0

Before loops or performance critical code, it’s a good idea to align arraysexplicitly (e.g., by putting them in the same Dataset or usingalign()) to avoid the overhead of repeated alignment with eachoperation. See Align and reindex for more details.

Note

There is no automatic alignment between arguments when performing in-placearithmetic operations such as +=. You will need to usemanual alignment. This ensures in-placearithmetic never needs to modify data types.

Coordinates

Although index coordinates are aligned, other coordinates are not, and if theirvalues conflict, they will be dropped. This is necessary, for example, becauseindexing turns 1D coordinates into scalar coordinates:

  1. In [71]: arr[0]
  2. Out[71]:
  3. <xarray.DataArray ()>
  4. array(0)
  5. Coordinates:
  6. x int64 0
  7.  
  8. In [72]: arr[1]
  9. Out[72]:
  10. <xarray.DataArray ()>
  11. array(1)
  12. Coordinates:
  13. x int64 1
  14.  
  15. # notice that the scalar coordinate 'x' is silently dropped
  16. In [73]: arr[1] - arr[0]
  17. Out[73]:
  18. <xarray.DataArray ()>
  19. array(1)

Still, xarray will persist other coordinates in arithmetic, as long as thereare no conflicting values:

  1. # only one argument has the 'x' coordinate
  2. In [74]: arr[0] + 1
  3. Out[74]:
  4. <xarray.DataArray ()>
  5. array(1)
  6. Coordinates:
  7. x int64 0
  8.  
  9. # both arguments have the same 'x' coordinate
  10. In [75]: arr[0] - arr[0]
  11. Out[75]:
  12. <xarray.DataArray ()>
  13. array(0)
  14. Coordinates:
  15. x int64 0

Math with datasets

Datasets support arithmetic operations by automatically looping over all datavariables:

  1. In [76]: ds = xr.Dataset({'x_and_y': (('x', 'y'), np.random.randn(3, 5)),
  2. ....: 'x_only': ('x', np.random.randn(3))},
  3. ....: coords=arr.coords)
  4. ....:
  5.  
  6. In [77]: ds > 0
  7. Out[77]:
  8. <xarray.Dataset>
  9. Dimensions: (x: 3, y: 5)
  10. Coordinates:
  11. * x (x) int64 0 1 2
  12. Dimensions without coordinates: y
  13. Data variables:
  14. x_and_y (x, y) bool True False False False True ... True True False False
  15. x_only (x) bool True False True

Datasets support most of the same methods found on data arrays:

  1. In [78]: ds.mean(dim='x')
  2. Out[78]:
  3. <xarray.Dataset>
  4. Dimensions: (y: 5)
  5. Dimensions without coordinates: y
  6. Data variables:
  7. x_and_y (y) float64 -0.06634 0.3027 -0.6106 -0.9014 -0.644
  8. x_only float64 0.138
  9.  
  10. In [79]: abs(ds)
  11. Out[79]:
  12. <xarray.Dataset>
  13. Dimensions: (x: 3, y: 5)
  14. Coordinates:
  15. * x (x) int64 0 1 2
  16. Dimensions without coordinates: y
  17. Data variables:
  18. x_and_y (x, y) float64 0.4691 0.2829 1.509 1.136 ... 0.7216 0.7068 1.04
  19. x_only (x) float64 0.2719 0.425 0.567

Datasets also support NumPy ufuncs (requires NumPy v1.13 or newer), oralternatively you can use apply() to apply a functionto each variable in a dataset:

  1. In [80]: np.sin(ds)
  2. Out[80]:
  3. <xarray.Dataset>
  4. Dimensions: (x: 3, y: 5)
  5. Coordinates:
  6. * x (x) int64 0 1 2
  7. Dimensions without coordinates: y
  8. Data variables:
  9. x_and_y (x, y) float64 0.4521 -0.2791 -0.9981 ... 0.6606 -0.6494 -0.8622
  10. x_only (x) float64 0.2685 -0.4123 0.5371
  11.  
  12. In [81]: ds.apply(np.sin)
  13. Out[81]:
  14. <xarray.Dataset>
  15. Dimensions: (x: 3, y: 5)
  16. Coordinates:
  17. * x (x) int64 0 1 2
  18. Dimensions without coordinates: y
  19. Data variables:
  20. x_and_y (x, y) float64 0.4521 -0.2791 -0.9981 ... 0.6606 -0.6494 -0.8622
  21. x_only (x) float64 0.2685 -0.4123 0.5371

Datasets also use looping over variables for broadcasting in binaryarithmetic. You can do arithmetic between any DataArray and a dataset:

  1. In [82]: ds + arr
  2. Out[82]:
  3. <xarray.Dataset>
  4. Dimensions: (x: 3, y: 5)
  5. Coordinates:
  6. * x (x) int64 0 1 2
  7. Dimensions without coordinates: y
  8. Data variables:
  9. x_and_y (x, y) float64 0.4691 -0.2829 -1.509 -1.136 ... 2.722 1.293 0.9604
  10. x_only (x) float64 0.2719 0.575 2.567

Arithmetic between two datasets matches data variables of the same name:

  1. In [83]: ds2 = xr.Dataset({'x_and_y': 0, 'x_only': 100})
  2.  
  3. In [84]: ds - ds2
  4. Out[84]:
  5. <xarray.Dataset>
  6. Dimensions: (x: 3, y: 5)
  7. Coordinates:
  8. * x (x) int64 0 1 2
  9. Dimensions without coordinates: y
  10. Data variables:
  11. x_and_y (x, y) float64 0.4691 -0.2829 -1.509 ... 0.7216 -0.7068 -1.04
  12. x_only (x) float64 -99.73 -100.4 -99.43

Similarly to index based alignment, the result has the intersection of allmatching data variables.

Wrapping custom computation

It doesn’t always make sense to do computation directly with xarray objects:

  • In the inner loop of performance limited code, using xarray can addconsiderable overhead compared to using NumPy or native Python types.This is particularly true when working with scalars or small arrays (lessthan ~1e6 elements). Keeping track of labels and ensuring their consistencyadds overhead, and xarray’s core itself is not especially fast, because it’swritten in Python rather than a compiled language like C. Also, xarray’shigh level label-based APIs removes low-level control over how operationsare implemented.

  • Even if speed doesn’t matter, it can be important to wrap existing code, orto support alternative interfaces that don’t use xarray objects.

For these reasons, it is often well-advised to write low-level routines thatwork with NumPy arrays, and to wrap these routines to work with xarray objects.However, adding support for labels on both Dataset andDataArray can be a bit of a chore.

To make this easier, xarray supplies the apply_ufunc() helperfunction, designed for wrapping functions that support broadcasting andvectorization on unlabeled arrays in the style of a NumPyuniversal function (“ufunc” for short).apply_ufunc takes care of everything needed for an idiomatic xarray wrapper,including alignment, broadcasting, looping over Dataset variables (ifneeded), and merging of coordinates. In fact, many internal xarrayfunctions/methods are written using apply_ufunc.

Simple functions that act independently on each value should work withoutany additional arguments:

  1. In [85]: squared_error = lambda x, y: (x - y) ** 2
  2.  
  3. In [86]: arr1 = xr.DataArray([0, 1, 2, 3], dims='x')
  4.  
  5. In [87]: xr.apply_ufunc(squared_error, arr1, 1)
  6. Out[87]:
  7. <xarray.DataArray (x: 4)>
  8. array([1, 0, 1, 4])
  9. Dimensions without coordinates: x

For using more complex operations that consider some array values collectively,it’s important to understand the idea of “core dimensions” from NumPy’sgeneralized ufuncs. Core dimensions are defined as dimensionsthat should not be broadcast over. Usually, they correspond to the fundamentaldimensions over which an operation is defined, e.g., the summed axis innp.sum. A good clue that core dimensions are needed is the presence of anaxis argument on the corresponding NumPy function.

With apply_ufunc, core dimensions are recognized by name, and then moved tothe last dimension of any input arguments before applying the given function.This means that for functions that accept an axis argument, you usually needto set axis=-1. As an example, here is how we would wrapnumpy.linalg.norm() to calculate the vector norm:

  1. def vector_norm(x, dim, ord=None):
  2. return xr.apply_ufunc(np.linalg.norm, x,
  3. input_core_dims=[[dim]],
  4. kwargs={'ord': ord, 'axis': -1})
  1. In [88]: vector_norm(arr1, dim='x')
  2. Out[88]:
  3. <xarray.DataArray ()>
  4. array(3.741657)

Because apply_ufunc follows a standard convention for ufuncs, it playsnicely with tools for building vectorized functions, likenumpy.broadcast_arrays() and numpy.vectorize(). For high performanceneeds, consider using Numba’s vectorize and guvectorize.

In addition to wrapping functions, apply_ufunc can automatically parallelizemany functions when using dask by setting dask='parallelized'. SeeAutomatic parallelization for details.

apply_ufunc() also supports some advanced options forcontrolling alignment of variables and the form of the result. See thedocstring for full details and more examples.