GNSSPropagator.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.propagation.analytical.gnss;

  18. import org.hipparchus.Field;
  19. import org.hipparchus.analysis.differentiation.Gradient;
  20. import org.hipparchus.analysis.differentiation.GradientField;
  21. import org.hipparchus.analysis.differentiation.UnivariateDerivative2;
  22. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  23. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  24. import org.hipparchus.linear.MatrixUtils;
  25. import org.hipparchus.linear.QRDecomposition;
  26. import org.hipparchus.linear.RealMatrix;
  27. import org.hipparchus.linear.RealVector;
  28. import org.hipparchus.util.FastMath;
  29. import org.hipparchus.util.FieldSinCos;
  30. import org.hipparchus.util.SinCos;
  31. import org.orekit.attitudes.Attitude;
  32. import org.orekit.attitudes.AttitudeProvider;
  33. import org.orekit.frames.Frame;
  34. import org.orekit.orbits.CartesianOrbit;
  35. import org.orekit.orbits.FieldKeplerianAnomalyUtility;
  36. import org.orekit.orbits.KeplerianAnomalyUtility;
  37. import org.orekit.orbits.KeplerianOrbit;
  38. import org.orekit.orbits.Orbit;
  39. import org.orekit.orbits.PositionAngleType;
  40. import org.orekit.propagation.AbstractMatricesHarvester;
  41. import org.orekit.propagation.SpacecraftState;
  42. import org.orekit.propagation.analytical.AbstractAnalyticalPropagator;
  43. import org.orekit.propagation.analytical.gnss.data.FieldGnssOrbitalElements;
  44. import org.orekit.propagation.analytical.gnss.data.GNSSOrbitalElements;
  45. import org.orekit.time.AbsoluteDate;
  46. import org.orekit.utils.DoubleArrayDictionary;
  47. import org.orekit.utils.FieldPVCoordinates;
  48. import org.orekit.utils.PVCoordinates;
  49. import org.orekit.utils.ParameterDriver;
  50. import org.orekit.utils.TimeSpanMap.Span;

  51. import java.util.ArrayList;
  52. import java.util.Collections;
  53. import java.util.List;

  54. /** Common handling of {@link AbstractAnalyticalPropagator} methods for GNSS propagators.
  55.  * <p>
  56.  * This class allows to provide easily a subset of {@link AbstractAnalyticalPropagator} methods
  57.  * for specific GNSS propagators.
  58.  * </p>
  59.  * @author Pascal Parraud
  60.  */
  61. public class GNSSPropagator extends AbstractAnalyticalPropagator {

  62.     /** Maximum number of iterations for internal loops.
  63.      * @since 13.0
  64.      */
  65.     private static final int MAX_ITER = 100;

  66.     /** Tolerance on position for rebuilding orbital elements from initial state.
  67.      * @since 13.0
  68.      */
  69.     private static final double TOL_P = 1.0e-6;

  70.     /** Tolerance on velocity for rebuilding orbital elements from initial state.
  71.      * @since 13.0
  72.      */
  73.     private static final double TOL_V = 1.0e-9;

  74.     /** Number of free parameters for orbital elements.
  75.      * @since 13.0
  76.      */
  77.     private static final int FREE_PARAMETERS = 6;

  78.     /** Convergence parameter.
  79.      * @since 13.0
  80.      */
  81.     private static final double EPS = 1.0e-12;

  82.     /** The GNSS propagation model used. */
  83.     private GNSSOrbitalElements<?> orbitalElements;

  84.     /** The ECI frame used for GNSS propagation. */
  85.     private final Frame eci;

  86.     /** The ECEF frame used for GNSS propagation. */
  87.     private final Frame ecef;

  88.     /**
  89.      * Build a new instance.
  90.      * @param orbitalElements GNSS orbital elements
  91.      * @param eci Earth Centered Inertial frame
  92.      * @param ecef Earth Centered Earth Fixed frame
  93.      * @param provider attitude provider
  94.      * @param mass satellite mass (kg)
  95.      */
  96.     GNSSPropagator(final GNSSOrbitalElements<?> orbitalElements, final Frame eci,
  97.                    final Frame ecef, final AttitudeProvider provider, final double mass) {
  98.         super(provider);
  99.         // Stores the GNSS orbital elements
  100.         this.orbitalElements = orbitalElements;
  101.         // Sets the Earth Centered Inertial frame
  102.         this.eci  = eci;
  103.         // Sets the Earth Centered Earth Fixed frame
  104.         this.ecef = ecef;
  105.         // Sets initial state
  106.         final Orbit orbit = propagateOrbit(orbitalElements.getDate());
  107.         final Attitude attitude = provider.getAttitude(orbit, orbit.getDate(), orbit.getFrame());

  108.         // calling the method from base class because the one overridden below recomputes the orbital elements
  109.         super.resetInitialState(new SpacecraftState(orbit, attitude).withMass(mass));

  110.     }

  111.     /**
  112.      * Build a new instance from an initial state.
  113.      * <p>
  114.      * The Keplerian elements already present in the {@code nonKeplerianElements} argument
  115.      * will be ignored as it is the {@code initialState} argument that will be used to
  116.      * build the complete orbital elements of the propagator
  117.      * </p>
  118.      * @param initialState         initial state
  119.      * @param nonKeplerianElements non-Keplerian orbital elements (the Keplerian orbital elements will be ignored)
  120.      * @param ecef                 Earth Centered Earth Fixed frame
  121.      * @param provider             attitude provider
  122.      * @param mass                 spacecraft mass
  123.      * @since 13.0
  124.      */
  125.     GNSSPropagator(final SpacecraftState initialState, final GNSSOrbitalElements<?> nonKeplerianElements,
  126.                    final Frame ecef, final AttitudeProvider provider, final double mass) {
  127.         this(buildOrbitalElements(initialState, nonKeplerianElements, ecef, provider, mass),
  128.              initialState.getFrame(), ecef, provider, initialState.getMass());
  129.     }

  130.     /**
  131.      * Gets the Earth Centered Inertial frame used to propagate the orbit.
  132.      *
  133.      * @return the ECI frame
  134.      */
  135.     public Frame getECI() {
  136.         return eci;
  137.     }

  138.     /**
  139.      * Gets the Earth Centered Earth Fixed frame used to propagate GNSS orbits according to the
  140.      * Interface Control Document.
  141.      *
  142.      * @return the ECEF frame
  143.      */
  144.     public Frame getECEF() {
  145.         return ecef;
  146.     }

  147.     /**
  148.      * Gets the Earth gravity coefficient used for GNSS propagation.
  149.      *
  150.      * @return the Earth gravity coefficient.
  151.      */
  152.     public double getMU() {
  153.         return orbitalElements.getMu();
  154.     }

  155.     /** Get the underlying GNSS propagation orbital elements.
  156.      * @return the underlying GNSS orbital elements
  157.      * @since 13.0
  158.      */
  159.     public GNSSOrbitalElements<?> getOrbitalElements() {
  160.         return orbitalElements;
  161.     }

  162.     /** {@inheritDoc}
  163.      * @since 13.0
  164.      */
  165.     @Override
  166.     protected AbstractMatricesHarvester createHarvester(final String stmName, final RealMatrix initialStm,
  167.                                                         final DoubleArrayDictionary initialJacobianColumns) {

  168.         // Create the harvester
  169.         final GnssHarvester harvester = new GnssHarvester(this, stmName, initialStm, initialJacobianColumns);

  170.         // Update the list of additional state provider
  171.         addAdditionalDataProvider(harvester);

  172.         // Return the configured harvester
  173.         return harvester;

  174.     }

  175.     /** {@inheritDoc}
  176.      * @since 13.0
  177.      */
  178.     @Override
  179.     protected List<String> getJacobiansColumnsNames() {
  180.         final List<String> columnsNames = new ArrayList<>();
  181.         for (final ParameterDriver driver : orbitalElements.getParametersDrivers()) {
  182.             if (driver.isSelected() && !columnsNames.contains(driver.getNamesSpanMap().getFirstSpan().getData())) {
  183.                 // As driver with same name should have same NamesSpanMap we only check if the first span is present,
  184.                 // if not we add all span names to columnsNames
  185.                 for (Span<String> span = driver.getNamesSpanMap().getFirstSpan(); span != null; span = span.next()) {
  186.                     columnsNames.add(span.getData());
  187.                 }
  188.             }
  189.         }
  190.         Collections.sort(columnsNames);
  191.         return columnsNames;
  192.     }

  193.     /** {@inheritDoc} */
  194.     @Override
  195.     public Orbit propagateOrbit(final AbsoluteDate date) {
  196.         // Gets the PVCoordinates in ECEF frame
  197.         final PVCoordinates pvaInECEF = propagateInEcef(date);
  198.         // Transforms the PVCoordinates to ECI frame
  199.         final PVCoordinates pvaInECI = ecef.getTransformTo(eci, date).transformPVCoordinates(pvaInECEF);
  200.         // Returns the Cartesian orbit
  201.         return new CartesianOrbit(pvaInECI, eci, date, getMU());
  202.     }

  203.     /**
  204.      * Gets the PVCoordinates of the GNSS SV in {@link #getECEF() ECEF frame}.
  205.      *
  206.      * <p>The algorithm uses automatic differentiation to compute velocity and
  207.      * acceleration.</p>
  208.      *
  209.      * @param date the computation date
  210.      * @return the GNSS SV PVCoordinates in {@link #getECEF() ECEF frame}
  211.      */
  212.     public PVCoordinates propagateInEcef(final AbsoluteDate date) {
  213.         // Duration from GNSS ephemeris Reference date
  214.         final UnivariateDerivative2 tk = new UnivariateDerivative2(getTk(date), 1.0, 0.0);
  215.         // Semi-major axis
  216.         final UnivariateDerivative2 ak = tk.multiply(orbitalElements.getADot()).add(orbitalElements.getSma());
  217.         // Mean motion
  218.         final UnivariateDerivative2 nA = tk.multiply(orbitalElements.getDeltaN0Dot() * 0.5).
  219.                                          add(orbitalElements.getDeltaN0()).
  220.                                          add(orbitalElements.getMeanMotion0());
  221.         // Mean anomaly
  222.         final UnivariateDerivative2 mk = tk.multiply(nA).add(orbitalElements.getM0());
  223.         // Eccentric Anomaly
  224.         final UnivariateDerivative2 e  = tk.newInstance(orbitalElements.getE());
  225.         final UnivariateDerivative2 ek = FieldKeplerianAnomalyUtility.ellipticMeanToEccentric(e, mk);
  226.         // True Anomaly
  227.         final UnivariateDerivative2 vk =  FieldKeplerianAnomalyUtility.ellipticEccentricToTrue(e, ek);
  228.         // Argument of Latitude
  229.         final UnivariateDerivative2 phik    = vk.add(orbitalElements.getPa());
  230.         final FieldSinCos<UnivariateDerivative2> cs2phi = FastMath.sinCos(phik.multiply(2));
  231.         // Argument of Latitude Correction
  232.         final UnivariateDerivative2 dphik = cs2phi.cos().multiply(orbitalElements.getCuc()).add(cs2phi.sin().multiply(orbitalElements.getCus()));
  233.         // Radius Correction
  234.         final UnivariateDerivative2 drk = cs2phi.cos().multiply(orbitalElements.getCrc()).add(cs2phi.sin().multiply(orbitalElements.getCrs()));
  235.         // Inclination Correction
  236.         final UnivariateDerivative2 dik = cs2phi.cos().multiply(orbitalElements.getCic()).add(cs2phi.sin().multiply(orbitalElements.getCis()));
  237.         // Corrected Argument of Latitude
  238.         final FieldSinCos<UnivariateDerivative2> csuk = FastMath.sinCos(phik.add(dphik));
  239.         // Corrected Radius
  240.         final UnivariateDerivative2 rk = ek.cos().multiply(-orbitalElements.getE()).add(1).multiply(ak).add(drk);
  241.         // Corrected Inclination
  242.         final UnivariateDerivative2 ik  = tk.multiply(orbitalElements.getIDot()).add(orbitalElements.getI0()).add(dik);
  243.         final FieldSinCos<UnivariateDerivative2> csik = FastMath.sinCos(ik);
  244.         // Positions in orbital plane
  245.         final UnivariateDerivative2 xk = csuk.cos().multiply(rk);
  246.         final UnivariateDerivative2 yk = csuk.sin().multiply(rk);
  247.         // Corrected longitude of ascending node
  248.         final double thetaDot = orbitalElements.getAngularVelocity();
  249.         final FieldSinCos<UnivariateDerivative2> csomk =
  250.             FastMath.sinCos(tk.multiply(orbitalElements.getOmegaDot() - thetaDot).
  251.                             add(orbitalElements.getOmega0() - thetaDot * orbitalElements.getTime()));
  252.         // returns the Earth-fixed coordinates
  253.         final FieldVector3D<UnivariateDerivative2> positionWithDerivatives =
  254.                         new FieldVector3D<>(xk.multiply(csomk.cos()).subtract(yk.multiply(csomk.sin()).multiply(csik.cos())),
  255.                                             xk.multiply(csomk.sin()).add(yk.multiply(csomk.cos()).multiply(csik.cos())),
  256.                                             yk.multiply(csik.sin()));
  257.         return new PVCoordinates(positionWithDerivatives);
  258.     }

  259.     /**
  260.      * Gets the duration from GNSS Reference epoch.
  261.      * <p>This takes the GNSS week roll-over into account.</p>
  262.      * @param date the considered date
  263.      * @return the duration from GNSS orbit Reference epoch (s)
  264.      */
  265.     private double getTk(final AbsoluteDate date) {
  266.         final double cycleDuration = orbitalElements.getCycleDuration();
  267.         // Time from ephemeris reference epoch
  268.         double tk = date.durationFrom(orbitalElements.getDate());
  269.         // Adjusts the time to take roll over week into account
  270.         while (tk > 0.5 * cycleDuration) {
  271.             tk -= cycleDuration;
  272.         }
  273.         while (tk < -0.5 * cycleDuration) {
  274.             tk += cycleDuration;
  275.         }
  276.         // Returns the time from ephemeris reference epoch
  277.         return tk;
  278.     }

  279.     /** {@inheritDoc} */
  280.     @Override
  281.     public Frame getFrame() {
  282.         return eci;
  283.     }

  284.     /** {@inheritDoc} */
  285.     @Override
  286.     protected double getMass(final AbsoluteDate date) {
  287.         return getInitialState().getMass();
  288.     }

  289.     /** {@inheritDoc} */
  290.     @Override
  291.     public void resetInitialState(final SpacecraftState state) {
  292.         orbitalElements = buildOrbitalElements(state, orbitalElements, ecef, getAttitudeProvider(), state.getMass());
  293.         final Orbit orbit = propagateOrbit(orbitalElements.getDate());
  294.         final Attitude attitude = getAttitudeProvider().getAttitude(orbit, orbit.getDate(), orbit.getFrame());
  295.         super.resetInitialState(new SpacecraftState(orbit, attitude).withMass(state.getMass()));
  296.     }

  297.     /** {@inheritDoc} */
  298.     @Override
  299.     protected void resetIntermediateState(final SpacecraftState state, final boolean forward) {
  300.         resetInitialState(state);
  301.     }

  302.     /**
  303.      * Build orbital elements from initial state.
  304.      * <p>
  305.      * This method is roughly the inverse of {@link #propagateInEcef(AbsoluteDate)}, except it starts from a state in
  306.      * inertial frame
  307.      * </p>
  308.      *
  309.      * @param initialState    initial state
  310.      * @param nonKeplerianElements non-Keplerian orbital elements (the Keplerian orbital elements will be overridden)
  311.      * @param ecef            Earth Centered Earth Fixed frame
  312.      * @param provider        attitude provider
  313.      * @param mass            satellite mass (kg)
  314.      * @return orbital elements that generate the {@code initialState} when used with a propagator
  315.      * @since 13.0
  316.      */
  317.     private static GNSSOrbitalElements<?> buildOrbitalElements(final SpacecraftState initialState,
  318.                                                                final GNSSOrbitalElements<?> nonKeplerianElements,
  319.                                                                final Frame ecef, final AttitudeProvider provider,
  320.                                                                final double mass) {

  321.         // get approximate initial orbit
  322.         final Frame frozenEcef = ecef.getFrozenFrame(initialState.getFrame(), initialState.getDate(), "frozen");
  323.         final KeplerianOrbit orbit = approximateInitialOrbit(initialState, nonKeplerianElements, frozenEcef);

  324.         // refine orbit using simple differential correction to reach target PV
  325.         final PVCoordinates targetPV = initialState.getPVCoordinates();
  326.         final FieldGnssOrbitalElements<Gradient, ?> gElements = convert(nonKeplerianElements, orbit);
  327.         for (int i = 0; i < MAX_ITER; ++i) {

  328.             // get position-velocity derivatives with respect to initial orbit
  329.             final FieldGnssPropagator<Gradient> gPropagator =
  330.                 new FieldGnssPropagator<>(gElements, initialState.getFrame(), ecef, provider,
  331.                                           gElements.getMu().newInstance(mass));
  332.             final FieldPVCoordinates<Gradient> gPV = gPropagator.getInitialState().getPVCoordinates();

  333.             // compute Jacobian matrix
  334.             final RealMatrix jacobian = MatrixUtils.createRealMatrix(FREE_PARAMETERS, FREE_PARAMETERS);
  335.             jacobian.setRow(0, gPV.getPosition().getX().getGradient());
  336.             jacobian.setRow(1, gPV.getPosition().getY().getGradient());
  337.             jacobian.setRow(2, gPV.getPosition().getZ().getGradient());
  338.             jacobian.setRow(3, gPV.getVelocity().getX().getGradient());
  339.             jacobian.setRow(4, gPV.getVelocity().getY().getGradient());
  340.             jacobian.setRow(5, gPV.getVelocity().getZ().getGradient());

  341.             // linear correction to get closer to target PV
  342.             final RealVector residuals = MatrixUtils.createRealVector(FREE_PARAMETERS);
  343.             residuals.setEntry(0, targetPV.getPosition().getX() - gPV.getPosition().getX().getValue());
  344.             residuals.setEntry(1, targetPV.getPosition().getY() - gPV.getPosition().getY().getValue());
  345.             residuals.setEntry(2, targetPV.getPosition().getZ() - gPV.getPosition().getZ().getValue());
  346.             residuals.setEntry(3, targetPV.getVelocity().getX() - gPV.getVelocity().getX().getValue());
  347.             residuals.setEntry(4, targetPV.getVelocity().getY() - gPV.getVelocity().getY().getValue());
  348.             residuals.setEntry(5, targetPV.getVelocity().getZ() - gPV.getVelocity().getZ().getValue());
  349.             final RealVector correction = new QRDecomposition(jacobian, EPS).getSolver().solve(residuals);

  350.             // update initial orbit
  351.             gElements.setSma(gElements.getSma().add(correction.getEntry(0)));
  352.             gElements.setE(gElements.getE().add(correction.getEntry(1)));
  353.             gElements.setI0(gElements.getI0().add(correction.getEntry(2)));
  354.             gElements.setPa(gElements.getPa().add(correction.getEntry(3)));
  355.             gElements.setOmega0(gElements.getOmega0().add(correction.getEntry(4)));
  356.             gElements.setM0(gElements.getM0().add(correction.getEntry(5)));

  357.             final double deltaP = FastMath.sqrt(residuals.getEntry(0) * residuals.getEntry(0) +
  358.                                                 residuals.getEntry(1) * residuals.getEntry(1) +
  359.                                                 residuals.getEntry(2) * residuals.getEntry(2));
  360.             final double deltaV = FastMath.sqrt(residuals.getEntry(3) * residuals.getEntry(3) +
  361.                                                 residuals.getEntry(4) * residuals.getEntry(4) +
  362.                                                 residuals.getEntry(5) * residuals.getEntry(5));

  363.             if (deltaP < TOL_P && deltaV < TOL_V) {
  364.                 break;
  365.             }

  366.         }

  367.         return gElements.toNonField();

  368.     }

  369.     /** Compute approximate initial orbit.
  370.      * @param initialState         initial state
  371.      * @param nonKeplerianElements non-Keplerian orbital elements (the Keplerian orbital elements will be ignored)
  372.      * @param frozenEcef           inertial frame aligned with Earth Centered Earth Fixed frame at orbit date
  373.      * @return approximate initial orbit that generate a state close to {@code initialState}
  374.      * @since 13.0
  375.      */
  376.     private static KeplerianOrbit approximateInitialOrbit(final SpacecraftState initialState,
  377.                                                           final GNSSOrbitalElements<?> nonKeplerianElements,
  378.                                                           final Frame frozenEcef) {

  379.         // rotate the state to a frame that is inertial but aligned with Earth frame,
  380.         // as analytical model is expressed in Earth frame
  381.         final PVCoordinates pv = initialState.getPVCoordinates(frozenEcef);
  382.         final Vector3D      p  = pv.getPosition();
  383.         final Vector3D      v  = pv.getVelocity();

  384.         // compute Keplerian orbital parameters
  385.         final double   rk  = p.getNorm();

  386.         // compute orbital plane orientation
  387.         final Vector3D normal = pv.getMomentum().normalize();
  388.         final double   cosIk  = normal.getZ();
  389.         final double   ik     = Vector3D.angle(normal, Vector3D.PLUS_K);

  390.         // compute position in orbital plane
  391.         final double   q   = FastMath.hypot(normal.getX(), normal.getY());
  392.         final double   cos = -normal.getY() / q;
  393.         final double   sin =  normal.getX() / q;
  394.         final double   xk  =  p.getX() * cos + p.getY() * sin;
  395.         final double   yk  = (p.getY() * cos - p.getX() * sin) / cosIk;

  396.         // corrected latitude argument
  397.         final double   uk  = FastMath.atan2(yk, xk);

  398.         // recover latitude argument before correction, using a fixed-point method
  399.         double phi = uk;
  400.         for (int i = 0; i < MAX_ITER; ++i) {
  401.             final double previous = phi;
  402.             final SinCos cs2Phi = FastMath.sinCos(2 * phi);
  403.             phi = uk - (cs2Phi.cos() * nonKeplerianElements.getCuc() + cs2Phi.sin() * nonKeplerianElements.getCus());
  404.             if (FastMath.abs(phi - previous) <= EPS) {
  405.                 break;
  406.             }
  407.         }
  408.         final SinCos cs2phi = FastMath.sinCos(2 * phi);

  409.         // recover plane orientation before correction
  410.         // here, we know that tk = 0 since our orbital elements will be at initial state date
  411.         final double i0  = ik - (cs2phi.cos() * nonKeplerianElements.getCic() + cs2phi.sin() * nonKeplerianElements.getCis());
  412.         final double om0 = FastMath.atan2(sin, cos) +
  413.                            nonKeplerianElements.getAngularVelocity() * nonKeplerianElements.getTime();

  414.         // recover eccentricity and anomaly
  415.         final double mu = initialState.getOrbit().getMu();
  416.         final double rV2OMu           = rk * v.getNormSq() / mu;
  417.         final double sma              = rk / (2 - rV2OMu);
  418.         final double eCosE            = rV2OMu - 1;
  419.         final double eSinE            = Vector3D.dotProduct(p, v) / FastMath.sqrt(mu * sma);
  420.         final double e                = FastMath.hypot(eCosE, eSinE);
  421.         final double eccentricAnomaly = FastMath.atan2(eSinE, eCosE);
  422.         final double aop              = phi - eccentricAnomaly;
  423.         final double meanAnomaly      = KeplerianAnomalyUtility.ellipticEccentricToMean(e, eccentricAnomaly);

  424.         return new KeplerianOrbit(sma, e, i0, aop, om0, meanAnomaly, PositionAngleType.MEAN,
  425.                                   PositionAngleType.MEAN, frozenEcef,
  426.                                   initialState.getDate(), mu);

  427.     }

  428.     /** Convert orbital elements to gradient.
  429.      * @param elements   primitive double elements
  430.      * @param orbit      Keplerian orbit
  431.      * @return converted elements, set up as gradient relative to Keplerian orbit
  432.      * @since 13.0
  433.      */
  434.     private static FieldGnssOrbitalElements<Gradient, ?> convert(final GNSSOrbitalElements<?> elements,
  435.                                                                  final KeplerianOrbit orbit) {

  436.         final Field<Gradient> field = GradientField.getField(FREE_PARAMETERS);
  437.         final FieldGnssOrbitalElements<Gradient, ?> gElements = elements.toField(field);

  438.         // Keplerian parameters
  439.         gElements.setSma(Gradient.variable(FREE_PARAMETERS, 0, orbit.getA()));
  440.         gElements.setE(Gradient.variable(FREE_PARAMETERS, 1, orbit.getE()));
  441.         gElements.setI0(Gradient.variable(FREE_PARAMETERS, 2, orbit.getI()));
  442.         gElements.setPa(Gradient.variable(FREE_PARAMETERS, 3, orbit.getPerigeeArgument()));
  443.         gElements.setOmega0(Gradient.variable(FREE_PARAMETERS, 4, orbit.getRightAscensionOfAscendingNode()));
  444.         gElements.setM0(Gradient.variable(FREE_PARAMETERS, 5, orbit.getMeanAnomaly()));

  445.         return gElements;

  446.     }

  447. }