33 #define omp_get_max_threads() 1
38 #include "muse_resampling.h"
39 #include "muse_instrument.h"
41 #include "muse_astro.h"
42 #include "muse_cplwrappers.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"
52 #include "muse_data_format_z.h"
120 const cpl_propertylist *aWCS)
122 cpl_ensure_code(aParams, CPL_ERROR_NULL_INPUT);
125 cpl_wcs_delete(aParams->
wcs);
127 return CPL_ERROR_NONE;
129 if (cpl_propertylist_has(aWCS,
"CTYPE3") &&
130 !strncmp(cpl_propertylist_get_string(aWCS,
"CTYPE3"),
"AWAV-LOG", 9)) {
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);
141 rc = cpl_error_get_code();
144 printf(
"CDi_j matrix:\n");
145 cpl_matrix_dump(cpl_wcs_get_cd(aParams->
wcs), stdout);
166 cpl_wcs_delete(aParams->
wcs);
183 muse_resampling_weight_function_linear(
double r)
185 return r == 0 ? FLT_MAX : 1. / r;
201 muse_resampling_weight_function_quadratic(
double r2)
203 return r2 == 0 ? FLT_MAX : 1. / r2;
221 muse_resampling_weight_function_renka(
double r,
double r_c)
225 }
else if (r >= r_c) {
228 double p = (r_c - r) / (r_c * r);
242 muse_resampling_weight_function_sinc(
double r)
244 return fabs(r) < DBL_EPSILON ? 1. : sin(CPL_MATH_PI * r) / (CPL_MATH_PI * r);
259 muse_resampling_weight_function_lanczos(
double dx,
double dy,
double dz,
unsigned int n)
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);
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)
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;
300 if (x <= 0 || y <= 0 || z <= 0) {
305 return (x > xin ? xin : x) * (y > yin ? yin : y) * (z > zin ? zin : z)
326 static cpl_error_code
330 cpl_ensure_code(aCube && aPixtable && aPixgrid, CPL_ERROR_NULL_INPUT);
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");
338 cpl_boolean loglambda = !strncmp(cpl_propertylist_get_string(aCube->
header,
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);
352 double xnorm = 1., ynorm = 1., znorm = 1. / kMuseSpectralSamplingA;
358 ptxoff = cpl_propertylist_get_double(aPixtable->
header,
"CRVAL1");
359 ptyoff = cpl_propertylist_get_double(aPixtable->
header,
"CRVAL2");
362 #ifdef ESO_ENABLE_DEBUG
364 if (getenv(
"MUSE_DEBUG_NEAREST")) {
365 debug = atoi(getenv(
"MUSE_DEBUG_NEAREST"));
370 #ifdef ESO_ENABLE_DEBUG
371 #pragma omp parallel for default(none) \
372 shared(aCube, aPixgrid, cd33, crpix3, crval3, data, debug, dq, lbda, \
373 loglambda, ptxoff, ptyoff, stat, stdout, wcs, xnorm, xpos, \
376 #pragma omp parallel for default(none) \
377 shared(aCube, aPixgrid, cd33, crpix3, crval3, data, dq, lbda, \
378 loglambda, ptxoff, ptyoff, stat, stdout, wcs, xnorm, xpos, \
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));
386 double lambda = (l + 1. - crpix3) * cd33 + crval3;
388 lambda = crval3 * exp((l + 1. - crpix3) * cd33 / crval3);
392 for (i = 0; i < aPixgrid->size_x; i++) {
394 for (j = 0; j < aPixgrid->size_y; j++) {
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]];
413 #ifdef ESO_ENABLE_DEBUG
415 printf(
"only: %f,%f\n", data[rows[0]], stat[rows[0]]);
419 }
else if (n_rows >= 2) {
421 cpl_size n, nbest = -1;
422 double dbest = FLT_MAX;
423 for (n = 0; n < n_rows; n++) {
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);
432 dx *= cos(y * CPL_MATH_RAD_DEG);
434 #ifdef ESO_ENABLE_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]],
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
452 printf(
"closest: %f,%f\n", data[rows[nbest]], stat[rows[nbest]]);
458 pdq[i + j * aPixgrid->size_x] = EURO3D_MISSDATA;
465 return CPL_ERROR_NONE;
490 static cpl_error_code
495 cpl_ensure_code(aCube && aPixtable && aPixgrid && aParams,
496 CPL_ERROR_NULL_INPUT);
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");
504 cpl_boolean loglambda = !strncmp(cpl_propertylist_get_string(aCube->
header,
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);
517 int *dq = cpl_table_get_data_int(aPixtable->
table, MUSE_PIXTABLE_DQ);
523 double xnorm = 1., ynorm = 1., znorm = 1. / kMuseSpectralSamplingA;
529 ptxoff = cpl_propertylist_get_double(aPixtable->
header,
"CRVAL1");
530 ptyoff = cpl_propertylist_get_double(aPixtable->
header,
"CRVAL2");
534 double zoutefac = exp(1.5 * cd33 / crval3) - exp(0.5 * cd33 / crval3);
536 double renka_rc = aParams->
rc
537 * sqrt((wcs->cd11*xnorm)*(wcs->cd11*xnorm) + (wcs->
cd22*ynorm)*(wcs->
cd22*ynorm)
538 + (cd33*znorm)*(cd33*znorm));
540 int ld = aParams->
ld;
543 cpl_msg_info(__func__,
"Overriding loop distance ld=%d", ld);
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);
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"));
559 if (getenv(
"MUSE_DEBUG_WEIGHTED_X")) {
560 debugx = atoi(getenv(
"MUSE_DEBUG_WEIGHTED_X"));
561 if (debugx < 1 || debugx > aPixgrid->size_x) {
565 if (getenv(
"MUSE_DEBUG_WEIGHTED_Y")) {
566 debugy = atoi(getenv(
"MUSE_DEBUG_WEIGHTED_Y"));
567 if (debugy < 1 || debugy > aPixgrid->size_y) {
571 if (getenv(
"MUSE_DEBUG_WEIGHTED_Z")) {
572 debugz = atoi(getenv(
"MUSE_DEBUG_WEIGHTED_Z"));
573 if (debugz < 1 || debugz > aPixgrid->size_z) {
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);
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,
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);
596 cpl_imagelist *wcube = NULL;
597 if (getenv(
"MUSE_DEBUG_WEIGHT_CUBE")) {
598 cpl_msg_debug(__func__,
"Weighted resampling: creating weight cube");
599 wcube = cpl_imagelist_new();
601 for (i = 0; i < aPixgrid->size_z; i++) {
602 cpl_image *image = cpl_image_new(aPixgrid->size_x, aPixgrid->size_y,
604 cpl_imagelist_set(wcube, image, i);
608 if (getenv(
"MUSE_DEBUG_WEIGHTED_GRID")) {
609 char *fn = getenv(
"MUSE_DEBUG_WEIGHTED_GRID");
610 FILE *grid = fopen(fn,
"w");
612 cpl_msg_info(__func__,
"Writing grid to \"%s\"", fn);
613 fprintf(grid,
"# i j l x y lambda\n");
615 for (l = 0; l < aPixgrid->size_z; l++) {
616 double lambda = (l + 1. - crpix3) * cd33 + crval3;
618 for (i = 0; i < aPixgrid->size_x; i++) {
620 for (j = 0; j < aPixgrid->size_y; j++) {
628 fprintf(grid,
"%03d %03d %04d %.10f %.10f %8.3f\n", i+1, j+1, l+1,
635 cpl_msg_warning(__func__,
"Writing grid to \"%s\" failed!", fn);
640 #ifdef ESO_ENABLE_DEBUG
641 #pragma omp parallel for default(none) \
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, \
648 #pragma omp parallel for default(none) \
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)
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));
659 double lambda = (l + 1. - crpix3) * cd33 + crval3;
662 lambda = crval3 * exp((l + 1. - crpix3) * cd33 / crval3);
663 zout2 = crval3 * exp((l - crpix3) * cd33 / crval3) * zoutefac;
665 float *pwcube = NULL;
667 pwcube = cpl_image_get_data_float(cpl_imagelist_get(wcube, l));
671 for (i = 0; i < aPixgrid->size_x; i++) {
673 for (j = 0; j < aPixgrid->size_y; j++) {
681 double sumdata = 0, sumstat = 0, sumweight = 0;
683 #ifdef ESO_ENABLE_DEBUG
684 cpl_size *pointlist = NULL;
685 double *pointweights = NULL;
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));
694 #ifdef ESO_ENABLE_DEBUG
696 printf(
"cell %d %d %d:\n", i, j, l);
701 for (i2 = i - ld; i2 <= i + ld; i2++) {
703 for (j2 = j - ld; j2 <= j + ld; j2++) {
705 for (l2 = l - ld; l2 <= l + ld; l2++) {
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);
716 for (n = 0; n < n_rows2; n++) {
718 #ifdef ESO_ENABLE_DEBUG
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);
729 double dx = fabs(x - (xpos[rows2[n]] + ptxoff)),
730 dy = fabs(y - (ypos[rows2[n]] + ptyoff)),
731 dlambda = fabs(lambda - lbda[rows2[n]]),
739 dx *= cos(y * CPL_MATH_RAD_DEG);
745 r2 = dx*dx + dy*dy + dlambda*dlambda;
749 weight = muse_resampling_weight_function_renka(sqrt(r2), renka_rc);
751 weight = muse_resampling_weight_function_drizzle(xsz, ysz, zsz,
755 weight = muse_resampling_weight_function_linear(sqrt(r2));
757 weight = muse_resampling_weight_function_quadratic(r2);
759 weight = muse_resampling_weight_function_lanczos(dx, dy, dlambda, ld);
764 weight *= wght[rows2[n]];
766 #ifdef ESO_ENABLE_DEBUG
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);
779 sumdata += data[rows2[n]] * weight;
780 sumstat += stat[rows2[n]] * weight*weight;
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));
793 pointlist[npoints-1] = rows2[n] + 1;
794 pointweights[npoints-1] = weight;
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,
806 weight, sumweight, sumdata, npoints);
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);
821 while (++m < npointlist && pointlist[m] != 0) {
823 printf(
"%12"CPL_SIZE_FORMAT
" %8.5f ", pointlist[m] - 1,
827 printf(
"sums: %g %g %g --> %g %g\n", sumdata, sumstat, sumweight,
828 sumdata / sumweight, sumstat / pow(sumweight, 2));
830 cpl_free(pointweights);
833 if (debug & 1 && npoints) {
834 printf(
" sumdata: %e %e (%d)", sumdata, sumweight, npoints);
841 if (!npoints || !isnormal(sumweight)) {
842 pdq[i + j * aPixgrid->size_x] = EURO3D_MISSDATA;
843 #ifdef ESO_ENABLE_DEBUG
845 printf(
" -> no points or weird weight\n");
852 sumdata /= sumweight;
853 sumstat /= sumweight*sumweight;
854 #ifdef ESO_ENABLE_DEBUG
856 printf(
" -> %e (%e)\n", sumdata, sumstat);
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;
866 pwcube[i + j * aPixgrid->size_x] = sumweight;
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());
881 cpl_msg_info(__func__,
"Saved weight cube as \"%s\"", fn);
883 cpl_imagelist_delete(wcube);
886 return CPL_ERROR_NONE;
900 static cpl_error_code
904 cpl_ensure_code(aPixtable && aParams, CPL_ERROR_NULL_INPUT);
905 const char func[] =
"muse_resampling_cube";
908 if (aParams->dlambda == 0.0) {
909 aParams->dlambda = kMuseSpectralSamplingA;
911 aParams->dlambda /= 1.6;
917 if (aParams->
dx == 0.0) {
920 if (aParams->dy == 0.0) {
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;
927 if (aParams->
dx == 0.0) {
928 aParams->
dx = kMuseSpaxelSizeX_WFM / 3600.;
930 aParams->
dx = kMuseSpaxelSizeX_NFM / 3600.;
934 aParams->
dx /= 3600.;
936 if (aParams->dy == 0.0) {
937 aParams->dy = kMuseSpaxelSizeY_WFM / 3600.;
939 aParams->dy = kMuseSpaxelSizeY_NFM / 3600.;
943 aParams->dy /= 3600.;
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;
965 static cpl_error_code
968 int *aX,
int *aY,
int *aZ)
970 cpl_ensure_code(aPixtable && aParams && aX && aY && aZ, CPL_ERROR_NULL_INPUT);
971 const char func[] =
"muse_resampling_cube";
972 double x1, y1, x2, y2;
985 *aX = lround(fabs(x2 - x1) / aParams->
dx) + 1;
986 *aY = lround(fabs(y2 - y1) / aParams->dy) + 1;
989 *aZ = (int)ceil((lmax - lmin) / aParams->dlambda) + 1;
991 *aZ = (int)ceil(lmin / aParams->dlambda * log(lmax / lmin)) + 1;
993 cpl_msg_info(func,
"Output cube size %d x %d x %d (fit to data)",
995 return CPL_ERROR_NONE;
1010 muse_resampling_override_size_int(
int *aV,
const char *aKeyword,
int aValue)
1012 if (aValue <= 0 || !aV) {
1015 const char func[] =
"muse_resampling_cube";
1016 cpl_msg_info(func,
"Overriding %s=%d (was %d)", aKeyword, aValue, *aV);
1035 static cpl_error_code
1036 muse_resampling_override_size(
int *aX,
int *aY,
int *aZ,
1039 cpl_ensure_code(aX && aY && aZ && aParams, CPL_ERROR_NULL_INPUT);
1040 if (!aParams->
wcs) {
1041 return CPL_ERROR_NONE;
1043 const char func[] =
"muse_resampling_cube";
1045 const cpl_array *dims = cpl_wcs_get_image_dims(aParams->
wcs);
1047 cpl_msg_debug(func,
"No dimensions to override were specified");
1048 return CPL_ERROR_NONE;
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));
1056 return CPL_ERROR_NONE;
1072 muse_resampling_override_wcs_double(
muse_datacube *aCube,
const char *aKeyword,
1073 double aValue,
int aError)
1075 if (aError || !aCube) {
1076 cpl_msg_debug(
"double",
"%s=%#g (%d)", aKeyword, aValue, aError);
1079 const char func[] =
"muse_resampling_cube";
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);
1086 cpl_propertylist_update_bool(aCube->
header,
"MUSE_RESAMPLING_WCS_OVERRIDDEN",
1103 static cpl_error_code
1107 cpl_ensure_code(aCube && aCube->
header && aParams, CPL_ERROR_NULL_INPUT);
1108 if (!aParams->
wcs) {
1109 return CPL_ERROR_NONE;
1111 const char func[] =
"muse_resampling_cube";
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);
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);
1121 cpl_msg_debug(func,
"No CRVALj to override were specified");
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);
1127 cpl_msg_debug(func,
"No CRPIXi to override were specified");
1130 int naxes = cpl_matrix_get_ncol(cd);
1131 if (cpl_matrix_get_determinant(cd) == 0.) {
1132 cpl_msg_warning(func,
"det(CDi_j) = 0, not overriding CDi_j!");
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!");
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));
1156 cpl_msg_debug(func,
"No CDi_j to override were specified");
1161 if (cpl_array_get_size(crval) > 2) {
1163 muse_resampling_override_wcs_double(aCube,
"CRVAL3",
1164 cpl_array_get_double(crval, 2, &err) * 1e10, err);
1167 muse_resampling_override_wcs_double(aCube,
"CRPIX3", cpl_array_get_double(crpix, 2, &err), err);
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));
1175 return CPL_ERROR_NONE;
1189 static cpl_error_code
1192 cpl_ensure_code(aCube && aPixtable, CPL_ERROR_NULL_INPUT);
1193 const char func[] =
"muse_resampling_cube";
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;
1204 double xoff1, yoff1, xoff2, yoff2;
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);
1228 cpl_propertylist_update_double(aCube->
header,
"CRPIX1", crpix1);
1229 cpl_propertylist_update_double(aCube->
header,
"CRPIX2", crpix2);
1231 return CPL_ERROR_NONE;
1236 const char *crtypestring[] = {
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), \
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) { \
1250 } else if (dq[rows[irow]]) { \
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]]); \
1258 #define MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END \
1259 if (aType == MUSE_RESAMPLING_CRSTATS_IRAF) { \
1261 level += data[rows[irow]]; \
1262 dev += stat[rows[irow]]; \
1265 if (npixels >= sxsize) { \
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); \
1272 vdata = cpl_image_get_data_double(sdata); \
1274 vdata[npixels] = data[rows[irow]]; \
1275 vdata[npixels + sxsize] = stat[rows[irow]]; \
1309 const char *
id =
"muse_resampling_cube";
1311 cpl_msg_warning(
id,
"Unknown type (%u) for cosmic-ray rejection statistics,"
1312 " resetting to \"%s\"", aType,
1316 cpl_msg_info(
id,
"Using %s statistics (%.3f sigma level) for cosmic-ray"
1317 " rejection", crtypestring[aType], aSigma);
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"));
1325 if (getenv(
"MUSE_DEBUG_CRREJECT_X")) {
1326 debugx = atoi(getenv(
"MUSE_DEBUG_CRREJECT_X"));
1327 if (debugx < 1 || debugx > aPixgrid->size_x) {
1331 if (getenv(
"MUSE_DEBUG_CRREJECT_Y")) {
1332 debugy = atoi(getenv(
"MUSE_DEBUG_CRREJECT_Y"));
1333 if (debugy < 1 || debugy > aPixgrid->size_y) {
1337 if (getenv(
"MUSE_DEBUG_CRREJECT_Z")) {
1338 debugz = atoi(getenv(
"MUSE_DEBUG_CRREJECT_Z"));
1339 if (debugz < 1 || debugz > aPixgrid->size_z) {
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);
1355 enum dirtype dir = DIR_NONE;
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");
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");
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");
1386 unsigned badinclude = EURO3D_COSMICRAY;
1388 #ifdef ESO_ENABLE_DEBUG
1389 #pragma omp parallel for default(none) \
1390 shared(aPixgrid, aPixtable, aSigma, badinclude, aType, crtypestring, \
1391 data, debug, debugx, debugy, debugz, dir, dq, id, stat, stdout)
1393 #pragma omp parallel for default(none) \
1394 shared(aPixgrid, aSigma, aType, badinclude, crtypestring, data, dir, \
1397 for (l = 0; l < aPixgrid->size_z; l++) {
1398 #define MUSE_CRREJECT_MAX_NPIX 1000
1399 int sxsize = MUSE_CRREJECT_MAX_NPIX;
1400 cpl_image *sdata = NULL;
1401 double *vdata = NULL;
1405 sdata = cpl_image_new(sxsize, 2, CPL_TYPE_DOUBLE);
1406 vdata = cpl_image_get_data_double(sdata);
1410 for (i = 0; i < aPixgrid->size_x; i++) {
1412 for (j = 0; j < aPixgrid->size_y; j++) {
1415 double level = 0., dev = 0;
1418 cpl_size npixels = 0;
1421 for (i2 = i - CR_LD; i2 <= i + CR_LD; i2++) {
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
1432 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END
1437 for (j2 = j - CR_LD; j2 <= j + CR_LD; j2++) {
1439 if (dir == DIR_X && j2 == j) {
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
1448 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END
1459 dev = sqrt(dev) / npixels;
1463 sflags = CPL_STATS_MEDIAN | CPL_STATS_MEDIAN_DEV;
1465 sflags = CPL_STATS_MEAN | CPL_STATS_STDEV;
1467 cpl_stats *sstats = cpl_stats_new_from_image_window(sdata, sflags,
1470 level = cpl_stats_get_median(sstats);
1471 dev = cpl_stats_get_median_dev(sstats);
1473 level = cpl_stats_get_mean(sstats);
1474 dev = cpl_stats_get_stdev(sstats);
1476 cpl_stats_delete(sstats);
1478 double limit = level + aSigma * dev;
1479 #ifdef ESO_ENABLE_DEBUG
1480 if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
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);
1486 cpl_stats *ssdata = cpl_stats_new_from_image_window(sdata, CPL_STATS_ALL,
1488 *ssstat = cpl_stats_new_from_image_window(sdata, CPL_STATS_ALL,
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);
1497 cpl_stats_delete(ssdata);
1498 cpl_stats_delete(ssstat);
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);
1521 #ifdef ESO_ENABLE_DEBUG
1528 cpl_image_delete(sdata);
1560 cpl_ensure(aParams, CPL_ERROR_NULL_INPUT, NULL);
1572 e3d->
header = cpl_propertylist_duplicate(cube->
header);
1573 cpl_propertylist_erase_regexp(e3d->
header,
1574 "^SIMPLE$|^BITPIX$|^NAXIS|^EURO3D$|^E3D_",
1576 cpl_propertylist_append_char(e3d->
header,
"EURO3D",
'T');
1577 cpl_propertylist_set_comment(e3d->
header,
"EURO3D",
1578 "file conforms to Euro3D standard");
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)) {
1588 cpl_errorstate_set(prestate);
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");
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");
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"));
1622 cpl_propertylist_erase_regexp(e3d->
header,
"^C...*3$|^CD3_.$", 0);
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);
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);
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));
1643 cpl_table_set_column_savetype(e3d->
dtable,
"SELECTED", CPL_TYPE_BOOL);
1648 cpl_table_set_column_unit(e3d->
dtable,
"XPOS",
"deg");
1649 cpl_table_set_column_unit(e3d->
dtable,
"YPOS",
"deg");
1653 unsigned euro3d_ignore = EURO3D_OUTSDRANGE | EURO3D_MISSDATA;
1654 cpl_vector *vusable = cpl_vector_new(nx * ny);
1656 for (i = 1; i <= nx; i++) {
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);
1670 int l, nusable = 0, nstart = -1;
1671 for (l = 0; l < nz; l++) {
1673 unsigned dq = cpl_image_get(cpl_imagelist_get(cube->
dq, l), i, j, &err);
1675 if (nstart < 0 && (dq & euro3d_ignore)) {
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),
1683 cpl_array_set_float(astat, nusable,
1684 sqrt(cpl_image_get(cpl_imagelist_get(cube->
stat, l),
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);
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);
1703 cpl_array_delete(adata);
1704 cpl_array_delete(adq);
1705 cpl_array_delete(astat);
1707 cpl_vector_set(vusable, itable, nusable);
1708 if (nstart != -1 && nusable > 0) {
1710 cpl_table_unselect_row(e3d->
dtable, itable);
1717 cpl_vector_set_size(vusable, itable);
1718 int maxusable = cpl_vector_get_max(vusable);
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),
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) {
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);
1737 cpl_table_erase_selected(e3d->
dtable);
1738 cpl_vector_delete(vusable);
1742 cpl_table_fill_column_window_int(e3d->
dtable,
"SELECTED", 0, itable,
1745 cpl_table_fill_column_window_int(e3d->
dtable,
"NSPAX", 0, itable, 1);
1747 cpl_table_fill_column_window_int(e3d->
dtable,
"GROUP_N", 0, itable, 1);
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");
1758 cpl_table_set_int(e3d->
gtable,
"GROUP_N", 0, 1);
1759 cpl_table_set_string(e3d->
gtable,
"G_SHAPE", 0,
"R");
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.) {
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);
1773 cpl_table_set_float(e3d->
gtable,
"G_POSWAV", 0,
1774 (kMuseNominalLambdaMin + kMuseNominalLambdaMax) / 2.);
1775 cpl_table_set_float(e3d->
gtable,
"G_AIRMAS", 0,
1777 cpl_table_set_float(e3d->
gtable,
"G_PARANG", 0,
1779 cpl_table_set_float(e3d->
gtable,
"G_PRESSU", 0,
1782 cpl_table_set_float(e3d->
gtable,
"G_TEMPER", 0,
1784 cpl_table_set_float(e3d->
gtable,
"G_HUMID", 0,
1847 cpl_ensure(aPixtable && aParams, CPL_ERROR_NULL_INPUT, NULL);
1849 CPL_ERROR_ILLEGAL_INPUT, NULL);
1857 muse_resampling_check_deltas(aPixtable, aParams);
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,
1867 double time = cpl_test_get_walltime();
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,
1882 cpl_propertylist_update_string(cube->
header,
"BUNIT",
1883 cpl_table_get_column_unit(aPixtable->
table,
1884 MUSE_PIXTABLE_DATA));
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);
1893 cpl_propertylist_copy_property_regexp(cube->
header, aPixtable->
header,
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.);
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.);
1907 cpl_propertylist_erase(cube->
header,
"WCSAXES");
1910 cpl_propertylist_append_double(cube->
header,
"CRPIX1", 1.);
1911 cpl_propertylist_append_double(cube->
header,
"CRVAL1",
1912 cpl_propertylist_get_float(aPixtable->
header,
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,
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);
1927 cpl_propertylist_append_double(cube->
header,
"CD1_2", 0.);
1928 cpl_propertylist_append_double(cube->
header,
"CD2_1", 0.);
1931 cpl_propertylist_append_string(cube->
header,
"CTYPE3",
"AWAV-LOG");
1933 cpl_propertylist_append_string(cube->
header,
"CTYPE3",
"AWAV");
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,
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.);
1948 muse_resampling_override_wcs(cube, aParams);
1950 muse_resampling_fit_data_to_grid(cube, aPixtable);
1954 cube->
data = cpl_imagelist_new();
1955 cube->
dq = cpl_imagelist_new();
1956 cube->
stat = cpl_imagelist_new();
1958 for (i = 0; i < zsize; i++) {
1959 cpl_image *image = cpl_image_new(xsize, ysize, CPL_TYPE_FLOAT);
1961 cpl_imagelist_set(cube->
data, image, i);
1963 cpl_imagelist_set(cube->
stat, cpl_image_duplicate(image), i);
1967 image = cpl_image_new(xsize, ysize, CPL_TYPE_INT);
1969 cpl_image_add_scalar(image, EURO3D_OUTSDRANGE);
1970 cpl_imagelist_set(cube->
dq, image, i);
1977 xsize, ysize, zsize);
1980 muse_resampling_crreject(aPixtable, pixgrid, aParams->
crtype,
1985 double timeinit = cpl_test_get_walltime(),
1986 cpuinit = cpl_test_get_cputime();
1989 cpl_error_code rc = CPL_ERROR_NONE;
1990 switch (aParams->
method) {
1992 cpl_msg_info(__func__,
"Starting resampling, using method \"nearest\"");
1993 rc = muse_resampling_cube_nearest(cube, aPixtable, pixgrid);
1996 cpl_msg_info(__func__,
"Starting resampling, using method \"renka\" "
1997 "(critical radius rc=%f, loop distance ld=%d)", aParams->
rc,
1999 rc = muse_resampling_cube_weighted(cube, aPixtable, pixgrid, aParams);
2004 cpl_msg_info(__func__,
"Starting resampling, using method \"%s\" (loop "
2012 rc = muse_resampling_cube_weighted(cube, aPixtable, pixgrid, aParams);
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);
2021 cpl_msg_debug(__func__,
"Method %d (no resampling)", aParams->
method);
2024 cpl_msg_error(__func__,
"Don't know this resampling method: %d",
2026 cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
2027 rc = CPL_ERROR_UNSUPPORTED_MODE;
2031 double timefini = cpl_test_get_walltime(),
2032 cpufini = cpl_test_get_cputime();
2036 *aPixgrid = pixgrid;
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());
2091 cpl_ensure(aPixtable && aPixtable->
table && aPixgrid && aParams &&
2092 aCube && aCube->
header, CPL_ERROR_NULL_INPUT, NULL);
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);
2109 int *dq = cpl_table_get_data_int(aPixtable->
table, MUSE_PIXTABLE_DQ);
2115 double xnorm = 1., ynorm = 1.,
2123 ptxoff = cpl_propertylist_get_double(aPixtable->
header,
"CRVAL1");
2124 ptyoff = cpl_propertylist_get_double(aPixtable->
header,
"CRVAL2");
2127 double renka_rc = aParams->
rc
2128 * sqrt(pow(wcs->cd11*xnorm, 2) + pow(wcs->
cd22*xnorm, 2));
2130 int ld = aParams->
ld;
2133 cpl_msg_info(__func__,
"Overriding loop distance ld=%d", ld);
2137 double xsz = aParams->
pfx / xnorm,
2138 ysz = aParams->pfy / ynorm,
2139 xout = fabs(wcs->cd11), yout = fabs(wcs->
cd22);
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);
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;
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) {
2168 while (l > 0 && fabs(fresp[l]) < DBL_EPSILON) {
2171 cpl_msg_debug(__func__,
"filter wavelength range %.1f..%.1fA", lmin, lmax);
2175 cpl_msg_debug(__func__,
"full wavelength range %.1f..%.1fA", lmin, lmax);
2179 #pragma omp parallel for default(none) \
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++) {
2185 for (j = 0; j < aPixgrid->size_y; j++) {
2193 double sumdata = 0, sumstat = 0, sumweight = 0;
2194 cpl_size npoints = 0;
2198 for (i2 = i - ld; i2 <= i + ld; i2++) {
2200 for (j2 = j - ld; j2 <= j + ld; j2++) {
2203 for (l2 = 0; l2 < aPixgrid->size_z; l2++) {
2207 for (n = 0; n < n_rows2; n++) {
2211 if (usevariance && !isnormal(stat[rows2[n]])) {
2214 if (lbda[rows2[n]] > lmax || lbda[rows2[n]] < lmin) {
2218 double dx = fabs(x - (xpos[rows2[n]] + ptxoff)),
2219 dy = fabs(y - (ypos[rows2[n]] + ptyoff)),
2224 dx *= cos(y * CPL_MATH_RAD_DEG);
2233 weight = muse_resampling_weight_function_renka(sqrt(r2), renka_rc);
2235 weight = muse_resampling_weight_function_drizzle(xsz, ysz, 1.,
2239 weight = muse_resampling_weight_function_linear(sqrt(r2));
2241 weight = muse_resampling_weight_function_quadratic(r2);
2243 weight = muse_resampling_weight_function_lanczos(dx, dy, 0., ld);
2248 weight *= wght[rows2[n]];
2257 weight /= stat[rows2[n]];
2259 sumweight += weight;
2260 sumdata += data[rows2[n]] * weight;
2261 sumstat += stat[rows2[n]] * weight*weight;
2271 if (!npoints || !isnormal(sumweight)) {
2272 pdq[i + j * aPixgrid->size_x] = EURO3D_MISSDATA;
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;
2310 static cpl_error_code
2314 cpl_ensure_code(aImage && aPixtable && aPixgrid, CPL_ERROR_NULL_INPUT);
2315 aImage->
data = cpl_image_new(aPixgrid->size_x, aPixgrid->size_z,
2318 double crval2 = cpl_propertylist_get_double(aImage->
header,
"CRVAL2"),
2319 cd22 = cpl_propertylist_get_double(aImage->
header,
"CD2_2");
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);
2325 #ifdef ESO_ENABLE_DEBUG
2327 if (getenv(
"MUSE_DEBUG_NEAREST")) {
2328 debug = atoi(getenv(
"MUSE_DEBUG_NEAREST"));
2333 for (i = 0; i < aPixgrid->size_x; i++) {
2335 for (j = 0; j < aPixgrid->size_z; j++) {
2341 imagedata[i + j * aPixgrid->size_x] = data[rows[0]];
2342 #ifdef ESO_ENABLE_DEBUG
2344 printf(
"only: %f\n", data[rows[0]]);
2348 }
else if (n_rows >= 2) {
2350 cpl_size n, nbest = -1;
2351 double dbest = FLT_MAX;
2352 for (n = 0; n < n_rows; n++) {
2355 double dlambda = fabs(j * cd22 + crval2 - lbda[rows[n]]);
2356 #ifdef ESO_ENABLE_DEBUG
2358 printf(
"%f, %f, d: %f best: %f (%f)\n",
2359 j * cd22 + crval2, lbda[rows[n]],
2360 dlambda, dbest, data[rows[n]]);
2363 if (dlambda < dbest) {
2368 imagedata[i + j * aPixgrid->size_x] = data[rows[nbest]];
2369 #ifdef ESO_ENABLE_DEBUG
2371 printf(
"closest: %f\n", data[rows[nbest]]);
2381 return CPL_ERROR_NONE;
2401 static cpl_error_code
2405 cpl_ensure_code(aImage && aPixtable && aPixgrid, CPL_ERROR_NULL_INPUT);
2406 aImage->
data = cpl_image_new(aPixgrid->size_x, aPixgrid->size_z,
2409 double crval2 = cpl_propertylist_get_double(aImage->
header,
"CRVAL2"),
2410 cd22 = cpl_propertylist_get_double(aImage->
header,
"CD2_2");
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);
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"));
2422 if (getenv(
"MUSE_DEBUG_WEIGHTED_X")) {
2423 debugx = atoi(getenv(
"MUSE_DEBUG_WEIGHTED_X"));
2424 if (debugx < 1 || debugx > aPixgrid->size_x) {
2428 if (getenv(
"MUSE_DEBUG_WEIGHTED_Z")) {
2429 debugz = atoi(getenv(
"MUSE_DEBUG_WEIGHTED_Z"));
2430 if (debugz < 1 || debugz > aPixgrid->size_z) {
2440 double renka_rc = 1.25;
2445 for (i = 0; i < aPixgrid->size_x; i++) {
2447 for (j = 0; j < aPixgrid->size_z; j++) {
2448 double sumdata = 0, sumweight = 0,
2449 lambda = j * cd22 + crval2;
2451 #ifdef ESO_ENABLE_DEBUG
2452 #define MAX_NPIX_2D 50
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));
2463 for (j2 = j - ld; j2 <= j + ld; j2++) {
2467 for (n = 0; n < n_rows; n++) {
2468 double dlambda = fabs(lambda - lbda[rows[n]]),
2469 weight = muse_resampling_weight_function_renka(dlambda,
2471 sumweight += weight;
2472 sumdata += data[rows[n]] * weight;
2474 #ifdef ESO_ENABLE_DEBUG
2475 if (debug & 128 && i+1 == debugx && j+1 == debugz &&
2476 npoints-1 < MAX_NPIX_2D) {
2478 pointlist[npoints-1] = rows[n] + 1;
2479 pointweights[npoints-1] = weight;
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",
2487 weight, sumweight, sumdata, npoints);
2493 #ifdef ESO_ENABLE_DEBUG
2494 if (debug & 128 && i+1 == debugx && j+1 == debugz) {
2495 printf(
"pixelnumber weight ");
2498 while (++m < MAX_NPIX_2D && pointlist[m] != 0) {
2500 printf(
"%12d %8.5f ", pointlist[m] - 1, pointweights[m]);
2504 cpl_free(pointlist);
2505 cpl_free(pointweights);
2517 imagedata[i + j * aPixgrid->size_x] = sumdata / sumweight;
2521 return CPL_ERROR_NONE;
2545 double aDX,
double aLambdaMin,
2546 double aLambdaMax,
double aDLambda)
2548 double dlambda = aDLambda;
2551 aLambdaMin, aLambdaMax,
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);
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);
2572 cpl_propertylist_append_string(image->
header,
"CUNIT1",
"pixel");
2573 cpl_propertylist_append_string(image->
header,
"CUNIT2",
"Angstrom");
2576 cpl_propertylist_append_string(image->
header,
"CTYPE1",
"PIXEL");
2579 cpl_propertylist_append_string(image->
header,
"CTYPE2",
"AWAV");
2581 cpl_propertylist_append_double(image->
header,
"CD1_2", 0.);
2582 cpl_propertylist_append_double(image->
header,
"CD2_1", 0.);
2585 cpl_error_code rc = CPL_ERROR_NONE;
2588 rc = muse_resampling_image_nearest(image, aPixtable, pixgrid);
2592 rc = muse_resampling_image_weighted(image, aPixtable, pixgrid);
2596 if (rc != CPL_ERROR_NONE) {
2597 cpl_msg_error(__func__,
"resampling failed: %s", cpl_error_get_message());
2648 double aDX,
double aDLambda)
2650 cpl_ensure(aPixtable, CPL_ERROR_NULL_INPUT, NULL);
2651 if (aDLambda == 0.0) {
2652 aDLambda = kMuseSpectralSamplingA;
2662 cpl_msg_info(__func__,
"Using nearest neighbor sampling (%d) along "
2663 "wavelengths.", aMethod);
2666 cpl_msg_info(__func__,
"Using renka-weighted interpolation (%d) along "
2667 "wavelengths.", aMethod);
2670 cpl_msg_error(__func__,
"Don't know this resampling method: %d", aMethod);
2671 cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
2684 return muse_resampling_image_selected(aPixtable, aMethod,
2685 aDX == 0.0 ? 1. : aDX,
2686 lmin, lmax, aDLambda);
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),
2714 #pragma omp parallel for default(none) \
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++) {
2719 cpl_msg_debug(__func__,
"slice pixel table %d", i_slice + 1);
2723 slice_image[i_slice] = NULL;
2726 slice_image[i_slice] = muse_resampling_image_selected(pt, aMethod, dx,
2727 lmin, lmax, aDLambda);
2732 for (i_slice = 0; i_slice < n_slices; i_slice++) {
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));
2738 if (slice == NULL) {
2743 image->
header = cpl_propertylist_duplicate(slice->
header);
2746 cpl_image_delete(image->
data);
2750 cpl_image_delete(image->
dq);
2755 cpl_image_delete(image->
stat);
2759 slice_image[i_slice] = NULL;
2783 cpl_ensure(aPixtable != NULL, CPL_ERROR_NULL_INPUT, NULL);
2788 cpl_size lsize = floor((lmax - lmin) / aBinwidth) + 2;
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");
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));
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");
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);
2817 if (cpl_table_has_column(aPixtable->
table, MUSE_PIXTABLE_WEIGHT)) {
2818 wght = cpl_table_get_data_float(aPixtable->
table, MUSE_PIXTABLE_WEIGHT);
2820 int *dq = cpl_table_get_data_int(aPixtable->
table, MUSE_PIXTABLE_DQ);
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
2829 for (isel = 0; isel < nsel; isel++) {
2830 cpl_size i_row = sel[isel];
2831 if ((dq[i_row] != EURO3D_GOODPIXEL)) {
2835 double l = (lbda[i_row] - lmin) / aBinwidth;
2837 cpl_size il = floor(l);
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;
2853 cpl_array_delete(asel);
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);
2863 cpl_table_erase_selected(spectrum);
2866 cpl_table_divide_columns(spectrum,
"data",
"weight");
2867 cpl_table_divide_columns(spectrum,
"stat",
"weight");
2868 cpl_table_erase_column(spectrum,
"weight");
muse_pixtable_wcs
State of the astrometric calibration of a MUSE pixel table.
Structure definition of a MUSE datacube.
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.
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.
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
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...
#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 %)
A structure containing a spatial two-axis WCS.
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...
muse_pixtable ** muse_pixtable_extracted_get_slices(muse_pixtable *aPixtable)
Extract one pixel table per IFU and slice.
cpl_image * data
the data extension
#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.
void muse_utils_memory_dump(const char *aMarker)
Display the current memory usage of the given program.
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)
cpl_image * stat
the statistics extension
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.
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
double muse_astro_parangle(const cpl_propertylist *aHeader)
Properly average parallactic angle values of one exposure.
Structure definition of MUSE three extension FITS file.
cpl_table * table
The pixel table.
cpl_propertylist * header
the FITS header
muse_resampling_crstats_type crtype
cpl_image * dq
the data quality extension
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.
muse_wcs * muse_wcs_new(cpl_propertylist *aHeader)
Create a new WCS structure from a given FITS header.
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.
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.
cpl_imagelist * data
the cube containing the actual data values
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.
Structure definition of a Euro3D datacube.
static const cpl_size * muse_pixgrid_get_rows(muse_pixgrid *aPixels, cpl_size aIndex)
Return a pointer to the rows stored in one pixel.
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
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.
cpl_table * gtable
the table containing the Euro3D groups
cpl_propertylist * header
the FITS header
double muse_astro_posangle(const cpl_propertylist *aHeader)
Derive the position angle of an observation from information in a FITS header.
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.
#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)
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...
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.
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.
#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.
muse_ins_mode muse_pfits_get_mode(const cpl_propertylist *aHeaders)
find out the observation mode
double muse_astro_airmass(cpl_propertylist *aHeader)
Derive the effective airmass of an observation from information in a FITS header. ...
cpl_imagelist * stat
the cube containing the data variance
cpl_propertylist * header
The FITS header.
#define MUSE_HDR_PT_XHI
FITS header keyword contains the upper limit of the data in x-ion.