Orbit.java

  1. /* Copyright 2002-2013 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.OrekitException;
  26. import org.orekit.errors.OrekitMessages;
  27. import org.orekit.frames.Frame;
  28. import org.orekit.frames.Transform;
  29. import org.orekit.time.AbsoluteDate;
  30. import org.orekit.time.TimeInterpolable;
  31. import org.orekit.time.TimeShiftable;
  32. import org.orekit.time.TimeStamped;
  33. import org.orekit.utils.PVCoordinates;
  34. import org.orekit.utils.PVCoordinatesProvider;

  35. /**
  36.  * This class handles orbital parameters.

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

  60.     /** Serializable UID. */
  61.     private static final long serialVersionUID = 438733454597999578L;

  62.     /** Frame in which are defined the orbital parameters. */
  63.     private final Frame frame;

  64.     /** Date of the orbital parameters. */
  65.     private final AbsoluteDate date;

  66.     /** Value of mu used to compute position and velocity (m<sup>3</sup>/s<sup>2</sup>). */
  67.     private final double mu;

  68.     /** Computed PVCoordinates. */
  69.     private transient PVCoordinates pvCoordinates;

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

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

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

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

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

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

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

  105.     /** Set the orbit from Cartesian parameters.
  106.      * @param pvCoordinates the position and velocity in the inertial frame
  107.      * @param frame the frame in which the {@link PVCoordinates} are defined
  108.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  109.      * @param date date of the orbital parameters
  110.      * @param mu central attraction coefficient (m^3/s^2)
  111.      * @exception IllegalArgumentException if frame is not a {@link
  112.      * Frame#isPseudoInertial pseudo-inertial frame}
  113.      */
  114.     protected Orbit(final PVCoordinates pvCoordinates, final Frame frame,
  115.                     final AbsoluteDate date, final double mu)
  116.         throws IllegalArgumentException {
  117.         ensurePseudoInertialFrame(frame);
  118.         this.date = date;
  119.         this.mu = mu;
  120.         this.pvCoordinates = pvCoordinates;
  121.         this.frame = frame;
  122.     }

  123.     /** Get the orbit type.
  124.      * @return orbit type
  125.      */
  126.     public abstract OrbitType getType();

  127.     /** Ensure the defining frame is a pseudo-inertial frame.
  128.      * @param frame frame to check
  129.      * @exception IllegalArgumentException if frame is not a {@link
  130.      * Frame#isPseudoInertial pseudo-inertial frame}
  131.      */
  132.     private static void ensurePseudoInertialFrame(final Frame frame)
  133.         throws IllegalArgumentException {
  134.         if (!frame.isPseudoInertial()) {
  135.             throw OrekitException.createIllegalArgumentException(
  136.                 OrekitMessages.NON_PSEUDO_INERTIAL_FRAME_NOT_SUITABLE_FOR_DEFINING_ORBITS,
  137.                 frame.getName());
  138.         }
  139.     }

  140.     /** Get the frame in which the orbital parameters are defined.
  141.      * @return frame in which the orbital parameters are defined
  142.      */
  143.     public Frame getFrame() {
  144.         return frame;
  145.     }

  146.     /** Get the semi-major axis.
  147.      * <p>Note that the semi-major axis is considered negative for hyperbolic orbits.</p>
  148.      * @return semi-major axis (m)
  149.      */
  150.     public abstract double getA();

  151.     /** Get the first component of the equinoctial eccentricity vector.
  152.      * @return first component of the equinoctial eccentricity vector
  153.      */
  154.     public abstract double getEquinoctialEx();

  155.     /** Get the second component of the equinoctial eccentricity vector.
  156.      * @return second component of the equinoctial eccentricity vector
  157.      */
  158.     public abstract double getEquinoctialEy();

  159.     /** Get the first component of the inclination vector.
  160.      * @return first component of the inclination vector
  161.      */
  162.     public abstract double getHx();

  163.     /** Get the second component of the inclination vector.
  164.      * @return second component of the inclination vector
  165.      */
  166.     public abstract double getHy();

  167.     /** Get the eccentric longitude argument.
  168.      * @return E + &omega; + &Omega; eccentric longitude argument (rad)
  169.      */
  170.     public abstract double getLE();

  171.     /** Get the true longitude argument.
  172.      * @return v + &omega; + &Omega; true longitude argument (rad)
  173.      */
  174.     public abstract double getLv();

  175.     /** Get the mean longitude argument.
  176.      * @return M + &omega; + &Omega; mean longitude argument (rad)
  177.      */
  178.     public abstract double getLM();

  179.     // Additional orbital elements

  180.     /** Get the eccentricity.
  181.      * @return eccentricity
  182.      */
  183.     public abstract double getE();

  184.     /** Get the inclination.
  185.      * @return inclination (rad)
  186.      */
  187.     public abstract double getI();

  188.     /** Get the central acceleration constant.
  189.      * @return central acceleration constant
  190.      */
  191.     public double getMu() {
  192.         return mu;
  193.     }

  194.     /** Get the keplerian period.
  195.      * <p>The keplerian period is computed directly from semi major axis
  196.      * and central acceleration constant.</p>
  197.      * @return keplerian period in seconds, or positive infinity for hyperbolic orbits
  198.      */
  199.     public double getKeplerianPeriod() {
  200.         final double a = getA();
  201.         return (a < 0) ? Double.POSITIVE_INFINITY : 2.0 * FastMath.PI * a * FastMath.sqrt(a / mu);
  202.     }

  203.     /** Get the keplerian mean motion.
  204.      * <p>The keplerian mean motion is computed directly from semi major axis
  205.      * and central acceleration constant.</p>
  206.      * @return keplerian mean motion in radians per second
  207.      */
  208.     public double getKeplerianMeanMotion() {
  209.         final double absA = FastMath.abs(getA());
  210.         return FastMath.sqrt(mu / absA) / absA;
  211.     }

  212.     /** Get the date of orbital parameters.
  213.      * @return date of the orbital parameters
  214.      */
  215.     public AbsoluteDate getDate() {
  216.         return date;
  217.     }

  218.     /** Get the {@link PVCoordinates} in a specified frame.
  219.      * @param outputFrame frame in which the position/velocity coordinates shall be computed
  220.      * @return pvCoordinates in the specified output frame
  221.      * @exception OrekitException if transformation between frames cannot be computed
  222.      * @see #getPVCoordinates()
  223.      */
  224.     public PVCoordinates getPVCoordinates(final Frame outputFrame)
  225.         throws OrekitException {
  226.         if (pvCoordinates == null) {
  227.             pvCoordinates = initPVCoordinates();
  228.         }

  229.         // If output frame requested is the same as definition frame,
  230.         // PV coordinates are returned directly
  231.         if (outputFrame == frame) {
  232.             return pvCoordinates;
  233.         }

  234.         // Else, PV coordinates are transformed to output frame
  235.         final Transform t = frame.getTransformTo(outputFrame, date);
  236.         return t.transformPVCoordinates(pvCoordinates);
  237.     }

  238.     /** {@inheritDoc} */
  239.     public PVCoordinates getPVCoordinates(final AbsoluteDate otherDate, final Frame otherFrame)
  240.         throws OrekitException {
  241.         return shiftedBy(otherDate.durationFrom(getDate())).getPVCoordinates(otherFrame);
  242.     }


  243.     /** Get the {@link PVCoordinates} in definition frame.
  244.      * @return pvCoordinates in the definition frame
  245.      * @see #getPVCoordinates(Frame)
  246.      */
  247.     public PVCoordinates getPVCoordinates() {
  248.         if (pvCoordinates == null) {
  249.             pvCoordinates = initPVCoordinates();
  250.         }
  251.         return pvCoordinates;
  252.     }

  253.     /** Compute the position/velocity coordinates from the canonical parameters.
  254.      * @return computed position/velocity coordinates
  255.      */
  256.     protected abstract PVCoordinates initPVCoordinates();

  257.     /** Get a time-shifted orbit.
  258.      * <p>
  259.      * The orbit can be slightly shifted to close dates. This shift is based on
  260.      * a simple keplerian model. It is <em>not</em> intended as a replacement
  261.      * for proper orbit and attitude propagation but should be sufficient for
  262.      * small time shifts or coarse accuracy.
  263.      * </p>
  264.      * @param dt time shift in seconds
  265.      * @return a new orbit, shifted with respect to the instance (which is immutable)
  266.      */
  267.     public abstract Orbit shiftedBy(final double dt);

  268.     /** Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
  269.      * <p>
  270.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  271.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  272.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  273.      * </p>
  274.      * @param type type of the position angle to use
  275.      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
  276.      * is larger than 6x6, only the 6x6 upper left corner will be modified
  277.      */
  278.     public void getJacobianWrtCartesian(final PositionAngle type, final double[][] jacobian) {

  279.         final double[][] cachedJacobian;
  280.         synchronized (this) {
  281.             switch (type) {
  282.             case MEAN :
  283.                 if (jacobianMeanWrtCartesian == null) {
  284.                     // first call, we need to compute the jacobian and cache it
  285.                     jacobianMeanWrtCartesian = computeJacobianMeanWrtCartesian();
  286.                 }
  287.                 cachedJacobian = jacobianMeanWrtCartesian;
  288.                 break;
  289.             case ECCENTRIC :
  290.                 if (jacobianEccentricWrtCartesian == null) {
  291.                     // first call, we need to compute the jacobian and cache it
  292.                     jacobianEccentricWrtCartesian = computeJacobianEccentricWrtCartesian();
  293.                 }
  294.                 cachedJacobian = jacobianEccentricWrtCartesian;
  295.                 break;
  296.             case TRUE :
  297.                 if (jacobianTrueWrtCartesian == null) {
  298.                     // first call, we need to compute the jacobian and cache it
  299.                     jacobianTrueWrtCartesian = computeJacobianTrueWrtCartesian();
  300.                 }
  301.                 cachedJacobian = jacobianTrueWrtCartesian;
  302.                 break;
  303.             default :
  304.                 throw OrekitException.createInternalError(null);
  305.             }
  306.         }

  307.         // fill the user provided array
  308.         for (int i = 0; i < cachedJacobian.length; ++i) {
  309.             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
  310.         }

  311.     }

  312.     /** Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.
  313.      * <p>
  314.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  315.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  316.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  317.      * </p>
  318.      * @param type type of the position angle to use
  319.      * @param jacobian placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix
  320.      * is larger than 6x6, only the 6x6 upper left corner will be modified
  321.      */
  322.     public void getJacobianWrtParameters(final PositionAngle type, final double[][] jacobian) {

  323.         final double[][] cachedJacobian;
  324.         synchronized (this) {
  325.             switch (type) {
  326.             case MEAN :
  327.                 if (jacobianWrtParametersMean == null) {
  328.                     // first call, we need to compute the jacobian and cache it
  329.                     jacobianWrtParametersMean = createInverseJacobian(type);
  330.                 }
  331.                 cachedJacobian = jacobianWrtParametersMean;
  332.                 break;
  333.             case ECCENTRIC :
  334.                 if (jacobianWrtParametersEccentric == null) {
  335.                     // first call, we need to compute the jacobian and cache it
  336.                     jacobianWrtParametersEccentric = createInverseJacobian(type);
  337.                 }
  338.                 cachedJacobian = jacobianWrtParametersEccentric;
  339.                 break;
  340.             case TRUE :
  341.                 if (jacobianWrtParametersTrue == null) {
  342.                     // first call, we need to compute the jacobian and cache it
  343.                     jacobianWrtParametersTrue = createInverseJacobian(type);
  344.                 }
  345.                 cachedJacobian = jacobianWrtParametersTrue;
  346.                 break;
  347.             default :
  348.                 throw OrekitException.createInternalError(null);
  349.             }
  350.         }

  351.         // fill the user-provided array
  352.         for (int i = 0; i < cachedJacobian.length; ++i) {
  353.             System.arraycopy(cachedJacobian[i], 0, jacobian[i], 0, cachedJacobian[i].length);
  354.         }

  355.     }

  356.     /** Create an inverse Jacobian.
  357.      * @param type type of the position angle to use
  358.      * @return inverse Jacobian
  359.      */
  360.     private double[][] createInverseJacobian(final PositionAngle type) {

  361.         // get the direct Jacobian
  362.         final double[][] directJacobian = new double[6][6];
  363.         getJacobianWrtCartesian(type, directJacobian);

  364.         // invert the direct Jacobian
  365.         final RealMatrix matrix = MatrixUtils.createRealMatrix(directJacobian);
  366.         final DecompositionSolver solver = new QRDecomposition(matrix).getSolver();
  367.         return solver.getInverse().getData();

  368.     }

  369.     /** Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
  370.      * <p>
  371.      * Element {@code jacobian[i][j]} is the derivative of parameter i of the orbit with
  372.      * respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
  373.      * whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
  374.      * </p>
  375.      * @return 6x6 Jacobian matrix
  376.      * @see #computeJacobianEccentricWrtCartesian()
  377.      * @see #computeJacobianTrueWrtCartesian()
  378.      */
  379.     protected abstract double[][] computeJacobianMeanWrtCartesian();

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

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

  402.     /** Add the contribution of the Keplerian motion to parameters derivatives
  403.      * <p>
  404.      * This method is used by integration-based propagators to evaluate the part of Keplerian
  405.      * motion to evolution of the orbital state.
  406.      * </p>
  407.      * @param type type of the position angle in the state
  408.      * @param gm attraction coefficient to use
  409.      * @param pDot array containing orbital state derivatives to update (the Keplerian
  410.      * part must be <em>added</em> to the array components, as the array may already
  411.      * contain some non-zero elements corresponding to non-Keplerian parts)
  412.      */
  413.     public abstract void addKeplerContribution(final PositionAngle type, final double gm, double[] pDot);

  414.         /** Fill a Jacobian half row with a single vector.
  415.      * @param a coefficient of the vector
  416.      * @param v vector
  417.      * @param row Jacobian matrix row
  418.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  419.      */
  420.     protected static void fillHalfRow(final double a, final Vector3D v, final double[] row, final int j) {
  421.         row[j]     = a * v.getX();
  422.         row[j + 1] = a * v.getY();
  423.         row[j + 2] = a * v.getZ();
  424.     }

  425.     /** Fill a Jacobian half row with a linear combination of vectors.
  426.      * @param a1 coefficient of the first vector
  427.      * @param v1 first vector
  428.      * @param a2 coefficient of the second vector
  429.      * @param v2 second vector
  430.      * @param row Jacobian matrix row
  431.      * @param j index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
  432.      */
  433.     protected static void fillHalfRow(final double a1, final Vector3D v1, final double a2, final Vector3D v2,
  434.                                       final double[] row, final int j) {
  435.         row[j]     = a1 * v1.getX() + a2 * v2.getX();
  436.         row[j + 1] = a1 * v1.getY() + a2 * v2.getY();
  437.         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ();
  438.     }

  439.     /** Fill a Jacobian half row with a linear combination of vectors.
  440.      * @param a1 coefficient of the first vector
  441.      * @param v1 first vector
  442.      * @param a2 coefficient of the second vector
  443.      * @param v2 second vector
  444.      * @param a3 coefficient of the third vector
  445.      * @param v3 third 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 a3, final Vector3D v3,
  451.                                       final double[] row, final int j) {
  452.         row[j]     = a1 * v1.getX() + a2 * v2.getX() + a3 * v3.getX();
  453.         row[j + 1] = a1 * v1.getY() + a2 * v2.getY() + a3 * v3.getY();
  454.         row[j + 2] = a1 * v1.getZ() + a2 * v2.getZ() + a3 * v3.getZ();
  455.     }

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

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

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

  521. }