{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "## Tutorial on how to use `timestaps` in Field construction" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "from parcels import Field\n", "from glob import glob\n", "import numpy as np" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Some NetCDF files, such as for example those from the [World Ocean Atlas](https://www.nodc.noaa.gov/OC5/woa18/), have time calendars that can't be parsed by `xarray`. These result in a `ValueError: unable to decode time units`, for example when the calendar is in 'months since' a particular date.\n", "\n", "In these cases, a workaround in Parcels is to use the `timestamps` argument in `Field` (or `FieldSet`) creation. Here, we show how this works for example temperature data from the World Ocean Atlas in the Pacific Ocean" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The following cell raises an error, since the calendar of the World Ocean Atlas data is in \"months since 1955-01-01 00:00:00\"" ] }, { "cell_type": "code", "execution_count": 2, "metadata": { "tags": [ "raises-exception" ] }, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "WARNING: File WOA_data/woa18_decav_t10_04.nc could not be decoded properly by xarray (version 0.20.1).\n", " It will be opened with no decoding. Filling values might be wrongly parsed.\n" ] }, { "ename": "RuntimeError", "evalue": "Xarray could not convert the calendar. If youre using from_netcdf, try using the timestamps keyword in the construction of your Field. See also the tutorial at https://nbviewer.jupyter.org/github/OceanParcels/parcels/blob/master/parcels/examples/tutorial_timestamps.ipynb", "output_type": "error", "traceback": [ "\u001b[0;31m---------------------------------------------------------------------------\u001b[0m", "\u001b[0;31mKeyError\u001b[0m Traceback (most recent call last)", "File \u001b[0;32m/opt/anaconda3/envs/py3_parcels/lib/python3.8/site-packages/xarray/coding/times.py:240\u001b[0m, in \u001b[0;36mdecode_cf_datetime\u001b[0;34m(num_dates, units, calendar, use_cftime)\u001b[0m\n\u001b[1;32m 239\u001b[0m \u001b[38;5;28;01mtry\u001b[39;00m:\n\u001b[0;32m--> 240\u001b[0m dates \u001b[38;5;241m=\u001b[39m \u001b[43m_decode_datetime_with_pandas\u001b[49m\u001b[43m(\u001b[49m\u001b[43mflat_num_dates\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43munits\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43mcalendar\u001b[49m\u001b[43m)\u001b[49m\n\u001b[1;32m 241\u001b[0m \u001b[38;5;28;01mexcept\u001b[39;00m (\u001b[38;5;167;01mKeyError\u001b[39;00m, OutOfBoundsDatetime, \u001b[38;5;167;01mOverflowError\u001b[39;00m):\n", "File \u001b[0;32m/opt/anaconda3/envs/py3_parcels/lib/python3.8/site-packages/xarray/coding/times.py:186\u001b[0m, in \u001b[0;36m_decode_datetime_with_pandas\u001b[0;34m(flat_num_dates, units, calendar)\u001b[0m\n\u001b[1;32m 185\u001b[0m delta, ref_date \u001b[38;5;241m=\u001b[39m _unpack_netcdf_time_units(units)\n\u001b[0;32m--> 186\u001b[0m delta \u001b[38;5;241m=\u001b[39m \u001b[43m_netcdf_to_numpy_timeunit\u001b[49m\u001b[43m(\u001b[49m\u001b[43mdelta\u001b[49m\u001b[43m)\u001b[49m\n\u001b[1;32m 187\u001b[0m \u001b[38;5;28;01mtry\u001b[39;00m:\n", "File \u001b[0;32m/opt/anaconda3/envs/py3_parcels/lib/python3.8/site-packages/xarray/coding/times.py:83\u001b[0m, in \u001b[0;36m_netcdf_to_numpy_timeunit\u001b[0;34m(units)\u001b[0m\n\u001b[1;32m 82\u001b[0m units \u001b[38;5;241m=\u001b[39m \u001b[38;5;124mf\u001b[39m\u001b[38;5;124m\"\u001b[39m\u001b[38;5;132;01m{\u001b[39;00munits\u001b[38;5;132;01m}\u001b[39;00m\u001b[38;5;124ms\u001b[39m\u001b[38;5;124m\"\u001b[39m\n\u001b[0;32m---> 83\u001b[0m \u001b[38;5;28;01mreturn\u001b[39;00m \u001b[43m{\u001b[49m\n\u001b[1;32m 84\u001b[0m \u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mnanoseconds\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m:\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mns\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m,\u001b[49m\n\u001b[1;32m 85\u001b[0m \u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mmicroseconds\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m:\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mus\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m,\u001b[49m\n\u001b[1;32m 86\u001b[0m \u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mmilliseconds\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m:\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mms\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m,\u001b[49m\n\u001b[1;32m 87\u001b[0m \u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mseconds\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m:\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43ms\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m,\u001b[49m\n\u001b[1;32m 88\u001b[0m \u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mminutes\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m:\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mm\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m,\u001b[49m\n\u001b[1;32m 89\u001b[0m \u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mhours\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m:\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mh\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m,\u001b[49m\n\u001b[1;32m 90\u001b[0m \u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mdays\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m:\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[38;5;124;43mD\u001b[39;49m\u001b[38;5;124;43m\"\u001b[39;49m\u001b[43m,\u001b[49m\n\u001b[1;32m 91\u001b[0m \u001b[43m\u001b[49m\u001b[43m}\u001b[49m\u001b[43m[\u001b[49m\u001b[43munits\u001b[49m\u001b[43m]\u001b[49m\n", "\u001b[0;31mKeyError\u001b[0m: 'months'", "\nDuring handling of the above exception, another exception occurred:\n", "\u001b[0;31mValueError\u001b[0m Traceback (most recent call last)", "File \u001b[0;32m/opt/anaconda3/envs/py3_parcels/lib/python3.8/site-packages/xarray/coding/times.py:153\u001b[0m, in \u001b[0;36m_decode_cf_datetime_dtype\u001b[0;34m(data, units, calendar, use_cftime)\u001b[0m\n\u001b[1;32m 152\u001b[0m \u001b[38;5;28;01mtry\u001b[39;00m:\n\u001b[0;32m--> 153\u001b[0m result \u001b[38;5;241m=\u001b[39m \u001b[43mdecode_cf_datetime\u001b[49m\u001b[43m(\u001b[49m\u001b[43mexample_value\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43munits\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43mcalendar\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43muse_cftime\u001b[49m\u001b[43m)\u001b[49m\n\u001b[1;32m 154\u001b[0m \u001b[38;5;28;01mexcept\u001b[39;00m \u001b[38;5;167;01mException\u001b[39;00m:\n", "File \u001b[0;32m/opt/anaconda3/envs/py3_parcels/lib/python3.8/site-packages/xarray/coding/times.py:242\u001b[0m, in \u001b[0;36mdecode_cf_datetime\u001b[0;34m(num_dates, units, calendar, use_cftime)\u001b[0m\n\u001b[1;32m 241\u001b[0m \u001b[38;5;28;01mexcept\u001b[39;00m (\u001b[38;5;167;01mKeyError\u001b[39;00m, OutOfBoundsDatetime, \u001b[38;5;167;01mOverflowError\u001b[39;00m):\n\u001b[0;32m--> 242\u001b[0m dates \u001b[38;5;241m=\u001b[39m \u001b[43m_decode_datetime_with_cftime\u001b[49m\u001b[43m(\u001b[49m\n\u001b[1;32m 243\u001b[0m \u001b[43m \u001b[49m\u001b[43mflat_num_dates\u001b[49m\u001b[38;5;241;43m.\u001b[39;49m\u001b[43mastype\u001b[49m\u001b[43m(\u001b[49m\u001b[38;5;28;43mfloat\u001b[39;49m\u001b[43m)\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43munits\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43mcalendar\u001b[49m\n\u001b[1;32m 244\u001b[0m \u001b[43m \u001b[49m\u001b[43m)\u001b[49m\n\u001b[1;32m 246\u001b[0m \u001b[38;5;28;01mif\u001b[39;00m (\n\u001b[1;32m 247\u001b[0m dates[np\u001b[38;5;241m.\u001b[39mnanargmin(num_dates)]\u001b[38;5;241m.\u001b[39myear \u001b[38;5;241m<\u001b[39m \u001b[38;5;241m1678\u001b[39m\n\u001b[1;32m 248\u001b[0m \u001b[38;5;129;01mor\u001b[39;00m dates[np\u001b[38;5;241m.\u001b[39mnanargmax(num_dates)]\u001b[38;5;241m.\u001b[39myear \u001b[38;5;241m>\u001b[39m\u001b[38;5;241m=\u001b[39m \u001b[38;5;241m2262\u001b[39m\n\u001b[1;32m 249\u001b[0m ):\n", "File \u001b[0;32m/opt/anaconda3/envs/py3_parcels/lib/python3.8/site-packages/xarray/coding/times.py:174\u001b[0m, in \u001b[0;36m_decode_datetime_with_cftime\u001b[0;34m(num_dates, units, calendar)\u001b[0m\n\u001b[1;32m 172\u001b[0m \u001b[38;5;28;01mraise\u001b[39;00m \u001b[38;5;167;01mModuleNotFoundError\u001b[39;00m(\u001b[38;5;124m\"\u001b[39m\u001b[38;5;124mNo module named \u001b[39m\u001b[38;5;124m'\u001b[39m\u001b[38;5;124mcftime\u001b[39m\u001b[38;5;124m'\u001b[39m\u001b[38;5;124m\"\u001b[39m)\n\u001b[1;32m 173\u001b[0m \u001b[38;5;28;01mreturn\u001b[39;00m np\u001b[38;5;241m.\u001b[39masarray(\n\u001b[0;32m--> 174\u001b[0m \u001b[43mcftime\u001b[49m\u001b[38;5;241;43m.\u001b[39;49m\u001b[43mnum2date\u001b[49m\u001b[43m(\u001b[49m\u001b[43mnum_dates\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43munits\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43mcalendar\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43monly_use_cftime_datetimes\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[38;5;28;43;01mTrue\u001b[39;49;00m\u001b[43m)\u001b[49m\n\u001b[1;32m 175\u001b[0m )\n", "File \u001b[0;32msrc/cftime/_cftime.pyx:498\u001b[0m, in \u001b[0;36mcftime._cftime.num2date\u001b[0;34m()\u001b[0m\n", "File \u001b[0;32msrc/cftime/_cftime.pyx:99\u001b[0m, in \u001b[0;36mcftime._cftime._dateparse\u001b[0;34m()\u001b[0m\n", "\u001b[0;31mValueError\u001b[0m: 'months since' units only allowed for '360_day' calendar", "\nDuring handling of the above exception, another exception occurred:\n", "\u001b[0;31mValueError\u001b[0m Traceback (most recent call last)", "File \u001b[0;32m~/Codes/parcels/parcels/tools/converters.py:233\u001b[0m, in \u001b[0;36mconvert_xarray_time_units\u001b[0;34m(ds, time)\u001b[0m\n\u001b[1;32m 232\u001b[0m \u001b[38;5;28;01mtry\u001b[39;00m:\n\u001b[0;32m--> 233\u001b[0m da2 \u001b[38;5;241m=\u001b[39m \u001b[43mxr\u001b[49m\u001b[38;5;241;43m.\u001b[39;49m\u001b[43mdecode_cf\u001b[49m\u001b[43m(\u001b[49m\u001b[43mda2\u001b[49m\u001b[43m)\u001b[49m\n\u001b[1;32m 234\u001b[0m \u001b[38;5;28;01mexcept\u001b[39;00m \u001b[38;5;167;01mValueError\u001b[39;00m:\n", "File \u001b[0;32m/opt/anaconda3/envs/py3_parcels/lib/python3.8/site-packages/xarray/conventions.py:646\u001b[0m, in \u001b[0;36mdecode_cf\u001b[0;34m(obj, concat_characters, mask_and_scale, decode_times, decode_coords, drop_variables, use_cftime, decode_timedelta)\u001b[0m\n\u001b[1;32m 644\u001b[0m \u001b[38;5;28;01mraise\u001b[39;00m \u001b[38;5;167;01mTypeError\u001b[39;00m(\u001b[38;5;124m\"\u001b[39m\u001b[38;5;124mcan only decode Dataset or DataStore objects\u001b[39m\u001b[38;5;124m\"\u001b[39m)\n\u001b[0;32m--> 646\u001b[0m \u001b[38;5;28mvars\u001b[39m, attrs, coord_names \u001b[38;5;241m=\u001b[39m \u001b[43mdecode_cf_variables\u001b[49m\u001b[43m(\u001b[49m\n\u001b[1;32m 647\u001b[0m \u001b[43m \u001b[49m\u001b[38;5;28;43mvars\u001b[39;49m\u001b[43m,\u001b[49m\n\u001b[1;32m 648\u001b[0m \u001b[43m \u001b[49m\u001b[43mattrs\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 649\u001b[0m \u001b[43m \u001b[49m\u001b[43mconcat_characters\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 650\u001b[0m \u001b[43m \u001b[49m\u001b[43mmask_and_scale\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 651\u001b[0m \u001b[43m \u001b[49m\u001b[43mdecode_times\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 652\u001b[0m \u001b[43m \u001b[49m\u001b[43mdecode_coords\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 653\u001b[0m \u001b[43m \u001b[49m\u001b[43mdrop_variables\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[43mdrop_variables\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 654\u001b[0m \u001b[43m \u001b[49m\u001b[43muse_cftime\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[43muse_cftime\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 655\u001b[0m \u001b[43m \u001b[49m\u001b[43mdecode_timedelta\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[43mdecode_timedelta\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 656\u001b[0m \u001b[43m\u001b[49m\u001b[43m)\u001b[49m\n\u001b[1;32m 657\u001b[0m ds \u001b[38;5;241m=\u001b[39m Dataset(\u001b[38;5;28mvars\u001b[39m, attrs\u001b[38;5;241m=\u001b[39mattrs)\n", "File \u001b[0;32m/opt/anaconda3/envs/py3_parcels/lib/python3.8/site-packages/xarray/conventions.py:512\u001b[0m, in \u001b[0;36mdecode_cf_variables\u001b[0;34m(variables, attributes, concat_characters, mask_and_scale, decode_times, decode_coords, drop_variables, use_cftime, decode_timedelta)\u001b[0m\n\u001b[1;32m 506\u001b[0m stack_char_dim \u001b[38;5;241m=\u001b[39m (\n\u001b[1;32m 507\u001b[0m concat_characters\n\u001b[1;32m 508\u001b[0m \u001b[38;5;129;01mand\u001b[39;00m v\u001b[38;5;241m.\u001b[39mdtype \u001b[38;5;241m==\u001b[39m \u001b[38;5;124m\"\u001b[39m\u001b[38;5;124mS1\u001b[39m\u001b[38;5;124m\"\u001b[39m\n\u001b[1;32m 509\u001b[0m \u001b[38;5;129;01mand\u001b[39;00m v\u001b[38;5;241m.\u001b[39mndim \u001b[38;5;241m>\u001b[39m \u001b[38;5;241m0\u001b[39m\n\u001b[1;32m 510\u001b[0m \u001b[38;5;129;01mand\u001b[39;00m stackable(v\u001b[38;5;241m.\u001b[39mdims[\u001b[38;5;241m-\u001b[39m\u001b[38;5;241m1\u001b[39m])\n\u001b[1;32m 511\u001b[0m )\n\u001b[0;32m--> 512\u001b[0m new_vars[k] \u001b[38;5;241m=\u001b[39m \u001b[43mdecode_cf_variable\u001b[49m\u001b[43m(\u001b[49m\n\u001b[1;32m 513\u001b[0m \u001b[43m \u001b[49m\u001b[43mk\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 514\u001b[0m \u001b[43m \u001b[49m\u001b[43mv\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 515\u001b[0m \u001b[43m \u001b[49m\u001b[43mconcat_characters\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[43mconcat_characters\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 516\u001b[0m \u001b[43m \u001b[49m\u001b[43mmask_and_scale\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[43mmask_and_scale\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 517\u001b[0m \u001b[43m \u001b[49m\u001b[43mdecode_times\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[43mdecode_times\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 518\u001b[0m \u001b[43m \u001b[49m\u001b[43mstack_char_dim\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[43mstack_char_dim\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 519\u001b[0m \u001b[43m \u001b[49m\u001b[43muse_cftime\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[43muse_cftime\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 520\u001b[0m \u001b[43m \u001b[49m\u001b[43mdecode_timedelta\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[43mdecode_timedelta\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 521\u001b[0m \u001b[43m\u001b[49m\u001b[43m)\u001b[49m\n\u001b[1;32m 522\u001b[0m \u001b[38;5;28;01mif\u001b[39;00m decode_coords \u001b[38;5;129;01min\u001b[39;00m [\u001b[38;5;28;01mTrue\u001b[39;00m, \u001b[38;5;124m\"\u001b[39m\u001b[38;5;124mcoordinates\u001b[39m\u001b[38;5;124m\"\u001b[39m, \u001b[38;5;124m\"\u001b[39m\u001b[38;5;124mall\u001b[39m\u001b[38;5;124m\"\u001b[39m]:\n", "File \u001b[0;32m/opt/anaconda3/envs/py3_parcels/lib/python3.8/site-packages/xarray/conventions.py:360\u001b[0m, in \u001b[0;36mdecode_cf_variable\u001b[0;34m(name, var, concat_characters, mask_and_scale, decode_times, decode_endianness, stack_char_dim, use_cftime, decode_timedelta)\u001b[0m\n\u001b[1;32m 359\u001b[0m \u001b[38;5;28;01mif\u001b[39;00m decode_times:\n\u001b[0;32m--> 360\u001b[0m var \u001b[38;5;241m=\u001b[39m \u001b[43mtimes\u001b[49m\u001b[38;5;241;43m.\u001b[39;49m\u001b[43mCFDatetimeCoder\u001b[49m\u001b[43m(\u001b[49m\u001b[43muse_cftime\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[43muse_cftime\u001b[49m\u001b[43m)\u001b[49m\u001b[38;5;241;43m.\u001b[39;49m\u001b[43mdecode\u001b[49m\u001b[43m(\u001b[49m\u001b[43mvar\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43mname\u001b[49m\u001b[38;5;241;43m=\u001b[39;49m\u001b[43mname\u001b[49m\u001b[43m)\u001b[49m\n\u001b[1;32m 362\u001b[0m dimensions, data, attributes, encoding \u001b[38;5;241m=\u001b[39m variables\u001b[38;5;241m.\u001b[39munpack_for_decoding(var)\n", "File \u001b[0;32m/opt/anaconda3/envs/py3_parcels/lib/python3.8/site-packages/xarray/coding/times.py:526\u001b[0m, in \u001b[0;36mCFDatetimeCoder.decode\u001b[0;34m(self, variable, name)\u001b[0m\n\u001b[1;32m 525\u001b[0m calendar \u001b[38;5;241m=\u001b[39m pop_to(attrs, encoding, \u001b[38;5;124m\"\u001b[39m\u001b[38;5;124mcalendar\u001b[39m\u001b[38;5;124m\"\u001b[39m)\n\u001b[0;32m--> 526\u001b[0m dtype \u001b[38;5;241m=\u001b[39m \u001b[43m_decode_cf_datetime_dtype\u001b[49m\u001b[43m(\u001b[49m\u001b[43mdata\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43munits\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43mcalendar\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;28;43mself\u001b[39;49m\u001b[38;5;241;43m.\u001b[39;49m\u001b[43muse_cftime\u001b[49m\u001b[43m)\u001b[49m\n\u001b[1;32m 527\u001b[0m transform \u001b[38;5;241m=\u001b[39m partial(\n\u001b[1;32m 528\u001b[0m decode_cf_datetime,\n\u001b[1;32m 529\u001b[0m units\u001b[38;5;241m=\u001b[39munits,\n\u001b[1;32m 530\u001b[0m calendar\u001b[38;5;241m=\u001b[39mcalendar,\n\u001b[1;32m 531\u001b[0m use_cftime\u001b[38;5;241m=\u001b[39m\u001b[38;5;28mself\u001b[39m\u001b[38;5;241m.\u001b[39muse_cftime,\n\u001b[1;32m 532\u001b[0m )\n", "File \u001b[0;32m/opt/anaconda3/envs/py3_parcels/lib/python3.8/site-packages/xarray/coding/times.py:163\u001b[0m, in \u001b[0;36m_decode_cf_datetime_dtype\u001b[0;34m(data, units, calendar, use_cftime)\u001b[0m\n\u001b[1;32m 158\u001b[0m msg \u001b[38;5;241m=\u001b[39m (\n\u001b[1;32m 159\u001b[0m \u001b[38;5;124mf\u001b[39m\u001b[38;5;124m\"\u001b[39m\u001b[38;5;124munable to decode time units \u001b[39m\u001b[38;5;132;01m{\u001b[39;00munits\u001b[38;5;132;01m!r}\u001b[39;00m\u001b[38;5;124m with \u001b[39m\u001b[38;5;132;01m{\u001b[39;00mcalendar_msg\u001b[38;5;132;01m!r}\u001b[39;00m\u001b[38;5;124m. Try \u001b[39m\u001b[38;5;124m\"\u001b[39m\n\u001b[1;32m 160\u001b[0m \u001b[38;5;124m\"\u001b[39m\u001b[38;5;124mopening your dataset with decode_times=False or installing cftime \u001b[39m\u001b[38;5;124m\"\u001b[39m\n\u001b[1;32m 161\u001b[0m \u001b[38;5;124m\"\u001b[39m\u001b[38;5;124mif it is not installed.\u001b[39m\u001b[38;5;124m\"\u001b[39m\n\u001b[1;32m 162\u001b[0m )\n\u001b[0;32m--> 163\u001b[0m \u001b[38;5;28;01mraise\u001b[39;00m \u001b[38;5;167;01mValueError\u001b[39;00m(msg)\n\u001b[1;32m 164\u001b[0m \u001b[38;5;28;01melse\u001b[39;00m:\n", "\u001b[0;31mValueError\u001b[0m: unable to decode time units 'months since 1955-01-01 00:00:00' with 'the default calendar'. Try opening your dataset with decode_times=False or installing cftime if it is not installed.", "\nDuring handling of the above exception, another exception occurred:\n", "\u001b[0;31mRuntimeError\u001b[0m Traceback (most recent call last)", "Input \u001b[0;32mIn [2]\u001b[0m, in \u001b[0;36m\u001b[0;34m()\u001b[0m\n\u001b[0;32m----> 1\u001b[0m tempfield \u001b[38;5;241m=\u001b[39m \u001b[43mField\u001b[49m\u001b[38;5;241;43m.\u001b[39;49m\u001b[43mfrom_netcdf\u001b[49m\u001b[43m(\u001b[49m\u001b[43mglob\u001b[49m\u001b[43m(\u001b[49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[38;5;124;43mWOA_data/woa18_decav_*_04.nc\u001b[39;49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[43m)\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[38;5;124;43mt_an\u001b[39;49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\n\u001b[1;32m 2\u001b[0m \u001b[43m \u001b[49m\u001b[43m{\u001b[49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[38;5;124;43mlon\u001b[39;49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[43m:\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[38;5;124;43mlon\u001b[39;49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[38;5;124;43mlat\u001b[39;49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[43m:\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[38;5;124;43mlat\u001b[39;49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[38;5;124;43mtime\u001b[39;49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[43m:\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[38;5;124;43mtime\u001b[39;49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[43m}\u001b[49m\u001b[43m)\u001b[49m\n", "File \u001b[0;32m~/Codes/parcels/parcels/field.py:390\u001b[0m, in \u001b[0;36mField.from_netcdf\u001b[0;34m(cls, filenames, variable, dimensions, indices, grid, mesh, timestamps, allow_time_extrapolation, time_periodic, deferred_load, **kwargs)\u001b[0m\n\u001b[1;32m 385\u001b[0m \u001b[38;5;28;01mraise\u001b[39;00m \u001b[38;5;167;01mRuntimeError\u001b[39;00m(\u001b[38;5;124m'\u001b[39m\u001b[38;5;124mMultiple files given but no time dimension specified\u001b[39m\u001b[38;5;124m'\u001b[39m)\n\u001b[1;32m 387\u001b[0m \u001b[38;5;28;01mif\u001b[39;00m grid \u001b[38;5;129;01mis\u001b[39;00m \u001b[38;5;28;01mNone\u001b[39;00m:\n\u001b[1;32m 388\u001b[0m \u001b[38;5;66;03m# Concatenate time variable to determine overall dimension\u001b[39;00m\n\u001b[1;32m 389\u001b[0m \u001b[38;5;66;03m# across multiple files\u001b[39;00m\n\u001b[0;32m--> 390\u001b[0m time, time_origin, timeslices, dataFiles \u001b[38;5;241m=\u001b[39m \u001b[38;5;28;43mcls\u001b[39;49m\u001b[38;5;241;43m.\u001b[39;49m\u001b[43mcollect_timeslices\u001b[49m\u001b[43m(\u001b[49m\u001b[43mtimestamps\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43mdata_filenames\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 391\u001b[0m \u001b[43m \u001b[49m\u001b[43m_grid_fb_class\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43mdimensions\u001b[49m\u001b[43m,\u001b[49m\n\u001b[1;32m 392\u001b[0m \u001b[43m \u001b[49m\u001b[43mindices\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43mnetcdf_engine\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[43mnetcdf_decodewarning\u001b[49m\u001b[43m)\u001b[49m\n\u001b[1;32m 393\u001b[0m grid \u001b[38;5;241m=\u001b[39m Grid\u001b[38;5;241m.\u001b[39mcreate_grid(lon, lat, depth, time, time_origin\u001b[38;5;241m=\u001b[39mtime_origin, mesh\u001b[38;5;241m=\u001b[39mmesh)\n\u001b[1;32m 394\u001b[0m grid\u001b[38;5;241m.\u001b[39mtimeslices \u001b[38;5;241m=\u001b[39m timeslices\n", "File \u001b[0;32m~/Codes/parcels/parcels/field.py:248\u001b[0m, in \u001b[0;36mField.collect_timeslices\u001b[0;34m(timestamps, data_filenames, _grid_fb_class, dimensions, indices, netcdf_engine, netcdf_decodewarning)\u001b[0m\n\u001b[1;32m 245\u001b[0m \u001b[38;5;28;01mfor\u001b[39;00m fname \u001b[38;5;129;01min\u001b[39;00m data_filenames:\n\u001b[1;32m 246\u001b[0m \u001b[38;5;28;01mwith\u001b[39;00m _grid_fb_class(fname, dimensions, indices, netcdf_engine\u001b[38;5;241m=\u001b[39mnetcdf_engine,\n\u001b[1;32m 247\u001b[0m netcdf_decodewarning\u001b[38;5;241m=\u001b[39mnetcdf_decodewarning) \u001b[38;5;28;01mas\u001b[39;00m filebuffer:\n\u001b[0;32m--> 248\u001b[0m ftime \u001b[38;5;241m=\u001b[39m \u001b[43mfilebuffer\u001b[49m\u001b[38;5;241;43m.\u001b[39;49m\u001b[43mtime\u001b[49m\n\u001b[1;32m 249\u001b[0m timeslices\u001b[38;5;241m.\u001b[39mappend(ftime)\n\u001b[1;32m 250\u001b[0m dataFiles\u001b[38;5;241m.\u001b[39mappend([fname] \u001b[38;5;241m*\u001b[39m \u001b[38;5;28mlen\u001b[39m(ftime))\n", "File \u001b[0;32m~/Codes/parcels/parcels/fieldfilebuffer.py:171\u001b[0m, in \u001b[0;36mNetcdfFileBuffer.time\u001b[0;34m(self)\u001b[0m\n\u001b[1;32m 169\u001b[0m \u001b[38;5;129m@property\u001b[39m\n\u001b[1;32m 170\u001b[0m \u001b[38;5;28;01mdef\u001b[39;00m \u001b[38;5;21mtime\u001b[39m(\u001b[38;5;28mself\u001b[39m):\n\u001b[0;32m--> 171\u001b[0m \u001b[38;5;28;01mreturn\u001b[39;00m \u001b[38;5;28;43mself\u001b[39;49m\u001b[38;5;241;43m.\u001b[39;49m\u001b[43mtime_access\u001b[49m\u001b[43m(\u001b[49m\u001b[43m)\u001b[49m\n", "File \u001b[0;32m~/Codes/parcels/parcels/fieldfilebuffer.py:181\u001b[0m, in \u001b[0;36mNetcdfFileBuffer.time_access\u001b[0;34m(self)\u001b[0m\n\u001b[1;32m 178\u001b[0m \u001b[38;5;28;01mreturn\u001b[39;00m np\u001b[38;5;241m.\u001b[39marray([\u001b[38;5;28;01mNone\u001b[39;00m])\n\u001b[1;32m 180\u001b[0m time_da \u001b[38;5;241m=\u001b[39m \u001b[38;5;28mself\u001b[39m\u001b[38;5;241m.\u001b[39mdataset[\u001b[38;5;28mself\u001b[39m\u001b[38;5;241m.\u001b[39mdimensions[\u001b[38;5;124m'\u001b[39m\u001b[38;5;124mtime\u001b[39m\u001b[38;5;124m'\u001b[39m]]\n\u001b[0;32m--> 181\u001b[0m \u001b[43mconvert_xarray_time_units\u001b[49m\u001b[43m(\u001b[49m\u001b[43mtime_da\u001b[49m\u001b[43m,\u001b[49m\u001b[43m \u001b[49m\u001b[38;5;28;43mself\u001b[39;49m\u001b[38;5;241;43m.\u001b[39;49m\u001b[43mdimensions\u001b[49m\u001b[43m[\u001b[49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[38;5;124;43mtime\u001b[39;49m\u001b[38;5;124;43m'\u001b[39;49m\u001b[43m]\u001b[49m\u001b[43m)\u001b[49m\n\u001b[1;32m 182\u001b[0m time \u001b[38;5;241m=\u001b[39m np\u001b[38;5;241m.\u001b[39marray([time_da[\u001b[38;5;28mself\u001b[39m\u001b[38;5;241m.\u001b[39mdimensions[\u001b[38;5;124m'\u001b[39m\u001b[38;5;124mtime\u001b[39m\u001b[38;5;124m'\u001b[39m]]\u001b[38;5;241m.\u001b[39mdata]) \u001b[38;5;28;01mif\u001b[39;00m \u001b[38;5;28mlen\u001b[39m(time_da\u001b[38;5;241m.\u001b[39mshape) \u001b[38;5;241m==\u001b[39m \u001b[38;5;241m0\u001b[39m \u001b[38;5;28;01melse\u001b[39;00m np\u001b[38;5;241m.\u001b[39marray(time_da[\u001b[38;5;28mself\u001b[39m\u001b[38;5;241m.\u001b[39mdimensions[\u001b[38;5;124m'\u001b[39m\u001b[38;5;124mtime\u001b[39m\u001b[38;5;124m'\u001b[39m]])\n\u001b[1;32m 183\u001b[0m \u001b[38;5;28;01mif\u001b[39;00m \u001b[38;5;28misinstance\u001b[39m(time[\u001b[38;5;241m0\u001b[39m], datetime\u001b[38;5;241m.\u001b[39mdatetime):\n", "File \u001b[0;32m~/Codes/parcels/parcels/tools/converters.py:235\u001b[0m, in \u001b[0;36mconvert_xarray_time_units\u001b[0;34m(ds, time)\u001b[0m\n\u001b[1;32m 233\u001b[0m da2 \u001b[38;5;241m=\u001b[39m xr\u001b[38;5;241m.\u001b[39mdecode_cf(da2)\n\u001b[1;32m 234\u001b[0m \u001b[38;5;28;01mexcept\u001b[39;00m \u001b[38;5;167;01mValueError\u001b[39;00m:\n\u001b[0;32m--> 235\u001b[0m \u001b[38;5;28;01mraise\u001b[39;00m \u001b[38;5;167;01mRuntimeError\u001b[39;00m(\u001b[38;5;124m'\u001b[39m\u001b[38;5;124mXarray could not convert the calendar. If you\u001b[39m\u001b[38;5;124m'\u001b[39m\u001b[38;5;124m'\u001b[39m\u001b[38;5;124mre using from_netcdf, \u001b[39m\u001b[38;5;124m'\u001b[39m\n\u001b[1;32m 236\u001b[0m \u001b[38;5;124m'\u001b[39m\u001b[38;5;124mtry using the timestamps keyword in the construction of your Field. \u001b[39m\u001b[38;5;124m'\u001b[39m\n\u001b[1;32m 237\u001b[0m \u001b[38;5;124m'\u001b[39m\u001b[38;5;124mSee also the tutorial at https://nbviewer.jupyter.org/github/OceanParcels/\u001b[39m\u001b[38;5;124m'\u001b[39m\n\u001b[1;32m 238\u001b[0m \u001b[38;5;124m'\u001b[39m\u001b[38;5;124mparcels/blob/master/parcels/examples/tutorial_timestamps.ipynb\u001b[39m\u001b[38;5;124m'\u001b[39m)\n\u001b[1;32m 239\u001b[0m ds[time] \u001b[38;5;241m=\u001b[39m da2[time]\n", "\u001b[0;31mRuntimeError\u001b[0m: Xarray could not convert the calendar. If youre using from_netcdf, try using the timestamps keyword in the construction of your Field. See also the tutorial at https://nbviewer.jupyter.org/github/OceanParcels/parcels/blob/master/parcels/examples/tutorial_timestamps.ipynb" ] } ], "source": [ "tempfield = Field.from_netcdf(glob('WOA_data/woa18_decav_*_04.nc'), 't_an', \n", " {'lon': 'lon', 'lat': 'lat', 'time': 'time'})" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "However, we can create our own numpy array of timestamps associated with each of the 12 snapshots in the netcdf file" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "timestamps = np.expand_dims(np.array([np.datetime64('2001-%.2d-15' %m) for m in range(1,13)]), axis=1)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "And then we can add the `timestamps` as an extra argument" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [], "source": [ "tempfield = Field.from_netcdf(glob('WOA_data/woa18_decav_*_04.nc'), 't_an',\n", " {'lon': 'lon', 'lat': 'lat', 'time': 'time'},\n", " netcdf_decodewarning=False,\n", " timestamps=timestamps)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note, by the way, that adding the `time_periodic=True` argument to `Field.from_netcdf()` will also mean that the climatology can be cycled for multiple years." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Furthermore, note that we used `netcdf_decodewarning=False` in the `FieldSet.from_netcdf()` call above. This is to silence an expected warning because the time dimension in the `coordinates.nc` file can't be decoded by `xarray`." ] } ], "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.8.13" } }, "nbformat": 4, "nbformat_minor": 2 }