Demonstrates plotting hourly tidal displacements for the Arctic Ocean
OTIS format tidal solutions provided by Ohio State University and ESR
Global Tide Model (GOT) solutions provided by Richard Ray at GSFC
Finite Element Solution (FES) provided by AVISO
calc_astrol_longitudes.py
: computes the basic astronomical mean longitudescalc_delta_time.py
: calculates difference between universal and dynamic timeconvert_ll_xy.py
: convert lat/lon points to and from projected coordinatesload_constituent.py
: loads parameters for a given tidal constituentload_nodal_corrections.py
: load the nodal corrections for tidal constituentsinfer_minor_corrections.py
: return corrections for minor constituentsmodel.py
: retrieves tide model parameters for named tide modelsread_tide_model.py
: extract tidal harmonic constants from OTIS tide modelsread_netcdf_model.py
: extract tidal harmonic constants from netcdf modelsread_GOT_model.py
: extract tidal harmonic constants from GSFC GOT modelsread_FES_model.py
: extract tidal harmonic constants from FES tide modelspredict_tide.py
: predict tidal elevation at a single time using harmonic constantsThis notebook uses Jupyter widgets to set parameters for calculating the tidal maps.
import os
import pyproj
import datetime
import numpy as np
import matplotlib
matplotlib.rcParams['axes.linewidth'] = 2.0
matplotlib.rcParams["animation.html"] = "jshtml"
import matplotlib.pyplot as plt
import matplotlib.animation as animation
import cartopy.crs as ccrs
import ipywidgets as widgets
from IPython.display import HTML
import pyTMD.time
import pyTMD.model
import pyTMD.tools
from pyTMD.calc_delta_time import calc_delta_time
from pyTMD.read_tide_model import extract_tidal_constants
from pyTMD.read_netcdf_model import extract_netcdf_constants
from pyTMD.read_GOT_model import extract_GOT_constants
from pyTMD.read_FES_model import extract_FES_constants
from pyTMD.infer_minor_corrections import infer_minor_corrections
from pyTMD.predict_tide import predict_tide
# autoreload
%load_ext autoreload
%autoreload 2
# available model list
model_list = sorted(pyTMD.model.global_ocean() + pyTMD.model.arctic_ocean())
# display widgets for setting directory and model
TMDwidgets = pyTMD.tools.widgets()
TMDwidgets.model.options = model_list
TMDwidgets.model.value = 'TPXO8-atlas'
widgets.VBox([
TMDwidgets.directory,
TMDwidgets.model,
TMDwidgets.compress,
])
# get model parameters
model = pyTMD.model(TMDwidgets.directory.value,
format=TMDwidgets.atlas.value,
compressed=TMDwidgets.compress.value
).elevation(TMDwidgets.model.value)
# create an image around the Arctic Ocean
# use NSIDC Polar Stereographic definitions
# https://nsidc.org/data/polar-stereo/ps_grids.html
xlimits = [-3850000,3750000]
ylimits = [-5350000,5850000]
spacing = [5e3,-5e3]
# x and y coordinates
x = np.arange(xlimits[0],xlimits[1]+spacing[0],spacing[0])
y = np.arange(ylimits[1],ylimits[0]+spacing[1],spacing[1])
xgrid,ygrid = np.meshgrid(x,y)
# x and y dimensions
nx = int((xlimits[1]-xlimits[0])/spacing[0])+1
ny = int((ylimits[0]-ylimits[1])/spacing[1])+1
# convert image coordinates from polar stereographic to latitude/longitude
crs1 = pyproj.CRS.from_string("epsg:{0:d}".format(3413))
crs2 = pyproj.CRS.from_string("epsg:{0:d}".format(4326))
transformer = pyproj.Transformer.from_crs(crs1, crs2, always_xy=True)
lon,lat = transformer.transform(xgrid.flatten(), ygrid.flatten())
# convert from calendar date to days relative to Jan 1, 1992 (48622 MJD)
YMD = TMDwidgets.datepick.value
tide_time = pyTMD.time.convert_calendar_dates(YMD.year, YMD.month,
YMD.day, hour=np.arange(24))
# delta time (TT - UT1) file
delta_file = pyTMD.utilities.get_data_path(['data','merged_deltat.data'])
# read tidal constants and interpolate to grid points
if model.format in ('OTIS','ATLAS','ESR'):
amp,ph,D,c = extract_tidal_constants(lon, lat, model.grid_file,
model.model_file, model.projection, type=model.type,
method='spline', grid=model.format)
DELTAT = np.zeros_like(tide_time)
elif (model.format == 'netcdf'):
amp,ph,D,c = extract_netcdf_constants(lon, lat, model.grid_file,
model.model_file, type=model.type, method='spline',
scale=model.scale, compressed=model.compressed)
DELTAT = np.zeros_like(tide_time)
elif (model.format == 'GOT'):
amp,ph,c = extract_GOT_constants(lon, lat, model.model_file,
method='spline', scale=model.scale,
compressed=model.compressed)
# interpolate delta times from calendar dates to tide time
DELTAT = calc_delta_time(delta_file, tide_time)
elif (model.format == 'FES'):
amp,ph = extract_FES_constants(lon, lat, model.model_file,
type=model.type, version=model.version, method='spline',
scale=model.scale, compressed=model.compressed)
c = model.constituents
# interpolate delta times from calendar dates to tide time
DELTAT = calc_delta_time(delta_file, tide_time)
# calculate complex phase in radians for Euler's
cph = -1j*ph*np.pi/180.0
# calculate constituent oscillation
hc = amp*np.exp(cph)
# allocate for tide map calculated every hour
tide_cm = np.ma.zeros((ny,nx,24))
for hour in range(24):
# predict tidal elevations at time and infer minor corrections
TIDE = predict_tide(tide_time[hour], hc, c, deltat=DELTAT[hour],
corrections=model.format)
MINOR = infer_minor_corrections(tide_time[hour], hc, c,
deltat=DELTAT[hour], corrections=model.format)
# add major and minor components and reform grid
# convert from meters to centimeters
tide_cm[:,:,hour] = 100.0*np.reshape((TIDE+MINOR),(ny,nx))
# output Arctic Ocean Tide Animation
projection = ccrs.Stereographic(central_longitude=-45.0,
central_latitude=+90.0,true_scale_latitude=+70.0)
fig, ax = plt.subplots(num=1, figsize=(8,9),
subplot_kw=dict(projection=projection))
# plot tide height
vmin,vmax = (np.min(tide_cm), np.max(tide_cm))
extent = (xlimits[0],xlimits[1],ylimits[0],ylimits[1])
im = ax.imshow(np.zeros((ny,nx)), interpolation='nearest',
vmin=vmin, vmax=vmax, transform=projection,
extent=extent, origin='upper', animated=True)
# add 50m resolution cartopy coastlines
ax.coastlines('50m')
# Add colorbar and adjust size
# pad = distance from main plot axis
# extend = add extension triangles to upper and lower bounds
# options: neither, both, min, max
# shrink = percent size of colorbar
# aspect = lengthXwidth aspect of colorbar
cbar = plt.colorbar(im, ax=ax, pad=0.025, extend='both',
extendfrac=0.0375, shrink=0.90, aspect=25.5, drawedges=False)
# rasterized colorbar to remove lines
cbar.solids.set_rasterized(True)
# Add label to the colorbar
cbar.ax.set_ylabel('{0} Tide Height'.format(model.name), fontsize=13)
cbar.ax.set_xlabel('cm', fontsize=13)
cbar.ax.xaxis.set_label_coords(0.50, 1.04)
# ticks lines all the way across
cbar.ax.tick_params(which='both', width=1, length=19,
labelsize=13, direction='in')
# add title (date and time)
ttl = ax.set_title(None, fontsize=13)
# set x and y limits
ax.set_xlim(xlimits)
ax.set_ylim(ylimits)
# stronger linewidth on frame
ax.spines['geo'].set_linewidth(2.0)
ax.spines['geo'].set_capstyle('projecting')
# adjust subplot within figure
fig.subplots_adjust(left=0.02,right=0.98,bottom=0.05,top=0.95)
# animate each map
def animate_maps(hour):
# set map data
im.set_data(tide_cm[:,:,hour])
# set title
args = (YMD.year,YMD.month,YMD.day,hour)
ttl.set_text('{0:4d}-{1:02d}-{2:02d}T{3:02d}:00:00'.format(*args))
# set animation
anim = animation.FuncAnimation(fig, animate_maps, frames=24)
%matplotlib inline
HTML(anim.to_jshtml())