Skip to content

OrbitTrajectory

OrbitTrajectory is a specialized trajectory container for orbital mechanics that tracks reference frames (ECI/ECEF) and orbital representations (Cartesian/Keplerian). Unlike Trajectory which store frame-agnostic data, OrbitTrajectory understands orbital mechanics and enables automatic conversions between reference frames and representations.

Use OrbitTrajectory when:

  • Working with orbital mechanics applications
  • Need to convert between ECI and ECEF frames
  • Need to convert between Cartesian and Keplerian representations
  • Want frame/representation metadata tracked automatically
  • Working with propagators that output orbital trajectories

OrbitTrajectory implements the OrbitalTrajectory trait in addition to Trajectory and Interpolatable, providing orbital-specific functionality on top of the standard trajectory interface.

Initialization

Empty Trajectory - Cartesian Representation

For cartesian representation, the frame can be ECI or ECEF. The AngleFormat must be None for Cartesian representations

import brahe as bh

bh.initialize_eop()

# Create trajectory in ECI frame, Cartesian representation
traj_eci = bh.OrbitTrajectory(
    6,  # State dimension (position + velocity)
    bh.OrbitFrame.ECI,
    bh.OrbitRepresentation.CARTESIAN,
    None,  # No angle format for Cartesian
)
print(f"Frame (str): {traj_eci.frame}")  # Output: ECI
print(
    f"Frame (repr): {repr(traj_eci.frame)}"
)  # Output: OrbitFrame(Earth-Centered Inertial)
print(f"Representation (str): {traj_eci.representation}")  # Output: Cartesian
print(
    f"Representation (repr): {repr(traj_eci.representation)}"
)  # Output: OrbitRepresentation(Cartesian)

# Create trajectory in ECEF frame, Cartesian representation
traj_ecef = bh.OrbitTrajectory(
    6, bh.OrbitFrame.ECEF, bh.OrbitRepresentation.CARTESIAN, None
)
print(f"Frame (str): {traj_ecef.frame}")  # Output: ECEF
print(
    f"Frame (repr): {repr(traj_ecef.frame)}"
)  # Output: OrbitFrame(Earth-Centered Earth-Fixed)
use brahe as bh;
use bh::trajectories::SOrbitTrajectory;
use bh::trajectories::traits::{OrbitFrame, OrbitRepresentation};

fn main() {
    bh::initialize_eop().unwrap();

    // Create trajectory in ECI frame, Cartesian representation
    let traj_eci = SOrbitTrajectory::new(
        OrbitFrame::ECI,
        OrbitRepresentation::Cartesian,
        None
    );
    println!("Frame (Display): {}", traj_eci.frame);
    println!("Frame (Debug): {:?}", traj_eci.frame);
    println!("Representation (Display): {}", traj_eci.representation);
    println!("Representation (Debug): {:?}", traj_eci.representation);

    // Create trajectory in ECEF frame, Cartesian representation
    let traj_ecef = SOrbitTrajectory::new(
        OrbitFrame::ECEF,
        OrbitRepresentation::Cartesian,
        None
    );
    println!("Frame (Display): {}", traj_ecef.frame);
    println!("Frame (Debug): {:?}", traj_ecef.frame);
}
Output
1
2
3
4
5
6
Frame (str): ECI
Frame (repr): OrbitFrame(Earth-Centered Inertial)
Representation (str): Cartesian
Representation (repr): OrbitRepresentation(Cartesian)
Frame (str): ECEF
Frame (repr): OrbitFrame(Earth-Centered Earth-Fixed)
1
2
3
4
5
6
Frame (Display): ECI
Frame (Debug): OrbitFrame(Earth-Centered Inertial)
Representation (Display): Cartesian
Representation (Debug): OrbitRepresentation(Cartesian)
Frame (Display): ECEF
Frame (Debug): OrbitFrame(Earth-Centered Earth-Fixed)

Empty Trajectory - Keplerian Elements

To create an empty trajectory in Keplerian representation you must specify the frame as ECI and provide an AngleFormat.

import brahe as bh

bh.initialize_eop()

# Create trajectory in ECI frame, Keplerian representation with radians
traj_kep_rad = bh.OrbitTrajectory(
    6,  # State dimension (6 orbital elements)
    bh.OrbitFrame.ECI,
    bh.OrbitRepresentation.KEPLERIAN,
    bh.AngleFormat.RADIANS,  # Required for Keplerian
)

# Create trajectory in ECI frame, Keplerian representation with degrees
traj_kep_deg = bh.OrbitTrajectory(
    6, bh.OrbitFrame.ECI, bh.OrbitRepresentation.KEPLERIAN, bh.AngleFormat.DEGREES
)
use brahe as bh;
use bh::trajectories::SOrbitTrajectory;
use bh::trajectories::traits::{OrbitFrame, OrbitRepresentation};
use bh::AngleFormat;

fn main() {
    bh::initialize_eop().unwrap();

    // Create trajectory in ECI frame, Keplerian representation with radians
    let _traj_kep_rad = SOrbitTrajectory::new(
        OrbitFrame::ECI,
        OrbitRepresentation::Keplerian,
        Some(AngleFormat::Radians)
    );

    // Create trajectory in ECI frame, Keplerian representation with degrees
    let _traj_kep_deg = SOrbitTrajectory::new(
        OrbitFrame::ECI,
        OrbitRepresentation::Keplerian,
        Some(AngleFormat::Degrees)
    );
}
Output


From Existing Data

You can also initialize an OrbitTrajectory from existing epoch and state data:

import brahe as bh
import numpy as np

bh.initialize_eop()

# Create epochs
epoch0 = bh.Epoch.from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, bh.TimeSystem.UTC)
epoch1 = epoch0 + 60.0
epoch2 = epoch0 + 120.0

# Create Cartesian states in ECI
state0 = np.array([bh.R_EARTH + 500e3, 0.0, 0.0, 0.0, 7600.0, 0.0])
state1 = np.array([bh.R_EARTH + 500e3, 456000.0, 0.0, -7600.0, 0.0, 0.0])
state2 = np.array([bh.R_EARTH + 500e3, 0.0, 0.0, 0.0, -7600.0, 0.0])

