Indexing and selecting data

xarray offers extremely flexible indexing routines that combine the best features of NumPy and pandas for data selection.

The most basic way to access elements of a DataArray object is to use Python’s [] syntax, such as array[i, j], where i and j are both integers. As xarray objects can store coordinates corresponding to each dimension of an array, label-based indexing similar to pandas.DataFrame.loc is also possible. In label-based indexing, the element position i is automatically looked-up from the coordinate values.

Dimensions of xarray objects have names, so you can also lookup the dimensions by name, instead of remembering their positional order.

Thus in total, xarray supports four different kinds of indexing, as described below and summarized in this table:

Dimension lookup

Index lookup

DataArray syntax

Dataset syntax

Positional

By integer

da[:, 0]

not available

Positional

By label

da.loc[:, 'IA']

not available

By name

By integer

da.isel(space=0) or
da[dict(space=0)]

ds.isel(space=0) or
ds[dict(space=0)]

By name

By label

da.sel(space='IA') or
da.loc[dict(space='IA')]

ds.sel(space='IA') or
ds.loc[dict(space='IA')]

More advanced indexing is also possible for all the methods by supplying DataArray objects as indexer. See Vectorized Indexing for the details.

Positional indexing

Indexing a DataArray directly works (mostly) just like it does for numpy arrays, except that the returned object is always another DataArray:

In [1]: da = xr.DataArray(
   ...:     np.random.rand(4, 3),
   ...:     [
   ...:         ("time", pd.date_range("2000-01-01", periods=4)),
   ...:         ("space", ["IA", "IL", "IN"]),
   ...:     ],
   ...: )
   ...: 

In [2]: da[:2]
Out[2]: 
<xarray.DataArray (time: 2, space: 3)>
array([[0.127, 0.967, 0.26 ],
       [0.897, 0.377, 0.336]])
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01 2000-01-02
  * space    (space) <U2 'IA' 'IL' 'IN'

In [3]: da[0, 0]
Out[3]: 
<xarray.DataArray ()>
array(0.127)
Coordinates:
    time     datetime64[ns] 2000-01-01
    space    <U2 'IA'

In [4]: da[:, [2, 1]]
Out[4]: 
<xarray.DataArray (time: 4, space: 2)>
array([[0.26 , 0.967],
       [0.336, 0.377],
       [0.123, 0.84 ],
       [0.448, 0.373]])
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01 2000-01-02 2000-01-03 2000-01-04
  * space    (space) <U2 'IN' 'IL'

Attributes are persisted in all indexing operations.

Warning

Positional indexing deviates from the NumPy when indexing with multiple arrays like da[[0, 1], [0, 1]], as described in Vectorized Indexing.

xarray also supports label-based indexing, just like pandas. Because we use a pandas.Index under the hood, label based indexing is very fast. To do label based indexing, use the loc attribute:

In [5]: da.loc["2000-01-01":"2000-01-02", "IA"]
Out[5]: 
<xarray.DataArray (time: 2)>
array([0.127, 0.897])
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01 2000-01-02
    space    <U2 'IA'

In this example, the selected is a subpart of the array in the range ‘2000-01-01’:’2000-01-02’ along the first coordinate time and with ‘IA’ value from the second coordinate space.

You can perform any of the label indexing operations supported by pandas, including indexing with individual, slices and arrays of labels, as well as indexing with boolean arrays. Like pandas, label based indexing in xarray is inclusive of both the start and stop bounds.

Setting values with label based indexing is also supported:

In [6]: da.loc["2000-01-01", ["IL", "IN"]] = -10

In [7]: da
Out[7]: 
<xarray.DataArray (time: 4, space: 3)>
array([[  0.127, -10.   , -10.   ],
       [  0.897,   0.377,   0.336],
       [  0.451,   0.84 ,   0.123],
       [  0.543,   0.373,   0.448]])
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01 2000-01-02 2000-01-03 2000-01-04
  * space    (space) <U2 'IA' 'IL' 'IN'

Indexing with dimension names

With the dimension names, we do not have to rely on dimension order and can use them explicitly to slice data. There are two ways to do this:

  1. Use a dictionary as the argument for array positional or label based array indexing:

    # index by integer array indices
    In [8]: da[dict(space=0, time=slice(None, 2))]
    Out[8]: 
    <xarray.DataArray (time: 2)>
    array([0.127, 0.897])
    Coordinates:
      * time     (time) datetime64[ns] 2000-01-01 2000-01-02
        space    <U2 'IA'
    
    # index by dimension coordinate labels
    In [9]: da.loc[dict(time=slice("2000-01-01", "2000-01-02"))]
    Out[9]: 
    <xarray.DataArray (time: 2, space: 3)>
    array([[  0.127, -10.   , -10.   ],
           [  0.897,   0.377,   0.336]])
    Coordinates:
      * time     (time) datetime64[ns] 2000-01-01 2000-01-02
      * space    (space) <U2 'IA' 'IL' 'IN'
    
  2. Use the sel() and isel() convenience methods:

    # index by integer array indices
    In [10]: da.isel(space=0, time=slice(None, 2))
    Out[10]: 
    <xarray.DataArray (time: 2)>
    array([0.127, 0.897])
    Coordinates:
      * time     (time) datetime64[ns] 2000-01-01 2000-01-02
        space    <U2 'IA'
    
    # index by dimension coordinate labels
    In [11]: da.sel(time=slice("2000-01-01", "2000-01-02"))
    Out[11]: 
    <xarray.DataArray (time: 2, space: 3)>
    array([[  0.127, -10.   , -10.   ],
           [  0.897,   0.377,   0.336]])
    Coordinates:
      * time     (time) datetime64[ns] 2000-01-01 2000-01-02
      * space    (space) <U2 'IA' 'IL' 'IN'
    

