boule.Sphere

class boule.Sphere(name, radius, geocentric_grav_const, angular_velocity, long_name=None, reference=None)[source]

Reference sphere (zero flattening ellipsoids)

Represents a rotating reference ellipsoid with zero flattening. It is defined by three parameters (radius, geocentric gravitational constant, and angular velocity) and offers other derived quantities.

All attributes of this class are read-only and cannot be changed after instantiation.

All parameters are in SI units.

Note

Must be used instead of boule.Ellipsoid to account for singularities due to zero flattening (and thus zero eccentricity) in normal gravity calculations.

Parameters
  • name (str) – A short name for the sphere, for example 'Moon'.

  • radius (float) – The radius of the sphere [meters].

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

  • angular_velocity (float) – The angular velocity of the rotating sphere (omega) [rad s^-1].

  • long_name (str or None) – A long name for the sphere, for example "Moon Reference System" (optional).

  • reference (str or None) – Citation for the sphere parameter values (optional).

Examples

We can define a sphere by specifying the 3 key numerical parameters:

>>> sphere = Sphere(
...     name="Moon",
...     long_name="That's no moon",
...     radius=1,
...     geocentric_grav_const=2,
...     angular_velocity=0.5,
... )
>>> print(sphere) 
Sphere(name='Moon', ...)
>>> print(sphere.long_name)
That's no moon

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

>>> print("{:.2f}".format(sphere.semimajor_axis))
1.00
>>> print("{:.2f}".format(sphere.semiminor_axis))
1.00
>>> print("{:.2f}".format(sphere.mean_radius))
1.00
>>> print("{:.2f}".format(sphere.gravity_equator))
1.75
>>> print("{:.2f}".format(sphere.gravity_pole))
2.00

Normal gravity (the magnitude of the gravity potential) can be calculated at any latitude and height. Note that this method returns values in mGal instead of m/s².

>>> print("{:.2f}".format(sphere.normal_gravity(latitude=0, height=0)))
175000.00
>>> print("{:.2f}".format(sphere.normal_gravity(latitude=90, height=0)))
200000.00

The flag si_units will return the Normal gravity in m/s².

>>> print("{:.2f}".format(sphere.normal_gravity(latitude=0, height=0, si_units=True)))
1.75
>>> print("{:.2f}".format(sphere.normal_gravity(latitude=90, height=0, si_units=True)))
2.00

The flattening and eccentricities will all be zero:

>>> print("{:.2f}".format(sphere.flattening))
0.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

Methods Summary

Sphere.geocentric_radius(latitude[, geodetic])

Distance from the center of the ellipsoid to its surface.

Sphere.geodetic_to_spherical(longitude, …)

Convert from geodetic to geocentric spherical coordinates.

Sphere.normal_gravity(latitude, height[, …])

Calculate normal gravity at any latitude and height

Sphere.prime_vertical_radius(sinlat)

Calculate the prime vertical radius for a given geodetic latitude

Sphere.spherical_to_geodetic(longitude, …)

Convert from geocentric spherical to geodetic coordinates.


Sphere.geocentric_radius(latitude, geodetic=True)

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 the radial coordinates to do that in the first place.

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.

Sphere.geodetic_to_spherical(longitude, latitude, height)

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.

Sphere.normal_gravity(latitude, height, si_units=False)[source]

Calculate normal gravity at any latitude and height

Computes the magnitude of the gradient of the gravity potential (gravitational + centrifugal; see [Heiskanen-Moritz]) generated by the sphere at the given latitude \(\theta\) and height \(h\):

\[\gamma(\theta, h) = \sqrt{\left( \frac{GM}{(R + h)^2} \right)^2 + \left(\omega^2 (R + h) - 2\frac{GM}{(R + h)^2} \right) \omega^2 (R + h) \cos^2 \theta}\]

in which \(R\) is the sphere radius, \(G\) is the gravitational constant, \(M\) is the mass of the sphere, and \(\omega\) is the angular velocity.

Note

A sphere under rotation is not in hydrostatic equilibrium. Therefore, it is not it’s own equipotential gravity surface (as is the case for the ellipsoid) and the normal gravity vector is not normal to the surface of the sphere.

Parameters
  • latitude (float or array) – The latitude where the normal gravity will be computed (in degrees). For a reference sphere there is no difference between geodetic and spherical latitudes.

  • height (float or array) – The height (above the surface of the sphere) of the computation point (in meters).

  • si_units (bool) – Return the value in mGal (False, default) or SI units (True)

Returns

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

Sphere.prime_vertical_radius(sinlat)

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 semi-major 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 semi-major axis

Sphere.spherical_to_geodetic(longitude, spherical_latitude, radius)

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.Sphere