# Create trajectory from data
epochs = [epoch0, epoch1, epoch2]
states = np.array([state0, state1, state2])  # Flattened array
traj = bh.OrbitTrajectory.from_orbital_data(
    epochs, states, bh.OrbitFrame.ECI, bh.OrbitRepresentation.CARTESIAN, None
)

print(f"Trajectory length: {len(traj)}")
use brahe as bh;
use bh::time::Epoch;
use bh::trajectories::SOrbitTrajectory;
use bh::trajectories::traits::{OrbitFrame, OrbitRepresentation, OrbitalTrajectory};
use bh::traits::Trajectory;
use bh::constants::R_EARTH;
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();

    // Create epochs
    let epoch0 = Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0,
        bh::time::TimeSystem::UTC);
    let epoch1 = epoch0 + 60.0;
    let epoch2 = epoch0 + 120.0;

    // Create Cartesian states
    let state0 = na::SVector::<f64, 6>::new(
        R_EARTH + 500e3, 0.0, 0.0, 0.0, 7600.0, 0.0
    );
    let state1 = na::SVector::<f64, 6>::new(
        R_EARTH + 500e3, 456000.0, 0.0, -7600.0, 0.0, 0.0
    );
    let state2 = na::SVector::<f64, 6>::new(
        R_EARTH + 500e3, 0.0, 0.0, 0.0, -7600.0, 0.0
    );

    // Create trajectory from data
    let epochs = vec![epoch0, epoch1, epoch2];
    let states = vec![state0, state1, state2];
    let traj = SOrbitTrajectory::from_orbital_data(
        epochs,
        states,
        OrbitFrame::ECI,
        OrbitRepresentation::Cartesian,
        None, // Angle Format
        None  // No covariances
    );

    println!("Trajectory length: {}", traj.len());
}
Output
Trajectory length: 3
Trajectory length: 3

From Propagator

The most common way to get an OrbitTrajectory from a propagator. All orbit propagators in Brahe have a *.trajectory attribute which is an OrbitTrajectory.

See the Propagators section for more details on propagators.

import brahe as bh
import numpy as np

bh.initialize_eop()

# Define orbital elements for a 500 km circular orbit
a = bh.R_EARTH + 500e3
e = 0.001
i = 97.8  # Sun-synchronous
raan = 15.0
argp = 30.0
M = 0.0
oe = np.array([a, e, i, raan, argp, M])

# Create epoch and propagator
epoch = bh.Epoch.from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, bh.TimeSystem.UTC)
propagator = bh.KeplerianPropagator.from_keplerian(
    epoch, oe, bh.AngleFormat.DEGREES, 60.0
)

# Propagate for several steps
propagator.propagate_steps(10)

# Access the trajectory
traj = propagator.trajectory
print(f"Trajectory length: {len(traj)}")  # Output: 11 (initial + 10 steps)
print(f"Frame: {traj.frame}")  # Output: OrbitFrame.ECI
print(f"Representation: {traj.representation}")  # Output: Keplerian
use brahe as bh;
use bh::time::Epoch;
use bh::traits::{Trajectory, SStatePropagator};
use bh::{KeplerianPropagator, R_EARTH, AngleFormat};
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();

    // Define orbital elements
    let oe = na::SVector::<f64, 6>::new(
        R_EARTH + 500e3, 0.001, 97.8_f64.to_radians(),
        15.0_f64.to_radians(), 30.0_f64.to_radians(), 0.0
    );

    // Create epoch and propagator
    let epoch = Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0,
        bh::time::TimeSystem::UTC);
    let mut propagator = KeplerianPropagator::from_keplerian(
        epoch, oe, AngleFormat::Radians, 60.0
    );

    // Propagate for several steps
    propagator.propagate_steps(10);

    // Access the trajectory
    let traj = &propagator.trajectory;
    println!("Trajectory length: {}", traj.len());  // Output: 11
    println!("Frame: {}", traj.frame);  // Output: ECI
    println!("Representation: {}", traj.representation);  // Output: Keplerian
}
Output
1
2
3
Trajectory length: 11
Frame: ECI
Representation: Keplerian
1
2
3
Trajectory length: 11
Frame: ECI
Representation: Keplerian

Frame Conversions

The key feature of OrbitTrajectory is automatic frame conversions of the trajectory data to different reference frames and representations. In particular, with a single method call you can convert between ECI and ECEF frames, and between Cartesian and Keplerian representations.

Converting ECI to ECEF

Convert a trajectory from Earth-Centered Inertial (ECI) to Earth-Centered Earth-Fixed (ECEF):

import brahe as bh
import numpy as np

bh.initialize_eop()

# Create trajectory in ECI frame
traj_eci = bh.OrbitTrajectory(
    6, bh.OrbitFrame.ECI, bh.OrbitRepresentation.CARTESIAN, None
)

# Add states in ECI
epoch0 = bh.Epoch.from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, bh.TimeSystem.UTC)
for i in range(5):
    epoch = epoch0 + i * 60.0
    # Define state at epoch
    state_eci = np.array([bh.R_EARTH + 500e3, i * 100e3, 0.0, 0.0, 7600.0, 0.0])
    traj_eci.add(epoch, state_eci)

print(f"Original frame: {traj_eci.frame}")
print(f"Original representation: {traj_eci.representation}")

# Convert all states in trajectory to ECEF
traj_ecef = traj_eci.to_ecef()

print(f"\nConverted frame: {traj_ecef.frame}")
print(f"Converted representation: {traj_ecef.representation}")
print(f"Same number of states: {len(traj_ecef)}")