The arguments to these methods can be any objects that could index the array along the dimension given by the keyword, e.g., labels for an individual value, Python slice objects or 1-dimensional arrays.

Note

We would love to be able to do indexing with labeled dimension names inside brackets, but unfortunately, Python does yet not support indexing with keyword arguments like da[space=0]

Nearest neighbor lookups

The label based selection methods sel(), reindex() and reindex_like() all support method and tolerance keyword argument. The method parameter allows for enabling nearest neighbor (inexact) lookups by use of the methods 'pad', 'backfill' or 'nearest':

In [12]: da = xr.DataArray([1, 2, 3], [("x", [0, 1, 2])])

In [13]: da.sel(x=[1.1, 1.9], method="nearest")
Out[13]: 
<xarray.DataArray (x: 2)>
array([2, 3])
Coordinates:
  * x        (x) int64 1 2

In [14]: da.sel(x=0.1, method="backfill")
Out[14]: 
<xarray.DataArray ()>
array(2)
Coordinates:
    x        int64 1

In [15]: da.reindex(x=[0.5, 1, 1.5, 2, 2.5], method="pad")
Out[15]: 
<xarray.DataArray (x: 5)>
array([1, 2, 2, 3, 3])
Coordinates:
  * x        (x) float64 0.5 1.0 1.5 2.0 2.5

Tolerance limits the maximum distance for valid matches with an inexact lookup:

In [16]: da.reindex(x=[1.1, 1.5], method="nearest", tolerance=0.2)
Out[16]: 
<xarray.DataArray (x: 2)>
array([ 2., nan])
Coordinates:
  * x        (x) float64 1.1 1.5

The method parameter is not yet supported if any of the arguments to .sel() is a slice object:

In [17]: da.sel(x=slice(1, 3), method="nearest")
NotImplementedError

However, you don’t need to use method to do inexact slicing. Slicing already returns all values inside the range (inclusive), as long as the index labels are monotonic increasing:

In [18]: da.sel(x=slice(0.9, 3.1))
Out[18]: 
<xarray.DataArray (x: 2)>
array([2, 3])
Coordinates:
  * x        (x) int64 1 2

Indexing axes with monotonic decreasing labels also works, as long as the slice or .loc arguments are also decreasing:

In [19]: reversed_da = da[::-1]

In [20]: reversed_da.loc[3.1:0.9]
Out[20]: 
<xarray.DataArray (x: 2)>
array([3, 2])
Coordinates:
  * x        (x) int64 2 1

Note

If you want to interpolate along coordinates rather than looking up the nearest neighbors, use interp() and interp_like(). See interpolation for the details.

Dataset indexing

We can also use these methods to index all variables in a dataset simultaneously, returning a new dataset:

In [21]: da = xr.DataArray(
   ....:     np.random.rand(4, 3),
   ....:     [
   ....:         ("time", pd.date_range("2000-01-01", periods=4)),
   ....:         ("space", ["IA", "IL", "IN"]),
   ....:     ],
   ....: )
   ....: 

In [22]: ds = da.to_dataset(name="foo")

In [23]: ds.isel(space=[0], time=[0])
Out[23]: 
<xarray.Dataset>
Dimensions:  (space: 1, time: 1)
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01
  * space    (space) <U2 'IA'
Data variables:
    foo      (time, space) float64 0.1294

In [24]: ds.sel(time="2000-01-01")
Out[24]: 
<xarray.Dataset>
Dimensions:  (space: 3)
Coordinates:
    time     datetime64[ns] 2000-01-01
  * space    (space) <U2 'IA' 'IL' 'IN'
Data variables:
    foo      (space) float64 0.1294 0.8599 0.8204

Positional indexing on a dataset is not supported because the ordering of dimensions in a dataset is somewhat ambiguous (it can vary between different arrays). However, you can do normal indexing with dimension names:

In [25]: ds[dict(space=[0], time=[0])]
Out[25]: 
<xarray.Dataset>
Dimensions:  (space: 1, time: 1)
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01
  * space    (space) <U2 'IA'
Data variables:
    foo      (time, space) float64 0.1294

In [26]: ds.loc[dict(time="2000-01-01")]
Out[26]: 
<xarray.Dataset>
Dimensions:  (space: 3)
Coordinates:
    time     datetime64[ns] 2000-01-01
  * space    (space) <U2 'IA' 'IL' 'IN'
Data variables:
    foo      (space) float64 0.1294 0.8599 0.8204

Using indexing to assign values to a subset of dataset (e.g., ds[dict(space=0)] = 1) is not yet supported.

Dropping labels and dimensions

The drop_sel() method returns a new object with the listed index labels along a dimension dropped:

In [27]: ds.drop_sel(space=["IN", "IL"])
Out[27]: 
<xarray.Dataset>
Dimensions:  (space: 1, time: 4)
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01 2000-01-02 2000-01-03 2000-01-04
  * space    (space) <U2 'IA'
Data variables:
    foo      (time, space) float64 0.1294 0.3521 0.5948 0.2355

drop_sel is both a Dataset and DataArray method.

Use drop_dims() to drop a full dimension from a Dataset. Any variables with these dimensions are also dropped:

