Skip to content

Coordinate Transformations

In dynamics we are primarily concerned with the object state, which encodes the location and often the rate of change of location of an object. This state is often expressed in a particular coordinate system, such as Cartesian coordinates (position and velocity) or Keplerian orbital elements. It also supports transforming to topocentric coordinates such as Geodetic coordinates. Brahe provides functions for transforming between different coordinate systems.

Keplerian and Cartesian Transformations

We can take a state in Keplerian orbital elements and transform it into Cartesian orbital elements using state_koe_to_eci. We can invert this transformation using state_eci_to_koe to get back the original Keplerian orbital elements.

import brahe as bh
import numpy as np

# Initialize EOP
bh.initialize_eop()

# Initialize a Keplerian state
# Define orbital elements [a, e, i, Ω, ω, M] in meters and degrees
# LEO satellite: 500 km altitude, 97.8° inclination (approx sun-synchronous)
oe_deg = np.array(
    [
        bh.R_EARTH + 500e3,  # Semi-major axis (m)
        0.01,  # Eccentricity
        97.8,  # Inclination (deg)
        15.0,  # Right ascension of ascending node (deg)
        30.0,  # Argument of periapsis (deg)
        45.0,  # Mean anomaly (deg)
    ]
)

# Convert orbital elements to Cartesian state using degrees
x_deg = bh.state_koe_to_eci(oe_deg, bh.AngleFormat.DEGREES)

# Convert back to degrees
oe_deg_2 = bh.state_eci_to_koe(x_deg, bh.AngleFormat.DEGREES)

print("Original Keplerian elements:")
for i, elem in enumerate(oe_deg):
    print(f"  [{i}]: {elem:.3f}")
print("Converted Cartesian state:")
for i, elem in enumerate(x_deg):
    print(f"  [{i}]: {elem:.3f}")
print("Back to Keplerian elements:")
for i, elem in enumerate(oe_deg_2):
    print(f"  [{i}]: {elem:.3f}")
use brahe as bh;
use nalgebra as na;

fn main() {
    // Initialize EOP
    bh::initialize_eop().unwrap();

    // Initialize a Keplerian state
    // Define orbital elements [a, e, i, Ω, ω, M] in meters and degrees
    // LEO satellite: 500 km altitude, 97.8° inclination (approx sun-synchronous)
    let oe_deg = na::vector![
        bh::R_EARTH + 500e3,  // Semi-major axis (m)
        0.01,                  // Eccentricity
        97.8,                  // Inclination (deg)
        15.0,                  // Right ascension of ascending node (deg)
        30.0,                  // Argument of periapsis (deg)
        45.0,                  // Mean anomaly (deg)
    ];

    // Convert orbital elements to Cartesian state using degrees
    let x_deg = bh::state_koe_to_eci(oe_deg, bh::AngleFormat::Degrees);

    // Convert back to degrees
    let oe_deg_2 = bh::state_eci_to_koe(x_deg, bh::AngleFormat::Degrees);

    println!("Original Keplerian elements:");
    for (i, elem) in oe_deg.iter().enumerate() {
        println!("  [{i}]: {elem:.3}");
    }
    println!("Converted Cartesian state:");
    for (i, elem) in x_deg.iter().enumerate() {
        println!("  [{i}]: {elem:.3}");
    }
    println!("Back to Keplerian elements:");
    for (i, elem) in oe_deg_2.iter().enumerate() {
        println!("  [{i}]: {elem:.3}");
    }
}
Output
Original Keplerian elements:
  [0]: 6878136.300
  [1]: 0.010
  [2]: 97.800
  [3]: 15.000
  [4]: 30.000
  [5]: 45.000
Converted Cartesian state:
  [0]: 1848964.106
  [1]: -434937.468
  [2]: 6560410.530
  [3]: -7098.380
  [4]: -2173.345
  [5]: 1913.333
