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 #define STORE_SLIT_WIDTH
00028 #define STORE_BIN_WIDTH
00029
00030
00031
00032
00033 #include <math.h>
00034 #include <string.h>
00035 #include <fenv.h>
00036
00037 #include <cpl.h>
00038
00039 #include "muse_instrument.h"
00040 #include "muse_lsf.h"
00041 #include "muse_optimize.h"
00042 #include "muse_pfits.h"
00043 #include "muse_pixtable.h"
00044 #include "muse_quality.h"
00045 #include "muse_resampling.h"
00046 #include "muse_tracing.h"
00047 #include "muse_utils.h"
00048 #include "muse_wavecalib.h"
00049
00050
00054
00055
00074 static cpl_size
00075 muse_lsf_check_arc_line(cpl_table *aPixtable, double aLambda) {
00076 cpl_size nrows = cpl_table_get_nrow(aPixtable);
00077
00078 if (nrows <= 100) {
00079 return 0;
00080 }
00081 cpl_errorstate prestate = cpl_errorstate_get();
00082
00083
00084 cpl_vector *pos = cpl_vector_new(nrows);
00085 cpl_vector *val = cpl_vector_new(nrows);
00086 cpl_vector *err = cpl_vector_new(nrows);
00087 double *d_pos = cpl_vector_get_data(pos);
00088 double *d_val = cpl_vector_get_data(val);
00089 double *d_err = cpl_vector_get_data(err);
00090 float *d_lambda = cpl_table_get_data_float(aPixtable, MUSE_PIXTABLE_LAMBDA);
00091 float *d_data = cpl_table_get_data_float(aPixtable, MUSE_PIXTABLE_DATA);
00092 float *d_stat = cpl_table_get_data_float(aPixtable, MUSE_PIXTABLE_STAT);
00093 cpl_size i, s = 1;
00094 double bglevel = 0;
00095 for (i = 0; i < nrows; i++) {
00096 if (fabs(d_lambda[i] - aLambda) > 5.5) {
00097 bglevel += d_data[i];
00098 s += 1;
00099 }
00100 d_pos[i] = d_lambda[i];
00101 d_val[i] = d_data[i];
00102 d_err[i] = sqrt(d_stat[i]);
00103 }
00104 bglevel /= s;
00105 double meansigma = sqrt(cpl_table_get_column_mean(aPixtable, MUSE_PIXTABLE_STAT));
00106 double sigma, area, mse;
00107 double xc = aLambda;
00108 cpl_fit_mode fit_pars = CPL_FIT_STDEV | CPL_FIT_AREA | CPL_FIT_CENTROID;
00109 cpl_vector_fit_gaussian(pos, NULL, val, NULL, fit_pars, &xc, &sigma,
00110 &area, &bglevel, &mse, NULL, NULL);
00111
00112
00113 for (i = 0; i < nrows; i++) {
00114 double e = (d_pos[i] - xc)/sigma;
00115 d_val[i] -= area / (CPL_MATH_SQRT2PI * sigma) * exp(-e*e/2) - bglevel;
00116 }
00117
00118
00119 double maxerror = 0.5 * sqrt(area) + 0.15 * area;
00120 cpl_size inv = 0;
00121 cpl_table_unselect_all(aPixtable);
00122 for (i = 0; i < nrows; i++) {
00123 if ((d_val[i] < -maxerror) || (d_val[i] > maxerror)) {
00124 inv++;
00125 cpl_table_select_row(aPixtable, i);
00126 }
00127 }
00128 cpl_table_erase_selected(aPixtable);
00129 nrows = cpl_table_get_nrow(aPixtable);
00130 cpl_vector_delete(pos);
00131 cpl_vector_delete(val);
00132 cpl_vector_delete(err);
00133
00134 cpl_table_fill_column_window_float(aPixtable, "line_lambda",
00135 0, nrows, xc);
00136 cpl_table_fill_column_window_float(aPixtable, "line_flux",
00137 0, nrows, area);
00138 cpl_table_fill_column_window_float(aPixtable, "line_background",
00139 0, nrows, bglevel);
00140
00141
00142 if (!cpl_errorstate_is_equal(prestate)) {
00143 cpl_errorstate_set(prestate);
00144 return 0;
00145 }
00146
00147 if (area/meansigma < 50) {
00148 return 0;
00149 }
00150
00151
00152 if ((sigma > 3) || (sigma < 0.7)) {
00153 uint32_t origin = (uint32_t)cpl_table_get_int(aPixtable,
00154 MUSE_PIXTABLE_ORIGIN,
00155 0, NULL);
00156 int i_ifu = muse_pixtable_origin_get_ifu(origin);
00157 int i_slice = muse_pixtable_origin_get_slice(origin);
00158 cpl_msg_debug(__func__, "Slice %2i.%02i: "
00159 "Ignoring line %.1f with implausible width %f (flux=%.0f)",
00160 i_ifu, i_slice, aLambda, sigma, area);
00161 return 0;
00162 }
00163 return nrows;
00164 }
00165
00185 static cpl_table *
00186 select_arc_line(muse_pixtable *aPixtable, double aLambda, double aWindow) {
00187 cpl_size low = muse_cpltable_find_sorted(aPixtable->table, MUSE_PIXTABLE_LAMBDA,
00188 aLambda - aWindow);
00189 cpl_size high = muse_cpltable_find_sorted(aPixtable->table, MUSE_PIXTABLE_LAMBDA,
00190 aLambda + aWindow);
00191 cpl_table *sel = cpl_table_extract(aPixtable->table, low + 1, high - low);
00192 cpl_table_select_all(sel);
00193 cpl_table_and_selected_int(sel, "dq", CPL_NOT_EQUAL_TO, 0);
00194 cpl_table_erase_selected(sel);
00195
00196 cpl_table_new_column(sel, "line_lambda", CPL_TYPE_FLOAT);
00197 cpl_table_new_column(sel, "line_flux", CPL_TYPE_FLOAT);
00198 cpl_table_new_column(sel, "line_background", CPL_TYPE_FLOAT);
00199
00200 cpl_size nrows = muse_lsf_check_arc_line(sel, aLambda);
00201
00202 if (nrows == 0) {
00203 cpl_table_delete(sel);
00204 return NULL;
00205 }
00206 return sel;
00207 }
00208
00209
00232 muse_pixtable *
00233 muse_lsf_create_arcpixtable(muse_imagelist *aImages,
00234 cpl_table *aTrace, cpl_table *aWave,
00235 cpl_table *aArcLines, int aQuality, double aWindow)
00236 {
00237 cpl_ensure(aImages, CPL_ERROR_NULL_INPUT, NULL);
00238 cpl_ensure(aTrace, CPL_ERROR_NULL_INPUT, NULL);
00239 cpl_ensure(aWave, CPL_ERROR_NULL_INPUT, NULL);
00240 cpl_ensure(aArcLines, CPL_ERROR_NULL_INPUT, NULL);
00241
00242
00243 muse_pixtable *pt = cpl_calloc(1, sizeof(muse_pixtable));
00244 pt->table = muse_cpltable_new(muse_pixtable_def, 0);
00245 cpl_table_new_column(pt->table, "line_lambda", CPL_TYPE_FLOAT);
00246 cpl_table_new_column(pt->table, "line_flux", CPL_TYPE_FLOAT);
00247 cpl_table_new_column(pt->table, "line_background", CPL_TYPE_FLOAT);
00248
00249 cpl_size i_image, n_images = muse_imagelist_get_size(aImages);
00250 #if 0
00251 #pragma omp parallel for default(none) \
00252 shared(i_image, n_images, aImages, aTrace, aWave, aArcLines, \
00253 aQuality, aWindow, pt)
00254 #endif
00255 for (i_image = 0; i_image < n_images; i_image++) {
00256 muse_image *image = muse_imagelist_get(aImages, i_image);
00257 muse_pixtable *pixtable = muse_pixtable_create(image, aTrace, aWave, NULL);
00258 if (pixtable == NULL) {
00259 continue;
00260 }
00261
00262 if (!pt->header) {
00263 pt->header = cpl_propertylist_duplicate(pixtable->header);
00264 }
00265 char *l = muse_utils_header_get_lamp_names(pixtable->header, '|');
00266 if (strlen(l) == 0) {
00267 cpl_msg_warning(__func__,
00268 "Ignoring frame without arc lamp switched on");
00269 cpl_free(l);
00270 cpl_free(pixtable);
00271 continue;
00272 }
00273 char *lamp = cpl_malloc(strlen(l)+3);
00274 sprintf(lamp, "|%s|", l);
00275 cpl_free(l);
00276 cpl_size i_line;
00277 cpl_size n_entries = 0;
00278 muse_pixtable **slice_pixtable = muse_pixtable_extracted_get_slices(pixtable);
00279 int n_slices = muse_pixtable_extracted_get_size(slice_pixtable);
00280 for (i_line = 0; i_line < cpl_table_get_nrow(aArcLines); i_line++) {
00281 if (cpl_table_get_int(aArcLines, "quality", i_line, NULL) < aQuality) {
00282 continue;
00283 }
00284 const char *n = muse_wave_lines_get_lampname(aArcLines, i_line);
00285 char *name = cpl_malloc(strlen(n) + 3);
00286 sprintf(name, "|%s|", n);
00287 if (strstr(name, lamp)) {
00288 double lambda = cpl_table_get_float(aArcLines, "lambda", i_line, NULL);
00289 int i_slice;
00290 #if 0
00291 #pragma omp parallel for default(none) \
00292 shared(i_slice, n_slices, n_entries, lambda, slice_pixtable, \
00293 aWindow, pt)
00294 #endif
00295 for (i_slice = 0; i_slice < n_slices; i_slice++) {
00296 cpl_table *sel = select_arc_line(slice_pixtable[i_slice],
00297 lambda, aWindow);
00298 if (sel != NULL) {
00299 #if 0
00300 #pragma omp critical(construct_pixtable)
00301 #endif
00302 cpl_table_insert(pt->table, sel, cpl_table_get_nrow(pt->table));
00303 #if 0
00304 #pragma omp atomic
00305 #endif
00306 n_entries += cpl_table_get_nrow(sel);
00307 cpl_table_delete(sel);
00308 }
00309 }
00310 }
00311 cpl_free(name);
00312 }
00313 muse_pixtable_extracted_delete(slice_pixtable);
00314 muse_pixtable_delete(pixtable);
00315 cpl_msg_info(__func__, "Using %"CPL_SIZE_FORMAT" entries with lamp %s",
00316 n_entries, lamp);
00317 cpl_free(lamp);
00318 }
00319 return pt;
00320 }
00321
00345 static cpl_polynomial *
00346 muse_lsf_fit_polynomial(cpl_table *aPixtable, int aOrderLsf, int aOrderLambda)
00347 {
00348 cpl_errorstate prestate = cpl_errorstate_get();
00349 cpl_size nrows = cpl_table_get_nrow(aPixtable);
00350 if (nrows <= (aOrderLsf + 1) * (aOrderLambda)) {
00351 return NULL;
00352 }
00353 double *p_dlambda = cpl_table_get_data_double(aPixtable, "d_lambda");
00354 cpl_ensure(p_dlambda != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
00355 float *p_line = cpl_table_get_data_float(aPixtable, "line_lambda");
00356 cpl_ensure(p_line != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
00357 float *p_data = cpl_table_get_data_float(aPixtable, MUSE_PIXTABLE_DATA);
00358 cpl_ensure(p_data != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
00359 float *p_flux = cpl_table_get_data_float(aPixtable, "line_flux");
00360 cpl_ensure(p_flux != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
00361 float *p_background = cpl_table_get_data_float(aPixtable, "line_background");
00362 cpl_ensure(p_background != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
00363
00364 cpl_matrix *coord = cpl_matrix_new(2, nrows);
00365 cpl_vector *data = cpl_vector_new(nrows);
00366 double *coord_ptr = cpl_matrix_get_data(coord);
00367 double *data_ptr = cpl_vector_get_data(data);
00368 cpl_size i_row;
00369
00370 for (i_row = 0; i_row < nrows; i_row++) {
00371 coord_ptr[i_row] = p_dlambda[i_row];
00372 coord_ptr[i_row + nrows] = p_line[i_row];
00373 data_ptr[i_row] = (p_data[i_row] - p_background[i_row]) / p_flux[i_row];
00374 }
00375
00376 cpl_polynomial *p = cpl_polynomial_new(2);
00377 const cpl_size maxdeg2d[] = { aOrderLsf, aOrderLambda };
00378 cpl_polynomial_fit(p, coord, NULL, data, NULL, CPL_TRUE, NULL, maxdeg2d);
00379 cpl_matrix_delete(coord);
00380 cpl_vector_delete(data);
00381
00382 if (!cpl_errorstate_is_equal(prestate)) {
00383 cpl_errorstate_set(prestate);
00384 return NULL;
00385 }
00386
00387 return p;
00388 }
00389
00417 cpl_error_code
00418 muse_lsf_fit_slice(const muse_pixtable *aPixtable, cpl_image *aLsfImage,
00419 muse_wcs *aWCS, double aLsfRegressionWindow)
00420 {
00421 cpl_ensure_code(aPixtable && aPixtable->table, CPL_ERROR_NULL_INPUT);
00422 cpl_ensure_code(aLsfImage, CPL_ERROR_NULL_INPUT);
00423 cpl_ensure_code(aWCS, CPL_ERROR_NULL_INPUT);
00424 cpl_ensure_code(aLsfRegressionWindow > 0, CPL_ERROR_ILLEGAL_INPUT);
00425
00426 cpl_size nrows = muse_pixtable_get_nrow(aPixtable);
00427 if (nrows == 0) {
00428 return CPL_ERROR_NONE;
00429 }
00430
00431 cpl_table *pixtable = cpl_table_duplicate(aPixtable->table);
00432 uint32_t origin = (uint32_t)cpl_table_get_int(pixtable,
00433 MUSE_PIXTABLE_ORIGIN, 0, NULL);
00434 int i_ifu = muse_pixtable_origin_get_ifu(origin);
00435 int i_slice = muse_pixtable_origin_get_slice(origin);
00436 cpl_msg_info(__func__, "processing slice %2i.%02i"
00437 " with %"CPL_SIZE_FORMAT" entries",
00438 i_ifu, i_slice, nrows);
00439
00440
00441
00442
00443
00444 cpl_table_cast_column(pixtable, MUSE_PIXTABLE_LAMBDA, "d_lambda", CPL_TYPE_DOUBLE);
00445 cpl_table_subtract_columns(pixtable, "d_lambda", "line_lambda");
00446 cpl_propertylist *order = cpl_propertylist_new();
00447 cpl_propertylist_append_bool(order, "d_lambda", CPL_FALSE);
00448 cpl_table_sort(pixtable, order);
00449 cpl_propertylist_delete(order);
00450
00451
00452
00453 cpl_size i_lsf;
00454 cpl_size n_lsf = cpl_image_get_size_x(aLsfImage);
00455 cpl_size n_lambda = cpl_image_get_size_y(aLsfImage);
00456 for (i_lsf = 1; i_lsf <= n_lsf; i_lsf++) {
00457 double x = aWCS->crval1 + (i_lsf - aWCS->crpix1) * aWCS->cd11;
00458 double x_ref = CPL_MAX(x, cpl_table_get_column_min(pixtable, "d_lambda")
00459 + aLsfRegressionWindow);
00460 x_ref = CPL_MIN(x, cpl_table_get_column_max(pixtable, "d_lambda")
00461 - aLsfRegressionWindow);
00462 cpl_size low = muse_cpltable_find_sorted(pixtable, "d_lambda",
00463 x_ref - aLsfRegressionWindow);
00464 cpl_size high = muse_cpltable_find_sorted(pixtable, "d_lambda",
00465 x_ref + aLsfRegressionWindow);
00466 cpl_table *selected = cpl_table_extract(pixtable, low+1, high-low);
00467 cpl_polynomial *p = muse_lsf_fit_polynomial(selected, 2, 3);
00468 if (p != NULL) {
00469 cpl_vector *c = cpl_vector_new(2);
00470 cpl_vector_set(c, 0, x);
00471 cpl_size i_lambda;
00472
00473 for (i_lambda = 1; i_lambda <= n_lambda; i_lambda++) {
00474 double y = aWCS->crval2 + (i_lambda - aWCS->crpix2) * aWCS->cd22;
00475 cpl_vector_set(c, 1, y);
00476 double flux = cpl_polynomial_eval(p, c);
00477 cpl_image_set(aLsfImage, i_lsf, i_lambda, flux);
00478 }
00479 cpl_vector_delete(c);
00480 } else {
00481 cpl_msg_warning(__func__,
00482 "Failed polynomial fit %2i.%02i %.2f;"
00483 " %"CPL_SIZE_FORMAT" entries",
00484 i_ifu, i_slice, x, cpl_table_get_nrow(selected));
00485 cpl_size i_lambda;
00486 for (i_lambda = 1; i_lambda <= n_lambda; i_lambda++) {
00487 cpl_image_reject(aLsfImage, i_lsf, i_lambda);
00488 }
00489 }
00490 cpl_polynomial_delete(p);
00491 cpl_table_delete(selected);
00492 }
00493 cpl_table_delete(pixtable);
00494
00495
00496 cpl_size i_lambda;
00497 for (i_lambda = 1; i_lambda <= n_lambda; i_lambda++) {
00498 double line_flux = cpl_image_get_flux_window(aLsfImage, 1, i_lambda,
00499 n_lsf, i_lambda) * aWCS->cd11;
00500 for (i_lsf = 1; i_lsf <= n_lsf; i_lsf++) {
00501 int invalid;
00502 double flux = cpl_image_get(aLsfImage, i_lsf, i_lambda, &invalid);
00503 if (!invalid) {
00504 cpl_image_set(aLsfImage, i_lsf, i_lambda, flux / line_flux);
00505 }
00506 }
00507 }
00508 return CPL_ERROR_NONE;
00509 }
00510
00519 muse_lsf_cube *
00520 muse_lsf_cube_new(double aLsfHalfRange, cpl_size aNLsf, cpl_size aNLambda,
00521 const cpl_propertylist *aHeader)
00522 {
00523 muse_lsf_cube *lsf = cpl_calloc(1, sizeof(muse_lsf_cube));
00524 lsf->header = cpl_propertylist_new();
00525 if (aHeader) {
00526
00527 const char *regex = MUSE_HDR_OVSC_REGEXP"|"MUSE_WCS_KEYS"|"MUSE_HDR_PT_REGEXP;
00528 cpl_propertylist_copy_property_regexp(lsf->header, aHeader, regex, 1);
00529 }
00530 lsf->img = cpl_imagelist_new();
00531 cpl_size i;
00532 for (i = 0; i < kMuseSlicesPerCCD; i++) {
00533 cpl_imagelist_set(lsf->img,
00534 cpl_image_new(aNLsf, aNLambda, CPL_TYPE_FLOAT), i);
00535 }
00536 double lsf_step = 2*aLsfHalfRange / (aNLsf - 1),
00537 lambda_step = (kMuseNominalLambdaMax - kMuseNominalLambdaMin)
00538 / (aNLambda - 1);
00539 lsf->wcs = cpl_calloc(1, sizeof(muse_wcs));
00540 lsf->wcs->cd11 = lsf_step;
00541 lsf->wcs->cd12 = 0.;
00542 lsf->wcs->cd21 = 0.;
00543 lsf->wcs->cd22 = lambda_step;
00544 lsf->wcs->crval1 = -aLsfHalfRange;
00545 lsf->wcs->crval2 = kMuseNominalLambdaMin;
00546 lsf->wcs->crpix1 = 1.;
00547 lsf->wcs->crpix2 = 1.;
00548 return lsf;
00549 }
00550
00555 void
00556 muse_lsf_cube_delete(muse_lsf_cube *aLsfCube)
00557 {
00558 if (!aLsfCube) {
00559 return;
00560 }
00561 cpl_propertylist_delete(aLsfCube->header);
00562 cpl_imagelist_delete(aLsfCube->img);
00563 cpl_free(aLsfCube->wcs);
00564 cpl_free(aLsfCube);
00565 }
00566
00572 cpl_error_code
00573 muse_lsf_cube_save(muse_lsf_cube *aLsfCube, const char *aFileName)
00574 {
00575 cpl_ensure_code(aLsfCube, CPL_ERROR_NULL_INPUT);
00576 cpl_error_code rc = cpl_propertylist_save(aLsfCube->header, aFileName,
00577 CPL_IO_CREATE);
00578 if (rc != CPL_ERROR_NONE) {
00579 return rc;
00580 }
00581
00582 cpl_propertylist *header = cpl_propertylist_new();
00583 cpl_propertylist_append_string(header, "EXTNAME", "LSF_PROFILE");
00584 cpl_propertylist_append_int(header, "WCSAXES", 2);
00585 cpl_propertylist_append_double(header, "CD1_1", aLsfCube->wcs->cd11);
00586 cpl_propertylist_append_double(header, "CD1_2", aLsfCube->wcs->cd12);
00587 cpl_propertylist_append_double(header, "CD2_1", aLsfCube->wcs->cd21);
00588 cpl_propertylist_append_double(header, "CD2_2", aLsfCube->wcs->cd22);
00589 cpl_propertylist_append_double(header, "CRPIX1", aLsfCube->wcs->crpix1);
00590 cpl_propertylist_append_double(header, "CRPIX2", aLsfCube->wcs->crpix2);
00591 cpl_propertylist_append_double(header, "CRVAL1", aLsfCube->wcs->crval1);
00592 cpl_propertylist_append_double(header, "CRVAL2", aLsfCube->wcs->crval2);
00593 cpl_propertylist_append_string(header, "CTYPE1", "PARAM");
00594 cpl_propertylist_append_string(header, "CTYPE2", "AWAV");
00595 cpl_propertylist_append_string(header, "CUNIT1", "Angstrom");
00596 cpl_propertylist_append_string(header, "CUNIT2", "Angstrom");
00597
00598 rc = cpl_imagelist_save(aLsfCube->img, aFileName,
00599 CPL_TYPE_FLOAT, header, CPL_IO_EXTEND);
00600 cpl_propertylist_delete(header);
00601 return rc;
00602 }
00603
00604
00616
00617 muse_lsf_cube **
00618 muse_lsf_cube_load_all(muse_processing *aProcessing)
00619 {
00620 cpl_ensure(aProcessing, CPL_ERROR_NULL_INPUT, NULL);
00621 muse_lsf_cube **lsf = cpl_calloc(kMuseNumIFUs, sizeof(muse_lsf_cube *));
00622 unsigned char ifu, nlsfs = 0;
00623 for (ifu = 1; ifu <= kMuseNumIFUs; ifu++) {
00624 cpl_frameset *frames = muse_frameset_find(aProcessing->inframes,
00625 MUSE_TAG_LSF_PROFILE, ifu, CPL_FALSE);
00626 cpl_errorstate es = cpl_errorstate_get();
00627 cpl_frame *frame = cpl_frameset_get_position(frames, 0);
00628 if (frame == NULL) {
00629 cpl_msg_warning(__func__, "No %s (cube format) specified for IFU %2hhu!",
00630 MUSE_TAG_LSF_PROFILE, ifu);
00631 cpl_errorstate_set(es);
00632 cpl_frameset_delete(frames);
00633 continue;
00634 }
00635 const char *fn = cpl_frame_get_filename(frame);
00636 lsf[ifu-1] = muse_lsf_cube_load(fn, ifu);
00637 if (lsf[ifu-1] == NULL) {
00638 cpl_msg_warning(__func__, "Could not load LSF (cube format) for IFU %2hhu "
00639 "from \"%s\"!", ifu, fn);
00640 cpl_frameset_delete(frames);
00641 muse_lsf_cube_delete_all(lsf);
00642 return NULL;
00643 }
00644 muse_processing_append_used(aProcessing, frame, CPL_FRAME_GROUP_CALIB, 1);
00645 cpl_frameset_delete(frames);
00646 nlsfs++;
00647 }
00648 if (!nlsfs) {
00649 cpl_msg_error(__func__, "Did not load any %ss (cube format)!",
00650 MUSE_TAG_LSF_PROFILE);
00651 muse_lsf_cube_delete_all(lsf);
00652 return NULL;
00653 }
00654 cpl_msg_info(__func__, "Successfully loaded %s%hhu %ss (cube format).",
00655 nlsfs == kMuseNumIFUs ? "all ": "", nlsfs, MUSE_TAG_LSF_PROFILE);
00656 return lsf;
00657 }
00658
00662 void
00663 muse_lsf_cube_delete_all(muse_lsf_cube **aLsfCube) {
00664 if (aLsfCube != NULL) {
00665 int ifu;
00666 for (ifu = 0; ifu < kMuseNumIFUs; ifu++) {
00667 muse_lsf_cube_delete(aLsfCube[ifu]);
00668 }
00669 cpl_free(aLsfCube);
00670 }
00671 }
00672
00680 muse_lsf_cube *
00681 muse_lsf_cube_load(const char *aFileName, unsigned char aIFU)
00682 {
00683 cpl_ensure(aFileName, CPL_ERROR_NULL_INPUT, NULL);
00684
00685
00686 int ext = cpl_fits_find_extension(aFileName, "LSF_PROFILE");
00687 char *extname = NULL;
00688 if (ext <= 0) {
00689 extname = cpl_sprintf("CHAN%02hhu.LSF_PROFILE", aIFU);
00690 ext = cpl_fits_find_extension(aFileName, extname);
00691 if (ext <= 0) {
00692 cpl_free(extname);
00693 cpl_error_set(__func__, CPL_ERROR_DATA_NOT_FOUND);
00694 return NULL;
00695 }
00696 }
00697 cpl_free(extname);
00698
00699 muse_lsf_cube *lsf = cpl_calloc(1, sizeof(muse_lsf_cube));
00700 lsf->header = cpl_propertylist_load(aFileName, 0);
00701 lsf->img = cpl_imagelist_load(aFileName, CPL_TYPE_DOUBLE, ext);
00702 lsf->wcs = cpl_calloc(1, sizeof(muse_wcs));
00703 if (lsf->img == NULL) {
00704 muse_lsf_cube_delete(lsf);
00705 return NULL;
00706 }
00707
00708 cpl_propertylist *header = cpl_propertylist_load(aFileName, ext);
00709 if (header == NULL) {
00710 muse_lsf_cube_delete(lsf);
00711 return NULL;
00712 }
00713 lsf->wcs->cd11 = muse_pfits_get_cd(header, 1, 1);
00714 lsf->wcs->cd12 = muse_pfits_get_cd(header, 1, 2);
00715 lsf->wcs->cd21 = muse_pfits_get_cd(header, 2, 1);
00716 lsf->wcs->cd22 = muse_pfits_get_cd(header, 2, 2);
00717 lsf->wcs->crpix1 = muse_pfits_get_crpix(header, 1);
00718 lsf->wcs->crpix2 = muse_pfits_get_crpix(header, 2);
00719 lsf->wcs->crval1 = muse_pfits_get_crval(header, 1);
00720 lsf->wcs->crval2 = muse_pfits_get_crval(header, 2);
00721
00722
00723
00724
00725
00726 cpl_propertylist_delete(header);
00727
00728 return lsf;
00729 }
00730
00741 cpl_image *
00742 muse_lsf_average_cube_all(muse_lsf_cube **aLsfCube, muse_pixtable *aPixtable)
00743 {
00744 cpl_ensure(aLsfCube, CPL_ERROR_NULL_INPUT, NULL);
00745 cpl_image *img = NULL;
00746 cpl_size wsum = 0;
00747 cpl_size w[kMuseNumIFUs][kMuseSlicesPerCCD];
00748 cpl_size ifu;
00749 for (ifu = 0; ifu < kMuseNumIFUs; ifu++) {
00750 cpl_size i_slice;
00751 for (i_slice = 0; i_slice < kMuseSlicesPerCCD; i_slice++) {
00752 w[ifu][i_slice] = (aPixtable == NULL)?1.0:0.0;
00753 }
00754 }
00755
00756 if (aPixtable != NULL) {
00757 cpl_size n_rows = muse_pixtable_get_nrow(aPixtable);
00758 cpl_size i_row;
00759 uint32_t *origin = (uint32_t*)
00760 cpl_table_get_data_int(aPixtable->table, MUSE_PIXTABLE_ORIGIN);
00761 for (i_row = 0; i_row < n_rows; i_row++) {
00762 int i_ifu = muse_pixtable_origin_get_ifu(origin[i_row]);
00763 int i_slice = muse_pixtable_origin_get_slice(origin[i_row]);
00764 w[i_ifu-1][i_slice-1]++;
00765 }
00766 }
00767
00768 for (ifu = 0; ifu < kMuseNumIFUs; ifu++) {
00769 if (aLsfCube[ifu] == NULL) {
00770 continue;
00771 }
00772 cpl_size i_slice;
00773 cpl_size n_slices = cpl_imagelist_get_size(aLsfCube[ifu]->img);
00774 for (i_slice = 0; i_slice < n_slices; i_slice++) {
00775 if (w[ifu][i_slice] > 0) {
00776 cpl_image *m = cpl_imagelist_get(aLsfCube[ifu]->img, i_slice);
00777 m = cpl_image_duplicate(m);
00778 cpl_image_multiply_scalar(m, w[ifu][i_slice]);
00779 wsum += w[ifu][i_slice];
00780 if (img == NULL) {
00781 img = m;
00782 } else {
00783 cpl_errorstate pre = cpl_errorstate_get();
00784 cpl_error_code rc = cpl_image_add(img, m);
00785 cpl_image_delete(m);
00786 if (rc != CPL_ERROR_NONE) {
00787 cpl_msg_warning(__func__, "Could not add cube of IFU %"
00788 CPL_SIZE_FORMAT": %s",
00789 ifu+1, cpl_error_get_message());
00790 cpl_errorstate_set(pre);
00791 break;
00792 }
00793 }
00794 }
00795 }
00796 }
00797 if ((img != NULL) && (wsum > 0)) {
00798 cpl_image_divide_scalar(img, wsum);
00799 return img;
00800 } else {
00801 cpl_image_delete(img);
00802 return NULL;
00803 }
00804 }
00805
00806
00807 muse_wcs *
00808 muse_lsf_cube_get_wcs_all(muse_lsf_cube **aLsfCube) {
00809 cpl_size ifu;
00810 for (ifu = 0; ifu < kMuseNumIFUs; ifu++) {
00811 if (aLsfCube[ifu] != NULL) {
00812 return aLsfCube[ifu]->wcs;
00813 }
00814 }
00815 return NULL;
00816 }
00817
00818
00829
00830 cpl_error_code
00831 muse_lsf_fold_rectangle(cpl_image *aLsfImage, const muse_wcs *aWCS,
00832 double aBinWidth)
00833 {
00834 cpl_ensure_code(aLsfImage && aWCS, CPL_ERROR_NULL_INPUT);
00835 cpl_ensure_code(aBinWidth > 0., CPL_ERROR_ILLEGAL_INPUT);
00836
00837 cpl_size ncol_kernel = aBinWidth / aWCS->cd11;
00838 ncol_kernel = ((ncol_kernel + 1) / 2) * 2 + 1;
00839 double r = ncol_kernel - aBinWidth / aWCS->cd11;
00840 cpl_matrix *kernel = cpl_matrix_new(1, ncol_kernel);
00841 cpl_matrix_fill(kernel, 1.0);
00842 cpl_matrix_set(kernel, 0, 0, 1 - r/2);
00843 cpl_matrix_set(kernel, 0, ncol_kernel-1, 1 - r/2);
00844 cpl_image *lsfImage0 = cpl_image_duplicate(aLsfImage);
00845 cpl_image_filter(aLsfImage, lsfImage0, kernel, CPL_FILTER_LINEAR,
00846 CPL_BORDER_FILTER);
00847 cpl_matrix_delete(kernel);
00848 cpl_image_delete(lsfImage0);
00849
00850 return CPL_ERROR_NONE;
00851 }
00852
00853
00866
00867 cpl_error_code
00868 muse_lsf_apply(const cpl_image *aLsfImage, const muse_wcs *aWCS,
00869 cpl_array *aVal, double aLambda)
00870 {
00871 cpl_ensure_code(aLsfImage, CPL_ERROR_NULL_INPUT);
00872 cpl_ensure_code(aWCS, CPL_ERROR_NULL_INPUT);
00873 cpl_ensure_code(aVal, CPL_ERROR_NULL_INPUT);
00874 cpl_size n_lsf = cpl_image_get_size_x(aLsfImage);
00875 cpl_size n_lambda = cpl_image_get_size_y(aLsfImage);
00876
00877
00878 double lambda_pix = (aLambda - aWCS->crval2) / aWCS->cd22 + aWCS->crpix2;
00879 if (lambda_pix < 1) {
00880 lambda_pix = 1;
00881 }
00882 if (lambda_pix > n_lambda) {
00883 lambda_pix = n_lambda;
00884 }
00885 cpl_size i_lambda = floor(lambda_pix);
00886 lambda_pix -= i_lambda;
00887
00888
00889
00890 cpl_array *lsf_lambda = cpl_array_new(n_lsf + 4, CPL_TYPE_DOUBLE);
00891 cpl_array *lsf_values = cpl_array_new(n_lsf + 4, CPL_TYPE_DOUBLE);
00892 cpl_size i;
00893 for (i = 1; i <= n_lsf; i++) {
00894 int res;
00895 double z = cpl_image_get(aLsfImage, i, i_lambda, &res);
00896 if (lambda_pix > 0) {
00897 z = z * (1 - lambda_pix)
00898 + cpl_image_get(aLsfImage, i, i_lambda+1, &res) * lambda_pix;
00899 }
00900 cpl_array_set(lsf_values, i+1, z);
00901 cpl_array_set(lsf_lambda, i+1, aWCS->crval1 + (i - aWCS->crpix1) * aWCS->cd11);
00902 }
00903
00904
00905 cpl_array_set(lsf_lambda, 0, -10000);
00906 cpl_array_set(lsf_values, 0, 0);
00907 cpl_array_set(lsf_lambda, 1, aWCS->crval1 - aWCS->crpix1 * aWCS->cd11);
00908 cpl_array_set(lsf_values, 1, 0);
00909 cpl_array_set(lsf_lambda, n_lsf+2, aWCS->crval1 + (n_lsf+1 - aWCS->crpix1) * aWCS->cd11);
00910 cpl_array_set(lsf_values, n_lsf+2, 0);
00911 cpl_array_set(lsf_lambda, n_lsf+3, 10000);
00912 cpl_array_set(lsf_values, n_lsf+3, 0);
00913
00914
00915 cpl_array *v = cpl_array_duplicate(lsf_values);
00916 cpl_array_multiply(v, lsf_lambda);
00917 double offset = cpl_array_get_mean(v)/cpl_array_get_mean(lsf_values);
00918 cpl_array_delete(v);
00919 cpl_array_subtract_scalar (lsf_lambda, offset);
00920
00921
00922 double norm = cpl_array_get_mean(lsf_values) * (n_lsf + 4) * aWCS->cd11;
00923 cpl_array_divide_scalar(lsf_values, norm);
00924
00925 cpl_array *values = muse_cplarray_interpolate_linear(aVal, lsf_lambda,
00926 lsf_values);
00927
00928
00929 memcpy(cpl_array_get_data_double(aVal), cpl_array_get_data_double(values),
00930 cpl_array_get_size(aVal) * sizeof(double));
00931
00932 cpl_array_delete(values);
00933 cpl_array_delete(lsf_lambda);
00934 cpl_array_delete(lsf_values);
00935 return CPL_ERROR_NONE;
00936 }
00937