In [28]: ds.drop_dims("time")
Out[28]: 
<xarray.Dataset>
Dimensions:  (space: 3)
Coordinates:
  * space    (space) <U2 'IA' 'IL' 'IN'
Data variables:
    *empty*

Masking with where

Indexing methods on xarray objects generally return a subset of the original data. However, it is sometimes useful to select an object with the same shape as the original data, but with some elements masked. To do this type of selection in xarray, use where():

In [29]: da = xr.DataArray(np.arange(16).reshape(4, 4), dims=["x", "y"])

In [30]: da.where(da.x + da.y < 4)
Out[30]: 
<xarray.DataArray (x: 4, y: 4)>
array([[ 0.,  1.,  2.,  3.],
       [ 4.,  5.,  6., nan],
       [ 8.,  9., nan, nan],
       [12., nan, nan, nan]])
Dimensions without coordinates: x, y

This is particularly useful for ragged indexing of multi-dimensional data, e.g., to apply a 2D mask to an image. Note that where follows all the usual xarray broadcasting and alignment rules for binary operations (e.g., +) between the object being indexed and the condition, as described in Computation:

In [31]: da.where(da.y < 2)
Out[31]: 
<xarray.DataArray (x: 4, y: 4)>
array([[ 0.,  1., nan, nan],
       [ 4.,  5., nan, nan],
       [ 8.,  9., nan, nan],
       [12., 13., nan, nan]])
Dimensions without coordinates: x, y

By default where maintains the original size of the data. For cases where the selected data size is much smaller than the original data, use of the option drop=True clips coordinate elements that are fully masked:

In [32]: da.where(da.y < 2, drop=True)
Out[32]: 
<xarray.DataArray (x: 4, y: 2)>
array([[ 0.,  1.],
       [ 4.,  5.],
       [ 8.,  9.],
       [12., 13.]])
Dimensions without coordinates: x, y

Selecting values with isin

To check whether elements of an xarray object contain a single object, you can compare with the equality operator == (e.g., arr == 3). To check multiple values, use isin():

In [33]: da = xr.DataArray([1, 2, 3, 4, 5], dims=["x"])

In [34]: da.isin([2, 4])
Out[34]: 
<xarray.DataArray (x: 5)>
array([False,  True, False,  True, False])
Dimensions without coordinates: x

isin() works particularly well with where() to support indexing by arrays that are not already labels of an array:

In [35]: lookup = xr.DataArray([-1, -2, -3, -4, -5], dims=["x"])

In [36]: da.where(lookup.isin([-2, -4]), drop=True)
Out[36]: 
<xarray.DataArray (x: 2)>
array([2., 4.])
Dimensions without coordinates: x

However, some caution is in order: when done repeatedly, this type of indexing is significantly slower than using sel().

Vectorized Indexing

Like numpy and pandas, xarray supports indexing many array elements at once in a vectorized manner.

If you only provide integers, slices, or unlabeled arrays (array without dimension names, such as np.ndarray, list, but not DataArray() or Variable()) indexing can be understood as orthogonally. Each indexer component selects independently along the corresponding dimension, similar to how vector indexing works in Fortran or MATLAB, or after using the numpy.ix_() helper:

In [37]: da = xr.DataArray(
   ....:     np.arange(12).reshape((3, 4)),
   ....:     dims=["x", "y"],
   ....:     coords={"x": [0, 1, 2], "y": ["a", "b", "c", "d"]},
   ....: )
   ....: 

In [38]: da
Out[38]: 
<xarray.DataArray (x: 3, y: 4)>
array([[ 0,  1,  2,  3],
       [ 4,  5,  6,  7],
       [ 8,  9, 10, 11]])
Coordinates:
  * x        (x) int64 0 1 2
  * y        (y) <U1 'a' 'b' 'c' 'd'

In [39]: da[[0, 2, 2], [1, 3]]
Out[39]: 
<xarray.DataArray (x: 3, y: 2)>
array([[ 1,  3],
       [ 9, 11],
       [ 9, 11]])
Coordinates:
  * x        (x) int64 0 2 2
  * y        (y) <U1 'b' 'd'

For more flexibility, you can supply DataArray() objects as indexers. Dimensions on resultant arrays are given by the ordered union of the indexers’ dimensions:

In [40]: ind_x = xr.DataArray([0, 1], dims=["x"])

In [41]: ind_y = xr.DataArray([0, 1], dims=["y"])

In [42]: da[ind_x, ind_y]  # orthogonal indexing
Out[42]: 
<xarray.DataArray (x: 2, y: 2)>
array([[0, 1],
       [4, 5]])
Coordinates:
  * x        (x) int64 0 1
  * y        (y) <U1 'a' 'b'

In [43]: da[ind_x, ind_x]  # vectorized indexing
Out[43]: 
<xarray.DataArray (x: 2)>
array([0, 5])
Coordinates:
  * x        (x) int64 0 1
    y        (x) <U1 'a' 'b'

Slices or sequences/arrays without named-dimensions are treated as if they have the same dimension which is indexed along:

# Because [0, 1] is used to index along dimension 'x',
# it is assumed to have dimension 'x'
In [44]: da[[0, 1], ind_x]
Out[44]: 
<xarray.DataArray (x: 2)>
array([0, 5])
Coordinates:
  * x        (x) int64 0 1
    y        (x) <U1 'a' 'b'

Furthermore, you can use multi-dimensional DataArray() as indexers, where the resultant array dimension is also determined by indexers’ dimension:

