29#include <cxmessages.h>
30#include <cxstrutils.h>
34#include <cpl_parameter.h>
40#include "gichebyshev.h"
67typedef struct GiBiasResults GiBiasResults;
75_giraffe_biasresults_clear(GiBiasResults *self)
85 cx_string_delete(self->limits);
90 cpl_matrix_delete(self->coeffs);
95 cpl_image_delete(self->model);
111_giraffe_method_string(cx_string *
string, GiBiasMethod method,
116 case GIBIAS_METHOD_UNIFORM:
117 cx_string_set(
string,
"UNIFORM");
120 case GIBIAS_METHOD_PLANE:
121 cx_string_set(
string,
"PLANE");
124 case GIBIAS_METHOD_CURVE:
125 cx_string_set(
string,
"CURVE");
128 case GIBIAS_METHOD_PROFILE:
129 cx_string_set(
string,
"PROFILE");
132 case GIBIAS_METHOD_MASTER:
133 cx_string_set(
string,
"MASTER");
136 case GIBIAS_METHOD_ZMASTER:
137 cx_string_set(
string,
"ZMASTER");
144 if (option != GIBIAS_OPTION_UNDEFINED) {
146 case GIBIAS_OPTION_PLANE:
147 cx_string_append(
string,
"+PLANE");
150 case GIBIAS_OPTION_CURVE:
151 cx_string_append(
string,
"+CURVE");
165_giraffe_stringify_coefficients(cx_string *
string, cpl_matrix *matrix)
170 cxchar *tmp = cx_line_alloc();
172 cxint nr = cpl_matrix_get_nrow(matrix);
173 cxint nc = cpl_matrix_get_ncol(matrix);
175 cxsize sz = cx_line_max();
177 cxdouble *data = cpl_matrix_get_data(matrix);
180 for (i = 0; i < nr; i++) {
181 for (j = 0; j < nc; j++) {
182 snprintf(tmp, sz,
"%g", data[i * nc + j]);
183 cx_string_append(
string, tmp);
185 if (i != nr - 1 || j < nc - 1) {
186 cx_string_append(
string,
":");
214_giraffe_compare_overscans(
const GiImage* image1,
const GiImage* image2)
217 cxint32 l1ovscx = -1;
218 cxint32 l1ovscy = -1;
219 cxint32 l1prscx = -1;
220 cxint32 l1prscy = -1;
221 cxint32 l2ovscx = -1;
222 cxint32 l2ovscy = -1;
223 cxint32 l2prscx = -1;
224 cxint32 l2prscy = -1;
226 cpl_propertylist *l1, *l2;
229 cx_assert(image1 != NULL && image2 != NULL);
234 if (cpl_propertylist_has(l1, GIALIAS_OVSCX)) {
235 l1ovscx = cpl_propertylist_get_int(l1, GIALIAS_OVSCX);
237 if (cpl_propertylist_has(l1, GIALIAS_OVSCY)) {
238 l1ovscy = cpl_propertylist_get_int(l1, GIALIAS_OVSCY);
240 if (cpl_propertylist_has(l1, GIALIAS_PRSCX)) {
241 l1prscx = cpl_propertylist_get_int(l1, GIALIAS_PRSCX);
243 if (cpl_propertylist_has(l1, GIALIAS_PRSCY)) {
244 l1prscy = cpl_propertylist_get_int(l1, GIALIAS_PRSCY);
247 if (cpl_propertylist_has(l2, GIALIAS_OVSCX)) {
248 l2ovscx = cpl_propertylist_get_int(l2, GIALIAS_OVSCX);
250 if (cpl_propertylist_has(l2, GIALIAS_OVSCY)) {
251 l2ovscy = cpl_propertylist_get_int(l2, GIALIAS_OVSCY);
253 if (cpl_propertylist_has(l2, GIALIAS_PRSCX)) {
254 l2prscx = cpl_propertylist_get_int(l2, GIALIAS_PRSCX);
256 if (cpl_propertylist_has(l2, GIALIAS_PRSCY)) {
257 l2prscy = cpl_propertylist_get_int(l2, GIALIAS_PRSCY);
260 if (l1ovscx != l2ovscx || l1ovscy != l2ovscy) {
264 if (l1prscx != l2prscx || l1prscy != l2prscy) {
286inline static cpl_matrix*
287_giraffe_bias_get_areas(
const cxchar*
string)
290 cxchar** regions = NULL;
292 cpl_matrix* areas = NULL;
295 cx_assert(
string != NULL);
297 regions = cx_strsplit(
string,
",", -1);
299 if (regions == NULL) {
306 const cxsize nvalues = 4;
311 while (regions[i] != NULL) {
316 areas = cpl_matrix_new(nregions, nvalues);
319 while (regions[i] != NULL) {
321 register cxsize j = 0;
323 cxchar** limits = cx_strsplit(regions[i],
":", nvalues);
326 if (limits == NULL) {
328 cpl_matrix_delete(areas);
335 for (j = 0; j < nvalues; ++j) {
343 if (limits[j] == NULL || *limits[j] ==
'\0') {
350 value = strtol(limits[j], &last, 10);
356 if ((errno == ERANGE &&
357 (value == LONG_MAX || value == LONG_MIN)) ||
358 (errno != 0 && value == 0)) {
375 cpl_matrix_set(areas, i, j, value);
384 cpl_matrix_delete(areas);
387 cx_strfreev(regions);
398 cx_strfreev(regions);
455_giraffe_bias_compute_mean(GiBiasResults* results,
const cpl_image* image,
456 const cpl_matrix* areas, cxdouble kappa, cxint numiter,
457 cxdouble maxfraction)
460 const cxchar*
const fctid =
"giraffe_bias_compute_mean";
463 const cxdouble* pdimg = NULL;
466 cxint img_dimx, img_dimy;
470 cxint x0, x1, x2, x3;
475 cxlong naccepted = 0L;
477 cxlong pixcount = 0L;
479 cxdouble currfraction = 2.;
481 cx_string* tmp = NULL;
483 cpl_matrix* matrix_zz1;
484 cpl_matrix* matrix_zz1diff;
491 if (results->limits == NULL) {
492 cpl_msg_info(fctid,
"Unable to store biaslimits return parameter, "
497 if (cpl_image_get_type(image) != CPL_TYPE_DOUBLE) {
498 cpl_msg_info(fctid,
"Only images allowed of type double, "
513 cpl_msg_info(fctid,
"Bias Areas: Missing bias areas, "
519 ba_num_rows = cpl_matrix_get_nrow(areas);
521 for (j = 0; j < ba_num_rows; j++) {
522 x3 = (cxint) cpl_matrix_get(areas, j, 3);
523 x2 = (cxint) cpl_matrix_get(areas, j, 2);
524 x1 = (cxint) cpl_matrix_get(areas, j, 1);
525 x0 = (cxint) cpl_matrix_get(areas, j, 0);
527 ntotal += (cxulong) ((x3 - x2 + 1) * (x1 - x0 + 1));
531 cpl_msg_info(fctid,
"Bias Areas: Inconsistent specification, "
536 matrix_zz1 = cpl_matrix_new(ntotal, 1);
537 matrix_zz1diff = cpl_matrix_new(ntotal, 1);
544 img_dimx = cpl_image_get_size_x(image);
545 img_dimy = cpl_image_get_size_y(image);
547 cx_string_set(results->limits,
"");
549 for (j = 0; j < ba_num_rows; j++) {
551 x3 = (cxint) cpl_matrix_get(areas, j, 3);
552 x2 = (cxint) cpl_matrix_get(areas, j, 2);
553 x1 = (cxint) cpl_matrix_get(areas, j, 1);
554 x0 = (cxint) cpl_matrix_get(areas, j, 0);
556 if ((x0 > img_dimx) || (x1 > img_dimx) ||
557 (x2 > img_dimy) || (x3 > img_dimy)) {
561 tmp = cx_string_new();
563 cx_string_sprintf(tmp,
"%d:%d:%d:%d;", x0, x1, x2, x3);
564 cx_string_append(results->limits, cx_string_get(tmp));
566 cx_string_delete(tmp);
569 pdimg = cpl_image_get_data_double_const(image);
571 for (l = x2; l < x3 + 1; l++) {
572 for (k = x0; k < x1 + 1; k++) {
573 cpl_matrix_set(matrix_zz1, n, 1, pdimg[k + l * img_dimx]);
580 cpl_msg_info(fctid,
"Bias Areas: Validation failed, aborting ...");
582 cpl_matrix_delete(matrix_zz1);
583 cpl_matrix_delete(matrix_zz1diff);
588 cpl_msg_info(fctid,
"Bias Areas: Using %s",
589 cx_string_get(results->limits));
600 cpl_msg_info(fctid,
"Sigma Clipping : Start");
604 while ((naccepted > 0) && (curriter < numiter) &&
605 (currfraction > maxfraction)) {
607 cxdouble ksigma = 0.;
609 results->mean = cpl_matrix_get_mean(matrix_zz1);
612 for (k = 0; k < cpl_matrix_get_nrow(matrix_zz1); k++) {
613 cpl_matrix_set(matrix_zz1diff, k, 0,
614 cpl_matrix_get(matrix_zz1, k, 0) - results->mean);
617 cpl_msg_info(fctid,
"bias[%d]: mean = %5g, sigma = %6.3g, "
618 "ratio = %6.3g, accepted = %ld\n", curriter,
619 results->mean, sigma, currfraction, naccepted);
621 ksigma = sigma * kappa;
624 for (l = 0; l < cpl_matrix_get_nrow(matrix_zz1); l++) {
625 if (fabs(cpl_matrix_get(matrix_zz1diff, l, 0)) > ksigma) {
629 cpl_matrix_set(matrix_zz1, pixcount, 0,
630 cpl_matrix_get(matrix_zz1, l, 0));
634 cpl_matrix_resize(matrix_zz1, 0, 0, 0, pixcount -
635 cpl_matrix_get_nrow(matrix_zz1));
637 cpl_matrix_resize(matrix_zz1diff, 0, 0, 0, pixcount -
638 cpl_matrix_get_nrow(matrix_zz1diff));
640 if (pixcount == naccepted) {
644 naccepted = pixcount;
646 currfraction = (cxdouble) naccepted / (cxdouble) ntotal;
650 cpl_msg_info(fctid,
"Sigma Clipping : End");
657 results->mean = cpl_matrix_get_mean(matrix_zz1);
660 results->sigma = results->rms / sqrt(cpl_matrix_get_nrow(matrix_zz1));
663 cpl_msg_info(fctid,
"Sigma Clipping Results : bias[%d]: mean = %5g, "
664 "sigma = %6.3g, ratio = %6.3g, accepted=%ld\n", curriter,
665 results->mean, results->rms, currfraction, naccepted);
672 if (matrix_zz1 != NULL) {
673 cpl_matrix_delete(matrix_zz1);
676 if (matrix_zz1diff != NULL) {
677 cpl_matrix_delete(matrix_zz1diff);
738_giraffe_bias_compute_plane(GiBiasResults* results,
const cpl_image* image,
739 const cpl_matrix* areas, cxdouble kappa,
740 cxint numiter, cxdouble maxfraction)
743 const cxchar*
const fctid =
"giraffe_bias_compute_plane";
754 cxsize naccepted = 0;
756 cxdouble fraction = 1.;
759 cpl_matrix* xbs = NULL;
760 cpl_matrix* ybs = NULL;
761 cpl_matrix* zbs = NULL;
762 cpl_matrix* coeffs = NULL;
765 cx_assert(results->limits != NULL);
766 cx_assert(results->coeffs == NULL);
768 cx_assert(areas != NULL);
770 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
777 nareas = cpl_matrix_get_nrow(areas);
779 for (j = 0; j < nareas; j++) {
781 cxint x3 = (cxint) cpl_matrix_get(areas, j, 3);
782 cxint x2 = (cxint) cpl_matrix_get(areas, j, 2);
783 cxint x1 = (cxint) cpl_matrix_get(areas, j, 1);
784 cxint x0 = (cxint) cpl_matrix_get(areas, j, 0);
786 ntotal += (cxsize) ((x3 - x2 + 1) * (x1 - x0 + 1));
792 cpl_msg_info(fctid,
"Bias Areas: Inconsistent specification, "
798 nx = cpl_image_get_size_x(image);
799 ny = cpl_image_get_size_y(image);
801 cpl_msg_info(fctid,
"Bias Areas: specified are %zu points in %dx%d "
802 "image", ntotal, nx, ny);
813 cx_string_set(results->limits,
"");
815 xbs = cpl_matrix_new(ntotal, 1);
816 ybs = cpl_matrix_new(ntotal, 1);
817 zbs = cpl_matrix_new(1, ntotal);
819 for (j = 0; j < nareas; ++j) {
821 const cxdouble* _img = cpl_image_get_data_double_const(image);
825 cx_string* tmp = NULL;
828 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
829 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
830 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
831 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
833 if ((x0 > nx) || (x1 > nx) || (x2 > ny) || (x3 > ny)) {
837 tmp = cx_string_new();
839 cx_string_sprintf(tmp,
"%d:%d:%d:%d;", x0, x1, x2, x3);
840 cx_string_append(results->limits, cx_string_get(tmp));
842 cx_string_delete(tmp);
845 for (k = x2; k < x3 + 1; ++k) {
847 register cxint l = 0;
850 for (l = x0; l < x1 + 1; ++l) {
852 cpl_matrix_set(xbs, n, 0, l);
853 cpl_matrix_set(ybs, n, 0, k);
854 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
863 cpl_matrix_set_size(xbs, n, 1);
864 cpl_matrix_set_size(ybs, n, 1);
865 cpl_matrix_set_size(zbs, 1, n);
869 cpl_msg_info(fctid,
"Bias Areas: Validation failed, aborting...");
871 cpl_matrix_delete(xbs);
872 cpl_matrix_delete(ybs);
873 cpl_matrix_delete(zbs);
881 cpl_msg_info(fctid,
"Bias Areas: Using %s [%zu pixels]",
882 cx_string_get(results->limits), ntotal);
889 cpl_msg_info(fctid,
"Sigma Clipping : Start");
893 while ((naccepted > 0) && (iteration < numiter) &&
894 (fraction > maxfraction)) {
898 cxdouble ksigma = 0.;
900 cpl_matrix* base = NULL;
901 cpl_matrix* fit = NULL;
904 base = cpl_matrix_new(3, naccepted);
908 cpl_msg_info(fctid,
"Sigma Clipping: Error creating design "
911 cpl_matrix_delete(zbs);
912 cpl_matrix_delete(ybs);
913 cpl_matrix_delete(xbs);
918 for (k = 0; k < naccepted; ++k) {
920 cpl_matrix_set(base, 0, k, 1.);
921 cpl_matrix_set(base, 1, k, cpl_matrix_get(xbs, k, 0));
922 cpl_matrix_set(base, 2, k, cpl_matrix_get(ybs, k, 0));
926 cpl_matrix_delete(coeffs);
931 if (coeffs == NULL) {
933 cpl_msg_info(fctid,
"Sigma Clipping : Error in least square "
934 "solution, aborting...");
936 cpl_matrix_delete(base);
939 cpl_matrix_delete(xbs);
940 cpl_matrix_delete(ybs);
941 cpl_matrix_delete(zbs);
952 fit = cpl_matrix_product_create(coeffs, base);
954 cpl_matrix_delete(base);
957 results->mean = cpl_matrix_get_mean(fit);
961 cpl_msg_info(fctid,
"Sigma Clipping : bias plane[%d]: %g + "
962 "%g * x + %g * y, sigma = %.5g, ratio = %.4g, "
963 "accepted = %zu\n", iteration,
964 cpl_matrix_get(coeffs, 0, 0),
965 cpl_matrix_get(coeffs, 0, 1),
966 cpl_matrix_get(coeffs, 0, 2),
967 sigma, fraction, naccepted);
974 ksigma = sigma * kappa;
978 for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
980 register cxdouble z = cpl_matrix_get(zbs, 0, j);
982 if (fabs(cpl_matrix_get(fit, 0, j) - z) > ksigma) {
986 cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
988 cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
990 cpl_matrix_set(zbs, 0, n, z);
995 cpl_matrix_set_size(xbs, n, 1);
996 cpl_matrix_set_size(ybs, n, 1);
997 cpl_matrix_set_size(zbs, 1, n);
999 cpl_matrix_delete(fit);
1002 if (n == naccepted) {
1008 fraction = (cxdouble)naccepted / (cxdouble)ntotal;
1013 cpl_msg_info(fctid,
"Sigma Clipping : End");
1020 results->coeffs = coeffs;
1021 results->rms = sigma;
1029 results->sigma = sigma / sqrt(cpl_matrix_get_ncol(zbs));
1031 cpl_msg_info(fctid,
"Sigma Clipping Results (%d/%zu, sigma = %g)",
1032 iteration, naccepted, results->rms);
1039 cpl_matrix_delete(xbs);
1042 cpl_matrix_delete(ybs);
1045 cpl_matrix_delete(zbs);
1048 return EXIT_SUCCESS;
1095_giraffe_bias_compute_curve(GiBiasResults* results,
const cpl_image* image,
1096 const cpl_matrix *areas, cxdouble kappa,
1097 cxint numiter, cxdouble maxfraction,
1098 cxdouble xdeg, cxdouble ydeg,
1099 cxdouble xstep, cxdouble ystep)
1102 const cxchar*
const fctid =
"giraffe_bias_compute_curve";
1108 cxint iteration = 0;
1112 cxsize naccepted = 0;
1114 cxdouble fraction = 1.;
1115 cxdouble sigma = 0.;
1117 cpl_matrix* xbs = NULL;
1118 cpl_matrix* ybs = NULL;
1119 cpl_matrix* zbs = NULL;
1121 GiChebyshev2D* fit = NULL;
1124 cx_assert(results != NULL);
1125 cx_assert(results->limits != NULL);
1126 cx_assert(results->coeffs == NULL);
1128 cx_assert(areas != NULL);
1130 cx_assert(image != NULL);
1131 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
1138 nareas = cpl_matrix_get_nrow(areas);
1140 for (j = 0; j < nareas; ++j) {
1142 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
1143 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
1144 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
1145 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
1147 ntotal += (cxulong) (ceil((1. + x3 - x2) / ystep) *
1148 ceil((1. + x1 - x0) / xstep));
1151 nx = cpl_image_get_size_x(image);
1152 ny = cpl_image_get_size_y(image);
1154 cpl_msg_info(fctid,
"Bias Areas: Found %zu points in %dx%d image",
1163 results->sigma = 0.;
1166 cx_string_set(results->limits,
"");
1168 xbs = cpl_matrix_new(ntotal, 1);
1169 ybs = cpl_matrix_new(ntotal, 1);
1170 zbs = cpl_matrix_new(1, ntotal);
1172 for (j = 0; j < nareas; ++j) {
1174 const cxdouble* _img = cpl_image_get_data_double_const(image);
1178 cx_string* tmp = NULL;
1181 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
1182 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
1183 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
1184 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
1186 if ((x0 > nx) || (x1 > nx) ||
1187 (x2 > ny) || (x3 > ny)) {
1191 tmp = cx_string_new();
1193 cx_string_sprintf(tmp,
"%d:%d:%d:%d;", x0, x1, x2, x3);
1194 cx_string_append(results->limits, cx_string_get(tmp));
1196 cx_string_delete(tmp);
1199 for (k = x2; k < x3 + 1; k += ystep) {
1201 register cxint l = 0;
1204 for (l = x0; l < x1 + 1; l += xstep) {
1206 cpl_matrix_set(xbs, n, 0, l);
1207 cpl_matrix_set(ybs, n, 0, k);
1208 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
1217 cpl_matrix_set_size(xbs, n, 1);
1218 cpl_matrix_set_size(ybs, n, 1);
1219 cpl_matrix_set_size(zbs, 1, n);
1223 cpl_msg_info(fctid,
"Bias Areas: Validation failed, aborting...");
1225 cpl_matrix_delete(xbs);
1226 cpl_matrix_delete(ybs);
1227 cpl_matrix_delete(zbs);
1235 cpl_msg_info(fctid,
"Bias Areas: Using %s [%zu pixels]",
1236 cx_string_get(results->limits), ntotal);
1243 cpl_msg_info(fctid,
"Sigma Clipping : Start");
1247 while ((naccepted > 0) && (iteration < numiter) &&
1248 (fraction > maxfraction)) {
1252 cxdouble ksigma = 0.;
1254 cpl_matrix* base = NULL;
1255 cpl_matrix* coeffs = NULL;
1256 cpl_matrix* _coeffs = NULL;
1257 cpl_matrix* _fit = NULL;
1260 base = giraffe_chebyshev_base2d(0., 0., nx, ny,
1261 xdeg, ydeg, xbs, ybs);
1265 cpl_msg_info(fctid,
"Sigma Clipping: Error creating design "
1268 cpl_matrix_delete(zbs);
1269 cpl_matrix_delete(ybs);
1270 cpl_matrix_delete(xbs);
1277 cpl_matrix_delete(base);
1280 if (_coeffs == NULL) {
1282 cpl_msg_info(fctid,
"Sigma Clipping : Error in least square "
1283 "solution, aborting...");
1285 cpl_matrix_delete(xbs);
1286 cpl_matrix_delete(ybs);
1287 cpl_matrix_delete(zbs);
1298 coeffs = cpl_matrix_wrap(xdeg, ydeg, cpl_matrix_get_data(_coeffs));
1301 giraffe_chebyshev2d_delete(fit);
1305 fit = giraffe_chebyshev2d_new(xdeg - 1, ydeg - 1);
1306 status = giraffe_chebyshev2d_set(fit, 0., nx, 0., ny,
1311 giraffe_chebyshev2d_delete(fit);
1314 cpl_matrix_unwrap(coeffs);
1317 cpl_matrix_delete(_coeffs);
1320 cpl_matrix_delete(xbs);
1321 cpl_matrix_delete(ybs);
1322 cpl_matrix_delete(zbs);
1328 cpl_matrix_unwrap(coeffs);
1331 cpl_matrix_delete(_coeffs);
1334 _fit = cpl_matrix_new(1, cpl_matrix_get_ncol(zbs));
1336 for (j = 0; j < cpl_matrix_get_ncol(_fit); ++j) {
1338 cxdouble x = cpl_matrix_get(xbs, n, 0);
1339 cxdouble y = cpl_matrix_get(ybs, n, 0);
1340 cxdouble z = giraffe_chebyshev2d_eval(fit, x, y);
1342 cpl_matrix_set(_fit, 0, j, z);
1346 results->mean = cpl_matrix_get_mean(_fit);
1350 cpl_msg_info(fctid,
"Sigma Clipping : bias surface[%d]: "
1351 "sigma = %8.5g, ratio = %7.4g, accepted = %zu\n",
1352 iteration, sigma, fraction, naccepted);
1359 ksigma = sigma * kappa;
1363 for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
1365 register cxdouble z = cpl_matrix_get(zbs, 0, j);
1367 if (fabs(cpl_matrix_get(_fit, 0, j) - z) > ksigma) {
1371 cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
1372 cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
1373 cpl_matrix_set(zbs, 0, n, z);
1378 cpl_matrix_set_size(xbs, n, 1);
1379 cpl_matrix_set_size(ybs, n, 1);
1380 cpl_matrix_set_size(zbs, 1, n);
1382 cpl_matrix_delete(_fit);
1386 if (n == naccepted) {
1392 fraction = (cxdouble)naccepted / (cxdouble)ntotal;
1397 cpl_msg_info(fctid,
"Sigma Clipping : End");
1404 results->coeffs = cpl_matrix_duplicate(giraffe_chebyshev2d_coeffs(fit));
1405 results->rms = sigma;
1413 results->sigma = sigma / sqrt(cpl_matrix_get_ncol(zbs));
1415 cpl_msg_info(fctid,
"Sigma Clipping Results (%d/%zu, sigma = %g)",
1416 iteration, naccepted, results->rms);
1423 giraffe_chebyshev2d_delete(fit);
1426 cpl_matrix_delete(xbs);
1429 cpl_matrix_delete(ybs);
1432 cpl_matrix_delete(zbs);
1435 return EXIT_SUCCESS;
1464_giraffe_bias_compute_profile(GiBiasResults* results,
const cpl_image* image,
1465 const cpl_matrix* areas, cxchar axis)
1468 const cxchar*
const fctid =
"_giraffe_bias_compute_profile";
1482 cxdouble sigma = 0.;
1484 cpl_matrix* _areas = NULL;
1486 cpl_image* profile = NULL;
1487 cpl_image* model = NULL;
1490 cx_assert(results != NULL);
1491 cx_assert(results->limits != NULL);
1492 cx_assert(results->model == NULL);
1494 cx_assert(areas != NULL);
1496 cx_assert(image != NULL);
1497 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
1499 cx_assert((axis ==
'x') || (axis ==
'y'));
1502 nx = cpl_image_get_size_x(image);
1503 ny = cpl_image_get_size_y(image);
1510 _areas = cpl_matrix_duplicate(areas);
1511 cx_assert(_areas != NULL);
1513 nareas = cpl_matrix_get_nrow(_areas);
1517 for (j = 0; j < nareas; ++j) {
1519 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1520 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1521 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1522 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1527 cpl_matrix_set(_areas, j, 0, x_0);
1532 cpl_matrix_set(_areas, j, 1, x_1);
1537 cpl_matrix_set(_areas, j, 2, y_0);
1540 if (y_1 != ny - 1) {
1542 cpl_matrix_set(_areas, j, 3, y_1);
1545 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1546 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
1556 for (j = 0; j < nareas; ++j) {
1558 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1559 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1560 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1561 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1566 cpl_matrix_set(_areas, j, 0, x_0);
1571 cpl_matrix_set(_areas, j, 1, x_1);
1576 cpl_matrix_set(_areas, j, 2, y_0);
1581 cpl_matrix_set(_areas, j, 3, y_1);
1584 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1585 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
1595 cpl_msg_info(fctid,
"Bias Areas: Found %zu points in %dx%d image",
1604 results->sigma = 0.;
1607 cx_string_set(results->limits,
"");
1609 profile = cpl_image_new(sx, sy, CPL_TYPE_DOUBLE);
1611 for (j = 0; j < nareas; ++j) {
1613 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1614 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1615 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1616 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1618 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1620 cx_string* tmp = cx_string_new();
1623 cx_string_sprintf(tmp,
"%d:%d:%d:%d;", x_0, x_1, y_0, y_1);
1624 cx_string_append(results->limits, cx_string_get(tmp));
1626 cx_string_delete(tmp);
1644 cxint sz = x_1 - x_0 + 1;
1646 cxdouble residual = 0.;
1648 const cxdouble* _img = cpl_image_get_data_double_const(image);
1650 cxdouble* _profile = cpl_image_get_data_double(profile);
1653 for (i = y_0; i < y_1 + 1; ++i) {
1655 register cxint k = i * nx;
1656 register cxint l = 0;
1658 cxdouble m = giraffe_array_mean(&_img[k + x_0], sz);
1660 _profile[i * sx + j] = m;
1662 for (l = x_0; l < x_1 + 1; ++l) {
1663 residual += (_img[k + l] - m) * (_img[k + l] - m);
1673 rms += residual / (sz - 1);
1674 sigma += residual / (sz * (sz - 1));
1680 cxint sz = y_1 - y_0 + 1;
1682 cxdouble residual = 0.;
1684 const cxdouble* _img = cpl_image_get_data_double_const(image);
1686 cxdouble* _profile = cpl_image_get_data_double(profile);
1687 cxdouble* buffer = cx_calloc(sz,
sizeof(cxdouble));
1690 for (i = x_0; i < x_1 + 1; ++i) {
1692 register cxint l = 0;
1694 register cxdouble m = 0.;
1697 for (l = y_0; l < y_1 + 1; ++l) {
1698 buffer[l - y_0] = _img[l * nx + i];
1701 m = giraffe_array_mean(buffer, sz);
1702 _profile[i * sx + j] = m;
1704 for (l = 0; l < sz; ++l) {
1705 residual += (buffer[l] - m) * (buffer[l] - m);
1715 rms += residual / (sz - 1);
1716 sigma += residual / (sz * (sz - 1));
1727 cpl_matrix_delete(_areas);
1735 model = cpl_image_collapse_create(profile, 1);
1737 cpl_image_delete(profile);
1740 cpl_image_divide_scalar(model, nvalid);
1741 results->model = model;
1745 results->mean = cpl_image_get_mean(results->model);
1748 results->rms = sqrt(rms / (sy * nvalid));
1749 results->sigma = sqrt(sigma / sy) / nvalid;
1752 return EXIT_SUCCESS;
1782_giraffe_bias_fit_profile(GiBiasResults* results,
const cpl_image* profile,
1788 cxint sy = cpl_image_get_size_y(profile);
1789 cxint nc = order + 1;
1791 cxdouble* _profile = (cxdouble*)cpl_image_get_data_double_const(profile);
1793 cpl_matrix* base = NULL;
1794 cpl_matrix* baseT = NULL;
1795 cpl_matrix* coeff = NULL;
1796 cpl_matrix* fit = NULL;
1797 cpl_matrix* x = cpl_matrix_new(sy, 1);
1798 cpl_matrix* y = cpl_matrix_wrap(sy, 1, _profile);
1799 cpl_matrix* covy = cpl_matrix_new(sy, sy);
1800 cpl_matrix* covc = cpl_matrix_new(nc, nc);
1803 if ((x == NULL) || (y == NULL) || (covy == NULL) || (covc == NULL)) {
1805 cpl_matrix_delete(covc);
1806 cpl_matrix_delete(covy);
1807 cpl_matrix_unwrap(y);
1808 cpl_matrix_delete(x);
1819 for (i = 0; i < sy; ++i)
1821 cpl_matrix_set(x, i, 0, i);
1824 baseT = giraffe_chebyshev_base1d(0., sy, nc, x);
1825 base = cpl_matrix_transpose_create(baseT);
1827 cpl_matrix_delete(baseT);
1830 cpl_matrix_delete(x);
1834 cpl_matrix_unwrap(y);
1838 cpl_matrix_fill_diagonal(covy, results->rms * results->rms, 0);
1848 cpl_matrix_delete(covy);
1851 if (coeff == NULL) {
1853 cpl_matrix_delete(base);
1856 cpl_matrix_unwrap(y);
1864 fit = cpl_matrix_product_create(base, coeff);
1866 cpl_matrix_delete(coeff);
1869 cpl_matrix_delete(base);
1874 cpl_matrix_unwrap(y);
1886 for (i = 0; i < sy; ++i)
1888 cpl_matrix_set(y, i, 0, cpl_matrix_get(fit, i, 0));
1891 cpl_matrix_delete(fit);
1894 cpl_matrix_unwrap(y);
1904 results->mean = cpl_image_get_mean(results->model);
1905 results->sigma = sqrt(cpl_matrix_get(covc, 0, 0));
1908 cpl_image_save(results->model,
"idiot.fits", CPL_BPP_IEEE_FLOAT,
1916 cpl_matrix_delete(covc);
1919 return EXIT_SUCCESS;
1953_giraffe_bias_compute(GiBiasResults* results, cpl_image* img,
1954 const cpl_matrix* biasareas,
1955 const cpl_image* master_bias,
1956 const cpl_image* bad_pixel_map,
1957 GiBiasMethod method, GiBiasOption option,
1958 GiBiasModel model, cxbool remove_bias, cxdouble mbias,
1959 cxdouble xdeg, cxdouble ydeg, cxdouble xstep,
1960 cxdouble ystep, cxdouble sigma, cxint numiter,
1961 cxdouble maxfraction)
1964 const cxchar*
const fctid =
"giraffe_bias_compute";
1972 cxint master_bias_dimx;
1973 cxint master_bias_dimy;
1974 cxint bad_pixel_dimx;
1975 cxint bad_pixel_dimy;
1977 cxint error_code = 0;
1979 cx_string* s = NULL;
1981 const cpl_matrix* coeff = NULL;
1984 cx_assert(results != NULL);
1985 cx_assert(results->limits == NULL);
1986 cx_assert(results->coeffs == NULL);
1987 cx_assert(results->model == NULL);
1989 cx_assert(img != NULL);
1990 cx_assert(cpl_image_get_type(img) == CPL_TYPE_DOUBLE);
1992 cx_assert(biasareas != NULL);
1999 img_dimx = cpl_image_get_size_x(img);
2000 img_dimy = cpl_image_get_size_y(img);
2002 img_centerx = img_dimx / 2;
2003 img_centery = img_dimy / 2;
2005 results->limits = cx_string_new();
2013 cpl_msg_info(fctid,
"Bias correction will be done.");
2016 cpl_msg_warning(fctid,
"No Bias correction will be done!");
2020 if (model == GIBIAS_MODEL_MEAN) {
2026 cpl_msg_info(fctid,
"Using bias model 'MEAN' ...");
2028 if ((option == GIBIAS_OPTION_PLANE) ||
2029 (method == GIBIAS_METHOD_PLANE)) {
2031 cpl_msg_error(fctid,
"Can not use MEAN model with PLANE method.");
2036 error_code = _giraffe_bias_compute_mean(results, img, biasareas,
2037 sigma, numiter, maxfraction);
2041 else if ((method == GIBIAS_METHOD_CURVE) ||
2042 ((method != GIBIAS_METHOD_PROFILE) &&
2043 (option == GIBIAS_OPTION_CURVE))) {
2050 cpl_msg_info(fctid,
"Using bias model 'CURVE' ...");
2052 error_code = _giraffe_bias_compute_curve(results, img, biasareas,
2053 sigma, numiter, maxfraction, xdeg, ydeg, xstep, ystep);
2055 if (error_code != EXIT_SUCCESS) {
2057 cpl_msg_error(fctid,
"Error during calculation of bias curve, "
2063 coeff = results->coeffs;
2072 cpl_msg_info(fctid,
"Using bias model 'PLANE (FITTED)' ...");
2074 error_code = _giraffe_bias_compute_plane(results, img, biasareas,
2075 sigma, numiter, maxfraction);
2077 if (error_code == EXIT_SUCCESS) {
2079 coeff = results->coeffs;
2081 results->mean = cpl_matrix_get(coeff, 0, 0) +
2082 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2083 cpl_matrix_get(coeff, 0, 2) * img_centery;
2088 cpl_msg_error(fctid,
"Error during calculation of bias plane, "
2101 s = cx_string_new();
2102 _giraffe_method_string(s, method, option);
2104 cpl_msg_info(fctid,
"Using bias method '%s'", cx_string_get(s));
2106 cx_string_delete(s);
2111 case GIBIAS_METHOD_UNIFORM:
2118 if (model == GIBIAS_MODEL_MEAN) {
2120 cpl_msg_info(fctid,
"bias value (areas mean) = %f, "
2121 "bias sigma = %f", results->mean, results->sigma);
2126 cpl_msg_info(fctid,
"bias value at (%d, %d) = %f, "
2127 "bias sigma = %f", img_centerx, img_centery,
2128 results->mean, results->sigma);
2132 if (remove_bias == TRUE) {
2134 cpl_image_subtract_scalar(img, results->mean);
2141 case GIBIAS_METHOD_PLANE:
2144 if (coeff == NULL) {
2146 error_code = _giraffe_bias_compute_plane(results, img,
2147 biasareas, sigma, numiter, maxfraction);
2149 if (error_code == EXIT_SUCCESS) {
2151 coeff = results->coeffs;
2153 results->mean = cpl_matrix_get(coeff, 0, 0) +
2154 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2155 cpl_matrix_get(coeff, 0, 2) * img_centery;
2160 cpl_msg_error(fctid,
"Error during calculation of bias "
2161 "plane, aborting...");
2168 cpl_msg_info(fctid,
"bias value at (%d, %d) = %f, bias sigma = %f, "
2169 "bias plane = %g + %g * x + %g * y", img_centerx,
2170 img_centery, results->mean, results->sigma,
2171 cpl_matrix_get(coeff, 0, 0),
2172 cpl_matrix_get(coeff, 0, 1),
2173 cpl_matrix_get(coeff, 0, 2));
2175 if (remove_bias == TRUE) {
2179 cxdouble* _img = cpl_image_get_data_double(img);
2182 cpl_image* bsmodel = cpl_image_new(img_dimx, img_dimy,
2184 cxdouble* _bsmodel = cpl_image_get_data_double(bsmodel);
2188 for (j = 0; j < img_dimy; j++) {
2190 cxsize jj = j * img_dimx;
2192 for (i = 0; i < img_dimx; i++) {
2195 _bsmodel[jj + i] = cpl_matrix_get(coeff, 0, 0)
2196 + cpl_matrix_get(coeff, 0, 1) * i
2197 + cpl_matrix_get(coeff, 0, 2) * j;
2200 _img[jj + i] -= cpl_matrix_get(coeff, 0, 0)
2201 + cpl_matrix_get(coeff, 0, 1) * i
2202 + cpl_matrix_get(coeff, 0, 2) * j;
2209 cpl_image_save(bsmodel,
"idiot.fits", CPL_BPP_IEEE_FLOAT,
2211 cpl_image_delete(bsmodel);
2220 case GIBIAS_METHOD_CURVE:
2223 if (coeff == NULL) {
2226 _giraffe_bias_compute_curve(results, img, biasareas,
2227 sigma, numiter, maxfraction,
2228 xdeg, ydeg, xstep, ystep);
2230 if (error_code != EXIT_SUCCESS) {
2232 cpl_msg_error(fctid,
"Error during calculation of bias "
2233 "surface curve, aborting...");
2238 coeff = results->coeffs;
2242 cpl_msg_info(fctid,
"bias value %f, bias sigma = %f\n",
2243 results->mean, results->sigma);
2246 if (remove_bias == TRUE) {
2248 cpl_image* bsmodel = NULL;
2250 cpl_matrix* x = cpl_matrix_new(img_dimx * img_dimy, 1);
2251 cpl_matrix* y = cpl_matrix_new(img_dimx * img_dimy, 1);
2252 cpl_matrix* fit = NULL;
2255 for (j = 0; j < img_dimy; ++j) {
2257 register cxsize jj = j * img_dimx;
2259 for (i = 0; i < img_dimx; ++i) {
2261 cpl_matrix_set(x, jj + i, 0, i);
2262 cpl_matrix_set(y, jj + i, 0, j);
2268 fit = giraffe_chebyshev_fit2d(0., 0., img_dimx, img_dimy,
2271 cpl_matrix_delete(y);
2274 cpl_matrix_delete(x);
2277 bsmodel = cpl_image_wrap_double(img_dimx, img_dimy,
2278 cpl_matrix_get_data(fit));
2281 cpl_image_save(bsmodel,
"idiot.fits", CPL_BPP_IEEE_FLOAT,
2285 cpl_image_subtract(img, bsmodel);
2287 cpl_image_unwrap(bsmodel);
2290 cpl_matrix_delete(fit);
2299 case GIBIAS_METHOD_PROFILE:
2302 error_code = _giraffe_bias_compute_profile(results, img,
2305 if (error_code != EXIT_SUCCESS) {
2307 cpl_msg_error(fctid,
"Error computing the bias area(s) "
2313 if (option == GIBIAS_OPTION_CURVE) {
2315 error_code = _giraffe_bias_fit_profile(results, results->model,
2317 if (error_code != EXIT_SUCCESS) {
2319 cpl_msg_error(fctid,
"Error fitting the bias profile");
2326 if (remove_bias == TRUE) {
2328 const cxdouble* _bias =
2329 cpl_image_get_data_double(results->model);
2331 cxdouble* _img = cpl_image_get_data_double(img);
2334 cx_assert(_bias != NULL);
2335 cx_assert(_img != NULL);
2337 for (j = 0; j < img_dimy; ++j) {
2339 cxsize jj = j * img_dimx;
2342 for (i = 0; i < img_dimx; ++i) {
2343 _img[jj + i] -= _bias[j];
2354 case GIBIAS_METHOD_MASTER:
2355 case GIBIAS_METHOD_ZMASTER:
2358 cxdouble biasdrift = 0.;
2360 if (master_bias == NULL) {
2362 cpl_msg_error(fctid,
"Master Bias Frame required with MASTER "
2368 master_bias_dimx = cpl_image_get_size_x(master_bias);
2369 master_bias_dimy = cpl_image_get_size_y(master_bias);
2371 if ((master_bias_dimx != img_dimx) ||
2372 (master_bias_dimy != img_dimy)) {
2374 cpl_msg_error(fctid,
"Invalid master bias! Size should "
2375 "be [%d, %d] but is [%d, %d].", img_dimx, img_dimy,
2376 master_bias_dimx, master_bias_dimy);
2381 if (coeff == NULL) {
2383 if (option == GIBIAS_OPTION_CURVE) {
2385 error_code = _giraffe_bias_compute_curve(results, img,
2386 biasareas, sigma, numiter, maxfraction, xdeg,
2387 ydeg, xstep, ystep);
2389 if (error_code != EXIT_SUCCESS) {
2391 cpl_msg_error(fctid,
"Error during calculation of "
2392 "bias surface curve, aborting...");
2397 coeff = results->coeffs;
2406 error_code = _giraffe_bias_compute_plane(results, img,
2407 biasareas, sigma, numiter, maxfraction);
2409 if (error_code == EXIT_SUCCESS) {
2411 coeff = results->coeffs;
2413 results->mean = cpl_matrix_get(coeff, 0, 0) +
2414 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2415 cpl_matrix_get(coeff, 0, 2) * img_centery;
2420 cpl_msg_error(fctid,
"Error during calculation of "
2421 "bias surface plane, aborting...");
2430 if (method == GIBIAS_METHOD_ZMASTER) {
2438 biasdrift = results->mean - mbias;
2440 cpl_msg_info(fctid,
"Using bias drift value %.4g", biasdrift);
2449 if (option == GIBIAS_OPTION_CURVE) {
2453 cpl_msg_info(fctid,
"Computed bias mean = %.4f, bias "
2454 "sigma = %.4e", results->mean, results->sigma);
2459 if (option == GIBIAS_OPTION_PLANE) {
2461 cpl_msg_info(fctid,
"Bias plane = %.4e + "
2462 "%.4e * x + %.4e * y",
2463 cpl_matrix_get(coeff, 0, 0),
2464 cpl_matrix_get(coeff, 0, 1),
2465 cpl_matrix_get(coeff, 0, 2));
2466 cpl_msg_info(fctid,
"Computed bias mean = %.4f, "
2467 "bias sigma = %.4e at (%d, %d)",
2468 results->mean, results->sigma,
2469 img_centerx, img_centery);
2474 cpl_msg_info(fctid,
"Computed bias mean = %.4f, bias "
2476 results->mean, results->sigma);
2487 if (remove_bias == TRUE) {
2504 if (bad_pixel_map == NULL) {
2506 cpl_image_subtract(img, master_bias);
2508 if (biasdrift != 0.) {
2509 cpl_image_subtract_scalar(img, biasdrift);
2515 bad_pixel_dimx = cpl_image_get_size_x(bad_pixel_map);
2516 bad_pixel_dimy = cpl_image_get_size_y(bad_pixel_map);
2518 if ((bad_pixel_dimx != img_dimx) ||
2519 (bad_pixel_dimy != img_dimy)) {
2521 cpl_msg_error(fctid,
"Invalid bad pixel map size "
2522 "should be [%d, %d] but is [%d, %d].",
2524 bad_pixel_dimx, bad_pixel_dimy);
2529 if (option == GIBIAS_OPTION_CURVE) {
2531 const cxint* _bpix =
2532 cpl_image_get_data_int_const(bad_pixel_map);
2534 const cxdouble* _mbias =
2535 cpl_image_get_data_double_const(master_bias);
2537 cxdouble* _img = cpl_image_get_data_double(img);
2540 cpl_matrix_new(img_dimx * img_dimy, 1);
2542 cpl_matrix_new(img_dimx * img_dimy, 1);
2543 cpl_matrix* fit = NULL;
2546 for (j = 0; j < img_dimy; ++j) {
2548 register cxsize jj = j * img_dimx;
2550 for (i = 0; i < img_dimx; ++i) {
2552 cpl_matrix_set(x, jj + i, 0, i);
2553 cpl_matrix_set(y, jj + i, 0, j);
2558 fit = giraffe_chebyshev_fit2d(0., 0.,
2562 cpl_matrix_delete(y);
2565 cpl_matrix_delete(x);
2569 for (j = 0; j < img_dimy; ++j) {
2571 register cxsize jj = j * img_dimx;
2573 for (i = 0; i < img_dimx; ++i) {
2575 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2577 _img[jj + i] -= (_mbias[jj + i] +
2582 _img[jj + i] -= cpl_matrix_get(fit, i, j);
2589 cpl_matrix_delete(fit);
2593 else if (option == GIBIAS_OPTION_PLANE) {
2595 const cxint* _bpix =
2596 cpl_image_get_data_int_const(bad_pixel_map);
2598 const cxdouble* _mbias =
2599 cpl_image_get_data_double_const(master_bias);
2601 cxdouble* _img = cpl_image_get_data_double(img);
2604 for (j = 0; j < img_dimy; j++) {
2606 cxsize jj = j * img_dimx;
2608 for (i = 0; i < img_dimx; i++) {
2610 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2612 _img[jj + i] -= (_mbias[jj + i] +
2619 cpl_matrix_get(coeff, 0, 0) +
2620 cpl_matrix_get(coeff, 0, 1) * j +
2621 cpl_matrix_get(coeff, 0, 2) * i;
2631 const cxint* _bpix =
2632 cpl_image_get_data_int_const(bad_pixel_map);
2634 const cxdouble* _mbias =
2635 cpl_image_get_data_double_const(master_bias);
2637 cxdouble* _img = cpl_image_get_data_double(img);
2640 for (j = 0; j < img_dimy; j++) {
2642 cxsize jj = j * img_dimx;
2644 for (i = 0; i < img_dimx; i++) {
2646 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2648 _img[jj + i] -= (_mbias[jj + i] +
2654 _img[jj + i] -= results->mean;
2675 cpl_msg_error(fctid,
"Invalid method, aborting...");
2684 cpl_msg_info(fctid,
"Resulting biaslimits : %s",
2685 cx_string_get(results->limits));
2713 const cxchar *
const _id =
"giraffe_get_raw_areas";
2723 cxint32 tprescx = 0;
2724 cxint32 tprescy = 0;
2725 cxint32 tovrscx = 0;
2726 cxint32 tovrscy = 0;
2731 cpl_propertylist *properties;
2736 cpl_error_set(_id, CPL_ERROR_NULL_INPUT);
2740 winx = cpl_propertylist_get_int(properties, GIALIAS_WINX);
2741 winy = cpl_propertylist_get_int(properties, GIALIAS_WINY);
2743 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
2744 tprescx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
2749 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
2750 tprescy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
2755 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
2756 tovrscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
2761 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
2762 tovrscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
2769 m_tmp = cpl_matrix_new(1, 4);
2773 cpl_matrix_set(m_tmp, row, 0, 0.0);
2774 cpl_matrix_set(m_tmp, row, 1, (cxdouble) prescx - 1);
2775 cpl_matrix_set(m_tmp, row, 2, 0.0);
2776 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2778 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2784 cpl_matrix_set(m_tmp, row, 0, (cxdouble) winx - ovrscx);
2785 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2786 cpl_matrix_set(m_tmp, row, 2, 0.0);
2787 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2789 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2793 if (!(prescy == 0 || prescx == 0 || ovrscx == 0)) {
2795 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2796 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2797 cpl_matrix_set(m_tmp, row, 2, 0.0);
2798 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2800 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2804 else if (!(prescy == 0 || prescx == 0)) {
2806 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2807 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2808 cpl_matrix_set(m_tmp, row, 2, 0.0);
2809 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2811 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2815 else if (!(prescy == 0 || ovrscx == 0)) {
2817 cpl_matrix_set(m_tmp, row, 0, 0.0);
2818 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2819 cpl_matrix_set(m_tmp, row, 2, 0.0);
2820 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2822 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2826 else if (prescy > 0) {
2828 cpl_matrix_set(m_tmp, row, 0, 0.0);
2829 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2830 cpl_matrix_set(m_tmp, row, 2, 0.0);
2831 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2833 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2837 if (!(ovrscy == 0 || prescx == 0 || ovrscx == 0)) {
2839 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2840 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2841 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2842 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2844 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2848 else if (!(ovrscy == 0 || prescx == 0)) {
2850 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2851 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2852 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2853 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2855 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2859 else if (!(ovrscy == 0 || ovrscx == 0)) {
2861 cpl_matrix_set(m_tmp, row, 0, 0.0);
2862 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2863 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2864 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2866 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2870 else if (ovrscy > 0) {
2872 cpl_matrix_set(m_tmp, row, 0, 0.0);
2873 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2874 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2875 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2877 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2881 cpl_matrix_resize(m_tmp, 0, -1, 0, 0);
2885 cpl_matrix_delete(m_tmp);
2915 const cxchar *fctid =
"giraffe_trim_raw_areas";
2918 cxint newlx, newly, newux, newuy;
2932 cpl_msg_error(fctid,
"Image does not contain any properties!");
2936 nx = cpl_image_get_size_x(_image);
2937 ny = cpl_image_get_size_y(_image);
2939 if (cpl_propertylist_has(properties, GIALIAS_NAXIS1)) {
2940 cxint _nx = cpl_propertylist_get_int(properties, GIALIAS_NAXIS1);
2943 cpl_msg_warning(fctid,
"Image size (%d) and image property `%s' "
2944 "(%d) are not consistent! Using image size ...",
2945 nx, GIALIAS_NAXIS1, _nx);
2949 cpl_msg_warning(fctid,
"Image does not contain any property `%s'. "
2950 "Using image size (%d)", GIALIAS_NAXIS1, nx);
2954 if (cpl_propertylist_has(properties, GIALIAS_NAXIS2)) {
2955 cxint _ny = cpl_propertylist_get_int(properties, GIALIAS_NAXIS2);
2958 cpl_msg_warning(fctid,
"Image size (%d) and image property `%s' "
2959 "(%d) are not consistent! Using image size ...",
2960 ny, GIALIAS_NAXIS2, _ny);
2964 cpl_msg_warning(fctid,
"Image does not contain any property `%s'. "
2965 "Using image size (%d)", GIALIAS_NAXIS2, ny);
2968 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
2969 ovscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
2972 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
2973 ovscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
2976 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
2977 prscx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
2980 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
2981 prscy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
2993 tmp = cpl_image_extract(_image, newlx, newly, newux, newuy);
2996 cpl_image_delete(tmp);
3000 nx = cpl_image_get_size_x(_image);
3001 ny = cpl_image_get_size_y(_image);
3008 cpl_propertylist_set_int(properties, GIALIAS_NAXIS1, nx);
3009 cpl_propertylist_set_int(properties, GIALIAS_NAXIS2, ny);
3011 if (cpl_propertylist_has(properties, GIALIAS_CRPIX1)) {
3012 cxdouble crpix1 = cpl_propertylist_get_double(properties,
3020 crpix1 += (cxdouble)prscx;
3021 cpl_propertylist_set_double(properties, GIALIAS_CRPIX1, crpix1);
3024 if (cpl_propertylist_has(properties, GIALIAS_CRPIX2)) {
3025 cxdouble crpix2 = cpl_propertylist_get_double(properties,
3028 crpix2 -= (cxdouble) prscy;
3029 cpl_propertylist_set_double(properties, GIALIAS_CRPIX2, crpix2);
3032 cpl_propertylist_erase(properties, GIALIAS_OVSCX);
3033 cpl_propertylist_erase(properties, GIALIAS_OVSCY);
3034 cpl_propertylist_erase(properties, GIALIAS_PRSCX);
3035 cpl_propertylist_erase(properties, GIALIAS_PRSCY);
3105 const GiImage* master_bias,
const GiImage* bad_pixels,
3106 const cpl_matrix* biaslimits,
const GiBiasConfig* config)
3109 const cxchar*
const fctid =
"giraffe_bias_remove";
3113 cxdouble bias_value = 0.;
3117 cpl_matrix* areas = NULL;
3124 cpl_propertylist* properties;
3126 GiBiasResults bias_results = {0., 0., 0., NULL, NULL, NULL};
3129 cx_assert(raw != NULL);
3130 cx_assert(config != NULL);
3131 cx_assert(biaslimits == NULL);
3133 if (result == NULL) {
3142 if (strncmp(config->areas,
"None", 4) == 0) {
3144 cpl_msg_info(fctid,
"No bias areas specified, using pre/overscan"
3145 "regions of the raw frame.");
3149 if (areas == NULL) {
3150 if (cpl_error_get_code() == CPL_ERROR_NULL_INPUT) {
3151 cpl_msg_error(fctid,
"Invalid raw image! Image does not "
3152 "contain any properties!");
3155 cpl_msg_error(fctid,
"Invalid raw image! Image does not "
3156 "contain or has invalid pre- and overscan "
3166 areas = _giraffe_bias_get_areas(config->areas);
3168 if (areas == NULL) {
3170 cpl_msg_error(fctid,
"Invalid bias area specification '%s'!",
3176 cpl_msg_info(fctid,
"Using bias area(s) '%s' for bias computation",
3192 if (master_bias != NULL) {
3193 cpl_propertylist *_properties;
3194 cxbool is_same_size = FALSE;
3196 is_same_size = _giraffe_compare_overscans(raw, master_bias);
3203 if (is_same_size==FALSE) {
3204 cpl_msg_info(fctid,
"PRE/OVRSCAN Regions do not match between "
3205 "RAW Image and Masterbias Image.");
3212 if (cpl_propertylist_has(_properties, GIALIAS_BIASVALUE)) {
3213 bias_value = cpl_propertylist_get_double(_properties,
3224 if (bad_pixels != NULL) {
3225 cxbool is_same_size = FALSE;
3227 is_same_size = _giraffe_compare_overscans(raw, bad_pixels);
3234 if (is_same_size == FALSE) {
3235 cpl_msg_info(fctid,
"PRE/OVRSCAN Regions do not match between "
3236 "RAW Image and Bad Pixel Image.");
3251 timage = cpl_image_duplicate(_raw);
3259 error_code = _giraffe_bias_compute(&bias_results, timage,
3260 areas, _master_bias,
3261 _bad_pixels, config->method,
3262 config->option, config->model,
3263 config->remove, bias_value,
3264 config->xdeg, config->ydeg,
3265 config->xstep, config->ystep,
3266 config->sigma, config->iterations,
3269 cpl_matrix_delete(areas);
3279 _giraffe_biasresults_clear(&bias_results);
3281 cpl_msg_error(fctid,
"Bias computation failed with error %d!",
3294 if (config->remove == TRUE) {
3297 cpl_image_delete(timage);
3304 cpl_image_delete(timage);
3318 if (config->remove == TRUE) {
3320 cpl_propertylist_set_int(properties, GIALIAS_BITPIX, -32);
3321 cpl_propertylist_set_double(properties, GIALIAS_BZERO, 0.);
3322 cpl_propertylist_set_double(properties, GIALIAS_BSCALE, 1.);
3324 if (cpl_propertylist_has(properties, GIALIAS_GIRFTYPE)) {
3325 cpl_propertylist_set_string(properties,
3326 GIALIAS_GIRFTYPE,
"BRMIMG");
3329 cpl_propertylist_append_string(properties,
3330 GIALIAS_GIRFTYPE,
"BRMIMG");
3332 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3333 "GIRAFFE bias subtracted frame");
3350 if (cpl_propertylist_has(properties, GIALIAS_CONAD)) {
3351 cxdouble conad = cpl_propertylist_get_double(properties, GIALIAS_CONAD);
3357 cpl_msg_info(fctid,
"Converting bias subtracted frame to "
3358 "electrons; conversion factor is %.2f e-/ADU",
3360 cpl_image_multiply_scalar(_result, conad);
3363 cpl_msg_warning(fctid,
"Invalid conversion factor %.2f "
3364 "e-/ADU! Bias subtracted image will not be "
3365 "converted.", conad);
3369 cpl_msg_info(fctid,
"Conversion factor not found. Bias subtracted "
3370 "image will remain in ADU");
3375 s = cx_string_new();
3377 _giraffe_method_string(s, config->method, config->option);
3378 cpl_propertylist_append_string(properties, GIALIAS_BIASMETHOD, cx_string_get(s));
3380 cx_string_delete(s);
3382 cpl_propertylist_append_double(properties, GIALIAS_BIASVALUE,
3384 cpl_propertylist_append_double(properties, GIALIAS_BIASERROR,
3385 bias_results.sigma);
3386 cpl_propertylist_append_double(properties, GIALIAS_BIASSIGMA,
3389 if (bias_results.coeffs) {
3390 cpl_propertylist_append_double(properties, GIALIAS_BCLIPSIGMA,
3392 cpl_propertylist_append_int(properties, GIALIAS_BCLIPNITER,
3393 config->iterations);
3394 cpl_propertylist_append_double(properties, GIALIAS_BCLIPMFRAC,
3398 if (bias_results.limits) {
3399 cpl_propertylist_append_string(properties, GIALIAS_BIASAREAS,
3400 cx_string_get(bias_results.limits));
3403 if (bias_results.coeffs) {
3404 s = cx_string_new();
3406 _giraffe_stringify_coefficients(s, bias_results.coeffs);
3407 cpl_propertylist_append_string(properties, GIALIAS_BIASPLANE,
3410 cx_string_delete(s);
3418 _giraffe_biasresults_clear(&bias_results);
3442 GiBiasConfig* config = NULL;
3449 config = cx_calloc(1,
sizeof *config);
3456 config->method = GIBIAS_METHOD_UNDEFINED;
3457 config->option = GIBIAS_OPTION_UNDEFINED;
3458 config->model = GIBIAS_MODEL_UNDEFINED;
3464 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.remove");
3465 config->remove = cpl_parameter_get_bool(p);
3467 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.method");
3468 s = cpl_parameter_get_string(p);
3470 if (!strcmp(s,
"UNIFORM")) {
3471 config->method = GIBIAS_METHOD_UNIFORM;
3474 if (!strcmp(s,
"PLANE")) {
3475 config->method = GIBIAS_METHOD_PLANE;
3478 if (!strcmp(s,
"CURVE")) {
3479 config->method = GIBIAS_METHOD_CURVE;
3482 if (!strcmp(s,
"PROFILE")) {
3483 config->method = GIBIAS_METHOD_PROFILE;
3486 if (!strcmp(s,
"MASTER")) {
3487 config->method = GIBIAS_METHOD_MASTER;
3490 if (!strcmp(s,
"ZMASTER")) {
3491 config->method = GIBIAS_METHOD_ZMASTER;
3494 if (!strcmp(s,
"PROFILE+CURVE")) {
3495 config->method = GIBIAS_METHOD_PROFILE;
3496 config->option = GIBIAS_OPTION_CURVE;
3499 if (!strcmp(s,
"MASTER+PLANE")) {
3500 config->method = GIBIAS_METHOD_MASTER;
3501 config->option = GIBIAS_OPTION_PLANE;
3504 if (!strcmp(s,
"ZMASTER+PLANE")) {
3505 config->method = GIBIAS_METHOD_ZMASTER;
3506 config->option = GIBIAS_OPTION_PLANE;
3509 if (!strcmp(s,
"MASTER+CURVE")) {
3510 config->method = GIBIAS_METHOD_MASTER;
3511 config->option = GIBIAS_OPTION_CURVE;
3514 if (!strcmp(s,
"ZMASTER+CURVE")) {
3515 config->method = GIBIAS_METHOD_ZMASTER;
3516 config->option = GIBIAS_OPTION_CURVE;
3519 cx_assert(config->method != GIBIAS_METHOD_UNDEFINED);
3521 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.areas");
3522 s = cpl_parameter_get_string(p);
3523 config->areas = cx_strdup(s);
3525 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.sigma");
3526 config->sigma = cpl_parameter_get_double(p);
3528 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.iterations");
3529 config->iterations = cpl_parameter_get_int(p);
3531 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.fraction");
3532 config->fraction = cpl_parameter_get_double(p);
3534 if (config->method == GIBIAS_METHOD_CURVE ||
3535 config->option == GIBIAS_OPTION_CURVE) {
3536 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.xorder");
3537 config->xdeg = 1 + cpl_parameter_get_int(p);
3539 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.yorder");
3540 config->ydeg = 1 + cpl_parameter_get_int(p);
3543 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.xstep");
3544 config->xstep = cpl_parameter_get_int(p);
3546 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.ystep");
3547 config->ystep = cpl_parameter_get_int(p);
3571 if (config->areas) {
3572 cx_free(config->areas);
3573 config->areas = NULL;
3605 p = cpl_parameter_new_value(
"giraffe.biasremoval.remove",
3607 "Enable bias removal",
3608 "giraffe.biasremoval",
3610 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"remove-bias");
3611 cpl_parameterlist_append(list, p);
3614 p = cpl_parameter_new_enum(
"giraffe.biasremoval.method",
3616 "Bias removal method",
3617 "giraffe.biasremoval",
3618 "PROFILE", 11,
"UNIFORM",
"PLANE",
"CURVE",
3619 "PROFILE",
"MASTER",
"ZMASTER",
"PROFILE+CURVE",
3620 "MASTER+PLANE",
"MASTER+CURVE",
"ZMASTER+PLANE",
3622 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-method");
3623 cpl_parameterlist_append(list, p);
3626 p = cpl_parameter_new_value(
"giraffe.biasremoval.areas",
3628 "Bias areas to use "
3629 "(Xl0:Xr0:Yl0:Yr0, ... ,Xln:Xrn:Yln:Yrn)",
3630 "giraffe.biasremoval",
3632 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-areas");
3633 cpl_parameterlist_append(list, p);
3636 p = cpl_parameter_new_value(
"giraffe.biasremoval.sigma",
3638 "Sigma Clipping: sigma threshold factor",
3639 "giraffe.biasremoval",
3641 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-sigma");
3642 cpl_parameterlist_append(list, p);
3645 p = cpl_parameter_new_value(
"giraffe.biasremoval.iterations",
3647 "Sigma Clipping: maximum number of "
3649 "giraffe.biasremoval",
3651 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-niter");
3652 cpl_parameterlist_append(list, p);
3655 p = cpl_parameter_new_value(
"giraffe.biasremoval.fraction",
3657 "Sigma Clipping: minimum fraction of points "
3658 "accepted/total [0.0..1.0]",
3659 "giraffe.biasremoval",
3661 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-mfrac");
3662 cpl_parameterlist_append(list, p);
3665 p = cpl_parameter_new_value(
"giraffe.biasremoval.xorder",
3667 "Order of X polynomial fit (CURVE method "
3669 "giraffe.biasremoval",
3671 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-xorder");
3672 cpl_parameterlist_append(list, p);
3674 p = cpl_parameter_new_value(
"giraffe.biasremoval.yorder",
3676 "Order of Y polynomial fit (CURVE method "
3678 "giraffe.biasremoval",
3680 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-yorder");
3681 cpl_parameterlist_append(list, p);
3684 p = cpl_parameter_new_value(
"giraffe.biasremoval.xstep",
3686 "Sampling step along X (CURVE method only)",
3687 "giraffe.biasremoval",
3689 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-xstep");
3690 cpl_parameterlist_append(list, p);
3693 p = cpl_parameter_new_value(
"giraffe.biasremoval.ystep",
3695 "Sampling step along Y (CURVE method only)",
3696 "giraffe.biasremoval",
3698 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-ystep");
3699 cpl_parameterlist_append(list, p);
cxint giraffe_trim_raw_areas(GiImage *image)
Remove pre- and overscan ares from an image.
GiBiasConfig * giraffe_bias_config_create(cpl_parameterlist *list)
Creates a setup structure for a bias removal task.
void giraffe_bias_config_add(cpl_parameterlist *list)
Adds parameters for the bias removal.
cxint giraffe_bias_remove(GiImage *result, const GiImage *raw, const GiImage *master_bias, const GiImage *bad_pixels, const cpl_matrix *biaslimits, const GiBiasConfig *config)
Removes the bias from an image.
cpl_matrix * giraffe_get_raw_areas(const GiImage *image)
Create bias areas from an image.
void giraffe_bias_config_destroy(GiBiasConfig *config)
Destroys a bias removal setup structure.
cpl_image * giraffe_image_get(const GiImage *self)
Gets the image data.
cpl_propertylist * giraffe_image_get_properties(const GiImage *self)
Get the properties of an image.
cxint giraffe_image_set(GiImage *self, cpl_image *image)
Sets the image data.
cxint giraffe_image_set_properties(GiImage *self, cpl_propertylist *properties)
Attaches a property list to an image.
cpl_matrix * giraffe_matrix_solve_cholesky(const cpl_matrix *A, const cpl_matrix *b, const cpl_matrix *Cb, cpl_matrix *Cx)
Solve a linear system using the Cholesky decomposition.
cxdouble giraffe_matrix_sigma_fit(const cpl_matrix *matrix, const cpl_matrix *matrix_fit)
Compute sigma of matrix fit.
cpl_matrix * giraffe_matrix_leastsq(const cpl_matrix *mA, const cpl_matrix *mB)
Computes the solution of an equation using a pseudo-inverse.
cxdouble giraffe_matrix_sigma_mean(const cpl_matrix *matrix, cxdouble mean)
Compute sigma of matrix elements, with a given mean value.