MUSE Pipeline Reference Manual  2.1.1
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-2015 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_phys.h"
41 #include "muse_quality.h"
42 #include "muse_utils.h"
43 #include "muse_wcs.h"
44 
45 /*----------------------------------------------------------------------------*
46  * Debugging/feature Macros *
47  *----------------------------------------------------------------------------*/
48 #define DEBUG_MUSE_DARCORRECT 0 /* if 1, output debugging in muse_dar_correct() */
49 #define DEBUG_MUSE_DARCHECK 0 /* if 1, output debug stuff in muse_dar_check(), *
50  * if > 1 also write files */
51 #define MUSE_DARCHECK_GRID_FITS 0 /* if 1, use grid-based fitting in muse_dar_check(); *
52  * it is more robust and slightly more accurate, but *
53  * about 10 times slower! */
54 
55 /*----------------------------------------------------------------------------*/
59 /*----------------------------------------------------------------------------*/
60 
63 /*----------------------------------------------------------------------------*/
137 /*----------------------------------------------------------------------------*/
138 cpl_error_code
139 muse_dar_correct(muse_pixtable *aPixtable, double aLambdaRef)
140 {
141  cpl_ensure_code(aPixtable && aPixtable->header, CPL_ERROR_NULL_INPUT);
142  if (aLambdaRef > 22000. || aLambdaRef < 3500.) {
143  cpl_msg_info(__func__, "Reference lambda %.3f Angstrom: skipping DAR "
144  "correction", aLambdaRef);
145  cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_NAME, -1.);
146  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_NAME,
147  MUSE_HDR_PT_DAR_C_SKIP);
148  return CPL_ERROR_NONE;
149  }
150  if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_NAME) &&
151  cpl_propertylist_get_double(aPixtable->header, MUSE_HDR_PT_DAR_NAME) > 0.) {
152  cpl_msg_info(__func__, "pixel table already corrected: skipping DAR "
153  "correction");
154  return CPL_ERROR_NONE;
155  }
156 
157  /* -------------------------------------------------------------- *
158  * get environmental and instrumental values from the FITS header *
159  * -------------------------------------------------------------- */
160  double airmass = muse_astro_airmass(aPixtable->header);
161  cpl_ensure_code(airmass >= 1., cpl_error_get_code());
162  /* simple zenith distance in radians */
163  double z = acos(1./airmass);
164 
165  /* compute the mean parallactic angle during exposure */
166  cpl_errorstate prestate = cpl_errorstate_get();
167  double parang = muse_astro_parangle(aPixtable->header);
168  cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
169 
170  /* compute the position angle on the sky from the angles we have */
171  prestate = cpl_errorstate_get();
172  double posang = muse_astro_posangle(aPixtable->header);
173  cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
174 
175  cpl_boolean isWFM = muse_pfits_get_mode(aPixtable->header)
176  < MUSE_MODE_NFM_AO_N;
177  if (!isWFM) {
178  cpl_msg_warning(__func__, "For NFM instrument mode DAR correction should "
179  "not be needed!");
180  }
181 
182  /* query values and set defaults in case we didn't get proper values */
183  prestate = cpl_errorstate_get();
184  double temp = muse_pfits_get_temp(aPixtable->header); /* temperature [Celsius] */
185  if (cpl_errorstate_is_equal(prestate)) {
186  temp += 273.15; /* T[K] */
187  } else {
188  cpl_msg_warning(__func__, "FITS header does not contain temperature: %s",
189  cpl_error_get_message());
190  temp = 11.5 + 273.15; /* T[K] */
191  }
192  prestate = cpl_errorstate_get();
193  double rhum = muse_pfits_get_rhum(aPixtable->header); /* relative humidity [%] */
194  if (!cpl_errorstate_is_equal(prestate)) {
195  cpl_msg_warning(__func__, "FITS header does not contain relative humidity: %s",
196  cpl_error_get_message());
197  rhum = 14.5;
198  }
199  rhum /= 100.; /* convert from % to fraction */
200  prestate = cpl_errorstate_get();
201  double pres = (muse_pfits_get_pres_start(aPixtable->header) /* pressure [mbar] */
202  + muse_pfits_get_pres_end(aPixtable->header)) / 2.;
203  if (!cpl_errorstate_is_equal(prestate)) {
204  cpl_msg_warning(__func__, "FITS header does not contain pressure values: %s",
205  cpl_error_get_message());
206  pres = 743.;
207  }
208 
209  /* -------------------------------------------------------------------------- *
210  * choose method and compute the refractive index at the reference wavelength *
211  * -------------------------------------------------------------------------- */
212  enum {
213  MUSE_DAR_METHOD_FILIPPENKO = 0,
214  MUSE_DAR_METHOD_OWENS,
215  MUSE_DAR_METHOD_EDLEN,
216  MUSE_DAR_METHOD_CIDDOR,
217  };
218  int method = MUSE_DAR_METHOD_FILIPPENKO;
219  if (getenv("MUSE_DAR_CORRECT_METHOD") &&
220  !strncmp(getenv("MUSE_DAR_CORRECT_METHOD"), "Owens", 6)) {
221  method = MUSE_DAR_METHOD_OWENS;
222  } else if (getenv("MUSE_DAR_CORRECT_METHOD") &&
223  !strncmp(getenv("MUSE_DAR_CORRECT_METHOD"), "Edlen", 6)) {
224  method = MUSE_DAR_METHOD_EDLEN;
225  } else if (getenv("MUSE_DAR_CORRECT_METHOD") &&
226  !strncmp(getenv("MUSE_DAR_CORRECT_METHOD"), "Ciddor", 7)) {
227  method = MUSE_DAR_METHOD_CIDDOR;
228  }
229  /* compute refractive index for reference wavelength in um and *
230  * output properties in "natural" (for the formulae) units */
231  double d1, d2,
232  fp, /* water vapor pressure in mmHg */
233  nr0, /* refractive index of air at reference wavelength */
234  p = pres * 100, /* atmopheric pressure [Pa] */
235  t = temp - 273.15, /* temperature [Celsius] */
236  pv, /* partial vapor pressure */
237  xv; /* mole fraction */
238  if (method == MUSE_DAR_METHOD_OWENS) {
239  muse_phys_nrindex_owens_coeffs(temp, rhum, pres, &d1, &d2);
240  nr0 = muse_phys_nrindex_owens(aLambdaRef / 10000., d1, d2);
241  cpl_msg_info(__func__, "DAR for T=%.2f K, RH=%.2f %%, pres=%.1f mbar, "
242  "at airmass=%.4f, using ref. wavelength %.6f um (Owens)",
243  temp, rhum*100., pres, 1./cos(z), aLambdaRef / 10000.);
244  } else if (method == MUSE_DAR_METHOD_EDLEN ||
245  method == MUSE_DAR_METHOD_CIDDOR) {
246  /* follow the detailed procedure from *
247  * http://emtoolbox.nist.gov/Wavelength/Documentation.asp#AppendixA */
248  /* Calculate saturation vapor pressure psv(t) over water,
249  * equations (A1) to (A7) */
250  const double T = temp, /* temperature [K] */
251  K1 = 1.16705214528E+03,
252  K2 = -7.24213167032E+05,
253  K3 = -1.70738469401E+01,
254  K4 = 1.20208247025E+04,
255  K5 = -3.23255503223E+06,
256  K6 = 1.49151086135E+01,
257  K7 = -4.82326573616E+03,
258  K8 = 4.05113405421E+05,
259  K9 = -2.38555575678E-01,
260  K10 = 6.50175348448E+02,
261  Omega = T + K9 / (T - K10),
262  A = Omega*Omega + K1 * Omega + K2,
263  B = K3 * Omega*Omega + K4 * Omega + K5,
264  C = K6 * Omega*Omega + K7 * Omega + K8,
265  X = -B + sqrt(B*B - 4 * A * C);
266  double psv_w = 1e6 * pow(2 * C / X, 4);
267  /* Calculate saturation vapor pressure psv(t) over ice, *
268  * equations (A8) tl (A13) */
269  const double A1 = -13.928169,
270  A2 = 34.7078238,
271  Theta = T / 273.16,
272  Y = A1 * (1 - pow(Theta, -1.5)) + A2 * (1 - pow(Theta, -1.25));
273  double psv_i = 611.657 * exp(Y);
274  /* determine Humidity */
275  /* for the Edlén equation: for relative humidity RH [%], calculate *
276  * partial pressure using psw_w for t > 0 [Celsius] and psw_i for t < 0 */
277  if (t > 0.) {
278  pv = rhum * psv_w; /* use saturation pressure over water */
279  } else {
280  pv = rhum * psv_i; /* use saturation pressure over ice */
281  }
282  /* for the Ciddor equation: express humidity as a mole fraction */
283  const double alpha = 1.00062,
284  beta = 3.14e-8,
285  gamma = 5.60e-7;
286  double f = alpha + beta * p + gamma * t*t; /* enhancement factor f(p,t) */
287  /* for relative humidity RH [%], calculate mole fraction xv using psw(t) */
288  if (t > 0.) {
289  xv = rhum * f * psv_w / p; /* use saturation pressure over water */
290  } else {
291  xv = rhum * f * psv_i / p; /* use saturation pressure over ice */
292  }
293  if (method == MUSE_DAR_METHOD_EDLEN) {
294  nr0 = muse_phys_nrindex_edlen(aLambdaRef / 10000., t, p, pv);
295  } else { /* MUSE_DAR_METHOD_CIDDOR */
296  nr0 = muse_phys_nrindex_ciddor(aLambdaRef / 10000., T, p, xv, 450.);
297  }
298  cpl_msg_info(__func__, "DAR for T=%.2f degC, RH=%.2f %%, p=%.1f Pa, "
299  "at airmass=%.4f, using ref. wavelength %.6f um (%s)",
300  t, rhum*100., p, 1./cos(z), aLambdaRef / 10000.,
301  method == MUSE_DAR_METHOD_EDLEN ? "Edlen" : "Ciddor");
302  } else {
303  /* use the Owens formula to derive saturation pressure, still needs T[K] */
305  /* using that, derive the water vapor pressure in mmHg */
306  fp = rhum * ps * MUSE_PHYS_hPa_TO_mmHg;
307  temp -= 273.15; /* temperature (again) in degrees Celsius */
308  pres *= MUSE_PHYS_hPa_TO_mmHg; /* need the pressure in mmHg */
309  nr0 = muse_phys_nrindex_filippenko(aLambdaRef / 10000., temp, fp, pres);
310  cpl_msg_info(__func__, "DAR for T=%.2f degC, fp=%.3f mmHg, P=%.1f mmHg, "
311  "at airmass=%.4f, using ref. wavelength %.6f um (Filippenko)",
312  temp, fp, pres, 1./cos(z), aLambdaRef / 10000.);
313  }
314 
315  /* ---------------------------------- *
316  * compute the direction of the shift *
317  * ---------------------------------- */
318  /* xshift is in E-W direction for posang = 0, yshift is N-S */
319  double xshift = -sin((parang + posang) * CPL_MATH_RAD_DEG),
320  yshift = cos((parang + posang) * CPL_MATH_RAD_DEG);
321 
322  /* apply correction for spaxel size */
323  muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
324  cpl_ensure_code(wcstype == MUSE_PIXTABLE_WCS_PIXEL ||
325  wcstype == MUSE_PIXTABLE_WCS_CELSPH, CPL_ERROR_INVALID_TYPE);
326  double fscale = 3600.; /* assume scale in arcsec */
327  if (wcstype == MUSE_PIXTABLE_WCS_CELSPH) {
328  fscale = 1. / 3600; /* scale in degrees */
329  double xscale, yscale;
330  muse_wcs_get_scales(aPixtable->header, &xscale, &yscale);
331  xshift /= -xscale; /* need to reverse the sign for RA */
332  yshift /= yscale;
333  } else { /* assume nominal spatial scale */
334  if (isWFM) {
335  xshift /= kMuseSpaxelSizeX_WFM;
336  yshift /= kMuseSpaxelSizeY_WFM;
337  } else {
338  xshift /= kMuseSpaxelSizeX_NFM;
339  yshift /= kMuseSpaxelSizeY_NFM;
340  }
341  }
342 
343  /* ------------------------------ *
344  * apply correction to all pixels *
345  * ------------------------------ */
346  /* diff. refr. base in arcsec converted from radians (Filippenko does *
347  * the conversion using x206265 which converts radians to arcsec) */
348  double dr0 = tan(z) * fscale * CPL_MATH_DEG_RAD;
349 #if DEBUG_MUSE_DARCORRECT
350  double wl;
351  for (wl = 4650.; wl <= 9300; wl += 50) {
352  double nr;
353  if (method == MUSE_DAR_METHOD_OWENS) {
354  nr = muse_phys_nrindex_owens(wl / 10000., d1, d2);
355  } else if (method == MUSE_DAR_METHOD_EDLEN) {
356  nr = muse_phys_nrindex_edlen(wl / 10000., t, p, pv);
357  } else if (method == MUSE_DAR_METHOD_CIDDOR) {
358  nr = muse_phys_nrindex_ciddor(wl / 10000., temp, p, xv, 450.);
359  } else {
360  nr = muse_phys_nrindex_filippenko(wl / 10000., temp, fp, pres);
361  }
362  double dr = dr0 * (nr0 - nr);
363  printf("%.1f Angstrom: %f ==> %f %f %s (%s)\n", wl, dr,
364  xshift * dr, yshift * dr,
365  wcstype == MUSE_PIXTABLE_WCS_CELSPH ? "deg" : "pix", __func__);
366  }
367  fflush(stdout);
368 #endif
369  float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
370  *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
371  *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA);
372  cpl_size i, nmax = muse_pixtable_get_nrow(aPixtable);
373  #pragma omp parallel for default(none) \
374  shared(d1, d2, dr0, fp, lbda, method, nmax, nr0, p, pres, pv, t, \
375  temp, xpos, ypos, xshift, xv, yshift)
376  for (i = 0; i < nmax; i++) {
377  double nr, lambda = lbda[i] * 1e-4;
378  if (method == MUSE_DAR_METHOD_OWENS) {
379  nr = muse_phys_nrindex_owens(lambda, d1, d2);
380  } else if (method == MUSE_DAR_METHOD_EDLEN) {
381  nr = muse_phys_nrindex_edlen(lambda, t, p, pv);
382  } else if (method == MUSE_DAR_METHOD_CIDDOR) {
383  nr = muse_phys_nrindex_ciddor(lambda, temp, p, xv, 450.);
384  } else {
385  nr = muse_phys_nrindex_filippenko(lambda, temp, fp, pres);
386  }
387  double dr = dr0 * (nr0 - nr);
388  /* correct coordinates directly in the pixel table */
389  xpos[i] += xshift * dr;
390  ypos[i] += yshift * dr;
391  } /* for i (all pixel table rows) */
392  /* add the header */
393  cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_NAME,
394  aLambdaRef);
395  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_NAME,
396  MUSE_HDR_PT_DAR_COMMENT);
397  /* need to recompute the (spatial) pixel table limits, but save the *
398  * previous contents because that's what we need for WCS calibration */
399  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XLO,
400  cpl_propertylist_get_float(aPixtable->header,
401  MUSE_HDR_PT_XLO));
402  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XHI,
403  cpl_propertylist_get_float(aPixtable->header,
404  MUSE_HDR_PT_XHI));
405  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YLO,
406  cpl_propertylist_get_float(aPixtable->header,
407  MUSE_HDR_PT_YLO));
408  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YHI,
409  cpl_propertylist_get_float(aPixtable->header,
410  MUSE_HDR_PT_YHI));
411  muse_pixtable_compute_limits(aPixtable);
412 
413  return CPL_ERROR_NONE;
414 } /* muse_dar_correct() */
415 
416 /*----------------------------------------------------------------------------*/
470 /*----------------------------------------------------------------------------*/
471 cpl_error_code
472 muse_dar_check(muse_pixtable *aPixtable, double *aMaxShift,
473  cpl_boolean aCorrect, muse_datacube **aCube)
474 {
475  cpl_ensure_code(aMaxShift, CPL_ERROR_NULL_INPUT);
476  *aMaxShift = -99.;
477  cpl_ensure_code(aPixtable, CPL_ERROR_NULL_INPUT);
478  muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
479  cpl_ensure_code(wcstype == MUSE_PIXTABLE_WCS_PIXEL, CPL_ERROR_INVALID_TYPE);
480  if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_CHECK)) {
481  cpl_msg_info(__func__, "pixel table already checked for DAR accuracy");
482  }
483  cpl_boolean corrected = CPL_FALSE;
484  if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_CORR)) {
485  cpl_msg_info(__func__, "pixel table already corrected for DAR residuals");
486  corrected = CPL_TRUE;
487  }
488  /* set the halfsize of the window to measure the centroids with respect to *
489  * the reference wavelength; +/-5 pix should suffice in case DAR was *
490  * previously corrected; in other cases, increase it to +/-15 pix */
491  int cenhalfsize = 5;
492  /* polynomial order for the fits, default 2, if not corrected yet, use 4; *
493  * here use the existence of the header keyword only as a first guess */
494  cpl_boolean residuals = cpl_propertylist_has(aPixtable->header,
496  int order = residuals || corrected ? 2 : 4;
497 
498  /* create cube, with coarse sampling in wavelength */
499  cpl_msg_info(__func__, "Intermediate resampling to %s atmospheric refraction"
500  "%s (using polynomial order %d):", aCorrect ? "correct" : "check",
501  residuals ? " residuals" : "", order);
502  cpl_msg_indent_more();
503  muse_resampling_params *params =
505  params->pfx = 1.; /* large pixfrac to not create sudden spatial shifts */
506  params->pfy = 1.;
507  params->pfl = 1.;
508  params->dlambda = 10.;
509  params->crsigma = 7.; /* moderate CR rejection */
510  muse_pixgrid *grid;
511  muse_datacube *cube = muse_resampling_cube(aPixtable, params, &grid);
513  cpl_msg_indent_less();
514  if (!cube) {
515  muse_pixgrid_delete(grid);
516  muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
517  return cpl_error_set_message(__func__, cpl_error_get_code(),
518  "Failure while creating cube!");
519  }
520  if (aCube) {
521  *aCube = cube;
522  }
523  double crpix3 = muse_pfits_get_crpix(cube->header, 3),
524  crval3 = muse_pfits_get_crval(cube->header, 3),
525  cdelt3 = muse_pfits_get_cd(cube->header, 3, 3);
526  cpl_errorstate state = cpl_errorstate_get(); /* header may be missing */
527  double lambdaref = cpl_propertylist_get_double(aPixtable->header,
529  lambda1 = (2 - crpix3) * cdelt3 + crval3, /* 2nd plane */
530  lambda3 = (grid->nz - 1 - crpix3) * cdelt3 + crval3, /* next to last */
531  lambda2 = (lambda1 + lambda3) / 2.; /* central plane */
532  if (!cpl_errorstate_is_equal(state) || lambdaref < 0) {
533  cpl_errorstate_set(state);
534  /* use reference wavelength in the center of the wavelength range */
535  lambdaref = lambda2;
536  if (!cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_CORR)) {
537  cenhalfsize = 15;
538  order = 4;
539  }
540  cpl_msg_warning(__func__, "Data is not DAR corrected, using reference "
541  "wavelength at %.3f Angstrom and polynomial order %d!",
542  lambdaref, order);
543  }
544  if (lambdaref > lambda3 || lambdaref < lambda1) { /* ref. outside cube! */
545  lambdaref = lambda2;
546  cpl_msg_warning(__func__, "Reference lambda outside cube, using reference "
547  "wavelength at %.3f Angstrom!", lambdaref);
548  }
549 
550  /* detect objects in the cube, using image planes near reference wavelength */
551  int iplane, irefplane = lround((lambdaref - crval3) / cdelt3 + crpix3) - 1;
552  /* medianing three images removes all cosmic rays but continuum objects stay *
553  * (need to use muse_images and the muse_combine function because *
554  * cpl_imagelist_collapse_median_create() disregards all bad pixels) */
556  unsigned int ilist = 0;
557  for (iplane = irefplane - 1; iplane <= irefplane + 1; iplane++) {
558  muse_image *image = muse_image_new();
559  image->data = cpl_image_duplicate(cpl_imagelist_get(cube->data, iplane));
560  image->dq = cpl_image_duplicate(cpl_imagelist_get(cube->dq, iplane));
561  image->stat = cpl_image_duplicate(cpl_imagelist_get(cube->stat, iplane));
562  muse_imagelist_set(list, image, ilist++);
563  } /* for iplane (planes around ref. wavelength) */
564  muse_image *median = muse_combine_median_create(list);
565  muse_imagelist_delete(list);
566  if (!median) {
567  if (!aCube) {
568  muse_datacube_delete(cube);
569  }
570  muse_pixgrid_delete(grid);
571  muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
572  return cpl_error_set_message(__func__, cpl_error_get_code(),
573  "Failure while creating detection image!");
574  }
575 #if DEBUG_MUSE_DARCHECK > 1
576  median->header = cpl_propertylist_new(); /* suppress errors while saving */
577  cpl_propertylist_append_string(median->header, "BUNIT",
578  cpl_propertylist_get_string(cube->header,
579  "BUNIT"));
580  muse_image_save(median, "IMAGE_darcheck_median.fits");
581 #endif
582  /* reject data in image to signify bad pixels to CPL routines */
584  /* use high sigmas for detection */
585  double dsigmas[] = { 50., 30., 10., 8., 6., 5. };
586  int ndsigmas = sizeof(dsigmas) / sizeof(double);
587  cpl_vector *vsigmas = cpl_vector_wrap(ndsigmas, dsigmas);
588  cpl_size isigma = -1;
589  state = cpl_errorstate_get();
590  cpl_apertures *apertures = cpl_apertures_extract(median->data, vsigmas, &isigma);
591  int nx = cpl_image_get_size_x(median->data),
592  ny = cpl_image_get_size_y(median->data);
593  muse_image_delete(median);
594  int napertures = apertures ? cpl_apertures_get_size(apertures) : 0;
595  if (napertures < 1) {
596  cpl_vector_unwrap(vsigmas);
597  cpl_apertures_delete(apertures);
598  cpl_errorstate_set(state); /* reset state to before aperture extraction */
599  if (!aCube) {
600  muse_datacube_delete(cube);
601  }
602  muse_pixgrid_delete(grid);
603  muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
604  return cpl_error_set_message(__func__, CPL_ERROR_DATA_NOT_FOUND, "No sources "
605  "found for DAR check down to %.1f sigma "
606  "limit on plane %d", dsigmas[ndsigmas - 1],
607  irefplane + 1);
608  }
609  cpl_msg_debug(__func__, "The %.1f sigma threshold was used to find %d source%s,"
610  " centered on plane %d", cpl_vector_get(vsigmas, isigma),
611  napertures, napertures == 1 ? "" : "s", iplane+1);
612  cpl_vector_unwrap(vsigmas);
613 #if DEBUG_MUSE_DARCHECK
614  cpl_apertures_dump(apertures, stdout);
615  fflush(stdout);
616 #if DEBUG_MUSE_DARCHECK > 1
617  muse_datacube_save(cube, "DATACUBE_darcheck.fits");
618 #endif
619 #endif
620 
621  /* restrict box sizes to be at most the size of the image plane */
622  int xhalfsize = 2*cenhalfsize > nx ? nx/2 : cenhalfsize,
623  yhalfsize = 2*cenhalfsize > ny ? ny/2 : cenhalfsize;
624 
625  /* loop through the wavelength planes and compute centroids for all apertures */
626  int l, nlambda = cpl_imagelist_get_size(cube->data);
627  cpl_vector *vlambda = cpl_vector_new(nlambda);
628  cpl_matrix *moffsets = cpl_matrix_new(2, nlambda);
629 #if !DEBUG_MUSE_DARCHECK
630  #pragma omp parallel for default(none) \
631  shared(apertures, aPixtable, cdelt3, crpix3, crval3, cube, grid, \
632  moffsets, napertures, nlambda, nx, ny, vlambda, xhalfsize, \
633  yhalfsize)
634 #endif
635  for (l = 0; l < nlambda; l++) {
636  double xoffset = 0., yoffset = 0.,
637  lambda = (l+1 - crpix3) * cdelt3 + crval3;
638  int n, noffsets = 0;
639  for (n = 1; n <= napertures; n++) {
640  /* use detection center to construct (larger) window for centroid */
641  double xc = cpl_apertures_get_centroid_x(apertures, n),
642  yc = cpl_apertures_get_centroid_y(apertures, n);
643  int x1 = lround(xc) - xhalfsize,
644  x2 = lround(xc) + xhalfsize,
645  y1 = lround(yc) - yhalfsize,
646  y2 = lround(yc) + yhalfsize;
647  /* force window to be inside the image */
648  if (x1 < 1) x1 = 1;
649  if (y1 < 1) y1 = 1;
650  if (x2 > nx) x2 = nx;
651  if (y2 > ny) y2 = ny;
652 
653 #define MUSE_DARCHECK_MAX_INTERNAL_SHIFT 5 /* accept max. 5 pix shift between *
654  * image centroid and fit center */
655 #if MUSE_DARCHECK_GRID_FITS /* carry out centroid and fit on the un-resampled *
656  * data, i.e on the pixel grid for better S/N */
657 #if DEBUG_MUSE_DARCHECK
658  const char *method = "grid/moffat";
659 #endif
660  /* transfer data from the pixel table into matrices */
661  float *cdata = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
662  *cstat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
663  int *cdq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
664 #define MUSE_DARCHECK_MAX_NPIX 7000 /* start and increase in fixed size bins */
665  cpl_matrix *pos = cpl_matrix_new(MUSE_DARCHECK_MAX_NPIX, 2);
666  cpl_vector *val = cpl_vector_new(MUSE_DARCHECK_MAX_NPIX),
667  *err = cpl_vector_new(MUSE_DARCHECK_MAX_NPIX);
668  int i, npix = 0, nsize = MUSE_DARCHECK_MAX_NPIX;
669  for (i = x1 - 1; i < x2; i++) {
670  int j;
671  for (j = y1 - 1; j < y2; j++) {
672  cpl_size idx = muse_pixgrid_get_index(grid, i, j, l, CPL_TRUE),
673  ipx, npx = muse_pixgrid_get_count(grid, idx);
674  const cpl_size *rows = muse_pixgrid_get_rows(grid, idx);
675  for (ipx = 0; ipx < npx; ipx++) {
676  if (cdq[rows[ipx]] != EURO3D_GOODPIXEL) {
677  continue;
678  }
679  if (npix >= nsize) { /* increase size of the buffers */
680  nsize += MUSE_DARCHECK_MAX_NPIX;
681  cpl_matrix_resize(pos, 0, nsize - cpl_matrix_get_nrow(pos), 0, 0);
682  cpl_vector_set_size(val, nsize);
683  cpl_vector_set_size(err, nsize);
684  }
685  cpl_matrix_set(pos, npix, 0, i+1);
686  cpl_matrix_set(pos, npix, 1, j+1);
687  cpl_vector_set(val, npix, cdata[rows[ipx]]);
688  cpl_vector_set(err, npix, sqrt(cstat[rows[ipx]]));
689  npix++;
690  } /* for ipx */
691  } /* for j */
692  } /* for i */
693  if (npix < 8) { /* need 8 pixels for 8 Moffat parameters */
694  cpl_matrix_delete(pos);
695  cpl_vector_delete(val);
696  cpl_vector_delete(err);
697  continue;
698  }
699  cpl_matrix_set_size(pos, npix, 2);
700  cpl_vector_set_size(val, npix);
701  cpl_vector_set_size(err, npix);
702 
703  /* compute centroid first, as a kind of first guess */
704  double xc1, yc1;
705  muse_utils_get_centroid(pos, val, NULL/*err*/, &xc1, &yc1,
707 
708  /* now use a 2D Moffat fit for a refined position */
709  cpl_array *pmoffat = cpl_array_new(8, CPL_TYPE_DOUBLE);
710  cpl_array_set(pmoffat, 2, xc1);
711  cpl_array_set(pmoffat, 3, yc1);
712  cpl_error_code rc = muse_utils_fit_moffat_2d(pos, val, err, pmoffat, NULL,
713  NULL, NULL, NULL);
714  double xc2 = cpl_array_get(pmoffat, 2, NULL),
715  yc2 = cpl_array_get(pmoffat, 3, NULL);
716  cpl_array_delete(pmoffat);
717  cpl_matrix_delete(pos);
718  cpl_vector_delete(val);
719  cpl_vector_delete(err);
720 
721  /* if the fit worked and both methods give roughly similar *
722  * positions we should have a valid estimate of the center */
723  if (rc == CPL_ERROR_NONE &&
724  fabs(xc1 - xc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT &&
725  fabs(yc1 - yc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT) {
726  xoffset += xc - xc2;
727  yoffset += yc - yc2;
728  noffsets++;
729  }
730 #else /* non-MUSE_DARCHECK_GRID_FITS follows, using wavelength plane images */
731 #if DEBUG_MUSE_DARCHECK
732  const char *method = "image/gauss";
733 #endif
734  /* add CPL bad pixel mask in image plane to get reliable results */
735  cpl_image *plane = cpl_imagelist_get(cube->data, l),
736  *plerr = cpl_image_power_create(cpl_imagelist_get(cube->stat, l), 0.5);
737  /* make sure to exclude bad pixels from the fits below */
738  muse_quality_image_reject_using_dq(plane, cpl_imagelist_get(cube->dq, l), plerr);
739  /* the error image may contain more bad pixels! */
740  cpl_image_reject_from_mask(plane, cpl_image_get_bpm(plerr));
741 
742  cpl_image *xtr = cpl_image_extract(plane, x1, y1, x2, y2);
743  int npix = (x2-x1+1) * (y2-y1+1) /* number of (good) pixels involved */
744  - cpl_image_count_rejected(xtr);
745  cpl_image_delete(xtr);
746  if (npix < 7) { /* need 7 pixels for 7 Gaussian parameters */
747  cpl_image_delete(plerr);
748  continue;
749  }
750 
751  /* compute centroid first, as a kind of first guess */
752  double xc1, yc1;
753  muse_utils_image_get_centroid_window(plane, x1, y1, x2, y2, &xc1, &yc1,
755 
756  /* do Gaussian fit, assuming that the object is roughly stellar */
757  cpl_array *pgauss = cpl_array_new(7, CPL_TYPE_DOUBLE);
758  cpl_array_set(pgauss, 3, xc1);
759  cpl_array_set(pgauss, 4, yc1);
760  cpl_error_code rc = cpl_fit_image_gaussian(plane, plerr,
761  lround(xc), lround(yc),
762  2*xhalfsize, 2*yhalfsize,
763  pgauss, NULL, NULL,
764  NULL, NULL, NULL,
765  NULL, NULL, NULL, NULL);
766  double xc2 = cpl_array_get(pgauss, 3, NULL),
767  yc2 = cpl_array_get(pgauss, 4, NULL);
768  cpl_array_delete(pgauss);
769  cpl_image_delete(plerr);
770 
771  /* if the fit worked and both methods give roughly similar *
772  * positions we should have a valid estimate of the center */
773  if (rc == CPL_ERROR_NONE &&
774  fabs(xc1 - xc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT &&
775  fabs(yc1 - yc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT) {
776  xoffset += xc - xc2;
777  yoffset += yc - yc2;
778  noffsets++;
779  }
780 #endif /* end of MUSE_DARCHECK_GRID_FITS */
781 #if DEBUG_MUSE_DARCHECK /* DEBUG output; compare direct centroid and fit position */
782  printf("%.1f Angstrom plane %d aper %d (%f,%f): %f %f (%s) %f %f (%d, %s)\n",
783  lambda, l+1, n, xc, yc, xc - xc2, yc - yc2, __func__,
784  xc2 - xc1, yc2 - yc1, npix, method);
785  fflush(stdout);
786 #endif
787  } /* for n (all apertures) */
788  /* record average offset in the matrix */
789  cpl_vector_set(vlambda, l, lambda);
790  cpl_matrix_set(moffsets, 0, l, xoffset / noffsets);
791  cpl_matrix_set(moffsets, 1, l, yoffset / noffsets);
792 #if DEBUG_MUSE_DARCHECK
793  printf("%.1f Angstrom plane %d all-apers: %f %f (%s)\n",
794  lambda, l+1, xoffset / noffsets, yoffset / noffsets, __func__);
795  fflush(stdout);
796 #endif
797  } /* for l (all wavelengths) */
798  cpl_apertures_delete(apertures);
799 
800  /* compute correction polynomials and maximum relative offset */
801  cpl_vector *vlam1 = cpl_vector_duplicate(vlambda),
802  *vlam2 = cpl_vector_duplicate(vlambda);
803  cpl_matrix *mlam1 = cpl_matrix_wrap(1, nlambda, cpl_vector_unwrap(vlam1)),
804  *mlam2 = cpl_matrix_wrap(1, nlambda, cpl_vector_unwrap(vlam2));
805  cpl_matrix *mxoff = cpl_matrix_extract_row(moffsets, 0),
806  *myoff = cpl_matrix_extract_row(moffsets, 1);
807  cpl_vector *vxoff = cpl_vector_wrap(nlambda, cpl_matrix_unwrap(mxoff)),
808  *vyoff = cpl_vector_wrap(nlambda, cpl_matrix_unwrap(myoff));
809  /* use iterative polynomial fits, separately in both spatial dimensions, *
810  * to throw out the most extreme outliers, and compute the correction */
811  double msex, msey, chisqx, chisqy;
812  cpl_polynomial *px = muse_utils_iterate_fit_polynomial(mlam1, vxoff,
813  NULL, NULL, order, 3.,
814  &msex, &chisqx),
815  *py = muse_utils_iterate_fit_polynomial(mlam2, vyoff,
816  NULL, NULL, order, 3.,
817  &msey, &chisqy);
818  cpl_vector_delete(vxoff);
819  cpl_vector_delete(vyoff);
820  cpl_matrix_delete(mlam1);
821  cpl_matrix_delete(mlam2);
822 #if DEBUG_MUSE_DARCHECK
823  printf("polynomial fit in x (order %d, rms=%f, chisq=%g):\n", order,
824  sqrt(msex), chisqx);
825  cpl_polynomial_dump(px, stdout);
826  printf("polynomial fit in y (order %d, rms=%f, chisq=%g):\n", order,
827  sqrt(msey), chisqy);
828  cpl_polynomial_dump(py, stdout);
829  fflush(stdout);
830 #endif
831 
832  /* shift to the zeropoint of the offsets at the reference wavelength; *
833  * this should be more stable than just taking the raw input value */
834  double xref = cpl_polynomial_eval_1d(px, lambdaref, NULL),
835  yref = cpl_polynomial_eval_1d(py, lambdaref, NULL);
836  cpl_size pows = 0; /* set the zero-order coefficient */
837  cpl_polynomial_set_coeff(px, &pows, cpl_polynomial_get_coeff(px, &pows) - xref);
838  cpl_polynomial_set_coeff(py, &pows, cpl_polynomial_get_coeff(py, &pows) - yref);
839 #if DEBUG_MUSE_DARCHECK
840  printf("polynomial in x (xref=%f at %f --> %f):\n", xref, lambdaref,
841  cpl_polynomial_eval_1d(px, lambdaref, NULL));
842  cpl_polynomial_dump(px, stdout);
843  printf("polynomial in y (yref=%f at %f --> %f):\n", yref, lambdaref,
844  cpl_polynomial_eval_1d(py, lambdaref, NULL));
845  cpl_polynomial_dump(py, stdout);
846  fflush(stdout);
847 #endif
848 
849  /* take the simple approach: search for the largest shift */
850  double xmin = DBL_MAX, xmax = -DBL_MAX,
851  ymin = DBL_MAX, ymax = -DBL_MAX;
852  for (l = 0; l < nlambda; l++) {
853  double lambda = cpl_vector_get(vlambda, l),
854  xshift = cpl_polynomial_eval_1d(px, lambda, NULL),
855  yshift = cpl_polynomial_eval_1d(py, lambda, NULL);
856  if (xshift < xmin) xmin = xshift;
857  if (xshift > xmax) xmax = xshift;
858  if (yshift < ymin) ymin = yshift;
859  if (yshift > ymax) ymax = yshift;
860 #if DEBUG_MUSE_DARCHECK
861  printf("%.1f Angstrom plane %d all-fitted: %f %f (%s)\n",
862  lambda, l+1, xshift, yshift, __func__);
863  fflush(stdout);
864 #endif
865  } /* for l (all wavelengths) */
866  cpl_vector_delete(vlambda);
867  cpl_matrix_delete(moffsets);
868 
869  /* we now have the extrema, compute the max shift per axis */
870  double maxdiffx = xmax - xmin,
871  maxdiffy = ymax - ymin;
872  /* correct from shift in pixels to shift in arcsec */
873  if (muse_pfits_get_mode(aPixtable->header) < MUSE_MODE_NFM_AO_N) {
874  maxdiffx *= kMuseSpaxelSizeX_WFM;
875  maxdiffy *= kMuseSpaxelSizeY_WFM;
876  } else {
877  maxdiffx *= kMuseSpaxelSizeX_NFM;
878  maxdiffy *= kMuseSpaxelSizeY_NFM;
879  }
880  *aMaxShift = fmax(maxdiffx, maxdiffy);
881  cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_CHECK,
882  *aMaxShift);
883  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_CHECK,
884  MUSE_HDR_PT_DAR_CHECK_C);
885 
886  muse_pixgrid_delete(grid);
887  if (!aCube) {
888  muse_datacube_delete(cube);
889  }
890 
891  /* carry out residual DAR correction, if requested */
892  if (aCorrect) {
893  cpl_msg_debug(__func__, "Correcting pixel table for the shift (max = %f "
894  "arcsec)", *aMaxShift);
895  float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
896  *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
897  *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA);
898  cpl_size i, nrow = muse_pixtable_get_nrow(aPixtable);
899  #pragma omp parallel for default(none) \
900  shared(lbda, nrow, px, py, xpos, ypos)
901  for (i = 0; i < nrow; i++) {
902  /* correct coordinates directly in the pixel table */
903  xpos[i] += cpl_polynomial_eval_1d(px, lbda[i], NULL);
904  ypos[i] += cpl_polynomial_eval_1d(py, lbda[i], NULL);
905  } /* for i (all pixel table rows) */
906 
907  /* add the header */
908  cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_CORR,
909  lambdaref);
910  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_CORR,
911  MUSE_HDR_PT_DAR_CORR_C);
912 
913  /* need to recompute the (spatial) pixel table limits, but save the *
914  * previous contents, if that was not already done */
915  if (!cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_PREDAR_XLO)) {
916  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XLO,
917  cpl_propertylist_get_float(aPixtable->header,
918  MUSE_HDR_PT_XLO));
919  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XHI,
920  cpl_propertylist_get_float(aPixtable->header,
921  MUSE_HDR_PT_XHI));
922  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YLO,
923  cpl_propertylist_get_float(aPixtable->header,
924  MUSE_HDR_PT_YLO));
925  cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YHI,
926  cpl_propertylist_get_float(aPixtable->header,
927  MUSE_HDR_PT_YHI));
928  }
929  muse_pixtable_compute_limits(aPixtable);
930  } /* if aCorrect */
931  cpl_polynomial_delete(px);
932  cpl_polynomial_delete(py);
933  /* reset cosmic ray statuses in aPixtable, since the CR rejection *
934  * done here might not be appropriate for the final datacube */
935  muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
936 
937  return CPL_ERROR_NONE;
938 } /* muse_dar_check() */
939 
#define MUSE_PIXTABLE_DQ
Definition: muse_pixtable.h:49
#define MUSE_PIXTABLE_XPOS
Definition: muse_pixtable.h:51
muse_pixtable_wcs
State of the astrometric calibration of a MUSE pixel table.
Structure definition of a MUSE datacube.
Definition: muse_datacube.h:48
#define MUSE_HDR_PT_PREDAR_YLO
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
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:1825
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
Definition: muse_image.c:85
double muse_pfits_get_cd(const cpl_propertylist *aHeaders, unsigned int aAxisI, unsigned int aAxisJ)
find out the WCS coordinate at the reference point
Definition: muse_pfits.c:446
#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:1024
double muse_pfits_get_crval(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS coordinate at the reference point
Definition: muse_pfits.c:423
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_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:139
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
muse_image * muse_combine_median_create(muse_imagelist *aImages)
Median combine a list of input images.
Definition: muse_combine.c:316
muse_datacube * muse_resampling_cube(muse_pixtable *aPixtable, muse_resampling_params *aParams, muse_pixgrid **aGrid)
Resample a pixel table onto a regular grid structure representing a FITS NAXIS=3 datacube.
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
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:472
The pixel grid.
Definition: muse_pixgrid.h:65
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
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
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. ...
#define MUSE_PIXTABLE_DATA
Definition: muse_pixtable.h:48
cpl_image * dq
the data quality extension
Definition: muse_image.h:56
double muse_pfits_get_crpix(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS reference point
Definition: muse_pfits.c:401
void muse_pixgrid_delete(muse_pixgrid *aGrid)
Delete a pixel grid and remove its memory.
Definition: muse_pixgrid.c:575
#define MUSE_HDR_PT_PREDAR_YHI
Structure definition of MUSE pixel table.
static const cpl_size * muse_pixgrid_get_rows(muse_pixgrid *aGrid, cpl_size aIndex)
Return a pointer to the rows stored in one pixel.
Definition: muse_pixgrid.h:171
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:1878
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:2234
#define MUSE_HDR_PT_PREDAR_XLO
cpl_imagelist * data
the cube containing the actual data values
Definition: muse_datacube.h:76
#define MUSE_HDR_PT_PREDAR_XHI
cpl_imagelist * dq
the optional cube containing the bad pixel status
Definition: muse_datacube.h:81
cpl_propertylist * header
the FITS header
Definition: muse_datacube.h:57
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:1060
#define MUSE_PIXTABLE_STAT
Definition: muse_pixtable.h:50
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:405
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:863
muse_imagelist * muse_imagelist_new(void)
Create a new (empty) MUSE image list.
#define MUSE_PIXTABLE_LAMBDA
Definition: muse_pixtable.h:53
static cpl_size muse_pixgrid_get_index(muse_pixgrid *aGrid, 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:100
#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:1006
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:1234
double muse_phys_nrindex_owens_saturation_pressure(double temp)
Compute the saturation pressure using the Owens calibration.
Definition: muse_phys.c:54
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
#define MUSE_PIXTABLE_YPOS
Definition: muse_pixtable.h:52
#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 *aGrid, cpl_size aIndex)
Return the number of rows stored in one pixel.
Definition: muse_pixgrid.h:140
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:1352
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:1328
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:86
cpl_propertylist * header
The FITS header.
#define MUSE_HDR_PT_XHI
FITS header keyword contains the upper limit of the data in x-ion.