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 tidal values following IERS Conventions
Calling Sequence
import pyTMD.predict
ht = pyTMD.predict.time_series(time, hc, con)
- pyTMD.predict.map(t: float | numpy.ndarray, hc: ndarray, constituents: list | numpy.ndarray, deltat: float | numpy.ndarray = 0.0, corrections: str = 'OTIS')[source]
Predict tides at a single time using harmonic constants [1]
- 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
- Returns
- ht: np.ndarray
tide values reconstructed using the nodal corrections
References
- 1
G. D. Egbert and S. Y. Erofeeva, “Efficient Inverse Modeling of Barotropic Ocean Tides,” Journal of Atmospheric and Oceanic Technology, 19(2), 183–204, (2002). doi: 10.1175/1520-0426(2002)019<0183:EIMOBO>2.0.CO;2
- pyTMD.predict.drift(t: float | numpy.ndarray, hc: ndarray, constituents: list | numpy.ndarray, deltat: float | numpy.ndarray = 0.0, corrections: str = 'OTIS')[source]
Predict tides at multiple times and locations using harmonic constants [1]
- 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
- Returns
- ht: np.ndarray
tidal time series reconstructed using the nodal corrections
References
- 1
G. D. Egbert and S. Y. Erofeeva, “Efficient Inverse Modeling of Barotropic Ocean Tides,” Journal of Atmospheric and Oceanic Technology, 19(2), 183–204, (2002). doi: 10.1175/1520-0426(2002)019<0183:EIMOBO>2.0.CO;2
- pyTMD.predict.time_series(t: float | numpy.ndarray, hc: ndarray, constituents: list | numpy.ndarray, deltat: float | numpy.ndarray = 0.0, corrections: str = 'OTIS')[source]
Predict tidal time series at a single location using harmonic constants [1]
- 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
- Returns
- ht: np.ndarray
tidal time series reconstructed using the nodal corrections
References
- 1
G. D. Egbert and S. Y. Erofeeva, “Efficient Inverse Modeling of Barotropic Ocean Tides,” Journal of Atmospheric and Oceanic Technology, 19(2), 183–204, (2002). doi: 10.1175/1520-0426(2002)019<0183:EIMOBO>2.0.CO;2
- pyTMD.predict.infer_minor(t: float | numpy.ndarray, zmajor: ndarray, constituents: list | numpy.ndarray, **kwargs)[source]
Infer the tidal values for minor constituents using their relation with major constituents [1] [2] [3] [4]
- 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
- Returns
- dh: np.ndarray
tidal time series for minor constituents
References
- 1
A. T. Doodson and H. D. Warburg, “Admiralty Manual of Tides”, HMSO, London, (1941).
- 2
P. Schureman, “Manual of Harmonic Analysis and Prediction of Tides,” US Coast and Geodetic Survey, Special Publication, 98, (1958).
- 3
M. G. G. Foreman and R. F. Henry, “The harmonic analysis of tidal model time series,” Advances in Water Resources, 12(3), 109–120, (1989). doi: 10.1016/0309-1708(89)90017-1
- 4
G. D. Egbert and S. Y. Erofeeva, “Efficient Inverse Modeling of Barotropic Ocean Tides,” Journal of Atmospheric and Oceanic Technology, 19(2), 183–204, (2002). doi: 10.1175/1520-0426(2002)019<0183:EIMOBO>2.0.CO;2
- pyTMD.predict.equilibrium_tide(t: ndarray, lat: ndarray)[source]
Compute the long-period equilibrium tides the summation of fifteen tidal spectral lines from Cartwright-Tayler-Edden tables [1] [2]
- Parameters
- t: np.ndarray
time (days relative to January 1, 1992)
- lat: np.ndarray
latitudes (degrees north)
- Returns
- lpet: np.ndarray
long-period equilibrium tide in meters
References
- 1
D. E. Cartwright and R. J. Tayler, “New Computations of the Tide-generating Potential,” Geophysical Journal of the Royal Astronomical Society, 23(1), 45–73. (1971). doi: 10.1111/j.1365-246X.1971.tb01803.x
- 2
D. E. Cartwright and A. C. Edden, “Corrected Tables of Tidal Harmonics,” Geophysical Journal of the Royal Astronomical Society, 33(3), 253–264, (1973). doi: 10.1111/j.1365-246X.1973.tb03420.x
- 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 [1] [2] [3] [4]
- Parameters
- t: np.ndarray
Time (days relative to January 1, 1992)
- XYZ: np.ndarray
Cartesian coordinates of the station (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
References
- 1
P. M. Mathews, B. A. Buffett, T. A. Herring and I. I Shapiro, “Forced nutations of the Earth: Influence of inner core dynamics: 1. Theory”, Journal of Geophysical Research: Solid Earth, 96(B5), 8219–8242, (1991). doi: 10.1029/90JB01955
- 2
P. M. Mathews, V. Dehant and J. M. Gipson, “Tidal station displacements”, Journal of Geophysical Research: Solid Earth, 102(B9), 20469–20477, (1997). doi: 10.1029/97JB01515
- 3
J. C. Ries, R. J. Eanes, C. K. Shum and M. M. Watkins, “Progress in the determination of the gravitational coefficient of the Earth”, Geophysical Research Letters, 19(6), 529–531, (1992). doi: 10.1029/92GL00259
- 4
J. M. Wahr, “Body tides on an elliptical, rotating, elastic and oceanless Earth”, Geophysical Journal of the Royal Astronomical Society, 64(3), 677–703, (1981). doi: 10.1111/j.1365-246X.1981.tb02690.x
- 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 station (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 station (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 station (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 station (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 station (meters)
- MJD: np.ndarray
Modified Julian Day (MJD)
- pyTMD.predict._free_to_mean(XYZ: ndarray, h2: float | numpy.ndarray, l2: float | numpy.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 station (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