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.analytical;
18  
19  import java.util.Arrays;
20  import java.util.List;
21  
22  import org.hipparchus.analysis.differentiation.Gradient;
23  import org.hipparchus.linear.MatrixUtils;
24  import org.hipparchus.linear.RealMatrix;
25  import org.orekit.orbits.FieldOrbit;
26  import org.orekit.orbits.OrbitType;
27  import org.orekit.orbits.PositionAngleType;
28  import org.orekit.propagation.AbstractMatricesHarvester;
29  import org.orekit.propagation.AdditionalDataProvider;
30  import org.orekit.propagation.FieldSpacecraftState;
31  import org.orekit.propagation.SpacecraftState;
32  import org.orekit.time.AbsoluteDate;
33  import org.orekit.time.FieldAbsoluteDate;
34  import org.orekit.utils.DoubleArrayDictionary;
35  import org.orekit.utils.FieldPVCoordinates;
36  import org.orekit.utils.ParameterDriver;
37  import org.orekit.utils.TimeSpanMap;
38  import org.orekit.utils.TimeSpanMap.Span;
39  
40  /**
41   * Base class harvester between two-dimensional Jacobian
42   * matrices and analytical orbit propagator.
43   * @author Thomas Paulet
44   * @author Bryan Cazabonne
45   * @since 11.1
46   */
47  public abstract class AbstractAnalyticalMatricesHarvester extends AbstractMatricesHarvester implements AdditionalDataProvider<double[]> {
48  
49      /** Columns names for parameters. */
50      private List<String> columnsNames;
51  
52      /** Epoch of the last computed state transition matrix. */
53      private AbsoluteDate epoch;
54  
55      /** Analytical derivatives that apply to State Transition Matrix. */
56      private final double[][] analyticalDerivativesStm;
57  
58      /** Analytical derivatives that apply to Jacobians columns. */
59      private final DoubleArrayDictionary analyticalDerivativesJacobianColumns;
60  
61      /** Propagator bound to this harvester. */
62      private final AbstractAnalyticalPropagator propagator;
63  
64      /** Simple constructor.
65       * <p>
66       * The arguments for initial matrices <em>must</em> be compatible with the
67       * {@link org.orekit.orbits.OrbitType orbit type}
68       * and {@link PositionAngleType position angle} that will be used by propagator
69       * </p>
70       * @param propagator propagator bound to this harvester
71       * @param stmName State Transition Matrix state name
72       * @param initialStm initial State Transition Matrix ∂Y/∂Y₀,
73       * if null (which is the most frequent case), assumed to be 6x6 identity
74       * @param initialJacobianColumns initial columns of the Jacobians matrix with respect to parameters,
75       * if null or if some selected parameters are missing from the dictionary, the corresponding
76       * initial column is assumed to be 0
77       */
78      protected AbstractAnalyticalMatricesHarvester(final AbstractAnalyticalPropagator propagator, final String stmName,
79                                                    final RealMatrix initialStm, final DoubleArrayDictionary initialJacobianColumns) {
80          super(stmName, initialStm, initialJacobianColumns);
81          this.propagator                           = propagator;
82          this.epoch                                = propagator.getInitialState().getDate();
83          this.columnsNames                         = null;
84          this.analyticalDerivativesStm             = getInitialStateTransitionMatrix().getData();
85          this.analyticalDerivativesJacobianColumns = new DoubleArrayDictionary();
86      }
87  
88      /** {@inheritDoc} */
89      @Override
90      public List<String> getJacobiansColumnsNames() {
91          return columnsNames == null ? propagator.getJacobiansColumnsNames() : columnsNames;
92      }
93  
94      /** {@inheritDoc} */
95      @Override
96      public void freezeColumnsNames() {
97          columnsNames = getJacobiansColumnsNames();
98      }
99  
100     /** {@inheritDoc} */
101     @Override
102     public String getName() {
103         return getStmName();
104     }
105 
106     /** {@inheritDoc} */
107     @Override
108     public double[] getAdditionalData(final SpacecraftState state) {
109         // Update the partial derivatives if needed
110         updateDerivativesIfNeeded(state);
111         // Return the state transition matrix in an array
112         return toArray(analyticalDerivativesStm);
113     }
114 
115     /** {@inheritDoc} */
116     @Override
117     public RealMatrix getStateTransitionMatrix(final SpacecraftState state) {
118         // Check if additional state is defined
119         if (!state.hasAdditionalData(getName())) {
120             return null;
121         }
122         // Return the state transition matrix
123         return toRealMatrix(state.getAdditionalState(getName()));
124     }
125 
126     /** {@inheritDoc} */
127     @Override
128     public RealMatrix getParametersJacobian(final SpacecraftState state) {
129         // Update the partial derivatives if needed
130         updateDerivativesIfNeeded(state);
131 
132         // Estimated parameters
133         final List<String> names = getJacobiansColumnsNames();
134         if (names == null || names.isEmpty()) {
135             return null;
136         }
137 
138         // Initialize Jacobian
139         final RealMatrix dYdP = MatrixUtils.createRealMatrix(STATE_DIMENSION, names.size());
140 
141         // Add the derivatives
142         for (int j = 0; j < names.size(); ++j) {
143             final double[] column = analyticalDerivativesJacobianColumns.get(names.get(j));
144             if (column != null) {
145                 for (int i = 0; i < STATE_DIMENSION; i++) {
146                     dYdP.addToEntry(i, j, column[i]);
147                 }
148             }
149         }
150 
151         // Return
152         return dYdP;
153     }
154 
155     /** {@inheritDoc} */
156     @Override
157     public void setReferenceState(final SpacecraftState reference) {
158 
159         // reset derivatives to zero
160         for (final double[] row : analyticalDerivativesStm) {
161             Arrays.fill(row, 0.0);
162         }
163         analyticalDerivativesJacobianColumns.clear();
164 
165         final AbstractAnalyticalGradientConverter converter           = getGradientConverter();
166         final FieldAbstractAnalyticalPropagator<Gradient> gPropagator = converter.getPropagator();
167 
168         // Compute Jacobian
169         final AbsoluteDate                   target     = reference.getDate();
170         final FieldAbsoluteDate<Gradient>    start      = gPropagator.getInitialState().getDate();
171         final double                         dt         = target.durationFrom(start.toAbsoluteDate());
172         final FieldSpacecraftState<Gradient> state      = gPropagator.getInitialState();
173         final Gradient[]                     parameters = converter.getParameters(state, converter);
174         final FieldOrbit<Gradient>           gOrbit     = gPropagator.propagateOrbit(start.shiftedBy(dt), parameters);
175         final FieldPVCoordinates<Gradient>   gPv        = gOrbit.getPVCoordinates();
176 
177         final double[] derivativesX  = gPv.getPosition().getX().getGradient();
178         final double[] derivativesY  = gPv.getPosition().getY().getGradient();
179         final double[] derivativesZ  = gPv.getPosition().getZ().getGradient();
180         final double[] derivativesVx = gPv.getVelocity().getX().getGradient();
181         final double[] derivativesVy = gPv.getVelocity().getY().getGradient();
182         final double[] derivativesVz = gPv.getVelocity().getZ().getGradient();
183 
184         // Update Jacobian with respect to state
185         addToRow(derivativesX,  0);
186         addToRow(derivativesY,  1);
187         addToRow(derivativesZ,  2);
188         addToRow(derivativesVx, 3);
189         addToRow(derivativesVy, 4);
190         addToRow(derivativesVz, 5);
191 
192         // Partial derivatives of the state with respect to propagation parameters
193         int paramsIndex = converter.getFreeStateParameters();
194         for (ParameterDriver driver : converter.getParametersDrivers()) {
195             if (driver.isSelected()) {
196 
197                 final TimeSpanMap<String> driverNameSpanMap = driver.getNamesSpanMap();
198                 // for each span (for each estimated value) corresponding name is added
199                 for (Span<String> span = driverNameSpanMap.getFirstSpan(); span != null; span = span.next()) {
200                     // get the partials derivatives for this driver
201                     DoubleArrayDictionary.Entry entry = analyticalDerivativesJacobianColumns.getEntry(span.getData());
202                     if (entry == null) {
203                         // create an entry filled with zeroes
204                         analyticalDerivativesJacobianColumns.put(span.getData(), new double[STATE_DIMENSION]);
205                         entry = analyticalDerivativesJacobianColumns.getEntry(span.getData());
206                     }
207 
208                     // add the contribution of the current force model
209                     entry.increment(new double[] {
210                         derivativesX[paramsIndex], derivativesY[paramsIndex], derivativesZ[paramsIndex],
211                         derivativesVx[paramsIndex], derivativesVy[paramsIndex], derivativesVz[paramsIndex]
212                     });
213                     ++paramsIndex;
214                 }
215             }
216         }
217 
218         // Update the epoch of the last computed partial derivatives
219         epoch = target;
220 
221     }
222 
223     /** Update the partial derivatives (if needed).
224      * @param state current spacecraft state
225      */
226     private void updateDerivativesIfNeeded(final SpacecraftState state) {
227         if (!state.getDate().isEqualTo(epoch)) {
228             setReferenceState(state);
229         }
230     }
231 
232     /** Fill State Transition Matrix rows.
233      * @param derivatives derivatives of a component
234      * @param index component index
235      */
236     private void addToRow(final double[] derivatives, final int index) {
237         for (int i = 0; i < 6; i++) {
238             analyticalDerivativesStm[index][i] += derivatives[i];
239         }
240     }
241 
242     /** Convert an array to a matrix (6x6 dimension).
243      * @param array input array
244      * @return the corresponding matrix
245      */
246     private RealMatrix toRealMatrix(final double[] array) {
247         final RealMatrix matrix = MatrixUtils.createRealMatrix(STATE_DIMENSION, STATE_DIMENSION);
248         int index = 0;
249         for (int i = 0; i < STATE_DIMENSION; ++i) {
250             for (int j = 0; j < STATE_DIMENSION; ++j) {
251                 matrix.setEntry(i, j, array[index++]);
252             }
253         }
254         return matrix;
255     }
256 
257     /** Set the STM data into an array.
258      * @param matrix STM matrix
259      * @return an array containing the STM data
260      */
261     private double[] toArray(final double[][] matrix) {
262         final double[] array = new double[STATE_DIMENSION * STATE_DIMENSION];
263         int index = 0;
264         for (int i = 0; i < STATE_DIMENSION; ++i) {
265             final double[] row = matrix[i];
266             for (int j = 0; j < STATE_DIMENSION; ++j) {
267                 array[index++] = row[j];
268             }
269         }
270         return array;
271     }
272 
273     /** {@inheritDoc} */
274     @Override
275     public OrbitType getOrbitType() {
276         // Set to CARTESIAN because analytical gradient converters uses Cartesian representation
277         return OrbitType.CARTESIAN;
278     }
279 
280     /** {@inheritDoc} */
281     @Override
282     public PositionAngleType getPositionAngleType() {
283         // Irrelevant: set a default value
284         return PositionAngleType.MEAN;
285     }
286 
287     /**
288      * Get the gradient converter related to the analytical orbit propagator.
289      * @return the gradient converter
290      */
291     public abstract AbstractAnalyticalGradientConverter getGradientConverter();
292 
293 }