MUSE Pipeline Reference Manual  2.1.1
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-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 <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 /*---------------------------------------------------------------------------*/
70 /*---------------------------------------------------------------------------*/
71 static double *
72 muse_datacube_collapse_filter_buffer(const muse_table *aFilter,
73  double aCRVAL, double aCRPIX, double aCDELT,
74  cpl_boolean aIsLog, int *aL1, int *aL2,
75  double *aFraction)
76 {
77  if (!aFilter || !aL1 || !aL2) {
78  return NULL;
79  }
80  if (!aFilter->table) {
81  return NULL;
82  }
83  double *fdata = (double *)cpl_calloc(*aL2, sizeof(double));
84 
85  int l;
86  for (l = 0; l < *aL2; l++) {
87  double lambda = (l + 1. - aCRPIX) * aCDELT + aCRVAL;
88  if (aIsLog) {
89  lambda = aCRVAL * exp((l + 1. - aCRPIX) * aCDELT / aCRVAL);
90  }
91  fdata[l] = muse_flux_response_interpolate(aFilter->table, lambda, NULL,
93  } /* for l (z / wavelengths) */
94 
95  /* try to find filter edges */
96  l = 0;
97  while (l < *aL2 && fabs(fdata[l]) < DBL_EPSILON) {
98  *aL1 = l++;
99  }
100  l = *aL2 - 1;
101  while (l > *aL1 && fabs(fdata[l]) < DBL_EPSILON) {
102  *aL2 = l--;
103  }
104  double lbda1 = (*aL1 + 1. - aCRPIX) * aCDELT + aCRVAL,
105  lbda2 = (*aL2 + 1. - aCRPIX) * aCDELT + aCRVAL;
106  if (aIsLog) {
107  lbda1 = aCRVAL * exp((*aL1 + 1. - aCRPIX) * aCDELT / aCRVAL);
108  lbda2 = aCRVAL * exp((*aL2 + 1. - aCRPIX) * aCDELT / aCRVAL);
109  }
110  double fraction = muse_utils_filter_fraction(aFilter, lbda1, lbda2);
111  cpl_msg_debug(__func__, "Filter wavelength range %.1f..%.1fA (cube planes "
112  "%d..%d), %.2f%% of filter area covered by data.", lbda1, lbda2,
113  *aL1, *aL2, fraction * 100.);
114  if (aFraction) {
115  *aFraction = fraction;
116  }
117 
118  return fdata;
119 } /* muse_datacube_collapse_filter_buffer() */
120 
121 /*---------------------------------------------------------------------------*/
143 /*---------------------------------------------------------------------------*/
144 muse_image *
146 {
147  cpl_ensure(aEuro3D && aEuro3D->dtable && aEuro3D->hdata, CPL_ERROR_NULL_INPUT,
148  NULL);
149 
150  /* guess dimensions of the output image from the Euro3D data */
151  muse_wcs *wcs = muse_wcs_new(aEuro3D->header);
152  /* similar test as in muse_pixtable_wcs_check() */
153  wcs->iscelsph = CPL_FALSE;
154  const char *unitx = cpl_table_get_column_unit(aEuro3D->dtable, "XPOS"),
155  *unity = cpl_table_get_column_unit(aEuro3D->dtable, "YPOS");
156  cpl_ensure(unitx && unity, CPL_ERROR_DATA_NOT_FOUND, NULL);
157  /* should be equal in the first three characters */
158  if (!strncmp(unitx, unity, 4) && !strncmp(unitx, "deg", 4)) {
159  wcs->iscelsph = CPL_TRUE;
160  }
161  /* extreme coordinates in pixel or celestial coordinates */
162  double xmin = cpl_table_get_column_min(aEuro3D->dtable, "XPOS"),
163  xmax = cpl_table_get_column_max(aEuro3D->dtable, "XPOS"),
164  ymin = cpl_table_get_column_min(aEuro3D->dtable, "YPOS"),
165  ymax = cpl_table_get_column_max(aEuro3D->dtable, "YPOS");
166  double x1, x2, y1, y2; /* pixel coordinates */
167  if (wcs->iscelsph) {
168  wcs->crval1 /= CPL_MATH_DEG_RAD; /* convert to radians */
169  wcs->crval2 /= CPL_MATH_DEG_RAD;
170  muse_wcs_pixel_from_celestial_fast(wcs, xmin / CPL_MATH_DEG_RAD,
171  ymin / CPL_MATH_DEG_RAD, &x1, &y1);
172  muse_wcs_pixel_from_celestial_fast(wcs, xmax / CPL_MATH_DEG_RAD,
173  ymax / CPL_MATH_DEG_RAD, &x2, &y2);
174  } else {
175  x1 = xmin;
176  x2 = xmax;
177  y1 = ymin;
178  y2 = ymax;
179  }
180  int zmin = cpl_table_get_column_min(aEuro3D->dtable, "SPEC_STA"),
181  zmax = cpl_table_get_column_max(aEuro3D->dtable, "SPEC_STA"),
182  nx = lround(fabs(x2 - x1)) + 1,
183  ny = lround(fabs(y2 - y1)) + 1,
184  nz = zmax - zmin + 1;
185 
186  /* count how many valid elements really are in the most shifted spectrum */
187  cpl_size zmaxpos = -1;
188  cpl_table_get_column_maxpos(aEuro3D->dtable, "SPEC_STA", &zmaxpos);
189  const cpl_array *amax = cpl_table_get_array(aEuro3D->dtable, "DATA_SPE",
190  zmaxpos);
191  int l, nsize = cpl_array_get_size(amax);
192  for (l = nsize - 1; l > 0; l--) {
193  /* find last valid element */
194  if (cpl_array_is_valid(amax, l) == 1) {
195  break;
196  }
197  }
198  nz += l + 1; /* add the number of valid elements in spectral direction */
199  int nspec = cpl_table_get_nrow(aEuro3D->dtable);
200  cpl_msg_debug(__func__, "Euro3D dimensions: %dx%dx%d (z = %d - %d, valid %d),"
201  " %d spectra", nx, ny, nz, zmax, zmin, l + 1, nspec);
202 
203  /* resample the filter curve on the wavelength grid */
204  double crvals = cpl_propertylist_get_double(aEuro3D->hdata, "CRVALS"),
205  cdelts = cpl_propertylist_get_double(aEuro3D->hdata, "CDELTS");
206  const char *ctypes = cpl_propertylist_get_string(aEuro3D->hdata, "CTYPES");
207  cpl_boolean loglambda = ctypes && (!strncmp(ctypes, "AWAV-LOG", 9) ||
208  !strncmp(ctypes, "WAVE-LOG", 9));
209  cpl_msg_debug(__func__, "spectral WCS: %f / %f %s", crvals, cdelts,
210  loglambda ? "log" : "linear");
211  int l1 = 0, l2 = nz; /* loop boundaries for filter operation */
212  double *fdata = NULL;
213  double ffraction = 1.; /* filter fraction */
214  if (aFilter) {
215  fdata = muse_datacube_collapse_filter_buffer(aFilter, crvals, zmin,
216  cdelts, loglambda, &l1, &l2,
217  &ffraction);
218  }
219 
220  /* create output muse_image, but only with two image extensions */
221  muse_image *image = muse_image_new();
222  image->header = cpl_propertylist_duplicate(aEuro3D->header);
223  /* the axis3 WCS was already erased in muse_resampling_euro3d() */
224  if (aFilter) {
225  muse_utils_filter_copy_properties(image->header, aFilter, ffraction);
226  }
227  image->data = cpl_image_new(nx, ny, CPL_TYPE_FLOAT);
228  float *outdata = cpl_image_get_data_float(image->data);
229  image->dq = cpl_image_new(nx, ny, CPL_TYPE_INT);
230  /* pre-fill with Euro3D status for missing data */
231  cpl_image_add_scalar(image->dq, EURO3D_MISSDATA);
232  int *outdq = cpl_image_get_data_int(image->dq);
233  /* image->stat remains NULL, variance information is lost by collapsing */
234 
235  /* check if we need to use the variance for weighting */
236  cpl_boolean usevariance = CPL_FALSE;
237  if (getenv("MUSE_COLLAPSE_USE_VARIANCE") &&
238  atoi(getenv("MUSE_COLLAPSE_USE_VARIANCE")) > 0) {
239  usevariance = CPL_TRUE;
240  }
241 
242  /* loop through all spaxels and over all wavelengths */
243  int k, noutside = 0;
244  #pragma omp parallel for default(none) private(l) /* as req. by Ralf */ \
245  shared(aEuro3D, fdata, l1, l2, noutside, nspec, nx, ny, outdata, \
246  outdq, usevariance, wcs)
247  for (k = 0; k < nspec; k++) {
248  int err;
249  double xpos = cpl_table_get(aEuro3D->dtable, "XPOS", k, &err),
250  ypos = cpl_table_get(aEuro3D->dtable, "YPOS", k, &err);
251  if (err) {
252  cpl_msg_warning(__func__, "spectrum %d in Euro3D table does not have "
253  "position information!", k + 1);
254  continue;
255  }
256  /* coordinate in the output image (indices starting at 0) */
257  double xpx, ypx;
258  if (wcs->iscelsph) {
259  muse_wcs_pixel_from_celestial_fast(wcs, xpos * CPL_MATH_RAD_DEG,
260  ypos * CPL_MATH_RAD_DEG, &xpx, &ypx);
261  } else {
262  muse_wcs_pixel_from_projplane_fast(wcs, xpos, ypos, &xpx, &ypx);
263  }
264  int i = lround(xpx) - 1,
265  j = lround(ypx) - 1;
266  if (i >= nx || j >= ny) {
267  /* should not happen, but skip this spectrum to avoid segfault */
268  noutside++;
269  continue;
270  }
271 
272  int nstart = cpl_table_get_int(aEuro3D->dtable, "SPEC_STA", k, &err);
273  const cpl_array *adata = cpl_table_get_array(aEuro3D->dtable, "DATA_SPE", k),
274  *adq = cpl_table_get_array(aEuro3D->dtable, "QUAL_SPE", k),
275  *astat = NULL;
276  if (usevariance) {
277  astat = cpl_table_get_array(aEuro3D->dtable, "STAT_SPE", k);
278  }
279 
280  double sumdata = 0., sumweight = 0.;
281  for (l = l1; l < l2; l++) {
282  /* array index for this wavelength */
283  int idx = l - nstart + 1;
284 
285  /* filter weight with fallback for operation without filter */
286  double fweight = fdata ? fdata[l] : 1.;
287  cpl_errorstate prestate = cpl_errorstate_get();
288  int dq = cpl_array_get_int(adq, idx, &err);
289  if (dq || err) { /* if pixel marked bad or inaccessible */
290  cpl_errorstate_set(prestate);
291  continue; /* ignore bad pixels of any kind */
292  }
293  /* make the variance information optional */
294  double variance = 1.;
295  if (usevariance) {
296  variance = astat ? cpl_array_get(astat, idx, &err) : 1.;
297  /* use ^2 as for Euro3D we store errors not variances directly */
298  variance *= variance;
299  if (err > 0) {
300  variance = DBL_MAX;
301  }
302  if (err || !isnormal(variance)) {
303  continue;
304  }
305  }
306  double data = cpl_array_get(adata, idx, &err);
307  /* weight by inverse variance */
308  sumdata += data * fweight / variance;
309  sumweight += fweight / variance;
310  } /* for l (z / wavelengths) */
311  if (isnormal(sumweight)) {
312  outdata[i + j*nx] = sumdata / sumweight;
313  outdq[i + j*nx] = EURO3D_GOODPIXEL;
314  } else { /* leave bad/missing pixels zero, but mark them bad */
315  outdq[i + j*nx] = EURO3D_MISSDATA;
316  }
317  } /* for k (all spaxels) */
318  cpl_free(wcs);
319  cpl_free(fdata); /* won't crash on NULL */
320  if (noutside > 0) {
321  cpl_msg_warning(__func__, "Skipped %d spaxels, due to their location "
322  "outside the recostructed image!", noutside);
323  }
324 
325  return image;
326 } /* muse_euro3dcube_collapse() */
327 
328 /*---------------------------------------------------------------------------*/
349 /*---------------------------------------------------------------------------*/
350 muse_image *
352 {
353  cpl_ensure(aCube && aCube->data && aCube->header,
354  CPL_ERROR_NULL_INPUT, NULL);
355  cpl_ensure(cpl_image_get_type(cpl_imagelist_get(aCube->data, 0)) == CPL_TYPE_FLOAT,
356  CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
357  if (aCube->dq) {
358  cpl_ensure(cpl_image_get_type(cpl_imagelist_get(aCube->dq, 0)) == CPL_TYPE_INT,
359  CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
360  }
361  if (aCube->stat) {
362  cpl_ensure(cpl_image_get_type(cpl_imagelist_get(aCube->stat, 0)) == CPL_TYPE_FLOAT,
363  CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
364  }
365 
366  /* find in/output dimensions */
367  int nx = cpl_image_get_size_x(cpl_imagelist_get(aCube->data, 0)),
368  ny = cpl_image_get_size_y(cpl_imagelist_get(aCube->data, 0)),
369  nz = cpl_imagelist_get_size(aCube->data);
370  /* resample the filter curve on the wavelength grid */
371  double crpix3 = muse_pfits_get_crpix(aCube->header, 3),
372  crval3 = muse_pfits_get_crval(aCube->header, 3),
373  cd33 = muse_pfits_get_cd(aCube->header, 3, 3);
374  const char *ctype3 = muse_pfits_get_ctype(aCube->header, 3);
375  cpl_boolean loglambda = ctype3 && (!strncmp(ctype3, "AWAV-LOG", 9) ||
376  !strncmp(ctype3, "WAVE-LOG", 9));
377  int l1 = 0, l2 = nz; /* loop boundaries for filter operation */
378  double *fdata = NULL;
379  double ffraction = 1.; /* filter fraction */
380  if (aFilter) {
381  fdata = muse_datacube_collapse_filter_buffer(aFilter, crval3, crpix3,
382  cd33, loglambda, &l1, &l2,
383  &ffraction);
384  }
385 
386  /* create output muse_image, but only with two image extensions */
387  muse_image *image = muse_image_new();
388  image->header = cpl_propertylist_duplicate(aCube->header);
389  cpl_propertylist_erase_regexp(image->header, "^C...*3$|^CD3_.$", 0);
390  if (aFilter) {
391  muse_utils_filter_copy_properties(image->header, aFilter, ffraction);
392  }
393  image->data = cpl_image_new(nx, ny, CPL_TYPE_FLOAT);
394  float *outdata = cpl_image_get_data_float(image->data);
395  image->dq = cpl_image_new(nx, ny, CPL_TYPE_INT);
396  int *outdq = cpl_image_get_data_int(image->dq);
397  /* image->stat remains NULL, variance information is lost by collapsing */
398 
399  /* check if we need to use the variance for weighting */
400  cpl_boolean usevariance = CPL_FALSE;
401  if (getenv("MUSE_COLLAPSE_USE_VARIANCE") &&
402  atoi(getenv("MUSE_COLLAPSE_USE_VARIANCE")) > 0) {
403  usevariance = CPL_TRUE;
404  }
405 
406  /* loop through all pixels in all wavelength planes */
407  int i;
408  #pragma omp parallel for default(none) /* as req. by Ralf */ \
409  shared(aCube, nx, ny, l1, l2, fdata, outdata, outdq, usevariance)
410  for (i = 0; i < nx; i++) {
411  int j;
412  for (j = 0; j < ny; j++) {
413  double sumdata = 0., sumweight = 0.;
414  int l;
415  for (l = l1; l < l2; l++) {
416  /* filter weight with fallback for operation without filter */
417  double fweight = fdata ? fdata[l] : 1.;
418  float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->data, l)),
419  *pstat = NULL;
420  if (!isfinite(pdata[i + j*nx])) {
421  continue;
422  }
423  if (aCube->dq) {
424  int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->dq, l));
425  if (pdq && pdq[i + j*nx]) {
426  continue; /* ignore bad pixels of any kind */
427  }
428  }
429  /* make the variance information optional */
430  double variance = 1.;
431  if (usevariance) {
432  pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->stat, l));
433  variance = pstat ? pstat[i + j*nx] : 1.;
434  if (!isnormal(variance)) {
435  continue;
436  }
437  }
438  /* weight by inverse variance */
439  sumdata += pdata[i + j*nx] * fweight / variance;
440  sumweight += fweight / variance;
441  } /* for l (z / wavelengths) */
442  if (isnormal(sumweight)) {
443  outdata[i + j*nx] = sumdata / sumweight;
444  outdq[i + j*nx] = EURO3D_GOODPIXEL;
445  } else { /* leave bad/missing pixels zero, but mark them bad */
446  outdq[i + j*nx] = EURO3D_MISSDATA;
447  }
448  } /* for j (y pixels) */
449  } /* for i (x pixels) */
450 
451  cpl_free(fdata); /* won't crash on NULL */
452 
453  return image;
454 } /* muse_datacube_collapse() */
455 
456 /*---------------------------------------------------------------------------*/
471 /*---------------------------------------------------------------------------*/
472 cpl_error_code
473 muse_datacube_save_recimages(const char *aFilename, muse_imagelist *aImages,
474  cpl_array *aNames)
475 {
476  cpl_ensure_code(aFilename, CPL_ERROR_NULL_INPUT);
477  cpl_error_code rc = CPL_ERROR_NONE;
478  if (!aImages || !aNames) {
479  return rc; /* this is a valid case, return without error */
480  }
481  unsigned int i, nimages = muse_imagelist_get_size(aImages);
482  for (i = 0; i < nimages; i++) {
483  muse_image *image = muse_imagelist_get(aImages, i);
484 
485  /* create small header with EXTNAME=filtername and WCS for images */
486  cpl_propertylist *header = cpl_propertylist_new();
487  cpl_errorstate es = cpl_errorstate_get();
488  const char *unit = muse_pfits_get_bunit(image->header),
489  *ucomment = cpl_propertylist_get_comment(image->header, "BUNIT");
490  if (!cpl_errorstate_is_equal(es) && !unit) {
491  cpl_errorstate_set(es); /* swallow errors, if there is no unit to propagate */
492  }
493  char dataext[KEYWORD_LENGTH], obj[KEYWORD_LENGTH],
494  *dqext = NULL, *statext = NULL;
495  snprintf(dataext, KEYWORD_LENGTH, "%s", cpl_array_get_string(aNames, i));
496  if (image->dq) {
497  dqext = cpl_sprintf("%s_%s", cpl_array_get_string(aNames, i), EXTNAME_DQ);
498  }
499  if (image->stat) {
500  statext = cpl_sprintf("%s_%s", cpl_array_get_string(aNames, i), EXTNAME_STAT);
501  }
502  snprintf(obj, KEYWORD_LENGTH, "%s", cpl_array_get_string(aNames, i));
503  cpl_propertylist_append_string(header, "EXTNAME", dataext);
504  cpl_propertylist_set_comment(header, "EXTNAME", "reconstructed image (data values)");
505  /* propagate the unit in the same way as in muse_image_save() */
506  if (unit) {
507  cpl_propertylist_append_string(header, "BUNIT", unit);
508  cpl_propertylist_set_comment(header, "BUNIT", ucomment);
509  }
510  /* update the OBJECT header with the filter name */
511  muse_utils_copy_modified_header(image->header, header, "OBJECT", obj);
512  cpl_propertylist_copy_property_regexp(header, image->header,
513  MUSE_WCS_KEYS"|"MUSE_HDR_FILTER_REGEXP,
514  0);
515  muse_utils_set_hduclass(header, "DATA", dataext, dqext, statext);
516  rc = cpl_image_save(image->data, aFilename, CPL_TYPE_UNSPECIFIED, header,
517  CPL_IO_EXTEND);
518 
519  if (image->dq) {
520  cpl_propertylist_update_string(header, "EXTNAME", dqext);
521  cpl_propertylist_set_comment(header, "EXTNAME", "reconstructed image (bad pixel status values)");
522  cpl_propertylist_erase(header, "BUNIT"); /* no unit for data quality */
523  snprintf(obj, KEYWORD_LENGTH, "%s, %s", cpl_array_get_string(aNames, i),
524  EXTNAME_DQ);
525  muse_utils_copy_modified_header(image->header, header, "OBJECT", obj);
526  muse_utils_set_hduclass(header, "QUALITY", dataext, dqext, statext);
527  rc = cpl_image_save(image->dq, aFilename, CPL_TYPE_INT, header,
528  CPL_IO_EXTEND);
529  }
530  if (image->stat) {
531  cpl_propertylist_update_string(header, "EXTNAME", statext);
532  cpl_propertylist_set_comment(header, "EXTNAME", "reconstructed image (variance)");
533  if (unit) {
534  char *ustat = cpl_sprintf("(%s)**2", unit); /* variance in squared units */
535  cpl_propertylist_append_string(header, "BUNIT", ustat);
536  cpl_free(ustat);
537  }
538  snprintf(obj, KEYWORD_LENGTH, "%s, %s", cpl_array_get_string(aNames, i),
539  EXTNAME_STAT);
540  muse_utils_copy_modified_header(image->header, header, "OBJECT", obj);
541  muse_utils_set_hduclass(header, "ERROR", dataext, dqext, statext);
542  rc = cpl_image_save(image->stat, aFilename, CPL_TYPE_UNSPECIFIED, header,
543  CPL_IO_EXTEND);
544  }
545 
546  cpl_propertylist_delete(header);
547  cpl_free(dqext);
548  cpl_free(statext);
549  } /* for i (all reconstructed images) */
550 
551  return rc;
552 } /* muse_datacube_save_recimages */
553 
554 /*---------------------------------------------------------------------------*/
569 /*---------------------------------------------------------------------------*/
570 cpl_error_code
571 muse_euro3dcube_save(muse_euro3dcube *aEuro3D, const char *aFilename)
572 {
573  cpl_ensure_code(aEuro3D && aFilename, CPL_ERROR_NULL_INPUT);
574 
575  /* save primary header and the data table in the first extension */
576  cpl_error_code rc = cpl_table_save(aEuro3D->dtable, aEuro3D->header,
577  aEuro3D->hdata, aFilename, CPL_IO_CREATE);
578  if (rc != CPL_ERROR_NONE) {
579  cpl_msg_warning(__func__, "failed to save data part of the Euro3D table: "
580  "%s", cpl_error_get_message());
581  }
582  /* save group header and data in the second extension */
583  rc = cpl_table_save(aEuro3D->gtable, NULL, aEuro3D->hgroup, aFilename,
584  CPL_IO_EXTEND);
585  if (rc != CPL_ERROR_NONE) {
586  cpl_msg_warning(__func__, "failed to save group part of the Euro3D table: "
587  "%s", cpl_error_get_message());
588  }
589 
590  rc = muse_datacube_save_recimages(aFilename, aEuro3D->recimages,
591  aEuro3D->recnames);
592  return rc;
593 } /* muse_euro3dcube_save() */
594 
595 /*---------------------------------------------------------------------------*/
605 /*---------------------------------------------------------------------------*/
606 void
608 {
609  /* if the euro3d object doesn't exists at all, we don't need to do anything */
610  if (!aEuro3D) {
611  return;
612  }
613 
614  /* checks for the existence of the sub-images *
615  * are done in the CPL functions */
616  cpl_table_delete(aEuro3D->dtable);
617  cpl_table_delete(aEuro3D->gtable);
618  cpl_propertylist_delete(aEuro3D->header);
619  cpl_propertylist_delete(aEuro3D->hdata);
620  cpl_propertylist_delete(aEuro3D->hgroup);
622  cpl_array_delete(aEuro3D->recnames);
623  cpl_free(aEuro3D);
624 } /* muse_euro3dcube_delete() */
625 
626 /*---------------------------------------------------------------------------*/
636 /*---------------------------------------------------------------------------*/
637 static cpl_error_code
638 muse_datacube_convert_dq_recimages(muse_datacube *aCube)
639 {
640  cpl_ensure_code(aCube, CPL_ERROR_NULL_INPUT);
641  if (!aCube->recimages) {
642  return CPL_ERROR_NONE; /* valid case, return without error */
643  }
644  unsigned int k, nimages = muse_imagelist_get_size(aCube->recimages);
645  for (k = 0; k < nimages; k++) {
646  muse_image *image = muse_imagelist_get(aCube->recimages, k);
647  if (!image->dq) { /* no point trying to convert... */
648  continue;
649  }
650  muse_image_dq_to_nan(image);
651  } /* for k (all reconstructed images) */
652 
653  return CPL_ERROR_NONE;
654 } /* muse_datacube_convert_dq_recimages() */
655 
656 /*---------------------------------------------------------------------------*/
669 /*---------------------------------------------------------------------------*/
670 cpl_error_code
672 {
673  cpl_ensure_code(aCube && aCube->data && aCube->stat && aCube->dq,
674  CPL_ERROR_NULL_INPUT);
675 
676  /* loop through all wavelength planes and then all pixels per plane */
677  int l, nx = cpl_image_get_size_x(cpl_imagelist_get(aCube->data, 0)),
678  ny = cpl_image_get_size_y(cpl_imagelist_get(aCube->data, 0)),
679  nz = cpl_imagelist_get_size(aCube->data);
680  #pragma omp parallel for default(none) /* as req. by Ralf */ \
681  shared(aCube, nx, ny, nz)
682  for (l = 0; l < nz; l++) {
683  int i;
684  for (i = 0; i < nx; i++) {
685  int j;
686  for (j = 0; j < ny; j++) {
687  float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->data, l)),
688  *pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->stat, l));
689  int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->dq, l));
690  if (pdq[i + j*nx] == EURO3D_GOODPIXEL) {
691  continue; /* nothing to do for good pixels */
692  }
693  /* set bad pixels of any type to not-a-number in both extensions */
694  pdata[i + j*nx] = NAN; /* supposed to be quiet NaN */
695  pstat[i + j*nx] = NAN;
696  } /* for j (y direction) */
697  } /* for i (x direction) */
698  } /* for l (wavelength planes) */
699 
700  /* deallocate DQ and set it to NULL to be able to check for it */
701  cpl_imagelist_delete(aCube->dq);
702  aCube->dq = NULL;
703 
704  /* do the same for the reconstructed images, if there are any */
705  muse_datacube_convert_dq_recimages(aCube);
706 
707  return CPL_ERROR_NONE;
708 } /* muse_datacube_convert_dq() */
709 
710 /*---------------------------------------------------------------------------*/
739 /*---------------------------------------------------------------------------*/
740 cpl_error_code
741 muse_datacube_save(muse_datacube *aCube, const char *aFilename)
742 {
743  cpl_ensure_code(aCube && aCube->header && aFilename, CPL_ERROR_NULL_INPUT);
744 
745  /* save headers into primary area */
746  cpl_propertylist *header = cpl_propertylist_new();
747  /* copy all headers, except the main WCS keys */
748  cpl_propertylist_copy_property_regexp(header, aCube->header,
749  MUSE_WCS_KEYS, 1);
750  cpl_error_code rc = cpl_propertylist_save(header, aFilename, CPL_IO_CREATE);
751  cpl_propertylist_delete(header);
752 
753  /* save data */
754  header = cpl_propertylist_new();
755  cpl_propertylist_append_string(header, "EXTNAME", EXTNAME_DATA);
756  cpl_propertylist_set_comment(header, "EXTNAME", EXTNAME_DATA_COMMENT);
757  muse_utils_copy_modified_header(aCube->header, header, "OBJECT",
758  EXTNAME_DATA);
759  cpl_propertylist_copy_property_regexp(header, aCube->header,
760  MUSE_WCS_KEYS"|^BUNIT", 0);
761  muse_utils_set_hduclass(header, "DATA", "DATA",
762  aCube->dq ? "DQ" : NULL, aCube->stat ? "STAT" : NULL);
763  rc = cpl_imagelist_save(aCube->data, aFilename, CPL_TYPE_FLOAT, header,
764  CPL_IO_EXTEND);
765  cpl_propertylist_delete(header);
766  /* save bad pixels, if available */
767  if (rc == CPL_ERROR_NONE && aCube->dq) {
768  header = cpl_propertylist_new();
769  cpl_propertylist_append_string(header, "EXTNAME", EXTNAME_DQ);
770  cpl_propertylist_set_comment(header, "EXTNAME", EXTNAME_DQ_COMMENT);
771  muse_utils_copy_modified_header(aCube->header, header, "OBJECT",
772  EXTNAME_DQ);
773  cpl_propertylist_copy_property_regexp(header, aCube->header,
774  MUSE_WCS_KEYS, 0);
775  muse_utils_set_hduclass(header, "QUALITY", "DATA", "DQ",
776  aCube->stat ? "STAT" : NULL);
777  rc = cpl_imagelist_save(aCube->dq, aFilename, CPL_TYPE_INT, header,
778  CPL_IO_EXTEND);
779  cpl_propertylist_delete(header);
780  }
781  /* save variance, if available */
782  if (rc == CPL_ERROR_NONE && aCube->stat) {
783  header = cpl_propertylist_new();
784  cpl_propertylist_append_string(header, "EXTNAME", EXTNAME_STAT);
785  cpl_propertylist_set_comment(header, "EXTNAME", EXTNAME_STAT_COMMENT);
786  /* add the correct BUNIT, depending on the data units */
787  if (!strncmp(muse_pfits_get_bunit(aCube->header),
788  kMuseFluxUnitString, strlen(kMuseFluxUnitString) + 1)) {
789  /* flux calibrated data */
790  cpl_propertylist_append_string(header, "BUNIT", kMuseFluxStatString);
791  }
792  muse_utils_copy_modified_header(aCube->header, header, "OBJECT",
793  EXTNAME_STAT);
794  cpl_propertylist_copy_property_regexp(header, aCube->header,
795  MUSE_WCS_KEYS, 0);
796  muse_utils_set_hduclass(header, "ERROR", "DATA", aCube->dq ? "DQ" : NULL,
797  "STAT");
798  rc = cpl_imagelist_save(aCube->stat, aFilename, CPL_TYPE_FLOAT, header,
799  CPL_IO_EXTEND);
800  cpl_propertylist_delete(header);
801  }
802 
803  rc = muse_datacube_save_recimages(aFilename, aCube->recimages, aCube->recnames);
804  return rc;
805 } /* muse_datacube_save() */
806 
807 /*---------------------------------------------------------------------------*/
808 /*
809  @private
810  @error{set CPL_ERROR_NULL_INPUT\, return NULL, aFilename is NULL}
811  @error{set CPL_ERROR_ILLEGAL_INPUT\, return NULL, aFilename does not exist}
812  @error{set CPL_ERROR_DATA_NOT_FOUND\, return NULL,
813  aFilename does not contain a DATA extension}
814  */
815 /*---------------------------------------------------------------------------*/
816 static cpl_propertylist *
817 muse_datacube_load_header(const char *aFilename)
818 {
819  cpl_ensure(aFilename, CPL_ERROR_NULL_INPUT, NULL);
820  int extdata = cpl_fits_find_extension(aFilename, "DATA");
821  cpl_ensure(extdata >= 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
822  cpl_ensure(extdata > 0, CPL_ERROR_DATA_NOT_FOUND, NULL);
823 
824  cpl_propertylist *header = cpl_propertylist_load(aFilename, 0);
825  cpl_propertylist *hdata = cpl_propertylist_load(aFilename, extdata);
826  cpl_propertylist_copy_property_regexp(header, hdata,
827  MUSE_WCS_KEYS"|BUNIT", 0);
828  cpl_propertylist_delete(hdata);
829  return header;
830 } /* muse_datacube_load_header() */
831 
832 /*---------------------------------------------------------------------------*/
862 /*---------------------------------------------------------------------------*/
864 muse_datacube_load(const char *aFilename)
865 {
866  cpl_ensure(aFilename, CPL_ERROR_NULL_INPUT, NULL);
867  muse_datacube *cube = cpl_calloc(1, sizeof(muse_datacube));
868  /* load primary header and merge in relevant entries from the DATA header */
869  cpl_errorstate state = cpl_errorstate_get();
870  cube->header = muse_datacube_load_header(aFilename);
871  if (!cpl_errorstate_is_equal(state) || !cube->header) {
872  cpl_msg_error(__func__, "Loading cube-like headers from \"%s\" failed!",
873  aFilename);
874  cpl_free(cube);
875  return NULL;
876  }
877  int ext = cpl_fits_find_extension(aFilename, "DATA");
878  cube->data = cpl_imagelist_load(aFilename, CPL_TYPE_FLOAT, ext);
879  /* DQ is usually not written to disk, but see if it's there and load it, if so */
880  ext = cpl_fits_find_extension(aFilename, "DQ");
881  if (ext > 0) {
882  cube->stat = cpl_imagelist_load(aFilename, CPL_TYPE_INT, ext);
883  }
884  ext = cpl_fits_find_extension(aFilename, "STAT");
885  if (ext > 0) {
886  cube->stat = cpl_imagelist_load(aFilename, CPL_TYPE_FLOAT, ext);
887  }
888  int next = cpl_fits_count_extensions(aFilename);
889  while (++ext <= next) {
890  /* there is (one more) collapsed image, load it */
891  muse_image *image = muse_image_new();
892  image->header = cpl_propertylist_load(aFilename, ext);
893  image->data = cpl_image_load(aFilename, CPL_TYPE_UNSPECIFIED, 0, ext);
894  /* search for an load DQ and STAT extensions of this image, if present */
895  const char *extname = muse_pfits_get_extname(image->header);
896  char *extname2 = cpl_sprintf("%s_DQ", extname);
897  int ext2 = cpl_fits_find_extension(aFilename, extname2);
898  cpl_free(extname2);
899  if (ext2 > 0) {
900  image->dq = cpl_image_load(aFilename, CPL_TYPE_INT, 0, ext2);
901  ext = ext2; /* skip this as a normal reconstructed image */
902  }
903  extname2 = cpl_sprintf("%s_STAT", extname);
904  ext2 = cpl_fits_find_extension(aFilename, extname2);
905  cpl_free(extname2);
906  if (ext2 > 0) {
907  image->stat = cpl_image_load(aFilename, CPL_TYPE_UNSPECIFIED, 0, ext2);
908  ext = ext2; /* skip this as a normal reconstructed image */
909  }
910 
911  if (!cube->recnames) {
912  cube->recnames = cpl_array_new(1, CPL_TYPE_STRING);
913  } else {
914  cpl_array_set_size(cube->recnames, cpl_array_get_size(cube->recnames) + 1);
915  }
916  /* append name of this reconstructed image at the end of the names array */
917  cpl_array_set_string(cube->recnames, cpl_array_get_size(cube->recnames) - 1,
918  extname);
919  if (!cube->recimages) {
920  cube->recimages = muse_imagelist_new();
921  }
922  /* append the image to the end of the reconstructed image list */
923  muse_imagelist_set(cube->recimages, image,
925  }
926  return cube;
927 } /* muse_datacube_load() */
928 
929 /*---------------------------------------------------------------------------*/
949 /*---------------------------------------------------------------------------*/
950 cpl_error_code
952 {
953  cpl_ensure_code(aCube && aAppend, CPL_ERROR_NULL_INPUT);
954 
955  cpl_size ndata = cpl_imagelist_get_size(aCube->data),
956  nstat = cpl_imagelist_get_size(aCube->stat);
957  cpl_ensure_code(ndata == nstat, CPL_ERROR_ILLEGAL_INPUT);
958  cpl_size i, n = cpl_imagelist_get_size(aAppend->data);
959  cpl_ensure_code(n == cpl_imagelist_get_size(aAppend->stat),
960  CPL_ERROR_ILLEGAL_INPUT);
961 
962  /* check sizes of the image planes */
963  cpl_size nx1 = cpl_image_get_size_x(cpl_imagelist_get(aCube->data, ndata - 1)),
964  ny1 = cpl_image_get_size_y(cpl_imagelist_get(aCube->data, ndata - 1)),
965  nx2 = cpl_image_get_size_x(cpl_imagelist_get(aAppend->data, 0)),
966  ny2 = cpl_image_get_size_y(cpl_imagelist_get(aAppend->data, 0));
967  cpl_ensure_code(nx1 == nx2 && ny1 == ny2, CPL_ERROR_ILLEGAL_INPUT);
968  const char *ctype1 = muse_pfits_get_ctype(aCube->header, 3),
969  *ctype2 = muse_pfits_get_ctype(aCube->header, 3);
970  cpl_ensure_code(ctype1 && ctype2, CPL_ERROR_UNSUPPORTED_MODE);
971  cpl_ensure_code((!strncmp(ctype1, "AWAV", 5) && !strncmp(ctype2, "AWAV", 5)) ||
972  (!strncmp(ctype1, "WAVE", 5) && !strncmp(ctype2, "WAVE", 5)),
973  CPL_ERROR_UNSUPPORTED_MODE);
974  /* compute wavelengths of last plane of first cube and first plane *
975  * of the cube to append, to check that they are adjacent */
976  double pix1 = muse_pfits_get_crpix(aCube->header, 3),
977  val1 = muse_pfits_get_crval(aCube->header, 3),
978  cd1 = muse_pfits_get_cd(aCube->header, 3, 3),
979  pix2 = muse_pfits_get_crpix(aAppend->header, 3),
980  val2 = muse_pfits_get_crval(aAppend->header, 3),
981  cd2 = muse_pfits_get_cd(aAppend->header, 3, 3),
982  lbda1 = ((ndata - 1) + 1. - pix1) * cd1 + val1,
983  lbda2 = (0 + 1. - pix2) * cd2 + val2;
984  /* AWAV-LOG and WAVE-LOG axes need more complex handling, so exclude them above */
985  cpl_msg_debug(__func__, "lambdas: %f %f (%f %f)", lbda1, lbda2, cd1, cd2);
986  cpl_ensure_code(fabs(cd1 - cd2) < DBL_EPSILON &&
987  fabs(lbda2 - cd2 - lbda1) < DBL_EPSILON,
988  CPL_ERROR_ILLEGAL_INPUT);
989 
990  /* reconstructed images are now obsolete, remove them */
991  cpl_array_delete(aCube->recnames);
992  aCube->recnames = NULL;
994  aCube->recimages = NULL;
995 
996  /* Check if there is a usable DQ component in both cubes. If so, *
997  * propagate them as the other two components, otherwise delete *
998  * the incoming DQ of the first (= output) cube. */
999  cpl_boolean usedq = CPL_FALSE;
1000  if (aCube->dq && cpl_imagelist_get_size(aCube->dq) == ndata &&
1001  aAppend->dq && cpl_imagelist_get_size(aAppend->dq) == n) {
1002  usedq = CPL_TRUE;
1003  } else {
1004  cpl_imagelist_delete(aCube->dq);
1005  aCube->dq = NULL;
1006  }
1007 
1008  /* append the CPL imagelists */
1009  for (i = 0; i < n; i++) {
1010  cpl_imagelist_set(aCube->data,
1011  cpl_image_duplicate(cpl_imagelist_get(aAppend->data, i)),
1012  cpl_imagelist_get_size(aCube->data));
1013  if (usedq) {
1014  cpl_imagelist_set(aCube->dq,
1015  cpl_image_duplicate(cpl_imagelist_get(aAppend->dq, i)),
1016  cpl_imagelist_get_size(aCube->dq));
1017  }
1018  cpl_imagelist_set(aCube->stat,
1019  cpl_image_duplicate(cpl_imagelist_get(aAppend->stat, i)),
1020  cpl_imagelist_get_size(aCube->stat));
1021  }
1022  return CPL_ERROR_NONE;
1023 } /* muse_datacube_concat() */
1024 
1025 /*---------------------------------------------------------------------------*/
1035 /*---------------------------------------------------------------------------*/
1036 void
1038 {
1039  /* if the cube does not exists at all, we don't need to do anything */
1040  if (!aCube) {
1041  return;
1042  }
1043 
1044  /* checks for the existence of the sub-images *
1045  * are done in the CPL functions */
1046  cpl_imagelist_delete(aCube->data);
1047  aCube->data = NULL;
1048  cpl_imagelist_delete(aCube->dq);
1049  aCube->dq = NULL;
1050  cpl_imagelist_delete(aCube->stat);
1051  aCube->stat = NULL;
1052 
1053  /* delete the FITS header, too */
1054  cpl_propertylist_delete(aCube->header);
1055  aCube->header = NULL;
1056 
1057  /* remove reconstructed image data, if present */
1059  cpl_array_delete(aCube->recnames);
1060  cpl_free(aCube);
1061 } /* muse_datacube_delete() */
1062 
Structure definition of a MUSE datacube.
Definition: muse_datacube.h:48
Structure definition for a collection of muse_images.
double muse_pfits_get_cd(const cpl_propertylist *aHeaders, unsigned int aAxisI, unsigned int aAxisJ)
find out the WCS coordinate at the reference point
Definition: muse_pfits.c:446
const char * muse_pfits_get_extname(const cpl_propertylist *aHeaders)
find out the extension name
Definition: muse_pfits.c:205
cpl_error_code muse_datacube_concat(muse_datacube *aCube, const muse_datacube *aAppend)
Concatenate one datacube at the end of another one.
double muse_pfits_get_crval(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS coordinate at the reference point
Definition: muse_pfits.c:423
A structure containing a spatial two-axis WCS.
Definition: muse_wcs.h:105
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...
cpl_error_code muse_utils_filter_copy_properties(cpl_propertylist *aHeader, const muse_table *aFilter, double aFraction)
Add/propagate filter properties to header of collapsed image.
Definition: muse_utils.c:934
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:245
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
const char * muse_pfits_get_bunit(const cpl_propertylist *aHeaders)
find out the unit string
Definition: muse_pfits.c:223
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:71
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:2611
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_pfits_get_crpix(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS reference point
Definition: muse_pfits.c:401
double muse_flux_response_interpolate(const cpl_table *aResponse, double aLambda, double *aError, muse_flux_interpolation_type aType)
Compute linearly interpolated response of some kind at given wavelength.
Definition: muse_flux.c:338
muse_wcs * muse_wcs_new(cpl_propertylist *aHeader)
Create a new WCS structure from a given FITS header.
Definition: muse_wcs.c:1870
cpl_error_code muse_datacube_convert_dq(muse_datacube *aCube)
Convert the DQ extension of a datacube to NANs in DATA and STAT.
#define MUSE_WCS_KEYS
regular expression for WCS properties
Definition: muse_wcs.h:48
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:176
cpl_table * table
The table.
Definition: muse_table.h:49
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:76
Structure definition of a Euro3D datacube.
Definition: muse_datacube.h:97
Structure to store a table together with a property list.
Definition: muse_table.h:43
double muse_utils_filter_fraction(const muse_table *aFilter, double aLambda1, double aLambda2)
Compute fraction of filter area covered by the data range.
Definition: muse_utils.c:893
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
const char * muse_pfits_get_ctype(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS axis type string
Definition: muse_pfits.c:469
cpl_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_imagelist * muse_imagelist_new(void)
Create a new (empty) MUSE image list.
cpl_boolean iscelsph
Definition: muse_wcs.h:112
muse_image * muse_datacube_collapse(muse_datacube *aCube, const muse_table *aFilter)
Integrate a FITS NAXIS=3 datacube along the wavelength direction.
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:107
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: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:2570
muse_image * muse_euro3dcube_collapse(muse_euro3dcube *aEuro3D, const muse_table *aFilter)
Integrate a Euro3D datacube along the wavelength direction.
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