boule.Ellipsoid#

class boule.Ellipsoid(name, semimajor_axis, flattening, geocentric_grav_const, angular_velocity, long_name=None, reference=None, comments=None)[source]#

A rotating oblate ellipsoid.

The ellipsoid is defined by four parameters: semimajor axis, flattening, geocentric gravitational constant, and angular velocity. It spins around its semiminor axis and has constant gravity potential at its surface. The internal density structure of the ellipsoid is unspecified but must be such that the constant potential condition is satisfied.

This class is read-only: Input parameters and attributes cannot be changed after instantiation.

Units: All input parameters and derived attributes are in SI units.

Parameters:
namestr

A short name for the ellipsoid, for example "WGS84".

semimajor_axisfloat

The semimajor axis of the ellipsoid. The equatorial (large) radius. Definition: \(a\). Units: \(m\).

flatteningfloat

The (first) flattening of the ellipsoid. Definition: \(f = (a - b)/a\). Units: adimensional.

geocentric_grav_constfloat

The geocentric gravitational constant. The product of the mass of the ellipsoid \(M\) and the gravitational constant \(G\). Definition: \(GM\). Units: \(m^3.s^{-2}\).

angular_velocityfloat

The angular velocity of the rotating ellipsoid. Definition: \(\omega\). Units: \(\\rad.s^{-1}\).

long_namestr or None

A long name for the ellipsoid, for example "World Geodetic System 1984" (optional).

referencestr or None

Citation for the ellipsoid parameter values (optional).

commentsstr or None

Additional comments regarding the ellipsoid (optional).

Notes

Caution

Use boule.Sphere if you desire zero flattening because there are singularities for this particular case in the normal gravity calculations.

Examples

We can define an ellipsoid by setting the 4 key numerical parameters and some metadata about where they came from:

>>> ellipsoid = Ellipsoid(
...     name="WGS84",
...     long_name="World Geodetic System 1984",
...     semimajor_axis=6378137,
...     flattening=1 / 298.257223563,
...     geocentric_grav_const=3986004.418e8,
...     angular_velocity=7292115e-11,
...     reference="Hofmann-Wellenhof & Moritz (2006)",
...     comments="This is the same as the boule WGS84 ellipsoid.",
... )
>>> print(ellipsoid) 
WGS84 - World Geodetic System 1984
Oblate ellipsoid:
  • Semimajor axis: 6378137 m
  • Flattening: 0.0033528106647474805
  • GM: 398600441800000.0 m³/s²
  • Angular velocity: 7.292115e-05 rad/s
Source:
  Hofmann-Wellenhof & Moritz (2006)
Comments:
  This is the same as the boule WGS84 ellipsoid.
>>> print(ellipsoid.long_name)
World Geodetic System 1984

The class then defines several derived attributes based on the input parameters:

>>> print(f"{ellipsoid.semiminor_axis:.4f} m")
6356752.3142 m
>>> print(f"{ellipsoid.linear_eccentricity:.8f} m")
521854.00842339 m
>>> print(f"{ellipsoid.first_eccentricity:.13e}")
8.1819190842621e-02
>>> print(f"{ellipsoid.second_eccentricity:.13e}")
8.2094437949696e-02
>>> print(f"{ellipsoid.mean_radius:.4f} m")
6370994.4018 m
>>> print(f"{ellipsoid.semiaxes_mean_radius:.4f} m")
6371008.7714 m
>>> print(f"{ellipsoid.volume_equivalent_radius:.4f} m")
6371000.7900 m
>>> print(f"{ellipsoid.mass:.10e} kg")
5.9721684941e+24 kg
>>> print(f"{ellipsoid.mean_density:.0f} kg/m³")
5513 kg/m³
>>> print(f"{ellipsoid.volume * 1e-9:.5e} km³")
1.08321e+12 km³
>>> print(f"{ellipsoid.area:.10e} m²")
5.1006562172e+14 m²
>>> print(f"{ellipsoid.area_equivalent_radius:0.4f} m")
6371007.1809 m
>>> print(f"{ellipsoid.gravity_equator:.10f} m/s²")
9.7803253359 m/s²
>>> print(f"{ellipsoid.gravity_pole:.10f} m/s²")
9.8321849379 m/s²
>>> print(f"{ellipsoid.reference_normal_gravity_potential:.3f} m²/s²")
62636851.715 m²/s²

