00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #ifdef HAVE_CONFIG_H
00023 #include <config.h>
00024 #endif
00025
00026
00027
00028
00029 #include <cpl.h>
00030 #include <math.h>
00031 #include <string.h>
00032
00033 #include "muse_flux.h"
00034 #include "muse_instrument.h"
00035
00036 #include "muse_astro.h"
00037 #include "muse_cplwrappers.h"
00038 #include "muse_pfits.h"
00039 #include "muse_quality.h"
00040 #include "muse_resampling.h"
00041 #include "muse_utils.h"
00042 #include "muse_wcs.h"
00043
00044
00048
00049
00052
00066
00067 muse_flux_object *
00068 muse_flux_object_new(void)
00069 {
00070 muse_flux_object *flux = cpl_calloc(1, sizeof(muse_flux_object));
00071
00072 flux->raref = NAN;
00073 flux->decref = NAN;
00074 return flux;
00075 }
00076
00077
00087
00088 void
00089 muse_flux_object_delete(muse_flux_object *aFluxObj)
00090 {
00091 if (!aFluxObj) {
00092 return;
00093 }
00094 muse_datacube_delete(aFluxObj->cube);
00095 aFluxObj->cube = NULL;
00096 muse_image_delete(aFluxObj->intimage);
00097 aFluxObj->intimage = NULL;
00098 cpl_table_delete(aFluxObj->sensitivity);
00099 aFluxObj->sensitivity = NULL;
00100 cpl_table_delete(aFluxObj->response);
00101 aFluxObj->response = NULL;
00102 cpl_table_delete(aFluxObj->telluric);
00103 aFluxObj->telluric = NULL;
00104 cpl_table_delete(aFluxObj->tellbands);
00105 aFluxObj->tellbands = NULL;
00106 cpl_free(aFluxObj);
00107 }
00108
00109
00117
00118 static double
00119 muse_flux_reference_table_sampling(cpl_table *aTable)
00120 {
00121 cpl_ensure(aTable, CPL_ERROR_NULL_INPUT, 0.);
00122 cpl_table_unselect_all(aTable);
00123 cpl_table_or_selected_double(aTable, "lambda", CPL_NOT_LESS_THAN,
00124 kMuseNominalLambdaMin);
00125 cpl_table_and_selected_double(aTable, "lambda", CPL_NOT_GREATER_THAN,
00126 kMuseNominalLambdaMax);
00127 cpl_size nsel = cpl_table_count_selected(aTable);
00128 cpl_array *asel = cpl_table_where_selected(aTable);
00129 cpl_size *sel = cpl_array_get_data_cplsize(asel);
00130 double lmin = cpl_table_get_double(aTable, "lambda", sel[0], NULL),
00131 lmax = cpl_table_get_double(aTable, "lambda", sel[nsel - 1], NULL);
00132 cpl_array_delete(asel);
00133 return (lmax - lmin) / nsel;
00134 }
00135
00136
00166
00167 cpl_error_code
00168 muse_flux_reference_table_check(cpl_table *aTable)
00169 {
00170 cpl_ensure_code(aTable, CPL_ERROR_NULL_INPUT);
00171
00172 const char *flam1 = "erg/s/cm**2/Angstrom",
00173 *flam2 = "erg/s/cm^2/Angstrom";
00174
00175 cpl_error_code rc = CPL_ERROR_NONE;
00176 cpl_errorstate prestate = cpl_errorstate_get();
00177
00178 if (cpl_table_has_column(aTable, "lambda") &&
00179 cpl_table_has_column(aTable, "flux") &&
00180 cpl_table_get_column_unit(aTable, "lambda") &&
00181 cpl_table_get_column_unit(aTable, "flux") &&
00182 !strncmp(cpl_table_get_column_unit(aTable, "lambda"), "Angstrom", 9) &&
00183 (!strncmp(cpl_table_get_column_unit(aTable, "flux"), flam1, strlen(flam1)) ||
00184 !strncmp(cpl_table_get_column_unit(aTable, "flux"), flam2, strlen(flam2)))) {
00185
00186
00187 if (cpl_table_get_column_type(aTable, "lambda") != CPL_TYPE_DOUBLE) {
00188 cpl_msg_debug(__func__, "Casting lambda column to double");
00189 cpl_table_cast_column(aTable, "lambda", NULL, CPL_TYPE_DOUBLE);
00190 }
00191 if (cpl_table_get_column_type(aTable, "flux") != CPL_TYPE_DOUBLE) {
00192 cpl_msg_debug(__func__, "Casting flux column to double");
00193 cpl_table_cast_column(aTable, "flux", NULL, CPL_TYPE_DOUBLE);
00194 }
00195
00196 if (cpl_table_has_column(aTable, "fluxerr")) {
00197 if (cpl_table_get_column_type(aTable, "fluxerr") != CPL_TYPE_DOUBLE) {
00198 cpl_msg_debug(__func__, "Casting fluxerr column to double");
00199 cpl_table_cast_column(aTable, "fluxerr", NULL, CPL_TYPE_DOUBLE);
00200 }
00201 const char *unit = cpl_table_get_column_unit(aTable, "fluxerr");
00202 if (!unit || (strncmp(unit, flam1, strlen(flam1)) &&
00203 strncmp(unit, flam2, strlen(flam2)))) {
00204 cpl_msg_debug(__func__, "Erasing fluxerr column because of unexpected "
00205 "unit (%s)", unit);
00206 cpl_table_erase_column(aTable, "fluxerr");
00207 }
00208 }
00209 cpl_msg_info(__func__, "Found MUSE format, average sampling %.3f Angstrom/bin"
00210 " over MUSE range", muse_flux_reference_table_sampling(aTable));
00211 } else if (cpl_table_has_column(aTable, "WAVELENGTH") &&
00212 cpl_table_has_column(aTable, "FLUX") &&
00213 cpl_table_get_column_unit(aTable, "WAVELENGTH") &&
00214 cpl_table_get_column_unit(aTable, "FLUX") &&
00215 !strncmp(cpl_table_get_column_unit(aTable, "WAVELENGTH"), "ANGSTROMS", 10) &&
00216 !strncmp(cpl_table_get_column_unit(aTable, "FLUX"), "FLAM", 5)) {
00217 #if 0
00218 printf("input HST CALSPEC table:\n");
00219 cpl_table_dump_structure(aTable, stdout);
00220 cpl_table_dump(aTable, cpl_table_get_nrow(aTable)/2, 3, stdout);
00221 fflush(stdout);
00222 #endif
00223
00224 cpl_table_cast_column(aTable, "WAVELENGTH", "lambda", CPL_TYPE_DOUBLE);
00225 cpl_table_cast_column(aTable, "FLUX", "flux", CPL_TYPE_DOUBLE);
00226 cpl_table_erase_column(aTable, "WAVELENGTH");
00227 cpl_table_erase_column(aTable, "FLUX");
00228 cpl_table_set_column_unit(aTable, "lambda", "Angstrom");
00229 cpl_table_set_column_unit(aTable, "flux", flam1);
00230
00231
00232 if (cpl_table_has_column(aTable, "STATERROR") &&
00233 cpl_table_has_column(aTable, "SYSERROR") &&
00234 cpl_table_get_column_unit(aTable, "STATERROR") &&
00235 cpl_table_get_column_unit(aTable, "SYSERROR") &&
00236 !strncmp(cpl_table_get_column_unit(aTable, "STATERROR"), "FLAM", 5) &&
00237 !strncmp(cpl_table_get_column_unit(aTable, "SYSERROR"), "FLAM", 5)) {
00238
00239
00240 cpl_table_cast_column(aTable, "STATERROR", "fluxerr", CPL_TYPE_DOUBLE);
00241 cpl_table_erase_column(aTable, "STATERROR");
00242 cpl_table_cast_column(aTable, "SYSERROR", NULL, CPL_TYPE_DOUBLE);
00243 cpl_table_power_column(aTable, "fluxerr", 2);
00244 cpl_table_power_column(aTable, "SYSERROR", 2);
00245 cpl_table_add_columns(aTable, "fluxerr", "SYSERROR");
00246 cpl_table_erase_column(aTable, "SYSERROR");
00247 cpl_table_power_column(aTable, "fluxerr", 0.5);
00248 cpl_table_set_column_unit(aTable, "fluxerr", flam1);
00249 }
00250
00251
00252
00253
00254
00255 if (cpl_table_has_column(aTable, "FWHM")) {
00256 cpl_table_erase_column(aTable, "FWHM");
00257 }
00258 if (cpl_table_has_column(aTable, "DATAQUAL")) {
00259 cpl_table_erase_column(aTable, "DATAQUAL");
00260 }
00261 if (cpl_table_has_column(aTable, "TOTEXP")) {
00262 cpl_table_erase_column(aTable, "TOTEXP");
00263 }
00264
00265 cpl_size irow, nrow = cpl_table_get_nrow(aTable);
00266 for (irow = 0; irow < nrow; irow++) {
00267 double lambda = cpl_table_get_double(aTable, "lambda", irow, NULL);
00268 cpl_table_set_double(aTable, "lambda", irow,
00269 muse_astro_wavelength_vacuum_to_air(lambda));
00270 }
00271 #if 0
00272 printf("converted HST CALSPEC table:\n");
00273 cpl_table_dump_structure(aTable, stdout);
00274 cpl_table_dump(aTable, cpl_table_get_nrow(aTable)/2, 3, stdout);
00275 fflush(stdout);
00276 #endif
00277 cpl_msg_info(__func__, "Found HST CALSPEC format on input, converted to "
00278 "MUSE format; average sampling %.3f Angstrom/bin over MUSE "
00279 "range (assumed vacuum wavelengths on input, converted to air).",
00280 muse_flux_reference_table_sampling(aTable));
00281 } else {
00282 cpl_msg_error(__func__, "Unknown format found!");
00283 #if 0
00284 cpl_table_dump_structure(aTable, stdout);
00285 #endif
00286 rc = CPL_ERROR_INCOMPATIBLE_INPUT;
00287 }
00288
00289
00290 if (!cpl_errorstate_is_equal(prestate)) {
00291 rc = cpl_error_get_code();
00292 }
00293 return rc;
00294 }
00295
00296
00334
00335 double
00336 muse_flux_response_interpolate(const cpl_table *aResponse, double aLambda,
00337 double *aError, muse_flux_interpolation_type aType)
00338 {
00339 double rv = 0.;
00340 if (aType == MUSE_FLUX_TELLURIC) {
00341 rv = 1.;
00342 }
00343 cpl_ensure(aResponse, CPL_ERROR_NULL_INPUT, rv);
00344 int size = cpl_table_get_nrow(aResponse);
00345 cpl_ensure(size > 0, cpl_error_get_code(), rv);
00346
00347
00348 const double *lbda = cpl_table_get_data_double_const(aResponse, "lambda"),
00349 *resp = NULL, *rerr = NULL;
00350 switch (aType) {
00351 case MUSE_FLUX_RESP_FILTER:
00352 resp = cpl_table_get_data_double_const(aResponse, "throughput");
00353 break;
00354 case MUSE_FLUX_RESP_FLUX:
00355 resp = cpl_table_get_data_double_const(aResponse, "response");
00356 if (aError) {
00357 rerr = cpl_table_get_data_double_const(aResponse, "resperr");
00358 }
00359 break;
00360 case MUSE_FLUX_RESP_STD_FLUX:
00361 resp = cpl_table_get_data_double_const(aResponse, "flux");
00362 if (aError) {
00363 rerr = cpl_table_get_data_double_const(aResponse, "fluxerr");
00364 }
00365 break;
00366 case MUSE_FLUX_RESP_EXTINCT:
00367 resp = cpl_table_get_data_double_const(aResponse, "extinction");
00368 break;
00369 case MUSE_FLUX_TELLURIC:
00370 resp = cpl_table_get_data_double_const(aResponse, "ftelluric");
00371 if (aError) {
00372 rerr = cpl_table_get_data_double_const(aResponse, "ftellerr");
00373 }
00374 break;
00375 default:
00376 cpl_error_set(__func__, CPL_ERROR_UNSUPPORTED_MODE);
00377 return rv;
00378 }
00379 cpl_ensure(lbda && resp, cpl_error_get_code(), rv);
00380 if (aError) {
00381 cpl_ensure(rerr, cpl_error_get_code(), rv);
00382 }
00383
00384
00385 if (aLambda < lbda[0]) {
00386 return rv;
00387 }
00388 if (aLambda > lbda[size-1]) {
00389 return rv;
00390 }
00391
00392
00393 double response = rv, resperror = 0.;
00394 int l = 0, r = size - 1,
00395 m = (l + r) / 2;
00396 while (CPL_TRUE) {
00397 if (aLambda >= lbda[m] && aLambda <= lbda[m+1]) {
00398
00399 double lquot = (aLambda - lbda[m]) / (lbda[m+1] - lbda[m]);
00400 response = resp[m] + (resp[m+1] - resp[m]) * lquot;
00401 if (rerr) {
00402
00403
00404
00405 resperror = sqrt(pow(rerr[m] * (1 - lquot), 2.)
00406 + pow(rerr[m+1] * lquot, 2.));
00407 }
00408 #if 0
00409 cpl_msg_debug(__func__, "Found at m=%d (%f: %f+/-%f) / "
00410 "m+1=%d (%f: %f+/-%f) -> %f: %f+/-%f",
00411 m, lbda[m], resp[m], rerr ? rerr[m] : 0.,
00412 m+1, lbda[m+1], resp[m+1], rerr ? rerr[m+1] : 0.,
00413 aLambda, response, resperror);
00414 #endif
00415 break;
00416 }
00417
00418 if (aLambda < lbda[m]) {
00419 r = m;
00420 }
00421 if (aLambda > lbda[m]) {
00422 l = m;
00423 }
00424 m = (l + r) / 2;
00425 }
00426
00427 #if 0
00428 cpl_msg_debug(__func__, "Response %g+/-%g at lambda=%fA", response, resperror,
00429 aLambda);
00430 #endif
00431 if (aError && rerr) {
00432 *aError = resperror;
00433 }
00434 return response;
00435 }
00436
00437
00452
00453 static double
00454 muse_flux_image_sky(cpl_image *aImage, double aX, double aY, double aHalfSize,
00455 unsigned int aDSky, float *aError)
00456 {
00457 if (aError) {
00458 *aError = FLT_MAX;
00459 }
00460
00461 int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
00462 y1 = aY - aHalfSize, y2 = aY + aHalfSize;
00463 unsigned char nskyarea = 0;
00464 double skylevel = 0., skyerror = 0.;
00465
00466 cpl_errorstate state = cpl_errorstate_get();
00467 cpl_stats_mode mode = CPL_STATS_MEDIAN | CPL_STATS_MEDIAN_DEV;
00468 cpl_stats *s = cpl_stats_new_from_image_window(aImage, mode,
00469 x1 - aDSky, y1, x1 - 1, y2);
00470 if (s) {
00471
00472
00473 nskyarea++;
00474 skylevel += cpl_stats_get_median(s);
00475 skyerror += pow(cpl_stats_get_median_dev(s), 2);
00476 cpl_stats_delete(s);
00477 }
00478
00479 s = cpl_stats_new_from_image_window(aImage,mode,
00480 x2 + 1, y1, x2 + aDSky, y2);
00481 if (s) {
00482 nskyarea++;
00483 skylevel += cpl_stats_get_median(s);
00484 skyerror += pow(cpl_stats_get_median_dev(s), 2);
00485 cpl_stats_delete(s);
00486 }
00487
00488 s = cpl_stats_new_from_image_window(aImage,mode,
00489 x1, y1 - aDSky, x2, y1 - 1);
00490 if (s) {
00491 nskyarea++;
00492 skylevel += cpl_stats_get_median(s);
00493 skyerror += pow(cpl_stats_get_median_dev(s), 2);
00494 cpl_stats_delete(s);
00495 }
00496
00497 s = cpl_stats_new_from_image_window(aImage,mode,
00498 x1, y2 + 1, x2, y2 + aDSky);
00499 if (s) {
00500 nskyarea++;
00501 skylevel += cpl_stats_get_median(s);
00502 skyerror += pow(cpl_stats_get_median_dev(s), 2);
00503 cpl_stats_delete(s);
00504 }
00505 if (nskyarea == 0) {
00506 return 0.;
00507 }
00508 skylevel /= nskyarea;
00509 skyerror = sqrt(skyerror) / nskyarea;
00510 if (!cpl_errorstate_is_equal(state)) {
00511 cpl_errorstate_set(state);
00512 }
00513 #if 0
00514 cpl_msg_debug(__func__, "skylevel = %f +/- %f (%u sky areas)",
00515 skylevel, skyerror, nskyarea);
00516 #endif
00517 if (aError) {
00518 *aError = skyerror;
00519 }
00520 return skylevel;
00521 }
00522
00523
00540
00541 static double
00542 muse_flux_image_gaussian(cpl_image *aImage, cpl_image *aImErr, double aX,
00543 double aY, double aHalfSize, unsigned int aDSky,
00544 unsigned int aMaxBad, float *aFErr)
00545 {
00546
00547
00548 UNUSED_ARGUMENT(aMaxBad);
00549
00550 if (aFErr) {
00551 *aFErr = FLT_MAX;
00552 }
00553
00554 cpl_array *params = cpl_array_new(7, CPL_TYPE_DOUBLE),
00555 *parerr = cpl_array_new(7, CPL_TYPE_DOUBLE);
00556
00557
00558
00559 cpl_errorstate state = cpl_errorstate_get();
00560 double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky, NULL);
00561 if (!cpl_errorstate_is_equal(state)) {
00562
00563
00564 cpl_errorstate_set(state);
00565 }
00566 cpl_array_set_double(params, 0, skylevel);
00567 cpl_array_set_double(params, 3, aX);
00568 cpl_array_set_double(params, 4, aY);
00569 double rms = 0, chisq = 0;
00570
00571 int nx = cpl_image_get_size_x(aImage),
00572 ny = cpl_image_get_size_y(aImage),
00573 xsize = fmin(aHalfSize, fmin(aX - 1., nx - aX)) * 2,
00574 ysize = fmin(aHalfSize, fmin(aY - 1., ny - aY)) * 2;
00575 cpl_error_code rc = cpl_fit_image_gaussian(aImage, aImErr, aX, aY, xsize, ysize,
00576 params, parerr, NULL, &rms, &chisq,
00577 NULL, NULL, NULL, NULL, NULL);
00578 if (rc != CPL_ERROR_NONE) {
00579 if (rc != CPL_ERROR_ILLEGAL_INPUT) {
00580 cpl_msg_debug(__func__, "rc = %d: %s", rc, cpl_error_get_message());
00581 }
00582 cpl_array_delete(params);
00583 cpl_array_delete(parerr);
00584 return 0;
00585 }
00586 double flux = cpl_array_get_double(params, 1, NULL),
00587 ferr = cpl_array_get_double(parerr, 1, NULL);
00588 #if 0
00589 double fwhmx = cpl_array_get_double(params, 5, NULL) * CPL_MATH_FWHM_SIG,
00590 fwhmy = cpl_array_get_double(params, 6, NULL) * CPL_MATH_FWHM_SIG;
00591 cpl_msg_debug(__func__, "%.3f,%.3f: %g+/-%g (bg: %g, FWHM: %.3f,%.3f, %g, %g)",
00592 cpl_array_get_double(params, 3, NULL), cpl_array_get_double(params, 4, NULL),
00593 flux, ferr, cpl_array_get_double(params, 0, NULL), fwhmx, fwhmy,
00594 rms, chisq);
00595 #endif
00596 #if 0
00597 cpl_msg_debug(__func__, "skylevel = %f", cpl_array_get_double(params, 0, NULL));
00598 cpl_msg_debug(__func__, "measured flux %f +/- %f", flux, ferr);
00599 #endif
00600 cpl_array_delete(params);
00601 cpl_array_delete(parerr);
00602 if (aFErr) {
00603 *aFErr = ferr;
00604 }
00605 return flux;
00606 }
00607
00608
00627
00628 static double
00629 muse_flux_image_moffat(cpl_image *aImage, cpl_image *aImErr, double aX,
00630 double aY, double aHalfSize, unsigned int aDSky,
00631 unsigned int aMaxBad, float *aFErr)
00632 {
00633 if (aFErr) {
00634 *aFErr = FLT_MAX;
00635 }
00636
00637 int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
00638 y1 = aY - aHalfSize, y2 = aY + aHalfSize,
00639 nx = cpl_image_get_size_x(aImage),
00640 ny = cpl_image_get_size_y(aImage);
00641 if (x1 < 1) {
00642 x1 = 1;
00643 }
00644 if (x2 > nx) {
00645 x2 = nx;
00646 }
00647 if (y1 < 1) {
00648 y1 = 1;
00649 }
00650 if (y2 > ny) {
00651 y2 = ny;
00652 }
00653 int npoints = (x2 - x1 + 1) * (y2 - y1 + 1);
00654
00655 cpl_matrix *pos = cpl_matrix_new(npoints, 2);
00656 cpl_vector *values = cpl_vector_new(npoints),
00657 *errors = cpl_vector_new(npoints);
00658 float *derr = cpl_image_get_data_float(aImErr);
00659 int i, idx = 0;
00660 for (i = x1; i <= x2; i++) {
00661 int j;
00662 for (j = y1; j <= y2; j++) {
00663 int err;
00664 double data = cpl_image_get(aImage, i, j, &err);
00665 if (err) {
00666 continue;
00667 }
00668 cpl_matrix_set(pos, idx, 0, i);
00669 cpl_matrix_set(pos, idx, 1, j);
00670 cpl_vector_set(values, idx, data);
00671 cpl_vector_set(errors, idx, derr[(i-1) + (j-1)*nx]);
00672 idx++;
00673 }
00674 }
00675
00676
00677 if (idx < 16 || (npoints - idx) > (int)aMaxBad) {
00678 cpl_matrix_delete(pos);
00679 cpl_vector_delete(values);
00680 cpl_vector_delete(errors);
00681 return 0;
00682 }
00683 cpl_matrix_set_size(pos, idx, 2);
00684 cpl_vector_set_size(values, idx);
00685 cpl_vector_set_size(errors, idx);
00686
00687 cpl_array *params = cpl_array_new(8, CPL_TYPE_DOUBLE),
00688 *parerr = cpl_array_new(8, CPL_TYPE_DOUBLE);
00689
00690 cpl_errorstate state = cpl_errorstate_get();
00691 double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky, NULL);
00692 if (!cpl_errorstate_is_equal(state)) {
00693
00694
00695 cpl_errorstate_set(state);
00696 }
00697 cpl_array_set_double(params, 0, skylevel);
00698 cpl_array_set_double(params, 2, aX);
00699 cpl_array_set_double(params, 3, aY);
00700 double rms = 0, chisq = 0;
00701 cpl_error_code rc = muse_utils_fit_moffat_2d(pos, values, errors,
00702 params, parerr, NULL,
00703 &rms, &chisq);
00704 cpl_matrix_delete(pos);
00705 cpl_vector_delete(values);
00706 cpl_vector_delete(errors);
00707 if (rc != CPL_ERROR_NONE) {
00708 if (rc != CPL_ERROR_ILLEGAL_INPUT) {
00709 cpl_msg_debug(__func__, "rc = %d: %s", rc, cpl_error_get_message());
00710 }
00711 cpl_array_delete(params);
00712 cpl_array_delete(parerr);
00713 return 0;
00714 }
00715
00716 double flux = cpl_array_get_double(params, 1, NULL);
00717 if (aFErr) {
00718 *aFErr = cpl_array_get_double(parerr, 1, NULL);
00719 }
00720 #if 0
00721 cpl_msg_debug(__func__, "skylevel = %f", cpl_array_get_double(params, 0, NULL));
00722 cpl_msg_debug(__func__, "measured flux %f +/- %f", flux, cpl_array_get_double(parerr, 1, NULL));
00723 #endif
00724 cpl_array_delete(params);
00725 cpl_array_delete(parerr);
00726 return flux;
00727 }
00728
00729
00748
00749 static double
00750 muse_flux_image_square(cpl_image *aImage, cpl_image *aImErr, double aX,
00751 double aY, double aHalfSize, unsigned int aDSky,
00752 unsigned int aMaxBad, float *aFErr)
00753 {
00754 if (aFErr) {
00755 *aFErr = FLT_MAX;
00756 }
00757 int x1 = aX - aHalfSize, x2 = aX + aHalfSize,
00758 y1 = aY - aHalfSize, y2 = aY + aHalfSize,
00759 nx = cpl_image_get_size_x(aImage),
00760 ny = cpl_image_get_size_y(aImage);
00761 if (x1 < 1) {
00762 x1 = 1;
00763 }
00764 if (x2 > nx) {
00765 x2 = nx;
00766 }
00767 if (y1 < 1) {
00768 y1 = 1;
00769 }
00770 if (y2 > ny) {
00771 y2 = ny;
00772 }
00773 float skyerror;
00774 cpl_errorstate state = cpl_errorstate_get();
00775 double skylevel = muse_flux_image_sky(aImage, aX, aY, aHalfSize, aDSky,
00776 &skyerror);
00777 if (!cpl_errorstate_is_equal(state)) {
00778
00779
00780 cpl_errorstate_set(state);
00781 cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
00782 return 0.;
00783 }
00784
00785
00786
00787 cpl_image *region = cpl_image_extract(aImage, x1, y1, x2, y2);
00788 #if 0
00789 cpl_msg_debug(__func__, "region [%d:%d,%d:%d] %"CPL_SIZE_FORMAT" bad pixels",
00790 x1, y1, x2, y2, cpl_image_count_rejected(region));
00791 #endif
00792 if (cpl_image_count_rejected(region) > aMaxBad) {
00793 cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
00794 cpl_image_delete(region);
00795 return 0.;
00796 }
00797 cpl_image *regerr = cpl_image_extract(aImErr, x1, y1, x2, y2);
00798 if (cpl_image_count_rejected(region) > 0) {
00799 cpl_detector_interpolate_rejected(region);
00800 cpl_detector_interpolate_rejected(regerr);
00801 }
00802
00803
00804
00805 int npoints = (x2 - x1 + 1) * (y2 - y1 + 1),
00806
00807
00808
00809 nsky = 2 * aDSky * (x2 - x1 + y2 - y1 + 2);
00810 double flux = cpl_image_get_flux(region) - skylevel * npoints,
00811 ferr = sqrt(cpl_image_get_sqflux(regerr)
00812 + npoints * skyerror*skyerror * (1. + (double)npoints / nsky));
00813 #if 0
00814 cpl_msg_debug(__func__, "measured flux %f +/- %f (%d object pixels, %d pixels"
00815 " with %f with sky %f)", flux, ferr, npoints, nsky,
00816 cpl_image_get_flux(region), cpl_image_get_sqflux(regerr));
00817 #endif
00818 cpl_image_delete(region);
00819 cpl_image_delete(regerr);
00820 if (aFErr) {
00821 *aFErr = ferr;
00822 }
00823 return flux;
00824 }
00825
00826
00846
00847 static double
00848 muse_flux_image_circle(cpl_image *aImage, cpl_image *aImErr, double aX,
00849 double aY, double aAper, double aAnnu, double aDAnnu,
00850 unsigned int aMaxBad, float *aFErr)
00851 {
00852 if (aFErr) {
00853 *aFErr = FLT_MAX;
00854 }
00855 double rmax = ceil(fmax(aAper, aAnnu + aDAnnu));
00856 int x1 = aX - rmax, x2 = aX + rmax,
00857 y1 = aY - rmax, y2 = aY + rmax,
00858 nx = cpl_image_get_size_x(aImage),
00859 ny = cpl_image_get_size_y(aImage);
00860 if (x1 < 1) {
00861 x1 = 1;
00862 }
00863 if (x2 > nx) {
00864 x2 = nx;
00865 }
00866 if (y1 < 1) {
00867 y1 = 1;
00868 }
00869 if (y2 > ny) {
00870 y2 = ny;
00871 }
00872
00873
00874 cpl_vector *vbg = cpl_vector_new((x2 - x1 + 1) * (y2 - y1 + 1)),
00875 *vbe = cpl_vector_new((x2 - x1 + 1) * (y2 - y1 + 1));
00876 unsigned int nbad = 0, nbg = 0;
00877 int i;
00878 for (i = x1; i <= x2; i++) {
00879 int j;
00880 for (j = y1; j <= y2; j++) {
00881 double r = sqrt(pow(aX - i, 2) + pow(aY - j, 2));
00882 if (r <= aAper) {
00883 nbad += cpl_image_is_rejected(aImage, i, j) == 1;
00884 }
00885 if (r < aAnnu || r > aAnnu + aDAnnu) {
00886 continue;
00887 }
00888 int err;
00889 double value = cpl_image_get(aImage, i, j, &err);
00890 if (err) {
00891 continue;
00892 }
00893 cpl_vector_set(vbg, nbg, value);
00894 cpl_vector_set(vbe, nbg, cpl_image_get(aImErr, i, j, &err));
00895 nbg++;
00896 }
00897 }
00898 if (nbg <= 0) {
00899 cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
00900 cpl_vector_delete(vbg);
00901 cpl_vector_delete(vbe);
00902 return 0.;
00903 }
00904 cpl_vector_set_size(vbg, nbg);
00905 cpl_vector_set_size(vbe, nbg);
00906 cpl_matrix *pos = cpl_matrix_new(1, nbg);
00907 double mse;
00908 cpl_polynomial *fit = muse_utils_iterate_fit_polynomial(pos, vbg, vbe, NULL,
00909 0, 3., &mse, NULL);
00910 #if 0
00911 unsigned int nrej = nbg - cpl_vector_get_size(vbg);
00912 #endif
00913 nbg = cpl_vector_get_size(vbg);
00914 cpl_size pows = 0;
00915 double smean = cpl_polynomial_get_coeff(fit, &pows),
00916 sstdev = sqrt(mse);
00917 cpl_polynomial_delete(fit);
00918 cpl_matrix_delete(pos);
00919 cpl_vector_delete(vbg);
00920 cpl_vector_delete(vbe);
00921 #if 0
00922 cpl_msg_debug(__func__, "sky: %d pixels (%d rejected), %f +/- %f; found %d "
00923 "bad pixels inside aperture", nbg, nrej, smean, sstdev, nbad);
00924 #endif
00925 if (nbad > aMaxBad) {
00926 cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
00927 return 0.;
00928 }
00929
00930
00931 if (nbad > 0) {
00932 cpl_detector_interpolate_rejected(aImage);
00933 cpl_detector_interpolate_rejected(aImErr);
00934 }
00935
00936
00937 double flux = 0.,
00938 ferr = 0.;
00939 unsigned int nobj = 0;
00940 for (i = x1; i <= x2; i++) {
00941 int j;
00942 for (j = y1; j <= y2; j++) {
00943 double r = sqrt(pow(aX - i, 2) + pow(aY - j, 2));
00944 if (r > aAper) {
00945 continue;
00946 }
00947 int err;
00948 double value = cpl_image_get(aImage, i, j, &err),
00949 error = cpl_image_get(aImErr, i, j, &err);
00950 flux += value;
00951 ferr += error*error;
00952 nobj++;
00953 }
00954 }
00955 flux -= smean * nobj;
00956
00957
00958
00959
00960 ferr = sqrt(ferr + nobj * sstdev*sstdev * (1. + (double)nobj / nbg));
00961 #if 0
00962 cpl_msg_debug(__func__, "flux: %d pixels (%d interpolated), %f +/- %f",
00963 nobj, nbad, flux, ferr);
00964 #endif
00965 if (aFErr) {
00966 *aFErr = ferr;
00967 }
00968 return flux;
00969 }
00970
00971
01005
01006 muse_image *
01007 muse_flux_integrate_cube(muse_datacube *aCube, cpl_apertures *aApertures,
01008 muse_flux_profile_type aProfile)
01009 {
01010 cpl_ensure(aCube && aApertures, CPL_ERROR_NULL_INPUT, NULL);
01011 switch (aProfile) {
01012 case MUSE_FLUX_PROFILE_GAUSSIAN:
01013 cpl_msg_info(__func__, "Gaussian profile fits for flux integration");
01014 break;
01015 case MUSE_FLUX_PROFILE_MOFFAT:
01016 cpl_msg_info(__func__, "Moffat profile fits for flux integration");
01017 break;
01018 case MUSE_FLUX_PROFILE_CIRCLE:
01019 cpl_msg_info(__func__, "Circular flux integration");
01020 break;
01021 case MUSE_FLUX_PROFILE_EQUAL_SQUARE:
01022 cpl_msg_info(__func__, "Simple square window flux integration");
01023 break;
01024 default:
01025 cpl_msg_error(__func__, "Unknown flux integration method!");
01026 cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
01027 return NULL;
01028 }
01029
01030
01031 int naper = cpl_apertures_get_size(aApertures),
01032 nlambda = cpl_imagelist_get_size(aCube->data),
01033 nplane = cpl_imagelist_get_size(aCube->data) / 2;
01034 cpl_image *cim = cpl_imagelist_get(aCube->data, nplane);
01035 muse_image *intimage = muse_image_new();
01036 intimage->data = cpl_image_new(nlambda, naper, CPL_TYPE_FLOAT);
01037 intimage->dq = cpl_image_new(nlambda, naper, CPL_TYPE_INT);
01038 intimage->stat = cpl_image_new(nlambda, naper, CPL_TYPE_FLOAT);
01039
01040 intimage->header = cpl_propertylist_new();
01041 cpl_propertylist_append_double(intimage->header, "CRVAL1",
01042 muse_pfits_get_crval(aCube->header, 3));
01043 cpl_propertylist_append_double(intimage->header, "CRPIX1",
01044 muse_pfits_get_crpix(aCube->header, 3));
01045 cpl_propertylist_append_double(intimage->header, "CD1_1",
01046 muse_pfits_get_cd(aCube->header, 3, 3));
01047 cpl_propertylist_append_string(intimage->header, "CTYPE1",
01048 muse_pfits_get_ctype(aCube->header, 3));
01049 cpl_propertylist_append_string(intimage->header, "CUNIT1",
01050 muse_pfits_get_cunit(aCube->header, 3));
01051
01052 cpl_propertylist_append_double(intimage->header, "CRVAL2", 1.);
01053 cpl_propertylist_append_double(intimage->header, "CRPIX2", 1.);
01054 cpl_propertylist_append_double(intimage->header, "CD2_2", 1.);
01055 cpl_propertylist_append_string(intimage->header, "CTYPE2", "PIXEL");
01056 cpl_propertylist_append_string(intimage->header, "CUNIT2", "pixel");
01057 cpl_propertylist_append_double(intimage->header, "CD1_2", 0.);
01058 cpl_propertylist_append_double(intimage->header, "CD2_1", 0.);
01059
01060 cpl_propertylist_append_string(intimage->header, "DATE-OBS",
01061 cpl_propertylist_get_string(aCube->header,
01062 "DATE-OBS"));
01063 cpl_propertylist_append_string(intimage->header, "BUNIT",
01064 muse_pfits_get_bunit(aCube->header));
01065 cpl_propertylist_append_double(intimage->header, "EXPTIME",
01066 muse_pfits_get_exptime(aCube->header));
01067 cpl_propertylist_append_string(intimage->header, "ESO INS MODE",
01068 cpl_propertylist_get_string(aCube->header,
01069 "ESO INS MODE"));
01070
01071
01072 cpl_errorstate ps = cpl_errorstate_get();
01073 double fwhm = (muse_pfits_get_fwhm_start(aCube->header)
01074 + muse_pfits_get_fwhm_end(aCube->header)) / 2.;
01075 if (muse_pfits_get_mode(aCube->header) < MUSE_MODE_NFM_AO_N) {
01076 fwhm /= (kMuseSpaxelSizeX_WFM + kMuseSpaxelSizeY_WFM) / 2.;
01077 } else {
01078 fwhm /= (kMuseSpaxelSizeX_NFM + kMuseSpaxelSizeY_NFM) / 2.;
01079 }
01080 if (!cpl_errorstate_is_equal(ps)) {
01081 double xc = cpl_apertures_get_centroid_x(aApertures, 1),
01082 yc = cpl_apertures_get_centroid_y(aApertures, 1),
01083 xfwhm, yfwhm;
01084 cpl_image_get_fwhm(cim, lround(xc), lround(yc), &xfwhm, &yfwhm);
01085 if (xfwhm > 0. && yfwhm > 0.) {
01086 fwhm = (xfwhm + yfwhm) / 2.;
01087 } else if (xfwhm > 0.) {
01088 fwhm = xfwhm;
01089 } else if (yfwhm > 0.) {
01090 fwhm = yfwhm;
01091 } else {
01092 fwhm = 5.;
01093 }
01094 cpl_errorstate_set(ps);
01095 cpl_msg_debug(__func__, "Using roughly estimated reference FWHM (%.3f pix) "
01096 "instead of DIMM seeing", fwhm);
01097 } else {
01098 cpl_msg_debug(__func__, "Using DIMM seeing of %.3f pix for reference FWHM",
01099 fwhm);
01100 }
01101
01102
01103
01104 cpl_image *sizes = cpl_image_new(nlambda, naper, CPL_TYPE_DOUBLE);
01105 double *psizes = cpl_image_get_data_double(sizes);
01106
01107 float *data = cpl_image_get_data_float(intimage->data),
01108 *stat = cpl_image_get_data_float(intimage->stat);
01109 int *dq = cpl_image_get_data_int(intimage->dq);
01110
01111 int l, ngood = 0, nillegal = 0, nbadbg = 0;
01112 #pragma omp parallel for default(none) \
01113 shared(aApertures, aCube, aProfile, data, dq, fwhm, naper, nbadbg, \
01114 ngood, nillegal, nlambda, psizes, stat)
01115 for (l = 0; l < nlambda; l++) {
01116 cpl_image *plane = cpl_imagelist_get(aCube->data, l),
01117 *pldq = aCube->dq ? cpl_imagelist_get(aCube->dq, l) : NULL,
01118 *plerr = cpl_image_duplicate(cpl_imagelist_get(aCube->stat, l));
01119 #if 0
01120 cpl_stats *stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
01121 cpl_msg_debug(__func__, "lambda = %d/%f %s", l + 1,
01122 (l + 1 - muse_pfits_get_crpix(aCube->header, 3))
01123 * muse_pfits_get_cd(aCube->header, 3, 3)
01124 + muse_pfits_get_crval(aCube->header, 3),
01125 muse_pfits_get_cunit(aCube->header, 3));
01126 cpl_msg_debug(__func__, "variance: %g...%g...%g", cpl_stats_get_min(stats),
01127 cpl_stats_get_mean(stats), cpl_stats_get_max(stats));
01128 cpl_stats_delete(stats);
01129 #endif
01130
01131 if (pldq) {
01132 muse_quality_image_reject_using_dq(plane, pldq, plerr);
01133 } else {
01134 cpl_image_reject_value(plane, CPL_VALUE_NAN);
01135 cpl_image_reject_value(plerr, CPL_VALUE_NAN);
01136 }
01137 #if 0
01138 stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
01139 cpl_msg_debug(__func__, "cut variance: %g...%g...%g (%"CPL_SIZE_FORMAT" bad"
01140 " pixel)", cpl_stats_get_min(stats), cpl_stats_get_mean(stats),
01141 cpl_stats_get_max(stats), cpl_image_count_rejected(plane));
01142 cpl_stats_delete(stats);
01143 #endif
01144
01145 cpl_image_power(plerr, 0.5);
01146 #if 0
01147 stats = cpl_stats_new_from_image(plerr, CPL_STATS_ALL);
01148 cpl_msg_debug(__func__, "errors: %g...%g...%g", cpl_stats_get_min(stats),
01149 cpl_stats_get_mean(stats), cpl_stats_get_max(stats));
01150 cpl_stats_delete(stats);
01151 #endif
01152 cpl_errorstate state = cpl_errorstate_get();
01153 int n;
01154 for (n = 1; n <= naper; n++) {
01155
01156 double xc = cpl_apertures_get_centroid_x(aApertures, n),
01157 yc = cpl_apertures_get_centroid_y(aApertures, n),
01158 size = sqrt(cpl_apertures_get_npix(aApertures, n)),
01159 xfwhm, yfwhm;
01160 cpl_errorstate prestate = cpl_errorstate_get();
01161 cpl_image_get_fwhm(plane, lround(xc), lround(yc), &xfwhm, &yfwhm);
01162 if (xfwhm < 0 || yfwhm < 0) {
01163 data[l + (n-1) * nlambda] = 0.;
01164 stat[l + (n-1) * nlambda] = FLT_MAX;
01165 cpl_errorstate_set(prestate);
01166 continue;
01167 }
01168
01169 double halfsize = fmax(1.5 * (xfwhm + yfwhm), 3. * fwhm);
01170 if (halfsize < size / 2) {
01171 halfsize = size / 2;
01172 }
01173 psizes[l + (n-1) * nlambda] = halfsize;
01174 #if 0
01175 cpl_msg_debug(__func__, "%.2f,%.2f FWHM %.2f %.2f size %.2f --> %.2f",
01176 xc, yc, xfwhm, yfwhm, size, halfsize * 2.);
01177 #endif
01178
01179 switch (aProfile) {
01180 case MUSE_FLUX_PROFILE_GAUSSIAN:
01181 data[l + (n-1) * nlambda] = muse_flux_image_gaussian(plane, plerr, xc, yc,
01182 halfsize, 5, 10,
01183 &stat[l + (n-1) * nlambda]);
01184 break;
01185 case MUSE_FLUX_PROFILE_MOFFAT:
01186 data[l + (n-1) * nlambda] = muse_flux_image_moffat(plane, plerr, xc, yc,
01187 halfsize, 5, 10,
01188 &stat[l + (n-1) * nlambda]);
01189 break;
01190 case MUSE_FLUX_PROFILE_CIRCLE: {
01191
01192
01193 double radius = 4./3. * halfsize,
01194 rannu = radius * 5. / 4.;
01195 psizes[l + (n-1) * nlambda] = radius;
01196 data[l + (n-1) * nlambda] = muse_flux_image_circle(plane, plerr, xc, yc,
01197 radius, rannu, 10, 10,
01198 &stat[l + (n-1) * nlambda]);
01199 break;
01200 }
01201 default:
01202 data[l + (n-1) * nlambda] = muse_flux_image_square(plane, plerr, xc, yc,
01203 halfsize, 5, 10,
01204 &stat[l + (n-1) * nlambda]);
01205 }
01206 if (data[l + (n-1) * nlambda] < 0 || !isfinite(data[l + (n-1) * nlambda])) {
01207 data[l + (n-1) * nlambda] = 0.;
01208 dq[l + (n-1) * nlambda] = EURO3D_MISSDATA;
01209 stat[l + (n-1) * nlambda] = FLT_MAX;
01210 }
01211 }
01212
01213
01214
01215 if (!cpl_errorstate_is_equal(state)) {
01216 if (cpl_error_get_code() == CPL_ERROR_ILLEGAL_INPUT) {
01217 cpl_errorstate_set(state);
01218 #pragma omp atomic
01219 nillegal++;
01220 } else if (cpl_error_get_code() == CPL_ERROR_DATA_NOT_FOUND) {
01221 cpl_errorstate_set(state);
01222 #pragma omp atomic
01223 nbadbg++;
01224 }
01225 } else {
01226 #pragma omp atomic
01227 ngood++;
01228 }
01229
01230 cpl_image_delete(plerr);
01231 }
01232
01233
01234 cpl_image_reject_value(sizes, CPL_VALUE_ZERO);
01235 int n;
01236 for (n = 1; n <= naper; n++) {
01237 if (aProfile == MUSE_FLUX_PROFILE_CIRCLE) {
01238 cpl_msg_info(__func__, "Radiuses used for circular flux integration for "
01239 "source %d: %f +/- %f (%f) %f..%f", n,
01240 cpl_image_get_mean_window(sizes, 1, n, nlambda, n),
01241 cpl_image_get_stdev_window(sizes, 1, n, nlambda, n),
01242 cpl_image_get_median_window(sizes, 1, n, nlambda, n),
01243 cpl_image_get_min_window(sizes, 1, n, nlambda, n),
01244 cpl_image_get_max_window(sizes, 1, n, nlambda, n));
01245 } else {
01246 cpl_msg_info(__func__, "Half-sizes used for flux integration for source "
01247 "%d: %f +/- %f (%f) %f..%f", n,
01248 cpl_image_get_mean_window(sizes, 1, n, nlambda, n),
01249 cpl_image_get_stdev_window(sizes, 1, n, nlambda, n),
01250 cpl_image_get_median_window(sizes, 1, n, nlambda, n),
01251 cpl_image_get_min_window(sizes, 1, n, nlambda, n),
01252 cpl_image_get_max_window(sizes, 1, n, nlambda, n));
01253 }
01254 }
01255 #if 0
01256 cpl_image_save(sizes, "sizes.fits", CPL_TYPE_UNSPECIFIED, NULL, CPL_IO_CREATE);
01257 #endif
01258 cpl_image_delete(sizes);
01259
01260
01261 cpl_propertylist_append_int(intimage->header, MUSE_HDR_FLUX_NOBJ, naper);
01262
01263
01264 cpl_propertylist *wcs1 = muse_wcs_create_default(),
01265 *wcs = muse_wcs_apply_cd(aCube->header, wcs1);
01266 cpl_propertylist_delete(wcs1);
01267
01268 double crpix1 = muse_pfits_get_crpix(aCube->header, 1)
01269 + (1. + cpl_image_get_size_x(cim)) / 2.,
01270 crpix2 = muse_pfits_get_crpix(aCube->header, 2)
01271 + (1. + cpl_image_get_size_y(cim)) / 2.;
01272 cpl_propertylist_update_double(wcs, "CRPIX1", crpix1);
01273 cpl_propertylist_update_double(wcs, "CRPIX2", crpix2);
01274 cpl_propertylist_update_double(wcs, "CRVAL1", muse_pfits_get_ra(aCube->header));
01275 cpl_propertylist_update_double(wcs, "CRVAL2", muse_pfits_get_dec(aCube->header));
01276 for (n = 1; n <= naper; n++) {
01277
01278 double xc = cpl_apertures_get_centroid_x(aApertures, n),
01279 yc = cpl_apertures_get_centroid_y(aApertures, n),
01280 ra, dec;
01281 muse_wcs_celestial_from_pixel(wcs, xc, yc, &ra, &dec);
01282 double flux = cpl_image_get_flux_window(intimage->data, 1, n, nlambda, n);
01283 cpl_msg_debug(__func__, "Source %02d: %.3f,%.3f pix, %f,%f deg, flux %e %s",
01284 n, xc, yc, ra, dec, flux, muse_pfits_get_bunit(intimage->header));
01285 char kw[KEYWORD_LENGTH];
01286 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_X, n);
01287 cpl_propertylist_append_float(intimage->header, kw, xc);
01288 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_Y, n);
01289 cpl_propertylist_append_float(intimage->header, kw, yc);
01290 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_RA, n);
01291 cpl_propertylist_append_double(intimage->header, kw, ra);
01292 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_DEC, n);
01293 cpl_propertylist_append_double(intimage->header, kw, dec);
01294 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_FLUX, n);
01295 cpl_propertylist_append_double(intimage->header, kw, flux);
01296 }
01297 cpl_propertylist_delete(wcs);
01298
01299 if (nillegal > 0 || nbadbg > 0) {
01300 cpl_msg_warning(__func__, "Successful fits in %d wavelength planes, but "
01301 "encountered %d \"Illegal input\" errors and %d bad "
01302 "background determinations", ngood, nillegal, nbadbg);
01303 } else {
01304 cpl_msg_info(__func__, "Successful fits in %d wavelength planes", ngood);
01305 }
01306
01307 return intimage;
01308 }
01309
01310
01342
01343 cpl_error_code
01344 muse_flux_integrate_std(muse_pixtable *aPixtable, muse_flux_profile_type aProfile,
01345 muse_flux_object *aFluxObj)
01346 {
01347 cpl_ensure_code(aPixtable && aFluxObj, CPL_ERROR_NULL_INPUT);
01348 switch (aProfile) {
01349 case MUSE_FLUX_PROFILE_GAUSSIAN:
01350 case MUSE_FLUX_PROFILE_MOFFAT:
01351 case MUSE_FLUX_PROFILE_CIRCLE:
01352 case MUSE_FLUX_PROFILE_EQUAL_SQUARE:
01353 break;
01354 default:
01355 return cpl_error_set(__func__, CPL_ERROR_ILLEGAL_INPUT);
01356 }
01357
01358 if (getenv("MUSE_DEBUG_FLUX") && atoi(getenv("MUSE_DEBUG_FLUX")) > 2) {
01359 const char *fn = "flux__pixtable.fits";
01360 cpl_msg_info(__func__, "Saving pixel table as \"%s\"", fn);
01361 muse_pixtable_save(aPixtable, fn);
01362 }
01363 muse_resampling_params *params =
01364 muse_resampling_params_new(MUSE_RESAMPLE_WEIGHTED_DRIZZLE);
01365 params->pfx = 1.;
01366 params->pfy = 1.;
01367 params->pfl = 1.;
01368
01369
01370 params->dlambda = kMuseSpectralSamplingA;
01371 params->crtype = MUSE_RESAMPLING_CRSTATS_MEDIAN;
01372 params->crsigma = 25.;
01373 muse_datacube *cube = muse_resampling_cube(aPixtable, params, NULL);
01374 if (cube) {
01375 aFluxObj->cube = cube;
01376 }
01377 muse_resampling_params_delete(params);
01378 if (getenv("MUSE_DEBUG_FLUX") && atoi(getenv("MUSE_DEBUG_FLUX")) >= 2) {
01379 const char *fn = "flux__cube.fits";
01380 cpl_msg_info(__func__, "Saving cube as \"%s\"", fn);
01381 muse_datacube_save(aFluxObj->cube, fn);
01382 }
01383 int nplane = cpl_imagelist_get_size(cube->data) / 2;
01384 cpl_image *cim = cpl_imagelist_get(cube->data, nplane),
01385 *cdq = cpl_imagelist_get(cube->dq, nplane);
01386
01387
01388 muse_quality_image_reject_using_dq(cim, cdq, NULL);
01389
01390 double dsigmas[] = { 50., 30., 10., 8., 6., 5. };
01391 cpl_vector *vsigmas = cpl_vector_wrap(sizeof(dsigmas) / sizeof(double),
01392 dsigmas);
01393 cpl_size isigma = -1;
01394 cpl_apertures *apertures = cpl_apertures_extract(cim, vsigmas, &isigma);
01395 int napertures = apertures ? cpl_apertures_get_size(apertures) : 0;
01396 if (napertures < 1) {
01397
01398 cpl_msg_error(__func__, "No sources for flux integration found down to %.1f"
01399 " sigma limit on plane %d",
01400 cpl_vector_get(vsigmas, cpl_vector_get_size(vsigmas) - 1),
01401 nplane + 1);
01402 cpl_vector_unwrap(vsigmas);
01403 cpl_apertures_delete(apertures);
01404 return cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
01405 }
01406 cpl_msg_debug(__func__, "The %.1f sigma threshold was used to find %d source%s"
01407 " on plane %d", cpl_vector_get(vsigmas, isigma), napertures,
01408 napertures == 1 ? "" : "s", nplane + 1);
01409 cpl_vector_unwrap(vsigmas);
01410 #if 0
01411 cpl_apertures_dump(apertures, stdout);
01412 fflush(stdout);
01413 #endif
01414
01415
01416 muse_image *intimage = muse_flux_integrate_cube(cube, apertures, aProfile);
01417 cpl_apertures_delete(apertures);
01418 aFluxObj->intimage = intimage;
01419
01420 return CPL_ERROR_NONE;
01421 }
01422
01423
01424
01425
01426
01427
01428
01429
01430 static const double kTelluricBands[][4] = {
01431 { 6273., 6320., 6213., 6380. },
01432 { 6864., 6967., 6750., 7130. },
01433 { 7164., 7325., 7070., 7580. },
01434 { 7590., 7700., 7470., 7830. },
01435 { 8131., 8345., 7900., 8600. },
01436 { 8952., 9028., 8850., 9082. },
01437 { 9274., 9770., 9080., 9263. },
01438 { -1., -1., -1., -1. }
01439 };
01440
01441
01445
01446 const muse_cpltable_def muse_response_tellbands_def[] = {
01447 { "lmin", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
01448 "lower limit of the telluric region", CPL_TRUE },
01449 { "lmax", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
01450 "upper limit of the telluric region", CPL_TRUE },
01451 { "bgmin", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
01452 "lower limit of the background region", CPL_TRUE },
01453 { "bgmax", CPL_TYPE_DOUBLE, "Angstrom", "%8.3f",
01454 "upper limit of the background region", CPL_TRUE },
01455 { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
01456 };
01457
01458
01465
01466 static void
01467 muse_flux_response_set_telluric_bands(muse_flux_object *aFluxObj,
01468 const cpl_table *aTellBands)
01469 {
01470 if (!aFluxObj) {
01471 cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
01472 return;
01473 }
01474
01475 if (aTellBands && muse_cpltable_check(aTellBands, muse_response_tellbands_def)
01476 == CPL_ERROR_NONE) {
01477 cpl_msg_debug(__func__, "using given table for telluric bands");
01478 aFluxObj->tellbands = cpl_table_duplicate(aTellBands);
01479 return;
01480 }
01481
01482 unsigned int ntell = sizeof(kTelluricBands) / sizeof(kTelluricBands[0]) - 1;
01483 cpl_msg_debug(__func__, "using builtin regions for telluric bands (%u "
01484 "entries)", ntell);
01485 aFluxObj->tellbands = muse_cpltable_new(muse_response_tellbands_def, ntell);
01486 cpl_table *tb = aFluxObj->tellbands;
01487 unsigned int k;
01488 for (k = 0; k < ntell; k++) {
01489 cpl_table_set_double(tb, "lmin", k, kTelluricBands[k][0]);
01490 cpl_table_set_double(tb, "lmax", k, kTelluricBands[k][1]);
01491 cpl_table_set_double(tb, "bgmin", k, kTelluricBands[k][2]);
01492 cpl_table_set_double(tb, "bgmax", k, kTelluricBands[k][3]);
01493 }
01494 if (getenv("MUSE_DEBUG_FLUX") && atoi(getenv("MUSE_DEBUG_FLUX")) >= 2) {
01495 const char *fn = "flux__tellregions.fits";
01496 cpl_msg_info(__func__, "Saving telluric bands table as \"%s\"", fn);
01497 cpl_table_save(tb, NULL, NULL, fn, CPL_IO_CREATE);
01498 }
01499 return;
01500 }
01501
01502
01515
01516 static void
01517 muse_flux_response_dump_sensitivity(muse_flux_object *aFluxObj,
01518 const char *aName)
01519 {
01520 char *dodebug = getenv("MUSE_DEBUG_FLUX");
01521 if (!dodebug || (dodebug && atoi(dodebug) <= 0)) {
01522 return;
01523 }
01524 char *fn = cpl_sprintf("flux__sens_%s.ascii", aName);
01525 FILE *fp = fopen(fn, "w");
01526 fprintf(fp, "#");
01527 cpl_table_dump(aFluxObj->sensitivity, 0,
01528 cpl_table_get_nrow(aFluxObj->sensitivity), fp);
01529 fclose(fp);
01530 cpl_msg_debug(__func__, "Written %"CPL_SIZE_FORMAT" datapoints to \"%s\"",
01531 cpl_table_get_nrow(aFluxObj->sensitivity), fn);
01532 cpl_free(fn);
01533 }
01534
01535
01554
01555 static void
01556 muse_flux_response_sensitivity(muse_flux_object *aFluxObj,
01557 unsigned int aStar, const cpl_table *aReference,
01558 double aAirmass, const cpl_table *aExtinct)
01559 {
01560 double crval = muse_pfits_get_crval(aFluxObj->intimage->header, 1),
01561 cdelt = muse_pfits_get_cd(aFluxObj->intimage->header, 1, 1),
01562 crpix = muse_pfits_get_crpix(aFluxObj->intimage->header, 1),
01563 exptime = muse_pfits_get_exptime(aFluxObj->intimage->header);
01564 int nlambda = cpl_image_get_size_x(aFluxObj->intimage->data);
01565
01566 aFluxObj->sensitivity = cpl_table_new(nlambda);
01567 cpl_table *sensitivity = aFluxObj->sensitivity;
01568 cpl_table_new_column(sensitivity, "lambda", CPL_TYPE_DOUBLE);
01569 cpl_table_new_column(sensitivity, "sens", CPL_TYPE_DOUBLE);
01570 cpl_table_new_column(sensitivity, "serr", CPL_TYPE_DOUBLE);
01571 cpl_table_new_column(sensitivity, "dq", CPL_TYPE_INT);
01572 cpl_table_set_column_format(sensitivity, "dq", "%u");
01573 float *data = cpl_image_get_data_float(aFluxObj->intimage->data),
01574 *stat = cpl_image_get_data_float(aFluxObj->intimage->stat);
01575 int l, idx = 0;
01576 for (l = 0; l < nlambda; l++) {
01577 if (data[l + aStar*nlambda] <= 0. ||
01578 stat[l + aStar*nlambda] <= 0. ||
01579 stat[l + aStar*nlambda] == FLT_MAX) {
01580 continue;
01581 }
01582
01583 double lambda = crval + cdelt * (l + 1 - crpix),
01584
01585 extinct = !aExtinct ? 0.
01586 : muse_flux_response_interpolate(aExtinct, lambda, NULL,
01587 MUSE_FLUX_RESP_EXTINCT);
01588 cpl_errorstate prestate = cpl_errorstate_get();
01589 double referr = 0.,
01590 ref = muse_flux_response_interpolate(aReference, lambda, &referr,
01591 MUSE_FLUX_RESP_STD_FLUX);
01592
01593 if (!cpl_errorstate_is_equal(prestate)) {
01594 cpl_errorstate_set(prestate);
01595 ref = muse_flux_response_interpolate(aReference, lambda, NULL,
01596 MUSE_FLUX_RESP_STD_FLUX);
01597 }
01598
01599
01600 double c = 2.5 * log10(data[l + aStar*nlambda]
01601 / exptime / cdelt / ref)
01602 + aAirmass * extinct,
01603 cerr = sqrt(pow(referr / ref, 2) + stat[l + aStar*nlambda]
01604 / pow(data[l + aStar*nlambda], 2))
01605 * 2.5 / CPL_MATH_LN10;
01606 cpl_table_set_double(sensitivity, "lambda", idx, lambda);
01607 cpl_table_set_double(sensitivity, "sens", idx, c);
01608 cpl_table_set_double(sensitivity, "serr", idx, cerr);
01609 cpl_table_set_int(sensitivity, "dq", idx, EURO3D_GOODPIXEL);
01610 idx++;
01611 }
01612
01613 cpl_table_set_size(sensitivity, idx);
01614 }
01615
01616
01629
01630 static void
01631 muse_flux_response_mark_questionable(muse_flux_object *aFluxObj)
01632 {
01633 if (!aFluxObj) {
01634 cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
01635 return;
01636 }
01637 if (!aFluxObj->tellbands) {
01638 cpl_error_set(__func__, CPL_ERROR_NULL_INPUT);
01639 return;
01640 }
01641 cpl_table *tsens = aFluxObj->sensitivity,
01642 *tb = aFluxObj->tellbands;
01643
01644
01645 cpl_boolean isnominal = muse_pfits_get_mode(aFluxObj->intimage->header)
01646 > MUSE_MODE_WFM_NONAO_X;
01647
01648
01649 int irow, nfluxes = cpl_table_get_nrow(tsens);
01650 for (irow = 0; irow < nfluxes; irow++) {
01651 double lambda = cpl_table_get_double(tsens, "lambda", irow, NULL);
01652 unsigned int dq = EURO3D_GOODPIXEL,
01653 k, nk = cpl_table_get_nrow(tb);
01654 for (k = 0; k < nk; k++) {
01655 double lmin = cpl_table_get_double(tb, "lmin", k, NULL),
01656 lmax = cpl_table_get_double(tb, "lmax", k, NULL);
01657 if (lambda >= lmin && lambda <= lmax) {
01658 dq |= EURO3D_TELLURIC;
01659 }
01660 }
01661 if (isnominal && lambda < kMuseNominalCutoff) {
01662 dq |= EURO3D_OUTSDRANGE;
01663 }
01664 cpl_table_set_int(tsens, "dq", irow, dq);
01665 }
01666 }
01667
01668
01691
01692 static cpl_polynomial *
01693 muse_flux_response_fit(muse_flux_object *aFluxObj,
01694 double aLambda1, double aLambda2,
01695 unsigned int aOrder, double aRSigma, double *aRMSE)
01696 {
01697 cpl_table *tsens = aFluxObj->sensitivity;
01698 cpl_table_select_all(tsens);
01699 cpl_table_and_selected_int(tsens, "dq", CPL_NOT_EQUAL_TO, EURO3D_GOODPIXEL);
01700 cpl_table_and_selected_int(tsens, "dq", CPL_NOT_EQUAL_TO, EURO3D_TELLCOR);
01701 cpl_table_or_selected_double(tsens, "lambda", CPL_LESS_THAN, aLambda1);
01702 cpl_table_or_selected_double(tsens, "lambda", CPL_GREATER_THAN, aLambda2);
01703
01704 cpl_table *tunwanted = cpl_table_extract_selected(tsens);
01705 cpl_table_erase_selected(tsens);
01706 muse_flux_response_dump_sensitivity(aFluxObj, "fitinput");
01707
01708
01709
01710 int nrow = cpl_table_get_nrow(tsens);
01711 cpl_matrix *lambdas = cpl_matrix_new(1, nrow);
01712 cpl_vector *sens = cpl_vector_new(nrow),
01713 *serr = cpl_vector_new(nrow);
01714 memcpy(cpl_matrix_get_data(lambdas),
01715 cpl_table_get_data_double_const(tsens, "lambda"), nrow*sizeof(double));
01716 memcpy(cpl_vector_get_data(sens),
01717 cpl_table_get_data_double_const(tsens, "sens"), nrow*sizeof(double));
01718 memcpy(cpl_vector_get_data(serr),
01719 cpl_table_get_data_double_const(tsens, "serr"), nrow*sizeof(double));
01720
01721
01722 double chisq, mse;
01723 cpl_polynomial *fit = muse_utils_iterate_fit_polynomial(lambdas, sens, serr,
01724 tsens, aOrder, aRSigma,
01725 &mse, &chisq);
01726 int nout = cpl_vector_get_size(sens);
01727 #if 0
01728 cpl_msg_debug(__func__, "transferred %d entries (%.3f...%.3f) for the "
01729 "order %u fit, %d entries are left, RMS %f", nrow, aLambda1,
01730 aLambda2, aOrder, nout, sqrt(mse));
01731 #endif
01732 cpl_matrix_delete(lambdas);
01733 cpl_vector_delete(sens);
01734 cpl_vector_delete(serr);
01735 if (aRMSE) {
01736 *aRMSE = mse / (nout - aOrder - 1);
01737 }
01738
01739
01740 cpl_table_insert(tsens, tunwanted, nout);
01741 cpl_table_delete(tunwanted);
01742 return fit;
01743 }
01744
01745
01764
01765 static void
01766 muse_flux_response_telluric(muse_flux_object *aFluxObj, double aAirmass)
01767 {
01768 cpl_table *tsens = aFluxObj->sensitivity,
01769 *tb = aFluxObj->tellbands;
01770 cpl_table_new_column(tsens, "sens_orig", CPL_TYPE_DOUBLE);
01771 cpl_table_new_column(tsens, "serr_orig", CPL_TYPE_DOUBLE);
01772 cpl_table_new_column(tsens, "tellcor", CPL_TYPE_DOUBLE);
01773 unsigned int k, nk = cpl_table_get_nrow(tb);
01774 for (k = 0; k < nk; k++) {
01775 double lmin = cpl_table_get_double(tb, "lmin", k, NULL),
01776 lmax = cpl_table_get_double(tb, "lmax", k, NULL),
01777 bgmin = cpl_table_get_double(tb, "bgmin", k, NULL),
01778 bgmax = cpl_table_get_double(tb, "bgmax", k, NULL),
01779 datamin = cpl_table_get_column_min(tsens, "lambda"),
01780 datamax = cpl_table_get_column_max(tsens, "lambda");
01781 cpl_boolean extrapolate = CPL_FALSE;
01782 if (bgmax < lmax || datamax < lmax) {
01783 extrapolate = CPL_TRUE;
01784 }
01785 if (bgmin > lmin || datamin > lmin) {
01786 extrapolate = CPL_TRUE;
01787 }
01788 if (datamin > lmax || datamax < lmin) {
01789
01790 cpl_msg_warning(__func__, "Telluric region %u (range %.2f...%.2f, "
01791 "reference region %.2f...%.2f) outside data range "
01792 "(%.2f..%.2f)!", k + 1, lmin, lmax, bgmin, bgmax,
01793 datamin, datamax);
01794 continue;
01795 }
01796
01797 unsigned int order = extrapolate ? 1 : 2;
01798
01799
01800 double rmse = 0.;
01801 cpl_errorstate state = cpl_errorstate_get();
01802 cpl_polynomial *fit = muse_flux_response_fit(aFluxObj, bgmin, bgmax,
01803 order, 3., &rmse);
01804 if (!cpl_errorstate_is_equal(state) || !fit) {
01805 cpl_msg_warning(__func__, "Telluric region %u (range %.2f...%.2f, "
01806 "reference region %.2f...%.2f) could not be fitted!",
01807 k + 1, lmin, lmax, bgmin, bgmax);
01808 cpl_errorstate_set(state);
01809 cpl_polynomial_delete(fit);
01810 continue;
01811 }
01812 cpl_msg_debug(__func__, "Telluric region %u: %.2f...%.2f, reference region "
01813 "%.2f...%.2f", k + 1, lmin, lmax, bgmin, bgmax);
01814 #if 0
01815 cpl_polynomial_dump(fit, stdout);
01816 fflush(stdout);
01817 #endif
01818 int irow, nrow = cpl_table_get_nrow(tsens);
01819 for (irow = 0; irow < nrow; irow++) {
01820 double lambda = cpl_table_get_double(tsens, "lambda", irow, NULL);
01821 if (lambda >= lmin && lambda <= lmax &&
01822 (unsigned int)cpl_table_get_int(tsens, "dq", irow, NULL)
01823 == EURO3D_TELLURIC) {
01824 double origval = cpl_table_get_double(tsens, "sens", irow, NULL),
01825 origerr = cpl_table_get_double(tsens, "serr", irow, NULL),
01826 interpval = cpl_polynomial_eval_1d(fit, lambda, NULL);
01827 cpl_table_set_int(tsens, "dq", irow, EURO3D_TELLCOR);
01828 cpl_table_set_double(tsens, "sens_orig", irow, origval);
01829 cpl_table_set_double(tsens, "sens", irow, interpval);
01830
01831
01832 cpl_table_set_double(tsens, "serr_orig", irow, origerr);
01833 cpl_table_set_double(tsens, "serr", irow, sqrt(origerr*origerr + rmse));
01834
01835 if (interpval > origval) {
01836
01837 double ftelluric = pow(10, -0.4 * (interpval - origval));
01838 cpl_table_set(tsens, "tellcor", irow, ftelluric);
01839 } else {
01840 cpl_table_set_double(tsens, "tellcor", irow, 1.);
01841 }
01842 }
01843 }
01844 cpl_polynomial_delete(fit);
01845 }
01846
01847 cpl_table_power_column(tsens, "tellcor", 1. / aAirmass);
01848 }
01849
01850
01867
01868 static void
01869 muse_flux_response_extrapolate(muse_flux_object *aFluxObj, double aDistance)
01870 {
01871 cpl_table *tsens = aFluxObj->sensitivity;
01872 cpl_propertylist *order = cpl_propertylist_new();
01873 cpl_propertylist_append_bool(order, "lambda", CPL_FALSE);
01874 cpl_table_sort(tsens, order);
01875
01876 int nrow = cpl_table_get_nrow(tsens);
01877
01878 double lambda1 = cpl_table_get_double(tsens, "lambda", 0, NULL),
01879 serr1 = cpl_table_get_double(tsens, "serr", 0, NULL);
01880 unsigned int dq = cpl_table_get_int(tsens, "dq", 0, NULL);
01881 int irow = 0;
01882 while (dq != EURO3D_GOODPIXEL && dq != EURO3D_TELLCOR) {
01883 lambda1 = cpl_table_get_double(tsens, "lambda", ++irow, NULL);
01884 serr1 = cpl_table_get_double(tsens, "serr", irow, NULL);
01885 dq = cpl_table_get_int(tsens, "dq", irow, NULL);
01886 }
01887 double lambda2 = cpl_table_get_double(tsens, "lambda", nrow - 1, NULL),
01888 serr2 = cpl_table_get_double(tsens, "serr", nrow - 1, NULL);
01889 dq = cpl_table_get_int(tsens, "dq", nrow - 1, NULL);
01890 irow = nrow - 1;
01891 while (dq != EURO3D_GOODPIXEL && dq != EURO3D_TELLCOR) {
01892 lambda2 = cpl_table_get_double(tsens, "lambda", --irow, NULL);
01893 serr2 = cpl_table_get_double(tsens, "serr", irow, NULL);
01894 dq = cpl_table_get_int(tsens, "dq", irow, NULL);
01895 }
01896 cpl_polynomial *fit1 = muse_flux_response_fit(aFluxObj, lambda1,
01897 lambda1 + aDistance, 1, 5., NULL),
01898 *fit2 = muse_flux_response_fit(aFluxObj, lambda2 - aDistance,
01899 lambda2, 1, 5., NULL);
01900 nrow = cpl_table_get_nrow(tsens);
01901 double d1 = (lambda1 - kMuseLambdaMinX) / 100.,
01902 d2 = (kMuseLambdaMaxX - lambda2) / 100.;
01903 cpl_table_set_size(tsens, nrow + 200);
01904 irow = nrow;
01905 double l;
01906 if (fit1) {
01907 for (l = kMuseLambdaMinX; l <= lambda1 - d1; l += d1) {
01908 double sens = cpl_polynomial_eval_1d(fit1, l, NULL),
01909
01910 serr = 2. * (lambda1 - l) / 50. * serr1 + serr1;
01911 if (sens <= 0) {
01912 cpl_table_set_invalid(tsens, "lambda", irow++);
01913 cpl_msg_debug(__func__, "invalid blueward extrapolation: %.3f %f +/- %f",
01914 l, sens, serr);
01915 continue;
01916 }
01917 cpl_table_set_double(tsens, "lambda", irow, l);
01918 cpl_table_set_double(tsens, "sens", irow, sens);
01919 cpl_table_set_double(tsens, "serr", irow, serr);
01920 cpl_table_set_int(tsens, "dq", irow++, (int)EURO3D_OUTSDRANGE);
01921 }
01922 cpl_msg_debug(__func__, "Extrapolated blue end: %.1f...%.1f Angstrom (using"
01923 " data from %.1f...%.1f Angstrom)", kMuseLambdaMinX,
01924 lambda1 - d1, lambda1, lambda1 + aDistance);
01925 }
01926 if (fit2) {
01927 for (l = lambda2 + d2; fit2 && l <= kMuseLambdaMaxX; l += d2) {
01928 double sens = cpl_polynomial_eval_1d(fit2, l, NULL),
01929 serr = 2. * (l - lambda2) / 50. * serr2 + serr2;
01930 if (sens <= 0) {
01931 cpl_table_set_invalid(tsens, "lambda", irow++);
01932 cpl_msg_debug(__func__, "invalid redward extrapolation: %.3f %f +/- %f",
01933 l, sens, serr);
01934 continue;
01935 }
01936 cpl_table_set_double(tsens, "lambda", irow, l);
01937 cpl_table_set_double(tsens, "sens", irow, sens);
01938 cpl_table_set_double(tsens, "serr", irow, serr);
01939 cpl_table_set_int(tsens, "dq", irow++, (int)EURO3D_OUTSDRANGE);
01940 }
01941 cpl_msg_debug(__func__, "Extrapolated red end: %.1f...%.1f Angstrom (using "
01942 "data from %.1f...%.1f Angstrom)", lambda2 + d2,
01943 kMuseLambdaMaxX, lambda2 - aDistance, lambda2);
01944 }
01945 #if 0
01946 cpl_polynomial_dump(fit1, stdout);
01947 cpl_polynomial_dump(fit2, stdout);
01948 fflush(stdout);
01949 #endif
01950 cpl_polynomial_delete(fit1);
01951 cpl_polynomial_delete(fit2);
01952
01953 cpl_table_select_all(tsens);
01954 cpl_table_and_selected_invalid(tsens, "sens");
01955 cpl_table_erase_selected(tsens);
01956
01957 cpl_table_sort(tsens, order);
01958 cpl_propertylist_delete(order);
01959 }
01960
01961
02013
02014 cpl_error_code
02015 muse_flux_response_compute(muse_flux_object *aFluxObj,
02016 muse_flux_selection_type aSelect, double aAirmass,
02017 const cpl_table *aReference,
02018 const cpl_table *aTellBands,
02019 const cpl_table *aExtinct)
02020 {
02021 cpl_ensure_code(aFluxObj && aFluxObj->intimage && aReference,
02022 CPL_ERROR_NULL_INPUT);
02023 cpl_ensure_code(aAirmass >= 1., CPL_ERROR_ILLEGAL_INPUT);
02024 if (!aExtinct) {
02025 cpl_msg_warning(__func__, "Extinction table not given!");
02026 }
02027 if (aSelect == MUSE_FLUX_SELECT_NEAREST &&
02028 (!isfinite(aFluxObj->raref) || !isfinite(aFluxObj->decref))) {
02029 cpl_msg_warning(__func__, "Reference position %f,%f contains infinite "
02030 "values, using flux to select star!", aFluxObj->raref,
02031 aFluxObj->decref);
02032 aSelect = MUSE_FLUX_SELECT_BRIGHTEST;
02033 }
02034 muse_flux_response_set_telluric_bands(aFluxObj, aTellBands);
02035
02036 int nobjects = cpl_image_get_size_y(aFluxObj->intimage->data);
02037 const char *bunit = muse_pfits_get_bunit(aFluxObj->intimage->header);
02038
02039 double flux = 0., dmin = DBL_MAX;
02040 int n, nstar = 1, nstardist = 1;
02041 for (n = 1; n <= nobjects; n++) {
02042 char kw[KEYWORD_LENGTH];
02043 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_RA, n);
02044 double ra = cpl_propertylist_get_double(aFluxObj->intimage->header, kw);
02045 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_DEC, n);
02046 double dec = cpl_propertylist_get_double(aFluxObj->intimage->header, kw),
02047 dthis = muse_astro_angular_distance(ra, dec, aFluxObj->raref,
02048 aFluxObj->decref);
02049 cpl_msg_debug(__func__, "distance(%d) = %f arcsec", n, dthis * 3600.);
02050 if (fabs(dthis) < dmin) {
02051 dmin = dthis;
02052 nstardist = n;
02053 }
02054 snprintf(kw, KEYWORD_LENGTH, MUSE_HDR_FLUX_OBJn_FLUX, n);
02055 double this = cpl_propertylist_get_double(aFluxObj->intimage->header, kw);
02056 cpl_msg_debug(__func__, "flux(%d) = %e %s", n, this, bunit);
02057 if (this > flux) {
02058 flux = this;
02059 nstar = n;
02060 }
02061 }
02062 int nselected;
02063 char *outstring = NULL;
02064 if (aSelect == MUSE_FLUX_SELECT_BRIGHTEST) {
02065 outstring = cpl_sprintf("Selected the brightest star (%d of %d; %.3e %s)"
02066 " as reference source", nstar, nobjects, flux, bunit);
02067 nselected = nstar;
02068 } else {
02069 outstring = cpl_sprintf("Selected the nearest star (%d of %d; %.2f arcsec) "
02070 "as reference source", nstar, nobjects, dmin*3600.);
02071 nselected = nstardist;
02072 }
02073 cpl_msg_info(__func__, "%s", outstring);
02074 cpl_free(outstring);
02075
02076 muse_flux_response_sensitivity(aFluxObj,
02077 nselected - 1, aReference,
02078 aAirmass, aExtinct);
02079 muse_flux_response_dump_sensitivity(aFluxObj, "initial");
02080
02081 muse_flux_response_mark_questionable(aFluxObj);
02082 muse_flux_response_dump_sensitivity(aFluxObj, "intermediate");
02083
02084 cpl_table_select_all(aFluxObj->sensitivity);
02085 cpl_table_and_selected_int(aFluxObj->sensitivity, "dq", CPL_EQUAL_TO,
02086 (int)EURO3D_OUTSDRANGE);
02087 cpl_table_erase_selected(aFluxObj->sensitivity);
02088 muse_flux_response_dump_sensitivity(aFluxObj, "intercut");
02089
02090
02091
02092 muse_flux_response_telluric(aFluxObj, aAirmass);
02093 muse_flux_response_dump_sensitivity(aFluxObj, "interpolated");
02094
02095 muse_flux_response_extrapolate(aFluxObj, 150.);
02096 muse_flux_response_dump_sensitivity(aFluxObj, "extrapolated");
02097
02098 return CPL_ERROR_NONE;
02099 }
02100
02101
02110
02111 const muse_cpltable_def muse_flux_responsetable_def[] = {
02112 { "lambda", CPL_TYPE_DOUBLE, "Angstrom", "%7.2f", "wavelength", CPL_TRUE },
02113 { "response", CPL_TYPE_DOUBLE,
02114 "2.5*log10((count/s/Angstrom)/(erg/s/cm**2/Angstrom))", "%.4e",
02115 "instrument response derived from standard star", CPL_TRUE },
02116 { "resperr", CPL_TYPE_DOUBLE,
02117 "2.5*log10((count/s/Angstrom)/(erg/s/cm**2/Angstrom))", "%.4e",
02118 "instrument response error derived from standard star", CPL_TRUE },
02119 { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
02120 };
02121
02122
02131
02132 const muse_cpltable_def muse_flux_tellurictable_def[] = {
02133 { "lambda", CPL_TYPE_DOUBLE, "Angstrom", "%7.2f", "wavelength", CPL_TRUE },
02134 { "ftelluric", CPL_TYPE_DOUBLE, "", "%.5f",
02135 "the telluric correction factor, normalized to an airmass of 1", CPL_TRUE },
02136 { "ftellerr", CPL_TYPE_DOUBLE, "", "%.5f",
02137 "the error of the telluric correction factor", CPL_TRUE },
02138 { NULL, 0, NULL, NULL, NULL, CPL_FALSE }
02139 };
02140
02141
02160
02161 static void
02162 muse_flux_get_response_table_smooth(cpl_table *aResp, double aHalfwidth,
02163 double aLambdaMin, double aLambdaMax,
02164 cpl_boolean aAverage)
02165 {
02166
02167 cpl_table_duplicate_column(aResp, "sens", aResp, "response");
02168 cpl_table_duplicate_column(aResp, "serr", aResp, "resperr");
02169
02170
02171 cpl_table_select_all(aResp);
02172 cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_LESS_THAN, aLambdaMin);
02173 cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_GREATER_THAN, aLambdaMax);
02174
02175 cpl_boolean sym = cpl_table_count_selected(aResp) < cpl_table_get_nrow(aResp);
02176 cpl_msg_debug(__func__, "%s smoothing response +/- %.3f Angstrom between %.3f "
02177 "and %.3f Angstrom", sym ? "symmetrical" : "", aHalfwidth,
02178 aLambdaMin, aLambdaMax);
02179
02180
02181 int i, n = cpl_table_get_nrow(aResp);
02182 for (i = 0; i < n; i++) {
02183 if (!cpl_table_is_selected(aResp, i)) {
02184 continue;
02185 }
02186 double lambda = cpl_table_get_double(aResp, "lambda", i, NULL);
02187 int j = i, j1 = i, j2 = i;
02188
02189 while (--j > 0 && cpl_table_is_selected(aResp, j) &&
02190 lambda - cpl_table_get_double(aResp, "lambda", j, NULL) <= aHalfwidth) {
02191 j1 = j;
02192 }
02193 j = i;
02194 while (++j < n && cpl_table_is_selected(aResp, j) &&
02195 cpl_table_get_double(aResp, "lambda", j, NULL) - lambda <= aHalfwidth) {
02196 j2 = j;
02197 }
02198 if (sym) {
02199 int jd1 = i - j1,
02200 jd2 = j2 - i;
02201 if (jd1 < jd2) {
02202 j2 = i + jd1;
02203 } else {
02204 j1 = i - jd2;
02205 }
02206 }
02207
02208 double *sens = cpl_table_get_data_double(aResp, "sens"),
02209 *serr = cpl_table_get_data_double(aResp, "serr");
02210 cpl_vector *v = cpl_vector_wrap(j2 - j1 + 1, sens + j1),
02211 *ve = cpl_vector_wrap(j2 - j1 + 1, serr + j1);
02212 if (aAverage) {
02213
02214 double mean = cpl_vector_get_mean(v),
02215 stdev = j2 == j1 ? 0. : cpl_vector_get_stdev(v),
02216 rerr = cpl_table_get_double(aResp, "resperr", i, NULL);
02217 cpl_table_set_double(aResp, "response", i, mean);
02218 cpl_table_set_double(aResp, "resperr", i, sqrt(rerr*rerr + stdev*stdev));
02219 } else {
02220
02221 double median = cpl_vector_get_median_const(v),
02222 mdev = muse_cplvector_get_adev_const(v, median),
02223 mederr = cpl_vector_get_median_const(ve);
02224 if (j2 == j1) {
02225 mdev = cpl_table_get_double(aResp, "serr", j1, NULL);
02226 }
02227 if (mdev < mederr) {
02228 mdev = mederr;
02229 }
02230 #if 0
02231 cpl_msg_debug(__func__, "%d %.3f %d...%d --> %f +/- %f", i, lambda, j1, j2,
02232 median, mdev);
02233 #endif
02234 cpl_table_set_double(aResp, "response", i, median);
02235 cpl_table_set_double(aResp, "resperr", i, mdev);
02236 }
02237 cpl_vector_unwrap(v);
02238 cpl_vector_unwrap(ve);
02239 }
02240
02241
02242 cpl_table_erase_column(aResp, "sens");
02243 cpl_table_erase_column(aResp, "serr");
02244 }
02245
02246
02263
02264 static unsigned int
02265 muse_flux_get_response_table_collect_points(const cpl_table *aTable,
02266 double aLambda, double aLDist,
02267 cpl_matrix *aPos, cpl_vector *aVal,
02268 cpl_vector *aErr)
02269 {
02270 unsigned int np = 0;
02271 int irow, nrow = cpl_table_get_nrow(aTable);
02272 for (irow = 0; irow < nrow; irow++) {
02273 double lambda = cpl_table_get(aTable, "lambda", irow, NULL);
02274 if (lambda < aLambda - aLDist || lambda > aLambda + aLDist) {
02275 continue;
02276 }
02277 cpl_matrix_set(aPos, 0, np, lambda);
02278 cpl_vector_set(aVal, np, cpl_table_get(aTable, "sens", irow, NULL));
02279 cpl_vector_set(aErr, np, cpl_table_get(aTable, "serr", irow, NULL));
02280 np++;
02281 }
02282 cpl_matrix_set_size(aPos, 1, np);
02283 cpl_vector_set_size(aVal, np);
02284 cpl_vector_set_size(aErr, np);
02285 return np;
02286 }
02287
02288
02313
02314 static void
02315 muse_flux_get_response_table_piecewise_poly(cpl_table *aResp, double aDMin,
02316 double aDMax, float aRSigma)
02317 {
02318
02319 cpl_table_duplicate_column(aResp, "sens", aResp, "response");
02320 cpl_table_duplicate_column(aResp, "serr", aResp, "resperr");
02321
02322
02323 unsigned int npold = 0, njumps = 0;
02324 double ldistold = -1,
02325 lambdaold = -1;
02326 cpl_array *jumppos = cpl_array_new(0, CPL_TYPE_DOUBLE),
02327 *jumplen = cpl_array_new(0, CPL_TYPE_DOUBLE);
02328
02329
02330 int irow, nrow = cpl_table_get_nrow(aResp);
02331 for (irow = 0; irow < nrow; irow++) {
02332 double lambda = cpl_table_get(aResp, "lambda", irow, NULL);
02333
02334
02335
02336 double ldist = aDMax;
02337 if ((lambda >= 5700 && lambda < 6200) || (lambda >= 6900 && lambda < 7200)) {
02338 ldist = aDMin;
02339 }
02340
02341 cpl_matrix *pos = cpl_matrix_new(1, nrow);
02342 cpl_vector *val = cpl_vector_new(nrow),
02343 *err = cpl_vector_new(nrow);
02344 unsigned int np = muse_flux_get_response_table_collect_points(aResp,
02345 lambda, ldist,
02346 pos, val, err);
02347 if (ldistold < 0) {
02348 npold = np;
02349 ldistold = ldist;
02350 lambdaold = lambda;
02351 }
02352 #if 0
02353 printf("%f Angstrom %u points (%u, %.3f):\n", lambda, np, npold, (double)np / npold - 1.);
02354
02355 fflush(stdout);
02356 #endif
02357
02358 if (np > 10 && fabs((double)np / npold - 1.) > 0.1) {
02359 cpl_msg_debug(__func__, "possible jump, changed at lambda %.3f (%u -> %u, "
02360 "%.3f -> %.3f)", lambda, npold, np, ldistold, ldist);
02361 cpl_array_set_size(jumppos, ++njumps);
02362 cpl_array_set_size(jumplen, njumps);
02363 cpl_array_set_double(jumppos, njumps - 1, (lambdaold + lambda) / 2.);
02364 cpl_array_set_double(jumplen, njumps - 1, (ldistold + ldist) / 2.);
02365 }
02366
02367
02368 unsigned int order = np > 3 ? 3 : np - 1;
02369 double mse;
02370
02371 cpl_polynomial *poly = muse_utils_iterate_fit_polynomial(pos, val, err,
02372 NULL, order, aRSigma,
02373 &mse, NULL);
02374 #if 0
02375 if (fabs(lambda - 40861.3) < 1.) {
02376 cpl_vector *res = cpl_vector_new(cpl_vector_get_size(val));
02377 cpl_vector_fill_polynomial_fit_residual(res, val, NULL, poly, pos, NULL);
02378 double rms = sqrt(cpl_vector_product(res, res) / cpl_vector_get_size(res));
02379 cpl_msg_debug(__func__, "lambda %f rms %f (%u/%d points)", lambda, rms,
02380 (unsigned)cpl_vector_get_size(val), np);
02381
02382 cpl_plot_vector(NULL, NULL, NULL, val);
02383 cpl_plot_vector(NULL, NULL, NULL, res);
02384 cpl_vector_delete(res);
02385 }
02386 #endif
02387 cpl_matrix_delete(pos);
02388 cpl_vector_delete(val);
02389 cpl_vector_delete(err);
02390 double resp = cpl_polynomial_eval_1d(poly, lambda, NULL);
02391 cpl_polynomial_delete(poly);
02392 cpl_table_set(aResp, "lambda", irow, lambda);
02393 cpl_table_set(aResp, "response", irow, resp);
02394 double serr = cpl_table_get(aResp, "serr", irow, NULL);
02395 cpl_table_set(aResp, "resperr", irow, sqrt(mse + serr*serr));
02396
02397 npold = np;
02398 ldistold = ldist;
02399 lambdaold = lambda;
02400 }
02401
02402
02403 cpl_table_erase_column(aResp, "sens");
02404 cpl_table_erase_column(aResp, "serr");
02405
02406 #if 0
02407 printf("jumppos (%u):\n", njumps);
02408 cpl_array_dump(jumppos, 0, 10000, stdout);
02409 printf("jumplen:\n");
02410 cpl_array_dump(jumplen, 0, 10000, stdout);
02411 fflush(stdout);
02412 #endif
02413
02414
02415
02416 unsigned int iarr;
02417 for (iarr = 0; iarr < njumps; iarr++) {
02418 double lambda = cpl_array_get_double(jumppos, iarr, NULL),
02419 ldist = cpl_array_get_double(jumplen, iarr, NULL) / 2;
02420
02421 cpl_table_select_all(aResp);
02422 cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_LESS_THAN, lambda - 5.);
02423 cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_GREATER_THAN, lambda + 5.);
02424 int nsel = cpl_table_count_selected(aResp);
02425 if (nsel <= 1) {
02426 cpl_msg_debug(__func__, "Only %d points near jump around %.1f Angstrom, "
02427 "not doing anything", nsel, lambda);
02428 continue;
02429 }
02430 cpl_table *xresp = cpl_table_extract_selected(aResp);
02431 double stdev = cpl_table_get_column_stdev(xresp, "response");
02432 cpl_table_dump(xresp, 0, nsel, stdout);
02433 fflush(stdout);
02434 cpl_table_delete(xresp);
02435 if (stdev < 0.001) {
02436 cpl_msg_debug(__func__, "%d points near jump around %.1f Angstrom, stdev "
02437 "only %f, not doing anything", nsel, lambda, stdev);
02438 continue;
02439 }
02440 cpl_msg_debug(__func__, "%d points near jump around %.1f Angstrom, stdev "
02441 "%f, erasing the region", nsel, lambda, stdev);
02442
02443
02444
02445 cpl_table_select_all(aResp);
02446 cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_LESS_THAN, lambda - ldist);
02447 cpl_table_and_selected_double(aResp, "lambda", CPL_NOT_GREATER_THAN, lambda + ldist);
02448 cpl_table_erase_selected(aResp);
02449 }
02450 cpl_array_delete(jumppos);
02451 cpl_array_delete(jumplen);
02452 }
02453
02454
02479
02480 cpl_error_code
02481 muse_flux_get_response_table(muse_flux_object *aFluxObj,
02482 muse_flux_smooth_type aSmooth)
02483 {
02484 cpl_ensure_code(aFluxObj && aFluxObj->sensitivity, CPL_ERROR_NULL_INPUT);
02485 cpl_ensure_code(aSmooth <= MUSE_FLUX_SMOOTH_PPOLY, CPL_ERROR_ILLEGAL_INPUT);
02486
02487 int nrow = cpl_table_get_nrow(aFluxObj->sensitivity);
02488 cpl_table *resp = muse_cpltable_new(muse_flux_responsetable_def, nrow);
02489
02490 const double *lambdas = cpl_table_get_data_double_const(aFluxObj->sensitivity,
02491 "lambda"),
02492 *sens = cpl_table_get_data_double_const(aFluxObj->sensitivity,
02493 "sens"),
02494 *serr = cpl_table_get_data_double_const(aFluxObj->sensitivity,
02495 "serr");
02496 cpl_table_copy_data_double(resp, "lambda", lambdas);
02497 cpl_table_copy_data_double(resp, "response", sens);
02498 cpl_table_copy_data_double(resp, "resperr", serr);
02499
02500
02501 if (aSmooth == MUSE_FLUX_SMOOTH_MEDIAN) {
02502 cpl_msg_info(__func__, "Smoothing response curve with median filter");
02503
02504 muse_flux_get_response_table_smooth(resp, 15., 0., 20000., CPL_FALSE);
02505 } else if (aSmooth == MUSE_FLUX_SMOOTH_PPOLY) {
02506 cpl_msg_info(__func__, "Smoothing response curve with piecewise polynomial");
02507
02508 muse_flux_get_response_table_piecewise_poly(resp, 150., 150., 3.);
02509
02510
02511 muse_flux_get_response_table_smooth(resp, 15., 0., 20000., CPL_TRUE);
02512 } else {
02513 cpl_msg_warning(__func__, "NOT smoothing the response curve at all!");
02514 }
02515
02516
02517 aFluxObj->response = resp;
02518 return CPL_ERROR_NONE;
02519 }
02520
02521
02540
02541 cpl_error_code
02542 muse_flux_get_telluric_table(muse_flux_object *aFluxObj)
02543 {
02544 cpl_ensure_code(aFluxObj && aFluxObj->sensitivity, CPL_ERROR_NULL_INPUT);
02545 cpl_table *tsens = aFluxObj->sensitivity;
02546 int nrow = cpl_table_get_nrow(tsens);
02547 cpl_table *tell = muse_cpltable_new(muse_flux_tellurictable_def, nrow);
02548
02549 cpl_table_fill_column_window_double(tell, "lambda", 0, nrow, 0);
02550 cpl_table_copy_data_double(tell, "lambda",
02551 cpl_table_get_data_double_const(tsens, "lambda"));
02552 cpl_table_fill_column_window_double(tell, "ftelluric", 0, nrow, 0);
02553 cpl_table_copy_data_double(tell, "ftelluric",
02554 cpl_table_get_data_double_const(tsens, "tellcor"));
02555
02556 #define TELL_MAX_ERR 0.1
02557 #define TELL_MIN_ERR 1e-4
02558 cpl_table_fill_column_window_double(tell, "ftellerr", 0, nrow, TELL_MAX_ERR);
02559
02560
02561
02562 cpl_table_duplicate_column(tell, "tellcor", tsens, "tellcor");
02563
02564 cpl_table_unselect_all(tell);
02565 int irow;
02566 for (irow = 0; irow < nrow; irow++) {
02567 int err;
02568 cpl_table_get_double(tell, "tellcor", irow, &err);
02569 if (err == 0) {
02570 continue;
02571 }
02572
02573 cpl_errorstate state = cpl_errorstate_get();
02574 double ftellcor = cpl_table_get_double(tell, "tellcor", irow - 1, &err);
02575 if (!cpl_errorstate_is_equal(state)) {
02576 cpl_errorstate_set(state);
02577 }
02578 if (err == 0 && ftellcor != 1.) {
02579 cpl_table_set_double(tell, "ftelluric", irow, 1.);
02580 continue;
02581 }
02582
02583 state = cpl_errorstate_get();
02584 ftellcor = cpl_table_get_double(tell, "tellcor", irow + 1, &err);
02585 if (!cpl_errorstate_is_equal(state)) {
02586 cpl_errorstate_set(state);
02587 }
02588 if (err == 0 && ftellcor != 1.) {
02589 cpl_table_set_double(tell, "ftelluric", irow, 1.);
02590 continue;
02591 }
02592 cpl_table_select_row(tell, irow);
02593 }
02594 cpl_table_erase_selected(tell);
02595 cpl_table_erase_column(tell, "tellcor");
02596
02597
02598 nrow = cpl_table_get_nrow(tell);
02599 for (irow = 0; irow < nrow; irow++) {
02600 int err;
02601 double dftell = 1. - cpl_table_get_double(tell, "ftelluric", irow, &err),
02602 ftellerr = cpl_table_get_double(tell, "ftellerr", irow, &err);
02603 if (ftellerr > dftell) {
02604 ftellerr = fmax(dftell, TELL_MIN_ERR);
02605 cpl_table_set_double(tell, "ftellerr", irow, ftellerr);
02606 }
02607 }
02608
02609 aFluxObj->telluric = tell;
02610 return CPL_ERROR_NONE;
02611 }
02612
02613
02650
02651 cpl_error_code
02652 muse_flux_calibrate(muse_pixtable *aPixtable, const cpl_table *aResponse,
02653 const cpl_table *aExtinction, const cpl_table *aTelluric)
02654 {
02655 cpl_ensure_code(aPixtable && aPixtable->header && aResponse,
02656 CPL_ERROR_NULL_INPUT);
02657 const char *unitdata = cpl_table_get_column_unit(aPixtable->table,
02658 MUSE_PIXTABLE_DATA);
02659 cpl_ensure_code(unitdata && !strncmp(unitdata, "count", 6),
02660 CPL_ERROR_INCOMPATIBLE_INPUT);
02661
02662 if (!aExtinction) {
02663 cpl_msg_warning(__func__, "%s missing!", MUSE_TAG_EXTINCT_TABLE);
02664 }
02665
02666 double exptime = muse_pfits_get_exptime(aPixtable->header);
02667 if (exptime <= 0.) {
02668 cpl_msg_warning(__func__, "Non-positive EXPTIME, not doing flux calibration!");
02669 return CPL_ERROR_ILLEGAL_INPUT;
02670 }
02671 double airmass = muse_astro_airmass(aPixtable->header);
02672 if (airmass < 0) {
02673 cpl_msg_warning(__func__, "Airmass unknown, not doing extinction "
02674 "correction: %s", cpl_error_get_message());
02675
02676 airmass = 0.;
02677 }
02678
02679 cpl_table *telluric = NULL;
02680 if (aTelluric) {
02681
02682 telluric = cpl_table_duplicate(aTelluric);
02683
02684
02685 cpl_table_power_column(telluric, "ftelluric", -airmass);
02686 }
02687
02688 cpl_msg_info(__func__, "Starting flux calibration (exptime=%.2fs, airmass=%.4f),"
02689 " %s telluric correction", exptime, airmass,
02690 aTelluric ? "with" : "without ("MUSE_TAG_STD_TELLURIC" not given)");
02691 float *lambda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA),
02692 *data = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
02693 *stat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
02694 cpl_size i, nrow = muse_pixtable_get_nrow(aPixtable);
02695 #pragma omp parallel for default(none) \
02696 shared(aExtinction, aResponse, airmass, data, exptime, lambda, nrow, \
02697 stat, telluric)
02698 for (i = 0; i < nrow; i++) {
02699
02700 double ddata = data[i], dstat = stat[i];
02701
02702
02703 if (aExtinction) {
02704 double fext = pow(10., 0.4 * airmass
02705 * muse_flux_response_interpolate(aExtinction,
02706 lambda[i], NULL,
02707 MUSE_FLUX_RESP_EXTINCT));
02708 #if 0
02709 printf("%f, data/stat = %f/%f -> ", fext, ddata, dstat);
02710 #endif
02711 ddata *= fext;
02712 dstat *= (fext * fext);
02713 #if 0
02714 printf(" --> %f/%f\n", ddata, dstat), fflush(NULL);
02715 #endif
02716 }
02717
02718
02719
02720
02721 double dlambda = kMuseSpectralSamplingA,
02722 dlamerr = 0.02,
02723
02724 rerr, resp = muse_flux_response_interpolate(aResponse, lambda[i],
02725 &rerr,
02726 MUSE_FLUX_RESP_FLUX);
02727
02728 resp = pow(10., 0.4 * resp);
02729
02730
02731 rerr = rerr * CPL_MATH_LN10 * resp / 2.5;
02732 #if 0
02733 printf("%f/%f/%f, %e/%e, data/stat = %e/%e -> ", lambda[i], dlambda, dlamerr, resp, rerr,
02734 ddata, dstat);
02735 #endif
02736 dstat = dstat * pow((1./(resp * exptime * dlambda)), 2)
02737 + pow(ddata * rerr / (resp*resp * exptime * dlambda), 2)
02738 + pow(ddata * dlamerr / (resp * exptime * dlambda*dlambda), 2);
02739 ddata /= (resp * exptime * dlambda);
02740 #if 0
02741 printf("%e/%e\n", ddata, dstat), fflush(NULL);
02742 #endif
02743
02744
02745
02746 ddata *= kMuseFluxUnitFactor;
02747 dstat *= kMuseFluxStatFactor;
02748
02749
02750
02751 if (lambda[i] < kTelluricBands[0][0] || !telluric) {
02752 data[i] = ddata;
02753 stat[i] = dstat;
02754 continue;
02755 }
02756 double terr, tell = muse_flux_response_interpolate(telluric, lambda[i],
02757 &terr,
02758 MUSE_FLUX_TELLURIC);
02759 data[i] = ddata * tell;
02760 stat[i] = tell*tell * dstat + ddata*ddata * terr*terr;
02761 }
02762 cpl_table_delete(telluric);
02763
02764
02765 cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_DATA,
02766 kMuseFluxUnitString);
02767 cpl_table_set_column_unit(aPixtable->table, MUSE_PIXTABLE_STAT,
02768 kMuseFluxStatString);
02769
02770
02771 cpl_propertylist_update_bool(aPixtable->header, MUSE_HDR_PT_FLUXCAL,
02772 CPL_TRUE);
02773 cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_FLUXCAL,
02774 MUSE_HDR_PT_FLUXCAL_COMMENT);
02775 return CPL_ERROR_NONE;
02776 }
02777