# Compare first states
_, state_eci = traj_eci.first()
_, state_ecef = traj_ecef.first()
print(
    f"\nFirst ECI state: [{state_eci[0]:.2f}, {state_eci[1]:.2f}, {state_eci[2]:.2f}] m"
)
print(
    f"First ECEF state: [{state_ecef[0]:.2f}, {state_ecef[1]:.2f}, {state_ecef[2]:.2f}] m"
)
use brahe as bh;
use bh::time::Epoch;
use bh::trajectories::SOrbitTrajectory;
use bh::trajectories::traits::{OrbitFrame, OrbitRepresentation, OrbitalTrajectory};
use bh::traits::Trajectory;
use bh::constants::R_EARTH;
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();

    // Create trajectory in ECI frame
    let mut traj_eci = SOrbitTrajectory::new(
        OrbitFrame::ECI,
        OrbitRepresentation::Cartesian,
        None
    );

    // Add states in ECI
    let epoch0 = Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0,
        bh::time::TimeSystem::UTC);
    for i in 0..5 {
        let epoch = epoch0 + (i as f64) * 60.0;
        let state_eci = na::SVector::<f64, 6>::new(
            R_EARTH + 500e3, (i as f64) * 100e3, 0.0, 0.0, 7600.0, 0.0
        );
        traj_eci.add(epoch, state_eci);
    }

    println!("Original frame: {:?}", traj_eci.frame);
    println!("Original representation: {:?}", traj_eci.representation);

    // Convert all states in trajectory to ECEF
    let traj_ecef = traj_eci.to_ecef();

    println!("\nConverted frame: {:?}", traj_ecef.frame);
    println!("Converted representation: {:?}", traj_ecef.representation);
    println!("Same number of states: {}", traj_ecef.len());

    // Compare first states
    let state_eci_first = traj_eci.state_at_idx(0).unwrap();
    let state_ecef_first = traj_ecef.state_at_idx(0).unwrap();
    println!("\nFirst state ECI: [{}, {}, {}] m",
        state_eci_first[0], state_eci_first[1], state_eci_first[2]
    );
    println!("First state ECEF: [{}, {}, {}] m",
        state_ecef_first[0], state_ecef_first[1], state_ecef_first[2]
    );
}
Output
1
2
3
4
5
6
7
8
9
Original frame: ECI
Original representation: Cartesian

Converted frame: ECEF
Converted representation: Cartesian
Same number of states: 5

First ECI state: [6878136.30, 0.00, 0.00] m
First ECEF state: [-1176064.06, -6776826.51, 15961.82] m
1
2
3
4
5
6
7
8
9
Original frame: OrbitFrame(Earth-Centered Inertial)
Original representation: OrbitRepresentation(Cartesian)

Converted frame: OrbitFrame(Earth-Centered Earth-Fixed)
Converted representation: OrbitRepresentation(Cartesian)
Same number of states: 5

First state ECI: [6878136.3, 0, 0] m
First state ECEF: [-1176064.0596141217, -6776826.507241378, 15961.823588606136] m

Converting ECEF to ECI

Convert from ECEF back to ECI:

import brahe as bh
import numpy as np

bh.initialize_eop()

# Create trajectory in ECEF frame
traj_ecef = bh.OrbitTrajectory(
    6, bh.OrbitFrame.ECEF, bh.OrbitRepresentation.CARTESIAN, None
)

# Add dummy states in ECEF
epoch0 = bh.Epoch.from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, bh.TimeSystem.UTC)
for i in range(3):
    epoch = epoch0 + i * 60.0
    # Define state at epoch
    state_ecef = np.array([bh.R_EARTH + 500e3, 0.0, 0.0, 0.0, 0.0, 7600.0])
    traj_ecef.add(epoch, state_ecef)

print(f"Original frame: {traj_ecef.frame}")  # Output: OrbitFrame.ECEF

# Convert to ECI
traj_eci = traj_ecef.to_eci()

print(f"Converted frame: {traj_eci.frame}")  # Output: OrbitFrame.ECI
print(f"Trajectory length: {len(traj_eci)}")  # Output: 3

# Iterate over converted states
for epoch, state_eci in traj_eci:
    pos_mag = np.linalg.norm(state_eci[0:3])
    vel_mag = np.linalg.norm(state_eci[3:6])
    print(f"Epoch: {epoch}")
    print(f"  Position magnitude: {pos_mag / 1e3:.2f} km")
    print(f"  Velocity magnitude: {vel_mag:.2f} m/s")
use brahe as bh;
use bh::time::Epoch;
use bh::trajectories::SOrbitTrajectory;
use bh::trajectories::traits::{OrbitFrame, OrbitRepresentation, OrbitalTrajectory};
use bh::traits::Trajectory;
use bh::constants::R_EARTH;
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();

    // Create trajectory in ECEF frame
    let mut traj_ecef = SOrbitTrajectory::new(
        OrbitFrame::ECEF,
        OrbitRepresentation::Cartesian,
        None
    );

    // Add dummy states in ECEF
    let epoch0 = Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0,
        bh::time::TimeSystem::UTC);
    for i in 0..3 {
        let epoch = epoch0 + (i as f64) * 60.0;
        let state_ecef = na::SVector::<f64, 6>::new(
            R_EARTH + 500e3, 0.0, 0.0, 0.0, 0.0, 7600.0
        );
        traj_ecef.add(epoch, state_ecef);
    }

    println!("Original frame: {:?}", traj_ecef.frame);

    // Convert to ECI
    let traj_eci = traj_ecef.to_eci();

    println!("Converted frame: {:?}", traj_eci.frame);
    println!("Trajectory length: {}", traj_eci.len());

    // Iterate over converted states
    for (epoch, state_eci) in &traj_eci {
        let pos_mag = state_eci.fixed_rows::<3>(0).norm();
        let vel_mag = state_eci.fixed_rows::<3>(3).norm();
        println!("Epoch: {}", epoch);
        println!("  Position magnitude: {:.2} km", pos_mag / 1e3);
        println!("  Velocity magnitude: {:.2} m/s", vel_mag);
    }
}
Output
Original frame: ECEF
Converted frame: ECI
Trajectory length: 3
Epoch: 2024-01-01 00:00:00.000 UTC
  Position magnitude: 6878.14 km
  Velocity magnitude: 7616.53 m/s
