HAWKI Pipeline Reference Manual  1.8.12
hawki_distortion.c
1 /* $Id: hawki_distortion.c,v 1.32 2011/02/23 11:49:37 cgarcia Exp $
2  *
3  * This file is part of the HAWKI Pipeline
4  * Copyright (C) 2002,2003 European Southern Observatory
5  *
6  * This program is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * This program is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this program; if not, write to the Free Software
18  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19  */
20 
21 /*
22  * $Author: cgarcia $
23  * $Date: 2011/02/23 11:49:37 $
24  * $Revision: 1.32 $
25  * $Name: hawki-1_8_12 $
26  */
27 
28 #ifdef HAVE_CONFIG_H
29 #include <config.h>
30 #endif
31 
32 //Minimization algorithm hard-coded constants
33 #define HAWKI_DISTORTION_MAX_ITER 10000
34 #define HAWKI_DISTORTION_TOLERANCE 0.001
35 #define HAWKI_DISTORTION_MAX_ITER2 100000
36 #define HAWKI_DISTORTION_TOLERANCE2 0.0001
37 
38 /*-----------------------------------------------------------------------------
39  Includes
40  -----------------------------------------------------------------------------*/
41 
42 #include <math.h>
43 #include <cpl.h>
44 #include <cxdeque.h>
45 #ifdef HAVE_LIBGSL
46 #include <gsl/gsl_multimin.h>
47 #endif
48 
49 #include "hawki_distortion.h"
50 #include "hawki_dfs.h"
51 #include "hawki_utils.h"
52 #include "hawki_load.h"
53 
54 
55 
56 /*----------------------------------------------------------------------------*/
60 /*----------------------------------------------------------------------------*/
61 
70 {
71  const cpl_table ** ref_catalogues;
72  const cpl_table * matching_sets;
73  cpl_bivector * offsets;
74  hawki_distortion * distortion;
75  int ncats;
76 };
77 
78 
79 //Private functions
80 
82 (const hawki_distortion * distortion,
83  double x_pos,
84  double y_pos,
85  double * x_dist,
86  double * y_dist);
87 
89 (const cpl_table ** ref_catalogues,
90  const cpl_bivector * cat_offsets,
91  const cpl_table * matching_sets,
92  int ncats,
93  hawki_distortion * distortion);
94 
95 #ifdef HAVE_LIBGSL
96 double hawki_distortion_gsl_obj_function
97 (const gsl_vector * dist_param,
98  void * args);
99 
100 int hawki_distortion_update_solution_from_param
101 (hawki_distortion * distortion,
102  const gsl_vector * dist_param);
103 
104 int hawki_distortion_update_offsets_from_param
105 (cpl_bivector * offsets,
106  const gsl_vector * dist_param);
107 
108 int hawki_distortion_update_param_from_solution
109 (gsl_vector * dist_param,
110  const hawki_distortion * distortion);
111 
112 int hawki_distortion_update_param_from_offsets
113 (gsl_vector * dist_param,
114  const cpl_bivector * offsets);
115 #endif
116 
118 (double * x_val, double * y_val, int *pos_flag,
119  int nvals, int *nvalid, double *var_x, double *var_y);
120 
121 
122 /*----------------------------------------------------------------------------*/
130 /*----------------------------------------------------------------------------*/
131 hawki_distortion * hawki_distortion_grid_new
132 (int detector_nx,
133  int detector_ny,
134  int grid_size)
135 {
136  hawki_distortion * distortion;
137 
138  //Allocate the structure
139  distortion = cpl_malloc(sizeof(hawki_distortion));
140 
141  //Allocate the images
142  distortion->dist_x = cpl_image_new
143  (grid_size, grid_size, CPL_TYPE_FLOAT);
144  distortion->dist_y = cpl_image_new
145  (grid_size, grid_size, CPL_TYPE_FLOAT);
146 
147  //Create the transformation between distortion images and the detector
148  distortion->x_cdelt = detector_nx / (double)grid_size;
149  distortion->y_cdelt = detector_ny / (double)grid_size;
150  distortion->x_crval = 0.5 + 0.5 * distortion->x_cdelt;
151  distortion->y_crval = 0.5 + 0.5 * distortion->y_cdelt;
152 
153  return distortion;
154 }
155 
156 /*----------------------------------------------------------------------------*/
161 /*----------------------------------------------------------------------------*/
163 (hawki_distortion * distortion)
164 {
165  if(distortion == NULL)
166  return;
167  cpl_image_delete(distortion->dist_x);
168  cpl_image_delete(distortion->dist_y);
169  cpl_free(distortion);
170 }
171 
172 /*----------------------------------------------------------------------------*/
180 /*----------------------------------------------------------------------------*/
181 hawki_distortion * hawki_distortion_load
182 (const cpl_frame * dist_x,
183  const cpl_frame * dist_y,
184  int idet)
185 {
186  const char * file_dist_x;
187  const char * file_dist_y;
188  hawki_distortion * distortion;
189  int iext;
190  cpl_propertylist * plist;
191 
192  //Allocate the structure
193  distortion = cpl_malloc(sizeof(hawki_distortion));
194 
195  //Read the images
196  file_dist_x = cpl_frame_get_filename(dist_x);
197  file_dist_y = cpl_frame_get_filename(dist_y);
198  distortion->dist_x = hawki_load_frame_detector
199  (dist_x, idet, CPL_TYPE_FLOAT);
200  distortion->dist_y = hawki_load_frame_detector
201  (dist_y, idet, CPL_TYPE_FLOAT);
202 
203  //Read the WCS keywords
204  iext = hawki_get_ext_from_detector(file_dist_x, idet);
205  plist = cpl_propertylist_load(file_dist_x, iext);
206  distortion->x_crval = cpl_propertylist_get_double(plist, "CRVAL1");
207  distortion->x_cdelt = cpl_propertylist_get_double(plist, "CDELT1");
208  distortion->y_crval = cpl_propertylist_get_double(plist, "CRVAL2");
209  distortion->y_cdelt = cpl_propertylist_get_double(plist, "CDELT2");
210  if(cpl_propertylist_get_double(plist, "CRPIX1") != 1 ||
211  cpl_propertylist_get_double(plist, "CRPIX2") != 1)
212  {
213  cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
214  __FILE__, __LINE__,"Wrong CRPIX? keywords");
215  cpl_image_delete(distortion->dist_x);
216  cpl_image_delete(distortion->dist_y);
217  cpl_propertylist_delete(plist);
218  cpl_free(distortion);
219  return NULL;
220  }
221  cpl_propertylist_delete(plist);
222  //Check that the keywords in X and Y are compatibles;
223  plist = cpl_propertylist_load(file_dist_y, iext);
224  if(distortion->x_crval != cpl_propertylist_get_double(plist, "CRVAL1") ||
225  distortion->x_cdelt != cpl_propertylist_get_double(plist, "CDELT1") ||
226  distortion->y_crval != cpl_propertylist_get_double(plist, "CRVAL2") ||
227  distortion->y_cdelt != cpl_propertylist_get_double(plist, "CDELT2") ||
228  cpl_propertylist_get_double(plist, "CRPIX1") != 1 ||
229  cpl_propertylist_get_double(plist, "CRPIX2") != 1)
230  {
231  cpl_error_set_message_macro(cpl_func, CPL_ERROR_ILLEGAL_INPUT, __FILE__,
232  __LINE__,"WCS keywords mismatch in X and Y distortions");
233  cpl_image_delete(distortion->dist_x);
234  cpl_image_delete(distortion->dist_y);
235  cpl_propertylist_delete(plist);
236  cpl_free(distortion);
237  return NULL;
238  }
239  cpl_propertylist_delete(plist);
240 
241  return distortion;
242 }
243 
244 /*----------------------------------------------------------------------------*/
250 /*----------------------------------------------------------------------------*/
252 (const hawki_distortion * distortion)
253 {
254  if(distortion == NULL)
255  {
256  cpl_error_set(__func__,CPL_ERROR_ILLEGAL_INPUT);
257  }
258  return cpl_image_get_size_x(distortion->dist_x);
259 }
260 
261 /*----------------------------------------------------------------------------*/
267 /*----------------------------------------------------------------------------*/
269 (const hawki_distortion * distortion)
270 {
271  if(distortion == NULL)
272  {
273  cpl_error_set(__func__,CPL_ERROR_ILLEGAL_INPUT);
274  }
275  return cpl_image_get_size_y(distortion->dist_x);
276 }
277 
278 /*----------------------------------------------------------------------------*/
285 /*----------------------------------------------------------------------------*/
287 (cpl_image ** ilist,
288  const cpl_frame * frame_dist_x,
289  const cpl_frame * frame_dist_y)
290 {
291  cpl_image * corr[HAWKI_NB_DETECTORS];
292  hawki_distortion * distortion;
293  cpl_image * dist_x;
294  cpl_image * dist_y;
295  int idet, j ;
296 
297  /* Test entries */
298  if (ilist == NULL) return -1 ;
299  if (frame_dist_x == NULL) return -1 ;
300  if (frame_dist_y == NULL) return -1 ;
301 
302  /* Loop on the 4 chips */
303  for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
304  {
305  int nx;
306  int ny;
307 
308  /* Get the image size */
309  nx = cpl_image_get_size_x(ilist[idet]);
310  ny = cpl_image_get_size_y(ilist[idet]);
311 
312  /* Load the distortion */
313  if ((distortion = hawki_distortion_load
314  (frame_dist_x, frame_dist_y, idet + 1)) == NULL)
315  {
316  cpl_msg_error(__func__, "Cannot load the distortion for chip %d",
317  idet+1) ;
318  for (j=0 ; j<idet ; j++) cpl_image_delete(corr[j]) ;
319  return -1 ;
320  }
321 
322  /* Create the offsets images */
323  dist_x = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
324  dist_y = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
325  if (hawki_distortion_create_maps_detector
326  (distortion, dist_x, dist_y))
327  {
328  cpl_msg_error(__func__, "Cannot create the distortion maps") ;
329  cpl_image_delete(dist_x);
330  cpl_image_delete(dist_y);
331  for (j=0 ; j<idet ; j++) cpl_image_delete(corr[j]) ;
332  return -1;
333  }
334 
335  /* Correct this image */
336  corr[idet] = hawki_distortion_correct_detector(ilist[idet], dist_x, dist_y);
337  if(corr[idet] == NULL)
338  {
339  cpl_msg_error(__func__, "Cannot correct the distortion") ;
340  hawki_distortion_delete(distortion);
341  cpl_image_delete(dist_x);
342  cpl_image_delete(dist_y);
343  for (j=0 ; j<idet; j++) cpl_image_delete(corr[j]) ;
344  return -1 ;
345  }
346  hawki_distortion_delete(distortion);
347  cpl_image_delete(dist_x) ;
348  cpl_image_delete(dist_y);
349  }
350 
351  /* Store the results */
352  for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
353  {
354  cpl_image_delete(ilist[idet]) ;
355  ilist[idet] = corr[idet] ;
356  }
357 
358  /* Return */
359  return 0 ;
360 }
361 
362 /*----------------------------------------------------------------------------*/
369 /*----------------------------------------------------------------------------*/
371 (cpl_image * image,
372  cpl_image * dist_x,
373  cpl_image * dist_y)
374 {
375  cpl_image * corr;
376  cpl_vector * profile ;
377 
378  /* Test entries */
379  if (image == NULL) return NULL;
380  if (dist_x == NULL) return NULL;
381  if (dist_y == NULL) return NULL;
382 
383  /* Create the output image */
384  corr = cpl_image_new(cpl_image_get_size_x(image),
385  cpl_image_get_size_y(image), CPL_TYPE_FLOAT) ;
386 
387  /* Create the interpolation profile */
388  profile = cpl_vector_new(CPL_KERNEL_DEF_SAMPLES) ;
389  cpl_vector_fill_kernel_profile(profile, CPL_KERNEL_DEFAULT,
390  CPL_KERNEL_DEF_WIDTH) ;
391 
392  /* Apply the distortion */
393  if (cpl_image_warp(corr, image, dist_x, dist_y, profile,
394  CPL_KERNEL_DEF_WIDTH, profile,
395  CPL_KERNEL_DEF_WIDTH) != CPL_ERROR_NONE)
396  {
397  cpl_msg_error(__func__, "Cannot warp the image") ;
398  cpl_image_delete(corr) ;
399  cpl_vector_delete(profile) ;
400  return NULL;
401  }
402  cpl_vector_delete(profile) ;
403 
404  /* Return */
405  return corr;
406 }
407 
408 /*----------------------------------------------------------------------------*/
421 /*----------------------------------------------------------------------------*/
423 (const hawki_distortion * distortion,
424  double x_pos,
425  double y_pos,
426  double * x_pos_distcorr,
427  double * y_pos_distcorr)
428 {
429  double x_dist;
430  double y_dist;
431 
432  if(distortion == NULL)
433  {
434  cpl_error_set("hawki_distortion_correct_coords", CPL_ERROR_ILLEGAL_INPUT);
435  return -1;
436  }
437 
439  (distortion, x_pos, y_pos, &x_dist, &y_dist);
440 
441  *x_pos_distcorr = x_pos - x_dist;
442  *y_pos_distcorr = y_pos - y_dist;
443 
444  return 0;
445 }
446 
447 
448 /*----------------------------------------------------------------------------*/
465 /*----------------------------------------------------------------------------*/
467 (const hawki_distortion * distortion,
468  double x_pos,
469  double y_pos,
470  double * x_pos_distinvcorr,
471  double * y_pos_distinvcorr)
472 {
473  double x_dist = 0;
474  double y_dist = 0;
475  int i;
476  int niter = 3;
477 
478  if(distortion == NULL)
479  {
480  cpl_error_set("hawki_distortion_inverse_correct_coords", CPL_ERROR_ILLEGAL_INPUT);
481  return -1;
482  }
483  for(i = 0; i < niter; ++i)
484  {
486  (distortion, x_pos + x_dist, y_pos + y_dist, &x_dist, &y_dist);
487  }
488 
489 
490  /* Apply the correction in the inverse direction */
491  *x_pos_distinvcorr = x_pos + x_dist;
492  *y_pos_distinvcorr = y_pos + y_dist;
493 
494  return 0;
495 }
496 
497 /*----------------------------------------------------------------------------*/
524 /*----------------------------------------------------------------------------*/
526 (const hawki_distortion * distortion,
527  double x_pos,
528  double y_pos,
529  double * x_dist,
530  double * y_dist)
531 {
532  int ix1;
533  int ix2;
534  int iy1;
535  int iy2;
536  int nx;
537  int ny;
538  double x1_pos;
539  double x2_pos;
540  double y1_pos;
541  double y2_pos;
542  double dx11;
543  double dx12;
544  double dx21;
545  double dx22;
546  double dy11;
547  double dy12;
548  double dy21;
549  double dy22;
550  int isnull;
551 
552  /* Get the size of the distortion images */
553  nx = cpl_image_get_size_x(distortion->dist_x);
554  ny = cpl_image_get_size_y(distortion->dist_x);
555 
556  //This uses bilinear interpolation
557  //Get lower left corner. This assumes CRPIX? =1 and ix, iy start with 0
558  ix1 = (int)floor((x_pos - distortion->x_crval)/distortion->x_cdelt);
559  iy1 = (int)floor((y_pos - distortion->y_crval)/distortion->y_cdelt);
560  //Handle extrapolation
561  if(ix1 < 0)
562  ix1 = 0;
563  if(ix1 >= nx - 1)
564  ix1 = nx - 2;
565  if(iy1 < 0)
566  iy1 = 0;
567  if(iy1 >= ny - 1)
568  iy1 = ny - 2;
569  //Get upper right corner
570  ix2 = ix1 + 1;
571  iy2 = iy1 + 1;
572 
573  //Get the position values
574  //This implies that CRPIX? = 1 and that ix, iy start at 0.
575  x1_pos = ix1 * distortion->x_cdelt + distortion->x_crval;
576  x2_pos = ix2 * distortion->x_cdelt + distortion->x_crval;
577  y1_pos = iy1 * distortion->y_cdelt + distortion->y_crval;
578  y2_pos = iy2 * distortion->y_cdelt + distortion->y_crval;
579 
580  //Get the values used to interpolate
581  //The +1 is because cpl_image_get uses FITS convention
582  dx11 = cpl_image_get(distortion->dist_x, ix1 + 1, iy1 + 1, &isnull);
583  dx21 = cpl_image_get(distortion->dist_x, ix2 + 1, iy1 + 1, &isnull);
584  dx12 = cpl_image_get(distortion->dist_x, ix1 + 1, iy2 + 1, &isnull);
585  dx22 = cpl_image_get(distortion->dist_x, ix2 + 1, iy2 + 1, &isnull);
586  dy11 = cpl_image_get(distortion->dist_y, ix1 + 1, iy1 + 1, &isnull);
587  dy21 = cpl_image_get(distortion->dist_y, ix2 + 1, iy1 + 1, &isnull);
588  dy12 = cpl_image_get(distortion->dist_y, ix1 + 1, iy2 + 1, &isnull);
589  dy22 = cpl_image_get(distortion->dist_y, ix2 + 1, iy2 + 1, &isnull);
590 
591 
592  //Compute the final values
593  *x_dist = dx11 * (x2_pos - x_pos) * (y2_pos - y_pos) +
594  dx21 * (x_pos - x1_pos) * (y2_pos - y_pos) +
595  dx12 * (x2_pos - x_pos) * (y_pos - y1_pos) +
596  dx22 * (x_pos - x1_pos) * (y_pos - y1_pos);
597  *x_dist = *x_dist / (x2_pos -x1_pos) / (y2_pos -y1_pos);
598  *y_dist = dy11 * (x2_pos - x_pos) * (y2_pos - y_pos) +
599  dy21 * (x_pos - x1_pos) * (y2_pos - y_pos) +
600  dy12 * (x2_pos - x_pos) * (y_pos - y1_pos) +
601  dy22 * (x_pos - x1_pos) * (y_pos - y1_pos);
602  *y_dist = *y_dist / (x2_pos -x1_pos) / (y2_pos -y1_pos);
603 
604  return 0;
605 
606 }
607 
608 /*----------------------------------------------------------------------------*/
615 /*----------------------------------------------------------------------------*/
617 (cpl_imagelist * ilist,
618  cpl_image ** dist_x,
619  cpl_image ** dist_y)
620 {
621  cpl_image * corr[HAWKI_NB_DETECTORS] ;
622  int i, j ;
623 
624  /* Test entries */
625  if (ilist == NULL) return -1 ;
626  if (dist_x == NULL) return -1 ;
627  if (dist_y == NULL) return -1 ;
628 
629  /* Loop on the 4 chips */
630  for (i=0 ; i<HAWKI_NB_DETECTORS ; i++)
631  {
632  cpl_image * cur_image;
633 
634  /* Get the current image */
635  cur_image = cpl_imagelist_get(ilist, i);
636 
637  /* Correct this image */
638  corr[i] = hawki_distortion_correct_detector(cur_image,dist_x[i],dist_y[i]);
639  if(corr[i] == NULL)
640  {
641  cpl_msg_error(__func__, "Cannot correct the distortion") ;
642  for (j=0 ; j<i ; j++) cpl_image_delete(corr[j]) ;
643  return -1 ;
644  }
645  }
646 
647  /* Store the results */
648  for (i=0 ; i<HAWKI_NB_DETECTORS ; i++)
649  cpl_imagelist_set(ilist, corr[i], i);
650 
651  /* Return */
652  return 0 ;
653 }
654 
657 /*----------------------------------------------------------------------------*/
665 /*----------------------------------------------------------------------------*/
666 int hawki_distortion_create_maps_detector
667 (const hawki_distortion * distortion,
668  cpl_image * dist_x,
669  cpl_image * dist_y)
670 {
671  double * pdx;
672  double * pdy;
673  int nx, ny;
674  int pos;
675  int i, j;
676 
677  /* Test entries */
678  if (distortion == NULL) return -1 ;
679  if (dist_x == NULL) return -1 ;
680  if (dist_y == NULL) return -1 ;
681 
682  /* Initialise */
683  nx = cpl_image_get_size_x(dist_x) ;
684  ny = cpl_image_get_size_y(dist_x) ;
685  if (cpl_image_get_size_x(dist_y) != nx) return -1 ;
686  if (cpl_image_get_size_y(dist_y) != ny) return -1 ;
687 
688  /* Access to the data */
689  pdx = cpl_image_get_data_double(dist_x) ;
690  pdy = cpl_image_get_data_double(dist_y) ;
691 
692  /* Loop on the pixels */
693  for (j=0 ; j<ny ; j++) {
694  for (i=0 ; i<nx ; i++) {
695  double x_dist;
696  double y_dist;
697 
698  pos = i+j*nx ;
699 
701  (distortion, (double)i, (double)j, &x_dist, &y_dist);
702 
703  pdx[pos] = x_dist;
704  pdy[pos] = y_dist;
705  }
706  }
707 
708  return 0 ;
709 }
710 
711 /*----------------------------------------------------------------------------*/
735 /*----------------------------------------------------------------------------*/
736 hawki_distortion * hawki_distortion_compute_solution
737 (const cpl_table ** ref_catalogues,
738  const cpl_bivector * initial_offsets,
739  const cpl_table * matching_sets,
740  int ncats,
741  int detector_nx,
742  int detector_ny,
743  int grid_size,
744  const hawki_distortion * initial_guess,
745  double * rms)
746 {
747 #ifdef HAVE_LIBGSL
748 
749  hawki_distortion * distortion;
750  cpl_bivector * offsets;
751  gsl_multimin_fminimizer * minimizer;
752  gsl_vector * step_size;
753  gsl_vector * init_param;
754  gsl_multimin_function minimize_function;
755 
756  int nfitparam = 0;
757  int iparam;
758  double tolerance = HAWKI_DISTORTION_TOLERANCE;
759  double tolerance2 = HAWKI_DISTORTION_TOLERANCE2;
760  int minimizer_status;
761  int iter = 0;
763 
764  /* Initialize the distortion solution */
765  if(initial_guess == NULL)
766  {
767  distortion = hawki_distortion_grid_new
768  (detector_nx, detector_ny, grid_size);
769  }
770  else
771  {
772  distortion = cpl_malloc(sizeof(hawki_distortion));
773  distortion->dist_x = cpl_image_duplicate(initial_guess->dist_x);
774  distortion->dist_y = cpl_image_duplicate(initial_guess->dist_y);
775  distortion->x_cdelt = initial_guess->x_cdelt;
776  distortion->x_crval = initial_guess->x_crval;
777  distortion->y_cdelt = initial_guess->y_cdelt;
778  distortion->y_crval = initial_guess->y_crval;
779  //We have to fit all the distortion coefficients plus
780  //the offsets of the catalogues
781  }
782  //We have to fit all the distortion coefficients plus
783  //the offsets of the catalogues
784  nfitparam = grid_size * grid_size * 2 + ncats * 2;
785  offsets = cpl_bivector_duplicate(initial_offsets);
786  if(cpl_table_get_nrow(matching_sets) * 2 < nfitparam)
787  {
788  cpl_msg_error(__func__,"Too few matches to compute distortion (< %d)",
789  nfitparam);
790  hawki_distortion_delete(distortion);
791  return NULL;
792  }
793 
794  /* Setup function to minimize */
795  args.ref_catalogues = ref_catalogues;
796  args.matching_sets = matching_sets;
797  args.distortion = distortion;
798  args.offsets = offsets;
799  args.ncats = ncats;
800 
801  minimize_function.f = hawki_distortion_gsl_obj_function;
802  minimize_function.n = nfitparam;
803  minimize_function.params = &args;
804 
805  /* Setup minimizer */
806  minimizer = gsl_multimin_fminimizer_alloc
807  (gsl_multimin_fminimizer_nmsimplex, nfitparam);
808  step_size = gsl_vector_alloc(nfitparam);
809  init_param = gsl_vector_alloc(nfitparam);
810  for(iparam = 0; iparam < grid_size * grid_size * 2; ++iparam)
811  gsl_vector_set(step_size, iparam, 5);
812  for(iparam = grid_size * grid_size * 2;
813  iparam < nfitparam; ++iparam)
814  gsl_vector_set(step_size, iparam, 1);
815  hawki_distortion_update_param_from_solution(init_param, distortion);
816  hawki_distortion_update_param_from_offsets(init_param, offsets);
817  gsl_multimin_fminimizer_set(minimizer, &minimize_function,
818  init_param, step_size);
819 
820  do
821  {
822  ++iter;
823  minimizer_status = gsl_multimin_fminimizer_iterate (minimizer);
824  if(minimizer_status)
825  break;
826  minimizer_status = gsl_multimin_test_size
827  (gsl_multimin_fminimizer_size(minimizer), tolerance);
828  cpl_msg_debug(__func__,"Iteration %d Minimum: %g",
829  iter, gsl_multimin_fminimizer_minimum(minimizer));
830  }
831  while (minimizer_status == GSL_CONTINUE && iter < HAWKI_DISTORTION_MAX_ITER);
832 
833  cpl_msg_warning(__func__, "rms before computing %f", hawki_distortion_compute_rms(ref_catalogues, offsets,
834  matching_sets, ncats, distortion));
835 
836 
837  //Do it again to avoid local minimum
838  gsl_multimin_fminimizer_set(minimizer, &minimize_function,
839  gsl_multimin_fminimizer_x(minimizer), step_size);
840  iter = 0;
841  do
842  {
843  ++iter;
844  minimizer_status = gsl_multimin_fminimizer_iterate (minimizer);
845  if(minimizer_status)
846  break;
847  minimizer_status = gsl_multimin_test_size
848  (gsl_multimin_fminimizer_size(minimizer), tolerance2);
849  cpl_msg_debug(__func__,"2nd run Iteration %d Minimum: %g",
850  iter, gsl_multimin_fminimizer_minimum(minimizer));
851  }
852  while (minimizer_status == GSL_CONTINUE && iter < HAWKI_DISTORTION_MAX_ITER2);
853 
854  /* Update the distortion solution */
855  hawki_distortion_update_solution_from_param
856  (distortion, gsl_multimin_fminimizer_x(minimizer));
857  hawki_distortion_update_offsets_from_param
858  (offsets, gsl_multimin_fminimizer_x(minimizer));
859 
860  *rms = hawki_distortion_compute_rms(ref_catalogues, offsets,
861  matching_sets, ncats, distortion);
862 
863  /* Free and return */
864  gsl_multimin_fminimizer_free(minimizer);
865  gsl_vector_free(init_param);
866  gsl_vector_free(step_size);
867  cpl_bivector_delete(offsets);
868 
869  return distortion;
870 #else
871  cpl_msg_error(__func__,"Not compiled with GSL support.");
872  return NULL;
873 #endif
874 }
875 
876 #ifdef HAVE_LIBGSL
877 /*----------------------------------------------------------------------------*/
885 /*----------------------------------------------------------------------------*/
886 double hawki_distortion_gsl_obj_function
887 (const gsl_vector * dist_param,
888  void * args)
889 {
890  struct _hawki_distortion_obj_function_args_ args_struct;
891  const cpl_table ** ref_catalogues;
892  const cpl_table * matching_sets;
893  hawki_distortion * distortion;
894  cpl_bivector * offsets;
895  int ncats;
896  double rms;
897  double objective_function;
898 
899 
900 
901  args_struct = *(struct _hawki_distortion_obj_function_args_ * )args;
902  ref_catalogues = args_struct.ref_catalogues;
903  matching_sets = args_struct.matching_sets;
904  distortion = args_struct.distortion;
905  offsets = args_struct.offsets;
906  ncats = args_struct.ncats;
907 
908  hawki_distortion_update_solution_from_param(distortion, dist_param);
909  hawki_distortion_update_offsets_from_param(offsets, dist_param);
910 
911  rms = hawki_distortion_compute_rms(ref_catalogues, offsets,
912  matching_sets, ncats, distortion);
913 
914  objective_function = rms;
915 
916 
917  cpl_msg_debug(__func__,"Objective function: %g", objective_function);
918  return objective_function;
919 }
920 #endif
921 
922 /*----------------------------------------------------------------------------*/
928 /*----------------------------------------------------------------------------*/
930 (const cpl_table ** ref_catalogues,
931  const cpl_bivector * cat_offsets,
932  const cpl_table * matching_sets,
933  int ncats,
934  hawki_distortion * distortion)
935 {
936  int icat;
937  int imatch;
938  int nmatch;
939  double rms = 0;
940  const double * x_cat_offsets;
941  const double * y_cat_offsets;
942  const double ** x_cat_cols;
943  const double ** y_cat_cols;
944  const cpl_array ** match_arrays;
945  double ** x_pos_values;
946  double ** y_pos_values;
947  int ** pos_flag;
948 
949 
950 
951  nmatch = cpl_table_get_nrow(matching_sets);
952 
953  x_cat_offsets = cpl_vector_get_data_const
954  (cpl_bivector_get_x_const(cat_offsets));
955  y_cat_offsets = cpl_vector_get_data_const
956  (cpl_bivector_get_y_const(cat_offsets));
957 
958  x_cat_cols = cpl_malloc(sizeof(double *) * ncats);
959  y_cat_cols = cpl_malloc(sizeof(double *) * ncats);
960  for(icat = 0; icat < ncats; ++icat)
961  {
962  x_cat_cols[icat] = cpl_table_get_data_double_const(ref_catalogues[icat],
963  HAWKI_COL_OBJ_POSX);
964  y_cat_cols[icat] = cpl_table_get_data_double_const(ref_catalogues[icat],
965  HAWKI_COL_OBJ_POSY);
966  }
967 
968  match_arrays = cpl_malloc(sizeof(cpl_array *) * nmatch);
969  x_pos_values = cpl_malloc(sizeof(double *) * nmatch);
970  y_pos_values = cpl_malloc(sizeof(double *) * nmatch);
971  pos_flag = cpl_malloc(sizeof(int *) * nmatch);
972  for(imatch = 0; imatch < nmatch; ++imatch)
973  {
974  match_arrays[imatch] = cpl_table_get_array(matching_sets,
975  HAWKI_COL_MATCHING_SETS, imatch);
976  x_pos_values[imatch] = cpl_malloc(sizeof(double) * ncats);
977  y_pos_values[imatch] = cpl_malloc(sizeof(double) * ncats);
978  pos_flag[imatch] = cpl_malloc(sizeof(int) * ncats);
979  }
980 
981 #ifdef _OPENMP
982 #pragma omp parallel for private(icat,imatch) reduction(+:rms)
983 #endif
984  for(imatch = 0; imatch < nmatch; ++imatch)
985  {
986  int nstddev;
987  double var_x;
988  double var_y;
989 
990  for(icat = 0; icat < ncats; ++icat)
991  {
992  int iobj;
993  double x_cat_offset;
994  double y_cat_offset;
995 
996  x_cat_offset = x_cat_offsets[icat];
997  y_cat_offset = y_cat_offsets[icat];
998 
999  if((iobj = cpl_array_get(match_arrays[imatch], icat, NULL)) != -1)
1000  {
1001  double x_cat;
1002  double y_cat;
1003  double x_dist_corr;
1004  double y_dist_corr;
1005  double x_dist;
1006  double y_dist;
1007  double x_glob;
1008  double y_glob;
1009 
1010  x_cat = x_cat_cols[icat][iobj];
1011  y_cat = y_cat_cols[icat][iobj];
1012 
1013 
1014  //These 4 lines of code are from hawki_distortion_correct_coords.
1015  //The are repeated here to avoid a cpl call, which is not thread-safe
1016  //Two checks to ensure thread-safety:
1017  //We have ensured outside the loop that distortion->dist_x and
1018  //distortion->dist_y are not null.
1019  //We have checked outside the loop the mask has all the points valid
1021  (distortion, x_cat, y_cat, &x_dist, &y_dist);
1022  x_dist_corr = x_cat - x_dist;
1023  y_dist_corr = y_cat - y_dist;
1024 
1025 
1026  x_glob = x_dist_corr + x_cat_offset;
1027  y_glob = y_dist_corr + y_cat_offset;
1028  x_pos_values[imatch][icat] = x_glob;
1029  y_pos_values[imatch][icat] = y_glob;
1030  pos_flag[imatch][icat] = 1;
1031  }
1032  else
1033  pos_flag[imatch][icat] = 0;
1034  }
1035 
1036  hawki_distortion_get_flag_vars(x_pos_values[imatch],
1037  y_pos_values[imatch], pos_flag[imatch],
1038  ncats, &nstddev, &var_x, &var_y);
1039 
1040  //The rms is counted as many times as this star is the list of catalogs.
1041  rms += sqrt(var_x + var_y) * nstddev;
1042 
1043  }
1044  cpl_free(x_cat_cols);
1045  cpl_free(y_cat_cols);
1046  for(imatch = 0; imatch < nmatch; ++imatch)
1047  {
1048  cpl_free(x_pos_values[imatch]);
1049  cpl_free(y_pos_values[imatch]);
1050  cpl_free(pos_flag[imatch]);
1051  }
1052  cpl_free(x_pos_values);
1053  cpl_free(y_pos_values);
1054  cpl_free(pos_flag);
1055  cpl_free(match_arrays);
1056 
1057  return rms;
1058 }
1059 
1060 #ifdef HAVE_LIBGSL
1061 /*----------------------------------------------------------------------------*/
1067 /*----------------------------------------------------------------------------*/
1068 int hawki_distortion_update_solution_from_param
1069 (hawki_distortion * distortion,
1070  const gsl_vector * dist_param)
1071 {
1072  int ipoint;
1073  int ix;
1074  int iy;
1075  int nx;
1076  int ny;
1077  double x_dist_mean;
1078  double y_dist_mean;
1079 
1080  nx = cpl_image_get_size_x(distortion->dist_x);
1081  ny = cpl_image_get_size_y(distortion->dist_x);
1082  for(ix = 0; ix < nx; ++ix)
1083  for(iy = 0; iy < ny; ++iy)
1084  {
1085  ipoint = ix + iy * nx;
1086  cpl_image_set(distortion->dist_x, ix+1, iy+1,
1087  gsl_vector_get(dist_param, ipoint * 2));
1088  cpl_image_set(distortion->dist_y, ix+1, iy+1,
1089  gsl_vector_get(dist_param, ipoint * 2 + 1));
1090  }
1091 
1092  /* Normalize to mean(distorsion) = 0 */
1093  x_dist_mean = cpl_image_get_mean(distortion->dist_x);
1094  y_dist_mean = cpl_image_get_mean(distortion->dist_y);
1095  cpl_image_subtract_scalar(distortion->dist_x, x_dist_mean);
1096  cpl_image_subtract_scalar(distortion->dist_y, y_dist_mean);
1097 
1098  return 0;
1099 }
1100 
1101 /*----------------------------------------------------------------------------*/
1107 /*----------------------------------------------------------------------------*/
1108 int hawki_distortion_update_offsets_from_param
1109 (cpl_bivector * offsets,
1110  const gsl_vector * dist_param)
1111 {
1112  int i;
1113  int ncats;
1114  int nparam;
1115 
1116  ncats = cpl_bivector_get_size(offsets);
1117  nparam = dist_param->size;
1118  for(i = 0; i < ncats; ++i)
1119  {
1120  cpl_vector_set(cpl_bivector_get_x(offsets), i,
1121  gsl_vector_get(dist_param, nparam - 2 * ncats + 2 * i));
1122  cpl_vector_set(cpl_bivector_get_y(offsets), i,
1123  gsl_vector_get(dist_param, nparam - 2 * ncats + 2 * i + 1));
1124  }
1125 
1126  return 0;
1127 }
1128 
1129 /*----------------------------------------------------------------------------*/
1135 /*----------------------------------------------------------------------------*/
1136 int hawki_distortion_update_param_from_solution
1137 (gsl_vector * dist_param,
1138  const hawki_distortion * distortion)
1139 {
1140  int ipoint;
1141  int ix;
1142  int iy;
1143  int nx;
1144  int ny;
1145  int isnull;
1146 
1147  nx = cpl_image_get_size_x(distortion->dist_x);
1148  ny = cpl_image_get_size_y(distortion->dist_y);
1149  for(ix = 0; ix < nx; ++ix)
1150  for(iy = 0; iy < ny; ++iy)
1151  {
1152  ipoint = ix + iy * nx;
1153  gsl_vector_set(dist_param, ipoint * 2,
1154  cpl_image_get(distortion->dist_x, ix+1, iy+1, &isnull));
1155  gsl_vector_set(dist_param, ipoint * 2 + 1,
1156  cpl_image_get(distortion->dist_y, ix+1, iy+1, &isnull));
1157  }
1158 
1159  return 0;
1160 }
1161 
1162 /*----------------------------------------------------------------------------*/
1167 /*----------------------------------------------------------------------------*/
1168 int hawki_distortion_update_param_from_offsets
1169 (gsl_vector * dist_param,
1170  const cpl_bivector * offsets)
1171 {
1172  int i;
1173  int ncats;
1174  int nparam;
1175 
1176  ncats = cpl_bivector_get_size(offsets);
1177  nparam = dist_param->size;
1178  for(i = 0; i < ncats; ++i)
1179  {
1180  gsl_vector_set(dist_param, nparam - 2 * ncats + 2 * i,
1181  cpl_vector_get(cpl_bivector_get_x_const(offsets), i));
1182  gsl_vector_set(dist_param, nparam - 2 * ncats + 2 * i + 1,
1183  cpl_vector_get(cpl_bivector_get_y_const(offsets), i));
1184  }
1185 
1186  return 0;
1187 }
1188 #endif
1189 
1190 /*----------------------------------------------------------------------------*/
1196 /*----------------------------------------------------------------------------*/
1198 (double * x_val, double * y_val, int *pos_flag,
1199  int nvals, int *nvalid, double *var_x, double *var_y)
1200 {
1201  double varsum_x = 0.0;
1202  double varsum_y = 0.0;
1203  double mean_x = 0.0;
1204  double mean_y = 0.0;
1205  int i;
1206 
1207  *nvalid = 0;
1208  for (i=0; i < nvals; i++)
1209  {
1210  if(pos_flag[i] == 1)
1211  {
1212  const double delta_x = (double)x_val[i] - mean_x;
1213  const double delta_y = (double)y_val[i] - mean_y;
1214 
1215  varsum_x += *nvalid * delta_x * delta_x / (*nvalid + 1.0);
1216  varsum_y += *nvalid * delta_y * delta_y / (*nvalid + 1.0);
1217  mean_x += delta_x / (*nvalid + 1.0);
1218  mean_y += delta_y / (*nvalid + 1.0);
1219  (*nvalid)++;
1220  }
1221  }
1222 
1223  /* Compute the bias-corrected standard deviation.
1224  - With the recurrence relation rounding can likely not cause
1225  the variance to become negative, but check just to be safe */
1226  *var_x = varsum_x / (double) (*nvalid - 1);
1227  *var_y = varsum_y / (double) (*nvalid - 1);
1228 }