MUSE Pipeline Reference Manual  1.0.2
muse_dar.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) 2001 Enrico Marchetti, European Southern Observatory
7  *
8  * This program is free software; you can redistribute it and/or modify
9  * it under the terms of the GNU General Public License as published by
10  * the Free Software Foundation; either version 2 of the License, or
11  * (at your option) any later version.
12  *
13  * This program is distributed in the hope that it will be useful,
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  * GNU General Public License for more details.
17  *
18  * You should have received a copy of the GNU General Public License
19  * along with this program; if not, write to the Free Software
20  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
21  */
22 
23 #ifdef HAVE_CONFIG_H
24 #include <config.h>
25 #endif
26 
27 /*----------------------------------------------------------------------------*
28  * Includes *
29  *----------------------------------------------------------------------------*/
30 #include <cpl.h>
31 #include <math.h>
32 #include <string.h>
33 
34 #include "muse_dar.h"
35 #include "muse_instrument.h"
36 
37 #include "muse_astro.h"
38 #include "muse_combine.h"
39 #include "muse_pfits.h"
40 #include "muse_quality.h"
41 #include "muse_utils.h"
42 #include "muse_wcs.h"
43 
44 /*----------------------------------------------------------------------------*
45  * Debugging/feature Macros *
46  *----------------------------------------------------------------------------*/
47 #define DEBUG_MUSE_DARCORRECT 0 /* if 1, output debugging in muse_dar_correct() */
48 #define DEBUG_MUSE_DARCHECK 0 /* if 1, output debug stuff in muse_dar_check(), *
49  * if > 1 also write files */
50 #define MUSE_DARCHECK_GRID_FITS 0 /* if 1, use grid-based fitting in muse_dar_check(); *
51  * it is more robust and slightly more accurate, but *
52  * about 10 times slower! */
53 
54 /*----------------------------------------------------------------------------*/
58 /*----------------------------------------------------------------------------*/
59 
62 /*----------------------------------------------------------------------------*/
71 /*----------------------------------------------------------------------------*/
72 static double
73 muse_dar_owens_saturation_pressure(double temp)
74 {
75  return -10474.0 + 116.43*temp - 0.43284*temp*temp + 0.00053840*pow(temp,3);
76 }
77 
78 /*----------------------------------------------------------------------------*/
88 /*----------------------------------------------------------------------------*/
89 static void
90 muse_dar_owens_coeffs(double temp, double rhum, double pres,
91  double *d1, double *d2)
92 {
93  /* pressure helpers: saturation pressure */
94  double ps = muse_dar_owens_saturation_pressure(temp),
95  p2 = rhum * ps, /* water vapor pressure */
96  p1 = pres - p2; /* dry air pressure */
97  /* helpers/coefficients */
98  *d1 = p1/temp * (1. + p1*(57.90e-8 - 9.3250e-4/temp + 0.25844/temp/temp));
99  *d2 = p2/temp * (1. + p2*(1. + 3.7e-4*p2)
100  * (-2.37321e-3 + 2.23366/temp - 710.792/temp/temp
101  + 7.75141e4/pow(temp,3)));
102 }
103 
104 /*----------------------------------------------------------------------------*/
114 /*----------------------------------------------------------------------------*/
115 static inline double
116 muse_dar_owens_nrindex(double l, double d1, double d2)
117 {
118  double lisq = 1./(l*l); /* inverse square of the wavelength */
119  return ((2371.34 + 683939.7/(130. - lisq) + 4547.3/(38.9 - lisq)) * d1
120  + (6487.31 + 58.058 * lisq - 0.71150 * lisq*lisq
121  + 0.08851 * lisq*lisq*lisq) * d2) * 1.0e-8;
122 }
123 
124 /*----------------------------------------------------------------------------*/
135 /*----------------------------------------------------------------------------*/
136 static inline double
137 muse_dar_filippenko_nrindex(double l, double T, double f, double P)
138 {
139  double lisq = 1./(l*l); /* inverse square of the wavelength */
140  /* 10^6 [n(lambda) - 1] at standard environmental conditions, Eq. (1) */
141  double nl = 64.328 + 29498.1 / (146 - lisq) + 255.4 / (41 - lisq);
142  /* correction for non-standard conditions, Eq. (2) */
143  nl *= P / 720.883 * (1 + (1.049 - 0.0157 * T) * 1e-6 * P) / (1 + 0.003661 * T);
144  /* correction for water vapor, Eq. (3) */
145  nl -= (0.0624 - 0.000680 * lisq) / (1 + 0.003661 * T) * f;
146  /* convert to refractive index n(lambda) */
147  nl = nl * 1e-6 + 1;
148  return nl;
149 }
150 
151 /*----------------------------------------------------------------------------*/
162 /*----------------------------------------------------------------------------*/
163 static inline double
164 muse_dar_edlen_nrindex(double l, double t, double p, double pv)
165 {
166  /* define constants */
167  const double A = 8342.54, B = 2406147,
168  C = 15998, D = 96095.43,
169  E = 0.601, F = 0.00972, G = 0.003661;
170  double S = 1./(l*l), /* supposed to be "laser vacuum wavelength"?! XXX */
171  /* intermediate result for p, pv, and t */
172  n_5 = 1 + 1e-8 * (A + B / (130 - S) + C / (38.9 - S)),
173  X = (1 + 1e-8 * (E - F * t) * p) / (1 + G * t),
174  n_tp = 1 + p * (n_5 - 1) * X / D,
175  n = n_tp - 10e-10 * (292.75 / (t + 273.15)) * (3.7345 - 0.0401 * S) * pv;
176  return n;
177 } /* muse_dar_edlen_nrindex() */
178 
179 /*----------------------------------------------------------------------------*/
195 /*----------------------------------------------------------------------------*/
196 static inline double
197 muse_dar_ciddor_nrindex(double l, double T, double p, double xv, double xCO2)
198 {
199  /* define constants */
200  const double w0 = 295.235, w1 = 2.6422, w2 = -0.03238, w3 = 0.004028,
201  k0 = 238.0185, k1 = 5792105, k2 = 57.362, k3 = 167917,
202  a0 = 1.58123e-6, a1 = -2.9331e-8, a2 = 1.1043e-10,
203  b0 = 5.707e-6, b1 = -2.051e-8,
204  c0 = 1.9898e-4, c1 = -2.376e-6,
205  d = 1.83e-11, e = -0.765e-8,
206  p_R1 = 101325, T_R1 = 288.15,
207  Z_a = 0.9995922115,
208  rho_vs = 0.00985938,
209  R = 8.314472, M_v = 0.018015;
210  double t = T - 273.15, /* temperature [Celsius] */
211  S = 1./(l*l), /* supposed to be "laser vacuum wavelength"?! XXX */
212  /* intermediate results that depend on S */
213  r_as = 1e-8 * ((k1 / (k0 - S)) + (k3 / (k2 - S))),
214  r_vs = 1.022e-8 * (w0 + w1 * S + w2 * S*S + w3 * S*S*S),
215  /* using the CO2 concentration, calculate intermediate values for CO2 */
216  M_a = 0.0289635 + 1.2011e-8 * (xCO2 - 400),
217  r_axs = r_as * (1 + 5.34e-7 * (xCO2 - 450)),
218  /* compressibility and density components */
219  Z_m = 1 - (p / T) * (a0 + a1 * t + a2 * t*t + (b0 + b1 * t) * xv
220  + (c0 + c1 * t) * xv*xv)
221  + (p/T)*(p/T) * (d + e * xv*xv),
222  rho_axs = p_R1 * M_a / (Z_a * R * T_R1),
223  rho_v = xv * p * M_v / (Z_m * R * T),
224  rho_a = (1 - xv) * p * M_a / (Z_m * R * T),
225  n = 1 + (rho_a / rho_axs) * r_axs + (rho_v / rho_vs) * r_vs;
226  return n;
227 } /* muse_dar_ciddor_nrindex() */
228 
229 /*----------------------------------------------------------------------------*/
304 /*----------------------------------------------------------------------------*/
305 cpl_error_code
306 muse_dar_correct(muse_pixtable *aPixtable, double aLambdaRef)
307 {
308  cpl_ensure_code(aPixtable && aPixtable->header, CPL_ERROR_NULL_INPUT);
309  if (aLambdaRef > 22000. || aLambdaRef < 3500.) {
310  cpl_msg_info(__func__, "Reference lambda %.3f Angstrom: skipping DAR "
311  "correction", aLambdaRef);
312  cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_NAME, -1.);
313  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_NAME,
314  MUSE_HDR_PT_DAR_C_SKIP);
315  return CPL_ERROR_NONE;
316  }
317  if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_NAME) &&
318  cpl_propertylist_get_double(aPixtable->header, MUSE_HDR_PT_DAR_NAME) > 0.) {
319  cpl_msg_info(__func__, "pixel table already corrected: skipping DAR "
320  "correction");
321  return CPL_ERROR_NONE;
322  }
323 
324  /* -------------------------------------------------------------- *
325  * get environmental and instrumental values from the FITS header *
326  * -------------------------------------------------------------- */
327  double airmass = muse_astro_airmass(aPixtable->header);
328  cpl_ensure_code(airmass >= 1., cpl_error_get_code());
329  /* simple zenith distance in radians */
330  double z = acos(1./airmass);
331 
332  /* compute the mean parallactic angle during exposure */
333  cpl_errorstate prestate = cpl_errorstate_get();
334  double parang = muse_astro_parangle(aPixtable->header);
335  cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
336 
337  /* compute the position angle on the sky from the angles we have */
338  prestate = cpl_errorstate_get();
339  double posang = muse_astro_posangle(aPixtable->header);
340  cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
341 
342  cpl_boolean isWFM = muse_pfits_get_mode(aPixtable->header)
343  <= MUSE_MODE_WFM_AO_N;
344  if (!isWFM) {
345  cpl_msg_warning(__func__, "For NFM instrument mode DAR correction should "
346  "not be needed!");
347  }
348 
349  /* query values and set defaults in case we didn't get proper values */
350  prestate = cpl_errorstate_get();
351  double temp = muse_pfits_get_temp(aPixtable->header); /* temperature [Celsius] */
352  if (cpl_errorstate_is_equal(prestate)) {
353  temp += 273.15; /* T[K] */
354  } else {
355  cpl_msg_warning(__func__, "FITS header does not contain temperature: %s",
356  cpl_error_get_message());
357  temp = 11.5 + 273.15; /* T[K] */
358  }
359  prestate = cpl_errorstate_get();
360  double rhum = muse_pfits_get_rhum(aPixtable->header); /* relative humidity [%] */
361  if (!cpl_errorstate_is_equal(prestate)) {
362  cpl_msg_warning(__func__, "FITS header does not contain relative humidity: %s",
363  cpl_error_get_message());
364  rhum = 14.5;
365  }
366  rhum /= 100.; /* convert from % to fraction */
367  prestate = cpl_errorstate_get();
368  double pres = (muse_pfits_get_pres_start(aPixtable->header) /* pressure [mbar] */
369  + muse_pfits_get_pres_end(aPixtable->header)) / 2.;
370  if (!cpl_errorstate_is_equal(prestate)) {
371  cpl_msg_warning(__func__, "FITS header does not contain pressure values: %s",
372  cpl_error_get_message());
373  pres = 743.;
374  }
375 
376  /* -------------------------------------------------------------------------- *
377  * choose method and compute the refractive index at the reference wavelength *
378  * -------------------------------------------------------------------------- */
379  enum {
380  MUSE_DAR_METHOD_FILIPPENKO = 0,
381  MUSE_DAR_METHOD_OWENS,
382  MUSE_DAR_METHOD_EDLEN,
383  MUSE_DAR_METHOD_CIDDOR,
384  };
385  int method = MUSE_DAR_METHOD_FILIPPENKO;
386  if (getenv("MUSE_DAR_CORRECT_METHOD") &&
387  !strncmp(getenv("MUSE_DAR_CORRECT_METHOD"), "Owens", 6)) {
388  method = MUSE_DAR_METHOD_OWENS;
389  } else if (getenv("MUSE_DAR_CORRECT_METHOD") &&
390  !strncmp(getenv("MUSE_DAR_CORRECT_METHOD"), "Edlen", 6)) {
391  method = MUSE_DAR_METHOD_EDLEN;
392  } else if (getenv("MUSE_DAR_CORRECT_METHOD") &&
393  !strncmp(getenv("MUSE_DAR_CORRECT_METHOD"), "Ciddor", 7)) {
394  method = MUSE_DAR_METHOD_CIDDOR;
395  }
396  /* compute refractive index for reference wavelength in um and *
397  * output properties in "natural" (for the formulae) units */
398  double d1, d2,
399  fp, /* water vapor pressure in mmHg */
400  nr0, /* refractive index of air at reference wavelength */
401  p = pres * 100, /* atmopheric pressure [Pa] */
402  t = temp - 273.15, /* temperature [Celsius] */
403  pv, /* partial vapor pressure */
404  xv; /* mole fraction */
405  if (method == MUSE_DAR_METHOD_OWENS) {
406  muse_dar_owens_coeffs(temp, rhum, pres, &d1, &d2);
407  nr0 = muse_dar_owens_nrindex(aLambdaRef / 10000., d1, d2);
408  cpl_msg_info(__func__, "DAR for T=%.2f K, RH=%.2f %%, pres=%.1f mbar, "
409  "at airmass=%.4f, using ref. wavelength %.6f um (Owens)",
410  temp, rhum*100., pres, 1./cos(z), aLambdaRef / 10000.);
411  } else if (method == MUSE_DAR_METHOD_EDLEN ||
412  method == MUSE_DAR_METHOD_CIDDOR) {
413  /* follow the detailed procedure from *
414  * http://emtoolbox.nist.gov/Wavelength/Documentation.asp#AppendixA */
415  /* Calculate saturation vapor pressure psv(t) over water,
416  * equations (A1) to (A7) */
417  const double T = temp, /* temperature [K] */
418  K1 = 1.16705214528E+03,
419  K2 = -7.24213167032E+05,
420  K3 = -1.70738469401E+01,
421  K4 = 1.20208247025E+04,
422  K5 = -3.23255503223E+06,
423  K6 = 1.49151086135E+01,
424  K7 = -4.82326573616E+03,
425  K8 = 4.05113405421E+05,
426  K9 = -2.38555575678E-01,
427  K10 = 6.50175348448E+02,
428  Omega = T + K9 / (T - K10),
429  A = Omega*Omega + K1 * Omega + K2,
430  B = K3 * Omega*Omega + K4 * Omega + K5,
431  C = K6 * Omega*Omega + K7 * Omega + K8,
432  X = -B + sqrt(B*B - 4 * A * C);
433  double psv_w = 1e6 * pow(2 * C / X, 4);
434  /* Calculate saturation vapor pressure psv(t) over ice, *
435  * equations (A8) tl (A13) */
436  const double A1 = -13.928169,
437  A2 = 34.7078238,
438  Theta = T / 273.16,
439  Y = A1 * (1 - pow(Theta, -1.5)) + A2 * (1 - pow(Theta, -1.25));
440  double psv_i = 611.657 * exp(Y);
441  /* determine Humidity */
442  /* for the EdlĂ©n equation: for relative humidity RH [%], calculate *
443  * partial pressure using psw_w for t > 0 [Celsius] and psw_i for t < 0 */
444  if (t > 0.) {
445  pv = rhum * psv_w; /* use saturation pressure over water */
446  } else {
447  pv = rhum * psv_i; /* use saturation pressure over ice */
448  }
449  /* for the Ciddor equation: express humidity as a mole fraction */
450  const double alpha = 1.00062,
451  beta = 3.14e-8,
452  gamma = 5.60e-7;
453  double f = alpha + beta * p + gamma * t*t; /* enhancement factor f(p,t) */
454  /* for relative humidity RH [%], calculate mole fraction xv using psw(t) */
455  if (t > 0.) {
456  xv = rhum * f * psv_w / p; /* use saturation pressure over water */
457  } else {
458  xv = rhum * f * psv_i / p; /* use saturation pressure over ice */
459  }
460  if (method == MUSE_DAR_METHOD_EDLEN) {
461  nr0 = muse_dar_edlen_nrindex(aLambdaRef / 10000., t, p, pv);
462  } else { /* MUSE_DAR_METHOD_CIDDOR */
463  nr0 = muse_dar_ciddor_nrindex(aLambdaRef / 10000., T, p, xv, 450.);
464  }
465  cpl_msg_info(__func__, "DAR for T=%.2f degC, RH=%.2f %%, p=%.1f Pa, "
466  "at airmass=%.4f, using ref. wavelength %.6f um (%s)",
467  t, rhum*100., p, 1./cos(z), aLambdaRef / 10000.,
468  method == MUSE_DAR_METHOD_EDLEN ? "Edlen" : "Ciddor");
469  } else {
470 #define CONVERT_hPa_TO_mmHg 0.75006158 /* conversion from hPa (or mbar) to *
471  * mmHg, needed for Filippenko */
472  /* use the Owens formula to derive saturation pressure, still needs T[K] */
473  double ps = muse_dar_owens_saturation_pressure(temp);
474  /* using that, derive the water vapor pressure in mmHg */
475  fp = rhum * ps * CONVERT_hPa_TO_mmHg;
476  temp -= 273.15; /* temperature (again) in degrees Celsius */
477  pres *= CONVERT_hPa_TO_mmHg; /* need the pressure in mmHg */
478  nr0 = muse_dar_filippenko_nrindex(aLambdaRef / 10000., temp, fp, pres);
479  cpl_msg_info(__func__, "DAR for T=%.2f degC, fp=%.3f mmHg, P=%.1f mmHg, "
480  "at airmass=%.4f, using ref. wavelength %.6f um (Filippenko)",
481  temp, fp, pres, 1./cos(z), aLambdaRef / 10000.);
482  }
483 
484  /* ---------------------------------- *
485  * compute the direction of the shift *
486  * ---------------------------------- */
487  /* xshift is in E-W direction for posang = 0, yshift is N-S */
488  double xshift = -sin((parang + posang) * CPL_MATH_RAD_DEG),
489  yshift = cos((parang + posang) * CPL_MATH_RAD_DEG);
490 
491  /* apply correction for spaxel size */
492  muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
493  cpl_ensure_code(wcstype == MUSE_PIXTABLE_WCS_PIXEL ||
494  wcstype == MUSE_PIXTABLE_WCS_CELSPH, CPL_ERROR_INVALID_TYPE);
495  double fscale = 3600.; /* assume scale in arcsec */
496  if (wcstype == MUSE_PIXTABLE_WCS_CELSPH) {
497  fscale = 1. / 3600; /* scale in degrees */
498  double xscale, yscale;
499  muse_wcs_get_scales(aPixtable->header, &xscale, &yscale);
500  xshift /= -xscale; /* need to reverse the sign for RA */
501  yshift /= yscale;
502  } else { /* assume nominal spatial scale */
503  if (isWFM) {
504  xshift /= kMuseSpaxelSizeX_WFM;
505  yshift /= kMuseSpaxelSizeY_WFM;
506  } else {
507  xshift /= kMuseSpaxelSizeX_NFM;
508  yshift /= kMuseSpaxelSizeY_NFM;
509  }
510  }
511 
512  /* ------------------------------ *
513  * apply correction to all pixels *
514  * ------------------------------ */
515  /* diff. refr. base in arcsec converted from radians (Filippenko does *
516  * the conversion using x206265 which converts radians to arcsec) */
517  double dr0 = tan(z) * fscale * CPL_MATH_DEG_RAD;
518 #if DEBUG_MUSE_DARCORRECT
519  double wl;
520  for (wl = 4650.; wl <= 9300; wl += 50) {
521  double nr;
522  if (method == MUSE_DAR_METHOD_OWENS) {
523  nr = muse_dar_owens_nrindex(wl / 10000., d1, d2);
524  } else if (method == MUSE_DAR_METHOD_EDLEN) {
525  nr = muse_dar_edlen_nrindex(wl / 10000., t, p, pv);
526  } else if (method == MUSE_DAR_METHOD_CIDDOR) {
527  nr = muse_dar_ciddor_nrindex(wl / 10000., temp, p, xv, 450.);
528  } else {
529  nr = muse_dar_filippenko_nrindex(wl / 10000., temp, fp, pres);
530  }
531  double dr = dr0 * (nr0 - nr);
532  printf("%.1f Angstrom: %f ==> %f %f %s (%s)\n", wl, dr,
533  xshift * dr, yshift * dr,
534  wcstype == MUSE_PIXTABLE_WCS_CELSPH ? "deg" : "pix", __func__);
535  }
536  fflush(stdout);
537 #endif
538  float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
539  *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
540  *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA);
541  cpl_size i, nmax = muse_pixtable_get_nrow(aPixtable);
542  #pragma omp parallel for default(none) \
543  shared(d1, d2, dr0, fp, lbda, method, nmax, nr0, p, pres, pv, t, \
544  temp, xpos, ypos, xshift, xv, yshift)
545  for (i = 0; i < nmax; i++) {
546  double nr, lambda = lbda[i] * 1e-4;
547  if (method == MUSE_DAR_METHOD_OWENS) {
548  nr = muse_dar_owens_nrindex(lambda, d1, d2);
549  } else if (method == MUSE_DAR_METHOD_EDLEN) {
550  nr = muse_dar_edlen_nrindex(lambda, t, p, pv);
551  } else if (method == MUSE_DAR_METHOD_CIDDOR) {
552  nr = muse_dar_ciddor_nrindex(lambda, temp, p, xv, 450.);
553  } else {
554  nr = muse_dar_filippenko_nrindex(lambda, temp, fp, pres);
555  }
556  double dr = dr0 * (nr0 - nr);
557  /* correct coordinates directly in the pixel table */
558  xpos[i] += xshift * dr;
559  ypos[i] += yshift * dr;
560  } /* for i (all pixel table rows) */
561  /* add the header */
562  cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_NAME,
563  aLambdaRef);
564  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_NAME,
565  MUSE_HDR_PT_DAR_COMMENT);
566  /* need to recompute the (spatial) pixel table limits, but save the *
567  * previous contents because that's what we need for WCS calibration */
568  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XLO,
569  cpl_propertylist_get_float(aPixtable->header,
570  MUSE_HDR_PT_XLO));
571  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XHI,
572  cpl_propertylist_get_float(aPixtable->header,
573  MUSE_HDR_PT_XHI));
574  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YLO,
575  cpl_propertylist_get_float(aPixtable->header,
576  MUSE_HDR_PT_YLO));
577  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YHI,
578  cpl_propertylist_get_float(aPixtable->header,
579  MUSE_HDR_PT_YHI));
580  muse_pixtable_compute_limits(aPixtable);
581 
582  return CPL_ERROR_NONE;
583 } /* muse_dar_correct() */
584 
585 /*----------------------------------------------------------------------------*/
639 /*----------------------------------------------------------------------------*/
640 cpl_error_code
641 muse_dar_check(muse_pixtable *aPixtable, double *aMaxShift,
642  cpl_boolean aCorrect, muse_datacube **aCube)
643 {
644  cpl_ensure_code(aMaxShift, CPL_ERROR_NULL_INPUT);
645  *aMaxShift = -99.;
646  cpl_ensure_code(aPixtable, CPL_ERROR_NULL_INPUT);
647  muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
648  cpl_ensure_code(wcstype == MUSE_PIXTABLE_WCS_PIXEL, CPL_ERROR_INVALID_TYPE);
649  if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_CHECK)) {
650  cpl_msg_info(__func__, "pixel table already checked for DAR accuracy");
651  }
652  if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_CORR)) {
653  cpl_msg_info(__func__, "pixel table already corrected for DAR residuals");
654  }
655  /* set the halfsize of the window to measure the centroids with respect to *
656  * the reference wavelength; +/-5 pix should suffice in case DAR was *
657  * previously corrected; in other cases, increase it to +/-15 pix */
658  int cenhalfsize = 5;
659  /* polynomial order for the fits, default 2, if not corrected yet, use 4 */
660  int order = 2;
661 
662  /* create cube, with coarse sampling in wavelength */
663  muse_resampling_params *params =
665  params->pfx = 1.; /* large pixfrac to not create sudden spatial shifts */
666  params->pfy = 1.;
667  params->pfl = 1.;
668  params->dlambda = 10.;
669  params->crsigma = 7.; /* moderate CR rejection */
670  muse_pixgrid *grid;
671  muse_datacube *cube = muse_resampling_cube(aPixtable, params, &grid);
673  if (!cube) {
674  muse_pixgrid_delete(grid);
675  muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
676  return cpl_error_set_message(__func__, cpl_error_get_code(),
677  "Failure while creating cube!");
678  }
679  if (aCube) {
680  *aCube = cube;
681  }
682  double crpix3 = cpl_propertylist_get_double(cube->header, "CRPIX3"),
683  crval3 = cpl_propertylist_get_double(cube->header, "CRVAL3"),
684  cdelt3 = cpl_propertylist_get_double(cube->header, "CD3_3");
685  cpl_errorstate state = cpl_errorstate_get(); /* header may be missing */
686  double lambdaref = cpl_propertylist_get_double(aPixtable->header,
688  lambda1 = (2 - crpix3) * cdelt3 + crval3, /* 2nd plane */
689  lambda3 = (grid->size_z - 1 - crpix3) * cdelt3 + crval3, /* next to last */
690  lambda2 = (lambda1 + lambda3) / 2.; /* central plane */
691  if (!cpl_errorstate_is_equal(state) || lambdaref < 0) {
692  cpl_errorstate_set(state);
693  /* use reference wavelength in the center of the wavelength range */
694  lambdaref = lambda2;
695  cpl_msg_warning(__func__, "Data is not DAR corrected, using reference "
696  "wavelength at %.3f Angstrom!", lambdaref);
697  if (!cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_CORR)) {
698  cenhalfsize = 15;
699  order = 4;
700  }
701  }
702  if (lambdaref > lambda3 || lambdaref < lambda1) { /* ref. outside cube! */
703  lambdaref = lambda2;
704  cpl_msg_warning(__func__, "Reference lambda outside cube, using reference "
705  "wavelength at %.3f Angstrom!", lambdaref);
706  }
707 
708  /* detect objects in the cube, using image planes near reference wavelength */
709  int iplane, irefplane = lround((lambdaref - crval3) / cdelt3 + crpix3) - 1;
710  /* medianing three images removes all cosmic rays but continuum objects stay *
711  * (need to use muse_images and the muse_combine function because *
712  * cpl_imagelist_collapse_median_create() disregards all bad pixels) */
714  unsigned int ilist = 0;
715  for (iplane = irefplane - 1; iplane <= irefplane + 1; iplane++) {
716  muse_image *image = muse_image_new();
717  image->data = cpl_image_duplicate(cpl_imagelist_get(cube->data, iplane));
718  image->dq = cpl_image_duplicate(cpl_imagelist_get(cube->dq, iplane));
719  image->stat = cpl_image_duplicate(cpl_imagelist_get(cube->stat, iplane));
720  muse_imagelist_set(list, image, ilist++);
721  } /* for iplane (planes around ref. wavelength) */
722  muse_image *median = muse_combine_median_create(list);
723  muse_imagelist_delete(list);
724  if (!median) {
725  if (!aCube) {
726  muse_datacube_delete(cube);
727  }
728  muse_pixgrid_delete(grid);
729  muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
730  return cpl_error_set_message(__func__, cpl_error_get_code(),
731  "Failure while creating detection image!");
732  }
733 #if DEBUG_MUSE_DARCHECK > 1
734  median->header = cpl_propertylist_new(); /* suppress errors while saving */
735  cpl_propertylist_append_string(median->header, "BUNIT",
736  cpl_propertylist_get_string(cube->header,
737  "BUNIT"));
738  muse_image_save(median, "IMAGE_darcheck_median.fits");
739 #endif
740  /* reject data in image to signify bad pixels to CPL routines */
742  /* use high sigmas for detection */
743  double dsigmas[] = { 50., 30., 10., 8., 6., 5. };
744  int ndsigmas = sizeof(dsigmas) / sizeof(double);
745  cpl_vector *vsigmas = cpl_vector_wrap(ndsigmas, dsigmas);
746  cpl_size isigma = -1;
747  state = cpl_errorstate_get();
748  cpl_apertures *apertures = cpl_apertures_extract(median->data, vsigmas, &isigma);
749  int nx = cpl_image_get_size_x(median->data),
750  ny = cpl_image_get_size_y(median->data);
751  muse_image_delete(median);
752  int napertures = apertures ? cpl_apertures_get_size(apertures) : 0;
753  if (napertures < 1) {
754  cpl_vector_unwrap(vsigmas);
755  cpl_apertures_delete(apertures);
756  cpl_errorstate_set(state); /* reset state to before aperture extraction */
757  if (!aCube) {
758  muse_datacube_delete(cube);
759  }
760  muse_pixgrid_delete(grid);
761  muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
762  return cpl_error_set_message(__func__, CPL_ERROR_DATA_NOT_FOUND, "No sources "
763  "found for DAR check down to %.1f sigma "
764  "limit on plane %d", dsigmas[ndsigmas - 1],
765  irefplane + 1);
766  }
767  cpl_msg_debug(__func__, "The %.1f sigma threshold was used to find %d source%s,"
768  " centered on plane %d", cpl_vector_get(vsigmas, isigma),
769  napertures, napertures == 1 ? "" : "s", iplane+1);
770  cpl_vector_unwrap(vsigmas);
771 #if DEBUG_MUSE_DARCHECK
772  cpl_apertures_dump(apertures, stdout);
773  fflush(stdout);
774 #if DEBUG_MUSE_DARCHECK > 1
775  muse_datacube_save(cube, "DATACUBE_darcheck.fits");
776 #endif
777 #endif
778 
779  /* restrict box sizes to be at most the size of the image plane */
780  int xhalfsize = 2*cenhalfsize > nx ? nx/2 : cenhalfsize,
781  yhalfsize = 2*cenhalfsize > ny ? ny/2 : cenhalfsize;
782 
783  /* loop through the wavelength planes and compute centroids for all apertures */
784  int l, nlambda = cpl_imagelist_get_size(cube->data);
785  cpl_vector *vlambda = cpl_vector_new(nlambda);
786  cpl_matrix *moffsets = cpl_matrix_new(2, nlambda);
787 #if !DEBUG_MUSE_DARCHECK
788  #pragma omp parallel for default(none) \
789  shared(apertures, aPixtable, cdelt3, crpix3, crval3, cube, grid, \
790  moffsets, napertures, nlambda, nx, ny, vlambda, xhalfsize, \
791  yhalfsize)
792 #endif
793  for (l = 0; l < nlambda; l++) {
794  double xoffset = 0., yoffset = 0.,
795  lambda = (l+1 - crpix3) * cdelt3 + crval3;
796  int n, noffsets = 0;
797  for (n = 1; n <= napertures; n++) {
798  /* use detection center to construct (larger) window for centroid */
799  double xc = cpl_apertures_get_centroid_x(apertures, n),
800  yc = cpl_apertures_get_centroid_y(apertures, n);
801  int x1 = lround(xc) - xhalfsize,
802  x2 = lround(xc) + xhalfsize,
803  y1 = lround(yc) - yhalfsize,
804  y2 = lround(yc) + yhalfsize;
805  /* force window to be inside the image */
806  if (x1 < 1) x1 = 1;
807  if (y1 < 1) y1 = 1;
808  if (x2 > nx) x2 = nx;
809  if (y2 > ny) y2 = ny;
810 
811 #define MUSE_DARCHECK_MAX_INTERNAL_SHIFT 5 /* accept max. 5 pix shift between *
812  * image centroid and fit center */
813 #if MUSE_DARCHECK_GRID_FITS /* carry out centroid and fit on the un-resampled *
814  * data, i.e on the pixel grid for better S/N */
815 #if DEBUG_MUSE_DARCHECK
816  const char *method = "grid/moffat";
817 #endif
818  /* transfer data from the pixel table into matrices */
819  float *cdata = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
820  *cstat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
821  int *cdq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
822 #define MUSE_DARCHECK_MAX_NPIX 7000 /* start and increase in fixed size bins */
823  cpl_matrix *pos = cpl_matrix_new(MUSE_DARCHECK_MAX_NPIX, 2);
824  cpl_vector *val = cpl_vector_new(MUSE_DARCHECK_MAX_NPIX),
825  *err = cpl_vector_new(MUSE_DARCHECK_MAX_NPIX);
826  int i, npix = 0, nsize = MUSE_DARCHECK_MAX_NPIX;
827  for (i = x1 - 1; i < x2; i++) {
828  int j;
829  for (j = y1 - 1; j < y2; j++) {
830  cpl_size idx = muse_pixgrid_get_index(grid, i, j, l, CPL_TRUE),
831  ipx, npx = muse_pixgrid_get_count(grid, idx);
832  const cpl_size *rows = muse_pixgrid_get_rows(grid, idx);
833  for (ipx = 0; ipx < npx; ipx++) {
834  if (cdq[rows[ipx]] != EURO3D_GOODPIXEL) {
835  continue;
836  }
837  if (npix >= nsize) { /* increase size of the buffers */
838  nsize += MUSE_DARCHECK_MAX_NPIX;
839  cpl_matrix_resize(pos, 0, nsize - cpl_matrix_get_nrow(pos), 0, 0);
840  cpl_vector_set_size(val, nsize);
841  cpl_vector_set_size(err, nsize);
842  }
843  cpl_matrix_set(pos, npix, 0, i+1);
844  cpl_matrix_set(pos, npix, 1, j+1);
845  cpl_vector_set(val, npix, cdata[rows[ipx]]);
846  cpl_vector_set(err, npix, sqrt(cstat[rows[ipx]]));
847  npix++;
848  } /* for ipx */
849  } /* for j */
850  } /* for i */
851  if (npix < 8) { /* need 8 pixels for 8 Moffat parameters */
852  cpl_matrix_delete(pos);
853  cpl_vector_delete(val);
854  cpl_vector_delete(err);
855  continue;
856  }
857  cpl_matrix_set_size(pos, npix, 2);
858  cpl_vector_set_size(val, npix);
859  cpl_vector_set_size(err, npix);
860 
861  /* compute centroid first, as a kind of first guess */
862  double xc1, yc1;
863  muse_utils_get_centroid(pos, val, NULL/*err*/, &xc1, &yc1,
865 
866  /* now use a 2D Moffat fit for a refined position */
867  cpl_array *pmoffat = cpl_array_new(8, CPL_TYPE_DOUBLE);
868  cpl_array_set(pmoffat, 2, xc1);
869  cpl_array_set(pmoffat, 3, yc1);
870  cpl_error_code rc = muse_utils_fit_moffat_2d(pos, val, err, pmoffat, NULL,
871  NULL, NULL, NULL);
872  double xc2 = cpl_array_get(pmoffat, 2, NULL),
873  yc2 = cpl_array_get(pmoffat, 3, NULL);
874  cpl_array_delete(pmoffat);
875  cpl_matrix_delete(pos);
876  cpl_vector_delete(val);
877  cpl_vector_delete(err);
878 
879  /* if the fit worked and both methods give roughly similar *
880  * positions we should have a valid estimate of the center */
881  if (rc == CPL_ERROR_NONE &&
882  fabs(xc1 - xc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT &&
883  fabs(yc1 - yc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT) {
884  xoffset += xc - xc2;
885  yoffset += yc - yc2;
886  noffsets++;
887  }
888 #else /* non-MUSE_DARCHECK_GRID_FITS follows, using wavelength plane images */
889 #if DEBUG_MUSE_DARCHECK
890  const char *method = "image/gauss";
891 #endif
892  /* add CPL bad pixel mask in image plane to get reliable results */
893  cpl_image *plane = cpl_imagelist_get(cube->data, l),
894  *plerr = cpl_image_power_create(cpl_imagelist_get(cube->stat, l), 0.5);
895  /* make sure to exclude bad pixels from the fits below */
896  muse_quality_image_reject_using_dq(plane, cpl_imagelist_get(cube->dq, l), plerr);
897  /* the error image may contain more bad pixels! */
898  cpl_image_reject_from_mask(plane, cpl_image_get_bpm(plerr));
899 
900  cpl_image *xtr = cpl_image_extract(plane, x1, y1, x2, y2);
901  int npix = (x2-x1+1) * (y2-y1+1) /* number of (good) pixels involved */
902  - cpl_image_count_rejected(xtr);
903  cpl_image_delete(xtr);
904  if (npix < 7) { /* need 7 pixels for 7 Gaussian parameters */
905  cpl_image_delete(plerr);
906  continue;
907  }
908 
909  /* compute centroid first, as a kind of first guess */
910  double xc1, yc1;
911  muse_utils_image_get_centroid_window(plane, x1, y1, x2, y2, &xc1, &yc1,
913 
914  /* do Gaussian fit, assuming that the object is roughly stellar */
915  cpl_array *pgauss = cpl_array_new(7, CPL_TYPE_DOUBLE);
916  cpl_array_set(pgauss, 3, xc1);
917  cpl_array_set(pgauss, 4, yc1);
918  cpl_error_code rc = cpl_fit_image_gaussian(plane, plerr,
919  lround(xc), lround(yc),
920  2*xhalfsize, 2*yhalfsize,
921  pgauss, NULL, NULL,
922  NULL, NULL, NULL,
923  NULL, NULL, NULL, NULL);
924  double xc2 = cpl_array_get(pgauss, 3, NULL),
925  yc2 = cpl_array_get(pgauss, 4, NULL);
926  cpl_array_delete(pgauss);
927  cpl_image_delete(plerr);
928 
929  /* if the fit worked and both methods give roughly similar *
930  * positions we should have a valid estimate of the center */
931  if (rc == CPL_ERROR_NONE &&
932  fabs(xc1 - xc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT &&
933  fabs(yc1 - yc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT) {
934  xoffset += xc - xc2;
935  yoffset += yc - yc2;
936  noffsets++;
937  }
938 #endif /* end of MUSE_DARCHECK_GRID_FITS */
939 #if DEBUG_MUSE_DARCHECK /* DEBUG output; compare direct centroid and fit position */
940  printf("%.1f Angstrom plane %d aper %d (%f,%f): %f %f (%s) %f %f (%d, %s)\n",
941  lambda, l+1, n, xc, yc, xc - xc2, yc - yc2, __func__,
942  xc2 - xc1, yc2 - yc1, npix, method);
943  fflush(stdout);
944 #endif
945  } /* for n (all apertures) */
946  /* record average offset in the matrix */
947  cpl_vector_set(vlambda, l, lambda);
948  cpl_matrix_set(moffsets, 0, l, xoffset / noffsets);
949  cpl_matrix_set(moffsets, 1, l, yoffset / noffsets);
950 #if DEBUG_MUSE_DARCHECK
951  printf("%.1f Angstrom plane %d all-apers: %f %f (%s)\n",
952  lambda, l+1, xoffset / noffsets, yoffset / noffsets, __func__);
953  fflush(stdout);
954 #endif
955  } /* for l (all wavelengths) */
956  cpl_apertures_delete(apertures);
957 
958  /* compute correction polynomials and maximum relative offset */
959  cpl_vector *vlam1 = cpl_vector_duplicate(vlambda),
960  *vlam2 = cpl_vector_duplicate(vlambda);
961  cpl_matrix *mlam1 = cpl_matrix_wrap(1, nlambda, cpl_vector_unwrap(vlam1)),
962  *mlam2 = cpl_matrix_wrap(1, nlambda, cpl_vector_unwrap(vlam2));
963  cpl_matrix *mxoff = cpl_matrix_extract_row(moffsets, 0),
964  *myoff = cpl_matrix_extract_row(moffsets, 1);
965  cpl_vector *vxoff = cpl_vector_wrap(nlambda, cpl_matrix_unwrap(mxoff)),
966  *vyoff = cpl_vector_wrap(nlambda, cpl_matrix_unwrap(myoff));
967  /* use iterative polynomial fits, separately in both spatial dimensions, *
968  * to throw out the most extreme outliers, and compute the correction */
969  double msex, msey, chisqx, chisqy;
970  cpl_polynomial *px = muse_utils_iterate_fit_polynomial(mlam1, vxoff,
971  NULL, NULL, order, 3.,
972  &msex, &chisqx),
973  *py = muse_utils_iterate_fit_polynomial(mlam2, vyoff,
974  NULL, NULL, order, 3.,
975  &msey, &chisqy);
976  cpl_vector_delete(vxoff);
977  cpl_vector_delete(vyoff);
978  cpl_matrix_delete(mlam1);
979  cpl_matrix_delete(mlam2);
980 #if DEBUG_MUSE_DARCHECK
981  printf("polynomial fit in x (order %d, rms=%f, chisq=%g):\n", order,
982  sqrt(msex), chisqx);
983  cpl_polynomial_dump(px, stdout);
984  printf("polynomial fit in y (order %d, rms=%f, chisq=%g):\n", order,
985  sqrt(msey), chisqy);
986  cpl_polynomial_dump(py, stdout);
987  fflush(stdout);
988 #endif
989 
990  /* shift to the zeropoint of the offsets at the reference wavelength; *
991  * this should be more stable than just taking the raw input value */
992  double xref = cpl_polynomial_eval_1d(px, lambdaref, NULL),
993  yref = cpl_polynomial_eval_1d(py, lambdaref, NULL);
994  cpl_size pows = 0; /* set the zero-order coefficient */
995  cpl_polynomial_set_coeff(px, &pows, cpl_polynomial_get_coeff(px, &pows) - xref);
996  cpl_polynomial_set_coeff(py, &pows, cpl_polynomial_get_coeff(py, &pows) - yref);
997 #if DEBUG_MUSE_DARCHECK
998  printf("polynomial in x (xref=%f at %f --> %f):\n", xref, lambdaref,
999  cpl_polynomial_eval_1d(px, lambdaref, NULL));
1000  cpl_polynomial_dump(px, stdout);
1001  printf("polynomial in y (yref=%f at %f --> %f):\n", yref, lambdaref,
1002  cpl_polynomial_eval_1d(py, lambdaref, NULL));
1003  cpl_polynomial_dump(py, stdout);
1004  fflush(stdout);
1005 #endif
1006 
1007  /* take the simple approach: search for the largest shift */
1008  double xmin = DBL_MAX, xmax = -DBL_MAX,
1009  ymin = DBL_MAX, ymax = -DBL_MAX;
1010  for (l = 0; l < nlambda; l++) {
1011  double lambda = cpl_vector_get(vlambda, l),
1012  xshift = cpl_polynomial_eval_1d(px, lambda, NULL),
1013  yshift = cpl_polynomial_eval_1d(py, lambda, NULL);
1014  if (xshift < xmin) xmin = xshift;
1015  if (xshift > xmax) xmax = xshift;
1016  if (yshift < ymin) ymin = yshift;
1017  if (yshift > ymax) ymax = yshift;
1018 #if DEBUG_MUSE_DARCHECK
1019  printf("%.1f Angstrom plane %d all-fitted: %f %f (%s)\n",
1020  lambda, l+1, xshift, yshift, __func__);
1021  fflush(stdout);
1022 #endif
1023  } /* for l (all wavelengths) */
1024  cpl_vector_delete(vlambda);
1025  cpl_matrix_delete(moffsets);
1026 
1027  /* we now have the extrema, compute the max shift per axis */
1028  double maxdiffx = xmax - xmin,
1029  maxdiffy = ymax - ymin;
1030  /* correct from shift in pixels to shift in arcsec */
1031  if (muse_pfits_get_mode(aPixtable->header) <= MUSE_MODE_WFM_AO_N) {
1032  maxdiffx *= kMuseSpaxelSizeX_WFM;
1033  maxdiffy *= kMuseSpaxelSizeY_WFM;
1034  } else {
1035  maxdiffx *= kMuseSpaxelSizeX_NFM;
1036  maxdiffy *= kMuseSpaxelSizeY_NFM;
1037  }
1038  *aMaxShift = fmax(maxdiffx, maxdiffy);
1039  cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_CHECK,
1040  *aMaxShift);
1041  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_CHECK,
1042  MUSE_HDR_PT_DAR_CHECK_C);
1043 
1044  muse_pixgrid_delete(grid);
1045  if (!aCube) {
1046  muse_datacube_delete(cube);
1047  }
1048 
1049  /* carry out residual DAR correction, if requested */
1050  if (aCorrect) {
1051  cpl_msg_debug(__func__, "Correcting pixel table for the shift (max = %f "
1052  "arcsec)", *aMaxShift);
1053  float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
1054  *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
1055  *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA);
1056  cpl_size i, nrow = muse_pixtable_get_nrow(aPixtable);
1057  #pragma omp parallel for default(none) \
1058  shared(lbda, nrow, px, py, xpos, ypos)
1059  for (i = 0; i < nrow; i++) {
1060  /* correct coordinates directly in the pixel table */
1061  xpos[i] += cpl_polynomial_eval_1d(px, lbda[i], NULL);
1062  ypos[i] += cpl_polynomial_eval_1d(py, lbda[i], NULL);
1063  } /* for i (all pixel table rows) */
1064 
1065  /* add the header */
1066  cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_CORR,
1067  lambdaref);
1068  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_CORR,
1069  MUSE_HDR_PT_DAR_CORR_C);
1070 
1071  /* need to recompute the (spatial) pixel table limits, but save the *
1072  * previous contents, if that was not already done */
1073  if (!cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_PREDAR_XLO)) {
1074  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XLO,
1075  cpl_propertylist_get_float(aPixtable->header,
1076  MUSE_HDR_PT_XLO));
1077  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XHI,
1078  cpl_propertylist_get_float(aPixtable->header,
1079  MUSE_HDR_PT_XHI));
1080  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YLO,
1081  cpl_propertylist_get_float(aPixtable->header,
1082  MUSE_HDR_PT_YLO));
1083  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YHI,
1084  cpl_propertylist_get_float(aPixtable->header,
1085  MUSE_HDR_PT_YHI));
1086  }
1087  muse_pixtable_compute_limits(aPixtable);
1088  } /* if aCorrect */
1089  cpl_polynomial_delete(px);
1090  cpl_polynomial_delete(py);
1091  /* reset cosmic ray statuses in aPixtable, since the CR rejection *
1092  * done here might not be appropriate for the final datacube */
1093  muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
1094 
1095  return CPL_ERROR_NONE;
1096 } /* muse_dar_check() */
1097 
muse_pixtable_wcs
State of the astrometric calibration of a MUSE pixel table.
Structure definition of a MUSE datacube.
Definition: muse_datacube.h:47
#define MUSE_HDR_PT_PREDAR_YLO
Structure definition for a collection of muse_images.
cpl_error_code muse_wcs_get_scales(cpl_propertylist *aHeader, double *aXScale, double *aYScale)
Compute the spatial scales (in degrees) from the FITS header WCS.
Definition: muse_wcs.c:1492
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
Definition: muse_image.c:85
#define MUSE_HDR_PT_XLO
FITS header keyword contains the lower limit of the data in x-direction.
double muse_pfits_get_rhum(const cpl_propertylist *aHeaders)
find out the relavtive humidity (in %)
Definition: muse_pfits.c:877
cpl_image * data
the data extension
Definition: muse_image.h:46
cpl_size muse_pixtable_get_nrow(const muse_pixtable *aPixtable)
get the number of rows within the pixel table
cpl_error_code muse_dar_correct(muse_pixtable *aPixtable, double aLambdaRef)
Correct the pixel coordinates of all pixels of a given pixel table for differential atmospheric refra...
Definition: muse_dar.c:306
muse_image * muse_combine_median_create(muse_imagelist *aImages)
Median combine a list of input images.
Definition: muse_combine.c:316
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:895
cpl_image * stat
the statistics extension
Definition: muse_image.h:64
void muse_datacube_delete(muse_datacube *aCube)
Deallocate memory associated to a muse_datacube object.
void muse_imagelist_delete(muse_imagelist *aList)
Free the memory of the MUSE image list.
cpl_error_code muse_dar_check(muse_pixtable *aPixtable, double *aMaxShift, cpl_boolean aCorrect, muse_datacube **aCube)
check that the continuum of objects in the frame is well-aligned perpendicular to the spatial axis...
Definition: muse_dar.c:641
The pixel grid.
Definition: muse_pixgrid.h:56
double muse_astro_parangle(const cpl_propertylist *aHeader)
Properly average parallactic angle values of one exposure.
Definition: muse_astro.c:449
Structure definition of MUSE three extension FITS file.
Definition: muse_image.h:40
cpl_table * table
The pixel table.
cpl_propertylist * header
the FITS header
Definition: muse_image.h:72
cpl_error_code muse_pixtable_reset_dq(muse_pixtable *aPixtable, unsigned int aDQ)
Reset a given bad pixel status (DQ flag) for all pixels in the table.
cpl_error_code muse_datacube_save(muse_datacube *aCube, const char *aFilename)
Save the three cube extensions and the FITS headers of a MUSE datacube to a file. ...
cpl_image * dq
the data quality extension
Definition: muse_image.h:56
#define MUSE_HDR_PT_PREDAR_YHI
Structure definition of MUSE pixel table.
muse_pixtable_wcs muse_pixtable_wcs_check(muse_pixtable *aPixtable)
Check the state of the world coordinate system of a pixel table.
#define MUSE_HDR_PT_YLO
FITS header keyword contains the lower limit of the data in y-direction.
cpl_error_code muse_utils_fit_moffat_2d(const cpl_matrix *aPositions, const cpl_vector *aValues, const cpl_vector *aErrors, cpl_array *aParams, cpl_array *aPErrors, const cpl_array *aPFlags, double *aRMS, double *aRedChisq)
Fit a 2D Moffat function to a given set of data.
Definition: muse_utils.c:1822
muse_resampling_params * muse_resampling_params_new(muse_resampling_type aMethod)
Create the resampling parameters structure.
cpl_polynomial * muse_utils_iterate_fit_polynomial(cpl_matrix *aPos, cpl_vector *aVal, cpl_vector *aErr, cpl_table *aExtra, const unsigned int aOrder, const double aRSigma, double *aMSE, double *aChiSq)
Iterate a polynomial fit.
Definition: muse_utils.c:2178
#define MUSE_HDR_PT_PREDAR_XLO
cpl_imagelist * data
the cube containing the actual data values
Definition: muse_datacube.h:75
#define MUSE_HDR_PT_PREDAR_XHI
static const cpl_size * muse_pixgrid_get_rows(muse_pixgrid *aPixels, cpl_size aIndex)
Return a pointer to the rows stored in one pixel.
Definition: muse_pixgrid.h:150
muse_datacube * muse_resampling_cube(muse_pixtable *aPixtable, muse_resampling_params *aParams, muse_pixgrid **aPixgrid)
Resample a pixel table onto a regular grid structure representing a FITS NAXIS=3 datacube.
cpl_imagelist * dq
the optional cube containing the bad pixel status
Definition: muse_datacube.h:80
void muse_pixgrid_delete(muse_pixgrid *aPixels)
Delete a pixgrid and remove its memory.
Definition: muse_pixgrid.c:398
cpl_propertylist * header
the FITS header
Definition: muse_datacube.h:56
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_pres_end(const cpl_propertylist *aHeaders)
find out the ambient pressure at end of exposure (in mbar)
Definition: muse_pfits.c:913
cpl_error_code muse_image_save(muse_image *aImage, const char *aFilename)
Save the three image extensions and the FITS headers of a MUSE image to a file.
Definition: muse_image.c:399
cpl_error_code muse_image_reject_from_dq(muse_image *aImage)
Reject pixels of a muse_image depending on its DQ data.
Definition: muse_image.c:856
muse_imagelist * muse_imagelist_new(void)
Create a new (empty) MUSE image list.
#define MUSE_HDR_PT_DAR_NAME
double muse_pfits_get_temp(const cpl_propertylist *aHeaders)
find out the ambient temperature (in degrees Celsius)
Definition: muse_pfits.c:859
cpl_error_code muse_utils_image_get_centroid_window(cpl_image *aImage, int aX1, int aY1, int aX2, int aY2, double *aX, double *aY, muse_utils_centroid_type aBgType)
Compute centroid over an image window, optionally marginalizing over the background.
Definition: muse_utils.c:1178
Resampling parameters.
void muse_resampling_params_delete(muse_resampling_params *aParams)
Delete a resampling parameters structure.
#define MUSE_HDR_PT_DAR_CORR
muse_image * muse_image_new(void)
Allocate memory for a new muse_image object.
Definition: muse_image.c:66
#define MUSE_HDR_PT_DAR_CHECK
static cpl_size muse_pixgrid_get_index(muse_pixgrid *aPixels, cpl_size aX, cpl_size aY, cpl_size aZ, cpl_boolean aAllowOutside)
Get the grid index determined from all three coordinates.
Definition: muse_pixgrid.h:89
#define MUSE_HDR_PT_YHI
FITS header keyword contains the upper limit of the data in y-direction.
cpl_error_code muse_imagelist_set(muse_imagelist *aList, muse_image *aImage, unsigned int aIdx)
Set the muse_image of given list index.
cpl_error_code muse_pixtable_compute_limits(muse_pixtable *aPixtable)
(Re-)Compute the limits of the coordinate columns of a pixel table.
static cpl_size muse_pixgrid_get_count(muse_pixgrid *aPixels, cpl_size aIndex)
Return the number of rows stored in one pixel.
Definition: muse_pixgrid.h:129
int muse_quality_image_reject_using_dq(cpl_image *aData, cpl_image *aDQ, cpl_image *aStat)
Reject pixels of one or two images on a DQ image.
Definition: muse_quality.c:628
muse_ins_mode muse_pfits_get_mode(const cpl_propertylist *aHeaders)
find out the observation mode
Definition: muse_pfits.c:1097
cpl_error_code muse_utils_get_centroid(const cpl_matrix *aPositions, const cpl_vector *aValues, const cpl_vector *aErrors, double *aX, double *aY, muse_utils_centroid_type aBgType)
Compute centroid of a two-dimensional dataset.
Definition: muse_utils.c:1272
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
cpl_imagelist * stat
the cube containing the data variance
Definition: muse_datacube.h:85
cpl_propertylist * header
The FITS header.
#define MUSE_HDR_PT_XHI
FITS header keyword contains the upper limit of the data in x-ion.