39 #ifndef CPL_SIZE_FORMAT
40 #define CPL_SIZE_FORMAT "d"
44 #if defined CPL_VERSION_CODE && CPL_VERSION_CODE <= CPL_VERSION(5, 92, 0)
45 #define cpl_apertures_get_pos_x cpl_apertures_get_max_x
50 #include "irplib_distortion.h"
52 #include "irplib_flat.h"
53 #include "irplib_utils.h"
54 #include "irplib_polynomial.h"
63 #define IRPLIB_MAX(A,B) ((A) > (B) ? (A) : (B))
64 #define IRPLIB_MIN(A,B) ((A) < (B) ? (A) : (B))
66 #define ARC_MINGOODPIX 100
67 #define ARC_MINARCLENFACT 2.0
68 #define ARC_MINNBARCS 4
69 #define ARC_RANGE_FACT 3.0
70 #define ARC_WINDOWSIZE 32
72 #define TRESH_MEDIAN_MIN 0.0
73 #define TRESH_SIGMA_MAX 200.0
85 static cpl_apertures * irplib_distortion_detect_arcs(cpl_image *,
86 cpl_image **,
int,
int,
double,
int,
int,
int,
int);
87 static cpl_error_code irplib_distortion_fill_border(cpl_image *,
int,
int,
89 static int irplib_distortion_threshold1d(cpl_image *,
double, cpl_image *,
91 static cpl_error_code irplib_distortion_purge_arcs(cpl_apertures **, cpl_image *,
92 const cpl_image *,
int,
int,
94 static cpl_error_code irplib_distortion_fill_arc_positions(cpl_bivector *,
98 const cpl_apertures *);
100 static double irplib_distortion_get_row_centroid(
const cpl_image *,
101 const cpl_image *,
int,
int);
103 static int irplib_distortion_sub_hor_lowpass(cpl_image *,
int);
104 static cpl_image * irplib_distortion_remove_ramp(
const cpl_image *);
106 static cpl_error_code irplib_image_filter_background_line(cpl_image *,
107 const cpl_image *,
int, cpl_boolean) ;
109 static cpl_error_code irplib_polynomial_fit_2d(cpl_polynomial *,
110 const cpl_bivector *,
111 const cpl_vector *,
int,
114 static cpl_matrix * irplib_matrix_product_normal_create(
const cpl_matrix *);
152 cpl_polynomial * irplib_distortion_estimate(
153 const cpl_image * org,
163 cpl_apertures ** arcs)
165 cpl_image * local_im;
166 cpl_image * filtered;
167 cpl_image * label_image;
168 double rightmost, leftmost;
170 cpl_vector * values_to_fit;
172 cpl_polynomial * poly2d;
174 const int nx = cpl_image_get_size_x(org);
175 const int ny = cpl_image_get_size_y(org);
176 const int min_arc_range = (int)(nx / ARC_RANGE_FACT);
180 cpl_ensure(org != NULL, CPL_ERROR_NULL_INPUT, NULL);
181 cpl_ensure(kappa >= 0.0, CPL_ERROR_ILLEGAL_INPUT, NULL);
182 cpl_ensure(max_arc_width > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
187 filtered = cpl_image_new(nx, ny, cpl_image_get_type(org));
189 irplib_image_filter_background_line(filtered, org, max_arc_width, CPL_TRUE);
192 local_im = irplib_distortion_remove_ramp(filtered);
193 cpl_image_delete(filtered);
198 cpl_error_ensure(local_im != NULL, cpl_error_get_code(),
199 return(NULL),
"Cannot clean the image");
202 *arcs = irplib_distortion_detect_arcs(local_im, &label_image, arc_sat,
203 max_arc_width, kappa, xmin, ymin,
206 cpl_image_delete(local_im);
207 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
208 "Cannot detect the arcs");
211 n_arcs = cpl_apertures_get_size(*arcs);
212 cpl_msg_info(cpl_func,
"%d detected arcs", n_arcs);
215 rightmost = leftmost = cpl_apertures_get_pos_x(*arcs, 1);
216 for (i=1; i<n_arcs; i++) {
217 if (cpl_apertures_get_pos_x(*arcs, i+1) < leftmost)
218 leftmost = cpl_apertures_get_pos_x(*arcs, i+1);
219 if (cpl_apertures_get_pos_x(*arcs, i+1) > rightmost)
220 rightmost = cpl_apertures_get_pos_x(*arcs, i+1);
222 if ((
int)(rightmost-leftmost) < min_arc_range) {
223 #if defined CPL_HAVE_VA_ARGS && CPL_HAVE_VA_ARGS != 0
224 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
225 "too narrow range (%g-%g)<%d",
226 rightmost, leftmost, min_arc_range);
228 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
231 cpl_apertures_delete(*arcs);
232 cpl_image_delete(local_im);
233 cpl_image_delete(label_image);
239 cpl_msg_info(cpl_func,
"Create deformation grid");
240 grid = cpl_bivector_new(n_arcs * ny);
241 values_to_fit = cpl_vector_new(n_arcs * ny);
243 if (irplib_distortion_fill_arc_positions(grid, values_to_fit, local_im,
244 label_image, *arcs)){
245 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
246 "cannot get arcs positions");
247 cpl_apertures_delete(*arcs);
248 cpl_image_delete(local_im);
249 cpl_image_delete(label_image);
253 cpl_image_delete(label_image);
254 cpl_image_delete(local_im);
257 poly2d = cpl_polynomial_new(2);
258 if (irplib_polynomial_fit_2d(poly2d, grid, values_to_fit, degree,
260 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
261 "cannot apply the 2d fit");
262 cpl_bivector_delete(grid);
263 cpl_vector_delete(values_to_fit);
264 cpl_apertures_delete(*arcs);
269 cpl_msg_info(cpl_func,
270 "Fitted a %d. degree 2D-polynomial to %"CPL_SIZE_FORMAT
" points "
271 "with mean-square error: %g", degree,
272 cpl_vector_get_size(values_to_fit), mse);
275 cpl_bivector_delete(grid);
276 cpl_vector_delete(values_to_fit);
299 static cpl_apertures * irplib_distortion_detect_arcs(
301 cpl_image ** label_im,
310 const int ny = cpl_image_get_size_y(im);
312 const int min_arclen = (int)(ny / ARC_MINARCLENFACT);
315 cpl_image * collapsed;
317 double threshold, fillval, median_val, sigma;
326 median_val = cpl_image_get_median_dev(im, &sigma);
327 fillval = median_val-sigma/2.0;
328 if (irplib_distortion_fill_border(im, xmin, ymin, xmax, ymax, fillval)) {
329 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
330 "cannot fill bad zones");
335 filt_im = cpl_image_duplicate(im);
336 if (irplib_distortion_sub_hor_lowpass(filt_im, ARC_WINDOWSIZE) == -1) {
337 cpl_image_delete(filt_im);
342 median_val = cpl_image_get_median_dev(filt_im, &sigma);
345 if (median_val < TRESH_MEDIAN_MIN) median_val = TRESH_MEDIAN_MIN;
346 if (sigma > TRESH_SIGMA_MAX) sigma = TRESH_SIGMA_MAX;
349 threshold = median_val + sigma * kappa;
352 collapsed = cpl_image_collapse_median_create(filt_im, 0, 0, 0);
355 if (irplib_distortion_threshold1d(filt_im, median_val, collapsed, 0.0)==-1) {
356 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
357 "cannot threshold the filtered image");
358 cpl_image_delete(filt_im);
359 cpl_image_delete(collapsed);
362 cpl_image_delete(collapsed);
365 bin_im = cpl_mask_threshold_image_create(filt_im, threshold,
367 cpl_image_delete(filt_im);
368 if (bin_im == NULL) {
369 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
370 "cannot binarise the image");
375 ngoodpix = cpl_mask_count(bin_im);
376 if (ngoodpix < ARC_MINGOODPIX) {
377 #if defined CPL_HAVE_VA_ARGS && CPL_HAVE_VA_ARGS != 0
378 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
379 "Too few (%d) white pixels", ngoodpix);
381 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
382 "Too few white pixels");
384 cpl_mask_delete(bin_im);
389 filter = cpl_mask_new(3, 3);
390 cpl_mask_not(filter);
391 cpl_mask_filter(bin_im, bin_im, filter, CPL_FILTER_OPENING,
393 cpl_mask_delete(filter);
396 *label_im = cpl_image_labelise_mask_create(bin_im, &nobj);
397 cpl_mask_delete(bin_im);
400 if ((det = cpl_apertures_new_from_image(im, *label_im)) == NULL) {
401 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
402 "Cannot compute arcs stats");
403 cpl_image_delete(*label_im);
409 if (irplib_distortion_purge_arcs(&det, *label_im, im, min_arclen,
410 max_arc_width, arc_sat)) {
411 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
412 "Cannot purge the arcs");
413 cpl_image_delete(*label_im);
415 cpl_apertures_delete(det);
418 if (cpl_apertures_get_size(det) < ARC_MINNBARCS) {
419 #if defined CPL_HAVE_VA_ARGS && CPL_HAVE_VA_ARGS != 0
420 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
421 "Not enough valid arcs (%"CPL_SIZE_FORMAT
" < %d)",
422 cpl_apertures_get_size(det), ARC_MINNBARCS);
424 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
425 "Not enough valid arcs, min="
426 IRPLIB_STRINGIFY(ARC_MINNBARCS));
428 cpl_image_delete(*label_im);
430 cpl_apertures_delete(det);
449 static cpl_error_code irplib_distortion_fill_border(cpl_image *
self,
456 const int nx = cpl_image_get_size_x(
self);
457 const int ny = cpl_image_get_size_y(
self);
458 float * pfi = cpl_image_get_data_float(
self);
459 const float fvalue = (float)fillval;
463 cpl_ensure_code(pfi != NULL, cpl_error_get_code());
466 xmin = IRPLIB_MIN(xmin, nx+1);
467 ymax = IRPLIB_MIN(ymax, ny);
470 xmax = IRPLIB_MAX(xmax, xmin - 1);
471 ymin = IRPLIB_MIN(ymin, ymax + 1);
475 for (j = 0; j < ymin-1; j++) {
476 for (i = 0; i < nx; i++) {
477 pfi[i+j*nx] = fvalue;
482 for (; j < ymax; j++) {
483 for (i = 0; i < xmin-1; i++) {
484 pfi[i+j*nx] = fvalue;
486 for (i = xmax; i < nx; i++) {
487 pfi[i+j*nx] = fvalue;
492 for (; j < ny; j++) {
493 for (i = 0; i < nx; i++) {
494 pfi[i+j*nx] = fvalue;
498 return CPL_ERROR_NONE;
501 static int irplib_distortion_threshold1d(
513 if (im == NULL)
return -1;
514 if (im1d == NULL)
return -1;
515 if (cpl_image_get_type(im) != CPL_TYPE_FLOAT)
return -1;
516 if (cpl_image_get_type(im1d) != CPL_TYPE_FLOAT)
return -1;
519 pim = cpl_image_get_data_float(im);
520 pim1d = cpl_image_get_data_float(im1d);
521 nx = cpl_image_get_size_x(im);
522 ny = cpl_image_get_size_y(im);
526 if (pim1d[i] < threshold) {
527 for (j=0; j<ny; j++) pim[i+j*nx] = (
float)newval;
534 static int irplib_distortion_sub_hor_lowpass(
540 cpl_vector * avglinehi;
541 cpl_vector * avglinelo;
544 int lopos, hipos, nx, ny;
548 if (im == NULL)
return -1;
549 if (filt_size <= 0)
return -1;
552 nx = cpl_image_get_size_x(im);
553 ny = cpl_image_get_size_y(im);
555 hipos = (int)(3*ny/4);
558 if ((linehi = cpl_vector_new_from_image_row(im, hipos)) == NULL) {
561 if ((linelo = cpl_vector_new_from_image_row(im, lopos)) == NULL) {
562 cpl_vector_delete(linehi);
567 if ((avglinehi = cpl_vector_filter_median_create(linehi,
568 filt_size)) == NULL) {
569 cpl_vector_delete(linehi);
570 cpl_vector_delete(linelo);
573 cpl_vector_delete(linehi);
575 if ((avglinelo = cpl_vector_filter_median_create(linelo,
576 filt_size)) == NULL) {
577 cpl_vector_delete(linelo);
578 cpl_vector_delete(avglinehi);
581 cpl_vector_delete(linelo);
584 cpl_vector_add(avglinehi, avglinelo);
585 cpl_vector_delete(avglinelo);
586 cpl_vector_divide_scalar(avglinehi, 2.0);
589 pavglinehi = cpl_vector_get_data(avglinehi);
590 pim = cpl_image_get_data_float(im);
591 for (i=0; i<nx; i++) {
592 for (j=0; j<ny; j++) {
593 pim[i+j*nx] -= pavglinehi[i];
596 cpl_vector_delete(avglinehi);
614 cpl_error_code irplib_distortion_purge_arcs(cpl_apertures **
self,
616 const cpl_image * arc_im,
621 const double ycenter = 0.5 * (1 + cpl_image_get_size_y(arc_im));
629 cpl_ensure_code(
self != NULL, CPL_ERROR_NULL_INPUT);
632 narcs = cpl_apertures_get_size(*
self);
634 cpl_ensure_code(narcs > 0, CPL_ERROR_DATA_NOT_FOUND);
635 cpl_ensure_code(cpl_image_get_type(lab_im) == CPL_TYPE_INT,
636 CPL_ERROR_ILLEGAL_INPUT);
639 relabel = cpl_calloc(narcs,
sizeof(
int));
642 for (i = 0; i < narcs; i++) {
645 + cpl_apertures_get_top(*
self, i+1)
646 - cpl_apertures_get_bottom(*
self, i+1);
648 if (cpl_apertures_get_top(*
self, i+1) < ycenter)
continue;
649 if (cpl_apertures_get_bottom(*
self, i+1) > ycenter)
continue;
651 if (arclen > min_arclen) {
652 const int arcwidth = 1
653 + cpl_apertures_get_right(*
self, i+1)
654 - cpl_apertures_get_left(*
self, i+1);
655 if (arcwidth < max_arcwidth) {
656 const int edge = cpl_apertures_get_left_y(*
self, i+1);
658 const double mean = cpl_apertures_get_mean(*
self, i+1);
659 if (mean < arc_sat) {
660 relabel[i] = ++nkeep;
662 if (nkeep == i+1) ifirst = nkeep;
671 int * plabim = cpl_image_get_data_int(lab_im);
672 const int npix = cpl_image_get_size_x(lab_im)
673 * cpl_image_get_size_y(lab_im);
677 #if defined CPL_HAVE_VA_ARGS && CPL_HAVE_VA_ARGS != 0
678 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
679 "All %d arc(s) are invalid", narcs);
681 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
682 "All arcs are invalid");
686 for (i = 0; i < npix; i++) {
687 const int label = plabim[i];
689 if (label < 0 || label > narcs)
break;
690 if (label >= ifirst) plabim[i] = relabel[label-1];
696 return cpl_error_set(cpl_func, plabim[i] < 0
697 ? CPL_ERROR_ILLEGAL_INPUT
698 : CPL_ERROR_INCOMPATIBLE_INPUT);
702 cpl_apertures_delete(*
self);
703 *
self = cpl_apertures_new_from_image(arc_im, lab_im);
709 cpl_msg_info(cpl_func,
"Purged %d of %d arcs (1st purged=%d)", narcs - nkeep,
713 cpl_ensure_code(*
self != NULL, cpl_error_get_code());
715 return CPL_ERROR_NONE;
734 static cpl_error_code
735 irplib_distortion_fill_arc_positions(cpl_bivector * grid,
736 cpl_vector * fitvalues,
737 const cpl_image * in,
738 const cpl_image * label_im,
739 const cpl_apertures * det)
741 const int narcs = cpl_apertures_get_size(det);
742 int nfitvals = cpl_vector_get_size(fitvalues);
743 const int nx = cpl_image_get_size_x(label_im);
744 const int ny = cpl_image_get_size_y(label_im);
745 cpl_image * filt_img;
747 cpl_vector * gridx = cpl_bivector_get_x(grid);
748 cpl_vector * gridy = cpl_bivector_get_y(grid);
749 cpl_polynomial* dist1d;
750 cpl_matrix * dist1dx = NULL;
751 cpl_vector * dist1dy = NULL;
758 cpl_ensure_code(nfitvals > 0, CPL_ERROR_DATA_NOT_FOUND);
759 cpl_ensure_code(narcs > 0, CPL_ERROR_DATA_NOT_FOUND);
760 cpl_ensure_code(cpl_image_get_type(label_im) == CPL_TYPE_INT,
761 CPL_ERROR_TYPE_MISMATCH);
764 if (nfitvals < narcs * ny) {
765 nfitvals = narcs * ny;
766 cpl_vector_set_size(fitvalues, nfitvals);
768 if (cpl_vector_get_size(gridx) < nfitvals ||
769 cpl_vector_get_size(gridy) < nfitvals) {
770 cpl_vector_set_size(gridx, nfitvals);
771 cpl_vector_set_size(gridy, nfitvals);
775 dgridx = cpl_vector_get_data(gridx);
776 dgridy = cpl_vector_get_data(gridy);
777 dfitv = cpl_vector_get_data(fitvalues);
780 kernel = cpl_mask_new(3, 3);
781 cpl_mask_not(kernel);
782 filt_img = cpl_image_new(nx, ny, cpl_image_get_type(in));
783 cpl_image_filter_mask(filt_img, in, kernel, CPL_FILTER_MEDIAN,
785 cpl_mask_delete(kernel);
787 dist1d = cpl_polynomial_new(1);
789 for (obj = 0; obj < narcs; obj++) {
791 const int * plabel_im = cpl_image_get_data_int_const(label_im);
792 const int ndist1d = cpl_apertures_get_top(det, obj+1)
793 - cpl_apertures_get_bottom(det, obj+1) + 1;
794 cpl_boolean sampsym = CPL_TRUE;
798 (void)cpl_matrix_unwrap(dist1dx);
799 (void)cpl_vector_unwrap(dist1dy);
800 dist1dx = cpl_matrix_wrap(1, ndist1d, dgridy + ndone);
801 dist1dy = cpl_vector_wrap(ndist1d, dfitv + ndone);
805 for (j = cpl_apertures_get_bottom(det, obj+1)-1;
806 j < cpl_apertures_get_top(det, obj+1); j++) {
808 for (i = 0; i < nx; i++) {
809 if (plabel_im[i + j * nx] == obj + 1)
break;
813 cpl_errorstate prestate = cpl_errorstate_get();
815 const double x_finepos
816 = irplib_distortion_get_row_centroid(filt_img, label_im,
818 if (!cpl_errorstate_is_equal(prestate)) {
819 irplib_error_recover(prestate,
"Could not find X-position "
820 "for line %d at y=%d (x=%d)",
822 }
else if (x_finepos >= 0.0) {
823 cpl_matrix_set(dist1dx, 0, k, 1.0 + j);
824 cpl_vector_set(dist1dy, k, 1.0 + x_finepos);
825 if (k > 0 && j != 1 + prevj) sampsym = CPL_FALSE;
832 double ref_xpos, grad;
833 cpl_error_code error;
834 const cpl_boolean did_drop = k != ndist1d;
835 const cpl_size mindeg = 0;
836 const cpl_size maxdeg = 2;
840 dist1dx = cpl_matrix_wrap(1, k, cpl_matrix_unwrap(dist1dx));
841 dist1dy = cpl_vector_wrap(k, cpl_vector_unwrap(dist1dy));
844 error = cpl_polynomial_fit(dist1d, dist1dx, &sampsym, dist1dy, NULL,
845 CPL_FALSE, &mindeg, &maxdeg);
847 cpl_msg_error(cpl_func,
"1D-fit failed");
851 ref_xpos = cpl_polynomial_eval_1d(dist1d, 0.5 * (ny + 1), &grad);
853 for (j = cpl_apertures_get_bottom(det, obj+1)-1;
854 j < cpl_apertures_get_top(det, obj+1); j++) {
855 const double xpos = cpl_polynomial_eval_1d(dist1d, j+1.0, NULL);
857 dfitv [ndone] = xpos;
858 dgridx[ndone] = ref_xpos;
862 dgridy[ndone] = 1.0 + j;
865 cpl_msg_info(cpl_func,
"Line %d has center gradient %g", obj+1,
870 cpl_image_delete(filt_img);
871 cpl_polynomial_delete(dist1d);
872 (void)cpl_matrix_unwrap(dist1dx);
873 (void)cpl_vector_unwrap(dist1dy);
875 cpl_msg_info(cpl_func,
"Found %d fitting points ("
876 "expected up to %d points)", ndone, nfitvals);
878 cpl_ensure_code(obj == narcs, cpl_error_get_code());
880 cpl_ensure_code(ndone > 0, CPL_ERROR_DATA_NOT_FOUND);
882 cpl_vector_set_size(fitvalues, ndone);
883 cpl_vector_set_size(gridx, ndone);
884 cpl_vector_set_size(gridy, ndone);
886 return CPL_ERROR_NONE;
900 static double irplib_distortion_get_row_centroid(
const cpl_image * im,
901 const cpl_image * label_im,
905 const int nx = cpl_image_get_size_x(im);
906 const int ny = cpl_image_get_size_y(im);
907 const int ynx = y * nx;
908 const float * pim = cpl_image_get_data_float_const(im);
909 const int * plabel_im = cpl_image_get_data_int_const(label_im);
918 cpl_ensure(pim != NULL, cpl_error_get_code(), -1.0);
919 cpl_ensure(plabel_im != NULL, cpl_error_get_code(), -2.0);
920 cpl_ensure(x >= 0, CPL_ERROR_ILLEGAL_INPUT, -3.0);
921 cpl_ensure(y >= 0, CPL_ERROR_ILLEGAL_INPUT, -4.0);
922 cpl_ensure(x < nx, CPL_ERROR_ILLEGAL_INPUT, -5.0);
923 cpl_ensure(y < ny, CPL_ERROR_ILLEGAL_INPUT, -6.0);
925 max = (double)pim[x + ynx];
926 objnum = plabel_im[x + ynx];
930 const double val = (double)pim[x + ynx];
936 if (firstpos < 0) firstpos = x;
949 }
while (x < nx && objnum == plabel_im[x + ynx]);
951 cpl_ensure(sum > 0.0, CPL_ERROR_DATA_NOT_FOUND, -7.0);
959 return (wsum < sum * firstpos || wsum > sum * lastpos)
960 ? maxpos : wsum / sum;
970 #define IS_NB_TESTPOINTS 8
971 #define IS_MIN_SLOPE 0.01
972 #define IS_MAX_SLOPE_DIF 0.075
973 #define IS_MAX_FIT_EDGE_DIF 0.05
974 #define IS_MIN_RAMP 10.0
975 #define IS_MAX_MNERR 13.0
976 #define IS_MAX_MNERR_DIF 8.0
977 #define IS_MAX_INTER_DIF 20.0
978 #define IS_SKIPZONE 2.5
979 #define SQR(x) ((x)*(x))
980 static cpl_image * irplib_distortion_remove_ramp(
const cpl_image * in)
983 const int nx = cpl_image_get_size_x(in);
984 const int ny = cpl_image_get_size_y(in);
985 const int yhi = (int)(ny/2);
986 const int ylo = yhi - 1;
988 cpl_vector * tmp_vector;
989 cpl_bivector * testpointlo;
990 double * testpointlo_x;
991 double * testpointlo_y;
992 cpl_bivector * testpointhi;
993 double * testpointhi_x;
994 double * testpointhi_y;
995 const int spacing = ny / (IS_SKIPZONE*IS_NB_TESTPOINTS);
996 double rampdif, fitslope;
1000 double * median_data;
1001 double medianerrlo, medianerrhi;
1008 cpl_ensure(cpl_image_get_type(in) == CPL_TYPE_FLOAT,
1009 CPL_ERROR_UNSUPPORTED_MODE, NULL);
1011 if (ny < IS_SKIPZONE * IS_NB_TESTPOINTS){
1012 #if defined CPL_HAVE_VA_ARGS && CPL_HAVE_VA_ARGS != 0
1013 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
1014 "image has %d lines, min="
1015 IRPLIB_STRINGIFY(IS_SKIPZONE)
"*"
1016 IRPLIB_STRINGIFY(IS_NB_TESTPOINTS), ny);
1018 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
1019 "image has too few lines, min="
1020 IRPLIB_STRINGIFY(IS_SKIPZONE)
"*"
1021 IRPLIB_STRINGIFY(IS_NB_TESTPOINTS));
1027 testpointhi = cpl_bivector_new(IS_NB_TESTPOINTS);
1028 testpointhi_x = cpl_bivector_get_x_data(testpointhi);
1029 testpointhi_y = cpl_bivector_get_y_data(testpointhi);
1030 testpointlo = cpl_bivector_new(IS_NB_TESTPOINTS);
1031 testpointlo_x = cpl_bivector_get_x_data(testpointlo);
1032 testpointlo_y = cpl_bivector_get_y_data(testpointlo);
1033 for (i=0; i<IS_NB_TESTPOINTS; i++) {
1034 y = yhi + i * spacing;
1035 tmp_vector = cpl_vector_new_from_image_row(in, y+1);
1036 testpointhi_x[i] = y - ny / 2;
1037 testpointhi_y[i] = cpl_vector_get_median_const(tmp_vector);
1038 cpl_vector_delete(tmp_vector);
1039 y = ylo - i * spacing;
1040 tmp_vector = cpl_vector_new_from_image_row(in, y+1);
1041 testpointlo_x[IS_NB_TESTPOINTS-i-1] = y;
1042 testpointlo_y[IS_NB_TESTPOINTS-i-1]=cpl_vector_get_median_const(tmp_vector);
1043 cpl_vector_delete(tmp_vector);
1048 testpointhi_y, IS_NB_TESTPOINTS);
1050 testpointlo_y, IS_NB_TESTPOINTS);
1053 median = cpl_vector_new(IS_NB_TESTPOINTS);
1054 median_data = cpl_vector_get_data(median);
1055 for (i=0; i<IS_NB_TESTPOINTS; i++) {
1056 median_data[i]=SQR(testpointhi_y[i]
1057 - pol_coefhi[0] - pol_coefhi[1] * testpointhi_x[i]);
1059 medianerrhi = cpl_vector_get_median(median);
1060 for (i=0; i<IS_NB_TESTPOINTS; i++) {
1061 median_data[i]=SQR(testpointlo_y[i]
1062 - pol_coeflo[0] - pol_coeflo[1] * testpointlo_x[i]);
1064 medianerrlo = cpl_vector_get_median(median);
1065 cpl_vector_delete(median);
1066 rampdif = testpointlo_y[IS_NB_TESTPOINTS-1] - testpointhi_y[0];
1067 slope = rampdif / (ny/2.0);
1068 fitslope = (pol_coefhi[1] + pol_coeflo[1]) / 2.0;
1070 cpl_bivector_delete(testpointlo);
1071 cpl_bivector_delete(testpointhi);
1074 if (fabs(rampdif)<IS_MIN_RAMP ||
1075 fabs(pol_coefhi[1]) < IS_MIN_SLOPE ||
1076 fabs(pol_coeflo[1]) < IS_MIN_SLOPE ||
1077 pol_coefhi[1]/pol_coeflo[1]<0.5 ||
1078 pol_coefhi[1]/pol_coeflo[1]>2.0 ||
1079 fabs(pol_coefhi[1]-pol_coeflo[1])>IS_MAX_SLOPE_DIF ||
1080 fabs(pol_coefhi[0]-pol_coeflo[0]) > IS_MAX_INTER_DIF ||
1081 medianerrlo> IS_MAX_MNERR ||
1082 medianerrhi> IS_MAX_MNERR ||
1083 fabs(medianerrlo-medianerrhi) >IS_MAX_MNERR_DIF ||
1084 fabs(slope-fitslope) > IS_MAX_FIT_EDGE_DIF ||
1085 slope/fitslope<0.5 ||
1086 slope/fitslope>2.0) ramp_present = 0;
1087 else ramp_present = 1;
1089 cpl_free(pol_coeflo);
1090 cpl_free(pol_coefhi);
1093 out = cpl_image_duplicate(in);
1094 pout = cpl_image_get_data_float(out);
1095 if (ramp_present == 1) {
1096 for (j=0; j<ny/2; j++) {
1097 val = slope * (j-ny/2);
1098 for (i=0; i<nx; i++)
1099 pout[i+j*nx] -= val;
1101 for (j=ny/2; j<ny; j++) {
1102 val = slope * (j-ny);
1103 for (i=0; i<nx; i++)
1104 pout[i+j*nx] -= val;
1127 static cpl_error_code irplib_image_filter_background_line(cpl_image *
self,
1128 const cpl_image * other,
1130 cpl_boolean vertical)
1132 const int nx = cpl_image_get_size_x(
self);
1133 const int ny = cpl_image_get_size_y(
self);
1134 const int msize = 1 + 2 * hsize;
1136 cpl_image * background;
1137 cpl_error_code error = CPL_ERROR_NONE;
1139 cpl_ensure_code(
self != NULL, CPL_ERROR_NULL_INPUT);
1140 cpl_ensure_code(hsize >= 0, CPL_ERROR_ILLEGAL_INPUT);
1142 if (other == NULL) other =
self;
1144 mask = vertical ? cpl_mask_new(msize, 1) : cpl_mask_new(1, msize);
1146 error |= cpl_mask_not(mask);
1148 background = cpl_image_new(nx, ny, cpl_image_get_type(other));
1150 error |= cpl_image_filter_mask(background, other, mask, CPL_FILTER_MEDIAN,
1152 cpl_mask_delete(mask);
1154 if (
self != other) {
1155 error |= cpl_image_copy(
self, other, 1, 1);
1158 error |= cpl_image_subtract(
self, background);
1159 cpl_image_delete(background);
1161 return error ? cpl_error_set_where(cpl_func) : CPL_ERROR_NONE;
1191 static cpl_matrix * irplib_matrix_product_normal_create(
const cpl_matrix *
self)
1195 cpl_matrix * product;
1196 const double * ai = cpl_matrix_get_data_const(
self);
1199 const int m = cpl_matrix_get_nrow(
self);
1200 const int n = cpl_matrix_get_ncol(
self);
1204 cpl_ensure(
self != NULL, CPL_ERROR_NULL_INPUT, NULL);
1210 product = cpl_matrix_new(m, m);
1211 bwrite = cpl_matrix_get_data(product);
1213 bwrite = (
double *) cpl_malloc(m * m *
sizeof(
double));
1214 product = cpl_matrix_wrap(m, m, bwrite);
1218 for (i = 0; i < m; i++, bwrite += m, ai += n) {
1220 for (j = i; j < m; j++, aj += n) {
1222 for (k = 0; k < n; k++) {
1223 sum += ai[k] * aj[k];
1248 static cpl_error_code irplib_polynomial_fit_2d(cpl_polynomial *
self,
1249 const cpl_bivector * xy_pos,
1250 const cpl_vector * values,
1251 int degree,
double fixy,
1255 const int np = cpl_bivector_get_size(xy_pos);
1257 const int nc1 = 1+degree;
1260 const int nc = nc1 * (1 + nc1) / 2 - nc1;
1265 #ifdef IRPLIB_DISTORTION_ASSERT
1266 const double * coeffs1d;
1275 cpl_error_code error;
1278 cpl_ensure_code(
self != NULL, CPL_ERROR_NULL_INPUT);
1279 cpl_ensure_code(cpl_polynomial_get_dimension(
self) == 2,
1280 CPL_ERROR_INVALID_TYPE);
1281 cpl_ensure_code(np > 0, cpl_error_get_code());
1282 cpl_ensure_code(values != NULL, CPL_ERROR_NULL_INPUT);
1284 cpl_ensure_code(cpl_vector_get_size(values) == np,
1285 CPL_ERROR_INCOMPATIBLE_INPUT);
1287 cpl_ensure_code(degree > 0, CPL_ERROR_ILLEGAL_INPUT);
1288 cpl_ensure_code(np >= nc, CPL_ERROR_DATA_NOT_FOUND);
1291 yhat = cpl_vector_duplicate(cpl_bivector_get_y_const(xy_pos));
1292 cpl_vector_subtract_scalar(yhat, fixy);
1295 xhat = cpl_vector_duplicate(cpl_bivector_get_x_const(xy_pos));
1296 zhat = cpl_vector_duplicate(values);
1297 cpl_vector_subtract(zhat, xhat);
1302 dmv = (
double*)cpl_malloc(nc*np*
sizeof(
double));
1303 mv = cpl_matrix_wrap(nc, np, dmv);
1306 for (i=0; i < np; i++) {
1307 const double x = cpl_vector_get(xhat, i);
1308 const double y = cpl_vector_get(yhat, i);
1312 for (degy = 1; degy <= degree; degy++) {
1314 for (degx = 0; degx <= degree-degy; degx++, j++) {
1315 dmv[np * j + i] = xvalue * yvalue;
1322 cpl_vector_delete(xhat);
1323 cpl_vector_delete(yhat);
1326 mb = cpl_matrix_wrap(np, 1, cpl_vector_get_data(zhat));
1329 mx = cpl_matrix_product_create(mv, mb);
1331 cpl_matrix_unwrap(mb);
1332 cpl_vector_delete(zhat);
1335 mh = irplib_matrix_product_normal_create(mv);
1336 cpl_matrix_delete(mv);
1339 error = cpl_matrix_decomp_chol(mh) || cpl_matrix_solve_chol(mh, mx);
1341 cpl_matrix_delete(mh);
1344 cpl_matrix_delete(mx);
1345 cpl_ensure_code(0, error);
1350 #ifdef IRPLIB_DISTORTION_ASSERT
1351 coeffs1d = cpl_matrix_get_data(mx);
1355 for (degy = 1; degy <= degree; degy++) {
1357 for (degx = 0; degx <= degree-degy; degx++, j++) {
1360 cpl_polynomial_set_coeff(
self, powers, cpl_matrix_get(mx, j, 0));
1365 cpl_matrix_delete(mx);
1370 cpl_polynomial_set_coeff(
self, powers, 1.0);
1373 cpl_polynomial_shift_1d(
self, 1, -fixy);
1377 const cpl_vector * x_pos = cpl_bivector_get_x_const(xy_pos);
1378 const cpl_vector * y_pos = cpl_bivector_get_y_const(xy_pos);
1379 cpl_vector * x_val = cpl_vector_new(2);
1383 for (i=0; i<np; i++) {
1384 cpl_vector_set(x_val, 0, cpl_vector_get(x_pos, i));
1385 cpl_vector_set(x_val, 1, cpl_vector_get(y_pos, i));
1387 residue = cpl_vector_get(values, i)
1388 - cpl_polynomial_eval(
self, x_val);
1389 *mse += residue * residue;
1391 cpl_vector_delete(x_val);
1396 return CPL_ERROR_NONE;