X-shooter Pipeline Reference Manual 3.8.15
xsh_baryvel.c
Go to the documentation of this file.
1/* *
2 * This file is part of the ESO UVES Pipeline *
3 * Copyright (C) 2004,2005 European Southern Observatory *
4 * *
5 * This library 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, 51 Franklin St, Fifth Floor, Boston, MA 02111-1307 USA *
18 * */
19
20/*
21 * $Author: amodigli $
22 * $Date: 2012-12-15 19:02:40 $
23 * $Revision: 1.5 $
24 * $Name: not supported by cvs2svn $
25 */
26
27#ifdef HAVE_CONFIG_H
28# include <config.h>
29#endif
30
31/*----------------------------------------------------------------------------*/
44/*----------------------------------------------------------------------------*/
47/*-----------------------------------------------------------------------------
48 Includes
49 -----------------------------------------------------------------------------*/
50
51#include <xsh_baryvel.h>
52#include <xsh_parameters.h>
53
54#include <xsh_pfits.h>
55#include <xsh_utils.h>
56#include <xsh_error.h>
57#include <xsh_msg.h>
58
59#include <irplib_sdp_spectrum.h>
60
61#include <cpl.h>
62
63#include <math.h>
64
65/*-----------------------------------------------------------------------------
66 Local functions
67 -----------------------------------------------------------------------------*/
68static void deg2dms(double in_val,
69 double *degs,
70 double *minutes,
71 double *seconds);
72
73static void deg2hms(double in_val,
74 double *hour,
75 double *min,
76 double *sec);
77
78static void compxy(double inputr[19], char inputc[4],
79 double outputr[4],
80 double utr, double mod_juldat);
81
82static void barvel(double DJE, double DEQ,
83 double DVELH[4], double DVELB[4]);
84
85
86/*----------------------------------------------------------------------------*/
93/*----------------------------------------------------------------------------*/
94void
95xsh_baryvel(const cpl_propertylist *raw_header,
96 double *bary_corr,
97 double *helio_corr)
98{
99
100 double outputr[4];
101
102//inputc(1:3) = "+++"
103 char inputc[] = "X+++"; /* 0th index not used */
104
105//define/local rneg/r/1/1 1.0
106 double rneg = 1.0;
107
108// write/keyw inputr/r/1/18 0.0 all
109 double inputr[19]; /* Do not use the zeroth element */
110
111
112/*
113 qc_ra = m$value({p1},O_POS(1))
114 qc_dec = m$value({p1},O_POS(2))
115 qc_geolat = m$value({p1},{h_geolat})
116 qc_geolon = m$value({p1},{h_geolon})
117 qc_obs_time = m$value({p1},O_TIME(7)) !using an image as input it take the
118 !date from the descriptor O_TIME(1,2,3)
119 !and the UT from O_TIME(5)
120*/
121 double qc_ra;
122 double qc_dec;
123 double qc_geolat;
124 double qc_geolon;
125
126 double utr;
127 double mod_juldat;
128
129 double ra_hour, ra_min, ra_sec;
130 double dec_deg, dec_min, dec_sec;
131 double lat_deg, lat_min, lat_sec;
132 double lon_deg, lon_min, lon_sec;
133
134 if(cpl_propertylist_has(raw_header,"RA")) {
135 check_msg( qc_ra = xsh_pfits_get_ra(raw_header), /* in degrees */
136 "Error getting object right ascension");
137 } else {
138 xsh_msg_warning("RA FITS keyword not present in header, barycor not computed");
139 return;
140 }
141
142 if(cpl_propertylist_has(raw_header,"DEC")) {
143 check_msg( qc_dec = xsh_pfits_get_dec(raw_header),
144 "Error getting object declination");
145 } else {
146 xsh_msg_warning("DEC FITS keyword not present in header, barycor not computed");
147 return;
148 }
149 check_msg( qc_geolat = xsh_pfits_get_geolat(raw_header),
150 "Error getting telescope latitude");
151 check_msg( qc_geolon = xsh_pfits_get_geolon(raw_header),
152 "Error getting telescope longitude");
153
154 /* double qc_obs_time = xsh_pfits_get_exptime(raw_header); Not used! */
155
156 check_msg( utr = xsh_pfits_get_utc(raw_header),
157 "Error reading UTC");
158 check_msg( mod_juldat = xsh_pfits_get_mjdobs(raw_header),
159 "Error julian date");
160
161 deg2hms(qc_ra, &ra_hour, &ra_min, &ra_sec);
162 deg2dms(qc_dec, &dec_deg, &dec_min, &dec_sec);
163 deg2dms(qc_geolat, &lat_deg, &lat_min, &lat_sec);
164 deg2dms(qc_geolon, &lon_deg, &lon_min, &lon_sec);
165
166// inputr(1) = m$value({p1},o_time(1))
167// inputr(2) = m$value({p1},o_time(2))
168// inputr(3) = m$value({p1},o_time(3))
169// inputr(4) = m$value({p1},o_time(5)) !UT in real hours
170// inputr[1] = year; not needed, pass mjd instead
171// inputr[2] = month;
172// inputr[3] = day;
173// inputr[4] = ut_hour; not needed, pass ut instead
174// inputr[5] = ut_min;
175// inputr[6] = ut_sec;
176
177// write/keyw inputr/r/7/3 {p4}
178 inputr[7] = lon_deg;
179 inputr[8] = lon_min;
180 inputr[9] = lon_sec;
181
182 //rneg = (inputr(7)*3600.)+(inputr(8)*60.)+inputr(9)
183 rneg = (inputr[7]*3600.)+(inputr[8]*60.)+inputr[9];
184 //inputc(1:1) = p4(1:1)
185 inputc[1] = (lon_deg >= 0) ? '+' : '-';
186 //if rneg .lt. 0.0 inputc(1:1) = "-"
187 if (rneg < 0) inputc[1] = '-';
188
189// write/keyw inputr/r/10/3 {p5},0,0
190 inputr[10] = lat_deg;
191 inputr[11] = lat_min;
192 inputr[12] = lat_sec;
193
194// rneg = (inputr(10)*3600.)+(inputr(11)*60.)+inputr(12)
195 rneg = (inputr[10]*3600.)+(inputr[11]*60.)+inputr[12];
196// inputc(2:2) = p5(1:1)
197 inputc[2] = (lat_deg >= 0) ? '+' : '-';
198// if rneg .lt. 0.0 inputc(2:2) = "-"
199 if (rneg < 0) inputc[2] = '-';
200
201// write/keyw inputr/r/13/3 {p2},0,0
202 inputr[13] = ra_hour;
203 inputr[14] = ra_min;
204 inputr[15] = ra_sec;
205
206// write/keyw inputr/r/16/3 {p3},0,0
207 inputr[16] = dec_deg;
208 inputr[17] = dec_min;
209 inputr[18] = dec_sec;
210
211// inputc(3:3) = p3(1:1)
212 inputc[3] = (dec_deg >= 0) ? '+' : '-';
213// rneg = (inputr(16)*3600.)+(inputr(17)*60.)+inputr(18)
214 rneg = (inputr[16]*3600.)+(inputr[17]*60.)+inputr[18];
215// if rneg .lt. 0.0 inputc(3:3) = "-"
216 if (rneg < 0) inputc[3] = '-';
217
218
219//C INPUTR/R/1/3 date: year,month,day
220//C INPUTR/R/4/3 universal time: hour,min,sec
221//C INPUTR/R/7/3 EAST longitude of observatory: degree,min,sec !! NOTE
222//C INPUTR/R/10/3 latitude of observatory: degree,min,sec
223//C INPUTR/R/13/3 right ascension: hour,min,sec
224//C INPUTR/R/16/3 declination: degree,min,sec
225
226 //write/keyw action BA !indicate barycorr stuff
227 //run MID_EXE:COMPXY !compute the corrections
228
229 compxy(inputr, inputc, outputr, utr, mod_juldat);
230
231// set/format f14.6,g24.12
232// xsh_msg_debug(" Barycentric correction time: {outputd(1)} day");
233// xsh_msg_debug(" Heliocentric correction time: {outputd(2)} day");
234// xsh_msg_debug(" ");
235 xsh_msg_debug(" Total barycentric RV correction: %f km/s", outputr[1]);
236 xsh_msg_debug(" Total heliocentric RV correction: %f km/s", outputr[2]);
237 xsh_msg_debug(" (incl. diurnal RV correction of %f km/s)", outputr[3]);
238// xsh_msg_debug(" ");
239// xsh_msg_debug("Descriptor O_TIME of image {p1} used for date and UT.");
240
241 *bary_corr = outputr[1];
242 *helio_corr = outputr[2];
243
244 cleanup:
245 return;
246}
247
248
249/*----------------------------------------------------------------------------*/
271/*----------------------------------------------------------------------------*/
272static void
273compxy(double inputr[19], char inputc[4],
274 double outputr[4],
275 double utr, double mod_juldat)
276{
277
278// INTEGER IAV,STAT,KUN(1),KNUL,N
279// INTEGER MADRID
280//
281// DOUBLE PRECISION UTR,STR,T0,DL,THETA0,PE,ST0HG,STG,GAST,R1
282 double STR;
283
284// double utr Not used. Use FITS header value instead
285 double t0, dl, theta0, pe, st0hg, stg;
286// DOUBLE PRECISION JD,JD0H,JD00,ZERO
287 double jd, jd0h;
288// DOUBLE PRECISION DCORB(3),DCORH(3),DVELB(3),DVELH(3)
289 double dvelb[4], dvelh[4];
290// DOUBLE PRECISION ALP,BCT,BEOV,BERV,DEL,EDV
291 double alp, del, beov, berv, EDV;
292// DOUBLE PRECISION HAR,HCT,HEOV,HERV,PHI,PI
293 double HAR, phi, heov, herv;
294// DOUBLE PRECISION EQX0,EQX1
295// DOUBLE PRECISION A0R,A1R,D0R,D1R
296// DOUBLE PRECISION DSMALL,DTEMP(3)
297//
298// REAL DATE0(3),DATE1(3),DATE00(3),A0(3),A1(3),D0(3),D1(3)
299// REAL DATE(3),UT(3),OLONG(3),ST(3)
300// double ut[4];
301// REAL OLAT(3),ALPHA(3),DELTA(3)
302// REAL RBUF(20)
303 double *rbuf;
304//
305// CHARACTER ACTIO*2,SIGNS*3,INPSGN*3
306 char inpsgn[4];
307//
308// COMMON /VMR/MADRID(1)
309//
310// DATA PI /3.1415926535897928D0/
311// DATA DSMALL /1.D-38/
312
313
314 double *olong, *olat, *alpha, *delta;
315
316//1000 SIGNS = '+++'
317 char signs[] = "+++";
318
319// CALL STKRDR('INPUTR',1,20,IAV,RBUF,KUN,KNUL,STAT)
320 rbuf = inputr;
321// CALL STKRDC('INPUTC',1,1,3,IAV,INPSGN,KUN,KNUL,STAT)
322 inpsgn[1] = inputc[1];
323 inpsgn[2] = inputc[2];
324 inpsgn[3] = inputc[3];
325
326
327// EQUIVALENCE (RBUF(1),DATE(1)),(RBUF(7),OLONG(1))
328// double *date = rbuf + 1 - 1; Not used, use the explicitly passed MJD instead
329 olong = rbuf + 7 - 1;
330// EQUIVALENCE (RBUF(10),OLAT(1)),(RBUF(13),ALPHA(1))
331 olat = rbuf + 10 - 1;
332 alpha = rbuf + 13 - 1;
333// EQUIVALENCE (RBUF(16),DELTA(1))
334 delta = rbuf + 16 - 1;
335
336
337
338// DO 1100 N=1,3
339// UT(N) = RBUF(N+3)
340//1100 CONTINUE
341// for (n = 1; n <= 3; n++)
342// {
343// ut[n] = rbuf[n+3];
344// }
345
346// ... convert UT to real hours, calculate Julian date
347
348// UTR = UT(1)+UT(2)/60.D0+UT(3)/3600.D0
349// utr = ut[1]+ut[2]/60. +ut[3]/3600.;
350
351 /* We know this one already but convert seconds -> hours */
352 utr /= 3600;
353
354// CALL JULDAT(DATE,UTR,JD)
355 jd = mod_juldat + 2400000.5;
356
357// ... likewise convert longitude and latitude of observatory to real hours
358// ... and degrees, respectively; take care of signs
359// ... NOTE: east longitude is assumed for input !!
360
361// IF ((OLONG(1).LT.0.0) .OR. (OLONG(2).LT.0.0) .OR.
362// + (OLONG(3).LT.0.0) .OR. (INPSGN(1:1).EQ.'-')) THEN
363 if (olong[1] < 0 || olong[2] < 0 ||
364 olong[3] < 0 || inpsgn[1] == '-') {
365// SIGNS(1:1) = '-'
366 signs[1] = '-';
367// OLONG(1) = ABS(OLONG(1))
368// OLONG(2) = ABS(OLONG(2))
369// OLONG(3) = ABS(OLONG(3))
370 olong[1] = fabs(olong[1]);
371 olong[2] = fabs(olong[2]);
372 olong[3] = fabs(olong[3]);
373// ENDIF
374 }
375
376// DL = OLONG(1)+OLONG(2)/60.D0+OLONG(3)/3600.D0
377 dl = olong[1]+olong[2]/60. +olong[3]/3600.;
378
379// IF (SIGNS(1:1).EQ.'-') DL = -DL ! negative longitude
380 if (signs[1] == '-') dl = -dl;
381
382// DL = -DL*24.D0/360.D0 ! convert back to west longitude
383 dl = -dl*24. /360.;
384
385// IF ((OLAT(1).LT.0.0) .OR. (OLAT(2).LT.0.0) .OR.
386// + (OLAT(3).LT.0.0) .OR. (INPSGN(2:2).EQ.'-')) THEN
387 if (olat[1] < 0 || olat[2] < 0 ||
388 olat[3] < 0 || inpsgn[2] == '-') {
389// SIGNS(2:2) = '-'
390 signs[2] = '-';
391
392// OLAT(1) = ABS(OLAT(1))
393// OLAT(2) = ABS(OLAT(2))
394// OLAT(3) = ABS(OLAT(3))
395 olat[1] = fabs(olat[1]);
396 olat[2] = fabs(olat[2]);
397 olat[3] = fabs(olat[3]);
398// ENDIF
399 }
400
401// PHI = OLAT(1)+OLAT(2)/60.D0+OLAT(3)/3600.D0
402 phi = olat[1]+olat[2]/60. +olat[3]/3600.;
403
404// IF (SIGNS(2:2).EQ.'-') PHI = -PHI ! negative latitude
405 if (signs[2] == '-') phi = -phi;
406
407// PHI = PHI*PI/180.D0
408 phi = phi*M_PI/180. ;
409
410// ... convert right ascension and declination to real radians
411
412// ALP = (ALPHA(1)*3600D0+ALPHA(2)*60D0+ALPHA(3))*PI /(12.D0*3600.D0)
413 alp = (alpha[1]*3600. +alpha[2]*60. +alpha[3])*M_PI/(12. *3600. );
414
415// IF ((DELTA(1).LT.0.0) .OR. (DELTA(2).LT.0.0) .OR.
416// + (DELTA(3).LT.0.0) .OR. (INPSGN(3:3).EQ.'-')) THEN
417 if (delta[1] < 0 || delta[2] < 0 ||
418 delta[3] < 0 || inpsgn[3] == '-') {
419// SIGNS(3:3) = '-'
420 signs[3] = '-';
421// DELTA(1) = ABS(DELTA(1))
422// DELTA(2) = ABS(DELTA(2))
423// DELTA(3) = ABS(DELTA(3))
424 delta[1] = fabs(delta[1]);
425 delta[2] = fabs(delta[2]);
426 delta[3] = fabs(delta[3]);
427// ENDIF
428 }
429
430// DEL = (DELTA(1)*3600.D0 + DELTA(2)*60.D0 + DELTA(3))
431// + * PI/(3600.D0*180.D0)
432 del = (delta[1]*3600.0 + delta[2]*60. + delta[3])
433 * M_PI/(3600. *180. );
434
435
436// IF (SIGNS(3:3).EQ.'-') DEL = -DEL ! negative declination
437 if (signs[3] == '-') del = - del;
438
439// ... calculate earth's orbital velocity in rectangular coordinates X,Y,Z
440// ... for both heliocentric and barycentric frames (DVELH, DVELB)
441// ... Note that setting the second argument of BARVEL to zero as done below
442// ... means that the input coordinates will not be corrected for precession.
443
444// CALL BARVEL(JD,0.0D0,DVELH,DVELB)
445 barvel(jd, 0.0, dvelh, dvelb);
446
447// ... with the rectangular velocity components known, the respective projections
448// ... HEOV and BEOV on a given line of sight (ALP,DEL) can be determined:
449
450// ... REFERENCE: THE ASTRONOMICAL ALMANAC 1982 PAGE:B17
451
452// BEOV=DVELB(1)*DCOS(ALP)*DCOS(DEL)+
453// 1 DVELB(2)*DSIN(ALP)*DCOS(DEL)+
454// 2 DVELB(3)*DSIN(DEL)
455 beov =
456 dvelb[1]*cos(alp)*cos(del)+
457 dvelb[2]*sin(alp)*cos(del)+
458 dvelb[3]*sin(del);
459
460// HEOV=DVELH(1)*DCOS(ALP)*DCOS(DEL)+
461// 1 DVELH(2)*DSIN(ALP)*DCOS(DEL)+
462// 2 DVELH(3)*DSIN(DEL)
463 heov =
464 dvelh[1]*cos(alp)*cos(del)+
465 dvelh[2]*sin(alp)*cos(del)+
466 dvelh[3]*sin(del);
467
468
469// ... For determination also of the contribution due to the diurnal rotation of
470// ... the earth (EDV), the hour angle (HAR) is needed at which the observation
471// ... was made which requires conversion of UT to sidereal time (ST).
472
473// ... Therefore, first compute ST at 0 hours UT (ST0HG)
474
475// ... REFERENCE : MEEUS J.,1980,ASTRONOMICAL FORMULAE FOR CALCULATORS
476
477// CALL JULDAT(DATE,ZERO,JD0H)
478 jd0h = jd - (utr/24.0);
479
480// T0=(JD0H-2415020.D0)/36525.D0
481 t0 = (jd0h-2415020. )/36525. ;
482
483// THETA0=0.276919398D0+100.0021359D0*T0+0.000001075D0*T0*T0
484 theta0 = 0.276919398 +100.0021359 *t0+0.000001075 *t0*t0 ;
485
486// PE=DINT(THETA0)
487 pe = (int) theta0;
488
489// THETA0=THETA0-PE
490 theta0 = theta0 - pe;
491
492// ST0HG=THETA0*24.D0
493 st0hg = theta0*24. ;
494
495// ... now do the conversion UT -> ST (MEAN SIDEREAL TIME)
496
497// ... REFERENCE : THE ASTRONOMICAL ALMANAC 1983, P B7
498// ... IN 1983: 1 MEAN SOLAR DAY = 1.00273790931 MEAN SIDEREAL DAYS
499// ... ST WITHOUT EQUATION OF EQUINOXES CORRECTION => ACCURACY +/- 1 SEC
500//
501// STG=ST0HG+UTR*1.00273790931D0
502 stg = st0hg+utr*1.00273790931 ;
503
504// IF (STG.LT.DL) STG=STG+24.D0
505 if (stg < dl) stg = stg +24. ;
506
507// STR=STG-DL
508 STR = stg-dl;
509
510// IF (STR.GE.24.D0) STR=STR-24.D0
511 if (STR >= 24. ) STR = STR-24. ;
512
513// STR = STR*PI/12.D0 ! ST in radians
514 STR = STR*M_PI/12. ;
515
516// HAR=STR-ALP ! hour angle of observation
517 HAR = STR-alp;
518
519// EDV=-0.4654D0*DSIN(HAR)*DCOS(DEL)*DCOS(PHI)
520 EDV = -0.4654 * sin(HAR)* cos(del)* cos(phi);
521
522// ... the total correction (in km/s) is the sum of orbital and diurnal components
523
524// HERV=HEOV+EDV
525 herv=heov+EDV;
526// BERV=BEOV+EDV
527 berv=beov+EDV;
528
529 /* The following is not needed. Do not translate */
530
531#if 0
532// ... Calculation of the barycentric and heliocentric correction times
533// ... (BCT and HCT) requires knowledge of the earth's position in its
534// ... orbit. Subroutine BARCOR returns the rectangular barycentric (DCORB)
535// ... and heliocentric (DCORH) coordinates.
536
537// CALL BARCOR(DCORH,DCORB)
538
539// ... from this, the correction times (in days) can be determined:
540// ... (REFERENCE: THE ASTRONOMICAL ALMANAC 1982 PAGE:B16)
541
542// BCT=+0.0057756D0*(DCORB(1)*DCOS(ALP)*DCOS(DEL)+
543// 1 DCORB(2)*DSIN(ALP)*DCOS(DEL)+
544// 2 DCORB(3)* DSIN(DEL))
545// HCT=+0.0057756D0*(DCORH(1)*DCOS(ALP)*DCOS(DEL)+
546// 1 DCORH(2)*DSIN(ALP)*DCOS(DEL)+
547// 2 DCORH(3)* DSIN(DEL))
548
549//... write results to keywords
550
551// CALL STKWRD('OUTPUTD',BCT,1,1,KUN,STAT) ! barycentric correction time
552// CALL STKWRD('OUTPUTD',HCT,2,1,KUN,STAT) ! heliocentric correction time
553#endif
554
555
556// RBUF(1) = BERV ! barocentric RV correction
557// RBUF(2) = HERV ! heliocentric RV correction
558// ... (note that EDV is already contained in both BERV and HERV)
559// RBUF(3) = EDV ! diurnal RV correction
560 rbuf[1] = berv;
561 rbuf[2] = herv;
562 rbuf[3] = EDV;
563
564// CALL STKWRR('OUTPUTR',RBUF,1,3,KUN,STAT)
565 outputr[1] = rbuf[1];
566 outputr[2] = rbuf[2];
567 outputr[3] = rbuf[3];
568// GOTO 9000
569 return;
570}
571
572/* @cond Convert FORTRAN indexing -> C indexing */
573#define DCFEL(x,y) dcfel[y][x]
574#define DCFEPS(x,y) dcfeps[y][x]
575#define CCSEL(x,y) ccsel[y][x]
576#define DCARGS(x,y) dcargs[y][x]
577#define CCAMPS(x,y) ccamps[y][x]
578#define CCSEC(x,y) ccsec[y][x]
579#define DCARGM(x,y) dcargm[y][x]
580#define CCAMPM(x,y) ccampm[y][x]
581#define DCEPS(x) dceps[x]
582#define FORBEL(x) forbel[x]
583#define SORBEL(x) sorbel[x]
584#define SN(x) sn[x]
585#define SINLP(x) sinlp[x]
586#define COSLP(x) coslp[x]
587#define CCPAMV(x) ccpamv[x]
588/* @endcond */
589/*----------------------------------------------------------------------------*/
602/*----------------------------------------------------------------------------*/
603
604// SUBROUTINE BARVEL(DJE,DEQ,DVELH,DVELB)
605
606static
607void barvel(double DJE, double DEQ,
608 double DVELH[4], double DVELB[4])
609{
610// DOUBLE PRECISION DJE,DEQ,DVELH(3),DVELB(3),SN(4)
611 double sn[5];
612// DOUBLE PRECISION DT,DTL,DCT0,DCJUL,DTSQ,DLOCAL,DC2PI,CC2PI
613 double DT,DTL,DTSQ,DLOCAL;
614// DOUBLE PRECISION DRD,DRLD,DCSLD,DC1
615 double DRD,DRLD;
616// DOUBLE PRECISION DXBD,DYBD,DZBD,DZHD,DXHD,DYHD
617 double DXBD,DYBD,DZBD,DZHD,DXHD,DYHD;
618// DOUBLE PRECISION DYAHD,DZAHD,DYABD,DZABD
619 double DYAHD,DZAHD,DYABD,DZABD;
620// DOUBLE PRECISION DML,DEPS,PHI,PHID,PSID,DPARAM,PARAM
621 double DML,DEPS,PHI,PHID,PSID,DPARAM,PARAM;
622// DOUBLE PRECISION CCFDI,CCKM,CCMLD,PLON,POMG,PECC
623 double PLON,POMG,PECC;
624// DOUBLE PRECISION PERTL,PERTLD,PERTRD,PERTP,PERTR,PERTPD
625 double PERTL,PERTLD,PERTRD,PERTP,PERTR,PERTPD;
626// DOUBLE PRECISION SINA,CCSGD,DC1MME,TL
627 double SINA,TL;
628// DOUBLE PRECISION CCSEC3,COSA,ESQ
629 double COSA,ESQ;
630// DOUBLE PRECISION DCFEL(3,8),DCEPS(3),CCSEL(3,17),DCARGS(2,15)
631// DOUBLE PRECISION CCAMPS(5,15),CCSEC(3,4),DCARGM(2,3)
632// DOUBLE PRECISION CCAMPM(4,3),CCPAMV(4)
633// DOUBLE PRECISION A,B,E,F,G,SINF,COSF,T,TSQ,TWOE,TWOG
634 double A,B,F,SINF,COSF,T,TSQ,TWOE,TWOG;
635//C
636// DOUBLE PRECISION DPREMA(3,3),DPSI,D1PDRO,DSINLS
637 double DPSI,D1PDRO,DSINLS;
638// DOUBLE PRECISION DCOSLS,DSINEP,DCOSEP
639 double DCOSLS,DSINEP,DCOSEP;
640// DOUBLE PRECISION FORBEL(7),SORBEL(17),SINLP(4),COSLP(4)
641 double forbel[8], sorbel[18], sinlp[5], coslp[5];
642// DOUBLE PRECISION SINLM,COSLM,SIGMA
643 double SINLM,COSLM,SIGMA;
644//C
645// INTEGER IDEQ,K,N
646// int IDEQ;
647 int K,N;
648//C
649// COMMON /BARXYZ/ DPREMA,DPSI,D1PDRO,DSINLS,DCOSLS,
650// + DSINEP,DCOSEP,FORBEL,SORBEL,SINLP,
651// + COSLP,SINLM,COSLM,SIGMA,IDEQ
652
653// EQUIVALENCE (SORBEL(1),E),(FORBEL(1),G)
654 double *E = sorbel + 1 - 1;
655 double *G = forbel + 1 - 1;
656//C
657// DATA DC2PI/6.2831853071796D0/,CC2PI/6.283185/,
658 double DC2PI = 6.2831853071796E0;
659 double CC2PI = 6.283185; /* ??? */
660
661// *DC1/1.0D0/,DCT0/2415020.0D0/,DCJUL/36525.0D0/
662 double DC1 = 1.0;
663 double DCT0 = 2415020.0E0;
664 double DCJUL = 36525.0E0;
665//C
666// DATA DCFEL/ 1.7400353D+00, 6.2833195099091D+02, 5.2796D-06,
667// * 6.2565836D+00, 6.2830194572674D+02,-2.6180D-06,
668// * 4.7199666D+00, 8.3997091449254D+03,-1.9780D-05,
669// * 1.9636505D-01, 8.4334662911720D+03,-5.6044D-05,
670// * 4.1547339D+00, 5.2993466764997D+01, 5.8845D-06,
671// * 4.6524223D+00, 2.1354275911213D+01, 5.6797D-06,
672// * 4.2620486D+00, 7.5025342197656D+00, 5.5317D-06,
673// * 1.4740694D+00, 3.8377331909193D+00, 5.6093D-06/
674
675 double dcfel[][4] = { {0, 0, 0, 0},
676 {0, 1.7400353E+00, 6.2833195099091E+02, 5.2796E-06},
677 {0, 6.2565836E+00, 6.2830194572674E+02,-2.6180E-06},
678 {0, 4.7199666E+00, 8.3997091449254E+03,-1.9780E-05},
679 {0, 1.9636505E-01, 8.4334662911720E+03,-5.6044E-05},
680 {0, 4.1547339E+00, 5.2993466764997E+01, 5.8845E-06},
681 {0, 4.6524223E+00, 2.1354275911213E+01, 5.6797E-06},
682 {0, 4.2620486E+00, 7.5025342197656E+00, 5.5317E-06},
683 {0, 1.4740694E+00, 3.8377331909193E+00, 5.6093E-06} };
684
685//C
686// DATA DCEPS/ 4.093198D-01,-2.271110D-04,-2.860401D-08/
687 double dceps[4] = {0, 4.093198E-01,-2.271110E-04,-2.860401E-08};
688
689//C
690// DATA CCSEL/ 1.675104D-02,-4.179579D-05,-1.260516D-07,
691// * 2.220221D-01, 2.809917D-02, 1.852532D-05,
692// * 1.589963D+00, 3.418075D-02, 1.430200D-05,
693// * 2.994089D+00, 2.590824D-02, 4.155840D-06,
694// * 8.155457D-01, 2.486352D-02, 6.836840D-06,
695// * 1.735614D+00, 1.763719D-02, 6.370440D-06,
696// * 1.968564D+00, 1.524020D-02,-2.517152D-06,
697// * 1.282417D+00, 8.703393D-03, 2.289292D-05,
698// * 2.280820D+00, 1.918010D-02, 4.484520D-06,
699// * 4.833473D-02, 1.641773D-04,-4.654200D-07,
700// * 5.589232D-02,-3.455092D-04,-7.388560D-07,
701// * 4.634443D-02,-2.658234D-05, 7.757000D-08,
702// * 8.997041D-03, 6.329728D-06,-1.939256D-09,
703// * 2.284178D-02,-9.941590D-05, 6.787400D-08,
704// * 4.350267D-02,-6.839749D-05,-2.714956D-07,
705// * 1.348204D-02, 1.091504D-05, 6.903760D-07,
706// * 3.106570D-02,-1.665665D-04,-1.590188D-07/
707
708 double ccsel[][4] = { {0, 0, 0, 0},
709 {0, 1.675104E-02, -4.179579E-05, -1.260516E-07},
710 {0, 2.220221E-01, 2.809917E-02, 1.852532E-05},
711 {0, 1.589963E+00, 3.418075E-02, 1.430200E-05},
712 {0, 2.994089E+00, 2.590824E-02, 4.155840E-06},
713 {0, 8.155457E-01, 2.486352E-02, 6.836840E-06},
714 {0, 1.735614E+00, 1.763719E-02, 6.370440E-06},
715 {0, 1.968564E+00, 1.524020E-02, -2.517152E-06},
716 {0, 1.282417E+00, 8.703393E-03, 2.289292E-05},
717 {0, 2.280820E+00, 1.918010E-02, 4.484520E-06},
718 {0, 4.833473E-02, 1.641773E-04, -4.654200E-07},
719 {0, 5.589232E-02, -3.455092E-04, -7.388560E-07},
720 {0, 4.634443E-02, -2.658234E-05, 7.757000E-08},
721 {0, 8.997041E-03, 6.329728E-06, -1.939256E-09},
722 {0, 2.284178E-02, -9.941590E-05, 6.787400E-08},
723 {0, 4.350267E-02, -6.839749E-05, -2.714956E-07},
724 {0, 1.348204E-02, 1.091504E-05, 6.903760E-07},
725 {0, 3.106570E-02, -1.665665E-04, -1.590188E-07} };
726
727
728
729// DATA DCARGS/ 5.0974222D+00,-7.8604195454652D+02,
730// * 3.9584962D+00,-5.7533848094674D+02,
731// * 1.6338070D+00,-1.1506769618935D+03,
732// * 2.5487111D+00,-3.9302097727326D+02,
733// * 4.9255514D+00,-5.8849265665348D+02,
734// * 1.3363463D+00,-5.5076098609303D+02,
735// * 1.6072053D+00,-5.2237501616674D+02,
736// * 1.3629480D+00,-1.1790629318198D+03,
737// * 5.5657014D+00,-1.0977134971135D+03,
738// * 5.0708205D+00,-1.5774000881978D+02,
739// * 3.9318944D+00, 5.2963464780000D+01,
740// * 4.8989497D+00, 3.9809289073258D+01,
741// * 1.3097446D+00, 7.7540959633708D+01,
742// * 3.5147141D+00, 7.9618578146517D+01,
743// * 3.5413158D+00,-5.4868336758022D+02/
744
745 double dcargs[][3] = { {0, 0, 0},
746 {0, 5.0974222E+00, -7.8604195454652E+02},
747 {0, 3.9584962E+00, -5.7533848094674E+02},
748 {0, 1.6338070E+00, -1.1506769618935E+03},
749 {0, 2.5487111E+00, -3.9302097727326E+02},
750 {0, 4.9255514E+00, -5.8849265665348E+02},
751 {0, 1.3363463E+00, -5.5076098609303E+02},
752 {0, 1.6072053E+00, -5.2237501616674E+02},
753 {0, 1.3629480E+00, -1.1790629318198E+03},
754 {0, 5.5657014E+00, -1.0977134971135E+03},
755 {0, 5.0708205E+00, -1.5774000881978E+02},
756 {0, 3.9318944E+00, 5.2963464780000E+01},
757 {0, 4.8989497E+00, 3.9809289073258E+01},
758 {0, 1.3097446E+00, 7.7540959633708E+01},
759 {0, 3.5147141E+00, 7.9618578146517E+01},
760 {0, 3.5413158E+00, -5.4868336758022E+02} };
761
762// DATA CCAMPS/
763// *-2.279594D-5, 1.407414D-5, 8.273188D-6, 1.340565D-5,-2.490817D-7,
764// *-3.494537D-5, 2.860401D-7, 1.289448D-7, 1.627237D-5,-1.823138D-7,
765// * 6.593466D-7, 1.322572D-5, 9.258695D-6,-4.674248D-7,-3.646275D-7,
766// * 1.140767D-5,-2.049792D-5,-4.747930D-6,-2.638763D-6,-1.245408D-7,
767// * 9.516893D-6,-2.748894D-6,-1.319381D-6,-4.549908D-6,-1.864821D-7,
768// * 7.310990D-6,-1.924710D-6,-8.772849D-7,-3.334143D-6,-1.745256D-7,
769// *-2.603449D-6, 7.359472D-6, 3.168357D-6, 1.119056D-6,-1.655307D-7,
770// *-3.228859D-6, 1.308997D-7, 1.013137D-7, 2.403899D-6,-3.736225D-7,
771// * 3.442177D-7, 2.671323D-6, 1.832858D-6,-2.394688D-7,-3.478444D-7,
772// * 8.702406D-6,-8.421214D-6,-1.372341D-6,-1.455234D-6,-4.998479D-8,
773// *-1.488378D-6,-1.251789D-5, 5.226868D-7,-2.049301D-7, 0.0D0,
774// *-8.043059D-6,-2.991300D-6, 1.473654D-7,-3.154542D-7, 0.0D0,
775// * 3.699128D-6,-3.316126D-6, 2.901257D-7, 3.407826D-7, 0.0D0,
776// * 2.550120D-6,-1.241123D-6, 9.901116D-8, 2.210482D-7, 0.0D0,
777// *-6.351059D-7, 2.341650D-6, 1.061492D-6, 2.878231D-7, 0.0D0/
778
779 double ccamps[][6] =
780 {{0, 0, 0, 0, 0, 0},
781 {0, -2.279594E-5, 1.407414E-5, 8.273188E-6, 1.340565E-5, -2.490817E-7},
782 {0, -3.494537E-5, 2.860401E-7, 1.289448E-7, 1.627237E-5, -1.823138E-7},
783 {0, 6.593466E-7, 1.322572E-5, 9.258695E-6, -4.674248E-7, -3.646275E-7},
784 {0, 1.140767E-5, -2.049792E-5, -4.747930E-6, -2.638763E-6, -1.245408E-7},
785 {0, 9.516893E-6, -2.748894E-6, -1.319381E-6, -4.549908E-6, -1.864821E-7},
786 {0, 7.310990E-6, -1.924710E-6, -8.772849E-7, -3.334143E-6, -1.745256E-7},
787 {0, -2.603449E-6, 7.359472E-6, 3.168357E-6, 1.119056E-6, -1.655307E-7},
788 {0, -3.228859E-6, 1.308997E-7, 1.013137E-7, 2.403899E-6, -3.736225E-7},
789 {0, 3.442177E-7, 2.671323E-6, 1.832858E-6, -2.394688E-7, -3.478444E-7},
790 {0, 8.702406E-6, -8.421214E-6, -1.372341E-6, -1.455234E-6, -4.998479E-8},
791 {0, -1.488378E-6, -1.251789E-5, 5.226868E-7, -2.049301E-7, 0.0E0},
792 {0, -8.043059E-6, -2.991300E-6, 1.473654E-7, -3.154542E-7, 0.0E0},
793 {0, 3.699128E-6, -3.316126E-6, 2.901257E-7, 3.407826E-7, 0.0E0},
794 {0, 2.550120E-6, -1.241123E-6, 9.901116E-8, 2.210482E-7, 0.0E0},
795 {0, -6.351059E-7, 2.341650E-6, 1.061492E-6, 2.878231E-7, 0.0E0}};
796
797
798// DATA CCSEC3/-7.757020D-08/
799 double CCSEC3 = -7.757020E-08;
800//C
801// DATA CCSEC/ 1.289600D-06, 5.550147D-01, 2.076942D+00,
802// * 3.102810D-05, 4.035027D+00, 3.525565D-01,
803// * 9.124190D-06, 9.990265D-01, 2.622706D+00,
804// * 9.793240D-07, 5.508259D+00, 1.559103D+01/
805
806 double ccsec[][4] = { {0, 0, 0, 0},
807 {0, 1.289600E-06, 5.550147E-01, 2.076942E+00},
808 {0, 3.102810E-05, 4.035027E+00, 3.525565E-01},
809 {0, 9.124190E-06, 9.990265E-01, 2.622706E+00},
810 {0, 9.793240E-07, 5.508259E+00, 1.559103E+01}};
811
812//C
813// DATA DCSLD/ 1.990987D-07/, CCSGD/ 1.990969D-07/
814 double DCSLD = 1.990987E-07, CCSGD = 1.990969E-07;
815//C
816// DATA CCKM/3.122140D-05/, CCMLD/2.661699D-06/, CCFDI/2.399485D-07/
817 double CCKM = 3.122140E-05, CCMLD = 2.661699E-06, CCFDI = 2.399485E-07;
818//C
819// DATA DCARGM/ 5.1679830D+00, 8.3286911095275D+03,
820// * 5.4913150D+00,-7.2140632838100D+03,
821// * 5.9598530D+00, 1.5542754389685D+04/
822
823 double dcargm[][3] = {{0, 0, 0},
824 {0, 5.1679830E+00, 8.3286911095275E+03},
825 {0, 5.4913150E+00, -7.2140632838100E+03},
826 {0, 5.9598530E+00, 1.5542754389685E+04}};
827//C
828// DATA CCAMPM/
829// * 1.097594D-01, 2.896773D-07, 5.450474D-02, 1.438491D-07,
830// * -2.223581D-02, 5.083103D-08, 1.002548D-02,-2.291823D-08,
831// * 1.148966D-02, 5.658888D-08, 8.249439D-03, 4.063015D-08/
832
833 double ccampm[][5] = {{0, 0, 0, 0, 0},
834 {0, 1.097594E-01, 2.896773E-07, 5.450474E-02, 1.438491E-07},
835 {0, -2.223581E-02, 5.083103E-08, 1.002548E-02, -2.291823E-08},
836 {0, 1.148966E-02, 5.658888E-08, 8.249439E-03, 4.063015E-08} };
837
838//C
839// DATA CCPAMV/8.326827D-11,1.843484D-11,1.988712D-12,1.881276D-12/,
840 double ccpamv[] = {0, 8.326827E-11, 1.843484E-11, 1.988712E-12, 1.881276E-12};
841// * DC1MME/0.99999696D0/
842 double DC1MME = 0.99999696E0;
843//C
844
845// IDEQ=DEQ
846// IDEQ=DEQ;
847
848// DT=(DJE-DCT0)/DCJUL
849 DT=(DJE-DCT0)/DCJUL;
850
851// T=DT
852 T=DT;
853
854// DTSQ=DT*DT
855 DTSQ=DT*DT;
856
857// TSQ=DTSQ
858 TSQ=DTSQ;
859
860 DML = 0; /* Suppress warning */
861// DO 100, K=1,8
862 for (K = 1; K <= 8; K++) {
863
864// DLOCAL=DMOD(DCFEL(1,K)+DT*DCFEL(2,K)+DTSQ*DCFEL(3,K),DC2PI)
865 DLOCAL=fmod(DCFEL(1,K)+DT*DCFEL(2,K)+DTSQ*DCFEL(3,K),DC2PI);
866
867// IF (K.EQ.1) DML=DLOCAL
868 if (K == 1) DML=DLOCAL;
869
870// IF (K.NE.1) FORBEL(K-1)=DLOCAL
871 if (K != 1) FORBEL(K-1)=DLOCAL;
872// 100 CONTINUE
873 }
874
875// DEPS=DMOD(DCEPS(1)+DT*DCEPS(2)+DTSQ*DCEPS(3), DC2PI)
876 DEPS=fmod(DCEPS(1)+DT*DCEPS(2)+DTSQ*DCEPS(3), DC2PI);
877
878// DO 200, K=1,17
879 for (K = 1; K <= 17; K++) {
880
881// SORBEL(K)=DMOD(CCSEL(1,K)+T*CCSEL(2,K)+TSQ*CCSEL(3,K),CC2PI)
882 SORBEL(K)=fmod(CCSEL(1,K)+T*CCSEL(2,K)+TSQ*CCSEL(3,K),CC2PI);
883
884// 200 CONTINUE
885 }
886
887// DO 300, K=1,4
888 for (K = 1; K <= 4; K++) {
889
890// A=DMOD(CCSEC(2,K)+T*CCSEC(3,K),CC2PI)
891 A=fmod(CCSEC(2,K)+T*CCSEC(3,K),CC2PI);
892
893// SN(K)=DSIN(A)
894 SN(K)=sin(A);
895// 300 CONTINUE
896 }
897
898// PERTL = CCSEC(1,1) *SN(1) +CCSEC(1,2)*SN(2)
899// * +(CCSEC(1,3)+T*CCSEC3)*SN(3) +CCSEC(1,4)*SN(4)
900
901 PERTL = CCSEC(1,1) *SN(1) +CCSEC(1,2)*SN(2)
902 +(CCSEC(1,3)+T*CCSEC3)*SN(3) +CCSEC(1,4)*SN(4);
903
904// PERTLD=0.0
905// PERTR =0.0
906// PERTRD=0.0
907 PERTLD=0.0;
908 PERTR =0.0;
909 PERTRD=0.0;
910
911// DO 400, K=1,15
912 for (K = 1; K <= 15; K++) {
913// A=DMOD(DCARGS(1,K)+DT*DCARGS(2,K), DC2PI)
914 A=fmod(DCARGS(1,K)+DT*DCARGS(2,K), DC2PI);
915
916// COSA=DCOS(A)
917 COSA=cos(A);
918
919// SINA=DSIN(A)
920 SINA=sin(A);
921
922// PERTL =PERTL+CCAMPS(1,K)*COSA+CCAMPS(2,K)*SINA
923 PERTL =PERTL+CCAMPS(1,K)*COSA+CCAMPS(2,K)*SINA;
924
925// PERTR =PERTR+CCAMPS(3,K)*COSA+CCAMPS(4,K)*SINA;
926 PERTR =PERTR+CCAMPS(3,K)*COSA+CCAMPS(4,K)*SINA;
927
928// IF (K.GE.11) GO TO 400
929 if (K >= 11) break;
930
931// PERTLD=PERTLD+(CCAMPS(2,K)*COSA-CCAMPS(1,K)*SINA)*CCAMPS(5,K)
932 PERTLD=PERTLD+(CCAMPS(2,K)*COSA-CCAMPS(1,K)*SINA)*CCAMPS(5,K);
933
934// PERTRD=PERTRD+(CCAMPS(4,K)*COSA-CCAMPS(3,K)*SINA)*CCAMPS(5,K)
935 PERTRD=PERTRD+(CCAMPS(4,K)*COSA-CCAMPS(3,K)*SINA)*CCAMPS(5,K);
936
937// 400 CONTINUE
938 }
939
940// ESQ=E*E
941 ESQ=E[1]*E[1];
942
943// DPARAM=DC1-ESQ
944 DPARAM=DC1-ESQ;
945
946// PARAM=DPARAM
947 PARAM=DPARAM;
948
949// TWOE=E+E
950 TWOE=E[1]+E[1];
951
952// TWOG=G+G
953 TWOG=G[1]+G[1];
954
955// PHI=TWOE*((1.0-ESQ*0.125D0)*DSIN(G)+E*0.625D0*DSIN(TWOG)
956// * +ESQ*0.5416667D0*DSIN(G+TWOG) )
957
958 PHI=TWOE*((1.0-ESQ*0.125 )*sin(G[1])+E[1]*0.625 *sin(TWOG)
959 +ESQ*0.5416667 *sin(G[1]+TWOG) ) ;
960
961 //F=G+PHI
962 F=G[1]+PHI;
963
964 //SINF=DSIN(F)
965 SINF=sin(F);
966
967 //COSF=DCOS(F)
968 COSF=cos(F);
969
970 //DPSI=DPARAM/(DC1+E*COSF)
971 DPSI=DPARAM/(DC1+E[1]*COSF);
972
973// PHID=TWOE*CCSGD*((1.0+ESQ*1.5D0)*COSF+E[1]*(1.25D0-SINF*SINF*0.5D0))
974 PHID=TWOE*CCSGD*((1.0+ESQ*1.5 )*COSF+E[1]*(1.25 -SINF*SINF*0.5 ));
975
976// PSID=CCSGD*E*SINF/SQRT(PARAM)
977 PSID=CCSGD*E[1]*SINF/sqrt(PARAM);
978
979// D1PDRO=(DC1+PERTR)
980 D1PDRO=(DC1+PERTR);
981
982// DRD=D1PDRO*(PSID+DPSI*PERTRD)
983 DRD=D1PDRO*(PSID+DPSI*PERTRD);
984
985// DRLD=D1PDRO*DPSI*(DCSLD+PHID+PERTLD)
986 DRLD=D1PDRO*DPSI*(DCSLD+PHID+PERTLD);
987
988// DTL=DMOD(DML+PHI+PERTL, DC2PI)
989 DTL=fmod(DML+PHI+PERTL, DC2PI);
990
991// DSINLS=DSIN(DTL)
992 DSINLS=sin(DTL);
993
994// DCOSLS=DCOS(DTL)
995 DCOSLS=cos(DTL);
996
997// DXHD = DRD*DCOSLS-DRLD*DSINLS
998 DXHD = DRD*DCOSLS-DRLD*DSINLS;
999
1000// DYHD = DRD*DSINLS+DRLD*DCOSLS
1001 DYHD = DRD*DSINLS+DRLD*DCOSLS;
1002
1003// PERTL =0.0
1004 PERTL =0.0;
1005// PERTLD=0.0
1006 PERTLD=0.0;
1007// PERTP =0.0
1008 PERTP =0.0;
1009// PERTPD=0.0
1010 PERTPD=0.0;
1011
1012 //DO 500 K=1,3
1013 for (K = 1; K <= 3; K++) {
1014 //A=DMOD(DCARGM(1,K)+DT*DCARGM(2,K), DC2PI)
1015 A=fmod(DCARGM(1,K)+DT*DCARGM(2,K), DC2PI);
1016
1017 //SINA =DSIN(A)
1018 SINA =sin(A);
1019
1020 //COSA =DCOS(A)
1021 COSA =cos(A);
1022
1023 //PERTL =PERTL +CCAMPM(1,K)*SINA
1024 PERTL =PERTL +CCAMPM(1,K)*SINA;
1025
1026 //PERTLD=PERTLD+CCAMPM(2,K)*COSA
1027 PERTLD=PERTLD+CCAMPM(2,K)*COSA;
1028
1029 //PERTP =PERTP +CCAMPM(3,K)*COSA
1030 PERTP =PERTP +CCAMPM(3,K)*COSA;
1031
1032 //PERTPD=PERTPD-CCAMPM(4,K)*SINA
1033 PERTPD=PERTPD-CCAMPM(4,K)*SINA;
1034
1035// 500 CONTINUE
1036 }
1037
1038 //TL=FORBEL(2)+PERTL
1039 TL=FORBEL(2)+PERTL;
1040
1041// SINLM=DSIN(TL)
1042 SINLM=sin(TL);
1043
1044// COSLM=DCOS(TL)
1045 COSLM=cos(TL);
1046
1047// SIGMA=CCKM/(1.0+PERTP)
1048 SIGMA=CCKM/(1.0+PERTP);
1049
1050// A=SIGMA*(CCMLD+PERTLD)
1051 A=SIGMA*(CCMLD+PERTLD);
1052
1053// B=SIGMA*PERTPD
1054 B=SIGMA*PERTPD;
1055
1056// DXHD=DXHD+A*SINLM+B*COSLM
1057 DXHD=DXHD+A*SINLM+B*COSLM;
1058
1059// DYHD=DYHD-A*COSLM+B*SINLM
1060 DYHD=DYHD-A*COSLM+B*SINLM;
1061
1062// DZHD= -SIGMA*CCFDI*DCOS(FORBEL(3))
1063 DZHD= -SIGMA*CCFDI* cos(FORBEL(3));
1064
1065// DXBD=DXHD*DC1MME
1066 DXBD=DXHD*DC1MME;
1067
1068// DYBD=DYHD*DC1MME
1069 DYBD=DYHD*DC1MME;
1070// DZBD=DZHD*DC1MME
1071 DZBD=DZHD*DC1MME;
1072
1073// DO 600 K=1,4
1074 for (K = 1; K <= 4; K++) {
1075
1076 //PLON=FORBEL(K+3)
1077 PLON=FORBEL(K+3);
1078
1079 //POMG=SORBEL(K+1)
1080 POMG=SORBEL(K+1);
1081
1082 //PECC=SORBEL(K+9)
1083 PECC=SORBEL(K+9);
1084
1085 //TL=DMOD(PLON+2.0*PECC*DSIN(PLON-POMG), CC2PI)
1086 TL=fmod(PLON+2.0*PECC* sin(PLON-POMG), CC2PI);
1087
1088 //SINLP(K)=DSIN(TL)
1089 SINLP(K)= sin(TL);
1090
1091 //COSLP(K)=DCOS(TL)
1092 COSLP(K)= cos(TL);
1093
1094 //DXBD=DXBD+CCPAMV(K)*(SINLP(K)+PECC*DSIN(POMG))
1095 DXBD=DXBD+CCPAMV(K)*(SINLP(K)+PECC*sin(POMG));
1096
1097 //DYBD=DYBD-CCPAMV(K)*(COSLP(K)+PECC*DCOS(POMG))
1098 DYBD=DYBD-CCPAMV(K)*(COSLP(K)+PECC*cos(POMG));
1099
1100 //DZBD=DZBD-CCPAMV(K)*SORBEL(K+13)*DCOS(PLON-SORBEL(K+5))
1101 DZBD=DZBD-CCPAMV(K)*SORBEL(K+13)*cos(PLON-SORBEL(K+5));
1102
1103// 600 CONTINUE
1104 }
1105
1106 //DCOSEP=DCOS(DEPS)
1107 DCOSEP=cos(DEPS);
1108 //DSINEP=DSIN(DEPS)
1109 DSINEP=sin(DEPS);
1110 //DYAHD=DCOSEP*DYHD-DSINEP*DZHD
1111 DYAHD=DCOSEP*DYHD-DSINEP*DZHD;
1112 //DZAHD=DSINEP*DYHD+DCOSEP*DZHD
1113 DZAHD=DSINEP*DYHD+DCOSEP*DZHD;
1114 //DYABD=DCOSEP*DYBD-DSINEP*DZBD
1115 DYABD=DCOSEP*DYBD-DSINEP*DZBD;
1116 //DZABD=DSINEP*DYBD+DCOSEP*DZBD
1117 DZABD=DSINEP*DYBD+DCOSEP*DZBD;
1118
1119 //DVELH(1)=DXHD
1120 DVELH[1]=DXHD;
1121 //DVELH(2)=DYAHD
1122 DVELH[2]=DYAHD;
1123 //DVELH(3)=DZAHD
1124 DVELH[3]=DZAHD;
1125
1126 //DVELB(1)=DXBD
1127 DVELB[1]=DXBD;
1128 //DVELB(2)=DYABD
1129 DVELB[2]=DYABD;
1130 //DVELB(3)=DZABD
1131 DVELB[3]=DZABD;
1132 //DO 800 N=1,3
1133 for (N = 1; N <= 3; N++) {
1134 //DVELH(N)=DVELH(N)*1.4959787D8
1135 DVELH[N]=DVELH[N]*1.4959787E8;
1136 //DVELB(N)=DVELB(N)*1.4959787D8
1137 DVELB[N]=DVELB[N]*1.4959787E8;
1138// 800 CONTINUE
1139 }
1140// RETURN
1141 return;
1142}
1143
1144/*----------------------------------------------------------------------------*/
1152/*----------------------------------------------------------------------------*/
1153
1154static void
1155deg2dms(double in_val,
1156 double *degs,
1157 double *minutes,
1158 double *seconds)
1159{
1160 deg2hms(in_val*15, degs, minutes, seconds);
1161}
1162
1166#define MIDAS_BUG 0
1167/*----------------------------------------------------------------------------*/
1175/*----------------------------------------------------------------------------*/
1176
1177static void
1178deg2hms(double in_val,
1179 double *hours,
1180 double *minutes,
1181 double *seconds)
1182{
1183// define/parameter p1 ? num "Enter value in deg units"
1184// define/local in_val/d/1/1 {p1}
1185//define/local out_val/c/1/80 " " all
1186//define/local hours/i/1/1 0
1187//define/local minutes/i/1/1 0
1188//define/local seconds/d/1/1 0
1189
1190//define/local tmp/d/1/1 0
1191 double tmp;
1192//define/local hold/c/1/80 " " all
1193//define/local sign/c/1/1 " "
1194
1195 char sign;
1196
1197//hold = "{in_val}"
1198//if m$index(hold,"-") .gt. 0 then
1199// in_val = m$abs(in_val)
1200// sign = "-"
1201//else
1202// sign = "+"
1203//endif
1204 if (in_val < 0) {
1205 in_val = fabs(in_val);
1206 sign = '-';
1207 }
1208 else {
1209 sign = '+';
1210 }
1211
1212//set/format i1
1214// tmp = in_val / 15
1215 tmp = in_val / 15;
1216
1217// hours = tmp !takes the integer part = hours
1218#if MIDAS_BUG
1219 *hours= xsh_round_double(tmp);
1220#else
1221 *hours= (int) tmp;
1222#endif
1223
1224// tmp = tmp - hours !takes the mantissa
1225 tmp = tmp - *hours;
1226// tmp = tmp * 60 !converts the mantissa in minutes
1227 tmp = tmp * 60;
1228
1229// minutes = tmp !takes the integer part = minutes
1230#if MIDAS_BUG
1231 *minutes= xsh_round_double(tmp);
1232#else
1233 *minutes= (int) tmp;
1234#endif
1235
1236// tmp = tmp - minutes !takes the mantissa
1237 tmp = tmp - *minutes;
1238
1239// seconds = tmp * 60 !converts the mantissa in seconds = seconds (with decimal)
1240 *seconds= tmp * 60;
1241
1242//out_val = "{sign}{hours},{minutes},{seconds}"
1243
1244 /* Rather than returning it explicitly, just attach sign to hours */
1245 if (sign == '-') *hours = -(*hours);
1246
1247 return;
1248}
1249
1250double xsh_baryvel_correct_value(double value, double correction) {
1251 double return_value = 1 + (correction / (CPL_PHYS_C * 1e-3)); // km/s
1252 return_value *= value ;
1253 return return_value ;
1254}
1255
1256/* PIPE-10061: Helper function to extract correction factor */
1258 cpl_propertylist* keys, xsh_bary_corr_param* bary_corr_param) {
1259
1260 double baryvalue = 0.0 ;
1261
1262 if (bary_corr_param->corr_method == BARY_CORR_BARY &&
1263 cpl_propertylist_has(keys, XSH_QC_VRAD_BARYCOR)) {
1264 baryvalue = cpl_propertylist_get_double(keys, XSH_QC_VRAD_BARYCOR);
1265 // TEST - make the correction huge
1266 // baryvalue = -100000. ;
1267 } else if (bary_corr_param->corr_method == BARY_CORR_HELIO &&
1268 cpl_propertylist_has(keys, XSH_QC_VRAD_HELICOR)) {
1269 baryvalue = cpl_propertylist_get_double(keys, XSH_QC_VRAD_HELICOR);
1270 // TEST - make the correction huge
1271 // baryvalue = 100000. ;
1272 }
1273
1274 return baryvalue ;
1275}
1276
1277/* PIPE-10061: This is a shameless copy of esotk_barycorr_adjust_header */
1278cpl_error_code
1279xsh_baryvel_adjust_header(cpl_propertylist * header,
1280 const char * keyword, double barycorr){
1281 /* Adjust IDP header keywords */
1282 if (cpl_propertylist_has(header, keyword)) {
1283 double value = cpl_propertylist_get_double(header, keyword);
1284 double value_c = xsh_baryvel_correct_value(value, barycorr);
1285 cpl_msg_info(__func__, "Updating %s: %f --> %f",
1286 keyword, value, value_c);
1287// getchar();
1288 cpl_propertylist_update_double(header, keyword, value_c);
1289 }
1290 return cpl_error_get_code();
1291}
1292
1293// FIXME use variables here once found to be working
1294cpl_error_code xsh_baryvel_set_specsys(cpl_propertylist * header,
1295 xsh_bary_corr_param* bary_param) {
1296 if (bary_param->corr_method == BARY_CORR_BARY) {
1297 cpl_propertylist_update_string(header, "SPECSYS",
1298 "BARYCORR");
1299 }
1300 else if (bary_param->corr_method == BARY_CORR_HELIO) {
1301 cpl_propertylist_update_string(header, "SPECSYS",
1302 "HELICORR");
1303 }
1304 else if (bary_param->corr_method == BARY_CORR_NONE) {
1305 cpl_propertylist_update_string(header, "SPECSYS",
1306 "TOPOCENT");
1307 }
1308}
1309
1310cpl_error_code xsh_baryvel_correct_header(cpl_propertylist * header,
1311 cpl_propertylist * bary_header,
1312 xsh_bary_corr_param * bary_param) {
1313 if (bary_param != NULL) {
1314 double baryvalue = 0.0 ;
1315 check(baryvalue = xsh_baryvel_get_correction_factor(bary_header,
1316 bary_param));
1317 if (fabs(baryvalue) > 1e-10) {
1318 check(xsh_baryvel_adjust_header(header, "CRVAL1",
1319 baryvalue));
1320 check(xsh_baryvel_adjust_header(header, "CDELT1",
1321 baryvalue));
1322 check(xsh_baryvel_set_specsys(header, bary_param));
1323 }
1324 }
1325
1326 cleanup:
1327
1328 return cpl_error_get_code();
1329}
1330#ifdef XSH_UNUSED
1331/* PIPE-10061: Implement a correction for irplib_sdp_spectrum objects */
1332cpl_error_code xsh_baryvel_correct_irplib_spectrum(irplib_sdp_spectrum* self,
1333 xsh_bary_corr_param* bary_corr_param,
1334 cpl_parameterlist* keys) {
1335 // Grab the data to alter
1336 // Note we don't need to explicitly write this data back, we're
1337 // modifying in place due to the pointer
1338 cpl_array* wave_data = irplib_sdp_spectrum_get_column_data(
1339 self, XSH_SDP_COLUMN_WAVE);
1340
1341 double baryvalue = 0.0;
1342 check(baryvalue = xsh_baryvel_get_correction_factor(keys, bary_corr_param));
1343
1344 if (fabs(baryvalue) > 1e-10) {
1345 check(cpl_array_multiply_scalar(wave_data,
1346 (1.0 + (baryvalue / CPL_PHYS_C * 1e-3))));
1347
1348
1349 // Alter other necessary header keywords
1350 // This replicates commands in xsh_sdp_spectrum_create
1351 cpl_msg_info(__func__, "Going to set WAVELMIN %f --> %f",
1352 irplib_sdp_spectrum_get_wavelmin(self),
1354 irplib_sdp_spectrum_get_wavelmin(self), baryvalue
1355 ));
1356// getchar();
1357 check(irplib_sdp_spectrum_set_wavelmin(
1358 self,
1360 irplib_sdp_spectrum_get_wavelmin(self), baryvalue
1361 )));
1362 check(irplib_sdp_spectrum_set_wavelmax(
1363 self,
1365 irplib_sdp_spectrum_get_wavelmax(self), baryvalue
1366 )));
1367 check(irplib_sdp_spectrum_set_specbin(
1368 self,
1370 irplib_sdp_spectrum_get_specbin(self), baryvalue
1371 )));
1372 check(irplib_sdp_spectrum_set_tdmin(
1373 self, cpl_array_get_min(wave_data)));
1374 check(irplib_sdp_spectrum_set_tdmax(
1375 self, cpl_array_get_max(wave_data)));
1376 check(irplib_sdp_spectrum_set_specval(
1377 self,
1378 0.5 * (cpl_array_get_min(wave_data) +
1379 cpl_array_get_max(wave_data))));
1380 check(irplib_sdp_spectrum_set_specbw(
1381 self,
1382 cpl_array_get_max(wave_data) - cpl_array_get_min(wave_data)));
1383
1384 // Adjust the SPECSYS header keyword
1385 if (bary_corr_param->corr_method == BARY_CORR_BARY) {
1386 irplib_sdp_spectrum_set_specsys(self, "BARYCORR");
1387 } else if (bary_corr_param->corr_method == BARY_CORR_HELIO) {
1388 irplib_sdp_spectrum_set_specsys(self, "HELICORR");
1389 }
1390 } else {
1391 cpl_msg_warning(__func__, "Unable to do barycentric correction on a "
1392 "file");
1393 }
1394
1395 cleanup:
1396 return cpl_error_get_code() ;
1397}
1398
1400#endif
1401#if 0 /* Not used / needed.
1402 We simply get the julian date from the input FITS header */
1403
1404// SUBROUTINE JULDAT(INDATE,UTR,JD)
1405//C++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1406//C
1407//C.IDENTIFICATION
1408//C FORTRAN subroutine JULDAT version 1.0 870102
1409//C original coding: D. Gillet ESO - Garching
1410//C variables renamed and restructured: D. Baade ST-ECF, Garching
1411//C
1412//C.KEYWORDS
1413//C geocentric Julian date
1414//C
1415//C.PURPOSE
1416//C calculate geocentric Julian date for any civil date (time in UT)
1417//C
1418//C.ALGORITHM
1419//C adapted from MEEUS J.,1980, ASTRONOMICAL FORMULAE FOR CALCULATORS
1420//C
1421//C.INPUT/OUTPUT
1422//C the following are passed from and to the calling program:
1423//C INDATE(3) : civil date as year,month,day OR year.fraction
1424//C UT : universal time expressed in real hours
1425//C JD : real geocentric Julian date
1426//C
1427//C.REVISIONS
1428//C made to accept also REAL dates D. Baade 910408
1429//C
1430//C-------------------------------------------------------------------------------
1431//C
1432
1433static void
1434juldat(double *INDATE,
1435 double UTR,
1436 double *JD)
1437{
1438// DOUBLE PRECISION YP,P,C,A,UT
1439 double UT;
1440
1441// DOUBLE PRECISION UTR,JD
1442//C
1443// INTEGER STAT,IA,IB,IC,ND,DATE(3)
1444 int DATE[4];
1445//C
1446// REAL INDATE(3),FRAC
1447//C
1448
1449// UT=UTR / 24.0D0
1450 UT=UTR / 24.0;
1451
1452// CHECK FORMAT OF DATE: may be either year,month,date OR year.fraction,0,0
1453// (Note that the fraction of the year must NOT include fractions of a day.)
1454// For all other formats exit and terminate also calling command sequence.
1455//
1456// IF ((INDATE(1)-INT(INDATE(1))).GT.1.0E-6) THEN
1457// IF ((INDATE(2).GT.1.0E-6).OR.(INDATE(3).GT.1.0E-6))
1458// + CALL STETER(1,'Error: Date was entered in wrong format.')
1459
1460// copy date input buffer copy to other buffer so that calling program
1461// does not notice any changes
1462
1463// FIRST CASE: format was year.fraction
1464
1465// DATE(1)=INT(INDATE(1))
1466// FRAC=INDATE(1)-DATE(1)
1467// DATE(2)=1
1468// DATE(3)=1
1469// ELSE
1470//
1471// SECOND CASE: format was year,month,day
1472//
1473
1474// DATE(1)=NINT(INDATE(1))
1475 DATE[1]=xsh_round_double(INDATE[1]);
1476
1477 //FRAC=0
1478 FRAC = 0;
1479
1480 //DATE(2)=NINT(INDATE(2))
1481 DATE[2]=xsh_round_double(INDATE[2]);
1482 //DATE(3)=NINT(INDATE(3))
1483 DATE[3]=xsh_round_double(INDATE[3]);
1484
1485 //IF ((DATE(2).EQ.0).AND.(DATE(3).EQ.0)) THEN
1486 if ((DATE[2] == 0) && (DATE[3] == 0)) {
1487 //DATE(2)=1
1488 DATE[2]=1;
1489 //DATE(3)=1
1490 DATE[3]=1;
1491// ENDIF
1492 }
1493
1494// IF ((DATE(2).LT.1).OR.(DATE(2).GT.12))
1495// + CALL STETER(1,'Error: such a month does not exist')
1496// IF ((DATE(3).LT.1).OR.(DATE(3).GT.31))
1497// + CALL STETER(1,'Error: such a day does not exist')
1498// ENDIF
1499
1500// from here on, the normal procedure applies which is based on the
1501// format year,month,day:
1502
1503 //IF (DATE(2) .GT. 2) THEN
1504 if (DATE[2] > 2) {
1505 //YP=DATE(1)
1506 YP=DATE[1];
1507 //P=DATE[2]
1508 P=DATE[2];
1509// ELSE
1510 } else {
1511 //YP=DATE(1)-1
1512 YP=DATE[1]-1;
1513 //P=DATE(2)+12.0
1514 P=DATE(2)+12.0;
1515// ENDIF
1516 }
1517
1518// C = DATE(1) + DATE(2)*1.D-2 + DATE(3)*1.D-4 + UT*1.D-6
1519 C = DATE[1] + DATE[2]*1.E-2 + DATE[3]*1.E-4 + UT*1.E-6;
1520
1521// IF (C .GE. 1582.1015D0) THEN
1522 if (C > 1582.1015E0) {
1523 //IA=IDINT(YP/100.D0)
1524 IA=(int) (YP/100.D0);
1525 //A=DBLE(IA)
1526 A=IA;
1527 //IB=2-IA+IDINT(A/4.D0)
1528 IB=2-IA+((int)(A/4.D0));
1529 //ELSE
1530 } else {
1531 //IB=0
1532 IB=0;
1533 //ENDIF
1534 }
1535
1536// JD = DINT(365.25D0*YP) + DINT(30.6001D0*(P+1.D0)) + DATE(3) + UT
1537// * + DBLE(IB) + 1720994.5D0
1538 *JD = ((int) (365.25E0*YP)) + ((int)(30.6001D0*(P+1.D0))) + DATE[3] + UT
1539 + IB + 1720994.5E0;
1540
1541// finally, take into account fraction of year (if any), respect leap
1542// year conventions
1543//
1544// IF (FRAC.GT.1.0E-6) THEN
1545 if (FRAC > 1.0E-6) {
1546 //ND=365
1547 ND=365;
1548
1549 //IF (C.GE.1582.1015D0) THEN
1550 IF (C >= 1582.1015E0) {
1551 //IC = MOD(DATE(1),4)
1552 IC = DATE[1] % 4;
1553 //IF (IC.EQ.0) THEN
1554 if (IC == 0) {
1555 //ND=366
1556 ND=366;
1557 //IC = MOD(DATE(1),100)
1558 IC = DATE[1] % 100;
1559 //IF (IC.EQ.0) THEN
1560 if (IC == 0) {
1561 //IC = MOD(DATE(1),400)
1562 IC = DATE[1] % 400;
1563 //IF (IC.NE.0) ND=365
1564 if (IC != 0) ND=365;
1565 //ENDIF
1566 }
1567 //ENDIF
1568 }
1569 //ENDIF
1570 }
1571
1572 //IF ( ABS(FRAC*ND-NINT(FRAC*ND)).GT.0.3) THEN
1573 if (fabs(FRAC*ND-xsh_round_double(FRAC*ND)) > 0.3) {
1574// CALL STTPUT
1575// + ('Warning: Fraction of year MAY not correspond to ',STAT)
1576// CALL STTPUT(' integer number of days.',STAT)
1577 xsh_msg_warning("Fraction of year MAY not correspond to "
1578 "integer number of days");
1579// ENDIF
1580 }
1581
1582// JD = JD+NINT(FRAC*ND)
1583 *JD = *JD+xsh_round_double(FRAC*ND);
1584// ENDIF
1585 }
1586
1587// RETURN
1588 return;
1589}
1590#endif
#define N
static void compxy(double inputr[19], char inputc[4], double outputr[4], double utr, double mod_juldat)
Compute velocity correction.
Definition: xsh_baryvel.c:273
double xsh_baryvel_get_correction_factor(cpl_propertylist *keys, xsh_bary_corr_param *bary_corr_param)
Definition: xsh_baryvel.c:1257
static void barvel(double DJE, double DEQ, double DVELH[4], double DVELB[4])
compute rectangular heliocentric and barycentric components of the earth's orbital velocity
Definition: xsh_baryvel.c:607
cpl_error_code xsh_baryvel_adjust_header(cpl_propertylist *header, const char *keyword, double barycorr)
Definition: xsh_baryvel.c:1279
cpl_error_code xsh_baryvel_correct_header(cpl_propertylist *header, cpl_propertylist *bary_header, xsh_bary_corr_param *bary_param)
Definition: xsh_baryvel.c:1310
static void deg2dms(double in_val, double *degs, double *minutes, double *seconds)
convert degrees -> degrees, minutes, seconds
Definition: xsh_baryvel.c:1155
static void deg2hms(double in_val, double *hour, double *min, double *sec)
convert hours -> degrees, minutes, seconds
Definition: xsh_baryvel.c:1178
double xsh_baryvel_correct_value(double value, double correction)
Definition: xsh_baryvel.c:1250
void xsh_baryvel(const cpl_propertylist *raw_header, double *bary_corr, double *helio_corr)
Compute velocity correction.
Definition: xsh_baryvel.c:95
cpl_error_code xsh_baryvel_set_specsys(cpl_propertylist *header, xsh_bary_corr_param *bary_param)
Definition: xsh_baryvel.c:1294
#define check(COMMAND)
Definition: xsh_error.h:71
#define check_msg(COMMAND,...)
Definition: xsh_error.h:62
#define xsh_msg_warning(...)
Print an warning message.
Definition: xsh_msg.h:88
#define xsh_msg_debug(...)
Print a debug message.
Definition: xsh_msg.h:99
double xsh_pfits_get_geolon(const cpl_propertylist *plist)
Find out the telescope longitude.
Definition: xsh_pfits.c:128
double xsh_pfits_get_utc(const cpl_propertylist *plist)
Find out the observation time.
Definition: xsh_pfits.c:142
double xsh_pfits_get_dec(const cpl_propertylist *plist)
Get the Right Ascension.
Definition: xsh_pfits.c:3500
double xsh_pfits_get_geolat(const cpl_propertylist *plist)
Find out the telescope latitude.
Definition: xsh_pfits.c:114
double xsh_pfits_get_ra(const cpl_propertylist *plist)
Get the Right Ascension.
Definition: xsh_pfits.c:3479
double xsh_pfits_get_mjdobs(const cpl_propertylist *plist)
Find out the modified julian observation date.
Definition: xsh_pfits.c:186
long xsh_round_double(double x)
Computes round(x)
Definition: xsh_utils.c:4383
enum bary_corr_method corr_method
#define C(i)
@ BARY_CORR_NONE
@ BARY_CORR_HELIO
@ BARY_CORR_BARY
#define XSH_SDP_COLUMN_WAVE
Definition: xsh_pfits.h:337
#define XSH_QC_VRAD_HELICOR
Definition: xsh_pfits_qc.h:46
#define XSH_QC_VRAD_BARYCOR
Definition: xsh_pfits_qc.h:44
#define M_PI
Definition: xsh_utils.h:43