EquinoctialOrbit.java

  1. /* Copyright 2002-2025 CS GROUP
  2.  * Licensed to CS GROUP (CS) under one or more
  3.  * contributor license agreements.  See the NOTICE file distributed with
  4.  * this work for additional information regarding copyright ownership.
  5.  * CS licenses this file to You under the Apache License, Version 2.0
  6.  * (the "License"); you may not use this file except in compliance with
  7.  * the License.  You may obtain a copy of the License at
  8.  *
  9.  *   http://www.apache.org/licenses/LICENSE-2.0
  10.  *
  11.  * Unless required by applicable law or agreed to in writing, software
  12.  * distributed under the License is distributed on an "AS IS" BASIS,
  13.  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
  14.  * See the License for the specific language governing permissions and
  15.  * limitations under the License.
  16.  */
  17. package org.orekit.orbits;

  18. import org.hipparchus.analysis.differentiation.UnivariateDerivative1;
  19. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  20. import org.hipparchus.util.FastMath;
  21. import org.hipparchus.util.SinCos;
  22. import org.orekit.errors.OrekitIllegalArgumentException;
  23. import org.orekit.errors.OrekitInternalError;
  24. import org.orekit.errors.OrekitMessages;
  25. import org.orekit.frames.Frame;
  26. import org.orekit.frames.KinematicTransform;
  27. import org.orekit.time.AbsoluteDate;
  28. import org.orekit.time.TimeOffset;
  29. import org.orekit.utils.PVCoordinates;
  30. import org.orekit.utils.TimeStampedPVCoordinates;


  31. /**
  32.  * This class handles equinoctial orbital parameters, which can support both
  33.  * circular and equatorial orbits.
  34.  * <p>
  35.  * The parameters used internally are the equinoctial elements which can be
  36.  * related to Keplerian elements as follows:
  37.  *   <pre>
  38.  *     a
  39.  *     ex = e cos(ω + Ω)
  40.  *     ey = e sin(ω + Ω)
  41.  *     hx = tan(i/2) cos(Ω)
  42.  *     hy = tan(i/2) sin(Ω)
  43.  *     lv = v + ω + Ω
  44.  *   </pre>
  45.  * where ω stands for the Perigee Argument and Ω stands for the
  46.  * Right Ascension of the Ascending Node.
  47.  *
  48.  * <p>
  49.  * The conversion equations from and to Keplerian elements given above hold only
  50.  * when both sides are unambiguously defined, i.e. when orbit is neither equatorial
  51.  * nor circular. When orbit is either equatorial or circular, the equinoctial
  52.  * parameters are still unambiguously defined whereas some Keplerian elements
  53.  * (more precisely ω and Ω) become ambiguous. For this reason, equinoctial
  54.  * parameters are the recommended way to represent orbits. Note however than
  55.  * the present implementation does not handle non-elliptical cases.
  56.  * </p>
  57.  * <p>
  58.  * The instance <code>EquinoctialOrbit</code> is guaranteed to be immutable.
  59.  * </p>
  60.  * @see    Orbit
  61.  * @see    KeplerianOrbit
  62.  * @see    CircularOrbit
  63.  * @see    CartesianOrbit
  64.  * @author Mathieu Rom&eacute;ro
  65.  * @author Luc Maisonobe
  66.  * @author Guylaine Prat
  67.  * @author Fabien Maussion
  68.  * @author V&eacute;ronique Pommier-Maurussane
  69.  */
  70. public class EquinoctialOrbit extends Orbit implements PositionAngleBased<EquinoctialOrbit> {

  71.     /** Semi-major axis (m). */
  72.     private final double a;

  73.     /** First component of the eccentricity vector. */
  74.     private final double ex;

  75.     /** Second component of the eccentricity vector. */
  76.     private final double ey;

  77.     /** First component of the inclination vector. */
  78.     private final double hx;

  79.     /** Second component of the inclination vector. */
  80.     private final double hy;

  81.     /** Cached longitude argument (rad). */
  82.     private final double cachedL;

  83.     /** Cache type of position angle (longitude argument). */
  84.     private final PositionAngleType cachedPositionAngleType;

  85.     /** Semi-major axis derivative (m/s). */
  86.     private final double aDot;

  87.     /** First component of the eccentricity vector derivative. */
  88.     private final double exDot;

  89.     /** Second component of the eccentricity vector derivative. */
  90.     private final double eyDot;

  91.     /** First component of the inclination vector derivative. */
  92.     private final double hxDot;

  93.     /** Second component of the inclination vector derivative. */
  94.     private final double hyDot;

  95.     /** Derivative of cached longitude argument (rad/s). */
  96.     private final double cachedLDot;

  97.     /** Partial Cartesian coordinates (position and velocity are valid, acceleration may be missing). */
  98.     private transient PVCoordinates partialPV;

  99.     /** Creates a new instance.
  100.      * @param a  semi-major axis (m)
  101.      * @param ex e cos(ω + Ω), first component of eccentricity vector
  102.      * @param ey e sin(ω + Ω), second component of eccentricity vector
  103.      * @param hx tan(i/2) cos(Ω), first component of inclination vector
  104.      * @param hy tan(i/2) sin(Ω), second component of inclination vector
  105.      * @param l  (M or E or v) + ω + Ω, mean, eccentric or true longitude argument (rad)
  106.      * @param type type of longitude argument
  107.      * @param cachedPositionAngleType type of cached longitude argument
  108.      * @param frame the frame in which the parameters are defined
  109.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  110.      * @param date date of the orbital parameters
  111.      * @param mu central attraction coefficient (m³/s²)
  112.      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
  113.      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
  114.      * @since 12.1
  115.      */
  116.     public EquinoctialOrbit(final double a, final double ex, final double ey,
  117.                             final double hx, final double hy, final double l,
  118.                             final PositionAngleType type, final PositionAngleType cachedPositionAngleType,
  119.                             final Frame frame, final AbsoluteDate date, final double mu)
  120.         throws IllegalArgumentException {
  121.         this(a, ex, ey, hx, hy, l,
  122.              0., 0., 0., 0., 0., computeKeplerianLDot(type, a, ex, ey, mu, l, type),
  123.              type, cachedPositionAngleType, frame, date, mu);
  124.     }

  125.     /** Creates a new instance without derivatives and with cached position angle same as value inputted.
  126.      * @param a  semi-major axis (m)
  127.      * @param ex e cos(ω + Ω), first component of eccentricity vector
  128.      * @param ey e sin(ω + Ω), second component of eccentricity vector
  129.      * @param hx tan(i/2) cos(Ω), first component of inclination vector
  130.      * @param hy tan(i/2) sin(Ω), second component of inclination vector
  131.      * @param l  (M or E or v) + ω + Ω, mean, eccentric or true longitude argument (rad)
  132.      * @param type type of longitude argument
  133.      * @param frame the frame in which the parameters are defined
  134.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  135.      * @param date date of the orbital parameters
  136.      * @param mu central attraction coefficient (m³/s²)
  137.      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
  138.      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
  139.      */
  140.     public EquinoctialOrbit(final double a, final double ex, final double ey,
  141.                             final double hx, final double hy, final double l,
  142.                             final PositionAngleType type,
  143.                             final Frame frame, final AbsoluteDate date, final double mu)
  144.             throws IllegalArgumentException {
  145.         this(a, ex, ey, hx, hy, l, type, type, frame, date, mu);
  146.     }

  147.     /** Creates a new instance.
  148.      * @param a  semi-major axis (m)
  149.      * @param ex e cos(ω + Ω), first component of eccentricity vector
  150.      * @param ey e sin(ω + Ω), second component of eccentricity vector
  151.      * @param hx tan(i/2) cos(Ω), first component of inclination vector
  152.      * @param hy tan(i/2) sin(Ω), second component of inclination vector
  153.      * @param l  (M or E or v) + ω + Ω, mean, eccentric or true longitude argument (rad)
  154.      * @param aDot  semi-major axis derivative (m/s)
  155.      * @param exDot d(e cos(ω + Ω))/dt, first component of eccentricity vector derivative
  156.      * @param eyDot d(e sin(ω + Ω))/dt, second component of eccentricity vector derivative
  157.      * @param hxDot d(tan(i/2) cos(Ω))/dt, first component of inclination vector derivative
  158.      * @param hyDot d(tan(i/2) sin(Ω))/dt, second component of inclination vector derivative
  159.      * @param lDot  d(M or E or v) + ω + Ω)/dr, mean, eccentric or true longitude argument  derivative (rad/s)
  160.      * @param type type of longitude argument
  161.      * @param cachedPositionAngleType of cached longitude argument
  162.      * @param frame the frame in which the parameters are defined
  163.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  164.      * @param date date of the orbital parameters
  165.      * @param mu central attraction coefficient (m³/s²)
  166.      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
  167.      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
  168.      * @since 12.1
  169.      */
  170.     public EquinoctialOrbit(final double a, final double ex, final double ey,
  171.                             final double hx, final double hy, final double l,
  172.                             final double aDot, final double exDot, final double eyDot,
  173.                             final double hxDot, final double hyDot, final double lDot,
  174.                             final PositionAngleType type, final PositionAngleType cachedPositionAngleType,
  175.                             final Frame frame, final AbsoluteDate date, final double mu)
  176.         throws IllegalArgumentException {
  177.         super(frame, date, mu);
  178.         if (ex * ex + ey * ey >= 1.0) {
  179.             throw new OrekitIllegalArgumentException(OrekitMessages.HYPERBOLIC_ORBIT_NOT_HANDLED_AS,
  180.                                                      getClass().getName());
  181.         }
  182.         this.cachedPositionAngleType = cachedPositionAngleType;
  183.         this.a     = a;
  184.         this.aDot  = aDot;
  185.         this.ex    = ex;
  186.         this.exDot = exDot;
  187.         this.ey    = ey;
  188.         this.eyDot = eyDot;
  189.         this.hx    = hx;
  190.         this.hxDot = hxDot;
  191.         this.hy    = hy;
  192.         this.hyDot = hyDot;

  193.         final UnivariateDerivative1 lUD = initializeCachedL(l, lDot, type);
  194.         this.cachedL = lUD.getValue();
  195.         this.cachedLDot = lUD.getFirstDerivative();

  196.         this.partialPV = null;

  197.     }

  198.     /** Creates a new instance with derivatives and with cached position angle same as value inputted.
  199.      * @param a  semi-major axis (m)
  200.      * @param ex e cos(ω + Ω), first component of eccentricity vector
  201.      * @param ey e sin(ω + Ω), second component of eccentricity vector
  202.      * @param hx tan(i/2) cos(Ω), first component of inclination vector
  203.      * @param hy tan(i/2) sin(Ω), second component of inclination vector
  204.      * @param l  (M or E or v) + ω + Ω, mean, eccentric or true longitude argument (rad)
  205.      * @param aDot  semi-major axis derivative (m/s)
  206.      * @param exDot d(e cos(ω + Ω))/dt, first component of eccentricity vector derivative
  207.      * @param eyDot d(e sin(ω + Ω))/dt, second component of eccentricity vector derivative
  208.      * @param hxDot d(tan(i/2) cos(Ω))/dt, first component of inclination vector derivative
  209.      * @param hyDot d(tan(i/2) sin(Ω))/dt, second component of inclination vector derivative
  210.      * @param lDot  d(M or E or v) + ω + Ω)/dr, mean, eccentric or true longitude argument  derivative (rad/s)
  211.      * @param type type of longitude argument
  212.      * @param frame the frame in which the parameters are defined
  213.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  214.      * @param date date of the orbital parameters
  215.      * @param mu central attraction coefficient (m³/s²)
  216.      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
  217.      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
  218.      */
  219.     public EquinoctialOrbit(final double a, final double ex, final double ey,
  220.                             final double hx, final double hy, final double l,
  221.                             final double aDot, final double exDot, final double eyDot,
  222.                             final double hxDot, final double hyDot, final double lDot,
  223.                             final PositionAngleType type,
  224.                             final Frame frame, final AbsoluteDate date, final double mu)
  225.             throws IllegalArgumentException {
  226.         this(a, ex, ey, hx, hy, l, aDot, exDot, eyDot, hxDot, hyDot, lDot, type, type, frame, date, mu);
  227.     }

  228.     /** Constructor from Cartesian parameters.
  229.      *
  230.      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
  231.      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
  232.      * use {@code mu} and the position to compute the acceleration, including
  233.      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
  234.      *
  235.      * @param pvCoordinates the position, velocity and acceleration
  236.      * @param frame the frame in which are defined the {@link PVCoordinates}
  237.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  238.      * @param mu central attraction coefficient (m³/s²)
  239.      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
  240.      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
  241.      */
  242.     public EquinoctialOrbit(final TimeStampedPVCoordinates pvCoordinates,
  243.                             final Frame frame, final double mu)
  244.         throws IllegalArgumentException {
  245.         super(pvCoordinates, frame, mu);

  246.         //  compute semi-major axis
  247.         final Vector3D pvP   = pvCoordinates.getPosition();
  248.         final Vector3D pvV   = pvCoordinates.getVelocity();
  249.         final Vector3D pvA   = pvCoordinates.getAcceleration();
  250.         final double r2      = pvP.getNormSq();
  251.         final double r       = FastMath.sqrt(r2);
  252.         final double V2      = pvV.getNormSq();
  253.         final double rV2OnMu = r * V2 / mu;

  254.         // compute semi-major axis
  255.         a = r / (2 - rV2OnMu);

  256.         if (!isElliptical()) {
  257.             throw new OrekitIllegalArgumentException(OrekitMessages.HYPERBOLIC_ORBIT_NOT_HANDLED_AS,
  258.                                                      getClass().getName());
  259.         }

  260.         // compute inclination vector
  261.         final Vector3D w = pvCoordinates.getMomentum().normalize();
  262.         final double d = 1.0 / (1 + w.getZ());
  263.         hx = -d * w.getY();
  264.         hy =  d * w.getX();

  265.         // compute true longitude argument
  266.         cachedPositionAngleType = PositionAngleType.TRUE;
  267.         final double cLv = (pvP.getX() - d * pvP.getZ() * w.getX()) / r;
  268.         final double sLv = (pvP.getY() - d * pvP.getZ() * w.getY()) / r;
  269.         cachedL = FastMath.atan2(sLv, cLv);

  270.         // compute eccentricity vector
  271.         final double eSE = Vector3D.dotProduct(pvP, pvV) / FastMath.sqrt(mu * a);
  272.         final double eCE = rV2OnMu - 1;
  273.         final double e2  = eCE * eCE + eSE * eSE;
  274.         final double f   = eCE - e2;
  275.         final double g   = FastMath.sqrt(1 - e2) * eSE;
  276.         ex = a * (f * cLv + g * sLv) / r;
  277.         ey = a * (f * sLv - g * cLv) / r;

  278.         partialPV = pvCoordinates;

  279.         if (hasNonKeplerianAcceleration(pvCoordinates, mu)) {
  280.             // we have a relevant acceleration, we can compute derivatives

  281.             final double[][] jacobian = new double[6][6];
  282.             getJacobianWrtCartesian(PositionAngleType.MEAN, jacobian);

  283.             final Vector3D keplerianAcceleration    = new Vector3D(-mu / (r * r2), pvP);
  284.             final Vector3D nonKeplerianAcceleration = pvA.subtract(keplerianAcceleration);
  285.             final double   aX                       = nonKeplerianAcceleration.getX();
  286.             final double   aY                       = nonKeplerianAcceleration.getY();
  287.             final double   aZ                       = nonKeplerianAcceleration.getZ();
  288.             aDot  = jacobian[0][3] * aX + jacobian[0][4] * aY + jacobian[0][5] * aZ;
  289.             exDot = jacobian[1][3] * aX + jacobian[1][4] * aY + jacobian[1][5] * aZ;
  290.             eyDot = jacobian[2][3] * aX + jacobian[2][4] * aY + jacobian[2][5] * aZ;
  291.             hxDot = jacobian[3][3] * aX + jacobian[3][4] * aY + jacobian[3][5] * aZ;
  292.             hyDot = jacobian[4][3] * aX + jacobian[4][4] * aY + jacobian[4][5] * aZ;

  293.             // in order to compute true longitude argument derivative, we must compute
  294.             // mean longitude argument derivative including Keplerian motion and convert to true longitude argument
  295.             final double lMDot = getKeplerianMeanMotion() +
  296.                                  jacobian[5][3] * aX + jacobian[5][4] * aY + jacobian[5][5] * aZ;
  297.             final UnivariateDerivative1 exUD = new UnivariateDerivative1(ex, exDot);
  298.             final UnivariateDerivative1 eyUD = new UnivariateDerivative1(ey, eyDot);
  299.             final UnivariateDerivative1 lMUD = new UnivariateDerivative1(getLM(), lMDot);
  300.             final UnivariateDerivative1 lvUD = FieldEquinoctialLongitudeArgumentUtility.meanToTrue(exUD, eyUD, lMUD);
  301.             cachedLDot = lvUD.getFirstDerivative();

  302.         } else {
  303.             // acceleration is either almost zero or NaN,
  304.             // we assume acceleration was not known
  305.             // we don't set up derivatives
  306.             aDot  = 0.;
  307.             exDot = 0.;
  308.             eyDot = 0.;
  309.             hxDot = 0.;
  310.             hyDot = 0.;
  311.             cachedLDot = computeKeplerianLDot(cachedPositionAngleType, a, ex, ey, mu, cachedL, cachedPositionAngleType);
  312.         }

  313.     }

  314.     /** Constructor from Cartesian parameters.
  315.      *
  316.      * <p> The acceleration provided in {@code pvCoordinates} is accessible using
  317.      * {@link #getPVCoordinates()} and {@link #getPVCoordinates(Frame)}. All other methods
  318.      * use {@code mu} and the position to compute the acceleration, including
  319.      * {@link #shiftedBy(double)} and {@link #getPVCoordinates(AbsoluteDate, Frame)}.
  320.      *
  321.      * @param pvCoordinates the position end velocity
  322.      * @param frame the frame in which are defined the {@link PVCoordinates}
  323.      * (<em>must</em> be a {@link Frame#isPseudoInertial pseudo-inertial frame})
  324.      * @param date date of the orbital parameters
  325.      * @param mu central attraction coefficient (m³/s²)
  326.      * @exception IllegalArgumentException if eccentricity is equal to 1 or larger or
  327.      * if frame is not a {@link Frame#isPseudoInertial pseudo-inertial frame}
  328.      */
  329.     public EquinoctialOrbit(final PVCoordinates pvCoordinates, final Frame frame,
  330.                             final AbsoluteDate date, final double mu)
  331.         throws IllegalArgumentException {
  332.         this(new TimeStampedPVCoordinates(date, pvCoordinates), frame, mu);
  333.     }

  334.     /** Constructor from any kind of orbital parameters.
  335.      * @param op orbital parameters to copy
  336.      */
  337.     public EquinoctialOrbit(final Orbit op) {
  338.         super(op.getFrame(), op.getDate(), op.getMu());
  339.         a         = op.getA();
  340.         aDot      = op.getADot();
  341.         ex        = op.getEquinoctialEx();
  342.         exDot     = op.getEquinoctialExDot();
  343.         ey        = op.getEquinoctialEy();
  344.         eyDot     = op.getEquinoctialEyDot();
  345.         hx        = op.getHx();
  346.         hxDot     = op.getHxDot();
  347.         hy        = op.getHy();
  348.         hyDot     = op.getHyDot();
  349.         cachedPositionAngleType = PositionAngleType.TRUE;
  350.         cachedL   = op.getLv();
  351.         cachedLDot = op.hasNonKeplerianAcceleration() ? op.getLvDot() :
  352.                 computeKeplerianLDot(cachedPositionAngleType, a, ex, ey, op.getMu(), cachedL, cachedPositionAngleType);
  353.         partialPV = null;
  354.     }

  355.     /** {@inheritDoc} */
  356.     @Override
  357.     public boolean hasNonKeplerianAcceleration() {
  358.         return aDot != 0. || exDot != 0. || eyDot != 0. || hxDot != 0. || hyDot != 0. ||
  359.                 FastMath.abs(cachedLDot - computeKeplerianLDot(cachedPositionAngleType, a, ex, ey, getMu(), cachedL, cachedPositionAngleType)) > TOLERANCE_POSITION_ANGLE_RATE;
  360.     }

  361.     /** {@inheritDoc} */
  362.     @Override
  363.     public OrbitType getType() {
  364.         return OrbitType.EQUINOCTIAL;
  365.     }

  366.     /** {@inheritDoc} */
  367.     @Override
  368.     public double getA() {
  369.         return a;
  370.     }

  371.     /** {@inheritDoc} */
  372.     @Override
  373.     public double getADot() {
  374.         return aDot;
  375.     }

  376.     /** {@inheritDoc} */
  377.     @Override
  378.     public double getEquinoctialEx() {
  379.         return ex;
  380.     }

  381.     /** {@inheritDoc} */
  382.     @Override
  383.     public double getEquinoctialExDot() {
  384.         return exDot;
  385.     }

  386.     /** {@inheritDoc} */
  387.     @Override
  388.     public double getEquinoctialEy() {
  389.         return ey;
  390.     }

  391.     /** {@inheritDoc} */
  392.     @Override
  393.     public double getEquinoctialEyDot() {
  394.         return eyDot;
  395.     }

  396.     /** {@inheritDoc} */
  397.     @Override
  398.     public double getHx() {
  399.         return hx;
  400.     }

  401.     /** {@inheritDoc} */
  402.     @Override
  403.     public double getHxDot() {
  404.         return hxDot;
  405.     }

  406.     /** {@inheritDoc} */
  407.     @Override
  408.     public double getHy() {
  409.         return hy;
  410.     }

  411.     /** {@inheritDoc} */
  412.     @Override
  413.     public double getHyDot() {
  414.         return hyDot;
  415.     }

  416.     /** {@inheritDoc} */
  417.     @Override
  418.     public double getLv() {
  419.         switch (cachedPositionAngleType) {
  420.             case TRUE:
  421.                 return cachedL;

  422.             case ECCENTRIC:
  423.                 return EquinoctialLongitudeArgumentUtility.eccentricToTrue(ex, ey, cachedL);

  424.             case MEAN:
  425.                 return EquinoctialLongitudeArgumentUtility.meanToTrue(ex, ey, cachedL);

  426.             default:
  427.                 throw new OrekitInternalError(null);
  428.         }
  429.     }

  430.     /** {@inheritDoc} */
  431.     @Override
  432.     public double getLvDot() {
  433.         switch (cachedPositionAngleType) {
  434.             case ECCENTRIC:
  435.                 final UnivariateDerivative1 lEUD = new UnivariateDerivative1(cachedL, cachedLDot);
  436.                 final UnivariateDerivative1 exUD     = new UnivariateDerivative1(ex,     exDot);
  437.                 final UnivariateDerivative1 eyUD     = new UnivariateDerivative1(ey,     eyDot);
  438.                 final UnivariateDerivative1 lvUD = FieldEquinoctialLongitudeArgumentUtility.eccentricToTrue(exUD, eyUD,
  439.                         lEUD);
  440.                 return lvUD.getFirstDerivative();

  441.             case TRUE:
  442.                 return cachedLDot;

  443.             case MEAN:
  444.                 final UnivariateDerivative1 lMUD = new UnivariateDerivative1(cachedL, cachedLDot);
  445.                 final UnivariateDerivative1 exUD2    = new UnivariateDerivative1(ex,     exDot);
  446.                 final UnivariateDerivative1 eyUD2    = new UnivariateDerivative1(ey,     eyDot);
  447.                 final UnivariateDerivative1 lvUD2 = FieldEquinoctialLongitudeArgumentUtility.meanToTrue(exUD2,
  448.                         eyUD2, lMUD);
  449.                 return lvUD2.getFirstDerivative();

  450.             default:
  451.                 throw new OrekitInternalError(null);
  452.         }
  453.     }

  454.     /** {@inheritDoc} */
  455.     @Override
  456.     public double getLE() {
  457.         switch (cachedPositionAngleType) {
  458.             case TRUE:
  459.                 return EquinoctialLongitudeArgumentUtility.trueToEccentric(ex, ey, cachedL);

  460.             case ECCENTRIC:
  461.                 return cachedL;

  462.             case MEAN:
  463.                 return EquinoctialLongitudeArgumentUtility.meanToEccentric(ex, ey, cachedL);

  464.             default:
  465.                 throw new OrekitInternalError(null);
  466.         }
  467.     }

  468.     /** {@inheritDoc} */
  469.     @Override
  470.     public double getLEDot() {
  471.         switch (cachedPositionAngleType) {
  472.             case TRUE:
  473.                 final UnivariateDerivative1 lvUD = new UnivariateDerivative1(cachedL, cachedLDot);
  474.                 final UnivariateDerivative1 exUD     = new UnivariateDerivative1(ex,     exDot);
  475.                 final UnivariateDerivative1 eyUD     = new UnivariateDerivative1(ey,     eyDot);
  476.                 final UnivariateDerivative1 lEUD = FieldEquinoctialLongitudeArgumentUtility.trueToEccentric(exUD, eyUD,
  477.                         lvUD);
  478.                 return lEUD.getFirstDerivative();

  479.             case ECCENTRIC:
  480.                 return cachedLDot;

  481.             case MEAN:
  482.                 final UnivariateDerivative1 lMUD = new UnivariateDerivative1(cachedL, cachedLDot);
  483.                 final UnivariateDerivative1 exUD2    = new UnivariateDerivative1(ex,     exDot);
  484.                 final UnivariateDerivative1 eyUD2    = new UnivariateDerivative1(ey,     eyDot);
  485.                 final UnivariateDerivative1 lEUD2 = FieldEquinoctialLongitudeArgumentUtility.meanToEccentric(exUD2,
  486.                         eyUD2, lMUD);
  487.                 return lEUD2.getFirstDerivative();

  488.             default:
  489.                 throw new OrekitInternalError(null);
  490.         }
  491.     }

  492.     /** {@inheritDoc} */
  493.     @Override
  494.     public double getLM() {
  495.         switch (cachedPositionAngleType) {
  496.             case TRUE:
  497.                 return EquinoctialLongitudeArgumentUtility.trueToMean(ex, ey, cachedL);

  498.             case MEAN:
  499.                 return cachedL;

  500.             case ECCENTRIC:
  501.                 return EquinoctialLongitudeArgumentUtility.eccentricToMean(ex, ey, cachedL);

  502.             default:
  503.                 throw new OrekitInternalError(null);
  504.         }
  505.     }

  506.     /** {@inheritDoc} */
  507.     @Override
  508.     public double getLMDot() {
  509.         switch (cachedPositionAngleType) {
  510.             case TRUE:
  511.                 final UnivariateDerivative1 lvUD = new UnivariateDerivative1(cachedL, cachedLDot);
  512.                 final UnivariateDerivative1 exUD     = new UnivariateDerivative1(ex,     exDot);
  513.                 final UnivariateDerivative1 eyUD     = new UnivariateDerivative1(ey,     eyDot);
  514.                 final UnivariateDerivative1 lMUD = FieldEquinoctialLongitudeArgumentUtility.trueToMean(exUD, eyUD, lvUD);
  515.                 return lMUD.getFirstDerivative();

  516.             case MEAN:
  517.                 return cachedLDot;

  518.             case ECCENTRIC:
  519.                 final UnivariateDerivative1 lEUD = new UnivariateDerivative1(cachedL, cachedLDot);
  520.                 final UnivariateDerivative1 exUD2    = new UnivariateDerivative1(ex,     exDot);
  521.                 final UnivariateDerivative1 eyUD2    = new UnivariateDerivative1(ey,     eyDot);
  522.                 final UnivariateDerivative1 lMUD2 = FieldEquinoctialLongitudeArgumentUtility.eccentricToMean(exUD2,
  523.                         eyUD2, lEUD);
  524.                 return lMUD2.getFirstDerivative();

  525.             default:
  526.                 throw new OrekitInternalError(null);
  527.         }
  528.     }

  529.     /** Get the longitude argument.
  530.      * @param type type of the angle
  531.      * @return longitude argument (rad)
  532.      */
  533.     public double getL(final PositionAngleType type) {
  534.         return (type == PositionAngleType.MEAN) ? getLM() :
  535.                                               ((type == PositionAngleType.ECCENTRIC) ? getLE() :
  536.                                                                                    getLv());
  537.     }

  538.     /** Get the longitude argument derivative.
  539.      * @param type type of the angle
  540.      * @return longitude argument derivative (rad/s)
  541.      */
  542.     public double getLDot(final PositionAngleType type) {
  543.         return (type == PositionAngleType.MEAN) ? getLMDot() :
  544.                                               ((type == PositionAngleType.ECCENTRIC) ? getLEDot() :
  545.                                                                                    getLvDot());
  546.     }

  547.     /** {@inheritDoc} */
  548.     @Override
  549.     public double getE() {
  550.         return FastMath.sqrt(ex * ex + ey * ey);
  551.     }

  552.     /** {@inheritDoc} */
  553.     @Override
  554.     public double getEDot() {
  555.         if (!hasNonKeplerianAcceleration()) {
  556.             return 0.;
  557.         }
  558.         return (ex * exDot + ey * eyDot) / FastMath.sqrt(ex * ex + ey * ey);
  559.     }

  560.     /** {@inheritDoc} */
  561.     @Override
  562.     public double getI() {
  563.         return 2 * FastMath.atan(FastMath.sqrt(hx * hx + hy * hy));
  564.     }

  565.     /** {@inheritDoc} */
  566.     @Override
  567.     public double getIDot() {
  568.         if (!hasNonKeplerianAcceleration()) {
  569.             return 0.;
  570.         }
  571.         final double h2 = hx * hx + hy * hy;
  572.         final double h  = FastMath.sqrt(h2);
  573.         return 2 * (hx * hxDot + hy * hyDot) / (h * (1 + h2));
  574.     }

  575.     /** Compute position and velocity but not acceleration.
  576.      */
  577.     private void computePVWithoutA() {

  578.         if (partialPV != null) {
  579.             // already computed
  580.             return;
  581.         }

  582.         // get equinoctial parameters
  583.         final double lE = getLE();

  584.         // inclination-related intermediate parameters
  585.         final double hx2   = hx * hx;
  586.         final double hy2   = hy * hy;
  587.         final double factH = 1. / (1 + hx2 + hy2);

  588.         // reference axes defining the orbital plane
  589.         final double ux = (1 + hx2 - hy2) * factH;
  590.         final double uy =  2 * hx * hy * factH;
  591.         final double uz = -2 * hy * factH;

  592.         final double vx = uy;
  593.         final double vy = (1 - hx2 + hy2) * factH;
  594.         final double vz =  2 * hx * factH;

  595.         // eccentricity-related intermediate parameters
  596.         final double exey = ex * ey;
  597.         final double ex2  = ex * ex;
  598.         final double ey2  = ey * ey;
  599.         final double e2   = ex2 + ey2;
  600.         final double eta  = 1 + FastMath.sqrt(1 - e2);
  601.         final double beta = 1. / eta;

  602.         // eccentric longitude argument
  603.         final SinCos scLe   = FastMath.sinCos(lE);
  604.         final double cLe    = scLe.cos();
  605.         final double sLe    = scLe.sin();
  606.         final double exCeyS = ex * cLe + ey * sLe;

  607.         // coordinates of position and velocity in the orbital plane
  608.         final double x      = a * ((1 - beta * ey2) * cLe + beta * exey * sLe - ex);
  609.         final double y      = a * ((1 - beta * ex2) * sLe + beta * exey * cLe - ey);

  610.         final double factor = FastMath.sqrt(getMu() / a) / (1 - exCeyS);
  611.         final double xdot   = factor * (-sLe + beta * ey * exCeyS);
  612.         final double ydot   = factor * ( cLe - beta * ex * exCeyS);

  613.         final Vector3D position =
  614.                         new Vector3D(x * ux + y * vx, x * uy + y * vy, x * uz + y * vz);
  615.         final Vector3D velocity =
  616.                         new Vector3D(xdot * ux + ydot * vx, xdot * uy + ydot * vy, xdot * uz + ydot * vz);
  617.         partialPV = new PVCoordinates(position, velocity);

  618.     }

  619.     /** Initialize cached argument of longitude with rate.
  620.      * @param l input argument of longitude
  621.      * @param lDot rate of input argument of longitude
  622.      * @param inputType position angle type passed as input
  623.      * @return argument of longitude to cache with rate
  624.      * @since 12.1
  625.      */
  626.     private UnivariateDerivative1 initializeCachedL(final double l, final double lDot,
  627.                                                     final PositionAngleType inputType) {
  628.         if (cachedPositionAngleType == inputType) {
  629.             return new UnivariateDerivative1(l, lDot);

  630.         } else {
  631.             final UnivariateDerivative1 exUD = new UnivariateDerivative1(ex, exDot);
  632.             final UnivariateDerivative1 eyUD = new UnivariateDerivative1(ey, eyDot);
  633.             final UnivariateDerivative1 lUD = new UnivariateDerivative1(l, lDot);

  634.             switch (cachedPositionAngleType) {

  635.                 case ECCENTRIC:
  636.                     if (inputType == PositionAngleType.MEAN) {
  637.                         return FieldEquinoctialLongitudeArgumentUtility.meanToEccentric(exUD, eyUD, lUD);
  638.                     } else {
  639.                         return FieldEquinoctialLongitudeArgumentUtility.trueToEccentric(exUD, eyUD, lUD);
  640.                     }

  641.                 case TRUE:
  642.                     if (inputType == PositionAngleType.MEAN) {
  643.                         return FieldEquinoctialLongitudeArgumentUtility.meanToTrue(exUD, eyUD, lUD);
  644.                     } else {
  645.                         return FieldEquinoctialLongitudeArgumentUtility.eccentricToTrue(exUD, eyUD, lUD);
  646.                     }

  647.                 case MEAN:
  648.                     if (inputType == PositionAngleType.TRUE) {
  649.                         return FieldEquinoctialLongitudeArgumentUtility.trueToMean(exUD, eyUD, lUD);
  650.                     } else {
  651.                         return FieldEquinoctialLongitudeArgumentUtility.eccentricToMean(exUD, eyUD, lUD);
  652.                     }

  653.                 default:
  654.                     throw new OrekitInternalError(null);

  655.             }

  656.         }

  657.     }

  658.     /** Compute non-Keplerian part of the acceleration from first time derivatives.
  659.      * @return non-Keplerian part of the acceleration
  660.      */
  661.     private Vector3D nonKeplerianAcceleration() {

  662.         final double[][] dCdP = new double[6][6];
  663.         getJacobianWrtParameters(PositionAngleType.MEAN, dCdP);

  664.         final double nonKeplerianMeanMotion = getLMDot() - getKeplerianMeanMotion();
  665.         final double nonKeplerianAx = dCdP[3][0] * aDot  + dCdP[3][1] * exDot + dCdP[3][2] * eyDot +
  666.                                       dCdP[3][3] * hxDot + dCdP[3][4] * hyDot + dCdP[3][5] * nonKeplerianMeanMotion;
  667.         final double nonKeplerianAy = dCdP[4][0] * aDot  + dCdP[4][1] * exDot + dCdP[4][2] * eyDot +
  668.                                       dCdP[4][3] * hxDot + dCdP[4][4] * hyDot + dCdP[4][5] * nonKeplerianMeanMotion;
  669.         final double nonKeplerianAz = dCdP[5][0] * aDot  + dCdP[5][1] * exDot + dCdP[5][2] * eyDot +
  670.                                       dCdP[5][3] * hxDot + dCdP[5][4] * hyDot + dCdP[5][5] * nonKeplerianMeanMotion;

  671.         return new Vector3D(nonKeplerianAx, nonKeplerianAy, nonKeplerianAz);

  672.     }

  673.     /** {@inheritDoc} */
  674.     @Override
  675.     protected Vector3D initPosition() {

  676.         // get equinoctial parameters
  677.         final double lE = getLE();

  678.         // inclination-related intermediate parameters
  679.         final double hx2   = hx * hx;
  680.         final double hy2   = hy * hy;
  681.         final double factH = 1. / (1 + hx2 + hy2);

  682.         // reference axes defining the orbital plane
  683.         final double ux = (1 + hx2 - hy2) * factH;
  684.         final double uy =  2 * hx * hy * factH;
  685.         final double uz = -2 * hy * factH;

  686.         final double vx = uy;
  687.         final double vy = (1 - hx2 + hy2) * factH;
  688.         final double vz =  2 * hx * factH;

  689.         // eccentricity-related intermediate parameters
  690.         final double exey = ex * ey;
  691.         final double ex2  = ex * ex;
  692.         final double ey2  = ey * ey;
  693.         final double e2   = ex2 + ey2;
  694.         final double eta  = 1 + FastMath.sqrt(1 - e2);
  695.         final double beta = 1. / eta;

  696.         // eccentric longitude argument
  697.         final SinCos scLe   = FastMath.sinCos(lE);
  698.         final double cLe    = scLe.cos();
  699.         final double sLe    = scLe.sin();

  700.         // coordinates of position and velocity in the orbital plane
  701.         final double x      = a * ((1 - beta * ey2) * cLe + beta * exey * sLe - ex);
  702.         final double y      = a * ((1 - beta * ex2) * sLe + beta * exey * cLe - ey);

  703.         return new Vector3D(x * ux + y * vx, x * uy + y * vy, x * uz + y * vz);

  704.     }

  705.     /** {@inheritDoc} */
  706.     @Override
  707.     protected TimeStampedPVCoordinates initPVCoordinates() {

  708.         // position and velocity
  709.         computePVWithoutA();

  710.         // acceleration
  711.         final double r2 = partialPV.getPosition().getNormSq();
  712.         final Vector3D keplerianAcceleration = new Vector3D(-getMu() / (r2 * FastMath.sqrt(r2)), partialPV.getPosition());
  713.         final Vector3D acceleration = hasNonKeplerianRates() ?
  714.                                       keplerianAcceleration.add(nonKeplerianAcceleration()) :
  715.                                       keplerianAcceleration;

  716.         return new TimeStampedPVCoordinates(getDate(), partialPV.getPosition(), partialPV.getVelocity(), acceleration);

  717.     }

  718.     /** {@inheritDoc} */
  719.     @Override
  720.     public EquinoctialOrbit inFrame(final Frame inertialFrame) {
  721.         final PVCoordinates pvCoordinates;
  722.         if (hasNonKeplerianAcceleration()) {
  723.             pvCoordinates = getPVCoordinates(inertialFrame);
  724.         } else {
  725.             final KinematicTransform transform = getFrame().getKinematicTransformTo(inertialFrame, getDate());
  726.             pvCoordinates = transform.transformOnlyPV(getPVCoordinates());
  727.         }
  728.         final EquinoctialOrbit equinoctialOrbit = new EquinoctialOrbit(pvCoordinates, inertialFrame, getDate(), getMu());
  729.         if (equinoctialOrbit.getCachedPositionAngleType() == getCachedPositionAngleType()) {
  730.             return equinoctialOrbit;
  731.         } else {
  732.             return equinoctialOrbit.withCachedPositionAngleType(getCachedPositionAngleType());
  733.         }
  734.     }

  735.     /** {@inheritDoc} */
  736.     @Override
  737.     public EquinoctialOrbit withCachedPositionAngleType(final PositionAngleType positionAngleType) {
  738.         return new EquinoctialOrbit(a, ex, ey, hx, hy, getL(positionAngleType), aDot, exDot, eyDot, hxDot, hyDot,
  739.                 getLDot(positionAngleType), positionAngleType, getFrame(), getDate(), getMu());
  740.     }

  741.     /** {@inheritDoc} */
  742.     @Override
  743.     public EquinoctialOrbit shiftedBy(final double dt) {
  744.         return shiftedBy(new TimeOffset(dt));
  745.     }

  746.     /** {@inheritDoc} */
  747.     @Override
  748.     public EquinoctialOrbit shiftedBy(final TimeOffset dt) {

  749.         final double dtS = dt.toDouble();

  750.         // use Keplerian-only motion
  751.         final EquinoctialOrbit keplerianShifted = new EquinoctialOrbit(a, ex, ey, hx, hy,
  752.                                                                        getLM() + getKeplerianMeanMotion() * dtS,
  753.                                                                        PositionAngleType.MEAN, cachedPositionAngleType,
  754.                                                                        getFrame(),
  755.                                                                        getDate().shiftedBy(dt), getMu());

  756.         if (hasNonKeplerianRates()) {

  757.             // extract non-Keplerian acceleration from first time derivatives
  758.             final Vector3D nonKeplerianAcceleration = nonKeplerianAcceleration();

  759.             // add quadratic effect of non-Keplerian acceleration to Keplerian-only shift
  760.             keplerianShifted.computePVWithoutA();
  761.             final Vector3D fixedP   = new Vector3D(1, keplerianShifted.partialPV.getPosition(),
  762.                                                    0.5 * dtS * dtS, nonKeplerianAcceleration);
  763.             final double   fixedR2 = fixedP.getNormSq();
  764.             final double   fixedR  = FastMath.sqrt(fixedR2);
  765.             final Vector3D fixedV  = new Vector3D(1, keplerianShifted.partialPV.getVelocity(),
  766.                                                   dtS, nonKeplerianAcceleration);
  767.             final Vector3D fixedA  = new Vector3D(-getMu() / (fixedR2 * fixedR), keplerianShifted.partialPV.getPosition(),
  768.                                                   1, nonKeplerianAcceleration);

  769.             // build a new orbit, taking non-Keplerian acceleration into account
  770.             return new EquinoctialOrbit(new TimeStampedPVCoordinates(keplerianShifted.getDate(),
  771.                                                                      fixedP, fixedV, fixedA),
  772.                                         keplerianShifted.getFrame(), keplerianShifted.getMu());

  773.         } else {
  774.             // Keplerian-only motion is all we can do
  775.             return keplerianShifted;
  776.         }

  777.     }

  778.     /** {@inheritDoc} */
  779.     @Override
  780.     protected double[][] computeJacobianMeanWrtCartesian() {

  781.         final double[][] jacobian = new double[6][6];

  782.         // compute various intermediate parameters
  783.         computePVWithoutA();
  784.         final Vector3D position = partialPV.getPosition();
  785.         final Vector3D velocity = partialPV.getVelocity();
  786.         final double r2         = position.getNormSq();
  787.         final double r          = FastMath.sqrt(r2);
  788.         final double r3         = r * r2;

  789.         final double mu         = getMu();
  790.         final double sqrtMuA    = FastMath.sqrt(a * mu);
  791.         final double a2         = a * a;

  792.         final double e2         = ex * ex + ey * ey;
  793.         final double oMe2       = 1 - e2;
  794.         final double epsilon    = FastMath.sqrt(oMe2);
  795.         final double beta       = 1 / (1 + epsilon);
  796.         final double ratio      = epsilon * beta;

  797.         final double hx2        = hx * hx;
  798.         final double hy2        = hy * hy;
  799.         final double hxhy       = hx * hy;

  800.         // precomputing equinoctial frame unit vectors (f, g, w)
  801.         final Vector3D f  = new Vector3D(1 - hy2 + hx2, 2 * hxhy, -2 * hy).normalize();
  802.         final Vector3D g  = new Vector3D(2 * hxhy, 1 + hy2 - hx2, 2 * hx).normalize();
  803.         final Vector3D w  = Vector3D.crossProduct(position, velocity).normalize();

  804.         // coordinates of the spacecraft in the equinoctial frame
  805.         final double x    = Vector3D.dotProduct(position, f);
  806.         final double y    = Vector3D.dotProduct(position, g);
  807.         final double xDot = Vector3D.dotProduct(velocity, f);
  808.         final double yDot = Vector3D.dotProduct(velocity, g);

  809.         // drDot / dEx = dXDot / dEx * f + dYDot / dEx * g
  810.         final double c1 = a / (sqrtMuA * epsilon);
  811.         final double c2 = a * sqrtMuA * beta / r3;
  812.         final double c3 = sqrtMuA / (r3 * epsilon);
  813.         final Vector3D drDotSdEx = new Vector3D( c1 * xDot * yDot - c2 * ey * x - c3 * x * y, f,
  814.                                                 -c1 * xDot * xDot - c2 * ey * y + c3 * x * x, g);

  815.         // drDot / dEy = dXDot / dEy * f + dYDot / dEy * g
  816.         final Vector3D drDotSdEy = new Vector3D( c1 * yDot * yDot + c2 * ex * x - c3 * y * y, f,
  817.                                                 -c1 * xDot * yDot + c2 * ex * y + c3 * x * y, g);

  818.         // da
  819.         final Vector3D vectorAR = new Vector3D(2 * a2 / r3, position);
  820.         final Vector3D vectorARDot = new Vector3D(2 * a2 / mu, velocity);
  821.         fillHalfRow(1, vectorAR,    jacobian[0], 0);
  822.         fillHalfRow(1, vectorARDot, jacobian[0], 3);

  823.         // dEx
  824.         final double d1 = -a * ratio / r3;
  825.         final double d2 = (hy * xDot - hx * yDot) / (sqrtMuA * epsilon);
  826.         final double d3 = (hx * y - hy * x) / sqrtMuA;
  827.         final Vector3D vectorExRDot =
  828.             new Vector3D((2 * x * yDot - xDot * y) / mu, g, -y * yDot / mu, f, -ey * d3 / epsilon, w);
  829.         fillHalfRow(ex * d1, position, -ey * d2, w, epsilon / sqrtMuA, drDotSdEy, jacobian[1], 0);
  830.         fillHalfRow(1, vectorExRDot, jacobian[1], 3);

  831.         // dEy
  832.         final Vector3D vectorEyRDot =
  833.             new Vector3D((2 * xDot * y - x * yDot) / mu, f, -x * xDot / mu, g, ex * d3 / epsilon, w);
  834.         fillHalfRow(ey * d1, position, ex * d2, w, -epsilon / sqrtMuA, drDotSdEx, jacobian[2], 0);
  835.         fillHalfRow(1, vectorEyRDot, jacobian[2], 3);

  836.         // dHx
  837.         final double h = (1 + hx2 + hy2) / (2 * sqrtMuA * epsilon);
  838.         fillHalfRow(-h * xDot, w, jacobian[3], 0);
  839.         fillHalfRow( h * x,    w, jacobian[3], 3);

  840.         // dHy
  841.         fillHalfRow(-h * yDot, w, jacobian[4], 0);
  842.         fillHalfRow( h * y,    w, jacobian[4], 3);

  843.         // dLambdaM
  844.         final double l = -ratio / sqrtMuA;
  845.         fillHalfRow(-1 / sqrtMuA, velocity, d2, w, l * ex, drDotSdEx, l * ey, drDotSdEy, jacobian[5], 0);
  846.         fillHalfRow(-2 / sqrtMuA, position, ex * beta, vectorEyRDot, -ey * beta, vectorExRDot, d3, w, jacobian[5], 3);

  847.         return jacobian;

  848.     }

  849.     /** {@inheritDoc} */
  850.     @Override
  851.     protected double[][] computeJacobianEccentricWrtCartesian() {

  852.         // start by computing the Jacobian with mean angle
  853.         final double[][] jacobian = computeJacobianMeanWrtCartesian();

  854.         // Differentiating the Kepler equation lM = lE - ex sin lE + ey cos lE leads to:
  855.         // dlM = (1 - ex cos lE - ey sin lE) dE - sin lE dex + cos lE dey
  856.         // which is inverted and rewritten as:
  857.         // dlE = a/r dlM + sin lE a/r dex - cos lE a/r dey
  858.         final SinCos scLe  = FastMath.sinCos(getLE());
  859.         final double cosLe = scLe.cos();
  860.         final double sinLe = scLe.sin();
  861.         final double aOr   = 1 / (1 - ex * cosLe - ey * sinLe);

  862.         // update longitude row
  863.         final double[] rowEx = jacobian[1];
  864.         final double[] rowEy = jacobian[2];
  865.         final double[] rowL  = jacobian[5];
  866.         for (int j = 0; j < 6; ++j) {
  867.             rowL[j] = aOr * (rowL[j] + sinLe * rowEx[j] - cosLe * rowEy[j]);
  868.         }

  869.         return jacobian;

  870.     }

  871.     /** {@inheritDoc} */
  872.     @Override
  873.     protected double[][] computeJacobianTrueWrtCartesian() {

  874.         // start by computing the Jacobian with eccentric angle
  875.         final double[][] jacobian = computeJacobianEccentricWrtCartesian();

  876.         // Differentiating the eccentric longitude equation
  877.         // tan((lv - lE)/2) = [ex sin lE - ey cos lE] / [sqrt(1-ex^2-ey^2) + 1 - ex cos lE - ey sin lE]
  878.         // leads to
  879.         // cT (dlv - dlE) = cE dlE + cX dex + cY dey
  880.         // with
  881.         // cT = [d^2 + (ex sin lE - ey cos lE)^2] / 2
  882.         // d  = 1 + sqrt(1-ex^2-ey^2) - ex cos lE - ey sin lE
  883.         // cE = (ex cos lE + ey sin lE) (sqrt(1-ex^2-ey^2) + 1) - ex^2 - ey^2
  884.         // cX =  sin lE (sqrt(1-ex^2-ey^2) + 1) - ey + ex (ex sin lE - ey cos lE) / sqrt(1-ex^2-ey^2)
  885.         // cY = -cos lE (sqrt(1-ex^2-ey^2) + 1) + ex + ey (ex sin lE - ey cos lE) / sqrt(1-ex^2-ey^2)
  886.         // which can be solved to find the differential of the true longitude
  887.         // dlv = (cT + cE) / cT dlE + cX / cT deX + cY / cT deX
  888.         final SinCos scLe      = FastMath.sinCos(getLE());
  889.         final double cosLe     = scLe.cos();
  890.         final double sinLe     = scLe.sin();
  891.         final double eSinE     = ex * sinLe - ey * cosLe;
  892.         final double ecosE     = ex * cosLe + ey * sinLe;
  893.         final double e2        = ex * ex + ey * ey;
  894.         final double epsilon   = FastMath.sqrt(1 - e2);
  895.         final double onePeps   = 1 + epsilon;
  896.         final double d         = onePeps - ecosE;
  897.         final double cT        = (d * d + eSinE * eSinE) / 2;
  898.         final double cE        = ecosE * onePeps - e2;
  899.         final double cX        = ex * eSinE / epsilon - ey + sinLe * onePeps;
  900.         final double cY        = ey * eSinE / epsilon + ex - cosLe * onePeps;
  901.         final double factorLe  = (cT + cE) / cT;
  902.         final double factorEx  = cX / cT;
  903.         final double factorEy  = cY / cT;

  904.         // update longitude row
  905.         final double[] rowEx = jacobian[1];
  906.         final double[] rowEy = jacobian[2];
  907.         final double[] rowL = jacobian[5];
  908.         for (int j = 0; j < 6; ++j) {
  909.             rowL[j] = factorLe * rowL[j] + factorEx * rowEx[j] + factorEy * rowEy[j];
  910.         }

  911.         return jacobian;

  912.     }

  913.     /** {@inheritDoc} */
  914.     @Override
  915.     public void addKeplerContribution(final PositionAngleType type, final double gm,
  916.                                       final double[] pDot) {
  917.         pDot[5] += computeKeplerianLDot(type, a, ex, ey, gm, cachedL, cachedPositionAngleType);
  918.     }

  919.     /**
  920.      * Compute rate of argument of longitude.
  921.      * @param type position angle type of rate
  922.      * @param a semi major axis
  923.      * @param ex ex
  924.      * @param ey ey
  925.      * @param mu mu
  926.      * @param l argument of longitude
  927.      * @param cachedType position angle type of passed l
  928.      * @return first-order time derivative for l
  929.      * @since 12.2
  930.      */
  931.     private static double computeKeplerianLDot(final PositionAngleType type, final double a, final double ex,
  932.                                                final double ey, final double mu,
  933.                                                final double l, final PositionAngleType cachedType) {
  934.         final double n  = FastMath.sqrt(mu / a) / a;
  935.         if (type == PositionAngleType.MEAN) {
  936.             return n;
  937.         }
  938.         final double oMe2;
  939.         final double ksi;
  940.         final SinCos sc;
  941.         if (type == PositionAngleType.ECCENTRIC) {
  942.             sc = FastMath.sinCos(EquinoctialLongitudeArgumentUtility.convertL(cachedType, l, ex, ey, type));
  943.             ksi  = 1. / (1 - ex * sc.cos() - ey * sc.sin());
  944.             return n * ksi;
  945.         } else { // TRUE
  946.             sc = FastMath.sinCos(EquinoctialLongitudeArgumentUtility.convertL(cachedType, l, ex, ey, type));
  947.             oMe2 = 1 - ex * ex - ey * ey;
  948.             ksi  = 1 + ex * sc.cos() + ey * sc.sin();
  949.             return n * ksi * ksi / (oMe2 * FastMath.sqrt(oMe2));
  950.         }
  951.     }

  952.     /**  Returns a string representation of this equinoctial parameters object.
  953.      * @return a string representation of this object
  954.      */
  955.     public String toString() {
  956.         return new StringBuilder().append("equinoctial parameters: ").append('{').
  957.                                   append("a: ").append(a).
  958.                                   append("; ex: ").append(ex).append("; ey: ").append(ey).
  959.                                   append("; hx: ").append(hx).append("; hy: ").append(hy).
  960.                                   append("; lv: ").append(FastMath.toDegrees(getLv())).
  961.                                   append(";}").toString();
  962.     }

  963.     /** {@inheritDoc} */
  964.     @Override
  965.     public PositionAngleType getCachedPositionAngleType() {
  966.         return cachedPositionAngleType;
  967.     }

  968.     /** {@inheritDoc} */
  969.     @Override
  970.     public boolean hasNonKeplerianRates() {
  971.         return hasNonKeplerianAcceleration();
  972.     }

  973.     /** {@inheritDoc} */
  974.     @Override
  975.     public EquinoctialOrbit withKeplerianRates() {
  976.         final PositionAngleType positionAngleType = getCachedPositionAngleType();
  977.         return new EquinoctialOrbit(getA(), getEquinoctialEx(), getEquinoctialEy(), getHx(), getHy(),
  978.                 getL(positionAngleType), positionAngleType, getFrame(), getDate(), getMu());
  979.     }

  980. }