GeographicLib 2.6
Loading...
Searching...
No Matches
GeodesicLine3.hpp
Go to the documentation of this file.
1/**
2 * \file GeodesicLine3.hpp
3 * \brief Header for GeographicLib::Triaxial::GeodesicLine3 class
4 *
5 * Copyright (c) Charles Karney (2024-2025) <karney@alum.mit.edu> and licensed
6 * under the MIT/X11 License. For more information, see
7 * https://geographiclib.sourceforge.io/
8 **********************************************************************/
9
10#if !defined(GEOGRAPHICLIB_GEODESICLINE3_HPP)
11#define GEOGRAPHICLIB_GEODESICLINE3_HPP 1
12
13#include <utility>
14#include <vector>
20
21#if defined(_MSC_VER)
22// Squelch warnings about dll vs fline + fline::fics + gline + gline::gics
23# pragma warning (push)
24# pragma warning (disable: 4251)
25#endif
26
27namespace GeographicLib {
28 namespace Triaxial {
29
30 /**
31 * \brief The direct geodesic problem for a triaxial ellipsoid
32 *
33 * This is an implementation of
34 * <a href="https://books.google.com/books?id=RbwGAAAAYAAJ&pg=PA309">
35 * Jacobi's method for finding geodesics on a triaxial ellipsoid</a> (1839).
36
37 * This class performs the quadrature necessary to evaluate to the integrals
38 * for the course and distance equations and it solves the 2 coupled
39 * equations to determine the latitude and longitude of the destination.
40 * This functionality is used by the Geodesic3 class to solve the inverse
41 * geodesic problem.
42 *
43 * Example of use:
44 * \include example-GeodesicLine3.cpp
45 *
46 * <a href="Geod3Solve.1.html">Geod3Solve</a> is a command-line utility
47 * providing access to the functionality of Geodesic3 and GeodesicLine3.
48 *
49 * The class experimental::TriaxialGeodesicODE solves the direct geodesic
50 * problem by integrating the corresponing ordinary differential equations.
51 * This is \e not part of %GeographicLib itself because it used Boost for
52 * solving the ODEs and Boost is not a requirement for using %GeographicLib.
53 **********************************************************************/
54 class GEOGRAPHICLIB_EXPORT GeodesicLine3 {
55 private:
56 /// \cond SKIP
57 friend class Geodesic3; // For access to fline, gline, etc.
58 /// \endcond
59 using real = Math::real;
60 using ang = Angle;
61 static constexpr int maxit_ = 300;
62
63 class hfun {
64 // This combines ffun abd gfun in order to minimize the duplication of
65 // code.
66 // Establish consistent notation for coordinates:
67 // theta = rotating coord, omega-90 or beta for circum- or transpolar
68 // phi = librating coord, beta or omega-90 for circum- or transpolar
69 // psi = rotating equivalent of phi
70 // zeta = a generic version of theta or psi
71 // u = possible change of vars for theta
72 // v = possible change of vars for psi
73 // w = a generic version of u or v
74 private:
75 real _kap, _kapp, _eps, _mu, _sqrtmu, _sqrtkap, _sqrtkapp;
76 bool _distp, _tx;
78 TrigfunExt _fun;
79 real _max;
80 bool _umb, _meridr, _meridl, _biaxr, _biaxl;
81 // The f functions
82 // mu > 0
83 static real fthtp(real c, real kap, real kapp, real eps, real mu);
84 static real fup(real cn, real kap, real kapp, real eps, real mu);
85 // mu == 0
86 static real dfp(real c, real kap, real kapp, real eps);
87 static real dfvp(real cn, real dn, real kap, real kapp, real eps);
88 // mu < 0
89 static real fpsip(real s, real c, real kap, real kapp,
90 real eps, real mu);
91 static real fvp(real dn, real kap, real kapp, real eps, real mu);
92 // biaxial variant for kap = 0, kapp = 1
93 static real fthtbiax(real tht, real eps, real mu);
94 // biaxial variant for kap = 1, kapp = 0, mu <= 0, !_tx
95 static real dfpsibiax(real s, real c, real eps, real mu);
96 // NOT USED
97 // biaxial variant for kap = 1, kapp = 0, mu <= 0, _tx
98 // static real dfvbiax(real dn, real eps, real mu);
99
100 // The g functions
101 // _mu > 0
102 static real gthtp(real c, real kap, real kapp, real eps, real mu);
103 static real gup(real cn, real dn, real kap, real kapp,
104 real eps, real mu);
105 // _mu == 0
106 static real g0p(real c, real kap, real kapp, real eps);
107 static real g0vp(real cn, real kap, real kapp, real eps);
108 // _mu < 0
109 static real gpsip(real s, real c, real kap, real kapp,
110 real eps, real mu);
111 static real gvp(real cn, real dn, real kap, real kapp,
112 real eps, real mu);
113 // biaxial variants for kap = 0, kapp = 1, mu >= 0
114 static real gthtbiax(real tht, real eps, real mu);
115 // biaxial variants for kap = 1, kapp = 0, mu <= 0
116 static real gpsibiax(real s, real c, real eps, real mu);
117 // NOT USED
118 // static real gvbiax(real cn, real dn, real eps, real mu);
119
120 // Return atan(m * tan(x)) keeping result continuous. Only defined for
121 // !signbit(m).
122 static real modang(real x, real m) {
123 return ang::radians(x).modang(m).radians();
124 }
125 real root(real z, real u0, int* countn, int* countb, real tol = 0) const;
126 public:
127 // Summary of f and g functions
128 // _meridr = kap == 0 && mu == 0 (biaxial meridian rotating coordinate)
129 // f = fthtbiax / sqrt(mu), fthtbiax = 1 inversion trivial
130 // g = gthtbiax, gthtbiax = 0
131 // _meridl = kapp == 0 && mu -= 0 (biaxial meridian librating coordinate)
132 // f = (atan(sqrt(-mu)*tan(psi)) - dfpsiobl) / sqrt(-mu)
133 // mu == 0 **
134 // atan(sqrt(-mu)*tan(psi)) -> round(psi/pi)*pi
135 // g = gpsibiax
136 // mu > 0 (rotating coordinate)
137 // f = ftht or fu
138 // g = gtht or gu **
139 // mu < 0 (librating coordinate)
140 // f = fpsi or fv
141 // g = gpsi or gv **
142 // _umb = !biaxial && mu == 0 (triaxial umbilic)
143 // f = (u + df)/sqrt(kap*kapp) or (u + dfv)/sqrt(kap*kapp)
144 // g = g0 or g0v
145 //
146 // Handing of f funtions needs to be handled specially for
147 // _merid[lr]
148 // leading order behavior is seperated out
149 // N.B. inverse of f is discontinuous for mu == 0
150 // functions which need special treatment
151 // operator(), deriv, inv, Slope, Max
152 // _umb
153 // leading order behavior is seperated out
154 // otherwise everything is handled by generic routines
155 //
156 // Handling of g function is always generic
157 // _umb
158 // functions which need special treatment
159 // Slope, fwd, rev, Max
160 //
161 hfun() {}
162 hfun(bool distp, real kap, real kapp, real eps, real mu,
163 const Geodesic3& tg);
164 real operator()(real u) const;
165 // THIS ISN"T USED
166 // ang operator()(const ang& ang) const;
167 real deriv(real u) const;
168 real df(real u) const { return _fun(u); }
169 real dfp(real u) const { return _fun.deriv(u); }
170
171 // Accurate (to tol) inverse by direct Newton (not using _finv)
172 real inv(real z, int* countn = nullptr, int* countb = nullptr) const;
173
174 real fwd(real zeta) const {
175 return _umb ? lam(zeta, _sqrtkapp) :
176 (_tx ? _ell.F(zeta) : zeta);
177 }
178 real rev(real w) const {
179 return _umb ? gd(w, _sqrtkapp) :
180 (_tx ? _ell.am(w) : w);
181 }
182 int NCoeffs() const { return _fun.NCoeffs(); }
183 bool txp() const { return _tx; }
184 real HalfPeriod() const {
185 return _umb ? Math::infinity() : (_tx ? _ell.K() : Math::pi()/2);
186 }
187 real Slope() const {
188 using std::sqrt;
189 return _umb ? 1 :
190 !_distp && _meridl ? 0 :
191 !_distp && _biaxl ? 1 - sqrt(-_mu) * _fun.Slope() :
192 _fun.Slope();
193 }
194 real Max() const { return _max; }
195 real MaxPlus() const {
196 using std::fmax;
197 return fmax(_max, HalfPeriod() * (Slope() == 0 ? 1 : Slope()) / 1000);
198 }
199 };
200
201 class fline {
202 private:
203 Geodesic3 _tg;
204 Geodesic3::gamblk _gm;
205 real _deltashift;
206 hfun _fpsi, _ftht;
207 public:
208 class fics {
209 // bundle of data setting the initial conditions for a geodesic
210 public:
211 // bool transpolar;
212 // alp1 is angle measured from line of const rotating coording
213 ang tht1, phi1, alp1, // rotating, librating starting point
214 psi1, // phi1 transformed to rotating angle psi
215 // Angles about which quantities oscillate
216 // circumpolar:
217 // omg0 not used, bet0 = cardinal(even), alp0 = cardinal(odd)
218 // transpolar:
219 // bet0 not used, omg0 = cardinal(even), alp0 = cardinal(even)
220 // umbilical
221 // bet0, omg0 = cardinal(even) = middle of starting segment
222 // !umbalt: alp0 = cardinal(odd)
223 // umbalt: alp0 = cardinal(even)
225 // Angle versions of u0, v0, delta, defer for now
226 // ang u0a, v0a, deltaa;
227 real u0, v0, delta; // starting point geodesic
228 int Nx, Ex; // Northgoing / eastgoing relative to tht
229 fics() {}
230 fics(const fline& f,
231 ang bet1, ang omg1, ang alp1);
232 void setquadrant(const fline& f, unsigned q);
233 void pos1(bool transpolar,
234 ang& bet10, ang& omg10, ang& alp10) const;
235 };
236 class disttx {
237 // bundle of data to pass along for distance
238 public:
240 int ind2;
241 };
242 fline() {}
243 fline(const Geodesic3& tg, bool neg = false);
244 fline(const Geodesic3& tg, Geodesic3::gamblk gm);
245 const hfun& fpsi() const { return _fpsi; }
246 const hfun& ftht() const { return _ftht; }
247 const hfun& fbet() const {
248 return !transpolar() ? fpsi() : ftht();
249 }
250 const hfun& fomg() const {
251 return transpolar() ? fpsi() : ftht();
252 }
253 const Geodesic3& tg() const { return _tg; }
254 const Ellipsoid3& t() const { return _tg.t(); }
255 real gamma() const { return _gm.gamma; }
256 real gammax() const { return _gm.gammax; }
257 real kx2() const { return _gm.kx2; }
258 real kxp2() const { return _gm.kxp2; }
259 real kx() const { return _gm.kx; }
260 real kxp() const { return _gm.kxp; }
261 real nu() const { return _gm.nu; }
262 real nup() const { return _gm.nup; }
263 real deltashift() const { return _deltashift; }
264 bool transpolar() const { return _gm.transpolar; }
265 const Geodesic3::gamblk& gm() const { return _gm; }
266 // Run fline to its first intersection with
267 // (for betp) bet2 and return omg2
268 // (for !betp) omg2 and return bet2
269 real Hybrid0(const fics& fic, ang bet2, ang omg2,
270 bool betp = true) const;
271 // Run fline to its first intersection with bet and return resulting
272 // bet2a, omg2a, alp2a (without angle normalization) and distance
273 // calculation object
274 disttx Hybrid(const fics& fic, ang betomg2,
275 ang& bet2a, ang& omg2a, ang& alp2a,
276 bool betp = true) const;
277 disttx ArcPos0(const fics& fic, ang tau12,
278 ang& bet2a, ang& omg2a, ang& alp2a,
279 bool betp = true) const;
280 };
281
282 class gline {
283 private:
284 Geodesic3 _tg;
285 Geodesic3::gamblk _gm;
286 hfun _gpsi, _gtht;
287 public:
288 real s0;
289 class gics {
290 // bundle of data setting the initial conditions for a distance calc
291 public:
292 real sig1, s13; // starting point
293 gics() {}
294 gics(const gline& g, const fline::fics& fic);
295 };
296 gline() {}
297 gline(const Geodesic3& tg, bool neg = false);
298 gline(const Geodesic3& tg, const Geodesic3::gamblk& gm);
299 real gamma() const { return _gm.gamma; }
300 real gammax() const { return _gm.gammax; }
301 real kx2() const { return _gm.kx2; }
302 real kxp2() const { return _gm.kxp2; }
303 real kx() const { return _gm.kx; }
304 real kxp() const { return _gm.kxp; }
305 bool transpolar() const { return _gm.transpolar; }
306 const hfun& gpsi() const { return _gpsi; }
307 const hfun& gtht() const { return _gtht; }
308 const hfun& gbet() const {
309 return !transpolar() ? gpsi() : gtht();
310 }
311 const hfun& gomg() const {
312 return transpolar() ? gpsi() : gtht();
313 }
314 const Geodesic3& tg() const { return _tg; }
315 const Ellipsoid3& t() const { return _tg.t(); }
316 const Geodesic3::gamblk& gm() const { return _gm; }
317 real dist(gics ic, fline::disttx d) const;
318 };
319
320 class zvals {
321 public:
322 real z, fz, gz;
323 zvals(real z0 = 0, real fz0 = 0, real gz0 = 0)
324 : z(z0), fz(fz0), gz(gz0) {}
325 bool operator<(const zvals& t) const { return z < t.z; }
326 bool operator==(const zvals& t) const { return z == t.z; }
327 };
328
329 class zset {
330 private:
331 std::vector<zvals> _s;
332 public:
333 zset(const zvals& a, const zvals& b)
334 : _s({a, b})
335 {
336 if (a == b)
337 // Allow coincident start and end values
338 _s.resize(1);
339 else if (!(a < b && a.fz <= b.fz && a.gz <= b.gz))
340 throw GeographicLib::GeographicErr("bad zset initializer");
341 }
342 int num() const { return int(_s.size()); }
343 const zvals& val(int i) const { return _s[i]; }
344 const zvals& min() const { return _s[0]; }
345 const zvals& max() const { return _s.back(); }
346 int insert(zvals& t, int flag = 0);
347 real bisect() const {
348 // return (min().z + max().z) / 2;
349 // return z in the middle of biggest gap
350 if (num() == 1)
351 return min().z;
352 real maxgap = -1; int maxind = 0;
353 for (int i = 0; i < num() - 1; ++i) {
354 real gap = _s[i+1].z - _s[i].z;
355 if (gap > maxgap) {
356 maxgap = gap;
357 maxind = i;
358 }
359 }
360 return (_s[maxind].z + _s[maxind+1].z) / 2;
361 }
362 };
363
364 Geodesic3 _tg;
365 fline _f;
366 fline::fics _fic;
367 gline _g;
368 gline::gics _gic;
369
370 static void solve2(real f0, real g0,
371 const hfun& fx, const hfun& fy,
372 const hfun& gx, const hfun& gy,
373 real& x, real& y,
374 int* countn = nullptr, int* countb = nullptr);
375 static void solve2u(real f0, real g0,
376 const hfun& fx, const hfun& fy,
377 const hfun& gx, const hfun& gy,
378 real& x, real& y,
379 int* countn = nullptr, int* countb = nullptr);
380 static void newt2(real f0, real g0,
381 const hfun& fx, const hfun& fy,
382 const hfun& gx, const hfun& gy,
383 real xa, real xb, real xscale,
384 real ya, real yb, real yscale,
385 real fscale, real gscale,
386 real& x, real& y,
387 int* countn = nullptr, int* countb = nullptr);
388 static void zsetsinsert(zset& xset, zset& yset,
389 zvals& xfg, zvals& yfg,
390 real f0, real g0);
391 static void zsetsdiag(const zset& xset, const zset& yset,
392 real f0, real g0);
393 static std::pair<real, real> zsetsbisect(const zset& xset, const zset& yset,
394 real f0, real g0, bool secant);
395 static real bigclamp(real x, real mult = 1) {
396 using std::fmax, std::fmin;
397 real z = mult * Geodesic3::BigValue();
398 return Math::clamp(x, -z, z);
399 }
400 static real lamang0(ang x, real mult = 1) {
401 // lam(x) when x is an ang -- no clamping
402 using std::asinh, std::fabs;
403 return asinh(mult * x.t());
404 }
405 static real lamang(ang x, real mult = 1) {
406 // lam(x) when x is an ang -- with clamping
407 // A consistent large value for x near pi/2.
408 return bigclamp(lamang0(x, mult));
409 }
410 static real lam(real x, real mult = 1) {
411 using std::tan, std::asinh, std::fabs, std::copysign;
412 // A consistent large value for x near pi/2. Also deals with the issue
413 // that tan(pi/2) may be negative, e.g., for long doubles.
414 return fabs(x) >= Math::pi()/2 ? copysign(Geodesic3::BigValue(), x) :
415 asinh(mult * tan(x));
416 }
417 static real gd(real x, real mult = 1) {
418 using std::atan, std::sinh;
419 return atan(sinh(x) / mult);
420 }
421 static ang anglam(real u, real mult = 1) {
422 using std::sinh;
423 return ang(sinh(u), mult, 0);
424 }
425 static real mcosh(real u, real mult = 1) {
426 using std::cosh, std::sinh, std::hypot;
427 return mult == 1 ? cosh(u) : hypot(sinh(u), mult) / mult;
428 }
429
430 const hfun& fbet() const { return _f.fbet(); }
431 const hfun& fomg() const { return _f.fomg(); }
432 const hfun& fpsi() const { return _f.fpsi(); }
433 const hfun& ftht() const { return _f.ftht(); }
434 const hfun& gbet() const { return _g.gbet(); }
435 const hfun& gomg() const { return _g.gomg(); }
436 const hfun& gpsi() const { return _g.gpsi(); }
437 const hfun& gtht() const { return _g.gtht(); }
438
439 // remainder with result in
440 // [-y/2, y/2) if alt = false (default)
441 // (-y/2, y/2] if alt = true
442 static std::pair<real, real> remx(real x, real y, bool alt = false) {
443 using std::remainder, std::rint;
444 real z = remainder(x, y);
445 if (alt) {
446 if (z == -y/2) z = y/2;
447 } else {
448 if (z == y/2) z = -y/2;
449 }
450 return std::pair<real, real>(z, rint((x - z) / y));
451 }
452 // remainder with x in
453 // [-pi/2, pi/2) if alt = false (default)
454 // (-pi/2, pi/2] if alt = true
455 // equivalent to remx(x.radians(), Math::pi(), alt)
456 static std::pair<real, real> remx(ang x, bool alt = false) {
457 using std::signbit;
458 real m = 0;
459 if (signbit(x.c())) {
460 x -= ang::cardinal(2);
461 ++m;
462 }
463 if (x.c() == 0) {
464 if (alt && signbit(x.s())) {
465 x += ang::cardinal(2);
466 --m;
467 } else if (!alt && !signbit(x.s())) {
468 x -= ang::cardinal(2);
469 ++m;
470 }
471 }
472 return std::pair<real, real>(x.radians0(), m + 2*x.n());
473 }
474 static bool biaxspecial(const Geodesic3& tg, real gamma) {
475 using std::fabs;
476 return Geodesic3::biaxp_ && (tg.t().k2() == 0 || tg.t().kp2() == 0) &&
477 gamma != 0 && fabs(gamma) < tg._ellipthresh;
478 }
479 // Private constructor to assemble the pieces of the class on exiting
480 // Geodesic3::Inverse.
481 GeodesicLine3(fline f, fline::fics fic, gline g, gline::gics gic);
482 // Private constructor to provide the umbilical fline object.
483 GeodesicLine3(const Geodesic3& tg);
484 // Find first crossing of bet = bet2. NaN returned if
485 // no crossing. Assume bet1 in [-90,0] (or maybe (-90,0]) and bet2 in
486 // [bet2, -bet2]. Special cases, bet2 = bet1...
487 // bet1 < 0, alp1 in [-90,90], omg2 = 0
488 // bet1 == 0, alp1 in (-90,90), omg2 = 0,
489 // alp1 = +/-90 omg2 = conj pt
490 void Hybrid(Angle betomg2,
491 Angle& bet2a, Angle& omg2a, Angle& alp2a,
492 real& s12, bool betp = true) const;
493 public:
494 /**
495 * Constructor for a geodesic line specified at point 1
496 *
497 * @param[in] tg the underlying Geodesic3 object.
498 * @param[in] bet1 the ellipsoidal latitude of point 1.
499 * @param[in] omg1 the ellipsoidal longitude of point 1.
500 * @param[in] alp1 the forward azimuth of the geodesic at point 1.
501 **********************************************************************/
502 GeodesicLine3(const Geodesic3& tg, Angle bet1, Angle omg1, Angle alp1);
503 /**
504 * Constructor for a geodesic line specified at point 1 in degrees
505 *
506 * @param[in] tg the underlying Geodesic3 object.
507 * @param[in] bet1 the ellipsoidal latitude of point 1.
508 * @param[in] omg1 the ellipsoidal longitude of point 1 (in degrees).
509 * @param[in] alp1 the forward azimuth of the geodesic at point 1 (in
510 * degrees).
511 **********************************************************************/
512 GeodesicLine3(const Geodesic3& tg, real bet1, real omg1, real alp1);
513 /**
514 * Find point 2 a given distance from point 1
515 *
516 * @param[in] s12 the distance from point 1 to point 2.
517 * @param[out] bet2 the ellipsoidal latitude of point 2.
518 * @param[out] omg2 the ellipsoidal longitude of point 2.
519 * @param[out] alp2 the forward azimuth of the geodesic at point 2.
520 **********************************************************************/
521 void Position(real s12, Angle& bet2, Angle& omg2, Angle& alp2) const;
522 /**
523 * Find point 2 a given distance from point 1 in degrees
524 *
525 * @param[in] s12 the distance from point 1 to point 2.
526 * @param[out] bet2 the ellipsoidal latitude of point 2 (in degrees).
527 * @param[out] omg2 the ellipsoidal longitude of point 2 (in degrees).
528 * @param[out] alp2 the forward azimuth of the geodesic at point 2 (in
529 * degrees).
530 * @param[in] unroll if true (the default) "unroll" the coordinates;
531 * otherwise reduce them to their conventional ranges.
532 **********************************************************************/
533 void Position(real s12, real& bet2, real& omg2, real& alp2,
534 bool unroll = true) const;
535 /**
536 * Define a reference point 3 on the geodesic line
537 *
538 * @param[in] s13 distance from point 1 to point 3.
539 **********************************************************************/
540 void SetDistance(real s13) { _gic.s13 = s13; }
541 /**
542 * @return \e s13 the distance from point 1 to reference point 3.
543 **********************************************************************/
544 real Distance() const { return _gic.s13; }
545 /**
546 * Return the coordinates of point 1
547 *
548 * @param[out] bet1 the ellipsoidal latitude of point 1.
549 * @param[out] omg1 the ellipsoidal longitude of point 1.
550 * @param[out] alp1 the forward azimuth of the geodesic at point 1.
551 **********************************************************************/
552 void pos1(Angle& bet1, Angle& omg1, Angle& alp1) const;
553 /**
554 * Return the coordinates of point 1 in degrees
555 *
556 * @param[out] bet1 the ellipsoidal latitude of point 1 (in degrees).
557 * @param[out] omg1 the ellipsoidal longitude of point 1 (in degrees).
558 * @param[out] alp1 the forward azimuth of the geodesic at point 1 (in
559 * degrees).
560 * @param[in] unroll if true (the default), return the coordinates used to
561 * specify this line; otherwise reduce the coordinates for point 1 to
562 * their conventional ranges.
563 **********************************************************************/
564 void pos1(real& bet1, real& omg1, real& alp1, bool unroll = true) const;
565 };
566
567 } // namespace Triaxial
568} // namespace GeographicLib
569
570#if defined(_MSC_VER)
571# pragma warning (pop)
572#endif
573
574#endif // GEOGRAPHICLIB_GEODESICLINE3_HPP
Header for the GeographicLib::AngleT class.
Header for GeographicLib::Constants class.
#define GEOGRAPHICLIB_EXPORT
Definition Constants.hpp:59
Header for GeographicLib::EllipticFunction class.
GeographicLib::Angle ang
GeographicLib::Math::real real
Header for GeographicLib::Triaxial::Geodesic3 class.
Header for GeographicLib::Trigfun class.
static AngleT cardinal(Math::real q)
Elliptic integrals and functions.
Math::real F(real phi) const
static T infinity()
Definition Math.cpp:308
static T pi()
Definition Math.hpp:187
The solution of the geodesic problem for a triaxial ellipsoid.
Definition Geodesic3.hpp:65
const Ellipsoid3 & t() const
void pos1(Angle &bet1, Angle &omg1, Angle &alp1) const
A function defined by its derivative and its inverse.
Definition Trigfun.hpp:498
real deriv(real x) const
Definition Trigfun.hpp:564
Namespace for operations on triaxial ellipsoids.
Namespace for GeographicLib.
AngleT< Math::real > Angle
Definition Angle.hpp:760