In [45]: ind = xr.DataArray([[0, 1], [0, 1]], dims=["a", "b"])

In [46]: da[ind]
Out[46]: 
<xarray.DataArray (a: 2, b: 2, y: 4)>
array([[[0, 1, 2, 3],
        [4, 5, 6, 7]],

       [[0, 1, 2, 3],
        [4, 5, 6, 7]]])
Coordinates:
    x        (a, b) int64 0 1 0 1
  * y        (y) <U1 'a' 'b' 'c' 'd'
Dimensions without coordinates: a, b

Similar to how NumPy’s advanced indexing works, vectorized indexing for xarray is based on our broadcasting rules. See Indexing rules for the complete specification.

Vectorized indexing also works with isel, loc, and sel:

In [47]: ind = xr.DataArray([[0, 1], [0, 1]], dims=["a", "b"])

In [48]: da.isel(y=ind)  # same as da[:, ind]
Out[48]: 
<xarray.DataArray (x: 3, a: 2, b: 2)>
array([[[0, 1],
        [0, 1]],

       [[4, 5],
        [4, 5]],

       [[8, 9],
        [8, 9]]])
Coordinates:
  * x        (x) int64 0 1 2
    y        (a, b) object 'a' 'b' 'a' 'b'
Dimensions without coordinates: a, b

In [49]: ind = xr.DataArray([["a", "b"], ["b", "a"]], dims=["a", "b"])

In [50]: da.loc[:, ind]  # same as da.sel(y=ind)
Out[50]: 
<xarray.DataArray (x: 3, a: 2, b: 2)>
array([[[0, 1],
        [1, 0]],

       [[4, 5],
        [5, 4]],

       [[8, 9],
        [9, 8]]])
Coordinates:
  * x        (x) int64 0 1 2
    y        (a, b) object 'a' 'b' 'b' 'a'
Dimensions without coordinates: a, b

These methods may also be applied to Dataset objects

In [51]: ds = da.to_dataset(name="bar")

In [52]: ds.isel(x=xr.DataArray([0, 1, 2], dims=["points"]))
Out[52]: 
<xarray.Dataset>
Dimensions:  (points: 3, y: 4)
Coordinates:
    x        (points) int64 0 1 2
  * y        (y) <U1 'a' 'b' 'c' 'd'
Dimensions without coordinates: points
Data variables:
    bar      (points, y) int64 0 1 2 3 4 5 6 7 8 9 10 11

Vectorized indexing may be used to extract information from the nearest grid cells of interest, for example, the nearest climate model grid cells to a collection specified weather station latitudes and longitudes.

In [53]: ds = xr.tutorial.open_dataset("air_temperature")

# Define target latitude and longitude (where weather stations might be)
In [54]: target_lon = xr.DataArray([200, 201, 202, 205], dims="points")

In [55]: target_lat = xr.DataArray([31, 41, 42, 42], dims="points")

# Retrieve data at the grid cells nearest to the target latitudes and longitudes
In [56]: da = ds["air"].sel(lon=target_lon, lat=target_lat, method="nearest")

In [57]: da
Out[57]: 
<xarray.DataArray 'air' (time: 2920, points: 4)>
array([[293.1    , 284.6    , 283.19998, 282.6    ],
       [293.19998, 283.29   , 282.79   , 283.5    ],
       [292.4    , 282.     , 280.79   , 282.4    ],
       ...,
       [288.88998, 282.49   , 281.29   , 280.99   ],
       [288.29   , 282.09   , 280.38998, 280.59   ],
       [289.49   , 282.09   , 280.59   , 280.99   ]], dtype=float32)
Coordinates:
    lat      (points) float32 30.0 40.0 42.5 42.5
    lon      (points) float32 200.0 200.0 202.5 205.0
  * time     (time) datetime64[ns] 2013-01-01 ... 2014-12-31T18:00:00
Dimensions without coordinates: points
Attributes:
    long_name:     4xDaily Air temperature at sigma level 995
    units:         degK
    precision:     2
    GRIB_id:       11
    GRIB_name:     TMP
    var_desc:      Air temperature
    dataset:       NMC Reanalysis
    level_desc:    Surface
    statistic:     Individual Obs
    parent_stat:   Other
    actual_range:  [185.16 322.1 ]

Tip

If you are lazily loading your data from disk, not every form of vectorized indexing is supported (or if supported, may not be supported efficiently). You may find increased performance by loading your data into memory first, e.g., with load().

Note

If an indexer is a DataArray(), its coordinates should not conflict with the selected subpart of the target array (except for the explicitly indexed dimensions with .loc/.sel). Otherwise, IndexError will be raised.

Assigning values with indexing

To select and assign values to a portion of a DataArray() you can use indexing with .loc :

In [58]: ds = xr.tutorial.open_dataset("air_temperature")

# add an empty 2D dataarray
In [59]: ds["empty"] = xr.full_like(ds.air.mean("time"), fill_value=0)

# modify one grid point using loc()
In [60]: ds["empty"].loc[dict(lon=260, lat=30)] = 100

# modify a 2D region using loc()
In [61]: lc = ds.coords["lon"]

In [62]: la = ds.coords["lat"]

In [63]: ds["empty"].loc[
   ....:     dict(lon=lc[(lc > 220) & (lc < 260)], lat=la[(la > 20) & (la < 60)])
   ....: ] = 100
   ....: 

or where():

