40 #include "irplib_utils.h"
42 #include "hawki_alloc.h"
43 #include "hawki_utils.h"
44 #include "hawki_calib.h"
45 #include "hawki_distortion.h"
46 #include "hawki_load.h"
47 #include "hawki_save.h"
48 #include "hawki_pfits.h"
49 #include "hawki_dfs.h"
50 #include "irplib_cat.h"
51 #include "irplib_stdstar.h"
52 #include "irplib_match_cats.h"
53 #include "hawki_match_cats.h"
62 int cpl_plugin_get_info(cpl_pluginlist * list);
64 static int hawki_cal_distortion_create(cpl_plugin *) ;
65 static int hawki_cal_distortion_exec(cpl_plugin *) ;
66 static int hawki_cal_distortion_destroy(cpl_plugin *) ;
67 static int hawki_cal_distortion(cpl_parameterlist * parlist,
68 cpl_frameset * frameset);
71 static int hawki_cal_distortion_get_apertures_from_raw_distor
72 (cpl_frameset * raw_target,
73 const cpl_frame * flat,
74 const cpl_frame * dark,
75 const cpl_frame * bpm,
76 cpl_image ** master_sky,
78 cpl_apertures *** apertures);
80 static int hawki_cal_distortion_load_master_calib
81 (
const cpl_frame * flat,
82 const cpl_frame * dark,
83 const cpl_frame * bpm,
84 cpl_imagelist ** flat_images,
85 cpl_imagelist ** dark_images,
86 cpl_imagelist ** bpm_images);
88 static cpl_image ** hawki_cal_distortion_get_master_sky
89 (cpl_frameset * raw_target,
90 const cpl_frame * flat,
91 const cpl_frame * dark,
92 const cpl_frame * bpm);
94 static int hawki_cal_distortion_subtract_sky
95 (cpl_imagelist * distor_corrected,
96 cpl_image * master_sky);
98 static hawki_distortion ** hawki_cal_distortion_compute_dist_solution
99 (cpl_apertures *** apertures,
101 cpl_bivector * offsets,
103 int * nmatched_pairs,
105 hawki_distortion ** distortion_guess);
107 static cpl_apertures * hawki_cal_distortion_get_image_apertures
111 static int hawki_cal_distortion_fill_obj_pos
112 (cpl_table * objects_positions,
113 cpl_apertures * apertures);
115 static int hawki_cal_distortion_add_offset_to_positions
116 (cpl_table ** objects_positions,
117 cpl_bivector * offsets);
119 static int hawki_cal_distortion_fit_first_order_solution
120 (hawki_distortion * distortion,
121 cpl_polynomial * fit2d_x,
122 cpl_polynomial * fit2d_y);
124 static cpl_propertylist ** hawki_cal_distortion_qc
125 (hawki_distortion ** distortion,
126 int * nmatched_pairs,
129 static int hawki_cal_distortion_save
130 (hawki_distortion ** distortion,
131 cpl_parameterlist * parlist,
132 cpl_propertylist ** qclists,
133 cpl_frameset * recipe_set);
135 static int hawki_cal_distortion_retrieve_input_param
136 (cpl_parameterlist * parlist);
150 } hawki_cal_distortion_config;
152 static char hawki_cal_distortion_description[] =
153 "hawki_cal_distortion -- HAWK-I distortion and astrometry autocalibration.\n\n"
154 "The input files must be tagged:\n"
155 "distortion_field.fits "HAWKI_CAL_DISTOR_RAW
"\n"
156 "sky_distortion.fits "HAWKI_CAL_DISTOR_SKY_RAW
"\n"
157 "flat-file.fits "HAWKI_CALPRO_FLAT
" (optional)\n"
158 "dark-file.fits "HAWKI_CALPRO_DARK
" (optional)\n"
159 "bpm-file.fits "HAWKI_CALPRO_BPM
" (optional)\n\n"
160 "The recipe creates as an output:\n"
161 "hawki_cal_distortion_distx.fits ("HAWKI_CALPRO_DISTORTION_X
") \n"
162 "hawki_cal_distortion_disty.fits ("HAWKI_CALPRO_DISTORTION_Y
") \n\n"
163 "The recipe performs the following steps:\n"
164 "-Basic calibration of astrometry fields\n"
165 "-Autocalibration of distortion, using method in A&A 454,1029 (2006)\n\n"
167 "esorex exits with an error code of 0 if the recipe completes successfully\n"
183 int cpl_plugin_get_info(cpl_pluginlist * list)
185 cpl_recipe * recipe = cpl_calloc(1,
sizeof(*recipe)) ;
186 cpl_plugin * plugin = &recipe->interface ;
188 cpl_plugin_init(plugin,
190 HAWKI_BINARY_VERSION,
191 CPL_PLUGIN_TYPE_RECIPE,
192 "hawki_cal_distortion",
193 "Distortion autocalibration",
194 hawki_cal_distortion_description,
195 "Cesar Enrique Garcia Dabo",
198 hawki_cal_distortion_create,
199 hawki_cal_distortion_exec,
200 hawki_cal_distortion_destroy);
202 cpl_pluginlist_append(list, plugin) ;
217 static int hawki_cal_distortion_create(cpl_plugin * plugin)
219 cpl_recipe * recipe ;
223 if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE)
224 recipe = (cpl_recipe *)plugin ;
228 recipe->parameters = cpl_parameterlist_new() ;
229 if (recipe->parameters == NULL)
234 p = cpl_parameter_new_value(
"hawki.hawki_cal_distortion.sigma_det",
235 CPL_TYPE_DOUBLE,
"detection level",
236 "hawki.hawki_cal_distortion", 6.) ;
237 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"sigma_det") ;
238 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV) ;
239 cpl_parameterlist_append(recipe->parameters, p) ;
242 p = cpl_parameter_new_value(
"hawki.hawki_cal_distortion.grid_points",
244 "number of points in distortion grid",
245 "hawki.hawki_cal_distortion", 9) ;
246 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"grid_points") ;
247 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV) ;
248 cpl_parameterlist_append(recipe->parameters, p) ;
251 p = cpl_parameter_new_value(
"hawki.hawki_cal_distortion.borders",
253 "number of pixels to trim at the borders",
254 "hawki.hawki_cal_distortion", 6) ;
255 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"borders") ;
256 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV) ;
257 cpl_parameterlist_append(recipe->parameters, p) ;
260 p = cpl_parameter_new_value(
"hawki.hawki_cal_distortion.subtract_linear",
262 "Subtract a linear term to the solution",
263 "hawki.hawki_cal_distortion", TRUE) ;
264 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"subtract_linear") ;
265 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV) ;
266 cpl_parameterlist_append(recipe->parameters, p) ;
279 static int hawki_cal_distortion_exec(cpl_plugin * plugin)
281 cpl_recipe * recipe ;
284 if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE)
285 recipe = (cpl_recipe *)plugin ;
291 return hawki_cal_distortion(recipe->parameters, recipe->frames) ;
301 static int hawki_cal_distortion_destroy(cpl_plugin * plugin)
303 cpl_recipe * recipe ;
306 if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE)
307 recipe = (cpl_recipe *)plugin ;
310 cpl_parameterlist_delete(recipe->parameters) ;
321 static int hawki_cal_distortion(cpl_parameterlist * parlist,
322 cpl_frameset * frameset)
324 const cpl_frame * flat = NULL;
325 const cpl_frame * dark = NULL;
326 const cpl_frame * bpm = NULL;
327 cpl_frameset * distorframes = NULL;
328 cpl_frameset * skyframes = NULL;
329 const cpl_frame * distorxguess = NULL;
330 const cpl_frame * distoryguess = NULL;
331 hawki_distortion ** distortionguess = NULL;
332 hawki_distortion ** distortion = NULL;
333 cpl_propertylist ** qclists = NULL;
334 cpl_image ** master_sky = NULL;
335 cpl_bivector * nominal_offsets = NULL;
336 cpl_apertures ** apertures[HAWKI_NB_DETECTORS];
337 int nmatched_pairs[HAWKI_NB_DETECTORS];
338 double rms[HAWKI_NB_DETECTORS];
346 if(hawki_cal_distortion_retrieve_input_param(parlist))
348 cpl_msg_error(__func__,
"Wrong parameters");
354 cpl_msg_error(__func__,
"Cannot identify RAW and CALIB frames") ;
359 cpl_msg_info(__func__,
"Identifying calibration data");
360 flat = cpl_frameset_find_const(frameset, HAWKI_CALPRO_FLAT);
361 dark = cpl_frameset_find_const(frameset, HAWKI_CALPRO_DARK);
362 bpm = cpl_frameset_find_const(frameset, HAWKI_CALPRO_BPM);
365 cpl_msg_info(__func__,
"Identifying distortion and sky data");
367 if (distorframes == NULL)
369 cpl_msg_error(__func__,
"Distortion images have to be provided (%s)",
370 HAWKI_CAL_DISTOR_RAW);
375 if (skyframes == NULL)
377 cpl_msg_error(__func__,
"Sky images have to be provided (%s)",
378 HAWKI_CAL_DISTOR_SKY_RAW);
379 cpl_frameset_delete(distorframes);
383 distorxguess = cpl_frameset_find_const(frameset, HAWKI_CALPRO_DISTORTION_X);
384 distoryguess = cpl_frameset_find_const(frameset, HAWKI_CALPRO_DISTORTION_Y);
385 if(distorxguess != NULL && distoryguess != NULL)
391 cpl_msg_info(__func__,
"Computing the master sky image");
392 master_sky = hawki_cal_distortion_get_master_sky(skyframes, flat, dark, bpm);
393 if(master_sky == NULL)
395 cpl_msg_error(__func__,
"Cannot get master sky image") ;
396 cpl_frameset_delete(distorframes);
397 cpl_frameset_delete(skyframes);
402 cpl_msg_info(__func__,
"Getting objects from distortion images");
403 if(hawki_cal_distortion_get_apertures_from_raw_distor
404 (distorframes, flat, dark, bpm, master_sky,
405 hawki_cal_distortion_config.sigma_det, apertures) == -1)
407 cpl_msg_error(__func__,
408 "Cannot get objects from distortion images");
409 for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
410 cpl_image_delete(master_sky[idet]);
411 cpl_free(master_sky);
412 cpl_frameset_delete(distorframes);
413 cpl_frameset_delete(skyframes);
416 for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
417 cpl_image_delete(master_sky[idet]);
418 cpl_free(master_sky);
421 cpl_msg_info(__func__,
"Getting the nominal offsets");
423 if (nominal_offsets == NULL)
425 cpl_msg_error(__func__,
"Cannot load the header offsets") ;
426 cpl_frameset_delete(distorframes);
427 cpl_frameset_delete(skyframes);
433 cpl_vector_multiply_scalar(cpl_bivector_get_x(nominal_offsets), -1.0);
434 cpl_vector_multiply_scalar(cpl_bivector_get_y(nominal_offsets), -1.0);
437 cpl_msg_indent_more();
438 for (ioff=0 ; ioff<cpl_bivector_get_size(nominal_offsets) ; ioff++)
440 cpl_msg_info(__func__,
"Telescope offsets (Frame %d): %g %g", ioff+1,
441 cpl_bivector_get_x_data(nominal_offsets)[ioff],
442 cpl_bivector_get_y_data(nominal_offsets)[ioff]);
444 cpl_msg_indent_less();
447 cpl_msg_info(__func__,
"Computing the distortion");
448 nframes = cpl_frameset_get_size(distorframes);
449 distortion = hawki_cal_distortion_compute_dist_solution
450 (apertures, nframes, nominal_offsets,
451 hawki_cal_distortion_config.grid_points,
454 cpl_bivector_delete(nominal_offsets);
455 if(distortion == NULL)
457 cpl_frameset_delete(distorframes);
458 cpl_frameset_delete(skyframes);
463 qclists = hawki_cal_distortion_qc(distortion, nmatched_pairs, rms);
466 cpl_msg_info(__func__,
"Saving products");
467 if(hawki_cal_distortion_save(distortion,
468 parlist, qclists, frameset) == -1)
470 cpl_msg_error(__func__,
"Could not save products");
471 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
472 cpl_propertylist_delete(qclists[idet]);
474 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
476 cpl_free(distortion);
477 cpl_frameset_delete(distorframes);
478 cpl_frameset_delete(skyframes);
483 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
484 cpl_propertylist_delete(qclists[idet]);
486 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
488 cpl_free(distortion);
489 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
491 for(iframe = 0 ; iframe < nframes; iframe++)
492 cpl_apertures_delete(apertures[idet][iframe]);
493 cpl_free(apertures[idet]);
495 cpl_frameset_delete(distorframes);
496 cpl_frameset_delete(skyframes);
500 if (cpl_error_get_code())
502 cpl_msg_error(__func__,
503 "HAWK-I pipeline could not recover from previous errors");
521 static int hawki_cal_distortion_load_master_calib
522 (
const cpl_frame * flat,
523 const cpl_frame * dark,
524 const cpl_frame * bpm,
525 cpl_imagelist ** flat_images,
526 cpl_imagelist ** dark_images,
527 cpl_imagelist ** bpm_images)
529 cpl_errorstate error_prevstate = cpl_errorstate_get();
537 cpl_msg_info(__func__,
"Loading the calibration data") ;
541 if(*flat_images == NULL)
543 cpl_msg_error(__func__,
"Error reading flat") ;
550 if(*dark_images == NULL)
552 cpl_msg_error(__func__,
"Error reading dark") ;
553 cpl_imagelist_delete(*flat_images);
560 if(*bpm_images == NULL)
562 cpl_msg_error(__func__,
"Error reading bpm") ;
563 cpl_imagelist_delete(*flat_images);
564 cpl_imagelist_delete(*dark_images);
569 if(!cpl_errorstate_is_equal(error_prevstate ))
571 cpl_msg_error(__func__,
"A problem happened loading calibration");
572 cpl_imagelist_delete(*flat_images);
573 cpl_imagelist_delete(*dark_images);
574 cpl_imagelist_delete(*bpm_images);
591 static int hawki_cal_distortion_get_apertures_from_raw_distor
592 (cpl_frameset * raw_distor,
593 const cpl_frame * flat,
594 const cpl_frame * dark,
595 const cpl_frame * bpm,
596 cpl_image ** master_sky,
598 cpl_apertures *** apertures)
600 cpl_imagelist * flat_images;
601 cpl_imagelist * dark_images;
602 cpl_imagelist * bpm_images;
603 cpl_propertylist * plist;
610 cpl_errorstate error_prevstate = cpl_errorstate_get();
613 cpl_msg_indent_more();
616 hawki_cal_distortion_load_master_calib
617 (flat, dark, bpm, &flat_images, &dark_images, &bpm_images);
622 if ((plist=cpl_propertylist_load
623 (cpl_frame_get_filename
624 (cpl_frameset_get_first_const(raw_distor)), 0)) == NULL)
626 cpl_msg_error(__func__,
"Cannot get header from frame");
627 cpl_imagelist_delete(flat_images);
628 cpl_imagelist_delete(dark_images);
632 cpl_imagelist_multiply_scalar(dark_images, science_dit);
633 cpl_propertylist_delete(plist);
637 nframes = cpl_frameset_get_size(raw_distor);
638 for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
640 cpl_imagelist * distor_serie;
641 cpl_imagelist * distor_serie_trimmed;
642 cpl_image * flat_det = NULL;
643 cpl_image * dark_det = NULL;
644 cpl_image * bpm_det = NULL;
646 cpl_msg_info(__func__,
"Working on detector %d", idet + 1);
647 cpl_msg_indent_more();
650 cpl_msg_info(__func__,
"Loading distortion images");
652 if(distor_serie== NULL)
654 cpl_msg_error(__func__,
"Error reading distortion images");
659 if(flat_images != NULL)
660 flat_det = cpl_imagelist_get(flat_images, idet);
661 if(dark_images != NULL)
662 dark_det = cpl_imagelist_get(dark_images, idet);
663 if(bpm_images != NULL)
664 bpm_det = cpl_imagelist_get(bpm_images, idet);
666 if(!cpl_errorstate_is_equal(error_prevstate ))
667 cpl_msg_warning(cpl_func,
"PPPPPP5");
670 cpl_msg_info(__func__,
"Applying basic calibration") ;
671 cpl_msg_indent_more();
673 (distor_serie, flat_det, dark_det, bpm_det) == -1)
675 cpl_msg_error(__func__,
"Cannot calibrate frame") ;
676 cpl_imagelist_delete(flat_images);
677 cpl_imagelist_delete(dark_images);
678 cpl_imagelist_delete(bpm_images);
679 cpl_imagelist_delete(distor_serie);
680 cpl_msg_indent_less() ;
683 cpl_msg_indent_less();
686 if (hawki_cal_distortion_config.borders > 0)
689 hawki_cal_distortion_config.borders);
690 cpl_imagelist_delete(distor_serie);
693 distor_serie_trimmed = distor_serie;
696 cpl_msg_info(__func__,
"Subtracting master sky") ;
697 if(hawki_cal_distortion_subtract_sky(distor_serie_trimmed, master_sky[idet]) == -1)
699 cpl_msg_error(__func__,
"Cannot subtract the sky") ;
700 cpl_imagelist_delete(flat_images);
701 cpl_imagelist_delete(dark_images);
702 cpl_imagelist_delete(bpm_images);
703 cpl_imagelist_delete(distor_serie_trimmed);
708 apertures[idet] = cpl_malloc(
sizeof(*(apertures[idet])) * nframes);
709 for(iframe = 0 ; iframe < nframes; iframe++)
711 cpl_image * this_image = cpl_imagelist_get(distor_serie_trimmed, iframe);
712 cpl_msg_info(__func__,
"Working with distortion image %d", iframe+1);
713 cpl_msg_indent_more();
714 apertures[idet][iframe] =
715 hawki_cal_distortion_get_image_apertures(this_image, sigma_det);
716 cpl_msg_indent_less();
720 cpl_imagelist_delete(distor_serie_trimmed);
721 cpl_msg_indent_less();
723 cpl_imagelist_delete(flat_images);
724 cpl_imagelist_delete(dark_images);
725 cpl_imagelist_delete(bpm_images);
728 if(!cpl_errorstate_is_equal(error_prevstate ))
730 cpl_msg_error(__func__,
"A problem happened detecting objects");
736 static cpl_image ** hawki_cal_distortion_get_master_sky
737 (cpl_frameset * raw_sky_frames,
738 const cpl_frame * flat,
739 const cpl_frame * dark,
740 const cpl_frame * bpm)
742 cpl_propertylist * plist;
746 cpl_imagelist * flat_images;
747 cpl_imagelist * dark_images;
748 cpl_imagelist * bpm_images;
749 cpl_image ** bkg_images = NULL;
751 cpl_errorstate error_prevstate = cpl_errorstate_get();
754 cpl_msg_indent_more();
757 hawki_cal_distortion_load_master_calib
758 (flat, dark, bpm, &flat_images, &dark_images, &bpm_images);
763 if ((plist=cpl_propertylist_load
764 (cpl_frame_get_filename
765 (cpl_frameset_get_first_const(raw_sky_frames)), 0)) == NULL)
767 cpl_msg_error(__func__,
"Cannot get header from frame");
768 cpl_imagelist_delete(flat_images);
769 cpl_imagelist_delete(dark_images);
770 cpl_imagelist_delete(bpm_images);
774 cpl_imagelist_multiply_scalar(dark_images, science_dit);
775 cpl_propertylist_delete(plist);
779 bkg_images = cpl_malloc(HAWKI_NB_DETECTORS *
sizeof(*bkg_images));
780 for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
782 cpl_imagelist * sky_serie;
783 cpl_imagelist * sky_serie_trimmed;
784 cpl_image * flat_det = NULL;
785 cpl_image * dark_det = NULL;
786 cpl_image * bpm_det = NULL;
792 cpl_msg_error(__func__,
"Error reading object image") ;
797 if(flat_images != NULL)
798 flat_det = cpl_imagelist_get(flat_images, idet);
799 if(dark_images != NULL)
800 dark_det = cpl_imagelist_get(dark_images, idet);
801 if(bpm_images != NULL)
802 bpm_det = cpl_imagelist_get(bpm_images, idet);
805 cpl_msg_info(__func__,
"Working on detector %d", idet + 1);
806 cpl_msg_indent_more();
808 (sky_serie, flat_det, dark_det, bpm_det) == -1)
810 cpl_msg_error(__func__,
"Cannot calibrate frame") ;
811 cpl_imagelist_delete(flat_images);
812 cpl_imagelist_delete(dark_images);
813 cpl_imagelist_delete(bpm_images);
814 cpl_imagelist_delete(sky_serie);
815 for(jdet = 0; jdet < idet; ++jdet)
816 cpl_image_delete(bkg_images[jdet]);
817 cpl_free(bkg_images);
818 cpl_msg_indent_less() ;
823 if (hawki_cal_distortion_config.borders > 0)
826 hawki_cal_distortion_config.borders);
827 cpl_imagelist_delete(sky_serie);
830 sky_serie_trimmed = sky_serie;
833 if ((bkg_images[idet] =
834 cpl_imagelist_collapse_median_create(sky_serie_trimmed)) == NULL)
836 cpl_msg_error(__func__,
"Cannot compute the median of obj images");
837 cpl_imagelist_delete(flat_images);
838 cpl_imagelist_delete(dark_images);
839 cpl_imagelist_delete(bpm_images);
840 cpl_imagelist_delete(sky_serie_trimmed);
841 for(jdet = 0; jdet < idet; ++jdet)
842 cpl_image_delete(bkg_images[jdet]);
843 cpl_free(bkg_images);
844 cpl_msg_indent_less();
847 cpl_imagelist_delete(sky_serie_trimmed);
848 cpl_msg_indent_less();
852 for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
853 cpl_image_subtract_scalar(bkg_images[idet],
854 cpl_image_get_median(bkg_images[idet]));
857 cpl_msg_indent_less();
858 cpl_imagelist_delete(flat_images);
859 cpl_imagelist_delete(dark_images);
860 cpl_imagelist_delete(bpm_images);
862 if(!cpl_errorstate_is_equal(error_prevstate ))
864 cpl_msg_error(__func__,
"A problem happened with basic calibration");
865 for(idet = 0; idet < HAWKI_NB_DETECTORS; ++idet)
866 cpl_image_delete(bkg_images[idet]);
867 cpl_free(bkg_images);
873 static int hawki_cal_distortion_subtract_sky
874 (cpl_imagelist * distor_corrected,
875 cpl_image * master_sky)
877 cpl_errorstate error_prevstate = cpl_errorstate_get();
881 ndistor = cpl_imagelist_get_size(distor_corrected);
882 for(idist = 0; idist < ndistor; ++idist)
884 cpl_image * target_image =
885 cpl_imagelist_get(distor_corrected, idist);
887 cpl_image_subtract_scalar
888 (target_image, cpl_image_get_median(target_image));
890 if (cpl_image_subtract
891 (target_image, master_sky)!=CPL_ERROR_NONE)
893 cpl_msg_error(cpl_func,
"Cannot apply the bkg to the images");
899 if(!cpl_errorstate_is_equal(error_prevstate ))
901 cpl_msg_error(__func__,
"A problem happened with sky subtraction");
919 static hawki_distortion ** hawki_cal_distortion_compute_dist_solution
920 (cpl_apertures *** apertures,
922 cpl_bivector * offsets,
924 int * nmatched_pairs,
926 hawki_distortion ** distortion_guess)
929 hawki_distortion ** distortion = NULL;
931 cpl_errorstate error_prevstate = cpl_errorstate_get();
934 distortion = cpl_malloc(HAWKI_NB_DETECTORS *
sizeof(*distortion));
937 cpl_msg_indent_more();
938 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
940 cpl_table ** obj_pos;
941 cpl_table ** obj_pos_offset;
944 hawki_distortion * dist_guess;
945 cpl_polynomial * fit2d_x = NULL;
946 cpl_polynomial * fit2d_y = NULL;
949 cpl_msg_info(__func__,
"Working on detector %d", idet+1);
950 cpl_msg_indent_more();
954 cpl_malloc(
sizeof(*obj_pos) * nframes);
956 cpl_malloc(
sizeof(*obj_pos_offset) * nframes);
957 for(iframe = 0; iframe < nframes; ++iframe)
959 obj_pos[iframe] = cpl_table_new(0);
961 (obj_pos[iframe], HAWKI_COL_OBJ_POSX, CPL_TYPE_DOUBLE);
963 (obj_pos[iframe], HAWKI_COL_OBJ_POSY, CPL_TYPE_DOUBLE);
967 for(iframe = 0 ; iframe < nframes ; ++iframe)
969 cpl_apertures * this_apertures;
972 this_apertures = apertures[idet][iframe];
975 hawki_cal_distortion_fill_obj_pos(obj_pos[iframe],
977 obj_pos_offset[iframe] = cpl_table_duplicate(obj_pos[iframe]);
981 hawki_cal_distortion_add_offset_to_positions
982 (obj_pos_offset, offsets);
985 cpl_msg_info(__func__,
"Matching all catalogs (may take a while)");
986 matches = irplib_match_cat_pairs(obj_pos_offset, nframes,
988 for(iframe = 0; iframe < nframes; ++iframe)
989 cpl_table_delete(obj_pos_offset[iframe]);
990 cpl_free(obj_pos_offset);
993 cpl_msg_error(__func__,
"Cannot match objects ");
994 for(iframe = 0; iframe < nframes; ++iframe)
995 cpl_table_delete(obj_pos[iframe]);
999 cpl_msg_info(__func__,
"Number of matching pairs %"CPL_SIZE_FORMAT,
1000 cpl_table_get_nrow(matches));
1001 nmatched_pairs[idet] = cpl_table_get_nrow(matches);
1004 cpl_msg_info(__func__,
"Computing distortion with the matched objects");
1005 cpl_msg_info(__func__,
" (This step will take a long time to run)");
1006 if(distortion_guess != NULL)
1007 dist_guess = distortion_guess[idet];
1010 distortion[idet] = hawki_distortion_compute_solution
1011 ((
const cpl_table **)obj_pos, offsets, matches,
1012 nframes, HAWKI_DET_NPIX_X , HAWKI_DET_NPIX_Y, grid_points,
1013 dist_guess, rms + idet);
1014 if(distortion[idet] == NULL)
1017 cpl_msg_error(__func__,
"Could not get the distortion");
1018 for(iframe = 0; iframe < nframes; ++iframe)
1019 cpl_table_delete(obj_pos[iframe]);
1021 for(jdet = 0; jdet < idet; ++jdet)
1023 cpl_table_delete(matches);
1028 if(hawki_cal_distortion_config.subtract_linear)
1030 cpl_msg_info(__func__,
"Subtracting first order polynomial");
1031 fit2d_x = cpl_polynomial_new(2);
1032 fit2d_y = cpl_polynomial_new(2);
1033 hawki_cal_distortion_fit_first_order_solution
1034 (distortion[idet], fit2d_x, fit2d_y);
1038 for(iframe = 0; iframe < nframes; ++iframe)
1039 cpl_table_delete(obj_pos[iframe]);
1041 if(hawki_cal_distortion_config.subtract_linear)
1043 cpl_polynomial_delete(fit2d_x);
1044 cpl_polynomial_delete(fit2d_y);
1046 cpl_table_delete(matches);
1047 cpl_msg_indent_less();
1050 if(!cpl_errorstate_is_equal(error_prevstate ))
1052 cpl_msg_error(__func__,
"A problem happened computing the distortion");
1053 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
1061 static cpl_apertures * hawki_cal_distortion_get_image_apertures
1065 cpl_apertures * apertures = NULL;
1066 cpl_mask * kernel = NULL;
1067 cpl_mask * object_mask = NULL;
1068 cpl_image * labels = NULL;
1075 median = cpl_image_get_median_dev(image, &med_dist);
1076 threshold = median + sigma_det * med_dist;
1077 cpl_msg_info(__func__,
"Detection threshold: %f", threshold);
1080 object_mask = cpl_mask_threshold_image_create
1081 (image, threshold, DBL_MAX);
1082 if (object_mask == NULL)
1086 kernel = cpl_mask_new(3,3);
1087 cpl_mask_not(kernel);
1089 if (cpl_mask_filter(object_mask, object_mask, kernel,
1090 CPL_FILTER_OPENING, CPL_BORDER_ZERO) != CPL_ERROR_NONE) {
1091 cpl_mask_delete(object_mask) ;
1092 cpl_mask_delete(kernel) ;
1095 cpl_mask_delete(kernel);
1098 labels = cpl_image_labelise_mask_create(object_mask, &nobj);
1101 cpl_mask_delete(object_mask);
1104 cpl_mask_delete(object_mask);
1105 cpl_msg_info(__func__,
"Number of objects detected: %"CPL_SIZE_FORMAT,
1109 apertures = cpl_apertures_new_from_image(image, labels);
1110 if (apertures == NULL)
1112 cpl_image_delete(labels);
1115 cpl_image_delete(labels);
1119 static int hawki_cal_distortion_fill_obj_pos
1120 (cpl_table * objects_positions,
1121 cpl_apertures * apertures)
1125 double border_off = 0;
1128 if(hawki_cal_distortion_config.borders > 0)
1129 border_off = hawki_cal_distortion_config.borders;
1131 nobjs = cpl_apertures_get_size(apertures);
1132 cpl_table_set_size(objects_positions, nobjs);
1134 for (iobj=0 ; iobj<nobjs ; iobj++)
1137 cpl_table_set_double(objects_positions, HAWKI_COL_OBJ_POSX, iobj,
1138 cpl_apertures_get_centroid_x(apertures,
1139 iobj+1) + border_off);
1140 cpl_table_set_double(objects_positions, HAWKI_COL_OBJ_POSY, iobj,
1141 cpl_apertures_get_centroid_y(apertures,
1142 iobj+1) + border_off);
1148 static int hawki_cal_distortion_add_offset_to_positions
1149 (cpl_table ** objects_positions,
1150 cpl_bivector * offsets)
1157 nframes = cpl_bivector_get_size(offsets);
1159 for(iframe = 0 ; iframe < nframes ; ++iframe)
1164 offset_x = cpl_bivector_get_x_data(offsets)[iframe];
1165 offset_y = cpl_bivector_get_y_data(offsets)[iframe];
1166 nobjs = cpl_table_get_nrow(objects_positions[iframe]);
1167 for (iobj=0 ; iobj<nobjs ; iobj++)
1169 cpl_table_set_double(objects_positions[iframe],
1170 HAWKI_COL_OBJ_POSX, iobj,
1171 cpl_table_get_double(objects_positions[iframe],
1172 HAWKI_COL_OBJ_POSX, iobj, &null) + offset_x);
1173 cpl_table_set_double(objects_positions[iframe],
1174 HAWKI_COL_OBJ_POSY, iobj,
1175 cpl_table_get_double(objects_positions[iframe],
1176 HAWKI_COL_OBJ_POSY, iobj, &null) + offset_y);
1183 static int hawki_cal_distortion_fit_first_order_solution
1184 (hawki_distortion * distortion,
1185 cpl_polynomial * fit2d_x,
1186 cpl_polynomial * fit2d_y)
1188 cpl_matrix * pixel_pos;
1189 cpl_vector * dist_x_val;
1190 cpl_vector * dist_y_val;
1196 const cpl_size mindeg2d[] = {0, 0};
1197 const cpl_size maxdeg2d[] = {1, 1};
1198 cpl_errorstate error_prevstate = cpl_errorstate_get();
1200 cpl_image * dist_x_plane;
1201 cpl_image * dist_y_plane;
1208 pixel_pos = cpl_matrix_new(2, nx * ny);
1209 dist_x_val = cpl_vector_new(nx*ny);
1210 dist_y_val = cpl_vector_new(nx*ny);
1211 for(i = 0; i < nx; ++i)
1212 for(j = 0; j < ny; ++j)
1214 cpl_matrix_set(pixel_pos, 0, i + nx * j, (
double)i);
1215 cpl_matrix_set(pixel_pos, 1, i + nx * j, (
double)j);
1216 cpl_vector_set(dist_x_val, i + nx * j,
1217 cpl_image_get(distortion->dist_x, i+1, j+1, &null));
1218 cpl_vector_set(dist_y_val, i + nx * j,
1219 cpl_image_get(distortion->dist_y, i+1, j+1, &null));
1223 cpl_polynomial_fit(fit2d_x, pixel_pos, NULL, dist_x_val,
1224 NULL, CPL_FALSE, mindeg2d, maxdeg2d);
1225 cpl_polynomial_fit(fit2d_y, pixel_pos, NULL, dist_y_val,
1226 NULL, CPL_FALSE, mindeg2d, maxdeg2d);
1228 cpl_polynomial_set_coeff(fit2d_x, mindeg2d, 0.);
1229 cpl_polynomial_set_coeff(fit2d_y, mindeg2d, 0.);
1232 pix = cpl_vector_new(2);
1233 dist_x_plane = cpl_image_new(nx,ny,cpl_image_get_type(distortion->dist_x));
1234 dist_y_plane = cpl_image_new(nx,ny,cpl_image_get_type(distortion->dist_y));
1235 for(i = 0; i < nx; ++i)
1236 for(j = 0; j < ny; ++j)
1240 cpl_vector_set(pix, 0, (
double)i);
1241 cpl_vector_set(pix, 1, (
double)j);
1242 fit_value_x = cpl_polynomial_eval(fit2d_x, pix);
1243 fit_value_y = cpl_polynomial_eval(fit2d_y, pix);
1244 cpl_image_set(dist_x_plane, i+1, j+1, fit_value_x);
1245 cpl_image_set(dist_y_plane, i+1, j+1, fit_value_y);
1247 cpl_image_subtract(distortion->dist_x, dist_x_plane);
1248 cpl_image_subtract(distortion->dist_y, dist_y_plane);
1251 dist_x_mean = cpl_image_get_mean(distortion->dist_x);
1252 dist_y_mean = cpl_image_get_mean(distortion->dist_y);
1253 cpl_msg_warning(__func__,
"Subtracting mean distortion in X %f",dist_x_mean);
1254 cpl_msg_warning(__func__,
"Subtracting mean distortion in Y %f",dist_y_mean);
1255 cpl_image_subtract_scalar(distortion->dist_x, dist_x_mean);
1256 cpl_image_subtract_scalar(distortion->dist_y, dist_y_mean);
1259 cpl_matrix_delete(pixel_pos);
1260 cpl_vector_delete(dist_x_val);
1261 cpl_vector_delete(dist_y_val);
1262 cpl_vector_delete(pix);
1263 cpl_image_delete(dist_x_plane);
1264 cpl_image_delete(dist_y_plane);
1266 if(!cpl_errorstate_is_equal(error_prevstate ))
1268 cpl_msg_error(__func__,
"A problem happened computing the linear term");
1269 cpl_msg_error(__func__,
"Error %s",cpl_error_get_message());
1282 static cpl_propertylist ** hawki_cal_distortion_qc
1283 (hawki_distortion ** distortion,
1284 int * nmatched_pairs,
1288 cpl_propertylist ** qclists;
1291 qclists = cpl_malloc(HAWKI_NB_DETECTORS *
sizeof(cpl_propertylist*)) ;
1294 for(idet = 0 ; idet < HAWKI_NB_DETECTORS ; ++idet)
1297 qclists[idet] = cpl_propertylist_new() ;
1299 cpl_propertylist_append_double
1300 (qclists[idet],
"ESO QC DIST NMATCHED", nmatched_pairs[idet]);
1302 cpl_propertylist_append_double
1303 (qclists[idet],
"ESO QC DIST TOTAL RMS", rms[idet]);
1332 static int hawki_cal_distortion_save
1333 (hawki_distortion ** distortion,
1334 cpl_parameterlist * parlist,
1335 cpl_propertylist ** qclists,
1336 cpl_frameset * recipe_set)
1338 const char * recipe_name =
"hawki_cal_distortion";
1344 (
const hawki_distortion **) distortion,
1347 (
const cpl_propertylist **)qclists,
1348 "hawki_cal_distortion_x.fits",
1349 "hawki_cal_distortion_y.fits");
1355 static int hawki_cal_distortion_retrieve_input_param
1356 (cpl_parameterlist * parlist)
1358 cpl_parameter * par ;
1361 par = cpl_parameterlist_find
1362 (parlist,
"hawki.hawki_cal_distortion.sigma_det");
1363 hawki_cal_distortion_config.sigma_det = cpl_parameter_get_double(par);
1364 par = cpl_parameterlist_find
1365 (parlist,
"hawki.hawki_cal_distortion.grid_points");
1366 hawki_cal_distortion_config.grid_points = cpl_parameter_get_int(par);
1367 par = cpl_parameterlist_find
1368 (parlist,
"hawki.hawki_cal_distortion.borders");
1369 hawki_cal_distortion_config.borders = cpl_parameter_get_int(par);
1370 par = cpl_parameterlist_find
1371 (parlist,
"hawki.hawki_cal_distortion.subtract_linear");
1372 hawki_cal_distortion_config.subtract_linear = cpl_parameter_get_bool(par);