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