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
46 #include <gsl/gsl_multimin.h>
49 #include "hawki_distortion.h"
50 #include "hawki_dfs.h"
51 #include "hawki_utils.h"
52 #include "hawki_load.h"
71 const cpl_table ** ref_catalogues;
72 const cpl_table * matching_sets;
73 cpl_bivector * offsets;
74 hawki_distortion * distortion;
82 (
const hawki_distortion * distortion,
89 (
const cpl_table ** ref_catalogues,
90 const cpl_bivector * cat_offsets,
91 const cpl_table * matching_sets,
93 hawki_distortion * distortion);
96 double hawki_distortion_gsl_obj_function
97 (
const gsl_vector * dist_param,
100 int hawki_distortion_update_solution_from_param
101 (hawki_distortion * distortion,
102 const gsl_vector * dist_param);
104 int hawki_distortion_update_offsets_from_param
105 (cpl_bivector * offsets,
106 const gsl_vector * dist_param);
108 int hawki_distortion_update_param_from_solution
109 (gsl_vector * dist_param,
110 const hawki_distortion * distortion);
112 int hawki_distortion_update_param_from_offsets
113 (gsl_vector * dist_param,
114 const cpl_bivector * offsets);
118 (
double * x_val,
double * y_val,
int *pos_flag,
119 int nvals,
int *nvalid,
double *var_x,
double *var_y);
136 hawki_distortion * distortion;
139 distortion = cpl_malloc(
sizeof(hawki_distortion));
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);
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;
163 (hawki_distortion * distortion)
165 if(distortion == NULL)
167 cpl_image_delete(distortion->dist_x);
168 cpl_image_delete(distortion->dist_y);
169 cpl_free(distortion);
182 (
const cpl_frame * dist_x,
183 const cpl_frame * dist_y,
186 const char * file_dist_x;
187 const char * file_dist_y;
188 hawki_distortion * distortion;
190 cpl_propertylist * plist;
193 distortion = cpl_malloc(
sizeof(hawki_distortion));
196 file_dist_x = cpl_frame_get_filename(dist_x);
197 file_dist_y = cpl_frame_get_filename(dist_y);
199 (dist_x, idet, CPL_TYPE_FLOAT);
201 (dist_y, idet, CPL_TYPE_FLOAT);
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)
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);
221 cpl_propertylist_delete(plist);
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)
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);
239 cpl_propertylist_delete(plist);
252 (
const hawki_distortion * distortion)
254 if(distortion == NULL)
256 cpl_error_set(__func__,CPL_ERROR_ILLEGAL_INPUT);
258 return cpl_image_get_size_x(distortion->dist_x);
269 (
const hawki_distortion * distortion)
271 if(distortion == NULL)
273 cpl_error_set(__func__,CPL_ERROR_ILLEGAL_INPUT);
275 return cpl_image_get_size_y(distortion->dist_x);
288 const cpl_frame * frame_dist_x,
289 const cpl_frame * frame_dist_y)
291 cpl_image * corr[HAWKI_NB_DETECTORS];
292 hawki_distortion * distortion;
298 if (ilist == NULL)
return -1 ;
299 if (frame_dist_x == NULL)
return -1 ;
300 if (frame_dist_y == NULL)
return -1 ;
303 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
309 nx = cpl_image_get_size_x(ilist[idet]);
310 ny = cpl_image_get_size_y(ilist[idet]);
314 (frame_dist_x, frame_dist_y, idet + 1)) == NULL)
316 cpl_msg_error(__func__,
"Cannot load the distortion for chip %d",
318 for (j=0 ; j<idet ; j++) cpl_image_delete(corr[j]) ;
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))
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]) ;
337 if(corr[idet] == NULL)
339 cpl_msg_error(__func__,
"Cannot correct the 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]) ;
347 cpl_image_delete(dist_x) ;
348 cpl_image_delete(dist_y);
352 for (idet=0 ; idet<HAWKI_NB_DETECTORS ; idet++)
354 cpl_image_delete(ilist[idet]) ;
355 ilist[idet] = corr[idet] ;
376 cpl_vector * profile ;
379 if (image == NULL)
return NULL;
380 if (dist_x == NULL)
return NULL;
381 if (dist_y == NULL)
return NULL;
384 corr = cpl_image_new(cpl_image_get_size_x(image),
385 cpl_image_get_size_y(image), CPL_TYPE_FLOAT) ;
388 profile = cpl_vector_new(CPL_KERNEL_DEF_SAMPLES) ;
389 cpl_vector_fill_kernel_profile(profile, CPL_KERNEL_DEFAULT,
390 CPL_KERNEL_DEF_WIDTH) ;
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)
397 cpl_msg_error(__func__,
"Cannot warp the image") ;
398 cpl_image_delete(corr) ;
399 cpl_vector_delete(profile) ;
402 cpl_vector_delete(profile) ;
423 (
const hawki_distortion * distortion,
426 double * x_pos_distcorr,
427 double * y_pos_distcorr)
432 if(distortion == NULL)
434 cpl_error_set(
"hawki_distortion_correct_coords", CPL_ERROR_ILLEGAL_INPUT);
439 (distortion, x_pos, y_pos, &x_dist, &y_dist);
441 *x_pos_distcorr = x_pos - x_dist;
442 *y_pos_distcorr = y_pos - y_dist;
467 (
const hawki_distortion * distortion,
470 double * x_pos_distinvcorr,
471 double * y_pos_distinvcorr)
478 if(distortion == NULL)
480 cpl_error_set(
"hawki_distortion_inverse_correct_coords", CPL_ERROR_ILLEGAL_INPUT);
483 for(i = 0; i < niter; ++i)
486 (distortion, x_pos + x_dist, y_pos + y_dist, &x_dist, &y_dist);
491 *x_pos_distinvcorr = x_pos + x_dist;
492 *y_pos_distinvcorr = y_pos + y_dist;
526 (
const hawki_distortion * distortion,
553 nx = cpl_image_get_size_x(distortion->dist_x);
554 ny = cpl_image_get_size_y(distortion->dist_x);
558 ix1 = (int)floor((x_pos - distortion->x_crval)/distortion->x_cdelt);
559 iy1 = (int)floor((y_pos - distortion->y_crval)/distortion->y_cdelt);
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;
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);
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);
617 (cpl_imagelist * ilist,
621 cpl_image * corr[HAWKI_NB_DETECTORS] ;
625 if (ilist == NULL)
return -1 ;
626 if (dist_x == NULL)
return -1 ;
627 if (dist_y == NULL)
return -1 ;
630 for (i=0 ; i<HAWKI_NB_DETECTORS ; i++)
632 cpl_image * cur_image;
635 cur_image = cpl_imagelist_get(ilist, i);
641 cpl_msg_error(__func__,
"Cannot correct the distortion") ;
642 for (j=0 ; j<i ; j++) cpl_image_delete(corr[j]) ;
648 for (i=0 ; i<HAWKI_NB_DETECTORS ; i++)
649 cpl_imagelist_set(ilist, corr[i], i);
666 int hawki_distortion_create_maps_detector
667 (
const hawki_distortion * distortion,
678 if (distortion == NULL)
return -1 ;
679 if (dist_x == NULL)
return -1 ;
680 if (dist_y == NULL)
return -1 ;
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 ;
689 pdx = cpl_image_get_data_double(dist_x) ;
690 pdy = cpl_image_get_data_double(dist_y) ;
693 for (j=0 ; j<ny ; j++) {
694 for (i=0 ; i<nx ; i++) {
701 (distortion, (
double)i, (
double)j, &x_dist, &y_dist);
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,
744 const hawki_distortion * initial_guess,
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;
758 double tolerance = HAWKI_DISTORTION_TOLERANCE;
759 double tolerance2 = HAWKI_DISTORTION_TOLERANCE2;
760 int minimizer_status;
765 if(initial_guess == NULL)
768 (detector_nx, detector_ny, grid_size);
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;
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)
788 cpl_msg_error(__func__,
"Too few matches to compute distortion (< %d)",
795 args.ref_catalogues = ref_catalogues;
796 args.matching_sets = matching_sets;
797 args.distortion = distortion;
798 args.offsets = offsets;
801 minimize_function.f = hawki_distortion_gsl_obj_function;
802 minimize_function.n = nfitparam;
803 minimize_function.params = &args;
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);
823 minimizer_status = gsl_multimin_fminimizer_iterate (minimizer);
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));
831 while (minimizer_status == GSL_CONTINUE && iter < HAWKI_DISTORTION_MAX_ITER);
834 matching_sets, ncats, distortion));
838 gsl_multimin_fminimizer_set(minimizer, &minimize_function,
839 gsl_multimin_fminimizer_x(minimizer), step_size);
844 minimizer_status = gsl_multimin_fminimizer_iterate (minimizer);
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));
852 while (minimizer_status == GSL_CONTINUE && iter < HAWKI_DISTORTION_MAX_ITER2);
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));
861 matching_sets, ncats, distortion);
864 gsl_multimin_fminimizer_free(minimizer);
865 gsl_vector_free(init_param);
866 gsl_vector_free(step_size);
867 cpl_bivector_delete(offsets);
871 cpl_msg_error(__func__,
"Not compiled with GSL support.");
886 double hawki_distortion_gsl_obj_function
887 (
const gsl_vector * dist_param,
891 const cpl_table ** ref_catalogues;
892 const cpl_table * matching_sets;
893 hawki_distortion * distortion;
894 cpl_bivector * offsets;
897 double objective_function;
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;
908 hawki_distortion_update_solution_from_param(distortion, dist_param);
909 hawki_distortion_update_offsets_from_param(offsets, dist_param);
912 matching_sets, ncats, distortion);
914 objective_function = rms;
917 cpl_msg_debug(__func__,
"Objective function: %g", objective_function);
918 return objective_function;
930 (
const cpl_table ** ref_catalogues,
931 const cpl_bivector * cat_offsets,
932 const cpl_table * matching_sets,
934 hawki_distortion * distortion)
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;
951 nmatch = cpl_table_get_nrow(matching_sets);
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));
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)
962 x_cat_cols[icat] = cpl_table_get_data_double_const(ref_catalogues[icat],
964 y_cat_cols[icat] = cpl_table_get_data_double_const(ref_catalogues[icat],
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)
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);
982 #pragma omp parallel for private(icat,imatch) reduction(+:rms)
984 for(imatch = 0; imatch < nmatch; ++imatch)
990 for(icat = 0; icat < ncats; ++icat)
996 x_cat_offset = x_cat_offsets[icat];
997 y_cat_offset = y_cat_offsets[icat];
999 if((iobj = cpl_array_get(match_arrays[imatch], icat, NULL)) != -1)
1010 x_cat = x_cat_cols[icat][iobj];
1011 y_cat = y_cat_cols[icat][iobj];
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;
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;
1033 pos_flag[imatch][icat] = 0;
1037 y_pos_values[imatch], pos_flag[imatch],
1038 ncats, &nstddev, &var_x, &var_y);
1041 rms += sqrt(var_x + var_y) * nstddev;
1044 cpl_free(x_cat_cols);
1045 cpl_free(y_cat_cols);
1046 for(imatch = 0; imatch < nmatch; ++imatch)
1048 cpl_free(x_pos_values[imatch]);
1049 cpl_free(y_pos_values[imatch]);
1050 cpl_free(pos_flag[imatch]);
1052 cpl_free(x_pos_values);
1053 cpl_free(y_pos_values);
1055 cpl_free(match_arrays);
1068 int hawki_distortion_update_solution_from_param
1069 (hawki_distortion * distortion,
1070 const gsl_vector * dist_param)
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)
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));
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);
1108 int hawki_distortion_update_offsets_from_param
1109 (cpl_bivector * offsets,
1110 const gsl_vector * dist_param)
1116 ncats = cpl_bivector_get_size(offsets);
1117 nparam = dist_param->size;
1118 for(i = 0; i < ncats; ++i)
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));
1136 int hawki_distortion_update_param_from_solution
1137 (gsl_vector * dist_param,
1138 const hawki_distortion * distortion)
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)
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));
1168 int hawki_distortion_update_param_from_offsets
1169 (gsl_vector * dist_param,
1170 const cpl_bivector * offsets)
1176 ncats = cpl_bivector_get_size(offsets);
1177 nparam = dist_param->size;
1178 for(i = 0; i < ncats; ++i)
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));
1198 (
double * x_val,
double * y_val,
int *pos_flag,
1199 int nvals,
int *nvalid,
double *var_x,
double *var_y)
1201 double varsum_x = 0.0;
1202 double varsum_y = 0.0;
1203 double mean_x = 0.0;
1204 double mean_y = 0.0;
1208 for (i=0; i < nvals; i++)
1210 if(pos_flag[i] == 1)
1212 const double delta_x = (double)x_val[i] - mean_x;
1213 const double delta_y = (double)y_val[i] - mean_y;
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);
1226 *var_x = varsum_x / (double) (*nvalid - 1);
1227 *var_y = varsum_y / (double) (*nvalid - 1);