MUSE Pipeline Reference Manual  0.18.1
muse_flux.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  *
5  * This file is part of the MUSE Instrument Pipeline
6  * Copyright (C) 2005-2011 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_flux.h"
35 #include "muse_instrument.h"
36 
37 #include "muse_astro.h"
38 #include "muse_cplwrappers.h"
39 #include "muse_pfits.h"
40 #include "muse_quality.h"
41 #include "muse_resampling.h"
42 #include "muse_utils.h"
43 
44 /*----------------------------------------------------------------------------*/
48 /*----------------------------------------------------------------------------*/
49 
52 /*----------------------------------------------------------------------------*/
64 /*----------------------------------------------------------------------------*/
67 {
68  muse_flux_object *flux = cpl_calloc(1, sizeof(muse_flux_object));
69  return flux;
70 }
71 
72 /*----------------------------------------------------------------------------*/
82 /*----------------------------------------------------------------------------*/
83 void
85 {
86  if (!aFluxObj) {
87  return;
88  }
89  muse_datacube_delete(aFluxObj->cube);
90  aFluxObj->cube = NULL;
91  muse_image_delete(aFluxObj->intimage);
92  aFluxObj->intimage = NULL;
93  cpl_table_delete(aFluxObj->sensitivity);
94  aFluxObj->sensitivity = NULL;
95  cpl_table_delete(aFluxObj->response);
96  aFluxObj->response = NULL;
97  cpl_table_delete(aFluxObj->telluric);
98  aFluxObj->telluric = NULL;
99  cpl_table_delete(aFluxObj->tellbands);
100  aFluxObj->tellbands = NULL;
101  cpl_free(aFluxObj);
102 }
103 
104 /*----------------------------------------------------------------------------*/
112 /*----------------------------------------------------------------------------*/
113 static double
115 {
116  cpl_ensure(aTable, CPL_ERROR_NULL_INPUT, 0.);
117  cpl_table_unselect_all(aTable);
118  cpl_table_or_selected_double(aTable, "lambda", CPL_NOT_LESS_THAN,
119  kMuseNominalLambdaMin);
120  cpl_table_and_selected_double(aTable, "lambda", CPL_NOT_GREATER_THAN,
121  kMuseNominalLambdaMax);
122  cpl_size nsel = cpl_table_count_selected(aTable);
123  cpl_array *asel = cpl_table_where_selected(aTable);
124  cpl_size *sel = cpl_array_get_data_cplsize(asel);
125  double lmin = cpl_table_get_double(aTable, "lambda", sel[0], NULL),
126  lmax = cpl_table_get_double(aTable, "lambda", sel[nsel - 1], NULL);
127  cpl_array_delete(asel);
128  return (lmax - lmin) / nsel;
129 } /* muse_flux_reference_table_sampling() */
130 
131 /*----------------------------------------------------------------------------*/
160 /*----------------------------------------------------------------------------*/
161 cpl_error_code
163 {
164  cpl_ensure_code(aTable, CPL_ERROR_NULL_INPUT);
165 
166  const char *flam1 = "erg/s/cm**2/Angstrom",
167  *flam2 = "erg/s/cm^2/Angstrom";
168 
169  cpl_error_code rc = CPL_ERROR_NONE;
170  cpl_errorstate prestate = cpl_errorstate_get();
171  /* check two different types of tables: MUSE specific or HST CALSPEC */
172  if (cpl_table_has_column(aTable, "lambda") &&
173  cpl_table_has_column(aTable, "flux") &&
174  cpl_table_get_column_unit(aTable, "lambda") &&
175  cpl_table_get_column_unit(aTable, "flux") &&
176  !strncmp(cpl_table_get_column_unit(aTable, "lambda"), "Angstrom", 9) &&
177  (!strncmp(cpl_table_get_column_unit(aTable, "flux"), flam1, strlen(flam1)) ||
178  !strncmp(cpl_table_get_column_unit(aTable, "flux"), flam2, strlen(flam2)))) {
179  /* normal case: MUSE STD_FLUX_TABLE as specified; still need to *
180  * check, if we need to convert the column types (could be e.g. float) */
181  if (cpl_table_get_column_type(aTable, "lambda") != CPL_TYPE_DOUBLE) {
182  cpl_msg_debug(__func__, "Casting lambda column to double");
183  cpl_table_cast_column(aTable, "lambda", NULL, CPL_TYPE_DOUBLE);
184  }
185  if (cpl_table_get_column_type(aTable, "flux") != CPL_TYPE_DOUBLE) {
186  cpl_msg_debug(__func__, "Casting flux column to double");
187  cpl_table_cast_column(aTable, "flux", NULL, CPL_TYPE_DOUBLE);
188  }
189  /* check optional column */
190  if (cpl_table_has_column(aTable, "fluxerr")) {
191  if (cpl_table_get_column_type(aTable, "fluxerr") != CPL_TYPE_DOUBLE) {
192  cpl_msg_debug(__func__, "Casting fluxerr column to double");
193  cpl_table_cast_column(aTable, "fluxerr", NULL, CPL_TYPE_DOUBLE);
194  }
195  const char *unit = cpl_table_get_column_unit(aTable, "fluxerr");
196  if (!unit || (strncmp(unit, flam1, strlen(flam1)) &&
197  strncmp(unit, flam2, strlen(flam2)))) {
198  cpl_msg_debug(__func__, "Erasing fluxerr column because of unexpected "
199  "unit (%s)", unit);
200  cpl_table_erase_column(aTable, "fluxerr"); /* wrong unit, erase */
201  }
202  } /* if has fluxerr */
203  cpl_msg_info(__func__, "Found MUSE format, average sampling %.3f Angstrom/bin"
204  " over MUSE range", muse_flux_reference_table_sampling(aTable));
205  } else if (cpl_table_has_column(aTable, "WAVELENGTH") &&
206  cpl_table_has_column(aTable, "FLUX") &&
207  cpl_table_get_column_unit(aTable, "WAVELENGTH") &&
208  cpl_table_get_column_unit(aTable, "FLUX") &&
209  !strncmp(cpl_table_get_column_unit(aTable, "WAVELENGTH"), "ANGSTROMS", 10) &&
210  !strncmp(cpl_table_get_column_unit(aTable, "FLUX"), "FLAM", 5)) {
211 #if 0
212  printf("input HST CALSPEC table:\n");
213  cpl_table_dump_structure(aTable, stdout);
214  cpl_table_dump(aTable, cpl_table_get_nrow(aTable)/2, 3, stdout);
215  fflush(stdout);
216 #endif
217  /* other allowed case: HST CALSPEC format */
218  cpl_table_cast_column(aTable, "WAVELENGTH", "lambda", CPL_TYPE_DOUBLE);
219  cpl_table_cast_column(aTable, "FLUX", "flux", CPL_TYPE_DOUBLE);
220  cpl_table_erase_column(aTable, "WAVELENGTH");
221  cpl_table_erase_column(aTable, "FLUX");
222  cpl_table_set_column_unit(aTable, "lambda", "Angstrom");
223  cpl_table_set_column_unit(aTable, "flux", flam1);
224  /* if the table comes with the typical STATERROR/SYSERROR separation, *
225  * convert them into a single combined fluxerr column */
226  if (cpl_table_has_column(aTable, "STATERROR") &&
227  cpl_table_has_column(aTable, "SYSERROR") &&
228  cpl_table_get_column_unit(aTable, "STATERROR") &&
229  cpl_table_get_column_unit(aTable, "SYSERROR") &&
230  !strncmp(cpl_table_get_column_unit(aTable, "STATERROR"), "FLAM", 5) &&
231  !strncmp(cpl_table_get_column_unit(aTable, "SYSERROR"), "FLAM", 5)) {
232  /* Cast to double before, not to lose precision, then compute *
233  * fluxerr = sqrt(STATERROR**2 + SYSERROR**2) */
234  cpl_table_cast_column(aTable, "STATERROR", "fluxerr", CPL_TYPE_DOUBLE);
235  cpl_table_erase_column(aTable, "STATERROR");
236  cpl_table_cast_column(aTable, "SYSERROR", NULL, CPL_TYPE_DOUBLE);
237  cpl_table_power_column(aTable, "fluxerr", 2);
238  cpl_table_power_column(aTable, "SYSERROR", 2);
239  cpl_table_add_columns(aTable, "fluxerr", "SYSERROR");
240  cpl_table_erase_column(aTable, "SYSERROR");
241  cpl_table_power_column(aTable, "fluxerr", 0.5);
242  cpl_table_set_column_unit(aTable, "fluxerr", flam1);
243  } /* if error columns */
244  /* XXX how to handle invalid entries in the STATERROR column *
245  * in telluric regions (e.g. in gd105_005.fits)? */
246  /* XXX how to handle DATAQUAL column? */
247 
248  /* erase further columns we don't need */
249  if (cpl_table_has_column(aTable, "FWHM")) {
250  cpl_table_erase_column(aTable, "FWHM");
251  }
252  if (cpl_table_has_column(aTable, "DATAQUAL")) {
253  cpl_table_erase_column(aTable, "DATAQUAL");
254  }
255  if (cpl_table_has_column(aTable, "TOTEXP")) {
256  cpl_table_erase_column(aTable, "TOTEXP");
257  }
258 #if 0
259  printf("converted HST CALSPEC table:\n");
260  cpl_table_dump_structure(aTable, stdout);
261  cpl_table_dump(aTable, cpl_table_get_nrow(aTable)/2, 3, stdout);
262  fflush(stdout);
263 #endif
264  cpl_msg_info(__func__, "Found HST CALSPEC format on input, converted to "
265  "MUSE format; average sampling %.3f Angstrom/bin over MUSE "
266  "range.", muse_flux_reference_table_sampling(aTable));
267  } else {
268  cpl_msg_error(__func__, "Unknown format found!");
269 #if 0
270  cpl_table_dump_structure(aTable, stdout);
271 #endif
272  rc = CPL_ERROR_INCOMPATIBLE_INPUT;
273  } /* else: no recognized format */
274 
275  /* check for errors in the above (casting!) before returning */
276  if (!cpl_errorstate_is_equal(prestate)) {
277  rc = cpl_error_get_code();
278  }
279  return rc;
280 } /* muse_flux_reference_table_check() */
281 
282 /*----------------------------------------------------------------------------*/
320 /*----------------------------------------------------------------------------*/
321 double
322 muse_flux_response_interpolate(const cpl_table *aResponse, double aLambda,
323  double *aError, muse_flux_interpolation_type aType)
324 {
325  double rv = 0.;
326  if (aType == MUSE_FLUX_TELLURIC) {
327  rv = 1.;
328  }
329  cpl_ensure(aResponse, CPL_ERROR_NULL_INPUT, rv);
330  int size = cpl_table_get_nrow(aResponse);
331  cpl_ensure(size > 0, cpl_error_get_code(), rv);
332 
333  /* access the correct table column(s) depending on the type */
334  const double *lbda = cpl_table_get_data_double_const(aResponse, "lambda"),
335  *resp = NULL, *rerr = NULL;
336  switch (aType) {
338  resp = cpl_table_get_data_double_const(aResponse, "throughput");
339  break;
340  case MUSE_FLUX_RESP_FLUX:
341  resp = cpl_table_get_data_double_const(aResponse, "response");
342  if (aError) {
343  rerr = cpl_table_get_data_double_const(aResponse, "resperr");
344  }
345  break;
347  resp = cpl_table_get_data_double_const(aResponse, "flux");
348  if (aError) {
349  rerr = cpl_table_get_data_double_const(aResponse, "fluxerr");
350  }
351  break;
353  resp = cpl_table_get_data_double_const(aResponse, "extinction");
354  break;
355  case MUSE_FLUX_TELLURIC:
356  resp = cpl_table_get_data_double_const(aResponse, "ftelluric");
357  if (aError) {
358  rerr = cpl_table_get_data_double_const(aResponse, "ftellerr");
359  }
360  break;
361  default:
362  cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
363  return rv;
364  } /* switch aType */
365  cpl_ensure(lbda && resp, cpl_error_get_code(), rv);
366  if (aError) {
367  cpl_ensure(rerr, cpl_error_get_code(), rv);
368  }
369 
370  /* outside wavelength range of table */
371  if (aLambda < lbda[0]) {
372  return rv;
373  }
374  if (aLambda > lbda[size-1]) {
375  return rv;
376  }
377 
378  /* binary search for the correct wavelength */
379  double response = rv, resperror = 0.;
380  int l = 0, r = size - 1, /* left and right end */
381  m = (l + r) / 2; /* middle index */
382  while (CPL_TRUE) {
383  if (aLambda >= lbda[m] && aLambda <= lbda[m+1]) {
384  /* found right interval, so interpolate */
385  double lquot = (aLambda - lbda[m]) / (lbda[m+1] - lbda[m]);
386  response = resp[m] + (resp[m+1] - resp[m]) * lquot;
387  if (rerr) { /* missing error information should be non-fatal */
388  /* checked again that the derivatives leading to this error estimate *
389  * are correct; apparently it's normal then, that the resulting *
390  * errors can be smaller than the two separate input errors */
391  resperror = sqrt(pow(rerr[m] * (1 - lquot), 2.)
392  + pow(rerr[m+1] * lquot, 2.));
393  }
394 #if 0
395  cpl_msg_debug(__func__, "Found at m=%d (%f: %f+/-%f) / "
396  "m+1=%d (%f: %f+/-%f) -> %f: %f+/-%f",
397  m, lbda[m], resp[m], rerr ? rerr[m] : 0.,
398  m+1, lbda[m+1], resp[m+1], rerr ? rerr[m+1] : 0.,
399  aLambda, response, resperror);
400 #endif
401  break;
402  }
403  /* create next interval */
404  if (aLambda < lbda[m]) {
405  r = m;
406  }
407  if (aLambda > lbda[m]) {
408  l = m;
409  }
410  m = (l + r) / 2;
411  } /* while */
412 
413 #if 0
414  cpl_msg_debug(__func__, "Response %g+/-%g at lambda=%fA", response, resperror,
415  aLambda);
416 #endif
417  if (aError && rerr) {
418  *aError = resperror;
419  }
420  return response;
421 } /* muse_flux_response_interpolate() */
422 
423 /*----------------------------------------------------------------------------*/
438 /*----------------------------------------------------------------------------*/
439 static double
440 muse_flux_image_sky(cpl_image *aImage, double aX, double aY, double aHalfSize,
441  unsigned int aDSky, float *aError)
442 {
443  if (aError) {
444  *aError = FLT_MAX;
445  }
446  /* coordinates of inner window */
447  int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
448  y1 = aY - aHalfSize, y2 = aY + aHalfSize;
449  unsigned char nskyarea = 0;
450  double skylevel = 0., skyerror = 0.;
451  /* left */
452  cpl_errorstate state = cpl_errorstate_get();
453  cpl_stats_mode mode = CPL_STATS_MEDIAN | CPL_STATS_MEDIAN_DEV;
454  cpl_stats *s = cpl_stats_new_from_image_window(aImage, mode,
455  x1 - aDSky, y1, x1 - 1, y2);
456  if (s) {
457  /* only if there was no error, the area was inside the image *
458  * boundaries and we can use it, otherwise the value is undefined */
459  nskyarea++;
460  skylevel += cpl_stats_get_median(s);
461  skyerror += pow(cpl_stats_get_median_dev(s), 2);
462  cpl_stats_delete(s);
463  }
464  /* right */
465  s = cpl_stats_new_from_image_window(aImage,mode,
466  x2 + 1, y1, x2 + aDSky, y2);
467  if (s) {
468  nskyarea++;
469  skylevel += cpl_stats_get_median(s);
470  skyerror += pow(cpl_stats_get_median_dev(s), 2);
471  cpl_stats_delete(s);
472  }
473  /* bottom */
474  s = cpl_stats_new_from_image_window(aImage,mode,
475  x1, y1 - aDSky, x2, y1 - 1);
476  if (s) {
477  nskyarea++;
478  skylevel += cpl_stats_get_median(s);
479  skyerror += pow(cpl_stats_get_median_dev(s), 2);
480  cpl_stats_delete(s);
481  }
482  /* top */
483  s = cpl_stats_new_from_image_window(aImage,mode,
484  x1, y2 + 1, x2, y2 + aDSky);
485  if (s) {
486  nskyarea++;
487  skylevel += cpl_stats_get_median(s);
488  skyerror += pow(cpl_stats_get_median_dev(s), 2);
489  cpl_stats_delete(s);
490  }
491  if (nskyarea == 0) {
492  return 0.;
493  }
494  skylevel /= nskyarea;
495  skyerror = sqrt(skyerror) / nskyarea;
496  if (!cpl_errorstate_is_equal(state)) { /* reset error code */
497  cpl_errorstate_set(state);
498  }
499 #if 0
500  cpl_msg_debug(__func__, "skylevel = %f +/- %f (%u sky areas)",
501  skylevel, skyerror, nskyarea);
502 #endif
503  if (aError) {
504  *aError = skyerror;
505  }
506  return skylevel;
507 } /* muse_flux_image_sky() */
508 
509 /*----------------------------------------------------------------------------*/
526 /*----------------------------------------------------------------------------*/
527 static double
528 muse_flux_image_gaussian(cpl_image *aImage, cpl_image *aImErr, double aX,
529  double aY, double aHalfSize, unsigned int aDSky,
530  unsigned int aMaxBad, float *aFErr)
531 {
532  /* there is no simple way to count bad pixels inside an *
533  * image window, so ignore this argument at the moment */
534  UNUSED_ARGUMENT(aMaxBad);
535 
536  if (aFErr) { /* set high variance for an error case */
537  *aFErr = FLT_MAX;
538  }
539 
540  cpl_array *params = cpl_array_new(7, CPL_TYPE_DOUBLE),
541  *parerr = cpl_array_new(7, CPL_TYPE_DOUBLE);
542  /* Set some first-guess parameters to help the fitting function. *
543  * Just set background and central position, cpl_fit_image_gaussian() *
544  * finds good defaults for everything else. */
545  cpl_errorstate state = cpl_errorstate_get();
546  double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky, NULL);
547  if (!cpl_errorstate_is_equal(state)) {
548  /* if background determination fails, a default of 0 should *
549  * be good enough, so that we can ignore this failure */
550  cpl_errorstate_set(state);
551  }
552  cpl_array_set_double(params, 0, skylevel);
553  cpl_array_set_double(params, 3, aX);
554  cpl_array_set_double(params, 4, aY);
555  double rms = 0, chisq = 0;
556  /* function wants full widths but at most out to the image boundary */
557  int nx = cpl_image_get_size_x(aImage),
558  ny = cpl_image_get_size_y(aImage),
559  xsize = fmin(aHalfSize, fmin(aX - 1., nx - aX)) * 2,
560  ysize = fmin(aHalfSize, fmin(aY - 1., ny - aY)) * 2;
561  cpl_error_code rc = cpl_fit_image_gaussian(aImage, aImErr, aX, aY, xsize, ysize,
562  params, parerr, NULL, &rms, &chisq,
563  NULL, NULL, NULL, NULL, NULL);
564  if (rc != CPL_ERROR_NONE) {
565  if (rc != CPL_ERROR_ILLEGAL_INPUT) {
566  cpl_msg_debug(__func__, "rc = %d: %s", rc, cpl_error_get_message());
567  }
568  cpl_array_delete(params);
569  cpl_array_delete(parerr);
570  return 0;
571  }
572  double flux = cpl_array_get_double(params, 1, NULL),
573  ferr = cpl_array_get_double(parerr, 1, NULL);
574 #if 0 /* DEBUG */
575  double fwhmx = cpl_array_get_double(params, 5, NULL) * CPL_MATH_FWHM_SIG,
576  fwhmy = cpl_array_get_double(params, 6, NULL) * CPL_MATH_FWHM_SIG;
577  cpl_msg_debug(__func__, "%.3f,%.3f: %g+/-%g (bg: %g, FWHM: %.3f,%.3f, %g, %g)",
578  cpl_array_get_double(params, 3, NULL), cpl_array_get_double(params, 4, NULL),
579  flux, ferr, cpl_array_get_double(params, 0, NULL), fwhmx, fwhmy,
580  rms, chisq);
581 #endif
582 #if 0 /* DEBUG */
583  cpl_msg_debug(__func__, "skylevel = %f", cpl_array_get_double(params, 0, NULL));
584  cpl_msg_debug(__func__, "measured flux %f +/- %f", flux, ferr);
585 #endif
586  cpl_array_delete(params);
587  cpl_array_delete(parerr);
588  if (aFErr) {
589  *aFErr = ferr;
590  }
591  return flux;
592 } /* muse_flux_image_gaussian() */
593 
594 /*----------------------------------------------------------------------------*/
613 /*----------------------------------------------------------------------------*/
614 static double
615 muse_flux_image_moffat(cpl_image *aImage, cpl_image *aImErr, double aX,
616  double aY, double aHalfSize, unsigned int aDSky,
617  unsigned int aMaxBad, float *aFErr)
618 {
619  if (aFErr) { /* set high variance for an error case */
620  *aFErr = FLT_MAX;
621  }
622  /* extract image regions around the fiducial peak into matrix and vectors */
623  int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
624  y1 = aY - aHalfSize, y2 = aY + aHalfSize,
625  nx = cpl_image_get_size_x(aImage),
626  ny = cpl_image_get_size_y(aImage);
627  if (x1 < 1) {
628  x1 = 1;
629  }
630  if (x2 > nx) {
631  x2 = nx;
632  }
633  if (y1 < 1) {
634  y1 = 1;
635  }
636  if (y2 > ny) {
637  y2 = ny;
638  }
639  int npoints = (x2 - x1 + 1) * (y2 - y1 + 1);
640 
641  cpl_matrix *pos = cpl_matrix_new(npoints, 2);
642  cpl_vector *values = cpl_vector_new(npoints),
643  *errors = cpl_vector_new(npoints);
644  float *derr = cpl_image_get_data_float(aImErr);
645  int i, idx = 0;
646  for (i = x1; i <= x2; i++) {
647  int j;
648  for (j = y1; j <= y2; j++) {
649  int err;
650  double data = cpl_image_get(aImage, i, j, &err);
651  if (err) { /* bad pixel or error */
652  continue;
653  }
654  cpl_matrix_set(pos, idx, 0, i);
655  cpl_matrix_set(pos, idx, 1, j);
656  cpl_vector_set(values, idx, data);
657  cpl_vector_set(errors, idx, derr[(i-1) + (j-1)*nx]);
658  idx++;
659  } /* for j (y pixels) */
660  } /* for i (x pixels) */
661  /* need at least something like 4x4 pixels for a solid fit; *
662  * also check missing entries against max number of bad pixels */
663  if (idx < 16 || (npoints - idx) > (int)aMaxBad) {
664  cpl_matrix_delete(pos);
665  cpl_vector_delete(values);
666  cpl_vector_delete(errors);
667  return 0;
668  }
669  cpl_matrix_set_size(pos, idx, 2);
670  cpl_vector_set_size(values, idx);
671  cpl_vector_set_size(errors, idx);
672 
673  cpl_array *params = cpl_array_new(8, CPL_TYPE_DOUBLE),
674  *parerr = cpl_array_new(8, CPL_TYPE_DOUBLE);
675  /* Set some first-guess parameters to help the fitting function. */
676  cpl_errorstate state = cpl_errorstate_get();
677  double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky, NULL);
678  if (!cpl_errorstate_is_equal(state)) {
679  /* if background determination fails, a default of 0 should *
680  * be good enough, so that we can ignore this failure */
681  cpl_errorstate_set(state);
682  }
683  cpl_array_set_double(params, 0, skylevel);
684  cpl_array_set_double(params, 2, aX);
685  cpl_array_set_double(params, 3, aY);
686  double rms = 0, chisq = 0;
687  cpl_error_code rc = muse_utils_fit_moffat_2d(pos, values, errors,
688  params, parerr, NULL,
689  &rms, &chisq);
690  cpl_matrix_delete(pos);
691  cpl_vector_delete(values);
692  cpl_vector_delete(errors);
693  if (rc != CPL_ERROR_NONE) {
694  if (rc != CPL_ERROR_ILLEGAL_INPUT) {
695  cpl_msg_debug(__func__, "rc = %d: %s", rc, cpl_error_get_message());
696  }
697  cpl_array_delete(params);
698  cpl_array_delete(parerr);
699  return 0;
700  }
701 
702  double flux = cpl_array_get_double(params, 1, NULL);
703  if (aFErr) {
704  *aFErr = cpl_array_get_double(parerr, 1, NULL);
705  }
706 #if 0 /* DEBUG */
707  cpl_msg_debug(__func__, "skylevel = %f", cpl_array_get_double(params, 0, NULL));
708  cpl_msg_debug(__func__, "measured flux %f +/- %f", flux, cpl_array_get_double(parerr, 1, NULL));
709 #endif
710  cpl_array_delete(params);
711  cpl_array_delete(parerr);
712  return flux;
713 } /* muse_flux_image_moffat() */
714 
715 /*----------------------------------------------------------------------------*/
734 /*----------------------------------------------------------------------------*/
735 static double
736 muse_flux_image_square(cpl_image *aImage, cpl_image *aImErr, double aX,
737  double aY, double aHalfSize, unsigned int aDSky,
738  unsigned int aMaxBad, float *aFErr)
739 {
740  if (aFErr) { /* set high variance for an error case */
741  *aFErr = FLT_MAX;
742  }
743  int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
744  y1 = aY - aHalfSize, y2 = aY + aHalfSize,
745  nx = cpl_image_get_size_x(aImage),
746  ny = cpl_image_get_size_y(aImage);
747  if (x1 < 1) {
748  x1 = 1;
749  }
750  if (x2 > nx) {
751  x2 = nx;
752  }
753  if (y1 < 1) {
754  y1 = 1;
755  }
756  if (y2 > ny) {
757  y2 = ny;
758  }
759  float skyerror;
760  cpl_errorstate state = cpl_errorstate_get();
761  double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky,
762  &skyerror);
763  if (!cpl_errorstate_is_equal(state)) {
764  /* background determination is critical for this method, *
765  * reset the error but return with zero anyway */
766  cpl_errorstate_set(state);
767  cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
768  return 0.; /* fail on missing background level */
769  }
770 
771  /* extract the measurement region; fail on too many bad pixels, *
772  * but interpolate, if there are only a few bad ones */
773  cpl_image *region = cpl_image_extract(aImage, x1, y1, x2, y2);
774 #if 0 /* DEBUG */
775  cpl_msg_debug(__func__, "region [%d:%d,%d:%d] %"CPL_SIZE_FORMAT" bad pixels",
776  x1, y1, x2, y2, cpl_image_count_rejected(region));
777 #endif
778  if (cpl_image_count_rejected(region) > aMaxBad) {
779  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
780  cpl_image_delete(region);
781  return 0.; /* fail on too many bad pixels */
782  }
783  cpl_image *regerr = cpl_image_extract(aImErr, x1, y1, x2, y2);
784  if (cpl_image_count_rejected(region) > 0) {
785  cpl_detector_interpolate_rejected(region);
786  cpl_detector_interpolate_rejected(regerr);
787  }
788 
789  /* integrated flux, subtracted by the sky over the size of the *
790  * aperture; an error modeled approximately after IRAF phot */
791  int npoints = (x2 - x1 + 1) * (y2 - y1 + 1),
792  /* number of sky pixels should be counted in muse_flux_image_sky(), *
793  * but it's really not so important to add another parameter, the *
794  * error that we compute here is probably too small to be useful... */
795  nsky = 2 * aDSky * (x2 - x1 + y2 - y1 + 2);
796  double flux = cpl_image_get_flux(region) - skylevel * npoints,
797  ferr = sqrt(cpl_image_get_sqflux(regerr)
798  + npoints * skyerror*skyerror * (1. + (double)npoints / nsky));
799 #if 0 /* DEBUG */
800  cpl_msg_debug(__func__, "measured flux %f +/- %f (%d object pixels, %d pixels"
801  " with %f with sky %f)", flux, ferr, npoints, nsky,
802  cpl_image_get_flux(region), cpl_image_get_sqflux(regerr));
803 #endif
804  cpl_image_delete(region);
805  cpl_image_delete(regerr);
806  if (aFErr) {
807  *aFErr = ferr;
808  }
809  return flux;
810 } /* muse_flux_image_square() */
811 
812 /*----------------------------------------------------------------------------*/
832 /*----------------------------------------------------------------------------*/
833 static double
834 muse_flux_image_circle(cpl_image *aImage, cpl_image *aImErr, double aX,
835  double aY, double aAper, double aAnnu, double aDAnnu,
836  unsigned int aMaxBad, float *aFErr)
837 {
838  if (aFErr) { /* set high variance, for error cases */
839  *aFErr = FLT_MAX;
840  }
841  double rmax = ceil(fmax(aAper, aAnnu + aDAnnu));
842  int x1 = aX - rmax, x2 = aX + rmax,
843  y1 = aY - rmax, y2 = aY + rmax,
844  nx = cpl_image_get_size_x(aImage),
845  ny = cpl_image_get_size_y(aImage);
846  if (x1 < 1) {
847  x1 = 1;
848  }
849  if (x2 > nx) {
850  x2 = nx;
851  }
852  if (y1 < 1) {
853  y1 = 1;
854  }
855  if (y2 > ny) {
856  y2 = ny;
857  }
858  /* first loop to collect the background and *
859  * count bad pixels inside the aperture */
860  cpl_vector *vbg = cpl_vector_new((x2 - x1 + 1) * (y2 - y1 + 1)),
861  *vbe = cpl_vector_new((x2 - x1 + 1) * (y2 - y1 + 1));
862  unsigned int nbad = 0, nbg = 0;
863  int i;
864  for (i = x1; i <= x2; i++) {
865  int j;
866  for (j = y1; j <= y2; j++) {
867  double r = sqrt(pow(aX - i, 2) + pow(aY - j, 2));
868  if (r <= aAper) {
869  nbad += cpl_image_is_rejected(aImage, i, j) == 1;
870  }
871  if (r < aAnnu || r > aAnnu + aDAnnu) {
872  continue;
873  }
874  int err;
875  double value = cpl_image_get(aImage, i, j, &err);
876  if (err) { /* exclude bad pixels */
877  continue;
878  }
879  cpl_vector_set(vbg, nbg, value);
880  cpl_vector_set(vbe, nbg, cpl_image_get(aImErr, i, j, &err));
881  nbg++;
882  } /* for j (vertical pixels) */
883  } /* for i (horizontal pixels) */
884  if (nbg <= 0) {
885  cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
886  cpl_vector_delete(vbg);
887  cpl_vector_delete(vbe);
888  return 0.; /* fail on missing background pixels */
889  }
890  cpl_vector_set_size(vbg, nbg);
891  cpl_vector_set_size(vbe, nbg);
892  cpl_matrix *pos = cpl_matrix_new(1, nbg); /* we don't care about positions... */
893  double mse;
894  cpl_polynomial *fit = muse_utils_iterate_fit_polynomial(pos, vbg, vbe, NULL,
895  0, 3., &mse, NULL);
896 #if 0 /* DEBUG */
897  unsigned int nrej = nbg - cpl_vector_get_size(vbg);
898 #endif
899  nbg = cpl_vector_get_size(vbg);
900  cpl_size pows = 0; /* get the zero-order coefficient */
901  double smean = cpl_polynomial_get_coeff(fit, &pows),
902  sstdev = sqrt(mse);
903  cpl_polynomial_delete(fit);
904  cpl_matrix_delete(pos);
905  cpl_vector_delete(vbg);
906  cpl_vector_delete(vbe);
907 #if 0 /* DEBUG */
908  cpl_msg_debug(__func__, "sky: %d pixels (%d rejected), %f +/- %f; found %d "
909  "bad pixels inside aperture", nbg, nrej, smean, sstdev, nbad);
910 #endif
911  if (nbad > aMaxBad) { /* too many bad pixels inside integration area? */
912  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
913  return 0.; /* fail on too many bad pixels */
914  }
915 
916  /* now replace the few bad pixels by interpolation */
917  if (nbad > 0) {
918  cpl_detector_interpolate_rejected(aImage);
919  cpl_detector_interpolate_rejected(aImErr);
920  }
921 
922  /* second loop to integrate the flux */
923  double flux = 0.,
924  ferr = 0.;
925  unsigned int nobj = 0;
926  for (i = x1; i <= x2; i++) {
927  int j;
928  for (j = y1; j <= y2; j++) {
929  double r = sqrt(pow(aX - i, 2) + pow(aY - j, 2));
930  if (r > aAper) {
931  continue;
932  }
933  int err;
934  double value = cpl_image_get(aImage, i, j, &err),
935  error = cpl_image_get(aImErr, i, j, &err);
936  flux += value;
937  ferr += error*error;
938  nobj++;
939  } /* for j (vertical pixels) */
940  } /* for i (horizontal pixels) */
941  flux -= smean * nobj;
942  /* Compute error like IRAF phot: *
943  * error = sqrt (flux / epadu + area * stdev**2 + *
944  * area**2 * stdev**2 / nsky) *
945  * We take our summed error instead of the error computed via the gain */
946  ferr = sqrt(ferr + nobj * sstdev*sstdev * (1. + (double)nobj / nbg));
947 #if 0 /* DEBUG */
948  cpl_msg_debug(__func__, "flux: %d pixels (%d interpolated), %f +/- %f",
949  nobj, nbad, flux, ferr);
950 #endif
951  if (aFErr) {
952  *aFErr = ferr;
953  }
954  return flux;
955 } /* muse_flux_image_circle() */
956 
957 /*----------------------------------------------------------------------------*/
1004 /*----------------------------------------------------------------------------*/
1005 cpl_error_code
1007  muse_flux_object *aFluxObj)
1008 {
1009  cpl_ensure_code(aPixtable && aFluxObj, CPL_ERROR_NULL_INPUT);
1010  switch (aProfile) {
1012  cpl_msg_info(__func__, "Gaussian profile fits for flux integration");
1013  break;
1015  cpl_msg_info(__func__, "Moffat profile fits for flux integration");
1016  break;
1018  cpl_msg_info(__func__, "Circular flux integration");
1019  break;
1021  cpl_msg_info(__func__, "Simple square window flux integration");
1022  break;
1023  default:
1024  return cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
1025  }
1026 
1027  if (getenv("MUSE_DEBUG_FLUX") && atoi(getenv("MUSE_DEBUG_FLUX")) > 2) {
1028  const char *fn = "flux__pixtable.fits";
1029  cpl_msg_info(__func__, "Saving pixel table as \"%s\"", fn);
1030  muse_pixtable_save(aPixtable, fn);
1031  }
1032  muse_resampling_params *params =
1034  params->pfx = 1.; /* large pixfrac to be sure to cover most gaps */
1035  params->pfy = 1.;
1036  params->pfl = 1.;
1037  /* resample at nominal resolution, the conversion to [1/Angstrom] *
1038  * is later done when computing the sensitivity function */
1039  params->dlambda = kMuseSpectralSamplingA;
1041  params->crsigma = 10.;
1042  muse_datacube *cube = muse_resampling_cube(aPixtable, params, NULL);
1043  if (cube) {
1044  aFluxObj->cube = cube;
1045  }
1047  if (getenv("MUSE_DEBUG_FLUX") && atoi(getenv("MUSE_DEBUG_FLUX")) >= 2) {
1048  const char *fn = "flux__cube.fits";
1049  cpl_msg_info(__func__, "Saving cube as \"%s\"", fn);
1050  muse_datacube_save(aFluxObj->cube, fn);
1051  }
1052  int nplane = cpl_imagelist_get_size(cube->data) / 2; /* central plane */
1053  cpl_image *cim = cpl_imagelist_get(cube->data, nplane);
1054  /* use high sigmas for detection */
1055  double dsigmas[] = { 50., 30., 10., 8., 6., 5. };
1056  cpl_vector *vsigmas = cpl_vector_wrap(sizeof(dsigmas) / sizeof(double),
1057  dsigmas);
1058  cpl_size isigma = -1;
1059  cpl_apertures *apertures = cpl_apertures_extract(cim, vsigmas, &isigma);
1060  int napertures = apertures ? cpl_apertures_get_size(apertures) : 0;
1061  if (napertures < 1) {
1062  cpl_vector_unwrap(vsigmas);
1063  cpl_apertures_delete(apertures);
1064  return cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
1065  }
1066  cpl_msg_debug(__func__, "The %.1f sigma threshold was used to find %d source%s",
1067  cpl_vector_get(vsigmas, isigma), napertures, napertures == 1 ? "" : "s");
1068  cpl_vector_unwrap(vsigmas);
1069 #if 1 /* DEBUG */
1070  cpl_apertures_dump(apertures, stdout);
1071  fflush(stdout);
1072 #endif
1073 
1074  /* construct image of number of wavelengths x number of stars */
1075  int nlambda = cpl_imagelist_get_size(cube->data);
1076  muse_image *intimage = muse_image_new();
1077  intimage->data = cpl_image_new(nlambda, napertures, CPL_TYPE_FLOAT);
1078  intimage->stat = cpl_image_new(nlambda, napertures, CPL_TYPE_FLOAT);
1079  /* copy wavelength WCS from 3rd axis of cube to x-axis of image */
1080  intimage->header = cpl_propertylist_new();
1081  cpl_propertylist_update_double(intimage->header, "CRVAL1",
1082  cpl_propertylist_get_double(cube->header,
1083  "CRVAL3"));
1084  cpl_propertylist_update_double(intimage->header, "CRPIX1",
1085  cpl_propertylist_get_double(cube->header,
1086  "CRPIX3"));
1087  cpl_propertylist_update_double(intimage->header, "CD1_1",
1088  cpl_propertylist_get_double(cube->header,
1089  "CD3_3"));
1090  cpl_propertylist_update_string(intimage->header, "CTYPE1",
1091  cpl_propertylist_get_string(cube->header,
1092  "CTYPE3"));
1093  cpl_propertylist_update_string(intimage->header, "CUNIT1",
1094  cpl_propertylist_get_string(cube->header,
1095  "CUNIT3"));
1096  /* fill the 2nd axis with standards */
1097  cpl_propertylist_update_double(intimage->header, "CRVAL2", 1.);
1098  cpl_propertylist_update_double(intimage->header, "CRPIX2", 1.);
1099  cpl_propertylist_update_double(intimage->header, "CD2_2", 1.);
1100  cpl_propertylist_update_string(intimage->header, "CTYPE2", "PIXEL");
1101  cpl_propertylist_update_string(intimage->header, "CUNIT2", "pixel");
1102  cpl_propertylist_update_double(intimage->header, "CD1_2", 0.);
1103  cpl_propertylist_update_double(intimage->header, "CD2_1", 0.);
1104  /* we need the data units, exposure time, and instrument mode, too */
1105  cpl_propertylist_append_string(intimage->header, "BUNIT",
1106  cpl_propertylist_get_string(cube->header,
1107  "BUNIT"));
1108  cpl_propertylist_update_double(intimage->header, "EXPTIME",
1110  cpl_propertylist_update_string(intimage->header, "ESO INS MODE",
1111  cpl_propertylist_get_string(cube->header,
1112  "ESO INS MODE"));
1113 
1114  /* get DIMM seeing from the headers, convert it to size in pixels */
1115  cpl_errorstate ps = cpl_errorstate_get();
1116  double fwhm = (muse_pfits_get_fwhm_start(cube->header)
1117  + muse_pfits_get_fwhm_end(cube->header)) / 2.;
1118  if (muse_pfits_get_mode(cube->header) < MUSE_MODE_NFM_AO_N) {
1119  fwhm /= (kMuseSpaxelSizeX_WFM + kMuseSpaxelSizeY_WFM) / 2.;
1120  } else { /* for NFM */
1121  fwhm /= (kMuseSpaxelSizeX_NFM + kMuseSpaxelSizeY_NFM) / 2.;
1122  }
1123  if (!cpl_errorstate_is_equal(ps)) { /* some headers are missing */
1124  double xc = cpl_apertures_get_centroid_x(apertures, 1),
1125  yc = cpl_apertures_get_centroid_y(apertures, 1),
1126  xfwhm, yfwhm;
1127  cpl_image_get_fwhm(cim, lround(xc), lround(yc), &xfwhm, &yfwhm);
1128  if (xfwhm > 0. && yfwhm > 0.) {
1129  fwhm = (xfwhm + yfwhm) / 2.;
1130  } else if (xfwhm > 0.) {
1131  fwhm = xfwhm;
1132  } else if (yfwhm > 0.) {
1133  fwhm = yfwhm;
1134  } else {
1135  fwhm = 5.; /* total failure to measure it, assume 1 arcsec seeing */
1136  }
1137  cpl_errorstate_set(ps);
1138  cpl_msg_warning(__func__, "Using roughly estimated FWHM (%.3f pix) instead "
1139  "of DIMM seeing", fwhm);
1140  } else {
1141  cpl_msg_debug(__func__, "Using DIMM seeing of %.3f pix", fwhm);
1142  }
1143 
1144  float *data = cpl_image_get_data_float(intimage->data),
1145  *stat = cpl_image_get_data_float(intimage->stat);
1146  /* loop over all wavelengths and measure the flux */
1147  int l, ngood = 0, nillegal = 0, nbadbg = 0; /* count good fits and errors */
1148  #pragma omp parallel for default(none) /* as req. by Ralf */ \
1149  shared(apertures, aProfile, cube, data, fwhm, napertures, nbadbg, \
1150  ngood, nillegal, nlambda, stat)
1151  for (l = 0; l < nlambda; l++) {
1152  cpl_image *plane = cpl_imagelist_get(cube->data, l),
1153  *pldq = cpl_imagelist_get(cube->dq, l),
1154  *plerr = cpl_image_duplicate(cpl_imagelist_get(cube->stat, l));
1155 #if 0 /* DEBUG */
1156  cpl_stats *stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1157  cpl_msg_debug(__func__, "lambda = %d/%f %s", l + 1,
1158  (l + 1 - cpl_propertylist_get_double(cube->header, "CRPIX3"))
1159  * cpl_propertylist_get_double(cube->header, "CD3_3")
1160  + cpl_propertylist_get_double(cube->header, "CRVAL3"),
1161  cpl_propertylist_get_string(cube->header, "CUNIT3"));
1162  cpl_msg_debug(__func__, "variance: %g...%g...%g", cpl_stats_get_min(stats),
1163  cpl_stats_get_mean(stats), cpl_stats_get_max(stats));
1164  cpl_stats_delete(stats);
1165 #endif
1166  /* make sure to exclude bad pixels from the fits below */
1167  muse_quality_image_reject_using_dq(plane, pldq, plerr);
1168 #if 0 /* DEBUG */
1169  stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1170  cpl_msg_debug(__func__, "cut variance: %g...%g...%g (%"CPL_SIZE_FORMAT" bad"
1171  " pixel)", cpl_stats_get_min(stats), cpl_stats_get_mean(stats),
1172  cpl_stats_get_max(stats), cpl_image_count_rejected(plane));
1173  cpl_stats_delete(stats);
1174 #endif
1175  /* convert variable to sigmas */
1176  cpl_image_power(plerr, 0.5);
1177 #if 0 /* DEBUG */
1178  stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1179  cpl_msg_debug(__func__, "errors: %g...%g...%g", cpl_stats_get_min(stats),
1180  cpl_stats_get_mean(stats), cpl_stats_get_max(stats));
1181  cpl_stats_delete(stats);
1182 #endif
1183  cpl_errorstate state = cpl_errorstate_get();
1184  int n;
1185  for (n = 1; n <= napertures; n++) {
1186  /* use detection aperture to construct much larger one for measurement */
1187  double xc = cpl_apertures_get_centroid_x(apertures, n),
1188  yc = cpl_apertures_get_centroid_y(apertures, n),
1189  size = sqrt(cpl_apertures_get_npix(apertures, n)),
1190  xfwhm, yfwhm;
1191  cpl_errorstate prestate = cpl_errorstate_get();
1192  cpl_image_get_fwhm(plane, lround(xc), lround(yc), &xfwhm, &yfwhm);
1193  if (xfwhm < 0 || yfwhm < 0) {
1194  data[l + (n-1) * nlambda] = 0.;
1195  stat[l + (n-1) * nlambda] = FLT_MAX;
1196  cpl_errorstate_set(prestate);
1197  continue;
1198  }
1199  /* half size for flux integration, at least 3 x FWHM */
1200  double halfsize = 1.5 * (xfwhm + yfwhm);
1201  if (halfsize < size / 2) { /* at least the size of the det. aperture */
1202  halfsize = size / 2;
1203  }
1204 #if 0 /* DEBUG */
1205  cpl_msg_debug(__func__, "%.2f,%.2f FWHM %.2f %.2f size %.2f --> %.2f",
1206  xc, yc, xfwhm, yfwhm, size, halfsize * 2.);
1207 #endif
1208 
1209  switch (aProfile) {
1211  data[l + (n-1) * nlambda] = muse_flux_image_gaussian(plane, plerr, xc, yc,
1212  halfsize, 5, 10,
1213  &stat[l + (n-1) * nlambda]);
1214  break;
1216  data[l + (n-1) * nlambda] = muse_flux_image_moffat(plane, plerr, xc, yc,
1217  halfsize, 5, 10,
1218  &stat[l + (n-1) * nlambda]);
1219  break;
1220  case MUSE_FLUX_PROFILE_CIRCLE: {
1221  /* the circular method needs larger region to properly integrate *
1222  * everything at least something like 4 x FWHM to be sure */
1223  double maxfwhm = fmax((xfwhm + yfwhm) / 2., fwhm),
1224  radius = halfsize < maxfwhm ? maxfwhm : halfsize,
1225  rannu = radius * 5. / 4.; /* background annulus a bit larger */
1226  data[l + (n-1) * nlambda] = muse_flux_image_circle(plane, plerr, xc, yc,
1227  radius, rannu, 10, 10,
1228  &stat[l + (n-1) * nlambda]);
1229  break;
1230  } /* case MUSE_FLUX_PROFILE_CIRCLE */
1231  default: /* MUSE_FLUX_PROFILE_EQUAL_SQUARE */
1232  data[l + (n-1) * nlambda] = muse_flux_image_square(plane, plerr, xc, yc,
1233  halfsize, 5, 10,
1234  &stat[l + (n-1) * nlambda]);
1235  } /* switch */
1236  if (data[l + (n-1) * nlambda] < 0 || !isfinite(data[l + (n-1) * nlambda])) {
1237  data[l + (n-1) * nlambda] = 0.;
1238  stat[l + (n-1) * nlambda] = FLT_MAX;
1239  }
1240  } /* for n (all apertures) */
1241 
1242  /* count "Illegal input" errors and for those reset the state, as there *
1243  * can be many of them (one for each incompletely filled wavelength plane */
1244  if (!cpl_errorstate_is_equal(state)) {
1245  if (cpl_error_get_code() == CPL_ERROR_ILLEGAL_INPUT) {
1246  cpl_errorstate_set(state);
1247  #pragma omp atomic
1248  nillegal++;
1249  } else if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
1250  cpl_errorstate_set(state);
1251  #pragma omp atomic
1252  nbadbg++;
1253  }
1254  } else {
1255  #pragma omp atomic
1256  ngood++;
1257  }
1258 
1259  cpl_image_delete(plerr);
1260  } /* for l (all wavelengths) */
1261  cpl_apertures_delete(apertures);
1262  if (getenv("MUSE_DEBUG_FLUX") && atoi(getenv("MUSE_DEBUG_FLUX")) >= 2) {
1263  const char *fn = "flux__intimage.fits";
1264  cpl_msg_info(__func__, "Saving flux image as \"%s\"", fn);
1265  intimage->dq = cpl_image_new(nlambda, napertures, CPL_TYPE_INT);
1266  muse_image_save(intimage, fn);
1267  }
1268  aFluxObj->intimage = intimage;
1269 
1270  if (nillegal > 0 || nbadbg > 0) {
1271  cpl_msg_warning(__func__, "Successful fits in %d wavelength planes, but "
1272  "encountered %d \"Illegal input\" errors and %d bad "
1273  "background determinations", ngood, nillegal, nbadbg);
1274  } else {
1275  cpl_msg_info(__func__, "Successful fits in %d wavelength planes", ngood);
1276  }
1277  return CPL_ERROR_NONE;
1278 } /* muse_flux_integrate_std() */
1279 
1280 /* prominent telluric absorption bands, together with clean regions: *
1281  * [0] lower limit of telluric region *
1282  * [1] upper limit of telluric region *
1283  * [2] lower limit of fit region (excludes telluric region) *
1284  * [3] upper limit of fit region (excludes telluric region) *
1285  * created from by-hand measurements on spectra from the ESO sky *
1286  * model, and from the Keck list of telluric lines */
1287 static const double kTelluricBands[][4] = {
1288  { 6273., 6320., 6213., 6380. }, /* now +/-60 Angstrom around... */
1289  { 6864., 6967., 6750., 7130. }, /* B-band */
1290  { 7164., 7325., 7070., 7580. },
1291  { 7590., 7700., 7470., 7830. }, /* A-band */
1292  { 8131., 8345., 7900., 8600. },
1293  { 8952., 9028., 8850., 9082. }, /* XXX very rough! */
1294  { 9274., 9770., 9080., 9263. }, /* XXX very very rough! */
1295  { -1., -1., -1., -1. }
1296 };
1297 
1298 /*----------------------------------------------------------------------------*/
1302 /*----------------------------------------------------------------------------*/
1304  { "lmin", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
1305  "lower limit of the telluric region", CPL_TRUE },
1306  { "lmax", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
1307  "upper limit of the telluric region", CPL_TRUE },
1308  { "bgmin", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
1309  "lower limit of the background region", CPL_TRUE },
1310  { "bgmax", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
1311  "upper limit of the background region", CPL_TRUE },
1312  { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
1313 };
1314 
1315 /*----------------------------------------------------------------------------*/
1322 /*----------------------------------------------------------------------------*/
1323 static void
1324 muse_flux_response_set_telluric_bands(muse_flux_object *aFluxObj,
1325  const cpl_table *aTellBands)
1326 {
1327  if (!aFluxObj) {
1328  cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
1329  return;
1330  }
1331  /* if a valid table was passed, just set that */
1332  if (aTellBands && muse_cpltable_check(aTellBands, muse_response_tellbands_def)
1333  == CPL_ERROR_NONE) {
1334  cpl_msg_debug(__func__, "using given table for telluric bands");
1335  aFluxObj->tellbands = cpl_table_duplicate(aTellBands);
1336  return;
1337  }
1338  /* create a table for the default regions */
1339  unsigned int ntell = sizeof(kTelluricBands) / sizeof(kTelluricBands[0]) - 1;
1340  cpl_msg_debug(__func__, "using builtin regions for telluric bands (%u "
1341  "entries)", ntell);
1342  aFluxObj->tellbands = muse_cpltable_new(muse_response_tellbands_def, ntell);
1343  cpl_table *tb = aFluxObj->tellbands; /* shortcut */
1344  unsigned int k;
1345  for (k = 0; k < ntell; k++) {
1346  cpl_table_set_double(tb, "lmin", k, kTelluricBands[k][0]);
1347  cpl_table_set_double(tb, "lmax", k, kTelluricBands[k][1]);
1348  cpl_table_set_double(tb, "bgmin", k, kTelluricBands[k][2]);
1349  cpl_table_set_double(tb, "bgmax", k, kTelluricBands[k][3]);
1350  } /* for k */
1351  if (getenv("MUSE_DEBUG_FLUX") && atoi(getenv("MUSE_DEBUG_FLUX")) >= 2) {
1352  const char *fn = "flux__tellregions.fits";
1353  cpl_msg_info(__func__, "Saving telluric bands table as \"%s\"", fn);
1354  cpl_table_save(tb, NULL, NULL, fn, CPL_IO_CREATE);
1355  }
1356  return;
1357 } /* muse_flux_response_set_telluric_bands() */
1358 
1359 /*----------------------------------------------------------------------------*/
1372 /*----------------------------------------------------------------------------*/
1373 static void
1374 muse_flux_response_dump_sensitivity(muse_flux_object *aFluxObj,
1375  const char *aName)
1376 {
1377  char *dodebug = getenv("MUSE_DEBUG_FLUX");
1378  if (!dodebug || (dodebug && atoi(dodebug) <= 0)) {
1379  return;
1380  }
1381  char *fn = cpl_sprintf("flux__sens_%s.ascii", aName);
1382  FILE *fp = fopen(fn, "w");
1383  fprintf(fp, "#"); /* prefix first line (table header) for easier plotting */
1384  cpl_table_dump(aFluxObj->sensitivity, 0,
1385  cpl_table_get_nrow(aFluxObj->sensitivity), fp);
1386  fclose(fp);
1387  cpl_msg_debug(__func__, "Written %"CPL_SIZE_FORMAT" datapoints to \"%s\"",
1388  cpl_table_get_nrow(aFluxObj->sensitivity), fn);
1389  cpl_free(fn);
1390 } /* muse_flux_response_dump_sensitivity() */
1391 
1392 /*----------------------------------------------------------------------------*/
1411 /*----------------------------------------------------------------------------*/
1412 static void
1413 muse_flux_response_sensitivity(muse_flux_object *aFluxObj,
1414  unsigned int aStar, const cpl_table *aReference,
1415  double aAirmass, const cpl_table *aExtinct)
1416 {
1417  double crval = cpl_propertylist_get_double(aFluxObj->intimage->header, "CRVAL1"),
1418  cdelt = cpl_propertylist_get_double(aFluxObj->intimage->header, "CD1_1"),
1419  crpix = cpl_propertylist_get_double(aFluxObj->intimage->header, "CRPIX1"),
1420  exptime = muse_pfits_get_exptime(aFluxObj->intimage->header);
1421  int nlambda = cpl_image_get_size_x(aFluxObj->intimage->data);
1422 
1423  aFluxObj->sensitivity = cpl_table_new(nlambda);
1424  cpl_table *sensitivity = aFluxObj->sensitivity;
1425  cpl_table_new_column(sensitivity, "lambda", CPL_TYPE_DOUBLE);
1426  cpl_table_new_column(sensitivity, "sens", CPL_TYPE_DOUBLE);
1427  cpl_table_new_column(sensitivity, "serr", CPL_TYPE_DOUBLE);
1428  cpl_table_new_column(sensitivity, "dq", CPL_TYPE_INT);
1429  cpl_table_set_column_format(sensitivity, "dq", "%u");
1430  float *data = cpl_image_get_data_float(aFluxObj->intimage->data),
1431  *stat = cpl_image_get_data_float(aFluxObj->intimage->stat);
1432  int l, idx = 0;
1433  for (l = 0; l < nlambda; l++) {
1434  if (data[l + aStar*nlambda] <= 0. ||
1435  stat[l + aStar*nlambda] <= 0. ||
1436  stat[l + aStar*nlambda] == FLT_MAX) { /* exclude bad fits */
1437  continue;
1438  }
1439 
1440  double lambda = crval + cdelt * (l + 1 - crpix),
1441  /* interpolate extinction curve at this wavelength */
1442  extinct = !aExtinct ? 0. /* no extinction term */
1443  : muse_flux_response_interpolate(aExtinct, lambda, NULL,
1445  cpl_errorstate prestate = cpl_errorstate_get();
1446  double referr = 0.,
1447  ref = muse_flux_response_interpolate(aReference, lambda, &referr,
1449  /* on error, try again without trying to find the fluxerr column */
1450  if (!cpl_errorstate_is_equal(prestate)) {
1451  cpl_errorstate_set(prestate); /* don't want to propagate this error outside */
1452  ref = muse_flux_response_interpolate(aReference, lambda, NULL,
1454  }
1455  /* calibration factor at this wavelength: ratio of observed count *
1456  * rate corrected for extinction to expected flux in magnitudes */
1457  double c = 2.5 * log10(data[l + aStar*nlambda]
1458  / exptime / cdelt / ref)
1459  + aAirmass * extinct,
1460  cerr = sqrt(pow(referr / ref, 2) + stat[l + aStar*nlambda]
1461  / pow(data[l + aStar*nlambda], 2))
1462  * 2.5 / CPL_MATH_LN10;
1463  cpl_table_set_double(sensitivity, "lambda", idx, lambda);
1464  cpl_table_set_double(sensitivity, "sens", idx, c);
1465  cpl_table_set_double(sensitivity, "serr", idx, cerr);
1466  cpl_table_set_int(sensitivity, "dq", idx, EURO3D_GOODPIXEL);
1467  idx++;
1468  } /* for l (all wavelengths) */
1469  /* cut the data to the used size */
1470  cpl_table_set_size(sensitivity, idx);
1471 } /* muse_flux_response_sensitivity() */
1472 
1473 /*----------------------------------------------------------------------------*/
1486 /*----------------------------------------------------------------------------*/
1487 static void
1488 muse_flux_response_mark_questionable(muse_flux_object *aFluxObj)
1489 {
1490  if (!aFluxObj) {
1491  cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
1492  return;
1493  }
1494  if (!aFluxObj->tellbands) {
1495  cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
1496  return;
1497  }
1498  cpl_table *tsens = aFluxObj->sensitivity,
1499  *tb = aFluxObj->tellbands;
1500 
1501  /* check for MUSE setup: is it nominal wavelength range? */
1502  cpl_boolean isnominal = muse_pfits_get_mode(aFluxObj->intimage->header)
1503  > MUSE_MODE_WFM_NONAO_X;
1504  /* exclude regions within the telluric absorption bands *
1505  * and outside wavelength range */
1506  int irow, nfluxes = cpl_table_get_nrow(tsens);
1507  for (irow = 0; irow < nfluxes; irow++) {
1508  double lambda = cpl_table_get_double(tsens, "lambda", irow, NULL);
1509  unsigned int dq = EURO3D_GOODPIXEL,
1510  k, nk = cpl_table_get_nrow(tb);
1511  for (k = 0; k < nk; k++) {
1512  double lmin = cpl_table_get_double(tb, "lmin", k, NULL),
1513  lmax = cpl_table_get_double(tb, "lmax", k, NULL);
1514  if (lambda >= lmin && lambda <= lmax) {
1515  dq |= EURO3D_TELLURIC;
1516  }
1517  } /* for k */
1518  if (isnominal && lambda < kMuseNominalCutoff) {
1519  dq |= EURO3D_OUTSDRANGE;
1520  }
1521  cpl_table_set_int(tsens, "dq", irow, dq);
1522  } /* for irow (all table) */
1523 } /* muse_flux_response_mark_questionable() */
1524 
1525 /*----------------------------------------------------------------------------*/
1548 /*----------------------------------------------------------------------------*/
1549 static cpl_polynomial *
1550 muse_flux_response_fit(muse_flux_object *aFluxObj,
1551  double aLambda1, double aLambda2,
1552  unsigned int aOrder, double aRSigma, double *aRMSE)
1553 {
1554  cpl_table *tsens = aFluxObj->sensitivity;
1555  cpl_table_select_all(tsens); /* default, but make sure it's true...*/
1556  cpl_table_and_selected_int(tsens, "dq", CPL_NOT_EQUAL_TO, EURO3D_GOODPIXEL);
1557  cpl_table_and_selected_int(tsens, "dq", CPL_NOT_EQUAL_TO, EURO3D_TELLCOR);
1558  cpl_table_or_selected_double(tsens, "lambda", CPL_LESS_THAN, aLambda1);
1559  cpl_table_or_selected_double(tsens, "lambda", CPL_GREATER_THAN, aLambda2);
1560  /* keep the "bad" ones around */
1561  cpl_table *tunwanted = cpl_table_extract_selected(tsens);
1562  cpl_table_erase_selected(tsens);
1563  muse_flux_response_dump_sensitivity(aFluxObj, "fitinput");
1564 
1565  /* convert sensitivity table to matrix (lambda) and vectors (sens *
1566  * and serr), exclude pixels marked as not EURO3D_GOODPIXEL */
1567  int nrow = cpl_table_get_nrow(tsens);
1568  cpl_matrix *lambdas = cpl_matrix_new(1, nrow);
1569  cpl_vector *sens = cpl_vector_new(nrow),
1570  *serr = cpl_vector_new(nrow);
1571  memcpy(cpl_matrix_get_data(lambdas),
1572  cpl_table_get_data_double_const(tsens, "lambda"), nrow*sizeof(double));
1573  memcpy(cpl_vector_get_data(sens),
1574  cpl_table_get_data_double_const(tsens, "sens"), nrow*sizeof(double));
1575  memcpy(cpl_vector_get_data(serr),
1576  cpl_table_get_data_double_const(tsens, "serr"), nrow*sizeof(double));
1577 
1578  /* do the fit */
1579  double chisq, mse;
1580  cpl_polynomial *fit = muse_utils_iterate_fit_polynomial(lambdas, sens, serr,
1581  tsens, aOrder, aRSigma,
1582  &mse, &chisq);
1583  int nout = cpl_vector_get_size(sens);
1584 #if 0
1585  cpl_msg_debug(__func__, "transferred %d entries (%.3f...%.3f) for the "
1586  "order %u fit, %d entries are left, RMS %f", nrow, aLambda1,
1587  aLambda2, aOrder, nout, sqrt(mse));
1588 #endif
1589  cpl_matrix_delete(lambdas);
1590  cpl_vector_delete(sens);
1591  cpl_vector_delete(serr);
1592  if (aRMSE) {
1593  *aRMSE = mse / (nout - aOrder - 1);
1594  }
1595 
1596  /* put "bad" entries back, at the end of the table */
1597  cpl_table_insert(tsens, tunwanted, nout);
1598  cpl_table_delete(tunwanted);
1599  return fit;
1600 } /* muse_flux_response_fit() */
1601 
1602 /*----------------------------------------------------------------------------*/
1621 /*----------------------------------------------------------------------------*/
1622 static void
1623 muse_flux_response_telluric(muse_flux_object *aFluxObj, double aAirmass)
1624 {
1625  cpl_table *tsens = aFluxObj->sensitivity,
1626  *tb = aFluxObj->tellbands;
1627  cpl_table_new_column(tsens, "sens_orig", CPL_TYPE_DOUBLE);
1628  cpl_table_new_column(tsens, "serr_orig", CPL_TYPE_DOUBLE);
1629  cpl_table_new_column(tsens, "tellcor", CPL_TYPE_DOUBLE);
1630  unsigned int k, nk = cpl_table_get_nrow(tb);
1631  for (k = 0; k < nk; k++) {
1632  unsigned int order = 2;
1633  double lmin = cpl_table_get_double(tb, "lmin", k, NULL),
1634  lmax = cpl_table_get_double(tb, "lmax", k, NULL),
1635  bgmin = cpl_table_get_double(tb, "bgmin", k, NULL),
1636  bgmax = cpl_table_get_double(tb, "bgmax", k, NULL);
1637  /* if we extrapolate (redward) then use a linear fit only */
1638  if (bgmax < lmax) {
1639  order = 1;
1640  }
1641  double rmse = 0.;
1642  cpl_polynomial *fit = muse_flux_response_fit(aFluxObj, bgmin, bgmax,
1643  order, 3., &rmse);
1644  cpl_msg_debug(__func__, "Telluric region %u: %.2f...%.2f, reference region "
1645  "%.2f...%.2f", k + 1, lmin, lmax, bgmin, bgmax);
1646 #if 0
1647  cpl_polynomial_dump(fit, stdout);
1648  fflush(stdout);
1649 #endif
1650  int irow, nrow = cpl_table_get_nrow(tsens);
1651  for (irow = 0; irow < nrow; irow++) {
1652  double lambda = cpl_table_get_double(tsens, "lambda", irow, NULL);
1653  if (lambda >= lmin && lambda <= lmax &&
1654  (unsigned int)cpl_table_get_int(tsens, "dq", irow, NULL)
1655  == EURO3D_TELLURIC) {
1656  double origval = cpl_table_get_double(tsens, "sens", irow, NULL),
1657  origerr = cpl_table_get_double(tsens, "serr", irow, NULL),
1658  interpval = cpl_polynomial_eval_1d(fit, lambda, NULL);
1659  cpl_table_set_int(tsens, "dq", irow, EURO3D_TELLCOR);
1660  cpl_table_set_double(tsens, "sens_orig", irow, origval);
1661  cpl_table_set_double(tsens, "sens", irow, interpval);
1662  /* correct the error bars of the fitted points *
1663  * by adding in quadrature the reduced MSE */
1664  cpl_table_set_double(tsens, "serr_orig", irow, origerr);
1665  cpl_table_set_double(tsens, "serr", irow, sqrt(origerr*origerr + rmse));
1666 
1667  if (interpval > origval) {
1668  /* compute the factor, as flux ratio */
1669  double ftelluric = pow(10, -0.4 * (interpval - origval));
1670  cpl_table_set(tsens, "tellcor", irow, ftelluric);
1671  } else {
1672  cpl_table_set_double(tsens, "tellcor", irow, 1.);
1673  }
1674  }
1675  } /* for irow */
1676  cpl_polynomial_delete(fit);
1677  } /* for k (telluric line regions) */
1678  /* correct for the airmass */
1679  cpl_table_power_column(tsens, "tellcor", 1. / aAirmass);
1680 } /* muse_flux_response_telluric() */
1681 
1682 /*----------------------------------------------------------------------------*/
1698 /*----------------------------------------------------------------------------*/
1699 static void
1700 muse_flux_response_extrapolate(muse_flux_object *aFluxObj)
1701 {
1702  cpl_table *tsens = aFluxObj->sensitivity;
1703  cpl_propertylist *order = cpl_propertylist_new();
1704  cpl_propertylist_append_bool(order, "lambda", CPL_FALSE);
1705  cpl_table_sort(tsens, order);
1706 
1707  int nrow = cpl_table_get_nrow(tsens);
1708  /* first and last good wavelength and corresponsing error */
1709  double lambda1 = cpl_table_get_double(tsens, "lambda", 0, NULL),
1710  serr1 = cpl_table_get_double(tsens, "serr", 0, NULL);
1711  unsigned int dq = cpl_table_get_int(tsens, "dq", 0, NULL);
1712  int irow = 0;
1713  while (dq != EURO3D_GOODPIXEL && dq != EURO3D_TELLCOR) {
1714  lambda1 = cpl_table_get_double(tsens, "lambda", ++irow, NULL);
1715  serr1 = cpl_table_get_double(tsens, "serr", irow, NULL);
1716  dq = cpl_table_get_int(tsens, "dq", irow, NULL);
1717  }
1718  double lambda2 = cpl_table_get_double(tsens, "lambda", nrow - 1, NULL),
1719  serr2 = cpl_table_get_double(tsens, "serr", nrow - 1, NULL);
1720  dq = cpl_table_get_int(tsens, "dq", nrow - 1, NULL);
1721  irow = nrow - 1;
1722  while (dq != EURO3D_GOODPIXEL && dq != EURO3D_TELLCOR) {
1723  lambda2 = cpl_table_get_double(tsens, "lambda", --irow, NULL);
1724  serr2 = cpl_table_get_double(tsens, "serr", irow, NULL);
1725  dq = cpl_table_get_int(tsens, "dq", irow, NULL);
1726  }
1727  const double ldist = 100.;
1728  cpl_polynomial *fit1 = muse_flux_response_fit(aFluxObj, lambda1, lambda1 + ldist,
1729  1, 5., NULL),
1730  *fit2 = muse_flux_response_fit(aFluxObj, lambda2 - ldist, lambda2,
1731  1, 5., NULL);
1732  nrow = cpl_table_get_nrow(tsens); /* fitting may have erased rows! */
1733  double d1 = (lambda1 - kMuseLambdaMinX) / 10., /* want 10 additional entries... */
1734  d2 = (kMuseLambdaMaxX - lambda2) / 10.; /* ... at each end */
1735  cpl_table_set_size(tsens, nrow + 20);
1736  irow = nrow;
1737  double l;
1738  for (l = kMuseLambdaMinX; l <= lambda1 - d1; l += d1) {
1739  double sens = cpl_polynomial_eval_1d(fit1, l, NULL),
1740  /* error that doubles every ldist Angstrom */
1741  serr = 2. * (lambda1 - l) / ldist * serr1 + serr1;
1742  if (sens <= 0) {
1743  cpl_table_set_invalid(tsens, "lambda", irow++);
1744  cpl_msg_debug(__func__, "invalid blueward extrapolation: %.3f %f +/- %f",
1745  l, sens, serr);
1746  continue;
1747  }
1748  cpl_table_set_double(tsens, "lambda", irow, l);
1749  cpl_table_set_double(tsens, "sens", irow, sens);
1750  cpl_table_set_double(tsens, "serr", irow, serr);
1751  cpl_table_set_int(tsens, "dq", irow++, (int)EURO3D_OUTSDRANGE);
1752  } /* for l */
1753  for (l = lambda2 + d2; l <= kMuseLambdaMaxX; l += d2) {
1754  double sens = cpl_polynomial_eval_1d(fit2, l, NULL),
1755  serr = 2. * (l - lambda2) / ldist * serr2 + serr2;
1756  if (sens <= 0) {
1757  cpl_table_set_invalid(tsens, "lambda", irow++);
1758  cpl_msg_debug(__func__, "invalid redward extrapolation: %.3f %f +/- %f",
1759  l, sens, serr);
1760  continue;
1761  }
1762  cpl_table_set_double(tsens, "lambda", irow, l);
1763  cpl_table_set_double(tsens, "sens", irow, sens);
1764  cpl_table_set_double(tsens, "serr", irow, serr);
1765  cpl_table_set_int(tsens, "dq", irow++, (int)EURO3D_OUTSDRANGE);
1766  } /* for l */
1767  cpl_msg_debug(__func__, "Extrapolated regions %.1f...%.1f Angstrom (using data "
1768  "from %.1f...%.1f A) and %.1f...%.1f Angstrom (using %.1f...%.1f A)",
1769  kMuseLambdaMinX, lambda1 - d1, lambda1, lambda1 + ldist,
1770  lambda2 + d2, kMuseLambdaMaxX, lambda2 - ldist, lambda2);
1771 #if 0
1772  cpl_polynomial_dump(fit1, stdout);
1773  cpl_polynomial_dump(fit2, stdout);
1774  fflush(stdout);
1775 #endif
1776  cpl_polynomial_delete(fit1);
1777  cpl_polynomial_delete(fit2);
1778  /* clean up invalid entries */
1779  cpl_table_select_all(tsens);
1780  cpl_table_and_selected_invalid(tsens, "sens");
1781  cpl_table_erase_selected(tsens);
1782  /* sort the resulting table again */
1783  cpl_table_sort(tsens, order);
1784  cpl_propertylist_delete(order);
1785 } /* muse_flux_response_extrapolate() */
1786 
1787 /*----------------------------------------------------------------------------*/
1831 /*----------------------------------------------------------------------------*/
1832 cpl_error_code
1834  const cpl_table *aReference,
1835  const cpl_table *aTellBands,
1836  const cpl_table *aExtinct)
1837 {
1838  cpl_ensure_code(aFluxObj && aFluxObj->intimage && aReference,
1839  CPL_ERROR_NULL_INPUT);
1840  cpl_ensure_code(aAirmass >= 1., CPL_ERROR_ILLEGAL_INPUT);
1841  if (!aExtinct) {
1842  cpl_msg_warning(__func__, "Extinction table not given!");
1843  }
1844  muse_flux_response_set_telluric_bands(aFluxObj, aTellBands);
1845 
1846  int nlambda = cpl_image_get_size_x(aFluxObj->intimage->data),
1847  nobjects = cpl_image_get_size_y(aFluxObj->intimage->data);
1848  /* XXX for now, find one star, the brightest, to use here */
1849  double flux = 0.;
1850  int n, nstar = 1;
1851  for (n = 1; n <= nobjects; n++) {
1852  double this = cpl_image_get_flux_window(aFluxObj->intimage->data, 1, n,
1853  nlambda, n);
1854  cpl_msg_debug(__func__, "flux(%d) = %e", n, this);
1855  if (this > flux) {
1856  flux = this;
1857  nstar = n;
1858  }
1859  } /* for n (all objects) */
1860  if (nobjects > 1) {
1861  cpl_msg_warning(__func__, "Assuming that the brightest star (%d of %d) is "
1862  "the reference source!", nstar, nobjects);
1863  }
1864  /* table of sensitivity function, its sigma, and a Euro3D-like quality flag */
1865  muse_flux_response_sensitivity(aFluxObj,
1866  nstar - 1, aReference,
1867  aAirmass, aExtinct);
1868  muse_flux_response_dump_sensitivity(aFluxObj, "initial");
1869 
1870  muse_flux_response_mark_questionable(aFluxObj);
1871  muse_flux_response_dump_sensitivity(aFluxObj, "intermediate");
1872  /* remove the ones outside the wavelength range */
1873  cpl_table_select_all(aFluxObj->sensitivity);
1874  cpl_table_and_selected_int(aFluxObj->sensitivity, "dq", CPL_EQUAL_TO,
1875  (int)EURO3D_OUTSDRANGE);
1876  cpl_table_erase_selected(aFluxObj->sensitivity);
1877  muse_flux_response_dump_sensitivity(aFluxObj, "intercut");
1878 
1879  /* interpolate telluric (dq == EURO3D_TELLURIC) regions *
1880  * and compute telluric correction factor */
1881  muse_flux_response_telluric(aFluxObj, aAirmass);
1882  muse_flux_response_dump_sensitivity(aFluxObj, "interpolated");
1883  /* extend the wavelength range using linear extrapolation */
1884  muse_flux_response_extrapolate(aFluxObj);
1885  muse_flux_response_dump_sensitivity(aFluxObj, "extrapolated");
1886 
1887  return CPL_ERROR_NONE;
1888 } /* muse_flux_response_compute() */
1889 
1890 /*----------------------------------------------------------------------------*/
1899 /*----------------------------------------------------------------------------*/
1901  { "lambda", CPL_TYPE_DOUBLE, "Angstrom", "%7.2f", "wavelength", CPL_TRUE },
1902  { "response", CPL_TYPE_DOUBLE,
1903  "2.5*log10((count/s/Angstrom)/(erg/s/cm**2/Angstrom))", "%.4e",
1904  "instrument response derived from standard star", CPL_TRUE },
1905  { "resperr", CPL_TYPE_DOUBLE,
1906  "2.5*log10((count/s/Angstrom)/(erg/s/cm**2/Angstrom))", "%.4e",
1907  "instrument response error derived from standard star", CPL_TRUE },
1908  { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
1909 };
1910 
1911 /*----------------------------------------------------------------------------*/
1920 /*----------------------------------------------------------------------------*/
1922  { "lambda", CPL_TYPE_DOUBLE, "Angstrom", "%7.2f", "wavelength", CPL_TRUE },
1923  { "ftelluric", CPL_TYPE_DOUBLE, "", "%.5f",
1924  "the telluric correction factor, normalized to an airmass of 1", CPL_TRUE },
1925  { "ftellerr", CPL_TYPE_DOUBLE, "", "%.5f",
1926  "the error of the telluric correction factor", CPL_TRUE },
1927  { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
1928 };
1929 
1930 /*----------------------------------------------------------------------------*/
1945 /*----------------------------------------------------------------------------*/
1946 static void
1948  cpl_table *aResp, double aHalfwidth)
1949 {
1950  /* sliding median to get the values, and its median deviation for the error */
1951  int i, n = cpl_table_get_nrow(aResp);
1952  for (i = 0; i < n; i++) {
1953  double lambda = cpl_table_get_double(aResp, "lambda", i, NULL);
1954  int j = i, j1 = i, j2 = i;
1955  /* search for the range */
1956  while (--j > 0 &&
1957  lambda - cpl_table_get_double(aResp, "lambda", j, NULL) <= aHalfwidth) {
1958  j1 = j;
1959  }
1960  j = i;
1961  while (++j < n &&
1962  cpl_table_get_double(aResp, "lambda", j, NULL) - lambda <= aHalfwidth) {
1963  j2 = j;
1964  }
1965  double *sens = cpl_table_get_data_double(aFluxObj->sensitivity, "sens"),
1966  *serr = cpl_table_get_data_double(aFluxObj->sensitivity, "serr");
1967  cpl_vector *v = cpl_vector_wrap(j2 - j1 + 1, sens + j1),
1968  *ve = cpl_vector_wrap(j2 - j1 + 1, serr + j1);
1969  double median = cpl_vector_get_median_const(v),
1970  mdev = muse_cplvector_get_adev_const(v, median),
1971  mederr = cpl_vector_get_median_const(ve);
1972  if (j2 == j1) { /* for single points, copy the error */
1973  mdev = cpl_table_get_double(aFluxObj->sensitivity, "serr", j1, NULL);
1974  }
1975  if (mdev < mederr) {
1976  mdev = mederr;
1977  }
1978 #if 0
1979  cpl_msg_debug(__func__, "%d %.3f %d...%d --> %f +/- %f", i, lambda, j1, j2,
1980  median, mdev);
1981 #endif
1982  cpl_vector_unwrap(v);
1983  cpl_vector_unwrap(ve);
1984  cpl_table_set_double(aResp, "response", i, median);
1985  cpl_table_set_double(aResp, "resperr", i, mdev);
1986  } /* for i (all aResp rows) */
1987 } /* muse_flux_get_response_table_smooth() */
1988 
1989 /*----------------------------------------------------------------------------*/
2007 /*----------------------------------------------------------------------------*/
2008 cpl_error_code
2010 {
2011  cpl_ensure_code(aFluxObj && aFluxObj->sensitivity, CPL_ERROR_NULL_INPUT);
2012  int nrow = cpl_table_get_nrow(aFluxObj->sensitivity);
2013  cpl_table *resp = muse_cpltable_new(muse_flux_responsetable_def, nrow);
2014  /* copy the lambda column from the sensitivity table */
2015  const double *lambdas = cpl_table_get_data_double_const(aFluxObj->sensitivity,
2016  "lambda");
2017  cpl_table_copy_data_double(resp, "lambda", lambdas);
2018  cpl_table_fill_column_window_double(resp, "response", 0, nrow, 0);
2019  /* first guess error: a constant value, the smoothing below will adapt it */
2020  cpl_table_fill_column_window_double(resp, "resperr", 0, nrow, 0.008);
2021 
2022  /* now fill it with real values, computed using a *
2023  * smoothing median over +/- 15 Angstrom width */
2024  muse_flux_get_response_table_smooth(aFluxObj, resp, 15.);
2025 
2026  aFluxObj->response = resp;
2027  return CPL_ERROR_NONE;
2028 } /* muse_flux_get_response_table() */
2029 
2030 /*----------------------------------------------------------------------------*/
2049 /*----------------------------------------------------------------------------*/
2050 cpl_error_code
2052 {
2053  cpl_ensure_code(aFluxObj && aFluxObj->sensitivity, CPL_ERROR_NULL_INPUT);
2054  cpl_table *tsens = aFluxObj->sensitivity;
2055  int nrow = cpl_table_get_nrow(tsens);
2056  cpl_table *tell = muse_cpltable_new(muse_flux_tellurictable_def, nrow);
2057  /* copy the lambda and tellcor columns from the sensitivity table */
2058  cpl_table_fill_column_window_double(tell, "lambda", 0, nrow, 0);
2059  cpl_table_copy_data_double(tell, "lambda",
2060  cpl_table_get_data_double_const(tsens, "lambda"));
2061  cpl_table_fill_column_window_double(tell, "ftelluric", 0, nrow, 0);
2062  cpl_table_copy_data_double(tell, "ftelluric",
2063  cpl_table_get_data_double_const(tsens, "tellcor"));
2064  /* XXX no (good) error estimates, for use constant error as starting point */
2065 #define TELL_MAX_ERR 0.1
2066 #define TELL_MIN_ERR 1e-4
2067  cpl_table_fill_column_window_double(tell, "ftellerr", 0, nrow, TELL_MAX_ERR);
2068 
2069  /* duplicate the tellcor column, to get the invalidity info; *
2070  * pad telluric correction factors with 1. in the new table */
2071  cpl_table_duplicate_column(tell, "tellcor", tsens, "tellcor");
2072  /* pad entries adjacent to the ones with a real telluric factor with 1 */
2073  cpl_table_unselect_all(tell);
2074  int irow;
2075  for (irow = 0; irow < nrow; irow++) {
2076  int err;
2077  cpl_table_get_double(tell, "tellcor", irow, &err); /* ignore value */
2078  if (err == 0) { /* has a valid entry -> do nothing */
2079  continue;
2080  }
2081  /* invalid entry, check previous one (again) */
2082  cpl_errorstate state = cpl_errorstate_get();
2083  double ftellcor = cpl_table_get_double(tell, "tellcor", irow - 1, &err);
2084  if (!cpl_errorstate_is_equal(state)) { /* recover from possible errors */
2085  cpl_errorstate_set(state);
2086  }
2087  if (err == 0 && ftellcor != 1.) { /* exist and is not 1 -> pad */
2088  cpl_table_set_double(tell, "ftelluric", irow, 1.);
2089  continue;
2090  }
2091  /* check the next one, too */
2092  state = cpl_errorstate_get();
2093  ftellcor = cpl_table_get_double(tell, "tellcor", irow + 1, &err);
2094  if (!cpl_errorstate_is_equal(state)) { /* recover from possible errors */
2095  cpl_errorstate_set(state);
2096  }
2097  if (err == 0 && ftellcor != 1.) { /* exist and is not 1 -> pad */
2098  cpl_table_set_double(tell, "ftelluric", irow, 1.);
2099  continue;
2100  }
2101  cpl_table_select_row(tell, irow); /* surrounded by invalid -> select */
2102  } /* for irow */
2103  cpl_table_erase_selected(tell); /* erase all the still invalid ones */
2104  cpl_table_erase_column(tell, "tellcor");
2105 
2106  /* next pass: adjust the error to be only about the distance between 1 and the value */
2107  nrow = cpl_table_get_nrow(tell);
2108  for (irow = 0; irow < nrow; irow++) {
2109  int err;
2110  double dftell = 1. - cpl_table_get_double(tell, "ftelluric", irow, &err),
2111  ftellerr = cpl_table_get_double(tell, "ftellerr", irow, &err);
2112  if (ftellerr > dftell) {
2113  ftellerr = fmax(dftell, TELL_MIN_ERR);
2114  cpl_table_set_double(tell, "ftellerr", irow, ftellerr);
2115  }
2116  } /* for irow */
2117 
2118  aFluxObj->telluric = tell;
2119  return CPL_ERROR_NONE;
2120 } /* muse_flux_get_telluric_table() */
2121 
2122 /*----------------------------------------------------------------------------*/
2159 /*----------------------------------------------------------------------------*/
2160 cpl_error_code
2161 muse_flux_calibrate(muse_pixtable *aPixtable, const cpl_table *aResponse,
2162  const cpl_table *aExtinction, const cpl_table *aTelluric)
2163 {
2164  cpl_ensure_code(aPixtable && aPixtable->header && aResponse,
2165  CPL_ERROR_NULL_INPUT);
2166  const char *unitdata = cpl_table_get_column_unit(aPixtable->table,
2167  MUSE_PIXTABLE_DATA);
2168  cpl_ensure_code(unitdata && !strncmp(unitdata, "count", 6),
2169  CPL_ERROR_INCOMPATIBLE_INPUT);
2170  /* warn for non-critical failure */
2171  if (!aExtinction) {
2172  cpl_msg_warning(__func__, "%s missing!", MUSE_TAG_EXTINCT_TABLE);
2173  }
2174 
2175  double exptime = muse_pfits_get_exptime(aPixtable->header);
2176  if (exptime <= 0.) {
2177  cpl_msg_warning(__func__, "Non-positive EXPTIME, not doing flux calibration!");
2178  return CPL_ERROR_ILLEGAL_INPUT;
2179  }
2180  double airmass = muse_astro_airmass(aPixtable->header);
2181  if (airmass < 0) {
2182  cpl_msg_warning(__func__, "Airmass unknown, not doing extinction "
2183  "correction: %s", cpl_error_get_message());
2184  /* reset to zero so that it has no effect */
2185  airmass = 0.;
2186  }
2187 
2188  cpl_table *telluric = NULL;
2189  if (aTelluric) {
2190  /* duplicate the telluric correction table to apply the airmass correction */
2191  telluric = cpl_table_duplicate(aTelluric);
2192  /* use negative exponent to be able to multiply (instead of divide) *
2193  * by the correction factor (should be slightly faster) */
2194  cpl_table_power_column(telluric, "ftelluric", -airmass);
2195  }
2196 
2197  cpl_msg_info(__func__, "Starting flux calibration (exptime=%.2fs, airmass=%.4f),"
2198  " %s telluric correction", exptime, airmass,
2199  aTelluric ? "with" : "without ("MUSE_TAG_STD_TELLURIC" not given)");
2200  float *lambda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
2201  *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
2202  *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
2203  cpl_size i, nrow = muse_pixtable_get_nrow(aPixtable);
2204  #pragma omp parallel for default(none) /* as req. by Ralf */ \
2205  shared(aExtinction, aResponse, airmass, data, exptime, lambda, nrow, \
2206  stat, telluric)
2207  for (i = 0; i < nrow; i++) {
2208  /* double values for intermediate results of this row */
2209  double ddata = data[i], dstat = stat[i];
2210 
2211  /* correct for extinction */
2212  if (aExtinction) {
2213  double fext = pow(10., 0.4 * airmass
2214  * muse_flux_response_interpolate(aExtinction,
2215  lambda[i], NULL,
2217 #if 0
2218  printf("%f, data/stat = %f/%f -> ", fext, ddata, dstat);
2219 #endif
2220  ddata *= fext;
2221  dstat *= (fext * fext);
2222 #if 0
2223  printf(" --> %f/%f\n", ddata, dstat), fflush(NULL);
2224 #endif
2225  }
2226 
2227  double dlambda = kMuseSpectralSamplingA,
2228 #if 0 /* XXX the difference in lambda coverage per pixel seems to be *
2229  * corrected for by the flat-field, so assuming a constant *
2230  * value for all pixels seems to be the correct thing to do */
2231  /* XXX Fit of *
2232  * g(x)=c0+c1*x+c2*x**2+c3*x**3+c4*x**4+c5*x**5 *
2233  * to data of lambda vs. dlambda of pixel 1900 to 2100 using a *
2234  * WAVE_MAP produced using INM data of 2010-03-12 results in the *
2235  * following polynomial coefficients: *
2236  * c0 = 1.26279 +/- 0.002745 (0.2173%) *
2237  * c1 = 7.12041e-05 +/- 2.087e-06 (2.931%) *
2238  * c2 = -2.97879e-08 +/- 6.246e-10 (2.097%) *
2239  * c3 = 4.653e-12 +/- 9.205e-14 (1.978%) *
2240  * c4 = -3.30195e-16 +/- 6.682e-18 (2.024%) *
2241  * c5 = 8.70389e-21 +/- 1.913e-22 (2.198%) *
2242  * So that's hardcoded here for now. */
2243  double dlambda = 1.26279
2244  + 7.12041e-05 * lambda[i]
2245  + -2.97879e-08 * pow(lambda[i], 2)
2246  + 4.653e-12 * pow(lambda[i], 3)
2247  + -3.30195e-16 * pow(lambda[i], 4)
2248  + 8.70389e-21 * pow(lambda[i], 5),
2249 #endif
2250  dlamerr = 0.02, /* XXX assume fixed error for the moment */
2251  /* resp/rerr get returned in mag units as in the table */
2252  rerr, resp = muse_flux_response_interpolate(aResponse, lambda[i],
2253  &rerr,
2255  /* convert from 2.5 log10(x) to non-log flux units */
2256  resp = pow(10., 0.4 * resp);
2257  /* magerr = 2.5 / log(10) * error / flux (see IRAF phot docs) *
2258  * ==> error = magerr / 2.5 * log(10) * flux */
2259  rerr = rerr * CPL_MATH_LN10 * resp / 2.5;
2260 #if 0
2261  printf("%f/%f/%f, %e/%e, data/stat = %e/%e -> ", lambda[i], dlambda, dlamerr, resp, rerr,
2262  ddata, dstat);
2263 #endif
2264  dstat = dstat * pow((1./(resp * exptime * dlambda)), 2)
2265  + pow(ddata * rerr / (resp*resp * exptime * dlambda), 2)
2266  + pow(ddata * dlamerr / (resp * exptime * dlambda*dlambda), 2);
2267  ddata /= (resp * exptime * dlambda);
2268 #if 0
2269  printf("%e/%e\n", ddata, dstat), fflush(NULL);
2270 #endif
2271 
2272  /* now convert to the float values to be stored in the pixel table, *
2273  * scaled by kMuseFluxUnitFactor to keep the floats from underflowing */
2274  ddata *= kMuseFluxUnitFactor;
2275  dstat *= kMuseFluxStatFactor;
2276 
2277  /* do the telluric correction, if the wavelength is redward of the start *
2278  * of the telluric regions, and if a telluric correction was given */
2279  if (lambda[i] < kTelluricBands[0][0] || !telluric) {
2280  data[i] = ddata;
2281  stat[i] = dstat;
2282  continue; /* skip telluric correction in the blue */
2283  }
2284  double terr, tell = muse_flux_response_interpolate(telluric, lambda[i],
2285  &terr,
2287  data[i] = ddata * tell;
2288  stat[i] = tell*tell * dstat + ddata*ddata * terr*terr;
2289  } /* for i (pixel table rows) */
2290  cpl_table_delete(telluric); /* NULL check done there... */
2291 
2292  /* now set the table column headers reflecting the flux units used */
2293  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_DATA,
2294  kMuseFluxUnitString);
2295  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_STAT,
2296  kMuseFluxStatString);
2297 
2298  /* add the status header */
2299  cpl_propertylist_update_bool(aPixtable->header, MUSE_HDR_PT_FLUXCAL,
2300  CPL_TRUE);
2301  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_FLUXCAL,
2302  MUSE_HDR_PT_FLUXCAL_COMMENT);
2303  return CPL_ERROR_NONE;
2304 } /* muse_flux_calibrate() */
2305 
Structure definition of a MUSE datacube.
Definition: muse_datacube.h:48
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
Definition: muse_image.c:84
cpl_size muse_pixtable_get_nrow(muse_pixtable *aPixtable)
get the number of rows within the pixel table
#define MUSE_HDR_PT_FLUXCAL
cpl_image * data
the data extension
Definition: muse_image.h:47
const muse_cpltable_def muse_flux_responsetable_def[]
MUSE response table definition.
Definition: muse_flux.c:1900
muse_flux_profile_type
Type of optimal profile to use.
Definition: muse_flux.h:97
cpl_image * stat
the statistics extension
Definition: muse_image.h:65
void muse_datacube_delete(muse_datacube *aCube)
Deallocate memory associated to a muse_datacube object.
cpl_error_code muse_flux_reference_table_check(cpl_table *aTable)
Check and/or adapt the standard flux reference table format.
Definition: muse_flux.c:162
const muse_cpltable_def muse_response_tellbands_def[]
Table definition for a telluric bands table.
Definition: muse_flux.c:1303
Structure definition of MUSE three extension FITS file.
Definition: muse_image.h:41
muse_flux_object * muse_flux_object_new(void)
Allocate memory for a new muse_flux_object object.
Definition: muse_flux.c:66
static double muse_flux_reference_table_sampling(cpl_table *aTable)
Compute average sampling for a MUSE-format flux reference table.
Definition: muse_flux.c:114
cpl_table * table
The pixel table.
cpl_propertylist * header
the FITS header
Definition: muse_image.h:73
cpl_error_code muse_cpltable_check(const cpl_table *aTable, const muse_cpltable_def *aDef)
Check whether the table contains the fields of the definition.
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. ...
muse_resampling_crstats_type crtype
cpl_image * dq
the data quality extension
Definition: muse_image.h:57
cpl_table * muse_cpltable_new(const muse_cpltable_def *aDef, cpl_size aLength)
Create an empty table according to the specified definition.
cpl_error_code muse_flux_get_telluric_table(muse_flux_object *aFluxObj)
Get the table of the telluric correction.
Definition: muse_flux.c:2051
double muse_flux_response_interpolate(const cpl_table *aResponse, double aLambda, double *aError, muse_flux_interpolation_type aType)
Compute linearly interpolated response of some kind at given wavelength.
Definition: muse_flux.c:322
double muse_pfits_get_fwhm_end(const cpl_propertylist *aHeaders)
find out the ambient seeing at end of exposure (in arcsec)
Definition: muse_pfits.c:858
cpl_error_code muse_flux_get_response_table(muse_flux_object *aFluxObj)
Get the table of the standard star response function.
Definition: muse_flux.c:2009
Structure definition of MUSE pixel table.
Flux object to store data needed while computing the flux calibration.
Definition: muse_flux.h:57
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:1473
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:1829
static void muse_flux_get_response_table_smooth(muse_flux_object *aFluxObj, cpl_table *aResp, double aHalfwidth)
Get the table of the standard star response function.
Definition: muse_flux.c:1947
double muse_pfits_get_fwhm_start(const cpl_propertylist *aHeaders)
find out the ambient seeing at start of exposure (in arcsec)
Definition: muse_pfits.c:840
cpl_error_code muse_flux_calibrate(muse_pixtable *aPixtable, const cpl_table *aResponse, const cpl_table *aExtinction, const cpl_table *aTelluric)
Convert the input pixel table from counts to fluxes.
Definition: muse_flux.c:2161
cpl_imagelist * data
the cube containing the actual data values
Definition: muse_datacube.h:76
muse_flux_interpolation_type
Type of table interpolation to use.
Definition: muse_flux.h:84
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:81
const muse_cpltable_def muse_flux_tellurictable_def[]
MUSE telluric correction table definition.
Definition: muse_flux.c:1921
cpl_error_code muse_pixtable_save(muse_pixtable *aPixtable, const char *aFilename)
Save a MUSE pixel table to a file on disk.
cpl_error_code muse_flux_integrate_std(muse_pixtable *aPixtable, muse_flux_profile_type aProfile, muse_flux_object *aFluxObj)
Integrate the flux of the standard star(s) in the field over all wavelengths.
Definition: muse_flux.c:1006
cpl_propertylist * header
the FITS header
Definition: muse_datacube.h:57
cpl_error_code muse_flux_response_compute(muse_flux_object *aFluxObj, double aAirmass, const cpl_table *aReference, const cpl_table *aTellBands, const cpl_table *aExtinct)
Compare measured flux distribution over wavelength with calibrated stellar fluxes and derive instrume...
Definition: muse_flux.c:1833
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:396
double muse_pfits_get_exptime(const cpl_propertylist *aHeaders)
find out the exposure time
Definition: muse_pfits.c:278
Resampling parameters.
Definition of a cpl table structure.
void muse_resampling_params_delete(muse_resampling_params *aParams)
Delete a resampling parameters structure.
double muse_cplvector_get_adev_const(const cpl_vector *aVector, double aCenter)
Compute the average absolute deviation of a (constant) vector.
muse_image * muse_image_new(void)
Allocate memory for a new muse_image object.
Definition: muse_image.c:65
void muse_flux_object_delete(muse_flux_object *aFluxObj)
Deallocate memory associated to a muse_flux_object.
Definition: muse_flux.c:84
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:623
muse_ins_mode muse_pfits_get_mode(const cpl_propertylist *aHeaders)
find out the observation mode
Definition: muse_pfits.c:1003
double muse_astro_airmass(cpl_propertylist *aHeader)
Derive the effective airmass of an observation from information in a FITS header. ...
Definition: muse_astro.c:322
cpl_imagelist * stat
the cube containing the data variance
Definition: muse_datacube.h:86
cpl_propertylist * header
The FITS header.