37 #include <cxstrutils.h>
41 #include <cpl_parameter.h>
47 #include "gichebyshev.h"
60 struct GiBiasResults {
74 typedef struct GiBiasResults GiBiasResults;
82 _giraffe_biasresults_clear(GiBiasResults *
self)
92 cx_string_delete(self->limits);
97 cpl_matrix_delete(self->coeffs);
102 cpl_image_delete(self->model);
118 _giraffe_method_string(cx_string *
string, GiBiasMethod method,
123 case GIBIAS_METHOD_UNIFORM:
124 cx_string_set(
string,
"UNIFORM");
127 case GIBIAS_METHOD_PLANE:
128 cx_string_set(
string,
"PLANE");
131 case GIBIAS_METHOD_CURVE:
132 cx_string_set(
string,
"CURVE");
135 case GIBIAS_METHOD_PROFILE:
136 cx_string_set(
string,
"PROFILE");
139 case GIBIAS_METHOD_MASTER:
140 cx_string_set(
string,
"MASTER");
143 case GIBIAS_METHOD_ZMASTER:
144 cx_string_set(
string,
"ZMASTER");
151 if (option != GIBIAS_OPTION_UNDEFINED) {
153 case GIBIAS_OPTION_PLANE:
154 cx_string_append(
string,
"+PLANE");
157 case GIBIAS_OPTION_CURVE:
158 cx_string_append(
string,
"+CURVE");
172 _giraffe_stringify_coefficients(cx_string *
string, cpl_matrix *matrix)
177 cxchar *tmp = cx_line_alloc();
179 cxint nr = cpl_matrix_get_nrow(matrix);
180 cxint nc = cpl_matrix_get_ncol(matrix);
182 cxsize sz = cx_line_max();
184 cxdouble *data = cpl_matrix_get_data(matrix);
187 for (i = 0; i < nr; i++) {
188 for (j = 0; j < nc; j++) {
189 snprintf(tmp, sz,
"%g", data[i * nc + j]);
190 cx_string_append(
string, tmp);
192 if (i != nr - 1 || j < nc - 1) {
193 cx_string_append(
string,
":");
221 _giraffe_compare_overscans(
const GiImage* image1,
const GiImage* image2)
224 cxint32 l1ovscx = -1;
225 cxint32 l1ovscy = -1;
226 cxint32 l1prscx = -1;
227 cxint32 l1prscy = -1;
228 cxint32 l2ovscx = -1;
229 cxint32 l2ovscy = -1;
230 cxint32 l2prscx = -1;
231 cxint32 l2prscy = -1;
233 cpl_propertylist *l1, *l2;
236 cx_assert(image1 != NULL && image2 != NULL);
241 if (cpl_propertylist_has(l1, GIALIAS_OVSCX)) {
242 l1ovscx = cpl_propertylist_get_int(l1, GIALIAS_OVSCX);
244 if (cpl_propertylist_has(l1, GIALIAS_OVSCY)) {
245 l1ovscy = cpl_propertylist_get_int(l1, GIALIAS_OVSCY);
247 if (cpl_propertylist_has(l1, GIALIAS_PRSCX)) {
248 l1prscx = cpl_propertylist_get_int(l1, GIALIAS_PRSCX);
250 if (cpl_propertylist_has(l1, GIALIAS_PRSCY)) {
251 l1prscy = cpl_propertylist_get_int(l1, GIALIAS_PRSCY);
254 if (cpl_propertylist_has(l2, GIALIAS_OVSCX)) {
255 l2ovscx = cpl_propertylist_get_int(l2, GIALIAS_OVSCX);
257 if (cpl_propertylist_has(l2, GIALIAS_OVSCY)) {
258 l2ovscy = cpl_propertylist_get_int(l2, GIALIAS_OVSCY);
260 if (cpl_propertylist_has(l2, GIALIAS_PRSCX)) {
261 l2prscx = cpl_propertylist_get_int(l2, GIALIAS_PRSCX);
263 if (cpl_propertylist_has(l2, GIALIAS_PRSCY)) {
264 l2prscy = cpl_propertylist_get_int(l2, GIALIAS_PRSCY);
267 if (l1ovscx != l2ovscx || l1ovscy != l2ovscy) {
271 if (l1prscx != l2prscx || l1prscy != l2prscy) {
293 inline static cpl_matrix*
294 _giraffe_bias_get_areas(
const cxchar*
string)
297 cxchar** regions = NULL;
299 cpl_matrix* areas = NULL;
302 cx_assert(
string != NULL);
304 regions = cx_strsplit(
string,
",", -1);
306 if (regions == NULL) {
313 const cxsize nvalues = 4;
318 while (regions[i] != NULL) {
323 areas = cpl_matrix_new(nregions, nvalues);
326 while (regions[i] != NULL) {
328 register cxsize j = 0;
330 cxchar** limits = cx_strsplit(regions[i],
":", nvalues);
333 if (limits == NULL) {
335 cpl_matrix_delete(areas);
342 for (j = 0; j < nvalues; ++j) {
350 if (limits[j] == NULL || *limits[j] ==
'\0') {
357 value = strtol(limits[j], &last, 10);
363 if ((errno == ERANGE &&
364 (value == LONG_MAX || value == LONG_MIN)) ||
365 (errno != 0 && value == 0)) {
382 cpl_matrix_set(areas, i, j, value);
391 cpl_matrix_delete(areas);
394 cx_strfreev(regions);
405 cx_strfreev(regions);
462 _giraffe_bias_compute_mean(GiBiasResults* results,
const cpl_image* image,
463 const cpl_matrix* areas, cxdouble kappa, cxint numiter,
464 cxdouble maxfraction)
467 const cxchar*
const fctid =
"giraffe_bias_compute_mean";
470 const cxdouble* pdimg = NULL;
473 cxint img_dimx, img_dimy;
477 cxint x0, x1, x2, x3;
482 cxlong naccepted = 0L;
484 cxlong pixcount = 0L;
486 cxdouble currfraction = 2.;
488 cx_string* tmp = NULL;
490 cpl_matrix* matrix_zz1;
491 cpl_matrix* matrix_zz1diff;
498 if (results->limits == NULL) {
499 cpl_msg_info(fctid,
"Unable to store biaslimits return parameter, "
504 if (cpl_image_get_type(image) != CPL_TYPE_DOUBLE) {
505 cpl_msg_info(fctid,
"Only images allowed of type double, "
520 cpl_msg_info(fctid,
"Bias Areas: Missing bias areas, "
525 ba_num_cols = cpl_matrix_get_ncol(areas);
526 ba_num_rows = cpl_matrix_get_nrow(areas);
528 for (j = 0; j < ba_num_rows; j++) {
529 x3 = (cxint) cpl_matrix_get(areas, j, 3);
530 x2 = (cxint) cpl_matrix_get(areas, j, 2);
531 x1 = (cxint) cpl_matrix_get(areas, j, 1);
532 x0 = (cxint) cpl_matrix_get(areas, j, 0);
534 ntotal += (cxulong) ((x3 - x2 + 1) * (x1 - x0 + 1));
538 cpl_msg_info(fctid,
"Bias Areas: Inconsistent specification, "
543 matrix_zz1 = cpl_matrix_new(ntotal, 1);
544 matrix_zz1diff = cpl_matrix_new(ntotal, 1);
551 img_dimx = cpl_image_get_size_x(image);
552 img_dimy = cpl_image_get_size_y(image);
554 cx_string_set(results->limits,
"");
556 for (j = 0; j < ba_num_rows; j++) {
558 x3 = (cxint) cpl_matrix_get(areas, j, 3);
559 x2 = (cxint) cpl_matrix_get(areas, j, 2);
560 x1 = (cxint) cpl_matrix_get(areas, j, 1);
561 x0 = (cxint) cpl_matrix_get(areas, j, 0);
563 if ((x0 > img_dimx) || (x1 > img_dimx) ||
564 (x2 > img_dimy) || (x3 > img_dimy)) {
568 tmp = cx_string_new();
570 cx_string_sprintf(tmp,
"%d:%d:%d:%d;", x0, x1, x2, x3);
571 cx_string_append(results->limits, cx_string_get(tmp));
573 cx_string_delete(tmp);
576 pdimg = cpl_image_get_data_double_const(image);
578 for (l = x2; l < x3 + 1; l++) {
579 for (k = x0; k < x1 + 1; k++) {
580 cpl_matrix_set(matrix_zz1, n, 1, pdimg[k + l * img_dimx]);
587 cpl_msg_info(fctid,
"Bias Areas: Validation failed, aborting ...");
589 cpl_matrix_delete(matrix_zz1);
590 cpl_matrix_delete(matrix_zz1diff);
595 cpl_msg_info(fctid,
"Bias Areas: Using %s",
596 cx_string_get(results->limits));
607 cpl_msg_info(fctid,
"Sigma Clipping : Start");
611 while ((naccepted > 0) && (curriter < numiter) &&
612 (currfraction > maxfraction)) {
614 cxdouble ksigma = 0.;
616 results->mean = cpl_matrix_get_mean(matrix_zz1);
619 for (k = 0; k < cpl_matrix_get_nrow(matrix_zz1); k++) {
620 cpl_matrix_set(matrix_zz1diff, k, 0,
621 cpl_matrix_get(matrix_zz1, k, 0) - results->mean);
624 cpl_msg_info(fctid,
"bias[%d]: mean = %5g, sigma = %6.3g, "
625 "ratio = %6.3g, accepted = %ld\n", curriter,
626 results->mean, sigma, currfraction, naccepted);
628 ksigma = sigma * kappa;
631 for (l = 0; l < cpl_matrix_get_nrow(matrix_zz1); l++) {
632 if (fabs(cpl_matrix_get(matrix_zz1diff, l, 0)) > ksigma) {
636 cpl_matrix_set(matrix_zz1, pixcount, 0,
637 cpl_matrix_get(matrix_zz1, l, 0));
641 cpl_matrix_resize(matrix_zz1, 0, 0, 0, pixcount -
642 cpl_matrix_get_nrow(matrix_zz1));
644 cpl_matrix_resize(matrix_zz1diff, 0, 0, 0, pixcount -
645 cpl_matrix_get_nrow(matrix_zz1diff));
647 if (pixcount == naccepted) {
651 naccepted = pixcount;
653 currfraction = (cxdouble) naccepted / (cxdouble) ntotal;
657 cpl_msg_info(fctid,
"Sigma Clipping : End");
664 results->mean = cpl_matrix_get_mean(matrix_zz1);
667 results->sigma = results->rms / sqrt(cpl_matrix_get_nrow(matrix_zz1));
670 cpl_msg_info(fctid,
"Sigma Clipping Results : bias[%d]: mean = %5g, "
671 "sigma = %6.3g, ratio = %6.3g, accepted=%ld\n", curriter,
672 results->mean, results->rms, currfraction, naccepted);
679 if (matrix_zz1 != NULL) {
680 cpl_matrix_delete(matrix_zz1);
683 if (matrix_zz1diff != NULL) {
684 cpl_matrix_delete(matrix_zz1diff);
745 _giraffe_bias_compute_plane(GiBiasResults* results,
const cpl_image* image,
746 const cpl_matrix* areas, cxdouble kappa,
747 cxint numiter, cxdouble maxfraction)
750 const cxchar*
const fctid =
"giraffe_bias_compute_plane";
761 cxsize naccepted = 0;
763 cxdouble fraction = 1.;
766 cpl_matrix* xbs = NULL;
767 cpl_matrix* ybs = NULL;
768 cpl_matrix* zbs = NULL;
769 cpl_matrix* coeffs = NULL;
772 cx_assert(results->limits != NULL);
773 cx_assert(results->coeffs == NULL);
775 cx_assert(areas != NULL);
777 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
784 nareas = cpl_matrix_get_nrow(areas);
786 for (j = 0; j < nareas; j++) {
788 cxint x3 = (cxint) cpl_matrix_get(areas, j, 3);
789 cxint x2 = (cxint) cpl_matrix_get(areas, j, 2);
790 cxint x1 = (cxint) cpl_matrix_get(areas, j, 1);
791 cxint x0 = (cxint) cpl_matrix_get(areas, j, 0);
793 ntotal += (cxsize) ((x3 - x2 + 1) * (x1 - x0 + 1));
799 cpl_msg_info(fctid,
"Bias Areas: Inconsistent specification, "
805 nx = cpl_image_get_size_x(image);
806 ny = cpl_image_get_size_y(image);
808 cpl_msg_info(fctid,
"Bias Areas: specified are %zu points in %dx%d "
809 "image", ntotal, nx, ny);
820 cx_string_set(results->limits,
"");
822 xbs = cpl_matrix_new(ntotal, 1);
823 ybs = cpl_matrix_new(ntotal, 1);
824 zbs = cpl_matrix_new(1, ntotal);
826 for (j = 0; j < nareas; ++j) {
828 const cxdouble* _img = cpl_image_get_data_double_const(image);
832 cx_string* tmp = NULL;
835 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
836 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
837 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
838 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
840 if ((x0 > nx) || (x1 > nx) || (x2 > ny) || (x3 > ny)) {
844 tmp = cx_string_new();
846 cx_string_sprintf(tmp,
"%d:%d:%d:%d;", x0, x1, x2, x3);
847 cx_string_append(results->limits, cx_string_get(tmp));
849 cx_string_delete(tmp);
852 for (k = x2; k < x3 + 1; ++k) {
854 register cxint l = 0;
857 for (l = x0; l < x1 + 1; ++l) {
859 cpl_matrix_set(xbs, n, 0, l);
860 cpl_matrix_set(ybs, n, 0, k);
861 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
870 cpl_matrix_set_size(xbs, n, 1);
871 cpl_matrix_set_size(ybs, n, 1);
872 cpl_matrix_set_size(zbs, 1, n);
876 cpl_msg_info(fctid,
"Bias Areas: Validation failed, aborting...");
878 cpl_matrix_delete(xbs);
879 cpl_matrix_delete(ybs);
880 cpl_matrix_delete(zbs);
888 cpl_msg_info(fctid,
"Bias Areas: Using %s [%zu pixels]",
889 cx_string_get(results->limits), ntotal);
896 cpl_msg_info(fctid,
"Sigma Clipping : Start");
900 while ((naccepted > 0) && (iteration < numiter) &&
901 (fraction > maxfraction)) {
905 cxdouble ksigma = 0.;
907 cpl_matrix* base = NULL;
908 cpl_matrix* fit = NULL;
911 base = cpl_matrix_new(3, naccepted);
915 cpl_msg_info(fctid,
"Sigma Clipping: Error creating design "
918 cpl_matrix_delete(zbs);
919 cpl_matrix_delete(ybs);
920 cpl_matrix_delete(xbs);
925 for (k = 0; k < naccepted; ++k) {
927 cpl_matrix_set(base, 0, k, 1.);
928 cpl_matrix_set(base, 1, k, cpl_matrix_get(xbs, k, 0));
929 cpl_matrix_set(base, 2, k, cpl_matrix_get(ybs, k, 0));
933 cpl_matrix_delete(coeffs);
938 if (coeffs == NULL) {
940 cpl_msg_info(fctid,
"Sigma Clipping : Error in least square "
941 "solution, aborting...");
943 cpl_matrix_delete(base);
946 cpl_matrix_delete(xbs);
947 cpl_matrix_delete(ybs);
948 cpl_matrix_delete(zbs);
959 fit = cpl_matrix_product_create(coeffs, base);
961 cpl_matrix_delete(base);
964 results->mean = cpl_matrix_get_mean(fit);
968 cpl_msg_info(fctid,
"Sigma Clipping : bias plane[%d]: %g + "
969 "%g * x + %g * y, sigma = %.5g, ratio = %.4g, "
970 "accepted = %zu\n", iteration,
971 cpl_matrix_get(coeffs, 0, 0),
972 cpl_matrix_get(coeffs, 0, 1),
973 cpl_matrix_get(coeffs, 0, 2),
974 sigma, fraction, naccepted);
981 ksigma = sigma * kappa;
985 for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
987 register cxdouble z = cpl_matrix_get(zbs, 0, j);
989 if (fabs(cpl_matrix_get(fit, 0, j) - z) > ksigma) {
993 cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
995 cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
997 cpl_matrix_set(zbs, 0, n, z);
1002 cpl_matrix_set_size(xbs, n, 1);
1003 cpl_matrix_set_size(ybs, n, 1);
1004 cpl_matrix_set_size(zbs, 1, n);
1006 cpl_matrix_delete(fit);
1009 if (n == naccepted) {
1015 fraction = (cxdouble)naccepted / (cxdouble)ntotal;
1020 cpl_msg_info(fctid,
"Sigma Clipping : End");
1027 results->coeffs = coeffs;
1028 results->rms = sigma;
1036 results->sigma = sigma / sqrt(cpl_matrix_get_ncol(zbs));
1038 cpl_msg_info(fctid,
"Sigma Clipping Results (%d/%zu, sigma = %g)",
1039 iteration, naccepted, results->rms);
1046 cpl_matrix_delete(xbs);
1049 cpl_matrix_delete(ybs);
1052 cpl_matrix_delete(zbs);
1055 return EXIT_SUCCESS;
1102 _giraffe_bias_compute_curve(GiBiasResults* results,
const cpl_image* image,
1103 const cpl_matrix *areas, cxdouble kappa,
1104 cxint numiter, cxdouble maxfraction,
1105 cxdouble xdeg, cxdouble ydeg,
1106 cxdouble xstep, cxdouble ystep)
1109 const cxchar*
const fctid =
"giraffe_bias_compute_curve";
1115 cxint iteration = 0;
1119 cxsize naccepted = 0;
1121 cxdouble fraction = 1.;
1122 cxdouble sigma = 0.;
1124 cpl_matrix* xbs = NULL;
1125 cpl_matrix* ybs = NULL;
1126 cpl_matrix* zbs = NULL;
1128 GiChebyshev2D* fit = NULL;
1131 cx_assert(results != NULL);
1132 cx_assert(results->limits != NULL);
1133 cx_assert(results->coeffs == NULL);
1135 cx_assert(areas != NULL);
1137 cx_assert(image != NULL);
1138 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
1145 nareas = cpl_matrix_get_nrow(areas);
1147 for (j = 0; j < nareas; ++j) {
1149 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
1150 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
1151 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
1152 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
1154 ntotal += (cxulong) (ceil((1. + x3 - x2) / ystep) *
1155 ceil((1. + x1 - x0) / xstep));
1158 nx = cpl_image_get_size_x(image);
1159 ny = cpl_image_get_size_y(image);
1161 cpl_msg_info(fctid,
"Bias Areas: Found %zu points in %dx%d image",
1170 results->sigma = 0.;
1173 cx_string_set(results->limits,
"");
1175 xbs = cpl_matrix_new(ntotal, 1);
1176 ybs = cpl_matrix_new(ntotal, 1);
1177 zbs = cpl_matrix_new(1, ntotal);
1179 for (j = 0; j < nareas; ++j) {
1181 const cxdouble* _img = cpl_image_get_data_double_const(image);
1185 cx_string* tmp = NULL;
1188 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
1189 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
1190 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
1191 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
1193 if ((x0 > nx) || (x1 > nx) ||
1194 (x2 > ny) || (x3 > ny)) {
1198 tmp = cx_string_new();
1200 cx_string_sprintf(tmp,
"%d:%d:%d:%d;", x0, x1, x2, x3);
1201 cx_string_append(results->limits, cx_string_get(tmp));
1203 cx_string_delete(tmp);
1206 for (k = x2; k < x3 + 1; k += ystep) {
1208 register cxint l = 0;
1211 for (l = x0; l < x1 + 1; l += xstep) {
1213 cpl_matrix_set(xbs, n, 0, l);
1214 cpl_matrix_set(ybs, n, 0, k);
1215 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
1224 cpl_matrix_set_size(xbs, n, 1);
1225 cpl_matrix_set_size(ybs, n, 1);
1226 cpl_matrix_set_size(zbs, 1, n);
1230 cpl_msg_info(fctid,
"Bias Areas: Validation failed, aborting...");
1232 cpl_matrix_delete(xbs);
1233 cpl_matrix_delete(ybs);
1234 cpl_matrix_delete(zbs);
1242 cpl_msg_info(fctid,
"Bias Areas: Using %s [%zu pixels]",
1243 cx_string_get(results->limits), ntotal);
1250 cpl_msg_info(fctid,
"Sigma Clipping : Start");
1254 while ((naccepted > 0) && (iteration < numiter) &&
1255 (fraction > maxfraction)) {
1259 cxdouble ksigma = 0.;
1261 cpl_matrix* base = NULL;
1262 cpl_matrix* coeffs = NULL;
1263 cpl_matrix* _coeffs = NULL;
1264 cpl_matrix* _fit = NULL;
1267 base = giraffe_chebyshev_base2d(0., 0., nx, ny,
1268 xdeg, ydeg, xbs, ybs);
1272 cpl_msg_info(fctid,
"Sigma Clipping: Error creating design "
1275 cpl_matrix_delete(zbs);
1276 cpl_matrix_delete(ybs);
1277 cpl_matrix_delete(xbs);
1284 cpl_matrix_delete(base);
1287 if (_coeffs == NULL) {
1289 cpl_msg_info(fctid,
"Sigma Clipping : Error in least square "
1290 "solution, aborting...");
1292 cpl_matrix_delete(xbs);
1293 cpl_matrix_delete(ybs);
1294 cpl_matrix_delete(zbs);
1305 coeffs = cpl_matrix_wrap(xdeg, ydeg, cpl_matrix_get_data(_coeffs));
1308 giraffe_chebyshev2d_delete(fit);
1312 fit = giraffe_chebyshev2d_new(xdeg - 1, ydeg - 1);
1313 status = giraffe_chebyshev2d_set(fit, 0., nx, 0., ny,
1318 giraffe_chebyshev2d_delete(fit);
1321 cpl_matrix_unwrap(coeffs);
1324 cpl_matrix_delete(_coeffs);
1327 cpl_matrix_delete(xbs);
1328 cpl_matrix_delete(ybs);
1329 cpl_matrix_delete(zbs);
1335 cpl_matrix_unwrap(coeffs);
1338 cpl_matrix_delete(_coeffs);
1341 _fit = cpl_matrix_new(1, cpl_matrix_get_ncol(zbs));
1343 for (j = 0; j < cpl_matrix_get_ncol(_fit); ++j) {
1345 cxdouble x = cpl_matrix_get(xbs, n, 0);
1346 cxdouble y = cpl_matrix_get(ybs, n, 0);
1347 cxdouble z = giraffe_chebyshev2d_eval(fit, x, y);
1349 cpl_matrix_set(_fit, 0, j, z);
1353 results->mean = cpl_matrix_get_mean(_fit);
1357 cpl_msg_info(fctid,
"Sigma Clipping : bias surface[%d]: "
1358 "sigma = %8.5g, ratio = %7.4g, accepted = %zu\n",
1359 iteration, sigma, fraction, naccepted);
1366 ksigma = sigma * kappa;
1370 for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
1372 register cxdouble z = cpl_matrix_get(zbs, 0, j);
1374 if (fabs(cpl_matrix_get(_fit, 0, j) - z) > ksigma) {
1378 cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
1379 cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
1380 cpl_matrix_set(zbs, 0, n, z);
1385 cpl_matrix_set_size(xbs, n, 1);
1386 cpl_matrix_set_size(ybs, n, 1);
1387 cpl_matrix_set_size(zbs, 1, n);
1389 cpl_matrix_delete(_fit);
1393 if (n == naccepted) {
1399 fraction = (cxdouble)naccepted / (cxdouble)ntotal;
1404 cpl_msg_info(fctid,
"Sigma Clipping : End");
1411 results->coeffs = cpl_matrix_duplicate(giraffe_chebyshev2d_coeffs(fit));
1412 results->rms = sigma;
1420 results->sigma = sigma / sqrt(cpl_matrix_get_ncol(zbs));
1422 cpl_msg_info(fctid,
"Sigma Clipping Results (%d/%zu, sigma = %g)",
1423 iteration, naccepted, results->rms);
1430 giraffe_chebyshev2d_delete(fit);
1433 cpl_matrix_delete(xbs);
1436 cpl_matrix_delete(ybs);
1439 cpl_matrix_delete(zbs);
1442 return EXIT_SUCCESS;
1471 _giraffe_bias_compute_profile(GiBiasResults* results,
const cpl_image* image,
1472 const cpl_matrix* areas, cxchar axis)
1475 const cxchar*
const fctid =
"_giraffe_bias_compute_profile";
1490 cxdouble sigma = 0.;
1492 cpl_matrix* _areas = NULL;
1494 cpl_image* profile = NULL;
1495 cpl_image* model = NULL;
1498 cx_assert(results != NULL);
1499 cx_assert(results->limits != NULL);
1500 cx_assert(results->model == NULL);
1502 cx_assert(areas != NULL);
1504 cx_assert(image != NULL);
1505 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
1507 cx_assert((axis ==
'x') || (axis ==
'y'));
1510 nx = cpl_image_get_size_x(image);
1511 ny = cpl_image_get_size_y(image);
1518 _areas = cpl_matrix_duplicate(areas);
1519 cx_assert(_areas != NULL);
1521 nareas = cpl_matrix_get_nrow(_areas);
1525 for (j = 0; j < nareas; ++j) {
1527 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1528 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1529 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1530 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1535 cpl_matrix_set(_areas, j, 0, x_0);
1540 cpl_matrix_set(_areas, j, 1, x_1);
1545 cpl_matrix_set(_areas, j, 2, y_0);
1548 if (y_1 != ny - 1) {
1550 cpl_matrix_set(_areas, j, 3, y_1);
1553 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1554 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
1564 for (j = 0; j < nareas; ++j) {
1566 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1567 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1568 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1569 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1574 cpl_matrix_set(_areas, j, 0, x_0);
1579 cpl_matrix_set(_areas, j, 1, x_1);
1584 cpl_matrix_set(_areas, j, 2, y_0);
1589 cpl_matrix_set(_areas, j, 3, y_1);
1592 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1593 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
1603 cpl_msg_info(fctid,
"Bias Areas: Found %zu points in %dx%d image",
1612 results->sigma = 0.;
1615 cx_string_set(results->limits,
"");
1617 profile = cpl_image_new(sx, sy, CPL_TYPE_DOUBLE);
1619 for (j = 0; j < nareas; ++j) {
1621 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1622 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1623 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1624 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1626 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1628 cx_string* tmp = cx_string_new();
1631 cx_string_sprintf(tmp,
"%d:%d:%d:%d;", x_0, x_1, y_0, y_1);
1632 cx_string_append(results->limits, cx_string_get(tmp));
1634 cx_string_delete(tmp);
1652 cxint sz = x_1 - x_0 + 1;
1654 cxdouble residual = 0.;
1656 const cxdouble* _img = cpl_image_get_data_double_const(image);
1658 cxdouble* _profile = cpl_image_get_data_double(profile);
1661 for (i = y_0; i < y_1 + 1; ++i) {
1663 register cxint k = i * nx;
1664 register cxint l = 0;
1666 cxdouble m = giraffe_array_mean(&_img[k + x_0], sz);
1668 _profile[i * sx + j] = m;
1670 for (l = x_0; l < x_1 + 1; ++l) {
1671 residual += (_img[k + l] - m) * (_img[k + l] - m);
1682 rms += residual / (sz - 1);
1683 sigma += residual / (sz * (sz - 1));
1689 cxint sz = y_1 - y_0 + 1;
1691 cxdouble residual = 0.;
1693 const cxdouble* _img = cpl_image_get_data_double_const(image);
1695 cxdouble* _profile = cpl_image_get_data_double(profile);
1696 cxdouble* buffer = cx_calloc(sz,
sizeof(cxdouble));
1699 for (i = x_0; i < x_1 + 1; ++i) {
1701 register cxint l = 0;
1703 register cxdouble m = 0.;
1706 for (l = y_0; l < y_1 + 1; ++l) {
1707 buffer[l - y_0] = _img[l * nx + i];
1710 m = giraffe_array_mean(buffer, sz);
1711 _profile[i * sx + j] = m;
1713 for (l = 0; l < sz; ++l) {
1714 residual += (buffer[l] - m) * (buffer[l] - m);
1725 rms += residual / (sz - 1);
1726 sigma += residual / (sz * (sz - 1));
1737 cpl_matrix_delete(_areas);
1745 model = cpl_image_collapse_create(profile, 1);
1747 cpl_image_delete(profile);
1750 cpl_image_divide_scalar(model, nvalid);
1751 results->model = model;
1755 results->mean = cpl_image_get_mean(results->model);
1758 results->rms = sqrt(rms / (sy * nvalid));
1759 results->sigma = sqrt(sigma / sy) / nvalid;
1762 return EXIT_SUCCESS;
1792 _giraffe_bias_fit_profile(GiBiasResults* results,
const cpl_image* profile,
1798 cxint sy = cpl_image_get_size_y(profile);
1799 cxint nc = order + 1;
1801 cxdouble* _profile = (cxdouble*)cpl_image_get_data_double_const(profile);
1803 cpl_matrix* base = NULL;
1804 cpl_matrix* baseT = NULL;
1805 cpl_matrix* coeff = NULL;
1806 cpl_matrix* fit = NULL;
1807 cpl_matrix* x = cpl_matrix_new(sy, 1);
1808 cpl_matrix* y = cpl_matrix_wrap(sy, 1, _profile);
1809 cpl_matrix* covy = cpl_matrix_new(sy, sy);
1810 cpl_matrix* covc = cpl_matrix_new(nc, nc);
1813 if ((x == NULL) || (y == NULL) || (covy == NULL) || (covc == NULL)) {
1815 cpl_matrix_delete(covc);
1816 cpl_matrix_delete(covy);
1817 cpl_matrix_unwrap(y);
1818 cpl_matrix_delete(x);
1829 for (i = 0; i < sy; ++i)
1831 cpl_matrix_set(x, i, 0, i);
1834 baseT = giraffe_chebyshev_base1d(0., sy, nc, x);
1835 base = cpl_matrix_transpose_create(baseT);
1837 cpl_matrix_delete(baseT);
1840 cpl_matrix_delete(x);
1844 cpl_matrix_unwrap(y);
1848 cpl_matrix_fill_diagonal(covy, results->rms * results->rms, 0);
1858 cpl_matrix_delete(covy);
1861 if (coeff == NULL) {
1863 cpl_matrix_delete(base);
1866 cpl_matrix_unwrap(y);
1874 fit = cpl_matrix_product_create(base, coeff);
1876 cpl_matrix_delete(coeff);
1879 cpl_matrix_delete(base);
1884 cpl_matrix_unwrap(y);
1896 for (i = 0; i < sy; ++i)
1898 cpl_matrix_set(y, i, 0, cpl_matrix_get(fit, i, 0));
1901 cpl_matrix_delete(fit);
1904 cpl_matrix_unwrap(y);
1914 results->mean = cpl_image_get_mean(results->model);
1915 results->sigma = sqrt(cpl_matrix_get(covc, 0, 0));
1918 cpl_image_save(results->model,
"idiot.fits", CPL_BPP_IEEE_FLOAT,
1926 cpl_matrix_delete(covc);
1929 return EXIT_SUCCESS;
1963 _giraffe_bias_compute(GiBiasResults* results, cpl_image* img,
1964 const cpl_matrix* biasareas,
1965 const cpl_image* master_bias,
1966 const cpl_image* bad_pixel_map,
1967 GiBiasMethod method, GiBiasOption option,
1968 GiBiasModel model, cxbool remove_bias, cxdouble mbias,
1969 cxdouble xdeg, cxdouble ydeg, cxdouble xstep,
1970 cxdouble ystep, cxdouble sigma, cxint numiter,
1971 cxdouble maxfraction)
1974 const cxchar*
const fctid =
"giraffe_bias_compute";
1982 cxint master_bias_dimx;
1983 cxint master_bias_dimy;
1984 cxint bad_pixel_dimx;
1985 cxint bad_pixel_dimy;
1987 cxint error_code = 0;
1989 cx_string* s = NULL;
1991 const cpl_matrix* coeff = NULL;
1994 cx_assert(results != NULL);
1995 cx_assert(results->limits == NULL);
1996 cx_assert(results->coeffs == NULL);
1997 cx_assert(results->model == NULL);
1999 cx_assert(img != NULL);
2000 cx_assert(cpl_image_get_type(img) == CPL_TYPE_DOUBLE);
2002 cx_assert(biasareas != NULL);
2009 img_dimx = cpl_image_get_size_x(img);
2010 img_dimy = cpl_image_get_size_y(img);
2012 img_centerx = img_dimx / 2;
2013 img_centery = img_dimy / 2;
2015 results->limits = cx_string_new();
2023 cpl_msg_info(fctid,
"Bias correction will be done.");
2026 cpl_msg_warning(fctid,
"No Bias correction will be done!");
2030 if (model == GIBIAS_MODEL_MEAN) {
2036 cpl_msg_info(fctid,
"Using bias model 'MEAN' ...");
2038 if ((option == GIBIAS_OPTION_PLANE) ||
2039 (method == GIBIAS_METHOD_PLANE)) {
2041 cpl_msg_error(fctid,
"Can not use MEAN model with PLANE method.");
2046 error_code = _giraffe_bias_compute_mean(results, img, biasareas,
2047 sigma, numiter, maxfraction);
2051 else if ((method == GIBIAS_METHOD_CURVE) ||
2052 ((method != GIBIAS_METHOD_PROFILE) &&
2053 (option == GIBIAS_OPTION_CURVE))) {
2060 cpl_msg_info(fctid,
"Using bias model 'CURVE' ...");
2062 error_code = _giraffe_bias_compute_curve(results, img, biasareas,
2063 sigma, numiter, maxfraction, xdeg, ydeg, xstep, ystep);
2065 if (error_code != EXIT_SUCCESS) {
2067 cpl_msg_error(fctid,
"Error during calculation of bias curve, "
2073 coeff = results->coeffs;
2082 cpl_msg_info(fctid,
"Using bias model 'PLANE (FITTED)' ...");
2084 error_code = _giraffe_bias_compute_plane(results, img, biasareas,
2085 sigma, numiter, maxfraction);
2087 if (error_code == EXIT_SUCCESS) {
2089 coeff = results->coeffs;
2091 results->mean = cpl_matrix_get(coeff, 0, 0) +
2092 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2093 cpl_matrix_get(coeff, 0, 2) * img_centery;
2098 cpl_msg_error(fctid,
"Error during calculation of bias plane, "
2111 s = cx_string_new();
2112 _giraffe_method_string(s, method, option);
2114 cpl_msg_info(fctid,
"Using bias method '%s'", cx_string_get(s));
2116 cx_string_delete(s);
2121 case GIBIAS_METHOD_UNIFORM:
2128 if (model == GIBIAS_MODEL_MEAN) {
2130 cpl_msg_info(fctid,
"bias value (areas mean) = %f, "
2131 "bias sigma = %f", results->mean, results->sigma);
2136 cpl_msg_info(fctid,
"bias value at (%d, %d) = %f, "
2137 "bias sigma = %f", img_centerx, img_centery,
2138 results->mean, results->sigma);
2142 if (remove_bias == TRUE) {
2144 cpl_image_subtract_scalar(img, results->mean);
2151 case GIBIAS_METHOD_PLANE:
2154 if (coeff == NULL) {
2156 error_code = _giraffe_bias_compute_plane(results, img,
2157 biasareas, sigma, numiter, maxfraction);
2159 if (error_code == EXIT_SUCCESS) {
2161 coeff = results->coeffs;
2163 results->mean = cpl_matrix_get(coeff, 0, 0) +
2164 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2165 cpl_matrix_get(coeff, 0, 2) * img_centery;
2170 cpl_msg_error(fctid,
"Error during calculation of bias "
2171 "plane, aborting...");
2178 cpl_msg_info(fctid,
"bias value at (%d, %d) = %f, bias sigma = %f, "
2179 "bias plane = %g + %g * x + %g * y", img_centerx,
2180 img_centery, results->mean, results->sigma,
2181 cpl_matrix_get(coeff, 0, 0),
2182 cpl_matrix_get(coeff, 0, 1),
2183 cpl_matrix_get(coeff, 0, 2));
2185 if (remove_bias == TRUE) {
2189 cxdouble* _img = cpl_image_get_data_double(img);
2192 cpl_image* bsmodel = cpl_image_new(img_dimx, img_dimy,
2194 cxdouble* _bsmodel = cpl_image_get_data_double(bsmodel);
2198 for (j = 0; j < img_dimy; j++) {
2200 cxsize jj = j * img_dimx;
2202 for (i = 0; i < img_dimx; i++) {
2205 _bsmodel[jj + i] = cpl_matrix_get(coeff, 0, 0)
2206 + cpl_matrix_get(coeff, 0, 1) * i
2207 + cpl_matrix_get(coeff, 0, 2) * j;
2210 _img[jj + i] -= cpl_matrix_get(coeff, 0, 0)
2211 + cpl_matrix_get(coeff, 0, 1) * i
2212 + cpl_matrix_get(coeff, 0, 2) * j;
2219 cpl_image_save(bsmodel,
"idiot.fits", CPL_BPP_IEEE_FLOAT,
2221 cpl_image_delete(bsmodel);
2230 case GIBIAS_METHOD_CURVE:
2233 if (coeff == NULL) {
2236 _giraffe_bias_compute_curve(results, img, biasareas,
2237 sigma, numiter, maxfraction,
2238 xdeg, ydeg, xstep, ystep);
2240 if (error_code != EXIT_SUCCESS) {
2242 cpl_msg_error(fctid,
"Error during calculation of bias "
2243 "surface curve, aborting...");
2248 coeff = results->coeffs;
2252 cpl_msg_info(fctid,
"bias value %f, bias sigma = %f\n",
2253 results->mean, results->sigma);
2256 if (remove_bias == TRUE) {
2258 cpl_image* bsmodel = NULL;
2260 cpl_matrix* x = cpl_matrix_new(img_dimx * img_dimy, 1);
2261 cpl_matrix* y = cpl_matrix_new(img_dimx * img_dimy, 1);
2262 cpl_matrix* fit = NULL;
2265 for (j = 0; j < img_dimy; ++j) {
2267 register cxsize jj = j * img_dimx;
2269 for (i = 0; i < img_dimx; ++i) {
2271 cpl_matrix_set(x, jj + i, 0, i);
2272 cpl_matrix_set(y, jj + i, 0, j);
2278 fit = giraffe_chebyshev_fit2d(0., 0., img_dimx, img_dimy,
2281 cpl_matrix_delete(y);
2284 cpl_matrix_delete(x);
2287 bsmodel = cpl_image_wrap_double(img_dimx, img_dimy,
2288 cpl_matrix_get_data(fit));
2291 cpl_image_save(bsmodel,
"idiot.fits", CPL_BPP_IEEE_FLOAT,
2295 cpl_image_subtract(img, bsmodel);
2297 cpl_image_unwrap(bsmodel);
2300 cpl_matrix_delete(fit);
2309 case GIBIAS_METHOD_PROFILE:
2312 error_code = _giraffe_bias_compute_profile(results, img,
2315 if (error_code != EXIT_SUCCESS) {
2317 cpl_msg_error(fctid,
"Error computing the bias area(s) "
2323 if (option == GIBIAS_OPTION_CURVE) {
2325 error_code = _giraffe_bias_fit_profile(results, results->model,
2327 if (error_code != EXIT_SUCCESS) {
2329 cpl_msg_error(fctid,
"Error fitting the bias profile");
2336 if (remove_bias == TRUE) {
2338 const cxdouble* _bias =
2339 cpl_image_get_data_double(results->model);
2341 cxdouble* _img = cpl_image_get_data_double(img);
2344 cx_assert(_bias != NULL);
2345 cx_assert(_img != NULL);
2347 for (j = 0; j < img_dimy; ++j) {
2349 cxsize jj = j * img_dimx;
2352 for (i = 0; i < img_dimx; ++i) {
2353 _img[jj + i] -= _bias[j];
2364 case GIBIAS_METHOD_MASTER:
2365 case GIBIAS_METHOD_ZMASTER:
2368 cxdouble biasdrift = 0.;
2370 if (master_bias == NULL) {
2372 cpl_msg_error(fctid,
"Master Bias Frame required with MASTER "
2378 master_bias_dimx = cpl_image_get_size_x(master_bias);
2379 master_bias_dimy = cpl_image_get_size_y(master_bias);
2381 if ((master_bias_dimx != img_dimx) ||
2382 (master_bias_dimy != img_dimy)) {
2384 cpl_msg_error(fctid,
"Invalid master bias! Size should "
2385 "be [%d, %d] but is [%d, %d].", img_dimx, img_dimy,
2386 master_bias_dimx, master_bias_dimy);
2391 if (coeff == NULL) {
2393 if (option == GIBIAS_OPTION_CURVE) {
2395 error_code = _giraffe_bias_compute_curve(results, img,
2396 biasareas, sigma, numiter, maxfraction, xdeg,
2397 ydeg, xstep, ystep);
2399 if (error_code != EXIT_SUCCESS) {
2401 cpl_msg_error(fctid,
"Error during calculation of "
2402 "bias surface curve, aborting...");
2407 coeff = results->coeffs;
2416 error_code = _giraffe_bias_compute_plane(results, img,
2417 biasareas, sigma, numiter, maxfraction);
2419 if (error_code == EXIT_SUCCESS) {
2421 coeff = results->coeffs;
2423 results->mean = cpl_matrix_get(coeff, 0, 0) +
2424 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2425 cpl_matrix_get(coeff, 0, 2) * img_centery;
2430 cpl_msg_error(fctid,
"Error during calculation of "
2431 "bias surface plane, aborting...");
2440 if (method == GIBIAS_METHOD_ZMASTER) {
2448 biasdrift = results->mean - mbias;
2450 cpl_msg_info(fctid,
"Using bias drift value %.4g", biasdrift);
2459 if (option == GIBIAS_OPTION_CURVE) {
2463 cpl_msg_info(fctid,
"Computed bias mean = %.4f, bias "
2464 "sigma = %.4e", results->mean, results->sigma);
2469 if (option == GIBIAS_OPTION_PLANE) {
2471 cpl_msg_info(fctid,
"Bias plane = %.4e + "
2472 "%.4e * x + %.4e * y",
2473 cpl_matrix_get(coeff, 0, 0),
2474 cpl_matrix_get(coeff, 0, 1),
2475 cpl_matrix_get(coeff, 0, 2));
2476 cpl_msg_info(fctid,
"Computed bias mean = %.4f, "
2477 "bias sigma = %.4e at (%d, %d)",
2478 results->mean, results->sigma,
2479 img_centerx, img_centery);
2484 cpl_msg_info(fctid,
"Computed bias mean = %.4f, bias "
2486 results->mean, results->sigma);
2497 if (remove_bias == TRUE) {
2514 if (bad_pixel_map == NULL) {
2516 cpl_image_subtract(img, master_bias);
2518 if (biasdrift != 0.) {
2519 cpl_image_subtract_scalar(img, biasdrift);
2525 bad_pixel_dimx = cpl_image_get_size_x(bad_pixel_map);
2526 bad_pixel_dimy = cpl_image_get_size_y(bad_pixel_map);
2528 if ((bad_pixel_dimx != img_dimx) ||
2529 (bad_pixel_dimy != img_dimy)) {
2531 cpl_msg_error(fctid,
"Invalid bad pixel map size "
2532 "should be [%d, %d] but is [%d, %d].",
2534 bad_pixel_dimx, bad_pixel_dimy);
2539 if (option == GIBIAS_OPTION_CURVE) {
2541 const cxint* _bpix =
2542 cpl_image_get_data_int_const(bad_pixel_map);
2544 const cxdouble* _mbias =
2545 cpl_image_get_data_double_const(master_bias);
2547 cxdouble* _img = cpl_image_get_data_double(img);
2550 cpl_matrix_new(img_dimx * img_dimy, 1);
2552 cpl_matrix_new(img_dimx * img_dimy, 1);
2553 cpl_matrix* fit = NULL;
2556 for (j = 0; j < img_dimy; ++j) {
2558 register cxsize jj = j * img_dimx;
2560 for (i = 0; i < img_dimx; ++i) {
2562 cpl_matrix_set(x, jj + i, 0, i);
2563 cpl_matrix_set(y, jj + i, 0, j);
2568 fit = giraffe_chebyshev_fit2d(0., 0.,
2572 cpl_matrix_delete(y);
2575 cpl_matrix_delete(x);
2579 for (j = 0; j < img_dimy; ++j) {
2581 register cxsize jj = j * img_dimx;
2583 for (i = 0; i < img_dimx; ++i) {
2585 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2587 _img[jj + i] -= (_mbias[jj + i] +
2592 _img[jj + i] -= cpl_matrix_get(fit, i, j);
2599 cpl_matrix_delete(fit);
2603 else if (option == GIBIAS_OPTION_PLANE) {
2605 const cxint* _bpix =
2606 cpl_image_get_data_int_const(bad_pixel_map);
2608 const cxdouble* _mbias =
2609 cpl_image_get_data_double_const(master_bias);
2611 cxdouble* _img = cpl_image_get_data_double(img);
2614 for (j = 0; j < img_dimy; j++) {
2616 cxsize jj = j * img_dimx;
2618 for (i = 0; i < img_dimx; i++) {
2620 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2622 _img[jj + i] -= (_mbias[jj + i] +
2629 cpl_matrix_get(coeff, 0, 0) +
2630 cpl_matrix_get(coeff, 0, 1) * j +
2631 cpl_matrix_get(coeff, 0, 2) * i;
2641 const cxint* _bpix =
2642 cpl_image_get_data_int_const(bad_pixel_map);
2644 const cxdouble* _mbias =
2645 cpl_image_get_data_double_const(master_bias);
2647 cxdouble* _img = cpl_image_get_data_double(img);
2650 for (j = 0; j < img_dimy; j++) {
2652 cxsize jj = j * img_dimx;
2654 for (i = 0; i < img_dimx; i++) {
2656 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2658 _img[jj + i] -= (_mbias[jj + i] +
2664 _img[jj + i] -= results->mean;
2685 cpl_msg_error(fctid,
"Invalid method, aborting...");
2694 cpl_msg_info(fctid,
"Resulting biaslimits : %s",
2695 cx_string_get(results->limits));
2723 const cxchar *
const _id =
"giraffe_get_raw_areas";
2733 cxint32 tprescx = 0;
2734 cxint32 tprescy = 0;
2735 cxint32 tovrscx = 0;
2736 cxint32 tovrscy = 0;
2741 cpl_propertylist *properties;
2746 cpl_error_set(_id, CPL_ERROR_NULL_INPUT);
2750 winx = cpl_propertylist_get_int(properties, GIALIAS_WINX);
2751 winy = cpl_propertylist_get_int(properties, GIALIAS_WINY);
2753 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
2754 tprescx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
2759 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
2760 tprescy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
2765 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
2766 tovrscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
2771 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
2772 tovrscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
2779 m_tmp = cpl_matrix_new(1, 4);
2783 cpl_matrix_set(m_tmp, row, 0, 0.0);
2784 cpl_matrix_set(m_tmp, row, 1, (cxdouble) prescx - 1);
2785 cpl_matrix_set(m_tmp, row, 2, 0.0);
2786 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2788 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2794 cpl_matrix_set(m_tmp, row, 0, (cxdouble) winx - ovrscx);
2795 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2796 cpl_matrix_set(m_tmp, row, 2, 0.0);
2797 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2799 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2803 if (!(prescy == 0 || prescx == 0 || ovrscx == 0)) {
2805 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2806 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2807 cpl_matrix_set(m_tmp, row, 2, 0.0);
2808 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2810 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2814 else if (!(prescy == 0 || prescx == 0)) {
2816 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2817 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2818 cpl_matrix_set(m_tmp, row, 2, 0.0);
2819 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2821 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2825 else if (!(prescy == 0 || ovrscx == 0)) {
2827 cpl_matrix_set(m_tmp, row, 0, 0.0);
2828 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2829 cpl_matrix_set(m_tmp, row, 2, 0.0);
2830 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2832 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2836 else if (prescy > 0) {
2838 cpl_matrix_set(m_tmp, row, 0, 0.0);
2839 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2840 cpl_matrix_set(m_tmp, row, 2, 0.0);
2841 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2843 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2847 if (!(ovrscy == 0 || prescx == 0 || ovrscx == 0)) {
2849 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2850 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2851 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2852 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2854 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2858 else if (!(ovrscy == 0 || prescx == 0)) {
2860 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2861 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2862 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2863 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2865 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2869 else if (!(ovrscy == 0 || ovrscx == 0)) {
2871 cpl_matrix_set(m_tmp, row, 0, 0.0);
2872 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2873 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2874 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2876 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2880 else if (ovrscy > 0) {
2882 cpl_matrix_set(m_tmp, row, 0, 0.0);
2883 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2884 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2885 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2887 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2891 cpl_matrix_resize(m_tmp, 0, -1, 0, 0);
2895 cpl_matrix_delete(m_tmp);
2925 const cxchar *fctid =
"giraffe_trim_raw_areas";
2928 cxint newlx, newly, newux, newuy;
2942 cpl_msg_error(fctid,
"Image does not contain any properties!");
2946 nx = cpl_image_get_size_x(_image);
2947 ny = cpl_image_get_size_y(_image);
2949 if (cpl_propertylist_has(properties, GIALIAS_NAXIS1)) {
2950 cxint _nx = cpl_propertylist_get_int(properties, GIALIAS_NAXIS1);
2953 cpl_msg_warning(fctid,
"Image size (%d) and image property `%s' "
2954 "(%d) are not consistent! Using image size ...",
2955 nx, GIALIAS_NAXIS1, _nx);
2959 cpl_msg_warning(fctid,
"Image does not contain any property `%s'. "
2960 "Using image size (%d)", GIALIAS_NAXIS1, nx);
2964 if (cpl_propertylist_has(properties, GIALIAS_NAXIS2)) {
2965 cxint _ny = cpl_propertylist_get_int(properties, GIALIAS_NAXIS2);
2968 cpl_msg_warning(fctid,
"Image size (%d) and image property `%s' "
2969 "(%d) are not consistent! Using image size ...",
2970 ny, GIALIAS_NAXIS2, _ny);
2974 cpl_msg_warning(fctid,
"Image does not contain any property `%s'. "
2975 "Using image size (%d)", GIALIAS_NAXIS2, ny);
2978 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
2979 ovscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
2982 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
2983 ovscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
2986 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
2987 prscx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
2990 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
2991 prscy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
3003 tmp = cpl_image_extract(_image, newlx, newly, newux, newuy);
3006 cpl_image_delete(tmp);
3010 nx = cpl_image_get_size_x(_image);
3011 ny = cpl_image_get_size_y(_image);
3018 cpl_propertylist_set_int(properties, GIALIAS_NAXIS1, nx);
3019 cpl_propertylist_set_int(properties, GIALIAS_NAXIS2, ny);
3021 if (cpl_propertylist_has(properties, GIALIAS_CRPIX1)) {
3022 cxdouble crpix1 = cpl_propertylist_get_double(properties,
3030 crpix1 += (cxdouble)prscx;
3031 cpl_propertylist_set_double(properties, GIALIAS_CRPIX1, crpix1);
3034 if (cpl_propertylist_has(properties, GIALIAS_CRPIX2)) {
3035 cxdouble crpix2 = cpl_propertylist_get_double(properties,
3038 crpix2 -= (cxdouble) prscy;
3039 cpl_propertylist_set_double(properties, GIALIAS_CRPIX2, crpix2);
3042 cpl_propertylist_erase(properties, GIALIAS_OVSCX);
3043 cpl_propertylist_erase(properties, GIALIAS_OVSCY);
3044 cpl_propertylist_erase(properties, GIALIAS_PRSCX);
3045 cpl_propertylist_erase(properties, GIALIAS_PRSCY);
3115 const GiImage* master_bias,
const GiImage* bad_pixels,
3116 const cpl_matrix* biaslimits,
const GiBiasConfig* config)
3119 const cxchar*
const fctid =
"giraffe_bias_remove";
3123 cxdouble bias_value = 0.;
3127 cpl_matrix* areas = NULL;
3134 cpl_propertylist* properties;
3136 GiBiasResults bias_results = {0., 0., 0., NULL, NULL, NULL};
3139 cx_assert(raw != NULL);
3140 cx_assert(config != NULL);
3141 cx_assert(biaslimits == NULL);
3143 if (result == NULL) {
3152 if (strncmp(config->areas,
"None", 4) == 0) {
3154 cpl_msg_info(fctid,
"No bias areas specified, using pre/overscan"
3155 "regions of the raw frame.");
3159 if (areas == NULL) {
3160 if (cpl_error_get_code() == CPL_ERROR_NULL_INPUT) {
3161 cpl_msg_error(fctid,
"Invalid raw image! Image does not "
3162 "contain any properties!");
3165 cpl_msg_error(fctid,
"Invalid raw image! Image does not "
3166 "contain or has invalid pre- and overscan "
3176 areas = _giraffe_bias_get_areas(config->areas);
3178 if (areas == NULL) {
3180 cpl_msg_error(fctid,
"Invalid bias area specification '%s'!",
3186 cpl_msg_info(fctid,
"Using bias area(s) '%s' for bias computation",
3202 if (master_bias != NULL) {
3203 cpl_propertylist *_properties;
3204 cxbool is_same_size = FALSE;
3206 is_same_size = _giraffe_compare_overscans(raw, master_bias);
3213 if (is_same_size==FALSE) {
3214 cpl_msg_info(fctid,
"PRE/OVRSCAN Regions do not match between "
3215 "RAW Image and Masterbias Image.");
3222 if (cpl_propertylist_has(_properties, GIALIAS_BIASVALUE)) {
3223 bias_value = cpl_propertylist_get_double(_properties,
3234 if (bad_pixels != NULL) {
3235 cxbool is_same_size = FALSE;
3237 is_same_size = _giraffe_compare_overscans(raw, bad_pixels);
3244 if (is_same_size == FALSE) {
3245 cpl_msg_info(fctid,
"PRE/OVRSCAN Regions do not match between "
3246 "RAW Image and Bad Pixel Image.");
3261 timage = cpl_image_duplicate(_raw);
3269 error_code = _giraffe_bias_compute(&bias_results, timage,
3270 areas, _master_bias,
3271 _bad_pixels, config->method,
3272 config->option, config->model,
3273 config->remove, bias_value,
3274 config->xdeg, config->ydeg,
3275 config->xstep, config->ystep,
3276 config->sigma, config->iterations,
3279 cpl_matrix_delete(areas);
3289 _giraffe_biasresults_clear(&bias_results);
3291 cpl_msg_error(fctid,
"Bias computation failed with error %d!",
3304 if (config->remove == TRUE) {
3307 cpl_image_delete(timage);
3314 cpl_image_delete(timage);
3328 if (config->remove == TRUE) {
3330 cpl_propertylist_set_int(properties, GIALIAS_BITPIX, -32);
3331 cpl_propertylist_set_double(properties, GIALIAS_BZERO, 0.);
3332 cpl_propertylist_set_double(properties, GIALIAS_BSCALE, 1.);
3334 if (cpl_propertylist_has(properties, GIALIAS_GIRFTYPE)) {
3335 cpl_propertylist_set_string(properties,
3336 GIALIAS_GIRFTYPE,
"BRMIMG");
3339 cpl_propertylist_append_string(properties,
3340 GIALIAS_GIRFTYPE,
"BRMIMG");
3342 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3343 "GIRAFFE bias subtracted frame");
3359 if (cpl_propertylist_has(properties, GIALIAS_CONAD)) {
3360 cxdouble conad = cpl_propertylist_get_double(properties, GIALIAS_CONAD);
3366 cpl_msg_info(fctid,
"Converting bias subtracted frame to "
3367 "electrons; conversion factor is %.2f e-/ADU",
3369 cpl_image_multiply_scalar(_result, conad);
3372 cpl_msg_warning(fctid,
"Invalid conversion factor %.2f "
3373 "e-/ADU! Bias subtracted image will not be "
3374 "converted.", conad);
3378 cpl_msg_info(fctid,
"Conversion factor not found. Bias subtracted "
3379 "image will remain in ADU");
3384 s = cx_string_new();
3386 _giraffe_method_string(s, config->method, config->option);
3387 cpl_propertylist_append_string(properties, GIALIAS_BIASMETHOD, cx_string_get(s));
3389 cx_string_delete(s);
3391 cpl_propertylist_append_double(properties, GIALIAS_BIASVALUE,
3393 cpl_propertylist_append_double(properties, GIALIAS_BIASERROR,
3394 bias_results.sigma);
3395 cpl_propertylist_append_double(properties, GIALIAS_BIASSIGMA,
3398 if (bias_results.coeffs) {
3399 cpl_propertylist_append_double(properties, GIALIAS_BCLIPSIGMA,
3401 cpl_propertylist_append_int(properties, GIALIAS_BCLIPNITER,
3402 config->iterations);
3403 cpl_propertylist_append_double(properties, GIALIAS_BCLIPMFRAC,
3407 if (bias_results.limits) {
3408 cpl_propertylist_append_string(properties, GIALIAS_BIASAREAS,
3409 cx_string_get(bias_results.limits));
3412 if (bias_results.coeffs) {
3413 s = cx_string_new();
3415 _giraffe_stringify_coefficients(s, bias_results.coeffs);
3416 cpl_propertylist_append_string(properties, GIALIAS_BIASPLANE,
3419 cx_string_delete(s);
3427 _giraffe_biasresults_clear(&bias_results);
3451 GiBiasConfig* config = NULL;
3458 config = cx_calloc(1,
sizeof *config);
3465 config->method = GIBIAS_METHOD_UNDEFINED;
3466 config->option = GIBIAS_OPTION_UNDEFINED;
3467 config->model = GIBIAS_MODEL_UNDEFINED;
3473 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.remove");
3474 config->remove = cpl_parameter_get_bool(p);
3476 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.method");
3477 s = cpl_parameter_get_string(p);
3479 if (!strcmp(s,
"UNIFORM")) {
3480 config->method = GIBIAS_METHOD_UNIFORM;
3483 if (!strcmp(s,
"PLANE")) {
3484 config->method = GIBIAS_METHOD_PLANE;
3487 if (!strcmp(s,
"CURVE")) {
3488 config->method = GIBIAS_METHOD_CURVE;
3491 if (!strcmp(s,
"PROFILE")) {
3492 config->method = GIBIAS_METHOD_PROFILE;
3495 if (!strcmp(s,
"MASTER")) {
3496 config->method = GIBIAS_METHOD_MASTER;
3499 if (!strcmp(s,
"ZMASTER")) {
3500 config->method = GIBIAS_METHOD_ZMASTER;
3503 if (!strcmp(s,
"PROFILE+CURVE")) {
3504 config->method = GIBIAS_METHOD_PROFILE;
3505 config->option = GIBIAS_OPTION_CURVE;
3508 if (!strcmp(s,
"MASTER+PLANE")) {
3509 config->method = GIBIAS_METHOD_MASTER;
3510 config->option = GIBIAS_OPTION_PLANE;
3513 if (!strcmp(s,
"ZMASTER+PLANE")) {
3514 config->method = GIBIAS_METHOD_ZMASTER;
3515 config->option = GIBIAS_OPTION_PLANE;
3518 if (!strcmp(s,
"MASTER+CURVE")) {
3519 config->method = GIBIAS_METHOD_MASTER;
3520 config->option = GIBIAS_OPTION_CURVE;
3523 if (!strcmp(s,
"ZMASTER+CURVE")) {
3524 config->method = GIBIAS_METHOD_ZMASTER;
3525 config->option = GIBIAS_OPTION_CURVE;
3528 cx_assert(config->method != GIBIAS_METHOD_UNDEFINED);
3530 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.areas");
3531 s = cpl_parameter_get_string(p);
3532 config->areas = cx_strdup(s);
3534 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.sigma");
3535 config->sigma = cpl_parameter_get_double(p);
3537 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.iterations");
3538 config->iterations = cpl_parameter_get_int(p);
3540 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.fraction");
3541 config->fraction = cpl_parameter_get_double(p);
3543 if (config->method == GIBIAS_METHOD_CURVE ||
3544 config->option == GIBIAS_OPTION_CURVE) {
3545 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.xorder");
3546 config->xdeg = 1 + cpl_parameter_get_int(p);
3548 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.yorder");
3549 config->ydeg = 1 + cpl_parameter_get_int(p);
3552 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.xstep");
3553 config->xstep = cpl_parameter_get_int(p);
3555 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.ystep");
3556 config->ystep = cpl_parameter_get_int(p);
3580 if (config->areas) {
3581 cx_free(config->areas);
3582 config->areas = NULL;
3614 p = cpl_parameter_new_value(
"giraffe.biasremoval.remove",
3616 "Enable bias removal",
3617 "giraffe.biasremoval",
3619 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"remove-bias");
3620 cpl_parameterlist_append(list, p);
3623 p = cpl_parameter_new_enum(
"giraffe.biasremoval.method",
3625 "Bias removal method",
3626 "giraffe.biasremoval",
3627 "PROFILE", 11,
"UNIFORM",
"PLANE",
"CURVE",
3628 "PROFILE",
"MASTER",
"ZMASTER",
"PROFILE+CURVE",
3629 "MASTER+PLANE",
"MASTER+CURVE",
"ZMASTER+PLANE",
3631 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-method");
3632 cpl_parameterlist_append(list, p);
3635 p = cpl_parameter_new_value(
"giraffe.biasremoval.areas",
3637 "Bias areas to use "
3638 "(Xl0:Xr0:Yl0:Yr0, ... ,Xln:Xrn:Yln:Yrn)",
3639 "giraffe.biasremoval",
3641 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-areas");
3642 cpl_parameterlist_append(list, p);
3645 p = cpl_parameter_new_value(
"giraffe.biasremoval.sigma",
3647 "Sigma Clipping: sigma threshold factor",
3648 "giraffe.biasremoval",
3650 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-sigma");
3651 cpl_parameterlist_append(list, p);
3654 p = cpl_parameter_new_value(
"giraffe.biasremoval.iterations",
3656 "Sigma Clipping: maximum number of "
3658 "giraffe.biasremoval",
3660 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-niter");
3661 cpl_parameterlist_append(list, p);
3664 p = cpl_parameter_new_value(
"giraffe.biasremoval.fraction",
3666 "Sigma Clipping: minimum fraction of points "
3667 "accepted/total [0.0..1.0]",
3668 "giraffe.biasremoval",
3670 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-mfrac");
3671 cpl_parameterlist_append(list, p);
3674 p = cpl_parameter_new_value(
"giraffe.biasremoval.xorder",
3676 "Order of X polynomial fit (CURVE method "
3678 "giraffe.biasremoval",
3680 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-xorder");
3681 cpl_parameterlist_append(list, p);
3683 p = cpl_parameter_new_value(
"giraffe.biasremoval.yorder",
3685 "Order of Y polynomial fit (CURVE method "
3687 "giraffe.biasremoval",
3689 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-yorder");
3690 cpl_parameterlist_append(list, p);
3693 p = cpl_parameter_new_value(
"giraffe.biasremoval.xstep",
3695 "Sampling step along X (CURVE method only)",
3696 "giraffe.biasremoval",
3698 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-xstep");
3699 cpl_parameterlist_append(list, p);
3702 p = cpl_parameter_new_value(
"giraffe.biasremoval.ystep",
3704 "Sampling step along Y (CURVE method only)",
3705 "giraffe.biasremoval",
3707 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-ystep");
3708 cpl_parameterlist_append(list, p);