Back to Keplerian elements:
  [0]: 6878136.300
  [1]: 0.010
  [2]: 97.800
  [3]: 15.000
  [4]: 30.000
  [5]: 45.000
Original Keplerian elements:
  [0]: 6878136.300
  [1]: 0.010
  [2]: 97.800
  [3]: 15.000
  [4]: 30.000
  [5]: 45.000
Converted Cartesian state:
  [0]: 1848964.106
  [1]: -434937.468
  [2]: 6560410.530
  [3]: -7098.380
  [4]: -2173.345
  [5]: 1913.333
Back to Keplerian elements:
  [0]: 6878136.300
  [1]: 0.010
  [2]: 97.800
  [3]: 15.000
  [4]: 30.000
  [5]: 45.000

Cartesian to Geodetic

Similarly we can convert from a Cartesian state in ECEF to geodetic coordinates (longitude, latitude, altitude) using state_ecef_to_geodetic.

import brahe as bh
import numpy as np

# Initialize EOP
bh.initialize_eop()

# Create an epoch
epc = bh.Epoch(2024, 1, 1, 0, 0, 0)

# Initialize a Keplerian state
# Define orbital elements [a, e, i, Ω, ω, M] in meters and degrees
# LEO satellite: 500 km altitude, 97.8° inclination (approx sun-synchronous)
oe_deg = np.array(
    [
        bh.R_EARTH + 500e3,  # Semi-major axis (m)
        0.01,  # Eccentricity
        97.8,  # Inclination (deg)
        15.0,  # Right ascension of ascending node (deg)
        30.0,  # Argument of periapsis (deg)
        45.0,  # Mean anomaly (deg)
    ]
)

# Convert orbital elements to Cartesian state using degrees
x_deg = bh.state_koe_to_eci(oe_deg, bh.AngleFormat.DEGREES)

# Covert ECI cartesian state to ECEF cartesian state
x_ecef = bh.state_eci_to_ecef(epc, x_deg)

# Convert ECEF cartesian state to geodetic coordinates
lat, lon, alt = bh.position_ecef_to_geodetic(x_ecef[:3], angle_format=bh.AngleFormat.DEGREES)
print(f"Geodetic coordinates (lat, lon, alt): {lat:.6f}°, {lon:.6f}°, {alt/1e3:.3f} km")
use brahe as bh;
use nalgebra as na;

fn main() {
    // Initialize EOP
    bh::initialize_eop().unwrap();

    // Create an epoch
    let epc = bh::Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, bh::TimeSystem::UTC);

    // Initialize a Keplerian state
    // Define orbital elements [a, e, i, Ω, ω, M] in meters and degrees
    // LEO satellite: 500 km altitude, 97.8° inclination (approx sun-synchronous)
    let oe_deg = na::vector![
        bh::R_EARTH + 500e3,  // Semi-major axis (m)
        0.01,                  // Eccentricity
        97.8,                  // Inclination (deg)
        15.0,                  // Right ascension of ascending node (deg)
        30.0,                  // Argument of periapsis (deg)
        45.0,                  // Mean anomaly (deg)
    ];

    // Convert orbital elements to Cartesian state using degrees
    let x_deg = bh::state_koe_to_eci(oe_deg, bh::AngleFormat::Degrees);

    // Convert ECI cartesian state to ECEF cartesian state
    let x_ecef = bh::state_eci_to_ecef(epc, x_deg);

    // Convert ECEF cartesian state to geodetic coordinates
    let geodetic = bh::position_ecef_to_geodetic(x_ecef.fixed_rows::<3>(0).into(), bh::AngleFormat::Degrees);
    println!("Geodetic coordinates (lat, lon, alt): {:.6}°, {:.6}°, {:.3} km", geodetic[0], geodetic[1], geodetic[2] / 1e3);
}
Output
Geodetic coordinates (lat, lon, alt): -113.194611°, 74.076494°, 471.474 km
Geodetic coordinates (lat, lon, alt): -113.194611°, 74.076494°, 471.474 km