nyx_space/cosmic/
bplane.rs

1/*
2    Nyx, blazing fast astrodynamics
3    Copyright (C) 2018-onwards Christopher Rabotin <christopher.rabotin@gmail.com>
4
5    This program is free software: you can redistribute it and/or modify
6    it under the terms of the GNU Affero General Public License as published
7    by the Free Software Foundation, either version 3 of the License, or
8    (at your option) any later version.
9
10    This program is distributed in the hope that it will be useful,
11    but WITHOUT ANY WARRANTY; without even the implied warranty of
12    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
13    GNU Affero General Public License for more details.
14
15    You should have received a copy of the GNU Affero General Public License
16    along with this program.  If not, see <https://www.gnu.org/licenses/>.
17*/
18
19use anise::prelude::{Frame, Orbit};
20
21use super::{AstroError, AstroPhysicsSnafu, OrbitDual, OrbitPartial};
22use crate::cosmic::NotHyperbolicSnafu;
23use crate::linalg::{Matrix2, Matrix3, Vector2, Vector3};
24use crate::md::objective::Objective;
25use crate::md::{AstroSnafu, StateParameter, TargetingError};
26use crate::time::{Duration, Epoch, Unit};
27use crate::utils::between_pm_180;
28use hyperdual::linalg::norm;
29use hyperdual::{Float, OHyperdual};
30
31use snafu::{ensure, ResultExt};
32use std::convert::From;
33use std::fmt;
34
35/// Stores a B-Plane
36#[derive(Copy, Clone, Debug)]
37
38pub struct BPlane {
39    /// The $B_T$ component, in kilometers
40    pub b_t: OrbitPartial,
41    /// The $B_R$ component, in kilometers
42    pub b_r: OrbitPartial,
43    /// The Linearized Time of Flight
44    pub ltof_s: OrbitPartial,
45    /// The B-Plane rotation matrix
46    pub str_dcm: Matrix3<f64>,
47    /// The frame in which this B Plane was computed
48    pub frame: Frame,
49    /// The time of computation
50    pub epoch: Epoch,
51}
52
53impl BPlane {
54    /// Returns a newly define B-Plane if the orbit is hyperbolic and already in Dual form
55    pub fn from_dual(orbit: OrbitDual) -> Result<Self, AstroError> {
56        ensure!(
57            orbit.ecc().context(AstroPhysicsSnafu)?.real() > 1.0,
58            NotHyperbolicSnafu
59        );
60
61        let one = OHyperdual::from(1.0);
62        let zero = OHyperdual::from(0.0);
63
64        let e_hat =
65            orbit.evec().context(AstroPhysicsSnafu)? / orbit.ecc().context(AstroPhysicsSnafu)?.dual;
66        let h_hat = orbit.hvec() / orbit.hmag().dual;
67        let n_hat = h_hat.cross(&e_hat);
68
69        // The reals implementation (which was initially validated) was:
70        // let s = e_hat / orbit.ecc() + (1.0 - (1.0 / orbit.ecc()).powi(2)).sqrt() * n_hat;
71        // let s_hat = s / s.norm();
72
73        let ecc = orbit.ecc().context(AstroPhysicsSnafu)?;
74
75        let incoming_asymptote_fact = (one - (one / ecc.dual).powi(2)).sqrt();
76
77        let s = Vector3::new(
78            e_hat[0] / ecc.dual + incoming_asymptote_fact * n_hat[0],
79            e_hat[1] / ecc.dual + incoming_asymptote_fact * n_hat[1],
80            e_hat[2] / ecc.dual + incoming_asymptote_fact * n_hat[2],
81        );
82
83        let s_hat = s / norm(&s); // Just to make sure to renormalize everything
84
85        // The reals implementation (which was initially validated) was:
86        // let b_vec = orbit.semi_minor_axis()
87        //     * ((1.0 - (1.0 / orbit.ecc()).powi(2)).sqrt() * e_hat
88        //         - (1.0 / orbit.ecc() * n_hat));
89        let semi_minor_axis = orbit.semi_minor_axis_km().context(AstroPhysicsSnafu)?;
90
91        let b_vec = Vector3::new(
92            semi_minor_axis.dual
93                * (incoming_asymptote_fact * e_hat[0] - ((one / ecc.dual) * n_hat[0])),
94            semi_minor_axis.dual
95                * (incoming_asymptote_fact * e_hat[1] - ((one / ecc.dual) * n_hat[1])),
96            semi_minor_axis.dual
97                * (incoming_asymptote_fact * e_hat[2] - ((one / ecc.dual) * n_hat[2])),
98        );
99
100        let t = s_hat.cross(&Vector3::new(zero, zero, one));
101        let t_hat = t / norm(&t);
102        let r_hat = s_hat.cross(&t_hat);
103
104        // Build the rotation matrix from inertial to B Plane
105        let str_rot = Matrix3::new(
106            s_hat[0].real(),
107            s_hat[1].real(),
108            s_hat[2].real(),
109            t_hat[0].real(),
110            t_hat[1].real(),
111            t_hat[2].real(),
112            r_hat[0].real(),
113            r_hat[1].real(),
114            r_hat[2].real(),
115        );
116
117        Ok(BPlane {
118            b_r: OrbitPartial {
119                dual: b_vec.dot(&r_hat),
120                param: StateParameter::BdotR,
121            },
122            b_t: OrbitPartial {
123                dual: b_vec.dot(&t_hat),
124                param: StateParameter::BdotT,
125            },
126            ltof_s: OrbitPartial {
127                dual: b_vec.dot(&s_hat) / orbit.vmag_km_s().dual,
128                param: StateParameter::BLTOF,
129            },
130            str_dcm: str_rot,
131            frame: orbit.frame,
132            epoch: orbit.dt,
133        })
134    }
135
136    /// Returns a newly defined B-Plane if the orbit is hyperbolic.
137    pub fn new(orbit: Orbit) -> Result<Self, AstroError> {
138        // Convert to OrbitDual so we can target it
139        Self::from_dual(OrbitDual::from(orbit))
140    }
141
142    /// Returns the DCM to convert to the B Plane from the inertial frame
143    pub fn inertial_to_bplane(&self) -> Matrix3<f64> {
144        self.str_dcm
145    }
146
147    /// Returns the Jacobian of the B plane (BR, BT, LTOF) with respect to the velocity
148    pub fn jacobian(&self) -> Matrix3<f64> {
149        Matrix3::new(
150            self.b_r.wtr_vx(),
151            self.b_r.wtr_vy(),
152            self.b_r.wtr_vz(),
153            self.b_t.wtr_vx(),
154            self.b_t.wtr_vy(),
155            self.b_t.wtr_vz(),
156            self.ltof_s.wtr_vx(),
157            self.ltof_s.wtr_vy(),
158            self.ltof_s.wtr_vz(),
159        )
160    }
161
162    /// Returns the Jacobian of the B plane (BR, BT) with respect to two of the velocity components
163    pub fn jacobian2(&self, invariant: StateParameter) -> Result<Matrix2<f64>, AstroError> {
164        match invariant {
165            StateParameter::VX => Ok(Matrix2::new(
166                self.b_r.wtr_vy(),
167                self.b_r.wtr_vz(),
168                self.b_t.wtr_vy(),
169                self.b_t.wtr_vz(),
170            )),
171            StateParameter::VY => Ok(Matrix2::new(
172                self.b_r.wtr_vx(),
173                self.b_r.wtr_vz(),
174                self.b_t.wtr_vx(),
175                self.b_t.wtr_vz(),
176            )),
177            StateParameter::VZ => Ok(Matrix2::new(
178                self.b_r.wtr_vx(),
179                self.b_r.wtr_vy(),
180                self.b_t.wtr_vx(),
181                self.b_t.wtr_vy(),
182            )),
183            _ => Err(AstroError::BPlaneInvariant),
184        }
185    }
186}
187
188impl BPlane {
189    pub fn b_dot_t(&self) -> f64 {
190        self.b_t.real()
191    }
192
193    pub fn b_dot_r(&self) -> f64 {
194        self.b_r.real()
195    }
196
197    pub fn ltof(&self) -> Duration {
198        self.ltof_s.real() * Unit::Second
199    }
200
201    /// Returns the B plane angle in degrees between -180 and 180
202    pub fn angle(&self) -> f64 {
203        between_pm_180(self.b_dot_r().atan2(self.b_dot_t()).to_degrees())
204    }
205
206    /// Returns the B plane vector magnitude, in kilometers
207    pub fn mag(&self) -> f64 {
208        (self.b_dot_t().powi(2) + self.b_dot_r().powi(2)).sqrt()
209    }
210}
211
212impl fmt::Display for BPlane {
213    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
214        write!(
215            f,
216            "[{}] {} B-Plane: B∙R = {:.3} km\tB∙T = {:.3} km\tAngle = {:.3} deg",
217            self.frame,
218            self.epoch,
219            self.b_dot_r(),
220            self.b_dot_t(),
221            self.angle()
222        )
223    }
224}
225
226#[derive(Copy, Clone, Debug)]
227pub struct BPlaneTarget {
228    /// The $B_T$ component, in kilometers
229    pub b_t_km: f64,
230    /// The $B_R$ component, in kilometers
231    pub b_r_km: f64,
232    /// The Linearized Time of Flight, in seconds
233    pub ltof_s: f64,
234
235    /// The tolerance on the $B_T$ component, in kilometers
236    pub tol_b_t_km: f64,
237    /// The tolerance on the $B_R$ component, in kilometers
238    pub tol_b_r_km: f64,
239    /// The tolerance on the Linearized Time of Flight, in seconds
240    pub tol_ltof_s: f64,
241}
242
243impl BPlaneTarget {
244    /// Initializes a new B Plane target with only the targets and the default tolerances.
245    /// Default tolerances are 1 millimeter in positions and 1 second in LTOF
246    pub fn from_targets(b_t_km: f64, b_r_km: f64, ltof: Duration) -> Self {
247        let tol_ltof: Duration = 6.0 * Unit::Hour;
248        Self {
249            b_t_km,
250            b_r_km,
251            ltof_s: ltof.to_seconds(),
252            tol_b_t_km: 1e-6,
253            tol_b_r_km: 1e-6,
254            tol_ltof_s: tol_ltof.to_seconds(),
255        }
256    }
257
258    /// Initializes a new B Plane target with only the B Plane targets (not LTOF constraint) and the default tolerances.
259    /// Default tolerances are 1 millimeter in positions. Here, the LTOF tolerance is set to 100 days.
260    pub fn from_bt_br(b_t_km: f64, b_r_km: f64) -> Self {
261        let ltof_tol: Duration = 100 * Unit::Day;
262        Self {
263            b_t_km,
264            b_r_km,
265            ltof_s: 0.0,
266            tol_b_t_km: 1e-6,
267            tol_b_r_km: 1e-6,
268            tol_ltof_s: ltof_tol.to_seconds(),
269        }
270    }
271
272    pub fn ltof_target_set(&self) -> bool {
273        self.ltof_s.abs() > 1e-10
274    }
275
276    pub fn to_objectives(self) -> [Objective; 2] {
277        self.to_objectives_with_tolerance(1.0)
278    }
279
280    pub fn to_objectives_with_tolerance(self, tol_km: f64) -> [Objective; 2] {
281        [
282            Objective::within_tolerance(StateParameter::BdotR, self.b_r_km, tol_km),
283            Objective::within_tolerance(StateParameter::BdotT, self.b_t_km, tol_km),
284        ]
285    }
286
287    /// Includes the linearized time of flight as an objective
288    pub fn to_all_objectives_with_tolerance(self, tol_km: f64) -> [Objective; 3] {
289        [
290            Objective::within_tolerance(StateParameter::BdotR, self.b_r_km, tol_km),
291            Objective::within_tolerance(StateParameter::BdotT, self.b_t_km, tol_km),
292            Objective::within_tolerance(StateParameter::BLTOF, self.ltof_s, self.tol_ltof_s * 1e5),
293        ]
294    }
295}
296
297impl fmt::Display for BPlaneTarget {
298    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
299        write!(
300            f,
301            "B-Plane target: B∙R = {:.3} km (+/- {:.1} m)\tB∙T = {:.3} km (+/- {:.1} m)",
302            self.b_r_km,
303            self.tol_b_r_km * 1e-3,
304            self.b_t_km,
305            self.tol_b_t_km * 1e-3,
306        )
307    }
308}
309
310/// Returns the Delta V (in km/s) needed to achieve the B Plane specified by B dot R and B dot T.
311/// If no LTOF target is set, this method will fix VX, VY and VZ successively and use the minimum of those as a seed for the LTOF variation finding.
312/// If the 3x3 search is worse than any of the 2x2s, then a 2x2 will be returned.
313/// This uses the hyperdual formulation of the Jacobian and will also vary the linearize time of flight (LTOF).
314pub fn try_achieve_b_plane(
315    orbit: Orbit,
316    target: BPlaneTarget,
317) -> Result<(Vector3<f64>, BPlane), TargetingError> {
318    let mut total_dv = Vector3::zeros();
319    let mut attempt_no = 0;
320    let max_iter = 10;
321
322    let mut real_orbit = orbit;
323    let mut prev_b_plane_err = f64::INFINITY;
324
325    if !target.ltof_target_set() {
326        // If no LTOF is targeted, we'll solve this with a least squared approach.
327        loop {
328            if attempt_no > max_iter {
329                return Err(TargetingError::TooManyIterations);
330            }
331
332            // Build current B Plane
333            let b_plane = BPlane::new(real_orbit).context(AstroSnafu)?;
334
335            // Check convergence
336            let br_err = target.b_r_km - b_plane.b_dot_r();
337            let bt_err = target.b_t_km - b_plane.b_dot_t();
338
339            if br_err.abs() < target.tol_b_r_km && bt_err.abs() < target.tol_b_t_km {
340                return Ok((total_dv, b_plane));
341            }
342
343            // Build the error vector
344            let b_plane_err = Vector2::new(br_err, bt_err);
345
346            if b_plane_err.norm() >= prev_b_plane_err {
347                // If the error is not going down, we'll raise an error
348                return Err(TargetingError::CorrectionIneffective {
349                    prev_val: prev_b_plane_err,
350                    cur_val: b_plane_err.norm(),
351                    action: "Delta-V correction is ineffective at reducing the B-Plane error",
352                });
353            }
354            prev_b_plane_err = b_plane_err.norm();
355
356            // Grab the first two rows of the Jacobian (discard the rest).
357            let full_jac = b_plane.jacobian();
358            let jac = full_jac.fixed_rows::<2>(0);
359            // Solve the Least Squares / compute the delta-v
360            let dv = jac.transpose() * (jac * jac.transpose()).try_inverse().unwrap() * b_plane_err;
361
362            total_dv[0] += dv[0];
363            total_dv[1] += dv[1];
364            total_dv[2] += dv[2];
365
366            // Rebuild a new orbit
367            real_orbit.velocity_km_s.x += dv[0];
368            real_orbit.velocity_km_s.y += dv[1];
369            real_orbit.velocity_km_s.z += dv[2];
370
371            attempt_no += 1;
372        }
373    } else {
374        // The LTOF targeting seems to break often, but it's still implemented
375        loop {
376            if attempt_no > max_iter {
377                return Err(TargetingError::TooManyIterations);
378            }
379
380            // Build current B Plane
381            let b_plane = BPlane::new(real_orbit).context(AstroSnafu)?;
382
383            // Check convergence
384            let br_err = target.b_r_km - b_plane.b_dot_r();
385            let bt_err = target.b_t_km - b_plane.b_dot_t();
386            let ltof_err = target.ltof_s - b_plane.ltof_s.real();
387
388            if br_err.abs() < target.tol_b_r_km
389                && bt_err.abs() < target.tol_b_t_km
390                && ltof_err.abs() < target.tol_ltof_s
391            {
392                return Ok((total_dv, b_plane));
393            }
394
395            // Build the error vector
396            let b_plane_err = Vector3::new(br_err, bt_err, ltof_err);
397
398            if b_plane_err.norm() >= prev_b_plane_err {
399                return Err(TargetingError::CorrectionIneffective {
400                    prev_val: prev_b_plane_err,
401                    cur_val: b_plane_err.norm(),
402                    action: "LTOF enabled correction is failing. Try to not set an LTOF target.",
403                });
404            }
405            prev_b_plane_err = b_plane_err.norm();
406
407            // Compute the delta-v
408            let dv = b_plane.jacobian() * b_plane_err;
409
410            total_dv[0] += dv[0];
411            total_dv[1] += dv[1];
412            total_dv[2] += dv[2];
413
414            // Rebuild a new orbit
415            real_orbit.velocity_km_s.x += dv[0];
416            real_orbit.velocity_km_s.y += dv[1];
417            real_orbit.velocity_km_s.z += dv[2];
418
419            attempt_no += 1;
420        }
421    }
422}