{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Working with Multidimensional Coordinates\n", "\n", "Author: [Ryan Abernathey](https://github.com/rabernat)\n", "\n", "Many datasets have _physical coordinates_ which differ from their _logical coordinates_. Xarray provides several ways to plot and analyze such datasets." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "ExecuteTime": { "end_time": "2018-11-28T20:49:56.068395Z", "start_time": "2018-11-28T20:49:56.035349Z" } }, "outputs": [], "source": [ "%matplotlib inline\n", "import numpy as np\n", "import pandas as pd\n", "import xarray as xr\n", "import cartopy.crs as ccrs\n", "from matplotlib import pyplot as plt" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "As an example, consider this dataset from the [xarray-data](https://github.com/pydata/xarray-data) repository." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "ExecuteTime": { "end_time": "2018-11-28T20:50:13.629720Z", "start_time": "2018-11-28T20:50:13.484542Z" } }, "outputs": [], "source": [ "ds = xr.tutorial.open_dataset(\"rasm\").load()\n", "ds" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In this example, the _logical coordinates_ are `x` and `y`, while the _physical coordinates_ are `xc` and `yc`, which represent the longitudes and latitudes of the data." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "ExecuteTime": { "end_time": "2018-11-28T20:50:15.836061Z", "start_time": "2018-11-28T20:50:15.768376Z" } }, "outputs": [], "source": [ "print(ds.xc.attrs)\n", "print(ds.yc.attrs)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Plotting ##\n", "\n", "Let's examine these coordinate variables by plotting them." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "ExecuteTime": { "end_time": "2018-11-28T20:50:17.928556Z", "start_time": "2018-11-28T20:50:17.031211Z" } }, "outputs": [], "source": [ "fig, (ax1, ax2) = plt.subplots(ncols=2, figsize=(14, 4))\n", "ds.xc.plot(ax=ax1)\n", "ds.yc.plot(ax=ax2)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note that the variables `xc` (longitude) and `yc` (latitude) are two-dimensional scalar fields.\n", "\n", "If we try to plot the data variable `Tair`, by default we get the logical coordinates." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "ExecuteTime": { "end_time": "2018-11-28T20:50:20.567749Z", "start_time": "2018-11-28T20:50:19.999393Z" } }, "outputs": [], "source": [ "ds.Tair[0].plot()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "In order to visualize the data on a conventional latitude-longitude grid, we can take advantage of xarray's ability to apply [cartopy](http://scitools.org.uk/cartopy/index.html) map projections." ] }, { "cell_type": "code", "execution_count": null, "metadata": { "ExecuteTime": { "end_time": "2018-11-28T20:50:31.131708Z", "start_time": "2018-11-28T20:50:30.444697Z" } }, "outputs": [], "source": [ "plt.figure(figsize=(14, 6))\n", "ax = plt.axes(projection=ccrs.PlateCarree())\n", "ax.set_global()\n", "ds.Tair[0].plot.pcolormesh(\n", " ax=ax, transform=ccrs.PlateCarree(), x=\"xc\", y=\"yc\", add_colorbar=False\n", ")\n", "ax.coastlines()\n", "ax.set_ylim([0, 90]);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Multidimensional Groupby ##\n", "\n", "The above example allowed us to visualize the data on a regular latitude-longitude grid. But what if we want to do a calculation that involves grouping over one of these physical coordinates (rather than the logical coordinates), for example, calculating the mean temperature at each latitude. This can be achieved using xarray's `groupby` function, which accepts multidimensional variables. By default, `groupby` will use every unique value in the variable, which is probably not what we want. Instead, we can use the `groupby_bins` function to specify the output coordinates of the group. " ] }, { "cell_type": "code", "execution_count": null, "metadata": { "ExecuteTime": { "end_time": "2018-11-28T20:50:43.670463Z", "start_time": "2018-11-28T20:50:43.245501Z" } }, "outputs": [], "source": [ "# define two-degree wide latitude bins\n", "lat_bins = np.arange(0, 91, 2)\n", "# define a label for each bin corresponding to the central latitude\n", "lat_center = np.arange(1, 90, 2)\n", "# group according to those bins and take the mean\n", "Tair_lat_mean = ds.Tair.groupby_bins(\"yc\", lat_bins, labels=lat_center).mean(\n", " dim=xr.ALL_DIMS\n", ")\n", "# plot the result\n", "Tair_lat_mean.plot()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The resulting coordinate for the `groupby_bins` operation got the `_bins` suffix appended: `yc_bins`. This help us distinguish it from the original multidimensional variable `yc`.\n", "\n", "**Note**: This group-by-latitude approach does not take into account the finite-size geometry of grid cells. It simply bins each value according to the coordinates at the cell center. Xarray has no understanding of grid cells and their geometry. More precise geographic regridding for xarray data is available via the [xesmf](https://xesmf.readthedocs.io) package." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [] } ], "metadata": { "anaconda-cloud": {}, "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.6.8" }, "toc": { "base_numbering": 1, "nav_menu": {}, "number_sections": true, "sideBar": true, "skip_h1_title": false, "title_cell": "Table of Contents", "title_sidebar": "Contents", "toc_cell": true, "toc_position": {}, "toc_section_display": true, "toc_window_display": true } }, "nbformat": 4, "nbformat_minor": 2 }