1   /* Copyright 2002-2015 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.apache.commons.math3.geometry.euclidean.threed.Vector3D;
22  import org.apache.commons.math3.linear.DecompositionSolver;
23  import org.apache.commons.math3.linear.MatrixUtils;
24  import org.apache.commons.math3.linear.QRDecomposition;
25  import org.apache.commons.math3.linear.RealMatrix;
26  import org.apache.commons.math3.util.FastMath;
27  import org.orekit.errors.OrekitException;
28  import org.orekit.errors.OrekitMessages;
29  import org.orekit.frames.Frame;
30  import org.orekit.frames.Transform;
31  import org.orekit.time.AbsoluteDate;
32  import org.orekit.time.TimeInterpolable;
33  import org.orekit.time.TimeShiftable;
34  import org.orekit.time.TimeStamped;
35  import org.orekit.utils.PVCoordinatesProvider;
36  import org.orekit.utils.TimeStampedPVCoordinates;
37  
38  /**
39   * This class handles orbital parameters.
40  
41   * <p>
42   * For user convenience, both the Cartesian and the equinoctial elements
43   * are provided by this class, regardless of the canonical representation
44   * implemented in the derived class (which may be classical keplerian
45   * elements for example).
46   * </p>
47   * <p>
48   * The parameters are defined in a frame specified by the user. It is important
49   * to make sure this frame is consistent: it probably is inertial and centered
50   * on the central body. This information is used for example by some
51   * force models.
52   * </p>
53   * <p>
54   * The object <code>OrbitalParameters</code> is guaranteed to be immutable.
55   * </p>
56   * @author Luc Maisonobe
57   * @author Guylaine Prat
58   * @author Fabien Maussion
59   * @author V&eacute;ronique Pommier-Maurussane
60   */
61  public abstract class Orbit
62      implements TimeStamped, TimeShiftable<Orbit>, TimeInterpolable<Orbit>,
63                 Serializable, PVCoordinatesProvider {
64  
65      /** Serializable UID. */
66      private static final long serialVersionUID = 438733454597999578L;
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 PVCoordinates. */
78      private transient TimeStampedPVCoordinates pvCoordinates;
79  
80      /** Jacobian of the orbital parameters with mean angle with respect to the Cartesian coordinates. */
81      private transient double[][] jacobianMeanWrtCartesian;
82  
83      /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with mean angle. */
84      private transient double[][] jacobianWrtParametersMean;
85  
86      /** Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian coordinates. */
87      private transient double[][] jacobianEccentricWrtCartesian;
88  
89      /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with eccentric angle. */
90      private transient double[][] jacobianWrtParametersEccentric;
91  
92      /** Jacobian of the orbital parameters with true angle with respect to the Cartesian coordinates. */
93      private transient double[][] jacobianTrueWrtCartesian;
94  
95      /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with true angle. */
96      private transient double[][] jacobianWrtParametersTrue;
97  
98      /** Default constructor.
99       * Build a new instance with arbitrary default elements.
100      * @param frame the frame in which the parameters are defined
101      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
102      * @param date date of the orbital parameters
103      * @param mu central attraction coefficient (m^3/s^2)
104      * @exception IllegalArgumentException if frame is not a {@link
105      * Frame#isPseudoInertial pseudo-inertial frame}
106      */
107     protected Orbit(final Frame frame, final AbsoluteDate date, final double mu)
108         throws IllegalArgumentException {
109         ensurePseudoInertialFrame(frame);
110         this.date                      = date;
111         this.mu                        = mu;
112         this.pvCoordinates             = null;
113         this.frame                     = frame;
114         jacobianMeanWrtCartesian       = null;
115         jacobianWrtParametersMean      = null;
116         jacobianEccentricWrtCartesian  = null;
117         jacobianWrtParametersEccentric = null;
118         jacobianTrueWrtCartesian       = null;
119         jacobianWrtParametersTrue      = null;
120     }
121 
122     /** Set the orbit from Cartesian parameters.
123      *
124      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
125      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
126      * use {@code mu} and the position to compute the acceleration, including
127      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
128      *
129      * @param pvCoordinates the position and velocity in the inertial frame
130      * @param frame the frame in which the {@link TimeStampedPVCoordinates} are defined
131      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
132      * @param mu central attraction coefficient (m^3/s^2)
133      * @exception IllegalArgumentException if frame is not a {@link
134      * Frame#isPseudoInertial pseudo-inertial frame}
135      */
136     protected Orbit(final TimeStampedPVCoordinates pvCoordinates, final Frame frame, final double mu)
137         throws IllegalArgumentException {
138         ensurePseudoInertialFrame(frame);
139         this.date = pvCoordinates.getDate();
140         this.mu = mu;
141         if (pvCoordinates.getAcceleration().getNormSq() == 0) {
142             // the acceleration was not provided,
143             // compute it from Newtonian attraction
144             final double r2 = pvCoordinates.getPosition().getNormSq();
145             final double r3 = r2 * FastMath.sqrt(r2);
146             this.pvCoordinates = new TimeStampedPVCoordinates(pvCoordinates.getDate(),
147                                                               pvCoordinates.getPosition(),
148                                                               pvCoordinates.getVelocity(),
149                                                               new Vector3D(-mu / r3, pvCoordinates.getPosition()));
150         } else {
151             this.pvCoordinates = pvCoordinates;
152         }
153         this.frame = frame;
154     }
155 
156     /** Get the orbit type.
157      * @return orbit type
158      */
159     public abstract OrbitType getType();
160 
161     /** Ensure the defining frame is a pseudo-inertial frame.
162      * @param frame frame to check
163      * @exception IllegalArgumentException if frame is not a {@link
164      * Frame#isPseudoInertial pseudo-inertial frame}
165      */
166     private static void ensurePseudoInertialFrame(final Frame frame)
167         throws IllegalArgumentException {
168         if (!frame.isPseudoInertial()) {
169             throw OrekitException.createIllegalArgumentException(
170                 OrekitMessages.NON_PSEUDO_INERTIAL_FRAME_NOT_SUITABLE_FOR_DEFINING_ORBITS,
171                 frame.getName());
172         }
173     }
174 
175     /** Get the frame in which the orbital parameters are defined.
176      * @return frame in which the orbital parameters are defined
177      */
178     public Frame getFrame() {
179         return frame;
180     }
181 
182     /** Get the semi-major axis.
183      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
184      * @return semi-major axis (m)
185      */
186     public abstract double getA();
187 
188     /** Get the first component of the equinoctial eccentricity vector.
189      * @return first component of the equinoctial eccentricity vector
190      */
191     public abstract double getEquinoctialEx();
192 
193     /** Get the second component of the equinoctial eccentricity vector.
194      * @return second component of the equinoctial eccentricity vector
195      */
196     public abstract double getEquinoctialEy();
197 
198     /** Get the first component of the inclination vector.
199      * @return first component of the inclination vector
200      */
201     public abstract double getHx();
202 
203     /** Get the second component of the inclination vector.
204      * @return second component of the inclination vector
205      */
206     public abstract double getHy();
207 
208     /** Get the eccentric longitude argument.
209      * @return E + ω + Ω eccentric longitude argument (rad)
210      */
211     public abstract double getLE();
212 
213     /** Get the true longitude argument.
214      * @return v + ω + Ω true longitude argument (rad)
215      */
216     public abstract double getLv();
217 
218     /** Get the mean longitude argument.
219      * @return M + ω + Ω mean longitude argument (rad)
220      */
221     public abstract double getLM();
222 
223     // Additional orbital elements
224 
225     /** Get the eccentricity.
226      * @return eccentricity
227      */
228     public abstract double getE();
229 
230     /** Get the inclination.
231      * @return inclination (rad)
232      */
233     public abstract double getI();
234 
235     /** Get the central acceleration constant.
236      * @return central acceleration constant
237      */
238     public double getMu() {
239         return mu;
240     }
241 
242     /** Get the keplerian period.
243      * <p>The keplerian period is computed directly from semi major axis
244      * and central acceleration constant.</p>
245      * @return keplerian period in seconds, or positive infinity for hyperbolic orbits
246      */
247     public double getKeplerianPeriod() {
248         final double a = getA();
249         return (a < 0) ? Double.POSITIVE_INFINITY : 2.0 * FastMath.PI * a * FastMath.sqrt(a / mu);
250     }
251 
252     /** Get the keplerian mean motion.
253      * <p>The keplerian mean motion is computed directly from semi major axis
254      * and central acceleration constant.</p>
255      * @return keplerian mean motion in radians per second
256      */
257     public double getKeplerianMeanMotion() {
258         final double absA = FastMath.abs(getA());
259         return FastMath.sqrt(mu / absA) / absA;
260     }
261 
262     /** Get the date of orbital parameters.
263      * @return date of the orbital parameters
264      */
265     public AbsoluteDate getDate() {
266         return date;
267     }
268 
269     /** Get the {@link TimeStampedPVCoordinates} in a specified frame.
270      * @param outputFrame frame in which the position/velocity coordinates shall be computed
271      * @return pvCoordinates in the specified output frame
272      * @exception OrekitException if transformation between frames cannot be computed
273      * @see #getPVCoordinates()
274      */
275     public TimeStampedPVCoordinates getPVCoordinates(final Frame outputFrame)
276         throws OrekitException {
277         if (pvCoordinates == null) {
278             pvCoordinates = initPVCoordinates();
279         }
280 
281         // If output frame requested is the same as definition frame,
282         // PV coordinates are returned directly
283         if (outputFrame == frame) {
284             return pvCoordinates;
285         }
286 
287         // Else, PV coordinates are transformed to output frame
288         final Transform t = frame.getTransformTo(outputFrame, date);
289         return t.transformPVCoordinates(pvCoordinates);
290     }
291 
292     /** {@inheritDoc} */
293     public TimeStampedPVCoordinates getPVCoordinates(final AbsoluteDate otherDate, final Frame otherFrame)
294         throws OrekitException {
295         return shiftedBy(otherDate.durationFrom(getDate())).getPVCoordinates(otherFrame);
296     }
297 
298 
299     /** Get the {@link TimeStampedPVCoordinates} in definition frame.
300      * @return pvCoordinates in the definition frame
301      * @see #getPVCoordinates(Frame)
302      */
303     public TimeStampedPVCoordinates getPVCoordinates() {
304         if (pvCoordinates == null) {
305             pvCoordinates = initPVCoordinates();
306         }
307         return pvCoordinates;
308     }
309 
310     /** Compute the position/velocity coordinates from the canonical parameters.
311      * @return computed position/velocity coordinates
312      */
313     protected abstract TimeStampedPVCoordinates initPVCoordinates();
314 
315     /** Get a time-shifted orbit.
316      * <p>
317      * The orbit can be slightly shifted to close dates. This shift is based on
318      * a simple keplerian model. It is <em>not</em> intended as a replacement
319      * for proper orbit and attitude propagation but should be sufficient for
320      * small time shifts or coarse accuracy.
321      * </p>
322      * @param dt time shift in seconds
323      * @return a new orbit, shifted with respect to the instance (which is immutable)
324      */
325     public abstract Orbit shiftedBy(final double dt);
326 
327     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
328      * <p>
329      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
330      * respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter
331      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
332      * </p>
333      * @param type type of the position angle to use
334      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
335      * is larger than 6x6, only the 6x6 upper left corner will be modified
336      */
337     public void getJacobianWrtCartesian(final PositionAngle type, final double[][] jacobian) {
338 
339         final double[][] cachedJacobian;
340         synchronized (this) {
341             switch (type) {
342             case MEAN :
343                 if (jacobianMeanWrtCartesian == null) {
344                     // first call, we need to compute the jacobian and cache it
345                     jacobianMeanWrtCartesian = computeJacobianMeanWrtCartesian();
346                 }
347                 cachedJacobian = jacobianMeanWrtCartesian;
348                 break;
349             case ECCENTRIC :
350                 if (jacobianEccentricWrtCartesian == null) {
351                     // first call, we need to compute the jacobian and cache it
352                     jacobianEccentricWrtCartesian = computeJacobianEccentricWrtCartesian();
353                 }
354                 cachedJacobian = jacobianEccentricWrtCartesian;
355                 break;
356             case TRUE :
357                 if (jacobianTrueWrtCartesian == null) {
358                     // first call, we need to compute the jacobian and cache it
359                     jacobianTrueWrtCartesian = computeJacobianTrueWrtCartesian();
360                 }
361                 cachedJacobian = jacobianTrueWrtCartesian;
362                 break;
363             default :
364                 throw OrekitException.createInternalError(null);
365             }
366         }
367 
368         // fill the user provided array
369         for (int i = 0; i < cachedJacobian.length; ++i) {
370             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
371         }
372 
373     }
374 
375     /** Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.
376      * <p>
377      * Element {@code jacobian[i][j]} is the derivative of Cartesian coordinate i of the orbit with
378      * respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate
379      * x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters.
380      * </p>
381      * @param type type of the position angle to use
382      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
383      * is larger than 6x6, only the 6x6 upper left corner will be modified
384      */
385     public void getJacobianWrtParameters(final PositionAngle type, final double[][] jacobian) {
386 
387         final double[][] cachedJacobian;
388         synchronized (this) {
389             switch (type) {
390             case MEAN :
391                 if (jacobianWrtParametersMean == null) {
392                     // first call, we need to compute the jacobian and cache it
393                     jacobianWrtParametersMean = createInverseJacobian(type);
394                 }
395                 cachedJacobian = jacobianWrtParametersMean;
396                 break;
397             case ECCENTRIC :
398                 if (jacobianWrtParametersEccentric == null) {
399                     // first call, we need to compute the jacobian and cache it
400                     jacobianWrtParametersEccentric = createInverseJacobian(type);
401                 }
402                 cachedJacobian = jacobianWrtParametersEccentric;
403                 break;
404             case TRUE :
405                 if (jacobianWrtParametersTrue == null) {
406                     // first call, we need to compute the jacobian and cache it
407                     jacobianWrtParametersTrue = createInverseJacobian(type);
408                 }
409                 cachedJacobian = jacobianWrtParametersTrue;
410                 break;
411             default :
412                 throw OrekitException.createInternalError(null);
413             }
414         }
415 
416         // fill the user-provided array
417         for (int i = 0; i < cachedJacobian.length; ++i) {
418             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
419         }
420 
421     }
422 
423     /** Create an inverse Jacobian.
424      * @param type type of the position angle to use
425      * @return inverse Jacobian
426      */
427     private double[][] createInverseJacobian(final PositionAngle type) {
428 
429         // get the direct Jacobian
430         final double[][] directJacobian = new double[6][6];
431         getJacobianWrtCartesian(type, directJacobian);
432 
433         // invert the direct Jacobian
434         final RealMatrix matrix = MatrixUtils.createRealMatrix(directJacobian);
435         final DecompositionSolver solver = new QRDecomposition(matrix).getSolver();
436         return solver.getInverse().getData();
437 
438     }
439 
440     /** Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
441      * <p>
442      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
443      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
444      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
445      * </p>
446      * @return 6x6 Jacobian matrix
447      * @see #computeJacobianEccentricWrtCartesian()
448      * @see #computeJacobianTrueWrtCartesian()
449      */
450     protected abstract double[][] computeJacobianMeanWrtCartesian();
451 
452     /** Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.
453      * <p>
454      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
455      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
456      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
457      * </p>
458      * @return 6x6 Jacobian matrix
459      * @see #computeJacobianMeanWrtCartesian()
460      * @see #computeJacobianTrueWrtCartesian()
461      */
462     protected abstract double[][] computeJacobianEccentricWrtCartesian();
463 
464     /** Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.
465      * <p>
466      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
467      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
468      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
469      * </p>
470      * @return 6x6 Jacobian matrix
471      * @see #computeJacobianMeanWrtCartesian()
472      * @see #computeJacobianEccentricWrtCartesian()
473      */
474     protected abstract double[][] computeJacobianTrueWrtCartesian();
475 
476     /** Add the contribution of the Keplerian motion to parameters derivatives
477      * <p>
478      * This method is used by integration-based propagators to evaluate the part of Keplerian
479      * motion to evolution of the orbital state.
480      * </p>
481      * @param type type of the position angle in the state
482      * @param gm attraction coefficient to use
483      * @param pDot array containing orbital state derivatives to update (the Keplerian
484      * part must be <em>added</em> to the array components, as the array may already
485      * contain some non-zero elements corresponding to non-Keplerian parts)
486      */
487     public abstract void addKeplerContribution(final PositionAngle type, final double gm, double[] pDot);
488 
489         /** Fill a Jacobian half row with a single vector.
490      * @param a coefficient of the vector
491      * @param v vector
492      * @param row Jacobian matrix row
493      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
494      */
495     protected static void fillHalfRow(final double a, final Vector3D v, final double[] row, final int j) {
496         row[j]     = a * v.getX();
497         row[j + 1] = a * v.getY();
498         row[j + 2] = a * v.getZ();
499     }
500 
501     /** Fill a Jacobian half row with a linear combination of vectors.
502      * @param a1 coefficient of the first vector
503      * @param v1 first vector
504      * @param a2 coefficient of the second vector
505      * @param v2 second vector
506      * @param row Jacobian matrix row
507      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
508      */
509     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
510                                       final double[] row, final int j) {
511         row[j]     = a1 * v1.getX() + a2 * v2.getX();
512         row[j + 1] = a1 * v1.getY() + a2 * v2.getY();
513         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ();
514     }
515 
516     /** Fill a Jacobian half row with a linear combination of vectors.
517      * @param a1 coefficient of the first vector
518      * @param v1 first vector
519      * @param a2 coefficient of the second vector
520      * @param v2 second vector
521      * @param a3 coefficient of the third vector
522      * @param v3 third vector
523      * @param row Jacobian matrix row
524      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
525      */
526     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
527                                       final double a3, final Vector3D v3,
528                                       final double[] row, final int j) {
529         row[j]     = a1 * v1.getX() + a2 * v2.getX() + a3 * v3.getX();
530         row[j + 1] = a1 * v1.getY() + a2 * v2.getY() + a3 * v3.getY();
531         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ() + a3 * v3.getZ();
532     }
533 
534     /** Fill a Jacobian half row with a linear combination of vectors.
535      * @param a1 coefficient of the first vector
536      * @param v1 first vector
537      * @param a2 coefficient of the second vector
538      * @param v2 second vector
539      * @param a3 coefficient of the third vector
540      * @param v3 third vector
541      * @param a4 coefficient of the fourth vector
542      * @param v4 fourth vector
543      * @param row Jacobian matrix row
544      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
545      */
546     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
547                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
548                                       final double[] row, final int j) {
549         row[j]     = a1 * v1.getX() + a2 * v2.getX() + a3 * v3.getX() + a4 * v4.getX();
550         row[j + 1] = a1 * v1.getY() + a2 * v2.getY() + a3 * v3.getY() + a4 * v4.getY();
551         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ() + a3 * v3.getZ() + a4 * v4.getZ();
552     }
553 
554     /** Fill a Jacobian half row with a linear combination of vectors.
555      * @param a1 coefficient of the first vector
556      * @param v1 first vector
557      * @param a2 coefficient of the second vector
558      * @param v2 second vector
559      * @param a3 coefficient of the third vector
560      * @param v3 third vector
561      * @param a4 coefficient of the fourth vector
562      * @param v4 fourth vector
563      * @param a5 coefficient of the fifth vector
564      * @param v5 fifth vector
565      * @param row Jacobian matrix row
566      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
567      */
568     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
569                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
570                                       final double a5, final Vector3D v5,
571                                       final double[] row, final int j) {
572         row[j]     = a1 * v1.getX() + a2 * v2.getX() + a3 * v3.getX() + a4 * v4.getX() + a5 * v5.getX();
573         row[j + 1] = a1 * v1.getY() + a2 * v2.getY() + a3 * v3.getY() + a4 * v4.getY() + a5 * v5.getY();
574         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ() + a3 * v3.getZ() + a4 * v4.getZ() + a5 * v5.getZ();
575     }
576 
577     /** Fill a Jacobian half row with a linear combination of vectors.
578      * @param a1 coefficient of the first vector
579      * @param v1 first vector
580      * @param a2 coefficient of the second vector
581      * @param v2 second vector
582      * @param a3 coefficient of the third vector
583      * @param v3 third vector
584      * @param a4 coefficient of the fourth vector
585      * @param v4 fourth vector
586      * @param a5 coefficient of the fifth vector
587      * @param v5 fifth vector
588      * @param a6 coefficient of the sixth vector
589      * @param v6 sixth vector
590      * @param row Jacobian matrix row
591      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
592      */
593     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
594                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
595                                       final double a5, final Vector3D v5, final double a6, final Vector3D v6,
596                                       final double[] row, final int j) {
597         row[j]     = a1 * v1.getX() + a2 * v2.getX() + a3 * v3.getX() + a4 * v4.getX() + a5 * v5.getX() + a6 * v6.getX();
598         row[j + 1] = a1 * v1.getY() + a2 * v2.getY() + a3 * v3.getY() + a4 * v4.getY() + a5 * v5.getY() + a6 * v6.getY();
599         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ() + a3 * v3.getZ() + a4 * v4.getZ() + a5 * v5.getZ() + a6 * v6.getZ();
600     }
601 
602 }