/[PAMELA software]/quicklook/OrbitalRate/src/sgp4.cpp
ViewVC logotype

Annotation of /quicklook/OrbitalRate/src/sgp4.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.1 - (hide annotations) (download)
Tue Dec 5 19:49:19 2006 UTC (17 years, 11 months ago) by pam-rm2
Branch: MAIN
CVS Tags: v2r02, v2r01, v2r00, HEAD
New version of OrbitalRate quicklook.  Initial import.
Nico

1 pam-rm2 1.1 //
2     // globals.cpp
3     //
4     #include <sgp4.h>
5    
6     //////////////////////////////////////////////////////////////////////////////
7     double sqr(const double x)
8     {
9     return (x * x);
10     }
11    
12     //////////////////////////////////////////////////////////////////////////////
13     double Fmod2p(const double arg)
14     {
15     double modu = fmod(arg, TWOPI);
16    
17     if (modu < 0.0)
18     modu += TWOPI;
19    
20     return modu;
21     }
22    
23     //////////////////////////////////////////////////////////////////////////////
24     // AcTan()
25     // ArcTangent of sin(x) / cos(x). The advantage of this function over arctan()
26     // is that it returns the correct quadrant of the angle.
27     double AcTan(const double sinx, const double cosx)
28     {
29     double ret;
30    
31     if (cosx == 0.0)
32     {
33     if (sinx > 0.0)
34     ret = PI / 2.0;
35     else
36     ret = 3.0 * PI / 2.0;
37     }
38     else
39     {
40     if (cosx > 0.0)
41     ret = atan(sinx / cosx);
42     else
43     ret = PI + atan(sinx / cosx);
44     }
45    
46     return ret;
47     }
48    
49     //////////////////////////////////////////////////////////////////////////////
50     double rad2deg(const double r)
51     {
52     const double DEG_PER_RAD = 180.0 / PI;
53     return r * DEG_PER_RAD;
54     }
55    
56     //////////////////////////////////////////////////////////////////////////////
57     double deg2rad(const double d)
58     {
59     const double RAD_PER_DEG = PI / 180.0;
60     return d * RAD_PER_DEG;
61     }
62    
63     //
64     // coord.cpp
65     //
66     // Copyright (c) 2003 Michael F. Henry
67     //
68    
69     //////////////////////////////////////////////////////////////////////
70     // cCoordGeo Class
71     //////////////////////////////////////////////////////////////////////
72    
73     cCoordGeo::cCoordGeo()
74     {
75     m_Lat = 0.0;
76     m_Lon = 0.0;
77     m_Alt = 0.0;
78     }
79    
80     //////////////////////////////////////////////////////////////////////
81     // cCoordTopo Class
82     //////////////////////////////////////////////////////////////////////
83    
84     cCoordTopo::cCoordTopo()
85     {
86     m_Az = 0.0;
87     m_El = 0.0;
88     m_Range = 0.0;
89     m_RangeRate = 0.0;
90    
91     }
92    
93    
94    
95     //
96     // cVector.cpp
97     //
98     // Copyright (c) 2001-2003 Michael F. Henry
99     //
100     //*****************************************************************************
101     // Multiply each component in the vector by 'factor'.
102     //*****************************************************************************
103     void cVector::Mul(double factor)
104     {
105     m_x *= factor;
106     m_y *= factor;
107     m_z *= factor;
108     m_w *= fabs(factor);
109     }
110    
111     //*****************************************************************************
112     // Subtract a vector from this one.
113     //*****************************************************************************
114     void cVector::Sub(const cVector& vec)
115     {
116     m_x -= vec.m_x;
117     m_y -= vec.m_y;
118     m_z -= vec.m_z;
119     m_w -= vec.m_w;
120     }
121    
122     //*****************************************************************************
123     // Calculate the angle between this vector and another
124     //*****************************************************************************
125     double cVector::Angle(const cVector& vec) const
126     {
127     return acos(Dot(vec) / (Magnitude() * vec.Magnitude()));
128     }
129    
130     //*****************************************************************************
131     //
132     //*****************************************************************************
133     double cVector::Magnitude() const
134     {
135     return sqrt((m_x * m_x) +
136     (m_y * m_y) +
137     (m_z * m_z));
138     }
139    
140     //*****************************************************************************
141     // Return the dot product
142     //*****************************************************************************
143     double cVector::Dot(const cVector& vec) const
144     {
145     return (m_x * vec.m_x) +
146     (m_y * vec.m_y) +
147     (m_z * vec.m_z);
148     }
149     //
150     // cJulian.cpp
151     //
152     // This class encapsulates Julian dates with the epoch of 12:00 noon (12:00 UT)
153     // on January 1, 4713 B.C. Some epoch dates:
154     // 01/01/1990 00:00 UTC - 2447892.5
155     // 01/01/1990 12:00 UTC - 2447893.0
156     // 01/01/2000 00:00 UTC - 2451544.5
157     // 01/01/2001 00:00 UTC - 2451910.5
158     //
159     // Note the Julian day begins at noon, which allows astronomers to have all
160     // the dates in a single observing session the same.
161     //
162     // References:
163     // "Astronomical Formulae for Calculators", Jean Meeus
164     // "Satellite Communications", Dennis Roddy, 2nd Edition, 1995.
165     //
166     // Copyright (c) 2003 Michael F. Henry
167     //
168     // mfh 12/24/2003
169     //
170    
171     //////////////////////////////////////////////////////////////////////////////
172     // Create a Julian date object from a time_t object. time_t objects store the
173     // number of seconds since midnight UTC January 1, 1970.
174     cJulian::cJulian(time_t time)
175     {
176     struct tm *ptm = gmtime(&time);
177     assert(ptm);
178    
179     int year = ptm->tm_year + 1900;
180     double day = ptm->tm_yday + 1 +
181     (ptm->tm_hour +
182     ((ptm->tm_min +
183     (ptm->tm_sec / 60.0)) / 60.0)) / 24.0;
184    
185     Initialize(year, day);
186     }
187    
188     //////////////////////////////////////////////////////////////////////////////
189     // Create a Julian date object from a year and day of year.
190     // Example parameters: year = 2001, day = 1.5 (Jan 1 12h)
191     cJulian::cJulian(int year, double day)
192     {
193     Initialize(year, day);
194     }
195    
196     //////////////////////////////////////////////////////////////////////////////
197     // Create a Julian date object.
198     cJulian::cJulian(int year, // i.e., 2004
199     int mon, // 1..12
200     int day, // 1..31
201     int hour, // 0..23
202     int min, // 0..59
203     double sec /* = 0.0 */) // 0..(59.999999...)
204    
205     {
206     // Calculate N, the day of the year (1..366)
207     int N;
208     int F1 = (int)((275.0 * mon) / 9.0);
209     int F2 = (int)((mon + 9.0) / 12.0);
210    
211     if (IsLeapYear(year))
212     {
213     // Leap year
214     N = F1 - F2 + day - 30;
215     }
216     else
217     {
218     // Common year
219     N = F1 - (2 * F2) + day - 30;
220     }
221    
222     double dblDay = N + (hour + (min + (sec / 60.0)) / 60.0) / 24.0;
223    
224     Initialize(year, dblDay);
225     }
226    
227     //////////////////////////////////////////////////////////////////////////////
228     void cJulian::Initialize(int year, double day)
229     {
230     // 1582 A.D.: 10 days removed from calendar
231     // 3000 A.D.: Arbitrary error checking limit
232     assert((year > 1582) && (year < 3000));
233     assert((day >= 0.0) && (day <= 366.5));
234    
235     // Now calculate Julian date
236    
237     year--;
238    
239     // Centuries are not leap years unless they divide by 400
240     int A = (year / 100);
241     int B = 2 - A + (A / 4);
242    
243     double NewYears = (int)(365.25 * year) +
244     (int)(30.6001 * 14) +
245     1720994.5 + B; // 1720994.5 = Oct 30, year -1
246    
247     m_Date = NewYears + day;
248     }
249    
250     //////////////////////////////////////////////////////////////////////////////
251     // getComponent()
252     // Return requested components of date.
253     // Year : Includes the century.
254     // Month: 1..12
255     // Day : 1..31 including fractional part
256     void cJulian::getComponent(int *pYear,
257     int *pMon /* = NULL */,
258     double *pDOM /* = NULL */) const
259     {
260     assert(pYear != NULL);
261    
262     double jdAdj = getDate() + 0.5;
263     int Z = (int)jdAdj; // integer part
264     double F = jdAdj - Z; // fractional part
265     double alpha = (int)((Z - 1867216.25) / 36524.25);
266     double A = Z + 1 + alpha - (int)(alpha / 4.0);
267     double B = A + 1524.0;
268     int C = (int)((B - 122.1) / 365.25);
269     int D = (int)(C * 365.25);
270     int E = (int)((B - D) / 30.6001);
271    
272     double DOM = B - D - (int)(E * 30.6001) + F;
273     int month = (E < 13.5) ? (E - 1) : (E - 13);
274     int year = (month > 2.5) ? (C - 4716) : (C - 4715);
275    
276     *pYear = year;
277    
278     if (pMon != NULL)
279     *pMon = month;
280    
281     if (pDOM != NULL)
282     *pDOM = DOM;
283     }
284    
285     //////////////////////////////////////////////////////////////////////////////
286     // toGMST()
287     // Calculate Greenwich Mean Sidereal Time for the Julian date. The return value
288     // is the angle, in radians, measuring eastward from the Vernal Equinox to the
289     // prime meridian. This angle is also referred to as "ThetaG" (Theta GMST).
290     //
291     // References:
292     // The 1992 Astronomical Almanac, page B6.
293     // Explanatory Supplement to the Astronomical Almanac, page 50.
294     // Orbital Coordinate Systems, Part III, Dr. T.S. Kelso, Satellite Times,
295     // Nov/Dec 1995
296     double cJulian::toGMST() const
297     {
298     const double UT = fmod(m_Date + 0.5, 1.0);
299     const double TU = (FromJan1_12h_2000() - UT) / 36525.0;
300    
301     double GMST = 24110.54841 + TU *
302     (8640184.812866 + TU * (0.093104 - TU * 6.2e-06));
303    
304     GMST = fmod(GMST + SEC_PER_DAY * OMEGA_E * UT, SEC_PER_DAY);
305    
306     if (GMST < 0.0)
307     GMST += SEC_PER_DAY; // "wrap" negative modulo value
308    
309     return (TWOPI * (GMST / SEC_PER_DAY));
310     }
311    
312     //////////////////////////////////////////////////////////////////////////////
313     // toLMST()
314     // Calculate Local Mean Sidereal Time for given longitude (for this date).
315     // The longitude is assumed to be in radians measured west from Greenwich.
316     // The return value is the angle, in radians, measuring eastward from the
317     // Vernal Equinox to the given longitude.
318     double cJulian::toLMST(double lon) const
319     {
320     return fmod(toGMST() + lon, TWOPI);
321     }
322    
323     //////////////////////////////////////////////////////////////////////////////
324     // toTime()
325     // Convert to type time_t
326     // Avoid using this function as it discards the fractional seconds of the
327     // time component.
328     time_t cJulian::toTime() const
329     {
330     int nYear;
331     int nMonth;
332     double dblDay;
333    
334     getComponent(&nYear, &nMonth, &dblDay);
335    
336     // dblDay is the fractional Julian Day (i.e., 29.5577).
337     // Save the whole number day in nDOM and convert dblDay to
338     // the fractional portion of day.
339     int nDOM = (int)dblDay;
340    
341     dblDay -= nDOM;
342    
343     const int SEC_PER_MIN = 60;
344     const int SEC_PER_HR = 60 * SEC_PER_MIN;
345     const int SEC_PER_DAY = 24 * SEC_PER_HR;
346    
347     int secs = (int)((dblDay * SEC_PER_DAY) + 0.5);
348    
349     // Create a "struct tm" type.
350     // NOTE:
351     // The "struct tm" type has a 1-second resolution. Any fractional
352     // component of the "seconds" time value is discarded.
353     struct tm tGMT;
354     memset(&tGMT, 0, sizeof(tGMT));
355    
356     tGMT.tm_year = nYear - 1900; // 2001 is 101
357     tGMT.tm_mon = nMonth - 1; // January is 0
358     tGMT.tm_mday = nDOM; // First day is 1
359     tGMT.tm_hour = secs / SEC_PER_HR;
360     tGMT.tm_min = (secs % SEC_PER_HR) / SEC_PER_MIN;
361     tGMT.tm_sec = (secs % SEC_PER_HR) % SEC_PER_MIN;
362     tGMT.tm_isdst = 0; // No conversion desired
363    
364     time_t tEpoch = mktime(&tGMT);
365    
366     if (tEpoch != -1)
367     {
368     // Valid time_t value returned from mktime().
369     // mktime() expects a local time which means that tEpoch now needs
370     // to be adjusted by the difference between this time zone and GMT.
371     tEpoch -= timezone;
372     }
373    
374     return tEpoch;
375     }
376     //
377     // cTle.cpp
378     // This class encapsulates a single set of standard NORAD two line elements.
379     //
380     // Copyright 1996-2005 Michael F. Henry
381     //
382     // Note: The column offsets are ZERO based.
383    
384     // Name
385     const int TLE_LEN_LINE_DATA = 69; const int TLE_LEN_LINE_NAME = 22;
386    
387     // Line 1
388     const int TLE1_COL_SATNUM = 2; const int TLE1_LEN_SATNUM = 5;
389     const int TLE1_COL_INTLDESC_A = 9; const int TLE1_LEN_INTLDESC_A = 2;
390     const int TLE1_COL_INTLDESC_B = 11; const int TLE1_LEN_INTLDESC_B = 3;
391     const int TLE1_COL_INTLDESC_C = 14; const int TLE1_LEN_INTLDESC_C = 3;
392     const int TLE1_COL_EPOCH_A = 18; const int TLE1_LEN_EPOCH_A = 2;
393     const int TLE1_COL_EPOCH_B = 20; const int TLE1_LEN_EPOCH_B = 12;
394     const int TLE1_COL_MEANMOTIONDT = 33; const int TLE1_LEN_MEANMOTIONDT = 10;
395     const int TLE1_COL_MEANMOTIONDT2 = 44; const int TLE1_LEN_MEANMOTIONDT2 = 8;
396     const int TLE1_COL_BSTAR = 53; const int TLE1_LEN_BSTAR = 8;
397     const int TLE1_COL_EPHEMTYPE = 62; const int TLE1_LEN_EPHEMTYPE = 1;
398     const int TLE1_COL_ELNUM = 64; const int TLE1_LEN_ELNUM = 4;
399    
400     // Line 2
401     const int TLE2_COL_SATNUM = 2; const int TLE2_LEN_SATNUM = 5;
402     const int TLE2_COL_INCLINATION = 8; const int TLE2_LEN_INCLINATION = 8;
403     const int TLE2_COL_RAASCENDNODE = 17; const int TLE2_LEN_RAASCENDNODE = 8;
404     const int TLE2_COL_ECCENTRICITY = 26; const int TLE2_LEN_ECCENTRICITY = 7;
405     const int TLE2_COL_ARGPERIGEE = 34; const int TLE2_LEN_ARGPERIGEE = 8;
406     const int TLE2_COL_MEANANOMALY = 43; const int TLE2_LEN_MEANANOMALY = 8;
407     const int TLE2_COL_MEANMOTION = 52; const int TLE2_LEN_MEANMOTION = 11;
408     const int TLE2_COL_REVATEPOCH = 63; const int TLE2_LEN_REVATEPOCH = 5;
409    
410     /////////////////////////////////////////////////////////////////////////////
411     cTle::cTle(string& strName, string& strLine1, string& strLine2)
412     {
413     m_strName = strName;
414     m_strLine1 = strLine1;
415     m_strLine2 = strLine2;
416    
417     Initialize();
418     }
419    
420     /////////////////////////////////////////////////////////////////////////////
421     cTle::cTle(const cTle &tle)
422     {
423     m_strName = tle.m_strName;
424     m_strLine1 = tle.m_strLine1;
425     m_strLine2 = tle.m_strLine2;
426    
427     for (int fld = FLD_FIRST; fld < FLD_LAST; fld++)
428     {
429     m_Field[fld] = tle.m_Field[fld];
430     }
431    
432     m_mapCache = tle.m_mapCache;
433     }
434    
435     /////////////////////////////////////////////////////////////////////////////
436     cTle::~cTle()
437     {
438     }
439    
440     /////////////////////////////////////////////////////////////////////////////
441     // getField()
442     // Return requested field as a double (function return value) or as a text
443     // string (*pstr) in the units requested (eUnit). Set 'bStrUnits' to true
444     // to have units appended to text string.
445     //
446     // Note: numeric return values are cached; asking for the same field more
447     // than once incurs minimal overhead.
448     double cTle::getField(eField fld,
449     eUnits units, /* = U_NATIVE */
450     string *pstr /* = NULL */,
451     bool bStrUnits /* = false */) const
452     {
453     assert((FLD_FIRST <= fld) && (fld < FLD_LAST));
454     assert((U_FIRST <= units) && (units < U_LAST));
455    
456     if (pstr)
457     {
458     // Return requested field in string form.
459     *pstr = m_Field[fld];
460    
461     if (bStrUnits)
462     *pstr += getUnits(fld);
463    
464     return 0.0;
465     }
466     else
467     {
468     // Return requested field in floating-point form.
469     // Return cache contents if it exists, else populate cache
470     FldKey key = Key(units, fld);
471    
472     if (m_mapCache.find(key) == m_mapCache.end())
473     {
474     // Value not in cache; add it
475     double valNative = atof(m_Field[fld].c_str());
476     double valConv = ConvertUnits(valNative, fld, units);
477     m_mapCache[key] = valConv;
478    
479     return valConv;
480     }
481     else
482     {
483     // return cached value
484     return m_mapCache[key];
485     }
486     }
487     }
488    
489     //////////////////////////////////////////////////////////////////////////////
490     // Convert the given field into the requested units. It is assumed that
491     // the value being converted is in the TLE format's "native" form.
492     double cTle::ConvertUnits(double valNative, // value to convert
493     eField fld, // what field the value is
494     eUnits units) // what units to convert to
495     {
496     switch (fld)
497     {
498     case FLD_I:
499     case FLD_RAAN:
500     case FLD_ARGPER:
501     case FLD_M:
502     {
503     // The native TLE format is DEGREES
504     if (units == U_RAD)
505     return valNative * RADS_PER_DEG;
506     }
507    
508     case FLD_NORADNUM:
509     case FLD_INTLDESC:
510     case FLD_SET:
511     case FLD_EPOCHYEAR:
512     case FLD_EPOCHDAY:
513     case FLD_ORBITNUM:
514     case FLD_E:
515     case FLD_MMOTION:
516     case FLD_MMOTIONDT:
517     case FLD_MMOTIONDT2:
518     case FLD_BSTAR:
519     case FLD_LAST:
520     { // do nothing
521    
522     }
523    
524     }
525    
526     return valNative; // return value in unconverted native format
527     }
528    
529     //////////////////////////////////////////////////////////////////////////////
530     string cTle::getUnits(eField fld) const
531     {
532     static const string strDegrees = " degrees";
533     static const string strRevsPerDay = " revs / day";
534     static const string strNull;
535    
536     switch (fld)
537     {
538     case FLD_I:
539     case FLD_RAAN:
540     case FLD_ARGPER:
541     case FLD_M:
542     return strDegrees;
543    
544     case FLD_MMOTION:
545     return strRevsPerDay;
546    
547     default:
548     return strNull;
549     }
550     }
551    
552     /////////////////////////////////////////////////////////////////////////////
553     // ExpToDecimal()
554     // Converts TLE-style exponential notation of the form [ |-]00000[+|-]0 to
555     // decimal notation. Assumes implied decimal point to the left of the first
556     // number in the string, i.e.,
557     // " 12345-3" = 0.00012345
558     // "-23429-5" = -0.0000023429
559     // " 40436+1" = 4.0436
560     string cTle::ExpToDecimal(const string &str)
561     {
562     const int COL_EXP_SIGN = 6;
563     const int LEN_EXP = 2;
564    
565     const int LEN_BUFREAL = 32; // max length of buffer to hold floating point
566     // representation of input string.
567     int nMan;
568     int nExp;
569    
570     // sscanf(%d) will read up to the exponent sign
571     sscanf(str.c_str(), "%d", &nMan);
572    
573     double dblMan = nMan;
574     bool bNeg = (nMan < 0);
575    
576     if (bNeg)
577     dblMan *= -1;
578    
579     // Move decimal place to left of first digit
580     while (dblMan >= 1.0)
581     dblMan /= 10.0;
582    
583     if (bNeg)
584     dblMan *= -1;
585    
586     // now read exponent
587     sscanf(str.substr(COL_EXP_SIGN, LEN_EXP).c_str(), "%d", &nExp);
588    
589     double dblVal = dblMan * pow(10.0, nExp);
590     char szVal[LEN_BUFREAL];
591    
592     snprintf(szVal, sizeof(szVal), "%.9f", dblVal);
593    
594     string strVal = szVal;
595    
596     return strVal;
597    
598     } // ExpToDecimal()
599    
600     /////////////////////////////////////////////////////////////////////////////
601     // Initialize()
602     // Initialize the string array.
603     void cTle::Initialize()
604     {
605     // Have we already been initialized?
606     if (m_Field[FLD_NORADNUM].size())
607     return;
608    
609     assert(!m_strName.empty());
610     assert(!m_strLine1.empty());
611     assert(!m_strLine2.empty());
612    
613     m_Field[FLD_NORADNUM] = m_strLine1.substr(TLE1_COL_SATNUM, TLE1_LEN_SATNUM);
614     m_Field[FLD_INTLDESC] = m_strLine1.substr(TLE1_COL_INTLDESC_A,
615     TLE1_LEN_INTLDESC_A +
616     TLE1_LEN_INTLDESC_B +
617     TLE1_LEN_INTLDESC_C);
618     m_Field[FLD_EPOCHYEAR] =
619     m_strLine1.substr(TLE1_COL_EPOCH_A, TLE1_LEN_EPOCH_A);
620    
621     m_Field[FLD_EPOCHDAY] =
622     m_strLine1.substr(TLE1_COL_EPOCH_B, TLE1_LEN_EPOCH_B);
623    
624     if (m_strLine1[TLE1_COL_MEANMOTIONDT] == '-')
625     {
626     // value is negative
627     m_Field[FLD_MMOTIONDT] = "-0";
628     }
629     else
630     m_Field[FLD_MMOTIONDT] = "0";
631    
632     m_Field[FLD_MMOTIONDT] += m_strLine1.substr(TLE1_COL_MEANMOTIONDT + 1,
633     TLE1_LEN_MEANMOTIONDT);
634    
635     // decimal point assumed; exponential notation
636     m_Field[FLD_MMOTIONDT2] = ExpToDecimal(
637     m_strLine1.substr(TLE1_COL_MEANMOTIONDT2,
638     TLE1_LEN_MEANMOTIONDT2));
639     // decimal point assumed; exponential notation
640     m_Field[FLD_BSTAR] = ExpToDecimal(m_strLine1.substr(TLE1_COL_BSTAR,
641     TLE1_LEN_BSTAR));
642     //TLE1_COL_EPHEMTYPE
643     //TLE1_LEN_EPHEMTYPE
644     m_Field[FLD_SET] = m_strLine1.substr(TLE1_COL_ELNUM, TLE1_LEN_ELNUM);
645    
646     TrimLeft(m_Field[FLD_SET]);
647    
648     //TLE2_COL_SATNUM
649     //TLE2_LEN_SATNUM
650    
651     m_Field[FLD_I] = m_strLine2.substr(TLE2_COL_INCLINATION,
652     TLE2_LEN_INCLINATION);
653     TrimLeft(m_Field[FLD_I]);
654    
655     m_Field[FLD_RAAN] = m_strLine2.substr(TLE2_COL_RAASCENDNODE,
656     TLE2_LEN_RAASCENDNODE);
657     TrimLeft(m_Field[FLD_RAAN]);
658    
659     // decimal point is assumed
660     m_Field[FLD_E] = "0.";
661     m_Field[FLD_E] += m_strLine2.substr(TLE2_COL_ECCENTRICITY,
662     TLE2_LEN_ECCENTRICITY);
663    
664     m_Field[FLD_ARGPER] = m_strLine2.substr(TLE2_COL_ARGPERIGEE,
665     TLE2_LEN_ARGPERIGEE);
666     TrimLeft(m_Field[FLD_ARGPER]);
667    
668     m_Field[FLD_M] = m_strLine2.substr(TLE2_COL_MEANANOMALY,
669     TLE2_LEN_MEANANOMALY);
670     TrimLeft(m_Field[FLD_M]);
671    
672     m_Field[FLD_MMOTION] = m_strLine2.substr(TLE2_COL_MEANMOTION,
673     TLE2_LEN_MEANMOTION);
674     TrimLeft(m_Field[FLD_MMOTION]);
675    
676     m_Field[FLD_ORBITNUM] = m_strLine2.substr(TLE2_COL_REVATEPOCH,
677     TLE2_LEN_REVATEPOCH);
678     TrimLeft(m_Field[FLD_ORBITNUM]);
679    
680     } // InitStrVars()
681    
682     /////////////////////////////////////////////////////////////////////////////
683     // IsTleFormat()
684     // Returns true if "str" is a valid data line of a two-line element set,
685     // else false.
686     //
687     // To be valid a line must:
688     // Have as the first character the line number
689     // Have as the second character a blank
690     // Be TLE_LEN_LINE_DATA characters long
691     // Have a valid checksum (note: no longer required as of 12/96)
692     //
693     bool cTle::IsValidLine(string& str, eTleLine line)
694     {
695     TrimLeft(str);
696     TrimRight(str);
697    
698     size_t nLen = str.size();
699    
700     if (nLen != (uint)TLE_LEN_LINE_DATA)
701     return false;
702    
703     // First char in string must be line number
704     if ((str[0] - '0') != line)
705     return false;
706    
707     // Second char in string must be blank
708     if (str[1] != ' ')
709     return false;
710    
711     /*
712     NOTE: 12/96
713     The requirement that the last char in the line data must be a valid
714     checksum is too restrictive.
715    
716     // Last char in string must be checksum
717     int nSum = CheckSum(str);
718    
719     if (nSum != (str[TLE_LEN_LINE_DATA - 1] - '0'))
720     return false;
721     */
722    
723     return true;
724    
725     } // IsTleFormat()
726    
727     /////////////////////////////////////////////////////////////////////////////
728     // CheckSum()
729     // Calculate the check sum for a given line of TLE data, the last character
730     // of which is the current checksum. (Although there is no check here,
731     // the current checksum should match the one we calculate.)
732     // The checksum algorithm:
733     // Each number in the data line is summed, modulo 10.
734     // Non-numeric characters are zero, except minus signs, which are 1.
735     //
736     int cTle::CheckSum(const string& str)
737     {
738     // The length is "- 1" because we don't include the current (existing)
739     // checksum character in the checksum calculation.
740     size_t len = str.size() - 1;
741     int xsum = 0;
742    
743     for (size_t i = 0; i < len; i++)
744     {
745     char ch = str[i];
746     if (isdigit(ch))
747     xsum += (ch - '0');
748     else if (ch == '-')
749     xsum++;
750     }
751    
752     return (xsum % 10);
753    
754     } // CheckSum()
755    
756     /////////////////////////////////////////////////////////////////////////////
757     void cTle::TrimLeft(string& s)
758     {
759     while (s[0] == ' ')
760     s.erase(0, 1);
761     }
762    
763     /////////////////////////////////////////////////////////////////////////////
764     void cTle::TrimRight(string& s)
765     {
766     while (s[s.size() - 1] == ' ')
767     s.erase(s.size() - 1);
768     }
769    
770     //
771     // cEci.cpp
772     //
773     // Copyright (c) 2002-2003 Michael F. Henry
774     //
775     //////////////////////////////////////////////////////////////////////
776     // cEci Class
777     //////////////////////////////////////////////////////////////////////
778     cEci::cEci(const cVector &pos,
779     const cVector &vel,
780     const cJulian &date,
781     bool IsAeUnits /* = true */)
782     {
783     m_pos = pos;
784     m_vel = vel;
785     m_date = date;
786     m_VecUnits = (IsAeUnits ? UNITS_AE : UNITS_NONE);
787     }
788    
789     //////////////////////////////////////////////////////////////////////
790     // cEci(cCoordGeo&, cJulian&)
791     // Calculate the ECI coordinates of the location "geo" at time "date".
792     // Assumes geo coordinates are km-based.
793     // Assumes the earth is an oblate spheroid as defined in WGS '72.
794     // Reference: The 1992 Astronomical Almanac, page K11
795     // Reference: www.celestrak.com (Dr. TS Kelso)
796     cEci::cEci(const cCoordGeo &geo, const cJulian &date)
797     {
798     m_VecUnits = UNITS_KM;
799    
800     double mfactor = TWOPI * (OMEGA_E / SEC_PER_DAY);
801     double lat = geo.m_Lat;
802     double lon = geo.m_Lon;
803     double alt = geo.m_Alt;
804    
805     // Calculate Local Mean Sidereal Time (theta)
806     double theta = date.toLMST(lon);
807     double c = 1.0 / sqrt(1.0 + F * (F - 2.0) * sqr(sin(lat)));
808     double s = sqr(1.0 - F) * c;
809     double achcp = (XKMPER_WGS72 * c + alt) * cos(lat);
810    
811     m_date = date;
812    
813     m_pos.m_x = achcp * cos(theta); // km
814     m_pos.m_y = achcp * sin(theta); // km
815     m_pos.m_z = (XKMPER_WGS72 * s + alt) * sin(lat); // km
816     m_pos.m_w = sqrt(sqr(m_pos.m_x) +
817     sqr(m_pos.m_y) +
818     sqr(m_pos.m_z)); // range, km
819    
820     m_vel.m_x = -mfactor * m_pos.m_y; // km / sec
821     m_vel.m_y = mfactor * m_pos.m_x;
822     m_vel.m_z = 0.0;
823     m_vel.m_w = sqrt(sqr(m_vel.m_x) + // range rate km/sec^2
824     sqr(m_vel.m_y));
825     }
826    
827     //////////////////////////////////////////////////////////////////////////////
828     // toGeo()
829     // Return the corresponding geodetic position (based on the current ECI
830     // coordinates/Julian date).
831     // Assumes the earth is an oblate spheroid as defined in WGS '72.
832     // Side effects: Converts the position and velocity vectors to km-based units.
833     // Reference: The 1992 Astronomical Almanac, page K12.
834     // Reference: www.celestrak.com (Dr. TS Kelso)
835     cCoordGeo cEci::toGeo()
836     {
837     ae2km(); // Vectors must be in kilometer-based units
838    
839     double theta = AcTan(m_pos.m_y, m_pos.m_x);
840     double lon = fmod(theta - m_date.toGMST(), TWOPI);
841    
842     if (lon < 0.0)
843     lon += TWOPI; // "wrap" negative modulo
844    
845     double r = sqrt(sqr(m_pos.m_x) + sqr(m_pos.m_y));
846     double e2 = F * (2.0 - F);
847     double lat = AcTan(m_pos.m_z, r);
848    
849     const double delta = 1.0e-07;
850     double phi;
851     double c;
852    
853     do
854     {
855     phi = lat;
856     c = 1.0 / sqrt(1.0 - e2 * sqr(sin(phi)));
857     lat = AcTan(m_pos.m_z + XKMPER_WGS72 * c * e2 * sin(phi), r);
858     }
859     while (fabs(lat - phi) > delta);
860    
861     double alt = r / cos(lat) - XKMPER_WGS72 * c;
862    
863     return cCoordGeo(lat, lon, alt); // radians, radians, kilometers
864     }
865    
866     //////////////////////////////////////////////////////////////////////////////
867     // ae2km()
868     // Convert the position and velocity vector units from AE-based units
869     // to kilometer based units.
870     void cEci::ae2km()
871     {
872     if (UnitsAreAe())
873     {
874     MulPos(XKMPER_WGS72 / AE); // km
875     MulVel((XKMPER_WGS72 / AE) * (MIN_PER_DAY / 86400)); // km/sec
876     m_VecUnits = UNITS_KM;
877     }
878     }
879     //
880     // cNoradBase.cpp
881     //
882     // Historical Note:
883     // The equations used here (and in derived classes) to determine satellite
884     // ECI coordinates/velocity come from the December, 1980 NORAD document
885     // "Space Track Report No. 3". The report details 6 orbital models and
886     // provides FORTRAN IV implementations of each. The classes here
887     // implement only two of the orbital models: SGP4 and SDP4. These two models,
888     // one for "near-earth" objects and one for "deep space" objects, are widely
889     // used in satellite tracking software and can produce very accurate results
890     // when used with current NORAD two-line element datum.
891     //
892     // The NORAD FORTRAN IV SGP4/SDP4 implementations were converted to Pascal by
893     // Dr. TS Kelso in 1995. In 1996 these routines were ported in a straight-
894     // forward manner to C++ by Varol Okan. The SGP4/SDP4 classes here were
895     // written by Michael F. Henry in 2002-03 and are a modern C++ re-write of
896     // the work done by Okan. In addition to introducing an object-oriented
897     // architecture, the last residues of the original FORTRAN code (such as
898     // labels and gotos) were eradicated.
899     //
900     // For excellent information on the underlying physics of orbits, visible
901     // satellite observations, current NORAD TLE data, and other related material,
902     // see http://www.celestrak.com which is maintained by Dr. TS Kelso.
903     //
904     // Copyright (c) 2003 Michael F. Henry
905     //
906     // mfh 12/07/2003
907     //
908     //////////////////////////////////////////////////////////////////////////////
909     cNoradBase::cNoradBase(const cOrbit &orbit) :
910     m_Orbit(orbit)
911     {
912     Initialize();
913     }
914    
915     cNoradBase& cNoradBase::operator=(const cNoradBase &b)
916     {
917     // m_Orbit is a "const" member var, so cast away its
918     // "const-ness" in order to complete the assigment.
919     *(const_cast<cOrbit*>(&m_Orbit)) = b.m_Orbit;
920    
921     return *this;
922     }
923    
924     //////////////////////////////////////////////////////////////////////////////
925     // Initialize()
926     // Perform the initialization of member variables, specifically the variables
927     // used by derived-class objects to calculate ECI coordinates.
928     void cNoradBase::Initialize()
929     {
930     // Initialize any variables which are time-independent when
931     // calculating the ECI coordinates of the satellite.
932     m_satInc = m_Orbit.Inclination();
933     m_satEcc = m_Orbit.Eccentricity();
934    
935     m_cosio = cos(m_satInc);
936     m_theta2 = m_cosio * m_cosio;
937     m_x3thm1 = 3.0 * m_theta2 - 1.0;
938     m_eosq = m_satEcc * m_satEcc;
939     m_betao2 = 1.0 - m_eosq;
940     m_betao = sqrt(m_betao2);
941    
942     // The "recovered" semi-minor axis and mean motion.
943     m_aodp = m_Orbit.SemiMinor();
944     m_xnodp = m_Orbit.mnMotionRec();
945    
946     // For perigee below 156 km, the values of S and QOMS2T are altered.
947     m_perigee = XKMPER_WGS72 * (m_aodp * (1.0 - m_satEcc) - AE);
948    
949     m_s4 = S;
950     m_qoms24 = QOMS2T;
951    
952     if (m_perigee < 156.0)
953     {
954     m_s4 = m_perigee - 78.0;
955    
956     if (m_perigee <= 98.0)
957     {
958     m_s4 = 20.0;
959     }
960    
961     m_qoms24 = pow((120.0 - m_s4) * AE / XKMPER_WGS72, 4.0);
962     m_s4 = m_s4 / XKMPER_WGS72 + AE;
963     }
964    
965     const double pinvsq = 1.0 / (m_aodp * m_aodp * m_betao2 * m_betao2);
966    
967     m_tsi = 1.0 / (m_aodp - m_s4);
968     m_eta = m_aodp * m_satEcc * m_tsi;
969     m_etasq = m_eta * m_eta;
970     m_eeta = m_satEcc * m_eta;
971    
972     const double psisq = fabs(1.0 - m_etasq);
973    
974     m_coef = m_qoms24 * pow(m_tsi,4.0);
975     m_coef1 = m_coef / pow(psisq,3.5);
976    
977     const double c2 = m_coef1 * m_xnodp *
978     (m_aodp * (1.0 + 1.5 * m_etasq + m_eeta * (4.0 + m_etasq)) +
979     0.75 * CK2 * m_tsi / psisq * m_x3thm1 *
980     (8.0 + 3.0 * m_etasq * (8.0 + m_etasq)));
981    
982     m_c1 = m_Orbit.BStar() * c2;
983     m_sinio = sin(m_satInc);
984    
985     const double a3ovk2 = -XJ3 / CK2 * pow(AE,3.0);
986    
987     m_c3 = m_coef * m_tsi * a3ovk2 * m_xnodp * AE * m_sinio / m_satEcc;
988     m_x1mth2 = 1.0 - m_theta2;
989     m_c4 = 2.0 * m_xnodp * m_coef1 * m_aodp * m_betao2 *
990     (m_eta * (2.0 + 0.5 * m_etasq) +
991     m_satEcc * (0.5 + 2.0 * m_etasq) -
992     2.0 * CK2 * m_tsi / (m_aodp * psisq) *
993     (-3.0 * m_x3thm1 * (1.0 - 2.0 * m_eeta + m_etasq * (1.5 - 0.5 * m_eeta)) +
994     0.75 * m_x1mth2 *
995     (2.0 * m_etasq - m_eeta * (1.0 + m_etasq)) *
996     cos(2.0 * m_Orbit.ArgPerigee())));
997    
998     const double theta4 = m_theta2 * m_theta2;
999     const double temp1 = 3.0 * CK2 * pinvsq * m_xnodp;
1000     const double temp2 = temp1 * CK2 * pinvsq;
1001     const double temp3 = 1.25 * CK4 * pinvsq * pinvsq * m_xnodp;
1002    
1003     m_xmdot = m_xnodp + 0.5 * temp1 * m_betao * m_x3thm1 +
1004     0.0625 * temp2 * m_betao *
1005     (13.0 - 78.0 * m_theta2 + 137.0 * theta4);
1006    
1007     const double x1m5th = 1.0 - 5.0 * m_theta2;
1008    
1009     m_omgdot = -0.5 * temp1 * x1m5th + 0.0625 * temp2 *
1010     (7.0 - 114.0 * m_theta2 + 395.0 * theta4) +
1011     temp3 * (3.0 - 36.0 * m_theta2 + 49.0 * theta4);
1012    
1013     const double xhdot1 = -temp1 * m_cosio;
1014    
1015     m_xnodot = xhdot1 + (0.5 * temp2 * (4.0 - 19.0 * m_theta2) +
1016     2.0 * temp3 * (3.0 - 7.0 * m_theta2)) * m_cosio;
1017     m_xnodcf = 3.5 * m_betao2 * xhdot1 * m_c1;
1018     m_t2cof = 1.5 * m_c1;
1019     m_xlcof = 0.125 * a3ovk2 * m_sinio *
1020     (3.0 + 5.0 * m_cosio) / (1.0 + m_cosio);
1021     m_aycof = 0.25 * a3ovk2 * m_sinio;
1022     m_x7thm1 = 7.0 * m_theta2 - 1.0;
1023     }
1024    
1025     //////////////////////////////////////////////////////////////////////////////
1026     bool cNoradBase::FinalPosition(double incl, double omega,
1027     double e, double a,
1028     double xl, double xnode,
1029     double xn, double tsince,
1030     cEci &eci)
1031     {
1032     if ((e * e) > 1.0)
1033     {
1034     // error in satellite data
1035     return false;
1036     }
1037    
1038     double beta = sqrt(1.0 - e * e);
1039    
1040     // Long period periodics
1041     double axn = e * cos(omega);
1042     double temp = 1.0 / (a * beta * beta);
1043     double xll = temp * m_xlcof * axn;
1044     double aynl = temp * m_aycof;
1045     double xlt = xl + xll;
1046     double ayn = e * sin(omega) + aynl;
1047    
1048     // Solve Kepler's Equation
1049    
1050     double capu = Fmod2p(xlt - xnode);
1051     double temp2 = capu;
1052     double temp3 = 0.0;
1053     double temp4 = 0.0;
1054     double temp5 = 0.0;
1055     double temp6 = 0.0;
1056     double sinepw = 0.0;
1057     double cosepw = 0.0;
1058     bool fDone = false;
1059    
1060     for (int i = 1; (i <= 10) && !fDone; i++)
1061     {
1062     sinepw = sin(temp2);
1063     cosepw = cos(temp2);
1064     temp3 = axn * sinepw;
1065     temp4 = ayn * cosepw;
1066     temp5 = axn * cosepw;
1067     temp6 = ayn * sinepw;
1068    
1069     double epw = (capu - temp4 + temp3 - temp2) /
1070     (1.0 - temp5 - temp6) + temp2;
1071    
1072     if (fabs(epw - temp2) <= E6A)
1073     fDone = true;
1074     else
1075     temp2 = epw;
1076     }
1077    
1078     // Short period preliminary quantities
1079     double ecose = temp5 + temp6;
1080     double esine = temp3 - temp4;
1081     double elsq = axn * axn + ayn * ayn;
1082     temp = 1.0 - elsq;
1083     double pl = a * temp;
1084     double r = a * (1.0 - ecose);
1085     double temp1 = 1.0 / r;
1086     double rdot = XKE * sqrt(a) * esine * temp1;
1087     double rfdot = XKE * sqrt(pl) * temp1;
1088     temp2 = a * temp1;
1089     double betal = sqrt(temp);
1090     temp3 = 1.0 / (1.0 + betal);
1091     double cosu = temp2 * (cosepw - axn + ayn * esine * temp3);
1092     double sinu = temp2 * (sinepw - ayn - axn * esine * temp3);
1093     double u = AcTan(sinu, cosu);
1094     double sin2u = 2.0 * sinu * cosu;
1095     double cos2u = 2.0 * cosu * cosu - 1.0;
1096    
1097     temp = 1.0 / pl;
1098     temp1 = CK2 * temp;
1099     temp2 = temp1 * temp;
1100    
1101     // Update for short periodics
1102     double rk = r * (1.0 - 1.5 * temp2 * betal * m_x3thm1) +
1103     0.5 * temp1 * m_x1mth2 * cos2u;
1104     double uk = u - 0.25 * temp2 * m_x7thm1 * sin2u;
1105     double xnodek = xnode + 1.5 * temp2 * m_cosio * sin2u;
1106     double xinck = incl + 1.5 * temp2 * m_cosio * m_sinio * cos2u;
1107     double rdotk = rdot - xn * temp1 * m_x1mth2 * sin2u;
1108     double rfdotk = rfdot + xn * temp1 * (m_x1mth2 * cos2u + 1.5 * m_x3thm1);
1109    
1110     // Orientation vectors
1111     double sinuk = sin(uk);
1112     double cosuk = cos(uk);
1113     double sinik = sin(xinck);
1114     double cosik = cos(xinck);
1115     double sinnok = sin(xnodek);
1116     double cosnok = cos(xnodek);
1117     double xmx = -sinnok * cosik;
1118     double xmy = cosnok * cosik;
1119     double ux = xmx * sinuk + cosnok * cosuk;
1120     double uy = xmy * sinuk + sinnok * cosuk;
1121     double uz = sinik * sinuk;
1122     double vx = xmx * cosuk - cosnok * sinuk;
1123     double vy = xmy * cosuk - sinnok * sinuk;
1124     double vz = sinik * cosuk;
1125    
1126     // Position
1127     double x = rk * ux;
1128     double y = rk * uy;
1129     double z = rk * uz;
1130    
1131     cVector vecPos(x, y, z);
1132    
1133     // Validate on altitude
1134     double altKm = (vecPos.Magnitude() * (XKMPER_WGS72 / AE));
1135    
1136     if ((altKm < XKMPER_WGS72) || (altKm > (2 * GEOSYNC_ALT)))
1137     return false;
1138    
1139     // Velocity
1140     double xdot = rdotk * ux + rfdotk * vx;
1141     double ydot = rdotk * uy + rfdotk * vy;
1142     double zdot = rdotk * uz + rfdotk * vz;
1143    
1144     cVector vecVel(xdot, ydot, zdot);
1145    
1146     cJulian gmt = m_Orbit.Epoch();
1147     gmt.addMin(tsince);
1148    
1149     eci = cEci(vecPos, vecVel, gmt);
1150    
1151     return true;
1152     }
1153    
1154     //
1155     // cNoradSGP4.cpp
1156     //
1157     // NORAD SGP4 implementation. See historical note in cNoradBase.cpp
1158     // Copyright (c) 2003 Michael F. Henry
1159     //
1160     // mfh 12/07/2003
1161     //
1162     //////////////////////////////////////////////////////////////////////////////
1163     cNoradSGP4::cNoradSGP4(const cOrbit &orbit) :
1164     cNoradBase(orbit)
1165     {
1166     m_c5 = 2.0 * m_coef1 * m_aodp * m_betao2 *
1167     (1.0 + 2.75 * (m_etasq + m_eeta) + m_eeta * m_etasq);
1168     m_omgcof = m_Orbit.BStar() * m_c3 * cos(m_Orbit.ArgPerigee());
1169     m_xmcof = -TWOTHRD * m_coef * m_Orbit.BStar() * AE / m_eeta;
1170     m_delmo = pow(1.0 + m_eta * cos(m_Orbit.mnAnomaly()), 3.0);
1171     m_sinmo = sin(m_Orbit.mnAnomaly());
1172     }
1173    
1174    
1175     //////////////////////////////////////////////////////////////////////////////
1176     // getPosition()
1177     // This procedure returns the ECI position and velocity for the satellite
1178     // in the orbit at the given number of minutes since the TLE epoch time
1179     // using the NORAD Simplified General Perturbation 4, near earth orbit
1180     // model.
1181     //
1182     // tsince - Time in minutes since the TLE epoch (GMT).
1183     // eci - ECI object to hold position information.
1184     // To convert the returned ECI position vector to km,
1185     // multiply each component by:
1186     // (XKMPER_WGS72 / AE).
1187     // To convert the returned ECI velocity vector to km/sec,
1188     // multiply each component by:
1189     // (XKMPER_WGS72 / AE) * (MIN_PER_DAY / 86400).
1190    
1191     bool cNoradSGP4::getPosition(double tsince, cEci &eci)
1192     {
1193     // For m_perigee less than 220 kilometers, the isimp flag is set and
1194     // the equations are truncated to linear variation in sqrt a and
1195     // quadratic variation in mean anomaly. Also, the m_c3 term, the
1196     // delta omega term, and the delta m term are dropped.
1197     bool isimp = false;
1198     if ((m_aodp * (1.0 - m_satEcc) / AE) < (220.0 / XKMPER_WGS72 + AE))
1199     {
1200     isimp = true;
1201     }
1202    
1203     double d2 = 0.0;
1204     double d3 = 0.0;
1205     double d4 = 0.0;
1206    
1207     double t3cof = 0.0;
1208     double t4cof = 0.0;
1209     double t5cof = 0.0;
1210    
1211     if (!isimp)
1212     {
1213     double c1sq = m_c1 * m_c1;
1214    
1215     d2 = 4.0 * m_aodp * m_tsi * c1sq;
1216    
1217     double temp = d2 * m_tsi * m_c1 / 3.0;
1218    
1219     d3 = (17.0 * m_aodp + m_s4) * temp;
1220     d4 = 0.5 * temp * m_aodp * m_tsi *
1221     (221.0 * m_aodp + 31.0 * m_s4) * m_c1;
1222     t3cof = d2 + 2.0 * c1sq;
1223     t4cof = 0.25 * (3.0 * d3 + m_c1 * (12.0 * d2 + 10.0 * c1sq));
1224     t5cof = 0.2 * (3.0 * d4 + 12.0 * m_c1 * d3 + 6.0 *
1225     d2 * d2 + 15.0 * c1sq * (2.0 * d2 + c1sq));
1226     }
1227    
1228     // Update for secular gravity and atmospheric drag.
1229     double xmdf = m_Orbit.mnAnomaly() + m_xmdot * tsince;
1230     double omgadf = m_Orbit.ArgPerigee() + m_omgdot * tsince;
1231     double xnoddf = m_Orbit.RAAN() + m_xnodot * tsince;
1232     double omega = omgadf;
1233     double xmp = xmdf;
1234     double tsq = tsince * tsince;
1235     double xnode = xnoddf + m_xnodcf * tsq;
1236     double tempa = 1.0 - m_c1 * tsince;
1237     double tempe = m_Orbit.BStar() * m_c4 * tsince;
1238     double templ = m_t2cof * tsq;
1239    
1240     if (!isimp)
1241     {
1242     double delomg = m_omgcof * tsince;
1243     double delm = m_xmcof * (pow(1.0 + m_eta * cos(xmdf), 3.0) - m_delmo);
1244     double temp = delomg + delm;
1245    
1246     xmp = xmdf + temp;
1247     omega = omgadf - temp;
1248    
1249     double tcube = tsq * tsince;
1250     double tfour = tsince * tcube;
1251    
1252     tempa = tempa - d2 * tsq - d3 * tcube - d4 * tfour;
1253     tempe = tempe + m_Orbit.BStar() * m_c5 * (sin(xmp) - m_sinmo);
1254     templ = templ + t3cof * tcube + tfour * (t4cof + tsince * t5cof);
1255     }
1256    
1257     double a = m_aodp * sqr(tempa);
1258     double e = m_satEcc - tempe;
1259    
1260    
1261     double xl = xmp + omega + xnode + m_xnodp * templ;
1262     double xn = XKE / pow(a, 1.5);
1263    
1264     return FinalPosition(m_satInc, omgadf, e, a, xl, xnode, xn, tsince, eci);
1265     }
1266    
1267     //
1268     // cNoradSDP4.cpp
1269     //
1270     // NORAD SDP4 implementation. See historical note in cNoradBase.cpp
1271     // Copyright (c) 2003 Michael F. Henry
1272     //
1273     // mfh 12/07/2003
1274     //
1275    
1276     const double zns = 1.19459E-5; const double c1ss = 2.9864797E-6;
1277     const double zes = 0.01675; const double znl = 1.5835218E-4;
1278     const double c1l = 4.7968065E-7; const double zel = 0.05490;
1279     const double zcosis = 0.91744867; const double zsinis = 0.39785416;
1280     const double zsings = -0.98088458; const double zcosgs = 0.1945905;
1281     const double q22 = 1.7891679E-6; const double q31 = 2.1460748E-6;
1282     const double q33 = 2.2123015E-7; const double g22 = 5.7686396;
1283     const double g32 = 0.95240898; const double g44 = 1.8014998;
1284     const double g52 = 1.0508330; const double g54 = 4.4108898;
1285     const double root22 = 1.7891679E-6; const double root32 = 3.7393792E-7;
1286     const double root44 = 7.3636953E-9; const double root52 = 1.1428639E-7;
1287     const double root54 = 2.1765803E-9; const double thdt = 4.3752691E-3;
1288    
1289     //////////////////////////////////////////////////////////////////////////////
1290     cNoradSDP4::cNoradSDP4(const cOrbit &orbit) :
1291     cNoradBase(orbit)
1292     {
1293     m_sing = sin(m_Orbit.ArgPerigee());
1294     m_cosg = cos(m_Orbit.ArgPerigee());
1295    
1296     dp_savtsn = 0.0;
1297     dp_zmos = 0.0;
1298     dp_se2 = 0.0;
1299     dp_se3 = 0.0;
1300     dp_si2 = 0.0;
1301     dp_si3 = 0.0;
1302     dp_sl2 = 0.0;
1303     dp_sl3 = 0.0;
1304     dp_sl4 = 0.0;
1305     dp_sghs = 0.0;
1306     dp_sgh2 = 0.0;
1307     dp_sgh3 = 0.0;
1308     dp_sgh4 = 0.0;
1309     dp_sh2 = 0.0;
1310     dp_sh3 = 0.0;
1311     dp_zmol = 0.0;
1312     dp_ee2 = 0.0;
1313     dp_e3 = 0.0;
1314     dp_xi2 = 0.0;
1315     dp_xi3 = 0.0;
1316     dp_xl2 = 0.0;
1317     dp_xl3 = 0.0;
1318     dp_xl4 = 0.0;
1319     dp_xgh2 = 0.0;
1320     dp_xgh3 = 0.0;
1321     dp_xgh4 = 0.0;
1322     dp_xh2 = 0.0;
1323     dp_xh3 = 0.0;
1324     dp_xqncl = 0.0;
1325     dp_thgr = 0.0;
1326     dp_omegaq = 0.0;
1327     dp_sse = 0.0;
1328     dp_ssi = 0.0;
1329     dp_ssl = 0.0;
1330     dp_ssh = 0.0;
1331     dp_ssg = 0.0;
1332     dp_d2201 = 0.0;
1333     dp_d2211 = 0.0;
1334     dp_d3210 = 0.0;
1335     dp_d3222 = 0.0;
1336     dp_d4410 = 0.0;
1337     dp_d4422 = 0.0;
1338     dp_d5220 = 0.0;
1339     dp_d5232 = 0.0;
1340     dp_d5421 = 0.0;
1341     dp_d5433 = 0.0;
1342     dp_xlamo = 0.0;
1343     dp_del1 = 0.0;
1344     dp_del2 = 0.0;
1345     dp_del3 = 0.0;
1346     dp_fasx2 = 0.0;
1347     dp_fasx4 = 0.0;
1348     dp_fasx6 = 0.0;
1349     dp_xfact = 0.0;
1350     dp_xli = 0.0;
1351     dp_xni = 0.0;
1352     dp_atime = 0.0;
1353     dp_stepp = 0.0;
1354     dp_stepn = 0.0;
1355     dp_step2 = 0.0;
1356    
1357     dp_iresfl = false;
1358     dp_isynfl = false;
1359    
1360     }
1361    
1362    
1363     /////////////////////////////////////////////////////////////////////////////
1364     bool cNoradSDP4::DeepInit(double *eosq, double *sinio, double *cosio,
1365     double *betao, double *aodp, double *theta2,
1366     double *sing, double *cosg, double *betao2,
1367     double *xmdot, double *omgdot, double *xnodott)
1368     {
1369     eqsq = *eosq;
1370     siniq = *sinio;
1371     cosiq = *cosio;
1372     rteqsq = *betao;
1373     ao = *aodp;
1374     cosq2 = *theta2;
1375     sinomo = *sing;
1376     cosomo = *cosg;
1377     bsq = *betao2;
1378     xlldot = *xmdot;
1379     omgdt = *omgdot;
1380     xnodot = *xnodott;
1381    
1382     // Deep space initialization
1383     cJulian jd = m_Orbit.Epoch();
1384    
1385     dp_thgr = jd.toGMST();
1386    
1387     double eq = m_Orbit.Eccentricity();
1388     double aqnv = 1.0 / ao;
1389    
1390     dp_xqncl = m_Orbit.Inclination();
1391    
1392     double xmao = m_Orbit.mnAnomaly();
1393     double xpidot = omgdt + xnodot;
1394     double sinq = sin(m_Orbit.RAAN());
1395     double cosq = cos(m_Orbit.RAAN());
1396    
1397     dp_omegaq = m_Orbit.ArgPerigee();
1398    
1399     // Initialize lunar solar terms
1400     double day = jd.FromJan1_12h_1900();
1401    
1402     if (day != dpi_day)
1403     {
1404     dpi_day = day;
1405     dpi_xnodce = 4.5236020 - 9.2422029E-4 * day;
1406     dpi_stem = sin(dpi_xnodce);
1407     dpi_ctem = cos(dpi_xnodce);
1408     dpi_zcosil = 0.91375164 - 0.03568096 * dpi_ctem;
1409     dpi_zsinil = sqrt(1.0 - dpi_zcosil * dpi_zcosil);
1410     dpi_zsinhl = 0.089683511 *dpi_stem / dpi_zsinil;
1411     dpi_zcoshl = sqrt(1.0 - dpi_zsinhl * dpi_zsinhl);
1412     dpi_c = 4.7199672 + 0.22997150 * day;
1413     dpi_gam = 5.8351514 + 0.0019443680 * day;
1414     dp_zmol = Fmod2p(dpi_c - dpi_gam);
1415     dpi_zx = 0.39785416 * dpi_stem / dpi_zsinil;
1416     dpi_zy = dpi_zcoshl * dpi_ctem + 0.91744867 * dpi_zsinhl * dpi_stem;
1417     dpi_zx = AcTan(dpi_zx,dpi_zy) + dpi_gam - dpi_xnodce;
1418     dpi_zcosgl = cos(dpi_zx);
1419     dpi_zsingl = sin(dpi_zx);
1420     dp_zmos = 6.2565837 + 0.017201977 * day;
1421     dp_zmos = Fmod2p(dp_zmos);
1422     }
1423    
1424     dp_savtsn = 1.0e20;
1425    
1426     double zcosg = zcosgs;
1427     double zsing = zsings;
1428     double zcosi = zcosis;
1429     double zsini = zsinis;
1430     double zcosh = cosq;
1431     double zsinh = sinq;
1432     double cc = c1ss;
1433     double zn = zns;
1434     double ze = zes;
1435     double zmo = dp_zmos;
1436     double xnoi = 1.0 / m_xnodp;
1437    
1438     double a1; double a3; double a7; double a8; double a9; double a10;
1439     double a2; double a4; double a5; double a6; double x1; double x2;
1440     double x3; double x4; double x5; double x6; double x7; double x8;
1441     double z31; double z32; double z33; double z1; double z2; double z3;
1442     double z11; double z12; double z13; double z21; double z22; double z23;
1443     double s3; double s2; double s4; double s1; double s5; double s6;
1444     double s7;
1445     double se = 0.0; double si = 0.0; double sl = 0.0;
1446     double sgh = 0.0; double sh = 0.0;
1447    
1448     // Apply the solar and lunar terms on the first pass, then re-apply the
1449     // solar terms again on the second pass.
1450    
1451     for (int pass = 1; pass <= 2; pass++)
1452     {
1453     // Do solar terms
1454     a1 = zcosg * zcosh + zsing * zcosi * zsinh;
1455     a3 = -zsing * zcosh + zcosg * zcosi * zsinh;
1456     a7 = -zcosg * zsinh + zsing * zcosi * zcosh;
1457     a8 = zsing * zsini;
1458     a9 = zsing * zsinh + zcosg * zcosi * zcosh;
1459     a10 = zcosg * zsini;
1460     a2 = cosiq * a7 + siniq * a8;
1461     a4 = cosiq * a9 + siniq * a10;
1462     a5 = -siniq * a7 + cosiq * a8;
1463     a6 = -siniq * a9 + cosiq * a10;
1464     x1 = a1 * cosomo + a2 * sinomo;
1465     x2 = a3 * cosomo + a4 * sinomo;
1466     x3 = -a1 * sinomo + a2 * cosomo;
1467     x4 = -a3 * sinomo + a4 * cosomo;
1468     x5 = a5 * sinomo;
1469     x6 = a6 * sinomo;
1470     x7 = a5 * cosomo;
1471     x8 = a6 * cosomo;
1472     z31 = 12.0 * x1 * x1 - 3.0 * x3 * x3;
1473     z32 = 24.0 * x1 * x2 - 6.0 * x3 * x4;
1474     z33 = 12.0 * x2 * x2 - 3.0 * x4 * x4;
1475     z1 = 3.0 * (a1 * a1 + a2 * a2) + z31 * eqsq;
1476     z2 = 6.0 * (a1 * a3 + a2 * a4) + z32 * eqsq;
1477     z3 = 3.0 * (a3 * a3 + a4 * a4) + z33 * eqsq;
1478     z11 = -6.0 * a1 * a5 + eqsq*(-24.0 * x1 * x7 - 6.0 * x3 * x5);
1479     z12 = -6.0 * (a1 * a6 + a3 * a5) +
1480     eqsq * (-24.0 * (x2 * x7 + x1 * x8) - 6.0 * (x3 * x6 + x4 * x5));
1481     z13 = -6.0 * a3 * a6 + eqsq * (-24.0 * x2 * x8 - 6.0 * x4 * x6);
1482     z21 = 6.0 * a2 * a5 + eqsq * (24.0 * x1 * x5 - 6.0 * x3 * x7);
1483     z22 = 6.0*(a4 * a5 + a2 * a6) +
1484     eqsq * (24.0 * (x2 * x5 + x1 * x6) - 6.0 * (x4 * x7 + x3 * x8));
1485     z23 = 6.0 * a4 * a6 + eqsq*(24.0 * x2 * x6 - 6.0 * x4 * x8);
1486     z1 = z1 + z1 + bsq * z31;
1487     z2 = z2 + z2 + bsq * z32;
1488     z3 = z3 + z3 + bsq * z33;
1489     s3 = cc * xnoi;
1490     s2 = -0.5 * s3/rteqsq;
1491     s4 = s3 * rteqsq;
1492     s1 = -15.0 * eq * s4;
1493     s5 = x1 * x3 + x2 * x4;
1494     s6 = x2 * x3 + x1 * x4;
1495     s7 = x2 * x4 - x1 * x3;
1496     se = s1 * zn * s5;
1497     si = s2 * zn * (z11 + z13);
1498     sl = -zn * s3 * (z1 + z3 - 14.0 - 6.0 * eqsq);
1499     sgh = s4 * zn * (z31 + z33 - 6.0);
1500     sh = -zn * s2 * (z21 + z23);
1501    
1502     if (dp_xqncl < 5.2359877E-2)
1503     sh = 0.0;
1504    
1505     dp_ee2 = 2.0 * s1 * s6;
1506     dp_e3 = 2.0 * s1 * s7;
1507     dp_xi2 = 2.0 * s2 * z12;
1508     dp_xi3 = 2.0 * s2 * (z13 - z11);
1509     dp_xl2 = -2.0 * s3 * z2;
1510     dp_xl3 = -2.0 * s3 * (z3 - z1);
1511     dp_xl4 = -2.0 * s3 * (-21.0 - 9.0 * eqsq) * ze;
1512     dp_xgh2 = 2.0 * s4 * z32;
1513     dp_xgh3 = 2.0 * s4 * (z33 - z31);
1514     dp_xgh4 = -18.0 * s4 * ze;
1515     dp_xh2 = -2.0 * s2 * z22;
1516     dp_xh3 = -2.0 * s2 * (z23 - z21);
1517    
1518     if (pass == 1)
1519     {
1520     // Do lunar terms
1521     dp_sse = se;
1522     dp_ssi = si;
1523     dp_ssl = sl;
1524     dp_ssh = sh / siniq;
1525     dp_ssg = sgh - cosiq * dp_ssh;
1526     dp_se2 = dp_ee2;
1527     dp_si2 = dp_xi2;
1528     dp_sl2 = dp_xl2;
1529     dp_sgh2 = dp_xgh2;
1530     dp_sh2 = dp_xh2;
1531     dp_se3 = dp_e3;
1532     dp_si3 = dp_xi3;
1533     dp_sl3 = dp_xl3;
1534     dp_sgh3 = dp_xgh3;
1535     dp_sh3 = dp_xh3;
1536     dp_sl4 = dp_xl4;
1537     dp_sgh4 = dp_xgh4;
1538     zcosg = dpi_zcosgl;
1539     zsing = dpi_zsingl;
1540     zcosi = dpi_zcosil;
1541     zsini = dpi_zsinil;
1542     zcosh = dpi_zcoshl * cosq + dpi_zsinhl * sinq;
1543     zsinh = sinq * dpi_zcoshl - cosq * dpi_zsinhl;
1544     zn = znl;
1545     cc = c1l;
1546     ze = zel;
1547     zmo = dp_zmol;
1548     }
1549     }
1550    
1551     dp_sse = dp_sse + se;
1552     dp_ssi = dp_ssi + si;
1553     dp_ssl = dp_ssl + sl;
1554     dp_ssg = dp_ssg + sgh - cosiq / siniq * sh;
1555     dp_ssh = dp_ssh + sh / siniq;
1556    
1557     // Geopotential resonance initialization for 12 hour orbits
1558     dp_iresfl = false;
1559     dp_isynfl = false;
1560    
1561     bool bInitOnExit = true;
1562     double g310;
1563     double f220;
1564     double bfact = 0.0;
1565    
1566     if ((m_xnodp >= 0.0052359877) || (m_xnodp <= 0.0034906585))
1567     {
1568     if ((m_xnodp < 8.26E-3) || (m_xnodp > 9.24E-3) || (eq < 0.5))
1569     {
1570     bInitOnExit = false;
1571     }
1572     else
1573     {
1574     dp_iresfl = true;
1575    
1576     double eoc = eq * eqsq;
1577     double g201 = -0.306 - (eq - 0.64) * 0.440;
1578    
1579     double g211; double g322;
1580    
1581     double g410; double g422;
1582     double g520;
1583    
1584     if (eq <= 0.65)
1585     {
1586     g211 = 3.616 - 13.247 * eq + 16.290 * eqsq;
1587     g310 = -19.302 + 117.390 * eq - 228.419 * eqsq + 156.591 * eoc;
1588     g322 = -18.9068 + 109.7927 * eq - 214.6334 * eqsq + 146.5816 * eoc;
1589     g410 = -41.122 + 242.694 * eq - 471.094 * eqsq + 313.953 * eoc;
1590     g422 = -146.407 + 841.880 * eq - 1629.014 * eqsq + 1083.435 * eoc;
1591     g520 = -532.114 + 3017.977 * eq - 5740.0 * eqsq + 3708.276 * eoc;
1592     }
1593     else
1594     {
1595     g211 = -72.099 + 331.819 * eq - 508.738 * eqsq + 266.724 * eoc;
1596     g310 = -346.844 + 1582.851 * eq - 2415.925 * eqsq + 1246.113 * eoc;
1597     g322 = -342.585 + 1554.908 * eq - 2366.899 * eqsq + 1215.972 * eoc;
1598     g410 = -1052.797 + 4758.686 * eq - 7193.992 * eqsq + 3651.957 * eoc;
1599     g422 = -3581.69 + 16178.11 * eq - 24462.77 * eqsq + 12422.52 * eoc;
1600    
1601     if (eq <= 0.715)
1602     g520 = 1464.74 - 4664.75 * eq + 3763.64 * eqsq;
1603     else
1604     g520 = -5149.66 + 29936.92 * eq - 54087.36 * eqsq + 31324.56 * eoc;
1605     }
1606    
1607     double g533;
1608     double g521;
1609     double g532;
1610    
1611     if (eq < 0.7)
1612     {
1613     g533 = -919.2277 + 4988.61 * eq - 9064.77 * eqsq + 5542.21 * eoc;
1614     g521 = -822.71072 + 4568.6173 * eq - 8491.4146 * eqsq + 5337.524 * eoc;
1615     g532 = -853.666 + 4690.25 * eq - 8624.77 * eqsq + 5341.4 * eoc;
1616     }
1617     else
1618     {
1619     g533 = -37995.78 + 161616.52 * eq - 229838.2 * eqsq + 109377.94 * eoc;
1620     g521 = -51752.104 + 218913.95 * eq - 309468.16 * eqsq + 146349.42 * eoc;
1621     g532 = -40023.88 + 170470.89 * eq - 242699.48 * eqsq + 115605.82 * eoc;
1622     }
1623    
1624     double sini2 = siniq * siniq;
1625     f220 = 0.75*(1.0 + 2.0 * cosiq + cosq2);
1626     double f221 = 1.5 * sini2;
1627     double f321 = 1.875 * siniq*(1.0 - 2.0 * cosiq - 3.0 * cosq2);
1628     double f322 = -1.875 * siniq*(1.0 + 2.0 * cosiq - 3.0 * cosq2);
1629     double f441 = 35.0 * sini2 * f220;
1630     double f442 = 39.3750 * sini2 * sini2;
1631     double f522 = 9.84375 * siniq*(sini2*(1.0 - 2.0 * cosiq - 5.0 * cosq2) +
1632     0.33333333*(-2.0 + 4.0 * cosiq + 6.0 * cosq2));
1633     double f523 = siniq*(4.92187512 * sini2*(-2.0 - 4.0 * cosiq + 10.0 * cosq2) +
1634     6.56250012 * (1.0 + 2.0 * cosiq - 3.0 * cosq2));
1635     double f542 = 29.53125 * siniq*(2.0 - 8.0 * cosiq + cosq2 * (-12.0 + 8.0 * cosiq + 10.0 * cosq2));
1636     double f543 = 29.53125 * siniq*(-2.0 - 8.0 * cosiq + cosq2 * (12.0 + 8.0 * cosiq - 10.0 * cosq2));
1637     double xno2 = m_xnodp * m_xnodp;
1638     double ainv2 = aqnv * aqnv;
1639     double temp1 = 3.0 * xno2 * ainv2;
1640     double temp = temp1 * root22;
1641    
1642     dp_d2201 = temp * f220 * g201;
1643     dp_d2211 = temp * f221 * g211;
1644     temp1 = temp1 * aqnv;
1645     temp = temp1 * root32;
1646     dp_d3210 = temp * f321 * g310;
1647     dp_d3222 = temp * f322 * g322;
1648     temp1 = temp1 * aqnv;
1649     temp = 2.0 * temp1 * root44;
1650     dp_d4410 = temp * f441 * g410;
1651     dp_d4422 = temp * f442 * g422;
1652     temp1 = temp1 * aqnv;
1653     temp = temp1 * root52;
1654     dp_d5220 = temp * f522 * g520;
1655     dp_d5232 = temp * f523 * g532;
1656     temp = 2.0 * temp1 * root54;
1657     dp_d5421 = temp * f542 * g521;
1658     dp_d5433 = temp * f543 * g533;
1659     dp_xlamo = xmao + m_Orbit.RAAN() + m_Orbit.RAAN() - dp_thgr - dp_thgr;
1660     bfact = xlldot + xnodot + xnodot - thdt - thdt;
1661     bfact = bfact + dp_ssl + dp_ssh + dp_ssh;
1662     }
1663     }
1664     else
1665     {
1666     // Synchronous resonance terms initialization
1667     dp_iresfl = true;
1668     dp_isynfl = true;
1669     double g200 = 1.0 + eqsq * (-2.5 + 0.8125 * eqsq);
1670     g310 = 1.0 + 2.0 * eqsq;
1671     double g300 = 1.0 + eqsq * (-6.0 + 6.60937 * eqsq);
1672     f220 = 0.75 * (1.0 + cosiq) * (1.0 + cosiq);
1673     double f311 = 0.9375 * siniq * siniq * (1.0 + 3 * cosiq) - 0.75 * (1.0 + cosiq);
1674     double f330 = 1.0 + cosiq;
1675     f330 = 1.875 * f330 * f330 * f330;
1676     dp_del1 = 3.0 * m_xnodp * m_xnodp * aqnv * aqnv;
1677     dp_del2 = 2.0 * dp_del1 * f220 * g200 * q22;
1678     dp_del3 = 3.0 * dp_del1 * f330 * g300 * q33 * aqnv;
1679     dp_del1 = dp_del1 * f311 * g310 * q31 * aqnv;
1680     dp_fasx2 = 0.13130908;
1681     dp_fasx4 = 2.8843198;
1682     dp_fasx6 = 0.37448087;
1683     dp_xlamo = xmao + m_Orbit.RAAN() + m_Orbit.ArgPerigee() - dp_thgr;
1684     bfact = xlldot + xpidot - thdt;
1685     bfact = bfact + dp_ssl + dp_ssg + dp_ssh;
1686     }
1687    
1688     if (bInitOnExit)
1689     {
1690     dp_xfact = bfact - m_xnodp;
1691    
1692     // Initialize integrator
1693     dp_xli = dp_xlamo;
1694     dp_xni = m_xnodp;
1695     dp_atime = 0.0;
1696     dp_stepp = 720.0;
1697     dp_stepn = -720.0;
1698     dp_step2 = 259200.0;
1699     }
1700    
1701     *eosq = eqsq;
1702     *sinio = siniq;
1703     *cosio = cosiq;
1704     *betao = rteqsq;
1705     *aodp = ao;
1706     *theta2 = cosq2;
1707     *sing = sinomo;
1708     *cosg = cosomo;
1709     *betao2 = bsq;
1710     *xmdot = xlldot;
1711     *omgdot = omgdt;
1712     *xnodott = xnodot;
1713    
1714     return true;
1715     }
1716    
1717     //////////////////////////////////////////////////////////////////////////////
1718     bool cNoradSDP4::DeepCalcDotTerms(double *pxndot, double *pxnddt, double *pxldot)
1719     {
1720     // Dot terms calculated
1721     if (dp_isynfl)
1722     {
1723     *pxndot = dp_del1 * sin(dp_xli - dp_fasx2) +
1724     dp_del2 * sin(2.0 * (dp_xli - dp_fasx4)) +
1725     dp_del3 * sin(3.0 * (dp_xli - dp_fasx6));
1726     *pxnddt = dp_del1 * cos(dp_xli - dp_fasx2) +
1727     2.0 * dp_del2 * cos(2.0 * (dp_xli - dp_fasx4)) +
1728     3.0 * dp_del3 * cos(3.0 * (dp_xli - dp_fasx6));
1729     }
1730     else
1731     {
1732     double xomi = dp_omegaq + omgdt * dp_atime;
1733     double x2omi = xomi + xomi;
1734     double x2li = dp_xli + dp_xli;
1735    
1736     *pxndot = dp_d2201 * sin(x2omi + dp_xli - g22) +
1737     dp_d2211 * sin(dp_xli - g22) +
1738     dp_d3210 * sin(xomi + dp_xli - g32) +
1739     dp_d3222 * sin(-xomi + dp_xli - g32) +
1740     dp_d4410 * sin(x2omi + x2li - g44) +
1741     dp_d4422 * sin(x2li - g44) +
1742     dp_d5220 * sin(xomi + dp_xli - g52) +
1743     dp_d5232 * sin(-xomi + dp_xli - g52) +
1744     dp_d5421 * sin(xomi + x2li - g54) +
1745     dp_d5433 * sin(-xomi + x2li - g54);
1746    
1747     *pxnddt = dp_d2201 * cos(x2omi + dp_xli - g22) +
1748     dp_d2211 * cos(dp_xli - g22) +
1749     dp_d3210 * cos(xomi + dp_xli - g32) +
1750     dp_d3222 * cos(-xomi + dp_xli - g32) +
1751     dp_d5220 * cos(xomi + dp_xli - g52) +
1752     dp_d5232 * cos(-xomi + dp_xli - g52) +
1753     2.0 * (dp_d4410 * cos(x2omi + x2li - g44) +
1754     dp_d4422 * cos(x2li - g44) +
1755     dp_d5421 * cos(xomi + x2li - g54) +
1756     dp_d5433 * cos(-xomi + x2li - g54));
1757     }
1758    
1759     *pxldot = dp_xni + dp_xfact;
1760     *pxnddt = (*pxnddt) * (*pxldot);
1761    
1762     return true;
1763     }
1764    
1765     //////////////////////////////////////////////////////////////////////////////
1766     void cNoradSDP4::DeepCalcIntegrator(double *pxndot, double *pxnddt,
1767     double *pxldot, const double &delt)
1768     {
1769     DeepCalcDotTerms(pxndot, pxnddt, pxldot);
1770    
1771     dp_xli = dp_xli + (*pxldot) * delt + (*pxndot) * dp_step2;
1772     dp_xni = dp_xni + (*pxndot) * delt + (*pxnddt) * dp_step2;
1773     dp_atime = dp_atime + delt;
1774     }
1775    
1776     //////////////////////////////////////////////////////////////////////////////
1777     bool cNoradSDP4::DeepSecular(double *xmdf, double *omgadf, double *xnode,
1778     double *emm, double *xincc, double *xnn,
1779     double *tsince)
1780     {
1781     xll = *xmdf;
1782     omgasm = *omgadf;
1783     xnodes = *xnode;
1784     xn = *xnn;
1785     t = *tsince;
1786    
1787     // Deep space secular effects
1788     xll = xll + dp_ssl * t;
1789     omgasm = omgasm + dp_ssg * t;
1790     xnodes = xnodes + dp_ssh * t;
1791     _em = m_Orbit.Eccentricity() + dp_sse * t;
1792     xinc = m_Orbit.Inclination() + dp_ssi * t;
1793    
1794     if (xinc < 0.0)
1795     {
1796     xinc = -xinc;
1797     xnodes = xnodes + PI;
1798     omgasm = omgasm - PI;
1799     }
1800    
1801     double xnddt = 0.0;
1802     double xndot = 0.0;
1803     double xldot = 0.0;
1804     double ft = 0.0;
1805     double delt = 0.0;
1806    
1807     bool fDone = false;
1808    
1809     if (dp_iresfl)
1810     {
1811     while (!fDone)
1812     {
1813     if ((dp_atime == 0.0) ||
1814     ((t >= 0.0) && (dp_atime < 0.0)) ||
1815     ((t < 0.0) && (dp_atime >= 0.0)))
1816     {
1817     if (t < 0)
1818     delt = dp_stepn;
1819     else
1820     delt = dp_stepp;
1821    
1822     // Epoch restart
1823     dp_atime = 0.0;
1824     dp_xni = m_xnodp;
1825     dp_xli = dp_xlamo;
1826    
1827     fDone = true;
1828     }
1829     else
1830     {
1831     if (fabs(t) < fabs(dp_atime))
1832     {
1833     delt = dp_stepp;
1834    
1835     if (t >= 0.0)
1836     delt = dp_stepn;
1837    
1838     DeepCalcIntegrator(&xndot, &xnddt, &xldot, delt);
1839     }
1840     else
1841     {
1842     delt = dp_stepn;
1843    
1844     delt = dp_stepp;
1845    
1846     fDone = true;
1847     }
1848     }
1849     }
1850    
1851     while (fabs(t - dp_atime) >= dp_stepp)
1852     {
1853     DeepCalcIntegrator(&xndot, &xnddt, &xldot, delt);
1854     }
1855    
1856     ft = t - dp_atime;
1857    
1858     DeepCalcDotTerms(&xndot, &xnddt, &xldot);
1859    
1860     xn = dp_xni + xndot * ft + xnddt * ft * ft * 0.5;
1861    
1862     double xl = dp_xli + xldot * ft + xndot * ft * ft * 0.5;
1863     double temp = -xnodes + dp_thgr + t * thdt;
1864    
1865     xll = xl - omgasm + temp;
1866    
1867     if (!dp_isynfl)
1868     xll = xl + temp + temp;
1869     }
1870    
1871     *xmdf = xll;
1872     *omgadf = omgasm;
1873     *xnode = xnodes;
1874     *emm = _em;
1875     *xincc = xinc;
1876     *xnn = xn;
1877     *tsince = t;
1878    
1879     return true;
1880     }
1881    
1882     //////////////////////////////////////////////////////////////////////////////
1883     bool cNoradSDP4::DeepPeriodics(double *e, double *xincc,
1884     double *omgadf, double *xnode,
1885     double *xmam)
1886     {
1887     _em = *e;
1888     xinc = *xincc;
1889     omgasm = *omgadf;
1890     xnodes = *xnode;
1891     xll = *xmam;
1892    
1893     // Lunar-solar periodics
1894     double sinis = sin(xinc);
1895     double cosis = cos(xinc);
1896    
1897     double sghs = 0.0;
1898     double shs = 0.0;
1899     double sh1 = 0.0;
1900     double pe = 0.0;
1901     double pinc = 0.0;
1902     double pl = 0.0;
1903     double sghl = 0.0;
1904    
1905     if (fabs(dp_savtsn - t) >= 30.0)
1906     {
1907     dp_savtsn = t;
1908    
1909     double zm = dp_zmos + zns * t;
1910     double zf = zm + 2.0 * zes * sin(zm);
1911     double sinzf = sin(zf);
1912     double f2 = 0.5 * sinzf * sinzf - 0.25;
1913     double f3 = -0.5 * sinzf * cos(zf);
1914     double ses = dp_se2 * f2 + dp_se3 * f3;
1915     double sis = dp_si2 * f2 + dp_si3 * f3;
1916     double sls = dp_sl2 * f2 + dp_sl3 * f3 + dp_sl4 * sinzf;
1917    
1918     sghs = dp_sgh2 * f2 + dp_sgh3 * f3 + dp_sgh4 * sinzf;
1919     shs = dp_sh2 * f2 + dp_sh3 * f3;
1920     zm = dp_zmol + znl * t;
1921     zf = zm + 2.0 * zel * sin(zm);
1922     sinzf = sin(zf);
1923     f2 = 0.5 * sinzf * sinzf - 0.25;
1924     f3 = -0.5 * sinzf * cos(zf);
1925    
1926     double sel = dp_ee2 * f2 + dp_e3 * f3;
1927     double sil = dp_xi2 * f2 + dp_xi3 * f3;
1928     double sll = dp_xl2 * f2 + dp_xl3 * f3 + dp_xl4 * sinzf;
1929    
1930     sghl = dp_xgh2 * f2 + dp_xgh3 * f3 + dp_xgh4 * sinzf;
1931     sh1 = dp_xh2 * f2 + dp_xh3 * f3;
1932     pe = ses + sel;
1933     pinc = sis + sil;
1934     pl = sls + sll;
1935     }
1936    
1937     double pgh = sghs + sghl;
1938     double ph = shs + sh1;
1939     xinc = xinc + pinc;
1940     _em = _em + pe;
1941    
1942     if (dp_xqncl >= 0.2)
1943     {
1944     // Apply periodics directly
1945     ph = ph / siniq;
1946     pgh = pgh - cosiq * ph;
1947     omgasm = omgasm + pgh;
1948     xnodes = xnodes + ph;
1949     xll = xll + pl;
1950     }
1951     else
1952     {
1953     // Apply periodics with Lyddane modification
1954     double sinok = sin(xnodes);
1955     double cosok = cos(xnodes);
1956     double alfdp = sinis * sinok;
1957     double betdp = sinis * cosok;
1958     double dalf = ph * cosok + pinc * cosis * sinok;
1959     double dbet = -ph * sinok + pinc * cosis * cosok;
1960    
1961     alfdp = alfdp + dalf;
1962     betdp = betdp + dbet;
1963    
1964     double xls = xll + omgasm + cosis * xnodes;
1965     double dls = pl + pgh - pinc * xnodes * sinis;
1966    
1967     xls = xls + dls;
1968     xnodes = AcTan(alfdp, betdp);
1969     xll = xll + pl;
1970     omgasm = xls - xll - cos(xinc) * xnodes;
1971     }
1972    
1973     *e = _em;
1974     *xincc = xinc;
1975     *omgadf = omgasm;
1976    
1977     *xnode = xnodes;
1978     *xmam = xll;
1979    
1980     return true;
1981     }
1982    
1983     //////////////////////////////////////////////////////////////////////////////
1984     // getPosition()
1985     // This procedure returns the ECI position and velocity for the satellite
1986     // in the orbit at the given number of minutes since the TLE epoch time
1987     // using the NORAD Simplified General Perturbation 4, "deep space" orbit
1988     // model.
1989     //
1990     // tsince - Time in minutes since the TLE epoch (GMT).
1991     // pECI - pointer to location to store the ECI data.
1992     // To convert the returned ECI position vector to km,
1993     // multiply each component by:
1994     // (XKMPER_WGS72 / AE).
1995     // To convert the returned ECI velocity vector to km/sec,
1996     // multiply each component by:
1997     // (XKMPER_WGS72 / AE) * (MIN_PER_DAY / 86400).
1998     bool cNoradSDP4::getPosition(double tsince, cEci &eci)
1999     {
2000     DeepInit(&m_eosq, &m_sinio, &m_cosio, &m_betao, &m_aodp, &m_theta2,
2001     &m_sing, &m_cosg, &m_betao2, &m_xmdot, &m_omgdot, &m_xnodot);
2002    
2003     // Update for secular gravity and atmospheric drag
2004     double xmdf = m_Orbit.mnAnomaly() + m_xmdot * tsince;
2005     double omgadf = m_Orbit.ArgPerigee() + m_omgdot * tsince;
2006     double xnoddf = m_Orbit.RAAN() + m_xnodot * tsince;
2007     double tsq = tsince * tsince;
2008     double xnode = xnoddf + m_xnodcf * tsq;
2009     double tempa = 1.0 - m_c1 * tsince;
2010     double tempe = m_Orbit.BStar() * m_c4 * tsince;
2011     double templ = m_t2cof * tsq;
2012     double xn = m_xnodp;
2013     double em;
2014     double xinc;
2015    
2016     DeepSecular(&xmdf, &omgadf, &xnode, &em, &xinc, &xn, &tsince);
2017    
2018     double a = pow(XKE / xn, TWOTHRD) * sqr(tempa);
2019     double e = em - tempe;
2020     double xmam = xmdf + m_xnodp * templ;
2021    
2022     DeepPeriodics(&e, &xinc, &omgadf, &xnode, &xmam);
2023    
2024     double xl = xmam + omgadf + xnode;
2025    
2026     xn = XKE / pow(a, 1.5);
2027    
2028     return FinalPosition(xinc, omgadf, e, a, xl, xnode, xn, tsince, eci);
2029     }
2030    
2031    
2032     // cOrbit.cpp
2033     //
2034     // Copyright (c) 2002-2003 Michael F. Henry
2035     //
2036     // mfh 11/15/2003
2037     //
2038     //////////////////////////////////////////////////////////////////////
2039     cOrbit::cOrbit(const cTle &tle) :
2040     m_tle(tle),
2041     m_pNoradModel(NULL)
2042     {
2043     m_tle.Initialize();
2044    
2045     int epochYear = (int)m_tle.getField(cTle::FLD_EPOCHYEAR);
2046     double epochDay = m_tle.getField(cTle::FLD_EPOCHDAY );
2047    
2048     if (epochYear < 57)
2049     epochYear += 2000;
2050     else
2051     epochYear += 1900;
2052    
2053     m_jdEpoch = cJulian(epochYear, epochDay);
2054    
2055     m_secPeriod = -1.0;
2056    
2057     // Recover the original mean motion and semimajor axis from the
2058     // input elements.
2059     double mm = mnMotion();
2060     double rpmin = mm * 2 * PI / MIN_PER_DAY; // rads per minute
2061    
2062     double a1 = pow(XKE / rpmin, TWOTHRD);
2063     double e = Eccentricity();
2064     double i = Inclination();
2065     double temp = (1.5 * CK2 * (3.0 * sqr(cos(i)) - 1.0) /
2066     pow(1.0 - e * e, 1.5));
2067     double delta1 = temp / (a1 * a1);
2068     double a0 = a1 *
2069     (1.0 - delta1 *
2070     ((1.0 / 3.0) + delta1 *
2071     (1.0 + 134.0 / 81.0 * delta1)));
2072    
2073     double delta0 = temp / (a0 * a0);
2074    
2075     m_mnMotionRec = rpmin / (1.0 + delta0);
2076     m_aeAxisSemiMinorRec = a0 / (1.0 - delta0);
2077     m_aeAxisSemiMajorRec = m_aeAxisSemiMinorRec / sqrt(1.0 - (e * e));
2078     m_kmPerigeeRec = XKMPER_WGS72 * (m_aeAxisSemiMajorRec * (1.0 - e) - AE);
2079     m_kmApogeeRec = XKMPER_WGS72 * (m_aeAxisSemiMajorRec * (1.0 + e) - AE);
2080    
2081     if (2.0 * PI / m_mnMotionRec >= 225.0)
2082     {
2083     // SDP4 - period >= 225 minutes.
2084     m_pNoradModel = new cNoradSDP4(*this);
2085     }
2086     else
2087     {
2088     // SGP4 - period < 225 minutes
2089     m_pNoradModel = new cNoradSGP4(*this);
2090     }
2091     }
2092    
2093     /////////////////////////////////////////////////////////////////////////////
2094     cOrbit::~cOrbit()
2095     {
2096     delete m_pNoradModel;
2097     }
2098    
2099     //////////////////////////////////////////////////////////////////////////////
2100     // Return the period in seconds
2101     double cOrbit::Period() const
2102     {
2103     if (m_secPeriod < 0.0)
2104     {
2105     // Calculate the period using the recovered mean motion.
2106     if (m_mnMotionRec == 0)
2107     m_secPeriod = 0.0;
2108     else
2109     m_secPeriod = (2 * PI) / m_mnMotionRec * 60.0;
2110     }
2111    
2112     return m_secPeriod;
2113     }
2114    
2115     //////////////////////////////////////////////////////////////////////////////
2116     // Returns elapsed number of seconds from epoch to given time.
2117     // Note: "Predicted" TLEs can have epochs in the future.
2118     double cOrbit::TPlusEpoch(const cJulian &gmt) const
2119     {
2120     return gmt.spanSec(Epoch());
2121     }
2122    
2123     //////////////////////////////////////////////////////////////////////////////
2124     // Returns the mean anomaly in radians at given GMT.
2125     // At epoch, the mean anomaly is given by the elements data.
2126     double cOrbit::mnAnomaly(cJulian gmt) const
2127     {
2128     double span = TPlusEpoch(gmt);
2129     double P = Period();
2130    
2131     assert(P != 0.0);
2132    
2133     return fmod(mnAnomaly() + (TWOPI * (span / P)), TWOPI);
2134     }
2135    
2136     //////////////////////////////////////////////////////////////////////////////
2137     // getPosition()
2138     // This procedure returns the ECI position and velocity for the satellite
2139     // at "tsince" minutes from the (GMT) TLE epoch. The vectors returned in
2140     // the ECI object are kilometer-based.
2141     // tsince - Time in minutes since the TLE epoch (GMT).
2142     bool cOrbit::getPosition(double tsince, cEci *pEci) const
2143     {
2144     bool rc;
2145    
2146     rc = m_pNoradModel->getPosition(tsince, *pEci);
2147    
2148     pEci->ae2km();
2149    
2150     return rc;
2151     }
2152    
2153     //////////////////////////////////////////////////////////////////////////////
2154     // SatName()
2155     // Return the name of the satellite. If requested, the NORAD number is
2156     // appended to the end of the name, i.e., "ISS (ZARYA) #25544".
2157     // The name of the satellite with the NORAD number appended is important
2158     // because many satellites, especially debris, have the same name and
2159     // would otherwise appear to be the same satellite in ouput data.
2160     string cOrbit::SatName(bool fAppendId /* = false */) const
2161     {
2162     string str = m_tle.getName();
2163    
2164     if (fAppendId)
2165     {
2166     string strId;
2167    
2168     m_tle.getField(cTle::FLD_NORADNUM, cTle::U_NATIVE, &strId);
2169     str = str + " #" + strId;
2170     }
2171    
2172     return str;
2173     }
2174    

  ViewVC Help
Powered by ViewVC 1.1.23