• Skip to content
  • Skip to link menu
KDE API Reference
  • KDE API Reference
  • kdeedu API Reference
  • KDE Home
  • Contact Us
 

marble

  • sources
  • kde-4.12
  • kdeedu
  • marble
  • src
  • plugins
  • render
  • satellites
  • mex
astrolib.cpp
Go to the documentation of this file.
1 //
2 // This file is part of the Marble Virtual Globe.
3 //
4 // This program is free software licensed under the GNU LGPL. You can
5 // find a copy of this license in LICENSE.txt in the top directory of
6 // the source code.
7 //
8 // Copyright 2012 Gerhard HOLTKAMP
9 //
10 
11 /* =========================================================================
12  Procedures needed for standard astronomy programs.
13  The majority of procedures in this unit are taken from Montenbruck,
14  Pfleger "Astronomie mit dem Personal Computer", Springer Verlag, 1989
15  as well as from the "Explanatory Supplement to the Astronomical Almanac"
16  University Science Books, Mill Valley, California, 1992
17  and modified correspondingly.
18 
19  Linux C++ version
20  Copyright : Gerhard HOLTKAMP 11-JAN-2012
21  ========================================================================= */
22 
23 #include <cmath>
24 using namespace std;
25 
26 #include "attlib.h"
27 #include "astrolib.h"
28 
29 double frac (double f)
30  { return fmod(f,1.0); }
31 
32 /*--------------- Function ddd ------------------------------------------*/
33 
34 double ddd (int d, int m, double s)
35  {
36  /*
37  Conversion from Degrees, Minutes, Seconds into decimal degrees
38 
39  INPUT :
40  d : degrees
41  m : minutes
42  s : seconds (and fractions)
43 
44  OUTPUT :
45  dd : decimal degrees
46  */
47  double dd;
48 
49  if ( (d < 0) || (m < 0) || (s < 0)) dd = -1.0;
50  else dd = 1.0;
51  dd = dd * (fabs(double(d)) + fabs(double(m))/60.0 + fabs(s)/3600.0);
52 
53  return dd;
54  }
55 
56 /*--------------- Function dms ------------------------------------------*/
57 
58 void dms (double dd, int &d, int &m, double &s)
59  /*
60  Conversion from decimal degrees into Degrees, Minutes, Seconds
61 
62  INPUT :
63  dd : decimal degrees
64 
65  OUTPUT :
66  d : degrees
67  m : minutes
68  s : seconds (and fractions)
69  */
70  {
71  double d1;
72 
73  d1 = fabs(dd);
74  // d = int(floor(d1));
75  d = int(d1);
76  d1 = (d1 - d) * 60.0;
77  // m = int(floor(d1));
78  m = int(d1);
79  s = (d1 - m) * 60.0;
80  if (dd < 0)
81  {
82  if (d != 0) d = -d;
83  else
84  {
85  if (m != 0) m = -m;
86  else s = -s;
87  };
88  };
89  }
90 
91  /*------------ FUNCTION mjd ----------------------------------------------*/
92 
93 double mjd (int day, int month, int year, double hour)
94 /*
95  Modified Julian Date ( MJD = Julian Date - 2400000.5)
96  valid for every date
97  Julian Calendar up to 4-OCT-1582,
98  Gregorian Calendar from 15-OCT-1582.
99  */
100  {
101  double a;
102  long int b, c;
103 
104  a = 10000.0 * year + 100.0 * month + day;
105  if (month <= 2)
106  {
107  month = month + 12;
108  year = year - 1;
109  };
110  if (a <= 15821004.1)
111  {
112  b = ((year+4716)/4) - 1181;
113  if (year < -4716)
114  {
115  c = year + 4717;
116  c = -c;
117  c = c / 4;
118  c = -c;
119  b = c - 1182;
120  };
121  }
122  else b = (year/400) - (year/100) + (year/4);
123  // { b = -2 + floor((year+4716)/4) - 1179;}
124  // else {b = floor(year/400) - floor(year/100) + floor(year/4);};
125 
126  a = 365.0 * year - 679004.0;
127  a = a + b + int(30.6001 * (month + 1)) + day + hour / 24.0;
128 
129  return a;
130  }
131 
132 /*---------------- Function julcent ---------------------------------------*/
133 
134 double julcent (double mjuld)
135  {
136  /*
137  Julian Centuries since J2000.0 of Modified Julian Date mjuld.
138  */
139  return (mjuld - 51544.5) / 36525.0;
140  }
141 
142 /*---------------- Function caldat -----------------------------------------*/
143 
144 void caldat (double mjd, int &day, int &month, int &year, double &hour)
145  /*
146  Calculate the calendar date from the Modified Julian Date
147 
148  INPUT :
149  mjd : Modified Julian Date (Julian Date - 2400000.5)
150 
151  OUTPUT :
152  day, month, year : corresponding date
153  hour : corresponding hours of the above date
154  */
155  {
156  long int b, c, d, e, f, jd0;
157 
158  jd0 = long(mjd + 2400001.0);
159  if (jd0 < 2299161) c = jd0 + 1524; /* Julian calendar */
160  else
161  { /* Gregorian calendar */
162  b = long (( jd0 - 1867216.25) / 36524.25);
163  c = jd0 + b - (b/4) + 1525;
164  };
165 
166  if (mjd < -2400001.0) // special case for year < -4712
167  {
168  if (mjd == floor(mjd)) jd0 = jd0 + 1;
169  c = long((-jd0 - 0.1)/ 365.25);
170  c = c + 1;
171  year = -4712 - c;
172  d = c / 4;
173  d = c * 365 + d; // number of days from JAN 1, -4712
174  f = d + jd0; // day of the year
175  if ((c % 4) == 0) e = 61;
176  else e = 60;
177  if (f == 0)
178  {
179  year = year - 1;
180  month = 12;
181  day = 31;
182  f = 500; // set as a marker
183  };
184  if (f < e)
185  {
186  if (f < 32)
187  {
188  month = 1;
189  day = f;
190  }
191  else
192  {
193  month = 2;
194  day = f - 31;
195  };
196  }
197  else
198  {
199  if (f < 500)
200  {
201  f = f - e;
202  month = long((f + 123.0) / 30.6001);
203  day = f - long(month * 30.6001) + 123;
204  month = month - 1;
205  };
206  };
207  }
208  else // normal case
209  {
210  d = long ((c - 122.1) / 365.25);
211  e = 365 * d + (d/4);
212  f = long ((c - e) / 30.6001);
213  day = c - e - long(30.6001 * f);
214  month = f - 1 - 12 * (f / 14);
215  year = d - 4715 - ((7 + month) / 10);
216  };
217 
218  hour = 24.0 * (mjd - floor(mjd));
219  }
220 
221 /*---------------- Function DefTdUt --------------------------------------*/
222 double DefTdUt (int yi)
223  {
224  /*
225  Get a suitable default for value of TDT - UT in year yi in seconds.
226  */
227 
228  const int td[9] = {55,55,56,56,57,58,58,59,60};
229  double t, result;
230  int yy, yr;
231 
232  yr = yi;
233 
234  if (yr > 1899)
235  {
236  if (yr >= 2000) // this corrects for observations until DEC 2006
237  {
238  if (yr > 2006) yr -= 1;
239  if (yr > 2006) yr -= 1; // no leap second at the end of 2007
240  if (yr > 2007) yr -= 1; // no leap second at the end of 2009
241  if (yr > 2007) yr -= 1; // no leap second at the end of 2010
242  yr -= 6;
243  if (yr < 1999) yr = 1999;
244  };
245  t = yr - 2000;
246  t = t / 100.0;
247  result = (((-339.84*t - 516.12)*t -160.22)*t + 92.23)*t + 71.28;
248  if (yr > 2013)
249  {
250  t = yr - 2013;
251  t = t / 100.0;
252  result = (27.5*t + 75.0)*t + 73.0;
253  }
254  else if (yr > 1994)
255  {
256  result -= 6.28;
257  }
258  else if (yr > 1985)
259  {
260  yy = yr - 1986;
261  result = td[yy];
262  };
263  };
264 
265  if (yr < 1900)
266  {
267  t = yr - 1800;
268  t = t / 100;
269  if (yr > 1649)
270  {
271  t = t - 0.19;
272  result = 5.156 + 13.3066 * t * t;
273  if (yr > 1864) result = -0.6*(yr - 1865) + 6;
274  if (yr > 1884) result = -0.2*(yr - 1885) - 6;
275  }
276  else
277  {
278  if (yr > 947) result = 25.5 * t*t;
279  else result = 1360 + (44.3*t + 320.0) * t;
280  };
281  };
282 
283 // round to full seconds
284 
285  if (result < 0)
286  {
287  t = -result;
288  t += 0.5;
289  result = - floor(t);
290  }
291  else result = floor(result + 0.5);
292 
293  result += 0.184; // because of TT = TAI + 32.184;
294 
295  return result;
296  }
297 
298 /*---------------- Function lsidtim -------------------------------------*/
299 
300 double lsidtim (double jd, double lambda, double ep2)
301  {
302  /* Calculate the Apparent Local Mean Sidereal Time (in decimal hours) for
303  the Modified Julian Date jd and geographic longitude lambda (in degrees)
304  and with the correction (AST - MST) eps (given in seconds)
305  */
306  double lmst;
307  double t, ut, gmst;
308  int mjd0;
309 
310  mjd0 = int(jd);
311  ut = (jd - mjd0)*24.0;
312  t = (mjd0 - 51544.5) / 36525.0;
313  gmst = 6.697374558 + 1.0027379093*ut
314  + (8640184.812866 + (0.093104 - 6.2e-6*t)*t)*t/3600.0;
315  lmst = 24.0 * frac((gmst + lambda/15.0) / 24.0);
316  lmst = lmst + ep2 / 3600.0;
317 
318  return lmst;
319  }
320 
321 /*---------------- Function eps -----------------------------------------*/
322 
323 double eps (double t) // obliquity of ecliptic
324  {
325  /*
326  Obliquety of ecliptic eps (in radians)
327  at time t (in Julian Centuries from J2000.0)
328  */
329  double tp;
330 
331  tp = 23.43929111 - (46.815+(0.00059-0.001813*t)*t)*t/3600.0;
332  tp = 1.74532925199e-2 * tp;
333 
334  return tp;
335  }
336 
337 /*---------------- Function eclequ -----------------------------------------*/
338 
339 Vec3 eclequ (double t, Vec3& r1) // ecliptic -> equatorial
340  {
341  /*
342  Convert position vector r1 from ecliptic into equatorial coordinates
343  at t (in Julian Centuries from J2000.0 )
344  */
345 
346  Mat3 m;
347  Vec3 r2;
348 
349  m = xrot (-eps(t));
350  r2 = mxvct (m, r1);
351  return r2;
352  }
353 
354 /*---------------- Function equecl -----------------------------------------*/
355 
356 Vec3 equecl (double t, Vec3& r1) // equatorial -> ecliptic
357  {
358  /*
359  Convert position vector r1 from equatorial into ecliptic coordinates
360  at t (in Julian Centuries from J2000.0)
361  */
362 
363  Mat3 m;
364  Vec3 r2;
365 
366  m = xrot (eps(t));
367  r2 = mxvct (m, r1);
368  return r2;
369  }
370 
371 /*---------------- Function pmatecl -----------------------------------------*/
372 
373 Mat3 pmatecl (double t1, double t2) // ecl. precession
374  {
375  /*
376  Calculate ecliptic precession matrix a from t1 to t2
377  (times in Julian Centuries from J2000.0)
378  */
379 
380  const double secrad = 4.8481368111e-6; // arcsec -> radians
381 
382  double ppi, pii, pa, dt;
383  Mat3 m1, m2, a;
384 
385  dt = t2 - t1;
386  ppi = 174.876383889 + (((3289.4789+0.60622*t1)*t1) +
387  ((-869.8089-0.50491*t1) + 0.03536*dt)*dt) / 3600.0;
388  ppi = ppi * 1.74532925199e-2;
389  pii = ((47.0029-(0.06603-0.000598*t1)*t1) +
390  ((-0.03302+0.000598*t1)+0.00006*dt)*dt)*dt * secrad;
391  pa = ((5029.0966+(2.22226-0.000042*t1)*t1) +
392  ((1.11113-0.000042*t1)-0.000006*dt)*dt)*dt * secrad;
393 
394  pa = ppi + pa;
395  m1 = zrot (-pa);
396  m2 = xrot (pii);
397  m1 = m1 * m2;
398  m2 = zrot (ppi);
399  a = m1 * m2;
400 
401  return a;
402  }
403 
404 /*---------------- Function pmatequ --------------------------------------*/
405 
406 Mat3 pmatequ (double t1, double t2) // equ. precession
407  {
408  /*
409  Calculate equatorial precession matrix a from t1 to t2
410  (times in Julian Centuries from J2000.0)
411  */
412  const double secrad = 4.8481368111e-6; // arcsec -> radians
413 
414  double dt, zeta, z, theta;
415  Mat3 m1, m2, a;
416 
417  dt = t2 - t1;
418  zeta = ((2306.2181+(1.39656-0.000139*t1)*t1) +
419  ((0.30188-0.000345*t1)+0.017998*dt)*dt)*dt * secrad;
420  z = zeta + ((0.7928+0.000411*t1)+0.000205*dt)*dt*dt * secrad;
421  theta = ((2004.3109-(0.8533+0.000217*t1)*t1) -
422  ((0.42665+0.000217*t1)+0.041833*dt)*dt)*dt * secrad;
423 
424  m1 = zrot (-z);
425  m2 = yrot (theta);
426  m1 = m1 * m2;
427  m2 = zrot (-zeta);
428  a = m1 * m2;
429 
430  return a;
431  }
432 
433 /*---------------- Function nutmat --------------------------------------*/
434 
435 Mat3 nutmat (double t, double& ep2, bool hpr)
436  {
437  /*
438  Calculate nutation matrix a from mean to true equatorial coordinates
439  at time t (in Julian Centuries from J2000.0)
440  Also calculates the correction ep2 for apparent sidereal time in sec
441 
442  if hpr is true a high precision is used, otherwise a low precision
443  (only the first 50 terms of the nutation theory are used)
444  */
445  const int ntb1 = 15;
446  const int tb1[ntb1][5] =
447  {
448  { 0, 0, 0, 0, 1}, // 1
449  { 0, 0, 0, 0, 2}, // 2
450  { 0, 0, 2,-2, 2}, // 9
451  { 0, 1, 0, 0, 0}, // 10
452  { 0, 1, 2,-2, 2}, // 11
453  { 0,-1, 2,-2, 2}, // 12
454  { 0, 0, 2,-2, 1}, // 13
455  { 0, 2, 0, 0, 0}, // 16
456  { 0, 2, 2,-2, 2}, // 18
457  { 0, 0, 2, 0, 2}, // 31
458  { 1, 0, 0, 0, 0}, // 32
459  { 0, 0, 2, 0, 1}, // 33
460  { 1, 0, 2, 0, 2}, // 34
461  { 1, 0, 0, 0, 1}, // 38
462  { -1, 0, 0, 0, 1}, // 39
463  };
464 
465  const int ntb2 = 35;
466  const int tb2[ntb2][5] =
467  {
468  { -2, 0, 2, 0, 1}, // 3
469  { 2, 0,-2, 0, 0}, // 4
470  { -2, 0, 2, 0, 2}, // 5
471  { 1,-1, 0,-1, 0}, // 6
472  { 0,-2, 2,-2, 1}, // 7
473  { 2, 0,-2, 0, 1}, // 8
474  { 2, 0, 0,-2, 0}, // 14
475  { 0, 0, 2,-2, 0}, // 15
476  { 0, 1, 0, 0, 1}, // 17
477  { 0,-1, 0, 0, 1}, // 19
478  { -2, 0, 0, 2, 1}, // 20
479  { 0,-1, 2,-2, 1}, // 21
480  { 2, 0, 0,-2, 1}, // 22
481  { 0, 1, 2,-2, 1}, // 23
482  { 1, 0, 0,-1, 0}, // 24
483  { 2, 1, 0,-2, 0}, // 25
484  { 0, 0,-2, 2, 1}, // 26
485  { 0, 1,-2, 2, 0}, // 27
486  { 0, 1, 0, 0, 2}, // 28
487  { -1, 0, 0, 1, 1}, // 29
488  { 0, 1, 2,-2, 0}, // 30
489  { 1, 0, 0,-2, 0}, // 35
490  { -1, 0, 2, 0, 2}, // 36
491  { 0, 0, 0, 2, 0}, // 37
492  { -1, 0, 2, 2, 2}, // 40
493  { 1, 0, 2, 0, 1}, // 41
494  { 0, 0, 2, 2, 2}, // 42
495  { 2, 0, 0, 0, 0}, // 43
496  { 1, 0, 2,-2, 2}, // 44
497  { 2, 0, 2, 0, 2}, // 45
498  { 0, 0, 2, 0, 0}, // 46
499  { -1, 0, 2, 0, 1}, // 47
500  { -1, 0, 0, 2, 1}, // 48
501  { 1, 0, 0,-2, 1}, // 49
502  { -1, 0, 2, 2, 1}, // 50
503  };
504 
505  const double tb3[ntb1][4] =
506  {
507  {-171996.0,-174.2, 92025.0, 8.9 }, // 1
508  { 2062.0, 0.2, -895.0, 0.5 }, // 2
509  { -13187.0, -1.6, 5736.0, -3.1 }, // 9
510  { 1426.0, -3.4, 54.0, -0.1 }, // 10
511  { -517.0, 1.2, 224.0, -0.6 }, // 11
512  { 217.0, -0.5, -95.0, 0.3 }, // 12
513  { 129.0, 0.1, -70.0, 0.0 }, // 13
514  { 17.0, -0.1, 0.0, 0.0 }, // 16
515  { -16.0, 0.1, 7.0, 0.0 }, // 18
516  { -2274.0, -0.2, 977.0, -0.5 }, // 31
517  { 712.0, 0.1, -7.0, 0.0 }, // 32
518  { -386.0, -0.4, 200.0, 0.0 }, // 33
519  { -301.0, 0.0, 129.0, -0.1 }, // 34
520  { 63.0, 0.1, -33.0, 0.0 }, // 38
521  { -58.0, -0.1, 32.0, 0.0 }, // 39
522  };
523 
524  const double tb4[ntb2][2] =
525  {
526  { 46.0, -24.0}, // 3
527  { 11.0, 0.0 }, // 4
528  { -3.0, 1.0 }, // 5
529  { -3.0, 0.0 }, // 6
530  { -2.0, 1.0 }, // 7
531  { 1.0, 0.0 }, // 8
532  { 48.0, 1.0 }, // 14
533  {-22.0, 0.0 }, // 15
534  {-15.0, 9.0 }, // 17
535  {-12.0, 6.0 }, // 19
536  { -6.0, 3.0 }, // 20
537  { -5.0, 3.0 }, // 21
538  { 4.0, -2.0 }, // 22
539  { 4.0, -2.0 }, // 23
540  { -4.0, 0.0 }, // 24
541  { 1.0, 0.0 }, // 25
542  { 1.0, 0.0 }, // 26
543  { -1.0, 0.0 }, // 27
544  { 1.0, 0.0 }, // 28
545  { 1.0, 0.0 }, // 29
546  { -1.0, 0.0 }, // 30
547  {-158.0,-1.0 }, // 35
548  {123.0,-53.0 }, // 36
549  { 63.0, -2.0 }, // 37
550  {-59.0, 26.0 }, // 40
551  {-51.0, 27.0 }, // 41
552  {-38.0, 16.0 }, // 42
553  { 29.0, -1.0 }, // 43
554  { 29.0,-12.0 }, // 44
555  {-31.0, 13.0 }, // 45
556  { 26.0, -1.0 }, // 46
557  { 21.0,-10.0 }, // 47
558  { 16.0, -8.0 }, // 48
559  {-13.0, 7.0 }, // 49
560  {-10.0, 5.0 }, // 50
561  };
562 
563  const double secrad = 4.8481368111e-6; // arcsec -> radians
564  const double p2 = 2.0 * M_PI;
565 
566  double ls, lm, d, f, n, dpsi, deps, ep0;
567  int j;
568  Mat3 m1, m2, a;
569 
570  if (hpr)
571  {
572  lm = 2.355548393544 +
573  (8328.691422883903 + (0.000151795164 + 0.000000310281*t)*t)*t;
574  ls = 6.240035939326 +
575  (628.301956024185 + (-0.000002797375 - 0.000000058178*t)*t)*t;
576  f = 1.627901933972 +
577  (8433.466158318464 + (-0.000064271750 + 0.000000053330*t)*t)*t;
578  d = 5.198469513580 +
579  (7771.377146170650 + (-0.000033408511 + 0.000000092115*t)*t)*t;
580  n = 2.182438624361 +
581  (-33.757045933754 + (0.000036142860 + 0.000000038785*t)*t)*t;
582 
583  lm = fmod(lm,p2);
584  ls = fmod(ls,p2);
585  f = fmod(f,p2);
586  d = fmod(d,p2);
587  n = fmod(n,p2);
588 
589  dpsi = 0.0;
590  deps = 0.0;
591  for(j=0; j<ntb1; ++j)
592  {
593  ep0 = tb1[j][0]*lm + tb1[j][1]*ls + tb1[j][2]*f + tb1[j][3]*d + tb1[j][4]*n;
594  dpsi = dpsi + (tb3[j][0]+tb3[j][1]*t) * sin(ep0);
595  deps = deps + (tb3[j][2]+tb3[j][3]*t) * cos(ep0);
596  };
597  for(j=0; j<ntb2; ++j)
598  {
599  ep0 = tb2[j][0]*lm + tb2[j][1]*ls + tb2[j][2]*f + tb2[j][3]*d + tb2[j][4]*n;
600  dpsi = dpsi + tb4[j][0] * sin(ep0);
601  deps = deps + tb4[j][1] * cos(ep0);
602  };
603  dpsi = 1.0e-4 * dpsi * secrad;
604  deps = 1.0e-4 * deps * secrad;
605  }
606 
607  else // low precision
608  {
609  ls = p2 * frac (0.993133+99.997306*t); // mean anomaly sun
610  d = p2 * frac (0.827362+1236.853087*t); // diff long. moon-sun
611  f = p2 * frac (0.259089+1342.227826*t); // dist. node
612  n = p2 * frac (0.347346 - 5.372447*t); // long. node
613 
614  dpsi = (-17.2*sin(n) - 1.319*sin(2*(f-d+n)) - 0.227*sin(2*(f+n))
615  + 0.206*sin(2*n) + 0.143*sin(ls)) * secrad;
616  deps = (+9.203*cos(n) + 0.574*cos(2*(f-d+n)) + 0.098*cos(2*(f+n))
617  -0.09*cos(2*n) ) * secrad;
618  };
619 
620  ep0 = eps (t);
621  ep2 = ep0 + deps;
622  m1 = xrot (ep0);
623  m2 = zrot (-dpsi);
624  m1 *= m2;
625  m2 = xrot (-ep2);
626  a = m2 * m1;
627  ep2 = dpsi * cos (ep2);
628  ep2 *= 13750.9870831; // convert radians into time-seconds
629 
630  return a;
631  }
632 
633 /*---------------- Function nutecl --------------------------------------*/
634 
635 Mat3 nutecl (double t, double& ep2) // nutation matrix (ecliptic)
636  {
637  /*
638  Calculate nutation matrix a from mean to true ecliptic coordinates
639  at time t (in Julian Centuries from J2000.0)
640  Also calculates the correction ep2 for apparent sidereal time in sec
641  */
642 
643  const double secrad = 4.8481368111e-6; // arcsec -> radians
644  const double p2 = 2.0 * M_PI;
645 
646  double ls, d, f, n, dpsi, deps, ep0;
647  Mat3 m1, m2, a;
648 
649  ls = p2 * frac (0.993133+99.997306*t); // mean anomaly sun
650  d = p2 * frac (0.827362+1236.853087*t); // diff long. moon-sun
651  f = p2 * frac (0.259089+1342.227826*t); // dist. node
652  n = p2 * frac (0.347346 - 5.372447*t); // long. node
653 
654  dpsi = (-17.2*sin(n) - 1.319*sin(2*(f-d+n)) - 0.227*sin(2*(f+n))
655  + 0.206*sin(2*n) + 0.143*sin(ls)) * secrad;
656 
657  deps = (+9.203*cos(n) + 0.574*cos(2*(f-d+n)) + 0.098*cos(2*(f+n))
658  -0.09*cos(2*n) ) * secrad;
659 
660  ep0 = eps (t);
661  ep2 = ep0 + deps;
662  m1 = xrot (-deps);
663  m2 = zrot (-dpsi);
664  a = m1 * m2;
665  ep2 = dpsi * cos (ep2);
666  ep2 *= 13750.9870831; // convert radians into time-seconds
667 
668  return a;
669  }
670 
671 /*---------------- Function pmatequ --------------------------------------*/
672 
673 Mat3 PoleMx (double xp, double yp)
674  {
675  /* Returns Polar Motion matrix.
676  xp and yp are the coordinates of the Celestial Ephemeris Pole
677  with respect to the IERS Reference Pole in arcsec as published
678  in the IERS Bulletin B.
679  */
680  const double arctrd = M_PI / (180.0*3600.0);
681  Mat3 res;
682  double xr, yr;
683 
684  xr = xp * arctrd;
685  yr = yp * arctrd;
686  res.assign(1.0,0.0,xr,0.0,1.0,-yr,-xr,yr,1.0);
687 
688  return res;
689  }
690 
691 /*---------------- Function aberrat --------------------------------------*/
692 
693 Vec3 aberrat (double t, Vec3& ve) // aberration
694  {
695  /*
696  Correct position vector ve for aberration into new position va
697  at time t (in Julian Centuries from J2000.0)
698  */
699  const double p2 = 2.0 * M_PI;
700 
701  double l, cl, d0;
702  Vec3 va;
703 
704  d0 = abs(ve);
705  l = p2 * frac(0.27908+100.00214*t);
706  cl = cos(l)*d0;
707  va[0] = ve[0] - 9.934e-5 * sin(l) * d0;
708  va[1] = ve[1] + 9.125e-5 * cl;
709  va[2] = ve[2] + 3.927e-5 * cl;
710 
711  return va;
712  }
713 
714 /*--------------------- Function GeoPos ----------------------------------*/
715 
716 Vec3 GeoPos (double jd, double ep2, double lat, double lng, double ht)
717  {
718  /* Return the geocentric vector (in the equatorial system) of the
719  geographic position given by the latitude lat and longitude lng
720  (in radians) and height ht (in m) at the MJD-time jd (UT).
721  ep2 : correction for apparent sidereal time in sec
722 
723  The length unit of the vector is in terms of the equatorial
724  Earth radius (6378.137 km)
725  */
726  double const e2 = 6.69438499959e-3;
727  double np, h, sp;
728 
729  Vec3 r;
730 
731  sp = sin(lat);
732  h = ht / 6378.137e3;
733  np = 1.0 / (sqrt(1.0 - e2*sp*sp));
734  r[2] = ((1.0 - e2)*np + h)*sp;
735 
736  sp = (np + h) * cos(lat);
737  np = lsidtim(jd, (lng*180.0/M_PI), ep2) * M_PI/12.0;
738 
739  r[0] = sp * cos(np);
740  r[1] = sp * sin(np);
741 
742  return r;
743  }
744 
745 Vec3 GeoPos (double jd, double ep2, double lat, double lng, double ht,
746  double xp, double yp)
747  {
748  /* Return the geocentric vector (in the equatorial system) of the
749  geographic position given by the latitude lat and longitude lng
750  (in radians) and height ht (in m) at the MJD-time jd (UT).
751 
752  ep2 : correction for apparent sidereal time in sec
753  xp, yp: coordinates of polar motion in arcsec.
754 
755  The length unit of the vector is in terms of the equatorial
756  Earth radius (6378.137 km)
757  */
758  double const e2 = 6.69438499959e-3;
759  double np, h, sp;
760  Mat3 prx;
761  Vec3 r;
762 
763  sp = sin(lat);
764  h = ht / 6378.137e3;
765  np = 1.0 / (sqrt(1.0 - e2*sp*sp));
766  r[2] = ((1.0 - e2)*np + h)*sp;
767  sp = (np + h) * cos(lat);
768  r[0] = sp * cos(lng);
769  r[1] = sp * sin(lng);
770 
771  // correct for polar motion
772  if ((xp != 0) || (yp != 0))
773  {
774  prx = mxtrn(PoleMx(xp, yp));
775  r = mxvct(prx, r);
776  };
777 
778  // convert into equatorial
779  np = lsidtim(jd, 0.0, ep2) * M_PI/12.0;
780  prx = zrot(-np);
781  r = mxvct(prx, r);
782 
783  return r;
784  }
785 
786 /*--------------------- Function EquHor ----------------------------------*/
787 
788 Vec3 EquHor (double jd, double ep2, double lat, double lng, Vec3 r)
789  {
790  /* convert vector r from the equatorial system into the horizontal system.
791  jd = MJD-time (UT)
792  ep2 : correction for apparent sidereal time in sec
793  lat, lng : geographic latitude and longitude (in radians)
794  */
795  double lst;
796  Vec3 s;
797  Mat3 mx;
798 
799  lst = lsidtim(jd, (lng*180.0/M_PI), ep2) * M_PI/12.0;
800  mx = zrot(lst);
801  s = mxvct(mx, r);
802  mx = yrot(M_PI/2.0 - lat);
803  s = mxvct(mx, s);
804 
805  return s;
806  }
807 
808 /*--------------------- Function HorEqu ----------------------------------*/
809 
810 
811 Vec3 HorEqu (double jd, double ep2, double lat, double lng, Vec3 r)
812  {
813  /* convert vector r from the horizontal system into the equatorial system.
814  jd = MJD-time (UT)
815  ep2 : correction for apparent sidereal time in sec
816  lat, lng : geographic latitude and longitude (in radians)
817  */
818  double lst;
819  Vec3 s;
820  Mat3 mx;
821 
822  mx = yrot(lat - M_PI/2.0);
823  s = mxvct(mx, r);
824  lst = -lsidtim(jd, (lng*180.0/M_PI), ep2) * M_PI/12.0;
825  mx = zrot(lst);
826  s = mxvct(mx, s);
827 
828  return s;
829  }
830 
831 /*--------------------- Function AppPos ----------------------------------*/
832 
833 void AppPos (double jd, double ep2, double lat, double lng, double ht,
834  int solsys, Vec3 r, double& azim, double& elev, double& dist)
835  {
836  /* get apparent position in the horizontal system
837  jd = MJD-time (UT)
838  ep2 : correction for apparent sidereal time in sec
839  lat, lng : geographic latitude and longitude (in radians)
840  ht : height above normal in meters.
841  solsys : = 1 if object is in solar system and parallax has to be
842  taken into account, 0 otherwise.
843  r = vector of celestial object. The unit of length of this vector
844  has to be in terms of the equatorial Earth radius (6378.14 km)
845  if solsys = 1, otherwise it's arbitrary.
846  azim : azimuth in radians (0 is to the North).
847  elev : elevation in radians
848  dist : distance (if solsys = 1; otherwise abs(r))
849  */
850  Vec3 s;
851 
852  // correct for topocentric position (parallax)
853  if (solsys) s = r - GeoPos(jd, ep2, lat, lng, ht);
854  else s = r;
855 
856  s = EquHor(jd, ep2, lat, lng, s);
857  s = carpol(s);
858  dist = s[0];
859  elev = s[2];
860  azim = M_PI - s[1];
861  }
862 
863 /*--------------------- Function AppRADec ----------------------------------*/
864 
865 void AppRADec (double jd, double ep2, double lat, double lng,
866  double azim, double elev, double& ra, double& dec)
867  {
868  /* get apparent position in the horizontal system
869  jd = MJD-time (UT)
870  ep2 : correction for apparent sidereal time in sec
871  lat, lng : geographic latitude and longitude (in radians)
872  azim : azimuth in radians (0 is to the North).
873  elev : elevation in radians
874  ra : Right Ascension (in radians)
875  dec : Declination (in radians)
876  */
877  Vec3 s;
878 
879  s[0] = 1.0;
880  s[1] = M_PI - azim;
881  s[2] = elev;
882  s = polcar(s);
883  s = HorEqu(jd, ep2, lat, lng, s);
884  s = carpol(s);
885  dec = s[2];
886  ra = s[1];
887  }
888 
889 /*------------------------- Refract ---------------------------------*/
890 
891 double Refract (double h, double p, double t)
892  {
893  /* Calculate atmospheric refraction with low precision
894  h = height (in radians) of object
895  p = presure (in millibars)
896  t = temperature (in degrees Celsius)
897 
898  RETURN: refraction angle in radians
899  */
900  double const raddeg = 180.0 / M_PI;
901  double r;
902 
903  r = h * raddeg;
904  r = (r + 7.31 / (r + 4.4)) / raddeg;
905  r = (0.28*p/(t+273.0)) * 0.0167 / tan(r);
906  r = r / raddeg;
907 
908  return r;
909  }
910 
911 /*---------------- Function eccanom --------------------------------------*/
912 
913 double eccanom (double man, double ecc)
914  {
915  /*
916  Solve Kepler equation for eccentric anomaly of elliptical orbits.
917  man : Mean Anomaly (in radians)
918  ecc : eccentricity
919  */
920  const double p2 = 2.0*M_PI;
921  const double eps = 1E-11;
922  const int maxit = 15;
923 
924  double m, e, f;
925  int i, mi;
926 
927  m=man/p2;
928  mi = int(m);
929  m=p2*(m-mi);
930  if (m < 0) m=m+p2;
931  if (ecc < 0.8) e=m;
932  else e=M_PI;
933  f = e - ecc*sin(e) - m;
934  i=0;
935 
936  while ((fabs(f) > eps) && (i < maxit))
937  {
938  e = e - f / (1.0 - ecc*cos(e));
939  f = e - ecc*sin(e) - m;
940  i = i + 1;
941  }
942 
943  return e;
944  }
945 
946 
947 /*---------------- Function hypanom --------------------------------------*/
948 
949 double hypanom (double mh, double ecc)
950  {
951  /*
952  Solve Kepler equation for eccentric anomaly of hyperbolic orbits.
953  mh : Mean Anomaly
954  ecc : eccentricity
955  */
956  const double eps = 1E-11;
957  const int maxit = 15;
958 
959  double h, f;
960  int i;
961 
962  h = log (2.0*fabs(mh)/ecc+1.8);
963  if (mh < 0.0) h = -h;
964  f = ecc * sinh(h) - h - mh;
965  i = 0;
966 
967  while ((fabs(f) > eps*(1.0+fabs(h+mh))) && (i < maxit))
968  {
969  h = h - f / (ecc*cosh(h) - 1.0);
970  f = ecc*sinh(h) - h - mh;
971  i = i + 1;
972  }
973 
974  return h;
975  }
976 
977 
978 /*------------------ Function ellip ------------------------------------*/
979 
980 void ellip (double gm, double t0, double t, double a, double ecc,
981  double m0, Vec3& r1, Vec3& v1)
982  {
983  /*
984  Calculate position r1 and velocity v1 of for elliptic orbit at time t.
985  gm : Gravitational Constant
986  t0 : epoch
987  a : semim-ajor axis
988  ecc : eccentricity
989  m0 : Mean Anomaly at epoch (in radians).
990  The units must be consistent with each other.
991  */
992  double m, e, fac, k, c, s;
993 
994  if (fabs(a) < 1e-60) a = 1e-60;
995  k = gm / a;
996  if (k >= 0) k = sqrt(k);
997  else k = 0; // just in case
998 
999  m = k * (t - t0) / a + m0;
1000  e = eccanom(m, ecc);
1001  fac = sqrt(1.0 - ecc*ecc);
1002  c = cos(e);
1003  s = sin(e);
1004  r1.assign (a*(c - ecc), a*fac*s, 0.0);
1005  m = 1.0 - ecc*c;
1006  v1.assign (-k*s/m, k*fac*c/m, 0.0);
1007  }
1008 
1009 /*---------------- Function hyperb --------------------------------------*/
1010 
1011 void hyperb (double gm, double t0, double t, double a, double ecc,
1012  Vec3& r1, Vec3& v1)
1013  {
1014  /*
1015  Calculate position r1 and velocity v1 of for hyperbolic orbit at time t.
1016  gm : Gravitational Constant
1017  t0 : time of perihelion passage
1018  a : semi-major axis
1019  ecc : eccentricity
1020  The units must be consistent with each other.
1021  */
1022  double mh, h, fac, k, c, s;
1023 
1024  a = fabs(a);
1025  if (a < 1e-60) a = 1e-60;
1026  k = gm / a;
1027  if (k >= 0) k = sqrt(k);
1028  else k = 0; // just in case
1029 
1030  mh = k * (t - t0) / a;
1031  h = hypanom (mh, ecc);
1032  fac = sqrt(ecc*ecc-1.0);
1033  c = cosh(h);
1034  s = sinh(h);
1035  r1.assign (a*(ecc-c), a*fac*s, 0.0);
1036  mh = ecc*c - 1.0;
1037  v1.assign (-k*s/mh, k*fac*c/mh, 0.0);
1038  }
1039 
1040 /*---------------- Function parab --------------------------------------*/
1041 
1042  void stumpff (double e2, double& c1, double& c2, double& c3)
1043  {
1044  /*
1045  Calculation of Stumpff functions c1=sin(e)/c,
1046  c2=(1-cos(e))/e^2 and c3=(e-sin(e))/e^3 for the
1047  argument e2=e^2 (e: eccentric anomaly)
1048  */
1049  const double eps=1E-12;
1050  double n, add;
1051 
1052  c1=0.0; c2=0.0; c3=0.0; add=1.0; n=1.0;
1053  do
1054  {
1055  c1=c1+add; add=add/(2.0*n);
1056  c2=c2+add; add=add/(2.0*n+1);
1057  c3=c3+add; add=-e2*add;
1058  n=n+1.0;
1059  } while (abs(add)>eps);
1060  }
1061 
1062 void parab (double gm, double t0, double t, double q, double ecc,
1063  Vec3& r1, Vec3& v1)
1064  {
1065  /*
1066  Calculate position r1 and velocity v1 of for parabolic
1067  (or near parabolic) orbit at time t using Stumpff's method.
1068  gm : Gravitational Constant
1069  t0 : time of perihelion passage
1070  q : perihelion distance
1071  ecc : eccentricity
1072  The units must be consistent with each other.
1073  */
1074  const double eps = 1E-9;
1075  const int maxit = 15;
1076 
1077  double e2, e20, fac, c1, c2, c3, k, tau, a, u, u2;
1078  double x, y;
1079  int i, fle;
1080 
1081  ecc = fabs(ecc);
1082  e2=0.0; fac=0.5*ecc; i=0;
1083 
1084  q = fabs(q);
1085  if (q < 1e-40) q = 1e-40;
1086  k = gm / (q*(1+ecc));
1087 
1088  if (k >= 0) k = sqrt(k);
1089  else k = 0; // just in case
1090 
1091  tau = gm / (q*q*q);
1092  if (tau >= 0) tau= 1.5*sqrt(tau)*(t-t0);
1093  else tau = 0;
1094 
1095  do
1096  {
1097  i = i + 1;
1098  e20 = e2;
1099  if (fac < 0.0) a = -1.0;
1100  else a = tau * sqrt(fac);
1101  a = sqrt(a*a+1.0)+a;
1102  if (a > 0.0) a = exp(log(a)/3.0);
1103  if (a == 0.0) u = 0.0;
1104  else u = a - 1.0/a;
1105  u2 = u*u;
1106  if (fac != 0) e2 = u2*(1.0-ecc)/fac;
1107  else e2 = 1;
1108  stumpff (e2, c1, c2, c3);
1109  fac = 3.0*ecc*c3;
1110  fle = (abs(e2-e20) < eps) || (i > maxit);
1111  } while (!fle);
1112 
1113  if (fac != 0)
1114  {
1115  tau = q*(1.0+u2*c2*ecc/fac);
1116  x = q*(1.0-u2*c2/fac);
1117  y = (1.0+ecc)/fac;
1118  if (y >= 0) y = q*sqrt(y)*u*c1;
1119  else y = 0;
1120  r1.assign (x, y, 0.0);
1121  v1.assign (-k*y/tau, k*(x/tau+ecc), 0.0);
1122  }
1123  else // just in case
1124  {
1125  r1.assign (0, 0, 0);
1126  v1.assign (0, 0, 0);
1127  };
1128  }
1129 
1130 /*---------------- Function kepler --------------------------------------*/
1131 
1132 void kepler (double gm, double t0, double t, double m0, double a, double ecc,
1133  double ran, double aper, double inc, Vec3& r1, Vec3& v1)
1134  {
1135  /*
1136  Calculate position r1 and velocity v1 of body at time t as an
1137  undisturbed 2-body Kepler problem.
1138 
1139  gm : Gravitational Constant
1140  t0 : time of perihelion passage (hyperbolic or near-parabolic)
1141  epoch of m0 for elliptical orbits.
1142  m0 : Mean Anomaly at epoch for elliptical orbits in degrees in the
1143  range 0 to 360. If m0 < 0 the calculations will be done as
1144  near-parabolic. Set m0 = 0 for hyperbolic orbits.
1145  a : semi-major axis for elliptical (positive) or hyperbolic orbits
1146  (negative). If m0 < 0, a signifies perihelion distance q instead
1147  ecc : eccentricity
1148  ran : Right Ascension of Ascending Node (in degrees)
1149  aper : Argument of Perihelion (in degrees)
1150  inc : Inclination (in degrees)
1151 
1152  The units must be consistent with each other.
1153 
1154  NOTE : If ecc = 1 the orbit will always be calculated as a parabolic one.
1155  If ecc < 1 (elliptical orbits) two possibilities exist:
1156  If m0 >= 0 the calculation will be done in the classical
1157  way with m0 signifying the mean anomaly at time t0 and
1158  a the semi-major axis.
1159  If m0 < 0 the orbit will be calculated as a near-parabolic
1160  orbit with a signifying the perihelion distance at time t0.
1161  (This latter method will be useful for high eccentricity
1162  comet orbits for which typically q is given instead of a)
1163  If ecc > 1 (hyperbolic orbits) two possibilities exist:
1164  If m0 >= 0 (the value is ignored for hyperbolic orbits)
1165  a classical hyperbolic calculation is performed with a
1166  signifying the semi-major axis (negative for hyperbolas).
1167  If m0 < 0 the orbit will be calculated as a near-parabolic
1168  orbit with a signifying the perihelion distance. (This can
1169  be useful for comet orbits which would rarely have an
1170  eccentricity >> 1)
1171  In either case t0 signifies the time of perihelion passage.
1172  */
1173  const double dgrd = M_PI / 180.0;
1174  enum {ell, par, hyp} kepc;
1175  Mat3 m1,m2;
1176 
1177  kepc = ell; // just to keep the compiler happy
1178 
1179  // convert into radians
1180  m0 *= dgrd;
1181  ran *= dgrd;
1182  aper *= dgrd;
1183  inc *= dgrd;
1184 
1185  // calculate the position and velocity within the orbit plane
1186  if (ecc == 1.0) kepc = par;
1187  if (ecc < 1.0)
1188  {
1189  if (m0 < 0.0) kepc = par;
1190  else kepc = ell;
1191  }
1192  if (ecc > 1.0)
1193  {
1194  if (m0 < 0.0) kepc = par;
1195  else kepc = hyp;
1196  }
1197 
1198  switch (kepc)
1199  {
1200  case ell : ellip (gm, t0, t, a, ecc, m0, r1, v1); break;
1201  case par : parab (gm, t0, t, a, ecc, r1, v1); break;
1202  case hyp : hyperb (gm, t0, t, a, ecc, r1, v1); break;
1203  }
1204 
1205  // convert into reference plane
1206  m1 = zrot (-aper);
1207  m2 = xrot (-inc);
1208  m1 *= m2;
1209  m2 = zrot (-ran);
1210  m2 = m2 * m1;
1211 
1212  r1 = mxvct (m2, r1);
1213  v1 = mxvct (m2, v1);
1214  }
1215 
1216 /*---------------- Function oscelm --------------------------------------*/
1217 
1218 void oscelm (double gm, double t, Vec3& r1, Vec3& v1,
1219  double& t0, double& m0, double& a, double& ecc,
1220  double& ran, double& aper, double& inc)
1221  {
1222  /*
1223  Get the osculating Kepler elements at epoch t from position r1 and
1224  velocity v1 of body.
1225 
1226  gm : Gravitational Constant. If set negative, its absolute value will
1227  be taken as gm and the perihelion distance will be returned
1228  instead of the semimajor axis.
1229  t0 : time of perihelion passage
1230  m0 : Mean Anomaly at epoch for elliptical orbits in degrees in the
1231  range 0 to 360. Set m0 = 0 for hyperbolic orbits.
1232  m0 will be set to -1 if perihelion distance is to be returned
1233  instead of a.
1234  a : semi-major axis for elliptical (positive) or hyperbolic orbits
1235  (negative). If m0 < 0, a signifies perihelion distance q instead.
1236  If gm is negative, the perihelion distance q will be returned
1237  ecc : eccentricity
1238  ran : Right Ascension of Ascending Node (in degrees)
1239  aper : Argument of Perihelion (in degrees)
1240  inc : Inclination (in degrees)
1241 
1242  The units must be consistent with each other.
1243  */
1244  const double rddg = 180.0/M_PI;
1245  const double p2 = 2.0*M_PI;
1246 
1247  Vec3 c;
1248  double cabs, u, cu, su, p, r;
1249  int pflg; // 1, if perihelion distance to be returned
1250 
1251  pflg = 0;
1252  if (gm < 0.0)
1253  {
1254  gm = -gm;
1255  pflg = 1;
1256  };
1257 
1258  if (gm < 1e-60) gm = 1e-60;
1259  c = r1 * v1;
1260  cabs = abs(c);
1261  if (fabs(cabs) < 1e-40) cabs = 1e-40;
1262  ran = atan20 (c[0], -c[1]);
1263  inc = c[2]/cabs;
1264  if (fabs(inc) <= 1.0) inc = acos (inc);
1265  else inc = 0;
1266 
1267  r = abs(r1);
1268  if (fabs(r) < 1e-40) r = 1e-40;
1269  su = sin(inc);
1270  if (su != 0.0) su = r1[2]/su;
1271  cu = r1[0]*cos(ran)+r1[1]*sin(ran);
1272  u = atan20 (su, cu); // argument of latitude
1273 
1274  a = abs(v1);
1275  a = 2.0/r - a*a/gm;
1276  if (fabs(a) < 1.0E-30) ecc = 1.0; // parabolic
1277  else
1278  {
1279  a = 1.0/a; // semimajor axis
1280  ecc = 0;
1281  };
1282 
1283  p = cabs*cabs/gm; // semilatum rectum
1284 
1285  if (ecc == 1.0)
1286  {
1287  p = p / 2.0; // perihelion distance
1288  a = 2*p;
1289  }
1290  else
1291  {
1292  ecc = 1.0 - p/a;
1293  if (ecc >= 0) ecc = sqrt(ecc);
1294  else ecc = 0;
1295  p = p / (1.0 + ecc); // perihelion distance
1296  }
1297 
1298  if (fabs(a) > 1e-60) cu = (1.0 - r/a);
1299  else cu = 0; // just in case
1300  su = dot(r1,v1) / sqrt(fabs(a)*gm);
1301  if (ecc < 1.0)
1302  {
1303  m0 = atan20(su,cu); // eccentric anomaly
1304  su = sin(m0);
1305  cu = cos(m0);
1306  aper = 1.0-ecc*ecc;
1307  if (aper >= 0) aper = atan20(sqrt(aper)*su,(cu-ecc)); // true anomaly
1308  m0 = m0 - ecc*su; // mean anomaly
1309  }
1310  else if (ecc > 1.0)
1311  {
1312  su = su / ecc;
1313  m0 = su + sqrt(su*su + 1.0);
1314  if (m0 >= 0) m0 = log(m0); // = asinh (su); hyperbolic anomaly
1315  aper = (ecc+1.0)/(ecc-1.0);
1316  if (aper >= 0) aper = 2.0 * atan(sqrt(aper)*tanh(m0/2.0));
1317  m0 = ecc*su-m0;
1318  }
1319 
1320  if (ecc != 1.0)
1321  {
1322  aper = u - aper; // argument of perihelion
1323  u = fabs(a);
1324  t0 = u/gm;
1325  if (t0 >= 0) t0 = t - u*sqrt(t0)*m0; // time of perihelion transit
1326  else t0 = t;
1327  }
1328  else // parabolic
1329  {
1330  pflg = 1;
1331 
1332  aper = 2.0 * atan(su); // true anomaly
1333  aper = u - aper; // argument of perihelion
1334  t0 = 2.0*p*p*p/gm;
1335  if (t0 >= 0) t0 = t - sqrt(t0) * (su*su*su/3.0 + su);
1336  else t0 = t;
1337  }
1338 
1339  if (m0 < 0.0) m0 = m0 + p2;
1340  if (ran < 0.0) ran = ran + p2;
1341  if (aper < 0.0) aper = aper + p2;
1342  m0 = rddg * m0;
1343  ran = rddg * ran;
1344  aper = rddg * aper;
1345  inc = rddg * inc;
1346 
1347  if (ecc > 1.0) m0 = 0.0;
1348  if (pflg)
1349  {
1350  a = p;
1351  m0 = -1.0;
1352  }
1353  }
1354 
1355 /*-------------------- QuickSun --------------------------------------*/
1356 
1357 Vec3 QuickSun (double t) // low precision position of the Sun at time t
1358  {
1359  /* Low precision position of the Sun at time t given in
1360  Julian centuries since J2000.
1361  Valid only between 1950 and 2050.
1362 
1363  Returns the position vector (in A.U.) in ecliptic coordinates
1364 
1365  */
1366  double n, g, l;
1367  Vec3 rs;
1368 
1369  n = 36525.0 * t;
1370  l = 280.460 + 0.9856474*n; // mean longitude
1371  g = 357.528 + 0.9856003*n; // mean anomaly
1372  l = 2.0*M_PI * frac(l/360.0);
1373  g = 2.0*M_PI * frac(g/360.0);
1374  l = l + (1.915*sin(g) + 0.020*sin(2.0*g)) * 1.74532925199e-2; // ecl.long.
1375  g = 1.00014 - 0.01671*cos(g) - 0.00014*cos(2.0*g); // radius
1376  rs[0] = g * cos(l);
1377  rs[1] = g * sin(l);
1378  rs[2] = 0;
1379 
1380  return rs;
1381  }
1382 
1383 /*-------------------- Class Sun200 --------------------------------------*/
1384  /*
1385  Ecliptic coordinates (in A.U.) and velocity (in A.U./day)
1386  of the sun for Equinox of Date given in Julian centuries since J2000.
1387  ======================================================================
1388  */
1389 Sun200::Sun200 ()
1390  { }
1391 
1392 Vec3 Sun200::position (double t) // position of the Sun at time t
1393  {
1394  Vec3 rs, vs;
1395 
1396  state (t, rs, vs);
1397 
1398  return rs;
1399  }
1400 
1401 void Sun200::state (double t, Vec3& rs, Vec3& vs)
1402  {
1403  /* State vector rs (position) and vs (velocity) of the Sun in
1404  ecliptic of date coordinates at time t (in Julian Centruries
1405  since J2000).
1406  */
1407  const double p2 = 2.0 * M_PI;
1408 
1409  double l, b, r;
1410  int i;
1411 
1412  tt = t;
1413 
1414  dl = 0.0; dr = 0.0; db = 0.0;
1415  m2 = p2 * frac(0.1387306 + 162.5485917*t);
1416  m3 = p2 * frac(0.9931266 + 99.9973604*t);
1417  m4 = p2 * frac(0.0543250 + 53.1666028*t);
1418  m5 = p2 * frac(0.0551750 + 8.4293972*t);
1419  m6 = p2 * frac(0.8816500 + 3.3938722*t);
1420  d = p2 * frac(0.8274 + 1236.8531*t);
1421  a= p2 * frac(0.3749 + 1325.5524*t);
1422  uu = p2 * frac(0.2591 + 1342.2278*t);
1423 
1424  c3[1] = 1.0; s3[1] = 0.0;
1425  c3[2] = cos(m3); s3[2] = sin(m3);
1426  c3[0] = c3[2]; s3[0] = -s3[2];
1427 
1428  for (i=3; i<9; ++i)
1429  addthe(c3[i-1],s3[i-1],c3[2],s3[2],c3[i],s3[i]);
1430  pertven(); pertmar(); pertjup(); pertsat(); pertmoo();
1431 
1432  dl = dl + 6.4 * sin(p2*(0.6983 + 0.0561*t))
1433  + 1.87 * sin(p2*(0.5764 + 0.4174*t))
1434  + 0.27 * sin(p2*(0.4189 + 0.3306*t))
1435  + 0.20 * sin(p2*(0.3581 + 2.4814*t));
1436  l = p2 * frac(0.7859453 + m3/p2 +((6191.2+1.1*t)*t+dl)/1296.0E3);
1437  r = 1.0001398 - 0.0000007*t + dr*1E-6;
1438  b = db * 4.8481368111E-6;
1439 
1440  cl= cos(l); sl=sin(l); cb=cos(b); sb=sin(b);
1441  rs[0]=r*cl*cb; rs[1]=r*sl*cb; rs[2]=r*sb;
1442 
1443  /* velocity calculation
1444  gms=2.9591202986e-4 AU^3/d^2
1445  e=0.0167086
1446  sqrt(gms/a)=1.7202085e-2 AU/d
1447  sqrt(gms/a)*sqrt(1-e^2)=1.71996836e-2 AU/d
1448  nu=M + 2e*sin(M) = M + 0.0334172 * sin(M) */
1449 
1450  uu = m3 + 0.0334172*sin(m3); // eccentric Anomaly E
1451  d = cos(uu);
1452  uu = sin(uu);
1453  a = 1.0 - 0.0167086*d;
1454  vs[0] = -1.7202085e-2*uu/a; // velocity in orbit plane
1455  vs[1] = 1.71996836e-2*d/a;
1456  uu = atan2 (0.9998604*uu, (d-0.0167086)); // true anomaly
1457  d = cos(uu);
1458  uu = sin(uu);
1459  dr = d*vs[0]+uu*vs[1];
1460  dl = (d*vs[1]-uu*vs[0]) / r;
1461 
1462  vs[0] = dr*cl*cb - dl*r*sl*cb;
1463  vs[1] = dr*sl*cb + dl*r*cl*cb;
1464  vs[2] = dr*sb;
1465  }
1466 
1467 void Sun200::addthe (double c1, double s1, double c2, double s2,
1468  double& cc, double& ss)
1469  {
1470  cc=c1*c2-s1*s2;
1471  ss=s1*c2+c1*s2;
1472  }
1473 
1474 void Sun200::term (int i1, int i, int it, double dlc, double dls, double drc,
1475  double drs, double dbc, double dbs)
1476  {
1477  if (it == 0) addthe (c3[i1+1],s3[i1+1],c[i+8],s[i+8],u,v);
1478  else
1479  {
1480  u=u*tt;
1481  v=v*tt;
1482  }
1483  dl = dl + dlc*u + dls*v;
1484  dr = dr + drc*u + drs*v;
1485  db = db + dbc*u + dbs*v;
1486  }
1487 
1488 void Sun200::pertven() // Kepler terms and perturbations by Venus
1489  {
1490  int i;
1491 
1492  c[8]=1.0; s[8]=0.0; c[7]=cos(m2); s[7]=-sin(m2);
1493  for (i=7; i>2; i--)
1494  addthe(c[i],s[i],c[7],s[7],c[i-1],s[i-1]);
1495 
1496  term (1, 0, 0, -0.22, 6892.76, -16707.37, -0.54, 0.0, 0.0);
1497  term (1, 0, 1, -0.06, -17.35, 42.04, -0.15, 0.0, 0.0);
1498  term (1, 0, 2, -0.01, -0.05, 0.13, -0.02, 0.0, 0.0);
1499  term (2, 0, 0, 0.0, 71.98, -139.57, 0.0, 0.0, 0.0);
1500  term (2, 0, 1, 0.0, -0.36, 0.7, 0.0, 0.0, 0.0);
1501  term (3, 0, 0, 0.0, 1.04, -1.75, 0.0, 0.0, 0.0);
1502  term (0, -1, 0, 0.03, -0.07, -0.16, -0.07, 0.02, -0.02);
1503  term (1, -1, 0, 2.35, -4.23, -4.75, -2.64, 0.0, 0.0);
1504  term (1, -2, 0, -0.1, 0.06, 0.12, 0.2, 0.02, 0.0);
1505  term (2, -1, 0, -0.06, -0.03, 0.2, -0.01, 0.01, -0.09);
1506  term (2, -2, 0, -4.7, 2.9, 8.28, 13.42, 0.01, -0.01);
1507  term (3, -2, 0, 1.8, -1.74, -1.44, -1.57, 0.04, -0.06);
1508  term (3, -3, 0, -0.67, 0.03, 0.11, 2.43, 0.01, 0.0);
1509  term (4, -2, 0, 0.03, -0.03, 0.1, 0.09, 0.01, -0.01);
1510  term (4, -3, 0, 1.51, -0.4, -0.88, -3.36, 0.18, -0.1);
1511  term (4, -4, 0, -0.19, -0.09, -0.38, 0.77, 0.0, 0.0);
1512  term (5, -3, 0, 0.76, -0.68, 0.3, 0.37, 0.01, 0.0);
1513  term (5, -4, 0, -0.14, -0.04, -0.11, 0.43, -0.03, 0.0);
1514  term (5, -5, 0, -0.05, -0.07, -0.31, 0.21, 0.0, 0.0);
1515  term (6, -4, 0, 0.15, -0.04, -0.06, -0.21, 0.01, 0.0);
1516  term (6, -5, 0, -0.03, -0.03, -0.09, 0.09, -0.01, 0.0);
1517  term (6, -6, 0, 0.0, -0.04, -0.18, 0.02, 0.0, 0.0);
1518  term (7, -5, 0, -0.12, -0.03, -0.08, 0.31, -0.02, -0.01);
1519  }
1520 
1521 void Sun200::pertmar() // Kepler terms and perturbations by Mars
1522  {
1523  int i;
1524 
1525  c[7] = cos(m4); s[7] = -sin(m4);
1526  for (i=7; i>0; i--)
1527  addthe(c[i],s[i],c[7],s[7],c[i-1],s[i-1]);
1528  term (1, -1, 0, -0.22, 0.17, -0.21, -0.27, 0.0, 0.0);
1529  term (1, -2, 0, -1.66, 0.62, 0.16, 0.28, 0.0, 0.0);
1530  term (2, -2, 0, 1.96, 0.57, -1.32, 4.55, 0.0, 0.01);
1531  term (2, -3, 0, 0.4, 0.15, -0.17, 0.46, 0.0, 0.0);
1532  term (2, -4, 0, 0.53, 0.26, 0.09, -0.22, 0.0, 0.0);
1533  term (3, -3, 0, 0.05, 0.12, -0.35, 0.15, 0.0, 0.0);
1534  term (3, -4, 0, -0.13, -0.48, 1.06, -0.29, 0.01, 0.0);
1535  term (3, -5, 0, -0.04, -0.2, 0.2, -0.04, 0.0, 0.0);
1536  term (4, -4, 0, 0.0, -0.03, 0.1, 0.04, 0.0, 0.0);
1537  term (4, -5, 0, 0.05, -0.07, 0.2, 0.14, 0.0, 00);
1538  term (4, -6, 0, -0.1, 0.11, -0.23, -0.22, 0.0, 0.0);
1539  term (5, -7, 0, -0.05, 0.0, 0.01, -0.14, 0.0, 0.0);
1540  term (5, -8, 0, 0.05, 0.01, -0.02, 0.1, 0.0, 0.0);
1541  }
1542 
1543 void Sun200::pertjup() // Kepler terms and perturbations by Jupiter
1544  {
1545  int i;
1546 
1547  c[7] = cos(m5); s[7] = -sin(m5);
1548  for (i=7; i>4; i--)
1549  addthe(c[i],s[i],c[7],s[7],c[i-1],s[i-1]);
1550  term (1, -1, 0, 0.01, 0.07, 0.18, -0.02, 0.0, -0.02);
1551  term (0, -1, 0, -0.31, 2.58, 0.52, 0.34, 0.02, 0.0);
1552  term (1, -1, 0, -7.21, -0.06, 0.13, -16.27, 0.0, -0.02);
1553  term (1, -2, 0, -0.54, -1.52, 3.09, -1.12, 0.01, -0.17);
1554  term (1, -3, 0, -0.03, -0.21, 0.38, -0.06, 0.0, -0.02);
1555  term (2, -1, 0, -0.16, 0.05, -0.18, -0.31, 0.01, 0.0);
1556  term (2, -2, 0, 0.14, -2.73, 9.23, 0.48, 0.0, 0.0);
1557  term (2, -3, 0, 0.07, -0.55, 1.83, 0.25, 0.01, 0.0);
1558  term (2, -4, 0, 0.02, -0.08, 0.25, 0.06, 0.0, 0.0);
1559  term (3, -2, 0, 0.01, -0.07, 0.16, 0.04, 0.0, 0.0);
1560  term (3, -3, 0, -0.16, -0.03, 0.08, -0.64, 0.0, 0.0);
1561  term (3, -4, 0, -0.04, -0.01, 0.03, -0.17, 0.0, 0.0);
1562  }
1563 
1564 void Sun200::pertsat() // Kepler terms and perturbations by Saturn
1565  {
1566  c[7] = cos(m6); s[7] = -sin(m6);
1567  addthe(c[7],s[7],c[7],s[7],c[6],s[6]);
1568  term (0, -1, 0, 0.0, 0.32, 0.01, 0.0, 0.0, 0.0);
1569  term (1, -1, 0, -0.08, -0.41, 0.97, -0.18, 0.0, -0.01);
1570  term (1, -2, 0, 0.04, 0.1, -0.23, 0.1, 0.0, 0.0);
1571  term (2, -2, 0, 0.04, 0.1, -0.35, 0.13, 0.0, 0.0);
1572  }
1573 
1574 void Sun200::pertmoo() // corrections for Earth-Moon center of gravity
1575  {
1576  dl = dl + 6.45*sin(d) - 0.42*sin(d-a) + 0.18*sin(d+a) + 0.17*sin(d-m3) - 0.06*sin(d+m3);
1577  dr = dr + 30.76*cos(d) - 3.06*cos(d-a) + 0.85*cos(d+a) - 0.58*cos(d+m3) + 0.57*cos(d-m3);
1578  db = db + 0.576*sin(uu);
1579  }
1580 
1581 /*--------------------- Class moon200 ----------------------------------*/
1582 
1583  /*
1584  Position vector (in Earth radii) of the Moon referred to Ecliptic
1585  for Equinox of Date.
1586  t is the time in Julian centuries since J2000.
1587  */
1588 
1589 Moon200::Moon200 ()
1590  { }
1591 
1592 Vec3 Moon200::position (double t) // position of the Moon at time t
1593  {
1594  Vec3 rm;
1595 
1596  const double arc = 206264.81; // 3600*180/pi = ''/r
1597  const double pi2 = M_PI * 2.0;
1598 
1599  double lambda, beta, r, fac;
1600 
1601  minit(t);
1602  solar1(); solar2(); solar3(); solarn(n); planetary(t);
1603 
1604  lambda = pi2 * frac ((l0+dlam/arc) / pi2);
1605  if (lambda < 0) lambda = lambda + pi2;
1606  s = f + ds / arc;
1607 
1608  fac = 1.000002708 + 139.978 * dgam;
1609  beta = (fac*(18518.511+1.189+gam1c)*sin(s)-6.24*sin(3*s)+n);
1610  beta = beta * 4.8481368111e-6;
1611  sinpi = sinpi * 0.999953253;
1612  r = arc / sinpi;
1613 
1614  rm.assign (r, lambda, beta);
1615  rm = polcar(rm);
1616 
1617  return rm;
1618  }
1619 
1620 
1621 void Moon200::addthe (double c1, double s1, double c2, double s2,
1622  double& c, double& s)
1623  {
1624  c=c1*c2-s1*s2;
1625  s=s1*c2+c1*s2;
1626  }
1627 
1628 double Moon200::sinus (double phi)
1629  {
1630  /* sin(phi) in units of 1r = 360ø */
1631  return sin (2.0*M_PI * frac(phi));
1632  }
1633 
1634 void Moon200::long_periodic (double t)
1635  {
1636  /* long periodic changes of mean elements */
1637 
1638  double s1, s2, s3, s4, s5, s6, s7;
1639 
1640  s1 = sinus (0.19833 + 0.05611*t);
1641  s2 = sinus (0.27869 + 0.04508*t);
1642  s3 = sinus (0.16827 - 0.36903*t);
1643  s4 = sinus (0.34734 - 5.37261*t);
1644  s5 = sinus (0.10498 - 5.37899*t);
1645  s6 = sinus (0.42681 - 0.41855*t);
1646  s7 = sinus (0.14943 - 5.37511*t);
1647  dl0 = 0.84*s1 + 0.31*s2 + 14.27*s3 + 7.26*s4 + 0.28*s5 + 0.24*s6;
1648  dl = 2.94*s1 + 0.31*s2 + 14.27*s3 + 9.34*s4 + 1.12*s5 + 0.83*s6;
1649  dls = -6.4*s1 - 1.89*s6;
1650  df = 0.21*s1+0.31*s2+14.27*s3-88.7*s4-15.3*s5+0.24*s6-1.86*s7;
1651  dd = dl0 - dls;
1652  dgam = -3332E-9 * sinus (0.59734 - 5.37261*t)
1653  -539E-9 * sinus (0.35498 - 5.37899*t)
1654  -64E-9 * sinus (0.39943 - 5.37511*t);
1655  }
1656 
1657 void Moon200::minit(double t)
1658  {
1659  /* calculate mean elements l (moon), F (node distance)
1660  l' (sun) and D (moon's elongation) */
1661 
1662  const double arc = 206264.81; // 3600*180/pi = ''/r
1663  const double pi2 = M_PI * 2.0;
1664 
1665  int i, j, max;
1666  double t2, arg, fac;
1667 
1668  max = 0; // just to keep the compiler happy
1669  arg = 0;
1670  fac = 0;
1671 
1672  t2 = t * t;
1673  dlam = 0; ds = 0; gam1c = 0; sinpi = 3422.7;
1674  long_periodic (t);
1675  l0 = pi2*frac(0.60643382+1336.85522467*t-0.00000313*t2) + dl0/arc;
1676  l = pi2*frac(0.37489701+1325.55240982*t+0.00002565*t2) + dl/arc;
1677  ls = pi2*frac(0.99312619+99.99735956*t-0.00000044*t2) + dls/arc;
1678  f = pi2*frac(0.25909118+1342.22782980*t-0.00000892*t2) + df/arc;
1679  d = pi2*frac(0.82736186+1236.85308708*t-0.00000397*t2) + dd/arc;
1680  for (i=0; i<4; ++i)
1681  {
1682  switch (i)
1683  {
1684  case 0: arg=l; max=4; fac=1.000002208; break;
1685  case 1: arg=ls; max=3; fac=0.997504612-0.002495388*t; break;
1686  case 2: arg=f; max=4; fac=1.000002708+139.978*dgam; break;
1687  case 3: arg=d; max=6; fac=1.0; break;
1688  }
1689  co[6][i] = 1.0; co[7][i] = cos (arg) * fac;
1690  si[6][i] = 0.0; si[7][i] = sin (arg) * fac;
1691  for (j = 2; j <= max; ++j)
1692  addthe(co[j+5][i],si[j+5][i],co[7][i],si[7][i],co[j+6][i],si[j+6][i]);
1693  for (j = 1; j <= max; ++j)
1694  {
1695  co[6-j][i]=co[j+6][i];
1696  si[6-j][i]=-si[j+6][i];
1697  }
1698  }
1699  }
1700 
1701 void Moon200::term (int p, int q, int r, int s, double& x, double& y)
1702  {
1703  // calculate x=cos(p*arg1+q*arg2...); y=sin(p*arg1+q*arg2...)
1704  int i[4];
1705  int k;
1706 
1707  i[0]=p; i[1]=q; i[2]=r; i[3]=s; x=1.0; y=0.0;
1708  for (k=0; k<4; k++)
1709  if (i[k]!=0) addthe(x,y,co[i[k]+6][k],si[i[k]+6][k],x,y);
1710  }
1711 
1712 void Moon200::addsol(double coeffl, double coeffs, double coeffg,
1713  double coeffp, int p, int q, int r, int s)
1714  {
1715  double x, y;
1716  term (p,q,r,s,x,y);
1717  dlam = dlam + coeffl*y; ds = ds + coeffs*y;
1718  gam1c = gam1c + coeffg*x; sinpi = sinpi + coeffp*x;
1719  }
1720 
1721 void Moon200::solar1()
1722  {
1723  addsol ( 13.902, 14.06, -0.001, 0.2607, 0, 0, 0, 4);
1724  addsol ( 0.403, -4.01, 0.394, 0.0023, 0, 0, 0, 3);
1725  addsol ( 2369.912, 2373.36, 0.601, 28.2333, 0, 0, 0, 2);
1726  addsol ( -125.154, -112.79, -0.725, -0.9781, 0, 0, 0, 1);
1727  addsol ( 1.979, 6.98, -0.445, 0.0433, 1, 0, 0, 4);
1728  addsol (191.953, 192.72, 0.029, 3.0861, 1, 0, 0, 2);
1729  addsol (-8.466, -13.51, 0.455, -0.1093, 1, 0, 0, 1);
1730  addsol (22639.5, 22609.07, 0.079, 186.5398, 1, 0, 0, 0);
1731  addsol (18.609, 3.59, -0.094, 0.0118, 1, 0, 0, -1);
1732  addsol (-4586.465, -4578.13, -0.077, 34.3117, 1, 0, 0, -2);
1733  addsol (3.215, 5.44, 0.192, -0.0386, 1, 0, 0, -3);
1734  addsol (-38.428, -38.64, 0.001, 0.6008, 1, 0, 0, -4);
1735  addsol (-0.393, -1.43, -0.092, 0.0086, 1, 0, 0, -6);
1736  addsol (-0.289, -1.59, 0.123, -0.0053, 0, 1, 0, 4);
1737  addsol (-24.420, -25.1, 0.04, -0.3, 0, 1, 0, 2);
1738  addsol (18.023, 17.93, 0.007, 0.1494, 0, 1, 0, 1);
1739  addsol (-668.146, -126.98, -1.302, -0.3997, 0, 1, 0, 0);
1740  addsol (0.56, 0.32, -0.001, -0.0037, 0, 1, 0, -1);
1741  addsol (-165.145, -165.06, 0.054, 1.9178, 0, 1, 0, -2);
1742  addsol (-1.877, -6.46, -0.416, 0.0339, 0, 1, 0, -4);
1743  addsol (0.213, 1.02, -0.074, 0.0054, 2, 0, 0, 4);
1744  addsol (14.387, 14.78, -0.017, 0.2833, 2, 0, 0, 2);
1745  addsol (-0.586, -1.2, 0.054, -0.01, 2, 0, 0, 1);
1746  addsol (769.016, 767.96, 0.107, 10.1657, 2, 0, 0, 0);
1747  addsol (1.75, 2.01, -0.018, 0.0155, 2, 0, 0, -1);
1748  addsol (-211.656, -152.53, 5.679, -0.3039, 2, 0, 0, -2);
1749  addsol (1.225, 0.91, -0.03, -0.0088, 2, 0, 0, -3);
1750  addsol (-30.773, -34.07, -0.308, 0.3722, 2, 0, 0, -4);
1751  addsol (-0.57, -1.4, -0.074, 0.0109, 2, 0, 0, -6);
1752  addsol (-2.921, -11.75, 0.787, -0.0484, 1, 1, 0, 2);
1753  addsol (1.267, 1.52, -0.022, 0.0164, 1, 1, 0, 1);
1754  addsol (-109.673, -115.18, 0.461, -0.949, 1, 1, 0, 0);
1755  addsol (-205.962, -182.36, 2.056, 1.4437, 1, 1, 0, -2);
1756  addsol (0.233, 0.36, 0.012, -0.0025, 1, 1, 0, -3);
1757  addsol (-4.391, -9.66, -0.471, 0.0673, 1, 1, 0, -4);
1758  }
1759 
1760 void Moon200::solar2()
1761  {
1762  addsol (0.283, 1.53, -0.111, 0.006, 1, -1, 0, 4);
1763  addsol (14.577, 31.7, -1.54, 0.2302, 1, -1, 0, 2);
1764  addsol (147.687, 138.76, 0.679, 1.1528, 1, -1, 0, 0);
1765  addsol (-1.089, 0.55, 0.021, 0.0, 1, -1, 0, -1);
1766  addsol (28.475, 23.59, -0.443, -0.2257, 1, -1, 0, -2);
1767  addsol (-0.276, -0.38, -0.006, -0.0036, 1, -1, 0, -3);
1768  addsol (0.636, 2.27, 0.146, -0.0102, 1, -1, 0, -4);
1769  addsol (-0.189, -1.68, 0.131, -0.0028, 0, 2, 0, 2);
1770  addsol (-7.486, -0.66, -0.037, -0.0086, 0, 2, 0, 0);
1771  addsol (-8.096, -16.35, -0.74, 0.0918, 0, 2, 0, -2);
1772  addsol (-5.741, -0.04, 0.0, -0.0009, 0, 0, 2, 2);
1773  addsol (0.255, 0.0, 0.0, 0.0, 0, 0, 2, 1);
1774  addsol (-411.608, -0.2, 0.0, -0.0124, 0, 0, 2, 0);
1775  addsol (0.584, 0.84, 0.0, 0.0071, 0, 0, 2, -1);
1776  addsol (-55.173, -52.14, 0.0, -0.1052, 0, 0, 2, -2);
1777  addsol (0.254, 0.25, 0.0, -0.0017, 0, 0, 2, -3);
1778  addsol (0.025, -1.67, 0.0, 0.0031, 0, 0, 2, -4);
1779  addsol (1.06, 2.96, -0.166, 0.0243, 3, 0, 0, 2);
1780  addsol (36.124, 50.64, -1.3, 0.6215, 3, 0, 0, 0);
1781  addsol (-13.193, -16.4, 0.258, -0.1187, 3, 0, 0, -2);
1782  addsol (-1.187, -0.74, 0.042, 0.0074, 3, 0, 0, -4);
1783  addsol (-0.293, -0.31, -0.002, 0.0046, 3, 0, 0, -6);
1784  addsol (-0.29, -1.45, 0.116, -0.0051, 2, 1, 0, 2);
1785  addsol (-7.649, -10.56, 0.259, -0.1038, 2, 1, 0, 0);
1786  addsol (-8.627, -7.59, 0.078, -0.0192, 2, 1, 0, -2);
1787  addsol (-2.74, -2.54, 0.022, 0.0324, 2, 1, 0, -4);
1788  addsol (1.181, 3.32, -0.212, 0.0213, 2, -1, 0, 2);
1789  addsol (9.703, 11.67, -0.151, 0.1268, 2, -1, 0, 0);
1790  addsol (-0.352, -0.37, 0.001, -0.0028, 2, -1, 0, -1);
1791  addsol (-2.494, -1.17, -0.003, -0.0017, 2, -1, 0, -2);
1792  addsol (0.36, 0.2, -0.012, -0.0043, 2, -1, 0, -4);
1793  addsol (-1.167, -1.25, 0.008, -0.0106, 1, 2, 0, 0);
1794  addsol (-7.412, -6.12, 0.117, 0.0484, 1, 2, 0, -2);
1795  addsol (-0.311, -0.65, -0.032, 0.0044, 1, 2, 0, -4);
1796  addsol (0.757, 1.82, -0.105, 0.0112, 1, -2, 0, 2);
1797  addsol (2.58, 2.32, 0.027, 0.0196, 1, -2, 0, 0);
1798  addsol (2.533, 2.4, -0.014, -0.0212, 1, -2, 0, -2);
1799  addsol (-0.344, -0.57, -0.025, 0.0036, 0, 3, 0, -2);
1800  addsol (-0.992, -0.02, 0.0, 0.0, 1, 0, 2, 2);
1801  addsol (-45.099, -0.02, 0.0, -0.0010, 1, 0, 2, 0);
1802  addsol (-0.179, -9.52, 0.0, -0.0833, 1, 0, 2, -2);
1803  addsol (-0.301, -0.33, 0.0, 0.0014, 1, 0, 2, -4);
1804  addsol (-6.382, -3.37, 0.0, -0.0481, 1, 0, -2, 2);
1805  addsol (39.528, 85.13, 0.0, -0.7136, 1, 0, -2, 0);
1806  addsol (9.366, 0.71, 0.0, -0.0112, 1, 0, -2, -2);
1807  addsol (0.202, 0.02, 0.0, 0.0, 1, 0, -2, -4);
1808  }
1809 
1810 void Moon200::solar3()
1811  {
1812  addsol (0.415, 0.1, 0.0, 0.0013, 0, 1, 2, 0);
1813  addsol (-2.152, -2.26, 0.0, -0.0066, 0, 1, 2, -2);
1814  addsol (-1.44, -1.3, 0.0, 0.0014, 0, 1, -2, 2);
1815  addsol (0.384, -0.04, 0.0, 0.0, 0, 1, -2, -2);
1816  addsol (1.938, 3.6, -0.145, 0.0401, 4, 0, 0, 0);
1817  addsol (-0.952, -1.58, 0.052, -0.0130, 4, 0, 0, -2);
1818  addsol (-0.551, -0.94, 0.032, -0.0097, 3, 1, 0, 0);
1819  addsol (-0.482, -0.57, 0.005, -0.0045, 3, 1, 0, -2);
1820  addsol (0.681, 0.96, -0.026, 0.0115, 3, -1, 0, 0);
1821  addsol (-0.297, -0.27, 0.002, -0.0009, 2, 2, 0, -2);
1822  addsol (0.254, 0.21, -0.003, 0.0, 2, -2, 0, -2);
1823  addsol (-0.25, -0.22, 0.004, 0.0014, 1, 3, 0, -2);
1824  addsol (-3.996, 0.0, 0.0, 0.0004, 2, 0, 2, 0);
1825  addsol (0.557, -0.75, 0.0, -0.009, 2, 0, 2, -2);
1826  addsol (-0.459, -0.38, 0.0, -0.0053, 2, 0, -2, 2);
1827  addsol (-1.298, 0.74, 0.0, 0.0004, 2, 0, -2, 0);
1828  addsol (0.538, 1.14, 0.0, -0.0141, 2, 0, -2, -2);
1829  addsol (0.263, 0.02, 0.0, 0.0, 1, 1, 2, 0);
1830  addsol (0.426, 0.07, 0.0, -0.0006, 1, 1, -2, -2);
1831  addsol (-0.304, 0.03, 0.0, 0.0003, 1, -1, 2, 0);
1832  addsol (-0.372, -0.19, 0.0, -0.0027, 1, -1, -2, 2);
1833  addsol (0.418, 0.0, 0.0, 0.0, 0, 0, 4, 0);
1834  addsol (-0.330, -0.04, 0.0, 0.0, 3, 0, 2, 0);
1835  }
1836 
1837 void Moon200::addn (double coeffn, int p, int q, int r, int s,
1838  double& n, double&x, double& y)
1839  {
1840  term (p,q,r,s,x,y);
1841  n=n+coeffn*y;
1842  }
1843 
1844 void Moon200::solarn (double& n)
1845  {
1846  // perturbation N of ecliptic latitude
1847  double x, y;
1848 
1849  n = 0.0;
1850  addn (-526.069, 0, 0, 1, -2, n, x, y);
1851  addn (-3.352, 0, 0, 1, -4, n, x, y);
1852  addn (44.297, 1, 0, 1, -2, n, x, y);
1853  addn (-6.0, 1, 0, 1, -4, n, x, y);
1854  addn (20.599, -1, 0, 1, 0, n, x, y);
1855  addn (-30.598, -1, 0, 1, -2, n, x, y);
1856  addn (-24.649, -2, 0, 1, 0, n, x, y);
1857  addn (-2.0, -2, 0, 1, -2, n, x, y);
1858  addn (-22.571, 0, 1, 1, -2, n, x, y);
1859  addn (10.985, 0, -1, 1, -2, n, x, y);
1860  }
1861 
1862 void Moon200::planetary (double t)
1863  {
1864  // perturbations of the ecliptic longitude by Venus and Jupiter
1865 
1866  dlam = dlam
1867  + 0.82*sinus(0.7736 -62.5512*t)+0.31*sinus(0.0466-125.1025*t)
1868  + 0.35*sinus(0.5785-25.1042*t)+0.66*sinus(0.4591+1335.8075*t)
1869  + 0.64*sinus(0.313-91.568*t)+1.14*sinus(0.148+1331.2898*t)
1870  + 0.21*sinus(0.5918+1056.5859*t)+0.44*sinus(0.5784+1322.8595*t)
1871  + 0.24*sinus(0.2275-5.7374*t)+0.28*sinus(0.2965+2.6929*t)
1872  + 0.33*sinus(0.3132+6.3368*t);
1873  }
1874 
1875 /*-------------------- Class Eclipse --------------------------------------*/
1876 /*
1877  Calculalations for Solar and Lunar Eclipses
1878  */
1879 
1880 Eclipse::Eclipse()
1881  {
1882  // some assignments just in case
1883  rs.assign(1,0,0);
1884  rm.assign(1,0,0);
1885  eshadow.assign (1,0,0);
1886  rint.assign (1,0,0);
1887  d_umbra = 0;
1888  ep2 = 0;
1889  }
1890 
1891 int Eclipse::solar (double jd, double tdut, double& phi, double& lamda)
1892  {
1893  /* Calculates various items about a solar eclipse.
1894  INPUT:
1895  jd = Modified Julian Date (UT)
1896  tdut = TDT - UT in sec
1897 
1898  OUTPUT:
1899  phi, lamda: Geographic latitude and longitude (in radians)
1900  of center of shadow if central eclipse
1901 
1902  RETURN:
1903  0 : No eclipse
1904  1 : Partial eclipse
1905  2 : Non-central annular eclipse
1906  3 : Non-central total eclipse
1907  4 : Annular eclipse
1908  5 : Total eclipse
1909  */
1910 
1911  const double flat = 0.996633; // flatting of the Earth
1912  const double ds = 218.245445; // diameter of Sun in Earth radii
1913  const double dm = 0.544986; // diameter of Moon in Earth radii
1914  double s0, s, dlt, r2, r0;
1915  int phase;
1916  Vec3 ve;
1917 
1918  // get the apparent equatorial coordinates of the sun and the moon
1919  equ_sun_moon(jd, tdut);
1920 
1921  rs[2]/=flat; // adjust for flatting of the Earth
1922  rm[2]/=flat;
1923  rint.assign ();
1924  lamda = 0;
1925  phi = 0;
1926 
1927  // intersect shadow axis with Earth
1928  eshadow = rm - rs;
1929  eshadow = vnorm(eshadow); // direction vector of shadow
1930 
1931  s0 = - dot(rm, eshadow); // distance Moon - fundamental plane
1932  r2 = dot (rm,rm);
1933  dlt = s0*s0 + 1.0 - r2;
1934  r0 = 1.0 - dlt;
1935  if (r0 > 0) r0 = sqrt (r0);
1936  else r0 = 0; // distance center of Earth - shadow axis
1937 
1938  r2 = abs(rs - rm);
1939  d_umbra = (ds - dm) * s0 / r2 - dm;// diameter of umbra at fundamental plane
1940  d_penumbra = (ds + dm) * s0 / r2 + dm;
1941 
1942  // get phase of eclipse
1943  if (r0 < 1.0)
1944  {
1945  if (dlt > 0) dlt = sqrt(dlt);
1946  else dlt = 0;
1947  s = s0 - dlt; // distance Moon - fundamental plane
1948  d_umbra = (ds - dm) * s / r2 - dm; // diameter of umbra at surface
1949  rint = rm + s * eshadow;
1950  rint[2] *= flat; // vector to intersection
1951  ve = carpol(rint);
1952  lamda = ve[1] - lsidtim(jd,0,ep2)*0.261799387799; // geographic coordinates
1953  if (lamda > M_PI) lamda -= 2.0*M_PI;
1954  if (lamda < (-M_PI)) lamda += 2.0*M_PI;
1955  phi = sqrt(rint[0]*rint[0] + rint[1]*rint[1])*0.993305615;
1956  phi = atan2(rint[2],phi);
1957 
1958  if (d_umbra > 0) phase = 4; // central annular eclipse
1959  else phase = 5; // central total eclipse
1960  }
1961  else
1962  {
1963  if (r0 < (1.0 + 0.5 * fabs(d_umbra)))
1964  {
1965  if (d_umbra > 0) phase = 2; // non-central annular eclipse
1966  else phase = 3; // non-central total eclipse
1967  }
1968  else
1969  {
1970  if (r0 < (1.0 + 0.5*d_penumbra)) phase = 1; // partial eclipse
1971  else phase = 0; // no eclipse
1972  }
1973  }
1974 
1975  rs[2]*=flat; // restore from flatting of the Earth
1976  rm[2]*=flat;
1977 
1978  return phase;
1979  }
1980 
1981 void Eclipse::maxpos (double jd, double tdut, double& phi, double& lamda)
1982  {
1983  /* Calculates the geographic position of the place of maximum eclipse
1984  at the time jd for a non-central eclipse. (NOTE that in case of
1985  a central eclipse the maximum position will be calculated by
1986  function solar. Do not use this function in that case as the
1987  result would be incorrect! No check is being done in this routine
1988  whether an eclipse is visible at all. Use maxpos only in case of a
1989  confirmed partial or non-central eclipse.
1990 
1991  INPUT:
1992  jd = Modified Julian Date (UT)
1993  tdut = TDT - UT in sec
1994 
1995  OUTPUT:
1996  phi, lamda: Geographic latitude and longitude (in radians)
1997  of maximum eclipse at the time
1998  */
1999 
2000  const double flat = 0.996633; // flatting of the Earth
2001  double s0;
2002  Vec3 ve;
2003 
2004  // get the apparent equatorial coordinates of the sun and the moon
2005  equ_sun_moon(jd, tdut);
2006 
2007  rs[2]/=flat; // adjust for flatting of the Earth
2008  rm[2]/=flat;
2009  rint.assign ();
2010  lamda = 0;
2011  phi = 0;
2012 
2013  // intersect shadow axis with Earth
2014  eshadow = rm - rs;
2015  eshadow = vnorm(eshadow); // direction vector of shadow
2016 
2017  s0 = - dot(rm, eshadow); // distance Moon - fundamental plane
2018  rint = rm + s0 * eshadow;
2019  rint = vnorm(rint); // normalize to 1 Earth radius
2020  rint[2] *= flat; // vector to position of maximum eclipse
2021  ve = carpol(rint);
2022 
2023  lamda = ve[1] - lsidtim(jd,0,ep2)*0.261799387799; // geographic coordinates
2024  if (lamda > M_PI) lamda -= 2.0*M_PI;
2025  if (lamda < (-M_PI)) lamda += 2.0*M_PI;
2026  phi = sqrt(rint[0]*rint[0] + rint[1]*rint[1])*0.993305615;
2027  phi = atan2(rint[2],phi);
2028 
2029  rs[2]*=flat; // restore from flatting of the Earth
2030  rm[2]*=flat;
2031  }
2032 
2033 void Eclipse::penumd (double jd, double tdut, Vec3& vrm, Vec3& ves,
2034  double& dpn, double& pang)
2035  {
2036  /* Calculates various items needed for finding out the northern and
2037  southern border of the penumbra during a solar eclipse.
2038 
2039  INPUT:
2040  jd = Modified Julian Date (UT)
2041  tdut = TDT - UT in sec
2042 
2043  OUTPUT:
2044  vrm = position vector of Moon adjusted for flattening
2045  ves = unit vector pointing into the direction of the center of the shadow
2046  dpn = diameter of the penumbra at the fundamental plane
2047  pang = angle of penumbral half-cone in radians
2048  */
2049 
2050  const double flat = 0.996633; // flatting of the Earth
2051  const double ds = 218.245445; // diameter of Sun in Earth radii
2052  const double dm = 0.544986; // diameter of Moon in Earth radii
2053  double s0, r2;
2054 
2055  // get the apparent equatorial coordinates of the sun and the moon
2056  equ_sun_moon(jd, tdut);
2057  rs[2]/=flat; // adjust for flatting of the Earth
2058  rm[2]/=flat;
2059 
2060  // intersect shadow axis with Earth
2061  eshadow = rm - rs;
2062  pang = abs(eshadow);
2063  eshadow = vnorm(eshadow); // direction vector of shadow
2064  ves = eshadow;
2065  vrm = rm;
2066 
2067  s0 = - dot(rm, eshadow); // distance Moon - fundamental plane
2068 
2069  r2 = abs(rs - rm);
2070  dpn = (ds + dm) * s0 / r2 + dm;
2071 
2072  // penumbral angle
2073  pang = asin((dm + ds) / (2.0 * pang));
2074 
2075  rs[2]*=flat; // restore from flatting of the Earth
2076  rm[2]*=flat;
2077  }
2078 
2079 void Eclipse::umbra (double jd, double tdut, Vec3& vrm, Vec3& ves,
2080  double& dpn, double& pang)
2081  {
2082  /* Calculates various items needed for finding out the northern and
2083  southern border of the umbra during a solar eclipse.
2084 
2085  INPUT:
2086  jd = Modified Julian Date (UT)
2087  tdut = TDT - UT in sec
2088 
2089  OUTPUT:
2090  vrm = position vector of Moon adjusted for flattening
2091  ves = unit vector pointing into the direction of the center of the shadow
2092  dpn = diameter of the umbra at the fundamental plane
2093  pang = angle of umbral half-cone in radians
2094  */
2095 
2096  const double flat = 0.996633; // flatting of the Earth
2097  const double ds = 218.245445; // diameter of Sun in Earth radii
2098  const double dm = 0.544986; // diameter of Moon in Earth radii
2099  double s0, r2;
2100 
2101  // get the apparent equatorial coordinates of the sun and the moon
2102  equ_sun_moon(jd, tdut);
2103  rs[2]/=flat; // adjust for flatting of the Earth
2104  rm[2]/=flat;
2105 
2106  // intersect shadow axis with Earth
2107  eshadow = rm - rs;
2108  pang = abs(eshadow);
2109  eshadow = vnorm(eshadow); // direction vector of shadow
2110  ves = eshadow;
2111  vrm = rm;
2112 
2113  s0 = - dot(rm, eshadow); // distance Moon - fundamental plane
2114 
2115  r2 = abs(rs - rm);
2116  dpn = (ds - dm) * s0 / r2 - dm;
2117 
2118  // umbral angle
2119  pang = asin((ds - dm) / (2.0 * pang));
2120 
2121  rs[2]*=flat; // restore from flatting of the Earth
2122  rm[2]*=flat;
2123  }
2124 
2125 void Eclipse::equ_sun_moon(double jd, double tdut)
2126  {
2127  /* Get the equatorial coordinates of the Sun and the Moon and
2128  store them in rs and rm.
2129  Also store the time t in Julian Centuries and the correction
2130  ep2 for the Apparent Sidereal Time.
2131  jd = Modified Julian Date (UT)
2132  tdut = TDT - UT in sec
2133  */
2134  double ae = 23454.77992; // 149597870.0/6378.14 = 1AE -> Earth Radii
2135  Mat3 mx;
2136 
2137  t = julcent (jd) + tdut / 3.15576e9; // =(86400.0 * 36525.0);
2138  rs = sun.position(t);
2139  rm = moon.position(t);
2140  rs = eclequ(t,rs);
2141  rm = eclequ(t,rm);
2142 
2143  // correct Moon coordinates for center of figure
2144  mx = zrot(-2.4240684e-6); // +0.5" in longitude
2145  rm = mxvct(mx,rm);
2146  mx = yrot(-1.2120342e-6); // -0.25" in latitude
2147  rm = mxvct(mx,rm);
2148  mx = nutmat(t,ep2); // nutation
2149  rs = mxvct(mx,rs); // apparent coordinates
2150  rs = aberrat(t,rs);
2151  rs *= ae;
2152  rm = mxvct(mx,rm);
2153  }
2154 
2155 double Eclipse::duration (double jd, double tdut, double& width)
2156  {
2157  /* Get the duration of a central eclipse at the center point
2158  for MJD jd and TDT-UT tdut in seconds.
2159  Also return the width of the umbra in km.
2160  A call to solar with the respective jd must have been done
2161  prior to calling this function !!!
2162  */
2163  const double omega = 4.3755e-3; // radial velocity of Earth in rad/min.
2164 
2165  double dur, lm, pa, umbold;
2166  Vec3 rold, eold, rsold, rmold;
2167  Mat3 mx;
2168 
2169  // save old values for jd
2170  rold = rint;
2171  eold = eshadow;
2172  umbold = d_umbra;
2173  rsold = rs;
2174  rmold = rm;
2175 
2176  dur = 0.1; // 0.1 min
2177  if (solar(jd+dur/1440.0,tdut, lm, pa) > 3)
2178  {
2179  mx = zrot(dur*omega);
2180  rint = mxvct(mx,rint);
2181  rint -= rold;
2182  pa = dot (rint, eold);
2183  lm = dot (rint, rint) - pa*pa;
2184  if (lm > 0) lm = sqrt(lm);
2185  else lm = 0;
2186  if (lm > 0) dur = fabs(umbold) / lm * dur * 60.0;
2187  else dur = 0;
2188  }
2189 
2190  else dur = -1;
2191 
2192  // restore old values for jd
2193  d_umbra = umbold;
2194  eshadow = eold;
2195  eold = rold*rint; // direction perpendicular of apparent shadow movement
2196  rint = rold;
2197  rs = rsold;
2198  rm = rmold;
2199 
2200  // get width of umbra at center location
2201  rold = vnorm (eold);
2202  pa = dot(rold,eshadow);
2203  if (pa > 1.0) pa = 1.0;
2204  if (pa < -1.0) pa = -1.0;
2205  pa = fabs(sin(acos(pa)));
2206  if (pa < 0.00001) pa = 0.00001;
2207  width = d_umbra / pa * 6378.14; // allow for projection
2208  width = fabs(width);
2209 
2210  return dur;
2211  }
2212 
2213 Vec3 Eclipse::GetRSun () // get Earth - Sun vector in Earth radii
2214  {
2215  return rs;
2216  }
2217 
2218 Vec3 Eclipse::GetRMoon () // get Earth - Moon vector in Earth radii
2219  {
2220  return rm;
2221  }
2222 
2223 double Eclipse::GetEp2 () // get the ep2 value
2224  {
2225  return ep2;
2226  }
2227 
2228 int Eclipse::lunar (double jd, double tdut)
2229  {
2230  /* check whether lunar eclipse is in progress at
2231  Modified Julian Date jd (UT).
2232  tdut = TDT - UT in sec
2233 
2234  RETURN:
2235  0 : No eclipse
2236  1 : Partial Penumbral Eclipse
2237  2 : Total Penumbral Eclipse
2238  3 : Partial Umbral Eclipse
2239  4 : Total Umbral Eclipse
2240  */
2241 
2242  const double dm = 0.544986; // diameter of Moon in Earth radii
2243  const double ds = 218.245445; // diameter of Sun in Earth radii
2244 
2245  double umbra, penumbra;
2246  double r2, s0, sep;
2247  int phase;
2248  Vec3 v1, v2;
2249 
2250  // get position of Sun and Moon
2251  equ_sun_moon(jd, tdut);
2252 
2253  // get radius of umbra and penumbra
2254  r2 = abs(rs);
2255  s0 = abs (rm);
2256  umbra = 1.02*fabs((ds - 2.0) * s0 / r2 - 2.0)*0.5; // radius of umbra
2257  penumbra = 1.02*fabs((ds + 2.0) * s0 / r2 + 2.0)*0.5;//radius of penumbra
2258  /* (the factor 1.02 allows for enlargment of shadow due to
2259  Earth's atmosphere) */
2260 
2261  // get angular separation of center of shadow and Moon
2262  r2 = abs(rm);
2263  sep = dot(rs,rm)/(abs(rs)*r2);
2264  if (fabs(sep) > 1.0) sep = 1.0;
2265  sep = acos(sep); // in radians
2266  sep = fabs(tan(sep)*r2); // distance of Moon and shadow in Earth radii
2267 
2268  // Now check the kind of eclipse
2269  if (sep < (umbra - dm/2.0)) phase = 4;
2270  else if (sep < (umbra + dm/2.0)) phase = 3;
2271  else if (sep < (penumbra - dm/2.0)) phase = 2;
2272  else if (sep < (penumbra + dm/2.0)) phase = 1;
2273  else phase = 0;
2274 
2275  return phase;
2276  }
2277 
2278 
stumpff
void stumpff(double e2, double &c1, double &c2, double &c3)
Definition: astrolib.cpp:1042
Sun200::state
void state(double t, Vec3 &rs, Vec3 &vs)
Definition: astrolib.cpp:1401
Refract
double Refract(double h, double p, double t)
Definition: astrolib.cpp:891
Mat3
Definition: attlib.h:62
GeoPos
Vec3 GeoPos(double jd, double ep2, double lat, double lng, double ht)
Definition: astrolib.cpp:716
Eclipse::duration
double duration(double jd, double tdut, double &width)
Definition: astrolib.cpp:2155
parab
void parab(double gm, double t0, double t, double q, double ecc, Vec3 &r1, Vec3 &v1)
Definition: astrolib.cpp:1062
Eclipse::Eclipse
Eclipse()
Definition: astrolib.cpp:1880
equecl
Vec3 equecl(double t, Vec3 &r1)
Definition: astrolib.cpp:356
eclequ
Vec3 eclequ(double t, Vec3 &r1)
Definition: astrolib.cpp:339
Vec3
Definition: attlib.h:26
Eclipse::GetRMoon
Vec3 GetRMoon()
Definition: astrolib.cpp:2218
Moon200::position
Vec3 position(double t)
Definition: astrolib.cpp:1592
Sun200::Sun200
Sun200()
Definition: astrolib.cpp:1389
Eclipse::solar
int solar(double jd, double tdut, double &phi, double &lamda)
Definition: astrolib.cpp:1891
yrot
Mat3 yrot(double a)
Definition: attlib.cpp:510
atan20
double atan20(double y, double x)
Definition: attlib.cpp:26
ddd
double ddd(int d, int m, double s)
Definition: astrolib.cpp:34
astrolib.h
lsidtim
double lsidtim(double jd, double lambda, double ep2)
Definition: astrolib.cpp:300
HorEqu
Vec3 HorEqu(double jd, double ep2, double lat, double lng, Vec3 r)
Definition: astrolib.cpp:811
eps
double eps(double t)
Definition: astrolib.cpp:323
Eclipse::lunar
int lunar(double jd, double tdut)
Definition: astrolib.cpp:2228
DefTdUt
double DefTdUt(int yi)
Definition: astrolib.cpp:222
PoleMx
Mat3 PoleMx(double xp, double yp)
Definition: astrolib.cpp:673
ellip
void ellip(double gm, double t0, double t, double a, double ecc, double m0, Vec3 &r1, Vec3 &v1)
Definition: astrolib.cpp:980
carpol
Vec3 carpol(const Vec3 &c)
Definition: attlib.cpp:185
vnorm
Vec3 vnorm(const Vec3 &c)
Definition: attlib.cpp:171
attlib.h
Eclipse::maxpos
void maxpos(double jd, double tdut, double &phi, double &lamda)
Definition: astrolib.cpp:1981
zrot
Mat3 zrot(double a)
Definition: attlib.cpp:531
pmatequ
Mat3 pmatequ(double t1, double t2)
Definition: astrolib.cpp:406
oscelm
void oscelm(double gm, double t, Vec3 &r1, Vec3 &v1, double &t0, double &m0, double &a, double &ecc, double &ran, double &aper, double &inc)
Definition: astrolib.cpp:1218
QuickSun
Vec3 QuickSun(double t)
Definition: astrolib.cpp:1357
Eclipse::GetRSun
Vec3 GetRSun()
Definition: astrolib.cpp:2213
frac
double frac(double f)
Definition: astrolib.cpp:29
Mat3::assign
void assign(double x11, double x12, double x13, double x21, double x22, double x23, double x31, double x32, double x33)
Definition: attlib.cpp:250
caldat
void caldat(double mjd, int &day, int &month, int &year, double &hour)
Definition: astrolib.cpp:144
nutecl
Mat3 nutecl(double t, double &ep2)
Definition: astrolib.cpp:635
pmatecl
Mat3 pmatecl(double t1, double t2)
Definition: astrolib.cpp:373
abs
double abs(const Vec3 &c)
Definition: attlib.cpp:101
nutmat
Mat3 nutmat(double t, double &ep2, bool hpr)
Definition: astrolib.cpp:435
mxtrn
Mat3 mxtrn(const Mat3 &m1)
Definition: attlib.cpp:380
kepler
void kepler(double gm, double t0, double t, double m0, double a, double ecc, double ran, double aper, double inc, Vec3 &r1, Vec3 &v1)
Definition: astrolib.cpp:1132
Eclipse::GetEp2
double GetEp2()
Definition: astrolib.cpp:2223
AppPos
void AppPos(double jd, double ep2, double lat, double lng, double ht, int solsys, Vec3 r, double &azim, double &elev, double &dist)
Definition: astrolib.cpp:833
hyperb
void hyperb(double gm, double t0, double t, double a, double ecc, Vec3 &r1, Vec3 &v1)
Definition: astrolib.cpp:1011
Sun200::position
Vec3 position(double t)
Definition: astrolib.cpp:1392
dms
void dms(double dd, int &d, int &m, double &s)
Definition: astrolib.cpp:58
eccanom
double eccanom(double man, double ecc)
Definition: astrolib.cpp:913
mjd
double mjd(int day, int month, int year, double hour)
Definition: astrolib.cpp:93
Vec3::assign
void assign(double x=0, double y=0, double z=0)
Definition: attlib.cpp:47
hypanom
double hypanom(double mh, double ecc)
Definition: astrolib.cpp:949
Eclipse::penumd
void penumd(double jd, double tdut, Vec3 &vrm, Vec3 &ves, double &dpn, double &pang)
Definition: astrolib.cpp:2033
M_PI
#define M_PI
Definition: GeoDataCoordinates.h:26
xrot
Mat3 xrot(double a)
Definition: attlib.cpp:489
julcent
double julcent(double mjuld)
Definition: astrolib.cpp:134
aberrat
Vec3 aberrat(double t, Vec3 &ve)
Definition: astrolib.cpp:693
EquHor
Vec3 EquHor(double jd, double ep2, double lat, double lng, Vec3 r)
Definition: astrolib.cpp:788
dot
double dot(const Vec3 &c1, const Vec3 &c2)
Definition: attlib.cpp:109
AppRADec
void AppRADec(double jd, double ep2, double lat, double lng, double azim, double elev, double &ra, double &dec)
Definition: astrolib.cpp:865
Moon200::Moon200
Moon200()
Definition: astrolib.cpp:1589
Eclipse::umbra
void umbra(double jd, double tdut, Vec3 &vrm, Vec3 &ves, double &dpn, double &pang)
Definition: astrolib.cpp:2079
polcar
Vec3 polcar(const Vec3 &c)
Definition: attlib.cpp:208
mxvct
Vec3 mxvct(const Mat3 &m1, Vec3 &v1)
Definition: attlib.cpp:474
This file is part of the KDE documentation.
Documentation copyright © 1996-2014 The KDE developers.
Generated on Tue Oct 14 2014 22:38:49 by doxygen 1.8.7 written by Dimitri van Heesch, © 1997-2006

KDE's Doxygen guidelines are available online.

marble

Skip menu "marble"
  • Main Page
  • Namespace List
  • Namespace Members
  • Alphabetical List
  • Class List
  • Class Hierarchy
  • Class Members
  • File List
  • File Members
  • Related Pages

kdeedu API Reference

Skip menu "kdeedu API Reference"
  • Analitza
  •     lib
  • kalgebra
  • kalzium
  •   libscience
  • kanagram
  • kig
  •   lib
  • klettres
  • kstars
  • libkdeedu
  •   keduvocdocument
  • marble
  • parley
  • rocs
  •   App
  •   RocsCore
  •   VisualEditor
  •   stepcore

Search



Report problems with this website to our bug tracking system.
Contact the specific authors with questions and comments about the page contents.

KDE® and the K Desktop Environment® logo are registered trademarks of KDE e.V. | Legal