spatial

  • Spatial transformation routines

  • Gravitational and ellipsoidal parameters [23, 42]

Source code

General Methods

pyTMD.spatial.data_type(x: ndarray, y: ndarray, t: ndarray) str[source]

Determines input data type based on variable dimensions

Parameters:
x: np.ndarray

x-dimension coordinates

y: np.ndarray

y-dimension coordinates

t: np.ndarray

time-dimension coordinates

Returns:
string denoting input data type
  • 'time series'

  • 'drift'

  • 'grid'

class pyTMD.spatial.datum(**kwargs)[source]

Class for gravitational and ellipsoidal parameters [23, 55]

Parameters:
ellipsoid: str, default ‘WGS84’

Reference ellipsoid name

  • 'CLK66': Clarke 1866

  • 'CLK80': Clarke 1880

  • 'GRS67': Geodetic Reference System 1967

  • 'GRS80': Geodetic Reference System 1980

  • 'HGH80': Hughes 1980 Ellipsoid

  • 'WGS60': World Geodetic System 1960

  • 'WGS66': World Geodetic System 1966

  • 'WGS72': World Geodetic System 1972

  • 'WGS84': World Geodetic System 1984

  • 'ATS77': Quasi-earth centred ellipsoid for ATS77

  • 'NAD27': North American Datum 1927

  • 'NAD83': North American Datum 1983

  • 'INTER': International

  • 'KRASS': Krassovsky (USSR)

  • 'SGS90': Soviet Geodetic System 1990

  • 'HLMRT': Helmert 1906 Ellipsoid

  • 'HOUGH': Hough 1960 Ellipsoid

  • 'AIRY': Airy (1830)

  • 'MAIRY': Modified Airy (1849)

  • 'MERIT': MERIT 1983 ellipsoid

  • 'TOPEX': TOPEX/POSEIDON ellipsoid

  • 'EGM96': EGM 1996 gravity model

  • 'IAG75': International Association of Geodesy (1975)

  • 'IAU64': International Astronomical Union (1964)

  • 'IAU76': International Astronomical Union (1976)

  • 'IERS89': IERS Numerical Standards (1989)

  • 'IERS': IERS Numerical Standards (2010)

units: str, default `MKS`

Output units

  • 'MKS': meters, kilograms, seconds

  • 'CGS': centimeters, grams, seconds

Attributes:
a_axis: float

Semi-major axis of the ellipsoid

flat: float

Flattening of the ellipsoid

omega: float

Angular velocity of the Earth

GM: float

Geocentric gravitational constant

property rad_e: float

Average radius of the Earth with same volume as ellipsoid

property b_axis: float

Semi-minor axis of the ellipsoid

property ratio: float

Ratio between ellipsoidal axes

property rad_p: float

Polar radius of curvature

property ecc: float

Linear eccentricity

property ecc1: float

First numerical eccentricity

property ecc2: float

Second numerical eccentricity

property m: float

m Parameter

property f2: float

f2 component

property f4: float

f4 component

property q: float

q Parameter

property q0: float

q0 Parameter

property J2: float

Oblateness coefficient

property C20: float

Normalized C20 harmonic

property gamma_a: float

Normal gravity at the equator

property gamma_b: float

Normal gravity at the pole

gamma_0(theta) float[source]

Normal gravity at colatitudes

Parameters:
theta: float

Colatitudes in radians

gamma_h(theta, height) float[source]

Normal gravity at colatitudes and heights

Parameters:
theta: float

Colatitudes in radians

height: float

Height above ellipsoid

property dk: float

Ratio between gravity at pole versus gravity at equator

property U0: float

Normal potential at the ellipsoid

property area: float

Surface area of the ellipsoid

property volume: float

Volume of the ellipsoid

property rho_e: float

Average density

pyTMD.spatial.convert_ellipsoid(lat1: ndarray, h1: ndarray, a1: float, f1: float, a2: float, f2: float, eps: float = 1e-12, itmax: int = 10)[source]

Convert latitudes and heights to a different ellipsoid using Newton-Raphson [33]

Parameters:
lat1: np.ndarray

latitude of input ellipsoid in degrees

h1: np.ndarray

height above input ellipsoid in meters

a1: float

semi-major axis of input ellipsoid

f1: float

flattening of input ellipsoid

a2: float

semi-major axis of output ellipsoid

f2: float

flattening of output ellipsoid

eps: float, default 1e-12

tolerance to prevent division by small numbers and to determine convergence

itmax: int, default 10

maximum number of iterations to use in Newton-Raphson