Use the class methods for calculating normal gravity and other geometric quantities.

Attributes:
area

The area of the ellipsoid.

area_equivalent_radius

The area equivalent radius of the ellipsoid.

eccentricity

Alias for the first eccentricity.

first_eccentricity

The (first) eccentricity of the ellipsoid.

gravity_equator

The normal gravity at the equator.

gravity_pole

The normal gravity at the pole.

linear_eccentricity

The linear eccentricity of the ellipsoid.

mass

The mass of the ellipsoid.

mean_density

The mean density of the ellipsoid.

mean_radius

The mean radius of the ellipsoid.

reference_normal_gravity_potential

The normal gravity potential on the surface of the ellipsoid.

second_eccentricity

The second eccentricity of the ellipsoid.

semiaxes_mean_radius

The arithmetic mean radius of the ellipsoid semi-axes [Moritz1988].

semimajor_axis_longitude

The semimajor axis longitude of the ellipsoid is equal to zero.

semimedium_axis

The semimedium axis of the ellipsoid is equal to its semimajor axis.

semiminor_axis

The semiminor (small/polar) axis of the ellipsoid.

thirdflattening

The third flattening of the ellipsoid (used in geodetic calculations).

volume

The volume bounded by the ellipsoid.

volume_equivalent_radius

The volume equivalent radius of the ellipsoid.

Methods

cartesian_to_ellipsoidal_harmonic(coordinates)

Convert from geocentric Cartesian to ellipsoidal harmonic coordinates.

cartesian_to_geodetic(coordinates)

Convert from geocentric Cartesian to geodetic coordinates.

cartesian_to_spherical(coordinates)

Convert from geocentric Cartesian coordinates to spherical coordinates.

centrifugal_potential(coordinates, *[, ...])

Centrifugal potential of the rotating ellipsoid.

ellipsoidal_harmonic_to_cartesian(coordinates)

Convert from ellipsoidal harmonic to geocentric Cartesian coordinates.

ellipsoidal_harmonic_to_geodetic(coordinates)

Convert from ellipsoidal-harmonic coordinates to geodetic coordinates.

ellipsoidal_harmonic_to_spherical(coordinates)

Convert from ellipsoidal harmonic to geocentric spherical coordinates.

geocentric_radius(latitude, *[, ...])

Radial distance from the center of the ellipsoid to its surface.

geodetic_to_cartesian(coordinates)

Convert from geodetic to geocentric Cartesian coordinates.

geodetic_to_ellipsoidal_harmonic(coordinates)

Convert from geodetic to ellipsoidal harmonic coordinates.

geodetic_to_spherical(coordinates)

Convert from geodetic to geocentric spherical coordinates.

normal_gravitational_potential(coordinates, *)

Calculate the normal gravitational potential of the ellipsoid.

normal_gravity(coordinates, *[, ...])

Calculate the normal gravity of the ellipsoid.

normal_gravity_potential(coordinates, *[, ...])

Calculate the normal gravity potential of the ellipsoid.

prime_vertical_radius(sinlat)

Calculate the prime vertical radius of curvature for a given geodetic latitude.

spherical_to_cartesian(coordinates)

Convert from spherical coordinates to geocentric Cartesian coordinates.

spherical_to_ellipsoidal_harmonic(coordinates)

Convert from geocentric spherical to ellipsoidal harmonic coordinates.

spherical_to_geodetic(coordinates)

Convert from geocentric spherical to geodetic coordinates.

Attributes#

Ellipsoid.area#

The area of the ellipsoid.

Definition: \(A = 2 \pi a^2 \left(1 + \dfrac{b^2}{e a^2} \text{arctanh}\,e \right)\).

Units: \(m^2\).

Ellipsoid.area_equivalent_radius#

The area equivalent radius of the ellipsoid.

Definition: \(R_2 = \sqrt{A / (4 \pi)}\).

Units: \(m\).

Ellipsoid.eccentricity#

Alias for the first eccentricity.

Ellipsoid.first_eccentricity#

The (first) eccentricity of the ellipsoid.

The ratio between the linear eccentricity and the semimajor axis.

Definition: \(e = \dfrac{\sqrt{a^2 - b^2}}{a} = \sqrt{2f - f^2}\).

Units: adimensional.