# modify one grid point using xr.where()
In [64]: ds["empty"] = xr.where(
   ....:     (ds.coords["lat"] == 20) & (ds.coords["lon"] == 260), 100, ds["empty"]
   ....: )
   ....: 

# or modify a 2D region using xr.where()
In [65]: mask = (
   ....:     (ds.coords["lat"] > 20)
   ....:     & (ds.coords["lat"] < 60)
   ....:     & (ds.coords["lon"] > 220)
   ....:     & (ds.coords["lon"] < 260)
   ....: )
   ....: 

In [66]: ds["empty"] = xr.where(mask, 100, ds["empty"])

Vectorized indexing can also be used to assign values to xarray object.

In [67]: da = xr.DataArray(
   ....:     np.arange(12).reshape((3, 4)),
   ....:     dims=["x", "y"],
   ....:     coords={"x": [0, 1, 2], "y": ["a", "b", "c", "d"]},
   ....: )
   ....: 

In [68]: da
Out[68]: 
<xarray.DataArray (x: 3, y: 4)>
array([[ 0,  1,  2,  3],
       [ 4,  5,  6,  7],
       [ 8,  9, 10, 11]])
Coordinates:
  * x        (x) int64 0 1 2
  * y        (y) <U1 'a' 'b' 'c' 'd'

In [69]: da[0] = -1  # assignment with broadcasting

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

In [71]: ind_x = xr.DataArray([0, 1], dims=["x"])

In [72]: ind_y = xr.DataArray([0, 1], dims=["y"])

In [73]: da[ind_x, ind_y] = -2  # assign -2 to (ix, iy) = (0, 0) and (1, 1)

In [74]: da
Out[74]: 
<xarray.DataArray (x: 3, y: 4)>
array([[-2, -2, -1, -1],
       [-2, -2,  6,  7],
       [ 8,  9, 10, 11]])
Coordinates:
  * x        (x) int64 0 1 2
  * y        (y) <U1 'a' 'b' 'c' 'd'

In [75]: da[ind_x, ind_y] += 100  # increment is also possible

In [76]: da
Out[76]: 
<xarray.DataArray (x: 3, y: 4)>
array([[98, 98, -1, -1],
       [98, 98,  6,  7],
       [ 8,  9, 10, 11]])
Coordinates:
  * x        (x) int64 0 1 2
  * y        (y) <U1 'a' 'b' 'c' 'd'

Like numpy.ndarray, value assignment sometimes works differently from what one may expect.

In [77]: da = xr.DataArray([0, 1, 2, 3], dims=["x"])

In [78]: ind = xr.DataArray([0, 0, 0], dims=["x"])

In [79]: da[ind] -= 1

In [80]: da
Out[80]: 
<xarray.DataArray (x: 4)>
array([-1,  1,  2,  3])
Dimensions without coordinates: x

Where the 0th element will be subtracted 1 only once. This is because v[0] = v[0] - 1 is called three times, rather than v[0] = v[0] - 1 - 1 - 1. See Assigning values to indexed arrays for the details.

Note

Dask array does not support value assignment (see Parallel computing with Dask for the details).

Note

Coordinates in both the left- and right-hand-side arrays should not conflict with each other. Otherwise, IndexError will be raised.

Warning

Do not try to assign values when using any of the indexing methods isel or sel:

# DO NOT do this
da.isel(space=0) = 0

Assigning values with the chained indexing using .sel or .isel fails silently.

In [81]: da = xr.DataArray([0, 1, 2, 3], dims=["x"])

# DO NOT do this
In [82]: da.isel(x=[0, 1, 2])[1] = -1

In [83]: da
Out[83]: 
<xarray.DataArray (x: 4)>
array([0, 1, 2, 3])
Dimensions without coordinates: x

More advanced indexing

The use of DataArray() objects as indexers enables very flexible indexing. The following is an example of the pointwise indexing:

In [84]: da = xr.DataArray(np.arange(56).reshape((7, 8)), dims=["x", "y"])

In [85]: da
Out[85]: 
<xarray.DataArray (x: 7, y: 8)>
array([[ 0,  1,  2, ...,  5,  6,  7],
       [ 8,  9, 10, ..., 13, 14, 15],
       [16, 17, 18, ..., 21, 22, 23],
       ...,
       [32, 33, 34, ..., 37, 38, 39],
       [40, 41, 42, ..., 45, 46, 47],
       [48, 49, 50, ..., 53, 54, 55]])
Dimensions without coordinates: x, y

In [86]: da.isel(x=xr.DataArray([0, 1, 6], dims="z"), y=xr.DataArray([0, 1, 0], dims="z"))
Out[86]: 
<xarray.DataArray (z: 3)>
array([ 0,  9, 48])
Dimensions without coordinates: z

where three elements at (ix, iy) = ((0, 0), (1, 1), (6, 0)) are selected and mapped along a new dimension z.

If you want to add a coordinate to the new dimension z, you can supply a DataArray with a coordinate,

In [87]: da.isel(
   ....:     x=xr.DataArray([0, 1, 6], dims="z", coords={"z": ["a", "b", "c"]}),
   ....:     y=xr.DataArray([0, 1, 0], dims="z"),
   ....: )
   ....: 
Out[87]: 
<xarray.DataArray (z: 3)>
array([ 0,  9, 48])
Coordinates:
  * z        (z) <U1 'a' 'b' 'c'

Analogously, label-based pointwise-indexing is also possible by the .sel method:

