1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17 package org.orekit.propagation.numerical.cr3bp;
18
19 import java.util.Arrays;
20
21 import org.hipparchus.analysis.differentiation.DerivativeStructure;
22 import org.hipparchus.linear.Array2DRowRealMatrix;
23 import org.hipparchus.linear.RealMatrix;
24 import org.orekit.bodies.CR3BPSystem;
25 import org.orekit.propagation.SpacecraftState;
26 import org.orekit.propagation.integration.AdditionalDerivativesProvider;
27 import org.orekit.time.AbsoluteDate;
28
29
30
31
32
33
34 @SuppressWarnings("deprecation")
35 public class STMEquations
36 implements AdditionalDerivativesProvider,
37 org.orekit.propagation.integration.AdditionalEquations {
38
39
40 private static final int DIM = 6;
41
42
43 private final CR3BPSystem syst;
44
45
46 private final String name;
47
48
49 private final double[][] jacobian = new double[DIM][DIM];
50
51
52
53
54 public STMEquations(final CR3BPSystem syst) {
55 this.syst = syst;
56 this.name = "stmEquations";
57
58
59 for (int j = 0; j < jacobian.length; ++j) {
60 Arrays.fill(jacobian[j], 0.0);
61 }
62
63 jacobian[0][3] = 1.0;
64 jacobian[1][4] = 1.0;
65 jacobian[2][5] = 1.0;
66 jacobian[3][4] = 2.0;
67 jacobian[4][3] = -2.0;
68 }
69
70
71
72
73
74 public SpacecraftState setInitialPhi(final SpacecraftState s) {
75 final int stateDimension = 36;
76 final double[] phi = new double[stateDimension];
77 for (int i = 0; i < stateDimension; i = i + 7) {
78 phi[i] = 1.0;
79 }
80 return s.addAdditionalState(name, phi);
81 }
82
83
84 public void init(final SpacecraftState initialState, final AbsoluteDate target) {
85
86 AdditionalDerivativesProvider.super.init(initialState, target);
87 }
88
89
90 public double[] computeDerivatives(final SpacecraftState s, final double[] pDot) {
91
92 System.arraycopy(derivatives(s), 0, pDot, 0, pDot.length);
93 return null;
94 }
95
96
97 public double[] derivatives(final SpacecraftState s) {
98
99
100 final double[] phi = s.getAdditionalState(getName());
101 final double[] dPhi = new double[phi.length];
102
103
104 final DerivativeStructure potential = new CR3BPForceModel(syst).getPotential(s);
105
106
107 final double[] dU = potential.getAllDerivatives();
108
109
110 final int idXX = potential.getFactory().getCompiler().getPartialDerivativeIndex(2, 0, 0);
111 final int idXY = potential.getFactory().getCompiler().getPartialDerivativeIndex(1, 1, 0);
112 final int idXZ = potential.getFactory().getCompiler().getPartialDerivativeIndex(1, 0, 1);
113 final int idYY = potential.getFactory().getCompiler().getPartialDerivativeIndex(0, 2, 0);
114 final int idYZ = potential.getFactory().getCompiler().getPartialDerivativeIndex(0, 1, 1);
115 final int idZZ = potential.getFactory().getCompiler().getPartialDerivativeIndex(0, 0, 2);
116
117
118 jacobian[3][0] = dU[idXX];
119 jacobian[4][1] = dU[idYY];
120 jacobian[5][2] = dU[idZZ];
121 jacobian[3][1] = dU[idXY];
122 jacobian[4][0] = jacobian[3][1];
123 jacobian[3][2] = dU[idXZ];
124 jacobian[5][0] = jacobian[3][2];
125 jacobian[4][2] = dU[idYZ];
126 jacobian[5][1] = jacobian[4][2];
127
128
129 for (int k = 0; k < DIM; k++) {
130 for (int l = 0; l < DIM; l++) {
131 for (int i = 0; i < DIM; i++) {
132 dPhi[DIM * k + l] =
133 dPhi[DIM * k + l] + jacobian[k][i] * phi[DIM * i + l];
134 }
135 }
136 }
137
138 return dPhi;
139
140 }
141
142
143 public String getName() {
144 return name;
145 }
146
147
148 @Override
149 public int getDimension() {
150 return DIM * DIM;
151 }
152
153
154
155
156
157 public RealMatrix getStateTransitionMatrix(final SpacecraftState s) {
158 final double[][] phi2dA = new double[DIM][DIM];
159 final double[] stm = s.getAdditionalState(getName());
160 for (int i = 0; i < DIM; i++) {
161 for (int j = 0; j < 6; j++) {
162 phi2dA[i][j] = stm[DIM * i + j];
163 }
164 }
165 return new Array2DRowRealMatrix(phi2dA, false);
166 }
167 }