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 java.util.Map;
20 import java.util.TreeMap;
21
22 import org.hipparchus.Field;
23 import org.hipparchus.RealFieldElement;
24 import org.hipparchus.analysis.polynomials.PolynomialFunction;
25 import org.hipparchus.analysis.polynomials.PolynomialsUtils;
26 import org.hipparchus.complex.Complex;
27 import org.hipparchus.exception.NullArgumentException;
28 import org.hipparchus.random.MersenneTwister;
29 import org.hipparchus.util.CombinatoricsUtils;
30 import org.hipparchus.util.Decimal64Field;
31 import org.hipparchus.util.FastMath;
32 import org.hipparchus.util.MathArrays;
33 import org.junit.Assert;
34 import org.junit.Test;
35 import org.orekit.errors.OrekitException;
36 import org.orekit.propagation.semianalytical.dsst.utilities.CoefficientsFactory.NSKey;
37
38 public class CoefficientFactoryTest {
39
40 private static final double eps0 = 0.;
41 private static final double eps10 = 1e-10;
42 private static final double eps12 = 1e-12;
43
44 /** Map of the Qns derivatives, for each (n, s) couple. */
45 private static Map<NSKey, PolynomialFunction> QNS_MAP = new TreeMap<NSKey, PolynomialFunction>();
46
47 @Test
48 public void testVns() {
49 final int order = 100;
50 TreeMap<NSKey, Double> Vns = CoefficientsFactory.computeVns(order);
51
52 // Odd terms are null
53 for (int i = 0; i < order; i++) {
54 for (int j = 0; j < i + 1; j++) {
55 if ((i - j) % 2 != 0) {
56 Assert.assertEquals(0d, Vns.get(new NSKey(i, j)), eps0);
57 }
58 }
59 }
60
61 // Check the first coefficients :
62 Assert.assertEquals(1, Vns.get(new NSKey(0, 0)), eps0);
63 Assert.assertEquals(0.5, Vns.get(new NSKey(1, 1)), eps0);
64 Assert.assertEquals(-0.5, Vns.get(new NSKey(2, 0)), eps0);
65 Assert.assertEquals(1 / 8d, Vns.get(new NSKey(2, 2)), eps0);
66 Assert.assertEquals(-1 / 8d, Vns.get(new NSKey(3, 1)), eps0);
67 Assert.assertEquals(1 / 48d, Vns.get(new NSKey(3, 3)), eps0);
68 Assert.assertEquals(3 / 8d, Vns.get(new NSKey(4, 0)), eps0);
69 Assert.assertEquals(-1 / 48d, Vns.get(new NSKey(4, 2)), eps0);
70 Assert.assertEquals(1 / 384d, Vns.get(new NSKey(4, 4)), eps0);
71 Assert.assertEquals(1 / 16d, Vns.get(new NSKey(5, 1)), eps0);
72 Assert.assertEquals(-1 / 384d, Vns.get(new NSKey(5, 3)), eps0);
73 Assert.assertEquals(1 / 3840d, Vns.get(new NSKey(5, 5)), eps0);
74 Assert.assertEquals(Vns.lastKey().getN(), order - 1);
75 Assert.assertEquals(Vns.lastKey().getS(), order - 1);
76 }
77
78 /**
79 * Test the direct computation method : the getVmns is using the Vns computation to compute the
80 * current element
81 */
82 @Test
83 public void testVmns() {
84 Assert.assertEquals(getVmns2(0, 0, 0), CoefficientsFactory.getVmns(0, 0, 0), eps0);
85 Assert.assertEquals(getVmns2(0, 1, 1), CoefficientsFactory.getVmns(0, 1, 1), eps0);
86 Assert.assertEquals(getVmns2(0, 2, 2), CoefficientsFactory.getVmns(0, 2, 2), eps0);
87 Assert.assertEquals(getVmns2(0, 3, 1), CoefficientsFactory.getVmns(0, 3, 1), eps0);
88 Assert.assertEquals(getVmns2(0, 3, 3), CoefficientsFactory.getVmns(0, 3, 3), eps0);
89 Assert.assertEquals(getVmns2(2, 2, 2), CoefficientsFactory.getVmns(2, 2, 2), eps0);
90 final double vmnsp = getVmns2(12, 26, 20);
91 Assert.assertEquals(vmnsp,
92 CoefficientsFactory.getVmns(12, 26, 20),
93 FastMath.abs(eps12 * vmnsp));
94 final double vmnsm = getVmns2(12, 27, -21);
95 Assert.assertEquals(vmnsm,
96 CoefficientsFactory.getVmns(12, 27, -21),
97 Math.abs(eps12 * vmnsm));
98 }
99
100 /** Error if m > n */
101 @Test(expected = OrekitException.class)
102 public void testVmnsError() {
103 // if m > n
104 CoefficientsFactory.getVmns(3, 2, 1);
105 }
106
107 /**
108 * Qns test based on two computation method. As methods are independent, if they give the same
109 * results, we assume them to be consistent.
110 */
111 @Test
112 public void testQns() {
113 Assert.assertEquals(1., getQnsPolynomialValue(0, 0, 0), 0.);
114 // Method comparison :
115 final int nmax = 10;
116 final int smax = 10;
117 final MersenneTwister random = new MersenneTwister(123456789);
118 for (int g = 0; g < 1000; g++) {
119 final double gamma = random.nextDouble();
120 double[][] qns = CoefficientsFactory.computeQns(gamma, nmax, smax);
121 for (int n = 0; n <= nmax; n++) {
122 final int sdim = FastMath.min(smax + 2, n);
123 for (int s = 0; s <= sdim; s++) {
124 final double qp = getQnsPolynomialValue(gamma, n, s);
125 Assert.assertEquals(qns[n][s], qp, FastMath.abs(eps10 * qns[n][s]));
126 }
127 }
128 }
129 }
130
131 @Test
132 public void testQnsField() {
133 doTestQnsField(Decimal64Field.getInstance());
134 }
135
136 /**
137 * Qns test based on two computation method. As methods are independent, if they give the same
138 * results, we assume them to be consistent.
139 */
140 private <T extends RealFieldElement<T>> void doTestQnsField(Field<T> field) {
141 final T zero = field.getZero();
142 Assert.assertEquals(1., getQnsPolynomialValue(0, 0, 0), 0.);
143 // Method comparison :
144 final int nmax = 10;
145 final int smax = 10;
146 final MersenneTwister random = new MersenneTwister(123456789);
147 for (int g = 0; g < 1000; g++) {
148 final T gamma = zero.add(random.nextDouble());
149 T[][] qns = CoefficientsFactory.computeQns(gamma, nmax, smax);
150 for (int n = 0; n <= nmax; n++) {
151 final int sdim = FastMath.min(smax + 2, n);
152 for (int s = 0; s <= sdim; s++) {
153 final T qp = getQnsPolynomialValue(gamma, n, s);
154 Assert.assertEquals(qns[n][s].getReal(), qp.getReal(), FastMath.abs(qns[n][s].multiply(eps10)).getReal());
155 }
156 }
157 }
158 }
159
160 /** Gs and Hs computation test based on 2 independent methods.
161 * If they give same results, we assume them to be consistent.
162 */
163 @Test
164 public void testGsHs() {
165 final int s = 50;
166 final MersenneTwister random = new MersenneTwister(123456789);
167 for (int i = 0; i < 10; i++) {
168 final double k = random.nextDouble();
169 final double h = random.nextDouble();
170 final double a = random.nextDouble();
171 final double b = random.nextDouble();
172 final double[][] GH = CoefficientsFactory.computeGsHs(k, h, a, b, s);
173 for (int j = 1; j < s; j++) {
174 final double[] GsHs = getGsHs(k, h, a, b, j);
175 Assert.assertEquals(GsHs[0], GH[0][j], FastMath.abs(eps12 * GsHs[0]));
176 Assert.assertEquals(GsHs[1], GH[1][j], FastMath.abs(eps12 * GsHs[1]));
177 }
178 }
179 }
180
181 @Test
182 public void testGsHsField() {
183 doTestGsHsField(Decimal64Field.getInstance());
184 }
185
186 /** Gs and Hs computation test based on 2 independent methods.
187 * If they give same results, we assume them to be consistent.
188 */
189 private <T extends RealFieldElement<T>> void doTestGsHsField(Field<T> field) {
190 final T zero = field.getZero();
191 final int s = 50;
192 final MersenneTwister random = new MersenneTwister(123456789);
193 for (int i = 0; i < 10; i++) {
194 final T k = zero.add(random.nextDouble());
195 final T h = zero.add(random.nextDouble());
196 final T a = zero.add(random.nextDouble());
197 final T b = zero.add(random.nextDouble());
198 final T[][] GH = CoefficientsFactory.computeGsHs(k, h, a, b, s, field);
199 for (int j = 1; j < s; j++) {
200 final T[] GsHs = getGsHs(k, h, a, b, j, field);
201 Assert.assertEquals(GsHs[0].getReal(), GH[0][j].getReal(), FastMath.abs(GsHs[0].multiply(eps12)).getReal());
202 Assert.assertEquals(GsHs[1].getReal(), GH[1][j].getReal(), FastMath.abs(GsHs[1].multiply(eps12)).getReal());
203 }
204 }
205 }
206
207 /**
208 * Direct computation for the Vmns coefficient from equation 2.7.1 - (6)
209 */
210 private static double getVmns2(final int m,
211 final int n,
212 final int s) {
213 double vmsn = 0d;
214 if ((n - s) % 2 == 0) {
215 final int coef = (s > 0 || s % 2 == 0) ? 1 : -1;
216 final int ss = (s > 0) ? s : -s;
217 final double num = FastMath.pow(-1, (n - ss) / 2) *
218 CombinatoricsUtils.factorialDouble(n + ss) *
219 CombinatoricsUtils.factorialDouble(n - ss);
220 final double den = FastMath.pow(2, n) *
221 CombinatoricsUtils.factorialDouble(n - m) *
222 CombinatoricsUtils.factorialDouble((n + ss) / 2) *
223 CombinatoricsUtils.factorialDouble((n - ss) / 2);
224 vmsn = coef * num / den;
225 }
226 return vmsn;
227 }
228
229 /** Get the Q<sub>ns</sub> value from 2.8.1-(4) evaluated in γ This method is using the
230 * Legendre polynomial to compute the Q<sub>ns</sub>'s one. This direct computation method
231 * allows to store the polynomials value in a static map. If the Q<sub>ns</sub> had been
232 * computed already, they just will be evaluated at γ
233 *
234 * @param gamma γ angle for which Q<sub>ns</sub> is evaluated
235 * @param n n value
236 * @param s s value
237 * @return the polynomial value evaluated at γ
238 */
239 private static double getQnsPolynomialValue(final double gamma, final int n, final int s) {
240 PolynomialFunction derivative;
241 if (QNS_MAP.containsKey(new NSKey(n, s))) {
242 derivative = QNS_MAP.get(new NSKey(n, s));
243 } else {
244 final PolynomialFunction legendre = PolynomialsUtils.createLegendrePolynomial(n);
245 derivative = legendre;
246 for (int i = 0; i < s; i++) {
247 derivative = (PolynomialFunction) derivative.polynomialDerivative();
248 }
249 QNS_MAP.put(new NSKey(n, s), derivative);
250 }
251 return derivative.value(gamma);
252 }
253
254 /** Get the Q<sub>ns</sub> value from 2.8.1-(4) evaluated in γ This method is using the
255 * Legendre polynomial to compute the Q<sub>ns</sub>'s one. This direct computation method
256 * allows to store the polynomials value in a static map. If the Q<sub>ns</sub> had been
257 * computed already, they just will be evaluated at γ
258 *
259 * @param gamma γ angle for which Q<sub>ns</sub> is evaluated
260 * @param n n value
261 * @param s s value
262 * @return the polynomial value evaluated at γ
263 */
264 private static <T extends RealFieldElement<T>> T getQnsPolynomialValue(final T gamma, final int n, final int s) {
265 PolynomialFunction derivative;
266 if (QNS_MAP.containsKey(new NSKey(n, s))) {
267 derivative = QNS_MAP.get(new NSKey(n, s));
268 } else {
269 final PolynomialFunction legendre = PolynomialsUtils.createLegendrePolynomial(n);
270 derivative = legendre;
271 for (int i = 0; i < s; i++) {
272 derivative = (PolynomialFunction) derivative.polynomialDerivative();
273 }
274 QNS_MAP.put(new NSKey(n, s), derivative);
275 }
276 return derivative.value(gamma);
277 }
278
279 /** Compute directly G<sub>s</sub> and H<sub>s</sub> coefficients from equation 3.1-(4).
280 * @param k x-component of the eccentricity vector
281 * @param h y-component of the eccentricity vector
282 * @param a 1st direction cosine
283 * @param b 2nd direction cosine
284 * @param s development order
285 * @return Array of G<sub>s</sub> and H<sub>s</sub> values for s.<br>
286 * The 1st element contains the G<sub>s</sub> value.
287 * The 2nd element contains the H<sub>s</sub> value.
288 */
289 private static double[] getGsHs(final double k, final double h,
290 final double a, final double b, final int s) {
291 final Complex as = new Complex(k, h).pow(s);
292 final Complex bs = new Complex(a, -b).pow(s);
293 final Complex asbs = as.multiply(bs);
294 return new double[] {asbs.getReal(), asbs.getImaginary()};
295 }
296
297 /** Compute directly G<sub>s</sub> and H<sub>s</sub> coefficients from equation 3.1-(4).
298 * @param k x-component of the eccentricity vector
299 * @param h y-component of the eccentricity vector
300 * @param a 1st direction cosine
301 * @param b 2nd direction cosine
302 * @param s development order
303 * @return Array of G<sub>s</sub> and H<sub>s</sub> values for s.<br>
304 * The 1st element contains the G<sub>s</sub> value.
305 * The 2nd element contains the H<sub>s</sub> value.
306 */
307 private static <T extends RealFieldElement<T>> T[] getGsHs(final T k, final T h,
308 final T a, final T b, final int s,
309 final Field<T> field) {
310 final FieldComplex<T> as = new FieldComplex<>(k, h).pow(s);
311 final FieldComplex<T> bs = new FieldComplex<>(a, b.negate()).pow(s);
312 final FieldComplex<T> asbs = as.multiply(bs);
313 final T[] values = MathArrays.buildArray(field, 2);
314 values[0] = asbs.getReal();
315 values[1] = asbs.getImaginary();
316 return values;
317 }
318
319 private static class FieldComplex <T extends RealFieldElement<T>> {
320
321 /** The imaginary part. */
322 private final T imaginary;
323
324 /** The real part. */
325 private final T real;
326
327 /**
328 * Create a complex number given the real and imaginary parts.
329 *
330 * @param real Real part.
331 * @param imaginary Imaginary part.
332 */
333 FieldComplex(final T real, final T imaginary) {
334 this.real = real;
335 this.imaginary = imaginary;
336 }
337
338 /**
339 * Access the real part.
340 *
341 * @return the real part.
342 */
343 public T getReal() {
344 return real;
345 }
346
347 /**
348 * Access the imaginary part.
349 *
350 * @return the imaginary part.
351 */
352 public T getImaginary() {
353 return imaginary;
354 }
355
356 /**
357 * Create a complex number given the real and imaginary parts.
358 *
359 * @param realPart Real part.
360 * @param imaginaryPart Imaginary part.
361 * @return a new complex number instance.
362 *
363 * @see #valueOf(double, double)
364 */
365 protected FieldComplex<T> createComplex(final T realPart, final T imaginaryPart) {
366 return new FieldComplex<>(realPart, imaginaryPart);
367 }
368
369 /**
370 * Returns a {@code Complex} whose value is {@code this * factor}.
371 * Implements preliminary checks for {@code NaN} and infinity followed by
372 * the definitional formula:
373 * <p>
374 * {@code (a + bi)(c + di) = (ac - bd) + (ad + bc)i}
375 * </p>
376 * <p>
377 * Returns finite values in components of the result per the definitional
378 * formula in all remaining cases.</p>
379 *
380 * @param factor value to be multiplied by this {@code Complex}.
381 * @return {@code this * factor}.
382 * @throws NullArgumentException if {@code factor} is {@code null}.
383 */
384 public FieldComplex<T> multiply(final FieldComplex<T> factor) throws NullArgumentException {
385 return createComplex(real.multiply(factor.real).subtract(imaginary.multiply(factor.imaginary)),
386 real.multiply(factor.imaginary).add(imaginary.multiply(factor.real)));
387 }
388
389 /**
390 * Returns a {@code Complex} whose value is {@code this * factor}, with {@code factor}
391 * interpreted as a integer number.
392 *
393 * @param factor value to be multiplied by this {@code Complex}.
394 * @return {@code this * factor}.
395 * @see #multiply(Complex)
396 */
397 public FieldComplex<T> multiply(final int factor) {
398 return createComplex(real.multiply(factor), imaginary.multiply(factor));
399 }
400
401 /**
402 * Returns of value of this complex number raised to the power of {@code x}.
403 *
404 * @param x exponent to which this {@code Complex} is to be raised.
405 * @return <code>this<sup>x</sup></code>.
406 * @see #pow(Complex)
407 */
408 public FieldComplex<T> pow(int x) {
409 return this.log().multiply(x).exp();
410 }
411
412 /**
413 * Compute the
414 * <a href="http://mathworld.wolfram.com/NaturalLogarithm.html" TARGET="_top">
415 * natural logarithm</a> of this complex number.
416 * Implements the formula:
417 * <pre>
418 * <code>
419 * log(a + bi) = ln(|a + bi|) + arg(a + bi)i
420 * </code>
421 * </pre>
422 * where ln on the right hand side is {@link FastMath#log},
423 * {@code |a + bi|} is the modulus, {@link Complex#abs}, and
424 * {@code arg(a + bi) = }{@link FastMath#atan2}(b, a).
425 * <p>
426 * Returns {@link Complex#NaN} if either real or imaginary part of the
427 * input argument is {@code NaN}.
428 * </p>
429 * Infinite (or critical) values in real or imaginary parts of the input may
430 * result in infinite or NaN values returned in parts of the result.
431 * <pre>
432 * Examples:
433 * <code>
434 * log(1 ± INFINITY i) = INFINITY ± (π/2)i
435 * log(INFINITY + i) = INFINITY + 0i
436 * log(-INFINITY + i) = INFINITY + πi
437 * log(INFINITY ± INFINITY i) = INFINITY ± (π/4)i
438 * log(-INFINITY ± INFINITY i) = INFINITY ± (3π/4)i
439 * log(0 + 0i) = -INFINITY + 0i
440 * </code>
441 * </pre>
442 *
443 * @return the value <code>ln this</code>, the natural logarithm
444 * of {@code this}.
445 */
446 public FieldComplex<T> log() {
447 return createComplex(FastMath.log(abs()),
448 FastMath.atan2(imaginary, real));
449 }
450
451 /**
452 * Compute the
453 * <a href="http://mathworld.wolfram.com/ExponentialFunction.html" TARGET="_top">
454 * exponential function</a> of this complex number.
455 * Implements the formula:
456 * <pre>
457 * <code>
458 * exp(a + bi) = exp(a)cos(b) + exp(a)sin(b)i
459 * </code>
460 * </pre>
461 * where the (real) functions on the right-hand side are
462 * {@link FastMath#exp}, {@link FastMath#cos}, and
463 * {@link FastMath#sin}.
464 * <p>
465 * Returns {@link Complex#NaN} if either real or imaginary part of the
466 * input argument is {@code NaN}.
467 * </p>
468 * Infinite values in real or imaginary parts of the input may result in
469 * infinite or NaN values returned in parts of the result.
470 * <pre>
471 * Examples:
472 * <code>
473 * exp(1 ± INFINITY i) = NaN + NaN i
474 * exp(INFINITY + i) = INFINITY + INFINITY i
475 * exp(-INFINITY + i) = 0 + 0i
476 * exp(±INFINITY ± INFINITY i) = NaN + NaN i
477 * </code>
478 * </pre>
479 *
480 * @return <code><i>e</i><sup>this</sup></code>.
481 */
482 public FieldComplex<T> exp() {
483 T expReal = FastMath.exp(real);
484 return createComplex(expReal.multiply(FastMath.cos(imaginary)),
485 expReal.multiply(FastMath.sin(imaginary)));
486 }
487
488 /**
489 * Return the absolute value of this complex number.
490 * Returns {@code NaN} if either real or imaginary part is {@code NaN}
491 * and {@code Double.POSITIVE_INFINITY} if neither part is {@code NaN},
492 * but at least one part is infinite.
493 *
494 * @return the absolute value.
495 */
496 public T abs() {
497 if (FastMath.abs(real).getReal() < FastMath.abs(imaginary).getReal()) {
498 if (imaginary.getReal() == 0.0) {
499 return FastMath.abs(real);
500 }
501 T q = real.divide(imaginary);
502 return FastMath.abs(imaginary).multiply(FastMath.sqrt(q.multiply(q).add(1.)));
503 } else {
504 if (real.getReal() == 0.0) {
505 return FastMath.abs(imaginary);
506 }
507 T q = imaginary.divide(real);
508 return FastMath.abs(real).multiply(FastMath.sqrt(q.multiply(q).add(1.)));
509 }
510 }
511
512 }
513 }