use crate::AngleFormat;
use crate::coordinates::{state_eci_to_koe, state_koe_to_eci};
use crate::math::SVector6;
use crate::relative_motion::{state_oe_to_roe, state_roe_to_oe};
pub fn state_eci_to_roe(
x_chief: SVector6,
x_deputy: SVector6,
angle_format: AngleFormat,
) -> SVector6 {
let oe_chief = state_eci_to_koe(x_chief, angle_format);
let oe_deputy = state_eci_to_koe(x_deputy, angle_format);
state_oe_to_roe(oe_chief, oe_deputy, angle_format)
}
pub fn state_roe_to_eci(x_chief: SVector6, roe: SVector6, angle_format: AngleFormat) -> SVector6 {
let oe_chief = state_eci_to_koe(x_chief, angle_format);
let oe_deputy = state_roe_to_oe(oe_chief, roe, angle_format);
state_koe_to_eci(oe_deputy, angle_format)
}
#[cfg(test)]
#[cfg_attr(coverage_nightly, coverage(off))]
mod tests {
use super::*;
use crate::constants::R_EARTH;
use crate::coordinates::state_koe_to_eci;
use crate::relative_motion::state_oe_to_roe;
use approx::assert_abs_diff_eq;
#[test]
fn test_state_eci_to_roe_degrees() {
let oe_chief = SVector6::new(R_EARTH + 700e3, 0.001, 97.8, 15.0, 30.0, 45.0);
let oe_deputy = SVector6::new(R_EARTH + 701e3, 0.0015, 97.85, 15.05, 30.05, 45.05);
let x_chief = state_koe_to_eci(oe_chief, AngleFormat::Degrees);
let x_deputy = state_koe_to_eci(oe_deputy, AngleFormat::Degrees);
let roe_from_eci = state_eci_to_roe(x_chief, x_deputy, AngleFormat::Degrees);
let roe_from_oe = state_oe_to_roe(oe_chief, oe_deputy, AngleFormat::Degrees);
assert_abs_diff_eq!(roe_from_eci[0], roe_from_oe[0], epsilon = 1e-10);
assert_abs_diff_eq!(roe_from_eci[1], roe_from_oe[1], epsilon = 1e-10);
assert_abs_diff_eq!(roe_from_eci[2], roe_from_oe[2], epsilon = 1e-10);
assert_abs_diff_eq!(roe_from_eci[3], roe_from_oe[3], epsilon = 1e-10);
assert_abs_diff_eq!(roe_from_eci[4], roe_from_oe[4], epsilon = 1e-10);
assert_abs_diff_eq!(roe_from_eci[5], roe_from_oe[5], epsilon = 1e-10);
}
#[test]
fn test_state_roe_to_eci_degrees() {
let oe_chief = SVector6::new(R_EARTH + 700e3, 0.001, 97.8, 15.0, 30.0, 45.0);
let x_chief = state_koe_to_eci(oe_chief, AngleFormat::Degrees);
let roe = SVector6::new(0.000142857, 0.05, 0.0005, -0.0003, 0.01, -0.02);
let x_deputy = state_roe_to_eci(x_chief, roe, AngleFormat::Degrees);
let pos_mag = (x_deputy[0].powi(2) + x_deputy[1].powi(2) + x_deputy[2].powi(2)).sqrt();
assert!(pos_mag > R_EARTH); assert!(pos_mag < R_EARTH + 2000e3);
let vel_mag = (x_deputy[3].powi(2) + x_deputy[4].powi(2) + x_deputy[5].powi(2)).sqrt();
assert!(vel_mag > 6000.0); assert!(vel_mag < 9000.0);
}
#[test]
fn test_state_eci_roe_roundtrip_degrees() {
let oe_chief = SVector6::new(R_EARTH + 700e3, 0.001, 97.8, 15.0, 30.0, 45.0);
let oe_deputy = SVector6::new(R_EARTH + 701e3, 0.0015, 97.85, 15.05, 30.05, 45.05);
let x_chief = state_koe_to_eci(oe_chief, AngleFormat::Degrees);
let x_deputy_orig = state_koe_to_eci(oe_deputy, AngleFormat::Degrees);
let roe = state_eci_to_roe(x_chief, x_deputy_orig, AngleFormat::Degrees);
let x_deputy_recovered = state_roe_to_eci(x_chief, roe, AngleFormat::Degrees);
assert_abs_diff_eq!(x_deputy_recovered[0], x_deputy_orig[0], epsilon = 1e-3);
assert_abs_diff_eq!(x_deputy_recovered[1], x_deputy_orig[1], epsilon = 1e-3);
assert_abs_diff_eq!(x_deputy_recovered[2], x_deputy_orig[2], epsilon = 1e-3);
assert_abs_diff_eq!(x_deputy_recovered[3], x_deputy_orig[3], epsilon = 1e-6);
assert_abs_diff_eq!(x_deputy_recovered[4], x_deputy_orig[4], epsilon = 1e-6);
assert_abs_diff_eq!(x_deputy_recovered[5], x_deputy_orig[5], epsilon = 1e-6);
}
#[test]
fn test_state_eci_to_roe_radians() {
use crate::constants::DEG2RAD;
let oe_chief = SVector6::new(
R_EARTH + 700e3,
0.001,
97.8 * DEG2RAD,
15.0 * DEG2RAD,
30.0 * DEG2RAD,
45.0 * DEG2RAD,
);
let oe_deputy = SVector6::new(
R_EARTH + 701e3,
0.0015,
97.85 * DEG2RAD,
15.05 * DEG2RAD,
30.05 * DEG2RAD,
45.05 * DEG2RAD,
);
let x_chief = state_koe_to_eci(oe_chief, AngleFormat::Radians);
let x_deputy = state_koe_to_eci(oe_deputy, AngleFormat::Radians);
let roe_from_eci = state_eci_to_roe(x_chief, x_deputy, AngleFormat::Radians);
let roe_from_oe = state_oe_to_roe(oe_chief, oe_deputy, AngleFormat::Radians);
assert_abs_diff_eq!(roe_from_eci[0], roe_from_oe[0], epsilon = 1e-10);
assert_abs_diff_eq!(roe_from_eci[1], roe_from_oe[1], epsilon = 1e-10);
assert_abs_diff_eq!(roe_from_eci[2], roe_from_oe[2], epsilon = 1e-10);
assert_abs_diff_eq!(roe_from_eci[3], roe_from_oe[3], epsilon = 1e-10);
assert_abs_diff_eq!(roe_from_eci[4], roe_from_oe[4], epsilon = 1e-10);
assert_abs_diff_eq!(roe_from_eci[5], roe_from_oe[5], epsilon = 1e-10);
}
#[test]
fn test_state_eci_roe_roundtrip_radians() {
use crate::constants::DEG2RAD;
let oe_chief = SVector6::new(
R_EARTH + 700e3,
0.001,
97.8 * DEG2RAD,
15.0 * DEG2RAD,
30.0 * DEG2RAD,
45.0 * DEG2RAD,
);
let oe_deputy = SVector6::new(
R_EARTH + 701e3,
0.0015,
97.85 * DEG2RAD,
15.05 * DEG2RAD,
30.05 * DEG2RAD,
45.05 * DEG2RAD,
);
let x_chief = state_koe_to_eci(oe_chief, AngleFormat::Radians);
let x_deputy_orig = state_koe_to_eci(oe_deputy, AngleFormat::Radians);
let roe = state_eci_to_roe(x_chief, x_deputy_orig, AngleFormat::Radians);
let x_deputy_recovered = state_roe_to_eci(x_chief, roe, AngleFormat::Radians);
assert_abs_diff_eq!(x_deputy_recovered[0], x_deputy_orig[0], epsilon = 1e-3);
assert_abs_diff_eq!(x_deputy_recovered[1], x_deputy_orig[1], epsilon = 1e-3);
assert_abs_diff_eq!(x_deputy_recovered[2], x_deputy_orig[2], epsilon = 1e-3);
assert_abs_diff_eq!(x_deputy_recovered[3], x_deputy_orig[3], epsilon = 1e-6);
assert_abs_diff_eq!(x_deputy_recovered[4], x_deputy_orig[4], epsilon = 1e-6);
assert_abs_diff_eq!(x_deputy_recovered[5], x_deputy_orig[5], epsilon = 1e-6);
}
}