Epoch: 2024-01-01 00:01:00.000 UTC
  Position magnitude: 6878.14 km
  Velocity magnitude: 7616.53 m/s
Epoch: 2024-01-01 00:02:00.000 UTC
  Position magnitude: 6878.14 km
  Velocity magnitude: 7616.53 m/s
Original frame: OrbitFrame(Earth-Centered Earth-Fixed)
Converted frame: OrbitFrame(Earth-Centered Inertial)
Trajectory length: 3
Epoch: 2024-01-01 00:00:00.000 UTC
  Position magnitude: 6878.14 km
  Velocity magnitude: 7616.53 m/s
Epoch: 2024-01-01 00:01:00.000 UTC
  Position magnitude: 6878.14 km
  Velocity magnitude: 7616.53 m/s
Epoch: 2024-01-01 00:02:00.000 UTC
  Position magnitude: 6878.14 km
  Velocity magnitude: 7616.53 m/s

Round-Trip Frame Conversion

Convert from ECI to ECEF and back to verify consistency:

import brahe as bh
import numpy as np

bh.initialize_eop()

# Create trajectory in ECI
traj_eci_original = bh.OrbitTrajectory(
    6, bh.OrbitFrame.ECI, bh.OrbitRepresentation.CARTESIAN, None
)

# Add a state
epoch = bh.Epoch.from_datetime(2024, 1, 1, 12, 0, 0.0, 0.0, bh.TimeSystem.UTC)
state_original = np.array([bh.R_EARTH + 500e3, 0.0, 0.0, 0.0, 7600.0, 0.0])
traj_eci_original.add(epoch, state_original)

# Convert to ECEF and back to ECI
traj_ecef = traj_eci_original.to_ecef()
traj_eci_roundtrip = traj_ecef.to_eci()

# Compare original and round-trip states
_, state_roundtrip = traj_eci_roundtrip.first()
diff = np.abs(state_original - state_roundtrip)

print(f"Position difference: {np.linalg.norm(diff[0:3]):.6e} m")
print(f"Velocity difference: {np.linalg.norm(diff[3:6]):.6e} m/s")
use brahe as bh;
use bh::time::Epoch;
use bh::trajectories::SOrbitTrajectory;
use bh::trajectories::traits::{OrbitFrame, OrbitRepresentation, OrbitalTrajectory};
use bh::traits::Trajectory;
use bh::constants::R_EARTH;
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();

    // Create trajectory in ECI
    let mut traj_eci_original = SOrbitTrajectory::new(
        OrbitFrame::ECI,
        OrbitRepresentation::Cartesian,
        None
    );

    // Add a state
    let epoch = Epoch::from_datetime(2024, 1, 1, 12, 0, 0.0, 0.0,
        bh::time::TimeSystem::UTC);
    let state_original = na::SVector::<f64, 6>::new(
        R_EARTH + 500e3, 0.0, 0.0, 0.0, 7600.0, 0.0
    );
    traj_eci_original.add(epoch, state_original);

    // Convert to ECEF and back to ECI
    let traj_ecef = traj_eci_original.to_ecef();
    let traj_eci_roundtrip = traj_ecef.to_eci();

    // Compare original and round-trip states
    let (_, state_roundtrip) = traj_eci_roundtrip.first().unwrap();
    let diff = state_original - state_roundtrip;

    println!("Position difference: {:.6e} m",
        diff.fixed_rows::<3>(0).norm());
    println!("Velocity difference: {:.6e} m/s",
        diff.fixed_rows::<3>(3).norm());
}
Output
Position difference: 1.760544e-11 m
Velocity difference: 1.867232e-12 m/s
Position difference: 1.760544e-11 m
Velocity difference: 1.867232e-12 m/s

Converting Cartesian to Keplerian

Convert from Cartesian position/velocity to Keplerian orbital elements:

import brahe as bh
import numpy as np

bh.initialize_eop()

# Create trajectory in ECI Cartesian
traj_cart = bh.OrbitTrajectory(
    6, bh.OrbitFrame.ECI, bh.OrbitRepresentation.CARTESIAN, None
)

# Add Cartesian states
epoch0 = bh.Epoch.from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, bh.TimeSystem.UTC)
for i in range(3):
    epoch = epoch0 + i * 300.0  # 5-minute intervals
    # Use orbital elements to create realistic Cartesian states
    oe = np.array([bh.R_EARTH + 500e3, 0.001, 97.8, 15.0, 30.0, i * 10.0])
    state_cart = bh.state_koe_to_eci(oe, bh.AngleFormat.DEGREES)
    traj_cart.add(epoch, state_cart)

print(f"Original representation: {traj_cart.representation}")

# Convert to Keplerian with degrees
traj_kep = traj_cart.to_keplerian(bh.AngleFormat.DEGREES)

print(f"Converted representation: {traj_kep.representation}")
print(f"Angle format: {traj_kep.angle_format}")

# Examine Keplerian elements
for epoch, oe in traj_kep:
    print(f"\nEpoch: {epoch}")
    print(f"  Semi-major axis: {oe[0] / 1e3:.2f} km")
    print(f"  Eccentricity: {oe[1]:.6f}")
    print(f"  Inclination: {oe[2]:.2f}°")
    print(f"  RAAN: {oe[3]:.2f}°")
    print(f"  Argument of perigee: {oe[4]:.2f}°")
    print(f"  Mean anomaly: {oe[5]:.2f}°")
