MUSE Pipeline Reference Manual  0.18.5
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  *
5  * This file is part of the MUSE Instrument Pipeline
6  * Copyright (C) 2005-2011 European Southern Observatory
7  *
8  * This program is free software; you can redistribute it and/or modify
9  * it under the terms of the GNU General Public License as published by
10  * the Free Software Foundation; either version 2 of the License, or
11  * (at your option) any later version.
12  *
13  * This program is distributed in the hope that it will be useful,
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  * GNU General Public License for more details.
17  *
18  * You should have received a copy of the GNU General Public License
19  * along with this program; if not, write to the Free Software
20  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
21  */
22 
23 #ifdef HAVE_CONFIG_H
24 #include <config.h>
25 #endif
26 
27 /*----------------------------------------------------------------------------*
28  * Includes *
29  *----------------------------------------------------------------------------*/
30 #include <cpl.h>
31 #include <math.h>
32 #include <string.h>
33 
34 #include "muse_datacube.h"
35 
36 #include "muse_flux.h"
37 #include "muse_pfits.h"
38 #include "muse_quality.h"
39 #include "muse_utils.h"
40 #include "muse_wcs.h"
41 
42 /*----------------------------------------------------------------------------*/
49 /*----------------------------------------------------------------------------*/
50 
53 /*---------------------------------------------------------------------------*/
70 /*---------------------------------------------------------------------------*/
71 static double *
72 muse_datacube_collapse_filter_buffer(cpl_table *aFilter,
73  double aCRVAL, double aCRPIX, double aCDELT,
74  cpl_boolean aIsLog, int *aL1, int *aL2)
75 {
76  if (!aFilter || !aL1 || !aL2) {
77  return NULL;
78  }
79  double *fdata = (double *)cpl_calloc(*aL2, sizeof(double));
80 
81  int l;
82  for (l = 0; l < *aL2; l++) {
83  double lambda = (l + 1. - aCRPIX) * aCDELT + aCRVAL;
84  if (aIsLog) {
85  lambda = aCRVAL * exp((l + 1. - aCRPIX) * aCDELT / aCRVAL);
86  }
87  fdata[l] = muse_flux_response_interpolate(aFilter, lambda, NULL,
89  } /* for l (z / wavelengths) */
90 
91  /* try to find filter edges */
92  l = 0;
93  while (l < *aL2 && fabs(fdata[l]) < DBL_EPSILON) {
94  *aL1 = l++;
95  }
96  l = *aL2 - 1;
97  while (l > *aL1 && fabs(fdata[l]) < DBL_EPSILON) {
98  *aL2 = l--;
99  }
100  cpl_msg_debug(__func__, "filter wavelength range %.1f..%.1fA (planes %d..%d)",
101  (*aL1 + 1. - aCRPIX) * aCDELT + aCRVAL,
102  (*aL2 + 1. - aCRPIX) * aCDELT + aCRVAL,
103  *aL1, *aL2);
104 
105  return fdata;
106 } /* muse_datacube_collapse_filter_buffer() */
107 
108 /*---------------------------------------------------------------------------*/
129 /*---------------------------------------------------------------------------*/
130 muse_image *
131 muse_euro3dcube_collapse(muse_euro3dcube *aEuro3D, cpl_table *aFilter)
132 {
133  cpl_ensure(aEuro3D && aEuro3D->dtable && aEuro3D->hdata, CPL_ERROR_NULL_INPUT,
134  NULL);
135 
136  /* guess dimensions of the output image from the Euro3D data */
137  muse_wcs *wcs = muse_wcs_new(aEuro3D->header);
138  /* similar test as in muse_pixtable_wcs_check() */
139  wcs->iscelsph = CPL_FALSE;
140  const char *unitx = cpl_table_get_column_unit(aEuro3D->dtable, "XPOS"),
141  *unity = cpl_table_get_column_unit(aEuro3D->dtable, "YPOS");
142  cpl_ensure(unitx && unity, CPL_ERROR_DATA_NOT_FOUND, NULL);
143  /* should be equal in the first three characters */
144  if (!strncmp(unitx, unity, 4) && !strncmp(unitx, "deg", 4)) {
145  wcs->iscelsph = CPL_TRUE;
146  }
147  /* extreme coordinates in pixel or celestial coordinates */
148  double xmin = cpl_table_get_column_min(aEuro3D->dtable, "XPOS"),
149  xmax = cpl_table_get_column_max(aEuro3D->dtable, "XPOS"),
150  ymin = cpl_table_get_column_min(aEuro3D->dtable, "YPOS"),
151  ymax = cpl_table_get_column_max(aEuro3D->dtable, "YPOS");
152  double x1, x2, y1, y2; /* pixel coordinates */
153  if (wcs->iscelsph) {
154  wcs->crval1 /= CPL_MATH_DEG_RAD; /* convert to radians */
155  wcs->crval2 /= CPL_MATH_DEG_RAD;
156  muse_wcs_pixel_from_celestial_fast(wcs, xmin / CPL_MATH_DEG_RAD,
157  ymin / CPL_MATH_DEG_RAD, &x1, &y1);
158  muse_wcs_pixel_from_celestial_fast(wcs, xmax / CPL_MATH_DEG_RAD,
159  ymax / CPL_MATH_DEG_RAD, &x2, &y2);
160  } else {
161  x1 = xmin;
162  x2 = xmax;
163  y1 = ymin;
164  y2 = ymax;
165  }
166  int zmin = cpl_table_get_column_min(aEuro3D->dtable, "SPEC_STA"),
167  zmax = cpl_table_get_column_max(aEuro3D->dtable, "SPEC_STA"),
168  nx = lround(fabs(x2 - x1)) + 1,
169  ny = lround(fabs(y2 - y1)) + 1,
170  nz = zmax - zmin + 1;
171 
172  /* count how many valid elements really are in the most shifted spectrum */
173  cpl_size zmaxpos = -1;
174  cpl_table_get_column_maxpos(aEuro3D->dtable, "SPEC_STA", &zmaxpos);
175  const cpl_array *amax = cpl_table_get_array(aEuro3D->dtable, "DATA_SPE",
176  zmaxpos);
177  int l, nsize = cpl_array_get_size(amax);
178  for (l = nsize - 1; l > 0; l--) {
179  /* find last valid element */
180  if (cpl_array_is_valid(amax, l) == 1) {
181  break;
182  }
183  }
184  nz += l + 1; /* add the number of valid elements in spectral direction */
185  int nspec = cpl_table_get_nrow(aEuro3D->dtable);
186  cpl_msg_debug(__func__, "Euro3D dimensions: %dx%dx%d (z = %d - %d, valid %d),"
187  " %d spectra", nx, ny, nz, zmax, zmin, l + 1, nspec);
188 
189  /* resample the filter curve on the wavelength grid */
190  double crvals = cpl_propertylist_get_double(aEuro3D->hdata, "CRVALS"),
191  cdelts = cpl_propertylist_get_double(aEuro3D->hdata, "CDELTS");
192  const char *ctypes = cpl_propertylist_get_string(aEuro3D->hdata, "CTYPES");
193  cpl_boolean loglambda = ctypes && !strncmp(ctypes, "AWAV-LOG", 9);
194  cpl_msg_debug(__func__, "spectral WCS: %f / %f %s", crvals, cdelts,
195  loglambda ? "log" : "linear");
196  int l1 = 0, l2 = nz; /* loop boundaries for filter operation */
197  double *fdata = muse_datacube_collapse_filter_buffer(aFilter, crvals, zmin,
198  cdelts, loglambda, &l1, &l2);
199 
200  /* create output muse_image, but only with two image extensions */
201  muse_image *image = muse_image_new();
202  image->header = cpl_propertylist_duplicate(aEuro3D->header);
203  /* the axis3 WCS was already erased in muse_resampling_euro3d() */
204  image->data = cpl_image_new(nx, ny, CPL_TYPE_FLOAT);
205  float *outdata = cpl_image_get_data_float(image->data);
206  image->dq = cpl_image_new(nx, ny, CPL_TYPE_INT);
207  /* pre-fill with Euro3D status for missing data */
208  cpl_image_add_scalar(image->dq, EURO3D_MISSDATA);
209  int *outdq = cpl_image_get_data_int(image->dq);
210  /* image->stat remains NULL, variance information is lost by collapsing */
211 
212  /* check if we need to use the variance for weighting */
213  cpl_boolean usevariance = CPL_FALSE;
214  if (getenv("MUSE_COLLAPSE_USE_VARIANCE") &&
215  atoi(getenv("MUSE_COLLAPSE_USE_VARIANCE")) > 0) {
216  usevariance = CPL_TRUE;
217  }
218 
219  /* loop through all spaxels and over all wavelengths */
220  int k, noutside = 0;
221  #pragma omp parallel for default(none) private(l) /* as req. by Ralf */ \
222  shared(aEuro3D, fdata, l1, l2, noutside, nspec, nx, ny, outdata, \
223  outdq, usevariance, wcs)
224  for (k = 0; k < nspec; k++) {
225  int err;
226  double xpos = cpl_table_get(aEuro3D->dtable, "XPOS", k, &err),
227  ypos = cpl_table_get(aEuro3D->dtable, "YPOS", k, &err);
228  if (err) {
229  cpl_msg_warning(__func__, "spectrum %d in Euro3D table does not have "
230  "position information!", k + 1);
231  continue;
232  }
233  /* coordinate in the output image (indices starting at 0) */
234  double xpx, ypx;
235  if (wcs->iscelsph) {
236  muse_wcs_pixel_from_celestial_fast(wcs, xpos * CPL_MATH_RAD_DEG,
237  ypos * CPL_MATH_RAD_DEG, &xpx, &ypx);
238  } else {
239  muse_wcs_pixel_from_projplane_fast(wcs, xpos, ypos, &xpx, &ypx);
240  }
241  int i = lround(xpx) - 1,
242  j = lround(ypx) - 1;
243  if (i >= nx || j >= ny) {
244  /* should not happen, but skip this spectrum to avoid segfault */
245  noutside++;
246  continue;
247  }
248 
249  int nstart = cpl_table_get_int(aEuro3D->dtable, "SPEC_STA", k, &err);
250  const cpl_array *adata = cpl_table_get_array(aEuro3D->dtable, "DATA_SPE", k),
251  *adq = cpl_table_get_array(aEuro3D->dtable, "QUAL_SPE", k),
252  *astat = NULL;
253  if (usevariance) {
254  astat = cpl_table_get_array(aEuro3D->dtable, "STAT_SPE", k);
255  }
256 
257  double sumdata = 0., sumweight = 0.;
258  for (l = l1; l < l2; l++) {
259  /* array index for this wavelength */
260  int idx = l - nstart + 1;
261 
262  /* filter weight with fallback for operation without filter */
263  double fweight = fdata ? fdata[l] : 1.;
264  cpl_errorstate prestate = cpl_errorstate_get();
265  int dq = cpl_array_get_int(adq, idx, &err);
266  if (dq || err) { /* if pixel marked bad or inaccessible */
267  cpl_errorstate_set(prestate);
268  continue; /* ignore bad pixels of any kind */
269  }
270  /* make the variance information optional */
271  double variance = 1.;
272  if (usevariance) {
273  variance = astat ? cpl_array_get(astat, idx, &err) : 1.;
274  /* use ^2 as for Euro3D we store errors not variances directly */
275  variance *= variance;
276  if (err > 0) {
277  variance = DBL_MAX;
278  }
279  if (err || !isnormal(variance)) {
280  continue;
281  }
282  }
283  double data = cpl_array_get(adata, idx, &err);
284  /* weight by inverse variance */
285  sumdata += data * fweight / variance;
286  sumweight += fweight / variance;
287  } /* for l (z / wavelengths) */
288  if (isnormal(sumweight)) {
289  outdata[i + j*nx] = sumdata / sumweight;
290  outdq[i + j*nx] = EURO3D_GOODPIXEL;
291  } else { /* leave bad/missing pixels zero, but mark them bad */
292  outdq[i + j*nx] = EURO3D_MISSDATA;
293  }
294  } /* for k (all spaxels) */
295  cpl_free(wcs);
296  cpl_free(fdata); /* won't crash on NULL */
297  if (noutside > 0) {
298  cpl_msg_warning(__func__, "Skipped %d spaxels, due to their location "
299  "outside the recostructed image!", noutside);
300  }
301 
302  return image;
303 } /* muse_euro3dcube_collapse() */
304 
305 /*---------------------------------------------------------------------------*/
324 /*---------------------------------------------------------------------------*/
325 muse_image *
326 muse_datacube_collapse(muse_datacube *aCube, cpl_table *aFilter)
327 {
328  cpl_ensure(aCube && aCube->data && aCube->header, CPL_ERROR_NULL_INPUT, NULL);
329 
330  /* find in/output dimensions */
331  int nx = cpl_image_get_size_x(cpl_imagelist_get(aCube->data, 0)),
332  ny = cpl_image_get_size_y(cpl_imagelist_get(aCube->data, 0)),
333  nz = cpl_imagelist_get_size(aCube->data);
334  /* resample the filter curve on the wavelength grid */
335  double crpix3 = cpl_propertylist_get_double(aCube->header, "CRPIX3"),
336  crval3 = cpl_propertylist_get_double(aCube->header, "CRVAL3"),
337  cd33 = cpl_propertylist_get_double(aCube->header, "CD3_3");
338  const char *ctype3 = cpl_propertylist_get_string(aCube->header, "CTYPE3");
339  cpl_boolean loglambda = ctype3 && !strncmp(ctype3, "AWAV-LOG", 9);
340  int l1 = 0, l2 = nz; /* loop boundaries for filter operation */
341  double *fdata = muse_datacube_collapse_filter_buffer(aFilter, crval3, crpix3,
342  cd33, loglambda, &l1, &l2);
343 
344  /* create output muse_image, but only with two image extensions */
345  muse_image *image = muse_image_new();
346  image->header = cpl_propertylist_duplicate(aCube->header);
347  cpl_propertylist_erase_regexp(image->header, "^C...*3$|^CD3_.$", 0);
348  image->data = cpl_image_new(nx, ny, CPL_TYPE_FLOAT);
349  float *outdata = cpl_image_get_data_float(image->data);
350  image->dq = cpl_image_new(nx, ny, CPL_TYPE_INT);
351  int *outdq = cpl_image_get_data_int(image->dq);
352  /* image->stat remains NULL, variance information is lost by collapsing */
353 
354  /* check if we need to use the variance for weighting */
355  cpl_boolean usevariance = CPL_FALSE;
356  if (getenv("MUSE_COLLAPSE_USE_VARIANCE") &&
357  atoi(getenv("MUSE_COLLAPSE_USE_VARIANCE")) > 0) {
358  usevariance = CPL_TRUE;
359  }
360 
361  /* loop through all pixels in all wavelength planes */
362  int i;
363  #pragma omp parallel for default(none) /* as req. by Ralf */ \
364  shared(aCube, nx, ny, l1, l2, fdata, outdata, outdq, usevariance)
365  for (i = 0; i < nx; i++) {
366  int j;
367  for (j = 0; j < ny; j++) {
368  double sumdata = 0., sumweight = 0.;
369  int l;
370  for (l = l1; l < l2; l++) {
371  /* filter weight with fallback for operation without filter */
372  double fweight = fdata ? fdata[l] : 1.;
373  float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->data, l)),
374  *pstat = NULL;
375  if (aCube->dq) {
376  int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->dq, l));
377  if (pdq && pdq[i + j*nx]) {
378  continue; /* ignore bad pixels of any kind */
379  }
380  }
381  /* make the variance information optional */
382  double variance = 1.;
383  if (usevariance) {
384  pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->stat, l));
385  variance = pstat ? pstat[i + j*nx] : 1.;
386  if (!isnormal(variance)) {
387  continue;
388  }
389  }
390  /* weight by inverse variance */
391  sumdata += pdata[i + j*nx] * fweight / variance;
392  sumweight += fweight / variance;
393  } /* for l (z / wavelengths) */
394  if (isnormal(sumweight)) {
395  outdata[i + j*nx] = sumdata / sumweight;
396  outdq[i + j*nx] = EURO3D_GOODPIXEL;
397  } else { /* leave bad/missing pixels zero, but mark them bad */
398  outdq[i + j*nx] = EURO3D_MISSDATA;
399  }
400  } /* for j (y pixels) */
401  } /* for i (x pixels) */
402 
403  cpl_free(fdata); /* won't crash on NULL */
404 
405  return image;
406 } /* muse_datacube_collapse() */
407 
408 /*---------------------------------------------------------------------------*/
420 /*---------------------------------------------------------------------------*/
421 static cpl_error_code
422 muse_datacube_save_recimages(const char *aFilename, muse_imagelist *aImages,
423  cpl_array *aNames)
424 {
425  cpl_ensure_code(aFilename, CPL_ERROR_NULL_INPUT);
426  cpl_error_code rc = CPL_ERROR_NONE;
427  if (!aImages || !aNames) {
428  return rc; /* this is a valid case, return without error */
429  }
430  unsigned int i, nimages = muse_imagelist_get_size(aImages);
431  for (i = 0; i < nimages; i++) {
432  muse_image *image = muse_imagelist_get(aImages, i);
433 
434  /* create small header with EXTNAME=filtername and WCS for images */
435  cpl_propertylist *header = cpl_propertylist_new();
436  char dataext[KEYWORD_LENGTH], obj[KEYWORD_LENGTH],
437  *dqext = NULL, *statext = NULL;
438  snprintf(dataext, KEYWORD_LENGTH, "%s", cpl_array_get_string(aNames, i));
439  if (image->dq) {
440  dqext = cpl_sprintf("%s_%s", cpl_array_get_string(aNames, i), EXTNAME_DQ);
441  }
442  if (image->stat) {
443  statext = cpl_sprintf("%s_%s", cpl_array_get_string(aNames, i), EXTNAME_STAT);
444  }
445  snprintf(obj, KEYWORD_LENGTH, "%s", cpl_array_get_string(aNames, i));
446  cpl_propertylist_append_string(header, "EXTNAME", dataext);
447  cpl_propertylist_set_comment(header, "EXTNAME", "reconstructed image (data values)");
448  muse_utils_copy_modified_header(image->header, header, "OBJECT", obj);
449  cpl_propertylist_update_string(header, MUSE_HDR_FILTER,
450  cpl_array_get_string(aNames, i));
451  cpl_propertylist_copy_property_regexp(header, image->header,
452  MUSE_WCS_KEYS, 0);
453  muse_utils_set_hduclass(header, "DATA", dataext, dqext, statext);
454  rc = cpl_image_save(image->data, aFilename, CPL_TYPE_UNSPECIFIED, header,
455  CPL_IO_EXTEND);
456 
457  if (image->dq) {
458  cpl_propertylist_update_string(header, "EXTNAME", dqext);
459  cpl_propertylist_set_comment(header, "EXTNAME", "reconstructed image (bad pixel status values)");
460  snprintf(obj, KEYWORD_LENGTH, "%s, %s", cpl_array_get_string(aNames, i),
461  EXTNAME_DQ);
462  muse_utils_copy_modified_header(image->header, header, "OBJECT", obj);
463  muse_utils_set_hduclass(header, "QUALITY", dataext, dqext, statext);
464  rc = cpl_image_save(image->dq, aFilename, CPL_TYPE_UNSPECIFIED, header,
465  CPL_IO_EXTEND);
466  }
467  if (image->stat) {
468  cpl_propertylist_update_string(header, "EXTNAME", statext);
469  cpl_propertylist_set_comment(header, "EXTNAME", "reconstructed image (variance)");
470  snprintf(obj, KEYWORD_LENGTH, "%s, %s", cpl_array_get_string(aNames, i),
471  EXTNAME_STAT);
472  muse_utils_copy_modified_header(image->header, header, "OBJECT", obj);
473  muse_utils_set_hduclass(header, "ERROR", dataext, dqext, statext);
474  rc = cpl_image_save(image->stat, aFilename, CPL_TYPE_UNSPECIFIED, header,
475  CPL_IO_EXTEND);
476  }
477 
478  cpl_propertylist_delete(header);
479  cpl_free(dqext);
480  cpl_free(statext);
481  } /* for i (all reconstructed images) */
482 
483  return rc;
484 } /* muse_datacube_save_recimages */
485 
486 /*---------------------------------------------------------------------------*/
501 /*---------------------------------------------------------------------------*/
502 cpl_error_code
503 muse_euro3dcube_save(muse_euro3dcube *aEuro3D, const char *aFilename)
504 {
505  cpl_ensure_code(aEuro3D && aFilename, CPL_ERROR_NULL_INPUT);
506 
507  /* save primary header and the data table in the first extension */
508  cpl_error_code rc = cpl_table_save(aEuro3D->dtable, aEuro3D->header,
509  aEuro3D->hdata, aFilename, CPL_IO_DEFAULT);
510  if (rc != CPL_ERROR_NONE) {
511  cpl_msg_warning(__func__, "failed to save data part of the Euro3D table: "
512  "%s", cpl_error_get_message());
513  }
514  /* save group header and data in the second extension */
515  rc = cpl_table_save(aEuro3D->gtable, NULL, aEuro3D->hgroup, aFilename,
516  CPL_IO_EXTEND);
517  if (rc != CPL_ERROR_NONE) {
518  cpl_msg_warning(__func__, "failed to save group part of the Euro3D table: "
519  "%s", cpl_error_get_message());
520  }
521 
522  rc = muse_datacube_save_recimages(aFilename, aEuro3D->recimages,
523  aEuro3D->recnames);
524  return rc;
525 } /* muse_euro3dcube_save() */
526 
527 /*---------------------------------------------------------------------------*/
537 /*---------------------------------------------------------------------------*/
538 void
540 {
541  /* if the euro3d object doesn't exists at all, we don't need to do anything */
542  if (!aEuro3D) {
543  return;
544  }
545 
546  /* checks for the existence of the sub-images *
547  * are done in the CPL functions */
548  cpl_table_delete(aEuro3D->dtable);
549  cpl_table_delete(aEuro3D->gtable);
550  cpl_propertylist_delete(aEuro3D->header);
551  cpl_propertylist_delete(aEuro3D->hdata);
552  cpl_propertylist_delete(aEuro3D->hgroup);
554  cpl_array_delete(aEuro3D->recnames);
555  cpl_free(aEuro3D);
556 } /* muse_euro3dcube_delete() */
557 
558 /*---------------------------------------------------------------------------*/
568 /*---------------------------------------------------------------------------*/
569 static cpl_error_code
570 muse_datacube_convert_dq_recimages(muse_datacube *aCube)
571 {
572  cpl_ensure_code(aCube, CPL_ERROR_NULL_INPUT);
573  if (!aCube->recimages) {
574  return CPL_ERROR_NONE; /* valid case, return without error */
575  }
576  unsigned int k, nimages = muse_imagelist_get_size(aCube->recimages);
577  for (k = 0; k < nimages; k++) {
578  muse_image *image = muse_imagelist_get(aCube->recimages, k);
579  /* XXX the following could be made part of muse_image in the future */
580  if (!image->dq) { /* no point trying to convert... */
581  continue;
582  }
583  int *pdq = cpl_image_get_data_int(image->dq);
584  float *pdata = cpl_image_get_data_float(image->data),
585  *pstat = NULL;
586  if (image->stat) {
587  pstat = cpl_image_get_data_float(image->stat);
588  }
589  int i, nx = cpl_image_get_size_x(image->data),
590  ny = cpl_image_get_size_y(image->data);
591  for (i = 0; i < nx; i++) {
592  int j;
593  for (j = 0; j < ny; j++) {
594  if (pdq[i + j*nx] == EURO3D_GOODPIXEL) {
595  continue; /* nothing to do for good pixels */
596  }
597  /* set bad pixels of any type to not-a-number in both extensions */
598  pdata[i + j*nx] = NAN; /* supposed to be quiet NaN */
599  if (pstat) {
600  pstat[i + j*nx] = NAN;
601  }
602  } /* for j (y direction) */
603  } /* for i (x direction) */
604 
605  /* deallocate DQ and set it to NULL to be able to check for it */
606  cpl_image_delete(image->dq);
607  image->dq = NULL;
608  } /* for k (all reconstructed images) */
609 
610  return CPL_ERROR_NONE;
611 } /* muse_datacube_convert_dq_recimages() */
612 
613 /*---------------------------------------------------------------------------*/
626 /*---------------------------------------------------------------------------*/
627 cpl_error_code
629 {
630  cpl_ensure_code(aCube && aCube->data && aCube->stat && aCube->dq,
631  CPL_ERROR_NULL_INPUT);
632 
633  /* loop through all wavelength planes and then all pixels per plane */
634  int l, nx = cpl_image_get_size_x(cpl_imagelist_get(aCube->data, 0)),
635  ny = cpl_image_get_size_y(cpl_imagelist_get(aCube->data, 0)),
636  nz = cpl_imagelist_get_size(aCube->data);
637  #pragma omp parallel for default(none) /* as req. by Ralf */ \
638  shared(aCube, nx, ny, nz)
639  for (l = 0; l < nz; l++) {
640  int i;
641  for (i = 0; i < nx; i++) {
642  int j;
643  for (j = 0; j < ny; j++) {
644  float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->data, l)),
645  *pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->stat, l));
646  int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->dq, l));
647  if (pdq[i + j*nx] == EURO3D_GOODPIXEL) {
648  continue; /* nothing to do for good pixels */
649  }
650  /* set bad pixels of any type to not-a-number in both extensions */
651  pdata[i + j*nx] = NAN; /* supposed to be quiet NaN */
652  pstat[i + j*nx] = NAN;
653  } /* for j (y direction) */
654  } /* for i (x direction) */
655  } /* for l (wavelength planes) */
656 
657  /* deallocate DQ and set it to NULL to be able to check for it */
658  cpl_imagelist_delete(aCube->dq);
659  aCube->dq = NULL;
660 
661  /* do the same for the reconstructed images, if there are any */
662  muse_datacube_convert_dq_recimages(aCube);
663 
664  return CPL_ERROR_NONE;
665 } /* muse_datacube_convert_dq() */
666 
667 /*---------------------------------------------------------------------------*/
696 /*---------------------------------------------------------------------------*/
697 cpl_error_code
698 muse_datacube_save(muse_datacube *aCube, const char *aFilename)
699 {
700  cpl_ensure_code(aCube && aCube->header && aFilename, CPL_ERROR_NULL_INPUT);
701 
702  /* save headers into primary area */
703  cpl_propertylist *header = cpl_propertylist_new();
704  /* copy all headers, except the main WCS keys */
705  cpl_propertylist_copy_property_regexp(header, aCube->header,
706  MUSE_WCS_KEYS, 1);
707  cpl_error_code rc = cpl_propertylist_save(header, aFilename, CPL_IO_CREATE);
708  cpl_propertylist_delete(header);
709 
710  /* save data */
711  header = cpl_propertylist_new();
712  cpl_propertylist_append_string(header, "EXTNAME", EXTNAME_DATA);
713  cpl_propertylist_set_comment(header, "EXTNAME", EXTNAME_DATA_COMMENT);
714  muse_utils_copy_modified_header(aCube->header, header, "OBJECT",
715  EXTNAME_DATA);
716  cpl_propertylist_copy_property_regexp(header, aCube->header,
717  MUSE_WCS_KEYS"|^BUNIT", 0);
718  muse_utils_set_hduclass(header, "DATA", "DATA",
719  aCube->dq ? "DQ" : NULL, aCube->stat ? "STAT" : NULL);
720  rc = cpl_imagelist_save(aCube->data, aFilename, CPL_TYPE_FLOAT, header,
721  CPL_IO_EXTEND);
722  cpl_propertylist_delete(header);
723  /* save bad pixels, if available */
724  if (rc == CPL_ERROR_NONE && aCube->dq) {
725  header = cpl_propertylist_new();
726  cpl_propertylist_append_string(header, "EXTNAME", EXTNAME_DQ);
727  cpl_propertylist_set_comment(header, "EXTNAME", EXTNAME_DQ_COMMENT);
728  muse_utils_copy_modified_header(aCube->header, header, "OBJECT",
729  EXTNAME_DQ);
730  cpl_propertylist_copy_property_regexp(header, aCube->header,
731  MUSE_WCS_KEYS, 0);
732  muse_utils_set_hduclass(header, "QUALITY", "DATA", "DQ",
733  aCube->stat ? "STAT" : NULL);
734  rc = cpl_imagelist_save(aCube->dq, aFilename, CPL_TYPE_INT, header,
735  CPL_IO_EXTEND);
736  cpl_propertylist_delete(header);
737  }
738  /* save variance, if available */
739  if (rc == CPL_ERROR_NONE && aCube->stat) {
740  header = cpl_propertylist_new();
741  cpl_propertylist_append_string(header, "EXTNAME", EXTNAME_STAT);
742  cpl_propertylist_set_comment(header, "EXTNAME", EXTNAME_STAT_COMMENT);
743  /* add the correct BUNIT, depending on the data units */
744  if (!strncmp(cpl_propertylist_get_string(aCube->header, "BUNIT"),
745  kMuseFluxUnitString, strlen(kMuseFluxUnitString) + 1)) {
746  /* flux calibrated data */
747  cpl_propertylist_append_string(header, "BUNIT", kMuseFluxStatString);
748  }
749  muse_utils_copy_modified_header(aCube->header, header, "OBJECT",
750  EXTNAME_STAT);
751  cpl_propertylist_copy_property_regexp(header, aCube->header,
752  MUSE_WCS_KEYS, 0);
753  muse_utils_set_hduclass(header, "ERROR", "DATA", aCube->dq ? "DQ" : NULL,
754  "STAT");
755  rc = cpl_imagelist_save(aCube->stat, aFilename, CPL_TYPE_FLOAT, header,
756  CPL_IO_EXTEND);
757  cpl_propertylist_delete(header);
758  }
759 
760  rc = muse_datacube_save_recimages(aFilename, aCube->recimages, aCube->recnames);
761  return rc;
762 } /* muse_datacube_save() */
763 
764 /*---------------------------------------------------------------------------*/
774 /*---------------------------------------------------------------------------*/
775 void
777 {
778  /* if the cube does not exists at all, we don't need to do anything */
779  if (!aCube) {
780  return;
781  }
782 
783  /* checks for the existence of the sub-images *
784  * are done in the CPL functions */
785  cpl_imagelist_delete(aCube->data);
786  aCube->data = NULL;
787  cpl_imagelist_delete(aCube->dq);
788  aCube->dq = NULL;
789  cpl_imagelist_delete(aCube->stat);
790  aCube->stat = NULL;
791 
792  /* delete the FITS header, too */
793  cpl_propertylist_delete(aCube->header);
794  aCube->header = NULL;
795 
796  /* remove reconstructed image data, if present */
798  cpl_array_delete(aCube->recnames);
799  cpl_free(aCube);
800 } /* muse_datacube_delete() */
801 
Structure definition of a MUSE datacube.
Definition: muse_datacube.h:48
Structure definition for a collection of muse_images.
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:47
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:65
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:41
cpl_array * recnames
the reconstructed image filter names
Definition: muse_datacube.h:71
cpl_propertylist * header
the FITS header
Definition: muse_image.h:73
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:2181
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:57
double muse_flux_response_interpolate(const cpl_table *aResponse, double aLambda, double *aError, muse_flux_interpolation_type aType)
Compute linearly interpolated response of some kind at given wavelength.
Definition: muse_flux.c:322
muse_wcs * muse_wcs_new(cpl_propertylist *aHeader)
Create a new WCS structure from a given FITS header.
Definition: muse_wcs.c:1536
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
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:76
Structure definition of a Euro3D datacube.
Definition: muse_datacube.h:97
cpl_imagelist * dq
the optional cube containing the bad pixel status
Definition: muse_datacube.h:81
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:57
muse_image * muse_euro3dcube_collapse(muse_euro3dcube *aEuro3D, cpl_table *aFilter)
Integrate a Euro3D datacube along the wavelength direction.
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:65
double crval2
Definition: muse_wcs.h:95
muse_imagelist * recimages
the reconstructed image data
Definition: muse_datacube.h:64
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:2140
cpl_imagelist * stat
the cube containing the data variance
Definition: muse_datacube.h:86