MUSE Pipeline Reference Manual  1.0.2
muse_datacube.c
1 /* -*- Mode: C; tab-width: 2; indent-tabs-mode: nil; c-basic-offset: 2 -*- */
2 /* vim:set sw=2 sts=2 et cin: */
3 /*
4  * This file is part of the MUSE Instrument Pipeline
5  * Copyright (C) 2005-2014 European Southern Observatory
6  *
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_datacube.h"
34 
35 #include "muse_flux.h"
36 #include "muse_pfits.h"
37 #include "muse_quality.h"
38 #include "muse_utils.h"
39 #include "muse_wcs.h"
40 
41 /*----------------------------------------------------------------------------*/
48 /*----------------------------------------------------------------------------*/
49 
52 /*---------------------------------------------------------------------------*/
69 /*---------------------------------------------------------------------------*/
70 static double *
71 muse_datacube_collapse_filter_buffer(cpl_table *aFilter,
72  double aCRVAL, double aCRPIX, double aCDELT,
73  cpl_boolean aIsLog, int *aL1, int *aL2)
74 {
75  if (!aFilter || !aL1 || !aL2) {
76  return NULL;
77  }
78  double *fdata = (double *)cpl_calloc(*aL2, sizeof(double));
79 
80  int l;
81  for (l = 0; l < *aL2; l++) {
82  double lambda = (l + 1. - aCRPIX) * aCDELT + aCRVAL;
83  if (aIsLog) {
84  lambda = aCRVAL * exp((l + 1. - aCRPIX) * aCDELT / aCRVAL);
85  }
86  fdata[l] = muse_flux_response_interpolate(aFilter, lambda, NULL,
88  } /* for l (z / wavelengths) */
89 
90  /* try to find filter edges */
91  l = 0;
92  while (l < *aL2 && fabs(fdata[l]) < DBL_EPSILON) {
93  *aL1 = l++;
94  }
95  l = *aL2 - 1;
96  while (l > *aL1 && fabs(fdata[l]) < DBL_EPSILON) {
97  *aL2 = l--;
98  }
99  cpl_msg_debug(__func__, "filter wavelength range %.1f..%.1fA (planes %d..%d)",
100  (*aL1 + 1. - aCRPIX) * aCDELT + aCRVAL,
101  (*aL2 + 1. - aCRPIX) * aCDELT + aCRVAL,
102  *aL1, *aL2);
103 
104  return fdata;
105 } /* muse_datacube_collapse_filter_buffer() */
106 
107 /*---------------------------------------------------------------------------*/
128 /*---------------------------------------------------------------------------*/
129 muse_image *
130 muse_euro3dcube_collapse(muse_euro3dcube *aEuro3D, cpl_table *aFilter)
131 {
132  cpl_ensure(aEuro3D && aEuro3D->dtable && aEuro3D->hdata, CPL_ERROR_NULL_INPUT,
133  NULL);
134 
135  /* guess dimensions of the output image from the Euro3D data */
136  muse_wcs *wcs = muse_wcs_new(aEuro3D->header);
137  /* similar test as in muse_pixtable_wcs_check() */
138  wcs->iscelsph = CPL_FALSE;
139  const char *unitx = cpl_table_get_column_unit(aEuro3D->dtable, "XPOS"),
140  *unity = cpl_table_get_column_unit(aEuro3D->dtable, "YPOS");
141  cpl_ensure(unitx && unity, CPL_ERROR_DATA_NOT_FOUND, NULL);
142  /* should be equal in the first three characters */
143  if (!strncmp(unitx, unity, 4) && !strncmp(unitx, "deg", 4)) {
144  wcs->iscelsph = CPL_TRUE;
145  }
146  /* extreme coordinates in pixel or celestial coordinates */
147  double xmin = cpl_table_get_column_min(aEuro3D->dtable, "XPOS"),
148  xmax = cpl_table_get_column_max(aEuro3D->dtable, "XPOS"),
149  ymin = cpl_table_get_column_min(aEuro3D->dtable, "YPOS"),
150  ymax = cpl_table_get_column_max(aEuro3D->dtable, "YPOS");
151  double x1, x2, y1, y2; /* pixel coordinates */
152  if (wcs->iscelsph) {
153  wcs->crval1 /= CPL_MATH_DEG_RAD; /* convert to radians */
154  wcs->crval2 /= CPL_MATH_DEG_RAD;
155  muse_wcs_pixel_from_celestial_fast(wcs, xmin / CPL_MATH_DEG_RAD,
156  ymin / CPL_MATH_DEG_RAD, &x1, &y1);
157  muse_wcs_pixel_from_celestial_fast(wcs, xmax / CPL_MATH_DEG_RAD,
158  ymax / CPL_MATH_DEG_RAD, &x2, &y2);
159  } else {
160  x1 = xmin;
161  x2 = xmax;
162  y1 = ymin;
163  y2 = ymax;
164  }
165  int zmin = cpl_table_get_column_min(aEuro3D->dtable, "SPEC_STA"),
166  zmax = cpl_table_get_column_max(aEuro3D->dtable, "SPEC_STA"),
167  nx = lround(fabs(x2 - x1)) + 1,
168  ny = lround(fabs(y2 - y1)) + 1,
169  nz = zmax - zmin + 1;
170 
171  /* count how many valid elements really are in the most shifted spectrum */
172  cpl_size zmaxpos = -1;
173  cpl_table_get_column_maxpos(aEuro3D->dtable, "SPEC_STA", &zmaxpos);
174  const cpl_array *amax = cpl_table_get_array(aEuro3D->dtable, "DATA_SPE",
175  zmaxpos);
176  int l, nsize = cpl_array_get_size(amax);
177  for (l = nsize - 1; l > 0; l--) {
178  /* find last valid element */
179  if (cpl_array_is_valid(amax, l) == 1) {
180  break;
181  }
182  }
183  nz += l + 1; /* add the number of valid elements in spectral direction */
184  int nspec = cpl_table_get_nrow(aEuro3D->dtable);
185  cpl_msg_debug(__func__, "Euro3D dimensions: %dx%dx%d (z = %d - %d, valid %d),"
186  " %d spectra", nx, ny, nz, zmax, zmin, l + 1, nspec);
187 
188  /* resample the filter curve on the wavelength grid */
189  double crvals = cpl_propertylist_get_double(aEuro3D->hdata, "CRVALS"),
190  cdelts = cpl_propertylist_get_double(aEuro3D->hdata, "CDELTS");
191  const char *ctypes = cpl_propertylist_get_string(aEuro3D->hdata, "CTYPES");
192  cpl_boolean loglambda = ctypes && !strncmp(ctypes, "AWAV-LOG", 9);
193  cpl_msg_debug(__func__, "spectral WCS: %f / %f %s", crvals, cdelts,
194  loglambda ? "log" : "linear");
195  int l1 = 0, l2 = nz; /* loop boundaries for filter operation */
196  double *fdata = muse_datacube_collapse_filter_buffer(aFilter, crvals, zmin,
197  cdelts, loglambda, &l1, &l2);
198 
199  /* create output muse_image, but only with two image extensions */
200  muse_image *image = muse_image_new();
201  image->header = cpl_propertylist_duplicate(aEuro3D->header);
202  /* the axis3 WCS was already erased in muse_resampling_euro3d() */
203  image->data = cpl_image_new(nx, ny, CPL_TYPE_FLOAT);
204  float *outdata = cpl_image_get_data_float(image->data);
205  image->dq = cpl_image_new(nx, ny, CPL_TYPE_INT);
206  /* pre-fill with Euro3D status for missing data */
207  cpl_image_add_scalar(image->dq, EURO3D_MISSDATA);
208  int *outdq = cpl_image_get_data_int(image->dq);
209  /* image->stat remains NULL, variance information is lost by collapsing */
210 
211  /* check if we need to use the variance for weighting */
212  cpl_boolean usevariance = CPL_FALSE;
213  if (getenv("MUSE_COLLAPSE_USE_VARIANCE") &&
214  atoi(getenv("MUSE_COLLAPSE_USE_VARIANCE")) > 0) {
215  usevariance = CPL_TRUE;
216  }
217 
218  /* loop through all spaxels and over all wavelengths */
219  int k, noutside = 0;
220  #pragma omp parallel for default(none) private(l) /* as req. by Ralf */ \
221  shared(aEuro3D, fdata, l1, l2, noutside, nspec, nx, ny, outdata, \
222  outdq, usevariance, wcs)
223  for (k = 0; k < nspec; k++) {
224  int err;
225  double xpos = cpl_table_get(aEuro3D->dtable, "XPOS", k, &err),
226  ypos = cpl_table_get(aEuro3D->dtable, "YPOS", k, &err);
227  if (err) {
228  cpl_msg_warning(__func__, "spectrum %d in Euro3D table does not have "
229  "position information!", k + 1);
230  continue;
231  }
232  /* coordinate in the output image (indices starting at 0) */
233  double xpx, ypx;
234  if (wcs->iscelsph) {
235  muse_wcs_pixel_from_celestial_fast(wcs, xpos * CPL_MATH_RAD_DEG,
236  ypos * CPL_MATH_RAD_DEG, &xpx, &ypx);
237  } else {
238  muse_wcs_pixel_from_projplane_fast(wcs, xpos, ypos, &xpx, &ypx);
239  }
240  int i = lround(xpx) - 1,
241  j = lround(ypx) - 1;
242  if (i >= nx || j >= ny) {
243  /* should not happen, but skip this spectrum to avoid segfault */
244  noutside++;
245  continue;
246  }
247 
248  int nstart = cpl_table_get_int(aEuro3D->dtable, "SPEC_STA", k, &err);
249  const cpl_array *adata = cpl_table_get_array(aEuro3D->dtable, "DATA_SPE", k),
250  *adq = cpl_table_get_array(aEuro3D->dtable, "QUAL_SPE", k),
251  *astat = NULL;
252  if (usevariance) {
253  astat = cpl_table_get_array(aEuro3D->dtable, "STAT_SPE", k);
254  }
255 
256  double sumdata = 0., sumweight = 0.;
257  for (l = l1; l < l2; l++) {
258  /* array index for this wavelength */
259  int idx = l - nstart + 1;
260 
261  /* filter weight with fallback for operation without filter */
262  double fweight = fdata ? fdata[l] : 1.;
263  cpl_errorstate prestate = cpl_errorstate_get();
264  int dq = cpl_array_get_int(adq, idx, &err);
265  if (dq || err) { /* if pixel marked bad or inaccessible */
266  cpl_errorstate_set(prestate);
267  continue; /* ignore bad pixels of any kind */
268  }
269  /* make the variance information optional */
270  double variance = 1.;
271  if (usevariance) {
272  variance = astat ? cpl_array_get(astat, idx, &err) : 1.;
273  /* use ^2 as for Euro3D we store errors not variances directly */
274  variance *= variance;
275  if (err > 0) {
276  variance = DBL_MAX;
277  }
278  if (err || !isnormal(variance)) {
279  continue;
280  }
281  }
282  double data = cpl_array_get(adata, idx, &err);
283  /* weight by inverse variance */
284  sumdata += data * fweight / variance;
285  sumweight += fweight / variance;
286  } /* for l (z / wavelengths) */
287  if (isnormal(sumweight)) {
288  outdata[i + j*nx] = sumdata / sumweight;
289  outdq[i + j*nx] = EURO3D_GOODPIXEL;
290  } else { /* leave bad/missing pixels zero, but mark them bad */
291  outdq[i + j*nx] = EURO3D_MISSDATA;
292  }
293  } /* for k (all spaxels) */
294  cpl_free(wcs);
295  cpl_free(fdata); /* won't crash on NULL */
296  if (noutside > 0) {
297  cpl_msg_warning(__func__, "Skipped %d spaxels, due to their location "
298  "outside the recostructed image!", noutside);
299  }
300 
301  return image;
302 } /* muse_euro3dcube_collapse() */
303 
304 /*---------------------------------------------------------------------------*/
323 /*---------------------------------------------------------------------------*/
324 muse_image *
325 muse_datacube_collapse(muse_datacube *aCube, cpl_table *aFilter)
326 {
327  cpl_ensure(aCube && aCube->data && aCube->header, CPL_ERROR_NULL_INPUT, NULL);
328 
329  /* find in/output dimensions */
330  int nx = cpl_image_get_size_x(cpl_imagelist_get(aCube->data, 0)),
331  ny = cpl_image_get_size_y(cpl_imagelist_get(aCube->data, 0)),
332  nz = cpl_imagelist_get_size(aCube->data);
333  /* resample the filter curve on the wavelength grid */
334  double crpix3 = cpl_propertylist_get_double(aCube->header, "CRPIX3"),
335  crval3 = cpl_propertylist_get_double(aCube->header, "CRVAL3"),
336  cd33 = cpl_propertylist_get_double(aCube->header, "CD3_3");
337  const char *ctype3 = cpl_propertylist_get_string(aCube->header, "CTYPE3");
338  cpl_boolean loglambda = ctype3 && !strncmp(ctype3, "AWAV-LOG", 9);
339  int l1 = 0, l2 = nz; /* loop boundaries for filter operation */
340  double *fdata = muse_datacube_collapse_filter_buffer(aFilter, crval3, crpix3,
341  cd33, loglambda, &l1, &l2);
342 
343  /* create output muse_image, but only with two image extensions */
344  muse_image *image = muse_image_new();
345  image->header = cpl_propertylist_duplicate(aCube->header);
346  cpl_propertylist_erase_regexp(image->header, "^C...*3$|^CD3_.$", 0);
347  image->data = cpl_image_new(nx, ny, CPL_TYPE_FLOAT);
348  float *outdata = cpl_image_get_data_float(image->data);
349  image->dq = cpl_image_new(nx, ny, CPL_TYPE_INT);
350  int *outdq = cpl_image_get_data_int(image->dq);
351  /* image->stat remains NULL, variance information is lost by collapsing */
352 
353  /* check if we need to use the variance for weighting */
354  cpl_boolean usevariance = CPL_FALSE;
355  if (getenv("MUSE_COLLAPSE_USE_VARIANCE") &&
356  atoi(getenv("MUSE_COLLAPSE_USE_VARIANCE")) > 0) {
357  usevariance = CPL_TRUE;
358  }
359 
360  /* loop through all pixels in all wavelength planes */
361  int i;
362  #pragma omp parallel for default(none) /* as req. by Ralf */ \
363  shared(aCube, nx, ny, l1, l2, fdata, outdata, outdq, usevariance)
364  for (i = 0; i < nx; i++) {
365  int j;
366  for (j = 0; j < ny; j++) {
367  double sumdata = 0., sumweight = 0.;
368  int l;
369  for (l = l1; l < l2; l++) {
370  /* filter weight with fallback for operation without filter */
371  double fweight = fdata ? fdata[l] : 1.;
372  float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->data, l)),
373  *pstat = NULL;
374  if (!isfinite(pdata[i + j*nx])) {
375  continue;
376  }
377  if (aCube->dq) {
378  int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->dq, l));
379  if (pdq && pdq[i + j*nx]) {
380  continue; /* ignore bad pixels of any kind */
381  }
382  }
383  /* make the variance information optional */
384  double variance = 1.;
385  if (usevariance) {
386  pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->stat, l));
387  variance = pstat ? pstat[i + j*nx] : 1.;
388  if (!isnormal(variance)) {
389  continue;
390  }
391  }
392  /* weight by inverse variance */
393  sumdata += pdata[i + j*nx] * fweight / variance;
394  sumweight += fweight / variance;
395  } /* for l (z / wavelengths) */
396  if (isnormal(sumweight)) {
397  outdata[i + j*nx] = sumdata / sumweight;
398  outdq[i + j*nx] = EURO3D_GOODPIXEL;
399  } else { /* leave bad/missing pixels zero, but mark them bad */
400  outdq[i + j*nx] = EURO3D_MISSDATA;
401  }
402  } /* for j (y pixels) */
403  } /* for i (x pixels) */
404 
405  cpl_free(fdata); /* won't crash on NULL */
406 
407  return image;
408 } /* muse_datacube_collapse() */
409 
410 /*---------------------------------------------------------------------------*/
421 /*---------------------------------------------------------------------------*/
422 cpl_error_code
423 muse_datacube_save_recimages(const char *aFilename, muse_imagelist *aImages,
424  cpl_array *aNames)
425 {
426  cpl_ensure_code(aFilename, CPL_ERROR_NULL_INPUT);
427  cpl_error_code rc = CPL_ERROR_NONE;
428  if (!aImages || !aNames) {
429  return rc; /* this is a valid case, return without error */
430  }
431  unsigned int i, nimages = muse_imagelist_get_size(aImages);
432  for (i = 0; i < nimages; i++) {
433  muse_image *image = muse_imagelist_get(aImages, i);
434 
435  /* create small header with EXTNAME=filtername and WCS for images */
436  cpl_propertylist *header = cpl_propertylist_new();
437  char dataext[KEYWORD_LENGTH], obj[KEYWORD_LENGTH],
438  *dqext = NULL, *statext = NULL;
439  snprintf(dataext, KEYWORD_LENGTH, "%s", cpl_array_get_string(aNames, i));
440  if (image->dq) {
441  dqext = cpl_sprintf("%s_%s", cpl_array_get_string(aNames, i), EXTNAME_DQ);
442  }
443  if (image->stat) {
444  statext = cpl_sprintf("%s_%s", cpl_array_get_string(aNames, i), EXTNAME_STAT);
445  }
446  snprintf(obj, KEYWORD_LENGTH, "%s", cpl_array_get_string(aNames, i));
447  cpl_propertylist_append_string(header, "EXTNAME", dataext);
448  cpl_propertylist_set_comment(header, "EXTNAME", "reconstructed image (data values)");
449  muse_utils_copy_modified_header(image->header, header, "OBJECT", obj);
450  cpl_propertylist_update_string(header, MUSE_HDR_FILTER,
451  cpl_array_get_string(aNames, i));
452  cpl_propertylist_copy_property_regexp(header, image->header,
453  MUSE_WCS_KEYS, 0);
454  muse_utils_set_hduclass(header, "DATA", dataext, dqext, statext);
455  rc = cpl_image_save(image->data, aFilename, CPL_TYPE_UNSPECIFIED, header,
456  CPL_IO_EXTEND);
457 
458  if (image->dq) {
459  cpl_propertylist_update_string(header, "EXTNAME", dqext);
460  cpl_propertylist_set_comment(header, "EXTNAME", "reconstructed image (bad pixel status values)");
461  snprintf(obj, KEYWORD_LENGTH, "%s, %s", cpl_array_get_string(aNames, i),
462  EXTNAME_DQ);
463  muse_utils_copy_modified_header(image->header, header, "OBJECT", obj);
464  muse_utils_set_hduclass(header, "QUALITY", dataext, dqext, statext);
465  rc = cpl_image_save(image->dq, aFilename, CPL_TYPE_UNSPECIFIED, header,
466  CPL_IO_EXTEND);
467  }
468  if (image->stat) {
469  cpl_propertylist_update_string(header, "EXTNAME", statext);
470  cpl_propertylist_set_comment(header, "EXTNAME", "reconstructed image (variance)");
471  snprintf(obj, KEYWORD_LENGTH, "%s, %s", cpl_array_get_string(aNames, i),
472  EXTNAME_STAT);
473  muse_utils_copy_modified_header(image->header, header, "OBJECT", obj);
474  muse_utils_set_hduclass(header, "ERROR", dataext, dqext, statext);
475  rc = cpl_image_save(image->stat, aFilename, CPL_TYPE_UNSPECIFIED, header,
476  CPL_IO_EXTEND);
477  }
478 
479  cpl_propertylist_delete(header);
480  cpl_free(dqext);
481  cpl_free(statext);
482  } /* for i (all reconstructed images) */
483 
484  return rc;
485 } /* muse_datacube_save_recimages */
486 
487 /*---------------------------------------------------------------------------*/
502 /*---------------------------------------------------------------------------*/
503 cpl_error_code
504 muse_euro3dcube_save(muse_euro3dcube *aEuro3D, const char *aFilename)
505 {
506  cpl_ensure_code(aEuro3D && aFilename, CPL_ERROR_NULL_INPUT);
507 
508  /* save primary header and the data table in the first extension */
509  cpl_error_code rc = cpl_table_save(aEuro3D->dtable, aEuro3D->header,
510  aEuro3D->hdata, aFilename, CPL_IO_DEFAULT);
511  if (rc != CPL_ERROR_NONE) {
512  cpl_msg_warning(__func__, "failed to save data part of the Euro3D table: "
513  "%s", cpl_error_get_message());
514  }
515  /* save group header and data in the second extension */
516  rc = cpl_table_save(aEuro3D->gtable, NULL, aEuro3D->hgroup, aFilename,
517  CPL_IO_EXTEND);
518  if (rc != CPL_ERROR_NONE) {
519  cpl_msg_warning(__func__, "failed to save group part of the Euro3D table: "
520  "%s", cpl_error_get_message());
521  }
522 
523  rc = muse_datacube_save_recimages(aFilename, aEuro3D->recimages,
524  aEuro3D->recnames);
525  return rc;
526 } /* muse_euro3dcube_save() */
527 
528 /*---------------------------------------------------------------------------*/
538 /*---------------------------------------------------------------------------*/
539 void
541 {
542  /* if the euro3d object doesn't exists at all, we don't need to do anything */
543  if (!aEuro3D) {
544  return;
545  }
546 
547  /* checks for the existence of the sub-images *
548  * are done in the CPL functions */
549  cpl_table_delete(aEuro3D->dtable);
550  cpl_table_delete(aEuro3D->gtable);
551  cpl_propertylist_delete(aEuro3D->header);
552  cpl_propertylist_delete(aEuro3D->hdata);
553  cpl_propertylist_delete(aEuro3D->hgroup);
555  cpl_array_delete(aEuro3D->recnames);
556  cpl_free(aEuro3D);
557 } /* muse_euro3dcube_delete() */
558 
559 /*---------------------------------------------------------------------------*/
569 /*---------------------------------------------------------------------------*/
570 static cpl_error_code
571 muse_datacube_convert_dq_recimages(muse_datacube *aCube)
572 {
573  cpl_ensure_code(aCube, CPL_ERROR_NULL_INPUT);
574  if (!aCube->recimages) {
575  return CPL_ERROR_NONE; /* valid case, return without error */
576  }
577  unsigned int k, nimages = muse_imagelist_get_size(aCube->recimages);
578  for (k = 0; k < nimages; k++) {
579  muse_image *image = muse_imagelist_get(aCube->recimages, k);
580  if (!image->dq) { /* no point trying to convert... */
581  continue;
582  }
583  muse_image_dq_to_nan(image);
584  } /* for k (all reconstructed images) */
585 
586  return CPL_ERROR_NONE;
587 } /* muse_datacube_convert_dq_recimages() */
588 
589 /*---------------------------------------------------------------------------*/
602 /*---------------------------------------------------------------------------*/
603 cpl_error_code
605 {
606  cpl_ensure_code(aCube && aCube->data && aCube->stat && aCube->dq,
607  CPL_ERROR_NULL_INPUT);
608 
609  /* loop through all wavelength planes and then all pixels per plane */
610  int l, nx = cpl_image_get_size_x(cpl_imagelist_get(aCube->data, 0)),
611  ny = cpl_image_get_size_y(cpl_imagelist_get(aCube->data, 0)),
612  nz = cpl_imagelist_get_size(aCube->data);
613  #pragma omp parallel for default(none) /* as req. by Ralf */ \
614  shared(aCube, nx, ny, nz)
615  for (l = 0; l < nz; l++) {
616  int i;
617  for (i = 0; i < nx; i++) {
618  int j;
619  for (j = 0; j < ny; j++) {
620  float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->data, l)),
621  *pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->stat, l));
622  int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->dq, l));
623  if (pdq[i + j*nx] == EURO3D_GOODPIXEL) {
624  continue; /* nothing to do for good pixels */
625  }
626  /* set bad pixels of any type to not-a-number in both extensions */
627  pdata[i + j*nx] = NAN; /* supposed to be quiet NaN */
628  pstat[i + j*nx] = NAN;
629  } /* for j (y direction) */
630  } /* for i (x direction) */
631  } /* for l (wavelength planes) */
632 
633  /* deallocate DQ and set it to NULL to be able to check for it */
634  cpl_imagelist_delete(aCube->dq);
635  aCube->dq = NULL;
636 
637  /* do the same for the reconstructed images, if there are any */
638  muse_datacube_convert_dq_recimages(aCube);
639 
640  return CPL_ERROR_NONE;
641 } /* muse_datacube_convert_dq() */
642 
643 /*---------------------------------------------------------------------------*/
672 /*---------------------------------------------------------------------------*/
673 cpl_error_code
674 muse_datacube_save(muse_datacube *aCube, const char *aFilename)
675 {
676  cpl_ensure_code(aCube && aCube->header && aFilename, CPL_ERROR_NULL_INPUT);
677 
678  /* save headers into primary area */
679  cpl_propertylist *header = cpl_propertylist_new();
680  /* copy all headers, except the main WCS keys */
681  cpl_propertylist_copy_property_regexp(header, aCube->header,
682  MUSE_WCS_KEYS, 1);
683  cpl_error_code rc = cpl_propertylist_save(header, aFilename, CPL_IO_CREATE);
684  cpl_propertylist_delete(header);
685 
686  /* save data */
687  header = cpl_propertylist_new();
688  cpl_propertylist_append_string(header, "EXTNAME", EXTNAME_DATA);
689  cpl_propertylist_set_comment(header, "EXTNAME", EXTNAME_DATA_COMMENT);
690  muse_utils_copy_modified_header(aCube->header, header, "OBJECT",
691  EXTNAME_DATA);
692  cpl_propertylist_copy_property_regexp(header, aCube->header,
693  MUSE_WCS_KEYS"|^BUNIT", 0);
694  muse_utils_set_hduclass(header, "DATA", "DATA",
695  aCube->dq ? "DQ" : NULL, aCube->stat ? "STAT" : NULL);
696  rc = cpl_imagelist_save(aCube->data, aFilename, CPL_TYPE_FLOAT, header,
697  CPL_IO_EXTEND);
698  cpl_propertylist_delete(header);
699  /* save bad pixels, if available */
700  if (rc == CPL_ERROR_NONE && aCube->dq) {
701  header = cpl_propertylist_new();
702  cpl_propertylist_append_string(header, "EXTNAME", EXTNAME_DQ);
703  cpl_propertylist_set_comment(header, "EXTNAME", EXTNAME_DQ_COMMENT);
704  muse_utils_copy_modified_header(aCube->header, header, "OBJECT",
705  EXTNAME_DQ);
706  cpl_propertylist_copy_property_regexp(header, aCube->header,
707  MUSE_WCS_KEYS, 0);
708  muse_utils_set_hduclass(header, "QUALITY", "DATA", "DQ",
709  aCube->stat ? "STAT" : NULL);
710  rc = cpl_imagelist_save(aCube->dq, aFilename, CPL_TYPE_INT, header,
711  CPL_IO_EXTEND);
712  cpl_propertylist_delete(header);
713  }
714  /* save variance, if available */
715  if (rc == CPL_ERROR_NONE && aCube->stat) {
716  header = cpl_propertylist_new();
717  cpl_propertylist_append_string(header, "EXTNAME", EXTNAME_STAT);
718  cpl_propertylist_set_comment(header, "EXTNAME", EXTNAME_STAT_COMMENT);
719  /* add the correct BUNIT, depending on the data units */
720  if (!strncmp(cpl_propertylist_get_string(aCube->header, "BUNIT"),
721  kMuseFluxUnitString, strlen(kMuseFluxUnitString) + 1)) {
722  /* flux calibrated data */
723  cpl_propertylist_append_string(header, "BUNIT", kMuseFluxStatString);
724  }
725  muse_utils_copy_modified_header(aCube->header, header, "OBJECT",
726  EXTNAME_STAT);
727  cpl_propertylist_copy_property_regexp(header, aCube->header,
728  MUSE_WCS_KEYS, 0);
729  muse_utils_set_hduclass(header, "ERROR", "DATA", aCube->dq ? "DQ" : NULL,
730  "STAT");
731  rc = cpl_imagelist_save(aCube->stat, aFilename, CPL_TYPE_FLOAT, header,
732  CPL_IO_EXTEND);
733  cpl_propertylist_delete(header);
734  }
735 
736  rc = muse_datacube_save_recimages(aFilename, aCube->recimages, aCube->recnames);
737  return rc;
738 } /* muse_datacube_save() */
739 
740 /*---------------------------------------------------------------------------*/
741 /*
742  @private
743  @error{set CPL_ERROR_NULL_INPUT\, return NULL, aFilename is NULL}
744  @error{set CPL_ERROR_ILLEGAL_INPUT\, return NULL, aFilename does not exist}
745  @error{set CPL_ERROR_DATA_NOT_FOUND\, return NULL,
746  aFilename does not contain a DATA extension}
747  */
748 /*---------------------------------------------------------------------------*/
749 static cpl_propertylist *
750 muse_datacube_load_header(const char *aFilename)
751 {
752  cpl_ensure(aFilename, CPL_ERROR_NULL_INPUT, NULL);
753  int extdata = cpl_fits_find_extension(aFilename, "DATA");
754  cpl_ensure(extdata >= 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
755  cpl_ensure(extdata > 0, CPL_ERROR_DATA_NOT_FOUND, NULL);
756 
757  cpl_propertylist *header = cpl_propertylist_load(aFilename, 0);
758  cpl_propertylist *hdata = cpl_propertylist_load(aFilename, extdata);
759  cpl_propertylist_copy_property_regexp(header, hdata,
760  MUSE_WCS_KEYS"|BUNIT", 0);
761  cpl_propertylist_delete(hdata);
762  return header;
763 } /* muse_datacube_load_header() */
764 
765 /*---------------------------------------------------------------------------*/
786 /*---------------------------------------------------------------------------*/
788 muse_datacube_load(const char *aFilename)
789 {
790  cpl_ensure(aFilename, CPL_ERROR_NULL_INPUT, NULL);
791  muse_datacube *cube = cpl_calloc(1, sizeof(muse_datacube));
792  /* load primary header and merge in relevant entries from the DATA header */
793  cpl_errorstate state = cpl_errorstate_get();
794  cube->header = muse_datacube_load_header(aFilename);
795  if (!cpl_errorstate_is_equal(state) || !cube->header) {
796  cpl_msg_error(__func__, "Loading cube-like headers from \"%s\" failed!",
797  aFilename);
798  cpl_free(cube);
799  return NULL;
800  }
801  int ext = cpl_fits_find_extension(aFilename, "DATA");
802  cube->data = cpl_imagelist_load(aFilename, CPL_TYPE_UNSPECIFIED, ext);
803  /* DQ is usually not written to disk, but see if it's there and load it, if so */
804  ext = cpl_fits_find_extension(aFilename, "DQ");
805  if (ext > 0) {
806  cube->stat = cpl_imagelist_load(aFilename, CPL_TYPE_UNSPECIFIED, ext);
807  }
808  ext = cpl_fits_find_extension(aFilename, "STAT");
809  if (ext > 0) {
810  cube->stat = cpl_imagelist_load(aFilename, CPL_TYPE_UNSPECIFIED, ext);
811  }
812  int next = cpl_fits_count_extensions(aFilename);
813  while (++ext <= next) {
814  /* there is (one more) collapsed image, load it */
815  muse_image *image = muse_image_new();
816  image->header = cpl_propertylist_load(aFilename, ext);
817  image->data = cpl_image_load(aFilename, CPL_TYPE_UNSPECIFIED, 0, ext);
818  /* XXX forget about image_DQ and image_STAT for the moment, they are unlikely... */
819 
820  if (!cube->recnames) {
821  cube->recnames = cpl_array_new(1, CPL_TYPE_STRING);
822  } else {
823  cpl_array_set_size(cube->recnames, cpl_array_get_size(cube->recnames) + 1);
824  }
825  /* append name of this reconstructed image at the end of the names array */
826  cpl_array_set_string(cube->recnames, cpl_array_get_size(cube->recnames) - 1,
828  if (!cube->recimages) {
829  cube->recimages = muse_imagelist_new();
830  }
831  /* append the image to the end of the reconstructed image list */
832  muse_imagelist_set(cube->recimages, image,
834  }
835  return cube;
836 } /* muse_datacube_load() */
837 
838 /*---------------------------------------------------------------------------*/
848 /*---------------------------------------------------------------------------*/
849 void
851 {
852  /* if the cube does not exists at all, we don't need to do anything */
853  if (!aCube) {
854  return;
855  }
856 
857  /* checks for the existence of the sub-images *
858  * are done in the CPL functions */
859  cpl_imagelist_delete(aCube->data);
860  aCube->data = NULL;
861  cpl_imagelist_delete(aCube->dq);
862  aCube->dq = NULL;
863  cpl_imagelist_delete(aCube->stat);
864  aCube->stat = NULL;
865 
866  /* delete the FITS header, too */
867  cpl_propertylist_delete(aCube->header);
868  aCube->header = NULL;
869 
870  /* remove reconstructed image data, if present */
872  cpl_array_delete(aCube->recnames);
873  cpl_free(aCube);
874 } /* muse_datacube_delete() */
875 
Structure definition of a MUSE datacube.
Definition: muse_datacube.h:47
Structure definition for a collection of muse_images.
const char * muse_pfits_get_extname(const cpl_propertylist *aHeaders)
find out the extension name
Definition: muse_pfits.c:188
A structure containing a spatial two-axis WCS.
Definition: muse_wcs.h:93
cpl_error_code muse_euro3dcube_save(muse_euro3dcube *aEuro3D, const char *aFilename)
Save a Euro3D cube object to a file.
cpl_image * data
the data extension
Definition: muse_image.h:46
muse_datacube * muse_datacube_load(const char *aFilename)
Load header, DATA and optionally STAT and DQ extensions as well as the reconstructed images of a MUSE...
static void muse_wcs_pixel_from_projplane_fast(muse_wcs *aWCS, double aX, double aY, double *aXOut, double *aYOut)
Convert from projection plane coordinates to pixel coordinates.
Definition: muse_wcs.h:233
cpl_image * stat
the statistics extension
Definition: muse_image.h:64
void muse_datacube_delete(muse_datacube *aCube)
Deallocate memory associated to a muse_datacube object.
void muse_imagelist_delete(muse_imagelist *aList)
Free the memory of the MUSE image list.
cpl_propertylist * hgroup
the group FITS header
cpl_array * recnames
the reconstructed image filter names
muse_image * muse_datacube_collapse(muse_datacube *aCube, cpl_table *aFilter)
Integrate a FITS NAXIS=3 datacube along the wavelength direction.
muse_imagelist * recimages
the reconstructed image data
cpl_propertylist * header
the primary FITS header
Structure definition of MUSE three extension FITS file.
Definition: muse_image.h:40
cpl_array * recnames
the reconstructed image filter names
Definition: muse_datacube.h:70
cpl_propertylist * header
the FITS header
Definition: muse_image.h:72
cpl_error_code muse_utils_set_hduclass(cpl_propertylist *aHeader, const char *aClass2, const char *aExtData, const char *aExtDQ, const char *aExtStat)
Set HDU headers for the ESO FITS data format.
Definition: muse_utils.c:2530
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. ...
unsigned int muse_imagelist_get_size(muse_imagelist *aList)
Return the number of stored images.
cpl_image * dq
the data quality extension
Definition: muse_image.h:56
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:336
muse_wcs * muse_wcs_new(cpl_propertylist *aHeader)
Create a new WCS structure from a given FITS header.
Definition: muse_wcs.c:1537
cpl_error_code muse_datacube_convert_dq(muse_datacube *aCube)
Convert the DQ extension of a datacube to NANs in DATA and STAT.
muse_image * muse_imagelist_get(muse_imagelist *aList, unsigned int aIdx)
Get the muse_image of given list index.
static void muse_wcs_pixel_from_celestial_fast(muse_wcs *aWCS, double aRA, double aDEC, double *aX, double *aY)
Convert from celestial spherical coordinates to pixel coordinates.
Definition: muse_wcs.h:164
cpl_error_code muse_datacube_save_recimages(const char *aFilename, muse_imagelist *aImages, cpl_array *aNames)
Save reconstructed images of a cube in extra extensions.
void muse_euro3dcube_delete(muse_euro3dcube *aEuro3D)
Deallocate memory associated to a muse_euro3dcube object.
cpl_imagelist * data
the cube containing the actual data values
Definition: muse_datacube.h:75
Structure definition of a Euro3D datacube.
Definition: muse_datacube.h:96
cpl_imagelist * dq
the optional cube containing the bad pixel status
Definition: muse_datacube.h:80
cpl_table * dtable
the table containing the actual Euro3D data
cpl_propertylist * hdata
the data FITS header
cpl_table * gtable
the table containing the Euro3D groups
cpl_propertylist * header
the FITS header
Definition: muse_datacube.h:56
muse_image * muse_euro3dcube_collapse(muse_euro3dcube *aEuro3D, cpl_table *aFilter)
Integrate a Euro3D datacube along the wavelength direction.
muse_imagelist * muse_imagelist_new(void)
Create a new (empty) MUSE image list.
cpl_boolean iscelsph
Definition: muse_wcs.h:100
muse_image * muse_image_new(void)
Allocate memory for a new muse_image object.
Definition: muse_image.c:66
double crval2
Definition: muse_wcs.h:95
cpl_error_code muse_imagelist_set(muse_imagelist *aList, muse_image *aImage, unsigned int aIdx)
Set the muse_image of given list index.
muse_imagelist * recimages
the reconstructed image data
Definition: muse_datacube.h:63
cpl_error_code muse_utils_copy_modified_header(cpl_propertylist *aH1, cpl_propertylist *aH2, const char *aKeyword, const char *aString)
Copy a modified header keyword from one header to another.
Definition: muse_utils.c:2489
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:897
cpl_imagelist * stat
the cube containing the data variance
Definition: muse_datacube.h:85