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