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