1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17 package org.orekit.forces.gravity;
18
19
20 import java.util.Collections;
21 import java.util.List;
22
23 import org.hipparchus.CalculusFieldElement;
24 import org.hipparchus.analysis.differentiation.DerivativeStructure;
25 import org.hipparchus.analysis.differentiation.Gradient;
26 import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
27 import org.hipparchus.geometry.euclidean.threed.SphericalCoordinates;
28 import org.hipparchus.geometry.euclidean.threed.Vector3D;
29 import org.hipparchus.linear.Array2DRowRealMatrix;
30 import org.hipparchus.linear.RealMatrix;
31 import org.hipparchus.util.FastMath;
32 import org.hipparchus.util.MathArrays;
33 import org.orekit.forces.ForceModel;
34 import org.orekit.forces.gravity.potential.NormalizedSphericalHarmonicsProvider;
35 import org.orekit.forces.gravity.potential.NormalizedSphericalHarmonicsProvider.NormalizedSphericalHarmonics;
36 import org.orekit.forces.gravity.potential.TideSystem;
37 import org.orekit.forces.gravity.potential.TideSystemProvider;
38 import org.orekit.frames.FieldStaticTransform;
39 import org.orekit.frames.Frame;
40 import org.orekit.frames.StaticTransform;
41 import org.orekit.propagation.FieldSpacecraftState;
42 import org.orekit.propagation.SpacecraftState;
43 import org.orekit.time.AbsoluteDate;
44 import org.orekit.time.FieldAbsoluteDate;
45 import org.orekit.utils.FieldPVCoordinates;
46 import org.orekit.utils.ParameterDriver;
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75 public class HolmesFeatherstoneAttractionModel implements ForceModel, TideSystemProvider {
76
77
78
79
80
81 private static final int SCALING = 930;
82
83
84
85
86
87
88
89 private static final double MU_SCALE = FastMath.scalb(1.0, 32);
90
91
92 private final ParameterDriver gmParameterDriver;
93
94
95 private final NormalizedSphericalHarmonicsProvider provider;
96
97
98 private final Frame bodyFrame;
99
100
101 private final double[] gnmOj;
102
103
104 private final double[] hnmOj;
105
106
107 private final double[] enm;
108
109
110 private final double[] sectorial;
111
112
113
114
115
116
117 public HolmesFeatherstoneAttractionModel(final Frame centralBodyFrame,
118 final NormalizedSphericalHarmonicsProvider provider) {
119
120 gmParameterDriver = new ParameterDriver(NewtonianAttraction.CENTRAL_ATTRACTION_COEFFICIENT,
121 provider.getMu(), MU_SCALE, 0.0, Double.POSITIVE_INFINITY);
122
123 this.provider = provider;
124 this.bodyFrame = centralBodyFrame;
125
126
127
128 final int degree = provider.getMaxDegree();
129 final int size = FastMath.max(0, degree * (degree + 1) / 2 - 1);
130 gnmOj = new double[size];
131 hnmOj = new double[size];
132 enm = new double[size];
133
134
135
136
137
138 int index = 0;
139 for (int m = degree; m >= 0; --m) {
140 final int j = (m == 0) ? 2 : 1;
141 for (int n = FastMath.max(2, m + 1); n <= degree; ++n) {
142 final double f = (n - m) * (n + m + 1);
143 gnmOj[index] = 2 * (m + 1) / FastMath.sqrt(j * f);
144 hnmOj[index] = FastMath.sqrt((n + m + 2) * (n - m - 1) / (j * f));
145 enm[index] = FastMath.sqrt(f / j);
146 ++index;
147 }
148 }
149
150
151 sectorial = new double[degree + 1];
152 sectorial[0] = FastMath.scalb(1.0, -SCALING);
153 if (degree > 0) {
154 sectorial[1] = FastMath.sqrt(3) * sectorial[0];
155 }
156 for (int m = 2; m < sectorial.length; ++m) {
157 sectorial[m] = FastMath.sqrt((2 * m + 1) / (2.0 * m)) * sectorial[m - 1];
158 }
159
160 }
161
162
163 @Override
164 public boolean dependsOnPositionOnly() {
165 return true;
166 }
167
168
169 public TideSystem getTideSystem() {
170 return provider.getTideSystem();
171 }
172
173
174
175
176
177
178
179 public double getMu() {
180 return gmParameterDriver.getValue();
181 }
182
183
184
185
186
187 public double getMu(final AbsoluteDate date) {
188 return gmParameterDriver.getValue(date);
189 }
190
191
192
193
194
195
196
197 public double value(final AbsoluteDate date, final Vector3D position,
198 final double mu) {
199 return mu / position.getNorm() + nonCentralPart(date, position, mu);
200 }
201
202
203
204
205
206
207
208 public double nonCentralPart(final AbsoluteDate date, final Vector3D position, final double mu) {
209
210 final int degree = provider.getMaxDegree();
211 final int order = provider.getMaxOrder();
212 final NormalizedSphericalHarmonics harmonics = provider.onDate(date);
213
214
215 double[] pnm0Plus2 = new double[degree + 1];
216 double[] pnm0Plus1 = new double[degree + 1];
217 double[] pnm0 = new double[degree + 1];
218
219
220 final double x = position.getX();
221 final double y = position.getY();
222 final double z = position.getZ();
223 final double x2 = x * x;
224 final double y2 = y * y;
225 final double z2 = z * z;
226 final double rho2 = x2 + y2;
227 final double r2 = rho2 + z2;
228 final double r = FastMath.sqrt(r2);
229 final double rho = FastMath.sqrt(rho2);
230 final double t = z / r;
231 final double u = rho / r;
232 final double tOu = z / rho;
233
234
235 final double[] aOrN = createDistancePowersArray(provider.getAe() / r);
236
237
238 final double[][] cosSinLambda = createCosSinArrays(x / rho, y / rho);
239
240
241 int index = 0;
242 double value = 0;
243 for (int m = degree; m >= 0; --m) {
244
245
246 index = computeTesseral(m, degree, index, t, u, tOu,
247 pnm0Plus2, pnm0Plus1, null, pnm0, null, null);
248
249 if (m <= order) {
250
251
252
253 double sumDegreeS = 0;
254 double sumDegreeC = 0;
255 for (int n = FastMath.max(2, m); n <= degree; ++n) {
256 sumDegreeS += pnm0[n] * aOrN[n] * harmonics.getNormalizedSnm(n, m);
257 sumDegreeC += pnm0[n] * aOrN[n] * harmonics.getNormalizedCnm(n, m);
258 }
259
260
261 value = value * u + cosSinLambda[1][m] * sumDegreeS + cosSinLambda[0][m] * sumDegreeC;
262
263 }
264
265
266 final double[] tmp = pnm0Plus2;
267 pnm0Plus2 = pnm0Plus1;
268 pnm0Plus1 = pnm0;
269 pnm0 = tmp;
270
271 }
272
273
274 value = FastMath.scalb(value, SCALING);
275
276
277 return mu * value / r;
278
279 }
280
281
282
283
284
285
286
287 public double[] gradient(final AbsoluteDate date, final Vector3D position, final double mu) {
288
289 final int degree = provider.getMaxDegree();
290 final int order = provider.getMaxOrder();
291 final NormalizedSphericalHarmonics harmonics = provider.onDate(date);
292
293
294 double[] pnm0Plus2 = new double[degree + 1];
295 double[] pnm0Plus1 = new double[degree + 1];
296 double[] pnm0 = new double[degree + 1];
297 final double[] pnm1 = new double[degree + 1];
298
299
300 final double x = position.getX();
301 final double y = position.getY();
302 final double z = position.getZ();
303 final double x2 = x * x;
304 final double y2 = y * y;
305 final double z2 = z * z;
306 final double r2 = x2 + y2 + z2;
307 final double r = FastMath.sqrt (r2);
308 final double rho2 = x2 + y2;
309 final double rho = FastMath.sqrt(rho2);
310 final double t = z / r;
311 final double u = rho / r;
312 final double tOu = z / rho;
313
314
315 final double[] aOrN = createDistancePowersArray(provider.getAe() / r);
316
317
318 final double[][] cosSinLambda = createCosSinArrays(x / rho, y / rho);
319
320
321 int index = 0;
322 double value = 0;
323 final double[] gradient = new double[3];
324 for (int m = degree; m >= 0; --m) {
325
326
327 index = computeTesseral(m, degree, index, t, u, tOu,
328 pnm0Plus2, pnm0Plus1, null, pnm0, pnm1, null);
329
330 if (m <= order) {
331
332
333
334 double sumDegreeS = 0;
335 double sumDegreeC = 0;
336 double dSumDegreeSdR = 0;
337 double dSumDegreeCdR = 0;
338 double dSumDegreeSdTheta = 0;
339 double dSumDegreeCdTheta = 0;
340 for (int n = FastMath.max(2, m); n <= degree; ++n) {
341 final double qSnm = aOrN[n] * harmonics.getNormalizedSnm(n, m);
342 final double qCnm = aOrN[n] * harmonics.getNormalizedCnm(n, m);
343 final double nOr = n / r;
344 final double s0 = pnm0[n] * qSnm;
345 final double c0 = pnm0[n] * qCnm;
346 final double s1 = pnm1[n] * qSnm;
347 final double c1 = pnm1[n] * qCnm;
348 sumDegreeS += s0;
349 sumDegreeC += c0;
350 dSumDegreeSdR -= nOr * s0;
351 dSumDegreeCdR -= nOr * c0;
352 dSumDegreeSdTheta += s1;
353 dSumDegreeCdTheta += c1;
354 }
355
356
357
358
359
360 final double sML = cosSinLambda[1][m];
361 final double cML = cosSinLambda[0][m];
362 value = value * u + sML * sumDegreeS + cML * sumDegreeC;
363 gradient[0] = gradient[0] * u + sML * dSumDegreeSdR + cML * dSumDegreeCdR;
364 gradient[1] = gradient[1] * u + m * (cML * sumDegreeS - sML * sumDegreeC);
365 gradient[2] = gradient[2] * u + sML * dSumDegreeSdTheta + cML * dSumDegreeCdTheta;
366
367 }
368
369
370 final double[] tmp = pnm0Plus2;
371 pnm0Plus2 = pnm0Plus1;
372 pnm0Plus1 = pnm0;
373 pnm0 = tmp;
374
375 }
376
377
378 value = FastMath.scalb(value, SCALING);
379 gradient[0] = FastMath.scalb(gradient[0], SCALING);
380 gradient[1] = FastMath.scalb(gradient[1], SCALING);
381 gradient[2] = FastMath.scalb(gradient[2], SCALING);
382
383
384 final double muOr = mu / r;
385 value *= muOr;
386 gradient[0] = muOr * gradient[0] - value / r;
387 gradient[1] *= muOr;
388 gradient[2] *= muOr;
389
390
391 return new SphericalCoordinates(position).toCartesianGradient(gradient);
392
393 }
394
395
396
397
398
399
400
401
402 public <T extends CalculusFieldElement<T>> T[] gradient(final FieldAbsoluteDate<T> date, final FieldVector3D<T> position,
403 final T mu) {
404
405 final int degree = provider.getMaxDegree();
406 final int order = provider.getMaxOrder();
407 final NormalizedSphericalHarmonics harmonics = provider.onDate(date.toAbsoluteDate());
408 final T zero = date.getField().getZero();
409
410 T[] pnm0Plus2 = MathArrays.buildArray(date.getField(), degree + 1);
411 T[] pnm0Plus1 = MathArrays.buildArray(date.getField(), degree + 1);
412 T[] pnm0 = MathArrays.buildArray(date.getField(), degree + 1);
413 final T[] pnm1 = MathArrays.buildArray(date.getField(), degree + 1);
414
415
416 final T x = position.getX();
417 final T y = position.getY();
418 final T z = position.getZ();
419 final T x2 = x.square();
420 final T y2 = y.square();
421 final T rho2 = x2.add(y2);
422 final T rho = rho2.sqrt();
423 final T z2 = z.square();
424 final T r2 = rho2.add(z2);
425 final T r = r2.sqrt();
426 final T t = z.divide(r);
427 final T u = rho.divide(r);
428 final T tOu = z.divide(rho);
429
430
431 final T[] aOrN = createDistancePowersArray(r.reciprocal().multiply(provider.getAe()));
432
433
434 final T[][] cosSinLambda = createCosSinArrays(x.divide(rho), y.divide(rho));
435
436 int index = 0;
437 T value = zero;
438 final T[] gradient = MathArrays.buildArray(zero.getField(), 3);
439 for (int m = degree; m >= 0; --m) {
440
441
442 index = computeTesseral(m, degree, index, t, u, tOu,
443 pnm0Plus2, pnm0Plus1, null, pnm0, pnm1, null);
444 if (m <= order) {
445
446
447
448 T sumDegreeS = zero;
449 T sumDegreeC = zero;
450 T dSumDegreeSdR = zero;
451 T dSumDegreeCdR = zero;
452 T dSumDegreeSdTheta = zero;
453 T dSumDegreeCdTheta = zero;
454 for (int n = FastMath.max(2, m); n <= degree; ++n) {
455 final T qSnm = aOrN[n].multiply(harmonics.getNormalizedSnm(n, m));
456 final T qCnm = aOrN[n].multiply(harmonics.getNormalizedCnm(n, m));
457 final T nOr = r.reciprocal().multiply(n);
458 final T s0 = pnm0[n].multiply(qSnm);
459 final T c0 = pnm0[n].multiply(qCnm);
460 final T s1 = pnm1[n].multiply(qSnm);
461 final T c1 = pnm1[n].multiply(qCnm);
462 sumDegreeS = sumDegreeS .add(s0);
463 sumDegreeC = sumDegreeC .add(c0);
464 dSumDegreeSdR = dSumDegreeSdR .subtract(nOr.multiply(s0));
465 dSumDegreeCdR = dSumDegreeCdR .subtract(nOr.multiply(c0));
466 dSumDegreeSdTheta = dSumDegreeSdTheta.add(s1);
467 dSumDegreeCdTheta = dSumDegreeCdTheta.add(c1);
468 }
469
470
471
472
473
474 final T sML = cosSinLambda[1][m];
475 final T cML = cosSinLambda[0][m];
476 value = value .multiply(u).add(sML.multiply(sumDegreeS )).add(cML.multiply(sumDegreeC));
477 gradient[0] = gradient[0].multiply(u).add(sML.multiply(dSumDegreeSdR)).add(cML.multiply(dSumDegreeCdR));
478 gradient[1] = gradient[1].multiply(u).add(cML.multiply(sumDegreeS).subtract(sML.multiply(sumDegreeC)).multiply(m));
479 gradient[2] = gradient[2].multiply(u).add(sML.multiply(dSumDegreeSdTheta)).add(cML.multiply(dSumDegreeCdTheta));
480 }
481
482 final T[] tmp = pnm0Plus2;
483 pnm0Plus2 = pnm0Plus1;
484 pnm0Plus1 = pnm0;
485 pnm0 = tmp;
486
487 }
488
489 value = value.scalb(SCALING);
490 gradient[0] = gradient[0].scalb(SCALING);
491 gradient[1] = gradient[1].scalb(SCALING);
492 gradient[2] = gradient[2].scalb(SCALING);
493
494
495 final T muOr = r.reciprocal().multiply(mu);
496 value = value.multiply(muOr);
497 gradient[0] = muOr.multiply(gradient[0]).subtract(value.divide(r));
498 gradient[1] = gradient[1].multiply(muOr);
499 gradient[2] = gradient[2].multiply(muOr);
500
501
502
503
504
505
506 final T xPos = position.getX();
507 final T yPos = position.getY();
508 final T zPos = position.getZ();
509 final T rho2Pos = x.square().add(y.square());
510 final T rhoPos = rho2.sqrt();
511 final T r2Pos = rho2.add(z.square());
512 final T rPos = r2Pos.sqrt();
513
514 final T[][] jacobianPos = MathArrays.buildArray(zero.getField(), 3, 3);
515
516
517 jacobianPos[0][0] = xPos.divide(rPos);
518 jacobianPos[0][1] = yPos.divide(rPos);
519 jacobianPos[0][2] = zPos.divide(rPos);
520
521
522 jacobianPos[1][0] = yPos.negate().divide(rho2Pos);
523 jacobianPos[1][1] = xPos.divide(rho2Pos);
524
525
526
527 final T rhoPosTimesR2Pos = rhoPos.multiply(r2Pos);
528 jacobianPos[2][0] = xPos.multiply(zPos).divide(rhoPosTimesR2Pos);
529 jacobianPos[2][1] = yPos.multiply(zPos).divide(rhoPosTimesR2Pos);
530 jacobianPos[2][2] = rhoPos.negate().divide(r2Pos);
531 final T[] cartGradPos = MathArrays.buildArray(zero.getField(), 3);
532 cartGradPos[0] = gradient[0].multiply(jacobianPos[0][0]).add(gradient[1].multiply(jacobianPos[1][0])).add(gradient[2].multiply(jacobianPos[2][0]));
533 cartGradPos[1] = gradient[0].multiply(jacobianPos[0][1]).add(gradient[1].multiply(jacobianPos[1][1])).add(gradient[2].multiply(jacobianPos[2][1]));
534 cartGradPos[2] = gradient[0].multiply(jacobianPos[0][2]) .add(gradient[2].multiply(jacobianPos[2][2]));
535 return cartGradPos;
536
537 }
538
539
540
541
542
543
544
545 private GradientHessian gradientHessian(final AbsoluteDate date, final Vector3D position, final double mu) {
546
547 final int degree = provider.getMaxDegree();
548 final int order = provider.getMaxOrder();
549 final NormalizedSphericalHarmonics harmonics = provider.onDate(date);
550
551
552 double[] pnm0Plus2 = new double[degree + 1];
553 double[] pnm0Plus1 = new double[degree + 1];
554 double[] pnm0 = new double[degree + 1];
555 double[] pnm1Plus1 = new double[degree + 1];
556 double[] pnm1 = new double[degree + 1];
557 final double[] pnm2 = new double[degree + 1];
558
559
560 final double x = position.getX();
561 final double y = position.getY();
562 final double z = position.getZ();
563 final double x2 = x * x;
564 final double y2 = y * y;
565 final double z2 = z * z;
566 final double rho2 = x2 + y2;
567 final double rho = FastMath.sqrt(rho2);
568 final double r2 = rho2 + z2;
569 final double r = FastMath.sqrt(r2);
570 final double t = z / r;
571 final double u = rho / r;
572 final double tOu = z / rho;
573
574
575 final double[] aOrN = createDistancePowersArray(provider.getAe() / r);
576
577
578 final double[][] cosSinLambda = createCosSinArrays(x / rho, y / rho);
579
580
581 int index = 0;
582 double value = 0;
583 final double[] gradient = new double[3];
584 final double[][] hessian = new double[3][3];
585 for (int m = degree; m >= 0; --m) {
586
587
588 index = computeTesseral(m, degree, index, t, u, tOu,
589 pnm0Plus2, pnm0Plus1, pnm1Plus1, pnm0, pnm1, pnm2);
590
591 if (m <= order) {
592
593
594
595 double sumDegreeS = 0;
596 double sumDegreeC = 0;
597 double dSumDegreeSdR = 0;
598 double dSumDegreeCdR = 0;
599 double dSumDegreeSdTheta = 0;
600 double dSumDegreeCdTheta = 0;
601 double d2SumDegreeSdRdR = 0;
602 double d2SumDegreeSdRdTheta = 0;
603 double d2SumDegreeSdThetadTheta = 0;
604 double d2SumDegreeCdRdR = 0;
605 double d2SumDegreeCdRdTheta = 0;
606 double d2SumDegreeCdThetadTheta = 0;
607 for (int n = FastMath.max(2, m); n <= degree; ++n) {
608 final double qSnm = aOrN[n] * harmonics.getNormalizedSnm(n, m);
609 final double qCnm = aOrN[n] * harmonics.getNormalizedCnm(n, m);
610 final double nOr = n / r;
611 final double nnP1Or2 = nOr * (n + 1) / r;
612 final double s0 = pnm0[n] * qSnm;
613 final double c0 = pnm0[n] * qCnm;
614 final double s1 = pnm1[n] * qSnm;
615 final double c1 = pnm1[n] * qCnm;
616 final double s2 = pnm2[n] * qSnm;
617 final double c2 = pnm2[n] * qCnm;
618 sumDegreeS += s0;
619 sumDegreeC += c0;
620 dSumDegreeSdR -= nOr * s0;
621 dSumDegreeCdR -= nOr * c0;
622 dSumDegreeSdTheta += s1;
623 dSumDegreeCdTheta += c1;
624 d2SumDegreeSdRdR += nnP1Or2 * s0;
625 d2SumDegreeSdRdTheta -= nOr * s1;
626 d2SumDegreeSdThetadTheta += s2;
627 d2SumDegreeCdRdR += nnP1Or2 * c0;
628 d2SumDegreeCdRdTheta -= nOr * c1;
629 d2SumDegreeCdThetadTheta += c2;
630 }
631
632
633 final double sML = cosSinLambda[1][m];
634 final double cML = cosSinLambda[0][m];
635 value = value * u + sML * sumDegreeS + cML * sumDegreeC;
636 gradient[0] = gradient[0] * u + sML * dSumDegreeSdR + cML * dSumDegreeCdR;
637 gradient[1] = gradient[1] * u + m * (cML * sumDegreeS - sML * sumDegreeC);
638 gradient[2] = gradient[2] * u + sML * dSumDegreeSdTheta + cML * dSumDegreeCdTheta;
639 hessian[0][0] = hessian[0][0] * u + sML * d2SumDegreeSdRdR + cML * d2SumDegreeCdRdR;
640 hessian[1][0] = hessian[1][0] * u + m * (cML * dSumDegreeSdR - sML * dSumDegreeCdR);
641 hessian[2][0] = hessian[2][0] * u + sML * d2SumDegreeSdRdTheta + cML * d2SumDegreeCdRdTheta;
642 hessian[1][1] = hessian[1][1] * u - m * m * (sML * sumDegreeS + cML * sumDegreeC);
643 hessian[2][1] = hessian[2][1] * u + m * (cML * dSumDegreeSdTheta - sML * dSumDegreeCdTheta);
644 hessian[2][2] = hessian[2][2] * u + sML * d2SumDegreeSdThetadTheta + cML * d2SumDegreeCdThetadTheta;
645
646 }
647
648
649 final double[] tmp0 = pnm0Plus2;
650 pnm0Plus2 = pnm0Plus1;
651 pnm0Plus1 = pnm0;
652 pnm0 = tmp0;
653 final double[] tmp1 = pnm1Plus1;
654 pnm1Plus1 = pnm1;
655 pnm1 = tmp1;
656
657 }
658
659
660 value = FastMath.scalb(value, SCALING);
661 for (int i = 0; i < 3; ++i) {
662 gradient[i] = FastMath.scalb(gradient[i], SCALING);
663 for (int j = 0; j <= i; ++j) {
664 hessian[i][j] = FastMath.scalb(hessian[i][j], SCALING);
665 }
666 }
667
668
669
670 final double muOr = mu / r;
671 value *= muOr;
672 gradient[0] = muOr * gradient[0] - value / r;
673 gradient[1] *= muOr;
674 gradient[2] *= muOr;
675 hessian[0][0] = muOr * hessian[0][0] - 2 * gradient[0] / r;
676 hessian[1][0] = muOr * hessian[1][0] - gradient[1] / r;
677 hessian[2][0] = muOr * hessian[2][0] - gradient[2] / r;
678 hessian[1][1] *= muOr;
679 hessian[2][1] *= muOr;
680 hessian[2][2] *= muOr;
681
682
683 final SphericalCoordinates sc = new SphericalCoordinates(position);
684 return new GradientHessian(sc.toCartesianGradient(gradient),
685 sc.toCartesianHessian(hessian, gradient));
686
687
688 }
689
690
691 private static class GradientHessian {
692
693
694 private final double[] gradient;
695
696
697 private final double[][] hessian;
698
699
700
701
702
703
704
705
706 GradientHessian(final double[] gradient, final double[][] hessian) {
707 this.gradient = gradient;
708 this.hessian = hessian;
709 }
710
711
712
713
714 public double[] getGradient() {
715 return gradient;
716 }
717
718
719
720
721 public double[][] getHessian() {
722 return hessian;
723 }
724
725 }
726
727
728
729
730
731 private double[] createDistancePowersArray(final double aOr) {
732
733
734 final double[] aOrN = new double[provider.getMaxDegree() + 1];
735 aOrN[0] = 1;
736 if (provider.getMaxDegree() > 0) {
737 aOrN[1] = aOr;
738 }
739
740
741 for (int n = 2; n < aOrN.length; ++n) {
742 final int p = n / 2;
743 final int q = n - p;
744 aOrN[n] = aOrN[p] * aOrN[q];
745 }
746
747 return aOrN;
748
749 }
750
751
752
753
754
755 private <T extends CalculusFieldElement<T>> T[] createDistancePowersArray(final T aOr) {
756
757
758 final T[] aOrN = MathArrays.buildArray(aOr.getField(), provider.getMaxDegree() + 1);
759 aOrN[0] = aOr.getField().getOne();
760 if (provider.getMaxDegree() > 0) {
761 aOrN[1] = aOr;
762 }
763
764
765 for (int n = 2; n < aOrN.length; ++n) {
766 final int p = n / 2;
767 final int q = n - p;
768 aOrN[n] = aOrN[p].multiply(aOrN[q]);
769 }
770
771 return aOrN;
772
773 }
774
775
776
777
778
779
780
781 private double[][] createCosSinArrays(final double cosLambda, final double sinLambda) {
782
783
784 final double[][] cosSin = new double[2][provider.getMaxOrder() + 1];
785 cosSin[0][0] = 1;
786 cosSin[1][0] = 0;
787 if (provider.getMaxOrder() > 0) {
788 cosSin[0][1] = cosLambda;
789 cosSin[1][1] = sinLambda;
790
791
792 for (int m = 2; m < cosSin[0].length; ++m) {
793
794
795
796
797
798
799 final int p = m / 2;
800 final int q = m - p;
801
802 cosSin[0][m] = cosSin[0][p] * cosSin[0][q] - cosSin[1][p] * cosSin[1][q];
803 cosSin[1][m] = cosSin[1][p] * cosSin[0][q] + cosSin[0][p] * cosSin[1][q];
804 }
805 }
806
807 return cosSin;
808
809 }
810
811
812
813
814
815
816
817
818 private <T extends CalculusFieldElement<T>> T[][] createCosSinArrays(final T cosLambda, final T sinLambda) {
819
820 final T one = cosLambda.getField().getOne();
821 final T zero = cosLambda.getField().getZero();
822
823 final T[][] cosSin = MathArrays.buildArray(one.getField(), 2, provider.getMaxOrder() + 1);
824 cosSin[0][0] = one;
825 cosSin[1][0] = zero;
826 if (provider.getMaxOrder() > 0) {
827 cosSin[0][1] = cosLambda;
828 cosSin[1][1] = sinLambda;
829
830
831 for (int m = 2; m < cosSin[0].length; ++m) {
832
833
834
835
836
837
838 final int p = m / 2;
839 final int q = m - p;
840
841 cosSin[0][m] = cosSin[0][p].multiply(cosSin[0][q]).subtract(cosSin[1][p].multiply(cosSin[1][q]));
842 cosSin[1][m] = cosSin[1][p].multiply(cosSin[0][q]).add(cosSin[0][p].multiply(cosSin[1][q]));
843
844 }
845 }
846
847 return cosSin;
848
849 }
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872 private int computeTesseral(final int m, final int degree, final int index,
873 final double t, final double u, final double tOu,
874 final double[] pnm0Plus2, final double[] pnm0Plus1, final double[] pnm1Plus1,
875 final double[] pnm0, final double[] pnm1, final double[] pnm2) {
876
877 final double u2 = u * u;
878
879
880 int n = FastMath.max(2, m);
881 if (n == m) {
882 pnm0[n] = sectorial[n];
883 ++n;
884 }
885
886
887 int localIndex = index;
888 while (n <= degree) {
889
890
891 pnm0[n] = gnmOj[localIndex] * t * pnm0Plus1[n] - hnmOj[localIndex] * u2 * pnm0Plus2[n];
892
893 ++localIndex;
894 ++n;
895
896 }
897
898 if (pnm1 != null) {
899
900
901 n = FastMath.max(2, m);
902 if (n == m) {
903 pnm1[n] = m * tOu * pnm0[n];
904 ++n;
905 }
906
907
908 localIndex = index;
909 while (n <= degree) {
910
911
912 pnm1[n] = m * tOu * pnm0[n] - enm[localIndex] * u * pnm0Plus1[n];
913
914 ++localIndex;
915 ++n;
916
917 }
918
919 if (pnm2 != null) {
920
921
922 n = FastMath.max(2, m);
923 if (n == m) {
924 pnm2[n] = m * (tOu * pnm1[n] - pnm0[n] / u2);
925 ++n;
926 }
927
928
929 localIndex = index;
930 while (n <= degree) {
931
932
933 pnm2[n] = m * (tOu * pnm1[n] - pnm0[n] / u2) - enm[localIndex] * u * pnm1Plus1[n];
934
935 ++localIndex;
936 ++n;
937
938 }
939
940 }
941
942 }
943
944 return localIndex;
945
946 }
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970 private <T extends CalculusFieldElement<T>> int computeTesseral(final int m, final int degree, final int index,
971 final T t, final T u, final T tOu,
972 final T[] pnm0Plus2, final T[] pnm0Plus1, final T[] pnm1Plus1,
973 final T[] pnm0, final T[] pnm1, final T[] pnm2) {
974
975 final T u2 = u.square();
976 final T zero = u.getField().getZero();
977
978 int n = FastMath.max(2, m);
979 if (n == m) {
980 pnm0[n] = zero.newInstance(sectorial[n]);
981 ++n;
982 }
983
984
985 int localIndex = index;
986 while (n <= degree) {
987
988
989 pnm0[n] = t.multiply(gnmOj[localIndex]).multiply(pnm0Plus1[n]).subtract(u2.multiply(pnm0Plus2[n]).multiply(hnmOj[localIndex]));
990 ++localIndex;
991 ++n;
992
993 }
994 if (pnm1 != null) {
995
996
997 n = FastMath.max(2, m);
998 if (n == m) {
999 pnm1[n] = tOu.multiply(m).multiply(pnm0[n]);
1000 ++n;
1001 }
1002
1003
1004 localIndex = index;
1005 while (n <= degree) {
1006
1007
1008 pnm1[n] = tOu.multiply(m).multiply(pnm0[n]).subtract(u.multiply(enm[localIndex]).multiply(pnm0Plus1[n]));
1009
1010 ++localIndex;
1011 ++n;
1012
1013 }
1014
1015 if (pnm2 != null) {
1016
1017
1018 n = FastMath.max(2, m);
1019 if (n == m) {
1020 pnm2[n] = tOu.multiply(pnm1[n]).subtract(pnm0[n].divide(u2)).multiply(m);
1021 ++n;
1022 }
1023
1024
1025 localIndex = index;
1026 while (n <= degree) {
1027
1028
1029 pnm2[n] = tOu.multiply(pnm1[n]).subtract(pnm0[n].divide(u2)).multiply(m).subtract(u.multiply(pnm1Plus1[n]).multiply(enm[localIndex]));
1030 ++localIndex;
1031 ++n;
1032
1033 }
1034
1035 }
1036
1037 }
1038 return localIndex;
1039
1040 }
1041
1042
1043 @Override
1044 public Vector3D acceleration(final SpacecraftState s, final double[] parameters) {
1045
1046 final double mu = parameters[0];
1047
1048
1049 final AbsoluteDate date = s.getDate();
1050 final StaticTransform fromBodyFrame = bodyFrame.getStaticTransformTo(s.getFrame(), date);
1051 final StaticTransform toBodyFrame = fromBodyFrame.getInverse();
1052 final Vector3D position = toBodyFrame.transformPosition(s.getPosition());
1053
1054
1055 return fromBodyFrame.transformVector(new Vector3D(gradient(date, position, mu)));
1056
1057 }
1058
1059
1060 public <T extends CalculusFieldElement<T>> FieldVector3D<T> acceleration(final FieldSpacecraftState<T> s,
1061 final T[] parameters) {
1062
1063 final T mu = parameters[0];
1064
1065
1066 if (isGradientStateDerivative(s)) {
1067 @SuppressWarnings("unchecked")
1068 final FieldVector3D<Gradient> p = (FieldVector3D<Gradient>) s.getPosition();
1069 @SuppressWarnings("unchecked")
1070 final FieldVector3D<T> a = (FieldVector3D<T>) accelerationWrtState(s.getDate().toAbsoluteDate(),
1071 s.getFrame(), p,
1072 (Gradient) mu);
1073 return a;
1074 } else if (isDSStateDerivative(s)) {
1075 @SuppressWarnings("unchecked")
1076 final FieldVector3D<DerivativeStructure> p = (FieldVector3D<DerivativeStructure>) s.getPosition();
1077 @SuppressWarnings("unchecked")
1078 final FieldVector3D<T> a = (FieldVector3D<T>) accelerationWrtState(s.getDate().toAbsoluteDate(),
1079 s.getFrame(), p,
1080 (DerivativeStructure) mu);
1081 return a;
1082 }
1083
1084
1085 final FieldAbsoluteDate<T> date = s.getDate();
1086 final FieldStaticTransform<T> fromBodyFrame = bodyFrame.getStaticTransformTo(s.getFrame(), date);
1087 final FieldStaticTransform<T> toBodyFrame = fromBodyFrame.getInverse();
1088 final FieldVector3D<T> position = toBodyFrame.transformPosition(s.getPosition());
1089
1090
1091 return fromBodyFrame.transformVector(new FieldVector3D<>(gradient(date, position, mu)));
1092
1093 }
1094
1095
1096
1097
1098
1099
1100
1101 private <T extends CalculusFieldElement<T>> boolean isDSStateDerivative(final FieldSpacecraftState<T> state) {
1102 try {
1103 final DerivativeStructure dsMass = (DerivativeStructure) state.getMass();
1104 final int o = dsMass.getOrder();
1105 final int p = dsMass.getFreeParameters();
1106 if (o != 1 || p < 3) {
1107 return false;
1108 }
1109 @SuppressWarnings("unchecked")
1110 final FieldPVCoordinates<DerivativeStructure> pv = (FieldPVCoordinates<DerivativeStructure>) state.getPVCoordinates();
1111 return isVariable(pv.getPosition().getX(), 0) &&
1112 isVariable(pv.getPosition().getY(), 1) &&
1113 isVariable(pv.getPosition().getZ(), 2);
1114 } catch (ClassCastException cce) {
1115 return false;
1116 }
1117 }
1118
1119
1120
1121
1122
1123
1124
1125 private <T extends CalculusFieldElement<T>> boolean isGradientStateDerivative(final FieldSpacecraftState<T> state) {
1126 try {
1127 final Gradient gMass = (Gradient) state.getMass();
1128 final int p = gMass.getFreeParameters();
1129 if (p < 3) {
1130 return false;
1131 }
1132 @SuppressWarnings("unchecked")
1133 final FieldPVCoordinates<Gradient> pv = (FieldPVCoordinates<Gradient>) state.getPVCoordinates();
1134 return isVariable(pv.getPosition().getX(), 0) &&
1135 isVariable(pv.getPosition().getY(), 1) &&
1136 isVariable(pv.getPosition().getZ(), 2);
1137 } catch (ClassCastException cce) {
1138 return false;
1139 }
1140 }
1141
1142
1143
1144
1145
1146
1147
1148 private boolean isVariable(final DerivativeStructure ds, final int index) {
1149 final double[] derivatives = ds.getAllDerivatives();
1150 boolean check = true;
1151 for (int i = 1; i < derivatives.length; ++i) {
1152 check &= derivatives[i] == ((index + 1 == i) ? 1.0 : 0.0);
1153 }
1154 return check;
1155 }
1156
1157
1158
1159
1160
1161
1162
1163 private boolean isVariable(final Gradient g, final int index) {
1164 final double[] derivatives = g.getGradient();
1165 boolean check = true;
1166 for (int i = 0; i < derivatives.length; ++i) {
1167 check &= derivatives[i] == ((index == i) ? 1.0 : 0.0);
1168 }
1169 return check;
1170 }
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199 private FieldVector3D<DerivativeStructure> accelerationWrtState(final AbsoluteDate date, final Frame frame,
1200 final FieldVector3D<DerivativeStructure> position,
1201 final DerivativeStructure mu) {
1202
1203
1204 final int freeParameters = mu.getFreeParameters();
1205
1206
1207 final StaticTransform fromBodyFrame = bodyFrame.getStaticTransformTo(frame, date);
1208 final StaticTransform toBodyFrame = fromBodyFrame.getInverse();
1209 final Vector3D positionBody = toBodyFrame.transformPosition(position.toVector3D());
1210
1211
1212 final GradientHessian gh = gradientHessian(date, positionBody, mu.getReal());
1213
1214
1215 final double[] gInertial = fromBodyFrame.transformVector(new Vector3D(gh.getGradient())).toArray();
1216
1217
1218 final RealMatrix hBody = new Array2DRowRealMatrix(gh.getHessian(), false);
1219 final RealMatrix rot = new Array2DRowRealMatrix(toBodyFrame.getRotation().getMatrix());
1220 final RealMatrix hInertial = rot.transposeMultiply(hBody).multiply(rot);
1221
1222
1223 final double[] derivatives = new double[freeParameters + 1];
1224 final DerivativeStructure[] accDer = new DerivativeStructure[3];
1225 for (int i = 0; i < 3; ++i) {
1226
1227
1228 derivatives[0] = gInertial[i];
1229
1230
1231 derivatives[1] = hInertial.getEntry(i, 0);
1232 derivatives[2] = hInertial.getEntry(i, 1);
1233 derivatives[3] = hInertial.getEntry(i, 2);
1234
1235
1236 if (derivatives.length > 4 && isVariable(mu, 3)) {
1237 derivatives[4] = gInertial[i] / mu.getReal();
1238 }
1239
1240 accDer[i] = position.getX().getFactory().build(derivatives);
1241
1242 }
1243
1244 return new FieldVector3D<>(accDer);
1245
1246 }
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275 private FieldVector3D<Gradient> accelerationWrtState(final AbsoluteDate date, final Frame frame,
1276 final FieldVector3D<Gradient> position,
1277 final Gradient mu) {
1278
1279
1280 final int freeParameters = mu.getFreeParameters();
1281
1282
1283 final StaticTransform fromBodyFrame = bodyFrame.getStaticTransformTo(frame, date);
1284 final StaticTransform toBodyFrame = fromBodyFrame.getInverse();
1285 final Vector3D positionBody = toBodyFrame.transformPosition(position.toVector3D());
1286
1287
1288 final GradientHessian gh = gradientHessian(date, positionBody, mu.getReal());
1289
1290
1291 final double[] gInertial = fromBodyFrame.transformVector(new Vector3D(gh.getGradient())).toArray();
1292
1293
1294 final RealMatrix hBody = new Array2DRowRealMatrix(gh.getHessian(), false);
1295 final RealMatrix rot = new Array2DRowRealMatrix(toBodyFrame.getRotation().getMatrix());
1296 final RealMatrix hInertial = rot.transposeMultiply(hBody).multiply(rot);
1297
1298
1299 final double[] derivatives = new double[freeParameters];
1300 final Gradient[] accDer = new Gradient[3];
1301 for (int i = 0; i < 3; ++i) {
1302
1303
1304 derivatives[0] = hInertial.getEntry(i, 0);
1305 derivatives[1] = hInertial.getEntry(i, 1);
1306 derivatives[2] = hInertial.getEntry(i, 2);
1307
1308
1309 if (derivatives.length > 3 && isVariable(mu, 3)) {
1310 derivatives[3] = gInertial[i] / mu.getReal();
1311 }
1312
1313 accDer[i] = new Gradient(gInertial[i], derivatives);
1314
1315 }
1316
1317 return new FieldVector3D<>(accDer);
1318
1319 }
1320
1321
1322 public List<ParameterDriver> getParametersDrivers() {
1323 return Collections.singletonList(gmParameterDriver);
1324 }
1325
1326 }