Ellipsoid.gravity_equator#

The normal gravity at the equator.

This is the norm of the gravity acceleration vector (gravitational + centrifugal accelerations) at the equator on the surface of the ellipsoid.

Units: \(m/s^2\).

Ellipsoid.gravity_pole#

The normal gravity at the pole.

This is the norm of the gravity acceleration vector (gravitational + centrifugal accelerations) at the poles on the surface of the ellipsoid.

Units: \(m/s^2\).

Ellipsoid.linear_eccentricity#

The linear eccentricity of the ellipsoid.

The distance between the ellipsoid’s center and one of its foci.

Definition: \(E = \sqrt{a^2 - b^2}\).

Units: \(m\).

Ellipsoid.mass#

The mass of the ellipsoid.

Definition: \(M = GM / G\).

Units: \(kg\).

Ellipsoid.mean_density#

The mean density of the ellipsoid.

Definition: \(\rho = M / V\).

Units: \(kg / m^3\).

Ellipsoid.mean_radius#

The mean radius of the ellipsoid.

This is equivalent to the degree 0 spherical harmonic coefficient of the ellipsoid shape.

Definition: \(R_0 = \dfrac{1}{4 \pi} {\displaystyle \int_0^{\pi} \int_0^{2 \pi}} r(\theta) \sin \theta \, d\theta \, d\lambda\)

in which \(r\) is the ellipsoid spherical radius, \(\theta\) is spherical latitude, and \(\lambda\) is spherical longitude.

Units: \(m\).

Ellipsoid.reference_normal_gravity_potential#

The normal gravity potential on the surface of the ellipsoid.

Definition: \(U_0 = \dfrac{GM}{E} \arctan{\dfrac{E}{b}} + \dfrac{1}{3} \omega^2 a^2\).

Units: \(m^2 / s^2\).

Ellipsoid.second_eccentricity#

The second eccentricity of the ellipsoid.

The ratio between the linear eccentricity and the semiminor axis.

Definition: \(e^\prime = \dfrac{\sqrt{a^2 - b^2}}{b} = \dfrac{\sqrt{2f - f^2}}{1 - f}\).

Units: adimensional.

Ellipsoid.semiaxes_mean_radius#

The arithmetic mean radius of the ellipsoid semi-axes [Moritz1988].

Definition: \(R_1 = (2a + b)/3\).

Units: \(m\).

Ellipsoid.semimajor_axis_longitude#

The semimajor axis longitude of the ellipsoid is equal to zero.

Definition: \(\lambda_a = 0\).

Units: \(m\).

Ellipsoid.semimedium_axis#

The semimedium axis of the ellipsoid is equal to its semimajor axis.

Units: \(m\).

Ellipsoid.semiminor_axis#

The semiminor (small/polar) axis of the ellipsoid.

Definition: \(b = a (1 - f)\).

Units: \(m\).

Ellipsoid.thirdflattening#

The third flattening of the ellipsoid (used in geodetic calculations).

Definition: \(f^{\prime\prime}= \dfrac{a -b}{a + b}\).

Units: adimensional.

Ellipsoid.volume#

The volume bounded by the ellipsoid.

Definition: \(V = \dfrac{4}{3} \pi a^2 b\).

Units: \(m^3\).

Ellipsoid.volume_equivalent_radius#

The volume equivalent radius of the ellipsoid.

Definition: \(R_3 = \left(\dfrac{3}{4 \pi} V \right)^{1/3}\).

Units: \(m\).

Methods#

Ellipsoid.cartesian_to_ellipsoidal_harmonic(coordinates)[source]#

Convert from geocentric Cartesian to ellipsoidal harmonic coordinates.

In the geocentric Cartesian system, z is aligned with the Earth’s axis of rotation and points towards the North pole, x and y are contained in the equatorial plane, x coincides with the Greenwich meridian, and y completes right-handed coordinate system.

The conversion is done by passing through the geodetic coordinate system.

Parameters:
coordinatestuple = (x, y, z)

The x, y, z coordinates of the points in a geocentric Cartesian coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. All coordinates must be in meters.

Returns:
converted_coordinatestuple = (longitude, latitude_reduced, u)

The converted longitude, reduced (or parametric) latitude, and the coordinate u (the semiminor axis of the ellipsoid that passes through the input coordinates) in a ellipsoidal harmonic coordinate system. The shape of each element will be compatible with the shape of the inputs. Longitude and latitude will be in degrees and u in meters.

