Skip to content

Orbital Properties

The orbits module also provides functions to help derive useful properties of an orbit.

Periapsis Properties

Two common properties of interest are the distance and velocity of an object at the periapsis of its orbit. The periapsis will be the point of the closest approach as well as the point of greatest speed throughout the orbit.

Info

The word periapsis is formed by combination of the Greek words "peri-", meaning around, about and "apsis" meaning "arch or vault". An apsis is the farthest or nearest point in the orbit of a planetary body about its primary body.

Thereforce periapsis of an orbit is the point of closest approach of the orbiting body with respect to its central body. The suffix can be modified to indicate the point of closest approach to a specific celestical body. The perigee is the point of cloest approach to an object orbiting Earth. The perihelion is the point of closest approach to the Sun.

The periapsis velocity is given by periapsis_velocity or, for an Earth-orbiting object, the function perigee_velocity can be used. perigee_velocity simplified the call by supplying the gravitational parameter of Earth to the function call. Periapsis velocity is calculated by $$ v_{p} = \sqrt{\frac{\mu}{a}}\sqrt{\frac{1+e}{1-e}} $$ where \(\mu\) is the gravitational constant of the central body, \(a\) is the semi-major axis of the orbit, and \(e\) is the eccentricity.

Another useful parameter of the periapsis is the distance of the object to the central body. Equation (2-75) from Vallado1 $$ r_p = \frac{a(1-e^2)}{1+e} = a(1-e) $$

use rastro::constants::{R_EARTH, GM_EARTH};
use rastro::orbits::{periapsis_velocity, perigee_velocity, periapsis_distance};

fn main() {
    let a = R_EARTH + 500.0e3;
    let e = 0.01;

    // Compute periapsis velocity
    let periapsis_velocity = periapsis_velocity(a, e, GM_EARTH);
    println!("Periapsis velocity: {:.3}", periapsis_velocity);
    // Periapsis velocity: 7689.119

    // Compute as a perigee velocity
    let perigee_velocity = perigee_velocity(a, e);
    println!("Perigee velocity:   {:.3}", perigee_velocity);
    assert_eq!(periapsis_velocity, perigee_velocity);
    // Perigee velocity:   7689.119


    // Compute periapsis distance
    let periapsis_distance = periapsis_distance(a, e);
    println!("Periapsis distance: {:.3}", periapsis_distance);
    // Periapsis distance: 6809354.937
}
import rastro

if __name__ == '__main__':
    a = rastro.R_EARTH + 500.0e3
    e = 0.01

    # Compute periapsis velocity
    periapsis_velocity = rastro.periapsis_velocity(a, e, rastro.GM_EARTH)
    print(f"Periapsis velocity: {periapsis_velocity:.3f}")
    # Periapsis velocity: 7689.119

    # Compute as a perigee velocity
    perigee_velocity = rastro.perigee_velocity(a, e)
    print(f"Perigee velocity:   {perigee_velocity:.3f}")
    assert periapsis_velocity == perigee_velocity
    # Perigee velocity:   7689.119

    # Compute periapsis distance
    periapsis_distance = rastro.periapsis_distance(a, e)
    print(f"Periapsis distance: {periapsis_distance:.3f}")
    # Periapsis distance: 6809354.937

Apoapsis Properties

The distance and velocity of an object at the apoapsis of an orbit. The apoapsis is the furthest point from the central body. It is also the point of lowest speed throughout the orbit.

Info

The word apoapsis is formed by combination of the Greek words "apo-", meaning away from, separate, or apart from and the word "apsis".

Thereforce apoapsis of an orbit is the farthest point of an orbiting body with respect to its central body. The suffix can be modified to indicate the farthest point with respect to a specific primary celestical body. The apogee is furthest point away for an object orbiting Earth. The aphelion is the furthest away from an object orbiting the Sun.

The apoapsis velocity is given by apoapsis_velocity or, for an Earth-orbiting object, the function apoapsis_velocity can be used. Apoapsis velocity is given by $$ v_{a} = \sqrt{\frac{\mu}{a}}\sqrt{\frac{1-e}{1+e}} $$ The distance of the object to the central body is given by $$ r_a = \frac{a(1-e^2)}{1+e} = a(1+e) $$

