{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "## Calculate ENSO Skill as a Function of Initial Month vs. Lead Time" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### In this example, we demonstrate: \n", "1. How to remotely access data from the North American Multi-model Ensemble (NMME) hindcast database and set it up to be used in `climpred`. \n", "2. How to calculate the Anomaly Correlation Coefficient (ACC) using monthly data\n", "3. How to calculate and plot historical forecast skill of the Nino3.4 index as function of initialization month and lead time." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### The North American Multi-model Ensemble (NMME)\n", "\n", "Further information on NMME is available from [Kirtman et al. 2014](https://journals.ametsoc.org/doi/full/10.1175/BAMS-D-12-00050.1) and the [NMME project website](https://www.cpc.ncep.noaa.gov/products/NMME/)\n", "\n", "The NMME public database is hosted on the International Research Institute for Climate and Society (IRI) data server http://iridl.ldeo.columbia.edu/SOURCES/.Models/.NMME/\n", "\n", "Since the NMME data server is accessed via this notebook, the time for the notebook to run may take a few minutes and vary depending on the speed that data is downloaded." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### Definitions\n", "\n", "**Anomalies**\n", ": Departure from normal, where normal is defined as the climatological value based on the average value for each month over all years.\n", "\n", "**Nino3.4**\n", ": An index used to represent the evolution of the El Nino-Southern Oscillation (ENSO). Calculated as the average sea surface temperature (SST) anomalies in the region 5S-5N; 190-240" ] }, { "cell_type": "code", "execution_count": 1, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:26:56.235431Z", "start_time": "2020-01-06T20:26:54.064651Z" } }, "outputs": [], "source": [ "import warnings\n", "\n", "import matplotlib.pyplot as plt\n", "import xarray as xr\n", "import pandas as pd\n", "import numpy as np\n", "from tqdm.auto import tqdm\n", "\n", "from climpred import HindcastEnsemble\n", "import climpred" ] }, { "cell_type": "code", "execution_count": 2, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:26:56.238990Z", "start_time": "2020-01-06T20:26:56.236743Z" } }, "outputs": [], "source": [ "warnings.filterwarnings(\"ignore\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Function to set 360 calendar to 360_day calendar and decond cf times" ] }, { "cell_type": "code", "execution_count": 3, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:26:56.595041Z", "start_time": "2020-01-06T20:26:56.592324Z" } }, "outputs": [], "source": [ "def decode_cf(ds, time_var):\n", " if ds[time_var].attrs['calendar'] == '360':\n", " ds[time_var].attrs['calendar'] = '360_day'\n", " ds = xr.decode_cf(ds, decode_times=True)\n", " return ds" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Load the monthly sea surface temperature (SST) hindcast data for the NCEP-CFSv2 model from the NMME data server. This is a large dataset, so we allow dask to chunk the data as it chooses. " ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:26:57.233967Z", "start_time": "2020-01-06T20:26:56.904799Z" } }, "outputs": [ { "data": { "text/html": [ "
<xarray.Dataset>\n",
       "Dimensions:  (L: 10, M: 24, S: 348, X: 360, Y: 181)\n",
       "Coordinates:\n",
       "  * S        (S) object 1982-01-01 00:00:00 ... 2010-12-01 00:00:00\n",
       "  * M        (M) float32 1.0 2.0 3.0 4.0 5.0 6.0 ... 20.0 21.0 22.0 23.0 24.0\n",
       "  * X        (X) float32 0.0 1.0 2.0 3.0 4.0 ... 355.0 356.0 357.0 358.0 359.0\n",
       "  * L        (L) float32 0.5 1.5 2.5 3.5 4.5 5.5 6.5 7.5 8.5 9.5\n",
       "  * Y        (Y) float32 -90.0 -89.0 -88.0 -87.0 -86.0 ... 87.0 88.0 89.0 90.0\n",
       "Data variables:\n",
       "    sst      (S, L, M, Y, X) float32 dask.array<chunksize=(6, 5, 8, 181, 360), meta=np.ndarray>\n",
       "Attributes:\n",
       "    Conventions:  IRIDL
" ], "text/plain": [ "\n", "Dimensions: (L: 10, M: 24, S: 348, X: 360, Y: 181)\n", "Coordinates:\n", " * S (S) object 1982-01-01 00:00:00 ... 2010-12-01 00:00:00\n", " * M (M) float32 1.0 2.0 3.0 4.0 5.0 6.0 ... 20.0 21.0 22.0 23.0 24.0\n", " * X (X) float32 0.0 1.0 2.0 3.0 4.0 ... 355.0 356.0 357.0 358.0 359.0\n", " * L (L) float32 0.5 1.5 2.5 3.5 4.5 5.5 6.5 7.5 8.5 9.5\n", " * Y (Y) float32 -90.0 -89.0 -88.0 -87.0 -86.0 ... 87.0 88.0 89.0 90.0\n", "Data variables:\n", " sst (S, L, M, Y, X) float32 dask.array\n", "Attributes:\n", " Conventions: IRIDL" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "url = 'http://iridl.ldeo.columbia.edu/SOURCES/.Models/.NMME/NCEP-CFSv2/.HINDCAST/.MONTHLY/.sst/dods'\n", "fcstds = decode_cf(xr.open_dataset(url, decode_times=False, \n", " chunks={'S': 'auto', 'L': 'auto', 'M':'auto'}),'S')\n", "\n", "fcstds" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The NMME data dimensions correspond to the following `climpred` dimension definitions: `X=lon`,`L=lead`,`Y=lat`,`M=member`, `S=init`. We will rename the dimensions to their `climpred` names." ] }, { "cell_type": "code", "execution_count": 5, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:26:58.652387Z", "start_time": "2020-01-06T20:26:58.649738Z" } }, "outputs": [], "source": [ "fcstds=fcstds.rename({'S': 'init','L': 'lead','M': 'member', 'X': 'lon', 'Y': 'lat'})" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Let's make sure that the `lead` dimension is set properly for `climpred`. NMME data stores `leads` as 0.5, 1.5, 2.5, etc, which correspond to 0, 1, 2, ... months since initialization. We will change the `lead` to be integers starting with zero. `climpred` also requires that `lead` dimension has an attribute called `units` indicating what time units the `lead` is assocated with. Options are: `years,seasons,months,weeks,pentads,days`. For the monthly NMME data, the `lead` `units` are `months`. " ] }, { "cell_type": "code", "execution_count": 6, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:26:58.932066Z", "start_time": "2020-01-06T20:26:58.928398Z" } }, "outputs": [], "source": [ "fcstds['lead']=(fcstds['lead']-0.5).astype('int')\n", "fcstds['lead'].attrs={'units': 'months'}" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now we need to make sure that the `init` dimension is set properly for `climpred`. For monthly data, the `init` dimension must be a `xr.cfdateTimeIndex` or a `pd.datetimeIndex`. We convert the `init` values to `pd.datatimeIndex`." ] }, { "cell_type": "code", "execution_count": 7, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:27:00.174981Z", "start_time": "2020-01-06T20:27:00.160531Z" } }, "outputs": [], "source": [ "fcstds['init']=pd.to_datetime(fcstds.init.values.astype(str))\n", "fcstds['init']=pd.to_datetime(fcstds['init'].dt.strftime('%Y%m01 00:00'))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Next, we want to get the verification SST data from the data server" ] }, { "cell_type": "code", "execution_count": 8, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:27:01.118878Z", "start_time": "2020-01-06T20:27:00.760475Z" } }, "outputs": [ { "data": { "text/html": [ "
<xarray.Dataset>\n",
       "Dimensions:  (T: 405, X: 360, Y: 181)\n",
       "Coordinates:\n",
       "  * Y        (Y) float32 -90.0 -89.0 -88.0 -87.0 -86.0 ... 87.0 88.0 89.0 90.0\n",
       "  * X        (X) float32 0.0 1.0 2.0 3.0 4.0 ... 355.0 356.0 357.0 358.0 359.0\n",
       "  * T        (T) object 1982-01-16 00:00:00 ... 2015-09-16 00:00:00\n",
       "Data variables:\n",
       "    sst      (T, Y, X) float32 ...\n",
       "Attributes:\n",
       "    Conventions:  IRIDL
" ], "text/plain": [ "\n", "Dimensions: (T: 405, X: 360, Y: 181)\n", "Coordinates:\n", " * Y (Y) float32 -90.0 -89.0 -88.0 -87.0 -86.0 ... 87.0 88.0 89.0 90.0\n", " * X (X) float32 0.0 1.0 2.0 3.0 4.0 ... 355.0 356.0 357.0 358.0 359.0\n", " * T (T) object 1982-01-16 00:00:00 ... 2015-09-16 00:00:00\n", "Data variables:\n", " sst (T, Y, X) float32 ...\n", "Attributes:\n", " Conventions: IRIDL" ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "obsurl='http://iridl.ldeo.columbia.edu/SOURCES/.Models/.NMME/.OIv2_SST/.sst/dods'\n", "verifds = decode_cf(xr.open_dataset(obsurl, decode_times=False),'T')\n", "verifds" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Rename the dimensions to correspond to climpred dimensions" ] }, { "cell_type": "code", "execution_count": 9, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:27:01.531745Z", "start_time": "2020-01-06T20:27:01.528852Z" } }, "outputs": [], "source": [ "verifds=verifds.rename({'T': 'time','X': 'lon', 'Y': 'lat'})" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Convert the `time` data to be of type `pd.datetimeIndex`" ] }, { "cell_type": "code", "execution_count": 10, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:27:02.082512Z", "start_time": "2020-01-06T20:27:02.072568Z" } }, "outputs": [], "source": [ "verifds['time']=pd.to_datetime(verifds.time.values.astype(str))\n", "verifds['time']=pd.to_datetime(verifds['time'].dt.strftime('%Y%m01 00:00'))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Subset the data to 1982-2010" ] }, { "cell_type": "code", "execution_count": 11, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:27:03.439898Z", "start_time": "2020-01-06T20:27:03.432591Z" } }, "outputs": [], "source": [ "fcstds=fcstds.sel(init=slice('1982-01-01','2010-12-01'))\n", "verifds=verifds.sel(time=slice('1982-01-01','2010-12-01'))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Calculate the Nino3.4 index for forecast and verification." ] }, { "cell_type": "code", "execution_count": 12, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:27:04.626000Z", "start_time": "2020-01-06T20:27:03.912618Z" } }, "outputs": [], "source": [ "fcstnino34=fcstds.sel(lat=slice(-5,5),lon=slice(190,240)).mean(['lat','lon'])\n", "verifnino34=verifds.sel(lat=slice(-5,5),lon=slice(190,240)).mean(['lat','lon'])\n", "\n", "fcstclimo = fcstnino34.groupby('init.month').mean('init')\n", "fcst = (fcstnino34.groupby('init.month') - fcstclimo)\n", "\n", "verifclimo = verifnino34.groupby('time.month').mean('time')\n", "verif = (verifnino34.groupby('time.month') - verifclimo)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Because will will calculate the anomaly correlation coefficient over all `time` for verification and `init` for the hindcasts, we need to rechunk the data so that these dimensions are in same chunk" ] }, { "cell_type": "code", "execution_count": 13, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:27:05.228399Z", "start_time": "2020-01-06T20:27:05.216120Z" } }, "outputs": [], "source": [ "fcst=fcst.chunk({'init':-1})\n", "verif=verif.chunk({'time':-1})" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Use the `climpred HindcastEnsemble` to calculate the anomaly correlation coefficient (ACC) as a function of initial month and `lead` " ] }, { "cell_type": "code", "execution_count": 14, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:31:40.476350Z", "start_time": "2020-01-06T20:27:06.112662Z" } }, "outputs": [ { "data": { "application/vnd.jupyter.widget-view+json": { "model_id": "a275d729d0f74b78b7613c296266f84b", "version_major": 2, "version_minor": 0 }, "text/plain": [ "HBox(children=(FloatProgress(value=0.0, max=12.0), HTML(value='')))" ] }, "metadata": {}, "output_type": "display_data" }, { "name": "stdout", "output_type": "stream", "text": [ "\n" ] } ], "source": [ "skill=np.zeros((fcst['lead'].size, 12))\n", "for im in tqdm(np.arange(0,12)):\n", " hindcast = HindcastEnsemble(fcst.sel(init=fcst['init.month']==im+1))\n", " hindcast = hindcast.add_observations(verif, 'observations')\n", " skillds = hindcast.verify(metric='acc')\n", " skill[:,im]=skillds['sst'].values" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Plot the ACC as function of Initial Month and lead-time" ] }, { "cell_type": "code", "execution_count": 15, "metadata": { "ExecuteTime": { "end_time": "2020-01-06T20:32:08.431825Z", "start_time": "2020-01-06T20:32:08.262267Z" } }, "outputs": [ { "data": { "text/plain": [ "Text(0, 0.5, 'Lead Time (Months)')" ] }, "execution_count": 15, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAW4AAAEWCAYAAABG030jAAAABHNCSVQICAgIfAhkiAAAAAlwSFlzAAALEgAACxIB0t1+/AAAADh0RVh0U29mdHdhcmUAbWF0cGxvdGxpYiB2ZXJzaW9uMy4xLjIsIGh0dHA6Ly9tYXRwbG90bGliLm9yZy8li6FKAAAgAElEQVR4nO3deZxcZZ3v8c+3OwkhhBAgwNUEJGpEGS9rRBR1QEQBEVyYK4iKuCAqi7sRN3RmvM64oqJ55SKbMiACo0FRQBS4DoLsEAhLJkJowxbWkL2T3/zxnMaiUl11Ol2nzznN9/16nVe6zjn1nF9VV3791HOeRRGBmZnVR0/ZAZiZ2dA4cZuZ1YwTt5lZzThxm5nVjBO3mVnNOHGbmdWME7c9J0g6UdKpZcdh1g1O3BUj6V5JD0napGHfByVd0fBYko6XNE/SMkl9kn4h6X9nx8+QtFrS0w3bLdmx7SVFw/57Jc3qENMeki6W9ISkxyT9RdJR2bG9Ja1rutZF2bHJkk6T9KCkpZLulvS5HO/BRpJ+Ium+7Hk3STqgzfnvy17TZ5r290naGyAivh4RH+x07Ryx7SjpekmPZ9vvJe2Y43kzJK2U9LMc5+6dvZ7Ptjg2TtJJku7Jfvf3Zu/x9g3nvEnSVdl794ikKyUdPNTXatXlxF1NY4AT2hw/OTt+PLAF8BLgl8CbG87594iY2LDt3FTG5IiYCBwOfFnS/q0uJOlVwB+AK4EXA1sCHwEaE+nipmu9Jdv/XWAi8DJgM+Bg4L87vHZIr/9+4B+z530JOK8xObXwGPA5SZNylD8ci4FDSe/7FGAucG6O550CXJfzGkeSXs+RLY6dT3of30V6b3YGbgD2BZB0KPAL4CxgGrAN8GXgLS3KsrqKCG8V2oB7gVmk/7iTs30fBK7Ifp4BrAX2aFPGGcC/DHJseyCAMQ37rgM+Pcj5fwJOaXOtvYG+QY7NA946yLHZwLea9v0K+OQg598KvGOQY+/L4rwI+ErD/j5g7+znk4CfNb0HRwKLgCXAFxqetxHwPVKSXpz9vFGL644BPgYs7/A7PQw4rzGGNudOAJZmz1kNzGw49gZgBbDtIM9V9no+U/bn2Fuxm2vc1XQ9cAXw6RbH9iUlyr8M9yJZk8tewD8AN7U4PgF4FamWtyGuAf5V0lGSZjQd+w/gnZKUXWtz4I20qL1K2ob0reL2Dtf7EvAJSVvkjO81wA6k9/TLkl6W7f8CsCewC6lGuwfwxaaYngBWAj8Avj7YBbJvAF8DPpUzpncAT5NqzZcA72049gbgLxFx/yDP3QHYlg3/fVlNOHFX15eB4yRt1bR/S+CBHM//dNYmPbCd2XR8CalWfyowKyIub1HG5qTPSKfrPb/pWv8n238ccDZwLHCHpAUNbdX/n1TrfW32+FDgzxGxuLFgSWOzMs6MiDvbBRERNwOXAh3b0TNfjYgVEXELcAspSQMcAXwtIh6OiEeArwLvabrWZFJTxbG0+KPX4J+Bn7RJts2OBH4eEWtJf9wOz94D6Py73zL7N8/nw2rMibuiImIe8GtSs0mjR4Hn5SjiWxExuWFrbi+dEhGbR8TLIuL78EzPi4EbjLOBx4F1Oa63uOla52WvYUWkm4K7k5LKecAvJG0REUGqXR+elfEuUoJ+hqQe4KekJoNjc7xmSH/wPiLpf+U498GGn5eT2uMBng/c13Dsvmzfs0TEMlKTz1mStm4+LmkXUi35u3kCl7QtsA9/fx9+BYzn7/cuOv3uH83+zfP5sBpz4q62rwAfAqY27LscmCZpZrcvliXZgRuMx0TEcuDPpK/vwy37KVKTwibA9Gz3OcChkl4AvBK4YOD8rAnlJ6Sba++IiDU5r3MncCFw4jDCXQy8oOHxdtm+VnpI7dJTWxzbm9SevkjSg6Smr3dIunGQst6TlXdRdv5CUuIeaC75PbCHpGmDPP8u0k3dYf++rNqcuCssIhYAPyf1HhnYdw/wI+CcrNvYOEnjJR3WqVvfBvos8D5Jn5G0JYCknSV17Ekh6UuSXjEQI6knzBOkBENE3AQ8QmquuSQinmh4+o9JvVHeEhErhhjzV4GjgMlDfN6Ac4AvStpK0hRSLf5n2WvaT9Kuknqz9uvvkL6ZzG9RzhzgRaS28l1ItfPfAG8a5LrvzWLfpWF7B/BmSVtGxO+By4D/lLS7pDGSNpV0jKT3Z99iPgl8KbuvMElSj6TXSJqzge+FVZATd/V9jVRLbXQ88ENSF7MnSF3s3kbqVTHgs3p23+olG3LxiLgaeH22LZT0GCkhXZzn6cDppPb0xcB+wJsj4umGc84hNSf8x8COrAb+YVLierDhNRyRM+a/kppYmt+3vP6FdIP4VuA24MZsH6Q/BucAT5Le9xcD+0fEyiz2EyX9NotjeUQ8OLCRbjquzNrNn0XSnqTa+SmNz4mIucAC/t6kdCjpvf95FsM8YCapNk5EnA+8E3g/6T1/KIv9Vxv4XlgFKf2RNjOzunCN28ysZgpL3Nkw3IclzWvYt4Wky7LhupdlfXfNzEatVrmw6bgkfT/rLnurpN06lVlkjfsMoHkY9Szg8oiYQeodUcTNNDOzKjmD9XNhowNII6JnAEeTbsy3VVjijoirSAM8Gh0CDAwEORN4a1HXNzOrgkFyYaNDgLMiuQaYLKltX/wx3Qwwh20i4gGAiHig1aCFAZKOJv31YZNe7f7SSRt1P5o02rrror+YG74q8vtRT73ei3Vri7up3juumPdi7cp1hZTbO7HA/8b9xcQc/YUUC8CNS1cuiYjmEcdD8mJtEstZ2/G8B1h1O2nqgwFzImKoXS+nkvrfD+jL9g06AnakE3du2YufAzBzi43j2je8uPsXGVPMf9D+R3ONFRmyMZv2FlIuAOOK+atQ1Hux6qni/udv8oLxhZT71N3LCyl30l4F3ip6ZFUhxfY/Vtzvb/xld97X+az2lrOWDz9rDFZrJ3H3yogY7mC4Vomobc1kpHuVPDTwFSD79+ERvr6ZWUciJcdOW5f0kSYHGzCNwUfqQnevnctc/j7H8JF4UICZVZBIzRGdti6ZC7w3612yJ/DkQJPyYAprKpF0DmmuhimS+kjzbnyDNCH+B0jzBv9TUdc3MxuObtVqB8mFYwEiYjZpJOyBpBGyy0nTNbRVWOKOiMMHObRvUdc0M+uWbiXuNrlw4HiQFuTIrbI3J83MyiJa3zGsCiduM7MWqjwfiBO3mVkLTtxmZjUy0Kukqqocm5lZKQb6cVeVE7eZWQtO3GZmNePEbWZWI24qMTOrGd+cNDOrIde4zcxqxiMnzcxqxG3cZmY15MRtZlYjrnGbmdVQlZNjlWMzMyuFa9xmZjXkxD1cY3vQ1O6vvB2PrO56mQBrV68rptxHiykXYN2atotKb7DejYr5+Pf0FtdZa+mCFYWUO+nlEwspd+2CZYWUC7BmeTGfuRWPrimk3G5xjdvMrIacuM3MakRAb9lBtOHEbWbWgmvcZmY148RtZlYjvjlpZlZDytNxqZjOWB05cZuZtdCjHFnZidvMrBoE9FR4XlcnbjOzFpSnxl0SJ24zs2bK2cZdEiduM7MWnLjNzGoktXG7qcTMrFYqXOF24jYza6WnwiNwnLjNzJpI4V4lZmZ1437cZmY1U+VeJaW04kj6hKTbJc2TdI6k7i9vY2Y2DCI6brnKkfaXdJekBZJmtTi+maSLJN2S5cWjOpU54olb0lTgeGBmRLycNF/5YSMdh5nZYES6Odlp61iO1AucAhwA7AgcLmnHptM+BtwRETsDewPfljSuXbll3TcdA2wsaQwwAVhcUhxmZutT6sfdacthD2BBRCyMiNXAucAhTecEsKkkAROBx4D+doWOeBt3RPxN0reARcAK4NKIuLT5PElHA0cDbDdhLPG3lV2PRTM26XqZAGMeLWYR4iIXWF31ZNvPyQbrHVdM3WD8FmMLKRegZ2wxjZvL7lleSLmbvKyYzzHA0zctLaTcKrcfD8gZ4xRJ1zc8nhMRcxoeTwXub3jcB7yyqYwfAnNJFdhNgXdGRNtVmstoKtmc9BdnOvB8YBNJ724+LyLmRMTMiJi51UZVXv3NzEYb5dyAJQN5KtvmtCiqWXNV/U3AzaR8uAvwQ0mT2sVXRlPJG4C/RsQjEbEGuBB4dQlxmJkNaqAvd7sthz5g24bH01i/afgo4MJIFgB/BV7artAyEvciYE9JE7I2nX2B+SXEYWY2KKnzlsN1wAxJ07MbjoeRmkUaLSLlQSRtA+wALGxXaBlt3NdKOh+4kdQAfxPQ/PXCzKw0EvT2DH/kZET0SzoWuITUg+60iLhd0jHZ8dnAPwNnSLqN1LTyuYhY0q7cUgbgRMRXgK+UcW0zszy6df80Ii4GLm7aN7vh58XAG4dSpkdOmpm1UOWeL07cZmbr8SRTZma14sWCzcxqyE0lZmZ1Ii9dZmZWO65xm5nViABVuJHbidvMrJlAXnPSzKxe3FRiZlYnEhpT3cztxG1m1oIqXOV24jYza5JuTpYdxeByJW5JWwN7kSb6XgHMA67vtEqDmVlt1bXGLWkfYBawBWn61YeB8cBbgRdl07N+OyKeKjpQM7MRU/NeJQcCH4qIRc0HsoV+DwL2Ay4oIDYzs9LUth93RHymzbF+4Jddj8jMrGyCngovdZu3jfsE4HRgKXAqsCswq9Xq7EWItbBu2dqul9tz19NdLxOgd9r4QsodU8B7MODJvmJuVxR1G2TtmuLmkZg0baNCyl29tJjf3xPXFddSOWm7Yj7L6wr8/TGvS+VUuMadtxXn/Vk79huBrUiLW36jsKjMzEokurbmZCHydgccCPFA4PSIuEVV7uRoZjYcUn3buBvcIOlSYDrweUmbAu4KaGajVp17lQz4ALALsDAilkvaktRcYmY2KlW5USFX4o6IdZIeAnbMugGamY1aEmgU9Cr5N+CdwB3AwK3xAK4qKC4zs1KNhjbutwI7RMSqIoMxM6uEmo+cHLAQGAs4cZvZc0Nd27gl/YDUJLIcuFnS5TQk74g4vtjwzMzKUeca9/XZvzcAc5uOVXcJZDOzYZCgp7emNe6IOBPSkPeIOLnxWDYM3sxsFKr2AJy8XwaObLHvfV2Mw8ysOpRzK0mnNu7DgXcB0yU1NpVsCjxaZGBmZmWqcxv31cADwBTg2w37lwK3FhWUmVnZqtxU0qmN+z7gPuBVIxOOmVkFlDz7Xye5vgxIerukeyQ9KekpSUslebkyMxuVBGiMOm5lyTsA59+Bt0TE/CKDMTOrBDEqFlJ4qJtJW9JkSedLulPSfEluijGzaunJseUgaX9Jd0laIGnWIOfsLelmSbdLurJTmXlr3NdL+jlpjcnGkZMX5nx+s5OB30XEoZLGARM2sBwzs+7rUo1bUi9wCmlR9T7gOklzI+KOhnMmAz8C9o+IRZK27lRu3sQ9iTTs/Y0N+wIYcuKWNAl4HVk/8IhYDaweajlmZoXqTnfAPYAFEbEQQNK5wCGkmVYHvAu4MCIWAUTEw50KzTsfdzcXTXgh8AhwuqSdScPpT4iIZY0nSToaOBpgu/GeAtzMRpAEY3Jl7imSrm94PCci5jQ8ngrc3/C4D3hlUxkvAcZKuoI0RubkiDir3UXzzsc9DfgBsBeppv0nUrLty/P8FtfcDTguIq6VdDIwC/hS40nZi58DsPtmG8e6/gKmRllRzOprPT39hZQ7dpPiZnYfO66Ycpc+XcwohrWPF7fi/br+FYWUO27TYiog69YWN23Q4hufLqTcaa/drJByuyrfR3dJRMxsc7xVe0vzL2wMsDuwL7Ax8GdJ10TE3cMLDU4nTTL1fNJfkIuyfRuiD+iLiGuzx+eTErmZWTUMtHF32jrrA7ZteDwNWNzinN9FxLKIWEJaoGbndoXmTdxbRcTpEdGfbWcAW+V87rNExIPA/ZJ2yHbty7Pbe8zMytedXiXXATMkTc86YhzG+jOt/gp4raQxkiaQmlLa9uLL+91tiaR3A+dkjw9neHOVHAecnb2QhXjhYTOrki71KomIfknHApcAvcBpEXG7pGOy47MjYr6k35GmEVkHnBoR89qVmzdxvx/4IfBdUvvM1dm+DRIRNwPt2oXMzMrVpfE3EXExcHHTvtlNj78JfDNvmXl7lSwCDs5bqJlZreXvVVKKTtO6fr/dcS9dZmajVnXzdsca9zHAPOA80p3Q6g7eNzPrpgrPVdIpcT8P+CfgnUA/8HPggoh4vOjAzMxKo2ovpNA2tIh4NLvruQ9piPpk4HZJ7xmJ4MzMStOdftyFyDtycjdSF8D9gN+ShqmbmY1Oor5t3JK+ChxE6gx+LvD5iChmPLeZWZXUtVcJaf6QhaThlzsDX1daz0dARMROxYZnZlaCOte4gekjEoWZWaWU24bdSafEvSgi2k49JkmdzjEzq5WK17g7hfZHScdJ2q5xp6Rxkl4v6UzgyOLCMzMrSY17lexPmpPkHEnTgSeA8aTJUi4FvpvNO2JmNnoIKHEV907aJu6IWElaC+1HksYCU4AVEfHESARnZlaaCjeV5F6SIyLWAA8UGIuZWTV0aVrXongxRzOzVpy4zcxqpsJNJblDk/QCSW/Ift5Y0qbFhWVmViLl6FFS4V4lAEj6EHA0sAXwItKCl7NJ60UWT9A7tvtv0uplxawU3lvEivTAysfXFFIuQFE98dcVVO7yVcWteL98VTHl9j5WzOdty62LSyAqqOh7Ll9aTMHd1FvdppK8Ne6PAXsBTwFExD3A1kUFZWZWKuXcSpK3jXtVRKzO5ilB0hjS2pNmZqNTUV83uiBvjftKSScCG0vaD/gFcFFxYZmZlazCNe68iXsW8AhwG/Bh0orFXywqKDOz0kmdt5LkXeV9HfD/ss3MbHSr+SRTAEg6SNJNkh6T9JSkpZKeKjo4M7PS1L07IPA94O3AbZ7C1cxGv3KbQjrJm7jvB+Y5aZvZc0Z183buxP1Z4GJJVwLPDE+IiO8UEpWZWZnEqKhx/yvwNGku7nHFhWNmVhEVvjmZN3FvERFvLDQSM7MqqfDsgHn/pvxekhO3mT03DDSVVLQf91DmKvmdpBXuDmhmzwUVztu5B+B4Clcze26p681JSS+NiDsl7dbqeETcWExYZmYlq27e7ljj/iRpHu5vtzgWwOu7HpGZWdlU7sjITjol7tMAImKfEYjFzKw6Kpy4O92c/FFRF5bUm81/8uuirmFmtsG6NK2rpP0l3SVpgaRZbc57haS1kg7tVGaZXcxPAOaXeH0zs9a61B1QUi9wCnAAsCNwuKQdBznv34BL8oTXqankhZLmDnYwIg7Oc5FmkqYBbyaNyPzkhpRhZlao7rSU7AEsiIiFAJLOBQ4B7mg67zjgAuAVeQrtlLgfofWNyeH6Hmn+k0G7GUo6mnRjlOepl7v+a0XXg9j2JXkHjg7NRpuPLaTc1fevLKRcgCeXFvPlq7eg73Rjeoub7+zJZcV8LtZGMW2max4oZhFigInjiym3yMWeuyZfG/cUSdc3PJ4TEXMaHk8lTdI3oA94ZWMBkqYCbyN19uhK4l4aEVfmKSgvSQcBD0fEDZL2Huy87MXPAfiHMRt5VkIzG1n5/s4uiYiZQyylOZ99D/hcRKxVzr7jnRL3vblKGZq9gIMlHUiatGqSpJ9FxLsLuJaZ2dCJbvUq6QO2bXg8DVjcdM5M4NwsaU8BDpTUHxG/HKzQtl9kI+LtGxZr2zI/HxHTImJ74DDgD07aZlYtOW5M5qsdXwfMkDRd0jhSznvWfcOImB4R22c58Xzgo+2SNuSfHdDM7LmlCxXuiOiXdCypt0gvcFpE3C7pmOz47A0pt9TEHRFXAFeUGYOZ2Xq611RCRFwMXNy0r2XCjoj35Smz01wlLecoabiI5yoxs9GprpNM8feugONJDei3kP4W7QRcC7ymuNDMzErUU90lcDrdnNwnm6fkPmC3iJgZEbsDuwILRiJAM7ORJ1BP560kedu4XxoRtw08iIh5knYpKCYzs3J1sY27CHkT93xJpwI/I3UefzeeZ8TMRrMat3EPOAr4CGliKICrgB8XEpGZWelUalNIJ3mXLlsJfDfbzMxGv7rXuCXNAP4vaVrCZ6adiYgXFhSXmVl5JOit7kRYeb8LnE5qGukH9gHOAn5aVFBmZqWr8DLveRP3xhFxOaCIuC8iTsLrTZrZaFbhxJ335uRKST3APdm4+78BWxcXlplZiUSlb07mjezjwATgeGB3UnfAI4sKysysXNkq7522kuTtVXIdgKSIiKOKDcnMrAJ6an5zUtKrJN1BNuhG0s6SClsB3sysVF1aLLgoeZtKvge8CXgUICJuAV5XVFBmZuUaBU0lABFxf9N6aMWtUGpmVrYK35zMm7jvl/RqILLld45nBOcqWbtOPFHAytuP3VDMX8zpfU8WUu5GmxW37sWqNcV8SItajX3ShOLqDU8tL6bcJ1cV8x4X9bsDWFFQzP3rqjsq8RkVHjmZ97dyDPAx0lLzfcAuwEeLCsrMrFQVb+PO26tkCXBE4z5JHye1fZuZjTKjY8h7K5/sWhRmZlVT9xr3IKrbAGRmNhwDTSUVNZzEXcxdJzOz0tV4Pm5JS2mdoAVsXEhEZmZVUNelyyJi05EKxMysUkZpU4mZ2egkVXquEiduM7NWXOM2M6uZnprenDQze24SVe7x7MRtZtZsFPfjNjMbveSbk2ZmNVLukPZOnLjNzFqp68hJM7PnLte4zczqpcJNJSP+XUDStpL+KGm+pNslnTDSMZiZtSdSeuy05ShJ2l/SXZIWSJrV4vgRkm7Ntqsl7dypzDJq3P3ApyLiRkmbAjdIuiwi7ighFjOz9YmuDMCR1AucAuxHWj3sOklzm/LdX4F/jIjHJR0AzAFe2a7cEa9xR8QDEXFj9vNS0tqVU0c6DjOz9pRj62gPYEFELIyI1cC5wCGNJ0TE1RHxePbwGmBap0JLbeOWtD2wK3Bti2NHA0cDbOOmeDMbUbnn454i6fqGx3MiYk7D46nA/Q2P+2hfm/4A8NtOFy0tI0qaCFwAfDwinmo+nr34OQAv0fhYWsBK1sWsxQ56eHwh5a57qLibJRuPW1dIuRM26i+k3ChwGY/NJ64ppNxlq4oZ0PFYge/FqoJiXltIqV2W7+bkkoiY2a6UFvta/sYk7UNK3K/pdNFSEreksaSkfXZEXFhGDGZm7XWlotQHbNvweBqweL0rSTsBpwIHRMSjnQoto1eJgJ8A8yPiOyN9fTOzzrKmkk5bZ9cBMyRNlzQOOAyY+6wrSdsBFwLviYi78xRaRo17L+A9wG2Sbs72nRgRF5cQi5nZ+gTqwsjJiOiXdCxwCdALnBYRt0s6Jjs+G/gysCXwo1Svpb9D88vIJ+6I+BNVHpJkZgZdG/KeVUovbto3u+HnDwIfHEqZ7q5hZrYez8dtZlY/FR7y7sRtZtaKZwc0M6uT3ANwSuHEbWbWTDhxm5nVj9u4zczqxTcnzczqZGA+7mpy4jYza8U1bjOzOhGomJkRu8GJ28ysFde4zczqxonbzKw+5AE4Zmb146YSM7O6ceI2M6sR9yoxM6sfN5UMzypgUdlBDMHq/mJuakwopNRk3YpiahfjCip3u4nFrB4Pxa0gv9mEYmJ+Yllx/43XW9W2S5YXVG53+eakmVl9CNe4zczqxd0Bzczqx4nbzKxOPDugmVn9uI3bzKxuXOM2M6sX17jNzGpEcuI2M6sdD3k3M6sb17jNzGrEA3DMzGrINW4zs3pxjdvMrE6Ea9xmZnUioMe9SszMasY1bjOzGqn2AJxSWt8l7S/pLkkLJM0qIwYzs/Z6cmyddcp3Sr6fHb9V0m55IhtRknqBU4ADgB2BwyXtONJxmJm1NTDsvd3WsYhc+e4AYEa2HQ38uFO5ZdS49wAWRMTCiFgNnAscUkIcZmaDyFZ577R1liffHQKcFck1wGRJz2tXaBlt3FOB+xse9wGvbD5J0tGkvz4Aqz7M3fNGILZumgIsKTuIIahXvE8DdYs5qVvMdYsXYIfhFnDDDfMvkWZOyXHqeEnXNzyeExFzGh7nyXetzpkKPDDYRctI3K2+X6y3rnb24ucASLo+ImYWHVg31S3musULjnkk1C1eSDEPt4yI2L8bsZAv3+XKiY3KaCrpA7ZteDwNWFxCHGZmRcuT74acE8tI3NcBMyRNlzQOOAyYW0IcZmZFy5Pv5gLvzXqX7Ak8GRGDNpNACU0lEdEv6VjgEqAXOC0ibu/wtDkdjldR3WKuW7zgmEdC3eKFCsU8WL6TdEx2fDZwMXAgsABYDhzVqVxFtG1KMTOziqnu9FdmZtaSE7eZWc1UOnHXbWi8pG0l/VHSfEm3Szqh7JjykNQr6SZJvy47ljwkTZZ0vqQ7s/f6VWXH1ImkT2SfiXmSzpE0vuyYmkk6TdLDkuY17NtC0mWS7sn+3bzMGJsNEvM3s8/GrZL+U9LkMmMsQmUTd02HxvcDn4qIlwF7Ah+rQcwAJwDzyw5iCE4GfhcRLwV2puKxS5oKHA/MjIiXk25SHVZuVC2dATT3X54FXB4RM4DLs8dVcgbrx3wZ8PKI2Am4G/j8SAdVtMombmo4ND4iHoiIG7Ofl5ISytRyo2pP0jTgzcCpZceSh6RJwOuAnwBExOqIeKLcqHIZA2wsaQwwgQqOXYiIq4DHmnYfApyZ/Xwm8NYRDaqDVjFHxKUR0Z89vIbUL3pUqXLiHmwYaC1I2h7YFbi23Eg6+h7wWWBd2YHk9ELgEeD0rHnnVEmblB1UOxHxN+BbwCLSMOYnI+LScqPKbZuBPsXZv1uXHM9QvR/4bdlBdFuVE/eQh4FWhaSJwAXAxyPiqbLjGYykg4CHI+KGsmMZgjHAbsCPI2JXYBnV+/r+LFm78CHAdOD5wCaS3l1uVKOfpC+Qmi/PLjuWbqty4q7l0HhJY0lJ++yIuLDseDrYCzhY0r2kpqjXS/pZuSF11Af0RcTAN5nzSYm8yt4A/DUiHomINcCFwKtLjimvhwZmqsv+fbjkeHKRdCRwEHBEjMLBKlVO3LUbGi9JpLbX+RHxnbLj6SQiPh8R0yJie9L7+4eIqHRNMCIeBO6XNDAD3L7AHSWGlMciYE9JE7LPyL5U/IZqg7nAkdnPRwK/KjGWXCTtD3wOOPQmypoAAAK4SURBVDgilpcdTxEqm7izmwsDQ0XnA+flGBpftr2A95Bqrjdn24FlBzUKHQecLelWYBfg6yXH01b27eB84EbgNtL/u8oMyx4g6Rzgz8AOkvokfQD4BrCfpHuA/bLHlTFIzD8ENgUuy/4Pzi41yAJ4yLuZWc1UtsZtZmatOXGbmdWME7eZWc04cZuZ1YwTt5lZzThx27BJejrHOacOTLgl6cSmY1dv6DUkhaSfNjweI+mRDZ3pMJt58KMNj/euy6yJ9tzhxG0jIiI+GBEDA2VObDo2nFGEy4CXS9o4e7wf8LdhlDcZ+GjHs8xK5MRtXZPVTq9omCv77GykINn+mZK+QZol72ZJZ2fHns7+nSjpckk3SrpNUt7ZIH9LmuEQ4HDgnIaYtpD0y2xu5msk7ZTtPymby/kKSQslHZ895RvAi7L4vpntm9jqNZmVxYnbum1X4OOkOdRfSBpN+oyImAWsiIhdIuKIpueuBN4WEbsB+wDfzpkkzwUOyxYn2Ilnz8j4VeCmbG7mE4GzGo69FHgTaQrhr2TzzMwC/juL7zN5XpPZSHPitm77S0T0RcQ64GZg+yE8V8DXs6HsvydN47tNpydFxK3ZdQ4nrZjd6DXAT7Pz/gBsKWmz7NhvImJVRCwhTZ402LWG85rMum5M2QHYqLOq4ee1DO0zdgSwFbB7RKzJZi3Mu8TXXNKc13sDWzbsbzc9cN5Yh/OazLrONW4rw5qsWaLZZqT5wddI2gd4wRDKPA34WkTc1rT/KtIfBCTtDSzpMEf6UtIERWaV5ZqDlWEOcKukG5vauc8GLpJ0PalJ4s68BUZEH2ktymYnkVbLuRVYzt+nKB2snEcl/Ve2+Oxvgd/kjcFspHh2QDOzmnFTiZlZzThxm5nVjBO3mVnNOHGbmdWME7eZWc04cZuZ1YwTt5lZzfwP7cuAvoku22MAAAAASUVORK5CYII=\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" } ], "source": [ "plt.pcolormesh(skill,cmap=plt.cm.YlOrRd,vmin=0.0,vmax=1.0)\n", "plt.colorbar()\n", "plt.title('NCEP-CFSv2 Nino3.4 ACC')\n", "plt.xlabel('Initial Month')\n", "plt.ylabel('Lead Time (Months)')" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### References\n", "\n", "1. Kirtman, B.P., D. Min, J.M. Infanti, J.L. Kinter, D.A. Paolino, Q. Zhang, H. van den Dool, S. Saha, M.P. Mendez, E. Becker, P. Peng, P. Tripp, J. Huang, D.G. DeWitt, M.K. Tippett, A.G. Barnston, S. Li, A. Rosati, S.D. Schubert, M. Rienecker, M. Suarez, Z.E. Li, J. Marshak, Y. Lim, J. Tribbia, K. Pegion, W.J. Merryfield, B. Denis, and E.F. Wood, 2014: The North American Multimodel Ensemble: Phase-1 Seasonal-to-Interannual Prediction; Phase-2 toward Developing Intraseasonal Prediction. Bull. Amer. Meteor. Soc., 95, 585–601, https://doi.org/10.1175/BAMS-D-12-00050.1" ] } ], "metadata": { "kernelspec": { "display_name": "Python (aoes)", "language": "python", "name": "aoes" }, "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.7" }, "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": false, "toc_position": {}, "toc_section_display": true, "toc_window_display": false }, "varInspector": { "cols": { "lenName": 16, "lenType": 16, "lenVar": 40 }, "kernels_config": { "python": { "delete_cmd_postfix": "", "delete_cmd_prefix": "del ", "library": "var_list.py", "varRefreshCmd": "print(var_dic_list())" }, "r": { "delete_cmd_postfix": ") ", "delete_cmd_prefix": "rm(", "library": "var_list.r", "varRefreshCmd": "cat(var_dic_list()) " } }, "types_to_exclude": [ "module", "function", "builtin_function_or_method", "instance", "_Feature" ], "window_display": false } }, "nbformat": 4, "nbformat_minor": 4 }