00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 #ifdef HAVE_CONFIG_H
00024 #include <config.h>
00025 #endif
00026
00027
00028
00029
00030 #include <cpl.h>
00031 #include <math.h>
00032 #include <string.h>
00033
00034 #include "muse_dar.h"
00035 #include "muse_instrument.h"
00036
00037 #include "muse_astro.h"
00038 #include "muse_combine.h"
00039 #include "muse_pfits.h"
00040 #include "muse_phys.h"
00041 #include "muse_quality.h"
00042 #include "muse_utils.h"
00043 #include "muse_wcs.h"
00044
00045
00046
00047
00048 #define DEBUG_MUSE_DARCORRECT 0
00049 #define DEBUG_MUSE_DARCHECK 0
00050
00051 #define MUSE_DARCHECK_GRID_FITS 0
00052
00053
00054
00055
00059
00060
00063
00137
00138 cpl_error_code
00139 muse_dar_correct(muse_pixtable *aPixtable, double aLambdaRef)
00140 {
00141 cpl_ensure_code(aPixtable && aPixtable->header, CPL_ERROR_NULL_INPUT);
00142 if (aLambdaRef > 22000. || aLambdaRef < 3500.) {
00143 cpl_msg_info(__func__, "Reference lambda %.3f Angstrom: skipping DAR "
00144 "correction", aLambdaRef);
00145 cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_NAME, -1.);
00146 cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_NAME,
00147 MUSE_HDR_PT_DAR_C_SKIP);
00148 return CPL_ERROR_NONE;
00149 }
00150 if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_NAME) &&
00151 cpl_propertylist_get_double(aPixtable->header, MUSE_HDR_PT_DAR_NAME) > 0.) {
00152 cpl_msg_info(__func__, "pixel table already corrected: skipping DAR "
00153 "correction");
00154 return CPL_ERROR_NONE;
00155 }
00156
00157
00158
00159
00160 double airmass = muse_astro_airmass(aPixtable->header);
00161 cpl_ensure_code(airmass >= 1., cpl_error_get_code());
00162
00163 double z = acos(1./airmass);
00164
00165
00166 cpl_errorstate prestate = cpl_errorstate_get();
00167 double parang = muse_astro_parangle(aPixtable->header);
00168 cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
00169
00170
00171 prestate = cpl_errorstate_get();
00172 double posang = muse_astro_posangle(aPixtable->header);
00173 cpl_ensure_code(cpl_errorstate_is_equal(prestate), cpl_error_get_code());
00174
00175 cpl_boolean isWFM = muse_pfits_get_mode(aPixtable->header)
00176 <= MUSE_MODE_WFM_AO_N;
00177 if (!isWFM) {
00178 cpl_msg_warning(__func__, "For NFM instrument mode DAR correction should "
00179 "not be needed!");
00180 }
00181
00182
00183 prestate = cpl_errorstate_get();
00184 double temp = muse_pfits_get_temp(aPixtable->header);
00185 if (cpl_errorstate_is_equal(prestate)) {
00186 temp += 273.15;
00187 } else {
00188 cpl_msg_warning(__func__, "FITS header does not contain temperature: %s",
00189 cpl_error_get_message());
00190 temp = 11.5 + 273.15;
00191 }
00192 prestate = cpl_errorstate_get();
00193 double rhum = muse_pfits_get_rhum(aPixtable->header);
00194 if (!cpl_errorstate_is_equal(prestate)) {
00195 cpl_msg_warning(__func__, "FITS header does not contain relative humidity: %s",
00196 cpl_error_get_message());
00197 rhum = 14.5;
00198 }
00199 rhum /= 100.;
00200 prestate = cpl_errorstate_get();
00201 double pres = (muse_pfits_get_pres_start(aPixtable->header)
00202 + muse_pfits_get_pres_end(aPixtable->header)) / 2.;
00203 if (!cpl_errorstate_is_equal(prestate)) {
00204 cpl_msg_warning(__func__, "FITS header does not contain pressure values: %s",
00205 cpl_error_get_message());
00206 pres = 743.;
00207 }
00208
00209
00210
00211
00212 enum {
00213 MUSE_DAR_METHOD_FILIPPENKO = 0,
00214 MUSE_DAR_METHOD_OWENS,
00215 MUSE_DAR_METHOD_EDLEN,
00216 MUSE_DAR_METHOD_CIDDOR,
00217 };
00218 int method = MUSE_DAR_METHOD_FILIPPENKO;
00219 if (getenv("MUSE_DAR_CORRECT_METHOD") &&
00220 !strncmp(getenv("MUSE_DAR_CORRECT_METHOD"), "Owens", 6)) {
00221 method = MUSE_DAR_METHOD_OWENS;
00222 } else if (getenv("MUSE_DAR_CORRECT_METHOD") &&
00223 !strncmp(getenv("MUSE_DAR_CORRECT_METHOD"), "Edlen", 6)) {
00224 method = MUSE_DAR_METHOD_EDLEN;
00225 } else if (getenv("MUSE_DAR_CORRECT_METHOD") &&
00226 !strncmp(getenv("MUSE_DAR_CORRECT_METHOD"), "Ciddor", 7)) {
00227 method = MUSE_DAR_METHOD_CIDDOR;
00228 }
00229
00230
00231 double d1, d2,
00232 fp,
00233 nr0,
00234 p = pres * 100,
00235 t = temp - 273.15,
00236 pv,
00237 xv;
00238 if (method == MUSE_DAR_METHOD_OWENS) {
00239 muse_phys_nrindex_owens_coeffs(temp, rhum, pres, &d1, &d2);
00240 nr0 = muse_phys_nrindex_owens(aLambdaRef / 10000., d1, d2);
00241 cpl_msg_info(__func__, "DAR for T=%.2f K, RH=%.2f %%, pres=%.1f mbar, "
00242 "at airmass=%.4f, using ref. wavelength %.6f um (Owens)",
00243 temp, rhum*100., pres, 1./cos(z), aLambdaRef / 10000.);
00244 } else if (method == MUSE_DAR_METHOD_EDLEN ||
00245 method == MUSE_DAR_METHOD_CIDDOR) {
00246
00247
00248
00249
00250 const double T = temp,
00251 K1 = 1.16705214528E+03,
00252 K2 = -7.24213167032E+05,
00253 K3 = -1.70738469401E+01,
00254 K4 = 1.20208247025E+04,
00255 K5 = -3.23255503223E+06,
00256 K6 = 1.49151086135E+01,
00257 K7 = -4.82326573616E+03,
00258 K8 = 4.05113405421E+05,
00259 K9 = -2.38555575678E-01,
00260 K10 = 6.50175348448E+02,
00261 Omega = T + K9 / (T - K10),
00262 A = Omega*Omega + K1 * Omega + K2,
00263 B = K3 * Omega*Omega + K4 * Omega + K5,
00264 C = K6 * Omega*Omega + K7 * Omega + K8,
00265 X = -B + sqrt(B*B - 4 * A * C);
00266 double psv_w = 1e6 * pow(2 * C / X, 4);
00267
00268
00269 const double A1 = -13.928169,
00270 A2 = 34.7078238,
00271 Theta = T / 273.16,
00272 Y = A1 * (1 - pow(Theta, -1.5)) + A2 * (1 - pow(Theta, -1.25));
00273 double psv_i = 611.657 * exp(Y);
00274
00275
00276
00277 if (t > 0.) {
00278 pv = rhum * psv_w;
00279 } else {
00280 pv = rhum * psv_i;
00281 }
00282
00283 const double alpha = 1.00062,
00284 beta = 3.14e-8,
00285 gamma = 5.60e-7;
00286 double f = alpha + beta * p + gamma * t*t;
00287
00288 if (t > 0.) {
00289 xv = rhum * f * psv_w / p;
00290 } else {
00291 xv = rhum * f * psv_i / p;
00292 }
00293 if (method == MUSE_DAR_METHOD_EDLEN) {
00294 nr0 = muse_phys_nrindex_edlen(aLambdaRef / 10000., t, p, pv);
00295 } else {
00296 nr0 = muse_phys_nrindex_ciddor(aLambdaRef / 10000., T, p, xv, 450.);
00297 }
00298 cpl_msg_info(__func__, "DAR for T=%.2f degC, RH=%.2f %%, p=%.1f Pa, "
00299 "at airmass=%.4f, using ref. wavelength %.6f um (%s)",
00300 t, rhum*100., p, 1./cos(z), aLambdaRef / 10000.,
00301 method == MUSE_DAR_METHOD_EDLEN ? "Edlen" : "Ciddor");
00302 } else {
00303
00304 double ps = muse_phys_nrindex_owens_saturation_pressure(temp);
00305
00306 fp = rhum * ps * MUSE_PHYS_hPa_TO_mmHg;
00307 temp -= 273.15;
00308 pres *= MUSE_PHYS_hPa_TO_mmHg;
00309 nr0 = muse_phys_nrindex_filippenko(aLambdaRef / 10000., temp, fp, pres);
00310 cpl_msg_info(__func__, "DAR for T=%.2f degC, fp=%.3f mmHg, P=%.1f mmHg, "
00311 "at airmass=%.4f, using ref. wavelength %.6f um (Filippenko)",
00312 temp, fp, pres, 1./cos(z), aLambdaRef / 10000.);
00313 }
00314
00315
00316
00317
00318
00319 double xshift = -sin((parang + posang) * CPL_MATH_RAD_DEG),
00320 yshift = cos((parang + posang) * CPL_MATH_RAD_DEG);
00321
00322
00323 muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
00324 cpl_ensure_code(wcstype == MUSE_PIXTABLE_WCS_PIXEL ||
00325 wcstype == MUSE_PIXTABLE_WCS_CELSPH, CPL_ERROR_INVALID_TYPE);
00326 double fscale = 3600.;
00327 if (wcstype == MUSE_PIXTABLE_WCS_CELSPH) {
00328 fscale = 1. / 3600;
00329 double xscale, yscale;
00330 muse_wcs_get_scales(aPixtable->header, &xscale, &yscale);
00331 xshift /= -xscale;
00332 yshift /= yscale;
00333 } else {
00334 if (isWFM) {
00335 xshift /= kMuseSpaxelSizeX_WFM;
00336 yshift /= kMuseSpaxelSizeY_WFM;
00337 } else {
00338 xshift /= kMuseSpaxelSizeX_NFM;
00339 yshift /= kMuseSpaxelSizeY_NFM;
00340 }
00341 }
00342
00343
00344
00345
00346
00347
00348 double dr0 = tan(z) * fscale * CPL_MATH_DEG_RAD;
00349 #if DEBUG_MUSE_DARCORRECT
00350 double wl;
00351 for (wl = 4650.; wl <= 9300; wl += 50) {
00352 double nr;
00353 if (method == MUSE_DAR_METHOD_OWENS) {
00354 nr = muse_phys_nrindex_owens(wl / 10000., d1, d2);
00355 } else if (method == MUSE_DAR_METHOD_EDLEN) {
00356 nr = muse_phys_nrindex_edlen(wl / 10000., t, p, pv);
00357 } else if (method == MUSE_DAR_METHOD_CIDDOR) {
00358 nr = muse_phys_nrindex_ciddor(wl / 10000., temp, p, xv, 450.);
00359 } else {
00360 nr = muse_phys_nrindex_filippenko(wl / 10000., temp, fp, pres);
00361 }
00362 double dr = dr0 * (nr0 - nr);
00363 printf("%.1f Angstrom: %f ==> %f %f %s (%s)\n", wl, dr,
00364 xshift * dr, yshift * dr,
00365 wcstype == MUSE_PIXTABLE_WCS_CELSPH ? "deg" : "pix", __func__);
00366 }
00367 fflush(stdout);
00368 #endif
00369 float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
00370 *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
00371 *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA);
00372 cpl_size i, nmax = muse_pixtable_get_nrow(aPixtable);
00373 #pragma omp parallel for default(none) \
00374 shared(d1, d2, dr0, fp, lbda, method, nmax, nr0, p, pres, pv, t, \
00375 temp, xpos, ypos, xshift, xv, yshift)
00376 for (i = 0; i < nmax; i++) {
00377 double nr, lambda = lbda[i] * 1e-4;
00378 if (method == MUSE_DAR_METHOD_OWENS) {
00379 nr = muse_phys_nrindex_owens(lambda, d1, d2);
00380 } else if (method == MUSE_DAR_METHOD_EDLEN) {
00381 nr = muse_phys_nrindex_edlen(lambda, t, p, pv);
00382 } else if (method == MUSE_DAR_METHOD_CIDDOR) {
00383 nr = muse_phys_nrindex_ciddor(lambda, temp, p, xv, 450.);
00384 } else {
00385 nr = muse_phys_nrindex_filippenko(lambda, temp, fp, pres);
00386 }
00387 double dr = dr0 * (nr0 - nr);
00388
00389 xpos[i] += xshift * dr;
00390 ypos[i] += yshift * dr;
00391 }
00392
00393 cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_NAME,
00394 aLambdaRef);
00395 cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_NAME,
00396 MUSE_HDR_PT_DAR_COMMENT);
00397
00398
00399 cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XLO,
00400 cpl_propertylist_get_float(aPixtable->header,
00401 MUSE_HDR_PT_XLO));
00402 cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XHI,
00403 cpl_propertylist_get_float(aPixtable->header,
00404 MUSE_HDR_PT_XHI));
00405 cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YLO,
00406 cpl_propertylist_get_float(aPixtable->header,
00407 MUSE_HDR_PT_YLO));
00408 cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YHI,
00409 cpl_propertylist_get_float(aPixtable->header,
00410 MUSE_HDR_PT_YHI));
00411 muse_pixtable_compute_limits(aPixtable);
00412
00413 return CPL_ERROR_NONE;
00414 }
00415
00416
00470
00471 cpl_error_code
00472 muse_dar_check(muse_pixtable *aPixtable, double *aMaxShift,
00473 cpl_boolean aCorrect, muse_datacube **aCube)
00474 {
00475 cpl_ensure_code(aMaxShift, CPL_ERROR_NULL_INPUT);
00476 *aMaxShift = -99.;
00477 cpl_ensure_code(aPixtable, CPL_ERROR_NULL_INPUT);
00478 muse_pixtable_wcs wcstype = muse_pixtable_wcs_check(aPixtable);
00479 cpl_ensure_code(wcstype == MUSE_PIXTABLE_WCS_PIXEL, CPL_ERROR_INVALID_TYPE);
00480 if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_CHECK)) {
00481 cpl_msg_info(__func__, "pixel table already checked for DAR accuracy");
00482 }
00483 if (cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_CORR)) {
00484 cpl_msg_info(__func__, "pixel table already corrected for DAR residuals");
00485 }
00486
00487
00488
00489 int cenhalfsize = 5;
00490
00491 int order = 2;
00492
00493
00494 muse_resampling_params *params =
00495 muse_resampling_params_new(MUSE_RESAMPLE_WEIGHTED_DRIZZLE);
00496 params->pfx = 1.;
00497 params->pfy = 1.;
00498 params->pfl = 1.;
00499 params->dlambda = 10.;
00500 params->crsigma = 7.;
00501 muse_pixgrid *grid;
00502 muse_datacube *cube = muse_resampling_cube(aPixtable, params, &grid);
00503 muse_resampling_params_delete(params);
00504 if (!cube) {
00505 muse_pixgrid_delete(grid);
00506 muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
00507 return cpl_error_set_message(__func__, cpl_error_get_code(),
00508 "Failure while creating cube!");
00509 }
00510 if (aCube) {
00511 *aCube = cube;
00512 }
00513 double crpix3 = muse_pfits_get_crpix(cube->header, 3),
00514 crval3 = muse_pfits_get_crval(cube->header, 3),
00515 cdelt3 = muse_pfits_get_cd(cube->header, 3, 3);
00516 cpl_errorstate state = cpl_errorstate_get();
00517 double lambdaref = cpl_propertylist_get_double(aPixtable->header,
00518 MUSE_HDR_PT_DAR_NAME),
00519 lambda1 = (2 - crpix3) * cdelt3 + crval3,
00520 lambda3 = (grid->size_z - 1 - crpix3) * cdelt3 + crval3,
00521 lambda2 = (lambda1 + lambda3) / 2.;
00522 if (!cpl_errorstate_is_equal(state) || lambdaref < 0) {
00523 cpl_errorstate_set(state);
00524
00525 lambdaref = lambda2;
00526 cpl_msg_warning(__func__, "Data is not DAR corrected, using reference "
00527 "wavelength at %.3f Angstrom!", lambdaref);
00528 if (!cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_DAR_CORR)) {
00529 cenhalfsize = 15;
00530 order = 4;
00531 }
00532 }
00533 if (lambdaref > lambda3 || lambdaref < lambda1) {
00534 lambdaref = lambda2;
00535 cpl_msg_warning(__func__, "Reference lambda outside cube, using reference "
00536 "wavelength at %.3f Angstrom!", lambdaref);
00537 }
00538
00539
00540 int iplane, irefplane = lround((lambdaref - crval3) / cdelt3 + crpix3) - 1;
00541
00542
00543
00544 muse_imagelist *list = muse_imagelist_new();
00545 unsigned int ilist = 0;
00546 for (iplane = irefplane - 1; iplane <= irefplane + 1; iplane++) {
00547 muse_image *image = muse_image_new();
00548 image->data = cpl_image_duplicate(cpl_imagelist_get(cube->data, iplane));
00549 image->dq = cpl_image_duplicate(cpl_imagelist_get(cube->dq, iplane));
00550 image->stat = cpl_image_duplicate(cpl_imagelist_get(cube->stat, iplane));
00551 muse_imagelist_set(list, image, ilist++);
00552 }
00553 muse_image *median = muse_combine_median_create(list);
00554 muse_imagelist_delete(list);
00555 if (!median) {
00556 if (!aCube) {
00557 muse_datacube_delete(cube);
00558 }
00559 muse_pixgrid_delete(grid);
00560 muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
00561 return cpl_error_set_message(__func__, cpl_error_get_code(),
00562 "Failure while creating detection image!");
00563 }
00564 #if DEBUG_MUSE_DARCHECK > 1
00565 median->header = cpl_propertylist_new();
00566 cpl_propertylist_append_string(median->header, "BUNIT",
00567 cpl_propertylist_get_string(cube->header,
00568 "BUNIT"));
00569 muse_image_save(median, "IMAGE_darcheck_median.fits");
00570 #endif
00571
00572 muse_image_reject_from_dq(median);
00573
00574 double dsigmas[] = { 50., 30., 10., 8., 6., 5. };
00575 int ndsigmas = sizeof(dsigmas) / sizeof(double);
00576 cpl_vector *vsigmas = cpl_vector_wrap(ndsigmas, dsigmas);
00577 cpl_size isigma = -1;
00578 state = cpl_errorstate_get();
00579 cpl_apertures *apertures = cpl_apertures_extract(median->data, vsigmas, &isigma);
00580 int nx = cpl_image_get_size_x(median->data),
00581 ny = cpl_image_get_size_y(median->data);
00582 muse_image_delete(median);
00583 int napertures = apertures ? cpl_apertures_get_size(apertures) : 0;
00584 if (napertures < 1) {
00585 cpl_vector_unwrap(vsigmas);
00586 cpl_apertures_delete(apertures);
00587 cpl_errorstate_set(state);
00588 if (!aCube) {
00589 muse_datacube_delete(cube);
00590 }
00591 muse_pixgrid_delete(grid);
00592 muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
00593 return cpl_error_set_message(__func__, CPL_ERROR_DATA_NOT_FOUND, "No sources "
00594 "found for DAR check down to %.1f sigma "
00595 "limit on plane %d", dsigmas[ndsigmas - 1],
00596 irefplane + 1);
00597 }
00598 cpl_msg_debug(__func__, "The %.1f sigma threshold was used to find %d source%s,"
00599 " centered on plane %d", cpl_vector_get(vsigmas, isigma),
00600 napertures, napertures == 1 ? "" : "s", iplane+1);
00601 cpl_vector_unwrap(vsigmas);
00602 #if DEBUG_MUSE_DARCHECK
00603 cpl_apertures_dump(apertures, stdout);
00604 fflush(stdout);
00605 #if DEBUG_MUSE_DARCHECK > 1
00606 muse_datacube_save(cube, "DATACUBE_darcheck.fits");
00607 #endif
00608 #endif
00609
00610
00611 int xhalfsize = 2*cenhalfsize > nx ? nx/2 : cenhalfsize,
00612 yhalfsize = 2*cenhalfsize > ny ? ny/2 : cenhalfsize;
00613
00614
00615 int l, nlambda = cpl_imagelist_get_size(cube->data);
00616 cpl_vector *vlambda = cpl_vector_new(nlambda);
00617 cpl_matrix *moffsets = cpl_matrix_new(2, nlambda);
00618 #if !DEBUG_MUSE_DARCHECK
00619 #pragma omp parallel for default(none) \
00620 shared(apertures, aPixtable, cdelt3, crpix3, crval3, cube, grid, \
00621 moffsets, napertures, nlambda, nx, ny, vlambda, xhalfsize, \
00622 yhalfsize)
00623 #endif
00624 for (l = 0; l < nlambda; l++) {
00625 double xoffset = 0., yoffset = 0.,
00626 lambda = (l+1 - crpix3) * cdelt3 + crval3;
00627 int n, noffsets = 0;
00628 for (n = 1; n <= napertures; n++) {
00629
00630 double xc = cpl_apertures_get_centroid_x(apertures, n),
00631 yc = cpl_apertures_get_centroid_y(apertures, n);
00632 int x1 = lround(xc) - xhalfsize,
00633 x2 = lround(xc) + xhalfsize,
00634 y1 = lround(yc) - yhalfsize,
00635 y2 = lround(yc) + yhalfsize;
00636
00637 if (x1 < 1) x1 = 1;
00638 if (y1 < 1) y1 = 1;
00639 if (x2 > nx) x2 = nx;
00640 if (y2 > ny) y2 = ny;
00641
00642 #define MUSE_DARCHECK_MAX_INTERNAL_SHIFT 5
00643
00644 #if MUSE_DARCHECK_GRID_FITS
00645
00646 #if DEBUG_MUSE_DARCHECK
00647 const char *method = "grid/moffat";
00648 #endif
00649
00650 float *cdata = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_DATA),
00651 *cstat = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_STAT);
00652 int *cdq = cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_DQ);
00653 #define MUSE_DARCHECK_MAX_NPIX 7000
00654 cpl_matrix *pos = cpl_matrix_new(MUSE_DARCHECK_MAX_NPIX, 2);
00655 cpl_vector *val = cpl_vector_new(MUSE_DARCHECK_MAX_NPIX),
00656 *err = cpl_vector_new(MUSE_DARCHECK_MAX_NPIX);
00657 int i, npix = 0, nsize = MUSE_DARCHECK_MAX_NPIX;
00658 for (i = x1 - 1; i < x2; i++) {
00659 int j;
00660 for (j = y1 - 1; j < y2; j++) {
00661 cpl_size idx = muse_pixgrid_get_index(grid, i, j, l, CPL_TRUE),
00662 ipx, npx = muse_pixgrid_get_count(grid, idx);
00663 const cpl_size *rows = muse_pixgrid_get_rows(grid, idx);
00664 for (ipx = 0; ipx < npx; ipx++) {
00665 if (cdq[rows[ipx]] != EURO3D_GOODPIXEL) {
00666 continue;
00667 }
00668 if (npix >= nsize) {
00669 nsize += MUSE_DARCHECK_MAX_NPIX;
00670 cpl_matrix_resize(pos, 0, nsize - cpl_matrix_get_nrow(pos), 0, 0);
00671 cpl_vector_set_size(val, nsize);
00672 cpl_vector_set_size(err, nsize);
00673 }
00674 cpl_matrix_set(pos, npix, 0, i+1);
00675 cpl_matrix_set(pos, npix, 1, j+1);
00676 cpl_vector_set(val, npix, cdata[rows[ipx]]);
00677 cpl_vector_set(err, npix, sqrt(cstat[rows[ipx]]));
00678 npix++;
00679 }
00680 }
00681 }
00682 if (npix < 8) {
00683 cpl_matrix_delete(pos);
00684 cpl_vector_delete(val);
00685 cpl_vector_delete(err);
00686 continue;
00687 }
00688 cpl_matrix_set_size(pos, npix, 2);
00689 cpl_vector_set_size(val, npix);
00690 cpl_vector_set_size(err, npix);
00691
00692
00693 double xc1, yc1;
00694 muse_utils_get_centroid(pos, val, NULL, &xc1, &yc1,
00695 MUSE_UTILS_CENTROID_MEDIAN);
00696
00697
00698 cpl_array *pmoffat = cpl_array_new(8, CPL_TYPE_DOUBLE);
00699 cpl_array_set(pmoffat, 2, xc1);
00700 cpl_array_set(pmoffat, 3, yc1);
00701 cpl_error_code rc = muse_utils_fit_moffat_2d(pos, val, err, pmoffat, NULL,
00702 NULL, NULL, NULL);
00703 double xc2 = cpl_array_get(pmoffat, 2, NULL),
00704 yc2 = cpl_array_get(pmoffat, 3, NULL);
00705 cpl_array_delete(pmoffat);
00706 cpl_matrix_delete(pos);
00707 cpl_vector_delete(val);
00708 cpl_vector_delete(err);
00709
00710
00711
00712 if (rc == CPL_ERROR_NONE &&
00713 fabs(xc1 - xc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT &&
00714 fabs(yc1 - yc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT) {
00715 xoffset += xc - xc2;
00716 yoffset += yc - yc2;
00717 noffsets++;
00718 }
00719 #else
00720 #if DEBUG_MUSE_DARCHECK
00721 const char *method = "image/gauss";
00722 #endif
00723
00724 cpl_image *plane = cpl_imagelist_get(cube->data, l),
00725 *plerr = cpl_image_power_create(cpl_imagelist_get(cube->stat, l), 0.5);
00726
00727 muse_quality_image_reject_using_dq(plane, cpl_imagelist_get(cube->dq, l), plerr);
00728
00729 cpl_image_reject_from_mask(plane, cpl_image_get_bpm(plerr));
00730
00731 cpl_image *xtr = cpl_image_extract(plane, x1, y1, x2, y2);
00732 int npix = (x2-x1+1) * (y2-y1+1)
00733 - cpl_image_count_rejected(xtr);
00734 cpl_image_delete(xtr);
00735 if (npix < 7) {
00736 cpl_image_delete(plerr);
00737 continue;
00738 }
00739
00740
00741 double xc1, yc1;
00742 muse_utils_image_get_centroid_window(plane, x1, y1, x2, y2, &xc1, &yc1,
00743 MUSE_UTILS_CENTROID_MEAN);
00744
00745
00746 cpl_array *pgauss = cpl_array_new(7, CPL_TYPE_DOUBLE);
00747 cpl_array_set(pgauss, 3, xc1);
00748 cpl_array_set(pgauss, 4, yc1);
00749 cpl_error_code rc = cpl_fit_image_gaussian(plane, plerr,
00750 lround(xc), lround(yc),
00751 2*xhalfsize, 2*yhalfsize,
00752 pgauss, NULL, NULL,
00753 NULL, NULL, NULL,
00754 NULL, NULL, NULL, NULL);
00755 double xc2 = cpl_array_get(pgauss, 3, NULL),
00756 yc2 = cpl_array_get(pgauss, 4, NULL);
00757 cpl_array_delete(pgauss);
00758 cpl_image_delete(plerr);
00759
00760
00761
00762 if (rc == CPL_ERROR_NONE &&
00763 fabs(xc1 - xc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT &&
00764 fabs(yc1 - yc2) < MUSE_DARCHECK_MAX_INTERNAL_SHIFT) {
00765 xoffset += xc - xc2;
00766 yoffset += yc - yc2;
00767 noffsets++;
00768 }
00769 #endif
00770 #if DEBUG_MUSE_DARCHECK
00771 printf("%.1f Angstrom plane %d aper %d (%f,%f): %f %f (%s) %f %f (%d, %s)\n",
00772 lambda, l+1, n, xc, yc, xc - xc2, yc - yc2, __func__,
00773 xc2 - xc1, yc2 - yc1, npix, method);
00774 fflush(stdout);
00775 #endif
00776 }
00777
00778 cpl_vector_set(vlambda, l, lambda);
00779 cpl_matrix_set(moffsets, 0, l, xoffset / noffsets);
00780 cpl_matrix_set(moffsets, 1, l, yoffset / noffsets);
00781 #if DEBUG_MUSE_DARCHECK
00782 printf("%.1f Angstrom plane %d all-apers: %f %f (%s)\n",
00783 lambda, l+1, xoffset / noffsets, yoffset / noffsets, __func__);
00784 fflush(stdout);
00785 #endif
00786 }
00787 cpl_apertures_delete(apertures);
00788
00789
00790 cpl_vector *vlam1 = cpl_vector_duplicate(vlambda),
00791 *vlam2 = cpl_vector_duplicate(vlambda);
00792 cpl_matrix *mlam1 = cpl_matrix_wrap(1, nlambda, cpl_vector_unwrap(vlam1)),
00793 *mlam2 = cpl_matrix_wrap(1, nlambda, cpl_vector_unwrap(vlam2));
00794 cpl_matrix *mxoff = cpl_matrix_extract_row(moffsets, 0),
00795 *myoff = cpl_matrix_extract_row(moffsets, 1);
00796 cpl_vector *vxoff = cpl_vector_wrap(nlambda, cpl_matrix_unwrap(mxoff)),
00797 *vyoff = cpl_vector_wrap(nlambda, cpl_matrix_unwrap(myoff));
00798
00799
00800 double msex, msey, chisqx, chisqy;
00801 cpl_polynomial *px = muse_utils_iterate_fit_polynomial(mlam1, vxoff,
00802 NULL, NULL, order, 3.,
00803 &msex, &chisqx),
00804 *py = muse_utils_iterate_fit_polynomial(mlam2, vyoff,
00805 NULL, NULL, order, 3.,
00806 &msey, &chisqy);
00807 cpl_vector_delete(vxoff);
00808 cpl_vector_delete(vyoff);
00809 cpl_matrix_delete(mlam1);
00810 cpl_matrix_delete(mlam2);
00811 #if DEBUG_MUSE_DARCHECK
00812 printf("polynomial fit in x (order %d, rms=%f, chisq=%g):\n", order,
00813 sqrt(msex), chisqx);
00814 cpl_polynomial_dump(px, stdout);
00815 printf("polynomial fit in y (order %d, rms=%f, chisq=%g):\n", order,
00816 sqrt(msey), chisqy);
00817 cpl_polynomial_dump(py, stdout);
00818 fflush(stdout);
00819 #endif
00820
00821
00822
00823 double xref = cpl_polynomial_eval_1d(px, lambdaref, NULL),
00824 yref = cpl_polynomial_eval_1d(py, lambdaref, NULL);
00825 cpl_size pows = 0;
00826 cpl_polynomial_set_coeff(px, &pows, cpl_polynomial_get_coeff(px, &pows) - xref);
00827 cpl_polynomial_set_coeff(py, &pows, cpl_polynomial_get_coeff(py, &pows) - yref);
00828 #if DEBUG_MUSE_DARCHECK
00829 printf("polynomial in x (xref=%f at %f --> %f):\n", xref, lambdaref,
00830 cpl_polynomial_eval_1d(px, lambdaref, NULL));
00831 cpl_polynomial_dump(px, stdout);
00832 printf("polynomial in y (yref=%f at %f --> %f):\n", yref, lambdaref,
00833 cpl_polynomial_eval_1d(py, lambdaref, NULL));
00834 cpl_polynomial_dump(py, stdout);
00835 fflush(stdout);
00836 #endif
00837
00838
00839 double xmin = DBL_MAX, xmax = -DBL_MAX,
00840 ymin = DBL_MAX, ymax = -DBL_MAX;
00841 for (l = 0; l < nlambda; l++) {
00842 double lambda = cpl_vector_get(vlambda, l),
00843 xshift = cpl_polynomial_eval_1d(px, lambda, NULL),
00844 yshift = cpl_polynomial_eval_1d(py, lambda, NULL);
00845 if (xshift < xmin) xmin = xshift;
00846 if (xshift > xmax) xmax = xshift;
00847 if (yshift < ymin) ymin = yshift;
00848 if (yshift > ymax) ymax = yshift;
00849 #if DEBUG_MUSE_DARCHECK
00850 printf("%.1f Angstrom plane %d all-fitted: %f %f (%s)\n",
00851 lambda, l+1, xshift, yshift, __func__);
00852 fflush(stdout);
00853 #endif
00854 }
00855 cpl_vector_delete(vlambda);
00856 cpl_matrix_delete(moffsets);
00857
00858
00859 double maxdiffx = xmax - xmin,
00860 maxdiffy = ymax - ymin;
00861
00862 if (muse_pfits_get_mode(aPixtable->header) <= MUSE_MODE_WFM_AO_N) {
00863 maxdiffx *= kMuseSpaxelSizeX_WFM;
00864 maxdiffy *= kMuseSpaxelSizeY_WFM;
00865 } else {
00866 maxdiffx *= kMuseSpaxelSizeX_NFM;
00867 maxdiffy *= kMuseSpaxelSizeY_NFM;
00868 }
00869 *aMaxShift = fmax(maxdiffx, maxdiffy);
00870 cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_CHECK,
00871 *aMaxShift);
00872 cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_CHECK,
00873 MUSE_HDR_PT_DAR_CHECK_C);
00874
00875 muse_pixgrid_delete(grid);
00876 if (!aCube) {
00877 muse_datacube_delete(cube);
00878 }
00879
00880
00881 if (aCorrect) {
00882 cpl_msg_debug(__func__, "Correcting pixel table for the shift (max = %f "
00883 "arcsec)", *aMaxShift);
00884 float *xpos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_XPOS),
00885 *ypos = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_YPOS),
00886 *lbda = cpl_table_get_data_float(aPixtable->table, MUSE_PIXTABLE_LAMBDA);
00887 cpl_size i, nrow = muse_pixtable_get_nrow(aPixtable);
00888 #pragma omp parallel for default(none) \
00889 shared(lbda, nrow, px, py, xpos, ypos)
00890 for (i = 0; i < nrow; i++) {
00891
00892 xpos[i] += cpl_polynomial_eval_1d(px, lbda[i], NULL);
00893 ypos[i] += cpl_polynomial_eval_1d(py, lbda[i], NULL);
00894 }
00895
00896
00897 cpl_propertylist_update_double(aPixtable->header, MUSE_HDR_PT_DAR_CORR,
00898 lambdaref);
00899 cpl_propertylist_set_comment(aPixtable->header, MUSE_HDR_PT_DAR_CORR,
00900 MUSE_HDR_PT_DAR_CORR_C);
00901
00902
00903
00904 if (!cpl_propertylist_has(aPixtable->header, MUSE_HDR_PT_PREDAR_XLO)) {
00905 cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XLO,
00906 cpl_propertylist_get_float(aPixtable->header,
00907 MUSE_HDR_PT_XLO));
00908 cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_XHI,
00909 cpl_propertylist_get_float(aPixtable->header,
00910 MUSE_HDR_PT_XHI));
00911 cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YLO,
00912 cpl_propertylist_get_float(aPixtable->header,
00913 MUSE_HDR_PT_YLO));
00914 cpl_propertylist_update_float(aPixtable->header, MUSE_HDR_PT_PREDAR_YHI,
00915 cpl_propertylist_get_float(aPixtable->header,
00916 MUSE_HDR_PT_YHI));
00917 }
00918 muse_pixtable_compute_limits(aPixtable);
00919 }
00920 cpl_polynomial_delete(px);
00921 cpl_polynomial_delete(py);
00922
00923
00924 muse_pixtable_reset_dq(aPixtable, EURO3D_COSMICRAY);
00925
00926 return CPL_ERROR_NONE;
00927 }
00928