File indexing completed on 2025-01-05 03:58:42
0001 // SPDX-License-Identifier: LGPL-2.1-or-later 0002 // 0003 // SPDX-FileCopyrightText: 2012 Gerhard HOLTKAMP 0004 // 0005 0006 /*************************************************************************** 0007 * Calculate Spacecraft around other planets * 0008 * * 0009 * * 0010 * Open Source Code. License: GNU LGPL Version 2+ * 0011 * * 0012 * Author: Gerhard HOLTKAMP, 26-AUG-2012 * 0013 ***************************************************************************/ 0014 0015 /*------------ include files and definitions -----------------------------*/ 0016 #include "planetarySats.h" 0017 #include <fstream> 0018 #include <iostream> 0019 #include <cstdlib> 0020 #include <cstring> 0021 #include <cmath> 0022 #include <ctime> 0023 using namespace std; 0024 0025 #include "astrolib.h" 0026 0027 // ################ Planetary Sats Class #################### 0028 0029 PlanetarySats::PlanetarySats() 0030 { 0031 plsatinit(); 0032 } 0033 0034 PlanetarySats::~PlanetarySats() 0035 { 0036 0037 } 0038 0039 double PlanetarySats::atan23 (double y, double x) 0040 { 0041 // redefine atan2 so that it doesn't crash when both x and y are 0 0042 double result; 0043 0044 if ((x == 0) && (y == 0)) result = 0; 0045 else result = atan2 (y, x); 0046 0047 return result; 0048 } 0049 0050 void PlanetarySats::plsatinit() 0051 { 0052 // initialize planetary sat data 0053 pls_moonflg = false; 0054 pls_day = 1; 0055 pls_month = 1; 0056 pls_year = 2012; 0057 pls_hour = 0; 0058 pls_minute = 0; 0059 pls_second = 0; 0060 pls_del_auto = 1; 0061 pls_step = 60.0; 0062 pls_delta_rt = 0.0; 0063 getTime(); 0064 getMars(); 0065 0066 strcpy (pls_satelmfl, "./planetarysats.txt"); 0067 0068 } 0069 0070 void PlanetarySats::getTime () // Get System Time and Date 0071 { 0072 time_t tt; 0073 int hh, mm, ss; 0074 double jd, hr; 0075 0076 tt = time(NULL); 0077 jd = 40587.0 + tt/86400.0; // seconds since 1-JAN-1970 0078 0079 jd = jd + pls_delta_rt / 24.0; 0080 caldat(jd, hh, mm, ss, hr); 0081 pls_year = ss; 0082 pls_month = mm; 0083 pls_day = hh; 0084 0085 dms(hr, hh, mm, jd); 0086 pls_hour = hh; 0087 pls_minute = mm; 0088 pls_second = int(jd); 0089 if (pls_del_auto) pls_del_tdut = DefTdUt(pls_year); 0090 setMJD(pls_year, pls_month, pls_day, hh, mm, jd); 0091 }; 0092 0093 void PlanetarySats::setStepWidth(double s) 0094 { 0095 // set the step width (in seconds) used for calculations 0096 if (s < 0.01) pls_step = 0.01; 0097 else pls_step = s; 0098 } 0099 0100 void PlanetarySats::setDeltaRT(double drt) 0101 { 0102 pls_delta_rt = drt; 0103 } 0104 0105 void PlanetarySats::setDeltaTAI_UTC(double d) 0106 { 0107 // c is the difference between TAI and UTC according to the IERS 0108 // we have to add 32.184 sec to get to the difference TT - UT 0109 // which is used in the calculations here 0110 0111 pls_del_tdut = d + 32.184; 0112 pls_del_auto = 0; 0113 } 0114 0115 void PlanetarySats::setAutoTAI_UTC() 0116 { 0117 // set the difference between TAI and UTC according to the IERS 0118 pls_del_auto = true; 0119 pls_del_tdut = DefTdUt(pls_year); 0120 } 0121 0122 void PlanetarySats::setMJD(int year, int month, int day, int hour, int min, double sec) 0123 { 0124 // set the (MJD-) time currently used for calculations to year, month, day, hour, min, sec 0125 double jd; 0126 0127 pls_year = year; 0128 pls_month = month; 0129 pls_day = day; 0130 pls_hour = hour; 0131 pls_minute = min; 0132 pls_second = sec; 0133 0134 jd = ddd(hour, min, sec); 0135 jd = mjd(day, month, year, jd); 0136 0137 pls_time = jd; 0138 0139 if (pls_del_auto) pls_del_tdut = DefTdUt(pls_year); 0140 0141 } 0142 0143 void PlanetarySats::getDatefromMJD(double mjd, int &year, int &month, int &day, int &hour, int &min, double &sec) 0144 { 0145 // convert times given in Modified Julian Date (MJD) into conventional date and time 0146 0147 double magn; 0148 0149 caldat((mjd), day, month, year, magn); 0150 dms (magn,hour,min,sec); 0151 if (sec>30.0) min++; 0152 if (min>59) 0153 { 0154 hour++; 0155 min=0; 0156 }; 0157 } 0158 0159 void PlanetarySats::setSatFile(char* fname) 0160 { 0161 strcpy (pls_satelmfl, fname); 0162 0163 } 0164 0165 void PlanetarySats::setStateVector(double mjd, double x, double y, double z, double vx, double vy, double vz) 0166 { 0167 pls_rep[0] = x; 0168 pls_rep[1] = y; 0169 pls_rep[2] = z; 0170 pls_vep[0] = vx; 0171 pls_vep[1] = vy; 0172 pls_vep[2] = vz; 0173 0174 int year = 0; 0175 int month = 0; 0176 int day = 0; 0177 int hour = 0; 0178 int min = 0; 0179 double sec = 0; 0180 getDatefromMJD(mjd, year, month, day, hour, min, sec); 0181 setMJD(year, month, day, hour, min, sec); 0182 pls_tepoch = pls_time; 0183 //pls_tepoch = pls_time + pls_del_tdut / 86400.0; // epoch in TT 0184 } 0185 0186 int PlanetarySats::getStateVector(int nsat) 0187 { 0188 // read the state vector from the planetary sat file 0189 // nsat = number of eligible sat to select (1 if first or only sat) 0190 // RETURN number of eligible sat selected, 0 if no suitable sats of file problems 0191 0192 int fsc, j, k, nst, utc; 0193 int yr, month, day, hour, min; 0194 double sec; 0195 bool searching; 0196 ifstream cfle; 0197 char satname[40]; // name of satellite 0198 char plntname[40]; // name of planet 0199 0200 strcpy(satname, ""); 0201 strcpy(plntname, ""); 0202 nst = 0; 0203 0204 searching = true; 0205 fsc = 0; 0206 0207 cfle.open(pls_satelmfl, ios::in); 0208 if (!cfle) 0209 { 0210 searching = false; 0211 cfle.close(); 0212 }; 0213 if(searching) 0214 { 0215 while (searching) 0216 { 0217 fsc = 1; 0218 0219 if (!cfle.getline(satname,40)) fsc = 0; 0220 else 0221 { 0222 k = strlen(satname); 0223 if ((k>1) && (satname[0] == '#')) 0224 { 0225 for (j=1; j<k; ++j) 0226 { 0227 pls_satname[j-1] = satname[j]; 0228 if(pls_satname[j-1] == '\n') pls_satname[j-1] = '\0'; 0229 }; 0230 pls_satname[k-1] = '\0'; 0231 } 0232 else fsc = 0; 0233 }; 0234 0235 if (cfle.eof()) 0236 { 0237 fsc = 0; 0238 searching = false; 0239 }; 0240 0241 if (fsc) 0242 { 0243 if (!cfle.getline(plntname, 40)) fsc = 0; 0244 else 0245 { 0246 k = strlen(plntname); 0247 if (k>0) 0248 { 0249 if(plntname[k-1] == '\n') plntname[k-1] = '\0'; 0250 if((k>1) && (plntname[k-2] == '\r')) plntname[k-2] = '\0'; 0251 }; 0252 }; 0253 }; 0254 0255 if (fsc) 0256 { 0257 cfle >> yr >> month >> day >> hour >> min >> sec >> utc; 0258 if (cfle.bad()) fsc = 0; 0259 if (cfle.eof()) 0260 { 0261 fsc = 0; 0262 searching = false; 0263 }; 0264 }; 0265 0266 if (fsc) 0267 { 0268 cfle >> pls_rep[0] >> pls_rep[1] >> pls_rep[2]; 0269 if (cfle.bad()) fsc = 0; 0270 if (cfle.eof()) 0271 { 0272 fsc = 0; 0273 searching = false; 0274 }; 0275 }; 0276 0277 if (fsc) 0278 { 0279 cfle >> pls_vep[0] >> pls_vep[1] >> pls_vep[2]; 0280 if (cfle.bad()) fsc = 0; 0281 }; 0282 0283 if (fsc) 0284 { 0285 if (strncmp(pls_plntname, plntname, 4) == 0) nst++; 0286 if (nst == nsat) 0287 { 0288 searching = false; 0289 0290 setMJD(yr, month, day, hour, min, sec); 0291 pls_tepoch = pls_time; 0292 if (utc) pls_tepoch = pls_tepoch + pls_del_tdut/86400.0; // epoch in TT 0293 }; 0294 }; 0295 }; 0296 cfle.close(); 0297 }; 0298 0299 if (fsc == 0) nst = 0; 0300 0301 return nst; 0302 } 0303 0304 void PlanetarySats::setPlanet(char* pname) 0305 { 0306 pls_moonflg = false; 0307 strcpy(pls_plntname, pname); 0308 if (strncmp("Mars", pname, 4) == 0) getMars(); 0309 if (strncmp("Venus", pname, 4) == 0) getVenus(); 0310 if (strncmp("Mercury", pname, 4) == 0) getMercury(); 0311 if (strncmp("Moon", pname, 4) == 0) getMoon(); 0312 } 0313 0314 void PlanetarySats::stateToKepler() 0315 { 0316 // convert state vector (mean equatorial J2000.0) into planetary Kepler elements 0317 0318 double t, dt, ag, gm, re, j2; 0319 double n, c, w, a, ecc, inc; 0320 Vec3 r1, v1; 0321 Mat3 mx; 0322 0323 dt = (pls_tepoch - 51544.5) / 36525.0; 0324 gm = pls_GM * 7.4649600000; // convert from m^3/s^2 into km^3/d^2 0325 re = pls_R0; 0326 j2 = pls_J2; 0327 0328 // convert into planet equatorial reference frame 0329 if (pls_moonflg) 0330 { 0331 mx = mxidn(); 0332 r1 = mxvct (mx, pls_rep); 0333 v1 = mxvct (mx, pls_vep); 0334 0335 } 0336 else 0337 { 0338 ag = (pls_axl0 + pls_axl1 * dt) * M_PI / 180.0; 0339 mx = zrot(ag + M_PI / 2.0); 0340 r1 = mxvct (mx, pls_rep); 0341 v1 = mxvct (mx, pls_vep); 0342 0343 ag = (pls_axb0 + pls_axb1 * dt) * M_PI / 180.0; 0344 mx = xrot(M_PI / 2.0 - ag); 0345 r1 = mxvct (mx, r1); 0346 v1 = mxvct (mx, v1); 0347 }; 0348 0349 v1 *= 86400.0; // convert into km / day 0350 0351 oscelm(gm, pls_tepoch, r1, v1, t, pls_m0, pls_a, pls_ecc, pls_ra, pls_per, pls_inc); 0352 0353 // now the mean motion 0354 a = pls_a; 0355 ecc = pls_ecc; 0356 inc = pls_inc; 0357 0358 // preliminary n 0359 if (a == 0) a = 1.0; // just in case 0360 if (a < 0) a = -a; // just in case 0361 n = sqrt (gm / (a*a*a)); 0362 0363 // correct for J2 term 0364 w = 1.0 - ecc*ecc; 0365 if (w > 0) 0366 { 0367 w = pow (w, -1.5); 0368 c = sin (inc*M_PI/180.0); 0369 n = n*(1.0 + 1.5*j2*re*re/(a*a) * w * (1.0 - 1.5*c*c)); 0370 } 0371 else n = 1.0; // do something to avoid a domain error 0372 0373 n = n / (2.0*M_PI); 0374 if (n > 1000.0) n = 1000.0; // avoid possible errors 0375 0376 pls_n0 = n; 0377 0378 } 0379 0380 void PlanetarySats::getKeplerElements(double &perc, double &apoc, double &inc, double &ecc, double &ra, double &tano, double &m0, double &a, double &n0) 0381 { 0382 // get Kepler elements of orbit with regard to the planetary equator and prime meridian. 0383 0384 double t, gm; 0385 Vec3 r1, v1; 0386 Mat3 mx; 0387 0388 0389 if (pls_moonflg) // for the Moon we normally work in J2000. Now get it into planetary 0390 { 0391 gm = pls_GM * 7.4649600000; // convert from m^3/s^2 into km^3/d^2 0392 0393 mx = getSelenographic(pls_tepoch); 0394 r1 = mxvct (mx, pls_rep); 0395 v1 = mxvct (mx, pls_vep); 0396 v1 *= 86400.0; // convert into km / day 0397 0398 oscelm(gm, pls_tepoch, r1, v1, t, m0, a, ecc, ra, tano, inc); 0399 0400 // now the mean motion 0401 if (a == 0) a = 1.0; // just in case 0402 if (a < 0) a = -a; // just in case 0403 n0 = sqrt (gm / (a*a*a)); 0404 n0 = n0 / (2.0*M_PI); 0405 } 0406 else 0407 { 0408 a = pls_a; 0409 n0 = pls_n0; 0410 m0 = pls_m0; 0411 tano = pls_per; 0412 ra = pls_ra; 0413 ecc = pls_ecc; 0414 inc = pls_inc; 0415 }; 0416 0417 perc = pls_a * (1.0 - pls_ecc) - pls_R0; 0418 apoc = pls_a * (1.0 + pls_ecc) - pls_R0; 0419 0420 } 0421 0422 int PlanetarySats::selectSat(char* sname) 0423 { 0424 // select specified satellite 0425 // RETURN 1 if successful, 0 if no suitable satellite found 0426 0427 int nst, res, sl; 0428 bool searching; 0429 0430 searching = true; 0431 nst = 1; 0432 sl = strlen(sname); 0433 0434 while (searching) 0435 { 0436 res = getStateVector(nst); 0437 if (res) 0438 { 0439 if (strncmp(pls_satname, sname, sl) == 0) searching = false; 0440 } 0441 else searching = false; 0442 nst++; 0443 }; 0444 0445 return res; 0446 } 0447 0448 void PlanetarySats::getSatName(char* sname) const 0449 { 0450 strcpy (sname, pls_satname); 0451 } 0452 0453 void PlanetarySats::currentPos() 0454 { 0455 getSatPos(pls_time); 0456 } 0457 0458 void PlanetarySats::nextStep() 0459 { 0460 pls_time = pls_time + pls_step / 86400.0; 0461 getSatPos(pls_time); 0462 } 0463 0464 double PlanetarySats::getLastMJD() const 0465 { 0466 return pls_time; 0467 } 0468 0469 void PlanetarySats::getPlanetographic(double &lng, double &lat, double &height) const 0470 { 0471 // planetographic coordinates from current state vector 0472 0473 lng = pls_lng; 0474 lat = pls_lat; 0475 height = pls_height; 0476 0477 } 0478 0479 void PlanetarySats::getFixedFrame(double &x, double &y, double &z, double &vx, double &vy, double &vz) 0480 { 0481 // last state vector coordinates in planetary fixed frame 0482 0483 x = pls_r[0]; 0484 y = pls_r[1]; 0485 z = pls_r[2]; 0486 vx = pls_v[0]; 0487 vy = pls_v[1]; 0488 vz = pls_v[2]; 0489 } 0490 0491 void PlanetarySats::getSatPos (double tutc) 0492 { 0493 // Get Position of Satellite at MJD-time t (UTC) 0494 0495 const double mp2 = 2.0*M_PI; 0496 0497 double t, dt, m0, ran, aper, inc, a, ecc, n0, re; 0498 double f, e, c, s, k, sh, j2, gm, fac, b1, b2, b3; 0499 int j; 0500 0501 Vec3 r1, v1, rg1, s2; 0502 Mat3 m1, m2; 0503 0504 // prepare orbit calculation 0505 0506 t = tutc + pls_del_tdut/86400.0; 0507 dt = t - pls_tepoch; 0508 0509 ecc = pls_ecc; 0510 if (ecc >= 1.0) ecc = 0.999; // to avoid crashes 0511 a = pls_a; 0512 n0 = mp2 * pls_n0; 0513 if (a < 1.0) a = 1.0; // avoid possible crashes later on 0514 re = pls_R0; 0515 f = pls_flat; 0516 j2 = pls_J2; 0517 gm = pls_GM * 7.4649600000; // convert from m^3/s^2 into km^3/d^2 0518 0519 // get current orbit elements ready 0520 aper = (re / a) / (1.0 - ecc*ecc); 0521 aper = 1.5 * j2 * aper*aper * n0; 0522 m0 = pls_inc*M_PI/180.0; 0523 ran = -aper * cos(m0) * dt; 0524 m0 = sin (m0); 0525 aper = aper * (2.0 - 2.5*m0*m0) * dt; 0526 ran = pls_ra * M_PI/180.0 + ran; 0527 aper = pls_per * M_PI/180.0 + aper; 0528 m0 = pls_m0*M_PI/180.0 + n0*dt; 0529 inc = pls_inc * M_PI/180.0; 0530 0531 // solve Kepler equation 0532 if (a < 1.0) a = 1.0; // avoid possible crashes later on 0533 e = eccanom(m0, ecc); 0534 fac = sqrt(1.0 - ecc*ecc); 0535 c = cos(e); 0536 s = sin(e); 0537 r1.assign (a*(c - ecc), a*fac*s, 0.0); 0538 0539 m0 = 1.0 - ecc*c; 0540 k = sqrt(gm/a); 0541 v1.assign (-k*s/m0, k*fac*c/m0, 0.0); 0542 0543 // convert into reference plane 0544 m1 = zrot (-aper); 0545 m2 = xrot (-inc); 0546 m1 *= m2; 0547 m2 = zrot (-ran); 0548 m2 = m2 * m1; 0549 r1 = mxvct (m2, r1); 0550 v1 = mxvct (m2, v1); 0551 v1 /= 86400.0; 0552 0553 // save state vector in planet-fixed frame 0554 0555 if (pls_moonflg) m1 = getSelenographic(t); 0556 else m1 = zrot((pls_W + pls_Wd*(t-51544.5))*(M_PI/180.0)); 0557 pls_r = mxvct(m1,r1); 0558 pls_v = mxvct(m1,v1); 0559 pls_r *= 1000.0; 0560 pls_v *= 1000.0; 0561 0562 // get the groundtrack coordinates 0563 0564 rg1 = mxvct(m1,r1); 0565 0566 s2 = carpol(rg1); 0567 pls_lat = s2[2]; // just preliminary 0568 pls_lng = s2[1]; // measured with the motion of rotation! 0569 if (pls_lng > mp2) pls_lng -= mp2; 0570 if (pls_lng < -M_PI) pls_lng += mp2; 0571 if (pls_lng > M_PI) pls_lng -= mp2; 0572 0573 // get height above ellipsoid and geodetic latitude 0574 if (abs(r1) > 0.1) 0575 { 0576 if (f == 0) pls_height = abs(r1) - re; 0577 else 0578 { 0579 c = f*(2.0 - f); 0580 s = r1[0]*r1[0] + r1[1]*r1[1]; 0581 sh = c*r1[2]; 0582 0583 for (j=0; j<4; ++j) 0584 { 0585 b1 = r1[2] + sh; 0586 b3 = sqrt(s+b1*b1); 0587 if (b3 < 1e-5) b3 = sin(pls_lat); // just in case 0588 else b3 = b1 / b3; 0589 b2 = re / sqrt(1.0 - c*b3*b3); 0590 sh = b2 * c * b3; 0591 }; 0592 0593 sh = r1[2] + sh; 0594 pls_lat = atan20(sh,sqrt(s)); 0595 sh = sqrt(s+sh*sh) - b2; 0596 pls_height = sh; 0597 } 0598 } 0599 else pls_height = 0; // this should never happen 0600 0601 pls_lat = pls_lat * 180.0 / M_PI; 0602 pls_lng = pls_lng * 180.0 / M_PI; 0603 0604 } 0605 0606 0607 void PlanetarySats::getMars() // Mars planetary constants 0608 { 0609 pls_J2 = 1.964e-3; 0610 pls_R0 = 3397.2; 0611 pls_flat = 0.00647630; 0612 pls_axl0 = 317.681; 0613 pls_axl1 = -0.108; 0614 pls_axb0 = 52.886; 0615 pls_axb1 = -0.061; 0616 pls_W = 176.868; 0617 pls_Wd = 350.8919830; 0618 pls_GM = 4.282828596416e+13; // 4.282837405582e+13 0619 } 0620 0621 void PlanetarySats::getVenus() // Venus planetary constants 0622 { 0623 pls_J2 = 0.027e-3; 0624 pls_R0 = 6051.9; 0625 pls_flat = 0.0; 0626 pls_axl0 = 272.72; 0627 pls_axl1 = 0.0; 0628 pls_axb0 = 67.15; 0629 pls_axb1 = 0.0; 0630 pls_W = 160.26; 0631 pls_Wd = -1.4813596; 0632 pls_GM = 3.24858761e+14; 0633 } 0634 0635 void PlanetarySats::getMercury() // Mercury planetary constants 0636 { 0637 pls_J2 = 0.0; 0638 pls_R0 = 2439.7; 0639 pls_flat = 0.0; 0640 pls_axl0 = 281.01; 0641 pls_axl1 = -0.033; 0642 pls_axb0 = 61.45; 0643 pls_axb1 = -0.005; 0644 pls_W = 329.71; 0645 pls_Wd = 6.1385025; 0646 pls_GM = 2.20320802e+13; 0647 } 0648 0649 void PlanetarySats::getMoon() // Moon planetary constants 0650 { 0651 pls_moonflg = true; 0652 pls_J2 = 0.2027e-3; 0653 pls_R0 = 1738.0; 0654 pls_flat = 0.0; 0655 pls_axl0 = 0.0; 0656 pls_axl1 = 0.0; 0657 pls_axb0 = 90.0; 0658 pls_axb1 = 0.0; 0659 pls_W = 0.0; 0660 pls_Wd = 13.17635898; 0661 pls_GM = 4.90279412e+12; 0662 } 0663 0664 Mat3 PlanetarySats::getSelenographic (double jd) 0665 { 0666 // Calculate the Matrix to transform from Mean of J2000 into selenographic 0667 // coordinates at MJD time jd. 0668 0669 double t, gam, gmp, l, omg, mln; 0670 double a, b, c, ic, gn, gp, omp; 0671 double const degrad = M_PI / 180.0; 0672 Mat3 m1, m2; 0673 0674 t = (jd - 15019.5) / 36525.0; 0675 gam = 281.2208333 + ((0.33333333e-5*t + 0.45277778e-3)*t + 1.7191750)*t; 0676 gmp = 334.3295556 + ((-0.125e-4*t - 0.10325e-1)*t + 4069.0340333)*t; 0677 l = 279.6966778 + (0.3025e-3*t + 36000.768925)*t; 0678 omg = 259.1832750 + ((0.22222222e-5*t + 0.20777778e-2)*t - 1934.1420083)*t; 0679 b = 23.45229444 + ((0.50277778e-6*t -0.16388889e-5)*t - 0.130125e-1)*t; 0680 mln = 270.4341639 + ((0.1888889e-5*t -0.11333e-2)*t + 481267.8831417)*t; 0681 ic = 1.535*degrad; 0682 gn = (l - gam)*degrad; 0683 gp = (mln - gmp)*degrad; 0684 omp = (gmp - omg)*degrad; 0685 a = -107.0*cos(gp) + 37.0*cos(gp+2.0*omp) -11.0*cos(2.0*(gp+omp)); 0686 a = a*0.000277778*degrad + ic; 0687 c = (-109.0*sin(gp) + 37.0*sin(gp+2.0*omp) - 11.0*sin(2.0*(gp+omp)))/sin(ic); 0688 c = (c*0.000277778 + omg)*degrad; 0689 gn = -12.0*sin(gp) + 59.0*sin(gn) + 18.0*sin(2.0*omp); 0690 gn = gn*0.000277778*degrad; // tau 0691 0692 b *= degrad; 0693 gam = cos(a)*cos(b) + sin(a)*sin(b)*cos(c); 0694 gmp = gam*gam; 0695 if(gmp > 1.0) gmp = 0; 0696 else gmp = sqrt(1.0 - gmp); 0697 gam = atan23(gmp,gam); // theta 0698 t = cos(a)*sin(b) - sin(a)*sin(b)*cos(c); 0699 l = -sin(a)*sin(c); 0700 0701 gmp = atan23(l,t); // phi 0702 t = sin(a)*cos(b) - cos(a)*sin(b)*cos(c); 0703 l = -sin(b)*sin(c); 0704 a = atan23(l,t); // delta 0705 c = a + mln*degrad + gn - c; // psi 0706 0707 // libration rotation matrix from Mean equator to true selenographic 0708 m1 = zrot(gmp); 0709 m2 = xrot(gam); 0710 m1 = m2 * m1; 0711 m2 = zrot(c); 0712 m1 = m2 * m1; 0713 0714 t = julcent(jd); 0715 m2 = pmatequ(0,t); // convert from mean of J2000 to mean of epoch 0716 m1 = m1 * m2; 0717 0718 return m1; 0719 } 0720