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