ApsideDetectionAdaptableIntervalFactory.java

  1. /* Copyright 2022-2024 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. import org.orekit.propagation.events.AdaptableInterval;

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

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

  40.     /**
  41.      * Method providing a candidate {@link AdaptableInterval} for arbitrary apside detection with forward propagation.
  42.      * It uses a Keplerian, eccentric approximation.
  43.      * @return adaptable interval for forward apside detection
  44.      */
  45.     public static AdaptableInterval getForwardApsideDetectionAdaptableInterval() {
  46.         return state -> {
  47.             final Orbit orbit = state.getOrbit();
  48.             final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
  49.             final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
  50.             final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
  51.             final double durationToNextPeriapsis = computeKeplerianDurationToNextPeriapsis(meanAnomaly, meanMotion);
  52.             final double durationToNextApoapsis = computeKeplerianDurationToNextApoapsis(meanAnomaly, meanMotion);
  53.             return FastMath.min(durationToNextPeriapsis, durationToNextApoapsis);
  54.         };
  55.     }

  56.     /**
  57.      * Method providing a candidate {@link AdaptableInterval} for arbitrary apside detection with backward propagation.
  58.      * It uses a Keplerian, eccentric approximation.
  59.      * @return adaptable interval for backward apside detection
  60.      */
  61.     public static AdaptableInterval getBackwardApsideDetectionAdaptableInterval() {
  62.         return state -> {
  63.             final Orbit orbit = state.getOrbit();
  64.             final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
  65.             final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
  66.             final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
  67.             final double durationFromPreviousPeriapsis = computeKeplerianDurationFromPreviousPeriapsis(meanAnomaly,
  68.                     meanMotion);
  69.             final double durationFromPreviousApoapsis = computeKeplerianDurationFromPreviousApoapsis(meanAnomaly,
  70.                     meanMotion);
  71.             return FastMath.min(durationFromPreviousApoapsis, durationFromPreviousPeriapsis);
  72.         };
  73.     }

  74.     /**
  75.      * Method providing a candidate {@link AdaptableInterval} for periapsis detection with forward propagation.
  76.      * It uses a Keplerian, eccentric approximation.
  77.      * @return adaptable interval for forward periaspsis detection
  78.      */
  79.     public static AdaptableInterval getForwardPeriapsisDetectionAdaptableInterval() {
  80.         return state -> {
  81.             final Orbit orbit = state.getOrbit();
  82.             final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
  83.             final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
  84.             final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
  85.             return computeKeplerianDurationToNextPeriapsis(meanAnomaly, meanMotion);
  86.         };
  87.     }

  88.     /**
  89.      * Method providing a candidate {@link AdaptableInterval} for periapsis detection with backward propagation.
  90.      * It uses a Keplerian, eccentric approximation.
  91.      * @return adaptable interval for backward periaspsis detection
  92.      */
  93.     public static AdaptableInterval getBackwardPeriapsisDetectionAdaptableInterval() {
  94.         return state -> {
  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.             return computeKeplerianDurationFromPreviousPeriapsis(meanAnomaly, meanMotion);
  100.         };
  101.     }

  102.     /**
  103.      * Method providing a candidate {@link AdaptableInterval} for apoapsis detection with forward propagation.
  104.      * It uses a Keplerian, eccentric approximation.
  105.      * @return adaptable interval for forward apoapsis detection
  106.      */
  107.     public static AdaptableInterval getForwardApoapsisDetectionAdaptableInterval() {
  108.         return state -> {
  109.             final Orbit orbit = state.getOrbit();
  110.             final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
  111.             final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
  112.             final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
  113.             return computeKeplerianDurationToNextApoapsis(meanAnomaly, meanMotion);
  114.         };
  115.     }

  116.     /**
  117.      * Method providing a candidate {@link AdaptableInterval} for apoapsis detection with backward propagation.
  118.      * It uses a Keplerian, eccentric approximation.
  119.      * @return adaptable interval for backward apoapsis detection
  120.      */
  121.     public static AdaptableInterval getBackwardApoapsisDetectionAdaptableInterval() {
  122.         return state -> {
  123.             final Orbit orbit = state.getOrbit();
  124.             final KeplerianOrbit keplerianOrbit = convertOrbitIntoKeplerianOne(orbit);
  125.             final double meanMotion = keplerianOrbit.getKeplerianMeanMotion();
  126.             final double meanAnomaly = keplerianOrbit.getMeanAnomaly();
  127.             return computeKeplerianDurationFromPreviousApoapsis(meanAnomaly, meanMotion);
  128.         };
  129.     }

  130.     /**
  131.      * Convert a generic {@link Orbit} into a {@link KeplerianOrbit}.
  132.      * @param orbit orbit to convert
  133.      * @return Keplerian orbit
  134.      */
  135.     private static KeplerianOrbit convertOrbitIntoKeplerianOne(final Orbit orbit) {
  136.         return (KeplerianOrbit) OrbitType.KEPLERIAN.convertType(orbit);
  137.     }

  138.     /**
  139.      * Method computing time to go until next periapsis, assuming Keplerian motion.
  140.      * @param meanAnomaly mean anomaly
  141.      * @param meanMotion Keplerian mean motion
  142.      * @return duration to next periapsis
  143.      */
  144.     private static double computeKeplerianDurationToNextPeriapsis(final double meanAnomaly,
  145.                                                                   final double meanMotion) {
  146.         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, FastMath.PI);
  147.         return (MathUtils.TWO_PI - normalizedMeanAnomaly) / meanMotion;
  148.     }

  149.     /**
  150.      * Method computing time elapsed since last periapsis, assuming Keplerian motion.
  151.      * @param meanAnomaly mean anomaly
  152.      * @param meanMotion Keplerian mean motion
  153.      * @return duration elapsed since last periapsis
  154.      */
  155.     public static double computeKeplerianDurationFromPreviousPeriapsis(final double meanAnomaly,
  156.                                                                        final double meanMotion) {
  157.         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, FastMath.PI);
  158.         return normalizedMeanAnomaly / meanMotion;
  159.     }

  160.     /**
  161.      * Method computing time to go until next apoapsis, assuming Keplerian motion.
  162.      * @param meanAnomaly mean anomaly
  163.      * @param meanMotion Keplerian mean motion
  164.      * @return duration to next apoapsis
  165.      */
  166.     private static double computeKeplerianDurationToNextApoapsis(final double meanAnomaly,
  167.                                                                  final double meanMotion) {
  168.         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, MathUtils.TWO_PI);
  169.         return (MathUtils.TWO_PI + FastMath.PI - normalizedMeanAnomaly) / meanMotion;
  170.     }

  171.     /**
  172.      * Method computing time elapsed since last apoapsis, assuming Keplerian motion.
  173.      * @param meanAnomaly mean anomaly
  174.      * @param meanMotion Keplerian mean motion
  175.      * @return duration elapsed since last apoapsis
  176.      */
  177.     public static double computeKeplerianDurationFromPreviousApoapsis(final double meanAnomaly,
  178.                                                                       final double meanMotion) {
  179.         final double normalizedMeanAnomaly = MathUtils.normalizeAngle(meanAnomaly, MathUtils.TWO_PI);
  180.         return (normalizedMeanAnomaly - FastMath.PI) / meanMotion;
  181.     }
  182. }