In [88]: da = xr.DataArray(
   ....:     np.random.rand(4, 3),
   ....:     [
   ....:         ("time", pd.date_range("2000-01-01", periods=4)),
   ....:         ("space", ["IA", "IL", "IN"]),
   ....:     ],
   ....: )
   ....: 

In [89]: times = xr.DataArray(
   ....:     pd.to_datetime(["2000-01-03", "2000-01-02", "2000-01-01"]), dims="new_time"
   ....: )
   ....: 

In [90]: da.sel(space=xr.DataArray(["IA", "IL", "IN"], dims=["new_time"]), time=times)
Out[90]: 
<xarray.DataArray (new_time: 3)>
array([0.92, 0.34, 0.59])
Coordinates:
    time      (new_time) datetime64[ns] 2000-01-03 2000-01-02 2000-01-01
    space     (new_time) <U2 'IA' 'IL' 'IN'
  * new_time  (new_time) datetime64[ns] 2000-01-03 2000-01-02 2000-01-01

Align and reindex

xarray’s reindex, reindex_like and align impose a DataArray or Dataset onto a new set of coordinates corresponding to dimensions. The original values are subset to the index labels still found in the new labels, and values corresponding to new labels not found in the original object are in-filled with NaN.

xarray operations that combine multiple objects generally automatically align their arguments to share the same indexes. However, manual alignment can be useful for greater control and for increased performance.

To reindex a particular dimension, use reindex():

In [91]: da.reindex(space=["IA", "CA"])
Out[91]: 
<xarray.DataArray (time: 4, space: 2)>
array([[0.574,   nan],
       [0.245,   nan],
       [0.92 ,   nan],
       [0.754,   nan]])
Coordinates:
  * space    (space) <U2 'IA' 'CA'
  * time     (time) datetime64[ns] 2000-01-01 2000-01-02 2000-01-03 2000-01-04

The reindex_like() method is a useful shortcut. To demonstrate, we will make a subset DataArray with new values:

In [92]: foo = da.rename("foo")

In [93]: baz = (10 * da[:2, :2]).rename("baz")

In [94]: baz
Out[94]: 
<xarray.DataArray 'baz' (time: 2, space: 2)>
array([[5.74 , 0.613],
       [2.453, 3.404]])
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01 2000-01-02
  * space    (space) <U2 'IA' 'IL'

Reindexing foo with baz selects out the first two values along each dimension:

In [95]: foo.reindex_like(baz)
Out[95]: 
<xarray.DataArray 'foo' (time: 2, space: 2)>
array([[0.574, 0.061],
       [0.245, 0.34 ]])
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01 2000-01-02
  * space    (space) object 'IA' 'IL'

The opposite operation asks us to reindex to a larger shape, so we fill in the missing values with NaN:

In [96]: baz.reindex_like(foo)
Out[96]: 
<xarray.DataArray 'baz' (time: 4, space: 3)>
array([[5.74 , 0.613,   nan],
       [2.453, 3.404,   nan],
       [  nan,   nan,   nan],
       [  nan,   nan,   nan]])
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01 2000-01-02 2000-01-03 2000-01-04
  * space    (space) object 'IA' 'IL' 'IN'

The align() function lets us perform more flexible database-like 'inner', 'outer', 'left' and 'right' joins:

In [97]: xr.align(foo, baz, join="inner")
Out[97]: 
(<xarray.DataArray 'foo' (time: 2, space: 2)>
 array([[0.574, 0.061],
        [0.245, 0.34 ]])
 Coordinates:
   * time     (time) datetime64[ns] 2000-01-01 2000-01-02
   * space    (space) <U2 'IA' 'IL',
 <xarray.DataArray 'baz' (time: 2, space: 2)>
 array([[5.74 , 0.613],
        [2.453, 3.404]])
 Coordinates:
   * time     (time) datetime64[ns] 2000-01-01 2000-01-02
   * space    (space) <U2 'IA' 'IL')

In [98]: xr.align(foo, baz, join="outer")
Out[98]: 
(<xarray.DataArray 'foo' (time: 4, space: 3)>
 array([[0.574, 0.061, 0.59 ],
        [0.245, 0.34 , 0.985],
        [0.92 , 0.038, 0.862],
        [0.754, 0.405, 0.344]])
 Coordinates:
   * time     (time) datetime64[ns] 2000-01-01 2000-01-02 2000-01-03 2000-01-04
   * space    (space) <U2 'IA' 'IL' 'IN',
 <xarray.DataArray 'baz' (time: 4, space: 3)>
 array([[5.74 , 0.613,   nan],
        [2.453, 3.404,   nan],
        [  nan,   nan,   nan],
        [  nan,   nan,   nan]])
 Coordinates:
   * time     (time) datetime64[ns] 2000-01-01 2000-01-02 2000-01-03 2000-01-04
   * space    (space) <U2 'IA' 'IL' 'IN')

Both reindex_like and align work interchangeably between DataArray and Dataset objects, and with any number of matching dimension names:

In [99]: ds
Out[99]: 
<xarray.Dataset>
Dimensions:  (lat: 25, lon: 53, time: 2920)
Coordinates:
  * lat      (lat) float32 75.0 72.5 70.0 67.5 65.0 ... 25.0 22.5 20.0 17.5 15.0
  * lon      (lon) float32 200.0 202.5 205.0 207.5 ... 322.5 325.0 327.5 330.0
  * time     (time) datetime64[ns] 2013-01-01 ... 2014-12-31T18:00:00
