NormalGravity.hpp 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400
  1. /**
  2. * \file NormalGravity.hpp
  3. * \brief Header for GeographicLib::NormalGravity class
  4. *
  5. * Copyright (c) Charles Karney (2011-2020) <charles@karney.com> and licensed
  6. * under the MIT/X11 License. For more information, see
  7. * https://geographiclib.sourceforge.io/
  8. **********************************************************************/
  9. #if !defined(GEOGRAPHICLIB_NORMALGRAVITY_HPP)
  10. #define GEOGRAPHICLIB_NORMALGRAVITY_HPP 1
  11. #include <GeographicLib/Constants.hpp>
  12. #include <GeographicLib/Geocentric.hpp>
  13. namespace GeographicLib {
  14. /**
  15. * \brief The normal gravity of the earth
  16. *
  17. * "Normal" gravity refers to an idealization of the earth which is modeled
  18. * as an rotating ellipsoid. The eccentricity of the ellipsoid, the rotation
  19. * speed, and the distribution of mass within the ellipsoid are such that the
  20. * ellipsoid is a "level ellipoid", a surface of constant potential
  21. * (gravitational plus centrifugal). The acceleration due to gravity is
  22. * therefore perpendicular to the surface of the ellipsoid.
  23. *
  24. * Because the distribution of mass within the ellipsoid is unspecified, only
  25. * the potential exterior to the ellipsoid is well defined. In this class,
  26. * the mass is assumed to be to concentrated on a "focal disc" of radius,
  27. * (<i>a</i><sup>2</sup> &minus; <i>b</i><sup>2</sup>)<sup>1/2</sup>, where
  28. * \e a is the equatorial radius of the ellipsoid and \e b is its polar
  29. * semi-axis. In the case of an oblate ellipsoid, the mass is concentrated
  30. * on a "focal rod" of length 2(<i>b</i><sup>2</sup> &minus;
  31. * <i>a</i><sup>2</sup>)<sup>1/2</sup>. As a result the potential is well
  32. * defined everywhere.
  33. *
  34. * There is a closed solution to this problem which is implemented here.
  35. * Series "approximations" are only used to evaluate certain combinations of
  36. * elementary functions where use of the closed expression results in a loss
  37. * of accuracy for small arguments due to cancellation of the leading terms.
  38. * However these series include sufficient terms to give full machine
  39. * precision.
  40. *
  41. * Although the formulation used in this class applies to ellipsoids with
  42. * arbitrary flattening, in practice, its use should be limited to about
  43. * <i>b</i>/\e a &isin; [0.01, 100] or \e f &isin; [&minus;99, 0.99].
  44. *
  45. * Definitions:
  46. * - <i>V</i><sub>0</sub>, the gravitational contribution to the normal
  47. * potential;
  48. * - &Phi;, the rotational contribution to the normal potential;
  49. * - \e U = <i>V</i><sub>0</sub> + &Phi;, the total potential;
  50. * - <b>&Gamma;</b> = &nabla;<i>V</i><sub>0</sub>, the acceleration due to
  51. * mass of the earth;
  52. * - <b>f</b> = &nabla;&Phi;, the centrifugal acceleration;
  53. * - <b>&gamma;</b> = &nabla;\e U = <b>&Gamma;</b> + <b>f</b>, the normal
  54. * acceleration;
  55. * - \e X, \e Y, \e Z, geocentric coordinates;
  56. * - \e x, \e y, \e z, local cartesian coordinates used to denote the east,
  57. * north and up directions.
  58. *
  59. * References:
  60. * - C. Somigliana, Teoria generale del campo gravitazionale dell'ellissoide
  61. * di rotazione, Mem. Soc. Astron. Ital, <b>4</b>, 541--599 (1929).
  62. * - W. A. Heiskanen and H. Moritz, Physical Geodesy (Freeman, San
  63. * Francisco, 1967), Secs. 1-19, 2-7, 2-8 (2-9, 2-10), 6-2 (6-3).
  64. * - B. Hofmann-Wellenhof, H. Moritz, Physical Geodesy (Second edition,
  65. * Springer, 2006) https://doi.org/10.1007/978-3-211-33545-1
  66. * - H. Moritz, Geodetic Reference System 1980, J. Geodesy 54(3), 395-405
  67. * (1980) https://doi.org/10.1007/BF02521480
  68. *
  69. * For more information on normal gravity see \ref normalgravity.
  70. *
  71. * Example of use:
  72. * \include example-NormalGravity.cpp
  73. **********************************************************************/
  74. class GEOGRAPHICLIB_EXPORT NormalGravity {
  75. private:
  76. static const int maxit_ = 20;
  77. typedef Math::real real;
  78. friend class GravityModel;
  79. real _a, _GM, _omega, _f, _J2, _omega2, _aomega2;
  80. real _e2, _ep2, _b, _E, _U0, _gammae, _gammap, _Q0, _k, _fstar;
  81. Geocentric _earth;
  82. static real atanzz(real x, bool alt) {
  83. // This routine obeys the identity
  84. // atanzz(x, alt) = atanzz(-x/(1+x), !alt)
  85. //
  86. // Require x >= -1. Best to call with alt, s.t. x >= 0; this results in
  87. // a call to atan, instead of asin, or to asinh, instead of atanh.
  88. using std::sqrt; using std::abs; using std::atan; using std::asin;
  89. using std::asinh; using std::atanh;
  90. real z = sqrt(abs(x));
  91. return x == 0 ? 1 :
  92. (alt ?
  93. (!(x < 0) ? asinh(z) : asin(z)) / sqrt(abs(x) / (1 + x)) :
  94. (!(x < 0) ? atan(z) : atanh(z)) / z);
  95. }
  96. static real atan7series(real x);
  97. static real atan5series(real x);
  98. static real Qf(real x, bool alt);
  99. static real Hf(real x, bool alt);
  100. static real QH3f(real x, bool alt);
  101. real Jn(int n) const;
  102. void Initialize(real a, real GM, real omega, real f_J2, bool geometricp);
  103. public:
  104. /** \name Setting up the normal gravity
  105. **********************************************************************/
  106. ///@{
  107. /**
  108. * Constructor for the normal gravity.
  109. *
  110. * @param[in] a equatorial radius (meters).
  111. * @param[in] GM mass constant of the ellipsoid
  112. * (meters<sup>3</sup>/seconds<sup>2</sup>); this is the product of \e G
  113. * the gravitational constant and \e M the mass of the earth (usually
  114. * including the mass of the earth's atmosphere).
  115. * @param[in] omega the angular velocity (rad s<sup>&minus;1</sup>).
  116. * @param[in] f_J2 either the flattening of the ellipsoid \e f or the
  117. * the dynamical form factor \e J2.
  118. * @param[out] geometricp if true (the default), then \e f_J2 denotes the
  119. * flattening, else it denotes the dynamical form factor \e J2.
  120. * @exception if \e a is not positive or if the other parameters do not
  121. * obey the restrictions given below.
  122. *
  123. * The shape of the ellipsoid can be given in one of two ways:
  124. * - geometrically (\e geomtricp = true), the ellipsoid is defined by the
  125. * flattening \e f = (\e a &minus; \e b) / \e a, where \e a and \e b are
  126. * the equatorial radius and the polar semi-axis. The parameters should
  127. * obey \e a &gt; 0, \e f &lt; 1. There are no restrictions on \e GM or
  128. * \e omega, in particular, \e GM need not be positive.
  129. * - physically (\e geometricp = false), the ellipsoid is defined by the
  130. * dynamical form factor <i>J</i><sub>2</sub> = (\e C &minus; \e A) /
  131. * <i>Ma</i><sup>2</sup>, where \e A and \e C are the equatorial and
  132. * polar moments of inertia and \e M is the mass of the earth. The
  133. * parameters should obey \e a &gt; 0, \e GM &gt; 0 and \e J2 &lt; 1/3
  134. * &minus; (<i>omega</i><sup>2</sup><i>a</i><sup>3</sup>/<i>GM</i>)
  135. * 8/(45&pi;). There is no restriction on \e omega.
  136. **********************************************************************/
  137. NormalGravity(real a, real GM, real omega, real f_J2,
  138. bool geometricp = true);
  139. /**
  140. * A default constructor for the normal gravity. This sets up an
  141. * uninitialized object and is used by GravityModel which constructs this
  142. * object before it has read in the parameters for the reference ellipsoid.
  143. **********************************************************************/
  144. NormalGravity() : _a(-1) {}
  145. ///@}
  146. /** \name Compute the gravity
  147. **********************************************************************/
  148. ///@{
  149. /**
  150. * Evaluate the gravity on the surface of the ellipsoid.
  151. *
  152. * @param[in] lat the geographic latitude (degrees).
  153. * @return &gamma; the acceleration due to gravity, positive downwards
  154. * (m s<sup>&minus;2</sup>).
  155. *
  156. * Due to the axial symmetry of the ellipsoid, the result is independent of
  157. * the value of the longitude. This acceleration is perpendicular to the
  158. * surface of the ellipsoid. It includes the effects of the earth's
  159. * rotation.
  160. **********************************************************************/
  161. Math::real SurfaceGravity(real lat) const;
  162. /**
  163. * Evaluate the gravity at an arbitrary point above (or below) the
  164. * ellipsoid.
  165. *
  166. * @param[in] lat the geographic latitude (degrees).
  167. * @param[in] h the height above the ellipsoid (meters).
  168. * @param[out] gammay the northerly component of the acceleration
  169. * (m s<sup>&minus;2</sup>).
  170. * @param[out] gammaz the upward component of the acceleration
  171. * (m s<sup>&minus;2</sup>); this is usually negative.
  172. * @return \e U the corresponding normal potential
  173. * (m<sup>2</sup> s<sup>&minus;2</sup>).
  174. *
  175. * Due to the axial symmetry of the ellipsoid, the result is independent of
  176. * the value of the longitude and the easterly component of the
  177. * acceleration vanishes, \e gammax = 0. The function includes the effects
  178. * of the earth's rotation. When \e h = 0, this function gives \e gammay =
  179. * 0 and the returned value matches that of NormalGravity::SurfaceGravity.
  180. **********************************************************************/
  181. Math::real Gravity(real lat, real h, real& gammay, real& gammaz)
  182. const;
  183. /**
  184. * Evaluate the components of the acceleration due to gravity and the
  185. * centrifugal acceleration in geocentric coordinates.
  186. *
  187. * @param[in] X geocentric coordinate of point (meters).
  188. * @param[in] Y geocentric coordinate of point (meters).
  189. * @param[in] Z geocentric coordinate of point (meters).
  190. * @param[out] gammaX the \e X component of the acceleration
  191. * (m s<sup>&minus;2</sup>).
  192. * @param[out] gammaY the \e Y component of the acceleration
  193. * (m s<sup>&minus;2</sup>).
  194. * @param[out] gammaZ the \e Z component of the acceleration
  195. * (m s<sup>&minus;2</sup>).
  196. * @return \e U = <i>V</i><sub>0</sub> + &Phi; the sum of the
  197. * gravitational and centrifugal potentials
  198. * (m<sup>2</sup> s<sup>&minus;2</sup>).
  199. *
  200. * The acceleration given by <b>&gamma;</b> = &nabla;\e U =
  201. * &nabla;<i>V</i><sub>0</sub> + &nabla;&Phi; = <b>&Gamma;</b> + <b>f</b>.
  202. **********************************************************************/
  203. Math::real U(real X, real Y, real Z,
  204. real& gammaX, real& gammaY, real& gammaZ) const;
  205. /**
  206. * Evaluate the components of the acceleration due to the gravitational
  207. * force in geocentric coordinates.
  208. *
  209. * @param[in] X geocentric coordinate of point (meters).
  210. * @param[in] Y geocentric coordinate of point (meters).
  211. * @param[in] Z geocentric coordinate of point (meters).
  212. * @param[out] GammaX the \e X component of the acceleration due to the
  213. * gravitational force (m s<sup>&minus;2</sup>).
  214. * @param[out] GammaY the \e Y component of the acceleration due to the
  215. * @param[out] GammaZ the \e Z component of the acceleration due to the
  216. * gravitational force (m s<sup>&minus;2</sup>).
  217. * @return <i>V</i><sub>0</sub> the gravitational potential
  218. * (m<sup>2</sup> s<sup>&minus;2</sup>).
  219. *
  220. * This function excludes the centrifugal acceleration and is appropriate
  221. * to use for space applications. In terrestrial applications, the
  222. * function NormalGravity::U (which includes this effect) should usually be
  223. * used.
  224. **********************************************************************/
  225. Math::real V0(real X, real Y, real Z,
  226. real& GammaX, real& GammaY, real& GammaZ) const;
  227. /**
  228. * Evaluate the centrifugal acceleration in geocentric coordinates.
  229. *
  230. * @param[in] X geocentric coordinate of point (meters).
  231. * @param[in] Y geocentric coordinate of point (meters).
  232. * @param[out] fX the \e X component of the centrifugal acceleration
  233. * (m s<sup>&minus;2</sup>).
  234. * @param[out] fY the \e Y component of the centrifugal acceleration
  235. * (m s<sup>&minus;2</sup>).
  236. * @return &Phi; the centrifugal potential (m<sup>2</sup>
  237. * s<sup>&minus;2</sup>).
  238. *
  239. * &Phi; is independent of \e Z, thus \e fZ = 0. This function
  240. * NormalGravity::U sums the results of NormalGravity::V0 and
  241. * NormalGravity::Phi.
  242. **********************************************************************/
  243. Math::real Phi(real X, real Y, real& fX, real& fY) const;
  244. ///@}
  245. /** \name Inspector functions
  246. **********************************************************************/
  247. ///@{
  248. /**
  249. * @return true if the object has been initialized.
  250. **********************************************************************/
  251. bool Init() const { return _a > 0; }
  252. /**
  253. * @return \e a the equatorial radius of the ellipsoid (meters). This is
  254. * the value used in the constructor.
  255. **********************************************************************/
  256. Math::real EquatorialRadius() const
  257. { return Init() ? _a : Math::NaN(); }
  258. /**
  259. * @return \e GM the mass constant of the ellipsoid
  260. * (m<sup>3</sup> s<sup>&minus;2</sup>). This is the value used in the
  261. * constructor.
  262. **********************************************************************/
  263. Math::real MassConstant() const
  264. { return Init() ? _GM : Math::NaN(); }
  265. /**
  266. * @return <i>J</i><sub><i>n</i></sub> the dynamical form factors of the
  267. * ellipsoid.
  268. *
  269. * If \e n = 2 (the default), this is the value of <i>J</i><sub>2</sub>
  270. * used in the constructor. Otherwise it is the zonal coefficient of the
  271. * Legendre harmonic sum of the normal gravitational potential. Note that
  272. * <i>J</i><sub><i>n</i></sub> = 0 if \e n is odd. In most gravity
  273. * applications, fully normalized Legendre functions are used and the
  274. * corresponding coefficient is <i>C</i><sub><i>n</i>0</sub> =
  275. * &minus;<i>J</i><sub><i>n</i></sub> / sqrt(2 \e n + 1).
  276. **********************************************************************/
  277. Math::real DynamicalFormFactor(int n = 2) const
  278. { return Init() ? ( n == 2 ? _J2 : Jn(n)) : Math::NaN(); }
  279. /**
  280. * @return &omega; the angular velocity of the ellipsoid (rad
  281. * s<sup>&minus;1</sup>). This is the value used in the constructor.
  282. **********************************************************************/
  283. Math::real AngularVelocity() const
  284. { return Init() ? _omega : Math::NaN(); }
  285. /**
  286. * @return <i>f</i> the flattening of the ellipsoid (\e a &minus; \e b)/\e
  287. * a.
  288. **********************************************************************/
  289. Math::real Flattening() const
  290. { return Init() ? _f : Math::NaN(); }
  291. /**
  292. * @return &gamma;<sub>e</sub> the normal gravity at equator (m
  293. * s<sup>&minus;2</sup>).
  294. **********************************************************************/
  295. Math::real EquatorialGravity() const
  296. { return Init() ? _gammae : Math::NaN(); }
  297. /**
  298. * @return &gamma;<sub>p</sub> the normal gravity at poles (m
  299. * s<sup>&minus;2</sup>).
  300. **********************************************************************/
  301. Math::real PolarGravity() const
  302. { return Init() ? _gammap : Math::NaN(); }
  303. /**
  304. * @return <i>f*</i> the gravity flattening (&gamma;<sub>p</sub> &minus;
  305. * &gamma;<sub>e</sub>) / &gamma;<sub>e</sub>.
  306. **********************************************************************/
  307. Math::real GravityFlattening() const
  308. { return Init() ? _fstar : Math::NaN(); }
  309. /**
  310. * @return <i>U</i><sub>0</sub> the constant normal potential for the
  311. * surface of the ellipsoid (m<sup>2</sup> s<sup>&minus;2</sup>).
  312. **********************************************************************/
  313. Math::real SurfacePotential() const
  314. { return Init() ? _U0 : Math::NaN(); }
  315. /**
  316. * @return the Geocentric object used by this instance.
  317. **********************************************************************/
  318. const Geocentric& Earth() const { return _earth; }
  319. /**
  320. * \deprecated An old name for EquatorialRadius().
  321. **********************************************************************/
  322. GEOGRAPHICLIB_DEPRECATED("Use EquatorialRadius()")
  323. Math::real MajorRadius() const { return EquatorialRadius(); }
  324. ///@}
  325. /**
  326. * A global instantiation of NormalGravity for the WGS84 ellipsoid.
  327. **********************************************************************/
  328. static const NormalGravity& WGS84();
  329. /**
  330. * A global instantiation of NormalGravity for the GRS80 ellipsoid.
  331. **********************************************************************/
  332. static const NormalGravity& GRS80();
  333. /**
  334. * Compute the flattening from the dynamical form factor.
  335. *
  336. * @param[in] a equatorial radius (meters).
  337. * @param[in] GM mass constant of the ellipsoid
  338. * (meters<sup>3</sup>/seconds<sup>2</sup>); this is the product of \e G
  339. * the gravitational constant and \e M the mass of the earth (usually
  340. * including the mass of the earth's atmosphere).
  341. * @param[in] omega the angular velocity (rad s<sup>&minus;1</sup>).
  342. * @param[in] J2 the dynamical form factor.
  343. * @return \e f the flattening of the ellipsoid.
  344. *
  345. * This routine requires \e a &gt; 0, \e GM &gt; 0, \e J2 &lt; 1/3 &minus;
  346. * <i>omega</i><sup>2</sup><i>a</i><sup>3</sup>/<i>GM</i> 8/(45&pi;). A
  347. * NaN is returned if these conditions do not hold. The restriction to
  348. * positive \e GM is made because for negative \e GM two solutions are
  349. * possible.
  350. **********************************************************************/
  351. static Math::real J2ToFlattening(real a, real GM, real omega, real J2);
  352. /**
  353. * Compute the dynamical form factor from the flattening.
  354. *
  355. * @param[in] a equatorial radius (meters).
  356. * @param[in] GM mass constant of the ellipsoid
  357. * (meters<sup>3</sup>/seconds<sup>2</sup>); this is the product of \e G
  358. * the gravitational constant and \e M the mass of the earth (usually
  359. * including the mass of the earth's atmosphere).
  360. * @param[in] omega the angular velocity (rad s<sup>&minus;1</sup>).
  361. * @param[in] f the flattening of the ellipsoid.
  362. * @return \e J2 the dynamical form factor.
  363. *
  364. * This routine requires \e a &gt; 0, \e GM &ne; 0, \e f &lt; 1. The
  365. * values of these parameters are not checked.
  366. **********************************************************************/
  367. static Math::real FlatteningToJ2(real a, real GM, real omega, real f);
  368. };
  369. } // namespace GeographicLib
  370. #endif // GEOGRAPHICLIB_NORMALGRAVITY_HPP