Vulcan
Aerospace Engineering Utilities Built on Janus
Loading...
Searching...
No Matches
AnomalyConversions.hpp
Go to the documentation of this file.
1// Vulcan Anomaly Conversions
2// Mean, Eccentric, and True anomaly solvers
3#pragma once
4
5#include <cmath>
6#include <janus/janus.hpp>
8
10
24template <typename Scalar>
25Scalar mean_to_eccentric(const Scalar &M, const Scalar &e, double tol = 1e-12,
26 int max_iter = 50) {
27 // Initial guess
28 Scalar E = M;
29
30 if constexpr (std::is_same_v<Scalar, double>) {
31 // Newton-Raphson for numeric
32 for (int iter = 0; iter < max_iter; ++iter) {
33 double f = E - e * std::sin(E) - M;
34 double f_prime = 1.0 - e * std::cos(E);
35 double dE = f / f_prime;
36 E -= dE;
37 if (std::abs(dE) < tol)
38 break;
39 }
40 } else {
41 // Fixed Newton steps for symbolic (enables autodiff)
42 for (int iter = 0; iter < 10; ++iter) {
43 Scalar f = E - e * janus::sin(E) - M;
44 Scalar f_prime = 1.0 - e * janus::cos(E);
45 E = E - f / f_prime;
46 }
47 }
48
49 return E;
50}
51
60template <typename Scalar>
61Scalar eccentric_to_true(const Scalar &E, const Scalar &e) {
62 const Scalar cos_E = janus::cos(E);
63 const Scalar sin_E = janus::sin(E);
64 const Scalar sqrt_term = janus::sqrt(1.0 - e * e);
65
66 return janus::atan2(sqrt_term * sin_E, cos_E - e);
67}
68
77template <typename Scalar>
78Scalar true_to_eccentric(const Scalar &nu, const Scalar &e) {
79 const Scalar cos_nu = janus::cos(nu);
80 const Scalar sin_nu = janus::sin(nu);
81 const Scalar sqrt_term = janus::sqrt(1.0 - e * e);
82
83 return janus::atan2(sqrt_term * sin_nu, e + cos_nu);
84}
85
96template <typename Scalar>
97Scalar eccentric_to_mean(const Scalar &E, const Scalar &e) {
98 return E - e * janus::sin(E);
99}
100
109template <typename Scalar>
110Scalar true_to_mean(const Scalar &nu, const Scalar &e) {
111 Scalar E = true_to_eccentric(nu, e);
112 return eccentric_to_mean(E, e);
113}
114
123template <typename Scalar>
124Scalar mean_to_true(const Scalar &M, const Scalar &e) {
125 Scalar E = mean_to_eccentric(M, e);
126 return eccentric_to_true(E, e);
127}
128
129} // namespace vulcan::orbital::anomaly
Definition AnomalyConversions.hpp:9
Scalar mean_to_true(const Scalar &M, const Scalar &e)
Convert mean anomaly to true anomaly.
Definition AnomalyConversions.hpp:124
Scalar mean_to_eccentric(const Scalar &M, const Scalar &e, double tol=1e-12, int max_iter=50)
Solve Kepler's equation: M = E - e*sin(E).
Definition AnomalyConversions.hpp:25
Scalar true_to_mean(const Scalar &nu, const Scalar &e)
Convert true anomaly to mean anomaly.
Definition AnomalyConversions.hpp:110
Scalar eccentric_to_mean(const Scalar &E, const Scalar &e)
Convert eccentric anomaly to mean anomaly.
Definition AnomalyConversions.hpp:97
Scalar true_to_eccentric(const Scalar &nu, const Scalar &e)
Convert true anomaly to eccentric anomaly.
Definition AnomalyConversions.hpp:78
Scalar eccentric_to_true(const Scalar &E, const Scalar &e)
Convert eccentric anomaly to true anomaly.
Definition AnomalyConversions.hpp:61