00001 /* * 00002 * This file is part of the ESO UVES Pipeline * 00003 * Copyright (C) 2004,2005 European Southern Observatory * 00004 * * 00005 * This library is free software; you can redistribute it and/or modify * 00006 * it under the terms of the GNU General Public License as published by * 00007 * the Free Software Foundation; either version 2 of the License, or * 00008 * (at your option) any later version. * 00009 * * 00010 * This program is distributed in the hope that it will be useful, * 00011 * but WITHOUT ANY WARRANTY; without even the implied warranty of * 00012 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * 00013 * GNU General Public License for more details. * 00014 * * 00015 * You should have received a copy of the GNU General Public License * 00016 * along with this program; if not, write to the Free Software * 00017 * Foundation, 51 Franklin St, Fifth Floor, Boston, MA 02111-1307 USA * 00018 * */ 00019 00020 /* 00021 * $Author: amodigli $ 00022 * $Date: 2007/06/06 08:17:33 $ 00023 * $Revision: 1.8 $ 00024 * $Name: uves-3_3_1 $ 00025 * $Log: uves_baryvel.c,v $ 00026 * Revision 1.8 2007/06/06 08:17:33 amodigli 00027 * replace tab with 4 spaces 00028 * 00029 * Revision 1.7 2007/04/24 12:50:29 jmlarsen 00030 * Replaced cpl_propertylist -> uves_propertylist which is much faster 00031 * 00032 * Revision 1.6 2007/03/15 12:33:16 jmlarsen 00033 * Removed redundant explicit array size 00034 * 00035 * Revision 1.5 2006/11/06 15:19:41 jmlarsen 00036 * Removed unused include directives 00037 * 00038 * Revision 1.4 2006/10/05 06:44:58 jmlarsen 00039 * Declared functions static 00040 * 00041 * Revision 1.3 2006/10/04 10:59:04 jmlarsen 00042 * Implemented QC.VRAD parameters 00043 * 00044 * Revision 1.2 2006/10/04 09:55:44 jmlarsen 00045 * Implemented 00046 * 00047 * Revision 1.4 2006/08/17 13:56:52 jmlarsen 00048 * Reduced max line length 00049 * 00050 * Revision 1.3 2005/12/19 16:17:56 jmlarsen 00051 * Replaced bool -> int 00052 * 00053 */ 00054 00055 #ifdef HAVE_CONFIG_H 00056 # include <config.h> 00057 #endif 00058 00059 /*----------------------------------------------------------------------------*/ 00072 /*----------------------------------------------------------------------------*/ 00075 /*----------------------------------------------------------------------------- 00076 Includes 00077 -----------------------------------------------------------------------------*/ 00078 00079 #include <uves_baryvel.h> 00080 00081 #include <uves_pfits.h> 00082 #include <uves_utils.h> 00083 #include <uves_error.h> 00084 #include <uves_msg.h> 00085 00086 #include <cpl.h> 00087 00088 #include <math.h> 00089 00090 /*----------------------------------------------------------------------------- 00091 Local functions 00092 -----------------------------------------------------------------------------*/ 00093 static void deg2dms(double in_val, 00094 double *degs, 00095 double *minutes, 00096 double *seconds); 00097 00098 static void deg2hms(double in_val, 00099 double *hour, 00100 double *min, 00101 double *sec); 00102 00103 static void compxy(double inputr[19], char inputc[4], 00104 double outputr[4], 00105 double utr, double mod_juldat); 00106 00107 static void barvel(double DJE, double DEQ, 00108 double DVELH[4], double DVELB[4]); 00109 00110 00111 /*----------------------------------------------------------------------------*/ 00118 /*----------------------------------------------------------------------------*/ 00119 void 00120 uves_baryvel(const uves_propertylist *raw_header, 00121 double *bary_corr, 00122 double *helio_corr) 00123 { 00124 00125 double outputr[4]; 00126 00127 //inputc(1:3) = "+++" 00128 char inputc[] = "X+++"; /* 0th index not used */ 00129 00130 //define/local rneg/r/1/1 1.0 00131 double rneg = 1.0; 00132 00133 // write/keyw inputr/r/1/18 0.0 all 00134 double inputr[19]; /* Do not use the zeroth element */ 00135 00136 00137 /* 00138 qc_ra = m$value({p1},O_POS(1)) 00139 qc_dec = m$value({p1},O_POS(2)) 00140 qc_geolat = m$value({p1},{h_geolat}) 00141 qc_geolon = m$value({p1},{h_geolon}) 00142 qc_obs_time = m$value({p1},O_TIME(7)) !using an image as input it take the 00143 !date from the descriptor O_TIME(1,2,3) 00144 !and the UT from O_TIME(5) 00145 */ 00146 double qc_ra; 00147 double qc_dec; 00148 double qc_geolat; 00149 double qc_geolon; 00150 00151 double utr; 00152 double mod_juldat; 00153 00154 double ra_hour, ra_min, ra_sec; 00155 double dec_deg, dec_min, dec_sec; 00156 double lat_deg, lat_min, lat_sec; 00157 double lon_deg, lon_min, lon_sec; 00158 00159 check( qc_ra = uves_pfits_get_ra(raw_header), /* in degrees */ 00160 "Error getting object right ascension"); 00161 check( qc_dec = uves_pfits_get_dec(raw_header), 00162 "Error getting object declination"); 00163 00164 check( qc_geolat = uves_pfits_get_geolat(raw_header), 00165 "Error getting telescope latitude"); 00166 check( qc_geolon = uves_pfits_get_geolon(raw_header), 00167 "Error getting telescope longitude"); 00168 00169 /* double qc_obs_time = uves_pfits_get_exptime(raw_header); Not used! */ 00170 00171 check( utr = uves_pfits_get_utc(raw_header), 00172 "Error reading UTC"); 00173 check( mod_juldat = uves_pfits_get_mjdobs(raw_header), 00174 "Error julian date"); 00175 00176 deg2hms(qc_ra, &ra_hour, &ra_min, &ra_sec); 00177 deg2dms(qc_dec, &dec_deg, &dec_min, &dec_sec); 00178 deg2dms(qc_geolat, &lat_deg, &lat_min, &lat_sec); 00179 deg2dms(qc_geolon, &lon_deg, &lon_min, &lon_sec); 00180 00181 // inputr(1) = m$value({p1},o_time(1)) 00182 // inputr(2) = m$value({p1},o_time(2)) 00183 // inputr(3) = m$value({p1},o_time(3)) 00184 // inputr(4) = m$value({p1},o_time(5)) !UT in real hours 00185 // inputr[1] = year; not needed, pass mjd instead 00186 // inputr[2] = month; 00187 // inputr[3] = day; 00188 // inputr[4] = ut_hour; not needed, pass ut instead 00189 // inputr[5] = ut_min; 00190 // inputr[6] = ut_sec; 00191 00192 // write/keyw inputr/r/7/3 {p4} 00193 inputr[7] = lon_deg; 00194 inputr[8] = lon_min; 00195 inputr[9] = lon_sec; 00196 00197 //rneg = (inputr(7)*3600.)+(inputr(8)*60.)+inputr(9) 00198 rneg = (inputr[7]*3600.)+(inputr[8]*60.)+inputr[9]; 00199 //inputc(1:1) = p4(1:1) 00200 inputc[1] = (lon_deg >= 0) ? '+' : '-'; 00201 //if rneg .lt. 0.0 inputc(1:1) = "-" 00202 if (rneg < 0) inputc[1] = '-'; 00203 00204 // write/keyw inputr/r/10/3 {p5},0,0 00205 inputr[10] = lat_deg; 00206 inputr[11] = lat_min; 00207 inputr[12] = lat_sec; 00208 00209 // rneg = (inputr(10)*3600.)+(inputr(11)*60.)+inputr(12) 00210 rneg = (inputr[10]*3600.)+(inputr[11]*60.)+inputr[12]; 00211 // inputc(2:2) = p5(1:1) 00212 inputc[2] = (lat_deg >= 0) ? '+' : '-'; 00213 // if rneg .lt. 0.0 inputc(2:2) = "-" 00214 if (rneg < 0) inputc[2] = '-'; 00215 00216 // write/keyw inputr/r/13/3 {p2},0,0 00217 inputr[13] = ra_hour; 00218 inputr[14] = ra_min; 00219 inputr[15] = ra_sec; 00220 00221 // write/keyw inputr/r/16/3 {p3},0,0 00222 inputr[16] = dec_deg; 00223 inputr[17] = dec_min; 00224 inputr[18] = dec_sec; 00225 00226 // inputc(3:3) = p3(1:1) 00227 inputc[3] = (dec_deg >= 0) ? '+' : '-'; 00228 // rneg = (inputr(16)*3600.)+(inputr(17)*60.)+inputr(18) 00229 rneg = (inputr[16]*3600.)+(inputr[17]*60.)+inputr[18]; 00230 // if rneg .lt. 0.0 inputc(3:3) = "-" 00231 if (rneg < 0) inputc[3] = '-'; 00232 00233 00234 //C INPUTR/R/1/3 date: year,month,day 00235 //C INPUTR/R/4/3 universal time: hour,min,sec 00236 //C INPUTR/R/7/3 EAST longitude of observatory: degree,min,sec !! NOTE 00237 //C INPUTR/R/10/3 latitude of observatory: degree,min,sec 00238 //C INPUTR/R/13/3 right ascension: hour,min,sec 00239 //C INPUTR/R/16/3 declination: degree,min,sec 00240 00241 //write/keyw action BA !indicate barycorr stuff 00242 //run MID_EXE:COMPXY !compute the corrections 00243 00244 compxy(inputr, inputc, outputr, utr, mod_juldat); 00245 00246 // set/format f14.6,g24.12 00247 // uves_msg_debug(" Barycentric correction time: {outputd(1)} day"); 00248 // uves_msg_debug(" Heliocentric correction time: {outputd(2)} day"); 00249 // uves_msg_debug(" "); 00250 uves_msg_debug(" Total barycentric RV correction: %f km/s", outputr[1]); 00251 uves_msg_debug(" Total heliocentric RV correction: %f km/s", outputr[2]); 00252 uves_msg_debug(" (incl. diurnal RV correction of %f km/s)", outputr[3]); 00253 // uves_msg_debug(" "); 00254 // uves_msg_debug("Descriptor O_TIME of image {p1} used for date and UT."); 00255 00256 *bary_corr = outputr[1]; 00257 *helio_corr = outputr[2]; 00258 00259 cleanup: 00260 return; 00261 } 00262 00263 00264 /*----------------------------------------------------------------------------*/ 00286 /*----------------------------------------------------------------------------*/ 00287 static void 00288 compxy(double inputr[19], char inputc[4], 00289 double outputr[4], 00290 double utr, double mod_juldat) 00291 { 00292 00293 // INTEGER IAV,STAT,KUN(1),KNUL,N 00294 // INTEGER MADRID 00295 // 00296 // DOUBLE PRECISION UTR,STR,T0,DL,THETA0,PE,ST0HG,STG,GAST,R1 00297 double STR; 00298 00299 // double utr Not used. Use FITS header value instead 00300 double t0, dl, theta0, pe, st0hg, stg; 00301 // DOUBLE PRECISION JD,JD0H,JD00,ZERO 00302 double jd, jd0h; 00303 // DOUBLE PRECISION DCORB(3),DCORH(3),DVELB(3),DVELH(3) 00304 double dvelb[4], dvelh[4]; 00305 // DOUBLE PRECISION ALP,BCT,BEOV,BERV,DEL,EDV 00306 double alp, del, beov, berv, EDV; 00307 // DOUBLE PRECISION HAR,HCT,HEOV,HERV,PHI,PI 00308 double HAR, phi, heov, herv; 00309 // DOUBLE PRECISION EQX0,EQX1 00310 // DOUBLE PRECISION A0R,A1R,D0R,D1R 00311 // DOUBLE PRECISION DSMALL,DTEMP(3) 00312 // 00313 // REAL DATE0(3),DATE1(3),DATE00(3),A0(3),A1(3),D0(3),D1(3) 00314 // REAL DATE(3),UT(3),OLONG(3),ST(3) 00315 // double ut[4]; 00316 // REAL OLAT(3),ALPHA(3),DELTA(3) 00317 // REAL RBUF(20) 00318 double *rbuf; 00319 // 00320 // CHARACTER ACTIO*2,SIGNS*3,INPSGN*3 00321 char inpsgn[4]; 00322 // 00323 // COMMON /VMR/MADRID(1) 00324 // 00325 // DATA PI /3.1415926535897928D0/ 00326 // DATA DSMALL /1.D-38/ 00327 00328 00329 double *olong, *olat, *alpha, *delta; 00330 00331 //1000 SIGNS = '+++' 00332 char signs[] = "+++"; 00333 00334 // CALL STKRDR('INPUTR',1,20,IAV,RBUF,KUN,KNUL,STAT) 00335 rbuf = inputr; 00336 // CALL STKRDC('INPUTC',1,1,3,IAV,INPSGN,KUN,KNUL,STAT) 00337 inpsgn[1] = inputc[1]; 00338 inpsgn[2] = inputc[2]; 00339 inpsgn[3] = inputc[3]; 00340 00341 00342 // EQUIVALENCE (RBUF(1),DATE(1)),(RBUF(7),OLONG(1)) 00343 // double *date = rbuf + 1 - 1; Not used, use the explicitly passed MJD instead 00344 olong = rbuf + 7 - 1; 00345 // EQUIVALENCE (RBUF(10),OLAT(1)),(RBUF(13),ALPHA(1)) 00346 olat = rbuf + 10 - 1; 00347 alpha = rbuf + 13 - 1; 00348 // EQUIVALENCE (RBUF(16),DELTA(1)) 00349 delta = rbuf + 16 - 1; 00350 00351 00352 00353 // DO 1100 N=1,3 00354 // UT(N) = RBUF(N+3) 00355 //1100 CONTINUE 00356 // for (n = 1; n <= 3; n++) 00357 // { 00358 // ut[n] = rbuf[n+3]; 00359 // } 00360 00361 // ... convert UT to real hours, calculate Julian date 00362 00363 // UTR = UT(1)+UT(2)/60.D0+UT(3)/3600.D0 00364 // utr = ut[1]+ut[2]/60. +ut[3]/3600.; 00365 00366 /* We know this one already but convert seconds -> hours */ 00367 utr /= 3600; 00368 00369 // CALL JULDAT(DATE,UTR,JD) 00370 jd = mod_juldat + 2400000.5; 00371 00372 // ... likewise convert longitude and latitude of observatory to real hours 00373 // ... and degrees, respectively; take care of signs 00374 // ... NOTE: east longitude is assumed for input !! 00375 00376 // IF ((OLONG(1).LT.0.0) .OR. (OLONG(2).LT.0.0) .OR. 00377 // + (OLONG(3).LT.0.0) .OR. (INPSGN(1:1).EQ.'-')) THEN 00378 if (olong[1] < 0 || olong[2] < 0 || 00379 olong[3] < 0 || inpsgn[1] == '-') { 00380 // SIGNS(1:1) = '-' 00381 signs[1] = '-'; 00382 // OLONG(1) = ABS(OLONG(1)) 00383 // OLONG(2) = ABS(OLONG(2)) 00384 // OLONG(3) = ABS(OLONG(3)) 00385 olong[1] = fabs(olong[1]); 00386 olong[2] = fabs(olong[2]); 00387 olong[3] = fabs(olong[3]); 00388 // ENDIF 00389 } 00390 00391 // DL = OLONG(1)+OLONG(2)/60.D0+OLONG(3)/3600.D0 00392 dl = olong[1]+olong[2]/60. +olong[3]/3600.; 00393 00394 // IF (SIGNS(1:1).EQ.'-') DL = -DL ! negative longitude 00395 if (signs[1] == '-') dl = -dl; 00396 00397 // DL = -DL*24.D0/360.D0 ! convert back to west longitude 00398 dl = -dl*24. /360.; 00399 00400 // IF ((OLAT(1).LT.0.0) .OR. (OLAT(2).LT.0.0) .OR. 00401 // + (OLAT(3).LT.0.0) .OR. (INPSGN(2:2).EQ.'-')) THEN 00402 if (olat[1] < 0 || olat[2] < 0 || 00403 olat[3] < 0 || inpsgn[2] == '-') { 00404 // SIGNS(2:2) = '-' 00405 signs[2] = '-'; 00406 00407 // OLAT(1) = ABS(OLAT(1)) 00408 // OLAT(2) = ABS(OLAT(2)) 00409 // OLAT(3) = ABS(OLAT(3)) 00410 olat[1] = fabs(olat[1]); 00411 olat[2] = fabs(olat[2]); 00412 olat[3] = fabs(olat[3]); 00413 // ENDIF 00414 } 00415 00416 // PHI = OLAT(1)+OLAT(2)/60.D0+OLAT(3)/3600.D0 00417 phi = olat[1]+olat[2]/60. +olat[3]/3600.; 00418 00419 // IF (SIGNS(2:2).EQ.'-') PHI = -PHI ! negative latitude 00420 if (signs[2] == '-') phi = -phi; 00421 00422 // PHI = PHI*PI/180.D0 00423 phi = phi*M_PI/180. ; 00424 00425 // ... convert right ascension and declination to real radians 00426 00427 // ALP = (ALPHA(1)*3600D0+ALPHA(2)*60D0+ALPHA(3))*PI /(12.D0*3600.D0) 00428 alp = (alpha[1]*3600. +alpha[2]*60. +alpha[3])*M_PI/(12. *3600. ); 00429 00430 // IF ((DELTA(1).LT.0.0) .OR. (DELTA(2).LT.0.0) .OR. 00431 // + (DELTA(3).LT.0.0) .OR. (INPSGN(3:3).EQ.'-')) THEN 00432 if (delta[1] < 0 || delta[2] < 0 || 00433 delta[3] < 0 || inpsgn[3] == '-') { 00434 // SIGNS(3:3) = '-' 00435 signs[3] = '-'; 00436 // DELTA(1) = ABS(DELTA(1)) 00437 // DELTA(2) = ABS(DELTA(2)) 00438 // DELTA(3) = ABS(DELTA(3)) 00439 delta[1] = fabs(delta[1]); 00440 delta[2] = fabs(delta[2]); 00441 delta[3] = fabs(delta[3]); 00442 // ENDIF 00443 } 00444 00445 // DEL = (DELTA(1)*3600.D0 + DELTA(2)*60.D0 + DELTA(3)) 00446 // + * PI/(3600.D0*180.D0) 00447 del = (delta[1]*3600.0 + delta[2]*60. + delta[3]) 00448 * M_PI/(3600. *180. ); 00449 00450 00451 // IF (SIGNS(3:3).EQ.'-') DEL = -DEL ! negative declination 00452 if (signs[3] == '-') del = - del; 00453 00454 // ... calculate earth's orbital velocity in rectangular coordinates X,Y,Z 00455 // ... for both heliocentric and barycentric frames (DVELH, DVELB) 00456 // ... Note that setting the second argument of BARVEL to zero as done below 00457 // ... means that the input coordinates will not be corrected for precession. 00458 00459 // CALL BARVEL(JD,0.0D0,DVELH,DVELB) 00460 barvel(jd, 0.0, dvelh, dvelb); 00461 00462 // ... with the rectangular velocity components known, the respective projections 00463 // ... HEOV and BEOV on a given line of sight (ALP,DEL) can be determined: 00464 00465 // ... REFERENCE: THE ASTRONOMICAL ALMANAC 1982 PAGE:B17 00466 00467 // BEOV=DVELB(1)*DCOS(ALP)*DCOS(DEL)+ 00468 // 1 DVELB(2)*DSIN(ALP)*DCOS(DEL)+ 00469 // 2 DVELB(3)*DSIN(DEL) 00470 beov = 00471 dvelb[1]*cos(alp)*cos(del)+ 00472 dvelb[2]*sin(alp)*cos(del)+ 00473 dvelb[3]*sin(del); 00474 00475 // HEOV=DVELH(1)*DCOS(ALP)*DCOS(DEL)+ 00476 // 1 DVELH(2)*DSIN(ALP)*DCOS(DEL)+ 00477 // 2 DVELH(3)*DSIN(DEL) 00478 heov = 00479 dvelh[1]*cos(alp)*cos(del)+ 00480 dvelh[2]*sin(alp)*cos(del)+ 00481 dvelh[3]*sin(del); 00482 00483 00484 // ... For determination also of the contribution due to the diurnal rotation of 00485 // ... the earth (EDV), the hour angle (HAR) is needed at which the observation 00486 // ... was made which requires conversion of UT to sidereal time (ST). 00487 00488 // ... Therefore, first compute ST at 0 hours UT (ST0HG) 00489 00490 // ... REFERENCE : MEEUS J.,1980,ASTRONOMICAL FORMULAE FOR CALCULATORS 00491 00492 // CALL JULDAT(DATE,ZERO,JD0H) 00493 jd0h = jd - (utr/24.0); 00494 00495 // T0=(JD0H-2415020.D0)/36525.D0 00496 t0 = (jd0h-2415020. )/36525. ; 00497 00498 // THETA0=0.276919398D0+100.0021359D0*T0+0.000001075D0*T0*T0 00499 theta0 = 0.276919398 +100.0021359 *t0+0.000001075 *t0*t0 ; 00500 00501 // PE=DINT(THETA0) 00502 pe = (int) theta0; 00503 00504 // THETA0=THETA0-PE 00505 theta0 = theta0 - pe; 00506 00507 // ST0HG=THETA0*24.D0 00508 st0hg = theta0*24. ; 00509 00510 // ... now do the conversion UT -> ST (MEAN SIDEREAL TIME) 00511 00512 // ... REFERENCE : THE ASTRONOMICAL ALMANAC 1983, P B7 00513 // ... IN 1983: 1 MEAN SOLAR DAY = 1.00273790931 MEAN SIDEREAL DAYS 00514 // ... ST WITHOUT EQUATION OF EQUINOXES CORRECTION => ACCURACY +/- 1 SEC 00515 // 00516 // STG=ST0HG+UTR*1.00273790931D0 00517 stg = st0hg+utr*1.00273790931 ; 00518 00519 // IF (STG.LT.DL) STG=STG+24.D0 00520 if (stg < dl) stg = stg +24. ; 00521 00522 // STR=STG-DL 00523 STR = stg-dl; 00524 00525 // IF (STR.GE.24.D0) STR=STR-24.D0 00526 if (STR >= 24. ) STR = STR-24. ; 00527 00528 // STR = STR*PI/12.D0 ! ST in radians 00529 STR = STR*M_PI/12. ; 00530 00531 // HAR=STR-ALP ! hour angle of observation 00532 HAR = STR-alp; 00533 00534 // EDV=-0.4654D0*DSIN(HAR)*DCOS(DEL)*DCOS(PHI) 00535 EDV = -0.4654 * sin(HAR)* cos(del)* cos(phi); 00536 00537 // ... the total correction (in km/s) is the sum of orbital and diurnal components 00538 00539 // HERV=HEOV+EDV 00540 herv=heov+EDV; 00541 // BERV=BEOV+EDV 00542 berv=beov+EDV; 00543 00544 /* The following is not needed. Do not translate */ 00545 00546 #if 0 00547 // ... Calculation of the barycentric and heliocentric correction times 00548 // ... (BCT and HCT) requires knowledge of the earth's position in its 00549 // ... orbit. Subroutine BARCOR returns the rectangular barycentric (DCORB) 00550 // ... and heliocentric (DCORH) coordinates. 00551 00552 // CALL BARCOR(DCORH,DCORB) 00553 00554 // ... from this, the correction times (in days) can be determined: 00555 // ... (REFERENCE: THE ASTRONOMICAL ALMANAC 1982 PAGE:B16) 00556 00557 // BCT=+0.0057756D0*(DCORB(1)*DCOS(ALP)*DCOS(DEL)+ 00558 // 1 DCORB(2)*DSIN(ALP)*DCOS(DEL)+ 00559 // 2 DCORB(3)* DSIN(DEL)) 00560 // HCT=+0.0057756D0*(DCORH(1)*DCOS(ALP)*DCOS(DEL)+ 00561 // 1 DCORH(2)*DSIN(ALP)*DCOS(DEL)+ 00562 // 2 DCORH(3)* DSIN(DEL)) 00563 00564 //... write results to keywords 00565 00566 // CALL STKWRD('OUTPUTD',BCT,1,1,KUN,STAT) ! barycentric correction time 00567 // CALL STKWRD('OUTPUTD',HCT,2,1,KUN,STAT) ! heliocentric correction time 00568 #endif 00569 00570 00571 // RBUF(1) = BERV ! barocentric RV correction 00572 // RBUF(2) = HERV ! heliocentric RV correction 00573 // ... (note that EDV is already contained in both BERV and HERV) 00574 // RBUF(3) = EDV ! diurnal RV correction 00575 rbuf[1] = berv; 00576 rbuf[2] = herv; 00577 rbuf[3] = EDV; 00578 00579 // CALL STKWRR('OUTPUTR',RBUF,1,3,KUN,STAT) 00580 outputr[1] = rbuf[1]; 00581 outputr[2] = rbuf[2]; 00582 outputr[3] = rbuf[3]; 00583 // GOTO 9000 00584 return; 00585 } 00586 00587 /* @cond Convert FORTRAN indexing -> C indexing */ 00588 #define DCFEL(x,y) dcfel[y][x] 00589 #define DCFEPS(x,y) dcfeps[y][x] 00590 #define CCSEL(x,y) ccsel[y][x] 00591 #define DCARGS(x,y) dcargs[y][x] 00592 #define CCAMPS(x,y) ccamps[y][x] 00593 #define CCSEC(x,y) ccsec[y][x] 00594 #define DCARGM(x,y) dcargm[y][x] 00595 #define CCAMPM(x,y) ccampm[y][x] 00596 #define DCEPS(x) dceps[x] 00597 #define FORBEL(x) forbel[x] 00598 #define SORBEL(x) sorbel[x] 00599 #define SN(x) sn[x] 00600 #define SINLP(x) sinlp[x] 00601 #define COSLP(x) coslp[x] 00602 #define CCPAMV(x) ccpamv[x] 00603 /* @endcond */ 00604 /*----------------------------------------------------------------------------*/ 00617 /*----------------------------------------------------------------------------*/ 00618 00619 // SUBROUTINE BARVEL(DJE,DEQ,DVELH,DVELB) 00620 00621 static 00622 void barvel(double DJE, double DEQ, 00623 double DVELH[4], double DVELB[4]) 00624 { 00625 // DOUBLE PRECISION DJE,DEQ,DVELH(3),DVELB(3),SN(4) 00626 double sn[5]; 00627 // DOUBLE PRECISION DT,DTL,DCT0,DCJUL,DTSQ,DLOCAL,DC2PI,CC2PI 00628 double DT,DTL,DTSQ,DLOCAL; 00629 // DOUBLE PRECISION DRD,DRLD,DCSLD,DC1 00630 double DRD,DRLD; 00631 // DOUBLE PRECISION DXBD,DYBD,DZBD,DZHD,DXHD,DYHD 00632 double DXBD,DYBD,DZBD,DZHD,DXHD,DYHD; 00633 // DOUBLE PRECISION DYAHD,DZAHD,DYABD,DZABD 00634 double DYAHD,DZAHD,DYABD,DZABD; 00635 // DOUBLE PRECISION DML,DEPS,PHI,PHID,PSID,DPARAM,PARAM 00636 double DML,DEPS,PHI,PHID,PSID,DPARAM,PARAM; 00637 // DOUBLE PRECISION CCFDI,CCKM,CCMLD,PLON,POMG,PECC 00638 double PLON,POMG,PECC; 00639 // DOUBLE PRECISION PERTL,PERTLD,PERTRD,PERTP,PERTR,PERTPD 00640 double PERTL,PERTLD,PERTRD,PERTP,PERTR,PERTPD; 00641 // DOUBLE PRECISION SINA,CCSGD,DC1MME,TL 00642 double SINA,TL; 00643 // DOUBLE PRECISION CCSEC3,COSA,ESQ 00644 double COSA,ESQ; 00645 // DOUBLE PRECISION DCFEL(3,8),DCEPS(3),CCSEL(3,17),DCARGS(2,15) 00646 // DOUBLE PRECISION CCAMPS(5,15),CCSEC(3,4),DCARGM(2,3) 00647 // DOUBLE PRECISION CCAMPM(4,3),CCPAMV(4) 00648 // DOUBLE PRECISION A,B,E,F,G,SINF,COSF,T,TSQ,TWOE,TWOG 00649 double A,B,F,SINF,COSF,T,TSQ,TWOE,TWOG; 00650 //C 00651 // DOUBLE PRECISION DPREMA(3,3),DPSI,D1PDRO,DSINLS 00652 double DPSI,D1PDRO,DSINLS; 00653 // DOUBLE PRECISION DCOSLS,DSINEP,DCOSEP 00654 double DCOSLS,DSINEP,DCOSEP; 00655 // DOUBLE PRECISION FORBEL(7),SORBEL(17),SINLP(4),COSLP(4) 00656 double forbel[8], sorbel[18], sinlp[5], coslp[5]; 00657 // DOUBLE PRECISION SINLM,COSLM,SIGMA 00658 double SINLM,COSLM,SIGMA; 00659 //C 00660 // INTEGER IDEQ,K,N 00661 int IDEQ,K,N; 00662 //C 00663 // COMMON /BARXYZ/ DPREMA,DPSI,D1PDRO,DSINLS,DCOSLS, 00664 // + DSINEP,DCOSEP,FORBEL,SORBEL,SINLP, 00665 // + COSLP,SINLM,COSLM,SIGMA,IDEQ 00666 00667 // EQUIVALENCE (SORBEL(1),E),(FORBEL(1),G) 00668 double *E = sorbel + 1 - 1; 00669 double *G = forbel + 1 - 1; 00670 //C 00671 // DATA DC2PI/6.2831853071796D0/,CC2PI/6.283185/, 00672 double DC2PI = 6.2831853071796E0; 00673 double CC2PI = 6.283185; /* ??? */ 00674 00675 // *DC1/1.0D0/,DCT0/2415020.0D0/,DCJUL/36525.0D0/ 00676 double DC1 = 1.0; 00677 double DCT0 = 2415020.0E0; 00678 double DCJUL = 36525.0E0; 00679 //C 00680 // DATA DCFEL/ 1.7400353D+00, 6.2833195099091D+02, 5.2796D-06, 00681 // * 6.2565836D+00, 6.2830194572674D+02,-2.6180D-06, 00682 // * 4.7199666D+00, 8.3997091449254D+03,-1.9780D-05, 00683 // * 1.9636505D-01, 8.4334662911720D+03,-5.6044D-05, 00684 // * 4.1547339D+00, 5.2993466764997D+01, 5.8845D-06, 00685 // * 4.6524223D+00, 2.1354275911213D+01, 5.6797D-06, 00686 // * 4.2620486D+00, 7.5025342197656D+00, 5.5317D-06, 00687 // * 1.4740694D+00, 3.8377331909193D+00, 5.6093D-06/ 00688 00689 double dcfel[][4] = { {0, 0, 0, 0}, 00690 {0, 1.7400353E+00, 6.2833195099091E+02, 5.2796E-06}, 00691 {0, 6.2565836E+00, 6.2830194572674E+02,-2.6180E-06}, 00692 {0, 4.7199666E+00, 8.3997091449254E+03,-1.9780E-05}, 00693 {0, 1.9636505E-01, 8.4334662911720E+03,-5.6044E-05}, 00694 {0, 4.1547339E+00, 5.2993466764997E+01, 5.8845E-06}, 00695 {0, 4.6524223E+00, 2.1354275911213E+01, 5.6797E-06}, 00696 {0, 4.2620486E+00, 7.5025342197656E+00, 5.5317E-06}, 00697 {0, 1.4740694E+00, 3.8377331909193E+00, 5.6093E-06} }; 00698 00699 //C 00700 // DATA DCEPS/ 4.093198D-01,-2.271110D-04,-2.860401D-08/ 00701 double dceps[4] = {0, 4.093198E-01,-2.271110E-04,-2.860401E-08}; 00702 00703 //C 00704 // DATA CCSEL/ 1.675104D-02,-4.179579D-05,-1.260516D-07, 00705 // * 2.220221D-01, 2.809917D-02, 1.852532D-05, 00706 // * 1.589963D+00, 3.418075D-02, 1.430200D-05, 00707 // * 2.994089D+00, 2.590824D-02, 4.155840D-06, 00708 // * 8.155457D-01, 2.486352D-02, 6.836840D-06, 00709 // * 1.735614D+00, 1.763719D-02, 6.370440D-06, 00710 // * 1.968564D+00, 1.524020D-02,-2.517152D-06, 00711 // * 1.282417D+00, 8.703393D-03, 2.289292D-05, 00712 // * 2.280820D+00, 1.918010D-02, 4.484520D-06, 00713 // * 4.833473D-02, 1.641773D-04,-4.654200D-07, 00714 // * 5.589232D-02,-3.455092D-04,-7.388560D-07, 00715 // * 4.634443D-02,-2.658234D-05, 7.757000D-08, 00716 // * 8.997041D-03, 6.329728D-06,-1.939256D-09, 00717 // * 2.284178D-02,-9.941590D-05, 6.787400D-08, 00718 // * 4.350267D-02,-6.839749D-05,-2.714956D-07, 00719 // * 1.348204D-02, 1.091504D-05, 6.903760D-07, 00720 // * 3.106570D-02,-1.665665D-04,-1.590188D-07/ 00721 00722 double ccsel[][4] = { {0, 0, 0, 0}, 00723 {0, 1.675104E-02, -4.179579E-05, -1.260516E-07}, 00724 {0, 2.220221E-01, 2.809917E-02, 1.852532E-05}, 00725 {0, 1.589963E+00, 3.418075E-02, 1.430200E-05}, 00726 {0, 2.994089E+00, 2.590824E-02, 4.155840E-06}, 00727 {0, 8.155457E-01, 2.486352E-02, 6.836840E-06}, 00728 {0, 1.735614E+00, 1.763719E-02, 6.370440E-06}, 00729 {0, 1.968564E+00, 1.524020E-02, -2.517152E-06}, 00730 {0, 1.282417E+00, 8.703393E-03, 2.289292E-05}, 00731 {0, 2.280820E+00, 1.918010E-02, 4.484520E-06}, 00732 {0, 4.833473E-02, 1.641773E-04, -4.654200E-07}, 00733 {0, 5.589232E-02, -3.455092E-04, -7.388560E-07}, 00734 {0, 4.634443E-02, -2.658234E-05, 7.757000E-08}, 00735 {0, 8.997041E-03, 6.329728E-06, -1.939256E-09}, 00736 {0, 2.284178E-02, -9.941590E-05, 6.787400E-08}, 00737 {0, 4.350267E-02, -6.839749E-05, -2.714956E-07}, 00738 {0, 1.348204E-02, 1.091504E-05, 6.903760E-07}, 00739 {0, 3.106570E-02, -1.665665E-04, -1.590188E-07} }; 00740 00741 00742 00743 // DATA DCARGS/ 5.0974222D+00,-7.8604195454652D+02, 00744 // * 3.9584962D+00,-5.7533848094674D+02, 00745 // * 1.6338070D+00,-1.1506769618935D+03, 00746 // * 2.5487111D+00,-3.9302097727326D+02, 00747 // * 4.9255514D+00,-5.8849265665348D+02, 00748 // * 1.3363463D+00,-5.5076098609303D+02, 00749 // * 1.6072053D+00,-5.2237501616674D+02, 00750 // * 1.3629480D+00,-1.1790629318198D+03, 00751 // * 5.5657014D+00,-1.0977134971135D+03, 00752 // * 5.0708205D+00,-1.5774000881978D+02, 00753 // * 3.9318944D+00, 5.2963464780000D+01, 00754 // * 4.8989497D+00, 3.9809289073258D+01, 00755 // * 1.3097446D+00, 7.7540959633708D+01, 00756 // * 3.5147141D+00, 7.9618578146517D+01, 00757 // * 3.5413158D+00,-5.4868336758022D+02/ 00758 00759 double dcargs[][3] = { {0, 0, 0}, 00760 {0, 5.0974222E+00, -7.8604195454652E+02}, 00761 {0, 3.9584962E+00, -5.7533848094674E+02}, 00762 {0, 1.6338070E+00, -1.1506769618935E+03}, 00763 {0, 2.5487111E+00, -3.9302097727326E+02}, 00764 {0, 4.9255514E+00, -5.8849265665348E+02}, 00765 {0, 1.3363463E+00, -5.5076098609303E+02}, 00766 {0, 1.6072053E+00, -5.2237501616674E+02}, 00767 {0, 1.3629480E+00, -1.1790629318198E+03}, 00768 {0, 5.5657014E+00, -1.0977134971135E+03}, 00769 {0, 5.0708205E+00, -1.5774000881978E+02}, 00770 {0, 3.9318944E+00, 5.2963464780000E+01}, 00771 {0, 4.8989497E+00, 3.9809289073258E+01}, 00772 {0, 1.3097446E+00, 7.7540959633708E+01}, 00773 {0, 3.5147141E+00, 7.9618578146517E+01}, 00774 {0, 3.5413158E+00, -5.4868336758022E+02} }; 00775 00776 // DATA CCAMPS/ 00777 // *-2.279594D-5, 1.407414D-5, 8.273188D-6, 1.340565D-5,-2.490817D-7, 00778 // *-3.494537D-5, 2.860401D-7, 1.289448D-7, 1.627237D-5,-1.823138D-7, 00779 // * 6.593466D-7, 1.322572D-5, 9.258695D-6,-4.674248D-7,-3.646275D-7, 00780 // * 1.140767D-5,-2.049792D-5,-4.747930D-6,-2.638763D-6,-1.245408D-7, 00781 // * 9.516893D-6,-2.748894D-6,-1.319381D-6,-4.549908D-6,-1.864821D-7, 00782 // * 7.310990D-6,-1.924710D-6,-8.772849D-7,-3.334143D-6,-1.745256D-7, 00783 // *-2.603449D-6, 7.359472D-6, 3.168357D-6, 1.119056D-6,-1.655307D-7, 00784 // *-3.228859D-6, 1.308997D-7, 1.013137D-7, 2.403899D-6,-3.736225D-7, 00785 // * 3.442177D-7, 2.671323D-6, 1.832858D-6,-2.394688D-7,-3.478444D-7, 00786 // * 8.702406D-6,-8.421214D-6,-1.372341D-6,-1.455234D-6,-4.998479D-8, 00787 // *-1.488378D-6,-1.251789D-5, 5.226868D-7,-2.049301D-7, 0.0D0, 00788 // *-8.043059D-6,-2.991300D-6, 1.473654D-7,-3.154542D-7, 0.0D0, 00789 // * 3.699128D-6,-3.316126D-6, 2.901257D-7, 3.407826D-7, 0.0D0, 00790 // * 2.550120D-6,-1.241123D-6, 9.901116D-8, 2.210482D-7, 0.0D0, 00791 // *-6.351059D-7, 2.341650D-6, 1.061492D-6, 2.878231D-7, 0.0D0/ 00792 00793 double ccamps[][6] = 00794 {{0, 0, 0, 0, 0, 0}, 00795 {0, -2.279594E-5, 1.407414E-5, 8.273188E-6, 1.340565E-5, -2.490817E-7}, 00796 {0, -3.494537E-5, 2.860401E-7, 1.289448E-7, 1.627237E-5, -1.823138E-7}, 00797 {0, 6.593466E-7, 1.322572E-5, 9.258695E-6, -4.674248E-7, -3.646275E-7}, 00798 {0, 1.140767E-5, -2.049792E-5, -4.747930E-6, -2.638763E-6, -1.245408E-7}, 00799 {0, 9.516893E-6, -2.748894E-6, -1.319381E-6, -4.549908E-6, -1.864821E-7}, 00800 {0, 7.310990E-6, -1.924710E-6, -8.772849E-7, -3.334143E-6, -1.745256E-7}, 00801 {0, -2.603449E-6, 7.359472E-6, 3.168357E-6, 1.119056E-6, -1.655307E-7}, 00802 {0, -3.228859E-6, 1.308997E-7, 1.013137E-7, 2.403899E-6, -3.736225E-7}, 00803 {0, 3.442177E-7, 2.671323E-6, 1.832858E-6, -2.394688E-7, -3.478444E-7}, 00804 {0, 8.702406E-6, -8.421214E-6, -1.372341E-6, -1.455234E-6, -4.998479E-8}, 00805 {0, -1.488378E-6, -1.251789E-5, 5.226868E-7, -2.049301E-7, 0.0E0}, 00806 {0, -8.043059E-6, -2.991300E-6, 1.473654E-7, -3.154542E-7, 0.0E0}, 00807 {0, 3.699128E-6, -3.316126E-6, 2.901257E-7, 3.407826E-7, 0.0E0}, 00808 {0, 2.550120E-6, -1.241123E-6, 9.901116E-8, 2.210482E-7, 0.0E0}, 00809 {0, -6.351059E-7, 2.341650E-6, 1.061492E-6, 2.878231E-7, 0.0E0}}; 00810 00811 00812 // DATA CCSEC3/-7.757020D-08/ 00813 double CCSEC3 = -7.757020E-08; 00814 //C 00815 // DATA CCSEC/ 1.289600D-06, 5.550147D-01, 2.076942D+00, 00816 // * 3.102810D-05, 4.035027D+00, 3.525565D-01, 00817 // * 9.124190D-06, 9.990265D-01, 2.622706D+00, 00818 // * 9.793240D-07, 5.508259D+00, 1.559103D+01/ 00819 00820 double ccsec[][4] = { {0, 0, 0, 0}, 00821 {0, 1.289600E-06, 5.550147E-01, 2.076942E+00}, 00822 {0, 3.102810E-05, 4.035027E+00, 3.525565E-01}, 00823 {0, 9.124190E-06, 9.990265E-01, 2.622706E+00}, 00824 {0, 9.793240E-07, 5.508259E+00, 1.559103E+01}}; 00825 00826 //C 00827 // DATA DCSLD/ 1.990987D-07/, CCSGD/ 1.990969D-07/ 00828 double DCSLD = 1.990987E-07, CCSGD = 1.990969E-07; 00829 //C 00830 // DATA CCKM/3.122140D-05/, CCMLD/2.661699D-06/, CCFDI/2.399485D-07/ 00831 double CCKM = 3.122140E-05, CCMLD = 2.661699E-06, CCFDI = 2.399485E-07; 00832 //C 00833 // DATA DCARGM/ 5.1679830D+00, 8.3286911095275D+03, 00834 // * 5.4913150D+00,-7.2140632838100D+03, 00835 // * 5.9598530D+00, 1.5542754389685D+04/ 00836 00837 double dcargm[][3] = {{0, 0, 0}, 00838 {0, 5.1679830E+00, 8.3286911095275E+03}, 00839 {0, 5.4913150E+00, -7.2140632838100E+03}, 00840 {0, 5.9598530E+00, 1.5542754389685E+04}}; 00841 //C 00842 // DATA CCAMPM/ 00843 // * 1.097594D-01, 2.896773D-07, 5.450474D-02, 1.438491D-07, 00844 // * -2.223581D-02, 5.083103D-08, 1.002548D-02,-2.291823D-08, 00845 // * 1.148966D-02, 5.658888D-08, 8.249439D-03, 4.063015D-08/ 00846 00847 double ccampm[][5] = {{0, 0, 0, 0, 0}, 00848 {0, 1.097594E-01, 2.896773E-07, 5.450474E-02, 1.438491E-07}, 00849 {0, -2.223581E-02, 5.083103E-08, 1.002548E-02, -2.291823E-08}, 00850 {0, 1.148966E-02, 5.658888E-08, 8.249439E-03, 4.063015E-08} }; 00851 00852 //C 00853 // DATA CCPAMV/8.326827D-11,1.843484D-11,1.988712D-12,1.881276D-12/, 00854 double ccpamv[] = {0, 8.326827E-11, 1.843484E-11, 1.988712E-12, 1.881276E-12}; 00855 // * DC1MME/0.99999696D0/ 00856 double DC1MME = 0.99999696E0; 00857 //C 00858 00859 // IDEQ=DEQ 00860 IDEQ=DEQ; 00861 00862 // DT=(DJE-DCT0)/DCJUL 00863 DT=(DJE-DCT0)/DCJUL; 00864 00865 // T=DT 00866 T=DT; 00867 00868 // DTSQ=DT*DT 00869 DTSQ=DT*DT; 00870 00871 // TSQ=DTSQ 00872 TSQ=DTSQ; 00873 00874 DML = 0; /* Suppress warning */ 00875 // DO 100, K=1,8 00876 for (K = 1; K <= 8; K++) { 00877 00878 // DLOCAL=DMOD(DCFEL(1,K)+DT*DCFEL(2,K)+DTSQ*DCFEL(3,K),DC2PI) 00879 DLOCAL=fmod(DCFEL(1,K)+DT*DCFEL(2,K)+DTSQ*DCFEL(3,K),DC2PI); 00880 00881 // IF (K.EQ.1) DML=DLOCAL 00882 if (K == 1) DML=DLOCAL; 00883 00884 // IF (K.NE.1) FORBEL(K-1)=DLOCAL 00885 if (K != 1) FORBEL(K-1)=DLOCAL; 00886 // 100 CONTINUE 00887 } 00888 00889 // DEPS=DMOD(DCEPS(1)+DT*DCEPS(2)+DTSQ*DCEPS(3), DC2PI) 00890 DEPS=fmod(DCEPS(1)+DT*DCEPS(2)+DTSQ*DCEPS(3), DC2PI); 00891 00892 // DO 200, K=1,17 00893 for (K = 1; K <= 17; K++) { 00894 00895 // SORBEL(K)=DMOD(CCSEL(1,K)+T*CCSEL(2,K)+TSQ*CCSEL(3,K),CC2PI) 00896 SORBEL(K)=fmod(CCSEL(1,K)+T*CCSEL(2,K)+TSQ*CCSEL(3,K),CC2PI); 00897 00898 // 200 CONTINUE 00899 } 00900 00901 // DO 300, K=1,4 00902 for (K = 1; K <= 4; K++) { 00903 00904 // A=DMOD(CCSEC(2,K)+T*CCSEC(3,K),CC2PI) 00905 A=fmod(CCSEC(2,K)+T*CCSEC(3,K),CC2PI); 00906 00907 // SN(K)=DSIN(A) 00908 SN(K)=sin(A); 00909 // 300 CONTINUE 00910 } 00911 00912 // PERTL = CCSEC(1,1) *SN(1) +CCSEC(1,2)*SN(2) 00913 // * +(CCSEC(1,3)+T*CCSEC3)*SN(3) +CCSEC(1,4)*SN(4) 00914 00915 PERTL = CCSEC(1,1) *SN(1) +CCSEC(1,2)*SN(2) 00916 +(CCSEC(1,3)+T*CCSEC3)*SN(3) +CCSEC(1,4)*SN(4); 00917 00918 // PERTLD=0.0 00919 // PERTR =0.0 00920 // PERTRD=0.0 00921 PERTLD=0.0; 00922 PERTR =0.0; 00923 PERTRD=0.0; 00924 00925 // DO 400, K=1,15 00926 for (K = 1; K <= 15; K++) { 00927 // A=DMOD(DCARGS(1,K)+DT*DCARGS(2,K), DC2PI) 00928 A=fmod(DCARGS(1,K)+DT*DCARGS(2,K), DC2PI); 00929 00930 // COSA=DCOS(A) 00931 COSA=cos(A); 00932 00933 // SINA=DSIN(A) 00934 SINA=sin(A); 00935 00936 // PERTL =PERTL+CCAMPS(1,K)*COSA+CCAMPS(2,K)*SINA 00937 PERTL =PERTL+CCAMPS(1,K)*COSA+CCAMPS(2,K)*SINA; 00938 00939 // PERTR =PERTR+CCAMPS(3,K)*COSA+CCAMPS(4,K)*SINA; 00940 PERTR =PERTR+CCAMPS(3,K)*COSA+CCAMPS(4,K)*SINA; 00941 00942 // IF (K.GE.11) GO TO 400 00943 if (K >= 11) break; 00944 00945 // PERTLD=PERTLD+(CCAMPS(2,K)*COSA-CCAMPS(1,K)*SINA)*CCAMPS(5,K) 00946 PERTLD=PERTLD+(CCAMPS(2,K)*COSA-CCAMPS(1,K)*SINA)*CCAMPS(5,K); 00947 00948 // PERTRD=PERTRD+(CCAMPS(4,K)*COSA-CCAMPS(3,K)*SINA)*CCAMPS(5,K) 00949 PERTRD=PERTRD+(CCAMPS(4,K)*COSA-CCAMPS(3,K)*SINA)*CCAMPS(5,K); 00950 00951 // 400 CONTINUE 00952 } 00953 00954 // ESQ=E*E 00955 ESQ=E[1]*E[1]; 00956 00957 // DPARAM=DC1-ESQ 00958 DPARAM=DC1-ESQ; 00959 00960 // PARAM=DPARAM 00961 PARAM=DPARAM; 00962 00963 // TWOE=E+E 00964 TWOE=E[1]+E[1]; 00965 00966 // TWOG=G+G 00967 TWOG=G[1]+G[1]; 00968 00969 // PHI=TWOE*((1.0-ESQ*0.125D0)*DSIN(G)+E*0.625D0*DSIN(TWOG) 00970 // * +ESQ*0.5416667D0*DSIN(G+TWOG) ) 00971 00972 PHI=TWOE*((1.0-ESQ*0.125 )*sin(G[1])+E[1]*0.625 *sin(TWOG) 00973 +ESQ*0.5416667 *sin(G[1]+TWOG) ) ; 00974 00975 //F=G+PHI 00976 F=G[1]+PHI; 00977 00978 //SINF=DSIN(F) 00979 SINF=sin(F); 00980 00981 //COSF=DCOS(F) 00982 COSF=cos(F); 00983 00984 //DPSI=DPARAM/(DC1+E*COSF) 00985 DPSI=DPARAM/(DC1+E[1]*COSF); 00986 00987 // PHID=TWOE*CCSGD*((1.0+ESQ*1.5D0)*COSF+E[1]*(1.25D0-SINF*SINF*0.5D0)) 00988 PHID=TWOE*CCSGD*((1.0+ESQ*1.5 )*COSF+E[1]*(1.25 -SINF*SINF*0.5 )); 00989 00990 // PSID=CCSGD*E*SINF/SQRT(PARAM) 00991 PSID=CCSGD*E[1]*SINF/sqrt(PARAM); 00992 00993 // D1PDRO=(DC1+PERTR) 00994 D1PDRO=(DC1+PERTR); 00995 00996 // DRD=D1PDRO*(PSID+DPSI*PERTRD) 00997 DRD=D1PDRO*(PSID+DPSI*PERTRD); 00998 00999 // DRLD=D1PDRO*DPSI*(DCSLD+PHID+PERTLD) 01000 DRLD=D1PDRO*DPSI*(DCSLD+PHID+PERTLD); 01001 01002 // DTL=DMOD(DML+PHI+PERTL, DC2PI) 01003 DTL=fmod(DML+PHI+PERTL, DC2PI); 01004 01005 // DSINLS=DSIN(DTL) 01006 DSINLS=sin(DTL); 01007 01008 // DCOSLS=DCOS(DTL) 01009 DCOSLS=cos(DTL); 01010 01011 // DXHD = DRD*DCOSLS-DRLD*DSINLS 01012 DXHD = DRD*DCOSLS-DRLD*DSINLS; 01013 01014 // DYHD = DRD*DSINLS+DRLD*DCOSLS 01015 DYHD = DRD*DSINLS+DRLD*DCOSLS; 01016 01017 // PERTL =0.0 01018 PERTL =0.0; 01019 // PERTLD=0.0 01020 PERTLD=0.0; 01021 // PERTP =0.0 01022 PERTP =0.0; 01023 // PERTPD=0.0 01024 PERTPD=0.0; 01025 01026 //DO 500 K=1,3 01027 for (K = 1; K <= 3; K++) { 01028 //A=DMOD(DCARGM(1,K)+DT*DCARGM(2,K), DC2PI) 01029 A=fmod(DCARGM(1,K)+DT*DCARGM(2,K), DC2PI); 01030 01031 //SINA =DSIN(A) 01032 SINA =sin(A); 01033 01034 //COSA =DCOS(A) 01035 COSA =cos(A); 01036 01037 //PERTL =PERTL +CCAMPM(1,K)*SINA 01038 PERTL =PERTL +CCAMPM(1,K)*SINA; 01039 01040 //PERTLD=PERTLD+CCAMPM(2,K)*COSA 01041 PERTLD=PERTLD+CCAMPM(2,K)*COSA; 01042 01043 //PERTP =PERTP +CCAMPM(3,K)*COSA 01044 PERTP =PERTP +CCAMPM(3,K)*COSA; 01045 01046 //PERTPD=PERTPD-CCAMPM(4,K)*SINA 01047 PERTPD=PERTPD-CCAMPM(4,K)*SINA; 01048 01049 // 500 CONTINUE 01050 } 01051 01052 //TL=FORBEL(2)+PERTL 01053 TL=FORBEL(2)+PERTL; 01054 01055 // SINLM=DSIN(TL) 01056 SINLM=sin(TL); 01057 01058 // COSLM=DCOS(TL) 01059 COSLM=cos(TL); 01060 01061 // SIGMA=CCKM/(1.0+PERTP) 01062 SIGMA=CCKM/(1.0+PERTP); 01063 01064 // A=SIGMA*(CCMLD+PERTLD) 01065 A=SIGMA*(CCMLD+PERTLD); 01066 01067 // B=SIGMA*PERTPD 01068 B=SIGMA*PERTPD; 01069 01070 // DXHD=DXHD+A*SINLM+B*COSLM 01071 DXHD=DXHD+A*SINLM+B*COSLM; 01072 01073 // DYHD=DYHD-A*COSLM+B*SINLM 01074 DYHD=DYHD-A*COSLM+B*SINLM; 01075 01076 // DZHD= -SIGMA*CCFDI*DCOS(FORBEL(3)) 01077 DZHD= -SIGMA*CCFDI* cos(FORBEL(3)); 01078 01079 // DXBD=DXHD*DC1MME 01080 DXBD=DXHD*DC1MME; 01081 01082 // DYBD=DYHD*DC1MME 01083 DYBD=DYHD*DC1MME; 01084 // DZBD=DZHD*DC1MME 01085 DZBD=DZHD*DC1MME; 01086 01087 // DO 600 K=1,4 01088 for (K = 1; K <= 4; K++) { 01089 01090 //PLON=FORBEL(K+3) 01091 PLON=FORBEL(K+3); 01092 01093 //POMG=SORBEL(K+1) 01094 POMG=SORBEL(K+1); 01095 01096 //PECC=SORBEL(K+9) 01097 PECC=SORBEL(K+9); 01098 01099 //TL=DMOD(PLON+2.0*PECC*DSIN(PLON-POMG), CC2PI) 01100 TL=fmod(PLON+2.0*PECC* sin(PLON-POMG), CC2PI); 01101 01102 //SINLP(K)=DSIN(TL) 01103 SINLP(K)= sin(TL); 01104 01105 //COSLP(K)=DCOS(TL) 01106 COSLP(K)= cos(TL); 01107 01108 //DXBD=DXBD+CCPAMV(K)*(SINLP(K)+PECC*DSIN(POMG)) 01109 DXBD=DXBD+CCPAMV(K)*(SINLP(K)+PECC*sin(POMG)); 01110 01111 //DYBD=DYBD-CCPAMV(K)*(COSLP(K)+PECC*DCOS(POMG)) 01112 DYBD=DYBD-CCPAMV(K)*(COSLP(K)+PECC*cos(POMG)); 01113 01114 //DZBD=DZBD-CCPAMV(K)*SORBEL(K+13)*DCOS(PLON-SORBEL(K+5)) 01115 DZBD=DZBD-CCPAMV(K)*SORBEL(K+13)*cos(PLON-SORBEL(K+5)); 01116 01117 // 600 CONTINUE 01118 } 01119 01120 //DCOSEP=DCOS(DEPS) 01121 DCOSEP=cos(DEPS); 01122 //DSINEP=DSIN(DEPS) 01123 DSINEP=sin(DEPS); 01124 //DYAHD=DCOSEP*DYHD-DSINEP*DZHD 01125 DYAHD=DCOSEP*DYHD-DSINEP*DZHD; 01126 //DZAHD=DSINEP*DYHD+DCOSEP*DZHD 01127 DZAHD=DSINEP*DYHD+DCOSEP*DZHD; 01128 //DYABD=DCOSEP*DYBD-DSINEP*DZBD 01129 DYABD=DCOSEP*DYBD-DSINEP*DZBD; 01130 //DZABD=DSINEP*DYBD+DCOSEP*DZBD 01131 DZABD=DSINEP*DYBD+DCOSEP*DZBD; 01132 01133 //DVELH(1)=DXHD 01134 DVELH[1]=DXHD; 01135 //DVELH(2)=DYAHD 01136 DVELH[2]=DYAHD; 01137 //DVELH(3)=DZAHD 01138 DVELH[3]=DZAHD; 01139 01140 //DVELB(1)=DXBD 01141 DVELB[1]=DXBD; 01142 //DVELB(2)=DYABD 01143 DVELB[2]=DYABD; 01144 //DVELB(3)=DZABD 01145 DVELB[3]=DZABD; 01146 //DO 800 N=1,3 01147 for (N = 1; N <= 3; N++) { 01148 //DVELH(N)=DVELH(N)*1.4959787D8 01149 DVELH[N]=DVELH[N]*1.4959787E8; 01150 //DVELB(N)=DVELB(N)*1.4959787D8 01151 DVELB[N]=DVELB[N]*1.4959787E8; 01152 // 800 CONTINUE 01153 } 01154 // RETURN 01155 return; 01156 } 01157 01158 /*----------------------------------------------------------------------------*/ 01166 /*----------------------------------------------------------------------------*/ 01167 01168 static void 01169 deg2dms(double in_val, 01170 double *degs, 01171 double *minutes, 01172 double *seconds) 01173 { 01174 deg2hms(in_val*15, degs, minutes, seconds); 01175 } 01176 01180 #define MIDAS_BUG 0 01181 /*----------------------------------------------------------------------------*/ 01189 /*----------------------------------------------------------------------------*/ 01190 01191 static void 01192 deg2hms(double in_val, 01193 double *hours, 01194 double *minutes, 01195 double *seconds) 01196 { 01197 // define/parameter p1 ? num "Enter value in deg units" 01198 // define/local in_val/d/1/1 {p1} 01199 //define/local out_val/c/1/80 " " all 01200 //define/local hours/i/1/1 0 01201 //define/local minutes/i/1/1 0 01202 //define/local seconds/d/1/1 0 01203 01204 //define/local tmp/d/1/1 0 01205 double tmp; 01206 //define/local hold/c/1/80 " " all 01207 //define/local sign/c/1/1 " " 01208 01209 char sign; 01210 01211 //hold = "{in_val}" 01212 //if m$index(hold,"-") .gt. 0 then 01213 // in_val = m$abs(in_val) 01214 // sign = "-" 01215 //else 01216 // sign = "+" 01217 //endif 01218 if (in_val < 0) { 01219 in_val = fabs(in_val); 01220 sign = '-'; 01221 } 01222 else { 01223 sign = '+'; 01224 } 01225 01226 //set/format i1 01228 // tmp = in_val / 15 01229 tmp = in_val / 15; 01230 01231 // hours = tmp !takes the integer part = hours 01232 #if MIDAS_BUG 01233 *hours= uves_round_double(tmp); 01234 #else 01235 *hours= (int) tmp; 01236 #endif 01237 01238 // tmp = tmp - hours !takes the mantissa 01239 tmp = tmp - *hours; 01240 // tmp = tmp * 60 !converts the mantissa in minutes 01241 tmp = tmp * 60; 01242 01243 // minutes = tmp !takes the integer part = minutes 01244 #if MIDAS_BUG 01245 *minutes= uves_round_double(tmp); 01246 #else 01247 *minutes= (int) tmp; 01248 #endif 01249 01250 // tmp = tmp - minutes !takes the mantissa 01251 tmp = tmp - *minutes; 01252 01253 // seconds = tmp * 60 !converts the mantissa in seconds = seconds (with decimal) 01254 *seconds= tmp * 60; 01255 01256 //out_val = "{sign}{hours},{minutes},{seconds}" 01257 01258 /* Rather than returning it explicitly, just attach sign to hours */ 01259 if (sign == '-') *hours = -(*hours); 01260 01261 return; 01262 } 01263 01266 #if 0 /* Not used / needed. 01267 We simply get the julian date from the input FITS header */ 01268 01269 // SUBROUTINE JULDAT(INDATE,UTR,JD) 01270 //C++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 01271 //C 01272 //C.IDENTIFICATION 01273 //C FORTRAN subroutine JULDAT version 1.0 870102 01274 //C original coding: D. Gillet ESO - Garching 01275 //C variables renamed and restructured: D. Baade ST-ECF, Garching 01276 //C 01277 //C.KEYWORDS 01278 //C geocentric Julian date 01279 //C 01280 //C.PURPOSE 01281 //C calculate geocentric Julian date for any civil date (time in UT) 01282 //C 01283 //C.ALGORITHM 01284 //C adapted from MEEUS J.,1980, ASTRONOMICAL FORMULAE FOR CALCULATORS 01285 //C 01286 //C.INPUT/OUTPUT 01287 //C the following are passed from and to the calling program: 01288 //C INDATE(3) : civil date as year,month,day OR year.fraction 01289 //C UT : universal time expressed in real hours 01290 //C JD : real geocentric Julian date 01291 //C 01292 //C.REVISIONS 01293 //C made to accept also REAL dates D. Baade 910408 01294 //C 01295 //C------------------------------------------------------------------------------- 01296 //C 01297 01298 static void 01299 juldat(double *INDATE, 01300 double UTR, 01301 double *JD) 01302 { 01303 // DOUBLE PRECISION YP,P,C,A,UT 01304 double UT; 01305 01306 // DOUBLE PRECISION UTR,JD 01307 //C 01308 // INTEGER STAT,IA,IB,IC,ND,DATE(3) 01309 int DATE[4]; 01310 //C 01311 // REAL INDATE(3),FRAC 01312 //C 01313 01314 // UT=UTR / 24.0D0 01315 UT=UTR / 24.0; 01316 01317 // CHECK FORMAT OF DATE: may be either year,month,date OR year.fraction,0,0 01318 // (Note that the fraction of the year must NOT include fractions of a day.) 01319 // For all other formats exit and terminate also calling command sequence. 01320 // 01321 // IF ((INDATE(1)-INT(INDATE(1))).GT.1.0E-6) THEN 01322 // IF ((INDATE(2).GT.1.0E-6).OR.(INDATE(3).GT.1.0E-6)) 01323 // + CALL STETER(1,'Error: Date was entered in wrong format.') 01324 01325 // copy date input buffer copy to other buffer so that calling program 01326 // does not notice any changes 01327 01328 // FIRST CASE: format was year.fraction 01329 01330 // DATE(1)=INT(INDATE(1)) 01331 // FRAC=INDATE(1)-DATE(1) 01332 // DATE(2)=1 01333 // DATE(3)=1 01334 // ELSE 01335 // 01336 // SECOND CASE: format was year,month,day 01337 // 01338 01339 // DATE(1)=NINT(INDATE(1)) 01340 DATE[1]=uves_round_double(INDATE[1]); 01341 01342 //FRAC=0 01343 FRAC = 0; 01344 01345 //DATE(2)=NINT(INDATE(2)) 01346 DATE[2]=uves_round_double(INDATE[2]); 01347 //DATE(3)=NINT(INDATE(3)) 01348 DATE[3]=uves_round_double(INDATE[3]); 01349 01350 //IF ((DATE(2).EQ.0).AND.(DATE(3).EQ.0)) THEN 01351 if ((DATE[2] == 0) && (DATE[3] == 0)) { 01352 //DATE(2)=1 01353 DATE[2]=1; 01354 //DATE(3)=1 01355 DATE[3]=1; 01356 // ENDIF 01357 } 01358 01359 // IF ((DATE(2).LT.1).OR.(DATE(2).GT.12)) 01360 // + CALL STETER(1,'Error: such a month does not exist') 01361 // IF ((DATE(3).LT.1).OR.(DATE(3).GT.31)) 01362 // + CALL STETER(1,'Error: such a day does not exist') 01363 // ENDIF 01364 01365 // from here on, the normal procedure applies which is based on the 01366 // format year,month,day: 01367 01368 //IF (DATE(2) .GT. 2) THEN 01369 if (DATE[2] > 2) { 01370 //YP=DATE(1) 01371 YP=DATE[1]; 01372 //P=DATE[2] 01373 P=DATE[2]; 01374 // ELSE 01375 } else { 01376 //YP=DATE(1)-1 01377 YP=DATE[1]-1; 01378 //P=DATE(2)+12.0 01379 P=DATE(2)+12.0; 01380 // ENDIF 01381 } 01382 01383 // C = DATE(1) + DATE(2)*1.D-2 + DATE(3)*1.D-4 + UT*1.D-6 01384 C = DATE[1] + DATE[2]*1.E-2 + DATE[3]*1.E-4 + UT*1.E-6; 01385 01386 // IF (C .GE. 1582.1015D0) THEN 01387 if (C > 1582.1015E0) { 01388 //IA=IDINT(YP/100.D0) 01389 IA=(int) (YP/100.D0); 01390 //A=DBLE(IA) 01391 A=IA; 01392 //IB=2-IA+IDINT(A/4.D0) 01393 IB=2-IA+((int)(A/4.D0)); 01394 //ELSE 01395 } else { 01396 //IB=0 01397 IB=0; 01398 //ENDIF 01399 } 01400 01401 // JD = DINT(365.25D0*YP) + DINT(30.6001D0*(P+1.D0)) + DATE(3) + UT 01402 // * + DBLE(IB) + 1720994.5D0 01403 *JD = ((int) (365.25E0*YP)) + ((int)(30.6001D0*(P+1.D0))) + DATE[3] + UT 01404 + IB + 1720994.5E0; 01405 01406 // finally, take into account fraction of year (if any), respect leap 01407 // year conventions 01408 // 01409 // IF (FRAC.GT.1.0E-6) THEN 01410 if (FRAC > 1.0E-6) { 01411 //ND=365 01412 ND=365; 01413 01414 //IF (C.GE.1582.1015D0) THEN 01415 IF (C >= 1582.1015E0) { 01416 //IC = MOD(DATE(1),4) 01417 IC = DATE[1] % 4; 01418 //IF (IC.EQ.0) THEN 01419 if (IC == 0) { 01420 //ND=366 01421 ND=366; 01422 //IC = MOD(DATE(1),100) 01423 IC = DATE[1] % 100; 01424 //IF (IC.EQ.0) THEN 01425 if (IC == 0) { 01426 //IC = MOD(DATE(1),400) 01427 IC = DATE[1] % 400; 01428 //IF (IC.NE.0) ND=365 01429 if (IC != 0) ND=365; 01430 //ENDIF 01431 } 01432 //ENDIF 01433 } 01434 //ENDIF 01435 } 01436 01437 //IF ( ABS(FRAC*ND-NINT(FRAC*ND)).GT.0.3) THEN 01438 if (fabs(FRAC*ND-uves_round_double(FRAC*ND)) > 0.3) { 01439 // CALL STTPUT 01440 // + ('Warning: Fraction of year MAY not correspond to ',STAT) 01441 // CALL STTPUT(' integer number of days.',STAT) 01442 uves_msg_warning("Fraction of year MAY not correspond to " 01443 "integer number of days"); 01444 // ENDIF 01445 } 01446 01447 // JD = JD+NINT(FRAC*ND) 01448 *JD = *JD+uves_round_double(FRAC*ND); 01449 // ENDIF 01450 } 01451 01452 // RETURN 01453 return; 01454 } 01455 #endif
1.4.6