1   /* Copyright 2022-2025 Romain Serra
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.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  
25  /**
26   * Factory class for {@link AdaptableInterval} suitable for apside detection on eccentric orbits.
27   * It requires {@link org.orekit.propagation.SpacecraftState} to be based on {@link Orbit} in order to work.
28   * @see AdaptableInterval
29   * @see org.orekit.propagation.events.ApsideDetector
30   * @see org.orekit.propagation.events.EventSlopeFilter
31   * @author Romain Serra
32   * @since 12.1
33   */
34  public class ApsideDetectionAdaptableIntervalFactory {
35  
36      /**
37       * Private constructor.
38       */
39      private ApsideDetectionAdaptableIntervalFactory() {
40          // factory class
41      }
42  
43      /**
44       * Method providing a candidate {@link AdaptableInterval} for arbitrary apside detection.
45       * It uses a Keplerian, eccentric approximation.
46       * @return adaptable interval for apside detection
47       */
48      public static AdaptableInterval getApsideDetectionAdaptableInterval() {
49          return (state, isForward) -> {
50              final Orbit orbit = state.getOrbit();
51              final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
52              final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
53              final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
54              if (isForward) {
55                  final double durationToNextPeriapsis = computeKeplerianDurationToNextPeriapsis(meanAnomaly, meanMotion);
56                  final double durationToNextApoapsis = computeKeplerianDurationToNextApoapsis(meanAnomaly, meanMotion);
57                  return FastMath.min(durationToNextPeriapsis, durationToNextApoapsis);
58              }
59              else {
60                  final double durationFromPreviousPeriapsis = computeKeplerianDurationFromPreviousPeriapsis(meanAnomaly,
61                          meanMotion);
62                  final double durationFromPreviousApoapsis = computeKeplerianDurationFromPreviousApoapsis(meanAnomaly,
63                          meanMotion);
64                  return FastMath.min(durationFromPreviousApoapsis, durationFromPreviousPeriapsis);
65              }
66          };
67      }
68  
69      /**
70       * Method providing a candidate {@link AdaptableInterval} for periapsis detection.
71       * It uses a Keplerian, eccentric approximation.
72       * @return adaptable interval for periaspsis detection
73       */
74      public static AdaptableInterval getPeriapsisDetectionAdaptableInterval() {
75          return (state, isForward) -> {
76              final Orbit orbit = state.getOrbit();
77              final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
78              final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
79              final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
80              if (isForward) {
81                  return computeKeplerianDurationToNextPeriapsis(meanAnomaly, meanMotion);
82              } else {
83                  return computeKeplerianDurationFromPreviousPeriapsis(meanAnomaly, meanMotion);
84              }
85          };
86      }
87  
88      /**
89       * Method providing a candidate {@link AdaptableInterval} for apoapsis detection.
90       * It uses a Keplerian, eccentric approximation.
91       * @return adaptable interval for apoapsis detection
92       */
93      public static AdaptableInterval getApoapsisDetectionAdaptableInterval() {
94          return (state, isForward) -> {
95              final Orbit orbit = state.getOrbit();
96              final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
97              final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
98              final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
99              if (isForward) {
100                 return computeKeplerianDurationToNextApoapsis(meanAnomaly, meanMotion);
101             }
102             else {
103                 return computeKeplerianDurationFromPreviousApoapsis(meanAnomaly, meanMotion);
104             }
105         };
106     }
107 
108     /**
109      * Convert a generic {@link Orbit} into a {@link KeplerianOrbit}.
110      * @param orbit orbit to convert
111      * @return Keplerian orbit
112      */
113     private static KeplerianOrbit convertOrbitIntoKeplerianOne(final Orbit orbit) {
114         return (KeplerianOrbit) OrbitType.KEPLERIAN.convertType(orbit);
115     }
116 
117     /**
118      * Method computing time to go until next periapsis, assuming Keplerian motion.
119      * @param meanAnomaly mean anomaly
120      * @param meanMotion Keplerian mean motion
121      * @return duration to next periapsis
122      */
123     private static double computeKeplerianDurationToNextPeriapsis(final double meanAnomaly,
124                                                                   final double meanMotion) {
125         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, FastMath.PI);
126         return (MathUtils.TWO_PI - normalizedMeanAnomaly) / meanMotion;
127     }
128 
129     /**
130      * Method computing time elapsed since last periapsis, assuming Keplerian motion.
131      * @param meanAnomaly mean anomaly
132      * @param meanMotion Keplerian mean motion
133      * @return duration elapsed since last periapsis
134      */
135     public static double computeKeplerianDurationFromPreviousPeriapsis(final double meanAnomaly,
136                                                                        final double meanMotion) {
137         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, FastMath.PI);
138         return normalizedMeanAnomaly / meanMotion;
139     }
140 
141     /**
142      * Method computing time to go until next apoapsis, assuming Keplerian motion.
143      * @param meanAnomaly mean anomaly
144      * @param meanMotion Keplerian mean motion
145      * @return duration to next apoapsis
146      */
147     private static double computeKeplerianDurationToNextApoapsis(final double meanAnomaly,
148                                                                  final double meanMotion) {
149         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, MathUtils.TWO_PI);
150         return (MathUtils.TWO_PI + FastMath.PI - normalizedMeanAnomaly) / meanMotion;
151     }
152 
153     /**
154      * Method computing time elapsed since last apoapsis, assuming Keplerian motion.
155      * @param meanAnomaly mean anomaly
156      * @param meanMotion Keplerian mean motion
157      * @return duration elapsed since last apoapsis
158      */
159     public static double computeKeplerianDurationFromPreviousApoapsis(final double meanAnomaly,
160                                                                       final double meanMotion) {
161         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, MathUtils.TWO_PI);
162         return (normalizedMeanAnomaly - FastMath.PI) / meanMotion;
163     }
164 }