MUSE Pipeline Reference Manual  2.1.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  * This file is part of the MUSE Instrument Pipeline
5  * Copyright (C) 2005-2017 European Southern Observatory
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software
19  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
20  */
21 
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
25 
26 /*----------------------------------------------------------------------------*
27  * Includes *
28  *----------------------------------------------------------------------------*/
29 #include <cpl.h>
30 #include <math.h>
31 #include <string.h>
32 
33 #include "muse_flux.h"
34 #include "muse_instrument.h"
35 
36 #include "muse_astro.h"
37 #include "muse_cplwrappers.h"
38 #include "muse_pfits.h"
39 #include "muse_quality.h"
40 #include "muse_resampling.h"
41 #include "muse_utils.h"
42 #include "muse_wcs.h"
43 
44 /*----------------------------------------------------------------------------*/
48 /*----------------------------------------------------------------------------*/
49 
52 /*----------------------------------------------------------------------------*/
66 /*----------------------------------------------------------------------------*/
69 {
70  muse_flux_object *flux = cpl_calloc(1, sizeof(muse_flux_object));
71  /* signify non-filled reference coordinates */
72  flux->raref = NAN;
73  flux->decref = NAN;
74  return flux;
75 }
76 
77 /*----------------------------------------------------------------------------*/
87 /*----------------------------------------------------------------------------*/
88 void
90 {
91  if (!aFluxObj) {
92  return;
93  }
94  muse_datacube_delete(aFluxObj->cube);
95  aFluxObj->cube = NULL;
96  muse_image_delete(aFluxObj->intimage);
97  aFluxObj->intimage = NULL;
98  cpl_table_delete(aFluxObj->reference);
99  aFluxObj->reference = NULL;
100  cpl_table_delete(aFluxObj->sensitivity);
101  aFluxObj->sensitivity = NULL;
102  muse_table_delete(aFluxObj->response);
103  aFluxObj->response = NULL;
104  muse_table_delete(aFluxObj->telluric);
105  aFluxObj->telluric = NULL;
106  cpl_table_delete(aFluxObj->tellbands);
107  aFluxObj->tellbands = NULL;
108  cpl_free(aFluxObj);
109 }
110 
111 /*----------------------------------------------------------------------------*/
119 /*----------------------------------------------------------------------------*/
120 static double
122 {
123  cpl_ensure(aTable, CPL_ERROR_NULL_INPUT, 0.);
124  cpl_table_unselect_all(aTable);
125  cpl_table_or_selected_double(aTable, "lambda", CPL_NOT_LESS_THAN,
126  kMuseNominalLambdaMin);
127  cpl_table_and_selected_double(aTable, "lambda", CPL_NOT_GREATER_THAN,
128  kMuseNominalLambdaMax);
129  cpl_size nsel = cpl_table_count_selected(aTable);
130  cpl_array *asel = cpl_table_where_selected(aTable);
131  cpl_size *sel = cpl_array_get_data_cplsize(asel);
132  double lmin = cpl_table_get_double(aTable, "lambda", sel[0], NULL),
133  lmax = cpl_table_get_double(aTable, "lambda", sel[nsel - 1], NULL);
134  cpl_array_delete(asel);
135  return (lmax - lmin) / nsel;
136 } /* muse_flux_reference_table_sampling() */
137 
138 /*----------------------------------------------------------------------------*/
168 /*----------------------------------------------------------------------------*/
169 cpl_error_code
171 {
172  cpl_ensure_code(aTable, CPL_ERROR_NULL_INPUT);
173 
174  const char *flam1 = "erg/s/cm**2/Angstrom",
175  *flam2 = "erg/s/cm^2/Angstrom";
176 
177  cpl_error_code rc = CPL_ERROR_NONE;
178  cpl_errorstate prestate = cpl_errorstate_get();
179  /* check two different types of tables: MUSE specific or HST CALSPEC */
180  if (cpl_table_has_column(aTable, "lambda") &&
181  cpl_table_has_column(aTable, "flux") &&
182  cpl_table_get_column_unit(aTable, "lambda") &&
183  cpl_table_get_column_unit(aTable, "flux") &&
184  !strncmp(cpl_table_get_column_unit(aTable, "lambda"), "Angstrom", 9) &&
185  (!strncmp(cpl_table_get_column_unit(aTable, "flux"), flam1, strlen(flam1)) ||
186  !strncmp(cpl_table_get_column_unit(aTable, "flux"), flam2, strlen(flam2)))) {
187  /* normal case: MUSE STD_FLUX_TABLE as specified; still need to *
188  * check, if we need to convert the column types (could be e.g. float) */
189  if (cpl_table_get_column_type(aTable, "lambda") != CPL_TYPE_DOUBLE) {
190  cpl_msg_debug(__func__, "Casting lambda column to double");
191  cpl_table_cast_column(aTable, "lambda", NULL, CPL_TYPE_DOUBLE);
192  }
193  if (cpl_table_get_column_type(aTable, "flux") != CPL_TYPE_DOUBLE) {
194  cpl_msg_debug(__func__, "Casting flux column to double");
195  cpl_table_cast_column(aTable, "flux", NULL, CPL_TYPE_DOUBLE);
196  }
197  /* check optional column */
198  if (cpl_table_has_column(aTable, "fluxerr")) {
199  if (cpl_table_get_column_type(aTable, "fluxerr") != CPL_TYPE_DOUBLE) {
200  cpl_msg_debug(__func__, "Casting fluxerr column to double");
201  cpl_table_cast_column(aTable, "fluxerr", NULL, CPL_TYPE_DOUBLE);
202  }
203  const char *unit = cpl_table_get_column_unit(aTable, "fluxerr");
204  if (!unit || (strncmp(unit, flam1, strlen(flam1)) &&
205  strncmp(unit, flam2, strlen(flam2)))) {
206  cpl_msg_debug(__func__, "Erasing fluxerr column because of unexpected "
207  "unit (%s)", unit);
208  cpl_table_erase_column(aTable, "fluxerr"); /* wrong unit, erase */
209  }
210  } /* if has fluxerr */
211  cpl_msg_info(__func__, "Found MUSE format, average sampling %.3f Angstrom/bin"
212  " over MUSE range", muse_flux_reference_table_sampling(aTable));
213  } else if (cpl_table_has_column(aTable, "WAVELENGTH") &&
214  cpl_table_has_column(aTable, "FLUX") &&
215  cpl_table_get_column_unit(aTable, "WAVELENGTH") &&
216  cpl_table_get_column_unit(aTable, "FLUX") &&
217  !strncmp(cpl_table_get_column_unit(aTable, "WAVELENGTH"), "ANGSTROMS", 10) &&
218  !strncmp(cpl_table_get_column_unit(aTable, "FLUX"), "FLAM", 5)) {
219 #if 0
220  printf("input HST CALSPEC table:\n");
221  cpl_table_dump_structure(aTable, stdout);
222  cpl_table_dump(aTable, cpl_table_get_nrow(aTable)/2, 3, stdout);
223  fflush(stdout);
224 #endif
225  /* other allowed case: HST CALSPEC format */
226  cpl_table_cast_column(aTable, "WAVELENGTH", "lambda", CPL_TYPE_DOUBLE);
227  cpl_table_cast_column(aTable, "FLUX", "flux", CPL_TYPE_DOUBLE);
228  cpl_table_erase_column(aTable, "WAVELENGTH");
229  cpl_table_erase_column(aTable, "FLUX");
230  cpl_table_set_column_unit(aTable, "lambda", "Angstrom");
231  cpl_table_set_column_unit(aTable, "flux", flam1);
232  /* if the table comes with the typical STATERROR/SYSERROR separation, *
233  * convert them into a single combined fluxerr column */
234  if (cpl_table_has_column(aTable, "STATERROR") &&
235  cpl_table_has_column(aTable, "SYSERROR") &&
236  cpl_table_get_column_unit(aTable, "STATERROR") &&
237  cpl_table_get_column_unit(aTable, "SYSERROR") &&
238  !strncmp(cpl_table_get_column_unit(aTable, "STATERROR"), "FLAM", 5) &&
239  !strncmp(cpl_table_get_column_unit(aTable, "SYSERROR"), "FLAM", 5)) {
240  /* Cast to double before, not to lose precision, then compute *
241  * fluxerr = sqrt(STATERROR**2 + SYSERROR**2) */
242  cpl_table_cast_column(aTable, "STATERROR", "fluxerr", CPL_TYPE_DOUBLE);
243  cpl_table_erase_column(aTable, "STATERROR");
244  cpl_table_cast_column(aTable, "SYSERROR", NULL, CPL_TYPE_DOUBLE);
245  cpl_table_power_column(aTable, "fluxerr", 2);
246  cpl_table_power_column(aTable, "SYSERROR", 2);
247  cpl_table_add_columns(aTable, "fluxerr", "SYSERROR");
248  cpl_table_erase_column(aTable, "SYSERROR");
249  cpl_table_power_column(aTable, "fluxerr", 0.5);
250  cpl_table_set_column_unit(aTable, "fluxerr", flam1);
251  } /* if error columns */
252  /* XXX how to handle invalid entries in the STATERROR column *
253  * in telluric regions (e.g. in gd105_005.fits)? */
254  /* XXX how to handle DATAQUAL column? */
255 
256  /* erase further columns we don't need */
257  if (cpl_table_has_column(aTable, "FWHM")) {
258  cpl_table_erase_column(aTable, "FWHM");
259  }
260  if (cpl_table_has_column(aTable, "DATAQUAL")) {
261  cpl_table_erase_column(aTable, "DATAQUAL");
262  }
263  if (cpl_table_has_column(aTable, "TOTEXP")) {
264  cpl_table_erase_column(aTable, "TOTEXP");
265  }
266  /* convert from vacuum to air wavelengths */
267  cpl_size irow, nrow = cpl_table_get_nrow(aTable);
268  for (irow = 0; irow < nrow; irow++) {
269  double lambda = cpl_table_get_double(aTable, "lambda", irow, NULL);
270  cpl_table_set_double(aTable, "lambda", irow,
272  } /* for irow (all table rows) */
273 #if 0
274  printf("converted HST CALSPEC table:\n");
275  cpl_table_dump_structure(aTable, stdout);
276  cpl_table_dump(aTable, cpl_table_get_nrow(aTable)/2, 3, stdout);
277  fflush(stdout);
278 #endif
279  cpl_msg_info(__func__, "Found HST CALSPEC format on input, converted to "
280  "MUSE format; average sampling %.3f Angstrom/bin over MUSE "
281  "range (assumed vacuum wavelengths on input, converted to air).",
283  } else {
284  cpl_msg_error(__func__, "Unknown format found!");
285 #if 0
286  cpl_table_dump_structure(aTable, stdout);
287 #endif
288  rc = CPL_ERROR_INCOMPATIBLE_INPUT;
289  } /* else: no recognized format */
290 
291  /* check for errors in the above (casting!) before returning */
292  if (!cpl_errorstate_is_equal(prestate)) {
293  rc = cpl_error_get_code();
294  }
295  return rc;
296 } /* muse_flux_reference_table_check() */
297 
298 /*----------------------------------------------------------------------------*/
336 /*----------------------------------------------------------------------------*/
337 double
338 muse_flux_response_interpolate(const cpl_table *aResponse, double aLambda,
339  double *aError, muse_flux_interpolation_type aType)
340 {
341  double rv = 0.;
342  if (aType == MUSE_FLUX_TELLURIC) {
343  rv = 1.;
344  }
345  cpl_ensure(aResponse, CPL_ERROR_NULL_INPUT, rv);
346  int size = cpl_table_get_nrow(aResponse);
347  cpl_ensure(size > 0, cpl_error_get_code(), rv);
348 
349  /* access the correct table column(s) depending on the type */
350  const double *lbda = cpl_table_get_data_double_const(aResponse, "lambda"),
351  *resp = NULL, *rerr = NULL;
352  switch (aType) {
354  resp = cpl_table_get_data_double_const(aResponse, "throughput");
355  break;
356  case MUSE_FLUX_RESP_FLUX:
357  resp = cpl_table_get_data_double_const(aResponse, "response");
358  if (aError) {
359  rerr = cpl_table_get_data_double_const(aResponse, "resperr");
360  }
361  break;
363  resp = cpl_table_get_data_double_const(aResponse, "flux");
364  if (aError) {
365  rerr = cpl_table_get_data_double_const(aResponse, "fluxerr");
366  }
367  break;
369  resp = cpl_table_get_data_double_const(aResponse, "extinction");
370  break;
371  case MUSE_FLUX_TELLURIC:
372  resp = cpl_table_get_data_double_const(aResponse, "ftelluric");
373  if (aError) {
374  rerr = cpl_table_get_data_double_const(aResponse, "ftellerr");
375  }
376  break;
377  default:
378  cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
379  return rv;
380  } /* switch aType */
381  cpl_ensure(lbda && resp, cpl_error_get_code(), rv);
382  if (aError) {
383  cpl_ensure(rerr, cpl_error_get_code(), rv);
384  }
385 
386  /* outside wavelength range of table */
387  if (aLambda < lbda[0]) {
388  return rv;
389  }
390  if (aLambda > lbda[size-1]) {
391  return rv;
392  }
393 
394  /* binary search for the correct wavelength */
395  double response = rv, resperror = 0.;
396  int l = 0, r = size - 1, /* left and right end */
397  m = (l + r) / 2; /* middle index */
398  while (CPL_TRUE) {
399  if (aLambda >= lbda[m] && aLambda <= lbda[m+1]) {
400  /* found right interval, so interpolate */
401  double lquot = (aLambda - lbda[m]) / (lbda[m+1] - lbda[m]);
402  response = resp[m] + (resp[m+1] - resp[m]) * lquot;
403  if (rerr) { /* missing error information should be non-fatal */
404  /* checked again that the derivatives leading to this error estimate *
405  * are correct; apparently it's normal then, that the resulting *
406  * errors can be smaller than the two separate input errors */
407  resperror = sqrt(pow(rerr[m] * (1 - lquot), 2.)
408  + pow(rerr[m+1] * lquot, 2.));
409  }
410 #if 0
411  cpl_msg_debug(__func__, "Found at m=%d (%f: %f+/-%f) / "
412  "m+1=%d (%f: %f+/-%f) -> %f: %f+/-%f",
413  m, lbda[m], resp[m], rerr ? rerr[m] : 0.,
414  m+1, lbda[m+1], resp[m+1], rerr ? rerr[m+1] : 0.,
415  aLambda, response, resperror);
416 #endif
417  break;
418  }
419  /* create next interval */
420  if (aLambda < lbda[m]) {
421  r = m;
422  }
423  if (aLambda > lbda[m]) {
424  l = m;
425  }
426  m = (l + r) / 2;
427  } /* while */
428 
429 #if 0
430  cpl_msg_debug(__func__, "Response %g+/-%g at lambda=%fA", response, resperror,
431  aLambda);
432 #endif
433  if (aError && rerr) {
434  *aError = resperror;
435  }
436  return response;
437 } /* muse_flux_response_interpolate() */
438 
439 /*----------------------------------------------------------------------------*/
454 /*----------------------------------------------------------------------------*/
455 static double
456 muse_flux_image_sky(cpl_image *aImage, double aX, double aY, double aHalfSize,
457  unsigned int aDSky, float *aError)
458 {
459  if (aError) {
460  *aError = FLT_MAX;
461  }
462  /* coordinates of inner window */
463  int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
464  y1 = aY - aHalfSize, y2 = aY + aHalfSize;
465  unsigned char nskyarea = 0;
466  double skylevel = 0., skyerror = 0.;
467  /* left */
468  cpl_errorstate state = cpl_errorstate_get();
469  cpl_stats_mode mode = CPL_STATS_MEDIAN | CPL_STATS_MEDIAN_DEV;
470  cpl_stats *s = cpl_stats_new_from_image_window(aImage, mode,
471  x1 - aDSky, y1, x1 - 1, y2);
472  if (s) {
473  /* only if there was no error, the area was inside the image *
474  * boundaries and we can use it, otherwise the value is undefined */
475  nskyarea++;
476  skylevel += cpl_stats_get_median(s);
477  skyerror += pow(cpl_stats_get_median_dev(s), 2);
478  cpl_stats_delete(s);
479  }
480  /* right */
481  s = cpl_stats_new_from_image_window(aImage,mode,
482  x2 + 1, y1, x2 + aDSky, y2);
483  if (s) {
484  nskyarea++;
485  skylevel += cpl_stats_get_median(s);
486  skyerror += pow(cpl_stats_get_median_dev(s), 2);
487  cpl_stats_delete(s);
488  }
489  /* bottom */
490  s = cpl_stats_new_from_image_window(aImage,mode,
491  x1, y1 - aDSky, x2, y1 - 1);
492  if (s) {
493  nskyarea++;
494  skylevel += cpl_stats_get_median(s);
495  skyerror += pow(cpl_stats_get_median_dev(s), 2);
496  cpl_stats_delete(s);
497  }
498  /* top */
499  s = cpl_stats_new_from_image_window(aImage,mode,
500  x1, y2 + 1, x2, y2 + aDSky);
501  if (s) {
502  nskyarea++;
503  skylevel += cpl_stats_get_median(s);
504  skyerror += pow(cpl_stats_get_median_dev(s), 2);
505  cpl_stats_delete(s);
506  }
507  if (nskyarea == 0) {
508  return 0.;
509  }
510  skylevel /= nskyarea;
511  skyerror = sqrt(skyerror) / nskyarea;
512  if (!cpl_errorstate_is_equal(state)) { /* reset error code */
513  cpl_errorstate_set(state);
514  }
515 #if 0
516  cpl_msg_debug(__func__, "skylevel = %f +/- %f (%u sky areas)",
517  skylevel, skyerror, nskyarea);
518 #endif
519  if (aError) {
520  *aError = skyerror;
521  }
522  return skylevel;
523 } /* muse_flux_image_sky() */
524 
525 /*----------------------------------------------------------------------------*/
542 /*----------------------------------------------------------------------------*/
543 static double
544 muse_flux_image_gaussian(cpl_image *aImage, cpl_image *aImErr, double aX,
545  double aY, double aHalfSize, unsigned int aDSky,
546  unsigned int aMaxBad, float *aFErr)
547 {
548  /* there is no simple way to count bad pixels inside an *
549  * image window, so ignore this argument at the moment */
550  UNUSED_ARGUMENT(aMaxBad);
551 
552  if (aFErr) { /* set high variance for an error case */
553  *aFErr = FLT_MAX;
554  }
555 
556  cpl_array *params = cpl_array_new(7, CPL_TYPE_DOUBLE),
557  *parerr = cpl_array_new(7, CPL_TYPE_DOUBLE);
558  /* Set some first-guess parameters to help the fitting function. *
559  * Just set background and central position, cpl_fit_image_gaussian() *
560  * finds good defaults for everything else. */
561  cpl_errorstate state = cpl_errorstate_get();
562  double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky, NULL);
563  if (!cpl_errorstate_is_equal(state)) {
564  /* if background determination fails, a default of 0 should *
565  * be good enough, so that we can ignore this failure */
566  cpl_errorstate_set(state);
567  }
568  cpl_array_set_double(params, 0, skylevel);
569  cpl_array_set_double(params, 3, aX);
570  cpl_array_set_double(params, 4, aY);
571  double rms = 0, chisq = 0;
572  /* function wants full widths but at most out to the image boundary */
573  int nx = cpl_image_get_size_x(aImage),
574  ny = cpl_image_get_size_y(aImage),
575  xsize = fmin(aHalfSize, fmin(aX - 1., nx - aX)) * 2,
576  ysize = fmin(aHalfSize, fmin(aY - 1., ny - aY)) * 2;
577  cpl_error_code rc = cpl_fit_image_gaussian(aImage, aImErr, aX, aY, xsize, ysize,
578  params, parerr, NULL, &rms, &chisq,
579  NULL, NULL, NULL, NULL, NULL);
580  if (rc != CPL_ERROR_NONE) {
581  if (rc != CPL_ERROR_ILLEGAL_INPUT) {
582  cpl_msg_debug(__func__, "rc = %d: %s", rc, cpl_error_get_message());
583  }
584  cpl_array_delete(params);
585  cpl_array_delete(parerr);
586  return 0;
587  }
588  double flux = cpl_array_get_double(params, 1, NULL),
589  ferr = cpl_array_get_double(parerr, 1, NULL);
590 #if 0 /* DEBUG */
591  double fwhmx = cpl_array_get_double(params, 5, NULL) * CPL_MATH_FWHM_SIG,
592  fwhmy = cpl_array_get_double(params, 6, NULL) * CPL_MATH_FWHM_SIG;
593  cpl_msg_debug(__func__, "%.3f,%.3f: %g+/-%g (bg: %g, FWHM: %.3f,%.3f, %g, %g)",
594  cpl_array_get_double(params, 3, NULL), cpl_array_get_double(params, 4, NULL),
595  flux, ferr, cpl_array_get_double(params, 0, NULL), fwhmx, fwhmy,
596  rms, chisq);
597 #endif
598 #if 0 /* DEBUG */
599  cpl_msg_debug(__func__, "skylevel = %f", cpl_array_get_double(params, 0, NULL));
600  cpl_msg_debug(__func__, "measured flux %f +/- %f", flux, ferr);
601 #endif
602  cpl_array_delete(params);
603  cpl_array_delete(parerr);
604  if (aFErr) {
605  *aFErr = ferr;
606  }
607  return flux;
608 } /* muse_flux_image_gaussian() */
609 
610 /*----------------------------------------------------------------------------*/
629 /*----------------------------------------------------------------------------*/
630 static double
631 muse_flux_image_moffat(cpl_image *aImage, cpl_image *aImErr, double aX,
632  double aY, double aHalfSize, unsigned int aDSky,
633  unsigned int aMaxBad, float *aFErr)
634 {
635  if (aFErr) { /* set high variance for an error case */
636  *aFErr = FLT_MAX;
637  }
638  /* extract image regions around the fiducial peak into matrix and vectors */
639  int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
640  y1 = aY - aHalfSize, y2 = aY + aHalfSize,
641  nx = cpl_image_get_size_x(aImage),
642  ny = cpl_image_get_size_y(aImage);
643  if (x1 < 1) {
644  x1 = 1;
645  }
646  if (x2 > nx) {
647  x2 = nx;
648  }
649  if (y1 < 1) {
650  y1 = 1;
651  }
652  if (y2 > ny) {
653  y2 = ny;
654  }
655  int npoints = (x2 - x1 + 1) * (y2 - y1 + 1);
656 
657  cpl_matrix *pos = cpl_matrix_new(npoints, 2);
658  cpl_vector *values = cpl_vector_new(npoints),
659  *errors = cpl_vector_new(npoints);
660  float *derr = cpl_image_get_data_float(aImErr);
661  int i, idx = 0;
662  for (i = x1; i <= x2; i++) {
663  int j;
664  for (j = y1; j <= y2; j++) {
665  int err;
666  double data = cpl_image_get(aImage, i, j, &err);
667  if (err) { /* bad pixel or error */
668  continue;
669  }
670  cpl_matrix_set(pos, idx, 0, i);
671  cpl_matrix_set(pos, idx, 1, j);
672  cpl_vector_set(values, idx, data);
673  cpl_vector_set(errors, idx, derr[(i-1) + (j-1)*nx]);
674  idx++;
675  } /* for j (y pixels) */
676  } /* for i (x pixels) */
677  /* need at least something like 4x4 pixels for a solid fit; *
678  * also check missing entries against max number of bad pixels */
679  if (idx < 16 || (npoints - idx) > (int)aMaxBad) {
680  cpl_matrix_delete(pos);
681  cpl_vector_delete(values);
682  cpl_vector_delete(errors);
683  return 0;
684  }
685  cpl_matrix_set_size(pos, idx, 2);
686  cpl_vector_set_size(values, idx);
687  cpl_vector_set_size(errors, idx);
688 
689  cpl_array *params = cpl_array_new(8, CPL_TYPE_DOUBLE),
690  *parerr = cpl_array_new(8, CPL_TYPE_DOUBLE);
691  /* Set some first-guess parameters to help the fitting function. */
692  cpl_errorstate state = cpl_errorstate_get();
693  double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky, NULL);
694  if (!cpl_errorstate_is_equal(state)) {
695  /* if background determination fails, a default of 0 should *
696  * be good enough, so that we can ignore this failure */
697  cpl_errorstate_set(state);
698  }
699  cpl_array_set_double(params, 0, skylevel);
700  cpl_array_set_double(params, 2, aX);
701  cpl_array_set_double(params, 3, aY);
702  double rms = 0, chisq = 0;
703  cpl_error_code rc = muse_utils_fit_moffat_2d(pos, values, errors,
704  params, parerr, NULL,
705  &rms, &chisq);
706  cpl_matrix_delete(pos);
707  cpl_vector_delete(values);
708  cpl_vector_delete(errors);
709  if (rc != CPL_ERROR_NONE) {
710  if (rc != CPL_ERROR_ILLEGAL_INPUT) {
711  cpl_msg_debug(__func__, "rc = %d: %s", rc, cpl_error_get_message());
712  }
713  cpl_array_delete(params);
714  cpl_array_delete(parerr);
715  return 0;
716  }
717 
718  double flux = cpl_array_get_double(params, 1, NULL);
719  if (aFErr) {
720  *aFErr = cpl_array_get_double(parerr, 1, NULL);
721  }
722 #if 0 /* DEBUG */
723  cpl_msg_debug(__func__, "skylevel = %f", cpl_array_get_double(params, 0, NULL));
724  cpl_msg_debug(__func__, "measured flux %f +/- %f", flux, cpl_array_get_double(parerr, 1, NULL));
725 #endif
726  cpl_array_delete(params);
727  cpl_array_delete(parerr);
728  return flux;
729 } /* muse_flux_image_moffat() */
730 
731 /*----------------------------------------------------------------------------*/
750 /*----------------------------------------------------------------------------*/
751 static double
752 muse_flux_image_square(cpl_image *aImage, cpl_image *aImErr, double aX,
753  double aY, double aHalfSize, unsigned int aDSky,
754  unsigned int aMaxBad, float *aFErr)
755 {
756  if (aFErr) { /* set high variance for an error case */
757  *aFErr = FLT_MAX;
758  }
759  int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
760  y1 = aY - aHalfSize, y2 = aY + aHalfSize,
761  nx = cpl_image_get_size_x(aImage),
762  ny = cpl_image_get_size_y(aImage);
763  if (x1 < 1) {
764  x1 = 1;
765  }
766  if (x2 > nx) {
767  x2 = nx;
768  }
769  if (y1 < 1) {
770  y1 = 1;
771  }
772  if (y2 > ny) {
773  y2 = ny;
774  }
775  float skyerror;
776  cpl_errorstate state = cpl_errorstate_get();
777  double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky,
778  &skyerror);
779  if (!cpl_errorstate_is_equal(state)) {
780  /* background determination is critical for this method, *
781  * reset the error but return with zero anyway */
782  cpl_errorstate_set(state);
783  cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
784  return 0.; /* fail on missing background level */
785  }
786 
787  /* extract the measurement region; fail on too many bad pixels, *
788  * but interpolate, if there are only a few bad ones */
789  cpl_image *region = cpl_image_extract(aImage, x1, y1, x2, y2);
790 #if 0 /* DEBUG */
791  cpl_msg_debug(__func__, "region [%d:%d,%d:%d] %"CPL_SIZE_FORMAT" bad pixels",
792  x1, y1, x2, y2, cpl_image_count_rejected(region));
793 #endif
794  if (cpl_image_count_rejected(region) > aMaxBad) {
795  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
796  cpl_image_delete(region);
797  return 0.; /* fail on too many bad pixels */
798  }
799  cpl_image *regerr = cpl_image_extract(aImErr, x1, y1, x2, y2);
800  if (cpl_image_count_rejected(region) > 0) {
801  cpl_detector_interpolate_rejected(region);
802  cpl_detector_interpolate_rejected(regerr);
803  }
804 
805  /* integrated flux, subtracted by the sky over the size of the *
806  * aperture; an error modeled approximately after IRAF phot */
807  int npoints = (x2 - x1 + 1) * (y2 - y1 + 1),
808  /* number of sky pixels should be counted in muse_flux_image_sky(), *
809  * but it's really not so important to add another parameter, the *
810  * error that we compute here is probably too small to be useful... */
811  nsky = 2 * aDSky * (x2 - x1 + y2 - y1 + 2);
812  double flux = cpl_image_get_flux(region) - skylevel * npoints,
813  ferr = sqrt(cpl_image_get_sqflux(regerr)
814  + npoints * skyerror*skyerror * (1. + (double)npoints / nsky));
815 #if 0 /* DEBUG */
816  cpl_msg_debug(__func__, "measured flux %f +/- %f (%d object pixels, %d pixels"
817  " with %f with sky %f)", flux, ferr, npoints, nsky,
818  cpl_image_get_flux(region), cpl_image_get_sqflux(regerr));
819 #endif
820  cpl_image_delete(region);
821  cpl_image_delete(regerr);
822  if (aFErr) {
823  *aFErr = ferr;
824  }
825  return flux;
826 } /* muse_flux_image_square() */
827 
828 /*----------------------------------------------------------------------------*/
848 /*----------------------------------------------------------------------------*/
849 static double
850 muse_flux_image_circle(cpl_image *aImage, cpl_image *aImErr, double aX,
851  double aY, double aAper, double aAnnu, double aDAnnu,
852  unsigned int aMaxBad, float *aFErr)
853 {
854  if (aFErr) { /* set high variance, for error cases */
855  *aFErr = FLT_MAX;
856  }
857  double rmax = ceil(fmax(aAper, aAnnu + aDAnnu));
858  int x1 = aX - rmax, x2 = aX + rmax,
859  y1 = aY - rmax, y2 = aY + rmax,
860  nx = cpl_image_get_size_x(aImage),
861  ny = cpl_image_get_size_y(aImage);
862  if (x1 < 1) {
863  x1 = 1;
864  }
865  if (x2 > nx) {
866  x2 = nx;
867  }
868  if (y1 < 1) {
869  y1 = 1;
870  }
871  if (y2 > ny) {
872  y2 = ny;
873  }
874  /* first loop to collect the background and *
875  * count bad pixels inside the aperture */
876  cpl_vector *vbg = cpl_vector_new((x2 - x1 + 1) * (y2 - y1 + 1)),
877  *vbe = cpl_vector_new((x2 - x1 + 1) * (y2 - y1 + 1));
878  unsigned int nbad = 0, nbg = 0;
879  int i;
880  for (i = x1; i <= x2; i++) {
881  int j;
882  for (j = y1; j <= y2; j++) {
883  double r = sqrt(pow(aX - i, 2) + pow(aY - j, 2));
884  if (r <= aAper) {
885  nbad += cpl_image_is_rejected(aImage, i, j) == 1;
886  }
887  if (r < aAnnu || r > aAnnu + aDAnnu) {
888  continue;
889  }
890  int err;
891  double value = cpl_image_get(aImage, i, j, &err);
892  if (err) { /* exclude bad pixels */
893  continue;
894  }
895  cpl_vector_set(vbg, nbg, value);
896  cpl_vector_set(vbe, nbg, cpl_image_get(aImErr, i, j, &err));
897  nbg++;
898  } /* for j (vertical pixels) */
899  } /* for i (horizontal pixels) */
900  if (nbg <= 0) {
901  cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
902  cpl_vector_delete(vbg);
903  cpl_vector_delete(vbe);
904  return 0.; /* fail on missing background pixels */
905  }
906  cpl_vector_set_size(vbg, nbg);
907  cpl_vector_set_size(vbe, nbg);
908  cpl_matrix *pos = cpl_matrix_new(1, nbg); /* we don't care about positions... */
909  double mse;
910  cpl_polynomial *fit = muse_utils_iterate_fit_polynomial(pos, vbg, vbe, NULL,
911  0, 3., &mse, NULL);
912 #if 0 /* DEBUG */
913  unsigned int nrej = nbg - cpl_vector_get_size(vbg);
914 #endif
915  nbg = cpl_vector_get_size(vbg);
916  cpl_size pows = 0; /* get the zero-order coefficient */
917  double smean = cpl_polynomial_get_coeff(fit, &pows),
918  sstdev = sqrt(mse);
919  cpl_polynomial_delete(fit);
920  cpl_matrix_delete(pos);
921  cpl_vector_delete(vbg);
922  cpl_vector_delete(vbe);
923 #if 0 /* DEBUG */
924  cpl_msg_debug(__func__, "sky: %d pixels (%d rejected), %f +/- %f; found %d "
925  "bad pixels inside aperture", nbg, nrej, smean, sstdev, nbad);
926 #endif
927  if (nbad > aMaxBad) { /* too many bad pixels inside integration area? */
928  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
929  return 0.; /* fail on too many bad pixels */
930  }
931 
932  /* now replace the few bad pixels by interpolation */
933  if (nbad > 0) {
934  cpl_detector_interpolate_rejected(aImage);
935  cpl_detector_interpolate_rejected(aImErr);
936  }
937 
938  /* second loop to integrate the flux */
939  double flux = 0.,
940  ferr = 0.;
941  unsigned int nobj = 0;
942  for (i = x1; i <= x2; i++) {
943  int j;
944  for (j = y1; j <= y2; j++) {
945  double r = sqrt(pow(aX - i, 2) + pow(aY - j, 2));
946  if (r > aAper) {
947  continue;
948  }
949  int err;
950  double value = cpl_image_get(aImage, i, j, &err),
951  error = cpl_image_get(aImErr, i, j, &err);
952  flux += value;
953  ferr += error*error;
954  nobj++;
955  } /* for j (vertical pixels) */
956  } /* for i (horizontal pixels) */
957  flux -= smean * nobj;
958  /* Compute error like IRAF phot: *
959  * error = sqrt (flux / epadu + area * stdev**2 + *
960  * area**2 * stdev**2 / nsky) *
961  * We take our summed error instead of the error computed via the gain */
962  ferr = sqrt(ferr + nobj * sstdev*sstdev * (1. + (double)nobj / nbg));
963 #if 0 /* DEBUG */
964  cpl_msg_debug(__func__, "flux: %d pixels (%d interpolated), %f +/- %f",
965  nobj, nbad, flux, ferr);
966 #endif
967  if (aFErr) {
968  *aFErr = ferr;
969  }
970  return flux;
971 } /* muse_flux_image_circle() */
972 
973 /*----------------------------------------------------------------------------*/
1007 /*----------------------------------------------------------------------------*/
1008 muse_image *
1009 muse_flux_integrate_cube(muse_datacube *aCube, cpl_apertures *aApertures,
1010  muse_flux_profile_type aProfile)
1011 {
1012  cpl_ensure(aCube && aApertures, CPL_ERROR_NULL_INPUT, NULL);
1013  switch (aProfile) {
1015  cpl_msg_info(__func__, "Gaussian profile fits for flux integration");
1016  break;
1018  cpl_msg_info(__func__, "Moffat profile fits for flux integration");
1019  break;
1021  cpl_msg_info(__func__, "Circular flux integration");
1022  break;
1024  cpl_msg_info(__func__, "Simple square window flux integration");
1025  break;
1026  default:
1027  cpl_msg_error(__func__, "Unknown flux integration method!");
1028  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
1029  return NULL;
1030  }
1031 
1032  /* construct image of number of wavelengths x number of stars */
1033  int naper = cpl_apertures_get_size(aApertures), /* can only be > 0 */
1034  nlambda = cpl_imagelist_get_size(aCube->data),
1035  nplane = cpl_imagelist_get_size(aCube->data) / 2; /* central plane */
1036  cpl_image *cim = cpl_imagelist_get(aCube->data, nplane);
1037  muse_image *intimage = muse_image_new();
1038  intimage->data = cpl_image_new(nlambda, naper, CPL_TYPE_FLOAT);
1039  intimage->dq = cpl_image_new(nlambda, naper, CPL_TYPE_INT);
1040  intimage->stat = cpl_image_new(nlambda, naper, CPL_TYPE_FLOAT);
1041  /* copy wavelength WCS from 3rd axis of cube to x-axis of image */
1042  intimage->header = cpl_propertylist_new();
1043  cpl_propertylist_append_double(intimage->header, "CRVAL1",
1044  muse_pfits_get_crval(aCube->header, 3));
1045  cpl_propertylist_append_double(intimage->header, "CRPIX1",
1046  muse_pfits_get_crpix(aCube->header, 3));
1047  cpl_propertylist_append_double(intimage->header, "CD1_1",
1048  muse_pfits_get_cd(aCube->header, 3, 3));
1049  cpl_propertylist_append_string(intimage->header, "CTYPE1",
1050  muse_pfits_get_ctype(aCube->header, 3));
1051  cpl_propertylist_append_string(intimage->header, "CUNIT1",
1052  muse_pfits_get_cunit(aCube->header, 3));
1053  /* fill the 2nd axis with standards */
1054  cpl_propertylist_append_double(intimage->header, "CRVAL2", 1.);
1055  cpl_propertylist_append_double(intimage->header, "CRPIX2", 1.);
1056  cpl_propertylist_append_double(intimage->header, "CD2_2", 1.);
1057  cpl_propertylist_append_string(intimage->header, "CTYPE2", "PIXEL");
1058  cpl_propertylist_append_string(intimage->header, "CUNIT2", "pixel");
1059  cpl_propertylist_append_double(intimage->header, "CD1_2", 0.);
1060  cpl_propertylist_append_double(intimage->header, "CD2_1", 0.);
1061  /* we need the date, data units, exposure time, and instrument mode, too */
1062  cpl_propertylist_append_string(intimage->header, "DATE-OBS",
1063  cpl_propertylist_get_string(aCube->header,
1064  "DATE-OBS"));
1065  cpl_propertylist_append_string(intimage->header, "BUNIT",
1066  muse_pfits_get_bunit(aCube->header));
1067  cpl_propertylist_append_double(intimage->header, "EXPTIME",
1068  muse_pfits_get_exptime(aCube->header));
1069  cpl_propertylist_append_string(intimage->header, "ESO INS MODE",
1070  cpl_propertylist_get_string(aCube->header,
1071  "ESO INS MODE"));
1072  if (cpl_propertylist_has(aCube->header, MUSE_HDR_FLUX_FFCORR)) {
1073  cpl_propertylist_append_bool(intimage->header, MUSE_HDR_FLUX_FFCORR, CPL_TRUE);
1074  cpl_propertylist_set_comment(intimage->header, MUSE_HDR_FLUX_FFCORR,
1075  MUSE_HDR_FLUX_FFCORR_C);
1076  }
1077 
1078  /* get DIMM seeing from the headers, convert it to size in pixels */
1079  cpl_errorstate ps = cpl_errorstate_get();
1080  double fwhm = (muse_pfits_get_fwhm_start(aCube->header)
1081  + muse_pfits_get_fwhm_end(aCube->header)) / 2.;
1082  if (muse_pfits_get_mode(aCube->header) < MUSE_MODE_NFM_AO_N) {
1083  fwhm /= (kMuseSpaxelSizeX_WFM + kMuseSpaxelSizeY_WFM) / 2.;
1084  } else { /* for NFM */
1085  fwhm /= (kMuseSpaxelSizeX_NFM + kMuseSpaxelSizeY_NFM) / 2.;
1086  }
1087  if (!cpl_errorstate_is_equal(ps)) { /* some headers are missing */
1088  double xc = cpl_apertures_get_centroid_x(aApertures, 1),
1089  yc = cpl_apertures_get_centroid_y(aApertures, 1),
1090  xfwhm, yfwhm;
1091  cpl_image_get_fwhm(cim, lround(xc), lround(yc), &xfwhm, &yfwhm);
1092  if (xfwhm > 0. && yfwhm > 0.) {
1093  fwhm = (xfwhm + yfwhm) / 2.;
1094  } else if (xfwhm > 0.) {
1095  fwhm = xfwhm;
1096  } else if (yfwhm > 0.) {
1097  fwhm = yfwhm;
1098  } else {
1099  fwhm = 5.; /* total failure to measure it, assume 1 arcsec seeing */
1100  }
1101  cpl_errorstate_set(ps);
1102  cpl_msg_debug(__func__, "Using roughly estimated reference FWHM (%.3f pix) "
1103  "instead of DIMM seeing", fwhm);
1104  } else {
1105  cpl_msg_debug(__func__, "Using DIMM seeing of %.3f pix for reference FWHM",
1106  fwhm);
1107  }
1108 
1109  /* track the sizes used for the fits/integrations in an image, *
1110  * record the half-sizes (or radiuses) in on row per object */
1111  cpl_image *sizes = cpl_image_new(nlambda, naper, CPL_TYPE_DOUBLE);
1112  double *psizes = cpl_image_get_data_double(sizes);
1113  /* access pointers for the flux-image */
1114  float *data = cpl_image_get_data_float(intimage->data),
1115  *stat = cpl_image_get_data_float(intimage->stat);
1116  int *dq = cpl_image_get_data_int(intimage->dq);
1117  /* loop over all wavelengths and measure the flux */
1118  int l, ngood = 0, nillegal = 0, nbadbg = 0; /* count good fits and errors */
1119  #pragma omp parallel for default(none) /* as req. by Ralf */ \
1120  shared(aApertures, aCube, aProfile, data, dq, fwhm, naper, nbadbg, \
1121  ngood, nillegal, nlambda, psizes, stat)
1122  for (l = 0; l < nlambda; l++) {
1123  cpl_image *plane = cpl_imagelist_get(aCube->data, l),
1124  *pldq = aCube->dq ? cpl_imagelist_get(aCube->dq, l) : NULL,
1125  *plerr = cpl_image_duplicate(cpl_imagelist_get(aCube->stat, l));
1126 #if 0 /* DEBUG */
1127  cpl_stats *stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1128  cpl_msg_debug(__func__, "lambda = %d/%f %s", l + 1,
1129  (l + 1 - muse_pfits_get_crpix(aCube->header, 3))
1130  * muse_pfits_get_cd(aCube->header, 3, 3)
1131  + muse_pfits_get_crval(aCube->header, 3),
1132  muse_pfits_get_cunit(aCube->header, 3));
1133  cpl_msg_debug(__func__, "variance: %g...%g...%g", cpl_stats_get_min(stats),
1134  cpl_stats_get_mean(stats), cpl_stats_get_max(stats));
1135  cpl_stats_delete(stats);
1136 #endif
1137  /* make sure to exclude bad pixels from the fits below */
1138  if (pldq) {
1139  muse_quality_image_reject_using_dq(plane, pldq, plerr);
1140  } else {
1141  cpl_image_reject_value(plane, CPL_VALUE_NAN); /* mark NANs */
1142  cpl_image_reject_value(plerr, CPL_VALUE_NAN);
1143  }
1144 #if 0 /* DEBUG */
1145  stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1146  cpl_msg_debug(__func__, "cut variance: %g...%g...%g (%"CPL_SIZE_FORMAT" bad"
1147  " pixel)", cpl_stats_get_min(stats), cpl_stats_get_mean(stats),
1148  cpl_stats_get_max(stats), cpl_image_count_rejected(plane));
1149  cpl_stats_delete(stats);
1150 #endif
1151  /* convert variable to sigmas */
1152  cpl_image_power(plerr, 0.5);
1153 #if 0 /* DEBUG */
1154  stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1155  cpl_msg_debug(__func__, "errors: %g...%g...%g", cpl_stats_get_min(stats),
1156  cpl_stats_get_mean(stats), cpl_stats_get_max(stats));
1157  cpl_stats_delete(stats);
1158 #endif
1159  cpl_errorstate state = cpl_errorstate_get();
1160  int n;
1161  for (n = 1; n <= naper; n++) {
1162  /* use detection aperture to construct much larger one for measurement */
1163  double xc = cpl_apertures_get_centroid_x(aApertures, n),
1164  yc = cpl_apertures_get_centroid_y(aApertures, n),
1165  size = sqrt(cpl_apertures_get_npix(aApertures, n)),
1166  xfwhm, yfwhm;
1167  cpl_errorstate prestate = cpl_errorstate_get();
1168  cpl_image_get_fwhm(plane, lround(xc), lround(yc), &xfwhm, &yfwhm);
1169  if (xfwhm < 0 || yfwhm < 0) {
1170  data[l + (n-1) * nlambda] = 0.;
1171  stat[l + (n-1) * nlambda] = FLT_MAX;
1172  cpl_errorstate_set(prestate);
1173  continue;
1174  }
1175  /* half size for flux integration, at least 3 x FWHM */
1176  double halfsize = fmax(1.5 * (xfwhm + yfwhm), 3. * fwhm);
1177  if (halfsize < size / 2) { /* at least the size of the det. aperture */
1178  halfsize = size / 2;
1179  }
1180  psizes[l + (n-1) * nlambda] = halfsize;
1181 #if 0 /* DEBUG */
1182  cpl_msg_debug(__func__, "%.2f,%.2f FWHM %.2f %.2f size %.2f --> %.2f",
1183  xc, yc, xfwhm, yfwhm, size, halfsize * 2.);
1184 #endif
1185 
1186  switch (aProfile) {
1188  data[l + (n-1) * nlambda] = muse_flux_image_gaussian(plane, plerr, xc, yc,
1189  halfsize, 5, 10,
1190  &stat[l + (n-1) * nlambda]);
1191  break;
1193  data[l + (n-1) * nlambda] = muse_flux_image_moffat(plane, plerr, xc, yc,
1194  halfsize, 5, 10,
1195  &stat[l + (n-1) * nlambda]);
1196  break;
1197  case MUSE_FLUX_PROFILE_CIRCLE: {
1198  /* the circular method needs larger region to properly integrate *
1199  * everything at least something like 4 x FWHM to be sure */
1200  double radius = 4./3. * halfsize,
1201  rannu = radius * 5. / 4.; /* background annulus a bit larger */
1202  psizes[l + (n-1) * nlambda] = radius;
1203  data[l + (n-1) * nlambda] = muse_flux_image_circle(plane, plerr, xc, yc,
1204  radius, rannu, 10, 10,
1205  &stat[l + (n-1) * nlambda]);
1206  break;
1207  } /* case MUSE_FLUX_PROFILE_CIRCLE */
1208  default: /* MUSE_FLUX_PROFILE_EQUAL_SQUARE */
1209  data[l + (n-1) * nlambda] = muse_flux_image_square(plane, plerr, xc, yc,
1210  halfsize, 5, 10,
1211  &stat[l + (n-1) * nlambda]);
1212  } /* switch */
1213  if (data[l + (n-1) * nlambda] < 0 || !isfinite(data[l + (n-1) * nlambda])) {
1214  data[l + (n-1) * nlambda] = 0.; /* should not contribute to flux */
1215  dq[l + (n-1) * nlambda] = EURO3D_MISSDATA; /* mark as bad in DQ extension */
1216  stat[l + (n-1) * nlambda] = FLT_MAX;
1217  }
1218  } /* for n (all apertures) */
1219 
1220  /* count "Illegal input" errors and for those reset the state, as there *
1221  * can be many of them (one for each incompletely filled wavelength plane */
1222  if (!cpl_errorstate_is_equal(state)) {
1223  if (cpl_error_get_code() == CPL_ERROR_ILLEGAL_INPUT) {
1224  cpl_errorstate_set(state);
1225  #pragma omp atomic
1226  nillegal++;
1227  } else if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
1228  cpl_errorstate_set(state);
1229  #pragma omp atomic
1230  nbadbg++;
1231  }
1232  } else {
1233  #pragma omp atomic
1234  ngood++;
1235  }
1236 
1237  cpl_image_delete(plerr);
1238  } /* for l (all wavelengths) */
1239 
1240  /* output statistics for the sizes used, mask out unset positions in the image */
1241  cpl_image_reject_value(sizes, CPL_VALUE_ZERO);
1242  int n;
1243  for (n = 1; n <= naper; n++) {
1244  if (aProfile == MUSE_FLUX_PROFILE_CIRCLE) {
1245  cpl_msg_info(__func__, "Radiuses used for circular flux integration for "
1246  "object %d: %f +/- %f (%f) %f..%f", n,
1247  cpl_image_get_mean_window(sizes, 1, n, nlambda, n),
1248  cpl_image_get_stdev_window(sizes, 1, n, nlambda, n),
1249  cpl_image_get_median_window(sizes, 1, n, nlambda, n),
1250  cpl_image_get_min_window(sizes, 1, n, nlambda, n),
1251  cpl_image_get_max_window(sizes, 1, n, nlambda, n));
1252  } else {
1253  cpl_msg_info(__func__, "Half-sizes used for flux integration for object "
1254  "%d: %f +/- %f (%f) %f..%f", n,
1255  cpl_image_get_mean_window(sizes, 1, n, nlambda, n),
1256  cpl_image_get_stdev_window(sizes, 1, n, nlambda, n),
1257  cpl_image_get_median_window(sizes, 1, n, nlambda, n),
1258  cpl_image_get_min_window(sizes, 1, n, nlambda, n),
1259  cpl_image_get_max_window(sizes, 1, n, nlambda, n));
1260  } /* else */
1261  } /* for n (all apertures) */
1262 #if 0
1263  cpl_image_save(sizes, "sizes.fits", CPL_TYPE_UNSPECIFIED, NULL, CPL_IO_CREATE);
1264 #endif
1265  cpl_image_delete(sizes);
1266 
1267  /* add headers about the objects to the integrated image */
1268  cpl_propertylist_append_int(intimage->header, MUSE_HDR_FLUX_NOBJ, naper);
1269  cpl_propertylist_set_comment(intimage->header, MUSE_HDR_FLUX_NOBJ,
1270  MUSE_HDR_FLUX_NOBJ_C);
1271  /* create a basic WCS, assuming nominal MUSE properties, for an *
1272  * estimate of the celestial position of all detected objects */
1273  cpl_propertylist *wcs1 = muse_wcs_create_default(),
1274  *wcs = muse_wcs_apply_cd(aCube->header, wcs1);
1275  cpl_propertylist_delete(wcs1);
1276  /* update WCS just as in muse_resampling_cube() */
1277  double crpix1 = muse_pfits_get_crpix(aCube->header, 1)
1278  + (1. + cpl_image_get_size_x(cim)) / 2.,
1279  crpix2 = muse_pfits_get_crpix(aCube->header, 2)
1280  + (1. + cpl_image_get_size_y(cim)) / 2.;
1281  cpl_propertylist_update_double(wcs, "CRPIX1", crpix1);
1282  cpl_propertylist_update_double(wcs, "CRPIX2", crpix2);
1283  cpl_propertylist_update_double(wcs, "CRVAL1", muse_pfits_get_ra(aCube->header));
1284  cpl_propertylist_update_double(wcs, "CRVAL2", muse_pfits_get_dec(aCube->header));
1285  for (n = 1; n <= naper; n++) {
1286  /* use detection aperture to construct much larger one for measurement */
1287  double xc = cpl_apertures_get_centroid_x(aApertures, n),
1288  yc = cpl_apertures_get_centroid_y(aApertures, n),
1289  ra, dec;
1290  muse_wcs_celestial_from_pixel(wcs, xc, yc, &ra, &dec);
1291  double flux = cpl_image_get_flux_window(intimage->data, 1, n, nlambda, n);
1292  cpl_msg_debug(__func__, "Object %02d: %.3f,%.3f pix, %f,%f deg, flux %e %s",
1293  n, xc, yc, ra, dec, flux, muse_pfits_get_bunit(intimage->header));
1294  char kw[KEYWORD_LENGTH];
1295  snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_X, n);
1296  cpl_propertylist_append_float(intimage->header, kw, xc);
1297  cpl_propertylist_set_comment(intimage->header, kw, MUSE_HDR_FLUX_OBJn_X_C);
1298  snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_Y, n);
1299  cpl_propertylist_append_float(intimage->header, kw, yc);
1300  cpl_propertylist_set_comment(intimage->header, kw, MUSE_HDR_FLUX_OBJn_Y_C);
1301  snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_RA, n);
1302  cpl_propertylist_append_double(intimage->header, kw, ra);
1303  cpl_propertylist_set_comment(intimage->header, kw, MUSE_HDR_FLUX_OBJn_RA_C);
1304  snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_DEC, n);
1305  cpl_propertylist_append_double(intimage->header, kw, dec);
1306  cpl_propertylist_set_comment(intimage->header, kw, MUSE_HDR_FLUX_OBJn_DEC_C);
1307  snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_FLUX, n);
1308  cpl_propertylist_append_double(intimage->header, kw, flux);
1309  cpl_propertylist_set_comment(intimage->header, kw, MUSE_HDR_FLUX_OBJn_FLUX_C);
1310  } /* for n (all apertures) */
1311  cpl_propertylist_delete(wcs);
1312 
1313  if (nillegal > 0 || nbadbg > 0) {
1314  cpl_msg_warning(__func__, "Successful fits in %d wavelength planes, but "
1315  "encountered %d \"Illegal input\" errors and %d bad "
1316  "background determinations", ngood, nillegal, nbadbg);
1317  } else {
1318  cpl_msg_info(__func__, "Successful fits in %d wavelength planes", ngood);
1319  }
1320 
1321  return intimage;
1322 } /* muse_flux_integrate_cube() */
1323 
1324 /*----------------------------------------------------------------------------*/
1360 /*----------------------------------------------------------------------------*/
1361 cpl_error_code
1363  muse_flux_object *aFluxObj)
1364 {
1365  cpl_ensure_code(aPixtable && aFluxObj, CPL_ERROR_NULL_INPUT);
1366  switch (aProfile) {
1371  break;
1372  default:
1373  return cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
1374  }
1375 
1376  /* if remove notch filter flags from NaD range, since we do want a *
1377  * (near-zero) throughput measurement inside that filtered range; *
1378  * but keep original DQ column around to reset it later */
1379  muse_ins_mode insmode = muse_pfits_get_mode(aPixtable->header);
1380  if (insmode >= MUSE_MODE_WFM_AO_E) {
1381  cpl_table_duplicate_column(aPixtable->table, MUSE_PIXTABLE_DQ"_NAD",
1382  aPixtable->table, MUSE_PIXTABLE_DQ);
1383  int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
1384  cpl_size irow, nrow = muse_pixtable_get_nrow(aPixtable);
1385  for (irow = 0; irow < nrow; irow++) {
1386  if (dq[irow] & EURO3D_NOTCH_NAD) {
1387  dq[irow] &= ~EURO3D_NOTCH_NAD; /* reset that bit */
1388  } /* if */
1389  } /* for */
1390  } /* if AO mode */
1391 
1392  if (getenv("MUSE_DEBUG_FLUX") && atoi(getenv("MUSE_DEBUG_FLUX")) > 2) {
1393  const char *fn = "flux__pixtable.fits";
1394  cpl_msg_info(__func__, "Saving pixel table as \"%s\"", fn);
1395  muse_pixtable_save(aPixtable, fn);
1396  }
1397  muse_resampling_params *params =
1399  params->pfx = 1.; /* large pixfrac to be sure to cover most gaps */
1400  params->pfy = 1.;
1401  params->pfl = 1.;
1402  /* resample at nominal resolution, the conversion to [1/Angstrom] *
1403  * is later done when computing the sensitivity function */
1404  params->dlambda = kMuseSpectralSamplingA;
1406  params->crsigma = 25.;
1407  muse_datacube *cube = muse_resampling_cube(aPixtable, params, NULL);
1408  if (cube) {
1409  aFluxObj->cube = cube;
1410  }
1412  /* reinstate the original DQ column */
1413  if (cpl_table_has_column(aPixtable->table, MUSE_PIXTABLE_DQ"_NAD")) {
1414  cpl_table_erase_column(aPixtable->table, MUSE_PIXTABLE_DQ);
1415  cpl_table_name_column(aPixtable->table, MUSE_PIXTABLE_DQ"_NAD",
1417  } /* if */
1418  if (getenv("MUSE_DEBUG_FLUX") && atoi(getenv("MUSE_DEBUG_FLUX")) >= 2) {
1419  const char *fn = "flux__cube.fits";
1420  cpl_msg_info(__func__, "Saving cube as \"%s\"", fn);
1421  muse_datacube_save(aFluxObj->cube, fn);
1422  }
1423  int nplane = cpl_imagelist_get_size(cube->data) / 2; /* central plane */
1424  cpl_image *cim = cpl_imagelist_get(cube->data, nplane),
1425  *cdq = cpl_imagelist_get(cube->dq, nplane);
1426  /* make sure to reject bad pixels before using *
1427  * object detection using image statistics */
1428  muse_quality_image_reject_using_dq(cim, cdq, NULL);
1429  /* use high sigmas for detection */
1430  double dsigmas[] = { 50., 30., 10., 8., 6., 5. };
1431  cpl_vector *vsigmas = cpl_vector_wrap(sizeof(dsigmas) / sizeof(double),
1432  dsigmas);
1433  cpl_size isigma = -1;
1434  cpl_apertures *apertures = cpl_apertures_extract(cim, vsigmas, &isigma);
1435  int napertures = apertures ? cpl_apertures_get_size(apertures) : 0;
1436  if (napertures < 1) {
1437  /* isigma is still -1 in this case, so take the last vector entry */
1438  cpl_msg_error(__func__, "No object for flux integration found down to %.1f"
1439  " sigma limit on plane %d",
1440  cpl_vector_get(vsigmas, cpl_vector_get_size(vsigmas) - 1),
1441  nplane + 1);
1442  cpl_vector_unwrap(vsigmas);
1443  cpl_apertures_delete(apertures);
1444  return cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
1445  }
1446  cpl_msg_debug(__func__, "The %.1f sigma threshold was used to find %d object%s"
1447  " on plane %d", cpl_vector_get(vsigmas, isigma), napertures,
1448  napertures == 1 ? "" : "s", nplane + 1);
1449  cpl_vector_unwrap(vsigmas);
1450 #if 0 /* DEBUG */
1451  cpl_apertures_dump(apertures, stdout);
1452  fflush(stdout);
1453 #endif
1454 
1455  /* now do the flux integration */
1456  muse_image *intimage = muse_flux_integrate_cube(cube, apertures, aProfile);
1457  cpl_apertures_delete(apertures);
1458  aFluxObj->intimage = intimage; /* save integrated fluxes into in/out struct */
1459 
1460  return CPL_ERROR_NONE;
1461 } /* muse_flux_integrate_std() */
1462 
1463 /* prominent telluric absorption bands, together with clean regions: *
1464  * [0] lower limit of telluric region *
1465  * [1] upper limit of telluric region *
1466  * [2] lower limit of fit region (excludes telluric region) *
1467  * [3] upper limit of fit region (excludes telluric region) *
1468  * created from by-hand measurements on spectra from the ESO sky *
1469  * model, and from the Keck list of telluric lines */
1470 static const double kTelluricBands[][4] = {
1471  { 6273., 6320., 6213., 6380. }, /* now +/-60 Angstrom around... */
1472  { 6864., 6967., 6750., 7130. }, /* B-band */
1473  { 7164., 7325., 7070., 7580. },
1474  { 7590., 7700., 7470., 7830. }, /* A-band */
1475  { 8131., 8345., 7900., 8600. },
1476  { 8952., 9028., 8850., 9082. }, /* very rough! */
1477  { 9274., 9770., 9080., 9263. }, /* very very rough! */
1478  { -1., -1., -1., -1. }
1479 };
1480 
1481 /*----------------------------------------------------------------------------*/
1485 /*----------------------------------------------------------------------------*/
1487  { "lmin", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
1488  "lower limit of the telluric region", CPL_TRUE },
1489  { "lmax", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
1490  "upper limit of the telluric region", CPL_TRUE },
1491  { "bgmin", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
1492  "lower limit of the background region", CPL_TRUE },
1493  { "bgmax", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
1494  "upper limit of the background region", CPL_TRUE },
1495  { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
1496 };
1497 
1498 /*----------------------------------------------------------------------------*/
1508 /*----------------------------------------------------------------------------*/
1509 static cpl_error_code
1510 muse_flux_response_set_telluric_bands(muse_flux_object *aFluxObj,
1511  const cpl_table *aTellBands)
1512 {
1513  cpl_ensure_code(aFluxObj, CPL_ERROR_NULL_INPUT);
1514 
1515  /* if a valid table was passed, just set that */
1516  if (aTellBands && muse_cpltable_check(aTellBands, muse_response_tellbands_def)
1517  == CPL_ERROR_NONE) {
1518  cpl_msg_debug(__func__, "using given table for telluric bands");
1519  aFluxObj->tellbands = cpl_table_duplicate(aTellBands);
1520  return CPL_ERROR_NONE;
1521  }
1522  /* create a table for the default regions */
1523  unsigned int ntell = sizeof(kTelluricBands) / sizeof(kTelluricBands[0]) - 1;
1524  cpl_msg_debug(__func__, "using builtin regions for telluric bands (%u "
1525  "entries)", ntell);
1526  aFluxObj->tellbands = muse_cpltable_new(muse_response_tellbands_def, ntell);
1527  cpl_table *tb = aFluxObj->tellbands; /* shortcut */
1528  unsigned int k;
1529  for (k = 0; k < ntell; k++) {
1530  cpl_table_set_double(tb, "lmin", k, kTelluricBands[k][0]);
1531  cpl_table_set_double(tb, "lmax", k, kTelluricBands[k][1]);
1532  cpl_table_set_double(tb, "bgmin", k, kTelluricBands[k][2]);
1533  cpl_table_set_double(tb, "bgmax", k, kTelluricBands[k][3]);
1534  } /* for k */
1535  if (getenv("MUSE_DEBUG_FLUX") && atoi(getenv("MUSE_DEBUG_FLUX")) >= 2) {
1536  const char *fn = "flux__tellregions.fits";
1537  cpl_msg_info(__func__, "Saving telluric bands table as \"%s\"", fn);
1538  cpl_table_save(tb, NULL, NULL, fn, CPL_IO_CREATE);
1539  }
1540  return CPL_ERROR_NONE;
1541 } /* muse_flux_response_set_telluric_bands() */
1542 
1543 /*----------------------------------------------------------------------------*/
1556 /*----------------------------------------------------------------------------*/
1557 static void
1558 muse_flux_response_dump_sensitivity(muse_flux_object *aFluxObj,
1559  const char *aName)
1560 {
1561  char *dodebug = getenv("MUSE_DEBUG_FLUX");
1562  if (!dodebug || (dodebug && atoi(dodebug) <= 0)) {
1563  return;
1564  }
1565  char *fn = cpl_sprintf("flux__sens_%s.ascii", aName);
1566  FILE *fp = fopen(fn, "w");
1567  fprintf(fp, "#"); /* prefix first line (table header) for easier plotting */
1568  cpl_table_dump(aFluxObj->sensitivity, 0,
1569  cpl_table_get_nrow(aFluxObj->sensitivity), fp);
1570  fclose(fp);
1571  cpl_msg_debug(__func__, "Written %"CPL_SIZE_FORMAT" datapoints to \"%s\"",
1572  cpl_table_get_nrow(aFluxObj->sensitivity), fn);
1573  cpl_free(fn);
1574 } /* muse_flux_response_dump_sensitivity() */
1575 
1576 /*----------------------------------------------------------------------------*/
1597 /*----------------------------------------------------------------------------*/
1598 static cpl_error_code
1599 muse_flux_response_sensitivity(muse_flux_object *aFluxObj,
1600  unsigned int aStar, const cpl_table *aReference,
1601  double aAirmass, const cpl_table *aExtinct)
1602 {
1603  cpl_ensure_code(aFluxObj && aFluxObj->intimage, CPL_ERROR_NULL_INPUT);
1604 
1605  double crval = muse_pfits_get_crval(aFluxObj->intimage->header, 1),
1606  cdelt = muse_pfits_get_cd(aFluxObj->intimage->header, 1, 1),
1607  crpix = muse_pfits_get_crpix(aFluxObj->intimage->header, 1),
1608  exptime = muse_pfits_get_exptime(aFluxObj->intimage->header);
1609  int nlambda = cpl_image_get_size_x(aFluxObj->intimage->data);
1610 
1611  aFluxObj->sensitivity = cpl_table_new(nlambda);
1612  cpl_table *sensitivity = aFluxObj->sensitivity;
1613  cpl_table_new_column(sensitivity, "lambda", CPL_TYPE_DOUBLE);
1614  cpl_table_new_column(sensitivity, "sens", CPL_TYPE_DOUBLE);
1615  cpl_table_new_column(sensitivity, "serr", CPL_TYPE_DOUBLE);
1616  cpl_table_new_column(sensitivity, "dq", CPL_TYPE_INT);
1617  cpl_table_set_column_format(sensitivity, "dq", "%u");
1618  float *data = cpl_image_get_data_float(aFluxObj->intimage->data),
1619  *stat = cpl_image_get_data_float(aFluxObj->intimage->stat);
1620  int l, idx = 0;
1621  for (l = 0; l < nlambda; l++) {
1622  if (data[l + aStar*nlambda] <= 0. ||
1623  stat[l + aStar*nlambda] <= 0. ||
1624  stat[l + aStar*nlambda] == FLT_MAX) { /* exclude bad fits */
1625  continue;
1626  }
1627 
1628  double lambda = crval + cdelt * (l + 1 - crpix),
1629  /* interpolate extinction curve at this wavelength */
1630  extinct = !aExtinct ? 0. /* no extinction term */
1631  : muse_flux_response_interpolate(aExtinct, lambda, NULL,
1633  cpl_errorstate prestate = cpl_errorstate_get();
1634  double referr = 0.,
1635  ref = muse_flux_response_interpolate(aReference, lambda, &referr,
1637  /* on error, try again without trying to find the fluxerr column */
1638  if (!cpl_errorstate_is_equal(prestate)) {
1639  cpl_errorstate_set(prestate); /* don't want to propagate this error outside */
1640  ref = muse_flux_response_interpolate(aReference, lambda, NULL,
1642  }
1643  /* calibration factor at this wavelength: ratio of observed count *
1644  * rate corrected for extinction to expected flux in magnitudes */
1645  double c = 2.5 * log10(data[l + aStar*nlambda]
1646  / exptime / cdelt / ref)
1647  + aAirmass * extinct,
1648  cerr = sqrt(pow(referr / ref, 2) + stat[l + aStar*nlambda]
1649  / pow(data[l + aStar*nlambda], 2))
1650  * 2.5 / CPL_MATH_LN10;
1651  cpl_table_set_double(sensitivity, "lambda", idx, lambda);
1652  cpl_table_set_double(sensitivity, "sens", idx, c);
1653  cpl_table_set_double(sensitivity, "serr", idx, cerr);
1654  cpl_table_set_int(sensitivity, "dq", idx, EURO3D_GOODPIXEL);
1655  idx++;
1656  } /* for l (all wavelengths) */
1657  /* cut the data to the used size */
1658  cpl_table_set_size(sensitivity, idx);
1659  return CPL_ERROR_NONE;
1660 } /* muse_flux_response_sensitivity() */
1661 
1662 /*----------------------------------------------------------------------------*/
1678 /*----------------------------------------------------------------------------*/
1679 static cpl_error_code
1680 muse_flux_response_mark_questionable(muse_flux_object *aFluxObj)
1681 {
1682  cpl_ensure_code(aFluxObj && aFluxObj->intimage, CPL_ERROR_NULL_INPUT);
1683 
1684  cpl_propertylist *header = aFluxObj->intimage->header; /* shortcut */
1685  cpl_table *tsens = aFluxObj->sensitivity,
1686  *tb = aFluxObj->tellbands;
1687 
1688  /* check for MUSE setup: is it nominal wavelength range? */
1689  cpl_boolean isnominal = muse_pfits_get_mode(header) == MUSE_MODE_WFM_NONAO_N,
1690  isAON = (muse_pfits_get_mode(header) == MUSE_MODE_WFM_AO_N)
1691  || (muse_pfits_get_mode(header) == MUSE_MODE_NFM_AO_N);
1692  /* exclude regions within the telluric absorption bands *
1693  * and outside wavelength range */
1694  int irow, nfluxes = cpl_table_get_nrow(tsens);
1695  for (irow = 0; irow < nfluxes; irow++) {
1696  double lambda = cpl_table_get_double(tsens, "lambda", irow, NULL);
1697  unsigned int dq = EURO3D_GOODPIXEL,
1698  k, nk = cpl_table_get_nrow(tb);
1699  for (k = 0; k < nk; k++) {
1700  double lmin = cpl_table_get_double(tb, "lmin", k, NULL),
1701  lmax = cpl_table_get_double(tb, "lmax", k, NULL);
1702  if (lambda >= lmin && lambda <= lmax) {
1703  dq |= EURO3D_TELLURIC;
1704  }
1705  } /* for k */
1706  if ((isnominal && lambda < kMuseNominalCutoff) ||
1707  (isAON && lambda < kMuseAOCutoff)) {
1708  dq |= EURO3D_OUTSDRANGE;
1709  }
1710  cpl_table_set_int(tsens, "dq", irow, dq);
1711  } /* for irow (all table) */
1712  return CPL_ERROR_NONE;
1713 } /* muse_flux_response_mark_questionable() */
1714 
1715 /*----------------------------------------------------------------------------*/
1739 /*----------------------------------------------------------------------------*/
1740 static cpl_polynomial *
1741 muse_flux_response_fit(muse_flux_object *aFluxObj,
1742  double aLambda1, double aLambda2,
1743  unsigned int aOrder, double aRSigma, double *aRMSE)
1744 {
1745  cpl_ensure(aFluxObj && aFluxObj->sensitivity, CPL_ERROR_NULL_INPUT, NULL);
1746 
1747  cpl_table *tsens = aFluxObj->sensitivity;
1748  cpl_table_select_all(tsens); /* default, but make sure it's true...*/
1749  cpl_table_and_selected_int(tsens, "dq", CPL_NOT_EQUAL_TO, EURO3D_GOODPIXEL);
1750  cpl_table_and_selected_int(tsens, "dq", CPL_NOT_EQUAL_TO, EURO3D_TELLCOR);
1751  cpl_table_or_selected_double(tsens, "lambda", CPL_LESS_THAN, aLambda1);
1752  cpl_table_or_selected_double(tsens, "lambda", CPL_GREATER_THAN, aLambda2);
1753  /* keep the "bad" ones around */
1754  cpl_table *tunwanted = cpl_table_extract_selected(tsens);
1755  cpl_table_erase_selected(tsens);
1756  muse_flux_response_dump_sensitivity(aFluxObj, "fitinput");
1757 
1758  /* convert sensitivity table to matrix (lambda) and vectors (sens *
1759  * and serr), exclude pixels marked as not EURO3D_GOODPIXEL */
1760  int nrow = cpl_table_get_nrow(tsens);
1761  cpl_matrix *lambdas = cpl_matrix_new(1, nrow);
1762  cpl_vector *sens = cpl_vector_new(nrow),
1763  *serr = cpl_vector_new(nrow);
1764  memcpy(cpl_matrix_get_data(lambdas),
1765  cpl_table_get_data_double_const(tsens, "lambda"), nrow*sizeof(double));
1766  memcpy(cpl_vector_get_data(sens),
1767  cpl_table_get_data_double_const(tsens, "sens"), nrow*sizeof(double));
1768  memcpy(cpl_vector_get_data(serr),
1769  cpl_table_get_data_double_const(tsens, "serr"), nrow*sizeof(double));
1770 
1771  /* do the fit */
1772  double chisq, mse;
1773  cpl_polynomial *fit = muse_utils_iterate_fit_polynomial(lambdas, sens, serr,
1774  tsens, aOrder, aRSigma,
1775  &mse, &chisq);
1776  int nout = cpl_vector_get_size(sens);
1777 #if 0
1778  cpl_msg_debug(__func__, "transferred %d entries (%.3f...%.3f) for the "
1779  "order %u fit, %d entries are left, RMS %f", nrow, aLambda1,
1780  aLambda2, aOrder, nout, sqrt(mse));
1781 #endif
1782  cpl_matrix_delete(lambdas);
1783  cpl_vector_delete(sens);
1784  cpl_vector_delete(serr);
1785  if (aRMSE) {
1786  *aRMSE = mse / (nout - aOrder - 1);
1787  }
1788 
1789  /* put "bad" entries back, at the end of the table */
1790  cpl_table_insert(tsens, tunwanted, nout);
1791  cpl_table_delete(tunwanted);
1792  return fit;
1793 } /* muse_flux_response_fit() */
1794 
1795 /*----------------------------------------------------------------------------*/
1815 /*----------------------------------------------------------------------------*/
1816 static cpl_error_code
1817 muse_flux_response_telluric(muse_flux_object *aFluxObj, double aAirmass)
1818 {
1819  cpl_ensure_code(aFluxObj, CPL_ERROR_NULL_INPUT);
1820 
1821  cpl_table *tsens = aFluxObj->sensitivity,
1822  *tb = aFluxObj->tellbands;
1823  cpl_table_new_column(tsens, "sens_orig", CPL_TYPE_DOUBLE);
1824  cpl_table_new_column(tsens, "serr_orig", CPL_TYPE_DOUBLE);
1825  cpl_table_new_column(tsens, "tellcor", CPL_TYPE_DOUBLE);
1826  unsigned int k, nk = cpl_table_get_nrow(tb);
1827  for (k = 0; k < nk; k++) {
1828  double lmin = cpl_table_get_double(tb, "lmin", k, NULL),
1829  lmax = cpl_table_get_double(tb, "lmax", k, NULL),
1830  bgmin = cpl_table_get_double(tb, "bgmin", k, NULL),
1831  bgmax = cpl_table_get_double(tb, "bgmax", k, NULL),
1832  datamin = cpl_table_get_column_min(tsens, "lambda"),
1833  datamax = cpl_table_get_column_max(tsens, "lambda");
1834  cpl_boolean extrapolate = CPL_FALSE;
1835  if (bgmax < lmax || datamax < lmax) {
1836  extrapolate = CPL_TRUE;
1837  }
1838  if (bgmin > lmin || datamin > lmin) {
1839  extrapolate = CPL_TRUE;
1840  }
1841  if (datamin > lmax || datamax < lmin) {
1842  /* no sense trying to even extrapolate linearly */
1843  cpl_msg_warning(__func__, "Telluric region %u (range %.2f...%.2f, "
1844  "reference region %.2f...%.2f) outside data range "
1845  "(%.2f..%.2f)!", k + 1, lmin, lmax, bgmin, bgmax,
1846  datamin, datamax);
1847  continue;
1848  }
1849  /* if we extrapolate (often redward) then use a linear fit only */
1850  unsigned int order = extrapolate ? 1 : 2;
1851  /* the telluric regions themselves are already marked, so *
1852  * they don't need to be touched again there before the fit */
1853  double rmse = 0.;
1854  cpl_errorstate state = cpl_errorstate_get();
1855  cpl_polynomial *fit = muse_flux_response_fit(aFluxObj, bgmin, bgmax,
1856  order, 3., &rmse);
1857  if (!cpl_errorstate_is_equal(state) || !fit) {
1858  cpl_msg_warning(__func__, "Telluric region %u (range %.2f...%.2f, "
1859  "reference region %.2f...%.2f) could not be fitted!",
1860  k + 1, lmin, lmax, bgmin, bgmax);
1861  cpl_errorstate_set(state); /* swallow the errors */
1862  cpl_polynomial_delete(fit); /* just in case the polynomial was created */
1863  continue;
1864  }
1865  cpl_msg_debug(__func__, "Telluric region %u: %.2f...%.2f, reference region "
1866  "%.2f...%.2f", k + 1, lmin, lmax, bgmin, bgmax);
1867 #if 0
1868  cpl_polynomial_dump(fit, stdout);
1869  fflush(stdout);
1870 #endif
1871  int irow, nrow = cpl_table_get_nrow(tsens);
1872  for (irow = 0; irow < nrow; irow++) {
1873  double lambda = cpl_table_get_double(tsens, "lambda", irow, NULL);
1874  if (lambda >= lmin && lambda <= lmax &&
1875  (unsigned int)cpl_table_get_int(tsens, "dq", irow, NULL)
1876  == EURO3D_TELLURIC) {
1877  double origval = cpl_table_get_double(tsens, "sens", irow, NULL),
1878  origerr = cpl_table_get_double(tsens, "serr", irow, NULL),
1879  interpval = cpl_polynomial_eval_1d(fit, lambda, NULL);
1880  cpl_table_set_int(tsens, "dq", irow, EURO3D_TELLCOR);
1881  cpl_table_set_double(tsens, "sens_orig", irow, origval);
1882  cpl_table_set_double(tsens, "sens", irow, interpval);
1883  /* correct the error bars of the fitted points *
1884  * by adding in quadrature the reduced MSE */
1885  cpl_table_set_double(tsens, "serr_orig", irow, origerr);
1886  cpl_table_set_double(tsens, "serr", irow, sqrt(origerr*origerr + rmse));
1887 
1888  if (interpval > origval) {
1889  /* compute the factor, as flux ratio */
1890  double ftelluric = pow(10, -0.4 * (interpval - origval));
1891  cpl_table_set(tsens, "tellcor", irow, ftelluric);
1892  } else {
1893  cpl_table_set_double(tsens, "tellcor", irow, 1.);
1894  }
1895  }
1896  } /* for irow */
1897  cpl_polynomial_delete(fit);
1898  } /* for k (telluric line regions) */
1899  /* correct for the airmass */
1900  cpl_table_power_column(tsens, "tellcor", 1. / aAirmass);
1901  return CPL_ERROR_NONE;
1902 } /* muse_flux_response_telluric() */
1903 
1904 /*----------------------------------------------------------------------------*/
1922 /*----------------------------------------------------------------------------*/
1923 static cpl_error_code
1924 muse_flux_response_extrapolate(muse_flux_object *aFluxObj, double aDistance)
1925 {
1926  cpl_ensure_code(aFluxObj, CPL_ERROR_NULL_INPUT);
1927 
1928  cpl_table *tsens = aFluxObj->sensitivity;
1929  cpl_propertylist *order = cpl_propertylist_new();
1930  cpl_propertylist_append_bool(order, "lambda", CPL_FALSE);
1931  cpl_table_sort(tsens, order);
1932 
1933  int nrow = cpl_table_get_nrow(tsens);
1934  /* first and last good wavelength and corresponsing error */
1935  double lambda1 = cpl_table_get_double(tsens, "lambda", 0, NULL),
1936  serr1 = cpl_table_get_double(tsens, "serr", 0, NULL);
1937  unsigned int dq = cpl_table_get_int(tsens, "dq", 0, NULL);
1938  int irow = 0;
1939  while (dq != EURO3D_GOODPIXEL && dq != EURO3D_TELLCOR) {
1940  lambda1 = cpl_table_get_double(tsens, "lambda", ++irow, NULL);
1941  serr1 = cpl_table_get_double(tsens, "serr", irow, NULL);
1942  dq = cpl_table_get_int(tsens, "dq", irow, NULL);
1943  }
1944  double lambda2 = cpl_table_get_double(tsens, "lambda", nrow - 1, NULL),
1945  serr2 = cpl_table_get_double(tsens, "serr", nrow - 1, NULL);
1946  dq = cpl_table_get_int(tsens, "dq", nrow - 1, NULL);
1947  irow = nrow - 1;
1948  while (dq != EURO3D_GOODPIXEL && dq != EURO3D_TELLCOR) {
1949  lambda2 = cpl_table_get_double(tsens, "lambda", --irow, NULL);
1950  serr2 = cpl_table_get_double(tsens, "serr", irow, NULL);
1951  dq = cpl_table_get_int(tsens, "dq", irow, NULL);
1952  }
1953  cpl_polynomial *fit1 = muse_flux_response_fit(aFluxObj, lambda1,
1954  lambda1 + aDistance, 1, 5., NULL),
1955  *fit2 = muse_flux_response_fit(aFluxObj, lambda2 - aDistance,
1956  lambda2, 1, 5., NULL);
1957  nrow = cpl_table_get_nrow(tsens); /* fitting may have erased rows! */
1958  double d1 = (lambda1 - kMuseLambdaMinX) / 100., /* want 10 additional entries... */
1959  d2 = (kMuseLambdaMaxX - lambda2) / 100.; /* ... at each end */
1960  cpl_table_set_size(tsens, nrow + 200);
1961  irow = nrow;
1962  double l;
1963  if (fit1) { /* extrapolate blue end */
1964  for (l = kMuseLambdaMinX; l <= lambda1 - d1; l += d1) {
1965  double sens = cpl_polynomial_eval_1d(fit1, l, NULL),
1966  /* error that doubles every 50 Angstrom */
1967  serr = 2. * (lambda1 - l) / 50. * serr1 + serr1;
1968  if (sens <= 0) {
1969  cpl_table_set_invalid(tsens, "lambda", irow++);
1970  cpl_msg_debug(__func__, "invalid blueward extrapolation: %.3f %f +/- %f",
1971  l, sens, serr);
1972  continue;
1973  }
1974  cpl_table_set_double(tsens, "lambda", irow, l);
1975  cpl_table_set_double(tsens, "sens", irow, sens);
1976  cpl_table_set_double(tsens, "serr", irow, serr);
1977  cpl_table_set_int(tsens, "dq", irow++, (int)EURO3D_OUTSDRANGE);
1978  } /* for l */
1979  cpl_msg_debug(__func__, "Extrapolated blue end: %.1f...%.1f Angstrom (using"
1980  " data from %.1f...%.1f Angstrom)", kMuseLambdaMinX,
1981  lambda1 - d1, lambda1, lambda1 + aDistance);
1982  } /* if blue end */
1983  if (fit2) { /* extrapolate red end */
1984  for (l = lambda2 + d2; fit2 && l <= kMuseLambdaMaxX; l += d2) {
1985  double sens = cpl_polynomial_eval_1d(fit2, l, NULL),
1986  serr = 2. * (l - lambda2) / 50. * serr2 + serr2;
1987  if (sens <= 0) {
1988  cpl_table_set_invalid(tsens, "lambda", irow++);
1989  cpl_msg_debug(__func__, "invalid redward extrapolation: %.3f %f +/- %f",
1990  l, sens, serr);
1991  continue;
1992  }
1993  cpl_table_set_double(tsens, "lambda", irow, l);
1994  cpl_table_set_double(tsens, "sens", irow, sens);
1995  cpl_table_set_double(tsens, "serr", irow, serr);
1996  cpl_table_set_int(tsens, "dq", irow++, (int)EURO3D_OUTSDRANGE);
1997  } /* for l */
1998  cpl_msg_debug(__func__, "Extrapolated red end: %.1f...%.1f Angstrom (using "
1999  "data from %.1f...%.1f Angstrom)", lambda2 + d2,
2000  kMuseLambdaMaxX, lambda2 - aDistance, lambda2);
2001  } /* if red end */
2002 #if 0
2003  cpl_polynomial_dump(fit1, stdout);
2004  cpl_polynomial_dump(fit2, stdout);
2005  fflush(stdout);
2006 #endif
2007  cpl_polynomial_delete(fit1);
2008  cpl_polynomial_delete(fit2);
2009  /* clean up invalid entries */
2010  cpl_table_select_all(tsens);
2011  cpl_table_and_selected_invalid(tsens, "sens");
2012  cpl_table_erase_selected(tsens);
2013  /* sort the resulting table again */
2014  cpl_table_sort(tsens, order);
2015  cpl_propertylist_delete(order);
2016  return CPL_ERROR_NONE;
2017 } /* muse_flux_response_extrapolate() */
2018 
2019 /*----------------------------------------------------------------------------*/
2068 /*----------------------------------------------------------------------------*/
2069 cpl_error_code
2071  muse_flux_selection_type aSelect, double aAirmass,
2072  const cpl_table *aReference,
2073  const cpl_table *aTellBands,
2074  const cpl_table *aExtinct)
2075 {
2076  cpl_ensure_code(aFluxObj && aFluxObj->intimage && aReference,
2077  CPL_ERROR_NULL_INPUT);
2078  cpl_ensure_code(aAirmass >= 1., CPL_ERROR_ILLEGAL_INPUT);
2079  if (!aExtinct) {
2080  cpl_msg_warning(__func__, "Extinction table not given!");
2081  }
2082  if (aSelect == MUSE_FLUX_SELECT_NEAREST &&
2083  (!isfinite(aFluxObj->raref) || !isfinite(aFluxObj->decref))) {
2084  cpl_msg_warning(__func__, "Reference position %f,%f contains infinite "
2085  "values, using flux to select star!", aFluxObj->raref,
2086  aFluxObj->decref);
2087  aSelect = MUSE_FLUX_SELECT_BRIGHTEST;
2088  }
2089  muse_flux_response_set_telluric_bands(aFluxObj, aTellBands);
2090 
2091  int nobjects = cpl_image_get_size_y(aFluxObj->intimage->data);
2092  const char *bunit = muse_pfits_get_bunit(aFluxObj->intimage->header);
2093  /* find brightest and nearest star */
2094  double flux = 0., dmin = DBL_MAX;
2095  int n, nstar = 1, nstardist = 1;
2096  for (n = 1; n <= nobjects; n++) {
2097  char kw[KEYWORD_LENGTH];
2098  snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_RA, n);
2099  double ra = cpl_propertylist_get_double(aFluxObj->intimage->header, kw);
2100  snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_DEC, n);
2101  double dec = cpl_propertylist_get_double(aFluxObj->intimage->header, kw),
2102  dthis = muse_astro_angular_distance(ra, dec, aFluxObj->raref,
2103  aFluxObj->decref);
2104  cpl_msg_debug(__func__, "distance(%d) = %f arcsec", n, dthis * 3600.);
2105  if (fabs(dthis) < dmin) {
2106  dmin = dthis;
2107  nstardist = n;
2108  }
2109  snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_FLUX, n);
2110  double this = cpl_propertylist_get_double(aFluxObj->intimage->header, kw);
2111  cpl_msg_debug(__func__, "flux(%d) = %e %s", n, this, bunit);
2112  if (this > flux) {
2113  flux = this;
2114  nstar = n;
2115  }
2116  } /* for n (all objects) */
2117  int nselected;
2118  char *outstring = NULL,
2119  *comment = NULL;
2120  if (aSelect == MUSE_FLUX_SELECT_BRIGHTEST) {
2121  outstring = cpl_sprintf("Selected the brightest star (%d of %d; %.3e %s)"
2122  " as reference object", nstar, nobjects, flux, bunit);
2123  comment = cpl_sprintf(MUSE_HDR_FLUX_NSEL_C, "brightest");
2124  nselected = nstar;
2125  } else {
2126  outstring = cpl_sprintf("Selected the nearest star (%d of %d; %.2f arcsec) "
2127  "as reference object", nstar, nobjects, dmin*3600.);
2128  comment = cpl_sprintf(MUSE_HDR_FLUX_NSEL_C, "nearest");
2129  nselected = nstardist;
2130  }
2131  cpl_msg_info(__func__, "%s", outstring);
2132  cpl_free(outstring);
2133 
2134  /* add selected object to the header of the intimage, *
2135  * use the comment to give the reason */
2136  cpl_propertylist_append_int(aFluxObj->intimage->header, MUSE_HDR_FLUX_NSEL,
2137  nselected);
2138  cpl_propertylist_set_comment(aFluxObj->intimage->header, MUSE_HDR_FLUX_NSEL,
2139  comment);
2140  cpl_free(comment);
2141 
2142  /* table of sensitivity function, its sigma, and a Euro3D-like quality flag */
2143  muse_flux_response_sensitivity(aFluxObj,
2144  nselected - 1, aReference,
2145  aAirmass, aExtinct);
2146  muse_flux_response_dump_sensitivity(aFluxObj, "initial");
2147 
2148  muse_flux_response_mark_questionable(aFluxObj);
2149  muse_flux_response_dump_sensitivity(aFluxObj, "intermediate");
2150  /* remove the ones outside the wavelength range */
2151  cpl_table_select_all(aFluxObj->sensitivity);
2152  cpl_table_and_selected_int(aFluxObj->sensitivity, "dq", CPL_EQUAL_TO,
2153  (int)EURO3D_OUTSDRANGE);
2154  cpl_table_erase_selected(aFluxObj->sensitivity);
2155  muse_flux_response_dump_sensitivity(aFluxObj, "intercut");
2156 
2157  /* interpolate telluric (dq == EURO3D_TELLURIC) regions *
2158  * and compute telluric correction factor */
2159  muse_flux_response_telluric(aFluxObj, aAirmass);
2160  muse_flux_response_dump_sensitivity(aFluxObj, "interpolated");
2161  /* extend the wavelength range using linear extrapolation */
2162  double width = 150.; /* use 150 Angstrom to fit line */
2163  if (cpl_propertylist_has(aFluxObj->intimage->header, MUSE_HDR_FLUX_FFCORR)) {
2164  /* The flat-field corrected curve has a sharp edge at the blue end so *
2165  * that we can only use 15 Angstrom max to fit the linear extrapolation *
2166  * function to the data; at the red end this change is irrelevant, since *
2167  * it's already smooth due to the correction of the telluric feature. */
2168  width = 15.;
2169  }
2170  muse_flux_response_extrapolate(aFluxObj, width);
2171  muse_flux_response_dump_sensitivity(aFluxObj, "extrapolated");
2172 
2173  /* store the used reference table in the flux-object, since we *
2174  * are about to return with an error code indicating success */
2175  aFluxObj->reference = cpl_table_duplicate(aReference);
2176 
2177  return CPL_ERROR_NONE;
2178 } /* muse_flux_response_compute() */
2179 
2180 /*----------------------------------------------------------------------------*/
2189 /*----------------------------------------------------------------------------*/
2191  { "lambda", CPL_TYPE_DOUBLE, "Angstrom", "%7.2f", "wavelength", CPL_TRUE },
2192  { "response", CPL_TYPE_DOUBLE,
2193  "2.5*log10((count/s/Angstrom)/(erg/s/cm**2/Angstrom))", "%.4e",
2194  "instrument response derived from standard star", CPL_TRUE },
2195  { "resperr", CPL_TYPE_DOUBLE,
2196  "2.5*log10((count/s/Angstrom)/(erg/s/cm**2/Angstrom))", "%.4e",
2197  "instrument response error derived from standard star", CPL_TRUE },
2198  { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
2199 };
2200 
2201 /*----------------------------------------------------------------------------*/
2210 /*----------------------------------------------------------------------------*/
2212  { "lambda", CPL_TYPE_DOUBLE, "Angstrom", "%7.2f", "wavelength", CPL_TRUE },
2213  { "ftelluric", CPL_TYPE_DOUBLE, "", "%.5f",
2214  "the telluric correction factor, normalized to an airmass of 1", CPL_TRUE },
2215  { "ftellerr", CPL_TYPE_DOUBLE, "", "%.5f",
2216  "the error of the telluric correction factor", CPL_TRUE },
2217  { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
2218 };
2219 
2220 /*----------------------------------------------------------------------------*/
2243 /*----------------------------------------------------------------------------*/
2244 static void
2245 muse_flux_get_response_table_smooth(cpl_table *aResp, double aHalfwidth,
2246  double aLambdaMin, double aLambdaMax,
2247  double aGapMin, double aGapMax,
2248  cpl_boolean aAverage)
2249 {
2250  cpl_msg_debug(__func__, "gap (%.3f..%.3f) wavelength range (%.3f..%.3f)",
2251  aGapMin, aGapMax, aLambdaMin, aLambdaMax);
2252  /* if a gap was given, call this function recursively, *
2253  * once for each spectral range outside the gap */
2254  if (aGapMin < aGapMax) { /* valid gap */
2255  if (aGapMin > aLambdaMin && aGapMax < aLambdaMax) {
2256  /* standard case for NaD notch: gap within the wavelength range *
2257  * need two calls, separately for both full ranges */
2258  muse_flux_get_response_table_smooth(aResp, aHalfwidth, aLambdaMin, aGapMin,
2259  0.1, -0.1, aAverage);
2260  muse_flux_get_response_table_smooth(aResp, aHalfwidth, aGapMax, aLambdaMax,
2261  0.1, -0.1, aAverage);
2262 #if 0 /* other cases that currently do not happen */
2263  } else if (aGapMin < aLambdaMin && aGapMax > aLambdaMax) {
2264  /* bad case: only NaD range, no smoothing possible */
2265  cpl_msg_warning(__func__, "No smoothing, gap (%.3f..%.3f) is larger than "
2266  "wavelength range (%.3f..%.3f).", aGapMin, aGapMax,
2267  aLambdaMin, aLambdaMax);
2268  } else if (aGapMin <= aLambdaMin && aGapMax < aLambdaMax) {
2269  /* 3rd case: gap at lower limit of wavelength range, only one call */
2270  muse_flux_get_response_table_smooth(aResp, aHalfwidth, aGapMax, aLambdaMax,
2271  0.1, -0.1, aAverage);
2272  } else if (aGapMin > aLambdaMin && aGapMax >= aLambdaMax) {
2273  /* 4th case: gap at upper limit of wavelength range , only one call */
2274  muse_flux_get_response_table_smooth(aResp, aHalfwidth, aLambdaMin, aGapMin,
2275  0.1, -0.1, aAverage);
2276 #endif
2277  }
2278  return;
2279  } /* if valid gap */
2280 
2281  /* duplicate the input columns, to not disturb the smoothing while running */
2282  cpl_table_duplicate_column(aResp, "sens", aResp, "response");
2283  cpl_table_duplicate_column(aResp, "serr", aResp, "resperr");
2284 
2285  /* select the rows which to use for the smoothing */
2286  cpl_table_select_all(aResp);
2287  cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_LESS_THAN, aLambdaMin);
2288  cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_GREATER_THAN, aLambdaMax);
2289 
2290  cpl_boolean sym = cpl_table_count_selected(aResp) < cpl_table_get_nrow(aResp);
2291  cpl_msg_debug(__func__, "%s smoothing response +/- %.3f Angstrom between %.3f "
2292  "and %.3f Angstrom", sym ? "symmetrical" : "", aHalfwidth,
2293  aLambdaMin, aLambdaMax);
2294 
2295  /* sliding median to get the values, and its median deviation for the error */
2296  int i, n = cpl_table_get_nrow(aResp);
2297  for (i = 0; i < n; i++) {
2298  if (!cpl_table_is_selected(aResp, i)) {
2299  continue;
2300  }
2301  double lambda = cpl_table_get_double(aResp, "lambda", i, NULL);
2302  int j = i, j1 = i, j2 = i;
2303  /* search for the range */
2304  while (--j > 0 && cpl_table_is_selected(aResp, j) &&
2305  lambda - cpl_table_get_double(aResp, "lambda", j, NULL) <= aHalfwidth) {
2306  j1 = j;
2307  }
2308  j = i;
2309  while (++j < n && cpl_table_is_selected(aResp, j) &&
2310  cpl_table_get_double(aResp, "lambda", j, NULL) - lambda <= aHalfwidth) {
2311  j2 = j;
2312  }
2313  if (sym) { /* adjust ranges to the smaller one for symmetrical smoothing */
2314  int jd1 = i - j1,
2315  jd2 = j2 - i;
2316  if (jd1 < jd2) {
2317  j2 = i + jd1;
2318  } else {
2319  j1 = i - jd2;
2320  }
2321  } /* if sym */
2322  double *sens = cpl_table_get_data_double(aResp, "sens"),
2323  *serr = cpl_table_get_data_double(aResp, "serr");
2324  cpl_vector *v = cpl_vector_wrap(j2 - j1 + 1, sens + j1),
2325  *ve = cpl_vector_wrap(j2 - j1 + 1, serr + j1);
2326  if (aAverage) {
2327  /* sliding average, use real stdev for >1 points */
2328  double mean = cpl_vector_get_mean(v),
2329  stdev = j2 == j1 ? 0. : cpl_vector_get_stdev(v),
2330  rerr = cpl_table_get_double(aResp, "resperr", i, NULL);
2331  cpl_table_set_double(aResp, "response", i, mean);
2332  cpl_table_set_double(aResp, "resperr", i, sqrt(rerr*rerr + stdev*stdev));
2333 #if 0
2334  cpl_msg_debug(__func__, "%d %.3f %d...%d --> %f +/- %f", i, lambda, j1, j2,
2335  mean, sqrt(rerr*rerr + stdev*stdev));
2336 #endif
2337  } else {
2338  /* sliding median */
2339  double median = cpl_vector_get_median_const(v),
2340  mdev = muse_cplvector_get_adev_const(v, median),
2341  mederr = cpl_vector_get_median_const(ve);
2342  if (j2 == j1) { /* for single points, copy the error */
2343  mdev = cpl_table_get_double(aResp, "serr", j1, NULL);
2344  }
2345  if (mdev < mederr) {
2346  mdev = mederr;
2347  }
2348 #if 0
2349  cpl_msg_debug(__func__, "%d %.3f %d...%d --> %f +/- %f", i, lambda, j1, j2,
2350  median, mdev);
2351 #endif
2352  cpl_table_set_double(aResp, "response", i, median);
2353  cpl_table_set_double(aResp, "resperr", i, mdev);
2354  }
2355  cpl_vector_unwrap(v);
2356  cpl_vector_unwrap(ve);
2357  } /* for i (all aResp rows) */
2358 
2359  /* erase the extra columns again */
2360  cpl_table_erase_column(aResp, "sens");
2361  cpl_table_erase_column(aResp, "serr");
2362 } /* muse_flux_get_response_table_smooth() */
2363 
2364 /*----------------------------------------------------------------------------*/
2379 /*----------------------------------------------------------------------------*/
2380 static unsigned int
2381 muse_flux_get_response_table_collect_points(const cpl_table *aTable,
2382  double aLambda, double aLDist,
2383  cpl_matrix *aPos, cpl_vector *aVal,
2384  cpl_vector *aErr)
2385 {
2386 
2387  unsigned int np = 0; /* counter for number of transfered points */
2388  int irow, nrow = cpl_table_get_nrow(aTable),
2389  nsel = cpl_table_count_selected(aTable);
2390  /* use the selected rows, if there are some (but not all are selected) */
2391  cpl_boolean selected = nsel > 0 && nsel != nrow;
2392  for (irow = 0; irow < nrow; irow++) {
2393  if (selected && !cpl_table_is_selected(aTable, irow)) {
2394  continue; /* ignore this row if it's not selected but we want those */
2395  }
2396  double lambda = cpl_table_get(aTable, "lambda", irow, NULL);
2397  if (lambda < aLambda - aLDist || lambda > aLambda + aLDist) {
2398  continue;
2399  }
2400  cpl_matrix_set(aPos, 0, np, lambda);
2401  cpl_vector_set(aVal, np, cpl_table_get(aTable, "sens", irow, NULL));
2402  cpl_vector_set(aErr, np, cpl_table_get(aTable, "serr", irow, NULL));
2403  np++;
2404  } /* for irow */
2405  cpl_matrix_set_size(aPos, 1, np);
2406  cpl_vector_set_size(aVal, np);
2407  cpl_vector_set_size(aErr, np);
2408  return np;
2409 } /* muse_flux_get_response_table_collect_points() */
2410 
2411 /*----------------------------------------------------------------------------*/
2433 /*----------------------------------------------------------------------------*/
2434 static void
2435 muse_flux_get_response_table_piecewise_poly(cpl_table *aResp, double aLambdaMin,
2436  double aLambdaMax, double aGapMin,
2437  double aGapMax, double aDStep,
2438  float aRSigma)
2439 {
2440  /* if a gap was given, call this function recursively, *
2441  * once for each spectral range outside the gap */
2442  cpl_msg_debug(__func__, "gap (%.3f..%.3f) wavelength range (%.3f..%.3f)",
2443  aGapMin, aGapMax, aLambdaMin, aLambdaMax);
2444  if (aGapMin < aGapMax) { /* valid gap */
2445  if (aGapMin > aLambdaMin && aGapMax < aLambdaMax) {
2446  /* standard case for NaD notch: gap within the wavelength range *
2447  * need two calls, separately for both full ranges */
2448  muse_flux_get_response_table_piecewise_poly(aResp, aLambdaMin, aGapMin,
2449  0.1, -0.1, aDStep, aRSigma);
2450  muse_flux_get_response_table_piecewise_poly(aResp, aGapMax, aLambdaMax,
2451  0.1, -0.1, aDStep, aRSigma);
2452 #if 0 /* other cases that currently do not happen */
2453  } else if (aGapMin < aLambdaMin && aGapMax > aLambdaMax) {
2454  /* bad case: only NaD range, no smoothing possible */
2455  cpl_msg_warning(__func__, "No piecewise polynomial smoothing, gap (%.3f.."
2456  "%.3f) is larger than wavelength range (%.3f..%.3f).",
2457  aGapMin, aGapMax, aLambdaMin, aLambdaMax);
2458  } else if (aGapMin <= aLambdaMin && aGapMax < aLambdaMax) {
2459  /* 3rd case: gap at lower limit of wavelength range, only one call */
2460  muse_flux_get_response_table_piecewise_poly(aResp, aGapMax, aLambdaMax,
2461  0.1, -0.1, aDStep, aRSigma);
2462  } else if (aGapMin > aLambdaMin && aGapMax >= aLambdaMax) {
2463  /* 4th case: gap at upper limit of wavelength range , only one call */
2464  muse_flux_get_response_table_piecewise_poly(aResp, aLambdaMin, aGapMin,
2465  0.1, -0.1, aDStep, aRSigma);
2466 #endif
2467  }
2468  return;
2469  } /* if valid gap */
2470 
2471  /* duplicate the input columns, to not disturb the smoothing while running */
2472  cpl_table_duplicate_column(aResp, "sens", aResp, "response");
2473  cpl_table_duplicate_column(aResp, "serr", aResp, "resperr");
2474 
2475  /* select the rows which to use for the smoothing */
2476  cpl_table_select_all(aResp);
2477  cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_LESS_THAN, aLambdaMin);
2478  cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_GREATER_THAN, aLambdaMax);
2479 
2480  /* variables to keep track of jumps that need to be fixed afterwards */
2481  unsigned int npold = 0, njumps = 0;
2482  double ldistold = -1,
2483  lambdaold = -1;
2484  cpl_array *jumppos = cpl_array_new(0, CPL_TYPE_DOUBLE),
2485  *jumplen = cpl_array_new(0, CPL_TYPE_DOUBLE);
2486 
2487  /* compute the piecewise cubic polynomial for the data around each input datapoint */
2488  int irow, nrow = cpl_table_get_nrow(aResp);
2489  for (irow = 0; irow < nrow; irow++) {
2490  if (!cpl_table_is_selected(aResp, irow)) {
2491  continue;
2492  }
2493 
2494  double lambda = cpl_table_get(aResp, "lambda", irow, NULL);
2495  /* set up breakpoints every aDStep Angstrom */
2496  double ldist = aDStep;
2497  /* collect all data ldist / 2 below and ldist / 2 above the knot wavelength */
2498  cpl_matrix *pos = cpl_matrix_new(1, nrow); /* start with far too long ones */
2499  cpl_vector *val = cpl_vector_new(nrow),
2500  *err = cpl_vector_new(nrow);
2501  unsigned int np = muse_flux_get_response_table_collect_points(aResp,
2502  lambda, ldist,
2503  pos, val, err);
2504  if (ldistold < 0) {
2505  npold = np;
2506  ldistold = ldist;
2507  lambdaold = lambda;
2508  }
2509 #if 0
2510  printf("%f Angstrom %u points (%u, %.3f, %f):\n", lambda, np, npold,
2511  (double)np / npold - 1., cpl_vector_get_mean(val));
2512  cpl_matrix_dump(pos, stdout);
2513  fflush(stdout);
2514 #endif
2515  /* the number of points changed more than 10% */
2516  if (np > 10 && fabs((double)np / npold - 1.) > 0.1) {
2517  cpl_msg_debug(__func__, "possible jump, changed at lambda %.3f (%u -> %u, "
2518  "%.3f -> %.3f)", lambda, npold, np, ldistold, ldist);
2519  cpl_array_set_size(jumppos, ++njumps);
2520  cpl_array_set_size(jumplen, njumps);
2521  cpl_array_set_double(jumppos, njumps - 1, (lambdaold + lambda) / 2.);
2522  cpl_array_set_double(jumplen, njumps - 1, (ldistold + ldist) / 2.);
2523  }
2524  /* we want to simulate a cubic spline, i.e. use order 3, *
2525  * but we have to do with the number of points we got */
2526  unsigned int order = np > 3 ? 3 : np - 1;
2527  double mse;
2528  /* fit the polynomial, but without rejection, i.e. use high rejection sigma */
2529  cpl_polynomial *poly = muse_utils_iterate_fit_polynomial(pos, val, err,
2530  NULL, order, aRSigma,
2531  &mse, NULL);
2532 #if 0
2533  if (fabs(lambda - 4861.3) < 1.) {
2534  cpl_vector *res = cpl_vector_new(cpl_vector_get_size(val));
2535  cpl_vector_fill_polynomial_fit_residual(res, val, NULL, poly, pos, NULL);
2536  double rms = sqrt(cpl_vector_product(res, res) / cpl_vector_get_size(res));
2537  cpl_msg_debug(__func__, "lambda %f rms %f (%u/%d points)", lambda, rms,
2538  (unsigned)cpl_vector_get_size(val), np);
2539  //const cpl_vector *v[] = { NULL, val, res };
2540  cpl_plot_vector(NULL, NULL, NULL, val);
2541  cpl_plot_vector(NULL, NULL, NULL, res);
2542  cpl_vector_delete(res);
2543  }
2544 #endif
2545  cpl_matrix_delete(pos);
2546  cpl_vector_delete(val);
2547  cpl_vector_delete(err);
2548  double resp = cpl_polynomial_eval_1d(poly, lambda, NULL);
2549  cpl_polynomial_delete(poly);
2550  cpl_table_set(aResp, "lambda", irow, lambda);
2551  cpl_table_set(aResp, "response", irow, resp);
2552  double serr = cpl_table_get(aResp, "serr", irow, NULL);
2553  cpl_table_set(aResp, "resperr", irow, sqrt(mse + serr*serr));
2554 
2555  npold = np;
2556  ldistold = ldist;
2557  lambdaold = lambda;
2558  } /* for i (all rows) */
2559 
2560  /* erase the extra columns again */
2561  cpl_table_erase_column(aResp, "sens");
2562  cpl_table_erase_column(aResp, "serr");
2563 
2564 #if 0
2565  printf("jumppos (%u):\n", njumps);
2566  cpl_array_dump(jumppos, 0, 10000, stdout);
2567  printf("jumplen:\n");
2568  cpl_array_dump(jumplen, 0, 10000, stdout);
2569  fflush(stdout);
2570 #endif
2571 
2572  /* this needs median smoothing afterwards as well *
2573  * to remove the jumps where the length changed */
2574  unsigned int iarr;
2575  for (iarr = 0; iarr < njumps; iarr++) {
2576  double lambda = cpl_array_get_double(jumppos, iarr, NULL),
2577  ldist = cpl_array_get_double(jumplen, iarr, NULL) / 2;
2578  /* check that the step in the +/- 5 Angstrom region actually worth worrying about */
2579  cpl_table_select_all(aResp);
2580  cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_LESS_THAN, lambda - 5.);
2581  cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_GREATER_THAN, lambda + 5.);
2582  int nsel = cpl_table_count_selected(aResp);
2583  if (nsel <= 1) {
2584  cpl_msg_debug(__func__, "Only %d points near jump around %.1f Angstrom, "
2585  "not doing anything", nsel, lambda);
2586  continue;
2587  }
2588  cpl_table *xresp = cpl_table_extract_selected(aResp);
2589  double stdev = cpl_table_get_column_stdev(xresp, "response");
2590 #if 0
2591  cpl_table_dump(xresp, 0, nsel, stdout);
2592  fflush(stdout);
2593 #endif
2594  cpl_table_delete(xresp);
2595  if (stdev < 0.001) {
2596  cpl_msg_debug(__func__, "%d points near jump around %.1f Angstrom, stdev "
2597  "only %f, not doing anything", nsel, lambda, stdev);
2598  continue;
2599  }
2600  cpl_msg_debug(__func__, "%d points near jump around %.1f Angstrom, stdev "
2601  "%f, erasing the region", nsel, lambda, stdev);
2602 
2603  /* erase the affected part of the response, linear *
2604  * interpolation when applying will take care of the gap */
2605  cpl_table_select_all(aResp);
2606  cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_LESS_THAN, lambda - ldist);
2607  cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_GREATER_THAN, lambda + ldist);
2608  cpl_table_erase_selected(aResp);
2609  } /* for iarr (all jump points) */
2610  cpl_array_delete(jumppos);
2611  cpl_array_delete(jumplen);
2612 } /* muse_flux_get_response_table_piecewise_poly() */
2613 
2614 /*----------------------------------------------------------------------------*/
2642 /*----------------------------------------------------------------------------*/
2643 cpl_error_code
2645  muse_flux_smooth_type aSmooth)
2646 {
2647  cpl_ensure_code(aFluxObj && aFluxObj->sensitivity, CPL_ERROR_NULL_INPUT);
2648  cpl_ensure_code(aSmooth <= MUSE_FLUX_SMOOTH_PPOLY, CPL_ERROR_ILLEGAL_INPUT);
2649 
2650  int nrow = cpl_table_get_nrow(aFluxObj->sensitivity);
2651  cpl_table *resp = muse_cpltable_new(muse_flux_responsetable_def, nrow);
2652  /* copy the relevant columns from the sensitivity table */
2653  const double *lambdas = cpl_table_get_data_double_const(aFluxObj->sensitivity,
2654  "lambda"),
2655  *sens = cpl_table_get_data_double_const(aFluxObj->sensitivity,
2656  "sens"),
2657  *serr = cpl_table_get_data_double_const(aFluxObj->sensitivity,
2658  "serr");
2659  cpl_table_copy_data_double(resp, "lambda", lambdas);
2660  cpl_table_copy_data_double(resp, "response", sens);
2661  cpl_table_copy_data_double(resp, "resperr", serr);
2662 
2663  /* Check, if we need to exclude the blue part. In nominal non-AO mode as *
2664  * well as in nominal AO mode, with flat-field spectrum correction, the *
2665  * sensitivity has a very narrow but real kink where the blue cutoff *
2666  * filter starts. We need to stop smoothing just redward of that. *
2667  * Also set up the smoothing gap here, if needed to skip the NaD range. */
2668  double blueend = kMuseLambdaMinX,
2669  gapmin = 0.1, gapmax = -0.1; /* setup no gap for a start */
2670  if (aFluxObj->cube && aFluxObj->cube->header &&
2671  cpl_propertylist_has(aFluxObj->cube->header, MUSE_HDR_FLUX_FFCORR)) {
2672  muse_ins_mode mode = muse_pfits_get_mode(aFluxObj->cube->header);
2673  if (mode == MUSE_MODE_WFM_NONAO_N) {
2674  /* normal blue cutoff filter, but no notch gap */
2675  blueend = kMuseNominalCutoffKink;
2676  } else if ((mode == MUSE_MODE_WFM_AO_N) || (mode == MUSE_MODE_NFM_AO_N)) {
2677  /* AO notch filter, with different blue cutoff and a gap */
2678  blueend = kMuseAOCutoffKink;
2679  gapmin = kMuseNa2LambdaMin;
2680  gapmax = kMuseNa2LambdaMax;
2681  } else if (mode == MUSE_MODE_WFM_AO_E) {
2682  /* AO-E notch filter, no blue cutoff but a different notch gap */
2683  gapmin = kMuseNaLambdaMin;
2684  gapmax = kMuseNaLambdaMax;
2685  } /* else if */
2686  } /* if */
2687 
2688  /* now smooth the response */
2689  if (aSmooth == MUSE_FLUX_SMOOTH_MEDIAN) {
2690  cpl_msg_info(__func__, "Smoothing response curve with median filter");
2691  /* keep the unsmoothed data in extra columns */
2692  cpl_table_duplicate_column(resp, "response_unsmoothed", resp, "response");
2693  cpl_table_duplicate_column(resp, "resperr_unsmoothed", resp, "resperr");
2694  /* use a sliding median over +/- 15 Angstrom width */
2695  muse_flux_get_response_table_smooth(resp, 15., blueend, kMuseLambdaMaxX,
2696  gapmin, gapmax, CPL_FALSE);
2697  } else if (aSmooth == MUSE_FLUX_SMOOTH_PPOLY) {
2698  cpl_msg_info(__func__, "Smoothing response curve with piecewise polynomial");
2699  /* keep the unsmoothed data in extra columns */
2700  cpl_table_duplicate_column(resp, "response_unsmoothed", resp, "response");
2701  cpl_table_duplicate_column(resp, "resperr_unsmoothed", resp, "resperr");
2702  /* use a piecewise polynomial, i.e. a simple, non-contiguous spline */
2703  muse_flux_get_response_table_piecewise_poly(resp, blueend, kMuseLambdaMaxX,
2704  gapmin, gapmax, 150., 3.);
2705  /* but this usually needs extra smoothing after *
2706  * the fact with a (symmetric) sliding average */
2707  muse_flux_get_response_table_smooth(resp, 15., blueend, kMuseLambdaMaxX,
2708  gapmin, gapmax, CPL_TRUE);
2709  } else {
2710  cpl_msg_warning(__func__, "NOT smoothing the response curve at all!");
2711  }
2712 
2713  /* set the table in the flux object */
2714  aFluxObj->response = muse_table_new();
2715  aFluxObj->response->table = resp;
2716  /* base the header on the cube properties, but remove some stuff */
2717  if (aFluxObj->cube) {
2718  aFluxObj->response->header = cpl_propertylist_duplicate(aFluxObj->cube->header);
2719  } else {
2720  /* at least create an empty FITS header for the response table */
2721  aFluxObj->response->header = cpl_propertylist_new();
2722  }
2723  cpl_propertylist_erase_regexp(aFluxObj->response->header, MUSE_WCS_KEYS
2724  "|^NAXIS|BUNIT", 0);
2725  /* the flat-spectrum correction status (MUSE_HDR_FLUX_FFCORR) should be *
2726  * set here since the cube creation propagates it from the pixel table, *
2727  * but there's no good way to check it here... */
2728  return CPL_ERROR_NONE;
2729 } /* muse_flux_get_response_table() */
2730 
2731 /*----------------------------------------------------------------------------*/
2750 /*----------------------------------------------------------------------------*/
2751 cpl_error_code
2753 {
2754  cpl_ensure_code(aFluxObj && aFluxObj->sensitivity, CPL_ERROR_NULL_INPUT);
2755  cpl_table *tsens = aFluxObj->sensitivity;
2756  int nrow = cpl_table_get_nrow(tsens);
2757  cpl_table *tell = muse_cpltable_new(muse_flux_tellurictable_def, nrow);
2758  /* copy the lambda and tellcor columns from the sensitivity table */
2759  cpl_table_fill_column_window_double(tell, "lambda", 0, nrow, 0);
2760  cpl_table_copy_data_double(tell, "lambda",
2761  cpl_table_get_data_double_const(tsens, "lambda"));
2762  cpl_table_fill_column_window_double(tell, "ftelluric", 0, nrow, 0);
2763  cpl_table_copy_data_double(tell, "ftelluric",
2764  cpl_table_get_data_double_const(tsens, "tellcor"));
2765  /* no (good) error estimates, use constant error as starting point */
2766 #define TELL_MAX_ERR 0.1
2767 #define TELL_MIN_ERR 1e-4
2768  cpl_table_fill_column_window_double(tell, "ftellerr", 0, nrow, TELL_MAX_ERR);
2769 
2770  /* duplicate the tellcor column, to get the invalidity info; *
2771  * pad telluric correction factors with 1. in the new table */
2772  cpl_table_duplicate_column(tell, "tellcor", tsens, "tellcor");
2773  /* pad entries adjacent to the ones with a real telluric factor with 1 */
2774  cpl_table_unselect_all(tell);
2775  int irow;
2776  for (irow = 0; irow < nrow; irow++) {
2777  int err;
2778  cpl_table_get_double(tell, "tellcor", irow, &err); /* ignore value */
2779  if (err == 0) { /* has a valid entry -> do nothing */
2780  continue;
2781  }
2782  /* invalid entry, check previous one (again) */
2783  cpl_errorstate state = cpl_errorstate_get();
2784  double ftellcor = cpl_table_get_double(tell, "tellcor", irow - 1, &err);
2785  if (!cpl_errorstate_is_equal(state)) { /* recover from possible errors */
2786  cpl_errorstate_set(state);
2787  }
2788  if (err == 0 && ftellcor != 1.) { /* exist and is not 1 -> pad */
2789  cpl_table_set_double(tell, "ftelluric", irow, 1.);
2790  continue;
2791  }
2792  /* check the next one, too */
2793  state = cpl_errorstate_get();
2794  ftellcor = cpl_table_get_double(tell, "tellcor", irow + 1, &err);
2795  if (!cpl_errorstate_is_equal(state)) { /* recover from possible errors */
2796  cpl_errorstate_set(state);
2797  }
2798  if (err == 0 && ftellcor != 1.) { /* exist and is not 1 -> pad */
2799  cpl_table_set_double(tell, "ftelluric", irow, 1.);
2800  continue;
2801  }
2802  cpl_table_select_row(tell, irow); /* surrounded by invalid -> select */
2803  } /* for irow */
2804  cpl_table_erase_selected(tell); /* erase all the still invalid ones */
2805  cpl_table_erase_column(tell, "tellcor");
2806 
2807  /* next pass: adjust the error to be only about the distance between 1 and the value */
2808  nrow = cpl_table_get_nrow(tell);
2809  for (irow = 0; irow < nrow; irow++) {
2810  int err;
2811  double dftell = 1. - cpl_table_get_double(tell, "ftelluric", irow, &err),
2812  ftellerr = cpl_table_get_double(tell, "ftellerr", irow, &err);
2813  if (ftellerr > dftell) {
2814  ftellerr = fmax(dftell, TELL_MIN_ERR);
2815  cpl_table_set_double(tell, "ftellerr", irow, ftellerr);
2816  }
2817  } /* for irow */
2818 
2819  aFluxObj->telluric = muse_table_new();
2820  aFluxObj->telluric->table = tell;
2821  /* base the header on the cube properties, but remove some stuff */
2822  aFluxObj->telluric->header = cpl_propertylist_duplicate(aFluxObj->cube->header);
2823  cpl_propertylist_erase_regexp(aFluxObj->telluric->header, MUSE_WCS_KEYS
2824  "|^NAXIS|BUNIT", 0);
2825  /* again, flat-spectrum correction status should be transferred by this */
2826  return CPL_ERROR_NONE;
2827 } /* muse_flux_get_telluric_table() */
2828 
2829 /*----------------------------------------------------------------------------*/
2856 /*----------------------------------------------------------------------------*/
2857 cpl_error_code
2859 {
2860  cpl_ensure_code(aFluxObj && (aFluxObj->sensitivity || aFluxObj->response),
2861  CPL_ERROR_NULL_INPUT);
2862  cpl_boolean hasresp = aFluxObj->response != NULL;
2863 
2864  /* if we have the headers, we can transfer OBS.TARG.NAME *
2865  * to a QC parameter to hold the star's name */
2866  if (hasresp && aFluxObj->response->header && aFluxObj->cube &&
2867  aFluxObj->cube->header) {
2868  cpl_errorstate state = cpl_errorstate_get();
2869  const char *name = muse_pfits_get_targname(aFluxObj->cube->header);
2870  if (!name) {
2871  cpl_msg_warning(__func__, "Unknown standard star in exposure (missing "
2872  "OBS.TARG.NAME)!");
2873  cpl_errorstate_set(state);
2874  name = "UNKNOWN";
2875  }
2876  cpl_propertylist_update_string(aFluxObj->response->header, QC_STD_NAME,
2877  name);
2878  } /* if reponse header */
2879 
2880  /* shortcut to the relevant table, *
2881  * modify it by appending a throughput column ("throughput") */
2882  cpl_table *thruput = hasresp ? aFluxObj->response->table
2883  : aFluxObj->sensitivity;
2884  cpl_ensure_code(thruput, CPL_ERROR_DATA_NOT_FOUND);
2885 
2886  cpl_msg_info(__func__, "Computing throughput using effective VLT area of %.1f"
2887  " cm**2, from %s curve", kVLTArea,
2888  hasresp ? "smoothed response" : "unsmoothed sensitivity");
2889 
2890  const double h = CPL_PHYS_H * 1e7, /* h in cgs units */
2891  c = CPL_PHYS_C * 1e10; /* c in Angstrom */
2892  if (!cpl_table_has_column(thruput, "throughput")) {
2893  cpl_table_new_column(thruput, "throughput", CPL_TYPE_DOUBLE);
2894  }
2895  /* easier done in a loop than with CPL table functions... */
2896  int i, n = cpl_table_get_nrow(thruput);
2897  for (i = 0; i < n; i++) {
2898  int err1, err2;
2899  double lambda = cpl_table_get(thruput, "lambda", i, &err1),
2900  sens = cpl_table_get(thruput, hasresp ? "response" : "sens", i, &err2);
2901  if (err1 || err2) {
2902  cpl_table_set_invalid(thruput, "throughput", i);
2903  continue;
2904  }
2905  double value = h * c / kVLTArea / lambda * pow(10., 0.4 * sens);
2906 #if 0
2907  cpl_msg_debug(__func__, "%.3f %.3f --> %e, %e ==> %.3f", lambda,
2908  sens, h * c / kVLTArea / lambda, pow(10., 0.4 * sens),
2909  value);
2910 #endif
2911  cpl_table_set_double(thruput, "throughput", i, value);
2912  } /* for i (all table rows) */
2913 
2914  /* average the throughput curve at +/- 100 Angstrom around five points */
2915  cpl_msg_indent_more();
2916  float lambda;
2917  for (lambda = 5000.; lambda < 9100.; lambda += 1000.) {
2918  cpl_table_unselect_all(thruput);
2919  cpl_table_or_selected_double(thruput, "lambda", CPL_NOT_LESS_THAN,
2920  lambda - 100.);
2921  cpl_table_and_selected_double(thruput, "lambda", CPL_NOT_GREATER_THAN,
2922  lambda + 100.);
2923  cpl_table *t = cpl_table_extract_selected(thruput);
2924  double mean = cpl_table_get_column_mean(t, "throughput"),
2925  stdev = cpl_table_get_column_stdev(t, "throughput");
2926  cpl_msg_info(__func__, "Throughput at %.0f +/- 100 Angstrom: %.4f +/- %.4f",
2927  lambda, mean, stdev);
2928  cpl_table_delete(t);
2929 
2930  /* save as QC parameters in the response table header (if created) */
2931  if (hasresp && aFluxObj->response->header) {
2932  char *kw = cpl_sprintf(QC_STD_THRU, lambda);
2933  cpl_propertylist_update_float(aFluxObj->response->header, kw, mean);
2934  cpl_free(kw);
2935  } /* if reponse header */
2936  } /* for lambda */
2937  cpl_msg_indent_less();
2938 
2939  return CPL_ERROR_NONE;
2940 } /* muse_flux_compute_qc() */
2941 
2942 /*----------------------------------------------------------------------------*/
2970 /*----------------------------------------------------------------------------*/
2971 cpl_error_code
2973  const char *aName)
2974 {
2975  cpl_ensure_code(aFluxObj && (aFluxObj->sensitivity || aFluxObj->response) &&
2976  aFluxObj->reference && aFilter && aFilter->table,
2977  CPL_ERROR_NULL_INPUT);
2978 
2979  /* use shortnames for the filters */
2980  char *p = aName ? strrchr(aName, '_') : NULL;
2981  const char *fshort = NULL;
2982  if (p) {
2983  fshort = p + 1;
2984  } else {
2985  fshort = "UNKNOWN";
2986  cpl_msg_warning(__func__, "%s filter given for QC zeropoint computation!",
2987  fshort);
2988  }
2989  char *kw = cpl_sprintf(QC_STD_ZP, fshort); /* FITS keyword name */
2990 
2991  /* shortcut to response curve, prefer the (smoothed!) final version */
2992  cpl_boolean hasresp = aFluxObj->response != NULL;
2993  cpl_table *thruput = hasresp ? aFluxObj->response->table
2994  : aFluxObj->sensitivity;
2995  const double h = CPL_PHYS_H * 1e7, /* h in cgs units */
2996  c = CPL_PHYS_C * 1e10; /* c in Angstrom */
2997 
2998  /* Now compute the weighted sum of the (smoothed) response curve *
2999  * over the filter function, which is the zeropoint in this filter. */
3000  double fref = 0., /* summed flux of the reference spectrum */
3001  fobs = 0.; /* summed flux of the observed spectrum */
3002  int i, n = cpl_table_get_nrow(thruput);
3003  for (i = 0; i < n; i++) {
3004  int err1, err2;
3005  double lambda = cpl_table_get(thruput, "lambda", i, &err1),
3006  sens = cpl_table_get(thruput, hasresp ? "response" : "sens", i, &err2),
3007  ref = muse_flux_response_interpolate(aFluxObj->reference, lambda,
3008  NULL, MUSE_FLUX_RESP_STD_FLUX);
3009  if (err1 || err2) {
3010  continue;
3011  }
3012  double fweight = muse_flux_response_interpolate(aFilter->table, lambda, NULL,
3014 
3015  /* recompute (smoothed) measured value at this wavelength, add it up */
3016  fobs += pow(10., 0.4 * sens) /* convert back to count (e-) */
3017  * ref /* take out the reference spectrum */
3018  * fweight /* weight by the filter function */
3019  * h * c / kVLTArea / lambda; /* convert to f-lambda */
3020  /* the weighted reference flux */
3021  fref += ref * fweight;
3022  } /* for i (all thruput table rows) */
3023 
3024  /* no normalization by the total filter-weight is necessary, *
3025  * as we used the same weights in numerator and denominator */
3026  double zp = -2.5 * log10(fobs / fref); /* zeropoint is in magnitudes */
3027  cpl_msg_indent_more();
3028  cpl_msg_info(__func__, "Zeropoint in filter %s: %.3f mag (throughput %.3f)",
3029  fshort, zp, pow(10., -0.4 * zp));
3030  cpl_msg_indent_less();
3031  if (aFluxObj->response->header) {
3032  /* finally add the QC keyword to the header */
3033  cpl_propertylist_update_float(aFluxObj->response->header, kw, zp);
3034  }
3035  cpl_free(kw);
3036 
3037  return CPL_ERROR_NONE;
3038 } /* muse_flux_compute_qc_zp() */
3039 
3040 /*----------------------------------------------------------------------------*/
3081 /*----------------------------------------------------------------------------*/
3082 cpl_error_code
3083 muse_flux_calibrate(muse_pixtable *aPixtable, const muse_table *aResponse,
3084  const cpl_table *aExtinction, const muse_table *aTelluric)
3085 {
3086  cpl_ensure_code(aPixtable && aPixtable->header && aResponse,
3087  CPL_ERROR_NULL_INPUT);
3088  const char *unitdata = cpl_table_get_column_unit(aPixtable->table,
3090  cpl_ensure_code(unitdata && !strncmp(unitdata, "count", 6),
3091  CPL_ERROR_INCOMPATIBLE_INPUT);
3092 
3093  /* the code below needs only the table-component of the *
3094  * response curve and the telluric correction spectrum */
3095  cpl_table *response = NULL,
3096  *telluric = NULL;
3097  /* get the tag of the pixel table for possible outputs */
3098  const char *catg = muse_pfits_get_pro_catg(aPixtable->header);
3099 #define FF_MSG "(flat-field spectrum %scorrected)"
3100  if (aResponse) {
3101  response = aResponse->table;
3102  /* cross-check the header of the response curve against *
3103  * the flat-field spectrum status of the pixel table */
3104  cpl_boolean respff = cpl_propertylist_has(aResponse->header,
3105  MUSE_HDR_FLUX_FFCORR),
3106  ptff = cpl_propertylist_has(aPixtable->header,
3108  if ((respff && !ptff) || (!respff && ptff)) {
3109  return cpl_error_set_message(__func__, CPL_ERROR_INCOMPATIBLE_INPUT,
3110  "Cannot apply this %s "FF_MSG" to this %s "FF_MSG,
3111  MUSE_TAG_STD_RESPONSE, !respff ? "un" : "",
3112  catg, !ptff ? "un" : "");
3113  } /* if unmatched ff status */
3114  } /* if aResponse */
3115  double airmass = muse_astro_airmass(aPixtable->header);
3116  if (airmass < 0) {
3117  cpl_msg_warning(__func__, "Airmass unknown, not doing extinction "
3118  "correction: %s", cpl_error_get_message());
3119  /* reset to zero so that it has no effect */
3120  airmass = 0.;
3121  }
3122  if (aTelluric) {
3123  /* duplicate the telluric correction table to apply the airmass correction */
3124  telluric = cpl_table_duplicate(aTelluric->table);
3125  /* use negative exponent to be able to multiply (instead of divide) *
3126  * by the correction factor (should be slightly faster) */
3127  cpl_table_power_column(telluric, "ftelluric", -airmass);
3128  /* not critical, but still cross-check flat-field spectrum status here */
3129  cpl_boolean tellff = cpl_propertylist_has(aTelluric->header,
3130  MUSE_HDR_FLUX_FFCORR),
3131  ptff = cpl_propertylist_has(aPixtable->header,
3133  if ((tellff && !ptff) || (!tellff && ptff)) {
3134  cpl_msg_warning(__func__, "Applying %s "FF_MSG" to %s "FF_MSG", this may "
3135  "not be correct!", MUSE_TAG_STD_TELLURIC,
3136  !tellff ? "un" : "", catg, !ptff ? "un" : "");
3137  } /* if unmatched ff status */
3138  } /* if aTelluric */
3139 
3140  /* warn for non-critical failure */
3141  if (!aExtinction) {
3142  cpl_msg_warning(__func__, "%s missing!", MUSE_TAG_EXTINCT_TABLE);
3143  }
3144 
3145  double exptime = muse_pfits_get_exptime(aPixtable->header);
3146  if (exptime <= 0.) {
3147  cpl_msg_warning(__func__, "Non-positive EXPTIME, not doing flux calibration!");
3148  cpl_table_delete(telluric);
3149  return CPL_ERROR_ILLEGAL_INPUT;
3150  }
3151 
3152  cpl_msg_info(__func__, "Starting flux calibration (exptime=%.2fs, airmass=%.4f),"
3153  " %s telluric correction", exptime, airmass,
3154  aTelluric ? "with" : "without ("MUSE_TAG_STD_TELLURIC" not given)");
3155  float *lambda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
3156  *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
3157  *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
3158  cpl_size i, nrow = muse_pixtable_get_nrow(aPixtable);
3159  #pragma omp parallel for default(none) /* as req. by Ralf */ \
3160  shared(aExtinction, airmass, data, exptime, lambda, nrow, response, \
3161  stat, telluric)
3162  for (i = 0; i < nrow; i++) {
3163  /* double values for intermediate results of this row */
3164  double ddata = data[i], dstat = stat[i];
3165 
3166  /* correct for extinction */
3167  if (aExtinction) {
3168  double fext = pow(10., 0.4 * airmass
3169  * muse_flux_response_interpolate(aExtinction,
3170  lambda[i], NULL,
3172 #if 0
3173  printf("%f, data/stat = %f/%f -> ", fext, ddata, dstat);
3174 #endif
3175  ddata *= fext;
3176  dstat *= (fext * fext);
3177 #if 0
3178  printf(" --> %f/%f\n", ddata, dstat), fflush(NULL);
3179 #endif
3180  }
3181 
3182  /* the difference in lambda coverage per pixel seems to be *
3183  * corrected for by the flat-field, so assuming a constant *
3184  * value for all pixels seems to be the correct thing to do */
3185  double dlambda = kMuseSpectralSamplingA,
3186  dlamerr = 0.02, /* assume typical fixed error for the moment */
3187  /* resp/rerr get returned in mag units as in the table */
3188  rerr, resp = muse_flux_response_interpolate(response, lambda[i],
3189  &rerr,
3191  /* convert from 2.5 log10(x) to non-log flux units */
3192  resp = pow(10., 0.4 * resp);
3193  /* magerr = 2.5 / log(10) * error / flux (see IRAF phot docs) *
3194  * ==> error = magerr / 2.5 * log(10) * flux */
3195  rerr = rerr * CPL_MATH_LN10 * resp / 2.5;
3196 #if 0
3197  printf("%f/%f/%f, %e/%e, data/stat = %e/%e -> ", lambda[i], dlambda, dlamerr, resp, rerr,
3198  ddata, dstat);
3199 #endif
3200  dstat = dstat * pow((1./(resp * exptime * dlambda)), 2)
3201  + pow(ddata * rerr / (resp*resp * exptime * dlambda), 2)
3202  + pow(ddata * dlamerr / (resp * exptime * dlambda*dlambda), 2);
3203  ddata /= (resp * exptime * dlambda);
3204 #if 0
3205  printf("%e/%e\n", ddata, dstat), fflush(NULL);
3206 #endif
3207 
3208  /* now convert to the float values to be stored in the pixel table, *
3209  * scaled by kMuseFluxUnitFactor to keep the floats from underflowing */
3210  ddata *= kMuseFluxUnitFactor;
3211  dstat *= kMuseFluxStatFactor;
3212 
3213  /* do the telluric correction, if the wavelength is redward of the start *
3214  * of the telluric regions, and if a telluric correction was given */
3215  if (lambda[i] < kTelluricBands[0][0] || !telluric) {
3216  data[i] = ddata;
3217  stat[i] = dstat;
3218  continue; /* skip telluric correction in the blue */
3219  }
3220  double terr, tell = muse_flux_response_interpolate(telluric, lambda[i],
3221  &terr,
3223  data[i] = ddata * tell;
3224  stat[i] = tell*tell * dstat + ddata*ddata * terr*terr;
3225  } /* for i (pixel table rows) */
3226  cpl_table_delete(telluric); /* NULL check done there... */
3227 
3228  /* now set the table column headers reflecting the flux units used */
3229  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_DATA,
3230  kMuseFluxUnitString);
3231  cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_STAT,
3232  kMuseFluxStatString);
3233 
3234  /* add the status header */
3235  cpl_propertylist_update_bool(aPixtable->header, MUSE_HDR_PT_FLUXCAL,
3236  CPL_TRUE);
3237  cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_FLUXCAL,
3238  MUSE_HDR_PT_FLUXCAL_COMMENT);
3239  return CPL_ERROR_NONE;
3240 } /* muse_flux_calibrate() */
3241 
#define MUSE_PIXTABLE_DQ
Definition: muse_pixtable.h:49
Structure definition of a MUSE datacube.
Definition: muse_datacube.h:48
cpl_error_code muse_flux_compute_qc(muse_flux_object *aFluxObj)
Compute QC parameters, related to on-sky throughput.
Definition: muse_flux.c:2858
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
muse_image * intimage
Definition: muse_flux.h:83
muse_flux_smooth_type
Type of response curve smoothing to use.
Definition: muse_flux.h:146
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_pfits_get_ra(const cpl_propertylist *aHeaders)
find out the right ascension
Definition: muse_pfits.c:267
#define MUSE_HDR_PT_FLUXCAL
cpl_image * data
the data extension
Definition: muse_image.h:46
const muse_cpltable_def muse_flux_responsetable_def[]
MUSE response table definition.
Definition: muse_flux.c:2190
cpl_size muse_pixtable_get_nrow(const muse_pixtable *aPixtable)
get the number of rows within the pixel table
muse_flux_profile_type
Type of optimal profile to use.
Definition: muse_flux.h:124
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.
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.
#define MUSE_HDR_PT_FFCORR
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:170
cpl_table * reference
Definition: muse_flux.h:85
const char * muse_pfits_get_bunit(const cpl_propertylist *aHeaders)
find out the unit string
Definition: muse_pfits.c:223
muse_table * muse_table_new(void)
Allocate memory for a new muse_table object.
Definition: muse_table.c:63
const muse_cpltable_def muse_response_tellbands_def[]
Table definition for a telluric bands table.
Definition: muse_flux.c:1486
cpl_table * sensitivity
Definition: muse_flux.h:88
Structure definition of MUSE three extension FITS file.
Definition: muse_image.h:40
cpl_error_code muse_flux_calibrate(muse_pixtable *aPixtable, const muse_table *aResponse, const cpl_table *aExtinction, const muse_table *aTelluric)
Convert the input pixel table from counts to fluxes.
Definition: muse_flux.c:3083
muse_flux_object * muse_flux_object_new(void)
Allocate memory for a new muse_flux_object object.
Definition: muse_flux.c:68
static double muse_flux_reference_table_sampling(cpl_table *aTable)
Compute average sampling for a MUSE-format flux reference table.
Definition: muse_flux.c:121
cpl_table * table
The pixel table.
const char * muse_pfits_get_cunit(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS axis unit string
Definition: muse_pfits.c:491
cpl_propertylist * header
the FITS header
Definition: muse_image.h:72
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. ...
#define MUSE_PIXTABLE_DATA
Definition: muse_pixtable.h:48
muse_table * telluric
Definition: muse_flux.h:92
muse_resampling_crstats_type crtype
cpl_image * dq
the data quality extension
Definition: muse_image.h:56
const char * muse_pfits_get_targname(const cpl_propertylist *aHeaders)
find out the ESO observation target name
Definition: muse_pfits.c:1677
cpl_table * muse_cpltable_new(const muse_cpltable_def *aDef, cpl_size aLength)
Create an empty table according to the specified definition.
double muse_pfits_get_crpix(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS reference point
Definition: muse_pfits.c:401
cpl_error_code muse_flux_get_telluric_table(muse_flux_object *aFluxObj)
Get the table of the telluric correction.
Definition: muse_flux.c:2752
cpl_error_code muse_flux_compute_qc_zp(muse_flux_object *aFluxObj, const muse_table *aFilter, const char *aName)
Compute QC zeropoint for given filter.
Definition: muse_flux.c:2972
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:338
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:1097
Structure definition of MUSE pixel table.
Flux object to store data needed while computing the flux calibration.
Definition: muse_flux.h:78
double muse_astro_wavelength_vacuum_to_air(double aVac)
Compute air wavelength for a given vacuum wavelength.
Definition: muse_astro.c:534
#define MUSE_WCS_KEYS
regular expression for WCS properties
Definition: muse_wcs.h:48
cpl_error_code muse_wcs_celestial_from_pixel(cpl_propertylist *aHeader, double aX, double aY, double *aRA, double *aDEC)
Convert from pixel coordinates to celestial spherical coordinates.
Definition: muse_wcs.c:1559
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
cpl_error_code muse_flux_response_compute(muse_flux_object *aFluxObj, muse_flux_selection_type aSelect, 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:2070
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
muse_table * response
Definition: muse_flux.h:90
cpl_table * table
The table.
Definition: muse_table.h:49
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:1078
muse_datacube * cube
Definition: muse_flux.h:80
cpl_imagelist * data
the cube containing the actual data values
Definition: muse_datacube.h:76
double muse_pfits_get_dec(const cpl_propertylist *aHeaders)
find out the declination
Definition: muse_pfits.c:285
muse_flux_interpolation_type
Type of table interpolation to use.
Definition: muse_flux.h:111
Structure to store a table together with a property list.
Definition: muse_table.h:43
cpl_error_code muse_flux_get_response_table(muse_flux_object *aFluxObj, muse_flux_smooth_type aSmooth)
Get the table of the standard star response function.
Definition: muse_flux.c:2644
void muse_table_delete(muse_table *aTable)
Deallocate memory associated to a muse_table object.
Definition: muse_table.c:80
cpl_imagelist * dq
the optional cube containing the bad pixel status
Definition: muse_datacube.h:81
cpl_table * tellbands
Definition: muse_flux.h:94
cpl_propertylist * header
the header
Definition: muse_table.h:56
const muse_cpltable_def muse_flux_tellurictable_def[]
MUSE telluric correction table definition.
Definition: muse_flux.c:2211
const char * muse_pfits_get_ctype(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS axis type string
Definition: muse_pfits.c:469
cpl_error_code muse_pixtable_save(muse_pixtable *aPixtable, const char *aFilename)
Save a MUSE pixel table to a file on disk.
muse_image * muse_flux_integrate_cube(muse_datacube *aCube, cpl_apertures *aApertures, muse_flux_profile_type aProfile)
Integrate the flux of the standard star(s) given a datacube.
Definition: muse_flux.c:1009
cpl_error_code muse_flux_integrate_std(muse_pixtable *aPixtable, muse_flux_profile_type aProfile, muse_flux_object *aFluxObj)
Reconstruct a cube, detect the standard star, and integrate its flux.
Definition: muse_flux.c:1362
cpl_propertylist * header
the FITS header
Definition: muse_datacube.h:57
#define MUSE_PIXTABLE_STAT
Definition: muse_pixtable.h:50
#define MUSE_PIXTABLE_LAMBDA
Definition: muse_pixtable.h:53
double muse_pfits_get_exptime(const cpl_propertylist *aHeaders)
find out the exposure time
Definition: muse_pfits.c:382
const char * muse_pfits_get_pro_catg(const cpl_propertylist *aHeaders)
find out the PRO category
Definition: muse_pfits.c:159
cpl_propertylist * muse_wcs_apply_cd(const cpl_propertylist *aExpHeader, const cpl_propertylist *aCalHeader)
Apply the CDi_j matrix of an astrometric solution to an observation.
Definition: muse_wcs.c:1237
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:66
muse_flux_selection_type
Type of star selection to use.
Definition: muse_flux.h:136
void muse_flux_object_delete(muse_flux_object *aFluxObj)
Deallocate memory associated to a muse_flux_object.
Definition: muse_flux.c:89
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
double muse_astro_angular_distance(double aRA1, double aDEC1, double aRA2, double aDEC2)
Compute angular distance in the sky between two positions.
Definition: muse_astro.c:496
muse_ins_mode muse_pfits_get_mode(const cpl_propertylist *aHeaders)
find out the observation mode
Definition: muse_pfits.c:1352
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.
cpl_propertylist * muse_wcs_create_default(void)
Create FITS headers containing a default (relative) WCS.
Definition: muse_wcs.c:1174