8#ifndef mx_astro_orbitUtils_hpp
9#define mx_astro_orbitUtils_hpp
11#include "../math/constants.hpp"
27template <
typename units>
29 typename units::realT m2,
30 typename units::realT a
44template <
typename units>
46 typename units::realT m2,
47 typename units::realT P
50 typedef typename units::realT realT;
52 constexpr realT two =
static_cast<realT
>( 2 );
61template <
typename realT>
76template <
typename realT>
84 realT tol = KEPLER_TOL,
85 long itmax = KEPLER_ITMAX
92 std::cerr <<
"e < 0 in rf_elements\n";
98 std::cerr <<
"e = 1 in rf_elements, which is not currently handled.\n";
106 f = 2.0 * std::atan( std::sqrt( ( e + 1.0 ) / ( e - 1.0 ) ) * std::tanh( E / 2.0 ) );
107 r = a * ( 1.0 - e * std::cosh( E ) );
118 f = 2.0 * std::atan( std::sqrt( ( 1.0 + e ) / ( 1.0 - e ) ) * std::tan( E / 2.0 ) );
119 r = a * ( 1.0 - e * std::cos( E ) );
135template <
typename realT>
141 return sin( f + w ) * sin( inc );
153template <
typename realT>
156 realT alf = std::acos( cos_alf );
165template <
typename units>
166typename units::realT
orbitRV(
typename units::realT m1,
167 typename units::realT m2,
168 typename units::realT inc,
169 typename units::realT a,
170 typename units::realT e,
171 typename units::realT w,
172 typename units::realT f
175 typedef typename units::realT realT;
177 realT msini = m2 * sin( inc );
179 realT rv = (msini)*std::sqrt( ( constants::G<units>() / ( m1 + m2 ) ) / ( a * ( 1 - e * e ) ) );
181 rv *= ( std::cos( w + f ) + e * std::cos( w ) );
200template <
typename vectorT,
typename realT>
230 realT _M, _r, _f, _E, _D;
233 realT cos_i = cos( i );
234 realT sin_i = sin( i );
235 realT cos_W = cos( W );
236 realT sin_W = sin( W );
238 realT cos_wf, sin_wf;
240 for(
size_t j = 0; j < N; j++ )
250 cos_wf = cos( w + _f );
251 sin_wf = sin( w + _f );
255 ( *x )[j] = _r * ( cos_W * cos_wf - sin_W * sin_wf * cos_i );
257 ( *y )[j] = _r * ( sin_W * cos_wf + cos_W * sin_wf * cos_i );
259 ( *z )[j] = _r * sin_wf * sin_i;
266 ( *rProj )[j] = sqrt( ( *x )[j] * ( *x )[j] + ( *y )[j] * ( *y )[j] );
270 realT _y = _r * ( sin_W * cos_wf + cos_W * sin_wf * cos_i );
271 ( *rProj )[j] = sqrt( ( *x )[j] * ( *x )[j] + _y * _y );
275 realT _x = _r * ( cos_W * cos_wf - sin_W * sin_wf * cos_i );
276 ( *rProj )[j] = sqrt( _x * _x + ( *y )[j] * ( *y )[j] );
280 realT _x = _r * ( cos_W * cos_wf - sin_W * sin_wf * cos_i );
281 realT _y = _r * ( sin_W * cos_wf + cos_W * sin_wf * cos_i );
282 ( *rProj )[j] = sqrt( _x * _x + _y * _y );
290 ( *cos_alf )[j] = sin_wf * sin_i;
293 ( *phi )[j] = orbitLambertPhase<realT>( sin_wf * sin_i );
309template <
typename vectorT,
typename arithT>
328 &x, &y, (vectorT *)0, &r, &rProj, (vectorT *)0, (vectorT *)0, &phi, t, N, a, P, e, t0, i, w, W );
constexpr floatT six_fifths()
Return 6/5 in the specified precision.
long solve_kepler(realT &E, realT &D, realT e, realT M, realT tol=KEPLER_TOL, long itmax=KEPLER_ITMAX)
Solve Kepler's equation for any e. Uses solve_kepler_danby if e < 1.0, hyperbolic_kepler otherwise.
realT orbitPhaseCosine(realT f, realT w, realT inc)
Get the orbital phase at true anomaly f. Calculates the cos(alfa) where alfa is the orbital phase.
realT orbitLambertPhase(realT cos_alf)
Get the lambertian phase function at orbital phase specified as cos(alfa), where alfa is the phase an...
units::realT orbitRV(typename units::realT m1, typename units::realT m2, typename units::realT inc, typename units::realT a, typename units::realT e, typename units::realT w, typename units::realT f)
Calculate the radial velocity of an orbiting body inclined relative to an observer.
units::realT orbitSemiMaj(typename units::realT m1, typename units::realT m2, typename units::realT P)
Calculate the semi-major axis of an orbit, given the masses and Period, using SI units.
int orbitCartesianWork(vectorT *x, vectorT *y, vectorT *z, vectorT *r, vectorT *rProj, vectorT *f, vectorT *cos_alf, vectorT *phi, const vectorT &t, const size_t N, realT a, realT P, realT e, realT t0, realT i, realT w, realT W)
Calculate various quantities of an orbit given the keplerian elements and a vector of times.
long orbitElements(realT &r, realT &f, realT E, realT D, realT e, realT M, realT a, realT tol=KEPLER_TOL, long itmax=KEPLER_ITMAX)
Calculate the separation and true anomaly for an orbit at the specified mean anomaly.
int orbitCartesian2DPhi(vectorT &x, vectorT &y, vectorT &r, vectorT &rProj, vectorT &phi, const vectorT &t, const size_t N, arithT a, arithT P, arithT e, arithT t0, arithT i, arithT w, arithT W)
units::realT orbitPeriod(typename units::realT m1, typename units::realT m2, typename units::realT a)
Calculate the period of an orbit, given the masses and semi-major axis.
realT orbitMeanAnol(realT t, realT t0, realT P)
Calculate the mean anomaly at time t, given the time of pericenter passage t0 and period P.
Declarations for the utilities related to the Kepler problem.