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  
19  import org.hipparchus.analysis.polynomials.SmoothStepFactory;
20  import org.hipparchus.linear.RealMatrix;
21  import org.orekit.frames.Frame;
22  import org.orekit.frames.LOFType;
23  import org.orekit.orbits.Orbit;
24  import org.orekit.orbits.OrbitType;
25  import org.orekit.orbits.PositionAngleType;
26  import org.orekit.time.AbsoluteDate;
27  import org.orekit.time.TimeInterpolator;
28  import org.orekit.time.TimeStampedPair;
29  
30  import java.util.List;
31  
32  /**
33   * State covariance blender.
34   * <p>
35   * Its purpose is to interpolate state covariance between tabulated state covariances by using the concept of blending,
36   * exposed in : "Efficient Covariance Interpolation using Blending of Approximate State Error Transitions" by Sergei
37   * Tanygin.
38   * <p>
39   * It propagates tabulated values to the interpolation date assuming a standard keplerian model and then blend each
40   * propagated covariances using a smoothstep function.
41   * <p>
42   * It gives accurate results as explained <a
43   * href="https://orekit.org/doc/technical-notes/Implementation_of_covariance_interpolation_in_Orekit.pdf">here</a>. In the
44   * very poorly tracked test case evolving in a highly dynamical environment mentioned in the linked thread, the user can
45   * expect at worst errors of less than 0.25% in position sigmas and less than 0.4% in velocity sigmas with steps of 40mn
46   * between tabulated values.
47   *
48   * @author Vincent Cucchietti
49   * @see org.hipparchus.analysis.polynomials.SmoothStepFactory
50   * @see org.hipparchus.analysis.polynomials.SmoothStepFactory.SmoothStepFunction
51   */
52  public class StateCovarianceBlender extends AbstractStateCovarianceInterpolator {
53  
54      /** Blending function. */
55      private final SmoothStepFactory.SmoothStepFunction blendingFunction;
56  
57      /**
58       * Constructor.
59       * <p>
60       * <b>BEWARE:</b> If the output local orbital frame is not considered pseudo-inertial, all the covariance components
61       * related to the velocity will be poorly interpolated. <b>Only the position covariance should be considered in this
62       * case.</b>
63       *
64       * @param blendingFunction blending function
65       * @param orbitInterpolator orbit interpolator
66       * @param outLOF local orbital frame
67       *
68       * @see Frame
69       * @see OrbitType
70       * @see PositionAngleType
71       */
72      public StateCovarianceBlender(final SmoothStepFactory.SmoothStepFunction blendingFunction,
73                                    final TimeInterpolator<Orbit> orbitInterpolator,
74                                    final LOFType outLOF) {
75          super(DEFAULT_INTERPOLATION_POINTS, 0., orbitInterpolator, outLOF);
76          this.blendingFunction = blendingFunction;
77      }
78  
79      /**
80       * Constructor.
81       *
82       * @param blendingFunction blending function
83       * @param orbitInterpolator orbit interpolator
84       * @param outFrame desired output covariance frame
85       * @param outPositionAngleType desired output position angle
86       * @param outOrbitType desired output orbit type
87       *
88       * @see Frame
89       * @see OrbitType
90       * @see PositionAngleType
91       */
92      public StateCovarianceBlender(final SmoothStepFactory.SmoothStepFunction blendingFunction,
93                                    final TimeInterpolator<Orbit> orbitInterpolator,
94                                    final Frame outFrame,
95                                    final OrbitType outOrbitType,
96                                    final PositionAngleType outPositionAngleType) {
97          super(DEFAULT_INTERPOLATION_POINTS, 0., orbitInterpolator, outFrame, outOrbitType, outPositionAngleType);
98          this.blendingFunction = blendingFunction;
99      }
100 
101     /** {@inheritDoc} */
102     @Override
103     protected StateCovariance computeInterpolatedCovarianceInOrbitFrame(
104             final List<TimeStampedPair<Orbit, StateCovariance>> uncertainStates,
105             final Orbit interpolatedOrbit) {
106 
107         // Necessarily only two sample for blending
108         final TimeStampedPair<Orbit, StateCovariance> previousUncertainState = uncertainStates.get(0);
109         final TimeStampedPair<Orbit, StateCovariance> nextUncertainState     = uncertainStates.get(1);
110 
111         // Get the interpolation date
112         final AbsoluteDate interpolationDate = interpolatedOrbit.getDate();
113 
114         // Propagate previous tabulated covariance to interpolation date
115         final StateCovariance forwardedCovariance =
116                 propagateCovarianceAnalytically(interpolationDate, interpolatedOrbit, previousUncertainState);
117 
118         // Propagate next tabulated covariance to interpolation date
119         final StateCovariance backwardedCovariance =
120                 propagateCovarianceAnalytically(interpolationDate, interpolatedOrbit, nextUncertainState);
121 
122         // Compute normalized time parameter
123         final double timeParameter =
124                 getTimeParameter(interpolationDate, previousUncertainState.getDate(), nextUncertainState.getDate());
125 
126         // Blend the covariance matrices
127         final double     blendingValue              = blendingFunction.value(timeParameter);
128         final RealMatrix forwardedCovarianceMatrix  = forwardedCovariance.getMatrix();
129         final RealMatrix backwardedCovarianceMatrix = backwardedCovariance.getMatrix();
130 
131         final RealMatrix blendedCovarianceMatrix =
132                 forwardedCovarianceMatrix.blendArithmeticallyWith(backwardedCovarianceMatrix, blendingValue);
133 
134         return new StateCovariance(blendedCovarianceMatrix, interpolationDate, interpolatedOrbit.getFrame(),
135                                    OrbitType.CARTESIAN, DEFAULT_POSITION_ANGLE);
136     }
137 
138     /**
139      * Propagate given state covariance to the interpolation date using keplerian motion.
140      *
141      * @param interpolationTime interpolation date at which we desire to interpolate the state covariance
142      * @param orbitAtInterpolatingTime orbit at interpolation date
143      * @param tabulatedUncertainState tabulated uncertain state
144      *
145      * @return state covariance at given interpolation date.
146      *
147      * @see StateCovariance
148      */
149     private StateCovariance propagateCovarianceAnalytically(final AbsoluteDate interpolationTime,
150                                                             final Orbit orbitAtInterpolatingTime,
151                                                             final TimeStampedPair<Orbit, StateCovariance> tabulatedUncertainState) {
152 
153         // Get orbit and covariance
154         final Orbit           tabulatedOrbit         = tabulatedUncertainState.getFirst();
155         final StateCovariance tabulatedCovariance    = tabulatedUncertainState.getSecond();
156         final Frame           interpolatedOrbitFrame = orbitAtInterpolatingTime.getFrame();
157 
158         // Express tabulated covariance in interpolated orbit frame for consistency
159         final StateCovariance tabulatedCovarianceInOrbitFrame =
160                 tabulatedCovariance.changeCovarianceFrame(tabulatedOrbit, interpolatedOrbitFrame);
161 
162         // First convert the covariance matrix to equinoctial elements to avoid singularities inherent to keplerian elements
163         final RealMatrix covarianceMatrixInEquinoctial =
164                 tabulatedCovarianceInOrbitFrame.changeCovarianceType(tabulatedOrbit, OrbitType.EQUINOCTIAL,
165                                                                      DEFAULT_POSITION_ANGLE).getMatrix();
166 
167         // Compute state error transition matrix in equinoctial elements (identical to the one in keplerian elements)
168         final RealMatrix stateErrorTransitionMatrixInEquinoctial =
169                 StateCovariance.getStm(tabulatedOrbit, interpolationTime.durationFrom(tabulatedOrbit.getDate()));
170 
171         // Propagate the covariance matrix using the previously computed state error transition matrix
172         final RealMatrix propagatedCovarianceMatrixInEquinoctial =
173                 stateErrorTransitionMatrixInEquinoctial.multiply(
174                         covarianceMatrixInEquinoctial.multiplyTransposed(stateErrorTransitionMatrixInEquinoctial));
175 
176         // Recreate a StateCovariance instance
177         final StateCovariance propagatedCovarianceInEquinoctial =
178                 new StateCovariance(propagatedCovarianceMatrixInEquinoctial, interpolationTime,
179                                     interpolatedOrbitFrame, OrbitType.EQUINOCTIAL, DEFAULT_POSITION_ANGLE);
180 
181         // Output propagated state covariance after converting back to cartesian elements
182         return propagatedCovarianceInEquinoctial.changeCovarianceType(orbitAtInterpolatingTime,
183                                                                       OrbitType.CARTESIAN, DEFAULT_POSITION_ANGLE);
184     }
185 }