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