Laplace’s tidal representation (pygeoid.tides.laplace)

class pygeoid.tides.laplace.LaplaceTidalEquation(bodies: list = ['moon', 'sun'], parts: list = ['zonal', 'tesseral', 'sectorial'], love: dict = {'h': <<class 'astropy.constants.constant.Constant'> name='Nominal degree 2 Love number h2' value=0.6078 uncertainty=0 unit='' reference='IERS Conventions(2010), IERS Technical Note 36, Verlagdes Bundesamts für Kartographie und Geodäsie, Frankfurt am Main, Germany.'>, 'k': <<class 'astropy.constants.constant.Constant'> name='Nominal degree 2 Love number k2' value=0.29525 uncertainty=0 unit='' reference='IERS Conventions(2010), IERS Technical Note 36, Verlagdes Bundesamts für Kartographie und Geodäsie, Frankfurt am Main, Germany.'>, 'l': <<class 'astropy.constants.constant.Constant'> name='Nominal degree 2 Shida number l2' value=0.0847 uncertainty=0 unit='' reference='IERS Conventions(2010), IERS Technical Note 36, Verlagdes Bundesamts für Kartographie und Geodäsie, Frankfurt am Main, Germany.'>}, gravimetric_factor: float = 1.1563, diminishing_factor: float = 0.6947, ephemeris: str = None)[source]

Class for Laplace’s tidal representation.

This class is a realision of the Laplace’s tidal representaion for the second-degree tidal potential:

W_2 = D (Z + T + S),

where D is called Doodson’s constant and Z, T and S are zonal, tesseral and sectorial parts of the tidal potential respectively. This representation is often applied for practical calculations of the total tidal potential and its derivatives, instead of more advanced and complicated harmonic decomposition.

If some tidal effect calculated from this class with all parts included (by default) and fot the elastic Earth is applied as a correction (effect taken with the negative sign), then the corrected value will be in the conventional tide free system, because a permanent part of the tide will also be removed. If mean or zero tide systems are desired, then permanent part of the tide should be restored in some way.

bodies

List of solar system bodies, the total tide from which will be taken into account. By default, only the tides from the Moon and the Sun is accounted since they are the largest.

Type:list, optional
parts

Which parts of Laplace’s tidal eqution are involved. There are ‘zonal’, ‘tesseral’ and ‘sectorial’ parts. By default, all parts are included.

Type:list, optional
love

Approximate Love numbers for elastic deformation calculations. This dictionary should contain ‘k’, ‘l’ and ‘h’ numbers of the second degree. Default values are taken from IERS Conventions 2010.

Type:dict, optional
gravimetric_factor

Gravimetric factor (delta). If None it will be calculated from Love numbers as delta = 1 + h - 3/2*k. Default value is approximate gravimetric factor from PREM model.

Type:float, optional
diminishing_factor

Diminishing factor (gamma). If None then it will be calculated from Love numbers as gamma = 1 + k - h. Default value is approximate diminishing factor from PREM model.

Type:float, optional
ephemeris

Ephemeris to use. If not given, use the one set with astropy.coordinates.solar_system_ephemeris.set (which is set to ‘builtin’ by default).

Type:str, optional

Notes

One can check which bodies are covered by a given ephemeris using:

>>> from astropy.coordinates import solar_system_ephemeris
>>> solar_system_ephemeris.bodies

References