Ellipsoid.cartesian_to_geodetic(coordinates)[source]#

Convert from geocentric Cartesian to geodetic coordinates.

The geodetic datum is defined by this ellipsoid. In the geocentric Cartesian system, z is aligned with the Earth’s axis of rotation and points towards the North pole, x and y are contained in the equatorial plane, x coincides with the Greenwich meridian, and y completes right-handed coordinate system.

Uses the non-iterative method of [Zhu1993] to make the conversion. This is accurate for points further than 43 km from the center of the ellipsoid.

Parameters:
coordinatestuple = (x, y, z)

The x, y, z coordinates of the points in a geocentric Cartesian coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. All coordinates must be in meters.

Returns:
converted_coordinatestuple = (longitude, latitude_geodetic, height)

The converted longitude, geodetic latitude, and geometric height in a geodetic coordinate system. The shape of each element will be compatible with the shape of the inputs. Longitude and latitude will be in degrees and height in meters.

References

[Zhu1993]

Ellipsoid.cartesian_to_spherical(coordinates)[source]#

Convert from geocentric Cartesian coordinates to spherical coordinates.

In the geocentric Cartesian system, z is aligned with the Earth’s axis of rotation and points towards the North pole, x and y are contained in the equatorial plane, x coincides with the Greenwich meridian, and y completes right-handed coordinate system.

Parameters:
coordinatestuple = (x, y, z)

The x, y, z coordinates of the points in a geocentric Cartesian coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. All coordinates must be in meters.

Returns:
converted_coordinatestuple = (longitude, latitude_spherical, radius)

The converted longitude, spherical latitude, and radius coordinates in a geocentric spherical coordinate system. The shape of each element will be compatible with the shape of the inputs. Longitude and latitude are in degrees and radius in meters

Ellipsoid.centrifugal_potential(coordinates, *, coordinate_system='geodetic')[source]#

Centrifugal potential of the rotating ellipsoid.

Calculate the centrifugal potential due to the rotation of the ellipsoid about its semiminor axis at the given points.

Parameters:
coordinatestuple = (coordinate1, coordinate2, coordinate3)

Tuple with 3 arrays containing the coordinates of the computation points. The meaning of the arrays is determined by the coordinate_system argument: longitude, geodetic latitude, and geometric height for a geodetic system; longitude, geocentric latitude, and radius for a geocentric spherical system; longitude, reduced latitude, and u for an ellipsoidal harmonic system; x, y, and z for a geocentric Cartesian system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitudes must be in degrees and height, radius, and u in meters. Since longitude is not used in computations (the potential is symmetric with longitude), it can be assigned None.

coordinate_systemstr

The coordinate system that will be assumed for the given coordinates. Should be one of: "geodetic" (default), "spherical", "cartesian", or "ellipsoidal harmonic".

Returns:
Phifloat or array

The centrifugal potential in m²/s².

Notes

The centrifugal potential \(\Phi\) at geodetic latitude \(\phi\) and height above the ellipsoid \(h\) (geometric height) is

\[\Phi(\phi, h) = \dfrac{1}{2} \omega^2 \left(N(\phi) + h\right)^2 \cos^2(\phi)\]

in which \(N(\phi)\) is the prime vertical radius of curvature of the ellipsoid and \(\omega\) is the angular velocity.

In geocentric Cartesian coordinates, the potential is

\[\Phi(x, y) = \dfrac{1}{2} \omega^2 \left(x^2 + y^2\right)^2\]

and in geocentric spherical coordinates, the potential is

\[\Phi(\theta, r) = \dfrac{1}{2} \omega^2 r^2 \cos^2\theta\]

in which \(\theta\) is the geocentric latitude and \(r\) is the radius.

For inputs in ellipsoidal harmonic coordinates, the coordinates will be converted to geodetic before calculation of the potential.

Ellipsoid.ellipsoidal_harmonic_to_cartesian(coordinates)[source]#

Convert from ellipsoidal harmonic to geocentric Cartesian coordinates.

In the geocentric Cartesian system, z is aligned with the Earth’s axis of rotation and points towards the North pole, x and y are contained in the equatorial plane, x coincides with the Greenwich meridian, and y completes right-handed coordinate system.

