compute
Calculates tidal elevations at points and times
Can use OTIS format tidal solutions provided by Ohio State University and ESR
Can use Global Tide Model (GOT) solutions provided by Richard Ray at GSFC
Can use Finite Element Solution (FES) models provided by AVISO
Calculates tidal currents at points and times
Can use OTIS format tidal solutions provided by Ohio State University and ESR
Can use Finite Element Solution (FES) models provided by AVISO
Calculates long-period equilibrium tides (LPET) at points and times
Uses the summation of fifteen tidal spectral lines from Cartwright and Edden, (1973)
Calculates radial pole load tides (LPT) at points and times
Following IERS Convention (2010) guidelines
Calculates radial ocean pole load tides (OPT) at points and times
Following IERS Convention (2010) guidelines
Calculates radial solid earth tides (SET) at points and times
Following IERS Convention (2010) guidelines
Calling Sequence
import pyTMD
tide_h = pyTMD.compute.tide_elevations(x, y, delta_time,
DIRECTORY=DIRECTORY, MODEL=MODEL, EPOCH=(2000,1,1,0,0,0),
EPSG=3031, TYPE='drift')
tide_uv = pyTMD.compute.tide_currents(x, y, delta_time,
DIRECTORY=DIRECTORY, MODEL=MODEL, EPOCH=(2000,1,1,0,0,0),
EPSG=3031, TYPE='drift')
- pyTMD.compute.corrections(x: ndarray, y: ndarray, delta_time: ndarray, CORRECTION: str = 'ocean', **kwargs)[source]
Wrapper function to compute tide corrections at points and times
- Parameters
- x: np.ndarray
x-coordinates in projection EPSG
- y: np.ndarray
y-coordinates in projection EPSG
- delta_time: np.ndarray
seconds since EPOCH or datetime array
- CORRECTION: str, default ‘ocean’
Correction type to compute
'ocean'
: ocean tide from model constituents'load'
: load tide from model constituents'LPET'
: long-period equilibrium tide'LPT'
: solid earth load pole tide'OPT'
: ocean pole tide'SET'
: solid earth tide
- **kwargs: dict
keyword arguments for correction functions
- Returns
- values: np.ndarray
tidal correction at coordinates and time in meters
- pyTMD.compute.tide_elevations(x: ndarray, y: ndarray, delta_time: ndarray, DIRECTORY: str | pathlib.Path | None = None, MODEL: str | None = None, ATLAS_FORMAT: str = 'netcdf', GZIP: bool = False, DEFINITION_FILE: str | pathlib.Path | None = None, EPSG: str | int = 3031, EPOCH: list | tuple = (2000, 1, 1, 0, 0, 0), TYPE: str | None = 'drift', TIME: str = 'UTC', METHOD: str = 'spline', EXTRAPOLATE: bool = False, CUTOFF: int | float = 10.0, INFER_MINOR: bool = True, APPLY_FLEXURE: bool = False, FILL_VALUE: float = nan, **kwargs)[source]
Compute ocean or load tides at points and times from model constituents
- Parameters
- x: np.ndarray
x-coordinates in projection EPSG
- y: np.ndarray
y-coordinates in projection EPSG
- delta_time: np.ndarray
seconds since EPOCH or datetime array
- DIRECTORY: str or NoneType, default None
working data directory for tide models
- MODEL: str or NoneType, default None
Tide model to use in correction
- ATLAS_FORMAT: str, default ‘netcdf’
ATLAS tide model format
'OTIS'
'netcdf'
- GZIP: bool, default False
Tide model files are gzip compressed
- DEFINITION_FILE: str or NoneType, default None
Tide model definition file for use
- EPSG: int, default: 3031 (Polar Stereographic South, WGS84)
Input coordinate system
- EPOCH: tuple, default (2000,1,1,0,0,0)
Time period for calculating delta times
- TYPE: str or NoneType, default ‘drift’
Input data type
None
: determined from input variable dimensions'drift'
: drift buoys or satellite/airborne altimetry'grid'
: spatial grids or images'time series'
: time series at a single point
- TIME: str, default ‘UTC’
Time type if need to compute leap seconds to convert to UTC
'GPS'
: leap seconds needed'LORAN'
: leap seconds needed (LORAN = GPS + 9 seconds)'TAI'
: leap seconds needed (TAI = GPS + 19 seconds)'UTC'
: no leap seconds needed'datetime'
: numpy datatime array in UTC
- METHOD: str
Interpolation method
`bilinear`
: quick bilinear interpolation`spline`
: scipy bivariate spline interpolation`linear`
,`nearest`
: scipy regular grid interpolations
- EXTRAPOLATE: bool, default False
Extrapolate with nearest-neighbors
- CUTOFF: int or float, default 10.0
Extrapolation cutoff in kilometers
Set to
np.inf
to extrapolate for all points- INFER_MINOR: bool, default True
Infer the height values for minor tidal constituents
- APPLY_FLEXURE: bool, default False
Apply ice flexure scaling factor to height values
Only valid for models containing flexure fields
- FILL_VALUE: float, default np.nan
Output invalid value
- Returns
- tide: np.ndarray
tidal elevation at coordinates and time in meters
- pyTMD.compute.tide_currents(x: ndarray, y: ndarray, delta_time: ndarray, DIRECTORY: str | pathlib.Path | None = None, MODEL: str | None = None, ATLAS_FORMAT: str = 'netcdf', GZIP: bool = False, DEFINITION_FILE: str | pathlib.Path | None = None, EPSG: str | int = 3031, EPOCH: list | tuple = (2000, 1, 1, 0, 0, 0), TYPE: str | None = 'drift', TIME: str = 'UTC', METHOD: str = 'spline', EXTRAPOLATE: bool = False, CUTOFF: int | float = 10.0, INFER_MINOR: bool = True, FILL_VALUE: float = nan, **kwargs)[source]
Compute ocean tide currents at points and times from model constituents
- Parameters
- x: np.ndarray
x-coordinates in projection EPSG
- y: np.ndarray
y-coordinates in projection EPSG
- delta_time: np.ndarray
seconds since EPOCH or datetime array
- DIRECTORY: str or NoneType, default None
working data directory for tide models
- MODEL: str or NoneType, default None
Tide model to use in correction
- ATLAS_FORMAT: str, default ‘netcdf’
ATLAS tide model format
'OTIS'
'netcdf'
- GZIP: bool, default False
Tide model files are gzip compressed
- DEFINITION_FILE: str or NoneType, default None
Tide model definition file for use
- EPSG: int, default: 3031 (Polar Stereographic South, WGS84)
Input coordinate system
- EPOCH: tuple, default (2000,1,1,0,0,0)
Time period for calculating delta times
- TYPE: str or NoneType, default ‘drift’
Input data type
None
: determined from input variable dimensions'drift'
: drift buoys or satellite/airborne altimetry'grid'
: spatial grids or images'time series'
: time series at a single point
- TIME: str, default ‘UTC’
Time type if need to compute leap seconds to convert to UTC
'GPS'
: leap seconds needed'LORAN'
: leap seconds needed (LORAN = GPS + 9 seconds)'TAI'
: leap seconds needed (TAI = GPS + 19 seconds)'UTC'
: no leap seconds needed'datetime'
: numpy datatime array in UTC
- METHOD: str
Interpolation method
`bilinear`
: quick bilinear interpolation`spline`
: scipy bivariate spline interpolation`linear`
,`nearest`
: scipy regular grid interpolations
- EXTRAPOLATE: bool, default False
Extrapolate with nearest-neighbors
- CUTOFF: int or float, default 10.0
Extrapolation cutoff in kilometers
Set to
np.inf
to extrapolate for all points- INFER_MINOR: bool, default True
Infer the height values for minor tidal constituents
- FILL_VALUE: float, default np.nan
Output invalid value
- Returns
- tide: dict
tidal currents at coordinates and times
- pyTMD.compute.LPET_elevations(x: ndarray, y: ndarray, delta_time: ndarray, EPSG: str | int = 3031, EPOCH: list | tuple = (2000, 1, 1, 0, 0, 0), TYPE: str | None = 'drift', TIME: str = 'UTC', **kwargs)[source]
Compute long-period equilibrium tidal elevations at points and times
- Parameters
- x: np.ndarray
x-coordinates in projection EPSG
- y: np.ndarray
y-coordinates in projection EPSG
- delta_time: np.ndarray
seconds since EPOCH or datetime array
- EPSG: int, default: 3031 (Polar Stereographic South, WGS84)
Input coordinate system
- EPOCH: tuple, default (2000,1,1,0,0,0)
Time period for calculating delta times
- TYPE: str or NoneType, default ‘drift’
Input data type
None
: determined from input variable dimensions'drift'
: drift buoys or satellite/airborne altimetry'grid'
: spatial grids or images'time series'
: time series at a single point
- TIME: str, default ‘UTC’
Time type if need to compute leap seconds to convert to UTC
'GPS'
: leap seconds needed'LORAN'
: leap seconds needed (LORAN = GPS + 9 seconds)'TAI'
: leap seconds needed (TAI = GPS + 19 seconds)'UTC'
: no leap seconds needed'datetime'
: numpy datatime array in UTC
- FILL_VALUE: float, default np.nan
Output invalid value
- Returns
- tide_lpe: np.ndarray
long-period equilibrium tide at coordinates and time in meters
- pyTMD.compute.LPT_displacements(x: ndarray, y: ndarray, delta_time: ndarray, h: float | numpy.ndarray = 0.0, EPSG: str | int = 3031, EPOCH: list | tuple = (2000, 1, 1, 0, 0, 0), TYPE: str | None = 'drift', TIME: str = 'UTC', ELLIPSOID: str = 'WGS84', CONVENTION: str = '2018', FILL_VALUE: float = nan, **kwargs)[source]
Compute radial load pole tide displacements at points and times following IERS Convention (2010) guidelines
- Parameters
- x: np.ndarray
x-coordinates in projection EPSG
- y: np.ndarray
y-coordinates in projection EPSG
- delta_time: np.ndarray
seconds since EPOCH or datetime array
- h: float or np.ndarray, default 0.0
height coordinates in meters
- EPSG: int, default: 3031 (Polar Stereographic South, WGS84)
Input coordinate system
- EPOCH: tuple, default (2000,1,1,0,0,0)
Time period for calculating delta times
- TYPE: str or NoneType, default ‘drift’
Input data type
None
: determined from input variable dimensions'drift'
: drift buoys or satellite/airborne altimetry'grid'
: spatial grids or images'time series'
: time series at a single point
- TIME: str, default ‘UTC’
Time type if need to compute leap seconds to convert to UTC
'GPS'
: leap seconds needed'LORAN'
: leap seconds needed (LORAN = GPS + 9 seconds)'TAI'
: leap seconds needed (TAI = GPS + 19 seconds)'UTC'
: no leap seconds needed'datetime'
: numpy datatime array in UTC
- ELLIPSOID: str, default ‘WGS84’
Ellipsoid for calculating Earth parameters
- CONVENTION: str, default ‘2018’
IERS Mean or Secular Pole Convention
'2003'
'2010'
'2015'
'2018'
- FILL_VALUE: float, default np.nan
Output invalid value
- Returns
- Srad: np.ndarray
solid earth pole tide at coordinates and time in meters
- pyTMD.compute.OPT_displacements(x: ndarray, y: ndarray, delta_time: ndarray, h: float | numpy.ndarray = 0.0, EPSG: str | int = 3031, EPOCH: list | tuple = (2000, 1, 1, 0, 0, 0), TYPE: str | None = 'drift', TIME: str = 'UTC', ELLIPSOID: str = 'WGS84', CONVENTION: str = '2018', METHOD: str = 'spline', FILL_VALUE: float = nan, **kwargs)[source]
Compute radial ocean pole tide displacements at points and times following IERS Convention (2010) guidelines
- Parameters
- x: np.ndarray
x-coordinates in projection EPSG
- y: np.ndarray
y-coordinates in projection EPSG
- delta_time: np.ndarray
seconds since EPOCH or datetime array
- h: float or np.ndarray, default 0.0
height coordinates in meters
- EPSG: int, default: 3031 (Polar Stereographic South, WGS84)
Input coordinate system
- EPOCH: tuple, default (2000,1,1,0,0,0)
Time period for calculating delta times
- TYPE: str or NoneType, default ‘drift’
Input data type
None
: determined from input variable dimensions'drift'
: drift buoys or satellite/airborne altimetry'grid'
: spatial grids or images'time series'
: time series at a single point
- TIME: str, default ‘UTC’
Time type if need to compute leap seconds to convert to UTC
'GPS'
: leap seconds needed'LORAN'
: leap seconds needed (LORAN = GPS + 9 seconds)'TAI'
: leap seconds needed (TAI = GPS + 19 seconds)'UTC'
: no leap seconds needed'datetime'
: numpy datatime array in UTC
- ELLIPSOID: str, default ‘WGS84’
Ellipsoid for calculating Earth parameters
- CONVENTION: str, default ‘2018’
IERS Mean or Secular Pole Convention
'2003'
'2010'
'2015'
'2018'
- METHOD: str
Interpolation method
`bilinear`
: quick bilinear interpolation`spline`
: scipy bivariate spline interpolation`linear`
,`nearest`
: scipy regular grid interpolations
- FILL_VALUE: float, default np.nan
Output invalid value
- Returns
- Urad: np.ndarray
ocean pole tide at coordinates and time in meters
- pyTMD.compute.SET_displacements(x: ndarray, y: ndarray, delta_time: ndarray, h: float | numpy.ndarray = 0.0, EPSG: str | int = 3031, EPOCH: list | tuple = (2000, 1, 1, 0, 0, 0), TYPE: str | None = 'drift', TIME: str = 'UTC', ELLIPSOID: str = 'WGS84', TIDE_SYSTEM='tide_free', EPHEMERIDES='approximate', **kwargs)[source]
Compute solid earth tidal elevations at points and times following IERS Convention (2010) guidelines
- Parameters
- x: np.ndarray
x-coordinates in projection EPSG
- y: np.ndarray
y-coordinates in projection EPSG
- delta_time: np.ndarray
seconds since EPOCH or datetime array
- h: float or np.ndarray, default 0.0
height coordinates in meters
- EPSG: int, default: 3031 (Polar Stereographic South, WGS84)
Input coordinate system
- EPOCH: tuple, default (2000,1,1,0,0,0)
Time period for calculating delta times
- TYPE: str or NoneType, default ‘drift’
Input data type
None
: determined from input variable dimensions'drift'
: drift buoys or satellite/airborne altimetry'grid'
: spatial grids or images'time series'
: time series at a single point
- TIME: str, default ‘UTC’
Time type if need to compute leap seconds to convert to UTC
'GPS'
: leap seconds needed'LORAN'
: leap seconds needed (LORAN = GPS + 9 seconds)'TAI'
: leap seconds needed (TAI = GPS + 19 seconds)'UTC'
: no leap seconds needed'datetime'
: numpy datatime array in UTC
- ELLIPSOID: str, default ‘WGS84’
Ellipsoid for calculating Earth parameters
- TIDE_SYSTEM: str, default ‘tide_free’
Permanent tide system for the output solid Earth tide
'tide_free'
: no permanent direct and indirect tidal potentials'mean_tide'
: permanent tidal potentials (direct and indirect)
- EPHEMERIDES: str, default ‘approximate’
Ephemerides for calculating Earth parameters
'approximate'
: approximate lunisolar parameters'JPL'
: computed from JPL ephmerides kernel
- Returns
- tide_se: np.ndarray
solid earth tide at coordinates and time in meters