33 #include "muse_flux.h"
34 #include "muse_instrument.h"
36 #include "muse_astro.h"
37 #include "muse_cplwrappers.h"
38 #include "muse_pfits.h"
39 #include "muse_quality.h"
40 #include "muse_resampling.h"
41 #include "muse_utils.h"
95 aFluxObj->
cube = NULL;
100 cpl_table_delete(aFluxObj->
response);
102 cpl_table_delete(aFluxObj->
telluric);
121 cpl_ensure(aTable, CPL_ERROR_NULL_INPUT, 0.);
122 cpl_table_unselect_all(aTable);
123 cpl_table_or_selected_double(aTable,
"lambda", CPL_NOT_LESS_THAN,
124 kMuseNominalLambdaMin);
125 cpl_table_and_selected_double(aTable,
"lambda", CPL_NOT_GREATER_THAN,
126 kMuseNominalLambdaMax);
127 cpl_size nsel = cpl_table_count_selected(aTable);
128 cpl_array *asel = cpl_table_where_selected(aTable);
129 cpl_size *sel = cpl_array_get_data_cplsize(asel);
130 double lmin = cpl_table_get_double(aTable,
"lambda", sel[0], NULL),
131 lmax = cpl_table_get_double(aTable,
"lambda", sel[nsel - 1], NULL);
132 cpl_array_delete(asel);
133 return (lmax - lmin) / nsel;
170 cpl_ensure_code(aTable, CPL_ERROR_NULL_INPUT);
172 const char *flam1 =
"erg/s/cm**2/Angstrom",
173 *flam2 =
"erg/s/cm^2/Angstrom";
175 cpl_error_code rc = CPL_ERROR_NONE;
176 cpl_errorstate prestate = cpl_errorstate_get();
178 if (cpl_table_has_column(aTable,
"lambda") &&
179 cpl_table_has_column(aTable,
"flux") &&
180 cpl_table_get_column_unit(aTable,
"lambda") &&
181 cpl_table_get_column_unit(aTable,
"flux") &&
182 !strncmp(cpl_table_get_column_unit(aTable,
"lambda"),
"Angstrom", 9) &&
183 (!strncmp(cpl_table_get_column_unit(aTable,
"flux"), flam1, strlen(flam1)) ||
184 !strncmp(cpl_table_get_column_unit(aTable,
"flux"), flam2, strlen(flam2)))) {
187 if (cpl_table_get_column_type(aTable,
"lambda") != CPL_TYPE_DOUBLE) {
188 cpl_msg_debug(__func__,
"Casting lambda column to double");
189 cpl_table_cast_column(aTable,
"lambda", NULL, CPL_TYPE_DOUBLE);
191 if (cpl_table_get_column_type(aTable,
"flux") != CPL_TYPE_DOUBLE) {
192 cpl_msg_debug(__func__,
"Casting flux column to double");
193 cpl_table_cast_column(aTable,
"flux", NULL, CPL_TYPE_DOUBLE);
196 if (cpl_table_has_column(aTable,
"fluxerr")) {
197 if (cpl_table_get_column_type(aTable,
"fluxerr") != CPL_TYPE_DOUBLE) {
198 cpl_msg_debug(__func__,
"Casting fluxerr column to double");
199 cpl_table_cast_column(aTable,
"fluxerr", NULL, CPL_TYPE_DOUBLE);
201 const char *unit = cpl_table_get_column_unit(aTable,
"fluxerr");
202 if (!unit || (strncmp(unit, flam1, strlen(flam1)) &&
203 strncmp(unit, flam2, strlen(flam2)))) {
204 cpl_msg_debug(__func__,
"Erasing fluxerr column because of unexpected "
206 cpl_table_erase_column(aTable,
"fluxerr");
209 cpl_msg_info(__func__,
"Found MUSE format, average sampling %.3f Angstrom/bin"
211 }
else if (cpl_table_has_column(aTable,
"WAVELENGTH") &&
212 cpl_table_has_column(aTable,
"FLUX") &&
213 cpl_table_get_column_unit(aTable,
"WAVELENGTH") &&
214 cpl_table_get_column_unit(aTable,
"FLUX") &&
215 !strncmp(cpl_table_get_column_unit(aTable,
"WAVELENGTH"),
"ANGSTROMS", 10) &&
216 !strncmp(cpl_table_get_column_unit(aTable,
"FLUX"),
"FLAM", 5)) {
218 printf(
"input HST CALSPEC table:\n");
219 cpl_table_dump_structure(aTable, stdout);
220 cpl_table_dump(aTable, cpl_table_get_nrow(aTable)/2, 3, stdout);
224 cpl_table_cast_column(aTable,
"WAVELENGTH",
"lambda", CPL_TYPE_DOUBLE);
225 cpl_table_cast_column(aTable,
"FLUX",
"flux", CPL_TYPE_DOUBLE);
226 cpl_table_erase_column(aTable,
"WAVELENGTH");
227 cpl_table_erase_column(aTable,
"FLUX");
228 cpl_table_set_column_unit(aTable,
"lambda",
"Angstrom");
229 cpl_table_set_column_unit(aTable,
"flux", flam1);
232 if (cpl_table_has_column(aTable,
"STATERROR") &&
233 cpl_table_has_column(aTable,
"SYSERROR") &&
234 cpl_table_get_column_unit(aTable,
"STATERROR") &&
235 cpl_table_get_column_unit(aTable,
"SYSERROR") &&
236 !strncmp(cpl_table_get_column_unit(aTable,
"STATERROR"),
"FLAM", 5) &&
237 !strncmp(cpl_table_get_column_unit(aTable,
"SYSERROR"),
"FLAM", 5)) {
240 cpl_table_cast_column(aTable,
"STATERROR",
"fluxerr", CPL_TYPE_DOUBLE);
241 cpl_table_erase_column(aTable,
"STATERROR");
242 cpl_table_cast_column(aTable,
"SYSERROR", NULL, CPL_TYPE_DOUBLE);
243 cpl_table_power_column(aTable,
"fluxerr", 2);
244 cpl_table_power_column(aTable,
"SYSERROR", 2);
245 cpl_table_add_columns(aTable,
"fluxerr",
"SYSERROR");
246 cpl_table_erase_column(aTable,
"SYSERROR");
247 cpl_table_power_column(aTable,
"fluxerr", 0.5);
248 cpl_table_set_column_unit(aTable,
"fluxerr", flam1);
255 if (cpl_table_has_column(aTable,
"FWHM")) {
256 cpl_table_erase_column(aTable,
"FWHM");
258 if (cpl_table_has_column(aTable,
"DATAQUAL")) {
259 cpl_table_erase_column(aTable,
"DATAQUAL");
261 if (cpl_table_has_column(aTable,
"TOTEXP")) {
262 cpl_table_erase_column(aTable,
"TOTEXP");
265 cpl_size irow, nrow = cpl_table_get_nrow(aTable);
266 for (irow = 0; irow < nrow; irow++) {
267 double lambda = cpl_table_get_double(aTable,
"lambda", irow, NULL);
268 cpl_table_set_double(aTable,
"lambda", irow,
272 printf(
"converted HST CALSPEC table:\n");
273 cpl_table_dump_structure(aTable, stdout);
274 cpl_table_dump(aTable, cpl_table_get_nrow(aTable)/2, 3, stdout);
277 cpl_msg_info(__func__,
"Found HST CALSPEC format on input, converted to "
278 "MUSE format; average sampling %.3f Angstrom/bin over MUSE "
279 "range (assumed vacuum wavelengths on input, converted to air).",
282 cpl_msg_error(__func__,
"Unknown format found!");
284 cpl_table_dump_structure(aTable, stdout);
286 rc = CPL_ERROR_INCOMPATIBLE_INPUT;
290 if (!cpl_errorstate_is_equal(prestate)) {
291 rc = cpl_error_get_code();
343 cpl_ensure(aResponse, CPL_ERROR_NULL_INPUT, rv);
344 int size = cpl_table_get_nrow(aResponse);
345 cpl_ensure(size > 0, cpl_error_get_code(), rv);
348 const double *lbda = cpl_table_get_data_double_const(aResponse,
"lambda"),
349 *resp = NULL, *rerr = NULL;
352 resp = cpl_table_get_data_double_const(aResponse,
"throughput");
355 resp = cpl_table_get_data_double_const(aResponse,
"response");
357 rerr = cpl_table_get_data_double_const(aResponse,
"resperr");
361 resp = cpl_table_get_data_double_const(aResponse,
"flux");
363 rerr = cpl_table_get_data_double_const(aResponse,
"fluxerr");
367 resp = cpl_table_get_data_double_const(aResponse,
"extinction");
370 resp = cpl_table_get_data_double_const(aResponse,
"ftelluric");
372 rerr = cpl_table_get_data_double_const(aResponse,
"ftellerr");
376 cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
379 cpl_ensure(lbda && resp, cpl_error_get_code(), rv);
381 cpl_ensure(rerr, cpl_error_get_code(), rv);
385 if (aLambda < lbda[0]) {
388 if (aLambda > lbda[size-1]) {
393 double response = rv, resperror = 0.;
394 int l = 0, r = size - 1,
397 if (aLambda >= lbda[m] && aLambda <= lbda[m+1]) {
399 double lquot = (aLambda - lbda[m]) / (lbda[m+1] - lbda[m]);
400 response = resp[m] + (resp[m+1] - resp[m]) * lquot;
405 resperror = sqrt(pow(rerr[m] * (1 - lquot), 2.)
406 + pow(rerr[m+1] * lquot, 2.));
409 cpl_msg_debug(__func__,
"Found at m=%d (%f: %f+/-%f) / "
410 "m+1=%d (%f: %f+/-%f) -> %f: %f+/-%f",
411 m, lbda[m], resp[m], rerr ? rerr[m] : 0.,
412 m+1, lbda[m+1], resp[m+1], rerr ? rerr[m+1] : 0.,
413 aLambda, response, resperror);
418 if (aLambda < lbda[m]) {
421 if (aLambda > lbda[m]) {
428 cpl_msg_debug(__func__,
"Response %g+/-%g at lambda=%fA", response, resperror,
431 if (aError && rerr) {
454 muse_flux_image_sky(cpl_image *aImage,
double aX,
double aY,
double aHalfSize,
455 unsigned int aDSky,
float *aError)
461 int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
462 y1 = aY - aHalfSize, y2 = aY + aHalfSize;
463 unsigned char nskyarea = 0;
464 double skylevel = 0., skyerror = 0.;
466 cpl_errorstate state = cpl_errorstate_get();
467 cpl_stats_mode mode = CPL_STATS_MEDIAN | CPL_STATS_MEDIAN_DEV;
468 cpl_stats *s = cpl_stats_new_from_image_window(aImage, mode,
469 x1 - aDSky, y1, x1 - 1, y2);
474 skylevel += cpl_stats_get_median(s);
475 skyerror += pow(cpl_stats_get_median_dev(s), 2);
479 s = cpl_stats_new_from_image_window(aImage,mode,
480 x2 + 1, y1, x2 + aDSky, y2);
483 skylevel += cpl_stats_get_median(s);
484 skyerror += pow(cpl_stats_get_median_dev(s), 2);
488 s = cpl_stats_new_from_image_window(aImage,mode,
489 x1, y1 - aDSky, x2, y1 - 1);
492 skylevel += cpl_stats_get_median(s);
493 skyerror += pow(cpl_stats_get_median_dev(s), 2);
497 s = cpl_stats_new_from_image_window(aImage,mode,
498 x1, y2 + 1, x2, y2 + aDSky);
501 skylevel += cpl_stats_get_median(s);
502 skyerror += pow(cpl_stats_get_median_dev(s), 2);
508 skylevel /= nskyarea;
509 skyerror = sqrt(skyerror) / nskyarea;
510 if (!cpl_errorstate_is_equal(state)) {
511 cpl_errorstate_set(state);
514 cpl_msg_debug(__func__,
"skylevel = %f +/- %f (%u sky areas)",
515 skylevel, skyerror, nskyarea);
542 muse_flux_image_gaussian(cpl_image *aImage, cpl_image *aImErr,
double aX,
543 double aY,
double aHalfSize,
unsigned int aDSky,
544 unsigned int aMaxBad,
float *aFErr)
548 UNUSED_ARGUMENT(aMaxBad);
554 cpl_array *params = cpl_array_new(7, CPL_TYPE_DOUBLE),
555 *parerr = cpl_array_new(7, CPL_TYPE_DOUBLE);
559 cpl_errorstate state = cpl_errorstate_get();
560 double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky, NULL);
561 if (!cpl_errorstate_is_equal(state)) {
564 cpl_errorstate_set(state);
566 cpl_array_set_double(params, 0, skylevel);
567 cpl_array_set_double(params, 3, aX);
568 cpl_array_set_double(params, 4, aY);
569 double rms = 0, chisq = 0;
571 int nx = cpl_image_get_size_x(aImage),
572 ny = cpl_image_get_size_y(aImage),
573 xsize = fmin(aHalfSize, fmin(aX - 1., nx - aX)) * 2,
574 ysize = fmin(aHalfSize, fmin(aY - 1., ny - aY)) * 2;
575 cpl_error_code rc = cpl_fit_image_gaussian(aImage, aImErr, aX, aY, xsize, ysize,
576 params, parerr, NULL, &rms, &chisq,
577 NULL, NULL, NULL, NULL, NULL);
578 if (rc != CPL_ERROR_NONE) {
579 if (rc != CPL_ERROR_ILLEGAL_INPUT) {
580 cpl_msg_debug(__func__,
"rc = %d: %s", rc, cpl_error_get_message());
582 cpl_array_delete(params);
583 cpl_array_delete(parerr);
586 double flux = cpl_array_get_double(params, 1, NULL),
587 ferr = cpl_array_get_double(parerr, 1, NULL);
589 double fwhmx = cpl_array_get_double(params, 5, NULL) * CPL_MATH_FWHM_SIG,
590 fwhmy = cpl_array_get_double(params, 6, NULL) * CPL_MATH_FWHM_SIG;
591 cpl_msg_debug(__func__,
"%.3f,%.3f: %g+/-%g (bg: %g, FWHM: %.3f,%.3f, %g, %g)",
592 cpl_array_get_double(params, 3, NULL), cpl_array_get_double(params, 4, NULL),
593 flux, ferr, cpl_array_get_double(params, 0, NULL), fwhmx, fwhmy,
597 cpl_msg_debug(__func__,
"skylevel = %f", cpl_array_get_double(params, 0, NULL));
598 cpl_msg_debug(__func__,
"measured flux %f +/- %f", flux, ferr);
600 cpl_array_delete(params);
601 cpl_array_delete(parerr);
629 muse_flux_image_moffat(cpl_image *aImage, cpl_image *aImErr,
double aX,
630 double aY,
double aHalfSize,
unsigned int aDSky,
631 unsigned int aMaxBad,
float *aFErr)
637 int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
638 y1 = aY - aHalfSize, y2 = aY + aHalfSize,
639 nx = cpl_image_get_size_x(aImage),
640 ny = cpl_image_get_size_y(aImage);
653 int npoints = (x2 - x1 + 1) * (y2 - y1 + 1);
655 cpl_matrix *pos = cpl_matrix_new(npoints, 2);
656 cpl_vector *values = cpl_vector_new(npoints),
657 *errors = cpl_vector_new(npoints);
658 float *derr = cpl_image_get_data_float(aImErr);
660 for (i = x1; i <= x2; i++) {
662 for (j = y1; j <= y2; j++) {
664 double data = cpl_image_get(aImage, i, j, &err);
668 cpl_matrix_set(pos, idx, 0, i);
669 cpl_matrix_set(pos, idx, 1, j);
670 cpl_vector_set(values, idx, data);
671 cpl_vector_set(errors, idx, derr[(i-1) + (j-1)*nx]);
677 if (idx < 16 || (npoints - idx) > (
int)aMaxBad) {
678 cpl_matrix_delete(pos);
679 cpl_vector_delete(values);
680 cpl_vector_delete(errors);
683 cpl_matrix_set_size(pos, idx, 2);
684 cpl_vector_set_size(values, idx);
685 cpl_vector_set_size(errors, idx);
687 cpl_array *params = cpl_array_new(8, CPL_TYPE_DOUBLE),
688 *parerr = cpl_array_new(8, CPL_TYPE_DOUBLE);
690 cpl_errorstate state = cpl_errorstate_get();
691 double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky, NULL);
692 if (!cpl_errorstate_is_equal(state)) {
695 cpl_errorstate_set(state);
697 cpl_array_set_double(params, 0, skylevel);
698 cpl_array_set_double(params, 2, aX);
699 cpl_array_set_double(params, 3, aY);
700 double rms = 0, chisq = 0;
702 params, parerr, NULL,
704 cpl_matrix_delete(pos);
705 cpl_vector_delete(values);
706 cpl_vector_delete(errors);
707 if (rc != CPL_ERROR_NONE) {
708 if (rc != CPL_ERROR_ILLEGAL_INPUT) {
709 cpl_msg_debug(__func__,
"rc = %d: %s", rc, cpl_error_get_message());
711 cpl_array_delete(params);
712 cpl_array_delete(parerr);
716 double flux = cpl_array_get_double(params, 1, NULL);
718 *aFErr = cpl_array_get_double(parerr, 1, NULL);
721 cpl_msg_debug(__func__,
"skylevel = %f", cpl_array_get_double(params, 0, NULL));
722 cpl_msg_debug(__func__,
"measured flux %f +/- %f", flux, cpl_array_get_double(parerr, 1, NULL));
724 cpl_array_delete(params);
725 cpl_array_delete(parerr);
750 muse_flux_image_square(cpl_image *aImage, cpl_image *aImErr,
double aX,
751 double aY,
double aHalfSize,
unsigned int aDSky,
752 unsigned int aMaxBad,
float *aFErr)
757 int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
758 y1 = aY - aHalfSize, y2 = aY + aHalfSize,
759 nx = cpl_image_get_size_x(aImage),
760 ny = cpl_image_get_size_y(aImage);
774 cpl_errorstate state = cpl_errorstate_get();
775 double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky,
777 if (!cpl_errorstate_is_equal(state)) {
780 cpl_errorstate_set(state);
781 cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
787 cpl_image *region = cpl_image_extract(aImage, x1, y1, x2, y2);
789 cpl_msg_debug(__func__,
"region [%d:%d,%d:%d] %"CPL_SIZE_FORMAT
" bad pixels",
790 x1, y1, x2, y2, cpl_image_count_rejected(region));
792 if (cpl_image_count_rejected(region) > aMaxBad) {
793 cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
794 cpl_image_delete(region);
797 cpl_image *regerr = cpl_image_extract(aImErr, x1, y1, x2, y2);
798 if (cpl_image_count_rejected(region) > 0) {
799 cpl_detector_interpolate_rejected(region);
800 cpl_detector_interpolate_rejected(regerr);
805 int npoints = (x2 - x1 + 1) * (y2 - y1 + 1),
809 nsky = 2 * aDSky * (x2 - x1 + y2 - y1 + 2);
810 double flux = cpl_image_get_flux(region) - skylevel * npoints,
811 ferr = sqrt(cpl_image_get_sqflux(regerr)
812 + npoints * skyerror*skyerror * (1. + (
double)npoints / nsky));
814 cpl_msg_debug(__func__,
"measured flux %f +/- %f (%d object pixels, %d pixels"
815 " with %f with sky %f)", flux, ferr, npoints, nsky,
816 cpl_image_get_flux(region), cpl_image_get_sqflux(regerr));
818 cpl_image_delete(region);
819 cpl_image_delete(regerr);
848 muse_flux_image_circle(cpl_image *aImage, cpl_image *aImErr,
double aX,
849 double aY,
double aAper,
double aAnnu,
double aDAnnu,
850 unsigned int aMaxBad,
float *aFErr)
855 double rmax = ceil(fmax(aAper, aAnnu + aDAnnu));
856 int x1 = aX - rmax, x2 = aX + rmax,
857 y1 = aY - rmax, y2 = aY + rmax,
858 nx = cpl_image_get_size_x(aImage),
859 ny = cpl_image_get_size_y(aImage);
874 cpl_vector *vbg = cpl_vector_new((x2 - x1 + 1) * (y2 - y1 + 1)),
875 *vbe = cpl_vector_new((x2 - x1 + 1) * (y2 - y1 + 1));
876 unsigned int nbad = 0, nbg = 0;
878 for (i = x1; i <= x2; i++) {
880 for (j = y1; j <= y2; j++) {
881 double r = sqrt(pow(aX - i, 2) + pow(aY - j, 2));
883 nbad += cpl_image_is_rejected(aImage, i, j) == 1;
885 if (r < aAnnu || r > aAnnu + aDAnnu) {
889 double value = cpl_image_get(aImage, i, j, &err);
893 cpl_vector_set(vbg, nbg, value);
894 cpl_vector_set(vbe, nbg, cpl_image_get(aImErr, i, j, &err));
899 cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
900 cpl_vector_delete(vbg);
901 cpl_vector_delete(vbe);
904 cpl_vector_set_size(vbg, nbg);
905 cpl_vector_set_size(vbe, nbg);
906 cpl_matrix *pos = cpl_matrix_new(1, nbg);
911 unsigned int nrej = nbg - cpl_vector_get_size(vbg);
913 nbg = cpl_vector_get_size(vbg);
915 double smean = cpl_polynomial_get_coeff(fit, &pows),
917 cpl_polynomial_delete(fit);
918 cpl_matrix_delete(pos);
919 cpl_vector_delete(vbg);
920 cpl_vector_delete(vbe);
922 cpl_msg_debug(__func__,
"sky: %d pixels (%d rejected), %f +/- %f; found %d "
923 "bad pixels inside aperture", nbg, nrej, smean, sstdev, nbad);
925 if (nbad > aMaxBad) {
926 cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
932 cpl_detector_interpolate_rejected(aImage);
933 cpl_detector_interpolate_rejected(aImErr);
939 unsigned int nobj = 0;
940 for (i = x1; i <= x2; i++) {
942 for (j = y1; j <= y2; j++) {
943 double r = sqrt(pow(aX - i, 2) + pow(aY - j, 2));
948 double value = cpl_image_get(aImage, i, j, &err),
949 error = cpl_image_get(aImErr, i, j, &err);
955 flux -= smean * nobj;
960 ferr = sqrt(ferr + nobj * sstdev*sstdev * (1. + (
double)nobj / nbg));
962 cpl_msg_debug(__func__,
"flux: %d pixels (%d interpolated), %f +/- %f",
963 nobj, nbad, flux, ferr);
1023 cpl_ensure_code(aPixtable && aFluxObj, CPL_ERROR_NULL_INPUT);
1026 cpl_msg_info(__func__,
"Gaussian profile fits for flux integration");
1029 cpl_msg_info(__func__,
"Moffat profile fits for flux integration");
1032 cpl_msg_info(__func__,
"Circular flux integration");
1035 cpl_msg_info(__func__,
"Simple square window flux integration");
1038 return cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
1041 if (getenv(
"MUSE_DEBUG_FLUX") && atoi(getenv(
"MUSE_DEBUG_FLUX")) > 2) {
1042 const char *fn =
"flux__pixtable.fits";
1043 cpl_msg_info(__func__,
"Saving pixel table as \"%s\"", fn);
1053 params->dlambda = kMuseSpectralSamplingA;
1058 aFluxObj->
cube = cube;
1061 if (getenv(
"MUSE_DEBUG_FLUX") && atoi(getenv(
"MUSE_DEBUG_FLUX")) >= 2) {
1062 const char *fn =
"flux__cube.fits";
1063 cpl_msg_info(__func__,
"Saving cube as \"%s\"", fn);
1066 int nplane = cpl_imagelist_get_size(cube->
data) / 2;
1067 cpl_image *cim = cpl_imagelist_get(cube->
data, nplane);
1069 double dsigmas[] = { 50., 30., 10., 8., 6., 5. };
1070 cpl_vector *vsigmas = cpl_vector_wrap(
sizeof(dsigmas) /
sizeof(
double),
1072 cpl_size isigma = -1;
1073 cpl_apertures *apertures = cpl_apertures_extract(cim, vsigmas, &isigma);
1074 int napertures = apertures ? cpl_apertures_get_size(apertures) : 0;
1075 if (napertures < 1) {
1076 cpl_vector_unwrap(vsigmas);
1077 cpl_apertures_delete(apertures);
1078 return cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
1080 cpl_msg_debug(__func__,
"The %.1f sigma threshold was used to find %d source%s",
1081 cpl_vector_get(vsigmas, isigma), napertures, napertures == 1 ?
"" :
"s");
1082 cpl_vector_unwrap(vsigmas);
1084 cpl_apertures_dump(apertures, stdout);
1089 int nlambda = cpl_imagelist_get_size(cube->
data);
1091 intimage->
data = cpl_image_new(nlambda, napertures, CPL_TYPE_FLOAT);
1092 intimage->
dq = cpl_image_new(nlambda, napertures, CPL_TYPE_INT);
1093 intimage->
stat = cpl_image_new(nlambda, napertures, CPL_TYPE_FLOAT);
1095 intimage->
header = cpl_propertylist_new();
1096 cpl_propertylist_append_double(intimage->
header,
"CRVAL1",
1097 cpl_propertylist_get_double(cube->
header,
1099 cpl_propertylist_append_double(intimage->
header,
"CRPIX1",
1100 cpl_propertylist_get_double(cube->
header,
1102 cpl_propertylist_append_double(intimage->
header,
"CD1_1",
1103 cpl_propertylist_get_double(cube->
header,
1105 cpl_propertylist_append_string(intimage->
header,
"CTYPE1",
1106 cpl_propertylist_get_string(cube->
header,
1108 cpl_propertylist_append_string(intimage->
header,
"CUNIT1",
1109 cpl_propertylist_get_string(cube->
header,
1112 cpl_propertylist_append_double(intimage->
header,
"CRVAL2", 1.);
1113 cpl_propertylist_append_double(intimage->
header,
"CRPIX2", 1.);
1114 cpl_propertylist_append_double(intimage->
header,
"CD2_2", 1.);
1115 cpl_propertylist_append_string(intimage->
header,
"CTYPE2",
"PIXEL");
1116 cpl_propertylist_append_string(intimage->
header,
"CUNIT2",
"pixel");
1117 cpl_propertylist_append_double(intimage->
header,
"CD1_2", 0.);
1118 cpl_propertylist_append_double(intimage->
header,
"CD2_1", 0.);
1120 cpl_propertylist_append_string(intimage->
header,
"DATE-OBS",
1121 cpl_propertylist_get_string(cube->
header,
1123 cpl_propertylist_append_string(intimage->
header,
"BUNIT",
1124 cpl_propertylist_get_string(cube->
header,
1126 cpl_propertylist_append_double(intimage->
header,
"EXPTIME",
1128 cpl_propertylist_append_string(intimage->
header,
"ESO INS MODE",
1129 cpl_propertylist_get_string(cube->
header,
1133 cpl_errorstate ps = cpl_errorstate_get();
1137 fwhm /= (kMuseSpaxelSizeX_WFM + kMuseSpaxelSizeY_WFM) / 2.;
1139 fwhm /= (kMuseSpaxelSizeX_NFM + kMuseSpaxelSizeY_NFM) / 2.;
1141 if (!cpl_errorstate_is_equal(ps)) {
1142 double xc = cpl_apertures_get_centroid_x(apertures, 1),
1143 yc = cpl_apertures_get_centroid_y(apertures, 1),
1145 cpl_image_get_fwhm(cim, lround(xc), lround(yc), &xfwhm, &yfwhm);
1146 if (xfwhm > 0. && yfwhm > 0.) {
1147 fwhm = (xfwhm + yfwhm) / 2.;
1148 }
else if (xfwhm > 0.) {
1150 }
else if (yfwhm > 0.) {
1155 cpl_errorstate_set(ps);
1156 cpl_msg_debug(__func__,
"Using roughly estimated reference FWHM (%.3f pix) "
1157 "instead of DIMM seeing", fwhm);
1159 cpl_msg_debug(__func__,
"Using DIMM seeing of %.3f pix for reference FWHM",
1165 cpl_image *sizes = cpl_image_new(nlambda, napertures, CPL_TYPE_DOUBLE);
1166 double *psizes = cpl_image_get_data_double(sizes);
1168 float *data = cpl_image_get_data_float(intimage->
data),
1169 *stat = cpl_image_get_data_float(intimage->
stat);
1170 int *dq = cpl_image_get_data_int(intimage->
dq);
1172 int l, ngood = 0, nillegal = 0, nbadbg = 0;
1173 #pragma omp parallel for default(none) \
1174 shared(apertures, aProfile, cube, data, dq, fwhm, napertures, nbadbg,\
1175 ngood, nillegal, nlambda, psizes, stat)
1176 for (l = 0; l < nlambda; l++) {
1177 cpl_image *plane = cpl_imagelist_get(cube->
data, l),
1178 *pldq = cpl_imagelist_get(cube->
dq, l),
1179 *plerr = cpl_image_duplicate(cpl_imagelist_get(cube->
stat, l));
1181 cpl_stats *stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1182 cpl_msg_debug(__func__,
"lambda = %d/%f %s", l + 1,
1183 (l + 1 - cpl_propertylist_get_double(cube->
header,
"CRPIX3"))
1184 * cpl_propertylist_get_double(cube->
header,
"CD3_3")
1185 + cpl_propertylist_get_double(cube->
header,
"CRVAL3"),
1186 cpl_propertylist_get_string(cube->
header,
"CUNIT3"));
1187 cpl_msg_debug(__func__,
"variance: %g...%g...%g", cpl_stats_get_min(stats),
1188 cpl_stats_get_mean(stats), cpl_stats_get_max(stats));
1189 cpl_stats_delete(stats);
1194 stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1195 cpl_msg_debug(__func__,
"cut variance: %g...%g...%g (%"CPL_SIZE_FORMAT
" bad"
1196 " pixel)", cpl_stats_get_min(stats), cpl_stats_get_mean(stats),
1197 cpl_stats_get_max(stats), cpl_image_count_rejected(plane));
1198 cpl_stats_delete(stats);
1201 cpl_image_power(plerr, 0.5);
1203 stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
1204 cpl_msg_debug(__func__,
"errors: %g...%g...%g", cpl_stats_get_min(stats),
1205 cpl_stats_get_mean(stats), cpl_stats_get_max(stats));
1206 cpl_stats_delete(stats);
1208 cpl_errorstate state = cpl_errorstate_get();
1210 for (n = 1; n <= napertures; n++) {
1212 double xc = cpl_apertures_get_centroid_x(apertures, n),
1213 yc = cpl_apertures_get_centroid_y(apertures, n),
1214 size = sqrt(cpl_apertures_get_npix(apertures, n)),
1216 cpl_errorstate prestate = cpl_errorstate_get();
1217 cpl_image_get_fwhm(plane, lround(xc), lround(yc), &xfwhm, &yfwhm);
1218 if (xfwhm < 0 || yfwhm < 0) {
1219 data[l + (n-1) * nlambda] = 0.;
1220 stat[l + (n-1) * nlambda] = FLT_MAX;
1221 cpl_errorstate_set(prestate);
1225 double halfsize = fmax(1.5 * (xfwhm + yfwhm), 3. * fwhm);
1226 if (halfsize < size / 2) {
1227 halfsize = size / 2;
1229 psizes[l + (n-1) * nlambda] = halfsize;
1231 cpl_msg_debug(__func__,
"%.2f,%.2f FWHM %.2f %.2f size %.2f --> %.2f",
1232 xc, yc, xfwhm, yfwhm, size, halfsize * 2.);
1237 data[l + (n-1) * nlambda] = muse_flux_image_gaussian(plane, plerr, xc, yc,
1239 &stat[l + (n-1) * nlambda]);
1242 data[l + (n-1) * nlambda] = muse_flux_image_moffat(plane, plerr, xc, yc,
1244 &stat[l + (n-1) * nlambda]);
1249 double radius = 4./3. * halfsize,
1250 rannu = radius * 5. / 4.;
1251 psizes[l + (n-1) * nlambda] = radius;
1252 data[l + (n-1) * nlambda] = muse_flux_image_circle(plane, plerr, xc, yc,
1253 radius, rannu, 10, 10,
1254 &stat[l + (n-1) * nlambda]);
1258 data[l + (n-1) * nlambda] = muse_flux_image_square(plane, plerr, xc, yc,
1260 &stat[l + (n-1) * nlambda]);
1262 if (data[l + (n-1) * nlambda] < 0 || !isfinite(data[l + (n-1) * nlambda])) {
1263 data[l + (n-1) * nlambda] = 0.;
1264 dq[l + (n-1) * nlambda] = EURO3D_MISSDATA;
1265 stat[l + (n-1) * nlambda] = FLT_MAX;
1271 if (!cpl_errorstate_is_equal(state)) {
1272 if (cpl_error_get_code() == CPL_ERROR_ILLEGAL_INPUT) {
1273 cpl_errorstate_set(state);
1276 }
else if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
1277 cpl_errorstate_set(state);
1286 cpl_image_delete(plerr);
1290 cpl_image_reject_value(sizes, CPL_VALUE_ZERO);
1292 for (n = 1; n <= napertures; n++) {
1294 cpl_msg_info(__func__,
"Radiuses used for circular flux integration for "
1295 "source %d: %f +/- %f (%f) %f..%f", n,
1296 cpl_image_get_mean_window(sizes, 1, n, nlambda, n),
1297 cpl_image_get_stdev_window(sizes, 1, n, nlambda, n),
1298 cpl_image_get_median_window(sizes, 1, n, nlambda, n),
1299 cpl_image_get_min_window(sizes, 1, n, nlambda, n),
1300 cpl_image_get_max_window(sizes, 1, n, nlambda, n));
1302 cpl_msg_info(__func__,
"Half-sizes used for flux integration for source "
1303 "%d: %f +/- %f (%f) %f..%f", n,
1304 cpl_image_get_mean_window(sizes, 1, n, nlambda, n),
1305 cpl_image_get_stdev_window(sizes, 1, n, nlambda, n),
1306 cpl_image_get_median_window(sizes, 1, n, nlambda, n),
1307 cpl_image_get_min_window(sizes, 1, n, nlambda, n),
1308 cpl_image_get_max_window(sizes, 1, n, nlambda, n));
1312 cpl_image_save(sizes,
"sizes.fits", CPL_TYPE_UNSPECIFIED, NULL, CPL_IO_CREATE);
1314 cpl_image_delete(sizes);
1317 cpl_propertylist_append_int(intimage->
header, MUSE_HDR_FLUX_NOBJ,
1323 cpl_propertylist_delete(wcs1);
1325 double crpix1 = cpl_propertylist_get_double(cube->
header,
"CRPIX1")
1326 + (1. + cpl_image_get_size_x(cim)) / 2.,
1327 crpix2 = cpl_propertylist_get_double(cube->
header,
"CRPIX2")
1328 + (1. + cpl_image_get_size_y(cim)) / 2.;
1329 cpl_propertylist_update_double(wcs,
"CRPIX1", crpix1);
1330 cpl_propertylist_update_double(wcs,
"CRPIX2", crpix2);
1333 for (n = 1; n <= napertures; n++) {
1335 double xc = cpl_apertures_get_centroid_x(apertures, n),
1336 yc = cpl_apertures_get_centroid_y(apertures, n),
1339 double flux = cpl_image_get_flux_window(intimage->
data, 1, n, nlambda, n);
1340 cpl_msg_debug(__func__,
"Source %02d: %.3f,%.3f pix, %f,%f deg, flux %e %s",
1341 n, xc, yc, ra, dec, flux,
1342 cpl_propertylist_get_string(intimage->
header,
"BUNIT"));
1343 char kw[KEYWORD_LENGTH];
1344 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_X, n);
1345 cpl_propertylist_append_float(intimage->
header, kw, xc);
1346 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_Y, n);
1347 cpl_propertylist_append_float(intimage->
header, kw, yc);
1348 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_RA, n);
1349 cpl_propertylist_append_double(intimage->
header, kw, ra);
1350 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_DEC, n);
1351 cpl_propertylist_append_double(intimage->
header, kw, dec);
1352 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_FLUX, n);
1353 cpl_propertylist_append_double(intimage->
header, kw, flux);
1355 cpl_propertylist_delete(wcs);
1356 cpl_apertures_delete(apertures);
1359 if (nillegal > 0 || nbadbg > 0) {
1360 cpl_msg_warning(__func__,
"Successful fits in %d wavelength planes, but "
1361 "encountered %d \"Illegal input\" errors and %d bad "
1362 "background determinations", ngood, nillegal, nbadbg);
1364 cpl_msg_info(__func__,
"Successful fits in %d wavelength planes", ngood);
1366 return CPL_ERROR_NONE;
1376 static const double kTelluricBands[][4] = {
1377 { 6273., 6320., 6213., 6380. },
1378 { 6864., 6967., 6750., 7130. },
1379 { 7164., 7325., 7070., 7580. },
1380 { 7590., 7700., 7470., 7830. },
1381 { 8131., 8345., 7900., 8600. },
1382 { 8952., 9028., 8850., 9082. },
1383 { 9274., 9770., 9080., 9263. },
1384 { -1., -1., -1., -1. }
1393 {
"lmin", CPL_TYPE_DOUBLE,
"Angstrom",
"%8.3f",
1394 "lower limit of the telluric region", CPL_TRUE },
1395 {
"lmax", CPL_TYPE_DOUBLE,
"Angstrom",
"%8.3f",
1396 "upper limit of the telluric region", CPL_TRUE },
1397 {
"bgmin", CPL_TYPE_DOUBLE,
"Angstrom",
"%8.3f",
1398 "lower limit of the background region", CPL_TRUE },
1399 {
"bgmax", CPL_TYPE_DOUBLE,
"Angstrom",
"%8.3f",
1400 "upper limit of the background region", CPL_TRUE },
1401 { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
1414 const cpl_table *aTellBands)
1417 cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
1422 == CPL_ERROR_NONE) {
1423 cpl_msg_debug(__func__,
"using given table for telluric bands");
1424 aFluxObj->
tellbands = cpl_table_duplicate(aTellBands);
1428 unsigned int ntell =
sizeof(kTelluricBands) /
sizeof(kTelluricBands[0]) - 1;
1429 cpl_msg_debug(__func__,
"using builtin regions for telluric bands (%u "
1434 for (k = 0; k < ntell; k++) {
1435 cpl_table_set_double(tb,
"lmin", k, kTelluricBands[k][0]);
1436 cpl_table_set_double(tb,
"lmax", k, kTelluricBands[k][1]);
1437 cpl_table_set_double(tb,
"bgmin", k, kTelluricBands[k][2]);
1438 cpl_table_set_double(tb,
"bgmax", k, kTelluricBands[k][3]);
1440 if (getenv(
"MUSE_DEBUG_FLUX") && atoi(getenv(
"MUSE_DEBUG_FLUX")) >= 2) {
1441 const char *fn =
"flux__tellregions.fits";
1442 cpl_msg_info(__func__,
"Saving telluric bands table as \"%s\"", fn);
1443 cpl_table_save(tb, NULL, NULL, fn, CPL_IO_CREATE);
1466 char *dodebug = getenv(
"MUSE_DEBUG_FLUX");
1467 if (!dodebug || (dodebug && atoi(dodebug) <= 0)) {
1470 char *fn = cpl_sprintf(
"flux__sens_%s.ascii", aName);
1471 FILE *fp = fopen(fn,
"w");
1476 cpl_msg_debug(__func__,
"Written %"CPL_SIZE_FORMAT
" datapoints to \"%s\"",
1503 unsigned int aStar,
const cpl_table *aReference,
1504 double aAirmass,
const cpl_table *aExtinct)
1506 double crval = cpl_propertylist_get_double(aFluxObj->
intimage->
header,
"CRVAL1"),
1507 cdelt = cpl_propertylist_get_double(aFluxObj->
intimage->
header,
"CD1_1"),
1508 crpix = cpl_propertylist_get_double(aFluxObj->
intimage->
header,
"CRPIX1"),
1510 int nlambda = cpl_image_get_size_x(aFluxObj->
intimage->
data);
1514 cpl_table_new_column(sensitivity,
"lambda", CPL_TYPE_DOUBLE);
1515 cpl_table_new_column(sensitivity,
"sens", CPL_TYPE_DOUBLE);
1516 cpl_table_new_column(sensitivity,
"serr", CPL_TYPE_DOUBLE);
1517 cpl_table_new_column(sensitivity,
"dq", CPL_TYPE_INT);
1518 cpl_table_set_column_format(sensitivity,
"dq",
"%u");
1519 float *data = cpl_image_get_data_float(aFluxObj->
intimage->
data),
1520 *stat = cpl_image_get_data_float(aFluxObj->
intimage->
stat);
1522 for (l = 0; l < nlambda; l++) {
1523 if (data[l + aStar*nlambda] <= 0. ||
1524 stat[l + aStar*nlambda] <= 0. ||
1525 stat[l + aStar*nlambda] == FLT_MAX) {
1529 double lambda = crval + cdelt * (l + 1 - crpix),
1531 extinct = !aExtinct ? 0.
1534 cpl_errorstate prestate = cpl_errorstate_get();
1539 if (!cpl_errorstate_is_equal(prestate)) {
1540 cpl_errorstate_set(prestate);
1546 double c = 2.5 * log10(data[l + aStar*nlambda]
1547 / exptime / cdelt / ref)
1548 + aAirmass * extinct,
1549 cerr = sqrt(pow(referr / ref, 2) + stat[l + aStar*nlambda]
1550 / pow(data[l + aStar*nlambda], 2))
1551 * 2.5 / CPL_MATH_LN10;
1552 cpl_table_set_double(sensitivity,
"lambda", idx, lambda);
1553 cpl_table_set_double(sensitivity,
"sens", idx, c);
1554 cpl_table_set_double(sensitivity,
"serr", idx, cerr);
1555 cpl_table_set_int(sensitivity,
"dq", idx, EURO3D_GOODPIXEL);
1559 cpl_table_set_size(sensitivity, idx);
1580 cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
1584 cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
1592 > MUSE_MODE_WFM_NONAO_X;
1595 int irow, nfluxes = cpl_table_get_nrow(tsens);
1596 for (irow = 0; irow < nfluxes; irow++) {
1597 double lambda = cpl_table_get_double(tsens,
"lambda", irow, NULL);
1598 unsigned int dq = EURO3D_GOODPIXEL,
1599 k, nk = cpl_table_get_nrow(tb);
1600 for (k = 0; k < nk; k++) {
1601 double lmin = cpl_table_get_double(tb,
"lmin", k, NULL),
1602 lmax = cpl_table_get_double(tb,
"lmax", k, NULL);
1603 if (lambda >= lmin && lambda <= lmax) {
1604 dq |= EURO3D_TELLURIC;
1607 if (isnominal && lambda < kMuseNominalCutoff) {
1608 dq |= EURO3D_OUTSDRANGE;
1610 cpl_table_set_int(tsens,
"dq", irow, dq);
1638 static cpl_polynomial *
1640 double aLambda1,
double aLambda2,
1641 unsigned int aOrder,
double aRSigma,
double *aRMSE)
1644 cpl_table_select_all(tsens);
1645 cpl_table_and_selected_int(tsens,
"dq", CPL_NOT_EQUAL_TO, EURO3D_GOODPIXEL);
1646 cpl_table_and_selected_int(tsens,
"dq", CPL_NOT_EQUAL_TO, EURO3D_TELLCOR);
1647 cpl_table_or_selected_double(tsens,
"lambda", CPL_LESS_THAN, aLambda1);
1648 cpl_table_or_selected_double(tsens,
"lambda", CPL_GREATER_THAN, aLambda2);
1650 cpl_table *tunwanted = cpl_table_extract_selected(tsens);
1651 cpl_table_erase_selected(tsens);
1652 muse_flux_response_dump_sensitivity(aFluxObj,
"fitinput");
1656 int nrow = cpl_table_get_nrow(tsens);
1657 cpl_matrix *lambdas = cpl_matrix_new(1, nrow);
1658 cpl_vector *sens = cpl_vector_new(nrow),
1659 *serr = cpl_vector_new(nrow);
1660 memcpy(cpl_matrix_get_data(lambdas),
1661 cpl_table_get_data_double_const(tsens,
"lambda"), nrow*
sizeof(
double));
1662 memcpy(cpl_vector_get_data(sens),
1663 cpl_table_get_data_double_const(tsens,
"sens"), nrow*
sizeof(
double));
1664 memcpy(cpl_vector_get_data(serr),
1665 cpl_table_get_data_double_const(tsens,
"serr"), nrow*
sizeof(
double));
1670 tsens, aOrder, aRSigma,
1672 int nout = cpl_vector_get_size(sens);
1674 cpl_msg_debug(__func__,
"transferred %d entries (%.3f...%.3f) for the "
1675 "order %u fit, %d entries are left, RMS %f", nrow, aLambda1,
1676 aLambda2, aOrder, nout, sqrt(mse));
1678 cpl_matrix_delete(lambdas);
1679 cpl_vector_delete(sens);
1680 cpl_vector_delete(serr);
1682 *aRMSE = mse / (nout - aOrder - 1);
1686 cpl_table_insert(tsens, tunwanted, nout);
1687 cpl_table_delete(tunwanted);
1716 cpl_table_new_column(tsens,
"sens_orig", CPL_TYPE_DOUBLE);
1717 cpl_table_new_column(tsens,
"serr_orig", CPL_TYPE_DOUBLE);
1718 cpl_table_new_column(tsens,
"tellcor", CPL_TYPE_DOUBLE);
1719 unsigned int k, nk = cpl_table_get_nrow(tb);
1720 for (k = 0; k < nk; k++) {
1721 double lmin = cpl_table_get_double(tb,
"lmin", k, NULL),
1722 lmax = cpl_table_get_double(tb,
"lmax", k, NULL),
1723 bgmin = cpl_table_get_double(tb,
"bgmin", k, NULL),
1724 bgmax = cpl_table_get_double(tb,
"bgmax", k, NULL),
1725 datamin = cpl_table_get_column_min(tsens,
"lambda"),
1726 datamax = cpl_table_get_column_max(tsens,
"lambda");
1727 cpl_boolean extrapolate = CPL_FALSE;
1728 if (bgmax < lmax || datamax < lmax) {
1729 extrapolate = CPL_TRUE;
1731 if (bgmin > lmin || datamin > lmin) {
1732 extrapolate = CPL_TRUE;
1734 if (datamin > lmax || datamax < lmin) {
1736 cpl_msg_warning(__func__,
"Telluric region %u (range %.2f...%.2f, "
1737 "reference region %.2f...%.2f) outside data range "
1738 "(%.2f..%.2f)!", k + 1, lmin, lmax, bgmin, bgmax,
1743 unsigned int order = extrapolate ? 1 : 2;
1747 cpl_errorstate state = cpl_errorstate_get();
1748 cpl_polynomial *fit = muse_flux_response_fit(aFluxObj, bgmin, bgmax,
1750 if (!cpl_errorstate_is_equal(state) || !fit) {
1751 cpl_msg_warning(__func__,
"Telluric region %u (range %.2f...%.2f, "
1752 "reference region %.2f...%.2f) could not be fitted!",
1753 k + 1, lmin, lmax, bgmin, bgmax);
1754 cpl_errorstate_set(state);
1755 cpl_polynomial_delete(fit);
1758 cpl_msg_debug(__func__,
"Telluric region %u: %.2f...%.2f, reference region "
1759 "%.2f...%.2f", k + 1, lmin, lmax, bgmin, bgmax);
1761 cpl_polynomial_dump(fit, stdout);
1764 int irow, nrow = cpl_table_get_nrow(tsens);
1765 for (irow = 0; irow < nrow; irow++) {
1766 double lambda = cpl_table_get_double(tsens,
"lambda", irow, NULL);
1767 if (lambda >= lmin && lambda <= lmax &&
1768 (
unsigned int)cpl_table_get_int(tsens,
"dq", irow, NULL)
1769 == EURO3D_TELLURIC) {
1770 double origval = cpl_table_get_double(tsens,
"sens", irow, NULL),
1771 origerr = cpl_table_get_double(tsens,
"serr", irow, NULL),
1772 interpval = cpl_polynomial_eval_1d(fit, lambda, NULL);
1773 cpl_table_set_int(tsens,
"dq", irow, EURO3D_TELLCOR);
1774 cpl_table_set_double(tsens,
"sens_orig", irow, origval);
1775 cpl_table_set_double(tsens,
"sens", irow, interpval);
1778 cpl_table_set_double(tsens,
"serr_orig", irow, origerr);
1779 cpl_table_set_double(tsens,
"serr", irow, sqrt(origerr*origerr + rmse));
1781 if (interpval > origval) {
1783 double ftelluric = pow(10, -0.4 * (interpval - origval));
1784 cpl_table_set(tsens,
"tellcor", irow, ftelluric);
1786 cpl_table_set_double(tsens,
"tellcor", irow, 1.);
1790 cpl_polynomial_delete(fit);
1793 cpl_table_power_column(tsens,
"tellcor", 1. / aAirmass);
1815 muse_flux_response_extrapolate(
muse_flux_object *aFluxObj,
double aDistance)
1818 cpl_propertylist *order = cpl_propertylist_new();
1819 cpl_propertylist_append_bool(order,
"lambda", CPL_FALSE);
1820 cpl_table_sort(tsens, order);
1822 int nrow = cpl_table_get_nrow(tsens);
1824 double lambda1 = cpl_table_get_double(tsens,
"lambda", 0, NULL),
1825 serr1 = cpl_table_get_double(tsens,
"serr", 0, NULL);
1826 unsigned int dq = cpl_table_get_int(tsens,
"dq", 0, NULL);
1828 while (dq != EURO3D_GOODPIXEL && dq != EURO3D_TELLCOR) {
1829 lambda1 = cpl_table_get_double(tsens,
"lambda", ++irow, NULL);
1830 serr1 = cpl_table_get_double(tsens,
"serr", irow, NULL);
1831 dq = cpl_table_get_int(tsens,
"dq", irow, NULL);
1833 double lambda2 = cpl_table_get_double(tsens,
"lambda", nrow - 1, NULL),
1834 serr2 = cpl_table_get_double(tsens,
"serr", nrow - 1, NULL);
1835 dq = cpl_table_get_int(tsens,
"dq", nrow - 1, NULL);
1837 while (dq != EURO3D_GOODPIXEL && dq != EURO3D_TELLCOR) {
1838 lambda2 = cpl_table_get_double(tsens,
"lambda", --irow, NULL);
1839 serr2 = cpl_table_get_double(tsens,
"serr", irow, NULL);
1840 dq = cpl_table_get_int(tsens,
"dq", irow, NULL);
1842 cpl_polynomial *fit1 = muse_flux_response_fit(aFluxObj, lambda1,
1843 lambda1 + aDistance, 1, 5., NULL),
1844 *fit2 = muse_flux_response_fit(aFluxObj, lambda2 - aDistance,
1845 lambda2, 1, 5., NULL);
1846 nrow = cpl_table_get_nrow(tsens);
1847 double d1 = (lambda1 - kMuseLambdaMinX) / 100.,
1848 d2 = (kMuseLambdaMaxX - lambda2) / 100.;
1849 cpl_table_set_size(tsens, nrow + 200);
1853 for (l = kMuseLambdaMinX; l <= lambda1 - d1; l += d1) {
1854 double sens = cpl_polynomial_eval_1d(fit1, l, NULL),
1856 serr = 2. * (lambda1 - l) / 50. * serr1 + serr1;
1858 cpl_table_set_invalid(tsens,
"lambda", irow++);
1859 cpl_msg_debug(__func__,
"invalid blueward extrapolation: %.3f %f +/- %f",
1863 cpl_table_set_double(tsens,
"lambda", irow, l);
1864 cpl_table_set_double(tsens,
"sens", irow, sens);
1865 cpl_table_set_double(tsens,
"serr", irow, serr);
1866 cpl_table_set_int(tsens,
"dq", irow++, (
int)EURO3D_OUTSDRANGE);
1868 cpl_msg_debug(__func__,
"Extrapolated blue end: %.1f...%.1f Angstrom (using"
1869 " data from %.1f...%.1f Angstrom)", kMuseLambdaMinX,
1870 lambda1 - d1, lambda1, lambda1 + aDistance);
1873 for (l = lambda2 + d2; fit2 && l <= kMuseLambdaMaxX; l += d2) {
1874 double sens = cpl_polynomial_eval_1d(fit2, l, NULL),
1875 serr = 2. * (l - lambda2) / 50. * serr2 + serr2;
1877 cpl_table_set_invalid(tsens,
"lambda", irow++);
1878 cpl_msg_debug(__func__,
"invalid redward extrapolation: %.3f %f +/- %f",
1882 cpl_table_set_double(tsens,
"lambda", irow, l);
1883 cpl_table_set_double(tsens,
"sens", irow, sens);
1884 cpl_table_set_double(tsens,
"serr", irow, serr);
1885 cpl_table_set_int(tsens,
"dq", irow++, (
int)EURO3D_OUTSDRANGE);
1887 cpl_msg_debug(__func__,
"Extrapolated red end: %.1f...%.1f Angstrom (using "
1888 "data from %.1f...%.1f Angstrom)", lambda2 + d2,
1889 kMuseLambdaMaxX, lambda2 - aDistance, lambda2);
1892 cpl_polynomial_dump(fit1, stdout);
1893 cpl_polynomial_dump(fit2, stdout);
1896 cpl_polynomial_delete(fit1);
1897 cpl_polynomial_delete(fit2);
1899 cpl_table_select_all(tsens);
1900 cpl_table_and_selected_invalid(tsens,
"sens");
1901 cpl_table_erase_selected(tsens);
1903 cpl_table_sort(tsens, order);
1904 cpl_propertylist_delete(order);
1963 const cpl_table *aReference,
1964 const cpl_table *aTellBands,
1965 const cpl_table *aExtinct)
1967 cpl_ensure_code(aFluxObj && aFluxObj->
intimage && aReference,
1968 CPL_ERROR_NULL_INPUT);
1969 cpl_ensure_code(aAirmass >= 1., CPL_ERROR_ILLEGAL_INPUT);
1971 cpl_msg_warning(__func__,
"Extinction table not given!");
1974 (!isfinite(aFluxObj->
raref) || !isfinite(aFluxObj->
decref))) {
1975 cpl_msg_warning(__func__,
"Reference position %f,%f contains infinite "
1976 "values, using flux to select star!", aFluxObj->
raref,
1980 muse_flux_response_set_telluric_bands(aFluxObj, aTellBands);
1982 int nobjects = cpl_image_get_size_y(aFluxObj->
intimage->
data);
1983 const char *bunit = cpl_propertylist_get_string(aFluxObj->
intimage->
header,
1986 double flux = 0., dmin = DBL_MAX;
1987 int n, nstar = 1, nstardist = 1;
1988 for (n = 1; n <= nobjects; n++) {
1989 char kw[KEYWORD_LENGTH];
1990 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_RA, n);
1991 double ra = cpl_propertylist_get_double(aFluxObj->
intimage->
header, kw);
1992 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_DEC, n);
1993 double dec = cpl_propertylist_get_double(aFluxObj->
intimage->
header, kw),
1996 cpl_msg_debug(__func__,
"distance(%d) = %f arcsec", n, dthis * 3600.);
1997 if (fabs(dthis) < dmin) {
2001 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_FLUX, n);
2002 double this = cpl_propertylist_get_double(aFluxObj->
intimage->
header, kw);
2003 cpl_msg_debug(__func__,
"flux(%d) = %e %s", n,
this, bunit);
2010 char *outstring = NULL;
2012 outstring = cpl_sprintf(
"Selected the brightest star (%d of %d; %.3e %s)"
2013 " as reference source", nstar, nobjects, flux, bunit);
2016 outstring = cpl_sprintf(
"Selected the nearest star (%d of %d; %.2f arcsec) "
2017 "as reference source", nstar, nobjects, dmin*3600.);
2018 nselected = nstardist;
2020 cpl_msg_info(__func__,
"%s", outstring);
2021 cpl_free(outstring);
2023 muse_flux_response_sensitivity(aFluxObj,
2024 nselected - 1, aReference,
2025 aAirmass, aExtinct);
2026 muse_flux_response_dump_sensitivity(aFluxObj,
"initial");
2028 muse_flux_response_mark_questionable(aFluxObj);
2029 muse_flux_response_dump_sensitivity(aFluxObj,
"intermediate");
2032 cpl_table_and_selected_int(aFluxObj->
sensitivity,
"dq", CPL_EQUAL_TO,
2033 (
int)EURO3D_OUTSDRANGE);
2035 muse_flux_response_dump_sensitivity(aFluxObj,
"intercut");
2039 muse_flux_response_telluric(aFluxObj, aAirmass);
2040 muse_flux_response_dump_sensitivity(aFluxObj,
"interpolated");
2042 muse_flux_response_extrapolate(aFluxObj, 150.);
2043 muse_flux_response_dump_sensitivity(aFluxObj,
"extrapolated");
2045 return CPL_ERROR_NONE;
2059 {
"lambda", CPL_TYPE_DOUBLE,
"Angstrom",
"%7.2f",
"wavelength", CPL_TRUE },
2060 {
"response", CPL_TYPE_DOUBLE,
2061 "2.5*log10((count/s/Angstrom)/(erg/s/cm**2/Angstrom))",
"%.4e",
2062 "instrument response derived from standard star", CPL_TRUE },
2063 {
"resperr", CPL_TYPE_DOUBLE,
2064 "2.5*log10((count/s/Angstrom)/(erg/s/cm**2/Angstrom))",
"%.4e",
2065 "instrument response error derived from standard star", CPL_TRUE },
2066 { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
2080 {
"lambda", CPL_TYPE_DOUBLE,
"Angstrom",
"%7.2f",
"wavelength", CPL_TRUE },
2081 {
"ftelluric", CPL_TYPE_DOUBLE,
"",
"%.5f",
2082 "the telluric correction factor, normalized to an airmass of 1", CPL_TRUE },
2083 {
"ftellerr", CPL_TYPE_DOUBLE,
"",
"%.5f",
2084 "the error of the telluric correction factor", CPL_TRUE },
2085 { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
2110 double aLambdaMin,
double aLambdaMax,
2111 cpl_boolean aAverage)
2114 cpl_table_duplicate_column(aResp,
"sens", aResp,
"response");
2115 cpl_table_duplicate_column(aResp,
"serr", aResp,
"resperr");
2118 cpl_table_select_all(aResp);
2119 cpl_table_and_selected_double(aResp,
"lambda", CPL_NOT_LESS_THAN, aLambdaMin);
2120 cpl_table_and_selected_double(aResp,
"lambda", CPL_NOT_GREATER_THAN, aLambdaMax);
2122 cpl_boolean sym = cpl_table_count_selected(aResp) < cpl_table_get_nrow(aResp);
2123 cpl_msg_debug(__func__,
"%s smoothing response +/- %.3f Angstrom between %.3f "
2124 "and %.3f Angstrom", sym ?
"symmetrical" :
"", aHalfwidth,
2125 aLambdaMin, aLambdaMax);
2128 int i, n = cpl_table_get_nrow(aResp);
2129 for (i = 0; i < n; i++) {
2130 if (!cpl_table_is_selected(aResp, i)) {
2133 double lambda = cpl_table_get_double(aResp,
"lambda", i, NULL);
2134 int j = i, j1 = i, j2 = i;
2136 while (--j > 0 && cpl_table_is_selected(aResp, j) &&
2137 lambda - cpl_table_get_double(aResp,
"lambda", j, NULL) <= aHalfwidth) {
2141 while (++j < n && cpl_table_is_selected(aResp, j) &&
2142 cpl_table_get_double(aResp,
"lambda", j, NULL) - lambda <= aHalfwidth) {
2155 double *sens = cpl_table_get_data_double(aResp,
"sens"),
2156 *serr = cpl_table_get_data_double(aResp,
"serr");
2157 cpl_vector *v = cpl_vector_wrap(j2 - j1 + 1, sens + j1),
2158 *ve = cpl_vector_wrap(j2 - j1 + 1, serr + j1);
2161 double mean = cpl_vector_get_mean(v),
2162 stdev = j2 == j1 ? 0. : cpl_vector_get_stdev(v),
2163 rerr = cpl_table_get_double(aResp,
"resperr", i, NULL);
2164 cpl_table_set_double(aResp,
"response", i, mean);
2165 cpl_table_set_double(aResp,
"resperr", i, sqrt(rerr*rerr + stdev*stdev));
2168 double median = cpl_vector_get_median_const(v),
2170 mederr = cpl_vector_get_median_const(ve);
2172 mdev = cpl_table_get_double(aResp,
"serr", j1, NULL);
2174 if (mdev < mederr) {
2178 cpl_msg_debug(__func__,
"%d %.3f %d...%d --> %f +/- %f", i, lambda, j1, j2,
2181 cpl_table_set_double(aResp,
"response", i, median);
2182 cpl_table_set_double(aResp,
"resperr", i, mdev);
2184 cpl_vector_unwrap(v);
2185 cpl_vector_unwrap(ve);
2189 cpl_table_erase_column(aResp,
"sens");
2190 cpl_table_erase_column(aResp,
"serr");
2212 muse_flux_get_response_table_collect_points(
const cpl_table *aTable,
2213 double aLambda,
double aLDist,
2214 cpl_matrix *aPos, cpl_vector *aVal,
2217 unsigned int np = 0;
2218 int irow, nrow = cpl_table_get_nrow(aTable);
2219 for (irow = 0; irow < nrow; irow++) {
2220 double lambda = cpl_table_get(aTable,
"lambda", irow, NULL);
2221 if (lambda < aLambda - aLDist || lambda > aLambda + aLDist) {
2224 cpl_matrix_set(aPos, 0, np, lambda);
2225 cpl_vector_set(aVal, np, cpl_table_get(aTable,
"sens", irow, NULL));
2226 cpl_vector_set(aErr, np, cpl_table_get(aTable,
"serr", irow, NULL));
2229 cpl_matrix_set_size(aPos, 1, np);
2230 cpl_vector_set_size(aVal, np);
2231 cpl_vector_set_size(aErr, np);
2262 muse_flux_get_response_table_piecewise_poly(cpl_table *aResp,
double aDMin,
2263 double aDMax,
float aRSigma)
2266 cpl_table_duplicate_column(aResp,
"sens", aResp,
"response");
2267 cpl_table_duplicate_column(aResp,
"serr", aResp,
"resperr");
2270 unsigned int npold = 0, njumps = 0;
2271 double ldistold = -1,
2273 cpl_array *jumppos = cpl_array_new(0, CPL_TYPE_DOUBLE),
2274 *jumplen = cpl_array_new(0, CPL_TYPE_DOUBLE);
2277 int irow, nrow = cpl_table_get_nrow(aResp);
2278 for (irow = 0; irow < nrow; irow++) {
2279 double lambda = cpl_table_get(aResp,
"lambda", irow, NULL);
2283 double ldist = aDMax;
2284 if ((lambda >= 5700 && lambda < 6200) || (lambda >= 6900 && lambda < 7200)) {
2288 cpl_matrix *pos = cpl_matrix_new(1, nrow);
2289 cpl_vector *val = cpl_vector_new(nrow),
2290 *err = cpl_vector_new(nrow);
2291 unsigned int np = muse_flux_get_response_table_collect_points(aResp,
2300 printf(
"%f Angstrom %u points (%u, %.3f):\n", lambda, np, npold, (
double)np / npold - 1.);
2305 if (np > 10 && fabs((
double)np / npold - 1.) > 0.1) {
2306 cpl_msg_debug(__func__,
"possible jump, changed at lambda %.3f (%u -> %u, "
2307 "%.3f -> %.3f)", lambda, npold, np, ldistold, ldist);
2308 cpl_array_set_size(jumppos, ++njumps);
2309 cpl_array_set_size(jumplen, njumps);
2310 cpl_array_set_double(jumppos, njumps - 1, (lambdaold + lambda) / 2.);
2311 cpl_array_set_double(jumplen, njumps - 1, (ldistold + ldist) / 2.);
2315 unsigned int order = np > 3 ? 3 : np - 1;
2319 NULL, order, aRSigma,
2322 if (fabs(lambda - 40861.3) < 1.) {
2323 cpl_vector *res = cpl_vector_new(cpl_vector_get_size(val));
2324 cpl_vector_fill_polynomial_fit_residual(res, val, NULL, poly, pos, NULL);
2325 double rms = sqrt(cpl_vector_product(res, res) / cpl_vector_get_size(res));
2326 cpl_msg_debug(__func__,
"lambda %f rms %f (%u/%d points)", lambda, rms,
2327 (
unsigned)cpl_vector_get_size(val), np);
2329 cpl_plot_vector(NULL, NULL, NULL, val);
2330 cpl_plot_vector(NULL, NULL, NULL, res);
2331 cpl_vector_delete(res);
2334 cpl_matrix_delete(pos);
2335 cpl_vector_delete(val);
2336 cpl_vector_delete(err);
2337 double resp = cpl_polynomial_eval_1d(poly, lambda, NULL);
2338 cpl_polynomial_delete(poly);
2339 cpl_table_set(aResp,
"lambda", irow, lambda);
2340 cpl_table_set(aResp,
"response", irow, resp);
2341 double serr = cpl_table_get(aResp,
"serr", irow, NULL);
2342 cpl_table_set(aResp,
"resperr", irow, sqrt(mse + serr*serr));
2350 cpl_table_erase_column(aResp,
"sens");
2351 cpl_table_erase_column(aResp,
"serr");
2354 printf(
"jumppos (%u):\n", njumps);
2355 cpl_array_dump(jumppos, 0, 10000, stdout);
2356 printf(
"jumplen:\n");
2357 cpl_array_dump(jumplen, 0, 10000, stdout);
2364 for (iarr = 0; iarr < njumps; iarr++) {
2365 double lambda = cpl_array_get_double(jumppos, iarr, NULL),
2366 ldist = cpl_array_get_double(jumplen, iarr, NULL) / 2;
2368 cpl_table_select_all(aResp);
2369 cpl_table_and_selected_double(aResp,
"lambda", CPL_NOT_LESS_THAN, lambda - 5.);
2370 cpl_table_and_selected_double(aResp,
"lambda", CPL_NOT_GREATER_THAN, lambda + 5.);
2371 int nsel = cpl_table_count_selected(aResp);
2373 cpl_msg_debug(__func__,
"Only %d points near jump around %.1f Angstrom, "
2374 "not doing anything", nsel, lambda);
2377 cpl_table *xresp = cpl_table_extract_selected(aResp);
2378 double stdev = cpl_table_get_column_stdev(xresp,
"response");
2379 cpl_table_dump(xresp, 0, nsel, stdout);
2381 cpl_table_delete(xresp);
2382 if (stdev < 0.001) {
2383 cpl_msg_debug(__func__,
"%d points near jump around %.1f Angstrom, stdev "
2384 "only %f, not doing anything", nsel, lambda, stdev);
2387 cpl_msg_debug(__func__,
"%d points near jump around %.1f Angstrom, stdev "
2388 "%f, erasing the region", nsel, lambda, stdev);
2392 cpl_table_select_all(aResp);
2393 cpl_table_and_selected_double(aResp,
"lambda", CPL_NOT_LESS_THAN, lambda - ldist);
2394 cpl_table_and_selected_double(aResp,
"lambda", CPL_NOT_GREATER_THAN, lambda + ldist);
2395 cpl_table_erase_selected(aResp);
2397 cpl_array_delete(jumppos);
2398 cpl_array_delete(jumplen);
2431 cpl_ensure_code(aFluxObj && aFluxObj->
sensitivity, CPL_ERROR_NULL_INPUT);
2434 int nrow = cpl_table_get_nrow(aFluxObj->
sensitivity);
2437 const double *lambdas = cpl_table_get_data_double_const(aFluxObj->
sensitivity,
2439 *sens = cpl_table_get_data_double_const(aFluxObj->
sensitivity,
2441 *serr = cpl_table_get_data_double_const(aFluxObj->
sensitivity,
2443 cpl_table_copy_data_double(resp,
"lambda", lambdas);
2444 cpl_table_copy_data_double(resp,
"response", sens);
2445 cpl_table_copy_data_double(resp,
"resperr", serr);
2449 cpl_msg_info(__func__,
"Smoothing response curve with median filter");
2453 cpl_msg_info(__func__,
"Smoothing response curve with piecewise polynomial");
2455 muse_flux_get_response_table_piecewise_poly(resp, 150., 150., 3.);
2460 cpl_msg_warning(__func__,
"NOT smoothing the response curve at all!");
2465 return CPL_ERROR_NONE;
2491 cpl_ensure_code(aFluxObj && aFluxObj->
sensitivity, CPL_ERROR_NULL_INPUT);
2493 int nrow = cpl_table_get_nrow(tsens);
2496 cpl_table_fill_column_window_double(tell,
"lambda", 0, nrow, 0);
2497 cpl_table_copy_data_double(tell,
"lambda",
2498 cpl_table_get_data_double_const(tsens,
"lambda"));
2499 cpl_table_fill_column_window_double(tell,
"ftelluric", 0, nrow, 0);
2500 cpl_table_copy_data_double(tell,
"ftelluric",
2501 cpl_table_get_data_double_const(tsens,
"tellcor"));
2503 #define TELL_MAX_ERR 0.1
2504 #define TELL_MIN_ERR 1e-4
2505 cpl_table_fill_column_window_double(tell,
"ftellerr", 0, nrow, TELL_MAX_ERR);
2509 cpl_table_duplicate_column(tell,
"tellcor", tsens,
"tellcor");
2511 cpl_table_unselect_all(tell);
2513 for (irow = 0; irow < nrow; irow++) {
2515 cpl_table_get_double(tell,
"tellcor", irow, &err);
2520 cpl_errorstate state = cpl_errorstate_get();
2521 double ftellcor = cpl_table_get_double(tell,
"tellcor", irow - 1, &err);
2522 if (!cpl_errorstate_is_equal(state)) {
2523 cpl_errorstate_set(state);
2525 if (err == 0 && ftellcor != 1.) {
2526 cpl_table_set_double(tell,
"ftelluric", irow, 1.);
2530 state = cpl_errorstate_get();
2531 ftellcor = cpl_table_get_double(tell,
"tellcor", irow + 1, &err);
2532 if (!cpl_errorstate_is_equal(state)) {
2533 cpl_errorstate_set(state);
2535 if (err == 0 && ftellcor != 1.) {
2536 cpl_table_set_double(tell,
"ftelluric", irow, 1.);
2539 cpl_table_select_row(tell, irow);
2541 cpl_table_erase_selected(tell);
2542 cpl_table_erase_column(tell,
"tellcor");
2545 nrow = cpl_table_get_nrow(tell);
2546 for (irow = 0; irow < nrow; irow++) {
2548 double dftell = 1. - cpl_table_get_double(tell,
"ftelluric", irow, &err),
2549 ftellerr = cpl_table_get_double(tell,
"ftellerr", irow, &err);
2550 if (ftellerr > dftell) {
2551 ftellerr = fmax(dftell, TELL_MIN_ERR);
2552 cpl_table_set_double(tell,
"ftellerr", irow, ftellerr);
2557 return CPL_ERROR_NONE;
2600 const cpl_table *aExtinction,
const cpl_table *aTelluric)
2602 cpl_ensure_code(aPixtable && aPixtable->
header && aResponse,
2603 CPL_ERROR_NULL_INPUT);
2604 const char *unitdata = cpl_table_get_column_unit(aPixtable->
table,
2605 MUSE_PIXTABLE_DATA);
2606 cpl_ensure_code(unitdata && !strncmp(unitdata,
"count", 6),
2607 CPL_ERROR_INCOMPATIBLE_INPUT);
2610 cpl_msg_warning(__func__,
"%s missing!", MUSE_TAG_EXTINCT_TABLE);
2614 if (exptime <= 0.) {
2615 cpl_msg_warning(__func__,
"Non-positive EXPTIME, not doing flux calibration!");
2616 return CPL_ERROR_ILLEGAL_INPUT;
2620 cpl_msg_warning(__func__,
"Airmass unknown, not doing extinction "
2621 "correction: %s", cpl_error_get_message());
2626 cpl_table *telluric = NULL;
2629 telluric = cpl_table_duplicate(aTelluric);
2632 cpl_table_power_column(telluric,
"ftelluric", -airmass);
2635 cpl_msg_info(__func__,
"Starting flux calibration (exptime=%.2fs, airmass=%.4f),"
2636 " %s telluric correction", exptime, airmass,
2637 aTelluric ?
"with" :
"without ("MUSE_TAG_STD_TELLURIC
" not given)");
2638 float *lambda = cpl_table_get_data_float(aPixtable->
table, MUSE_PIXTABLE_LAMBDA),
2639 *data = cpl_table_get_data_float(aPixtable->
table, MUSE_PIXTABLE_DATA),
2640 *stat = cpl_table_get_data_float(aPixtable->
table, MUSE_PIXTABLE_STAT);
2642 #pragma omp parallel for default(none) \
2643 shared(aExtinction, aResponse, airmass, data, exptime, lambda, nrow, \
2645 for (i = 0; i < nrow; i++) {
2647 double ddata = data[i], dstat = stat[i];
2651 double fext = pow(10., 0.4 * airmass
2656 printf(
"%f, data/stat = %f/%f -> ", fext, ddata, dstat);
2659 dstat *= (fext * fext);
2661 printf(
" --> %f/%f\n", ddata, dstat), fflush(NULL);
2668 double dlambda = kMuseSpectralSamplingA,
2675 resp = pow(10., 0.4 * resp);
2678 rerr = rerr * CPL_MATH_LN10 * resp / 2.5;
2680 printf(
"%f/%f/%f, %e/%e, data/stat = %e/%e -> ", lambda[i], dlambda, dlamerr, resp, rerr,
2683 dstat = dstat * pow((1./(resp * exptime * dlambda)), 2)
2684 + pow(ddata * rerr / (resp*resp * exptime * dlambda), 2)
2685 + pow(ddata * dlamerr / (resp * exptime * dlambda*dlambda), 2);
2686 ddata /= (resp * exptime * dlambda);
2688 printf(
"%e/%e\n", ddata, dstat), fflush(NULL);
2693 ddata *= kMuseFluxUnitFactor;
2694 dstat *= kMuseFluxStatFactor;
2698 if (lambda[i] < kTelluricBands[0][0] || !telluric) {
2706 data[i] = ddata * tell;
2707 stat[i] = tell*tell * dstat + ddata*ddata * terr*terr;
2709 cpl_table_delete(telluric);
2712 cpl_table_set_column_unit(aPixtable->
table, MUSE_PIXTABLE_DATA,
2713 kMuseFluxUnitString);
2714 cpl_table_set_column_unit(aPixtable->
table, MUSE_PIXTABLE_STAT,
2715 kMuseFluxStatString);
2721 MUSE_HDR_PT_FLUXCAL_COMMENT);
2722 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.
muse_flux_smooth_type
Type of response curve smoothing to use.
double muse_pfits_get_ra(const cpl_propertylist *aHeaders)
find out the right ascension
#define MUSE_HDR_PT_FLUXCAL
cpl_image * data
the data extension
const muse_cpltable_def muse_flux_responsetable_def[]
MUSE response table definition.
cpl_size muse_pixtable_get_nrow(const muse_pixtable *aPixtable)
get the number of rows within the pixel table
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)
Structure definition of MUSE pixel table.
Flux object to store data needed while computing the flux calibration.
double muse_astro_wavelength_vacuum_to_air(double aVac)
Compute air wavelength for a given vacuum wavelength.
cpl_error_code muse_wcs_celestial_from_pixel(cpl_propertylist *aHeader, double aX, double aY, double *aRA, double *aDEC)
Convert from pixel coordinates to celestial spherical coordinates.
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.
cpl_error_code muse_flux_response_compute(muse_flux_object *aFluxObj, muse_flux_selection_type aSelect, 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...
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.
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
double muse_pfits_get_dec(const cpl_propertylist *aHeaders)
find out the declination
muse_flux_interpolation_type
Type of table interpolation to use.
cpl_error_code muse_flux_get_response_table(muse_flux_object *aFluxObj, muse_flux_smooth_type aSmooth)
Get the table of the standard star response function.
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
double muse_pfits_get_exptime(const cpl_propertylist *aHeaders)
find out the exposure time
cpl_propertylist * muse_wcs_apply_cd(const cpl_propertylist *aExpHeader, const cpl_propertylist *aCalHeader)
Apply the CDi_j matrix of an astrometric solution to an observation.
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.
muse_flux_selection_type
Type of star selection to use.
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.
double muse_astro_angular_distance(double aRA1, double aDEC1, double aRA2, double aDEC2)
Compute angular distance in the sky between two positions.
muse_ins_mode muse_pfits_get_mode(const cpl_propertylist *aHeaders)
find out the observation mode
static void muse_flux_get_response_table_smooth(cpl_table *aResp, double aHalfwidth, double aLambdaMin, double aLambdaMax, cpl_boolean aAverage)
Get the table of the standard star response function.
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.
cpl_propertylist * muse_wcs_create_default(void)
Create FITS headers containing a default (relative) WCS.