{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Estimating Equilibrium Climate Sensitivity (ECS) in CMIP6 models\n", "\n", "**Authors:** [Henry Drake](https://eapsweb.mit.edu/people/hdrake) and [Ryan Abernathey](https://ocean-transport.github.io/)\n", "\n", "**Adapted from** https://github.com/hdrake/cmip6-temperature-demo/blob/master/notebooks/01_calculate_ECS_Gregory_method.ipynb \n", "\n", "---\n", "\n", "*Definition:* Equilibrium Climate Sensitivity is defined as change in global-mean near-surface air temperature (GMST) change due to an instantaneous doubling of CO$_{2}$ concentrations and once the coupled ocean-atmosphere-sea ice system has acheived a statistical equilibrium (i.e. at the top-of-atmosphere, incoming solar shortwave radiation is balanced by reflected solar shortwave and outgoing thermal longwave radiation).\n", "\n", "This notebook uses the [\"Gregory method\"](https://agupubs.onlinelibrary.wiley.com/doi/epdf/10.1029/2003GL018747) to approximate the ECS of CMIP6 models based on the first 150 years after an abrupt quadrupling of CO$_{2}$ concentrations. The \"Gregory Method\" extrapolates the quasi-linear relationship between GMST and radiative imbalance at the top-of-atmosphere to estimate how much warming would occur if the system were in radiative balance at the top-of-atmosphere, which is by definition the equilibrium response. In particular, we extrapolate the linear relationship that occurs between 100 and 150 years after the abrupt quadrupling. Since the radiative forcing due to CO$_{2}$ is a logarithmic function of the CO$_{2}$ concentration, the GMST change from a first doubling is roughly the same as for a second doubling (to first order, we can assume feedbacks as constant), which means that the GMST change due to a quadrupling of CO$_{2}$ is roughly $\\Delta T_{4 \\times \\text{CO}_{2}} = 2 \\times \\text{ECS}$. See also [Mauritsen et al. 2019](https://agupubs.onlinelibrary.wiley.com/doi/epdf/10.1029/2018MS001400) for a detailed application of the Gregory Method (with modifications) for the case of one specific CMIP6 model, the MPI-M Earth System Model.\n", "\n", "For another take on applying the Gregory method to estimate ECS, see [Angeline Pendergrass' code](https://github.com/apendergrass/cmip6-ecs)." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Python packages" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "from matplotlib import pyplot as plt\n", "import numpy as np\n", "import pandas as pd\n", "import xarray as xr\n", "import xesmf as xe\n", "import cartopy\n", "import dask\n", "from tqdm.autonotebook import tqdm # Fancy progress bars for our loops!\n", "import intake\n", "import fsspec\n", "\n", "%matplotlib inline\n", "plt.rcParams['figure.figsize'] = 12, 6\n", "%config InlineBackend.figure_format = 'retina' " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Compute Cluster\n", "\n", "Here we use a dask cluster to parallelize our analysis. The cluster scales up and down adaptively." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "from dask_gateway import Gateway\n", "from dask.distributed import Client\n", "\n", "gateway = Gateway()\n", "cluster = gateway.new_cluster()\n", "cluster.adapt(minimum=1, maximum=20)\n", "client = Client(cluster)\n", "cluster" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Data catalogs\n", "\n", "This notebook uses [`intake-esm`](https://intake-esm.readthedocs.io/en/latest/) to ingest and organize climate model output from the fresh-off-the-supercomputers Phase 6 of the Coupled Model Intercomparison Project (CMIP6). \n", "\n", "The file `https://storage.googleapis.com/cmip6/cmip6-zarr-consolidated-stores.csv` in google cloud storage contains thousands of lines of metadata, each describing an individual climate model experiment's simulated data.\n", "\n", "For example, the first line in the csv file contains the precipitation rate (`variable_id = 'pr'`), as a function of latitude, longitude, and time, in an individual climate model experiment with the BCC-ESM1 model (`source_id = 'BCC-ESM1'`) developed by the Beijing Climate Center (`institution_id = 'BCC'`). The model is *forced* by the forcing experiment SSP370 (`experiment_id = 'ssp370'`), which stands for the Shared Socio-Economic Pathway 3 that results in a change in radiative forcing of $\\Delta F = 7.0$ W/m$^{2}$ from pre-industrial to 2100. This simulation was run as part of the `AerChemMIP` activity, which is a spin-off of the CMIP activity that focuses specifically on how aerosol chemistry affects climate." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "df = pd.read_csv('https://storage.googleapis.com/cmip6/cmip6-zarr-consolidated-stores.csv')\n", "df.head()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The file `pangeo-cmip6.json` describes the structure of the CMIP6 metadata and is formatted so as to be read in by the `intake.open_esm_datastore` method, which categorizes all of the data pointers into a tiered collection. For example, this collection contains the simulated data from 28691 individual experiments, representing 48 different models from 23 different scientific institutions. There are 190 different climate variables (e.g. sea surface temperature, sea ice concentration, atmospheric winds, dissolved organic carbon in the ocean, etc.) available for 29 different forcing experiments.\n", "\n", "### Use Intake ESM\n", "\n", "[Intake-esm](https://intake-esm.readthedocs.io/) is a new package designed to make working with these data archives a bit simpler." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "col = intake.open_esm_datastore(\"https://storage.googleapis.com/cmip6/pangeo-cmip6.json\")\n", "col" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Here, we show the various forcing experiments that climate modellers ran in these simulations." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "df['experiment_id'].unique()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Analysis of Climate Model Output Data" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Loading data\n", "\n", "[intake-esm](https://intake-esm.readthedocs.io/) enables loading data directly into an [xarray.DataArray](http://xarray.pydata.org/en/stable/api.html#dataset), a metadata-aware extension of numpy arrays. `xarray` objects leverage [dask](https://dask.org/) to only read data into memory as needed for any specific operation (i.e. lazy evaluation). Think of `xarray` Datasets as ways of conveniently organizing large arrays of floating point numbers (e.g. climate model data) on an n-dimensional discrete grid, with important metadata such as units, variable, names, etc.\n", "\n", "Note that data on the cloud are in [zarr](https://zarr.readthedocs.io/en/stable/) format, an extension of the metadata-aware format [netcdf](https://www.unidata.ucar.edu/software/netcdf/) commonly used in geosciences.\n", "\n", "`intake-esm` has rules for aggegating datasets; these rules are defined in the collection-specification file.\n", "\n", "#### Choice of simulated forcing experiments\n", "\n", "Here, we choose the `piControl` experiment (in which CO2 concentrations are held fixed at a pre-industrial level of ~300 ppm) and `abrupt-4xCO2` experiment (in which CO2 concentrations are instantaneously quadrupled - or doubled twice - from a pre-industrial controrl state). Since the radiative forcing of CO2 is roughly a logarithmic function of CO2 concentrations, the ECS is roughly independent of the initial CO2 concentration. Thus, if one doubling of CO2 results in $ECS$ of warming, then two doublings (or, a quadrupling) results in $2 \\times ECS$ of warming.\n", "\n", "Ideally, we would choose the `abrupt-2xCO2` forcing experiment, but more 4xCO2 data are currently avaiable in Google Cloud Storage, making for a better example.\n", "\n", "### Prepare Data" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "query = dict(\n", " experiment_id=['abrupt-4xCO2','piControl'], # pick the `abrupt-4xCO2` and `piControl` forcing experiments\n", " table_id='Amon', # choose to look at atmospheric variables (A) saved at monthly resolution (mon)\n", " variable_id=['tas', 'rsut','rsdt','rlut'], # choose to look at near-surface air temperature (tas) as our variable\n", " member_id = 'r1i1p1f1', # arbitrarily pick one realization for each model (i.e. just one set of initial conditions)\n", ")\n", "\n", "col_subset = col.search(require_all_on=[\"source_id\"], **query)\n", "col_subset.df.groupby(\"source_id\")[\n", " [\"experiment_id\", \"variable_id\", \"table_id\"]\n", "].nunique()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The following functions help us load and homogenize the data.\n", "We use some [dask.delayed](https://docs.dask.org/en/latest/delayed.html) programming to open the datasets in parallel." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "def drop_all_bounds(ds):\n", " \"\"\"Drop coordinates like 'time_bounds' from datasets,\n", " which can lead to issues when merging.\"\"\"\n", " drop_vars = [vname for vname in ds.coords\n", " if (('_bounds') in vname ) or ('_bnds') in vname]\n", " return ds.drop(drop_vars)\n", "\n", "def open_dsets(df):\n", " \"\"\"Open datasets from cloud storage and return xarray dataset.\"\"\"\n", " dsets = [xr.open_zarr(fsspec.get_mapper(ds_url), consolidated=True)\n", " .pipe(drop_all_bounds)\n", " for ds_url in df.zstore]\n", " try:\n", " ds = xr.merge(dsets, join='exact')\n", " return ds\n", " except ValueError:\n", " return None\n", "\n", "def open_delayed(df):\n", " \"\"\"A dask.delayed wrapper around `open_dsets`.\n", " Allows us to open many datasets in parallel.\"\"\"\n", " return dask.delayed(open_dsets)(df)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Create a nested dictionary of models and experiments. It will be structured like this:\n", "\n", " {'CESM2':\n", " {\n", " 'piControl': ,\n", " 'abrupt-4xCO2': \n", " },\n", " ...\n", " }" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "from collections import defaultdict\n", "\n", "dsets = defaultdict(dict) \n", "for group, df in col_subset.df.groupby(by=['source_id', 'experiment_id']):\n", " dsets[group[0]][group[1]] = open_delayed(df)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Open one of the datasets directly, just to show what it looks like:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "%time open_dsets(df)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now use dask to do this in parallel on all of the datasets:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "dsets_ = dask.compute(dict(dsets))[0]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Reduce Data via Global Mean\n", "\n", "We don't want to load all of the raw model data into memory right away. Instead, we want to reduce the data by taking the global mean. We need to remember to weight this global mean by a factor proportional to `cos(lat)`." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "def get_lat_name(ds):\n", " \"\"\"Figure out what is the latitude coordinate for each dataset.\"\"\"\n", " for lat_name in ['lat', 'latitude']:\n", " if lat_name in ds.coords:\n", " return lat_name\n", " raise RuntimeError(\"Couldn't find a latitude coordinate\")\n", "\n", "def global_mean(ds):\n", " \"\"\"Return global mean of a whole dataset.\"\"\"\n", " lat = ds[get_lat_name(ds)]\n", " weight = np.cos(np.deg2rad(lat))\n", " weight /= weight.mean()\n", " other_dims = set(ds.dims) - {'time'}\n", " return (ds * weight).mean(other_dims)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We now apply this function, plus resampling to annual mean data, to all of the datasets.\n", "We also concatenate the experiments together into a single Dataset for each model.\n", "This is the most complex cell in the notebook. A lot is happening here." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "expts = ['piControl', 'abrupt-4xCO2']\n", "expt_da = xr.DataArray(expts, dims='experiment_id',\n", " coords={'experiment_id': expts})\n", "\n", "dsets_aligned = {}\n", "\n", "for k, v in tqdm(dsets_.items()):\n", " expt_dsets = v.values()\n", " if any([d is None for d in expt_dsets]):\n", " print(f\"Missing experiment for {k}\")\n", " continue\n", " \n", " for ds in expt_dsets:\n", " ds.coords['year'] = ds.time.dt.year - ds.time.dt.year[0]\n", " \n", " # workaround for\n", " # https://github.com/pydata/xarray/issues/2237#issuecomment-620961663\n", " dsets_ann_mean = [v[expt].pipe(global_mean)\n", " .swap_dims({'time': 'year'})\n", " .drop('time')\n", " .coarsen(year=12).mean()\n", " for expt in expts]\n", " \n", " # align everything with the 4xCO2 experiment\n", " dsets_aligned[k] = xr.concat(dsets_ann_mean, join='right',\n", " dim=expt_da)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Do the Computation\n", "\n", "Up to this point, no computations have actually happened. Everything has been \"lazy\". Now we trigger the computation to actual occur and load the global / annual mean data into memory." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "dsets_aligned_ = dask.compute(dsets_aligned)[0]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now we concatenate across models to produce one big dataset with all the required variables." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "source_ids = list(dsets_aligned_.keys())\n", "source_da = xr.DataArray(source_ids, dims='source_id',\n", " coords={'source_id': source_ids})\n", "\n", "big_ds = xr.concat([ds.reset_coords(drop=True)\n", " for ds in dsets_aligned_.values()],\n", " dim=source_da)\n", "big_ds" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Calculated Derived Variables\n", "\n", "We need to calculate the net radiative imbalance, plus the anomaly of the abrupt 4xCO2 run compared to the piControl run." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "big_ds['imbalance'] = big_ds['rsdt'] - big_ds['rsut'] - big_ds['rlut']\n", "\n", "ds_mean = big_ds[['tas', 'imbalance']].sel(experiment_id='piControl').mean(dim='year')\n", "ds_anom = big_ds[['tas', 'imbalance']] - ds_mean\n", "\n", "# add some metadata\n", "ds_anom.tas.attrs['long_name'] = 'Global Mean Surface Temp Anom'\n", "ds_anom.tas.attrs['units'] = 'K'\n", "ds_anom.imbalance.attrs['long_name'] = 'Global Mean Radiative Imbalance'\n", "ds_anom.imbalance.attrs['units'] = 'W m$^{-2}$'\n", "\n", "ds_anom" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Plot Timeseries\n", "\n", "Here we plot the global mean surface temperature for each model:" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "ds_anom.tas.plot.line(col='source_id', x='year', col_wrap=5)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We can see that the models cover different time intervals. Let's limit the rest of our analysis to the first 150 years." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "# limit to the gregory 150-year period\n", "first_150_years = slice(0, 149)\n", "ds_anom.tas.sel(year=first_150_years).plot.line(col='source_id', x='year', col_wrap=5)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Same thing for radiative imbalance." ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "ds_anom.imbalance.sel(year=first_150_years).plot.line(col='source_id', x='year', col_wrap=5)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Calculate ECS" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "ds_abrupt = ds_anom.sel(year=first_150_years, experiment_id='abrupt-4xCO2').reset_coords(drop=True)\n", "\n", "def calc_ecs(ds):\n", " # Some sources don't have all 150 years, drop those missing values.\n", " a, b = np.polyfit(ds.tas.dropna(\"year\"),\n", " ds.imbalance.dropna(\"year\"), 1)\n", " ecs = -0.5 * (b/a)\n", " return xr.DataArray(ecs)\n", "\n", "ds_abrupt['ecs'] = ds_abrupt.groupby('source_id').apply(calc_ecs)\n", "ds_abrupt" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Reproduce a plot similar to [Mark Zelinka's](https://twitter.com/mzelinka/status/1255534531144085513):" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "fg = ds_abrupt.plot.scatter(x='tas', y='imbalance', col='source_id', col_wrap=5)\n", "\n", "def calc_and_plot_ecs(x, y, **kwargs):\n", " x = x[~np.isnan(x)]\n", " y = y[~np.isnan(y)]\n", " a, b = np.polyfit(x, y, 1)\n", " ecs = -0.5 * b/a\n", " plt.autoscale(False)\n", " plt.plot([0, 10], np.polyval([a, b], [0, 10]), 'k')\n", " plt.text(4, 6, f'ECS = {ecs:3.2f}', fontdict={'weight': 'bold', 'size': 12})\n", " plt.grid()\n", " \n", "fg.map(calc_and_plot_ecs, 'tas', 'imbalance')" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "ds_abrupt.ecs.plot.hist();" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "ds_abrupt.ecs.to_dataframe().sort_values('ecs').plot(kind='bar')" ] } ], "metadata": { "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.7.7" } }, "nbformat": 4, "nbformat_minor": 4 }