MUSE Pipeline Reference Manual  0.18.5
muse_resampling.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 #ifndef _OPENMP
33 #define omp_get_max_threads() 1
34 #else
35 #include <omp.h>
36 #endif
37 
38 #include "muse_resampling.h"
39 #include "muse_instrument.h"
40 
41 #include "muse_astro.h"
42 #include "muse_cplwrappers.h"
43 #include "muse_dar.h"
44 #include "muse_dfs.h"
45 #include "muse_flux.h"
46 #include "muse_pfits.h"
47 #include "muse_pixgrid.h"
48 #include "muse_pixtable.h"
49 #include "muse_quality.h"
50 #include "muse_utils.h"
51 #include "muse_wcs.h"
52 #include "muse_data_format_z.h"
53 
54 /*----------------------------------------------------------------------------*/
58 /*----------------------------------------------------------------------------*/
59 
62 /*----------------------------------------------------------------------------*
63  * Resampling, collapsing, and saving in 3D *
64  *----------------------------------------------------------------------------*/
65 
66 /*---------------------------------------------------------------------------*/
79 /*---------------------------------------------------------------------------*/
82 {
83  cpl_ensure(aMethod <= MUSE_RESAMPLE_NONE, CPL_ERROR_ILLEGAL_INPUT, NULL);
84  muse_resampling_params *params = cpl_calloc(1, sizeof(muse_resampling_params));
85  params->method = aMethod;
86  /* leave crtype zero == MUSE_RESAMPLING_CRSTATS_IRAF by default */
87  /* leave crsigma zero, to not do CR rejection by default */
88  params->ld = 1;
89  params->pfx = 0.6;
90  params->pfy = 0.6;
91  params->pfl = 0.6;
92  params->rc = 1.25;
93  /* dx, dy, and dlambda, and wcs stay zero to trigger defaults */
94  /* leave tlambda zero == MUSE_RESAMPLING_DISP_AWAV by default */
95  return params;
96 } /* muse_resampling_params_new() */
97 
98 /*---------------------------------------------------------------------------*/
117 /*---------------------------------------------------------------------------*/
118 cpl_error_code
120  const cpl_propertylist *aWCS)
121 {
122  cpl_ensure_code(aParams, CPL_ERROR_NULL_INPUT);
123  if (!aWCS) {
125  cpl_wcs_delete(aParams->wcs);
126  aParams->wcs = NULL;
127  return CPL_ERROR_NONE;
128  }
129  if (cpl_propertylist_has(aWCS, "CTYPE3") &&
130  !strncmp(cpl_propertylist_get_string(aWCS, "CTYPE3"), "AWAV-LOG", 9)) {
132  } else { /* make sure it wasn't set to rubbish elsewhere */
134  }
135  cpl_errorstate state = cpl_errorstate_get();
136  cpl_error_code rc = CPL_ERROR_NONE;
137  aParams->wcs = cpl_wcs_new_from_propertylist(aWCS);
138  if (!cpl_errorstate_is_equal(state)) {
139  cpl_wcs_delete(aParams->wcs);
140  aParams->wcs = NULL;
141  rc = cpl_error_get_code();
142  }
143 #if 0
144  printf("CDi_j matrix:\n");
145  cpl_matrix_dump(cpl_wcs_get_cd(aParams->wcs), stdout);
146  fflush(stdout);
147 #endif
148  return rc;
149 } /* muse_resampling_params_set_wcs() */
150 
151 /*---------------------------------------------------------------------------*/
159 /*---------------------------------------------------------------------------*/
160 void
162 {
163  if (!aParams) {
164  return;
165  }
166  cpl_wcs_delete(aParams->wcs);
167  cpl_free(aParams);
168 } /* muse_resampling_params_delete() */
169 
170 /*---------------------------------------------------------------------------*/
181 /*---------------------------------------------------------------------------*/
182 static inline double
183 muse_resampling_weight_function_linear(double r)
184 {
185  return r == 0 ? FLT_MAX : 1. / r;
186 }
187 
188 /*---------------------------------------------------------------------------*/
199 /*---------------------------------------------------------------------------*/
200 static inline double
201 muse_resampling_weight_function_quadratic(double r2)
202 {
203  return r2 == 0 ? FLT_MAX : 1. / r2;
204 }
205 
206 /*---------------------------------------------------------------------------*/
219 /*---------------------------------------------------------------------------*/
220 static inline double
221 muse_resampling_weight_function_renka(double r, double r_c)
222 {
223  if (r == 0) {
224  return FLT_MAX;
225  } else if (r >= r_c) {
226  return DBL_MIN;
227  } else {
228  double p = (r_c - r) / (r_c * r);
229  return p*p;
230  }
231 }
232 
233 /*---------------------------------------------------------------------------*/
240 /*---------------------------------------------------------------------------*/
241 static inline double
242 muse_resampling_weight_function_sinc(double r)
243 {
244  return fabs(r) < DBL_EPSILON ? 1. : sin(CPL_MATH_PI * r) / (CPL_MATH_PI * r);
245 }
246 
247 /*---------------------------------------------------------------------------*/
257 /*---------------------------------------------------------------------------*/
258 static inline double
259 muse_resampling_weight_function_lanczos(double dx, double dy, double dz, unsigned int n)
260 {
261  return (fabs(dx) >= n || fabs(dy) >= n || fabs(dz) > n) ? 0.
262  : muse_resampling_weight_function_sinc(dx) * muse_resampling_weight_function_sinc(dx / n)
263  * muse_resampling_weight_function_sinc(dy) * muse_resampling_weight_function_sinc(dy / n)
264  * muse_resampling_weight_function_sinc(dz) * muse_resampling_weight_function_sinc(dz / n);
265 }
266 
267 /*---------------------------------------------------------------------------*/
286 /*---------------------------------------------------------------------------*/
287 static inline double
288 muse_resampling_weight_function_drizzle(double xin, double yin, double zin,
289  double xout, double yout, double zout,
290  double dx, double dy, double dz)
291 {
292  /* compute the three terms in the numerator: if the offset *
293  * plus the output halfsize is less than the input halfsize, *
294  * then that side is fully contained in the input pixel */
295  double x = (dx + xout / 2.) <= xin / 2. ? xout : (xin + xout) / 2. - dx,
296  y = (dy + yout / 2.) <= yin / 2. ? yout : (yin + yout) / 2. - dy,
297  z = (dz + zout / 2.) <= zin / 2. ? zout : (zin + zout) / 2. - dz;
298  /* any negative value means that the input pixel is completely *
299  * outside the target voxel, so it doesn't contribute */
300  if (x <= 0 || y <= 0 || z <= 0) {
301  return 0.;
302  }
303  /* any value > input size means this dimension of the input pixel *
304  * is completely inside the target voxel, so that's the limit! */
305  return (x > xin ? xin : x) * (y > yin ? yin : y) * (z > zin ? zin : z)
306  / (xin * yin * zin);
307 } /* muse_resampling_weight_function_drizzle() */
308 
309 /*---------------------------------------------------------------------------*/
325 /*---------------------------------------------------------------------------*/
326 static cpl_error_code
327 muse_resampling_cube_nearest(muse_datacube *aCube, muse_pixtable *aPixtable,
328  muse_pixgrid *aPixgrid)
329 {
330  cpl_ensure_code(aCube && aPixtable && aPixgrid, CPL_ERROR_NULL_INPUT);
331  double ptxoff = 0., /* zero by default ... */
332  ptyoff = 0., /* for pixel coordinates */
333  crval3 = cpl_propertylist_get_double(aCube->header, "CRVAL3"),
334  crpix3 = cpl_propertylist_get_double(aCube->header, "CRPIX3"),
335  cd33 = cpl_propertylist_get_double(aCube->header, "CD3_3");
336  muse_wcs *wcs = muse_wcs_new(aCube->header);
338  cpl_boolean loglambda = !strncmp(cpl_propertylist_get_string(aCube->header,
339  "CTYPE3"),
340  "AWAV-LOG", 9);
341  float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
342  *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
343  *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
344  *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
345  *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
346  int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
347 
348  /* If our data was astrometrically calibrated, we need to scale the *
349  * data units to the pixel size in all three dimensions so that the *
350  * radius computation works again. *
351  * Otherwise dx~5.6e-5deg won't contribute to the weighting at all. */
352  double xnorm = 1., ynorm = 1., znorm = 1. / kMuseSpectralSamplingA;
353  if (wcs->iscelsph) {
354  muse_wcs_get_scales(aPixtable->header, &xnorm, &ynorm);
355  xnorm = 1. / xnorm;
356  ynorm = 1. / ynorm;
357  /* need to use the real coordinate offset for celestial spherical */
358  ptxoff = cpl_propertylist_get_double(aPixtable->header, "CRVAL1");
359  ptyoff = cpl_propertylist_get_double(aPixtable->header, "CRVAL2");
360  }
361 
362 #ifdef ESO_ENABLE_DEBUG
363  int debug = 0;
364  if (getenv("MUSE_DEBUG_NEAREST")) {
365  debug = atoi(getenv("MUSE_DEBUG_NEAREST"));
366  }
367 #endif
368 
369  int l;
370 #ifdef ESO_ENABLE_DEBUG
371  #pragma omp parallel for default(none) /* as req. by Ralf */ \
372  shared(aCube, aPixgrid, cd33, crpix3, crval3, data, debug, dq, lbda, \
373  loglambda, ptxoff, ptyoff, stat, stdout, wcs, xnorm, xpos, \
374  ynorm, ypos, znorm)
375 #else
376  #pragma omp parallel for default(none) /* as req. by Ralf */ \
377  shared(aCube, aPixgrid, cd33, crpix3, crval3, data, dq, lbda, \
378  loglambda, ptxoff, ptyoff, stat, stdout, wcs, xnorm, xpos, \
379  ynorm, ypos, znorm)
380 #endif
381  for (l = 0; l < aPixgrid->size_z; l++) {
382  float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->data, l)),
383  *pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->stat, l));
384  int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->dq, l));
385  /* wavelength of center of current grid cell (l is index starting at 0) */
386  double lambda = (l + 1. - crpix3) * cd33 + crval3;
387  if (loglambda) {
388  lambda = crval3 * exp((l + 1. - crpix3) * cd33 / crval3);
389  }
390 
391  int i;
392  for (i = 0; i < aPixgrid->size_x; i++) {
393  int j;
394  for (j = 0; j < aPixgrid->size_y; j++) {
395  cpl_size idx = muse_pixgrid_get_index(aPixgrid, i, j, l, CPL_FALSE),
396  n_rows = muse_pixgrid_get_count(aPixgrid, idx);
397  const cpl_size *rows = muse_pixgrid_get_rows(aPixgrid, idx);
398 
399  /* x and y position of center of current grid cell (i, j start at 0) */
400  double x, y;
401  if (wcs->iscelsph) {
402  muse_wcs_celestial_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
403  } else {
404  muse_wcs_projplane_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
405  }
406 
407  if (n_rows == 1) {
408  /* if there is only one pixel in the cell, just use it */
409  pdata[i + j * aPixgrid->size_x] = data[rows[0]];
410  pstat[i + j * aPixgrid->size_x] = stat[rows[0]];
411  pdq[i + j * aPixgrid->size_x] = dq[rows[0]];
412 
413 #ifdef ESO_ENABLE_DEBUG
414  if (debug) {
415  printf("only: %f,%f\n", data[rows[0]], stat[rows[0]]);
416  fflush(stdout);
417  }
418 #endif
419  } else if (n_rows >= 2) {
420  /* loop through all available values and take the closest one */
421  cpl_size n, nbest = -1;
422  double dbest = FLT_MAX; /* some unlikely large value to start with*/
423  for (n = 0; n < n_rows; n++) {
424  /* the differences for the cell center and the current pixel */
425  double dx = fabs(x - xpos[rows[n]] + ptxoff) * xnorm,
426  dy = fabs(y - ypos[rows[n]] + ptyoff) * ynorm,
427  dlambda = fabs(lambda - lbda[rows[n]]) * znorm,
428  dthis = sqrt(dx*dx + dy*dy + dlambda*dlambda);
429  if (wcs->iscelsph) {
430  /* Not strictly necessary for NN, but still scale the RA *
431  * distance properly, see muse_resampling_cube_weighted(). */
432  dx *= cos(y * CPL_MATH_RAD_DEG);
433  }
434 #ifdef ESO_ENABLE_DEBUG
435  if (debug) {
436  printf("%f %f %f, %f %f %f, d: %f %f %f -> %f best: %f (%f,%f)\n",
437  x, y, lambda, xpos[rows[n]] + ptxoff, ypos[rows[n]] + ptyoff,
438  lbda[rows[n]], dx, dy, dlambda, dthis, dbest, data[rows[n]],
439  stat[rows[n]]);
440  }
441 #endif
442  if (dthis < dbest) {
443  nbest = n;
444  dbest = dthis;
445  }
446  }
447  pdata[i + j * aPixgrid->size_x] = data[rows[nbest]];
448  pstat[i + j * aPixgrid->size_x] = stat[rows[nbest]];
449  pdq[i + j * aPixgrid->size_x] = dq[rows[nbest]];
450 #ifdef ESO_ENABLE_DEBUG
451  if (debug) {
452  printf("closest: %f,%f\n", data[rows[nbest]], stat[rows[nbest]]);
453  fflush(stdout);
454  }
455 #endif
456  } else {
457  /* npix == 0: do nothing, pixel stays zero */
458  pdq[i + j * aPixgrid->size_x] = EURO3D_MISSDATA;
459  }
460  } /* for j (y direction) */
461  } /* for i (x direction) */
462  } /* for l (wavelength planes) */
463  cpl_free(wcs);
464 
465  return CPL_ERROR_NONE;
466 } /* muse_resampling_cube_nearest() */
467 
468 /*---------------------------------------------------------------------------*/
489 /*---------------------------------------------------------------------------*/
490 static cpl_error_code
491 muse_resampling_cube_weighted(muse_datacube *aCube, muse_pixtable *aPixtable,
492  muse_pixgrid *aPixgrid,
493  muse_resampling_params *aParams)
494 {
495  cpl_ensure_code(aCube && aPixtable && aPixgrid && aParams,
496  CPL_ERROR_NULL_INPUT);
497  double ptxoff = 0., /* zero by default ... */
498  ptyoff = 0., /* for pixel coordinates */
499  crval3 = cpl_propertylist_get_double(aCube->header, "CRVAL3"),
500  crpix3 = cpl_propertylist_get_double(aCube->header, "CRPIX3"),
501  cd33 = cpl_propertylist_get_double(aCube->header, "CD3_3");
502  muse_wcs *wcs = muse_wcs_new(aCube->header);
504  cpl_boolean loglambda = !strncmp(cpl_propertylist_get_string(aCube->header,
505  "CTYPE3"),
506  "AWAV-LOG", 9);
507  cpl_errorstate prestate = cpl_errorstate_get();
508  float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
509  *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
510  *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
511  *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
512  *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT),
513  *wght = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_WEIGHT);
514  if (!cpl_errorstate_is_equal(prestate)) {
515  cpl_errorstate_set(prestate); /* recover from missing weight column */
516  }
517  int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
518 
519  /* If our data was astrometrically calibrated, we need to scale the *
520  * data units to the pixel size in all three dimensions so that the *
521  * radius computation works again. *
522  * Otherwise dx~5.6e-5deg won't contribute to the weighting at all. */
523  double xnorm = 1., ynorm = 1., znorm = 1. / kMuseSpectralSamplingA;
524  if (wcs->iscelsph) {
525  muse_wcs_get_scales(aPixtable->header, &xnorm, &ynorm);
526  xnorm = 1. / xnorm;
527  ynorm = 1. / ynorm;
528  /* need to use the real coordinate offset for celestial spherical */
529  ptxoff = cpl_propertylist_get_double(aPixtable->header, "CRVAL1");
530  ptyoff = cpl_propertylist_get_double(aPixtable->header, "CRVAL2");
531  }
532  /* precomputed factor for output voxel size calculation in *
533  * wavelength direction, only needed for log-lambda axis */
534  double zoutefac = exp(1.5 * cd33 / crval3) - exp(0.5 * cd33 / crval3);
535  /* scale the input critical radius by the voxel radius */
536  double renka_rc = aParams->rc /* XXX beware of rotation! */
537  * sqrt((wcs->cd11*xnorm)*(wcs->cd11*xnorm) + (wcs->cd22*ynorm)*(wcs->cd22*ynorm)
538  + (cd33*znorm)*(cd33*znorm));
539  /* loop distance (to take into account surrounding pixels) verification */
540  int ld = aParams->ld;
541  if (ld <= 0) {
542  ld = 1;
543  cpl_msg_info(__func__, "Overriding loop distance ld=%d", ld);
544  }
545 
546  /* pixel sizes in all three directions, scaled by pixfrac, and *
547  * output pixel sizes (absolute values), as needed for drizzle */
548  double xsz = aParams->pfx / xnorm,
549  ysz = aParams->pfy / ynorm,
550  zsz = aParams->pfl / znorm,
551  xout = fabs(wcs->cd11), yout = fabs(wcs->cd22), zout = fabs(cd33);
552 
553 #ifdef ESO_ENABLE_DEBUG
554  int debug = 0, debugx = 0, debugy = 0, debugz = 0;
555  if (getenv("MUSE_DEBUG_WEIGHTED")) {
556  debug = atoi(getenv("MUSE_DEBUG_WEIGHTED"));
557  }
558  if (debug & 2) { /* need coordinates */
559  if (getenv("MUSE_DEBUG_WEIGHTED_X")) {
560  debugx = atoi(getenv("MUSE_DEBUG_WEIGHTED_X"));
561  if (debugx < 1 || debugx > aPixgrid->size_x) {
562  debugx = 0;
563  }
564  }
565  if (getenv("MUSE_DEBUG_WEIGHTED_Y")) {
566  debugy = atoi(getenv("MUSE_DEBUG_WEIGHTED_Y"));
567  if (debugy < 1 || debugy > aPixgrid->size_y) {
568  debugy = 0;
569  }
570  }
571  if (getenv("MUSE_DEBUG_WEIGHTED_Z")) {
572  debugz = atoi(getenv("MUSE_DEBUG_WEIGHTED_Z"));
573  if (debugz < 1 || debugz > aPixgrid->size_z) {
574  debugz = 0;
575  }
576  }
577  }
578  if (debug & 8) {
579  printf("parameters:\n cd=%e %e %e\n"
580  " corrected crpix3=%e\n norm=%e %e %e\n",
581  wcs->cd11, wcs->cd22, cd33, crpix3, xnorm, ynorm, znorm);
582  if (aParams->method == MUSE_RESAMPLE_WEIGHTED_DRIZZLE) {
583  printf(" drop sizes %e %e %e (pixfrac %f,%f,%f)\n"
584  " output sizes %e %e %e\n",
585  xsz, ysz, zsz, aParams->pfx, aParams->pfy, aParams->pfl,
586  xout, yout, zout);
587  } else {
588  printf(" resulting renka_rc: %e %e %e / %e %e %e --> %e\n",
589  pow(wcs->cd11, 2), pow(wcs->cd22, 2), cd33*cd33,
590  pow(wcs->cd11*xnorm, 2), pow(wcs->cd22*ynorm, 2),
591  pow(cd33*znorm, 2), renka_rc);
592  }
593  fflush(stdout);
594  }
595 #endif
596  cpl_imagelist *wcube = NULL;
597  if (getenv("MUSE_DEBUG_WEIGHT_CUBE")) { /* create a weight cube */
598  cpl_msg_debug(__func__, "Weighted resampling: creating weight cube");
599  wcube = cpl_imagelist_new();
600  int i;
601  for (i = 0; i < aPixgrid->size_z; i++) {
602  cpl_image *image = cpl_image_new(aPixgrid->size_x, aPixgrid->size_y,
603  CPL_TYPE_FLOAT);
604  cpl_imagelist_set(wcube, image, i);
605  } /* for i (all wavelength planes) */
606  } /* if weight cube */
607 
608  if (getenv("MUSE_DEBUG_WEIGHTED_GRID")) {
609  char *fn = getenv("MUSE_DEBUG_WEIGHTED_GRID");
610  FILE *grid = fopen(fn, "w");
611  if (grid) {
612  cpl_msg_info(__func__, "Writing grid to \"%s\"", fn);
613  fprintf(grid, "# i j l x y lambda\n");
614  int l;
615  for (l = 0; l < aPixgrid->size_z; l++) {
616  double lambda = (l + 1. - crpix3) * cd33 + crval3;
617  int i;
618  for (i = 0; i < aPixgrid->size_x; i++) {
619  int j;
620  for (j = 0; j < aPixgrid->size_y; j++) {
621  /* x and y position of center of current grid cell (i, j start at 0) */
622  double x, y;
623  if (wcs->iscelsph) {
624  muse_wcs_celestial_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
625  } else {
626  muse_wcs_projplane_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
627  }
628  fprintf(grid, "%03d %03d %04d %.10f %.10f %8.3f\n", i+1, j+1, l+1,
629  x, y, lambda);
630  } /* for j (y direction) */
631  } /* for i (x direction) */
632  } /* for l (wavelength planes) */
633  fclose(grid);
634  } else {
635  cpl_msg_warning(__func__, "Writing grid to \"%s\" failed!", fn);
636  }
637  } /* if grid output */
638 
639  int l;
640 #ifdef ESO_ENABLE_DEBUG
641  #pragma omp parallel for default(none) /* as req. by Ralf */ \
642  shared(aCube, aParams, aPixgrid, aPixtable, cd33, crpix3, crval3, \
643  data, debug, debugx, debugy, debugz, dq, lbda, ld, loglambda, \
644  ptxoff, ptyoff, renka_rc, stat, stdout, wcs, wcube, wght, \
645  xnorm, xout, xpos, xsz, ynorm, yout, ypos, ysz, znorm, zout, \
646  zoutefac, zsz)
647 #else
648  #pragma omp parallel for default(none) /* as req. by Ralf */ \
649  shared(aCube, aParams, aPixgrid, aPixtable, cd33, crpix3, crval3, \
650  data, dq, lbda, ld, loglambda, ptxoff, ptyoff, renka_rc, stat,\
651  stdout, wcs, wcube, wght, xnorm, xout, xpos, xsz, ynorm, yout,\
652  ypos, ysz, znorm, zout, zoutefac, zsz)
653 #endif
654  for (l = 0; l < aPixgrid->size_z; l++) {
655  float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->data, l)),
656  *pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->stat, l));
657  int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->dq, l));
658  /* wavelength of center of current grid cell (l is index starting at 0) */
659  double lambda = (l + 1. - crpix3) * cd33 + crval3;
660  double zout2 = zout; /* correct the output pixel size for log-lambda */
661  if (loglambda) {
662  lambda = crval3 * exp((l + 1. - crpix3) * cd33 / crval3);
663  zout2 = crval3 * exp((l - crpix3) * cd33 / crval3) * zoutefac;
664  }
665  float *pwcube = NULL;
666  if (wcube) { /* weight cube */
667  pwcube = cpl_image_get_data_float(cpl_imagelist_get(wcube, l));
668  } /* if weight cube */
669 
670  int i;
671  for (i = 0; i < aPixgrid->size_x; i++) {
672  int j;
673  for (j = 0; j < aPixgrid->size_y; j++) {
674  /* x and y position of center of current grid cell (i, j start at 0) */
675  double x, y;
676  if (wcs->iscelsph) {
677  muse_wcs_celestial_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
678  } else {
679  muse_wcs_projplane_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
680  }
681  double sumdata = 0, sumstat = 0, sumweight = 0;
682  int npoints = 0;
683 #ifdef ESO_ENABLE_DEBUG
684  cpl_size *pointlist = NULL;
685  double *pointweights = NULL;
686  int npointlist = 0;
687  if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
688  pointlist = cpl_calloc(100, sizeof(cpl_size));
689  pointweights = cpl_malloc(100 * sizeof(double));
690  npointlist = 100;
691  }
692 #endif
693 
694 #ifdef ESO_ENABLE_DEBUG
695  if (debug & 1) {
696  printf("cell %d %d %d:\n", i, j, l);
697  }
698 #endif
699  /* loop through surrounding cells and their contained pixels */
700  int i2;
701  for (i2 = i - ld; i2 <= i + ld; i2++) {
702  int j2;
703  for (j2 = j - ld; j2 <= j + ld; j2++) {
704  int l2;
705  for (l2 = l - ld; l2 <= l + ld; l2++) {
706  cpl_size idx2 = muse_pixgrid_get_index(aPixgrid, i2, j2, l2, CPL_FALSE),
707  n_rows2 = muse_pixgrid_get_count(aPixgrid, idx2);
708 #ifdef ESO_ENABLE_DEBUG
709  if (debug & 8 && n_rows2 < 1) {
710  printf("%d %d %d / %d %d %d (%"CPL_SIZE_FORMAT"): no rows!\n",
711  i+1, j+1, l+1, i2+1, j2+1, l2+1, idx2);
712  }
713 #endif
714  const cpl_size *rows2 = muse_pixgrid_get_rows(aPixgrid, idx2);
715  cpl_size n;
716  for (n = 0; n < n_rows2; n++) {
717  if (dq[rows2[n]]) { /* exclude all bad pixels */
718 #ifdef ESO_ENABLE_DEBUG
719  if (debug & 8) {
720  printf("%d %d %d / %d %d %d (%"CPL_SIZE_FORMAT", "
721  "%"CPL_SIZE_FORMAT"): bad!\n",
722  i+1, j+1, l+1, i2+1, j2+1, l2+1, idx2, n);
723  fflush(stdout);
724  }
725 #endif
726  continue;
727  }
728 
729  double dx = fabs(x - (xpos[rows2[n]] + ptxoff)),
730  dy = fabs(y - (ypos[rows2[n]] + ptyoff)),
731  dlambda = fabs(lambda - lbda[rows2[n]]),
732  r2 = 0;
733  if (wcs->iscelsph) {
734  /* Since the distances of RA in degrees get larger the *
735  * closer we get to the celestial pole, we have to *
736  * compensate for that by multiplying the distance in *
737  * RA by cos(delta), to make it comparable to the *
738  * distances in pixels for the differnt kernels below. */
739  dx *= cos(y * CPL_MATH_RAD_DEG);
740  }
741  if (aParams->method != MUSE_RESAMPLE_WEIGHTED_DRIZZLE) {
742  dx *= xnorm;
743  dy *= ynorm;
744  dlambda *= znorm;
745  r2 = dx*dx + dy*dy + dlambda*dlambda;
746  }
747  double weight = 0.;
748  if (aParams->method == MUSE_RESAMPLE_WEIGHTED_RENKA) {
749  weight = muse_resampling_weight_function_renka(sqrt(r2), renka_rc);
750  } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_DRIZZLE) {
751  weight = muse_resampling_weight_function_drizzle(xsz, ysz, zsz,
752  xout, yout, zout2,
753  dx, dy, dlambda);
754  } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_LINEAR) {
755  weight = muse_resampling_weight_function_linear(sqrt(r2));
756  } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_QUADRATIC) {
757  weight = muse_resampling_weight_function_quadratic(r2);
758  } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_LANCZOS) {
759  weight = muse_resampling_weight_function_lanczos(dx, dy, dlambda, ld);
760  }
761 
762  if (wght) { /* the pixel table does contain weights */
763  /* apply it on top of the weight computed here */
764  weight *= wght[rows2[n]];
765  }
766 #ifdef ESO_ENABLE_DEBUG
767  if (debug & 8) {
768  printf("%d %d %d / %d %d %d (%"CPL_SIZE_FORMAT", %"CPL_SIZE_FORMAT"):"
769  " x %e %e %e %e y %e %e %e %e l %e %e %e %e --> %e %e %e\n",
770  i+1, j+1, l+1, i2+1, j2+1, l2+1, idx2, n,
771  x, xpos[rows2[n]]+ptxoff, fabs(x - (xpos[rows2[n]]+ptxoff)), dx,
772  y, ypos[rows2[n]]+ptyoff, fabs(y - (ypos[rows2[n]]+ptyoff)), dy,
773  lambda, lbda[rows2[n]], fabs(lambda - lbda[rows2[n]]), dlambda,
774  r2, sqrt(r2), weight);
775  fflush(stdout);
776  }
777 #endif
778  sumweight += weight;
779  sumdata += data[rows2[n]] * weight;
780  sumstat += stat[rows2[n]] * weight*weight;
781  npoints++;
782 #ifdef ESO_ENABLE_DEBUG
783  if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
784  if (npoints > npointlist) {
785  pointlist = cpl_realloc(pointlist,
786  (npointlist + 100) * sizeof(cpl_size));
787  memset(pointlist + npointlist, 0, 100 * sizeof(cpl_size));
788  pointweights = cpl_realloc(pointweights,
789  (npointlist + 100) * sizeof(double));
790  npointlist += 100;
791  }
792  /* store row number instead of index, because we cannot be zero: */
793  pointlist[npoints-1] = rows2[n] + 1;
794  pointweights[npoints-1] = weight;
795  }
796 
797  if (debug & 4) {
798  cpl_size idx = muse_pixgrid_get_index(aPixgrid, i, j, l, CPL_FALSE),
799  count = muse_pixgrid_get_count(aPixgrid, idx);
800  if (count) {
801  printf(" pixel %4d,%4d,%4d (%8"CPL_SIZE_FORMAT"): "
802  "%2"CPL_SIZE_FORMAT" %2"CPL_SIZE_FORMAT" %f %f %f, "
803  " %e -> %e ==> %e %e (%d)\n", i+1, j+1, l+1, idx,
804  n, count, dx, dy, dlambda,
805  data[muse_pixgrid_get_rows(aPixgrid, idx)[n]],
806  weight, sumweight, sumdata, npoints);
807  }
808  }
809 #endif
810  } /* for n (all pixels in grid cell) */
811  } /* for l2 (lambda direction) */
812  } /* for j2 (y direction) */
813  } /* for i2 (x direction) */
814 
815 #ifdef ESO_ENABLE_DEBUG
816  if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
817  printf("cell center (%d, %d, %d): %14.7e %14.7e %9.3f\npixelnumber "
818  "weight ", debugx, debugy, debugz, x, y, lambda);
819  muse_pixtable_dump(aPixtable, 0, 0, 2);
820  int m = -1;
821  while (++m < npointlist && pointlist[m] != 0) {
822  /* access row using row index again instead of stored row number: */
823  printf("%12"CPL_SIZE_FORMAT" %8.5f ", pointlist[m] - 1,
824  pointweights[m]);
825  muse_pixtable_dump(aPixtable, pointlist[m] - 1, 1, 0);
826  }
827  printf("sums: %g %g %g --> %g %g\n", sumdata, sumstat, sumweight,
828  sumdata / sumweight, sumstat / pow(sumweight, 2));
829  cpl_free(pointlist);
830  cpl_free(pointweights);
831  }
832 
833  if (debug & 1 && npoints) {
834  printf(" sumdata: %e %e (%d)", sumdata, sumweight, npoints);
835  }
836 #endif
837 
838  /* if no points were found, we cannot divide by the summed weight *
839  * and don't need to set the output pixel value (it's 0 already), *
840  * only set the relevant Euro3D bad pixel flag */
841  if (!npoints || !isnormal(sumweight)) {
842  pdq[i + j * aPixgrid->size_x] = EURO3D_MISSDATA;
843 #ifdef ESO_ENABLE_DEBUG
844  if (debug & 1) {
845  printf(" -> no points or weird weight\n");
846  }
847 #endif
848  continue;
849  }
850 
851  /* divide results by weight of summed pixels */
852  sumdata /= sumweight;
853  sumstat /= sumweight*sumweight;
854 #ifdef ESO_ENABLE_DEBUG
855  if (debug & 1) {
856  printf(" -> %e (%e)\n", sumdata, sumstat);
857  }
858  if (debug) {
859  fflush(stdout); /* flush the output in any debug case */
860  }
861 #endif
862  pdata[i + j * aPixgrid->size_x] = sumdata;
863  pstat[i + j * aPixgrid->size_x] = sumstat;
864  pdq[i + j * aPixgrid->size_x] = EURO3D_GOODPIXEL; /* now we can mark it as good */
865  if (pwcube) {
866  pwcube[i + j * aPixgrid->size_x] = sumweight;
867  } /* if weight cube */
868  } /* for j (y direction) */
869  } /* for i (x direction) */
870  } /* for l (wavelength planes) */
871  cpl_free(wcs);
872 
873  if (wcube) { /* weight cube */
874  const char *fn = getenv("MUSE_DEBUG_WEIGHT_CUBE");
875  cpl_error_code rc = cpl_imagelist_save(wcube, fn, CPL_TYPE_UNSPECIFIED,
876  NULL, CPL_IO_CREATE);
877  if (rc != CPL_ERROR_NONE) {
878  cpl_msg_warning(__func__, "Failure to save weight cube as \"%s\": %s", fn,
879  cpl_error_get_message());
880  } else {
881  cpl_msg_info(__func__, "Saved weight cube as \"%s\"", fn);
882  }
883  cpl_imagelist_delete(wcube);
884  } /* if weight cube */
885 
886  return CPL_ERROR_NONE;
887 } /* muse_resampling_cube_weighted() */
888 
889 /*---------------------------------------------------------------------------*/
899 /*---------------------------------------------------------------------------*/
900 static cpl_error_code
901 muse_resampling_check_deltas(muse_pixtable *aPixtable,
902  muse_resampling_params *aParams)
903 {
904  cpl_ensure_code(aPixtable && aParams, CPL_ERROR_NULL_INPUT);
905  const char func[] = "muse_resampling_cube"; /* pretend to be in that fct... */
906 
907  /* wavelength direction */
908  if (aParams->dlambda == 0.0) {
909  aParams->dlambda = kMuseSpectralSamplingA;
910  if (aParams->tlambda == MUSE_RESAMPLING_DISP_AWAV_LOG) {
911  aParams->dlambda /= 1.6; /* XXX seems to be a reasonable value... */
912  }
913  }
914 
916  /* if pixel table in pixel units, take a shortcut */
917  if (aParams->dx == 0.0) {
918  aParams->dx = 1.0;
919  }
920  if (aParams->dy == 0.0) {
921  aParams->dy = 1.0;
922  }
923  cpl_msg_debug(func, "steps from parameters: %.2f pix, %.2f pix, %.3f Angstrom",
924  aParams->dx, aParams->dy, aParams->dlambda);
925  return CPL_ERROR_NONE;
926  }
927  if (aParams->dx == 0.0) {
928  aParams->dx = kMuseSpaxelSizeX_WFM / 3600.;
929  if (muse_pfits_get_mode(aPixtable->header) == MUSE_MODE_NFM_AO_N) {
930  aParams->dx = kMuseSpaxelSizeX_NFM / 3600.;
931  }
932  } else {
933  /* convert from arcsec to degrees */
934  aParams->dx /= 3600.;
935  }
936  if (aParams->dy == 0.0) {
937  aParams->dy = kMuseSpaxelSizeY_WFM / 3600.;
938  if (muse_pfits_get_mode(aPixtable->header) == MUSE_MODE_NFM_AO_N) {
939  aParams->dy = kMuseSpaxelSizeY_NFM / 3600.;
940  }
941  } else {
942  /* anything else should be interpreted as arcsec, but we need deg */
943  aParams->dy /= 3600.;
944  }
945  cpl_msg_debug(func, "steps from parameters: %f arcsec, %f arcsec, %.3f Angstrom",
946  aParams->dx * 3600., aParams->dy * 3600., aParams->dlambda);
947  return CPL_ERROR_NONE;
948 } /* muse_resampling_check_deltas() */
949 
950 /*---------------------------------------------------------------------------*/
964 /*---------------------------------------------------------------------------*/
965 static cpl_error_code
966 muse_resampling_compute_size(muse_pixtable *aPixtable,
967  muse_resampling_params *aParams,
968  int *aX, int *aY, int *aZ)
969 {
970  cpl_ensure_code(aPixtable && aParams && aX && aY && aZ, CPL_ERROR_NULL_INPUT);
971  const char func[] = "muse_resampling_cube"; /* pretend to be in that fct... */
972  double x1, y1, x2, y2;
973  float xmin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XLO),
974  xmax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XHI),
975  ymin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YLO),
976  ymax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YHI);
977  muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
978  if (wcstype == MUSE_PIXTABLE_WCS_CELSPH) {
979  muse_wcs_projplane_from_celestial(aPixtable->header, xmin, ymin, &x1, &y1);
980  muse_wcs_projplane_from_celestial(aPixtable->header, xmax, ymax, &x2, &y2);
981  } else {
982  muse_wcs_projplane_from_pixel(aPixtable->header, xmin, ymin, &x1, &y1);
983  muse_wcs_projplane_from_pixel(aPixtable->header, xmax, ymax, &x2, &y2);
984  }
985  *aX = lround(fabs(x2 - x1) / aParams->dx) + 1;
986  *aY = lround(fabs(y2 - y1) / aParams->dy) + 1;
987  float lmin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LLO),
988  lmax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LHI);
989  *aZ = (int)ceil((lmax - lmin) / aParams->dlambda) + 1;
990  if (aParams->tlambda == MUSE_RESAMPLING_DISP_AWAV_LOG) {
991  *aZ = (int)ceil(lmin / aParams->dlambda * log(lmax / lmin)) + 1;
992  }
993  cpl_msg_info(func, "Output cube size %d x %d x %d (fit to data)",
994  *aX, *aY, *aZ);
995  return CPL_ERROR_NONE;
996 } /* muse_resampling_compute_size() */
997 
998 /*---------------------------------------------------------------------------*/
1008 /*---------------------------------------------------------------------------*/
1009 static void
1010 muse_resampling_override_size_int(int *aV, const char *aKeyword, int aValue)
1011 {
1012  if (aValue <= 0 || !aV) {
1013  return;
1014  }
1015  const char func[] = "muse_resampling_cube"; /* pretend to be in that fct... */
1016  cpl_msg_info(func, "Overriding %s=%d (was %d)", aKeyword, aValue, *aV);
1017  *aV = aValue;
1018 } /* muse_resampling_override_size_int() */
1019 
1020 /*---------------------------------------------------------------------------*/
1034 /*---------------------------------------------------------------------------*/
1035 static cpl_error_code
1036 muse_resampling_override_size(int *aX, int *aY, int *aZ,
1037  muse_resampling_params *aParams)
1038 {
1039  cpl_ensure_code(aX && aY && aZ && aParams, CPL_ERROR_NULL_INPUT);
1040  if (!aParams->wcs) { /* quietly return */
1041  return CPL_ERROR_NONE;
1042  }
1043  const char func[] = "muse_resampling_cube"; /* pretend to be in that fct... */
1044  /* cube size overrides */
1045  const cpl_array *dims = cpl_wcs_get_image_dims(aParams->wcs);
1046  if (!dims) {
1047  cpl_msg_debug(func, "No dimensions to override were specified");
1048  return CPL_ERROR_NONE;
1049  }
1050  muse_resampling_override_size_int(aX, "NAXIS1", cpl_array_get_int(dims, 0, NULL));
1051  muse_resampling_override_size_int(aY, "NAXIS2", cpl_array_get_int(dims, 1, NULL));
1052  if (cpl_wcs_get_image_naxis(aParams->wcs) >= 3) {
1053  muse_resampling_override_size_int(aZ, "NAXIS3",
1054  cpl_array_get_int(dims, 2, NULL));
1055  }
1056  return CPL_ERROR_NONE;
1057 } /* muse_resampling_override_size() */
1058 
1059 /*---------------------------------------------------------------------------*/
1070 /*---------------------------------------------------------------------------*/
1071 static void
1072 muse_resampling_override_wcs_double(muse_datacube *aCube, const char *aKeyword,
1073  double aValue, int aError)
1074 {
1075  if (aError || !aCube) {
1076  cpl_msg_debug("double", "%s=%#g (%d)", aKeyword, aValue, aError);
1077  return;
1078  }
1079  const char func[] = "muse_resampling_cube"; /* pretend to be in that fct... */
1080  double old = cpl_propertylist_has(aCube->header, aKeyword)
1081  ? cpl_propertylist_get_double(aCube->header, aKeyword) : NAN;
1082  cpl_msg_info(func, "Overriding %s=%#g (was %#g)", aKeyword, aValue, old);
1083  cpl_propertylist_update_double(aCube->header, aKeyword, aValue);
1084  /* Leave the marked that something was overridden, will *
1085  * be evaluated and removed in muse_pixgrid_create()! */
1086  cpl_propertylist_update_bool(aCube->header, "MUSE_RESAMPLING_WCS_OVERRIDDEN",
1087  CPL_TRUE);
1088 } /* muse_resampling_override_wcs_double() */
1089 
1090 /*---------------------------------------------------------------------------*/
1102 /*---------------------------------------------------------------------------*/
1103 static cpl_error_code
1104 muse_resampling_override_wcs(muse_datacube *aCube,
1105  muse_resampling_params *aParams)
1106 {
1107  cpl_ensure_code(aCube && aCube->header && aParams, CPL_ERROR_NULL_INPUT);
1108  if (!aParams->wcs) { /* quietly return */
1109  return CPL_ERROR_NONE;
1110  }
1111  const char func[] = "muse_resampling_cube"; /* pretend to be in that fct... */
1112  const cpl_array *crval = cpl_wcs_get_crval(aParams->wcs),
1113  *crpix = cpl_wcs_get_crpix(aParams->wcs);
1114  const cpl_matrix *cd = cpl_wcs_get_cd(aParams->wcs);
1115  int err = 0;
1116  /* spatial axes overrides */
1117  if (crval) {
1118  muse_resampling_override_wcs_double(aCube, "CRVAL1", cpl_array_get_double(crval, 0, &err), err);
1119  muse_resampling_override_wcs_double(aCube, "CRVAL2", cpl_array_get_double(crval, 1, &err), err);
1120  } else {
1121  cpl_msg_debug(func, "No CRVALj to override were specified");
1122  }
1123  if (crpix) {
1124  muse_resampling_override_wcs_double(aCube, "CRPIX1", cpl_array_get_double(crpix, 0, &err), err);
1125  muse_resampling_override_wcs_double(aCube, "CRPIX2", cpl_array_get_double(crpix, 1, &err), err);
1126  } else {
1127  cpl_msg_debug(func, "No CRPIXi to override were specified");
1128  }
1129  if (cd) {
1130  int naxes = cpl_matrix_get_ncol(cd); /* nrow should be the same */
1131  if (cpl_matrix_get_determinant(cd) == 0.) {
1132  cpl_msg_warning(func, "det(CDi_j) = 0, not overriding CDi_j!");
1133  cd = NULL; /* don't override dispersion direction, either */
1134  } else if (naxes > 2 && (cpl_matrix_get(cd, 0, 2) != 0. ||
1135  cpl_matrix_get(cd, 2, 0) != 0. ||
1136  cpl_matrix_get(cd, 2, 1) != 0. ||
1137  cpl_matrix_get(cd, 1, 2) != 0.)) {
1138  cpl_msg_warning(func, "Axis 3 (dispersion) is not cleanly separated from "
1139  "axes 1 and 2, not overriding CDi_j!");
1140  cd = NULL; /* don't override dispersion direction, either */
1141  } else {
1142  cpl_errorstate es = cpl_errorstate_get();
1143  muse_resampling_override_wcs_double(aCube, "CD1_1", cpl_matrix_get(cd, 0, 0),
1144  !cpl_errorstate_is_equal(es));
1145  es = cpl_errorstate_get();
1146  muse_resampling_override_wcs_double(aCube, "CD1_2", cpl_matrix_get(cd, 0, 1),
1147  !cpl_errorstate_is_equal(es));
1148  es = cpl_errorstate_get();
1149  muse_resampling_override_wcs_double(aCube, "CD2_1", cpl_matrix_get(cd, 1, 0),
1150  !cpl_errorstate_is_equal(es));
1151  es = cpl_errorstate_get();
1152  muse_resampling_override_wcs_double(aCube, "CD2_2", cpl_matrix_get(cd, 1, 1),
1153  !cpl_errorstate_is_equal(es));
1154  }
1155  } else {
1156  cpl_msg_debug(func, "No CDi_j to override were specified");
1157  }
1158  /* wavelength axis overrides; apparently, values originally in Angstrom *
1159  * values are here in meters (no way to check this here, since we cannot *
1160  * access aParams->wcs->wcsptr!), so need to convert them back */
1161  if (cpl_array_get_size(crval) > 2) {
1162  if (crval) {
1163  muse_resampling_override_wcs_double(aCube, "CRVAL3",
1164  cpl_array_get_double(crval, 2, &err) * 1e10, err);
1165  }
1166  if (crpix) {
1167  muse_resampling_override_wcs_double(aCube, "CRPIX3", cpl_array_get_double(crpix, 2, &err), err);
1168  }
1169  if (cd) {
1170  cpl_errorstate es = cpl_errorstate_get();
1171  muse_resampling_override_wcs_double(aCube, "CD3_3", cpl_matrix_get(cd, 2, 2) * 1e10,
1172  !cpl_errorstate_is_equal(es));
1173  }
1174  } /* if 3D */
1175  return CPL_ERROR_NONE;
1176 } /* muse_resampling_override_wcs() */
1177 
1178 /*---------------------------------------------------------------------------*/
1188 /*---------------------------------------------------------------------------*/
1189 static cpl_error_code
1190 muse_resampling_fit_data_to_grid(muse_datacube *aCube, muse_pixtable *aPixtable)
1191 {
1192  cpl_ensure_code(aCube && aPixtable, CPL_ERROR_NULL_INPUT);
1193  const char func[] = "muse_resampling_cube"; /* pretend to be in that fct... */
1194  /* Check if we need to move the spatial zeropoint to fit into the *
1195  * computed grid at all. This is only the case, if the grid was *
1196  * computed automatically, not overridden. */
1197  if (cpl_propertylist_has(aCube->header, "MUSE_RESAMPLING_WCS_OVERRIDDEN") &&
1198  cpl_propertylist_get_bool(aCube->header, "MUSE_RESAMPLING_WCS_OVERRIDDEN")) {
1199  cpl_msg_warning(func, "Output WCS was forced, won't update CRPIX[12]!");
1200  cpl_propertylist_erase(aCube->header, "MUSE_RESAMPLING_WCS_OVERRIDDEN");
1201  return CPL_ERROR_NONE;
1202  }
1203  /* do determine offset of coordinate to the first (1,1) pixel */
1204  double xoff1, yoff1, xoff2, yoff2;
1205  float xlo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XLO),
1206  xhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XHI),
1207  ylo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YLO),
1208  yhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YHI);
1209  muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
1210  if (wcstype == MUSE_PIXTABLE_WCS_CELSPH) {
1211  muse_wcs_pixel_from_celestial(aCube->header, xlo, ylo, &xoff1, &yoff1);
1212  muse_wcs_pixel_from_celestial(aCube->header, xhi, yhi, &xoff2, &yoff2);
1213  } else {
1214  muse_wcs_pixel_from_projplane(aCube->header, xlo, ylo, &xoff1, &yoff1);
1215  muse_wcs_pixel_from_projplane(aCube->header, xhi, yhi, &xoff2, &yoff2);
1216  }
1217  /* the actual minimum offset we need depends on the axis-direction, *
1218  * i.e. the sign of CDi_j, so we take the minimum of the result here */
1219  double xoff = fmin(xoff1, xoff2) - 1,
1220  yoff = fmin(yoff1, yoff2) - 1;
1221  if (xoff != 0. || yoff != 0.) {
1222  double crpix1 = cpl_propertylist_get_double(aCube->header, "CRPIX1"),
1223  crpix2 = cpl_propertylist_get_double(aCube->header, "CRPIX2");
1224  cpl_msg_debug(func, "Updating CRPIX[12]=%#g,%#g (were %#g,%#g)",
1225  crpix1 - xoff, crpix2 - yoff, crpix1, crpix2);
1226  crpix1 -= xoff;
1227  crpix2 -= yoff;
1228  cpl_propertylist_update_double(aCube->header, "CRPIX1", crpix1);
1229  cpl_propertylist_update_double(aCube->header, "CRPIX2", crpix2);
1230  } /* if offsets */
1231  return CPL_ERROR_NONE;
1232 } /* muse_resampling_fit_data_to_grid() */
1233 
1234 /* corresponds to muse_resampling_crstats_type, *
1235  * keep in sync with those values! */
1236 const char *crtypestring[] = {
1237  "IRAF-like",
1238  "average",
1239  "median"
1240 };
1241 
1242 #define MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_START \
1243  cpl_size idx = muse_pixgrid_get_index(aPixgrid, i2, j2, l2, CPL_FALSE),\
1244  nrow = muse_pixgrid_get_count(aPixgrid, idx), \
1245  irow; \
1246  const cpl_size *rows = muse_pixgrid_get_rows(aPixgrid, idx); \
1247  for (irow = 0; irow < nrow; irow++) { \
1248  if (!dq[rows[irow]] || dq[rows[irow]] & badinclude) { \
1249  /* do nothing */ \
1250  } else if (dq[rows[irow]]) { \
1251  continue; /* exclude all other bad pixels */ \
1252  }
1253 #define MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_DEBUG \
1254  if (debug & 4 && i+1 == debugx && j+1 == debugy && l+1 == debugz) { \
1255  printf("%s: data / stat (%"CPL_SIZE_FORMAT") = %e / %e\n", \
1256  __func__, npixels+1, data[rows[irow]], stat[rows[irow]]); \
1257  }
1258 #define MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END \
1259  if (aType == MUSE_RESAMPLING_CRSTATS_IRAF) { \
1260  /* add the point from the pixel table to the vector */ \
1261  level += data[rows[irow]]; \
1262  dev += stat[rows[irow]]; \
1263  } else { \
1264  /* add the point to the npixels'th column of the image */ \
1265  if (npixels >= sxsize) { \
1266  /* is there no cpl_image_resize()?! */ \
1267  sxsize += MUSE_CRREJECT_MAX_NPIX; \
1268  cpl_image *tmp = cpl_image_new(sxsize, 2, CPL_TYPE_DOUBLE); \
1269  cpl_image_copy(tmp, sdata, 1, 1); \
1270  cpl_image_delete(sdata); \
1271  sdata = tmp; \
1272  vdata = cpl_image_get_data_double(sdata); \
1273  } \
1274  vdata[npixels] = data[rows[irow]]; \
1275  vdata[npixels + sxsize] = stat[rows[irow]]; \
1276  } \
1277  /* increase the index for all methods */ \
1278  npixels++; \
1279  } /* for irow (all pixels in grid cell) */
1280 
1281 /*---------------------------------------------------------------------------*/
1304 /*---------------------------------------------------------------------------*/
1305 static void
1306 muse_resampling_crreject(muse_pixtable *aPixtable, muse_pixgrid *aPixgrid,
1307  muse_resampling_crstats_type aType, double aSigma)
1308 {
1309  const char *id = "muse_resampling_cube"; /* pretend to be in that function */
1310  if (aType > MUSE_RESAMPLING_CRSTATS_MEDIAN) {
1311  cpl_msg_warning(id, "Unknown type (%u) for cosmic-ray rejection statistics,"
1312  " resetting to \"%s\"", aType,
1313  crtypestring[MUSE_RESAMPLING_CRSTATS_IRAF]);
1315  }
1316  cpl_msg_info(id, "Using %s statistics (%.3f sigma level) for cosmic-ray"
1317  " rejection", crtypestring[aType], aSigma);
1318 
1319 #ifdef ESO_ENABLE_DEBUG
1320  int debug = 0, debugx = 0, debugy = 0, debugz = 0;
1321  if (getenv("MUSE_DEBUG_CRREJECT")) {
1322  debug = atoi(getenv("MUSE_DEBUG_CRREJECT"));
1323  }
1324  if (debug & 2) { /* need coordinates */
1325  if (getenv("MUSE_DEBUG_CRREJECT_X")) {
1326  debugx = atoi(getenv("MUSE_DEBUG_CRREJECT_X"));
1327  if (debugx < 1 || debugx > aPixgrid->size_x) {
1328  debugx = 0;
1329  }
1330  }
1331  if (getenv("MUSE_DEBUG_CRREJECT_Y")) {
1332  debugy = atoi(getenv("MUSE_DEBUG_CRREJECT_Y"));
1333  if (debugy < 1 || debugy > aPixgrid->size_y) {
1334  debugy = 0;
1335  }
1336  }
1337  if (getenv("MUSE_DEBUG_CRREJECT_Z")) {
1338  debugz = atoi(getenv("MUSE_DEBUG_CRREJECT_Z"));
1339  if (debugz < 1 || debugz > aPixgrid->size_z) {
1340  debugz = 0;
1341  }
1342  }
1343  }
1344 #endif
1345 
1346  float *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
1347  *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
1348  int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
1349 
1350  enum dirtype {
1351  DIR_X = 0,
1352  DIR_Y = 1,
1353  DIR_NONE = 2
1354  };
1355  enum dirtype dir = DIR_NONE;
1356  /* XXX also exclude multiple exposures with different POSANG! */
1357  cpl_boolean haswcs = muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_CELSPH;
1358  double posang = muse_astro_posangle(aPixtable->header);
1359  const double palimit = 5.;
1360  if (!haswcs || (fabs(posang) < palimit || fabs(fabs(posang) - 180.) < palimit ||
1361  fabs(fabs(posang) - 360.) < palimit)) {
1362  cpl_msg_debug(id, "CR rejection: posang = %f deg --> DIR_X "
1363  "(%s / %s / %s / %s)", posang, haswcs ? "yes": "no",
1364  fabs(posang) < palimit ? "true" : "false",
1365  fabs(fabs(posang) - 180.) < palimit ? "true" : "false",
1366  fabs(fabs(posang) - 360.) < palimit ? "true" : "false");
1367  dir = DIR_X;
1368  } else if (fabs(fabs(posang) - 90.) < palimit ||
1369  fabs(fabs(posang) - 270.) < palimit) {
1370  cpl_msg_debug(id, "CR rejection: posang = %f deg --> DIR_Y "
1371  "(%s / %s / %s)", posang, haswcs ? "yes": "no",
1372  fabs(fabs(posang) - 90.) < palimit ? "true" : "false",
1373  fabs(fabs(posang) - 270.) < palimit ? "true" : "false");
1374  dir = DIR_Y;
1375  } else {
1376  cpl_msg_debug(id, "CR rejection: posang = %f deg --> DIR_NONE "
1377  "(%s / %s / %s / %s / %s / %s)", posang, haswcs ? "yes": "no",
1378  fabs(posang) < palimit ? "true" : "false",
1379  fabs(fabs(posang) - 90.) < palimit ? "true" : "false",
1380  fabs(fabs(posang) - 180.) < palimit ? "true" : "false",
1381  fabs(fabs(posang) - 270.) < palimit ? "true" : "false",
1382  fabs(fabs(posang) - 360.) < palimit ? "true" : "false");
1383  }
1384 
1385  /* when computing surrounding statistics, include these types of bad pixels */
1386  unsigned badinclude = EURO3D_COSMICRAY;
1387  int l;
1388 #ifdef ESO_ENABLE_DEBUG
1389  #pragma omp parallel for default(none) /* as req. by Ralf */ \
1390  shared(aPixgrid, aPixtable, aSigma, badinclude, aType, crtypestring, \
1391  data, debug, debugx, debugy, debugz, dir, dq, id, stat, stdout)
1392 #else
1393  #pragma omp parallel for default(none) /* as req. by Ralf */ \
1394  shared(aPixgrid, aSigma, aType, badinclude, crtypestring, data, dir, \
1395  dq, id, stat)
1396 #endif
1397  for (l = 0; l < aPixgrid->size_z; l++) {
1398 #define MUSE_CRREJECT_MAX_NPIX 1000 /* XXX compute size dynamically somehow */
1399  int sxsize = MUSE_CRREJECT_MAX_NPIX;
1400  cpl_image *sdata = NULL;
1401  double *vdata = NULL;
1402  if (aType > MUSE_RESAMPLING_CRSTATS_IRAF) {
1403  /* create image for statistics computations */
1404  /* XXX would memset'ting the buffer not be faster? */
1405  sdata = cpl_image_new(sxsize, 2, CPL_TYPE_DOUBLE);
1406  vdata = cpl_image_get_data_double(sdata);
1407  }
1408 
1409  int i;
1410  for (i = 0; i < aPixgrid->size_x; i++) {
1411  int j;
1412  for (j = 0; j < aPixgrid->size_y; j++) {
1413  /* in the IRAF-like case, compute local statistics, i.e. set *
1414  * the variables to zero (mean, and stdev from the variance) */
1415  double level = 0., dev = 0;
1416  /* loop through surrounding cells: compute local statistics */
1417 #define CR_LD 1 /* loop distance for CR rejection statistics */
1418  cpl_size npixels = 0; /* count pixels again, without bad ones */
1419  int i2, j2, l2;
1420  if (dir == DIR_Y) {
1421  for (i2 = i - CR_LD; i2 <= i + CR_LD; i2++) {
1422  int nwidth = CR_LD;
1423  if (i2 == i) {
1424  nwidth = 0; /* don't use neighbors in the same slice */
1425  }
1426  for (j2 = j - nwidth; j2 <= j + nwidth; j2++) {
1427  for (l2 = l - nwidth; l2 <= l + nwidth; l2++) {
1428  MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_START
1429 #ifdef ESO_ENABLE_DEBUG
1430  MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_DEBUG
1431 #endif
1432  MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END
1433  } /* for l2 (wavelength/z direction) */
1434  } /* for i2 (x direction) */
1435  } /* for j2 (y direction) */
1436  } else {
1437  for (j2 = j - CR_LD; j2 <= j + CR_LD; j2++) {
1438  int nwidth = CR_LD;
1439  if (dir == DIR_X && j2 == j) {
1440  nwidth = 0; /* don't use neighbors in the same slice */
1441  }
1442  for (i2 = i - nwidth; i2 <= i + nwidth; i2++) {
1443  for (l2 = l - nwidth; l2 <= l + nwidth; l2++) {
1444  MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_START
1445 #ifdef ESO_ENABLE_DEBUG
1446  MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_DEBUG
1447 #endif
1448  MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END
1449  } /* for l2 (wavelength/z direction) */
1450  } /* for i2 (x direction) */
1451  } /* for j2 (y direction) */
1452  } /* else */
1453  /* need at least 3 values to do something sensible */
1454  if (npixels < 3) {
1455  continue;
1456  }
1457  if (aType == MUSE_RESAMPLING_CRSTATS_IRAF) {
1458  level /= npixels; /* local mean value */
1459  dev = sqrt(dev) / npixels; /* local mean sigma value */
1460  } else {
1461  unsigned sflags;
1462  if (aType == MUSE_RESAMPLING_CRSTATS_MEDIAN) {
1463  sflags = CPL_STATS_MEDIAN | CPL_STATS_MEDIAN_DEV;
1464  } else {
1465  sflags = CPL_STATS_MEAN | CPL_STATS_STDEV;
1466  }
1467  cpl_stats *sstats = cpl_stats_new_from_image_window(sdata, sflags,
1468  1, 1, npixels, 1);
1469  if (aType == MUSE_RESAMPLING_CRSTATS_MEDIAN) {
1470  level = cpl_stats_get_median(sstats); /* local median value */
1471  dev = cpl_stats_get_median_dev(sstats); /* local median deviation */
1472  } else {
1473  level = cpl_stats_get_mean(sstats); /* local mean value */
1474  dev = cpl_stats_get_stdev(sstats); /* local standard deviation */
1475  }
1476  cpl_stats_delete(sstats);
1477  }
1478  double limit = level + aSigma * dev;
1479 #ifdef ESO_ENABLE_DEBUG
1480  if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
1481  if (aType == MUSE_RESAMPLING_CRSTATS_IRAF) {
1482  printf("%s: %03d,%03d,%04d: %.3f+/-%.3f (stats), %.3f (limit) (%"
1483  CPL_SIZE_FORMAT" values)\n", __func__, i+1, j+1, l+1, level,
1484  dev, limit, npixels);
1485  } else {
1486  cpl_stats *ssdata = cpl_stats_new_from_image_window(sdata, CPL_STATS_ALL,
1487  1, 1, npixels, 1),
1488  *ssstat = cpl_stats_new_from_image_window(sdata, CPL_STATS_ALL,
1489  1, 2, npixels, 2);
1490  printf("%s: %03d,%03d,%04d: %e +/- %e (%s), %e (limit) (%"CPL_SIZE_FORMAT
1491  " values); data stats:\n", __func__, i+1, j+1, l+1,
1492  level, dev, crtypestring[aType], limit, npixels);
1493  cpl_stats_dump(ssdata, CPL_STATS_ALL, stdout);
1494  printf("%s: variance stats:\n", __func__);
1495  cpl_stats_dump(ssstat, CPL_STATS_ALL, stdout);
1496  fflush(stdout);
1497  cpl_stats_delete(ssdata);
1498  cpl_stats_delete(ssstat);
1499  }
1500  }
1501 #endif
1502 
1503  /* loop through pixels in the central cell: mark cosmic rays */
1504  cpl_size idx = muse_pixgrid_get_index(aPixgrid, i, j, l, CPL_FALSE),
1505  nrow = muse_pixgrid_get_count(aPixgrid, idx),
1506  irow;
1507  const cpl_size *rows = muse_pixgrid_get_rows(aPixgrid, idx);
1508  for (irow = 0; irow < nrow; irow++) {
1509  if (data[rows[irow]] > limit) {
1510  dq[rows[irow]] |= EURO3D_COSMICRAY;
1511 #ifdef ESO_ENABLE_DEBUG
1512  if (debug & 1 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
1513  printf("%s: %03d,%03d,%04d: rejected row %"CPL_SIZE_FORMAT" (%"
1514  CPL_SIZE_FORMAT" of %"CPL_SIZE_FORMAT" in this gridcell):\t",
1515  __func__, i+1, j+1, l+1, rows[irow], irow+1, nrow);
1516  muse_pixtable_dump(aPixtable, rows[irow], 1, 0);
1517  }
1518 #endif
1519  }
1520  } /* for irow (all pixels in grid cell) */
1521 #ifdef ESO_ENABLE_DEBUG
1522  if (debug) {
1523  fflush(stdout);
1524  }
1525 #endif
1526  } /* for j (y direction) */
1527  } /* for i (x direction) */
1528  cpl_image_delete(sdata);
1529  } /* for l (wavelength planes) */
1530 } /* muse_resampling_crreject() */
1531 
1532 /*---------------------------------------------------------------------------*/
1555 /*---------------------------------------------------------------------------*/
1558  muse_resampling_params *aParams)
1559 {
1560  cpl_ensure(aParams, CPL_ERROR_NULL_INPUT, NULL);
1561  cpl_ensure(aParams->method < MUSE_RESAMPLE_NONE, CPL_ERROR_ILLEGAL_INPUT, NULL);
1562  /* all other error checking is already done in muse_resampling_cube() */
1563 
1564  /* call the existing resampling routine for datacubes */
1565  muse_datacube *cube = muse_resampling_cube(aPixtable, aParams, NULL);
1566  if (!cube) {
1567  return NULL;
1568  }
1569 
1570  /* create Euro3D object, copy header and set typical Euro3D primary headers */
1571  muse_euro3dcube *e3d = cpl_calloc(1, sizeof(muse_euro3dcube));
1572  e3d->header = cpl_propertylist_duplicate(cube->header);
1573  cpl_propertylist_erase_regexp(e3d->header,
1574  "^SIMPLE$|^BITPIX$|^NAXIS|^EURO3D$|^E3D_",
1575  0);
1576  cpl_propertylist_append_char(e3d->header, "EURO3D", 'T');
1577  cpl_propertylist_set_comment(e3d->header, "EURO3D",
1578  "file conforms to Euro3D standard");
1579  /* E3D_VERS is supposed to be a string! */
1580  cpl_propertylist_append_string(e3d->header, "E3D_VERS", "1.0");
1581  cpl_propertylist_set_comment(e3d->header, "E3D_VERS",
1582  "version number of the Euro3D format");
1583  cpl_errorstate prestate = cpl_errorstate_get();
1584  double darlambdaref = cpl_propertylist_get_double(e3d->header,
1586  if (!cpl_errorstate_is_equal(prestate)) {
1587  darlambdaref = -1.; /* set to negative value to be sure on error */
1588  cpl_errorstate_set(prestate);
1589  }
1590  if (darlambdaref > 0.) {
1591  cpl_propertylist_append_char(e3d->header, "E3D_ADC", 'T');
1592  cpl_propertylist_set_comment(e3d->header, "E3D_ADC",
1593  "data was corrected for atmospheric dispersion");
1594  } else {
1595  cpl_propertylist_append_char(e3d->header, "E3D_ADC", 'F');
1596  cpl_propertylist_set_comment(e3d->header, "E3D_ADC",
1597  "data not corrected for atmospheric dispersion");
1598  }
1599 
1600  /* fill WCS info in data header */
1601  e3d->hdata = cpl_propertylist_new();
1602  cpl_propertylist_append_string(e3d->hdata, "EXTNAME", "E3D_DATA");
1603  cpl_propertylist_set_comment(e3d->hdata, "EXTNAME",
1604  "This is the Euro3D data table extension");
1605  cpl_propertylist_append_string(e3d->hdata, "CTYPES",
1606  cpl_propertylist_get_string(e3d->header, "CTYPE3"));
1607  cpl_propertylist_set_comment(e3d->hdata, "CTYPES",
1608  cpl_propertylist_get_comment(e3d->header, "CTYPE3"));
1609  cpl_propertylist_append_string(e3d->hdata, "CUNITS",
1610  cpl_propertylist_get_string(e3d->header, "CUNIT3"));
1611  cpl_propertylist_set_comment(e3d->hdata, "CUNITS",
1612  cpl_propertylist_get_comment(e3d->header, "CUNIT3"));
1613  cpl_propertylist_append_double(e3d->hdata, "CRVALS",
1614  cpl_propertylist_get_double(e3d->header, "CRVAL3"));
1615  cpl_propertylist_set_comment(e3d->hdata, "CRVALS",
1616  cpl_propertylist_get_comment(e3d->header, "CRVAL3"));
1617  cpl_propertylist_append_double(e3d->hdata, "CDELTS",
1618  cpl_propertylist_get_double(e3d->header, "CD3_3"));
1619  cpl_propertylist_set_comment(e3d->hdata, "CDELTS",
1620  cpl_propertylist_get_comment(e3d->header, "CD3_3"));
1621  /* remove WCS for spatial axis from main header */
1622  cpl_propertylist_erase_regexp(e3d->header, "^C...*3$|^CD3_.$", 0);
1623 
1624  /* create and fill Euro3D table structure */
1625  int nx = cpl_image_get_size_x(cpl_imagelist_get(cube->data, 0)),
1626  ny = cpl_image_get_size_y(cpl_imagelist_get(cube->data, 0)),
1627  nz = cpl_imagelist_get_size(cube->data);
1629  /* set array columns to initial (expected) depth */
1630  /* XXX this is very expensive for big tables, we should make *
1631  * muse_cpltable_new() take an extra argument for the depth */
1632  cpl_table_set_column_depth(e3d->dtable, "DATA_SPE", nz);
1633  cpl_table_set_column_depth(e3d->dtable, "QUAL_SPE", nz);
1634  cpl_table_set_column_depth(e3d->dtable, "STAT_SPE", nz);
1635  /* set column units depending on data unit type in the pixel table */
1636  cpl_table_set_column_unit(e3d->dtable, "DATA_SPE",
1637  cpl_table_get_column_unit(aPixtable->table,
1638  MUSE_PIXTABLE_DATA));
1639  cpl_table_set_column_unit(e3d->dtable, "STAT_SPE",
1640  cpl_table_get_column_unit(aPixtable->table,
1641  MUSE_PIXTABLE_STAT));
1642  /* set column save types as needed */
1643  cpl_table_set_column_savetype(e3d->dtable, "SELECTED", CPL_TYPE_BOOL);
1644 
1645  muse_wcs *wcs = muse_wcs_new(cube->header);
1647  if (wcs->iscelsph) { /* if not, then the preset unit of "pix" is valid */
1648  cpl_table_set_column_unit(e3d->dtable, "XPOS", "deg");
1649  cpl_table_set_column_unit(e3d->dtable, "YPOS", "deg");
1650  }
1651 
1652  /* Loop over all spaxels in the cube and transfer the data */
1653  unsigned euro3d_ignore = EURO3D_OUTSDRANGE | EURO3D_MISSDATA;
1654  cpl_vector *vusable = cpl_vector_new(nx * ny); /* for statistics */
1655  int i, itable = 0;
1656  for (i = 1; i <= nx; i++) {
1657  int j;
1658  for (j = 1; j <= ny; j++) {
1659  cpl_array *adata = cpl_array_new(nz, CPL_TYPE_FLOAT),
1660  *adq = cpl_array_new(nz, CPL_TYPE_INT),
1661  *astat = cpl_array_new(nz, CPL_TYPE_FLOAT);
1662  /* WCS coordinates */
1663  double x, y;
1664  if (wcs->iscelsph) {
1665  muse_wcs_celestial_from_pixel_fast(wcs, i, j, &x, &y);
1666  } else {
1667  muse_wcs_projplane_from_pixel_fast(wcs, i, j, &x, &y);
1668  }
1669 
1670  int l, nusable = 0, nstart = -1; /* numbers needed for Euro3D format */
1671  for (l = 0; l < nz; l++) {
1672  int err;
1673  unsigned dq = cpl_image_get(cpl_imagelist_get(cube->dq, l), i, j, &err);
1674  /* check if we can ignore this as bad pixel at the start */
1675  if (nstart < 0 && (dq & euro3d_ignore)) {
1676  continue;
1677  }
1678  cpl_array_set_int(adq, nusable, dq);
1679  cpl_array_set_float(adata, nusable,
1680  cpl_image_get(cpl_imagelist_get(cube->data, l),
1681  i, j, &err));
1682  /* take the square root so that we get errors and not variances */
1683  cpl_array_set_float(astat, nusable,
1684  sqrt(cpl_image_get(cpl_imagelist_get(cube->stat, l),
1685  i, j, &err)));
1686  nusable++;
1687  if (nstart < 0) {
1688  nstart = l + 1;
1689  }
1690  } /* for l (z / wavelengths) */
1691 
1692  cpl_table_set_int(e3d->dtable, "SPEC_ID", itable, itable + 1);
1693  cpl_table_set_int(e3d->dtable, "SPEC_LEN", itable, nusable);
1694  cpl_table_set_int(e3d->dtable, "SPEC_STA", itable, nstart);
1695  cpl_table_set_double(e3d->dtable, "XPOS", itable, x);
1696  cpl_table_set_double(e3d->dtable, "YPOS", itable, y);
1697  /* not sure what to use for SPAX_ID, just leave it unset for now */
1698  /*cpl_table_set_string(e3d->dtable, "SPAX_ID", itable, ""); */
1699  cpl_table_set_array(e3d->dtable, "DATA_SPE", itable, adata);
1700  cpl_table_set_array(e3d->dtable, "QUAL_SPE", itable, adq);
1701  cpl_table_set_array(e3d->dtable, "STAT_SPE", itable, astat);
1702 
1703  cpl_array_delete(adata);
1704  cpl_array_delete(adq);
1705  cpl_array_delete(astat);
1706 
1707  cpl_vector_set(vusable, itable, nusable);
1708  if (nstart != -1 && nusable > 0) {
1709  /* good spectrum, unselect to not remove it below */
1710  cpl_table_unselect_row(e3d->dtable, itable);
1711  }
1712  itable++; /* move to next table row */
1713  } /* for j (y spaxels) */
1714  } /* for i (x spaxels) */
1715  cpl_free(wcs);
1716 
1717  cpl_vector_set_size(vusable, itable);
1718  int maxusable = cpl_vector_get_max(vusable);
1719 #if 0
1720  cpl_msg_debug(__func__, "filled %d of %d spaxels, usable are "
1721  "%f..%f(%f)+/-%f..%d pixels per spectrum", itable, nx * ny,
1722  cpl_vector_get_min(vusable), cpl_vector_get_mean(vusable),
1723  cpl_vector_get_median(vusable), cpl_vector_get_stdev(vusable),
1724  maxusable);
1725 #endif
1726  cpl_msg_debug(__func__, "filled %"CPL_SIZE_FORMAT" of %d spaxels, usable "
1727  "are max. %d of %d pixels per spectrum%s",
1728  itable - cpl_table_count_selected(e3d->dtable), nx * ny,
1729  maxusable, nz, maxusable == nz ? "" : ", cutting table");
1730  if (maxusable != nz) {
1731  /* set array columns to real depth (maximum usable spectral pixels found) */
1732  cpl_table_set_column_depth(e3d->dtable, "DATA_SPE", maxusable);
1733  cpl_table_set_column_depth(e3d->dtable, "QUAL_SPE", maxusable);
1734  cpl_table_set_column_depth(e3d->dtable, "STAT_SPE", maxusable);
1735  }
1736  /* resize table to number of used spaxels */
1737  cpl_table_erase_selected(e3d->dtable); /* erase empty spectra */
1738  cpl_vector_delete(vusable);
1739 
1740  /* fill columns that are the same for all spaxels */
1741  /* by default, all spectra are supposed to be selected */
1742  cpl_table_fill_column_window_int(e3d->dtable, "SELECTED", 0, itable,
1743  1 /* TRUE */);
1744  /* assume that we only used one resolution element */
1745  cpl_table_fill_column_window_int(e3d->dtable, "NSPAX", 0, itable, 1);
1746  /* we have only one group for MUSE */
1747  cpl_table_fill_column_window_int(e3d->dtable, "GROUP_N", 0, itable, 1);
1748 
1749  /* done converting the data, free the datacube */
1750  muse_datacube_delete(cube);
1751 
1752  /* fill info in group table */
1753  e3d->hgroup = cpl_propertylist_new();
1754  cpl_propertylist_append_string(e3d->hgroup, "EXTNAME", "E3D_GRP");
1755  cpl_propertylist_set_comment(e3d->hgroup, "EXTNAME",
1756  "This is the Euro3D group table extension");
1757  e3d->gtable = muse_cpltable_new(muse_euro3dcube_e3d_grp_def, 1); /* single group */
1758  cpl_table_set_int(e3d->gtable, "GROUP_N", 0, 1);
1759  cpl_table_set_string(e3d->gtable, "G_SHAPE", 0, "R"/*ECTANG" ULAR*/);
1760  cpl_table_set_float(e3d->gtable, "G_SIZE1", 0, aParams->dx);
1761  cpl_table_set_float(e3d->gtable, "G_ANGLE", 0, 0.);
1762  cpl_table_set_float(e3d->gtable, "G_SIZE2", 0, aParams->dy);
1763  if (darlambdaref > 0.) {
1764  /* set all properties to NaN as required by the Euro3D specs */
1765  cpl_table_set_float(e3d->gtable, "G_POSWAV", 0, NAN);
1766  cpl_table_set_float(e3d->gtable, "G_AIRMAS", 0, NAN);
1767  cpl_table_set_float(e3d->gtable, "G_PARANG", 0, NAN);
1768  cpl_table_set_float(e3d->gtable, "G_PRESSU", 0, NAN);
1769  cpl_table_set_float(e3d->gtable, "G_TEMPER", 0, NAN);
1770  cpl_table_set_float(e3d->gtable, "G_HUMID", 0, NAN);
1771  } else {
1772  /* because we then don't know any better, use the central MUSE wavelength */
1773  cpl_table_set_float(e3d->gtable, "G_POSWAV", 0,
1774  (kMuseNominalLambdaMin + kMuseNominalLambdaMax) / 2.);
1775  cpl_table_set_float(e3d->gtable, "G_AIRMAS", 0,
1776  muse_astro_airmass(e3d->header));
1777  cpl_table_set_float(e3d->gtable, "G_PARANG", 0,
1778  muse_astro_parangle(e3d->header));
1779  cpl_table_set_float(e3d->gtable, "G_PRESSU", 0,
1781  + muse_pfits_get_pres_start(e3d->header)) / 2.);
1782  cpl_table_set_float(e3d->gtable, "G_TEMPER", 0,
1783  muse_pfits_get_temp(e3d->header) + 273.15);
1784  cpl_table_set_float(e3d->gtable, "G_HUMID", 0,
1785  muse_pfits_get_rhum(e3d->header));
1786  }
1787 
1788  return e3d;
1789 } /* muse_resampling_euro3d() */
1790 
1791 /*---------------------------------------------------------------------------*/
1842 /*---------------------------------------------------------------------------*/
1843 muse_datacube *
1845  muse_pixgrid **aPixgrid)
1846 {
1847  cpl_ensure(aPixtable && aParams, CPL_ERROR_NULL_INPUT, NULL);
1848  cpl_ensure(muse_pixtable_get_type(aPixtable) == MUSE_PIXTABLE_TYPE_FULL,
1849  CPL_ERROR_ILLEGAL_INPUT, NULL);
1850  muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
1851  cpl_ensure(wcstype == MUSE_PIXTABLE_WCS_CELSPH ||
1852  wcstype == MUSE_PIXTABLE_WCS_PIXEL, CPL_ERROR_UNSUPPORTED_MODE,
1853  NULL);
1854 
1855  /* compute or set the size of the output grid depending on *
1856  * the inputs and the data available in the pixel table */
1857  muse_resampling_check_deltas(aPixtable, aParams);
1858  /* compute output sizes; wavelength is different in that it is *
1859  * more useful to contain partly empty areas within the field *
1860  * for the extreme ends, so use ceil() */
1861  int xsize = 0, ysize = 0, zsize = 0;
1862  muse_resampling_compute_size(aPixtable, aParams, &xsize, &ysize, &zsize);
1863  muse_resampling_override_size(&xsize, &ysize, &zsize, aParams);
1864  cpl_ensure(xsize > 0 && ysize > 0 && zsize > 0, CPL_ERROR_ILLEGAL_OUTPUT,
1865  NULL);
1866 
1867  double time = cpl_test_get_walltime();
1868 
1869  /* create the structure for the output datacube */
1870  muse_datacube *cube = cpl_calloc(1, sizeof(muse_datacube));
1871 
1872  /* copy and clean up incoming FITS headers */
1873  cube->header = cpl_propertylist_duplicate(aPixtable->header);
1874  cpl_propertylist_erase_regexp(cube->header,
1875  "^SIMPLE$|^BITPIX$|^NAXIS|^EXTEND$|^XTENSION$|"
1876  "^DATASUM$|^DATAMIN$|^DATAMAX$|^DATAMD5$|"
1877  "^PCOUNT$|^GCOUNT$|^HDUVERS$|^BLANK$|"
1878  "^BZERO$|^BSCALE$|^CHECKSUM$|^INHERIT$|"
1879  "^EXTNAME$|"MUSE_WCS_KEYS"|"MUSE_HDR_PT_REGEXP,
1880  0);
1881  /* set data unit depending on data unit type in the pixel table */
1882  cpl_propertylist_update_string(cube->header, "BUNIT",
1883  cpl_table_get_column_unit(aPixtable->table,
1884  MUSE_PIXTABLE_DATA));
1885  /* set NAXIS for later handling of the WCS */
1886  cpl_propertylist_update_int(cube->header, "NAXIS", 3);
1887  cpl_propertylist_update_int(cube->header, "NAXIS1", xsize);
1888  cpl_propertylist_update_int(cube->header, "NAXIS2", ysize);
1889  cpl_propertylist_update_int(cube->header, "NAXIS3", zsize);
1890  /* if pixel table was astrometrically calibrated, use its WCS headers *
1891  * Axis 1: x or RA, axis 2: y or DEC, axis 3: lambda */
1892  if (wcstype == MUSE_PIXTABLE_WCS_CELSPH) {
1893  cpl_propertylist_copy_property_regexp(cube->header, aPixtable->header,
1894  MUSE_WCS_KEYS, 0);
1895  cpl_propertylist_update_double(cube->header, "CD1_1", -aParams->dx);
1896  cpl_propertylist_update_double(cube->header, "CD2_2", aParams->dy);
1897  cpl_propertylist_update_double(cube->header, "CD1_2", 0.);
1898  cpl_propertylist_update_double(cube->header, "CD2_1", 0.);
1899  /* correct CRPIX by the central pixel coordinate */
1900  cpl_propertylist_update_double(cube->header, "CRPIX1",
1901  cpl_propertylist_get_double(cube->header, "CRPIX1")
1902  + (1. + xsize) / 2.);
1903  cpl_propertylist_update_double(cube->header, "CRPIX2",
1904  cpl_propertylist_get_double(cube->header, "CRPIX2")
1905  + (1. + ysize) / 2.);
1906  /* the pixel table had WCSAXES=2 which should not be propagated to the cube! */
1907  cpl_propertylist_erase(cube->header, "WCSAXES");
1908  } else {
1909  /* add basic WCS headers that are also required by the functions below. */
1910  cpl_propertylist_append_double(cube->header, "CRPIX1", 1.);
1911  cpl_propertylist_append_double(cube->header, "CRVAL1",
1912  cpl_propertylist_get_float(aPixtable->header,
1913  MUSE_HDR_PT_XLO));
1914  /* for linear axes CTYPE doesn't matter, as long as it's not in 4-3 form, *
1915  * see Greisen & Calabretta 2002 A&A 395, 1061 */
1916  cpl_propertylist_append_string(cube->header, "CTYPE1", "PIXEL");
1917  cpl_propertylist_append_string(cube->header, "CUNIT1", "pixel");
1918  cpl_propertylist_append_double(cube->header, "CD1_1", aParams->dx);
1919  cpl_propertylist_append_double(cube->header, "CRPIX2", 1.);
1920  cpl_propertylist_append_double(cube->header, "CRVAL2",
1921  cpl_propertylist_get_float(aPixtable->header,
1922  MUSE_HDR_PT_YLO));
1923  cpl_propertylist_append_string(cube->header, "CTYPE2", "PIXEL");
1924  cpl_propertylist_append_string(cube->header, "CUNIT2", "pixel");
1925  cpl_propertylist_append_double(cube->header, "CD2_2", aParams->dy);
1926  /* fill in empty cross-terms of the CDi_j matrix */
1927  cpl_propertylist_append_double(cube->header, "CD1_2", 0.);
1928  cpl_propertylist_append_double(cube->header, "CD2_1", 0.);
1929  }
1930  if (aParams->tlambda == MUSE_RESAMPLING_DISP_AWAV_LOG) {
1931  cpl_propertylist_append_string(cube->header, "CTYPE3", "AWAV-LOG");
1932  } else {
1933  cpl_propertylist_append_string(cube->header, "CTYPE3", "AWAV");
1934  }
1935  cpl_propertylist_append_string(cube->header, "CUNIT3", "Angstrom");
1936  cpl_propertylist_append_double(cube->header, "CD3_3", aParams->dlambda);
1937  cpl_propertylist_append_double(cube->header, "CRPIX3", 1.);
1938  cpl_propertylist_append_double(cube->header, "CRVAL3",
1939  cpl_propertylist_get_float(aPixtable->header,
1940  MUSE_HDR_PT_LLO));
1941  /* fill in empty cross-terms of the CDi_j matrix */
1942  cpl_propertylist_append_double(cube->header, "CD1_3", 0.);
1943  cpl_propertylist_append_double(cube->header, "CD2_3", 0.);
1944  cpl_propertylist_append_double(cube->header, "CD3_1", 0.);
1945  cpl_propertylist_append_double(cube->header, "CD3_2", 0.);
1946  /* after everything has been transferred from the pixel table headers, *
1947  * see if we want to override anything from an OUTPUT_WCS setup */
1948  muse_resampling_override_wcs(cube, aParams);
1949  /* check, if we need to move the zeropoint to fit the data into the grid */
1950  muse_resampling_fit_data_to_grid(cube, aPixtable);
1951 
1952  if (aParams->method < MUSE_RESAMPLE_NONE) {
1953  /* fill the cube for the data */
1954  cube->data = cpl_imagelist_new();
1955  cube->dq = cpl_imagelist_new();
1956  cube->stat = cpl_imagelist_new();
1957  int i;
1958  for (i = 0; i < zsize; i++) {
1959  cpl_image *image = cpl_image_new(xsize, ysize, CPL_TYPE_FLOAT);
1960  /* set as last image in the list to extend the list size by one */
1961  cpl_imagelist_set(cube->data, image, i);
1962  /* can use the same (empty) image for the variance for now */
1963  cpl_imagelist_set(cube->stat, cpl_image_duplicate(image), i);
1964 
1965  /* the other two are part of their lists now, we can reuse *
1966  * the variable without leaking */
1967  image = cpl_image_new(xsize, ysize, CPL_TYPE_INT);
1968  /* pre-fill with Euro3D status for outside data range */
1969  cpl_image_add_scalar(image, EURO3D_OUTSDRANGE);
1970  cpl_imagelist_set(cube->dq, image, i);
1971  } /* for i (all wavelength planes) */
1972  } /* if method not MUSE_RESAMPLE_NONE */
1973 
1974  muse_utils_memory_dump("muse_resampling_cube() before pixgrid");
1975  /* convert the pixel table into a pixel grid */
1976  muse_pixgrid *pixgrid = muse_pixgrid_create(aPixtable, cube->header,
1977  xsize, ysize, zsize);
1978 
1979  if (aParams->crsigma > 0) {
1980  muse_resampling_crreject(aPixtable, pixgrid, aParams->crtype,
1981  aParams->crsigma);
1982  }
1983  muse_utils_memory_dump("muse_resampling_cube() after pixgrid and crreject");
1984 
1985  double timeinit = cpl_test_get_walltime(),
1986  cpuinit = cpl_test_get_cputime();
1987 
1988  /* do the resampling */
1989  cpl_error_code rc = CPL_ERROR_NONE;
1990  switch (aParams->method) {
1991  case MUSE_RESAMPLE_NEAREST:
1992  cpl_msg_info(__func__, "Starting resampling, using method \"nearest\"");
1993  rc = muse_resampling_cube_nearest(cube, aPixtable, pixgrid);
1994  break;
1996  cpl_msg_info(__func__, "Starting resampling, using method \"renka\" "
1997  "(critical radius rc=%f, loop distance ld=%d)", aParams->rc,
1998  aParams->ld);
1999  rc = muse_resampling_cube_weighted(cube, aPixtable, pixgrid, aParams);
2000  break;
2004  cpl_msg_info(__func__, "Starting resampling, using method \"%s\" (loop "
2005  "distance ld=%d)",
2007  ? "linear"
2009  ? "quadratic"
2010  : "lanczos"),
2011  aParams->ld);
2012  rc = muse_resampling_cube_weighted(cube, aPixtable, pixgrid, aParams);
2013  break;
2015  cpl_msg_info(__func__, "Starting resampling, using method \"drizzle\" "
2016  "(pixfrac f=%.3f,%.3f,%.3f, loop distance ld=%d)",
2017  aParams->pfx, aParams->pfy, aParams->pfl, aParams->ld);
2018  rc = muse_resampling_cube_weighted(cube, aPixtable, pixgrid, aParams);
2019  break;
2020  case MUSE_RESAMPLE_NONE:
2021  cpl_msg_debug(__func__, "Method %d (no resampling)", aParams->method);
2022  break;
2023  default:
2024  cpl_msg_error(__func__, "Don't know this resampling method: %d",
2025  aParams->method);
2026  cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
2027  rc = CPL_ERROR_UNSUPPORTED_MODE;
2028  }
2029  muse_utils_memory_dump("muse_resampling_cube() after resampling");
2030 
2031  double timefini = cpl_test_get_walltime(),
2032  cpufini = cpl_test_get_cputime();
2033 
2034  /* now that we have resampled we can either remove the pixel grid or save it */
2035  if (aPixgrid) {
2036  *aPixgrid = pixgrid;
2037  } else {
2038  muse_pixgrid_delete(pixgrid);
2039  }
2040 
2041  cpl_msg_debug(__func__, "resampling took %.3fs (wall-clock) and %.3fs "
2042  "(%.3fs CPU, %d CPUs) for muse_resampling_cube*() alone",
2043  timefini - time, timefini - timeinit, cpufini - cpuinit,
2044  omp_get_max_threads());
2045  if (rc != CPL_ERROR_NONE) {
2046  cpl_msg_error(__func__, "resampling failed: %s", cpl_error_get_message());
2047  muse_datacube_delete(cube);
2048  return NULL;
2049  }
2050 
2051  muse_utils_memory_dump("muse_resampling_cube() end");
2052  return cube;
2053 } /* muse_resampling_cube() */
2054 
2055 /*---------------------------------------------------------------------------*/
2085 /*---------------------------------------------------------------------------*/
2086 muse_image *
2088  muse_datacube *aCube, cpl_table *aFilter,
2089  muse_resampling_params *aParams)
2090 {
2091  cpl_ensure(aPixtable && aPixtable->table && aPixgrid && aParams &&
2092  aCube && aCube->header, CPL_ERROR_NULL_INPUT, NULL);
2093  cpl_ensure(aParams->method < MUSE_RESAMPLE_NONE &&
2094  aParams->method > MUSE_RESAMPLE_NEAREST, CPL_ERROR_ILLEGAL_INPUT,
2095  NULL);
2096 
2097  muse_wcs *wcs = muse_wcs_new(aCube->header);
2099  cpl_errorstate prestate = cpl_errorstate_get();
2100  float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
2101  *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
2102  *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
2103  *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
2104  *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT),
2105  *wght = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_WEIGHT);
2106  if (!cpl_errorstate_is_equal(prestate)) {
2107  cpl_errorstate_set(prestate); /* recover from missing weight column */
2108  }
2109  int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
2110 
2111  /* If our data was astrometrically calibrated, we need to scale the *
2112  * data units to the pixel size in all three dimensions so that the *
2113  * radius computation works again. *
2114  * Otherwise dx~5.6e-5deg won't contribute to the weighting at all. */
2115  double xnorm = 1., ynorm = 1.,
2116  ptxoff = 0., /* zero by default ... */
2117  ptyoff = 0.; /* for pixel coordinates */
2118  if (wcs->iscelsph) {
2119  muse_wcs_get_scales(aPixtable->header, &xnorm, &ynorm);
2120  xnorm = 1. / xnorm;
2121  ynorm = 1. / ynorm;
2122  /* need to use the real coordinate offset for celestial spherical */
2123  ptxoff = cpl_propertylist_get_double(aPixtable->header, "CRVAL1");
2124  ptyoff = cpl_propertylist_get_double(aPixtable->header, "CRVAL2");
2125  }
2126  /* scale the input critical radius by the voxel radius */
2127  double renka_rc = aParams->rc /* XXX beware of rotation! */
2128  * sqrt(pow(wcs->cd11*xnorm, 2) + pow(wcs->cd22*xnorm, 2));
2129  /* loop distance (to take into account surrounding pixels) verification */
2130  int ld = aParams->ld;
2131  if (ld <= 0) {
2132  ld = 1;
2133  cpl_msg_info(__func__, "Overriding loop distance ld=%d", ld);
2134  }
2135  /* pixel sizes in the spatial directions, scaled by pixfrac, and *
2136  * output pixel sizes (absolute values), as needed for drizzle */
2137  double xsz = aParams->pfx / xnorm,
2138  ysz = aParams->pfy / ynorm,
2139  xout = fabs(wcs->cd11), yout = fabs(wcs->cd22);
2140 
2141  muse_image *image = muse_image_new();
2142  image->data = cpl_image_new(aPixgrid->size_x, aPixgrid->size_y, CPL_TYPE_FLOAT);
2143  image->dq = cpl_image_new(aPixgrid->size_x, aPixgrid->size_y, CPL_TYPE_INT);
2144  image->stat = cpl_image_new(aPixgrid->size_x, aPixgrid->size_y, CPL_TYPE_FLOAT);
2145  image->header = cpl_propertylist_duplicate(aCube->header);
2146  cpl_propertylist_erase_regexp(image->header, "^C...*3$|^CD3_.$", 0);
2147  float *pdata = cpl_image_get_data_float(image->data),
2148  *pstat = cpl_image_get_data_float(image->stat);
2149  int *pdq = cpl_image_get_data_int(image->dq);
2150 
2151  /* check if we need to use the variance for weighting */
2152  cpl_boolean usevariance = CPL_FALSE;
2153  if (getenv("MUSE_COLLAPSE_USE_VARIANCE") &&
2154  atoi(getenv("MUSE_COLLAPSE_USE_VARIANCE")) > 0) {
2155  usevariance = CPL_TRUE;
2156  }
2157 
2158  /* find filter range */
2159  double lmin, lmax;
2160  if (aFilter) {
2161  const double *flbda = cpl_table_get_data_double_const(aFilter, "lambda"),
2162  *fresp = cpl_table_get_data_double_const(aFilter, "throughput");
2163  int l = 0, nl = cpl_table_get_nrow(aFilter);
2164  while (l < nl && fabs(fresp[l]) < DBL_EPSILON) {
2165  lmin = flbda[l++];
2166  }
2167  l = nl - 1;
2168  while (l > 0 && fabs(fresp[l]) < DBL_EPSILON) {
2169  lmax = flbda[l--];
2170  }
2171  cpl_msg_debug(__func__, "filter wavelength range %.1f..%.1fA", lmin, lmax);
2172  } else {
2173  lmin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LLO);
2174  lmax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LHI);
2175  cpl_msg_debug(__func__, "full wavelength range %.1f..%.1fA", lmin, lmax);
2176  }
2177 
2178  int i;
2179  #pragma omp parallel for default(none) /* as req. by Ralf */ \
2180  shared(aFilter, aParams, aPixgrid, data, dq, lbda, ld, lmax, lmin, \
2181  pdata, pdq, pstat, ptxoff, ptyoff, renka_rc, stat, usevariance,\
2182  wcs, wght, xnorm, xout, xpos, xsz, ynorm, yout, ypos, ysz)
2183  for (i = 0; i < aPixgrid->size_x; i++) {
2184  int j;
2185  for (j = 0; j < aPixgrid->size_y; j++) {
2186  /* x and y position of center of current grid cell (i, j start at 0) */
2187  double x, y;
2188  if (wcs->iscelsph) {
2189  muse_wcs_celestial_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
2190  } else {
2191  muse_wcs_projplane_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
2192  }
2193  double sumdata = 0, sumstat = 0, sumweight = 0;
2194  cpl_size npoints = 0;
2195 
2196  /* loop through surrounding cells and their contained pixels */
2197  int i2;
2198  for (i2 = i - ld; i2 <= i + ld; i2++) {
2199  int j2;
2200  for (j2 = j - ld; j2 <= j + ld; j2++) {
2201  /* loop over full wavelength range! */
2202  int l2;
2203  for (l2 = 0; l2 < aPixgrid->size_z; l2++) {
2204  cpl_size idx2 = muse_pixgrid_get_index(aPixgrid, i2, j2, l2, CPL_FALSE),
2205  n, n_rows2 = muse_pixgrid_get_count(aPixgrid, idx2);
2206  const cpl_size *rows2 = muse_pixgrid_get_rows(aPixgrid, idx2);
2207  for (n = 0; n < n_rows2; n++) {
2208  if (dq[rows2[n]]) { /* exclude all bad pixels */
2209  continue;
2210  }
2211  if (usevariance && !isnormal(stat[rows2[n]])) {
2212  continue; /* weighting by inv. variance doesn`t work with extremes */
2213  }
2214  if (lbda[rows2[n]] > lmax || lbda[rows2[n]] < lmin) {
2215  continue; /* exclude pixels outside the filter range */
2216  }
2217 
2218  double dx = fabs(x - (xpos[rows2[n]] + ptxoff)),
2219  dy = fabs(y - (ypos[rows2[n]] + ptyoff)),
2220  r2 = 0;
2221  if (wcs->iscelsph) {
2222  /* Make the RA distance comparable to pixel *
2223  * sizes, see muse_resampling_cube_weighted(). */
2224  dx *= cos(y * CPL_MATH_RAD_DEG);
2225  }
2226  if (aParams->method != MUSE_RESAMPLE_WEIGHTED_DRIZZLE) {
2227  dx *= xnorm;
2228  dy *= ynorm;
2229  r2 = dx*dx + dy*dy;
2230  }
2231  double weight = 0.;
2232  if (aParams->method == MUSE_RESAMPLE_WEIGHTED_RENKA) {
2233  weight = muse_resampling_weight_function_renka(sqrt(r2), renka_rc);
2234  } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_DRIZZLE) {
2235  weight = muse_resampling_weight_function_drizzle(xsz, ysz, 1.,
2236  xout, yout, 1.,
2237  dx, dy, 0.);
2238  } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_LINEAR) {
2239  weight = muse_resampling_weight_function_linear(sqrt(r2));
2240  } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_QUADRATIC) {
2241  weight = muse_resampling_weight_function_quadratic(r2);
2242  } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_LANCZOS) {
2243  weight = muse_resampling_weight_function_lanczos(dx, dy, 0., ld);
2244  }
2245 
2246  if (wght) { /* the pixel table does contain weights */
2247  /* apply it on top of the weight computed here */
2248  weight *= wght[rows2[n]];
2249  }
2250  /* apply the filter curve on top of the weights, if a filter *
2251  * was given; otherwise the weight stays unchanged */
2252  if (aFilter) {
2253  weight *= muse_flux_response_interpolate(aFilter, lbda[rows2[n]],
2254  NULL, MUSE_FLUX_RESP_FILTER);
2255  }
2256  if (usevariance) { /* weight by inverse variance */
2257  weight /= stat[rows2[n]];
2258  }
2259  sumweight += weight;
2260  sumdata += data[rows2[n]] * weight;
2261  sumstat += stat[rows2[n]] * weight*weight;
2262  npoints++;
2263  } /* for n (all pixels in grid cell) */
2264  } /* for l2 (lambda direction) */
2265  } /* for j2 (y direction) */
2266  } /* for i2 (x direction) */
2267 
2268  /* if no points were found, we cannot divide by the summed weight *
2269  * and don't need to set the output pixel value (it's 0 already), *
2270  * only set the relevant Euro3D bad pixel flag */
2271  if (!npoints || !isnormal(sumweight)) {
2272  pdq[i + j * aPixgrid->size_x] = EURO3D_MISSDATA;
2273  continue;
2274  }
2275 
2276  /* divide results by weight of summed pixels */
2277  sumdata /= sumweight;
2278  sumstat /= sumweight*sumweight;
2279  pdata[i + j * aPixgrid->size_x] = sumdata;
2280  pstat[i + j * aPixgrid->size_x] = sumstat;
2281  pdq[i + j * aPixgrid->size_x] = EURO3D_GOODPIXEL; /* now we can mark it as good */
2282  } /* for j (y direction) */
2283  } /* for i (x direction) */
2284  cpl_free(wcs);
2285 
2286  return image;
2287 } /* muse_resampling_collapse_pixgrid() */
2288 
2289 
2290 /*----------------------------------------------------------------------------*
2291  * Resampling and collapsing in 2D *
2292  *----------------------------------------------------------------------------*/
2293 
2294 /*---------------------------------------------------------------------------*/
2309 /*---------------------------------------------------------------------------*/
2310 static cpl_error_code
2311 muse_resampling_image_nearest(muse_image *aImage, muse_pixtable *aPixtable,
2312  muse_pixgrid *aPixgrid)
2313 {
2314  cpl_ensure_code(aImage && aPixtable && aPixgrid, CPL_ERROR_NULL_INPUT);
2315  aImage->data = cpl_image_new(aPixgrid->size_x, aPixgrid->size_z,
2316  CPL_TYPE_FLOAT);
2317 
2318  double crval2 = cpl_propertylist_get_double(aImage->header, "CRVAL2"),
2319  cd22 = cpl_propertylist_get_double(aImage->header, "CD2_2");
2320  /* get all (relevant) table columns for easy pointer access */
2321  float *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
2322  *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA);
2323  float *imagedata = cpl_image_get_data_float(aImage->data);
2324 
2325 #ifdef ESO_ENABLE_DEBUG
2326  int debug = 0;
2327  if (getenv("MUSE_DEBUG_NEAREST")) {
2328  debug = atoi(getenv("MUSE_DEBUG_NEAREST"));
2329  }
2330 #endif
2331 
2332  int i;
2333  for (i = 0; i < aPixgrid->size_x; i++) {
2334  int j;
2335  for (j = 0; j < aPixgrid->size_z; j++) {
2336  cpl_size idx = muse_pixgrid_get_index(aPixgrid, i, 0, j, CPL_FALSE),
2337  n_rows = muse_pixgrid_get_count(aPixgrid, idx);
2338  const cpl_size *rows = muse_pixgrid_get_rows(aPixgrid, idx);
2339  if (n_rows == 1) {
2340  /* if there is only one pixel in the cell, just use it */
2341  imagedata[i + j * aPixgrid->size_x] = data[rows[0]];
2342 #ifdef ESO_ENABLE_DEBUG
2343  if (debug) {
2344  printf("only: %f\n", data[rows[0]]);
2345  fflush(stdout);
2346  }
2347 #endif
2348  } else if (n_rows >= 2) {
2349  /* loop through all available values and take the closest one */
2350  cpl_size n, nbest = -1;
2351  double dbest = FLT_MAX; /* some unlikely large value to start with*/
2352  for (n = 0; n < n_rows; n++) {
2353  /* the differences for the cell center and the current pixel; *
2354  * for simplicitly, just compare the difference in lambda */
2355  double dlambda = fabs(j * cd22 + crval2 - lbda[rows[n]]);
2356 #ifdef ESO_ENABLE_DEBUG
2357  if (debug) {
2358  printf("%f, %f, d: %f best: %f (%f)\n",
2359  j * cd22 + crval2, lbda[rows[n]],
2360  dlambda, dbest, data[rows[n]]);
2361  }
2362 #endif
2363  if (dlambda < dbest) {
2364  nbest = n;
2365  dbest = dlambda;
2366  }
2367  }
2368  imagedata[i + j * aPixgrid->size_x] = data[rows[nbest]];
2369 #ifdef ESO_ENABLE_DEBUG
2370  if (debug) {
2371  printf("closest: %f\n", data[rows[nbest]]);
2372  fflush(stdout);
2373  }
2374 #endif
2375  } else {
2376  /* npix == 0: do nothing, pixel stays zero */
2377  }
2378  } /* for j (y direction = wavelength planes) */
2379  } /* for i (x direction) */
2380 
2381  return CPL_ERROR_NONE;
2382 } /* muse_resampling_image_nearest() */
2383 
2384 /*---------------------------------------------------------------------------*/
2400 /*---------------------------------------------------------------------------*/
2401 static cpl_error_code
2402 muse_resampling_image_weighted(muse_image *aImage, muse_pixtable *aPixtable,
2403  muse_pixgrid *aPixgrid)
2404 {
2405  cpl_ensure_code(aImage && aPixtable && aPixgrid, CPL_ERROR_NULL_INPUT);
2406  aImage->data = cpl_image_new(aPixgrid->size_x, aPixgrid->size_z,
2407  CPL_TYPE_FLOAT);
2408 
2409  double crval2 = cpl_propertylist_get_double(aImage->header, "CRVAL2"),
2410  cd22 = cpl_propertylist_get_double(aImage->header, "CD2_2");
2411  /* get all (relevant) table columns for easy pointer access */
2412  float *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
2413  *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA);
2414  float *imagedata = cpl_image_get_data_float(aImage->data);
2415 
2416 #ifdef ESO_ENABLE_DEBUG
2417  int debug = 0, debugx = 0, debugz = 0;
2418  if (getenv("MUSE_DEBUG_WEIGHTED")) {
2419  debug = atoi(getenv("MUSE_DEBUG_WEIGHTED"));
2420  }
2421  if (debug & 128) { /* need coordinates */
2422  if (getenv("MUSE_DEBUG_WEIGHTED_X")) {
2423  debugx = atoi(getenv("MUSE_DEBUG_WEIGHTED_X"));
2424  if (debugx < 1 || debugx > aPixgrid->size_x) {
2425  debugx = 0;
2426  }
2427  }
2428  if (getenv("MUSE_DEBUG_WEIGHTED_Z")) {
2429  debugz = atoi(getenv("MUSE_DEBUG_WEIGHTED_Z"));
2430  if (debugz < 1 || debugz > aPixgrid->size_z) {
2431  debugz = 0;
2432  }
2433  }
2434  }
2435 #endif
2436 
2437  int ld = 1; /* loop distance, to take into account surrounding pixels */
2438  /* Renka critical radius; only 1.25 to 1.5 makes sense in *
2439  * this 1D interpolation, hardcode to 1.25 for the moment */
2440  double renka_rc = 1.25;
2441  /* scale critical radius by the scale factor in wavelength */
2442  renka_rc *= cd22;
2443 
2444  int i;
2445  for (i = 0; i < aPixgrid->size_x; i++) {
2446  int j;
2447  for (j = 0; j < aPixgrid->size_z; j++) {
2448  double sumdata = 0, sumweight = 0,
2449  lambda = j * cd22 + crval2;
2450  int npoints = 0;
2451 #ifdef ESO_ENABLE_DEBUG
2452 #define MAX_NPIX_2D 50 /* allow for quite a few pixels */
2453  int *pointlist = NULL;
2454  double *pointweights = NULL;
2455  if (debug & 128 && i+1 == debugx && j+1 == debugz) {
2456  pointlist = cpl_calloc(MAX_NPIX_2D, sizeof(int));
2457  pointweights = cpl_malloc(MAX_NPIX_2D * sizeof(double));
2458  }
2459 #endif
2460 
2461  /* go through this and the surrounding cells and their contained pixels */
2462  int j2;
2463  for (j2 = j - ld; j2 <= j + ld; j2++) {
2464  cpl_size idx = muse_pixgrid_get_index(aPixgrid, i, 0, j2, CPL_FALSE),
2465  n, n_rows = muse_pixgrid_get_count(aPixgrid, idx);
2466  const cpl_size *rows = muse_pixgrid_get_rows(aPixgrid, idx);
2467  for (n = 0; n < n_rows; n++) {
2468  double dlambda = fabs(lambda - lbda[rows[n]]),
2469  weight = muse_resampling_weight_function_renka(dlambda,
2470  renka_rc);
2471  sumweight += weight;
2472  sumdata += data[rows[n]] * weight;
2473  npoints++;
2474 #ifdef ESO_ENABLE_DEBUG
2475  if (debug & 128 && i+1 == debugx && j+1 == debugz &&
2476  npoints-1 < MAX_NPIX_2D) {
2477  /* store row number instead of index, because we cannot be zero: */
2478  pointlist[npoints-1] = rows[n] + 1;
2479  pointweights[npoints-1] = weight;
2480  }
2481 
2482  if (debug & 256) {
2483  printf(" pixel %4d,%4d (%8"CPL_SIZE_FORMAT"): %2"CPL_SIZE_FORMAT
2484  " %2"CPL_SIZE_FORMAT" %f, %f -> %f ==> %f %f (%d)\n",
2485  i, j2, idx, n, muse_pixgrid_get_count(aPixgrid, idx), dlambda,
2486  data[muse_pixgrid_get_rows(aPixgrid, idx)[n]],
2487  weight, sumweight, sumdata, npoints);
2488  }
2489 #endif
2490  } /* for n (all pixels in grid cell) */
2491  } /* for j2 (x direction) */
2492 
2493 #ifdef ESO_ENABLE_DEBUG
2494  if (debug & 128 && i+1 == debugx && j+1 == debugz) {
2495  printf("pixelnumber weight ");
2496  muse_pixtable_dump(aPixtable, 0, 0, 2);
2497  int m = -1;
2498  while (++m < MAX_NPIX_2D && pointlist[m] != 0) {
2499  /* access row using index again: */
2500  printf("%12d %8.5f ", pointlist[m] - 1, pointweights[m]);
2501  muse_pixtable_dump(aPixtable, pointlist[m] - 1, 1, 0);
2502  }
2503  fflush(stdout);
2504  cpl_free(pointlist);
2505  cpl_free(pointweights);
2506  }
2507  if (debug) {
2508  fflush(stdout); /* flush the output in any debug case */
2509  }
2510 #endif
2511  /* if no points were found, we cannot divide by the summed weight *
2512  * and don't need to set the output pixel value (it's 0 already) */
2513  if (!npoints) {
2514  continue;
2515  }
2516  /* divide results by weight of summed pixels */
2517  imagedata[i + j * aPixgrid->size_x] = sumdata / sumweight;
2518  } /* for j (y direction = wavelength planes) */
2519  } /* for i (x direction) */
2520 
2521  return CPL_ERROR_NONE;
2522 } /* muse_resampling_image_weighted() */
2523 
2524 /*---------------------------------------------------------------------------*/
2541 /*---------------------------------------------------------------------------*/
2542 static muse_image *
2543 muse_resampling_image_selected(muse_pixtable *aPixtable,
2544  muse_resampling_type aMethod,
2545  double aDX, double aLambdaMin,
2546  double aLambdaMax, double aDLambda)
2547 {
2548  double dlambda = aDLambda;
2549  float xmin = 0.;
2550  muse_pixgrid *pixgrid = muse_pixgrid_2d_create(aPixtable->table, aDX,
2551  aLambdaMin, aLambdaMax,
2552  dlambda, &xmin);
2553 
2554  /* create output image and add FITS header to it (we already need *
2555  * CRVAL2 and CD2_2 in the sub-functions called from here) */
2556  muse_image *image = muse_image_new();
2557  image->header = cpl_propertylist_new();
2558  const char *unit = cpl_table_get_column_unit(aPixtable->table, "data");
2559  cpl_propertylist_append_string(image->header, "BUNIT", unit);
2560  /* copy input header and append the necessary WCS keywords *
2561  * WCSAXES must preceed all other keywords, but defaults to NAXIS, *
2562  * so we can ignore this here for simplicity */
2563  cpl_propertylist_copy_property_regexp(image->header, aPixtable->header,
2564  MUSE_HDR_PT_REGEXP"|"MUSE_WCS_KEYS, 1);
2565  cpl_propertylist_append_double(image->header, "CRPIX1", 1.);
2566  cpl_propertylist_append_double(image->header, "CRPIX2", 1.);
2567  cpl_propertylist_append_double(image->header, "CRVAL1", 1.);
2568  cpl_propertylist_append_double(image->header, "CRVAL2", aLambdaMin);
2569  cpl_propertylist_append_double(image->header, "CD1_1", 1.);
2570  cpl_propertylist_append_double(image->header, "CD2_2", dlambda);
2571  /* units verified using FITS standard v3.0 */
2572  cpl_propertylist_append_string(image->header, "CUNIT1", "pixel");
2573  cpl_propertylist_append_string(image->header, "CUNIT2", "Angstrom");
2574  /* for linear axes CTYPE doesn't matter, as long as it's not in 4-3 form, *
2575  * see Greisen & Calabretta 2002 A&A 395, 1061 */
2576  cpl_propertylist_append_string(image->header, "CTYPE1", "PIXEL");
2577  /* although DS9 seems to like "LAMBDA" better, the type (for the wavelength *
2578  * direction) was verified using Greisen et al. 2006 A&A 446, 747 */
2579  cpl_propertylist_append_string(image->header, "CTYPE2", "AWAV");
2580  /* fill in empty cross-terms of the CDi_j matrix */
2581  cpl_propertylist_append_double(image->header, "CD1_2", 0.);
2582  cpl_propertylist_append_double(image->header, "CD2_1", 0.);
2583 
2584  /* do the resampling */
2585  cpl_error_code rc = CPL_ERROR_NONE;
2586  switch (aMethod) {
2587  case MUSE_RESAMPLE_NEAREST:
2588  rc = muse_resampling_image_nearest(image, aPixtable, pixgrid);
2589  break;
2590  default: /* MUSE_RESAMPLE_WEIGHTED_RENKA */
2591  /* the caller is responsible for checking for other methods */
2592  rc = muse_resampling_image_weighted(image, aPixtable, pixgrid);
2593  }
2594 
2595  muse_pixgrid_delete(pixgrid);
2596  if (rc != CPL_ERROR_NONE) {
2597  cpl_msg_error(__func__, "resampling failed: %s", cpl_error_get_message());
2598  muse_image_delete(image);
2599  return NULL;
2600  }
2601 
2602  return image;
2603 } /* muse_resampling_image_selected() */
2604 
2605 /*---------------------------------------------------------------------------*/
2645 /*---------------------------------------------------------------------------*/
2646 muse_image *
2648  double aDX, double aDLambda)
2649 {
2650  cpl_ensure(aPixtable, CPL_ERROR_NULL_INPUT, NULL);
2651  if (aDLambda == 0.0) {
2652  aDLambda = kMuseSpectralSamplingA;
2653  }
2654  muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
2655  cpl_ensure(wcstype == MUSE_PIXTABLE_WCS_CELSPH ||
2656  wcstype == MUSE_PIXTABLE_WCS_PIXEL, CPL_ERROR_UNSUPPORTED_MODE,
2657  NULL);
2658 
2659  /* check method and give info message */
2660  switch (aMethod) {
2661  case MUSE_RESAMPLE_NEAREST:
2662  cpl_msg_info(__func__, "Using nearest neighbor sampling (%d) along "
2663  "wavelengths.", aMethod);
2664  break;
2666  cpl_msg_info(__func__, "Using renka-weighted interpolation (%d) along "
2667  "wavelengths.", aMethod);
2668  break;
2669  default:
2670  cpl_msg_error(__func__, "Don't know this resampling method: %d", aMethod);
2671  cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
2672  return NULL;
2673  }
2674 
2675  /* compute the size of the output grid depending on the inputs and the *
2676  * data available in the pixel table */
2677  float lmin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LLO),
2678  lmax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LHI);
2679 
2681  /* Pretend that we are dealing with only one slice. *
2682  * We can do that because the pixel table was created in such *
2683  * a way that no slice overlaps in x-direction with another. */
2684  return muse_resampling_image_selected(aPixtable, aMethod,
2685  aDX == 0.0 ? 1. : aDX,
2686  lmin, lmax, aDLambda);
2687  }
2688 
2689  /* For pixel tables with full geometry information, resample each *
2690  * slice of each IFU into an image one by one and stack them next *
2691  * to each other (horizontally) in a large output image. */
2692  muse_pixtable **slice_pixtable = muse_pixtable_extracted_get_slices(aPixtable);
2693  int n_slices = muse_pixtable_extracted_get_size(slice_pixtable);
2694 
2695  double dx = aDX;
2696  if (dx == 0.0) {
2698  dx = 1.0;
2699  } else {
2700  double xsc, ysc;
2701  muse_wcs_get_scales(aPixtable->header, &xsc, &ysc);
2702  /* XXX need such an extra 1.2 factor to get good results: */
2703  dx = 1.20 * xsc;
2704  }
2705  }
2706  cpl_msg_debug(__func__, "Resampling %d slices to a 2D image, using bins of"
2707  " %e %s x %.3f Angstrom", n_slices, dx,
2708  cpl_table_get_column_unit(aPixtable->table, MUSE_PIXTABLE_XPOS),
2709  aDLambda);
2710 
2711  /* build one stacked image per slice */
2712  muse_image *slice_image[n_slices];
2713  int i_slice;
2714  #pragma omp parallel for default(none) /* as req. by Ralf */ \
2715  shared(aDLambda, aMethod, dx, i_slice, lmax, lmin, n_slices, \
2716  slice_image, slice_pixtable)
2717  for (i_slice = 0; i_slice < n_slices; i_slice++) {
2718 #if 0
2719  cpl_msg_debug(__func__, "slice pixel table %d", i_slice + 1);
2720 #endif
2721  muse_pixtable *pt = slice_pixtable[i_slice];
2722  if (muse_pixtable_get_nrow(pt) < 1) {
2723  slice_image[i_slice] = NULL;
2724  continue;
2725  }
2726  slice_image[i_slice] = muse_resampling_image_selected(pt, aMethod, dx,
2727  lmin, lmax, aDLambda);
2728  } /* for i_slice */
2729 
2730  /* concatenate all images into one large output image */
2731  muse_image *image = muse_image_new();
2732  for (i_slice = 0; i_slice < n_slices; i_slice++) {
2733  muse_image *slice = slice_image[i_slice];
2734 #if 0
2735  cpl_msg_debug(__func__, "slice %d: %ldx%ld", i_slice + 1,
2736  cpl_image_get_size_x(slice->data), cpl_image_get_size_y(slice->data));
2737 #endif
2738  if (slice == NULL) {
2739  continue;
2740  }
2741  if (!image->header) {
2742  /* use FITS header from the first slice */
2743  image->header = cpl_propertylist_duplicate(slice->header);
2744  }
2745  cpl_image *data = muse_cplimage_concat_x(image->data, slice->data);
2746  cpl_image_delete(image->data);
2747  image->data = data;
2748  if (slice->dq) {
2749  cpl_image *dq = muse_cplimage_concat_x(image->dq, slice->dq);
2750  cpl_image_delete(image->dq);
2751  image->dq = dq;
2752  }
2753  if (slice->stat) {
2754  cpl_image *stat = muse_cplimage_concat_x(image->stat, slice->stat);
2755  cpl_image_delete(image->stat);
2756  image->stat = stat;
2757  }
2758  muse_image_delete(slice);
2759  slice_image[i_slice] = NULL;
2760  } /* for i_slice */
2761  muse_pixtable_extracted_delete(slice_pixtable);
2762 
2763  return image;
2764 } /* muse_resampling_image() */
2765 
2766 /*---------------------------------------------------------------------------*/
2779 /*---------------------------------------------------------------------------*/
2780 cpl_table *
2781 muse_resampling_spectrum(muse_pixtable *aPixtable, double aBinwidth)
2782 {
2783  cpl_ensure(aPixtable != NULL, CPL_ERROR_NULL_INPUT, NULL);
2784 
2785  /* determine spectral range */
2786  double lmin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LLO);
2787  double lmax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LHI);
2788  cpl_size lsize = floor((lmax - lmin) / aBinwidth) + 2;
2789 
2790  /* create the empty spectral table */
2791  cpl_table *spectrum = muse_cpltable_new(muse_dataspectrum_def, lsize);
2792  cpl_table_fill_column_window(spectrum, "lambda", 0, lsize, 0);
2793  cpl_table_fill_column_window(spectrum, "data", 0, lsize, 0);
2794  cpl_table_fill_column_window(spectrum, "stat", 0, lsize, 0);
2795  cpl_table_fill_column_window(spectrum, "dq", 0, lsize, 0);
2796  double *spec_data = cpl_table_get_data_double(spectrum, "data");
2797  double *spec_stat = cpl_table_get_data_double(spectrum, "stat");
2798  double *spec_lbda = cpl_table_get_data_double(spectrum, "lambda");
2799  /* inherit spectral data units from the pixel table */
2800  cpl_table_set_column_unit(spectrum, "data",
2801  cpl_table_get_column_unit(aPixtable->table,
2802  MUSE_PIXTABLE_DATA));
2803  cpl_table_set_column_unit(spectrum, "stat",
2804  cpl_table_get_column_unit(aPixtable->table,
2805  MUSE_PIXTABLE_STAT));
2806 
2807  /* add a temporary column for the integrated weight of each bin */
2808  cpl_table_new_column(spectrum, "weight", CPL_TYPE_DOUBLE);
2809  cpl_table_fill_column_window(spectrum, "weight", 0, lsize, 0);
2810  double *spec_weight = cpl_table_get_data_double(spectrum, "weight");
2811 
2812  /* accessors for the pixel table */
2813  float *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA);
2814  float *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA);
2815  float *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
2816  float *wght = NULL;
2817  if (cpl_table_has_column(aPixtable->table, MUSE_PIXTABLE_WEIGHT)) {
2818  wght = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_WEIGHT);
2819  }
2820  int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
2821 
2822  /* loop through all selected pixel table rows and sort them *
2823  * into the output bins, assigning them linear weights */
2824  cpl_array *asel = cpl_table_where_selected(aPixtable->table);
2825  const cpl_size *sel = cpl_array_get_data_cplsize_const(asel);
2826  cpl_size isel, nsel = cpl_array_get_size(asel);
2827  cpl_msg_debug(__func__, "Resample spectrum from %" CPL_SIZE_FORMAT
2828  " pixels", nsel);
2829  for (isel = 0; isel < nsel; isel++) {
2830  cpl_size i_row = sel[isel];
2831  if ((dq[i_row] != EURO3D_GOODPIXEL)) {
2832  continue;
2833  }
2834  /* compute target bins and weights */
2835  double l = (lbda[i_row] - lmin) / aBinwidth;
2836  if (l < 0) l = 0;
2837  cpl_size il = floor(l);
2838  l -= il;
2839  double w0 = 1-l;
2840  double w1 = l;
2841  if (wght) {
2842  w0 *= wght[i_row];
2843  w1 *= wght[i_row];
2844  }
2845  /* weight the data in the two involved bins */
2846  spec_data[il] += w0 * data[i_row];
2847  spec_data[il+1] += w1 * data[i_row];
2848  spec_stat[il] += w0 * stat[i_row];
2849  spec_stat[il+1] += w1 * stat[i_row];
2850  spec_weight[il] += w0;
2851  spec_weight[il+1] += w1;
2852  } /* for isel (selected pixel table rows) */
2853  cpl_array_delete(asel);
2854 
2855  /* assign wavelengths and erase empty spectral bins */
2856  cpl_size ispec;
2857  for (ispec = 0; ispec < lsize; ispec++) {
2858  if (spec_weight[ispec] > 0) {
2859  spec_lbda[ispec] = ispec * aBinwidth + lmin;
2860  cpl_table_unselect_row(spectrum, ispec);
2861  }
2862  } /* for ispec (all spectrum points) */
2863  cpl_table_erase_selected(spectrum);
2864 
2865  /* apply the weights and delete the temporary weight column */
2866  cpl_table_divide_columns(spectrum, "data", "weight");
2867  cpl_table_divide_columns(spectrum, "stat", "weight");
2868  cpl_table_erase_column(spectrum, "weight");
2869  return spectrum;
2870 } /* muse_resampling_spectrum() */
2871 
muse_pixtable_wcs
State of the astrometric calibration of a MUSE pixel table.
Structure definition of a MUSE datacube.
Definition: muse_datacube.h:48
muse_image * muse_resampling_image(muse_pixtable *aPixtable, muse_resampling_type aMethod, double aDX, double aDLambda)
Resample a pixel table onto a two-dimensional regular grid.
cpl_error_code muse_wcs_get_scales(cpl_propertylist *aHeader, double *aXScale, double *aYScale)
Compute the spatial scales (in degrees) from the FITS header WCS.
Definition: muse_wcs.c:1491
cpl_error_code muse_wcs_pixel_from_projplane(cpl_propertylist *aHeader, double aX, double aY, double *aXOut, double *aYOut)
Convert from projection plane coordinates to pixel coordinates.
Definition: muse_wcs.c:1405
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
Definition: muse_image.c:84
void muse_pixtable_extracted_delete(muse_pixtable **aPixtables)
Delete a pixel table array.
muse_pixgrid * muse_pixgrid_2d_create(cpl_table *aTable, double aDX, double aZMin, double aZMax, double aDZ, float *aXMin)
Convert selected rows of a pixel table into 2D pixgrid, linking the grid points to entries (=rows) in...
Definition: muse_pixgrid.c:330
#define MUSE_HDR_PT_XLO
FITS header keyword contains the lower limit of the data in x-direction.
cpl_size muse_pixtable_extracted_get_size(muse_pixtable **aPixtables)
Get the size of an array of extracted pixel tables.
cpl_size muse_pixtable_get_nrow(muse_pixtable *aPixtable)
get the number of rows within the pixel table
double muse_pfits_get_rhum(const cpl_propertylist *aHeaders)
find out the relavtive humidity (in %)
Definition: muse_pfits.c:786
A structure containing a spatial two-axis WCS.
Definition: muse_wcs.h:93
double cd22
Definition: muse_wcs.h:96
cpl_error_code muse_resampling_params_set_wcs(muse_resampling_params *aParams, const cpl_propertylist *aWCS)
Set an output WCS (and wavelength scale) in the resampling parameters.
muse_pixgrid * muse_pixgrid_create(muse_pixtable *aPixtable, cpl_propertylist *aHeader, cpl_size aXSize, cpl_size aYSize, cpl_size aZSize)
Convert selected rows of a pixel table into pixel grid, linking the grid points to entries (=rows) in...
Definition: muse_pixgrid.c:168
muse_pixtable ** muse_pixtable_extracted_get_slices(muse_pixtable *aPixtable)
Extract one pixel table per IFU and slice.
cpl_image * data
the data extension
Definition: muse_image.h:47
#define MUSE_HDR_PT_LHI
FITS header keyword contains the upper limit of the data in spectral direction.
int muse_pixtable_get_type(muse_pixtable *aPixtable)
Determine the type of pixel table.
const muse_cpltable_def muse_dataspectrum_def[]
void muse_utils_memory_dump(const char *aMarker)
Display the current memory usage of the given program.
Definition: muse_utils.c:2272
cpl_error_code muse_pixtable_dump(muse_pixtable *aPixtable, cpl_size aStart, cpl_size aCount, unsigned char aDisplayHeader)
Dump a MUSE pixel table to the screen, resolving the origin column.
double muse_pfits_get_pres_start(const cpl_propertylist *aHeaders)
find out the ambient pressure at start of exposure (in mbar)
Definition: muse_pfits.c:804
cpl_image * stat
the statistics extension
Definition: muse_image.h:65
cpl_error_code muse_wcs_projplane_from_celestial(cpl_propertylist *aHeader, double aRA, double aDEC, double *aX, double *aY)
Convert from celestial spherical coordinates to projection plane coordinates.
Definition: muse_wcs.c:1324
void muse_datacube_delete(muse_datacube *aCube)
Deallocate memory associated to a muse_datacube object.
cpl_propertylist * hgroup
the group FITS header
cpl_propertylist * header
the primary FITS header
The pixel grid.
Definition: muse_pixgrid.h:56
double muse_astro_parangle(const cpl_propertylist *aHeader)
Properly average parallactic angle values of one exposure.
Definition: muse_astro.c:440
Structure definition of MUSE three extension FITS file.
Definition: muse_image.h:41
cpl_table * table
The pixel table.
cpl_propertylist * header
the FITS header
Definition: muse_image.h:73
muse_resampling_crstats_type crtype
cpl_image * dq
the data quality extension
Definition: muse_image.h:57
cpl_table * muse_cpltable_new(const muse_cpltable_def *aDef, cpl_size aLength)
Create an empty table according to the specified definition.
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
Structure definition of MUSE pixel table.
muse_pixtable_wcs muse_pixtable_wcs_check(muse_pixtable *aPixtable)
Check the state of the world coordinate system of a pixel table.
#define MUSE_HDR_PT_YLO
FITS header keyword contains the lower limit of the data in y-direction.
static void muse_wcs_celestial_from_pixel_fast(muse_wcs *aWCS, double aX, double aY, double *aRA, double *aDEC)
Convert from pixel coordinates to celestial spherical coordinates.
Definition: muse_wcs.h:126
muse_resampling_crstats_type
Cosmic ray rejection statistics type.
muse_resampling_params * muse_resampling_params_new(muse_resampling_type aMethod)
Create the resampling parameters structure.
cpl_error_code muse_wcs_projplane_from_pixel(cpl_propertylist *aHeader, double aX, double aY, double *aXOut, double *aYOut)
Convert from pixel coordinates to projection plane coordinates.
Definition: muse_wcs.c:1372
cpl_imagelist * data
the cube containing the actual data values
Definition: muse_datacube.h:76
cpl_error_code muse_wcs_pixel_from_celestial(cpl_propertylist *aHeader, double aRA, double aDEC, double *aX, double *aY)
Convert from celestial spherical coordinates to pixel coordinates.
Definition: muse_wcs.c:1269
Structure definition of a Euro3D datacube.
Definition: muse_datacube.h:97
static const cpl_size * muse_pixgrid_get_rows(muse_pixgrid *aPixels, cpl_size aIndex)
Return a pointer to the rows stored in one pixel.
Definition: muse_pixgrid.h:150
muse_datacube * muse_resampling_cube(muse_pixtable *aPixtable, muse_resampling_params *aParams, muse_pixgrid **aPixgrid)
Resample a pixel table onto a regular grid structure representing a FITS NAXIS=3 datacube.
cpl_imagelist * dq
the optional cube containing the bad pixel status
Definition: muse_datacube.h:81
muse_resampling_dispersion_type tlambda
#define MUSE_HDR_PT_LLO
FITS header keyword contains the lower limit of the data in spectral direction.
cpl_table * dtable
the table containing the actual Euro3D data
cpl_image * muse_cplimage_concat_x(const cpl_image *aImage1, const cpl_image *aImage2)
Concatenate two images in x direction.
cpl_propertylist * hdata
the data FITS header
muse_image * muse_resampling_collapse_pixgrid(muse_pixtable *aPixtable, muse_pixgrid *aPixgrid, muse_datacube *aCube, cpl_table *aFilter, muse_resampling_params *aParams)
Integrate a pixel table / pixel grid along the wavelength direction.
void muse_pixgrid_delete(muse_pixgrid *aPixels)
Delete a pixgrid and remove its memory.
Definition: muse_pixgrid.c:398
cpl_table * gtable
the table containing the Euro3D groups
cpl_propertylist * header
the FITS header
Definition: muse_datacube.h:57
const muse_cpltable_def muse_euro3dcube_e3d_data_def[]
double muse_astro_posangle(const cpl_propertylist *aHeader)
Derive the position angle of an observation from information in a FITS header.
Definition: muse_astro.c:404
static void muse_wcs_projplane_from_pixel_fast(muse_wcs *aWCS, double aX, double aY, double *aXOut, double *aYOut)
Convert from pixel coordinates to projection plane coordinates.
Definition: muse_wcs.h:204
cpl_boolean iscelsph
Definition: muse_wcs.h:100
#define MUSE_HDR_PT_DAR_NAME
muse_resampling_type
Resampling types.
double muse_pfits_get_temp(const cpl_propertylist *aHeaders)
find out the ambient temperature (in degrees Celsius)
Definition: muse_pfits.c:768
muse_euro3dcube * muse_resampling_euro3d(muse_pixtable *aPixtable, muse_resampling_params *aParams)
Resample a pixel table onto a regular grid structure representing a Euro3D format file...
Resampling parameters.
void muse_resampling_params_delete(muse_resampling_params *aParams)
Delete a resampling parameters structure.
muse_resampling_type method
muse_image * muse_image_new(void)
Allocate memory for a new muse_image object.
Definition: muse_image.c:65
cpl_table * muse_resampling_spectrum(muse_pixtable *aPixtable, double aBinwidth)
Resample the selected pixels of a pixel table into a spectrum.
static cpl_size muse_pixgrid_get_index(muse_pixgrid *aPixels, cpl_size aX, cpl_size aY, cpl_size aZ, cpl_boolean aAllowOutside)
Get the grid index determined from all three coordinates.
Definition: muse_pixgrid.h:89
#define MUSE_HDR_PT_YHI
FITS header keyword contains the upper limit of the data in y-direction.
static cpl_size muse_pixgrid_get_count(muse_pixgrid *aPixels, cpl_size aIndex)
Return the number of rows stored in one pixel.
Definition: muse_pixgrid.h:129
const muse_cpltable_def muse_euro3dcube_e3d_grp_def[]
muse_ins_mode muse_pfits_get_mode(const cpl_propertylist *aHeaders)
find out the observation mode
Definition: muse_pfits.c:1003
double muse_astro_airmass(cpl_propertylist *aHeader)
Derive the effective airmass of an observation from information in a FITS header. ...
Definition: muse_astro.c:322
cpl_imagelist * stat
the cube containing the data variance
Definition: muse_datacube.h:86
cpl_propertylist * header
The FITS header.
#define MUSE_HDR_PT_XHI
FITS header keyword contains the upper limit of the data in x-ion.