{ "cells": [ { "attachments": {}, "cell_type": "markdown", "metadata": {}, "source": [ "Plot Tide Forecasts\n", "===================\n", "\n", "Plots the daily tidal displacements for a given location\n", "\n", "OTIS format tidal solutions provided by Ohio State University and ESR \n", "- http://volkov.oce.orst.edu/tides/region.html \n", "- https://www.esr.org/research/polar-tide-models/list-of-polar-tide-models/\n", "- ftp://ftp.esr.org/pub/datasets/tmd/ \n", "\n", "Global Tide Model (GOT) solutions provided by Richard Ray at GSFC \n", "\n", "Finite Element Solution (FES) provided by AVISO \n", "- https://www.aviso.altimetry.fr/en/data/products/auxiliary-products/global-tide-fes.html\n", "\n", "#### Python Dependencies\n", " - [numpy: Scientific Computing Tools For Python](https://www.numpy.org) \n", " - [scipy: Scientific Tools for Python](https://www.scipy.org/) \n", " - [pyproj: Python interface to PROJ library](https://pypi.org/project/pyproj/) \n", " - [netCDF4: Python interface to the netCDF C library](https://unidata.github.io/netcdf4-python/) \n", " - [matplotlib: Python 2D plotting library](https://matplotlib.org/) \n", " - [ipyleaflet: Jupyter / Leaflet bridge enabling interactive maps](https://github.com/jupyter-widgets/ipyleaflet) \n", "\n", "#### Program Dependencies\n", "\n", "- `calc_astrol_longitudes.py`: computes the basic astronomical mean longitudes \n", "- `convert_ll_xy.py`: convert lat/lon points to and from projected coordinates \n", "- `load_constituent.py`: loads parameters for a given tidal constituent \n", "- `load_nodal_corrections.py`: load the nodal corrections for tidal constituents \n", "- `io.model.py`: retrieves tide model parameters for named tide models\n", "- `io.OTIS.py`: extract tidal harmonic constants from OTIS tide models \n", "- `io.ATLAS.py`: extract tidal harmonic constants from ATLAS netCDF4 tide models \n", "- `io.GOT.py`: extract tidal harmonic constants from GOT tide models \n", "- `io.FES.py`: extract tidal harmonic constants from FES tide models \n", "- `io.constituents.py`: basic tide model constituent class \n", "- `predict.py`: predict tidal values using harmonic constants \n", "- `time.py`: utilities for calculating time operations\n", "- `utilities.py`: download and management utilities for files\n", "\n", "This notebook uses Jupyter widgets to set parameters for calculating the tidal maps. " ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Load modules" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "from __future__ import print_function\n", "\n", "import os\n", "import datetime\n", "import numpy as np\n", "import matplotlib.pyplot as plt\n", "import ipywidgets as widgets\n", "import IPython.display\n", "\n", "# import tide programs\n", "import pyTMD.io\n", "import pyTMD.time\n", "import pyTMD.predict\n", "import pyTMD.tools\n", "import pyTMD.utilities\n", "\n", "# autoreload\n", "%load_ext autoreload\n", "%autoreload 2" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "# available model list\n", "model_list = sorted(pyTMD.io.model.ocean_elevation())\n", "# display widgets for setting directory and model\n", "TMDwidgets = pyTMD.tools.widgets()\n", "TMDwidgets.model.options = model_list\n", "TMDwidgets.model.value = 'GOT4.10'\n", "widgets.VBox([\n", " TMDwidgets.directory,\n", " TMDwidgets.model,\n", " TMDwidgets.atlas,\n", " TMDwidgets.compress,\n", " TMDwidgets.datepick\n", "])" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "# default coordinates to use\n", "LAT,LON = (32.86710263,-117.25750387)\n", "m = pyTMD.tools.leaflet(center=(LAT,LON), zoom=12,\n", " zoom_control=True, marker_control=True)\n", "# show map\n", "m.map" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "%matplotlib widget\n", "\n", "# get model parameters\n", "model = pyTMD.io.model(TMDwidgets.directory.value,\n", " format=TMDwidgets.atlas.value,\n", " compressed=TMDwidgets.compress.value\n", " ).elevation(TMDwidgets.model.value)\n", "\n", "# convert from calendar date to days relative to Jan 1, 1992 (48622 MJD)\n", "YMD = TMDwidgets.datepick.value\n", "# calculate a weeks forecast every minute\n", "minutes = np.arange(7*1440)\n", "# convert time from MJD to days relative to Jan 1, 1992 (48622 MJD)\n", "tide_time = pyTMD.time.convert_calendar_dates(YMD.year, YMD.month,\n", " YMD.day, minute=minutes)\n", "hours = minutes/60.0\n", "\n", "# create plot with tidal displacements, high and low tides and dates\n", "fig,ax1 = plt.subplots(num=1)\n", "xmax = np.ceil(hours[-1]).astype('i')\n", "l1, = ax1.plot([], [], 'k')\n", "l2, = ax1.plot([], [], 'r*')\n", "l3, = ax1.plot([], [], 'b*')\n", "for h in range(24,xmax,24):\n", " ax1.axvline(h,color='gray',lw=0.5,ls='dashed',dashes=(11,5))\n", "ax1.set_xlim(0,xmax)\n", "ax1.set_ylabel(f'{model.name} Tidal Displacement [cm]')\n", "args = (YMD.year,YMD.month,YMD.day)\n", "ax1.set_xlabel('Time from {0:4d}-{1:02d}-{2:02d} UTC [Hours]'.format(*args))\n", "ttl = ax1.set_title(None)\n", "fig.subplots_adjust(left=0.10,right=0.98,bottom=0.10,top=0.95)\n", "\n", "# delta time (TT - UT1) file\n", "delta_file = pyTMD.utilities.get_data_path(['data','merged_deltat.data'])\n", "\n", "# read tidal constants and interpolate to leaflet points\n", "if model.format in ('OTIS','ATLAS','ESR'):\n", " constituents = pyTMD.io.OTIS.read_constants(\n", " model.grid_file, model.model_file,\n", " model.projection, type=model.type,\n", " grid=model.format)\n", " c = constituents.fields\n", " DELTAT = np.zeros_like(tide_time)\n", "elif (model.format == 'netcdf'):\n", " constituents = pyTMD.io.ATLAS.read_constants(\n", " model.grid_file, model.model_file,\n", " type=model.type, compressed=model.compressed)\n", " c = constituents.fields\n", " DELTAT = np.zeros_like(tide_time)\n", "elif (model.format == 'GOT'):\n", " constituents = pyTMD.io.GOT.read_constants(\n", " model.model_file, compressed=model.compressed)\n", " c = constituents.fields\n", " # interpolate delta times from calendar dates to tide time\n", " DELTAT = pyTMD.time.interpolate_delta_time(delta_file, tide_time)\n", "elif (model.format == 'FES'):\n", " constituents = pyTMD.io.FES.read_constants(model.model_file,\n", " type=model.type, version=model.version,\n", " compressed=model.compressed)\n", " c = model.constituents\n", " # interpolate delta times from calendar dates to tide time\n", " DELTAT = pyTMD.time.interpolate_delta_time(delta_file, tide_time)\n", "\n", "# update the tide prediction and plot\n", "def update_tide_prediction(*args):\n", " # leaflet location\n", " LAT,LON = np.copy(m.marker.location)\n", " # verify longitudes\n", " LON = m.wrap_longitudes(LON)\n", " if model.format in ('OTIS','ATLAS','ESR'):\n", " amp,ph,D = pyTMD.io.OTIS.interpolate_constants(\n", " np.atleast_1d(LON), np.atleast_1d(LAT),\n", " constituents, model.projection, type=model.type,\n", " method='spline', extrapolate=True)\n", " elif (model.format == 'netcdf'):\n", " amp,ph,D = pyTMD.io.ATLAS.interpolate_constants(\n", " np.atleast_1d(LON), np.atleast_1d(LAT),\n", " constituents, type=model.type, scale=model.scale,\n", " method='spline', extrapolate=True)\n", " elif (model.format == 'GOT'):\n", " amp,ph = pyTMD.io.GOT.interpolate_constants(\n", " np.atleast_1d(LON), np.atleast_1d(LAT),\n", " constituents, scale=model.scale,\n", " method='spline', extrapolate=True)\n", " elif (model.format == 'FES'):\n", " amp,ph = pyTMD.io.FES.interpolate_constants(\n", " np.atleast_1d(LON), np.atleast_1d(LAT),\n", " constituents, scale=model.scale,\n", " method='spline', extrapolate=True)\n", " # calculate complex phase in radians for Euler's\n", " cph = -1j*ph*np.pi/180.0\n", " # calculate constituent oscillation\n", " hc = amp*np.exp(cph)\n", " # predict tidal elevations at time 1 and infer minor corrections\n", " TIDE = pyTMD.predict.time_series(tide_time, hc, c,\n", " deltat=DELTAT, corrections=model.format)\n", " MINOR = pyTMD.predict.infer_minor(tide_time, hc, c,\n", " deltat=DELTAT, corrections=model.format)\n", " TIDE.data[:] += MINOR.data[:]\n", " # convert to centimeters\n", " TIDE.data[:] *= 100.0\n", "\n", " # differentiate to calculate high and low tides\n", " diff = np.zeros_like(tide_time, dtype=np.float64)\n", " # forward differentiation for starting point\n", " diff[0] = TIDE.data[1] - TIDE.data[0]\n", " # backward differentiation for end point\n", " diff[-1] = TIDE.data[-1] - TIDE.data[-2]\n", " # centered differentiation for all others\n", " diff[1:-1] = (TIDE.data[2:] - TIDE.data[0:-2])/2.0\n", " # indices of high and low tides\n", " htindex, = np.nonzero((np.sign(diff[0:-1]) >= 0) & (np.sign(diff[1:]) < 0))\n", " ltindex, = np.nonzero((np.sign(diff[0:-1]) <= 0) & (np.sign(diff[1:]) > 0))\n", " # update plot data\n", " l1.set_data(hours, TIDE.data)\n", " l2.set_data(hours[htindex], TIDE.data[htindex])\n", " l3.set_data(hours[ltindex], TIDE.data[ltindex])\n", " # update plot title\n", " ttl.set_text(u'{0:0.6f}\\u00b0N {1:0.6f}\\u00b0W'.format(LAT,LON))\n", " ax1.relim()\n", " ax1.autoscale_view()\n", " fig.canvas.draw()\n", "\n", "# run tide prediction at initial location\n", "update_tide_prediction()\n", "# watch marker location for changes\n", "m.marker_text.observe(update_tide_prediction)" ] } ], "metadata": { "interpreter": { "hash": "31f2aee4e71d21fbe5cf8b01ff0e069b9275f58929596ceb00d14d90e3e16cd6" }, "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.10.6" } }, "nbformat": 4, "nbformat_minor": 4 }