Returns:
lat2: np.ndarray

latitude of output ellipsoid in degrees

h2: np.ndarray

height above output ellipsoid in meters

pyTMD.spatial.compute_delta_h(lat: ndarray, a1: float, f1: float, a2: float, f2: float)[source]

Compute difference in elevation for two ellipsoids at a given latitude using a simplified empirical relation [33]

Parameters:
lat: np.ndarray

latitudes (degrees north)

a1: float

semi-major axis of input ellipsoid

f1: float

flattening of input ellipsoid

a2: float

semi-major axis of output ellipsoid

f2: float

flattening of output ellipsoid

Returns:
delta_h: np.ndarray

difference in elevation for two ellipsoids

pyTMD.spatial.wrap_longitudes(lon: float | ndarray)[source]

Wraps longitudes to range from -180 to +180

Parameters:
lon: float or np.ndarray

longitude (degrees east)

pyTMD.spatial.to_dms(d: ndarray)[source]

Convert decimal degrees to degrees, minutes and seconds

Parameters:
d: np.ndarray

decimal degrees

Returns:
degree: np.ndarray

degrees

minute: np.ndarray

minutes (arcminutes)

second: np.ndarray

seconds (arcseconds)

pyTMD.spatial.from_dms(degree: ndarray, minute: ndarray, second: ndarray)[source]

Convert degrees, minutes and seconds to decimal degrees

Parameters:
degree: np.ndarray

degrees

minute: np.ndarray

minutes (arcminutes)

second: np.ndarray

seconds (arcseconds)

Returns:
d: np.ndarray

decimal degrees

pyTMD.spatial.to_cartesian(lon: ndarray, lat: ndarray, h: float | ndarray = 0.0, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805)[source]

Converts geodetic coordinates to Cartesian coordinates

Parameters:
lon: np.ndarray

longitude (degrees east)

lat: np.ndarray

latitude (degrees north)

h: float or np.ndarray, default 0.0

height above ellipsoid (or sphere)

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

for spherical coordinates set to radius of the Earth

flat: float, default 1.0/298.257223563

ellipsoidal flattening

for spherical coordinates set to 0

pyTMD.spatial.to_sphere(x: ndarray, y: ndarray, z: ndarray)[source]

Convert from cartesian coordinates to spherical coordinates

Parameters:
x, np.ndarray

cartesian x-coordinates

y, np.ndarray

cartesian y-coordinates

z, np.ndarray

cartesian z-coordinates

pyTMD.spatial.to_geodetic(x: ndarray, y: ndarray, z: ndarray, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805, method: str = 'bowring', eps: float = np.float64(2.220446049250313e-16), iterations: int = 10)[source]

Convert from cartesian coordinates to geodetic coordinates using either iterative or closed-form methods

Parameters:
x, np.ndarray

cartesian x-coordinates

y, np.ndarray

cartesian y-coordinates

z, np.ndarray

cartesian z-coordinates

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

flat: float, default 1.0/298.257223563

ellipsoidal flattening

method: str, default ‘bowring’

method to use for conversion

  • 'moritz': iterative solution

  • 'bowring': iterative solution

  • 'zhu': closed-form solution

eps: float, default np.finfo(np.float64).eps

tolerance for iterative methods

iterations: int, default 10

maximum number of iterations

pyTMD.spatial._moritz_iterative(x: ndarray, y: ndarray, z: ndarray, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805, eps: float = np.float64(2.220446049250313e-16), iterations: int = 10)[source]

Convert from cartesian coordinates to geodetic coordinates using the iterative solution of Hofmann-Wellenhof and Moritz [23]

Parameters:
x, np.ndarray

cartesian x-coordinates

y, np.ndarray

cartesian y-coordinates

z, np.ndarray

cartesian z-coordinates

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

flat: float, default 1.0/298.257223563

ellipsoidal flattening

eps: float, default np.finfo(np.float64).eps

tolerance for iterative method

iterations: int, default 10

maximum number of iterations

pyTMD.spatial._bowring_iterative(x: ndarray, y: ndarray, z: ndarray, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805, eps: float = np.float64(2.220446049250313e-16), iterations: int = 10)[source]

Convert from cartesian coordinates to geodetic coordinates using the iterative solution of Bowring [4], Bowring [5]

Parameters:
x, np.ndarray

cartesian x-coordinates

y, np.ndarray

cartesian y-coordinates

z, np.ndarray

cartesian z-coordinates

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

flat: float, default 1.0/298.257223563

ellipsoidal flattening

eps: float, default np.finfo(np.float64).eps

tolerance for iterative method

