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