predict
Predict tidal values using harmonic constants
At a single time (
map
) such as for imageryMultiple times and locations (
drift
) such as for airborne and satellite altimetryTime series at a location (
time_series
) such as to compare with tide gauges
Predicts tidal values from minor constituents inferred using major constituents
Predicts long-period equilibrium ocean tides
Predicts solid earth tides
Calling Sequence
import pyTMD.predict
ht = pyTMD.predict.time_series(time, hc, con)
- pyTMD.predict.map(t: float | ndarray, hc: ndarray, constituents: list | ndarray, deltat: float | ndarray = 0.0, corrections: str = 'OTIS', **kwargs)[source]
Predict tides at a single time using harmonic constants [19]
- Parameters:
- t: float or np.ndarray
days relative to 1992-01-01T00:00:00
- hc: np.ndarray
harmonic constant vector
- constituents: list or np.ndarray
tidal constituent IDs
- deltat: float or np.ndarray, default 0.0
time correction for converting to Ephemeris Time (days)
- corrections: str, default ‘OTIS’
use nodal corrections from OTIS/ATLAS or GOT/FES models
- **kwargs: dict
keyword arguments for nodal corrections functions
- Returns:
- ht: np.ndarray
tide values reconstructed using the nodal corrections
- pyTMD.predict.drift(t: float | ndarray, hc: ndarray, constituents: list | ndarray, deltat: float | ndarray = 0.0, corrections: str = 'OTIS', **kwargs)[source]
Predict tides at multiple times and locations using harmonic constants [19]
- Parameters:
- t: float or np.ndarray
days relative to 1992-01-01T00:00:00
- hc: np.ndarray
harmonic constant vector
- constituents: list or np.ndarray
tidal constituent IDs
- deltat: float or np.ndarray, default 0.0
time correction for converting to Ephemeris Time (days)
- corrections: str, default ‘OTIS’
use nodal corrections from OTIS/ATLAS or GOT/FES models
- **kwargs: dict
keyword arguments for nodal corrections functions
- Returns:
- ht: np.ndarray
tidal time series reconstructed using the nodal corrections
- pyTMD.predict.time_series(t: float | ndarray, hc: ndarray, constituents: list | ndarray, deltat: float | ndarray = 0.0, corrections: str = 'OTIS', **kwargs)[source]
Predict tidal time series at a single location using harmonic constants [19]
- Parameters:
- t: float or np.ndarray
days relative to 1992-01-01T00:00:00
- hc: np.ndarray
harmonic constant vector
- constituents: list or np.ndarray
tidal constituent IDs
- deltat: float or np.ndarray, default 0.0
time correction for converting to Ephemeris Time (days)
- corrections: str, default ‘OTIS’
use nodal corrections from OTIS/ATLAS or GOT/FES models
- **kwargs: dict
keyword arguments for nodal corrections functions
- Returns:
- ht: np.ndarray
tidal time series reconstructed using the nodal corrections
- pyTMD.predict.infer_minor(t: float | ndarray, zmajor: ndarray, constituents: list | ndarray, **kwargs)[source]
Infer the tidal values for minor constituents using their relation with major constituents [17, 19, 20, 50]
- Parameters:
- t: float or np.ndarray
days relative to 1992-01-01T00:00:00
- zmajor: np.ndarray
Complex HC for given constituents/points
- constituents: list
tidal constituent IDs
- deltat: float or np.ndarray, default 0.0
time correction for converting to Ephemeris Time (days)
- corrections: str, default ‘OTIS’
use nodal corrections from OTIS/ATLAS or GOT/FES models
- raise_exception: bool, default False
Raise a
ValueError
if major constituents are not found- minor: list or None, default None
tidal constituent IDs of minor constituents for inference
- raise_exception: bool, default False
Raise a
ValueError
if major constituents are not found
- Returns:
- dh: np.ndarray
tidal time series for minor constituents
- pyTMD.predict._infer_short_period(t: float | ndarray, zmajor: ndarray, constituents: list | ndarray, **kwargs)[source]
Infer the tidal values for short-period minor constituents using their relation with major constituents [19, 45]
- Parameters:
- t: float or np.ndarray
days relative to 1992-01-01T00:00:00
- zmajor: np.ndarray
Complex HC for given constituents/points
- constituents: list
tidal constituent IDs
- deltat: float or np.ndarray, default 0.0
time correction for converting to Ephemeris Time (days)
- corrections: str, default ‘OTIS’
use nodal corrections from OTIS/ATLAS or GOT/FES models
- minor: list or None, default None
tidal constituent IDs of minor constituents for inference
- raise_exception: bool, default False
Raise a
ValueError
if major constituents are not found
- Returns:
- dh: np.ndarray
tidal time series for minor constituents
- pyTMD.predict._infer_semi_diurnal(t: float | ndarray, zmajor: ndarray, constituents: list | ndarray, **kwargs)[source]
Infer the tidal values for semi-diurnal minor constituents using their relation with major constituents [10, 35, 45]
- Parameters:
- t: float or np.ndarray
days relative to 1992-01-01T00:00:00
- zmajor: np.ndarray
Complex HC for given constituents/points
- constituents: list
tidal constituent IDs
- deltat: float or np.ndarray, default 0.0
time correction for converting to Ephemeris Time (days)
- minor: list or None, default None
tidal constituent IDs of minor constituents for inference
- method: str, default ‘linear’
method for interpolating between major constituents
‘linear’: linear interpolation
‘admittance’: Munk-Cartwright admittance interpolation
- raise_exception: bool, default False
Raise a
ValueError
if major constituents are not found
- Returns:
- dh: np.ndarray
tidal time series for minor constituents
- pyTMD.predict._infer_diurnal(t: float | ndarray, zmajor: ndarray, constituents: list | ndarray, **kwargs)[source]
Infer the tidal values for diurnal minor constituents using their relation with major constituents taking into account resonance due to free core nutation [9, 35, 46, 59]
- Parameters:
- t: float or np.ndarray
days relative to 1992-01-01T00:00:00
- zmajor: np.ndarray
Complex HC for given constituents/points
- constituents: list
tidal constituent IDs
- deltat: float or np.ndarray, default 0.0
time correction for converting to Ephemeris Time (days)
- minor: list or None, default None
tidal constituent IDs of minor constituents for inference
- method: str, default ‘linear’
method for interpolating between major constituents
‘linear’: linear interpolation
‘admittance’: Munk-Cartwright admittance interpolation
- raise_exception: bool, default False
Raise a
ValueError
if major constituents are not found
- Returns:
- dh: np.ndarray
tidal time series for minor constituents
- pyTMD.predict._infer_long_period(t: float | ndarray, zmajor: ndarray, constituents: list | ndarray, **kwargs)[source]
Infer the tidal values for long-period minor constituents using their relation with major constituents [9, 45, 47]
- Parameters:
- t: float or np.ndarray
days relative to 1992-01-01T00:00:00
- zmajor: np.ndarray
Complex HC for given constituents/points
- constituents: list
tidal constituent IDs
- deltat: float or np.ndarray, default 0.0
time correction for converting to Ephemeris Time (days)
- minor: list or None, default None
tidal constituent IDs of minor constituents for inference
- raise_exception: bool, default False
Raise a
ValueError
if major constituents are not found
- Returns:
- dh: np.ndarray
tidal time series for minor constituents
- pyTMD.predict.equilibrium_tide(t: ndarray, lat: ndarray, **kwargs)[source]
Compute the long-period equilibrium tides the summation of fifteen tidal spectral lines from Cartwright-Tayler-Edden tables [9, 10]
- Parameters:
- t: np.ndarray
time (days relative to January 1, 1992)
- lat: np.ndarray
latitude (degrees north)
- deltat: float or np.ndarray, default 0.0
time correction for converting to Ephemeris Time (days)
- corrections: str, default ‘OTIS’
use nodal corrections from OTIS/ATLAS or GOT/FES models
- constituents: list
long-period tidal constituent IDs
- Returns:
- lpet: np.ndarray
long-period equilibrium tide in meters
- pyTMD.predict.load_pole_tide(t: ndarray, XYZ: ndarray, deltat: float = 0.0, gamma_0: float = 9.80665, omega: float = 7.2921151467e-05, h2: float = 0.6207, l2: float = 0.0836, convention: str = '2018')[source]
Estimate load pole tide displacements in Cartesian coordinates [42]
- Parameters:
- t: np.ndarray
Time (days relative to January 1, 1992)
- XYZ: np.ndarray
Cartesian coordinates of the prediction points (meters)
- deltat: float or np.ndarray, default 0.0
time correction for converting to Ephemeris Time (days)
- gamma_0: float, default 9.80665
Normal gravity (m/s^2)
- omega: float, default 7.2921151467e-5
Earth’s rotation rate (radians/second)
- h2: float, default 0.6207
Degree-2 Love number of vertical displacement
- l2: float, default 0.0836
Degree-2 Love (Shida) number of horizontal displacement
- convention: str, default ‘2018’
IERS Mean or Secular Pole Convention
'2003'
'2010'
'2015'
'2018'
- Returns:
- dxt: np.ndarray
Load pole tide displacements in meters in Cartesian coordinates
- pyTMD.predict.ocean_pole_tide(t: ndarray, XYZ: ndarray, UXYZ: ndarray, deltat: float = 0.0, gamma_0: float = 9.780325, a_axis: float = 6378136.3, GM: float = 398600441800000.0, omega: float = 7.2921151467e-05, rho_w: float = 1025.0, g2: complex = 0.687 + 0.0036j, convention: str = '2018')[source]
Estimate ocean pole tide displacements in Cartesian coordinates [13, 14, 42]
- Parameters:
- t: np.ndarray
Time (days relative to January 1, 1992)
- XYZ: np.ndarray
Cartesian coordinates of the prediction points (meters)
- UXYZ: np.ndarray
Ocean pole tide values from Desai (2002)
- deltat: float or np.ndarray, default 0.0
time correction for converting to Ephemeris Time (days)
- a_axis: float, default 6378136.3
Semi-major axis of the Earth (meters)
- gamma_0: float, default 9.780325
Normal gravity (m/s^2)
- GM: float, default 3.986004418e14
Earth’s gravitational constant [m^3/s^2]
- omega: float, default 7.2921151467e-5
Earth’s rotation rate (radians/second)
- rho_w: float, default 1025.0
Density of sea water [kg/m^3]
- g2: complex, default 0.6870 + 0.0036j
Degree-2 Love number differential (1 + k2 - h2)
- convention: str, default ‘2018’
IERS Mean or Secular Pole Convention
'2003'
'2010'
'2015'
'2018'
- Returns:
- dxt: np.ndarray
Load pole tide displacements in meters in Cartesian coordinates
- pyTMD.predict.solid_earth_tide(t: ndarray, XYZ: ndarray, SXYZ: ndarray, LXYZ: ndarray, a_axis: float = 6378136.6, tide_system: str = 'tide_free', **kwargs)[source]
Compute the solid Earth tides due to the gravitational attraction of the moon and sun [30, 32, 49, 57]
- Parameters:
- t: np.ndarray
Time (days relative to January 1, 1992)
- XYZ: np.ndarray
Cartesian coordinates of the prediction points (meters)
- SXYZ: np.ndarray
Earth-centered Earth-fixed coordinates of the sun (meters)
- LXYZ: np.ndarray
Earth-centered Earth-fixed coordinates of the moon (meters)
- a_axis: float, default 6378136.3
Semi-major axis of the Earth (meters)
- 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)
- Returns:
- dxt: np.ndarray
Solid Earth tide in meters in Cartesian coordinates
- pyTMD.predict._out_of_phase_diurnal(XYZ: ndarray, SXYZ: ndarray, LXYZ: ndarray, F2_solar: ndarray, F2_lunar: ndarray)[source]
Computes the out-of-phase corrections induced by mantle anelasticity in the diurnal band
- Parameters:
- XYZ: np.ndarray
Cartesian coordinates of the prediction points (meters)
- SXYZ: np.ndarray
Earth-centered Earth-fixed coordinates of the sun (meters)
- LXYZ: np.ndarray
Earth-centered Earth-fixed coordinates of the moon (meters)
- F2_solar: np.ndarray
Factors for the sun
- F2_lunar: np.ndarray
Factors for the moon
- pyTMD.predict._out_of_phase_semidiurnal(XYZ: ndarray, SXYZ: ndarray, LXYZ: ndarray, F2_solar: ndarray, F2_lunar: ndarray)[source]
Computes the out-of-phase corrections induced by mantle anelasticity in the semi-diurnal band
- Parameters:
- XYZ: np.ndarray
Cartesian coordinates of the prediction points (meters)
- SXYZ: np.ndarray
Earth-centered Earth-fixed coordinates of the sun (meters)
- LXYZ: np.ndarray
Earth-centered Earth-fixed coordinates of the moon (meters)
- F2_solar: np.ndarray
Factors for the sun
- F2_lunar: np.ndarray
Factors for the moon
- pyTMD.predict._latitude_dependence(XYZ: ndarray, SXYZ: ndarray, LXYZ: ndarray, F2_solar: ndarray, F2_lunar: ndarray)[source]
Computes the corrections induced by the latitude of the dependence given by L1
- Parameters:
- XYZ: np.ndarray
Cartesian coordinates of the prediction points (meters)
- SXYZ: np.ndarray
Earth-centered Earth-fixed coordinates of the sun (meters)
- LXYZ: np.ndarray
Earth-centered Earth-fixed coordinates of the moon (meters)
- F2_solar: np.ndarray
Factors for the sun
- F2_lunar: np.ndarray
Factors for the moon
- pyTMD.predict._frequency_dependence_diurnal(XYZ: ndarray, MJD: ndarray)[source]
Computes the in-phase and out-of-phase corrections induced by mantle anelasticity in the diurnal band
- Parameters:
- XYZ: np.ndarray
Cartesian coordinates of the prediction points (meters)
- MJD: np.ndarray
Modified Julian Day (MJD)
- pyTMD.predict._frequency_dependence_long_period(XYZ: ndarray, MJD: ndarray)[source]
Computes the in-phase and out-of-phase corrections induced by mantle anelasticity in the long-period band
- Parameters:
- XYZ: np.ndarray
Cartesian coordinates of the prediction points (meters)
- MJD: np.ndarray
Modified Julian Day (MJD)
- pyTMD.predict._free_to_mean(XYZ: ndarray, h2: float | ndarray, l2: float | ndarray)[source]
Calculate offsets for converting the permanent tide from a tide-free to a mean-tide state
- Parameters:
- XYZ: np.ndarray
Cartesian coordinates of the prediction points (meters)
- h2: float or np.ndarray
Degree-2 Love number of vertical displacement
- l2: float or np.ndarray
Degree-2 Love (Shida) number of horizontal displacement