Warning

The apoapsis position and velocity are valid for elliptic and circular orbits. For parabolic and hyperbolic orbits these two quantities are undefined.

use rastro::constants::{R_EARTH, GM_EARTH};
use rastro::orbits::{apoapsis_velocity, apogee_velocity, apoapsis_distance};

fn main() {
    let a = R_EARTH + 500.0e3;
    let e = 0.01;

    // Compute periapsis velocity
    let apoapsis_velocity = apoapsis_velocity(a, e, GM_EARTH);
    println!("Apoapsis velocity: {:.3}", apoapsis_velocity);
    // Apoapsis velocity: 7536.859

    // Compute as a perigee velocity
    let apogee_velocity = apogee_velocity(a, e);
    println!("Apogee velocity:   {:.3}", apogee_velocity);
    assert_eq!(apoapsis_velocity, apogee_velocity);
    // Apogee velocity:   7536.859


    // Compute periapsis distance
    let apoapsis_distance = apoapsis_distance(a, e);
    println!("Apoapsis distance: {:.3}", apoapsis_distance);
    // Apoapsis distance: 6946917.663
}
import rastro

if __name__ == '__main__':
    a = rastro.R_EARTH + 500.0e3
    e = 0.01

    # Compute periapsis velocity
    apoapsis_velocity = rastro.apoapsis_velocity(a, e, rastro.GM_EARTH)
    print(f"Apoapsis velocity: {apoapsis_velocity:.3f}")
    # Apoapsis velocity: 7536.859

    # Compute as a apogee velocity
    apogee_velocity = rastro.apogee_velocity(a, e)
    print(f"Apogee velocity:   {apogee_velocity:.3f}")
    assert apoapsis_velocity == apogee_velocity
    # Perigee velocity:   7536.859

    # Compute apoapsis distance
    apoapsis_distance = rastro.apoapsis_distance(a, e)
    print(f"Apoapsis distance: {apoapsis_distance:.3f}")
    # Periapsis distance: 6946917.663

Mean Motion and Semi-major Axis

A satellite's average angular rate over one orbit is it's mean motion \(n\), which can be calculated from the semi-major axis \(a\) and central body's gravitational parameter \(\mu\) $$ n = \sqrt{\frac{\mu}{a^3}} $$ mean_motion will calculate the mean motion of an object for an Earth-orbiting object while, mean_motion_general can calculate it for any arbitrary body when provided the graviational parameter of that body.

Since for some orbital data exchange formats, an object's orbit is characterized in terms of its mean motion instead of semi-major axis, RAstro provides semimajor_axis and semimajor_axis_general to invert above equation and recover the semi-major axis $$ a = \sqrt[3]{\frac{\mu}{n^2}} $$

Orbital Period

The orbital period \(T\) of a satellite is the time it takes for an orbit to progress through a full revolution. It is given by $$ T = 2\pi\sqrt{\frac{a^3}{\mu}} $$ orbital_period will return the orbital period for an Earth-orbiting object and orbital_period_general will do the same for any arbitrary body for which the gravitational constant is known.

The plot below shows how, for a circular orbit, as the semi-major axis increases, the orbital period increases while the speed decreases.

Plot Source
fig_orbital_period.py
# 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 numpy as np
import plotly.graph_objects as go
import plotly.io as pio
from   plotly.subplots import make_subplots
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 = make_subplots(specs=[[{"secondary_y": True}]])
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, 40],
    showline=True, linewidth=2, linecolor='Grey'
)
fig.update_yaxes(
    showgrid=True, gridwidth=1, gridcolor='LightGrey', range=[0, 10],
    showline=True, linewidth=2, linecolor='Grey'
)
fig.update_yaxes(
    showgrid=False, range=[0, 30], secondary_y=True
)
fig.update_xaxes(
    tickmode='linear', tick0=0, dtick=1, title_text=r"Satellite Altitude [1000 km]"
)
fig.update_yaxes(
    tickmode='linear', tick0=0, dtick=1, title_text=r"Velocity [km/s]"
)
fig.update_yaxes(
    tickmode='linear', tick0=0, dtick=5, title_text=r"Orbital Period [hours]", secondary_y=True
)

