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