00001 #include <math.h>
00002 #include <stdio.h>
00003 #include <sys/stat.h>
00004 #include "timeio.h"
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018 #define EPS 3.e-17
00019 static double ecanom (double e, double manom) {
00020 int i;
00021 double ecm_new, ecm_old;
00022
00023 ecm_old = manom + e*sin(manom);
00024 i = 0;
00025 for (;;) {
00026 ++i;
00027 ecm_new = ecm_old - (ecm_old - manom - e * sin( ecm_old)) /
00028 (1 - e * cos( ecm_old));
00029 if (fabs((ecm_new - ecm_old)/(ecm_new + ecm_old)) < EPS ||
00030 fabs((ecm_new - ecm_old)) < EPS ) break;
00031 ecm_old = ecm_new;
00032 if (i > 16) {
00033 fprintf (stderr, "ecanom: excessive iterations: %d\nold,new", i);
00034 fprintf (stderr, " %g %g %g\n", ecm_old, ecm_new, ecm_new-ecm_old);
00035 break;
00036 }
00037 }
00038 return (ecm_new);
00039 }
00040
00041 #define EPHEM_OFFSET (-11865398400.0)
00042
00043 #define EPHEM_SIZE 35
00044
00045 #define EPHEM_T 1
00046 #define EPHEM_TCENT 2
00047 #define EPHEM_JD 3
00048 #define EPHEM_ST 4
00049 #define EPHEM_RA 5
00050 #define EPHEM_DEC 6
00051 #define EPHEM_CLONG 7
00052 #define EPHEM_L0 8
00053 #define EPHEM_B0 9
00054 #define EPHEM_P 10
00055
00056 #define EPHEM_DIST 11
00057 #define EPHEM_RSUN 12
00058 #define EPHEM_VEARTH 13
00059 #define EPHEM_EOT 14
00060 #define EPHEM_B0V 15
00061 #define EPHEM_PHI 16
00062 #define EPHEM_SINP 17
00063 #define EPHEM_COSP 18
00064 #define EPHEM_PAROT 19
00065
00066 #define EPHEM_RHO 20
00067 #define EPHEM_L 21
00068 #define EPHEM_B 22
00069 #define EPHEM_R 23
00070 #define EPHEM_SECZ 24
00071 #define EPHEM_SINL 25
00072 #define EPHEM_SINB 26
00073
00074 #define EPHEM_IA 27
00075 #define EPHEM_LD 28
00076 #define EPHEM_AIRM 29
00077 #define EPHEM_DV 30
00078
00079 #define EPHEM_OBS_LON 31
00080 #define EPHEM_OBS_LAT 32
00081
00082
00083 static void solephem (double time, double *ephem) {
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131 double dwork1, dwork2, dsd;
00132 long year,modsd,ly,sd,doy;
00133 double tenths,tenths_sid_mean,tt;
00134 double days_1900,t,mean_sid_time,hams,rams,ra,lo,eot,pi,two_pi,
00135 clong,omega,xx,xs,xc,cosclo,b0,l0,l00,perigee,e,rday,
00136 clong_omega,p,r,rsun,sinincl,cosincl,tanincl,
00137 sinb0,dec,g,ge,v,eps,y;
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147 pi=3.14159265358979324;
00148 two_pi=6.28318530717995865;
00149 rday=two_pi/864000;
00150
00151
00152
00153 sinincl = 0.126198969135829761;
00154 cosincl = 0.992004949679714999;
00155 tanincl = 0.127216068001046929;
00156
00157
00158
00159
00160
00161 tenths = modf((double)time/86400.0, &dsd);
00162 sd = dsd;
00163 tenths *= 864000.0;
00164
00165 modf((double)sd/36524.0, &dwork1);
00166 modf((double)sd/146096.0, &dwork2);
00167 modsd = (long)(sd + dwork1 - dwork2);
00168 ly=modsd/1461;
00169 dwork1 = (modsd+1-1461*ly)/365.3;
00170 if(dwork1 < 0)
00171 {
00172 modf(fabs(dwork1), &dwork1);
00173 dwork1 = -1.0 - dwork1;
00174 }
00175 else
00176 modf(dwork1, &dwork1);
00177 year=4*ly+dwork1;
00178 modf((double)(year*365.25), &dwork1);
00179 doy= modsd + 1 - dwork1;
00180 year=1601+year;
00181
00182
00183
00184
00185
00186
00187 days_1900=dsd-109206.0-0.5;
00188
00189 t=days_1900/36525.0;
00190 tt = t + tenths/(864000.0*36525.0);
00191
00192 e = 0.01675104 - tt*(4.18e-5 + 1.26e-7*tt);
00193
00194 eps = 0.40931975520272993 - tt*(2.271109689158e-4 +
00195 tt*(2.86040072e-8 - tt*8.7751276e-9));
00196
00197 perigee=4.90822946686888692 + tt*( 3.000526416797e-2 +
00198 tt*( 7.9024630021e-6 + tt*5.81776417e-8));
00199
00200 lo = 4.8815279341118791 + tt*(628.33195094248177 + 5.2796209873e-6*tt);
00201
00202
00203 y = tan( (double)(eps/2.0) );
00204 y = y*y;
00205
00206 tenths_sid_mean=239258.36+86401845.42*t+
00207 0.929*t*t+1.0027379093*tenths;
00208
00209 mean_sid_time = two_pi*modf(tenths_sid_mean/864000.0, &dwork2);
00210
00211 hams=(tenths-432000.0)*rday; if( hams < 0 ) hams=hams+two_pi;
00212
00213 rams=mean_sid_time-hams; if( rams < 0 ) rams=rams+two_pi;
00214
00215 g = lo - perigee;
00216
00217 ge = ecanom(e, g);
00218
00219 v = 2.0*atan(sqrt((1.0+e)/(1.0-e))*tan(0.5*ge));
00220
00221 clong = v + perigee;
00222 if(clong > two_pi)clong -= two_pi;
00223 else if (clong < 0)clong += two_pi;
00224
00225 ra = atan2( (1.0-y)*sin(clong),(1+y)*cos(clong));
00226 if(ra < 0)ra += two_pi;
00227
00228 eot = rams - ra;
00229 if(eot > pi)eot -= two_pi;
00230 else if(eot < -pi)eot += two_pi;
00231
00232 dec = asin(sin(eps)*sin(clong));
00233
00234 omega=1.2857258823024895+2.436188747575e-4*(time-7857648000.0)/3.15576e7;
00235
00236 clong_omega=clong-omega;
00237 if( clong_omega < 0 )
00238 clong_omega= clong_omega+two_pi;
00239 else if( clong_omega >= two_pi)
00240 clong_omega= clong_omega-two_pi;
00241 xs=sin(clong_omega); xc= cos(clong_omega);
00242 cosclo = xc= cos(clong_omega);
00243
00244 sinb0= xs*sinincl;
00245 b0 = asin(sinb0);
00246
00247
00248
00249 l00 = atan2(-xs*cosincl,-xc)+5.1097298952004202 + 0.247564432906997103
00250 * (2430000.5 - 2415020.0 - days_1900 - (tenths/864000.0));
00251 l00 = modf(l00/two_pi, &dwork1);
00252 l00 *= two_pi;
00253 if(l00 < 0)l00 += two_pi;
00254 xx = xc * sinincl;
00255
00256
00257
00258
00259
00260
00261
00262 l0= (time - 7983921600.0)/(25.38*86400.0)
00263
00264 + 2.5 - clong_omega/two_pi - (year-1854);
00265 if( doy>90 && clong>omega ) l0=l0-1;
00266
00267 dwork1 = l0;
00268 modf(l0, &l0);
00269 l0 = l0 + (two_pi - l00)/two_pi;
00270
00271 if (dwork1 - l0 > 0.9) l0 = l0 + 1.0;
00272 else if (l0 - dwork1 > 0.9) l0 = l0 - 1.0;
00273
00274
00275
00276 p=atan(-tan(eps)*cos(clong))+atan(-tanincl*xc);
00277 xc=cos(clong-perigee);
00278 r=(1.0-e*e)/(1.0+e*xc);
00279 rsun=959.63/r;
00280
00281 ephem[EPHEM_T] = time;
00282 ephem[EPHEM_TCENT] = t;
00283 ephem[EPHEM_JD] = 2415020.0 + days_1900;
00284 ephem[EPHEM_ST] = mean_sid_time;
00285 ephem[EPHEM_RA] = ra;
00286 ephem[EPHEM_DEC] = dec;
00287 ephem[EPHEM_CLONG] = clong;
00288 ephem[EPHEM_L0] = l0*360;
00289 ephem[EPHEM_B0] = b0;
00290 ephem[EPHEM_P] = p;
00291 ephem[EPHEM_DIST] = r;
00292 ephem[EPHEM_RSUN] = rsun;
00293 ephem[EPHEM_VEARTH] = 29785.9*(1+e*xc);
00294
00295
00296 ephem[EPHEM_EOT] = eot;
00297
00298
00299
00300
00301
00302
00303
00304 ephem[EPHEM_B0V] = xx * ephem[EPHEM_VEARTH];
00305 ephem[EPHEM_PHI] = atan (0.4336*cos(clong)+p);
00306 ephem[EPHEM_SINP] = sin (p);
00307 ephem[EPHEM_COSP] = cos (p);
00308 ephem[EPHEM_PAROT] = atan (-tanincl*cosclo);
00309 return;
00310 }
00311
00312 static void calc_sun_ephemeris (TIME time, double ephem[], double obs_lon, double obs_lat) {
00313 solephem (time - EPHEM_OFFSET, ephem);
00314 ephem[EPHEM_OBS_LON] = obs_lon * M_PI / 180.0;
00315 ephem[EPHEM_OBS_LAT] = obs_lat * M_PI / 180.0;
00316 }
00317
00318 static int sds_flip (void *data, long length, int size) {
00319 long i,j;
00320 char *c, t;
00321 unsigned short *s;
00322 unsigned int *l;
00323
00324 if (data == NULL) return 1;
00325
00326 if (size == 2) {
00327 s = (unsigned short *)data;
00328 for (i=0; i<length; i++, s++) *s = (*s<<8) | (*s>>8);
00329 } else if (size == 4) {
00330 l = (unsigned int *)data;
00331 for (i=0; i<length; i++, l++)
00332 *l = (*l<<24) | ((*l&0xff00)<<8) | ((*l&0xff0000)>>8) | (*l>>24);
00333 } else {
00334 c = (char *)data;
00335 for (i=0; i<length; i++, c+=size)
00336 for (j=0; j<(size/2); j++) {
00337 t = *(c+j);
00338 *(c+j) = *(c+size-j-1);
00339 *(c+size-j-1) = t;
00340 }
00341 }
00342
00343 return 0;
00344 }
00345
00346 #define SUMMARY_FILE ("/home/soi/CM/tables/ephemeris/summary")
00347 #define TABLE_DT (600.0)
00348 #define SUMMARY_REC (6)
00349
00350 int soho_ephemeris (TIME obs_time, double *r, double *b, double *l,
00351 double *vr, double *vn, double *vw, TIME *table_mod_time) {
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363 static FILE *sumfile = NULL;
00364 static int firstcall = 1;
00365 static TIME tab_start, tab_stop;
00366 struct stat stat_buf;
00367 double obs_vel0[SUMMARY_REC], obs_vel1[SUMMARY_REC];
00368 double dt, f0, f1;
00369 int gap = 0, offset, reclen = SUMMARY_REC * sizeof (double);
00370
00371 if (firstcall) {
00372 stat (SUMMARY_FILE, &stat_buf);
00373 *table_mod_time = stat_buf.st_mtime + UNIX_EPOCH;
00374 sumfile = fopen (SUMMARY_FILE, "r");
00375 if (sumfile) {
00376 fread (&tab_start, sizeof (TIME), 1, sumfile);
00377 fread (&tab_stop, sizeof (TIME), 1, sumfile);
00378 #if __BYTE_ORDER == __LITTLE_ENDIAN
00379 sds_flip ((void *)&tab_start, 1, sizeof (TIME));
00380 sds_flip ((void *)&tab_stop, 1, sizeof (TIME));
00381 #endif
00382 }
00383 firstcall = 0;
00384 }
00385
00386 if (!sumfile || obs_time < tab_start || obs_time > (tab_stop + TABLE_DT)) {
00387 double p1, p2;
00388 double ephem[EPHEM_SIZE];
00389 calc_sun_ephemeris (obs_time, ephem, 0.0, 0.0);
00390 *r = 0.99 * ephem[EPHEM_DIST];
00391 *l = 360.0 - fmod (ephem[EPHEM_L0], 360.0);
00392 *b = ephem[EPHEM_B0] * 180.0 / M_PI;
00393 *vw = ephem[EPHEM_VEARTH];
00394 *vn = ephem[EPHEM_B0V];
00395 calc_sun_ephemeris (obs_time - 1.0e5, ephem, 0.0, 0.0);
00396 p1 = ephem[EPHEM_DIST];
00397 calc_sun_ephemeris (obs_time + 1.0e5, ephem, 0.0, 0.0);
00398 p2 = ephem[EPHEM_DIST];
00399 *vr = (1.496e11)*(p2 - p1)/2.0e5;
00400 return (sumfile ? -1 : -2);
00401 }
00402 if (obs_time >= tab_stop) {
00403 gap = 1;
00404 fseek (sumfile, -reclen, SEEK_END);
00405 fread (obs_vel0, sizeof (double), SUMMARY_REC, sumfile);
00406 #if __BYTE_ORDER == __LITTLE_ENDIAN
00407 sds_flip ((void *)obs_vel0, SUMMARY_REC, sizeof (double));
00408 #endif
00409 *vw = obs_vel0[4];
00410 *vn = obs_vel0[5];
00411 *r = obs_vel0[0];
00412 *vr = obs_vel0[3];
00413 *b = obs_vel0[2];
00414 *l = obs_vel0[1];
00415 return gap;
00416 }
00417 offset = reclen * floor ((obs_time - tab_start) / TABLE_DT) +
00418 2 * sizeof (TIME);
00419 fseek (sumfile, offset, SEEK_SET);
00420 fread (obs_vel0, sizeof (double), SUMMARY_REC, sumfile);
00421 fread (obs_vel1, sizeof (double), SUMMARY_REC, sumfile);
00422 #if __BYTE_ORDER == __LITTLE_ENDIAN
00423 sds_flip ((void *)obs_vel0, SUMMARY_REC, sizeof (double));
00424 sds_flip ((void *)obs_vel1, SUMMARY_REC, sizeof (double));
00425 #endif
00426 dt = fmod ((obs_time - tab_start), TABLE_DT) / TABLE_DT;
00427 f1 = dt; f0 = 1.0 - f1;
00428 *r = f0 * obs_vel0[0] + f1 * obs_vel1[0];
00429 while (obs_vel1[1] > obs_vel0[1]) obs_vel1[1] -= 360.0;
00430 *l = f0 * obs_vel0[1] + f1 * obs_vel1[1];
00431 while (*l < 0.0) *l += 360.0;
00432 *b = f0 * obs_vel0[2] + f1 * obs_vel1[2];
00433 *vr = f0 * obs_vel0[3] + f1 * obs_vel1[3];
00434 *vw = f0 * obs_vel0[4] + f1 * obs_vel1[4];
00435 *vn = f0 * obs_vel0[5] + f1 * obs_vel1[5];
00436 if (*r < 0.1) gap = 2;
00437 return gap;
00438 }