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
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.
importbraheasbhbh.initialize_eop()# Create trajectory in ECI frame, Keplerian representation with radianstraj_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 degreestraj_kep_deg=bh.OrbitTrajectory(6,bh.OrbitFrame.ECI,bh.OrbitRepresentation.KEPLERIAN,bh.AngleFormat.DEGREES)
usebraheasbh;usebh::trajectories::SOrbitTrajectory;usebh::trajectories::traits::{OrbitFrame,OrbitRepresentation};usebh::AngleFormat;fnmain(){bh::initialize_eop().unwrap();// Create trajectory in ECI frame, Keplerian representation with radianslet_traj_kep_rad=SOrbitTrajectory::new(OrbitFrame::ECI,OrbitRepresentation::Keplerian,Some(AngleFormat::Radians));// Create trajectory in ECI frame, Keplerian representation with degreeslet_traj_kep_deg=SOrbitTrajectory::new(OrbitFrame::ECI,OrbitRepresentation::Keplerian,Some(AngleFormat::Degrees));}
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.
importbraheasbhimportnumpyasnpbh.initialize_eop()# Define orbital elements for a 500 km circular orbita=bh.R_EARTH+500e3e=0.001i=97.8# Sun-synchronousraan=15.0argp=30.0M=0.0oe=np.array([a,e,i,raan,argp,M])# Create epoch and propagatorepoch=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 stepspropagator.propagate_steps(10)# Access the trajectorytraj=propagator.trajectoryprint(f"Trajectory length: {len(traj)}")# Output: 11 (initial + 10 steps)print(f"Frame: {traj.frame}")# Output: OrbitFrame.ECIprint(f"Representation: {traj.representation}")# Output: Keplerian
usebraheasbh;usebh::time::Epoch;usebh::traits::{Trajectory,SStatePropagator};usebh::{KeplerianPropagator,R_EARTH,AngleFormat};usenalgebraasna;fnmain(){bh::initialize_eop().unwrap();// Define orbital elementsletoe=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 propagatorletepoch=Epoch::from_datetime(2024,1,1,0,0,0.0,0.0,bh::time::TimeSystem::UTC);letmutpropagator=KeplerianPropagator::from_keplerian(epoch,oe,AngleFormat::Radians,60.0);// Propagate for several stepspropagator.propagate_steps(10);// Access the trajectorylettraj=&propagator.trajectory;println!("Trajectory length: {}",traj.len());// Output: 11println!("Frame: {}",traj.frame);// Output: ECIprintln!("Representation: {}",traj.representation);// Output: Keplerian}
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.
importbraheasbhimportnumpyasnpbh.initialize_eop()# Create trajectory in ECI frametraj_eci=bh.OrbitTrajectory(6,bh.OrbitFrame.ECI,bh.OrbitRepresentation.CARTESIAN,None)# Add states in ECIepoch0=bh.Epoch.from_datetime(2024,1,1,0,0,0.0,0.0,bh.TimeSystem.UTC)foriinrange(5):epoch=epoch0+i*60.0# Define state at epochstate_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 ECEFtraj_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")# Output:# 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
usebraheasbh;usebh::time::Epoch;usebh::trajectories::SOrbitTrajectory;usebh::trajectories::traits::{OrbitFrame,OrbitRepresentation,OrbitalTrajectory};usebh::traits::Trajectory;usebh::constants::R_EARTH;usenalgebraasna;fnmain(){bh::initialize_eop().unwrap();// Create trajectory in ECI frameletmuttraj_eci=SOrbitTrajectory::new(OrbitFrame::ECI,OrbitRepresentation::Cartesian,None);// Add states in ECIletepoch0=Epoch::from_datetime(2024,1,1,0,0,0.0,0.0,bh::time::TimeSystem::UTC);foriin0..5{letepoch=epoch0+(iasf64)*60.0;letstate_eci=na::SVector::<f64,6>::new(R_EARTH+500e3,(iasf64)*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 ECEFlettraj_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 statesletstate_eci_first=traj_eci.state_at_idx(0).unwrap();letstate_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:// 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.507241379, 15961.82358860613] m
importbraheasbhimportnumpyasnpbh.initialize_eop()# Create trajectory in ECItraj_eci_original=bh.OrbitTrajectory(6,bh.OrbitFrame.ECI,bh.OrbitRepresentation.CARTESIAN,None)# Add a stateepoch=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 ECItraj_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")# Expected: Very small differences (numerical precision)# Output:# Position difference: 2.499882e-10 m# Velocity difference: 1.829382e-12 m/s
usebraheasbh;usebh::time::Epoch;usebh::trajectories::SOrbitTrajectory;usebh::trajectories::traits::{OrbitFrame,OrbitRepresentation,OrbitalTrajectory};usebh::traits::Trajectory;usebh::constants::R_EARTH;usenalgebraasna;fnmain(){bh::initialize_eop().unwrap();// Create trajectory in ECIletmuttraj_eci_original=SOrbitTrajectory::new(OrbitFrame::ECI,OrbitRepresentation::Cartesian,None);// Add a stateletepoch=Epoch::from_datetime(2024,1,1,12,0,0.0,0.0,bh::time::TimeSystem::UTC);letstate_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 ECIlettraj_ecef=traj_eci_original.to_ecef();lettraj_eci_roundtrip=traj_ecef.to_eci();// Compare original and round-trip stateslet(_,state_roundtrip)=traj_eci_roundtrip.first().unwrap();letdiff=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());// Expected: Very small differences (numerical precision)}// Output:// Position difference: 2.499882e-10 m// Velocity difference: 1.829382e-12 m/s
importbraheasbhimportnumpyasnpbh.initialize_eop()# Create trajectory in ECI Cartesiantraj_cart=bh.OrbitTrajectory(6,bh.OrbitFrame.ECI,bh.OrbitRepresentation.CARTESIAN,None)# Add a stateepoch=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 radianstraj_kep_rad=traj_cart.to_keplerian(bh.AngleFormat.RADIANS)_,oe_rad=traj_kep_rad.first()# Convert to Keplerian with degreestraj_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}°")# Output:# Radians version:# Inclination: 0.900000 rad = 51.57°# Degrees version:# Inclination: 51.57°
usebraheasbh;usebh::time::Epoch;usebh::trajectories::SOrbitTrajectory;usebh::trajectories::traits::{OrbitFrame,OrbitRepresentation,OrbitalTrajectory};usebh::traits::Trajectory;usebh::{state_koe_to_eci,R_EARTH,AngleFormat};usenalgebraasna;fnmain(){bh::initialize_eop().unwrap();// Create trajectory in ECI Cartesianletmuttraj_cart=SOrbitTrajectory::new(OrbitFrame::ECI,OrbitRepresentation::Cartesian,None);// Add a stateletepoch=Epoch::from_datetime(2024,1,1,0,0,0.0,0.0,bh::time::TimeSystem::UTC);letoe=na::SVector::<f64,6>::new(R_EARTH+500e3,0.001,0.9,1.0,0.5,0.0);letstate_cart=state_koe_to_eci(oe,AngleFormat::Radians);traj_cart.add(epoch,state_cart);// Convert to Keplerian with radianslettraj_kep_rad=traj_cart.to_keplerian(AngleFormat::Radians);let(_,oe_rad)=traj_kep_rad.first().unwrap();// Convert to Keplerian with degreeslettraj_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:// Radians version:// Inclination: 0.900000 rad = 51.57°// Degrees version:// Inclination: 51.57°
importbraheasbhimportnumpyasnpbh.initialize_eop()# Create trajectorytraj=bh.OrbitTrajectory(6,bh.OrbitFrame.ECI,bh.OrbitRepresentation.CARTESIAN,None)# Add statesepoch0=bh.Epoch.from_datetime(2024,1,1,0,0,0.0,0.0,bh.TimeSystem.UTC)foriinrange(10):epoch=epoch0+i*60.0oe=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 propertiesprint(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 timeinterp_epoch=epoch0+45.0interp_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 statesfori,(epoch,state)inenumerate(traj):ifi<2:# Just show first twoprint(f"State {i}: Epoch={epoch}, Position magnitude={np.linalg.norm(state[0:3])/1e3:.2f} km")# 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
importbraheasbhimportnumpyasnpbh.initialize_eop()# 1. Define orbit and create propagatora=bh.R_EARTH+500e3# 500 km altitudee=0.001# Nearly circulari=97.8# Sun-synchronousraan=15.0argp=30.0M=0.0oe=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 periodperiod=bh.orbital_period(a)end_epoch=epoch+periodpropagator.propagate_to(end_epoch)# 3. Get trajectory in ECI Cartesiantraj_eci=propagator.trajectoryprint(f"Propagated {len(traj_eci)} states over {traj_eci.timespan()/60:.1f} minutes")# 4. Convert to ECEF to analyze ground tracktraj_ecef=traj_eci.to_ecef()print("\nGround track in ECEF frame:")fori,(epoch,state_ecef)inenumerate(traj_ecef):ifi%10==0:# Sample every 10 states# Convert ECEF to geodetic for latitude/longitudelat,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 evolutiontraj_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}°")
usebraheasbh;usebh::time::Epoch;usebh::traits::{Trajectory,SStatePropagator};usebh::{KeplerianPropagator,orbital_period,position_ecef_to_geodetic,R_EARTH,AngleFormat};usenalgebraasna;fnmain(){bh::initialize_eop().unwrap();// 1. Define orbit and create propagatorletoe=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);letepoch=Epoch::from_datetime(2024,1,1,0,0,0.0,0.0,bh::time::TimeSystem::UTC);letmutpropagator=KeplerianPropagator::from_keplerian(epoch,oe,AngleFormat::Radians,60.0);// 2. Propagate for one orbit periodletperiod=orbital_period(R_EARTH+500e3);letend_epoch=epoch+period;propagator.propagate_to(end_epoch);// 3. Get trajectory in ECI Cartesianlettraj_eci=&propagator.trajectory;println!("Propagated {} states over {:.1} minutes",traj_eci.len(),traj_eci.timespan().unwrap()/60.0);// 4. Convert to ECEFlettraj_ecef=traj_eci.to_ecef();println!("\nGround track in ECEF frame:");for(i,(epoch,state_ecef))intraj_ecef.into_iter().enumerate(){ifi%10==0{letpos_ecef:na::Vector3<f64>=state_ecef.fixed_rows::<3>(0).into();letlla=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 Keplerianlettraj_kep=traj_eci.to_keplerian(AngleFormat::Radians);letfirst_oe=traj_kep.state_at_idx(0).unwrap();letlast_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());}