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  
19  import org.hipparchus.util.FastMath;
20  import org.orekit.annotation.DefaultDataContext;
21  import org.orekit.data.DataContext;
22  import org.orekit.frames.Frame;
23  import org.orekit.orbits.Orbit;
24  import org.orekit.orbits.OrbitType;
25  import org.orekit.orbits.PositionAngleType;
26  import org.orekit.propagation.analytical.tle.TLE;
27  import org.orekit.propagation.analytical.tle.TLEPropagator;
28  import org.orekit.propagation.conversion.averaging.elements.AveragedKeplerianWithMeanAngle;
29  import org.orekit.time.AbsoluteDate;
30  import org.orekit.time.TimeScale;
31  import org.orekit.time.UTCScale;
32  
33  /**
34   * Class representing an averaged orbital state as in the TLE-related theory.
35   * Note it is the averaged mean motion that is written in a Two-Line Element and that, for now,
36   * conversions back and forth to averaged semi-major axis are approximated with the osculating ones.
37   *
38   * @author Romain Serra
39   * @see AveragedOrbitalState
40   * @see TLEPropagator
41   * @since 12.1
42   */
43  public class SGP4OrbitalState extends AbstractAveragedOrbitalState {
44  
45      /** B-star used internally and set to zero. Should not impact calculations. */
46      private static final double B_STAR = 0.;
47  
48      /** Averaged Keplerian elements. */
49      private final AveragedKeplerianWithMeanAngle averagedElements;
50      /** UTC time scale. */
51      private final UTCScale utc;
52  
53      /**
54       * Constructor.
55       * @param date epoch
56       * @param elements averaged orbital elements
57       * @param dataContext data context
58       */
59      public SGP4OrbitalState(final AbsoluteDate date,
60                              final AveragedKeplerianWithMeanAngle elements,
61                              final DataContext dataContext) {
62          this(date, elements, dataContext.getFrames().getTEME(),
63                  dataContext.getTimeScales().getUTC());
64      }
65  
66      /**
67       * Constructor with default data context.
68       * @param date epoch
69       * @param elements averaged orbital elements
70       */
71      @DefaultDataContext
72      public SGP4OrbitalState(final AbsoluteDate date,
73                              final AveragedKeplerianWithMeanAngle elements) {
74          this(date, elements, DataContext.getDefault());
75      }
76  
77      /**
78       * Private constructor.
79       * @param date epoch
80       * @param elements averaged orbital elements
81       * @param teme TEME frame
82       * @param utc UTC time scale
83       */
84      private SGP4OrbitalState(final AbsoluteDate date,
85                               final AveragedKeplerianWithMeanAngle elements,
86                               final Frame teme, final TimeScale utc) {
87          super(date, teme);
88          this.averagedElements = elements;
89          this.utc = (UTCScale) utc;
90      }
91  
92      /**
93       * Static constructor. Input frame is implicitly assumed to be TEME (it is not checked).
94       * @param tle TLE
95       * @param teme TEME frame (not checked)
96       * @return TLE-based averaged orbital state
97       */
98      public static SGP4OrbitalState of(final TLE tle, final Frame teme) {
99          final double semiMajorAxis = computeSemiMajorAxis(tle);
100         final AveragedKeplerianWithMeanAngle elements = new AveragedKeplerianWithMeanAngle(
101                 semiMajorAxis, tle.getE(), tle.getI(), tle.getPerigeeArgument(), tle.getRaan(),
102                 tle.getMeanAnomaly());
103         return new SGP4OrbitalState(tle.getDate(), elements, teme, tle.getUtc());
104     }
105 
106     /** {@inheritDoc} */
107     @Override
108     public double getMu() {
109         return getTleMu();
110     }
111 
112     /**
113      * Getter for TLE's Earth gravitational constant.
114      * @return mu.
115      */
116     private static double getTleMu() {
117         return TLEPropagator.getMU();
118     }
119 
120     /** {@inheritDoc} */
121     @Override
122     public OrbitType getOrbitType() {
123         return OrbitType.KEPLERIAN;
124     }
125 
126     /** {@inheritDoc} */
127     @Override
128     public PositionAngleType getPositionAngleType() {
129         return PositionAngleType.MEAN;
130     }
131 
132     /** {@inheritDoc} */
133     @Override
134     public AveragedKeplerianWithMeanAngle getAveragedElements() {
135         return averagedElements;
136     }
137 
138     /** {@inheritDoc} */
139     @Override
140     public Orbit toOsculatingOrbit() {
141         final TLEPropagator propagator = createPropagator();
142         return propagator.getInitialState().getOrbit();
143     }
144 
145     /**
146      * Create TLE propagator.
147      * @return propagator using relevant theory
148      */
149     private TLEPropagator createPropagator() {
150         final TLE tle = createTLE();
151         return TLEPropagator.selectExtrapolator(tle, getFrame());
152     }
153 
154     /**
155      * Create Orekit TLE.
156      * @return TLE
157      */
158     private TLE createTLE() {
159         final double averagedMeanMotion = computeMeanMotion(getAveragedElements()
160                 .getAveragedSemiMajorAxis());
161         final double averagedEccentricity = getAveragedElements().getAveragedEccentricity();
162         final double averagedInclination = getAveragedElements().getAveragedInclination();
163         final double averagedRAAN = getAveragedElements().getAveragedRightAscensionOfTheAscendingNode();
164         final double averagedPerigeeArgument = getAveragedElements().getAveragedPerigeeArgument();
165         final double averagedMeanAnomaly = getAveragedElements().getAveragedMeanAnomaly();
166         return new TLE(0, (char) 0, 2000, 1, "1", 0, 0,
167                 getDate(), averagedMeanMotion, 0., 0., averagedEccentricity,
168                 averagedInclination, averagedPerigeeArgument, averagedRAAN, averagedMeanAnomaly,
169                 1, B_STAR, utc);
170     }
171 
172     /**
173      * Convert averaged semi-major axis to averaged mean motion. Uses an approximate transformation
174      * (same as osculating).
175      * @param averagedSemiMajorAxis semi-major axis
176      * @return mean motion
177      */
178     private static double computeMeanMotion(final double averagedSemiMajorAxis) {
179         final double cubedSemiMajorAxis = averagedSemiMajorAxis * averagedSemiMajorAxis *
180                 averagedSemiMajorAxis;
181         return FastMath.sqrt(getTleMu() / cubedSemiMajorAxis);
182     }
183 
184     /**
185      * Compute semi-major axis from Two-Line Elements. Uses an approximate transformation
186      * (same as osculating).
187      * @param tle TLE
188      * @return semi-major axis
189      */
190     private static double computeSemiMajorAxis(final TLE tle) {
191         final double meanMotion = tle.getMeanMotion();
192         return FastMath.cbrt(getTleMu() / (meanMotion * meanMotion));
193     }
194 }