mxlib
c++ tools for analyzing astronomical data and other tasks by Jared R. Males. [git repo]
Loading...
Searching...
No Matches
orbitUtils.hpp
Go to the documentation of this file.
1/** \file orbitUtils.hpp
2 * \author Jared R. Males
3 * \brief Utilities related to orbits.
4 * \ingroup astrofiles
5 *
6 */
7
8#ifndef mx_astro_orbitUtils_hpp
9#define mx_astro_orbitUtils_hpp
10
11#include "../math/constants.hpp"
12#include "kepler.hpp"
13
14namespace mx
15{
16namespace astro
17{
18
19/// Calculate the period of an orbit, given the masses and semi-major axis.
20/**
21 * \tparam units specifies the unit system and precision, see \ref astrounits.
22 *
23 * \returns the period in the specified units
24 *
25 * \ingroup orbits
26 */
27template <typename units>
28typename units::realT orbitPeriod( typename units::realT m1, ///< [in] mass of one object.
29 typename units::realT m2, ///< [in] mass of the other object (can be 0).
30 typename units::realT a ///< [in] the semi-major axis of the orbit.
31)
32{
33 return math::two_pi<typename units::realT>() * sqrt( ( a * a * a ) / ( constants::G<units>() * ( m1 + m2 ) ) );
34}
35
36/// Calculate the semi-major axis of an orbit, given the masses and Period, using SI units.
37/**
38 * \tparam units specifies the unit system and precision, see \ref astrounits.
39 *
40 * \returns the semi-major axis in the specified units
41 *
42 * \ingroup orbits
43 */
44template <typename units>
45typename units::realT orbitSemiMaj( typename units::realT m1, ///< [in] mass of one object.
46 typename units::realT m2, ///< [in] mass of the other object (can be 0).
47 typename units::realT P ///< [in] the period of the orbit.
48)
49{
50 typedef typename units::realT realT;
51
52 constexpr realT two = static_cast<realT>( 2 ); // Just for convenience.
53
54 return pow( pow( P / ( math::two_pi<realT>() ), two ) * constants::G<units>() * ( m1 + m2 ), math::third<realT>() );
55}
56
57/// Calculate the mean anomaly at time t, given the time of pericenter passage t0 and period P.
58/**
59 * \ingroup orbits
60 */
61template <typename realT>
62realT orbitMeanAnol( realT t, realT t0, realT P )
63{
64 return math::two_pi<realT>() * ( t - t0 ) / P;
65}
66
67/// Calculate the separation and true anomaly for an orbit at the specified mean anomaly.
68/**
69 * \tparam realT is the type used for arithmetic
70 *
71 * \returns the numer of iterations in solving Kepler's equation, if successful.
72 * \returns -1 if itmax is exceeded.
73 *
74 * \ingroup orbits
75 */
76template <typename realT>
77long orbitElements( realT &r, ///< [out] is the orbital separation
78 realT &f, ///< [out] is the true anomaly
79 realT E, ///< [out] is the eccentric anomaly, provided since it is free
80 realT D, ///< [out] is the error after the last iteration of the Kepler solution.
81 realT e, ///< [in] is the eccentricity
82 realT M, ///< [in] is the mean anomaly
83 realT a, ///< [in] is the semi-major axis
84 realT tol = KEPLER_TOL, ///< [in] is the desired tolerance
85 long itmax = KEPLER_ITMAX ///< [in] is the maximum number of iterations
86)
87{
88 long its;
89
90 if( e < 0.0 )
91 {
92 std::cerr << "e < 0 in rf_elements\n";
93 return -1;
94 }
95
96 if( e == 1.0 )
97 {
98 std::cerr << "e = 1 in rf_elements, which is not currently handled.\n";
99 return -1;
100 }
101
102 its = solve_kepler( E, D, e, M, tol, itmax );
103
104 if( e > 1.0 )
105 {
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 ) );
108 }
109 else if( e < 1.0 )
110 {
111 if( e == 0.0 )
112 {
113 f = M;
114 r = a;
115 }
116 else
117 {
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 ) );
120 }
121 }
122
123 return its;
124}
125
126/// Get the orbital phase at true anomaly f. Calculates the cos(alfa) where alfa is the orbital phase.
127/** Only calculates cos(alfa) to save an operation when possible.
128 *
129 * \tparam realT is the type used for arithmetic
130 *
131 * \returns the cosine of the phase angle
132 *
133 * \ingroup orbits
134 */
135template <typename realT>
136realT orbitPhaseCosine( realT f, ///< [in] is the true anomoaly at which to calculate alfa
137 realT w, ///< [in] the argument of pericenter of the orbit
138 realT inc ///< [in] is the inclinatin of the orbit
139)
140{
141 return sin( f + w ) * sin( inc );
142}
143
144/// Get the lambertian phase function at orbital phase specified as cos(alfa), where alfa is the phase angle.
145/** Uses cos(alfa) to save an operation after calculating orbit phase.
146 *
147 * \tparam realT is the type used for arithmetic
148 *
149 * \returns the lambert phase function at alpha
150 *
151 * \ingroup orbits
152 */
153template <typename realT>
154realT orbitLambertPhase( realT cos_alf /**< [in] the cosine of the phase angle*/ )
155{
156 realT alf = std::acos( cos_alf );
157 return ( sin( alf ) + ( math::pi<realT>() - alf ) * cos_alf ) / math::pi<realT>();
158}
159
160/// Calculate the radial velocity of an orbiting body inclined relative to an observer.
161/** \returns the radial velocity of the body with mass m1 as it orbits the barycenter of its orbit w.r.t. m2.
162 *
163 * \ingroup orbits
164 */
165template <typename units>
166typename units::realT orbitRV( typename units::realT m1, ///< [in] the mass of the body for which the RV is computed
167 typename units::realT m2, ///< [in] the mass of the other body in the orbiting pair
168 typename units::realT inc, ///< [in] the inclination of the orbit of m2.
169 typename units::realT a, ///< [in] the semi-major axis of the orbit
170 typename units::realT e, ///< [in] the eccentricity of the orbit
171 typename units::realT w, ///< [in] the argument of pericenter
172 typename units::realT f ///< [in] the true anomaly at which the RV is to be calculated.
173)
174{
175 typedef typename units::realT realT;
176
177 realT msini = m2 * sin( inc );
178
179 realT rv = (msini)*std::sqrt( ( constants::G<units>() / ( m1 + m2 ) ) / ( a * ( 1 - e * e ) ) );
180
181 rv *= ( std::cos( w + f ) + e * std::cos( w ) );
182
183 return rv;
184}
185
186/// Calculate various quantities of an orbit given the keplerian elements and a vector of times.
187/** Only those quantities with non-null pointers are actually calculated.
188 *
189 * The orbital parameters with dimensions, (a, P, and t0) must be in consistent units.
190 *
191 * \tparam vectorT a type whose elements are accessed with the [] operator. No other requirements are placed on this
192 * type, so it could be a native pointer (i.e. double *). \tparam realT the type in which to perform calculations. Does
193 * not need to match the storage type of vectorT.
194 *
195 * \retval -1 on error calculating the orbit (from rf_elements)
196 * \retval 0 on success.
197 *
198 * \ingroup orbits
199 */
200template <typename vectorT, typename realT>
202 vectorT *x, ///< [out] [optional] If not NULL, will be the projected x positions of the orbit. Must be at least as
203 ///< long as t.
204 vectorT *y, ///< [out] [optional] If not NULL, will be the projected y positions of the orbit. Must be at least as
205 ///< long as t.
206 vectorT *z, ///< [out] [optional] If not NULL, will be the projected z positions of the orbit. Must be at least as
207 ///< long as t.
208 vectorT *r, ///< [out] [optional] If not NULL, will be the separation of the orbit. Must be at least as long as t.
209 vectorT *rProj, ///< [out] [optional] If not NULL, will be the projected separation of the orbit. Must be at least
210 ///< as long as t.
211 vectorT
212 *f, ///< [out] [optional] If not NULL, will be the true anomaly of the orbit. Must be at least as long as t.
213 vectorT *cos_alf, ///< [out] [optional] If not NULL, will be the phase angle of the orbit. Must be at least as long
214 ///< as t.
215 vectorT
216 *phi, ///< [out] [optional] If not NULL, will be the Lambertian phase function. Must be at least as long as t.
217 const vectorT &t, ///< [in] the times at which to calculate the projected positions
218 const size_t N, ///< [in] the number of points contained in t, and allocated in nx, ny, nz. Avoids requiring a
219 ///< size() member (i.e. if vectorT is native pointers).
220 realT a, ///< [in] the semi-major axis of the orbit.
221 realT P, ///< [in] the orbital period.
222 realT e, ///< [in] the eccentricity of the orbit.
223 realT t0, ///< [in] the time of pericenter passage of the orbit.
224 realT i, ///< [in] the inclination of the orbit [radians].
225 realT w, ///< [in] the argument of pericenter of the orbit [radians].
226 realT W ///< [in] the longitude of the ascending node of the orbit [radians].
227)
228{
229 long rv;
230 realT _M, _r, _f, _E, _D;
231
232 // Just do these calcs once
233 realT cos_i = cos( i );
234 realT sin_i = sin( i );
235 realT cos_W = cos( W );
236 realT sin_W = sin( W );
237
238 realT cos_wf, sin_wf;
239
240 for( size_t j = 0; j < N; j++ )
241 {
242 _M = orbitMeanAnol( t[j], t0, P );
243
244 rv = orbitElements( _r, _f, _E, _D, e, _M, a );
245
246 if( rv < 0 )
247 return -1;
248
249 // one calc
250 cos_wf = cos( w + _f );
251 sin_wf = sin( w + _f );
252
253 // Assumption is that these will be optimized away if NULL
254 if( x )
255 ( *x )[j] = _r * ( cos_W * cos_wf - sin_W * sin_wf * cos_i );
256 if( y )
257 ( *y )[j] = _r * ( sin_W * cos_wf + cos_W * sin_wf * cos_i );
258 if( z )
259 ( *z )[j] = _r * sin_wf * sin_i;
260 if( r )
261 ( *r )[j] = _r;
262 if( rProj )
263 {
264 if( x && y )
265 {
266 ( *rProj )[j] = sqrt( ( *x )[j] * ( *x )[j] + ( *y )[j] * ( *y )[j] );
267 }
268 else if( x && !y )
269 {
270 realT _y = _r * ( sin_W * cos_wf + cos_W * sin_wf * cos_i );
271 ( *rProj )[j] = sqrt( ( *x )[j] * ( *x )[j] + _y * _y );
272 }
273 else if( !x && y )
274 {
275 realT _x = _r * ( cos_W * cos_wf - sin_W * sin_wf * cos_i );
276 ( *rProj )[j] = sqrt( _x * _x + ( *y )[j] * ( *y )[j] );
277 }
278 else
279 {
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 );
283 }
284 }
285
286 if( f )
287 ( *f )[j] = _f;
288
289 if( cos_alf )
290 ( *cos_alf )[j] = sin_wf * sin_i;
291
292 if( phi )
293 ( *phi )[j] = orbitLambertPhase<realT>( sin_wf * sin_i );
294 }
295 return 0;
296}
297
298/// Calculate the cartesian x-y position and the Lambert phase function of an orbit given Keplerian elements and a
299/// vector of times.
300/**
301 * \tparam vectorT a type whose elements are accessed with the [] operator.
302 * \tparam arithT the type in which to perform arithmetic
303 *
304 * \retval -1 on error calculating the orbit (from orbitCartesianWork)
305 * \retval 0 on success.
306 *
307 * \ingroup orbits
308 */
309template <typename vectorT, typename arithT>
311 vectorT &x, ///< [out] the projected x positions of the orbit. Must be at least as long as t.
312 vectorT &y, ///< [out] the projected y positions of the orbit. Must be at least as long as t.
313 vectorT &r, ///< [out] the separation of the orbit. Must be at least as long as t.
314 vectorT &rProj, ///< [out] the projected separation of the orbit. Must be at least as long as t.
315 vectorT &phi, ///< [out] the Lambert phase function. Must be at least as long as t.
316 const vectorT &t, ///< [in] the times at which to calculate the projected positions
317 const size_t N, ///< [in] the number of points contained in N, and allocated in nx, ny, nz
318 arithT a, ///< [in] the semi-major axis of the orbit
319 arithT P, ///< [in] the orbital period
320 arithT e, ///< [in] the eccentricity of the orbit
321 arithT t0, ///< [in] the time of pericenter passage of the orbit
322 arithT i, ///< [in] the inclination of the orbit
323 arithT w, ///< [in] the argument of pericenter of the orbit
324 arithT W ///< [in] the longitude of the ascending node of the orbit.
325)
326{
327 return orbitCartesianWork(
328 &x, &y, (vectorT *)0, &r, &rProj, (vectorT *)0, (vectorT *)0, &phi, t, N, a, P, e, t0, i, w, W );
329}
330
331} // namespace astro
332} // namespace mx
333
334#endif // mx_astro_orbitUtils_hpp
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.
Definition kepler.hpp:156
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.
The mxlib c++ namespace.
Definition mxError.hpp:106