/* compute sky circumstances of an object in heliocentric hyperbolic orbit. */ static int obj_parabolic (Now *np, Obj *op) { double lsn, rsn; /* true geoc lng of sun; dist from sn to earth*/ double lam; /* geocentric ecliptic longitude */ double bet; /* geocentric ecliptic latitude */ double mag; /* magnitude */ double inc, om, Om; double lpd, psi, rp, rho; double dt; int pass; /* find solar ecliptical longitude and distance to sun from earth */ sunpos (mjed, &lsn, &rsn, 0); /* two passes to correct lam and bet for light travel time. */ dt = 0.0; for (pass = 0; pass < 2; pass++) { reduce_elements (op->p_epoch, mjd-dt, degrad(op->p_inc), degrad(op->p_om), degrad(op->p_Om), &inc, &om, &Om); comet (mjed-dt, op->p_ep, inc, om, op->p_qp, Om, &lpd, &psi, &rp, &rho, &lam, &bet); dt = rho*LTAU/3600.0/24.0; /* light travel time, in days / AU */ } /* fill in all of op->s_* stuff except s_size and s_mag */ cir_sky (np, lpd, psi, rp, &rho, lam, bet, lsn, rsn, op); /* compute magnitude and size */ gk_mag (op->p_g, op->p_k, rp, rho, &mag); set_smag (op, mag); op->s_size = (float)(op->p_size / rho); return (0); }
/* low precision ecliptic coordinates of Pluto from mean orbit. * Only for sake of completeness outside available perturbation theories. */ static void pluto_ell ( double mj, /* epoch */ double *ret) /* ecliptic coordinates {l,b,r} at equinox of date */ { /* mean orbital elements of Pluto. * The origin of these is somewhat obscure. */ double a = 39.543, /* semimajor axis, au */ e = 0.2490, /* excentricity */ inc0 = 17.140, /* inclination, deg */ Om0 = 110.307, /* long asc node, deg */ omeg0 = 113.768, /* arg of perihel, deg */ mjp = 2448045.539 - MJD0, /* epoch of perihel */ mjeq = J2000, /* equinox of elements */ n = 144.9600/36525.; /* daily motion, deg */ double inc, Om, omeg; /* orbital elements at epoch of date */ double ma, ea, nu; /* mean, excentric and true anomaly */ double lo, slo, clo; /* longitude in orbit from asc node */ reduce_elements(mjeq, mj, degrad(inc0), degrad(omeg0), degrad(Om0), &inc, &omeg, &Om); ma = degrad((mj - mjp) * n); anomaly(ma, e, &nu, &ea); ret[2] = a * (1.0 - e*cos(ea)); /* r */ lo = omeg + nu; slo = sin(lo); clo = cos(lo); ret[1] = asin(slo * sin(inc)); /* b */ ret[0] = atan2(slo * cos(inc), clo) + Om; /* l */ }
/* given a Julian date, return the lunar libration in lat and long, in rads. */ void llibration (double JD, double *llatp, double *llonp) { double lg, lt; /* arc seconds */ lg = gplan (JD, &liblon); lt = gplan (JD, &liblat); *llonp = degrad (lg/3600.0); *llatp = degrad (lt/3600.0); }
/* compute position of secondary component of a BINARYSTAR */ static int obj_2binary (Now *np, Obj *op) { if (op->b_nbp > 0) { /* we just have discrete pa/sep, project each from primary */ int i; for (i = 0; i < op->b_nbp; i++) { BinPos *bp = &op->b_bp[i]; bp->bp_dec = op->s_dec + bp->bp_sep*cos(bp->bp_pa); bp->bp_ra = op->s_ra + bp->bp_sep*sin(bp->bp_pa)/cos(op->s_dec); } } else { BinOrbit *bp = &op->b_bo; double t, theta, rho; mjd_year (mjd, &t); binaryStarOrbit (t, bp->bo_T, bp->bo_e, bp->bo_o, bp->bo_O, bp->bo_i, bp->bo_a, bp->bo_P, &theta, &rho); bp->bo_pa = (float)theta; bp->bo_sep = (float)rho; rho = degrad(rho/3600.); /* arc secs to rads */ bp->bo_dec = op->s_dec + rho*cos(theta); bp->bo_ra = op->s_ra + rho*sin(theta)/cos(op->s_dec); } return (0); }
/* set svis according to whether moon is in sun light */ static void moonSVis( Obj *sop, /* SUN */ Obj *jop, /* jupiter */ MoonData md[J_NMOONS]) { double esd = sop->s_edist; double eod = jop->s_edist; double sod = jop->s_sdist; double soa = degrad(jop->s_elong); double esa = asin(esd*sin(soa)/sod); double h = sod*jop->s_hlat; double nod = h*(1./eod - 1./sod); double sca = cos(esa), ssa = sin(esa); int i; for (i = 1; i < J_NMOONS; i++) { MoonData *mdp = &md[i]; double xp = sca*mdp->x + ssa*mdp->z; double yp = mdp->y; double zp = -ssa*mdp->x + sca*mdp->z; double ca = cos(nod), sa = sin(nod); double xpp = xp; double ypp = ca*yp - sa*zp; double zpp = sa*yp + ca*zp; int outside = xpp*xpp + ypp*ypp > 1.0; int infront = zpp > 0.0; mdp->svis = outside || infront; } }
/* from W. M. Smart */ static void binaryStarOrbit ( double t, /* desired ephemeris epoch, year */ double T, /* epoch of periastron, year */ double e, /* eccentricity */ double o, /* argument of periastron, degrees */ double O, /* ascending node, degrees */ double i, /* inclination, degrees */ double a, /* semi major axis, arcsecs */ double P, /* period, years */ double *thetap, /* position angle, rads E of N */ double *rhop) /* separation, arcsecs */ { double M, E, cosE, nu, cosnu, r, rho, theta; /* find mean anomaly, insure 0..2*PI */ M = 2*PI/P*(t-T); range (&M, 2*PI); /* solve for eccentric anomaly */ E = solveKepler (M, e); cosE = cos(E); /* find true anomaly and separation */ cosnu = (cosE - e)/(1.0 - e*cosE); r = a*(1.0 - e*e)/(1.0 + e*cosnu); nu = acos(cosnu); if (E > PI) nu = -nu; /* project onto sky */ theta = atan(tan(nu+degrad(o))*cos(degrad(i))) + degrad(O); rho = r*cos(nu+degrad(o))/cos(theta-degrad(O)); if (rho < 0) { theta += PI; rho = -rho; } range (&theta, 2*PI); *thetap = theta; *rhop = rho; }
static void transform_point( double dy, double rot, double pt[ 3 ] ) { double r, theta; pt[ 1 ] += dy; if ( pt[ 0 ] == 0.0 && pt[ 2 ] == 0.0 ) return; r = sqrt( pt[ 0 ] * pt[ 0 ] + pt[ 2 ] * pt[ 2 ] ); theta = PI / 2.0 - atan2( pt[ 2 ], pt[ 0 ] ) - degrad( rot ); pt[ 0 ] = r * sin( theta ); pt[ 2 ] = r * cos( theta ); }
static void unrefractLT15 (double pr, double tr, double aa, double *ta) { double aadeg = raddeg(aa); double r, a, b; a = ((2e-5*aadeg+1.96e-2)*aadeg+1.594e-1)*pr; b = (273+tr)*((8.45e-2*aadeg+5.05e-1)*aadeg+1); r = degrad(a/b); *ta = (aa < 0 && r < 0) ? aa : aa - r; /* 0 below ~5 degs */ }
/* given the modified Julian date, mj, find the mean obliquity of the * ecliptic, *eps, in radians. * * IAU expression (see e.g. Astron. Almanac 1984); stern */ void obliquity (double mj, double *eps) { static double lastmj = -16347, lasteps; if (mj != lastmj) { double t = (mj - J2000)/36525.; /* centuries from J2000 */ lasteps = degrad(23.4392911 + /* 23^ 26' 21".448 */ t * (-46.8150 + t * ( -0.00059 + t * ( 0.001813 )))/3600.0); lastmj = mj; } *eps = lasteps; }
/* find moon's circumstances now. */ static int moon_cir (Now *np, Obj *op) { double lsn, rsn; /* true geoc lng of sun; dist from sn to earth*/ double lam; /* geocentric ecliptic longitude */ double bet; /* geocentric ecliptic latitude */ double edistau; /* earth-moon dist, in au */ double el; /* elongation, rads east */ double ms; /* sun's mean anomaly */ double md; /* moon's mean anomaly */ double i; moon (mjed, &lam, &bet, &edistau, &ms, &md); /* mean ecliptic & EOD*/ sunpos (mjed, &lsn, &rsn, NULL); /* mean ecliptic & EOD*/ op->s_hlong = (float)lam; /* save geo in helio fields */ op->s_hlat = (float)bet; /* find angular separation from sun */ elongation (lam, bet, lsn, &el); op->s_elong = (float)raddeg(el); /* want degrees */ /* solve triangle of earth, sun, and elongation for moon-sun dist */ op->s_sdist = (float) sqrt (edistau*edistau + rsn*rsn - 2.0*edistau*rsn*cos(el)); /* TODO: improve mag; this is based on a flat moon model. */ i = -12.7 + 2.5*(log10(PI) - log10(PI/2*(1+1.e-6-cos(el)))) + 5*log10(edistau/.0025) /* dist */; set_smag (op, i); /* find phase -- allow for projection effects */ i = 0.1468*sin(el)*(1 - 0.0549*sin(md))/(1 - 0.0167*sin(ms)); op->s_phase = (float)((1+cos(PI-el-degrad(i)))/2*100); /* fill moon's ra/dec, alt/az in op and update for topo dist */ cir_pos (np, bet, lam, &edistau, op); op->s_edist = (float)edistau; op->s_size = (float)(3600*2.0*raddeg(asin(MRAD/MAU/edistau))); /* moon angular dia, seconds */ return (0); }
/* apply relativistic light bending correction to ra/dec; stern * * The algorithm is from: * Mean and apparent place computations in the new IAU * system. III - Apparent, topocentric, and astrometric * places of planets and stars * KAPLAN, G. H.; HUGHES, J. A.; SEIDELMANN, P. K.; * SMITH, C. A.; YALLOP, B. D. * Astronomical Journal (ISSN 0004-6256), vol. 97, April 1989, p. 1197-1210. * * This article is a very good collection of formulea for geocentric and * topocentric place calculation in general. The apparent and * astrometric place calculation in this file currently does not follow * the strict algorithm from this paper and hence is not fully correct. * The entire calculation is currently based on the rotating EOD frame and * not the "inertial" J2000 frame. */ static void deflect ( double mjd1, /* equinox */ double lpd, double psi, /* heliocentric ecliptical long / lat */ double rsn, double lsn, /* distance and longitude of sun */ double rho, /* geocentric distance */ double *ra, double *dec)/* geocentric equatoreal */ { double hra, hdec; /* object heliocentric equatoreal */ double el; /* HELIOCENTRIC elongation object--earth */ double g1, g2; /* relativistic weights */ double u[3]; /* object geocentric cartesian */ double q[3]; /* object heliocentric cartesian unit vect */ double e[3]; /* earth heliocentric cartesian unit vect */ double qe, uq, eu; /* scalar products */ int i; /* counter */ #define G 1.32712438e20 /* heliocentric grav const; in m^3*s^-2 */ #define c 299792458.0 /* speed of light in m/s */ elongation(lpd, psi, lsn-PI, &el); el = fabs(el); /* only continue if object is within about 10 deg around the sun, * not obscured by the sun's disc (radius 0.25 deg) and farther away * than the sun. * * precise geocentric deflection is: g1 * tan(el/2) * radially outwards from sun; the vector munching below * just applys this component-wise * Note: el = HELIOCENTRIC elongation. * g1 is always about 0.004 arc seconds * g2 varies from 0 (highest contribution) to 2 */ if (el<degrad(170) || el>degrad(179.75) || rho<rsn) return; /* get cartesian vectors */ sphcart(*ra, *dec, rho, u, u+1, u+2); ecl_eq(mjd1, psi, lpd, &hra, &hdec); sphcart(hra, hdec, 1.0, q, q+1, q+2); ecl_eq(mjd1, 0.0, lsn-PI, &hra, &hdec); sphcart(hra, hdec, 1.0, e, e+1, e+2); /* evaluate scalar products */ qe = uq = eu = 0.0; for(i=0; i<=2; ++i) { qe += q[i]*e[i]; uq += u[i]*q[i]; eu += e[i]*u[i]; } g1 = 2*G/(c*c*MAU)/rsn; g2 = 1 + qe; /* now deflect geocentric vector */ g1 /= g2; for(i=0; i<=2; ++i) u[i] += g1*(uq*e[i] - eu*q[i]); /* back to spherical */ cartsph(u[0], u[1], u[2], ra, dec, &rho); /* rho thrown away */ }
/* compute sky circumstances of an object in heliocentric hyperbolic orbit. */ static int obj_hyperbolic (Now *np, Obj *op) { double lsn, rsn; /* true geoc lng of sun; dist from sn to earth*/ double dt; /* light travel time to object */ double lg; /* helio long of earth */ double nu; /* true anomaly and eccentric anomaly */ double rp=0; /* distance from the sun */ double lo, slo, clo; /* angle from ascending node */ double inc; /* inclination */ double psi=0; /* heliocentric latitude */ double spsi=0, cpsi=0; /* trig of heliocentric latitude */ double lpd; /* heliocentric longitude */ double rho=0; /* distance from the Earth */ double om; /* arg of perihelion */ double Om; /* long of ascending node. */ double lam; /* geocentric ecliptic longitude */ double bet; /* geocentric ecliptic latitude */ double e; /* fast eccentricity */ double ll=0, sll, cll; /* helio angle between object and earth */ double mag; /* magnitude */ double a; /* mean distance */ double tp; /* time from perihelion (days) */ double rpd=0; double y; int pass; /* find solar ecliptical longitude and distance to sun from earth */ sunpos (mjed, &lsn, &rsn, 0); lg = lsn + PI; e = op->h_e; a = op->h_qp/(e - 1.0); /* correct for light time by computing position at time mjd, then * again at mjd-dt, where * dt = time it takes light to travel earth-object distance. */ dt = 0; for (pass = 0; pass < 2; pass++) { reduce_elements (op->h_epoch, mjd-dt, degrad(op->h_inc), degrad (op->h_om), degrad (op->h_Om), &inc, &om, &Om); tp = mjed - dt - op->h_ep; if (vrc (&nu, &rp, tp, op->h_e, op->h_qp) < 0) op->o_flags |= NOCIRCUM; nu = degrad(nu); lo = nu + om; slo = sin(lo); clo = cos(lo); spsi = slo*sin(inc); y = slo*cos(inc); psi = asin(spsi); lpd = atan(y/clo)+Om; if (clo<0) lpd += PI; range (&lpd, 2*PI); cpsi = cos(psi); rpd = rp*cpsi; ll = lpd-lg; rho = sqrt(rsn*rsn+rp*rp-2*rsn*rp*cpsi*cos(ll)); dt = rho*5.775518e-3; /* light travel time, in days */ } /* compute sin and cos of ll */ sll = sin(ll); cll = cos(ll); /* find geocentric ecliptic longitude and latitude */ if (rpd < rsn) lam = atan(-1*rpd*sll/(rsn-rpd*cll))+lg+PI; else lam = atan(rsn*sll/(rpd-rsn*cll))+lpd; range (&lam, 2*PI); bet = atan(rpd*spsi*sin(lam-lpd)/(cpsi*rsn*sll)); /* fill in all of op->s_* stuff except s_size and s_mag */ cir_sky (np, lpd, psi, rp, &rho, lam, bet, lsn, rsn, op); /* compute magnitude and size */ gk_mag (op->h_g, op->h_k, rp, rho, &mag); set_smag (op, mag); op->s_size = (float)(op->h_size / rho); return (0); }
/* compute location of GRS and Galilean moons. * if md == NULL, just to cml. * from "Astronomical Formulae for Calculators", 2nd ed, by Jean Meeus, * Willmann-Bell, Richmond, Va., U.S.A. (c) 1982, chapters 35 and 36. */ void meeus_jupiter( double d, double *cmlI, double *cmlII, /* central meridian longitude, rads */ MoonData md[J_NMOONS]) /* fill in md[1..NM-1].x/y/z for each moon. * N.B. md[0].ra/dec must already be set */ { #define dsin(x) sin(degrad(x)) #define dcos(x) cos(degrad(x)) double A, B, Del, J, K, M, N, R, V; double cor_u1, cor_u2, cor_u3, cor_u4; double solc, tmp, G, H, psi, r, r1, r2, r3, r4; double u1, u2, u3, u4; double lam, Ds; double z1, z2, z3, z4; double De, dsinDe; double theta, phi; double tvc, pvc; double salpha, calpha; int i; V = 134.63 + 0.00111587 * d; M = (358.47583 + 0.98560003*d); N = (225.32833 + 0.0830853*d) + 0.33 * dsin (V); J = 221.647 + 0.9025179*d - 0.33 * dsin(V); A = 1.916*dsin(M) + 0.02*dsin(2*M); B = 5.552*dsin(N) + 0.167*dsin(2*N); K = (J+A-B); R = 1.00014 - 0.01672 * dcos(M) - 0.00014 * dcos(2*M); r = 5.20867 - 0.25192 * dcos(N) - 0.00610 * dcos(2*N); Del = sqrt (R*R + r*r - 2*R*r*dcos(K)); psi = raddeg (asin (R/Del*dsin(K))); *cmlI = degrad(268.28 + 877.8169088*(d - Del/173) + psi - B); range (cmlI, 2*PI); *cmlII = degrad(290.28 + 870.1869088*(d - Del/173) + psi - B); range (cmlII, 2*PI); /* that's it if don't want moon info too */ if (!md) return; solc = (d - Del/173.); /* speed of light correction */ tmp = psi - B; u1 = 84.5506 + 203.4058630 * solc + tmp; u2 = 41.5015 + 101.2916323 * solc + tmp; u3 = 109.9770 + 50.2345169 * solc + tmp; u4 = 176.3586 + 21.4879802 * solc + tmp; G = 187.3 + 50.310674 * solc; H = 311.1 + 21.569229 * solc; cor_u1 = 0.472 * dsin (2*(u1-u2)); cor_u2 = 1.073 * dsin (2*(u2-u3)); cor_u3 = 0.174 * dsin (G); cor_u4 = 0.845 * dsin (H); r1 = 5.9061 - 0.0244 * dcos (2*(u1-u2)); r2 = 9.3972 - 0.0889 * dcos (2*(u2-u3)); r3 = 14.9894 - 0.0227 * dcos (G); r4 = 26.3649 - 0.1944 * dcos (H); md[1].x = -r1 * dsin (u1+cor_u1); md[2].x = -r2 * dsin (u2+cor_u2); md[3].x = -r3 * dsin (u3+cor_u3); md[4].x = -r4 * dsin (u4+cor_u4); lam = 238.05 + 0.083091*d + 0.33*dsin(V) + B; Ds = 3.07*dsin(lam + 44.5); De = Ds - 2.15*dsin(psi)*dcos(lam+24.) - 1.31*(r-Del)/Del*dsin(lam-99.4); dsinDe = dsin(De); z1 = r1 * dcos(u1+cor_u1); z2 = r2 * dcos(u2+cor_u2); z3 = r3 * dcos(u3+cor_u3); z4 = r4 * dcos(u4+cor_u4); md[1].y = z1*dsinDe; md[2].y = z2*dsinDe; md[3].y = z3*dsinDe; md[4].y = z4*dsinDe; /* compute sky transformation angle as triple vector product */ tvc = PI/2.0 - md[0].dec; pvc = md[0].ra; theta = PI/2.0 - POLE_DEC; phi = POLE_RA; salpha = -sin(tvc)*sin(theta)*(cos(pvc)*sin(phi) - sin(pvc)*cos(phi)); calpha = sqrt (1.0 - salpha*salpha); for (i = 0; i < J_NMOONS; i++) { double tx = md[i].x*calpha + md[i].y*salpha; double ty = -md[i].x*salpha + md[i].y*calpha; md[i].x = tx; md[i].y = ty; } md[1].z = z1; md[2].z = z2; md[3].z = z3; md[4].z = z4; }
/****************************************************************** * adapted from BdL FORTRAN Code; stern * * Reference : Bureau des Longitudes - PBGF9502 * * Object : calculate a VSOP87 position for a given time. * * Input : * * mj modified julian date, counted from J1900.0 * time scale : dynamical time TDB. * * obj object number as in astro.h, NB: not for pluto * * prec relative precision * * if prec is equal to 0 then the precision is the precision * p0 of the complete solution VSOP87. * Mercury p0 = 0.6 10**-8 * Venus p0 = 2.5 10**-8 * Earth p0 = 2.5 10**-8 * Mars p0 = 10.0 10**-8 * Jupiter p0 = 35.0 10**-8 * Saturn p0 = 70.0 10**-8 * Uranus p0 = 8.0 10**-8 * Neptune p0 = 42.0 10**-8 * * if prec is not equal to 0, let us say in between p0 and * 10**-3, the precision is : * for the positions : * - prec*a0 au for the distances. * - prec rad for the other variables. * for the velocities : * - prec*a0 au/day for the distances. * - prec rad/day for the other variables. * a0 is the semi-major axis of the body. * * Output : * * ret[6] array of the results (double). * * for spherical coordinates : * 1: longitude (rd) * 2: latitude (rd) * 3: radius (au) * #if VSOP_GETRATE: * 4: longitude velocity (rad/day) * 5: latitude velocity (rad/day) * 6: radius velocity (au/day) * * return: error index (int) * 0: no error. * 2: object out of range [MERCURY .. NEPTUNE, SUN] * 3: precision out of range [0.0 .. 1e-3] ******************************************************************/ int vsop87 (double mj, int obj, double prec, double *ret) { static double (*vx_map[])[3] = { /* data tables */ vx_mercury, vx_venus, vx_mars, vx_jupiter, vx_saturn, vx_uranus, vx_neptune, 0, vx_earth, }; static int (*vn_map[])[3] = { /* indexes */ vn_mercury, vn_venus, vn_mars, vn_jupiter, vn_saturn, vn_uranus, vn_neptune, 0, vn_earth, }; static double a0[] = { /* semimajor axes; for precision ctrl only */ 0.39, 0.72, 1.5, 5.2, 9.6, 19.2, 30.1, 39.5, 1.0, }; double (*vx_obj)[3] = vx_map[obj]; /* VSOP87 data and indexes */ int (*vn_obj)[3] = vn_map[obj]; double t[VSOP_MAXALPHA+1]; /* powers of time */ double t_abs[VSOP_MAXALPHA+1]; /* powers of abs(time) */ double q; /* aux for precision control */ int i, cooidx, alpha; /* misc indexes */ if (obj == PLUTO || obj > SUN) return (2); if (prec < 0.0 || prec > 1e-3) return(3); /* zero result array */ for (i = 0; i < 6; ++i) ret[i] = 0.0; /* time and its powers */ t[0] = 1.0; t[1] = (mj - J2000)/VSOP_A1000; for (i = 2; i <= VSOP_MAXALPHA; ++i) t[i] = t[i-1] * t[1]; t_abs[0] = 1.0; for (i = 1; i <= VSOP_MAXALPHA; ++i) t_abs[i] = fabs(t[i]); /* precision control */ q = -log10(prec + 1e-35) - 2; /* decades below 1e-2 */ q = VSOP_ASCALE * prec / 10.0 / q; /* reduce threshold progressively * for higher precision */ /* do the term summation; first the spatial dimensions */ for (cooidx = 0; cooidx < 3; ++cooidx) { /* then the powers of time */ for (alpha = 0; vn_obj[alpha+1][cooidx] ; ++alpha) { double p, term, termdot; /* precision threshold */ p= alpha ? q/(t_abs[alpha] + alpha*t_abs[alpha-1]*1e-4 + 1e-35) : q; #if VSOP_SPHERICAL if (cooidx == 2) /* scale by semimajor axis for radius */ #endif p *= a0[obj]; term = termdot = 0.0; for (i = vn_obj[alpha][cooidx]; i < vn_obj[alpha+1][cooidx]; ++i) { double a, b, c, arg; a = vx_obj[i][0]; if (a < p) continue; /* ignore small terms */ b = vx_obj[i][1]; c = vx_obj[i][2]; arg = b + c * t[1]; term += a * cos(arg); #if VSOP_GETRATE termdot += -c * a * sin(arg); #endif } ret[cooidx] += t[alpha] * term; #if VSOP_GETRATE ret[cooidx + 3] += t[alpha] * termdot + ((alpha > 0) ? alpha * t[alpha - 1] * term : 0.0); #endif } /* alpha */ } /* cooidx */ for (i = 0; i < 6; ++i) ret[i] /= VSOP_ASCALE; #if VSOP_SPHERICAL /* reduce longitude to 0..2pi */ ret[0] -= floor(ret[0]/(2.*PI)) * (2.*PI); #endif #if VSOP_GETRATE /* convert millenium rate to day rate */ for (i = 3; i < 6; ++i) ret[i] /= VSOP_A1000; #endif #if VSOP_SPHERICAL /* reduction from dynamical equinox of VSOP87 to FK5; */ if (prec < 5e-7) { /* 5e-7 rad = 0.1 arc seconds */ double L1, c1, s1; L1 = ret[0] - degrad(13.97 * t[1] - 0.031 * t[2]); c1 = cos(L1); s1 = sin(L1); ret[0] += degrad(-0.09033 + 0.03916 * (c1 + s1) * tan(ret[1]))/3600.0; ret[1] += degrad(0.03916 * (c1 - s1))/3600.0; } #endif return (0); }
/* get jupiter info in md[0], moon info in md[1..J_NMOONS-1]. * if !dir always use meeus model. * if !jop caller just wants md[] for names * N.B. we assume sop and jop are updated. */ void jupiter_data ( double Mjd, /* mjd */ char dir[], /* dir in which to look for helper files */ Obj *sop, /* Sun */ Obj *jop, /* jupiter */ double *sizep, /* jup angular diam, rads */ double *cmlI, double *cmlII, /* central meridian longitude, rads */ double *polera, double *poledec, /* pole location */ MoonData md[J_NMOONS]) /* return info */ { double JD; /* always copy back at least for name */ memcpy (md, jmd, sizeof(jmd)); /* pole */ if (polera) *polera = POLE_RA; if (poledec) *poledec = POLE_DEC; /* nothing else if repeat call or just want names */ if (Mjd == mdmjd || !jop) { if (jop) { *sizep = sizemjd; *cmlI = cmlImjd; *cmlII = cmlIImjd; } return; } JD = Mjd + MJD0; /* planet in [0] */ md[0].ra = jop->s_ra; md[0].dec = jop->s_dec; md[0].mag = get_mag(jop); md[0].x = 0; md[0].y = 0; md[0].z = 0; md[0].evis = 1; md[0].svis = 1; /* size is straight from jop */ *sizep = degrad(jop->s_size/3600.0); /* mags from JPL ephemeris */ md[1].mag = 5.7; md[2].mag = 5.8; md[3].mag = 5.3; md[4].mag = 6.7; /* get moon data from BDL if possible, else Meeus' model. * always use Meeus for cml */ if (use_bdl (JD, dir, md) == 0) meeus_jupiter (Mjd, cmlI, cmlII, NULL); else meeus_jupiter (Mjd, cmlI, cmlII, md); /* set visibilities */ moonSVis (sop, jop, md); moonPShad (sop, jop, md); moonEVis (md); moonTrans (md); /* fill in moon ra and dec */ moonradec (*sizep, md); /* save */ mdmjd = Mjd; sizemjd = *sizep; cmlImjd = *cmlI; cmlIImjd = *cmlII; memcpy (jmd, md, sizeof(jmd)); }
int main(int argc, char *argv[]) { Now *np; Now now; ObjF *op = NULL; char msg[256]; char *wfsroot, *filename, ra_hms[50], dec_dms[50]; double ra, dec, fov, mag, dra, ddec, dist, fo_ra, fo_dec; int nop = 0; int i; now.n_mjd = 0.0; now.n_lat = degrad(31.689); now.n_lng = degrad(-110.885); now.n_tz = 7.0; now.n_temp = 20.0; now.n_pressure = 700.0; now.n_elev = 2580.0/ERAD; now.n_dip = 0.0; now.n_epoch = EOD; sprintf(now.n_tznm, "UTC-7"); np = &now; time_fromsys(np); ra = hrrad(atof(argv[1])); dec = degrad(atof(argv[2])); fov = degrad(atof(argv[3])); mag = atof(argv[4]); /* precess(mjd, J2000, &ra, &dec); */ if (getenv("WFSROOT")) { wfsroot = getenv("WFSROOT"); filename = (char *) malloc(strlen(wfsroot) + 20); strcpy(filename, wfsroot); strcat(filename, "/wfscat/tycho.xe2"); } else { filename = "/mmt/shwfs/wfscat/tycho.xe2"; } nop = xe2fetch(filename, np, ra, dec, fov, mag, &op, msg); for (i=0; i<nop; i++) { dist = raddeg( acos( sin(op[i].fo_dec)*sin(dec) + cos(op[i].fo_dec)*cos(dec)*cos(ra-op[i].fo_ra) ) ); dist *= 60.0; fo_ra = radhr(op[i].fo_ra); fo_dec = raddeg(op[i].fo_dec); sprintf(ra_hms, "%02d:%02d:%05.2f", (int) fo_ra, (int) ((fo_ra - (int)fo_ra)*60.0), ( (fo_ra - (int)fo_ra)*60.0 - (int) ((fo_ra - (int)fo_ra)*60.0) )*60.0); if (fo_dec < 0.0 && (int)fo_dec >= 0) { sprintf(dec_dms, "-%02d:%02d:%06.3f", (int)fo_dec, abs( (int) ((fo_dec - (int)fo_dec)*60.0) ), fabs( ( (fo_dec - (int)fo_dec)*60.0 - (int) ((fo_dec - (int)fo_dec)*60.0) )*60.0) ); } else { sprintf(dec_dms, "%+03d:%02d:%06.3f", (int)fo_dec, abs( (int) ((fo_dec - (int)fo_dec)*60.0) ), fabs( ( (fo_dec - (int)fo_dec)*60.0 - (int) ((fo_dec - (int)fo_dec)*60.0) )*60.0) ); } printf("%15s %4.1f mag %2s %c %s %s %+07.3f %+07.2f %9.2f\n", op[i].co_name, op[i].co_mag/MAGSCALE, op[i].fo_spect, op[i].fo_class, ra_hms, dec_dms, 0.1*op[i].fo_pma/cos(op[i].fo_dec), 0.1*op[i].fo_pmd, dist); } return(0); }
/* compute sky circumstances of an object in heliocentric elliptic orbit at *np. */ static int obj_elliptical (Now *np, Obj *op) { double lsn, rsn; /* true geoc lng of sun; dist from sn to earth*/ double dt; /* light travel time to object */ double lg; /* helio long of earth */ double nu; /* true anomaly */ double rp=0; /* distance from the sun */ double lo, slo, clo; /* angle from ascending node */ double inc; /* inclination */ double psi=0; /* heliocentric latitude */ double spsi=0, cpsi=0; /* trig of heliocentric latitude */ double lpd; /* heliocentric longitude */ double rho=0; /* distance from the Earth */ double om; /* arg of perihelion */ double Om; /* long of ascending node. */ double lam; /* geocentric ecliptic longitude */ double bet; /* geocentric ecliptic latitude */ double ll=0, sll, cll; /* helio angle between object and earth */ double mag; /* magnitude */ double e_n; /* mean daily motion */ double tp; /* time from perihelion (days) */ double rpd=0; double y; int pass; /* find location of earth from sun now */ sunpos (mjed, &lsn, &rsn, 0); lg = lsn + PI; /* mean daily motion is derived fro mean distance */ e_n = 0.9856076686/pow((double)op->e_a, 1.5); /* correct for light time by computing position at time mjd, then * again at mjd-dt, where * dt = time it takes light to travel earth-object distance. */ dt = 0; for (pass = 0; pass < 2; pass++) { reduce_elements (op->e_epoch, mjd-dt, degrad(op->e_inc), degrad (op->e_om), degrad (op->e_Om), &inc, &om, &Om); tp = mjed - dt - (op->e_cepoch - op->e_M/e_n); if (vrc (&nu, &rp, tp, op->e_e, op->e_a*(1-op->e_e)) < 0) op->o_flags |= NOCIRCUM; nu = degrad(nu); lo = nu + om; slo = sin(lo); clo = cos(lo); spsi = slo*sin(inc); y = slo*cos(inc); psi = asin(spsi); lpd = atan(y/clo)+Om; if (clo<0) lpd += PI; range (&lpd, 2*PI); cpsi = cos(psi); rpd = rp*cpsi; ll = lpd-lg; rho = sqrt(rsn*rsn+rp*rp-2*rsn*rp*cpsi*cos(ll)); printf("\nrho from sqrt(): %f\n", rho); dt = rho*LTAU/3600.0/24.0; /* light travel time, in days / AU */ } /* compute sin and cos of ll */ sll = sin(ll); cll = cos(ll); /* find geocentric ecliptic longitude and latitude */ if (rpd < rsn) lam = atan(-1*rpd*sll/(rsn-rpd*cll))+lg+PI; else lam = atan(rsn*sll/(rpd-rsn*cll))+lpd; range (&lam, 2*PI); bet = atan(rpd*spsi*sin(lam-lpd)/(cpsi*rsn*sll)); /* fill in all of op->s_* stuff except s_size and s_mag */ rho; cir_sky (np, lpd, psi, rp, &rho, lam, bet, lsn, rsn, op); /* compute magnitude and size */ if (op->e_mag.whichm == MAG_HG) { /* the H and G parameters from the Astro. Almanac. */ hg_mag (op->e_mag.m1, op->e_mag.m2, rp, rho, rsn, &mag); if (op->e_size) op->s_size = (float)(op->e_size / rho); else op->s_size = (float)(h_albsize (op->e_mag.m1)/rho); } else { /* the g/k model of comets */ gk_mag (op->e_mag.m1, op->e_mag.m2, rp, rho, &mag); op->s_size = (float)(op->e_size / rho); } set_smag (op, mag); return (0); }
/* given the modified JD, mj, find the nutation in obliquity, *deps, and * the nutation in longitude, *dpsi, each in radians. */ void nutation ( double mj, double *deps, /* on input: precision parameter in arc seconds */ double *dpsi) { static double lastmj = -10000, lastdeps, lastdpsi; double T, T2, T3, T10; /* jul cent since J2000 */ double prec; /* series precis in arc sec */ int i, isecul; /* index in term table */ static double delcache[5][2*NUT_MAXMUL+1]; /* cache for multiples of delaunay args * [M',M,F,D,Om][-min*x, .. , 0, .., max*x] * make static to have unfilled fields cleared on init */ if (mj == lastmj) { *deps = lastdeps; *dpsi = lastdpsi; return; } prec = 0.0; #if 0 /* this is if deps should contain a precision value */ prec =* deps; if (prec < 0.0 || prec > 1.0) /* accept only sane value */ prec = 1.0; #endif /* augment for abundance of small terms */ prec *= NUT_SCALE/10; T = (mj - J2000)/36525.; T2 = T * T; T3 = T2 * T; T10 = T/10.; /* calculate delaunay args and place in cache */ for (i = 0; i < 5; ++i) { double x; short j; x = delaunay[i][0] + delaunay[i][1] * T + delaunay[i][2] * T2 + delaunay[i][3] * T3; /* convert to radians */ x /= SECPERCIRC; x -= floor(x); x *= 2.*PI; /* fill cache table */ for (j = 0; j <= 2*NUT_MAXMUL; ++j) delcache[i][j] = (j - NUT_MAXMUL) * x; } /* find dpsi and deps */ lastdpsi = lastdeps = 0.; for (i = isecul = 0; i < NUT_SERIES ; ++i) { double arg = 0., ampsin, ampcos; short j; if (ampconst[i][0] || ampconst[i][1]) { /* take non-secular terms from simple array */ ampsin = ampconst[i][0]; ampcos = ampconst[i][1]; } else { /* secular terms from different array */ ampsin = ampsecul[isecul][1] + ampsecul[isecul][2] * T10; ampcos = ampsecul[isecul][3] + ampsecul[isecul][4] * T10; ++isecul; } for (j = 0; j < 5; ++j) arg += delcache[j][NUT_MAXMUL + multarg[i][j]]; if (fabs(ampsin) >= prec) lastdpsi += ampsin * sin(arg); if (fabs(ampcos) >= prec) lastdeps += ampcos * cos(arg); } /* convert to radians. */ lastdpsi = degrad(lastdpsi/3600./NUT_SCALE); lastdeps = degrad(lastdeps/3600./NUT_SCALE); lastmj = mj; *deps = lastdeps; *dpsi = lastdpsi; }