Uses the conversion specified in [HofmannWellenhofMoritz2006].

Parameters:
coordinatestuple = (longitude, latitude_reduced, u)

Longitude, reduced (or parametric) latitude, and u (the semiminor axis of the ellipsoid that passes through the input coordinates) coordinates of the points in a ellipsoidal harmonic coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitude must be in degrees and u in meters.

Returns:
converted_coordinatestuple = (x, y, z)

The converted x, y, z coordinates in a geocentric Cartesian coordinate system. The shape of each element will be compatible with the shape of the inputs. All coordinates are in meters.

References

[HofmannWellenhofMoritz2006]

Ellipsoid.ellipsoidal_harmonic_to_geodetic(coordinates)[source]#

Convert from ellipsoidal-harmonic coordinates to geodetic coordinates.

The geodetic datum is defined by this ellipsoid.

Parameters:
coordinatestuple = (longitude, latitude_reduced, u)

Longitude, reduced (or parametric) latitude, and u (the semiminor axis of the ellipsoid that passes through the input coordinates) coordinates of the points in a ellipsoidal harmonic coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitude must be in degrees and u in meters. Since longitude is not affected by conversions, it can be assigned None.

Returns:
converted_coordinatestuple = (longitude, latitude_geodetic, height)

The converted longitude, geodetic latitude, and geometric height in a geodetic coordinate system. The shape of each element will be compatible with the shape of the inputs. If the input longitude is None, the output will also be None. Longitude and latitude will be in degrees and height in meters.

Ellipsoid.ellipsoidal_harmonic_to_spherical(coordinates)[source]#

Convert from ellipsoidal harmonic to geocentric spherical coordinates.

The conversion is performed by passing through the geodetic system. Because of this, there is a loss of accuracy when doing the round-trip between spherical and ellipsoidal harmonic coordinates. The latitude conversion is good to approximately 0.001 degrees and radius to approximately 1 meter.

Parameters:
coordinatestuple = (longitude, latitude_reduced, u)

Longitude, reduced (or parametric) latitude, and u (the semiminor axis of the ellipsoid that passes through the input coordinates) coordinates of the points in a ellipsoidal harmonic coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitude must be in degrees and u in meters. Since longitude is the same in both systems and is not used in the computations, it can be passed as None if desired.

Returns:
converted_coordinatestuple = (longitude, latitude_spherical, radius)

The converted longitude, geocentric spherical latitude, and radius in a geocentric spherical coordinate system. The shape of each element will be compatible with the shape of the inputs. If the input longitude is None, the output will also be None. Longitude and latitude will be in degrees and radius in meters.

Ellipsoid.geocentric_radius(latitude, *, coordinate_system='geodetic')[source]#

Radial distance from the center of the ellipsoid to its surface.

Can be calculated from either geocentric geodetic or geocentric spherical latitudes.

Parameters:
latitudefloat or array

Latitude coordinates on geodetic coordinate system in degrees.

coordinate_systemstr

The coordinate system that will be assumed for the given latitude. Should be one of: "geodetic" (default) for geodetic latitudes or "spherical" for geocentric spherical latitudes.

Returns:
geocentric_radiusfloat or array

The geocentric radius for the given latitude(s) in the same units as the ellipsoid axis.

Notes

The geocentric surface radius \(R\) is a function of the geocentric geodetic latitude \(\phi\) and the semimajor and semiminor axis, \(a\) and \(b\) [1]:

\[R(\phi) = \sqrt{ \dfrac{ (a^2\cos\phi)^2 + (b^2\sin\phi)^2 }{ (a\cos\phi)^2 + (b\sin\phi)^2 } }\]

Alternatively, the geocentric surface radius can also be calculated using the geocentric spherical latitude \(\theta\) by passing geodetic=False:

\[R(\theta) = \sqrt{ \dfrac{ 1 }{ (\frac{\cos\theta}{a})^2 + (\frac{\sin\theta}{b})^2 } }\]

This can be useful if you already have the geocentric spherical latitudes and need the geocentric radius of the ellipsoid (for example, in spherical harmonic synthesis). In these cases, the coordinate conversion route is not possible since we need the radial coordinates to do that in the first place.

References

Ellipsoid.geodetic_to_cartesian(coordinates)[source]#

Convert from geodetic to geocentric Cartesian coordinates.

