GIRAFFE Pipeline Reference Manual

girvcorrection.c

00001 /* $Id: girvcorrection.c,v 1.2 2008/03/17 13:56:17 rpalsa Exp $
00002  *
00003  * This file is part of the GIRAFFE Pipeline
00004  * Copyright (C) 2002-2006 European Southern Observatory
00005  *
00006  * This program is free software; you can redistribute it and/or modify
00007  * it under the terms of the GNU General Public License as published by
00008  * the Free Software Foundation; either version 2 of the License, or
00009  * (at your option) any later version.
00010  *
00011  * This program is distributed in the hope that it will be useful,
00012  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00013  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014  * GNU General Public License for more details.
00015  *
00016  * You should have received a copy of the GNU General Public License
00017  * along with this program; if not, write to the Free Software
00018  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00019  */
00020 
00021 /*
00022  * $Author: rpalsa $
00023  * $Date: 2008/03/17 13:56:17 $
00024  * $Revision: 1.2 $
00025  * $Name: giraffe-2_9 $
00026  */
00027 
00028 #ifdef HAVE_CONFIG_H
00029 #  include <config.h>
00030 #endif
00031 
00032 
00033 #include <math.h>
00034 #include <string.h>
00035 
00036 #include <cxtypes.h>
00037 
00038 #include <cpl_matrix.h>
00039 
00040 #include "girvcorrection.h"
00041 
00042 
00043 
00052 /*
00053  * Constants
00054  */
00055 
00056 static const cxdouble dct0   = 2415020.0;
00057 static const cxdouble dcjul  = 36525.0;
00058 static const cxdouble dc1900 = 1900.0;
00059 static const cxdouble dctrop = 365.24219572;
00060 static const cxdouble dcbes  = 0.313;
00061 
00062 static const cxdouble RV_DPI  =
00063     3.1415926535897932384626433832795028841971693993751;   /* pi      */
00064 
00065 static const cxdouble RV_D2PI =
00066     6.2831853071795864769252867665590057683943387987502;   /* 2pi     */
00067 
00068 static const cxdouble RV_D4PI =
00069     12.566370614359172953850573533118011536788677597500;   /* 4pi     */
00070 
00071 static const cxdouble RV_DPIBY2 =
00072     1.5707963267948966192313216916397514420985846996876;   /* pi/2   */
00073 
00074 static const cxdouble RV_DD2R =
00075     0.017453292519943295769236907684886127134428718885417; /* pi/180  */
00076 
00077 static const cxdouble RV_DAS2R =
00078     4.8481368110953599358991410235794797595635330237270e-6; /* pi/(180*3600) */
00079 
00080 
00081 /*
00082  * Slalib functions
00083  */
00084 
00085 /*
00086  *  - - - - - - - - - -
00087  *   s l a D e u l e r
00088  *  - - - - - - - - - -
00089  *
00090  *  Form a rotation matrix from the Euler angles - three successive
00091  *  rotations about specified Cartesian axes.
00092  *
00093  *  (double precision)
00094  *
00095  *  Given:
00096  *    *order char     specifies about which axes the rotations occur
00097  *    phi    double   1st rotation (radians)
00098  *    theta  double   2nd rotation (   "   )
00099  *    psi    double   3rd rotation (   "   )
00100  *
00101  *  Returned:
00102  *    Matrix   rmat  rotation matrix
00103  *
00104  *  A rotation is positive when the reference frame rotates
00105  *  anticlockwise as seen looking towards the origin from the
00106  *  positive region of the specified axis.
00107  *
00108  *  The characters of order define which axes the three successive
00109  *  rotations are about.  A typical value is 'zxz', indicating that
00110  *  rmat is to become the direction cosine matrix corresponding to
00111  *  rotations of the reference frame through phi radians about the
00112  *  old z-axis, followed by theta radians about the resulting x-axis,
00113  *  then psi radians about the resulting z-axis.
00114  *
00115  *  The axis names can be any of the following, in any order or
00116  *  combination:  x, y, z, uppercase or lowercase, 1, 2, 3.  Normal
00117  *  axis labelling/numbering conventions apply;  the xyz (=123)
00118  *  triad is right-handed.  Thus, the 'zxz' example given above
00119  *  could be written 'zxz' or '313' (or even 'zxz' or '3xz').  Order
00120  *  is terminated by length or by the first unrecognized character.
00121  *
00122  *  Fewer than three rotations are acceptable, in which case the later
00123  *  angle arguments are ignored.  Zero rotations leaves rmat set to the
00124  *  identity matrix.
00125  *
00126  *  Last revision:   9 December 1996
00127  *
00128  *  Copyright P.T.Wallace.  All rights reserved.
00129  */
00130 
00131 inline static void
00132 slaDeuler(const cxchar* order, cxdouble phi, cxdouble theta, cxdouble psi,
00133           cpl_matrix* rmat)
00134 {
00135     register cxint j, i, l, n, k;
00136     cxdouble result[3][3], rotn[3][3], angle, s, c , w, wm[3][3];
00137     cxchar axis;
00138 
00139     /* Initialize result matrix */
00140     for ( j = 0; j < 3; j++ ) {
00141         for ( i = 0; i < 3; i++ ) {
00142             result[i][j] = ( i == j ) ? 1.0 : 0.0;
00143         }
00144     }
00145 
00146     /* Establish length of axis string */
00147     l = strlen ( order );
00148 
00149     /* Look at each character of axis string until finished */
00150     for ( n = 0; n < 3; n++ ) {
00151         if ( n <= l ) {
00152 
00153             /* Initialize rotation matrix for the current rotation */
00154             for ( j = 0; j < 3; j++ ) {
00155                 for ( i = 0; i < 3; i++ ) {
00156                     rotn[i][j] = ( i == j ) ? 1.0 : 0.0;
00157                 }
00158             }
00159 
00160             /* Pick up the appropriate Euler angle and take sine & cosine */
00161             switch ( n ) {
00162             case 0 :
00163                 angle = phi;
00164                 break;
00165             case 1 :
00166                 angle = theta;
00167                 break;
00168             default:
00169                 angle = psi;
00170                 break;
00171             }
00172             s = sin ( angle );
00173             c = cos ( angle );
00174 
00175             /* Identify the axis */
00176             axis =  order[n];
00177             if ( ( axis == 'X' ) || ( axis == 'x' ) || ( axis == '1' ) ) {
00178 
00179                 /* Matrix for x-rotation */
00180                 rotn[1][1] = c;
00181                 rotn[1][2] = s;
00182                 rotn[2][1] = -s;
00183                 rotn[2][2] = c;
00184             }
00185             else if ( ( axis == 'Y' ) || ( axis == 'y' ) || ( axis == '2' ) ) {
00186 
00187                 /* Matrix for y-rotation */
00188                 rotn[0][0] = c;
00189                 rotn[0][2] = -s;
00190                 rotn[2][0] = s;
00191                 rotn[2][2] = c;
00192             }
00193             else if ( ( axis == 'Z' ) || ( axis == 'z' ) || ( axis == '3' ) ) {
00194 
00195                 /* Matrix for z-rotation */
00196                 rotn[0][0] = c;
00197                 rotn[0][1] = s;
00198                 rotn[1][0] = -s;
00199                 rotn[1][1] = c;
00200             }
00201             else {
00202 
00203                 /* Unrecognized character - fake end of string */
00204                 l = 0;
00205             }
00206 
00207             /* Apply the current rotation (matrix rotn x matrix result) */
00208             for ( i = 0; i < 3; i++ ) {
00209                 for ( j = 0; j < 3; j++ ) {
00210                     w = 0.0;
00211                     for ( k = 0; k < 3; k++ ) {
00212                         w += rotn[i][k] * result[k][j];
00213                     }
00214                     wm[i][j] = w;
00215                 }
00216             }
00217             for ( j = 0; j < 3; j++ ) {
00218                 for ( i= 0; i < 3; i++ ) {
00219                     result[i][j] = wm[i][j];
00220                 }
00221             }
00222         }
00223     }
00224 
00225     /* Copy the result */
00226     for ( j = 0; j < 3; j++ ) {
00227         for ( i = 0; i < 3; i++ ) {
00228             cpl_matrix_set(rmat, i, j, result[i][j]);
00229         }
00230     }
00231 
00232     return;
00233 
00234 }
00235 
00236 
00237 /*
00238  * @brief
00239  *  Calculates the matrix of general precession from ep0 to ep1
00240  *
00241  * @param ep0  Initial epoch of mean equator and equinox e.g., 1950.0
00242  * @param ep1  Final epoch of mean equator and equinox
00243  *
00244  * @return
00245  *   The 3 x 3 matrix of general precession.
00246  *
00247  * Form the matrix of precession between two epochs, using the
00248  * model of Simon et al (1994), which is suitable for long
00249  * periods of time.
00250  *
00251  * @notes
00252  *   -# The epochs are TDB (loosely ET) Julian epochs.
00253  *   -# The matrix is in the sense   v(ep1)  =  rmatp * v(ep0)
00254  *   -# The absolute accuracy of the model is limited by the
00255  *      uncertainty in the general precession, about 0.3 arcsec per
00256  *      1000 years.  The remainder of the formulation provides a
00257  *      precision of 1 mas over the interval from 1000AD to 3000AD,
00258  *      0.1 arcsec from 1000BC to 5000AD and 1 arcsec from
00259  *      4000BC to 8000AD.
00260  *
00261  * Reference:
00262  *  Simon, J. L., et al., 1994. Astron.Astrophys., 282, 663-683.
00263  *
00264  * Author:
00265  *  P. T. Wallace. 23 August 1994
00266  */
00267 
00268 inline static cpl_matrix*
00269 slaPrecession (cxdouble ep0, cxdouble ep1)
00270 {
00271 
00272     cxdouble t = 0.;
00273     cxdouble w = 0.;
00274     cxdouble z = 0.;
00275     cxdouble t0 = 0.;
00276     cxdouble zeta = 0.;
00277     cxdouble theta = 0.;
00278     cxdouble tas2r = 0.;
00279 
00280     cpl_matrix* mprec = NULL;
00281 
00282 
00283    /*
00284     * Interval between basic epoch J2000.0 and beginning epoch (1000JY)
00285     */
00286 
00287     t0 = ( ep0 - 2000.0 ) / 1000.0;
00288 
00289 
00290    /*
00291     * Interval over which precession required (1000JY)
00292     */
00293 
00294     t =  ( ep1 - ep0 ) / 1000.0;
00295 
00296 
00297    /*
00298     * Euler angles
00299     */
00300 
00301     tas2r = t * RV_DAS2R;
00302     w = 23060.9097 +
00303         (139.7459 +
00304          (-0.0038 +
00305           (-0.5918 +
00306            (-0.0037 +
00307             0.0007 * t0) * t0) * t0) * t0) * t0;
00308 
00309     zeta = (w +
00310             (30.2226 +
00311              (-0.2523 +
00312               (-0.3840 +
00313                (-0.0014 +
00314                 0.0007 * t0) * t0) * t0) * t0 +
00315              (18.0183 +
00316               (-0.1326 +
00317                (0.0006 +
00318                 0.0005 * t0) * t0) * t0 +
00319               (-0.0583 +
00320                (-0.0001 +
00321                 0.0007 * t0) * t0 +
00322                (-0.0285 +
00323                 -0.0002 * t) * t) * t) * t) * t) * tas2r;
00324 
00325     z = (w +
00326          (109.5270 +
00327           (0.2446 +
00328            (-1.3913 +
00329             (-0.0134 +
00330              0.0026 * t0) * t0) * t0) * t0 +
00331           (18.2667 +
00332            (-1.1400 +
00333             (-0.0173 +
00334              0.0044 * t0) * t0) * t0 +
00335            (-0.2821 +
00336             (-0.0093 +
00337              0.0032 * t0) * t0 +
00338             (-0.0301 +
00339              0.0006 * t0 +
00340              -0.0001 * t) * t) * t) * t) * t) * tas2r;
00341 
00342     theta = (20042.0207 +
00343              (-85.3131 +
00344               (-0.2111 +
00345                (0.3642 +
00346                 (0.0008 +
00347                  -0.0005 * t0) * t0) * t0) * t0) * t0 +
00348              (-42.6566 +
00349               (-0.2111 +
00350                (0.5463 +
00351                 (0.0017 +
00352                  -0.0012 * t0) * t0) * t0) * t0 +
00353               (-41.8238 +
00354                (0.0359 +
00355                 (0.0027 +
00356                  -0.0001 * t0) * t0) * t0 +
00357                (-0.0731 +
00358                 (0.0019 +
00359                  0.0009 * t0) * t0 +
00360                 (-0.0127 +
00361                  0.0011 * t0 + 0.0004 * t) * t) * t) * t) * t) * tas2r;
00362 
00363 
00364    /*
00365     * Rotation matrix
00366     */
00367 
00368     mprec = cpl_matrix_new(3, 3);
00369     slaDeuler("ZYZ", -zeta, theta, -z, mprec);
00370 
00371     return mprec;
00372 
00373 }
00374 
00375 
00376 
00377 /*
00378  * @brief
00379  *  Compute the mean local sideral time
00380  *
00381  * @param djd    julian date
00382  * @param dlong  observer's longitude (radians)
00383  *
00384  * @return
00385  *   Mean local sideral time (radians).
00386  *
00387  * @note
00388  *  Constants taken from the american ephemeris and nautical almanac, 1980)
00389  *  Translated from the FORTRAN version written by G. Torres (1989)
00390  */
00391 
00392 inline static cxdouble
00393 sideral_time(cxdouble djd, cxdouble dlong)
00394 {
00395 
00396     /*
00397      * Constants D1,D2,D3 for calculating Greenwich
00398      * Mean Sideral Time at 0 UT
00399      */
00400 
00401     const cxdouble d1 = 1.739935934667999;
00402     const cxdouble d2 = 6.283319509909095e02;
00403     const cxdouble d3 = 6.755878646261384e-06;
00404 
00405     const cxdouble df = 1.00273790934;
00406 
00407     cxdouble dut  = 0.;
00408     cxdouble dt   = 0.;
00409     cxdouble dst0 = 0.;
00410     cxdouble dst  = 0.;
00411     cxdouble djd0 = floor(djd) + 0.5;
00412 
00413     if (djd0 > djd) {
00414         djd0 -= 1.;
00415     }
00416 
00417     dut = (djd - djd0) * RV_D2PI;
00418 
00419     dt = (djd0 - dct0) / dcjul;
00420     dst0 = d1 + d2 * dt + d3 * dt * dt;
00421     dst0 = fmod(dst0, RV_D2PI);
00422     dst = df * dut + dst0 - dlong;
00423     dst = fmod(dst + RV_D4PI, RV_D2PI);
00424 
00425     return dst;
00426 
00427 }
00428 
00429 
00430 /*
00431  * @brief
00432  *  Compute correction from topocentric radial velocity to geocentric.
00433  *
00434  * @param dlat  observer's geodetic latitude (radians)
00435  * @param dalt  observer's elevation above sea level (meters)
00436  * @param dec   star's declination (radians) for mean
00437  *              equator and equinox of date
00438  * @param dha   star's hour angle (radians)
00439  *
00440  * @return
00441  *  Geocentric correction in km/s.
00442  *
00443  * Calculates the correction required to transform the topocentric radial
00444  * velocity of a given star to geocentric. The maximum error of this
00445  * routine is not expected to be larger than 0.1 m/s.
00446  *
00447  * @note
00448  *  vr = r.w.cos(dec).sin(hour angle), where r = geocentric radius at
00449  *  observer's latitude and elevation, and w = 2.pi/t, t = length of sideral
00450  *  day 23 hours 56 minutes and 4 seconds in seconds.
00451  *  The hour angle is positive east of the meridian.
00452  *  Other relevant equations from E. W. Woolard & G. M. Clemence (1966),
00453  *  Spherical Astronomy,  p.45 and p.48
00454  *
00455  * Author:
00456  *  - G. Torres (1989)
00457  *  - G. Simond (C translation)
00458  */
00459 
00460 inline static cxdouble
00461 geo_correction(cxdouble dlat, cxdouble dalt, cxdouble dec, cxdouble dha)
00462 {
00463 
00464     /*
00465      * Earth's equatorial radius (km)  (Geod. ref. sys., 1980)
00466      */
00467 
00468     const cxdouble da = 6378.137;
00469 
00470     /*
00471      * Polar flattening (Geod. ref. sys., 1980)
00472      */
00473 
00474     const cxdouble df =  1./298.257222;
00475 
00476     /*
00477      * Angular rotation rate (2.pi/t)
00478      */
00479 
00480     const cxdouble dw = RV_D2PI/86164.;
00481 
00482 
00483     const cxdouble de2 = df * (2.0 - df);
00484     const cxdouble dsdlats = sin (dlat) * sin (dlat);
00485 
00486     cxdouble d1    = 0.;
00487     cxdouble d2    = 0.;
00488     cxdouble dr0   = 0.;
00489     cxdouble dlatg = 0.;
00490     cxdouble drh   = 0.;
00491     cxdouble dvelg = 0.;
00492 
00493 
00494     /*
00495      * Calculate geocentric radius dr0 at sea level (km)
00496      */
00497 
00498     d1 = 1.0 - de2 * (2.0 - de2) * dsdlats;
00499     d2 = 1.0 - de2 * dsdlats;
00500     dr0 = da * sqrt(d1 / d2);
00501 
00502 
00503     /*
00504      * Calculate geocentric latitude dlatg
00505      */
00506 
00507     d1 = de2 * sin(2.0 * dlat);
00508     d2 = 2.0 * d2;
00509     dlatg = dlat - atan(d1 / d2);
00510 
00511 
00512     /*
00513      * Calculate geocentric radius drh at elevation dalt (km)
00514      */
00515 
00516     drh = dr0 * cos(dlatg) + (dalt / 1000.) * cos(dlat);
00517 
00518 
00519     /*
00520      * Projected component to star at declination = dec and
00521      * at hour angle = dha, in units of km/s
00522      */
00523 
00524     dvelg = dw * drh * cos(dec) * sin(dha);
00525 
00526     return dvelg;
00527 
00528 }
00529 
00530 
00531 /*
00532  * @brief
00533  *  Compute the heliocentric and barycentric velocity components of the earth.
00534  *
00535  *  @param dje  Julian ephermeris date
00536  *  @param deq  Epoch of mean equator and mean equinox of @em hvel
00537  *              and @em bvel. If @em deq is @c 0, both vectors refer
00538  *              to the mean equator and equinox of @em dje
00539  *
00540  * @return
00541  *  The components (dx/dt, dy/dt, dz/dt, k=1,2,3  A.U./s) of the
00542  *  heliocentric and barycentric velocity are returned in the
00543  *  vectors @em hvel and @em bvel respectively.
00544  *
00545  * Calculates the heliocentric and barycentric velocity components
00546  * of the earth.  The largest deviations from the JPL-de96 ephemeris
00547  * are 42 cm/s for both heliocentric and barycentric velocity
00548  * components.
00549  *
00550  * @note
00551  *  vr = r.w.cos(dec).sin(hour angle), where r = geocentric radius at
00552  *  observer's latitude and elevation, and w = 2.pi/t, t = length of sideral
00553  *  day (sec). The hour angle is positive east of the meridian.
00554  *  Other relevant equations from E. W. Woolard  & G. M. Clemence (1966),
00555  *  Spherical Astronomy,  p.45 and p.48
00556  *
00557  * Authors:
00558  *  - P. Stumpff (ibm-version 1979): astron. astrophys. suppl. ser. 41,
00559  *    1 (1980)
00560  *  - M. H.Slovak (vax 11/780 implementation 1986)
00561  *  - G. Torres (1989)
00562  *  - G. Simond (C translation)
00563  */
00564 
00565 inline static void
00566 earth_velocity(cxdouble dje, cxdouble deq, cxdouble* const hvel,
00567                cxdouble* const bvel)
00568 {
00569 
00570 
00571     /*
00572      * Constants dcfel(i, k) of fast-changing elements
00573      *
00574      *   i = 1          i = 2                 i = 3
00575      */
00576 
00577     const cxdouble dcfel[][3] = {
00578         {1.7400353e+00, 6.2833195099091e+02,  5.2796e-06},
00579         {6.2565836e+00, 6.2830194572674e+02, -2.6180e-06},
00580         {4.7199666e+00, 8.3997091449254e+03, -1.9780e-05},
00581         {1.9636505e-01, 8.4334662911720e+03, -5.6044e-05},
00582         {4.1547339e+00, 5.2993466764997e+01,  5.8845e-06},
00583         {4.6524223e+00, 2.1354275911213e+01,  5.6797e-06},
00584         {4.2620486e+00, 7.5025342197656e+00,  5.5317e-06},
00585         {1.4740694e+00, 3.8377331909193e+00,  5.6093e-06}
00586     };
00587 
00588 
00589     /*
00590      * Constants dceps and ccsel(i,k) of slowly changing elements
00591      *
00592      *  i = 1          i = 2          i = 3
00593      */
00594 
00595     const cxdouble dceps[3] = {
00596          4.093198e-01,
00597         -2.271110e-04,
00598         -2.860401e-08
00599     };
00600 
00601     const cxdouble ccsel[][3] = {
00602         {1.675104e-02, -4.179579e-05, -1.260516e-07},
00603         {2.220221e-01,  2.809917e-02,  1.852532e-05},
00604         {1.589963e+00,  3.418075e-02,  1.430200e-05},
00605         {2.994089e+00,  2.590824e-02,  4.155840e-06},
00606         {8.155457e-01,  2.486352e-02,  6.836840e-06},
00607         {1.735614e+00,  1.763719e-02,  6.370440e-06},
00608         {1.968564e+00,  1.524020e-02, -2.517152e-06},
00609         {1.282417e+00,  8.703393e-03,  2.289292e-05},
00610         {2.280820e+00,  1.918010e-02,  4.484520e-06},
00611         {4.833473e-02,  1.641773e-04, -4.654200e-07},
00612         {5.589232e-02, -3.455092e-04, -7.388560e-07},
00613         {4.634443e-02, -2.658234e-05,  7.757000e-08},
00614         {8.997041e-03,  6.329728e-06, -1.939256e-09},
00615         {2.284178e-02, -9.941590e-05,  6.787400e-08},
00616         {4.350267e-02, -6.839749e-05, -2.714956e-07},
00617         {1.348204e-02,  1.091504e-05,  6.903760e-07},
00618         {3.106570e-02, -1.665665e-04, -1.590188e-07}
00619     };
00620 
00621 
00622     /*
00623      * Constants of the arguments of the short-period perturbations by
00624      * the planets: dcargs(i, k)
00625      *
00626      *   i = 1        i = 2
00627      */
00628 
00629     const cxdouble dcargs[][2] = {
00630         {5.0974222e+00, -7.8604195454652e+02},
00631         {3.9584962e+00, -5.7533848094674e+02},
00632         {1.6338070e+00, -1.1506769618935e+03},
00633         {2.5487111e+00, -3.9302097727326e+02},
00634         {4.9255514e+00, -5.8849265665348e+02},
00635         {1.3363463e+00, -5.5076098609303e+02},
00636         {1.6072053e+00, -5.2237501616674e+02},
00637         {1.3629480e+00, -1.1790629318198e+03},
00638         {5.5657014e+00, -1.0977134971135e+03},
00639         {5.0708205e+00, -1.5774000881978e+02},
00640         {3.9318944e+00,  5.2963464780000e+01},
00641         {4.8989497e+00,  3.9809289073258e+01},
00642         {1.3097446e+00,  7.7540959633708e+01},
00643         {3.5147141e+00,  7.9618578146517e+01},
00644         {3.5413158e+00, -5.4868336758022e+02}
00645     };
00646 
00647 
00648     /*
00649      * Amplitudes ccamps(n, k) of the short-period perturbations
00650      *
00651      *   n = 1          n = 2         n = 3         n = 4        n = 5
00652      */
00653 
00654     const cxdouble ccamps[][5] = {
00655         {-2.279594e-5,  1.407414e-5,  8.273188e-6,  1.340565e-5, -2.490817e-7},
00656         {-3.494537e-5,  2.860401e-7,  1.289448e-7,  1.627237e-5, -1.823138e-7},
00657         { 6.593466e-7,  1.322572e-5,  9.258695e-6, -4.674248e-7, -3.646275e-7},
00658         { 1.140767e-5, -2.049792e-5, -4.747930e-6, -2.638763e-6, -1.245408e-7},
00659         { 9.516893e-6, -2.748894e-6, -1.319381e-6, -4.549908e-6, -1.864821e-7},
00660         { 7.310990e-6, -1.924710e-6, -8.772849e-7, -3.334143e-6, -1.745256e-7},
00661         {-2.603449e-6,  7.359472e-6,  3.168357e-6,  1.119056e-6, -1.655307e-7},
00662         {-3.228859e-6,  1.308997e-7,  1.013137e-7,  2.403899e-6, -3.736225e-7},
00663         { 3.442177e-7,  2.671323e-6,  1.832858e-6, -2.394688e-7, -3.478444e-7},
00664         { 8.702406e-6, -8.421214e-6, -1.372341e-6, -1.455234e-6, -4.998479e-8},
00665         {-1.488378e-6, -1.251789e-5,  5.226868e-7, -2.049301e-7,  0.0e0},
00666         {-8.043059e-6, -2.991300e-6,  1.473654e-7, -3.154542e-7,  0.0e0},
00667         { 3.699128e-6, -3.316126e-6,  2.901257e-7,  3.407826e-7,  0.0e0},
00668         { 2.550120e-6, -1.241123e-6,  9.901116e-8,  2.210482e-7,  0.0e0},
00669         {-6.351059e-7,  2.341650e-6,  1.061492e-6,  2.878231e-7,  0.0e0}
00670     };
00671 
00672 
00673     /*
00674      * Constants of the secular perturbations in longitude
00675      * ccsec3 and ccsec(n,k)
00676      *   n = 1         n = 2         n = 3
00677      */
00678 
00679     const cxdouble ccsec3 = -7.757020e-08;
00680 
00681     const cxdouble ccsec[][3] = {
00682         {1.289600e-06, 5.550147e-01, 2.076942e+00},
00683         {3.102810e-05, 4.035027e+00, 3.525565e-01},
00684         {9.124190e-06, 9.990265e-01, 2.622706e+00},
00685         {9.793240e-07, 5.508259e+00, 1.559103e+01}
00686     };
00687 
00688 
00689     /*
00690      * Sideral rate dcsld in longitude, rate ccsgd in mean anomaly
00691      */
00692 
00693     const cxdouble dcsld = 1.990987e-07;
00694     const cxdouble ccsgd = 1.990969e-07;
00695 
00696 
00697     /*
00698      * Some constants used in the calculation of the
00699      * lunar contribution
00700      */
00701 
00702     const cxdouble cckm  = 3.122140e-05;
00703     const cxdouble ccmld = 2.661699e-06;
00704     const cxdouble ccfdi = 2.399485e-07;
00705 
00706 
00707     /*
00708      * Constants dcargm(i,k) of the arguments of the
00709      * perturbations of the motion of the moon
00710      *
00711      *   i = 1          i = 2
00712      */
00713 
00714     const cxdouble dcargm[][2] = {
00715         {5.1679830e+00,  8.3286911095275e+03},
00716         {5.4913150e+00, -7.2140632838100e+03},
00717         {5.9598530e+00,  1.5542754389685e+04}
00718     };
00719 
00720 
00721     /*
00722      * Amplitudes ccampm(n,k) of the perturbations of the moon
00723      *    n = 1         n = 2         n = 3         n = 4
00724      */
00725 
00726     const cxdouble ccampm[][4] = {
00727         { 1.097594e-01, 2.896773e-07, 5.450474e-02,  1.438491e-07},
00728         {-2.223581e-02, 5.083103e-08, 1.002548e-02, -2.291823e-08},
00729         { 1.148966e-02, 5.658888e-08, 8.249439e-03,  4.063015e-08}
00730     };
00731 
00732 
00733     /*
00734      * ccpamv = a * m * dl/dt (planets)
00735      */
00736 
00737     const cxdouble ccpamv[4] = {
00738         8.326827e-11,
00739         1.843484e-11,
00740         1.988712e-12,
00741         1.881276e-12
00742     };
00743 
00744 
00745     /*
00746      * dc1mme = 1 - mass(earth + moon)
00747      */
00748 
00749     const cxdouble dc1mme = 0.99999696;
00750 
00751 
00752     register cxint k = 0;
00753     register cxint n = 0;
00754 
00755     cxint ideq = 0;
00756 
00757     cxdouble a      = 0.;
00758     cxdouble b      = 0.;
00759     cxdouble f      = 0.;
00760     cxdouble dt     = 0.;
00761     cxdouble t      = 0.;
00762     cxdouble tl     = 0.;
00763     cxdouble dtsq   = 0.;
00764     cxdouble tsq    = 0.;
00765     cxdouble dml    = 0.;
00766     cxdouble dlocal = 0.;
00767     cxdouble deps   = 0.;
00768     cxdouble pertl  = 0.;
00769     cxdouble pertld = 0.;
00770     cxdouble pertr  = 0.;
00771     cxdouble pertrd = 0.;
00772     cxdouble pertp  = 0.;
00773     cxdouble pertpd = 0.;
00774     cxdouble sina   = 0.;
00775     cxdouble cosa   = 0.;
00776     cxdouble twoe   = 0.;
00777     cxdouble twog   = 0.;
00778     cxdouble param  = 0.;
00779     cxdouble dparam = 0.;
00780     cxdouble dpsi   = 0.;
00781     cxdouble phi    = 0.;
00782     cxdouble phid   = 0.;
00783     cxdouble psid   = 0.;
00784     cxdouble sin_f  = 0.;
00785     cxdouble cos_f  = 0.;
00786     cxdouble esq    = 0.;
00787     cxdouble d1pdro = 0.;
00788     cxdouble drd    = 0.;
00789     cxdouble drld   = 0.;
00790     cxdouble dsinls = 0.;
00791     cxdouble dcosls = 0.;
00792     cxdouble dtl    = 0.;
00793     cxdouble dxhd   = 0.;
00794     cxdouble dyhd   = 0.;
00795     cxdouble dzhd   = 0.;
00796     cxdouble sinlm  = 0.;
00797     cxdouble coslm  = 0.;
00798     cxdouble sigma  = 0.;
00799     cxdouble plon   = 0.;
00800     cxdouble pomg   = 0.;
00801     cxdouble dxbd   = 0.;
00802     cxdouble dybd   = 0.;
00803     cxdouble dzbd   = 0.;
00804     cxdouble pecc   = 0.;
00805     cxdouble dcosep = 0.;
00806     cxdouble dsinep = 0.;
00807     cxdouble dyahd  = 0.;
00808     cxdouble dzahd  = 0.;
00809     cxdouble dyabd  = 0.;
00810     cxdouble dzabd  = 0.;
00811     cxdouble sn[4]      = {0., 0., 0., 0.};
00812     cxdouble sinlp[4]   = {0., 0., 0., 0.};
00813     cxdouble coslp[4]   = {0., 0., 0., 0.};
00814     cxdouble forbel[7]  = {0., 0., 0., 0., 0., 0., 0.};
00815     cxdouble sorbel[17];
00816 
00817 
00818     memset(sorbel, 0, sizeof sorbel);
00819 
00820 
00821     /*
00822      * Control-parameter ideq, and time-arguments
00823      */
00824 
00825     ideq = (cxint)deq;
00826     dt = (dje - dct0) / dcjul;
00827     t = dt;
00828     dtsq = dt * dt;
00829     tsq = dtsq;
00830 
00831 
00832     /*
00833      * Values of all elements for the instant dje
00834      */
00835 
00836     for (k = 0; k < 8; k++) {
00837 
00838         dlocal = fmod(dcfel[k][0] + dt * dcfel[k][1] + dtsq * dcfel[k][2],
00839                          RV_D2PI);
00840 
00841         if (k == 0) {
00842             dml = dlocal;
00843         }
00844 
00845         if (k != 0) {
00846             forbel[k - 1] = dlocal;
00847         }
00848 
00849     }
00850 
00851     deps = fmod(dceps[0] + dt * dceps[1] + dtsq * dceps[2], RV_D2PI);
00852 
00853     for (k = 0; k < 17; k++) {
00854 
00855         sorbel[k] = fmod(ccsel[k][0] + t * ccsel[k][1] + tsq * ccsel[k][2],
00856                          RV_D2PI);
00857 
00858     }
00859 
00860 
00861     /*
00862      * Secular perturbations in longitude
00863      */
00864 
00865     for (k = 0; k < 4; k++) {
00866 
00867         a = fmod(ccsec[k][1] + t * ccsec[k][2], RV_D2PI);
00868         sn[k] = sin(a);
00869 
00870     }
00871 
00872 
00873     /*
00874      * Periodic perturbations of the earth-moon barycenter
00875      */
00876 
00877     pertl = ccsec[0][0] * sn[0] + ccsec[1][0] * sn[1] +
00878         (ccsec[2][0] + t * ccsec3) * sn[2] + ccsec[3][0] * sn[3];
00879 
00880     pertld = 0.;
00881     pertr  = 0.;
00882     pertrd = 0.;
00883 
00884     for (k = 0; k < 15; k++) {
00885 
00886         a = fmod(dcargs[k][0] + dt * dcargs[k][1], RV_D2PI);
00887         cosa = cos (a);
00888         sina = sin (a);
00889         pertl += (ccamps[k][0] * cosa + ccamps[k][1] * sina);
00890         pertr += (ccamps[k][2] * cosa + ccamps[k][3] * sina);
00891 
00892         if (k >= 10) {
00893             continue;
00894         }
00895 
00896         pertld += ((ccamps[k][1] * cosa - ccamps[k][0] * sina) * ccamps[k][4]);
00897         pertrd += ((ccamps[k][3] * cosa - ccamps[k][2] * sina) * ccamps[k][4]);
00898 
00899     }
00900 
00901 
00902     /*
00903      * Elliptic part of the motion of the earth-moon barycenter
00904      */
00905 
00906     esq = sorbel[0] * sorbel[0];
00907     dparam = 1. - esq;
00908     param = dparam;
00909     twoe = sorbel[0] + sorbel[0];
00910     twog = forbel[0] + forbel[0];
00911     phi = twoe * ((1.0 - esq * (1.0 / 8.0)) * sin (forbel[0]) +
00912                   sorbel[0] * (5.0 / 8.0) * sin (twog) +
00913                   esq * 0.5416667 * sin (forbel[0] + twog));
00914     f = forbel[0] + phi;
00915     sin_f = sin(f);
00916     cos_f = cos(f);
00917     dpsi = dparam / (1. + sorbel[0] * cos_f);
00918     phid = twoe * ccsgd * ((1.0 + esq * 1.50) * cos_f +
00919                            sorbel[0] * (1.250 - sin_f * sin_f * 0.50));
00920     psid = ccsgd * sorbel[0] * sin_f / sqrt(param);
00921 
00922 
00923     /*
00924      * Perturbed heliocentric motion of the earth-moon barycenter
00925      */
00926 
00927     d1pdro = 1. + pertr;
00928     drd = d1pdro * (psid + dpsi * pertrd);
00929     drld = d1pdro * dpsi * (dcsld + phid + pertld);
00930     dtl = fmod(dml + phi + pertl, RV_D2PI);
00931     dsinls = sin(dtl);
00932     dcosls = cos(dtl);
00933     dxhd = drd * dcosls - drld * dsinls;
00934     dyhd = drd * dsinls + drld * dcosls;
00935 
00936 
00937     /*
00938      * Influence of eccentricity, evection and variation of
00939      * the geocentric motion of the moon
00940      */
00941 
00942     pertl  = 0.;
00943     pertld = 0.;
00944     pertp  = 0.;
00945     pertpd = 0.;
00946 
00947     for (k = 0; k < 3; k++) {
00948 
00949         a = fmod(dcargm[k][0] + dt * dcargm[k][1], RV_D2PI);
00950         sina = sin(a);
00951         cosa = cos(a);
00952         pertl += ccampm[k][0] * sina;
00953         pertld += ccampm[k][1] * cosa;
00954         pertp += ccampm[k][2] * cosa;
00955         pertpd -= ccampm[k][3] * sina;
00956 
00957     }
00958 
00959 
00960     /*
00961      * Heliocentric motion of the earth
00962      */
00963 
00964     tl = forbel[1] + pertl;
00965     sinlm = sin(tl);
00966     coslm = cos(tl);
00967     sigma = cckm / (1. + pertp);
00968     a = sigma * (ccmld + pertld);
00969     b = sigma * pertpd;
00970     dxhd = dxhd + a * sinlm + b * coslm;
00971     dyhd = dyhd - a * coslm + b * sinlm;
00972     dzhd = -sigma * ccfdi * cos(forbel[2]);
00973 
00974 
00975     /*
00976      * Barycentric motion of the earth
00977      */
00978 
00979     dxbd = dxhd * dc1mme;
00980     dybd = dyhd * dc1mme;
00981     dzbd = dzhd * dc1mme;
00982 
00983     for (k = 0; k < 4; k++) {
00984 
00985         plon = forbel[k + 3];
00986         pomg = sorbel[k + 1];
00987         pecc = sorbel[k + 9];
00988         tl = fmod(plon + 2.0 * pecc * sin (plon - pomg), RV_D2PI);
00989         sinlp[k] = sin(tl);
00990         coslp[k] = cos(tl);
00991         dxbd = dxbd + ccpamv[k] * (sinlp[k] + pecc * sin(pomg));
00992         dybd = dybd - ccpamv[k] * (coslp[k] + pecc * cos(pomg));
00993         dzbd = dzbd - ccpamv[k] * sorbel[k + 13] * cos(plon - sorbel[k + 5]);
00994 
00995     }
00996 
00997 
00998     /*
00999      * Transition to mean equator of date
01000      */
01001 
01002     dcosep = cos(deps);
01003     dsinep = sin(deps);
01004     dyahd = dcosep * dyhd - dsinep * dzhd;
01005     dzahd = dsinep * dyhd + dcosep * dzhd;
01006     dyabd = dcosep * dybd - dsinep * dzbd;
01007     dzabd = dsinep * dybd + dcosep * dzbd;
01008 
01009     if (ideq == 0) {
01010 
01011         hvel[0] = dxhd;
01012         hvel[1] = dyahd;
01013         hvel[2] = dzahd;
01014 
01015         bvel[0] = dxbd;
01016         bvel[1] = dyabd;
01017         bvel[2] = dzabd;
01018 
01019     }
01020     else {
01021 
01022         /*
01023          * General precession from epoch dje to deq
01024          */
01025 
01026         cxdouble deqdat = (dje - dct0 - dcbes) / dctrop + dc1900;
01027 
01028         cpl_matrix* prec = slaPrecession(deqdat, deq);
01029 
01030 
01031         for (n = 0; n < 3; n++) {
01032 
01033             hvel[n] =
01034                 dxhd  * cpl_matrix_get(prec, 0, n) +
01035                 dyahd * cpl_matrix_get(prec, 1, n) +
01036                 dzahd * cpl_matrix_get(prec, 2, n);
01037 
01038             bvel[n] =
01039                 dxbd  * cpl_matrix_get(prec, 0, n) +
01040                 dyabd * cpl_matrix_get(prec, 1, n) +
01041                 dzabd * cpl_matrix_get(prec, 2, n);
01042         }
01043 
01044         cpl_matrix_delete(prec);
01045 
01046     }
01047 
01048     return;
01049 
01050 }
01051 
01052 
01053 /*
01054  * Public functions
01055  */
01056 
01090 void
01091 giraffe_rvcorrection_compute(GiRvCorrection* rv,
01092                              cxdouble jdate, cxdouble longitude,
01093                              cxdouble latitude, cxdouble elevation,
01094                              cxdouble ra, cxdouble dec,
01095                              cxdouble equinox)
01096 {
01097 
01098     cxint i = 0;
01099 
01100     const cxdouble aukm = 1.4959787e08;   /* 1 UA = 149 597 870 km */
01101 
01102     cxdouble eqt    = 0.;
01103     cxdouble ha     = 0.;
01104     cxdouble ra2    = 0.;
01105     cxdouble dec2   = 0.;
01106     cxdouble dc[3]  = {0., 0., 0.};
01107     cxdouble dcc[3] = {0., 0., 0.};
01108     cxdouble hv[3]  = {0., 0., 0.};
01109     cxdouble bv[3]  = {0., 0., 0.};
01110     cxdouble _long  = longitude * RV_DD2R;
01111     cxdouble _lat   = latitude * RV_DD2R;
01112     cxdouble _ra    = ra * 15.0 * RV_DD2R;
01113     cxdouble _dec   = dec * RV_DD2R;
01114     cxdouble st     = sideral_time(jdate, _long);
01115 
01116     cpl_matrix* precession = NULL;
01117 
01118 
01119     /*
01120      * Precess r.a. and dec. to mean equator and equinox of date
01121      */
01122 
01123     eqt = (jdate - dct0 - dcbes) / dctrop + dc1900;
01124 
01125     dc[0] = cos(_ra) * cos(_dec);
01126     dc[1] = sin(_ra) * cos(_dec);
01127     dc[2] = sin(_dec);
01128 
01129     precession = slaPrecession(equinox, eqt);
01130 
01131     for (i = 0; i < 3; ++i) {
01132 
01133         dcc[i] =
01134             dc[0] * cpl_matrix_get(precession, i, 0) +
01135             dc[1] * cpl_matrix_get(precession, i, 1) +
01136             dc[2] * cpl_matrix_get(precession, i, 2);
01137 
01138     }
01139 
01140     cpl_matrix_delete(precession);
01141     precession = NULL;
01142 
01143 
01144     if (dcc[0] != 0.) {
01145 
01146         cxdouble darg = dcc[1] / dcc[0];
01147 
01148 
01149         ra2 = atan(darg);
01150 
01151         if (dcc[0] < 0.) {
01152             ra2 += RV_DPI;
01153         }
01154         else {
01155             if (dcc[1] < 0.) {
01156                 ra2 += RV_D2PI;
01157             }
01158         }
01159 
01160     }
01161     else {
01162 
01163         if (dcc[1] > 0.) {
01164             ra2 = RV_DPIBY2;
01165         }
01166         else {
01167             ra2 = 1.5 * RV_DPI;
01168         }
01169 
01170     }
01171 
01172     dec2 = asin(dcc[2]);
01173 
01174 
01175     /*
01176      * Calculate hour angle = local sideral time - r.a.
01177      */
01178 
01179     ha = st - ra2;
01180 
01181 
01182     /*
01183      * Calculate observer's geocentric velocity
01184      * (elevation assumed to be zero)
01185      */
01186 
01187     rv->gc = geo_correction(_lat, elevation, dec2, -ha);
01188 
01189 
01190     /*
01191      * Calculate components of earth's barycentric velocity,
01192      * bvel(i), i=1,2,3  in units of a.u./s
01193      */
01194 
01195     earth_velocity (jdate, eqt, hv, bv);
01196 
01197 
01198     /*
01199      * Project barycentric velocity to the direction of the
01200      * star, and convert to km/s
01201      */
01202 
01203     rv->bc = 0.;
01204     rv->hc = 0.;
01205 
01206     for (i = 0; i < 3; ++i) {
01207         rv->bc += bv[i] * dcc[i] * aukm;
01208         rv->hc += hv[i] * dcc[i] * aukm;
01209     }
01210 
01211     return;
01212 
01213 }

This file is part of the GIRAFFE Pipeline Reference Manual 2.9.0.
Documentation copyright © 2002-2006 European Southern Observatory.
Generated on Thu Jan 26 14:20:28 2012 by doxygen 1.6.3 written by Dimitri van Heesch, © 1997-2004