ApsideDetectionAdaptableIntervalFactory.java

  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. import org.hipparchus.util.FastMath;
  19. import org.hipparchus.util.MathUtils;
  20. import org.orekit.orbits.KeplerianOrbit;
  21. import org.orekit.orbits.Orbit;
  22. import org.orekit.orbits.OrbitType;

  23. /**
  24.  * Factory class for {@link AdaptableInterval} suitable for apside detection on eccentric orbits.
  25.  * It requires {@link org.orekit.propagation.SpacecraftState} to be based on {@link Orbit} in order to work.
  26.  * @see AdaptableInterval
  27.  * @see org.orekit.propagation.events.ApsideDetector
  28.  * @see org.orekit.propagation.events.EventSlopeFilter
  29.  * @author Romain Serra
  30.  * @since 12.1
  31.  */
  32. public class ApsideDetectionAdaptableIntervalFactory {

  33.     /**
  34.      * Private constructor.
  35.      */
  36.     private ApsideDetectionAdaptableIntervalFactory() {
  37.         // factory class
  38.     }

  39.     /**
  40.      * Method providing a candidate {@link AdaptableInterval} for arbitrary apside detection.
  41.      * It uses a Keplerian, eccentric approximation.
  42.      * @return adaptable interval for apside detection
  43.      */
  44.     public static AdaptableInterval getApsideDetectionAdaptableInterval() {
  45.         return (state, isForward) -> {
  46.             final Orbit orbit = state.getOrbit();
  47.             final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
  48.             final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
  49.             final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
  50.             if (isForward) {
  51.                 final double durationToNextPeriapsis = computeKeplerianDurationToNextPeriapsis(meanAnomaly, meanMotion);
  52.                 final double durationToNextApoapsis = computeKeplerianDurationToNextApoapsis(meanAnomaly, meanMotion);
  53.                 return FastMath.min(durationToNextPeriapsis, durationToNextApoapsis);
  54.             }
  55.             else {
  56.                 final double durationFromPreviousPeriapsis = computeKeplerianDurationFromPreviousPeriapsis(meanAnomaly,
  57.                         meanMotion);
  58.                 final double durationFromPreviousApoapsis = computeKeplerianDurationFromPreviousApoapsis(meanAnomaly,
  59.                         meanMotion);
  60.                 return FastMath.min(durationFromPreviousApoapsis, durationFromPreviousPeriapsis);
  61.             }
  62.         };
  63.     }

  64.     /**
  65.      * Method providing a candidate {@link AdaptableInterval} for periapsis detection.
  66.      * It uses a Keplerian, eccentric approximation.
  67.      * @return adaptable interval for periaspsis detection
  68.      */
  69.     public static AdaptableInterval getPeriapsisDetectionAdaptableInterval() {
  70.         return (state, isForward) -> {
  71.             final Orbit orbit = state.getOrbit();
  72.             final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
  73.             final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
  74.             final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
  75.             if (isForward) {
  76.                 return computeKeplerianDurationToNextPeriapsis(meanAnomaly, meanMotion);
  77.             } else {
  78.                 return computeKeplerianDurationFromPreviousPeriapsis(meanAnomaly, meanMotion);
  79.             }
  80.         };
  81.     }

  82.     /**
  83.      * Method providing a candidate {@link AdaptableInterval} for apoapsis detection.
  84.      * It uses a Keplerian, eccentric approximation.
  85.      * @return adaptable interval for apoapsis detection
  86.      */
  87.     public static AdaptableInterval getApoapsisDetectionAdaptableInterval() {
  88.         return (state, isForward) -> {
  89.             final Orbit orbit = state.getOrbit();
  90.             final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
  91.             final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
  92.             final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
  93.             if (isForward) {
  94.                 return computeKeplerianDurationToNextApoapsis(meanAnomaly, meanMotion);
  95.             }
  96.             else {
  97.                 return computeKeplerianDurationFromPreviousApoapsis(meanAnomaly, meanMotion);
  98.             }
  99.         };
  100.     }

  101.     /**
  102.      * Convert a generic {@link Orbit} into a {@link KeplerianOrbit}.
  103.      * @param orbit orbit to convert
  104.      * @return Keplerian orbit
  105.      */
  106.     private static KeplerianOrbit convertOrbitIntoKeplerianOne(final Orbit orbit) {
  107.         return (KeplerianOrbit) OrbitType.KEPLERIAN.convertType(orbit);
  108.     }

  109.     /**
  110.      * Method computing time to go until next periapsis, assuming Keplerian motion.
  111.      * @param meanAnomaly mean anomaly
  112.      * @param meanMotion Keplerian mean motion
  113.      * @return duration to next periapsis
  114.      */
  115.     private static double computeKeplerianDurationToNextPeriapsis(final double meanAnomaly,
  116.                                                                   final double meanMotion) {
  117.         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, FastMath.PI);
  118.         return (MathUtils.TWO_PI - normalizedMeanAnomaly) / meanMotion;
  119.     }

  120.     /**
  121.      * Method computing time elapsed since last periapsis, assuming Keplerian motion.
  122.      * @param meanAnomaly mean anomaly
  123.      * @param meanMotion Keplerian mean motion
  124.      * @return duration elapsed since last periapsis
  125.      */
  126.     public static double computeKeplerianDurationFromPreviousPeriapsis(final double meanAnomaly,
  127.                                                                        final double meanMotion) {
  128.         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, FastMath.PI);
  129.         return normalizedMeanAnomaly / meanMotion;
  130.     }

  131.     /**
  132.      * Method computing time to go until next apoapsis, assuming Keplerian motion.
  133.      * @param meanAnomaly mean anomaly
  134.      * @param meanMotion Keplerian mean motion
  135.      * @return duration to next apoapsis
  136.      */
  137.     private static double computeKeplerianDurationToNextApoapsis(final double meanAnomaly,
  138.                                                                  final double meanMotion) {
  139.         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, MathUtils.TWO_PI);
  140.         return (MathUtils.TWO_PI + FastMath.PI - normalizedMeanAnomaly) / meanMotion;
  141.     }

  142.     /**
  143.      * Method computing time elapsed since last apoapsis, assuming Keplerian motion.
  144.      * @param meanAnomaly mean anomaly
  145.      * @param meanMotion Keplerian mean motion
  146.      * @return duration elapsed since last apoapsis
  147.      */
  148.     public static double computeKeplerianDurationFromPreviousApoapsis(final double meanAnomaly,
  149.                                                                       final double meanMotion) {
  150.         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, MathUtils.TWO_PI);
  151.         return (normalizedMeanAnomaly - FastMath.PI) / meanMotion;
  152.     }
  153. }