Data variables:
    air      (time, lat, lon) float32 241.2 242.5 243.5 ... 296.5 296.2 295.7
    empty    (lat, lon) float32 0.0 0.0 0.0 0.0 0.0 0.0 ... 0.0 0.0 0.0 0.0 0.0
Attributes:
    Conventions:  COARDS
    title:        4x daily NMC reanalysis (1948)
    description:  Data is from NMC initialized reanalysis\n(4x/day).  These a...
    platform:     Model
    references:   http://www.esrl.noaa.gov/psd/data/gridded/data.ncep.reanaly...

In [100]: ds.reindex_like(baz)
Out[100]: 
<xarray.Dataset>
Dimensions:  (lat: 25, lon: 53, time: 2)
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01 2000-01-02
  * lat      (lat) float32 75.0 72.5 70.0 67.5 65.0 ... 25.0 22.5 20.0 17.5 15.0
  * lon      (lon) float32 200.0 202.5 205.0 207.5 ... 322.5 325.0 327.5 330.0
Data variables:
    air      (time, lat, lon) float32 nan nan nan nan nan ... nan nan nan nan
    empty    (lat, lon) float32 0.0 0.0 0.0 0.0 0.0 0.0 ... 0.0 0.0 0.0 0.0 0.0
Attributes:
    Conventions:  COARDS
    title:        4x daily NMC reanalysis (1948)
    description:  Data is from NMC initialized reanalysis\n(4x/day).  These a...
    platform:     Model
    references:   http://www.esrl.noaa.gov/psd/data/gridded/data.ncep.reanaly...

In [101]: other = xr.DataArray(["a", "b", "c"], dims="other")

# this is a no-op, because there are no shared dimension names
In [102]: ds.reindex_like(other)
Out[102]: 
<xarray.Dataset>
Dimensions:  (lat: 25, lon: 53, time: 2920)
Coordinates:
  * lat      (lat) float32 75.0 72.5 70.0 67.5 65.0 ... 25.0 22.5 20.0 17.5 15.0
  * lon      (lon) float32 200.0 202.5 205.0 207.5 ... 322.5 325.0 327.5 330.0
  * time     (time) datetime64[ns] 2013-01-01 ... 2014-12-31T18:00:00
Data variables:
    air      (time, lat, lon) float32 241.2 242.5 243.5 ... 296.5 296.2 295.7
    empty    (lat, lon) float32 0.0 0.0 0.0 0.0 0.0 0.0 ... 0.0 0.0 0.0 0.0 0.0
Attributes:
    Conventions:  COARDS
    title:        4x daily NMC reanalysis (1948)
    description:  Data is from NMC initialized reanalysis\n(4x/day).  These a...
    platform:     Model
    references:   http://www.esrl.noaa.gov/psd/data/gridded/data.ncep.reanaly...

Missing coordinate labels

Coordinate labels for each dimension are optional (as of xarray v0.9). Label based indexing with .sel and .loc uses standard positional, integer-based indexing as a fallback for dimensions without a coordinate label:

In [103]: da = xr.DataArray([1, 2, 3], dims="x")

In [104]: da.sel(x=[0, -1])
Out[104]: 
<xarray.DataArray (x: 2)>
array([1, 3])
Dimensions without coordinates: x

Alignment between xarray objects where one or both do not have coordinate labels succeeds only if all dimensions of the same name have the same length. Otherwise, it raises an informative error:

In [105]: xr.align(da, da[:2])
ValueError: arguments without labels along dimension 'x' cannot be aligned because they have different dimension sizes: {2, 3}

Underlying Indexes

xarray uses the pandas.Index internally to perform indexing operations. If you need to access the underlying indexes, they are available through the indexes attribute.

In [106]: da = xr.DataArray(
   .....:     np.random.rand(4, 3),
   .....:     [
   .....:         ("time", pd.date_range("2000-01-01", periods=4)),
   .....:         ("space", ["IA", "IL", "IN"]),
   .....:     ],
   .....: )
   .....: 

In [107]: da
Out[107]: 
<xarray.DataArray (time: 4, space: 3)>
array([[0.171, 0.395, 0.642],
       [0.275, 0.462, 0.871],
       [0.401, 0.611, 0.118],
       [0.702, 0.414, 0.342]])
Coordinates:
  * time     (time) datetime64[ns] 2000-01-01 2000-01-02 2000-01-03 2000-01-04
  * space    (space) <U2 'IA' 'IL' 'IN'

In [108]: da.indexes
Out[108]: 
time: DatetimeIndex(['2000-01-01', '2000-01-02', '2000-01-03', '2000-01-04'], dtype='datetime64[ns]', name='time', freq='D')
space: Index(['IA', 'IL', 'IN'], dtype='object', name='space')

In [109]: da.indexes["time"]
Out[109]: DatetimeIndex(['2000-01-01', '2000-01-02', '2000-01-03', '2000-01-04'], dtype='datetime64[ns]', name='time', freq='D')

Use get_index() to get an index for a dimension, falling back to a default pandas.RangeIndex if it has no coordinate labels:

In [110]: da = xr.DataArray([1, 2, 3], dims="x")

In [111]: da
Out[111]: 
<xarray.DataArray (x: 3)>
array([1, 2, 3])
Dimensions without coordinates: x

In [112]: da.get_index("x")
Out[112]: RangeIndex(start=0, stop=3, step=1, name='x')

Copies vs. Views

Whether array indexing returns a view or a copy of the underlying data depends on the nature of the labels.

