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