use brahe as bh;
use bh::time::Epoch;
use bh::trajectories::SOrbitTrajectory;
use bh::trajectories::traits::{OrbitFrame, OrbitRepresentation, OrbitalTrajectory};
use bh::traits::Trajectory;
use bh::{state_koe_to_eci, R_EARTH, AngleFormat};
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();

    // Create trajectory in ECI Cartesian
    let mut traj_cart = SOrbitTrajectory::new(
        OrbitFrame::ECI,
        OrbitRepresentation::Cartesian,
        None
    );

    // Add Cartesian states
    let epoch0 = Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0,
        bh::time::TimeSystem::UTC);
    for i in 0..3 {
        let epoch = epoch0 + (i as f64) * 300.0;
        let oe = na::SVector::<f64, 6>::new(
            R_EARTH + 500e3, 0.001, 97.8, 15.0, 30.0, (i as f64) * 10.0
        );
        let state_cart = state_koe_to_eci(oe, AngleFormat::Degrees);
        traj_cart.add(epoch, state_cart);
    }

    println!("Original representation: {:?}", traj_cart.representation);

    // Convert to Keplerian with radians
    let traj_kep = traj_cart.to_keplerian(AngleFormat::Degrees);

    println!("Converted representation: {:?}", traj_kep.representation);
    println!("Angle format: {:?}", traj_kep.angle_format);

    // Examine Keplerian elements
    for (epoch, oe) in &traj_kep {
        println!("\nEpoch: {}", epoch);
        println!("  Semi-major axis: {:.2} km", oe[0] / 1e3);
        println!("  Eccentricity: {:.6}", oe[1]);
        println!("  Inclination: {:.2}°", oe[2]);
        println!("  RAAN: {:.2}°", oe[3]);
        println!("  Argument of perigee: {:.2}°", oe[4]);
        println!("  Mean Anomaly: {:.2}°", oe[5]);
    }
}
Output
Original representation: Cartesian
Converted representation: Keplerian
Angle format: Degrees

Epoch: 2024-01-01 00:00:00.000 UTC
  Semi-major axis: 6878.14 km
  Eccentricity: 0.001000
  Inclination: 97.80°
  RAAN: 15.00°
  Argument of perigee: 30.00°
  Mean anomaly: 0.00°

Epoch: 2024-01-01 00:05:00.000 UTC
  Semi-major axis: 6878.14 km
  Eccentricity: 0.001000
  Inclination: 97.80°
  RAAN: 15.00°
  Argument of perigee: 30.00°
  Mean anomaly: 10.00°

Epoch: 2024-01-01 00:10:00.000 UTC
  Semi-major axis: 6878.14 km
  Eccentricity: 0.001000
  Inclination: 97.80°
  RAAN: 15.00°
  Argument of perigee: 30.00°
  Mean anomaly: 20.00°
Original representation: OrbitRepresentation(Cartesian)
Converted representation: OrbitRepresentation(Keplerian)
Angle format: Some(Degrees)

Epoch: 2024-01-01 00:00:00.000 UTC
  Semi-major axis: 6878.14 km
  Eccentricity: 0.001000
  Inclination: 97.80°
  RAAN: 15.00°
  Argument of perigee: 30.00°
  Mean Anomaly: 0.00°

Epoch: 2024-01-01 00:05:00.000 UTC
  Semi-major axis: 6878.14 km
  Eccentricity: 0.001000
  Inclination: 97.80°
  RAAN: 15.00°
  Argument of perigee: 30.00°
  Mean Anomaly: 10.00°

Epoch: 2024-01-01 00:10:00.000 UTC
  Semi-major axis: 6878.14 km
  Eccentricity: 0.001000
  Inclination: 97.80°
  RAAN: 15.00°
  Argument of perigee: 30.00°
  Mean Anomaly: 20.00°

Converting with Different Angle Formats

Convert to Keplerian with different angle formats:

import brahe as bh
import numpy as np

bh.initialize_eop()

# Create trajectory in ECI Cartesian
traj_cart = bh.OrbitTrajectory(
    6, bh.OrbitFrame.ECI, bh.OrbitRepresentation.CARTESIAN, None
)

# Add a state
epoch = bh.Epoch.from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, bh.TimeSystem.UTC)
oe = np.array([bh.R_EARTH + 500e3, 0.001, 0.9, 1.0, 0.5, 0.0])
state_cart = bh.state_koe_to_eci(oe, bh.AngleFormat.RADIANS)
traj_cart.add(epoch, state_cart)

# Convert to Keplerian with radians
traj_kep_rad = traj_cart.to_keplerian(bh.AngleFormat.RADIANS)
_, oe_rad = traj_kep_rad.first()

# Convert to Keplerian with degrees
traj_kep_deg = traj_cart.to_keplerian(bh.AngleFormat.DEGREES)
_, oe_deg = traj_kep_deg.first()

print("Radians version:")
print(f"  Inclination: {oe_rad[2]:.6f} rad = {np.degrees(oe_rad[2]):.2f}°")

print("\nDegrees version:")
print(f"  Inclination: {oe_deg[2]:.2f}°")
use brahe as bh;
use bh::time::Epoch;
use bh::trajectories::SOrbitTrajectory;
use bh::trajectories::traits::{OrbitFrame, OrbitRepresentation, OrbitalTrajectory};
use bh::traits::Trajectory;
use bh::{state_koe_to_eci, R_EARTH, AngleFormat};
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();

    // Create trajectory in ECI Cartesian
    let mut traj_cart = SOrbitTrajectory::new(
        OrbitFrame::ECI,
        OrbitRepresentation::Cartesian,
        None
    );

    // Add a state
    let epoch = Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0,
        bh::time::TimeSystem::UTC);
    let oe = na::SVector::<f64, 6>::new(
        R_EARTH + 500e3, 0.001, 0.9, 1.0, 0.5, 0.0
    );
    let state_cart = state_koe_to_eci(oe, AngleFormat::Radians);
    traj_cart.add(epoch, state_cart);

    // Convert to Keplerian with radians
    let traj_kep_rad = traj_cart.to_keplerian(AngleFormat::Radians);
    let (_, oe_rad) = traj_kep_rad.first().unwrap();

    // Convert to Keplerian with degrees
    let traj_kep_deg = traj_cart.to_keplerian(AngleFormat::Degrees);
    let (_, oe_deg) = traj_kep_deg.first().unwrap();

    println!("Radians version:");
    println!("  Inclination: {:.6} rad = {:.2}°", oe_rad[2], oe_rad[2].to_degrees());

    println!("\nDegrees version:");
    println!("  Inclination: {:.2}°", oe_deg[2]);
}
Output
1
2
3
4
5
Radians version:
  Inclination: 0.900000 rad = 51.57°

