StateCovariance.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;

  18. import org.hipparchus.analysis.differentiation.Gradient;
  19. import org.hipparchus.analysis.differentiation.GradientField;
  20. import org.hipparchus.linear.Array2DRowRealMatrix;
  21. import org.hipparchus.linear.MatrixUtils;
  22. import org.hipparchus.linear.RealMatrix;
  23. import org.hipparchus.util.MathArrays;
  24. import org.orekit.errors.OrekitException;
  25. import org.orekit.errors.OrekitMessages;
  26. import org.orekit.frames.Frame;
  27. import org.orekit.frames.KinematicTransform;
  28. import org.orekit.frames.LOF;
  29. import org.orekit.orbits.FieldOrbit;
  30. import org.orekit.orbits.Orbit;
  31. import org.orekit.orbits.OrbitType;
  32. import org.orekit.orbits.PositionAngleType;
  33. import org.orekit.time.AbsoluteDate;
  34. import org.orekit.time.FieldAbsoluteDate;
  35. import org.orekit.time.TimeStamped;
  36. import org.orekit.utils.CartesianCovarianceUtils;

  37. /** This class is the representation of a covariance matrix at a given date.
  38.  * <p>
  39.  * Currently, the covariance only represents the orbital elements.
  40.  * <p>
  41.  * It is possible to change the covariance frame by using the
  42.  * {@link #changeCovarianceFrame(Orbit, Frame)} or {@link #changeCovarianceFrame(Orbit, LOF)} method.
  43.  * These methods are based on Equations (18) and (20) of <i>Covariance Transformations for Satellite
  44.  * Flight Dynamics Operations</i> by David A. SVallado.
  45.  * <p>
  46.  * Finally, covariance orbit type can be changed using the
  47.  * {@link #changeCovarianceType(Orbit, OrbitType, PositionAngleType)} method.
  48.  *
  49.  * @author Bryan Cazabonne
  50.  * @author Vincent Cucchietti
  51.  * @since 11.3
  52.  */
  53. public class StateCovariance implements TimeStamped {

  54.     /** State dimension. */
  55.     public static final int STATE_DIMENSION = 6;

  56.     /** Default position angle for covariance expressed in Cartesian elements. */
  57.     private static final PositionAngleType DEFAULT_POSITION_ANGLE = PositionAngleType.TRUE;

  58.     /** Orbital covariance [6x6]. */
  59.     private final RealMatrix orbitalCovariance;

  60.     /** Covariance frame (can be null if LOF is defined). */
  61.     private final Frame frame;

  62.     /** Covariance LOF type (can be null if frame is defined). */
  63.     private final LOF lof;

  64.     /** Covariance epoch. */
  65.     private final AbsoluteDate epoch;

  66.     /** Covariance orbit type. */
  67.     private final OrbitType orbitType;

  68.     /** Covariance position angle type (not used if orbitType is CARTESIAN). */
  69.     private final PositionAngleType angleType;

  70.     /**
  71.      * Constructor.
  72.      * @param orbitalCovariance 6x6 orbital parameters covariance
  73.      * @param epoch epoch of the covariance
  74.      * @param lof covariance LOF type
  75.      */
  76.     public StateCovariance(final RealMatrix orbitalCovariance, final AbsoluteDate epoch, final LOF lof) {
  77.         this(orbitalCovariance, epoch, null, lof, OrbitType.CARTESIAN, DEFAULT_POSITION_ANGLE);
  78.     }

  79.     /**
  80.      * Constructor.
  81.      * @param orbitalCovariance 6x6 orbital parameters covariance
  82.      * @param epoch epoch of the covariance
  83.      * @param covarianceFrame covariance frame (inertial or Earth fixed)
  84.      * @param orbitType orbit type of the covariance (CARTESIAN if covarianceFrame is not pseudo-inertial)
  85.      * @param angleType position angle type of the covariance (not used if orbitType is CARTESIAN)
  86.      */
  87.     public StateCovariance(final RealMatrix orbitalCovariance, final AbsoluteDate epoch,
  88.                            final Frame covarianceFrame,
  89.                            final OrbitType orbitType, final PositionAngleType angleType) {
  90.         this(orbitalCovariance, epoch, covarianceFrame, null, orbitType, angleType);
  91.     }

  92.     /**
  93.      * Private constructor.
  94.      * @param orbitalCovariance 6x6 orbital parameters covariance
  95.      * @param epoch epoch of the covariance
  96.      * @param covarianceFrame covariance frame (inertial or Earth fixed)
  97.      * @param lof covariance LOF type
  98.      * @param orbitType orbit type of the covariance
  99.      * @param angleType position angle type of the covariance (not used if orbitType is CARTESIAN)
  100.      */
  101.     private StateCovariance(final RealMatrix orbitalCovariance, final AbsoluteDate epoch,
  102.                             final Frame covarianceFrame, final LOF lof,
  103.                             final OrbitType orbitType, final PositionAngleType angleType) {

  104.         checkFrameAndOrbitTypeConsistency(covarianceFrame, orbitType);

  105.         this.orbitalCovariance = orbitalCovariance;
  106.         this.epoch = epoch;
  107.         this.frame     = covarianceFrame;
  108.         this.lof       = lof;
  109.         this.orbitType = orbitType;
  110.         this.angleType = angleType;

  111.     }

  112.     /**
  113.      * Check constructor's inputs consistency.
  114.      *
  115.      * @param covarianceFrame covariance frame (inertial or Earth fixed)
  116.      * @param inputType orbit type of the covariance
  117.      *
  118.      * @throws OrekitException if input frame is not pseudo-inertial AND the orbit type is not Cartesian
  119.      */
  120.     public static void checkFrameAndOrbitTypeConsistency(final Frame covarianceFrame, final OrbitType inputType) {

  121.         // State covariance expressed in a celestial body frame
  122.         if (covarianceFrame != null) {

  123.             // Input frame is not pseudo-inertial
  124.             if (!covarianceFrame.isPseudoInertial() && inputType != OrbitType.CARTESIAN) {
  125.                 throw new OrekitException(OrekitMessages.WRONG_ORBIT_PARAMETERS_TYPE,
  126.                                           inputType.name(),
  127.                                           OrbitType.CARTESIAN.name());
  128.             }
  129.         }
  130.     }

  131.     /**
  132.      * Checks if input/output orbit and angle types are identical.
  133.      *
  134.      * @param inOrbitType input orbit type
  135.      * @param inAngleType input angle type
  136.      * @param outOrbitType output orbit type
  137.      * @param outAngleType output angle type
  138.      * @return flag defining if input/output orbit and angle types are identical
  139.      */
  140.     public static boolean inputAndOutputAreIdentical(final OrbitType inOrbitType, final PositionAngleType inAngleType,
  141.                                                      final OrbitType outOrbitType, final PositionAngleType outAngleType) {
  142.         return inOrbitType == outOrbitType && inAngleType == outAngleType;
  143.     }

  144.     /**
  145.      * Checks if input and output orbit types are both {@code OrbitType.CARTESIAN}.
  146.      *
  147.      * @param inOrbitType input orbit type
  148.      * @param outOrbitType output orbit type
  149.      * @return flag defining if input and output orbit types are both {@code OrbitType.CARTESIAN}
  150.      */
  151.     public static boolean inputAndOutputOrbitTypesAreCartesian(final OrbitType inOrbitType, final OrbitType outOrbitType) {
  152.         return inOrbitType == OrbitType.CARTESIAN && outOrbitType == OrbitType.CARTESIAN;
  153.     }

  154.     /** {@inheritDoc}. */
  155.     @Override
  156.     public AbsoluteDate getDate() {
  157.         return epoch;
  158.     }

  159.     /**
  160.      * Get the covariance matrix.
  161.      * @return the covariance matrix
  162.      */
  163.     public RealMatrix getMatrix() {
  164.         return orbitalCovariance;
  165.     }

  166.     /**
  167.      * Get the covariance orbit type.
  168.      * @return the covariance orbit type
  169.      */
  170.     public OrbitType getOrbitType() {
  171.         return orbitType;
  172.     }

  173.     /**
  174.      * Get the covariance angle type.
  175.      * @return the covariance angle type
  176.      */
  177.     public PositionAngleType getPositionAngleType() {
  178.         return angleType;
  179.     }

  180.     /**
  181.      * Get the covariance frame.
  182.      * @return the covariance frame (can be null)
  183.      * @see #getLOF()
  184.      */
  185.     public Frame getFrame() {
  186.         return frame;
  187.     }

  188.     /**
  189.      * Get the covariance LOF type.
  190.      * @return the covariance LOF type (can be null)
  191.      * @see #getFrame()
  192.      */
  193.     public LOF getLOF() {
  194.         return lof;
  195.     }

  196.     /**
  197.      * Get the covariance matrix in another orbit type.
  198.      * <p>
  199.      * The covariance orbit type <b>cannot</b> be changed if the covariance
  200.      * matrix is expressed in a {@link LOF local orbital frame} or a
  201.      * non-pseudo inertial frame.
  202.      * <p>
  203.      * As this type change uses the jacobian matrix of the transformation, it introduces a linear approximation.
  204.      * Hence, the current covariance matrix <b>will not exactly match</b> the new linearized case and the
  205.      * distribution will not follow a generalized Gaussian distribution anymore.
  206.      * <p>
  207.      * This is based on equation (1) to (6) from "Vallado, D. A. (2004). Covariance transformations for satellite flight
  208.      * dynamics operations."
  209.      *
  210.      * @param orbit orbit to which the covariance matrix should correspond
  211.      * @param outOrbitType target orbit type of the state covariance matrix
  212.      * @param outAngleType target position angle type of the state covariance matrix
  213.      * @return a new covariance state, expressed in the target orbit type with the target position angle
  214.      * @see #changeCovarianceFrame(Orbit, Frame)
  215.      */
  216.     public StateCovariance changeCovarianceType(final Orbit orbit, final OrbitType outOrbitType,
  217.                                                 final PositionAngleType outAngleType) {

  218.         // Handle case where the covariance is already expressed in the output type
  219.         if (outOrbitType == orbitType && (outOrbitType == OrbitType.CARTESIAN || outAngleType == angleType)) {
  220.             if (lof == null) {
  221.                 return new StateCovariance(orbitalCovariance, epoch, frame, orbitType, angleType);
  222.             }
  223.             else {
  224.                 return new StateCovariance(orbitalCovariance, epoch, lof);
  225.             }
  226.         }

  227.         // Check if the covariance is expressed in a celestial body frame
  228.         if (frame != null) {

  229.             // Check if the covariance is defined in an inertial frame
  230.             if (frame.isPseudoInertial()) {
  231.                 return changeTypeAndCreate(orbit, epoch, frame, orbitType, angleType, outOrbitType, outAngleType,
  232.                                            orbitalCovariance);
  233.             }

  234.             // The covariance is not defined in an inertial frame. The orbit type cannot be changed
  235.             throw new OrekitException(OrekitMessages.CANNOT_CHANGE_COVARIANCE_TYPE_IF_DEFINED_IN_NON_INERTIAL_FRAME);

  236.         }

  237.         // The covariance is not expressed in a celestial body frame. The orbit type cannot be changed
  238.         throw new OrekitException(OrekitMessages.CANNOT_CHANGE_COVARIANCE_TYPE_IF_DEFINED_IN_LOF);

  239.     }

  240.     /**
  241.      * Get the covariance in a given local orbital frame.
  242.      * <p>
  243.      * Changing the covariance frame is a linear process, this method does not introduce approximation unless a change
  244.      * in covariance orbit type is required.
  245.      * <p>
  246.      * This is based on equation (18) to (20) "from Vallado, D. A. (2004). Covariance transformations for satellite
  247.      * flight dynamics operations."
  248.      *
  249.      * @param orbit orbit to which the covariance matrix should correspond
  250.      * @param lofOut output local orbital frame
  251.      * @return a new covariance state, expressed in the output local orbital frame
  252.      */
  253.     public StateCovariance changeCovarianceFrame(final Orbit orbit, final LOF lofOut) {

  254.         // Verify current covariance frame
  255.         if (lof != null) {

  256.             // Check specific case where output covariance will be the same
  257.             if (lofOut == lof) {
  258.                 return new StateCovariance(orbitalCovariance, epoch, lof);
  259.             }

  260.             // Change the covariance local orbital frame
  261.             return changeFrameAndCreate(orbit, epoch, lof, lofOut, orbitalCovariance);

  262.         } else {

  263.             // Covariance is expressed in celestial body frame
  264.             return changeFrameAndCreate(orbit, epoch, frame, lofOut, orbitalCovariance, orbitType, angleType);

  265.         }

  266.     }

  267.     /**
  268.      * Get the covariance in the output frame.
  269.      * <p>
  270.      * Changing the covariance frame is a linear process, this method does not introduce approximation unless a change
  271.      * in covariance orbit type is required.
  272.      * <p>
  273.      * This is based on equation (18) to (20) "from Vallado, D. A. (2004). Covariance transformations for satellite
  274.      * flight dynamics operations."
  275.      *
  276.      * @param orbit orbit to which the covariance matrix should correspond
  277.      * @param frameOut output frame
  278.      * @return a new covariance state, expressed in the output frame
  279.      */
  280.     public StateCovariance changeCovarianceFrame(final Orbit orbit, final Frame frameOut) {

  281.         // Verify current covariance frame
  282.         if (lof != null) {

  283.             // Covariance is expressed in local orbital frame
  284.             return changeFrameAndCreate(orbit, epoch, lof, frameOut, orbitalCovariance);

  285.         } else {

  286.             // Check specific case where output covariance will be the same
  287.             if (frame == frameOut) {
  288.                 return new StateCovariance(orbitalCovariance, epoch, frame, orbitType, angleType);
  289.             }

  290.             // Change covariance frame
  291.             return changeFrameAndCreate(orbit, epoch, frame, frameOut, orbitalCovariance, orbitType, angleType);

  292.         }

  293.     }

  294.     /**
  295.      * Get a time-shifted covariance matrix.
  296.      * <p>
  297.      * The shifting model is a linearized, Keplerian one. In other words, it is based on a state transition matrix that
  298.      * is computed assuming Keplerian motion.
  299.      * <p>
  300.      * Shifting is <em>not</em> intended as a replacement for proper covariance propagation, but should be sufficient
  301.      * for small time shifts or coarse accuracy.
  302.      *
  303.      * @param orbit orbit to which the covariance matrix should correspond
  304.      * @param dt time shift in seconds
  305.      * @return a new covariance state, shifted with respect to the instance
  306.      */
  307.     public StateCovariance shiftedBy(final Orbit orbit, final double dt) {

  308.         // Shifted orbit
  309.         final Orbit shifted = orbit.shiftedBy(dt);

  310.         // State covariance expressed in celestial body frame
  311.         if (frame != null) {

  312.             // State covariance expressed in a pseudo-inertial frame
  313.             if (frame.isPseudoInertial()) {

  314.                 // Shift covariance by applying the STM
  315.                 final RealMatrix stm = getKeplerianStm(orbit, dt);
  316.                 return new StateCovariance(stm.multiply(getMatrix().multiplyTransposed(stm)), getDate().shiftedBy(dt),
  317.                         frame, orbitType, getPositionAngleType());

  318.             }

  319.             // State covariance expressed in a non pseudo-inertial frame
  320.             else {

  321.                 // Convert state covariance to orbit pseudo-inertial frame
  322.                 final StateCovariance inOrbitFrame = changeCovarianceFrame(orbit, orbit.getFrame());

  323.                 // Shift the state covariance
  324.                 final StateCovariance shiftedCovariance = inOrbitFrame.shiftedBy(orbit, dt);

  325.                 // Restore the initial covariance frame
  326.                 return shiftedCovariance.changeCovarianceFrame(shifted, frame);
  327.             }
  328.         }

  329.         // State covariance expressed in a commonly used local orbital frame (LOF)
  330.         else {

  331.             // Convert state covariance to orbit pseudo-inertial frame
  332.             final StateCovariance inOrbitFrame = changeCovarianceFrame(orbit, orbit.getFrame());

  333.             // Shift the state covariance
  334.             final StateCovariance shiftedCovariance = inOrbitFrame.shiftedBy(orbit, dt);

  335.             // Restore the initial covariance frame
  336.             return shiftedCovariance.changeCovarianceFrame(shifted, lof);
  337.         }

  338.     }

  339.     /**
  340.      * Create a covariance matrix in another orbit type.
  341.      * <p>
  342.      * As this type change uses the jacobian matrix of the transformation, it introduces a linear approximation.
  343.      * Hence, the input covariance matrix <b>will not exactly match</b> the new linearized case and the
  344.      * distribution will not follow a generalized Gaussian distribution anymore.
  345.      * <p>
  346.      * This is based on equation (1) to (6) from "Vallado, D. A. (2004). Covariance transformations for satellite flight
  347.      * dynamics operations."
  348.      *
  349.      * @param orbit orbit to which the covariance matrix should correspond
  350.      * @param date covariance epoch
  351.      * @param covFrame covariance frame
  352.      * @param inOrbitType initial orbit type of the state covariance matrix
  353.      * @param inAngleType initial position angle type of the state covariance matrix
  354.      * @param outOrbitType target orbit type of the state covariance matrix
  355.      * @param outAngleType target position angle type of the state covariance matrix
  356.      * @param inputCov input covariance
  357.      * @return the covariance expressed in the target orbit type with the target position angle
  358.      */
  359.     private static StateCovariance changeTypeAndCreate(final Orbit orbit,
  360.                                                        final AbsoluteDate date,
  361.                                                        final Frame covFrame,
  362.                                                        final OrbitType inOrbitType, final PositionAngleType inAngleType,
  363.                                                        final OrbitType outOrbitType, final PositionAngleType outAngleType,
  364.                                                        final RealMatrix inputCov) {

  365.         // Check if type change is really necessary, if not then return input covariance
  366.         if (inputAndOutputOrbitTypesAreCartesian(inOrbitType, outOrbitType) ||
  367.             inputAndOutputAreIdentical(inOrbitType, inAngleType, outOrbitType, outAngleType)) {
  368.             return new StateCovariance(inputCov, date, covFrame, inOrbitType, inAngleType);
  369.         }

  370.         // Notations:
  371.         // I: Input orbit type
  372.         // O: Output orbit type
  373.         // C: Cartesian parameters

  374.         // Compute dOutputdCartesian
  375.         final double[][] aOC               = new double[STATE_DIMENSION][STATE_DIMENSION];
  376.         final Orbit      orbitInOutputType = outOrbitType.convertType(orbit);
  377.         orbitInOutputType.getJacobianWrtCartesian(outAngleType, aOC);
  378.         final RealMatrix dOdC = new Array2DRowRealMatrix(aOC, false);

  379.         // Compute dCartesiandInput
  380.         final double[][] aCI              = new double[STATE_DIMENSION][STATE_DIMENSION];
  381.         final Orbit      orbitInInputType = inOrbitType.convertType(orbit);
  382.         orbitInInputType.getJacobianWrtParameters(inAngleType, aCI);
  383.         final RealMatrix dCdI = new Array2DRowRealMatrix(aCI, false);

  384.         // Compute dOutputdInput
  385.         final RealMatrix dOdI = dOdC.multiply(dCdI);

  386.         // Conversion of the state covariance matrix in target orbit type
  387.         final RealMatrix outputCov = dOdI.multiply(inputCov.multiplyTransposed(dOdI));

  388.         // Return the converted covariance
  389.         return new StateCovariance(outputCov, date, covFrame, outOrbitType, outAngleType);

  390.     }

  391.     /**
  392.      * Create a covariance matrix from a {@link LOF local orbital frame} to another
  393.      * {@link LOF local orbital frame}.
  394.      * <p>
  395.      * Changing the covariance frame is a linear process, this method does not introduce approximation.
  396.      * <p>
  397.      * The transformation is based on Equation (18) to (20) of "Covariance Transformations for Satellite Flight Dynamics
  398.      * Operations" by David A. Vallado.
  399.      *
  400.      * @param orbit orbit to which the covariance matrix should correspond
  401.      * @param date covariance epoch
  402.      * @param lofIn the local orbital frame in which the input covariance matrix is expressed
  403.      * @param lofOut the target local orbital frame
  404.      * @param inputCartesianCov input covariance {@code CARTESIAN})
  405.      * @return the covariance matrix expressed in the target commonly used local orbital frame in Cartesian elements
  406.      */
  407.     private static StateCovariance changeFrameAndCreate(final Orbit orbit, final AbsoluteDate date,
  408.                                                         final LOF lofIn, final LOF lofOut,
  409.                                                         final RealMatrix inputCartesianCov) {

  410.         // Builds the matrix to perform covariance transformation
  411.         final RealMatrix jacobianFromLofInToLofOut =
  412.                 getJacobian(LOF.transformFromLOFInToLOFOut(lofIn, lofOut, date, orbit.getPVCoordinates()));

  413.         // Get the Cartesian covariance matrix converted to frameOut
  414.         final RealMatrix cartesianCovarianceOut =
  415.                 jacobianFromLofInToLofOut.multiply(inputCartesianCov.multiplyTransposed(jacobianFromLofInToLofOut));

  416.         // Output converted covariance
  417.         return new StateCovariance(cartesianCovarianceOut, date, lofOut);

  418.     }

  419.     /**
  420.      * Convert the covariance matrix from a {@link Frame frame} to a {@link LOF local orbital frame}.
  421.      * <p>
  422.      * Changing the covariance frame is a linear process, this method does not introduce approximation unless a change
  423.      * in covariance orbit type is required.
  424.      * <p>
  425.      * The transformation is based on Equation (18) to (20) of "Covariance Transformations for Satellite Flight Dynamics
  426.      * Operations" by David A. Vallado.
  427.      * <p>
  428.      * As the frame transformation must be performed with the covariance expressed in Cartesian elements, both the orbit
  429.      * and position angle types of the input covariance must be provided.
  430.      * <p>
  431.      * <b>The output covariance matrix will provided in Cartesian orbit type and not converted back to
  432.      * its original expression (if input different from Cartesian elements).</b>
  433.      *
  434.      * @param orbit orbit to which the covariance matrix should correspond
  435.      * @param date covariance epoch
  436.      * @param frameIn the frame in which the input covariance matrix is expressed. In case the frame is <b>not</b>
  437.      * pseudo-inertial, the input covariance matrix is expected to be expressed in <b>Cartesian elements</b>.
  438.      * @param lofOut the target local orbital frame
  439.      * @param inputCov input covariance
  440.      * @param covOrbitType orbit type of the covariance matrix (used if frameIn is pseudo-inertial)
  441.      * @param covAngleType position angle type of the covariance matrix (used if frameIn is pseudo-inertial) (not used
  442.      * if covOrbitType equals {@code CARTESIAN})
  443.      * @return the covariance matrix expressed in the target local orbital frame in Cartesian elements
  444.      * @throws OrekitException if input frame is <b>not</b> pseudo-inertial <b>and</b> the input covariance is
  445.      * <b>not</b> expressed in Cartesian elements.
  446.      */
  447.     private static StateCovariance changeFrameAndCreate(final Orbit orbit, final AbsoluteDate date,
  448.                                                         final Frame frameIn, final LOF lofOut,
  449.                                                         final RealMatrix inputCov,
  450.                                                         final OrbitType covOrbitType,
  451.                                                         final PositionAngleType covAngleType) {

  452.         // Input frame is inertial
  453.         if (frameIn.isPseudoInertial()) {

  454.             // Convert input matrix to Cartesian parameters in input frame
  455.             final RealMatrix cartesianCovarianceIn =
  456.                     changeTypeAndCreate(orbit, date, frameIn, covOrbitType, covAngleType,
  457.                                         OrbitType.CARTESIAN, PositionAngleType.MEAN,
  458.                                         inputCov).getMatrix();

  459.             // Builds the matrix to perform covariance transformation
  460.             final RealMatrix jacobianFromFrameInToLofOut =
  461.                             getJacobian(lofOut.transformFromInertial(date, orbit.getPVCoordinates(frameIn)));

  462.             // Get the Cartesian covariance matrix converted to frameOut
  463.             final RealMatrix cartesianCovarianceOut =
  464.                     jacobianFromFrameInToLofOut.multiply(cartesianCovarianceIn.multiplyTransposed(jacobianFromFrameInToLofOut));

  465.             // Return converted covariance matrix expressed in Cartesian elements
  466.             return new StateCovariance(cartesianCovarianceOut, date, lofOut);

  467.         }

  468.         // Input frame is not inertial so the covariance matrix is expected to be in Cartesian elements
  469.         else {

  470.             // Method checkInputConsistency already verify that input covariance is defined in CARTESIAN type
  471.             final Frame orbitInertialFrame = orbit.getFrame();

  472.             // Compute rotation matrix from frameIn to orbit inertial frame
  473.             final RealMatrix cartesianCovarianceInOrbitFrame =
  474.                    changeFrameAndCreate(orbit, date, frameIn, orbitInertialFrame, inputCov,
  475.                                          OrbitType.CARTESIAN, PositionAngleType.MEAN).getMatrix();

  476.             // Convert from orbit inertial frame to lofOut
  477.             return changeFrameAndCreate(orbit, date, orbitInertialFrame, lofOut, cartesianCovarianceInOrbitFrame,
  478.                                         OrbitType.CARTESIAN, PositionAngleType.MEAN);

  479.         }

  480.     }

  481.     /**
  482.      * Convert the covariance matrix from a {@link LOF  local orbital frame} to a {@link Frame frame}.
  483.      * <p>
  484.      * Changing the covariance frame is a linear process, this method does not introduce approximation.
  485.      * <p>
  486.      * The transformation is based on Equation (18) to (20) of "Covariance Transformations for Satellite Flight Dynamics
  487.      * Operations" by David A. Vallado.
  488.      * <p>
  489.      * The <u>input</u> covariance matrix is necessarily provided in <b>Cartesian orbit type</b>.
  490.      * <p>
  491.      * The <u>output</u> covariance matrix will be expressed in <b>Cartesian elements</b>.
  492.      *
  493.      * @param orbit orbit to which the covariance matrix should correspond
  494.      * @param date covariance epoch
  495.      * @param lofIn the local orbital frame in which the input covariance matrix is expressed
  496.      * @param frameOut the target frame
  497.      * @param inputCartesianCov input covariance ({@code CARTESIAN})
  498.      * @return the covariance matrix expressed in the target frame in Cartesian elements
  499.      */
  500.     private static StateCovariance changeFrameAndCreate(final Orbit orbit, final AbsoluteDate date,
  501.                                                         final LOF lofIn, final Frame frameOut,
  502.                                                         final RealMatrix inputCartesianCov) {

  503.         // Output frame is pseudo-inertial
  504.         if (frameOut.isPseudoInertial()) {

  505.             // Builds the matrix to perform covariance transformation
  506.             final RealMatrix jacobianFromLofInToFrameOut =
  507.                     getJacobian(lofIn.transformFromInertial(date, orbit.getPVCoordinates(frameOut)).getInverse());

  508.             // Transform covariance
  509.             final RealMatrix transformedCovariance =
  510.                     jacobianFromLofInToFrameOut.multiply(inputCartesianCov.multiplyTransposed(jacobianFromLofInToFrameOut));

  511.             // Get the Cartesian covariance matrix converted to frameOut
  512.             return new StateCovariance(transformedCovariance, date, frameOut, OrbitType.CARTESIAN,
  513.                                        DEFAULT_POSITION_ANGLE);

  514.         }

  515.         // Output frame is not pseudo-inertial
  516.         else {

  517.             // Builds the matrix to perform covariance transformation
  518.             final RealMatrix jacobianFromLofInToOrbitFrame =
  519.                     getJacobian(lofIn.transformFromInertial(date, orbit.getPVCoordinates()).getInverse());

  520.             // Get the Cartesian covariance matrix converted to orbit inertial frame
  521.             final RealMatrix cartesianCovarianceInOrbitFrame = jacobianFromLofInToOrbitFrame.multiply(
  522.                     inputCartesianCov.multiplyTransposed(jacobianFromLofInToOrbitFrame));

  523.             // Get the Cartesian covariance matrix converted to frameOut
  524.             return changeFrameAndCreate(orbit, date, orbit.getFrame(), frameOut, cartesianCovarianceInOrbitFrame,
  525.                                         OrbitType.CARTESIAN, PositionAngleType.MEAN);
  526.         }

  527.     }

  528.     /**
  529.      * Get the covariance matrix in another frame.
  530.      * <p>
  531.      * The transformation is based on Equation (18) of "Covariance Transformations for Satellite Flight Dynamics
  532.      * Operations" by David A. Vallado.
  533.      * <p>
  534.      * Changing the covariance frame is a linear process, this method does not introduce approximation unless a change
  535.      * in covariance orbit type is required.
  536.      * <p>
  537.      * As the frame transformation must be performed with the covariance expressed in Cartesian elements, both the orbit
  538.      * and position angle types of the input covariance must be provided.
  539.      * <p>
  540.      * In case the <u>input</u> frame is <b>not</b> pseudo-inertial, the <u>input</u> covariance matrix is necessarily
  541.      * expressed in <b>Cartesian elements</b>.
  542.      * <p>
  543.      * In case the <u>output</u> frame is <b>not</b> pseudo-inertial, the <u>output</u> covariance matrix will be
  544.      * expressed in <b>Cartesian elements</b>.
  545.      *
  546.      * @param orbit orbit to which the covariance matrix should correspond
  547.      * @param date covariance epoch
  548.      * @param frameIn the frame in which the input covariance matrix is expressed
  549.      * @param frameOut the target frame
  550.      * @param inputCov input covariance
  551.      * @param covOrbitType orbit type of the covariance matrix (used if frameIn is pseudo-inertial)
  552.      * @param covAngleType position angle type of the covariance matrix (used if frameIn is pseudo-inertial) (<b>not</b>
  553.      * used if covOrbitType equals {@code CARTESIAN})
  554.      * @return the covariance matrix expressed in the target frame
  555.      * @throws OrekitException if input frame is <b>not</b> pseudo-inertial <b>and</b> the input covariance is
  556.      * <b>not</b> expressed in Cartesian elements.
  557.      */
  558.     private static StateCovariance changeFrameAndCreate(final Orbit orbit, final AbsoluteDate date,
  559.                                                         final Frame frameIn, final Frame frameOut,
  560.                                                         final RealMatrix inputCov,
  561.                                                         final OrbitType covOrbitType,
  562.                                                         final PositionAngleType covAngleType) {

  563.         // Input frame pseudo-inertial
  564.         if (frameIn.isPseudoInertial()) {

  565.             // Convert input matrix to Cartesian parameters in input frame
  566.             final RealMatrix cartesianCovarianceIn =
  567.                     changeTypeAndCreate(orbit, date, frameIn, covOrbitType, covAngleType,
  568.                                         OrbitType.CARTESIAN, PositionAngleType.MEAN,
  569.                                         inputCov).getMatrix();

  570.             // Get the Cartesian covariance matrix converted to frameOut
  571.             final RealMatrix cartesianCovarianceOut = CartesianCovarianceUtils.changeReferenceFrame(frameIn,
  572.                     cartesianCovarianceIn, date, frameOut);

  573.             // Output frame is pseudo-inertial
  574.             if (frameOut.isPseudoInertial()) {

  575.                 // Convert output Cartesian matrix to initial orbit type and position angle
  576.                 return changeTypeAndCreate(orbit, date, frameOut, OrbitType.CARTESIAN, PositionAngleType.MEAN,
  577.                                            covOrbitType, covAngleType, cartesianCovarianceOut);

  578.             }

  579.             // Output frame is not pseudo-inertial
  580.             else {

  581.                 // Output Cartesian matrix
  582.                 return new StateCovariance(cartesianCovarianceOut, date, frameOut, OrbitType.CARTESIAN,
  583.                                            DEFAULT_POSITION_ANGLE);

  584.             }
  585.         }

  586.         // Input frame is not pseudo-inertial
  587.         else {

  588.             // Method checkInputConsistency already verify that input covariance is defined in CARTESIAN type

  589.             // Convert covariance matrix to frameOut
  590.             final RealMatrix covInFrameOut = CartesianCovarianceUtils.changeReferenceFrame(frameIn,
  591.                     inputCov, date, frameOut);

  592.             // Output the Cartesian covariance matrix converted to frameOut
  593.             return new StateCovariance(covInFrameOut, date, frameOut, OrbitType.CARTESIAN, DEFAULT_POSITION_ANGLE);

  594.         }

  595.     }

  596.     /**
  597.      * Builds the matrix to perform covariance frame transformation.
  598.      * @param transform input transformation
  599.      * @return the matrix to perform the covariance frame transformation
  600.      */
  601.     private static RealMatrix getJacobian(final KinematicTransform transform) {
  602.         return MatrixUtils.createRealMatrix(transform.getPVJacobian());
  603.     }

  604.     /**
  605.      * Get the state transition matrix considering Keplerian contribution only and assuming equinoctial elements with mean anomaly.
  606.      *
  607.      * @param initialOrbit orbit to which the initial covariance matrix should correspond
  608.      * @param dt time difference between the two orbits
  609.      * @return the state transition matrix used to shift the covariance matrix
  610.      * @deprecated since 13.1. If you must, do:
  611.      * final RealMatrix stm = MatrixUtils.createRealIdentityMatrix(STATE_DIMENSION);
  612.      * double contribution = initialOrbit.getMeanAnomalyDotWrtA() * dt;
  613.      * stm.setEntry(5, 0, contribution);
  614.      */
  615.     @Deprecated
  616.     public static RealMatrix getStm(final Orbit initialOrbit, final double dt) {
  617.         // initialize the STM
  618.         final RealMatrix stm = MatrixUtils.createRealIdentityMatrix(STATE_DIMENSION);

  619.         // State transition matrix using Keplerian contribution only
  620.         final double contribution = initialOrbit.getMeanAnomalyDotWrtA() * dt;
  621.         stm.setEntry(5, 0, contribution);

  622.         // Return
  623.         return stm;
  624.     }

  625.     /**
  626.      * Get the state transition matrix of Keplerian motion.
  627.      * The coordinates used are those of the covariance itself.
  628.      *
  629.      * @param initialOrbit orbit to which the initial covariance matrix should correspond
  630.      * @param dt time difference between the two orbits
  631.      * @return the state transition matrix used to shift the covariance matrix
  632.      * @since 13.1
  633.      */
  634.     RealMatrix getKeplerianStm(final Orbit initialOrbit, final double dt) {

  635.         // compute derivatives of Keplerian motion
  636.         final GradientField field = GradientField.getField(STATE_DIMENSION);
  637.         final double[] stateVector = new double[STATE_DIMENSION];
  638.         orbitType.mapOrbitToArray(initialOrbit, getPositionAngleType(), stateVector, null);
  639.         final Gradient[] fieldGradientStateVector = MathArrays.buildArray(field, STATE_DIMENSION);
  640.         for (int i = 0; i < STATE_DIMENSION; i++) {
  641.             fieldGradientStateVector[i] = Gradient.variable(STATE_DIMENSION, i, stateVector[i]);
  642.         }
  643.         final FieldAbsoluteDate<Gradient> fieldGradientDate = new FieldAbsoluteDate<>(field, initialOrbit.getDate());
  644.         final FieldOrbit<Gradient> fieldOrbit = orbitType.mapArrayToOrbit(fieldGradientStateVector, null, getPositionAngleType(),
  645.                 fieldGradientDate, field.getOne().newInstance(initialOrbit.getMu()), initialOrbit.getFrame());
  646.         final FieldOrbit<Gradient> shiftedOrbit = fieldOrbit.shiftedBy(dt);  // automatic differentiation
  647.         final Gradient[] gradient = MathArrays.buildArray(field, STATE_DIMENSION);
  648.         orbitType.mapOrbitToArray(shiftedOrbit, getPositionAngleType(), gradient, null);

  649.         // Return state transition matrix
  650.         final RealMatrix stm = MatrixUtils.createRealIdentityMatrix(STATE_DIMENSION);
  651.         for (int i = 0; i < gradient.length; i++) {
  652.             stm.setRow(i, gradient[i].getGradient());
  653.         }
  654.         return stm;

  655.     }

  656. }