boule.Ellipsoid

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

Reference oblate ellipsoid.

The ellipsoid is oblate and spins around it’s minor axis. It is defined by four parameters and offers other derived quantities as read-only properties. In fact, all attributes of this class are read-only and cannot be changed after instantiation.

All ellipsoid parameters are in SI units.

Parameters
  • name (str) – A short name for the ellipsoid, for example 'WGS84'.

  • semimajor_axis (float) – The semi-major axis of the ellipsoid (equatorial radius), usually represented by “a” [meters].

  • flattening (float) – The flattening of the ellipsoid (f) [adimensional].

  • geocentric_grav_const (float) – The geocentric gravitational constant (GM) [m^3 s^-2].

  • angular_velocity (float) – The angular velocity of the rotating ellipsoid (omega) [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).

Examples

We can define a reference unit sphere by using 0 as the flattening:

>>> sphere = Ellipsoid(
...     name="sphere",
...     long_name="Unit sphere",
...     semimajor_axis=1,
...     flattening=0,
...     geocentric_grav_const=1,
...     angular_velocity=0
... )
>>> print(sphere) 
Ellipsoid(name='sphere', ...)
>>> print(sphere.long_name)
Unit sphere
>>> print("{:.2f}".format(sphere.semiminor_axis))
1.00
>>> print("{:.2f}".format(sphere.mean_radius))
1.00
>>> print("{:.2f}".format(sphere.linear_eccentricity))
0.00
>>> print("{:.2f}".format(sphere.first_eccentricity))
0.00
>>> print("{:.2f}".format(sphere.second_eccentricity))
0.00

Attributes

Ellipsoid.emm

Auxiliary quantity \(m = \omega^2 a^2 b / (GM)\)

Ellipsoid.first_eccentricity

The first eccentricity [adimensional]

Ellipsoid.gravity_equator

The norm of the gravity vector at the equator on the ellipsoid [m/s^2]

Ellipsoid.gravity_pole

The norm of the gravity vector at the poles on the ellipsoid [m/s^2]

Ellipsoid.linear_eccentricity

The linear eccentricity [meters]

Ellipsoid.mean_radius

The arithmetic mean radius [meters]

Ellipsoid.second_eccentricity

The second eccentricity [adimensional]

Ellipsoid.semiminor_axis

The small (polar) axis of the ellipsoid [meters]

Methods

Ellipsoid.geocentric_radius(latitude[, geodetic])

Distance from the center of the ellipsoid to its surface.

Ellipsoid.geodetic_to_spherical(longitude, …)

Convert from geodetic to geocentric spherical coordinates.

Ellipsoid.normal_gravity(latitude, height)

Calculate normal gravity at any latitude and height.

Ellipsoid.prime_vertical_radius(sinlat)

Calculate the prime vertical radius for a given geodetic latitude

Ellipsoid.spherical_to_geodetic(longitude, …)

Convert from geocentric spherical to geodetic coordinates.


Ellipsoid.geocentric_radius(latitude, geodetic=True)[source]

Distance from the center of the ellipsoid to its surface.

The geocentric radius and is a function of the geodetic latitude \(\phi\) and the semi-major and semi-minor axis, a and b:

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

See https://en.wikipedia.org/wiki/Earth_radius#Geocentric_radius

The same could be achieved with boule.Ellipsoid.geodetic_to_spherical by passing any value for the longitudes and heights equal to zero. This method provides a simpler and possibly faster alternative.

Alternatively, the geocentric radius can also be expressed in terms of the geocentric (spherical) latitude \(\theta\):

\[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 latitudes and need the geocentric radius of the ellipsoid (for example, in spherical harmonic analysis). In these cases, the coordinate conversion route is not possible since we need a radius to do that in the first place.

Boule generally tries to avoid doing coordinate conversions inside functions in favour of the user handling the conversions prior to input. This simplifies the code and makes sure that users know precisely which conversions are taking place. This method is an exception since a coordinate conversion route would not be possible.

Note

No elevation is taken into account (the height is zero). If you need the geocentric radius at a height other than zero, use boule.Ellipsoid.geodetic_to_spherical instead.

Parameters
  • latitude (float or array) – Latitude coordinates on geodetic coordinate system in degrees.

  • geodetic (bool) – If True (default), will assume that latitudes are geodetic latitudes. Otherwise, will that they are geocentric spherical latitudes.

Returns

geocentric_radius (float or array) – The geocentric radius for the given latitude(s) in the same units as the ellipsoid axis.

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)[source]

Calculate normal gravity at any latitude and height.

Computes the magnitude of the gradient of the gravity potential (gravitational + centrifugal) generated by the ellipsoid at the given latitude and (geometric) height. Uses of a closed form expression of [LiGotze2001].

Parameters
  • latitude (float or array) – The (geodetic) latitude where the normal gravity will be computed (in degrees).

  • height (float or array) – The ellipsoidal (geometric) height of computation point (in meters).

Returns

gamma (float or array) – The normal gravity in mGal.

Ellipsoid.prime_vertical_radius(sinlat)[source]

Calculate the prime vertical radius for a given geodetic latitude

The prime vertical radius is defined as:

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

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

This function receives the sine of the latitude as input to avoid repeated computations of trigonometric functions.

Parameters

sinlat (float or array-like) – Sine of the latitude angle.

Returns

prime_vertical_radius (float or array-like) – Prime vertical radius given in the same units as the semimajor axis

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.

Examples using boule.Ellipsoid