Degrees version:
  Inclination: 51.57°
1
2
3
4
5
Radians version:
  Inclination: 0.900000 rad = 51.57°

Degrees version:
  Inclination: 51.57°

Combined Frame and Representation Conversions

Every conversion method returns a new OrbitTrajectory instance, so you can chain conversions together if desired:

import brahe as bh
import numpy as np

bh.initialize_eop()

# Start with ECI Cartesian trajectory
traj_eci_cart = bh.OrbitTrajectory(
    6, bh.OrbitFrame.ECI, bh.OrbitRepresentation.CARTESIAN, None
)

# Add states
epoch = bh.Epoch.from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, bh.TimeSystem.UTC)
oe = np.array([bh.R_EARTH + 500e3, 0.001, 0.9, 1.0, 0.5, 0.0])
state_cart = bh.state_koe_to_eci(oe, bh.AngleFormat.RADIANS)
traj_eci_cart.add(epoch, state_cart)

print("Original:")
print(f"  Frame: {traj_eci_cart.frame}")
print(f"  Representation: {traj_eci_cart.representation}")

# Convert to ECEF frame (stays Cartesian)
traj_ecef_cart = traj_eci_cart.to_ecef()
print("\nAfter to_ecef():")
print(f"  Frame: {traj_ecef_cart.frame}")
print(f"  Representation: {traj_ecef_cart.representation}")

# Convert back to ECI
traj_eci_cart2 = traj_ecef_cart.to_eci()
print("\nAfter to_eci():")
print(f"  Frame: {traj_eci_cart2.frame}")
print(f"  Representation: {traj_eci_cart2.representation}")

# Convert to Keplerian (in ECI frame)
traj_eci_kep = traj_eci_cart2.to_keplerian(bh.AngleFormat.DEGREES)
print("\nAfter to_keplerian():")
print(f"  Frame: {traj_eci_kep.frame}")
print(f"  Representation: {traj_eci_kep.representation}")
print(f"  Angle format: {traj_eci_kep.angle_format}")
use brahe as bh;
use bh::time::Epoch;
use bh::trajectories::SOrbitTrajectory;
use bh::trajectories::traits::{OrbitFrame, OrbitRepresentation, OrbitalTrajectory};
use bh::traits::Trajectory;
use bh::{state_koe_to_eci, R_EARTH, AngleFormat};
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();

    // Start with ECI Cartesian trajectory
    let mut traj_eci_cart = SOrbitTrajectory::new(
        OrbitFrame::ECI,
        OrbitRepresentation::Cartesian,
        None
    );

    // Add states
    let epoch = Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0,
        bh::time::TimeSystem::UTC);
    let oe = na::SVector::<f64, 6>::new(
        R_EARTH + 500e3, 0.001, 0.9, 1.0, 0.5, 0.0
    );
    let state_cart = state_koe_to_eci(oe, AngleFormat::Radians);
    traj_eci_cart.add(epoch, state_cart);

    println!("Original:");
    println!("  Frame: {:?}", traj_eci_cart.frame);
    println!("  Representation: {:?}", traj_eci_cart.representation);

    // Convert to ECEF frame (stays Cartesian)
    let traj_ecef_cart = traj_eci_cart.to_ecef();
    println!("\nAfter to_ecef():");
    println!("  Frame: {:?}", traj_ecef_cart.frame);
    println!("  Representation: {:?}", traj_ecef_cart.representation);

    // Convert back to ECI
    let traj_eci_cart2 = traj_ecef_cart.to_eci();
    println!("\nAfter to_eci():");
    println!("  Frame: {:?}", traj_eci_cart2.frame);
    println!("  Representation: {:?}", traj_eci_cart2.representation);

    // Convert to Keplerian (in ECI frame)
    let traj_eci_kep = traj_eci_cart2.to_keplerian(AngleFormat::Radians);
    println!("\nAfter to_keplerian():");
    println!("  Frame: {:?}", traj_eci_kep.frame);
    println!("  Representation: {:?}", traj_eci_kep.representation);
}
Output
Original:
  Frame: ECI
  Representation: Cartesian

After to_ecef():
  Frame: ECEF
  Representation: Cartesian

After to_eci():
  Frame: ECI
  Representation: Cartesian

After to_keplerian():
  Frame: ECI
  Representation: Keplerian
  Angle format: Degrees
Original:
  Frame: OrbitFrame(Earth-Centered Inertial)
  Representation: OrbitRepresentation(Cartesian)

After to_ecef():
  Frame: OrbitFrame(Earth-Centered Earth-Fixed)
  Representation: OrbitRepresentation(Cartesian)

After to_eci():
  Frame: OrbitFrame(Earth-Centered Inertial)
  Representation: OrbitRepresentation(Cartesian)

After to_keplerian():
  Frame: OrbitFrame(Earth-Centered Inertial)
  Representation: OrbitRepresentation(Keplerian)

Standard Trajectory Operations

OrbitTrajectory supports all standard trajectory operations since it implements the Trajectory and Interpolatable traits:

import brahe as bh
import numpy as np

bh.initialize_eop()

# Create trajectory
traj = bh.OrbitTrajectory(6, bh.OrbitFrame.ECI, bh.OrbitRepresentation.CARTESIAN, None)

# Add states
epoch0 = bh.Epoch.from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, bh.TimeSystem.UTC)
for i in range(10):
    epoch = epoch0 + i * 60.0
    oe = np.array([bh.R_EARTH + 500e3, 0.001, 0.9, 1.0, 0.5, i * 0.1])
    state = bh.state_koe_to_eci(oe, bh.AngleFormat.RADIANS)
    traj.add(epoch, state)

# Query properties
print(f"Length: {len(traj)}")
print(f"Timespan: {traj.timespan():.1f} seconds")
print(f"Start epoch: {traj.start_epoch()}")
print(f"End epoch: {traj.end_epoch()}")

