1 /* Copyright 2002-2026 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.orekit.errors.OrekitException;
21 import org.orekit.frames.Frame;
22 import org.orekit.orbits.AbstractOrbitInterpolator;
23 import org.orekit.orbits.CartesianOrbit;
24 import org.orekit.orbits.Orbit;
25 import org.orekit.propagation.analytical.AbstractAnalyticalPropagator;
26 import org.orekit.time.AbsoluteDate;
27 import org.orekit.utils.PVCoordinates;
28
29 import java.util.List;
30
31 /**
32 * Orbit blender.
33 * <p>
34 * Its purpose is to interpolate orbit state between tabulated orbit states using the concept of blending, exposed in :
35 * "Efficient Covariance Interpolation using Blending of Approximate State Error Transitions" by Sergei Tanygin, and applying
36 * it to orbit states instead of covariances.
37 * <p>
38 * It propagates tabulated values to the interpolating time using given propagator and then blend each propagated
39 * states using a smoothstep function. It gives especially good results as explained
40 * <a href="https://orekit.org/doc/technical-notes/Implementation_of_covariance_interpolation_in_Orekit.pdf">here</a>
41 * compared to Hermite interpolation when time steps between tabulated values get significant (In LEO, > 10 mn for
42 * example).
43 * <p>
44 * <b>In most cases, an analytical propagator would be used to quickly fill the gap between tabulated values and recreate
45 * a dense ephemeris</b>.
46 * <p>
47 * However, a fully configured and accurate numerical propagator can be used to recreate an even
48 * more precise ephemeris in case the initial tabulated values were obtained from an external source.
49 * <p>
50 * Note that in the current implementation, the returned blended orbit is necessarily Cartesian.
51 *
52 * @author Vincent Cucchietti
53 * @see org.hipparchus.analysis.polynomials.SmoothStepFactory
54 * @see org.hipparchus.analysis.polynomials.SmoothStepFactory.SmoothStepFunction
55 * @see Propagator
56 * @see AbstractAnalyticalPropagator
57 *
58 * @since 12.0
59 */
60 public class OrbitBlender extends AbstractOrbitInterpolator {
61
62 /** Propagator used to propagate tabulated orbits to interpolating time. */
63 private final Propagator blendingPropagator;
64
65 /** Blending function. */
66 private final SmoothStepFactory.SmoothStepFunction blendingFunction;
67
68 /**
69 * Default constructor.
70 * <p>
71 * <b>In most cases, an analytical propagator would be used to quickly fill the gap between tabulated values and recreate
72 * a dense ephemeris</b>.
73 * <p>
74 * However, a fully configured and accurate numerical propagator can be used to recreate an even
75 * more precise ephemeris in case the initial tabulated values were obtained from an external source.
76 *
77 * @param blendingFunction
78 * {@link org.hipparchus.analysis.polynomials.SmoothStepFactory.SmoothStepFunction smoothstep function} used for
79 * blending
80 * @param blendingPropagator propagator used to propagate tabulated orbits to interpolating time
81 * @param outputInertialFrame output inertial frame
82 *
83 * @throws OrekitException if output frame is not inertial
84 * @see org.hipparchus.analysis.polynomials.SmoothStepFactory.SmoothStepFunction
85 */
86 public OrbitBlender(final SmoothStepFactory.SmoothStepFunction blendingFunction,
87 final Propagator blendingPropagator,
88 final Frame outputInertialFrame) {
89 super(DEFAULT_INTERPOLATION_POINTS, 0., outputInertialFrame);
90 this.blendingFunction = blendingFunction;
91 this.blendingPropagator = blendingPropagator;
92 }
93
94 /** {@inheritDoc} */
95 @Override
96 protected Orbit interpolate(final InterpolationData interpolationData) {
97
98 // Get first and last entry
99 final List<Orbit> neighborList = interpolationData.getNeighborList();
100 final Orbit previousOrbit = neighborList.getFirst();
101 final Orbit nextOrbit = neighborList.get(1);
102
103 // Propagate orbits
104 final AbsoluteDate interpolationDate = interpolationData.getInterpolationDate();
105 final Orbit forwardedOrbit = propagateOrbit(previousOrbit, interpolationDate);
106 final Orbit backwardedOrbit = propagateOrbit(nextOrbit, interpolationDate);
107
108 // Extract position-velocity-acceleration coordinates
109 final PVCoordinates forwardedPV = forwardedOrbit.getPVCoordinates(getOutputInertialFrame());
110 final PVCoordinates backwardedPV = backwardedOrbit.getPVCoordinates(getOutputInertialFrame());
111
112 // Blend PV coordinates
113 final double timeParameter = getTimeParameter(interpolationDate, previousOrbit.getDate(), nextOrbit.getDate());
114 final double blendingValue = blendingFunction.value(timeParameter);
115
116 final PVCoordinates blendedPV = forwardedPV.blendArithmeticallyWith(backwardedPV, blendingValue);
117
118 // Output new blended instance
119 return new CartesianOrbit(blendedPV, getOutputInertialFrame(), interpolationDate, previousOrbit.getMu());
120 }
121
122 /**
123 * Propagate orbit using predefined {@link Propagator propagator}.
124 *
125 * @param tabulatedOrbit tabulated orbit to propagate
126 * @param propagationDate propagation date
127 *
128 * @return orbit propagated to propagation date
129 */
130 private Orbit propagateOrbit(final Orbit tabulatedOrbit,
131 final AbsoluteDate propagationDate) {
132 blendingPropagator.resetInitialState(new SpacecraftState(tabulatedOrbit));
133 return blendingPropagator.propagate(propagationDate).getOrbit();
134 }
135 }