For positional (integer) indexing, xarray follows the same rules as NumPy:

  • Positional indexing with only integers and slices returns a view.

  • Positional indexing with arrays or lists returns a copy.

The rules for label based indexing are more complex:

  • Label-based indexing with only slices returns a view.

  • Label-based indexing with arrays returns a copy.

  • Label-based indexing with scalars returns a view or a copy, depending upon if the corresponding positional indexer can be represented as an integer or a slice object. The exact rules are determined by pandas.

Whether data is a copy or a view is more predictable in xarray than in pandas, so unlike pandas, xarray does not produce SettingWithCopy warnings. However, you should still avoid assignment with chained indexing.

Multi-level indexing

Just like pandas, advanced indexing on multi-level indexes is possible with loc and sel. You can slice a multi-index by providing multiple indexers, i.e., a tuple of slices, labels, list of labels, or any selector allowed by pandas:

In [113]: midx = pd.MultiIndex.from_product([list("abc"), [0, 1]], names=("one", "two"))

In [114]: mda = xr.DataArray(np.random.rand(6, 3), [("x", midx), ("y", range(3))])

In [115]: mda
Out[115]: 
<xarray.DataArray (x: 6, y: 3)>
array([[0.596, 0.2  , 0.1  ],
       [0.735, 0.017, 0.481],
       [0.096, 0.497, 0.839],
       [0.897, 0.733, 0.759],
       [0.561, 0.471, 0.139],
       [0.094, 0.942, 0.134]])
Coordinates:
  * x        (x) MultiIndex
  - one      (x) object 'a' 'a' 'b' 'b' 'c' 'c'
  - two      (x) int64 0 1 0 1 0 1
  * y        (y) int64 0 1 2

In [116]: mda.sel(x=(list("ab"), [0]))
Out[116]: 
<xarray.DataArray (x: 2, y: 3)>
array([[0.596, 0.2  , 0.1  ],
       [0.096, 0.497, 0.839]])
Coordinates:
  * x        (x) MultiIndex
  - one      (x) object 'a' 'b'
  - two      (x) int64 0 0
  * y        (y) int64 0 1 2

You can also select multiple elements by providing a list of labels or tuples or a slice of tuples:

In [117]: mda.sel(x=[("a", 0), ("b", 1)])
Out[117]: 
<xarray.DataArray (x: 2, y: 3)>
array([[0.596, 0.2  , 0.1  ],
       [0.897, 0.733, 0.759]])
Coordinates:
  * x        (x) MultiIndex
  - one      (x) object 'a' 'b'
  - two      (x) int64 0 1
  * y        (y) int64 0 1 2

Additionally, xarray supports dictionaries:

In [118]: mda.sel(x={"one": "a", "two": 0})
Out[118]: 
<xarray.DataArray (y: 3)>
array([0.596, 0.2  , 0.1  ])
Coordinates:
    x        object ('a', 0)
  * y        (y) int64 0 1 2

For convenience, sel also accepts multi-index levels directly as keyword arguments:

In [119]: mda.sel(one="a", two=0)
Out[119]: 
<xarray.DataArray (y: 3)>
array([0.596, 0.2  , 0.1  ])
Coordinates:
    x        object ('a', 0)
  * y        (y) int64 0 1 2

Note that using sel it is not possible to mix a dimension indexer with level indexers for that dimension (e.g., mda.sel(x={'one': 'a'}, two=0) will raise a ValueError).

Like pandas, xarray handles partial selection on multi-index (level drop). As shown below, it also renames the dimension / coordinate when the multi-index is reduced to a single index.

In [120]: mda.loc[{"one": "a"}, ...]
Out[120]: 
<xarray.DataArray (two: 2, y: 3)>
array([[0.596, 0.2  , 0.1  ],
       [0.735, 0.017, 0.481]])
Coordinates:
  * two      (two) int64 0 1
  * y        (y) int64 0 1 2

Unlike pandas, xarray does not guess whether you provide index levels or dimensions when using loc in some ambiguous cases. For example, for mda.loc[{'one': 'a', 'two': 0}] and mda.loc['a', 0] xarray always interprets (‘one’, ‘two’) and (‘a’, 0) as the names and labels of the 1st and 2nd dimension, respectively. You must specify all dimensions or use the ellipsis in the loc specifier, e.g. in the example above, mda.loc[{'one': 'a', 'two': 0}, :] or mda.loc[('a', 0), ...].

Indexing rules

Here we describe the full rules xarray uses for vectorized indexing. Note that this is for the purposes of explanation: for the sake of efficiency and to support various backends, the actual implementation is different.

  1. (Only for label based indexing.) Look up positional indexes along each dimension from the corresponding pandas.Index.

  2. A full slice object : is inserted for each dimension without an indexer.

  3. slice objects are converted into arrays, given by np.arange(*slice.indices(...)).

  4. Assume dimension names for array indexers without dimensions, such as np.ndarray and list, from the dimensions to be indexed along. For example, v.isel(x=[0, 1]) is understood as v.isel(x=xr.DataArray([0, 1], dims=['x'])).

  5. For each variable in a Dataset or DataArray (the array and its coordinates):

    1. Broadcast all relevant indexers based on their dimension names (see Broadcasting by dimension name for full details).

    2. Index the underling array by the broadcast indexers, using NumPy’s advanced indexing rules.

  6. If any indexer DataArray has coordinates and no coordinate with the same name exists, attach them to the indexed object.

Note

Only 1-dimensional boolean arrays can be used as indexers.