GIRAFFE Pipeline Reference Manual

girvcorrection.c
1 /* $Id$
2  *
3  * This file is part of the GIRAFFE Pipeline
4  * Copyright (C) 2002-2006 European Southern Observatory
5  *
6  * This program is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * This program is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this program; if not, write to the Free Software
18  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
19  */
20 
21 /*
22  * $Author$
23  * $Date$
24  * $Revision$
25  * $Name$
26  */
27 
28 #ifdef HAVE_CONFIG_H
29 # include <config.h>
30 #endif
31 
32 
33 #include <math.h>
34 #include <string.h>
35 
36 #include <cxtypes.h>
37 
38 #include <cpl_matrix.h>
39 
40 #include "girvcorrection.h"
41 
42 
43 
52 /*
53  * Constants
54  */
55 
56 static const cxdouble dct0 = 2415020.0;
57 static const cxdouble dcjul = 36525.0;
58 static const cxdouble dc1900 = 1900.0;
59 static const cxdouble dctrop = 365.24219572;
60 static const cxdouble dcbes = 0.313;
61 
62 static const cxdouble RV_DPI =
63  3.1415926535897932384626433832795028841971693993751; /* pi */
64 
65 static const cxdouble RV_D2PI =
66  6.2831853071795864769252867665590057683943387987502; /* 2pi */
67 
68 static const cxdouble RV_D4PI =
69  12.566370614359172953850573533118011536788677597500; /* 4pi */
70 
71 static const cxdouble RV_DPIBY2 =
72  1.5707963267948966192313216916397514420985846996876; /* pi/2 */
73 
74 static const cxdouble RV_DD2R =
75  0.017453292519943295769236907684886127134428718885417; /* pi/180 */
76 
77 static const cxdouble RV_DAS2R =
78  4.8481368110953599358991410235794797595635330237270e-6; /* pi/(180*3600) */
79 
80 
81 /*
82  * Slalib functions
83  */
84 
85 /*
86  * - - - - - - - - - -
87  * s l a D e u l e r
88  * - - - - - - - - - -
89  *
90  * Form a rotation matrix from the Euler angles - three successive
91  * rotations about specified Cartesian axes.
92  *
93  * (double precision)
94  *
95  * Given:
96  * *order char specifies about which axes the rotations occur
97  * phi double 1st rotation (radians)
98  * theta double 2nd rotation ( " )
99  * psi double 3rd rotation ( " )
100  *
101  * Returned:
102  * Matrix rmat rotation matrix
103  *
104  * A rotation is positive when the reference frame rotates
105  * anticlockwise as seen looking towards the origin from the
106  * positive region of the specified axis.
107  *
108  * The characters of order define which axes the three successive
109  * rotations are about. A typical value is 'zxz', indicating that
110  * rmat is to become the direction cosine matrix corresponding to
111  * rotations of the reference frame through phi radians about the
112  * old z-axis, followed by theta radians about the resulting x-axis,
113  * then psi radians about the resulting z-axis.
114  *
115  * The axis names can be any of the following, in any order or
116  * combination: x, y, z, uppercase or lowercase, 1, 2, 3. Normal
117  * axis labelling/numbering conventions apply; the xyz (=123)
118  * triad is right-handed. Thus, the 'zxz' example given above
119  * could be written 'zxz' or '313' (or even 'zxz' or '3xz'). Order
120  * is terminated by length or by the first unrecognized character.
121  *
122  * Fewer than three rotations are acceptable, in which case the later
123  * angle arguments are ignored. Zero rotations leaves rmat set to the
124  * identity matrix.
125  *
126  * Last revision: 9 December 1996
127  *
128  * Copyright P.T.Wallace. All rights reserved.
129  */
130 
131 inline static void
132 slaDeuler(const cxchar* order, cxdouble phi, cxdouble theta, cxdouble psi,
133  cpl_matrix* rmat)
134 {
135  register cxint j, i, l, n, k;
136  cxdouble result[3][3], rotn[3][3], angle, s, c , w, wm[3][3];
137  cxchar axis;
138 
139  /* Initialize result matrix */
140  for ( j = 0; j < 3; j++ ) {
141  for ( i = 0; i < 3; i++ ) {
142  result[i][j] = ( i == j ) ? 1.0 : 0.0;
143  }
144  }
145 
146  /* Establish length of axis string */
147  l = strlen ( order );
148 
149  /* Look at each character of axis string until finished */
150  for ( n = 0; n < 3; n++ ) {
151  if ( n <= l ) {
152 
153  /* Initialize rotation matrix for the current rotation */
154  for ( j = 0; j < 3; j++ ) {
155  for ( i = 0; i < 3; i++ ) {
156  rotn[i][j] = ( i == j ) ? 1.0 : 0.0;
157  }
158  }
159 
160  /* Pick up the appropriate Euler angle and take sine & cosine */
161  switch ( n ) {
162  case 0 :
163  angle = phi;
164  break;
165  case 1 :
166  angle = theta;
167  break;
168  default:
169  angle = psi;
170  break;
171  }
172  s = sin ( angle );
173  c = cos ( angle );
174 
175  /* Identify the axis */
176  axis = order[n];
177  if ( ( axis == 'X' ) || ( axis == 'x' ) || ( axis == '1' ) ) {
178 
179  /* Matrix for x-rotation */
180  rotn[1][1] = c;
181  rotn[1][2] = s;
182  rotn[2][1] = -s;
183  rotn[2][2] = c;
184  }
185  else if ( ( axis == 'Y' ) || ( axis == 'y' ) || ( axis == '2' ) ) {
186 
187  /* Matrix for y-rotation */
188  rotn[0][0] = c;
189  rotn[0][2] = -s;
190  rotn[2][0] = s;
191  rotn[2][2] = c;
192  }
193  else if ( ( axis == 'Z' ) || ( axis == 'z' ) || ( axis == '3' ) ) {
194 
195  /* Matrix for z-rotation */
196  rotn[0][0] = c;
197  rotn[0][1] = s;
198  rotn[1][0] = -s;
199  rotn[1][1] = c;
200  }
201  else {
202 
203  /* Unrecognized character - fake end of string */
204  l = 0;
205  }
206 
207  /* Apply the current rotation (matrix rotn x matrix result) */
208  for ( i = 0; i < 3; i++ ) {
209  for ( j = 0; j < 3; j++ ) {
210  w = 0.0;
211  for ( k = 0; k < 3; k++ ) {
212  w += rotn[i][k] * result[k][j];
213  }
214  wm[i][j] = w;
215  }
216  }
217  for ( j = 0; j < 3; j++ ) {
218  for ( i= 0; i < 3; i++ ) {
219  result[i][j] = wm[i][j];
220  }
221  }
222  }
223  }
224 
225  /* Copy the result */
226  for ( j = 0; j < 3; j++ ) {
227  for ( i = 0; i < 3; i++ ) {
228  cpl_matrix_set(rmat, i, j, result[i][j]);
229  }
230  }
231 
232  return;
233 
234 }
235 
236 
237 /*
238  * @brief
239  * Calculates the matrix of general precession from ep0 to ep1
240  *
241  * @param ep0 Initial epoch of mean equator and equinox e.g., 1950.0
242  * @param ep1 Final epoch of mean equator and equinox
243  *
244  * @return
245  * The 3 x 3 matrix of general precession.
246  *
247  * Form the matrix of precession between two epochs, using the
248  * model of Simon et al (1994), which is suitable for long
249  * periods of time.
250  *
251  * @notes
252  * -# The epochs are TDB (loosely ET) Julian epochs.
253  * -# The matrix is in the sense v(ep1) = rmatp * v(ep0)
254  * -# The absolute accuracy of the model is limited by the
255  * uncertainty in the general precession, about 0.3 arcsec per
256  * 1000 years. The remainder of the formulation provides a
257  * precision of 1 mas over the interval from 1000AD to 3000AD,
258  * 0.1 arcsec from 1000BC to 5000AD and 1 arcsec from
259  * 4000BC to 8000AD.
260  *
261  * Reference:
262  * Simon, J. L., et al., 1994. Astron.Astrophys., 282, 663-683.
263  *
264  * Author:
265  * P. T. Wallace. 23 August 1994
266  */
267 
268 inline static cpl_matrix*
269 slaPrecession (cxdouble ep0, cxdouble ep1)
270 {
271 
272  cxdouble t = 0.;
273  cxdouble w = 0.;
274  cxdouble z = 0.;
275  cxdouble t0 = 0.;
276  cxdouble zeta = 0.;
277  cxdouble theta = 0.;
278  cxdouble tas2r = 0.;
279 
280  cpl_matrix* mprec = NULL;
281 
282 
283  /*
284  * Interval between basic epoch J2000.0 and beginning epoch (1000JY)
285  */
286 
287  t0 = ( ep0 - 2000.0 ) / 1000.0;
288 
289 
290  /*
291  * Interval over which precession required (1000JY)
292  */
293 
294  t = ( ep1 - ep0 ) / 1000.0;
295 
296 
297  /*
298  * Euler angles
299  */
300 
301  tas2r = t * RV_DAS2R;
302  w = 23060.9097 +
303  (139.7459 +
304  (-0.0038 +
305  (-0.5918 +
306  (-0.0037 +
307  0.0007 * t0) * t0) * t0) * t0) * t0;
308 
309  zeta = (w +
310  (30.2226 +
311  (-0.2523 +
312  (-0.3840 +
313  (-0.0014 +
314  0.0007 * t0) * t0) * t0) * t0 +
315  (18.0183 +
316  (-0.1326 +
317  (0.0006 +
318  0.0005 * t0) * t0) * t0 +
319  (-0.0583 +
320  (-0.0001 +
321  0.0007 * t0) * t0 +
322  (-0.0285 +
323  -0.0002 * t) * t) * t) * t) * t) * tas2r;
324 
325  z = (w +
326  (109.5270 +
327  (0.2446 +
328  (-1.3913 +
329  (-0.0134 +
330  0.0026 * t0) * t0) * t0) * t0 +
331  (18.2667 +
332  (-1.1400 +
333  (-0.0173 +
334  0.0044 * t0) * t0) * t0 +
335  (-0.2821 +
336  (-0.0093 +
337  0.0032 * t0) * t0 +
338  (-0.0301 +
339  0.0006 * t0 +
340  -0.0001 * t) * t) * t) * t) * t) * tas2r;
341 
342  theta = (20042.0207 +
343  (-85.3131 +
344  (-0.2111 +
345  (0.3642 +
346  (0.0008 +
347  -0.0005 * t0) * t0) * t0) * t0) * t0 +
348  (-42.6566 +
349  (-0.2111 +
350  (0.5463 +
351  (0.0017 +
352  -0.0012 * t0) * t0) * t0) * t0 +
353  (-41.8238 +
354  (0.0359 +
355  (0.0027 +
356  -0.0001 * t0) * t0) * t0 +
357  (-0.0731 +
358  (0.0019 +
359  0.0009 * t0) * t0 +
360  (-0.0127 +
361  0.0011 * t0 + 0.0004 * t) * t) * t) * t) * t) * tas2r;
362 
363 
364  /*
365  * Rotation matrix
366  */
367 
368  mprec = cpl_matrix_new(3, 3);
369  slaDeuler("ZYZ", -zeta, theta, -z, mprec);
370 
371  return mprec;
372 
373 }
374 
375 
376 
377 /*
378  * @brief
379  * Compute the mean local sideral time
380  *
381  * @param djd julian date
382  * @param dlong observer's longitude (radians)
383  *
384  * @return
385  * Mean local sideral time (radians).
386  *
387  * @note
388  * Constants taken from the american ephemeris and nautical almanac, 1980)
389  * Translated from the FORTRAN version written by G. Torres (1989)
390  */
391 
392 inline static cxdouble
393 sideral_time(cxdouble djd, cxdouble dlong)
394 {
395 
396  /*
397  * Constants D1,D2,D3 for calculating Greenwich
398  * Mean Sideral Time at 0 UT
399  */
400 
401  const cxdouble d1 = 1.739935934667999;
402  const cxdouble d2 = 6.283319509909095e02;
403  const cxdouble d3 = 6.755878646261384e-06;
404 
405  const cxdouble df = 1.00273790934;
406 
407  cxdouble dut = 0.;
408  cxdouble dt = 0.;
409  cxdouble dst0 = 0.;
410  cxdouble dst = 0.;
411  cxdouble djd0 = floor(djd) + 0.5;
412 
413  if (djd0 > djd) {
414  djd0 -= 1.;
415  }
416 
417  dut = (djd - djd0) * RV_D2PI;
418 
419  dt = (djd0 - dct0) / dcjul;
420  dst0 = d1 + d2 * dt + d3 * dt * dt;
421  dst0 = fmod(dst0, RV_D2PI);
422  dst = df * dut + dst0 - dlong;
423  dst = fmod(dst + RV_D4PI, RV_D2PI);
424 
425  return dst;
426 
427 }
428 
429 
430 /*
431  * @brief
432  * Compute correction from topocentric radial velocity to geocentric.
433  *
434  * @param dlat observer's geodetic latitude (radians)
435  * @param dalt observer's elevation above sea level (meters)
436  * @param dec star's declination (radians) for mean
437  * equator and equinox of date
438  * @param dha star's hour angle (radians)
439  *
440  * @return
441  * Geocentric correction in km/s.
442  *
443  * Calculates the correction required to transform the topocentric radial
444  * velocity of a given star to geocentric. The maximum error of this
445  * routine is not expected to be larger than 0.1 m/s.
446  *
447  * @note
448  * vr = r.w.cos(dec).sin(hour angle), where r = geocentric radius at
449  * observer's latitude and elevation, and w = 2.pi/t, t = length of sideral
450  * day 23 hours 56 minutes and 4 seconds in seconds.
451  * The hour angle is positive east of the meridian.
452  * Other relevant equations from E. W. Woolard & G. M. Clemence (1966),
453  * Spherical Astronomy, p.45 and p.48
454  *
455  * Author:
456  * - G. Torres (1989)
457  * - G. Simond (C translation)
458  */
459 
460 inline static cxdouble
461 geo_correction(cxdouble dlat, cxdouble dalt, cxdouble dec, cxdouble dha)
462 {
463 
464  /*
465  * Earth's equatorial radius (km) (Geod. ref. sys., 1980)
466  */
467 
468  const cxdouble da = 6378.137;
469 
470  /*
471  * Polar flattening (Geod. ref. sys., 1980)
472  */
473 
474  const cxdouble df = 1./298.257222;
475 
476  /*
477  * Angular rotation rate (2.pi/t)
478  */
479 
480  const cxdouble dw = RV_D2PI/86164.;
481 
482 
483  const cxdouble de2 = df * (2.0 - df);
484  const cxdouble dsdlats = sin (dlat) * sin (dlat);
485 
486  cxdouble d1 = 0.;
487  cxdouble d2 = 0.;
488  cxdouble dr0 = 0.;
489  cxdouble dlatg = 0.;
490  cxdouble drh = 0.;
491  cxdouble dvelg = 0.;
492 
493 
494  /*
495  * Calculate geocentric radius dr0 at sea level (km)
496  */
497 
498  d1 = 1.0 - de2 * (2.0 - de2) * dsdlats;
499  d2 = 1.0 - de2 * dsdlats;
500  dr0 = da * sqrt(d1 / d2);
501 
502 
503  /*
504  * Calculate geocentric latitude dlatg
505  */
506 
507  d1 = de2 * sin(2.0 * dlat);
508  d2 = 2.0 * d2;
509  dlatg = dlat - atan(d1 / d2);
510 
511 
512  /*
513  * Calculate geocentric radius drh at elevation dalt (km)
514  */
515 
516  drh = dr0 * cos(dlatg) + (dalt / 1000.) * cos(dlat);
517 
518 
519  /*
520  * Projected component to star at declination = dec and
521  * at hour angle = dha, in units of km/s
522  */
523 
524  dvelg = dw * drh * cos(dec) * sin(dha);
525 
526  return dvelg;
527 
528 }
529 
530 
531 /*
532  * @brief
533  * Compute the heliocentric and barycentric velocity components of the earth.
534  *
535  * @param dje Julian ephermeris date
536  * @param deq Epoch of mean equator and mean equinox of @em hvel
537  * and @em bvel. If @em deq is @c 0, both vectors refer
538  * to the mean equator and equinox of @em dje
539  *
540  * @return
541  * The components (dx/dt, dy/dt, dz/dt, k=1,2,3 A.U./s) of the
542  * heliocentric and barycentric velocity are returned in the
543  * vectors @em hvel and @em bvel respectively.
544  *
545  * Calculates the heliocentric and barycentric velocity components
546  * of the earth. The largest deviations from the JPL-de96 ephemeris
547  * are 42 cm/s for both heliocentric and barycentric velocity
548  * components.
549  *
550  * @note
551  * vr = r.w.cos(dec).sin(hour angle), where r = geocentric radius at
552  * observer's latitude and elevation, and w = 2.pi/t, t = length of sideral
553  * day (sec). The hour angle is positive east of the meridian.
554  * Other relevant equations from E. W. Woolard & G. M. Clemence (1966),
555  * Spherical Astronomy, p.45 and p.48
556  *
557  * Authors:
558  * - P. Stumpff (ibm-version 1979): astron. astrophys. suppl. ser. 41,
559  * 1 (1980)
560  * - M. H.Slovak (vax 11/780 implementation 1986)
561  * - G. Torres (1989)
562  * - G. Simond (C translation)
563  */
564 
565 inline static void
566 earth_velocity(cxdouble dje, cxdouble deq, cxdouble* const hvel,
567  cxdouble* const bvel)
568 {
569 
570 
571  /*
572  * Constants dcfel(i, k) of fast-changing elements
573  *
574  * i = 1 i = 2 i = 3
575  */
576 
577  const cxdouble dcfel[][3] = {
578  {1.7400353e+00, 6.2833195099091e+02, 5.2796e-06},
579  {6.2565836e+00, 6.2830194572674e+02, -2.6180e-06},
580  {4.7199666e+00, 8.3997091449254e+03, -1.9780e-05},
581  {1.9636505e-01, 8.4334662911720e+03, -5.6044e-05},
582  {4.1547339e+00, 5.2993466764997e+01, 5.8845e-06},
583  {4.6524223e+00, 2.1354275911213e+01, 5.6797e-06},
584  {4.2620486e+00, 7.5025342197656e+00, 5.5317e-06},
585  {1.4740694e+00, 3.8377331909193e+00, 5.6093e-06}
586  };
587 
588 
589  /*
590  * Constants dceps and ccsel(i,k) of slowly changing elements
591  *
592  * i = 1 i = 2 i = 3
593  */
594 
595  const cxdouble dceps[3] = {
596  4.093198e-01,
597  -2.271110e-04,
598  -2.860401e-08
599  };
600 
601  const cxdouble ccsel[][3] = {
602  {1.675104e-02, -4.179579e-05, -1.260516e-07},
603  {2.220221e-01, 2.809917e-02, 1.852532e-05},
604  {1.589963e+00, 3.418075e-02, 1.430200e-05},
605  {2.994089e+00, 2.590824e-02, 4.155840e-06},
606  {8.155457e-01, 2.486352e-02, 6.836840e-06},
607  {1.735614e+00, 1.763719e-02, 6.370440e-06},
608  {1.968564e+00, 1.524020e-02, -2.517152e-06},
609  {1.282417e+00, 8.703393e-03, 2.289292e-05},
610  {2.280820e+00, 1.918010e-02, 4.484520e-06},
611  {4.833473e-02, 1.641773e-04, -4.654200e-07},
612  {5.589232e-02, -3.455092e-04, -7.388560e-07},
613  {4.634443e-02, -2.658234e-05, 7.757000e-08},
614  {8.997041e-03, 6.329728e-06, -1.939256e-09},
615  {2.284178e-02, -9.941590e-05, 6.787400e-08},
616  {4.350267e-02, -6.839749e-05, -2.714956e-07},
617  {1.348204e-02, 1.091504e-05, 6.903760e-07},
618  {3.106570e-02, -1.665665e-04, -1.590188e-07}
619  };
620 
621 
622  /*
623  * Constants of the arguments of the short-period perturbations by
624  * the planets: dcargs(i, k)
625  *
626  * i = 1 i = 2
627  */
628 
629  const cxdouble dcargs[][2] = {
630  {5.0974222e+00, -7.8604195454652e+02},
631  {3.9584962e+00, -5.7533848094674e+02},
632  {1.6338070e+00, -1.1506769618935e+03},
633  {2.5487111e+00, -3.9302097727326e+02},
634  {4.9255514e+00, -5.8849265665348e+02},
635  {1.3363463e+00, -5.5076098609303e+02},
636  {1.6072053e+00, -5.2237501616674e+02},
637  {1.3629480e+00, -1.1790629318198e+03},
638  {5.5657014e+00, -1.0977134971135e+03},
639  {5.0708205e+00, -1.5774000881978e+02},
640  {3.9318944e+00, 5.2963464780000e+01},
641  {4.8989497e+00, 3.9809289073258e+01},
642  {1.3097446e+00, 7.7540959633708e+01},
643  {3.5147141e+00, 7.9618578146517e+01},
644  {3.5413158e+00, -5.4868336758022e+02}
645  };
646 
647 
648  /*
649  * Amplitudes ccamps(n, k) of the short-period perturbations
650  *
651  * n = 1 n = 2 n = 3 n = 4 n = 5
652  */
653 
654  const cxdouble ccamps[][5] = {
655  {-2.279594e-5, 1.407414e-5, 8.273188e-6, 1.340565e-5, -2.490817e-7},
656  {-3.494537e-5, 2.860401e-7, 1.289448e-7, 1.627237e-5, -1.823138e-7},
657  { 6.593466e-7, 1.322572e-5, 9.258695e-6, -4.674248e-7, -3.646275e-7},
658  { 1.140767e-5, -2.049792e-5, -4.747930e-6, -2.638763e-6, -1.245408e-7},
659  { 9.516893e-6, -2.748894e-6, -1.319381e-6, -4.549908e-6, -1.864821e-7},
660  { 7.310990e-6, -1.924710e-6, -8.772849e-7, -3.334143e-6, -1.745256e-7},
661  {-2.603449e-6, 7.359472e-6, 3.168357e-6, 1.119056e-6, -1.655307e-7},
662  {-3.228859e-6, 1.308997e-7, 1.013137e-7, 2.403899e-6, -3.736225e-7},
663  { 3.442177e-7, 2.671323e-6, 1.832858e-6, -2.394688e-7, -3.478444e-7},
664  { 8.702406e-6, -8.421214e-6, -1.372341e-6, -1.455234e-6, -4.998479e-8},
665  {-1.488378e-6, -1.251789e-5, 5.226868e-7, -2.049301e-7, 0.0e0},
666  {-8.043059e-6, -2.991300e-6, 1.473654e-7, -3.154542e-7, 0.0e0},
667  { 3.699128e-6, -3.316126e-6, 2.901257e-7, 3.407826e-7, 0.0e0},
668  { 2.550120e-6, -1.241123e-6, 9.901116e-8, 2.210482e-7, 0.0e0},
669  {-6.351059e-7, 2.341650e-6, 1.061492e-6, 2.878231e-7, 0.0e0}
670  };
671 
672 
673  /*
674  * Constants of the secular perturbations in longitude
675  * ccsec3 and ccsec(n,k)
676  * n = 1 n = 2 n = 3
677  */
678 
679  const cxdouble ccsec3 = -7.757020e-08;
680 
681  const cxdouble ccsec[][3] = {
682  {1.289600e-06, 5.550147e-01, 2.076942e+00},
683  {3.102810e-05, 4.035027e+00, 3.525565e-01},
684  {9.124190e-06, 9.990265e-01, 2.622706e+00},
685  {9.793240e-07, 5.508259e+00, 1.559103e+01}
686  };
687 
688 
689  /*
690  * Sideral rate dcsld in longitude, rate ccsgd in mean anomaly
691  */
692 
693  const cxdouble dcsld = 1.990987e-07;
694  const cxdouble ccsgd = 1.990969e-07;
695 
696 
697  /*
698  * Some constants used in the calculation of the
699  * lunar contribution
700  */
701 
702  const cxdouble cckm = 3.122140e-05;
703  const cxdouble ccmld = 2.661699e-06;
704  const cxdouble ccfdi = 2.399485e-07;
705 
706 
707  /*
708  * Constants dcargm(i,k) of the arguments of the
709  * perturbations of the motion of the moon
710  *
711  * i = 1 i = 2
712  */
713 
714  const cxdouble dcargm[][2] = {
715  {5.1679830e+00, 8.3286911095275e+03},
716  {5.4913150e+00, -7.2140632838100e+03},
717  {5.9598530e+00, 1.5542754389685e+04}
718  };
719 
720 
721  /*
722  * Amplitudes ccampm(n,k) of the perturbations of the moon
723  * n = 1 n = 2 n = 3 n = 4
724  */
725 
726  const cxdouble ccampm[][4] = {
727  { 1.097594e-01, 2.896773e-07, 5.450474e-02, 1.438491e-07},
728  {-2.223581e-02, 5.083103e-08, 1.002548e-02, -2.291823e-08},
729  { 1.148966e-02, 5.658888e-08, 8.249439e-03, 4.063015e-08}
730  };
731 
732 
733  /*
734  * ccpamv = a * m * dl/dt (planets)
735  */
736 
737  const cxdouble ccpamv[4] = {
738  8.326827e-11,
739  1.843484e-11,
740  1.988712e-12,
741  1.881276e-12
742  };
743 
744 
745  /*
746  * dc1mme = 1 - mass(earth + moon)
747  */
748 
749  const cxdouble dc1mme = 0.99999696;
750 
751 
752  register cxint k = 0;
753  register cxint n = 0;
754 
755  cxint ideq = 0;
756 
757  cxdouble a = 0.;
758  cxdouble b = 0.;
759  cxdouble f = 0.;
760  cxdouble dt = 0.;
761  cxdouble t = 0.;
762  cxdouble tl = 0.;
763  cxdouble dtsq = 0.;
764  cxdouble tsq = 0.;
765  cxdouble dml = 0.;
766  cxdouble dlocal = 0.;
767  cxdouble deps = 0.;
768  cxdouble pertl = 0.;
769  cxdouble pertld = 0.;
770  cxdouble pertr = 0.;
771  cxdouble pertrd = 0.;
772  cxdouble pertp = 0.;
773  cxdouble pertpd = 0.;
774  cxdouble sina = 0.;
775  cxdouble cosa = 0.;
776  cxdouble twoe = 0.;
777  cxdouble twog = 0.;
778  cxdouble param = 0.;
779  cxdouble dparam = 0.;
780  cxdouble dpsi = 0.;
781  cxdouble phi = 0.;
782  cxdouble phid = 0.;
783  cxdouble psid = 0.;
784  cxdouble sin_f = 0.;
785  cxdouble cos_f = 0.;
786  cxdouble esq = 0.;
787  cxdouble d1pdro = 0.;
788  cxdouble drd = 0.;
789  cxdouble drld = 0.;
790  cxdouble dsinls = 0.;
791  cxdouble dcosls = 0.;
792  cxdouble dtl = 0.;
793  cxdouble dxhd = 0.;
794  cxdouble dyhd = 0.;
795  cxdouble dzhd = 0.;
796  cxdouble sinlm = 0.;
797  cxdouble coslm = 0.;
798  cxdouble sigma = 0.;
799  cxdouble plon = 0.;
800  cxdouble pomg = 0.;
801  cxdouble dxbd = 0.;
802  cxdouble dybd = 0.;
803  cxdouble dzbd = 0.;
804  cxdouble pecc = 0.;
805  cxdouble dcosep = 0.;
806  cxdouble dsinep = 0.;
807  cxdouble dyahd = 0.;
808  cxdouble dzahd = 0.;
809  cxdouble dyabd = 0.;
810  cxdouble dzabd = 0.;
811  cxdouble sn[4] = {0., 0., 0., 0.};
812  cxdouble sinlp[4] = {0., 0., 0., 0.};
813  cxdouble coslp[4] = {0., 0., 0., 0.};
814  cxdouble forbel[7] = {0., 0., 0., 0., 0., 0., 0.};
815  cxdouble sorbel[17];
816 
817 
818  memset(sorbel, 0, sizeof sorbel);
819 
820 
821  /*
822  * Control-parameter ideq, and time-arguments
823  */
824 
825  ideq = (cxint)deq;
826  dt = (dje - dct0) / dcjul;
827  t = dt;
828  dtsq = dt * dt;
829  tsq = dtsq;
830 
831 
832  /*
833  * Values of all elements for the instant dje
834  */
835 
836  for (k = 0; k < 8; k++) {
837 
838  dlocal = fmod(dcfel[k][0] + dt * dcfel[k][1] + dtsq * dcfel[k][2],
839  RV_D2PI);
840 
841  if (k == 0) {
842  dml = dlocal;
843  }
844 
845  if (k != 0) {
846  forbel[k - 1] = dlocal;
847  }
848 
849  }
850 
851  deps = fmod(dceps[0] + dt * dceps[1] + dtsq * dceps[2], RV_D2PI);
852 
853  for (k = 0; k < 17; k++) {
854 
855  sorbel[k] = fmod(ccsel[k][0] + t * ccsel[k][1] + tsq * ccsel[k][2],
856  RV_D2PI);
857 
858  }
859 
860 
861  /*
862  * Secular perturbations in longitude
863  */
864 
865  for (k = 0; k < 4; k++) {
866 
867  a = fmod(ccsec[k][1] + t * ccsec[k][2], RV_D2PI);
868  sn[k] = sin(a);
869 
870  }
871 
872 
873  /*
874  * Periodic perturbations of the earth-moon barycenter
875  */
876 
877  pertl = ccsec[0][0] * sn[0] + ccsec[1][0] * sn[1] +
878  (ccsec[2][0] + t * ccsec3) * sn[2] + ccsec[3][0] * sn[3];
879 
880  pertld = 0.;
881  pertr = 0.;
882  pertrd = 0.;
883 
884  for (k = 0; k < 15; k++) {
885 
886  a = fmod(dcargs[k][0] + dt * dcargs[k][1], RV_D2PI);
887  cosa = cos (a);
888  sina = sin (a);
889  pertl += (ccamps[k][0] * cosa + ccamps[k][1] * sina);
890  pertr += (ccamps[k][2] * cosa + ccamps[k][3] * sina);
891 
892  if (k >= 10) {
893  continue;
894  }
895 
896  pertld += ((ccamps[k][1] * cosa - ccamps[k][0] * sina) * ccamps[k][4]);
897  pertrd += ((ccamps[k][3] * cosa - ccamps[k][2] * sina) * ccamps[k][4]);
898 
899  }
900 
901 
902  /*
903  * Elliptic part of the motion of the earth-moon barycenter
904  */
905 
906  esq = sorbel[0] * sorbel[0];
907  dparam = 1. - esq;
908  param = dparam;
909  twoe = sorbel[0] + sorbel[0];
910  twog = forbel[0] + forbel[0];
911  phi = twoe * ((1.0 - esq * (1.0 / 8.0)) * sin (forbel[0]) +
912  sorbel[0] * (5.0 / 8.0) * sin (twog) +
913  esq * 0.5416667 * sin (forbel[0] + twog));
914  f = forbel[0] + phi;
915  sin_f = sin(f);
916  cos_f = cos(f);
917  dpsi = dparam / (1. + sorbel[0] * cos_f);
918  phid = twoe * ccsgd * ((1.0 + esq * 1.50) * cos_f +
919  sorbel[0] * (1.250 - sin_f * sin_f * 0.50));
920  psid = ccsgd * sorbel[0] * sin_f / sqrt(param);
921 
922 
923  /*
924  * Perturbed heliocentric motion of the earth-moon barycenter
925  */
926 
927  d1pdro = 1. + pertr;
928  drd = d1pdro * (psid + dpsi * pertrd);
929  drld = d1pdro * dpsi * (dcsld + phid + pertld);
930  dtl = fmod(dml + phi + pertl, RV_D2PI);
931  dsinls = sin(dtl);
932  dcosls = cos(dtl);
933  dxhd = drd * dcosls - drld * dsinls;
934  dyhd = drd * dsinls + drld * dcosls;
935 
936 
937  /*
938  * Influence of eccentricity, evection and variation of
939  * the geocentric motion of the moon
940  */
941 
942  pertl = 0.;
943  pertld = 0.;
944  pertp = 0.;
945  pertpd = 0.;
946 
947  for (k = 0; k < 3; k++) {
948 
949  a = fmod(dcargm[k][0] + dt * dcargm[k][1], RV_D2PI);
950  sina = sin(a);
951  cosa = cos(a);
952  pertl += ccampm[k][0] * sina;
953  pertld += ccampm[k][1] * cosa;
954  pertp += ccampm[k][2] * cosa;
955  pertpd -= ccampm[k][3] * sina;
956 
957  }
958 
959 
960  /*
961  * Heliocentric motion of the earth
962  */
963 
964  tl = forbel[1] + pertl;
965  sinlm = sin(tl);
966  coslm = cos(tl);
967  sigma = cckm / (1. + pertp);
968  a = sigma * (ccmld + pertld);
969  b = sigma * pertpd;
970  dxhd = dxhd + a * sinlm + b * coslm;
971  dyhd = dyhd - a * coslm + b * sinlm;
972  dzhd = -sigma * ccfdi * cos(forbel[2]);
973 
974 
975  /*
976  * Barycentric motion of the earth
977  */
978 
979  dxbd = dxhd * dc1mme;
980  dybd = dyhd * dc1mme;
981  dzbd = dzhd * dc1mme;
982 
983  for (k = 0; k < 4; k++) {
984 
985  plon = forbel[k + 3];
986  pomg = sorbel[k + 1];
987  pecc = sorbel[k + 9];
988  tl = fmod(plon + 2.0 * pecc * sin (plon - pomg), RV_D2PI);
989  sinlp[k] = sin(tl);
990  coslp[k] = cos(tl);
991  dxbd = dxbd + ccpamv[k] * (sinlp[k] + pecc * sin(pomg));
992  dybd = dybd - ccpamv[k] * (coslp[k] + pecc * cos(pomg));
993  dzbd = dzbd - ccpamv[k] * sorbel[k + 13] * cos(plon - sorbel[k + 5]);
994 
995  }
996 
997 
998  /*
999  * Transition to mean equator of date
1000  */
1001 
1002  dcosep = cos(deps);
1003  dsinep = sin(deps);
1004  dyahd = dcosep * dyhd - dsinep * dzhd;
1005  dzahd = dsinep * dyhd + dcosep * dzhd;
1006  dyabd = dcosep * dybd - dsinep * dzbd;
1007  dzabd = dsinep * dybd + dcosep * dzbd;
1008 
1009  if (ideq == 0) {
1010 
1011  hvel[0] = dxhd;
1012  hvel[1] = dyahd;
1013  hvel[2] = dzahd;
1014 
1015  bvel[0] = dxbd;
1016  bvel[1] = dyabd;
1017  bvel[2] = dzabd;
1018 
1019  }
1020  else {
1021 
1022  /*
1023  * General precession from epoch dje to deq
1024  */
1025 
1026  cxdouble deqdat = (dje - dct0 - dcbes) / dctrop + dc1900;
1027 
1028  cpl_matrix* prec = slaPrecession(deqdat, deq);
1029 
1030 
1031  for (n = 0; n < 3; n++) {
1032 
1033  hvel[n] =
1034  dxhd * cpl_matrix_get(prec, 0, n) +
1035  dyahd * cpl_matrix_get(prec, 1, n) +
1036  dzahd * cpl_matrix_get(prec, 2, n);
1037 
1038  bvel[n] =
1039  dxbd * cpl_matrix_get(prec, 0, n) +
1040  dyabd * cpl_matrix_get(prec, 1, n) +
1041  dzabd * cpl_matrix_get(prec, 2, n);
1042  }
1043 
1044  cpl_matrix_delete(prec);
1045 
1046  }
1047 
1048  return;
1049 
1050 }
1051 
1052 
1053 /*
1054  * Public functions
1055  */
1056 
1090 void
1091 giraffe_rvcorrection_compute(GiRvCorrection* rv,
1092  cxdouble jdate, cxdouble longitude,
1093  cxdouble latitude, cxdouble elevation,
1094  cxdouble ra, cxdouble dec,
1095  cxdouble equinox)
1096 {
1097 
1098  cxint i = 0;
1099 
1100  const cxdouble aukm = 1.4959787e08; /* 1 UA = 149 597 870 km */
1101 
1102  cxdouble eqt = 0.;
1103  cxdouble ha = 0.;
1104  cxdouble ra2 = 0.;
1105  cxdouble dec2 = 0.;
1106  cxdouble dc[3] = {0., 0., 0.};
1107  cxdouble dcc[3] = {0., 0., 0.};
1108  cxdouble hv[3] = {0., 0., 0.};
1109  cxdouble bv[3] = {0., 0., 0.};
1110  cxdouble _long = longitude * RV_DD2R;
1111  cxdouble _lat = latitude * RV_DD2R;
1112  cxdouble _ra = ra * 15.0 * RV_DD2R;
1113  cxdouble _dec = dec * RV_DD2R;
1114  cxdouble st = sideral_time(jdate, _long);
1115 
1116  cpl_matrix* precession = NULL;
1117 
1118 
1119  /*
1120  * Precess r.a. and dec. to mean equator and equinox of date
1121  */
1122 
1123  eqt = (jdate - dct0 - dcbes) / dctrop + dc1900;
1124 
1125  dc[0] = cos(_ra) * cos(_dec);
1126  dc[1] = sin(_ra) * cos(_dec);
1127  dc[2] = sin(_dec);
1128 
1129  precession = slaPrecession(equinox, eqt);
1130 
1131  for (i = 0; i < 3; ++i) {
1132 
1133  dcc[i] =
1134  dc[0] * cpl_matrix_get(precession, i, 0) +
1135  dc[1] * cpl_matrix_get(precession, i, 1) +
1136  dc[2] * cpl_matrix_get(precession, i, 2);
1137 
1138  }
1139 
1140  cpl_matrix_delete(precession);
1141  precession = NULL;
1142 
1143 
1144  if (dcc[0] != 0.) {
1145 
1146  cxdouble darg = dcc[1] / dcc[0];
1147 
1148 
1149  ra2 = atan(darg);
1150 
1151  if (dcc[0] < 0.) {
1152  ra2 += RV_DPI;
1153  }
1154  else {
1155  if (dcc[1] < 0.) {
1156  ra2 += RV_D2PI;
1157  }
1158  }
1159 
1160  }
1161  else {
1162 
1163  if (dcc[1] > 0.) {
1164  ra2 = RV_DPIBY2;
1165  }
1166  else {
1167  ra2 = 1.5 * RV_DPI;
1168  }
1169 
1170  }
1171 
1172  dec2 = asin(dcc[2]);
1173 
1174 
1175  /*
1176  * Calculate hour angle = local sideral time - r.a.
1177  */
1178 
1179  ha = st - ra2;
1180 
1181 
1182  /*
1183  * Calculate observer's geocentric velocity
1184  * (elevation assumed to be zero)
1185  */
1186 
1187  rv->gc = geo_correction(_lat, elevation, dec2, -ha);
1188 
1189 
1190  /*
1191  * Calculate components of earth's barycentric velocity,
1192  * bvel(i), i=1,2,3 in units of a.u./s
1193  */
1194 
1195  earth_velocity (jdate, eqt, hv, bv);
1196 
1197 
1198  /*
1199  * Project barycentric velocity to the direction of the
1200  * star, and convert to km/s
1201  */
1202 
1203  rv->bc = 0.;
1204  rv->hc = 0.;
1205 
1206  for (i = 0; i < 3; ++i) {
1207  rv->bc += bv[i] * dcc[i] * aukm;
1208  rv->hc += hv[i] * dcc[i] * aukm;
1209  }
1210 
1211  return;
1212 
1213 }

This file is part of the GIRAFFE Pipeline Reference Manual 2.12.
Documentation copyright © 2002-2006 European Southern Observatory.
Generated on Mon Mar 24 2014 11:43:53 by doxygen 1.8.2 written by Dimitri van Heesch, © 1997-2004