StateCovariance.java

  1. /* Copyright 2002-2022 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.linear.Array2DRowRealMatrix;
  19. import org.hipparchus.linear.MatrixUtils;
  20. import org.hipparchus.linear.RealMatrix;
  21. import org.hipparchus.util.FastMath;
  22. import org.orekit.errors.OrekitException;
  23. import org.orekit.errors.OrekitMessages;
  24. import org.orekit.frames.Frame;
  25. import org.orekit.frames.LOFType;
  26. import org.orekit.frames.Transform;
  27. import org.orekit.orbits.Orbit;
  28. import org.orekit.orbits.OrbitType;
  29. import org.orekit.orbits.PositionAngle;
  30. import org.orekit.time.AbsoluteDate;
  31. import org.orekit.time.TimeStamped;
  32. import org.orekit.utils.CartesianDerivativesFilter;

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

  50.     /** State dimension. */
  51.     private static final int STATE_DIMENSION = 6;

  52.     /** Default position angle for covariance expressed in Cartesian elements. */
  53.     private static final PositionAngle DEFAULT_POSITION_ANGLE = PositionAngle.TRUE;

  54.     /** Orbital covariance [6x6]. */
  55.     private final RealMatrix orbitalCovariance;

  56.     /** Covariance frame (can be null if lofType is defined). */
  57.     private final Frame frame;

  58.     /** Covariance LOF type (can be null if frame is defined). */
  59.     private final LOFType lofType;

  60.     /** Covariance epoch. */
  61.     private final AbsoluteDate epoch;

  62.     /** Covariance orbit type. */
  63.     private final OrbitType orbitType;

  64.     /** Covariance position angle type (not used if orbitType is CARTESIAN). */
  65.     private final PositionAngle angleType;

  66.     /**
  67.      * Constructor.
  68.      * @param orbitalCovariance 6x6 orbital parameters covariance
  69.      * @param epoch epoch of the covariance
  70.      * @param lofType covariance LOF type
  71.      */
  72.     public StateCovariance(final RealMatrix orbitalCovariance, final AbsoluteDate epoch, final LOFType lofType) {
  73.         this(orbitalCovariance, epoch, null, lofType, OrbitType.CARTESIAN, DEFAULT_POSITION_ANGLE);
  74.     }

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

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

  100.         checkInputConsistency(covarianceFrame, orbitType);

  101.         this.orbitalCovariance = orbitalCovariance;
  102.         this.epoch = epoch;
  103.         this.frame = covarianceFrame;
  104.         this.lofType = lofType;
  105.         this.orbitType = orbitType;
  106.         this.angleType = angleType;

  107.     }

  108.     /** Check constructor's inputs consistency.
  109.      *
  110.      * @param covarianceFrame covariance frame (inertial or Earth fixed)
  111.      * @param inputType orbit type of the covariance
  112.      */
  113.     private void checkInputConsistency(final Frame covarianceFrame, final OrbitType inputType) {

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

  116.             // Input frame is pseudo-inertial
  117.             if (!covarianceFrame.isPseudoInertial() && inputType != OrbitType.CARTESIAN) {
  118.                 throw new OrekitException(OrekitMessages.WRONG_ORBIT_PARAMETERS_TYPE,
  119.                                           inputType.name(),
  120.                                           OrbitType.CARTESIAN.name());
  121.             }

  122.         }
  123.     }

  124.     /** {@inheritDoc}. */
  125.     @Override
  126.     public AbsoluteDate getDate() {
  127.         return epoch;
  128.     }

  129.     /**
  130.      * Get the covariance matrix.
  131.      * @return the covariance matrix
  132.      */
  133.     public RealMatrix getMatrix() {
  134.         return orbitalCovariance;
  135.     }

  136.     /**
  137.      * Get the covariance orbit type.
  138.      * @return the covariance orbit type
  139.      */
  140.     public OrbitType getOrbitType() {
  141.         return orbitType;
  142.     }

  143.     /**
  144.      * Get the covariance angle type.
  145.      * @return the covariance angle type
  146.      */
  147.     public PositionAngle getPositionAngle() {
  148.         return angleType;
  149.     }

  150.     /**
  151.      * Get the covariance frame.
  152.      * @return the covariance frame (can be null)
  153.      * @see #getLOFType()
  154.      */
  155.     public Frame getFrame() {
  156.         return frame;
  157.     }

  158.     /**
  159.      * Get the covariance LOF type.
  160.      * @return the covariance LOF type (can be null)
  161.      * @see #getFrame()
  162.      */
  163.     public LOFType getLOFType() {
  164.         return lofType;
  165.     }

  166.     /**
  167.      * Get the covariance matrix in another orbit type.
  168.      * <p>
  169.      * The covariance orbit type <b>cannot</b> be changed if the covariance
  170.      * matrix is expressed in a {@link LOFType local orbital frame} or a
  171.      * non-pseudo inertial frame.
  172.      * <p>
  173.      * As this type change uses the jacobian matrix of the transformation, it introduces a linear approximation.
  174.      * Hence, the current covariance matrix <b>will not exactly match</b> the new linearized case and the
  175.      * distribution will not follow a generalized Gaussian distribution anymore.
  176.      * <p>
  177.      * This is based on equation (1) to (6) from "Vallado, D. A. (2004). Covariance transformations for satellite flight
  178.      * dynamics operations."
  179.      *
  180.      * @param orbit orbit to which the covariance matrix should correspond
  181.      * @param outOrbitType target orbit type of the state covariance matrix
  182.      * @param outAngleType target position angle type of the state covariance matrix
  183.      * @return a new covariance state, expressed in the target orbit type with the target position angle
  184.      * @see #changeCovarianceFrame(Orbit, Frame)
  185.      */
  186.     public StateCovariance changeCovarianceType(final Orbit orbit, final OrbitType outOrbitType,
  187.                                                 final PositionAngle outAngleType) {

  188.         // Check if the covariance expressed in a celestial body frame
  189.         if (frame != null) {

  190.             // Check if the covarianc is defined in inertial frame
  191.             if (frame.isPseudoInertial()) {
  192.                 return changeTypeAndCreate(orbit, epoch, frame, orbitType, angleType, outOrbitType, outAngleType,
  193.                                            orbitalCovariance);
  194.             }

  195.             // The covariance is not defined in inertial frame. The orbit type cannot be changes
  196.             throw new OrekitException(OrekitMessages.CANNOT_CHANGE_COVARIANCE_TYPE_IF_DEFINED_IN_NON_INERTIAL_FRAME);

  197.         }

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

  200.     }

  201.     /**
  202.      * Get the covariance in a given local orbital frame.
  203.      * <p>
  204.      * Changing the covariance frame is a linear process, this method does not introduce approximation unless a change
  205.      * in covariance orbit type is required.
  206.      * <p>
  207.      * This is based on equation (18) to (20) "from Vallado, D. A. (2004). Covariance transformations for satellite
  208.      * flight dynamics operations."
  209.      *
  210.      * @param orbit orbit to which the covariance matrix should correspond
  211.      * @param lofOut output local orbital frame
  212.      * @return a new covariance state, expressed in the output local orbital frame
  213.      */
  214.     public StateCovariance changeCovarianceFrame(final Orbit orbit, final LOFType lofOut) {

  215.         // Verify current covariance frame
  216.         if (lofType != null) {

  217.             // Change the covariance local orbital frame
  218.             return changeFrameAndCreate(orbit, epoch, lofType, lofOut, orbitalCovariance);

  219.         } else {

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

  222.         }

  223.     }

  224.     /**
  225.      * Get the covariance in the output frame.
  226.      * <p>
  227.      * Changing the covariance frame is a linear process, this method does not introduce approximation unless a change
  228.      * in covariance orbit type is required.
  229.      * <p>
  230.      * This is based on equation (18) to (20) "from Vallado, D. A. (2004). Covariance transformations for satellite
  231.      * flight dynamics operations."
  232.      *
  233.      * @param orbit orbit to which the covariance matrix should correspond
  234.      * @param frameOut output frame
  235.      * @return a new covariance state, expressed in the output frame
  236.      */
  237.     public StateCovariance changeCovarianceFrame(final Orbit orbit, final Frame frameOut) {

  238.         // Verify current covariance frame
  239.         if (lofType != null) {

  240.             // Covariance is expressed in local orbital frame
  241.             return changeFrameAndCreate(orbit, epoch, lofType, frameOut, orbitalCovariance);

  242.         } else {

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

  245.         }

  246.     }

  247.     /**
  248.      * Get a time-shifted covariance matrix.
  249.      * <p>
  250.      * The shifting model is a linearized, Keplerian one. In other words, it is based on a state transition matrix that
  251.      * is computed assuming Keplerian motion.
  252.      * <p>
  253.      * Shifting is <em>not</em> intended as a replacement for proper covariance propagation, but should be sufficient
  254.      * for small time shifts or coarse accuracy.
  255.      *
  256.      * @param orbit orbit to which the covariance matrix should correspond
  257.      * @param dt time shift in seconds
  258.      * @return a new covariance state, shifted with respect to the instance
  259.      */
  260.     public StateCovariance shiftedBy(final Orbit orbit, final double dt) {

  261.         // Shifted orbit
  262.         final Orbit shifted = orbit.shiftedBy(dt);

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

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

  267.                 // Compute STM
  268.                 final RealMatrix stm = getStm(orbit, dt);

  269.                 // Convert covariance in STM type (i.e., Equinoctial elements)
  270.                 final StateCovariance inStmType = changeTypeAndCreate(orbit, epoch, frame, orbitType, angleType,
  271.                                                                       OrbitType.EQUINOCTIAL, PositionAngle.MEAN,
  272.                                                                       orbitalCovariance);

  273.                 // Shift covariance by applying the STM
  274.                 final RealMatrix shiftedCov = stm.multiply(inStmType.getMatrix().multiplyTransposed(stm));

  275.                 // Restore the initial covariance type
  276.                 return changeTypeAndCreate(shifted, shifted.getDate(), frame,
  277.                                            OrbitType.EQUINOCTIAL, PositionAngle.MEAN,
  278.                                            orbitType, angleType, shiftedCov);
  279.             }

  280.             // State covariance expressed in a non pseudo-inertial frame
  281.             else {

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

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

  286.                 // Restore the initial covariance frame
  287.                 return shiftedCovariance.changeCovarianceFrame(shifted, frame);
  288.             }
  289.         }

  290.         // State covariance expressed in a commonly used local orbital frame (LOFType)
  291.         else {

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

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

  296.             // Restore the initial covariance frame
  297.             return shiftedCovariance.changeCovarianceFrame(shifted, lofType);
  298.         }

  299.     }

  300.     /**
  301.      * Create a covariance matrix in another orbit type.
  302.      * <p>
  303.      * As this type change uses the jacobian matrix of the transformation, it introduces a linear approximation.
  304.      * Hence, the input covariance matrix <b>will not exactly match</b> the new linearized case and the
  305.      * distribution will not follow a generalized Gaussian distribution anymore.
  306.      * <p>
  307.      * This is based on equation (1) to (6) from "Vallado, D. A. (2004). Covariance transformations for satellite flight
  308.      * dynamics operations."
  309.      *
  310.      * @param orbit orbit to which the covariance matrix should correspond
  311.      * @param date covariance epoch
  312.      * @param covFrame covariance frame
  313.      * @param inOrbitType initial orbit type of the state covariance matrix
  314.      * @param inAngleType initial position angle type of the state covariance matrix
  315.      * @param outOrbitType target orbit type of the state covariance matrix
  316.      * @param outAngleType target position angle type of the state covariance matrix
  317.      * @param inputCov input covariance
  318.      * @return the covariance expressed in the target orbit type with the target position angle
  319.      */
  320.     private static StateCovariance changeTypeAndCreate(final Orbit orbit, final AbsoluteDate date,
  321.                                                        final Frame covFrame,
  322.                                                        final OrbitType inOrbitType, final PositionAngle inAngleType,
  323.                                                        final OrbitType outOrbitType, final PositionAngle outAngleType,
  324.                                                        final RealMatrix inputCov) {

  325.         // Notations:
  326.         // I: Input orbit type
  327.         // O: Output orbit type
  328.         // C: Cartesian parameters

  329.         // Compute dOutputdCartesian
  330.         final double[][] aOC               = new double[STATE_DIMENSION][STATE_DIMENSION];
  331.         final Orbit      orbitInOutputType = outOrbitType.convertType(orbit);
  332.         orbitInOutputType.getJacobianWrtCartesian(outAngleType, aOC);
  333.         final RealMatrix dOdC = new Array2DRowRealMatrix(aOC, false);

  334.         // Compute dCartesiandInput
  335.         final double[][] aCI              = new double[STATE_DIMENSION][STATE_DIMENSION];
  336.         final Orbit      orbitInInputType = inOrbitType.convertType(orbit);
  337.         orbitInInputType.getJacobianWrtParameters(inAngleType, aCI);
  338.         final RealMatrix dCdI = new Array2DRowRealMatrix(aCI, false);

  339.         // Compute dOutputdInput
  340.         final RealMatrix dOdI = dOdC.multiply(dCdI);

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

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

  345.     }

  346.     /**
  347.      * Create a covariance matrix from a {@link LOFType local orbital frame} to another
  348.      * {@link LOFType local orbital frame}.
  349.      * <p>
  350.      * Changing the covariance frame is a linear process, this method does not introduce approximation.
  351.      * <p>
  352.      * The transformation is based on Equation (18) to (20) of "Covariance Transformations for Satellite Flight Dynamics
  353.      * Operations" by David A. Vallado.
  354.      *
  355.      * @param orbit orbit to which the covariance matrix should correspond
  356.      * @param date covariance epoch
  357.      * @param lofIn the local orbital frame in which the input covariance matrix is expressed
  358.      * @param lofOut the target local orbital frame
  359.      * @param inputCartesianCov input covariance {@code CARTESIAN})
  360.      * @return the covariance matrix expressed in the target commonly used local orbital frame in Cartesian elements
  361.      */
  362.     private static StateCovariance changeFrameAndCreate(final Orbit orbit, final AbsoluteDate date,
  363.                                                         final LOFType lofIn, final LOFType lofOut,
  364.                                                         final RealMatrix inputCartesianCov) {

  365.         // Builds the matrix to perform covariance transformation
  366.         final RealMatrix jacobianFromLofInToLofOut =
  367.                         getJacobian(LOFType.transformFromLOFInToLOFOut(lofIn, lofOut, date, orbit.getPVCoordinates()));

  368.         // Get the Cartesian covariance matrix converted to frameOut
  369.         final RealMatrix cartesianCovarianceOut =
  370.                 jacobianFromLofInToLofOut.multiply(inputCartesianCov.multiplyTransposed(jacobianFromLofInToLofOut));

  371.         // Output converted covariance
  372.         return new StateCovariance(cartesianCovarianceOut, date, lofOut);

  373.     }

  374.     /**
  375.      * Convert the covariance matrix from a {@link Frame frame} to a {@link LOFType local orbital frame}.
  376.      * <p>
  377.      * Changing the covariance frame is a linear process, this method does not introduce approximation unless a change
  378.      * in covariance orbit type is required.
  379.      * <p>
  380.      * The transformation is based on Equation (18) to (20) of "Covariance Transformations for Satellite Flight Dynamics
  381.      * Operations" by David A. Vallado.
  382.      * <p>
  383.      * As the frame transformation must be performed with the covariance expressed in Cartesian elements, both the orbit
  384.      * and position angle types of the input covariance must be provided.
  385.      * <p>
  386.      * <b>The output covariance matrix will provided in Cartesian orbit type and not converted back to
  387.      * its original expression (if input different from Cartesian elements).</b>
  388.      *
  389.      * @param orbit orbit to which the covariance matrix should correspond
  390.      * @param date covariance epoch
  391.      * @param frameIn the frame in which the input covariance matrix is expressed. In case the frame is <b>not</b>
  392.      * pseudo-inertial, the input covariance matrix is expected to be expressed in <b>Cartesian elements</b>.
  393.      * @param lofOut the target local orbital frame
  394.      * @param inputCov input covariance
  395.      * @param covOrbitType orbit type of the covariance matrix (used if frameIn is pseudo-inertial)
  396.      * @param covAngleType position angle type of the covariance matrix (used if frameIn is pseudo-inertial) (not used
  397.      * if covOrbitType equals {@code CARTESIAN})
  398.      * @return the covariance matrix expressed in the target local orbital frame in Cartesian elements
  399.      * @throws OrekitException if input frame is <b>not</b> pseudo-inertial <b>and</b> the input covariance is
  400.      * <b>not</b> expressed in Cartesian elements.
  401.      */
  402.     private static StateCovariance changeFrameAndCreate(final Orbit orbit, final AbsoluteDate date,
  403.                                                         final Frame frameIn, final LOFType lofOut,
  404.                                                         final RealMatrix inputCov,
  405.                                                         final OrbitType covOrbitType,
  406.                                                         final PositionAngle covAngleType) {

  407.         // Input frame is inertial
  408.         if (frameIn.isPseudoInertial()) {

  409.             // Convert input matrix to Cartesian parameters in input frame
  410.             final RealMatrix cartesianCovarianceIn =
  411.                     changeTypeAndCreate(orbit, date, frameIn, covOrbitType, covAngleType,
  412.                                         OrbitType.CARTESIAN, PositionAngle.MEAN,
  413.                                         inputCov).getMatrix();

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

  417.             // Get the Cartesian covariance matrix converted to frameOut
  418.             final RealMatrix cartesianCovarianceOut =
  419.                     jacobianFromFrameInToLofOut.multiply(cartesianCovarianceIn.multiplyTransposed(jacobianFromFrameInToLofOut));

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

  422.         }

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

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

  427.             // Compute rotation matrix from frameIn to orbit inertial frame
  428.             final RealMatrix cartesianCovarianceInOrbitFrame =
  429.                    changeFrameAndCreate(orbit, date, frameIn, orbitInertialFrame, inputCov,
  430.                                          OrbitType.CARTESIAN, PositionAngle.MEAN).getMatrix();

  431.             // Convert from orbit inertial frame to lofOut
  432.             return changeFrameAndCreate(orbit, date, orbitInertialFrame, lofOut, cartesianCovarianceInOrbitFrame,
  433.                                         OrbitType.CARTESIAN, PositionAngle.MEAN);

  434.         }

  435.     }

  436.     /**
  437.      * Convert the covariance matrix from a {@link LOFType  local orbital frame} to a {@link Frame frame}.
  438.      * <p>
  439.      * Changing the covariance frame is a linear process, this method does not introduce approximation.
  440.      * <p>
  441.      * The transformation is based on Equation (18) to (20) of "Covariance Transformations for Satellite Flight Dynamics
  442.      * Operations" by David A. Vallado.
  443.      * <p>
  444.      * The <u>input</u> covariance matrix is necessarily provided in <b>Cartesian orbit type</b>.
  445.      * <p>
  446.      * The <u>output</u> covariance matrix will be expressed in <b>Cartesian elements</b>.
  447.      *
  448.      * @param orbit orbit to which the covariance matrix should correspond
  449.      * @param date covariance epoch
  450.      * @param lofIn the local orbital frame in which the input covariance matrix is expressed
  451.      * @param frameOut the target frame
  452.      * @param inputCartesianCov input covariance ({@code CARTESIAN})
  453.      * @return the covariance matrix expressed in the target frame in Cartesian elements
  454.      */
  455.     private static StateCovariance changeFrameAndCreate(final Orbit orbit, final AbsoluteDate date,
  456.                                                         final LOFType lofIn, final Frame frameOut,
  457.                                                         final RealMatrix inputCartesianCov) {

  458.         // Output frame is pseudo-inertial
  459.         if (frameOut.isPseudoInertial()) {

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

  463.             // Transform covariance
  464.             final RealMatrix transformedCovariance =
  465.                     jacobianFromLofInToFrameOut.multiply(inputCartesianCov.multiplyTransposed(jacobianFromLofInToFrameOut));

  466.             // Get the Cartesian covariance matrix converted to frameOut
  467.             return new StateCovariance(transformedCovariance, date, frameOut, OrbitType.CARTESIAN,
  468.                                        DEFAULT_POSITION_ANGLE);

  469.         }

  470.         // Output frame is not pseudo-inertial
  471.         else {

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

  475.             // Get the Cartesian covariance matrix converted to orbit inertial frame
  476.             final RealMatrix cartesianCovarianceInOrbitFrame = jacobianFromLofInToOrbitFrame.multiply(
  477.                     inputCartesianCov.multiplyTransposed(jacobianFromLofInToOrbitFrame));

  478.             // Get the Cartesian covariance matrix converted to frameOut
  479.             return changeFrameAndCreate(orbit, date, orbit.getFrame(), frameOut, cartesianCovarianceInOrbitFrame,
  480.                                         OrbitType.CARTESIAN, PositionAngle.MEAN);
  481.         }

  482.     }

  483.     /**
  484.      * Get the covariance matrix in another frame.
  485.      * <p>
  486.      * The transformation is based on Equation (18) of "Covariance Transformations for Satellite Flight Dynamics
  487.      * Operations" by David A. Vallado.
  488.      * <p>
  489.      * Changing the covariance frame is a linear process, this method does not introduce approximation unless a change
  490.      * in covariance orbit type is required.
  491.      * <p>
  492.      * As the frame transformation must be performed with the covariance expressed in Cartesian elements, both the orbit
  493.      * and position angle types of the input covariance must be provided.
  494.      * <p>
  495.      * In case the <u>input</u> frame is <b>not</b> pseudo-inertial, the <u>input</u> covariance matrix is necessarily
  496.      * expressed in <b>Cartesian elements</b>.
  497.      * <p>
  498.      * In case the <u>output</u> frame is <b>not</b> pseudo-inertial, the <u>output</u> covariance matrix will be
  499.      * expressed in <b>Cartesian elements</b>.
  500.      *
  501.      * @param orbit orbit to which the covariance matrix should correspond
  502.      * @param date covariance epoch
  503.      * @param frameIn the frame in which the input covariance matrix is expressed
  504.      * @param frameOut the target frame
  505.      * @param inputCov input covariance
  506.      * @param covOrbitType orbit type of the covariance matrix (used if frameIn is pseudo-inertial)
  507.      * @param covAngleType position angle type of the covariance matrix (used if frameIn is pseudo-inertial) (<b>not</b>
  508.      * used if covOrbitType equals {@code CARTESIAN})
  509.      * @return the covariance matrix expressed in the target frame
  510.      * @throws OrekitException if input frame is <b>not</b> pseudo-inertial <b>and</b> the input covariance is
  511.      * <b>not</b> expressed in Cartesian elements.
  512.      */
  513.     private static StateCovariance changeFrameAndCreate(final Orbit orbit, final AbsoluteDate date,
  514.                                                         final Frame frameIn, final Frame frameOut,
  515.                                                         final RealMatrix inputCov,
  516.                                                         final OrbitType covOrbitType,
  517.                                                         final PositionAngle covAngleType) {

  518.         // Get the transform from the covariance frame to the output frame
  519.         final Transform inToOut = frameIn.getTransformTo(frameOut, orbit.getDate());

  520.         // Matrix to perform the covariance transformation
  521.         final RealMatrix j = getJacobian(inToOut);

  522.         // Input frame pseudo-inertial
  523.         if (frameIn.isPseudoInertial()) {

  524.             // Convert input matrix to Cartesian parameters in input frame
  525.             final RealMatrix cartesianCovarianceIn =
  526.                     changeTypeAndCreate(orbit, date, frameIn, covOrbitType, covAngleType,
  527.                                         OrbitType.CARTESIAN, PositionAngle.MEAN,
  528.                                         inputCov).getMatrix();

  529.             // Get the Cartesian covariance matrix converted to frameOut
  530.             final RealMatrix cartesianCovarianceOut = j.multiply(cartesianCovarianceIn.multiplyTransposed(j));

  531.             // Output frame is pseudo-inertial
  532.             if (frameOut.isPseudoInertial()) {

  533.                 // Convert output Cartesian matrix to initial orbit type and position angle
  534.                 return changeTypeAndCreate(orbit, date, frameOut, OrbitType.CARTESIAN, PositionAngle.MEAN,
  535.                                            covOrbitType, covAngleType, cartesianCovarianceOut);

  536.             }

  537.             // Output frame is not pseudo-inertial
  538.             else {

  539.                 // Output Cartesian matrix
  540.                 return new StateCovariance(cartesianCovarianceOut, date, frameOut, OrbitType.CARTESIAN,
  541.                                            DEFAULT_POSITION_ANGLE);

  542.             }
  543.         }

  544.         // Input frame is not pseudo-inertial
  545.         else {

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

  547.             // Convert covariance matrix to frameOut
  548.             final RealMatrix covInFrameOut = j.multiply(inputCov.multiplyTransposed(j));

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

  551.         }

  552.     }

  553.     /**
  554.      * Builds the matrix to perform covariance frame transformation.
  555.      * @param transform input transformation
  556.      * @return the matrix to perform the covariance frame transformation
  557.      */
  558.     private static RealMatrix getJacobian(final Transform transform) {
  559.         // Get the Jacobian of the transform
  560.         final double[][] jacobian = new double[STATE_DIMENSION][STATE_DIMENSION];
  561.         transform.getJacobian(CartesianDerivativesFilter.USE_PV, jacobian);
  562.         // Return
  563.         return new Array2DRowRealMatrix(jacobian, false);

  564.     }

  565.     /**
  566.      * Get the state transition matrix considering Keplerian contribution only.
  567.      *
  568.      * @param initialOrbit orbit to which the initial covariance matrix should correspond
  569.      * @param dt time difference between the two orbits
  570.      * @return the state transition matrix used to shift the covariance matrix
  571.      */
  572.     private RealMatrix getStm(final Orbit initialOrbit, final double dt) {

  573.         // initialize the STM
  574.         final RealMatrix stm = MatrixUtils.createRealIdentityMatrix(STATE_DIMENSION);

  575.         // State transition matrix using Keplerian contribution only
  576.         final double mu           = initialOrbit.getMu();
  577.         final double sma          = initialOrbit.getA();
  578.         final double contribution = -1.5 * dt * FastMath.sqrt(mu / FastMath.pow(sma, 5));
  579.         stm.setEntry(5, 0, contribution);

  580.         // Return
  581.         return stm;

  582.     }

  583. }