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_quadrants.h" 50 #include "muse_quality.h" 51 #include "muse_utils.h" 53 #include "muse_data_format_z.h" 121 cpl_ensure_code(aParams && aString, CPL_ERROR_NULL_INPUT);
123 int npf = cpl_array_get_size(pf);
124 cpl_error_code rc = CPL_ERROR_NONE;
127 aParams->
pfx = aParams->pfy = aParams->pfl = atof(cpl_array_get_string(pf, 0));
130 aParams->
pfx = aParams->pfy = atof(cpl_array_get_string(pf, 0));
131 aParams->pfl = atof(cpl_array_get_string(pf, 1));
134 aParams->
pfx = atof(cpl_array_get_string(pf, 0));
135 aParams->pfy = atof(cpl_array_get_string(pf, 1));
136 aParams->pfl = atof(cpl_array_get_string(pf, 2));
139 cpl_msg_warning(__func__,
"%d instead of 1-3 values (\"%s\") were given as " 140 "pixfrac, values remain unchanged (%.2f, %.2f, %.2f)!",
141 npf, aString, aParams->
pfx, aParams->pfy, aParams->pfl);
142 rc = CPL_ERROR_ILLEGAL_INPUT;
144 cpl_array_delete(pf);
170 const cpl_propertylist *aWCS)
172 cpl_ensure_code(aParams, CPL_ERROR_NULL_INPUT);
175 cpl_wcs_delete(aParams->
wcs);
177 return CPL_ERROR_NONE;
181 if (cpl_propertylist_has(aWCS,
"CTYPE3")) {
192 cpl_errorstate state = cpl_errorstate_get();
193 cpl_error_code rc = CPL_ERROR_NONE;
194 aParams->
wcs = cpl_wcs_new_from_propertylist(aWCS);
195 if (!cpl_errorstate_is_equal(state)) {
196 cpl_wcs_delete(aParams->
wcs);
198 rc = cpl_error_get_code();
201 printf(
"CDi_j matrix:\n");
202 cpl_matrix_dump(cpl_wcs_get_cd(aParams->
wcs), stdout);
223 cpl_wcs_delete(aParams->
wcs);
240 muse_resampling_weight_function_linear(
double r)
242 return r == 0 ? FLT_MAX : 1. / r;
258 muse_resampling_weight_function_quadratic(
double r2)
260 return r2 == 0 ? FLT_MAX : 1. / r2;
278 muse_resampling_weight_function_renka(
double r,
double r_c)
282 }
else if (r >= r_c) {
285 double p = (r_c - r) / (r_c * r);
299 muse_resampling_weight_function_sinc(
double r)
301 return fabs(r) < DBL_EPSILON ? 1. : sin(CPL_MATH_PI * r) / (CPL_MATH_PI * r);
316 muse_resampling_weight_function_lanczos(
double dx,
double dy,
double dz,
unsigned int n)
318 return (fabs(dx) >= n || fabs(dy) >= n || fabs(dz) > n) ? 0.
319 : muse_resampling_weight_function_sinc(dx) * muse_resampling_weight_function_sinc(dx / n)
320 * muse_resampling_weight_function_sinc(dy) * muse_resampling_weight_function_sinc(dy / n)
321 * muse_resampling_weight_function_sinc(dz) * muse_resampling_weight_function_sinc(dz / n);
345 muse_resampling_weight_function_drizzle(
double xin,
double yin,
double zin,
346 double xout,
double yout,
double zout,
347 double dx,
double dy,
double dz)
352 double x = (dx + xout / 2.) <= xin / 2. ? xout : (xin + xout) / 2. - dx,
353 y = (dy + yout / 2.) <= yin / 2. ? yout : (yin + yout) / 2. - dy,
354 z = (dz + zout / 2.) <= zin / 2. ? zout : (zin + zout) / 2. - dz;
357 if (x <= 0 || y <= 0 || z <= 0) {
362 return (x > xin ? xin : x) * (y > yin ? yin : y) * (z > zin ? zin : z)
383 static cpl_error_code
387 cpl_ensure_code(aCube && aPixtable && aGrid, CPL_ERROR_NULL_INPUT);
396 cpl_boolean loglambda = ctype3 && (!strncmp(ctype3,
"AWAV-LOG", 9) ||
397 !strncmp(ctype3,
"WAVE-LOG", 9));
409 double xnorm = 1., ynorm = 1., znorm = 1. / kMuseSpectralSamplingA;
419 #ifdef ESO_ENABLE_DEBUG 421 if (getenv(
"MUSE_DEBUG_NEAREST")) {
422 debug = atoi(getenv(
"MUSE_DEBUG_NEAREST"));
427 #ifdef ESO_ENABLE_DEBUG 428 #pragma omp parallel for default(none) \ 429 shared(aCube, aGrid, cd33, crpix3, crval3, data, debug, dq, lbda, \ 430 loglambda, ptxoff, ptyoff, stat, stdout, wcs, xnorm, xpos, \ 433 #pragma omp parallel for default(none) \ 434 shared(aCube, aGrid, cd33, crpix3, crval3, data, dq, lbda, \ 435 loglambda, ptxoff, ptyoff, stat, stdout, wcs, xnorm, xpos, \ 438 for (l = 0; l < aGrid->nz; l++) {
439 float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->
data, l)),
440 *pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->
stat, l));
441 int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->
dq, l));
443 double lambda = (l + 1. - crpix3) * cd33 + crval3;
445 lambda = crval3 * exp((l + 1. - crpix3) * cd33 / crval3);
449 for (i = 0; i < aGrid->nx; i++) {
451 for (j = 0; j < aGrid->ny; j++) {
466 pdata[i + j * aGrid->nx] = data[rows[0]];
467 pstat[i + j * aGrid->nx] = stat[rows[0]];
468 pdq[i + j * aGrid->nx] = dq[rows[0]];
470 #ifdef ESO_ENABLE_DEBUG 472 printf(
"only: %f,%f\n", data[rows[0]], stat[rows[0]]);
476 }
else if (n_rows >= 2) {
478 cpl_size n, nbest = -1;
479 double dbest = FLT_MAX;
480 for (n = 0; n < n_rows; n++) {
482 double dx = fabs(x - xpos[rows[n]] + ptxoff) * xnorm,
483 dy = fabs(y - ypos[rows[n]] + ptyoff) * ynorm,
484 dlambda = fabs(lambda - lbda[rows[n]]) * znorm,
485 dthis = sqrt(dx*dx + dy*dy + dlambda*dlambda);
489 dx *= cos(y * CPL_MATH_RAD_DEG);
491 #ifdef ESO_ENABLE_DEBUG 493 printf(
"%f %f %f, %f %f %f, d: %f %f %f -> %f best: %f (%f,%f)\n",
494 x, y, lambda, xpos[rows[n]] + ptxoff, ypos[rows[n]] + ptyoff,
495 lbda[rows[n]], dx, dy, dlambda, dthis, dbest, data[rows[n]],
504 pdata[i + j * aGrid->nx] = data[rows[nbest]];
505 pstat[i + j * aGrid->nx] = stat[rows[nbest]];
506 pdq[i + j * aGrid->nx] = dq[rows[nbest]];
507 #ifdef ESO_ENABLE_DEBUG 509 printf(
"closest: %f,%f\n", data[rows[nbest]], stat[rows[nbest]]);
515 pdq[i + j * aGrid->nx] = EURO3D_MISSDATA;
522 return CPL_ERROR_NONE;
547 static cpl_error_code
552 cpl_ensure_code(aCube && aPixtable && aGrid && aParams,
553 CPL_ERROR_NULL_INPUT);
562 cpl_boolean loglambda = ctype3 && (!strncmp(ctype3,
"AWAV-LOG", 9) ||
563 !strncmp(ctype3,
"WAVE-LOG", 9));
564 cpl_errorstate prestate = cpl_errorstate_get();
571 if (!cpl_errorstate_is_equal(prestate)) {
572 cpl_errorstate_set(prestate);
581 double xnorm = 1., ynorm = 1., znorm = 1. / kMuseSpectralSamplingA;
592 double zoutefac = exp(1.5 * cd33 / crval3) - exp(0.5 * cd33 / crval3);
594 double renka_rc = aParams->
rc 595 * sqrt((wcs->cd11*xnorm)*(wcs->cd11*xnorm) + (wcs->
cd22*ynorm)*(wcs->
cd22*ynorm)
596 + (cd33*znorm)*(cd33*znorm));
598 int ld = aParams->
ld;
601 cpl_msg_info(__func__,
"Overriding loop distance ld=%d", ld);
606 double xsz = aParams->
pfx / xnorm,
607 ysz = aParams->pfy / ynorm,
608 zsz = aParams->pfl / znorm,
609 xout = fabs(wcs->cd11), yout = fabs(wcs->
cd22), zout = fabs(cd33);
612 char *envmaskcube = getenv(
"MUSE_BADPIX_FROM_MASKCUBE"),
615 cpl_table *bpt = NULL;
616 cpl_imagelist *maskcube = NULL;
620 int nfiles = cpl_array_get_size(files);
622 maskname = cpl_strdup(cpl_array_get_string(files, 0));
623 badpix = cpl_strdup(cpl_array_get_string(files, 1));
624 }
else if (nfiles == 1) {
625 maskname = cpl_strdup(cpl_array_get_string(files, 0));
626 badpix = cpl_sprintf(
"BADPIX_FROM_MASKCUBE.fits");
627 cpl_msg_warning(__func__,
"MUSE_BADPIX_FROM_MASKCUBE with only one " 628 "argument, using \"%s\" forthe bad pixel table", badpix);
630 cpl_msg_warning(__func__,
"MUSE_BADPIX_FROM_MASKCUBE without arguments, " 631 "not creating a mask table!");
634 maskcube = cpl_imagelist_load(maskname, CPL_TYPE_INT, 0);
638 int xsize = cpl_image_get_size_x(cpl_imagelist_get(maskcube, 0)),
639 ysize = cpl_image_get_size_y(cpl_imagelist_get(maskcube, 0)),
640 zsize = cpl_imagelist_get_size(maskcube);
641 if (zsize != aGrid->nz) {
642 cpl_msg_warning(__func__,
"Invalid z size of mask cube (%d vs %d)",
643 zsize, (
int)aGrid->nz);
644 cpl_imagelist_delete(maskcube);
646 }
else if (xsize != aGrid->nx) {
647 cpl_msg_warning(__func__,
"Invalid x size of mask cube (%d vs %d)",
648 xsize, (
int)aGrid->nx);
649 cpl_imagelist_delete(maskcube);
651 }
else if (ysize != aGrid->ny) {
652 cpl_msg_warning(__func__,
"Invalid y size of mask cube (%d vs %d)",
653 ysize, (
int)aGrid->ny);
654 cpl_imagelist_delete(maskcube);
657 cpl_msg_info(__func__,
"Mask cube with valid size found (%d x %d x %d)",
658 xsize, ysize, zsize);
667 cpl_table_new_column(bpt,
"pweight", CPL_TYPE_DOUBLE);
668 cpl_table_new_column(bpt,
"ifu", CPL_TYPE_INT);
670 cpl_msg_warning(__func__,
"No mask cube found, not creating output bad " 673 cpl_array_delete(files);
676 #ifdef ESO_ENABLE_DEBUG 677 int debug = 0, debugx = 0, debugy = 0, debugz = 0;
678 if (getenv(
"MUSE_DEBUG_WEIGHTED")) {
679 debug = atoi(getenv(
"MUSE_DEBUG_WEIGHTED"));
682 if (getenv(
"MUSE_DEBUG_WEIGHTED_X")) {
683 debugx = atoi(getenv(
"MUSE_DEBUG_WEIGHTED_X"));
684 if (debugx < 1 || debugx > aGrid->nx) {
688 if (getenv(
"MUSE_DEBUG_WEIGHTED_Y")) {
689 debugy = atoi(getenv(
"MUSE_DEBUG_WEIGHTED_Y"));
690 if (debugy < 1 || debugy > aGrid->ny) {
694 if (getenv(
"MUSE_DEBUG_WEIGHTED_Z")) {
695 debugz = atoi(getenv(
"MUSE_DEBUG_WEIGHTED_Z"));
696 if (debugz < 1 || debugz > aGrid->nz) {
702 printf(
"parameters:\n cd=%e %e %e\n" 703 " corrected crpix3=%e\n norm=%e %e %e\n",
704 wcs->cd11, wcs->
cd22, cd33, crpix3, xnorm, ynorm, znorm);
706 printf(
" drop sizes %e %e %e (pixfrac %f,%f,%f)\n" 707 " output sizes %e %e %e\n",
708 xsz, ysz, zsz, aParams->
pfx, aParams->pfy, aParams->pfl,
711 printf(
" resulting renka_rc: %e %e %e / %e %e %e --> %e\n",
712 pow(wcs->cd11, 2), pow(wcs->
cd22, 2), cd33*cd33,
713 pow(wcs->cd11*xnorm, 2), pow(wcs->
cd22*ynorm, 2),
714 pow(cd33*znorm, 2), renka_rc);
719 cpl_imagelist *wcube = NULL;
720 if (getenv(
"MUSE_DEBUG_WEIGHT_CUBE")) {
721 cpl_msg_debug(__func__,
"Weighted resampling: creating weight cube");
722 wcube = cpl_imagelist_new();
724 for (i = 0; i < aGrid->nz; i++) {
725 cpl_image *image = cpl_image_new(aGrid->nx, aGrid->ny,
727 cpl_imagelist_set(wcube, image, i);
731 if (getenv(
"MUSE_DEBUG_WEIGHTED_GRID")) {
732 char *fn = getenv(
"MUSE_DEBUG_WEIGHTED_GRID");
733 FILE *grid = fopen(fn,
"w");
735 cpl_msg_info(__func__,
"Writing grid to \"%s\"", fn);
736 fprintf(grid,
"# i j l x y lambda\n");
738 for (l = 0; l < aGrid->nz; l++) {
739 double lambda = (l + 1. - crpix3) * cd33 + crval3;
741 for (i = 0; i < aGrid->nx; i++) {
743 for (j = 0; j < aGrid->ny; j++) {
751 fprintf(grid,
"%03d %03d %04d %.10f %.10f %8.3f\n", i+1, j+1, l+1,
758 cpl_msg_warning(__func__,
"Writing grid to \"%s\" failed!", fn);
763 #ifdef ESO_ENABLE_DEBUG 764 #pragma omp parallel for default(none) \ 765 shared(aCube, aParams, aGrid, aPixtable, bpt, cd33, crpix3, \ 766 crval3, data, debug, debugx, debugy, debugz, dq, lbda, ld, \ 767 loglambda, maskcube, orgn, ptxoff, ptyoff, renka_rc, stat, \ 768 stdout, wcs, wcube, wght, xnorm, xout, xpos, xsz, ynorm, \ 769 yout, ypos, ysz, znorm, zout, zoutefac, zsz) 771 #pragma omp parallel for default(none) \ 772 shared(aCube, aParams, aGrid, aPixtable, bpt, cd33, crpix3, \ 773 crval3, data, dq, lbda, ld, loglambda, maskcube, orgn, ptxoff,\ 774 ptyoff, renka_rc, stat, stdout, wcs, wcube, wght, xnorm, xout,\ 775 xpos, xsz, ynorm, yout, ypos, ysz, znorm, zout, zoutefac, zsz) 777 for (l = 0; l < aGrid->nz; l++) {
778 float *pdata = cpl_image_get_data_float(cpl_imagelist_get(aCube->
data, l)),
779 *pstat = cpl_image_get_data_float(cpl_imagelist_get(aCube->
stat, l));
780 int *pdq = cpl_image_get_data_int(cpl_imagelist_get(aCube->
dq, l));
782 double lambda = (l + 1. - crpix3) * cd33 + crval3;
785 lambda = crval3 * exp((l + 1. - crpix3) * cd33 / crval3);
786 zout2 = crval3 * exp((l - crpix3) * cd33 / crval3) * zoutefac;
788 float *pwcube = NULL;
790 pwcube = cpl_image_get_data_float(cpl_imagelist_get(wcube, l));
793 cpl_image *mask = NULL;
795 mask = cpl_imagelist_get(maskcube, l);
799 for (i = 0; i < aGrid->nx; i++) {
801 for (j = 0; j < aGrid->ny; j++) {
809 double sumdata = 0, sumstat = 0, sumweight = 0;
811 #ifdef ESO_ENABLE_DEBUG 812 cpl_size *pointlist = NULL;
813 double *pointweights = NULL;
815 if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
816 pointlist = cpl_calloc(100,
sizeof(cpl_size));
817 pointweights = cpl_malloc(100 *
sizeof(
double));
822 cpl_size *plist = NULL;
823 double *pweights = NULL;
826 if (mask && cpl_image_get(mask, i+1, j+1, &err) != 0) {
827 plist = cpl_calloc(100,
sizeof(cpl_size));
828 pweights = cpl_malloc(100 *
sizeof(
double));
832 #ifdef ESO_ENABLE_DEBUG 834 printf(
"cell %d %d %d:\n", i, j, l);
839 for (i2 = i - ld; i2 <= i + ld; i2++) {
841 for (j2 = j - ld; j2 <= j + ld; j2++) {
843 for (l2 = l - ld; l2 <= l + ld; l2++) {
846 #ifdef ESO_ENABLE_DEBUG 847 if (debug & 8 && n_rows2 < 1) {
848 printf(
"%d %d %d / %d %d %d (%"CPL_SIZE_FORMAT
"): no rows!\n",
849 i+1, j+1, l+1, i2+1, j2+1, l2+1, idx2);
854 for (n = 0; n < n_rows2; n++) {
856 #ifdef ESO_ENABLE_DEBUG 858 printf(
"%d %d %d / %d %d %d (%"CPL_SIZE_FORMAT
", " 859 "%"CPL_SIZE_FORMAT
"): bad!\n",
860 i+1, j+1, l+1, i2+1, j2+1, l2+1, idx2, n);
867 double dx = fabs(x - (xpos[rows2[n]] + ptxoff)),
868 dy = fabs(y - (ypos[rows2[n]] + ptyoff)),
869 dlambda = fabs(lambda - lbda[rows2[n]]),
877 dx *= cos(y * CPL_MATH_RAD_DEG);
883 r2 = dx*dx + dy*dy + dlambda*dlambda;
887 weight = muse_resampling_weight_function_renka(sqrt(r2), renka_rc);
889 weight = muse_resampling_weight_function_drizzle(xsz, ysz, zsz,
893 weight = muse_resampling_weight_function_linear(sqrt(r2));
895 weight = muse_resampling_weight_function_quadratic(r2);
897 weight = muse_resampling_weight_function_lanczos(dx, dy, dlambda, ld);
902 weight *= wght[rows2[n]];
904 #ifdef ESO_ENABLE_DEBUG 906 printf(
"%d %d %d / %d %d %d (%"CPL_SIZE_FORMAT
", %"CPL_SIZE_FORMAT
"):" 907 " x %e %e %e %e y %e %e %e %e l %e %e %e %e --> %e %e %e\n",
908 i+1, j+1, l+1, i2+1, j2+1, l2+1, idx2, n,
909 x, xpos[rows2[n]]+ptxoff, fabs(x - (xpos[rows2[n]]+ptxoff)), dx,
910 y, ypos[rows2[n]]+ptyoff, fabs(y - (ypos[rows2[n]]+ptyoff)), dy,
911 lambda, lbda[rows2[n]], fabs(lambda - lbda[rows2[n]]), dlambda,
912 r2, sqrt(r2), weight);
917 sumdata += data[rows2[n]] * weight;
918 sumstat += stat[rows2[n]] * weight*weight;
920 #ifdef ESO_ENABLE_DEBUG 921 if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
922 if (npoints > npointlist) {
923 pointlist = cpl_realloc(pointlist,
924 (npointlist + 100) *
sizeof(cpl_size));
925 memset(pointlist + npointlist, 0, 100 *
sizeof(cpl_size));
926 pointweights = cpl_realloc(pointweights,
927 (npointlist + 100) *
sizeof(
double));
931 pointlist[npoints-1] = rows2[n] + 1;
932 pointweights[npoints-1] = weight;
939 printf(
" pixel %4d,%4d,%4d (%8"CPL_SIZE_FORMAT
"): " 940 "%2"CPL_SIZE_FORMAT
" %2"CPL_SIZE_FORMAT
" %f %f %f, " 941 " %e -> %e ==> %e %e (%d)\n", i+1, j+1, l+1, idx,
942 n, count, dx, dy, dlambda,
944 weight, sumweight, sumdata, npoints);
949 if (plist && pweights) {
950 if (npoints > nplist) {
951 plist = cpl_realloc(plist, (nplist + 100) *
sizeof(cpl_size));
952 memset(plist + nplist, 0, 100 *
sizeof(cpl_size));
953 pweights = cpl_realloc(pweights, (nplist + 100) *
sizeof(
double));
957 plist[npoints-1] = rows2[n] + 1;
958 pweights[npoints-1] = weight;
965 #ifdef ESO_ENABLE_DEBUG 966 if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
967 printf(
"cell center (%d, %d, %d): %14.7e %14.7e %9.3f\npixelnumber " 968 "weight ", debugx, debugy, debugz, x, y, lambda);
971 while (++m < npointlist && pointlist[m] != 0) {
973 printf(
"%12"CPL_SIZE_FORMAT
" %8.5f ", pointlist[m] - 1,
977 printf(
"sums: %g %g %g --> %g %g\n", sumdata, sumstat, sumweight,
978 sumdata / sumweight, sumstat / pow(sumweight, 2));
980 cpl_free(pointweights);
983 if (debug & 1 && npoints) {
984 printf(
" sumdata: %e %e (%d)", sumdata, sumweight, npoints);
988 if (plist && pweights) {
990 #pragma omp critical(muse_bpt_size) 993 ibad = cpl_table_get_nrow(bpt);
995 cpl_table_set_size(bpt, ibad + nplist);
998 while (++m < nplist && plist[m] != 0) {
1000 cpl_size ipt = plist[m] - 1;
1005 #pragma omp critical(muse_bpt_size) 1008 cpl_table_set_int(bpt, MUSE_BADPIX_X, ibad, xraw);
1009 cpl_table_set_int(bpt, MUSE_BADPIX_Y, ibad, yraw);
1010 cpl_table_set_int(bpt, MUSE_BADPIX_DQ, ibad, EURO3D_BADOTHER);
1011 cpl_table_set_double(bpt,
"pweight", ibad, pweights[m]);
1012 cpl_table_set_int(bpt,
"ifu", ibad, ifu);
1015 #ifdef ESO_ENABLE_DEBUG 1016 if (debug & (i+1 == debugx && j+1 == debugy && l+1 == debugz)) {
1017 printf(
"%12"CPL_SIZE_FORMAT
" %8.5f %02hhu %4d %4d\n",
1018 ipt, pweights[m], ifu, xraw, yraw);
1032 if (!npoints || !isnormal(sumweight)) {
1033 pdq[i + j * aGrid->nx] = EURO3D_MISSDATA;
1034 #ifdef ESO_ENABLE_DEBUG 1036 printf(
" -> no points or weird weight\n");
1043 sumdata /= sumweight;
1044 sumstat /= sumweight*sumweight;
1045 #ifdef ESO_ENABLE_DEBUG 1047 printf(
" -> %e (%e)\n", sumdata, sumstat);
1053 pdata[i + j * aGrid->nx] = sumdata;
1054 pstat[i + j * aGrid->nx] = sumstat;
1055 pdq[i + j * aGrid->nx] = EURO3D_GOODPIXEL;
1057 pwcube[i + j * aGrid->nx] = sumweight;
1065 const char *fn = getenv(
"MUSE_DEBUG_WEIGHT_CUBE");
1066 cpl_error_code rc = cpl_imagelist_save(wcube, fn, CPL_TYPE_UNSPECIFIED,
1067 NULL, CPL_IO_CREATE);
1068 if (rc != CPL_ERROR_NONE) {
1069 cpl_msg_warning(__func__,
"Failure to save weight cube as \"%s\": %s", fn,
1070 cpl_error_get_message());
1072 cpl_msg_info(__func__,
"Saved weight cube as \"%s\"", fn);
1074 cpl_imagelist_delete(wcube);
1081 cpl_boolean isfirst = CPL_TRUE;
1082 for (nifu = 1; nifu <= kMuseNumIFUs; nifu++) {
1084 cpl_table_unselect_all(bpt);
1085 cpl_table_or_selected_int(bpt,
"ifu", CPL_EQUAL_TO, nifu);
1086 cpl_table_and_selected_double(bpt,
"pweight", CPL_GREATER_THAN, 0.);
1089 cpl_table *bptifu = cpl_table_extract_selected(bpt);
1090 cpl_size nbad = cpl_table_get_nrow(bptifu);
1092 cpl_msg_info(__func__,
"No bad pixels from mask in IFU %02hhu, " 1093 "nothing saved.", nifu);
1094 cpl_table_delete(bptifu);
1099 cpl_propertylist *sorting = cpl_propertylist_new();
1100 cpl_propertylist_append_bool(sorting, MUSE_BADPIX_X, CPL_FALSE);
1101 cpl_propertylist_append_bool(sorting, MUSE_BADPIX_Y, CPL_FALSE);
1102 cpl_table_sort(bptifu, sorting);
1103 cpl_propertylist_delete(sorting);
1106 cpl_table_unselect_all(bptifu);
1108 xprev = cpl_table_get_int(bptifu, MUSE_BADPIX_X, 0, NULL),
1109 yprev = cpl_table_get_int(bptifu, MUSE_BADPIX_Y, 0, NULL);
1110 for (irow = 1; irow < nbad; irow++) {
1111 int x = cpl_table_get_int(bptifu, MUSE_BADPIX_X, irow, NULL),
1112 y = cpl_table_get_int(bptifu, MUSE_BADPIX_Y, irow, NULL);
1113 if (x == xprev && y == yprev) {
1114 cpl_table_select_row(bptifu, irow);
1121 cpl_table_erase_selected(bptifu);
1124 cpl_table_erase_column(bptifu,
"pweight");
1125 cpl_table_erase_column(bptifu,
"ifu");
1129 cpl_propertylist *header = cpl_propertylist_new(),
1131 char *chan = cpl_sprintf(
"CHAN%02hhu", nifu);
1132 cpl_propertylist_append_string(header,
"EXTNAME", chan);
1134 pheader = cpl_propertylist_new();
1135 cpl_propertylist_append_string(pheader,
"INSTRUME",
"MUSE");
1136 cpl_propertylist_append_double(pheader,
"MJD-OBS",
1138 cpl_propertylist_append_string(pheader,
"DATE-OBS",
1140 cpl_propertylist_append_string(pheader,
"ESO DRS MUSE MASKCUBE",
1143 cpl_table_save(bptifu, pheader, header, badpix,
1144 isfirst ? CPL_IO_CREATE : CPL_IO_EXTEND);
1146 cpl_propertylist_delete(header);
1147 cpl_propertylist_delete(pheader);
1148 cpl_table_delete(bptifu);
1149 cpl_msg_info(__func__,
"Bad pixel table from mask cube for IFU %02hhu " 1150 "saved to \"%s\" (%"CPL_SIZE_FORMAT
" pixels).", nifu,
1152 isfirst = CPL_FALSE;
1154 cpl_table_delete(bpt);
1158 cpl_frameset *fset = cpl_frameset_new();
1159 cpl_frame *frame = cpl_frame_new();
1160 cpl_frame_set_filename(frame, badpix);
1161 cpl_frame_set_tag(frame, MUSE_TAG_BADPIX_TABLE);
1162 cpl_frame_set_type(frame, CPL_FRAME_TYPE_TABLE);
1163 cpl_frame_set_group(frame, CPL_FRAME_GROUP_PRODUCT);
1164 cpl_frame_set_level(frame, CPL_FRAME_LEVEL_FINAL);
1165 cpl_frameset_insert(fset, frame);
1166 cpl_dfs_sign_products(fset, CPL_DFS_SIGNATURE_DATAMD5 | CPL_DFS_SIGNATURE_CHECKSUM);
1167 cpl_frameset_delete(fset);
1171 cpl_imagelist_delete(maskcube);
1173 return CPL_ERROR_NONE;
1187 static cpl_error_code
1191 cpl_ensure_code(aPixtable && aParams, CPL_ERROR_NULL_INPUT);
1192 const char func[] =
"muse_resampling_cube";
1195 if (aParams->dlambda == 0.0) {
1196 aParams->dlambda = kMuseSpectralSamplingA;
1199 aParams->dlambda /= 1.6;
1205 if (aParams->
dx == 0.0) {
1208 if (aParams->dy == 0.0) {
1211 cpl_msg_debug(func,
"steps from parameters: %.2f pix, %.2f pix, %.3f Angstrom",
1212 aParams->
dx, aParams->dy, aParams->dlambda);
1213 return CPL_ERROR_NONE;
1215 if (aParams->
dx == 0.0) {
1216 aParams->
dx = kMuseSpaxelSizeX_WFM / 3600.;
1218 aParams->
dx = kMuseSpaxelSizeX_NFM / 3600.;
1222 aParams->
dx /= 3600.;
1224 if (aParams->dy == 0.0) {
1225 aParams->dy = kMuseSpaxelSizeY_WFM / 3600.;
1227 aParams->dy = kMuseSpaxelSizeY_NFM / 3600.;
1231 aParams->dy /= 3600.;
1233 cpl_msg_debug(func,
"steps from parameters: %f arcsec, %f arcsec, %.3f Angstrom",
1234 aParams->
dx * 3600., aParams->dy * 3600., aParams->dlambda);
1235 return CPL_ERROR_NONE;
1253 static cpl_error_code
1256 int *aX,
int *aY,
int *aZ)
1258 cpl_ensure_code(aPixtable && aParams && aX && aY && aZ, CPL_ERROR_NULL_INPUT);
1259 const char func[] =
"muse_resampling_cube";
1260 double x1, y1, x2, y2;
1273 *aX = lround(fabs(x2 - x1) / aParams->
dx) + 1;
1274 *aY = lround(fabs(y2 - y1) / aParams->dy) + 1;
1277 *aZ = (int)ceil((lmax - lmin) / aParams->dlambda) + 1;
1280 *aZ = (int)ceil(lmin / aParams->dlambda * log(lmax / lmin)) + 1;
1282 cpl_msg_info(func,
"Output cube size %d x %d x %d (fit to data)",
1284 return CPL_ERROR_NONE;
1299 muse_resampling_override_size_int(
int *aV,
const char *aKeyword,
int aValue)
1301 if (aValue <= 0 || !aV) {
1304 const char func[] =
"muse_resampling_cube";
1305 cpl_msg_info(func,
"Overriding %s=%d (was %d)", aKeyword, aValue, *aV);
1324 static cpl_error_code
1325 muse_resampling_override_size(
int *aX,
int *aY,
int *aZ,
1328 cpl_ensure_code(aX && aY && aZ && aParams, CPL_ERROR_NULL_INPUT);
1329 if (!aParams->
wcs) {
1330 return CPL_ERROR_NONE;
1332 const char func[] =
"muse_resampling_cube";
1334 const cpl_array *dims = cpl_wcs_get_image_dims(aParams->
wcs);
1336 cpl_msg_debug(func,
"No dimensions to override were specified");
1337 return CPL_ERROR_NONE;
1339 muse_resampling_override_size_int(aX,
"NAXIS1", cpl_array_get_int(dims, 0, NULL));
1340 muse_resampling_override_size_int(aY,
"NAXIS2", cpl_array_get_int(dims, 1, NULL));
1341 if (cpl_wcs_get_image_naxis(aParams->
wcs) >= 3) {
1342 muse_resampling_override_size_int(aZ,
"NAXIS3",
1343 cpl_array_get_int(dims, 2, NULL));
1345 return CPL_ERROR_NONE;
1361 muse_resampling_override_wcs_double(
muse_datacube *aCube,
const char *aKeyword,
1362 double aValue,
int aError)
1364 if (aError || !aCube) {
1365 cpl_msg_debug(
"double",
"%s=%#g (%d)", aKeyword, aValue, aError);
1368 const char func[] =
"muse_resampling_cube";
1369 double old = cpl_propertylist_has(aCube->
header, aKeyword)
1370 ? cpl_propertylist_get_double(aCube->
header, aKeyword) : NAN;
1371 cpl_msg_info(func,
"Overriding %s=%#g (was %#g)", aKeyword, aValue, old);
1372 cpl_propertylist_update_double(aCube->
header, aKeyword, aValue);
1375 cpl_propertylist_update_bool(aCube->
header,
"MUSE_RESAMPLING_WCS_OVERRIDDEN",
1392 static cpl_error_code
1396 cpl_ensure_code(aCube && aCube->
header && aParams, CPL_ERROR_NULL_INPUT);
1397 if (!aParams->
wcs) {
1398 return CPL_ERROR_NONE;
1400 const char func[] =
"muse_resampling_cube";
1401 const cpl_array *crval = cpl_wcs_get_crval(aParams->
wcs),
1402 *crpix = cpl_wcs_get_crpix(aParams->
wcs);
1403 const cpl_matrix *cd = cpl_wcs_get_cd(aParams->
wcs);
1407 muse_resampling_override_wcs_double(aCube,
"CRVAL1", cpl_array_get_double(crval, 0, &err), err);
1408 muse_resampling_override_wcs_double(aCube,
"CRVAL2", cpl_array_get_double(crval, 1, &err), err);
1410 cpl_msg_debug(func,
"No CRVALj to override were specified");
1413 muse_resampling_override_wcs_double(aCube,
"CRPIX1", cpl_array_get_double(crpix, 0, &err), err);
1414 muse_resampling_override_wcs_double(aCube,
"CRPIX2", cpl_array_get_double(crpix, 1, &err), err);
1416 cpl_msg_debug(func,
"No CRPIXi to override were specified");
1419 int naxes = cpl_matrix_get_ncol(cd);
1420 if (cpl_matrix_get_determinant(cd) == 0.) {
1421 cpl_msg_warning(func,
"det(CDi_j) = 0, not overriding CDi_j!");
1423 }
else if (naxes > 2 && (cpl_matrix_get(cd, 0, 2) != 0. ||
1424 cpl_matrix_get(cd, 2, 0) != 0. ||
1425 cpl_matrix_get(cd, 2, 1) != 0. ||
1426 cpl_matrix_get(cd, 1, 2) != 0.)) {
1427 cpl_msg_warning(func,
"Axis 3 (dispersion) is not cleanly separated from " 1428 "axes 1 and 2, not overriding CDi_j!");
1431 cpl_errorstate es = cpl_errorstate_get();
1432 muse_resampling_override_wcs_double(aCube,
"CD1_1", cpl_matrix_get(cd, 0, 0),
1433 !cpl_errorstate_is_equal(es));
1434 es = cpl_errorstate_get();
1435 muse_resampling_override_wcs_double(aCube,
"CD1_2", cpl_matrix_get(cd, 0, 1),
1436 !cpl_errorstate_is_equal(es));
1437 es = cpl_errorstate_get();
1438 muse_resampling_override_wcs_double(aCube,
"CD2_1", cpl_matrix_get(cd, 1, 0),
1439 !cpl_errorstate_is_equal(es));
1440 es = cpl_errorstate_get();
1441 muse_resampling_override_wcs_double(aCube,
"CD2_2", cpl_matrix_get(cd, 1, 1),
1442 !cpl_errorstate_is_equal(es));
1445 cpl_msg_debug(func,
"No CDi_j to override were specified");
1450 if (cpl_array_get_size(crval) > 2) {
1452 muse_resampling_override_wcs_double(aCube,
"CRVAL3",
1453 cpl_array_get_double(crval, 2, &err) * 1e10, err);
1456 muse_resampling_override_wcs_double(aCube,
"CRPIX3", cpl_array_get_double(crpix, 2, &err), err);
1459 cpl_errorstate es = cpl_errorstate_get();
1460 muse_resampling_override_wcs_double(aCube,
"CD3_3", cpl_matrix_get(cd, 2, 2) * 1e10,
1461 !cpl_errorstate_is_equal(es));
1464 return CPL_ERROR_NONE;
1478 static cpl_error_code
1481 cpl_ensure_code(aCube && aPixtable, CPL_ERROR_NULL_INPUT);
1482 const char func[] =
"muse_resampling_cube";
1486 if (cpl_propertylist_has(aCube->
header,
"MUSE_RESAMPLING_WCS_OVERRIDDEN") &&
1487 cpl_propertylist_get_bool(aCube->
header,
"MUSE_RESAMPLING_WCS_OVERRIDDEN")) {
1488 cpl_msg_debug(func,
"Output WCS was forced, won't update CRPIX[12]!");
1489 cpl_propertylist_erase(aCube->
header,
"MUSE_RESAMPLING_WCS_OVERRIDDEN");
1490 return CPL_ERROR_NONE;
1493 double xoff1, yoff1, xoff2, yoff2;
1508 double xoff = fmin(xoff1, xoff2) - 1,
1509 yoff = fmin(yoff1, yoff2) - 1;
1510 if (xoff != 0. || yoff != 0.) {
1513 cpl_msg_debug(func,
"Updating CRPIX[12]=%#g,%#g (were %#g,%#g)",
1514 crpix1 - xoff, crpix2 - yoff, crpix1, crpix2);
1517 cpl_propertylist_update_double(aCube->
header,
"CRPIX1", crpix1);
1518 cpl_propertylist_update_double(aCube->
header,
"CRPIX2", crpix2);
1520 return CPL_ERROR_NONE;
1525 const char *crtypestring[] = {
1531 #define MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_START \ 1532 cpl_size idx = muse_pixgrid_get_index(aGrid, i2, j2, l2, CPL_FALSE), \ 1533 nrow = muse_pixgrid_get_count(aGrid, idx), \ 1535 const cpl_size *rows = muse_pixgrid_get_rows(aGrid, idx); \ 1536 for (irow = 0; irow < nrow; irow++) { \ 1537 if (muse_quality_is_usable(dq[rows[irow]], badinclude)) { \ 1539 } else if (dq[rows[irow]]) { \ 1542 #define MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_DEBUG \ 1543 if (debug & 4 && i+1 == debugx && j+1 == debugy && l+1 == debugz) { \ 1544 printf("%s: data / stat (%"CPL_SIZE_FORMAT") = %e / %e\n", \ 1545 __func__, npixels+1, data[rows[irow]], stat[rows[irow]]); \ 1547 #define MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END \ 1548 if (aType == MUSE_RESAMPLING_CRSTATS_IRAF) { \ 1550 level += data[rows[irow]]; \ 1551 dev += stat[rows[irow]]; \ 1554 if (npixels >= sxsize) { \ 1556 sxsize += MUSE_CRREJECT_MAX_NPIX; \ 1557 cpl_image *tmp = cpl_image_new(sxsize, 2, CPL_TYPE_DOUBLE); \ 1558 cpl_image_copy(tmp, sdata, 1, 1); \ 1559 cpl_image_delete(sdata); \ 1561 vdata = cpl_image_get_data_double(sdata); \ 1563 vdata[npixels] = data[rows[irow]]; \ 1564 vdata[npixels + sxsize] = stat[rows[irow]]; \ 1598 const char *
id =
"muse_resampling_cube";
1600 cpl_msg_warning(
id,
"Unknown type (%u) for cosmic-ray rejection statistics," 1601 " resetting to \"%s\"", aType,
1605 cpl_msg_info(
id,
"Using %s statistics (%.3f sigma level) for cosmic-ray" 1606 " rejection", crtypestring[aType], aSigma);
1608 #ifdef ESO_ENABLE_DEBUG 1609 int debug = 0, debugx = 0, debugy = 0, debugz = 0;
1610 if (getenv(
"MUSE_DEBUG_CRREJECT")) {
1611 debug = atoi(getenv(
"MUSE_DEBUG_CRREJECT"));
1614 if (getenv(
"MUSE_DEBUG_CRREJECT_X")) {
1615 debugx = atoi(getenv(
"MUSE_DEBUG_CRREJECT_X"));
1616 if (debugx < 1 || debugx > aGrid->nx) {
1620 if (getenv(
"MUSE_DEBUG_CRREJECT_Y")) {
1621 debugy = atoi(getenv(
"MUSE_DEBUG_CRREJECT_Y"));
1622 if (debugy < 1 || debugy > aGrid->ny) {
1626 if (getenv(
"MUSE_DEBUG_CRREJECT_Z")) {
1627 debugz = atoi(getenv(
"MUSE_DEBUG_CRREJECT_Z"));
1628 if (debugz < 1 || debugz > aGrid->nz) {
1640 enum dirtype dir = DIR_NONE;
1644 const double palimit = 5.;
1645 if (!haswcs || (fabs(posang) < palimit || fabs(fabs(posang) - 180.) < palimit ||
1646 fabs(fabs(posang) - 360.) < palimit)) {
1647 cpl_msg_debug(
id,
"CR rejection: posang = %f deg --> DIR_X " 1648 "(%s / %s / %s / %s)", posang, haswcs ?
"yes":
"no",
1649 fabs(posang) < palimit ?
"true" :
"false",
1650 fabs(fabs(posang) - 180.) < palimit ?
"true" :
"false",
1651 fabs(fabs(posang) - 360.) < palimit ?
"true" :
"false");
1653 }
else if (fabs(fabs(posang) - 90.) < palimit ||
1654 fabs(fabs(posang) - 270.) < palimit) {
1655 cpl_msg_debug(
id,
"CR rejection: posang = %f deg --> DIR_Y " 1656 "(%s / %s / %s)", posang, haswcs ?
"yes":
"no",
1657 fabs(fabs(posang) - 90.) < palimit ?
"true" :
"false",
1658 fabs(fabs(posang) - 270.) < palimit ?
"true" :
"false");
1661 cpl_msg_debug(
id,
"CR rejection: posang = %f deg --> DIR_NONE " 1662 "(%s / %s / %s / %s / %s / %s)", posang, haswcs ?
"yes":
"no",
1663 fabs(posang) < palimit ?
"true" :
"false",
1664 fabs(fabs(posang) - 90.) < palimit ?
"true" :
"false",
1665 fabs(fabs(posang) - 180.) < palimit ?
"true" :
"false",
1666 fabs(fabs(posang) - 270.) < palimit ?
"true" :
"false",
1667 fabs(fabs(posang) - 360.) < palimit ?
"true" :
"false");
1676 uint32_t badinclude = EURO3D_COSMICRAY;
1678 #ifdef ESO_ENABLE_DEBUG 1679 #pragma omp parallel for default(none) \ 1680 shared(aGrid, aPixtable, aSigma, badinclude, aType, crtypestring, \ 1681 data, debug, debugx, debugy, debugz, dir, dq, id, stat, stdout) 1683 #pragma omp parallel for default(none) \ 1684 shared(aGrid, aSigma, aType, badinclude, crtypestring, data, dir, \ 1687 for (l = 0; l < aGrid->nz; l++) {
1688 #define MUSE_CRREJECT_MAX_NPIX 1000 1689 int sxsize = MUSE_CRREJECT_MAX_NPIX;
1690 cpl_image *sdata = NULL;
1691 double *vdata = NULL;
1695 sdata = cpl_image_new(sxsize, 2, CPL_TYPE_DOUBLE);
1696 vdata = cpl_image_get_data_double(sdata);
1700 for (i = 0; i < aGrid->nx; i++) {
1702 for (j = 0; j < aGrid->ny; j++) {
1705 double level = 0., dev = 0;
1708 cpl_size npixels = 0;
1711 for (i2 = i - CR_LD; i2 <= i + CR_LD; i2++) {
1716 for (j2 = j - nwidth; j2 <= j + nwidth; j2++) {
1717 for (l2 = l - nwidth; l2 <= l + nwidth; l2++) {
1718 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_START
1719 #ifdef ESO_ENABLE_DEBUG 1720 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_DEBUG
1722 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END
1727 for (j2 = j - CR_LD; j2 <= j + CR_LD; j2++) {
1729 if (dir == DIR_X && j2 == j) {
1732 for (i2 = i - nwidth; i2 <= i + nwidth; i2++) {
1733 for (l2 = l - nwidth; l2 <= l + nwidth; l2++) {
1734 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_START
1735 #ifdef ESO_ENABLE_DEBUG 1736 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_DEBUG
1738 MUSE_RESAMPLING_CRREJECT_COMPUTE_STATS_END
1749 dev = sqrt(dev) / npixels;
1753 sflags = CPL_STATS_MEDIAN | CPL_STATS_MAD;
1755 sflags = CPL_STATS_MEAN | CPL_STATS_STDEV;
1757 cpl_stats *sstats = cpl_stats_new_from_image_window(sdata, sflags,
1760 level = cpl_stats_get_median(sstats);
1761 dev = cpl_stats_get_mad(sstats);
1763 level = cpl_stats_get_mean(sstats);
1764 dev = cpl_stats_get_stdev(sstats);
1766 cpl_stats_delete(sstats);
1768 double limit = level + aSigma * dev;
1769 #ifdef ESO_ENABLE_DEBUG 1770 if (debug & 2 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
1772 printf(
"%s: %03d,%03d,%04d: %.3f+/-%.3f (stats), %.3f (limit) (%" 1773 CPL_SIZE_FORMAT
" values)\n", __func__, i+1, j+1, l+1, level,
1774 dev, limit, npixels);
1776 cpl_stats *ssdata = cpl_stats_new_from_image_window(sdata, CPL_STATS_ALL,
1778 *ssstat = cpl_stats_new_from_image_window(sdata, CPL_STATS_ALL,
1780 printf(
"%s: %03d,%03d,%04d: %e +/- %e (%s), %e (limit) (%"CPL_SIZE_FORMAT
1781 " values); data stats:\n", __func__, i+1, j+1, l+1,
1782 level, dev, crtypestring[aType], limit, npixels);
1783 cpl_stats_dump(ssdata, CPL_STATS_ALL, stdout);
1784 printf(
"%s: variance stats:\n", __func__);
1785 cpl_stats_dump(ssstat, CPL_STATS_ALL, stdout);
1787 cpl_stats_delete(ssdata);
1788 cpl_stats_delete(ssstat);
1798 for (irow = 0; irow < nrow; irow++) {
1799 if (data[rows[irow]] > limit) {
1800 dq[rows[irow]] |= EURO3D_COSMICRAY;
1801 #ifdef ESO_ENABLE_DEBUG 1802 if (debug & 1 && i+1 == debugx && j+1 == debugy && l+1 == debugz) {
1803 printf(
"%s: %03d,%03d,%04d: rejected row %"CPL_SIZE_FORMAT
" (%" 1804 CPL_SIZE_FORMAT
" of %"CPL_SIZE_FORMAT
" in this gridcell):\t",
1805 __func__, i+1, j+1, l+1, rows[irow], irow+1, nrow);
1811 #ifdef ESO_ENABLE_DEBUG 1818 cpl_image_delete(sdata);
1850 cpl_ensure(aParams, CPL_ERROR_NULL_INPUT, NULL);
1862 e3d->
header = cpl_propertylist_duplicate(cube->
header);
1863 cpl_propertylist_erase_regexp(e3d->
header,
1864 "^SIMPLE$|^BITPIX$|^NAXIS|^EURO3D$|^E3D_",
1866 cpl_propertylist_append_char(e3d->
header,
"EURO3D",
'T');
1867 cpl_propertylist_set_comment(e3d->
header,
"EURO3D",
1868 "file conforms to Euro3D standard");
1870 cpl_propertylist_append_string(e3d->
header,
"E3D_VERS",
"1.0");
1871 cpl_propertylist_set_comment(e3d->
header,
"E3D_VERS",
1872 "version number of the Euro3D format");
1873 cpl_errorstate prestate = cpl_errorstate_get();
1874 double darlambdaref = cpl_propertylist_get_double(e3d->
header,
1876 if (!cpl_errorstate_is_equal(prestate)) {
1878 cpl_errorstate_set(prestate);
1880 if (darlambdaref > 0.) {
1881 cpl_propertylist_append_char(e3d->
header,
"E3D_ADC",
'T');
1882 cpl_propertylist_set_comment(e3d->
header,
"E3D_ADC",
1883 "data was corrected for atmospheric dispersion");
1885 cpl_propertylist_append_char(e3d->
header,
"E3D_ADC",
'F');
1886 cpl_propertylist_set_comment(e3d->
header,
"E3D_ADC",
1887 "data not corrected for atmospheric dispersion");
1891 e3d->
hdata = cpl_propertylist_new();
1892 cpl_propertylist_append_string(e3d->
hdata,
"EXTNAME",
"E3D_DATA");
1893 cpl_propertylist_set_comment(e3d->
hdata,
"EXTNAME",
1894 "This is the Euro3D data table extension");
1895 cpl_propertylist_append_string(e3d->
hdata,
"CTYPES",
1897 cpl_propertylist_set_comment(e3d->
hdata,
"CTYPES",
1898 cpl_propertylist_get_comment(e3d->
header,
"CTYPE3"));
1899 cpl_propertylist_append_string(e3d->
hdata,
"CUNITS",
1901 cpl_propertylist_set_comment(e3d->
hdata,
"CUNITS",
1902 cpl_propertylist_get_comment(e3d->
header,
"CUNIT3"));
1903 cpl_propertylist_append_double(e3d->
hdata,
"CRVALS",
1905 cpl_propertylist_set_comment(e3d->
hdata,
"CRVALS",
1906 cpl_propertylist_get_comment(e3d->
header,
"CRVAL3"));
1907 cpl_propertylist_append_double(e3d->
hdata,
"CDELTS",
1909 cpl_propertylist_set_comment(e3d->
hdata,
"CDELTS",
1910 cpl_propertylist_get_comment(e3d->
header,
"CD3_3"));
1912 cpl_propertylist_erase_regexp(e3d->
header,
"^C...*3$|^CD3_.$", 0);
1915 int nx = cpl_image_get_size_x(cpl_imagelist_get(cube->
data, 0)),
1916 ny = cpl_image_get_size_y(cpl_imagelist_get(cube->
data, 0)),
1917 nz = cpl_imagelist_get_size(cube->
data);
1922 cpl_table_set_column_depth(e3d->
dtable,
"DATA_SPE", nz);
1923 cpl_table_set_column_depth(e3d->
dtable,
"QUAL_SPE", nz);
1924 cpl_table_set_column_depth(e3d->
dtable,
"STAT_SPE", nz);
1926 cpl_table_set_column_unit(e3d->
dtable,
"DATA_SPE",
1927 cpl_table_get_column_unit(aPixtable->
table,
1929 cpl_table_set_column_unit(e3d->
dtable,
"STAT_SPE",
1930 cpl_table_get_column_unit(aPixtable->
table,
1933 cpl_table_set_column_savetype(e3d->
dtable,
"SELECTED", CPL_TYPE_BOOL);
1938 cpl_table_set_column_unit(e3d->
dtable,
"XPOS",
"deg");
1939 cpl_table_set_column_unit(e3d->
dtable,
"YPOS",
"deg");
1943 unsigned euro3d_ignore = EURO3D_OUTSDRANGE | EURO3D_MISSDATA;
1944 cpl_vector *vusable = cpl_vector_new(nx * ny);
1946 for (i = 1; i <= nx; i++) {
1948 for (j = 1; j <= ny; j++) {
1949 cpl_array *adata = cpl_array_new(nz, CPL_TYPE_FLOAT),
1950 *adq = cpl_array_new(nz, CPL_TYPE_INT),
1951 *astat = cpl_array_new(nz, CPL_TYPE_FLOAT);
1960 int l, nusable = 0, nstart = -1;
1961 for (l = 0; l < nz; l++) {
1963 unsigned dq = cpl_image_get(cpl_imagelist_get(cube->
dq, l), i, j, &err);
1965 if (nstart < 0 && (dq & euro3d_ignore)) {
1968 cpl_array_set_int(adq, nusable, dq);
1969 cpl_array_set_float(adata, nusable,
1970 cpl_image_get(cpl_imagelist_get(cube->
data, l),
1973 cpl_array_set_float(astat, nusable,
1974 sqrt(cpl_image_get(cpl_imagelist_get(cube->
stat, l),
1982 cpl_table_set_int(e3d->
dtable,
"SPEC_ID", itable, itable + 1);
1983 cpl_table_set_int(e3d->
dtable,
"SPEC_LEN", itable, nusable);
1984 cpl_table_set_int(e3d->
dtable,
"SPEC_STA", itable, nstart);
1985 cpl_table_set_double(e3d->
dtable,
"XPOS", itable, x);
1986 cpl_table_set_double(e3d->
dtable,
"YPOS", itable, y);
1989 cpl_table_set_array(e3d->
dtable,
"DATA_SPE", itable, adata);
1990 cpl_table_set_array(e3d->
dtable,
"QUAL_SPE", itable, adq);
1991 cpl_table_set_array(e3d->
dtable,
"STAT_SPE", itable, astat);
1993 cpl_array_delete(adata);
1994 cpl_array_delete(adq);
1995 cpl_array_delete(astat);
1997 cpl_vector_set(vusable, itable, nusable);
1998 if (nstart != -1 && nusable > 0) {
2000 cpl_table_unselect_row(e3d->
dtable, itable);
2007 cpl_vector_set_size(vusable, itable);
2008 int maxusable = cpl_vector_get_max(vusable);
2010 cpl_msg_debug(__func__,
"filled %d of %d spaxels, usable are " 2011 "%f..%f(%f)+/-%f..%d pixels per spectrum", itable, nx * ny,
2012 cpl_vector_get_min(vusable), cpl_vector_get_mean(vusable),
2013 cpl_vector_get_median(vusable), cpl_vector_get_stdev(vusable),
2016 cpl_msg_debug(__func__,
"filled %"CPL_SIZE_FORMAT
" of %d spaxels, usable " 2017 "are max. %d of %d pixels per spectrum%s",
2018 itable - cpl_table_count_selected(e3d->
dtable), nx * ny,
2019 maxusable, nz, maxusable == nz ?
"" :
", cutting table");
2020 if (maxusable != nz) {
2022 cpl_table_set_column_depth(e3d->
dtable,
"DATA_SPE", maxusable);
2023 cpl_table_set_column_depth(e3d->
dtable,
"QUAL_SPE", maxusable);
2024 cpl_table_set_column_depth(e3d->
dtable,
"STAT_SPE", maxusable);
2027 cpl_table_erase_selected(e3d->
dtable);
2028 cpl_vector_delete(vusable);
2032 cpl_table_fill_column_window_int(e3d->
dtable,
"SELECTED", 0, itable,
2035 cpl_table_fill_column_window_int(e3d->
dtable,
"NSPAX", 0, itable, 1);
2037 cpl_table_fill_column_window_int(e3d->
dtable,
"GROUP_N", 0, itable, 1);
2043 e3d->
hgroup = cpl_propertylist_new();
2044 cpl_propertylist_append_string(e3d->
hgroup,
"EXTNAME",
"E3D_GRP");
2045 cpl_propertylist_set_comment(e3d->
hgroup,
"EXTNAME",
2046 "This is the Euro3D group table extension");
2048 cpl_table_set_int(e3d->
gtable,
"GROUP_N", 0, 1);
2049 cpl_table_set_string(e3d->
gtable,
"G_SHAPE", 0,
"R");
2050 cpl_table_set_float(e3d->
gtable,
"G_SIZE1", 0, aParams->
dx);
2051 cpl_table_set_float(e3d->
gtable,
"G_ANGLE", 0, 0.);
2052 cpl_table_set_float(e3d->
gtable,
"G_SIZE2", 0, aParams->dy);
2053 if (darlambdaref > 0.) {
2055 cpl_table_set_float(e3d->
gtable,
"G_POSWAV", 0, NAN);
2056 cpl_table_set_float(e3d->
gtable,
"G_AIRMAS", 0, NAN);
2057 cpl_table_set_float(e3d->
gtable,
"G_PARANG", 0, NAN);
2058 cpl_table_set_float(e3d->
gtable,
"G_PRESSU", 0, NAN);
2059 cpl_table_set_float(e3d->
gtable,
"G_TEMPER", 0, NAN);
2060 cpl_table_set_float(e3d->
gtable,
"G_HUMID", 0, NAN);
2063 cpl_table_set_float(e3d->
gtable,
"G_POSWAV", 0,
2064 (kMuseNominalLambdaMin + kMuseNominalLambdaMax) / 2.);
2065 cpl_table_set_float(e3d->
gtable,
"G_AIRMAS", 0,
2067 cpl_table_set_float(e3d->
gtable,
"G_PARANG", 0,
2069 cpl_table_set_float(e3d->
gtable,
"G_PRESSU", 0,
2072 cpl_table_set_float(e3d->
gtable,
"G_TEMPER", 0,
2074 cpl_table_set_float(e3d->
gtable,
"G_HUMID", 0,
2137 cpl_ensure(aPixtable && aParams, CPL_ERROR_NULL_INPUT, NULL);
2139 CPL_ERROR_ILLEGAL_INPUT, NULL);
2147 muse_resampling_check_deltas(aPixtable, aParams);
2151 int xsize = 0, ysize = 0, zsize = 0;
2152 muse_resampling_compute_size(aPixtable, aParams, &xsize, &ysize, &zsize);
2153 muse_resampling_override_size(&xsize, &ysize, &zsize, aParams);
2154 cpl_ensure(xsize > 0 && ysize > 0 && zsize > 0, CPL_ERROR_ILLEGAL_OUTPUT,
2157 double time = cpl_test_get_walltime();
2163 cube->
header = cpl_propertylist_duplicate(aPixtable->
header);
2164 cpl_propertylist_erase_regexp(cube->
header,
2165 "^SIMPLE$|^BITPIX$|^NAXIS|^EXTEND$|^XTENSION$|" 2166 "^DATASUM$|^DATAMIN$|^DATAMAX$|^DATAMD5$|" 2167 "^PCOUNT$|^GCOUNT$|^HDUVERS$|^BLANK$|" 2168 "^BZERO$|^BSCALE$|^CHECKSUM$|^INHERIT$|" 2172 cpl_propertylist_update_string(cube->
header,
"BUNIT",
2173 cpl_table_get_column_unit(aPixtable->
table,
2176 cpl_propertylist_update_int(cube->
header,
"NAXIS", 3);
2177 cpl_propertylist_update_int(cube->
header,
"NAXIS1", xsize);
2178 cpl_propertylist_update_int(cube->
header,
"NAXIS2", ysize);
2179 cpl_propertylist_update_int(cube->
header,
"NAXIS3", zsize);
2183 cpl_propertylist_copy_property_regexp(cube->
header, aPixtable->
header,
2185 cpl_propertylist_update_double(cube->
header,
"CD1_1", -aParams->
dx);
2186 cpl_propertylist_update_double(cube->
header,
"CD2_2", aParams->dy);
2187 cpl_propertylist_update_double(cube->
header,
"CD1_2", 0.);
2188 cpl_propertylist_update_double(cube->
header,
"CD2_1", 0.);
2190 cpl_propertylist_update_double(cube->
header,
"CRPIX1",
2192 + (1. + xsize) / 2.);
2193 cpl_propertylist_update_double(cube->
header,
"CRPIX2",
2195 + (1. + ysize) / 2.);
2197 cpl_propertylist_erase(cube->
header,
"WCSAXES");
2200 cpl_propertylist_append_double(cube->
header,
"CRPIX1", 1.);
2201 cpl_propertylist_append_double(cube->
header,
"CRVAL1",
2202 cpl_propertylist_get_float(aPixtable->
header,
2206 cpl_propertylist_append_string(cube->
header,
"CTYPE1",
"PIXEL");
2207 cpl_propertylist_append_string(cube->
header,
"CUNIT1",
"pixel");
2208 cpl_propertylist_append_double(cube->
header,
"CD1_1", aParams->
dx);
2209 cpl_propertylist_append_double(cube->
header,
"CRPIX2", 1.);
2210 cpl_propertylist_append_double(cube->
header,
"CRVAL2",
2211 cpl_propertylist_get_float(aPixtable->
header,
2213 cpl_propertylist_append_string(cube->
header,
"CTYPE2",
"PIXEL");
2214 cpl_propertylist_append_string(cube->
header,
"CUNIT2",
"pixel");
2215 cpl_propertylist_append_double(cube->
header,
"CD2_2", aParams->dy);
2217 cpl_propertylist_append_double(cube->
header,
"CD1_2", 0.);
2218 cpl_propertylist_append_double(cube->
header,
"CD2_1", 0.);
2222 cpl_propertylist_append_string(cube->
header,
"CTYPE3",
"AWAV");
2225 cpl_propertylist_append_string(cube->
header,
"CTYPE3",
"AWAV-LOG");
2228 cpl_propertylist_append_string(cube->
header,
"CTYPE3",
"WAVE");
2231 cpl_propertylist_append_string(cube->
header,
"CTYPE3",
"WAVE-LOG");
2234 cpl_propertylist_append_string(cube->
header,
"CTYPE3",
"UNKNOWN");
2236 cpl_propertylist_append_string(cube->
header,
"CUNIT3",
"Angstrom");
2237 cpl_propertylist_append_double(cube->
header,
"CD3_3", aParams->dlambda);
2238 cpl_propertylist_append_double(cube->
header,
"CRPIX3", 1.);
2239 cpl_propertylist_append_double(cube->
header,
"CRVAL3",
2240 cpl_propertylist_get_float(aPixtable->
header,
2243 cpl_propertylist_append_double(cube->
header,
"CD1_3", 0.);
2244 cpl_propertylist_append_double(cube->
header,
"CD2_3", 0.);
2245 cpl_propertylist_append_double(cube->
header,
"CD3_1", 0.);
2246 cpl_propertylist_append_double(cube->
header,
"CD3_2", 0.);
2249 cpl_propertylist_append_bool(cube->
header, MUSE_HDR_FLUX_FFCORR, CPL_TRUE);
2250 cpl_propertylist_set_comment(cube->
header, MUSE_HDR_FLUX_FFCORR,
2251 MUSE_HDR_FLUX_FFCORR_C);
2255 muse_resampling_override_wcs(cube, aParams);
2257 muse_resampling_fit_data_to_grid(cube, aPixtable);
2261 cube->
data = cpl_imagelist_new();
2262 cube->
dq = cpl_imagelist_new();
2263 cube->
stat = cpl_imagelist_new();
2265 for (i = 0; i < zsize; i++) {
2266 cpl_image *image = cpl_image_new(xsize, ysize, CPL_TYPE_FLOAT);
2268 cpl_imagelist_set(cube->
data, image, i);
2270 cpl_imagelist_set(cube->
stat, cpl_image_duplicate(image), i);
2274 image = cpl_image_new(xsize, ysize, CPL_TYPE_INT);
2276 cpl_image_add_scalar(image, EURO3D_OUTSDRANGE);
2277 cpl_imagelist_set(cube->
dq, image, i);
2284 xsize, ysize, zsize);
2287 muse_resampling_crreject(aPixtable, grid, aParams->
crtype,
2292 double timeinit = cpl_test_get_walltime(),
2293 cpuinit = cpl_test_get_cputime();
2296 cpl_error_code rc = CPL_ERROR_NONE;
2297 switch (aParams->
method) {
2299 cpl_msg_info(__func__,
"Starting resampling, using method \"nearest\"");
2300 rc = muse_resampling_cube_nearest(cube, aPixtable, grid);
2303 cpl_msg_info(__func__,
"Starting resampling, using method \"renka\" " 2304 "(critical radius rc=%f, loop distance ld=%d)", aParams->
rc,
2306 rc = muse_resampling_cube_weighted(cube, aPixtable, grid, aParams);
2311 cpl_msg_info(__func__,
"Starting resampling, using method \"%s\" (loop " 2319 rc = muse_resampling_cube_weighted(cube, aPixtable, grid, aParams);
2322 cpl_msg_info(__func__,
"Starting resampling, using method \"drizzle\" " 2323 "(pixfrac f=%.3f,%.3f,%.3f, loop distance ld=%d)",
2324 aParams->
pfx, aParams->pfy, aParams->pfl, aParams->
ld);
2325 rc = muse_resampling_cube_weighted(cube, aPixtable, grid, aParams);
2328 cpl_msg_debug(__func__,
"Method %d (no resampling)", aParams->
method);
2331 cpl_msg_error(__func__,
"Don't know this resampling method: %d",
2333 cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
2334 rc = CPL_ERROR_UNSUPPORTED_MODE;
2338 double timefini = cpl_test_get_walltime(),
2339 cpufini = cpl_test_get_cputime();
2348 cpl_msg_debug(__func__,
"resampling took %.3fs (wall-clock) and %.3fs " 2349 "(%.3fs CPU, %d CPUs) for muse_resampling_cube*() alone",
2350 timefini - time, timefini - timeinit, cpufini - cpuinit,
2351 omp_get_max_threads());
2352 if (rc != CPL_ERROR_NONE) {
2353 cpl_msg_error(__func__,
"resampling failed: %s", cpl_error_get_message());
2398 cpl_ensure(aPixtable && aPixtable->
table && aGrid && aParams &&
2399 aCube && aCube->
header, CPL_ERROR_NULL_INPUT, NULL);
2406 cpl_errorstate prestate = cpl_errorstate_get();
2413 if (!cpl_errorstate_is_equal(prestate)) {
2414 cpl_errorstate_set(prestate);
2422 double xnorm = 1., ynorm = 1.,
2434 double renka_rc = aParams->
rc 2435 * sqrt(pow(wcs->cd11*xnorm, 2) + pow(wcs->
cd22*xnorm, 2));
2437 int ld = aParams->
ld;
2440 cpl_msg_info(__func__,
"Overriding loop distance ld=%d", ld);
2444 double xsz = aParams->
pfx / xnorm,
2445 ysz = aParams->pfy / ynorm,
2446 xout = fabs(wcs->cd11), yout = fabs(wcs->
cd22);
2449 image->
data = cpl_image_new(aGrid->nx, aGrid->ny, CPL_TYPE_FLOAT);
2450 image->
dq = cpl_image_new(aGrid->nx, aGrid->ny, CPL_TYPE_INT);
2451 image->
stat = cpl_image_new(aGrid->nx, aGrid->ny, CPL_TYPE_FLOAT);
2452 image->
header = cpl_propertylist_duplicate(aCube->
header);
2453 cpl_propertylist_erase_regexp(image->
header,
"^C...*3$|^CD3_.$", 0);
2454 float *pdata = cpl_image_get_data_float(image->
data),
2455 *pstat = cpl_image_get_data_float(image->
stat);
2456 int *pdq = cpl_image_get_data_int(image->
dq);
2459 cpl_boolean usevariance = CPL_FALSE;
2460 if (getenv(
"MUSE_COLLAPSE_USE_VARIANCE") &&
2461 atoi(getenv(
"MUSE_COLLAPSE_USE_VARIANCE")) > 0) {
2462 usevariance = CPL_TRUE;
2466 cpl_table *filter = aFilter && aFilter->
table ? aFilter->
table : NULL;
2469 const double *flbda = cpl_table_get_data_double_const(filter,
"lambda"),
2470 *fresp = cpl_table_get_data_double_const(filter,
"throughput");
2471 int l = 0, nl = cpl_table_get_nrow(filter);
2472 while (l < nl && fabs(fresp[l]) < DBL_EPSILON) {
2476 while (l > 0 && fabs(fresp[l]) < DBL_EPSILON) {
2479 cpl_msg_debug(__func__,
"filter wavelength range %.1f..%.1fA", lmin, lmax);
2486 cpl_msg_debug(__func__,
"full wavelength range %.1f..%.1fA", lmin, lmax);
2490 #pragma omp parallel for default(none) \ 2491 shared(aParams, aGrid, data, dq, filter, lbda, ld, lmax, lmin, pdata,\ 2492 pdq, pstat, ptxoff, ptyoff, renka_rc, stat, usevariance, wcs, \ 2493 wght, xnorm, xout, xpos, xsz, ynorm, yout, ypos, ysz) 2494 for (i = 0; i < aGrid->nx; i++) {
2496 for (j = 0; j < aGrid->ny; j++) {
2504 double sumdata = 0, sumstat = 0, sumweight = 0;
2505 cpl_size npoints = 0;
2509 for (i2 = i - ld; i2 <= i + ld; i2++) {
2511 for (j2 = j - ld; j2 <= j + ld; j2++) {
2514 for (l2 = 0; l2 < aGrid->nz; l2++) {
2518 for (n = 0; n < n_rows2; n++) {
2522 if (usevariance && !isnormal(stat[rows2[n]])) {
2525 if (lbda[rows2[n]] > lmax || lbda[rows2[n]] < lmin) {
2529 double dx = fabs(x - (xpos[rows2[n]] + ptxoff)),
2530 dy = fabs(y - (ypos[rows2[n]] + ptyoff)),
2535 dx *= cos(y * CPL_MATH_RAD_DEG);
2544 weight = muse_resampling_weight_function_renka(sqrt(r2), renka_rc);
2546 weight = muse_resampling_weight_function_drizzle(xsz, ysz, 1.,
2550 weight = muse_resampling_weight_function_linear(sqrt(r2));
2552 weight = muse_resampling_weight_function_quadratic(r2);
2554 weight = muse_resampling_weight_function_lanczos(dx, dy, 0., ld);
2559 weight *= wght[rows2[n]];
2568 weight /= stat[rows2[n]];
2570 sumweight += weight;
2571 sumdata += data[rows2[n]] * weight;
2572 sumstat += stat[rows2[n]] * weight*weight;
2582 if (!npoints || !isnormal(sumweight)) {
2583 pdq[i + j * aGrid->nx] = EURO3D_MISSDATA;
2588 sumdata /= sumweight;
2589 sumstat /= sumweight*sumweight;
2590 pdata[i + j * aGrid->nx] = sumdata;
2591 pstat[i + j * aGrid->nx] = sumstat;
2592 pdq[i + j * aGrid->nx] = EURO3D_GOODPIXEL;
2622 static cpl_error_code
2626 cpl_ensure_code(aImage && aPixtable && aGrid, CPL_ERROR_NULL_INPUT);
2627 aImage->
data = cpl_image_new(aGrid->nx, aGrid->nz,
2636 float *imagedata = cpl_image_get_data_float(aImage->
data);
2638 #ifdef ESO_ENABLE_DEBUG 2640 if (getenv(
"MUSE_DEBUG_NEAREST")) {
2641 debug = atoi(getenv(
"MUSE_DEBUG_NEAREST"));
2646 for (i = 0; i < aGrid->nx; i++) {
2648 for (j = 0; j < aGrid->nz; j++) {
2652 if (n_rows == 1 && !dq[rows[0]]) {
2654 imagedata[i + j * aGrid->nx] = data[rows[0]];
2655 #ifdef ESO_ENABLE_DEBUG 2657 printf(
"only: %f\n", data[rows[0]]);
2661 }
else if (n_rows >= 2) {
2663 cpl_size n, nbest = -1;
2664 double dbest = FLT_MAX;
2665 for (n = 0; n < n_rows; n++) {
2671 double dlambda = fabs(j * cd22 + crval2 - lbda[rows[n]]);
2672 #ifdef ESO_ENABLE_DEBUG 2674 printf(
"%f, %f, d: %f best: %f (%f)\n",
2675 j * cd22 + crval2, lbda[rows[n]],
2676 dlambda, dbest, data[rows[n]]);
2679 if (dlambda < dbest) {
2684 imagedata[i + j * aGrid->nx] = data[rows[nbest]];
2685 #ifdef ESO_ENABLE_DEBUG 2687 printf(
"closest: %f\n", data[rows[nbest]]);
2697 return CPL_ERROR_NONE;
2718 static cpl_error_code
2722 cpl_ensure_code(aImage && aPixtable && aGrid, CPL_ERROR_NULL_INPUT);
2723 aImage->
data = cpl_image_new(aGrid->nx, aGrid->nz,
2732 float *imagedata = cpl_image_get_data_float(aImage->
data);
2734 #ifdef ESO_ENABLE_DEBUG 2735 int debug = 0, debugx = 0, debugz = 0;
2736 if (getenv(
"MUSE_DEBUG_WEIGHTED")) {
2737 debug = atoi(getenv(
"MUSE_DEBUG_WEIGHTED"));
2740 if (getenv(
"MUSE_DEBUG_WEIGHTED_X")) {
2741 debugx = atoi(getenv(
"MUSE_DEBUG_WEIGHTED_X"));
2742 if (debugx < 1 || debugx > aGrid->nx) {
2746 if (getenv(
"MUSE_DEBUG_WEIGHTED_Z")) {
2747 debugz = atoi(getenv(
"MUSE_DEBUG_WEIGHTED_Z"));
2748 if (debugz < 1 || debugz > aGrid->nz) {
2758 double renka_rc = 1.25;
2763 for (i = 0; i < aGrid->nx; i++) {
2765 for (j = 0; j < aGrid->nz; j++) {
2766 double sumdata = 0, sumweight = 0,
2767 lambda = j * cd22 + crval2;
2769 #ifdef ESO_ENABLE_DEBUG 2770 #define MAX_NPIX_2D 50 2771 int *pointlist = NULL;
2772 double *pointweights = NULL;
2773 if (debug & 128 && i+1 == debugx && j+1 == debugz) {
2774 pointlist = cpl_calloc(MAX_NPIX_2D,
sizeof(
int));
2775 pointweights = cpl_malloc(MAX_NPIX_2D *
sizeof(
double));
2781 for (j2 = j - ld; j2 <= j + ld; j2++) {
2785 for (n = 0; n < n_rows; n++) {
2789 double dlambda = fabs(lambda - lbda[rows[n]]),
2790 weight = muse_resampling_weight_function_renka(dlambda,
2792 sumweight += weight;
2793 sumdata += data[rows[n]] * weight;
2795 #ifdef ESO_ENABLE_DEBUG 2796 if (debug & 128 && i+1 == debugx && j+1 == debugz &&
2797 npoints-1 < MAX_NPIX_2D) {
2799 pointlist[npoints-1] = rows[n] + 1;
2800 pointweights[npoints-1] = weight;
2804 printf(
" pixel %4d,%4d (%8"CPL_SIZE_FORMAT
"): %2"CPL_SIZE_FORMAT
2805 " %2"CPL_SIZE_FORMAT
" %f, %f -> %f ==> %f %f (%d)\n",
2808 weight, sumweight, sumdata, npoints);
2814 #ifdef ESO_ENABLE_DEBUG 2815 if (debug & 128 && i+1 == debugx && j+1 == debugz) {
2816 printf(
"pixelnumber weight ");
2819 while (++m < MAX_NPIX_2D && pointlist[m] != 0) {
2821 printf(
"%12d %8.5f ", pointlist[m] - 1, pointweights[m]);
2825 cpl_free(pointlist);
2826 cpl_free(pointweights);
2838 imagedata[i + j * aGrid->nx] = sumdata / sumweight;
2842 return CPL_ERROR_NONE;
2866 double aDX,
double aLambdaMin,
2867 double aLambdaMax,
double aDLambda)
2869 double dlambda = aDLambda;
2872 aLambdaMin, aLambdaMax,
2878 image->
header = cpl_propertylist_new();
2879 const char *unit = cpl_table_get_column_unit(aPixtable->
table,
"data");
2880 cpl_propertylist_append_string(image->
header,
"BUNIT", unit);
2884 cpl_propertylist_copy_property_regexp(image->
header, aPixtable->
header,
2886 cpl_propertylist_append_double(image->
header,
"CRPIX1", 1.);
2887 cpl_propertylist_append_double(image->
header,
"CRPIX2", 1.);
2888 cpl_propertylist_append_double(image->
header,
"CRVAL1", 1.);
2889 cpl_propertylist_append_double(image->
header,
"CRVAL2", aLambdaMin);
2890 cpl_propertylist_append_double(image->
header,
"CD1_1", 1.);
2891 cpl_propertylist_append_double(image->
header,
"CD2_2", dlambda);
2893 cpl_propertylist_append_string(image->
header,
"CUNIT1",
"pixel");
2894 cpl_propertylist_append_string(image->
header,
"CUNIT2",
"Angstrom");
2897 cpl_propertylist_append_string(image->
header,
"CTYPE1",
"PIXEL");
2900 cpl_propertylist_append_string(image->
header,
"CTYPE2",
"AWAV");
2902 cpl_propertylist_append_double(image->
header,
"CD1_2", 0.);
2903 cpl_propertylist_append_double(image->
header,
"CD2_1", 0.);
2906 cpl_error_code rc = CPL_ERROR_NONE;
2909 rc = muse_resampling_image_nearest(image, aPixtable, grid);
2913 rc = muse_resampling_image_weighted(image, aPixtable, grid);
2917 if (rc != CPL_ERROR_NONE) {
2918 cpl_msg_error(__func__,
"resampling failed: %s", cpl_error_get_message());
2969 double aDX,
double aDLambda)
2971 cpl_ensure(aPixtable, CPL_ERROR_NULL_INPUT, NULL);
2972 if (aDLambda == 0.0) {
2973 aDLambda = kMuseSpectralSamplingA;
2983 cpl_msg_info(__func__,
"Using nearest neighbor sampling (%d) along " 2984 "wavelengths.", aMethod);
2987 cpl_msg_info(__func__,
"Using renka-weighted interpolation (%d) along " 2988 "wavelengths.", aMethod);
2991 cpl_msg_error(__func__,
"Don't know this resampling method: %d", aMethod);
2992 cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
3005 return muse_resampling_image_selected(aPixtable, aMethod,
3006 aDX == 0.0 ? 1. : aDX,
3007 lmin, lmax, aDLambda);
3027 cpl_msg_debug(__func__,
"Resampling %d slices to a 2D image, using bins of" 3028 " %e %s x %.3f Angstrom", n_slices, dx,
3035 #pragma omp parallel for default(none) \ 3036 shared(aDLambda, aMethod, dx, lmax, lmin, n_slices, \ 3037 slice_image, slice_pixtable) 3038 for (i_slice = 0; i_slice < n_slices; i_slice++) {
3040 cpl_msg_debug(__func__,
"slice pixel table %d", i_slice + 1);
3044 slice_image[i_slice] = NULL;
3047 slice_image[i_slice] = muse_resampling_image_selected(pt, aMethod, dx,
3048 lmin, lmax, aDLambda);
3053 for (i_slice = 0; i_slice < n_slices; i_slice++) {
3056 cpl_msg_debug(__func__,
"slice %d: %ldx%ld", i_slice + 1,
3057 cpl_image_get_size_x(slice->
data), cpl_image_get_size_y(slice->
data));
3059 if (slice == NULL) {
3064 image->
header = cpl_propertylist_duplicate(slice->
header);
3067 cpl_image_delete(image->
data);
3071 cpl_image_delete(image->
dq);
3076 cpl_image_delete(image->
stat);
3080 slice_image[i_slice] = NULL;
3086 cpl_propertylist_append_bool(image->
header, MUSE_HDR_FLUX_FFCORR, CPL_TRUE);
3087 cpl_propertylist_set_comment(image->
header, MUSE_HDR_FLUX_FFCORR,
3088 MUSE_HDR_FLUX_FFCORR_C);
3114 cpl_ensure(aPixtable && aPixtable->
header && aPixtable->
table,
3115 CPL_ERROR_NULL_INPUT, NULL);
3117 == CPL_ERROR_NONE, CPL_ERROR_ILLEGAL_INPUT, NULL);
3122 cpl_size lsize = floor((lmax - lmin) / aBinwidth) + 2;
3126 cpl_table_fill_column_window(spectrum,
"lambda", 0, lsize, 0);
3127 cpl_table_fill_column_window(spectrum,
"data", 0, lsize, 0);
3128 cpl_table_fill_column_window(spectrum,
"stat", 0, lsize, 0);
3129 cpl_table_fill_column_window(spectrum,
"dq", 0, lsize, 0);
3130 double *spec_data = cpl_table_get_data_double(spectrum,
"data");
3131 double *spec_stat = cpl_table_get_data_double(spectrum,
"stat");
3132 double *spec_lbda = cpl_table_get_data_double(spectrum,
"lambda");
3134 cpl_table_set_column_unit(spectrum,
"data",
3135 cpl_table_get_column_unit(aPixtable->
table,
3137 cpl_table_set_column_unit(spectrum,
"stat",
3138 cpl_table_get_column_unit(aPixtable->
table,
3142 cpl_table_new_column(spectrum,
"weight", CPL_TYPE_DOUBLE);
3143 cpl_table_fill_column_window(spectrum,
"weight", 0, lsize, 0);
3144 double *spec_weight = cpl_table_get_data_double(spectrum,
"weight");
3158 cpl_array *asel = cpl_table_where_selected(aPixtable->
table);
3159 const cpl_size *sel = cpl_array_get_data_cplsize_const(asel);
3160 cpl_size isel, nsel = cpl_array_get_size(asel);
3161 for (isel = 0; isel < nsel; isel++) {
3162 cpl_size i_row = sel[isel];
3163 if ((dq[i_row] != EURO3D_GOODPIXEL)) {
3168 if (! isfinite(data[i_row])) {
3172 double l = (lbda[i_row] - lmin) / aBinwidth;
3174 cpl_size il = floor(l);
3183 spec_data[il] += w0 * data[i_row];
3184 spec_data[il+1] += w1 * data[i_row];
3185 spec_stat[il] += w0 * stat[i_row];
3186 spec_stat[il+1] += w1 * stat[i_row];
3187 spec_weight[il] += w0;
3188 spec_weight[il+1] += w1;
3190 cpl_array_delete(asel);
3194 for (ispec = 0; ispec < lsize; ispec++) {
3195 if (spec_weight[ispec] > 0) {
3196 spec_lbda[ispec] = ispec * aBinwidth + lmin;
3197 cpl_table_unselect_row(spectrum, ispec);
3200 cpl_table_erase_selected(spectrum);
3203 cpl_table_divide_columns(spectrum,
"data",
"weight");
3204 cpl_table_divide_columns(spectrum,
"stat",
"weight");
3205 cpl_table_erase_column(spectrum,
"weight");
3243 float aLo,
float aHi,
unsigned char aIter)
3245 cpl_ensure(aPixtable && aPixtable->
header && aPixtable->
table,
3246 CPL_ERROR_NULL_INPUT, NULL);
3248 == CPL_ERROR_NONE, CPL_ERROR_ILLEGAL_INPUT, NULL);
3261 cpl_array *asel = cpl_table_where_selected(aPixtable->
table);
3262 const cpl_size *sel = cpl_array_get_data_cplsize_const(asel);
3263 cpl_size isel, nsel = cpl_array_get_size(asel);
3266 cpl_size nhi = 0, nlo = 0;
3268 for (i = 1; i <= aIter; i++) {
3269 cpl_size ispec, nbins = cpl_table_get_nrow(spectrum);
3270 double *sdata = cpl_table_get_data_double(spectrum,
"data"),
3271 *sstat = cpl_table_get_data_double(spectrum,
"stat"),
3272 *sstddev = cpl_malloc(nbins *
sizeof(
double));
3274 for (ispec = 0; ispec < nbins; ispec++) {
3275 sstddev[ispec] = sqrt(sstat[ispec]);
3279 for (isel = 0; isel < nsel; isel++) {
3280 cpl_size irow = sel[isel];
3281 if ((dq[irow] != EURO3D_GOODPIXEL)) {
3284 double l = lbda[irow];
3287 if ((ispec < nbins - 1) && (sdata[ispec] < sdata[ispec + 1])) {
3291 if (aHi > 0. && data[irow] > sdata[ispec] + aHi * sstddev[ispec]) {
3292 dq[irow] = EURO3D_OUTLIER;
3296 if (aLo > 0. && data[irow] < sdata[ispec] - aLo * sstddev[ispec]) {
3297 dq[irow] = EURO3D_OUTLIER;
3305 cpl_msg_debug(__func__,
"%"CPL_SIZE_FORMAT
" of %"CPL_SIZE_FORMAT
" pixels " 3306 "are outliers (%"CPL_SIZE_FORMAT
" low and %"CPL_SIZE_FORMAT
3307 " high, after %hhu iteration%s)", nlo + nhi, nsel, nlo, nhi,
3308 i, i == 1 ?
"" :
"s");
3310 cpl_table_delete(spectrum);
3313 cpl_array_delete(asel);
#define MUSE_PIXTABLE_XPOS
muse_image * muse_resampling_collapse_pixgrid(muse_pixtable *aPixtable, muse_pixgrid *aGrid, muse_datacube *aCube, const muse_table *aFilter, muse_resampling_params *aParams)
Integrate a pixel table / pixel grid along the wavelength direction.
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.
double muse_pfits_get_cd(const cpl_propertylist *aHeaders, unsigned int aAxisI, unsigned int aAxisJ)
find out the WCS coordinate at the reference point
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 pixel grid, linking the grid points to entries (=rows)...
#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.
double muse_pfits_get_rhum(const cpl_propertylist *aHeaders)
find out the relavtive humidity (in %)
double muse_pfits_get_crval(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS coordinate at the reference point
A structure containing a spatial two-axis WCS.
cpl_table * muse_resampling_spectrum_iterate(muse_pixtable *aPixtable, double aBinwidth, float aLo, float aHi, unsigned char aIter)
Iteratively resample selected pixels of a pixel table into spectrum.
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
cpl_size muse_pixtable_get_nrow(const muse_pixtable *aPixtable)
get the number of rows within the pixel table
#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.
cpl_error_code muse_utils_filter_copy_properties(cpl_propertylist *aHeader, const muse_table *aFilter, double aFraction)
Add/propagate filter properties to header of collapsed image.
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.
muse_datacube * muse_resampling_cube(muse_pixtable *aPixtable, muse_resampling_params *aParams, muse_pixgrid **aGrid)
Resample a pixel table onto a regular grid structure representing a FITS NAXIS=3 datacube.
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.
const char * muse_pfits_get_dateobs(const cpl_propertylist *aHeaders)
find out the date of observations
void muse_datacube_delete(muse_datacube *aCube)
Deallocate memory associated to a muse_datacube object.
cpl_propertylist * hgroup
the group FITS header
#define MUSE_HDR_PT_FFCORR
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
const char * muse_pfits_get_cunit(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS axis unit string
cpl_error_code muse_pixtable_reset_dq(muse_pixtable *aPixtable, unsigned int aDQ)
Reset a given bad pixel status (DQ flag) for all pixels in the table.
cpl_error_code muse_cpltable_check(const cpl_table *aTable, const muse_cpltable_def *aDef)
Check whether the table contains the fields of the definition.
#define MUSE_PIXTABLE_DATA
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_pfits_get_crpix(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS reference point
void muse_pixgrid_delete(muse_pixgrid *aGrid)
Delete a pixel grid and remove its memory.
unsigned int muse_pixtable_origin_get_x(uint32_t aOrigin, muse_pixtable *aPixtable, cpl_size aRow)
Get the horizontal coordinate from the encoded 32bit origin number.
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.
cpl_array * muse_cplarray_new_from_delimited_string(const char *aString, const char *aDelim)
Convert a delimited string into an array of strings.
#define MUSE_PIXTABLE_WEIGHT
Structure definition of MUSE pixel table.
static const cpl_size * muse_pixgrid_get_rows(muse_pixgrid *aGrid, cpl_size aIndex)
Return a pointer to the rows stored in one pixel.
#define MUSE_WCS_KEYS
regular expression for WCS properties
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_table * table
The table.
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.
Structure to store a table together with a property list.
#define MUSE_PIXTABLE_ORIGIN
double muse_utils_filter_fraction(const muse_table *aFilter, double aLambda1, double aLambda2)
Compute fraction of filter area covered by the data range.
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_size muse_cpltable_find_sorted(const cpl_table *aTable, const char *aColumn, double aValue)
Find a row in a table.
cpl_table * dtable
the table containing the actual Euro3D data
const char * muse_pfits_get_ctype(const cpl_propertylist *aHeaders, unsigned int aAxis)
find out the WCS axis type string
cpl_image * muse_cplimage_concat_x(const cpl_image *aImage1, const cpl_image *aImage2)
Concatenate two images in x direction.
cpl_error_code muse_resampling_params_set_pixfrac(muse_resampling_params *aParams, const char *aString)
Set resampling pixfrac given a string that can contain up to three floating-point values...
cpl_propertylist * hdata
the data FITS header
unsigned short muse_pixtable_origin_get_ifu(uint32_t aOrigin)
Get the IFU number from the encoded 32bit origin number.
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.
#define MUSE_PIXTABLE_STAT
double muse_pfits_get_mjdobs(const cpl_propertylist *aHeaders)
find out the Julian Date of the observation
const muse_cpltable_def muse_pixtable_def[]
MUSE pixel table definition.
#define MUSE_PIXTABLE_LAMBDA
static cpl_size muse_pixgrid_get_index(muse_pixgrid *aGrid, cpl_size aX, cpl_size aY, cpl_size aZ, cpl_boolean aAllowOutside)
Get the grid index determined from all three coordinates.
cpl_error_code muse_quadrants_coords_to_raw(cpl_propertylist *aHeader, int *aX, int *aY)
Convert coordinates of a trimmed image to raw-image coordinates.
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.
unsigned int muse_pixtable_origin_get_y(uint32_t aOrigin)
Get the vertical coordinate from the encoded 32bit origin number.
#define MUSE_PIXTABLE_YPOS
#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 *aGrid, 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.