## Generate and plot data

# Generate range of true anomalies from 0 to 41,000 km altitude in 1,000 km increments
alt = np.arange(0, 41000*1e3, 500*1e3)

# Compute velocity over altitude
vp = [rastro.perigee_velocity(rastro.R_EARTH + a, 0.0)/1e3 for a in alt]
fig.add_trace(go.Scatter(x=alt/1e6, y=vp, name="Velocity"), secondary_y=False)


# Compute orbital period over altitude
period = [rastro.orbital_period(rastro.R_EARTH + a)/3600 for a in alt]
fig.add_trace(go.Scatter(x=alt/1e6, y=period, name="Orbital Period"), secondary_y=True)

pio.write_html(fig, file=OUTFILE, include_plotlyjs='cdn', full_html=False, auto_play=False)

Sun-Synchronous Inclination

A Sun-syncrhonous orbit is one whose nodal precession rate matches the average rate of the Earth's rotation about the Sun. That is, it is an orbit where the right ascension of the ascending node changes at the same rate as the Sun moves relative to the Earth. Sun-synchronous orbits are often highly relevant an useful for optical Earth observation satellites that desire to have consistent illumination conditions. A Sun-synchronous orbit is guaranteed to cross the equator at the same local time (e.g. 2pm) at each crossing event.

Due to Earth's oblateness known as the \(J_2\) zonal harmonic, frequently referred as simply \(J_2\), there is a constant, secular drift of the right ascension of all Earth orbiting objects. Since the \(J_2\) harmonic is the second largest after that of point-mass gravity, it is a dominant effect on the orbit

\[ \dot{\Omega} = -\frac{3nR^2_EJ_2}{2a^2(1-e^2)^2} \]

when combined with the required nodal precession rate for sun synchronicity

\[ \dot{\Omega} = \frac{360 \; \text{deg}}{1 \; \text{year}}\frac{1 \; \text{year}}{365.2421897 \; \text{day}} = 0.98564736 \frac{\text{deg}}{\text{day}} \]

can be rearranged to solve for the inclination as a function of semi-major axis and eccentricity

\[ i = \arccos{\Big(-\frac{2a^{\frac{7}{2}}\dot{\Omega}_{ss}(1-e^2)^2}{3R^2_EJ_2\sqrt{\mu}}\Big)} \]

The function sun_syncrhonous_inclination calculates this inclination.

The figure below shows how the inclination required to maintain a sun-synchronous orbit varies with altitude. Rocket launches to Sun-synchronous orbits will commonly have a fixed altitude (around 500 to 1000 kilometers for many Earth observation missions) and zero eccentricity. The launch provider will select the inclination the rocket is then determined by the above equation to provide the desired effect.

Plot Source
fig_sun_synchronous_inclination.py
# 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 numpy as np
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=[300, 1000],
    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_xaxes(
    tickmode='linear', tick0=300, dtick=100, title_text=r"Satellite Altitude [km]"
)
fig.update_yaxes(
    tickmode='linear', title_text=r"Inclination [deg]"
)

## Generate and plot data

# Generate range of true anomalies
alt = np.arange(300e3, 1000e3, 1e3)

# Compute and plot eccentric anomaly for range of true anomalies
for e in [0.0, 0.1, 0.3, 0.5]:
    # Take output mod 360 to wrap from 0 to 2pi
    ssi = [rastro.sun_synchronous_inclination(rastro.R_EARTH + a, e, True) for a in alt]
    fig.add_trace(go.Scatter(x=alt/1e3, y=ssi, name=f"e = {e:.1f}"))

pio.write_html(fig, file=OUTFILE, include_plotlyjs='cdn', full_html=False, auto_play=False)

  1. D. Vallado, Fundamentals of Astrodynamics and Applications (4th Ed.), 2010
    https://celestrak.com/software/vallado-sw.php 

Back to top