iterations: int, default 10

maximum number of iterations

pyTMD.spatial._zhu_closed_form(x: ndarray, y: ndarray, z: ndarray, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805)[source]

Convert from cartesian coordinates to geodetic coordinates using the closed-form solution of Zhu [61]

Parameters:
x, np.ndarray

cartesian x-coordinates

y, np.ndarray

cartesian y-coordinates

z, np.ndarray

cartesian z-coordinates

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

flat: float, default 1.0/298.257223563

ellipsoidal flattening

pyTMD.spatial.to_ENU(x: ndarray, y: ndarray, z: ndarray, lon0: float | ndarray = 0.0, lat0: float | ndarray = 0.0, h0: float | ndarray = 0.0, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805)[source]

Convert from Earth-Centered Earth-Fixed (ECEF) cartesian coordinates to East-North-Up coordinates (ENU)

Parameters:
x, np.ndarray

cartesian x-coordinates

y, np.ndarray

cartesian y-coordinates

z, np.ndarray

cartesian z-coordinates

lon0: float or np.ndarray, default 0.0

reference longitude (degrees east)

lat0: float or np.ndarray, default 0.0

reference latitude (degrees north)

h0: float or np.ndarray, default 0.0

reference height (meters)

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

flat: float, default 1.0/298.257223563

ellipsoidal flattening

Returns:
E: np.ndarray

east coordinates

N: np.ndarray

north coordinates

U: np.ndarray

up coordinates

pyTMD.spatial.from_ENU(E: ndarray, N: ndarray, U: ndarray, lon0: float | ndarray = 0.0, lat0: float | ndarray = 0.0, h0: float | ndarray = 0.0, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805)[source]

Convert from East-North-Up coordinates (ENU) to Earth-Centered Earth-Fixed (ECEF) cartesian coordinates

Parameters:
E, np.ndarray

east coordinates

N, np.ndarray

north coordinates

U, np.ndarray

up coordinates

lon0: float or np.ndarray, default 0.0

reference longitude (degrees east)

lat0: float or np.ndarray, default 0.0

reference latitude (degrees north)

h0: float or np.ndarray, default 0.0

reference height (meters)

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

flat: float, default 1.0/298.257223563

ellipsoidal flattening

Returns:
x, float

cartesian x-coordinates

y, float

cartesian y-coordinates

z, float

cartesian z-coordinates

pyTMD.spatial.to_horizontal(E: ndarray, N: ndarray, U: ndarray)[source]

Convert from East-North-Up coordinates (ENU) to a celestial horizontal coordinate system (alt-az)

Parameters:
E: np.ndarray

east coordinates

N: np.ndarray

north coordinates

U: np.ndarray

up coordinates

Returns:
alpha: np.ndarray

altitude (elevation) angle in degrees

phi: np.ndarray

azimuth angle in degrees

D: np.ndarray

distance from observer to object in meters

pyTMD.spatial.to_zenith(x: ndarray, y: ndarray, z: ndarray, lon0: float | ndarray = 0.0, lat0: float | ndarray = 0.0, h0: float | ndarray = 0.0, a_axis: float = 6378137.0, flat: float = 0.0033528106647474805)[source]

Calculate zenith angle of an object from Earth-Centered Earth-Fixed (ECEF) cartesian coordinates

Parameters:
x, np.ndarray

cartesian x-coordinates

y, np.ndarray

cartesian y-coordinates

z, np.ndarray

cartesian z-coordinates

lon0: float or np.ndarray, default 0.0

reference longitude (degrees east)

lat0: float or np.ndarray, default 0.0

reference latitude (degrees north)

h0: float or np.ndarray, default 0.0

reference height (meters)

a_axis: float, default 6378137.0

semimajor axis of the ellipsoid

flat: float, default 1.0/298.257223563

ellipsoidal flattening

Returns:
zenith: np.ndarray

zenith angle of object in degrees

pyTMD.spatial.scale_factors(lat: ndarray, flat: float = 0.0033528106647474805, reference_latitude: float = 70.0, metric: str = 'area')[source]

Calculates scaling factors to account for polar stereographic distortion including special case of at the exact pole [52]

Parameters:
lat: np.ndarray

latitude (degrees north)

flat: float, default 1.0/298.257223563

ellipsoidal flattening

reference_latitude: float, default 70.0

reference latitude (true scale latitude)

metric: str, default ‘area’

metric to calculate scaling factors

  • 'distance': scale factors for distance

  • 'area': scale factors for area

Returns:
scale: np.ndarray

scaling factors at input latitudes