00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 #ifdef HAVE_CONFIG_H
00029 #include <config.h>
00030 #endif
00031
00032
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
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
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
00138 distortion = cpl_malloc(sizeof(hawki_distortion));
00139
00140
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
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
00192 distortion = cpl_malloc(sizeof(hawki_distortion));
00193
00194
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
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
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
00297 if (ilist == NULL) return -1 ;
00298 if (frame_dist_x == NULL) return -1 ;
00299 if (frame_dist_y == NULL) return -1 ;
00300
00301
00302 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
00303 {
00304 int nx;
00305 int ny;
00306
00307
00308 nx = cpl_image_get_size_x(ilist[idet]);
00309 ny = cpl_image_get_size_y(ilist[idet]);
00310
00311
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
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
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
00351 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
00352 {
00353 cpl_image_delete(ilist[idet]) ;
00354 ilist[idet] = corr[idet] ;
00355 }
00356
00357
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
00378 if (image == NULL) return NULL;
00379 if (dist_x == NULL) return NULL;
00380 if (dist_y == NULL) return NULL;
00381
00382
00383 corr = cpl_image_new(cpl_image_get_size_x(image),
00384 cpl_image_get_size_y(image), CPL_TYPE_FLOAT) ;
00385
00386
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
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
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
00490 *x_pos_distinvcorr = x_pos + x_dist;
00491 *y_pos_distinvcorr = y_pos + y_dist;
00492
00493 return 0;
00494 }
00495
00496
00510
00511 int hawki_distortion_interpolate_distortion
00512 (const hawki_distortion * distortion,
00513 double x_pos,
00514 double y_pos,
00515 double * x_dist,
00516 double * y_dist)
00517 {
00518 int ix1;
00519 int ix2;
00520 int iy1;
00521 int iy2;
00522 int nx;
00523 int ny;
00524 double x1_pos;
00525 double x2_pos;
00526 double y1_pos;
00527 double y2_pos;
00528 double dx11;
00529 double dx12;
00530 double dx21;
00531 double dx22;
00532 double dy11;
00533 double dy12;
00534 double dy21;
00535 double dy22;
00536 int isnull;
00537
00538
00539 nx = cpl_image_get_size_x(distortion->dist_x);
00540 ny = cpl_image_get_size_y(distortion->dist_x);
00541
00542
00543
00544 ix1 = (int)floor((x_pos - distortion->x_crval)/distortion->x_cdelt);
00545 iy1 = (int)floor((y_pos - distortion->y_crval)/distortion->y_cdelt);
00546
00547 if(ix1 < 0)
00548 ix1 = 0;
00549 if(ix1 >= nx - 1)
00550 ix1 = nx - 2;
00551 if(iy1 < 0)
00552 iy1 = 0;
00553 if(iy1 >= ny - 1)
00554 iy1 = ny - 2;
00555
00556 ix2 = ix1 + 1;
00557 iy2 = iy1 + 1;
00558
00559
00560
00561 x1_pos = ix1 * distortion->x_cdelt + distortion->x_crval;
00562 x2_pos = ix2 * distortion->x_cdelt + distortion->x_crval;
00563 y1_pos = iy1 * distortion->y_cdelt + distortion->y_crval;
00564 y2_pos = iy2 * distortion->y_cdelt + distortion->y_crval;
00565
00566
00567
00568 dx11 = cpl_image_get(distortion->dist_x, ix1 + 1, iy1 + 1, &isnull);
00569 dx21 = cpl_image_get(distortion->dist_x, ix2 + 1, iy1 + 1, &isnull);
00570 dx12 = cpl_image_get(distortion->dist_x, ix1 + 1, iy2 + 1, &isnull);
00571 dx22 = cpl_image_get(distortion->dist_x, ix2 + 1, iy2 + 1, &isnull);
00572 dy11 = cpl_image_get(distortion->dist_y, ix1 + 1, iy1 + 1, &isnull);
00573 dy21 = cpl_image_get(distortion->dist_y, ix2 + 1, iy1 + 1, &isnull);
00574 dy12 = cpl_image_get(distortion->dist_y, ix1 + 1, iy2 + 1, &isnull);
00575 dy22 = cpl_image_get(distortion->dist_y, ix2 + 1, iy2 + 1, &isnull);
00576
00577
00578
00579 *x_dist = dx11 * (x2_pos - x_pos) * (y2_pos - y_pos) +
00580 dx21 * (x_pos - x1_pos) * (y2_pos - y_pos) +
00581 dx12 * (x2_pos - x_pos) * (y_pos - y1_pos) +
00582 dx22 * (x_pos - x1_pos) * (y_pos - y1_pos);
00583 *x_dist = *x_dist / (x2_pos -x1_pos) / (y2_pos -y1_pos);
00584 *y_dist = dy11 * (x2_pos - x_pos) * (y2_pos - y_pos) +
00585 dy21 * (x_pos - x1_pos) * (y2_pos - y_pos) +
00586 dy12 * (x2_pos - x_pos) * (y_pos - y1_pos) +
00587 dy22 * (x_pos - x1_pos) * (y_pos - y1_pos);
00588 *y_dist = *y_dist / (x2_pos -x1_pos) / (y2_pos -y1_pos);
00589
00590 return 0;
00591
00592 }
00593
00594
00601
00602 int hawki_distortion_apply_maps
00603 (cpl_imagelist * ilist,
00604 cpl_image ** dist_x,
00605 cpl_image ** dist_y)
00606 {
00607 cpl_image * corr[HAWKI_NB_DETECTORS] ;
00608 int i, j ;
00609
00610
00611 if (ilist == NULL) return -1 ;
00612 if (dist_x == NULL) return -1 ;
00613 if (dist_y == NULL) return -1 ;
00614
00615
00616 for (i=0 ; i<HAWKI_NB_DETECTORS ; i++)
00617 {
00618 cpl_image * cur_image;
00619
00620
00621 cur_image = cpl_imagelist_get(ilist, i);
00622
00623
00624 corr[i] = hawki_distortion_correct_detector(cur_image,dist_x[i],dist_y[i]);
00625 if(corr[i] == NULL)
00626 {
00627 cpl_msg_error(__func__, "Cannot correct the distortion") ;
00628 for (j=0 ; j<i ; j++) cpl_image_delete(corr[j]) ;
00629 return -1 ;
00630 }
00631 }
00632
00633
00634 for (i=0 ; i<HAWKI_NB_DETECTORS ; i++)
00635 cpl_imagelist_set(ilist, corr[i], i);
00636
00637
00638 return 0 ;
00639 }
00640
00643
00651
00652 int hawki_distortion_create_maps_detector
00653 (const hawki_distortion * distortion,
00654 cpl_image * dist_x,
00655 cpl_image * dist_y)
00656 {
00657 double * pdx;
00658 double * pdy;
00659 int nx, ny;
00660 int pos;
00661 int i, j;
00662
00663
00664 if (distortion == NULL) return -1 ;
00665 if (dist_x == NULL) return -1 ;
00666 if (dist_y == NULL) return -1 ;
00667
00668
00669 nx = cpl_image_get_size_x(dist_x) ;
00670 ny = cpl_image_get_size_y(dist_x) ;
00671 if (cpl_image_get_size_x(dist_y) != nx) return -1 ;
00672 if (cpl_image_get_size_y(dist_y) != ny) return -1 ;
00673
00674
00675 pdx = cpl_image_get_data_double(dist_x) ;
00676 pdy = cpl_image_get_data_double(dist_y) ;
00677
00678
00679 for (j=0 ; j<ny ; j++) {
00680 for (i=0 ; i<nx ; i++) {
00681 double x_dist;
00682 double y_dist;
00683
00684 pos = i+j*nx ;
00685
00686 hawki_distortion_interpolate_distortion
00687 (distortion, (double)i, (double)j, &x_dist, &y_dist);
00688
00689 pdx[pos] = x_dist;
00690 pdy[pos] = y_dist;
00691 }
00692 }
00693
00694 return 0 ;
00695 }
00696
00697
00721
00722 hawki_distortion * hawki_distortion_compute_solution
00723 (const cpl_table ** ref_catalogues,
00724 const cpl_bivector * initial_offsets,
00725 const cpl_table * matching_sets,
00726 int ncats,
00727 int detector_nx,
00728 int detector_ny,
00729 int grid_size,
00730 const hawki_distortion * initial_guess,
00731 double * rms)
00732 {
00733 #ifdef HAVE_LIBGSL
00734
00735 hawki_distortion * distortion;
00736 cpl_bivector * offsets;
00737 gsl_multimin_fminimizer * minimizer;
00738 gsl_vector * step_size;
00739 gsl_vector * init_param;
00740 gsl_multimin_function minimize_function;
00741
00742 int nfitparam = 0;
00743 int iparam;
00744 double tolerance = HAWKI_DISTORTION_TOLERANCE;
00745 double tolerance2 = HAWKI_DISTORTION_TOLERANCE2;
00746 int minimizer_status;
00747 int iter = 0;
00748 struct _hawki_distortion_obj_function_args_ args;
00749
00750
00751 if(initial_guess == NULL)
00752 {
00753 distortion = hawki_distortion_grid_new
00754 (detector_nx, detector_ny, grid_size);
00755 }
00756 else
00757 {
00758 distortion = cpl_malloc(sizeof(hawki_distortion));
00759 distortion->dist_x = cpl_image_duplicate(initial_guess->dist_x);
00760 distortion->dist_y = cpl_image_duplicate(initial_guess->dist_y);
00761 distortion->x_cdelt = initial_guess->x_cdelt;
00762 distortion->x_crval = initial_guess->x_crval;
00763 distortion->y_cdelt = initial_guess->y_cdelt;
00764 distortion->y_crval = initial_guess->y_crval;
00765
00766
00767 }
00768
00769
00770 nfitparam = grid_size * grid_size * 2 + ncats * 2;
00771 offsets = cpl_bivector_duplicate(initial_offsets);
00772 if(cpl_table_get_nrow(matching_sets) * 2 < nfitparam)
00773 {
00774 cpl_msg_error(__func__,"Too few matches to compute distortion (< %d)",
00775 nfitparam);
00776 hawki_distortion_delete(distortion);
00777 return NULL;
00778 }
00779
00780
00781 args.ref_catalogues = ref_catalogues;
00782 args.matching_sets = matching_sets;
00783 args.distortion = distortion;
00784 args.offsets = offsets;
00785 args.ncats = ncats;
00786
00787 minimize_function.f = hawki_distortion_gsl_obj_function;
00788 minimize_function.n = nfitparam;
00789 minimize_function.params = &args;
00790
00791
00792 minimizer = gsl_multimin_fminimizer_alloc
00793 (gsl_multimin_fminimizer_nmsimplex, nfitparam);
00794 step_size = gsl_vector_alloc(nfitparam);
00795 init_param = gsl_vector_alloc(nfitparam);
00796 for(iparam = 0; iparam < grid_size * grid_size * 2; ++iparam)
00797 gsl_vector_set(step_size, iparam, 5);
00798 for(iparam = grid_size * grid_size * 2;
00799 iparam < nfitparam; ++iparam)
00800 gsl_vector_set(step_size, iparam, 1);
00801 hawki_distortion_update_param_from_solution(init_param, distortion);
00802 hawki_distortion_update_param_from_offsets(init_param, offsets);
00803 gsl_multimin_fminimizer_set(minimizer, &minimize_function,
00804 init_param, step_size);
00805
00806 do
00807 {
00808 ++iter;
00809 minimizer_status = gsl_multimin_fminimizer_iterate (minimizer);
00810 if(minimizer_status)
00811 break;
00812 minimizer_status = gsl_multimin_test_size
00813 (gsl_multimin_fminimizer_size(minimizer), tolerance);
00814 cpl_msg_debug(__func__,"Iteration %d Minimum: %g",
00815 iter, gsl_multimin_fminimizer_minimum(minimizer));
00816 }
00817 while (minimizer_status == GSL_CONTINUE && iter < HAWKI_DISTORTION_MAX_ITER);
00818
00819 cpl_msg_warning(__func__, "rms before computing %f", hawki_distortion_compute_rms(ref_catalogues, offsets,
00820 matching_sets, ncats, distortion));
00821
00822
00823
00824 gsl_multimin_fminimizer_set(minimizer, &minimize_function,
00825 gsl_multimin_fminimizer_x(minimizer), step_size);
00826 iter = 0;
00827 do
00828 {
00829 ++iter;
00830 minimizer_status = gsl_multimin_fminimizer_iterate (minimizer);
00831 if(minimizer_status)
00832 break;
00833 minimizer_status = gsl_multimin_test_size
00834 (gsl_multimin_fminimizer_size(minimizer), tolerance2);
00835 cpl_msg_debug(__func__,"2nd run Iteration %d Minimum: %g",
00836 iter, gsl_multimin_fminimizer_minimum(minimizer));
00837 }
00838 while (minimizer_status == GSL_CONTINUE && iter < HAWKI_DISTORTION_MAX_ITER2);
00839
00840
00841 hawki_distortion_update_solution_from_param
00842 (distortion, gsl_multimin_fminimizer_x(minimizer));
00843 hawki_distortion_update_offsets_from_param
00844 (offsets, gsl_multimin_fminimizer_x(minimizer));
00845
00846 *rms = hawki_distortion_compute_rms(ref_catalogues, offsets,
00847 matching_sets, ncats, distortion);
00848
00849
00850 gsl_multimin_fminimizer_free(minimizer);
00851 gsl_vector_free(init_param);
00852 gsl_vector_free(step_size);
00853 cpl_bivector_delete(offsets);
00854
00855 return distortion;
00856 #else
00857 cpl_msg_error(__func__,"Not compiled with GSL support.");
00858 return NULL;
00859 #endif
00860 }
00861
00862 #ifdef HAVE_LIBGSL
00863
00871
00872 double hawki_distortion_gsl_obj_function
00873 (const gsl_vector * dist_param,
00874 void * args)
00875 {
00876 struct _hawki_distortion_obj_function_args_ args_struct;
00877 const cpl_table ** ref_catalogues;
00878 const cpl_table * matching_sets;
00879 hawki_distortion * distortion;
00880 cpl_bivector * offsets;
00881 int ncats;
00882 double rms;
00883 double objective_function;
00884
00885
00886
00887 args_struct = *(struct _hawki_distortion_obj_function_args_ * )args;
00888 ref_catalogues = args_struct.ref_catalogues;
00889 matching_sets = args_struct.matching_sets;
00890 distortion = args_struct.distortion;
00891 offsets = args_struct.offsets;
00892 ncats = args_struct.ncats;
00893
00894 hawki_distortion_update_solution_from_param(distortion, dist_param);
00895 hawki_distortion_update_offsets_from_param(offsets, dist_param);
00896
00897 rms = hawki_distortion_compute_rms(ref_catalogues, offsets,
00898 matching_sets, ncats, distortion);
00899
00900 objective_function = rms;
00901
00902
00903 cpl_msg_debug(__func__,"Objective function: %g", objective_function);
00904 return objective_function;
00905 }
00906 #endif
00907
00908
00914
00915 double hawki_distortion_compute_rms
00916 (const cpl_table ** ref_catalogues,
00917 const cpl_bivector * cat_offsets,
00918 const cpl_table * matching_sets,
00919 int ncats,
00920 hawki_distortion * distortion)
00921 {
00922 int imatch;
00923 int nmatch;
00924 double rms = 0;
00925
00926
00927 nmatch = cpl_table_get_nrow(matching_sets);
00928 for(imatch = 0; imatch < nmatch; ++imatch)
00929 {
00930 cx_deque * x_pos_values;
00931 cx_deque * y_pos_values;
00932 const cpl_array * match;
00933 int icat;
00934 double rms_x;
00935 double rms_y;
00936
00937 match = cpl_table_get_array(matching_sets,
00938 HAWKI_COL_MATCHING_SETS, imatch);
00939
00940 x_pos_values = cx_deque_new();
00941 y_pos_values = cx_deque_new();
00942
00943 cpl_msg_debug(__func__,"Computing the stddev on match %d", imatch);
00944
00945 for(icat = 0; icat < ncats; ++icat)
00946 {
00947 int iobj;
00948 double x_cat_offset;
00949 double y_cat_offset;
00950
00951 x_cat_offset = cpl_vector_get_data_const
00952 (cpl_bivector_get_x_const(cat_offsets))[icat];
00953 y_cat_offset = cpl_vector_get_data_const
00954 (cpl_bivector_get_y_const(cat_offsets))[icat];
00955
00956 if((iobj = cpl_array_get(match, icat, NULL)) != -1)
00957 {
00958 double x_cat;
00959 double y_cat;
00960 double x_dist_corr;
00961 double y_dist_corr;
00962 double * x_glob = cx_malloc(sizeof(double));
00963 double * y_glob = cx_malloc(sizeof(double));
00964
00965 x_cat = cpl_table_get_double(ref_catalogues[icat],
00966 HAWKI_COL_OBJ_POSX, iobj, NULL);
00967 y_cat = cpl_table_get_double(ref_catalogues[icat],
00968 HAWKI_COL_OBJ_POSY, iobj, NULL);
00969
00970 hawki_distortion_correct_coords
00971 (distortion,
00972 x_cat, y_cat,
00973 &x_dist_corr, &y_dist_corr);
00974 *x_glob = x_dist_corr + x_cat_offset;
00975 *y_glob = y_dist_corr + y_cat_offset;
00976 cx_deque_push_back(x_pos_values, (void *)x_glob);
00977 cx_deque_push_back(y_pos_values, (void *)y_glob);
00978
00979 cpl_msg_debug(__func__,"In cat %d the pos is %f %f. Glob %f %f",
00980 icat, x_cat, y_cat, *x_glob, *y_glob);
00981 }
00982 }
00983
00984 rms_x = hawki_distortion_get_deque_stddev(x_pos_values);
00985 rms_y = hawki_distortion_get_deque_stddev(y_pos_values);
00986
00987 rms += sqrt(rms_x * rms_x + rms_y * rms_y) * cx_deque_size(x_pos_values);
00988
00989 cx_deque_destroy(x_pos_values, cx_free);
00990 cx_deque_destroy(y_pos_values, cx_free);
00991 }
00992
00993 return rms;
00994 }
00995
00996 #ifdef HAVE_LIBGSL
00997
01003
01004 int hawki_distortion_update_solution_from_param
01005 (hawki_distortion * distortion,
01006 const gsl_vector * dist_param)
01007 {
01008 int ipoint;
01009 int ix;
01010 int iy;
01011 int nx;
01012 int ny;
01013 double x_dist_mean;
01014 double y_dist_mean;
01015
01016 nx = cpl_image_get_size_x(distortion->dist_x);
01017 ny = cpl_image_get_size_y(distortion->dist_x);
01018 for(ix = 0; ix < nx; ++ix)
01019 for(iy = 0; iy < ny; ++iy)
01020 {
01021 ipoint = ix + iy * nx;
01022 cpl_image_set(distortion->dist_x, ix+1, iy+1,
01023 gsl_vector_get(dist_param, ipoint * 2));
01024 cpl_image_set(distortion->dist_y, ix+1, iy+1,
01025 gsl_vector_get(dist_param, ipoint * 2 + 1));
01026 }
01027
01028
01029 x_dist_mean = cpl_image_get_mean(distortion->dist_x);
01030 y_dist_mean = cpl_image_get_mean(distortion->dist_y);
01031 cpl_image_subtract_scalar(distortion->dist_x, x_dist_mean);
01032 cpl_image_subtract_scalar(distortion->dist_y, y_dist_mean);
01033
01034 return 0;
01035 }
01036
01037
01043
01044 int hawki_distortion_update_offsets_from_param
01045 (cpl_bivector * offsets,
01046 const gsl_vector * dist_param)
01047 {
01048 int i;
01049 int ncats;
01050 int nparam;
01051
01052 ncats = cpl_bivector_get_size(offsets);
01053 nparam = dist_param->size;
01054 for(i = 0; i < ncats; ++i)
01055 {
01056 cpl_vector_set(cpl_bivector_get_x(offsets), i,
01057 gsl_vector_get(dist_param, nparam - 2 * ncats + 2 * i));
01058 cpl_vector_set(cpl_bivector_get_y(offsets), i,
01059 gsl_vector_get(dist_param, nparam - 2 * ncats + 2 * i + 1));
01060 }
01061
01062 return 0;
01063 }
01064
01065
01071
01072 int hawki_distortion_update_param_from_solution
01073 (gsl_vector * dist_param,
01074 const hawki_distortion * distortion)
01075 {
01076 int ipoint;
01077 int ix;
01078 int iy;
01079 int nx;
01080 int ny;
01081 int isnull;
01082
01083 nx = cpl_image_get_size_x(distortion->dist_x);
01084 ny = cpl_image_get_size_y(distortion->dist_y);
01085 for(ix = 0; ix < nx; ++ix)
01086 for(iy = 0; iy < ny; ++iy)
01087 {
01088 ipoint = ix + iy * nx;
01089 gsl_vector_set(dist_param, ipoint * 2,
01090 cpl_image_get(distortion->dist_x, ix+1, iy+1, &isnull));
01091 gsl_vector_set(dist_param, ipoint * 2 + 1,
01092 cpl_image_get(distortion->dist_y, ix+1, iy+1, &isnull));
01093 }
01094
01095 return 0;
01096 }
01097
01098
01103
01104 int hawki_distortion_update_param_from_offsets
01105 (gsl_vector * dist_param,
01106 const cpl_bivector * offsets)
01107 {
01108 int i;
01109 int ncats;
01110 int nparam;
01111
01112 ncats = cpl_bivector_get_size(offsets);
01113 nparam = dist_param->size;
01114 for(i = 0; i < ncats; ++i)
01115 {
01116 gsl_vector_set(dist_param, nparam - 2 * ncats + 2 * i,
01117 cpl_vector_get(cpl_bivector_get_x_const(offsets), i));
01118 gsl_vector_set(dist_param, nparam - 2 * ncats + 2 * i + 1,
01119 cpl_vector_get(cpl_bivector_get_y_const(offsets), i));
01120 }
01121
01122 return 0;
01123 }
01124 #endif
01125
01126
01131
01132 double hawki_distortion_get_deque_stddev
01133 (cx_deque * deque)
01134 {
01135
01136 double varsum = 0.0;
01137 double mean = 0.0;
01138 size_t i;
01139 cx_deque_iterator it;
01140
01141 cpl_ensure(deque != NULL, CPL_ERROR_NULL_INPUT,-1);
01142 cpl_ensure(cx_deque_size(deque) > 1, CPL_ERROR_ILLEGAL_INPUT,-2);
01143
01144 for (i=0, it = cx_deque_begin(deque); i < cx_deque_size(deque);
01145 i++, it = cx_deque_next(deque, it)) {
01146 const double delta = *(double *)cx_deque_get(deque, it) - mean;
01147
01148 varsum += i * delta * delta / (i + 1.0);
01149 mean += delta / (i + 1.0);
01150 }
01151
01152 return sqrt(varsum / (double) (cx_deque_size(deque)-1));
01153
01154 }