1 /* Copyright 2002-2020 CS GROUP
2 * Licensed to CS GROUP (CS) under one or more
3 * contributor license agreements. See the NOTICE file distributed with
4 * this work for additional information regarding copyright ownership.
5 * CS licenses this file to You under the Apache License, Version 2.0
6 * (the "License"); you may not use this file except in compliance with
7 * the License. You may obtain a copy of the License at
8 *
9 * http://www.apache.org/licenses/LICENSE-2.0
10 *
11 * Unless required by applicable law or agreed to in writing, software
12 * distributed under the License is distributed on an "AS IS" BASIS,
13 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 * See the License for the specific language governing permissions and
15 * limitations under the License.
16 */
17 package org.orekit.propagation.semianalytical.dsst.utilities;
18
19 import org.hipparchus.geometry.euclidean.threed.Vector3D;
20 import org.hipparchus.util.FastMath;
21 import org.hipparchus.util.MathUtils;
22 import org.orekit.frames.Frame;
23 import org.orekit.orbits.Orbit;
24 import org.orekit.time.AbsoluteDate;
25
26
27 /** Container class for common parameters used by all DSST forces.
28 * <p>
29 * Most of them are defined in Danielson paper at § 2.1.
30 * </p>
31 * @author Pascal Parraud
32 */
33 public class AuxiliaryElements {
34
35 /** Orbit date. */
36 private final AbsoluteDate date;
37
38 /** Orbit frame. */
39 private final Frame frame;
40
41 /** Eccentricity. */
42 private final double ecc;
43
44 /** Keplerian mean motion. */
45 private final double n;
46
47 /** Keplerian period. */
48 private final double period;
49
50 /** Semi-major axis. */
51 private final double sma;
52
53 /** x component of eccentricity vector. */
54 private final double k;
55
56 /** y component of eccentricity vector. */
57 private final double h;
58
59 /** x component of inclination vector. */
60 private final double q;
61
62 /** y component of inclination vector. */
63 private final double p;
64
65 /** Mean longitude. */
66 private final double lm;
67
68 /** True longitude. */
69 private final double lv;
70
71 /** Eccentric longitude. */
72 private final double le;
73
74 /** Retrograde factor I.
75 * <p>
76 * DSST model needs equinoctial orbit as internal representation.
77 * Classical equinoctial elements have discontinuities when inclination
78 * is close to zero. In this representation, I = +1. <br>
79 * To avoid this discontinuity, another representation exists and equinoctial
80 * elements can be expressed in a different way, called "retrograde" orbit.
81 * This implies I = -1. <br>
82 * As Orekit doesn't implement the retrograde orbit, I is always set to +1.
83 * But for the sake of consistency with the theory, the retrograde factor
84 * has been kept in the formulas.
85 * </p>
86 */
87 private final int I;
88
89 /** B = sqrt(1 - h² - k²). */
90 private final double B;
91
92 /** C = 1 + p² + q². */
93 private final double C;
94
95 /** Equinoctial frame f vector. */
96 private final Vector3D f;
97
98 /** Equinoctial frame g vector. */
99 private final Vector3D g;
100
101 /** Equinoctial frame w vector. */
102 private final Vector3D w;
103
104 /** Direction cosine α. */
105 private final double alpha;
106
107 /** Direction cosine β. */
108 private final double beta;
109
110 /** Direction cosine γ. */
111 private final double gamma;
112
113 /** Simple constructor.
114 * @param orbit related mean orbit for auxiliary elements
115 * @param retrogradeFactor retrograde factor I [Eq. 2.1.2-(2)]
116 */
117 public AuxiliaryElements(final Orbit orbit, final int retrogradeFactor) {
118
119 // Date of the orbit
120 date = orbit.getDate();
121
122 // Orbit definition frame
123 frame = orbit.getFrame();
124
125 // Eccentricity
126 ecc = orbit.getE();
127
128 // Keplerian mean motion
129 n = orbit.getKeplerianMeanMotion();
130
131 // Keplerian period
132 period = orbit.getKeplerianPeriod();
133
134 // Equinoctial elements [Eq. 2.1.2-(1)]
135 sma = orbit.getA();
136 k = orbit.getEquinoctialEx();
137 h = orbit.getEquinoctialEy();
138 q = orbit.getHx();
139 p = orbit.getHy();
140 lm = MathUtils.normalizeAngle(orbit.getLM(), FastMath.PI);
141 lv = MathUtils.normalizeAngle(orbit.getLv(), FastMath.PI);
142 le = MathUtils.normalizeAngle(orbit.getLE(), FastMath.PI);
143
144 // Retrograde factor [Eq. 2.1.2-(2)]
145 I = retrogradeFactor;
146
147 final double k2 = k * k;
148 final double h2 = h * h;
149 final double q2 = q * q;
150 final double p2 = p * p;
151
152 // A, B, C parameters [Eq. 2.1.6-(1)]
153 B = FastMath.sqrt(1 - k2 - h2);
154 C = 1 + q2 + p2;
155
156 // Equinoctial reference frame [Eq. 2.1.4-(1)]
157 final double ooC = 1. / C;
158 final double px2 = 2. * p;
159 final double qx2 = 2. * q;
160 final double pq2 = px2 * q;
161 f = new Vector3D(ooC, new Vector3D(1. - p2 + q2, pq2, -px2 * I));
162 g = new Vector3D(ooC, new Vector3D(pq2 * I, (1. + p2 - q2) * I, qx2));
163 w = new Vector3D(ooC, new Vector3D(px2, -qx2, (1. - p2 - q2) * I));
164
165 // Direction cosines for central body [Eq. 2.1.9-(1)]
166 alpha = f.getZ();
167 beta = g.getZ();
168 gamma = w.getZ();
169 }
170
171 /** Get the date of the orbit.
172 * @return the date
173 */
174 public AbsoluteDate getDate() {
175 return date;
176 }
177
178 /** Get the definition frame of the orbit.
179 * @return the definition frame
180 */
181 public Frame getFrame() {
182 return frame;
183 }
184
185 /** Get the eccentricity.
186 * @return ecc
187 */
188 public double getEcc() {
189 return ecc;
190 }
191
192 /** Get the Keplerian mean motion.
193 * @return n
194 */
195 public double getMeanMotion() {
196 return n;
197 }
198
199 /** Get the Keplerian period.
200 * @return period
201 */
202 public double getKeplerianPeriod() {
203 return period;
204 }
205
206 /** Get the semi-major axis.
207 * @return the semi-major axis a
208 */
209 public double getSma() {
210 return sma;
211 }
212
213 /** Get the x component of eccentricity vector.
214 * <p>
215 * This element called k in DSST corresponds to ex
216 * for the {@link org.orekit.orbits.EquinoctialOrbit}
217 * </p>
218 * @return k
219 */
220 public double getK() {
221 return k;
222 }
223
224 /** Get the y component of eccentricity vector.
225 * <p>
226 * This element called h in DSST corresponds to ey
227 * for the {@link org.orekit.orbits.EquinoctialOrbit}
228 * </p>
229 * @return h
230 */
231 public double getH() {
232 return h;
233 }
234
235 /** Get the x component of inclination vector.
236 * <p>
237 * This element called q in DSST corresponds to hx
238 * for the {@link org.orekit.orbits.EquinoctialOrbit}
239 * </p>
240 * @return q
241 */
242 public double getQ() {
243 return q;
244 }
245
246 /** Get the y component of inclination vector.
247 * <p>
248 * This element called p in DSST corresponds to hy
249 * for the {@link org.orekit.orbits.EquinoctialOrbit}
250 * </p>
251 * @return p
252 */
253 public double getP() {
254 return p;
255 }
256
257 /** Get the mean longitude.
258 * @return lm
259 */
260 public double getLM() {
261 return lm;
262 }
263
264 /** Get the true longitude.
265 * @return lv
266 */
267 public double getLv() {
268 return lv;
269 }
270
271 /** Get the eccentric longitude.
272 * @return lf
273 */
274 public double getLf() {
275 return le;
276 }
277
278 /** Get the retrograde factor.
279 * @return the retrograde factor I
280 */
281 public int getRetrogradeFactor() {
282 return I;
283 }
284
285 /** Get B = sqrt(1 - e²).
286 * @return B
287 */
288 public double getB() {
289 return B;
290 }
291
292 /** Get C = 1 + p² + q².
293 * @return C
294 */
295 public double getC() {
296 return C;
297 }
298
299 /** Get equinoctial frame vector f.
300 * @return f vector
301 */
302 public Vector3D getVectorF() {
303 return f;
304 }
305
306 /** Get equinoctial frame vector g.
307 * @return g vector
308 */
309 public Vector3D getVectorG() {
310 return g;
311 }
312
313 /** Get equinoctial frame vector w.
314 * @return w vector
315 */
316 public Vector3D getVectorW() {
317 return w;
318 }
319
320 /** Get direction cosine α for central body.
321 * @return α
322 */
323 public double getAlpha() {
324 return alpha;
325 }
326
327 /** Get direction cosine β for central body.
328 * @return β
329 */
330 public double getBeta() {
331 return beta;
332 }
333
334 /** Get direction cosine γ for central body.
335 * @return γ
336 */
337 public double getGamma() {
338 return gamma;
339 }
340
341 }