{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Visualization Gallery\n", "\n", "This notebook shows common visualization issues encountered in xarray." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "import cartopy.crs as ccrs\n", "import matplotlib.pyplot as plt\n", "import xarray as xr\n", "\n", "%matplotlib inline" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Load example dataset:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "ds = xr.tutorial.load_dataset(\"air_temperature\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Multiple plots and map projections\n", "\n", "Control the map projection parameters on multiple axes\n", "\n", "This example illustrates how to plot multiple maps and control their extent\n", "and aspect ratio.\n", "\n", "For more details see [this discussion](https://github.com/pydata/xarray/issues/1397#issuecomment-299190567) on github." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "air = ds.air.isel(time=[0, 724]) - 273.15\n", "\n", "# This is the map projection we want to plot *onto*\n", "map_proj = ccrs.LambertConformal(central_longitude=-95, central_latitude=45)\n", "\n", "p = air.plot(\n", " transform=ccrs.PlateCarree(), # the data's projection\n", " col=\"time\",\n", " col_wrap=1, # multiplot settings\n", " aspect=ds.dims[\"lon\"] / ds.dims[\"lat\"], # for a sensible figsize\n", " subplot_kws={\"projection\": map_proj},\n", ") # the plot's projection\n", "\n", "# We have to set the map's options on all axes\n", "for ax in p.axes.flat:\n", " ax.coastlines()\n", " ax.set_extent([-160, -30, 5, 75])" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Centered colormaps\n", "\n", "Xarray's automatic colormaps choice" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "air = ds.air.isel(time=0)\n", "\n", "f, ((ax1, ax2), (ax3, ax4)) = plt.subplots(2, 2, figsize=(8, 6))\n", "\n", "# The first plot (in kelvins) chooses \"viridis\" and uses the data's min/max\n", "air.plot(ax=ax1, cbar_kwargs={\"label\": \"K\"})\n", "ax1.set_title(\"Kelvins: default\")\n", "ax2.set_xlabel(\"\")\n", "\n", "# The second plot (in celsius) now chooses \"BuRd\" and centers min/max around 0\n", "airc = air - 273.15\n", "airc.plot(ax=ax2, cbar_kwargs={\"label\": \"°C\"})\n", "ax2.set_title(\"Celsius: default\")\n", "ax2.set_xlabel(\"\")\n", "ax2.set_ylabel(\"\")\n", "\n", "# The center doesn't have to be 0\n", "air.plot(ax=ax3, center=273.15, cbar_kwargs={\"label\": \"K\"})\n", "ax3.set_title(\"Kelvins: center=273.15\")\n", "\n", "# Or it can be ignored\n", "airc.plot(ax=ax4, center=False, cbar_kwargs={\"label\": \"°C\"})\n", "ax4.set_title(\"Celsius: center=False\")\n", "ax4.set_ylabel(\"\")\n", "\n", "# Make it nice\n", "plt.tight_layout()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Control the plot's colorbar\n", "\n", "Use ``cbar_kwargs`` keyword to specify the number of ticks.\n", "The ``spacing`` kwarg can be used to draw proportional ticks." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "air2d = ds.air.isel(time=500)\n", "\n", "# Prepare the figure\n", "f, (ax1, ax2, ax3) = plt.subplots(1, 3, figsize=(14, 4))\n", "\n", "# Irregular levels to illustrate the use of a proportional colorbar\n", "levels = [245, 250, 255, 260, 265, 270, 275, 280, 285, 290, 310, 340]\n", "\n", "# Plot data\n", "air2d.plot(ax=ax1, levels=levels)\n", "air2d.plot(ax=ax2, levels=levels, cbar_kwargs={\"ticks\": levels})\n", "air2d.plot(\n", " ax=ax3, levels=levels, cbar_kwargs={\"ticks\": levels, \"spacing\": \"proportional\"}\n", ")\n", "\n", "# Show plots\n", "plt.tight_layout()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Multiple lines from a 2d DataArray\n", "\n", "Use ``xarray.plot.line`` on a 2d DataArray to plot selections as\n", "multiple lines.\n", "\n", "See ``plotting.multiplelines`` for more details." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "air = ds.air - 273.15 # to celsius\n", "\n", "# Prepare the figure\n", "f, (ax1, ax2) = plt.subplots(1, 2, figsize=(8, 4), sharey=True)\n", "\n", "# Selected latitude indices\n", "isel_lats = [10, 15, 20]\n", "\n", "# Temperature vs longitude plot - illustrates the \"hue\" kwarg\n", "air.isel(time=0, lat=isel_lats).plot.line(ax=ax1, hue=\"lat\")\n", "ax1.set_ylabel(\"°C\")\n", "\n", "# Temperature vs time plot - illustrates the \"x\" and \"add_legend\" kwargs\n", "air.isel(lon=30, lat=isel_lats).plot.line(ax=ax2, x=\"time\", add_legend=False)\n", "ax2.set_ylabel(\"\")\n", "\n", "# Show\n", "plt.tight_layout()" ] }, { "cell_type": "markdown", "metadata": { "jp-MarkdownHeadingCollapsed": true, "tags": [] }, "source": [ "## `imshow()` and rasterio map projections\n", "\n", "\n", "Using rasterio's projection information for more accurate plots.\n", "\n", "This example extends `recipes.rasterio` and plots the image in the\n", "original map projection instead of relying on pcolormesh and a map\n", "transformation." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "da = xr.tutorial.open_rasterio(\"RGB.byte\")\n", "\n", "# The data is in UTM projection. We have to set it manually until\n", "# https://github.com/SciTools/cartopy/issues/813 is implemented\n", "crs = ccrs.UTM(\"18\")\n", "\n", "# Plot on a map\n", "ax = plt.subplot(projection=crs)\n", "da.plot.imshow(ax=ax, rgb=\"band\", transform=crs)\n", "ax.coastlines(\"10m\", color=\"r\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Parsing rasterio geocoordinates\n", "\n", "Converting a projection's cartesian coordinates into 2D longitudes and\n", "latitudes.\n", "\n", "These new coordinates might be handy for plotting and indexing, but it should\n", "be kept in mind that a grid which is regular in projection coordinates will\n", "likely be irregular in lon/lat. It is often recommended to work in the data's\n", "original map projection (see `recipes.rasterio_rgb`)." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "from pyproj import Transformer\n", "import numpy as np\n", "\n", "da = xr.tutorial.open_rasterio(\"RGB.byte\")\n", "\n", "x, y = np.meshgrid(da[\"x\"], da[\"y\"])\n", "transformer = Transformer.from_crs(da.crs, \"EPSG:4326\", always_xy=True)\n", "lon, lat = transformer.transform(x, y)\n", "da.coords[\"lon\"] = ((\"y\", \"x\"), lon)\n", "da.coords[\"lat\"] = ((\"y\", \"x\"), lat)\n", "\n", "# Compute a greyscale out of the rgb image\n", "greyscale = da.mean(dim=\"band\")\n", "\n", "# Plot on a map\n", "ax = plt.subplot(projection=ccrs.PlateCarree())\n", "greyscale.plot(\n", " ax=ax,\n", " x=\"lon\",\n", " y=\"lat\",\n", " transform=ccrs.PlateCarree(),\n", " cmap=\"Greys_r\",\n", " shading=\"auto\",\n", " add_colorbar=False,\n", ")\n", "ax.coastlines(\"10m\", color=\"r\")" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3 (ipykernel)", "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.9.7" } }, "nbformat": 4, "nbformat_minor": 4 }