The geodetic datum is defined by this ellipsoid. In the geocentric Cartesian system, z is aligned with the Earth’s axis of rotation and points towards the North pole, x and y are contained in the equatorial plane, x coincides with the Greenwich meridian, and y completes right-handed coordinate system.

Parameters:
coordinatestuple = (longitude, latitude_geodetic, height)

Longitude, latitude, and geometric height coordinates of the points in a geodetic coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitude must be in degrees and height in meters.

Returns:
converted_coordinatestuple = (x, y, z)

The converted x, y, z coordinates in a geocentric Cartesian coordinate system. The shape of each element will be compatible with the shape of the inputs. All coordinates are in meters.

Ellipsoid.geodetic_to_ellipsoidal_harmonic(coordinates)[source]#

Convert from geodetic to ellipsoidal harmonic coordinates.

The geodetic datum is defined by this ellipsoid, and the coordinates are converted following [Lakshmanan1991] and [LiGotze2001].

Parameters:
coordinatestuple = (longitude, latitude_geodetic, height)

Longitude, latitude, and geometric height coordinates of the points in a geodetic coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitude must be in degrees and height in meters. Since longitude is not affected by conversions, it can be assigned None.

Returns:
converted_coordinatestuple = (longitude, latitude_reduced, u)

The converted longitude, reduced (or parametric) latitude, and the coordinate u (the semiminor axis of the ellipsoid that passes through the input coordinates) in a ellipsoidal harmonic coordinate system. The shape of each element will be compatible with the shape of the inputs. If the input longitude is None, the output will also be None. Longitude and latitude will be in degrees and u in meters.

Ellipsoid.geodetic_to_spherical(coordinates)[source]#

Convert from geodetic to geocentric spherical coordinates.

The geodetic datum is defined by this ellipsoid. The coordinates are converted following [Vermeille2002].

Parameters:
coordinatestuple = (longitude, latitude_geodetic, height)

Longitude, latitude, and geometric height coordinates of the points in a geodetic coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitude must be in degrees and height in meters. Since longitude is not affected by conversions, it can be assigned None.

Returns:
converted_coordinatestuple = (longitude, latitude_spherical, radius)

The converted longitude, geocentric spherical latitude, and radius in a geocentric spherical coordinate system. The shape of each element will be compatible with the shape of the inputs. If the input longitude is None, the output will also be None. Longitude and latitude will be in degrees and radius in meters.

Ellipsoid.normal_gravitational_potential(coordinates, *, coordinate_system='geodetic')[source]#

Calculate the normal gravitational potential of the ellipsoid.

Computes the gravitational potential generated by the ellipsoid at the given points. Does not include the centrifugal potential. See normal_gravity_potential for a version that includes the centrifugal component.

Caution

These expressions are only valid for heights on or above the surface of the ellipsoid.

Parameters:
coordinatestuple = (coordinate1, coordinate2, coordinate3)

Tuple with 3 arrays containing the coordinates of the computation points. The meaning of the arrays is determined by the coordinate_system argument: longitude, geodetic latitude, and geometric height for a geodetic system; longitude, geocentric latitude, and radius for a geocentric spherical system; longitude, reduced latitude, and u for an ellipsoidal harmonic system; x, y, and z for a geocentric Cartesian system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitudes must be in degrees and height, radius, and u in meters. Since longitude is not used in computations (the potential is symmetric with longitude), it can be assigned None.

coordinate_systemstr

The coordinate system that will be assumed for the given coordinates. Should be one of: "geodetic" (default), "spherical", "cartesian", or "ellipsoidal harmonic".

Returns:
Vfloat or array

The normal gravitational potential in m²/s².

Notes

Computes the gravitational potential generated by the ellipsoid at the given geodetic latitude \(\phi\) and height above the ellipsoid \(h\) (geometric height).

\[V(\beta, u) = \dfrac{GM}{E} \arctan{\dfrac{E}{u}} + \dfrac{1}{3} \omega^2 a^2 \dfrac{q}{q_0} P_2 (\sin \beta)\]

in which \(V\) is the gravitational potential of the ellipsoid (no centrifugal term), \(GM\) is the geocentric gravitational constant, \(E\) is the linear eccentricity, \(\omega\) is the angular rotation rate, \(q\) and \(q_0\) are auxiliary functions, \(P_2\) is the degree 2 unnormalized Legendre Polynomial, and \(u\) and \(\beta\) are ellipsoidal-harmonic coordinates corresponding to the input geodetic latitude and ellipsoidal height. See eq. 2-124 of [HofmannWellenhofMoritz2006].

