MUSE Pipeline Reference Manual  2.1.1
muse_astro.c
1 /* -*- Mode: C; tab-width: 2; indent-tabs-mode: nil; c-basic-offset: 2 -*- */
2 /* vim:set sw=2 sts=2 et cin: */
3 /*
4  * This file is part of the MUSE Instrument Pipeline
5  * Copyright (C) 2005-2014 European Southern Observatory
6  * (C) 2000-2002 European Southern Observatory
7  * (C) 2002-2006 European Southern Observatory
8  *
9  * This program is free software; you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation; either version 2 of the License, or
12  * (at your option) any later version.
13  *
14  * This program is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the GNU General Public License
20  * along with this program; if not, write to the Free Software
21  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
22  */
23 
24 #ifdef HAVE_CONFIG_H
25 #include <config.h>
26 #endif
27 
28 /*----------------------------------------------------------------------------*
29  * Includes *
30  *----------------------------------------------------------------------------*/
31 #include <cpl.h>
32 #include <math.h>
33 #include <string.h>
34 #include <ctype.h>
35 
36 #include "muse_astro.h"
37 
38 #include "muse_pfits.h"
39 
40 /*----------------------------------------------------------------------------*
41  * Defines *
42  *----------------------------------------------------------------------------*/
43 
44 #define MUSE_ASTRO_AIRMASS_APPROX 2 /* airmass approximation to use: *
45  * 0 Young & Irvine (1967) *
46  * 1 Young (1994) *
47  * 2 Hardie (1962) */
48 
49 /* Constants used by the radial-velocity correction (for precession etc.) */
50 static const double kT0 = 2415020.; /* t0: Julian date of Dec 31st, 1899, 12 UT */
51 static const double kEpoch1900 = 1900.; /* epoch B1900 */
52 static const double kJulianCentury = 36525.; /* Julian century [days] */
53 static const double kTropYear = 365.242198781; /* tropical year of B1900 [days] */
54 static const double kBesselOffset = 0.31352; /* offset for Besselian epoch [days] */
55 static const double kAU = 149597870.7; /* 1 au = 149 597 870.7 km (IAU 2012) */
56 
57 /*----------------------------------------------------------------------------*/
66 /*----------------------------------------------------------------------------*/
67 
70 /*----------------------------------------------------------------------------*/
86 /*----------------------------------------------------------------------------*/
87 static double
88 muse_astro_get_zenith_distance(double aHourAngle, double aDelta,
89  double aLatitude)
90 {
91  double p0 = sin(aLatitude) * sin(aDelta),
92  p1 = cos(aLatitude) * cos(aDelta),
93  z = p0 + cos(aHourAngle) * p1;
94  return fabs(z) < FLT_EPSILON ? 0. : z;
95 }
96 
97 #if MUSE_ASTRO_AIRMASS_APPROX == 0
98 /*----------------------------------------------------------------------------*/
110 /*----------------------------------------------------------------------------*/
111 static double
112 muse_astro_get_airmass_youngirvine(double aSecZ)
113 {
114  return aSecZ * (1. - 0.0012 * (pow(aSecZ, 2) - 1.));
115 }
116 #endif /* MUSE_ASTRO_AIRMASS_APPROX == 0 */
117 
118 #if MUSE_ASTRO_AIRMASS_APPROX == 1
119 /*----------------------------------------------------------------------------*/
133 /*----------------------------------------------------------------------------*/
134 static double
135 muse_astro_get_airmass_young(double aCosZt)
136 {
137  return (1.002432 * aCosZt*aCosZt + 0.148386 * aCosZt + 0.0096467)
138  / (aCosZt*aCosZt*aCosZt + 0.149864 * aCosZt*aCosZt + 0.0102963 * aCosZt
139  + 0.000303978);
140 }
141 #endif /* MUSE_ASTRO_AIRMASS_APPROX == 1 */
142 
143 #if MUSE_ASTRO_AIRMASS_APPROX == 2
144 /*----------------------------------------------------------------------------*/
157 /*----------------------------------------------------------------------------*/
158 static double
159 muse_astro_get_airmass_hardie(double aSecZ)
160 {
161  double secm1 = aSecZ - 1;
162  return aSecZ - 0.0018167 * secm1 - 0.002875 * secm1*secm1
163  - 0.0008083 * secm1*secm1*secm1;
164 }
165 #endif /* MUSE_ASTRO_AIRMASS_APPROX == 2 */
166 
167 /*----------------------------------------------------------------------------*/
205 /*----------------------------------------------------------------------------*/
206 double
207 muse_astro_compute_airmass(double aRA, double aDEC, double aLST,
208  double aExptime, double aLatitude)
209 {
210  cpl_ensure(aRA >= 0. && aRA < 360. && aDEC >= -90. && aDEC <= 90. &&
211  aLST >= 0. && aLST < 86400. && aLatitude >= -90. && aLatitude <= 90.,
212  CPL_ERROR_ILLEGAL_INPUT, -1.);
213  cpl_ensure(aExptime >= 0., CPL_ERROR_ILLEGAL_INPUT, -1.);
214 
215  /* Compute hour angle of the observation in degrees. */
216  double HA = aLST * 15./3600. - aRA;
217  /* Range adjustments. Angle between line of sight and the meridian is needed. */
218  if (HA < -180.) {
219  HA += 360.;
220  }
221  if (HA > 180.) {
222  HA -= 360.;
223  }
224 
225  /* Convert angles from degrees to radians. */
226  double delta = aDEC * CPL_MATH_RAD_DEG,
227  latitude = aLatitude * CPL_MATH_RAD_DEG,
228  hourangle = HA * CPL_MATH_RAD_DEG;
229 
230  /* Calculate airmass of the observation using the approximation given *
231  * by Young (1994). For non-zero exposure times these airmass values are *
232  * averaged using the weights given by Stetson. */
233  double cosz = muse_astro_get_zenith_distance(hourangle, delta, latitude);
234 #if MUSE_ASTRO_AIRMASS_APPROX == 2
235  double z = acos(cosz) * CPL_MATH_DEG_RAD;
236  const double zlimit = 80.;
237  if (z > zlimit) {
238  cpl_msg_warning(__func__, "Zenith angle %f > %f!", z, zlimit);
239  }
240 #endif
241  if (cosz == 0. || fabs(1. / cosz) < FLT_EPSILON ||
242  acos(cosz) > CPL_MATH_PI_2) {
243  cpl_msg_error(__func__, "Airmass computation unsuccessful. Object is below "
244  "the horizon at start (z = %f).", acos(cosz) * CPL_MATH_DEG_RAD);
245  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_OUTPUT);
246  return -1.;
247  }
248 
249  double airmass = 0.;
250 #if MUSE_ASTRO_AIRMASS_APPROX == 0
251  airmass = muse_astro_get_airmass_youngirvine(1. / cosz);
252 #endif
253 #if MUSE_ASTRO_AIRMASS_APPROX == 1
254  airmass = muse_astro_get_airmass_young(cosz);
255 #endif
256 #if MUSE_ASTRO_AIRMASS_APPROX == 2
257  airmass = muse_astro_get_airmass_hardie(1. / cosz);
258 #endif
259 #if MUSE_ASTRO_AIRMASS_APPROX > 2
260 #error set MUSE_ASTRO_AIRMASS_APPROX to 0, 1, 2
261 #endif
262 
263  /* if the exposure time is larger than zero, weight in airmass *
264  * at mid and end exposure according to Stetson's weights (which *
265  * are also used in IRAF/astcalc for the eairmass function */
266  if (aExptime > 0.) {
267  const double weights[] = {1./6., 2./3., 1./6.};
268  const int nweights = sizeof(weights) / sizeof(double);
269 
270  double timeStep = aExptime / (nweights - 1) * 15./3600. * CPL_MATH_RAD_DEG;
271  airmass *= weights[0];
272 
273  int i;
274  for (i = 1; i < nweights; i++) {
275  cosz = muse_astro_get_zenith_distance(hourangle + i * timeStep, delta, latitude);
276 #if MUSE_ASTRO_AIRMASS_APPROX == 2
277  z = acos(cosz) * CPL_MATH_DEG_RAD;
278  if (z > zlimit) {
279  cpl_msg_warning(__func__, "Zenith angle %f > %f!", z, zlimit);
280  }
281 #endif
282  if (cosz == 0. || fabs(1. / cosz) < FLT_EPSILON ||
283  acos(cosz) > CPL_MATH_PI_2) {
284  cpl_msg_error(__func__, "Airmass computation unsuccessful at timeStep. "
285  "Object is below the horizon at %s exposure (z=%f).",
286  i == 1 ? "mid" : "end", acos(cosz) * CPL_MATH_DEG_RAD);
287  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_OUTPUT);
288  return -1.;
289  }
290 
291 #if MUSE_ASTRO_AIRMASS_APPROX == 0
292  airmass += weights[i] * muse_astro_get_airmass_youngirvine(1. / cosz);
293 #endif
294 #if MUSE_ASTRO_AIRMASS_APPROX == 1
295  airmass += weights[i] * muse_astro_get_airmass_young(cosz);
296 #endif
297 #if MUSE_ASTRO_AIRMASS_APPROX == 2
298  airmass += weights[i] * muse_astro_get_airmass_hardie(1. / cosz);
299 #endif
300  } /* for i (weights) */
301  } /* if aExptime */
302 
303 #if MUSE_ASTRO_AIRMASS_APPROX == 0
304  /* Accuracy limit for airmass approximation of Young & Irvine */
305  const double airmasslimit = 4.;
306  if (airmass > airmasslimit) {
307  cpl_msg_warning(__func__, "Airmass larger than %f", airmasslimit);
308  }
309 #endif
310 
311  return airmass;
312 } /* muse_astro_compute_airmass() */
313 
314 /*----------------------------------------------------------------------------*/
329 /*----------------------------------------------------------------------------*/
330 double
331 muse_astro_airmass(cpl_propertylist *aHeader)
332 {
333  cpl_ensure(aHeader, CPL_ERROR_NULL_INPUT, -1.);
334 
335  cpl_errorstate prestate = cpl_errorstate_get();
336  double airm1 = muse_pfits_get_airmass_start(aHeader),
337  airm2 = muse_pfits_get_airmass_end(aHeader);
338  cpl_errorstate_set(prestate); /* we don't really care if these are missing */
339 
340  /* Try to read all necessary FITS headers. Reset them to something invalid *
341  * if they cannot be read, so that we get a real error message from *
342  * muse_astro_compute_airmass(), just for GEOLAT we don't care, because *
343  * muse_pfits knows a good default. */
344  prestate = cpl_errorstate_get();
345  double ra = muse_pfits_get_ra(aHeader);
346  if (!cpl_errorstate_is_equal(prestate)) {
347  ra = -1000.;
348  }
349  prestate = cpl_errorstate_get();
350  double dec = muse_pfits_get_dec(aHeader);
351  if (!cpl_errorstate_is_equal(prestate)) {
352  dec = -1000.;
353  }
354  prestate = cpl_errorstate_get();
355  double lst = muse_pfits_get_lst(aHeader);
356  if (!cpl_errorstate_is_equal(prestate)) {
357  lst = -1000.;
358  }
359  prestate = cpl_errorstate_get();
360  double exptime = muse_pfits_get_exptime(aHeader);
361  if (!cpl_errorstate_is_equal(prestate)) {
362  exptime = -1.;
363  }
364 
365  double airmass = muse_astro_compute_airmass(ra, dec, lst, exptime,
366  muse_pfits_get_geolat(aHeader));
367 
368  /* airmass computation failed, use simple average */
369  if (airmass < 0. && airm1 != 0. && airm2 != 0.) {
370  /* comparison doesn't make sense in this case, warn and return immediately */
371  double average = (airm1 + airm2) / 2.;
372  cpl_msg_warning(__func__, "airmass computation unsuccessful (%s), using simple "
373  "average of start and end values (%f)", cpl_error_get_message(),
374  average);
375  return average;
376  }
377 
378  /* compare computed and header values */
379  cpl_msg_debug(__func__, "airmass=%f (header %f, %f)", airmass, airm1, airm2);
380  if (airm1 != 0. && airm2 != 0.) {
381  cpl_boolean check = airmass > fmin(airm1, airm2) - 0.005
382  && airmass < fmax(airm1, airm2) + 0.005;
383  if (!check) {
384  cpl_msg_warning(__func__, "Computed airmass %.3f is NOT in the range "
385  "recorded in the FITS header (%f, %f)", airmass,
386  airm1, airm2);
387  }
388  }
389 
390  return airmass;
391 } /* muse_astro_airmass() */
392 
393 /*----------------------------------------------------------------------------*/
411 /*----------------------------------------------------------------------------*/
412 double
413 muse_astro_posangle(const cpl_propertylist *aHeader)
414 {
415  cpl_ensure(aHeader, CPL_ERROR_NULL_INPUT, 0.);
416  double posang = muse_pfits_get_drot_posang(aHeader);
417  const char *mode = muse_pfits_get_drot_mode(aHeader);
418  if (mode && !strncmp(mode, "SKY", 4)) {
419  posang *= -1.;
420  } else if (mode && strncmp(mode, "STAT", 5)) { /* not STAT, which we ignore */
421  cpl_msg_warning(__func__, "Derotator mode is neither SKY nor STAT! "
422  "Effective position angle may be wrong!");
423  } else if (!mode) { /* MODE missing */
424  cpl_msg_warning(__func__, "Derotator mode is not given! Effective position "
425  "angle may be wrong!");
426  }
427  return posang;
428 } /* muse_astro_posangle() */
429 
430 /*----------------------------------------------------------------------------*/
447 /*----------------------------------------------------------------------------*/
448 double
449 muse_astro_parangle(const cpl_propertylist *aHeader)
450 {
451  cpl_ensure(aHeader, CPL_ERROR_NULL_INPUT, 0.);
452  cpl_errorstate es = cpl_errorstate_get();
453  double p1 = muse_pfits_get_parang_start(aHeader),
454  p2 = muse_pfits_get_parang_end(aHeader);
455  if (!cpl_errorstate_is_equal(es)) {
456  cpl_msg_warning(__func__, "One or both TEL.PARANG keywords are missing!");
457  }
458 
459  double parang = (p1 + p2) / 2.;
460  if (fabs(p1 - p2) < 90.) {
461  /* are not going through the meridian with a 360 deg flip in *
462  * PARANG, so it should be safe to just average the two values */
463  return parang;
464  }
465  /* Flip in PARANG while observing this object, special (pretty *
466  * convoluted) handling needed: Both absolute values should be close to *
467  * 180., so compute the absolute distance from 180, the average of the *
468  * distance to 180 with the respective sign, and subtract the average *
469  * of those 180. The final sign is taken from the value whose distance *
470  * to the from respectively signed 180 is larger. */
471  double d1 = copysign(180. - fabs(p1), p1),
472  d2 = copysign(180. - fabs(p2), p2);
473  parang = 180. - fabs((d1 + d2) / 2.);
474  if (fabs(d1) > fabs(d2)) {
475  parang = copysign(parang, p1);
476  } else {
477  parang = copysign(parang, p2);
478  } /* else */
479  return parang;
480 } /* muse_astro_parangle() */
481 
482 /*----------------------------------------------------------------------------*/
494 /*----------------------------------------------------------------------------*/
495 double
496 muse_astro_angular_distance(double aRA1, double aDEC1, double aRA2, double aDEC2)
497 {
498  /* first convert all inputs to radians */
499  double ra1 = aRA1 * CPL_MATH_RAD_DEG,
500  dec1 = aDEC1 * CPL_MATH_RAD_DEG,
501  ra2 = aRA2 * CPL_MATH_RAD_DEG,
502  dec2 = aDEC2 * CPL_MATH_RAD_DEG;
503  /* compute formula components */
504  double dra = ra2 - ra1,
505  cosdra = cos(dra),
506  sindec1 = sin(dec1),
507  cosdec1 = cos(dec1),
508  sindec2 = sin(dec2),
509  cosdec2 = cos(dec2);
510  /* compute the two terms in the nominator */
511  double t1 = cosdec2 * sin(dra),
512  t2 = cosdec1 * sindec2 - sindec1 * cosdec2 * cosdra,
513  dsigma = atan2(sqrt(t1*t1 + t2*t2),
514  sindec1 * sindec2 + cosdec1 * cosdec2 * cosdra);
515  return dsigma * CPL_MATH_DEG_RAD;
516 } /* muse_astro_angular_distance() */
517 
518 /*----------------------------------------------------------------------------*/
532 /*----------------------------------------------------------------------------*/
533 double
535 {
536  return aVac / (1. + 2.735182e-4 + 131.4182 / pow(aVac, 2)
537  + 2.76249e8 / pow(aVac, 4));
538 } /* muse_astro_wavelength_vacuum_to_air() */
539 
540 /*----------------------------------------------------------------------------*/
567 /*----------------------------------------------------------------------------*/
568 static cpl_matrix *
569 muse_astro_create_rotation(const char *aMode, double aPhi, double aTheta,
570  double aPsi)
571 {
572  int nrotations = strlen(aMode);
573  nrotations = nrotations <= 3 ? nrotations : 3;
574  cpl_matrix *self = cpl_matrix_new(3, 3);
575 
576  /* Initialize the rotation matrix as identity */
577  cpl_matrix_fill_diagonal(self, 1., 0);
578 
579  if (!nrotations) {
580  return self;
581  }
582 
583  /* Create the general rotation matrix by successively left *
584  * multiplying the rotation matrices about the given axes. */
585  cpl_matrix *R = cpl_matrix_new(3, 3);
586  double angle[3] = { aPhi, aTheta, aPsi },
587  *_R = cpl_matrix_get_data(R);
588  int i;
589  for (i = 0; i < nrotations; ++i) {
590  double sa = sin(angle[i]),
591  ca = cos(angle[i]);
592  cpl_matrix_fill(R, 0.);
593  switch (tolower(aMode[i])) {
594  case 'x':
595  _R[0] = 1.;
596  _R[4] = ca;
597  _R[5] = sa;
598  _R[7] = -sa;
599  _R[8] = ca;
600  break;
601  case 'y':
602  _R[0] = ca;
603  _R[2] = -sa;
604  _R[4] = 1.;
605  _R[6] = sa;
606  _R[8] = ca;
607  break;
608  case 'z':
609  _R[0] = ca;
610  _R[1] = sa;
611  _R[3] = -sa;
612  _R[4] = ca;
613  _R[8] = 1.;
614  break;
615  default:
616  /* Invalid axis designation */
617  cpl_matrix_delete(R);
618  cpl_matrix_delete(self);
619  return NULL;
620  } /* switch */
621  cpl_matrix *mt = cpl_matrix_product_create(R, self);
622  cpl_matrix_delete(self);
623  self = mt;
624  } /* for i */
625  cpl_matrix_delete(R);
626  return self;
627 } /* muse_astro_create_rotation() */
628 
629 /*----------------------------------------------------------------------------*/
645 /*----------------------------------------------------------------------------*/
646 static cpl_matrix *
647 muse_astro_precession_matrix(double aEpoch1, double aEpoch2)
648 {
649  /* shortcuts to be in Simon et al. notation */
650  double sigma0 = 2000., /* the basic epoch of J2000 */
651  sigmaF = aEpoch1, /* the arbitrary fixed epoch */
652  sigmaD = aEpoch2; /* the epoch of the date */
653  /* compute the time differences, in units of 1000 years */
654  double T = (sigmaF - sigma0) / 1000., /* JD(sigmaF) - JD(sigma0) / 365250 */
655  t = (sigmaD - sigmaF) / 1000.; /* JD(sigmaD) - JD(sigmaF) / 365250 */
656  /* the angles, in units of arcsec */
657  double thetaA = (20042.0207 - 85.3131 * T - 0.2111 * T*T
658  + 0.3642 * T*T*T + 0.0008 * T*T*T*T - 0.0005 * T*T*T*T*T) * t
659  + (-42.6566 - 0.2111 * T + 0.5463 * T*T
660  + 0.0017 * T*T*T - 0.0012 * T*T*T*T) * t*t
661  + (-41.8238 + 0.0359 * T + 0.0027 *T*T - 0.0001 * T*T*T) * t*t*t
662  + (-0.0731 + 0.0019 * T + 0.0009 *T*T) * t*t*t*t
663  + (-0.0127 + 0.0011 * T) * t*t*t*t*t
664  + (0.0004) * t*t*t*t*t*t,
665  zetaA = (23060.9097 + 139.7495 * T - 0.0038 * T*T - 0.5918 * T*T*T
666  - 0.0037 * T*T*T*T + 0.0007 * T*T*T*T*T) * t
667  + (30.2226 - 0.2523 * T - 0.3840 * T*T - 0.0014 *T*T*T
668  + 0.0007 * T*T*T*T) * t*t
669  + (18.0183 - 0.1326 * T + 0.0006 * T*T + 0.0005 * T*T*T) * t*t*t
670  + (-0.0583 - 0.0001 * T + 0.0007 * T*T) * t*t*t*t
671  + (-0.0285) * t*t*t*t*t
672  + (-0.0002) * t*t*t*t*t*t,
673  zA = (23060.9097 +139.7495 * T -0.0038 * T*T -0.5918 * T*T*T
674  - 0.0037 * T*T*T*T + 0.0007 * T*T*T*T*T) * t
675  + (109.5270 + 0.2446 * T - 1.3913 * T*T - 0.0134 * T*T*T
676  + 0.0026 * T*T*T*T) * t*t
677  + (18.2667 - 1.1400 * T - 0.0173 * T*T + 0.0044 * T*T*T) * t*t*t
678  + (-0.2821 - 0.0093 * T + 0.0032 * T*T) * t*t*t*t
679  + (-0.0301 + 0.0006 * T) * t*t*t*t*t
680  + (-0.0001) * t*t*t*t*t*t;
681  /* convert the three angles to radians */
682  thetaA *= CPL_MATH_RAD_DEG / 3600.;
683  zetaA *= CPL_MATH_RAD_DEG / 3600.;
684  zA *= CPL_MATH_RAD_DEG / 3600.;
685 
686  /* create the matrix */
687  return muse_astro_create_rotation("zyz", -zetaA, thetaA, -zA);
688 } /* muse_astro_precession_matrix() */
689 
690 /*----------------------------------------------------------------------------*/
703 /*----------------------------------------------------------------------------*/
704 static double
705 muse_astro_sidereal_time(double aJD, double aLong)
706 {
707  /* Constants D1,D2,D3 for calculating Greenwich *
708  * Mean Sidereal Time at 0 UT */
709  const double d1 = 1.739935934667999,
710  d2 = 6.283319509909095e02,
711  d3 = 6.755878646261384e-06,
712  df = 1.00273790934;
713  double djd0 = floor(aJD) + 0.5;
714  if (djd0 > aJD) {
715  djd0 -= 1.;
716  }
717  double dut = (aJD - djd0) * CPL_MATH_2PI,
718  dt = (djd0 - kT0) / kJulianCentury,
719  dst0 = d1 + d2 * dt + d3 * dt * dt;
720  dst0 = fmod(dst0, CPL_MATH_2PI);
721  double dst = df * dut + dst0 - aLong;
722  dst = fmod(dst + 2.*CPL_MATH_2PI, CPL_MATH_2PI);
723  return dst;
724 } /* muse_astro_sidereal_time() */
725 
726 /*----------------------------------------------------------------------------*/
754 /*----------------------------------------------------------------------------*/
755 static double
756 muse_astro_geo_correction(double aLat, double aElev, double aDEC, double aHA)
757 {
758  /* Earth's equatorial radius (km) (Geod. ref. sys., 1980) */
759  const double da = 6378.137;
760 
761  /* Polar flattening (Geod. ref. sys., 1980) */
762  const double df = 1. / 298.257222;
763 
764  /* Angular rotation rate (2.pi/t) */
765  const double dw = CPL_MATH_2PI / 86164.;
766  const double de2 = df * (2.0 - df),
767  dsdlats = sin(aLat) * sin(aLat);
768 
769  /* Calculate geocentric radius dr0 at sea level (km) */
770  double d1 = 1.0 - de2 * (2.0 - de2) * dsdlats,
771  d2 = 1.0 - de2 * dsdlats,
772  dr0 = da * sqrt(d1 / d2);
773 
774  /* Calculate geocentric latitude dlatg */
775  d1 = de2 * sin(2.0 * aLat);
776  d2 = 2.0 * d2;
777  double dlatg = aLat - atan(d1 / d2);
778 
779  /* Calculate geocentric radius drh at elevation aElev (km) */
780  double drh = dr0 * cos(dlatg) + (aElev / 1000.) * cos(aLat);
781 
782  /* Projected component to star at declination = aDEC and *
783  * at hour angle = aHA, in units of km/s */
784  return dw * drh * cos(aDEC) * sin(aHA);
785 } /* muse_astro_geo_correction() */
786 
787 /*----------------------------------------------------------------------------*/
814 /*----------------------------------------------------------------------------*/
815 static void
816 muse_astro_earth_velocity(double aDJE, double aDEqu, double* const aHVel,
817  double* const aBVel)
818 {
819  /* Constants dcfel(i, k) of fast-changing elements *
820  * i = 1 i = 2 i = 3 */
821  const double dcfel[][3] = {
822  { 1.7400353e+00, 6.2833195099091e+02, 5.2796e-06 },
823  { 6.2565836e+00, 6.2830194572674e+02, -2.6180e-06 },
824  { 4.7199666e+00, 8.3997091449254e+03, -1.9780e-05 },
825  { 1.9636505e-01, 8.4334662911720e+03, -5.6044e-05 },
826  { 4.1547339e+00, 5.2993466764997e+01, 5.8845e-06 },
827  { 4.6524223e+00, 2.1354275911213e+01, 5.6797e-06 },
828  { 4.2620486e+00, 7.5025342197656e+00, 5.5317e-06 },
829  { 1.4740694e+00, 3.8377331909193e+00, 5.6093e-06 }
830  };
831 
832  /* Constants dceps and ccsel(i,k) of slowly changing elements *
833  * i = 1 i = 2 i = 3 */
834  const double dceps[3] = {
835  4.093198e-01, -2.271110e-04, -2.860401e-08
836  };
837  const double ccsel[][3] = {
838  { 1.675104e-02, -4.179579e-05, -1.260516e-07 },
839  { 2.220221e-01, 2.809917e-02, 1.852532e-05 },
840  { 1.589963e+00, 3.418075e-02, 1.430200e-05 },
841  { 2.994089e+00, 2.590824e-02, 4.155840e-06 },
842  { 8.155457e-01, 2.486352e-02, 6.836840e-06 },
843  { 1.735614e+00, 1.763719e-02, 6.370440e-06 },
844  { 1.968564e+00, 1.524020e-02, -2.517152e-06 },
845  { 1.282417e+00, 8.703393e-03, 2.289292e-05 },
846  { 2.280820e+00, 1.918010e-02, 4.484520e-06 },
847  { 4.833473e-02, 1.641773e-04, -4.654200e-07 },
848  { 5.589232e-02, -3.455092e-04, -7.388560e-07 },
849  { 4.634443e-02, -2.658234e-05, 7.757000e-08 },
850  { 8.997041e-03, 6.329728e-06, -1.939256e-09 },
851  { 2.284178e-02, -9.941590e-05, 6.787400e-08 },
852  { 4.350267e-02, -6.839749e-05, -2.714956e-07 },
853  { 1.348204e-02, 1.091504e-05, 6.903760e-07 },
854  { 3.106570e-02, -1.665665e-04, -1.590188e-07 }
855  };
856 
857  /* Constants of the arguments of the short-period perturbations by *
858  * the planets: dcargs(i, k) *
859  * i = 1 i = 2 */
860  const double dcargs[][2] = {
861  { 5.0974222e+00, -7.8604195454652e+02 },
862  { 3.9584962e+00, -5.7533848094674e+02 },
863  { 1.6338070e+00, -1.1506769618935e+03 },
864  { 2.5487111e+00, -3.9302097727326e+02 },
865  { 4.9255514e+00, -5.8849265665348e+02 },
866  { 1.3363463e+00, -5.5076098609303e+02 },
867  { 1.6072053e+00, -5.2237501616674e+02 },
868  { 1.3629480e+00, -1.1790629318198e+03 },
869  { 5.5657014e+00, -1.0977134971135e+03 },
870  { 5.0708205e+00, -1.5774000881978e+02 },
871  { 3.9318944e+00, 5.2963464780000e+01 },
872  { 4.8989497e+00, 3.9809289073258e+01 },
873  { 1.3097446e+00, 7.7540959633708e+01 },
874  { 3.5147141e+00, 7.9618578146517e+01 },
875  { 3.5413158e+00, -5.4868336758022e+02 }
876  };
877 
878  /* Amplitudes ccamps(n, k) of the short-period perturbations *
879  * n = 1 n = 2 n = 3 n = 4 n = 5 */
880  const double ccamps[][5] = {
881  { -2.279594e-5, 1.407414e-5, 8.273188e-6, 1.340565e-5, -2.490817e-7 },
882  { -3.494537e-5, 2.860401e-7, 1.289448e-7, 1.627237e-5, -1.823138e-7 },
883  { 6.593466e-7, 1.322572e-5, 9.258695e-6, -4.674248e-7, -3.646275e-7 },
884  { 1.140767e-5, -2.049792e-5, -4.747930e-6, -2.638763e-6, -1.245408e-7 },
885  { 9.516893e-6, -2.748894e-6, -1.319381e-6, -4.549908e-6, -1.864821e-7 },
886  { 7.310990e-6, -1.924710e-6, -8.772849e-7, -3.334143e-6, -1.745256e-7 },
887  { -2.603449e-6, 7.359472e-6, 3.168357e-6, 1.119056e-6, -1.655307e-7 },
888  { -3.228859e-6, 1.308997e-7, 1.013137e-7, 2.403899e-6, -3.736225e-7 },
889  { 3.442177e-7, 2.671323e-6, 1.832858e-6, -2.394688e-7, -3.478444e-7 },
890  { 8.702406e-6, -8.421214e-6, -1.372341e-6, -1.455234e-6, -4.998479e-8 },
891  { -1.488378e-6, -1.251789e-5, 5.226868e-7, -2.049301e-7, 0.0e0 },
892  { -8.043059e-6, -2.991300e-6, 1.473654e-7, -3.154542e-7, 0.0e0 },
893  { 3.699128e-6, -3.316126e-6, 2.901257e-7, 3.407826e-7, 0.0e0 },
894  { 2.550120e-6, -1.241123e-6, 9.901116e-8, 2.210482e-7, 0.0e0 },
895  { -6.351059e-7, 2.341650e-6, 1.061492e-6, 2.878231e-7, 0.0e0 }
896  };
897 
898  /* Constants of the secular perturbations in longitude *
899  * ccsec3 and ccsec(n,k) *
900  * n = 1 n = 2 n = 3 */
901  const double ccsec3 = -7.757020e-08,
902  ccsec[][3] = {
903  { 1.289600e-06, 5.550147e-01, 2.076942e+00 },
904  { 3.102810e-05, 4.035027e+00, 3.525565e-01 },
905  { 9.124190e-06, 9.990265e-01, 2.622706e+00 },
906  { 9.793240e-07, 5.508259e+00, 1.559103e+01 }
907  };
908 
909  /* Sidereal rate dcsld in longitude, rate ccsgd in mean anomaly */
910  const double dcsld = 1.990987e-07,
911  ccsgd = 1.990969e-07;
912 
913  /* Some constants used in the calculation of the *
914  * lunar contribution */
915  const double cckm = 3.122140e-05,
916  ccmld = 2.661699e-06,
917  ccfdi = 2.399485e-07;
918 
919  /* Constants dcargm(i,k) of the arguments of the *
920  * perturbations of the motion of the moon *
921  * i = 1 i = 2 */
922  const double dcargm[][2] = {
923  { 5.1679830e+00, 8.3286911095275e+03 },
924  { 5.4913150e+00, -7.2140632838100e+03 },
925  { 5.9598530e+00, 1.5542754389685e+04 }
926  };
927 
928  /* Amplitudes ccampm(n,k) of the perturbations of the moon *
929  * n = 1 n = 2 n = 3 n = 4 */
930  const double ccampm[][4] = {
931  { 1.097594e-01, 2.896773e-07, 5.450474e-02, 1.438491e-07 },
932  {-2.223581e-02, 5.083103e-08, 1.002548e-02, -2.291823e-08 },
933  { 1.148966e-02, 5.658888e-08, 8.249439e-03, 4.063015e-08 }
934  };
935 
936  /* ccpamv = a * m * dl/dt (planets) */
937  const double ccpamv[4] = { 8.326827e-11, 1.843484e-11,
938  1.988712e-12, 1.881276e-12 };
939 
940  /* dc1mme = 1 - mass(earth + moon) */
941  const double dc1mme = 0.99999696;
942 
943  /* Control-parameter ideq, and time-arguments */
944  int ideq = (int)aDEqu;
945  double dt = (aDJE - kT0) / kJulianCentury,
946  t = dt,
947  dtsq = dt * dt,
948  tsq = dtsq;
949 
950  /* Values of all elements for the instant aDJE */
951  double dml = 0.,
952  forbel[7] = { 0., 0., 0., 0., 0., 0., 0. };
953  int k;
954  for (k = 0; k < 8; k++) {
955  double dlocal = fmod(dcfel[k][0] + dt * dcfel[k][1] + dtsq * dcfel[k][2],
956  CPL_MATH_2PI);
957  if (k == 0) {
958  dml = dlocal;
959  }
960  if (k != 0) {
961  forbel[k - 1] = dlocal;
962  }
963  } /* for k */
964  double deps = fmod(dceps[0] + dt * dceps[1] + dtsq * dceps[2], CPL_MATH_2PI),
965  sorbel[17];
966  for (k = 0; k < 17; k++) {
967  sorbel[k] = fmod(ccsel[k][0] + t * ccsel[k][1] + tsq * ccsel[k][2],
968  CPL_MATH_2PI);
969  } /* for k */
970 
971  /* Secular perturbations in longitude */
972  double sn[4];
973  for (k = 0; k < 4; k++) {
974  double a = fmod(ccsec[k][1] + t * ccsec[k][2], CPL_MATH_2PI);
975  sn[k] = sin(a);
976  } /* for k */
977 
978  /* Periodic perturbations of the earth-moon barycenter */
979  double pertl = ccsec[0][0] * sn[0] + ccsec[1][0] * sn[1]
980  + (ccsec[2][0] + t * ccsec3) * sn[2] + ccsec[3][0] * sn[3],
981  pertld = 0.,
982  pertr = 0.,
983  pertrd = 0.;
984  for (k = 0; k < 15; k++) {
985  double a = fmod(dcargs[k][0] + dt * dcargs[k][1], CPL_MATH_2PI);
986  pertl += (ccamps[k][0] * cos(a) + ccamps[k][1] * sin(a));
987  pertr += (ccamps[k][2] * cos(a) + ccamps[k][3] * sin(a));
988  if (k >= 10) {
989  continue;
990  }
991  pertld += ((ccamps[k][1] * cos(a) - ccamps[k][0] * sin(a)) * ccamps[k][4]);
992  pertrd += ((ccamps[k][3] * cos(a) - ccamps[k][2] * sin(a)) * ccamps[k][4]);
993  } /* for k */
994 
995  /* Elliptic part of the motion of the earth-moon barycenter */
996  double esq = sorbel[0] * sorbel[0],
997  dparam = 1. - esq,
998  param = dparam,
999  twoe = sorbel[0] + sorbel[0],
1000  twog = forbel[0] + forbel[0],
1001  phi = twoe * ((1.0 - esq * (1.0 / 8.0)) * sin(forbel[0])
1002  + sorbel[0] * (5.0 / 8.0) * sin(twog)
1003  + esq * 0.5416667 * sin(forbel[0] + twog)),
1004  f = forbel[0] + phi,
1005  dpsi = dparam / (1. + sorbel[0] * cos(f)),
1006  phid = twoe * ccsgd * ((1.0 + esq * 1.50) * cos(f)
1007  + sorbel[0] * (1.250 - sin(f) * sin(f) * 0.50)),
1008  psid = ccsgd * sorbel[0] * sin(f) / sqrt(param);
1009 
1010  /* Perturbed heliocentric motion of the earth-moon barycenter */
1011  double d1pdro = 1. + pertr,
1012  drd = d1pdro * (psid + dpsi * pertrd),
1013  drld = d1pdro * dpsi * (dcsld + phid + pertld),
1014  dtl = fmod(dml + phi + pertl, CPL_MATH_2PI),
1015  dxhd = drd * cos(dtl) - drld * sin(dtl),
1016  dyhd = drd * sin(dtl) + drld * cos(dtl);
1017 
1018  /* Influence of eccentricity, evection and variation of *
1019  * the geocentric motion of the moon */
1020  pertl = 0.;
1021  pertld = 0.;
1022  double pertp = 0.,
1023  pertpd = 0.;
1024  for (k = 0; k < 3; k++) {
1025  double a = fmod(dcargm[k][0] + dt * dcargm[k][1], CPL_MATH_2PI);
1026  pertl += ccampm[k][0] * sin(a);
1027  pertld += ccampm[k][1] * cos(a);
1028  pertp += ccampm[k][2] * cos(a);
1029  pertpd -= ccampm[k][3] * sin(a);
1030  } /* for k */
1031 
1032  /* Heliocentric motion of the earth */
1033  double tl = forbel[1] + pertl,
1034  sigma = cckm / (1. + pertp),
1035  a = sigma * (ccmld + pertld),
1036  b = sigma * pertpd;
1037  dxhd = dxhd + a * sin(tl) + b * cos(tl);
1038  dyhd = dyhd - a * cos(tl) + b * sin(tl);
1039  double dzhd = -sigma * ccfdi * cos(forbel[2]);
1040 
1041  /* Barycentric motion of the earth */
1042  double dxbd = dxhd * dc1mme,
1043  dybd = dyhd * dc1mme,
1044  dzbd = dzhd * dc1mme;
1045  for (k = 0; k < 4; k++) {
1046  double plon = forbel[k + 3],
1047  pomg = sorbel[k + 1],
1048  pecc = sorbel[k + 9];
1049  tl = fmod(plon + 2.0 * pecc * sin(plon - pomg), CPL_MATH_2PI);
1050  dxbd = dxbd + ccpamv[k] * (sin(tl) + pecc * sin(pomg));
1051  dybd = dybd - ccpamv[k] * (cos(tl) + pecc * cos(pomg));
1052  dzbd = dzbd - ccpamv[k] * sorbel[k + 13] * cos(plon - sorbel[k + 5]);
1053  } /* for k */
1054 
1055  /* Transition to mean equator of date */
1056  double dyahd = cos(deps) * dyhd - sin(deps) * dzhd,
1057  dzahd = sin(deps) * dyhd + cos(deps) * dzhd,
1058  dyabd = cos(deps) * dybd - sin(deps) * dzbd,
1059  dzabd = sin(deps) * dybd + cos(deps) * dzbd;
1060  if (ideq == 0) {
1061  aHVel[0] = dxhd;
1062  aHVel[1] = dyahd;
1063  aHVel[2] = dzahd;
1064  aBVel[0] = dxbd;
1065  aBVel[1] = dyabd;
1066  aBVel[2] = dzabd;
1067  } else {
1068  /* General precession from epoch aDJE to aDEqu */
1069  double deqdat = (aDJE - kT0 - kBesselOffset) / kTropYear + kEpoch1900;
1070  cpl_matrix* prec = muse_astro_precession_matrix(deqdat, aDEqu);
1071  int n;
1072  for (n = 0; n < 3; n++) {
1073  aHVel[n] = dxhd * cpl_matrix_get(prec, 0, n)
1074  + dyahd * cpl_matrix_get(prec, 1, n)
1075  + dzahd * cpl_matrix_get(prec, 2, n);
1076  aBVel[n] = dxbd * cpl_matrix_get(prec, 0, n)
1077  + dyabd * cpl_matrix_get(prec, 1, n)
1078  + dzabd * cpl_matrix_get(prec, 2, n);
1079  }
1080  cpl_matrix_delete(prec);
1081  } /* else */
1082 } /* muse_astro_earth_velocity() */
1083 
1084 /*----------------------------------------------------------------------------*/
1109 /*----------------------------------------------------------------------------*/
1110 static void
1111 muse_rvcorrection_compute(muse_astro_rvcorr *aRV, double aJD, double aLong,
1112  double aLat, double aElev, double aRA, double aDEC,
1113  double aEqui)
1114 {
1115  /* Precess RA and DEC to mean equator and equinox of date */
1116  double eqt = (aJD - kT0 - kBesselOffset) / kTropYear + kEpoch1900,
1117  dc[3] = { cos(aRA * 15.0 * CPL_MATH_RAD_DEG) * cos(aDEC * CPL_MATH_RAD_DEG),
1118  sin(aRA * 15.0 * CPL_MATH_RAD_DEG) * cos(aDEC * CPL_MATH_RAD_DEG),
1119  sin(aDEC * CPL_MATH_RAD_DEG) };
1120  cpl_matrix *precession = muse_astro_precession_matrix(aEqui, eqt);
1121  double dcc[3];
1122  int i;
1123  for (i = 0; i < 3; ++i) {
1124  dcc[i] = dc[0] * cpl_matrix_get(precession, i, 0)
1125  + dc[1] * cpl_matrix_get(precession, i, 1)
1126  + dc[2] * cpl_matrix_get(precession, i, 2);
1127  }
1128  cpl_matrix_delete(precession);
1129  precession = NULL;
1130 
1131  double ra2 = 0.,
1132  dec2 = asin(dcc[2]);
1133  if (dcc[0] != 0.) {
1134  ra2 = atan(dcc[1] / dcc[0]);
1135  if (dcc[0] < 0.) {
1136  ra2 += CPL_MATH_PI;
1137  } else {
1138  if (dcc[1] < 0.) {
1139  ra2 += CPL_MATH_2PI;
1140  } /* if */
1141  } /* else */
1142  } else {
1143  if (dcc[1] > 0.) {
1144  ra2 = CPL_MATH_PI_2;
1145  } else {
1146  ra2 = 1.5 * CPL_MATH_PI;
1147  } /* else */
1148  } /* else */
1149 
1150  /* Calculate hour angle = local sidereal time - RA *
1151  * convert input values to units of radians */
1152  double st = muse_astro_sidereal_time(aJD, aLong * CPL_MATH_RAD_DEG),
1153  ha = st - ra2;
1154 
1155  /* Calculate observer's geocentric velocity *
1156  * (elevation assumed to be zero) */
1157  aRV->geo = muse_astro_geo_correction(aLat * CPL_MATH_RAD_DEG, aElev, dec2, -ha);
1158 
1159  /* Calculate components of earth's barycentric velocity, *
1160  * bvel(i), i=1,2,3 in units of au/s */
1161  double hv[3] = { 0., 0., 0. },
1162  bv[3] = { 0., 0., 0. };
1163  muse_astro_earth_velocity(aJD, eqt, hv, bv);
1164 
1165  /* Project barycentric velocity to the direction of the *
1166  * star, and convert to km/s */
1167  aRV->bary = 0.;
1168  aRV->helio = 0.;
1169  for (i = 0; i < 3; ++i) {
1170  aRV->bary += bv[i] * dcc[i] * kAU;
1171  aRV->helio += hv[i] * dcc[i] * kAU;
1172  } /* for i */
1173 } /* muse_rvcorrection_compute() */
1174 
1175 /*----------------------------------------------------------------------------*/
1195 /*----------------------------------------------------------------------------*/
1197 muse_astro_rvcorr_compute(const cpl_propertylist *aHeader)
1198 {
1199  muse_astro_rvcorr rvcorr = { 0., 0., 0. };
1200  cpl_ensure(aHeader, CPL_ERROR_NULL_INPUT, rvcorr);
1201 
1202  /* get all the necessary properties from the header, watching the state */
1203  cpl_errorstate state = cpl_errorstate_get();
1204  double exptime = muse_pfits_get_exptime(aHeader),
1205  jd = muse_pfits_get_mjdobs(aHeader),
1206  equinox = muse_pfits_get_equinox(aHeader),
1207  ra = muse_pfits_get_ra(aHeader) / 15., /* need RA in hours */
1208  dec = muse_pfits_get_dec(aHeader);
1209  /* Compute julian date of mid exposure. 2400000.5 is the offset between *
1210  * JD and MJD and corresponds to November 17th 1858 0:00 UT. */
1211  jd += 2400000.5 + 0.5 * exptime / (24. * 3600.);
1212  if (!cpl_errorstate_is_equal(state)) {
1213  cpl_error_set_message(__func__, CPL_ERROR_DATA_NOT_FOUND, "Could not find all"
1214  " properties necessary for radial velocity computation!");
1215  return rvcorr;
1216  }
1217  /* the following functions provide good defaults... */
1218  double lon = muse_pfits_get_geolon(aHeader),
1219  lat = muse_pfits_get_geolat(aHeader),
1220  elev = muse_pfits_get_geoelev(aHeader);
1221  if (!cpl_errorstate_is_equal(state)) { /* ... so ignore any errors */
1222  cpl_errorstate_set(state);
1223  }
1224 
1225  /* now call the function that originates from the GIRAFFE pipeline */
1226  muse_rvcorrection_compute(&rvcorr, jd, lon, lat, elev, ra, dec, equinox);
1227  /* add geocentric velocity to get total corrective velocities */
1228  rvcorr.bary = rvcorr.bary + rvcorr.geo;
1229  rvcorr.helio = rvcorr.helio + rvcorr.geo;
1230  return rvcorr;
1231 } /* muse_astro_rvcorr_compute() */
1232 
double muse_pfits_get_ra(const cpl_propertylist *aHeaders)
find out the right ascension
Definition: muse_pfits.c:267
double muse_pfits_get_airmass_start(const cpl_propertylist *aHeaders)
find out the airmass at start of exposure
Definition: muse_pfits.c:970
Structure to store bary-, helio-, and geocentric velocity corrections.
Definition: muse_astro.h:43
double muse_astro_parangle(const cpl_propertylist *aHeader)
Properly average parallactic angle values of one exposure.
Definition: muse_astro.c:449
double muse_pfits_get_drot_posang(const cpl_propertylist *aHeaders)
find out the MUSE derotator position angle (in degrees)
Definition: muse_pfits.c:1296
double muse_pfits_get_equinox(const cpl_propertylist *aHeaders)
find out the equinox
Definition: muse_pfits.c:306
double muse_pfits_get_geolat(const cpl_propertylist *aHeaders)
find out the telescope&#39;s latitude
Definition: muse_pfits.c:882
double muse_astro_wavelength_vacuum_to_air(double aVac)
Compute air wavelength for a given vacuum wavelength.
Definition: muse_astro.c:534
const char * muse_pfits_get_drot_mode(const cpl_propertylist *aHeaders)
find out the MUSE derotator mode
Definition: muse_pfits.c:1278
double muse_astro_compute_airmass(double aRA, double aDEC, double aLST, double aExptime, double aLatitude)
Compute the effective airmass of an observation.
Definition: muse_astro.c:207
double muse_pfits_get_dec(const cpl_propertylist *aHeaders)
find out the declination
Definition: muse_pfits.c:285
muse_astro_rvcorr muse_astro_rvcorr_compute(const cpl_propertylist *aHeader)
Compute radial velocity corrections given the header of an exposure.
Definition: muse_astro.c:1197
double muse_pfits_get_parang_start(const cpl_propertylist *aHeaders)
find out the parallactic angle at start of exposure (in degrees)
Definition: muse_pfits.c:1134
double muse_astro_posangle(const cpl_propertylist *aHeader)
Derive the position angle of an observation from information in a FITS header.
Definition: muse_astro.c:413
double muse_pfits_get_airmass_end(const cpl_propertylist *aHeaders)
find out the airmass at end of exposure
Definition: muse_pfits.c:988
double muse_pfits_get_mjdobs(const cpl_propertylist *aHeaders)
find out the Julian Date of the observation
Definition: muse_pfits.c:346
double muse_pfits_get_exptime(const cpl_propertylist *aHeaders)
find out the exposure time
Definition: muse_pfits.c:382
double muse_pfits_get_geoelev(const cpl_propertylist *aHeaders)
find out the telescope&#39;s elevation
Definition: muse_pfits.c:928
double muse_pfits_get_lst(const cpl_propertylist *aHeaders)
find out the local siderial time
Definition: muse_pfits.c:328
double muse_pfits_get_geolon(const cpl_propertylist *aHeaders)
find out the telescope&#39;s longitude
Definition: muse_pfits.c:906
double muse_astro_angular_distance(double aRA1, double aDEC1, double aRA2, double aDEC2)
Compute angular distance in the sky between two positions.
Definition: muse_astro.c:496
double muse_pfits_get_parang_end(const cpl_propertylist *aHeaders)
find out the parallactic angle at end of exposure (in degrees)
Definition: muse_pfits.c:1152
double muse_astro_airmass(cpl_propertylist *aHeader)
Derive the effective airmass of an observation from information in a FITS header. ...
Definition: muse_astro.c:331