1 /* Copyright 2002-2025 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.orbits;
18
19 import org.hipparchus.geometry.euclidean.threed.Vector3D;
20 import org.hipparchus.linear.DecompositionSolver;
21 import org.hipparchus.linear.MatrixUtils;
22 import org.hipparchus.linear.QRDecomposition;
23 import org.hipparchus.linear.RealMatrix;
24 import org.hipparchus.linear.RealVector;
25 import org.hipparchus.util.FastMath;
26 import org.hipparchus.util.MathArrays;
27 import org.orekit.errors.OrekitIllegalArgumentException;
28 import org.orekit.errors.OrekitInternalError;
29 import org.orekit.errors.OrekitMessages;
30 import org.orekit.frames.Frame;
31 import org.orekit.frames.StaticTransform;
32 import org.orekit.frames.Transform;
33 import org.orekit.time.AbsoluteDate;
34 import org.orekit.time.TimeOffset;
35 import org.orekit.utils.PVCoordinates;
36 import org.orekit.utils.ShiftablePVCoordinatesHolder;
37 import org.orekit.utils.TimeStampedPVCoordinates;
38
39 /**
40 * This class handles orbital parameters.
41
42 * <p>
43 * For user convenience, both the Cartesian and the equinoctial elements
44 * are provided by this class, regardless of the canonical representation
45 * implemented in the derived class (which may be classical Keplerian
46 * elements for example).
47 * </p>
48 * <p>
49 * The parameters are defined in a frame specified by the user. It is important
50 * to make sure this frame is consistent: it probably is inertial and centered
51 * on the central body. This information is used for example by some
52 * force models.
53 * </p>
54 * <p>
55 * Instance of this class are guaranteed to be immutable.
56 * </p>
57 * @author Luc Maisonobe
58 * @author Guylaine Prat
59 * @author Fabien Maussion
60 * @author Véronique Pommier-Maurussane
61 */
62 public abstract class Orbit
63 implements ShiftablePVCoordinatesHolder<Orbit> {
64
65 /** Absolute tolerance when checking if the rate of the position angle is Keplerian or not. */
66 protected static final double TOLERANCE_POSITION_ANGLE_RATE = 1e-15;
67
68 /** Frame in which are defined the orbital parameters. */
69 private final Frame frame;
70
71 /** Date of the orbital parameters. */
72 private final AbsoluteDate date;
73
74 /** Value of mu used to compute position and velocity (m³/s²). */
75 private final double mu;
76
77 /** Computed position.
78 * @since 12.0
79 */
80 private Vector3D position;
81
82 /** Computed PVCoordinates. */
83 private TimeStampedPVCoordinates pvCoordinates;
84
85 /** Jacobian of the orbital parameters with mean angle with respect to the Cartesian coordinates. */
86 private double[][] jacobianMeanWrtCartesian;
87
88 /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with mean angle. */
89 private double[][] jacobianWrtParametersMean;
90
91 /** Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian coordinates. */
92 private double[][] jacobianEccentricWrtCartesian;
93
94 /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with eccentric angle. */
95 private double[][] jacobianWrtParametersEccentric;
96
97 /** Jacobian of the orbital parameters with true angle with respect to the Cartesian coordinates. */
98 private double[][] jacobianTrueWrtCartesian;
99
100 /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with true angle. */
101 private double[][] jacobianWrtParametersTrue;
102
103 /** Default constructor.
104 * Build a new instance with arbitrary default elements.
105 * @param frame the frame in which the parameters are defined
106 * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
107 * @param date date of the orbital parameters
108 * @param mu central attraction coefficient (m^3/s^2)
109 * @exception IllegalArgumentException if frame is not a {@link
110 * Frame#isPseudoInertial pseudo-inertial frame}
111 */
112 protected Orbit(final Frame frame, final AbsoluteDate date, final double mu)
113 throws IllegalArgumentException {
114 ensurePseudoInertialFrame(frame);
115 this.date = date;
116 this.mu = mu;
117 this.pvCoordinates = null;
118 this.frame = frame;
119 jacobianMeanWrtCartesian = null;
120 jacobianWrtParametersMean = null;
121 jacobianEccentricWrtCartesian = null;
122 jacobianWrtParametersEccentric = null;
123 jacobianTrueWrtCartesian = null;
124 jacobianWrtParametersTrue = null;
125 }
126
127 /** Set the orbit from Cartesian parameters.
128 *
129 * <p> The acceleration provided in {@code pvCoordinates} is accessible using
130 * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
131 * use {@code mu} and the position to compute the acceleration, including
132 * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
133 *
134 * @param pvCoordinates the position and velocity in the inertial frame
135 * @param frame the frame in which the {@link TimeStampedPVCoordinates} are defined
136 * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
137 * @param mu central attraction coefficient (m^3/s^2)
138 * @exception IllegalArgumentException if frame is not a {@link
139 * Frame#isPseudoInertial pseudo-inertial frame}
140 */
141 protected Orbit(final TimeStampedPVCoordinates pvCoordinates, final Frame frame, final double mu)
142 throws IllegalArgumentException {
143 ensurePseudoInertialFrame(frame);
144 this.date = pvCoordinates.getDate();
145 this.mu = mu;
146 if (pvCoordinates.getAcceleration().getNormSq() == 0) {
147 // the acceleration was not provided,
148 // compute it from Newtonian attraction
149 final double r2 = pvCoordinates.getPosition().getNormSq();
150 final double r3 = r2 * FastMath.sqrt(r2);
151 this.pvCoordinates = new TimeStampedPVCoordinates(pvCoordinates.getDate(),
152 pvCoordinates.getPosition(),
153 pvCoordinates.getVelocity(),
154 new Vector3D(-mu / r3, pvCoordinates.getPosition()));
155 } else {
156 this.pvCoordinates = pvCoordinates;
157 }
158 this.frame = frame;
159 }
160
161 /** Compute non-Keplerian part of the acceleration from first time derivatives.
162 * @return non-Keplerian part of the acceleration
163 * @since 13.1
164 */
165 protected Vector3D nonKeplerianAcceleration() {
166
167 final double[][] dPdC = new double[6][6];
168 final PositionAngleType positionAngleType = PositionAngleType.MEAN;
169 getJacobianWrtCartesian(positionAngleType, dPdC);
170 final RealMatrix subMatrix = MatrixUtils.createRealMatrix(dPdC);
171
172 final DecompositionSolver solver = getDecompositionSolver(subMatrix);
173
174 final double[] derivatives = new double[6];
175 getType().mapOrbitToArray(this, positionAngleType, new double[6], derivatives);
176 derivatives[5] -= getKeplerianMeanMotion();
177
178 final RealVector solution = solver.solve(MatrixUtils.createRealVector(derivatives));
179 return new Vector3D(solution.getEntry(3), solution.getEntry(4), solution.getEntry(5));
180
181 }
182
183 /** Check if Cartesian coordinates include non-Keplerian acceleration.
184 * @param pva Cartesian coordinates
185 * @param mu central attraction coefficient
186 * @return true if Cartesian coordinates include non-Keplerian acceleration
187 */
188 protected static boolean hasNonKeplerianAcceleration(final PVCoordinates pva, final double mu) {
189
190 final Vector3D a = pva.getAcceleration();
191 if (a == null) {
192 return false;
193 }
194
195 final Vector3D p = pva.getPosition();
196 final double r2 = p.getNormSq();
197
198 // Check if acceleration is relatively close to 0 compared to the Keplerian acceleration
199 final double tolerance = mu * 1e-9;
200 final Vector3D aTimesR2 = a.scalarMultiply(r2);
201 if (aTimesR2.getNorm() < tolerance) {
202 return false;
203 }
204
205 if ((aTimesR2.add(p.normalize().scalarMultiply(mu))).getNorm() > tolerance) {
206 // we have a relevant acceleration, we can compute derivatives
207 return true;
208 } else {
209 // the provided acceleration is either too small to be reliable (probably even 0), or NaN
210 return false;
211 }
212 }
213
214 /** Returns true if and only if the orbit is elliptical i.e. has a non-negative semi-major axis.
215 * @return true if getA() is strictly greater than 0
216 * @since 12.0
217 */
218 public boolean isElliptical() {
219 return getA() > 0.;
220 }
221
222 /** Get the orbit type.
223 * @return orbit type
224 */
225 public abstract OrbitType getType();
226
227 /** Ensure the defining frame is a pseudo-inertial frame.
228 * @param frame frame to check
229 * @exception IllegalArgumentException if frame is not a {@link
230 * Frame#isPseudoInertial pseudo-inertial frame}
231 */
232 private static void ensurePseudoInertialFrame(final Frame frame)
233 throws IllegalArgumentException {
234 if (!frame.isPseudoInertial()) {
235 throw new OrekitIllegalArgumentException(OrekitMessages.NON_PSEUDO_INERTIAL_FRAME,
236 frame.getName());
237 }
238 }
239
240 /** Get the frame in which the orbital parameters are defined.
241 * @return frame in which the orbital parameters are defined
242 */
243 public Frame getFrame() {
244 return frame;
245 }
246
247 /** Get the semi-major axis.
248 * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
249 * @return semi-major axis (m)
250 */
251 public abstract double getA();
252
253 /** Get the semi-major axis derivative.
254 * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
255 * <p>
256 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
257 * </p>
258 * @return semi-major axis derivative (m/s)
259 * @since 9.0
260 */
261 public abstract double getADot();
262
263 /** Get the first component of the equinoctial eccentricity vector.
264 * @return first component of the equinoctial eccentricity vector
265 */
266 public abstract double getEquinoctialEx();
267
268 /** Get the first component of the equinoctial eccentricity vector derivative.
269 * <p>
270 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
271 * </p>
272 * @return first component of the equinoctial eccentricity vector derivative
273 * @since 9.0
274 */
275 public abstract double getEquinoctialExDot();
276
277 /** Get the second component of the equinoctial eccentricity vector.
278 * @return second component of the equinoctial eccentricity vector
279 */
280 public abstract double getEquinoctialEy();
281
282 /** Get the second component of the equinoctial eccentricity vector derivative.
283 * <p>
284 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
285 * </p>
286 * @return second component of the equinoctial eccentricity vector derivative
287 * @since 9.0
288 */
289 public abstract double getEquinoctialEyDot();
290
291 /** Get the first component of the inclination vector.
292 * @return first component of the inclination vector
293 */
294 public abstract double getHx();
295
296 /** Get the first component of the inclination vector derivative.
297 * <p>
298 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
299 * </p>
300 * @return first component of the inclination vector derivative
301
302 * @since 9.0
303 */
304 public abstract double getHxDot();
305
306 /** Get the second component of the inclination vector.
307 * @return second component of the inclination vector
308 */
309 public abstract double getHy();
310
311 /** Get the second component of the inclination vector derivative.
312 * <p>
313 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
314 * </p>
315 * @return second component of the inclination vector derivative
316 * @since 9.0
317 */
318 public abstract double getHyDot();
319
320 /** Get the eccentric longitude argument.
321 * @return E + ω + Ω eccentric longitude argument (rad)
322 */
323 public abstract double getLE();
324
325 /** Get the eccentric longitude argument derivative.
326 * <p>
327 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
328 * </p>
329 * @return d(E + ω + Ω)/dt eccentric longitude argument derivative (rad/s)
330 * @since 9.0
331 */
332 public abstract double getLEDot();
333
334 /** Get the true longitude argument.
335 * @return v + ω + Ω true longitude argument (rad)
336 */
337 public abstract double getLv();
338
339 /** Get the true longitude argument derivative.
340 * <p>
341 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
342 * </p>
343 * @return d(v + ω + Ω)/dt true longitude argument derivative (rad/s)
344 * @since 9.0
345 */
346 public abstract double getLvDot();
347
348 /** Get the mean longitude argument.
349 * @return M + ω + Ω mean longitude argument (rad)
350 */
351 public abstract double getLM();
352
353 /** Get the mean longitude argument derivative.
354 * <p>
355 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
356 * </p>
357 * @return d(M + ω + Ω)/dt mean longitude argument derivative (rad/s)
358 * @since 9.0
359 */
360 public abstract double getLMDot();
361
362 // Additional orbital elements
363
364 /** Get the eccentricity.
365 * @return eccentricity
366 */
367 public abstract double getE();
368
369 /** Get the eccentricity derivative.
370 * <p>
371 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
372 * </p>
373 * @return eccentricity derivative
374 * @since 9.0
375 */
376 public abstract double getEDot();
377
378 /** Get the inclination.
379 * @return inclination (rad)
380 */
381 public abstract double getI();
382
383 /** Get the inclination derivative.
384 * <p>
385 * If the orbit was created without derivatives, the value returned is {@link Double#NaN}.
386 * </p>
387 * @return inclination derivative (rad/s)
388 * @since 9.0
389 */
390 public abstract double getIDot();
391
392 /** Check if orbit includes non-Keplerian rates.
393 * @return true if orbit includes non-Keplerian derivatives
394 * @see #getADot()
395 * @see #getEquinoctialExDot()
396 * @see #getEquinoctialEyDot()
397 * @see #getHxDot()
398 * @see #getHyDot()
399 * @see #getLEDot()
400 * @see #getLvDot()
401 * @see #getLMDot()
402 * @see #getEDot()
403 * @see #getIDot()
404 * @since 13.0
405 */
406 public boolean hasNonKeplerianAcceleration() {
407 return hasNonKeplerianAcceleration(getPVCoordinates(), getMu());
408 }
409
410 /** Get the central acceleration constant.
411 * @return central acceleration constant
412 */
413 public double getMu() {
414 return mu;
415 }
416
417 /** Get the Keplerian period.
418 * <p>The Keplerian period is computed directly from semi major axis
419 * and central acceleration constant.</p>
420 * @return Keplerian period in seconds, or positive infinity for hyperbolic orbits
421 */
422 public double getKeplerianPeriod() {
423 final double a = getA();
424 return isElliptical() ? 2.0 * FastMath.PI * a * FastMath.sqrt(a / mu) : Double.POSITIVE_INFINITY;
425 }
426
427 /** Get the Keplerian mean motion.
428 * <p>The Keplerian mean motion is computed directly from semi major axis
429 * and central acceleration constant.</p>
430 * @return Keplerian mean motion in radians per second
431 */
432 public double getKeplerianMeanMotion() {
433 final double absA = FastMath.abs(getA());
434 return FastMath.sqrt(mu / absA) / absA;
435 }
436
437 /** Get the derivative of the mean anomaly with respect to the semi major axis.
438 * @return derivative of the mean anomaly with respect to the semi major axis
439 */
440 public double getMeanAnomalyDotWrtA() {
441 return -1.5 * getKeplerianMeanMotion() / getA();
442 }
443
444 /** Get the date of orbital parameters.
445 * @return date of the orbital parameters
446 */
447 public AbsoluteDate getDate() {
448 return date;
449 }
450
451 /** Get the {@link TimeStampedPVCoordinates} in a specified frame.
452 * @param outputFrame frame in which the position/velocity coordinates shall be computed
453 * @return pvCoordinates in the specified output frame
454 * @see #getPVCoordinates()
455 */
456 public TimeStampedPVCoordinates getPVCoordinates(final Frame outputFrame) {
457 if (pvCoordinates == null) {
458 pvCoordinates = initPVCoordinates();
459 }
460
461 // If output frame requested is the same as definition frame,
462 // PV coordinates are returned directly
463 if (outputFrame == frame) {
464 return pvCoordinates;
465 }
466
467 // Else, PV coordinates are transformed to output frame
468 final Transform t = frame.getTransformTo(outputFrame, date);
469 return t.transformPVCoordinates(pvCoordinates);
470 }
471
472 /** Get the position in a specified frame.
473 * @param outputFrame frame in which the position coordinates shall be computed
474 * @return position in the specified output frame
475 * @see #getPosition()
476 * @since 12.0
477 */
478 public Vector3D getPosition(final Frame outputFrame) {
479 if (position == null) {
480 position = initPosition();
481 }
482
483 // If output frame requested is the same as definition frame,
484 // Position vector is returned directly
485 if (outputFrame == frame) {
486 return position;
487 }
488
489 // Else, position vector is transformed to output frame
490 final StaticTransform t = frame.getStaticTransformTo(outputFrame, date);
491 return t.transformPosition(position);
492
493 }
494
495 /** Get the position in definition frame.
496 * @return position in the definition frame
497 * @see #getPVCoordinates()
498 * @since 12.0
499 */
500 public Vector3D getPosition() {
501 if (position == null) {
502 position = initPosition();
503 }
504 return position;
505 }
506
507 /** Get the {@link TimeStampedPVCoordinates} in definition frame.
508 * @return pvCoordinates in the definition frame
509 * @see #getPVCoordinates(Frame)
510 */
511 public TimeStampedPVCoordinates getPVCoordinates() {
512 if (pvCoordinates == null) {
513 pvCoordinates = initPVCoordinates();
514 position = pvCoordinates.getPosition();
515 }
516 return pvCoordinates;
517 }
518
519 /** Compute the position coordinates from the canonical parameters.
520 * @return computed position coordinates
521 * @since 12.0
522 */
523 protected abstract Vector3D initPosition();
524
525 /** Compute the position/velocity coordinates from the canonical parameters.
526 * @return computed position/velocity coordinates
527 */
528 protected abstract TimeStampedPVCoordinates initPVCoordinates();
529
530 /**
531 * Create a new object representing the same physical orbital state, but attached to a different reference frame.
532 * If the new frame is not inertial, an exception will be thrown.
533 *
534 * @param inertialFrame reference frame of output orbit
535 * @return orbit with different frame
536 * @since 13.0
537 */
538 public abstract Orbit inFrame(Frame inertialFrame);
539
540 /** Get a time-shifted orbit.
541 * <p>
542 * The orbit can be slightly shifted to close dates. The shifting model is a
543 * Keplerian one if no derivatives are available in the orbit, or Keplerian
544 * plus quadratic effect of the non-Keplerian acceleration if derivatives are
545 * available. Shifting is <em>not</em> intended as a replacement for proper
546 * orbit propagation but should be sufficient for small time shifts or coarse
547 * accuracy.
548 * </p>
549 * @param dt time shift in seconds
550 * @return a new orbit, shifted with respect to the instance (which is immutable)
551 */
552 @Override
553 public abstract Orbit shiftedBy(double dt);
554
555 /** Get a time-shifted orbit.
556 * <p>
557 * The orbit can be slightly shifted to close dates. The shifting model is a
558 * Keplerian one if no derivatives are available in the orbit, or Keplerian
559 * plus quadratic effect of the non-Keplerian acceleration if derivatives are
560 * available. Shifting is <em>not</em> intended as a replacement for proper
561 * orbit propagation but should be sufficient for small time shifts or coarse
562 * accuracy.
563 * </p>
564 * @param dt time shift
565 * @return a new orbit, shifted with respect to the instance (which is immutable)
566 */
567 @Override
568 public abstract Orbit shiftedBy(TimeOffset dt);
569
570 /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
571 * <p>
572 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
573 * respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter
574 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
575 * </p>
576 * @param type type of the position angle to use
577 * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
578 * is larger than 6x6, only the 6x6 upper left corner will be modified
579 */
580 public void getJacobianWrtCartesian(final PositionAngleType type, final double[][] jacobian) {
581
582 final double[][] cachedJacobian;
583 synchronized (this) {
584 switch (type) {
585 case MEAN :
586 if (jacobianMeanWrtCartesian == null) {
587 // first call, we need to compute the Jacobian and cache it
588 jacobianMeanWrtCartesian = computeJacobianMeanWrtCartesian();
589 }
590 cachedJacobian = jacobianMeanWrtCartesian;
591 break;
592 case ECCENTRIC :
593 if (jacobianEccentricWrtCartesian == null) {
594 // first call, we need to compute the Jacobian and cache it
595 jacobianEccentricWrtCartesian = computeJacobianEccentricWrtCartesian();
596 }
597 cachedJacobian = jacobianEccentricWrtCartesian;
598 break;
599 case TRUE :
600 if (jacobianTrueWrtCartesian == null) {
601 // first call, we need to compute the Jacobian and cache it
602 jacobianTrueWrtCartesian = computeJacobianTrueWrtCartesian();
603 }
604 cachedJacobian = jacobianTrueWrtCartesian;
605 break;
606 default :
607 throw new OrekitInternalError(null);
608 }
609 }
610
611 // fill the user provided array
612 for (int i = 0; i < cachedJacobian.length; ++i) {
613 System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
614 }
615
616 }
617
618 /** Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.
619 * <p>
620 * Element {@code jacobian[i][j]} is the derivative of Cartesian coordinate i of the orbit with
621 * respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate
622 * x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters.
623 * </p>
624 * @param type type of the position angle to use
625 * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
626 * is larger than 6x6, only the 6x6 upper left corner will be modified
627 */
628 public void getJacobianWrtParameters(final PositionAngleType type, final double[][] jacobian) {
629
630 final double[][] cachedJacobian;
631 synchronized (this) {
632 switch (type) {
633 case MEAN :
634 if (jacobianWrtParametersMean == null) {
635 // first call, we need to compute the Jacobian and cache it
636 jacobianWrtParametersMean = createInverseJacobian(type);
637 }
638 cachedJacobian = jacobianWrtParametersMean;
639 break;
640 case ECCENTRIC :
641 if (jacobianWrtParametersEccentric == null) {
642 // first call, we need to compute the Jacobian and cache it
643 jacobianWrtParametersEccentric = createInverseJacobian(type);
644 }
645 cachedJacobian = jacobianWrtParametersEccentric;
646 break;
647 case TRUE :
648 if (jacobianWrtParametersTrue == null) {
649 // first call, we need to compute the Jacobian and cache it
650 jacobianWrtParametersTrue = createInverseJacobian(type);
651 }
652 cachedJacobian = jacobianWrtParametersTrue;
653 break;
654 default :
655 throw new OrekitInternalError(null);
656 }
657 }
658
659 // fill the user-provided array
660 for (int i = 0; i < cachedJacobian.length; ++i) {
661 System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
662 }
663
664 }
665
666 /** Create an inverse Jacobian.
667 * @param type type of the position angle to use
668 * @return inverse Jacobian
669 */
670 private double[][] createInverseJacobian(final PositionAngleType type) {
671
672 // get the direct Jacobian
673 final double[][] directJacobian = new double[6][6];
674 getJacobianWrtCartesian(type, directJacobian);
675
676 // invert the direct Jacobian
677 final RealMatrix matrix = MatrixUtils.createRealMatrix(directJacobian);
678 final DecompositionSolver solver = getDecompositionSolver(matrix);
679 return solver.getInverse().getData();
680
681 }
682
683 /**
684 * Method to build a matrix decomposition solver.
685 * @param realMatrix matrix
686 * @return solver
687 * @since 13.1
688 */
689 protected DecompositionSolver getDecompositionSolver(final RealMatrix realMatrix) {
690 return new QRDecomposition(realMatrix).getSolver();
691 }
692
693 /** Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
694 * <p>
695 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
696 * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
697 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
698 * </p>
699 * <p>
700 * The array returned by this method will not be modified.
701 * </p>
702 * @return 6x6 Jacobian matrix
703 * @see #computeJacobianEccentricWrtCartesian()
704 * @see #computeJacobianTrueWrtCartesian()
705 */
706 protected abstract double[][] computeJacobianMeanWrtCartesian();
707
708 /** Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.
709 * <p>
710 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
711 * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
712 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
713 * </p>
714 * <p>
715 * The array returned by this method will not be modified.
716 * </p>
717 * @return 6x6 Jacobian matrix
718 * @see #computeJacobianMeanWrtCartesian()
719 * @see #computeJacobianTrueWrtCartesian()
720 */
721 protected abstract double[][] computeJacobianEccentricWrtCartesian();
722
723 /** Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.
724 * <p>
725 * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
726 * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
727 * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
728 * </p>
729 * <p>
730 * The array returned by this method will not be modified.
731 * </p>
732 * @return 6x6 Jacobian matrix
733 * @see #computeJacobianMeanWrtCartesian()
734 * @see #computeJacobianEccentricWrtCartesian()
735 */
736 protected abstract double[][] computeJacobianTrueWrtCartesian();
737
738 /** Add the contribution of the Keplerian motion to parameters derivatives
739 * <p>
740 * This method is used by integration-based propagators to evaluate the part of Keplerian
741 * motion to evolution of the orbital state.
742 * </p>
743 * @param type type of the position angle in the state
744 * @param gm attraction coefficient to use
745 * @param pDot array containing orbital state derivatives to update (the Keplerian
746 * part must be <em>added</em> to the array components, as the array may already
747 * contain some non-zero elements corresponding to non-Keplerian parts)
748 */
749 public abstract void addKeplerContribution(PositionAngleType type, double gm, double[] pDot);
750
751 /** Fill a Jacobian half row with a single vector.
752 * @param a coefficient of the vector
753 * @param v vector
754 * @param row Jacobian matrix row
755 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
756 */
757 protected static void fillHalfRow(final double a, final Vector3D v, final double[] row, final int j) {
758 row[j] = a * v.getX();
759 row[j + 1] = a * v.getY();
760 row[j + 2] = a * v.getZ();
761 }
762
763 /** Fill a Jacobian half row with a linear combination of vectors.
764 * @param a1 coefficient of the first vector
765 * @param v1 first vector
766 * @param a2 coefficient of the second vector
767 * @param v2 second vector
768 * @param row Jacobian matrix row
769 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
770 */
771 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
772 final double[] row, final int j) {
773 row[j] = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX());
774 row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY());
775 row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ());
776 }
777
778 /** Fill a Jacobian half row with a linear combination of vectors.
779 * @param a1 coefficient of the first vector
780 * @param v1 first vector
781 * @param a2 coefficient of the second vector
782 * @param v2 second vector
783 * @param a3 coefficient of the third vector
784 * @param v3 third vector
785 * @param row Jacobian matrix row
786 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
787 */
788 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
789 final double a3, final Vector3D v3,
790 final double[] row, final int j) {
791 row[j] = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX());
792 row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY());
793 row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ());
794 }
795
796 /** Fill a Jacobian half row with a linear combination of vectors.
797 * @param a1 coefficient of the first vector
798 * @param v1 first vector
799 * @param a2 coefficient of the second vector
800 * @param v2 second vector
801 * @param a3 coefficient of the third vector
802 * @param v3 third vector
803 * @param a4 coefficient of the fourth vector
804 * @param v4 fourth vector
805 * @param row Jacobian matrix row
806 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
807 */
808 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
809 final double a3, final Vector3D v3, final double a4, final Vector3D v4,
810 final double[] row, final int j) {
811 row[j] = MathArrays.linearCombination(a1, v1.getX(), a2, v2.getX(), a3, v3.getX(), a4, v4.getX());
812 row[j + 1] = MathArrays.linearCombination(a1, v1.getY(), a2, v2.getY(), a3, v3.getY(), a4, v4.getY());
813 row[j + 2] = MathArrays.linearCombination(a1, v1.getZ(), a2, v2.getZ(), a3, v3.getZ(), a4, v4.getZ());
814 }
815
816 /** Fill a Jacobian half row with a linear combination of vectors.
817 * @param a1 coefficient of the first vector
818 * @param v1 first vector
819 * @param a2 coefficient of the second vector
820 * @param v2 second vector
821 * @param a3 coefficient of the third vector
822 * @param v3 third vector
823 * @param a4 coefficient of the fourth vector
824 * @param v4 fourth vector
825 * @param a5 coefficient of the fifth vector
826 * @param v5 fifth vector
827 * @param row Jacobian matrix row
828 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
829 */
830 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
831 final double a3, final Vector3D v3, final double a4, final Vector3D v4,
832 final double a5, final Vector3D v5,
833 final double[] row, final int j) {
834 final double[] a = new double[] {
835 a1, a2, a3, a4, a5
836 };
837 row[j] = MathArrays.linearCombination(a, new double[] {
838 v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX()
839 });
840 row[j + 1] = MathArrays.linearCombination(a, new double[] {
841 v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY()
842 });
843 row[j + 2] = MathArrays.linearCombination(a, new double[] {
844 v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ()
845 });
846 }
847
848 /** Fill a Jacobian half row with a linear combination of vectors.
849 * @param a1 coefficient of the first vector
850 * @param v1 first vector
851 * @param a2 coefficient of the second vector
852 * @param v2 second vector
853 * @param a3 coefficient of the third vector
854 * @param v3 third vector
855 * @param a4 coefficient of the fourth vector
856 * @param v4 fourth vector
857 * @param a5 coefficient of the fifth vector
858 * @param v5 fifth vector
859 * @param a6 coefficient of the sixth vector
860 * @param v6 sixth vector
861 * @param row Jacobian matrix row
862 * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
863 */
864 protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
865 final double a3, final Vector3D v3, final double a4, final Vector3D v4,
866 final double a5, final Vector3D v5, final double a6, final Vector3D v6,
867 final double[] row, final int j) {
868 final double[] a = new double[] {
869 a1, a2, a3, a4, a5, a6
870 };
871 row[j] = MathArrays.linearCombination(a, new double[] {
872 v1.getX(), v2.getX(), v3.getX(), v4.getX(), v5.getX(), v6.getX()
873 });
874 row[j + 1] = MathArrays.linearCombination(a, new double[] {
875 v1.getY(), v2.getY(), v3.getY(), v4.getY(), v5.getY(), v6.getY()
876 });
877 row[j + 2] = MathArrays.linearCombination(a, new double[] {
878 v1.getZ(), v2.getZ(), v3.getZ(), v4.getZ(), v5.getZ(), v6.getZ()
879 });
880 }
881
882 }