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.numerical;
18  
19  import org.hipparchus.analysis.differentiation.Gradient;
20  import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
21  import org.hipparchus.linear.Array2DRowRealMatrix;
22  import org.hipparchus.linear.DecompositionSolver;
23  import org.hipparchus.linear.MatrixUtils;
24  import org.hipparchus.linear.QRDecomposition;
25  import org.hipparchus.linear.RealMatrix;
26  import org.orekit.errors.OrekitException;
27  import org.orekit.errors.OrekitMessages;
28  import org.orekit.forces.ForceModel;
29  import org.orekit.forces.gravity.ThirdBodyAttractionEpoch;
30  import org.orekit.propagation.FieldSpacecraftState;
31  import org.orekit.propagation.SpacecraftState;
32  import org.orekit.propagation.integration.AdditionalDerivativesProvider;
33  import org.orekit.propagation.integration.CombinedDerivatives;
34  import org.orekit.utils.ParameterDriver;
35  import org.orekit.utils.ParameterDriversList;
36  import org.orekit.utils.TimeSpanMap.Span;
37  
38  import java.util.IdentityHashMap;
39  import java.util.Map;
40  
41  /** Computes derivatives of the acceleration, including ThirdBodyAttraction.
42   *
43   * {@link AdditionalDerivativesProvider Provider} computing the partial derivatives
44   * of the state (orbit) with respect to initial state and force models parameters.
45   * <p>
46   * This set of equations are automatically added to a {@link NumericalPropagator numerical propagator}
47   * in order to compute partial derivatives of the orbit along with the orbit itself. This is
48   * useful for example in orbit determination applications.
49   * </p>
50   * <p>
51   * The partial derivatives with respect to initial state can be either dimension 6
52   * (orbit only) or 7 (orbit and mass).
53   * </p>
54   * <p>
55   * The partial derivatives with respect to force models parameters has a dimension
56   * equal to the number of selected parameters. Parameters selection is implemented at
57   * {@link ForceModel force models} level. Users must retrieve a {@link ParameterDriver
58   * parameter driver} using {@link ForceModel#getParameterDriver(String)} and then
59   * select it by calling {@link ParameterDriver#setSelected(boolean) setSelected(true)}.
60   * </p>
61   * <p>
62   * If several force models provide different {@link ParameterDriver drivers} for the
63   * same parameter name, selecting any of these drivers has the side effect of
64   * selecting all the drivers for this shared parameter. In this case, the partial
65   * derivatives will be the sum of the partial derivatives contributed by the
66   * corresponding force models. This case typically arises for central attraction
67   * coefficient, which has an influence on {@link org.orekit.forces.gravity.NewtonianAttraction
68   * Newtonian attraction}, {@link org.orekit.forces.gravity.HolmesFeatherstoneAttractionModel
69   * gravity field}, and {@link org.orekit.forces.gravity.Relativity relativity}.
70   * </p>
71   * @author V&eacute;ronique Pommier-Maurussane
72   * @author Luc Maisonobe
73   * @since 10.2
74   */
75  public class EpochDerivativesEquations
76      implements AdditionalDerivativesProvider  {
77  
78      /** State dimension, fixed to 6. */
79      public static final int STATE_DIMENSION = 6;
80  
81      /** Propagator computing state evolution. */
82      private final NumericalPropagator propagator;
83  
84      /** Selected parameters for Jacobian computation. */
85      private ParameterDriversList selected;
86  
87      /** Parameters map. */
88      private Map<String, Integer> map;
89  
90      /** Name. */
91      private final String name;
92  
93      /** Simple constructor.
94       * <p>
95       * Upon construction, this set of equations is <em>automatically</em> added to
96       * the propagator by calling its {@link
97       * NumericalPropagator#addAdditionalDerivativesProvider(AdditionalDerivativesProvider)} method. So
98       * there is no need to call this method explicitly for these equations.
99       * </p>
100      * @param name name of the partial derivatives equations
101      * @param propagator the propagator that will handle the orbit propagation
102      */
103     public EpochDerivativesEquations(final String name, final NumericalPropagator propagator) {
104         this.name                   = name;
105         this.selected               = null;
106         this.map                    = null;
107         this.propagator             = propagator;
108         propagator.addAdditionalDerivativesProvider(this);
109     }
110 
111     /** {@inheritDoc} */
112     public String getName() {
113         return name;
114     }
115 
116     /** {@inheritDoc} */
117     @Override
118     public int getDimension() {
119         freezeParametersSelection();
120         return 6 * (6 + selected.getNbParams() + 1);
121     }
122 
123     /** Freeze the selected parameters from the force models.
124      */
125     private void freezeParametersSelection() {
126         if (selected == null) {
127 
128             // first pass: gather all parameters, binding similar names together
129             selected = new ParameterDriversList();
130             for (final ForceModel provider : propagator.getAllForceModels()) {
131                 for (final ParameterDriver driver : provider.getParametersDrivers()) {
132                     selected.add(driver);
133                 }
134             }
135 
136             // second pass: now that shared parameter names are bound together,
137             // their selections status have been synchronized, we can filter them
138             selected.filter(true);
139 
140             // third pass: sort parameters lexicographically
141             selected.sort();
142 
143             // fourth pass: set up a map between parameters drivers and matrices columns
144             map = new IdentityHashMap<>();
145             int parameterIndex = 0;
146             int previousParameterIndex = 0;
147             for (final ParameterDriver selectedDriver : selected.getDrivers()) {
148                 for (final ForceModel provider : propagator.getAllForceModels()) {
149                     for (final ParameterDriver driver : provider.getParametersDrivers()) {
150                         if (driver.getName().equals(selectedDriver.getName())) {
151                             previousParameterIndex = parameterIndex;
152                             for (Span<String> span = driver.getNamesSpanMap().getFirstSpan(); span != null; span = span.next()) {
153                                 map.put(span.getData(), previousParameterIndex++);
154                             }
155                         }
156                     }
157                 }
158                 parameterIndex = previousParameterIndex;
159             }
160 
161         }
162     }
163 
164     /** Set the initial value of the Jacobian with respect to state and parameter.
165      * <p>
166      * This method is equivalent to call {@link #setInitialJacobians(SpacecraftState,
167      * double[][], double[][])} with dYdY0 set to the identity matrix and dYdP set
168      * to a zero matrix.
169      * </p>
170      * <p>
171      * The force models parameters for which partial derivatives are desired,
172      * <em>must</em> have been {@link ParameterDriver#setSelected(boolean) selected}
173      * before this method is called, so proper matrices dimensions are used.
174      * </p>
175      * @param s0 initial state
176      * @return state with initial Jacobians added
177      */
178     public SpacecraftState setInitialJacobians(final SpacecraftState s0) {
179         freezeParametersSelection();
180         final int epochStateDimension = 6;
181         final double[][] dYdY0 = new double[epochStateDimension][epochStateDimension];
182         final double[][] dYdP  = new double[epochStateDimension][selected.getNbValuesToEstimate() + 6];
183         for (int i = 0; i < epochStateDimension; ++i) {
184             dYdY0[i][i] = 1.0;
185         }
186         return setInitialJacobians(s0, dYdY0, dYdP);
187     }
188 
189     /** Set the initial value of the Jacobian with respect to state and parameter.
190      * <p>
191      * The returned state must be added to the propagator (it is not done
192      * automatically, as the user may need to add more states to it).
193      * </p>
194      * <p>
195      * The force models parameters for which partial derivatives are desired,
196      * <em>must</em> have been {@link ParameterDriver#setSelected(boolean) selected}
197      * before this method is called, and the {@code dY1dP} matrix dimension <em>must</em>
198      * be consistent with the selection.
199      * </p>
200      * @param s1 current state
201      * @param dY1dY0 Jacobian of current state at time t₁ with respect
202      * to state at some previous time t₀ (must be 6x6)
203      * @param dY1dP Jacobian of current state at time t₁ with respect
204      * to parameters (may be null if no parameters are selected)
205      * @return state with initial Jacobians added
206      */
207     public SpacecraftState setInitialJacobians(final SpacecraftState s1,
208                                                final double[][] dY1dY0, final double[][] dY1dP) {
209 
210         freezeParametersSelection();
211 
212         // Check dimensions
213         final int stateDimEpoch = dY1dY0.length;
214         if (stateDimEpoch != 6 || stateDimEpoch != dY1dY0[0].length) {
215             throw new OrekitException(OrekitMessages.STATE_JACOBIAN_NOT_6X6,
216                                       stateDimEpoch, dY1dY0[0].length);
217         }
218         if (dY1dP != null && stateDimEpoch != dY1dP.length) {
219             throw new OrekitException(OrekitMessages.STATE_AND_PARAMETERS_JACOBIANS_ROWS_MISMATCH,
220                                       stateDimEpoch, dY1dP.length);
221         }
222 
223         // store the matrices as a single dimension array
224         final double[] p = new double[STATE_DIMENSION * (STATE_DIMENSION + selected.getNbValuesToEstimate()) + 6];
225         setInitialJacobians(s1, dY1dY0, dY1dP, p);
226 
227         // set value in propagator
228         return s1.addAdditionalData(name, p);
229 
230     }
231 
232     /** Set the Jacobian with respect to state into a one-dimensional additional state array.
233      * <p>
234      * This method converts the Jacobians to Cartesian parameters and put the converted data
235      * in the one-dimensional {@code p} array.
236      * </p>
237      * @param state spacecraft state
238      * @param dY1dY0 Jacobian of current state at time t₁
239      * with respect to state at some previous time t₀
240      * @param dY1dP Jacobian of current state at time t₁
241      * with respect to parameters (may be null if there are no parameters)
242      * @param p placeholder where to put the one-dimensional additional state
243      */
244     public void setInitialJacobians(final SpacecraftState state, final double[][] dY1dY0,
245                                     final double[][] dY1dP, final double[] p) {
246 
247         // set up a converter
248         final RealMatrix dY1dC1 = MatrixUtils.createRealIdentityMatrix(STATE_DIMENSION);
249         final DecompositionSolver solver = new QRDecomposition(dY1dC1).getSolver();
250 
251         // convert the provided state Jacobian
252         final RealMatrix dC1dY0 = solver.solve(new Array2DRowRealMatrix(dY1dY0, false));
253 
254         // map the converted state Jacobian to one-dimensional array
255         int index = 0;
256         for (int i = 0; i < STATE_DIMENSION; ++i) {
257             for (int j = 0; j < STATE_DIMENSION; ++j) {
258                 p[index++] = dC1dY0.getEntry(i, j);
259             }
260         }
261 
262         if (selected.getNbValuesToEstimate() != 0) {
263             // convert the provided state Jacobian
264             final RealMatrix dC1dP = solver.solve(new Array2DRowRealMatrix(dY1dP, false));
265 
266             // map the converted parameters Jacobian to one-dimensional array
267             for (int i = 0; i < STATE_DIMENSION; ++i) {
268                 for (int j = 0; j < selected.getNbValuesToEstimate(); ++j) {
269                     p[index++] = dC1dP.getEntry(i, j);
270                 }
271             }
272         }
273 
274     }
275 
276     /** {@inheritDoc} */
277     public CombinedDerivatives combinedDerivatives(final SpacecraftState s) {
278 
279         // initialize acceleration Jacobians to zero
280         final int paramDimEpoch = selected.getNbValuesToEstimate() + 1; // added epoch
281         final int dimEpoch      = 3;
282         final double[][] dAccdParam = new double[dimEpoch][paramDimEpoch];
283         final double[][] dAccdPos   = new double[dimEpoch][dimEpoch];
284         final double[][] dAccdVel   = new double[dimEpoch][dimEpoch];
285 
286         final NumericalGradientConverter fullConverter    = new NumericalGradientConverter(s, 6, propagator.getAttitudeProvider());
287         final NumericalGradientConverter posOnlyConverter = new NumericalGradientConverter(s, 3, propagator.getAttitudeProvider());
288 
289         // compute acceleration Jacobians, finishing with the largest force: Newtonian attraction
290         for (final ForceModel forceModel : propagator.getAllForceModels()) {
291             final NumericalGradientConverter converter = forceModel.dependsOnPositionOnly() ? posOnlyConverter : fullConverter;
292             final FieldSpacecraftState<Gradient> dsState = converter.getState(forceModel);
293             final Gradient[] parameters = converter.getParametersAtStateDate(dsState, forceModel);
294 
295             final FieldVector3D<Gradient> acceleration = forceModel.acceleration(dsState, parameters);
296             final double[] derivativesX = acceleration.getX().getGradient();
297             final double[] derivativesY = acceleration.getY().getGradient();
298             final double[] derivativesZ = acceleration.getZ().getGradient();
299 
300             // update Jacobians with respect to state
301             addToRow(derivativesX, 0, converter.getFreeStateParameters(), dAccdPos, dAccdVel);
302             addToRow(derivativesY, 1, converter.getFreeStateParameters(), dAccdPos, dAccdVel);
303             addToRow(derivativesZ, 2, converter.getFreeStateParameters(), dAccdPos, dAccdVel);
304 
305             int index = converter.getFreeStateParameters();
306             for (ParameterDriver driver : forceModel.getParametersDrivers()) {
307                 if (driver.isSelected()) {
308                     for (Span<String> span = driver.getNamesSpanMap().getFirstSpan(); span != null; span = span.next()) {
309                         final int parameterIndex = map.get(span.getData());
310                         dAccdParam[0][parameterIndex] += derivativesX[index];
311                         dAccdParam[1][parameterIndex] += derivativesY[index];
312                         dAccdParam[2][parameterIndex] += derivativesZ[index];
313                         ++index;
314                     }
315                 }
316             }
317 
318             // Add the derivatives of the acceleration w.r.t. the Epoch
319             if (forceModel instanceof ThirdBodyAttractionEpoch) {
320                 final double[] parametersValues = new double[] {parameters[0].getValue()};
321                 final double[] derivatives = ((ThirdBodyAttractionEpoch) forceModel).getDerivativesToEpoch(s, parametersValues);
322                 dAccdParam[0][paramDimEpoch - 1] += derivatives[0];
323                 dAccdParam[1][paramDimEpoch - 1] += derivatives[1];
324                 dAccdParam[2][paramDimEpoch - 1] += derivatives[2];
325             }
326 
327         }
328 
329         // the variational equations of the complete state Jacobian matrix have the following form:
330 
331         // [        |        ]   [                 |                  ]   [     |     ]
332         // [  Adot  |  Bdot  ]   [  dVel/dPos = 0  |  dVel/dVel = Id  ]   [  A  |  B  ]
333         // [        |        ]   [                 |                  ]   [     |     ]
334         // ---------+---------   ------------------+------------------- * ------+------
335         // [        |        ]   [                 |                  ]   [     |     ]
336         // [  Cdot  |  Ddot  ] = [    dAcc/dPos    |     dAcc/dVel    ]   [  C  |  D  ]
337         // [        |        ]   [                 |                  ]   [     |     ]
338 
339         // The A, B, C and D sub-matrices and their derivatives (Adot ...) are 3x3 matrices
340 
341         // The expanded multiplication above can be rewritten to take into account
342         // the fixed values found in the sub-matrices in the left factor. This leads to:
343 
344         //     [ Adot ] = [ C ]
345         //     [ Bdot ] = [ D ]
346         //     [ Cdot ] = [ dAcc/dPos ] * [ A ] + [ dAcc/dVel ] * [ C ]
347         //     [ Ddot ] = [ dAcc/dPos ] * [ B ] + [ dAcc/dVel ] * [ D ]
348 
349         // The following loops compute these expressions taking care of the mapping of the
350         // (A, B, C, D) matrices into the single dimension array p and of the mapping of the
351         // (Adot, Bdot, Cdot, Ddot) matrices into the single dimension array pDot.
352 
353         // copy C and E into Adot and Bdot
354         final int stateDim = 6;
355         final double[] p = s.getAdditionalState(getName());
356         final double[] pDot = new double[p.length];
357         System.arraycopy(p, dimEpoch * stateDim, pDot, 0, dimEpoch * stateDim);
358 
359         // compute Cdot and Ddot
360         for (int i = 0; i < dimEpoch; ++i) {
361             final double[] dAdPi = dAccdPos[i];
362             final double[] dAdVi = dAccdVel[i];
363             for (int j = 0; j < stateDim; ++j) {
364                 pDot[(dimEpoch + i) * stateDim + j] =
365                     dAdPi[0] * p[j]                + dAdPi[1] * p[j +     stateDim] + dAdPi[2] * p[j + 2 * stateDim] +
366                     dAdVi[0] * p[j + 3 * stateDim] + dAdVi[1] * p[j + 4 * stateDim] + dAdVi[2] * p[j + 5 * stateDim];
367             }
368         }
369 
370         for (int k = 0; k < paramDimEpoch; ++k) {
371             // the variational equations of the parameters Jacobian matrix are computed
372             // one column at a time, they have the following form:
373             // [      ]   [                 |                  ]   [   ]   [                  ]
374             // [ Edot ]   [  dVel/dPos = 0  |  dVel/dVel = Id  ]   [ E ]   [  dVel/dParam = 0 ]
375             // [      ]   [                 |                  ]   [   ]   [                  ]
376             // --------   ------------------+------------------- * ----- + --------------------
377             // [      ]   [                 |                  ]   [   ]   [                  ]
378             // [ Fdot ] = [    dAcc/dPos    |     dAcc/dVel    ]   [ F ]   [    dAcc/dParam   ]
379             // [      ]   [                 |                  ]   [   ]   [                  ]
380 
381             // The E and F sub-columns and their derivatives (Edot, Fdot) are 3 elements columns.
382 
383             // The expanded multiplication and addition above can be rewritten to take into
384             // account the fixed values found in the sub-matrices in the left factor. This leads to:
385 
386             //     [ Edot ] = [ F ]
387             //     [ Fdot ] = [ dAcc/dPos ] * [ E ] + [ dAcc/dVel ] * [ F ] + [ dAcc/dParam ]
388 
389             // The following loops compute these expressions taking care of the mapping of the
390             // (E, F) columns into the single dimension array p and of the mapping of the
391             // (Edot, Fdot) columns into the single dimension array pDot.
392 
393             // copy F into Edot
394             final int columnTop = stateDim * stateDim + k;
395             pDot[columnTop]                     = p[columnTop + 3 * paramDimEpoch];
396             pDot[columnTop +     paramDimEpoch] = p[columnTop + 4 * paramDimEpoch];
397             pDot[columnTop + 2 * paramDimEpoch] = p[columnTop + 5 * paramDimEpoch];
398 
399             // compute Fdot
400             for (int i = 0; i < dimEpoch; ++i) {
401                 final double[] dAdP = dAccdPos[i];
402                 final double[] dAdV = dAccdVel[i];
403                 pDot[columnTop + (dimEpoch + i) * paramDimEpoch] =
404                     dAccdParam[i][k] +
405                     dAdP[0] * p[columnTop]                     + dAdP[1] * p[columnTop +     paramDimEpoch] + dAdP[2] * p[columnTop + 2 * paramDimEpoch] +
406                     dAdV[0] * p[columnTop + 3 * paramDimEpoch] + dAdV[1] * p[columnTop + 4 * paramDimEpoch] + dAdV[2] * p[columnTop + 5 * paramDimEpoch];
407             }
408 
409         }
410 
411         return new CombinedDerivatives(pDot, null);
412 
413     }
414 
415     /** Fill Jacobians rows.
416      * @param derivatives derivatives of a component of acceleration (along either x, y or z)
417      * @param index component index (0 for x, 1 for y, 2 for z)
418      * @param freeStateParameters number of free parameters, either 3 (position),
419      * 6 (position-velocity) or 7 (position-velocity-mass)
420      * @param dAccdPos Jacobian of acceleration with respect to spacecraft position
421      * @param dAccdVel Jacobian of acceleration with respect to spacecraft velocity
422      */
423     private void addToRow(final double[] derivatives, final int index, final int freeStateParameters,
424                           final double[][] dAccdPos, final double[][] dAccdVel) {
425 
426         for (int i = 0; i < 3; ++i) {
427             dAccdPos[index][i] += derivatives[i];
428         }
429         if (freeStateParameters > 3) {
430             for (int i = 0; i < 3; ++i) {
431                 dAccdVel[index][i] += derivatives[i + 3];
432             }
433         }
434 
435     }
436 
437 }
438