1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17 package org.orekit.propagation.events.intervals;
18
19 import org.hipparchus.util.FastMath;
20 import org.hipparchus.util.MathUtils;
21 import org.orekit.orbits.KeplerianOrbit;
22 import org.orekit.orbits.Orbit;
23 import org.orekit.orbits.OrbitType;
24 import org.orekit.propagation.events.AdaptableInterval;
25
26
27
28
29
30
31
32
33
34
35 public class ApsideDetectionAdaptableIntervalFactory {
36
37
38
39
40 private ApsideDetectionAdaptableIntervalFactory() {
41
42 }
43
44
45
46
47
48
49 public static AdaptableInterval getForwardApsideDetectionAdaptableInterval() {
50 return state -> {
51 final Orbit orbit = state.getOrbit();
52 final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
53 final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
54 final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
55 final double durationToNextPeriapsis = computeKeplerianDurationToNextPeriapsis(meanAnomaly, meanMotion);
56 final double durationToNextApoapsis = computeKeplerianDurationToNextApoapsis(meanAnomaly, meanMotion);
57 return FastMath.min(durationToNextPeriapsis, durationToNextApoapsis);
58 };
59 }
60
61
62
63
64
65
66 public static AdaptableInterval getBackwardApsideDetectionAdaptableInterval() {
67 return state -> {
68 final Orbit orbit = state.getOrbit();
69 final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
70 final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
71 final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
72 final double durationFromPreviousPeriapsis = computeKeplerianDurationFromPreviousPeriapsis(meanAnomaly,
73 meanMotion);
74 final double durationFromPreviousApoapsis = computeKeplerianDurationFromPreviousApoapsis(meanAnomaly,
75 meanMotion);
76 return FastMath.min(durationFromPreviousApoapsis, durationFromPreviousPeriapsis);
77 };
78 }
79
80
81
82
83
84
85 public static AdaptableInterval getForwardPeriapsisDetectionAdaptableInterval() {
86 return state -> {
87 final Orbit orbit = state.getOrbit();
88 final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
89 final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
90 final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
91 return computeKeplerianDurationToNextPeriapsis(meanAnomaly, meanMotion);
92 };
93 }
94
95
96
97
98
99
100 public static AdaptableInterval getBackwardPeriapsisDetectionAdaptableInterval() {
101 return state -> {
102 final Orbit orbit = state.getOrbit();
103 final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
104 final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
105 final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
106 return computeKeplerianDurationFromPreviousPeriapsis(meanAnomaly, meanMotion);
107 };
108 }
109
110
111
112
113
114
115 public static AdaptableInterval getForwardApoapsisDetectionAdaptableInterval() {
116 return state -> {
117 final Orbit orbit = state.getOrbit();
118 final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
119 final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
120 final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
121 return computeKeplerianDurationToNextApoapsis(meanAnomaly, meanMotion);
122 };
123 }
124
125
126
127
128
129
130 public static AdaptableInterval getBackwardApoapsisDetectionAdaptableInterval() {
131 return state -> {
132 final Orbit orbit = state.getOrbit();
133 final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
134 final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
135 final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
136 return computeKeplerianDurationFromPreviousApoapsis(meanAnomaly, meanMotion);
137 };
138 }
139
140
141
142
143
144
145 private static KeplerianOrbit convertOrbitIntoKeplerianOne(final Orbit orbit) {
146 return (KeplerianOrbit) OrbitType.KEPLERIAN.convertType(orbit);
147 }
148
149
150
151
152
153
154
155 private static double computeKeplerianDurationToNextPeriapsis(final double meanAnomaly,
156 final double meanMotion) {
157 final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, FastMath.PI);
158 return (MathUtils.TWO_PI - normalizedMeanAnomaly) / meanMotion;
159 }
160
161
162
163
164
165
166
167 public static double computeKeplerianDurationFromPreviousPeriapsis(final double meanAnomaly,
168 final double meanMotion) {
169 final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, FastMath.PI);
170 return normalizedMeanAnomaly / meanMotion;
171 }
172
173
174
175
176
177
178
179 private static double computeKeplerianDurationToNextApoapsis(final double meanAnomaly,
180 final double meanMotion) {
181 final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, MathUtils.TWO_PI);
182 return (MathUtils.TWO_PI + FastMath.PI - normalizedMeanAnomaly) / meanMotion;
183 }
184
185
186
187
188
189
190
191 public static double computeKeplerianDurationFromPreviousApoapsis(final double meanAnomaly,
192 final double meanMotion) {
193 final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, MathUtils.TWO_PI);
194 return (normalizedMeanAnomaly - FastMath.PI) / meanMotion;
195 }
196 }