1   /* Copyright 2002-2025 CS GROUP
2    * Licensed to CS GROUP (CS) under one or more
3    * contributor license agreements.  See the NOTICE file distributed with
4    * this work for additional information regarding copyright ownership.
5    * CS licenses this file to You under the Apache License, Version 2.0
6    * (the "License"); you may not use this file except in compliance with
7    * the License.  You may obtain a copy of the License at
8    *
9    *   http://www.apache.org/licenses/LICENSE-2.0
10   *
11   * Unless required by applicable law or agreed to in writing, software
12   * distributed under the License is distributed on an "AS IS" BASIS,
13   * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14   * See the License for the specific language governing permissions and
15   * limitations under the License.
16   */
17  package org.orekit.orbits;
18  
19  import org.hipparchus.analysis.polynomials.SmoothStepFactory;
20  import org.orekit.errors.OrekitException;
21  import org.orekit.frames.Frame;
22  import org.orekit.propagation.Propagator;
23  import org.orekit.propagation.SpacecraftState;
24  import org.orekit.propagation.analytical.AbstractAnalyticalPropagator;
25  import org.orekit.time.AbsoluteDate;
26  import org.orekit.utils.PVCoordinates;
27  
28  import java.util.List;
29  
30  /**
31   * Orbit blender.
32   * <p>
33   * Its purpose is to interpolate orbit state between tabulated orbit states using the concept of blending, exposed in :
34   * "Efficient Covariance Interpolation using Blending of Approximate State Error Transitions" by Sergei Tanygin, and applying
35   * it to orbit states instead of covariances.
36   * <p>
37   * It propagates tabulated values to the interpolating time using given propagator and then blend each propagated
38   * states using a smoothstep function. It gives especially good results as explained
39   * <a href="https://orekit.org/doc/technical-notes/Implementation_of_covariance_interpolation_in_Orekit.pdf">here</a>
40   * compared to Hermite interpolation when time steps between tabulated values get significant (In LEO, &gt; 10 mn for
41   * example).
42   * <p>
43   * <b>In most cases, an analytical propagator would be used to quickly fill the gap between tabulated values and recreate
44   * a dense ephemeris</b>.
45   * <p>
46   * However, a fully configured and accurate numerical propagator can be used to recreate an even
47   * more precise ephemeris in case the initial tabulated values were obtained from an external source.
48   * <p>
49   * Note that in the current implementation, the returned blended orbit is necessarily Cartesian.
50   *
51   * @author Vincent Cucchietti
52   * @see org.hipparchus.analysis.polynomials.SmoothStepFactory
53   * @see org.hipparchus.analysis.polynomials.SmoothStepFactory.SmoothStepFunction
54   * @see Propagator
55   * @see AbstractAnalyticalPropagator
56   *
57   * @since 12.0
58   */
59  public class OrbitBlender extends AbstractOrbitInterpolator {
60  
61      /** Propagator used to propagate tabulated orbits to interpolating time. */
62      private final Propagator blendingPropagator;
63  
64      /** Blending function. */
65      private final SmoothStepFactory.SmoothStepFunction blendingFunction;
66  
67      /**
68       * Default constructor.
69       * <p>
70       * <b>In most cases, an analytical propagator would be used to quickly fill the gap between tabulated values and recreate
71       * a dense ephemeris</b>.
72       * <p>
73       * However, a fully configured and accurate numerical propagator can be used to recreate an even
74       * more precise ephemeris in case the initial tabulated values were obtained from an external source.
75       *
76       * @param blendingFunction
77       * {@link org.hipparchus.analysis.polynomials.SmoothStepFactory.SmoothStepFunction smoothstep function} used for
78       * blending
79       * @param blendingPropagator propagator used to propagate tabulated orbits to interpolating time
80       * @param outputInertialFrame output inertial frame
81       *
82       * @throws OrekitException if output frame is not inertial
83       * @see org.hipparchus.analysis.polynomials.SmoothStepFactory.SmoothStepFunction
84       */
85      public OrbitBlender(final SmoothStepFactory.SmoothStepFunction blendingFunction,
86                          final Propagator blendingPropagator,
87                          final Frame outputInertialFrame) {
88          super(DEFAULT_INTERPOLATION_POINTS, 0., outputInertialFrame);
89          this.blendingFunction   = blendingFunction;
90          this.blendingPropagator = blendingPropagator;
91      }
92  
93      /** {@inheritDoc} */
94      @Override
95      protected Orbit interpolate(final InterpolationData interpolationData) {
96  
97          // Get first and last entry
98          final List<Orbit> neighborList  = interpolationData.getNeighborList();
99          final Orbit       previousOrbit = neighborList.get(0);
100         final Orbit       nextOrbit     = neighborList.get(1);
101 
102         // Propagate orbits
103         final AbsoluteDate interpolationDate = interpolationData.getInterpolationDate();
104         final Orbit forwardedOrbit  = propagateOrbit(previousOrbit, interpolationDate);
105         final Orbit backwardedOrbit = propagateOrbit(nextOrbit, interpolationDate);
106 
107         // Extract position-velocity-acceleration coordinates
108         final PVCoordinates forwardedPV  = forwardedOrbit.getPVCoordinates(getOutputInertialFrame());
109         final PVCoordinates backwardedPV = backwardedOrbit.getPVCoordinates(getOutputInertialFrame());
110 
111         // Blend PV coordinates
112         final double timeParameter = getTimeParameter(interpolationDate, previousOrbit.getDate(), nextOrbit.getDate());
113         final double blendingValue = blendingFunction.value(timeParameter);
114 
115         final PVCoordinates blendedPV = forwardedPV.blendArithmeticallyWith(backwardedPV, blendingValue);
116 
117         // Output new blended instance
118         return new CartesianOrbit(blendedPV, getOutputInertialFrame(), interpolationDate, previousOrbit.getMu());
119     }
120 
121     /**
122      * Propagate orbit using predefined {@link Propagator propagator}.
123      *
124      * @param tabulatedOrbit tabulated orbit to propagate
125      * @param propagationDate propagation date
126      *
127      * @return orbit propagated to propagation date
128      */
129     private Orbit propagateOrbit(final Orbit tabulatedOrbit,
130                                  final AbsoluteDate propagationDate) {
131         blendingPropagator.resetInitialState(new SpacecraftState(tabulatedOrbit));
132         return blendingPropagator.propagate(propagationDate).getOrbit();
133     }
134 }