Skip to main content

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