HAWKI Pipeline Reference Manual 1.8.12
hawki_distortion.c
00001 /* $Id: hawki_distortion.c,v 1.32 2011/02/23 11:49:37 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/02/23 11:49:37 $
00024  * $Revision: 1.32 $
00025  * $Name: hawki-1_8_12 $
00026  */
00027 
00028 #ifdef HAVE_CONFIG_H
00029 #include <config.h>
00030 #endif
00031 
00032 //Minimization algorithm hard-coded constants
00033 #define HAWKI_DISTORTION_MAX_ITER 10000
00034 #define HAWKI_DISTORTION_TOLERANCE 0.001
00035 #define HAWKI_DISTORTION_MAX_ITER2 100000
00036 #define HAWKI_DISTORTION_TOLERANCE2 0.0001
00037 
00038 /*-----------------------------------------------------------------------------
00039                                    Includes
00040  -----------------------------------------------------------------------------*/
00041 
00042 #include <math.h>
00043 #include <cpl.h>
00044 #include <cxdeque.h>
00045 #ifdef HAVE_LIBGSL
00046 #include <gsl/gsl_multimin.h>
00047 #endif
00048 
00049 #include "hawki_distortion.h"
00050 #include "hawki_dfs.h"
00051 #include "hawki_utils.h"
00052 #include "hawki_load.h"
00053 
00054 
00055 
00056 /*----------------------------------------------------------------------------*/
00060 /*----------------------------------------------------------------------------*/
00061 
00069 struct _hawki_distortion_obj_function_args_
00070 {
00071     const cpl_table   ** ref_catalogues;
00072     const cpl_table    * matching_sets;
00073     cpl_bivector       * offsets;
00074     hawki_distortion   * distortion;
00075     int                  ncats;
00076 };
00077     
00078 
00079 //Private functions
00080 
00081 int hawki_distortion_interpolate_distortion
00082 (const hawki_distortion * distortion, 
00083  double                   x_pos,
00084  double                   y_pos,
00085  double                 * x_dist, 
00086  double                 * y_dist);
00087 
00088 double hawki_distortion_compute_rms
00089 (const cpl_table    ** ref_catalogues,
00090  const cpl_bivector  * cat_offsets,
00091  const cpl_table     * matching_sets,
00092  int                   ncats,
00093  hawki_distortion    * distortion);
00094 
00095 #ifdef HAVE_LIBGSL
00096 double hawki_distortion_gsl_obj_function
00097 (const gsl_vector * dist_param,
00098  void             * args);
00099 
00100 int hawki_distortion_update_solution_from_param
00101 (hawki_distortion * distortion,
00102  const gsl_vector * dist_param);
00103 
00104 int hawki_distortion_update_offsets_from_param
00105 (cpl_bivector      * offsets,
00106  const gsl_vector  * dist_param);
00107 
00108 int hawki_distortion_update_param_from_solution
00109 (gsl_vector             * dist_param,
00110  const hawki_distortion * distortion);
00111 
00112 int hawki_distortion_update_param_from_offsets
00113 (gsl_vector              * dist_param,
00114  const cpl_bivector      * offsets);
00115 #endif
00116 
00117 void hawki_distortion_get_flag_vars
00118 (double * x_val, double * y_val, int *pos_flag,
00119  int nvals, int *nvalid, double *var_x, double *var_y);
00120 
00121 
00122 /*----------------------------------------------------------------------------*/
00130 /*----------------------------------------------------------------------------*/
00131 hawki_distortion * hawki_distortion_grid_new
00132 (int detector_nx, 
00133  int detector_ny, 
00134  int grid_size)
00135 {
00136     hawki_distortion * distortion;
00137     
00138     //Allocate the structure
00139     distortion = cpl_malloc(sizeof(hawki_distortion));
00140     
00141     //Allocate the images
00142     distortion->dist_x = cpl_image_new
00143         (grid_size, grid_size, CPL_TYPE_FLOAT);
00144     distortion->dist_y = cpl_image_new
00145         (grid_size, grid_size, CPL_TYPE_FLOAT);
00146     
00147     //Create the transformation between distortion images and the detector
00148     distortion->x_cdelt = detector_nx / (double)grid_size; 
00149     distortion->y_cdelt = detector_ny / (double)grid_size;
00150     distortion->x_crval  = 0.5 + 0.5 * distortion->x_cdelt;
00151     distortion->y_crval  = 0.5 + 0.5 * distortion->y_cdelt;
00152     
00153     return distortion;
00154 }
00155 
00156 /*----------------------------------------------------------------------------*/
00161 /*----------------------------------------------------------------------------*/
00162 void hawki_distortion_delete
00163 (hawki_distortion * distortion)
00164 {
00165     if(distortion == NULL)
00166         return;
00167     cpl_image_delete(distortion->dist_x);
00168     cpl_image_delete(distortion->dist_y);
00169     cpl_free(distortion);
00170 }
00171 
00172 /*----------------------------------------------------------------------------*/
00180 /*----------------------------------------------------------------------------*/
00181 hawki_distortion * hawki_distortion_load
00182 (const cpl_frame * dist_x,
00183  const cpl_frame * dist_y,
00184  int               idet)
00185 {
00186     const char       * file_dist_x;
00187     const char       * file_dist_y;
00188     hawki_distortion * distortion;
00189     int                iext;
00190     cpl_propertylist * plist;
00191     
00192     //Allocate the structure
00193     distortion = cpl_malloc(sizeof(hawki_distortion));
00194     
00195     //Read the images
00196     file_dist_x = cpl_frame_get_filename(dist_x);
00197     file_dist_y = cpl_frame_get_filename(dist_y);
00198     distortion->dist_x = hawki_load_frame_detector
00199         (dist_x, idet, CPL_TYPE_FLOAT);
00200     distortion->dist_y = hawki_load_frame_detector
00201         (dist_y, idet, CPL_TYPE_FLOAT);
00202     
00203     //Read the WCS keywords
00204     iext = hawki_get_ext_from_detector(file_dist_x, idet);
00205     plist = cpl_propertylist_load(file_dist_x, iext);
00206     distortion->x_crval = cpl_propertylist_get_double(plist, "CRVAL1");
00207     distortion->x_cdelt = cpl_propertylist_get_double(plist, "CDELT1");
00208     distortion->y_crval = cpl_propertylist_get_double(plist, "CRVAL2");
00209     distortion->y_cdelt = cpl_propertylist_get_double(plist, "CDELT2");
00210     if(cpl_propertylist_get_double(plist, "CRPIX1") != 1 ||
00211        cpl_propertylist_get_double(plist, "CRPIX2") != 1)
00212     {
00213         cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
00214                                     __FILE__, __LINE__,"Wrong CRPIX? keywords");
00215         cpl_image_delete(distortion->dist_x);
00216         cpl_image_delete(distortion->dist_y);
00217         cpl_propertylist_delete(plist);
00218         cpl_free(distortion);
00219         return NULL;
00220     }
00221     cpl_propertylist_delete(plist);
00222     //Check that the keywords in X and Y are compatibles;
00223     plist = cpl_propertylist_load(file_dist_y, iext);
00224     if(distortion->x_crval != cpl_propertylist_get_double(plist, "CRVAL1") ||
00225        distortion->x_cdelt != cpl_propertylist_get_double(plist, "CDELT1") ||
00226        distortion->y_crval != cpl_propertylist_get_double(plist, "CRVAL2") ||
00227        distortion->y_cdelt != cpl_propertylist_get_double(plist, "CDELT2") ||
00228        cpl_propertylist_get_double(plist, "CRPIX1") != 1 ||
00229        cpl_propertylist_get_double(plist, "CRPIX2") != 1)
00230     {
00231         cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT, __FILE__,
00232                 __LINE__,"WCS keywords mismatch in X and Y distortions");
00233         cpl_image_delete(distortion->dist_x);
00234         cpl_image_delete(distortion->dist_y);
00235         cpl_propertylist_delete(plist);
00236         cpl_free(distortion);
00237         return NULL;
00238     }
00239     cpl_propertylist_delete(plist);
00240     
00241     return distortion;
00242 }
00243 
00244 /*----------------------------------------------------------------------------*/
00250 /*----------------------------------------------------------------------------*/
00251 int hawki_distortion_get_size_x
00252 (const hawki_distortion * distortion)
00253 {
00254     if(distortion == NULL)
00255     {
00256         cpl_error_set(__func__,CPL_ERROR_ILLEGAL_INPUT);
00257     }
00258     return cpl_image_get_size_x(distortion->dist_x);
00259 }
00260 
00261 /*----------------------------------------------------------------------------*/
00267 /*----------------------------------------------------------------------------*/
00268 int hawki_distortion_get_size_y
00269 (const hawki_distortion * distortion)
00270 {
00271     if(distortion == NULL)
00272     {
00273         cpl_error_set(__func__,CPL_ERROR_ILLEGAL_INPUT);
00274     }
00275     return cpl_image_get_size_y(distortion->dist_x);
00276 }
00277 
00278 /*----------------------------------------------------------------------------*/
00285 /*----------------------------------------------------------------------------*/
00286 int hawki_distortion_correct_alldetectors
00287 (cpl_image        ** ilist,
00288  const cpl_frame   * frame_dist_x,
00289  const cpl_frame   * frame_dist_y)
00290 {
00291     cpl_image        *   corr[HAWKI_NB_DETECTORS];
00292     hawki_distortion *   distortion;
00293     cpl_image        *   dist_x;
00294     cpl_image        *   dist_y;
00295     int                  idet, j ;
00296 
00297     /* Test entries */
00298     if (ilist == NULL) return -1 ;
00299     if (frame_dist_x == NULL) return -1 ;
00300     if (frame_dist_y == NULL) return -1 ;
00301 
00302     /* Loop on the 4 chips */
00303     for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++) 
00304     {
00305         int nx;
00306         int ny;
00307         
00308         /* Get the image size */
00309         nx = cpl_image_get_size_x(ilist[idet]);
00310         ny = cpl_image_get_size_y(ilist[idet]);
00311 
00312         /* Load the distortion */
00313         if ((distortion = hawki_distortion_load
00314                  (frame_dist_x, frame_dist_y, idet + 1)) == NULL) 
00315         {
00316             cpl_msg_error(__func__, "Cannot load the distortion for chip %d", 
00317                     idet+1) ;
00318             for (j=0 ; j<idet ; j++) cpl_image_delete(corr[j]) ;
00319             return -1 ;
00320         }
00321 
00322         /* Create the offsets images */
00323         dist_x = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
00324         dist_y = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
00325         if (hawki_distortion_create_maps_detector
00326                 (distortion, dist_x, dist_y))
00327         {
00328             cpl_msg_error(__func__, "Cannot create the distortion maps") ;
00329             cpl_image_delete(dist_x);
00330             cpl_image_delete(dist_y);
00331             for (j=0 ; j<idet ; j++) cpl_image_delete(corr[j]) ;
00332             return -1;
00333         }
00334 
00335         /* Correct this image */
00336         corr[idet] = hawki_distortion_correct_detector(ilist[idet], dist_x, dist_y);
00337         if(corr[idet] == NULL)
00338         {
00339             cpl_msg_error(__func__, "Cannot correct the distortion") ;
00340             hawki_distortion_delete(distortion);
00341             cpl_image_delete(dist_x);
00342             cpl_image_delete(dist_y);
00343             for (j=0 ; j<idet; j++) cpl_image_delete(corr[j]) ;
00344             return -1 ;
00345         }
00346         hawki_distortion_delete(distortion);
00347         cpl_image_delete(dist_x) ;
00348         cpl_image_delete(dist_y);
00349     }
00350 
00351     /* Store the results */
00352     for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
00353     {
00354         cpl_image_delete(ilist[idet]) ;
00355         ilist[idet] = corr[idet] ;
00356     }
00357     
00358     /* Return */
00359     return 0 ;
00360 }
00361 
00362 /*----------------------------------------------------------------------------*/
00369 /*----------------------------------------------------------------------------*/
00370 cpl_image *  hawki_distortion_correct_detector
00371 (cpl_image       *  image,
00372  cpl_image       *  dist_x,
00373  cpl_image       *  dist_y)
00374 {
00375     cpl_image       *   corr;
00376     cpl_vector      *   profile ;
00377     
00378     /* Test entries */
00379     if (image  == NULL) return NULL;
00380     if (dist_x == NULL) return NULL;
00381     if (dist_y == NULL) return NULL;
00382 
00383     /* Create the output image */
00384     corr = cpl_image_new(cpl_image_get_size_x(image),
00385                          cpl_image_get_size_y(image), CPL_TYPE_FLOAT) ;
00386 
00387     /* Create the interpolation profile */
00388     profile = cpl_vector_new(CPL_KERNEL_DEF_SAMPLES) ;
00389     cpl_vector_fill_kernel_profile(profile, CPL_KERNEL_DEFAULT,
00390                                    CPL_KERNEL_DEF_WIDTH) ;
00391 
00392     /* Apply the distortion */
00393     if (cpl_image_warp(corr, image, dist_x, dist_y, profile, 
00394                        CPL_KERNEL_DEF_WIDTH, profile, 
00395                        CPL_KERNEL_DEF_WIDTH) != CPL_ERROR_NONE) 
00396     {
00397         cpl_msg_error(__func__, "Cannot warp the image") ;
00398         cpl_image_delete(corr) ;
00399         cpl_vector_delete(profile) ;
00400         return NULL;
00401      }
00402     cpl_vector_delete(profile) ;
00403 
00404     /* Return */
00405     return corr;
00406 }
00407 
00408 /*----------------------------------------------------------------------------*/
00421 /*----------------------------------------------------------------------------*/
00422 int hawki_distortion_correct_coords
00423 (const hawki_distortion * distortion, 
00424  double                   x_pos,
00425  double                   y_pos,
00426  double                 * x_pos_distcorr, 
00427  double                 * y_pos_distcorr)
00428 {
00429     double x_dist;
00430     double y_dist;
00431     
00432     if(distortion == NULL)
00433     {
00434         cpl_error_set("hawki_distortion_correct_coords", CPL_ERROR_ILLEGAL_INPUT);
00435         return -1;
00436     }
00437     
00438     hawki_distortion_interpolate_distortion
00439         (distortion, x_pos, y_pos, &x_dist, &y_dist);
00440 
00441     *x_pos_distcorr = x_pos - x_dist;
00442     *y_pos_distcorr = y_pos - y_dist;
00443     
00444     return 0;
00445 }
00446 
00447     
00448 /*----------------------------------------------------------------------------*/
00465 /*----------------------------------------------------------------------------*/
00466 int hawki_distortion_inverse_correct_coords
00467 (const hawki_distortion * distortion, 
00468  double                   x_pos,
00469  double                   y_pos,
00470  double                 * x_pos_distinvcorr, 
00471  double                 * y_pos_distinvcorr)
00472 {
00473     double x_dist = 0;
00474     double y_dist = 0;
00475     int    i;
00476     int    niter = 3;
00477     
00478     if(distortion == NULL)
00479     {
00480         cpl_error_set("hawki_distortion_inverse_correct_coords", CPL_ERROR_ILLEGAL_INPUT);
00481         return -1;
00482     }
00483     for(i = 0; i < niter; ++i)
00484     {
00485         hawki_distortion_interpolate_distortion
00486             (distortion, x_pos + x_dist, y_pos + y_dist, &x_dist, &y_dist);
00487     }
00488 
00489     
00490     /* Apply the correction in the inverse direction */
00491     *x_pos_distinvcorr = x_pos + x_dist;
00492     *y_pos_distinvcorr = y_pos + y_dist;
00493     
00494     return 0;
00495 }
00496 
00497 /*----------------------------------------------------------------------------*/
00524 /*----------------------------------------------------------------------------*/
00525 int hawki_distortion_interpolate_distortion
00526 (const hawki_distortion * distortion, 
00527  double                   x_pos,
00528  double                   y_pos,
00529  double                 * x_dist, 
00530  double                 * y_dist)
00531 {
00532     int             ix1;
00533     int             ix2;
00534     int             iy1;
00535     int             iy2;
00536     int             nx;
00537     int             ny;
00538     double          x1_pos;
00539     double          x2_pos;
00540     double          y1_pos;
00541     double          y2_pos;
00542     double          dx11;
00543     double          dx12;
00544     double          dx21;
00545     double          dx22;
00546     double          dy11;
00547     double          dy12;
00548     double          dy21;
00549     double          dy22;
00550     int             isnull;
00551 
00552     /* Get the size of the distortion images */
00553     nx = cpl_image_get_size_x(distortion->dist_x);
00554     ny = cpl_image_get_size_y(distortion->dist_x);
00555 
00556     //This uses bilinear interpolation
00557     //Get lower left corner. This assumes CRPIX? =1 and ix, iy start with 0
00558     ix1 = (int)floor((x_pos - distortion->x_crval)/distortion->x_cdelt);
00559     iy1 = (int)floor((y_pos - distortion->y_crval)/distortion->y_cdelt);
00560     //Handle extrapolation
00561     if(ix1 < 0)
00562         ix1 = 0;
00563     if(ix1 >= nx - 1)
00564         ix1 = nx - 2;
00565     if(iy1 < 0)
00566         iy1 = 0;
00567     if(iy1 >= ny - 1)
00568         iy1 = ny - 2;
00569     //Get upper right corner
00570     ix2 = ix1 + 1;
00571     iy2 = iy1 + 1;
00572     
00573     //Get the position values
00574     //This implies that CRPIX? = 1 and that ix, iy start at 0.
00575     x1_pos = ix1 * distortion->x_cdelt + distortion->x_crval; 
00576     x2_pos = ix2 * distortion->x_cdelt + distortion->x_crval; 
00577     y1_pos = iy1 * distortion->y_cdelt + distortion->y_crval; 
00578     y2_pos = iy2 * distortion->y_cdelt + distortion->y_crval; 
00579     
00580     //Get the values used to interpolate
00581     //The +1 is because cpl_image_get uses FITS convention
00582     dx11 = cpl_image_get(distortion->dist_x, ix1 + 1, iy1 + 1, &isnull);  
00583     dx21 = cpl_image_get(distortion->dist_x, ix2 + 1, iy1 + 1, &isnull);  
00584     dx12 = cpl_image_get(distortion->dist_x, ix1 + 1, iy2 + 1, &isnull);  
00585     dx22 = cpl_image_get(distortion->dist_x, ix2 + 1, iy2 + 1, &isnull);  
00586     dy11 = cpl_image_get(distortion->dist_y, ix1 + 1, iy1 + 1, &isnull);  
00587     dy21 = cpl_image_get(distortion->dist_y, ix2 + 1, iy1 + 1, &isnull);  
00588     dy12 = cpl_image_get(distortion->dist_y, ix1 + 1, iy2 + 1, &isnull);  
00589     dy22 = cpl_image_get(distortion->dist_y, ix2 + 1, iy2 + 1, &isnull);  
00590     
00591     
00592     //Compute the final values
00593     *x_dist = dx11 * (x2_pos - x_pos) * (y2_pos - y_pos) +
00594               dx21 * (x_pos - x1_pos) * (y2_pos - y_pos) +
00595               dx12 * (x2_pos - x_pos) * (y_pos - y1_pos) +
00596               dx22 * (x_pos - x1_pos) * (y_pos - y1_pos);
00597     *x_dist = *x_dist / (x2_pos -x1_pos) / (y2_pos -y1_pos);
00598     *y_dist = dy11 * (x2_pos - x_pos) * (y2_pos - y_pos) +
00599               dy21 * (x_pos - x1_pos) * (y2_pos - y_pos) +
00600               dy12 * (x2_pos - x_pos) * (y_pos - y1_pos) +
00601               dy22 * (x_pos - x1_pos) * (y_pos - y1_pos);
00602     *y_dist = *y_dist / (x2_pos -x1_pos) / (y2_pos -y1_pos);
00603 
00604     return 0;
00605 
00606 }
00607 
00608 /*----------------------------------------------------------------------------*/
00615 /*----------------------------------------------------------------------------*/
00616 int hawki_distortion_apply_maps
00617 (cpl_imagelist   *   ilist, 
00618  cpl_image       **  dist_x,
00619  cpl_image       **  dist_y)
00620 {
00621     cpl_image       *   corr[HAWKI_NB_DETECTORS] ;
00622     int                 i, j ;
00623 
00624     /* Test entries */
00625     if (ilist == NULL) return -1 ;
00626     if (dist_x == NULL) return -1 ;
00627     if (dist_y == NULL) return -1 ;
00628 
00629     /* Loop on the 4 chips */
00630     for (i=0 ; i<HAWKI_NB_DETECTORS ; i++)
00631     {
00632         cpl_image * cur_image;
00633 
00634         /* Get the current image */
00635         cur_image = cpl_imagelist_get(ilist, i);
00636         
00637         /* Correct this image */
00638         corr[i] = hawki_distortion_correct_detector(cur_image,dist_x[i],dist_y[i]);
00639         if(corr[i] == NULL)
00640         {
00641             cpl_msg_error(__func__, "Cannot correct the distortion") ;
00642             for (j=0 ; j<i ; j++) cpl_image_delete(corr[j]) ;
00643             return -1 ;
00644         }
00645     }
00646 
00647     /* Store the results */
00648     for (i=0 ; i<HAWKI_NB_DETECTORS ; i++) 
00649         cpl_imagelist_set(ilist, corr[i], i);
00650     
00651     /* Return */
00652     return 0 ;
00653 }
00654 
00657 /*----------------------------------------------------------------------------*/
00665 /*----------------------------------------------------------------------------*/
00666 int hawki_distortion_create_maps_detector
00667 (const hawki_distortion * distortion,
00668  cpl_image              * dist_x,
00669  cpl_image              * dist_y)
00670 {
00671     double              *   pdx;
00672     double              *   pdy;
00673     int                     nx, ny;
00674     int                     pos;
00675     int                     i, j;
00676 
00677     /* Test entries */
00678     if (distortion == NULL) return -1 ;
00679     if (dist_x == NULL) return -1 ;
00680     if (dist_y == NULL) return -1 ;
00681 
00682     /* Initialise */
00683     nx = cpl_image_get_size_x(dist_x) ;
00684     ny = cpl_image_get_size_y(dist_x) ;
00685     if (cpl_image_get_size_x(dist_y) != nx) return -1 ;
00686     if (cpl_image_get_size_y(dist_y) != ny) return -1 ;
00687   
00688     /* Access to the data */
00689     pdx = cpl_image_get_data_double(dist_x) ;
00690     pdy = cpl_image_get_data_double(dist_y) ;
00691 
00692     /* Loop on the pixels */
00693     for (j=0 ; j<ny ; j++) {
00694         for (i=0 ; i<nx ; i++) {
00695             double x_dist;
00696             double y_dist;
00697 
00698             pos = i+j*nx ;
00699 
00700             hawki_distortion_interpolate_distortion
00701                 (distortion, (double)i, (double)j, &x_dist, &y_dist);
00702             
00703             pdx[pos] = x_dist;
00704             pdy[pos] = y_dist;
00705         }
00706     }
00707 
00708     return 0 ; 
00709 }
00710 
00711 /*----------------------------------------------------------------------------*/
00735 /*----------------------------------------------------------------------------*/
00736 hawki_distortion *  hawki_distortion_compute_solution
00737 (const cpl_table       ** ref_catalogues,
00738  const cpl_bivector     * initial_offsets,
00739  const cpl_table        * matching_sets,
00740  int                      ncats,
00741  int                      detector_nx,
00742  int                      detector_ny,
00743  int                      grid_size,
00744  const hawki_distortion * initial_guess,
00745  double                 * rms)
00746 {
00747 #ifdef HAVE_LIBGSL
00748     
00749     hawki_distortion                          * distortion;
00750     cpl_bivector                              * offsets;
00751     gsl_multimin_fminimizer                   * minimizer;
00752     gsl_vector                                * step_size;
00753     gsl_vector                                * init_param;
00754     gsl_multimin_function                       minimize_function;
00755     
00756     int                                         nfitparam = 0;
00757     int                                         iparam;
00758     double                                      tolerance = HAWKI_DISTORTION_TOLERANCE;
00759     double                                      tolerance2 = HAWKI_DISTORTION_TOLERANCE2;
00760     int                                         minimizer_status;
00761     int                                         iter = 0;
00762     struct _hawki_distortion_obj_function_args_ args;
00763 
00764     /* Initialize the distortion solution */
00765     if(initial_guess == NULL)
00766     {
00767         distortion = hawki_distortion_grid_new
00768             (detector_nx, detector_ny, grid_size);
00769     }
00770     else
00771     {
00772         distortion = cpl_malloc(sizeof(hawki_distortion));
00773         distortion->dist_x  = cpl_image_duplicate(initial_guess->dist_x);
00774         distortion->dist_y  = cpl_image_duplicate(initial_guess->dist_y);
00775         distortion->x_cdelt = initial_guess->x_cdelt;
00776         distortion->x_crval = initial_guess->x_crval;
00777         distortion->y_cdelt = initial_guess->y_cdelt;
00778         distortion->y_crval = initial_guess->y_crval;
00779         //We have to fit all the distortion coefficients plus
00780         //the offsets of the catalogues
00781     }
00782     //We have to fit all the distortion coefficients plus
00783     //the offsets of the catalogues
00784     nfitparam = grid_size * grid_size * 2 + ncats * 2;
00785     offsets = cpl_bivector_duplicate(initial_offsets);
00786     if(cpl_table_get_nrow(matching_sets) * 2 < nfitparam)
00787     {
00788         cpl_msg_error(__func__,"Too few matches to compute distortion (< %d)",
00789                       nfitparam);
00790         hawki_distortion_delete(distortion);
00791         return NULL;
00792     }
00793 
00794     /* Setup function to minimize */
00795     args.ref_catalogues = ref_catalogues;
00796     args.matching_sets  = matching_sets;
00797     args.distortion     = distortion;
00798     args.offsets        = offsets;
00799     args.ncats          = ncats;
00800     
00801     minimize_function.f      = hawki_distortion_gsl_obj_function;
00802     minimize_function.n      = nfitparam;
00803     minimize_function.params = &args;
00804     
00805     /* Setup minimizer */
00806     minimizer = gsl_multimin_fminimizer_alloc
00807         (gsl_multimin_fminimizer_nmsimplex, nfitparam);
00808     step_size = gsl_vector_alloc(nfitparam);
00809     init_param = gsl_vector_alloc(nfitparam);
00810     for(iparam = 0; iparam < grid_size * grid_size * 2; ++iparam)
00811         gsl_vector_set(step_size, iparam, 5);
00812     for(iparam = grid_size * grid_size * 2;
00813         iparam < nfitparam; ++iparam)
00814         gsl_vector_set(step_size, iparam, 1);
00815     hawki_distortion_update_param_from_solution(init_param, distortion);
00816     hawki_distortion_update_param_from_offsets(init_param, offsets);
00817     gsl_multimin_fminimizer_set(minimizer, &minimize_function,
00818                                 init_param, step_size);
00819 
00820     do
00821     {
00822         ++iter;
00823         minimizer_status = gsl_multimin_fminimizer_iterate (minimizer);
00824         if(minimizer_status)
00825             break;
00826         minimizer_status = gsl_multimin_test_size
00827             (gsl_multimin_fminimizer_size(minimizer), tolerance);
00828         cpl_msg_debug(__func__,"Iteration %d Minimum: %g",
00829                       iter, gsl_multimin_fminimizer_minimum(minimizer));
00830     }
00831     while (minimizer_status == GSL_CONTINUE && iter < HAWKI_DISTORTION_MAX_ITER);
00832 
00833     cpl_msg_warning(__func__, "rms before computing %f", hawki_distortion_compute_rms(ref_catalogues, offsets, 
00834                                         matching_sets, ncats, distortion));
00835 
00836     
00837     //Do it again to avoid local minimum
00838     gsl_multimin_fminimizer_set(minimizer, &minimize_function,
00839                                 gsl_multimin_fminimizer_x(minimizer), step_size);
00840     iter = 0;
00841     do
00842     {
00843         ++iter;
00844         minimizer_status = gsl_multimin_fminimizer_iterate (minimizer);
00845         if(minimizer_status)
00846             break;
00847         minimizer_status = gsl_multimin_test_size
00848             (gsl_multimin_fminimizer_size(minimizer), tolerance2);
00849         cpl_msg_debug(__func__,"2nd run Iteration %d Minimum: %g",
00850                       iter, gsl_multimin_fminimizer_minimum(minimizer));
00851     }
00852     while (minimizer_status == GSL_CONTINUE && iter < HAWKI_DISTORTION_MAX_ITER2);
00853 
00854     /* Update the distortion solution */
00855     hawki_distortion_update_solution_from_param
00856         (distortion, gsl_multimin_fminimizer_x(minimizer));
00857     hawki_distortion_update_offsets_from_param
00858         (offsets, gsl_multimin_fminimizer_x(minimizer));
00859     
00860     *rms = hawki_distortion_compute_rms(ref_catalogues, offsets, 
00861                                         matching_sets, ncats, distortion);
00862     
00863     /* Free and return */
00864     gsl_multimin_fminimizer_free(minimizer);
00865     gsl_vector_free(init_param);
00866     gsl_vector_free(step_size);
00867     cpl_bivector_delete(offsets);
00868 
00869     return distortion;
00870 #else
00871     cpl_msg_error(__func__,"Not compiled with GSL support.");
00872     return NULL;
00873 #endif
00874 }
00875 
00876 #ifdef HAVE_LIBGSL
00877 /*----------------------------------------------------------------------------*/
00885 /*----------------------------------------------------------------------------*/
00886 double hawki_distortion_gsl_obj_function
00887 (const gsl_vector * dist_param,
00888  void             * args)
00889 {
00890     struct _hawki_distortion_obj_function_args_  args_struct;
00891     const cpl_table                           ** ref_catalogues;
00892     const cpl_table                            * matching_sets;
00893     hawki_distortion                           * distortion;
00894     cpl_bivector                               * offsets;
00895     int                                          ncats;
00896     double                                       rms;
00897     double                                       objective_function;
00898 
00899     
00900     
00901     args_struct    = *(struct _hawki_distortion_obj_function_args_ * )args;
00902     ref_catalogues = args_struct.ref_catalogues; 
00903     matching_sets  = args_struct.matching_sets; 
00904     distortion     = args_struct.distortion;
00905     offsets        = args_struct.offsets;
00906     ncats          = args_struct.ncats; 
00907     
00908     hawki_distortion_update_solution_from_param(distortion, dist_param);
00909     hawki_distortion_update_offsets_from_param(offsets, dist_param);
00910     
00911     rms = hawki_distortion_compute_rms(ref_catalogues, offsets, 
00912                                        matching_sets, ncats, distortion);
00913     
00914     objective_function = rms; 
00915     
00916     
00917     cpl_msg_debug(__func__,"Objective function: %g", objective_function);
00918     return objective_function;
00919 }
00920 #endif
00921 
00922 /*----------------------------------------------------------------------------*/
00928 /*----------------------------------------------------------------------------*/
00929 double hawki_distortion_compute_rms
00930 (const cpl_table   ** ref_catalogues,
00931  const cpl_bivector * cat_offsets,
00932  const cpl_table    * matching_sets,
00933  int                  ncats,
00934  hawki_distortion   * distortion)
00935 {
00936     int                icat;
00937     int                imatch;
00938     int                nmatch;
00939     double             rms = 0;
00940     const double     * x_cat_offsets;
00941     const double     * y_cat_offsets;
00942     const double    ** x_cat_cols; 
00943     const double    ** y_cat_cols;
00944     const cpl_array ** match_arrays;
00945     double          ** x_pos_values;
00946     double          ** y_pos_values;
00947     int             ** pos_flag;
00948 
00949     
00950     
00951     nmatch = cpl_table_get_nrow(matching_sets);
00952 
00953     x_cat_offsets = cpl_vector_get_data_const
00954         (cpl_bivector_get_x_const(cat_offsets));
00955     y_cat_offsets = cpl_vector_get_data_const
00956         (cpl_bivector_get_y_const(cat_offsets));
00957 
00958     x_cat_cols = cpl_malloc(sizeof(double *) * ncats);
00959     y_cat_cols = cpl_malloc(sizeof(double *) * ncats);
00960     for(icat = 0; icat < ncats; ++icat)
00961     {
00962         x_cat_cols[icat] = cpl_table_get_data_double_const(ref_catalogues[icat],
00963                                                            HAWKI_COL_OBJ_POSX);
00964         y_cat_cols[icat] = cpl_table_get_data_double_const(ref_catalogues[icat],
00965                                                            HAWKI_COL_OBJ_POSY);
00966     }
00967 
00968     match_arrays = cpl_malloc(sizeof(cpl_array *) * nmatch);
00969     x_pos_values = cpl_malloc(sizeof(double *) * nmatch);
00970     y_pos_values = cpl_malloc(sizeof(double *) * nmatch);
00971     pos_flag     = cpl_malloc(sizeof(int *) * nmatch);
00972     for(imatch = 0; imatch < nmatch; ++imatch)
00973     {
00974         match_arrays[imatch] = cpl_table_get_array(matching_sets, 
00975                 HAWKI_COL_MATCHING_SETS, imatch);
00976         x_pos_values[imatch] = cpl_malloc(sizeof(double) * ncats);
00977         y_pos_values[imatch] = cpl_malloc(sizeof(double) * ncats);
00978         pos_flag[imatch]     = cpl_malloc(sizeof(int) * ncats);
00979     }
00980 
00981 #ifdef _OPENMP    
00982 #pragma omp parallel for  private(icat,imatch) reduction(+:rms)
00983 #endif    
00984     for(imatch = 0; imatch < nmatch; ++imatch)
00985     {
00986         int                 nstddev;
00987         double              var_x;
00988         double              var_y;
00989                 
00990         for(icat = 0; icat < ncats; ++icat)
00991         {
00992             int    iobj;
00993             double x_cat_offset;
00994             double y_cat_offset;
00995             
00996             x_cat_offset = x_cat_offsets[icat];
00997             y_cat_offset = y_cat_offsets[icat];            
00998 
00999             if((iobj = cpl_array_get(match_arrays[imatch], icat, NULL)) != -1)
01000             {
01001                 double x_cat;
01002                 double y_cat;
01003                 double x_dist_corr;
01004                 double y_dist_corr;
01005                 double x_dist;
01006                 double y_dist;
01007                 double x_glob;
01008                 double y_glob;
01009 
01010                 x_cat = x_cat_cols[icat][iobj];
01011                 y_cat = y_cat_cols[icat][iobj];
01012 
01013 
01014                 //These 4 lines of code are from hawki_distortion_correct_coords.
01015                 //The are repeated here to avoid a cpl call, which is not thread-safe
01016                 //Two checks to ensure thread-safety:
01017                 //We have ensured outside the loop that distortion->dist_x and
01018                 //distortion->dist_y are not null. 
01019                 //We have checked outside the loop the mask has all the points valid
01020                 hawki_distortion_interpolate_distortion
01021                     (distortion, x_cat, y_cat, &x_dist, &y_dist);
01022                 x_dist_corr = x_cat - x_dist;
01023                 y_dist_corr = y_cat - y_dist;
01024 
01025 
01026                 x_glob = x_dist_corr + x_cat_offset;
01027                 y_glob = y_dist_corr + y_cat_offset;
01028                 x_pos_values[imatch][icat] = x_glob;
01029                 y_pos_values[imatch][icat] = y_glob;
01030                 pos_flag[imatch][icat] = 1;
01031            }
01032            else
01033                pos_flag[imatch][icat] = 0;
01034         }
01035 
01036         hawki_distortion_get_flag_vars(x_pos_values[imatch], 
01037                                        y_pos_values[imatch], pos_flag[imatch],
01038                                        ncats, &nstddev, &var_x, &var_y);
01039         
01040         //The rms is counted as many times as this star is the list of catalogs.
01041         rms += sqrt(var_x + var_y) * nstddev;
01042         
01043     }
01044     cpl_free(x_cat_cols);
01045     cpl_free(y_cat_cols);
01046     for(imatch = 0; imatch < nmatch; ++imatch)
01047     {
01048         cpl_free(x_pos_values[imatch]);
01049         cpl_free(y_pos_values[imatch]);
01050         cpl_free(pos_flag[imatch]);
01051     }
01052     cpl_free(x_pos_values);
01053     cpl_free(y_pos_values);
01054     cpl_free(pos_flag);
01055     cpl_free(match_arrays);
01056     
01057     return rms;
01058 }
01059 
01060 #ifdef HAVE_LIBGSL
01061 /*----------------------------------------------------------------------------*/
01067 /*----------------------------------------------------------------------------*/
01068 int hawki_distortion_update_solution_from_param
01069 (hawki_distortion  * distortion,
01070  const gsl_vector  * dist_param)
01071 {
01072     int     ipoint;
01073     int     ix;
01074     int     iy;
01075     int     nx;
01076     int     ny;
01077     double  x_dist_mean; 
01078     double  y_dist_mean; 
01079     
01080     nx  = cpl_image_get_size_x(distortion->dist_x);
01081     ny  = cpl_image_get_size_y(distortion->dist_x);
01082     for(ix = 0; ix < nx; ++ix)
01083         for(iy = 0; iy < ny; ++iy)
01084         {
01085             ipoint = ix + iy * nx;
01086             cpl_image_set(distortion->dist_x, ix+1, iy+1, 
01087                           gsl_vector_get(dist_param, ipoint * 2));
01088             cpl_image_set(distortion->dist_y, ix+1, iy+1, 
01089                           gsl_vector_get(dist_param, ipoint * 2 + 1));
01090         }
01091     
01092     /* Normalize to mean(distorsion) = 0 */
01093     x_dist_mean = cpl_image_get_mean(distortion->dist_x);
01094     y_dist_mean = cpl_image_get_mean(distortion->dist_y);
01095     cpl_image_subtract_scalar(distortion->dist_x, x_dist_mean);
01096     cpl_image_subtract_scalar(distortion->dist_y, y_dist_mean);
01097 
01098     return 0;
01099 }
01100 
01101 /*----------------------------------------------------------------------------*/
01107 /*----------------------------------------------------------------------------*/
01108 int hawki_distortion_update_offsets_from_param
01109 (cpl_bivector      * offsets,
01110  const gsl_vector  * dist_param)
01111 {
01112     int     i;
01113     int     ncats;
01114     int     nparam;
01115     
01116     ncats  = cpl_bivector_get_size(offsets);
01117     nparam = dist_param->size;
01118     for(i = 0; i < ncats; ++i)
01119     {
01120         cpl_vector_set(cpl_bivector_get_x(offsets), i,  
01121                        gsl_vector_get(dist_param, nparam - 2 * ncats + 2 * i));
01122         cpl_vector_set(cpl_bivector_get_y(offsets), i,  
01123                        gsl_vector_get(dist_param, nparam - 2 * ncats + 2 * i + 1));
01124     }
01125     
01126     return 0;
01127 }
01128 
01129 /*----------------------------------------------------------------------------*/
01135 /*----------------------------------------------------------------------------*/
01136 int hawki_distortion_update_param_from_solution
01137 (gsl_vector              * dist_param,
01138  const hawki_distortion  * distortion)
01139 {
01140     int  ipoint;
01141     int  ix;
01142     int  iy;
01143     int  nx;
01144     int  ny;
01145     int  isnull;
01146     
01147     nx  = cpl_image_get_size_x(distortion->dist_x);
01148     ny  = cpl_image_get_size_y(distortion->dist_y);
01149     for(ix = 0; ix < nx; ++ix)
01150         for(iy = 0; iy < ny; ++iy)
01151         {
01152             ipoint = ix + iy * nx;
01153             gsl_vector_set(dist_param, ipoint * 2, 
01154                            cpl_image_get(distortion->dist_x, ix+1, iy+1, &isnull));
01155             gsl_vector_set(dist_param, ipoint * 2 + 1, 
01156                            cpl_image_get(distortion->dist_y, ix+1, iy+1, &isnull));
01157     }
01158     
01159     return 0;
01160 }
01161 
01162 /*----------------------------------------------------------------------------*/
01167 /*----------------------------------------------------------------------------*/
01168 int hawki_distortion_update_param_from_offsets
01169 (gsl_vector              * dist_param,
01170  const cpl_bivector      * offsets)
01171 {
01172     int     i;
01173     int     ncats;
01174     int     nparam;
01175     
01176     ncats  = cpl_bivector_get_size(offsets);
01177     nparam = dist_param->size;
01178     for(i = 0; i < ncats; ++i)
01179     {
01180         gsl_vector_set(dist_param, nparam - 2 * ncats + 2 * i, 
01181                        cpl_vector_get(cpl_bivector_get_x_const(offsets), i));
01182         gsl_vector_set(dist_param, nparam - 2 * ncats + 2 * i + 1, 
01183                        cpl_vector_get(cpl_bivector_get_y_const(offsets), i));
01184     }
01185 
01186     return 0;
01187 }
01188 #endif
01189 
01190 /*----------------------------------------------------------------------------*/
01196 /*----------------------------------------------------------------------------*/
01197 void hawki_distortion_get_flag_vars
01198 (double * x_val, double * y_val, int *pos_flag,
01199  int nvals, int *nvalid, double *var_x, double *var_y)
01200 {
01201     double varsum_x = 0.0;
01202     double varsum_y = 0.0;
01203     double mean_x = 0.0;
01204     double mean_y = 0.0;
01205     int i;
01206 
01207     *nvalid = 0;
01208     for (i=0; i < nvals; i++)
01209     {
01210         if(pos_flag[i] == 1)
01211         {
01212             const double delta_x = (double)x_val[i] - mean_x;
01213             const double delta_y = (double)y_val[i] - mean_y;
01214 
01215             varsum_x += *nvalid * delta_x * delta_x / (*nvalid + 1.0);
01216             varsum_y += *nvalid * delta_y * delta_y / (*nvalid + 1.0);
01217             mean_x   += delta_x / (*nvalid + 1.0);
01218             mean_y   += delta_y / (*nvalid + 1.0);
01219             (*nvalid)++;
01220         }
01221     }
01222 
01223     /* Compute the bias-corrected standard deviation.
01224        - With the recurrence relation rounding can likely not cause
01225        the variance to become negative, but check just to be safe */
01226     *var_x = varsum_x / (double) (*nvalid - 1);
01227     *var_y = varsum_y / (double) (*nvalid - 1);
01228 }