Orbit.java

  1. /* Copyright 2002-2016 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. import java.io.Serializable;

  19. import org.apache.commons.math3.geometry.euclidean.threed.Vector3D;
  20. import org.apache.commons.math3.linear.DecompositionSolver;
  21. import org.apache.commons.math3.linear.MatrixUtils;
  22. import org.apache.commons.math3.linear.QRDecomposition;
  23. import org.apache.commons.math3.linear.RealMatrix;
  24. import org.apache.commons.math3.util.FastMath;
  25. import org.orekit.errors.OrekitIllegalArgumentException;
  26. import org.orekit.errors.OrekitException;
  27. import org.orekit.errors.OrekitInternalError;
  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.  * This class handles orbital parameters.

  39.  * <p>
  40.  * For user convenience, both the Cartesian and the equinoctial elements
  41.  * are provided by this class, regardless of the canonical representation
  42.  * implemented in the derived class (which may be classical keplerian
  43.  * elements for example).
  44.  * </p>
  45.  * <p>
  46.  * The parameters are defined in a frame specified by the user. It is important
  47.  * to make sure this frame is consistent: it probably is inertial and centered
  48.  * on the central body. This information is used for example by some
  49.  * force models.
  50.  * </p>
  51.  * <p>
  52.  * The object <code>OrbitalParameters</code> is guaranteed to be immutable.
  53.  * </p>
  54.  * @author Luc Maisonobe
  55.  * @author Guylaine Prat
  56.  * @author Fabien Maussion
  57.  * @author V&eacute;ronique Pommier-Maurussane
  58.  */
  59. public abstract class Orbit
  60.     implements TimeStamped, TimeShiftable<Orbit>, TimeInterpolable<Orbit>,
  61.                Serializable, PVCoordinatesProvider {

  62.     /** Serializable UID. */
  63.     private static final long serialVersionUID = 438733454597999578L;

  64.     /** Frame in which are defined the orbital parameters. */
  65.     private final Frame frame;

  66.     /** Date of the orbital parameters. */
  67.     private final AbsoluteDate date;

  68.     /** Value of mu used to compute position and velocity (m³/s²). */
  69.     private final double mu;

  70.     /** Computed PVCoordinates. */
  71.     private transient TimeStampedPVCoordinates pvCoordinates;

  72.     /** Jacobian of the orbital parameters with mean angle with respect to the Cartesian coordinates. */
  73.     private transient double[][] jacobianMeanWrtCartesian;

  74.     /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with mean angle. */
  75.     private transient double[][] jacobianWrtParametersMean;

  76.     /** Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian coordinates. */
  77.     private transient double[][] jacobianEccentricWrtCartesian;

  78.     /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with eccentric angle. */
  79.     private transient double[][] jacobianWrtParametersEccentric;

  80.     /** Jacobian of the orbital parameters with true angle with respect to the Cartesian coordinates. */
  81.     private transient double[][] jacobianTrueWrtCartesian;

  82.     /** Jacobian of the Cartesian coordinates with respect to the orbital parameters with true angle. */
  83.     private transient double[][] jacobianWrtParametersTrue;

  84.     /** Default constructor.
  85.      * Build a new instance with arbitrary default elements.
  86.      * @param frame the frame in which the parameters are defined
  87.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  88.      * @param date date of the orbital parameters
  89.      * @param mu central attraction coefficient (m^3/s^2)
  90.      * @exception IllegalArgumentException if frame is not a {@link
  91.      * Frame#isPseudoInertial pseudo-inertial frame}
  92.      */
  93.     protected Orbit(final Frame frame, final AbsoluteDate date, final double mu)
  94.         throws IllegalArgumentException {
  95.         ensurePseudoInertialFrame(frame);
  96.         this.date                      = date;
  97.         this.mu                        = mu;
  98.         this.pvCoordinates             = null;
  99.         this.frame                     = frame;
  100.         jacobianMeanWrtCartesian       = null;
  101.         jacobianWrtParametersMean      = null;
  102.         jacobianEccentricWrtCartesian  = null;
  103.         jacobianWrtParametersEccentric = null;
  104.         jacobianTrueWrtCartesian       = null;
  105.         jacobianWrtParametersTrue      = null;
  106.     }

  107.     /** Set the orbit from Cartesian parameters.
  108.      *
  109.      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
  110.      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
  111.      * use {@code mu} and the position to compute the acceleration, including
  112.      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
  113.      *
  114.      * @param pvCoordinates the position and velocity in the inertial frame
  115.      * @param frame the frame in which the {@link TimeStampedPVCoordinates} are defined
  116.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  117.      * @param mu central attraction coefficient (m^3/s^2)
  118.      * @exception IllegalArgumentException if frame is not a {@link
  119.      * Frame#isPseudoInertial pseudo-inertial frame}
  120.      */
  121.     protected Orbit(final TimeStampedPVCoordinates pvCoordinates, final Frame frame, final double mu)
  122.         throws IllegalArgumentException {
  123.         ensurePseudoInertialFrame(frame);
  124.         this.date = pvCoordinates.getDate();
  125.         this.mu = mu;
  126.         if (pvCoordinates.getAcceleration().getNormSq() == 0) {
  127.             // the acceleration was not provided,
  128.             // compute it from Newtonian attraction
  129.             final double r2 = pvCoordinates.getPosition().getNormSq();
  130.             final double r3 = r2 * FastMath.sqrt(r2);
  131.             this.pvCoordinates = new TimeStampedPVCoordinates(pvCoordinates.getDate(),
  132.                                                               pvCoordinates.getPosition(),
  133.                                                               pvCoordinates.getVelocity(),
  134.                                                               new Vector3D(-mu / r3, pvCoordinates.getPosition()));
  135.         } else {
  136.             this.pvCoordinates = pvCoordinates;
  137.         }
  138.         this.frame = frame;
  139.     }

  140.     /** Get the orbit type.
  141.      * @return orbit type
  142.      */
  143.     public abstract OrbitType getType();

  144.     /** Ensure the defining frame is a pseudo-inertial frame.
  145.      * @param frame frame to check
  146.      * @exception IllegalArgumentException if frame is not a {@link
  147.      * Frame#isPseudoInertial pseudo-inertial frame}
  148.      */
  149.     private static void ensurePseudoInertialFrame(final Frame frame)
  150.         throws IllegalArgumentException {
  151.         if (!frame.isPseudoInertial()) {
  152.             throw new OrekitIllegalArgumentException(OrekitMessages.NON_PSEUDO_INERTIAL_FRAME,
  153.                                                      frame.getName());
  154.         }
  155.     }

  156.     /** Get the frame in which the orbital parameters are defined.
  157.      * @return frame in which the orbital parameters are defined
  158.      */
  159.     public Frame getFrame() {
  160.         return frame;
  161.     }

  162.     /** Get the semi-major axis.
  163.      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
  164.      * @return semi-major axis (m)
  165.      */
  166.     public abstract double getA();

  167.     /** Get the first component of the equinoctial eccentricity vector.
  168.      * @return first component of the equinoctial eccentricity vector
  169.      */
  170.     public abstract double getEquinoctialEx();

  171.     /** Get the second component of the equinoctial eccentricity vector.
  172.      * @return second component of the equinoctial eccentricity vector
  173.      */
  174.     public abstract double getEquinoctialEy();

  175.     /** Get the first component of the inclination vector.
  176.      * @return first component of the inclination vector
  177.      */
  178.     public abstract double getHx();

  179.     /** Get the second component of the inclination vector.
  180.      * @return second component of the inclination vector
  181.      */
  182.     public abstract double getHy();

  183.     /** Get the eccentric longitude argument.
  184.      * @return E + ω + Ω eccentric longitude argument (rad)
  185.      */
  186.     public abstract double getLE();

  187.     /** Get the true longitude argument.
  188.      * @return v + ω + Ω true longitude argument (rad)
  189.      */
  190.     public abstract double getLv();

  191.     /** Get the mean longitude argument.
  192.      * @return M + ω + Ω mean longitude argument (rad)
  193.      */
  194.     public abstract double getLM();

  195.     // Additional orbital elements

  196.     /** Get the eccentricity.
  197.      * @return eccentricity
  198.      */
  199.     public abstract double getE();

  200.     /** Get the inclination.
  201.      * @return inclination (rad)
  202.      */
  203.     public abstract double getI();

  204.     /** Get the central acceleration constant.
  205.      * @return central acceleration constant
  206.      */
  207.     public double getMu() {
  208.         return mu;
  209.     }

  210.     /** Get the keplerian period.
  211.      * <p>The keplerian period is computed directly from semi major axis
  212.      * and central acceleration constant.</p>
  213.      * @return keplerian period in seconds, or positive infinity for hyperbolic orbits
  214.      */
  215.     public double getKeplerianPeriod() {
  216.         final double a = getA();
  217.         return (a < 0) ? Double.POSITIVE_INFINITY : 2.0 * FastMath.PI * a * FastMath.sqrt(a / mu);
  218.     }

  219.     /** Get the keplerian mean motion.
  220.      * <p>The keplerian mean motion is computed directly from semi major axis
  221.      * and central acceleration constant.</p>
  222.      * @return keplerian mean motion in radians per second
  223.      */
  224.     public double getKeplerianMeanMotion() {
  225.         final double absA = FastMath.abs(getA());
  226.         return FastMath.sqrt(mu / absA) / absA;
  227.     }

  228.     /** Get the date of orbital parameters.
  229.      * @return date of the orbital parameters
  230.      */
  231.     public AbsoluteDate getDate() {
  232.         return date;
  233.     }

  234.     /** Get the {@link TimeStampedPVCoordinates} in a specified frame.
  235.      * @param outputFrame frame in which the position/velocity coordinates shall be computed
  236.      * @return pvCoordinates in the specified output frame
  237.      * @exception OrekitException if transformation between frames cannot be computed
  238.      * @see #getPVCoordinates()
  239.      */
  240.     public TimeStampedPVCoordinates getPVCoordinates(final Frame outputFrame)
  241.         throws OrekitException {
  242.         if (pvCoordinates == null) {
  243.             pvCoordinates = initPVCoordinates();
  244.         }

  245.         // If output frame requested is the same as definition frame,
  246.         // PV coordinates are returned directly
  247.         if (outputFrame == frame) {
  248.             return pvCoordinates;
  249.         }

  250.         // Else, PV coordinates are transformed to output frame
  251.         final Transform t = frame.getTransformTo(outputFrame, date);
  252.         return t.transformPVCoordinates(pvCoordinates);
  253.     }

  254.     /** {@inheritDoc} */
  255.     public TimeStampedPVCoordinates getPVCoordinates(final AbsoluteDate otherDate, final Frame otherFrame)
  256.         throws OrekitException {
  257.         return shiftedBy(otherDate.durationFrom(getDate())).getPVCoordinates(otherFrame);
  258.     }


  259.     /** Get the {@link TimeStampedPVCoordinates} in definition frame.
  260.      * @return pvCoordinates in the definition frame
  261.      * @see #getPVCoordinates(Frame)
  262.      */
  263.     public TimeStampedPVCoordinates getPVCoordinates() {
  264.         if (pvCoordinates == null) {
  265.             pvCoordinates = initPVCoordinates();
  266.         }
  267.         return pvCoordinates;
  268.     }

  269.     /** Compute the position/velocity coordinates from the canonical parameters.
  270.      * @return computed position/velocity coordinates
  271.      */
  272.     protected abstract TimeStampedPVCoordinates initPVCoordinates();

  273.     /** Get a time-shifted orbit.
  274.      * <p>
  275.      * The orbit can be slightly shifted to close dates. This shift is based on
  276.      * a simple keplerian model. It is <em>not</em> intended as a replacement
  277.      * for proper orbit and attitude propagation but should be sufficient for
  278.      * small time shifts or coarse accuracy.
  279.      * </p>
  280.      * @param dt time shift in seconds
  281.      * @return a new orbit, shifted with respect to the instance (which is immutable)
  282.      */
  283.     public abstract Orbit shiftedBy(final double dt);

  284.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  285.      * <p>
  286.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  287.      * respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter
  288.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  289.      * </p>
  290.      * @param type type of the position angle to use
  291.      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
  292.      * is larger than 6x6, only the 6x6 upper left corner will be modified
  293.      */
  294.     public void getJacobianWrtCartesian(final PositionAngle type, final double[][] jacobian) {

  295.         final double[][] cachedJacobian;
  296.         synchronized (this) {
  297.             switch (type) {
  298.                 case MEAN :
  299.                     if (jacobianMeanWrtCartesian == null) {
  300.                         // first call, we need to compute the jacobian and cache it
  301.                         jacobianMeanWrtCartesian = computeJacobianMeanWrtCartesian();
  302.                     }
  303.                     cachedJacobian = jacobianMeanWrtCartesian;
  304.                     break;
  305.                 case ECCENTRIC :
  306.                     if (jacobianEccentricWrtCartesian == null) {
  307.                         // first call, we need to compute the jacobian and cache it
  308.                         jacobianEccentricWrtCartesian = computeJacobianEccentricWrtCartesian();
  309.                     }
  310.                     cachedJacobian = jacobianEccentricWrtCartesian;
  311.                     break;
  312.                 case TRUE :
  313.                     if (jacobianTrueWrtCartesian == null) {
  314.                         // first call, we need to compute the jacobian and cache it
  315.                         jacobianTrueWrtCartesian = computeJacobianTrueWrtCartesian();
  316.                     }
  317.                     cachedJacobian = jacobianTrueWrtCartesian;
  318.                     break;
  319.                 default :
  320.                     throw new OrekitInternalError(null);
  321.             }
  322.         }

  323.         // fill the user provided array
  324.         for (int i = 0; i < cachedJacobian.length; ++i) {
  325.             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
  326.         }

  327.     }

  328.     /** Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.
  329.      * <p>
  330.      * Element {@code jacobian[i][j]} is the derivative of Cartesian coordinate i of the orbit with
  331.      * respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate
  332.      * x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters.
  333.      * </p>
  334.      * @param type type of the position angle to use
  335.      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
  336.      * is larger than 6x6, only the 6x6 upper left corner will be modified
  337.      */
  338.     public void getJacobianWrtParameters(final PositionAngle type, final double[][] jacobian) {

  339.         final double[][] cachedJacobian;
  340.         synchronized (this) {
  341.             switch (type) {
  342.                 case MEAN :
  343.                     if (jacobianWrtParametersMean == null) {
  344.                         // first call, we need to compute the jacobian and cache it
  345.                         jacobianWrtParametersMean = createInverseJacobian(type);
  346.                     }
  347.                     cachedJacobian = jacobianWrtParametersMean;
  348.                     break;
  349.                 case ECCENTRIC :
  350.                     if (jacobianWrtParametersEccentric == null) {
  351.                         // first call, we need to compute the jacobian and cache it
  352.                         jacobianWrtParametersEccentric = createInverseJacobian(type);
  353.                     }
  354.                     cachedJacobian = jacobianWrtParametersEccentric;
  355.                     break;
  356.                 case TRUE :
  357.                     if (jacobianWrtParametersTrue == null) {
  358.                         // first call, we need to compute the jacobian and cache it
  359.                         jacobianWrtParametersTrue = createInverseJacobian(type);
  360.                     }
  361.                     cachedJacobian = jacobianWrtParametersTrue;
  362.                     break;
  363.                 default :
  364.                     throw new OrekitInternalError(null);
  365.             }
  366.         }

  367.         // fill the user-provided array
  368.         for (int i = 0; i < cachedJacobian.length; ++i) {
  369.             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
  370.         }

  371.     }

  372.     /** Create an inverse Jacobian.
  373.      * @param type type of the position angle to use
  374.      * @return inverse Jacobian
  375.      */
  376.     private double[][] createInverseJacobian(final PositionAngle type) {

  377.         // get the direct Jacobian
  378.         final double[][] directJacobian = new double[6][6];
  379.         getJacobianWrtCartesian(type, directJacobian);

  380.         // invert the direct Jacobian
  381.         final RealMatrix matrix = MatrixUtils.createRealMatrix(directJacobian);
  382.         final DecompositionSolver solver = new QRDecomposition(matrix).getSolver();
  383.         return solver.getInverse().getData();

  384.     }

  385.     /** Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
  386.      * <p>
  387.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  388.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  389.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  390.      * </p>
  391.      * @return 6x6 Jacobian matrix
  392.      * @see #computeJacobianEccentricWrtCartesian()
  393.      * @see #computeJacobianTrueWrtCartesian()
  394.      */
  395.     protected abstract double[][] computeJacobianMeanWrtCartesian();

  396.     /** Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.
  397.      * <p>
  398.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  399.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  400.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  401.      * </p>
  402.      * @return 6x6 Jacobian matrix
  403.      * @see #computeJacobianMeanWrtCartesian()
  404.      * @see #computeJacobianTrueWrtCartesian()
  405.      */
  406.     protected abstract double[][] computeJacobianEccentricWrtCartesian();

  407.     /** Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.
  408.      * <p>
  409.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  410.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  411.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  412.      * </p>
  413.      * @return 6x6 Jacobian matrix
  414.      * @see #computeJacobianMeanWrtCartesian()
  415.      * @see #computeJacobianEccentricWrtCartesian()
  416.      */
  417.     protected abstract double[][] computeJacobianTrueWrtCartesian();

  418.     /** Add the contribution of the Keplerian motion to parameters derivatives
  419.      * <p>
  420.      * This method is used by integration-based propagators to evaluate the part of Keplerian
  421.      * motion to evolution of the orbital state.
  422.      * </p>
  423.      * @param type type of the position angle in the state
  424.      * @param gm attraction coefficient to use
  425.      * @param pDot array containing orbital state derivatives to update (the Keplerian
  426.      * part must be <em>added</em> to the array components, as the array may already
  427.      * contain some non-zero elements corresponding to non-Keplerian parts)
  428.      */
  429.     public abstract void addKeplerContribution(final PositionAngle type, final double gm, double[] pDot);

  430.         /** Fill a Jacobian half row with a single vector.
  431.      * @param a coefficient of the vector
  432.      * @param v vector
  433.      * @param row Jacobian matrix row
  434.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  435.      */
  436.     protected static void fillHalfRow(final double a, final Vector3D v, final double[] row, final int j) {
  437.         row[j]     = a * v.getX();
  438.         row[j + 1] = a * v.getY();
  439.         row[j + 2] = a * v.getZ();
  440.     }

  441.     /** Fill a Jacobian half row with a linear combination of vectors.
  442.      * @param a1 coefficient of the first vector
  443.      * @param v1 first vector
  444.      * @param a2 coefficient of the second vector
  445.      * @param v2 second vector
  446.      * @param row Jacobian matrix row
  447.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  448.      */
  449.     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
  450.                                       final double[] row, final int j) {
  451.         row[j]     = a1 * v1.getX() + a2 * v2.getX();
  452.         row[j + 1] = a1 * v1.getY() + a2 * v2.getY();
  453.         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ();
  454.     }

  455.     /** Fill a Jacobian half row with a linear combination of vectors.
  456.      * @param a1 coefficient of the first vector
  457.      * @param v1 first vector
  458.      * @param a2 coefficient of the second vector
  459.      * @param v2 second vector
  460.      * @param a3 coefficient of the third vector
  461.      * @param v3 third vector
  462.      * @param row Jacobian matrix row
  463.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  464.      */
  465.     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
  466.                                       final double a3, final Vector3D v3,
  467.                                       final double[] row, final int j) {
  468.         row[j]     = a1 * v1.getX() + a2 * v2.getX() + a3 * v3.getX();
  469.         row[j + 1] = a1 * v1.getY() + a2 * v2.getY() + a3 * v3.getY();
  470.         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ() + a3 * v3.getZ();
  471.     }

  472.     /** Fill a Jacobian half row with a linear combination of vectors.
  473.      * @param a1 coefficient of the first vector
  474.      * @param v1 first vector
  475.      * @param a2 coefficient of the second vector
  476.      * @param v2 second vector
  477.      * @param a3 coefficient of the third vector
  478.      * @param v3 third vector
  479.      * @param a4 coefficient of the fourth vector
  480.      * @param v4 fourth vector
  481.      * @param row Jacobian matrix row
  482.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  483.      */
  484.     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
  485.                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
  486.                                       final double[] row, final int j) {
  487.         row[j]     = a1 * v1.getX() + a2 * v2.getX() + a3 * v3.getX() + a4 * v4.getX();
  488.         row[j + 1] = a1 * v1.getY() + a2 * v2.getY() + a3 * v3.getY() + a4 * v4.getY();
  489.         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ() + a3 * v3.getZ() + a4 * v4.getZ();
  490.     }

  491.     /** Fill a Jacobian half row with a linear combination of vectors.
  492.      * @param a1 coefficient of the first vector
  493.      * @param v1 first vector
  494.      * @param a2 coefficient of the second vector
  495.      * @param v2 second vector
  496.      * @param a3 coefficient of the third vector
  497.      * @param v3 third vector
  498.      * @param a4 coefficient of the fourth vector
  499.      * @param v4 fourth vector
  500.      * @param a5 coefficient of the fifth vector
  501.      * @param v5 fifth vector
  502.      * @param row Jacobian matrix row
  503.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  504.      */
  505.     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
  506.                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
  507.                                       final double a5, final Vector3D v5,
  508.                                       final double[] row, final int j) {
  509.         row[j]     = a1 * v1.getX() + a2 * v2.getX() + a3 * v3.getX() + a4 * v4.getX() + a5 * v5.getX();
  510.         row[j + 1] = a1 * v1.getY() + a2 * v2.getY() + a3 * v3.getY() + a4 * v4.getY() + a5 * v5.getY();
  511.         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ() + a3 * v3.getZ() + a4 * v4.getZ() + a5 * v5.getZ();
  512.     }

  513.     /** Fill a Jacobian half row with a linear combination of vectors.
  514.      * @param a1 coefficient of the first vector
  515.      * @param v1 first vector
  516.      * @param a2 coefficient of the second vector
  517.      * @param v2 second vector
  518.      * @param a3 coefficient of the third vector
  519.      * @param v3 third vector
  520.      * @param a4 coefficient of the fourth vector
  521.      * @param v4 fourth vector
  522.      * @param a5 coefficient of the fifth vector
  523.      * @param v5 fifth vector
  524.      * @param a6 coefficient of the sixth vector
  525.      * @param v6 sixth vector
  526.      * @param row Jacobian matrix row
  527.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  528.      */
  529.     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
  530.                                       final double a3, final Vector3D v3, final double a4, final Vector3D v4,
  531.                                       final double a5, final Vector3D v5, final double a6, final Vector3D v6,
  532.                                       final double[] row, final int j) {
  533.         row[j]     = a1 * v1.getX() + a2 * v2.getX() + a3 * v3.getX() + a4 * v4.getX() + a5 * v5.getX() + a6 * v6.getX();
  534.         row[j + 1] = a1 * v1.getY() + a2 * v2.getY() + a3 * v3.getY() + a4 * v4.getY() + a5 * v5.getY() + a6 * v6.getY();
  535.         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ() + a3 * v3.getZ() + a4 * v4.getZ() + a5 * v5.getZ() + a6 * v6.getZ();
  536.     }

  537. }