nyx_space/tools/lambert/
godding.rs1use super::{
20 LambertError, LambertInput, LambertSolution, TransferKind, LAMBERT_EPSILON,
21 LAMBERT_EPSILON_RAD, LAMBERT_EPSILON_TIME, MAX_ITERATIONS,
22};
23
24use core::f64::consts::PI;
25
26pub fn gooding(input: LambertInput, kind: TransferKind) -> Result<LambertSolution, LambertError> {
45 let r_init = input.initial_state.radius_km;
46 let r_final = input.final_state.radius_km;
47 let tof_s = (input.final_state.epoch - input.initial_state.epoch).to_seconds();
48 let mu_km3_s2 = input.mu_km2_s3();
49
50 let r_init_norm = r_init.norm();
51 let r_final_norm = r_final.norm();
52 let r_norm_product = r_init_norm * r_final_norm;
53 let cos_dnu = r_init.dot(&r_final) / r_norm_product;
54
55 let dm = kind.direction_of_motion(&r_final, &r_init)?;
56
57 let nu_init = r_init[1].atan2(r_init[0]);
58 let nu_final = r_final[1].atan2(r_final[0]);
59
60 let a = dm * (r_norm_product * (1.0 + cos_dnu)).sqrt();
61
62 if nu_final - nu_init < LAMBERT_EPSILON_RAD && a.abs() < LAMBERT_EPSILON {
63 return Err(LambertError::TargetsTooClose);
64 }
65
66 let mut phi_upper = 4.0 * PI.powi(2);
67 let mut phi_lower = -4.0 * PI.powi(2);
68 let mut phi = 0.0;
69
70 let mut c2: f64 = 1.0 / 2.0;
71 let mut c3: f64 = 1.0 / 6.0;
72 let mut iter: usize = 0;
73 let mut cur_tof: f64 = 0.0;
74 let mut y = 0.0;
75
76 while (cur_tof - tof_s).abs() > LAMBERT_EPSILON_TIME {
77 if iter > MAX_ITERATIONS {
78 return Err(LambertError::SolverMaxIter {
79 maxiter: MAX_ITERATIONS,
80 });
81 }
82 iter += 1;
83
84 y = r_init_norm + r_final_norm + a * (phi * c3 - 1.0) / c2.sqrt();
85 if a > 0.0 && y < 0.0 {
86 for _ in 0..500 {
87 phi += 0.1;
88 y = r_init_norm + r_final_norm + a * (phi * c3 - 1.0) / c2.sqrt();
89 if y >= 0.0 {
90 break;
91 }
92 }
93 if y < 0.0 {
94 return Err(LambertError::NotReasonablePhi);
95 }
96 }
97
98 let chi = (y / c2).sqrt();
99 cur_tof = (chi.powi(3) * c3 + a * y.sqrt()) / mu_km3_s2.sqrt();
100
101 if cur_tof < tof_s {
102 phi_lower = phi;
103 } else {
104 phi_upper = phi;
105 }
106
107 phi = (phi_upper + phi_lower) / 2.0;
108
109 if phi > LAMBERT_EPSILON {
110 let sqrt_phi = phi.sqrt();
111 let (s_sphi, c_sphi) = sqrt_phi.sin_cos();
112 c2 = (1.0 - c_sphi) / phi;
113 c3 = (sqrt_phi - s_sphi) / phi.powi(3).sqrt();
114 } else if phi < -LAMBERT_EPSILON {
115 let sqrt_phi = (-phi).sqrt();
116 c2 = (1.0 - sqrt_phi.cosh()) / phi;
117 c3 = (sqrt_phi.sinh() - sqrt_phi) / (-phi).powi(3).sqrt();
118 } else {
119 c2 = 0.5;
120 c3 = 1.0 / 6.0;
121 }
122 }
123
124 let f = 1.0 - y / r_init_norm;
125 let g_dot = 1.0 - y / r_final_norm;
126 let g = a * (y / mu_km3_s2).sqrt();
127
128 Ok(LambertSolution {
129 v_init_km_s: (r_final - f * r_init) / g,
130 v_final_km_s: (1.0 / g) * (g_dot * r_final - r_init),
131 phi_rad: phi,
132 input,
133 })
134}
135
136#[cfg(test)]
137mod ut_lambert_gooding {
138
139 use super::{gooding, TransferKind};
140 use crate::{linalg::Vector3, tools::lambert::LambertInput};
141 use anise::{frames::Frame, prelude::Orbit};
142 use hifitime::{Epoch, Unit};
143 #[test]
144 fn test_lambert_vallado_shortway() {
145 let frame = Frame {
146 ephemeris_id: 301,
147 orientation_id: 0,
148 mu_km3_s2: Some(3.98600433e5),
149 shape: None,
150 };
151 let initial_state = Orbit {
152 radius_km: Vector3::new(15945.34, 0.0, 0.0),
153 velocity_km_s: Vector3::zeros(),
154 epoch: Epoch::from_gregorian_utc_at_midnight(2025, 1, 1),
155 frame,
156 };
157 let final_state = Orbit {
158 radius_km: Vector3::new(12214.83899, 10249.46731, 0.0),
159 velocity_km_s: Vector3::zeros(),
160 epoch: Epoch::from_gregorian_utc_at_midnight(2025, 1, 1) + Unit::Minute * 76.0,
161 frame,
162 };
163
164 let input = LambertInput::from_planetary_states(initial_state, final_state).unwrap();
165
166 let exp_vi = Vector3::new(2.058913, 2.915965, 0.0);
167 let exp_vf = Vector3::new(-3.451565, 0.910315, 0.0);
168
169 let sol = gooding(input, TransferKind::ShortWay).unwrap();
170
171 assert!((sol.v_init_km_s - exp_vi).norm() < 1e-6);
172 assert!((sol.v_final_km_s - exp_vf).norm() < 1e-6);
173
174 println!("{sol:?}");
175 println!("{}\n{}", sol.transfer_orbit(), sol.arrival_orbit());
176 println!(
177 "v_inf+ = {} km/s -- c3 = {} km^2/s^2",
178 sol.v_inf_depart_km_s().norm(),
179 sol.c3_km2_s2()
180 );
181 println!("v_inf- = {} km/s", sol.v_inf_arrive_km_s().norm());
182 }
183
184 #[test]
185 fn test_lambert_vallado_longway() {
186 let frame = Frame {
187 ephemeris_id: 301,
188 orientation_id: 0,
189 mu_km3_s2: Some(3.98600433e5),
190 shape: None,
191 };
192 let initial_state = Orbit {
193 radius_km: Vector3::new(15945.34, 0.0, 0.0),
194 velocity_km_s: Vector3::zeros(),
195 epoch: Epoch::from_gregorian_utc_at_midnight(2025, 1, 1),
196 frame,
197 };
198 let final_state = Orbit {
199 radius_km: Vector3::new(12214.83899, 10249.46731, 0.0),
200 velocity_km_s: Vector3::zeros(),
201 epoch: Epoch::from_gregorian_utc_at_midnight(2025, 1, 1) + Unit::Minute * 76.0,
202 frame,
203 };
204
205 let input = LambertInput::from_planetary_states(initial_state, final_state).unwrap();
206 let exp_vi = Vector3::new(-3.811158, -2.003854, 0.0);
207 let exp_vf = Vector3::new(4.207569, 0.914724, 0.0);
208
209 let sol = gooding(input, TransferKind::LongWay).unwrap();
210
211 assert!((sol.v_init_km_s - exp_vi).norm() < 1e-6);
212 assert!((sol.v_final_km_s - exp_vf).norm() < 1e-6);
213 }
214}