34 #include "muse_flux.h"
35 #include "muse_instrument.h"
37 #include "muse_astro.h"
38 #include "muse_cplwrappers.h"
39 #include "muse_pfits.h"
40 #include "muse_quality.h"
41 #include "muse_resampling.h"
42 #include "muse_utils.h"
90 aFluxObj->cube = NULL;
92 aFluxObj->intimage = NULL;
93 cpl_table_delete(aFluxObj->sensitivity);
94 aFluxObj->sensitivity = NULL;
95 cpl_table_delete(aFluxObj->response);
96 aFluxObj->response = NULL;
97 cpl_table_delete(aFluxObj->telluric);
98 aFluxObj->telluric = NULL;
99 cpl_table_delete(aFluxObj->tellbands);
100 aFluxObj->tellbands = NULL;
116 cpl_ensure(aTable, CPL_ERROR_NULL_INPUT, 0.);
117 cpl_table_unselect_all(aTable);
118 cpl_table_or_selected_double(aTable,
"lambda", CPL_NOT_LESS_THAN,
119 kMuseNominalLambdaMin);
120 cpl_table_and_selected_double(aTable,
"lambda", CPL_NOT_GREATER_THAN,
121 kMuseNominalLambdaMax);
122 cpl_size nsel = cpl_table_count_selected(aTable);
123 cpl_array *asel = cpl_table_where_selected(aTable);
124 cpl_size *sel = cpl_array_get_data_cplsize(asel);
125 double lmin = cpl_table_get_double(aTable,
"lambda", sel[0], NULL),
126 lmax = cpl_table_get_double(aTable,
"lambda", sel[nsel - 1], NULL);
127 cpl_array_delete(asel);
128 return (lmax - lmin) / nsel;
164 cpl_ensure_code(aTable, CPL_ERROR_NULL_INPUT);
166 const char *flam1 =
"erg/s/cm**2/Angstrom",
167 *flam2 =
"erg/s/cm^2/Angstrom";
169 cpl_error_code rc = CPL_ERROR_NONE;
170 cpl_errorstate prestate = cpl_errorstate_get();
172 if (cpl_table_has_column(aTable,
"lambda") &&
173 cpl_table_has_column(aTable,
"flux") &&
174 cpl_table_get_column_unit(aTable,
"lambda") &&
175 cpl_table_get_column_unit(aTable,
"flux") &&
176 !strncmp(cpl_table_get_column_unit(aTable,
"lambda"),
"Angstrom", 9) &&
177 (!strncmp(cpl_table_get_column_unit(aTable,
"flux"), flam1, strlen(flam1)) ||
178 !strncmp(cpl_table_get_column_unit(aTable,
"flux"), flam2, strlen(flam2)))) {
181 if (cpl_table_get_column_type(aTable,
"lambda") != CPL_TYPE_DOUBLE) {
182 cpl_msg_debug(__func__,
"Casting lambda column to double");
183 cpl_table_cast_column(aTable,
"lambda", NULL, CPL_TYPE_DOUBLE);
185 if (cpl_table_get_column_type(aTable,
"flux") != CPL_TYPE_DOUBLE) {
186 cpl_msg_debug(__func__,
"Casting flux column to double");
187 cpl_table_cast_column(aTable,
"flux", NULL, CPL_TYPE_DOUBLE);
190 if (cpl_table_has_column(aTable,
"fluxerr")) {
191 if (cpl_table_get_column_type(aTable,
"fluxerr") != CPL_TYPE_DOUBLE) {
192 cpl_msg_debug(__func__,
"Casting fluxerr column to double");
193 cpl_table_cast_column(aTable,
"fluxerr", NULL, CPL_TYPE_DOUBLE);
195 const char *unit = cpl_table_get_column_unit(aTable,
"fluxerr");
196 if (!unit || (strncmp(unit, flam1, strlen(flam1)) &&
197 strncmp(unit, flam2, strlen(flam2)))) {
198 cpl_msg_debug(__func__,
"Erasing fluxerr column because of unexpected "
200 cpl_table_erase_column(aTable,
"fluxerr");
203 cpl_msg_info(__func__,
"Found MUSE format, average sampling %.3f Angstrom/bin"
205 }
else if (cpl_table_has_column(aTable,
"WAVELENGTH") &&
206 cpl_table_has_column(aTable,
"FLUX") &&
207 cpl_table_get_column_unit(aTable,
"WAVELENGTH") &&
208 cpl_table_get_column_unit(aTable,
"FLUX") &&
209 !strncmp(cpl_table_get_column_unit(aTable,
"WAVELENGTH"),
"ANGSTROMS", 10) &&
210 !strncmp(cpl_table_get_column_unit(aTable,
"FLUX"),
"FLAM", 5)) {
212 printf(
"input HST CALSPEC table:\n");
213 cpl_table_dump_structure(aTable, stdout);
214 cpl_table_dump(aTable, cpl_table_get_nrow(aTable)/2, 3, stdout);
218 cpl_table_cast_column(aTable,
"WAVELENGTH",
"lambda", CPL_TYPE_DOUBLE);
219 cpl_table_cast_column(aTable,
"FLUX",
"flux", CPL_TYPE_DOUBLE);
220 cpl_table_erase_column(aTable,
"WAVELENGTH");
221 cpl_table_erase_column(aTable,
"FLUX");
222 cpl_table_set_column_unit(aTable,
"lambda",
"Angstrom");
223 cpl_table_set_column_unit(aTable,
"flux", flam1);
226 if (cpl_table_has_column(aTable,
"STATERROR") &&
227 cpl_table_has_column(aTable,
"SYSERROR") &&
228 cpl_table_get_column_unit(aTable,
"STATERROR") &&
229 cpl_table_get_column_unit(aTable,
"SYSERROR") &&
230 !strncmp(cpl_table_get_column_unit(aTable,
"STATERROR"),
"FLAM", 5) &&
231 !strncmp(cpl_table_get_column_unit(aTable,
"SYSERROR"),
"FLAM", 5)) {
234 cpl_table_cast_column(aTable,
"STATERROR",
"fluxerr", CPL_TYPE_DOUBLE);
235 cpl_table_erase_column(aTable,
"STATERROR");
236 cpl_table_cast_column(aTable,
"SYSERROR", NULL, CPL_TYPE_DOUBLE);
237 cpl_table_power_column(aTable,
"fluxerr", 2);
238 cpl_table_power_column(aTable,
"SYSERROR", 2);
239 cpl_table_add_columns(aTable,
"fluxerr",
"SYSERROR");
240 cpl_table_erase_column(aTable,
"SYSERROR");
241 cpl_table_power_column(aTable,
"fluxerr", 0.5);
242 cpl_table_set_column_unit(aTable,
"fluxerr", flam1);
249 if (cpl_table_has_column(aTable,
"FWHM")) {
250 cpl_table_erase_column(aTable,
"FWHM");
252 if (cpl_table_has_column(aTable,
"DATAQUAL")) {
253 cpl_table_erase_column(aTable,
"DATAQUAL");
255 if (cpl_table_has_column(aTable,
"TOTEXP")) {
256 cpl_table_erase_column(aTable,
"TOTEXP");
259 printf(
"converted HST CALSPEC table:\n");
260 cpl_table_dump_structure(aTable, stdout);
261 cpl_table_dump(aTable, cpl_table_get_nrow(aTable)/2, 3, stdout);
264 cpl_msg_info(__func__,
"Found HST CALSPEC format on input, converted to "
265 "MUSE format; average sampling %.3f Angstrom/bin over MUSE "
268 cpl_msg_error(__func__,
"Unknown format found!");
270 cpl_table_dump_structure(aTable, stdout);
272 rc = CPL_ERROR_INCOMPATIBLE_INPUT;
276 if (!cpl_errorstate_is_equal(prestate)) {
277 rc = cpl_error_get_code();
329 cpl_ensure(aResponse, CPL_ERROR_NULL_INPUT, rv);
330 int size = cpl_table_get_nrow(aResponse);
331 cpl_ensure(size > 0, cpl_error_get_code(), rv);
334 const double *lbda = cpl_table_get_data_double_const(aResponse,
"lambda"),
335 *resp = NULL, *rerr = NULL;
338 resp = cpl_table_get_data_double_const(aResponse,
"throughput");
341 resp = cpl_table_get_data_double_const(aResponse,
"response");
343 rerr = cpl_table_get_data_double_const(aResponse,
"resperr");
347 resp = cpl_table_get_data_double_const(aResponse,
"flux");
349 rerr = cpl_table_get_data_double_const(aResponse,
"fluxerr");
353 resp = cpl_table_get_data_double_const(aResponse,
"extinction");
356 resp = cpl_table_get_data_double_const(aResponse,
"ftelluric");
358 rerr = cpl_table_get_data_double_const(aResponse,
"ftellerr");
362 cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
365 cpl_ensure(lbda && resp, cpl_error_get_code(), rv);
367 cpl_ensure(rerr, cpl_error_get_code(), rv);
371 if (aLambda < lbda[0]) {
374 if (aLambda > lbda[size-1]) {
379 double response = rv, resperror = 0.;
380 int l = 0, r = size - 1,
383 if (aLambda >= lbda[m] && aLambda <= lbda[m+1]) {
385 double lquot = (aLambda - lbda[m]) / (lbda[m+1] - lbda[m]);
386 response = resp[m] + (resp[m+1] - resp[m]) * lquot;
391 resperror = sqrt(pow(rerr[m] * (1 - lquot), 2.)
392 + pow(rerr[m+1] * lquot, 2.));
395 cpl_msg_debug(__func__,
"Found at m=%d (%f: %f+/-%f) / "
396 "m+1=%d (%f: %f+/-%f) -> %f: %f+/-%f",
397 m, lbda[m], resp[m], rerr ? rerr[m] : 0.,
398 m+1, lbda[m+1], resp[m+1], rerr ? rerr[m+1] : 0.,
399 aLambda, response, resperror);
404 if (aLambda < lbda[m]) {
407 if (aLambda > lbda[m]) {
414 cpl_msg_debug(__func__,
"Response %g+/-%g at lambda=%fA", response, resperror,
417 if (aError && rerr) {
440 muse_flux_image_sky(cpl_image *aImage,
double aX,
double aY,
double aHalfSize,
441 unsigned int aDSky,
float *aError)
447 int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
448 y1 = aY - aHalfSize, y2 = aY + aHalfSize;
449 unsigned char nskyarea = 0;
450 double skylevel = 0., skyerror = 0.;
452 cpl_errorstate state = cpl_errorstate_get();
453 cpl_stats_mode mode = CPL_STATS_MEDIAN | CPL_STATS_MEDIAN_DEV;
454 cpl_stats *s = cpl_stats_new_from_image_window(aImage, mode,
455 x1 - aDSky, y1, x1 - 1, y2);
460 skylevel += cpl_stats_get_median(s);
461 skyerror += pow(cpl_stats_get_median_dev(s), 2);
465 s = cpl_stats_new_from_image_window(aImage,mode,
466 x2 + 1, y1, x2 + aDSky, y2);
469 skylevel += cpl_stats_get_median(s);
470 skyerror += pow(cpl_stats_get_median_dev(s), 2);
474 s = cpl_stats_new_from_image_window(aImage,mode,
475 x1, y1 - aDSky, x2, y1 - 1);
478 skylevel += cpl_stats_get_median(s);
479 skyerror += pow(cpl_stats_get_median_dev(s), 2);
483 s = cpl_stats_new_from_image_window(aImage,mode,
484 x1, y2 + 1, x2, y2 + aDSky);
487 skylevel += cpl_stats_get_median(s);
488 skyerror += pow(cpl_stats_get_median_dev(s), 2);
494 skylevel /= nskyarea;
495 skyerror = sqrt(skyerror) / nskyarea;
496 if (!cpl_errorstate_is_equal(state)) {
497 cpl_errorstate_set(state);
500 cpl_msg_debug(__func__,
"skylevel = %f +/- %f (%u sky areas)",
501 skylevel, skyerror, nskyarea);
528 muse_flux_image_gaussian(cpl_image *aImage, cpl_image *aImErr,
double aX,
529 double aY,
double aHalfSize,
unsigned int aDSky,
530 unsigned int aMaxBad,
float *aFErr)
534 UNUSED_ARGUMENT(aMaxBad);
540 cpl_array *params = cpl_array_new(7, CPL_TYPE_DOUBLE),
541 *parerr = cpl_array_new(7, CPL_TYPE_DOUBLE);
545 cpl_errorstate state = cpl_errorstate_get();
546 double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky, NULL);
547 if (!cpl_errorstate_is_equal(state)) {
550 cpl_errorstate_set(state);
552 cpl_array_set_double(params, 0, skylevel);
553 cpl_array_set_double(params, 3, aX);
554 cpl_array_set_double(params, 4, aY);
555 double rms = 0, chisq = 0;
557 int nx = cpl_image_get_size_x(aImage),
558 ny = cpl_image_get_size_y(aImage),
559 xsize = fmin(aHalfSize, fmin(aX - 1., nx - aX)) * 2,
560 ysize = fmin(aHalfSize, fmin(aY - 1., ny - aY)) * 2;
561 cpl_error_code rc = cpl_fit_image_gaussian(aImage, aImErr, aX, aY, xsize, ysize,
562 params, parerr, NULL, &rms, &chisq,
563 NULL, NULL, NULL, NULL, NULL);
564 if (rc != CPL_ERROR_NONE) {
565 if (rc != CPL_ERROR_ILLEGAL_INPUT) {
566 cpl_msg_debug(__func__,
"rc = %d: %s", rc, cpl_error_get_message());
568 cpl_array_delete(params);
569 cpl_array_delete(parerr);
572 double flux = cpl_array_get_double(params, 1, NULL),
573 ferr = cpl_array_get_double(parerr, 1, NULL);
575 double fwhmx = cpl_array_get_double(params, 5, NULL) * CPL_MATH_FWHM_SIG,
576 fwhmy = cpl_array_get_double(params, 6, NULL) * CPL_MATH_FWHM_SIG;
577 cpl_msg_debug(__func__,
"%.3f,%.3f: %g+/-%g (bg: %g, FWHM: %.3f,%.3f, %g, %g)",
578 cpl_array_get_double(params, 3, NULL), cpl_array_get_double(params, 4, NULL),
579 flux, ferr, cpl_array_get_double(params, 0, NULL), fwhmx, fwhmy,
583 cpl_msg_debug(__func__,
"skylevel = %f", cpl_array_get_double(params, 0, NULL));
584 cpl_msg_debug(__func__,
"measured flux %f +/- %f", flux, ferr);
586 cpl_array_delete(params);
587 cpl_array_delete(parerr);
615 muse_flux_image_moffat(cpl_image *aImage, cpl_image *aImErr,
double aX,
616 double aY,
double aHalfSize,
unsigned int aDSky,
617 unsigned int aMaxBad,
float *aFErr)
623 int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
624 y1 = aY - aHalfSize, y2 = aY + aHalfSize,
625 nx = cpl_image_get_size_x(aImage),
626 ny = cpl_image_get_size_y(aImage);
639 int npoints = (x2 - x1 + 1) * (y2 - y1 + 1);
641 cpl_matrix *pos = cpl_matrix_new(npoints, 2);
642 cpl_vector *values = cpl_vector_new(npoints),
643 *errors = cpl_vector_new(npoints);
644 float *derr = cpl_image_get_data_float(aImErr);
646 for (i = x1; i <= x2; i++) {
648 for (j = y1; j <= y2; j++) {
650 double data = cpl_image_get(aImage, i, j, &err);
654 cpl_matrix_set(pos, idx, 0, i);
655 cpl_matrix_set(pos, idx, 1, j);
656 cpl_vector_set(values, idx, data);
657 cpl_vector_set(errors, idx, derr[(i-1) + (j-1)*nx]);
663 if (idx < 16 || (npoints - idx) > (
int)aMaxBad) {
664 cpl_matrix_delete(pos);
665 cpl_vector_delete(values);
666 cpl_vector_delete(errors);
669 cpl_matrix_set_size(pos, idx, 2);
670 cpl_vector_set_size(values, idx);
671 cpl_vector_set_size(errors, idx);
673 cpl_array *params = cpl_array_new(8, CPL_TYPE_DOUBLE),
674 *parerr = cpl_array_new(8, CPL_TYPE_DOUBLE);
676 cpl_errorstate state = cpl_errorstate_get();
677 double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky, NULL);
678 if (!cpl_errorstate_is_equal(state)) {
681 cpl_errorstate_set(state);
683 cpl_array_set_double(params, 0, skylevel);
684 cpl_array_set_double(params, 2, aX);
685 cpl_array_set_double(params, 3, aY);
686 double rms = 0, chisq = 0;
688 params, parerr, NULL,
690 cpl_matrix_delete(pos);
691 cpl_vector_delete(values);
692 cpl_vector_delete(errors);
693 if (rc != CPL_ERROR_NONE) {
694 if (rc != CPL_ERROR_ILLEGAL_INPUT) {
695 cpl_msg_debug(__func__,
"rc = %d: %s", rc, cpl_error_get_message());
697 cpl_array_delete(params);
698 cpl_array_delete(parerr);
702 double flux = cpl_array_get_double(params, 1, NULL);
704 *aFErr = cpl_array_get_double(parerr, 1, NULL);
707 cpl_msg_debug(__func__,
"skylevel = %f", cpl_array_get_double(params, 0, NULL));
708 cpl_msg_debug(__func__,
"measured flux %f +/- %f", flux, cpl_array_get_double(parerr, 1, NULL));
710 cpl_array_delete(params);
711 cpl_array_delete(parerr);
736 muse_flux_image_square(cpl_image *aImage, cpl_image *aImErr,
double aX,
737 double aY,
double aHalfSize,
unsigned int aDSky,
738 unsigned int aMaxBad,
float *aFErr)
743 int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
744 y1 = aY - aHalfSize, y2 = aY + aHalfSize,
745 nx = cpl_image_get_size_x(aImage),
746 ny = cpl_image_get_size_y(aImage);
760 cpl_errorstate state = cpl_errorstate_get();
761 double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky,
763 if (!cpl_errorstate_is_equal(state)) {
766 cpl_errorstate_set(state);
767 cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
773 cpl_image *region = cpl_image_extract(aImage, x1, y1, x2, y2);
775 cpl_msg_debug(__func__,
"region [%d:%d,%d:%d] %"CPL_SIZE_FORMAT
" bad pixels",
776 x1, y1, x2, y2, cpl_image_count_rejected(region));
778 if (cpl_image_count_rejected(region) > aMaxBad) {
779 cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
780 cpl_image_delete(region);
783 cpl_image *regerr = cpl_image_extract(aImErr, x1, y1, x2, y2);
784 if (cpl_image_count_rejected(region) > 0) {
785 cpl_detector_interpolate_rejected(region);
786 cpl_detector_interpolate_rejected(regerr);
791 int npoints = (x2 - x1 + 1) * (y2 - y1 + 1),
795 nsky = 2 * aDSky * (x2 - x1 + y2 - y1 + 2);
796 double flux = cpl_image_get_flux(region) - skylevel * npoints,
797 ferr = sqrt(cpl_image_get_sqflux(regerr)
798 + npoints * skyerror*skyerror * (1. + (
double)npoints / nsky));
800 cpl_msg_debug(__func__,
"measured flux %f +/- %f (%d object pixels, %d pixels"
801 " with %f with sky %f)", flux, ferr, npoints, nsky,
802 cpl_image_get_flux(region), cpl_image_get_sqflux(regerr));
804 cpl_image_delete(region);
805 cpl_image_delete(regerr);
834 muse_flux_image_circle(cpl_image *aImage, cpl_image *aImErr,
double aX,
835 double aY,
double aAper,
double aAnnu,
double aDAnnu,
836 unsigned int aMaxBad,
float *aFErr)
841 double rmax = ceil(fmax(aAper, aAnnu + aDAnnu));
842 int x1 = aX - rmax, x2 = aX + rmax,
843 y1 = aY - rmax, y2 = aY + rmax,
844 nx = cpl_image_get_size_x(aImage),
845 ny = cpl_image_get_size_y(aImage);
860 cpl_vector *vbg = cpl_vector_new((x2 - x1 + 1) * (y2 - y1 + 1)),
861 *vbe = cpl_vector_new((x2 - x1 + 1) * (y2 - y1 + 1));
862 unsigned int nbad = 0, nbg = 0;
864 for (i = x1; i <= x2; i++) {
866 for (j = y1; j <= y2; j++) {
867 double r = sqrt(pow(aX - i, 2) + pow(aY - j, 2));
869 nbad += cpl_image_is_rejected(aImage, i, j) == 1;
871 if (r < aAnnu || r > aAnnu + aDAnnu) {
875 double value = cpl_image_get(aImage, i, j, &err);
879 cpl_vector_set(vbg, nbg, value);
880 cpl_vector_set(vbe, nbg, cpl_image_get(aImErr, i, j, &err));
885 cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
886 cpl_vector_delete(vbg);
887 cpl_vector_delete(vbe);
890 cpl_vector_set_size(vbg, nbg);
891 cpl_vector_set_size(vbe, nbg);
892 cpl_matrix *pos = cpl_matrix_new(1, nbg);
897 unsigned int nrej = nbg - cpl_vector_get_size(vbg);
899 nbg = cpl_vector_get_size(vbg);
901 double smean = cpl_polynomial_get_coeff(fit, &pows),
903 cpl_polynomial_delete(fit);
904 cpl_matrix_delete(pos);
905 cpl_vector_delete(vbg);
906 cpl_vector_delete(vbe);
908 cpl_msg_debug(__func__,
"sky: %d pixels (%d rejected), %f +/- %f; found %d "
909 "bad pixels inside aperture", nbg, nrej, smean, sstdev, nbad);
911 if (nbad > aMaxBad) {
912 cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
918 cpl_detector_interpolate_rejected(aImage);
919 cpl_detector_interpolate_rejected(aImErr);
925 unsigned int nobj = 0;
926 for (i = x1; i <= x2; i++) {
928 for (j = y1; j <= y2; j++) {
929 double r = sqrt(pow(aX - i, 2) + pow(aY - j, 2));
934 double value = cpl_image_get(aImage, i, j, &err),
935 error = cpl_image_get(aImErr, i, j, &err);
941 flux -= smean * nobj;
946 ferr = sqrt(ferr + nobj * sstdev*sstdev * (1. + (
double)nobj / nbg));
948 cpl_msg_debug(__func__,
"flux: %d pixels (%d interpolated), %f +/- %f",
949 nobj, nbad, flux, ferr);
1009 cpl_ensure_code(aPixtable && aFluxObj, CPL_ERROR_NULL_INPUT);
1012 cpl_msg_info(__func__,
"Gaussian profile fits for flux integration");
1015 cpl_msg_info(__func__,
"Moffat profile fits for flux integration");
1018 cpl_msg_info(__func__,
"Circular flux integration");
1021 cpl_msg_info(__func__,
"Simple square window flux integration");
1024 return cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
1027 if (getenv(
"MUSE_DEBUG_FLUX") && atoi(getenv(
"MUSE_DEBUG_FLUX")) > 2) {
1028 const char *fn =
"flux__pixtable.fits";
1029 cpl_msg_info(__func__,
"Saving pixel table as \"%s\"", fn);
1039 params->dlambda = kMuseSpectralSamplingA;
1044 aFluxObj->cube = cube;
1047 if (getenv(
"MUSE_DEBUG_FLUX") && atoi(getenv(
"MUSE_DEBUG_FLUX")) >= 2) {
1048 const char *fn =
"flux__cube.fits";
1049 cpl_msg_info(__func__,
"Saving cube as \"%s\"", fn);
1052 int nplane = cpl_imagelist_get_size(cube->
data) / 2;
1053 cpl_image *cim = cpl_imagelist_get(cube->
data, nplane);
1055 double dsigmas[] = { 50., 30., 10., 8., 6., 5. };
1056 cpl_vector *vsigmas = cpl_vector_wrap(
sizeof(dsigmas) /
sizeof(
double),
1058 cpl_size isigma = -1;
1059 cpl_apertures *apertures = cpl_apertures_extract(cim, vsigmas, &isigma);
1060 int napertures = apertures ? cpl_apertures_get_size(apertures) : 0;
1061 if (napertures < 1) {
1062 cpl_vector_unwrap(vsigmas);
1063 cpl_apertures_delete(apertures);
1064 return cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
1066 cpl_msg_debug(__func__,
"The %.1f sigma threshold was used to find %d source%s",
1067 cpl_vector_get(vsigmas, isigma), napertures, napertures == 1 ?
"" :
"s");
1068 cpl_vector_unwrap(vsigmas);
1070 cpl_apertures_dump(apertures, stdout);
1075 int nlambda = cpl_imagelist_get_size(cube->
data);
1077 intimage->
data = cpl_image_new(nlambda, napertures, CPL_TYPE_FLOAT);
1078 intimage->
stat = cpl_image_new(nlambda, napertures, CPL_TYPE_FLOAT);
1080 intimage->
header = cpl_propertylist_new();
1081 cpl_propertylist_update_double(intimage->
header,
"CRVAL1",
1082 cpl_propertylist_get_double(cube->
header,
1084 cpl_propertylist_update_double(intimage->
header,
"CRPIX1",
1085 cpl_propertylist_get_double(cube->
header,
1087 cpl_propertylist_update_double(intimage->
header,
"CD1_1",
1088 cpl_propertylist_get_double(cube->
header,
1090 cpl_propertylist_update_string(intimage->
header,
"CTYPE1",
1091 cpl_propertylist_get_string(cube->
header,
1093 cpl_propertylist_update_string(intimage->
header,
"CUNIT1",
1094 cpl_propertylist_get_string(cube->
header,
1097 cpl_propertylist_update_double(intimage->
header,
"CRVAL2", 1.);
1098 cpl_propertylist_update_double(intimage->
header,
"CRPIX2", 1.);
1099 cpl_propertylist_update_double(intimage->
header,
"CD2_2", 1.);
1100 cpl_propertylist_update_string(intimage->
header,
"CTYPE2",
"PIXEL");
1101 cpl_propertylist_update_string(intimage->
header,
"CUNIT2",
"pixel");
1102 cpl_propertylist_update_double(intimage->
header,
"CD1_2", 0.);
1103 cpl_propertylist_update_double(intimage->
header,
"CD2_1", 0.);
1105 cpl_propertylist_append_string(intimage->
header,
"BUNIT",
1106 cpl_propertylist_get_string(cube->
header,
1108 cpl_propertylist_update_double(intimage->
header,
"EXPTIME",
1110 cpl_propertylist_update_string(intimage->
header,
"ESO INS MODE",
1111 cpl_propertylist_get_string(cube->
header,
1115 cpl_errorstate ps = cpl_errorstate_get();
1119 fwhm /= (kMuseSpaxelSizeX_WFM + kMuseSpaxelSizeY_WFM) / 2.;
1121 fwhm /= (kMuseSpaxelSizeX_NFM + kMuseSpaxelSizeY_NFM) / 2.;
1123 if (!cpl_errorstate_is_equal(ps)) {
1124 double xc = cpl_apertures_get_centroid_x(apertures, 1),
1125 yc = cpl_apertures_get_centroid_y(apertures, 1),
1127 cpl_image_get_fwhm(cim, lround(xc), lround(yc), &xfwhm, &yfwhm);
1128 if (xfwhm > 0. && yfwhm > 0.) {
1129 fwhm = (xfwhm + yfwhm) / 2.;
1130 }
else if (xfwhm > 0.) {
1132 }
else if (yfwhm > 0.) {
1137 cpl_errorstate_set(ps);
1138 cpl_msg_warning(__func__,
"Using roughly estimated FWHM (%.3f pix) instead "
1139 "of DIMM seeing", fwhm);
1141 cpl_msg_debug(__func__,
"Using DIMM seeing of %.3f pix", fwhm);
1144 float *data = cpl_image_get_data_float(intimage->
data),
1145 *stat = cpl_image_get_data_float(intimage->
stat);
1147 int l, ngood = 0, nillegal = 0, nbadbg = 0;
1148 #pragma omp parallel for default(none) \
1149 shared(apertures, aProfile, cube, data, fwhm, napertures, nbadbg, \
1150 ngood, nillegal, nlambda, stat)
1151 for (l = 0; l < nlambda; l++) {
1152 cpl_image *plane = cpl_imagelist_get(cube->
data, l),
1153 *pldq = cpl_imagelist_get(cube->
dq, l),
1154 *plerr = cpl_image_duplicate(cpl_imagelist_get(cube->
stat, l));
1156 cpl_stats *stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1157 cpl_msg_debug(__func__,
"lambda = %d/%f %s", l + 1,
1158 (l + 1 - cpl_propertylist_get_double(cube->
header,
"CRPIX3"))
1159 * cpl_propertylist_get_double(cube->
header,
"CD3_3")
1160 + cpl_propertylist_get_double(cube->
header,
"CRVAL3"),
1161 cpl_propertylist_get_string(cube->
header,
"CUNIT3"));
1162 cpl_msg_debug(__func__,
"variance: %g...%g...%g", cpl_stats_get_min(stats),
1163 cpl_stats_get_mean(stats), cpl_stats_get_max(stats));
1164 cpl_stats_delete(stats);
1169 stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1170 cpl_msg_debug(__func__,
"cut variance: %g...%g...%g (%"CPL_SIZE_FORMAT
" bad"
1171 " pixel)", cpl_stats_get_min(stats), cpl_stats_get_mean(stats),
1172 cpl_stats_get_max(stats), cpl_image_count_rejected(plane));
1173 cpl_stats_delete(stats);
1176 cpl_image_power(plerr, 0.5);
1178 stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1179 cpl_msg_debug(__func__,
"errors: %g...%g...%g", cpl_stats_get_min(stats),
1180 cpl_stats_get_mean(stats), cpl_stats_get_max(stats));
1181 cpl_stats_delete(stats);
1183 cpl_errorstate state = cpl_errorstate_get();
1185 for (n = 1; n <= napertures; n++) {
1187 double xc = cpl_apertures_get_centroid_x(apertures, n),
1188 yc = cpl_apertures_get_centroid_y(apertures, n),
1189 size = sqrt(cpl_apertures_get_npix(apertures, n)),
1191 cpl_errorstate prestate = cpl_errorstate_get();
1192 cpl_image_get_fwhm(plane, lround(xc), lround(yc), &xfwhm, &yfwhm);
1193 if (xfwhm < 0 || yfwhm < 0) {
1194 data[l + (n-1) * nlambda] = 0.;
1195 stat[l + (n-1) * nlambda] = FLT_MAX;
1196 cpl_errorstate_set(prestate);
1200 double halfsize = 1.5 * (xfwhm + yfwhm);
1201 if (halfsize < size / 2) {
1202 halfsize = size / 2;
1205 cpl_msg_debug(__func__,
"%.2f,%.2f FWHM %.2f %.2f size %.2f --> %.2f",
1206 xc, yc, xfwhm, yfwhm, size, halfsize * 2.);
1211 data[l + (n-1) * nlambda] = muse_flux_image_gaussian(plane, plerr, xc, yc,
1213 &stat[l + (n-1) * nlambda]);
1216 data[l + (n-1) * nlambda] = muse_flux_image_moffat(plane, plerr, xc, yc,
1218 &stat[l + (n-1) * nlambda]);
1223 double maxfwhm = fmax((xfwhm + yfwhm) / 2., fwhm),
1224 radius = halfsize < maxfwhm ? maxfwhm : halfsize,
1225 rannu = radius * 5. / 4.;
1226 data[l + (n-1) * nlambda] = muse_flux_image_circle(plane, plerr, xc, yc,
1227 radius, rannu, 10, 10,
1228 &stat[l + (n-1) * nlambda]);
1232 data[l + (n-1) * nlambda] = muse_flux_image_square(plane, plerr, xc, yc,
1234 &stat[l + (n-1) * nlambda]);
1236 if (data[l + (n-1) * nlambda] < 0 || !isfinite(data[l + (n-1) * nlambda])) {
1237 data[l + (n-1) * nlambda] = 0.;
1238 stat[l + (n-1) * nlambda] = FLT_MAX;
1244 if (!cpl_errorstate_is_equal(state)) {
1245 if (cpl_error_get_code() == CPL_ERROR_ILLEGAL_INPUT) {
1246 cpl_errorstate_set(state);
1249 }
else if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
1250 cpl_errorstate_set(state);
1259 cpl_image_delete(plerr);
1261 cpl_apertures_delete(apertures);
1262 if (getenv(
"MUSE_DEBUG_FLUX") && atoi(getenv(
"MUSE_DEBUG_FLUX")) >= 2) {
1263 const char *fn =
"flux__intimage.fits";
1264 cpl_msg_info(__func__,
"Saving flux image as \"%s\"", fn);
1265 intimage->
dq = cpl_image_new(nlambda, napertures, CPL_TYPE_INT);
1268 aFluxObj->intimage = intimage;
1270 if (nillegal > 0 || nbadbg > 0) {
1271 cpl_msg_warning(__func__,
"Successful fits in %d wavelength planes, but "
1272 "encountered %d \"Illegal input\" errors and %d bad "
1273 "background determinations", ngood, nillegal, nbadbg);
1275 cpl_msg_info(__func__,
"Successful fits in %d wavelength planes", ngood);
1277 return CPL_ERROR_NONE;
1287 static const double kTelluricBands[][4] = {
1288 { 6273., 6320., 6213., 6380. },
1289 { 6864., 6967., 6750., 7130. },
1290 { 7164., 7325., 7070., 7580. },
1291 { 7590., 7700., 7470., 7830. },
1292 { 8131., 8345., 7900., 8600. },
1293 { 8952., 9028., 8850., 9082. },
1294 { 9274., 9770., 9080., 9263. },
1295 { -1., -1., -1., -1. }
1304 {
"lmin", CPL_TYPE_DOUBLE,
"Angstrom",
"%8.3f",
1305 "lower limit of the telluric region", CPL_TRUE },
1306 {
"lmax", CPL_TYPE_DOUBLE,
"Angstrom",
"%8.3f",
1307 "upper limit of the telluric region", CPL_TRUE },
1308 {
"bgmin", CPL_TYPE_DOUBLE,
"Angstrom",
"%8.3f",
1309 "lower limit of the background region", CPL_TRUE },
1310 {
"bgmax", CPL_TYPE_DOUBLE,
"Angstrom",
"%8.3f",
1311 "upper limit of the background region", CPL_TRUE },
1312 { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
1325 const cpl_table *aTellBands)
1328 cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
1333 == CPL_ERROR_NONE) {
1334 cpl_msg_debug(__func__,
"using given table for telluric bands");
1335 aFluxObj->tellbands = cpl_table_duplicate(aTellBands);
1339 unsigned int ntell =
sizeof(kTelluricBands) /
sizeof(kTelluricBands[0]) - 1;
1340 cpl_msg_debug(__func__,
"using builtin regions for telluric bands (%u "
1343 cpl_table *tb = aFluxObj->tellbands;
1345 for (k = 0; k < ntell; k++) {
1346 cpl_table_set_double(tb,
"lmin", k, kTelluricBands[k][0]);
1347 cpl_table_set_double(tb,
"lmax", k, kTelluricBands[k][1]);
1348 cpl_table_set_double(tb,
"bgmin", k, kTelluricBands[k][2]);
1349 cpl_table_set_double(tb,
"bgmax", k, kTelluricBands[k][3]);
1351 if (getenv(
"MUSE_DEBUG_FLUX") && atoi(getenv(
"MUSE_DEBUG_FLUX")) >= 2) {
1352 const char *fn =
"flux__tellregions.fits";
1353 cpl_msg_info(__func__,
"Saving telluric bands table as \"%s\"", fn);
1354 cpl_table_save(tb, NULL, NULL, fn, CPL_IO_CREATE);
1377 char *dodebug = getenv(
"MUSE_DEBUG_FLUX");
1378 if (!dodebug || (dodebug && atoi(dodebug) <= 0)) {
1381 char *fn = cpl_sprintf(
"flux__sens_%s.ascii", aName);
1382 FILE *fp = fopen(fn,
"w");
1384 cpl_table_dump(aFluxObj->sensitivity, 0,
1385 cpl_table_get_nrow(aFluxObj->sensitivity), fp);
1387 cpl_msg_debug(__func__,
"Written %"CPL_SIZE_FORMAT
" datapoints to \"%s\"",
1388 cpl_table_get_nrow(aFluxObj->sensitivity), fn);
1414 unsigned int aStar,
const cpl_table *aReference,
1415 double aAirmass,
const cpl_table *aExtinct)
1417 double crval = cpl_propertylist_get_double(aFluxObj->intimage->
header,
"CRVAL1"),
1418 cdelt = cpl_propertylist_get_double(aFluxObj->intimage->
header,
"CD1_1"),
1419 crpix = cpl_propertylist_get_double(aFluxObj->intimage->
header,
"CRPIX1"),
1421 int nlambda = cpl_image_get_size_x(aFluxObj->intimage->
data);
1423 aFluxObj->sensitivity = cpl_table_new(nlambda);
1424 cpl_table *sensitivity = aFluxObj->sensitivity;
1425 cpl_table_new_column(sensitivity,
"lambda", CPL_TYPE_DOUBLE);
1426 cpl_table_new_column(sensitivity,
"sens", CPL_TYPE_DOUBLE);
1427 cpl_table_new_column(sensitivity,
"serr", CPL_TYPE_DOUBLE);
1428 cpl_table_new_column(sensitivity,
"dq", CPL_TYPE_INT);
1429 cpl_table_set_column_format(sensitivity,
"dq",
"%u");
1430 float *data = cpl_image_get_data_float(aFluxObj->intimage->
data),
1431 *stat = cpl_image_get_data_float(aFluxObj->intimage->
stat);
1433 for (l = 0; l < nlambda; l++) {
1434 if (data[l + aStar*nlambda] <= 0. ||
1435 stat[l + aStar*nlambda] <= 0. ||
1436 stat[l + aStar*nlambda] == FLT_MAX) {
1440 double lambda = crval + cdelt * (l + 1 - crpix),
1442 extinct = !aExtinct ? 0.
1445 cpl_errorstate prestate = cpl_errorstate_get();
1450 if (!cpl_errorstate_is_equal(prestate)) {
1451 cpl_errorstate_set(prestate);
1457 double c = 2.5 * log10(data[l + aStar*nlambda]
1458 / exptime / cdelt / ref)
1459 + aAirmass * extinct,
1460 cerr = sqrt(pow(referr / ref, 2) + stat[l + aStar*nlambda]
1461 / pow(data[l + aStar*nlambda], 2))
1462 * 2.5 / CPL_MATH_LN10;
1463 cpl_table_set_double(sensitivity,
"lambda", idx, lambda);
1464 cpl_table_set_double(sensitivity,
"sens", idx, c);
1465 cpl_table_set_double(sensitivity,
"serr", idx, cerr);
1466 cpl_table_set_int(sensitivity,
"dq", idx, EURO3D_GOODPIXEL);
1470 cpl_table_set_size(sensitivity, idx);
1491 cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
1494 if (!aFluxObj->tellbands) {
1495 cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
1498 cpl_table *tsens = aFluxObj->sensitivity,
1499 *tb = aFluxObj->tellbands;
1503 > MUSE_MODE_WFM_NONAO_X;
1506 int irow, nfluxes = cpl_table_get_nrow(tsens);
1507 for (irow = 0; irow < nfluxes; irow++) {
1508 double lambda = cpl_table_get_double(tsens,
"lambda", irow, NULL);
1509 unsigned int dq = EURO3D_GOODPIXEL,
1510 k, nk = cpl_table_get_nrow(tb);
1511 for (k = 0; k < nk; k++) {
1512 double lmin = cpl_table_get_double(tb,
"lmin", k, NULL),
1513 lmax = cpl_table_get_double(tb,
"lmax", k, NULL);
1514 if (lambda >= lmin && lambda <= lmax) {
1515 dq |= EURO3D_TELLURIC;
1518 if (isnominal && lambda < kMuseNominalCutoff) {
1519 dq |= EURO3D_OUTSDRANGE;
1521 cpl_table_set_int(tsens,
"dq", irow, dq);
1549 static cpl_polynomial *
1551 double aLambda1,
double aLambda2,
1552 unsigned int aOrder,
double aRSigma,
double *aRMSE)
1554 cpl_table *tsens = aFluxObj->sensitivity;
1555 cpl_table_select_all(tsens);
1556 cpl_table_and_selected_int(tsens,
"dq", CPL_NOT_EQUAL_TO, EURO3D_GOODPIXEL);
1557 cpl_table_and_selected_int(tsens,
"dq", CPL_NOT_EQUAL_TO, EURO3D_TELLCOR);
1558 cpl_table_or_selected_double(tsens,
"lambda", CPL_LESS_THAN, aLambda1);
1559 cpl_table_or_selected_double(tsens,
"lambda", CPL_GREATER_THAN, aLambda2);
1561 cpl_table *tunwanted = cpl_table_extract_selected(tsens);
1562 cpl_table_erase_selected(tsens);
1563 muse_flux_response_dump_sensitivity(aFluxObj,
"fitinput");
1567 int nrow = cpl_table_get_nrow(tsens);
1568 cpl_matrix *lambdas = cpl_matrix_new(1, nrow);
1569 cpl_vector *sens = cpl_vector_new(nrow),
1570 *serr = cpl_vector_new(nrow);
1571 memcpy(cpl_matrix_get_data(lambdas),
1572 cpl_table_get_data_double_const(tsens,
"lambda"), nrow*
sizeof(
double));
1573 memcpy(cpl_vector_get_data(sens),
1574 cpl_table_get_data_double_const(tsens,
"sens"), nrow*
sizeof(
double));
1575 memcpy(cpl_vector_get_data(serr),
1576 cpl_table_get_data_double_const(tsens,
"serr"), nrow*
sizeof(
double));
1581 tsens, aOrder, aRSigma,
1583 int nout = cpl_vector_get_size(sens);
1585 cpl_msg_debug(__func__,
"transferred %d entries (%.3f...%.3f) for the "
1586 "order %u fit, %d entries are left, RMS %f", nrow, aLambda1,
1587 aLambda2, aOrder, nout, sqrt(mse));
1589 cpl_matrix_delete(lambdas);
1590 cpl_vector_delete(sens);
1591 cpl_vector_delete(serr);
1593 *aRMSE = mse / (nout - aOrder - 1);
1597 cpl_table_insert(tsens, tunwanted, nout);
1598 cpl_table_delete(tunwanted);
1625 cpl_table *tsens = aFluxObj->sensitivity,
1626 *tb = aFluxObj->tellbands;
1627 cpl_table_new_column(tsens,
"sens_orig", CPL_TYPE_DOUBLE);
1628 cpl_table_new_column(tsens,
"serr_orig", CPL_TYPE_DOUBLE);
1629 cpl_table_new_column(tsens,
"tellcor", CPL_TYPE_DOUBLE);
1630 unsigned int k, nk = cpl_table_get_nrow(tb);
1631 for (k = 0; k < nk; k++) {
1632 unsigned int order = 2;
1633 double lmin = cpl_table_get_double(tb,
"lmin", k, NULL),
1634 lmax = cpl_table_get_double(tb,
"lmax", k, NULL),
1635 bgmin = cpl_table_get_double(tb,
"bgmin", k, NULL),
1636 bgmax = cpl_table_get_double(tb,
"bgmax", k, NULL);
1642 cpl_polynomial *fit = muse_flux_response_fit(aFluxObj, bgmin, bgmax,
1644 cpl_msg_debug(__func__,
"Telluric region %u: %.2f...%.2f, reference region "
1645 "%.2f...%.2f", k + 1, lmin, lmax, bgmin, bgmax);
1647 cpl_polynomial_dump(fit, stdout);
1650 int irow, nrow = cpl_table_get_nrow(tsens);
1651 for (irow = 0; irow < nrow; irow++) {
1652 double lambda = cpl_table_get_double(tsens,
"lambda", irow, NULL);
1653 if (lambda >= lmin && lambda <= lmax &&
1654 (
unsigned int)cpl_table_get_int(tsens,
"dq", irow, NULL)
1655 == EURO3D_TELLURIC) {
1656 double origval = cpl_table_get_double(tsens,
"sens", irow, NULL),
1657 origerr = cpl_table_get_double(tsens,
"serr", irow, NULL),
1658 interpval = cpl_polynomial_eval_1d(fit, lambda, NULL);
1659 cpl_table_set_int(tsens,
"dq", irow, EURO3D_TELLCOR);
1660 cpl_table_set_double(tsens,
"sens_orig", irow, origval);
1661 cpl_table_set_double(tsens,
"sens", irow, interpval);
1664 cpl_table_set_double(tsens,
"serr_orig", irow, origerr);
1665 cpl_table_set_double(tsens,
"serr", irow, sqrt(origerr*origerr + rmse));
1667 if (interpval > origval) {
1669 double ftelluric = pow(10, -0.4 * (interpval - origval));
1670 cpl_table_set(tsens,
"tellcor", irow, ftelluric);
1672 cpl_table_set_double(tsens,
"tellcor", irow, 1.);
1676 cpl_polynomial_delete(fit);
1679 cpl_table_power_column(tsens,
"tellcor", 1. / aAirmass);
1702 cpl_table *tsens = aFluxObj->sensitivity;
1703 cpl_propertylist *order = cpl_propertylist_new();
1704 cpl_propertylist_append_bool(order,
"lambda", CPL_FALSE);
1705 cpl_table_sort(tsens, order);
1707 int nrow = cpl_table_get_nrow(tsens);
1709 double lambda1 = cpl_table_get_double(tsens,
"lambda", 0, NULL),
1710 serr1 = cpl_table_get_double(tsens,
"serr", 0, NULL);
1711 unsigned int dq = cpl_table_get_int(tsens,
"dq", 0, NULL);
1713 while (dq != EURO3D_GOODPIXEL && dq != EURO3D_TELLCOR) {
1714 lambda1 = cpl_table_get_double(tsens,
"lambda", ++irow, NULL);
1715 serr1 = cpl_table_get_double(tsens,
"serr", irow, NULL);
1716 dq = cpl_table_get_int(tsens,
"dq", irow, NULL);
1718 double lambda2 = cpl_table_get_double(tsens,
"lambda", nrow - 1, NULL),
1719 serr2 = cpl_table_get_double(tsens,
"serr", nrow - 1, NULL);
1720 dq = cpl_table_get_int(tsens,
"dq", nrow - 1, NULL);
1722 while (dq != EURO3D_GOODPIXEL && dq != EURO3D_TELLCOR) {
1723 lambda2 = cpl_table_get_double(tsens,
"lambda", --irow, NULL);
1724 serr2 = cpl_table_get_double(tsens,
"serr", irow, NULL);
1725 dq = cpl_table_get_int(tsens,
"dq", irow, NULL);
1727 const double ldist = 100.;
1728 cpl_polynomial *fit1 = muse_flux_response_fit(aFluxObj, lambda1, lambda1 + ldist,
1730 *fit2 = muse_flux_response_fit(aFluxObj, lambda2 - ldist, lambda2,
1732 nrow = cpl_table_get_nrow(tsens);
1733 double d1 = (lambda1 - kMuseLambdaMinX) / 10.,
1734 d2 = (kMuseLambdaMaxX - lambda2) / 10.;
1735 cpl_table_set_size(tsens, nrow + 20);
1738 for (l = kMuseLambdaMinX; l <= lambda1 - d1; l += d1) {
1739 double sens = cpl_polynomial_eval_1d(fit1, l, NULL),
1741 serr = 2. * (lambda1 - l) / ldist * serr1 + serr1;
1743 cpl_table_set_invalid(tsens,
"lambda", irow++);
1744 cpl_msg_debug(__func__,
"invalid blueward extrapolation: %.3f %f +/- %f",
1748 cpl_table_set_double(tsens,
"lambda", irow, l);
1749 cpl_table_set_double(tsens,
"sens", irow, sens);
1750 cpl_table_set_double(tsens,
"serr", irow, serr);
1751 cpl_table_set_int(tsens,
"dq", irow++, (
int)EURO3D_OUTSDRANGE);
1753 for (l = lambda2 + d2; l <= kMuseLambdaMaxX; l += d2) {
1754 double sens = cpl_polynomial_eval_1d(fit2, l, NULL),
1755 serr = 2. * (l - lambda2) / ldist * serr2 + serr2;
1757 cpl_table_set_invalid(tsens,
"lambda", irow++);
1758 cpl_msg_debug(__func__,
"invalid redward extrapolation: %.3f %f +/- %f",
1762 cpl_table_set_double(tsens,
"lambda", irow, l);
1763 cpl_table_set_double(tsens,
"sens", irow, sens);
1764 cpl_table_set_double(tsens,
"serr", irow, serr);
1765 cpl_table_set_int(tsens,
"dq", irow++, (
int)EURO3D_OUTSDRANGE);
1767 cpl_msg_debug(__func__,
"Extrapolated regions %.1f...%.1f Angstrom (using data "
1768 "from %.1f...%.1f A) and %.1f...%.1f Angstrom (using %.1f...%.1f A)",
1769 kMuseLambdaMinX, lambda1 - d1, lambda1, lambda1 + ldist,
1770 lambda2 + d2, kMuseLambdaMaxX, lambda2 - ldist, lambda2);
1772 cpl_polynomial_dump(fit1, stdout);
1773 cpl_polynomial_dump(fit2, stdout);
1776 cpl_polynomial_delete(fit1);
1777 cpl_polynomial_delete(fit2);
1779 cpl_table_select_all(tsens);
1780 cpl_table_and_selected_invalid(tsens,
"sens");
1781 cpl_table_erase_selected(tsens);
1783 cpl_table_sort(tsens, order);
1784 cpl_propertylist_delete(order);
1834 const cpl_table *aReference,
1835 const cpl_table *aTellBands,
1836 const cpl_table *aExtinct)
1838 cpl_ensure_code(aFluxObj && aFluxObj->intimage && aReference,
1839 CPL_ERROR_NULL_INPUT);
1840 cpl_ensure_code(aAirmass >= 1., CPL_ERROR_ILLEGAL_INPUT);
1842 cpl_msg_warning(__func__,
"Extinction table not given!");
1844 muse_flux_response_set_telluric_bands(aFluxObj, aTellBands);
1846 int nlambda = cpl_image_get_size_x(aFluxObj->intimage->
data),
1847 nobjects = cpl_image_get_size_y(aFluxObj->intimage->
data);
1851 for (n = 1; n <= nobjects; n++) {
1852 double this = cpl_image_get_flux_window(aFluxObj->intimage->
data, 1, n,
1854 cpl_msg_debug(__func__,
"flux(%d) = %e", n,
this);
1861 cpl_msg_warning(__func__,
"Assuming that the brightest star (%d of %d) is "
1862 "the reference source!", nstar, nobjects);
1865 muse_flux_response_sensitivity(aFluxObj,
1866 nstar - 1, aReference,
1867 aAirmass, aExtinct);
1868 muse_flux_response_dump_sensitivity(aFluxObj,
"initial");
1870 muse_flux_response_mark_questionable(aFluxObj);
1871 muse_flux_response_dump_sensitivity(aFluxObj,
"intermediate");
1873 cpl_table_select_all(aFluxObj->sensitivity);
1874 cpl_table_and_selected_int(aFluxObj->sensitivity,
"dq", CPL_EQUAL_TO,
1875 (
int)EURO3D_OUTSDRANGE);
1876 cpl_table_erase_selected(aFluxObj->sensitivity);
1877 muse_flux_response_dump_sensitivity(aFluxObj,
"intercut");
1881 muse_flux_response_telluric(aFluxObj, aAirmass);
1882 muse_flux_response_dump_sensitivity(aFluxObj,
"interpolated");
1884 muse_flux_response_extrapolate(aFluxObj);
1885 muse_flux_response_dump_sensitivity(aFluxObj,
"extrapolated");
1887 return CPL_ERROR_NONE;
1901 {
"lambda", CPL_TYPE_DOUBLE,
"Angstrom",
"%7.2f",
"wavelength", CPL_TRUE },
1902 {
"response", CPL_TYPE_DOUBLE,
1903 "2.5*log10((count/s/Angstrom)/(erg/s/cm**2/Angstrom))",
"%.4e",
1904 "instrument response derived from standard star", CPL_TRUE },
1905 {
"resperr", CPL_TYPE_DOUBLE,
1906 "2.5*log10((count/s/Angstrom)/(erg/s/cm**2/Angstrom))",
"%.4e",
1907 "instrument response error derived from standard star", CPL_TRUE },
1908 { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
1922 {
"lambda", CPL_TYPE_DOUBLE,
"Angstrom",
"%7.2f",
"wavelength", CPL_TRUE },
1923 {
"ftelluric", CPL_TYPE_DOUBLE,
"",
"%.5f",
1924 "the telluric correction factor, normalized to an airmass of 1", CPL_TRUE },
1925 {
"ftellerr", CPL_TYPE_DOUBLE,
"",
"%.5f",
1926 "the error of the telluric correction factor", CPL_TRUE },
1927 { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
1948 cpl_table *aResp,
double aHalfwidth)
1951 int i, n = cpl_table_get_nrow(aResp);
1952 for (i = 0; i < n; i++) {
1953 double lambda = cpl_table_get_double(aResp,
"lambda", i, NULL);
1954 int j = i, j1 = i, j2 = i;
1957 lambda - cpl_table_get_double(aResp,
"lambda", j, NULL) <= aHalfwidth) {
1962 cpl_table_get_double(aResp,
"lambda", j, NULL) - lambda <= aHalfwidth) {
1965 double *sens = cpl_table_get_data_double(aFluxObj->sensitivity,
"sens"),
1966 *serr = cpl_table_get_data_double(aFluxObj->sensitivity,
"serr");
1967 cpl_vector *v = cpl_vector_wrap(j2 - j1 + 1, sens + j1),
1968 *ve = cpl_vector_wrap(j2 - j1 + 1, serr + j1);
1969 double median = cpl_vector_get_median_const(v),
1971 mederr = cpl_vector_get_median_const(ve);
1973 mdev = cpl_table_get_double(aFluxObj->sensitivity,
"serr", j1, NULL);
1975 if (mdev < mederr) {
1979 cpl_msg_debug(__func__,
"%d %.3f %d...%d --> %f +/- %f", i, lambda, j1, j2,
1982 cpl_vector_unwrap(v);
1983 cpl_vector_unwrap(ve);
1984 cpl_table_set_double(aResp,
"response", i, median);
1985 cpl_table_set_double(aResp,
"resperr", i, mdev);
2011 cpl_ensure_code(aFluxObj && aFluxObj->sensitivity, CPL_ERROR_NULL_INPUT);
2012 int nrow = cpl_table_get_nrow(aFluxObj->sensitivity);
2015 const double *lambdas = cpl_table_get_data_double_const(aFluxObj->sensitivity,
2017 cpl_table_copy_data_double(resp,
"lambda", lambdas);
2018 cpl_table_fill_column_window_double(resp,
"response", 0, nrow, 0);
2020 cpl_table_fill_column_window_double(resp,
"resperr", 0, nrow, 0.008);
2026 aFluxObj->response = resp;
2027 return CPL_ERROR_NONE;
2053 cpl_ensure_code(aFluxObj && aFluxObj->sensitivity, CPL_ERROR_NULL_INPUT);
2054 cpl_table *tsens = aFluxObj->sensitivity;
2055 int nrow = cpl_table_get_nrow(tsens);
2058 cpl_table_fill_column_window_double(tell,
"lambda", 0, nrow, 0);
2059 cpl_table_copy_data_double(tell,
"lambda",
2060 cpl_table_get_data_double_const(tsens,
"lambda"));
2061 cpl_table_fill_column_window_double(tell,
"ftelluric", 0, nrow, 0);
2062 cpl_table_copy_data_double(tell,
"ftelluric",
2063 cpl_table_get_data_double_const(tsens,
"tellcor"));
2065 #define TELL_MAX_ERR 0.1
2066 #define TELL_MIN_ERR 1e-4
2067 cpl_table_fill_column_window_double(tell,
"ftellerr", 0, nrow, TELL_MAX_ERR);
2071 cpl_table_duplicate_column(tell,
"tellcor", tsens,
"tellcor");
2073 cpl_table_unselect_all(tell);
2075 for (irow = 0; irow < nrow; irow++) {
2077 cpl_table_get_double(tell,
"tellcor", irow, &err);
2082 cpl_errorstate state = cpl_errorstate_get();
2083 double ftellcor = cpl_table_get_double(tell,
"tellcor", irow - 1, &err);
2084 if (!cpl_errorstate_is_equal(state)) {
2085 cpl_errorstate_set(state);
2087 if (err == 0 && ftellcor != 1.) {
2088 cpl_table_set_double(tell,
"ftelluric", irow, 1.);
2092 state = cpl_errorstate_get();
2093 ftellcor = cpl_table_get_double(tell,
"tellcor", irow + 1, &err);
2094 if (!cpl_errorstate_is_equal(state)) {
2095 cpl_errorstate_set(state);
2097 if (err == 0 && ftellcor != 1.) {
2098 cpl_table_set_double(tell,
"ftelluric", irow, 1.);
2101 cpl_table_select_row(tell, irow);
2103 cpl_table_erase_selected(tell);
2104 cpl_table_erase_column(tell,
"tellcor");
2107 nrow = cpl_table_get_nrow(tell);
2108 for (irow = 0; irow < nrow; irow++) {
2110 double dftell = 1. - cpl_table_get_double(tell,
"ftelluric", irow, &err),
2111 ftellerr = cpl_table_get_double(tell,
"ftellerr", irow, &err);
2112 if (ftellerr > dftell) {
2113 ftellerr = fmax(dftell, TELL_MIN_ERR);
2114 cpl_table_set_double(tell,
"ftellerr", irow, ftellerr);
2118 aFluxObj->telluric = tell;
2119 return CPL_ERROR_NONE;
2162 const cpl_table *aExtinction,
const cpl_table *aTelluric)
2164 cpl_ensure_code(aPixtable && aPixtable->
header && aResponse,
2165 CPL_ERROR_NULL_INPUT);
2166 const char *unitdata = cpl_table_get_column_unit(aPixtable->
table,
2167 MUSE_PIXTABLE_DATA);
2168 cpl_ensure_code(unitdata && !strncmp(unitdata,
"count", 6),
2169 CPL_ERROR_INCOMPATIBLE_INPUT);
2172 cpl_msg_warning(__func__,
"%s missing!", MUSE_TAG_EXTINCT_TABLE);
2176 if (exptime <= 0.) {
2177 cpl_msg_warning(__func__,
"Non-positive EXPTIME, not doing flux calibration!");
2178 return CPL_ERROR_ILLEGAL_INPUT;
2182 cpl_msg_warning(__func__,
"Airmass unknown, not doing extinction "
2183 "correction: %s", cpl_error_get_message());
2188 cpl_table *telluric = NULL;
2191 telluric = cpl_table_duplicate(aTelluric);
2194 cpl_table_power_column(telluric,
"ftelluric", -airmass);
2197 cpl_msg_info(__func__,
"Starting flux calibration (exptime=%.2fs, airmass=%.4f),"
2198 " %s telluric correction", exptime, airmass,
2199 aTelluric ?
"with" :
"without ("MUSE_TAG_STD_TELLURIC
" not given)");
2200 float *lambda = cpl_table_get_data_float(aPixtable->
table, MUSE_PIXTABLE_LAMBDA),
2201 *data = cpl_table_get_data_float(aPixtable->
table, MUSE_PIXTABLE_DATA),
2202 *stat = cpl_table_get_data_float(aPixtable->
table, MUSE_PIXTABLE_STAT);
2204 #pragma omp parallel for default(none) \
2205 shared(aExtinction, aResponse, airmass, data, exptime, lambda, nrow, \
2207 for (i = 0; i < nrow; i++) {
2209 double ddata = data[i], dstat = stat[i];
2213 double fext = pow(10., 0.4 * airmass
2218 printf(
"%f, data/stat = %f/%f -> ", fext, ddata, dstat);
2221 dstat *= (fext * fext);
2223 printf(
" --> %f/%f\n", ddata, dstat), fflush(NULL);
2227 double dlambda = kMuseSpectralSamplingA,
2243 double dlambda = 1.26279
2244 + 7.12041e-05 * lambda[i]
2245 + -2.97879e-08 * pow(lambda[i], 2)
2246 + 4.653e-12 * pow(lambda[i], 3)
2247 + -3.30195e-16 * pow(lambda[i], 4)
2248 + 8.70389e-21 * pow(lambda[i], 5),
2256 resp = pow(10., 0.4 * resp);
2259 rerr = rerr * CPL_MATH_LN10 * resp / 2.5;
2261 printf(
"%f/%f/%f, %e/%e, data/stat = %e/%e -> ", lambda[i], dlambda, dlamerr, resp, rerr,
2264 dstat = dstat * pow((1./(resp * exptime * dlambda)), 2)
2265 + pow(ddata * rerr / (resp*resp * exptime * dlambda), 2)
2266 + pow(ddata * dlamerr / (resp * exptime * dlambda*dlambda), 2);
2267 ddata /= (resp * exptime * dlambda);
2269 printf(
"%e/%e\n", ddata, dstat), fflush(NULL);
2274 ddata *= kMuseFluxUnitFactor;
2275 dstat *= kMuseFluxStatFactor;
2279 if (lambda[i] < kTelluricBands[0][0] || !telluric) {
2287 data[i] = ddata * tell;
2288 stat[i] = tell*tell * dstat + ddata*ddata * terr*terr;
2290 cpl_table_delete(telluric);
2293 cpl_table_set_column_unit(aPixtable->
table, MUSE_PIXTABLE_DATA,
2294 kMuseFluxUnitString);
2295 cpl_table_set_column_unit(aPixtable->
table, MUSE_PIXTABLE_STAT,
2296 kMuseFluxStatString);
2302 MUSE_HDR_PT_FLUXCAL_COMMENT);
2303 return CPL_ERROR_NONE;
Structure definition of a MUSE datacube.
void muse_image_delete(muse_image *aImage)
Deallocate memory associated to a muse_image object.
cpl_size muse_pixtable_get_nrow(muse_pixtable *aPixtable)
get the number of rows within the pixel table
#define MUSE_HDR_PT_FLUXCAL
cpl_image * data
the data extension
const muse_cpltable_def muse_flux_responsetable_def[]
MUSE response table definition.
muse_flux_profile_type
Type of optimal profile to use.
cpl_image * stat
the statistics extension
void muse_datacube_delete(muse_datacube *aCube)
Deallocate memory associated to a muse_datacube object.
cpl_error_code muse_flux_reference_table_check(cpl_table *aTable)
Check and/or adapt the standard flux reference table format.
const muse_cpltable_def muse_response_tellbands_def[]
Table definition for a telluric bands table.
Structure definition of MUSE three extension FITS file.
muse_flux_object * muse_flux_object_new(void)
Allocate memory for a new muse_flux_object object.
static double muse_flux_reference_table_sampling(cpl_table *aTable)
Compute average sampling for a MUSE-format flux reference table.
cpl_table * table
The pixel table.
cpl_propertylist * header
the FITS header
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.
cpl_error_code muse_datacube_save(muse_datacube *aCube, const char *aFilename)
Save the three cube extensions and the FITS headers of a MUSE datacube to a file. ...
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.
cpl_error_code muse_flux_get_telluric_table(muse_flux_object *aFluxObj)
Get the table of the telluric correction.
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.
double muse_pfits_get_fwhm_end(const cpl_propertylist *aHeaders)
find out the ambient seeing at end of exposure (in arcsec)
cpl_error_code muse_flux_get_response_table(muse_flux_object *aFluxObj)
Get the table of the standard star response function.
Structure definition of MUSE pixel table.
Flux object to store data needed while computing the flux calibration.
cpl_error_code muse_utils_fit_moffat_2d(const cpl_matrix *aPositions, const cpl_vector *aValues, const cpl_vector *aErrors, cpl_array *aParams, cpl_array *aPErrors, const cpl_array *aPFlags, double *aRMS, double *aRedChisq)
Fit a 2D Moffat function to a given set of data.
muse_resampling_params * muse_resampling_params_new(muse_resampling_type aMethod)
Create the resampling parameters structure.
cpl_polynomial * muse_utils_iterate_fit_polynomial(cpl_matrix *aPos, cpl_vector *aVal, cpl_vector *aErr, cpl_table *aExtra, const unsigned int aOrder, const double aRSigma, double *aMSE, double *aChiSq)
Iterate a polynomial fit.
static void muse_flux_get_response_table_smooth(muse_flux_object *aFluxObj, cpl_table *aResp, double aHalfwidth)
Get the table of the standard star response function.
double muse_pfits_get_fwhm_start(const cpl_propertylist *aHeaders)
find out the ambient seeing at start of exposure (in arcsec)
cpl_error_code muse_flux_calibrate(muse_pixtable *aPixtable, const cpl_table *aResponse, const cpl_table *aExtinction, const cpl_table *aTelluric)
Convert the input pixel table from counts to fluxes.
cpl_imagelist * data
the cube containing the actual data values
muse_flux_interpolation_type
Type of table interpolation to use.
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
const muse_cpltable_def muse_flux_tellurictable_def[]
MUSE telluric correction table definition.
cpl_error_code muse_pixtable_save(muse_pixtable *aPixtable, const char *aFilename)
Save a MUSE pixel table to a file on disk.
cpl_error_code muse_flux_integrate_std(muse_pixtable *aPixtable, muse_flux_profile_type aProfile, muse_flux_object *aFluxObj)
Integrate the flux of the standard star(s) in the field over all wavelengths.
cpl_propertylist * header
the FITS header
cpl_error_code muse_flux_response_compute(muse_flux_object *aFluxObj, double aAirmass, const cpl_table *aReference, const cpl_table *aTellBands, const cpl_table *aExtinct)
Compare measured flux distribution over wavelength with calibrated stellar fluxes and derive instrume...
cpl_error_code muse_image_save(muse_image *aImage, const char *aFilename)
Save the three image extensions and the FITS headers of a MUSE image to a file.
double muse_pfits_get_exptime(const cpl_propertylist *aHeaders)
find out the exposure time
Definition of a cpl table structure.
void muse_resampling_params_delete(muse_resampling_params *aParams)
Delete a resampling parameters structure.
double muse_cplvector_get_adev_const(const cpl_vector *aVector, double aCenter)
Compute the average absolute deviation of a (constant) vector.
muse_image * muse_image_new(void)
Allocate memory for a new muse_image object.
void muse_flux_object_delete(muse_flux_object *aFluxObj)
Deallocate memory associated to a muse_flux_object.
int muse_quality_image_reject_using_dq(cpl_image *aData, cpl_image *aDQ, cpl_image *aStat)
Reject pixels of one or two images on a DQ image.
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.