# Interpolate at intermediate time
interp_epoch = epoch0 + 45.0
interp_state = traj.interpolate(interp_epoch)
print(f"\nInterpolated state at {interp_epoch}:")
print(f"  Position (km): {interp_state[0:3] / 1e3}")
print(f"  Velocity (m/s): {interp_state[3:6]}")

# Iterate over states
for i, (epoch, state) in enumerate(traj):
    if i < 2:  # Just show first two
        print(
            f"State {i}: Epoch={epoch}, Position magnitude={np.linalg.norm(state[0:3]) / 1e3:.2f} km"
        )
use brahe as bh;
use bh::time::Epoch;
use bh::trajectories::SOrbitTrajectory;
use bh::trajectories::traits::{OrbitFrame, OrbitRepresentation};
use bh::traits::{Trajectory, InterpolatableTrajectory};
use bh::{state_koe_to_eci, R_EARTH, AngleFormat};
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();

    // Create trajectory
    let mut traj = SOrbitTrajectory::new(
        OrbitFrame::ECI,
        OrbitRepresentation::Cartesian,
        None
    );

    // Add states
    let epoch0 = Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0,
        bh::time::TimeSystem::UTC);
    for i in 0..10 {
        let epoch = epoch0 + (i as f64) * 60.0;
        let oe = na::SVector::<f64, 6>::new(
            R_EARTH + 500e3, 0.001, 0.9, 1.0, 0.5, (i as f64) * 0.1
        );
        let state = state_koe_to_eci(oe, AngleFormat::Radians);
        traj.add(epoch, state);
    }

    // Query properties
    println!("Length: {}", traj.len());
    println!("Timespan: {:.1} seconds", traj.timespan().unwrap());
    println!("Start epoch: {}", traj.start_epoch().unwrap());
    println!("End epoch: {}", traj.end_epoch().unwrap());

    // Interpolate at intermediate time
    let interp_epoch = epoch0 + 45.0;
    let interp_state = traj.interpolate(&interp_epoch).unwrap();
    println!("\nInterpolated state at {}:", interp_epoch);
    println!("  Position (km): [{}, {}, {}] km",
        interp_state[0] / 1e3, interp_state[1] / 1e3, interp_state[2] / 1e3
    );
    println!("  Velocity (m/s): [{}, {}, {}] m/s",
        interp_state[3], interp_state[4], interp_state[5]
    );

    // Iterate over states
    for (i, (epoch, state)) in traj.into_iter().enumerate().take(2) {
        let pos_mag = state.fixed_rows::<3>(0).norm();
        println!("State {}: Epoch={}, Position magnitude={:.2} km",
            i, epoch, pos_mag / 1e3);
    }
}
Output
Length: 10
Timespan: 540.0 seconds
Start epoch: 2024-01-01 00:00:00.000 UTC
End epoch: 2024-01-01 00:09:00.000 UTC

Interpolated state at 2024-01-01 00:00:45.000 UTC:
  Position (km): [1159.01597302 6101.29789026 2925.16369358]
  Velocity (m/s): [-5578.86734152 -1338.77483001  5004.22925364]
State 0: Epoch=2024-01-01 00:00:00.000 UTC, Position magnitude=6871.26 km
State 1: Epoch=2024-01-01 00:01:00.000 UTC, Position magnitude=6871.29 km
Length: 10
Timespan: 540.0 seconds
Start epoch: 2024-01-01 00:00:00.000 UTC
End epoch: 2024-01-01 00:09:00.000 UTC

Interpolated state at 2024-01-01 00:00:45.000 UTC:
  Position (km): [1159.0159730226278, 6101.297890257402, 2925.16369357997] km
  Velocity (m/s): [-5578.867341523014, -1338.7748300095711, 5004.22925363932] m/s
State 0: Epoch=2024-01-01 00:00:00.000 UTC, Position magnitude=6871.26 km
State 1: Epoch=2024-01-01 00:01:00.000 UTC, Position magnitude=6871.29 km

Practical Workflow Example

A complete example showing propagation, frame conversion, and analysis:

import brahe as bh
import numpy as np

bh.initialize_eop()

# 1. Define orbit and create propagator
a = bh.R_EARTH + 500e3  # 500 km altitude
e = 0.001  # Nearly circular
i = 97.8  # Sun-synchronous
raan = 15.0
argp = 30.0
M = 0.0
oe = np.array([a, e, i, raan, argp, M])

epoch = bh.Epoch.from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0, bh.TimeSystem.UTC)
propagator = bh.KeplerianPropagator.from_keplerian(
    epoch, oe, bh.AngleFormat.DEGREES, 60.0
)

# 2. Propagate for one orbit period
period = bh.orbital_period(a)
end_epoch = epoch + period
propagator.propagate_to(end_epoch)

# 3. Get trajectory in ECI Cartesian
traj_eci = propagator.trajectory
print(f"Propagated {len(traj_eci)} states over {traj_eci.timespan() / 60:.1f} minutes")

# 4. Convert to ECEF to analyze ground track
traj_ecef = traj_eci.to_ecef()
print("\nGround track in ECEF frame:")
for i, (epoch, state_ecef) in enumerate(traj_ecef):
    if i % 10 == 0:  # Sample every 10 states
        # Convert ECEF to geodetic for latitude/longitude
        lat, lon, alt = bh.position_ecef_to_geodetic(
            state_ecef[0:3], bh.AngleFormat.DEGREES
        )
        print(f"  {epoch}: Lat={lat:6.2f}°, Lon={lon:7.2f}°, Alt={alt / 1e3:6.2f} km")

# 5. Convert to Keplerian to analyze orbital evolution
traj_kep = traj_eci.to_keplerian(bh.AngleFormat.DEGREES)
first_oe = traj_kep.state_at_idx(0)
last_oe = traj_kep.state_at_idx(len(traj_kep) - 1)

