{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "You can run this notebook in a [live session](https://binder.pangeo.io/v2/gh/pangeo-data/climpred/main?urlpath=lab/tree/docs/source/examples/subseasonal/seasonal-enso-subx-example.ipynb) [binder badge](https://binder.pangeo.io/v2/gh/pangeo-data/climpred/main?urlpath=lab/tree/docs/source/examples/subseasonal/seasonal-enso-subx-example.ipynb) or view it [on Github](https://github.com/pangeo-data/climpred/blob/main/docs/source/examples/subseasonal/seasonal-enso-subx-example.ipynb)." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Calculate Seasonal ENSO Skill of the NMME model NCEP-CFSv2" ] }, { "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 seasonal data" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### The North American Multi-model Ensemble (NMME)\n", "\n", "Further information on NMME is available from {cite:t}`Kirtman2014` 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/](http://iridl.ldeo.columbia.edu/SOURCES/.Models/.NMME/)." ] }, { "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": {}, "outputs": [ { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 1;\n", " var nbb_unformatted_code = \"# linting\\n%load_ext nb_black\\n%load_ext lab_black\";\n", " var nbb_formatted_code = \"# linting\\n%load_ext nb_black\\n%load_ext lab_black\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "# linting\n", "%load_ext nb_black\n", "%load_ext lab_black" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 2;\n", " var nbb_unformatted_code = \"import matplotlib.pyplot as plt\\nimport xarray as xr\\nimport pandas as pd\\nimport numpy as np\\n\\nfrom climpred import HindcastEnsemble\\nimport climpred\";\n", " var nbb_formatted_code = \"import matplotlib.pyplot as plt\\nimport xarray as xr\\nimport pandas as pd\\nimport numpy as np\\n\\nfrom climpred import HindcastEnsemble\\nimport climpred\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "import matplotlib.pyplot as plt\n", "import xarray as xr\n", "import pandas as pd\n", "import numpy as np\n", "\n", "from climpred import HindcastEnsemble\n", "import climpred" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 3;\n", " var nbb_unformatted_code = \"def decode_cf(ds, time_var):\\n \\\"\\\"\\\"Decodes time dimension to CFTime standards.\\\"\\\"\\\"\\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\";\n", " var nbb_formatted_code = \"def decode_cf(ds, time_var):\\n \\\"\\\"\\\"Decodes time dimension to CFTime standards.\\\"\\\"\\\"\\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\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "def decode_cf(ds, time_var):\n", " \"\"\"Decodes time dimension to CFTime standards.\"\"\"\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": [ "The original monthly sea surface temperature (SST) hindcast data for the NCEP-CFSv2 model from IRIDL is a large dataset (~20GB). However, we can leverage `ingrid` on the IRIDL server and open server-side preprocessed data via `OpenDAP` into `xarray`. Averaging over longitude `X` and latitude `Y` and ensemble member `M` reduces the download size to just a few kB.\n", "\n", "- [http://iridl.ldeo.columbia.edu/dochelp/topics/DODS/fnlist.html](http://iridl.ldeo.columbia.edu/dochelp/topics/DODS/fnlist.html)\n", "- [https://iridl.ldeo.columbia.edu/dochelp/Documentation/funcindex.html?Set-Language=en](https://iridl.ldeo.columbia.edu/dochelp/Documentation/funcindex.html?Set-Language=en)\n", "\n", "We take the ensemble mean here, since we are just using deterministic metrics in this example. `climpred` will automatically take the mean over the ensemble while evaluating metrics with `comparison=\"e2o\"`, but this should make things more efficient so it doesn't have to be done multiple times." ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "CPU times: user 62.9 ms, sys: 17.9 ms, total: 80.8 ms\n", "Wall time: 800 ms\n" ] }, { "data": { "text/html": [ "
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
<xarray.Dataset>\n",
       "Dimensions:  (S: 348, L: 10)\n",
       "Coordinates:\n",
       "  * S        (S) object 1982-01-01 00:00:00 ... 2010-12-01 00:00:00\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",
       "Data variables:\n",
       "    sst      (S, L) float64 26.36 26.63 27.15 27.96 ... 27.93 27.18 26.31 25.96\n",
       "Attributes:\n",
       "    Conventions:  IRIDL
" ], "text/plain": [ "\n", "Dimensions: (S: 348, L: 10)\n", "Coordinates:\n", " * S (S) object 1982-01-01 00:00:00 ... 2010-12-01 00:00:00\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", "Data variables:\n", " sst (S, L) float64 26.36 26.63 27.15 27.96 ... 27.93 27.18 26.31 25.96\n", "Attributes:\n", " Conventions: IRIDL" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" }, { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 4;\n", " var nbb_unformatted_code = \"%%time\\n# Get NMME data for NCEP-CFSv2, SST\\nurl = 'http://iridl.ldeo.columbia.edu/SOURCES/.Models/.NMME/NCEP-CFSv2/.HINDCAST/.MONTHLY/.sst/X/190/240/RANGEEDGES/Y/-5/5/RANGEEDGES/[X%20Y%20M]average/dods'\\nfcstds = xr.open_dataset(url, decode_times=False)\\nfcstds = decode_cf(fcstds, 'S').compute()\\nfcstds\";\n", " var nbb_formatted_code = \"%%time\\n# Get NMME data for NCEP-CFSv2, SST\\nurl = 'http://iridl.ldeo.columbia.edu/SOURCES/.Models/.NMME/NCEP-CFSv2/.HINDCAST/.MONTHLY/.sst/X/190/240/RANGEEDGES/Y/-5/5/RANGEEDGES/[X%20Y%20M]average/dods'\\nfcstds = xr.open_dataset(url, decode_times=False)\\nfcstds = decode_cf(fcstds, 'S').compute()\\nfcstds\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "%%time\n", "# Get NMME data for NCEP-CFSv2, SST\n", "url = 'http://iridl.ldeo.columbia.edu/SOURCES/.Models/.NMME/NCEP-CFSv2/.HINDCAST/.MONTHLY/.sst/X/190/240/RANGEEDGES/Y/-5/5/RANGEEDGES/[X%20Y%20M]average/dods'\n", "fcstds = xr.open_dataset(url, decode_times=False)\n", "fcstds = decode_cf(fcstds, 'S').compute()\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`. `climpred` renames the dimensions based on their `attrs` `standard_name` when creating {py:class}`.HindcastEnsemble`, but we need to first adapt the coordinates manually." ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 5;\n", " var nbb_unformatted_code = \"fcstds = fcstds.rename({\\\"S\\\": \\\"init\\\", \\\"L\\\": \\\"lead\\\"})\";\n", " var nbb_formatted_code = \"fcstds = fcstds.rename({\\\"S\\\": \\\"init\\\", \\\"L\\\": \\\"lead\\\"})\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "fcstds = fcstds.rename({\"S\": \"init\", \"L\": \"lead\"})" ] }, { "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." ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 6;\n", " var nbb_unformatted_code = \"fcstds[\\\"lead\\\"] = (fcstds[\\\"lead\\\"] - 0.5).astype(\\\"int\\\")\";\n", " var nbb_formatted_code = \"fcstds[\\\"lead\\\"] = (fcstds[\\\"lead\\\"] - 0.5).astype(\\\"int\\\")\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "fcstds[\"lead\"] = (fcstds[\"lead\"] - 0.5).astype(\"int\")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Next, we want to get the verification SST data from the data server" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "CPU times: user 11.3 ms, sys: 9.96 ms, total: 21.2 ms\n", "Wall time: 530 ms\n" ] }, { "data": { "text/html": [ "
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
<xarray.Dataset>\n",
       "Dimensions:  (T: 405)\n",
       "Coordinates:\n",
       "  * T        (T) object 1982-01-16 00:00:00 ... 2015-09-16 00:00:00\n",
       "Data variables:\n",
       "    sst      (T) float64 26.81 26.78 27.26 28.06 ... 28.96 28.83 28.89 28.99\n",
       "Attributes:\n",
       "    Conventions:  IRIDL
" ], "text/plain": [ "\n", "Dimensions: (T: 405)\n", "Coordinates:\n", " * T (T) object 1982-01-16 00:00:00 ... 2015-09-16 00:00:00\n", "Data variables:\n", " sst (T) float64 26.81 26.78 27.26 28.06 ... 28.96 28.83 28.89 28.99\n", "Attributes:\n", " Conventions: IRIDL" ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" }, { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 7;\n", " var nbb_unformatted_code = \"%%time\\nobsurl='http://iridl.ldeo.columbia.edu/SOURCES/.Models/.NMME/.OIv2_SST/.sst/X/190/240/RANGEEDGES/Y/-5/5/RANGEEDGES/[X%20Y]average/dods'\\nverifds = decode_cf(xr.open_dataset(obsurl, decode_times=False),'T').compute()\\nverifds\";\n", " var nbb_formatted_code = \"%%time\\nobsurl='http://iridl.ldeo.columbia.edu/SOURCES/.Models/.NMME/.OIv2_SST/.sst/X/190/240/RANGEEDGES/Y/-5/5/RANGEEDGES/[X%20Y]average/dods'\\nverifds = decode_cf(xr.open_dataset(obsurl, decode_times=False),'T').compute()\\nverifds\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "%%time\n", "obsurl='http://iridl.ldeo.columbia.edu/SOURCES/.Models/.NMME/.OIv2_SST/.sst/X/190/240/RANGEEDGES/Y/-5/5/RANGEEDGES/[X%20Y]average/dods'\n", "verifds = decode_cf(xr.open_dataset(obsurl, decode_times=False),'T').compute()\n", "verifds" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Rename the dimensions to correspond to `climpred` dimensions." ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 8;\n", " var nbb_unformatted_code = \"verifds = verifds.rename({\\\"T\\\": \\\"time\\\"})\";\n", " var nbb_formatted_code = \"verifds = verifds.rename({\\\"T\\\": \\\"time\\\"})\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "verifds = verifds.rename({\"T\": \"time\"})" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Convert the `time` data to be on the first of the month and in the same calendar format as the forecast output. The time dimension is natively decoded to start in February, even though it starts in January. It is also labelled as mid-month, and we need to label it as the start of the month to ensure that the dates align properly." ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 9;\n", " var nbb_unformatted_code = \"verifds[\\\"time\\\"] = xr.cftime_range(\\n start=\\\"1982-01\\\", periods=verifds[\\\"time\\\"].size, freq=\\\"MS\\\", calendar=\\\"360_day\\\"\\n)\";\n", " var nbb_formatted_code = \"verifds[\\\"time\\\"] = xr.cftime_range(\\n start=\\\"1982-01\\\", periods=verifds[\\\"time\\\"].size, freq=\\\"MS\\\", calendar=\\\"360_day\\\"\\n)\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "verifds[\"time\"] = xr.cftime_range(\n", " start=\"1982-01\", periods=verifds[\"time\"].size, freq=\"MS\", calendar=\"360_day\"\n", ")" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Subset the data to 1982-2010" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 10;\n", " var nbb_unformatted_code = \"verifds = verifds.sel(time=slice(\\\"1982-01-01\\\", \\\"2010-12-01\\\"))\\nfcstds = fcstds.sel(init=slice(\\\"1982-01-01\\\", \\\"2010-12-01\\\"))\";\n", " var nbb_formatted_code = \"verifds = verifds.sel(time=slice(\\\"1982-01-01\\\", \\\"2010-12-01\\\"))\\nfcstds = fcstds.sel(init=slice(\\\"1982-01-01\\\", \\\"2010-12-01\\\"))\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "verifds = verifds.sel(time=slice(\"1982-01-01\", \"2010-12-01\"))\n", "fcstds = fcstds.sel(init=slice(\"1982-01-01\", \"2010-12-01\"))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Make Seasonal Averages with `center=True` and drop `NaN`s." ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 11;\n", " var nbb_unformatted_code = \"f = (\\n fcstds.rolling(lead=3, center=True)\\n .mean()\\n .dropna(\\\"lead\\\")\\n .isel(lead=slice(None, None, 3))\\n .assign_coords(lead=np.arange(3))\\n)\\no = verifds.rolling(time=3, center=True).mean().dropna(dim=\\\"time\\\")\";\n", " var nbb_formatted_code = \"f = (\\n fcstds.rolling(lead=3, center=True)\\n .mean()\\n .dropna(\\\"lead\\\")\\n .isel(lead=slice(None, None, 3))\\n .assign_coords(lead=np.arange(3))\\n)\\no = verifds.rolling(time=3, center=True).mean().dropna(dim=\\\"time\\\")\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "f = (\n", " fcstds.rolling(lead=3, center=True)\n", " .mean()\n", " .dropna(\"lead\")\n", " .isel(lead=slice(None, None, 3))\n", " .assign_coords(lead=np.arange(3))\n", ")\n", "o = verifds.rolling(time=3, center=True).mean().dropna(dim=\"time\")" ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [ { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 12;\n", " var nbb_unformatted_code = \"f[\\\"lead\\\"].attrs = {\\\"units\\\": \\\"seasons\\\"}\";\n", " var nbb_formatted_code = \"f[\\\"lead\\\"].attrs = {\\\"units\\\": \\\"seasons\\\"}\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "f[\"lead\"].attrs = {\"units\": \"seasons\"}" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Create a {py:class}`.HindcastEnsemble` object and {py:meth}`.HindcastEnsemble.verify` forecast skill." ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [ { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 13;\n", " var nbb_unformatted_code = \"hindcast = HindcastEnsemble(f).add_observations(o)\\nhindcast = hindcast.remove_seasonality(\\\"month\\\") # remove seasonal cycle\";\n", " var nbb_formatted_code = \"hindcast = HindcastEnsemble(f).add_observations(o)\\nhindcast = hindcast.remove_seasonality(\\\"month\\\") # remove seasonal cycle\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "hindcast = HindcastEnsemble(f).add_observations(o)\n", "hindcast = hindcast.remove_seasonality(\"month\") # remove seasonal cycle" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Compute the Anomaly Correlation Coefficient (ACC) {py:func}`.climpred.metrics._pearson_r` for 0, 1, and 2 season lead-times." ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [ { "data": { "text/html": [ "
\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "
<xarray.Dataset>\n",
       "Dimensions:  (lead: 3)\n",
       "Coordinates:\n",
       "  * lead     (lead) int64 0 1 2\n",
       "    skill    <U11 'initialized'\n",
       "Data variables:\n",
       "    sst      (lead) float64 0.8466 0.7671 0.6785\n",
       "Attributes:\n",
       "    prediction_skill_software:     climpred https://climpred.readthedocs.io/\n",
       "    skill_calculated_by_function:  HindcastEnsemble.verify()\n",
       "    number_of_initializations:     348\n",
       "    alignment:                     maximize\n",
       "    metric:                        pearson_r\n",
       "    comparison:                    e2o\n",
       "    dim:                           init\n",
       "    reference:                     []
" ], "text/plain": [ "\n", "Dimensions: (lead: 3)\n", "Coordinates:\n", " * lead (lead) int64 0 1 2\n", " skill " ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "skillds = hindcast.verify(\n", " metric=\"acc\", comparison=\"e2o\", dim=\"init\", alignment=\"maximize\"\n", ")\n", "skillds" ] }, { "cell_type": "code", "execution_count": 15, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "Text(0, 0.5, 'ACC')" ] }, "execution_count": 15, "metadata": {}, "output_type": "execute_result" }, { "data": { "image/png": "iVBORw0KGgoAAAANSUhEUgAAAYIAAAETCAYAAAA7wAFvAAAAOXRFWHRTb2Z0d2FyZQBNYXRwbG90bGliIHZlcnNpb24zLjUuMCwgaHR0cHM6Ly9tYXRwbG90bGliLm9yZy8/fFQqAAAACXBIWXMAAAsTAAALEwEAmpwYAAAZt0lEQVR4nO3df5xddX3n8debQGSR32R0IT9IxFCMaEFmQ1GkUGAbRBNWUBPtAl0h9mGDdv1RY7URKduH4hYrbmyNSrFUCAHRHSUaW35U2xXMUH4GCMwjBpKAMkFEsJYQeO8f5wwcbu5MZkLO3GTO+/l4nEfu+Z7vPedz753c9z2/ZZuIiGiuXTpdQEREdFaCICKi4RIEERENlyCIiGi4BEFERMMlCCIiGi5BENFhkp6S9KpO1xHNlSCIYZO0VtKjkl5eaTtH0k2VcUn6gKS7Jf1a0npJV0t6XTn9Mkmbyi+/geGOctpUSa60r5W0cCs1vVtSb9n/EUnfk3RsOe18Sc+0LOtPy2lzJN0u6VeSNkq6QdK0YbwHp0r6F0m/lPQzSV+VtNdLec9s72l7zdaWPYzaLpK0rnxND0r6s0H6nVm+z+cMY56XSdos6cA2035f0g8lPSmpX9I/S5pdmX6gpK+Vn8uTku6T9OnqexE7hgRBjNQ44INDTP9COf0DwP7AocC3gVMrfS4qv/wGht9umce+tvcE5gGLJM1qtyBJHwL+GvhL4JXAFOBLwJxKt6talnWRpFcDfw98GNgHmAYsBp7d6qsv+l8IHAS8BpgIfG4rz9nae7a9fA04zPbewBuB90h6e7WDpP2APwNWbW1m5Rf26cATwB+0TDsDuJrifZxE8f4vAt5WTt8f+DHwn4BjbO8FnAzsCxyyza8w6mE7Q4ZhDcBaYCHwC4ova4BzgJvKx9MpvkxnDjGPy4ALB5k2FTCwa6VtJfCRNn33AZ4C3jHEss4H/qFN+xnA7YM85yDgN8D+lbYjgY3Abm36vx24a1vfs3LcwKsr789i4DrgSeAW4JBK3zeW78kT5b9vHGS5E4G7gD9taf9b4P3ATcA5W/m8zwTWUYTY3ZV2AQ8BHx3iuReWy9+l03+3GbY+ZI0gRqqX4kvkI22mnQist/2Tl7qQchPTm4DXAre16XIMsDvwrW2Y/b8Bh0n6vKQTJO05MMH2wxS/ZE+v9H83cI3tZ9rM6zi2/ut6qPesnbnAp4H9gD7gf8Hzv7KvAy4BDgAuBq6TdMDAEyUtlPQUsB54OXBFZdpMoJsiDIbjLOBKYCnF+3VU2f5bwGTgmiGeexJwre3nhrms6KAEQWyLRcB5krpa2g8AHhnG8z9SbmMfGL7eMn0jxS/orwILbV/fZh4HABttb97Kst7ZsqyDXGyPP57iV/MyYGO5LXwgEK6g2CyFJFF8MV/ROmNJJ1N8WS4axmse7D1r51u2f1K+tm8AR5TtpwIP2L7c9mbbVwL3UW6OAbD9GWAv4A3A5RRrDkgaR7HZbMFwvpwlTQFOAK6w/XPgeoo1BCjeexj6sx7u30LsABIEMWK27wa+S7HJo+oxYIudim38b9v7VoazWqZPsL2f7dfYvgSg3Ak8sMP3PeWyJkjadSvLWtayrIfL13Cz7Xfa7gLeTPHL/hPlc74JHFPuID0OeA74UXWmkn6HIhzOsH3/1l7wEO9ZOz+rPP53YCCgDgIebOn7IEWgVZdl27dRbOL6dNn8fuBO2zcPY/kA/x241/bt5fg3gHdL2o3ivYehP+vh/i3EDiBBENvqU8C5vPhL6HpgkqTu7b0w26f4hR2+36DYfPM0cNp2mPdK4Frg8HL8ceAHwLsoNgsttf38ZXolHQn0AP9jkLWVwbR7z0biYeDglrYpwIZB+u/KCztmTwT+W3mk088o9jX8laT/M8hzzwReVel/MTABeAuwmmLfwemDPBfgn8rl5TtmJ5APKbaJ7T7gKoqjgwbaHqDY/HClpOMljZe0u6S5WzsMdBuW/wTF5pbFkk6TtIek3SSdIumioZ4r6VhJ50p6RTl+GDAbqP5avoLiy/AMXryd/XDg+8B5tr8zwpq3eM9GaDlwqIpDZneV9C5gBvBdSbtIep+k/cr9KzOBP6YIZ4CzKY5yOqIceinWFj7RsgwkHUMRIDMr/Q+nfE/KUPwQ8OeS/lDS3uXyj5W0pJzNxcDewNclHVzOd6KkiyW9fhtff9Sl03urM+w8A8URMCdVxicD/8GLj4ARxVEmqyg2a2yg+PJ7bTn9MmATxRE/A8PGctpUWo4aGkZN76H4Uvs1xSaV6yiPpGHwo4YOB74D/Lxc/lrgs1SOCqI47PFJYFXLc/+OYlNRtf5VL/E9az1q6MLKtOMpdsAPjB8L3Eqx7f9W4NiyfReKgPpFWdP9FIeJapC6bmKQo4YodiZ/s037TIq1sP3L8VkUm8yeAvrLeZ5a6X8QcGn5uTxJsT/jU8Aenf5bzvDiQeUHFhERDZVNQxERDZcgiIhouARBRETDJQgiIhouQRAR0XBbOytzhzNhwgRPnTq102VEROxUbr311o0uzqTfwk4XBFOnTqW3t7fTZURE7FQktV6e5HnZNBQR0XAJgoiIhksQREQ0XIIgIqLhEgQREQ2XIIiIaLgEQUREwyUIIiIabqc7oWy0TV14XadLqNXaz5za6RIiosOyRhAR0XAJgoiIhksQREQ0XIIgIqLhEgQREQ2XIIiIaLgEQUREw9UaBJJmSVotqU/SwjbTp0i6UdJtku6U9JY664mIiC3VFgSSxgGLgVOAGcA8STNaun0SWGb7SGAu8KW66omIiPbqXCOYCfTZXmN7E7AUmNPSx8De5eN9gIdrrCciItqo8xITE4F1lfH1wNEtfc4HfiDpPODlwEk11hMREW10emfxPOAy25OAtwCXS9qiJknzJfVK6u3v7x/1IiMixrI6g2ADMLkyPqlsq3ovsAzA9o+B3YEJrTOyvcR2t+3urq6umsqNiGimOoNgJTBd0jRJ4yl2Bve09HkIOBFA0msogiA/+SMiRlFtQWB7M7AAWAHcS3F00CpJF0iaXXb7MHCupDuAK4GzbbuumiIiYku13o/A9nJgeUvbosrje4A31VlDREQMrdM7iyMiosMSBBERDZcgiIhouNyzOMa0sXzP6dxvOraXrBFERDRcgiAiouESBBERDZcgiIhouARBRETDJQgiIhouQRAR0XAJgoiIhksQREQ0XIIgIqLhEgQREQ1XaxBImiVptaQ+SQvbTP+8pNvL4X5Jv6yznoiI2FJtF52TNA5YDJwMrAdWSuopb0YDgO3/Wel/HnBkXfVERER7da4RzAT6bK+xvQlYCswZov88ittVRkTEKKozCCYC6yrj68u2LUg6GJgG3FBjPRER0caOsrN4LnCN7WfbTZQ0X1KvpN7+/v5RLi0iYmyrMwg2AJMr45PKtnbmMsRmIdtLbHfb7u7q6tqOJUZERJ1BsBKYLmmapPEUX/Y9rZ0kHQbsB/y4xloiImIQtQWB7c3AAmAFcC+wzPYqSRdIml3pOhdYatt11RIREYOr9Z7FtpcDy1vaFrWMn19nDRERMbTcvD4idkhTF17X6RJqtfYzp3a6hOftKEcNRUREhyQIIiIaLkEQEdFwCYKIiIZLEERENFyCICKi4RIEERENlyCIiGi4BEFERMMlCCIiGi5BEBHRcAmCiIiGSxBERDRcgiAiouESBBERDVdrEEiaJWm1pD5JCwfp805J90haJemKOuuJiIgt1XZjGknjgMXAycB6YKWkHtv3VPpMBz4OvMn245JeUVc9ERHRXp1rBDOBPttrbG8ClgJzWvqcCyy2/TiA7UdrrCciItqoMwgmAusq4+vLtqpDgUMl/aukmyXNqrGeiIhoo9P3LN4VmA4cD0wCfijpdbZ/We0kaT4wH2DKlCmjXGJExNhW5xrBBmByZXxS2Va1Huix/YztnwL3UwTDi9heYrvbdndXV1dtBUdENFGdQbASmC5pmqTxwFygp6XPtynWBpA0gWJT0Zoaa4qIiBa1BYHtzcACYAVwL7DM9ipJF0iaXXZbATwm6R7gRuCjth+rq6aIiNhSrfsIbC8Hlre0Lao8NvChcoiIiA7ImcUREQ2XIIiIaLgEQUREwyUIIiIaLkEQEdFwCYKIiIZLEERENFyCICKi4RIEERENlyCIiGi4BEFERMMlCCIiGi5BEBHRcAmCiIiGSxBERDRcrUEgaZak1ZL6JC1sM/1sSf2Sbi+Hc+qsJyIitlTbjWkkjQMWAydT3Jt4paQe2/e0dL3K9oK66oiIiKHVuUYwE+izvcb2JmApMKfG5UVExDaoMwgmAusq4+vLtlanS7pT0jWSJtdYT0REtNHpncXfAabafj3wj8DX23WSNF9Sr6Te/v7+US0wImKsqzMINgDVX/iTyrbn2X7M9tPl6FeBo9rNyPYS2922u7u6umopNiKiqeoMgpXAdEnTJI0H5gI91Q6SDqyMzgburbGeiIhoo7ajhmxvlrQAWAGMAy61vUrSBUCv7R7gA5JmA5uBXwBn11VPRES0V1sQANheDixvaVtUefxx4ON11hAREUPr9M7iiIjosEGDQNLnJL2vTfv7JH2m3rIiImK0DLVG8HvAkjbtXwHeWk85EREx2oYKgpfZdmuj7ecA1VdSRESMpqGC4DeSprc2lm2/qa+kiIgYTUMdNbQI+J6kC4Fby7ZuiqN8/qTmuiIiYpQMGgS2vyfpNOCjwHll893A6bbvGoXaIiJiFAwaBJJ2B35u+6yW9i5Ju9v+j9qri4iI2g21j+AS4M1t2o8FPl9PORERMdqGCoKjbF/b2mj7W8Bx9ZUUERGjaagg2GMbnxcRETuRob7QH5U0s7WxbMtNASIixoihDh/9KLBM0mW8+PDRMykuKR0REWPAoGsEtn8CHE1xFvHZwMDRQ2dRhEFERIwBQ16G2vbPgU9JegMwjyIEjgO+OQq1RUTEKBjqPIJDKb785wEbgasA2T5hlGqLiIhRMNTO4vsorkD6VtvH2v4i8OxIZi5plqTVkvokLRyi3+mSLKl7JPOPiIiXbqggeDvwCHCjpK9IOpERXHVU0jhgMXAKMAOYJ2lGm357AR8EbhlJ4RERsX0MtbP427bnAocBN1JcaO4Vkv5G0n8dxrxnAn2219jeBCwF5rTp9xfAZ4FcsiIiogO2emKY7V/bvsL224BJwG3Ax4Yx74nAusr4+rLteeVO6Mm2rxt+yRERsT2N6Axh24/bXmL7xJe6YEm7ABcDHx5G3/mSeiX19vfnXLaIiO2pzktFbAAmV8YnlW0D9gIOB26StBb4HaCn3Q7jMny6bXd3dXXVWHJERPPUGQQrgemSpkkaT3E2cs/ARNtP2J5ge6rtqcDNwGzbvTXWFBERLWoLAtubgQXACuBeYJntVZIukDS7ruVGRMTIDHlm8UtlezmwvKVt0SB9j6+zloiIaC+Xk46IaLgEQUREwyUIIiIaLkEQEdFwCYKIiIZLEERENFyCICKi4RIEERENlyCIiGi4BEFERMMlCCIiGi5BEBHRcAmCiIiGSxBERDRcgiAiouFqDQJJsyStltQnaWGb6X8k6S5Jt0v6F0kz6qwnIiK2VFsQSBoHLAZOAWYA89p80V9h+3W2jwAuoriZfUREjKI61whmAn2219jeBCwF5lQ72P5VZfTlgGusJyIi2qjzVpUTgXWV8fXA0a2dJP0x8CFgPPB7NdYTERFtdHxnse3Ftg8BPgZ8sl0fSfMl9Urq7e/vH90CIyLGuDqDYAMwuTI+qWwbzFLgtHYTbC+x3W27u6ura/tVGBERtQbBSmC6pGmSxgNzgZ5qB0nTK6OnAg/UWE9ERLRR2z4C25slLQBWAOOAS22vknQB0Gu7B1gg6STgGeBx4Ky66omIiPbq3FmM7eXA8pa2RZXHH6xz+RERsXUd31kcERGdlSCIiGi4BEFERMMlCCIiGi5BEBHRcAmCiIiGSxBERDRcgiAiouESBBERDZcgiIhouARBRETDJQgiIhouQRAR0XAJgoiIhksQREQ0XIIgIqLhag0CSbMkrZbUJ2lhm+kfknSPpDslXS/p4DrriYiILdUWBJLGAYuBU4AZwDxJM1q63QZ02349cA1wUV31REREe3WuEcwE+myvsb0JWArMqXawfaPtfy9HbwYm1VhPRES0UWcQTATWVcbXl22DeS/wvRrriYiINmq9ef1wSfoDoBv43UGmzwfmA0yZMmUUK4uIGPvqXCPYAEyujE8q215E0knAJ4DZtp9uNyPbS2x32+7u6uqqpdiIiKaqMwhWAtMlTZM0HpgL9FQ7SDoS+DJFCDxaYy0RETGI2oLA9mZgAbACuBdYZnuVpAskzS67fQ7YE7ha0u2SegaZXURE1KTWfQS2lwPLW9oWVR6fVOfyIyJi63JmcUREwyUIIiIaLkEQEdFwCYKIiIZLEERENFyCICKi4RIEERENlyCIiGi4BEFERMMlCCIiGi5BEBHRcAmCiIiGSxBERDRcgiAiouESBBERDVdrEEiaJWm1pD5JC9tMP07Sv0naLOmMOmuJiIj2agsCSeOAxcApwAxgnqQZLd0eAs4GrqirjoiIGFqddyibCfTZXgMgaSkwB7hnoIPtteW052qsIyIihlDnpqGJwLrK+PqyLSIidiA7xc5iSfMl9Urq7e/v73Q5ERFjSp1BsAGYXBmfVLaNmO0ltrttd3d1dW2X4iIiolBnEKwEpkuaJmk8MBfoqXF5ERGxDWoLAtubgQXACuBeYJntVZIukDQbQNJ/kbQeeAfwZUmr6qonIiLaq/OoIWwvB5a3tC2qPF5JsckoIiI6ZKfYWRwREfVJEERENFyCICKi4RIEERENlyCIiGi4BEFERMMlCCIiGi5BEBHRcAmCiIiGSxBERDRcgiAiouESBBERDZcgiIhouARBRETDJQgiIhouQRAR0XC1BoGkWZJWS+qTtLDN9JdJuqqcfoukqXXWExERW6otCCSNAxYDpwAzgHmSZrR0ey/wuO1XA58HPltXPRER0V6dawQzgT7ba2xvApYCc1r6zAG+Xj6+BjhRkmqsKSIiWtR5z+KJwLrK+Hrg6MH62N4s6QngAGBjtZOk+cD8cvQpSatrqXjHMIGW118nZR1se8pnt3Mb65/fwYNNqPXm9duL7SXAkk7XMRok9dru7nQdMXL57HZuTf786tw0tAGYXBmfVLa17SNpV2Af4LEaa4qIiBZ1BsFKYLqkaZLGA3OBnpY+PcBZ5eMzgBtsu8aaIiKiRW2bhspt/guAFcA44FLbqyRdAPTa7gG+BlwuqQ/4BUVYNF0jNoGNUfnsdm6N/fyUH+AREc2WM4sjIhouQRAR0XAJgoiIhtspziMYyyQdRnGG9cSyaQPQY/vezlUVMfaV//cmArfYfqrSPsv29ztX2ejLGkEHSfoYxaU3BPykHARc2e4ifbHzkPSHna4hBifpA8D/Bc4D7pZUvfzNX3amqs7JUUMdJOl+4LW2n2lpHw+ssj29M5XFSyXpIdtTOl1HtCfpLuAY20+VVz2+Brjc9hck3Wb7yM5WOLqyaaizngMOAh5saT+wnBY7MEl3DjYJeOVo1hIjtsvA5iDbayUdD1wj6WCKz69REgSd9SfA9ZIe4IUL9E0BXg0s6FRRMWyvBH4feLylXcD/G/1yYgR+LukI27cDlGsGbwUuBV7X0co6IEHQQba/L+lQikt2V3cWr7T9bOcqi2H6LrDnwJdJlaSbRr2aGIkzgc3VBtubgTMlfbkzJXVO9hFERDRcjhqKiGi4BEFERMMlCCIiGi5BEDs9SU9tvdeI57lW0oQ27ZJ0g6S9y/FPSFol6U5Jt0tqvR3rqJL0T5L262QNsfNJEESMzFuAO2z/StIxwFuBN9h+PXASL75PdydcDry/wzXETiZBEGOSpEMkfV/SrZJ+VF5XBklvk3SLpNvKX8+vLNsPkPSD8tf9Vxn8pKL3UFyaAIoT/zbafhrA9kbbD5fzO0rSP5fLXyHpwLL9XEkrJd0h6ZuS9ijb3yHp7rL9h2Xb7pL+TtJdZb0nlO1nS7q2fH0PSLqoUl8PMG87vpXRBLYzZNipB+CpNm3XA9PLx0dT3AYVYD9eOGz6HOCvyseXAIvKx6cCBia0me+DwF7l4z2B24H7gS8Bv1u270ZxQllXOf4uijv0ARxQmdeFwHnl47uAieXjfct/P1x53mHAQ8DuwNnAGop7fO9e1jS5Mt8HqsvJkGFrQ04oizFH0p7AG4Grped/2L+s/HcScFX5C3088NOy/Tjg7QC2r5PUerbwgP1tP1n2e0rSUcCbgRPK+S4EeoHDgX8slz8OeKR8/uGSLgT2pQiSFWX7vwKXSVoGXFu2HQt8sVzWfZIeBA4tp11v+4ny9d4DHMwLm6Uepbh0yWNbfbMiyJnFMTbtAvzS9hFtpn0RuNh2T3l9mfNHOO/Nknax/RyAizPAbwJuKi9kdhZwK8VFA49p8/zLgNNs3yHpbOD4cj5/VO5oPhW4tQyYoTxdefwsL/6/vDvwmxG+rmiw7COIMcf2r4CfSnoHPH+kz2+Xk/ehuIwHFF/aA34IvLvsfwrFJqR2VgOvKvv9lqTqFWKPoNhMsxroKncmI2k3Sa8t++wFPCJpN4r9DZR9DrF9i+1FQD8wGfjRQJ/yUiRTynkPSsUqyH8G1g7VL6IqawQxFuwhaX1l/GKKL9C/kfRJim32S4E7KNYAri43/dwATCuf82mK+0Csoti+/9Agy7qO4ld8H8WmnS9K2pfiujV9wHzbmySdAVwiaR+K/2d/DawC/hy4heLL/haKYAD4XBkqoti/cQdwX/ka7irnf7btpyubu9o5CrjZxXVzIoYl1xqKGIFy38Lf2z6507W0I+kLFHe4u77TtcTOI5uGIkbA9iPAVwZOKNsB3Z0QiJHKGkFERMNljSAiouESBBERDZcgiIhouARBRETDJQgiIhru/wPqUjo6SM80UQAAAABJRU5ErkJggg==\n", "text/plain": [ "
" ] }, "metadata": { "needs_background": "light" }, "output_type": "display_data" }, { "data": { "application/javascript": [ "\n", " setTimeout(function() {\n", " var nbb_cell_id = 15;\n", " var nbb_unformatted_code = \"skillds.sst.to_series().plot(kind=\\\"bar\\\")\\nplt.title(\\\"NCEP-CFSv2 Nino34 ACC\\\")\\nplt.xlabel(\\\"Lead (Season)\\\")\\nplt.ylabel(\\\"ACC\\\")\";\n", " var nbb_formatted_code = \"skillds.sst.to_series().plot(kind=\\\"bar\\\")\\nplt.title(\\\"NCEP-CFSv2 Nino34 ACC\\\")\\nplt.xlabel(\\\"Lead (Season)\\\")\\nplt.ylabel(\\\"ACC\\\")\";\n", " var nbb_cells = Jupyter.notebook.get_cells();\n", " for (var i = 0; i < nbb_cells.length; ++i) {\n", " if (nbb_cells[i].input_prompt_number == nbb_cell_id) {\n", " if (nbb_cells[i].get_text() == nbb_unformatted_code) {\n", " nbb_cells[i].set_text(nbb_formatted_code);\n", " }\n", " break;\n", " }\n", " }\n", " }, 500);\n", " " ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "skillds.sst.to_series().plot(kind=\"bar\")\n", "plt.title(\"NCEP-CFSv2 Nino34 ACC\")\n", "plt.xlabel(\"Lead (Season)\")\n", "plt.ylabel(\"ACC\")\n", "plt.show()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "### References\n", "\n", "```{eval-rst}\n", ".. bibliography::\n", " :filter: docname in docnames\n", "```" ] } ], "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.10" }, "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 }