sinfo_baryvel.c

00001 /*                                                                            *
00002  *   This file is part of the ESO SINFONI 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: 2012/03/02 08:42:20 $
00023  * $Revision: 1.3 $
00024  * $Name: HEAD $
00025  * $Log: sinfo_baryvel.c,v $
00026  * Revision 1.3  2012/03/02 08:42:20  amodigli
00027  * fixed some typos on doxygen
00028  *
00029  * Revision 1.2  2009/04/28 11:42:18  amodigli
00030  * now return cpl_error_code
00031  *
00032  * Revision 1.1  2009/01/02 08:27:58  amodigli
00033  * added to repository
00034  *
00035  * Revision 1.8  2007/06/06 08:17:33  amodigli
00036  * replace tab with 4 spaces
00037  *
00038  */
00039 
00040 #ifdef HAVE_CONFIG_H
00041 #  include <config.h>
00042 #endif
00043 
00044 
00046 /*---------------------------------------------------------------------------*/
00060 /*---------------------------------------------------------------------------*/
00061 
00062 /*----------------------------------------------------------------------------
00063   Includes
00064   ---------------------------------------------------------------------------*/
00065 
00066 #include <sinfo_baryvel.h>
00067 
00068 #include <sinfo_pfits.h>
00069 #include <sinfo_utils.h>
00070 #include <sinfo_error.h>
00071 #include <sinfo_msg.h>
00072 #include <sinfo_functions.h>
00073 
00074 #include <cpl.h>
00075 
00076 #include <math.h>
00077 
00078 #define H_GEOLAT "ESO TEL GEOLAT"
00079 #define H_GEOLON "ESO TEL GEOLON"
00080 #define H_UTC "UTC"
00081 
00082 /*-----------------------------------------------------------------------------
00083   Local functions
00084   ---------------------------------------------------------------------------*/
00085 
00086 static double sinfo_pfits_get_geolat(const cpl_propertylist * plist);
00087 static double sinfo_pfits_get_geolon(const cpl_propertylist * plist);
00088 static double sinfo_pfits_get_utc(const cpl_propertylist * plist);
00089 
00090 
00091 
00092 static void deg2dms(double in_val, 
00093             double *degs,
00094             double *minutes,
00095             double *seconds);
00096 
00097 static void deg2hms(double in_val, 
00098             double *hour,
00099             double *min,
00100             double *sec);
00101 
00102 static void compxy(double inputr[19], char inputc[4],
00103            double outputr[4],
00104            double utr, double mod_juldat);
00105 
00106 static void barvel(double DJE, double DEQ,
00107            double DVELH[4], double DVELB[4]);
00108 
00109 
00110 
00111 
00112 /*--------------------------------------------------------------------------*/
00113 
00114 /*--------------------------------------------------------------------------*/
00120 /*--------------------------------------------------------------------------*/
00121 static double sinfo_pfits_get_geolat(const cpl_propertylist * plist)
00122 {
00123     double returnvalue = 0;
00124     
00125     check(returnvalue=cpl_propertylist_get_double(plist, H_GEOLAT), 
00126        "Error reading keyword '%s'", H_GEOLAT);
00127     
00128   cleanup:
00129     return returnvalue;
00130 }
00131 
00132 /*--------------------------------------------------------------------------*/
00138 /*--------------------------------------------------------------------------*/
00139 static double sinfo_pfits_get_geolon(const cpl_propertylist * plist)
00140 {
00141     double returnvalue = 0;
00142 
00143     check(returnvalue=cpl_propertylist_get_double(plist, H_GEOLON), 
00144        "Error reading keyword '%s'", H_GEOLON);
00145       
00146   cleanup:
00147     return returnvalue;
00148 }
00149 
00150 
00151 
00152 
00153 /*---------------------------------------------------------------------------*/
00159 /*---------------------------------------------------------------------------*/
00160 static double sinfo_pfits_get_utc(const cpl_propertylist * plist)
00161 {
00162     double returnvalue = 0;
00163 
00164     check(returnvalue=cpl_propertylist_get_double(plist, H_UTC), 
00165        "Error reading keyword '%s'", H_UTC);
00166      
00167   cleanup:
00168     return returnvalue;
00169 }
00170 
00171 
00172 
00173 #if 0   /* Not used / needed.
00174        We simply get the julian date from the input FITS header */
00175 
00176 //      SUBROUTINE JULDAT(INDATE,UTR,JD)
00177 //C++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
00178 //C
00179 //C.IDENTIFICATION
00180 //C  FORTRAN subroutine                    JULDAT     version 1.0       870102
00181 //C  original coding:                      D. Gillet        ESO - Garching
00182 //C  variables renamed and restructured:   D. Baade         ST-ECF, Garching
00183 //C
00184 //C.KEYWORDS
00185 //C  geocentric Julian date
00186 //C
00187 //C.PURPOSE
00188 //C  calculate geocentric Julian date for any civil date (time in UT)
00189 //C
00190 //C.ALGORITHM
00191 //C adapted from MEEUS J.,1980, ASTRONOMICAL FORMULAE FOR CALCULATORS
00192 //C
00193 //C.INPUT/OUTPUT
00194 //C the following are passed from and to the calling program:
00195 //C  INDATE(3)    :         civil date as year,month,day OR year.fraction
00196 //C  UT           :         universal time expressed in real hours
00197 //C  JD           :         real geocentric Julian date
00198 //C
00199 //C.REVISIONS
00200 //C made to accept also REAL dates         D. Baade             910408
00201 //C
00202 //C---------------------------------------------------------------------------
00203 //C
00204 
00205 static void 
00206 juldat(double *INDATE,
00207        double UTR,
00208        double *JD)
00209 {
00210   double UT;
00211 
00212   int DATE[4];
00213 
00214   UT=UTR / 24.0;
00215 
00216   /*
00217     CHECK FORMAT OF DATE: may be either year,month,date OR year.fraction,0,0 
00218     (Note that the fraction of the year must NOT include fractions of a day.)
00219     For all other formats exit and terminate also calling command sequence.
00220   
00221     IF ((INDATE(1)-INT(INDATE(1))).GT.1.0E-6) THEN 
00222     IF ((INDATE(2).GT.1.0E-6).OR.(INDATE(3).GT.1.0E-6)) 
00223     +       CALL   STETER(1,'Error: Date was entered in wrong format.')
00224 
00225     copy date input buffer copy to other buffer so that calling program 
00226     does not notice any changes
00227 
00228     FIRST CASE: format was year.fraction
00229 
00230     DATE(1)=INT(INDATE(1))
00231     FRAC=INDATE(1)-DATE(1)
00232     DATE(2)=1
00233     DATE(3)=1
00234     ELSE
00235   
00236     SECOND CASE: format was year,month,day
00237   */
00238 
00239   DATE[1]=sinfo_round_double(INDATE[1]);
00240 
00241   FRAC = 0;
00242 
00243   DATE[2]=sinfo_round_double(INDATE[2]);
00244 
00245   DATE[3]=sinfo_round_double(INDATE[3]);
00246 
00247   if ((DATE[2] == 0) &&  (DATE[3] == 0)) {
00248 
00249     DATE[2]=1;
00250 
00251     DATE[3]=1;
00252 
00253   }
00254 
00255   /*
00256     from here on, the normal procedure applies which is based on the 
00257     format year,month,day:
00258   */
00259   if (DATE[2] > 2) {
00260     YP=DATE[1];
00261     P=DATE[2];
00262   } else {
00263     YP=DATE[1]-1;
00264     P=DATE(2)+12.0;
00265   }
00266 
00267   C = DATE[1] + DATE[2]*1.E-2 + DATE[3]*1.E-4 + UT*1.E-6;
00268 
00269   if (C  >   1582.1015E0) {
00270     IA=(int) (YP/100.D0);
00271     A=IA;
00272     IB=2-IA+((int)(A/4.D0));
00273   } else {
00274     IB=0;
00275   }
00276 
00277   *JD = ((int) (365.25E0*YP)) + ((int)(30.6001D0*(P+1.D0))) + DATE[3] + UT
00278     + IB + 1720994.5E0;
00279 
00280   /*
00281     finally, take into account fraction of year (if any), respect leap
00282     year conventions
00283   */
00284   if (FRAC > 1.0E-6) {
00285     ND=365;
00286 
00287     IF (C >= 1582.1015E0) {
00288       IC = DATE[1] % 4;
00289       if (IC == 0) {
00290         ND=366;
00291         IC = DATE[1] % 100;
00292         if (IC == 0) {
00293       IC = DATE[1] % 400;
00294       if (IC != 0) ND=365;
00295         }
00296       }
00297     }
00298 
00299     if (fabs(FRAC*ND-sinfo_round_double(FRAC*ND)) > 0.3) {
00300       sinfo_msg_warning("Fraction of year MAY not correspond to "
00301             "integer number of days");
00302     }
00303 
00304     *JD = *JD+sinfo_round_double(FRAC*ND);
00305   }
00306 
00307   return;
00308 }
00309 
00310 #endif
00311 
00315 #define MIDAS_BUG 0
00316 /*---------------------------------------------------------------------------*/
00324 /*---------------------------------------------------------------------------*/
00325 
00326 static void
00327 deg2hms(double in_val, 
00328     double *hours,
00329     double *minutes,
00330     double *seconds)
00331 {
00332   double tmp;
00333   char sign;
00334   if (in_val < 0) {
00335     in_val = fabs(in_val);
00336     sign = '-';
00337   }
00338   else {
00339     sign = '+';
00340   }
00341 
00342   tmp   = in_val / 15;
00343 
00344   /* takes the integer part = hours */
00345 #if MIDAS_BUG
00346   *hours= sinfo_round_double(tmp);
00347 #else
00348   *hours= (int) tmp;
00349 #endif
00350 
00351   /* takes the mantissa */
00352   tmp   = tmp - *hours;
00353   /* converts the mantissa in minutes */
00354   tmp   = tmp * 60;
00355 
00356   /* takes the integer part = minutes */
00357 #if MIDAS_BUG
00358   *minutes= sinfo_round_double(tmp);
00359 #else
00360   *minutes= (int) tmp;
00361 #endif
00362 
00363   /* takes the mantissa */
00364   tmp   = tmp - *minutes;
00365 
00366   /* converts the mantissa in seconds = seconds (with decimal) */
00367   *seconds= tmp * 60;
00368 
00369   /* Rather than returning it explicitly, just  attach sign to hours */
00370   if (sign == '-') *hours = -(*hours);
00371 
00372   return;
00373 }
00374 
00375 /*---------------------------------------------------------------------------*/
00383 /*---------------------------------------------------------------------------*/
00384 
00385 static void
00386 deg2dms(double in_val, 
00387     double *degs,
00388     double *minutes,
00389     double *seconds)
00390 {
00391   deg2hms(in_val*15, degs, minutes, seconds);
00392 }
00393 
00394 
00395 
00396 
00397 
00398 /* @cond Convert FORTRAN indexing -> C indexing */
00399 #define DCFEL(x,y)  dcfel[y][x]
00400 #define DCFEPS(x,y) dcfeps[y][x]
00401 #define CCSEL(x,y)  ccsel[y][x]
00402 #define DCARGS(x,y) dcargs[y][x]
00403 #define CCAMPS(x,y) ccamps[y][x]
00404 #define CCSEC(x,y)  ccsec[y][x]
00405 #define DCARGM(x,y) dcargm[y][x]
00406 #define CCAMPM(x,y) ccampm[y][x]
00407 #define DCEPS(x)    dceps[x]
00408 #define FORBEL(x)   forbel[x]
00409 #define SORBEL(x)   sorbel[x]
00410 #define SN(x)       sn[x]
00411 #define SINLP(x)    sinlp[x]
00412 #define COSLP(x)    coslp[x]
00413 #define CCPAMV(x)   ccpamv[x]
00414 /* @endcond */
00415 /*---------------------------------------------------------------------------*/
00428 /*---------------------------------------------------------------------------*/
00429 
00430 
00431 static 
00432 void barvel(double DJE, double DEQ,
00433         double DVELH[4], double DVELB[4])
00434 {
00435   double sn[5];
00436   double DT,DTL,DTSQ,DLOCAL;
00437   double DRD,DRLD;
00438   double DXBD,DYBD,DZBD,DZHD,DXHD,DYHD;
00439   double DYAHD,DZAHD,DYABD,DZABD;
00440   double DML,DEPS,PHI,PHID,PSID,DPARAM,PARAM;
00441   double PLON,POMG,PECC;
00442   double PERTL,PERTLD,PERTRD,PERTP,PERTR,PERTPD;
00443   double SINA,TL;
00444   double COSA,ESQ;
00445   double A,B,F,SINF,COSF,T,TSQ,TWOE,TWOG;
00446 
00447   double DPSI,D1PDRO,DSINLS;
00448   double DCOSLS,DSINEP,DCOSEP;
00449   double forbel[8], sorbel[18], sinlp[5], coslp[5];
00450   double SINLM,COSLM,SIGMA;
00451   int IDEQ,K,N;
00452 
00453   double *E = sorbel + 1 - 1;
00454   double *G = forbel + 1 - 1;
00455   double DC2PI = 6.2831853071796E0;
00456   double CC2PI = 6.283185;             /* ??? */
00457 
00458   double DC1 = 1.0;
00459   double DCT0 = 2415020.0E0;
00460   double DCJUL = 36525.0E0;
00461 
00462   double dcfel[][4] = { {0, 0, 0, 0},
00463             {0, 1.7400353E+00, 6.2833195099091E+02, 5.2796E-06},
00464             {0, 6.2565836E+00, 6.2830194572674E+02,-2.6180E-06},
00465             {0, 4.7199666E+00, 8.3997091449254E+03,-1.9780E-05},
00466             {0, 1.9636505E-01, 8.4334662911720E+03,-5.6044E-05},
00467             {0, 4.1547339E+00, 5.2993466764997E+01, 5.8845E-06},
00468             {0, 4.6524223E+00, 2.1354275911213E+01, 5.6797E-06},
00469             {0, 4.2620486E+00, 7.5025342197656E+00, 5.5317E-06},
00470             {0, 1.4740694E+00, 3.8377331909193E+00, 5.6093E-06} };
00471     
00472   double dceps[4] = {0, 4.093198E-01,-2.271110E-04,-2.860401E-08};
00473 
00474   double ccsel[][4] = { {0, 0, 0, 0},
00475             {0, 1.675104E-02, -4.179579E-05, -1.260516E-07},
00476             {0, 2.220221E-01,  2.809917E-02,  1.852532E-05},
00477             {0, 1.589963E+00,  3.418075E-02,  1.430200E-05},
00478             {0, 2.994089E+00,  2.590824E-02,  4.155840E-06},
00479             {0, 8.155457E-01,  2.486352E-02,  6.836840E-06},
00480             {0, 1.735614E+00,  1.763719E-02,  6.370440E-06},
00481             {0, 1.968564E+00,  1.524020E-02, -2.517152E-06},
00482             {0, 1.282417E+00,  8.703393E-03,  2.289292E-05},
00483             {0, 2.280820E+00,  1.918010E-02,  4.484520E-06},
00484             {0, 4.833473E-02,  1.641773E-04, -4.654200E-07},
00485             {0, 5.589232E-02, -3.455092E-04, -7.388560E-07},
00486             {0, 4.634443E-02, -2.658234E-05,  7.757000E-08},
00487             {0, 8.997041E-03,  6.329728E-06, -1.939256E-09},
00488             {0, 2.284178E-02, -9.941590E-05,  6.787400E-08},
00489             {0, 4.350267E-02, -6.839749E-05, -2.714956E-07},
00490             {0, 1.348204E-02,  1.091504E-05,  6.903760E-07},
00491             {0, 3.106570E-02, -1.665665E-04, -1.590188E-07} };
00492 
00493 
00494   double dcargs[][3] = { {0, 0, 0},
00495              {0, 5.0974222E+00, -7.8604195454652E+02},
00496              {0, 3.9584962E+00, -5.7533848094674E+02},
00497              {0, 1.6338070E+00, -1.1506769618935E+03},
00498              {0, 2.5487111E+00, -3.9302097727326E+02},
00499              {0, 4.9255514E+00, -5.8849265665348E+02},
00500              {0, 1.3363463E+00, -5.5076098609303E+02},
00501              {0, 1.6072053E+00, -5.2237501616674E+02},
00502              {0, 1.3629480E+00, -1.1790629318198E+03},
00503              {0, 5.5657014E+00, -1.0977134971135E+03},
00504              {0, 5.0708205E+00, -1.5774000881978E+02},
00505              {0, 3.9318944E+00,  5.2963464780000E+01},
00506              {0, 4.8989497E+00,  3.9809289073258E+01},
00507              {0, 1.3097446E+00,  7.7540959633708E+01},
00508              {0, 3.5147141E+00,  7.9618578146517E+01},
00509              {0, 3.5413158E+00, -5.4868336758022E+02} };
00510 
00511 
00512   double ccamps[][6] = 
00513     {{0, 0, 0, 0, 0, 0},
00514      {0, -2.279594E-5,  1.407414E-5,  8.273188E-6,  1.340565E-5, -2.490817E-7},
00515      {0, -3.494537E-5,  2.860401E-7,  1.289448E-7,  1.627237E-5, -1.823138E-7},
00516      {0,  6.593466E-7,  1.322572E-5,  9.258695E-6, -4.674248E-7, -3.646275E-7},
00517      {0,  1.140767E-5, -2.049792E-5, -4.747930E-6, -2.638763E-6, -1.245408E-7},
00518      {0,  9.516893E-6, -2.748894E-6, -1.319381E-6, -4.549908E-6, -1.864821E-7},
00519      {0,  7.310990E-6, -1.924710E-6, -8.772849E-7, -3.334143E-6, -1.745256E-7},
00520      {0, -2.603449E-6,  7.359472E-6,  3.168357E-6,  1.119056E-6, -1.655307E-7},
00521      {0, -3.228859E-6,  1.308997E-7,  1.013137E-7,  2.403899E-6, -3.736225E-7},
00522      {0,  3.442177E-7,  2.671323E-6,  1.832858E-6, -2.394688E-7, -3.478444E-7},
00523      {0,  8.702406E-6, -8.421214E-6, -1.372341E-6, -1.455234E-6, -4.998479E-8},
00524      {0, -1.488378E-6, -1.251789E-5,  5.226868E-7, -2.049301E-7,  0.0E0},
00525      {0, -8.043059E-6, -2.991300E-6,  1.473654E-7, -3.154542E-7,  0.0E0},
00526      {0,  3.699128E-6, -3.316126E-6,  2.901257E-7,  3.407826E-7,  0.0E0},
00527      {0,  2.550120E-6, -1.241123E-6,  9.901116E-8,  2.210482E-7,  0.0E0},
00528      {0, -6.351059E-7,  2.341650E-6,  1.061492E-6,  2.878231E-7,  0.0E0}};
00529 
00530 
00531 
00532   double CCSEC3 = -7.757020E-08;
00533 
00534   double ccsec[][4] = { {0, 0, 0, 0},
00535             {0, 1.289600E-06,  5.550147E-01,  2.076942E+00},
00536             {0, 3.102810E-05,  4.035027E+00,  3.525565E-01},
00537             {0, 9.124190E-06,  9.990265E-01,  2.622706E+00},
00538             {0, 9.793240E-07,  5.508259E+00,  1.559103E+01}};
00539 
00540   double DCSLD =  1.990987E-07, CCSGD = 1.990969E-07;
00541 
00542   double CCKM = 3.122140E-05, CCMLD = 2.661699E-06, CCFDI = 2.399485E-07;
00543 
00544   double dcargm[][3] = {{0, 0, 0},
00545             {0, 5.1679830E+00,  8.3286911095275E+03},
00546             {0, 5.4913150E+00, -7.2140632838100E+03},
00547             {0, 5.9598530E+00,  1.5542754389685E+04}};
00548 
00549   double ccampm[][5] = {{0, 0, 0, 0, 0},
00550             {0,  1.097594E-01,  2.896773E-07,  5.450474E-02,  1.438491E-07},
00551             {0, -2.223581E-02,  5.083103E-08,  1.002548E-02, -2.291823E-08},
00552             {0,  1.148966E-02,  5.658888E-08,  8.249439E-03,  4.063015E-08} };
00553 
00554   double ccpamv[] = {0, 8.326827E-11, 1.843484E-11, 1.988712E-12, 1.881276E-12};
00555 
00556   double DC1MME = 0.99999696E0;
00557 
00558   IDEQ=DEQ;
00559 
00560   DT=(DJE-DCT0)/DCJUL;
00561 
00562   T=DT;
00563 
00564   DTSQ=DT*DT;
00565 
00566   TSQ=DTSQ;
00567 
00568   DML = 0;  /* Suppress warning */
00569   for (K = 1; K <= 8; K++) {
00570 
00571     DLOCAL=fmod(DCFEL(1,K)+DT*DCFEL(2,K)+DTSQ*DCFEL(3,K),DC2PI);
00572 
00573     if (K == 1)  DML=DLOCAL;
00574 
00575     if (K != 1)  FORBEL(K-1)=DLOCAL;
00576   }
00577 
00578   DEPS=fmod(DCEPS(1)+DT*DCEPS(2)+DTSQ*DCEPS(3), DC2PI);
00579 
00580   for (K = 1; K <= 17; K++) {
00581 
00582     SORBEL(K)=fmod(CCSEL(1,K)+T*CCSEL(2,K)+TSQ*CCSEL(3,K),CC2PI);
00583 
00584   }
00585 
00586   for (K = 1; K <= 4; K++) {
00587 
00588     A=fmod(CCSEC(2,K)+T*CCSEC(3,K),CC2PI);
00589 
00590     SN(K)=sin(A);
00591 
00592   }
00593 
00594   PERTL =  CCSEC(1,1)          *SN(1) +CCSEC(1,2)*SN(2)
00595     +(CCSEC(1,3)+T*CCSEC3)*SN(3) +CCSEC(1,4)*SN(4);
00596 
00597   PERTLD=0.0;
00598   PERTR =0.0;
00599   PERTRD=0.0;
00600 
00601   for (K = 1; K <= 15; K++) {
00602 
00603     A=fmod(DCARGS(1,K)+DT*DCARGS(2,K), DC2PI);
00604 
00605     COSA=cos(A);
00606 
00607     SINA=sin(A);
00608 
00609     PERTL =PERTL+CCAMPS(1,K)*COSA+CCAMPS(2,K)*SINA;
00610 
00611     PERTR =PERTR+CCAMPS(3,K)*COSA+CCAMPS(4,K)*SINA;
00612 
00613     if (K >= 11) break;
00614 
00615     PERTLD=PERTLD+(CCAMPS(2,K)*COSA-CCAMPS(1,K)*SINA)*CCAMPS(5,K);
00616 
00617     PERTRD=PERTRD+(CCAMPS(4,K)*COSA-CCAMPS(3,K)*SINA)*CCAMPS(5,K);
00618 
00619   }
00620 
00621 
00622   ESQ=E[1]*E[1];
00623 
00624   DPARAM=DC1-ESQ;
00625 
00626   PARAM=DPARAM;
00627 
00628   TWOE=E[1]+E[1];
00629 
00630   TWOG=G[1]+G[1];
00631 
00632   PHI=TWOE*((1.0-ESQ*0.125  )*sin(G[1])+E[1]*0.625  *sin(TWOG)
00633         +ESQ*0.5416667  *sin(G[1]+TWOG) ) ;
00634     
00635   F=G[1]+PHI;
00636 
00637   SINF=sin(F);
00638 
00639   COSF=cos(F);
00640 
00641   DPSI=DPARAM/(DC1+E[1]*COSF);
00642 
00643   PHID=TWOE*CCSGD*((1.0+ESQ*1.5  )*COSF+E[1]*(1.25  -SINF*SINF*0.5  ));
00644 
00645   PSID=CCSGD*E[1]*SINF/sqrt(PARAM);
00646 
00647   D1PDRO=(DC1+PERTR);
00648 
00649   DRD=D1PDRO*(PSID+DPSI*PERTRD);
00650 
00651   DRLD=D1PDRO*DPSI*(DCSLD+PHID+PERTLD);
00652 
00653   DTL=fmod(DML+PHI+PERTL, DC2PI);
00654 
00655   DSINLS=sin(DTL);
00656 
00657   DCOSLS=cos(DTL);
00658 
00659   DXHD = DRD*DCOSLS-DRLD*DSINLS;
00660 
00661   DYHD = DRD*DSINLS+DRLD*DCOSLS;
00662 
00663   PERTL =0.0;
00664 
00665   PERTLD=0.0;
00666 
00667   PERTP =0.0;
00668 
00669   PERTPD=0.0;
00670 
00671   for (K = 1; K <= 3; K++) {
00672     A=fmod(DCARGM(1,K)+DT*DCARGM(2,K), DC2PI);
00673 
00674     SINA  =sin(A);
00675 
00676     COSA  =cos(A);
00677 
00678     PERTL =PERTL +CCAMPM(1,K)*SINA;
00679 
00680     PERTLD=PERTLD+CCAMPM(2,K)*COSA;
00681 
00682     PERTP =PERTP +CCAMPM(3,K)*COSA;
00683 
00684     PERTPD=PERTPD-CCAMPM(4,K)*SINA;
00685   }
00686     
00687   TL=FORBEL(2)+PERTL;
00688 
00689   SINLM=sin(TL);
00690 
00691   COSLM=cos(TL);
00692 
00693   SIGMA=CCKM/(1.0+PERTP);
00694 
00695   A=SIGMA*(CCMLD+PERTLD);
00696 
00697   B=SIGMA*PERTPD;
00698 
00699   DXHD=DXHD+A*SINLM+B*COSLM;
00700 
00701   DYHD=DYHD-A*COSLM+B*SINLM;
00702 
00703   DZHD=    -SIGMA*CCFDI* cos(FORBEL(3));
00704 
00705   DXBD=DXHD*DC1MME;
00706 
00707   DYBD=DYHD*DC1MME;
00708 
00709   DZBD=DZHD*DC1MME;
00710 
00711   for (K = 1; K <= 4; K++) {
00712 
00713     PLON=FORBEL(K+3);
00714 
00715     POMG=SORBEL(K+1);
00716 
00717     PECC=SORBEL(K+9);
00718 
00719     TL=fmod(PLON+2.0*PECC* sin(PLON-POMG), CC2PI);
00720 
00721     SINLP(K)= sin(TL);
00722 
00723     COSLP(K)= cos(TL);
00724 
00725     DXBD=DXBD+CCPAMV(K)*(SINLP(K)+PECC*sin(POMG));
00726 
00727     DYBD=DYBD-CCPAMV(K)*(COSLP(K)+PECC*cos(POMG));
00728 
00729     DZBD=DZBD-CCPAMV(K)*SORBEL(K+13)*cos(PLON-SORBEL(K+5));
00730 
00731   }
00732     
00733   DCOSEP=cos(DEPS);
00734   DSINEP=sin(DEPS);
00735   DYAHD=DCOSEP*DYHD-DSINEP*DZHD;
00736   DZAHD=DSINEP*DYHD+DCOSEP*DZHD;
00737   DYABD=DCOSEP*DYBD-DSINEP*DZBD;
00738   DZABD=DSINEP*DYBD+DCOSEP*DZBD;
00739 
00740   DVELH[1]=DXHD;
00741   DVELH[2]=DYAHD;
00742   DVELH[3]=DZAHD;
00743 
00744   DVELB[1]=DXBD;
00745   DVELB[2]=DYABD;
00746   DVELB[3]=DZABD;
00747 
00748   for (N = 1; N <= 3; N++) {
00749     DVELH[N]=DVELH[N]*1.4959787E8;
00750     DVELB[N]=DVELB[N]*1.4959787E8;
00751   }
00752   return;
00753 }
00754 
00755 
00756 
00757 
00758 /*--------------------------------------------------------------------------*/
00780 /*--------------------------------------------------------------------------*/
00781 static void
00782 compxy(double inputr[19], char inputc[4],
00783        double outputr[4],
00784        double utr, double mod_juldat)
00785 {
00786   double STR;
00787   double t0, dl, theta0, pe, st0hg, stg;
00788   double jd, jd0h;
00789   double dvelb[4], dvelh[4];
00790   double alp, del, beov, berv, EDV;
00791   double HAR, phi, heov, herv;
00792   double *rbuf;
00793   char inpsgn[4];
00794   double *olong, *olat, *alpha, *delta;
00795   char signs[] = "+++";
00796   rbuf = inputr;
00797   inpsgn[1] = inputc[1];
00798   inpsgn[2] = inputc[2];
00799   inpsgn[3] = inputc[3];
00800   olong = rbuf + 7 - 1;
00801   olat  = rbuf + 10 - 1;
00802   alpha = rbuf + 13 - 1;
00803   delta = rbuf + 16 - 1;
00804   // ... convert UT to real hours, calculate Julian date
00805   /* We know this one already but convert seconds -> hours */
00806   utr /= 3600;
00807 
00808 
00809   jd = mod_juldat + 2400000.5;
00810   
00811   // ... likewise convert longitude and latitude of observatory to real hours
00812   // ... and degrees, respectively; take care of signs
00813   // ... NOTE: east longitude is assumed for input !!
00814 
00815   if (olong[1] < 0 || olong[2] < 0 ||
00816       olong[3] < 0 || inpsgn[1] == '-') {
00817     signs[1] = '-';
00818     olong[1] = fabs(olong[1]);
00819     olong[2] = fabs(olong[2]);
00820     olong[3] = fabs(olong[3]);
00821   }
00822   dl = olong[1]+olong[2]/60.  +olong[3]/3600.;
00823   if (signs[1]   == '-') dl = -dl;
00824   dl = -dl*24.  /360.;
00825 
00826   if (olat[1] < 0 || olat[2] < 0 ||
00827       olat[3] < 0 || inpsgn[2] == '-') {
00828     signs[2] = '-';
00829  
00830     olat[1] = fabs(olat[1]);
00831     olat[2] = fabs(olat[2]);
00832     olat[3] = fabs(olat[3]);
00833 
00834   }
00835 
00836   phi = olat[1]+olat[2]/60.  +olat[3]/3600.;
00837 
00838   if (signs[2]   == '-') phi = -phi;
00839 
00840   phi = phi*M_PI/180. ;
00841 
00842   // ... convert right ascension and declination to real radians
00843 
00844   alp = (alpha[1]*3600. +alpha[2]*60. +alpha[3])*M_PI/(12.  *3600.  );
00845 
00846   if (delta[1] < 0 || delta[2] < 0 ||
00847       delta[3] < 0 || inpsgn[3] == '-') {
00848 
00849     signs[3] = '-';
00850 
00851     delta[1] = fabs(delta[1]);
00852     delta[2] = fabs(delta[2]);
00853     delta[3] = fabs(delta[3]);
00854 
00855   }
00856 
00857   del = (delta[1]*3600.0  + delta[2]*60.   + delta[3])
00858     * M_PI/(3600. *180. );
00859 
00860 
00861 
00862   if (signs[3]   == '-') del = - del;
00863 
00864   // ... calculate earth's orbital velocity in rectangular coordinates X,Y,Z
00865   // ... for both heliocentric and barycentric frames (DVELH, DVELB)
00866   // ... Note that setting the second argument of BARVEL to zero as done below
00867   // ... means that the input coordinates will not be corrected for precession.
00868 
00869 
00870   barvel(jd, 0.0, dvelh, dvelb);
00871 
00872   // ... with the rectangular velocity components known, the respective projections
00873   // ... HEOV and BEOV on a given line of sight (ALP,DEL) can be determined:
00874 
00875   // ... REFERENCE: THE ASTRONOMICAL ALMANAC 1982 PAGE:B17
00876 
00877   beov =
00878     dvelb[1]*cos(alp)*cos(del)+
00879     dvelb[2]*sin(alp)*cos(del)+
00880     dvelb[3]*sin(del);
00881       
00882   heov =
00883     dvelh[1]*cos(alp)*cos(del)+
00884     dvelh[2]*sin(alp)*cos(del)+
00885     dvelh[3]*sin(del);
00886       
00887 
00888   // ... For determination also of the contribution due to the diurnal rotation of
00889   // ... the earth (EDV), the hour angle (HAR) is needed at which the observation
00890   // ... was made which requires conversion of UT to sidereal time (ST).
00891 
00892   // ... Therefore, first compute ST at 0 hours UT (ST0HG)
00893 
00894   // ... REFERENCE : MEEUS J.,1980,ASTRONOMICAL FORMULAE FOR CALCULATORS
00895 
00896 
00897   jd0h = jd - (utr/24.0);
00898       
00899   t0 = (jd0h-2415020.  )/36525. ;
00900       
00901 
00902   theta0 = 0.276919398  +100.0021359  *t0+0.000001075  *t0*t0 ;
00903 
00904   pe = (int) theta0;
00905 
00906   theta0 = theta0 - pe;
00907 
00908   st0hg = theta0*24. ;
00909 
00910   // ... now do the conversion UT -> ST (MEAN SIDEREAL TIME)
00911 
00912   // ... REFERENCE : THE ASTRONOMICAL ALMANAC 1983, P B7
00913   // ... IN 1983: 1 MEAN SOLAR DAY = 1.00273790931 MEAN SIDEREAL DAYS
00914   // ... ST WITHOUT EQUATION OF EQUINOXES CORRECTION => ACCURACY +/- 1 SEC
00915   //
00916   stg = st0hg+utr*1.00273790931 ;
00917       
00918   if (stg < dl) stg = stg +24. ;
00919 
00920   STR = stg-dl;
00921 
00922 
00923   if (STR >= 24. ) STR = STR-24. ;
00924 
00925   STR = STR*M_PI/12. ;
00926 
00927   HAR = STR-alp;
00928       
00929 
00930   EDV = -0.4654  * sin(HAR)* cos(del)* cos(phi);
00931 
00932   // ... the total correction (in km/s) is the sum of orbital and diurnal components
00933 
00934 
00935   herv=heov+EDV;
00936   berv=beov+EDV;
00937 
00938   /* The following is not needed. Do not translate */
00939 
00940 #if 0
00941   // ... Calculation of the barycentric and heliocentric correction times
00942   // ... (BCT and HCT) requires knowledge of the earth's position in its
00943   // ... orbit. Subroutine BARCOR returns the rectangular barycentric (DCORB)
00944   // ... and heliocentric (DCORH) coordinates.
00945 
00946   //      CALL BARCOR(DCORH,DCORB)
00947 
00948   // ... from this, the correction times (in days) can be determined:
00949   // ... (REFERENCE: THE ASTRONOMICAL ALMANAC 1982 PAGE:B16)
00950 
00951   //      BCT=+0.0057756D0*(DCORB(1)*DCOS(ALP)*DCOS(DEL)+
00952   //     1                DCORB(2)*DSIN(ALP)*DCOS(DEL)+
00953   //     2                DCORB(3)*          DSIN(DEL))
00954   //      HCT=+0.0057756D0*(DCORH(1)*DCOS(ALP)*DCOS(DEL)+
00955   //     1                DCORH(2)*DSIN(ALP)*DCOS(DEL)+
00956   //     2                DCORH(3)*          DSIN(DEL))
00957 
00958   //... write results to keywords
00959 
00960   //      CALL STKWRD('OUTPUTD',BCT,1,1,KUN,STAT)    ! barycentric correction time
00961   //      CALL STKWRD('OUTPUTD',HCT,2,1,KUN,STAT)    ! heliocentric correction time
00962 #endif
00963 
00964   rbuf[1] = berv;   /* barocentric RV correction */
00965   rbuf[2] = herv;   /* heliocentric RV correction */
00966   rbuf[3] = EDV;    /* diurnal RV correction */
00967 
00968 
00969   outputr[1] = rbuf[1];
00970   outputr[2] = rbuf[2];
00971   outputr[3] = rbuf[3];
00972 
00973   return;
00974 }
00975 
00976 
00977 
00978 /*----------------------------------------------------------------------------*/
00985 /*----------------------------------------------------------------------------*/
00986 cpl_error_code
00987 sinfo_baryvel(const cpl_propertylist *raw_header,
00988          double *bary_corr,
00989          double *helio_corr)
00990 {
00991 
00992     double outputr[4];
00993 
00994     char inputc[] = "X+++";       /* 0th index not used */
00995 
00996     double rneg = 1.0;
00997 
00998     double inputr[19];                  /* Do not use the zeroth element */
00999 
01000 
01001 /*
01002   qc_ra       = m$value({p1},O_POS(1))
01003   qc_dec      = m$value({p1},O_POS(2))
01004   qc_geolat   = m$value({p1},{h_geolat})
01005   qc_geolon   = m$value({p1},{h_geolon})
01006   qc_obs_time = m$value({p1},O_TIME(7))  !using an image as input it take the
01007                                          !date from the descriptor O_TIME(1,2,3)
01008                                          !and the UT from O_TIME(5)
01009 */
01010     double qc_ra;
01011     double qc_dec;
01012     double qc_geolat;
01013     double qc_geolon;
01014 
01015     double utr;
01016     double mod_juldat;
01017 
01018     double ra_hour, ra_min, ra_sec;
01019     double dec_deg, dec_min, dec_sec;
01020     double lat_deg, lat_min, lat_sec;
01021     double lon_deg, lon_min, lon_sec;
01022 
01023     check( qc_ra       = sinfo_pfits_get_ra(raw_header),  /* in degrees */
01024        "Error getting object right ascension");
01025     check( qc_dec      = sinfo_pfits_get_dec(raw_header),
01026        "Error getting object declination");
01027 
01028     check( qc_geolat   = sinfo_pfits_get_geolat(raw_header),
01029        "Error getting telescope latitude");
01030     check( qc_geolon   = sinfo_pfits_get_geolon(raw_header),
01031        "Error getting telescope longitude");
01032 
01033     /* double qc_obs_time = sinfo_pfits_get_exptime(raw_header);   Not used! */
01034 
01035     check( utr         = sinfo_pfits_get_utc(raw_header),
01036        "Error reading UTC");
01037     check( mod_juldat  = sinfo_pfits_get_mjdobs(raw_header),
01038        "Error julian date");
01039 
01040     deg2hms(qc_ra,     &ra_hour, &ra_min, &ra_sec);
01041     deg2dms(qc_dec,    &dec_deg, &dec_min, &dec_sec);
01042     deg2dms(qc_geolat, &lat_deg, &lat_min, &lat_sec);
01043     deg2dms(qc_geolon, &lon_deg, &lon_min, &lon_sec);
01044 
01045 
01046     inputr[7] = lon_deg;
01047     inputr[8] = lon_min;
01048     inputr[9] = lon_sec;
01049 
01050 
01051     rneg = (inputr[7]*3600.)+(inputr[8]*60.)+inputr[9];
01052 
01053     inputc[1] = (lon_deg >= 0) ? '+' : '-';
01054 
01055     if (rneg < 0) inputc[1] = '-';
01056 
01057 
01058     inputr[10] = lat_deg;
01059     inputr[11] = lat_min;
01060     inputr[12] = lat_sec;
01061 
01062 
01063     rneg = (inputr[10]*3600.)+(inputr[11]*60.)+inputr[12];
01064 
01065     inputc[2] = (lat_deg >= 0) ? '+' : '-';
01066 
01067     if (rneg < 0) inputc[2] = '-';
01068 
01069 
01070     inputr[13] = ra_hour;
01071     inputr[14] = ra_min;
01072     inputr[15] = ra_sec;
01073 
01074 
01075     inputr[16] = dec_deg;
01076     inputr[17] = dec_min;
01077     inputr[18] = dec_sec;
01078 
01079 
01080     inputc[3] = (dec_deg >= 0) ? '+' : '-';
01081 
01082     rneg = (inputr[16]*3600.)+(inputr[17]*60.)+inputr[18];
01083 
01084     if (rneg < 0) inputc[3] = '-';
01085     
01086 
01087 //C  INPUTR/R/1/3    date: year,month,day
01088 //C  INPUTR/R/4/3    universal time: hour,min,sec
01089 //C  INPUTR/R/7/3    EAST longitude of observatory: degree,min,sec  !! NOTE
01090 //C  INPUTR/R/10/3   latitude of observatory: degree,min,sec
01091 //C  INPUTR/R/13/3   right ascension: hour,min,sec
01092 //C  INPUTR/R/16/3   declination: degree,min,sec
01093 
01094     /* compute the corrections */
01095     compxy(inputr, inputc, outputr, utr, mod_juldat);
01096 
01097    sinfo_msg_debug("        Total barycentric RV correction:  %f km/s", outputr[1]);
01098    sinfo_msg_debug("        Total heliocentric RV correction: %f km/s", outputr[2]);
01099    sinfo_msg_debug("          (incl. diurnal RV correction of %f km/s)", outputr[3]);
01100 
01101 
01102    *bary_corr = outputr[1];
01103    *helio_corr = outputr[2];
01104 
01105   cleanup:
01106     if (cpl_error_get_code() != CPL_ERROR_NONE) {
01107        sinfo_check_rec_status(0);
01108     }
01109     return cpl_error_get_code();
01110 }

Generated on 3 Mar 2013 for SINFONI Pipeline Reference Manual by  doxygen 1.6.1