1use core::f64::consts::PI;
20
21use crate::errors::LambertError;
22
23use super::{LambertInput, LambertSolution, TransferKind, LAMBERT_EPSILON, MAX_ITERATIONS};
24
25pub fn izzo(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 ri_norm = r_init.norm();
51 let rf_norm = r_final.norm();
52
53 let chord = r_final - r_init;
54 let c_norm = chord.norm();
55
56 let s = (ri_norm + rf_norm + c_norm) * 0.5;
58
59 let i_r1 = r_init / ri_norm;
61 let i_r2 = r_final / rf_norm;
62
63 let mut i_h = i_r1.cross(&i_r2);
64 i_h /= i_h.norm(); if i_h.z < 0.0 {
70 i_h = -i_h;
71 }
72
73 let lambda = 1.0 - c_norm / s;
75 let mut m_lambda = lambda.sqrt();
76
77 let revs = match kind {
78 TransferKind::NRevs(n) => n.into(),
79 _ => 0_u32,
80 };
81
82 let dm = if matches!(kind, TransferKind::NRevs(_)) {
87 1.0
88 } else {
89 kind.direction_of_motion(&r_final, &r_init)?
90 };
91 let retrograde = dm < 0.0;
92
93 let (mut i_t1, mut i_t2) = if retrograde {
94 m_lambda = -m_lambda;
95 (i_r1.cross(&i_h), i_r2.cross(&i_h))
96 } else {
97 (i_h.cross(&i_r1), i_h.cross(&i_r2))
98 };
99 i_t1 /= i_t1.norm();
101 i_t2 /= i_t2.norm();
102
103 let t = (2.0 * mu_km3_s2 / s.powi(3)).sqrt() * tof_s;
105
106 let (x, y) = find_xy(
108 m_lambda,
109 t,
110 revs,
111 MAX_ITERATIONS,
112 LAMBERT_EPSILON,
113 LAMBERT_EPSILON,
114 true,
115 )?;
116 let gamma = (mu_km3_s2 * s / 2.0).sqrt();
118 let rho = (ri_norm - rf_norm) / c_norm;
119 let sigma = (1.0 - rho.powi(2)).sqrt();
120
121 let (v_r1, v_r2, v_t1, v_t2) = reconstruct(x, y, ri_norm, rf_norm, m_lambda, gamma, rho, sigma);
124
125 let v_init = v_r1 * (r_init / ri_norm) + v_t1 * i_t1;
127 let v_final = v_r2 * (r_final / rf_norm) + v_t2 * i_t2;
128
129 Ok(LambertSolution {
130 v_init_km_s: v_init,
131 v_final_km_s: v_final,
132 phi_rad: 0.0,
133 input,
134 })
135}
136
137pub fn reconstruct(
140 x: f64,
141 y: f64,
142 r1: f64,
143 r2: f64,
144 ll: f64,
145 gamma: f64,
146 rho: f64,
147 sigma: f64,
148) -> (f64, f64, f64, f64) {
149 let v_r1 = gamma * ((ll * y - x) - rho * (ll * y + x)) / r1;
150 let v_r2 = -gamma * ((ll * y - x) + rho * (ll * y + x)) / r2;
151 let v_t1 = gamma * sigma * (y + ll * x) / r1;
152 let v_t2 = gamma * sigma * (y + ll * x) / r2;
153
154 (v_r1, v_r2, v_t1, v_t2)
155}
156
157pub fn find_xy(
163 ll: f64,
164 t: f64,
165 m: u32,
166 maxiter: usize,
167 atol: f64,
168 rtol: f64,
169 low_path: bool,
170) -> Result<(f64, f64), LambertError> {
171 assert!(ll.abs() < 1.0, "|ll| must be less than 1");
174
175 let mut m_max = (t / PI).floor() as u32;
176 let t_00 = ll.acos() + ll * (1.0 - ll.powi(2)).sqrt(); if t < t_00 + (m_max as f64) * PI && m_max > 0 {
180 let (_, t_min) = compute_t_min(ll, m_max, maxiter, atol, rtol)?;
181 if t < t_min {
182 m_max -= 1;
183 }
184 }
185
186 if m > m_max {
188 return Err(LambertError::MultiRevNotFeasible { m, m_max });
189 }
190
191 let x_0 = initial_guess(t, ll, m, low_path);
193
194 let x = householder(x_0, t, ll, m, atol, rtol, maxiter)?;
197 let y = compute_y(x, ll);
198
199 Ok((x, y))
200}
201
202fn compute_y(x: f64, ll: f64) -> f64 {
203 (1.0 - ll.powi(2) * (1.0 - x.powi(2))).sqrt()
204}
205
206fn compute_psi(x: f64, y: f64, ll: f64) -> f64 {
210 if x > 1.0 {
211 ((y - x * ll) * (x.powi(2) - 1.0).sqrt()).asinh()
213 } else if x < 1.0 {
214 (x * y + ll * (1.0 - x.powi(2))).acos()
217 } else {
218 0.0
220 }
221}
222
223fn tof_equation(x: f64, t0: f64, ll: f64, m: u32) -> f64 {
225 let y = compute_y(x, ll);
226 tof_equation_y(x, y, t0, ll, m)
227}
228
229fn tof_equation_y(x: f64, y: f64, t0: f64, ll: f64, m: u32) -> f64 {
231 let t_ = if m == 0 && x.powi(2) > 0.6 && x.powi(2) < 1.4 {
232 let eta = y - ll * x;
233 let s_1 = (1.0 - ll - x * eta) * 0.5;
234 let q = 4.0 / 3.0 * hyp2f1b(s_1);
235 (eta.powi(3) * q + 4.0 * ll * eta) * 0.5
236 } else {
237 let psi = compute_psi(x, y, ll);
238 let m_float = m as f64;
239 let den = 1.0 - x.powi(2);
240
241 ((psi + m_float * PI) / den.abs().sqrt() - x + ll * y) / den
244 };
245
246 t_ - t0
247}
248
249fn tof_equation_p(x: f64, y: f64, t: f64, ll: f64) -> f64 {
251 (3.0 * t * x - 2.0 + 2.0 * ll.powi(3) * x / y) / (1.0 - x.powi(2))
252}
253
254fn tof_equation_p2(x: f64, y: f64, t: f64, dt: f64, ll: f64) -> f64 {
256 (3.0 * t + 5.0 * x * dt + 2.0 * (1.0 - ll.powi(2)) * ll.powi(3) / y.powi(3)) / (1.0 - x.powi(2))
257}
258
259fn tof_equation_p3(x: f64, y: f64, _t: f64, dt: f64, ddt: f64, ll: f64) -> f64 {
261 (7.0 * x * ddt + 8.0 * dt - 6.0 * (1.0 - ll.powi(2)) * ll.powi(5) * x / y.powi(5))
262 / (1.0 - x.powi(2))
263}
264
265fn compute_t_min(
268 ll: f64,
269 m: u32,
270 maxiter: usize,
271 atol: f64,
272 rtol: f64,
273) -> Result<(f64, f64), LambertError> {
274 if (ll - 1.0).abs() < 1e-9 {
276 let x_t_min = 0.0;
277 let t_min = tof_equation(x_t_min, 0.0, ll, m);
278 Ok((x_t_min, t_min))
279 } else if m == 0 {
280 let x_t_min = f64::INFINITY;
281 let t_min = 0.0;
282 Ok((x_t_min, t_min))
283 } else {
284 let x_i = 0.1;
286 let t_i = tof_equation(x_i, 0.0, ll, m);
287 let x_t_min = halley(x_i, t_i, ll, atol, rtol, maxiter)?;
288 let t_min = tof_equation(x_t_min, 0.0, ll, m);
289 Ok((x_t_min, t_min))
290 }
291}
292
293fn initial_guess(tof_s: f64, ll: f64, m: u32, low_path: bool) -> f64 {
310 if m == 0 {
311 let t_0 = ll.acos() + ll * (1.0 - ll.powi(2)).sqrt(); let t_1 = 2.0 * (1.0 - ll.powi(3)) / 3.0; if tof_s >= t_0 {
317 (t_0 / tof_s).powf(2.0 / 3.0) - 1.0
318 } else if tof_s < t_1 {
319 5.0 / 2.0 * t_1 / tof_s * (t_1 - tof_s) / (1.0 - ll.powi(5)) + 1.0
320 } else {
321 (2.0f64.ln() * (tof_s / t_0).ln() / (t_1 / t_0).ln()).exp() - 1.0
324 }
325 } else {
326 let m_float = m as f64;
328
329 let term_l = ((m_float * PI + PI) / (8.0 * tof_s)).powf(2.0 / 3.0);
332 let x_0l = (term_l - 1.0) / (term_l + 1.0);
333
334 let term_r = ((8.0 * tof_s) / (m_float * PI)).powf(2.0 / 3.0);
335 let x_0r = (term_r - 1.0) / (term_r + 1.0);
336
337 if low_path {
339 x_0l.max(x_0r)
340 } else {
341 x_0l.min(x_0r)
342 }
343 }
344}
345
346pub fn hyp2f1b(x: f64) -> f64 {
348 if x >= 1.0 {
349 return f64::INFINITY;
350 }
351
352 let mut res = 1.0;
353 let mut term = 1.0;
354 let mut ii = 0.0_f64; loop {
356 term *= (3.0 + ii) * (1.0 + ii) / (2.5 + ii) * x / (ii + 1.0);
357 let res_old = res;
358 res += term;
359
360 if res == res_old {
362 return res;
363 }
364 ii += 1.0;
365 }
366}
367
368pub fn halley(
373 mut p0: f64,
374 t0: f64,
375 ll: f64,
376 atol: f64,
377 rtol: f64,
378 maxiter: usize,
379) -> Result<f64, LambertError> {
380 for _ in 1..=maxiter {
381 let y = compute_y(p0, ll);
382 let fder = tof_equation_p(p0, y, t0, ll);
383 let fder2 = tof_equation_p2(p0, y, t0, fder, ll);
384
385 if fder2.abs() < 1e-14 {
387 return Err(LambertError::TargetsTooClose);
388 }
389
390 let fder3 = tof_equation_p3(p0, y, t0, fder, fder2, ll);
391
392 let p = p0 - 2.0 * fder * fder2 / (2.0 * fder2.powi(2) - fder * fder3);
394
395 if (p - p0).abs() < rtol * p0.abs() + atol {
396 return Ok(p);
397 }
398 p0 = p;
399 }
400
401 Err(LambertError::SolverMaxIter { maxiter })
402}
403
404pub fn householder(
407 mut p0: f64,
408 t0: f64,
409 ll: f64,
410 m: u32,
411 atol: f64,
412 rtol: f64,
413 maxiter: usize,
414) -> Result<f64, LambertError> {
415 for _ in 1..=maxiter {
416 let y = compute_y(p0, ll);
417 let fval = tof_equation_y(p0, y, t0, ll, m);
418 let t = fval + t0;
419 let fder = tof_equation_p(p0, y, t, ll);
420 let fder2 = tof_equation_p2(p0, y, t, fder, ll);
421 let fder3 = tof_equation_p3(p0, y, t, fder, fder2, ll);
422
423 let num = fder.powi(2) - fval * fder2 / 2.0;
425 let den = fder * (fder.powi(2) - fval * fder2) + fder3 * fval.powi(2) / 6.0;
426
427 if den.abs() < 1e-14 {
428 return Err(LambertError::TargetsTooClose);
429 }
430
431 let p = p0 - fval * (num / den);
432
433 if (p - p0).abs() < rtol * p0.abs() + atol {
434 return Ok(p);
435 }
436 p0 = p;
437 }
438 Err(LambertError::SolverMaxIter { maxiter })
439}
440
441#[cfg(test)]
442mod ut_lambert_izzo {
443 use super::{izzo, LambertInput, TransferKind};
444 use crate::linalg::Vector3;
445 use anise::prelude::{Frame, Orbit};
446 use hifitime::{Epoch, Unit};
447
448 #[test]
449 fn test_lambert_izzo_shortway() {
450 let frame = Frame {
453 ephemeris_id: 301,
454 orientation_id: 0,
455 mu_km3_s2: Some(3.98600433e5),
456 shape: None,
457 frozen_epoch: None,
458 force_inertial: false,
459 };
460 let initial_state = Orbit {
461 radius_km: Vector3::new(15945.34, 0.0, 0.0),
462 velocity_km_s: Vector3::zeros(),
463 epoch: Epoch::from_gregorian_utc_at_midnight(2025, 1, 1),
464 frame,
465 };
466 let final_state = Orbit {
467 radius_km: Vector3::new(12214.83899, 10249.46731, 0.0),
468 velocity_km_s: Vector3::zeros(),
469 epoch: Epoch::from_gregorian_utc_at_midnight(2025, 1, 1) + Unit::Minute * 76.0,
470 frame,
471 };
472
473 let input = LambertInput::from_planetary_states(initial_state, final_state).unwrap();
474
475 let exp_vi = Vector3::new(2.058913, 2.915965, 0.0);
476 let exp_vf = Vector3::new(-3.451565, 0.910315, 0.0);
477
478 let sol = izzo(input, TransferKind::ShortWay).unwrap();
479
480 println!("{sol:?}\t{exp_vi}\t{exp_vf}");
481
482 assert!((sol.v_init_km_s - exp_vi).norm() < 1e-5);
483 assert!((sol.v_final_km_s - exp_vf).norm() < 1e-5);
484
485 assert!(
486 sol.v_inf_outgoing_declination_deg().abs() < 1e-8,
487 "outgoing decl should be zero"
488 );
489 assert!(
490 (sol.v_inf_outgoing_right_ascension_deg() - 54.7747288).abs() < 1e-6,
491 "wrong outgoing RA"
492 );
493 }
494
495 #[test]
496 fn test_lambert_izzo_longway() {
497 let frame = Frame {
499 ephemeris_id: 301,
500 orientation_id: 0,
501 mu_km3_s2: Some(3.98600433e5),
502 shape: None,
503 frozen_epoch: None,
504 force_inertial: false,
505 };
506 let initial_state = Orbit {
507 radius_km: Vector3::new(15945.34, 0.0, 0.0),
508 velocity_km_s: Vector3::zeros(),
509 epoch: Epoch::from_gregorian_utc_at_midnight(2025, 1, 1),
510 frame,
511 };
512 let final_state = Orbit {
513 radius_km: Vector3::new(12214.83899, 10249.46731, 0.0),
514 velocity_km_s: Vector3::zeros(),
515 epoch: Epoch::from_gregorian_utc_at_midnight(2025, 1, 1) + Unit::Minute * 76.0,
516 frame,
517 };
518
519 let input = LambertInput::from_planetary_states(initial_state, final_state).unwrap();
520
521 let exp_vi = Vector3::new(-3.811158, -2.003854, 0.0);
522 let exp_vf = Vector3::new(4.207569, 0.914724, 0.0);
523
524 let sol = izzo(input, TransferKind::LongWay).unwrap();
525
526 assert!((sol.v_init_km_s - exp_vi).norm() < 1e-5);
527 assert!((sol.v_final_km_s - exp_vf).norm() < 1e-5);
528
529 assert!(
530 sol.v_inf_outgoing_declination_deg().abs() < 1e-8,
531 "outgoing decl should be zero"
532 );
533 assert!(
534 (sol.v_inf_outgoing_right_ascension_deg() + 152.2652291).abs() < 1e-6,
535 "wrong outgoing RA"
536 );
537 }
538
539 #[test]
540 fn test_lambert_izzo_direction_toggle_ih_negative() {
541 let frame = Frame {
545 ephemeris_id: 301,
546 orientation_id: 0,
547 mu_km3_s2: Some(3.98600433e5),
548 shape: None,
549 frozen_epoch: None,
550 force_inertial: false,
551 };
552 let t0 = Epoch::from_gregorian_utc_at_midnight(2025, 1, 1);
553
554 let r1 = Vector3::new(7000.0, 0.0, 0.0);
556 let r2 = Vector3::new(0.0, -7000.0, -500.0);
557
558 let initial = Orbit {
559 radius_km: r1,
560 velocity_km_s: Vector3::zeros(),
561 epoch: t0,
562 frame,
563 };
564 let final_ = Orbit {
565 radius_km: r2,
566 velocity_km_s: Vector3::zeros(),
567 epoch: t0 + Unit::Minute * 60.0,
568 frame,
569 };
570 let input = LambertInput::from_planetary_states(initial, final_).unwrap();
571
572 let short = izzo(input, TransferKind::ShortWay).unwrap();
573 let long = izzo(input, TransferKind::LongWay).unwrap();
574
575 assert!(
576 (short.v_init_km_s - long.v_init_km_s).norm() > 1e-3,
577 "direction toggle ignored — short and long-way collapsed: short={:?}, long={:?}",
578 short.v_init_km_s,
579 long.v_init_km_s,
580 );
581 }
582
583 #[test]
584 fn test_lambert_izzo_auto_matches_longway_for_long_arc() {
585 let frame = Frame {
586 ephemeris_id: 301,
587 orientation_id: 0,
588 mu_km3_s2: Some(3.98600433e5),
589 shape: None,
590 frozen_epoch: None,
591 force_inertial: false,
592 };
593 let t0 = Epoch::from_gregorian_utc_at_midnight(2025, 1, 1);
594
595 let initial = Orbit {
597 radius_km: Vector3::new(8000.0, 0.0, 0.0),
598 velocity_km_s: Vector3::zeros(),
599 epoch: t0,
600 frame,
601 };
602 let final_ = Orbit {
603 radius_km: Vector3::new(0.0, -8000.0, 0.0),
604 velocity_km_s: Vector3::zeros(),
605 epoch: t0 + Unit::Minute * 60.0,
606 frame,
607 };
608 let input = LambertInput::from_planetary_states(initial, final_).unwrap();
609
610 let auto = izzo(input, TransferKind::Auto).unwrap();
611 let long = izzo(input, TransferKind::LongWay).unwrap();
612 let short = izzo(input, TransferKind::ShortWay).unwrap();
613
614 assert!(
615 (auto.v_init_km_s - long.v_init_km_s).norm() < 1e-10,
616 "Auto should match LongWay for long-arc geometry",
617 );
618 assert!(
619 (auto.v_init_km_s - short.v_init_km_s).norm() > 1e-3,
620 "Auto should differ from ShortWay for long-arc geometry",
621 );
622 }
623
624 #[test]
625 fn test_lambert_izzo_auto_matches_shortway_for_short_arc() {
626 let frame = Frame {
627 ephemeris_id: 301,
628 orientation_id: 0,
629 mu_km3_s2: Some(3.98600433e5),
630 shape: None,
631 frozen_epoch: None,
632 force_inertial: false,
633 };
634 let t0 = Epoch::from_gregorian_utc_at_midnight(2025, 1, 1);
635
636 let initial = Orbit {
638 radius_km: Vector3::new(8000.0, 0.0, 0.0),
639 velocity_km_s: Vector3::zeros(),
640 epoch: t0,
641 frame,
642 };
643 let final_ = Orbit {
644 radius_km: Vector3::new(0.0, 8000.0, 0.0),
645 velocity_km_s: Vector3::zeros(),
646 epoch: t0 + Unit::Minute * 60.0,
647 frame,
648 };
649 let input = LambertInput::from_planetary_states(initial, final_).unwrap();
650
651 let auto = izzo(input, TransferKind::Auto).unwrap();
652 let short = izzo(input, TransferKind::ShortWay).unwrap();
653 let long = izzo(input, TransferKind::LongWay).unwrap();
654
655 assert!(
656 (auto.v_init_km_s - short.v_init_km_s).norm() < 1e-10,
657 "Auto should match ShortWay for short-arc geometry",
658 );
659 assert!(
660 (auto.v_init_km_s - long.v_init_km_s).norm() > 1e-3,
661 "Auto should differ from LongWay for short-arc geometry",
662 );
663 }
664}