MUSE Pipeline Reference Manual  2.1.1
muse_phys.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) 2015 European Southern Observatory
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software
19  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
20  */
21 
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
25 
26 /*----------------------------------------------------------------------------*
27  * Includes *
28  *----------------------------------------------------------------------------*/
29 #include <math.h>
30 #include <string.h>
31 
32 #include "muse_phys.h"
33 
34 #include "muse_pfits.h"
35 
36 /*----------------------------------------------------------------------------*/
40 /*----------------------------------------------------------------------------*/
41 
44 /*----------------------------------------------------------------------------*/
52 /*----------------------------------------------------------------------------*/
53 double
55 {
56  return -10474.0 + 116.43*temp - 0.43284*temp*temp + 0.00053840*pow(temp,3);
57 }
58 
59 /*----------------------------------------------------------------------------*/
68 /*----------------------------------------------------------------------------*/
69 void
70 muse_phys_nrindex_owens_coeffs(double temp, double rhum, double pres,
71  double *d1, double *d2)
72 {
73  /* pressure helpers: saturation pressure */
75  p2 = rhum * ps, /* water vapor pressure */
76  p1 = pres - p2; /* dry air pressure */
77  /* helpers/coefficients */
78  *d1 = p1/temp * (1. + p1*(57.90e-8 - 9.3250e-4/temp + 0.25844/temp/temp));
79  *d2 = p2/temp * (1. + p2*(1. + 3.7e-4*p2)
80  * (-2.37321e-3 + 2.23366/temp - 710.792/temp/temp
81  + 7.75141e4/pow(temp,3)));
82 }
83 
84 /*----------------------------------------------------------------------------*/
93 /*----------------------------------------------------------------------------*/
94 double
95 muse_phys_nrindex_owens(double l, double d1, double d2)
96 {
97  double lisq = 1./(l*l); /* inverse square of the wavelength */
98  return 1. + ((2371.34 + 683939.7/(130. - lisq) + 4547.3/(38.9 - lisq)) * d1
99  + (6487.31 + 58.058 * lisq - 0.71150 * lisq*lisq
100  + 0.08851 * lisq*lisq*lisq) * d2) * 1.0e-8;
101 }
102 
103 /*----------------------------------------------------------------------------*/
113 /*----------------------------------------------------------------------------*/
114 double
115 muse_phys_nrindex_edlen(double l, double t, double p, double pv)
116 {
117  /* define constants */
118  const double A = 8342.54, B = 2406147,
119  C = 15998, D = 96095.43,
120  E = 0.601, F = 0.00972, G = 0.003661;
121  double S = 1./(l*l), /* supposed to be "laser vacuum wavelength", but for *
122  * MUSE the approximation is OK if we use air as input */
123  /* intermediate result for p, pv, and t */
124  n_5 = 1 + 1e-8 * (A + B / (130 - S) + C / (38.9 - S)),
125  X = (1 + 1e-8 * (E - F * t) * p) / (1 + G * t),
126  n_tp = 1 + p * (n_5 - 1) * X / D,
127  n = n_tp - 10e-10 * (292.75 / (t + 273.15)) * (3.7345 - 0.0401 * S) * pv;
128  return n;
129 } /* muse_phys_nrindex_edlen() */
130 
131 /*----------------------------------------------------------------------------*/
146 /*----------------------------------------------------------------------------*/
147 double
148 muse_phys_nrindex_ciddor(double l, double T, double p, double xv, double xCO2)
149 {
150  /* define constants */
151  const double w0 = 295.235, w1 = 2.6422, w2 = -0.03238, w3 = 0.004028,
152  k0 = 238.0185, k1 = 5792105, k2 = 57.362, k3 = 167917,
153  a0 = 1.58123e-6, a1 = -2.9331e-8, a2 = 1.1043e-10,
154  b0 = 5.707e-6, b1 = -2.051e-8,
155  c0 = 1.9898e-4, c1 = -2.376e-6,
156  d = 1.83e-11, e = -0.765e-8,
157  p_R1 = 101325, T_R1 = 288.15,
158  Z_a = 0.9995922115,
159  rho_vs = 0.00985938,
160  R = 8.314472, M_v = 0.018015;
161  double t = T - 273.15, /* temperature [Celsius] */
162  S = 1./(l*l), /* supposed to be "laser vacuum wavelength", but for *
163  * MUSE the approximation is OK if we use air as input */
164  /* intermediate results that depend on S */
165  r_as = 1e-8 * ((k1 / (k0 - S)) + (k3 / (k2 - S))),
166  r_vs = 1.022e-8 * (w0 + w1 * S + w2 * S*S + w3 * S*S*S),
167  /* using the CO2 concentration, calculate intermediate values for CO2 */
168  M_a = 0.0289635 + 1.2011e-8 * (xCO2 - 400),
169  r_axs = r_as * (1 + 5.34e-7 * (xCO2 - 450)),
170  /* compressibility and density components */
171  Z_m = 1 - (p / T) * (a0 + a1 * t + a2 * t*t + (b0 + b1 * t) * xv
172  + (c0 + c1 * t) * xv*xv)
173  + (p/T)*(p/T) * (d + e * xv*xv),
174  rho_axs = p_R1 * M_a / (Z_a * R * T_R1),
175  rho_v = xv * p * M_v / (Z_m * R * T),
176  rho_a = (1 - xv) * p * M_a / (Z_m * R * T),
177  n = 1 + (rho_a / rho_axs) * r_axs + (rho_v / rho_vs) * r_vs;
178  return n;
179 } /* muse_phys_nrindex_ciddor() */
180 
181 /*----------------------------------------------------------------------------*/
191 /*----------------------------------------------------------------------------*/
192 double
193 muse_phys_nrindex_filippenko(double l, double T, double f, double P)
194 {
195  double lisq = 1./(l*l); /* inverse square of the wavelength */
196  /* 10^6 [n(lambda) - 1] at standard environmental conditions, Eq. (1) */
197  double nl = 64.328 + 29498.1 / (146 - lisq) + 255.4 / (41 - lisq);
198  /* correction for non-standard conditions, Eq. (2) */
199  nl *= P / 720.883 * (1 + (1.049 - 0.0157 * T) * 1e-6 * P) / (1 + 0.003661 * T);
200  /* correction for water vapor, Eq. (3) */
201  nl -= (0.0624 - 0.000680 * lisq) / (1 + 0.003661 * T) * f;
202  /* convert to refractive index n(lambda) */
203  nl = nl * 1e-6 + 1;
204  return nl;
205 }
206 
207 /*----------------------------------------------------------------------------*/
259 /*----------------------------------------------------------------------------*/
260 cpl_error_code
261 muse_phys_air_to_vacuum(muse_pixtable *aPixtable, unsigned int aFlags)
262 {
263  cpl_ensure_code(aPixtable && aPixtable->header && aPixtable->table,
264  CPL_ERROR_NULL_INPUT);
265  cpl_ensure_code(muse_cpltable_check(aPixtable->table, muse_pixtable_def)
266  == CPL_ERROR_NONE, CPL_ERROR_DATA_NOT_FOUND);
267  cpl_ensure_code(aFlags != 0, CPL_ERROR_UNSUPPORTED_MODE);
268 
269  /* check to see if the incoming table is in air wavelengths; *
270  * start by assuming that it is */
271  if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_SPEC_TYPE)) {
272  const char *spectype = cpl_propertylist_get_string(aPixtable->header,
274  if (spectype && !strncmp(spectype, "WAVE", 4)) {
275  cpl_msg_warning(__func__, "Pixel table has spectral type \"%s\", not in "
276  "air wavelengths!", spectype);
277  return cpl_error_set(__func__, CPL_ERROR_TYPE_MISMATCH);
278  }
279  if (spectype && strncmp(spectype, "AWAV", 4)) {
280  cpl_msg_warning(__func__, "Pixel table has unknown spectral type \"%s\", "
281  "not in air wavelengths, not doing conversion to vacuum!",
282  spectype);
283  return cpl_error_set(__func__, CPL_ERROR_TYPE_MISMATCH);
284  }
285  } /* if incoming spec type */
286 
287  /* set up International Standard Atmosphere by default */
288  double rhum = 0., /* dry aig */
289  temp = 15. + 273.15, /* [K] */
290  pres = 1013.25; /* [hPa] */
291  /* check if the 2nd bit is set, if yes query measured air parameters */
292  cpl_boolean realair = (aFlags & 0x3) >> 1;
293  if (realair) {
294  /* query measured values from pixel table header */
295  cpl_errorstate prestate = cpl_errorstate_get();
296  temp = muse_pfits_get_temp(aPixtable->header); /* temperature [Celsius] */
297  if (!cpl_errorstate_is_equal(prestate)) {
298  cpl_msg_warning(__func__, "Pixel table header does not contain temperature"
299  ", no conversion to vacuum: %s", cpl_error_get_message());
300  return CPL_ERROR_DATA_NOT_FOUND;
301  }
302  temp += 273.15; /* T[K] */
303 
304  prestate = cpl_errorstate_get();
305  rhum = muse_pfits_get_rhum(aPixtable->header); /* relative humidity [%] */
306  if (!cpl_errorstate_is_equal(prestate)) {
307  cpl_msg_warning(__func__, "Pixel table header does not contain relative "
308  "humidity, no conversion to vacuum: %s",
309  cpl_error_get_message());
310  return CPL_ERROR_DATA_NOT_FOUND;
311  }
312  rhum /= 100.; /* convert from % to fraction */
313 
314  prestate = cpl_errorstate_get();
315  pres = (muse_pfits_get_pres_start(aPixtable->header) /* pressure [mbar] */
316  + muse_pfits_get_pres_end(aPixtable->header)) / 2.;
317  if (!cpl_errorstate_is_equal(prestate)) {
318  cpl_msg_warning(__func__, "Pixel table header does not contain pressure "
319  "values, no conversion to vacuum: %s",
320  cpl_error_get_message());
321  return CPL_ERROR_DATA_NOT_FOUND;
322  }
323  } /* if MUSE_PHYS_CONVERT_MEASURED */
324 
325  /* check the method to use */
326  unsigned int method = aFlags & 0xC; /* cannot produce invalid method */
327 
328  /* compute refractive index for reference wavelength in um and *
329  * output properties in "natural" (for the formulae) units */
330  double d1, d2,
331  fp, /* water vapor pressure in mmHg */
332  p = pres * 100, /* atmopheric pressure [Pa] */
333  t = temp - 273.15, /* temperature [Celsius] */
334  pv, /* partial vapor pressure */
335  xv; /* mole fraction */
336  if (method == MUSE_PHYS_METHOD_OWENS) {
337  muse_phys_nrindex_owens_coeffs(temp, rhum, pres, &d1, &d2);
338  cpl_msg_info(__func__, "Air to vacuum conversion for T=%.2f K, RH=%.2f %%, "
339  "pres=%.1f mbar (%s, Owens)", temp, rhum*100., pres,
340  realair ? "measured parameters" : "standard air");
341  } else if (method == MUSE_PHYS_METHOD_EDLEN ||
342  method == MUSE_PHYS_METHOD_CIDDOR) {
343  /* follow the detailed procedure from *
344  * http://emtoolbox.nist.gov/Wavelength/Documentation.asp#AppendixA */
345  /* Calculate saturation vapor pressure psv(t) over water, *
346  * equations (A1) to (A7) */
347  const double T = temp, /* temperature [K] */
348  K1 = 1.16705214528E+03,
349  K2 = -7.24213167032E+05,
350  K3 = -1.70738469401E+01,
351  K4 = 1.20208247025E+04,
352  K5 = -3.23255503223E+06,
353  K6 = 1.49151086135E+01,
354  K7 = -4.82326573616E+03,
355  K8 = 4.05113405421E+05,
356  K9 = -2.38555575678E-01,
357  K10 = 6.50175348448E+02,
358  Omega = T + K9 / (T - K10),
359  A = Omega*Omega + K1 * Omega + K2,
360  B = K3 * Omega*Omega + K4 * Omega + K5,
361  C = K6 * Omega*Omega + K7 * Omega + K8,
362  X = -B + sqrt(B*B - 4 * A * C);
363  double psv_w = 1e6 * pow(2 * C / X, 4);
364  /* Calculate saturation vapor pressure psv(t) over ice, *
365  * equations (A8) tl (A13) */
366  const double A1 = -13.928169,
367  A2 = 34.7078238,
368  Theta = T / 273.16,
369  Y = A1 * (1 - pow(Theta, -1.5)) + A2 * (1 - pow(Theta, -1.25));
370  double psv_i = 611.657 * exp(Y);
371  /* determine Humidity */
372  /* for the Edlén equation: for relative humidity RH [%], calculate *
373  * partial pressure using psw_w for t > 0 [Celsius] and psw_i for t < 0 */
374  if (t > 0.) {
375  pv = rhum * psv_w; /* use saturation pressure over water */
376  } else {
377  pv = rhum * psv_i; /* use saturation pressure over ice */
378  }
379  /* for the Ciddor equation: express humidity as a mole fraction */
380  const double alpha = 1.00062,
381  beta = 3.14e-8,
382  gamma = 5.60e-7;
383  double f = alpha + beta * p + gamma * t*t; /* enhancement factor f(p,t) */
384  /* for relative humidity RH [%], calculate mole fraction xv using psw(t) */
385  if (t > 0.) {
386  xv = rhum * f * psv_w / p; /* use saturation pressure over water */
387  } else {
388  xv = rhum * f * psv_i / p; /* use saturation pressure over ice */
389  }
390  cpl_msg_info(__func__, "Air to vacuum conversion for T=%.2f degC, RH=%.2f "
391  "%%, p=%.1f Pa (%s, %s)", t, rhum*100., p,
392  realair ? "measured parameters" : "standard air",
393  method == MUSE_PHYS_METHOD_EDLEN ? "Edlen" : "Ciddor");
394  } else {
395  /* use the Owens formula to derive saturation pressure, still needs T[K] */
397  /* using that, derive the water vapor pressure in mmHg */
398  fp = rhum * ps * MUSE_PHYS_hPa_TO_mmHg;
399  temp -= 273.15; /* temperature (again) in degrees Celsius */
400  pres *= MUSE_PHYS_hPa_TO_mmHg; /* need the pressure in mmHg */
401  cpl_msg_info(__func__, "Air to vacuum conversion for T=%.2f degC, fp=%.3f "
402  "mmHg, P=%.1f mmHg (%s, Filippenko)", temp, fp, pres,
403  realair ? "measured parameters" : "standard air");
404  } /* else: Filippenko */
405 
406  float *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA);
407  cpl_size i, nmax = muse_pixtable_get_nrow(aPixtable);
408  #pragma omp parallel for default(none) \
409  shared(d1, d2, fp, lbda, method, nmax, p, pres, pv, t, temp, xv)
410  for (i = 0; i < nmax; i++) {
411  double nr, /* the refractive index that we want */
412  lambda = lbda[i] * 1e-4; /* the wavelength in um */
413  if (method == MUSE_PHYS_METHOD_OWENS) {
414  nr = muse_phys_nrindex_owens(lambda, d1, d2);
415  } else if (method == MUSE_PHYS_METHOD_EDLEN) {
416  nr = muse_phys_nrindex_edlen(lambda, t, p, pv);
417  } else if (method == MUSE_PHYS_METHOD_CIDDOR) {
418  /* The standard CO2 content at ground level is 0.0322% (see US 1976 *
419  * atmosphere), but Birch et al. apparently give 450., so use that. */
420  nr = muse_phys_nrindex_ciddor(lambda, temp, p, xv, 450.);
421  } else {
422  nr = muse_phys_nrindex_filippenko(lambda, temp, fp, pres);
423  }
424  /* directly apply the correction to the pixel table */
425  lbda[i] *= nr; /* lambda_vac = nrindex * lambda_air */
426  } /* for i (all pixel table rows) */
427 
428  /* need to recompute the (spectral) pixel table limits now */
429  muse_pixtable_compute_limits(aPixtable);
430 
431  /* add the header with type code as defined in FITS for vacuum */
432  cpl_propertylist_update_string(aPixtable->header, MUSE_HDR_PT_SPEC_TYPE,
433  "WAVE");
434  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_SPEC_TYPE,
435  MUSE_HDR_PT_SPEC_TYPE_COMMENT);
436 
437  return CPL_ERROR_NONE;
438 } /* muse_phys_air_to_vacuum() */
439 
double muse_phys_nrindex_edlen(double l, double t, double p, double pv)
Compute the refractive index for the given wavelength following Edlén.
Definition: muse_phys.c:115
double muse_pfits_get_rhum(const cpl_propertylist *aHeaders)
find out the relavtive humidity (in %)
Definition: muse_pfits.c:1024
double muse_phys_nrindex_filippenko(double l, double T, double f, double P)
Compute the refractive index for the given wavelength following Filippenko.
Definition: muse_phys.c:193
cpl_size muse_pixtable_get_nrow(const muse_pixtable *aPixtable)
get the number of rows within the pixel table
void muse_phys_nrindex_owens_coeffs(double temp, double rhum, double pres, double *d1, double *d2)
Compute the two coefficients for the Owens method.
Definition: muse_phys.c:70
double muse_pfits_get_pres_start(const cpl_propertylist *aHeaders)
find out the ambient pressure at start of exposure (in mbar)
Definition: muse_pfits.c:1042
double muse_phys_nrindex_ciddor(double l, double T, double p, double xv, double xCO2)
Compute the refractive index for the given wavelength following Ciddor.
Definition: muse_phys.c:148
cpl_table * table
The pixel table.
double muse_phys_nrindex_owens(double l, double d1, double d2)
Compute the refractive index for the given wavelength following Owens.
Definition: muse_phys.c:95
#define MUSE_HDR_PT_SPEC_TYPE
cpl_error_code muse_cpltable_check(const cpl_table *aTable, const muse_cpltable_def *aDef)
Check whether the table contains the fields of the definition.
Structure definition of MUSE pixel table.
double muse_pfits_get_pres_end(const cpl_propertylist *aHeaders)
find out the ambient pressure at end of exposure (in mbar)
Definition: muse_pfits.c:1060
cpl_error_code muse_phys_air_to_vacuum(muse_pixtable *aPixtable, unsigned int aFlags)
Convert a pixel table from air to vacuum wavelengths.
Definition: muse_phys.c:261
const muse_cpltable_def muse_pixtable_def[]
MUSE pixel table definition.
#define MUSE_PIXTABLE_LAMBDA
Definition: muse_pixtable.h:53
double muse_pfits_get_temp(const cpl_propertylist *aHeaders)
find out the ambient temperature (in degrees Celsius)
Definition: muse_pfits.c:1006
double muse_phys_nrindex_owens_saturation_pressure(double temp)
Compute the saturation pressure using the Owens calibration.
Definition: muse_phys.c:54
cpl_error_code muse_pixtable_compute_limits(muse_pixtable *aPixtable)
(Re-)Compute the limits of the coordinate columns of a pixel table.
cpl_propertylist * header
The FITS header.