1   /* Copyright 2002-2025 CS GROUP
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.analytical.tle;
18  
19  import org.hipparchus.CalculusFieldElement;
20  import org.hipparchus.util.FastMath;
21  import org.hipparchus.util.FieldSinCos;
22  import org.orekit.annotation.DefaultDataContext;
23  import org.orekit.attitudes.AttitudeProvider;
24  import org.orekit.data.DataContext;
25  import org.orekit.frames.Frame;
26  
27  /** This class contains methods to compute propagated coordinates with the SGP4 model.
28   * <p>
29   * The user should not bother in this class since it is handled internaly by the
30   * {@link TLEPropagator}.
31   * </p>
32   * <p>This implementation is largely inspired from the paper and source code <a
33   * href="https://www.celestrak.com/publications/AIAA/2006-6753/">Revisiting Spacetrack
34   * Report #3</a> and is fully compliant with its results and tests cases.</p>
35   * @author Felix R. Hoots, Ronald L. Roehrich, December 1980 (original fortran)
36   * @author David A. Vallado, Paul Crawford, Richard Hujsak, T.S. Kelso (C++ translation and improvements)
37   * @author Fabien Maussion (java translation)
38   * @author Thomas Paulet (field translation)
39   * @since 11.0
40   * @param <T> type of the field elements
41   */
42  public class FieldSGP4<T extends CalculusFieldElement<T>> extends FieldTLEPropagator<T> {
43  
44      /** If perige is less than 220 km, some calculus are avoided. */
45      private boolean lessThan220;
46  
47      /** (1 + eta * cos(M0))³. */
48      private T delM0;
49  
50      // CHECKSTYLE: stop JavadocVariable check
51      private T d2;
52      private T d3;
53      private T d4;
54      private T t3cof;
55      private T t4cof;
56      private T t5cof;
57      private T sinM0;
58      private T omgcof;
59      private T xmcof;
60      private T c5;
61      // CHECKSTYLE: resume JavadocVariable check
62  
63      /** Constructor for a unique initial TLE.
64       *
65       * <p>This constructor uses the {@link DataContext#getDefault() default data context}.
66       *
67       * @param initialTLE the TLE to propagate.
68       * @param attitudeProvider provider for attitude computation
69       * @param mass spacecraft mass (kg)
70       * @param parameters SGP4 and SDP4 model parameters
71       * @see #FieldSGP4(FieldTLE, AttitudeProvider, CalculusFieldElement, Frame, CalculusFieldElement[])
72       */
73      @DefaultDataContext
74      public FieldSGP4(final FieldTLE<T> initialTLE, final AttitudeProvider attitudeProvider,
75                  final T mass, final T[] parameters) {
76          this(initialTLE, attitudeProvider, mass,
77                  DataContext.getDefault().getFrames().getTEME(), parameters);
78      }
79  
80      /** Constructor for a unique initial TLE.
81       * @param initialTLE the TLE to propagate.
82       * @param attitudeProvider provider for attitude computation
83       * @param mass spacecraft mass (kg)
84       * @param teme the TEME frame to use for propagation.
85       * @param parameters SGP4 and SDP4 model parameters
86       */
87      public FieldSGP4(final FieldTLE<T> initialTLE,
88                  final AttitudeProvider attitudeProvider,
89                  final T mass,
90                  final Frame teme,
91                  final T[] parameters) {
92          super(initialTLE, attitudeProvider, mass, teme, parameters);
93      }
94  
95      /** Initialization proper to each propagator (SGP or SDP).
96       * @param parameters model parameters
97       */
98      protected void sxpInitialize(final T[] parameters) {
99  
100         final T bStar = parameters[0];
101         // For perigee less than 220 kilometers, the equations are truncated to
102         // linear variation in sqrt a and quadratic variation in mean anomaly.
103         // Also, the c3 term, the delta omega term, and the delta m term are dropped.
104         lessThan220 = perige.getReal() < 220;
105         if (!lessThan220) {
106             final FieldSinCos<T> scM0 = FastMath.sinCos(tle.getMeanAnomaly());
107             final T c1sq = c1.square();
108             delM0 = eta.multiply(scM0.cos()).add(1.0);
109             delM0 = delM0.multiply(delM0).multiply(delM0);
110             d2 = a0dp.multiply(tsi).multiply(c1sq).multiply(4.0);
111             final T temp = d2.multiply(tsi).multiply(c1).divide(3.0);
112             d3 = a0dp.multiply(17.0).add(s4).multiply(temp);
113             d4 = temp.multiply(0.5).multiply(a0dp).multiply(tsi).multiply(a0dp.multiply(221.0).add(s4.multiply(31.0))).multiply(c1);
114             t3cof = d2.add(c1sq.multiply(2));
115             t4cof = d3.multiply(3.0).add(c1.multiply(d2.multiply(12.0).add(c1sq.multiply(10)))).multiply(0.25);
116             t5cof = d4.multiply(3.0).add(c1.multiply(12.0).multiply(d3)).add(
117                     d2.multiply(d2).multiply(6.0)).add(c1sq.multiply(15.0).multiply(d2.multiply(2).add(c1sq))).multiply(0.2);
118             sinM0 = scM0.sin();
119             if (tle.getE().getReal() < 1e-4) {
120                 omgcof = c1sq.getField().getZero();
121                 xmcof = c1sq.getField().getZero();
122             } else  {
123                 final T c3 = coef.multiply(tsi).multiply(xn0dp).multiply(TLEConstants.A3OVK2 * TLEConstants.NORMALIZED_EQUATORIAL_RADIUS).multiply(sini0.divide(tle.getE()));
124                 xmcof = coef.multiply(bStar).divide(eeta).multiply(-TLEConstants.TWO_THIRD * TLEConstants.NORMALIZED_EQUATORIAL_RADIUS);
125                 omgcof = bStar.multiply(c3).multiply(FastMath.cos(tle.getPerigeeArgument()));
126             }
127         }
128 
129         c5 = coef1.multiply(2).multiply(a0dp).multiply(beta02).multiply(etasq.add(eeta).multiply(2.75).add(eeta.multiply(etasq)).add(1));
130         // initialized
131     }
132 
133     /** Propagation proper to each propagator (SGP or SDP).
134      * @param tSince the offset from initial epoch (min)
135      * @param parameters model parameters
136      */
137     protected void sxpPropagate(final T tSince, final T[] parameters) {
138 
139         // Update for secular gravity and atmospheric drag.
140         final T bStar = parameters[0];
141         final T xmdf = tle.getMeanAnomaly().add(xmdot.multiply(tSince));
142         final T omgadf = tle.getPerigeeArgument().add(omgdot.multiply(tSince));
143         final T xn0ddf = tle.getRaan().add(xnodot.multiply(tSince));
144         omega = omgadf;
145         T xmp = xmdf;
146         final T tsq = tSince.square();
147         xnode = xn0ddf.add(xnodcf.multiply(tsq));
148         T tempa = c1.multiply(tSince).negate().add(1.0);
149         T tempe = bStar.multiply(c4).multiply(tSince);
150         T templ = t2cof.multiply(tsq);
151 
152         if (!lessThan220) {
153             final T delomg = omgcof.multiply(tSince);
154             T delm = eta.multiply(FastMath.cos(xmdf)).add(1.0);
155             delm = xmcof.multiply(delm.square().multiply(delm).subtract(delM0));
156             final T temp = delomg.add(delm);
157             xmp = xmdf.add(temp);
158             omega = omgadf.subtract(temp);
159             final T tcube = tsq.multiply(tSince);
160             final T tfour = tSince.multiply(tcube);
161             tempa = tempa.subtract(d2.multiply(tsq)).subtract(d3.multiply(tcube)).subtract(d4.multiply(tfour));
162             tempe = tempe.add(bStar.multiply(c5).multiply(FastMath.sin(xmp).subtract(sinM0)));
163             templ = templ.add(t3cof.multiply(tcube)).add(tfour.multiply(t4cof.add(tSince.multiply(t5cof))));
164         }
165 
166         a = a0dp.multiply(tempa).multiply(tempa);
167         e = tle.getE().subtract(tempe);
168 
169         // A highly arbitrary lower limit on e,  of 1e-6:
170         if (e.getReal() < 1e-6) {
171             e = e.getField().getZero().newInstance(1e-6);
172         }
173 
174         xl = xmp.add(omega).add(xnode).add(xn0dp.multiply(templ));
175 
176         i = tle.getI();
177 
178     }
179 
180 }