True, Eccentric, and Mean Anomaly¶
This section deals with the conversion between true, eccentric, and mean anomaly.
True anomaly, frequently denoted \(\nu\), is the angular parameter that defines the position of an object moving along a Keplerian orbit. It is the angle between the eccentricity vector (vector pointing from the main pericenter to the periapsis) and the current position of the body in the orbital plane itself.
The eccentric anomaly, \(E\), is another angular parameter that defines the position of an object moving along a Keplerian orbit if viewed from the center of the ellipse.
Finally, the mean anomaly, \(M\), defines the fraction of an orbital period that has elapsed since the orbiting object has passed its periapsis. It is the angle from the pericenter an object moving on a fictitious circular orbit with the same semi-major axis would have progressed through in the same time as the body on the true elliptical orbit.
Conversion between all types of angular anomaly is possible. However, there is no known direct conversion between true and mean anomaly. Conversion between the two is accomplished by transformation through eccentric anomaly.
True and Eccentric Anomaly Conversions¶
To convert from true anomaly to eccentric anomaly, you can use the function
anomaly_eccentric_to_true
. To perform the reverse conversion use
anomaly_true_to_eccentric
.
Eccentric anomaly can be converted to true anomaly by using equations derived using equations
from Vallado1. Equation (2-12)
$$
\sin{\nu} = \frac{\sin{E}\sqrt{1-e^2}}{1 - e\cos{E}}
$$
can be divided by Equation (2-10)
$$
\cos{\nu} = \frac{\cos{E}-e}{1 - e\cos{E}}
$$
And rearranged to get
$$
\nu = \arctan{\frac{\sin{E}\sqrt{1-e^2}}{\cos{E}-e}}
$$
Which is what is implemented by anomaly_eccentric_to_true
. Similarly, Equations (2-9) from
Vallado can be rearranged to get
$$
E = \arctan{\frac{\sin{\nu}\sqrt{1-e^2}}{\cos{\nu}+e}}
$$
which allows for conversion from true anomaly to eccentric anomaly and is implemented in
anomaly_true_to_eccentric
.
use approx::assert_abs_diff_eq;
use rastro::orbits::{anomaly_true_to_eccentric, anomaly_eccentric_to_true};
fn main() {
let nu = 45.0; // Starting true anomaly
let e = 0.01; // Eccentricity
// Convert to eccentric anomaly
let ecc_anomaly = anomaly_true_to_eccentric(nu, e, true);
// Convert back from eccentric to true anomaly
let nu_2 = anomaly_eccentric_to_true(ecc_anomaly, e, true);
// Confirm equality to within tolerance
assert_abs_diff_eq!(nu, nu_2, epsilon=1e-14);
}
import pytest
import rastro
if __name__ == '__main__':
nu = 45.0 # Starting true anomaly
e = 0.01 # Eccentricity
# Convert to eccentric anomaly
ecc_anomaly = rastro.anomaly_true_to_eccentric(nu, e, True)
# Convert back from eccentric to true anomaly
nu_2 = rastro.anomaly_eccentric_to_true(ecc_anomaly, e, True)
# Confirm equality to within tolerance
assert nu == pytest.approx(nu_2, abs=1e-14)
Plot Source
# Generate plot of eccentric anomaly versus true anomaly for a range of eccentricies.
# Highlights the effect of eccentricity on the difference of the two.
import os
import pathlib
import plotly.graph_objects as go
import plotly.io as pio
import rastro
## Define Constants
SCRIPT_NAME = pathlib.Path(__file__).stem
OUTDIR = os.getenv("RASTRO_FIGURE_OUTPUT_DIR") # Build Environment Variable
OUTFILE = f"{OUTDIR}/{SCRIPT_NAME}.html"
## Create figure
fig = go.Figure()
fig.update_layout(dict(paper_bgcolor='rgba(0,0,0,0)', plot_bgcolor='rgba(0,0,0,0)'))
fig.update_xaxes(
showgrid=True, gridwidth=1, gridcolor='LightGrey', range=[0, 360],
showline=True, linewidth=2, linecolor='Grey'
)
fig.update_yaxes(
showgrid=True, gridwidth=1, gridcolor='LightGrey', range=[0, 360],
showline=True, linewidth=2, linecolor='Grey'
)
fig.update_layout(
xaxis=dict(tickmode='linear', tick0=0, dtick=30,
title_text=r"True Anomaly (deg)"),
yaxis=dict(tickmode='linear', tick0=0, dtick=30,
title_text=r"Eccentric Anomaly (deg)")
)
## Generate and plot data
# Generate range of true anomalies
nu = [x for x in range(0, 360)]
# Compute and plot eccentric anomaly for range of true anomalies
for e in [0.0, 0.1, 0.3, 0.5, 0.7, 0.9]:
# Take output mod 360 to wrap from 0 to 2pi
ecc = [rastro.anomaly_true_to_eccentric(x, e, True) % 360 for x in nu]
fig.add_trace(go.Scatter(x=nu, y=ecc, name=f"e = {e:.1f}"))
pio.write_html(fig, file=OUTFILE, include_plotlyjs='cdn', full_html=False, auto_play=False)
Eccentric and Mean Anomaly Conversions¶
To convert from true anomaly to eccentric anomaly, you can use the function
anomaly_eccentric_to_mean
. To perform the reverse conversion use
anomaly_mean_to_eccentric
.
Conversion from eccentric anomaly to mean anomaly is accomplished by application of Kepler's
equation
$$
M = E - e\sin{E}
$$
which is implemented in anomaly_eccentric_to_mean
.
Converting back from mean anomaly to eccentric anomaly is more challenging.
There is no known closed-form solution to convert from mean anomaly to eccentric anomaly.
Instead, we introduce the auxiliary equation
$$
f(E) = E - e\sin(E) - M
$$
And treat the problem as numerically solving for the root of \(f\) for a given \(M\). This iteration
can be accomplished using Newton's method. Starting from an initial guess \(E_0\) the value of
\(E_*\) can be iteratively updated using
$$
E_{i+1} = \frac{f(E_i)}{f^\prime(E_i)}= E_i - \frac{E_i - e\sin{E_i} - M}{1 - e\cos{E_i}}
$$
The algorithm is run until a coverage threshold of
$$
|E_{i+1} - E_i| \leq \Delta_{\text{tol}}
$$
is reached. The threshold set as 100 times floating-point machine precision 100 * f64::epsilon
.
This conversion is provided by anomaly_mean_to_eccentric
.
Warning
Because this is a numerical method, convergence is not guaranteed. There is an upper
limit of 10 iterations to reach convergence. Since convergence may not occur the output of
the function is a Result
, forcing the user to explicitly handle the case where the algorithm
does not converage.
Since Python lacks Rust's same error handling mechanisms, non-convergence will result in a runtime error.
use approx::assert_abs_diff_eq;
use rastro::orbits::{anomaly_eccentric_to_mean, anomaly_mean_to_eccentric};
fn main() {
let ecc = 45.0; // Starting true anomaly
let e = 0.01; // Eccentricity
// Convert to eccentric anomaly
let mean_anomaly = anomaly_eccentric_to_mean(ecc, e, true);
// Convert back from eccentric to true anomaly
let ecc_2 = anomaly_mean_to_eccentric(mean_anomaly, e, true).unwrap();
// Confirm equality to within tolerance
assert_abs_diff_eq!(ecc, ecc_2, epsilon=1e-14);
}
import pytest
import rastro
if __name__ == '__main__':
ecc = 45.0 # Starting eccentric anomaly
e = 0.01 # Eccentricity
# Convert to mean anomaly
mean_anomaly = rastro.anomaly_eccentric_to_mean(ecc, e, True)
# Convert back from mean to eccentric anomaly
ecc_2 = rastro.anomaly_mean_to_eccentric(mean_anomaly, e, True)
# Confirm equality to within tolerance
assert ecc == pytest.approx(ecc_2, abs=1e-14)
Plot Source
# Generate plot of mean anomaly versus eccentric anomaly for a range of eccentricies.
# Highlights the effect of eccentricity on the difference of the two.
import os
import pathlib
import plotly.graph_objects as go
import plotly.io as pio
import rastro
## Define Constants
SCRIPT_NAME = pathlib.Path(__file__).stem
OUTDIR = os.getenv("RASTRO_FIGURE_OUTPUT_DIR") # Build Environment Variable
OUTFILE = f"{OUTDIR}/{SCRIPT_NAME}.html"
## Create figure
fig = go.Figure()
fig.update_layout(dict(paper_bgcolor='rgba(0,0,0,0)', plot_bgcolor='rgba(0,0,0,0)'))
fig.update_xaxes(
showgrid=True, gridwidth=1, gridcolor='LightGrey', range=[0, 360],
showline=True, linewidth=2, linecolor='Grey'
)
fig.update_yaxes(
showgrid=True, gridwidth=1, gridcolor='LightGrey', range=[0, 360],
showline=True, linewidth=2, linecolor='Grey'
)
fig.update_layout(
xaxis=dict(tickmode='linear', tick0=0, dtick=30,
title_text=r"True Anomaly (deg)"),
yaxis=dict(tickmode='linear', tick0=0, dtick=30,
title_text=r"Eccentric Anomaly (deg)")
)
## Generate and plot data
# Generate range of true anomalies
ecc = [x for x in range(0, 360)]
# Compute and plot eccentric anomaly for range of true anomalies
for e in [0.0, 0.1, 0.3, 0.5, 0.7, 0.9]:
# Take output mod 360 to wrap from 0 to 2pi
mean = [rastro.anomaly_eccentric_to_mean(x, e, True) % 360 for x in ecc]
fig.add_trace(go.Scatter(x=ecc, y=mean, name=f"e = {e:.1f}"))
pio.write_html(fig, file=OUTFILE, include_plotlyjs='cdn', full_html=False, auto_play=False)
True and Mean Anomaly Conversions¶
Methods to convert from true anomaly to mean anomaly are
provided for convenience. These methods simply wrap successive calls to two
anomaly_true_to_mean
. To perform the reverse conversion use
anomaly_mean_to_true
.
use approx::assert_abs_diff_eq;
use rastro::orbits::{anomaly_eccentric_to_mean, anomaly_mean_to_eccentric};
fn main() {
let nu = 45.0; // Starting true anomaly
let e = 0.01; // Eccentricity
// Convert to eccentric anomaly
let mean_anomaly = anomaly_eccentric_to_mean(nu, e, true);
// Convert back from eccentric to true anomaly
let nu_2 = anomaly_mean_to_eccentric(mean_anomaly, e, true).unwrap();
// Confirm equality to within tolerance
assert_abs_diff_eq!(nu, nu_2, epsilon=1e-14);
}
import pytest
import rastro
if __name__ == '__main__':
nu = 45.0 # Starting true anomaly
e = 0.01 # Eccentricity
# Convert to mean anomaly
mean_anomaly = rastro.anomaly_true_to_mean(nu, e, True)
# Convert back from eccentric to true anomaly
nu_2 = rastro.anomaly_mean_to_true(mean_anomaly, e, True)
# Confirm equality to within tolerance
assert nu == pytest.approx(nu_2, abs=1e-14)
Plot Source
# Generate plot of mean anomaly versus true anomaly for a range of eccentricies.
# Highlights the effect of eccentricity on the difference of the two.
import os
import pathlib
import plotly.graph_objects as go
import plotly.io as pio
import rastro
## Define Constants
SCRIPT_NAME = pathlib.Path(__file__).stem
OUTDIR = os.getenv("RASTRO_FIGURE_OUTPUT_DIR") # Build Environment Variable
OUTFILE = f"{OUTDIR}/{SCRIPT_NAME}.html"
## Create figure
fig = go.Figure()
fig.update_layout(dict(paper_bgcolor='rgba(0,0,0,0)', plot_bgcolor='rgba(0,0,0,0)'))
fig.update_xaxes(
showgrid=True, gridwidth=1, gridcolor='LightGrey', range=[0, 360],
showline=True, linewidth=2, linecolor='Grey'
)
fig.update_yaxes(
showgrid=True, gridwidth=1, gridcolor='LightGrey', range=[0, 360],
showline=True, linewidth=2, linecolor='Grey'
)
fig.update_layout(
xaxis=dict(tickmode='linear', tick0=0, dtick=30,
title_text=r"True Anomaly (deg)"),
yaxis=dict(tickmode='linear', tick0=0, dtick=30,
title_text=r"Eccentric Anomaly (deg)")
)
## Generate and plot data
# Generate range of true anomalies
nu = [x for x in range(0, 360)]
# Compute and plot eccentric anomaly for range of true anomalies
for e in [0.0, 0.1, 0.3, 0.5, 0.7, 0.9]:
# Take output mod 360 to wrap from 0 to 2pi
mean = [rastro.anomaly_true_to_mean(x, e, True) % 360 for x in nu]
fig.add_trace(go.Scatter(x=nu, y=mean, name=f"e = {e:.1f}"))
pio.write_html(fig, file=OUTFILE, include_plotlyjs='cdn', full_html=False, auto_play=False)
-
D. Vallado, Fundamentals of Astrodynamics and Applications (4th Ed.), 2010
https://celestrak.com/software/vallado-sw.php ↩