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.
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.
importbraheasbhimportnumpyasnp# Initialize EOPbh.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,# Eccentricity97.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 degreesx_deg=bh.state_koe_to_eci(oe_deg,bh.AngleFormat.DEGREES)# Convert back to degreesoe_deg_2=bh.state_eci_to_koe(x_deg,bh.AngleFormat.DEGREES)print("Original Keplerian elements:")fori,eleminenumerate(oe_deg):print(f" [{i}]: {elem:.3f}")print("Converted Cartesian state:")fori,eleminenumerate(x_deg):print(f" [{i}]: {elem:.3f}")print("Back to Keplerian elements:")fori,eleminenumerate(oe_deg_2):print(f" [{i}]: {elem:.3f}")
usebraheasbh;usenalgebraasna;fnmain(){// Initialize EOPbh::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)letoe_deg=na::vector![bh::R_EARTH+500e3,// Semi-major axis (m)0.01,// Eccentricity97.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 degreesletx_deg=bh::state_koe_to_eci(oe_deg,bh::AngleFormat::Degrees);// Convert back to degreesletoe_deg_2=bh::state_eci_to_koe(x_deg,bh::AngleFormat::Degrees);println!("Original Keplerian elements:");for(i,elem)inoe_deg.iter().enumerate(){println!(" [{i}]: {elem:.3}");}println!("Converted Cartesian state:");for(i,elem)inx_deg.iter().enumerate(){println!(" [{i}]: {elem:.3}");}println!("Back to Keplerian elements:");for(i,elem)inoe_deg_2.iter().enumerate(){println!(" [{i}]: {elem:.3}");}}
importbraheasbhimportnumpyasnp# Initialize EOPbh.initialize_eop()# Create an epochepc=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,# Eccentricity97.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 degreesx_deg=bh.state_koe_to_eci(oe_deg,bh.AngleFormat.DEGREES)# Covert ECI cartesian state to ECEF cartesian statex_ecef=bh.state_eci_to_ecef(epc,x_deg)# Convert ECEF cartesian state to geodetic coordinateslat,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")
usebraheasbh;usenalgebraasna;fnmain(){// Initialize EOPbh::initialize_eop().unwrap();// Create an epochletepc=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)letoe_deg=na::vector![bh::R_EARTH+500e3,// Semi-major axis (m)0.01,// Eccentricity97.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 degreesletx_deg=bh::state_koe_to_eci(oe_deg,bh::AngleFormat::Degrees);// Convert ECI cartesian state to ECEF cartesian stateletx_ecef=bh::state_eci_to_ecef(epc,x_deg);// Convert ECEF cartesian state to geodetic coordinatesletgeodetic=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);}