Assumes that the internal density distribution of the ellipsoid is such that the gravity potential is constant at its surface.

Note

Since the calculations happen in ellipsoidal harmonic coordinates, passing inputs in any other coordinate system will require conversion, which may slow down computations if done in a loop.

Ellipsoid.normal_gravity(coordinates, *, coordinate_system='geodetic', si_units=False)[source]#

Calculate the normal gravity of the ellipsoid.

Computes the magnitude of the gradient of the gravity potential generated by this ellipsoid at any point outside the ellipsoid. Based on the closed-form expressions by [Lakshmanan1991] and corrected by [LiGotze2001]. This means that the free-air correction is not necessary to calculate normal gravity at the observation points.

Caution

These expressions are only valid for heights on or above the surface of the ellipsoid.

Parameters:
coordinatestuple = (coordinate1, coordinate2, coordinate3)

Tuple with 3 arrays containing the coordinates of the computation points. The meaning of the arrays is determined by the coordinate_system argument: longitude, geodetic latitude, and geometric height for a geodetic system; longitude, geocentric latitude, and radius for a geocentric spherical system; longitude, reduced latitude, and u for an ellipsoidal harmonic system; x, y, and z for a geocentric Cartesian system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitudes must be in degrees and height, radius, and u in meters. Since longitude is not used in computations (the potential is symmetric with longitude), it can be assigned None.

coordinate_systemstr

The coordinate system that will be assumed for the given coordinates. Should be one of: "geodetic" (default), "spherical", "cartesian", or "ellipsoidal harmonic".

si_unitsbool

Return the value in mGal (False, default) or m/s² (True).

Returns:
gammafloat or array

The normal gravity in mGal or m/s².

Notes

Normal gravity is defined as the magnitude of the gradient of the gravity potential generated by a reference ellipsoid at the given geodetic latitude \(\phi\) and height above the ellipsoid \(h\) (geometric height).

\[\gamma(\phi, h) = \|\vec{\nabla}U(\phi, h)\|\]

in which \(U = V + \Phi\) is the gravity potential of the ellipsoid, \(V\) is the gravitational potential of the ellipsoid, and \(\Phi\) is the centrifugal potential.

The equations used here assume that the internal density distribution of the ellipsoid is such that the gravity potential is constant at its surface. The specific internal density distribution is undefined.

Note

Since the calculations happen in ellipsoidal harmonic coordinates, passing inputs in any other coordinate system will require conversion, which may slow down computations if done in a loop.

Ellipsoid.normal_gravity_potential(coordinates, *, coordinate_system='geodetic')[source]#

Calculate the normal gravity potential of the ellipsoid.

Computes the gravity potential (gravitational + centrifugal) generated by the ellipsoid at the given points. See normal_gravitational_potential for a version that is purely gravitational.

Caution

These expressions are only valid for heights on or above the surface of the ellipsoid.

Parameters:
coordinatestuple = (coordinate1, coordinate2, coordinate3)

Tuple with 3 arrays containing the coordinates of the computation points. The meaning of the arrays is determined by the coordinate_system argument: longitude, geodetic latitude, and geometric height for a geodetic system; longitude, geocentric latitude, and radius for a geocentric spherical system; longitude, reduced latitude, and u for an ellipsoidal harmonic system; x, y, and z for a geocentric Cartesian system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitudes must be in degrees and height, radius, and u in meters. Since longitude is not used in computations (the potential is symmetric with longitude), it can be assigned None.

coordinate_systemstr

The coordinate system that will be assumed for the given coordinates. Should be one of: "geodetic" (default), "spherical", "cartesian", or "ellipsoidal harmonic".

Returns:
Vfloat or array

The normal gravity potential in m²/s².

Notes

Computes the gravity potential generated by the ellipsoid at the given geodetic latitude \(\phi\) and height above the ellipsoid \(h\) (geometric height).

\[U(\beta, u) = \dfrac{GM}{E} \arctan{\dfrac{E}{u}} + \dfrac{1}{2} \omega^2 a^2 \dfrac{q}{q_0} \left(\sin^2 \beta - \dfrac{1}{3}\right) + \dfrac{1}{2} \omega^2 \left(u^2 + E^2\right) \cos^2 \beta\]

