MUSE Pipeline Reference Manual  2.1.1
muse_idp.c
1 /* -*- Mode: C; tab-width: 2; indent-tabs-mode: nil; c-basic-offset: 2 -*- */
2 /* vim:set sw=2 sts=2 et cin: */
3 /*
4  * This file is part of the MUSE Instrument Pipeline
5  * Copyright (C) 2015 European Southern Observatory
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software
19  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
20  */
21 
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
25 
26 /*----------------------------------------------------------------------------*
27  * Includes *
28  *----------------------------------------------------------------------------*/
29 #include <string.h>
30 
31 #include "muse_dfs.h"
32 #include "muse_pfits.h"
33 #include "muse_pixtable.h"
34 #include "muse_cplwrappers.h"
35 #include "muse_wcs.h"
36 #include "muse_flux.h"
37 #include "muse_utils.h"
38 #include "muse_idp.h"
39 
40 /*----------------------------------------------------------------------------*
41  * Defines *
42  *----------------------------------------------------------------------------*/
43 
44 #define KEY_ASSON "ASSON"
45 #define KEY_ASSON_COMMENT "Associated file name"
46 #define KEY_DEC "DEC"
47 #define KEY_DEC_COMMENT "[deg] Image center (J2000)"
48 #define KEY_EXPTIME "EXPTIME"
49 #define KEY_EXPTIME_COMMENT "[s] Total integration time per pixel"
50 #define KEY_FLUXCAL "FLUXCAL"
51 #define KEY_FLUXCAL_COMMENT "Type of flux calibration (ABSOLUTE or UNCALIBRATED)"
52 #define KEY_FLUXCAL_VALUE_FALSE "UNCALIBRATED"
53 #define KEY_FLUXCAL_VALUE_TRUE "ABSOLUTE"
54 #define KEY_MJDOBS "MJD-OBS"
55 #define KEY_MJDOBS_COMMENT "[d] Start of observations (days)"
56 #define KEY_MJDEND "MJD-END"
57 #define KEY_MJDEND_COMMENT "[d] End of observations (days)"
58 #define KEY_OBID "OBID"
59 #define KEY_OBID_COMMENT "Observation block ID"
60 #define KEY_OBSTECH "OBSTECH"
61 #define KEY_OBSTECH_COMMENT "Technique for observation"
62 #define KEY_PROCSOFT "PROCSOFT"
63 #define KEY_PROCSOFT_COMMENT "ESO pipeline version"
64 #define KEY_PRODCATG "PRODCATG"
65 #define KEY_PRODCATG_COMMENT "Data product category"
66 #define KEY_PRODCATG_VALUE_IFS_CUBE "SCIENCE.CUBE.IFS"
67 #define KEY_PRODCATG_VALUE_FOV_IMAGE "ANCILLARY.IMAGE"
68 #define KEY_PROG_ID "PROG_ID"
69 #define KEY_PROG_ID_COMMENT "ESO programme identification"
70 #define KEY_PROG_ID_VALUE_MULTIPLE "MULTI"
71 #define KEY_PROGID "PROGID"
72 #define KEY_PROGID_COMMENT KEY_PROG_ID_COMMENT
73 #define KEY_PROV "PROV"
74 #define KEY_PROV_COMMENT "Originating raw science file"
75 #define KEY_RA "RA"
76 #define KEY_RA_COMMENT "[deg] Image center (J2000)"
77 #define KEY_REFERENC "REFERENC"
78 #define KEY_REFERENC_COMMENT "Reference publication"
79 #define KEY_TEXPTIME "TEXPTIME"
80 #define KEY_TEXPTIME_COMMENT "[s] Total integration time of all exposures"
81 #define KEY_WAVELMIN "WAVELMIN"
82 #define KEY_WAVELMIN_COMMENT "[nm] Minimum wavelength"
83 #define KEY_WAVELMAX "WAVELMAX"
84 #define KEY_WAVELMAX_COMMENT "[nm] Maximum wavelength"
85 #define KEY_SKY_RES "SKY_RES"
86 #define KEY_SKY_RES_COMMENT "[arcsec] FWHM effective spatial resolution (measured)"
87 #define KEY_SKY_RERR "SKY_RERR"
88 #define KEY_SKY_RERR_COMMENT "[arcsec] Error of SKY_RES (estimated)"
89 #define KEY_SPEC_RES "SPEC_RES"
90 #define KEY_SPEC_RES_COMMENT "Spectral resolving power at central wavelength"
91 #define KEY_SPEC_ERR "CRDER3"
92 #define KEY_SPEC_ERR_COMMENT "[Angstrom] Random error in spectral coordinate"
93 #define KEY_PIXNOISE "PIXNOISE"
94 #define KEY_PIXNOISE_COMMENT "[erg/s/cm**2/Angstrom] pixel-to-pixel noise"
95 #define KEY_ABMAGLIM "ABMAGLIM"
96 #define KEY_ABMAGLIM_COMMENT "5-sigma magnitude limit for point sources"
97 #define KEY_NCOMBINE "NCOMBINE"
98 #define KEY_NCOMBINE_COMMENT "No. of combined raw science data files"
99 
100 /* A regular expression matching SDP IFS standard keywords to be cleared. */
101 #define CLEAN_KEYS_REGEXP \
102  "^(" KEY_MJDEND "|" \
103  KEY_PROCSOFT "|" \
104  KEY_PRODCATG "|" \
105  KEY_PROG_ID "|" \
106  KEY_PROGID "[0-9]+|" \
107  KEY_OBID "[0-9]+|" \
108  KEY_OBSTECH "|" \
109  KEY_FLUXCAL "|" \
110  KEY_TEXPTIME "|" \
111  KEY_WAVELMIN "|" \
112  KEY_WAVELMAX "|" \
113  KEY_SKY_RES "|" \
114  KEY_SKY_RERR "|" \
115  KEY_SPEC_RES "|" \
116  KEY_PIXNOISE "|" \
117  KEY_ABMAGLIM "|" \
118  KEY_REFERENC "|" \
119  KEY_NCOMBINE "|" \
120  KEY_PROV "[0-9]+|" \
121  KEY_ASSON "[0-9]+" ")$"
122 
123 /*----------------------------------------------------------------------------*/
127 /*----------------------------------------------------------------------------*/
128 
131 /* Lookup table for spectral resolution vs. wavelength */
132 
133 typedef struct {
134  /* Wavelength in Angstrom at which the spectral resolution was measured. */
135  double lambda;
136  /* Measured spectral resolution. */
137  double R;
138 } muse_specres_lut_entry;
139 
140 /* The last entry in the data tables must have a zero wavelength entry! */
141 
142 static const muse_specres_lut_entry
143 muse_specres_data_wfm[] =
144 {
145  {4650.0000, 1674.77},
146  {4810.3448, 1769.45},
147  {4970.6897, 1864.50},
148  {5131.0345, 1959.52},
149  {5291.3793, 2054.19},
150  {5451.7242, 2148.34},
151  {5612.0690, 2241.81},
152  {5772.4138, 2334.46},
153  {5932.7587, 2426.18},
154  {6093.1035, 2516.84},
155  {6253.4483, 2606.28},
156  {6413.7932, 2694.36},
157  {6574.1380, 2780.89},
158  {6734.4828, 2865.70},
159  {6894.8277, 2948.59},
160  {7055.1725, 3029.32},
161  {7215.5173, 3107.70},
162  {7375.8622, 3183.46},
163  {7536.2070, 3256.36},
164  {7696.5518, 3326.12},
165  {7856.8967, 3392.45},
166  {8017.2415, 3455.01},
167  {8177.5863, 3513.44},
168  {8337.9312, 3567.37},
169  {8498.2760, 3616.34},
170  {8658.6208, 3659.87},
171  {8818.9657, 3697.42},
172  {8979.3105, 3728.38},
173  {9139.6553, 3752.11},
174  {9300.0002, 3767.93},
175  { 0.0000 , 0.00}
176 };
177 
178 /* XXX: Dummy values from the User Manual for the NFM. Update for production. */
179 static const muse_specres_lut_entry
180 muse_specres_data_nfm[] =
181 {
182  {4800.0000, 1740.00},
183  {9300.0000, 3450.00},
184  { 0.0000 , 0.00}
185 };
186 
187 /* Instrument properties and reference values */
188 
189 /* Number of detector readout ports */
190 static const unsigned int kMuseNumQuadrants = 4;
191 
192 /* Nominal field size [arcsec] */
193 static const double kMuseFovSizeWFM = 60.00;
194 static const double kMuseFovSizeNFM = 7.43;
195 
196 /* Nominal pixel scale wfm, nfm */
197 static const double kMusePixscaleWFM = 0.2;
198 static const double kMusePixscaleNFM = 0.025;
199 
200 /* Reference wavelength [Angstrom] */
201 static const double kMuseLambdaRef = 7000.;
202 
203 /* Instrument response at reference wavelength */
204 static const double kMuseResponseRef = 41.2;
205 
206 /* Typical error estimated from unresolved sources [arcsec] */
207 static const double kMuseExtraErrorIQE = 0.3;
208 
209 /* Typical values [arcsec] for the effective spatial resolution */
210 /* and its uncertainty, given by the median and the standard */
211 /* deviation of previous measurements. */
212 static const double kMuseSkyRes = 0.853823;
213 static const double kMuseSkyResError = 0.495547;
214 
215 
216 /* Unit conversion factors */
217 
218 static const double day2sec = 86400.0; /* Seconds per day */
219 static const double m2nm = 1.e9; /* meter to nanometer */
220 static const double nm2Angstrom = 10.; /* nanometer to Angstrom */
221 static const double deg2as = 3600.; /* degrees to arcseconds */
222 
223 
224 /*----------------------------------------------------------------------------*
225  * Functions *
226  *----------------------------------------------------------------------------*/
227 
228 /*----------------------------------------------------------------------------*/
237 /*----------------------------------------------------------------------------*/
238 static cpl_boolean
239 muse_pfits_get_fluxcal(const cpl_propertylist *aHeaders)
240 {
241  cpl_ensure(aHeaders, CPL_ERROR_NULL_INPUT, CPL_FALSE);
242  cpl_errorstate prestate = cpl_errorstate_get();
243  cpl_boolean value = cpl_propertylist_get_bool(aHeaders,
245  cpl_errorstate_set(prestate);
246  return value;
247 }
248 
249 /*----------------------------------------------------------------------------*/
266 /*----------------------------------------------------------------------------*/
267 static int
268 muse_idp_get_scales(const cpl_propertylist *aHeader, double *xscale, double *yscale)
269 {
270 
271  int flag = 0;
272 
273  double _xscale = 0.;
274  double _yscale = 0.;
275 
276  cpl_errorstate status = cpl_errorstate_get();
277 
278  if (muse_wcs_get_scales((cpl_propertylist *)aHeader, &_xscale, &_yscale) != CPL_ERROR_NONE) {
279  const char *insmode = muse_pfits_get_insmode(aHeader);
280  if (insmode) {
281  if (strncmp(insmode, "NFM", 3) == 0) {
282  _xscale = kMusePixscaleNFM;
283  _yscale = kMusePixscaleNFM;
284  } else {
285  _xscale = kMusePixscaleWFM;
286  _yscale = kMusePixscaleWFM;
287  }
288  flag = 1;
289  } else {
290  *xscale = _xscale;
291  *yscale = _yscale;
292  flag = -1;
293  }
294  } else {
295  _xscale *= deg2as;
296  _yscale *= deg2as;
297  }
298 
299  cpl_errorstate_set(status);
300 
301  if (flag == -1) {
302  cpl_msg_debug(__func__, "Keyword \"INS.MODE\" is missing!");
303  cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
304  } else {
305  *xscale = _xscale;
306  *yscale = _yscale;
307  }
308 
309  return flag;
310 }
311 
312 
313 /*----------------------------------------------------------------------------*/
326 /*----------------------------------------------------------------------------*/
327 static double
328 muse_idp_compute_specres(double aLambda, const muse_specres_lut_entry *aSpecresLut)
329 {
330 
331  /* Locate the appropriate wavelength bin in the lut and compute the *
332  * spectral by linear interpolation. */
333 
334  cpl_size i = 0;
335  while ((aSpecresLut[i].lambda > 0.) && (aLambda > aSpecresLut[i].lambda)) {
336  ++i;
337  }
338 
339  if (i == 0) {
340  return aSpecresLut[0].R;
341  }
342  else if (aSpecresLut[i].lambda == 0.) {
343  return aSpecresLut[i - 1].R;
344  }
345 
346  double t = (aLambda - aSpecresLut[i - 1].lambda) /
347  (aSpecresLut[i].lambda - aSpecresLut[i - 1].lambda);
348 
349  return (1 - t) * aSpecresLut[i - 1].R + t * aSpecresLut[i].R;
350 }
351 
352 /*----------------------------------------------------------------------------*/
366 /*----------------------------------------------------------------------------*/
367 static double
368 muse_idp_compute_iqe(double aLambda, double aSeeing, double aAirmass)
369 {
370 
371  const double rad2as = 180. / CPL_MATH_PI * 3600.; /* radians to arcsec */
372 
373 
374  /* Wavefront outer scale at Paranal [meters] */
375  const double L0 = 23.;
376 
377  /* UT M1 diameter [meters] from UT M1 vignetted area */
378  const double D = 8.1175365721399437;
379 
380  /* Kolb factor (cf. ESO Technical Report 12) */
381  const double Fkolb = 1. / (1. + 300. * D / L0) - 1.;
382 
383 
384  /* Scaled wavelength: lambda over lambda_ref */
385  double lambda = aLambda / 5000.;
386 
387  /* Fried parameter */
388  double r0 = 0.976 * 5.e-7 / aSeeing * rad2as *
389  pow(lambda, 1.2) * pow(aAirmass, -0.6);
390 
391  double iqe = aSeeing * pow(lambda, -0.2) * pow(aAirmass, 0.6);
392  iqe *= sqrt(1. + Fkolb * 2.183 * pow(r0 / L0, 0.356));
393 
394  return iqe;
395 }
396 
397 /*----------------------------------------------------------------------------*/
412 /*----------------------------------------------------------------------------*/
413 static double
414 muse_idp_compute_pixnoise(const muse_datacube *aCube,
415  double aVarianceMin, double aVarianceMax,
416  double aVarianceBin)
417 {
418  cpl_ensure(aCube != NULL, CPL_ERROR_NULL_INPUT, 0.);
419  cpl_ensure(aCube->stat, CPL_ERROR_DATA_NOT_FOUND, 0.);
420  cpl_ensure((aVarianceMin >= 0.) && (aVarianceMax > aVarianceMin),
421  CPL_ERROR_ILLEGAL_INPUT, 0.);
422  cpl_ensure(aVarianceBin > 0., CPL_ERROR_ILLEGAL_INPUT, 0.);
423 
424  cpl_size nbins = (cpl_size)((aVarianceMax - aVarianceMin) / aVarianceBin) + 1;
425  if (nbins < 2) {
426  cpl_error_set(__func__, CPL_ERROR_INCOMPATIBLE_INPUT);
427  return 0.;
428  }
429 
430  const cpl_imagelist *data = aCube->stat;
431  cpl_matrix *histogram = cpl_matrix_new(1, nbins);
432  double *hdata = cpl_matrix_get_data(histogram);
433 
434  cpl_size iplane;
435  for (iplane = 0; iplane < cpl_imagelist_get_size(data); ++iplane) {
436  const cpl_image *plane = cpl_imagelist_get_const(data, iplane);
437  const float *_data = cpl_image_get_data_float_const(plane);
438  cpl_size ndata = cpl_image_get_size_x(plane) * cpl_image_get_size_y(plane);
439  cpl_size idata;
440  for (idata = 0; idata < ndata; ++idata) {
441  register double v = _data[idata];
442  if ((isnan(v) == 0) && (v >= aVarianceMin) && (v <= aVarianceMax)) {
443  cpl_size ibin = (cpl_size)((v + aVarianceMin) / aVarianceBin);
444  hdata[ibin] += 1.;
445  }
446  }
447  }
448 
449  /* Get position of the histogram maximum. Since the histogram is a 1D *
450  * vector, the row position is just a dummy. */
451  cpl_size irow = 0.;
452  cpl_size icol = 0.;
453  cpl_matrix_get_maxpos(histogram, &irow, &icol);
454 
455  /* Use center of the histogram bin */
456  double mode = aVarianceMin + (icol + 0.5) * aVarianceBin;
457  cpl_matrix_delete(histogram);
458 
459  return mode;
460 }
461 
462 /*----------------------------------------------------------------------------*/
472 /*----------------------------------------------------------------------------*/
473 static double
474 muse_idp_compute_abmaglimit(const muse_datacube *aCube, double aWlenMin,
475  double aWlenMax, double aFwhm)
476 {
477  cpl_ensure(aCube != NULL, CPL_ERROR_NULL_INPUT, 0.);
478  cpl_ensure(aWlenMax > aWlenMin, CPL_ERROR_ILLEGAL_INPUT, 0.);
479  cpl_ensure(aFwhm > 0., CPL_ERROR_ILLEGAL_INPUT, 0.);
480 
481  /* Minimum number of valid data pixel needed to perform the computation *
482  * of the magnitude limit. This is given by using the standard deviation *
483  * in the computation of the limiting magnitude */
484  const cpl_size min_valid = 2;
485 
486  /* Create white light image and encode all invalid pixels as NaN for *
487  * the following analysis. */
488  muse_image *fov = muse_datacube_collapse((muse_datacube *)aCube, NULL);
490 
491  cpl_size nx = cpl_image_get_size_x(fov->data);
492  cpl_size ny = cpl_image_get_size_y(fov->data);
493 
494  /* Get median pixel value ignoring bad pixel values (NaN) */
495  cpl_image_reject_value(fov->data, CPL_VALUE_NOTFINITE);
496 #undef USE_HISTOGRAM
497 #ifndef USE_HISTOGRAM
498  double median = cpl_image_get_median(fov->data);
499 #else
500  double sigma = 0;
501  double median = cpl_image_get_mad(fov->data, &sigma);
502 #endif
503 
504 
505  /* Create cleaned FOV image from the original image using only *
506  * pixels which have less intensity than the median. All invalid *
507  * pixels are recorded in a mask image. To prepare for the following *
508  * convolution step the image is also cropped to an even number of *
509  * pixels along both axes (if needed). */
510  cpl_size _nx = (nx % 2) == 0 ? nx : nx - 1;
511  cpl_size _ny = (ny % 2) == 0 ? ny : ny - 1;
512  cpl_size npixel = _nx * _ny;
513 
514  cpl_image *fov_clean = cpl_image_new(_nx, _ny, CPL_TYPE_DOUBLE);
515  cpl_image_fill_window(fov_clean, 1, 1, _nx, _ny, median);
516 
517  cpl_mask *fov_invalid = cpl_mask_new(_nx, _ny);
518  cpl_binary *_fov_invalid = cpl_mask_get_data(fov_invalid);
519 
520  const float *_fov_data = cpl_image_get_data_float(fov->data);
521  double *_fov_clean = cpl_image_get_data_double(fov_clean);
522 
523  cpl_size iy;
524  for (iy = 0; iy < _ny; ++iy) {
525  cpl_size ix;
526  for (ix = 0; ix < _nx; ++ix) {
527  cpl_size _ipixel = nx * iy + ix;
528  if ((isnan(_fov_data[_ipixel]) == 0) &&
529  (_fov_data[_ipixel] < median)) {
530  _fov_clean[_nx * iy + ix] = _fov_data[_ipixel];
531  } else {
532  _fov_invalid[_nx * iy + ix] = CPL_BINARY_1;
533  }
534 
535  }
536  }
537  muse_image_delete(fov);
538 
539  /* Convolve the valid image with a Gaussian normalized to 1 at the peak. *
540  * muse_matrix_new_gaussian_2d() creates a 2d Gaussian normalized such *
541  * the integral (the total sum of the elements) is 1. Thus it has to be *
542  * renormalized. The half width (radius) of the convolution box is set to *
543  * 3 sigma. */
544  const double csigma = aFwhm / (2. * sqrt(2. * log(2.)));
545  const double radius = CPL_MAX(3. * csigma, 2.001);
546  const int nhalf = (int)(radius + 0.5);
547  cpl_matrix *ckernel = muse_matrix_new_gaussian_2d(nhalf, nhalf, csigma);
548  cpl_matrix_divide_scalar(ckernel, cpl_matrix_get_max(ckernel));
549 
550  cpl_image *fov_smooth = muse_convolve_image(fov_clean, ckernel);
551 
552  cpl_matrix_delete(ckernel);
553  cpl_image_delete(fov_clean);
554 
555  /* Perform a morphological dilation on the mask to create an additional *
556  * border of masked pixels around the already masked areas. This should *
557  * invalidated pixels in the smoothed image which are affected by border *
558  * effects during the convolution. */
559 
560  cpl_size nmask = 2 * nhalf + 1;
561  cpl_mask *mkernel = cpl_mask_new(nmask, nmask);
562  cpl_binary *_mkernel = cpl_mask_get_data(mkernel);
563  double rsquare = radius * radius;
564 
565  for (iy = 0; iy < nmask; ++iy) {
566  cpl_size ix;
567  double y = (double)iy - radius;
568  for (ix = 0; ix < nmask; ++ix) {
569  double x = (double)ix - radius;
570  _mkernel[nmask * iy + ix] = ((x * x + y * y) > rsquare) ?
571  CPL_BINARY_0 : CPL_BINARY_1;
572  }
573  }
574 
575  cpl_mask *fov_mask = cpl_mask_new(_nx, _ny);
576  cpl_mask_filter(fov_mask, fov_invalid, mkernel,
577  CPL_FILTER_DILATION, CPL_BORDER_ZERO);
578  cpl_mask_delete(mkernel);
579 
580  /* The filtering of the mask leaves an nhalf pixel wide border on each *
581  * side of the mask image which are not flagged. These are set as masked *
582  * here. *
583  */
584  cpl_binary *_fov_mask = cpl_mask_get_data(fov_mask);
585 
586  cpl_size ipixel;
587  for (ipixel = 0; ipixel < npixel; ++ipixel) {
588  cpl_size _x = ipixel % _nx;
589  cpl_size _y = ipixel / _nx;
590  if (((_x < nhalf) ||(_x > (_nx - nhalf - 1))) ||
591  ((_y < nhalf) ||(_y > (_ny - nhalf - 1)))) {
592  _fov_mask[_y * _nx + _x] = CPL_BINARY_1;
593  }
594  }
595 
596  /* Check whether enough valid pixels are left for the following analysis. *
597  * If not enough valid pixels are found revert the dilated mask to the *
598  * original and use this instead. */
599  if (npixel - cpl_mask_count(fov_mask) < min_valid) {
600  cpl_msg_debug(__func__, "Number of valid pixels in the dilated mask is"
601  "less than 2! Falling back to the original mask!");
602  cpl_mask_delete(fov_mask);
603  fov_mask = fov_invalid;
604  fov_invalid = NULL;
605  } else {
606  cpl_mask_delete(fov_invalid);
607  }
608 
609 #ifndef USE_HISTOGRAM
610  const double *_fov_smooth = cpl_image_get_data_double_const(fov_smooth);
611 
612  double *_fov_valid = cpl_calloc(npixel, sizeof *_fov_valid);
613 
614  cpl_size nvalid = 0;
615  for (ipixel = 0; ipixel < npixel; ++ipixel) {
616  if (_fov_mask[ipixel] == CPL_BINARY_0) {
617  _fov_valid[nvalid] = _fov_smooth[ipixel];
618  ++nvalid;
619  }
620  }
621  cpl_mask_delete(fov_mask);
622  cpl_image_delete(fov_smooth);
623 
624  cpl_vector *fov_valid = cpl_vector_wrap(nvalid, _fov_valid);
625  double sdev = cpl_vector_get_stdev(fov_valid);
626 
627  cpl_vector_unwrap(fov_valid);
628  cpl_free(_fov_valid);
629 
630  /* Compute the AB limiting magnitude: Numerical values are from *
631  * converting f(nu) to f(lambda) and from the unit conversion [Jy] *
632  * to [erg/cm**2/s/Angstrom]: *
633  * f(nu) = lambda**2 / c * f(lambda), *
634  * f(nu)/Jy = 3.34e4 * (lambda/Angstrom)**2 * *
635  * f(lambda)/(erg/cm**2/s/Angstrom) *
636  * and the monochromatic AB magnitude *
637  * m_AB = -2.5 * log10(f(nu) / (3631 Jy)) */
638 
639  double zeropoint = -2.5 * log10(5. * (3.34e4 * aWlenMin * aWlenMax / 3631.));
640  double abmaglimit = -2.5 * log10(2. * sdev / kMuseFluxUnitFactor) + zeropoint;
641 #else
642  /* Create histogram from the smoothed image of valid pixels. The bin size *
643  * is chosen to be 0.1 times the median absolute deviation of the pixel *
644  * values of the original image */
645  const double hstep = 0.1 * sigma;
646  double hmin = cpl_image_get_min(fov_smooth);
647  double hmax = cpl_image_get_max(fov_smooth);
648 
649  cpl_size nbins = (hmax - hmin) / hstep + 1;
650  cpl_vector *bins = cpl_vector_new(nbins);
651 
652  /* To prepare for the following fit of the histogram peak use the bin *
653  * center as the value of the histogram bin */
654  cpl_size ibin;
655  for (ibin = 0; ibin < nbins; ++ ibin) {
656  cpl_vector_set(bins, ibin, hmin + (ibin + 0.5) * hstep);
657  }
658 
659  cpl_vector *histogram = cpl_vector_new(nbins);
660  cpl_vector_fill(histogram, 0.);
661  double *hdata = cpl_vector_get_data(histogram);
662 
663  const double *_fov_smooth = cpl_image_get_data_double_const(fov_smooth);
664 
665  for (ipixel = 0; ipixel < npixel; ++ipixel) {
666  if ((_fov_mask[ipixel] == CPL_BINARY_0) &&
667  (_fov_smooth[ipixel] >= hmin) && (_fov_smooth[ipixel] <= hmax)) {
668  cpl_size jbin = (cpl_size)((_fov_smooth[ipixel] + fabs(hmin)) / hstep);
669  hdata[jbin] += 1.;
670  }
671  }
672 
673  cpl_mask_delete(fov_mask);
674  cpl_image_delete(fov_smooth);
675 
676  /* Fit a Gaussian to the histogram peak */
677  double center = 0.;
678  double sdev = 0.;
679  double area = 0.;
680  double bckgrnd = 0.;
681 
682  cpl_fit_mode mode = CPL_FIT_CENTROID | CPL_FIT_STDEV | CPL_FIT_AREA;
683  cpl_errorstate status = cpl_errorstate_get();
684  cpl_error_code ecode = cpl_vector_fit_gaussian(bins, NULL, histogram, NULL,
685  mode,
686  &center, &sdev, &area, &bckgrnd,
687  NULL, NULL, NULL);
688  cpl_vector_delete(histogram);
689  cpl_vector_delete(bins);
690 
691  double abmaglimit = 0.;
692  if (!cpl_errorstate_is_equal(status) && (ecode != CPL_ERROR_CONTINUE)) {
693  cpl_msg_debug(__func__, "Fit of a Gaussian to the histogram peak failed!");
694  cpl_errorstate_set(status);
695  } else {
696  if (ecode == CPL_ERROR_CONTINUE) {
697  cpl_msg_debug(__func__, "Fit of a Gaussian to the histogram peak did "
698  "not converge!");
699  cpl_errorstate_set(status);
700  }
701 
702  /* Compute the AB limiting magnitude: Numerical values are from *
703  * converting f(nu) to f(lambda) and from the unit conversion [Jy] *
704  * to [erg/cm**2/s/Angstrom]: *
705  * f(nu) = lambda**2 / c * f(lambda), *
706  * f(nu)/Jy = 3.34e4 * (lambda/Angstrom)**2 * *
707  * f(lambda)/(erg/cm**2/s/Angstrom) *
708  * and the monochromatic AB magnitude *
709  * m_AB = -2.5 * log10(f(nu) / (3631 Jy)) */
710 
711  double zeropoint = -2.5 *
712  log10(5. * (3.34e4 * aWlenMin * aWlenMax / 3631.));
713  abmaglimit = -2.5 * log10(2. * sdev / kMuseFluxUnitFactor) + zeropoint;
714  }
715 #endif
716 
717  return abmaglimit;
718 }
719 
720 
721 
722 muse_idp_properties *
723 muse_idp_properties_new(void) {
724  muse_idp_properties *properties = cpl_calloc(1, sizeof *properties);
725  return properties;
726 }
727 
728 void
729 muse_idp_properties_delete(muse_idp_properties *aProperties)
730 {
731  if (aProperties) {
732  cpl_array_delete(aProperties->obid);
733  cpl_array_delete(aProperties->progid);
734  cpl_propertylist_delete(aProperties->prov);
735  cpl_array_delete(aProperties->asson);
736 
737  /* TODO: Remove deallocator call for assoc once it has been *
738  * removed from the interface. */
739  cpl_array_delete(aProperties->assoc);
740  cpl_free((char *)aProperties->prodcatg);
741  cpl_free((char *)aProperties->procsoft);
742  cpl_free((char *)aProperties->obstech);
743  cpl_free((char *)aProperties->referenc);
744  }
745  cpl_free(aProperties);
746  return;
747 }
748 
749 
750 muse_idp_properties *
751 muse_idp_properties_collect(muse_processing *aProcessing,
752  const muse_datacube *aCube,
753  const char *aTag)
754 {
755  cpl_ensure(aProcessing, CPL_ERROR_NULL_INPUT, NULL);
756  cpl_table *exposures = muse_processing_sort_exposures(aProcessing);
757  cpl_ensure(exposures, CPL_ERROR_DATA_NOT_FOUND, NULL);
758  cpl_size nexposures = cpl_table_get_nrow(exposures);
759 
760  muse_idp_properties *properties = muse_idp_properties_new();
761  properties->obid = cpl_array_new(nexposures, CPL_TYPE_LONG);
762  properties->progid = cpl_array_new(nexposures, CPL_TYPE_STRING);
763  properties->prov = cpl_propertylist_new();
764  properties->fluxcal = CPL_TRUE;
765  properties->wlerror = 0.026; /* Fallback random error estimate of wavelength */
766 
767  /* Collect data from the raw data files */
768 
769  cpl_errorstate status = cpl_errorstate_get();
770 
771  cpl_size kexposure = 0;
772  cpl_size iexposure;
773  double fwhmmin = DBL_MAX;
774  double fwhmmax = -DBL_MAX;
775 
776  for (iexposure = 0; iexposure < nexposures; ++iexposure) {
777  const char *filename = NULL;
778  if (cpl_table_has_column(exposures, "00")) {
779  filename = cpl_table_get_string(exposures, "00", iexposure);
780  }
781  if (!filename) {
782  int ifu=1;
783  while (!filename && (ifu <= kMuseNumIFUs)) {
784  char *column = cpl_sprintf("%02d", ifu);
785  filename = cpl_table_get_string(exposures, column, iexposure);
786  cpl_free(column);
787  ++ifu;
788  }
789  }
790  if (!filename) {
791  cpl_msg_error(__func__, "No pixel table is available for exposure %"
792  CPL_SIZE_FORMAT " of %" CPL_SIZE_FORMAT "!",
793  iexposure + 1, nexposures);
794  cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
795  break;
796  } else {
797 
798  cpl_msg_debug(__func__, "Collecting IDP data from exposure '%s'",
799  filename);
800  cpl_propertylist *_header = cpl_propertylist_load(filename, 0);
801 
802  const char *progid = muse_pfits_get_progid(_header);
803  long obid = muse_pfits_get_obsid(_header);
804  double mjd = muse_pfits_get_mjdobs(_header);
805  double exptime = muse_pfits_get_exptime(_header);
806  double endtime = mjd + exptime / day2sec;
807 
808  /* Estimate the number of illuminated pixels from the nominal size *
809  * of the field-of-view in pixels. */
810  double npixel = 0.;
811 
812  cpl_errorstate _status = cpl_errorstate_get();
813 
814  const char *_insmode = muse_pfits_get_insmode(_header);
815  if (_insmode) {
816  if (strncmp(_insmode, "NFM", 3) == 0) {
817  npixel = (kMuseFovSizeNFM * kMuseFovSizeNFM) /
818  (kMusePixscaleNFM * kMusePixscaleNFM);
819  } else {
820  npixel = (kMuseFovSizeWFM * kMuseFovSizeWFM) /
821  (kMusePixscaleWFM * kMusePixscaleWFM);
822  }
823  }
824 
825  if (!cpl_errorstate_is_equal(_status)) {
826  if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
827  cpl_msg_debug(__func__, "Keyword \"INS.MODE\" is missing in '%s'!",
828  filename);
829  }
830  cpl_msg_error(__func__, "Cannot determine pixel scales of exposure '%s'",
831  filename);
832  cpl_propertylist_delete(_header);
833  break;
834  }
835 
836  properties->exptime += exptime * npixel;
837  properties->texptime += exptime;
838  properties->mjd_end = CPL_MAX(endtime, properties->mjd_end);
839 
840  /* Get image quality (at the reference wavelength) and the associated *
841  * error for this exposure */
842  double zstart = muse_pfits_get_airmass_start(_header);
843  double zend = muse_pfits_get_airmass_end(_header);
844 
845  _status = cpl_errorstate_get();
846  double fwhmstart = muse_pfits_get_fwhm_start(_header);
847  double fwhmend = muse_pfits_get_fwhm_end(_header);
848 
849  if (!cpl_errorstate_is_equal(_status)) {
850  cpl_msg_debug(__func__, "Missing or invalid seeing information "
851  "(\"TEL.AMBI.FWHM.START\", \"TEL.AMBI.FWHM.END\") in "
852  "exposure '%s'", filename);
853  cpl_errorstate_set(_status);
854  }
855 
856  double iqe = kMuseSkyRes;
857  double fwhmIA = muse_pfits_get_ia_fwhmlin(_header);
858 
859  if (cpl_errorstate_is_equal(_status)) {
860  iqe = muse_idp_compute_iqe(kMuseLambdaRef, fwhmIA,
861  0.5 * (zstart + zend));
862  }
863  else {
864  /* Ignore the error and use the default constant instead. */
865  /* Set fwhmstart and fwhmend to 0. This forces fwhmmin to */
866  /* zero and the uncertainty of SKY_RES to become the */
867  /* default constant. */
868  cpl_msg_debug(__func__, "Missing or invalid seeing information "
869  "(\"TEL.IA.FWHMLIN\") in exposure '%s'", filename);
870  cpl_errorstate_set(_status);
871  fwhmstart = 0.;
872  fwhmend = 0.;
873  }
874 
875  properties->skyres += exptime * iqe;
876 
877  fwhmmin = CPL_MIN(fwhmmin, CPL_MIN(fwhmstart, fwhmend));
878  fwhmmax = CPL_MAX(fwhmmax, CPL_MAX(fwhmstart, fwhmend));
879 
880  /* Observation ID, program ID, etc. */
881  cpl_array_set_long(properties->obid, kexposure, obid);
882  if (progid) {
883  cpl_array_set_string(properties->progid, kexposure, progid);
884  }
885 
886  if ((muse_pfits_get_fluxcal(_header) == CPL_FALSE) &&
887  (cpl_frameset_find(aProcessing->usedframes, MUSE_TAG_STD_RESPONSE) == NULL)) {
888  properties->fluxcal = CPL_FALSE;
889  }
890 
891  unsigned int nraw = cpl_propertylist_get_size(properties->prov);
892  const char *prov = muse_pfits_get_ancestor(_header);
893  if (!prov) {
894  prov = muse_pfits_get_arcfile(_header);
895  if (!prov) {
896  prov = muse_pfits_get_origfile(_header);
897  }
898  }
899 
900  if (prov) {
901  char *name = cpl_sprintf("PROV%-u", ++nraw);
902  cpl_propertylist_append_string(properties->prov, name, prov);
903  cpl_free(name);
904  } else {
905  unsigned int iraw = 1;
906  prov = muse_pfits_get_raw_filename(_header, iraw);
907  while (prov) {
908  char *name = cpl_sprintf("PROV%-u", ++nraw);
909  cpl_propertylist_append_string(properties->prov, name, prov);
910  prov = muse_pfits_get_raw_filename(_header, ++iraw);
911  cpl_free(name);
912  }
913  }
914 
915  cpl_propertylist_delete(_header);
916  ++kexposure;
917  }
918  }
919  cpl_table_delete(exposures);
920 
921  if (!cpl_errorstate_is_equal(status)) {
922  muse_idp_properties_delete(properties);
923  return NULL;
924  }
925 
926  if (kexposure != nexposures) {
927  cpl_msg_warning(__func__, "Could not collect raw data information of all "
928  "exposures. Got %" CPL_SIZE_FORMAT " of %" CPL_SIZE_FORMAT
929  "!", kexposure, nexposures);
930  }
931 
932  properties->ncombine = kexposure;
933  cpl_array_set_size(properties->obid, kexposure);
934  cpl_array_set_size(properties->progid, kexposure);
935 
936 
937  /* Collect ancillary products data */
938  cpl_msg_debug(__func__, "Collecting IDP ancillary products information");
939 
940  cpl_size nframes = cpl_frameset_get_size(aProcessing->outframes);
941  properties->asson = cpl_array_new(nframes, CPL_TYPE_STRING);
942 
943  cpl_size iancillary = 0;
944  cpl_size nptotal = 0;
945  cpl_size iframe;
946  for (iframe = 0; iframe < nframes; ++iframe) {
947  const cpl_frame *frame = cpl_frameset_get_position(aProcessing->outframes,
948  iframe);
949  const char *tag = cpl_frame_get_tag(frame);
950 
951  if (strcmp(tag, MUSE_TAG_IMAGE_FOV) == 0) {
952  cpl_array_set_string(properties->asson, iancillary,
953  cpl_frame_get_filename(frame));
954 
955  /* Estimate the number of illuminated pixels of the combined *
956  * field-of-view. The maximum number of valid pixels across all *
957  * the created field-of-view images is used as this estimate. */
958 
959  muse_image *_fov = muse_fov_load(cpl_frame_get_filename(frame));
960 
961  if (_fov) {
962  float *pixels = cpl_image_get_data_float(_fov->data);
963  cpl_size npix = cpl_image_get_size_x(_fov->data) * cpl_image_get_size_y(_fov->data) ;
964 
965  cpl_size _nptotal = 0;
966  cpl_size ipix;
967  for (ipix = 0; ipix < npix; ++ipix) {
968  if (isnan(pixels[ipix]) == 0) {
969  ++_nptotal;
970  }
971  }
972  nptotal = CPL_MAX(nptotal, _nptotal);
973 
974  muse_image_delete(_fov);
975  }
976  ++iancillary;
977  }
978  }
979  cpl_array_set_size(properties->asson, iancillary);
980 
981  /* Estimate exposure time per pixel. If the number of illuminated pixels *
982  * cannot be computed, or if there are no "valid" pixels the exposure *
983  * time is forced to zero. */
984 
985  if (nptotal == 0) {
986  properties->exptime = 0.;
987  } else {
988  properties->exptime /= nptotal;
989  }
990 
991  /* Compute product image quality estimate and its uncertainty */
992  properties->skyres /= properties->texptime;
993  if ((fwhmmin > 0.) && (fwhmmax > 0.)) {
994  properties->skyrerr = 0.5 * fabs(fwhmmax - fwhmmin);
995  properties->skyrerr = sqrt(properties->skyrerr * properties->skyrerr +
996  kMuseExtraErrorIQE * kMuseExtraErrorIQE);
997  } else {
998  cpl_msg_warning(__func__, "Got invalid seeing information! Using an "
999  "estimate for \"SKY_RERR\"!");
1000  properties->skyrerr = kMuseSkyResError;
1001  }
1002 
1003 
1004  /* Collect data from the product header. */
1005 
1006  cpl_msg_debug(__func__, "Collecting IDP data from products");
1007 
1008  const cpl_propertylist *header = aCube->header;
1009 
1010  /* Get FOV center coordinates and the limits of wavelength range. */
1011  cpl_wcs *wcs = cpl_wcs_new_from_propertylist(header);
1012  int naxes = cpl_wcs_get_image_naxis(wcs);
1013  const cpl_array *naxis = cpl_wcs_get_image_dims(wcs);
1014 
1015  cpl_matrix *position = cpl_matrix_new(2, naxes);
1016  cpl_matrix_fill_column(position, 0.5 * cpl_array_get_int(naxis, 0, NULL), 0);
1017  cpl_matrix_fill_column(position, 0.5 * cpl_array_get_int(naxis, 1, NULL), 1);
1018  cpl_matrix_set(position, 0, 2, 1.);
1019  cpl_matrix_set(position, 1, 2, cpl_array_get_int(naxis, 2, NULL));
1020 
1021  cpl_matrix *world;
1022  cpl_array *flags;
1023  cpl_error_code _status = cpl_wcs_convert(wcs, position, &world, &flags,
1024  CPL_WCS_PHYS2WORLD);
1025 
1026  cpl_array_delete(flags);
1027  cpl_matrix_delete(position);
1028 
1029  if (_status != CPL_ERROR_NONE) {
1030  cpl_matrix_delete(world);
1031  cpl_wcs_delete(wcs);
1032  muse_idp_properties_delete(properties);
1033  return NULL;
1034  }
1035 
1036  properties->fovcenter[0] = cpl_matrix_get(world, 0, 0);
1037  properties->fovcenter[1] = cpl_matrix_get(world, 0, 1);
1038 
1039  // XXX: Check it there is anything extra that needs to be done for a
1040  // AWAV-LOG axis
1041  properties->wlenrange[0] = cpl_matrix_get(world, 0, 2) * m2nm;
1042  properties->wlenrange[1] = cpl_matrix_get(world, 1, 2) * m2nm;
1043  cpl_matrix_delete(world);
1044  cpl_wcs_delete(wcs);
1045 
1046  /* Estimate the spectral resolution at the wavelength range center */
1047  const muse_specres_lut_entry *specres_data = muse_specres_data_wfm;
1048 
1049  const char *insmode = muse_pfits_get_insmode(header);
1050  if (strncmp(insmode, "NFM", 3) == 0) {
1051  specres_data = muse_specres_data_nfm;
1052  }
1053 
1054  double wlencenter = 0.5 * (properties->wlenrange[0] + properties->wlenrange[1]);
1055  properties->specres = muse_idp_compute_specres(wlencenter * nm2Angstrom,
1056  specres_data);
1057 
1058  /* Estimate the pixel-to-pixel noise from the statistical error of the *
1059  * data cube data. */
1060 
1061  status = cpl_errorstate_get();
1062 
1063  /* The arguments of the following call refer to variances and not sigmas. *
1064  * Thus, it considers uncertainties in the range -5 to 5 sigma, and a *
1065  * histogram bin size of 0.2 sigma. */
1066 
1067  double pixnoise = muse_idp_compute_pixnoise(aCube, 0, 25., 0.04);
1068  if (!cpl_errorstate_is_equal(status)) {
1069  double ron = 0.;
1070  double gain = 0.;
1071 
1072  unsigned int iquadrant;
1073  for (iquadrant = 1; iquadrant <= kMuseNumQuadrants; ++iquadrant) {
1074  ron += muse_pfits_get_ron(header, iquadrant);
1075  gain += muse_pfits_get_gain(header, iquadrant);
1076  }
1077 
1078  /* Here kMuseNumQuadrants^2 results from computing the mean ron and gain */
1079  properties->pixnoise = ron * gain * pow(10., -0.4 * kMuseResponseRef) /
1080  (double)(kMuseNumQuadrants * kMuseNumQuadrants);
1081  cpl_errorstate_set(status);
1082  } else {
1083  const char *bunit = muse_pfits_get_bunit(header);
1084 
1085  properties->pixnoise = sqrt(pixnoise);
1086 
1087  if (bunit &&
1088  !strncmp(bunit, kMuseFluxUnitString, strlen(kMuseFluxUnitString) + 1)) {
1089  properties->pixnoise /= kMuseFluxUnitFactor;
1090  } else {
1091  if (!bunit) {
1092  cpl_msg_warning(__func__,
1093  "No BUNIT data units are present in the data cube!");
1094  } else {
1095  cpl_msg_warning(__func__,
1096  "Data cube has unsupported data unit in BUNIT!");
1097  }
1098  cpl_msg_warning(__func__, "No scale factor has been applied to PIXNOISE!");
1099  }
1100  }
1101 
1102  /* Compute limiting magnitude */
1103  double xscale = 0.;
1104  double yscale = 0.;
1105 
1106  if (muse_idp_get_scales(header, &xscale, &yscale) == -1) {
1107  cpl_msg_error(__func__, "Cannot determine pixel scales of the data cube");
1108  muse_idp_properties_delete(properties);
1109  return NULL;
1110  }
1111 
1112  status = cpl_errorstate_get();
1113  double abmaglimit =
1114  muse_idp_compute_abmaglimit(aCube,
1115  properties->wlenrange[0] * nm2Angstrom,
1116  properties->wlenrange[1] * nm2Angstrom,
1117  properties->skyres / CPL_MAX(xscale, yscale));
1118  if (!cpl_errorstate_is_equal(status)) {
1119  properties->abmaglimit = -99.;
1120  cpl_errorstate_set(status);
1121  } else {
1122  properties->abmaglimit = abmaglimit;
1123  }
1124 
1125  /* Get the product category from the product header, unless it is *
1126  * overridden an explicitly given tag. */
1127  if (!aTag) {
1128  aTag = muse_pfits_get_pro_catg(header);
1129  }
1130  if (!aTag || (strcmp(aTag, MUSE_TAG_CUBE_FINAL) != 0)) {
1131  cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
1132  muse_idp_properties_delete(properties);
1133  return NULL;
1134  } else {
1135  properties->prodcatg = cpl_strdup(KEY_PRODCATG_VALUE_IFS_CUBE);
1136  }
1137 
1138  /* Get the most recently applied processing recipe of the product */
1139  unsigned int recid = 0;
1140  const char *procsoft = NULL;
1141  const char *pipeid;
1142  while ((pipeid = muse_pfits_get_pipe_id(header, ++recid)) != NULL) {
1143  procsoft = pipeid;
1144  }
1145 
1146  if (!procsoft) {
1147  properties->procsoft = cpl_strdup(PACKAGE "/" PACKAGE_VERSION);
1148  } else {
1149  properties->procsoft = cpl_strdup(procsoft);
1150  }
1151 
1152  /* The following IDP properties are simply hardcoded here since *
1153  * not supposed to change. */
1154  properties->obstech = cpl_strdup("IFU");
1155  properties->referenc = NULL;
1156 
1157  return properties;
1158 }
1159 
1160 cpl_error_code
1161 muse_idp_properties_update(cpl_propertylist *aHeader,
1162  const muse_idp_properties *aProperties)
1163 {
1164  cpl_ensure(aHeader && aProperties, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1165  cpl_ensure(cpl_array_get_size(aProperties->obid) == aProperties->ncombine,
1166  CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_ILLEGAL_INPUT);
1167  cpl_ensure(cpl_array_get_size(aProperties->progid) == aProperties->ncombine,
1168  CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_ILLEGAL_INPUT);
1169  cpl_ensure(cpl_propertylist_get_size(aProperties->prov) >= aProperties->ncombine,
1170  CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_ILLEGAL_INPUT);
1171 
1172  cpl_propertylist_erase_regexp(aHeader, CLEAN_KEYS_REGEXP, 0);
1173 
1174  cpl_propertylist_update_double(aHeader, KEY_RA, aProperties->fovcenter[0]);
1175  cpl_propertylist_set_comment(aHeader, KEY_RA, KEY_RA_COMMENT);
1176  cpl_propertylist_update_double(aHeader, KEY_DEC, aProperties->fovcenter[1]);
1177  cpl_propertylist_set_comment(aHeader, KEY_DEC, KEY_DEC_COMMENT);
1178  cpl_propertylist_update_double(aHeader, KEY_EXPTIME, aProperties->exptime);
1179  cpl_propertylist_set_comment(aHeader, KEY_EXPTIME, KEY_EXPTIME_COMMENT);
1180  cpl_propertylist_insert_after_double(aHeader, KEY_EXPTIME,
1181  KEY_TEXPTIME, aProperties->texptime);
1182  cpl_propertylist_set_comment(aHeader, KEY_TEXPTIME, KEY_TEXPTIME_COMMENT);
1183 
1184  cpl_propertylist_insert_after_int(aHeader, KEY_TEXPTIME,
1185  KEY_NCOMBINE, aProperties->ncombine);
1186  cpl_propertylist_set_comment(aHeader, KEY_NCOMBINE, KEY_NCOMBINE_COMMENT);
1187 
1188  cpl_propertylist_set_comment(aHeader, KEY_MJDOBS, KEY_MJDOBS_COMMENT);
1189  cpl_propertylist_insert_after_double(aHeader, KEY_MJDOBS,
1190  KEY_MJDEND, aProperties->mjd_end);
1191  cpl_propertylist_set_comment(aHeader, KEY_MJDEND, KEY_MJDEND_COMMENT);
1192 
1193  cpl_array *obids = cpl_array_duplicate(aProperties->obid);
1194  muse_cplarray_sort(obids, CPL_TRUE);
1195 
1196  /* Add observation block IDs, skipping duplicates */
1197  long obid = cpl_array_get_long(obids, 0, NULL);
1198  cpl_propertylist_update_long(aHeader, KEY_OBID "1", obid);
1199  cpl_propertylist_set_comment(aHeader, KEY_OBID "1", KEY_OBID_COMMENT);
1200 
1201  if (aProperties->ncombine > 1) {
1202  unsigned int ik = 1;
1203  cpl_size idx;
1204  for (idx = 1; idx < aProperties->ncombine; ++idx) {
1205  long _obid = cpl_array_get_long(obids, idx, NULL);
1206  if (_obid != obid) {
1207  char *key = cpl_sprintf(KEY_OBID "%-u", ++ik);
1208  cpl_propertylist_update_long(aHeader, key, _obid);
1209  cpl_propertylist_set_comment(aHeader, key, KEY_OBID_COMMENT);
1210  cpl_free(key);
1211  obid = _obid;
1212  }
1213  }
1214  }
1215  cpl_array_delete(obids);
1216 
1217  /* Add the ESO program ID(s). If the product is made from a single program *
1218  * PROG_ID is the ID string, otherwise it set to the special value MULTI, *
1219  * followed by a list of all IDs, skipping duplicate ID values */
1220  cpl_array *progids = cpl_array_duplicate(aProperties->progid);
1221  muse_cplarray_sort(progids, CPL_TRUE);
1222  const char *progid = cpl_array_get_string(progids, 0);
1223 
1224  if (aProperties->ncombine > 1) {
1225  unsigned int nprogid = 1;
1226  cpl_size idx;
1227  for (idx = 1; idx < aProperties->ncombine; ++idx) {
1228  const char *_progid = cpl_array_get_string(progids, idx);
1229  if (strcmp(_progid, progid) != 0) {
1230  progid = _progid;
1231  ++nprogid;
1232  }
1233  }
1234 
1235  progid = cpl_array_get_string(progids, 0);
1236  if (nprogid == 1) {
1237  cpl_propertylist_update_string(aHeader, KEY_PROG_ID, progid);
1238  } else {
1239  cpl_propertylist_update_string(aHeader, KEY_PROG_ID,
1240  KEY_PROG_ID_VALUE_MULTIPLE);
1241  cpl_propertylist_update_string(aHeader, KEY_PROGID "1", progid);
1242  cpl_propertylist_set_comment(aHeader, KEY_PROGID "1", KEY_PROGID_COMMENT);
1243 
1244  unsigned int ik = 1;
1245  for (idx = 1; idx < aProperties->ncombine; ++idx) {
1246  const char *_progid = cpl_array_get_string(progids, idx);
1247  if (strcmp(_progid, progid) != 0) {
1248  char *key = cpl_sprintf(KEY_PROGID "%-u", ++ik);
1249  cpl_propertylist_update_string(aHeader, key, _progid);
1250  cpl_propertylist_set_comment(aHeader, key, KEY_PROGID_COMMENT);
1251  cpl_free(key);
1252  progid = _progid;
1253  }
1254  }
1255  }
1256  cpl_propertylist_set_comment(aHeader, KEY_PROG_ID, KEY_PROG_ID_COMMENT);
1257 
1258  } else {
1259  cpl_propertylist_update_string(aHeader, KEY_PROG_ID, progid);
1260  cpl_propertylist_set_comment(aHeader, KEY_PROG_ID, KEY_PROG_ID_COMMENT);
1261  }
1262  cpl_array_delete(progids);
1263 
1264  /* Add raw data information */
1265  cpl_propertylist_append(aHeader, aProperties->prov);
1266 
1267  /* Add ancillary products file information */
1268  cpl_size idx;
1269  for (idx = 0; idx < cpl_array_get_size(aProperties->asson); ++idx) {
1270  char *name = cpl_sprintf(KEY_ASSON "%-" CPL_SIZE_FORMAT, idx + 1);
1271  cpl_propertylist_update_string(aHeader, name,
1272  cpl_array_get_string(aProperties->asson, idx));
1273  cpl_free(name);
1274  }
1275 
1276  cpl_propertylist_update_string(aHeader, KEY_PRODCATG, aProperties->prodcatg);
1277  cpl_propertylist_set_comment(aHeader, KEY_PRODCATG, KEY_PRODCATG_COMMENT);
1278 
1279  cpl_propertylist_update_string(aHeader, KEY_PROCSOFT, aProperties->procsoft);
1280  cpl_propertylist_set_comment(aHeader, KEY_PROCSOFT, KEY_PROCSOFT_COMMENT);
1281 
1282  cpl_propertylist_update_string(aHeader, KEY_OBSTECH, aProperties->obstech);
1283  cpl_propertylist_set_comment(aHeader, KEY_OBSTECH, KEY_OBSTECH_COMMENT);
1284 
1285  const char *fluxcal = KEY_FLUXCAL_VALUE_TRUE;
1286  if (aProperties->fluxcal == CPL_FALSE) {
1287  fluxcal = KEY_FLUXCAL_VALUE_FALSE;
1288  }
1289  cpl_propertylist_update_string(aHeader, KEY_FLUXCAL, fluxcal);
1290  cpl_propertylist_set_comment(aHeader, KEY_FLUXCAL, KEY_FLUXCAL_COMMENT);
1291 
1292  cpl_propertylist_insert_after_double(aHeader, KEY_FLUXCAL, KEY_WAVELMIN,
1293  aProperties->wlenrange[0]);
1294  cpl_propertylist_set_comment(aHeader, KEY_WAVELMIN, KEY_WAVELMIN_COMMENT);
1295  cpl_propertylist_insert_after_double(aHeader, KEY_WAVELMIN, KEY_WAVELMAX,
1296  aProperties->wlenrange[1]);
1297  cpl_propertylist_set_comment(aHeader, KEY_WAVELMAX, KEY_WAVELMAX_COMMENT);
1298  cpl_propertylist_insert_after_double(aHeader, KEY_WAVELMAX, KEY_SPEC_RES,
1299  aProperties->specres);
1300  cpl_propertylist_set_comment(aHeader, KEY_SPEC_RES, KEY_SPEC_RES_COMMENT);
1301  cpl_propertylist_insert_after_double(aHeader, KEY_SPEC_RES, KEY_SKY_RES,
1302  aProperties->skyres);
1303  cpl_propertylist_set_comment(aHeader, KEY_SKY_RES, KEY_SKY_RES_COMMENT);
1304  cpl_propertylist_insert_after_double(aHeader, KEY_SKY_RES, KEY_SKY_RERR,
1305  aProperties->skyrerr);
1306  cpl_propertylist_set_comment(aHeader, KEY_SKY_RERR, KEY_SKY_RERR_COMMENT);
1307  cpl_propertylist_insert_after_double(aHeader, KEY_SKY_RERR, KEY_PIXNOISE,
1308  aProperties->pixnoise);
1309  cpl_propertylist_set_comment(aHeader, KEY_PIXNOISE, KEY_PIXNOISE_COMMENT);
1310  cpl_propertylist_insert_after_double(aHeader, KEY_PIXNOISE, KEY_ABMAGLIM,
1311  aProperties->abmaglimit);
1312  cpl_propertylist_set_comment(aHeader, KEY_ABMAGLIM, KEY_ABMAGLIM_COMMENT);
1313 
1314  const char *reference = "";
1315  if (aProperties->referenc) {
1316  reference = aProperties->referenc;
1317  }
1318  cpl_propertylist_update_string(aHeader, KEY_REFERENC, reference);
1319  cpl_propertylist_set_comment(aHeader, KEY_REFERENC, KEY_REFERENC_COMMENT);
1320 
1321  cpl_propertylist_insert_after_double(aHeader, "CRVAL3", KEY_SPEC_ERR,
1322  aProperties->wlerror);
1323  cpl_propertylist_set_comment(aHeader, KEY_SPEC_ERR, KEY_SPEC_ERR_COMMENT);
1324 
1325  /* Fallback required for IDPs */
1326 
1327  if (!cpl_propertylist_has(aHeader, "CSYER1")) {
1328  cpl_propertylist_update_double(aHeader, "CSYER1", -1.);
1329  cpl_propertylist_set_comment(aHeader, "CSYER1",
1330  "[deg] Systematic error in coordinate");
1331  }
1332  if (!cpl_propertylist_has(aHeader, "CSYER2")) {
1333  cpl_propertylist_update_double(aHeader, "CSYER2", -1.);
1334  cpl_propertylist_set_comment(aHeader, "CSYER2",
1335  "[deg] Systematic error in coordinate");
1336  }
1337 
1338  return CPL_ERROR_NONE;
1339 }
1340 
1341 cpl_error_code
1342 muse_idp_properties_update_fov(muse_image *aFov)
1343 {
1344  cpl_error_code status = CPL_ERROR_NONE;
1345 
1346  cpl_ensure(aFov, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1347  cpl_ensure(aFov->header, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_ILLEGAL_INPUT);
1348 
1349  status = cpl_propertylist_update_string(aFov->header, KEY_PRODCATG,
1350  KEY_PRODCATG_VALUE_FOV_IMAGE);
1351  cpl_propertylist_set_comment(aFov->header, KEY_PRODCATG, KEY_PRODCATG_COMMENT);
1352 
1353  return status;
1354 }
1355 
Structure definition of a MUSE datacube.
Definition: muse_datacube.h:48
cpl_error_code muse_wcs_get_scales(cpl_propertylist *aHeader, double *aXScale, double *aYScale)
Compute the spatial scales (in degrees) from the FITS header WCS.
Definition: muse_wcs.c:1825
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
Definition: muse_image.c:85
const char * muse_pfits_get_insmode(const cpl_propertylist *aHeaders)
find out the observation mode
Definition: muse_pfits.c:1383
#define MUSE_HDR_PT_FLUXCAL
double muse_pfits_get_gain(const cpl_propertylist *aHeaders, unsigned char aQuadrant)
find the detector gain (in units of count/adu)
Definition: muse_pfits.c:686
cpl_image * data
the data extension
Definition: muse_image.h:46
const char * muse_pfits_get_ancestor(const cpl_propertylist *aHeaders)
find out the ancestor of a file.
Definition: muse_pfits.c:1739
const char * muse_pfits_get_raw_filename(const cpl_propertylist *aHeaders, unsigned int idx)
find out the i-th raw file name.
Definition: muse_pfits.c:1718
double muse_pfits_get_airmass_start(const cpl_propertylist *aHeaders)
find out the airmass at start of exposure
Definition: muse_pfits.c:970
cpl_error_code muse_cplarray_sort(cpl_array *aArray, cpl_boolean aOrder)
Sort float, int or double array by quicksort.
muse_image * muse_fov_load(const char *aFilename)
Load a FOV image into a MUSE image.
Definition: muse_utils.c:2864
const char * muse_pfits_get_bunit(const cpl_propertylist *aHeaders)
find out the unit string
Definition: muse_pfits.c:223
Structure definition of MUSE three extension FITS file.
Definition: muse_image.h:40
cpl_propertylist * header
the FITS header
Definition: muse_image.h:72
cpl_frameset * usedframes
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
cpl_frameset * outframes
const char * muse_pfits_get_arcfile(const cpl_propertylist *aHeaders)
find out the arcfile
Definition: muse_pfits.c:54
double muse_pfits_get_ron(const cpl_propertylist *aHeaders, unsigned char aQuadrant)
find the detector read-out noise
Definition: muse_pfits.c:660
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
const char * muse_pfits_get_progid(const cpl_propertylist *aHeaders)
find out the ESO program identification
Definition: muse_pfits.c:1659
long muse_pfits_get_obsid(const cpl_propertylist *aHeaders)
find out the observation block id
Definition: muse_pfits.c:1641
cpl_propertylist * header
the FITS header
Definition: muse_datacube.h:57
double muse_pfits_get_airmass_end(const cpl_propertylist *aHeaders)
find out the airmass at end of exposure
Definition: muse_pfits.c:988
double muse_pfits_get_mjdobs(const cpl_propertylist *aHeaders)
find out the Julian Date of the observation
Definition: muse_pfits.c:346
cpl_image * muse_convolve_image(const cpl_image *aImage, const cpl_matrix *aKernel)
Compute the convolution of an image with a kernel.
Definition: muse_utils.c:2748
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
muse_image * muse_datacube_collapse(muse_datacube *aCube, const muse_table *aFilter)
Integrate a FITS NAXIS=3 datacube along the wavelength direction.
cpl_matrix * muse_matrix_new_gaussian_2d(int aXHalfwidth, int aYHalfwidth, double aSigma)
Create a matrix that contains a normalized 2D Gaussian.
Definition: muse_utils.c:1099
const char * muse_pfits_get_pipe_id(const cpl_propertylist *aHeaders, unsigned int idx)
find out the value of the flux calibration flag
Definition: muse_pfits.c:1696
cpl_table * muse_processing_sort_exposures(muse_processing *aProcessing)
Sort input frames (containing lists of pixel table filenames) into different exposures.
double muse_pfits_get_ia_fwhmlin(const cpl_propertylist *aHeaders)
find out the image analysis FWHM from a linear fit (in arcsec)
Definition: muse_pfits.c:1260
const char * muse_pfits_get_origfile(const cpl_propertylist *aHeaders)
find out the origfile
Definition: muse_pfits.c:71
cpl_error_code muse_image_dq_to_nan(muse_image *aImage)
Convert pixels flagged in the DQ extension to NANs in DATA (and STAT, if present).
Definition: muse_image.c:904
cpl_imagelist * stat
the cube containing the data variance
Definition: muse_datacube.h:86