Orbital Properties¶
The orbits module provides functions to compute essential properties of satellite orbits, including orbital period, mean motion, periapsis/apoapsis characteristics, and specialized orbits like sun-synchronous configurations. These properties are fundamental for mission design, orbit determination, and trajectory analysis.
For complete API documentation, see Orbits API Reference.
Orbital Period¶
The orbital period \(T\) of a satellite is the time it takes to complete one full revolution around the central body. It is related to the semi-major axis \(a\) and gravitational parameter \(\mu\) by:
The orbital_period function computes the period for Earth-orbiting objects, while orbital_period_general accepts an explicit gravitational parameter for any celestial body.
import brahe as bh
bh.initialize_eop()
# Define orbit parameters
a = bh.R_EARTH + 500.0e3  # Semi-major axis (m) - LEO orbit at 500 km altitude
# Compute orbital period for Earth orbit (uses GM_EARTH internally)
period_earth = bh.orbital_period(a)
print(f"Orbital period (Earth): {period_earth:.3f} s")
print(f"Orbital period (Earth): {period_earth / 60:.3f} min")
# Compute orbital period for general body (explicit GM)
period_general = bh.orbital_period_general(a, bh.GM_EARTH)
print(f"Orbital period (general): {period_general:.3f} s")
# Verify they match
print(f"Difference: {abs(period_earth - period_general):.2e} s")
# Example with approximate GEO altitude
a_geo = bh.R_EARTH + 35786e3
period_geo = bh.orbital_period(a_geo)
print(f"\nGEO orbital period: {period_geo / 3600:.3f} hours")
# Expected output:
# Orbital period (Earth): 5676.977 s
# Orbital period (Earth): 94.616 min
# Orbital period (general): 5676.977 s
# Difference: 0.00e0 s
# GEO orbital period: 23.934 hours
use brahe as bh;
fn main() {
    bh::initialize_eop().unwrap();
    // Define orbit parameters
    let a = bh::constants::R_EARTH + 500.0e3; // Semi-major axis (m) - LEO orbit at 500 km altitude
    // Compute orbital period for Earth orbit (uses GM_EARTH internally)
    let period_earth = bh::orbits::orbital_period(a);
    println!("Orbital period (Earth): {:.3} s", period_earth);
    println!("Orbital period (Earth): {:.3} min", period_earth / 60.0);
    // Compute orbital period for general body (explicit GM)
    let period_general = bh::orbits::orbital_period_general(a, bh::constants::GM_EARTH);
    println!("Orbital period (general): {:.3} s", period_general);
    // Verify they match
    println!("Difference: {:.2e} s", (period_earth - period_general).abs());
    // Example with approximate GEO altitude
    let a_geo = bh::constants::R_EARTH + 35786e3;
    let period_geo = bh::orbits::orbital_period(a_geo);
    println!("\nGEO orbital period: {:.3} hours", period_geo / 3600.0);
    // Expected output:
    // Orbital period (Earth): 5676.977 s
    // Orbital period (Earth): 94.616 min
    // Orbital period (general): 5676.977 s
    // Difference: 0.00e0 s
    // GEO orbital period: 23.934 hours
}
The plot below shows how orbital period and velocity vary with altitude for circular Earth orbits:
Plot Source
import os
import pathlib
import sys
import numpy as np
import plotly.graph_objects as go
from plotly.subplots import make_subplots
import brahe as bh
# Add plots directory to path for importing brahe_theme
sys.path.insert(0, str(pathlib.Path(__file__).parent))
from brahe_theme import get_theme_colors, save_themed_html
# Configuration
SCRIPT_NAME = pathlib.Path(__file__).stem
OUTDIR = pathlib.Path(os.getenv("BRAHE_FIGURE_OUTPUT_DIR", "./docs/figures/"))
# Ensure output directory exists
os.makedirs(OUTDIR, exist_ok=True)
# Generate data
# Generate range of altitudes from 0 to 40,000 km in 500 km increments
alt = np.arange(0, 41000 * 1e3, 500 * 1e3)
# Compute velocity over altitude (km/s)
vp = [bh.perigee_velocity(bh.R_EARTH + a, 0.0) / 1e3 for a in alt]
# Compute orbital period over altitude (hours)
period = [bh.orbital_period(bh.R_EARTH + a) / 3600 for a in alt]
# Create figure with theme support
def create_figure(theme):
    """Create figure with theme-specific colors."""
    colors = get_theme_colors(theme)
    # Create subplot with secondary y-axis
    fig = make_subplots(specs=[[{"secondary_y": True}]])
    # Add velocity trace (primary y-axis)
    fig.add_trace(
        go.Scatter(
            x=alt / 1e6,
            y=vp,
            mode="lines",
            line=dict(color=colors["primary"], width=2),
            name="Velocity",
            showlegend=True,
        ),
        secondary_y=False,
    )
    # Add orbital period trace (secondary y-axis)
    fig.add_trace(
        go.Scatter(
            x=alt / 1e6,
            y=period,
            mode="lines",
            line=dict(color=colors["secondary"], width=2),
            name="Orbital Period",
            showlegend=True,
        ),
        secondary_y=True,
    )
    # Configure primary x-axis
    fig.update_xaxes(
        tickmode="linear",
        tick0=0,
        dtick=5,
        title_text="Satellite Altitude [1000 km]",
        range=[0, 40],
        showgrid=False,
    )
    # Configure primary y-axis (velocity)
    fig.update_yaxes(
        tickmode="linear",
        tick0=0,
        dtick=1,
        title_text="Velocity [km/s]",
        range=[0, 10],
        showgrid=False,
        secondary_y=False,
    )
    # Configure secondary y-axis (period)
    fig.update_yaxes(
        tickmode="linear",
        tick0=0,
        dtick=5,
        title_text="Orbital Period [hours]",
        range=[0, 30],
        showgrid=False,
        secondary_y=True,
    )
    return fig
