boule.Ellipsoid
Contents
boule.Ellipsoid#
- class boule.Ellipsoid(name, semimajor_axis, flattening, geocentric_grav_const, angular_velocity, long_name=None, reference=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 it’s 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
name (str) – A short name for the ellipsoid, for example
"WGS84"
.semimajor_axis (float) – The semimajor axis of the ellipsoid. The equatorial (large) radius. Definition: \(a\). Units: \(m\).
flattening (float) – The (first) flattening of the ellipsoid. Definition: \(f = (a - b)/a\). Units: adimensional.
geocentric_grav_const (float) – 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_velocity (float) – The angular velocity of the rotating ellipsoid. Definition: \(\omega\). Units: \(\\rad.s^{-1}\).
long_name (str or None) – A long name for the ellipsoid, for example
"World Geodetic System 1984"
(optional).reference (str or None) – Citation for the ellipsoid parameter values (optional).
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, B., & Moritz, H. (2006). Physical Geodesy " ... "(2nd, corr. ed. 2006 edition ed.). Wien ; New York: Springer." ... ), ... ) >>> print(ellipsoid) Ellipsoid(name='WGS84', ...) >>> 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") 6371008.7714 m >>> print(f"{ellipsoid.volume * 1e-9:.5e} km³") 1.08321e+12 km³ >>> print(f"{ellipsoid.gravity_equator:.10f} m/s²") 9.7803253359 m/s² >>> print(f"{ellipsoid.gravity_pole:.10f} m/s²") 9.8321849379 m/s²
Use the class methods for calculating normal gravity and other geometric quantities.
Attributes#
- 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 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 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: \(c = \sqrt{a^2 - b^2}\). Units: \(m\).
- Ellipsoid.mean_radius#
The arithmetic mean radius of the ellipsoid [Moritz1988]. Definition: \(R_1 = (2a + b)/3\). Units: \(m\).
- 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.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 c\). Units: \(m^3\).
Methods#
List of methods
|
Radial distance from the center of the ellipsoid to its surface. |
|
Convert from geodetic to geocentric spherical coordinates. |
|
Normal gravity of the ellipsoid at the given latitude and height. |
|
The prime vertical radius of curvature for a given geodetic latitude. |
|
Convert from geocentric spherical to geodetic coordinates. |
Methods documentation
- Ellipsoid.geocentric_radius(latitude, geodetic=True)[source]#
Radial distance from the center of the ellipsoid to its surface.
Can be calculated from either geocentric geodetic or geocentric spherical latitudes.
- Parameters
- Returns
geocentric_radius (float or array) – The geocentric radius for the given latitude(s) in the same units as the ellipsoid axis.
Tip
No elevation is taken into account. If you need the geocentric radius at a height other than zero, use
pymap3d.geodetic2spherical
instead.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_spherical(longitude, latitude, height)[source]#
Convert from geodetic to geocentric spherical coordinates.
The geodetic datum is defined by this ellipsoid. The coordinates are converted following [Vermeille2002].
- Parameters
longitude (array) – Longitude coordinates on geodetic coordinate system in degrees.
latitude (array) – Latitude coordinates on geodetic coordinate system in degrees.
height (array) – Ellipsoidal heights in meters.
- Returns
longitude (array) – Longitude coordinates on geocentric spherical coordinate system in degrees. The longitude coordinates are not modified during this conversion.
spherical_latitude (array) – Converted latitude coordinates on geocentric spherical coordinate system in degrees.
radius (array) – Converted spherical radius coordinates in meters.
- Ellipsoid.normal_gravity(latitude, height, si_units=False)[source]#
Normal gravity of the ellipsoid at the given latitude and height.
Computes the magnitude of the gradient of the gravity potential (gravitational + centrifugal; see [HofmannWellenhofMoritz2006]) generated by the 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.
Assumes that the internal density distribution of the ellipsoid is such that the gravity potential is constant at its surface.
Based on closed-form expressions by [Lakshmanan1991] and corrected by [LiGotze2001] which don’t require the free-air correction.
Caution
These expressions are only valid for heights on or above the surface of the ellipsoid.
- Parameters
- Returns
gamma (float or array) – The normal gravity in mGal or m/s².
- Ellipsoid.prime_vertical_radius(sinlat)[source]#
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
sinlat (float or array-like) – Sine of the geocentric geodetic latitude.
- Returns
prime_vertical_radius (float 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_geodetic(longitude, spherical_latitude, radius)[source]#
Convert from geocentric spherical to geodetic coordinates.
The geodetic datum is defined by this ellipsoid. The coordinates are converted following [Vermeille2002].
- Parameters
longitude (array) – Longitude coordinates on geocentric spherical coordinate system in degrees.
spherical_latitude (array) – Latitude coordinates on geocentric spherical coordinate system in degrees.
radius (array) – Spherical radius coordinates in meters.
- Returns
longitude (array) – Longitude coordinates on geodetic coordinate system in degrees. The longitude coordinates are not modified during this conversion.
latitude (array) – Converted latitude coordinates on geodetic coordinate system in degrees.
height (array) – Converted ellipsoidal height coordinates in meters.