The functions described here convert between Keplerian orbital elements and Cartesian state vectors. While these transformations are part of the "coordinates" module, they specifically deal with orbital mechanics - converting between two different coordinate representations of a satellite's orbit.
Understanding both representations is essential: Keplerian elements provide intuitive orbital parameters like size, shape, and orientation, while Cartesian states are necessary for numerical orbit propagation and applying perturbations.
Keplerian elements describe an orbit using six classical parameters:
\(a\): Semi-major axis (meters) - defines the orbit's size
\(e\): Eccentricity (dimensionless) - defines the orbit's shape (0 = circular, 0 < e < 1 = elliptical)
\(i\): Inclination (radians or degrees) - tilt of orbital plane relative to equator
\(\Omega\): Right ascension of ascending node (radians or degrees) - where orbit crosses equator going north
\(\omega\): Argument of periapsis (radians or degrees) - where orbit is closest to Earth
\(M\): Mean anomaly (radians or degrees) - position of satellite along orbit
In brahe, the combined vector has ordering [a, e, i, Ω, ω, M]
Info
Brahe uses mean anomaly as the default anomaly representation for Keplerian elements. Other anomaly types (eccentric, true) can be converted using the anomaly conversion functions in the Orbits module.
Cartesian states represent position and velocity in three-dimensional space:
Position: \([p_x, p_y, p_z]\) in meters
Velocity: \([v_x, v_y, v_z]\) in meters per second
In brahe, the state vector is combined as [p_x, p_y, p_z, v_x, v_y, v_z]
Cartesian states are typically expressed in an inertial reference frame like Earth-Centered Inertial (ECI), where the axes are fixed with respect to the stars rather than rotating with Earth.
Info
All position and velocity components in Cartesian states are in SI base units (meters and meters per second).
They must be in SI base units for inputs and are always returned in SI base units.
importbraheasbhimportnumpyasnpbh.initialize_eop()# Define orbital elements [a, e, i, Ω, ω, M] in meters and degrees# LEO satellite: 500 km altitude, 97.8° inclination (~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 degreesstate=bh.state_koe_to_eci(oe_deg,bh.AngleFormat.DEGREES)print("Cartesian state [x, y, z, vx, vy, vz] (m, m/s):")print(f"Position: [{state[0]:.3f}, {state[1]:.3f}, {state[2]:.3f}]")print(f"Velocity: [{state[3]:.6f}, {state[4]:.6f}, {state[5]:.6f}]")# Cartesian state (m, m/s):# Position: [1848964.106, -434937.468, 6560410.530]# Velocity: [-7098.379734, -2173.344867, 1913.333385]
usebraheasbh;usenalgebraasna;fnmain(){bh::initialize_eop().unwrap();// Define orbital elements [a, e, i, Ω, ω, M] in meters and degrees// LEO satellite: 500 km altitude, 97.8° inclination (~sun-synchronous)letoe_deg=na::SVector::<f64,6>::new(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 degreesletstate=bh::state_koe_to_eci(oe_deg,bh::AngleFormat::Degrees);println!("Cartesian state [x, y, z, vx, vy, vz] (m, m/s):");println!("Position: [{:.3}, {:.3}, {:.3}]",state[0],state[1],state[2]);println!("Velocity: [{:.6}, {:.6}, {:.6}]",state[3],state[4],state[5]);// Cartesian state (m, m/s):// Position: [1848964.106, -434937.468, 6560410.530]// Velocity: [-7098.379734, -2173.344867, 1913.333385]}
importbraheasbhimportnumpyasnpfrommathimportpibh.initialize_eop()# Define orbital elements [a, e, i, Ω, ω, M] in meters and degrees# LEO satellite: 500 km altitude, 97.8° inclination (~sun-synchronous)oe_deg=np.array([bh.R_EARTH+500e3,# Semi-major axis (m)0.01,# Eccentricitypi/4,# Inclination (rad)pi/8,# Right ascension of ascending node (rad)pi/2,# Argument of periapsis (rad)3*pi/4,# Mean anomaly (rad)])# Convert orbital elements to Cartesian state using degreesstate=bh.state_koe_to_eci(oe_deg,bh.AngleFormat.RADIANS)print("Cartesian state [x, y, z, vx, vy, vz] (m, m/s):")print(f"Position: [{state[0]:.3f}, {state[1]:.3f}, {state[2]:.3f}]")print(f"Velocity: [{state[3]:.6f}, {state[4]:.6f}, {state[5]:.6f}]")# Cartesian state (m, m/s):# Position: [-3117582.037, -5092452.343, -3511765.495]# Velocity: [6408.435846, -1407.501408, -3752.763969]
usebraheasbh;usenalgebraasna;usestd::f64::consts::PI;fnmain(){bh::initialize_eop().unwrap();// Define orbital elements [a, e, i, Ω, ω, M] in meters and radians// LEO satellite: 500 km altitude, 45° inclinationletoe_rad=na::SVector::<f64,6>::new(bh::R_EARTH+500e3,// Semi-major axis (m)0.01,// EccentricityPI/4.0,// Inclination (rad)PI/8.0,// Right ascension of ascending node (rad)PI/2.0,// Argument of periapsis (rad)3.0*PI/4.0// Mean anomaly (rad));// Convert orbital elements to Cartesian state using radiansletstate=bh::state_koe_to_eci(oe_rad,bh::AngleFormat::Radians);println!("Cartesian state [x, y, z, vx, vy, vz] (m, m/s):");println!("Position: [{:.3}, {:.3}, {:.3}]",state[0],state[1],state[2]);println!("Velocity: [{:.6}, {:.6}, {:.6}]",state[3],state[4],state[5]);// Cartesian state (m, m/s):// Position: [-3117582.037, -5092452.343, -3511765.495]// Velocity: [6408.435846, -1407.501408, -3752.763969]}
Info
The AngleFormat parameter only affects the three angular elements (i, Ω, ω, M). Semi-major axis is always in meters, and eccentricity is always dimensionless.
importbraheasbhimportnumpyasnpbh.initialize_eop()# Define Cartesian state vector [px, py, pz, vx, vy, vz] in meters and meters per secondstate=np.array([1848964.106,-434937.468,6560410.530,-7098.379734,-2173.344867,1913.333385])# Convert orbital elements to Cartesian state using degreesoe_deg=bh.state_eci_to_koe(state,bh.AngleFormat.DEGREES)print("Osculating state [a, e, i, Ω, ω, M] (deg):")print(f"Semi-major axis (m): {oe_deg[0]:.3f}")print(f"Eccentricity: {oe_deg[1]:.6f}")print(f"Inclination (deg): {oe_deg[2]:.6f}")print(f"RA of ascending node (deg): {oe_deg[3]:.6f}")print(f"Argument of periapsis (deg): {oe_deg[4]:.6f}")print(f"Mean anomaly (deg): {oe_deg[5]:.6f}")# Osculating state (deg):# Semi-major axis (m): 6878136.299# Eccentricity: 0.010000# Inclination (deg): 97.800000# RA of ascending node (deg): 15.000000# Argument of periapsis (deg): 30.000000# Mean anomaly (deg): 45.000000# You can also convert using radiansoe_rad=bh.state_eci_to_koe(state,bh.AngleFormat.RADIANS)print("\nOsculating state [a, e, i, Ω, ω, M] (rad):")print(f"Semi-major axis (m): {oe_rad[0]:.3f}")print(f"Eccentricity: {oe_rad[1]:.6f}")print(f"Inclination (rad): {oe_rad[2]:.6f}")print(f"RA of ascending node (rad): {oe_rad[3]:.6f}")print(f"Argument of periapsis (rad): {oe_rad[4]:.6f}")print(f"Mean anomaly (rad): {oe_rad[5]:.6f}")# Osculating state (rad):# Semi-major axis (m): 6878136.299# Eccentricity: 0.010000# Inclination (rad): 1.706932# RA of ascending node (rad): 0.261799# Argument of periapsis (rad): 0.523599# Mean anomaly (rad): 0.785398
usebraheasbh;usenalgebraasna;fnmain(){bh::initialize_eop().unwrap();// Define Cartesian state vector [px, py, pz, vx, vy, vz] in meters and meters per secondletstate=na::SVector::<f64,6>::new(1848964.106,-434937.468,6560410.530,-7098.379734,-2173.344867,1913.333385);// Convert Cartesian state to orbital elements using degreesletoe_deg=bh::state_eci_to_koe(state,bh::AngleFormat::Degrees);println!("Osculating state [a, e, i, Ω, ω, M] (deg):");println!("Semi-major axis (m): {:.3}",oe_deg[0]);println!("Eccentricity: {:.6}",oe_deg[1]);println!("Inclination (deg): {:.6}",oe_deg[2]);println!("RA of ascending node (deg): {:.6}",oe_deg[3]);println!("Argument of periapsis (deg): {:.6}",oe_deg[4]);println!("Mean anomaly (deg): {:.6}",oe_deg[5]);// Osculating state (deg):// Semi-major axis (m): 6878136.299// Eccentricity: 0.010000// Inclination (deg): 97.800000// RA of ascending node (deg): 15.000000// Argument of periapsis (deg): 30.000000// Mean anomaly (deg): 45.000000// You can also convert using radiansletoe_rad=bh::state_eci_to_koe(state,bh::AngleFormat::Radians);println!("\nOsculating state [a, e, i, Ω, ω, M] (rad):");println!("Semi-major axis (m): {:.3}",oe_rad[0]);println!("Eccentricity: {:.6}",oe_rad[1]);println!("Inclination (rad): {:.6}",oe_rad[2]);println!("RA of ascending node (rad): {:.6}",oe_rad[3]);println!("Argument of periapsis (rad): {:.6}",oe_rad[4]);println!("Mean anomaly (rad): {:.6}",oe_rad[5]);// Osculating state (rad):// Semi-major axis (m): 6878136.299// Eccentricity: 0.010000// Inclination (rad): 1.706932// RA of ascending node (rad): 0.261799// Argument of periapsis (rad): 0.523599// Mean anomaly (rad): 0.785398}