HAWKI Pipeline Reference Manual  1.8.12
hawki_cal_distortion.c
1 /* $Id: hawki_cal_distortion.c,v 1.17 2013/02/01 17:18:05 cgarcia Exp $
2  *
3  * This file is part of the HAWKI Pipeline
4  * Copyright (C) 2002,2003 European Southern Observatory
5  *
6  * This program is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * This program is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this program; if not, write to the Free Software
18  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19  */
20 
21 /*
22  * $Author: cgarcia $
23  * $Date: 2013/02/01 17:18:05 $
24  * $Revision: 1.17 $
25  * $Name: hawki-1_8_12 $
26  */
27 
28 #ifdef HAVE_CONFIG_H
29 #include <config.h>
30 #endif
31 
32 /*-----------------------------------------------------------------------------
33  Includes
34  -----------------------------------------------------------------------------*/
35 
36 #include <math.h>
37 #include <cpl.h>
38 #include <string.h>
39 
40 #include "irplib_utils.h"
41 
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"
54 
55 /*-----------------------------------------------------------------------------
56  Functions prototypes
57  -----------------------------------------------------------------------------*/
58 
59 #ifdef __cplusplus
60 extern "C"
61 #endif
62 int cpl_plugin_get_info(cpl_pluginlist * list);
63 
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);
69 
70 
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,
77  double sigma_det,
78  cpl_apertures *** apertures);
79 
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);
87 
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);
93 
94 static int hawki_cal_distortion_subtract_sky
95 (cpl_imagelist * distor_corrected,
96  cpl_image * master_sky);
97 
98 static hawki_distortion ** hawki_cal_distortion_compute_dist_solution
99 (cpl_apertures *** apertures,
100  int nframes,
101  cpl_bivector * offsets,
102  int grid_points,
103  int * nmatched_pairs,
104  double * rms,
105  hawki_distortion ** distortion_guess);
106 
107 static cpl_apertures * hawki_cal_distortion_get_image_apertures
108 (cpl_image * image,
109  double sigma_det);
110 
111 static int hawki_cal_distortion_fill_obj_pos
112 (cpl_table * objects_positions,
113  cpl_apertures * apertures);
114 
115 static int hawki_cal_distortion_add_offset_to_positions
116 (cpl_table ** objects_positions,
117  cpl_bivector * offsets);
118 
119 static int hawki_cal_distortion_fit_first_order_solution
120 (hawki_distortion * distortion,
121  cpl_polynomial * fit2d_x,
122  cpl_polynomial * fit2d_y);
123 
124 static cpl_propertylist ** hawki_cal_distortion_qc
125 (hawki_distortion ** distortion,
126  int * nmatched_pairs,
127  double * rms);
128 
129 static int hawki_cal_distortion_save
130 (hawki_distortion ** distortion,
131  cpl_parameterlist * parlist,
132  cpl_propertylist ** qclists,
133  cpl_frameset * recipe_set);
134 
135 static int hawki_cal_distortion_retrieve_input_param
136 (cpl_parameterlist * parlist);
137 
138 
139 /*-----------------------------------------------------------------------------
140  Static variables
141  -----------------------------------------------------------------------------*/
142 
143 static struct
144 {
145  /* Inputs */
146  double sigma_det;
147  int grid_points;
148  int borders;
149  int subtract_linear;
150 } hawki_cal_distortion_config;
151 
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"
166 "Return code:\n"
167 "esorex exits with an error code of 0 if the recipe completes successfully\n"
168 "or 1 otherwise";
169 
170 /*-----------------------------------------------------------------------------
171  Functions code
172  -----------------------------------------------------------------------------*/
173 
174 /*----------------------------------------------------------------------------*/
182 /*----------------------------------------------------------------------------*/
183 int cpl_plugin_get_info(cpl_pluginlist * list)
184 {
185  cpl_recipe * recipe = cpl_calloc(1, sizeof(*recipe)) ;
186  cpl_plugin * plugin = &recipe->interface ;
187 
188  cpl_plugin_init(plugin,
189  CPL_PLUGIN_API,
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",
196  PACKAGE_BUGREPORT,
198  hawki_cal_distortion_create,
199  hawki_cal_distortion_exec,
200  hawki_cal_distortion_destroy);
201 
202  cpl_pluginlist_append(list, plugin) ;
203 
204  return 0;
205 }
206 
207 /*----------------------------------------------------------------------------*/
216 /*----------------------------------------------------------------------------*/
217 static int hawki_cal_distortion_create(cpl_plugin * plugin)
218 {
219  cpl_recipe * recipe ;
220  cpl_parameter * p ;
221 
222  /* Get the recipe out of the plugin */
223  if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE)
224  recipe = (cpl_recipe *)plugin ;
225  else return -1 ;
226 
227  /* Create the parameters list in the cpl_recipe object */
228  recipe->parameters = cpl_parameterlist_new() ;
229  if (recipe->parameters == NULL)
230  return 1;
231 
232  /* Fill the parameters list */
233  /* --sigma_det */
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) ;
240 
241  /* --grid_points */
242  p = cpl_parameter_new_value("hawki.hawki_cal_distortion.grid_points",
243  CPL_TYPE_INT,
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) ;
249 
250  /* --borders */
251  p = cpl_parameter_new_value("hawki.hawki_cal_distortion.borders",
252  CPL_TYPE_INT,
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) ;
258 
259  /* --subtract_linear */
260  p = cpl_parameter_new_value("hawki.hawki_cal_distortion.subtract_linear",
261  CPL_TYPE_BOOL,
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) ;
267 
268  /* Return */
269  return 0;
270 }
271 
272 /*----------------------------------------------------------------------------*/
278 /*----------------------------------------------------------------------------*/
279 static int hawki_cal_distortion_exec(cpl_plugin * plugin)
280 {
281  cpl_recipe * recipe ;
282 
283  /* Get the recipe out of the plugin */
284  if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE)
285  recipe = (cpl_recipe *)plugin ;
286  else return -1 ;
287 
288  /* Issue a banner */
290 
291  return hawki_cal_distortion(recipe->parameters, recipe->frames) ;
292 }
293 
294 /*----------------------------------------------------------------------------*/
300 /*----------------------------------------------------------------------------*/
301 static int hawki_cal_distortion_destroy(cpl_plugin * plugin)
302 {
303  cpl_recipe * recipe ;
304 
305  /* Get the recipe out of the plugin */
306  if (cpl_plugin_get_type(plugin) == CPL_PLUGIN_TYPE_RECIPE)
307  recipe = (cpl_recipe *)plugin ;
308  else return -1 ;
309 
310  cpl_parameterlist_delete(recipe->parameters) ;
311  return 0 ;
312 }
313 
314 /*----------------------------------------------------------------------------*/
320 /*----------------------------------------------------------------------------*/
321 static int hawki_cal_distortion(cpl_parameterlist * parlist,
322  cpl_frameset * frameset)
323 {
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];
339  int idet;
340  int ioff;
341  int iframe;
342  int nframes;
343 
344 
345  /* Retrieve input parameters */
346  if(hawki_cal_distortion_retrieve_input_param(parlist))
347  {
348  cpl_msg_error(__func__, "Wrong parameters");
349  return -1;
350  }
351 
352  /* Identify the RAW and CALIB frames in the input frameset */
353  if (hawki_dfs_set_groups(frameset)) {
354  cpl_msg_error(__func__, "Cannot identify RAW and CALIB frames") ;
355  return -1 ;
356  }
357 
358  /* Retrieve calibration data */
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);
363 
364  /* Retrieve raw frames */
365  cpl_msg_info(__func__, "Identifying distortion and sky data");
366  distorframes = hawki_extract_frameset(frameset, HAWKI_CAL_DISTOR_RAW) ;
367  if (distorframes == NULL)
368  {
369  cpl_msg_error(__func__, "Distortion images have to be provided (%s)",
370  HAWKI_CAL_DISTOR_RAW);
371  return -1 ;
372  }
373  /* Retrieve sky frames */
374  skyframes = hawki_extract_frameset(frameset, HAWKI_CAL_DISTOR_SKY_RAW) ;
375  if (skyframes == NULL)
376  {
377  cpl_msg_error(__func__, "Sky images have to be provided (%s)",
378  HAWKI_CAL_DISTOR_SKY_RAW);
379  cpl_frameset_delete(distorframes);
380  return -1 ;
381  }
382  /* Retrieve the distortion first guess (if provided) */
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)
386  {
387  //distortionguess = hawki_distortion_load(distorxtguess)
388  }
389 
390  /* Get the master sky frame */
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)
394  {
395  cpl_msg_error(__func__, "Cannot get master sky image") ;
396  cpl_frameset_delete(distorframes);
397  cpl_frameset_delete(skyframes);
398  return -1;
399  }
400 
401  /* Aperture detection, basic reduction and sky subtraction of distortion images */
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)
406  {
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);
414  return -1 ;
415  }
416  for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
417  cpl_image_delete(master_sky[idet]);
418  cpl_free(master_sky);
419 
420  /* Get the nominal offsets from the header */
421  cpl_msg_info(__func__,"Getting the nominal offsets");
422  nominal_offsets = hawki_get_header_tel_offsets(distorframes);
423  if (nominal_offsets == NULL)
424  {
425  cpl_msg_error(__func__, "Cannot load the header offsets") ;
426  cpl_frameset_delete(distorframes);
427  cpl_frameset_delete(skyframes);
428  return -1;
429  }
430 
431  /* Get the oposite offsets. This is to change from
432  * telescope convention to cpl convention */
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);
435 
436  /* Print the header offsets */
437  cpl_msg_indent_more();
438  for (ioff=0 ; ioff<cpl_bivector_get_size(nominal_offsets) ; ioff++)
439  {
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]);
443  }
444  cpl_msg_indent_less();
445 
446  /* Get the distortion solution, the real stuff */
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,
452  nmatched_pairs, rms,
453  distortionguess);
454  cpl_bivector_delete(nominal_offsets);
455  if(distortion == NULL)
456  {
457  cpl_frameset_delete(distorframes);
458  cpl_frameset_delete(skyframes);
459  return -1;
460  }
461 
462  /* Get some QC */
463  qclists = hawki_cal_distortion_qc(distortion, nmatched_pairs, rms);
464 
465  /* Save the products */
466  cpl_msg_info(__func__,"Saving products");
467  if(hawki_cal_distortion_save(distortion,
468  parlist, qclists, frameset) == -1)
469  {
470  cpl_msg_error(__func__,"Could not save products");
471  for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
472  cpl_propertylist_delete(qclists[idet]);
473  cpl_free(qclists);
474  for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
475  hawki_distortion_delete(distortion[idet]);
476  cpl_free(distortion);
477  cpl_frameset_delete(distorframes);
478  cpl_frameset_delete(skyframes);
479  return -1;
480  }
481 
482  /* Free and return */
483  for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
484  cpl_propertylist_delete(qclists[idet]);
485  cpl_free(qclists);
486  for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
487  hawki_distortion_delete(distortion[idet]);
488  cpl_free(distortion);
489  for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
490  {
491  for(iframe = 0 ; iframe < nframes; iframe++)
492  cpl_apertures_delete(apertures[idet][iframe]);
493  cpl_free(apertures[idet]);
494  }
495  cpl_frameset_delete(distorframes);
496  cpl_frameset_delete(skyframes);
497 
498 
499  /* Return */
500  if (cpl_error_get_code())
501  {
502  cpl_msg_error(__func__,
503  "HAWK-I pipeline could not recover from previous errors");
504  return -1 ;
505  }
506  else return 0 ;
507 }
508 
509 /*----------------------------------------------------------------------------*/
520 /*----------------------------------------------------------------------------*/
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)
528 {
529  cpl_errorstate error_prevstate = cpl_errorstate_get();
530 
531  /* Initializing the pointers */
532  *flat_images = NULL;
533  *dark_images = NULL;
534  *bpm_images = NULL;
535 
536  /* Loading the calibration files */
537  cpl_msg_info(__func__, "Loading the calibration data") ;
538  if(flat != NULL)
539  {
540  *flat_images = hawki_load_frame(flat, CPL_TYPE_FLOAT);
541  if(*flat_images == NULL)
542  {
543  cpl_msg_error(__func__, "Error reading flat") ;
544  return -1;
545  }
546  }
547  if(dark != NULL)
548  {
549  *dark_images = hawki_load_frame(dark, CPL_TYPE_FLOAT);
550  if(*dark_images == NULL)
551  {
552  cpl_msg_error(__func__, "Error reading dark") ;
553  cpl_imagelist_delete(*flat_images);
554  return -1;
555  }
556  }
557  if(bpm != NULL)
558  {
559  *bpm_images = hawki_load_frame(bpm, CPL_TYPE_INT);
560  if(*bpm_images == NULL)
561  {
562  cpl_msg_error(__func__, "Error reading bpm") ;
563  cpl_imagelist_delete(*flat_images);
564  cpl_imagelist_delete(*dark_images);
565  return -1;
566  }
567  }
568 
569  if(!cpl_errorstate_is_equal(error_prevstate ))
570  {
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);
575  return -1;
576  }
577  return 0;
578 }
579 
580 /*----------------------------------------------------------------------------*/
590 /*----------------------------------------------------------------------------*/
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,
597  double sigma_det,
598  cpl_apertures *** apertures)
599 {
600  cpl_imagelist * flat_images;
601  cpl_imagelist * dark_images;
602  cpl_imagelist * bpm_images;
603  cpl_propertylist * plist;
604 
605  double science_dit;
606  int iframe;
607  int nframes;
608  int idet;
609 
610  cpl_errorstate error_prevstate = cpl_errorstate_get();
611 
612  /* Indentation */
613  cpl_msg_indent_more();
614 
615  /* Loading calibrations */
616  hawki_cal_distortion_load_master_calib
617  (flat, dark, bpm, &flat_images, &dark_images, &bpm_images);
618 
619  /* Multiply the dark image by the science exposure time */
620  if(dark != NULL)
621  {
622  if ((plist=cpl_propertylist_load
623  (cpl_frame_get_filename
624  (cpl_frameset_get_first_const(raw_distor)), 0)) == NULL)
625  {
626  cpl_msg_error(__func__, "Cannot get header from frame");
627  cpl_imagelist_delete(flat_images);
628  cpl_imagelist_delete(dark_images);
629  return -1;
630  }
631  science_dit = hawki_pfits_get_dit(plist);
632  cpl_imagelist_multiply_scalar(dark_images, science_dit);
633  cpl_propertylist_delete(plist);
634  }
635 
636  /* Loop on detectors */
637  nframes = cpl_frameset_get_size(raw_distor);
638  for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
639  {
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;
645 
646  cpl_msg_info(__func__, "Working on detector %d", idet + 1);
647  cpl_msg_indent_more();
648 
649  /* Loading the distortion images for one detector */
650  cpl_msg_info(__func__, "Loading distortion images");
651  distor_serie = hawki_load_detector(raw_distor, idet + 1, CPL_TYPE_FLOAT);
652  if(distor_serie== NULL)
653  {
654  cpl_msg_error(__func__, "Error reading distortion images");
655  return -1;
656  }
657 
658  /* Getting the calibs */
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);
665 
666  if(!cpl_errorstate_is_equal(error_prevstate ))
667  cpl_msg_warning(cpl_func,"PPPPPP5");
668 
669  /* Applying the calibrations */
670  cpl_msg_info(__func__, "Applying basic calibration") ;
671  cpl_msg_indent_more();
673  (distor_serie, flat_det, dark_det, bpm_det) == -1)
674  {
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() ;
681  return -1;
682  }
683  cpl_msg_indent_less();
684 
685  /* Discard the pixels on the sides */
686  if (hawki_cal_distortion_config.borders > 0)
687  {
688  distor_serie_trimmed = hawki_trim_detector_calib(distor_serie,
689  hawki_cal_distortion_config.borders);
690  cpl_imagelist_delete(distor_serie);
691  }
692  else
693  distor_serie_trimmed = distor_serie;
694 
695  /* Subtract sky */
696  cpl_msg_info(__func__, "Subtracting master sky") ;
697  if(hawki_cal_distortion_subtract_sky(distor_serie_trimmed, master_sky[idet]) == -1)
698  {
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);
704  return -1;
705  }
706 
707  /* Creating apertures */
708  apertures[idet] = cpl_malloc(sizeof(*(apertures[idet])) * nframes);
709  for(iframe = 0 ; iframe < nframes; iframe++)
710  {
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();
717  }
718 
719  /* Delete the list of images */
720  cpl_imagelist_delete(distor_serie_trimmed);
721  cpl_msg_indent_less();
722  }
723  cpl_imagelist_delete(flat_images);
724  cpl_imagelist_delete(dark_images);
725  cpl_imagelist_delete(bpm_images);
726 
727 
728  if(!cpl_errorstate_is_equal(error_prevstate ))
729  {
730  cpl_msg_error(__func__, "A problem happened detecting objects");
731  return -1 ;
732  }
733  return 0;
734 }
735 
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)
741 {
742  cpl_propertylist * plist;
743  double science_dit;
744  int idet;
745  int jdet;
746  cpl_imagelist * flat_images;
747  cpl_imagelist * dark_images;
748  cpl_imagelist * bpm_images;
749  cpl_image ** bkg_images = NULL;
750 
751  cpl_errorstate error_prevstate = cpl_errorstate_get();
752 
753  /* Indentation */
754  cpl_msg_indent_more();
755 
756  /* Reading calibrations */
757  hawki_cal_distortion_load_master_calib
758  (flat, dark, bpm, &flat_images, &dark_images, &bpm_images);
759 
760  /* Multiply the dark image by the science exposure time */
761  if(dark != NULL)
762  {
763  if ((plist=cpl_propertylist_load
764  (cpl_frame_get_filename
765  (cpl_frameset_get_first_const(raw_sky_frames)), 0)) == NULL)
766  {
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);
771  return NULL;
772  }
773  science_dit = hawki_pfits_get_dit(plist);
774  cpl_imagelist_multiply_scalar(dark_images, science_dit);
775  cpl_propertylist_delete(plist);
776  }
777 
778  /* Compute the sky median */
779  bkg_images = cpl_malloc(HAWKI_NB_DETECTORS * sizeof(*bkg_images));
780  for(idet = 0; idet < HAWKI_NB_DETECTORS ; ++idet)
781  {
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;
787 
788  /* Loading the sky images for one detector */
789  sky_serie = hawki_load_detector(raw_sky_frames, idet + 1, CPL_TYPE_FLOAT);
790  if(sky_serie== NULL)
791  {
792  cpl_msg_error(__func__, "Error reading object image") ;
793  return NULL;
794  }
795 
796  /* Getting the calibs */
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);
803 
804  /* Applying the calibrations */
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)
809  {
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() ;
819  return NULL;
820  }
821 
822  /* Discard the pixels on the sides */
823  if (hawki_cal_distortion_config.borders > 0)
824  {
825  sky_serie_trimmed = hawki_trim_detector_calib(sky_serie,
826  hawki_cal_distortion_config.borders);
827  cpl_imagelist_delete(sky_serie);
828  }
829  else
830  sky_serie_trimmed = sky_serie;
831 
832  /* Averaging */
833  if ((bkg_images[idet] =
834  cpl_imagelist_collapse_median_create(sky_serie_trimmed)) == NULL)
835  {
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();
845  return NULL;
846  }
847  cpl_imagelist_delete(sky_serie_trimmed);
848  cpl_msg_indent_less();
849  }
850 
851  /* Subtract the median of the frame */
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]));
855 
856  /* Cleaning up */
857  cpl_msg_indent_less();
858  cpl_imagelist_delete(flat_images);
859  cpl_imagelist_delete(dark_images);
860  cpl_imagelist_delete(bpm_images);
861 
862  if(!cpl_errorstate_is_equal(error_prevstate ))
863  {
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);
868  return NULL ;
869  }
870  return bkg_images;
871 }
872 
873 static int hawki_cal_distortion_subtract_sky
874 (cpl_imagelist * distor_corrected,
875  cpl_image * master_sky)
876 {
877  cpl_errorstate error_prevstate = cpl_errorstate_get();
878 
879  /* Subtract the background to each object frame */
880  int idist, ndistor;
881  ndistor = cpl_imagelist_get_size(distor_corrected);
882  for(idist = 0; idist < ndistor; ++idist)
883  {
884  cpl_image * target_image =
885  cpl_imagelist_get(distor_corrected, idist);
886  /* First subtract the median of the image */
887  cpl_image_subtract_scalar
888  (target_image, cpl_image_get_median(target_image));
889 
890  if (cpl_image_subtract
891  (target_image, master_sky)!=CPL_ERROR_NONE)
892  {
893  cpl_msg_error(cpl_func,"Cannot apply the bkg to the images");
894  return -1 ;
895  }
896  }
897 
898  /* Free and return */
899  if(!cpl_errorstate_is_equal(error_prevstate ))
900  {
901  cpl_msg_error(__func__, "A problem happened with sky subtraction");
902  return -1;
903  }
904  return 0;
905 }
906 
907 
908 /*----------------------------------------------------------------------------*/
918 /*----------------------------------------------------------------------------*/
919 static hawki_distortion ** hawki_cal_distortion_compute_dist_solution
920 (cpl_apertures *** apertures,
921  int nframes,
922  cpl_bivector * offsets,
923  int grid_points,
924  int * nmatched_pairs,
925  double * rms,
926  hawki_distortion ** distortion_guess)
927 {
928  int idet;
929  hawki_distortion ** distortion = NULL;
930 
931  cpl_errorstate error_prevstate = cpl_errorstate_get();
932 
933  /* Allocate the distortion */
934  distortion = cpl_malloc(HAWKI_NB_DETECTORS * sizeof(*distortion));
935 
936  /* Loop on the detectors */
937  cpl_msg_indent_more();
938  for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
939  {
940  cpl_table ** obj_pos;
941  cpl_table ** obj_pos_offset;
942  int iframe;
943  cpl_table * matches;
944  hawki_distortion * dist_guess;
945  cpl_polynomial * fit2d_x = NULL;
946  cpl_polynomial * fit2d_y = NULL;
947 
948 
949  cpl_msg_info(__func__, "Working on detector %d", idet+1);
950  cpl_msg_indent_more();
951 
952  /* Initialize the objects positions */
953  obj_pos =
954  cpl_malloc(sizeof(*obj_pos) * nframes);
955  obj_pos_offset =
956  cpl_malloc(sizeof(*obj_pos_offset) * nframes);
957  for(iframe = 0; iframe < nframes; ++iframe)
958  {
959  obj_pos[iframe] = cpl_table_new(0);
960  cpl_table_new_column
961  (obj_pos[iframe], HAWKI_COL_OBJ_POSX, CPL_TYPE_DOUBLE);
962  cpl_table_new_column
963  (obj_pos[iframe], HAWKI_COL_OBJ_POSY, CPL_TYPE_DOUBLE);
964  }
965 
966  /* Loop on images to fill object_positions */
967  for(iframe = 0 ; iframe < nframes ; ++iframe)
968  {
969  cpl_apertures * this_apertures;
970 
971  /* Create the detected apertures list */
972  this_apertures = apertures[idet][iframe];
973 
974  /* Fill the objects position table */
975  hawki_cal_distortion_fill_obj_pos(obj_pos[iframe],
976  this_apertures);
977  obj_pos_offset[iframe] = cpl_table_duplicate(obj_pos[iframe]);
978  }
979 
980  /* Get the objects positions with offsets */
981  hawki_cal_distortion_add_offset_to_positions
982  (obj_pos_offset, offsets);
983 
984  /* Get the all the matching pairs */
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);
991  if(matches == NULL)
992  {
993  cpl_msg_error(__func__, "Cannot match objects ");
994  for(iframe = 0; iframe < nframes; ++iframe)
995  cpl_table_delete(obj_pos[iframe]);
996  cpl_free(obj_pos);
997  return NULL;
998  }
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);
1002 
1003  /* Compute the distortion */
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];
1008  else
1009  dist_guess = NULL;
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)
1015  {
1016  int jdet;
1017  cpl_msg_error(__func__,"Could not get the distortion");
1018  for(iframe = 0; iframe < nframes; ++iframe)
1019  cpl_table_delete(obj_pos[iframe]);
1020  cpl_free(obj_pos);
1021  for(jdet = 0; jdet < idet; ++jdet)
1022  hawki_distortion_delete(distortion[idet]);
1023  cpl_table_delete(matches);
1024  return NULL;
1025  }
1026 
1027  /* Removing the first order polinomial to the distortion */
1028  if(hawki_cal_distortion_config.subtract_linear)
1029  {
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);
1035  }
1036 
1037  /* Free */
1038  for(iframe = 0; iframe < nframes; ++iframe)
1039  cpl_table_delete(obj_pos[iframe]);
1040  cpl_free(obj_pos);
1041  if(hawki_cal_distortion_config.subtract_linear)
1042  {
1043  cpl_polynomial_delete(fit2d_x);
1044  cpl_polynomial_delete(fit2d_y);
1045  }
1046  cpl_table_delete(matches);
1047  cpl_msg_indent_less();
1048  }
1049 
1050  if(!cpl_errorstate_is_equal(error_prevstate ))
1051  {
1052  cpl_msg_error(__func__, "A problem happened computing the distortion");
1053  for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
1054  hawki_distortion_delete(distortion[idet]);
1055  return NULL ;
1056  }
1057  /* Free and return */
1058  return distortion;
1059 }
1060 
1061 static cpl_apertures * hawki_cal_distortion_get_image_apertures
1062 (cpl_image * image,
1063  double sigma_det)
1064 {
1065  cpl_apertures * apertures = NULL;
1066  cpl_mask * kernel = NULL;
1067  cpl_mask * object_mask = NULL;
1068  cpl_image * labels = NULL;
1069  cpl_size nobj;
1070  double median;
1071  double med_dist;
1072  double threshold;
1073 
1074  /* Get the threshold */
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);
1078 
1079  /* Create the mask */
1080  object_mask = cpl_mask_threshold_image_create
1081  (image, threshold, DBL_MAX);
1082  if (object_mask == NULL)
1083  return NULL;
1084 
1085  /* Apply morphological opening to remove single pixel detections */
1086  kernel = cpl_mask_new(3,3);
1087  cpl_mask_not(kernel);
1088 
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) ;
1093  return NULL;
1094  }
1095  cpl_mask_delete(kernel);
1096 
1097  /* Labelise the different detected apertures */
1098  labels = cpl_image_labelise_mask_create(object_mask, &nobj);
1099  if (labels == NULL)
1100  {
1101  cpl_mask_delete(object_mask);
1102  return NULL;
1103  }
1104  cpl_mask_delete(object_mask);
1105  cpl_msg_info(__func__, "Number of objects detected: %"CPL_SIZE_FORMAT,
1106  nobj);
1107 
1108  /* Create the detected apertures list */
1109  apertures = cpl_apertures_new_from_image(image, labels);
1110  if (apertures == NULL)
1111  {
1112  cpl_image_delete(labels);
1113  return NULL;
1114  }
1115  cpl_image_delete(labels);
1116  return apertures;
1117 }
1118 
1119 static int hawki_cal_distortion_fill_obj_pos
1120 (cpl_table * objects_positions,
1121  cpl_apertures * apertures)
1122 {
1123  cpl_size nobjs;
1124  cpl_size iobj;
1125  double border_off = 0;
1126 
1127  /* Take into account that the images have been trimmed */
1128  if(hawki_cal_distortion_config.borders > 0)
1129  border_off = hawki_cal_distortion_config.borders;
1130 
1131  nobjs = cpl_apertures_get_size(apertures);
1132  cpl_table_set_size(objects_positions, nobjs);
1133 
1134  for (iobj=0 ; iobj<nobjs ; iobj++)
1135  {
1136  /* Fill with the already known information */
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);
1143  }
1144 
1145  return 0;
1146 }
1147 
1148 static int hawki_cal_distortion_add_offset_to_positions
1149 (cpl_table ** objects_positions,
1150  cpl_bivector * offsets)
1151 {
1152  int nframes;
1153  int iframe;
1154  cpl_size nobjs;
1155  cpl_size iobj;
1156 
1157  nframes = cpl_bivector_get_size(offsets);
1158 
1159  for(iframe = 0 ; iframe < nframes ; ++iframe)
1160  {
1161  double offset_x;
1162  double offset_y;
1163  int null;
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++)
1168  {
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);
1177  }
1178  }
1179 
1180  return 0;
1181 }
1182 
1183 static int hawki_cal_distortion_fit_first_order_solution
1184 (hawki_distortion * distortion,
1185  cpl_polynomial * fit2d_x,
1186  cpl_polynomial * fit2d_y)
1187 {
1188  cpl_matrix * pixel_pos;
1189  cpl_vector * dist_x_val;
1190  cpl_vector * dist_y_val;
1191  int nx;
1192  int ny;
1193  int i;
1194  int j;
1195  int null;
1196  const cpl_size mindeg2d[] = {0, 0};
1197  const cpl_size maxdeg2d[] = {1, 1};
1198  cpl_errorstate error_prevstate = cpl_errorstate_get();
1199  cpl_vector * pix;
1200  cpl_image * dist_x_plane;
1201  cpl_image * dist_y_plane;
1202  double dist_x_mean;
1203  double dist_y_mean;
1204 
1205  /* Fill the bivector with pixel positions in X,Y */
1206  nx = hawki_distortion_get_size_x(distortion);
1207  ny = hawki_distortion_get_size_y(distortion);
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)
1213  {
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));
1220  }
1221 
1222  /* Fit the polynomial */
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);
1227  /* Removing the constant term */
1228  cpl_polynomial_set_coeff(fit2d_x, mindeg2d, 0.);
1229  cpl_polynomial_set_coeff(fit2d_y, mindeg2d, 0.);
1230 
1231  /* Subtract the linear term */
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)
1237  {
1238  double fit_value_x;
1239  double fit_value_y;
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);
1246  }
1247  cpl_image_subtract(distortion->dist_x, dist_x_plane);
1248  cpl_image_subtract(distortion->dist_y, dist_y_plane);
1249 
1250  /* Subtract the mean distortion, again */
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);
1257 
1258  /* Free and return */
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);
1265 
1266  if(!cpl_errorstate_is_equal(error_prevstate ))
1267  {
1268  cpl_msg_error(__func__, "A problem happened computing the linear term");
1269  cpl_msg_error(__func__,"Error %s",cpl_error_get_message());
1270  //cpl_msg_error(__func__,"Where %s",cpl_error_get_where());
1271  return -1;
1272  }
1273  return 0;
1274 }
1275 
1276 /*----------------------------------------------------------------------------*/
1281 /*----------------------------------------------------------------------------*/
1282 static cpl_propertylist ** hawki_cal_distortion_qc
1283 (hawki_distortion ** distortion,
1284  int * nmatched_pairs,
1285  double * rms)
1286 {
1287  int idet;
1288  cpl_propertylist ** qclists;
1289 
1290  /* Allocate the qclists */
1291  qclists = cpl_malloc(HAWKI_NB_DETECTORS * sizeof(cpl_propertylist*)) ;
1292 
1293  /* Loop on the detectors to get the mean zpoint */
1294  for(idet = 0 ; idet < HAWKI_NB_DETECTORS ; ++idet)
1295  {
1296  /* Allocate this qclist */
1297  qclists[idet] = cpl_propertylist_new() ;
1298 
1299  cpl_propertylist_append_double
1300  (qclists[idet], "ESO QC DIST NMATCHED", nmatched_pairs[idet]);
1301 
1302  cpl_propertylist_append_double
1303  (qclists[idet], "ESO QC DIST TOTAL RMS", rms[idet]);
1304 
1305  /* Getting the jacobian of the distortion map */
1306  /* The jacobian has to be definitive positive in all the detector to
1307  * be have a biyective function invertible anywhere:
1308  * http://en.wikipedia.org/wiki/Jacobian_matrix_and_determinant#Jacobian_determinant
1309  * http://en.wikipedia.org/wiki/Inverse_function#Inverses_and_derivatives
1310  * This should be a QC check.
1311  */
1312 
1313 
1314  //cpl_propertylist_append_double
1315  //(qclists[idet], "ESO QC DIST JACOBIAN_1_1", jacobian[1][1]);
1316  }
1317 
1318  return qclists;
1319 }
1320 
1321 /*----------------------------------------------------------------------------*/
1331 /*----------------------------------------------------------------------------*/
1332 static int hawki_cal_distortion_save
1333 (hawki_distortion ** distortion,
1334  cpl_parameterlist * parlist,
1335  cpl_propertylist ** qclists,
1336  cpl_frameset * recipe_set)
1337 {
1338  const char * recipe_name = "hawki_cal_distortion";
1339 
1340  /* Write the distortion in both axes */
1341  hawki_distortion_save(recipe_set,
1342  parlist,
1343  recipe_set,
1344  (const hawki_distortion **) distortion,
1345  recipe_name,
1346  NULL,
1347  (const cpl_propertylist **)qclists,
1348  "hawki_cal_distortion_x.fits",
1349  "hawki_cal_distortion_y.fits");
1350 
1351  /* Free and return */
1352  return 0;
1353 }
1354 
1355 static int hawki_cal_distortion_retrieve_input_param
1356 (cpl_parameterlist * parlist)
1357 {
1358  cpl_parameter * par ;
1359 
1360  par = NULL ;
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);
1373 
1374 
1375  return 0;
1376 }