Coverage for pygeodesy/ellipsoids.py: 96%
766 statements
« prev ^ index » next coverage.py v7.10.7, created at 2026-03-25 15:01 -0400
« prev ^ index » next coverage.py v7.10.7, created at 2026-03-25 15:01 -0400
2# -*- coding: utf-8 -*-
4u'''Ellipsoidal and spherical earth models.
6Classes L{a_f2Tuple}, L{Ellipsoid} and L{Ellipsoid2}, an L{Ellipsoids} registry and
72 dozen functions to convert I{equatorial} radius, I{polar} radius, I{eccentricities},
8I{flattenings} and I{inverse flattening}.
10See module L{datums} for L{Datum} and L{Transform} information and other details.
12Following is the list of predefined L{Ellipsoid}s, all instantiated lazily.
14@var Ellipsoids.Airy1830: Ellipsoid(name='Airy1830', a=6377563.396, f=0.00334085, f_=299.3249646, b=6356256.90923729)
15@var Ellipsoids.AiryModified: Ellipsoid(name='AiryModified', a=6377340.189, f=0.00334085, f_=299.3249646, b=6356034.44793853)
16@var Ellipsoids.ATS1977: Ellipsoid(name='ATS1977', a=6378135, f=0.00335281, f_=298.257, b=6356750.30492159)
17@var Ellipsoids.Australia1966: Ellipsoid(name='Australia1966', a=6378160, f=0.00335289, f_=298.25, b=6356774.71919531)
18@var Ellipsoids.Bessel1841: Ellipsoid(name='Bessel1841', a=6377397.155, f=0.00334277, f_=299.1528128, b=6356078.962818)
19@var Ellipsoids.BesselModified: Ellipsoid(name='BesselModified', a=6377492.018, f=0.00334277, f_=299.1528128, b=6356173.5087127)
20@var Ellipsoids.CGCS2000: Ellipsoid(name='CGCS2000', a=6378137, f=0.00335281, f_=298.2572221, b=6356752.31414036)
21@var Ellipsoids.Clarke1866: Ellipsoid(name='Clarke1866', a=6378206.4, f=0.00339008, f_=294.97869821, b=6356583.8)
22@var Ellipsoids.Clarke1880: Ellipsoid(name='Clarke1880', a=6378249.145, f=0.00340756, f_=293.465, b=6356514.86954978)
23@var Ellipsoids.Clarke1880IGN: Ellipsoid(name='Clarke1880IGN', a=6378249.2, f=0.00340755, f_=293.46602129, b=6356515)
24@var Ellipsoids.Clarke1880Mod: Ellipsoid(name='Clarke1880Mod', a=6378249.145, f=0.00340755, f_=293.46630766, b=6356514.96639549)
25@var Ellipsoids.CPM1799: Ellipsoid(name='CPM1799', a=6375738.7, f=0.00299052, f_=334.39, b=6356671.92557493)
26@var Ellipsoids.Delambre1810: Ellipsoid(name='Delambre1810', a=6376428, f=0.00321027, f_=311.5, b=6355957.92616372)
27@var Ellipsoids.Engelis1985: Ellipsoid(name='Engelis1985', a=6378136.05, f=0.00335282, f_=298.2566, b=6356751.32272154)
28@var Ellipsoids.Everest1969: Ellipsoid(name='Everest1969', a=6377295.664, f=0.00332445, f_=300.8017, b=6356094.667915)
29@var Ellipsoids.Everest1975: Ellipsoid(name='Everest1975', a=6377299.151, f=0.00332445, f_=300.8017255, b=6356098.14512013)
30@var Ellipsoids.Fisher1968: Ellipsoid(name='Fisher1968', a=6378150, f=0.00335233, f_=298.3, b=6356768.33724438)
31@var Ellipsoids.GEM10C: Ellipsoid(name='GEM10C', a=6378137, f=0.00335281, f_=298.2572236, b=6356752.31424783)
32@var Ellipsoids.GPES: Ellipsoid(name='GPES', a=6378135, f=0, f_=0, b=6378135)
33@var Ellipsoids.GRS67: Ellipsoid(name='GRS67', a=6378160, f=0.00335292, f_=298.24716743, b=6356774.51609071)
34@var Ellipsoids.GRS80: Ellipsoid(name='GRS80', a=6378137, f=0.00335281, f_=298.2572221, b=6356752.31414035)
35@var Ellipsoids.Helmert1906: Ellipsoid(name='Helmert1906', a=6378200, f=0.00335233, f_=298.3, b=6356818.16962789)
36@var Ellipsoids.IAU76: Ellipsoid(name='IAU76', a=6378140, f=0.00335281, f_=298.257, b=6356755.28815753)
37@var Ellipsoids.IERS1989: Ellipsoid(name='IERS1989', a=6378136, f=0.00335281, f_=298.257, b=6356751.30156878)
38@var Ellipsoids.IERS1992TOPEX: Ellipsoid(name='IERS1992TOPEX', a=6378136.3, f=0.00335281, f_=298.25722356, b=6356751.61659215)
39@var Ellipsoids.IERS2003: Ellipsoid(name='IERS2003', a=6378136.6, f=0.00335282, f_=298.25642, b=6356751.85797165)
40@var Ellipsoids.Intl1924: Ellipsoid(name='Intl1924', a=6378388, f=0.003367, f_=297, b=6356911.94612795)
41@var Ellipsoids.Intl1967: Ellipsoid(name='Intl1967', a=6378157.5, f=0.0033529, f_=298.24961539, b=6356772.2)
42@var Ellipsoids.Krassovski1940: Ellipsoid(name='Krassovski1940', a=6378245, f=0.00335233, f_=298.3, b=6356863.01877305)
43@var Ellipsoids.Krassowsky1940: Ellipsoid(name='Krassowsky1940', a=6378245, f=0.00335233, f_=298.3, b=6356863.01877305)
44@var Ellipsoids.Maupertuis1738: Ellipsoid(name='Maupertuis1738', a=6397300, f=0.0052356, f_=191, b=6363806.28272251)
45@var Ellipsoids.Mercury1960: Ellipsoid(name='Mercury1960', a=6378166, f=0.00335233, f_=298.3, b=6356784.28360711)
46@var Ellipsoids.Mercury1968Mod: Ellipsoid(name='Mercury1968Mod', a=6378150, f=0.00335233, f_=298.3, b=6356768.33724438)
47@var Ellipsoids.NWL1965: Ellipsoid(name='NWL1965', a=6378145, f=0.00335289, f_=298.25, b=6356759.76948868)
48@var Ellipsoids.OSU86F: Ellipsoid(name='OSU86F', a=6378136.2, f=0.00335281, f_=298.2572236, b=6356751.51693008)
49@var Ellipsoids.OSU91A: Ellipsoid(name='OSU91A', a=6378136.3, f=0.00335281, f_=298.2572236, b=6356751.6165948)
50@var Ellipsoids.Plessis1817: Ellipsoid(name='Plessis1817', a=6376523, f=0.00324002, f_=308.64, b=6355862.93325557)
51@var Ellipsoids.PZ90: Ellipsoid(name='PZ90', a=6378136, f=0.0033528, f_=298.2578393, b=6356751.36174571)
52@var Ellipsoids.SGS85: Ellipsoid(name='SGS85', a=6378136, f=0.00335281, f_=298.257, b=6356751.30156878)
53@var Ellipsoids.SoAmerican1969: Ellipsoid(name='SoAmerican1969', a=6378160, f=0.00335289, f_=298.25, b=6356774.71919531)
54@var Ellipsoids.Sphere: Ellipsoid(name='Sphere', a=6371008.771415, f=0, f_=0, b=6371008.771415)
55@var Ellipsoids.SphereAuthalic: Ellipsoid(name='SphereAuthalic', a=6371000, f=0, f_=0, b=6371000)
56@var Ellipsoids.SpherePopular: Ellipsoid(name='SpherePopular', a=6378137, f=0, f_=0, b=6378137)
57@var Ellipsoids.Struve1860: Ellipsoid(name='Struve1860', a=6378298.3, f=0.00339294, f_=294.73, b=6356657.14266956)
58@var Ellipsoids.WGS60: Ellipsoid(name='WGS60', a=6378165, f=0.00335233, f_=298.3, b=6356783.28695944)
59@var Ellipsoids.WGS66: Ellipsoid(name='WGS66', a=6378145, f=0.00335289, f_=298.25, b=6356759.76948868)
60@var Ellipsoids.WGS72: Ellipsoid(name='WGS72', a=6378135, f=0.00335278, f_=298.26, b=6356750.52001609)
61@var Ellipsoids.WGS84: Ellipsoid(name='WGS84', a=6378137, f=0.00335281, f_=298.25722356, b=6356752.31424518)
62@var Ellipsoids.WGS84_NGS: Ellipsoid(name='WGS84_NGS', a=6378137, f=0.00335281, f_=298.2572221, b=6356752.31414035)
63'''
64# make sure int/int division yields float quotient, see .basics
65from __future__ import division as _; del _ # noqa: E702 ;
67# from pygeodesy.albers import AlbersEqualAreaCylindrical # _MODS
68from pygeodesy.basics import copysign0, isbool, _isin, isint, typename
69from pygeodesy.constants import EPS, EPS_2, EPS0, EPS02, EPS1, INF, NINF, \
70 _EPSqrt, PI_2, PI_3, PI2, PI4, R_M, R_MA, R_FM, \
71 _EPStol as _TOL, _floatuple as _T, _isfinite, \
72 _over, _SQRT2, _0_0s, _0_0, _0_5, _1_0, _1_EPS, \
73 _2_0, _4_0, _90_0, _0_25, _3_0 # PYCHOK used!
74# from pygeodesy.ellipses import Ellipse # _MODS
75from pygeodesy.errors import _AssertionError, IntersectionError, _ValueError, _xattr, _xkwds_not
76from pygeodesy.fmath import cbrt, cbrt2, fdot, fmean_, Fhorner, fpowers, hypot, hypot_, \
77 hypot1, hypot2, sqrt3, Fsum
78# from pygeodesy.fsums import Fsum # from .fmath
79# from pygeodesy.internals import typename # from .basics
80from pygeodesy.interns import NN, _a_, _Airy1830_, _AiryModified_, _b_, _Bessel1841_, \
81 _Clarke1866_, _Clarke1880IGN_, _DMAIN_, _DOT_, _f_, _GRS80_, \
82 _Intl1924_, _incompatible_, _invalid_, _Krassovski1940_, \
83 _Krassowsky1940_, _lat_, _meridional_, _negative_, _not_finite_, \
84 _prime_vertical_, _Sphere_, _SPACE_, _vs_, _WGS72_, _WGS84_
85# from pygeodesy.lazily import _ALL_LAZY, _ALL_MODS as _MODS # from .named
86from pygeodesy.named import _lazyNamedEnumItem as _lazy, _name__, _NamedEnum, \
87 _NamedEnumItem, _NamedTuple, _Pass, _ALL_LAZY, _MODS
88from pygeodesy.namedTuples import Circle4Tuple, Distance2Tuple, Vector3Tuple, Vector4Tuple
89from pygeodesy.props import deprecated_Property_RO, Property_RO, property_doc_, \
90 deprecated_property_RO, property_RO, property_ROver
91from pygeodesy.streprs import Fmt, fstr, instr, strs, unstr
92# from pygeodesy.triaxials.triaxial5 import _hartzell3, _plumbTo3 # _MODS
93from pygeodesy.units import Azimuth, Bearing, Distance, Float, Float_, Lamd, Lat, _Lat0, \
94 Meter, Meter2, Meter3, Phi, Phid, Radius, Radius_, Scalar
95from pygeodesy.utily import atan1, atan1d, atan2b, degrees90, m2radians, radians2m, sincos2d
97from math import asinh, atan, atanh, cos, degrees, exp, fabs, radians, sin, sinh, sqrt, tan # as _tan
99__all__ = _ALL_LAZY.ellipsoids
100__version__ = '26.03.25'
102_f_0_0 = Float(f =_0_0) # zero flattening
103_f__0_0 = Float(f_=_0_0) # zero inverse flattening
104# see U{WGS84_f<https://GeographicLib.SourceForge.io/C++/doc/classGeographicLib_1_1Constants.html>}
105# _f_WGS84 = Float(f =_1_0 / (298257223563 / 1000000000)) # 0.003_352_810_664_747_480_5
106_f__WGS84 = Float(f_=_1_0 / (1000000000 / 298257223563)) # 298.257223562_999_97 vs 298.257223563
109def _aux(lat, inverse, auxLat, clip=90):
110 '''Return a named auxiliary latitude in C{degrees}.
111 '''
112 return Lat(lat, clip=clip, name=_lat_ if inverse else typename(auxLat))
115def _sin2cos2(rad):
116 '''(INTERNAL) Return 2-tuple C{(sin(B{phi})**2, cos(B{phi})**2)}.
117 '''
118 if rad:
119 s2 = sin(rad)**2
120 if s2 > EPS:
121 c2 = _1_0 - s2
122 if c2 > EPS:
123 if c2 < EPS1:
124 return s2, c2
125 else:
126 return _1_0, _0_0 # rad == PI_2
127 return _0_0, _1_0 # rad == 0
130class a_f2Tuple(_NamedTuple):
131 '''2-Tuple C{(a, f)} specifying an ellipsoid by I{equatorial}
132 radius C{a} in C{meter} and scalar I{flattening} C{f}.
134 @see: Class L{Ellipsoid2}.
135 '''
136 _Names_ = (_a_, _f_) # name 'f' not 'f_'
137 _Units_ = (_Pass, _Pass)
139 def __new__(cls, a, f, **name):
140 '''New L{a_f2Tuple} ellipsoid specification.
142 @arg a: Equatorial radius (C{scalar} > 0).
143 @arg f: Flattening (C{scalar} < 1, negative for I{prolate}).
144 @kwarg name: Optional C{B{name}=NN} (C{str}).
146 @return: An L{a_f2Tuple}C{(a, f)}.
148 @raise UnitError: Invalid B{C{a}} or B{C{f}}.
150 @note: C{abs(B{f}) < }L{EPS<pygeodesy.constants.EPS>} is
151 forced to C{B{f}=0}, I{spherical}.
153 @note: Negative C{B{f}} produces a I{prolate} ellipsoid.
154 '''
155 a = Radius_(a=a) # low=EPS, high=None
156 f = Float_( f=f, low=None, high=EPS1)
157 if fabs(f) < EPS: # force spherical
158 f = _f_0_0
159 return _NamedTuple.__new__(cls, a, f, **name)
161 @Property_RO
162 def b(self):
163 '''Get the I{polar} radius (C{meter}), M{a * (1 - f)}.
164 '''
165 return a_f2b(self.a, self.f) # PYCHOK .a and .f
167 def ellipsoid(self, **name):
168 '''Return an L{Ellipsoid} for this 2-tuple C{(a, f)}.
170 @kwarg name: Optional C{B{name}=NN} (C{str}).
172 @raise NameError: A registered C{ellipsoid} with the
173 same B{C{name}} already exists.
174 '''
175 return Ellipsoid(self.a, f=self.f, name=self._name__(name)) # PYCHOK .a and .f
177 @Property_RO
178 def f_(self):
179 '''Get the I{inverse} flattening (C{scalar}), M{1 / f} == M{a / (a - b)}.
180 '''
181 return f2f_(self.f) # PYCHOK .f
184class Curvature2Tuple(_NamedTuple):
185 '''2-Tuple C{(meridional, prime_vertical)} of radii of curvature, both in
186 C{meter}, conventionally.
187 '''
188 _Names_ = (_meridional_, _prime_vertical_)
189 _Units_ = ( Meter, Meter)
191 @property_RO
192 def transverse(self):
193 '''Get this I{prime_vertical}, aka I{transverse} radius of curvature.
194 '''
195 return self.prime_vertical
198class Ellipsoid(_NamedEnumItem):
199 '''Ellipsoid with I{equatorial} and I{polar} radii, I{flattening}, I{inverse
200 flattening} and other, often used, I{cached} attributes, supporting
201 I{oblate} and I{prolate} ellipsoidal and I{spherical} earth models.
202 '''
203 _a = 0 # equatorial radius, semi-axis (C{meter})
204 _b = 0 # polar radius, semi-axis (C{meter}): a * (f - 1) / f
205 _f = 0 # (1st) flattening: (a - b) / a
206 _f_ = 0 # inverse flattening: 1 / f = a / (a - b)
208 _geodsolve = NN # means, use PYGEODESY_GEODSOLVE
209 _KsOrder = 8 # Krüger series order (4, 6 or 8)
210 _rhumbsolve = NN # means, use PYGEODESY_RHUMBSOLVE
212 def __init__(self, a, b=None, f_=None, f=None, **name):
213 '''New L{Ellipsoid} from the I{equatorial} radius I{and} either
214 the I{polar} radius or I{inverse flattening} or I{flattening}.
216 @arg a: Equatorial radius, semi-axis (C{meter}).
217 @arg b: Optional polar radius, semi-axis (C{meter}).
218 @arg f_: Inverse flattening: M{a / (a - b)} (C{float} >>> 1.0).
219 @arg f: Flattening: M{(a - b) / a} (C{scalar}, near zero for
220 spherical).
221 @kwarg name: Optional, unique C{B{name}=NN} (C{str}).
223 @raise NameError: Ellipsoid with the same B{C{name}} already exists.
225 @raise ValueError: Invalid B{C{a}}, B{C{b}}, B{C{f_}} or B{C{f}} or
226 B{C{f_}} and B{C{f}} are incompatible.
228 @note: M{abs(f_) > 1 / EPS} or M{abs(1 / f_) < EPS} is forced
229 to M{1 / f_ = 0}, spherical.
230 '''
231 ff_ = f, f_ # assertion below
232 n = _name__(**name) if name else NN
233 try:
234 a = Radius_(a=a) # low=EPS
235 if not _isfinite(a):
236 raise ValueError(_SPACE_(_a_, _not_finite_))
238 if b: # not _isin(b, None, _0_0)
239 b = Radius_(b=b) # low=EPS
240 f = a_b2f(a, b) if f is None else Float(f=f)
241 f_ = f2f_(f) if f_ is None else Float(f_=f_)
242 elif f is not None:
243 f = Float(f=f)
244 b = a_f2b(a, f)
245 f_ = f2f_(f) if f_ is None else Float(f_=f_)
246 elif f_:
247 f_ = Float(f_=f_)
248 b = a_f_2b(a, f_) # a * (f_ - 1) / f_
249 f = f_2f(f_)
250 else: # only a, spherical
251 f_ = f = 0
252 b = a # superfluous
254 if not (f < _1_0): # sanity check, see .ecef.Ecef.__init__
255 raise ValueError(_SPACE_(_f_, _invalid_))
256 if not _isfinite(b):
257 raise ValueError(_SPACE_(_b_, _not_finite_))
259 if fabs(f) < EPS or a == b or not f_: # spherical
260 b = a
261 f = _f_0_0
262 f_ = _f__0_0
264 except (TypeError, ValueError) as x:
265 d = _xkwds_not(None, b=b, f_=f_, f=f)
266 t = instr(self, a=a, name=n, **d)
267 raise _ValueError(t, cause=x)
269 self._a = a
270 self._b = b
271 self._f = f
272 self._f_ = f_
274 self._register(Ellipsoids, n)
276 if f and f_: # see test/testEllipsoidal
277 d = dict(eps=_TOL)
278 if None in ff_: # both f_ and f given
279 d.update(Error=_ValueError, txt=_incompatible_)
280 self._assert(_1_0 / f, f_=f_, **d)
281 self._assert(_1_0 / f_, f =f, **d)
282 self._assert(self.b2_a2, e21=self.e21, eps=EPS)
284 def __eq__(self, other):
285 '''Compare this and an other ellipsoid.
287 @arg other: The other ellipsoid (L{Ellipsoid} or L{Ellipsoid2}).
289 @return: C{True} if equal, C{False} otherwise.
290 '''
291 return self is other or (isinstance(other, Ellipsoid) and
292 self.a == other.a and
293 (self.f == other.f or self.b == other.b))
295 def __hash__(self):
296 return self._hash # memoized
298 @Property_RO
299 def a(self):
300 '''Get the I{equatorial} radius, semi-axis (C{meter}).
301 '''
302 return self._a
304 equatoradius = a # = Requatorial
306 @Property_RO
307 def a2(self):
308 '''Get the I{equatorial} radius I{squared} (C{meter} I{squared}), M{a**2}.
309 '''
310 return Meter2(a2=self.a**2)
312 @Property_RO
313 def a2_(self):
314 '''Get the inverse of the I{equatorial} radius I{squared} (C{meter} I{squared}), M{1 / a**2}.
315 '''
316 return Float(a2_=_1_0 / self.a2)
318 @Property_RO
319 def A(self):
320 '''Get the UTM I{meridional (or rectifying)} radius (C{meter}).
322 @note: C{A * PI / 2} ≈= L{L<Ellipsoid.L>}, the I{quarter meridian}.
324 @see: I{Meridian arc unit} U{Q<https://StudyLib.net/doc/7443565/>},
325 the mean, meridional length I{per radian}.
326 '''
327 A, n = self.a, self.n
328 if n:
329 d = (n + _1_0) * 1048576 / A
330 if d: # use 6 n**2 terms, half-way between the _KsOrder's 4, 6, 8
331 # <https://GeographicLib.SourceForge.io/C++/doc/tmseries30.html>
332 # <https://GeographicLib.SourceForge.io/C++/doc/transversemercator.html> and
333 # <https://www.MyGeodesy.id.AU/documents/Karney-Krueger%20equations.pdf> (3)
334 # A *= fhorner(n**2, 1048576, 262144, 16384, 4096, 1600, 784, 441) / 1048576) / (1 + n)
335 A = Radius(A=Fhorner(n**2, 1048576, 262144, 16384, 4096, 1600, 784, 441).fover(d))
336 return A
337# # Moritz, H. <https://Geodesy.Geology.Ohio-State.edu/course/refpapers/00740128.pdf>
338# # q = 4 / self.rocPolar
339# # Q = (1 - 3 / 4 * e'2 + 45 / 64 * e'4 - 175 / 256 * e'6 + 11025 / 16384 * e'8) * rocPolar
340# # = (4 + e'2 * (-3 + e'2 * (45 / 16 + e'2 * (-175 / 64 + e'2 * 11025 / 4096)))) / q
341# return Radius(Q=Fhorner(self.e22, 4, -3, 45 / 16, -175 / 64, 11025 / 4096).fover(q))
343 @Property_RO
344 def a_b(self):
345 '''Get the ratio I{equatorial} over I{polar} radius (C{float}), M{a / b} == M{1 / (1 - f)}.
346 '''
347 return Float(a_b=self.a / self.b if self.f else _1_0)
349 @Property_RO
350 def a2_b(self):
351 '''Get the I{polar} meridional (or polar) radius of curvature (C{meter}), M{a**2 / b}.
353 @see: U{Radii of Curvature
354 <https://WikiPedia.org/wiki/Earth_radius#Radii_of_curvature>}
355 and U{Moritz, H. (1980), Geodetic Reference System 1980
356 <https://WikiPedia.org/wiki/Earth_radius#cite_note-Moritz-2>}.
358 @note: Symbol C{c} is used by IUGG and IERS for the U{polar radius of curvature
359 <https://WikiPedia.org/wiki/Earth_radius#Radii_of_curvature>}, see L{c2}
360 and L{R2} or L{Rauthalic}.
361 '''
362 return Radius(a2_b=self.a2 / self.b if self.f else self.a) # = rocPolar
364 @Property_RO
365 def a2_b2(self):
366 '''Get the ratio I{equatorial} over I{polar} radius I{squared} (C{float}),
367 M{(a / b)**2} == M{1 / (1 - e**2)} == M{1 / (1 - e2)} == M{1 / e21}.
368 '''
369 return Float(a2_b2=self.a_b**2 if self.f else _1_0)
371 @Property_RO
372 def a_f(self):
373 '''Get the I{equatorial} radius and I{flattening} (L{a_f2Tuple}), see method C{toEllipsoid2}.
374 '''
375 return a_f2Tuple(self.a, self.f, name=self.name)
377 a_f2 = a_f # synonym
379 @Property_RO
380 def _albersCyl(self):
381 '''(INTERNAL) Helper for C{auxAuthalic}.
382 '''
383 return _MODS.albers.AlbersEqualAreaCylindrical(datum=self, name=self.name)
385 @Property_RO
386 def AlphaKs(self):
387 '''Get the I{Krüger} U{Alpha series coefficients<https://GeographicLib.SourceForge.io/C++/doc/tmseries30.html>} (C{KsOrder}C{-tuple}).
388 '''
389 return self._Kseries( # XXX int/int quotients may require from __future__ import division as _; del _ # noqa: E702 ;
390 # n n**2 n**3 n**4 n**5 n**6 n**7 n**8
391 _T(1/2, -2/3, 5/16, 41/180, -127/288, 7891/37800, 72161/387072, -18975107/50803200),
392 _T(13/48, -3/5, 557/1440, 281/630, -1983433/1935360, 13769/28800, 148003883/174182400), # PYCHOK unaligned
393 _T(61/240, -103/140, 15061/26880, 167603/181440, -67102379/29030400, 79682431/79833600), # PYCHOK unaligned
394 _T(49561/161280, -179/168, 6601661/7257600, 97445/49896, -40176129013/7664025600), # PYCHOK unaligned
395 _T(34729/80640, -3418889/1995840, 14644087/9123840, 2605413599/622702080), # PYCHOK unaligned
396 _T(212378941/319334400, -30705481/10378368, 175214326799/58118860800), # PYCHOK unaligned
397 _T(1522256789/1383782400, -16759934899/3113510400), # PYCHOK unaligned
398 _T(1424729850961/743921418240)) # PYCHOK unaligned
400 @Property_RO
401 def area(self):
402 '''Get the ellipsoid's surface area (C{meter} I{squared}), M{4 * PI * c2}.
404 @see: Properties L{areax}, L{c2} and L{R2} and functions
405 L{ellipsoidalExact.areaOf} and L{ellipsoidalKarney.areaOf}.
406 '''
407 return Meter2(area=self.c2 * PI4)
409 @Property_RO
410 def areax(self):
411 '''Get the ellipsoid's surface area (C{meter} I{squared}), M{4 * PI * c2x}, more
412 accurate for very I{oblate} ellipsoids.
414 @see: Properties L{area}, L{c2x} and L{R2x}, class L{GeodesicExact} and
415 functions L{ellipsoidalExact.areaOf} and L{ellipsoidalKarney.areaOf}.
416 '''
417 return Meter2(areax=self.c2x * PI4)
419 def _assert(self, val, eps=_TOL, f0=_0_0, Error=_AssertionError, txt=NN, **name_value):
420 '''(INTERNAL) Assert a C{name=value} vs C{val}.
421 '''
422 for n, v in name_value.items():
423 if fabs(v - val) > eps: # PYCHOK no cover
424 t = (v, _vs_, val)
425 t = _SPACE_.join(strs(t, prec=12, fmt=Fmt.g))
426 t = Fmt.EQUAL(self._DOT_(n), t)
427 raise Error(t, txt=txt or Fmt.exceeds_eps(eps))
428 return Float(v if self.f else f0, name=n)
429 raise Error(unstr(self._DOT_(typename(self._assert)), val,
430 eps=eps, f0=f0, **name_value))
432 def auxAuthalic(self, lat, inverse=False):
433 '''Compute the I{authalic} auxiliary latitude (Xi) or the I{inverse} thereof.
435 @arg lat: The geodetic (or I{authalic}) latitude (C{degrees90}).
436 @kwarg inverse: If C{True}, B{C{lat}} is the I{authalic} and
437 return the geodetic latitude (C{bool}).
439 @return: The I{authalic} (or geodetic) latitude in C{degrees90}.
441 @see: U{Inverse-/AuthalicLatitude<https://GeographicLib.SourceForge.io/
442 C++/doc/classGeographicLib_1_1Ellipsoid.html>}, U{Authalic latitude
443 <https://WikiPedia.org/wiki/Latitude#Authalic_latitude>}, and
444 U{Snyder<https://Pubs.USGS.gov/pp/1395/report.pdf>}, p 16.
445 '''
446 if self.f:
447 _f = self._albersCyl._tanf if inverse else \
448 self._albersCyl._txif # PYCHOK attr
449 lat = atan1d(_f(tan(Phid(lat)))) # PYCHOK attr
450 return _aux(lat, inverse, Ellipsoid.auxAuthalic)
452 def auxConformal(self, lat, inverse=False):
453 '''Compute the I{conformal} auxiliary latitude (Chi) or the I{inverse} thereof.
455 @arg lat: The geodetic (or I{conformal}) latitude (C{degrees90}).
456 @kwarg inverse: If C{True}, B{C{lat}} is the I{conformal} and
457 return the geodetic latitude (C{bool}).
459 @return: The I{conformal} (or geodetic) latitude in C{degrees90}.
461 @see: U{Inverse-/ConformalLatitude<https://GeographicLib.SourceForge.io/
462 C++/doc/classGeographicLib_1_1Ellipsoid.html>}, U{Conformal latitude
463 <https://WikiPedia.org/wiki/Latitude#Conformal_latitude>}, and
464 U{Snyder<https://Pubs.USGS.gov/pp/1395/report.pdf>}, pp 15-16.
465 '''
466 if self.f:
467 _f = self.es_tauf if inverse else self.es_taupf # PYCHOK attr
468 lat = atan1d(_f(tan(Phid(lat)))) # PYCHOK attr
469 return _aux(lat, inverse, Ellipsoid.auxConformal)
471 def auxGeocentric(self, lat, inverse=False, height=_0_0):
472 '''Compute the I{geocentric} auxiliary latitude (Theta) or the I{inverse} thereof.
474 @arg lat: The geodetic (or I{geocentric}) latitude (C{degrees90}).
475 @kwarg inverse: If C{True}, B{C{lat}} is the geocentric and
476 return the I{geocentric} latitude (C{bool}).
477 @kwarg height: Optional, ellipsoidal height (C{meter}).
479 @return: The I{geocentric} (or geodetic) latitude in C{degrees90}.
481 @see: U{Inverse-/GeocentricLatitude<https://GeographicLib.SourceForge.io/
482 C++/doc/classGeographicLib_1_1Ellipsoid.html>}, U{Geocentric latitude
483 <https://WikiPedia.org/wiki/Latitude#Geocentric_latitude>}, and
484 U{Snyder<https://Pubs.USGS.gov/pp/1395/report.pdf>}, pp 17-18.
485 '''
486 if self.f: # and lat
487 t = tan(Phid(lat))
488 f = self.b2_a2
489 if height:
490 if inverse:
491 lat = atan1d(t * f) # geodetic n
492 d, f = f, _1_0
493 else:
494 d = _1_0
495 n = self.rocPrimeVertical(lat)
496 f = _over(n * f + height, n * d + height)
497 elif inverse:
498 f = self.a2_b2
499 lat = atan1d(t * f)
500 return _aux(lat, inverse, Ellipsoid.auxGeocentric)
502 def auxIsometric(self, lat, inverse=False):
503 '''Compute the I{isometric} auxiliary latitude (Psi) or the I{inverse} thereof.
505 @arg lat: The geodetic (or I{isometric}) latitude (C{degrees}).
506 @kwarg inverse: If C{True}, B{C{lat}} is the I{isometric} and
507 return the geodetic latitude (C{bool}).
509 @return: The I{isometric} (or geodetic) latitude in C{degrees}.
511 @note: The I{isometric} latitude for geodetic C{+/-90} is far
512 outside the C{[-90..+90]} range but the inverse thereof
513 is the original geodetic latitude.
515 @see: U{Inverse-/IsometricLatitude<https://GeographicLib.SourceForge.io/
516 C++/doc/classGeographicLib_1_1Ellipsoid.html>}, U{Isometric latitude
517 <https://WikiPedia.org/wiki/Latitude#Isometric_latitude>}, and
518 U{Snyder<https://Pubs.USGS.gov/pp/1395/report.pdf>}, pp 15-16.
519 '''
520 if self.f:
521 r = Phid(lat, clip=0)
522 lat = degrees(atan1(self.es_tauf(sinh(r))) if inverse else
523 asinh(self.es_taupf(tan(r))))
524 # clip=0, since auxIsometric(+/-90) is far outside [-90..+90]
525 return _aux(lat, inverse, Ellipsoid.auxIsometric, clip=0)
527 def auxParametric(self, lat, inverse=False):
528 '''Compute the I{parametric} auxiliary latitude (Beta) or the I{inverse} thereof.
530 @arg lat: The geodetic (or I{parametric}) latitude (C{degrees90}).
531 @kwarg inverse: If C{True}, B{C{lat}} is the I{parametric} and
532 return the geodetic latitude (C{bool}).
534 @return: The I{parametric} (or geodetic) latitude in C{degrees90}.
536 @see: U{Inverse-/ParametricLatitude<https://GeographicLib.SourceForge.io/
537 C++/doc/classGeographicLib_1_1Ellipsoid.html>}, U{Parametric latitude
538 <https://WikiPedia.org/wiki/Latitude#Parametric_(or_reduced)_latitude>},
539 and U{Snyder<https://Pubs.USGS.gov/pp/1395/report.pdf>}, p 18.
540 '''
541 if self.f:
542 lat = self._Betad(Lat(lat), inverse=inverse)
543 return _aux(lat, inverse, Ellipsoid.auxParametric)
545 auxReduced = auxParametric # synonymous
547 def auxRectifying(self, lat, inverse=False):
548 '''Compute the I{rectifying} auxiliary latitude (Mu) or the I{inverse} thereof.
550 @arg lat: The geodetic (or I{rectifying}) latitude (C{degrees90}).
551 @kwarg inverse: If C{True}, B{C{lat}} is the I{rectifying} and
552 return the geodetic latitude (C{bool}).
554 @return: The I{rectifying} (or geodetic) latitude in C{degrees90}.
556 @see: U{Inverse-/RectifyingLatitude<https://GeographicLib.SourceForge.io/
557 C++/doc/classGeographicLib_1_1Ellipsoid.html>}, U{Rectifying latitude
558 <https://WikiPedia.org/wiki/Latitude#Rectifying_latitude>}, and
559 U{Snyder<https://Pubs.USGS.gov/pp/1395/report.pdf>}, pp 16-17.
560 '''
561 if self.f:
562 lat = Lat(lat)
563 if 0 < fabs(lat) < _90_0:
564 if inverse:
565 e = self._elliptic_e22
566 d = degrees90(e.fEinv(e.cE * lat / _90_0))
567 lat = self.auxParametric(d, inverse=True)
568 else:
569 lat = _over(self.Llat(lat), self.L) * _90_0
570 return _aux(lat, inverse, Ellipsoid.auxRectifying)
572 @Property_RO
573 def b(self):
574 '''Get the I{polar} radius, semi-axis (C{meter}).
575 '''
576 return self._b
578 polaradius = b # = Rpolar
580 @Property_RO
581 def b_a(self):
582 '''Get the ratio I{polar} over I{equatorial} radius (C{float}), M{b / a == f1 == 1 - f}.
584 @see: Property L{f1}.
585 '''
586 return self._assert(self.b / self.a, b_a=self.f1, f0=_1_0)
588 @Property_RO
589 def b2(self):
590 '''Get the I{polar} radius I{squared} (C{float}), M{b**2}.
591 '''
592 return Meter2(b2=self.b**2)
594 @Property_RO
595 def b2_a(self):
596 '''Get the I{equatorial} meridional radius of curvature (C{meter}), M{b**2 / a}, see C{rocMeridional}C{(0)}.
598 @see: U{Radii of Curvature<https://WikiPedia.org/wiki/Earth_radius#Radii_of_curvature>}.
599 '''
600 return Radius(b2_a=_over(self.b2, self.a) if self.f else self.b)
602 @Property_RO
603 def b2_a2(self):
604 '''Get the ratio I{polar} over I{equatorial} radius I{squared} (C{float}), M{(b / a)**2}
605 == M{(1 - f)**2} == M{1 - e**2} == C{e21}.
606 '''
607 return Float(b2_a2=self.b_a**2 if self.f else _1_0)
609 def _Betad(self, lat, inverse=False):
610 '''(INTERNAL) Get the I{parametric (or reduced) auxiliary latitude} or inverse thereof.
611 '''
612 s, c = sincos2d(lat) # like Karney's tand(lat)
613 s *= self.a_b if inverse else self.b_a
614 return atan1d(s, c) # (s * a, c * b) if inverse else (s * b, c * a)
616 @Property_RO
617 def BetaKs(self):
618 '''Get the I{Krüger} U{Beta series coefficients<https://GeographicLib.SourceForge.io/C++/doc/tmseries30.html>} (C{KsOrder}C{-tuple}).
619 '''
620 return self._Kseries( # XXX int/int quotients may require from __future__ import division as _; del _ # noqa: E702 ;
621 # n n**2 n**3 n**4 n**5 n**6 n**7 n**8
622 _T(1/2, -2/3, 37/96, -1/360, -81/512, 96199/604800, -5406467/38707200, 7944359/67737600),
623 _T(1/48, 1/15, -437/1440, 46/105, -1118711/3870720, 51841/1209600, 24749483/348364800), # PYCHOK unaligned
624 _T(17/480, -37/840, -209/4480, 5569/90720, 9261899/58060800, -6457463/17740800), # PYCHOK unaligned
625 _T(4397/161280, -11/504, -830251/7257600, 466511/2494800, 324154477/7664025600), # PYCHOK unaligned
626 _T(4583/161280, -108847/3991680, -8005831/63866880, 22894433/124540416), # PYCHOK unaligned
627 _T(20648693/638668800, -16363163/518918400, -2204645983/12915302400), # PYCHOK unaligne
628 _T(219941297/5535129600, -497323811/12454041600), # PYCHOK unaligned
629 _T(191773887257/3719607091200)) # PYCHOK unaligned
631 @deprecated_Property_RO
632 def c(self): # PYCHOK no cover
633 '''DEPRECATED, use property C{R2} or C{Rauthalic}.'''
634 return self.R2
636 @Property_RO
637 def c2(self):
638 '''Get the I{authalic} earth radius I{squared} (C{meter} I{squared}).
640 @see: Properties L{c2x}, L{area}, L{R2}, L{Rauthalic}, I{Karney's} U{equation (60)
641 <https://Link.Springer.com/article/10.1007%2Fs00190-012-0578-z>} and C++ U{Ellipsoid.Area
642 <https://GeographicLib.SourceForge.io/C++/doc/classGeographicLib_1_1Ellipsoid.html>},
643 U{Authalic radius<https://WikiPedia.org/wiki/Earth_radius#Authalic_radius>}, U{Surface area
644 <https://WikiPedia.org/wiki/Ellipsoid>} and U{surface area
645 <https://www.Numericana.com/answer/geometry.htm#oblate>}.
646 '''
647 return self._c2f(False)
649 @Property_RO
650 def c2x(self):
651 '''Get the I{authalic} earth radius I{squared} (C{meter} I{squared}), more accurate for very I{oblate}
652 ellipsoids.
654 @see: Properties L{c2}, L{areax}, L{R2x}, L{Rauthalicx}, class L{GeodesicExact} and I{Karney}'s comments at C++
655 attribute U{GeodesicExact._c2<https://GeographicLib.SourceForge.io/C++/doc/GeodesicExact_8cpp_source.html>}.
656 '''
657 return self._c2f(True)
659 def _c2f(self, c2x):
660 '''(INTERNAL) Helper for C{.c2} and C{.c2x}.
661 '''
662 f, c2 = self.f, self.b2
663 if f:
664 e = self.e
665 if e > EPS0:
666 if f > 0: # .isOblate
667 c2 *= (asinh(sqrt(self.e22abs)) if c2x else atanh(e)) / e
668 elif f < 0: # .isProlate
669 c2 *= atan1(e) / e # XXX asin?
670 c2 = Meter2(c2=(self.a2 + c2) * _0_5)
671 return c2
673 def circle4(self, lat):
674 '''Get the equatorial or a parallel I{circle of latitude}.
676 @arg lat: Geodetic latitude (C{degrees90} or C{str}).
678 @return: A L{Circle4Tuple}C{(radius, height, lat, beta)}.
680 @raise RangeError: Latitude B{C{lat}} outside valid range, only if
681 L{rangerrors<pygeodesy.rangerrors>} is C{True}.
683 @raise TypeError: Invalid B{C{lat}}.
685 @raise ValueError: Invalid B{C{lat}}.
687 @see: Definition of U{I{p} and I{z} under B{Parametric (or reduced) latitude}
688 <https://WikiPedia.org/wiki/Latitude>}, I{Karney's} C++ U{CircleRadius and CircleHeight
689 <https://GeographicLib.SourceForge.io/C++/doc/classGeographicLib_1_1Ellipsoid.html>}
690 and method C{Triaxial.ellipse5}.
691 '''
692 lat = _Lat0(lat)
693 if lat:
694 B = lat # beta
695 if fabs(lat) < _90_0:
696 if self.f:
697 B = self._Betad(lat)
698 z, r = sincos2d(B)
699 r *= self.a
700 z *= self.b
701 else: # near-polar
702 r, z = _0_0, copysign0(self.b, lat)
703 else: # equatorial
704 r = self.a
705 z = lat = B = _0_0
706 return Circle4Tuple(r, z, lat, B)
708 def degrees2m(self, deg, lat=0):
709 '''Convert an arc in degrees to a distance along the equator or
710 along a parallel of (geodetic) latitude.
712 @arg deg: The angle (C{degrees}).
713 @kwarg lat: Parallel latitude (C{degrees90}, C{str}).
715 @return: Distance (C{meter}, same units as the equatorial
716 and polar radii) or C{0} for near-polar B{C{lat}}.
718 @raise RangeError: Latitude B{C{lat}} outside valid range and
719 L{rangerrors<pygeodesy.rangerrors>} is C{True}.
721 @raise ValueError: Invalid B{C{deg}} or B{C{lat}}.
722 '''
723 return self.radians2m(radians(deg), lat=lat)
725 def distance2(self, lat0, lon0, lat1, lon1):
726 '''I{Approximate} the distance and (initial) bearing between
727 two points based on the U{local, flat earth approximation
728 <https://www.EdWilliams.org/avform.htm#flat>} aka U{Hubeny
729 <https://www.OVG.AT/de/vgi/files/pdf/3781/>} formula.
731 I{Suitable only for distances of several hundred Km or Miles
732 and only between points not near-polar}.
734 @arg lat0: From latitude (C{degrees}).
735 @arg lon0: From longitude (C{degrees}).
736 @arg lat1: To latitude (C{degrees}).
737 @arg lon1: To longitude (C{degrees}).
739 @return: A L{Distance2Tuple}C{(distance, initial)} with C{distance}
740 in same units as this ellipsoid's axes.
742 @note: The meridional and prime_vertical radii of curvature are
743 taken and scaled I{at the initial latitude}, see C{roc2}.
745 @see: Function L{pygeodesy.flatLocal}/L{pygeodesy.hubeny}.
746 '''
747 phi0 = Phid(lat0=lat0)
748 m, n = self.roc2_(phi0, scaled=True)
749 m *= Phid(lat1=lat1) - phi0
750 n *= Lamd(lon1=lon1) - Lamd(lon0=lon0)
751 return Distance2Tuple(hypot(m, n), atan2b(n, m))
753 @Property_RO
754 def e(self):
755 '''Get the I{unsigned, (1st) eccentricity} (C{float}), M{sqrt(1 - (b / a)**2))}, see C{a_b2e}.
757 @see: Property L{es}.
758 '''
759 return Float(e=sqrt(self.e2abs) if self.e2 else _0_0)
761 @deprecated_Property_RO
762 def e12(self): # see property ._e12
763 '''DEPRECATED, use property C{e21}.'''
764 return self.e21
766# @Property_RO
767# def _e12(self): # see property ._elliptic_e12
768# # (INTERNAL) until e12 above can be replaced with e21.
769# return self.e2 / (_1_0 - self.e2) # see I{Karney}'s Ellipsoid._e12 = e2 / (1 - e2)
771 @Property_RO
772 def e2(self):
773 '''Get the I{signed, (1st) eccentricity squared} (C{float}), M{f * (2 - f)
774 == 1 - (b / a)**2}, see C{a_b2e2}.
775 '''
776 return self._assert(a_b2e2(self.a, self.b), e2=f2e2(self.f))
778 @Property_RO
779 def e2abs(self):
780 '''Get the I{unsigned, (1st) eccentricity squared} (C{float}).
781 '''
782 return fabs(self.e2)
784 @Property_RO
785 def e21(self):
786 '''Get 1 less I{1st eccentricity squared} (C{float}), M{1 - e**2}
787 == M{1 - e2} == M{(1 - f)**2} == M{b**2 / a**2}, see C{b2_a2}.
788 '''
789 return self._assert((_1_0 - self.f)**2, e21=_1_0 - self.e2, f0=_1_0)
791# _e2m = e21 # see I{Karney}'s Ellipsoid._e2m = 1 - _e2
792 _1_e21 = a2_b2 # == M{1 / e21} == M{1 / (1 - e**2)}
794 @Property_RO
795 def e22(self):
796 '''Get the I{signed, 2nd eccentricity squared} (C{float}), M{e2 / (1 - e2)
797 == e2 / (1 - f)**2 == (a / b)**2 - 1}, see C{a_b2e22}.
798 '''
799 return self._assert(a_b2e22(self.a, self.b), e22=f2e22(self.f))
801 @Property_RO
802 def e22abs(self):
803 '''Get the I{unsigned, 2nd eccentricity squared} (C{float}).
804 '''
805 return fabs(self.e22)
807 @Property_RO
808 def e32(self):
809 '''Get the I{signed, 3rd eccentricity squared} (C{float}), M{e2 / (2 - e2)
810 == (a**2 - b**2) / (a**2 + b**2)}, see C{a_b2e32}.
811 '''
812 return self._assert(a_b2e32(self.a, self.b), e32=f2e32(self.f))
814 @Property_RO
815 def e32abs(self):
816 '''Get the I{unsigned, 3rd eccentricity squared} (C{float}).
817 '''
818 return fabs(self.e32)
820 @Property_RO
821 def e4(self):
822 '''Get the I{unsignd, (1st) eccentricity} to 4th power (C{float}), M{e**4 == e2**2}.
823 '''
824 return Float(e4=self.e2**2 if self.e2 else _0_0)
826 eccentricity = e # eccentricity
827# eccentricity2 = e2 # eccentricity squared
828 eccentricity1st2 = e2 # first eccentricity squared, signed
829 eccentricity2nd2 = e22 # second eccentricity squared, signed
830 eccentricity3rd2 = e32 # third eccentricity squared, signed
832 def ecef(self, Ecef=None):
833 '''Return U{ECEF<https://WikiPedia.org/wiki/ECEF>} converter.
835 @kwarg Ecef: ECEF class to use, default L{EcefKarney}.
837 @return: An ECEF converter for this C{ellipsoid}.
839 @raise TypeError: Invalid B{C{Ecef}}.
841 @see: Module L{pygeodesy.ecef}.
842 '''
843 return _MODS.ecef._4Ecef(self, Ecef)
845 @Property_RO
846 def _elliptic_e12(self): # see I{Karney}'s Ellipsoid._e12
847 '''(INTERNAL) Elliptic helper for C{Rhumb}.
848 '''
849 e12 = _over(self.e2, self.e2 - _1_0) # NOT DEPRECATED .e12!
850 return _MODS.elliptic.Elliptic(e12)
852 @Property_RO
853 def _elliptic_e22(self): # aka ._elliptic_ep2
854 '''(INTERNAL) Elliptic helper for C{auxRectifying}, C{L}, C{Llat}.
855 '''
856 return _MODS.elliptic.Elliptic(-self.e22abs) # complex
858 equatoradius = a # Requatorial
860 @Property_RO
861 def equatorimeter(self):
862 '''Get the ellipsoid's I{equatorial} perimeter (C{meter}).
863 '''
864 return Meter(equatorimeter=self.a * PI2)
866 def e2s(self, s):
867 '''Compute norm M{sqrt(1 - e2 * s**2)}.
869 @arg s: Sine value (C{scalar}).
871 @return: Norm (C{float}).
873 @raise ValueError: Invalid B{C{s}}.
874 '''
875 return sqrt(self.e2s2(s)) if self.e2 else _1_0
877 def e2s2(self, s):
878 '''Compute M{1 - e2 * s**2}.
880 @arg s: Sine value (C{scalar}).
882 @return: Result (C{float}).
884 @raise ValueError: Invalid B{C{s}}.
885 '''
886 r, e2 = _1_0, self.e2
887 if e2: # and s
888 try:
889 r -= e2 * Scalar(s=s)**2
890 if r < 0:
891 raise ValueError(_negative_)
892 except (TypeError, ValueError) as x:
893 t = self._DOT_(typename(Ellipsoid.e2s2))
894 raise _ValueError(t, s, cause=x)
895 return r
897 @Property_RO
898 def es(self):
899 '''Get the I{signed (1st) eccentricity} (C{float}).
901 @see: Property L{e}.
902 '''
903 # note, self.e is always non-negative
904 return Float(es=copysign0(self.e, self.f)) # see .ups
906 def es_atanh(self, x):
907 '''Compute M{es * atanh(es * x)} or M{-es * atan(es * x)}
908 for I{oblate} respectively I{prolate} ellipsoids where
909 I{es} is the I{signed} (1st) eccentricity.
911 @raise ValueError: Invalid B{C{x}}.
913 @see: Function U{Math::eatanhe<https://GeographicLib.SourceForge.io/
914 C++/doc/classGeographicLib_1_1Math.html>}.
915 '''
916 return self._es_atanh(Scalar(x=x)) if self.f else _0_0
918 def _es_atanh(self, x): # see .albers._atanhee, .AuxLat._atanhee
919 '''(INTERNAL) Helper for .es_atanh, ._es_taupf2 and ._exp_es_atanh.
920 '''
921 es = self.es # signOf(es) == signOf(f)
922 return es * (atanh(es * x) if es > 0 else # .isOblate
923 (-atan(es * x) if es < 0 else # .isProlate
924 _0_0)) # .isSpherical
926 @Property_RO
927 def es_c(self):
928 '''Get M{(1 - f) * exp(es_atanh(1))} (C{float}), M{b_a * exp(es_atanh(1))}.
929 '''
930 return Float(es_c=(self._exp_es_atanh_1 * self.b_a) if self.f else _1_0)
932 def es_tauf(self, taup):
933 '''Compute I{Karney}'s U{equations (19), (20) and (21)
934 <https://ArXiv.org/abs/1002.1417>}.
936 @see: I{Karney}'s C++ method U{Math::tauf<https://GeographicLib.
937 SourceForge.io/C++/doc/classGeographicLib_1_1Math.html>} and
938 and I{Veness}' JavaScript method U{toLatLon<https://www.
939 Movable-Type.co.UK/scripts/latlong-utm-mgrs.html>}.
940 '''
941 t = Scalar(taup=taup)
942 if self.f: # .isEllipsoidal
943 a = fabs(t)
944 T = (self._exp_es_atanh_1 if a > 70 else self._1_e21) * t
945 if fabs(T * _EPSqrt) < _2_0: # handles +/- INF and NAN
946 s = (a * _TOL) if a > _1_0 else _TOL
947 for T, _, d in self._es_tauf3(t, T): # max 2
948 if fabs(d) < s:
949 break
950 t = Scalar(tauf=T)
951 return t
953 def _es_tauf3(self, taup, T, N=9): # in .utm.Utm._toLLEB
954 '''(INTERNAL) Yield a 3-tuple C{(τi, iteration, delta)} for at most
955 B{C{N}} Newton iterations, converging rapidly except when C{delta}
956 toggles on +/-1.12e-16 or +/-4.47e-16, see C{.utm.Utm._toLLEB}.
957 '''
958 e = self._1_e21
959 _F2_ = Fsum(T).fsum2f_ # τ0
960 _tf2 = self._es_taupf2
961 for i in range(1, N + 1):
962 a, h = _tf2(T)
963 # = (taup - a) / hypot1(a) / ((e + T**2) / h)
964 d = _over((taup - a) * (T**2 + e), hypot1(a) * h)
965 T, d = _F2_(d) # τi, (τi - τi-1)
966 yield T, i, d
968 def es_taupf(self, tau):
969 '''Compute I{Karney}'s U{equations (7), (8) and (9)
970 <https://ArXiv.org/abs/1002.1417>}.
972 @see: I{Karney}'s C++ method U{Math::taupf<https://GeographicLib.
973 SourceForge.io/C++/doc/classGeographicLib_1_1Math.html>}.
974 '''
975 t = Scalar(tau=tau)
976 if self.f: # .isEllipsoidal
977 t, _ = self._es_taupf2(t)
978 t = Scalar(taupf=t)
979 return t
981 def _es_taupf2(self, tau):
982 '''(INTERNAL) Return 2-tuple C{(es_taupf(tau), hypot1(tau))}.
983 '''
984 if _isfinite(tau):
985 h = hypot1(tau)
986 s = sinh(self._es_atanh(tau / h))
987 a = hypot1(s) * tau - h * s
988 else:
989 a, h = tau, INF
990 return a, h
992 @Property_RO
993 def _exp_es_atanh_1(self):
994 '''(INTERNAL) Helper for .es_c and .es_tauf.
995 '''
996 return exp(self._es_atanh(_1_0)) if self.es else _1_0
998 @Property_RO
999 def f(self):
1000 '''Get the I{flattening} (C{scalar}), M{(a - b) / a}, C{0} for spherical, negative for prolate.
1001 '''
1002 return self._f
1004 @Property_RO
1005 def f_(self):
1006 '''Get the I{inverse flattening} (C{scalar}), M{1 / f} == M{a / (a - b)}, C{0} for spherical, see C{a_b2f_}.
1007 '''
1008 return self._f_
1010 @Property_RO
1011 def f1(self):
1012 '''Get the I{1 - flattening} (C{float}), M{f1 == 1 - f == b / a}.
1014 @see: Property L{b_a}.
1015 '''
1016 return Float(f1=_1_0 - self.f)
1018 @Property_RO
1019 def f2(self):
1020 '''Get the I{2nd flattening} (C{float}), M{(a - b) / b == f / (1 - f)}, C{0} for spherical, see C{a_b2f2}.
1021 '''
1022 return self._assert(self.a_b - _1_0, f2=f2f2(self.f))
1024 @deprecated_Property_RO
1025 def geodesic(self):
1026 '''DEPRECATED, use property C{geodesicw}.'''
1027 return self.geodesicw
1029 def geodesic_(self, exact=True):
1030 '''Get the an I{exact} C{Geodesic...} instance for this ellipsoid.
1032 @kwarg exact: If C{bool} return L{GeodesicExact}C{(exact=B{exact}, ...)},
1033 otherwise a L{Geodesic}, L{GeodesicExact} or L{GeodesicSolve}
1034 instance for I{this} ellipsoid.
1036 @return: The C{exact} geodesic (C{Geodesic...}).
1038 @raise TypeError: Invalid B{C{exact}}.
1040 @raise ValueError: Incompatible B{C{exact}} ellipsoid.
1041 '''
1042 if isbool(exact): # for consistenccy with C{.rhumb_}
1043 g = _MODS.geodesicx.GeodesicExact(self, C4order=30 if exact else 24,
1044 name=self.name)
1045 else:
1046 g = exact
1047 E = _xattr(g, ellipsoid=None)
1048 if not (E is self and isinstance(g, self._Geodesics)):
1049 raise _ValueError(exact=g, ellipsoid=E, txt_not_=self.name)
1050 return g
1052 @property_ROver
1053 def _Geodesics(self):
1054 '''(INTERNAL) Get all C{Geodesic...} classes, I{once}.
1055 '''
1056 t = (_MODS.geodesicx.GeodesicExact,
1057 _MODS.geodsolve.GeodesicSolve)
1058 try:
1059 t += (_MODS.geodesicw.Geodesic,
1060 _MODS.geodesicw._wrapped.Geodesic)
1061 except ImportError:
1062 pass
1063 return t # overwrite property_ROver
1065 @property_RO
1066 def geodesicw(self):
1067 '''Get this ellipsoid's I{wrapped} U{geodesicw.Geodesic
1068 <https://GeographicLib.SourceForge.io/Python/doc/code.html>}, provided
1069 I{Karney}'s U{geographiclib<https://PyPI.org/project/geographiclib>}
1070 package is installed.
1071 '''
1072 # if not self.isEllipsoidal:
1073 # raise _IsnotError(_ellipsoidal_, ellipsoid=self)
1074 return _MODS.geodesicw.Geodesic(self)
1076 @property_RO
1077 def geodesicx(self):
1078 '''Get this ellipsoid's I{exact} L{GeodesicExact}.
1079 '''
1080 # if not self.isEllipsoidal:
1081 # raise _IsnotError(_ellipsoidal_, ellipsoid=self)
1082 return _MODS.geodesicx.GeodesicExact(self, name=self.name)
1084 @property
1085 def geodsolve(self):
1086 '''Get this ellipsoid's L{GeodesicSolve}, the I{wrapper} around utility
1087 U{GeodSolve<https://GeographicLib.SourceForge.io/C++/doc/GeodSolve.1.html>},
1088 provided the path to the C{GeodSolve} executable is specified with env
1089 variable C{PYGEODESY_GEODSOLVE} or re-/set with this property..
1090 '''
1091 # if not self.isEllipsoidal:
1092 # raise _IsnotError(_ellipsoidal_, ellipsoid=self)
1093 return _MODS.geodsolve.GeodesicSolve(self, path=self._geodsolve, name=self.name)
1095 @geodsolve.setter # PYCHOK setter!
1096 def geodsolve(self, path):
1097 '''Re-/set the (fully qualified) path to the U{GeodSolve
1098 <https://GeographicLib.SourceForge.io/C++/doc/GeodSolve.1.html>} executable,
1099 overriding env variable C{PYGEODESY_GEODSOLVE} (C{str}).
1100 '''
1101 self._geodsolve = path
1103 def hartzell4(self, pov, los=None):
1104 '''Compute the intersection of this ellipsoid's surface and a Line-Of-Sight
1105 from a Point-Of-View in space.
1107 @arg pov: Point-Of-View outside this ellipsoid (C{Cartesian}, L{Ecef9Tuple}
1108 or L{Vector3d}).
1109 @kwarg los: Line-Of-Sight, I{direction} to this ellipsoid (L{Los}, L{Vector3d})
1110 or C{True} for the I{normal, perpendicular, plumb} to the surface
1111 of this ellipsoid or C{False} or C{None} to point to its center.
1113 @return: L{Vector4Tuple}C{(x, y, z, h)} with the cartesian coordinates C{x},
1114 C{y} and C{z} of the projection on or the intersection with this
1115 ellipsoid and the I{distance} C{h} from B{C{pov}} to C{(x, y, z)}
1116 along B{C{los}}, all in C{meter}, conventionally.
1118 @raise IntersectionError: Null B{C{pov}} or B{C{los}} vector, or B{C{pov}}
1119 is inside this ellipsoid or B{C{los}} points
1120 outside this ellipsoid or in opposite direction.
1122 @raise TypeError: Invalid B{C{pov}} or B{C{los}}.
1124 @see: U{I{Satellite Line-of-Sight Intersection with Earth}<https://StephenHartzell.
1125 Medium.com/satellite-line-of-sight-intersection-with-earth-d786b4a6a9b6>} and
1126 methods L{Ellipsoid.height4} and L{Triaxial.hartzell4}.
1127 '''
1128 try:
1129 m = _MODS._triaxials_triaxial5
1130 v, d, i = m._hartzell3(pov, los, self._triaxial)
1131 except Exception as x:
1132 raise IntersectionError(pov=pov, los=los, cause=x)
1133 return Vector4Tuple(v.x, v.y, v.z, d, iteration=i, name__=self.hartzell4)
1135 @Property_RO
1136 def _hash(self):
1137 return hash((self.a, self.f))
1139 def height4(self, xyz, normal=True):
1140 '''Compute the projection on and the height of a cartesian above or below
1141 this ellipsoid's surface.
1143 @arg xyz: The cartesian (C{Cartesian}, L{Ecef9Tuple}, L{Vector3d},
1144 L{Vector3Tuple} or L{Vector4Tuple}).
1145 @kwarg normal: If C{True}, the projection is perpendicular to (the nearest
1146 point on) this ellipsoid's surface, otherwise the C{radial}
1147 line to this ellipsoid's center (C{bool}).
1149 @return: L{Vector4Tuple}C{(x, y, z, h)} with the cartesian coordinates C{x},
1150 C{y} and C{z} of the projection on and the height C{h} above or
1151 below this ellipsoid's surface, all in C{meter}, conventionally.
1153 @raise ValueError: Null B{C{xyz}}.
1155 @raise TypeError: Non-cartesian B{C{xyz}}.
1157 @see: U{Distance to<https://StackOverflow.com/questions/22959698/distance-from-given-point-to-given-ellipse>}
1158 and U{intersection with<https://MathWorld.Wolfram.com/Ellipse-LineIntersection.html>} an ellipse and
1159 methods L{Ellipsoid.hartzell4} and L{Triaxial.height4}.
1160 '''
1161 v = _MODS.vector3d._otherV3d(xyz=xyz)
1162 r = v.length
1164 a, b, i = self.a, self.b, None
1165 if r < EPS0: # EPS
1166 v = v.times(_0_0)
1167 h = -a
1169 elif self.isSpherical:
1170 v = v.times(a / r)
1171 h = r - a
1173 elif normal: # perpendicular to ellipsoid
1174 x, y = hypot(v.x, v.y), fabs(v.z)
1175 if x < EPS0: # PYCHOK no cover
1176 z = copysign0(b, v.z)
1177 v = Vector3Tuple(v.x, v.y, z)
1178 h = y - b # polar
1179 elif y < EPS0: # PYCHOK no cover
1180 t = a / r
1181 v = v.times_(t, t, 0) # force z=0.0
1182 h = x - a # equatorial
1183 else: # normal in 1st quadrant
1184 m = _MODS._triaxials_triaxial5
1185 x, y, i = m._plumbTo3(x, y, self)
1186 t, v = v, v.times_(x, x, y)
1187 h = t.minus(v).length
1189 else: # radial to ellipsoid's center
1190 h = hypot_(a * v.z, b * v.x, b * v.y)
1191 t = (a * b / h) if h > EPS0 else _0_0 # EPS
1192 v = v.times(t)
1193 h = r * (_1_0 - t)
1195 return Vector4Tuple(v.x, v.y, v.z, h, iteration=i, name__=self.height4)
1197 def _heightB(self, sa, ca, z, p): # in ecef.EcefSudano, ecec.EcefVeness
1198 '''(INTERNAL) Height above ellipsoid (Bowring eqn 7) at C{lat}.
1199 '''
1200 # sa, ca = sincos2d(lat)
1201 # p = hypot(x, y) # distance to polar axis
1203 # r = a / self.e2s(sa) # length of normal terminated by polar axis
1204 # h = p * ca + z * sa - (a * a / r)
1205 return (p * ca + fabs(z * sa) - self.a * self.e2s(sa)) if sa else (p - self.a)
1207 @Property_RO
1208 def _heightMax(self):
1209 '''(INTERNAL) Get the height limit (C{meter}, conventionally).
1210 '''
1211 return self.a / EPS_2 # self.a * _2_EPS, about 12M lightyears
1213 def _hubeny_2(self, phi2, phi1, lam21, scaled=True, squared=True):
1214 '''(INTERNAL) like function C{pygeodesy.flatLocal_}/C{pygeodesy.hubeny_},
1215 returning the I{angular} distance in C{radians squared} or C{radians}
1216 '''
1217 m, n = self.roc2_((phi2 + phi1) * _0_5, scaled=scaled)
1218 h, r = (hypot2, self.a2_) if squared else (hypot, _1_0 / self.a)
1219 return h(m * (phi2 - phi1), n * lam21) * r
1221 @Property_RO
1222 def isEllipsoidal(self):
1223 '''Is this model I{ellipsoidal} (C{bool})?
1224 '''
1225 return self.f != 0
1227 @Property_RO
1228 def isOblate(self):
1229 '''Is this ellipsoid I{oblate} (C{bool})? I{Prolate} or
1230 spherical otherwise.
1231 '''
1232 return self.f > 0
1234 @Property_RO
1235 def isProlate(self):
1236 '''Is this ellipsoid I{prolate} (C{bool})? I{Oblate} or
1237 spherical otherwise.
1238 '''
1239 return self.f < 0
1241 @Property_RO
1242 def isSpherical(self):
1243 '''Is this ellipsoid I{spherical} (C{bool})?
1244 '''
1245 return self.f == 0
1247 def _Kseries(self, *AB8Ks):
1248 '''(INTERNAL) Compute the 4-, 6- or 8-th order I{Krüger} Alpha
1249 or Beta series coefficients per I{Karney}'s U{equations (35)
1250 and (36)<https://ArXiv.org/pdf/1002.1417v3.pdf>}.
1252 @arg AB8Ks: 8-Tuple of 8-th order I{Krüger} Alpha or Beta series
1253 coefficient tuples.
1255 @return: I{Krüger} series coefficients (L{KsOrder}C{-tuple}).
1257 @see: I{Karney}'s 30-th order U{TMseries30
1258 <https://GeographicLib.SourceForge.io/C++/doc/tmseries30.html>}.
1259 '''
1260 k = self.KsOrder
1261 if self.n:
1262 ns = fpowers(self.n, k)
1263 ks = tuple(fdot(AB8Ks[i][:k-i], *ns[i:]) for i in range(k))
1264 else:
1265 ks = _0_0s(k)
1266 return ks
1268 @property_doc_(''' the I{Krüger} series' order (C{int}), see properties C{AlphaKs}, C{BetaKs}.''')
1269 def KsOrder(self):
1270 '''Get the I{Krüger} series' order (C{int} 4, 6 or 8).
1271 '''
1272 return self._KsOrder
1274 @KsOrder.setter # PYCHOK setter!
1275 def KsOrder(self, order):
1276 '''Set the I{Krüger} series' order (C{int} 4, 6 or 8).
1278 @raise ValueError: Invalid B{C{order}}.
1279 '''
1280 if not (isint(order) and _isin(order, 4, 6, 8)):
1281 raise _ValueError(order=order)
1282 if self._KsOrder != order:
1283 Ellipsoid.AlphaKs._update(self)
1284 Ellipsoid.BetaKs._update(self)
1285 self._KsOrder = order
1287 @Property_RO
1288 def L(self):
1289 '''Get the I{quarter meridian} C{L}, aka the C{polar distance}
1290 along a meridian between the equator and a pole (C{meter}),
1291 M{b * Elliptic(-e2 / (1 - e2)).cE} or M{b * PI / 2}.
1292 '''
1293 r = self._elliptic_e22.cE if self.f else PI_2
1294 return Distance(L=self.b * r)
1296 def Llat(self, lat):
1297 '''Return the I{meridional length}, the distance along a meridian
1298 between the equator and a (geodetic) latitude, see C{L}.
1300 @arg lat: Geodetic latitude (C{degrees90}).
1302 @return: The meridional length at B{C{lat}}, negative on southern
1303 hemisphere (C{meter}).
1304 '''
1305 r = self._elliptic_e22.fEd(self.auxParametric(lat)) if self.f else Phid(lat)
1306 return Distance(Llat=self.b * r)
1308 Lmeridian = Llat # meridional distance
1310 @property_RO
1311 def _Lpd(self):
1312 '''Get the I{quarter meridian} per degree (C{meter}), M{self.L / 90}.
1313 '''
1314 return Meter(_Lpd=self.L / _90_0)
1316 @property_RO
1317 def _Lpr(self):
1318 '''Get the I{quarter meridian} per radian (C{meter}), M{self.L / PI_2}.
1319 '''
1320 return Meter(_Lpr=self.L / PI_2)
1322 @deprecated_Property_RO
1323 def majoradius(self): # PYCHOK no cover
1324 '''DEPRECATED, use property C{a} or C{Requatorial}.'''
1325 return self.a
1327 def m2degrees(self, distance, lat=0):
1328 '''Convert a distance to an arc in degrees along the equator or
1329 along a parallel of (geodetic) latitude.
1331 @arg distance: Distance (C{meter}).
1332 @kwarg lat: Parallel latitude (C{degrees90}, C{str}).
1334 @return: Angle (C{degrees}) or C{0} for near-polar B{C{lat}}.
1336 @raise RangeError: Latitude B{C{lat}} outside valid range and
1337 L{rangerrors<pygeodesy.rangerrors>} is C{True}.
1339 @raise ValueError: Invalid B{C{distance}} or B{C{lat}}.
1340 '''
1341 return degrees(self.m2radians(distance, lat=lat))
1343 def m2radians(self, distance, lat=0):
1344 '''Convert a distance to an angle along the equator or along
1345 a parallel of (geodetic) latitude.
1347 @arg distance: Distance (C{meter}).
1348 @kwarg lat: Parallel latitude (C{degrees90}, C{str}).
1350 @return: Angle (C{radians}) or C{INF} for near-polar B{C{lat}}.
1352 @raise RangeError: Latitude B{C{lat}} outside valid range and
1353 L{rangerrors<pygeodesy.rangerrors>} is C{True}.
1355 @raise ValueError: Invalid B{C{distance}} or B{C{lat}}.
1356 '''
1357 r = self.circle4(lat).radius if lat else self.a
1358 return m2radians(distance, radius=r, lat=0)
1360 @deprecated_Property_RO
1361 def minoradius(self): # PYCHOK no cover
1362 '''DEPRECATED, use property C{b}, C{polaradius} or C{Rpolar}.'''
1363 return self.b
1365 @Property_RO
1366 def n(self):
1367 '''Get the I{3rd flattening} (C{float}), M{f / (2 - f) == (a - b) / (a + b)}, see C{a_b2n}.
1368 '''
1369 return self._assert(a_b2n(self.a, self.b), n=f2n(self.f))
1371 flattening = f
1372 flattening1st = f
1373 flattening2nd = f2
1374 flattening3rd = n
1376 polaradius = b # Rpolar
1378 @property_RO
1379 def polarimeter(self):
1380 '''Get the ellipsoid's I{polar}, meridional perimeter (C{meter}).
1381 '''
1382 return Meter(polarimeter=self.L * _4_0)
1384# Q = A # I{meridian arc unit} C{Q}, the mean, meridional length I{per radian}
1386 @deprecated_Property_RO
1387 def quarteradius(self): # PYCHOK no cover
1388 '''DEPRECATED, use property C{L} or method C{Llat}.'''
1389 return self.L
1391 @Property_RO
1392 def R1(self):
1393 '''Get the I{mean} earth radius per I{IUGG} (C{meter}), M{(2 * a + b) / 3 == a * (1 - f / 3)}.
1395 @see: U{Earth radius<https://WikiPedia.org/wiki/Earth_radius>}
1396 and method C{Rgeometric}.
1397 '''
1398 return Radius(R1=fmean_(self.a, self.a, self.b) if self.f else self.a)
1400 Rmean = R1
1402 @Property_RO
1403 def R2(self):
1404 '''Get the I{authalic} earth radius (C{meter}), M{sqrt(c2)}.
1406 @see: C{R2x}, C{c2}, C{area} and U{Earth radius
1407 <https://WikiPedia.org/wiki/Earth_radius>}.
1408 '''
1409 return Radius(R2=sqrt(self.c2) if self.f else self.a)
1410# # Moritz, H. <https://Geodesy.Geology.Ohio-State.edu/course/refpapers/00740128.pdf>
1411# # R2 = (1 - 2/3 * e'2 + 26/45 * e'4 - 100/189 * e'6 + 7034/14175 * e'8) * rocPolar
1412# # = (3 + e'2 * (-2 + e'2 * (26/15 + e'2 * (-100/63 + e'2 * 7034/4725)))) * rocPolar / 3
1413# return Fhorner(self.e22, 3, -2, 26 / 15, -100 / 63, 7034 / 4725).fover(3 / self.rocPolar)
1415 Rauthalic = R2
1417 @Property_RO
1418 def R2x(self):
1419 '''Get the I{authalic} earth radius (C{meter}), M{sqrt(c2x)}.
1421 @see: C{R2}, C{c2x} and C{areax}.
1422 '''
1423 return Radius(R2x=sqrt(self.c2x) if self.f else self.a)
1425 Rauthalicx = R2x
1427 @Property_RO
1428 def R3(self):
1429 '''Get the I{volumetric} earth radius (C{meter}), M{(a * a * b)**(1/3)}.
1431 @see: U{Earth radius<https://WikiPedia.org/wiki/Earth_radius>} and C{volume}.
1432 '''
1433 r = (cbrt(self.b_a) * self.a) if self.f else self.a
1434 return Radius(R3=r)
1436 Rvolumetric = R3
1438 def radians2m(self, rad, lat=0):
1439 '''Convert an angle to the distance along the equator or along
1440 a parallel of (geodetic) latitude.
1442 @arg rad: The angle (C{radians}).
1443 @kwarg lat: Parallel latitude (C{degrees90}, C{str}).
1445 @return: Distance (C{meter}, same units as the equatorial
1446 and polar radii) or C{0} for near-polar B{C{lat}}.
1448 @raise RangeError: Latitude B{C{lat}} outside valid range and
1449 L{rangerrors<pygeodesy.rangerrors>} is C{True}.
1451 @raise ValueError: Invalid B{C{rad}} or B{C{lat}}.
1452 '''
1453 r = self.circle4(lat).radius if lat else self.a
1454 return radians2m(rad, radius=r, lat=0)
1456 @Property_RO
1457 def Rbiaxial(self):
1458 '''Get the I{biaxial, quadratic} mean earth radius (C{meter}), M{sqrt((a**2 + b**2) / 2)}.
1460 @see: C{Rtriaxial}
1461 '''
1462 a, b = self.a, self.b
1463 if b != a:
1464 b = hypot(a, b) / _SQRT2
1465 return Radius(Rbiaxial=b)
1467 Requatorial = a # for consistent naming
1469 def Rgeocentric(self, lat):
1470 '''Compute the I{geocentric} earth radius of (geodetic) latitude.
1472 @arg lat: Latitude (C{degrees90}).
1474 @return: Geocentric earth radius (C{meter}).
1476 @raise ValueError: Invalid B{C{lat}}.
1478 @see: U{Geocentric Radius
1479 <https://WikiPedia.org/wiki/Earth_radius#Geocentric_radius>}
1480 '''
1481 r, p = self.a, Phid(lat)
1482 if p and self.f:
1483 if fabs(p) < PI_2:
1484 s2, c2 = _sin2cos2(p)
1485 # R == sqrt((a2**2 * c2 + b2**2 * s2) / (a2 * c2 + b2 * s2))
1486 # == sqrt(a2**2 * (c2 + (b2 / a2)**2 * s2) / (a2 * (c2 + b2 / a2 * s2)))
1487 # == sqrt(a2 * (c2 + (b2 / a2)**2 * s2) / (c2 + (b2 / a2) * s2))
1488 # == a * sqrt((c2 + b2_a2 * b2_a2 * s2) / (c2 + b2_a2 * s2))
1489 s2 *= self.b2_a2
1490 r *= sqrt((c2 + self.b2_a2 * s2) / (c2 + s2))
1491 else:
1492 r = self.b
1493 return Radius(Rgeocentric=r)
1495 @Property_RO
1496 def Rgeometric(self):
1497 '''Get the I{geometric} mean earth radius (C{meter}), M{sqrt(a * b)}.
1499 @see: C{R1}.
1500 '''
1501 g = sqrt(self.a * self.b) if self.f else self.a
1502 return Radius(Rgeometric=g)
1504 def rhumb_(self, exact=True):
1505 '''Get the an I{exact} C{Rhumb...} instance for this ellipsoid.
1507 @kwarg exact: If C{bool} or C{None} return L{Rhumb}C{(exact=B{exact}, ...)},
1508 otherwise a L{Rhumb}, L{RhumbAux} or L{RhumbSolve} instance
1509 for I{this} ellipsoid.
1511 @return: The C{exact} rhumb (C{Rhumb...}).
1513 @raise TypeError: Invalid B{C{exact}}.
1515 @raise ValueError: Incompatible B{C{exact}} ellipsoid.
1516 '''
1517 if isbool(exact): # use Rhumb for backward compatibility
1518 r = _MODS.rhumb.ekx.Rhumb(self, exact=exact, name=self.name)
1519 else:
1520 r = exact
1521 E = _xattr(r, ellipsoid=None)
1522 if not (E is self and isinstance(r, self._Rhumbs)):
1523 raise _ValueError(exact=r, ellipsosid=E, txt_not_=self.name)
1524 return r
1526 @property_RO
1527 def rhumbaux(self):
1528 '''Get this ellipsoid's I{Auxiliary} C{rhumb.RhumbAux}.
1529 '''
1530 # if not self.isEllipsoidal:
1531 # raise _IsnotError(_ellipsoidal_, ellipsoid=self)
1532 return _MODS.rhumb.aux_.RhumbAux(self, name=self.name)
1534 @property_RO
1535 def rhumbekx(self):
1536 '''Get this ellipsoid's I{Elliptic, Krüger} C{rhumb.Rhumb}.
1537 '''
1538 # if not self.isEllipsoidal:
1539 # raise _IsnotError(_ellipsoidal_, ellipsoid=self)
1540 return _MODS.rhumb.ekx.Rhumb(self, name=self.name)
1542 @property_ROver
1543 def _Rhumbs(self):
1544 '''(INTERNAL) Get all C{Rhumb...} classes, I{once}.
1545 '''
1546 r = _MODS.rhumb
1547 return (r.aux_.RhumbAux, # overwrite property_ROver
1548 r.ekx.Rhumb, r.solve.RhumbSolve)
1550 @property
1551 def rhumbsolve(self):
1552 '''Get this ellipsoid's L{RhumbSolve}, the I{wrapper} around utility
1553 U{RhumbSolve<https://GeographicLib.SourceForge.io/C++/doc/GeodSolve.1.html>},
1554 provided the path to the C{RhumbSolve} executable is specified with env
1555 variable C{PYGEODESY_RHUMBSOLVE} or re-/set with this property.
1556 '''
1557 # if not self.isEllipsoidal:
1558 # raise _IsnotError(_ellipsoidal_, ellipsoid=self)
1559 return _MODS.rhumb.solve.RhumbSolve(self, path=self._rhumbsolve, name=self.name)
1561 @rhumbsolve.setter # PYCHOK setter!
1562 def rhumbsolve(self, path):
1563 '''Re-/set the (fully qualified) path to the U{RhumbSolve
1564 <https://GeographicLib.SourceForge.io/C++/doc/GeodSolve.1.html>} executable,
1565 overriding env variable C{PYGEODESY_RHUMBSOLVE} (C{str}).
1566 '''
1567 self._rhumbsolve = path
1569 @deprecated_property_RO
1570 def rhumbx(self):
1571 '''DEPRECATED on 2023.11.28, use property C{rhumbekx}.'''
1572 return self.rhumbekx
1574 def Rlat(self, lat):
1575 '''I{Average} the earth radius between C{equatoradius} at C{0} and
1576 C{polaradius} at C{+/-90} degrees latitude.
1578 @arg lat: Latitude (C{degrees90}, C{str} or C{Ang}).
1580 @return: Averaged earth radius (C{meter}).
1582 @raise RangeError: Latitude B{C{lat}} outside valid range, only if
1583 L{rangerrors<pygeodesy.rangerrors>} is C{True}.
1585 @raise TypeError: Invalid B{C{lat}}.
1587 @raise ValueError: Invalid B{C{lat}}.
1588 '''
1589 # r = a - (a - b) * |lat| / 90
1590 r = self.a
1591 if self.f: # .isEllipsoidal
1592 lat = _Lat0(lat)
1593 if lat:
1594 r -= (r - self.b) * fabs(lat) / _90_0
1595 r = Radius(Rlat=r)
1596 return r
1598 Rpolar = b # for consistent naming
1600 def roc1_(self, sa, ca=None):
1601 '''Compute the I{prime-vertical}, I{normal} radius of curvature
1602 of (geodetic) latitude, I{unscaled}.
1604 @arg sa: Sine of the latitude (C{float}, [-1.0..+1.0]).
1605 @kwarg ca: Optional cosine of the latitude (C{float}, [-1.0..+1.0])
1606 to use an alternate formula.
1608 @return: The prime-vertical radius of curvature (C{float}).
1610 @note: The delta between both formulae with C{Ellipsoids.WGS84}
1611 is less than 2 nanometer over the entire latitude range.
1613 @see: Method L{roc2_} and class L{EcefYou}.
1614 '''
1615 if sa and self.f: # .isEllipsoidal
1616 if ca is None:
1617 r = self.e2s2(sa) # see .roc2_ and _EcefBase._forward
1618 n = sqrt(self.a2 / r) if r > EPS02 else _0_0
1619 elif ca: # derived from EcefYou.forward
1620 h = hypot(ca, self.b_a * sa)
1621 n = self.a / h
1622 else:
1623 n = self.a2_b / fabs(sa)
1624 else:
1625 n = self.a
1626 return n
1628 def roc2(self, lat, scaled=False):
1629 '''Compute the I{meridional} and I{prime-vertical}, I{normal}
1630 radii of curvature of (geodetic) latitude.
1632 @arg lat: Latitude (C{degrees90}).
1633 @kwarg scaled: Scale prime_vertical by C{cos(radians(B{lat}))} (C{bool}).
1635 @return: An L{Curvature2Tuple}C{(meridional, prime_vertical)} with
1636 the radii of curvature.
1638 @raise ValueError: Invalid B{C{lat}}.
1640 @see: Methods L{roc2_} and L{roc1_}, U{Local, flat earth approximation
1641 <https://www.EdWilliams.org/avform.htm#flat>} and meridional and
1642 prime vertical U{Radii of Curvature<https://WikiPedia.org/wiki/
1643 Earth_radius#Radii_of_curvature>}.
1644 '''
1645 return self.roc2_(Phid(lat), scaled=scaled)
1647 def roc2_(self, phi, scaled=False):
1648 '''Compute the I{meridional} and I{prime-vertical}, I{normal} radii of
1649 curvature of (geodetic) latitude.
1651 @arg phi: Latitude (C{radians}).
1652 @kwarg scaled: Scale prime_vertical by C{cos(B{phi})} (C{bool}).
1654 @return: An L{Curvature2Tuple}C{(meridional, prime_vertical)} with the
1655 radii of curvature.
1657 @raise ValueError: Invalid B{C{phi}}.
1659 @see: Methods L{roc2} and L{roc1_}, property L{rocEquatorial2}, U{Local,
1660 flat earth approximation<https://www.EdWilliams.org/avform.htm#flat>}
1661 and the meridional and prime vertical U{Radii of Curvature
1662 <https://WikiPedia.org/wiki/Earth_radius#Radii_of_curvature>}.
1663 '''
1664 p = fabs(Phi(phi))
1665 if self.f:
1666 r = self.e2s2(sin(p))
1667 if r > EPS02:
1668 n = sqrt(self.a2 / r)
1669 m = n * self.e21 / r
1670 else:
1671 m = n = _0_0
1672 else:
1673 m = n = self.a
1674 if scaled and p:
1675 n *= cos(p) if p < PI_2 else _0_0
1676 return Curvature2Tuple(m, n)
1678 def rocAzimuth(self, lat, azimuth):
1679 '''Compute the I{directional} radius of curvature of (geodetic) latitude
1680 and C{azimuth} compass direction.
1682 @see: Method L{rocBearing<Ellipsoid.rocBearing>} for details, using C{azimuth} for C{bearing}.
1683 '''
1684 return Radius(rocAzimuth=self._rocDirectional(lat, Azimuth(azimuth)))
1686 def rocBearing(self, lat, bearing):
1687 '''Compute the I{directional} radius of curvature of (geodetic) latitude
1688 and C{bearing} compass direction.
1690 @arg lat: Latitude (C{degrees90}).
1691 @arg bearing: Direction (compass C{degrees360}).
1693 @return: Directional radius of curvature (C{meter}).
1695 @raise RangeError: Latitude B{C{lat}} outside valid range and
1696 L{rangerrors<pygeodesy.rangerrors>} is C{True}.
1698 @raise ValueError: Invalid B{C{lat}} or B{C{bearing}}.
1700 @see: U{Radii of Curvature<https://WikiPedia.org/wiki/Earth_radius#Radii_of_curvature>}
1701 '''
1702 return Radius(rocBearing=self._rocDirectional(lat, Bearing(bearing)))
1704 def _rocDirectional(self, lat, deg):
1705 '''(INTERNAL) Helper for C{rocAzimuth} and C{rocBearing}.
1706 '''
1707 if self.f:
1708 s2, c2 = _sin2cos2(radians(deg))
1709 m, n = self.roc2_(Phid(lat))
1710 if n < m: # == n / (c2 * n / m + s2)
1711 c2 *= n / m
1712 elif m < n: # == m / (c2 + s2 * m / n)
1713 s2 *= m / n
1714 n = m
1715 r = _over(n, c2 + s2) # == 1 / (c2 / m + s2 / n)
1716 else:
1717 r = self.b # == self.a
1718 return r
1720 @Property_RO
1721 def rocEquatorial2(self):
1722 '''Get the I{meridional} and I{prime-vertical}, I{normal} radii of curvature
1723 at the equator as L{Curvature2Tuple}C{(meridional, prime_vertical)}.
1725 @see: Methods L{rocMeridional} and L{rocPrimeVertical}, properties L{b2_a},
1726 L{a2_b}, C{rocPolar} and polar and equatorial U{Radii of Curvature
1727 <https://WikiPedia.org/wiki/Earth_radius#Radii_of_curvature>}.
1728 '''
1729 m = self.b2_a if self.f else self.a
1730 return Curvature2Tuple(m, self.a)
1732 def rocGauss(self, lat):
1733 '''Compute the I{Gaussian} radius of curvature of (geodetic) latitude.
1735 @arg lat: Latitude (C{degrees90}).
1737 @return: Gaussian radius of curvature (C{meter}).
1739 @raise ValueError: Invalid B{C{lat}}.
1741 @see: Non-directional U{Radii of Curvature<https://WikiPedia.org/wiki/
1742 Earth_radius#Radii_of_curvature>}
1743 '''
1744 # using ...
1745 # m, n = self.roc2_(Phid(lat))
1746 # return sqrt(m * n)
1747 # ... requires 1 or 2 sqrt
1748 g = self.b
1749 if self.f:
1750 s2, c2 = _sin2cos2(Phid(lat))
1751 g = _over(g, c2 + self.b2_a2 * s2)
1752 return Radius(rocGauss=g)
1754 def rocMean(self, lat):
1755 '''Compute the I{mean} radius of curvature of (geodetic) latitude.
1757 @arg lat: Latitude (C{degrees90}).
1759 @return: Mean radius of curvature (C{meter}).
1761 @raise ValueError: Invalid B{C{lat}}.
1763 @see: Non-directional U{Radii of Curvature<https://WikiPedia.org/wiki/
1764 Earth_radius#Radii_of_curvature>}
1765 '''
1766 if self.f:
1767 m, n = self.roc2_(Phid(lat))
1768 m *= _over(n * _2_0, m + n) # == 2 / (1 / m + 1 / n)
1769 else:
1770 m = self.a
1771 return Radius(rocMean=m)
1773 def rocMeridional(self, lat):
1774 '''Compute the I{meridional} radius of curvature of (geodetic) latitude.
1776 @arg lat: Latitude (C{degrees90}).
1778 @return: Meridional radius of curvature (C{meter}).
1780 @raise ValueError: Invalid B{C{lat}}.
1782 @see: Methods L{roc2} and L{roc2_}, U{Local, flat earth approximation
1783 <https://www.EdWilliams.org/avform.htm#flat>} and U{Radii of
1784 Curvature<https://WikiPedia.org/wiki/Earth_radius#Radii_of_curvature>}.
1785 '''
1786 r = self.roc2_(Phid(lat)) if lat else self.rocEquatorial2
1787 return Radius(rocMeridional=r.meridional)
1789 rocPolar = a2_b # synonymous
1791 def rocPrimeVertical(self, lat):
1792 '''Compute the I{prime-vertical}, I{normal} radius of curvature of
1793 (geodetic) latitude, aka the I{transverse} radius of curvature.
1795 @arg lat: Latitude (C{degrees90}).
1797 @return: Prime-vertical radius of curvature (C{meter}).
1799 @raise ValueError: Invalid B{C{lat}}.
1801 @see: Methods L{roc2}, L{roc2_} and L{roc1_}, U{Local, flat earth
1802 approximation<https://www.EdWilliams.org/avform.htm#flat>} and
1803 U{Radii of Curvature<https://WikiPedia.org/wiki/
1804 Earth_radius#Radii_of_curvature>}.
1805 '''
1806 r = self.roc1_(sin(Phid(lat))) if lat else self.a
1807 return Radius(rocPrimeVertical=r)
1809 rocTransverse = rocPrimeVertical # synonymous
1811 @deprecated_Property_RO
1812 def Rquadratic(self): # PYCHOK no cover
1813 '''DEPRECATED, use property C{Rbiaxial} or C{Rtriaxial}.'''
1814 return self.Rbiaxial
1816 @deprecated_Property_RO
1817 def Rr(self): # PYCHOK no cover
1818 '''DEPRECATED, use property C{Rrectifying}.'''
1819 return self.Rrectifying
1821 @Property_RO
1822 def Rrectifying(self):
1823 '''Get the I{rectifying} earth radius (C{meter}), M{((a**(3/2) + b**(3/2)) / 2)**(2/3)}.
1825 @see: U{Earth radius<https://WikiPedia.org/wiki/Earth_radius>}.
1826 '''
1827 r = self.a
1828 if self.f:
1829 r *= cbrt2((sqrt3(self.b_a) + _1_0) * _0_5)
1830 return Radius(Rrectifying=r)
1832 @deprecated_Property_RO
1833 def Rs(self): # PYCHOK no cover
1834 '''DEPRECATED, use property C{Rgeometric}.'''
1835 return self.Rgeometric
1837 @Property_RO
1838 def Rtriaxial(self):
1839 '''Get the I{triaxial, quadratic} mean earth radius (C{meter}), M{sqrt((3 * a**2 + b**2) / 4)}.
1841 @see: C{Rbiaxial}
1842 '''
1843 q, b = self.a, self.b
1844 if b < q:
1845 q *= sqrt((self.b2_a2 + _3_0) * _0_25)
1846 elif b > q:
1847 q = sqrt((self.a2_b2 * _3_0 + _1_0) * _0_25) * b
1848 return Radius(Rtriaxial=q)
1850 def toEllipse(self, **name):
1851 '''Get this ellipsoid's meridional ellipse L{Ellipse<pygeodesy.Ellipse>}.
1853 @kwarg name: Optional, unique C{B{name}=NN} (C{str}).
1854 '''
1855 return _MODS.ellipses.Ellipse(self.a, self.b, **name)
1857 def toEllipsoid2(self, **name):
1858 '''Get a copy of this ellipsoid as an L{Ellipsoid2}.
1860 @kwarg name: Optional, unique C{B{name}=NN} (C{str}).
1862 @see: Property C{a_f}.
1863 '''
1864 return Ellipsoid2(self, None, **name)
1866 def toStr(self, prec=8, terse=4, **sep_name): # PYCHOK expected
1867 '''Return this ellipsoid as a text string.
1869 @kwarg prec: Number of decimal digits, unstripped (C{int}).
1870 @kwarg terse: Limit the number of items (C{int}, 0...18),
1871 use C{B{terse}=0} or C{=None} for all.
1872 @kwarg sep_name: Optional C{B{name}=NN} (C{str}) or C{None}
1873 to exclude this ellipsoid's name and separator
1874 C{B{sep}=", "} to join the items (C{str}).
1876 @return: This C{Ellipsoid}'s attributes (C{str}).
1877 '''
1878 E = Ellipsoid
1879 t = (E.a, E.f, E.f_, E.b, E.f2, E.n, E.e,
1880 E.e2, E.e21, E.e22, E.e32,
1881 E.A, E.L, E.R1, E.R2, E.R3,
1882 E.Rbiaxial, E.Rtriaxial)
1883 if terse:
1884 t = t[:terse]
1885 return self._instr(prec=prec, props=t, **sep_name)
1887 def toTriaxial(self, **name):
1888 '''Convert this ellipsoid to a L{Triaxial_}.
1890 @kwarg name: Optional C{B{name}=NN} (C{str}).
1892 @return: A L{Triaxial_} or L{Triaxial} with the C{X} axis
1893 pointing east and C{Z} pointing north.
1895 @see: Method L{Triaxial_.toEllipsoid}.
1896 '''
1897 T = self._triaxial
1898 return T.copy(**name) if name else T
1900 @property_RO
1901 def _triaxial(self):
1902 '''(INTERNAL) Get this ellipsoid's un-/ordered C{Triaxial/_}.
1903 '''
1904 a, b, t = self.a, self.b, _MODS.triaxials
1905 T = t.Triaxial if a > b else t.Triaxial_
1906 return T(a, a, b, name=self.name)
1908 @Property_RO
1909 def volume(self):
1910 '''Get the ellipsoid's I{volume} (C{meter**3}), M{4 / 3 * PI * R3**3}.
1912 @see: C{R3}.
1913 '''
1914 return Meter3(volume=self.a2 * self.b * PI_3 * _4_0)
1917class Ellipsoid2(Ellipsoid):
1918 '''An L{Ellipsoid} specified by I{equatorial} radius and I{flattening}.
1919 '''
1920 def __init__(self, a, f=None, **name):
1921 '''New L{Ellipsoid2}.
1923 @arg a: Equatorial radius, semi-axis (C{meter}) or a previous
1924 L{Ellipsoid} instance.
1925 @arg f: Flattening: (C{float} < 1.0, negative for I{prolate}),
1926 if B{C{a}} is in C{meter}.
1927 @kwarg name: Optional, unique C{B{name}=NN} (C{str}).
1929 @raise NameError: Ellipsoid with that B{C{name}} already exists.
1931 @raise ValueError: Invalid B{C{a}} or B{C{f}}.
1933 @note: C{abs(B{f}) < EPS} is forced to C{B{f}=0}, I{spherical}.
1934 Negative C{B{f}} produces a I{prolate} ellipsoid.
1935 '''
1936 if f is None and isinstance(a, Ellipsoid):
1937 Ellipsoid.__init__(self, a.a, f =a.f,
1938 b=a.b, f_=a.f_, **name)
1939 else:
1940 Ellipsoid.__init__(self, a, f=f, **name)
1943def _ispherical_a_b(a, b):
1944 '''(INTERNAL) C{True} for spherical or invalid C{a} or C{b}.
1945 '''
1946 return a < EPS0 or b < EPS0 or fabs(a - b) < EPS0
1949def _ispherical_f(f):
1950 '''(INTERNAL) C{True} for spherical or invalid C{f}.
1951 '''
1952 return f > EPS1 or fabs(f) < EPS
1955def _ispherical_f_(f_):
1956 '''(INTERNAL) C{True} for spherical or invalid C{f_}.
1957 '''
1958 f_ = fabs(f_)
1959 return f_ < EPS or f_ > _1_EPS
1962def a_b2e(a, b):
1963 '''Return C{e}, the I{1st eccentricity} for a given I{equatorial} and I{polar} radius.
1965 @arg a: Equatorial radius (C{scalar} > 0).
1966 @arg b: Polar radius (C{scalar} > 0).
1968 @return: The I{unsigned}, (1st) eccentricity (C{float} or C{0}), M{sqrt(1 - (b / a)**2)}.
1970 @note: The result is always I{non-negative} and C{0} for I{near-spherical} ellipsoids.
1971 '''
1972 e2 = _a2b2e2(a, b, b2=False)
1973 return Float(e=sqrt(fabs(e2)) if e2 else _0_0) # == sqrt(fabs((a - b) * (a + b))) / a
1976def a_b2e2(a, b):
1977 '''Return C{e2}, the I{1st eccentricity squared} for a given I{equatorial} and I{polar} radius.
1979 @arg a: Equatorial radius (C{scalar} > 0).
1980 @arg b: Polar radius (C{scalar} > 0).
1982 @return: The I{signed}, (1st) eccentricity I{squared} (C{float} or C{0}), M{1 - (b / a)**2}.
1984 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
1985 for I{near-spherical} ellipsoids.
1986 '''
1987 return Float(e2=_a2b2e2(a, b, b2=False))
1990def a_b2e22(a, b):
1991 '''Return C{e22}, the I{2nd eccentricity squared} for a given I{equatorial} and I{polar} radius.
1993 @arg a: Equatorial radius (C{scalar} > 0).
1994 @arg b: Polar radius (C{scalar} > 0).
1996 @return: The I{signed}, 2nd eccentricity I{squared} (C{float} or C{0}), M{(a / b)**2 - 1}.
1998 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
1999 for I{near-spherical} ellipsoids.
2000 '''
2001 return Float(e22=_a2b2e2(a, b, a2=False))
2004def a_b2e32(a, b):
2005 '''Return C{e32}, the I{3rd eccentricity squared} for a given I{equatorial} and I{polar} radius.
2007 @arg a: Equatorial radius (C{scalar} > 0).
2008 @arg b: Polar radius (C{scalar} > 0).
2010 @return: The I{signed}, 3rd eccentricity I{squared} (C{float} or C{0}),
2011 M{(a**2 - b**2) / (a**2 + b**2)}.
2013 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2014 for I{near-spherical} ellipsoids.
2015 '''
2016 return Float(e32=_a2b2e2(a, b))
2019def _a2b2e2(a, b, a2=True, b2=True):
2020 '''(INTERNAL) Helper for C{a_b2e}, C{a_b2e2}, C{a_b2e22} and C{a_b2e32}.
2021 '''
2022 if _ispherical_a_b(a, b):
2023 e2 = _0_0
2024 else: # a > 0, b > 0
2025 a, b = (_1_0, b / a) if a > b else (a / b, _1_0)
2026 a2b2 = float(a - b) * (a + b)
2027 e2 = _over(a2b2, (a**2 if a2 else _0_0) +
2028 (b**2 if b2 else _0_0)) if a2b2 else _0_0
2029 return e2
2032def a_b2f(a, b):
2033 '''Return C{f}, the I{flattening} for a given I{equatorial} and I{polar} radius.
2035 @arg a: Equatorial radius (C{scalar} > 0).
2036 @arg b: Polar radius (C{scalar} > 0).
2038 @return: The flattening (C{scalar} or C{0}), M{(a - b) / a}.
2040 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2041 for I{near-spherical} ellipsoids.
2042 '''
2043 f = 0 if _ispherical_a_b(a, b) else _over(float(a - b), a)
2044 return _f_0_0 if _ispherical_f(f) else Float(f=f)
2047def a_b2f_(a, b):
2048 '''Return C{f_}, the I{inverse flattening} for a given I{equatorial} and I{polar} radius.
2050 @arg a: Equatorial radius (C{scalar} > 0).
2051 @arg b: Polar radius (C{scalar} > 0).
2053 @return: The inverse flattening (C{scalar} or C{0}), M{a / (a - b)}.
2055 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2056 for I{near-spherical} ellipsoids.
2057 '''
2058 f_ = 0 if _ispherical_a_b(a, b) else _over(a, float(a - b))
2059 return _f__0_0 if _ispherical_f_(f_) else Float(f_=f_)
2062def a_b2f2(a, b):
2063 '''Return C{f2}, the I{2nd flattening} for a given I{equatorial} and I{polar} radius.
2065 @arg a: Equatorial radius (C{scalar} > 0).
2066 @arg b: Polar radius (C{scalar} > 0).
2068 @return: The I{signed}, 2nd flattening (C{scalar} or C{0}), M{(a - b) / b}.
2070 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2071 for I{near-spherical} ellipsoids.
2072 '''
2073 t = 0 if _ispherical_a_b(a, b) else float(a - b)
2074 return Float(f2=_0_0 if fabs(t) < EPS0 else _over(t, b))
2077def a_b2n(a, b):
2078 '''Return C{n}, the I{3rd flattening} for a given I{equatorial} and I{polar} radius.
2080 @arg a: Equatorial radius (C{scalar} > 0).
2081 @arg b: Polar radius (C{scalar} > 0).
2083 @return: The I{signed}, 3rd flattening (C{scalar} or C{0}), M{(a - b) / (a + b)}.
2085 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2086 for I{near-spherical} ellipsoids.
2087 '''
2088 t = 0 if _ispherical_a_b(a, b) else float(a - b)
2089 return Float(n=_0_0 if fabs(t) < EPS0 else _over(t, a + b))
2092def a_f2b(a, f):
2093 '''Return C{b}, the I{polar} radius for a given I{equatorial} radius and I{flattening}.
2095 @arg a: Equatorial radius (C{scalar} > 0).
2096 @arg f: Flattening (C{scalar} < 1, negative for I{prolate}).
2098 @return: The polar radius (C{float}), M{a * (1 - f)}.
2099 '''
2100 b = a if _ispherical_f(f) else (a * (_1_0 - f))
2101 return Radius_(b=a if _ispherical_a_b(a, b) else b)
2104def a_f_2b(a, f_):
2105 '''Return C{b}, the I{polar} radius for a given I{equatorial} radius and I{inverse flattening}.
2107 @arg a: Equatorial radius (C{scalar} > 0).
2108 @arg f_: Inverse flattening (C{scalar} >>> 1).
2110 @return: The polar radius (C{float}), M{a * (f_ - 1) / f_}.
2111 '''
2112 b = a if _ispherical_f_(f_) else _over(a * (f_ - _1_0), f_)
2113 return Radius_(b=a if _ispherical_a_b(a, b) else b)
2116def b_f2a(b, f):
2117 '''Return C{a}, the I{equatorial} radius for a given I{polar} radius and I{flattening}.
2119 @arg b: Polar radius (C{scalar} > 0).
2120 @arg f: Flattening (C{scalar} < 1, negative for I{prolate}).
2122 @return: The equatorial radius (C{float}), M{b / (1 - f)}.
2123 '''
2124 t = _1_0 - f
2125 a = b if fabs(t) < EPS0 else _over(b, t)
2126 return Radius_(a=b if _ispherical_a_b(a, b) else a)
2129def b_f_2a(b, f_):
2130 '''Return C{a}, the I{equatorial} radius for a given I{polar} radius and I{inverse flattening}.
2132 @arg b: Polar radius (C{scalar} > 0).
2133 @arg f_: Inverse flattening (C{scalar} >>> 1).
2135 @return: The equatorial radius (C{float}), M{b * f_ / (f_ - 1)}.
2136 '''
2137 t = f_ - _1_0
2138 a = b if _ispherical_f_(f_) or fabs(t) < EPS0 \
2139 or fabs(t - f_) < EPS0 else _over(b * f_, t)
2140 return Radius_(a=b if _ispherical_a_b(a, b) else a)
2143def e2f(e):
2144 '''Return C{f}, the I{flattening} for a given I{1st eccentricity}.
2146 @arg e: The (1st) eccentricity (0 <= C{float} < 1)
2148 @return: The flattening (C{scalar} or C{0}).
2150 @see: Function L{e22f}.
2151 '''
2152 return e22f(e**2)
2155def e22f(e2):
2156 '''Return C{f}, the I{flattening} for a given I{1st eccentricity squared}.
2158 @arg e2: The (1st) eccentricity I{squared}, I{signed} (L{NINF} < C{float} < 1)
2160 @return: The flattening (C{float} or C{0}), M{e2 / (sqrt(1 - e2) + 1)}.
2161 '''
2162 return Float(f=_over(e2, sqrt(_1_0 - e2) + _1_0)) if e2 and e2 < _1_0 else _f_0_0
2165def f2e2(f):
2166 '''Return C{e2}, the I{1st eccentricity squared} for a given I{flattening}.
2168 @arg f: Flattening (C{scalar} < 1, negative for I{prolate}).
2170 @return: The I{signed}, (1st) eccentricity I{squared} (C{float} < 1), M{f * (2 - f)}.
2172 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2173 for I{near-spherical} ellipsoids.
2175 @see: U{Eccentricity conversions<https://GeographicLib.SourceForge.io/
2176 C++/doc/classGeographicLib_1_1Ellipsoid.html>} and U{Flattening
2177 <https://WikiPedia.org/wiki/Flattening>}.
2178 '''
2179 return Float(e2=_0_0 if _ispherical_f(f) else (f * (_2_0 - f)))
2182def f2e22(f):
2183 '''Return C{e22}, the I{2nd eccentricity squared} for a given I{flattening}.
2185 @arg f: Flattening (C{scalar} < 1, negative for I{prolate}).
2187 @return: The I{signed}, 2nd eccentricity I{squared} (C{float} > -1 or C{INF}),
2188 M{f * (2 - f) / (1 - f)**2}.
2190 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2191 for near-spherical ellipsoids.
2193 @see: U{Eccentricity conversions<https://GeographicLib.SourceForge.io/
2194 C++/doc/classGeographicLib_1_1Ellipsoid.html>}.
2195 '''
2196 # e2 / (1 - e2) == f * (2 - f) / (1 - f)**2
2197 t = (_1_0 - f)**2
2198 return Float(e22=INF if t < EPS0 else _over(f2e2(f), t)) # PYCHOK type
2201def f2e32(f):
2202 '''Return C{e32}, the I{3rd eccentricity squared} for a given I{flattening}.
2204 @arg f: Flattening (C{scalar} < 1, negative for I{prolate}).
2206 @return: The I{signed}, 3rd eccentricity I{squared} (C{float}),
2207 M{f * (2 - f) / (1 + (1 - f)**2)}.
2209 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2210 for I{near-spherical} ellipsoids.
2212 @see: U{Eccentricity conversions<https://GeographicLib.SourceForge.io/
2213 C++/doc/classGeographicLib_1_1Ellipsoid.html>}.
2214 '''
2215 # e2 / (2 - e2) == f * (2 - f) / (1 + (1 - f)**2)
2216 e2 = f2e2(f)
2217 return Float(e32=_over(e2, _2_0 - e2))
2220def f_2f(f_):
2221 '''Return C{f}, the I{flattening} for a given I{inverse flattening}.
2223 @arg f_: Inverse flattening (C{scalar} >>> 1).
2225 @return: The flattening (C{scalar} or C{0}), M{1 / f_}.
2227 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2228 for I{near-spherical} ellipsoids.
2229 '''
2230 f = 0 if _ispherical_f_(f_) else _over(_1_0, f_)
2231 return _f_0_0 if _ispherical_f(f) else Float(f=f) # PYCHOK type
2234def f2f_(f):
2235 '''Return C{f_}, the I{inverse flattening} for a given I{flattening}.
2237 @arg f: Flattening (C{scalar} < 1, negative for I{prolate}).
2239 @return: The inverse flattening (C{scalar} or C{0}), M{1 / f}.
2241 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2242 for I{near-spherical} ellipsoids.
2243 '''
2244 f_ = 0 if _ispherical_f(f) else _over(_1_0, f)
2245 return _f__0_0 if _ispherical_f_(f_) else Float(f_=f_) # PYCHOK type
2248def f2f2(f):
2249 '''Return C{f2}, the I{2nd flattening} for a given I{flattening}.
2251 @arg f: Flattening (C{scalar} < 1, negative for I{prolate}).
2253 @return: The I{signed}, 2nd flattening (C{scalar} or C{INF}), M{f / (1 - f)}.
2255 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2256 for I{near-spherical} ellipsoids.
2258 @see: U{Eccentricity conversions<https://GeographicLib.SourceForge.io/
2259 C++/doc/classGeographicLib_1_1Ellipsoid.html>} and U{Flattening
2260 <https://WikiPedia.org/wiki/Flattening>}.
2261 '''
2262 t = _1_0 - f
2263 return Float(f2=_0_0 if _ispherical_f(f) else
2264 (INF if fabs(t) < EPS else _over(f, t))) # PYCHOK type
2267def f2n(f):
2268 '''Return C{n}, the I{3rd flattening} for a given I{flattening}.
2270 @arg f: Flattening (C{scalar} < 1, negative for I{prolate}).
2272 @return: The I{signed}, 3rd flattening (-1 <= C{float} < 1),
2273 M{f / (2 - f)}.
2275 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2276 for I{near-spherical} ellipsoids.
2278 @see: U{Eccentricity conversions<https://GeographicLib.SourceForge.io/
2279 C++/doc/classGeographicLib_1_1Ellipsoid.html>} and U{Flattening
2280 <https://WikiPedia.org/wiki/Flattening>}.
2281 '''
2282 return Float(n=_0_0 if _ispherical_f(f) else _over(f, float(_2_0 - f)))
2285def n2e2(n):
2286 '''Return C{e2}, the I{1st eccentricity squared} for a given I{3rd flattening}.
2288 @arg n: The 3rd flattening (-1 <= C{scalar} < 1).
2290 @return: The I{signed}, (1st) eccentricity I{squared} (C{float} or NINF),
2291 M{4 * n / (1 + n)**2}.
2293 @note: The result is positive for I{oblate}, negative for I{prolate} or C{0}
2294 for I{near-spherical} ellipsoids.
2296 @see: U{Flattening<https://WikiPedia.org/wiki/Flattening>}.
2297 '''
2298 t = (n + _1_0)**2
2299 return Float(e2=_0_0 if fabs(n) < EPS0 else
2300 (NINF if t < EPS0 else _over(_4_0 * n, t)))
2303def n2f(n):
2304 '''Return C{f}, the I{flattening} for a given I{3rd flattening}.
2306 @arg n: The 3rd flattening (-1 <= C{scalar} < 1).
2308 @return: The flattening (C{scalar} or NINF), M{2 * n / (1 + n)}.
2310 @see: U{Eccentricity conversions<https://GeographicLib.SourceForge.io/
2311 C++/doc/classGeographicLib_1_1Ellipsoid.html>} and U{Flattening
2312 <https://WikiPedia.org/wiki/Flattening>}.
2313 '''
2314 t = n + _1_0
2315 f = 0 if fabs(n) < EPS0 else (NINF if t < EPS0 else _over(_2_0 * n, t))
2316 return _f_0_0 if _ispherical_f(f) else Float(f=f)
2319def n2f_(n):
2320 '''Return C{f_}, the I{inverse flattening} for a given I{3rd flattening}.
2322 @arg n: The 3rd flattening (-1 <= C{scalar} < 1).
2324 @return: The inverse flattening (C{scalar} or C{0}), M{1 / f}.
2326 @see: L{n2f} and L{f2f_}.
2327 '''
2328 return f2f_(n2f(n))
2331class Ellipsoids(_NamedEnum):
2332 '''(INTERNAL) L{Ellipsoid} registry, I{must} be a sub-class
2333 to accommodate the L{_LazyNamedEnumItem} properties.
2334 '''
2335 def _Lazy(self, a, b, f_, **kwds):
2336 '''(INTERNAL) Instantiate the L{Ellipsoid}.
2337 '''
2338 return Ellipsoid(a, b=b, f_=f_, **kwds)
2340Ellipsoids = Ellipsoids(Ellipsoid) # PYCHOK singleton
2341'''Some pre-defined L{Ellipsoid}s, all I{lazily} instantiated.'''
2342# <https://www.GNU.org/software/gama/manual/html_node/Supported-ellipsoids.html>
2343# <https://GSSC.ESA.int/navipedia/index.php/Reference_Frames_in_GNSS>
2344# <https://kb.OSU.edu/dspace/handle/1811/77986>
2345# <https://www.IBM.com/docs/en/db2/11.5?topic=systems-supported-spheroids>
2346# <https://w3.Energistics.org/archive/Epicentre/Epicentre_v3.0/DataModel/LogicalDictionary/StandardValues/ellipsoid.html>
2347# <https://GitHub.com/locationtech/proj4j/blob/master/src/main/java/org/locationtech/proj4j/datum/Ellipsoid.java>
2348Ellipsoids._assert( # <https://WikiPedia.org/wiki/Earth_ellipsoid>
2349 Airy1830 = _lazy(_Airy1830_, *_T(6377563.396, _0_0, 299.3249646)), # b=6356256.909
2350 AiryModified = _lazy(_AiryModified_, *_T(6377340.189, _0_0, 299.3249646)), # b=6356034.448
2351# APL4_9 = _lazy('APL4_9', *_T(6378137.0, _0_0, 298.24985392)), # Appl. Phys. Lab. 1965
2352# ANS = _lazy('ANS', *_T(6378160.0, _0_0, 298.25)), # Australian Nat. Spheroid
2353# AN_SA96 = _lazy('AN_SA96', *_T(6378160.0, _0_0, 298.24985392)), # Australian Nat. South America
2354 Australia1966 = _lazy('Australia1966', *_T(6378160.0, _0_0, 298.25)), # b=6356774.7192
2355 ATS1977 = _lazy('ATS1977', *_T(6378135.0, _0_0, 298.257)), # "Average Terrestrial System"
2356 Bessel1841 = _lazy(_Bessel1841_, *_T(6377397.155, 6356078.962818, 299.152812797)),
2357 BesselModified = _lazy('BesselModified', *_T(6377492.018, _0_0, 299.1528128)),
2358# BesselNamibia = _lazy('BesselNamibia', *_T(6377483.865, _0_0, 299.1528128)),
2359 CGCS2000 = _lazy('CGCS2000', *_T(R_MA, _0_0, 298.257222101)), # BeiDou Coord System (BDC)
2360# Clarke1858 = _lazy('Clarke1858', *_T(6378293.639, _0_0, 294.260676369)),
2361 Clarke1866 = _lazy(_Clarke1866_, *_T(6378206.4, 6356583.8, 294.978698214)),
2362 Clarke1880 = _lazy('Clarke1880', *_T(6378249.145, 6356514.86954978, 293.465)),
2363 Clarke1880IGN = _lazy(_Clarke1880IGN_, *_T(6378249.2, 6356515.0, 293.466021294)),
2364 Clarke1880Mod = _lazy('Clarke1880Mod', *_T(6378249.145, 6356514.96639549, 293.466307656)), # aka Clarke1880Arc
2365 CPM1799 = _lazy('CPM1799', *_T(6375738.7, 6356671.92557493, 334.39)), # Comm. des Poids et Mesures
2366 Delambre1810 = _lazy('Delambre1810', *_T(6376428.0, 6355957.92616372, 311.5)), # Belgium
2367 Engelis1985 = _lazy('Engelis1985', *_T(6378136.05, 6356751.32272154, 298.2566)),
2368# Everest1830 = _lazy('Everest1830', *_T(6377276.345, _0_0, 300.801699997)),
2369# Everest1948 = _lazy('Everest1948', *_T(6377304.063, _0_0, 300.801699997)),
2370# Everest1956 = _lazy('Everest1956', *_T(6377301.243, _0_0, 300.801699997)),
2371 Everest1969 = _lazy('Everest1969', *_T(6377295.664, 6356094.667915, 300.801699997)),
2372 Everest1975 = _lazy('Everest1975', *_T(6377299.151, 6356098.14512013, 300.8017255)),
2373 Fisher1968 = _lazy('Fisher1968', *_T(6378150.0, 6356768.33724438, 298.3)),
2374# Fisher1968Mod = _lazy('Fisher1968Mod', *_T(6378155.0, _0_0, 298.3)),
2375 GEM10C = _lazy('GEM10C', *_T(R_MA, 6356752.31424783, 298.2572236)),
2376 GPES = _lazy('GPES', *_T(6378135.0, 6356750.0, _0_0)), # "Gen. Purpose Earth Spheroid"
2377 GRS67 = _lazy('GRS67', *_T(6378160.0, _0_0, 298.247167427)), # Lucerne b=6356774.516
2378# GRS67Truncated = _lazy('GRS67Truncated', *_T(6378160.0, _0_0, 298.25)),
2379 GRS80 = _lazy(_GRS80_, *_T(R_MA, 6356752.314140347, 298.25722210088)), # IUGG, ITRS, ETRS89
2380# Hayford1924 = _lazy('Hayford1924', *_T(6378388.0, 6356911.94612795, None)), # aka Intl1924 f_=297
2381 Helmert1906 = _lazy('Helmert1906', *_T(6378200.0, 6356818.16962789, 298.3)),
2382# Hough1960 = _lazy('Hough1960', *_T(6378270.0, _0_0, 297.0)),
2383 IAU76 = _lazy('IAU76', *_T(6378140.0, _0_0, 298.257)), # Int'l Astronomical Union
2384 IERS1989 = _lazy('IERS1989', *_T(6378136.0, _0_0, 298.257)), # b=6356751.302
2385 IERS1992TOPEX = _lazy('IERS1992TOPEX', *_T(6378136.3, 6356751.61659215, 298.257223563)), # IERS/TOPEX/Poseidon/McCarthy
2386 IERS2003 = _lazy('IERS2003', *_T(6378136.6, 6356751.85797165, 298.25642)),
2387 Intl1924 = _lazy(_Intl1924_, *_T(6378388.0, _0_0, 297.0)), # aka Hayford b=6356911.9462795
2388 Intl1967 = _lazy('Intl1967', *_T(6378157.5, 6356772.2, 298.24961539)), # New Int'l
2389 Krassovski1940 = _lazy(_Krassovski1940_, *_T(6378245.0, 6356863.01877305, 298.3)), # spelling
2390 Krassowsky1940 = _lazy(_Krassowsky1940_, *_T(6378245.0, 6356863.01877305, 298.3)), # spelling
2391# Kaula = _lazy('Kaula', *_T(6378163.0, _0_0, 298.24)), # Kaula 1961
2392# Lerch = _lazy('Lerch', *_T(6378139.0, _0_0, 298.257)), # Lerch 1979
2393 Maupertuis1738 = _lazy('Maupertuis1738', *_T(6397300.0, 6363806.28272251, 191.0)), # France
2394 Mercury1960 = _lazy('Mercury1960', *_T(6378166.0, 6356784.28360711, 298.3)),
2395 Mercury1968Mod = _lazy('Mercury1968Mod', *_T(6378150.0, 6356768.33724438, 298.3)),
2396# MERIT = _lazy('MERIT', *_T(6378137.0, _0_0, 298.257)), # MERIT 1983
2397# NWL10D = _lazy('NWL10D', *_T(6378135.0, _0_0, 298.26)), # Naval Weapons Lab.
2398 NWL1965 = _lazy('NWL1965', *_T(6378145.0, 6356759.76948868, 298.25)), # Naval Weapons Lab.
2399# NWL9D = _lazy('NWL9D', *_T(6378145.0, 6356759.76948868, 298.25)), # NWL1965
2400 OSU86F = _lazy('OSU86F', *_T(6378136.2, 6356751.51693008, 298.2572236)),
2401 OSU91A = _lazy('OSU91A', *_T(6378136.3, 6356751.6165948, 298.2572236)),
2402# Plessis1817 = _lazy('Plessis1817', *_T(6397523.0, 6355863.0, 153.56512242)), # XXX incorrect?
2403 Plessis1817 = _lazy('Plessis1817', *_T(6376523.0, 6355862.93325557, 308.64)), # XXX IGN France 1972
2404# Prolate = _lazy('Prolate', *_T(6356752.3, R_MA, _0_0)),
2405 PZ90 = _lazy('PZ90', *_T(6378136.0, _0_0, 298.257839303)), # GLOSNASS PZ-90 and PZ-90.11
2406# SEAsia = _lazy('SEAsia', *_T(6378155.0, _0_0, 298.3)), # SouthEast Asia
2407 SGS85 = _lazy('SGS85', *_T(6378136.0, 6356751.30156878, 298.257)), # Soviet Geodetic System
2408 SoAmerican1969 = _lazy('SoAmerican1969', *_T(6378160.0, 6356774.71919531, 298.25)), # South American
2409 Sphere = _lazy(_Sphere_, *_T(R_M, R_M, _0_0)), # pseudo
2410 SphereAuthalic = _lazy('SphereAuthalic', *_T(R_FM, R_FM, _0_0)), # pseudo
2411 SpherePopular = _lazy('SpherePopular', *_T(R_MA, R_MA, _0_0)), # EPSG:3857 Spheroid
2412 Struve1860 = _lazy('Struve1860', *_T(6378298.3, 6356657.14266956, 294.73)),
2413# Walbeck = _lazy('Walbeck', *_T(6376896.0, _0_0, 302.78)),
2414# WarOffice = _lazy('WarOffice', *_T(6378300.0, _0_0, 296.0)),
2415 WGS60 = _lazy('WGS60', *_T(6378165.0, 6356783.28695944, 298.3)),
2416 WGS66 = _lazy('WGS66', *_T(6378145.0, 6356759.76948868, 298.25)),
2417 WGS72 = _lazy(_WGS72_, *_T(6378135.0, _0_0, 298.26)), # b=6356750.52
2418 WGS84 = _lazy(_WGS84_, *_T(R_MA, _0_0, _f__WGS84)), # b=6356752.3142451793
2419# U{NOAA/NOS/NGS/inverse<https://GitHub.com/noaa-ngs/inverse/blob/main/invers3d.f>}
2420 WGS84_NGS = _lazy('WGS84_NGS', *_T(R_MA, _0_0, 298.257222100882711243162836600094))
2421)
2423_EWGS84 = Ellipsoids.WGS84 # (INTERNAL) shared
2425if __name__ == _DMAIN_:
2427 from pygeodesy import nameof, printf
2429 for E in (_EWGS84, Ellipsoids.GRS80, # NAD83,
2430 Ellipsoids.Sphere, Ellipsoids.SpherePopular,
2431 Ellipsoid(_EWGS84.b, _EWGS84.a, name='_Prolate')):
2432 e = f2n(E.f) - E.n
2433 printf('# %s: %s', _DOT_('Ellipsoids', E.name), E.toStr(prec=10, terse=0), nl=1)
2434 printf('# e=%s, f_=%s, f=%s, n=%s (%s)', fstr(E.e, prec=13, fmt=Fmt.e),
2435 fstr(E.f_, prec=13, fmt=Fmt.e),
2436 fstr(E.f, prec=13, fmt=Fmt.e),
2437 fstr(E.n, prec=13, fmt=Fmt.e),
2438 fstr(e, prec=9, fmt=Fmt.e))
2439 printf('# %s %s', Ellipsoid.AlphaKs.name, fstr(E.AlphaKs, prec=20))
2440 printf('# %s %s', Ellipsoid.BetaKs.name, fstr(E.BetaKs, prec=20))
2441 printf('# %s %s', nameof(Ellipsoid.KsOrder), E.KsOrder) # property
2443 from pygeodesy.internals import _pregistry
2444 # __doc__ of this file, force all into registry
2445 _pregistry(Ellipsoids)
2447# % python3.13 -m pygeodesy.ellipsoids
2449# Ellipsoids.WGS84: name='WGS84', a=6378137, f=0.0033528107, f_=298.257223563, b=6356752.3142451793, f2=0.0033640898, n=0.0016792204, e=0.0818191908, e2=0.00669438, e21=0.99330562, e22=0.0067394967, e32=0.0033584313, A=6367449.1458234144, L=10001965.7293127235, R1=6371008.7714150595, R2=6371007.1809184738, R3=6371000.7900091587, Rbiaxial=6367453.6345163295, Rtriaxial=6372797.5559594007
2450# e=8.1819190842622e-02, f_=2.98257223563e+02, f=3.3528106647475e-03, n=1.6792203863837e-03 (0.0e+00)
2451# AlphaKs 0.00083773182062446994, 0.00000076085277735725, 0.00000000119764550324, 0.00000000000242917068, 0.00000000000000571182, 0.0000000000000000148, 0.00000000000000000004, 0.0
2452# BetaKs 0.00083773216405794875, 0.0000000590587015222, 0.00000000016734826653, 0.00000000000021647981, 0.00000000000000037879, 0.00000000000000000072, 0.0, 0.0
2453# KsOrder 8
2455# Ellipsoids.GRS80: name='GRS80', a=6378137, f=0.0033528107, f_=298.2572221009, b=6356752.3141403468, f2=0.0033640898, n=0.0016792204, e=0.081819191, e2=0.00669438, e21=0.99330562, e22=0.0067394968, e32=0.0033584313, A=6367449.1457710434, L=10001965.7292304561, R1=6371008.7713801153, R2=6371007.1808835147, R3=6371000.7899741363, Rbiaxial=6367453.6344640013, Rtriaxial=6372797.5559332585
2456# e=8.1819191042833e-02, f_=2.9825722210088e+02, f=3.3528106811837e-03, n=1.6792203946295e-03 (0.0e+00)
2457# AlphaKs 0.00083773182472890429, 0.00000076085278481561, 0.00000000119764552086, 0.00000000000242917073, 0.00000000000000571182, 0.0000000000000000148, 0.00000000000000000004, 0.0
2458# BetaKs 0.0008377321681623882, 0.00000005905870210374, 0.000000000167348269, 0.00000000000021647982, 0.00000000000000037879, 0.00000000000000000072, 0.0, 0.0
2459# KsOrder 8
2461# Ellipsoids.Sphere: name='Sphere', a=6371008.7714149999, f=0, f_=0, b=6371008.7714149999, f2=0, n=0, e=0, e2=0, e21=1, e22=0, e32=0, A=6371008.7714149999, L=10007557.1761167478, R1=6371008.7714149999, R2=6371008.7714149999, R3=6371008.7714149999, Rbiaxial=6371008.7714149999, Rtriaxial=6371008.7714149999
2462# e=0.0e+00, f_=0.0e+00, f=0.0e+00, n=0.0e+00 (0.0e+00)
2463# AlphaKs 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0
2464# BetaKs 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0
2465# KsOrder 8
2467# Ellipsoids.SpherePopular: name='SpherePopular', a=6378137, f=0, f_=0, b=6378137, f2=0, n=0, e=0, e2=0, e21=1, e22=0, e32=0, A=6378137, L=10018754.171394622, R1=6378137, R2=6378137, R3=6378137, Rbiaxial=6378137, Rtriaxial=6378137
2468# e=0.0e+00, f_=0.0e+00, f=0.0e+00, n=0.0e+00 (0.0e+00)
2469# AlphaKs 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0
2470# BetaKs 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0
2471# KsOrder 8
2473# Ellipsoids._Prolate: name='_Prolate', a=6356752.3142451793, f=-0.0033640898, f_=-297.257223563, b=6378137, f2=-0.0033528107, n=-0.0016792204, e=0.0820944379, e2=-0.0067394967, e21=1.0067394967, e22=-0.00669438, e32=-0.0033584313, A=6367449.1458234144, L=10035500.5204500332, R1=6363880.5428301189, R2=6363878.9413582645, R3=6363872.5644020075, Rbiaxial=6367453.6345163295, Rtriaxial=6362105.2243882557
2474# e=8.2094437949696e-02, f_=-2.97257223563e+02, f=-3.3640898209765e-03, n=-1.6792203863837e-03 (0.0e+00)
2475# AlphaKs -0.00084149152514366627, 0.00000076653480614871, -0.00000000120934503389, 0.0000000000024576225, -0.00000000000000578863, 0.00000000000000001502, -0.00000000000000000004, 0.0
2476# BetaKs -0.00084149187224351817, 0.00000005842735196773, -0.0000000001680487236, 0.00000000000021706261, -0.00000000000000038002, 0.00000000000000000073, -0.0, 0.0
2477# KsOrder 8
2479# **) MIT License
2480#
2481# Copyright (C) 2016-2026 -- mrJean1 at Gmail -- All Rights Reserved.
2482#
2483# Permission is hereby granted, free of charge, to any person obtaining a
2484# copy of this software and associated documentation files (the "Software"),
2485# to deal in the Software without restriction, including without limitation
2486# the rights to use, copy, modify, merge, publish, distribute, sublicense,
2487# and/or sell copies of the Software, and to permit persons to whom the
2488# Software is furnished to do so, subject to the following conditions:
2489#
2490# The above copyright notice and this permission notice shall be included
2491# in all copies or substantial portions of the Software.
2492#
2493# THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
2494# OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
2495# FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
2496# THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR
2497# OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE,
2498# ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR
2499# OTHER DEALINGS IN THE SOFTWARE.