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.semianalytical.dsst;
18  
19  import org.hipparchus.analysis.differentiation.Gradient;
20  import org.hipparchus.exception.LocalizedCoreFormats;
21  import org.hipparchus.linear.MatrixUtils;
22  import org.hipparchus.linear.RealMatrix;
23  import org.orekit.attitudes.AttitudeProvider;
24  import org.orekit.errors.OrekitException;
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.integration.CombinedDerivatives;
30  import org.orekit.propagation.semianalytical.dsst.forces.DSSTForceModel;
31  import org.orekit.propagation.semianalytical.dsst.utilities.FieldAuxiliaryElements;
32  import org.orekit.time.AbsoluteDate;
33  import org.orekit.utils.DoubleArrayDictionary;
34  import org.orekit.utils.ParameterDriver;
35  import org.orekit.utils.TimeSpanMap.Span;
36  
37  import java.util.HashMap;
38  import java.util.List;
39  import java.util.Map;
40  
41  /** Generator for State Transition Matrix.
42   * @author Luc Maisonobe
43   * @since 11.1
44   */
45  class DSSTStateTransitionMatrixGenerator implements AdditionalDerivativesProvider {
46  
47      /** Space dimension. */
48      private static final int SPACE_DIMENSION = 3;
49  
50      /** Retrograde factor I.
51       *  <p>
52       *  DSST model needs equinoctial orbit as internal representation.
53       *  Classical equinoctial elements have discontinuities when inclination
54       *  is close to zero. In this representation, I = +1. <br>
55       *  To avoid this discontinuity, another representation exists and equinoctial
56       *  elements can be expressed in a different way, called "retrograde" orbit.
57       *  This implies I = -1. <br>
58       *  As Orekit doesn't implement the retrograde orbit, I is always set to +1.
59       *  But for the sake of consistency with the theory, the retrograde factor
60       *  has been kept in the formulas.
61       *  </p>
62       */
63      private static final int I = 1;
64  
65      /** State dimension. */
66      public static final int STATE_DIMENSION = 2 * SPACE_DIMENSION;
67  
68      /** Name of the Cartesian STM additional state. */
69      private final String stmName;
70  
71      /** Force models used in propagation. */
72      private final List<DSSTForceModel> forceModels;
73  
74      /** Attitude provider used in propagation. */
75      private final AttitudeProvider attitudeProvider;
76  
77      /** Observers for partial derivatives. */
78      private final Map<String, DSSTPartialsObserver> partialsObservers;
79  
80      /** Mean or osculating. */
81      private final PropagationType propagationType;
82  
83      /** Simple constructor.
84       * @param stmName name of the Cartesian STM additional state
85       * @param forceModels force models used in propagation
86       * @param attitudeProvider attitude provider used in propagation
87       * @param propagationType mean or osculating.
88       */
89      DSSTStateTransitionMatrixGenerator(final String stmName,
90                                         final List<DSSTForceModel> forceModels,
91                                         final AttitudeProvider attitudeProvider,
92                                         final PropagationType propagationType) {
93          this.stmName           = stmName;
94          this.forceModels       = forceModels;
95          this.attitudeProvider  = attitudeProvider;
96          this.propagationType   = propagationType;
97          this.partialsObservers = new HashMap<>();
98      }
99  
100     /** Register an observer for partial derivatives.
101      * <p>
102      * The observer {@link DSSTPartialsObserver#partialsComputed(SpacecraftState, RealMatrix, double[])} partialsComputed}
103      * method will be called when partial derivatives are computed, as a side effect of
104      * calling {@link #computePartials(SpacecraftState)} (SpacecraftState)}
105      * </p>
106      * @param name name of the parameter driver this observer is interested in (may be null)
107      * @param observer observer to register
108      */
109     void addObserver(final String name, final DSSTPartialsObserver observer) {
110         partialsObservers.put(name, observer);
111     }
112 
113     /** {@inheritDoc} */
114     @Override
115     public String getName() {
116         return stmName;
117     }
118 
119     /** {@inheritDoc} */
120     @Override
121     public int getDimension() {
122         return STATE_DIMENSION * STATE_DIMENSION;
123     }
124 
125     @Override
126     @SuppressWarnings("unchecked")
127     public void init(final SpacecraftState initialState, final AbsoluteDate target) {
128         // initialize short period terms.
129         // the propagator will have called the non-field method
130         // so just call the field method here
131         // This should be a Field copy of the code in DSSTPropagator.beforeIntegration(...)
132         // but with just the short period initialization calls
133         // See also how the field state is set up in computePartials(...)
134 
135         final DSSTGradientConverter converter =
136                 new DSSTGradientConverter(initialState, attitudeProvider);
137 
138         // check if only mean elements must be used
139         final PropagationType type = propagationType;
140 
141         // initialize all perturbing forces
142         for (final DSSTForceModel forceModel : forceModels) {
143             final FieldSpacecraftState<Gradient> dsState = converter.getState(forceModel);
144             final Gradient[] parameters = converter.getParametersAtStateDate(dsState, forceModel);
145             final FieldAuxiliaryElements<Gradient> auxiliaryElements = new FieldAuxiliaryElements<>(dsState.getOrbit(), I);
146             forceModel.initializeShortPeriodTerms(auxiliaryElements, type, parameters);
147         }
148 
149         // if required, insert the special short periodics step handler
150         if (type == PropagationType.OSCULATING) {
151             // Compute short periodic coefficients for this point
152             for (DSSTForceModel forceModel : forceModels) {
153                 final FieldSpacecraftState<Gradient> dsState = converter.getState(forceModel);
154                 final Gradient[] parameters = converter.getParametersAtStateDate(dsState, forceModel);
155                 forceModel.updateShortPeriodTerms(parameters, dsState);
156             }
157         }
158 
159     }
160 
161     /** {@inheritDoc} */
162     @Override
163     public boolean yields(final SpacecraftState state) {
164         return !state.hasAdditionalData(getName());
165     }
166 
167     /** Set the initial value of the State Transition Matrix.
168      * <p>
169      * The returned state must be added to the propagator.
170      * </p>
171      * @param state initial state
172      * @param dYdY0 initial State Transition Matrix ∂Y/∂Y₀,
173      * if null (which is the most frequent case), assumed to be 6x6 identity
174      * @return state with initial STM (converted to Cartesian ∂C/∂Y₀) added
175      */
176     SpacecraftState setInitialStateTransitionMatrix(final SpacecraftState state, final RealMatrix dYdY0) {
177 
178         if (dYdY0 != null) {
179             if (dYdY0.getRowDimension() != STATE_DIMENSION ||
180                             dYdY0.getColumnDimension() != STATE_DIMENSION) {
181                 throw new OrekitException(LocalizedCoreFormats.DIMENSIONS_MISMATCH_2x2,
182                                           dYdY0.getRowDimension(), dYdY0.getColumnDimension(),
183                                           STATE_DIMENSION, STATE_DIMENSION);
184             }
185         }
186 
187         // flatten matrix
188         final double[] flat = new double[STATE_DIMENSION * STATE_DIMENSION];
189         int k = 0;
190         for (int i = 0; i < STATE_DIMENSION; ++i) {
191             for (int j = 0; j < STATE_DIMENSION; ++j) {
192                 flat[k++] = dYdY0.getEntry(i, j);
193             }
194         }
195 
196         // set additional state
197         return state.addAdditionalData(stmName, flat);
198 
199     }
200 
201     /** {@inheritDoc} */
202     public CombinedDerivatives combinedDerivatives(final SpacecraftState state) {
203 
204         final double[] p = state.getAdditionalState(getName());
205         final double[] res = new double[p.length];
206 
207         // perform matrix multiplication with matrices flatten
208         final RealMatrix factor = computePartials(state);
209         int index = 0;
210         for (int i = 0; i < STATE_DIMENSION; ++i) {
211             for (int j = 0; j < STATE_DIMENSION; ++j) {
212                 double sum = 0;
213                 for (int k = 0; k < STATE_DIMENSION; ++k) {
214                     sum += factor.getEntry(i, k) * p[j + k * STATE_DIMENSION];
215                 }
216                 res[index++] = sum;
217             }
218         }
219 
220         return new CombinedDerivatives(res, null);
221 
222     }
223 
224     /** Compute the various partial derivatives.
225      * @param state current spacecraft state
226      * @return factor matrix
227      */
228     private RealMatrix computePartials(final SpacecraftState state) {
229 
230         // set up containers for partial derivatives
231         final RealMatrix            factor               = MatrixUtils.createRealMatrix(STATE_DIMENSION, STATE_DIMENSION);
232         final DoubleArrayDictionary meanElementsPartials = new DoubleArrayDictionary();
233         final DSSTGradientConverter converter            = new DSSTGradientConverter(state, attitudeProvider);
234 
235         // Compute Jacobian
236         for (final DSSTForceModel forceModel : forceModels) {
237 
238             final FieldSpacecraftState<Gradient> dsState = converter.getState(forceModel);
239             final Gradient[] parameters = converter.getParametersAtStateDate(dsState, forceModel);
240             final FieldAuxiliaryElements<Gradient> auxiliaryElements = new FieldAuxiliaryElements<>(dsState.getOrbit(), I);
241 
242             final Gradient[] meanElementRate = forceModel.getMeanElementRate(dsState, auxiliaryElements, parameters);
243             final double[] derivativesA  = meanElementRate[0].getGradient();
244             final double[] derivativesEx = meanElementRate[1].getGradient();
245             final double[] derivativesEy = meanElementRate[2].getGradient();
246             final double[] derivativesHx = meanElementRate[3].getGradient();
247             final double[] derivativesHy = meanElementRate[4].getGradient();
248             final double[] derivativesL  = meanElementRate[5].getGradient();
249 
250             // update Jacobian with respect to state
251             addToRow(derivativesA,  0, factor);
252             addToRow(derivativesEx, 1, factor);
253             addToRow(derivativesEy, 2, factor);
254             addToRow(derivativesHx, 3, factor);
255             addToRow(derivativesHy, 4, factor);
256             addToRow(derivativesL,  5, factor);
257 
258             // partials derivatives with respect to parameters
259             int paramsIndex = converter.getFreeStateParameters();
260             for (ParameterDriver driver : forceModel.getParametersDrivers()) {
261                 if (driver.isSelected()) {
262                     // for each span (for each estimated value) corresponding name is added
263                     for (Span<String> span = driver.getNamesSpanMap().getFirstSpan(); span != null; span = span.next()) {
264 
265                         // get the partials derivatives for this driver
266                         DoubleArrayDictionary.Entry entry = meanElementsPartials.getEntry(span.getData());
267                         if (entry == null) {
268                             // create an entry filled with zeroes
269                             meanElementsPartials.put(span.getData(), new double[STATE_DIMENSION]);
270                             entry = meanElementsPartials.getEntry(span.getData());
271                         }
272                         // add the contribution of the current force model
273                         entry.increment(new double[] {
274                             derivativesA[paramsIndex], derivativesEx[paramsIndex], derivativesEy[paramsIndex],
275                             derivativesHx[paramsIndex], derivativesHy[paramsIndex], derivativesL[paramsIndex]
276                         });
277                         ++paramsIndex;
278                     }
279 
280                 }
281             }
282 
283         }
284 
285         // notify observers
286         for (Map.Entry<String, DSSTPartialsObserver> observersEntry : partialsObservers.entrySet()) {
287             final DoubleArrayDictionary.Entry entry = meanElementsPartials.getEntry(observersEntry.getKey());
288             observersEntry.getValue().partialsComputed(state, factor, entry == null ? new double[STATE_DIMENSION] : entry.getValue());
289         }
290 
291         return factor;
292 
293     }
294 
295     /** Fill Jacobians rows.
296      * @param derivatives derivatives of a component
297      * @param index component index (0 for a, 1 for ex, 2 for ey, 3 for hx, 4 for hy, 5 for l)
298      * @param factor Jacobian of mean elements rate with respect to mean elements
299      */
300     private void addToRow(final double[] derivatives, final int index, final RealMatrix factor) {
301         for (int i = 0; i < 6; i++) {
302             factor.addToEntry(index, i, derivatives[i]);
303         }
304     }
305 
306     /** Interface for observing partials derivatives. */
307     @FunctionalInterface
308     public interface DSSTPartialsObserver {
309 
310         /** Callback called when partial derivatives have been computed.
311          * @param state current spacecraft state
312          * @param factor factor matrix
313          * @param meanElementsPartials partials derivatives of mean elements rates with respect to the parameter driver
314          * that was registered (zero if no parameters were not selected or parameter is unknown)
315          */
316         void partialsComputed(SpacecraftState state, RealMatrix factor, double[] meanElementsPartials);
317 
318     }
319 
320 }
321