Working with Multidimensional Coordinates

Author: Ryan Abernathey

Many datasets have physical coordinates which differ from theirlogical coordinates. Xarray provides several ways to plot and analyzesuch datasets.

  1. In [1]: import numpy as np
  2.  
  3. In [2]: import pandas as pd
  4.  
  5. In [3]: import xarray as xr
  6.  
  7. In [4]: import netCDF4
  8.  
  9. In [5]: import cartopy.crs as ccrs
  10.  
  11. In [6]: import matplotlib.pyplot as plt

As an example, consider this dataset from thexarray-data repository.

  1. In [7]: ds = xr.tutorial.open_dataset('rasm').load()
  2.  
  3. In [8]: ds
  4. Out[8]:
  5. <xarray.Dataset>
  6. Dimensions: (time: 36, x: 275, y: 205)
  7. Coordinates:
  8. * time (time) object 1980-09-16 12:00:00 ... 1983-08-17 00:00:00
  9. xc (y, x) float64 189.2 189.4 189.6 189.7 ... 17.65 17.4 17.15 16.91
  10. yc (y, x) float64 16.53 16.78 17.02 17.27 ... 28.26 28.01 27.76 27.51
  11. Dimensions without coordinates: x, y
  12. Data variables:
  13. Tair (time, y, x) float64 nan nan nan nan nan ... 29.8 28.66 28.19 28.21
  14. Attributes:
  15. title: /workspace/jhamman/processed/R1002RBRxaaa01a/l...
  16. institution: U.W.
  17. source: RACM R1002RBRxaaa01a
  18. output_frequency: daily
  19. output_mode: averaged
  20. convention: CF-1.4
  21. references: Based on the initial model of Liang et al., 19...
  22. comment: Output from the Variable Infiltration Capacity...
  23. nco_openmp_thread_number: 1
  24. NCO: "4.6.0"
  25. history: Tue Dec 27 14:15:22 2016: ncatted -a dimension...

In this example, the logical coordinates are x and y, whilethe physical coordinates are xc and yc, which represent thelatitudes and longitude of the data.

  1. In [9]: ds.xc.attrs
  2. Out[9]:
  3. OrderedDict([('long_name', 'longitude of grid cell center'),
  4. ('units', 'degrees_east'),
  5. ('bounds', 'xv')])
  6.  
  7. In [10]: ds.yc.attrs
  8. Out[10]:
  9. OrderedDict([('long_name', 'latitude of grid cell center'),
  10. ('units', 'degrees_north'),
  11. ('bounds', 'yv')])

Plotting

Let’s examine these coordinate variables by plotting them.

  1. In [11]: fig, (ax1, ax2) = plt.subplots(ncols=2, figsize=(9,3))
  2.  
  3. In [12]: ds.xc.plot(ax=ax1);
  4.  
  5. In [13]: ds.yc.plot(ax=ax2);

../_images/xarray_multidimensional_coords_8_2.pngNote that the variables xc (longitude) and yc (latitude) aretwo-dimensional scalar fields.

If we try to plot the data variable Tair, by default we get thelogical coordinates.

  1. In [14]: ds.Tair[0].plot();

../_images/xarray_multidimensional_coords_10_1.pngIn order to visualize the data on a conventional latitude-longitudegrid, we can take advantage of xarray’s ability to applycartopy map projections.

  1. In [15]: plt.figure(figsize=(7,2));
  2.  
  3. In [16]: ax = plt.axes(projection=ccrs.PlateCarree());
  4.  
  5. In [17]: ds.Tair[0].plot.pcolormesh(ax=ax, transform=ccrs.PlateCarree(),
  6. ....: x='xc', y='yc', add_colorbar=False);
  7. ....:
  8.  
  9. In [18]: ax.coastlines();

../_images/xarray_multidimensional_coords_12_0.png

Multidimensional Groupby

The above example allowed us to visualize the data on a regularlatitude-longitude grid. But what if we want to do a calculation thatinvolves grouping over one of these physical coordinates (rather thanthe logical coordinates), for example, calculating the mean temperatureat each latitude. This can be achieved using xarray’s groupbyfunction, which accepts multidimensional variables. By default,groupby will use every unique value in the variable, which isprobably not what we want. Instead, we can use the groupby_binsfunction to specify the output coordinates of the group.

  1. # define two-degree wide latitude bins
  2. In [19]: lat_bins = np.arange(0, 91, 2)
  3.  
  4. # define a label for each bin corresponding to the central latitude
  5. In [20]: lat_center = np.arange(1, 90, 2)
  6.  
  7. # group according to those bins and take the mean
  8. In [21]: Tair_lat_mean = (ds.Tair.groupby_bins('xc', lat_bins, labels=lat_center)
  9. ....: .mean(xr.ALL_DIMS))
  10. ....:
  11.  
  12. # plot the result
  13. In [22]: Tair_lat_mean.plot();

../_images/xarray_multidimensional_coords_14_1.pngNote that the resulting coordinate for the groupby_bins operationgot the _bins suffix appended: xc_bins. This help us distinguishit from the original multidimensional variable xc.