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.io.Serializable;
21
22 import org.hipparchus.analysis.differentiation.DerivativeStructure;
23 import org.hipparchus.geometry.euclidean.threed.FieldRotation;
24 import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
25 import org.hipparchus.geometry.euclidean.threed.SphericalCoordinates;
26 import org.hipparchus.geometry.euclidean.threed.Vector3D;
27 import org.hipparchus.linear.Array2DRowRealMatrix;
28 import org.hipparchus.linear.RealMatrix;
29 import org.hipparchus.util.FastMath;
30 import org.orekit.errors.OrekitException;
31 import org.orekit.errors.OrekitInternalError;
32 import org.orekit.forces.AbstractForceModel;
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.Frame;
38 import org.orekit.frames.Transform;
39 import org.orekit.propagation.SpacecraftState;
40 import org.orekit.propagation.events.EventDetector;
41 import org.orekit.propagation.numerical.TimeDerivativesEquations;
42 import org.orekit.time.AbsoluteDate;
43 import org.orekit.utils.ParameterDriver;
44 import org.orekit.utils.ParameterObserver;
45
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 extends AbstractForceModel implements 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[] parametersDrivers;
92
93
94 private final NormalizedSphericalHarmonicsProvider provider;
95
96
97 private double mu;
98
99
100 private final Frame bodyFrame;
101
102
103 private final double[] gnmOj;
104
105
106 private final double[] hnmOj;
107
108
109 private final double[] enm;
110
111
112 private final double[] sectorial;
113
114
115
116
117
118
119 public HolmesFeatherstoneAttractionModel(final Frame centralBodyFrame,
120 final NormalizedSphericalHarmonicsProvider provider) {
121
122 this.parametersDrivers = new ParameterDriver[1];
123 try {
124 parametersDrivers[0] = new ParameterDriver(NewtonianAttraction.CENTRAL_ATTRACTION_COEFFICIENT,
125 provider.getMu(), MU_SCALE, 0.0, Double.POSITIVE_INFINITY);
126 parametersDrivers[0].addObserver(new ParameterObserver() {
127
128 @Override
129 public void valueChanged(final double previousValue, final ParameterDriver driver) {
130 HolmesFeatherstoneAttractionModel.this.mu = driver.getValue();
131 }
132 });
133 } catch (OrekitException oe) {
134
135 throw new OrekitInternalError(oe);
136 };
137
138 this.provider = provider;
139 this.mu = provider.getMu();
140 this.bodyFrame = centralBodyFrame;
141
142
143
144 final int degree = provider.getMaxDegree();
145 final int size = FastMath.max(0, degree * (degree + 1) / 2 - 1);
146 gnmOj = new double[size];
147 hnmOj = new double[size];
148 enm = new double[size];
149
150
151
152
153
154 int index = 0;
155 for (int m = degree; m >= 0; --m) {
156 final int j = (m == 0) ? 2 : 1;
157 for (int n = FastMath.max(2, m + 1); n <= degree; ++n) {
158 final double f = (n - m) * (n + m + 1);
159 gnmOj[index] = 2 * (m + 1) / FastMath.sqrt(j * f);
160 hnmOj[index] = FastMath.sqrt((n + m + 2) * (n - m - 1) / (j * f));
161 enm[index] = FastMath.sqrt(f / j);
162 ++index;
163 }
164 }
165
166
167 sectorial = new double[degree + 1];
168 sectorial[0] = FastMath.scalb(1.0, -SCALING);
169 sectorial[1] = FastMath.sqrt(3) * sectorial[0];
170 for (int m = 2; m < sectorial.length; ++m) {
171 sectorial[m] = FastMath.sqrt((2 * m + 1) / (2.0 * m)) * sectorial[m - 1];
172 }
173
174 }
175
176
177 public TideSystem getTideSystem() {
178 return provider.getTideSystem();
179 }
180
181
182
183
184
185
186
187 public double value(final AbsoluteDate date, final Vector3D position)
188 throws OrekitException {
189 return mu / position.getNorm() + nonCentralPart(date, position);
190 }
191
192
193
194
195
196
197
198 public double nonCentralPart(final AbsoluteDate date, final Vector3D position)
199 throws OrekitException {
200
201 final int degree = provider.getMaxDegree();
202 final int order = provider.getMaxOrder();
203 final NormalizedSphericalHarmonics harmonics = provider.onDate(date);
204
205
206 double[] pnm0Plus2 = new double[degree + 1];
207 double[] pnm0Plus1 = new double[degree + 1];
208 double[] pnm0 = new double[degree + 1];
209
210
211 final double x = position.getX();
212 final double y = position.getY();
213 final double z = position.getZ();
214 final double x2 = x * x;
215 final double y2 = y * y;
216 final double z2 = z * z;
217 final double r2 = x2 + y2 + z2;
218 final double r = FastMath.sqrt (r2);
219 final double rho = FastMath.sqrt(x2 + y2);
220 final double t = z / r;
221 final double u = rho / r;
222 final double tOu = z / rho;
223
224
225 final double[] aOrN = createDistancePowersArray(provider.getAe() / r);
226
227
228 final double[][] cosSinLambda = createCosSinArrays(position.getX() / rho, position.getY() / rho);
229
230
231 int index = 0;
232 double value = 0;
233 for (int m = degree; m >= 0; --m) {
234
235
236 index = computeTesseral(m, degree, index, t, u, tOu,
237 pnm0Plus2, pnm0Plus1, null, pnm0, null, null);
238
239 if (m <= order) {
240
241
242
243 double sumDegreeS = 0;
244 double sumDegreeC = 0;
245 for (int n = FastMath.max(2, m); n <= degree; ++n) {
246 sumDegreeS += pnm0[n] * aOrN[n] * harmonics.getNormalizedSnm(n, m);
247 sumDegreeC += pnm0[n] * aOrN[n] * harmonics.getNormalizedCnm(n, m);
248 }
249
250
251 value = value * u + cosSinLambda[1][m] * sumDegreeS + cosSinLambda[0][m] * sumDegreeC;
252
253 }
254
255
256 final double[] tmp = pnm0Plus2;
257 pnm0Plus2 = pnm0Plus1;
258 pnm0Plus1 = pnm0;
259 pnm0 = tmp;
260
261 }
262
263
264 value = FastMath.scalb(value, SCALING);
265
266
267 return mu * value / r;
268
269 }
270
271
272
273
274
275
276
277 public double[] gradient(final AbsoluteDate date, final Vector3D position)
278 throws OrekitException {
279
280 final int degree = provider.getMaxDegree();
281 final int order = provider.getMaxOrder();
282 final NormalizedSphericalHarmonics harmonics = provider.onDate(date);
283
284
285 double[] pnm0Plus2 = new double[degree + 1];
286 double[] pnm0Plus1 = new double[degree + 1];
287 double[] pnm0 = new double[degree + 1];
288 final double[] pnm1 = new double[degree + 1];
289
290
291 final double x = position.getX();
292 final double y = position.getY();
293 final double z = position.getZ();
294 final double x2 = x * x;
295 final double y2 = y * y;
296 final double z2 = z * z;
297 final double r2 = x2 + y2 + z2;
298 final double r = FastMath.sqrt (r2);
299 final double rho2 = x2 + y2;
300 final double rho = FastMath.sqrt(rho2);
301 final double t = z / r;
302 final double u = rho / r;
303 final double tOu = z / rho;
304
305
306 final double[] aOrN = createDistancePowersArray(provider.getAe() / r);
307
308
309 final double[][] cosSinLambda = createCosSinArrays(position.getX() / rho, position.getY() / rho);
310
311
312 int index = 0;
313 double value = 0;
314 final double[] gradient = new double[3];
315 for (int m = degree; m >= 0; --m) {
316
317
318 index = computeTesseral(m, degree, index, t, u, tOu,
319 pnm0Plus2, pnm0Plus1, null, pnm0, pnm1, null);
320
321 if (m <= order) {
322
323
324
325 double sumDegreeS = 0;
326 double sumDegreeC = 0;
327 double dSumDegreeSdR = 0;
328 double dSumDegreeCdR = 0;
329 double dSumDegreeSdTheta = 0;
330 double dSumDegreeCdTheta = 0;
331 for (int n = FastMath.max(2, m); n <= degree; ++n) {
332 final double qSnm = aOrN[n] * harmonics.getNormalizedSnm(n, m);
333 final double qCnm = aOrN[n] * harmonics.getNormalizedCnm(n, m);
334 final double nOr = n / r;
335 final double s0 = pnm0[n] * qSnm;
336 final double c0 = pnm0[n] * qCnm;
337 final double s1 = pnm1[n] * qSnm;
338 final double c1 = pnm1[n] * qCnm;
339 sumDegreeS += s0;
340 sumDegreeC += c0;
341 dSumDegreeSdR -= nOr * s0;
342 dSumDegreeCdR -= nOr * c0;
343 dSumDegreeSdTheta += s1;
344 dSumDegreeCdTheta += c1;
345 }
346
347
348
349
350
351 final double sML = cosSinLambda[1][m];
352 final double cML = cosSinLambda[0][m];
353 value = value * u + sML * sumDegreeS + cML * sumDegreeC;
354 gradient[0] = gradient[0] * u + sML * dSumDegreeSdR + cML * dSumDegreeCdR;
355 gradient[1] = gradient[1] * u + m * (cML * sumDegreeS - sML * sumDegreeC);
356 gradient[2] = gradient[2] * u + sML * dSumDegreeSdTheta + cML * dSumDegreeCdTheta;
357
358 }
359
360
361 final double[] tmp = pnm0Plus2;
362 pnm0Plus2 = pnm0Plus1;
363 pnm0Plus1 = pnm0;
364 pnm0 = tmp;
365
366 }
367
368
369 value = FastMath.scalb(value, SCALING);
370 gradient[0] = FastMath.scalb(gradient[0], SCALING);
371 gradient[1] = FastMath.scalb(gradient[1], SCALING);
372 gradient[2] = FastMath.scalb(gradient[2], SCALING);
373
374
375 final double muOr = mu / r;
376 value *= muOr;
377 gradient[0] = muOr * gradient[0] - value / r;
378 gradient[1] *= muOr;
379 gradient[2] *= muOr;
380
381
382 return new SphericalCoordinates(position).toCartesianGradient(gradient);
383
384 }
385
386
387
388
389
390
391
392 public GradientHessian gradientHessian(final AbsoluteDate date, final Vector3D position)
393 throws OrekitException {
394
395 final int degree = provider.getMaxDegree();
396 final int order = provider.getMaxOrder();
397 final NormalizedSphericalHarmonics harmonics = provider.onDate(date);
398
399
400 double[] pnm0Plus2 = new double[degree + 1];
401 double[] pnm0Plus1 = new double[degree + 1];
402 double[] pnm0 = new double[degree + 1];
403 double[] pnm1Plus1 = new double[degree + 1];
404 double[] pnm1 = new double[degree + 1];
405 final double[] pnm2 = new double[degree + 1];
406
407
408 final double x = position.getX();
409 final double y = position.getY();
410 final double z = position.getZ();
411 final double x2 = x * x;
412 final double y2 = y * y;
413 final double z2 = z * z;
414 final double r2 = x2 + y2 + z2;
415 final double r = FastMath.sqrt (r2);
416 final double rho2 = x2 + y2;
417 final double rho = FastMath.sqrt(rho2);
418 final double t = z / r;
419 final double u = rho / r;
420 final double tOu = z / rho;
421
422
423 final double[] aOrN = createDistancePowersArray(provider.getAe() / r);
424
425
426 final double[][] cosSinLambda = createCosSinArrays(position.getX() / rho, position.getY() / rho);
427
428
429 int index = 0;
430 double value = 0;
431 final double[] gradient = new double[3];
432 final double[][] hessian = new double[3][3];
433 for (int m = degree; m >= 0; --m) {
434
435
436 index = computeTesseral(m, degree, index, t, u, tOu,
437 pnm0Plus2, pnm0Plus1, pnm1Plus1, pnm0, pnm1, pnm2);
438
439 if (m <= order) {
440
441
442
443 double sumDegreeS = 0;
444 double sumDegreeC = 0;
445 double dSumDegreeSdR = 0;
446 double dSumDegreeCdR = 0;
447 double dSumDegreeSdTheta = 0;
448 double dSumDegreeCdTheta = 0;
449 double d2SumDegreeSdRdR = 0;
450 double d2SumDegreeSdRdTheta = 0;
451 double d2SumDegreeSdThetadTheta = 0;
452 double d2SumDegreeCdRdR = 0;
453 double d2SumDegreeCdRdTheta = 0;
454 double d2SumDegreeCdThetadTheta = 0;
455 for (int n = FastMath.max(2, m); n <= degree; ++n) {
456 final double qSnm = aOrN[n] * harmonics.getNormalizedSnm(n, m);
457 final double qCnm = aOrN[n] * harmonics.getNormalizedCnm(n, m);
458 final double nOr = n / r;
459 final double nnP1Or2 = nOr * (n + 1) / r;
460 final double s0 = pnm0[n] * qSnm;
461 final double c0 = pnm0[n] * qCnm;
462 final double s1 = pnm1[n] * qSnm;
463 final double c1 = pnm1[n] * qCnm;
464 final double s2 = pnm2[n] * qSnm;
465 final double c2 = pnm2[n] * qCnm;
466 sumDegreeS += s0;
467 sumDegreeC += c0;
468 dSumDegreeSdR -= nOr * s0;
469 dSumDegreeCdR -= nOr * c0;
470 dSumDegreeSdTheta += s1;
471 dSumDegreeCdTheta += c1;
472 d2SumDegreeSdRdR += nnP1Or2 * s0;
473 d2SumDegreeSdRdTheta -= nOr * s1;
474 d2SumDegreeSdThetadTheta += s2;
475 d2SumDegreeCdRdR += nnP1Or2 * c0;
476 d2SumDegreeCdRdTheta -= nOr * c1;
477 d2SumDegreeCdThetadTheta += c2;
478 }
479
480
481 final double sML = cosSinLambda[1][m];
482 final double cML = cosSinLambda[0][m];
483 value = value * u + sML * sumDegreeS + cML * sumDegreeC;
484 gradient[0] = gradient[0] * u + sML * dSumDegreeSdR + cML * dSumDegreeCdR;
485 gradient[1] = gradient[1] * u + m * (cML * sumDegreeS - sML * sumDegreeC);
486 gradient[2] = gradient[2] * u + sML * dSumDegreeSdTheta + cML * dSumDegreeCdTheta;
487 hessian[0][0] = hessian[0][0] * u + sML * d2SumDegreeSdRdR + cML * d2SumDegreeCdRdR;
488 hessian[1][0] = hessian[1][0] * u + m * (cML * dSumDegreeSdR - sML * dSumDegreeCdR);
489 hessian[2][0] = hessian[2][0] * u + sML * d2SumDegreeSdRdTheta + cML * d2SumDegreeCdRdTheta;
490 hessian[1][1] = hessian[1][1] * u - m * m * (sML * sumDegreeS + cML * sumDegreeC);
491 hessian[2][1] = hessian[2][1] * u + m * (cML * dSumDegreeSdTheta - sML * dSumDegreeCdTheta);
492 hessian[2][2] = hessian[2][2] * u + sML * d2SumDegreeSdThetadTheta + cML * d2SumDegreeCdThetadTheta;
493
494 }
495
496
497 final double[] tmp0 = pnm0Plus2;
498 pnm0Plus2 = pnm0Plus1;
499 pnm0Plus1 = pnm0;
500 pnm0 = tmp0;
501 final double[] tmp1 = pnm1Plus1;
502 pnm1Plus1 = pnm1;
503 pnm1 = tmp1;
504
505 }
506
507
508 value = FastMath.scalb(value, SCALING);
509 for (int i = 0; i < 3; ++i) {
510 gradient[i] = FastMath.scalb(gradient[i], SCALING);
511 for (int j = 0; j <= i; ++j) {
512 hessian[i][j] = FastMath.scalb(hessian[i][j], SCALING);
513 }
514 }
515
516
517
518 final double muOr = mu / r;
519 value *= muOr;
520 gradient[0] = muOr * gradient[0] - value / r;
521 gradient[1] *= muOr;
522 gradient[2] *= muOr;
523 hessian[0][0] = muOr * hessian[0][0] - 2 * gradient[0] / r;
524 hessian[1][0] = muOr * hessian[1][0] - gradient[1] / r;
525 hessian[2][0] = muOr * hessian[2][0] - gradient[2] / r;
526 hessian[1][1] *= muOr;
527 hessian[2][1] *= muOr;
528 hessian[2][2] *= muOr;
529
530
531 final SphericalCoordinates sc = new SphericalCoordinates(position);
532 return new GradientHessian(sc.toCartesianGradient(gradient),
533 sc.toCartesianHessian(hessian, gradient));
534
535
536 }
537
538
539 public static class GradientHessian implements Serializable {
540
541
542 private static final long serialVersionUID = 20130219L;
543
544
545 private final double[] gradient;
546
547
548 private final double[][] hessian;
549
550
551
552
553
554
555
556
557 public GradientHessian(final double[] gradient, final double[][] hessian) {
558 this.gradient = gradient;
559 this.hessian = hessian;
560 }
561
562
563
564
565 public double[] getGradient() {
566 return gradient;
567 }
568
569
570
571
572 public double[][] getHessian() {
573 return hessian;
574 }
575
576 }
577
578
579
580
581
582 private double[] createDistancePowersArray(final double aOr) {
583
584
585 final double[] aOrN = new double[provider.getMaxDegree() + 1];
586 aOrN[0] = 1;
587 aOrN[1] = aOr;
588
589
590 for (int n = 2; n < aOrN.length; ++n) {
591 final int p = n / 2;
592 final int q = n - p;
593 aOrN[n] = aOrN[p] * aOrN[q];
594 }
595
596 return aOrN;
597
598 }
599
600
601
602
603
604
605
606 private double[][] createCosSinArrays(final double cosLambda, final double sinLambda) {
607
608
609 final double[][] cosSin = new double[2][provider.getMaxOrder() + 1];
610 cosSin[0][0] = 1;
611 cosSin[1][0] = 0;
612 if (provider.getMaxOrder() > 0) {
613 cosSin[0][1] = cosLambda;
614 cosSin[1][1] = sinLambda;
615
616
617 for (int m = 2; m < cosSin[0].length; ++m) {
618
619
620
621
622
623
624 final int p = m / 2;
625 final int q = m - p;
626
627 cosSin[0][m] = cosSin[0][p] * cosSin[0][q] - cosSin[1][p] * cosSin[1][q];
628 cosSin[1][m] = cosSin[1][p] * cosSin[0][q] + cosSin[0][p] * cosSin[1][q];
629
630 }
631 }
632
633 return cosSin;
634
635 }
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658 private int computeTesseral(final int m, final int degree, final int index,
659 final double t, final double u, final double tOu,
660 final double[] pnm0Plus2, final double[] pnm0Plus1, final double[] pnm1Plus1,
661 final double[] pnm0, final double[] pnm1, final double[] pnm2) {
662
663 final double u2 = u * u;
664
665
666 int n = FastMath.max(2, m);
667 if (n == m) {
668 pnm0[n] = sectorial[n];
669 ++n;
670 }
671
672
673 int localIndex = index;
674 while (n <= degree) {
675
676
677 pnm0[n] = gnmOj[localIndex] * t * pnm0Plus1[n] - hnmOj[localIndex] * u2 * pnm0Plus2[n];
678
679 ++localIndex;
680 ++n;
681
682 }
683
684 if (pnm1 != null) {
685
686
687 n = FastMath.max(2, m);
688 if (n == m) {
689 pnm1[n] = m * tOu * pnm0[n];
690 ++n;
691 }
692
693
694 localIndex = index;
695 while (n <= degree) {
696
697
698 pnm1[n] = m * tOu * pnm0[n] - enm[localIndex] * u * pnm0Plus1[n];
699
700 ++localIndex;
701 ++n;
702
703 }
704
705 if (pnm2 != null) {
706
707
708 n = FastMath.max(2, m);
709 if (n == m) {
710 pnm2[n] = m * (tOu * pnm1[n] - pnm0[n] / u2);
711 ++n;
712 }
713
714
715 localIndex = index;
716 while (n <= degree) {
717
718
719 pnm2[n] = m * (tOu * pnm1[n] - pnm0[n] / u2) - enm[localIndex] * u * pnm1Plus1[n];
720
721 ++localIndex;
722 ++n;
723
724 }
725
726 }
727
728 }
729
730 return localIndex;
731
732 }
733
734
735 public void addContribution(final SpacecraftState s, final TimeDerivativesEquations adder)
736 throws OrekitException {
737
738
739 final AbsoluteDate date = s.getDate();
740 final Transform fromBodyFrame = bodyFrame.getTransformTo(s.getFrame(), date);
741 final Transform toBodyFrame = fromBodyFrame.getInverse();
742 final Vector3D position = toBodyFrame.transformPosition(s.getPVCoordinates().getPosition());
743
744
745 final Vector3D gInertial = fromBodyFrame.transformVector(new Vector3D(gradient(date, position)));
746
747 adder.addXYZAcceleration(gInertial.getX(), gInertial.getY(), gInertial.getZ());
748
749 }
750
751
752 public EventDetector[] getEventsDetectors() {
753 return new EventDetector[0];
754 }
755
756
757 public FieldVector3D<DerivativeStructure> accelerationDerivatives(final AbsoluteDate date, final Frame frame,
758 final FieldVector3D<DerivativeStructure> position, final FieldVector3D<DerivativeStructure> velocity,
759 final FieldRotation<DerivativeStructure> rotation, final DerivativeStructure mass)
760 throws OrekitException {
761
762
763 final Transform fromBodyFrame = bodyFrame.getTransformTo(frame, date);
764 final Transform toBodyFrame = fromBodyFrame.getInverse();
765 final Vector3D positionBody = toBodyFrame.transformPosition(position.toVector3D());
766
767
768 final GradientHessian gh = gradientHessian(date, positionBody);
769
770
771 final double[] gInertial = fromBodyFrame.transformVector(new Vector3D(gh.getGradient())).toArray();
772
773
774 final RealMatrix hBody = new Array2DRowRealMatrix(gh.getHessian(), false);
775 final RealMatrix rot = new Array2DRowRealMatrix(toBodyFrame.getRotation().getMatrix());
776 final RealMatrix hInertial = rot.transpose().multiply(hBody).multiply(rot);
777
778
779 final int parameters = mass.getFreeParameters();
780 final int order = mass.getOrder();
781 final double[] derivatives = new double[1 + parameters];
782 final DerivativeStructure[] accDer = new DerivativeStructure[3];
783 for (int i = 0; i < 3; ++i) {
784
785
786 derivatives[0] = gInertial[i];
787
788
789 derivatives[1] = hInertial.getEntry(i, 0);
790 derivatives[2] = hInertial.getEntry(i, 1);
791 derivatives[3] = hInertial.getEntry(i, 2);
792
793
794
795 accDer[i] = new DerivativeStructure(parameters, order, derivatives);
796
797 }
798
799 return new FieldVector3D<DerivativeStructure>(accDer);
800
801 }
802
803
804 public FieldVector3D<DerivativeStructure> accelerationDerivatives(final SpacecraftState s, final String paramName)
805 throws OrekitException, IllegalArgumentException {
806
807 complainIfNotSupported(paramName);
808
809
810 final AbsoluteDate date = s.getDate();
811 final Transform fromBodyFrame = bodyFrame.getTransformTo(s.getFrame(), date);
812 final Transform toBodyFrame = fromBodyFrame.getInverse();
813 final Vector3D position = toBodyFrame.transformPosition(s.getPVCoordinates().getPosition());
814
815
816 final Vector3D gInertial = fromBodyFrame.transformVector(new Vector3D(gradient(date, position)));
817
818 return new FieldVector3D<DerivativeStructure>(new DerivativeStructure(1, 1, gInertial.getX(), gInertial.getX() / mu),
819 new DerivativeStructure(1, 1, gInertial.getY(), gInertial.getY() / mu),
820 new DerivativeStructure(1, 1, gInertial.getZ(), gInertial.getZ() / mu));
821
822 }
823
824
825 public ParameterDriver[] getParametersDrivers() {
826 return parametersDrivers.clone();
827 }
828
829 }