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 #include <cxmemory.h>
38 #include <cxstrutils.h>
39 
40 #include <cpl_matrix.h>
41 
42 #include "girvcorrection.h"
43 
44 
45 
54 /*
55  * Constants
56  */
57 
58 static const cxdouble dct0 = 2415020.0;
59 static const cxdouble dcjul = 36525.0;
60 static const cxdouble dc1900 = 1900.0;
61 static const cxdouble dctrop = 365.24219572;
62 static const cxdouble dcbes = 0.313;
63 
64 static const cxdouble RV_DPI =
65  3.1415926535897932384626433832795028841971693993751; /* pi */
66 
67 static const cxdouble RV_D2PI =
68  6.2831853071795864769252867665590057683943387987502; /* 2pi */
69 
70 static const cxdouble RV_D4PI =
71  12.566370614359172953850573533118011536788677597500; /* 4pi */
72 
73 static const cxdouble RV_DPIBY2 =
74  1.5707963267948966192313216916397514420985846996876; /* pi/2 */
75 
76 static const cxdouble RV_DD2R =
77  0.017453292519943295769236907684886127134428718885417; /* pi/180 */
78 
79 static const cxdouble RV_DAS2R =
80  4.8481368110953599358991410235794797595635330237270e-6; /* pi/(180*3600) */
81 
82 
83 /*
84  * @brief
85  * Generate a rotation matrix from a set of Euler angles
86  *
87  * @param mode List of axis to rotate about (intrinsic rotations).
88  * @param phi Angle of the first rotation.
89  * @param theta Angle of the second rotation.
90  * @param psi Angle of the third rotation.
91  *
92  * @return
93  * The function returns the created 3x3 rotation matrix, or @c NULL in case
94  * an error occurs.
95  *
96  * The function creates a rotation matrix from the Euler angles @em phi,
97  * @em theta and @em psi, by applying three successive rotations of the
98  * respective angle about the cartesian axes given in the list @em mode.
99  *
100  * The rotations are intrinsic rotations and are carried out in the order
101  * the axes are given by @em mode. The axes designations which may appear
102  * in the list mode are 'x', 'y', and 'z' for the first, second, and third
103  * axis respectively.
104  *
105  * The angles are defined according to the right hand rule, i.e. they are
106  * positive if they represent a rotation which appears counter-clockwise
107  * when observed from a point on the positive part of the rotation axis, and
108  * they are negative for clockwise rotations.
109  */
110 
111 inline static cpl_matrix *
112 _giraffe_create_rotation(const cxchar *mode, cxdouble phi, cxdouble theta,
113  cxdouble psi)
114 {
115 
116  cxint nrotations = strlen(mode);
117 
118  cpl_matrix *self = cpl_matrix_new(3, 3);
119 
120 
121  nrotations = nrotations <= 3 ? nrotations : 3;
122 
123 
124  /*
125  * Initialize the rotation matrix as identity
126  */
127 
128  cpl_matrix_fill_diagonal(self, 1., 0);
129 
130 
131  /*
132  * Create the general rotation matrix by successively left
133  * multiplying the rotation matrices about the given axes.
134  */
135 
136  if (nrotations > 0) {
137 
138  const cxchar *_mode = cx_strlower(cx_strdup(mode));
139 
140  register cxint i = 0;
141 
142  cpl_matrix *R = cpl_matrix_new(3, 3);
143 
144  cxdouble angle[3] = {phi, theta, psi};
145  cxdouble *_R = cpl_matrix_get_data(R);
146 
147 
148  for (i = 0; i < nrotations; ++i) {
149 
150  cxdouble sa = sin(angle[i]);
151  cxdouble ca = cos(angle[i]);
152 
153  cpl_matrix *mt = NULL;
154 
155 
156  cpl_matrix_fill(R, 0.);
157 
158  switch (_mode[i]) {
159  case 'x':
160  {
161  _R[0] = 1.;
162  _R[4] = ca;
163  _R[5] = sa;
164  _R[7] = -sa;
165  _R[8] = ca;
166  break;
167  }
168 
169  case 'y':
170  {
171  _R[0] = ca;
172  _R[2] = -sa;
173  _R[4] = 1.;
174  _R[6] = sa;
175  _R[8] = ca;
176  break;
177  }
178 
179  case 'z':
180  {
181  _R[0] = ca;
182  _R[1] = sa;
183  _R[3] = -sa;
184  _R[4] = ca;
185  _R[8] = 1.;
186  break;
187  }
188 
189  default:
190  {
191  /* Invalid axis designation */
192 
193  cpl_matrix_delete(R);
194  cpl_matrix_delete(self);
195 
196  cx_free((cxchar*)_mode);
197 
198  return NULL;
199  break;
200  }
201 
202  }
203 
204  mt = cpl_matrix_product_create(R, self);
205 
206  cpl_matrix_delete(self);
207  self = mt;
208 
209  }
210 
211  cpl_matrix_delete(R);
212  cx_free((cxchar *)_mode);
213 
214  }
215 
216  return self;
217 
218 }
219 
220 
243 inline static cpl_matrix *
244 _giraffe_precession_matrix(cxdouble epoch0, cxdouble epoch1)
245 {
246 
247  /*
248  * Variable names follow the notation of Simon et al.
249  */
250 
251  const cxdouble sigma0 = 2000.; /* The basic epoch of J2000 */
252  const cxdouble sigmaF = epoch0; /* The arbitrary fixed epoch */
253  const cxdouble sigmaD = epoch1; /* The epoch of the date */
254 
255  /*
256  * Time differences between the given epochs, in units of 1000 years:
257  *
258  * JD(sigmaF) - JD(sigma0) / 365250.
259  * JD(sigmaD) - JD(sigmaF) / 365250.
260  */
261 
262  const cxdouble T = (sigmaF - sigma0) / 1000.;
263  const cxdouble t = (sigmaD - sigmaF) / 1000.;
264 
265 
266  cpl_matrix *mprc = NULL;
267 
268  cxdouble thetaA = RV_DAS2R;
269  cxdouble zetaA = RV_DAS2R;
270  cxdouble zA = RV_DAS2R;
271 
272  cxdouble T2 = T * T;
273  cxdouble T3 = T2 * T;
274  cxdouble T4 = T3 * T;
275  cxdouble T5 = T4 * T;
276 
277  cxdouble t2 = t * t;
278  cxdouble t3 = t2 * t;
279  cxdouble t4 = t3 * t;
280  cxdouble t5 = t4 * t;
281  cxdouble t6 = t5 * t;
282 
283 
284  /*
285  * Compute the Euler angles (in units of radians)
286  *
287  * The constants in the following are the constants given in the example
288  * code prcupd.f of Simon et al. (1994) converted to units of arcsec
289  * (the constants in the example are given in units of 1/10000 arcsec).
290  */
291 
292 
293  thetaA *= (20042.0207 - 85.3131 * T - 0.2111 * T2 + 0.3642 * T3 +
294  0.0008 * T4 -0.0005 * T5) * t +
295  (-42.6566 - 0.2111 * T + 0.5463 * T2 + 0.0017 * T3 -
296  0.0012 * T4) * t2 +
297  (-41.8238 + 0.0359 * T + 0.0027 * T2 - 0.0001 * T3) * t3 +
298  (-0.0731 + 0.0019 * T + 0.0009 * T2) * t4 +
299  (-0.0127 + 0.0011 * T) * t5 +
300  (0.0004) * t6;
301 
302  zetaA *= (23060.9097 + 139.7495 * T - 0.0038 * T2 - 0.5918 * T3 -
303  0.0037 * T4 + 0.0007 * T5) * t +
304  (30.2226 - 0.2523 * T - 0.3840 * T2 - 0.0014 * T3 +
305  0.0007 * T4) * t2 +
306  (18.0183 - 0.1326 * T + 0.0006 * T2 + 0.0005 * T3) * t3 +
307  (-0.0583 - 0.0001 * T + 0.0007 * T2) * t4 +
308  (-0.0285) * t5 +
309  (-0.0002) * t6;
310 
311  zA *= (23060.9097 + 139.7495 * T - 0.0038 * T2 - 0.5918 * T3 -
312  0.0037 * T4 + 0.0007 * T5) * t +
313  (109.5270 + 0.2446 * T - 1.3913 * T2 - 0.0134 * T3 +
314  0.0026 * T4) * t2 +
315  (18.2667 - 1.1400 * T - 0.0173 * T2 + 0.0044 * T3) * t3 +
316  (-0.2821 - 0.0093 * T + 0.0032 * T2) * t4 +
317  (-0.0301 + 0.0006 * T) * t5 +
318  (-0.0001) * t6;
319 
320 
321  /*
322  * Create precession matrix
323  */
324 
325  mprc = _giraffe_create_rotation("zyz", -zetaA, thetaA, -zA);
326 
327  return mprc;
328 
329 }
330 
331 
332 /*
333  * @brief
334  * Compute the mean local sideral time
335  *
336  * @param djd julian date
337  * @param dlong observer's longitude (radians)
338  *
339  * @return
340  * Mean local sideral time (radians).
341  *
342  * @note
343  * Constants taken from the american ephemeris and nautical almanac, 1980)
344  * Translated from the FORTRAN version written by G. Torres (1989)
345  */
346 
347 inline static cxdouble
348 sideral_time(cxdouble djd, cxdouble dlong)
349 {
350 
351  /*
352  * Constants D1,D2,D3 for calculating Greenwich
353  * Mean Sideral Time at 0 UT
354  */
355 
356  const cxdouble d1 = 1.739935934667999;
357  const cxdouble d2 = 6.283319509909095e02;
358  const cxdouble d3 = 6.755878646261384e-06;
359 
360  const cxdouble df = 1.00273790934;
361 
362  cxdouble dut = 0.;
363  cxdouble dt = 0.;
364  cxdouble dst0 = 0.;
365  cxdouble dst = 0.;
366  cxdouble djd0 = floor(djd) + 0.5;
367 
368  if (djd0 > djd) {
369  djd0 -= 1.;
370  }
371 
372  dut = (djd - djd0) * RV_D2PI;
373 
374  dt = (djd0 - dct0) / dcjul;
375  dst0 = d1 + d2 * dt + d3 * dt * dt;
376  dst0 = fmod(dst0, RV_D2PI);
377  dst = df * dut + dst0 - dlong;
378  dst = fmod(dst + RV_D4PI, RV_D2PI);
379 
380  return dst;
381 
382 }
383 
384 
385 /*
386  * @brief
387  * Compute correction from topocentric radial velocity to geocentric.
388  *
389  * @param dlat observer's geodetic latitude (radians)
390  * @param dalt observer's elevation above sea level (meters)
391  * @param dec star's declination (radians) for mean
392  * equator and equinox of date
393  * @param dha star's hour angle (radians)
394  *
395  * @return
396  * Geocentric correction in km/s.
397  *
398  * Calculates the correction required to transform the topocentric radial
399  * velocity of a given star to geocentric. The maximum error of this
400  * routine is not expected to be larger than 0.1 m/s.
401  *
402  * @note
403  * vr = r.w.cos(dec).sin(hour angle), where r = geocentric radius at
404  * observer's latitude and elevation, and w = 2.pi/t, t = length of sideral
405  * day 23 hours 56 minutes and 4 seconds in seconds.
406  * The hour angle is positive east of the meridian.
407  * Other relevant equations from E. W. Woolard & G. M. Clemence (1966),
408  * Spherical Astronomy, p.45 and p.48
409  *
410  * Author:
411  * - G. Torres (1989)
412  * - G. Simond (C translation)
413  */
414 
415 inline static cxdouble
416 geo_correction(cxdouble dlat, cxdouble dalt, cxdouble dec, cxdouble dha)
417 {
418 
419  /*
420  * Earth's equatorial radius (km) (Geod. ref. sys., 1980)
421  */
422 
423  const cxdouble da = 6378.137;
424 
425  /*
426  * Polar flattening (Geod. ref. sys., 1980)
427  */
428 
429  const cxdouble df = 1./298.257222;
430 
431  /*
432  * Angular rotation rate (2.pi/t)
433  */
434 
435  const cxdouble dw = RV_D2PI/86164.;
436 
437 
438  const cxdouble de2 = df * (2.0 - df);
439  const cxdouble dsdlats = sin (dlat) * sin (dlat);
440 
441  cxdouble d1 = 0.;
442  cxdouble d2 = 0.;
443  cxdouble dr0 = 0.;
444  cxdouble dlatg = 0.;
445  cxdouble drh = 0.;
446  cxdouble dvelg = 0.;
447 
448 
449  /*
450  * Calculate geocentric radius dr0 at sea level (km)
451  */
452 
453  d1 = 1.0 - de2 * (2.0 - de2) * dsdlats;
454  d2 = 1.0 - de2 * dsdlats;
455  dr0 = da * sqrt(d1 / d2);
456 
457 
458  /*
459  * Calculate geocentric latitude dlatg
460  */
461 
462  d1 = de2 * sin(2.0 * dlat);
463  d2 = 2.0 * d2;
464  dlatg = dlat - atan(d1 / d2);
465 
466 
467  /*
468  * Calculate geocentric radius drh at elevation dalt (km)
469  */
470 
471  drh = dr0 * cos(dlatg) + (dalt / 1000.) * cos(dlat);
472 
473 
474  /*
475  * Projected component to star at declination = dec and
476  * at hour angle = dha, in units of km/s
477  */
478 
479  dvelg = dw * drh * cos(dec) * sin(dha);
480 
481  return dvelg;
482 
483 }
484 
485 
486 /*
487  * @brief
488  * Compute the heliocentric and barycentric velocity components of the earth.
489  *
490  * @param dje Julian ephermeris date
491  * @param deq Epoch of mean equator and mean equinox of @em hvel
492  * and @em bvel. If @em deq is @c 0, both vectors refer
493  * to the mean equator and equinox of @em dje
494  *
495  * @return
496  * The components (dx/dt, dy/dt, dz/dt, k=1,2,3 A.U./s) of the
497  * heliocentric and barycentric velocity are returned in the
498  * vectors @em hvel and @em bvel respectively.
499  *
500  * Calculates the heliocentric and barycentric velocity components
501  * of the earth. The largest deviations from the JPL-de96 ephemeris
502  * are 42 cm/s for both heliocentric and barycentric velocity
503  * components.
504  *
505  * @note
506  * vr = r.w.cos(dec).sin(hour angle), where r = geocentric radius at
507  * observer's latitude and elevation, and w = 2.pi/t, t = length of sideral
508  * day (sec). The hour angle is positive east of the meridian.
509  * Other relevant equations from E. W. Woolard & G. M. Clemence (1966),
510  * Spherical Astronomy, p.45 and p.48
511  *
512  * Authors:
513  * - P. Stumpff (ibm-version 1979): astron. astrophys. suppl. ser. 41,
514  * 1 (1980)
515  * - M. H.Slovak (vax 11/780 implementation 1986)
516  * - G. Torres (1989)
517  * - G. Simond (C translation)
518  */
519 
520 inline static void
521 earth_velocity(cxdouble dje, cxdouble deq, cxdouble* const hvel,
522  cxdouble* const bvel)
523 {
524 
525 
526  /*
527  * Constants dcfel(i, k) of fast-changing elements
528  *
529  * i = 1 i = 2 i = 3
530  */
531 
532  const cxdouble dcfel[][3] = {
533  {1.7400353e+00, 6.2833195099091e+02, 5.2796e-06},
534  {6.2565836e+00, 6.2830194572674e+02, -2.6180e-06},
535  {4.7199666e+00, 8.3997091449254e+03, -1.9780e-05},
536  {1.9636505e-01, 8.4334662911720e+03, -5.6044e-05},
537  {4.1547339e+00, 5.2993466764997e+01, 5.8845e-06},
538  {4.6524223e+00, 2.1354275911213e+01, 5.6797e-06},
539  {4.2620486e+00, 7.5025342197656e+00, 5.5317e-06},
540  {1.4740694e+00, 3.8377331909193e+00, 5.6093e-06}
541  };
542 
543 
544  /*
545  * Constants dceps and ccsel(i,k) of slowly changing elements
546  *
547  * i = 1 i = 2 i = 3
548  */
549 
550  const cxdouble dceps[3] = {
551  4.093198e-01,
552  -2.271110e-04,
553  -2.860401e-08
554  };
555 
556  const cxdouble ccsel[][3] = {
557  {1.675104e-02, -4.179579e-05, -1.260516e-07},
558  {2.220221e-01, 2.809917e-02, 1.852532e-05},
559  {1.589963e+00, 3.418075e-02, 1.430200e-05},
560  {2.994089e+00, 2.590824e-02, 4.155840e-06},
561  {8.155457e-01, 2.486352e-02, 6.836840e-06},
562  {1.735614e+00, 1.763719e-02, 6.370440e-06},
563  {1.968564e+00, 1.524020e-02, -2.517152e-06},
564  {1.282417e+00, 8.703393e-03, 2.289292e-05},
565  {2.280820e+00, 1.918010e-02, 4.484520e-06},
566  {4.833473e-02, 1.641773e-04, -4.654200e-07},
567  {5.589232e-02, -3.455092e-04, -7.388560e-07},
568  {4.634443e-02, -2.658234e-05, 7.757000e-08},
569  {8.997041e-03, 6.329728e-06, -1.939256e-09},
570  {2.284178e-02, -9.941590e-05, 6.787400e-08},
571  {4.350267e-02, -6.839749e-05, -2.714956e-07},
572  {1.348204e-02, 1.091504e-05, 6.903760e-07},
573  {3.106570e-02, -1.665665e-04, -1.590188e-07}
574  };
575 
576 
577  /*
578  * Constants of the arguments of the short-period perturbations by
579  * the planets: dcargs(i, k)
580  *
581  * i = 1 i = 2
582  */
583 
584  const cxdouble dcargs[][2] = {
585  {5.0974222e+00, -7.8604195454652e+02},
586  {3.9584962e+00, -5.7533848094674e+02},
587  {1.6338070e+00, -1.1506769618935e+03},
588  {2.5487111e+00, -3.9302097727326e+02},
589  {4.9255514e+00, -5.8849265665348e+02},
590  {1.3363463e+00, -5.5076098609303e+02},
591  {1.6072053e+00, -5.2237501616674e+02},
592  {1.3629480e+00, -1.1790629318198e+03},
593  {5.5657014e+00, -1.0977134971135e+03},
594  {5.0708205e+00, -1.5774000881978e+02},
595  {3.9318944e+00, 5.2963464780000e+01},
596  {4.8989497e+00, 3.9809289073258e+01},
597  {1.3097446e+00, 7.7540959633708e+01},
598  {3.5147141e+00, 7.9618578146517e+01},
599  {3.5413158e+00, -5.4868336758022e+02}
600  };
601 
602 
603  /*
604  * Amplitudes ccamps(n, k) of the short-period perturbations
605  *
606  * n = 1 n = 2 n = 3 n = 4 n = 5
607  */
608 
609  const cxdouble ccamps[][5] = {
610  {-2.279594e-5, 1.407414e-5, 8.273188e-6, 1.340565e-5, -2.490817e-7},
611  {-3.494537e-5, 2.860401e-7, 1.289448e-7, 1.627237e-5, -1.823138e-7},
612  { 6.593466e-7, 1.322572e-5, 9.258695e-6, -4.674248e-7, -3.646275e-7},
613  { 1.140767e-5, -2.049792e-5, -4.747930e-6, -2.638763e-6, -1.245408e-7},
614  { 9.516893e-6, -2.748894e-6, -1.319381e-6, -4.549908e-6, -1.864821e-7},
615  { 7.310990e-6, -1.924710e-6, -8.772849e-7, -3.334143e-6, -1.745256e-7},
616  {-2.603449e-6, 7.359472e-6, 3.168357e-6, 1.119056e-6, -1.655307e-7},
617  {-3.228859e-6, 1.308997e-7, 1.013137e-7, 2.403899e-6, -3.736225e-7},
618  { 3.442177e-7, 2.671323e-6, 1.832858e-6, -2.394688e-7, -3.478444e-7},
619  { 8.702406e-6, -8.421214e-6, -1.372341e-6, -1.455234e-6, -4.998479e-8},
620  {-1.488378e-6, -1.251789e-5, 5.226868e-7, -2.049301e-7, 0.0e0},
621  {-8.043059e-6, -2.991300e-6, 1.473654e-7, -3.154542e-7, 0.0e0},
622  { 3.699128e-6, -3.316126e-6, 2.901257e-7, 3.407826e-7, 0.0e0},
623  { 2.550120e-6, -1.241123e-6, 9.901116e-8, 2.210482e-7, 0.0e0},
624  {-6.351059e-7, 2.341650e-6, 1.061492e-6, 2.878231e-7, 0.0e0}
625  };
626 
627 
628  /*
629  * Constants of the secular perturbations in longitude
630  * ccsec3 and ccsec(n,k)
631  * n = 1 n = 2 n = 3
632  */
633 
634  const cxdouble ccsec3 = -7.757020e-08;
635 
636  const cxdouble ccsec[][3] = {
637  {1.289600e-06, 5.550147e-01, 2.076942e+00},
638  {3.102810e-05, 4.035027e+00, 3.525565e-01},
639  {9.124190e-06, 9.990265e-01, 2.622706e+00},
640  {9.793240e-07, 5.508259e+00, 1.559103e+01}
641  };
642 
643 
644  /*
645  * Sideral rate dcsld in longitude, rate ccsgd in mean anomaly
646  */
647 
648  const cxdouble dcsld = 1.990987e-07;
649  const cxdouble ccsgd = 1.990969e-07;
650 
651 
652  /*
653  * Some constants used in the calculation of the
654  * lunar contribution
655  */
656 
657  const cxdouble cckm = 3.122140e-05;
658  const cxdouble ccmld = 2.661699e-06;
659  const cxdouble ccfdi = 2.399485e-07;
660 
661 
662  /*
663  * Constants dcargm(i,k) of the arguments of the
664  * perturbations of the motion of the moon
665  *
666  * i = 1 i = 2
667  */
668 
669  const cxdouble dcargm[][2] = {
670  {5.1679830e+00, 8.3286911095275e+03},
671  {5.4913150e+00, -7.2140632838100e+03},
672  {5.9598530e+00, 1.5542754389685e+04}
673  };
674 
675 
676  /*
677  * Amplitudes ccampm(n,k) of the perturbations of the moon
678  * n = 1 n = 2 n = 3 n = 4
679  */
680 
681  const cxdouble ccampm[][4] = {
682  { 1.097594e-01, 2.896773e-07, 5.450474e-02, 1.438491e-07},
683  {-2.223581e-02, 5.083103e-08, 1.002548e-02, -2.291823e-08},
684  { 1.148966e-02, 5.658888e-08, 8.249439e-03, 4.063015e-08}
685  };
686 
687 
688  /*
689  * ccpamv = a * m * dl/dt (planets)
690  */
691 
692  const cxdouble ccpamv[4] = {
693  8.326827e-11,
694  1.843484e-11,
695  1.988712e-12,
696  1.881276e-12
697  };
698 
699 
700  /*
701  * dc1mme = 1 - mass(earth + moon)
702  */
703 
704  const cxdouble dc1mme = 0.99999696;
705 
706 
707  register cxint k = 0;
708  register cxint n = 0;
709 
710  cxint ideq = 0;
711 
712  cxdouble a = 0.;
713  cxdouble b = 0.;
714  cxdouble f = 0.;
715  cxdouble dt = 0.;
716  cxdouble t = 0.;
717  cxdouble tl = 0.;
718  cxdouble dtsq = 0.;
719  cxdouble tsq = 0.;
720  cxdouble dml = 0.;
721  cxdouble dlocal = 0.;
722  cxdouble deps = 0.;
723  cxdouble pertl = 0.;
724  cxdouble pertld = 0.;
725  cxdouble pertr = 0.;
726  cxdouble pertrd = 0.;
727  cxdouble pertp = 0.;
728  cxdouble pertpd = 0.;
729  cxdouble sina = 0.;
730  cxdouble cosa = 0.;
731  cxdouble twoe = 0.;
732  cxdouble twog = 0.;
733  cxdouble param = 0.;
734  cxdouble dparam = 0.;
735  cxdouble dpsi = 0.;
736  cxdouble phi = 0.;
737  cxdouble phid = 0.;
738  cxdouble psid = 0.;
739  cxdouble sin_f = 0.;
740  cxdouble cos_f = 0.;
741  cxdouble esq = 0.;
742  cxdouble d1pdro = 0.;
743  cxdouble drd = 0.;
744  cxdouble drld = 0.;
745  cxdouble dsinls = 0.;
746  cxdouble dcosls = 0.;
747  cxdouble dtl = 0.;
748  cxdouble dxhd = 0.;
749  cxdouble dyhd = 0.;
750  cxdouble dzhd = 0.;
751  cxdouble sinlm = 0.;
752  cxdouble coslm = 0.;
753  cxdouble sigma = 0.;
754  cxdouble plon = 0.;
755  cxdouble pomg = 0.;
756  cxdouble dxbd = 0.;
757  cxdouble dybd = 0.;
758  cxdouble dzbd = 0.;
759  cxdouble pecc = 0.;
760  cxdouble dcosep = 0.;
761  cxdouble dsinep = 0.;
762  cxdouble dyahd = 0.;
763  cxdouble dzahd = 0.;
764  cxdouble dyabd = 0.;
765  cxdouble dzabd = 0.;
766  cxdouble sn[4] = {0., 0., 0., 0.};
767  cxdouble sinlp[4] = {0., 0., 0., 0.};
768  cxdouble coslp[4] = {0., 0., 0., 0.};
769  cxdouble forbel[7] = {0., 0., 0., 0., 0., 0., 0.};
770  cxdouble sorbel[17];
771 
772 
773  memset(sorbel, 0, sizeof sorbel);
774 
775 
776  /*
777  * Control-parameter ideq, and time-arguments
778  */
779 
780  ideq = (cxint)deq;
781  dt = (dje - dct0) / dcjul;
782  t = dt;
783  dtsq = dt * dt;
784  tsq = dtsq;
785 
786 
787  /*
788  * Values of all elements for the instant dje
789  */
790 
791  for (k = 0; k < 8; k++) {
792 
793  dlocal = fmod(dcfel[k][0] + dt * dcfel[k][1] + dtsq * dcfel[k][2],
794  RV_D2PI);
795 
796  if (k == 0) {
797  dml = dlocal;
798  }
799 
800  if (k != 0) {
801  forbel[k - 1] = dlocal;
802  }
803 
804  }
805 
806  deps = fmod(dceps[0] + dt * dceps[1] + dtsq * dceps[2], RV_D2PI);
807 
808  for (k = 0; k < 17; k++) {
809 
810  sorbel[k] = fmod(ccsel[k][0] + t * ccsel[k][1] + tsq * ccsel[k][2],
811  RV_D2PI);
812 
813  }
814 
815 
816  /*
817  * Secular perturbations in longitude
818  */
819 
820  for (k = 0; k < 4; k++) {
821 
822  a = fmod(ccsec[k][1] + t * ccsec[k][2], RV_D2PI);
823  sn[k] = sin(a);
824 
825  }
826 
827 
828  /*
829  * Periodic perturbations of the earth-moon barycenter
830  */
831 
832  pertl = ccsec[0][0] * sn[0] + ccsec[1][0] * sn[1] +
833  (ccsec[2][0] + t * ccsec3) * sn[2] + ccsec[3][0] * sn[3];
834 
835  pertld = 0.;
836  pertr = 0.;
837  pertrd = 0.;
838 
839  for (k = 0; k < 15; k++) {
840 
841  a = fmod(dcargs[k][0] + dt * dcargs[k][1], RV_D2PI);
842  cosa = cos (a);
843  sina = sin (a);
844  pertl += (ccamps[k][0] * cosa + ccamps[k][1] * sina);
845  pertr += (ccamps[k][2] * cosa + ccamps[k][3] * sina);
846 
847  if (k >= 10) {
848  continue;
849  }
850 
851  pertld += ((ccamps[k][1] * cosa - ccamps[k][0] * sina) * ccamps[k][4]);
852  pertrd += ((ccamps[k][3] * cosa - ccamps[k][2] * sina) * ccamps[k][4]);
853 
854  }
855 
856 
857  /*
858  * Elliptic part of the motion of the earth-moon barycenter
859  */
860 
861  esq = sorbel[0] * sorbel[0];
862  dparam = 1. - esq;
863  param = dparam;
864  twoe = sorbel[0] + sorbel[0];
865  twog = forbel[0] + forbel[0];
866  phi = twoe * ((1.0 - esq * (1.0 / 8.0)) * sin (forbel[0]) +
867  sorbel[0] * (5.0 / 8.0) * sin (twog) +
868  esq * 0.5416667 * sin (forbel[0] + twog));
869  f = forbel[0] + phi;
870  sin_f = sin(f);
871  cos_f = cos(f);
872  dpsi = dparam / (1. + sorbel[0] * cos_f);
873  phid = twoe * ccsgd * ((1.0 + esq * 1.50) * cos_f +
874  sorbel[0] * (1.250 - sin_f * sin_f * 0.50));
875  psid = ccsgd * sorbel[0] * sin_f / sqrt(param);
876 
877 
878  /*
879  * Perturbed heliocentric motion of the earth-moon barycenter
880  */
881 
882  d1pdro = 1. + pertr;
883  drd = d1pdro * (psid + dpsi * pertrd);
884  drld = d1pdro * dpsi * (dcsld + phid + pertld);
885  dtl = fmod(dml + phi + pertl, RV_D2PI);
886  dsinls = sin(dtl);
887  dcosls = cos(dtl);
888  dxhd = drd * dcosls - drld * dsinls;
889  dyhd = drd * dsinls + drld * dcosls;
890 
891 
892  /*
893  * Influence of eccentricity, evection and variation of
894  * the geocentric motion of the moon
895  */
896 
897  pertl = 0.;
898  pertld = 0.;
899  pertp = 0.;
900  pertpd = 0.;
901 
902  for (k = 0; k < 3; k++) {
903 
904  a = fmod(dcargm[k][0] + dt * dcargm[k][1], RV_D2PI);
905  sina = sin(a);
906  cosa = cos(a);
907  pertl += ccampm[k][0] * sina;
908  pertld += ccampm[k][1] * cosa;
909  pertp += ccampm[k][2] * cosa;
910  pertpd -= ccampm[k][3] * sina;
911 
912  }
913 
914 
915  /*
916  * Heliocentric motion of the earth
917  */
918 
919  tl = forbel[1] + pertl;
920  sinlm = sin(tl);
921  coslm = cos(tl);
922  sigma = cckm / (1. + pertp);
923  a = sigma * (ccmld + pertld);
924  b = sigma * pertpd;
925  dxhd = dxhd + a * sinlm + b * coslm;
926  dyhd = dyhd - a * coslm + b * sinlm;
927  dzhd = -sigma * ccfdi * cos(forbel[2]);
928 
929 
930  /*
931  * Barycentric motion of the earth
932  */
933 
934  dxbd = dxhd * dc1mme;
935  dybd = dyhd * dc1mme;
936  dzbd = dzhd * dc1mme;
937 
938  for (k = 0; k < 4; k++) {
939 
940  plon = forbel[k + 3];
941  pomg = sorbel[k + 1];
942  pecc = sorbel[k + 9];
943  tl = fmod(plon + 2.0 * pecc * sin (plon - pomg), RV_D2PI);
944  sinlp[k] = sin(tl);
945  coslp[k] = cos(tl);
946  dxbd = dxbd + ccpamv[k] * (sinlp[k] + pecc * sin(pomg));
947  dybd = dybd - ccpamv[k] * (coslp[k] + pecc * cos(pomg));
948  dzbd = dzbd - ccpamv[k] * sorbel[k + 13] * cos(plon - sorbel[k + 5]);
949 
950  }
951 
952 
953  /*
954  * Transition to mean equator of date
955  */
956 
957  dcosep = cos(deps);
958  dsinep = sin(deps);
959  dyahd = dcosep * dyhd - dsinep * dzhd;
960  dzahd = dsinep * dyhd + dcosep * dzhd;
961  dyabd = dcosep * dybd - dsinep * dzbd;
962  dzabd = dsinep * dybd + dcosep * dzbd;
963 
964  if (ideq == 0) {
965 
966  hvel[0] = dxhd;
967  hvel[1] = dyahd;
968  hvel[2] = dzahd;
969 
970  bvel[0] = dxbd;
971  bvel[1] = dyabd;
972  bvel[2] = dzabd;
973 
974  }
975  else {
976 
977  /*
978  * General precession from epoch dje to deq
979  */
980 
981  cxdouble deqdat = (dje - dct0 - dcbes) / dctrop + dc1900;
982 
983  cpl_matrix* prec = _giraffe_precession_matrix(deqdat, deq);
984 
985 
986  for (n = 0; n < 3; n++) {
987 
988  hvel[n] =
989  dxhd * cpl_matrix_get(prec, 0, n) +
990  dyahd * cpl_matrix_get(prec, 1, n) +
991  dzahd * cpl_matrix_get(prec, 2, n);
992 
993  bvel[n] =
994  dxbd * cpl_matrix_get(prec, 0, n) +
995  dyabd * cpl_matrix_get(prec, 1, n) +
996  dzabd * cpl_matrix_get(prec, 2, n);
997  }
998 
999  cpl_matrix_delete(prec);
1000 
1001  }
1002 
1003  return;
1004 
1005 }
1006 
1007 
1008 /*
1009  * Public functions
1010  */
1011 
1045 void
1046 giraffe_rvcorrection_compute(GiRvCorrection* rv,
1047  cxdouble jdate, cxdouble longitude,
1048  cxdouble latitude, cxdouble elevation,
1049  cxdouble ra, cxdouble dec,
1050  cxdouble equinox)
1051 {
1052 
1053  cxint i = 0;
1054 
1055  const cxdouble aukm = 1.4959787e08; /* 1 UA = 149 597 870 km */
1056 
1057  cxdouble eqt = 0.;
1058  cxdouble ha = 0.;
1059  cxdouble ra2 = 0.;
1060  cxdouble dec2 = 0.;
1061  cxdouble dc[3] = {0., 0., 0.};
1062  cxdouble dcc[3] = {0., 0., 0.};
1063  cxdouble hv[3] = {0., 0., 0.};
1064  cxdouble bv[3] = {0., 0., 0.};
1065  cxdouble _long = longitude * RV_DD2R;
1066  cxdouble _lat = latitude * RV_DD2R;
1067  cxdouble _ra = ra * 15.0 * RV_DD2R;
1068  cxdouble _dec = dec * RV_DD2R;
1069  cxdouble st = sideral_time(jdate, _long);
1070 
1071  cpl_matrix* precession = NULL;
1072 
1073 
1074  /*
1075  * Precess r.a. and dec. to mean equator and equinox of date
1076  */
1077 
1078  eqt = (jdate - dct0 - dcbes) / dctrop + dc1900;
1079 
1080  dc[0] = cos(_ra) * cos(_dec);
1081  dc[1] = sin(_ra) * cos(_dec);
1082  dc[2] = sin(_dec);
1083 
1084  precession = _giraffe_precession_matrix(equinox, eqt);
1085 
1086  for (i = 0; i < 3; ++i) {
1087 
1088  dcc[i] =
1089  dc[0] * cpl_matrix_get(precession, i, 0) +
1090  dc[1] * cpl_matrix_get(precession, i, 1) +
1091  dc[2] * cpl_matrix_get(precession, i, 2);
1092 
1093  }
1094 
1095  cpl_matrix_delete(precession);
1096  precession = NULL;
1097 
1098 
1099  if (dcc[0] != 0.) {
1100 
1101  cxdouble darg = dcc[1] / dcc[0];
1102 
1103 
1104  ra2 = atan(darg);
1105 
1106  if (dcc[0] < 0.) {
1107  ra2 += RV_DPI;
1108  }
1109  else {
1110  if (dcc[1] < 0.) {
1111  ra2 += RV_D2PI;
1112  }
1113  }
1114 
1115  }
1116  else {
1117 
1118  if (dcc[1] > 0.) {
1119  ra2 = RV_DPIBY2;
1120  }
1121  else {
1122  ra2 = 1.5 * RV_DPI;
1123  }
1124 
1125  }
1126 
1127  dec2 = asin(dcc[2]);
1128 
1129 
1130  /*
1131  * Calculate hour angle = local sideral time - r.a.
1132  */
1133 
1134  ha = st - ra2;
1135 
1136 
1137  /*
1138  * Calculate observer's geocentric velocity
1139  * (elevation assumed to be zero)
1140  */
1141 
1142  rv->gc = geo_correction(_lat, elevation, dec2, -ha);
1143 
1144 
1145  /*
1146  * Calculate components of earth's barycentric velocity,
1147  * bvel(i), i=1,2,3 in units of a.u./s
1148  */
1149 
1150  earth_velocity (jdate, eqt, hv, bv);
1151 
1152 
1153  /*
1154  * Project barycentric velocity to the direction of the
1155  * star, and convert to km/s
1156  */
1157 
1158  rv->bc = 0.;
1159  rv->hc = 0.;
1160 
1161  for (i = 0; i < 3; ++i) {
1162  rv->bc += bv[i] * dcc[i] * aukm;
1163  rv->hc += hv[i] * dcc[i] * aukm;
1164  }
1165 
1166  return;
1167 
1168 }
void giraffe_rvcorrection_compute(GiRvCorrection *rv, cxdouble jdate, cxdouble longitude, cxdouble latitude, cxdouble elevation, cxdouble ra, cxdouble dec, cxdouble equinox)
Compute heliocentric, barycentric and geocentric correction.

This file is part of the GIRAFFE Pipeline Reference Manual 2.14.
Documentation copyright © 2002-2006 European Southern Observatory.
Generated on Wed Mar 11 2015 13:19:42 by doxygen 1.8.9.1 written by Dimitri van Heesch, © 1997-2004