Skip to main content

nyx_space/dynamics/
drag.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::almanac::Almanac;
20use anise::constants::frames::IAU_EARTH_FRAME;
21use snafu::ResultExt;
22
23use super::{
24    DynamicsAlmanacSnafu, DynamicsAstroSnafu, DynamicsError, DynamicsPlanetarySnafu, ForceModel,
25};
26use crate::cosmic::{AstroError, AstroPhysicsSnafu, Frame, Spacecraft};
27use crate::linalg::{Matrix4x3, Vector3};
28use serde::{Deserialize, Serialize};
29use serde_dhall::StaticType;
30use std::fmt;
31use std::sync::Arc;
32
33/// Density in kg/m^3 and altitudes in meters, not kilometers!
34#[derive(Clone, Copy, Debug, Serialize, Deserialize, StaticType)]
35pub enum AtmDensity {
36    Constant(f64),
37    Exponential { rho0: f64, r0: f64, ref_alt_m: f64 },
38    StdAtm { max_alt_m: f64 },
39}
40
41/// `ConstantDrag` implements a constant drag model as defined in Vallado, 4th ed., page 551, with an important caveat.
42///
43/// **WARNING:** This basic model assumes that the velocity of the spacecraft is identical to the velocity of the upper atmosphere.
44/// This is a **bad** assumption and **should not** be used for high fidelity simulations.
45/// This will be resolved after https://github.com/nyx-space/nyx/issues/317 is implemented.
46#[derive(Clone)]
47pub struct ConstantDrag {
48    /// atmospheric density in kg/m^3
49    pub rho: f64,
50    /// Frame causing the drag
51    pub drag_frame: Frame,
52    /// Set to true to estimate the coefficient of drag
53    pub estimate: bool,
54}
55
56impl fmt::Display for ConstantDrag {
57    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
58        write!(
59            f,
60            "\tConstant Drag rho = {} kg/m^3 in frame {}",
61            self.rho, self.drag_frame
62        )
63    }
64}
65
66impl ForceModel for ConstantDrag {
67    fn estimation_index(&self) -> Option<usize> {
68        if self.estimate {
69            Some(7)
70        } else {
71            None
72        }
73    }
74
75    fn eom(&self, ctx: &Spacecraft, almanac: Arc<Almanac>) -> Result<Vector3<f64>, DynamicsError> {
76        let osc = almanac
77            .transform_to(ctx.orbit, self.drag_frame, None)
78            .context(DynamicsAlmanacSnafu {
79                action: "transforming into drag frame",
80            })?;
81
82        let velocity = osc.velocity_km_s;
83        // Note the 1e3 factor to convert drag units from ((kg * km^2 * s^-2) / m^1) to (kg * km * s^-2)
84        Ok(-0.5
85            * 1e3
86            * self.rho
87            * ctx.drag.coeff_drag
88            * ctx.drag.area_m2
89            * velocity.norm()
90            * velocity)
91    }
92
93    fn dual_eom(
94        &self,
95        _osc_ctx: &Spacecraft,
96        _almanac: Arc<Almanac>,
97    ) -> Result<(Vector3<f64>, Matrix4x3<f64>), DynamicsError> {
98        Err(DynamicsError::DynamicsAstro {
99            source: AstroError::PartialsUndefined,
100        })
101    }
102}
103
104/// `Drag` implements all three drag models.
105#[derive(Copy, Clone, Debug, Serialize, Deserialize, StaticType)]
106pub struct Drag {
107    /// Density computation method
108    pub density: AtmDensity,
109    /// Frame to compute the drag in
110    pub drag_frame: Frame,
111    /// Set to true to estimate the coefficient of drag
112    pub estimate: bool,
113}
114
115impl Drag {
116    /// Common exponential drag model for the Earth
117    pub fn earth_exp(almanac: Arc<Almanac>) -> Result<Arc<Self>, DynamicsError> {
118        Ok(Arc::new(Self {
119            density: AtmDensity::Exponential {
120                rho0: 3.614e-13,
121                r0: 700_000.0,
122                ref_alt_m: 88_667.0,
123            },
124            drag_frame: almanac.frame_info(IAU_EARTH_FRAME).context({
125                DynamicsPlanetarySnafu {
126                    action: "planetary data from third body not loaded",
127                }
128            })?,
129            estimate: false,
130        }))
131    }
132
133    /// Drag model which uses the standard atmosphere 1976 model for atmospheric density
134    pub fn std_atm1976(almanac: Arc<Almanac>) -> Result<Arc<Self>, DynamicsError> {
135        Ok(Arc::new(Self {
136            density: AtmDensity::StdAtm {
137                max_alt_m: 1_000_000.0,
138            },
139            drag_frame: almanac.frame_info(IAU_EARTH_FRAME).context({
140                DynamicsPlanetarySnafu {
141                    action: "planetary data from third body not loaded",
142                }
143            })?,
144            estimate: false,
145        }))
146    }
147}
148
149impl fmt::Display for Drag {
150    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
151        write!(
152            f,
153            "\tDrag density {:?} in frame {}",
154            self.density, self.drag_frame
155        )
156    }
157}
158
159impl ForceModel for Drag {
160    fn estimation_index(&self) -> Option<usize> {
161        if self.estimate {
162            Some(7)
163        } else {
164            None
165        }
166    }
167
168    fn eom(&self, ctx: &Spacecraft, almanac: Arc<Almanac>) -> Result<Vector3<f64>, DynamicsError> {
169        let integration_frame = ctx.orbit.frame;
170
171        let osc_drag_frame = almanac
172            .transform_to(ctx.orbit, self.drag_frame, None)
173            .context(DynamicsAlmanacSnafu {
174                action: "transforming into drag frame",
175            })?;
176
177        match self.density {
178            AtmDensity::Constant(rho) => {
179                let velocity = osc_drag_frame.velocity_km_s;
180                // Note the 1e3 factor to convert drag units from ((kg * km^2 * s^-2) / m^1) to (kg * km * s^-2)
181                Ok(-0.5
182                    * 1e3
183                    * rho
184                    * ctx.drag.coeff_drag
185                    * ctx.drag.area_m2
186                    * velocity.norm()
187                    * velocity)
188            }
189
190            AtmDensity::Exponential {
191                rho0,
192                r0,
193                ref_alt_m,
194            } => {
195                // Compute rho in the drag frame.
196                let rho = rho0
197                    * (-(osc_drag_frame.rmag_km()
198                        - (r0
199                            + self
200                                .drag_frame
201                                .mean_equatorial_radius_km()
202                                .context(AstroPhysicsSnafu)
203                                .context(DynamicsAstroSnafu)?))
204                        / ref_alt_m)
205                        .exp();
206
207                // TODO: Drag modeling will be improved in https://github.com/nyx-space/nyx/issues/317
208                // The frame will be double checked in this PR as well.
209                let velocity_integr_frame = almanac
210                    .transform_to(osc_drag_frame, integration_frame, None)
211                    .context(DynamicsAlmanacSnafu {
212                        action: "rotating into the integration frame",
213                    })?
214                    .velocity_km_s;
215
216                let velocity = velocity_integr_frame - osc_drag_frame.velocity_km_s;
217                // Note the 1e3 factor to convert drag units from ((kg * km^2 * s^-2) / m^1) to (kg * km * s^-2)
218                Ok(-0.5
219                    * 1e3
220                    * rho
221                    * ctx.drag.coeff_drag
222                    * ctx.drag.area_m2
223                    * velocity.norm()
224                    * velocity)
225            }
226
227            AtmDensity::StdAtm { max_alt_m } => {
228                let altitude_km = osc_drag_frame.rmag_km()
229                    - self
230                        .drag_frame
231                        .mean_equatorial_radius_km()
232                        .context(AstroPhysicsSnafu)
233                        .context(DynamicsAstroSnafu)?;
234                let rho = if altitude_km > max_alt_m / 1_000.0 {
235                    // Use a constant density
236                    10.0_f64.powf((-7e-5) * altitude_km - 14.464)
237                } else {
238                    // Code from AVS/Schaub's Basilisk
239                    // Calculating the density based on a scaled 6th order polynomial fit to the log of density
240                    let scale = (altitude_km - 526.8000) / 292.8563;
241                    let logdensity =
242                        0.34047 * scale.powi(6) - 0.5889 * scale.powi(5) - 0.5269 * scale.powi(4)
243                            + 1.0036 * scale.powi(3)
244                            + 0.60713 * scale.powi(2)
245                            - 2.3024 * scale
246                            - 12.575;
247
248                    /* Calculating density by raising 10 to the log of density */
249                    10.0_f64.powf(logdensity)
250                };
251
252                let velocity_integr_frame = almanac
253                    .transform_to(osc_drag_frame, integration_frame, None)
254                    .context(DynamicsAlmanacSnafu {
255                        action: "rotating into the integration frame",
256                    })?
257                    .velocity_km_s;
258
259                let velocity = velocity_integr_frame - osc_drag_frame.velocity_km_s;
260                // Note the 1e3 factor to convert drag units from ((kg * km^2 * s^-2) / m^1) to (kg * km * s^-2)
261                Ok(-0.5
262                    * 1e3
263                    * rho
264                    * ctx.drag.coeff_drag
265                    * ctx.drag.area_m2
266                    * velocity.norm()
267                    * velocity)
268            }
269        }
270    }
271
272    fn dual_eom(
273        &self,
274        _osc_ctx: &Spacecraft,
275        _almanac: Arc<Almanac>,
276    ) -> Result<(Vector3<f64>, Matrix4x3<f64>), DynamicsError> {
277        Err(DynamicsError::DynamicsAstro {
278            source: AstroError::PartialsUndefined,
279        })
280    }
281}