00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #ifdef HAVE_CONFIG_H
00023 #include <config.h>
00024 #endif
00025
00026
00027
00028
00029 #include <cpl.h>
00030 #include <math.h>
00031 #include <string.h>
00032 #ifndef _OPENMP
00033 #define omp_get_max_threads() 1
00034 #else
00035 #include <omp.h>
00036 #endif
00037
00038 #include "muse_resampling.h"
00039 #include "muse_instrument.h"
00040
00041 #include "muse_astro.h"
00042 #include "muse_cplwrappers.h"
00043 #include "muse_dar.h"
00044 #include "muse_dfs.h"
00045 #include "muse_flux.h"
00046 #include "muse_pfits.h"
00047 #include "muse_pixgrid.h"
00048 #include "muse_pixtable.h"
00049 #include "muse_quality.h"
00050 #include "muse_utils.h"
00051 #include "muse_wcs.h"
00052 #include "muse_data_format_z.h"
00053
00054
00058
00059
00062
00063
00064
00065
00066
00079
00080 muse_resampling_params *
00081 muse_resampling_params_new(muse_resampling_type aMethod)
00082 {
00083 cpl_ensure(aMethod <= MUSE_RESAMPLE_NONE, CPL_ERROR_ILLEGAL_INPUT, NULL);
00084 muse_resampling_params *params = cpl_calloc(1, sizeof(muse_resampling_params));
00085 params->method = aMethod;
00086
00087
00088 params->ld = 1;
00089 params->pfx = 0.6;
00090 params->pfy = 0.6;
00091 params->pfl = 0.6;
00092 params->rc = 1.25;
00093
00094
00095 return params;
00096 }
00097
00098
00117
00118 cpl_error_code
00119 muse_resampling_params_set_wcs(muse_resampling_params *aParams,
00120 const cpl_propertylist *aWCS)
00121 {
00122 cpl_ensure_code(aParams, CPL_ERROR_NULL_INPUT);
00123 if (!aWCS) {
00124 aParams->tlambda = MUSE_RESAMPLING_DISP_AWAV;
00125 cpl_wcs_delete(aParams->wcs);
00126 aParams->wcs = NULL;
00127 return CPL_ERROR_NONE;
00128 }
00129
00130 aParams->tlambda = MUSE_RESAMPLING_DISP_AWAV;
00131 if (cpl_propertylist_has(aWCS, "CTYPE3")) {
00132 if (!strncmp(muse_pfits_get_ctype(aWCS, 3), "AWAV-LOG", 9)) {
00133 aParams->tlambda = MUSE_RESAMPLING_DISP_AWAV_LOG;
00134 } else if (!strncmp(muse_pfits_get_ctype(aWCS, 3), "WAVE", 5)) {
00135 aParams->tlambda = MUSE_RESAMPLING_DISP_WAVE;
00136 } else if (!strncmp(muse_pfits_get_ctype(aWCS, 3), "WAVE-LOG", 9)) {
00137 aParams->tlambda = MUSE_RESAMPLING_DISP_WAVE_LOG;
00138 } else {
00139
00140 }
00141 }
00142 cpl_errorstate state = cpl_errorstate_get();
00143 cpl_error_code rc = CPL_ERROR_NONE;
00144 aParams->wcs = cpl_wcs_new_from_propertylist(aWCS);
00145 if (!cpl_errorstate_is_equal(state)) {
00146 cpl_wcs_delete(aParams->wcs);
00147 aParams->wcs = NULL;
00148 rc = cpl_error_get_code();
00149 }
00150 #if 0
00151 printf("CDi_j matrix:\n");
00152 cpl_matrix_dump(cpl_wcs_get_cd(aParams->wcs), stdout);
00153 fflush(stdout);
00154 #endif
00155 return rc;
00156 }
00157
00158
00166
00167 void
00168 muse_resampling_params_delete(muse_resampling_params *aParams)
00169 {
00170 if (!aParams) {
00171 return;
00172 }
00173 cpl_wcs_delete(aParams->wcs);
00174 cpl_free(aParams);
00175 }
00176
00177
00188
00189 static inline double
00190 muse_resampling_weight_function_linear(double r)
00191 {
00192 return r == 0 ? FLT_MAX : 1. / r;
00193 }
00194
00195
00206
00207 static inline double
00208 muse_resampling_weight_function_quadratic(double r2)
00209 {
00210 return r2 == 0 ? FLT_MAX : 1. / r2;
00211 }
00212
00213
00226
00227 static inline double
00228 muse_resampling_weight_function_renka(double r, double r_c)
00229 {
00230 if (r == 0) {
00231 return FLT_MAX;
00232 } else if (r >= r_c) {
00233 return DBL_MIN;
00234 } else {
00235 double p = (r_c - r) / (r_c * r);
00236 return p*p;
00237 }
00238 }
00239
00240
00247
00248 static inline double
00249 muse_resampling_weight_function_sinc(double r)
00250 {
00251 return fabs(r) < DBL_EPSILON ? 1. : sin(CPL_MATH_PI * r) / (CPL_MATH_PI * r);
00252 }
00253
00254
00264
00265 static inline double
00266 muse_resampling_weight_function_lanczos(double dx, double dy, double dz, unsigned int n)
00267 {
00268 return (fabs(dx) >= n || fabs(dy) >= n || fabs(dz) > n) ? 0.
00269 : muse_resampling_weight_function_sinc(dx) * muse_resampling_weight_function_sinc(dx / n)
00270 * muse_resampling_weight_function_sinc(dy) * muse_resampling_weight_function_sinc(dy / n)
00271 * muse_resampling_weight_function_sinc(dz) * muse_resampling_weight_function_sinc(dz / n);
00272 }
00273
00274
00293
00294 static inline double
00295 muse_resampling_weight_function_drizzle(double xin, double yin, double zin,
00296 double xout, double yout, double zout,
00297 double dx, double dy, double dz)
00298 {
00299
00300
00301
00302 double x = (dx + xout / 2.) <= xin / 2. ? xout : (xin + xout) / 2. - dx,
00303 y = (dy + yout / 2.) <= yin / 2. ? yout : (yin + yout) / 2. - dy,
00304 z = (dz + zout / 2.) <= zin / 2. ? zout : (zin + zout) / 2. - dz;
00305
00306
00307 if (x <= 0 || y <= 0 || z <= 0) {
00308 return 0.;
00309 }
00310
00311
00312 return (x > xin ? xin : x) * (y > yin ? yin : y) * (z > zin ? zin : z)
00313 / (xin * yin * zin);
00314 }
00315
00316
00332
00333 static cpl_error_code
00334 muse_resampling_cube_nearest(muse_datacube *aCube, muse_pixtable *aPixtable,
00335 muse_pixgrid *aPixgrid)
00336 {
00337 cpl_ensure_code(aCube && aPixtable && aPixgrid, CPL_ERROR_NULL_INPUT);
00338 double ptxoff = 0.,
00339 ptyoff = 0.,
00340 crval3 = muse_pfits_get_crval(aCube->header, 3),
00341 crpix3 = muse_pfits_get_crpix(aCube->header, 3),
00342 cd33 = muse_pfits_get_cd(aCube->header, 3, 3);
00343 const char *ctype3 = muse_pfits_get_ctype(aCube->header, 3);
00344 muse_wcs *wcs = muse_wcs_new(aCube->header);
00345 wcs->iscelsph = muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_CELSPH;
00346 cpl_boolean loglambda = ctype3 && (!strncmp(ctype3, "AWAV-LOG", 9) ||
00347 !strncmp(ctype3, "WAVE-LOG", 9));
00348 float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
00349 *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
00350 *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
00351 *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
00352 *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
00353 int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
00354
00355
00356
00357
00358
00359 double xnorm = 1., ynorm = 1., znorm = 1. / kMuseSpectralSamplingA;
00360 if (wcs->iscelsph) {
00361 muse_wcs_get_scales(aPixtable->header, &xnorm, &ynorm);
00362 xnorm = 1. / xnorm;
00363 ynorm = 1. / ynorm;
00364
00365 ptxoff = muse_pfits_get_crval(aPixtable->header, 1);
00366 ptyoff = muse_pfits_get_crval(aPixtable->header, 2);
00367 }
00368
00369 #ifdef ESO_ENABLE_DEBUG
00370 int debug = 0;
00371 if (getenv("MUSE_DEBUG_NEAREST")) {
00372 debug = atoi(getenv("MUSE_DEBUG_NEAREST"));
00373 }
00374 #endif
00375
00376 int l;
00377 #ifdef ESO_ENABLE_DEBUG
00378 #pragma omp parallel for default(none) \
00379 shared(aCube, aPixgrid, cd33, crpix3, crval3, data, debug, dq, lbda, \
00380 loglambda, ptxoff, ptyoff, stat, stdout, wcs, xnorm, xpos, \
00381 ynorm, ypos, znorm)
00382 #else
00383 #pragma omp parallel for default(none) \
00384 shared(aCube, aPixgrid, cd33, crpix3, crval3, data, dq, lbda, \
00385 loglambda, ptxoff, ptyoff, stat, stdout, wcs, xnorm, xpos, \
00386 ynorm, ypos, znorm)
00387 #endif
00388 for (l = 0; l < aPixgrid->size_z; l++) {
00389 float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->data, l)),
00390 *pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->stat, l));
00391 int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->dq, l));
00392
00393 double lambda = (l + 1. - crpix3) * cd33 + crval3;
00394 if (loglambda) {
00395 lambda = crval3 * exp((l + 1. - crpix3) * cd33 / crval3);
00396 }
00397
00398 int i;
00399 for (i = 0; i < aPixgrid->size_x; i++) {
00400 int j;
00401 for (j = 0; j < aPixgrid->size_y; j++) {
00402 cpl_size idx = muse_pixgrid_get_index(aPixgrid, i, j, l, CPL_FALSE),
00403 n_rows = muse_pixgrid_get_count(aPixgrid, idx);
00404 const cpl_size *rows = muse_pixgrid_get_rows(aPixgrid, idx);
00405
00406
00407 double x, y;
00408 if (wcs->iscelsph) {
00409 muse_wcs_celestial_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
00410 } else {
00411 muse_wcs_projplane_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
00412 }
00413
00414 if (n_rows == 1) {
00415
00416 pdata[i + j * aPixgrid->size_x] = data[rows[0]];
00417 pstat[i + j * aPixgrid->size_x] = stat[rows[0]];
00418 pdq[i + j * aPixgrid->size_x] = dq[rows[0]];
00419
00420 #ifdef ESO_ENABLE_DEBUG
00421 if (debug) {
00422 printf("only: %f,%f\n", data[rows[0]], stat[rows[0]]);
00423 fflush(stdout);
00424 }
00425 #endif
00426 } else if (n_rows >= 2) {
00427
00428 cpl_size n, nbest = -1;
00429 double dbest = FLT_MAX;
00430 for (n = 0; n < n_rows; n++) {
00431
00432 double dx = fabs(x - xpos[rows[n]] + ptxoff) * xnorm,
00433 dy = fabs(y - ypos[rows[n]] + ptyoff) * ynorm,
00434 dlambda = fabs(lambda - lbda[rows[n]]) * znorm,
00435 dthis = sqrt(dx*dx + dy*dy + dlambda*dlambda);
00436 if (wcs->iscelsph) {
00437
00438
00439 dx *= cos(y * CPL_MATH_RAD_DEG);
00440 }
00441 #ifdef ESO_ENABLE_DEBUG
00442 if (debug) {
00443 printf("%f %f %f, %f %f %f, d: %f %f %f -> %f best: %f (%f,%f)\n",
00444 x, y, lambda, xpos[rows[n]] + ptxoff, ypos[rows[n]] + ptyoff,
00445 lbda[rows[n]], dx, dy, dlambda, dthis, dbest, data[rows[n]],
00446 stat[rows[n]]);
00447 }
00448 #endif
00449 if (dthis < dbest) {
00450 nbest = n;
00451 dbest = dthis;
00452 }
00453 }
00454 pdata[i + j * aPixgrid->size_x] = data[rows[nbest]];
00455 pstat[i + j * aPixgrid->size_x] = stat[rows[nbest]];
00456 pdq[i + j * aPixgrid->size_x] = dq[rows[nbest]];
00457 #ifdef ESO_ENABLE_DEBUG
00458 if (debug) {
00459 printf("closest: %f,%f\n", data[rows[nbest]], stat[rows[nbest]]);
00460 fflush(stdout);
00461 }
00462 #endif
00463 } else {
00464
00465 pdq[i + j * aPixgrid->size_x] = EURO3D_MISSDATA;
00466 }
00467 }
00468 }
00469 }
00470 cpl_free(wcs);
00471
00472 return CPL_ERROR_NONE;
00473 }
00474
00475
00496
00497 static cpl_error_code
00498 muse_resampling_cube_weighted(muse_datacube *aCube, muse_pixtable *aPixtable,
00499 muse_pixgrid *aPixgrid,
00500 muse_resampling_params *aParams)
00501 {
00502 cpl_ensure_code(aCube && aPixtable && aPixgrid && aParams,
00503 CPL_ERROR_NULL_INPUT);
00504 double ptxoff = 0.,
00505 ptyoff = 0.,
00506 crval3 = muse_pfits_get_crval(aCube->header, 3),
00507 crpix3 = muse_pfits_get_crpix(aCube->header, 3),
00508 cd33 = muse_pfits_get_cd(aCube->header, 3, 3);
00509 const char *ctype3 = muse_pfits_get_ctype(aCube->header, 3);
00510 muse_wcs *wcs = muse_wcs_new(aCube->header);
00511 wcs->iscelsph = muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_CELSPH;
00512 cpl_boolean loglambda = ctype3 && (!strncmp(ctype3, "AWAV-LOG", 9) ||
00513 !strncmp(ctype3, "WAVE-LOG", 9));
00514 cpl_errorstate prestate = cpl_errorstate_get();
00515 float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
00516 *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
00517 *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
00518 *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
00519 *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT),
00520 *wght = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_WEIGHT);
00521 if (!cpl_errorstate_is_equal(prestate)) {
00522 cpl_errorstate_set(prestate);
00523 }
00524 int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
00525
00526
00527
00528
00529
00530 double xnorm = 1., ynorm = 1., znorm = 1. / kMuseSpectralSamplingA;
00531 if (wcs->iscelsph) {
00532 muse_wcs_get_scales(aPixtable->header, &xnorm, &ynorm);
00533 xnorm = 1. / xnorm;
00534 ynorm = 1. / ynorm;
00535
00536 ptxoff = muse_pfits_get_crval(aPixtable->header, 1);
00537 ptyoff = muse_pfits_get_crval(aPixtable->header, 2);
00538 }
00539
00540
00541 double zoutefac = exp(1.5 * cd33 / crval3) - exp(0.5 * cd33 / crval3);
00542
00543 double renka_rc = aParams->rc
00544 * sqrt((wcs->cd11*xnorm)*(wcs->cd11*xnorm) + (wcs->cd22*ynorm)*(wcs->cd22*ynorm)
00545 + (cd33*znorm)*(cd33*znorm));
00546
00547 int ld = aParams->ld;
00548 if (ld <= 0) {
00549 ld = 1;
00550 cpl_msg_info(__func__, "Overriding loop distance ld=%d", ld);
00551 }
00552
00553
00554
00555 double xsz = aParams->pfx / xnorm,
00556 ysz = aParams->pfy / ynorm,
00557 zsz = aParams->pfl / znorm,
00558 xout = fabs(wcs->cd11), yout = fabs(wcs->cd22), zout = fabs(cd33);
00559
00560 #ifdef ESO_ENABLE_DEBUG
00561 int debug = 0, debugx = 0, debugy = 0, debugz = 0;
00562 if (getenv("MUSE_DEBUG_WEIGHTED")) {
00563 debug = atoi(getenv("MUSE_DEBUG_WEIGHTED"));
00564 }
00565 if (debug & 2) {
00566 if (getenv("MUSE_DEBUG_WEIGHTED_X")) {
00567 debugx = atoi(getenv("MUSE_DEBUG_WEIGHTED_X"));
00568 if (debugx < 1 || debugx > aPixgrid->size_x) {
00569 debugx = 0;
00570 }
00571 }
00572 if (getenv("MUSE_DEBUG_WEIGHTED_Y")) {
00573 debugy = atoi(getenv("MUSE_DEBUG_WEIGHTED_Y"));
00574 if (debugy < 1 || debugy > aPixgrid->size_y) {
00575 debugy = 0;
00576 }
00577 }
00578 if (getenv("MUSE_DEBUG_WEIGHTED_Z")) {
00579 debugz = atoi(getenv("MUSE_DEBUG_WEIGHTED_Z"));
00580 if (debugz < 1 || debugz > aPixgrid->size_z) {
00581 debugz = 0;
00582 }
00583 }
00584 }
00585 if (debug & 8) {
00586 printf("parameters:\n cd=%e %e %e\n"
00587 " corrected crpix3=%e\n norm=%e %e %e\n",
00588 wcs->cd11, wcs->cd22, cd33, crpix3, xnorm, ynorm, znorm);
00589 if (aParams->method == MUSE_RESAMPLE_WEIGHTED_DRIZZLE) {
00590 printf(" drop sizes %e %e %e (pixfrac %f,%f,%f)\n"
00591 " output sizes %e %e %e\n",
00592 xsz, ysz, zsz, aParams->pfx, aParams->pfy, aParams->pfl,
00593 xout, yout, zout);
00594 } else {
00595 printf(" resulting renka_rc: %e %e %e / %e %e %e --> %e\n",
00596 pow(wcs->cd11, 2), pow(wcs->cd22, 2), cd33*cd33,
00597 pow(wcs->cd11*xnorm, 2), pow(wcs->cd22*ynorm, 2),
00598 pow(cd33*znorm, 2), renka_rc);
00599 }
00600 fflush(stdout);
00601 }
00602 #endif
00603 cpl_imagelist *wcube = NULL;
00604 if (getenv("MUSE_DEBUG_WEIGHT_CUBE")) {
00605 cpl_msg_debug(__func__, "Weighted resampling: creating weight cube");
00606 wcube = cpl_imagelist_new();
00607 int i;
00608 for (i = 0; i < aPixgrid->size_z; i++) {
00609 cpl_image *image = cpl_image_new(aPixgrid->size_x, aPixgrid->size_y,
00610 CPL_TYPE_FLOAT);
00611 cpl_imagelist_set(wcube, image, i);
00612 }
00613 }
00614
00615 if (getenv("MUSE_DEBUG_WEIGHTED_GRID")) {
00616 char *fn = getenv("MUSE_DEBUG_WEIGHTED_GRID");
00617 FILE *grid = fopen(fn, "w");
00618 if (grid) {
00619 cpl_msg_info(__func__, "Writing grid to \"%s\"", fn);
00620 fprintf(grid, "# i j l x y lambda\n");
00621 int l;
00622 for (l = 0; l < aPixgrid->size_z; l++) {
00623 double lambda = (l + 1. - crpix3) * cd33 + crval3;
00624 int i;
00625 for (i = 0; i < aPixgrid->size_x; i++) {
00626 int j;
00627 for (j = 0; j < aPixgrid->size_y; j++) {
00628
00629 double x, y;
00630 if (wcs->iscelsph) {
00631 muse_wcs_celestial_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
00632 } else {
00633 muse_wcs_projplane_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
00634 }
00635 fprintf(grid, "%03d %03d %04d %.10f %.10f %8.3f\n", i+1, j+1, l+1,
00636 x, y, lambda);
00637 }
00638 }
00639 }
00640 fclose(grid);
00641 } else {
00642 cpl_msg_warning(__func__, "Writing grid to \"%s\" failed!", fn);
00643 }
00644 }
00645
00646 int l;
00647 #ifdef ESO_ENABLE_DEBUG
00648 #pragma omp parallel for default(none) \
00649 shared(aCube, aParams, aPixgrid, aPixtable, cd33, crpix3, crval3, \
00650 data, debug, debugx, debugy, debugz, dq, lbda, ld, loglambda, \
00651 ptxoff, ptyoff, renka_rc, stat, stdout, wcs, wcube, wght, \
00652 xnorm, xout, xpos, xsz, ynorm, yout, ypos, ysz, znorm, zout, \
00653 zoutefac, zsz)
00654 #else
00655 #pragma omp parallel for default(none) \
00656 shared(aCube, aParams, aPixgrid, aPixtable, cd33, crpix3, crval3, \
00657 data, dq, lbda, ld, loglambda, ptxoff, ptyoff, renka_rc, stat,\
00658 stdout, wcs, wcube, wght, xnorm, xout, xpos, xsz, ynorm, yout,\
00659 ypos, ysz, znorm, zout, zoutefac, zsz)
00660 #endif
00661 for (l = 0; l < aPixgrid->size_z; l++) {
00662 float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->data, l)),
00663 *pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->stat, l));
00664 int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->dq, l));
00665
00666 double lambda = (l + 1. - crpix3) * cd33 + crval3;
00667 double zout2 = zout;
00668 if (loglambda) {
00669 lambda = crval3 * exp((l + 1. - crpix3) * cd33 / crval3);
00670 zout2 = crval3 * exp((l - crpix3) * cd33 / crval3) * zoutefac;
00671 }
00672 float *pwcube = NULL;
00673 if (wcube) {
00674 pwcube = cpl_image_get_data_float(cpl_imagelist_get(wcube, l));
00675 }
00676
00677 int i;
00678 for (i = 0; i < aPixgrid->size_x; i++) {
00679 int j;
00680 for (j = 0; j < aPixgrid->size_y; j++) {
00681
00682 double x, y;
00683 if (wcs->iscelsph) {
00684 muse_wcs_celestial_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
00685 } else {
00686 muse_wcs_projplane_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
00687 }
00688 double sumdata = 0, sumstat = 0, sumweight = 0;
00689 int npoints = 0;
00690 #ifdef ESO_ENABLE_DEBUG
00691 cpl_size *pointlist = NULL;
00692 double *pointweights = NULL;
00693 int npointlist = 0;
00694 if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
00695 pointlist = cpl_calloc(100, sizeof(cpl_size));
00696 pointweights = cpl_malloc(100 * sizeof(double));
00697 npointlist = 100;
00698 }
00699 #endif
00700
00701 #ifdef ESO_ENABLE_DEBUG
00702 if (debug & 1) {
00703 printf("cell %d %d %d:\n", i, j, l);
00704 }
00705 #endif
00706
00707 int i2;
00708 for (i2 = i - ld; i2 <= i + ld; i2++) {
00709 int j2;
00710 for (j2 = j - ld; j2 <= j + ld; j2++) {
00711 int l2;
00712 for (l2 = l - ld; l2 <= l + ld; l2++) {
00713 cpl_size idx2 = muse_pixgrid_get_index(aPixgrid, i2, j2, l2, CPL_FALSE),
00714 n_rows2 = muse_pixgrid_get_count(aPixgrid, idx2);
00715 #ifdef ESO_ENABLE_DEBUG
00716 if (debug & 8 && n_rows2 < 1) {
00717 printf("%d %d %d / %d %d %d (%"CPL_SIZE_FORMAT"): no rows!\n",
00718 i+1, j+1, l+1, i2+1, j2+1, l2+1, idx2);
00719 }
00720 #endif
00721 const cpl_size *rows2 = muse_pixgrid_get_rows(aPixgrid, idx2);
00722 cpl_size n;
00723 for (n = 0; n < n_rows2; n++) {
00724 if (dq[rows2[n]]) {
00725 #ifdef ESO_ENABLE_DEBUG
00726 if (debug & 8) {
00727 printf("%d %d %d / %d %d %d (%"CPL_SIZE_FORMAT", "
00728 "%"CPL_SIZE_FORMAT"): bad!\n",
00729 i+1, j+1, l+1, i2+1, j2+1, l2+1, idx2, n);
00730 fflush(stdout);
00731 }
00732 #endif
00733 continue;
00734 }
00735
00736 double dx = fabs(x - (xpos[rows2[n]] + ptxoff)),
00737 dy = fabs(y - (ypos[rows2[n]] + ptyoff)),
00738 dlambda = fabs(lambda - lbda[rows2[n]]),
00739 r2 = 0;
00740 if (wcs->iscelsph) {
00741
00742
00743
00744
00745
00746 dx *= cos(y * CPL_MATH_RAD_DEG);
00747 }
00748 if (aParams->method != MUSE_RESAMPLE_WEIGHTED_DRIZZLE) {
00749 dx *= xnorm;
00750 dy *= ynorm;
00751 dlambda *= znorm;
00752 r2 = dx*dx + dy*dy + dlambda*dlambda;
00753 }
00754 double weight = 0.;
00755 if (aParams->method == MUSE_RESAMPLE_WEIGHTED_RENKA) {
00756 weight = muse_resampling_weight_function_renka(sqrt(r2), renka_rc);
00757 } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_DRIZZLE) {
00758 weight = muse_resampling_weight_function_drizzle(xsz, ysz, zsz,
00759 xout, yout, zout2,
00760 dx, dy, dlambda);
00761 } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_LINEAR) {
00762 weight = muse_resampling_weight_function_linear(sqrt(r2));
00763 } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_QUADRATIC) {
00764 weight = muse_resampling_weight_function_quadratic(r2);
00765 } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_LANCZOS) {
00766 weight = muse_resampling_weight_function_lanczos(dx, dy, dlambda, ld);
00767 }
00768
00769 if (wght) {
00770
00771 weight *= wght[rows2[n]];
00772 }
00773 #ifdef ESO_ENABLE_DEBUG
00774 if (debug & 8) {
00775 printf("%d %d %d / %d %d %d (%"CPL_SIZE_FORMAT", %"CPL_SIZE_FORMAT"):"
00776 " x %e %e %e %e y %e %e %e %e l %e %e %e %e --> %e %e %e\n",
00777 i+1, j+1, l+1, i2+1, j2+1, l2+1, idx2, n,
00778 x, xpos[rows2[n]]+ptxoff, fabs(x - (xpos[rows2[n]]+ptxoff)), dx,
00779 y, ypos[rows2[n]]+ptyoff, fabs(y - (ypos[rows2[n]]+ptyoff)), dy,
00780 lambda, lbda[rows2[n]], fabs(lambda - lbda[rows2[n]]), dlambda,
00781 r2, sqrt(r2), weight);
00782 fflush(stdout);
00783 }
00784 #endif
00785 sumweight += weight;
00786 sumdata += data[rows2[n]] * weight;
00787 sumstat += stat[rows2[n]] * weight*weight;
00788 npoints++;
00789 #ifdef ESO_ENABLE_DEBUG
00790 if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
00791 if (npoints > npointlist) {
00792 pointlist = cpl_realloc(pointlist,
00793 (npointlist + 100) * sizeof(cpl_size));
00794 memset(pointlist + npointlist, 0, 100 * sizeof(cpl_size));
00795 pointweights = cpl_realloc(pointweights,
00796 (npointlist + 100) * sizeof(double));
00797 npointlist += 100;
00798 }
00799
00800 pointlist[npoints-1] = rows2[n] + 1;
00801 pointweights[npoints-1] = weight;
00802 }
00803
00804 if (debug & 4) {
00805 cpl_size idx = muse_pixgrid_get_index(aPixgrid, i, j, l, CPL_FALSE),
00806 count = muse_pixgrid_get_count(aPixgrid, idx);
00807 if (count) {
00808 printf(" pixel %4d,%4d,%4d (%8"CPL_SIZE_FORMAT"): "
00809 "%2"CPL_SIZE_FORMAT" %2"CPL_SIZE_FORMAT" %f %f %f, "
00810 " %e -> %e ==> %e %e (%d)\n", i+1, j+1, l+1, idx,
00811 n, count, dx, dy, dlambda,
00812 data[muse_pixgrid_get_rows(aPixgrid, idx)[n]],
00813 weight, sumweight, sumdata, npoints);
00814 }
00815 }
00816 #endif
00817 }
00818 }
00819 }
00820 }
00821
00822 #ifdef ESO_ENABLE_DEBUG
00823 if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
00824 printf("cell center (%d, %d, %d): %14.7e %14.7e %9.3f\npixelnumber "
00825 "weight ", debugx, debugy, debugz, x, y, lambda);
00826 muse_pixtable_dump(aPixtable, 0, 0, 2);
00827 int m = -1;
00828 while (++m < npointlist && pointlist[m] != 0) {
00829
00830 printf("%12"CPL_SIZE_FORMAT" %8.5f ", pointlist[m] - 1,
00831 pointweights[m]);
00832 muse_pixtable_dump(aPixtable, pointlist[m] - 1, 1, 0);
00833 }
00834 printf("sums: %g %g %g --> %g %g\n", sumdata, sumstat, sumweight,
00835 sumdata / sumweight, sumstat / pow(sumweight, 2));
00836 cpl_free(pointlist);
00837 cpl_free(pointweights);
00838 }
00839
00840 if (debug & 1 && npoints) {
00841 printf(" sumdata: %e %e (%d)", sumdata, sumweight, npoints);
00842 }
00843 #endif
00844
00845
00846
00847
00848 if (!npoints || !isnormal(sumweight)) {
00849 pdq[i + j * aPixgrid->size_x] = EURO3D_MISSDATA;
00850 #ifdef ESO_ENABLE_DEBUG
00851 if (debug & 1) {
00852 printf(" -> no points or weird weight\n");
00853 }
00854 #endif
00855 continue;
00856 }
00857
00858
00859 sumdata /= sumweight;
00860 sumstat /= sumweight*sumweight;
00861 #ifdef ESO_ENABLE_DEBUG
00862 if (debug & 1) {
00863 printf(" -> %e (%e)\n", sumdata, sumstat);
00864 }
00865 if (debug) {
00866 fflush(stdout);
00867 }
00868 #endif
00869 pdata[i + j * aPixgrid->size_x] = sumdata;
00870 pstat[i + j * aPixgrid->size_x] = sumstat;
00871 pdq[i + j * aPixgrid->size_x] = EURO3D_GOODPIXEL;
00872 if (pwcube) {
00873 pwcube[i + j * aPixgrid->size_x] = sumweight;
00874 }
00875 }
00876 }
00877 }
00878 cpl_free(wcs);
00879
00880 if (wcube) {
00881 const char *fn = getenv("MUSE_DEBUG_WEIGHT_CUBE");
00882 cpl_error_code rc = cpl_imagelist_save(wcube, fn, CPL_TYPE_UNSPECIFIED,
00883 NULL, CPL_IO_CREATE);
00884 if (rc != CPL_ERROR_NONE) {
00885 cpl_msg_warning(__func__, "Failure to save weight cube as \"%s\": %s", fn,
00886 cpl_error_get_message());
00887 } else {
00888 cpl_msg_info(__func__, "Saved weight cube as \"%s\"", fn);
00889 }
00890 cpl_imagelist_delete(wcube);
00891 }
00892
00893 return CPL_ERROR_NONE;
00894 }
00895
00896
00906
00907 static cpl_error_code
00908 muse_resampling_check_deltas(muse_pixtable *aPixtable,
00909 muse_resampling_params *aParams)
00910 {
00911 cpl_ensure_code(aPixtable && aParams, CPL_ERROR_NULL_INPUT);
00912 const char func[] = "muse_resampling_cube";
00913
00914
00915 if (aParams->dlambda == 0.0) {
00916 aParams->dlambda = kMuseSpectralSamplingA;
00917 if (aParams->tlambda == MUSE_RESAMPLING_DISP_AWAV_LOG ||
00918 aParams->tlambda == MUSE_RESAMPLING_DISP_WAVE_LOG) {
00919 aParams->dlambda /= 1.6;
00920 }
00921 }
00922
00923 if (muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_PIXEL) {
00924
00925 if (aParams->dx == 0.0) {
00926 aParams->dx = 1.0;
00927 }
00928 if (aParams->dy == 0.0) {
00929 aParams->dy = 1.0;
00930 }
00931 cpl_msg_debug(func, "steps from parameters: %.2f pix, %.2f pix, %.3f Angstrom",
00932 aParams->dx, aParams->dy, aParams->dlambda);
00933 return CPL_ERROR_NONE;
00934 }
00935 if (aParams->dx == 0.0) {
00936 aParams->dx = kMuseSpaxelSizeX_WFM / 3600.;
00937 if (muse_pfits_get_mode(aPixtable->header) == MUSE_MODE_NFM_AO_N) {
00938 aParams->dx = kMuseSpaxelSizeX_NFM / 3600.;
00939 }
00940 } else {
00941
00942 aParams->dx /= 3600.;
00943 }
00944 if (aParams->dy == 0.0) {
00945 aParams->dy = kMuseSpaxelSizeY_WFM / 3600.;
00946 if (muse_pfits_get_mode(aPixtable->header) == MUSE_MODE_NFM_AO_N) {
00947 aParams->dy = kMuseSpaxelSizeY_NFM / 3600.;
00948 }
00949 } else {
00950
00951 aParams->dy /= 3600.;
00952 }
00953 cpl_msg_debug(func, "steps from parameters: %f arcsec, %f arcsec, %.3f Angstrom",
00954 aParams->dx * 3600., aParams->dy * 3600., aParams->dlambda);
00955 return CPL_ERROR_NONE;
00956 }
00957
00958
00972
00973 static cpl_error_code
00974 muse_resampling_compute_size(muse_pixtable *aPixtable,
00975 muse_resampling_params *aParams,
00976 int *aX, int *aY, int *aZ)
00977 {
00978 cpl_ensure_code(aPixtable && aParams && aX && aY && aZ, CPL_ERROR_NULL_INPUT);
00979 const char func[] = "muse_resampling_cube";
00980 double x1, y1, x2, y2;
00981 float xmin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XLO),
00982 xmax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XHI),
00983 ymin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YLO),
00984 ymax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YHI);
00985 muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
00986 if (wcstype == MUSE_PIXTABLE_WCS_CELSPH) {
00987 muse_wcs_projplane_from_celestial(aPixtable->header, xmin, ymin, &x1, &y1);
00988 muse_wcs_projplane_from_celestial(aPixtable->header, xmax, ymax, &x2, &y2);
00989 } else {
00990 muse_wcs_projplane_from_pixel(aPixtable->header, xmin, ymin, &x1, &y1);
00991 muse_wcs_projplane_from_pixel(aPixtable->header, xmax, ymax, &x2, &y2);
00992 }
00993 *aX = lround(fabs(x2 - x1) / aParams->dx) + 1;
00994 *aY = lround(fabs(y2 - y1) / aParams->dy) + 1;
00995 float lmin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LLO),
00996 lmax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LHI);
00997 *aZ = (int)ceil((lmax - lmin) / aParams->dlambda) + 1;
00998 if (aParams->tlambda == MUSE_RESAMPLING_DISP_AWAV_LOG ||
00999 aParams->tlambda == MUSE_RESAMPLING_DISP_WAVE_LOG) {
01000 *aZ = (int)ceil(lmin / aParams->dlambda * log(lmax / lmin)) + 1;
01001 }
01002 cpl_msg_info(func, "Output cube size %d x %d x %d (fit to data)",
01003 *aX, *aY, *aZ);
01004 return CPL_ERROR_NONE;
01005 }
01006
01007
01017
01018 static void
01019 muse_resampling_override_size_int(int *aV, const char *aKeyword, int aValue)
01020 {
01021 if (aValue <= 0 || !aV) {
01022 return;
01023 }
01024 const char func[] = "muse_resampling_cube";
01025 cpl_msg_info(func, "Overriding %s=%d (was %d)", aKeyword, aValue, *aV);
01026 *aV = aValue;
01027 }
01028
01029
01043
01044 static cpl_error_code
01045 muse_resampling_override_size(int *aX, int *aY, int *aZ,
01046 muse_resampling_params *aParams)
01047 {
01048 cpl_ensure_code(aX && aY && aZ && aParams, CPL_ERROR_NULL_INPUT);
01049 if (!aParams->wcs) {
01050 return CPL_ERROR_NONE;
01051 }
01052 const char func[] = "muse_resampling_cube";
01053
01054 const cpl_array *dims = cpl_wcs_get_image_dims(aParams->wcs);
01055 if (!dims) {
01056 cpl_msg_debug(func, "No dimensions to override were specified");
01057 return CPL_ERROR_NONE;
01058 }
01059 muse_resampling_override_size_int(aX, "NAXIS1", cpl_array_get_int(dims, 0, NULL));
01060 muse_resampling_override_size_int(aY, "NAXIS2", cpl_array_get_int(dims, 1, NULL));
01061 if (cpl_wcs_get_image_naxis(aParams->wcs) >= 3) {
01062 muse_resampling_override_size_int(aZ, "NAXIS3",
01063 cpl_array_get_int(dims, 2, NULL));
01064 }
01065 return CPL_ERROR_NONE;
01066 }
01067
01068
01079
01080 static void
01081 muse_resampling_override_wcs_double(muse_datacube *aCube, const char *aKeyword,
01082 double aValue, int aError)
01083 {
01084 if (aError || !aCube) {
01085 cpl_msg_debug("double", "%s=%#g (%d)", aKeyword, aValue, aError);
01086 return;
01087 }
01088 const char func[] = "muse_resampling_cube";
01089 double old = cpl_propertylist_has(aCube->header, aKeyword)
01090 ? cpl_propertylist_get_double(aCube->header, aKeyword) : NAN;
01091 cpl_msg_info(func, "Overriding %s=%#g (was %#g)", aKeyword, aValue, old);
01092 cpl_propertylist_update_double(aCube->header, aKeyword, aValue);
01093
01094
01095 cpl_propertylist_update_bool(aCube->header, "MUSE_RESAMPLING_WCS_OVERRIDDEN",
01096 CPL_TRUE);
01097 }
01098
01099
01111
01112 static cpl_error_code
01113 muse_resampling_override_wcs(muse_datacube *aCube,
01114 muse_resampling_params *aParams)
01115 {
01116 cpl_ensure_code(aCube && aCube->header && aParams, CPL_ERROR_NULL_INPUT);
01117 if (!aParams->wcs) {
01118 return CPL_ERROR_NONE;
01119 }
01120 const char func[] = "muse_resampling_cube";
01121 const cpl_array *crval = cpl_wcs_get_crval(aParams->wcs),
01122 *crpix = cpl_wcs_get_crpix(aParams->wcs);
01123 const cpl_matrix *cd = cpl_wcs_get_cd(aParams->wcs);
01124 int err = 0;
01125
01126 if (crval) {
01127 muse_resampling_override_wcs_double(aCube, "CRVAL1", cpl_array_get_double(crval, 0, &err), err);
01128 muse_resampling_override_wcs_double(aCube, "CRVAL2", cpl_array_get_double(crval, 1, &err), err);
01129 } else {
01130 cpl_msg_debug(func, "No CRVALj to override were specified");
01131 }
01132 if (crpix) {
01133 muse_resampling_override_wcs_double(aCube, "CRPIX1", cpl_array_get_double(crpix, 0, &err), err);
01134 muse_resampling_override_wcs_double(aCube, "CRPIX2", cpl_array_get_double(crpix, 1, &err), err);
01135 } else {
01136 cpl_msg_debug(func, "No CRPIXi to override were specified");
01137 }
01138 if (cd) {
01139 int naxes = cpl_matrix_get_ncol(cd);
01140 if (cpl_matrix_get_determinant(cd) == 0.) {
01141 cpl_msg_warning(func, "det(CDi_j) = 0, not overriding CDi_j!");
01142 cd = NULL;
01143 } else if (naxes > 2 && (cpl_matrix_get(cd, 0, 2) != 0. ||
01144 cpl_matrix_get(cd, 2, 0) != 0. ||
01145 cpl_matrix_get(cd, 2, 1) != 0. ||
01146 cpl_matrix_get(cd, 1, 2) != 0.)) {
01147 cpl_msg_warning(func, "Axis 3 (dispersion) is not cleanly separated from "
01148 "axes 1 and 2, not overriding CDi_j!");
01149 cd = NULL;
01150 } else {
01151 cpl_errorstate es = cpl_errorstate_get();
01152 muse_resampling_override_wcs_double(aCube, "CD1_1", cpl_matrix_get(cd, 0, 0),
01153 !cpl_errorstate_is_equal(es));
01154 es = cpl_errorstate_get();
01155 muse_resampling_override_wcs_double(aCube, "CD1_2", cpl_matrix_get(cd, 0, 1),
01156 !cpl_errorstate_is_equal(es));
01157 es = cpl_errorstate_get();
01158 muse_resampling_override_wcs_double(aCube, "CD2_1", cpl_matrix_get(cd, 1, 0),
01159 !cpl_errorstate_is_equal(es));
01160 es = cpl_errorstate_get();
01161 muse_resampling_override_wcs_double(aCube, "CD2_2", cpl_matrix_get(cd, 1, 1),
01162 !cpl_errorstate_is_equal(es));
01163 }
01164 } else {
01165 cpl_msg_debug(func, "No CDi_j to override were specified");
01166 }
01167
01168
01169
01170 if (cpl_array_get_size(crval) > 2) {
01171 if (crval) {
01172 muse_resampling_override_wcs_double(aCube, "CRVAL3",
01173 cpl_array_get_double(crval, 2, &err) * 1e10, err);
01174 }
01175 if (crpix) {
01176 muse_resampling_override_wcs_double(aCube, "CRPIX3", cpl_array_get_double(crpix, 2, &err), err);
01177 }
01178 if (cd) {
01179 cpl_errorstate es = cpl_errorstate_get();
01180 muse_resampling_override_wcs_double(aCube, "CD3_3", cpl_matrix_get(cd, 2, 2) * 1e10,
01181 !cpl_errorstate_is_equal(es));
01182 }
01183 }
01184 return CPL_ERROR_NONE;
01185 }
01186
01187
01197
01198 static cpl_error_code
01199 muse_resampling_fit_data_to_grid(muse_datacube *aCube, muse_pixtable *aPixtable)
01200 {
01201 cpl_ensure_code(aCube && aPixtable, CPL_ERROR_NULL_INPUT);
01202 const char func[] = "muse_resampling_cube";
01203
01204
01205
01206 if (cpl_propertylist_has(aCube->header, "MUSE_RESAMPLING_WCS_OVERRIDDEN") &&
01207 cpl_propertylist_get_bool(aCube->header, "MUSE_RESAMPLING_WCS_OVERRIDDEN")) {
01208 cpl_msg_debug(func, "Output WCS was forced, won't update CRPIX[12]!");
01209 cpl_propertylist_erase(aCube->header, "MUSE_RESAMPLING_WCS_OVERRIDDEN");
01210 return CPL_ERROR_NONE;
01211 }
01212
01213 double xoff1, yoff1, xoff2, yoff2;
01214 float xlo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XLO),
01215 xhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_XHI),
01216 ylo = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YLO),
01217 yhi = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_YHI);
01218 muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
01219 if (wcstype == MUSE_PIXTABLE_WCS_CELSPH) {
01220 muse_wcs_pixel_from_celestial(aCube->header, xlo, ylo, &xoff1, &yoff1);
01221 muse_wcs_pixel_from_celestial(aCube->header, xhi, yhi, &xoff2, &yoff2);
01222 } else {
01223 muse_wcs_pixel_from_projplane(aCube->header, xlo, ylo, &xoff1, &yoff1);
01224 muse_wcs_pixel_from_projplane(aCube->header, xhi, yhi, &xoff2, &yoff2);
01225 }
01226
01227
01228 double xoff = fmin(xoff1, xoff2) - 1,
01229 yoff = fmin(yoff1, yoff2) - 1;
01230 if (xoff != 0. || yoff != 0.) {
01231 double crpix1 = muse_pfits_get_crpix(aCube->header, 1),
01232 crpix2 = muse_pfits_get_crpix(aCube->header, 2);
01233 cpl_msg_debug(func, "Updating CRPIX[12]=%#g,%#g (were %#g,%#g)",
01234 crpix1 - xoff, crpix2 - yoff, crpix1, crpix2);
01235 crpix1 -= xoff;
01236 crpix2 -= yoff;
01237 cpl_propertylist_update_double(aCube->header, "CRPIX1", crpix1);
01238 cpl_propertylist_update_double(aCube->header, "CRPIX2", crpix2);
01239 }
01240 return CPL_ERROR_NONE;
01241 }
01242
01243
01244
01245 const char *crtypestring[] = {
01246 "IRAF-like",
01247 "average",
01248 "median"
01249 };
01250
01251 #define MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_START \
01252 cpl_size idx = muse_pixgrid_get_index(aPixgrid, i2, j2, l2, CPL_FALSE),\
01253 nrow = muse_pixgrid_get_count(aPixgrid, idx), \
01254 irow; \
01255 const cpl_size *rows = muse_pixgrid_get_rows(aPixgrid, idx); \
01256 for (irow = 0; irow < nrow; irow++) { \
01257 if (muse_quality_is_usable(dq[rows[irow]], badinclude)) { \
01258 \
01259 } else if (dq[rows[irow]]) { \
01260 continue; \
01261 }
01262 #define MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_DEBUG \
01263 if (debug & 4 && i+1 == debugx && j+1 == debugy && l+1 == debugz) { \
01264 printf("%s: data / stat (%"CPL_SIZE_FORMAT") = %e / %e\n", \
01265 __func__, npixels+1, data[rows[irow]], stat[rows[irow]]); \
01266 }
01267 #define MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END \
01268 if (aType == MUSE_RESAMPLING_CRSTATS_IRAF) { \
01269 \
01270 level += data[rows[irow]]; \
01271 dev += stat[rows[irow]]; \
01272 } else { \
01273 \
01274 if (npixels >= sxsize) { \
01275 \
01276 sxsize += MUSE_CRREJECT_MAX_NPIX; \
01277 cpl_image *tmp = cpl_image_new(sxsize, 2, CPL_TYPE_DOUBLE); \
01278 cpl_image_copy(tmp, sdata, 1, 1); \
01279 cpl_image_delete(sdata); \
01280 sdata = tmp; \
01281 vdata = cpl_image_get_data_double(sdata); \
01282 } \
01283 vdata[npixels] = data[rows[irow]]; \
01284 vdata[npixels + sxsize] = stat[rows[irow]]; \
01285 } \
01286 \
01287 npixels++; \
01288 }
01289
01290
01313
01314 static void
01315 muse_resampling_crreject(muse_pixtable *aPixtable, muse_pixgrid *aPixgrid,
01316 muse_resampling_crstats_type aType, double aSigma)
01317 {
01318 const char *id = "muse_resampling_cube";
01319 if (aType > MUSE_RESAMPLING_CRSTATS_MEDIAN) {
01320 cpl_msg_warning(id, "Unknown type (%u) for cosmic-ray rejection statistics,"
01321 " resetting to \"%s\"", aType,
01322 crtypestring[MUSE_RESAMPLING_CRSTATS_MEDIAN]);
01323 aType = MUSE_RESAMPLING_CRSTATS_MEDIAN;
01324 }
01325 cpl_msg_info(id, "Using %s statistics (%.3f sigma level) for cosmic-ray"
01326 " rejection", crtypestring[aType], aSigma);
01327
01328 #ifdef ESO_ENABLE_DEBUG
01329 int debug = 0, debugx = 0, debugy = 0, debugz = 0;
01330 if (getenv("MUSE_DEBUG_CRREJECT")) {
01331 debug = atoi(getenv("MUSE_DEBUG_CRREJECT"));
01332 }
01333 if (debug & 2) {
01334 if (getenv("MUSE_DEBUG_CRREJECT_X")) {
01335 debugx = atoi(getenv("MUSE_DEBUG_CRREJECT_X"));
01336 if (debugx < 1 || debugx > aPixgrid->size_x) {
01337 debugx = 0;
01338 }
01339 }
01340 if (getenv("MUSE_DEBUG_CRREJECT_Y")) {
01341 debugy = atoi(getenv("MUSE_DEBUG_CRREJECT_Y"));
01342 if (debugy < 1 || debugy > aPixgrid->size_y) {
01343 debugy = 0;
01344 }
01345 }
01346 if (getenv("MUSE_DEBUG_CRREJECT_Z")) {
01347 debugz = atoi(getenv("MUSE_DEBUG_CRREJECT_Z"));
01348 if (debugz < 1 || debugz > aPixgrid->size_z) {
01349 debugz = 0;
01350 }
01351 }
01352 }
01353 #endif
01354
01355 enum dirtype {
01356 DIR_X = 0,
01357 DIR_Y = 1,
01358 DIR_NONE = 2
01359 };
01360 enum dirtype dir = DIR_NONE;
01361
01362 cpl_boolean haswcs = muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_CELSPH;
01363 double posang = muse_astro_posangle(aPixtable->header);
01364 const double palimit = 5.;
01365 if (!haswcs || (fabs(posang) < palimit || fabs(fabs(posang) - 180.) < palimit ||
01366 fabs(fabs(posang) - 360.) < palimit)) {
01367 cpl_msg_debug(id, "CR rejection: posang = %f deg --> DIR_X "
01368 "(%s / %s / %s / %s)", posang, haswcs ? "yes": "no",
01369 fabs(posang) < palimit ? "true" : "false",
01370 fabs(fabs(posang) - 180.) < palimit ? "true" : "false",
01371 fabs(fabs(posang) - 360.) < palimit ? "true" : "false");
01372 dir = DIR_X;
01373 } else if (fabs(fabs(posang) - 90.) < palimit ||
01374 fabs(fabs(posang) - 270.) < palimit) {
01375 cpl_msg_debug(id, "CR rejection: posang = %f deg --> DIR_Y "
01376 "(%s / %s / %s)", posang, haswcs ? "yes": "no",
01377 fabs(fabs(posang) - 90.) < palimit ? "true" : "false",
01378 fabs(fabs(posang) - 270.) < palimit ? "true" : "false");
01379 dir = DIR_Y;
01380 } else {
01381 cpl_msg_debug(id, "CR rejection: posang = %f deg --> DIR_NONE "
01382 "(%s / %s / %s / %s / %s / %s)", posang, haswcs ? "yes": "no",
01383 fabs(posang) < palimit ? "true" : "false",
01384 fabs(fabs(posang) - 90.) < palimit ? "true" : "false",
01385 fabs(fabs(posang) - 180.) < palimit ? "true" : "false",
01386 fabs(fabs(posang) - 270.) < palimit ? "true" : "false",
01387 fabs(fabs(posang) - 360.) < palimit ? "true" : "false");
01388 }
01389
01390
01391 float *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
01392 *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
01393 int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
01394
01395
01396 uint32_t badinclude = EURO3D_COSMICRAY;
01397 int l;
01398 #ifdef ESO_ENABLE_DEBUG
01399 #pragma omp parallel for default(none) \
01400 shared(aPixgrid, aPixtable, aSigma, badinclude, aType, crtypestring, \
01401 data, debug, debugx, debugy, debugz, dir, dq, id, stat, stdout)
01402 #else
01403 #pragma omp parallel for default(none) \
01404 shared(aPixgrid, aSigma, aType, badinclude, crtypestring, data, dir, \
01405 dq, id, stat)
01406 #endif
01407 for (l = 0; l < aPixgrid->size_z; l++) {
01408 #define MUSE_CRREJECT_MAX_NPIX 1000
01409 int sxsize = MUSE_CRREJECT_MAX_NPIX;
01410 cpl_image *sdata = NULL;
01411 double *vdata = NULL;
01412 if (aType > MUSE_RESAMPLING_CRSTATS_IRAF) {
01413
01414
01415 sdata = cpl_image_new(sxsize, 2, CPL_TYPE_DOUBLE);
01416 vdata = cpl_image_get_data_double(sdata);
01417 }
01418
01419 int i;
01420 for (i = 0; i < aPixgrid->size_x; i++) {
01421 int j;
01422 for (j = 0; j < aPixgrid->size_y; j++) {
01423
01424
01425 double level = 0., dev = 0;
01426
01427 #define CR_LD 1
01428 cpl_size npixels = 0;
01429 int i2, j2, l2;
01430 if (dir == DIR_Y) {
01431 for (i2 = i - CR_LD; i2 <= i + CR_LD; i2++) {
01432 int nwidth = CR_LD;
01433 if (i2 == i) {
01434 nwidth = 0;
01435 }
01436 for (j2 = j - nwidth; j2 <= j + nwidth; j2++) {
01437 for (l2 = l - nwidth; l2 <= l + nwidth; l2++) {
01438 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_START
01439 #ifdef ESO_ENABLE_DEBUG
01440 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_DEBUG
01441 #endif
01442 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END
01443 }
01444 }
01445 }
01446 } else {
01447 for (j2 = j - CR_LD; j2 <= j + CR_LD; j2++) {
01448 int nwidth = CR_LD;
01449 if (dir == DIR_X && j2 == j) {
01450 nwidth = 0;
01451 }
01452 for (i2 = i - nwidth; i2 <= i + nwidth; i2++) {
01453 for (l2 = l - nwidth; l2 <= l + nwidth; l2++) {
01454 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_START
01455 #ifdef ESO_ENABLE_DEBUG
01456 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_DEBUG
01457 #endif
01458 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END
01459 }
01460 }
01461 }
01462 }
01463
01464 if (npixels < 3) {
01465 continue;
01466 }
01467 if (aType == MUSE_RESAMPLING_CRSTATS_IRAF) {
01468 level /= npixels;
01469 dev = sqrt(dev) / npixels;
01470 } else {
01471 unsigned sflags;
01472 if (aType == MUSE_RESAMPLING_CRSTATS_MEDIAN) {
01473 sflags = CPL_STATS_MEDIAN | CPL_STATS_MAD;
01474 } else {
01475 sflags = CPL_STATS_MEAN | CPL_STATS_STDEV;
01476 }
01477 cpl_stats *sstats = cpl_stats_new_from_image_window(sdata, sflags,
01478 1, 1, npixels, 1);
01479 if (aType == MUSE_RESAMPLING_CRSTATS_MEDIAN) {
01480 level = cpl_stats_get_median(sstats);
01481 dev = cpl_stats_get_mad(sstats);
01482 } else {
01483 level = cpl_stats_get_mean(sstats);
01484 dev = cpl_stats_get_stdev(sstats);
01485 }
01486 cpl_stats_delete(sstats);
01487 }
01488 double limit = level + aSigma * dev;
01489 #ifdef ESO_ENABLE_DEBUG
01490 if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
01491 if (aType == MUSE_RESAMPLING_CRSTATS_IRAF) {
01492 printf("%s: %03d,%03d,%04d: %.3f+/-%.3f (stats), %.3f (limit) (%"
01493 CPL_SIZE_FORMAT" values)\n", __func__, i+1, j+1, l+1, level,
01494 dev, limit, npixels);
01495 } else {
01496 cpl_stats *ssdata = cpl_stats_new_from_image_window(sdata, CPL_STATS_ALL,
01497 1, 1, npixels, 1),
01498 *ssstat = cpl_stats_new_from_image_window(sdata, CPL_STATS_ALL,
01499 1, 2, npixels, 2);
01500 printf("%s: %03d,%03d,%04d: %e +/- %e (%s), %e (limit) (%"CPL_SIZE_FORMAT
01501 " values); data stats:\n", __func__, i+1, j+1, l+1,
01502 level, dev, crtypestring[aType], limit, npixels);
01503 cpl_stats_dump(ssdata, CPL_STATS_ALL, stdout);
01504 printf("%s: variance stats:\n", __func__);
01505 cpl_stats_dump(ssstat, CPL_STATS_ALL, stdout);
01506 fflush(stdout);
01507 cpl_stats_delete(ssdata);
01508 cpl_stats_delete(ssstat);
01509 }
01510 }
01511 #endif
01512
01513
01514 cpl_size idx = muse_pixgrid_get_index(aPixgrid, i, j, l, CPL_FALSE),
01515 nrow = muse_pixgrid_get_count(aPixgrid, idx),
01516 irow;
01517 const cpl_size *rows = muse_pixgrid_get_rows(aPixgrid, idx);
01518 for (irow = 0; irow < nrow; irow++) {
01519 if (data[rows[irow]] > limit) {
01520 dq[rows[irow]] |= EURO3D_COSMICRAY;
01521 #ifdef ESO_ENABLE_DEBUG
01522 if (debug & 1 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
01523 printf("%s: %03d,%03d,%04d: rejected row %"CPL_SIZE_FORMAT" (%"
01524 CPL_SIZE_FORMAT" of %"CPL_SIZE_FORMAT" in this gridcell):\t",
01525 __func__, i+1, j+1, l+1, rows[irow], irow+1, nrow);
01526 muse_pixtable_dump(aPixtable, rows[irow], 1, 0);
01527 }
01528 #endif
01529 }
01530 }
01531 #ifdef ESO_ENABLE_DEBUG
01532 if (debug) {
01533 fflush(stdout);
01534 }
01535 #endif
01536 }
01537 }
01538 cpl_image_delete(sdata);
01539 }
01540 }
01541
01542
01565
01566 muse_euro3dcube *
01567 muse_resampling_euro3d(muse_pixtable *aPixtable,
01568 muse_resampling_params *aParams)
01569 {
01570 cpl_ensure(aParams, CPL_ERROR_NULL_INPUT, NULL);
01571 cpl_ensure(aParams->method < MUSE_RESAMPLE_NONE, CPL_ERROR_ILLEGAL_INPUT, NULL);
01572
01573
01574
01575 muse_datacube *cube = muse_resampling_cube(aPixtable, aParams, NULL);
01576 if (!cube) {
01577 return NULL;
01578 }
01579
01580
01581 muse_euro3dcube *e3d = cpl_calloc(1, sizeof(muse_euro3dcube));
01582 e3d->header = cpl_propertylist_duplicate(cube->header);
01583 cpl_propertylist_erase_regexp(e3d->header,
01584 "^SIMPLE$|^BITPIX$|^NAXIS|^EURO3D$|^E3D_",
01585 0);
01586 cpl_propertylist_append_char(e3d->header, "EURO3D", 'T');
01587 cpl_propertylist_set_comment(e3d->header, "EURO3D",
01588 "file conforms to Euro3D standard");
01589
01590 cpl_propertylist_append_string(e3d->header, "E3D_VERS", "1.0");
01591 cpl_propertylist_set_comment(e3d->header, "E3D_VERS",
01592 "version number of the Euro3D format");
01593 cpl_errorstate prestate = cpl_errorstate_get();
01594 double darlambdaref = cpl_propertylist_get_double(e3d->header,
01595 MUSE_HDR_PT_DAR_NAME);
01596 if (!cpl_errorstate_is_equal(prestate)) {
01597 darlambdaref = -1.;
01598 cpl_errorstate_set(prestate);
01599 }
01600 if (darlambdaref > 0.) {
01601 cpl_propertylist_append_char(e3d->header, "E3D_ADC", 'T');
01602 cpl_propertylist_set_comment(e3d->header, "E3D_ADC",
01603 "data was corrected for atmospheric dispersion");
01604 } else {
01605 cpl_propertylist_append_char(e3d->header, "E3D_ADC", 'F');
01606 cpl_propertylist_set_comment(e3d->header, "E3D_ADC",
01607 "data not corrected for atmospheric dispersion");
01608 }
01609
01610
01611 e3d->hdata = cpl_propertylist_new();
01612 cpl_propertylist_append_string(e3d->hdata, "EXTNAME", "E3D_DATA");
01613 cpl_propertylist_set_comment(e3d->hdata, "EXTNAME",
01614 "This is the Euro3D data table extension");
01615 cpl_propertylist_append_string(e3d->hdata, "CTYPES",
01616 muse_pfits_get_ctype(e3d->header, 3));
01617 cpl_propertylist_set_comment(e3d->hdata, "CTYPES",
01618 cpl_propertylist_get_comment(e3d->header, "CTYPE3"));
01619 cpl_propertylist_append_string(e3d->hdata, "CUNITS",
01620 muse_pfits_get_cunit(e3d->header, 3));
01621 cpl_propertylist_set_comment(e3d->hdata, "CUNITS",
01622 cpl_propertylist_get_comment(e3d->header, "CUNIT3"));
01623 cpl_propertylist_append_double(e3d->hdata, "CRVALS",
01624 muse_pfits_get_crval(e3d->header, 3));
01625 cpl_propertylist_set_comment(e3d->hdata, "CRVALS",
01626 cpl_propertylist_get_comment(e3d->header, "CRVAL3"));
01627 cpl_propertylist_append_double(e3d->hdata, "CDELTS",
01628 muse_pfits_get_cd(e3d->header, 3, 3));
01629 cpl_propertylist_set_comment(e3d->hdata, "CDELTS",
01630 cpl_propertylist_get_comment(e3d->header, "CD3_3"));
01631
01632 cpl_propertylist_erase_regexp(e3d->header, "^C...*3$|^CD3_.$", 0);
01633
01634
01635 int nx = cpl_image_get_size_x(cpl_imagelist_get(cube->data, 0)),
01636 ny = cpl_image_get_size_y(cpl_imagelist_get(cube->data, 0)),
01637 nz = cpl_imagelist_get_size(cube->data);
01638 e3d->dtable = muse_cpltable_new(muse_euro3dcube_e3d_data_def, nx * ny);
01639
01640
01641
01642 cpl_table_set_column_depth(e3d->dtable, "DATA_SPE", nz);
01643 cpl_table_set_column_depth(e3d->dtable, "QUAL_SPE", nz);
01644 cpl_table_set_column_depth(e3d->dtable, "STAT_SPE", nz);
01645
01646 cpl_table_set_column_unit(e3d->dtable, "DATA_SPE",
01647 cpl_table_get_column_unit(aPixtable->table,
01648 MUSE_PIXTABLE_DATA));
01649 cpl_table_set_column_unit(e3d->dtable, "STAT_SPE",
01650 cpl_table_get_column_unit(aPixtable->table,
01651 MUSE_PIXTABLE_STAT));
01652
01653 cpl_table_set_column_savetype(e3d->dtable, "SELECTED", CPL_TYPE_BOOL);
01654
01655 muse_wcs *wcs = muse_wcs_new(cube->header);
01656 wcs->iscelsph = muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_CELSPH;
01657 if (wcs->iscelsph) {
01658 cpl_table_set_column_unit(e3d->dtable, "XPOS", "deg");
01659 cpl_table_set_column_unit(e3d->dtable, "YPOS", "deg");
01660 }
01661
01662
01663 unsigned euro3d_ignore = EURO3D_OUTSDRANGE | EURO3D_MISSDATA;
01664 cpl_vector *vusable = cpl_vector_new(nx * ny);
01665 int i, itable = 0;
01666 for (i = 1; i <= nx; i++) {
01667 int j;
01668 for (j = 1; j <= ny; j++) {
01669 cpl_array *adata = cpl_array_new(nz, CPL_TYPE_FLOAT),
01670 *adq = cpl_array_new(nz, CPL_TYPE_INT),
01671 *astat = cpl_array_new(nz, CPL_TYPE_FLOAT);
01672
01673 double x, y;
01674 if (wcs->iscelsph) {
01675 muse_wcs_celestial_from_pixel_fast(wcs, i, j, &x, &y);
01676 } else {
01677 muse_wcs_projplane_from_pixel_fast(wcs, i, j, &x, &y);
01678 }
01679
01680 int l, nusable = 0, nstart = -1;
01681 for (l = 0; l < nz; l++) {
01682 int err;
01683 unsigned dq = cpl_image_get(cpl_imagelist_get(cube->dq, l), i, j, &err);
01684
01685 if (nstart < 0 && (dq & euro3d_ignore)) {
01686 continue;
01687 }
01688 cpl_array_set_int(adq, nusable, dq);
01689 cpl_array_set_float(adata, nusable,
01690 cpl_image_get(cpl_imagelist_get(cube->data, l),
01691 i, j, &err));
01692
01693 cpl_array_set_float(astat, nusable,
01694 sqrt(cpl_image_get(cpl_imagelist_get(cube->stat, l),
01695 i, j, &err)));
01696 nusable++;
01697 if (nstart < 0) {
01698 nstart = l + 1;
01699 }
01700 }
01701
01702 cpl_table_set_int(e3d->dtable, "SPEC_ID", itable, itable + 1);
01703 cpl_table_set_int(e3d->dtable, "SPEC_LEN", itable, nusable);
01704 cpl_table_set_int(e3d->dtable, "SPEC_STA", itable, nstart);
01705 cpl_table_set_double(e3d->dtable, "XPOS", itable, x);
01706 cpl_table_set_double(e3d->dtable, "YPOS", itable, y);
01707
01708
01709 cpl_table_set_array(e3d->dtable, "DATA_SPE", itable, adata);
01710 cpl_table_set_array(e3d->dtable, "QUAL_SPE", itable, adq);
01711 cpl_table_set_array(e3d->dtable, "STAT_SPE", itable, astat);
01712
01713 cpl_array_delete(adata);
01714 cpl_array_delete(adq);
01715 cpl_array_delete(astat);
01716
01717 cpl_vector_set(vusable, itable, nusable);
01718 if (nstart != -1 && nusable > 0) {
01719
01720 cpl_table_unselect_row(e3d->dtable, itable);
01721 }
01722 itable++;
01723 }
01724 }
01725 cpl_free(wcs);
01726
01727 cpl_vector_set_size(vusable, itable);
01728 int maxusable = cpl_vector_get_max(vusable);
01729 #if 0
01730 cpl_msg_debug(__func__, "filled %d of %d spaxels, usable are "
01731 "%f..%f(%f)+/-%f..%d pixels per spectrum", itable, nx * ny,
01732 cpl_vector_get_min(vusable), cpl_vector_get_mean(vusable),
01733 cpl_vector_get_median(vusable), cpl_vector_get_stdev(vusable),
01734 maxusable);
01735 #endif
01736 cpl_msg_debug(__func__, "filled %"CPL_SIZE_FORMAT" of %d spaxels, usable "
01737 "are max. %d of %d pixels per spectrum%s",
01738 itable - cpl_table_count_selected(e3d->dtable), nx * ny,
01739 maxusable, nz, maxusable == nz ? "" : ", cutting table");
01740 if (maxusable != nz) {
01741
01742 cpl_table_set_column_depth(e3d->dtable, "DATA_SPE", maxusable);
01743 cpl_table_set_column_depth(e3d->dtable, "QUAL_SPE", maxusable);
01744 cpl_table_set_column_depth(e3d->dtable, "STAT_SPE", maxusable);
01745 }
01746
01747 cpl_table_erase_selected(e3d->dtable);
01748 cpl_vector_delete(vusable);
01749
01750
01751
01752 cpl_table_fill_column_window_int(e3d->dtable, "SELECTED", 0, itable,
01753 1 );
01754
01755 cpl_table_fill_column_window_int(e3d->dtable, "NSPAX", 0, itable, 1);
01756
01757 cpl_table_fill_column_window_int(e3d->dtable, "GROUP_N", 0, itable, 1);
01758
01759
01760 muse_datacube_delete(cube);
01761
01762
01763 e3d->hgroup = cpl_propertylist_new();
01764 cpl_propertylist_append_string(e3d->hgroup, "EXTNAME", "E3D_GRP");
01765 cpl_propertylist_set_comment(e3d->hgroup, "EXTNAME",
01766 "This is the Euro3D group table extension");
01767 e3d->gtable = muse_cpltable_new(muse_euro3dcube_e3d_grp_def, 1);
01768 cpl_table_set_int(e3d->gtable, "GROUP_N", 0, 1);
01769 cpl_table_set_string(e3d->gtable, "G_SHAPE", 0, "R");
01770 cpl_table_set_float(e3d->gtable, "G_SIZE1", 0, aParams->dx);
01771 cpl_table_set_float(e3d->gtable, "G_ANGLE", 0, 0.);
01772 cpl_table_set_float(e3d->gtable, "G_SIZE2", 0, aParams->dy);
01773 if (darlambdaref > 0.) {
01774
01775 cpl_table_set_float(e3d->gtable, "G_POSWAV", 0, NAN);
01776 cpl_table_set_float(e3d->gtable, "G_AIRMAS", 0, NAN);
01777 cpl_table_set_float(e3d->gtable, "G_PARANG", 0, NAN);
01778 cpl_table_set_float(e3d->gtable, "G_PRESSU", 0, NAN);
01779 cpl_table_set_float(e3d->gtable, "G_TEMPER", 0, NAN);
01780 cpl_table_set_float(e3d->gtable, "G_HUMID", 0, NAN);
01781 } else {
01782
01783 cpl_table_set_float(e3d->gtable, "G_POSWAV", 0,
01784 (kMuseNominalLambdaMin + kMuseNominalLambdaMax) / 2.);
01785 cpl_table_set_float(e3d->gtable, "G_AIRMAS", 0,
01786 muse_astro_airmass(e3d->header));
01787 cpl_table_set_float(e3d->gtable, "G_PARANG", 0,
01788 muse_astro_parangle(e3d->header));
01789 cpl_table_set_float(e3d->gtable, "G_PRESSU", 0,
01790 (muse_pfits_get_pres_start(e3d->header)
01791 + muse_pfits_get_pres_start(e3d->header)) / 2.);
01792 cpl_table_set_float(e3d->gtable, "G_TEMPER", 0,
01793 muse_pfits_get_temp(e3d->header) + 273.15);
01794 cpl_table_set_float(e3d->gtable, "G_HUMID", 0,
01795 muse_pfits_get_rhum(e3d->header));
01796 }
01797
01798 return e3d;
01799 }
01800
01801
01852
01853 muse_datacube *
01854 muse_resampling_cube(muse_pixtable *aPixtable, muse_resampling_params *aParams,
01855 muse_pixgrid **aPixgrid)
01856 {
01857 cpl_ensure(aPixtable && aParams, CPL_ERROR_NULL_INPUT, NULL);
01858 cpl_ensure(muse_pixtable_get_type(aPixtable) == MUSE_PIXTABLE_TYPE_FULL,
01859 CPL_ERROR_ILLEGAL_INPUT, NULL);
01860 muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
01861 cpl_ensure(wcstype == MUSE_PIXTABLE_WCS_CELSPH ||
01862 wcstype == MUSE_PIXTABLE_WCS_PIXEL, CPL_ERROR_UNSUPPORTED_MODE,
01863 NULL);
01864
01865
01866
01867 muse_resampling_check_deltas(aPixtable, aParams);
01868
01869
01870
01871 int xsize = 0, ysize = 0, zsize = 0;
01872 muse_resampling_compute_size(aPixtable, aParams, &xsize, &ysize, &zsize);
01873 muse_resampling_override_size(&xsize, &ysize, &zsize, aParams);
01874 cpl_ensure(xsize > 0 && ysize > 0 && zsize > 0, CPL_ERROR_ILLEGAL_OUTPUT,
01875 NULL);
01876
01877 double time = cpl_test_get_walltime();
01878
01879
01880 muse_datacube *cube = cpl_calloc(1, sizeof(muse_datacube));
01881
01882
01883 cube->header = cpl_propertylist_duplicate(aPixtable->header);
01884 cpl_propertylist_erase_regexp(cube->header,
01885 "^SIMPLE$|^BITPIX$|^NAXIS|^EXTEND$|^XTENSION$|"
01886 "^DATASUM$|^DATAMIN$|^DATAMAX$|^DATAMD5$|"
01887 "^PCOUNT$|^GCOUNT$|^HDUVERS$|^BLANK$|"
01888 "^BZERO$|^BSCALE$|^CHECKSUM$|^INHERIT$|"
01889 "^EXTNAME$|"MUSE_WCS_KEYS"|"MUSE_HDR_PT_REGEXP,
01890 0);
01891
01892 cpl_propertylist_update_string(cube->header, "BUNIT",
01893 cpl_table_get_column_unit(aPixtable->table,
01894 MUSE_PIXTABLE_DATA));
01895
01896 cpl_propertylist_update_int(cube->header, "NAXIS", 3);
01897 cpl_propertylist_update_int(cube->header, "NAXIS1", xsize);
01898 cpl_propertylist_update_int(cube->header, "NAXIS2", ysize);
01899 cpl_propertylist_update_int(cube->header, "NAXIS3", zsize);
01900
01901
01902 if (wcstype == MUSE_PIXTABLE_WCS_CELSPH) {
01903 cpl_propertylist_copy_property_regexp(cube->header, aPixtable->header,
01904 MUSE_WCS_KEYS, 0);
01905 cpl_propertylist_update_double(cube->header, "CD1_1", -aParams->dx);
01906 cpl_propertylist_update_double(cube->header, "CD2_2", aParams->dy);
01907 cpl_propertylist_update_double(cube->header, "CD1_2", 0.);
01908 cpl_propertylist_update_double(cube->header, "CD2_1", 0.);
01909
01910 cpl_propertylist_update_double(cube->header, "CRPIX1",
01911 muse_pfits_get_crpix(cube->header, 1)
01912 + (1. + xsize) / 2.);
01913 cpl_propertylist_update_double(cube->header, "CRPIX2",
01914 muse_pfits_get_crpix(cube->header, 2)
01915 + (1. + ysize) / 2.);
01916
01917 cpl_propertylist_erase(cube->header, "WCSAXES");
01918 } else {
01919
01920 cpl_propertylist_append_double(cube->header, "CRPIX1", 1.);
01921 cpl_propertylist_append_double(cube->header, "CRVAL1",
01922 cpl_propertylist_get_float(aPixtable->header,
01923 MUSE_HDR_PT_XLO));
01924
01925
01926 cpl_propertylist_append_string(cube->header, "CTYPE1", "PIXEL");
01927 cpl_propertylist_append_string(cube->header, "CUNIT1", "pixel");
01928 cpl_propertylist_append_double(cube->header, "CD1_1", aParams->dx);
01929 cpl_propertylist_append_double(cube->header, "CRPIX2", 1.);
01930 cpl_propertylist_append_double(cube->header, "CRVAL2",
01931 cpl_propertylist_get_float(aPixtable->header,
01932 MUSE_HDR_PT_YLO));
01933 cpl_propertylist_append_string(cube->header, "CTYPE2", "PIXEL");
01934 cpl_propertylist_append_string(cube->header, "CUNIT2", "pixel");
01935 cpl_propertylist_append_double(cube->header, "CD2_2", aParams->dy);
01936
01937 cpl_propertylist_append_double(cube->header, "CD1_2", 0.);
01938 cpl_propertylist_append_double(cube->header, "CD2_1", 0.);
01939 }
01940 switch (aParams->tlambda) {
01941 case MUSE_RESAMPLING_DISP_AWAV:
01942 cpl_propertylist_append_string(cube->header, "CTYPE3", "AWAV");
01943 break;
01944 case MUSE_RESAMPLING_DISP_AWAV_LOG:
01945 cpl_propertylist_append_string(cube->header, "CTYPE3", "AWAV-LOG");
01946 break;
01947 case MUSE_RESAMPLING_DISP_WAVE:
01948 cpl_propertylist_append_string(cube->header, "CTYPE3", "WAVE");
01949 break;
01950 case MUSE_RESAMPLING_DISP_WAVE_LOG:
01951 cpl_propertylist_append_string(cube->header, "CTYPE3", "WAVE-LOG");
01952 break;
01953 default:
01954 cpl_propertylist_append_string(cube->header, "CTYPE3", "UNKNOWN");
01955 }
01956 cpl_propertylist_append_string(cube->header, "CUNIT3", "Angstrom");
01957 cpl_propertylist_append_double(cube->header, "CD3_3", aParams->dlambda);
01958 cpl_propertylist_append_double(cube->header, "CRPIX3", 1.);
01959 cpl_propertylist_append_double(cube->header, "CRVAL3",
01960 cpl_propertylist_get_float(aPixtable->header,
01961 MUSE_HDR_PT_LLO));
01962
01963 cpl_propertylist_append_double(cube->header, "CD1_3", 0.);
01964 cpl_propertylist_append_double(cube->header, "CD2_3", 0.);
01965 cpl_propertylist_append_double(cube->header, "CD3_1", 0.);
01966 cpl_propertylist_append_double(cube->header, "CD3_2", 0.);
01967
01968
01969 muse_resampling_override_wcs(cube, aParams);
01970
01971 muse_resampling_fit_data_to_grid(cube, aPixtable);
01972
01973 if (aParams->method < MUSE_RESAMPLE_NONE) {
01974
01975 cube->data = cpl_imagelist_new();
01976 cube->dq = cpl_imagelist_new();
01977 cube->stat = cpl_imagelist_new();
01978 int i;
01979 for (i = 0; i < zsize; i++) {
01980 cpl_image *image = cpl_image_new(xsize, ysize, CPL_TYPE_FLOAT);
01981
01982 cpl_imagelist_set(cube->data, image, i);
01983
01984 cpl_imagelist_set(cube->stat, cpl_image_duplicate(image), i);
01985
01986
01987
01988 image = cpl_image_new(xsize, ysize, CPL_TYPE_INT);
01989
01990 cpl_image_add_scalar(image, EURO3D_OUTSDRANGE);
01991 cpl_imagelist_set(cube->dq, image, i);
01992 }
01993 }
01994
01995 muse_utils_memory_dump("muse_resampling_cube() before pixgrid");
01996
01997 muse_pixgrid *pixgrid = muse_pixgrid_create(aPixtable, cube->header,
01998 xsize, ysize, zsize);
01999
02000 if (aParams->crsigma > 0) {
02001 muse_resampling_crreject(aPixtable, pixgrid, aParams->crtype,
02002 aParams->crsigma);
02003 }
02004 muse_utils_memory_dump("muse_resampling_cube() after pixgrid and crreject");
02005
02006 double timeinit = cpl_test_get_walltime(),
02007 cpuinit = cpl_test_get_cputime();
02008
02009
02010 cpl_error_code rc = CPL_ERROR_NONE;
02011 switch (aParams->method) {
02012 case MUSE_RESAMPLE_NEAREST:
02013 cpl_msg_info(__func__, "Starting resampling, using method \"nearest\"");
02014 rc = muse_resampling_cube_nearest(cube, aPixtable, pixgrid);
02015 break;
02016 case MUSE_RESAMPLE_WEIGHTED_RENKA:
02017 cpl_msg_info(__func__, "Starting resampling, using method \"renka\" "
02018 "(critical radius rc=%f, loop distance ld=%d)", aParams->rc,
02019 aParams->ld);
02020 rc = muse_resampling_cube_weighted(cube, aPixtable, pixgrid, aParams);
02021 break;
02022 case MUSE_RESAMPLE_WEIGHTED_LINEAR:
02023 case MUSE_RESAMPLE_WEIGHTED_QUADRATIC:
02024 case MUSE_RESAMPLE_WEIGHTED_LANCZOS:
02025 cpl_msg_info(__func__, "Starting resampling, using method \"%s\" (loop "
02026 "distance ld=%d)",
02027 aParams->method == MUSE_RESAMPLE_WEIGHTED_LINEAR
02028 ? "linear"
02029 : (aParams->method == MUSE_RESAMPLE_WEIGHTED_QUADRATIC
02030 ? "quadratic"
02031 : "lanczos"),
02032 aParams->ld);
02033 rc = muse_resampling_cube_weighted(cube, aPixtable, pixgrid, aParams);
02034 break;
02035 case MUSE_RESAMPLE_WEIGHTED_DRIZZLE:
02036 cpl_msg_info(__func__, "Starting resampling, using method \"drizzle\" "
02037 "(pixfrac f=%.3f,%.3f,%.3f, loop distance ld=%d)",
02038 aParams->pfx, aParams->pfy, aParams->pfl, aParams->ld);
02039 rc = muse_resampling_cube_weighted(cube, aPixtable, pixgrid, aParams);
02040 break;
02041 case MUSE_RESAMPLE_NONE:
02042 cpl_msg_debug(__func__, "Method %d (no resampling)", aParams->method);
02043 break;
02044 default:
02045 cpl_msg_error(__func__, "Don't know this resampling method: %d",
02046 aParams->method);
02047 cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
02048 rc = CPL_ERROR_UNSUPPORTED_MODE;
02049 }
02050 muse_utils_memory_dump("muse_resampling_cube() after resampling");
02051
02052 double timefini = cpl_test_get_walltime(),
02053 cpufini = cpl_test_get_cputime();
02054
02055
02056 if (aPixgrid) {
02057 *aPixgrid = pixgrid;
02058 } else {
02059 muse_pixgrid_delete(pixgrid);
02060 }
02061
02062 cpl_msg_debug(__func__, "resampling took %.3fs (wall-clock) and %.3fs "
02063 "(%.3fs CPU, %d CPUs) for muse_resampling_cube*() alone",
02064 timefini - time, timefini - timeinit, cpufini - cpuinit,
02065 omp_get_max_threads());
02066 if (rc != CPL_ERROR_NONE) {
02067 cpl_msg_error(__func__, "resampling failed: %s", cpl_error_get_message());
02068 muse_datacube_delete(cube);
02069 return NULL;
02070 }
02071
02072 muse_utils_memory_dump("muse_resampling_cube() end");
02073 return cube;
02074 }
02075
02076
02106
02107 muse_image *
02108 muse_resampling_collapse_pixgrid(muse_pixtable *aPixtable, muse_pixgrid *aPixgrid,
02109 muse_datacube *aCube, cpl_table *aFilter,
02110 muse_resampling_params *aParams)
02111 {
02112 cpl_ensure(aPixtable && aPixtable->table && aPixgrid && aParams &&
02113 aCube && aCube->header, CPL_ERROR_NULL_INPUT, NULL);
02114 cpl_ensure(aParams->method < MUSE_RESAMPLE_NONE &&
02115 aParams->method > MUSE_RESAMPLE_NEAREST, CPL_ERROR_ILLEGAL_INPUT,
02116 NULL);
02117
02118 muse_wcs *wcs = muse_wcs_new(aCube->header);
02119 wcs->iscelsph = muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_CELSPH;
02120 cpl_errorstate prestate = cpl_errorstate_get();
02121 float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
02122 *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
02123 *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
02124 *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
02125 *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT),
02126 *wght = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_WEIGHT);
02127 if (!cpl_errorstate_is_equal(prestate)) {
02128 cpl_errorstate_set(prestate);
02129 }
02130 int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
02131
02132
02133
02134
02135
02136 double xnorm = 1., ynorm = 1.,
02137 ptxoff = 0.,
02138 ptyoff = 0.;
02139 if (wcs->iscelsph) {
02140 muse_wcs_get_scales(aPixtable->header, &xnorm, &ynorm);
02141 xnorm = 1. / xnorm;
02142 ynorm = 1. / ynorm;
02143
02144 ptxoff = muse_pfits_get_crval(aPixtable->header, 1);
02145 ptyoff = muse_pfits_get_crval(aPixtable->header, 2);
02146 }
02147
02148 double renka_rc = aParams->rc
02149 * sqrt(pow(wcs->cd11*xnorm, 2) + pow(wcs->cd22*xnorm, 2));
02150
02151 int ld = aParams->ld;
02152 if (ld <= 0) {
02153 ld = 1;
02154 cpl_msg_info(__func__, "Overriding loop distance ld=%d", ld);
02155 }
02156
02157
02158 double xsz = aParams->pfx / xnorm,
02159 ysz = aParams->pfy / ynorm,
02160 xout = fabs(wcs->cd11), yout = fabs(wcs->cd22);
02161
02162 muse_image *image = muse_image_new();
02163 image->data = cpl_image_new(aPixgrid->size_x, aPixgrid->size_y, CPL_TYPE_FLOAT);
02164 image->dq = cpl_image_new(aPixgrid->size_x, aPixgrid->size_y, CPL_TYPE_INT);
02165 image->stat = cpl_image_new(aPixgrid->size_x, aPixgrid->size_y, CPL_TYPE_FLOAT);
02166 image->header = cpl_propertylist_duplicate(aCube->header);
02167 cpl_propertylist_erase_regexp(image->header, "^C...*3$|^CD3_.$", 0);
02168 float *pdata = cpl_image_get_data_float(image->data),
02169 *pstat = cpl_image_get_data_float(image->stat);
02170 int *pdq = cpl_image_get_data_int(image->dq);
02171
02172
02173 cpl_boolean usevariance = CPL_FALSE;
02174 if (getenv("MUSE_COLLAPSE_USE_VARIANCE") &&
02175 atoi(getenv("MUSE_COLLAPSE_USE_VARIANCE")) > 0) {
02176 usevariance = CPL_TRUE;
02177 }
02178
02179
02180 double lmin, lmax;
02181 if (aFilter) {
02182 const double *flbda = cpl_table_get_data_double_const(aFilter, "lambda"),
02183 *fresp = cpl_table_get_data_double_const(aFilter, "throughput");
02184 int l = 0, nl = cpl_table_get_nrow(aFilter);
02185 while (l < nl && fabs(fresp[l]) < DBL_EPSILON) {
02186 lmin = flbda[l++];
02187 }
02188 l = nl - 1;
02189 while (l > 0 && fabs(fresp[l]) < DBL_EPSILON) {
02190 lmax = flbda[l--];
02191 }
02192 cpl_msg_debug(__func__, "filter wavelength range %.1f..%.1fA", lmin, lmax);
02193 } else {
02194 lmin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LLO);
02195 lmax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LHI);
02196 cpl_msg_debug(__func__, "full wavelength range %.1f..%.1fA", lmin, lmax);
02197 }
02198
02199 int i;
02200 #pragma omp parallel for default(none) \
02201 shared(aFilter, aParams, aPixgrid, data, dq, lbda, ld, lmax, lmin, \
02202 pdata, pdq, pstat, ptxoff, ptyoff, renka_rc, stat, usevariance,\
02203 wcs, wght, xnorm, xout, xpos, xsz, ynorm, yout, ypos, ysz)
02204 for (i = 0; i < aPixgrid->size_x; i++) {
02205 int j;
02206 for (j = 0; j < aPixgrid->size_y; j++) {
02207
02208 double x, y;
02209 if (wcs->iscelsph) {
02210 muse_wcs_celestial_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
02211 } else {
02212 muse_wcs_projplane_from_pixel_fast(wcs, i + 1, j + 1, &x, &y);
02213 }
02214 double sumdata = 0, sumstat = 0, sumweight = 0;
02215 cpl_size npoints = 0;
02216
02217
02218 int i2;
02219 for (i2 = i - ld; i2 <= i + ld; i2++) {
02220 int j2;
02221 for (j2 = j - ld; j2 <= j + ld; j2++) {
02222
02223 int l2;
02224 for (l2 = 0; l2 < aPixgrid->size_z; l2++) {
02225 cpl_size idx2 = muse_pixgrid_get_index(aPixgrid, i2, j2, l2, CPL_FALSE),
02226 n, n_rows2 = muse_pixgrid_get_count(aPixgrid, idx2);
02227 const cpl_size *rows2 = muse_pixgrid_get_rows(aPixgrid, idx2);
02228 for (n = 0; n < n_rows2; n++) {
02229 if (dq[rows2[n]]) {
02230 continue;
02231 }
02232 if (usevariance && !isnormal(stat[rows2[n]])) {
02233 continue;
02234 }
02235 if (lbda[rows2[n]] > lmax || lbda[rows2[n]] < lmin) {
02236 continue;
02237 }
02238
02239 double dx = fabs(x - (xpos[rows2[n]] + ptxoff)),
02240 dy = fabs(y - (ypos[rows2[n]] + ptyoff)),
02241 r2 = 0;
02242 if (wcs->iscelsph) {
02243
02244
02245 dx *= cos(y * CPL_MATH_RAD_DEG);
02246 }
02247 if (aParams->method != MUSE_RESAMPLE_WEIGHTED_DRIZZLE) {
02248 dx *= xnorm;
02249 dy *= ynorm;
02250 r2 = dx*dx + dy*dy;
02251 }
02252 double weight = 0.;
02253 if (aParams->method == MUSE_RESAMPLE_WEIGHTED_RENKA) {
02254 weight = muse_resampling_weight_function_renka(sqrt(r2), renka_rc);
02255 } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_DRIZZLE) {
02256 weight = muse_resampling_weight_function_drizzle(xsz, ysz, 1.,
02257 xout, yout, 1.,
02258 dx, dy, 0.);
02259 } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_LINEAR) {
02260 weight = muse_resampling_weight_function_linear(sqrt(r2));
02261 } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_QUADRATIC) {
02262 weight = muse_resampling_weight_function_quadratic(r2);
02263 } else if (aParams->method == MUSE_RESAMPLE_WEIGHTED_LANCZOS) {
02264 weight = muse_resampling_weight_function_lanczos(dx, dy, 0., ld);
02265 }
02266
02267 if (wght) {
02268
02269 weight *= wght[rows2[n]];
02270 }
02271
02272
02273 if (aFilter) {
02274 weight *= muse_flux_response_interpolate(aFilter, lbda[rows2[n]],
02275 NULL, MUSE_FLUX_RESP_FILTER);
02276 }
02277 if (usevariance) {
02278 weight /= stat[rows2[n]];
02279 }
02280 sumweight += weight;
02281 sumdata += data[rows2[n]] * weight;
02282 sumstat += stat[rows2[n]] * weight*weight;
02283 npoints++;
02284 }
02285 }
02286 }
02287 }
02288
02289
02290
02291
02292 if (!npoints || !isnormal(sumweight)) {
02293 pdq[i + j * aPixgrid->size_x] = EURO3D_MISSDATA;
02294 continue;
02295 }
02296
02297
02298 sumdata /= sumweight;
02299 sumstat /= sumweight*sumweight;
02300 pdata[i + j * aPixgrid->size_x] = sumdata;
02301 pstat[i + j * aPixgrid->size_x] = sumstat;
02302 pdq[i + j * aPixgrid->size_x] = EURO3D_GOODPIXEL;
02303 }
02304 }
02305 cpl_free(wcs);
02306
02307 return image;
02308 }
02309
02310
02311
02312
02313
02314
02315
02330
02331 static cpl_error_code
02332 muse_resampling_image_nearest(muse_image *aImage, muse_pixtable *aPixtable,
02333 muse_pixgrid *aPixgrid)
02334 {
02335 cpl_ensure_code(aImage && aPixtable && aPixgrid, CPL_ERROR_NULL_INPUT);
02336 aImage->data = cpl_image_new(aPixgrid->size_x, aPixgrid->size_z,
02337 CPL_TYPE_FLOAT);
02338
02339 double crval2 = muse_pfits_get_crval(aImage->header, 2),
02340 cd22 = muse_pfits_get_cd(aImage->header, 2, 2);
02341
02342 float *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
02343 *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA);
02344 float *imagedata = cpl_image_get_data_float(aImage->data);
02345
02346 #ifdef ESO_ENABLE_DEBUG
02347 int debug = 0;
02348 if (getenv("MUSE_DEBUG_NEAREST")) {
02349 debug = atoi(getenv("MUSE_DEBUG_NEAREST"));
02350 }
02351 #endif
02352
02353 int i;
02354 for (i = 0; i < aPixgrid->size_x; i++) {
02355 int j;
02356 for (j = 0; j < aPixgrid->size_z; j++) {
02357 cpl_size idx = muse_pixgrid_get_index(aPixgrid, i, 0, j, CPL_FALSE),
02358 n_rows = muse_pixgrid_get_count(aPixgrid, idx);
02359 const cpl_size *rows = muse_pixgrid_get_rows(aPixgrid, idx);
02360 if (n_rows == 1) {
02361
02362 imagedata[i + j * aPixgrid->size_x] = data[rows[0]];
02363 #ifdef ESO_ENABLE_DEBUG
02364 if (debug) {
02365 printf("only: %f\n", data[rows[0]]);
02366 fflush(stdout);
02367 }
02368 #endif
02369 } else if (n_rows >= 2) {
02370
02371 cpl_size n, nbest = -1;
02372 double dbest = FLT_MAX;
02373 for (n = 0; n < n_rows; n++) {
02374
02375
02376 double dlambda = fabs(j * cd22 + crval2 - lbda[rows[n]]);
02377 #ifdef ESO_ENABLE_DEBUG
02378 if (debug) {
02379 printf("%f, %f, d: %f best: %f (%f)\n",
02380 j * cd22 + crval2, lbda[rows[n]],
02381 dlambda, dbest, data[rows[n]]);
02382 }
02383 #endif
02384 if (dlambda < dbest) {
02385 nbest = n;
02386 dbest = dlambda;
02387 }
02388 }
02389 imagedata[i + j * aPixgrid->size_x] = data[rows[nbest]];
02390 #ifdef ESO_ENABLE_DEBUG
02391 if (debug) {
02392 printf("closest: %f\n", data[rows[nbest]]);
02393 fflush(stdout);
02394 }
02395 #endif
02396 } else {
02397
02398 }
02399 }
02400 }
02401
02402 return CPL_ERROR_NONE;
02403 }
02404
02405
02421
02422 static cpl_error_code
02423 muse_resampling_image_weighted(muse_image *aImage, muse_pixtable *aPixtable,
02424 muse_pixgrid *aPixgrid)
02425 {
02426 cpl_ensure_code(aImage && aPixtable && aPixgrid, CPL_ERROR_NULL_INPUT);
02427 aImage->data = cpl_image_new(aPixgrid->size_x, aPixgrid->size_z,
02428 CPL_TYPE_FLOAT);
02429
02430 double crval2 = muse_pfits_get_crval(aImage->header, 2),
02431 cd22 = muse_pfits_get_cd(aImage->header, 2, 2);
02432
02433 float *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
02434 *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA);
02435 float *imagedata = cpl_image_get_data_float(aImage->data);
02436
02437 #ifdef ESO_ENABLE_DEBUG
02438 int debug = 0, debugx = 0, debugz = 0;
02439 if (getenv("MUSE_DEBUG_WEIGHTED")) {
02440 debug = atoi(getenv("MUSE_DEBUG_WEIGHTED"));
02441 }
02442 if (debug & 128) {
02443 if (getenv("MUSE_DEBUG_WEIGHTED_X")) {
02444 debugx = atoi(getenv("MUSE_DEBUG_WEIGHTED_X"));
02445 if (debugx < 1 || debugx > aPixgrid->size_x) {
02446 debugx = 0;
02447 }
02448 }
02449 if (getenv("MUSE_DEBUG_WEIGHTED_Z")) {
02450 debugz = atoi(getenv("MUSE_DEBUG_WEIGHTED_Z"));
02451 if (debugz < 1 || debugz > aPixgrid->size_z) {
02452 debugz = 0;
02453 }
02454 }
02455 }
02456 #endif
02457
02458 int ld = 1;
02459
02460
02461 double renka_rc = 1.25;
02462
02463 renka_rc *= cd22;
02464
02465 int i;
02466 for (i = 0; i < aPixgrid->size_x; i++) {
02467 int j;
02468 for (j = 0; j < aPixgrid->size_z; j++) {
02469 double sumdata = 0, sumweight = 0,
02470 lambda = j * cd22 + crval2;
02471 int npoints = 0;
02472 #ifdef ESO_ENABLE_DEBUG
02473 #define MAX_NPIX_2D 50
02474 int *pointlist = NULL;
02475 double *pointweights = NULL;
02476 if (debug & 128 && i+1 == debugx && j+1 == debugz) {
02477 pointlist = cpl_calloc(MAX_NPIX_2D, sizeof(int));
02478 pointweights = cpl_malloc(MAX_NPIX_2D * sizeof(double));
02479 }
02480 #endif
02481
02482
02483 int j2;
02484 for (j2 = j - ld; j2 <= j + ld; j2++) {
02485 cpl_size idx = muse_pixgrid_get_index(aPixgrid, i, 0, j2, CPL_FALSE),
02486 n, n_rows = muse_pixgrid_get_count(aPixgrid, idx);
02487 const cpl_size *rows = muse_pixgrid_get_rows(aPixgrid, idx);
02488 for (n = 0; n < n_rows; n++) {
02489 double dlambda = fabs(lambda - lbda[rows[n]]),
02490 weight = muse_resampling_weight_function_renka(dlambda,
02491 renka_rc);
02492 sumweight += weight;
02493 sumdata += data[rows[n]] * weight;
02494 npoints++;
02495 #ifdef ESO_ENABLE_DEBUG
02496 if (debug & 128 && i+1 == debugx && j+1 == debugz &&
02497 npoints-1 < MAX_NPIX_2D) {
02498
02499 pointlist[npoints-1] = rows[n] + 1;
02500 pointweights[npoints-1] = weight;
02501 }
02502
02503 if (debug & 256) {
02504 printf(" pixel %4d,%4d (%8"CPL_SIZE_FORMAT"): %2"CPL_SIZE_FORMAT
02505 " %2"CPL_SIZE_FORMAT" %f, %f -> %f ==> %f %f (%d)\n",
02506 i, j2, idx, n, muse_pixgrid_get_count(aPixgrid, idx), dlambda,
02507 data[muse_pixgrid_get_rows(aPixgrid, idx)[n]],
02508 weight, sumweight, sumdata, npoints);
02509 }
02510 #endif
02511 }
02512 }
02513
02514 #ifdef ESO_ENABLE_DEBUG
02515 if (debug & 128 && i+1 == debugx && j+1 == debugz) {
02516 printf("pixelnumber weight ");
02517 muse_pixtable_dump(aPixtable, 0, 0, 2);
02518 int m = -1;
02519 while (++m < MAX_NPIX_2D && pointlist[m] != 0) {
02520
02521 printf("%12d %8.5f ", pointlist[m] - 1, pointweights[m]);
02522 muse_pixtable_dump(aPixtable, pointlist[m] - 1, 1, 0);
02523 }
02524 fflush(stdout);
02525 cpl_free(pointlist);
02526 cpl_free(pointweights);
02527 }
02528 if (debug) {
02529 fflush(stdout);
02530 }
02531 #endif
02532
02533
02534 if (!npoints) {
02535 continue;
02536 }
02537
02538 imagedata[i + j * aPixgrid->size_x] = sumdata / sumweight;
02539 }
02540 }
02541
02542 return CPL_ERROR_NONE;
02543 }
02544
02545
02562
02563 static muse_image *
02564 muse_resampling_image_selected(muse_pixtable *aPixtable,
02565 muse_resampling_type aMethod,
02566 double aDX, double aLambdaMin,
02567 double aLambdaMax, double aDLambda)
02568 {
02569 double dlambda = aDLambda;
02570 float xmin = 0.;
02571 muse_pixgrid *pixgrid = muse_pixgrid_2d_create(aPixtable->table, aDX,
02572 aLambdaMin, aLambdaMax,
02573 dlambda, &xmin);
02574
02575
02576
02577 muse_image *image = muse_image_new();
02578 image->header = cpl_propertylist_new();
02579 const char *unit = cpl_table_get_column_unit(aPixtable->table, "data");
02580 cpl_propertylist_append_string(image->header, "BUNIT", unit);
02581
02582
02583
02584 cpl_propertylist_copy_property_regexp(image->header, aPixtable->header,
02585 MUSE_HDR_PT_REGEXP"|"MUSE_WCS_KEYS, 1);
02586 cpl_propertylist_append_double(image->header, "CRPIX1", 1.);
02587 cpl_propertylist_append_double(image->header, "CRPIX2", 1.);
02588 cpl_propertylist_append_double(image->header, "CRVAL1", 1.);
02589 cpl_propertylist_append_double(image->header, "CRVAL2", aLambdaMin);
02590 cpl_propertylist_append_double(image->header, "CD1_1", 1.);
02591 cpl_propertylist_append_double(image->header, "CD2_2", dlambda);
02592
02593 cpl_propertylist_append_string(image->header, "CUNIT1", "pixel");
02594 cpl_propertylist_append_string(image->header, "CUNIT2", "Angstrom");
02595
02596
02597 cpl_propertylist_append_string(image->header, "CTYPE1", "PIXEL");
02598
02599
02600 cpl_propertylist_append_string(image->header, "CTYPE2", "AWAV");
02601
02602 cpl_propertylist_append_double(image->header, "CD1_2", 0.);
02603 cpl_propertylist_append_double(image->header, "CD2_1", 0.);
02604
02605
02606 cpl_error_code rc = CPL_ERROR_NONE;
02607 switch (aMethod) {
02608 case MUSE_RESAMPLE_NEAREST:
02609 rc = muse_resampling_image_nearest(image, aPixtable, pixgrid);
02610 break;
02611 default:
02612
02613 rc = muse_resampling_image_weighted(image, aPixtable, pixgrid);
02614 }
02615
02616 muse_pixgrid_delete(pixgrid);
02617 if (rc != CPL_ERROR_NONE) {
02618 cpl_msg_error(__func__, "resampling failed: %s", cpl_error_get_message());
02619 muse_image_delete(image);
02620 return NULL;
02621 }
02622
02623 return image;
02624 }
02625
02626
02666
02667 muse_image *
02668 muse_resampling_image(muse_pixtable *aPixtable, muse_resampling_type aMethod,
02669 double aDX, double aDLambda)
02670 {
02671 cpl_ensure(aPixtable, CPL_ERROR_NULL_INPUT, NULL);
02672 if (aDLambda == 0.0) {
02673 aDLambda = kMuseSpectralSamplingA;
02674 }
02675 muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
02676 cpl_ensure(wcstype == MUSE_PIXTABLE_WCS_CELSPH ||
02677 wcstype == MUSE_PIXTABLE_WCS_PIXEL, CPL_ERROR_UNSUPPORTED_MODE,
02678 NULL);
02679
02680
02681 switch (aMethod) {
02682 case MUSE_RESAMPLE_NEAREST:
02683 cpl_msg_info(__func__, "Using nearest neighbor sampling (%d) along "
02684 "wavelengths.", aMethod);
02685 break;
02686 case MUSE_RESAMPLE_WEIGHTED_RENKA:
02687 cpl_msg_info(__func__, "Using renka-weighted interpolation (%d) along "
02688 "wavelengths.", aMethod);
02689 break;
02690 default:
02691 cpl_msg_error(__func__, "Don't know this resampling method: %d", aMethod);
02692 cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
02693 return NULL;
02694 }
02695
02696
02697
02698 float lmin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LLO),
02699 lmax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LHI);
02700
02701 if (muse_pixtable_get_type(aPixtable) == MUSE_PIXTABLE_TYPE_SIMPLE) {
02702
02703
02704
02705 return muse_resampling_image_selected(aPixtable, aMethod,
02706 aDX == 0.0 ? 1. : aDX,
02707 lmin, lmax, aDLambda);
02708 }
02709
02710
02711
02712
02713 muse_pixtable **slice_pixtable = muse_pixtable_extracted_get_slices(aPixtable);
02714 int n_slices = muse_pixtable_extracted_get_size(slice_pixtable);
02715
02716 double dx = aDX;
02717 if (dx == 0.0) {
02718 if (muse_pixtable_wcs_check(aPixtable) == MUSE_PIXTABLE_WCS_PIXEL) {
02719 dx = 1.0;
02720 } else {
02721 double xsc, ysc;
02722 muse_wcs_get_scales(aPixtable->header, &xsc, &ysc);
02723
02724 dx = 1.20 * xsc;
02725 }
02726 }
02727 cpl_msg_debug(__func__, "Resampling %d slices to a 2D image, using bins of"
02728 " %e %s x %.3f Angstrom", n_slices, dx,
02729 cpl_table_get_column_unit(aPixtable->table, MUSE_PIXTABLE_XPOS),
02730 aDLambda);
02731
02732
02733 muse_image *slice_image[n_slices];
02734 int i_slice;
02735 #pragma omp parallel for default(none) \
02736 shared(aDLambda, aMethod, dx, lmax, lmin, n_slices, \
02737 slice_image, slice_pixtable)
02738 for (i_slice = 0; i_slice < n_slices; i_slice++) {
02739 #if 0
02740 cpl_msg_debug(__func__, "slice pixel table %d", i_slice + 1);
02741 #endif
02742 muse_pixtable *pt = slice_pixtable[i_slice];
02743 if (muse_pixtable_get_nrow(pt) < 1) {
02744 slice_image[i_slice] = NULL;
02745 continue;
02746 }
02747 slice_image[i_slice] = muse_resampling_image_selected(pt, aMethod, dx,
02748 lmin, lmax, aDLambda);
02749 }
02750
02751
02752 muse_image *image = muse_image_new();
02753 for (i_slice = 0; i_slice < n_slices; i_slice++) {
02754 muse_image *slice = slice_image[i_slice];
02755 #if 0
02756 cpl_msg_debug(__func__, "slice %d: %ldx%ld", i_slice + 1,
02757 cpl_image_get_size_x(slice->data), cpl_image_get_size_y(slice->data));
02758 #endif
02759 if (slice == NULL) {
02760 continue;
02761 }
02762 if (!image->header) {
02763
02764 image->header = cpl_propertylist_duplicate(slice->header);
02765 }
02766 cpl_image *data = muse_cplimage_concat_x(image->data, slice->data);
02767 cpl_image_delete(image->data);
02768 image->data = data;
02769 if (slice->dq) {
02770 cpl_image *dq = muse_cplimage_concat_x(image->dq, slice->dq);
02771 cpl_image_delete(image->dq);
02772 image->dq = dq;
02773 }
02774 if (slice->stat) {
02775 cpl_image *stat = muse_cplimage_concat_x(image->stat, slice->stat);
02776 cpl_image_delete(image->stat);
02777 image->stat = stat;
02778 }
02779 muse_image_delete(slice);
02780 slice_image[i_slice] = NULL;
02781 }
02782 muse_pixtable_extracted_delete(slice_pixtable);
02783
02784 return image;
02785 }
02786
02787
02803
02804 cpl_table *
02805 muse_resampling_spectrum(muse_pixtable *aPixtable, double aBinwidth)
02806 {
02807 cpl_ensure(aPixtable && aPixtable->header && aPixtable->table,
02808 CPL_ERROR_NULL_INPUT, NULL);
02809 cpl_ensure(muse_cpltable_check(aPixtable->table, muse_pixtable_def)
02810 == CPL_ERROR_NONE, CPL_ERROR_ILLEGAL_INPUT, NULL);
02811
02812
02813 double lmin = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LLO);
02814 double lmax = cpl_propertylist_get_float(aPixtable->header, MUSE_HDR_PT_LHI);
02815 cpl_size lsize = floor((lmax - lmin) / aBinwidth) + 2;
02816
02817
02818 cpl_table *spectrum = muse_cpltable_new(muse_dataspectrum_def, lsize);
02819 cpl_table_fill_column_window(spectrum, "lambda", 0, lsize, 0);
02820 cpl_table_fill_column_window(spectrum, "data", 0, lsize, 0);
02821 cpl_table_fill_column_window(spectrum, "stat", 0, lsize, 0);
02822 cpl_table_fill_column_window(spectrum, "dq", 0, lsize, 0);
02823 double *spec_data = cpl_table_get_data_double(spectrum, "data");
02824 double *spec_stat = cpl_table_get_data_double(spectrum, "stat");
02825 double *spec_lbda = cpl_table_get_data_double(spectrum, "lambda");
02826
02827 cpl_table_set_column_unit(spectrum, "data",
02828 cpl_table_get_column_unit(aPixtable->table,
02829 MUSE_PIXTABLE_DATA));
02830 cpl_table_set_column_unit(spectrum, "stat",
02831 cpl_table_get_column_unit(aPixtable->table,
02832 MUSE_PIXTABLE_STAT));
02833
02834
02835 cpl_table_new_column(spectrum, "weight", CPL_TYPE_DOUBLE);
02836 cpl_table_fill_column_window(spectrum, "weight", 0, lsize, 0);
02837 double *spec_weight = cpl_table_get_data_double(spectrum, "weight");
02838
02839
02840 float *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA);
02841 float *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA);
02842 float *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
02843 float *wght = NULL;
02844 if (cpl_table_has_column(aPixtable->table, MUSE_PIXTABLE_WEIGHT)) {
02845 wght = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_WEIGHT);
02846 }
02847 int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
02848
02849
02850
02851 cpl_array *asel = cpl_table_where_selected(aPixtable->table);
02852 const cpl_size *sel = cpl_array_get_data_cplsize_const(asel);
02853 cpl_size isel, nsel = cpl_array_get_size(asel);
02854 for (isel = 0; isel < nsel; isel++) {
02855 cpl_size i_row = sel[isel];
02856 if ((dq[i_row] != EURO3D_GOODPIXEL)) {
02857 continue;
02858 }
02859
02860
02861 if (! isfinite(data[i_row])) {
02862 continue;
02863 }
02864
02865 double l = (lbda[i_row] - lmin) / aBinwidth;
02866 if (l < 0) l = 0;
02867 cpl_size il = floor(l);
02868 l -= il;
02869 double w0 = 1-l;
02870 double w1 = l;
02871 if (wght) {
02872 w0 *= wght[i_row];
02873 w1 *= wght[i_row];
02874 }
02875
02876 spec_data[il] += w0 * data[i_row];
02877 spec_data[il+1] += w1 * data[i_row];
02878 spec_stat[il] += w0 * stat[i_row];
02879 spec_stat[il+1] += w1 * stat[i_row];
02880 spec_weight[il] += w0;
02881 spec_weight[il+1] += w1;
02882 }
02883 cpl_array_delete(asel);
02884
02885
02886 cpl_size ispec;
02887 for (ispec = 0; ispec < lsize; ispec++) {
02888 if (spec_weight[ispec] > 0) {
02889 spec_lbda[ispec] = ispec * aBinwidth + lmin;
02890 cpl_table_unselect_row(spectrum, ispec);
02891 }
02892 }
02893 cpl_table_erase_selected(spectrum);
02894
02895
02896 cpl_table_divide_columns(spectrum, "data", "weight");
02897 cpl_table_divide_columns(spectrum, "stat", "weight");
02898 cpl_table_erase_column(spectrum, "weight");
02899 return spectrum;
02900 }
02901
02902
02933
02934 cpl_table *
02935 muse_resampling_spectrum_iterate(muse_pixtable *aPixtable, double aBinwidth,
02936 float aLo, float aHi, unsigned char aIter)
02937 {
02938 cpl_ensure(aPixtable && aPixtable->header && aPixtable->table,
02939 CPL_ERROR_NULL_INPUT, NULL);
02940 cpl_ensure(muse_cpltable_check(aPixtable->table, muse_pixtable_def)
02941 == CPL_ERROR_NONE, CPL_ERROR_ILLEGAL_INPUT, NULL);
02942
02943
02944 cpl_table *spectrum = muse_resampling_spectrum(aPixtable, aBinwidth);
02945 if (aIter == 0) {
02946 return spectrum;
02947 }
02948
02949
02950 float *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
02951 *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA);
02952 int *dq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
02953
02954 cpl_array *asel = cpl_table_where_selected(aPixtable->table);
02955 const cpl_size *sel = cpl_array_get_data_cplsize_const(asel);
02956 cpl_size isel, nsel = cpl_array_get_size(asel);
02957
02958
02959 cpl_size nhi = 0, nlo = 0;
02960 unsigned char i;
02961 for (i = 1; i <= aIter; i++) {
02962 cpl_size ispec, nbins = cpl_table_get_nrow(spectrum);
02963 double *sdata = cpl_table_get_data_double(spectrum, "data"),
02964 *sstat = cpl_table_get_data_double(spectrum, "stat"),
02965 *sstddev = cpl_malloc(nbins * sizeof(double));
02966
02967 for (ispec = 0; ispec < nbins; ispec++) {
02968 sstddev[ispec] = sqrt(sstat[ispec]);
02969 }
02970
02971
02972 for (isel = 0; isel < nsel; isel++) {
02973 cpl_size irow = sel[isel];
02974 if ((dq[irow] != EURO3D_GOODPIXEL)) {
02975 continue;
02976 }
02977 double l = lbda[irow];
02978 ispec = muse_cpltable_find_sorted(spectrum, "lambda", l);
02979
02980 if ((ispec < nbins - 1) && (sdata[ispec] < sdata[ispec + 1])) {
02981 ispec++;
02982 }
02983
02984 if (aHi > 0. && data[irow] > sdata[ispec] + aHi * sstddev[ispec]) {
02985 dq[irow] = EURO3D_OUTLIER;
02986 nhi++;
02987 }
02988
02989 if (aLo > 0. && data[irow] < sdata[ispec] - aLo * sstddev[ispec]) {
02990 dq[irow] = EURO3D_OUTLIER;
02991 nlo++;
02992 }
02993 }
02994 cpl_free(sstddev);
02995
02996
02997
02998 cpl_msg_debug(__func__, "%"CPL_SIZE_FORMAT" of %"CPL_SIZE_FORMAT" pixels "
02999 "are outliers (%"CPL_SIZE_FORMAT" low and %"CPL_SIZE_FORMAT
03000 " high, after %hhu iteration%s)", nlo + nhi, nsel, nlo, nhi,
03001 i, i == 1 ? "" : "s");
03002
03003 cpl_table_delete(spectrum);
03004 spectrum = muse_resampling_spectrum(aPixtable, aBinwidth);
03005 }
03006 cpl_array_delete(asel);
03007
03008
03009 muse_pixtable_reset_dq(aPixtable, EURO3D_OUTLIER);
03010
03011 return spectrum;
03012 }
03013