00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 #ifdef HAVE_CONFIG_H
00029 #include <config.h>
00030 #endif
00031
00032
00033
00034
00035
00036 #include <math.h>
00037 #include <cpl.h>
00038 #include <string.h>
00039
00040 #include "irplib_utils.h"
00041
00042 #include "hawki_alloc.h"
00043 #include "hawki_utils.h"
00044 #include "hawki_calib.h"
00045 #include "hawki_distortion.h"
00046 #include "hawki_load.h"
00047 #include "hawki_save.h"
00048 #include "hawki_pfits.h"
00049 #include "hawki_dfs.h"
00050 #include "irplib_cat.h"
00051 #include "irplib_stdstar.h"
00052 #include "irplib_match_cats.h"
00053 #include "hawki_match_cats.h"
00054
00055
00056
00057
00058
00059 static int hawki_cal_distortion_create(cpl_plugin *) ;
00060 static int hawki_cal_distortion_exec(cpl_plugin *) ;
00061 static int hawki_cal_distortion_destroy(cpl_plugin *) ;
00062 static int hawki_cal_distortion(cpl_parameterlist * parlist,
00063 cpl_frameset * frameset);
00064
00065
00066 static int hawki_cal_distortion_get_apertures_from_raw_distor
00067 (cpl_frameset * raw_target,
00068 const cpl_frame * flat,
00069 const cpl_frame * dark,
00070 const cpl_frame * bpm,
00071 cpl_image ** master_sky,
00072 double sigma_det,
00073 cpl_apertures *** apertures);
00074
00075 static int hawki_cal_distortion_load_master_calib
00076 (const cpl_frame * flat,
00077 const cpl_frame * dark,
00078 const cpl_frame * bpm,
00079 cpl_imagelist ** flat_images,
00080 cpl_imagelist ** dark_images,
00081 cpl_imagelist ** bpm_images);
00082
00083 static cpl_image ** hawki_cal_distortion_get_master_sky
00084 (cpl_frameset * raw_target,
00085 const cpl_frame * flat,
00086 const cpl_frame * dark,
00087 const cpl_frame * bpm);
00088
00089 static int hawki_cal_distortion_subtract_sky
00090 (cpl_imagelist * distor_corrected,
00091 cpl_image * master_sky);
00092
00093 static hawki_distortion ** hawki_cal_distortion_compute_dist_solution
00094 (cpl_apertures *** apertures,
00095 int nframes,
00096 cpl_bivector * offsets,
00097 int grid_points,
00098 int * nmatched_pairs,
00099 double * rms,
00100 hawki_distortion ** distortion_guess);
00101
00102 static cpl_apertures * hawki_cal_distortion_get_image_apertures
00103 (cpl_image * image,
00104 double sigma_det);
00105
00106 static int hawki_cal_distortion_fill_obj_pos
00107 (cpl_table * objects_positions,
00108 cpl_apertures * apertures);
00109
00110 static int hawki_cal_distortion_add_offset_to_positions
00111 (cpl_table ** objects_positions,
00112 cpl_bivector * offsets);
00113
00114 static int hawki_cal_distortion_fit_first_order_solution
00115 (hawki_distortion * distortion,
00116 cpl_polynomial * fit2d_x,
00117 cpl_polynomial * fit2d_y);
00118
00119 static cpl_propertylist ** hawki_cal_distortion_qc
00120 (hawki_distortion ** distortion,
00121 int * nmatched_pairs,
00122 double * rms);
00123
00124 static int hawki_cal_distortion_save
00125 (hawki_distortion ** distortion,
00126 cpl_parameterlist * parlist,
00127 cpl_propertylist ** qclists,
00128 cpl_frameset * recipe_set);
00129
00130 static int hawki_cal_distortion_retrieve_input_param
00131 (cpl_parameterlist * parlist);
00132
00133
00134
00135
00136
00137
00138 static struct
00139 {
00140
00141 double sigma_det;
00142 int grid_points;
00143 int borders;
00144 int subtract_linear;
00145 } hawki_cal_distortion_config;
00146
00147 static char hawki_cal_distortion_description[] =
00148 "hawki_cal_distortion -- HAWK-I distortion and astrometry autocalibration.\n\n"
00149 "The input files must be tagged:\n"
00150 "distortion_field.fits "HAWKI_CAL_DISTOR_RAW"\n"
00151 "sky_distortion.fits "HAWKI_CAL_DISTOR_SKY_RAW"\n"
00152 "flat-file.fits "HAWKI_CALPRO_FLAT" (optional)\n"
00153 "dark-file.fits "HAWKI_CALPRO_DARK" (optional)\n"
00154 "bpm-file.fits "HAWKI_CALPRO_BPM" (optional)\n\n"
00155 "The recipe creates as an output:\n"
00156 "hawki_cal_distortion_distx.fits ("HAWKI_CALPRO_DISTORTION_X") \n"
00157 "hawki_cal_distortion_disty.fits ("HAWKI_CALPRO_DISTORTION_Y") \n\n"
00158 "The recipe performs the following steps:\n"
00159 "-Basic calibration of astrometry fields\n"
00160 "-Autocalibration of distortion, using method in A&A 454,1029 (2006)\n\n"
00161 "Return code:\n"
00162 "esorex exits with an error code of 0 if the recipe completes successfully\n"
00163 "or 1 otherwise";
00164
00165
00166
00167
00168
00169
00177
00178 int cpl_plugin_get_info(cpl_pluginlist * list)
00179 {
00180 cpl_recipe * recipe = cpl_calloc(1, sizeof(*recipe)) ;
00181 cpl_plugin * plugin = &recipe->interface ;
00182
00183 cpl_plugin_init(plugin,
00184 CPL_PLUGIN_API,
00185 HAWKI_BINARY_VERSION,
00186 CPL_PLUGIN_TYPE_RECIPE,
00187 "hawki_cal_distortion",
00188 "Distortion autocalibration",
00189 hawki_cal_distortion_description,
00190 "Cesar Enrique Garcia Dabo",
00191 PACKAGE_BUGREPORT,
00192 hawki_get_license(),
00193 hawki_cal_distortion_create,
00194 hawki_cal_distortion_exec,
00195 hawki_cal_distortion_destroy);
00196
00197 cpl_pluginlist_append(list, plugin) ;
00198
00199 return 0;
00200 }
00201
00202
00211
00212 static int hawki_cal_distortion_create(cpl_plugin * plugin)
00213 {
00214 cpl_recipe * recipe ;
00215 cpl_parameter * p ;
00216
00217
00218 if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE)
00219 recipe = (cpl_recipe *)plugin ;
00220 else return -1 ;
00221
00222
00223 recipe->parameters = cpl_parameterlist_new() ;
00224
00225
00226
00227 p = cpl_parameter_new_value("hawki.hawki_cal_distortion.sigma_det",
00228 CPL_TYPE_DOUBLE, "detection level",
00229 "hawki.hawki_cal_distortion", 6.) ;
00230 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "sigma_det") ;
00231 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV) ;
00232 cpl_parameterlist_append(recipe->parameters, p) ;
00233
00234
00235 p = cpl_parameter_new_value("hawki.hawki_cal_distortion.grid_points",
00236 CPL_TYPE_INT,
00237 "number of points in distortion grid",
00238 "hawki.hawki_cal_distortion", 9) ;
00239 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "grid_points") ;
00240 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV) ;
00241 cpl_parameterlist_append(recipe->parameters, p) ;
00242
00243
00244 p = cpl_parameter_new_value("hawki.hawki_cal_distortion.borders",
00245 CPL_TYPE_INT,
00246 "number of pixels to trim at the borders",
00247 "hawki.hawki_cal_distortion", 6) ;
00248 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "borders") ;
00249 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV) ;
00250 cpl_parameterlist_append(recipe->parameters, p) ;
00251
00252
00253 p = cpl_parameter_new_value("hawki.hawki_cal_distortion.subtract_linear",
00254 CPL_TYPE_BOOL,
00255 "Subtract a linear term to the solution",
00256 "hawki.hawki_cal_distortion", TRUE) ;
00257 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "subtract_linear") ;
00258 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV) ;
00259 cpl_parameterlist_append(recipe->parameters, p) ;
00260
00261
00262 return 0;
00263 }
00264
00265
00271
00272 static int hawki_cal_distortion_exec(cpl_plugin * plugin)
00273 {
00274 cpl_recipe * recipe ;
00275
00276
00277 if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE)
00278 recipe = (cpl_recipe *)plugin ;
00279 else return -1 ;
00280
00281
00282 hawki_print_banner();
00283
00284 return hawki_cal_distortion(recipe->parameters, recipe->frames) ;
00285 }
00286
00287
00293
00294 static int hawki_cal_distortion_destroy(cpl_plugin * plugin)
00295 {
00296 cpl_recipe * recipe ;
00297
00298
00299 if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE)
00300 recipe = (cpl_recipe *)plugin ;
00301 else return -1 ;
00302
00303 cpl_parameterlist_delete(recipe->parameters) ;
00304 return 0 ;
00305 }
00306
00307
00313
00314 static int hawki_cal_distortion(cpl_parameterlist * parlist,
00315 cpl_frameset * frameset)
00316 {
00317 const cpl_frame * flat = NULL;
00318 const cpl_frame * dark = NULL;
00319 const cpl_frame * bpm = NULL;
00320 cpl_frameset * distorframes = NULL;
00321 cpl_frameset * skyframes = NULL;
00322 const cpl_frame * distorxguess = NULL;
00323 const cpl_frame * distoryguess = NULL;
00324 hawki_distortion ** distortionguess = NULL;
00325 hawki_distortion ** distortion = NULL;
00326 cpl_propertylist ** qclists = NULL;
00327 cpl_image ** master_sky = NULL;
00328 cpl_bivector * nominal_offsets = NULL;
00329 cpl_apertures ** apertures[HAWKI_NB_DETECTORS];
00330 int nmatched_pairs[HAWKI_NB_DETECTORS];
00331 double rms[HAWKI_NB_DETECTORS];
00332 int idet;
00333 int ioff;
00334 int iframe;
00335 int nframes;
00336
00337
00338
00339 if(hawki_cal_distortion_retrieve_input_param(parlist))
00340 {
00341 cpl_msg_error(__func__, "Wrong parameters");
00342 return -1;
00343 }
00344
00345
00346 if (hawki_dfs_set_groups(frameset)) {
00347 cpl_msg_error(__func__, "Cannot identify RAW and CALIB frames") ;
00348 return -1 ;
00349 }
00350
00351
00352 cpl_msg_info(__func__, "Identifying calibration data");
00353 flat = cpl_frameset_find_const(frameset, HAWKI_CALPRO_FLAT);
00354 dark = cpl_frameset_find_const(frameset, HAWKI_CALPRO_DARK);
00355 bpm = cpl_frameset_find_const(frameset, HAWKI_CALPRO_BPM);
00356
00357
00358 cpl_msg_info(__func__, "Identifying distortion and sky data");
00359 distorframes = hawki_extract_frameset(frameset, HAWKI_CAL_DISTOR_RAW) ;
00360 if (distorframes == NULL)
00361 {
00362 cpl_msg_error(__func__, "Distortion images have to be provided (%s)",
00363 HAWKI_CAL_DISTOR_RAW);
00364 return -1 ;
00365 }
00366
00367 skyframes = hawki_extract_frameset(frameset, HAWKI_CAL_DISTOR_SKY_RAW) ;
00368 if (skyframes == NULL)
00369 {
00370 cpl_msg_error(__func__, "Sky images have to be provided (%s)",
00371 HAWKI_CAL_DISTOR_SKY_RAW);
00372 cpl_frameset_delete(distorframes);
00373 return -1 ;
00374 }
00375
00376 distorxguess = cpl_frameset_find_const(frameset, HAWKI_CALPRO_DISTORTION_X);
00377 distoryguess = cpl_frameset_find_const(frameset, HAWKI_CALPRO_DISTORTION_Y);
00378 if(distorxguess != NULL && distoryguess != NULL)
00379 {
00380
00381 }
00382
00383
00384 cpl_msg_info(__func__, "Computing the master sky image");
00385 master_sky = hawki_cal_distortion_get_master_sky(skyframes, flat, dark, bpm);
00386 if(master_sky == NULL)
00387 {
00388 cpl_msg_error(__func__, "Cannot get master sky image") ;
00389 cpl_frameset_delete(distorframes);
00390 cpl_frameset_delete(skyframes);
00391 return -1;
00392 }
00393
00394
00395 cpl_msg_info(__func__, "Getting objects from distortion images");
00396 if(hawki_cal_distortion_get_apertures_from_raw_distor
00397 (distorframes, flat, dark, bpm, master_sky,
00398 hawki_cal_distortion_config.sigma_det, apertures) == -1)
00399 {
00400 cpl_msg_error(__func__,
00401 "Cannot get objects from distortion images");
00402 for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
00403 cpl_image_delete(master_sky[idet]);
00404 cpl_free(master_sky);
00405 cpl_frameset_delete(distorframes);
00406 cpl_frameset_delete(skyframes);
00407 return -1 ;
00408 }
00409 for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
00410 cpl_image_delete(master_sky[idet]);
00411 cpl_free(master_sky);
00412
00413
00414 cpl_msg_info(__func__,"Getting the nominal offsets");
00415 nominal_offsets = hawki_get_header_tel_offsets(distorframes);
00416 if (nominal_offsets == NULL)
00417 {
00418 cpl_msg_error(__func__, "Cannot load the header offsets") ;
00419 cpl_frameset_delete(distorframes);
00420 cpl_frameset_delete(skyframes);
00421 return -1;
00422 }
00423
00424
00425
00426 cpl_vector_multiply_scalar(cpl_bivector_get_x(nominal_offsets), -1.0);
00427 cpl_vector_multiply_scalar(cpl_bivector_get_y(nominal_offsets), -1.0);
00428
00429
00430 cpl_msg_indent_more();
00431 for (ioff=0 ; ioff<cpl_bivector_get_size(nominal_offsets) ; ioff++)
00432 {
00433 cpl_msg_info(__func__, "Telescope offsets (Frame %d): %g %g", ioff+1,
00434 cpl_bivector_get_x_data(nominal_offsets)[ioff],
00435 cpl_bivector_get_y_data(nominal_offsets)[ioff]);
00436 }
00437 cpl_msg_indent_less();
00438
00439
00440 cpl_msg_info(__func__, "Computing the distortion");
00441 nframes = cpl_frameset_get_size(distorframes);
00442 distortion = hawki_cal_distortion_compute_dist_solution
00443 (apertures, nframes, nominal_offsets,
00444 hawki_cal_distortion_config.grid_points,
00445 nmatched_pairs, rms,
00446 distortionguess);
00447 cpl_bivector_delete(nominal_offsets);
00448 if(distortion == NULL)
00449 {
00450 cpl_frameset_delete(distorframes);
00451 cpl_frameset_delete(skyframes);
00452 return -1;
00453 }
00454
00455
00456 qclists = hawki_cal_distortion_qc(distortion, nmatched_pairs, rms);
00457
00458
00459 cpl_msg_info(__func__,"Saving products");
00460 if(hawki_cal_distortion_save(distortion,
00461 parlist, qclists, frameset) == -1)
00462 {
00463 cpl_msg_error(__func__,"Could not save products");
00464 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
00465 cpl_propertylist_delete(qclists[idet]);
00466 cpl_free(qclists);
00467 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
00468 hawki_distortion_delete(distortion[idet]);
00469 cpl_free(distortion);
00470 cpl_frameset_delete(distorframes);
00471 cpl_frameset_delete(skyframes);
00472 return -1;
00473 }
00474
00475
00476 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
00477 cpl_propertylist_delete(qclists[idet]);
00478 cpl_free(qclists);
00479 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
00480 hawki_distortion_delete(distortion[idet]);
00481 cpl_free(distortion);
00482 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
00483 {
00484 for(iframe = 0 ; iframe < nframes; iframe++)
00485 cpl_apertures_delete(apertures[idet][iframe]);
00486 cpl_free(apertures[idet]);
00487 }
00488 cpl_frameset_delete(distorframes);
00489 cpl_frameset_delete(skyframes);
00490
00491
00492
00493 if (cpl_error_get_code())
00494 {
00495 cpl_msg_error(__func__,
00496 "HAWK-I pipeline could not recover from previous errors");
00497 return -1 ;
00498 }
00499 else return 0 ;
00500 }
00501
00502
00513
00514 static int hawki_cal_distortion_load_master_calib
00515 (const cpl_frame * flat,
00516 const cpl_frame * dark,
00517 const cpl_frame * bpm,
00518 cpl_imagelist ** flat_images,
00519 cpl_imagelist ** dark_images,
00520 cpl_imagelist ** bpm_images)
00521 {
00522 cpl_errorstate error_prevstate = cpl_errorstate_get();
00523
00524
00525 *flat_images = NULL;
00526 *dark_images = NULL;
00527 *bpm_images = NULL;
00528
00529
00530 cpl_msg_info(__func__, "Loading the calibration data") ;
00531 if(flat != NULL)
00532 {
00533 *flat_images = hawki_load_frame(flat, CPL_TYPE_FLOAT);
00534 if(flat_images == NULL)
00535 {
00536 cpl_msg_error(__func__, "Error reading flat") ;
00537 return -1;
00538 }
00539 }
00540 if(dark != NULL)
00541 {
00542 *dark_images = hawki_load_frame(dark, CPL_TYPE_FLOAT);
00543 if(dark_images == NULL)
00544 {
00545 cpl_msg_error(__func__, "Error reading dark") ;
00546 cpl_imagelist_delete(*flat_images);
00547 return -1;
00548 }
00549 }
00550 if(bpm != NULL)
00551 {
00552 *bpm_images = hawki_load_frame(bpm, CPL_TYPE_INT);
00553 if(bpm_images == NULL)
00554 {
00555 cpl_msg_error(__func__, "Error reading bpm") ;
00556 cpl_imagelist_delete(*flat_images);
00557 cpl_imagelist_delete(*dark_images);
00558 return -1;
00559 }
00560 }
00561
00562 if(!cpl_errorstate_is_equal(error_prevstate ))
00563 {
00564 cpl_msg_error(__func__, "A problem happened loading calibration");
00565 cpl_imagelist_delete(*flat_images);
00566 cpl_imagelist_delete(*dark_images);
00567 cpl_imagelist_delete(*bpm_images);
00568 return -1;
00569 }
00570 return 0;
00571 }
00572
00573
00583
00584 static int hawki_cal_distortion_get_apertures_from_raw_distor
00585 (cpl_frameset * raw_distor,
00586 const cpl_frame * flat,
00587 const cpl_frame * dark,
00588 const cpl_frame * bpm,
00589 cpl_image ** master_sky,
00590 double sigma_det,
00591 cpl_apertures *** apertures)
00592 {
00593 cpl_imagelist * flat_images;
00594 cpl_imagelist * dark_images;
00595 cpl_imagelist * bpm_images;
00596 cpl_propertylist * plist;
00597
00598 double science_dit;
00599 int iframe;
00600 int nframes;
00601 int idet;
00602
00603 cpl_errorstate error_prevstate = cpl_errorstate_get();
00604
00605
00606 cpl_msg_indent_more();
00607
00608
00609 hawki_cal_distortion_load_master_calib
00610 (flat, dark, bpm, &flat_images, &dark_images, &bpm_images);
00611
00612
00613 if(dark != NULL)
00614 {
00615 if ((plist=cpl_propertylist_load
00616 (cpl_frame_get_filename
00617 (cpl_frameset_get_first_const(raw_distor)), 0)) == NULL)
00618 {
00619 cpl_msg_error(__func__, "Cannot get header from frame");
00620 cpl_imagelist_delete(flat_images);
00621 cpl_imagelist_delete(dark_images);
00622 return -1;
00623 }
00624 science_dit = hawki_pfits_get_dit(plist);
00625 cpl_imagelist_multiply_scalar(dark_images, science_dit);
00626 cpl_propertylist_delete(plist);
00627 }
00628
00629
00630 nframes = cpl_frameset_get_size(raw_distor);
00631 for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
00632 {
00633 cpl_imagelist * distor_serie;
00634 cpl_imagelist * distor_serie_trimmed;
00635
00636 cpl_msg_info(__func__, "Working on detector %d", idet + 1);
00637 cpl_msg_indent_more();
00638
00639
00640 cpl_msg_info(__func__, "Loading distortion images");
00641 distor_serie = hawki_load_detector(raw_distor, idet + 1, CPL_TYPE_FLOAT);
00642 if(distor_serie== NULL)
00643 {
00644 cpl_msg_error(__func__, "Error reading distortion images");
00645 return -1;
00646 }
00647
00648
00649 cpl_msg_info(__func__, "Applying basic calibration") ;
00650 cpl_msg_indent_more();
00651 if (hawki_flat_dark_bpm_detector_calib
00652 (distor_serie,
00653 cpl_imagelist_get(flat_images, idet),
00654 cpl_imagelist_get(dark_images, idet),
00655 cpl_imagelist_get(bpm_images, idet)) == -1)
00656 {
00657 cpl_msg_error(__func__, "Cannot calibrate frame") ;
00658 cpl_imagelist_delete(flat_images);
00659 cpl_imagelist_delete(dark_images);
00660 cpl_imagelist_delete(bpm_images);
00661 cpl_imagelist_delete(distor_serie);
00662 cpl_msg_indent_less() ;
00663 return -1;
00664 }
00665 cpl_msg_indent_less();
00666
00667
00668 if (hawki_cal_distortion_config.borders > 0)
00669 {
00670 distor_serie_trimmed = hawki_trim_detector_calib(distor_serie,
00671 hawki_cal_distortion_config.borders);
00672 cpl_imagelist_delete(distor_serie);
00673 }
00674 else
00675 distor_serie_trimmed = distor_serie;
00676
00677
00678 cpl_msg_info(__func__, "Subtracting master sky") ;
00679 if(hawki_cal_distortion_subtract_sky(distor_serie_trimmed, master_sky[idet]) == -1)
00680 {
00681 cpl_msg_error(__func__, "Cannot subtract the sky") ;
00682 cpl_imagelist_delete(flat_images);
00683 cpl_imagelist_delete(dark_images);
00684 cpl_imagelist_delete(bpm_images);
00685 cpl_imagelist_delete(distor_serie_trimmed);
00686 return -1;
00687 }
00688
00689
00690 apertures[idet] = cpl_malloc(sizeof(*(apertures[idet])) * nframes);
00691 for(iframe = 0 ; iframe < nframes; iframe++)
00692 {
00693 cpl_image * this_image = cpl_imagelist_get(distor_serie_trimmed, iframe);
00694 cpl_msg_info(__func__,"Working with distortion image %d", iframe+1);
00695 cpl_msg_indent_more();
00696 apertures[idet][iframe] =
00697 hawki_cal_distortion_get_image_apertures(this_image, sigma_det);
00698 cpl_msg_indent_less();
00699 }
00700
00701
00702 cpl_imagelist_delete(distor_serie_trimmed);
00703 cpl_msg_indent_less();
00704 }
00705 cpl_imagelist_delete(flat_images);
00706 cpl_imagelist_delete(dark_images);
00707 cpl_imagelist_delete(bpm_images);
00708
00709
00710 if(!cpl_errorstate_is_equal(error_prevstate ))
00711 {
00712 cpl_msg_error(__func__, "A problem happened with basic calibration");
00713 return -1 ;
00714 }
00715 return 0;
00716 }
00717
00718 static cpl_image ** hawki_cal_distortion_get_master_sky
00719 (cpl_frameset * raw_sky_frames,
00720 const cpl_frame * flat,
00721 const cpl_frame * dark,
00722 const cpl_frame * bpm)
00723 {
00724 cpl_propertylist * plist;
00725 double science_dit;
00726 int idet;
00727 int jdet;
00728 cpl_imagelist * flat_images;
00729 cpl_imagelist * dark_images;
00730 cpl_imagelist * bpm_images;
00731 cpl_image ** bkg_images = NULL;
00732
00733 cpl_errorstate error_prevstate = cpl_errorstate_get();
00734
00735
00736 cpl_msg_indent_more();
00737
00738
00739 hawki_cal_distortion_load_master_calib
00740 (flat, dark, bpm, &flat_images, &dark_images, &bpm_images);
00741
00742
00743 if(dark != NULL)
00744 {
00745 if ((plist=cpl_propertylist_load
00746 (cpl_frame_get_filename
00747 (cpl_frameset_get_first_const(raw_sky_frames)), 0)) == NULL)
00748 {
00749 cpl_msg_error(__func__, "Cannot get header from frame");
00750 cpl_imagelist_delete(flat_images);
00751 cpl_imagelist_delete(dark_images);
00752 cpl_imagelist_delete(bpm_images);
00753 return NULL;
00754 }
00755 science_dit = hawki_pfits_get_dit(plist);
00756 cpl_imagelist_multiply_scalar(dark_images, science_dit);
00757 cpl_propertylist_delete(plist);
00758 }
00759
00760
00761 bkg_images = cpl_malloc(HAWKI_NB_DETECTORS * sizeof(*bkg_images));
00762 for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
00763 {
00764 cpl_imagelist * sky_serie;
00765 cpl_imagelist * sky_serie_trimmed;
00766
00767
00768 sky_serie = hawki_load_detector(raw_sky_frames, idet + 1, CPL_TYPE_FLOAT);
00769 if(sky_serie== NULL)
00770 {
00771 cpl_msg_error(__func__, "Error reading object image") ;
00772 return NULL;
00773 }
00774
00775
00776 cpl_msg_info(__func__, "Working on detector %d", idet + 1);
00777 cpl_msg_indent_more();
00778 if (hawki_flat_dark_bpm_detector_calib
00779 (sky_serie,
00780 cpl_imagelist_get(flat_images, idet),
00781 cpl_imagelist_get(dark_images, idet),
00782 cpl_imagelist_get(bpm_images, idet)) == -1)
00783 {
00784 cpl_msg_error(__func__, "Cannot calibrate frame") ;
00785 cpl_imagelist_delete(flat_images);
00786 cpl_imagelist_delete(dark_images);
00787 cpl_imagelist_delete(bpm_images);
00788 cpl_imagelist_delete(sky_serie);
00789 for(jdet = 0; jdet < idet; ++jdet)
00790 cpl_image_delete(bkg_images[jdet]);
00791 cpl_free(bkg_images);
00792 cpl_msg_indent_less() ;
00793 return NULL;
00794 }
00795
00796
00797 if (hawki_cal_distortion_config.borders > 0)
00798 {
00799 sky_serie_trimmed = hawki_trim_detector_calib(sky_serie,
00800 hawki_cal_distortion_config.borders);
00801 cpl_imagelist_delete(sky_serie);
00802 }
00803 else
00804 sky_serie_trimmed = sky_serie;
00805
00806
00807 if ((bkg_images[idet] =
00808 cpl_imagelist_collapse_median_create(sky_serie_trimmed)) == NULL)
00809 {
00810 cpl_msg_error(__func__, "Cannot compute the median of obj images");
00811 cpl_imagelist_delete(flat_images);
00812 cpl_imagelist_delete(dark_images);
00813 cpl_imagelist_delete(bpm_images);
00814 cpl_imagelist_delete(sky_serie_trimmed);
00815 for(jdet = 0; jdet < idet; ++jdet)
00816 cpl_image_delete(bkg_images[jdet]);
00817 cpl_free(bkg_images);
00818 cpl_msg_indent_less();
00819 return NULL;
00820 }
00821 cpl_imagelist_delete(sky_serie_trimmed);
00822 cpl_msg_indent_less();
00823 }
00824
00825
00826 for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
00827 cpl_image_subtract_scalar(bkg_images[idet],
00828 cpl_image_get_median(bkg_images[idet]));
00829
00830
00831 cpl_msg_indent_less();
00832 cpl_imagelist_delete(flat_images);
00833 cpl_imagelist_delete(dark_images);
00834 cpl_imagelist_delete(bpm_images);
00835
00836 if(!cpl_errorstate_is_equal(error_prevstate ))
00837 {
00838 cpl_msg_error(__func__, "A problem happened with basic calibration");
00839 for(idet = 0; idet < HAWKI_NB_DETECTORS; ++idet)
00840 cpl_image_delete(bkg_images[idet]);
00841 cpl_free(bkg_images);
00842 return NULL ;
00843 }
00844 return bkg_images;
00845 }
00846
00847 static int hawki_cal_distortion_subtract_sky
00848 (cpl_imagelist * distor_corrected,
00849 cpl_image * master_sky)
00850 {
00851 cpl_errorstate error_prevstate = cpl_errorstate_get();
00852
00853
00854 int idist, ndistor;
00855 ndistor = cpl_imagelist_get_size(distor_corrected);
00856 for(idist = 0; idist < ndistor; ++idist)
00857 {
00858 cpl_image * target_image =
00859 cpl_imagelist_get(distor_corrected, idist);
00860
00861 cpl_image_subtract_scalar
00862 (target_image, cpl_image_get_median(target_image));
00863
00864 if (cpl_image_subtract
00865 (target_image, master_sky)!=CPL_ERROR_NONE)
00866 {
00867 cpl_msg_error(cpl_func,"Cannot apply the bkg to the images");
00868 return -1 ;
00869 }
00870 }
00871
00872
00873 if(!cpl_errorstate_is_equal(error_prevstate ))
00874 {
00875 cpl_msg_error(__func__, "A problem happened with sky subtraction");
00876 return -1;
00877 }
00878 return 0;
00879 }
00880
00881
00882
00892
00893 static hawki_distortion ** hawki_cal_distortion_compute_dist_solution
00894 (cpl_apertures *** apertures,
00895 int nframes,
00896 cpl_bivector * offsets,
00897 int grid_points,
00898 int * nmatched_pairs,
00899 double * rms,
00900 hawki_distortion ** distortion_guess)
00901 {
00902 int idet;
00903 hawki_distortion ** distortion = NULL;
00904
00905 cpl_errorstate error_prevstate = cpl_errorstate_get();
00906
00907
00908 distortion = cpl_malloc(HAWKI_NB_DETECTORS * sizeof(*distortion));
00909
00910
00911 cpl_msg_indent_more();
00912 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
00913 {
00914 cpl_table ** obj_pos;
00915 cpl_table ** obj_pos_offset;
00916 int iframe;
00917 cpl_table * matches;
00918 hawki_distortion * dist_guess;
00919 cpl_polynomial * fit2d_x = NULL;
00920 cpl_polynomial * fit2d_y = NULL;
00921
00922
00923 cpl_msg_info(__func__, "Working on detector %d", idet+1);
00924 cpl_msg_indent_more();
00925
00926
00927 obj_pos =
00928 cpl_malloc(sizeof(*obj_pos) * nframes);
00929 obj_pos_offset =
00930 cpl_malloc(sizeof(*obj_pos_offset) * nframes);
00931 for(iframe = 0; iframe < nframes; ++iframe)
00932 {
00933 obj_pos[iframe] = cpl_table_new(0);
00934 cpl_table_new_column
00935 (obj_pos[iframe], HAWKI_COL_OBJ_POSX, CPL_TYPE_DOUBLE);
00936 cpl_table_new_column
00937 (obj_pos[iframe], HAWKI_COL_OBJ_POSY, CPL_TYPE_DOUBLE);
00938 }
00939
00940
00941 for(iframe = 0 ; iframe < nframes ; ++iframe)
00942 {
00943 cpl_apertures * this_apertures;
00944
00945
00946 this_apertures = apertures[idet][iframe];
00947
00948
00949 hawki_cal_distortion_fill_obj_pos(obj_pos[iframe],
00950 this_apertures);
00951 obj_pos_offset[iframe] = cpl_table_duplicate(obj_pos[iframe]);
00952 }
00953
00954
00955 hawki_cal_distortion_add_offset_to_positions
00956 (obj_pos_offset, offsets);
00957
00958
00959
00960
00961
00962 cpl_msg_info(__func__, "Matching all catalogs (may take a while)");
00963 matches = irplib_match_cat_pairs(obj_pos_offset, nframes,
00964 hawki_match_condition_5_pix);
00965 for(iframe = 0; iframe < nframes; ++iframe)
00966 cpl_table_delete(obj_pos_offset[iframe]);
00967 cpl_free(obj_pos_offset);
00968 if(matches == NULL)
00969 {
00970 cpl_msg_error(__func__, "Cannot match objects ");
00971 for(iframe = 0; iframe < nframes; ++iframe)
00972 cpl_table_delete(obj_pos[iframe]);
00973 cpl_free(obj_pos);
00974 return NULL;
00975 }
00976 cpl_msg_info(__func__,"Number of matching pairs %d",
00977 cpl_table_get_nrow(matches));
00978 nmatched_pairs[idet] = cpl_table_get_nrow(matches);
00979
00980
00981
00982 cpl_msg_info(__func__, "Computing distortion with the matched objects");
00983 cpl_msg_info(__func__, " (This step will take a long time to run)");
00984 if(distortion_guess != NULL)
00985 dist_guess = distortion_guess[idet];
00986 else
00987 dist_guess = NULL;
00988 distortion[idet] = hawki_distortion_compute_solution
00989 ((const cpl_table **)obj_pos, offsets, matches,
00990 nframes, HAWKI_DET_NPIX_X , HAWKI_DET_NPIX_Y, grid_points,
00991 dist_guess, rms + idet);
00992 if(distortion[idet] == NULL)
00993 {
00994 int jdet;
00995 cpl_msg_error(__func__,"Could not get the distortion");
00996 for(iframe = 0; iframe < nframes; ++iframe)
00997 cpl_table_delete(obj_pos[iframe]);
00998 cpl_free(obj_pos);
00999 for(jdet = 0; jdet < idet; ++jdet)
01000 hawki_distortion_delete(distortion[idet]);
01001 cpl_table_delete(matches);
01002 return NULL;
01003 }
01004
01005
01006 if(hawki_cal_distortion_config.subtract_linear)
01007 {
01008 cpl_msg_info(__func__,"Subtracting first order polynomial");
01009 fit2d_x = cpl_polynomial_new(2);
01010 fit2d_y = cpl_polynomial_new(2);
01011 hawki_cal_distortion_fit_first_order_solution
01012 (distortion[idet], fit2d_x, fit2d_y);
01013 }
01014
01015
01016 for(iframe = 0; iframe < nframes; ++iframe)
01017 cpl_table_delete(obj_pos[iframe]);
01018 cpl_free(obj_pos);
01019 if(hawki_cal_distortion_config.subtract_linear)
01020 {
01021 cpl_polynomial_delete(fit2d_x);
01022 cpl_polynomial_delete(fit2d_y);
01023 }
01024 cpl_table_delete(matches);
01025 cpl_msg_indent_less();
01026 }
01027
01028 if(!cpl_errorstate_is_equal(error_prevstate ))
01029 {
01030 cpl_msg_error(__func__, "A problem happened computing the distortion");
01031 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
01032 hawki_distortion_delete(distortion[idet]);
01033 return NULL ;
01034 }
01035
01036 return distortion;
01037 }
01038
01039 static cpl_apertures * hawki_cal_distortion_get_image_apertures
01040 (cpl_image * image,
01041 double sigma_det)
01042 {
01043 cpl_apertures * apertures = NULL;
01044 cpl_matrix * kernel = NULL;
01045 cpl_mask * object_mask = NULL;
01046 cpl_image * labels = NULL;
01047 int nobj;
01048 double median;
01049 double med_dist;
01050 double threshold;
01051
01052
01053 median = cpl_image_get_median_dev(image, &med_dist);
01054 threshold = median + sigma_det * med_dist;
01055 cpl_msg_info(__func__,"Detection threshold: %f", threshold);
01056
01057
01058 object_mask = cpl_mask_threshold_image_create
01059 (image, threshold, DBL_MAX);
01060 if (object_mask == NULL)
01061 return NULL;
01062
01063
01064 kernel = cpl_matrix_new(3, 3);
01065 cpl_matrix_fill(kernel, 1.00) ;
01066 if (cpl_mask_opening(object_mask, kernel) != CPL_ERROR_NONE) {
01067 cpl_mask_delete(object_mask) ;
01068 cpl_matrix_delete(kernel) ;
01069 return NULL;
01070 }
01071 cpl_matrix_delete(kernel);
01072
01073
01074 labels = cpl_image_labelise_mask_create(object_mask, &nobj);
01075 if (labels == NULL)
01076 {
01077 cpl_mask_delete(object_mask);
01078 return NULL;
01079 }
01080 cpl_mask_delete(object_mask);
01081 cpl_msg_info(__func__, "Number of objects detected: %d", nobj) ;
01082
01083
01084 apertures = cpl_apertures_new_from_image(image, labels);
01085 if (apertures == NULL)
01086 {
01087 cpl_image_delete(labels);
01088 return NULL;
01089 }
01090 cpl_image_delete(labels);
01091 return apertures;
01092 }
01093
01094 static int hawki_cal_distortion_fill_obj_pos
01095 (cpl_table * objects_positions,
01096 cpl_apertures * apertures)
01097 {
01098 int nobjs;
01099 int iobj;
01100 double border_off = 0;
01101
01102
01103 if(hawki_cal_distortion_config.borders > 0)
01104 border_off = hawki_cal_distortion_config.borders;
01105
01106 nobjs = cpl_apertures_get_size(apertures);
01107 cpl_table_set_size(objects_positions, nobjs);
01108
01109 for (iobj=0 ; iobj<nobjs ; iobj++)
01110 {
01111
01112 cpl_table_set_double(objects_positions, HAWKI_COL_OBJ_POSX, iobj,
01113 cpl_apertures_get_centroid_x(apertures,
01114 iobj+1) + border_off);
01115 cpl_table_set_double(objects_positions, HAWKI_COL_OBJ_POSY, iobj,
01116 cpl_apertures_get_centroid_y(apertures,
01117 iobj+1) + border_off);
01118 }
01119
01120 return 0;
01121 }
01122
01123 static int hawki_cal_distortion_add_offset_to_positions
01124 (cpl_table ** objects_positions,
01125 cpl_bivector * offsets)
01126 {
01127 int nframes;
01128 int iframe;
01129 int nobjs;
01130 int iobj;
01131
01132 nframes = cpl_bivector_get_size(offsets);
01133
01134 for(iframe = 0 ; iframe < nframes ; ++iframe)
01135 {
01136 double offset_x;
01137 double offset_y;
01138 int null;
01139 offset_x = cpl_bivector_get_x_data(offsets)[iframe];
01140 offset_y = cpl_bivector_get_y_data(offsets)[iframe];
01141 nobjs = cpl_table_get_nrow(objects_positions[iframe]);
01142 for (iobj=0 ; iobj<nobjs ; iobj++)
01143 {
01144 cpl_table_set_double(objects_positions[iframe],
01145 HAWKI_COL_OBJ_POSX, iobj,
01146 cpl_table_get_double(objects_positions[iframe],
01147 HAWKI_COL_OBJ_POSX, iobj, &null) + offset_x);
01148 cpl_table_set_double(objects_positions[iframe],
01149 HAWKI_COL_OBJ_POSY, iobj,
01150 cpl_table_get_double(objects_positions[iframe],
01151 HAWKI_COL_OBJ_POSY, iobj, &null) + offset_y);
01152 }
01153 }
01154
01155 return 0;
01156 }
01157
01158 static int hawki_cal_distortion_fit_first_order_solution
01159 (hawki_distortion * distortion,
01160 cpl_polynomial * fit2d_x,
01161 cpl_polynomial * fit2d_y)
01162 {
01163 cpl_matrix * pixel_pos;
01164 cpl_vector * dist_x_val;
01165 cpl_vector * dist_y_val;
01166 int nx;
01167 int ny;
01168 int i;
01169 int j;
01170 int null;
01171 const int mindeg2d[] = {0, 0};
01172 const int maxdeg2d[] = {1, 1};
01173 cpl_errorstate error_prevstate = cpl_errorstate_get();
01174 cpl_vector * pix;
01175 cpl_image * dist_x_plane;
01176 cpl_image * dist_y_plane;
01177 double dist_x_mean;
01178 double dist_y_mean;
01179
01180
01181 nx = hawki_distortion_get_size_x(distortion);
01182 ny = hawki_distortion_get_size_y(distortion);
01183 pixel_pos = cpl_matrix_new(2, nx * ny);
01184 dist_x_val = cpl_vector_new(nx*ny);
01185 dist_y_val = cpl_vector_new(nx*ny);
01186 for(i = 0; i < nx; ++i)
01187 for(j = 0; j < ny; ++j)
01188 {
01189 cpl_matrix_set(pixel_pos, 0, i + nx * j, (double)i);
01190 cpl_matrix_set(pixel_pos, 1, i + nx * j, (double)j);
01191 cpl_vector_set(dist_x_val, i + nx * j,
01192 cpl_image_get(distortion->dist_x, i+1, j+1, &null));
01193 cpl_vector_set(dist_y_val, i + nx * j,
01194 cpl_image_get(distortion->dist_y, i+1, j+1, &null));
01195 }
01196
01197
01198 cpl_polynomial_fit(fit2d_x, pixel_pos, NULL, dist_x_val,
01199 NULL, CPL_FALSE, mindeg2d, maxdeg2d);
01200 cpl_polynomial_fit(fit2d_y, pixel_pos, NULL, dist_y_val,
01201 NULL, CPL_FALSE, mindeg2d, maxdeg2d);
01202
01203 cpl_polynomial_set_coeff(fit2d_x, mindeg2d, 0.);
01204 cpl_polynomial_set_coeff(fit2d_y, mindeg2d, 0.);
01205
01206
01207 pix = cpl_vector_new(2);
01208 dist_x_plane = cpl_image_new(nx,ny,cpl_image_get_type(distortion->dist_x));
01209 dist_y_plane = cpl_image_new(nx,ny,cpl_image_get_type(distortion->dist_y));
01210 for(i = 0; i < nx; ++i)
01211 for(j = 0; j < ny; ++j)
01212 {
01213 double fit_value_x;
01214 double fit_value_y;
01215 cpl_vector_set(pix, 0, (double)i);
01216 cpl_vector_set(pix, 1, (double)j);
01217 fit_value_x = cpl_polynomial_eval(fit2d_x, pix);
01218 fit_value_y = cpl_polynomial_eval(fit2d_y, pix);
01219 cpl_image_set(dist_x_plane, i+1, j+1, fit_value_x);
01220 cpl_image_set(dist_y_plane, i+1, j+1, fit_value_y);
01221 }
01222 cpl_image_subtract(distortion->dist_x, dist_x_plane);
01223 cpl_image_subtract(distortion->dist_y, dist_y_plane);
01224
01225
01226 dist_x_mean = cpl_image_get_mean(distortion->dist_x);
01227 dist_y_mean = cpl_image_get_mean(distortion->dist_y);
01228 cpl_msg_warning(__func__,"Subtracting mean distortion in X %f",dist_x_mean);
01229 cpl_msg_warning(__func__,"Subtracting mean distortion in Y %f",dist_y_mean);
01230 cpl_image_subtract_scalar(distortion->dist_x, dist_x_mean);
01231 cpl_image_subtract_scalar(distortion->dist_y, dist_y_mean);
01232
01233
01234 cpl_matrix_delete(pixel_pos);
01235 cpl_vector_delete(dist_x_val);
01236 cpl_vector_delete(dist_y_val);
01237 cpl_vector_delete(pix);
01238 cpl_image_delete(dist_x_plane);
01239 cpl_image_delete(dist_y_plane);
01240
01241 if(!cpl_errorstate_is_equal(error_prevstate ))
01242 {
01243 cpl_msg_error(__func__, "A problem happened computing the linear term");
01244 cpl_msg_error(__func__,"Error %s",cpl_error_get_message());
01245
01246 return -1;
01247 }
01248 return 0;
01249 }
01250
01251
01256
01257 static cpl_propertylist ** hawki_cal_distortion_qc
01258 (hawki_distortion ** distortion,
01259 int * nmatched_pairs,
01260 double * rms)
01261 {
01262 int idet;
01263 cpl_propertylist ** qclists;
01264
01265
01266 qclists = cpl_malloc(HAWKI_NB_DETECTORS * sizeof(cpl_propertylist*)) ;
01267
01268
01269 for(idet = 0 ; idet < HAWKI_NB_DETECTORS ; ++idet)
01270 {
01271
01272 qclists[idet] = cpl_propertylist_new() ;
01273
01274 cpl_propertylist_append_double
01275 (qclists[idet], "ESO QC DIST NMATCHED", nmatched_pairs[idet]);
01276
01277 cpl_propertylist_append_double
01278 (qclists[idet], "ESO QC DIST TOTAL RMS", rms[idet]);
01279
01280
01281
01282
01283
01284
01285
01286
01287
01288
01289
01290
01291 }
01292
01293 return qclists;
01294 }
01295
01296
01306
01307 static int hawki_cal_distortion_save
01308 (hawki_distortion ** distortion,
01309 cpl_parameterlist * parlist,
01310 cpl_propertylist ** qclists,
01311 cpl_frameset * recipe_set)
01312 {
01313 const char * recipe_name = "hawki_cal_distortion";
01314
01315
01316 hawki_distortion_save(recipe_set,
01317 parlist,
01318 recipe_set,
01319 (const hawki_distortion **) distortion,
01320 recipe_name,
01321 NULL,
01322 (const cpl_propertylist **)qclists,
01323 "hawki_cal_distortion_x.fits",
01324 "hawki_cal_distortion_y.fits");
01325
01326
01327 return 0;
01328 }
01329
01330 static int hawki_cal_distortion_retrieve_input_param
01331 (cpl_parameterlist * parlist)
01332 {
01333 cpl_parameter * par ;
01334
01335 par = NULL ;
01336 par = cpl_parameterlist_find
01337 (parlist, "hawki.hawki_cal_distortion.sigma_det");
01338 hawki_cal_distortion_config.sigma_det = cpl_parameter_get_double(par);
01339 par = cpl_parameterlist_find
01340 (parlist, "hawki.hawki_cal_distortion.grid_points");
01341 hawki_cal_distortion_config.grid_points = cpl_parameter_get_int(par);
01342 par = cpl_parameterlist_find
01343 (parlist, "hawki.hawki_cal_distortion.borders");
01344 hawki_cal_distortion_config.borders = cpl_parameter_get_int(par);
01345 par = cpl_parameterlist_find
01346 (parlist, "hawki.hawki_cal_distortion.subtract_linear");
01347 hawki_cal_distortion_config.subtract_linear = cpl_parameter_get_bool(par);
01348
01349
01350 return 0;
01351 }