1   /* Copyright 2002-2022 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.semianalytical.dsst;
18  
19  import java.util.IdentityHashMap;
20  import java.util.Map;
21  
22  import org.hipparchus.analysis.differentiation.Gradient;
23  import org.orekit.errors.OrekitException;
24  import org.orekit.errors.OrekitMessages;
25  import org.orekit.propagation.FieldSpacecraftState;
26  import org.orekit.propagation.PropagationType;
27  import org.orekit.propagation.SpacecraftState;
28  import org.orekit.propagation.integration.AdditionalDerivativesProvider;
29  import org.orekit.propagation.semianalytical.dsst.forces.DSSTForceModel;
30  import org.orekit.propagation.semianalytical.dsst.utilities.FieldAuxiliaryElements;
31  import org.orekit.time.AbsoluteDate;
32  import org.orekit.utils.ParameterDriver;
33  import org.orekit.utils.ParameterDriversList;
34  
35  /** {@link AdditionalDerivativesProvider derivatives provider} computing the partial derivatives
36   * of the state (orbit) with respect to initial state and force models parameters.
37   * <p>
38   * This set of equations are automatically added to a {@link DSSTPropagator DSST propagator}
39   * in order to compute partial derivatives of the orbit along with the orbit itself. This is
40   * useful for example in orbit determination applications.
41   * </p>
42   * <p>
43   * The partial derivatives with respect to initial state are dimension 6 (orbit only).
44   * </p>
45   * <p>
46   * The partial derivatives with respect to force models parameters has a dimension
47   * equal to the number of selected parameters. Parameters selection is implemented at
48   * {@link DSSTForceModel DSST force models} level. Users must retrieve a {@link ParameterDriver
49   * parameter driver} by looping on all drivers using {@link DSSTForceModel#getParametersDrivers()}
50   * and then select it by calling {@link ParameterDriver#setSelected(boolean) setSelected(true)}.
51   * </p>
52   * @author Bryan Cazabonne
53   * @since 10.0
54   * @deprecated as of 11.1, replaced by {@link
55   * org.orekit.propagation.Propagator#setupMatricesComputation(String,
56   * org.hipparchus.linear.RealMatrix, org.orekit.utils.DoubleArrayDictionary)}
57   */
58  @Deprecated
59  public class DSSTPartialDerivativesEquations
60      implements AdditionalDerivativesProvider,
61                 org.orekit.propagation.integration.AdditionalEquations {
62  
63      /** Retrograde factor I.
64       *  <p>
65       *  DSST model needs equinoctial orbit as internal representation.
66       *  Classical equinoctial elements have discontinuities when inclination
67       *  is close to zero. In this representation, I = +1. <br>
68       *  To avoid this discontinuity, another representation exists and equinoctial
69       *  elements can be expressed in a different way, called "retrograde" orbit.
70       *  This implies I = -1. <br>
71       *  As Orekit doesn't implement the retrograde orbit, I is always set to +1.
72       *  But for the sake of consistency with the theory, the retrograde factor
73       *  has been kept in the formulas.
74       *  </p>
75       */
76      private static final int I = 1;
77  
78      /** Propagator computing state evolution. */
79      private final DSSTPropagator propagator;
80  
81      /** Selected parameters for Jacobian computation. */
82      private ParameterDriversList selected;
83  
84      /** Parameters map. */
85      private Map<ParameterDriver, Integer> map;
86  
87      /** Name. */
88      private final String name;
89  
90      /** Flag for Jacobian matrices initialization. */
91      private boolean initialized;
92  
93      /** Type of the orbit used for the propagation.*/
94      private PropagationType propagationType;
95  
96      /** Simple constructor.
97       * <p>
98       * Upon construction, this set of equations is <em>automatically</em> added to
99       * the propagator by calling its {@link
100      * DSSTPropagator#addAdditionalDerivativesProvider(AdditionalDerivativesProvider)} method. So
101      * there is no need to call this method explicitly for these equations.
102      * </p>
103      * @param name name of the partial derivatives equations
104      * @param propagator the propagator that will handle the orbit propagation
105      * @param propagationType type of the orbit used for the propagation (mean or osculating)
106      */
107     public DSSTPartialDerivativesEquations(final String name,
108                                            final DSSTPropagator propagator,
109                                            final PropagationType propagationType) {
110         this.name                   = name;
111         this.selected               = null;
112         this.map                    = null;
113         this.propagator             = propagator;
114         this.initialized            = false;
115         this.propagationType        = propagationType;
116         propagator.addAdditionalDerivativesProvider(this);
117     }
118 
119     /** {@inheritDoc} */
120     public String getName() {
121         return name;
122     }
123 
124     /** {@inheritDoc} */
125     @Override
126     public int getDimension() {
127         freezeParametersSelection();
128         return 6 * (6 + selected.getNbParams());
129     }
130 
131     /** Freeze the selected parameters from the force models.
132      */
133     private void freezeParametersSelection() {
134         if (selected == null) {
135 
136             // first pass: gather all parameters, binding similar names together
137             selected = new ParameterDriversList();
138             for (final DSSTForceModel provider : propagator.getAllForceModels()) {
139                 for (final ParameterDriver driver : provider.getParametersDrivers()) {
140                     selected.add(driver);
141                 }
142             }
143 
144             // second pass: now that shared parameter names are bound together,
145             // their selections status have been synchronized, we can filter them
146             selected.filter(true);
147 
148             // third pass: sort parameters lexicographically
149             selected.sort();
150 
151             // fourth pass: set up a map between parameters drivers and matrices columns
152             map = new IdentityHashMap<ParameterDriver, Integer>();
153             int parameterIndex = 0;
154             for (final ParameterDriver selectedDriver : selected.getDrivers()) {
155                 for (final DSSTForceModel provider : propagator.getAllForceModels()) {
156                     for (final ParameterDriver driver : provider.getParametersDrivers()) {
157                         if (driver.getName().equals(selectedDriver.getName())) {
158                             map.put(driver, parameterIndex);
159                         }
160                     }
161                 }
162                 ++parameterIndex;
163             }
164 
165         }
166     }
167 
168     /** Set the initial value of the Jacobian with respect to state and parameter.
169      * <p>
170      * This method is equivalent to call {@link #setInitialJacobians(SpacecraftState,
171      * double[][], double[][])} with dYdY0 set to the identity matrix and dYdP set
172      * to a zero matrix.
173      * </p>
174      * <p>
175      * The force models parameters for which partial derivatives are desired,
176      * <em>must</em> have been {@link ParameterDriver#setSelected(boolean) selected}
177      * before this method is called, so proper matrices dimensions are used.
178      * </p>
179      * @param s0 initial state
180      * @return state with initial Jacobians added
181      */
182     public SpacecraftState setInitialJacobians(final SpacecraftState s0) {
183         freezeParametersSelection();
184         final int stateDimension = 6;
185         final double[][] dYdY0 = new double[stateDimension][stateDimension];
186         final double[][] dYdP  = new double[stateDimension][selected.getNbParams()];
187         for (int i = 0; i < stateDimension; ++i) {
188             dYdY0[i][i] = 1.0;
189         }
190         return setInitialJacobians(s0, dYdY0, dYdP);
191     }
192 
193     /** Set the initial value of the Jacobian with respect to state and parameter.
194      * <p>
195      * The returned state must be added to the propagator (it is not done
196      * automatically, as the user may need to add more states to it).
197      * </p>
198      * <p>
199      * The force models parameters for which partial derivatives are desired,
200      * <em>must</em> have been {@link ParameterDriver#setSelected(boolean) selected}
201      * before this method is called, and the {@code dY1dP} matrix dimension <em>must</em>
202      * be consistent with the selection.
203      * </p>
204      * @param s1 current state
205      * @param dY1dY0 Jacobian of current state at time t₁ with respect
206      * to state at some previous time t₀ (must be 6x6)
207      * @param dY1dP Jacobian of current state at time t₁ with respect
208      * to parameters (may be null if no parameters are selected)
209      * @return state with initial Jacobians added
210      */
211     public SpacecraftState setInitialJacobians(final SpacecraftState s1,
212                                                final double[][] dY1dY0, final double[][] dY1dP) {
213 
214         freezeParametersSelection();
215 
216         // Check dimensions
217         final int stateDim = dY1dY0.length;
218         if (stateDim != 6 || stateDim != dY1dY0[0].length) {
219             throw new OrekitException(OrekitMessages.STATE_JACOBIAN_NOT_6X6,
220                                       stateDim, dY1dY0[0].length);
221         }
222         if (dY1dP != null && stateDim != dY1dP.length) {
223             throw new OrekitException(OrekitMessages.STATE_AND_PARAMETERS_JACOBIANS_ROWS_MISMATCH,
224                                       stateDim, dY1dP.length);
225         }
226         if (dY1dP == null && selected.getNbParams() != 0 ||
227             dY1dP != null && selected.getNbParams() != dY1dP[0].length) {
228             throw new OrekitException(new OrekitException(OrekitMessages.INITIAL_MATRIX_AND_PARAMETERS_NUMBER_MISMATCH,
229                                                           dY1dP == null ? 0 : dY1dP[0].length, selected.getNbParams()));
230         }
231 
232         // store the matrices as a single dimension array
233         initialized = true;
234         final DSSTJacobiansMapper mapper = getMapper();
235         final double[] p = new double[mapper.getAdditionalStateDimension()];
236         mapper.setInitialJacobians(s1, dY1dY0, dY1dP, p);
237 
238         // set value in propagator
239         return s1.addAdditionalState(name, p);
240 
241     }
242 
243     /** Get a mapper between two-dimensional Jacobians and one-dimensional additional state.
244      * @return a mapper between two-dimensional Jacobians and one-dimensional additional state,
245      * with the same name as the instance
246      * @see #setInitialJacobians(SpacecraftState)
247      * @see #setInitialJacobians(SpacecraftState, double[][], double[][])
248      */
249     public DSSTJacobiansMapper getMapper() {
250         if (!initialized) {
251             throw new OrekitException(OrekitMessages.STATE_JACOBIAN_NOT_INITIALIZED);
252         }
253         return new DSSTJacobiansMapper(name, selected, propagator, map, propagationType);
254     }
255 
256     /** {@inheritDoc} */
257     public void init(final SpacecraftState initialState, final AbsoluteDate target) {
258         // FIXME: remove in 12.0 when AdditionalEquations is removed
259         AdditionalDerivativesProvider.super.init(initialState, target);
260     }
261 
262     /** {@inheritDoc} */
263     public double[] computeDerivatives(final SpacecraftState s, final double[] pDot) {
264         // FIXME: remove in 12.0 when AdditionalEquations is removed
265         System.arraycopy(derivatives(s), 0, pDot, 0, pDot.length);
266         return null;
267     }
268 
269     /** {@inheritDoc} */
270     public double[] derivatives(final SpacecraftState s) {
271 
272         // initialize Jacobians to zero
273         final int paramDim = selected.getNbParams();
274         final int dim = 6;
275         final double[][] dMeanElementRatedParam   = new double[dim][paramDim];
276         final double[][] dMeanElementRatedElement = new double[dim][dim];
277         final DSSTGradientConverter converter = new DSSTGradientConverter(s, propagator.getAttitudeProvider());
278 
279         // Compute Jacobian
280         for (final DSSTForceModel forceModel : propagator.getAllForceModels()) {
281 
282             final FieldSpacecraftState<Gradient> dsState = converter.getState(forceModel);
283             final Gradient[] parameters = converter.getParameters(dsState, forceModel);
284             final FieldAuxiliaryElements<Gradient> auxiliaryElements = new FieldAuxiliaryElements<>(dsState.getOrbit(), I);
285 
286             // "field" initialization of the force model if it was not done before
287             forceModel.initializeShortPeriodTerms(auxiliaryElements, propagationType, parameters);
288             final Gradient[] meanElementRate = forceModel.getMeanElementRate(dsState, auxiliaryElements, parameters);
289             final double[] derivativesA  = meanElementRate[0].getGradient();
290             final double[] derivativesEx = meanElementRate[1].getGradient();
291             final double[] derivativesEy = meanElementRate[2].getGradient();
292             final double[] derivativesHx = meanElementRate[3].getGradient();
293             final double[] derivativesHy = meanElementRate[4].getGradient();
294             final double[] derivativesL  = meanElementRate[5].getGradient();
295 
296             // update Jacobian with respect to state
297             addToRow(derivativesA,  0, dMeanElementRatedElement);
298             addToRow(derivativesEx, 1, dMeanElementRatedElement);
299             addToRow(derivativesEy, 2, dMeanElementRatedElement);
300             addToRow(derivativesHx, 3, dMeanElementRatedElement);
301             addToRow(derivativesHy, 4, dMeanElementRatedElement);
302             addToRow(derivativesL,  5, dMeanElementRatedElement);
303 
304             int index = converter.getFreeStateParameters();
305             for (ParameterDriver driver : forceModel.getParametersDrivers()) {
306                 if (driver.isSelected()) {
307                     final int parameterIndex = map.get(driver);
308                     dMeanElementRatedParam[0][parameterIndex] += derivativesA[index];
309                     dMeanElementRatedParam[1][parameterIndex] += derivativesEx[index];
310                     dMeanElementRatedParam[2][parameterIndex] += derivativesEy[index];
311                     dMeanElementRatedParam[3][parameterIndex] += derivativesHx[index];
312                     dMeanElementRatedParam[4][parameterIndex] += derivativesHy[index];
313                     dMeanElementRatedParam[5][parameterIndex] += derivativesL[index];
314                     ++index;
315                 }
316             }
317 
318         }
319 
320         // The variational equations of the complete state Jacobian matrix have the following form:
321 
322         //                     [ Adot ] = [ dMeanElementRatedElement ] * [ A ]
323 
324         // The A matrix and its derivative (Adot) are 6 * 6 matrices
325 
326         // The following loops compute these expression taking care of the mapping of the
327         // A matrix into the single dimension array p and of the mapping of the
328         // Adot matrix into the single dimension array pDot.
329 
330         final double[] p = s.getAdditionalState(getName());
331         final double[] pDot = new double[p.length];
332 
333         for (int i = 0; i < dim; i++) {
334             final double[] dMeanElementRatedElementi = dMeanElementRatedElement[i];
335             for (int j = 0; j < dim; j++) {
336                 pDot[j + dim * i] =
337                     dMeanElementRatedElementi[0] * p[j]           + dMeanElementRatedElementi[1] * p[j +     dim] + dMeanElementRatedElementi[2] * p[j + 2 * dim] +
338                     dMeanElementRatedElementi[3] * p[j + 3 * dim] + dMeanElementRatedElementi[4] * p[j + 4 * dim] + dMeanElementRatedElementi[5] * p[j + 5 * dim];
339             }
340         }
341 
342         final int columnTop = dim * dim;
343         for (int k = 0; k < paramDim; k++) {
344             // the variational equations of the parameters Jacobian matrix are computed
345             // one column at a time, they have the following form:
346 
347             //             [ Bdot ] = [ dMeanElementRatedElement ] * [ B ] + [ dMeanElementRatedParam ]
348 
349             // The B sub-columns and its derivative (Bdot) are 6 elements columns.
350 
351             // The following loops compute this expression taking care of the mapping of the
352             // B columns into the single dimension array p and of the mapping of the
353             // Bdot columns into the single dimension array pDot.
354 
355             for (int i = 0; i < dim; ++i) {
356                 final double[] dMeanElementRatedElementi = dMeanElementRatedElement[i];
357                 pDot[columnTop + (i + dim * k)] =
358                     dMeanElementRatedParam[i][k] +
359                     dMeanElementRatedElementi[0] * p[columnTop + k]                + dMeanElementRatedElementi[1] * p[columnTop + k +     paramDim] + dMeanElementRatedElementi[2] * p[columnTop + k + 2 * paramDim] +
360                     dMeanElementRatedElementi[3] * p[columnTop + k + 3 * paramDim] + dMeanElementRatedElementi[4] * p[columnTop + k + 4 * paramDim] + dMeanElementRatedElementi[5] * p[columnTop + k + 5 * paramDim];
361             }
362         }
363 
364         return pDot;
365 
366     }
367 
368     /** Fill Jacobians rows.
369      * @param derivatives derivatives of a component
370      * @param index component index (0 for a, 1 for ex, 2 for ey, 3 for hx, 4 for hy, 5 for l)
371      * @param dMeanElementRatedElement Jacobian of mean elements rate with respect to mean elements
372      */
373     private void addToRow(final double[] derivatives, final int index,
374                           final double[][] dMeanElementRatedElement) {
375 
376         for (int i = 0; i < 6; i++) {
377             dMeanElementRatedElement[index][i] += derivatives[i];
378         }
379 
380     }
381 
382 }