# Generate and save both themed versions
light_path, dark_path = save_themed_html(create_figure, OUTDIR / SCRIPT_NAME)
print(f"✓ Generated {light_path}")
print(f"✓ Generated {dark_path}")
From State Vector¶
When orbital elements are unknown but you have a Cartesian state vector, orbital_period_from_state computes the period directly from position and velocity:
import brahe as bh
import numpy as np
bh.initialize_eop()
# Define orbital elements for a LEO satellite
a = bh.R_EARTH + 500.0e3  # Semi-major axis (m)
e = 0.01  # Eccentricity
i = 97.8  # Inclination (degrees)
raan = 15.0  # Right ascension of ascending node (degrees)
argp = 30.0  # Argument of periapsis (degrees)
nu = 45.0  # True anomaly (degrees)
# Convert to Cartesian state
oe = np.array([a, e, i, raan, argp, nu])
state_eci = bh.state_osculating_to_cartesian(oe, bh.AngleFormat.DEGREES)
print("ECI State (position in km, velocity in km/s):")
print(
    f"  r = [{state_eci[0] / 1e3:.3f}, {state_eci[1] / 1e3:.3f}, {state_eci[2] / 1e3:.3f}] km"
)
print(
    f"  v = [{state_eci[3] / 1e3:.3f}, {state_eci[4] / 1e3:.3f}, {state_eci[5] / 1e3:.3f}] km/s"
)
# Compute orbital period from state vector
period = bh.orbital_period_from_state(state_eci, bh.GM_EARTH)
print(f"\nOrbital period from state: {period:.3f} s")
print(f"Orbital period from state: {period / 60:.3f} min")
# Verify against period computed from semi-major axis
period_from_sma = bh.orbital_period(a)
print(f"\nOrbital period from SMA: {period_from_sma:.3f} s")
print(f"Difference: {abs(period - period_from_sma):.2e} s")
# Expected output:
# ECI State (position in km, velocity in km/s):
#   r = [1848.964, -434.937, 6560.411] km
#   v = [-7.098, -2.173, 1.913] km/s
# Orbital period from state: 5676.977 s
# Orbital period from state: 94.616 min
# Orbital period from SMA: 5676.977 s
# Difference: 3.64e-12 s
use brahe as bh;
use nalgebra as na;
fn main() {
    bh::initialize_eop().unwrap();
    // Define orbital elements for a LEO satellite
    let a = bh::constants::R_EARTH + 500.0e3; // Semi-major axis (m)
    let e = 0.01; // Eccentricity
    let i = 97.8; // Inclination (degrees)
    let raan = 15.0; // Right ascension of ascending node (degrees)
    let argp = 30.0; // Argument of periapsis (degrees)
    let nu = 45.0; // True anomaly (degrees)
    // Convert to Cartesian state
    let oe = na::SVector::<f64, 6>::new(a, e, i, raan, argp, nu);
    let state_eci = bh::state_osculating_to_cartesian(oe, bh::constants::AngleFormat::Degrees);
    println!("ECI State (position in km, velocity in km/s):");
    println!("  r = [{:.3}, {:.3}, {:.3}] km", state_eci[0]/1e3, state_eci[1]/1e3, state_eci[2]/1e3);
    println!("  v = [{:.3}, {:.3}, {:.3}] km/s", state_eci[3]/1e3, state_eci[4]/1e3, state_eci[5]/1e3);
    // Compute orbital period from state vector
    let period = bh::orbits::orbital_period_from_state(&state_eci, bh::constants::GM_EARTH);
    println!("\nOrbital period from state: {:.3} s", period);
    println!("Orbital period from state: {:.3} min", period / 60.0);
    // Verify against period computed from semi-major axis
    let period_from_sma = bh::orbits::orbital_period(a);
    println!("\nOrbital period from SMA: {:.3} s", period_from_sma);
    println!("Difference: {:.2e} s", (period - period_from_sma).abs());
    // Expected output:
    // ECI State (position in km, velocity in km/s):
    // r = [1848.964, -434.937, 6560.411] km
    // v = [-7.098, -2.173, 1.913] km/s
    // Orbital period from state: 5676.977 s
    // Orbital period from state: 94.616 min
    // Orbital period from SMA: 5676.977 s
    // Difference: 3.64e-12 s
}
Semi-major Axis from Period¶
The inverse relationship allows computing semi-major axis when orbital period is known (useful for mission design):
import brahe as bh
bh.initialize_eop()
# Example 1: LEO satellite with 98-minute period
period_leo = 98 * 60  # 98 minutes in seconds
a_leo = bh.semimajor_axis_from_orbital_period(period_leo)
altitude_leo = a_leo - bh.R_EARTH
print("LEO Satellite (98 min period):")
print(f"  Semi-major axis: {a_leo:.3f} m")
print(f"  Altitude: {altitude_leo / 1e3:.3f} km")
# Example 2: Geosynchronous orbit (24-hour period)
period_geo = 24 * 3600  # 24 hours in seconds
a_geo = bh.semimajor_axis_from_orbital_period(period_geo)
altitude_geo = a_geo - bh.R_EARTH
print("\nGeosynchronous Orbit (24 hour period):")
print(f"  Semi-major axis: {a_geo:.3f} m")
print(f"  Altitude: {altitude_geo / 1e3:.3f} km")
# Example 3: Using general function for Moon orbit
period_moon = 27.3 * 24 * 3600  # 27.3 days in seconds
a_moon = bh.semimajor_axis_from_orbital_period_general(period_moon, bh.GM_EARTH)
print("\nMoon's orbit (27.3 day period):")
print(f"  Semi-major axis: {a_moon / 1e3:.3f} km")
# Verify round-trip conversion
period_verify = bh.orbital_period(a_leo)
print("\nRound-trip verification:")
print(f"  Original period: {period_leo:.3f} s")
print(f"  Computed period: {period_verify:.3f} s")
print(f"  Difference: {abs(period_leo - period_verify):.2e} s")
# Expected output:
# LEO Satellite (98 min period):
#   Semi-major axis: 7041160.278 m
#   Altitude: 663.024 km
# Geosynchronous Orbit (24 hour period):
#   Semi-major axis: 42241095.664 m
#   Altitude: 35862.959 km
# Moon's orbit (27.3 day period):
#   Semi-major axis: 382980.745 km
# Round-trip verification:
#   Original period: 5880.000 s
#   Computed period: 5880.000 s
#   Difference: 8.19e-12 s
use brahe as bh;
fn main() {
    bh::initialize_eop().unwrap();
    // Example 1: LEO satellite with 98-minute period
    let period_leo = 98.0 * 60.0; // 98 minutes in seconds
    let a_leo = bh::orbits::semimajor_axis_from_orbital_period(period_leo);
    let altitude_leo = a_leo - bh::constants::R_EARTH;
    println!("LEO Satellite (98 min period):");
    println!("  Semi-major axis: {:.3} m", a_leo);
    println!("  Altitude: {:.3} km", altitude_leo / 1e3);
    // Example 2: Geosynchronous orbit (24-hour period)
    let period_geo = 24.0 * 3600.0; // 24 hours in seconds
    let a_geo = bh::orbits::semimajor_axis_from_orbital_period(period_geo);
    let altitude_geo = a_geo - bh::constants::R_EARTH;
    println!("\nGeosynchronous Orbit (24 hour period):");
    println!("  Semi-major axis: {:.3} m", a_geo);
    println!("  Altitude: {:.3} km", altitude_geo / 1e3);
    // Example 3: Using general function for Moon orbit
    let period_moon = 27.3 * 24.0 * 3600.0; // 27.3 days in seconds
    let a_moon = bh::orbits::semimajor_axis_from_orbital_period_general(period_moon, bh::constants::GM_EARTH);
    println!("\nMoon's orbit (27.3 day period):");
    println!("  Semi-major axis: {:.3} km", a_moon / 1e3);
    // Verify round-trip conversion
    let period_verify = bh::orbits::orbital_period(a_leo);
    println!("\nRound-trip verification:");
    println!("  Original period: {:.3} s", period_leo);
    println!("  Computed period: {:.3} s", period_verify);
    println!("  Difference: {:.2e} s", (period_leo - period_verify).abs());
    // Expected output:
    // LEO Satellite (98 min period):
    //   Semi-major axis: 7041160.278 m
    //   Altitude: 663.024 km
    // Geosynchronous Orbit (24 hour period):
    //   Semi-major axis: 42241095.664 m
    //   Altitude: 35862.959 km
    // Moon's orbit (27.3 day period):
    //   Semi-major axis: 382980.745 km
    // Round-trip verification:
    //   Original period: 5880.000 s
    //   Computed period: 5880.000 s
    //   Difference: 8.19e-12 s
}
Mean Motion¶
A satellite's average angular rate over one orbit is its mean motion \(n\), calculated from the semi-major axis and gravitational parameter:
The mean_motion function computes this for Earth-orbiting objects, while mean_motion_general works for any celestial body. Both functions support output in radians or degrees per second via the angle_format parameter.
import brahe as bh
bh.initialize_eop()
# Define orbit parameters
a_leo = bh.R_EARTH + 500.0e3  # LEO satellite at 500 km altitude
a_geo = bh.R_EARTH + 35786e3  # GEO satellite
# Compute mean motion in radians/s (Earth-specific)
n_leo_rad = bh.mean_motion(a_leo, bh.AngleFormat.RADIANS)
n_geo_rad = bh.mean_motion(a_geo, bh.AngleFormat.RADIANS)
print("Mean Motion in radians/second:")
print(f"  LEO (500 km): {n_leo_rad:.6f} rad/s")
print(f"  GEO:          {n_geo_rad:.6f} rad/s")
# Compute mean motion in degrees/s
n_leo_deg = bh.mean_motion(a_leo, bh.AngleFormat.DEGREES)
n_geo_deg = bh.mean_motion(a_geo, bh.AngleFormat.DEGREES)
print("\nMean Motion in degrees/second:")
print(f"  LEO (500 km): {n_leo_deg:.6f} deg/s")
print(f"  GEO:          {n_geo_deg:.6f} deg/s")
# Convert to degrees/day (common unit for TLEs)
print("\nMean Motion in degrees/day:")
print(f"  LEO (500 km): {n_leo_deg * 86400:.3f} deg/day")
print(f"  GEO:          {n_geo_deg * 86400:.3f} deg/day")
# Verify using general function
n_leo_general = bh.mean_motion_general(a_leo, bh.GM_EARTH, bh.AngleFormat.RADIANS)
print(f"\nVerification (general function): {n_leo_general:.6f} rad/s")
print(f"Difference: {abs(n_leo_rad - n_leo_general):.2e} rad/s")
# Expected output:
# Mean Motion in radians/second:
#   LEO (500 km): 0.001107 rad/s
#   GEO:          0.000073 rad/s
# Mean Motion in degrees/second:
#   LEO (500 km): 0.063414 deg/s
#   GEO:          0.004178 deg/s
# Mean Motion in degrees/day:
#   LEO (500 km): 5478.972 deg/day
#   GEO:          360.986 deg/day
# Verification (general function): 0.001107 rad/s
# Difference: 0.00e+00 rad/s
use brahe as bh;
fn main() {
    bh::initialize_eop().unwrap();
    // Define orbit parameters
    let a_leo = bh::constants::R_EARTH + 500.0e3; // LEO satellite at 500 km altitude
    let a_geo = bh::constants::R_EARTH + 35786e3; // GEO satellite
    // Compute mean motion in radians/s (Earth-specific)
    let n_leo_rad = bh::orbits::mean_motion(a_leo, bh::constants::AngleFormat::Radians);
    let n_geo_rad = bh::orbits::mean_motion(a_geo, bh::constants::AngleFormat::Radians);
    println!("Mean Motion in radians/second:");
    println!("  LEO (500 km): {:.6} rad/s", n_leo_rad);
    println!("  GEO:          {:.6} rad/s", n_geo_rad);
    // Compute mean motion in degrees/s
    let n_leo_deg = bh::orbits::mean_motion(a_leo, bh::constants::AngleFormat::Degrees);
    let n_geo_deg = bh::orbits::mean_motion(a_geo, bh::constants::AngleFormat::Degrees);
    println!("\nMean Motion in degrees/second:");
    println!("  LEO (500 km): {:.6} deg/s", n_leo_deg);
    println!("  GEO:          {:.6} deg/s", n_geo_deg);
    // Convert to degrees/day (common unit for TLEs)
    println!("\nMean Motion in degrees/day:");
    println!("  LEO (500 km): {:.3} deg/day", n_leo_deg * 86400.0);
    println!("  GEO:          {:.3} deg/day", n_geo_deg * 86400.0);
    // Verify using general function
    let n_leo_general = bh::orbits::mean_motion_general(a_leo, bh::constants::GM_EARTH, bh::constants::AngleFormat::Radians);
    println!("\nVerification (general function): {:.6} rad/s", n_leo_general);
    println!("Difference: {:.2e} rad/s", (n_leo_rad - n_leo_general).abs());
    // Expected output:
    // Mean Motion in radians/second:
    //   LEO (500 km): 0.001107 rad/s
    //   GEO:          0.000073 rad/s
    // Mean Motion in degrees/second:
    //   LEO (500 km): 0.063414 deg/s
    //   GEO:          0.004178 deg/s
    // Mean Motion in degrees/day:
    //   LEO (500 km): 5478.972 deg/day
    //   GEO:          360.986 deg/day
    // Verification (general function): 0.001107 rad/s
    // Difference: 0.00e+00 rad/s
}
Semi-major Axis from Mean Motion¶
Since orbital data formats like TLEs specify mean motion instead of semi-major axis, the inverse computation is essential:
import brahe as bh
bh.initialize_eop()
# Example 1: ISS-like orbit with ~15.5 revolutions per day
n_iss = 15.5 * 360.0 / 86400.0  # Convert revs/day to deg/s
a_iss = bh.semimajor_axis(n_iss, bh.AngleFormat.DEGREES)
altitude_iss = a_iss - bh.R_EARTH
print("ISS-like Orbit (15.5 revs/day):")
print(f"  Mean motion: {n_iss:.6f} deg/s")
print(f"  Semi-major axis: {a_iss:.3f} m")
print(f"  Altitude: {altitude_iss / 1e3:.3f} km")
# Example 2: Geosynchronous orbit (1 revolution per day)
n_geo = 1.0 * 360.0 / 86400.0  # 1 rev/day in deg/s
a_geo = bh.semimajor_axis(n_geo, bh.AngleFormat.DEGREES)
altitude_geo = a_geo - bh.R_EARTH
print("\nGeosynchronous Orbit (1 rev/day):")
print(f"  Mean motion: {n_geo:.6f} deg/s")
print(f"  Semi-major axis: {a_geo:.3f} m")
print(f"  Altitude: {altitude_geo / 1e3:.3f} km")
# Example 3: Using radians
n_leo_rad = 0.001  # rad/s
a_leo = bh.semimajor_axis(n_leo_rad, bh.AngleFormat.RADIANS)
print("\nLEO from radians/s:")
print(f"  Mean motion: {n_leo_rad:.6f} rad/s")
print(f"  Semi-major axis: {a_leo:.3f} m")
print(f"  Altitude: {(a_leo - bh.R_EARTH) / 1e3:.3f} km")
# Verify round-trip conversion
n_verify = bh.mean_motion(a_iss, bh.AngleFormat.DEGREES)
print("\nRound-trip verification:")
print(f"  Original mean motion: {n_iss:.6f} deg/s")
print(f"  Computed mean motion: {n_verify:.6f} deg/s")
print(f"  Difference: {abs(n_iss - n_verify):.2e} deg/s")
# Expected output:
# ISS-like Orbit (15.5 revs/day):
#   Mean motion: 0.064583 deg/s
#   Semi-major axis: 6794863.068 m
#   Altitude: 416.727 km
# Geosynchronous Orbit (1 rev/day):
#   Mean motion: 0.004167 deg/s
#   Semi-major axis: 42241095.664 m
#   Altitude: 35862.959 km
# LEO from radians/s:
#   Mean motion: 0.001000 rad/s
#   Semi-major axis: 7359459.593 m
#   Altitude: 981.323 km
# Round-trip verification:
#   Original mean motion: 0.064583 deg/s
#   Computed mean motion: 0.064583 deg/s
#   Difference: 9.71e-17 deg/s
use brahe as bh;
fn main() {
    bh::initialize_eop().unwrap();
    // Example 1: ISS-like orbit with ~15.5 revolutions per day
    let n_iss = 15.5 * 360.0 / 86400.0; // Convert revs/day to deg/s
    let a_iss = bh::orbits::semimajor_axis(n_iss, bh::constants::AngleFormat::Degrees);
    let altitude_iss = a_iss - bh::constants::R_EARTH;
    println!("ISS-like Orbit (15.5 revs/day):");
    println!("  Mean motion: {:.6} deg/s", n_iss);
    println!("  Semi-major axis: {:.3} m", a_iss);
    println!("  Altitude: {:.3} km", altitude_iss / 1e3);
    // Example 2: Geosynchronous orbit (1 revolution per day)
    let n_geo = 1.0 * 360.0 / 86400.0; // 1 rev/day in deg/s
    let a_geo = bh::orbits::semimajor_axis(n_geo, bh::constants::AngleFormat::Degrees);
    let altitude_geo = a_geo - bh::constants::R_EARTH;
    println!("\nGeosynchronous Orbit (1 rev/day):");
    println!("  Mean motion: {:.6} deg/s", n_geo);
    println!("  Semi-major axis: {:.3} m", a_geo);
    println!("  Altitude: {:.3} km", altitude_geo / 1e3);
    // Example 3: Using radians
    let n_leo_rad = 0.001; // rad/s
    let a_leo = bh::orbits::semimajor_axis(n_leo_rad, bh::constants::AngleFormat::Radians);
    println!("\nLEO from radians/s:");
    println!("  Mean motion: {:.6} rad/s", n_leo_rad);
    println!("  Semi-major axis: {:.3} m", a_leo);
    println!("  Altitude: {:.3} km", (a_leo - bh::constants::R_EARTH) / 1e3);
    // Verify round-trip conversion
    let n_verify = bh::orbits::mean_motion(a_iss, bh::constants::AngleFormat::Degrees);
    println!("\nRound-trip verification:");
    println!("  Original mean motion: {:.6} deg/s", n_iss);
    println!("  Computed mean motion: {:.6} deg/s", n_verify);
    println!("  Difference: {:.2e} deg/s", (n_iss - n_verify).abs());
    // Expected output:
    // ISS-like Orbit (15.5 revs/day):
    //   Mean motion: 0.064583 deg/s
    //   Semi-major axis: 6794863.068 m
    //   Altitude: 416.727 km
    // Geosynchronous Orbit (1 rev/day):
    //   Mean motion: 0.004167 deg/s
    //   Semi-major axis: 42241095.664 m
    //   Altitude: 35862.959 km
    // LEO from radians/s:
    //   Mean motion: 0.001000 rad/s
    //   Semi-major axis: 7359459.593 m
    //   Altitude: 981.323 km
    // Round-trip verification:
    //   Original mean motion: 0.064583 deg/s
    //   Computed mean motion: 0.064583 deg/s
    //   Difference: 9.71e-17 deg/s
}
Periapsis Properties¶
The periapsis is the point of closest approach to the central body, where orbital velocity is greatest.
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.
Therefore periapsis is the point of closest approach of the orbiting body with respect to its central body. The suffix can be modified to indicate the closest approach to a specific celestial body: perigee for Earth, perihelion for the Sun.
Brahe provides functions to compute periapsis velocity, distance, and altitude based on orbital elements.
Velocity¶
The periapsis velocity is given by:
where \(\mu\) is the gravitational parameter, \(a\) is the semi-major axis, and \(e\) is the eccentricity.
Distance¶
The periapsis distance from the center of the central body is (from Vallado1 Equation 2-75):
Altitude¶
The periapsis altitude is the height above the surface of the central body:
where \(R_{body}\) is the radius of the central body. For Earth orbits, the perigee_altitude function provides a convenient wrapper using \(R_{\oplus}\).
Code Example¶
import brahe as bh
bh.initialize_eop()
# Define orbit parameters
a = bh.R_EARTH + 500.0e3  # Semi-major axis (m)
e = 0.01  # Eccentricity
# Compute periapsis velocity (generic)
periapsis_velocity = bh.periapsis_velocity(a, e, bh.GM_EARTH)
print(f"Periapsis velocity: {periapsis_velocity:.3f} m/s")
# Compute as a perigee velocity (Earth-specific)
perigee_velocity = bh.perigee_velocity(a, e)
print(f"Perigee velocity:   {perigee_velocity:.3f} m/s")
# Compute periapsis distance
periapsis_distance = bh.periapsis_distance(a, e)
print(f"Periapsis distance: {periapsis_distance / 1e3:.3f} km")
# Compute periapsis altitude (generic)
periapsis_altitude = bh.periapsis_altitude(a, e, bh.R_EARTH)
print(f"Periapsis altitude: {periapsis_altitude / 1e3:.3f} km")
# Compute as a perigee altitude (Earth-specific)
perigee_altitude = bh.perigee_altitude(a, e)
print(f"Perigee altitude:   {perigee_altitude / 1e3:.3f} km")
# Expected output:
# Periapsis velocity: 7689.119 m/s
# Perigee velocity:   7689.119 m/s
# Periapsis distance: 6809.355 km
# Periapsis altitude: 431.219 km
# Perigee altitude:   431.219 km
use brahe as bh;
fn main() {
    bh::initialize_eop().unwrap();
    // Define orbit parameters
    let a = bh::constants::R_EARTH + 500.0e3; // Semi-major axis (m)
    let e = 0.01; // Eccentricity
    // Compute periapsis velocity (generic)
    let periapsis_velocity = bh::orbits::periapsis_velocity(a, e, bh::constants::GM_EARTH);
    println!("Periapsis velocity: {:.3} m/s", periapsis_velocity);
    // Compute as a perigee velocity (Earth-specific)
    let perigee_velocity = bh::orbits::perigee_velocity(a, e);
    println!("Perigee velocity:   {:.3} m/s", perigee_velocity);
    // Compute periapsis distance
    let periapsis_distance = bh::orbits::periapsis_distance(a, e);
    println!("Periapsis distance: {:.3} km", periapsis_distance / 1e3);
    // Compute periapsis altitude (generic)
    let periapsis_altitude = bh::orbits::periapsis_altitude(a, e, bh::constants::R_EARTH);
    println!("Periapsis altitude: {:.3} km", periapsis_altitude / 1e3);
    // Compute as a perigee altitude (Earth-specific)
    let perigee_altitude = bh::orbits::perigee_altitude(a, e);
    println!("Perigee altitude:   {:.3} km", perigee_altitude / 1e3);
    // Expected output:
    // Periapsis velocity: 7689.119 m/s
    // Perigee velocity:   7689.119 m/s
    // Periapsis distance: 6809.355 km
    // Periapsis altitude: 431.219 km
    // Perigee altitude:   431.219 km
}
Apoapsis Properties¶
The apoapsis is the farthest point from the central body, where orbital velocity is lowest.
Info
The word apoapsis is formed by combination of the Greek words "apo-" (meaning away from, separate, or apart from) and "apsis".
Therefore apoapsis is the farthest point of an orbiting body with respect to its central body. The suffix can be modified to indicate the farthest point from a specific celestial body: apogee for Earth, aphelion for the Sun.
Brahe provides functions to compute apoapsis velocity, distance, and altitude based on orbital elements.
Warning
Apoapsis position, velocity, and altitude are only defined for elliptic and circular orbits. For parabolic and hyperbolic orbits, these quantities are undefined.
Velocity¶
The apoapsis velocity is given by:
Distance¶
The apoapsis distance from the center of the central body is:
Altitude¶
The apoapsis altitude is the height above the surface of the central body:
where \(R_{body}\) is the radius of the central body. For Earth orbits, the apogee_altitude function provides a convenient wrapper using \(R_{\oplus}\).
Code Example¶
import brahe as bh
bh.initialize_eop()
# Define orbit parameters
a = bh.R_EARTH + 500.0e3  # Semi-major axis (m)
e = 0.01  # Eccentricity
# Compute apoapsis velocity (generic)
apoapsis_velocity = bh.apoapsis_velocity(a, e, bh.GM_EARTH)
print(f"Apoapsis velocity: {apoapsis_velocity:.3f} m/s")
# Compute as an apogee velocity (Earth-specific)
apogee_velocity = bh.apogee_velocity(a, e)
print(f"Apogee velocity:   {apogee_velocity:.3f} m/s")
# Compute apoapsis distance
apoapsis_distance = bh.apoapsis_distance(a, e)
print(f"Apoapsis distance: {apoapsis_distance / 1e3:.3f} km")
# Compute apoapsis altitude (generic)
apoapsis_altitude = bh.apoapsis_altitude(a, e, bh.R_EARTH)
print(f"Apoapsis altitude: {apoapsis_altitude / 1e3:.3f} km")
# Compute as an apogee altitude (Earth-specific)
apogee_altitude = bh.apogee_altitude(a, e)
print(f"Apogee altitude:   {apogee_altitude / 1e3:.3f} km")
# Expected output:
# Apoapsis velocity: 7536.859 m/s
# Apogee velocity:   7536.859 m/s
# Apoapsis distance: 6946.918 km
# Apoapsis altitude: 568.781 km
# Apogee altitude:   568.781 km
use brahe as bh;
fn main() {
    bh::initialize_eop().unwrap();
    // Define orbit parameters
    let a = bh::constants::R_EARTH + 500.0e3; // Semi-major axis (m)
    let e = 0.01; // Eccentricity
    // Compute apoapsis velocity (generic)
    let apoapsis_velocity = bh::orbits::apoapsis_velocity(a, e, bh::constants::GM_EARTH);
    println!("Apoapsis velocity: {:.3} m/s", apoapsis_velocity);
    // Compute as an apogee velocity (Earth-specific)
    let apogee_velocity = bh::orbits::apogee_velocity(a, e);
    println!("Apogee velocity:   {:.3} m/s", apogee_velocity);
    // Compute apoapsis distance
    let apoapsis_distance = bh::orbits::apoapsis_distance(a, e);
    println!("Apoapsis distance: {:.3} km", apoapsis_distance / 1e3);
    // Compute apoapsis altitude (generic)
    let apoapsis_altitude = bh::orbits::apoapsis_altitude(a, e, bh::constants::R_EARTH);
    println!("Apoapsis altitude: {:.3} km", apoapsis_altitude / 1e3);
    // Compute as an apogee altitude (Earth-specific)
    let apogee_altitude = bh::orbits::apogee_altitude(a, e);
    println!("Apogee altitude:   {:.3} km", apogee_altitude / 1e3);
    // Expected output:
    // Apoapsis velocity: 7536.859 m/s
    // Apogee velocity:   7536.859 m/s
    // Apoapsis distance: 6946.918 km
    // Apoapsis altitude: 568.781 km
    // Apogee altitude:   568.781 km
}
Sun-Synchronous Inclination¶
A sun-synchronous orbit maintains a constant angle relative to the Sun by matching its nodal precession rate to Earth's annual revolution. The right ascension of the ascending node (\(\Omega\)) advances at the same rate as the Sun's apparent motion: approximately 0.9856°/day. This configuration is highly valuable for Earth observation satellites requiring consistent illumination conditions—a sun-synchronous satellite crosses the equator at the same local time on each pass (e.g., always at 2 PM).
Earth's oblateness, characterized by the \(J_2\) zonal harmonic, causes secular drift in \(\Omega\):
For sun-synchronicity, this must equal:
Solving for inclination as a function of semi-major axis and eccentricity:
The sun_synchronous_inclination function computes this required inclination:
"""
import brahe as bh
bh.initialize_eop()
# Example 1: Typical sun-synchronous LEO at 800 km altitude
a_leo = bh.R_EARTH + 800e3  # Semi-major axis
e_leo = 0.0  # Circular orbit
inc_leo_deg = bh.sun_synchronous_inclination(a_leo, e_leo, bh.AngleFormat.DEGREES)
inc_leo_rad = bh.sun_synchronous_inclination(a_leo, e_leo, bh.AngleFormat.RADIANS)
print("Sun-synchronous LEO (800 km, circular):")
print(f"  Inclination: {inc_leo_deg:.3f} degrees")
print(f"  Inclination: {inc_leo_rad:.6f} radians")
# Example 2: Different altitudes
altitudes = [500, 600, 700, 800, 900, 1000]  # km
print("\nSun-synchronous inclination vs altitude (circular orbits):")
for alt_km in altitudes:
    a = bh.R_EARTH + alt_km * 1e3
    inc = bh.sun_synchronous_inclination(a, 0.0, bh.AngleFormat.DEGREES)
    print(f"  {alt_km:4d} km: {inc:.3f} deg")
# Example 3: Effect of eccentricity
a_fixed = bh.R_EARTH + 700e3
eccentricities = [0.0, 0.001, 0.005, 0.01, 0.02]
print("\nSun-synchronous inclination vs eccentricity (700 km orbit):")
for e in eccentricities:
    inc = bh.sun_synchronous_inclination(a_fixed, e, bh.AngleFormat.DEGREES)
    print(f"  e = {e:.3f}: {inc:.3f} deg")
# Example 4: Practical mission example (Landsat-like)
a_landsat = bh.R_EARTH + 705e3
e_landsat = 0.0001
inc_landsat = bh.sun_synchronous_inclination(
    a_landsat, e_landsat, bh.AngleFormat.DEGREES
)
print("\nLandsat-like orbit (705 km, nearly circular):")
print(f"  Inclination: {inc_landsat:.3f} deg")
print(f"  Period: {bh.orbital_period(a_landsat) / 60:.3f} min")
# Expected output:
# Sun-synchronous LEO (800 km, circular):
#   Inclination: 98.603 degrees
#   Inclination: 1.720948 radians
# Sun-synchronous inclination vs altitude (circular orbits):
#    500 km: 97.402 deg
#    600 km: 97.788 deg
#    700 km: 98.188 deg
#    800 km: 98.603 deg
#    900 km: 99.033 deg
#   1000 km: 99.479 deg
# Sun-synchronous inclination vs eccentricity (700 km orbit):
#   e = 0.000: 98.188 deg
#   e = 0.001: 98.188 deg
#   e = 0.005: 98.187 deg
#   e = 0.010: 98.186 deg
#   e = 0.020: 98.181 deg
# Landsat-like orbit (705 km, nearly circular):
#   Inclination: 98.208 deg
#   Period: 98.878 min
#[allow(unused_imports)]
use brahe as bh;
fn main() {
    bh::initialize_eop().unwrap();
    // Example 1: Typical sun-synchronous LEO at 800 km altitude
    let a_leo = bh::constants::R_EARTH + 800e3; // Semi-major axis
    let e_leo = 0.0; // Circular orbit
    let inc_leo_deg = bh::orbits::sun_synchronous_inclination(a_leo, e_leo, bh::constants::AngleFormat::Degrees);
    let inc_leo_rad = bh::orbits::sun_synchronous_inclination(a_leo, e_leo, bh::constants::AngleFormat::Radians);
    println!("Sun-synchronous LEO (800 km, circular):");
    println!("  Inclination: {:.3} degrees", inc_leo_deg);
    println!("  Inclination: {:.6} radians", inc_leo_rad);
    // Example 2: Different altitudes
    let altitudes = [500.0, 600.0, 700.0, 800.0, 900.0, 1000.0]; // km
    println!("\nSun-synchronous inclination vs altitude (circular orbits):");
    for alt_km in altitudes.iter() {
        let a = bh::constants::R_EARTH + alt_km * 1e3;
        let inc = bh::orbits::sun_synchronous_inclination(a, 0.0, bh::constants::AngleFormat::Degrees);
        println!("  {:4} km: {:.3} deg", *alt_km as i32, inc);
    }
    // Example 3: Effect of eccentricity
    let a_fixed = bh::constants::R_EARTH + 700e3;
    let eccentricities = [0.0, 0.001, 0.005, 0.01, 0.02];
    println!("\nSun-synchronous inclination vs eccentricity (700 km orbit):");
    for e in eccentricities.iter() {
        let inc = bh::orbits::sun_synchronous_inclination(a_fixed, *e, bh::constants::AngleFormat::Degrees);
        println!("  e = {:.3}: {:.3} deg", e, inc);
    }
    // Example 4: Practical mission example (Landsat-like)
    let a_landsat = bh::constants::R_EARTH + 705e3;
    let e_landsat = 0.0001;
    let inc_landsat = bh::orbits::sun_synchronous_inclination(a_landsat, e_landsat, bh::constants::AngleFormat::Degrees);
    println!("\nLandsat-like orbit (705 km, nearly circular):");
    println!("  Inclination: {:.3} deg", inc_landsat);
    println!("  Period: {:.3} min", bh::orbits::orbital_period(a_landsat) / 60.0);
    // Expected output:
    // Sun-synchronous LEO (800 km, circular):
    //   Inclination: 98.603 degrees
    //   Inclination: 1.720948 radians
    // Sun-synchronous inclination vs altitude (circular orbits):
    //    500 km: 97.402 deg
    //    600 km: 97.788 deg
    //    700 km: 98.188 deg
    //    800 km: 98.603 deg
    //    900 km: 99.033 deg
    //   1000 km: 99.479 deg
    // Sun-synchronous inclination vs eccentricity (700 km orbit):
    //   e = 0.000: 98.188 deg
    //   e = 0.001: 98.188 deg
    //   e = 0.005: 98.187 deg
    //   e = 0.010: 98.186 deg
    //   e = 0.020: 98.181 deg
    // Landsat-like orbit (705 km, nearly circular):
    //   Inclination: 98.208 deg
    //   Period: 98.878 min
}
The plot below shows how the required inclination varies with altitude for sun-synchronous orbits:
Plot Source
import os
import pathlib
import sys
import numpy as np
import plotly.graph_objects as go
import brahe as bh
# Add plots directory to path for importing brahe_theme
sys.path.insert(0, str(pathlib.Path(__file__).parent))
from brahe_theme import get_theme_colors, save_themed_html
# Configuration
SCRIPT_NAME = pathlib.Path(__file__).stem
OUTDIR = pathlib.Path(os.getenv("BRAHE_FIGURE_OUTPUT_DIR", "./docs/figures/"))
# Ensure output directory exists
os.makedirs(OUTDIR, exist_ok=True)
# Generate data
# Generate range of altitudes from 300 to 1000 km in 1 km increments
alt = np.arange(300e3, 1000e3, 1e3)
# Compute sun-synchronous inclination for range of eccentricities
eccentricities = [0.0, 0.1, 0.3, 0.5]
ssi_data = {}
for e in eccentricities:
    ssi_data[e] = [
        bh.sun_synchronous_inclination(bh.R_EARTH + a, e, bh.AngleFormat.DEGREES)
        for a in alt
    ]
# Create figure with theme support
def create_figure(theme):
    """Create figure with theme-specific colors."""
    colors = get_theme_colors(theme)
    fig = go.Figure()
    # Color palette for different eccentricities
    color_palette = [
        colors["primary"],
        colors["secondary"],
        colors["accent"],
        colors["error"],
    ]
    # Add traces for each eccentricity
    for i, e in enumerate(eccentricities):
        fig.add_trace(
            go.Scatter(
                x=alt / 1e3,
                y=ssi_data[e],
                mode="lines",
                line=dict(color=color_palette[i % len(color_palette)], width=2),
                name=f"e = {e:.1f}",
                showlegend=True,
            )
        )
    # Configure axes
    fig.update_xaxes(
        tickmode="linear",
        tick0=300,
        dtick=100,
        title_text="Satellite Altitude [km]",
        range=[300, 1000],
        showgrid=False,
    )
    fig.update_yaxes(
        tickmode="linear",
        title_text="Inclination [deg]",
        showgrid=False,
    )
    return fig
# Generate and save both themed versions
light_path, dark_path = save_themed_html(create_figure, OUTDIR / SCRIPT_NAME)
print(f"✓ Generated {light_path}")
print(f"✓ Generated {dark_path}")
Most sun-synchronous Earth observation missions operate at altitudes between 500-1000 km with near-zero eccentricity. The launch provider selects the precise inclination based on the above equation to achieve the desired sun-synchronous behavior.
See Also¶
- Orbits API Reference - Complete Python API documentation
- Anomaly Conversions - Converting between true, eccentric, and mean anomaly
- 
D. Vallado, Fundamentals of Astrodynamics and Applications (4th Ed.), 2010 ↩