GIRAFFE Pipeline Reference Manual

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