8 #ifndef __mx_astro_constants_hpp__
9 #define __mx_astro_constants_hpp__
31 template<
typename realT>
34 return static_cast<realT
>(4.8481368111333441675396429478852851658848753880815e-06);
50 template<
typename units>
51 constexpr
typename units::realT
G()
53 return static_cast<typename units::realT
>(6.67408e-11) * (units::length*units::length*units::length)/(units::mass * units::time*units::time);
59 template<
typename units>
60 constexpr
typename units::realT
c()
62 return static_cast<typename units::realT
>(299792458)* (units::length/units::time);
70 template<
typename units>
71 constexpr
typename units::realT
k()
73 return static_cast<typename units::realT
>(1.38064852e-23) * units::energy / units::temperature;
81 template<
typename units>
82 constexpr
typename units::realT
sigma()
84 return static_cast<typename units::realT
>(5.670367e-8) * (units::energy/units::time) / (units::length*units::length) / (units::temperature*units::temperature*units::temperature*units::temperature);
90 template<
typename units>
91 constexpr
typename units::realT
h()
93 return static_cast<typename units::realT
>(6.626070040e-34) * units::energy * units::time;
99 template<
typename units>
100 constexpr
typename units::realT
day()
102 return static_cast<typename units::realT
>(86400.0)*units::time;
108 template<
typename units>
109 constexpr
typename units::realT
year()
111 return (
static_cast<typename units::realT
>(365.25)*
static_cast<typename units::realT
>(86400.0))*units::time;
117 template<
typename units>
118 constexpr
typename units::realT
au()
120 return static_cast<typename units::realT
>(149597870700.0)*units::length;
124 template<
typename units>
127 return au<units>()/tan_arcsec<typename units::realT>();
136 template<
typename units>
139 return static_cast<typename units::realT
>(6.957e8)*units::length;
147 template<
typename units>
150 return static_cast<typename units::realT
>(1361)*units::energy/units::time/(units::length*units::length);
158 template<
typename units>
161 return static_cast<typename units::realT
>(3.828e26)*units::energy/units::time;
169 template<
typename units>
172 return static_cast<typename units::realT
>(5772.0)*units::temperature;
180 template<
typename units>
181 constexpr
typename units::realT
GMSun()
183 return static_cast<typename units::realT
>(1.3271244e20)*(units::length*units::length*units::length)/(units::time*units::time);
192 template<
typename units>
195 return GMSun<units>()/G<units>();
203 template<
typename units>
206 return static_cast<typename units::realT
>(6.3781e6)*units::length;
214 template<
typename units>
217 return static_cast<typename units::realT
>(3.986004e14)*(units::length*units::length*units::length)/(units::time*units::time);
226 template<
typename units>
229 return GMEarth<units>()/G<units>();
237 template<
typename units>
240 return static_cast<typename units::realT
>(7.1492e7)*units::length;
248 template<
typename units>
251 return static_cast<typename units::realT
>(1.2668653e17)*(units::length*units::length*units::length)/(units::time*units::time);
260 template<
typename units>
263 return GMJupiter<units>()/G<units>();
276 template<
typename units>
279 return static_cast<typename units::realT
>(6.023657330e6);
287 template<
typename units>
290 return massSun<units>()/mrMercury<units>();
300 template<
typename units>
301 constexpr
typename units::realT radMercury()
303 return static_cast<typename units::realT
>(2.4397e6)*units::length;
314 template<
typename units>
315 constexpr
typename units::realT mrVenus()
317 return static_cast<typename units::realT
>(4.08523719e5);
325 template<
typename units>
326 constexpr
typename units::realT massVenus()
328 return massSun<units>()/mrVenus<units>();
338 template<
typename units>
341 return static_cast<typename units::realT
>(6.0518e6)*units::length;
352 template<
typename units>
355 return static_cast<typename units::realT
>( 3.09870359e6);
363 template<
typename units>
366 return massSun<units>()/mrMars<units>();
376 template<
typename units>
377 constexpr
typename units::realT radMars()
379 return static_cast<typename units::realT
>(3.39619e6)*units::length;
390 template<
typename units>
391 constexpr
typename units::realT mrSaturn()
393 return static_cast<typename units::realT
>( 3.4979018e3);
401 template<
typename units>
402 constexpr
typename units::realT massSaturn()
404 return massSun<units>()/mrSaturn<units>();
414 template<
typename units>
417 return static_cast<typename units::realT
>(6.0268e7)*units::length;
429 template<
typename units>
432 return static_cast<typename units::realT
>( 2.2902951e4);
440 template<
typename units>
443 return massSun<units>()/mrUranus<units>();
453 template<
typename units>
454 constexpr
typename units::realT radUranus()
456 return static_cast<typename units::realT
>(2.5559e7)*units::length;
467 template<
typename units>
468 constexpr
typename units::realT mrNeptune()
470 return static_cast<typename units::realT
>( 1.941226e4);
478 template<
typename units>
479 constexpr
typename units::realT massNeptune()
481 return massSun<units>()/mrNeptune<units>();
491 template<
typename units>
494 return static_cast<typename units::realT
>(2.4764e7)*units::length;
constexpr units::realT G()
Newton's Gravitational Constant.
constexpr units::realT radSaturn()
Radius of Mars.
constexpr units::realT radNeptune()
Radius of Uranus.
constexpr units::realT year()
Length of year.
constexpr units::realT GMJupiter()
Jupiter Mass Parameter.
constexpr units::realT radEarth()
Radius of Earth (nominal equatorial)
constexpr units::realT c()
The speed of light.
constexpr units::realT parsec()
The parsec.
constexpr units::realT massUranus()
Mass of Uranus.
constexpr units::realT radVenus()
Radius of Mercury.
constexpr units::realT mrMercury()
Mass ratio of Mercury to the Sun.
constexpr units::realT radSun()
Radius of the Sun.
constexpr units::realT massMars()
Mass of Mars.
constexpr units::realT GMEarth()
Earth Mass Parameter.
constexpr units::realT mrUranus()
Mass ratio of Uranus to the Sun.
constexpr units::realT massSun()
Solar Mass.
constexpr units::realT radJupiter()
Radius of Jupiter (nominal equatorial)
constexpr units::realT massMercury()
Mass of Mercury.
constexpr units::realT lumSun()
Luminosity of the Sun.
constexpr units::realT h()
Planck Constant.
constexpr units::realT GMSun()
Solar Mass Parameter.
constexpr units::realT sigma()
Stefan-Boltzmann Constant.
constexpr units::realT day()
Length of day.
constexpr units::realT mrMars()
Mass ratio of Mars to the Sun.
constexpr units::realT massEarth()
Earth Mass.
constexpr units::realT au()
Astronomical Unit.
constexpr units::realT TeffSun()
Effective Temperature of the Sun.
constexpr units::realT massJupiter()
Jupiter Mass.
constexpr units::realT solarIrrad()
Solar Irradiance at 1 au.
constexpr units::realT k()
Boltzmann Constant.
constexpr realT tan_arcsec()
Tangent of 1.0 arcsec.