[1]Melchior, P. (1983), The Tides of the Planet Earth. 2nd edition, Pergamon Press.
potential(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m")) -> Unit("m2 / s2")[source]

Return tidal potential.

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
Returns:

potential – Tidal potential.

Return type:

Quantity

potential_lat_derivative(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m")) -> Unit("m2 / s2")[source]

Return tidal potential latitude derivative.

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
Returns:

potential_dlat – Tidal potential latitude derivative.

Return type:

Quantity

potential_lon_derivative(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m")) -> Unit("m2 / s2")[source]

Return tidal potential longitude derivative.

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
Returns:

potential_dlon – Tidal potential longitude derivative.

Return type:

Quantity

potential_r_derivative(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m")) -> Unit("m / s2")[source]

Return tidal potential radial derivative.

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
Returns:

potential_dr – Tidal potential radial derivative.

Return type:

Quantity

gradient(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m")) -> Unit("m / s2")[source]

Return tidal potential gradient (tidal acceleration).

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
Returns:

grad – Tidal potential radial derivative.

Return type:

Quantity

displacement(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m"), gravity: Unit("m / s2") = <<class 'astropy.constants.codata2018.CODATA2018'> name='Standard acceleration of gravity' value=9.80665 uncertainty=0.0 unit='m / s2' reference='CODATA 2018'>, elastic: bool = True) -> Unit("m")[source]

Return direct tidal deformation for elastic or equilibrium Earth.

Direct effect is an equipotential surface deformation caused by tidal potential. Equilibrium tidal deformation is caused only by theoretical tide, elastic deformation takes into account the elastic properties of the Earth.

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
  • gravity (Quantity) – Mean global gravity of the Earth.
  • elastic (bool, optional) – If True, elastic deformation (with approximate Love numbers) will be returned, otherwise equilibrium deformation. Default is True.
Returns:

  • u_lat (~astropy.units.Quantity) – Deformation in the north-south direction.
  • u_lon (~astropy.units.Quantity) – Deformation in the east-west direction.
  • u_r (~astropy.units.Quantity) – Deformation in the radial direction.

deformation_potential(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m")) -> Unit("m2 / s2")[source]

Return deformation tidal potential.

An incremental (deforamtion) potential is produced by the re-distribution of masses which is factorized by a number, the second Love number k, on the earth’s surface.

This is also called the indirect effect.

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
Returns:

delta_potential – Deformation potential.

Return type:

Quantity

gravity(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m"), elastic: bool = True) -> Unit("m / s2")[source]

Return tidal gravity variation.

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
  • elastic (bool, optional) – If True then the Earth is elastic (deformable) and gravity change is multiplied by the gravimetric factor specified in the class instance.
Returns:

delta_g – Tidal gravity variation.

Return type:

Quantity

tilt(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m"), azimuth: Unit("deg") = None, gravity: Unit("m / s2") = <<class 'astropy.constants.codata2018.CODATA2018'> name='Standard acceleration of gravity' value=9.80665 uncertainty=0.0 unit='m / s2' reference='CODATA 2018'>, elastic: bool = True) -> Unit(dimensionless)[source]

Return tidal tilt.

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
  • azimuth (Quantity, optional) – If given then tilt will be returned in that direction. If None then full tilt angle will be returned. Default is None.
  • elastic (bool, optional) – If True then the Earth is elastic (deformable) and tilt is multiplied by the combination of Love numbers (h - (1 + k)).
  • gravity (Quantity, optional) – Mean global gravity of the Earth.
Returns:

tilt – Tidal tilt.

Return type:

Quantity

geoidal_height(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m"), gravity: Unit("m / s2") = <<class 'astropy.constants.codata2018.CODATA2018'> name='Standard acceleration of gravity' value=9.80665 uncertainty=0.0 unit='m / s2' reference='CODATA 2018'>) -> Unit("m")[source]

Return tidal variation of the geoidal height.

Geoidal heights are affected by direct and indirect effects.

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
  • gravity (Quantity, optional) – Mean global gravity of the Earth.
Returns:

delta_geoid – Geoidal height tidal variation.

Return type:

Quantity

physical_height(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m"), gravity: Unit("m / s2") = <<class 'astropy.constants.codata2018.CODATA2018'> name='Standard acceleration of gravity' value=9.80665 uncertainty=0.0 unit='m / s2' reference='CODATA 2018'>) -> Unit("m")[source]

Return tidal variation of the (absolute) physical heights.

Tidal variations of physically defined heights like orthometric, normal, or normal-orthometric heights can be derived as the difference between ellipsoidal and geoidal heights.

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
  • gravity (Quantity, optional) – Mean global gravity of the Earth.
Returns:

delta_h – Physical heights tidal variation.

Return type:

Quantity

physical_height_difference(time: astropy.time.core.Time, lat: Unit("deg"), lon: Unit("deg"), r: Unit("m"), azimuth: Unit("deg"), length: Unit("m"), gravity: Unit("m / s2") = <<class 'astropy.constants.codata2018.CODATA2018'> name='Standard acceleration of gravity' value=9.80665 uncertainty=0.0 unit='m / s2' reference='CODATA 2018'>) -> Unit("m")[source]

Return tidal variation of the physical heights difference.

Parameters:
  • time (Time) – Time of observation.
  • lat (Quantity) – Geocentric (spherical) latitude.
  • lon (Quantity) – Geocentric (spherical) longitude.
  • r (Quantity) – Geocentric radius (radial distance).
  • azimuth (Quantity) – Azimuth of the levelling line.
  • length (Quantity) – Length of the levelling line.
  • gravity (Quantity, optional) – Mean global gravity of the Earth.
Returns:

delta_h – Physical heights difference tidal variation.

Return type:

Quantity