print("\nOrbital element evolution:")
print(f"  Semi-major axis: {first_oe[0] / 1e3:.2f} km → {last_oe[0] / 1e3:.2f} km")
print(f"  Eccentricity: {first_oe[1]:.6f}{last_oe[1]:.6f}")
print(f"  Inclination: {first_oe[2]:.2f}° → {last_oe[2]:.2f}°")
print(f"  True anomaly: {first_oe[5]:.2f}° → {last_oe[5]:.2f}°")
use brahe as bh;
use bh::time::Epoch;
use bh::traits::{Trajectory, SStatePropagator};
use bh::{KeplerianPropagator, orbital_period, position_ecef_to_geodetic, R_EARTH, AngleFormat};
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();

    // 1. Define orbit and create propagator
    let oe = na::SVector::<f64, 6>::new(
        R_EARTH + 500e3, 0.001, 97.8_f64.to_radians(),
        15.0_f64.to_radians(), 30.0_f64.to_radians(), 0.0
    );

    let epoch = Epoch::from_datetime(2024, 1, 1, 0, 0, 0.0, 0.0,
        bh::time::TimeSystem::UTC);
    let mut propagator = KeplerianPropagator::from_keplerian(
        epoch, oe, AngleFormat::Radians, 60.0
    );

    // 2. Propagate for one orbit period
    let period = orbital_period(R_EARTH + 500e3);
    let end_epoch = epoch + period;
    propagator.propagate_to(end_epoch);

    // 3. Get trajectory in ECI Cartesian
    let traj_eci = &propagator.trajectory;
    println!("Propagated {} states over {:.1} minutes",
        traj_eci.len(), traj_eci.timespan().unwrap() / 60.0);

    // 4. Convert to ECEF
    let traj_ecef = traj_eci.to_ecef();
    println!("\nGround track in ECEF frame:");
    for (i, (epoch, state_ecef)) in traj_ecef.into_iter().enumerate() {
        if i % 10 == 0 {
            let pos_ecef: na::Vector3<f64> = state_ecef.fixed_rows::<3>(0).into();
            let lla = position_ecef_to_geodetic(pos_ecef, AngleFormat::Degrees);
            println!("  {}: Lat={:6.2}°, Lon={:7.2}°, Alt={:6.2} km",
                epoch, lla[0], lla[1], lla[2] / 1e3);
        }
    }

    // 5. Convert to Keplerian
    let traj_kep = traj_eci.to_keplerian(AngleFormat::Radians);
    let first_oe = traj_kep.state_at_idx(0).unwrap();
    let last_oe = traj_kep.state_at_idx(traj_kep.len() - 1).unwrap();

    println!("\nOrbital element evolution:");
    println!("  Semi-major axis: {:.2} km → {:.2} km",
        first_oe[0] / 1e3, last_oe[0] / 1e3);
    println!("  Eccentricity: {:.6} → {:.6}",
        first_oe[1], last_oe[1]);
    println!("  Inclination: {:.2}° → {:.2}°",
        first_oe[2].to_degrees(), last_oe[2].to_degrees());
}
Output
Propagated 96 states over 94.6 minutes

Ground track in ECEF frame:
  2024-01-01 00:00:00.000 UTC: Lat=-89.31°, Lon=  29.98°, Alt=498.43 km
  2024-01-01 00:10:00.000 UTC: Lat=-106.05°, Lon=  67.10°, Alt=512.71 km
  2024-01-01 00:20:00.000 UTC: Lat=114.91°, Lon=  72.06°, Alt=517.69 km
  2024-01-01 00:30:00.000 UTC: Lat= 93.18°, Lon=  35.41°, Alt=509.96 km
  2024-01-01 00:40:00.000 UTC: Lat= 84.82°, Lon=  -2.37°, Alt=506.12 km
  2024-01-01 00:50:00.000 UTC: Lat= 76.09°, Lon= -40.08°, Alt=515.59 km
  2024-01-01 01:00:00.000 UTC: Lat= 46.92°, Lon= -76.10°, Alt=524.72 km
  2024-01-01 01:10:00.000 UTC: Lat=-87.12°, Lon= -62.75°, Alt=517.32 km
  2024-01-01 01:20:00.000 UTC: Lat=-101.18°, Lon= -25.46°, Alt=500.05 km
  2024-01-01 01:30:00.000 UTC: Lat=-109.11°, Lon=  12.49°, Alt=494.44 km

Orbital element evolution:
  Semi-major axis: 6878.14 km → 6878.14 km
  Eccentricity: 0.001000 → 0.001000
  Inclination: 97.80° → 97.80°
  True anomaly: 0.00° → 360.00°
Propagated 96 states over 94.6 minutes

Ground track in ECEF frame:
  2024-01-01 00:00:00.000 UTC: Lat=-89.31°, Lon=  29.98°, Alt=498.43 km
  2024-01-01 00:10:00.000 UTC: Lat=-106.05°, Lon=  67.10°, Alt=512.71 km
  2024-01-01 00:20:00.000 UTC: Lat=114.91°, Lon=  72.06°, Alt=517.69 km
  2024-01-01 00:30:00.000 UTC: Lat= 93.18°, Lon=  35.41°, Alt=509.96 km
  2024-01-01 00:40:00.000 UTC: Lat= 84.82°, Lon=  -2.37°, Alt=506.12 km
  2024-01-01 00:50:00.000 UTC: Lat= 76.09°, Lon= -40.08°, Alt=515.59 km
  2024-01-01 01:00:00.000 UTC: Lat= 46.92°, Lon= -76.10°, Alt=524.72 km
  2024-01-01 01:10:00.000 UTC: Lat=-87.12°, Lon= -62.75°, Alt=517.32 km
  2024-01-01 01:20:00.000 UTC: Lat=-101.18°, Lon= -25.46°, Alt=500.05 km
  2024-01-01 01:30:00.000 UTC: Lat=-109.11°, Lon=  12.49°, Alt=494.44 km

Orbital element evolution:
  Semi-major axis: 6878.14 km → 6878.14 km
  Eccentricity: 0.001000 → 0.001000
  Inclination: 97.80° → 97.80°

See Also