in which \(U\) is the gravity potential of the ellipsoid, \(GM\) is the geocentric gravitational constant, \(E\) is the linear eccentricity, \(\omega\) is the angular rotation rate, \(q\) and \(q_0\) are auxiliary functions, and \(u\) and \(\beta\) are ellipsoidal-harmonic coordinates corresponding to the input geodetic latitude and ellipsoidal height. See eq. 2-126 of [HofmannWellenhofMoritz2006].

Assumes that the internal density distribution of the ellipsoid is such that the gravity potential is constant at its surface.

Note

Since the calculations happen in ellipsoidal harmonic coordinates, passing inputs in any other coordinate system will require conversion, which may slow down computations if done in a loop.

Ellipsoid.prime_vertical_radius(sinlat)[source]#

Calculate the prime vertical radius of curvature for a given geodetic latitude.

Note

This function receives the sine of the latitude as input to avoid repeated computations of trigonometric functions in methods/functions that rely on it.

Parameters:
sinlatfloat or array_like

Sine of the geocentric geodetic latitude.

Returns:
prime_vertical_radiusfloat or array_like

Prime vertical radius given in the same units as the semi-major axis

Notes

The prime vertical radius of curvature \(N\) is defined as [2]:

\[N(\phi) = \frac{a}{\sqrt{1 - e^2 \sin^2(\phi)}}\]

Where \(a\) is the semimajor axis and \(e\) is the first eccentricity.

References

Ellipsoid.spherical_to_cartesian(coordinates)[source]#

Convert from spherical coordinates to geocentric Cartesian coordinates.

In the geocentric Cartesian system, z is aligned with the Earth’s axis of rotation and points towards the North pole, x and y are contained in the equatorial plane, x coincides with the Greenwich meridian, and y completes right-handed coordinate system.

Parameters:
coordinatestuple = (longitude, latitude_spherical, height)

Longitude, latitude, and radius coordinates of the points in a geocentric spherical coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitude must be in degrees and radius in meters.

Returns:
converted_coordinatestuple = (x, y, z)

The converted x, y, z coordinates in a geocentric Cartesian coordinate system. The shape of each element will be compatible with the shape of the inputs. All coordinates are in meters.

Ellipsoid.spherical_to_ellipsoidal_harmonic(coordinates)[source]#

Convert from geocentric spherical to ellipsoidal harmonic coordinates.

The conversion is performed by passing through the geodetic system. Because of this, there is a loss of accuracy when doing the round-trip between spherical and ellipsoidal harmonic coordinates. The latitude conversion is good to approximately 0.001 degrees and radius to approximately 1 meter.

Parameters:
coordinatestuple = (longitude, latitude_spherical, height)

Longitude, latitude, and radius coordinates of the points in a geocentric spherical coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitude must be in degrees and radius in meters. Since longitude is not affected by conversions, it can be assigned None.

Returns:
converted_coordinatestuple = (longitude, latitude_reduced, u)

The converted longitude, reduced (or parametric) latitude, and the coordinate u (the semiminor axis of the ellipsoid that passes through the input coordinates) in a ellipsoidal harmonic coordinate system. The shape of each element will be compatible with the shape of the inputs. Longitude and latitude will be in degrees and u in meters. If the input longitude is None, the output will also be None.

Ellipsoid.spherical_to_geodetic(coordinates)[source]#

Convert from geocentric spherical to geodetic coordinates.

The geodetic datum is defined by this ellipsoid. The coordinates are converted following [Vermeille2002].

Parameters:
coordinatestuple = (longitude, latitude_spherical, height)

Longitude, latitude, and radius coordinates of the points in a geocentric spherical coordinate system. Each element can be a single number or an array. The shape of the arrays must be compatible. Longitude and latitude must be in degrees and radius in meters. Since longitude is not affected by conversions, it can be assigned None.

Returns:
converted_coordinatestuple = (longitude, latitude_geodetic, height)

The converted longitude, geodetic latitude, and geometric height in a geodetic coordinate system. The shape of each element will be compatible with the shape of the inputs. If the input longitude is None, the output will also be None. Longitude and latitude will be in degrees and height in meters.