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.linear.RealMatrix;
20 import org.orekit.propagation.SpacecraftState;
21 import org.orekit.propagation.integration.AdditionalDerivativesProvider;
22 import org.orekit.propagation.integration.CombinedDerivatives;
23
24 /** Generator for one column of a Jacobian matrix.
25 * <p>
26 * This generator is based on variational equations, so
27 * it implements {@link AdditionalDerivativesProvider} and
28 * computes only the derivative of the Jacobian column, to
29 * be integrated by the propagator alongside the primary state.
30 * </p>
31 * @author Luc Maisonobe
32 * @since 11.1
33 */
34 class DSSTIntegrableJacobianColumnGenerator
35 implements AdditionalDerivativesProvider, DSSTStateTransitionMatrixGenerator.DSSTPartialsObserver {
36
37 /** Name of the state for State Transition Matrix. */
38 private final String stmName;
39
40 /** Name of the parameter corresponding to the column. */
41 private final String columnName;
42
43 /** Last value computed for the partial derivatives. */
44 private final double[] pDot;
45
46 /** Simple constructor.
47 * <p>
48 * The generator for State Transition Matrix <em>must</em> be registered as
49 * an integrable generator to the same propagator as the instance, as it
50 * must be scheduled to update the state before the instance
51 * </p>
52 * @param stmGenerator generator for State Transition Matrix
53 * @param columnName name of the parameter corresponding to the column
54 */
55 DSSTIntegrableJacobianColumnGenerator(final DSSTStateTransitionMatrixGenerator stmGenerator, final String columnName) {
56 this.stmName = stmGenerator.getName();
57 this.columnName = columnName;
58 this.pDot = new double[getDimension()];
59 stmGenerator.addObserver(columnName, this);
60 }
61
62 /** {@inheritDoc} */
63 @Override
64 public String getName() {
65 return columnName;
66 }
67
68 /** Get the dimension of the generated column.
69 * @return dimension of the generated column
70 */
71 public int getDimension() {
72 return 6;
73 }
74
75 /** {@inheritDoc}
76 * <p>
77 * The column derivative can be computed only if the State Transition Matrix derivatives
78 * are available, as it implies the STM generator has already been run.
79 * </p>
80 */
81 @Override
82 public boolean yields(final SpacecraftState state) {
83 return !state.hasAdditionalStateDerivative(stmName);
84 }
85
86 /** {@inheritDoc} */
87 @Override
88 public void partialsComputed(final SpacecraftState state, final RealMatrix factor, final double[] meanElementsPartials) {
89 // retrieve current Jacobian column
90 final double[] p = state.getAdditionalState(getName());
91
92 // compute time derivative of the Jacobian column
93 System.arraycopy(factor.operate(p), 0, pDot, 0, pDot.length);
94 for (int i = 0; i < pDot.length; ++i) {
95 pDot[i] += meanElementsPartials[i];
96 }
97
98 }
99
100 /** {@inheritDoc} */
101 @Override
102 public CombinedDerivatives combinedDerivatives(final SpacecraftState s) {
103 return new CombinedDerivatives(pDot, null);
104 }
105
106 }
107