SGP4OrbitalState.java

  1. /* Copyright 2020-2025 Exotrail
  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.conversion.averaging;

  18. import org.hipparchus.util.FastMath;
  19. import org.orekit.annotation.DefaultDataContext;
  20. import org.orekit.data.DataContext;
  21. import org.orekit.frames.Frame;
  22. import org.orekit.orbits.Orbit;
  23. import org.orekit.orbits.OrbitType;
  24. import org.orekit.orbits.PositionAngleType;
  25. import org.orekit.propagation.analytical.tle.TLE;
  26. import org.orekit.propagation.analytical.tle.TLEPropagator;
  27. import org.orekit.propagation.conversion.averaging.elements.AveragedKeplerianWithMeanAngle;
  28. import org.orekit.time.AbsoluteDate;
  29. import org.orekit.time.TimeScale;
  30. import org.orekit.time.UTCScale;

  31. /**
  32.  * Class representing an averaged orbital state as in the TLE-related theory.
  33.  * Note it is the averaged mean motion that is written in a Two-Line Element and that, for now,
  34.  * conversions back and forth to averaged semi-major axis are approximated with the osculating ones.
  35.  *
  36.  * @author Romain Serra
  37.  * @see AveragedOrbitalState
  38.  * @see TLEPropagator
  39.  * @since 12.1
  40.  */
  41. public class SGP4OrbitalState extends AbstractAveragedOrbitalState {

  42.     /** B-star used internally and set to zero. Should not impact calculations. */
  43.     private static final double B_STAR = 0.;

  44.     /** Averaged Keplerian elements. */
  45.     private final AveragedKeplerianWithMeanAngle averagedElements;
  46.     /** UTC time scale. */
  47.     private final UTCScale utc;

  48.     /**
  49.      * Constructor.
  50.      * @param date epoch
  51.      * @param elements averaged orbital elements
  52.      * @param dataContext data context
  53.      */
  54.     public SGP4OrbitalState(final AbsoluteDate date,
  55.                             final AveragedKeplerianWithMeanAngle elements,
  56.                             final DataContext dataContext) {
  57.         this(date, elements, dataContext.getFrames().getTEME(),
  58.                 dataContext.getTimeScales().getUTC());
  59.     }

  60.     /**
  61.      * Constructor with default data context.
  62.      * @param date epoch
  63.      * @param elements averaged orbital elements
  64.      */
  65.     @DefaultDataContext
  66.     public SGP4OrbitalState(final AbsoluteDate date,
  67.                             final AveragedKeplerianWithMeanAngle elements) {
  68.         this(date, elements, DataContext.getDefault());
  69.     }

  70.     /**
  71.      * Private constructor.
  72.      * @param date epoch
  73.      * @param elements averaged orbital elements
  74.      * @param teme TEME frame
  75.      * @param utc UTC time scale
  76.      */
  77.     private SGP4OrbitalState(final AbsoluteDate date,
  78.                              final AveragedKeplerianWithMeanAngle elements,
  79.                              final Frame teme, final TimeScale utc) {
  80.         super(date, teme);
  81.         this.averagedElements = elements;
  82.         this.utc = (UTCScale) utc;
  83.     }

  84.     /**
  85.      * Static constructor. Input frame is implicitly assumed to be TEME (it is not checked).
  86.      * @param tle TLE
  87.      * @param teme TEME frame (not checked)
  88.      * @return TLE-based averaged orbital state
  89.      */
  90.     public static SGP4OrbitalState of(final TLE tle, final Frame teme) {
  91.         final double semiMajorAxis = computeSemiMajorAxis(tle);
  92.         final AveragedKeplerianWithMeanAngle elements = new AveragedKeplerianWithMeanAngle(
  93.                 semiMajorAxis, tle.getE(), tle.getI(), tle.getPerigeeArgument(), tle.getRaan(),
  94.                 tle.getMeanAnomaly());
  95.         return new SGP4OrbitalState(tle.getDate(), elements, teme, tle.getUtc());
  96.     }

  97.     /** {@inheritDoc} */
  98.     @Override
  99.     public double getMu() {
  100.         return getTleMu();
  101.     }

  102.     /**
  103.      * Getter for TLE's Earth gravitational constant.
  104.      * @return mu.
  105.      */
  106.     private static double getTleMu() {
  107.         return TLEPropagator.getMU();
  108.     }

  109.     /** {@inheritDoc} */
  110.     @Override
  111.     public OrbitType getOrbitType() {
  112.         return OrbitType.KEPLERIAN;
  113.     }

  114.     /** {@inheritDoc} */
  115.     @Override
  116.     public PositionAngleType getPositionAngleType() {
  117.         return PositionAngleType.MEAN;
  118.     }

  119.     /** {@inheritDoc} */
  120.     @Override
  121.     public AveragedKeplerianWithMeanAngle getAveragedElements() {
  122.         return averagedElements;
  123.     }

  124.     /** {@inheritDoc} */
  125.     @Override
  126.     public Orbit toOsculatingOrbit() {
  127.         final TLEPropagator propagator = createPropagator();
  128.         return propagator.getInitialState().getOrbit();
  129.     }

  130.     /**
  131.      * Create TLE propagator.
  132.      * @return propagator using relevant theory
  133.      */
  134.     private TLEPropagator createPropagator() {
  135.         final TLE tle = createTLE();
  136.         return TLEPropagator.selectExtrapolator(tle, getFrame());
  137.     }

  138.     /**
  139.      * Create Orekit TLE.
  140.      * @return TLE
  141.      */
  142.     private TLE createTLE() {
  143.         final double averagedMeanMotion = computeMeanMotion(getAveragedElements()
  144.                 .getAveragedSemiMajorAxis());
  145.         final double averagedEccentricity = getAveragedElements().getAveragedEccentricity();
  146.         final double averagedInclination = getAveragedElements().getAveragedInclination();
  147.         final double averagedRAAN = getAveragedElements().getAveragedRightAscensionOfTheAscendingNode();
  148.         final double averagedPerigeeArgument = getAveragedElements().getAveragedPerigeeArgument();
  149.         final double averagedMeanAnomaly = getAveragedElements().getAveragedMeanAnomaly();
  150.         return new TLE(0, (char) 0, 2000, 1, "1", 0, 0,
  151.                 getDate(), averagedMeanMotion, 0., 0., averagedEccentricity,
  152.                 averagedInclination, averagedPerigeeArgument, averagedRAAN, averagedMeanAnomaly,
  153.                 1, B_STAR, utc);
  154.     }

  155.     /**
  156.      * Convert averaged semi-major axis to averaged mean motion. Uses an approximate transformation
  157.      * (same as osculating).
  158.      * @param averagedSemiMajorAxis semi-major axis
  159.      * @return mean motion
  160.      */
  161.     private static double computeMeanMotion(final double averagedSemiMajorAxis) {
  162.         final double cubedSemiMajorAxis = averagedSemiMajorAxis * averagedSemiMajorAxis *
  163.                 averagedSemiMajorAxis;
  164.         return FastMath.sqrt(getTleMu() / cubedSemiMajorAxis);
  165.     }

  166.     /**
  167.      * Compute semi-major axis from Two-Line Elements. Uses an approximate transformation
  168.      * (same as osculating).
  169.      * @param tle TLE
  170.      * @return semi-major axis
  171.      */
  172.     private static double computeSemiMajorAxis(final TLE tle) {
  173.         final double meanMotion = tle.getMeanMotion();
  174.         return FastMath.cbrt(getTleMu() / (meanMotion * meanMotion));
  175.     }
  176. }