hawki_distortion.c

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

Generated on Thu Feb 17 17:13:07 2011 for HAWKI Pipeline Reference Manual by  doxygen 1.4.7