29 #include <cxstrutils.h>
33 #include <cpl_parameter.h>
39 #include "gichebyshev.h"
52 struct GiBiasResults {
66 typedef struct GiBiasResults GiBiasResults;
74 _giraffe_biasresults_clear(GiBiasResults *
self)
84 cx_string_delete(self->limits);
89 cpl_matrix_delete(self->coeffs);
94 cpl_image_delete(self->model);
110 _giraffe_method_string(cx_string *
string, GiBiasMethod method,
115 case GIBIAS_METHOD_UNIFORM:
116 cx_string_set(
string,
"UNIFORM");
119 case GIBIAS_METHOD_PLANE:
120 cx_string_set(
string,
"PLANE");
123 case GIBIAS_METHOD_CURVE:
124 cx_string_set(
string,
"CURVE");
127 case GIBIAS_METHOD_PROFILE:
128 cx_string_set(
string,
"PROFILE");
131 case GIBIAS_METHOD_MASTER:
132 cx_string_set(
string,
"MASTER");
135 case GIBIAS_METHOD_ZMASTER:
136 cx_string_set(
string,
"ZMASTER");
143 if (option != GIBIAS_OPTION_UNDEFINED) {
145 case GIBIAS_OPTION_PLANE:
146 cx_string_append(
string,
"+PLANE");
149 case GIBIAS_OPTION_CURVE:
150 cx_string_append(
string,
"+CURVE");
164 _giraffe_stringify_coefficients(cx_string *
string, cpl_matrix *matrix)
169 cxchar *tmp = cx_line_alloc();
171 cxint nr = cpl_matrix_get_nrow(matrix);
172 cxint nc = cpl_matrix_get_ncol(matrix);
174 cxsize sz = cx_line_max();
176 cxdouble *data = cpl_matrix_get_data(matrix);
179 for (i = 0; i < nr; i++) {
180 for (j = 0; j < nc; j++) {
181 snprintf(tmp, sz,
"%g", data[i * nc + j]);
182 cx_string_append(
string, tmp);
184 if (i != nr - 1 || j < nc - 1) {
185 cx_string_append(
string,
":");
213 _giraffe_compare_overscans(
const GiImage* image1,
const GiImage* image2)
216 cxint32 l1ovscx = -1;
217 cxint32 l1ovscy = -1;
218 cxint32 l1prscx = -1;
219 cxint32 l1prscy = -1;
220 cxint32 l2ovscx = -1;
221 cxint32 l2ovscy = -1;
222 cxint32 l2prscx = -1;
223 cxint32 l2prscy = -1;
225 cpl_propertylist *l1, *l2;
228 cx_assert(image1 != NULL && image2 != NULL);
233 if (cpl_propertylist_has(l1, GIALIAS_OVSCX)) {
234 l1ovscx = cpl_propertylist_get_int(l1, GIALIAS_OVSCX);
236 if (cpl_propertylist_has(l1, GIALIAS_OVSCY)) {
237 l1ovscy = cpl_propertylist_get_int(l1, GIALIAS_OVSCY);
239 if (cpl_propertylist_has(l1, GIALIAS_PRSCX)) {
240 l1prscx = cpl_propertylist_get_int(l1, GIALIAS_PRSCX);
242 if (cpl_propertylist_has(l1, GIALIAS_PRSCY)) {
243 l1prscy = cpl_propertylist_get_int(l1, GIALIAS_PRSCY);
246 if (cpl_propertylist_has(l2, GIALIAS_OVSCX)) {
247 l2ovscx = cpl_propertylist_get_int(l2, GIALIAS_OVSCX);
249 if (cpl_propertylist_has(l2, GIALIAS_OVSCY)) {
250 l2ovscy = cpl_propertylist_get_int(l2, GIALIAS_OVSCY);
252 if (cpl_propertylist_has(l2, GIALIAS_PRSCX)) {
253 l2prscx = cpl_propertylist_get_int(l2, GIALIAS_PRSCX);
255 if (cpl_propertylist_has(l2, GIALIAS_PRSCY)) {
256 l2prscy = cpl_propertylist_get_int(l2, GIALIAS_PRSCY);
259 if (l1ovscx != l2ovscx || l1ovscy != l2ovscy) {
263 if (l1prscx != l2prscx || l1prscy != l2prscy) {
285 inline static cpl_matrix*
286 _giraffe_bias_get_areas(
const cxchar*
string)
289 cxchar** regions = NULL;
291 cpl_matrix* areas = NULL;
294 cx_assert(
string != NULL);
296 regions = cx_strsplit(
string,
",", -1);
298 if (regions == NULL) {
305 const cxsize nvalues = 4;
310 while (regions[i] != NULL) {
315 areas = cpl_matrix_new(nregions, nvalues);
318 while (regions[i] != NULL) {
320 register cxsize j = 0;
322 cxchar** limits = cx_strsplit(regions[i],
":", nvalues);
325 if (limits == NULL) {
327 cpl_matrix_delete(areas);
334 for (j = 0; j < nvalues; ++j) {
342 if (limits[j] == NULL || *limits[j] ==
'\0') {
349 value = strtol(limits[j], &last, 10);
355 if ((errno == ERANGE &&
356 (value == LONG_MAX || value == LONG_MIN)) ||
357 (errno != 0 && value == 0)) {
374 cpl_matrix_set(areas, i, j, value);
383 cpl_matrix_delete(areas);
386 cx_strfreev(regions);
397 cx_strfreev(regions);
454 _giraffe_bias_compute_mean(GiBiasResults* results,
const cpl_image* image,
455 const cpl_matrix* areas, cxdouble kappa, cxint numiter,
456 cxdouble maxfraction)
459 const cxchar*
const fctid =
"giraffe_bias_compute_mean";
462 const cxdouble* pdimg = NULL;
465 cxint img_dimx, img_dimy;
469 cxint x0, x1, x2, x3;
474 cxlong naccepted = 0L;
476 cxlong pixcount = 0L;
478 cxdouble currfraction = 2.;
480 cx_string* tmp = NULL;
482 cpl_matrix* matrix_zz1;
483 cpl_matrix* matrix_zz1diff;
490 if (results->limits == NULL) {
491 cpl_msg_info(fctid,
"Unable to store biaslimits return parameter, "
496 if (cpl_image_get_type(image) != CPL_TYPE_DOUBLE) {
497 cpl_msg_info(fctid,
"Only images allowed of type double, "
512 cpl_msg_info(fctid,
"Bias Areas: Missing bias areas, "
518 ba_num_rows = cpl_matrix_get_nrow(areas);
520 for (j = 0; j < ba_num_rows; j++) {
521 x3 = (cxint) cpl_matrix_get(areas, j, 3);
522 x2 = (cxint) cpl_matrix_get(areas, j, 2);
523 x1 = (cxint) cpl_matrix_get(areas, j, 1);
524 x0 = (cxint) cpl_matrix_get(areas, j, 0);
526 ntotal += (cxulong) ((x3 - x2 + 1) * (x1 - x0 + 1));
530 cpl_msg_info(fctid,
"Bias Areas: Inconsistent specification, "
535 matrix_zz1 = cpl_matrix_new(ntotal, 1);
536 matrix_zz1diff = cpl_matrix_new(ntotal, 1);
543 img_dimx = cpl_image_get_size_x(image);
544 img_dimy = cpl_image_get_size_y(image);
546 cx_string_set(results->limits,
"");
548 for (j = 0; j < ba_num_rows; j++) {
550 x3 = (cxint) cpl_matrix_get(areas, j, 3);
551 x2 = (cxint) cpl_matrix_get(areas, j, 2);
552 x1 = (cxint) cpl_matrix_get(areas, j, 1);
553 x0 = (cxint) cpl_matrix_get(areas, j, 0);
555 if ((x0 > img_dimx) || (x1 > img_dimx) ||
556 (x2 > img_dimy) || (x3 > img_dimy)) {
560 tmp = cx_string_new();
562 cx_string_sprintf(tmp,
"%d:%d:%d:%d;", x0, x1, x2, x3);
563 cx_string_append(results->limits, cx_string_get(tmp));
565 cx_string_delete(tmp);
568 pdimg = cpl_image_get_data_double_const(image);
570 for (l = x2; l < x3 + 1; l++) {
571 for (k = x0; k < x1 + 1; k++) {
572 cpl_matrix_set(matrix_zz1, n, 1, pdimg[k + l * img_dimx]);
579 cpl_msg_info(fctid,
"Bias Areas: Validation failed, aborting ...");
581 cpl_matrix_delete(matrix_zz1);
582 cpl_matrix_delete(matrix_zz1diff);
587 cpl_msg_info(fctid,
"Bias Areas: Using %s",
588 cx_string_get(results->limits));
599 cpl_msg_info(fctid,
"Sigma Clipping : Start");
603 while ((naccepted > 0) && (curriter < numiter) &&
604 (currfraction > maxfraction)) {
606 cxdouble ksigma = 0.;
608 results->mean = cpl_matrix_get_mean(matrix_zz1);
611 for (k = 0; k < cpl_matrix_get_nrow(matrix_zz1); k++) {
612 cpl_matrix_set(matrix_zz1diff, k, 0,
613 cpl_matrix_get(matrix_zz1, k, 0) - results->mean);
616 cpl_msg_info(fctid,
"bias[%d]: mean = %5g, sigma = %6.3g, "
617 "ratio = %6.3g, accepted = %ld\n", curriter,
618 results->mean, sigma, currfraction, naccepted);
620 ksigma = sigma * kappa;
623 for (l = 0; l < cpl_matrix_get_nrow(matrix_zz1); l++) {
624 if (fabs(cpl_matrix_get(matrix_zz1diff, l, 0)) > ksigma) {
628 cpl_matrix_set(matrix_zz1, pixcount, 0,
629 cpl_matrix_get(matrix_zz1, l, 0));
633 cpl_matrix_resize(matrix_zz1, 0, 0, 0, pixcount -
634 cpl_matrix_get_nrow(matrix_zz1));
636 cpl_matrix_resize(matrix_zz1diff, 0, 0, 0, pixcount -
637 cpl_matrix_get_nrow(matrix_zz1diff));
639 if (pixcount == naccepted) {
643 naccepted = pixcount;
645 currfraction = (cxdouble) naccepted / (cxdouble) ntotal;
649 cpl_msg_info(fctid,
"Sigma Clipping : End");
656 results->mean = cpl_matrix_get_mean(matrix_zz1);
659 results->sigma = results->rms / sqrt(cpl_matrix_get_nrow(matrix_zz1));
662 cpl_msg_info(fctid,
"Sigma Clipping Results : bias[%d]: mean = %5g, "
663 "sigma = %6.3g, ratio = %6.3g, accepted=%ld\n", curriter,
664 results->mean, results->rms, currfraction, naccepted);
671 if (matrix_zz1 != NULL) {
672 cpl_matrix_delete(matrix_zz1);
675 if (matrix_zz1diff != NULL) {
676 cpl_matrix_delete(matrix_zz1diff);
737 _giraffe_bias_compute_plane(GiBiasResults* results,
const cpl_image* image,
738 const cpl_matrix* areas, cxdouble kappa,
739 cxint numiter, cxdouble maxfraction)
742 const cxchar*
const fctid =
"giraffe_bias_compute_plane";
753 cxsize naccepted = 0;
755 cxdouble fraction = 1.;
758 cpl_matrix* xbs = NULL;
759 cpl_matrix* ybs = NULL;
760 cpl_matrix* zbs = NULL;
761 cpl_matrix* coeffs = NULL;
764 cx_assert(results->limits != NULL);
765 cx_assert(results->coeffs == NULL);
767 cx_assert(areas != NULL);
769 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
776 nareas = cpl_matrix_get_nrow(areas);
778 for (j = 0; j < nareas; j++) {
780 cxint x3 = (cxint) cpl_matrix_get(areas, j, 3);
781 cxint x2 = (cxint) cpl_matrix_get(areas, j, 2);
782 cxint x1 = (cxint) cpl_matrix_get(areas, j, 1);
783 cxint x0 = (cxint) cpl_matrix_get(areas, j, 0);
785 ntotal += (cxsize) ((x3 - x2 + 1) * (x1 - x0 + 1));
791 cpl_msg_info(fctid,
"Bias Areas: Inconsistent specification, "
797 nx = cpl_image_get_size_x(image);
798 ny = cpl_image_get_size_y(image);
800 cpl_msg_info(fctid,
"Bias Areas: specified are %zu points in %dx%d "
801 "image", ntotal, nx, ny);
812 cx_string_set(results->limits,
"");
814 xbs = cpl_matrix_new(ntotal, 1);
815 ybs = cpl_matrix_new(ntotal, 1);
816 zbs = cpl_matrix_new(1, ntotal);
818 for (j = 0; j < nareas; ++j) {
820 const cxdouble* _img = cpl_image_get_data_double_const(image);
824 cx_string* tmp = NULL;
827 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
828 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
829 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
830 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
832 if ((x0 > nx) || (x1 > nx) || (x2 > ny) || (x3 > ny)) {
836 tmp = cx_string_new();
838 cx_string_sprintf(tmp,
"%d:%d:%d:%d;", x0, x1, x2, x3);
839 cx_string_append(results->limits, cx_string_get(tmp));
841 cx_string_delete(tmp);
844 for (k = x2; k < x3 + 1; ++k) {
846 register cxint l = 0;
849 for (l = x0; l < x1 + 1; ++l) {
851 cpl_matrix_set(xbs, n, 0, l);
852 cpl_matrix_set(ybs, n, 0, k);
853 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
862 cpl_matrix_set_size(xbs, n, 1);
863 cpl_matrix_set_size(ybs, n, 1);
864 cpl_matrix_set_size(zbs, 1, n);
868 cpl_msg_info(fctid,
"Bias Areas: Validation failed, aborting...");
870 cpl_matrix_delete(xbs);
871 cpl_matrix_delete(ybs);
872 cpl_matrix_delete(zbs);
880 cpl_msg_info(fctid,
"Bias Areas: Using %s [%zu pixels]",
881 cx_string_get(results->limits), ntotal);
888 cpl_msg_info(fctid,
"Sigma Clipping : Start");
892 while ((naccepted > 0) && (iteration < numiter) &&
893 (fraction > maxfraction)) {
897 cxdouble ksigma = 0.;
899 cpl_matrix* base = NULL;
900 cpl_matrix* fit = NULL;
903 base = cpl_matrix_new(3, naccepted);
907 cpl_msg_info(fctid,
"Sigma Clipping: Error creating design "
910 cpl_matrix_delete(zbs);
911 cpl_matrix_delete(ybs);
912 cpl_matrix_delete(xbs);
917 for (k = 0; k < naccepted; ++k) {
919 cpl_matrix_set(base, 0, k, 1.);
920 cpl_matrix_set(base, 1, k, cpl_matrix_get(xbs, k, 0));
921 cpl_matrix_set(base, 2, k, cpl_matrix_get(ybs, k, 0));
925 cpl_matrix_delete(coeffs);
930 if (coeffs == NULL) {
932 cpl_msg_info(fctid,
"Sigma Clipping : Error in least square "
933 "solution, aborting...");
935 cpl_matrix_delete(base);
938 cpl_matrix_delete(xbs);
939 cpl_matrix_delete(ybs);
940 cpl_matrix_delete(zbs);
951 fit = cpl_matrix_product_create(coeffs, base);
953 cpl_matrix_delete(base);
956 results->mean = cpl_matrix_get_mean(fit);
960 cpl_msg_info(fctid,
"Sigma Clipping : bias plane[%d]: %g + "
961 "%g * x + %g * y, sigma = %.5g, ratio = %.4g, "
962 "accepted = %zu\n", iteration,
963 cpl_matrix_get(coeffs, 0, 0),
964 cpl_matrix_get(coeffs, 0, 1),
965 cpl_matrix_get(coeffs, 0, 2),
966 sigma, fraction, naccepted);
973 ksigma = sigma * kappa;
977 for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
979 register cxdouble z = cpl_matrix_get(zbs, 0, j);
981 if (fabs(cpl_matrix_get(fit, 0, j) - z) > ksigma) {
985 cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
987 cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
989 cpl_matrix_set(zbs, 0, n, z);
994 cpl_matrix_set_size(xbs, n, 1);
995 cpl_matrix_set_size(ybs, n, 1);
996 cpl_matrix_set_size(zbs, 1, n);
998 cpl_matrix_delete(fit);
1001 if (n == naccepted) {
1007 fraction = (cxdouble)naccepted / (cxdouble)ntotal;
1012 cpl_msg_info(fctid,
"Sigma Clipping : End");
1019 results->coeffs = coeffs;
1020 results->rms = sigma;
1028 results->sigma = sigma / sqrt(cpl_matrix_get_ncol(zbs));
1030 cpl_msg_info(fctid,
"Sigma Clipping Results (%d/%zu, sigma = %g)",
1031 iteration, naccepted, results->rms);
1038 cpl_matrix_delete(xbs);
1041 cpl_matrix_delete(ybs);
1044 cpl_matrix_delete(zbs);
1047 return EXIT_SUCCESS;
1094 _giraffe_bias_compute_curve(GiBiasResults* results,
const cpl_image* image,
1095 const cpl_matrix *areas, cxdouble kappa,
1096 cxint numiter, cxdouble maxfraction,
1097 cxdouble xdeg, cxdouble ydeg,
1098 cxdouble xstep, cxdouble ystep)
1101 const cxchar*
const fctid =
"giraffe_bias_compute_curve";
1107 cxint iteration = 0;
1111 cxsize naccepted = 0;
1113 cxdouble fraction = 1.;
1114 cxdouble sigma = 0.;
1116 cpl_matrix* xbs = NULL;
1117 cpl_matrix* ybs = NULL;
1118 cpl_matrix* zbs = NULL;
1120 GiChebyshev2D* fit = NULL;
1123 cx_assert(results != NULL);
1124 cx_assert(results->limits != NULL);
1125 cx_assert(results->coeffs == NULL);
1127 cx_assert(areas != NULL);
1129 cx_assert(image != NULL);
1130 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
1137 nareas = cpl_matrix_get_nrow(areas);
1139 for (j = 0; j < nareas; ++j) {
1141 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
1142 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
1143 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
1144 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
1146 ntotal += (cxulong) (ceil((1. + x3 - x2) / ystep) *
1147 ceil((1. + x1 - x0) / xstep));
1150 nx = cpl_image_get_size_x(image);
1151 ny = cpl_image_get_size_y(image);
1153 cpl_msg_info(fctid,
"Bias Areas: Found %zu points in %dx%d image",
1162 results->sigma = 0.;
1165 cx_string_set(results->limits,
"");
1167 xbs = cpl_matrix_new(ntotal, 1);
1168 ybs = cpl_matrix_new(ntotal, 1);
1169 zbs = cpl_matrix_new(1, ntotal);
1171 for (j = 0; j < nareas; ++j) {
1173 const cxdouble* _img = cpl_image_get_data_double_const(image);
1177 cx_string* tmp = NULL;
1180 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
1181 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
1182 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
1183 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
1185 if ((x0 > nx) || (x1 > nx) ||
1186 (x2 > ny) || (x3 > ny)) {
1190 tmp = cx_string_new();
1192 cx_string_sprintf(tmp,
"%d:%d:%d:%d;", x0, x1, x2, x3);
1193 cx_string_append(results->limits, cx_string_get(tmp));
1195 cx_string_delete(tmp);
1198 for (k = x2; k < x3 + 1; k += ystep) {
1200 register cxint l = 0;
1203 for (l = x0; l < x1 + 1; l += xstep) {
1205 cpl_matrix_set(xbs, n, 0, l);
1206 cpl_matrix_set(ybs, n, 0, k);
1207 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
1216 cpl_matrix_set_size(xbs, n, 1);
1217 cpl_matrix_set_size(ybs, n, 1);
1218 cpl_matrix_set_size(zbs, 1, n);
1222 cpl_msg_info(fctid,
"Bias Areas: Validation failed, aborting...");
1224 cpl_matrix_delete(xbs);
1225 cpl_matrix_delete(ybs);
1226 cpl_matrix_delete(zbs);
1234 cpl_msg_info(fctid,
"Bias Areas: Using %s [%zu pixels]",
1235 cx_string_get(results->limits), ntotal);
1242 cpl_msg_info(fctid,
"Sigma Clipping : Start");
1246 while ((naccepted > 0) && (iteration < numiter) &&
1247 (fraction > maxfraction)) {
1251 cxdouble ksigma = 0.;
1253 cpl_matrix* base = NULL;
1254 cpl_matrix* coeffs = NULL;
1255 cpl_matrix* _coeffs = NULL;
1256 cpl_matrix* _fit = NULL;
1259 base = giraffe_chebyshev_base2d(0., 0., nx, ny,
1260 xdeg, ydeg, xbs, ybs);
1264 cpl_msg_info(fctid,
"Sigma Clipping: Error creating design "
1267 cpl_matrix_delete(zbs);
1268 cpl_matrix_delete(ybs);
1269 cpl_matrix_delete(xbs);
1276 cpl_matrix_delete(base);
1279 if (_coeffs == NULL) {
1281 cpl_msg_info(fctid,
"Sigma Clipping : Error in least square "
1282 "solution, aborting...");
1284 cpl_matrix_delete(xbs);
1285 cpl_matrix_delete(ybs);
1286 cpl_matrix_delete(zbs);
1297 coeffs = cpl_matrix_wrap(xdeg, ydeg, cpl_matrix_get_data(_coeffs));
1300 giraffe_chebyshev2d_delete(fit);
1304 fit = giraffe_chebyshev2d_new(xdeg - 1, ydeg - 1);
1305 status = giraffe_chebyshev2d_set(fit, 0., nx, 0., ny,
1310 giraffe_chebyshev2d_delete(fit);
1313 cpl_matrix_unwrap(coeffs);
1316 cpl_matrix_delete(_coeffs);
1319 cpl_matrix_delete(xbs);
1320 cpl_matrix_delete(ybs);
1321 cpl_matrix_delete(zbs);
1327 cpl_matrix_unwrap(coeffs);
1330 cpl_matrix_delete(_coeffs);
1333 _fit = cpl_matrix_new(1, cpl_matrix_get_ncol(zbs));
1335 for (j = 0; j < cpl_matrix_get_ncol(_fit); ++j) {
1337 cxdouble x = cpl_matrix_get(xbs, n, 0);
1338 cxdouble y = cpl_matrix_get(ybs, n, 0);
1339 cxdouble z = giraffe_chebyshev2d_eval(fit, x, y);
1341 cpl_matrix_set(_fit, 0, j, z);
1345 results->mean = cpl_matrix_get_mean(_fit);
1349 cpl_msg_info(fctid,
"Sigma Clipping : bias surface[%d]: "
1350 "sigma = %8.5g, ratio = %7.4g, accepted = %zu\n",
1351 iteration, sigma, fraction, naccepted);
1358 ksigma = sigma * kappa;
1362 for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
1364 register cxdouble z = cpl_matrix_get(zbs, 0, j);
1366 if (fabs(cpl_matrix_get(_fit, 0, j) - z) > ksigma) {
1370 cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
1371 cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
1372 cpl_matrix_set(zbs, 0, n, z);
1377 cpl_matrix_set_size(xbs, n, 1);
1378 cpl_matrix_set_size(ybs, n, 1);
1379 cpl_matrix_set_size(zbs, 1, n);
1381 cpl_matrix_delete(_fit);
1385 if (n == naccepted) {
1391 fraction = (cxdouble)naccepted / (cxdouble)ntotal;
1396 cpl_msg_info(fctid,
"Sigma Clipping : End");
1403 results->coeffs = cpl_matrix_duplicate(giraffe_chebyshev2d_coeffs(fit));
1404 results->rms = sigma;
1412 results->sigma = sigma / sqrt(cpl_matrix_get_ncol(zbs));
1414 cpl_msg_info(fctid,
"Sigma Clipping Results (%d/%zu, sigma = %g)",
1415 iteration, naccepted, results->rms);
1422 giraffe_chebyshev2d_delete(fit);
1425 cpl_matrix_delete(xbs);
1428 cpl_matrix_delete(ybs);
1431 cpl_matrix_delete(zbs);
1434 return EXIT_SUCCESS;
1463 _giraffe_bias_compute_profile(GiBiasResults* results,
const cpl_image* image,
1464 const cpl_matrix* areas, cxchar axis)
1467 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);
1674 rms += residual / (sz - 1);
1675 sigma += residual / (sz * (sz - 1));
1681 cxint sz = y_1 - y_0 + 1;
1683 cxdouble residual = 0.;
1685 const cxdouble* _img = cpl_image_get_data_double_const(image);
1687 cxdouble* _profile = cpl_image_get_data_double(profile);
1688 cxdouble* buffer = cx_calloc(sz,
sizeof(cxdouble));
1691 for (i = x_0; i < x_1 + 1; ++i) {
1693 register cxint l = 0;
1695 register cxdouble m = 0.;
1698 for (l = y_0; l < y_1 + 1; ++l) {
1699 buffer[l - y_0] = _img[l * nx + i];
1702 m = giraffe_array_mean(buffer, sz);
1703 _profile[i * sx + j] = m;
1705 for (l = 0; l < sz; ++l) {
1706 residual += (buffer[l] - m) * (buffer[l] - m);
1717 rms += residual / (sz - 1);
1718 sigma += residual / (sz * (sz - 1));
1729 cpl_matrix_delete(_areas);
1737 model = cpl_image_collapse_create(profile, 1);
1739 cpl_image_delete(profile);
1742 cpl_image_divide_scalar(model, nvalid);
1743 results->model = model;
1747 results->mean = cpl_image_get_mean(results->model);
1750 results->rms = sqrt(rms / (sy * nvalid));
1751 results->sigma = sqrt(sigma / sy) / nvalid;
1754 return EXIT_SUCCESS;
1784 _giraffe_bias_fit_profile(GiBiasResults* results,
const cpl_image* profile,
1790 cxint sy = cpl_image_get_size_y(profile);
1791 cxint nc = order + 1;
1793 cxdouble* _profile = (cxdouble*)cpl_image_get_data_double_const(profile);
1795 cpl_matrix* base = NULL;
1796 cpl_matrix* baseT = NULL;
1797 cpl_matrix* coeff = NULL;
1798 cpl_matrix* fit = NULL;
1799 cpl_matrix* x = cpl_matrix_new(sy, 1);
1800 cpl_matrix* y = cpl_matrix_wrap(sy, 1, _profile);
1801 cpl_matrix* covy = cpl_matrix_new(sy, sy);
1802 cpl_matrix* covc = cpl_matrix_new(nc, nc);
1805 if ((x == NULL) || (y == NULL) || (covy == NULL) || (covc == NULL)) {
1807 cpl_matrix_delete(covc);
1808 cpl_matrix_delete(covy);
1809 cpl_matrix_unwrap(y);
1810 cpl_matrix_delete(x);
1821 for (i = 0; i < sy; ++i)
1823 cpl_matrix_set(x, i, 0, i);
1826 baseT = giraffe_chebyshev_base1d(0., sy, nc, x);
1827 base = cpl_matrix_transpose_create(baseT);
1829 cpl_matrix_delete(baseT);
1832 cpl_matrix_delete(x);
1836 cpl_matrix_unwrap(y);
1840 cpl_matrix_fill_diagonal(covy, results->rms * results->rms, 0);
1850 cpl_matrix_delete(covy);
1853 if (coeff == NULL) {
1855 cpl_matrix_delete(base);
1858 cpl_matrix_unwrap(y);
1866 fit = cpl_matrix_product_create(base, coeff);
1868 cpl_matrix_delete(coeff);
1871 cpl_matrix_delete(base);
1876 cpl_matrix_unwrap(y);
1888 for (i = 0; i < sy; ++i)
1890 cpl_matrix_set(y, i, 0, cpl_matrix_get(fit, i, 0));
1893 cpl_matrix_delete(fit);
1896 cpl_matrix_unwrap(y);
1906 results->mean = cpl_image_get_mean(results->model);
1907 results->sigma = sqrt(cpl_matrix_get(covc, 0, 0));
1910 cpl_image_save(results->model,
"idiot.fits", CPL_BPP_IEEE_FLOAT,
1918 cpl_matrix_delete(covc);
1921 return EXIT_SUCCESS;
1955 _giraffe_bias_compute(GiBiasResults* results, cpl_image* img,
1956 const cpl_matrix* biasareas,
1957 const cpl_image* master_bias,
1958 const cpl_image* bad_pixel_map,
1959 GiBiasMethod method, GiBiasOption option,
1960 GiBiasModel model, cxbool remove_bias, cxdouble mbias,
1961 cxdouble xdeg, cxdouble ydeg, cxdouble xstep,
1962 cxdouble ystep, cxdouble sigma, cxint numiter,
1963 cxdouble maxfraction)
1966 const cxchar*
const fctid =
"giraffe_bias_compute";
1974 cxint master_bias_dimx;
1975 cxint master_bias_dimy;
1976 cxint bad_pixel_dimx;
1977 cxint bad_pixel_dimy;
1979 cxint error_code = 0;
1981 cx_string* s = NULL;
1983 const cpl_matrix* coeff = NULL;
1986 cx_assert(results != NULL);
1987 cx_assert(results->limits == NULL);
1988 cx_assert(results->coeffs == NULL);
1989 cx_assert(results->model == NULL);
1991 cx_assert(img != NULL);
1992 cx_assert(cpl_image_get_type(img) == CPL_TYPE_DOUBLE);
1994 cx_assert(biasareas != NULL);
2001 img_dimx = cpl_image_get_size_x(img);
2002 img_dimy = cpl_image_get_size_y(img);
2004 img_centerx = img_dimx / 2;
2005 img_centery = img_dimy / 2;
2007 results->limits = cx_string_new();
2015 cpl_msg_info(fctid,
"Bias correction will be done.");
2018 cpl_msg_warning(fctid,
"No Bias correction will be done!");
2022 if (model == GIBIAS_MODEL_MEAN) {
2028 cpl_msg_info(fctid,
"Using bias model 'MEAN' ...");
2030 if ((option == GIBIAS_OPTION_PLANE) ||
2031 (method == GIBIAS_METHOD_PLANE)) {
2033 cpl_msg_error(fctid,
"Can not use MEAN model with PLANE method.");
2038 error_code = _giraffe_bias_compute_mean(results, img, biasareas,
2039 sigma, numiter, maxfraction);
2043 else if ((method == GIBIAS_METHOD_CURVE) ||
2044 ((method != GIBIAS_METHOD_PROFILE) &&
2045 (option == GIBIAS_OPTION_CURVE))) {
2052 cpl_msg_info(fctid,
"Using bias model 'CURVE' ...");
2054 error_code = _giraffe_bias_compute_curve(results, img, biasareas,
2055 sigma, numiter, maxfraction, xdeg, ydeg, xstep, ystep);
2057 if (error_code != EXIT_SUCCESS) {
2059 cpl_msg_error(fctid,
"Error during calculation of bias curve, "
2065 coeff = results->coeffs;
2074 cpl_msg_info(fctid,
"Using bias model 'PLANE (FITTED)' ...");
2076 error_code = _giraffe_bias_compute_plane(results, img, biasareas,
2077 sigma, numiter, maxfraction);
2079 if (error_code == EXIT_SUCCESS) {
2081 coeff = results->coeffs;
2083 results->mean = cpl_matrix_get(coeff, 0, 0) +
2084 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2085 cpl_matrix_get(coeff, 0, 2) * img_centery;
2090 cpl_msg_error(fctid,
"Error during calculation of bias plane, "
2103 s = cx_string_new();
2104 _giraffe_method_string(s, method, option);
2106 cpl_msg_info(fctid,
"Using bias method '%s'", cx_string_get(s));
2108 cx_string_delete(s);
2113 case GIBIAS_METHOD_UNIFORM:
2120 if (model == GIBIAS_MODEL_MEAN) {
2122 cpl_msg_info(fctid,
"bias value (areas mean) = %f, "
2123 "bias sigma = %f", results->mean, results->sigma);
2128 cpl_msg_info(fctid,
"bias value at (%d, %d) = %f, "
2129 "bias sigma = %f", img_centerx, img_centery,
2130 results->mean, results->sigma);
2134 if (remove_bias == TRUE) {
2136 cpl_image_subtract_scalar(img, results->mean);
2143 case GIBIAS_METHOD_PLANE:
2146 if (coeff == NULL) {
2148 error_code = _giraffe_bias_compute_plane(results, img,
2149 biasareas, sigma, numiter, maxfraction);
2151 if (error_code == EXIT_SUCCESS) {
2153 coeff = results->coeffs;
2155 results->mean = cpl_matrix_get(coeff, 0, 0) +
2156 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2157 cpl_matrix_get(coeff, 0, 2) * img_centery;
2162 cpl_msg_error(fctid,
"Error during calculation of bias "
2163 "plane, aborting...");
2170 cpl_msg_info(fctid,
"bias value at (%d, %d) = %f, bias sigma = %f, "
2171 "bias plane = %g + %g * x + %g * y", img_centerx,
2172 img_centery, results->mean, results->sigma,
2173 cpl_matrix_get(coeff, 0, 0),
2174 cpl_matrix_get(coeff, 0, 1),
2175 cpl_matrix_get(coeff, 0, 2));
2177 if (remove_bias == TRUE) {
2181 cxdouble* _img = cpl_image_get_data_double(img);
2184 cpl_image* bsmodel = cpl_image_new(img_dimx, img_dimy,
2186 cxdouble* _bsmodel = cpl_image_get_data_double(bsmodel);
2190 for (j = 0; j < img_dimy; j++) {
2192 cxsize jj = j * img_dimx;
2194 for (i = 0; i < img_dimx; i++) {
2197 _bsmodel[jj + i] = cpl_matrix_get(coeff, 0, 0)
2198 + cpl_matrix_get(coeff, 0, 1) * i
2199 + cpl_matrix_get(coeff, 0, 2) * j;
2202 _img[jj + i] -= cpl_matrix_get(coeff, 0, 0)
2203 + cpl_matrix_get(coeff, 0, 1) * i
2204 + cpl_matrix_get(coeff, 0, 2) * j;
2211 cpl_image_save(bsmodel,
"idiot.fits", CPL_BPP_IEEE_FLOAT,
2213 cpl_image_delete(bsmodel);
2222 case GIBIAS_METHOD_CURVE:
2225 if (coeff == NULL) {
2228 _giraffe_bias_compute_curve(results, img, biasareas,
2229 sigma, numiter, maxfraction,
2230 xdeg, ydeg, xstep, ystep);
2232 if (error_code != EXIT_SUCCESS) {
2234 cpl_msg_error(fctid,
"Error during calculation of bias "
2235 "surface curve, aborting...");
2240 coeff = results->coeffs;
2244 cpl_msg_info(fctid,
"bias value %f, bias sigma = %f\n",
2245 results->mean, results->sigma);
2248 if (remove_bias == TRUE) {
2250 cpl_image* bsmodel = NULL;
2252 cpl_matrix* x = cpl_matrix_new(img_dimx * img_dimy, 1);
2253 cpl_matrix* y = cpl_matrix_new(img_dimx * img_dimy, 1);
2254 cpl_matrix* fit = NULL;
2257 for (j = 0; j < img_dimy; ++j) {
2259 register cxsize jj = j * img_dimx;
2261 for (i = 0; i < img_dimx; ++i) {
2263 cpl_matrix_set(x, jj + i, 0, i);
2264 cpl_matrix_set(y, jj + i, 0, j);
2270 fit = giraffe_chebyshev_fit2d(0., 0., img_dimx, img_dimy,
2273 cpl_matrix_delete(y);
2276 cpl_matrix_delete(x);
2279 bsmodel = cpl_image_wrap_double(img_dimx, img_dimy,
2280 cpl_matrix_get_data(fit));
2283 cpl_image_save(bsmodel,
"idiot.fits", CPL_BPP_IEEE_FLOAT,
2287 cpl_image_subtract(img, bsmodel);
2289 cpl_image_unwrap(bsmodel);
2292 cpl_matrix_delete(fit);
2301 case GIBIAS_METHOD_PROFILE:
2304 error_code = _giraffe_bias_compute_profile(results, img,
2307 if (error_code != EXIT_SUCCESS) {
2309 cpl_msg_error(fctid,
"Error computing the bias area(s) "
2315 if (option == GIBIAS_OPTION_CURVE) {
2317 error_code = _giraffe_bias_fit_profile(results, results->model,
2319 if (error_code != EXIT_SUCCESS) {
2321 cpl_msg_error(fctid,
"Error fitting the bias profile");
2328 if (remove_bias == TRUE) {
2330 const cxdouble* _bias =
2331 cpl_image_get_data_double(results->model);
2333 cxdouble* _img = cpl_image_get_data_double(img);
2336 cx_assert(_bias != NULL);
2337 cx_assert(_img != NULL);
2339 for (j = 0; j < img_dimy; ++j) {
2341 cxsize jj = j * img_dimx;
2344 for (i = 0; i < img_dimx; ++i) {
2345 _img[jj + i] -= _bias[j];
2356 case GIBIAS_METHOD_MASTER:
2357 case GIBIAS_METHOD_ZMASTER:
2360 cxdouble biasdrift = 0.;
2362 if (master_bias == NULL) {
2364 cpl_msg_error(fctid,
"Master Bias Frame required with MASTER "
2370 master_bias_dimx = cpl_image_get_size_x(master_bias);
2371 master_bias_dimy = cpl_image_get_size_y(master_bias);
2373 if ((master_bias_dimx != img_dimx) ||
2374 (master_bias_dimy != img_dimy)) {
2376 cpl_msg_error(fctid,
"Invalid master bias! Size should "
2377 "be [%d, %d] but is [%d, %d].", img_dimx, img_dimy,
2378 master_bias_dimx, master_bias_dimy);
2383 if (coeff == NULL) {
2385 if (option == GIBIAS_OPTION_CURVE) {
2387 error_code = _giraffe_bias_compute_curve(results, img,
2388 biasareas, sigma, numiter, maxfraction, xdeg,
2389 ydeg, xstep, ystep);
2391 if (error_code != EXIT_SUCCESS) {
2393 cpl_msg_error(fctid,
"Error during calculation of "
2394 "bias surface curve, aborting...");
2399 coeff = results->coeffs;
2408 error_code = _giraffe_bias_compute_plane(results, img,
2409 biasareas, sigma, numiter, maxfraction);
2411 if (error_code == EXIT_SUCCESS) {
2413 coeff = results->coeffs;
2415 results->mean = cpl_matrix_get(coeff, 0, 0) +
2416 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2417 cpl_matrix_get(coeff, 0, 2) * img_centery;
2422 cpl_msg_error(fctid,
"Error during calculation of "
2423 "bias surface plane, aborting...");
2432 if (method == GIBIAS_METHOD_ZMASTER) {
2440 biasdrift = results->mean - mbias;
2442 cpl_msg_info(fctid,
"Using bias drift value %.4g", biasdrift);
2451 if (option == GIBIAS_OPTION_CURVE) {
2455 cpl_msg_info(fctid,
"Computed bias mean = %.4f, bias "
2456 "sigma = %.4e", results->mean, results->sigma);
2461 if (option == GIBIAS_OPTION_PLANE) {
2463 cpl_msg_info(fctid,
"Bias plane = %.4e + "
2464 "%.4e * x + %.4e * y",
2465 cpl_matrix_get(coeff, 0, 0),
2466 cpl_matrix_get(coeff, 0, 1),
2467 cpl_matrix_get(coeff, 0, 2));
2468 cpl_msg_info(fctid,
"Computed bias mean = %.4f, "
2469 "bias sigma = %.4e at (%d, %d)",
2470 results->mean, results->sigma,
2471 img_centerx, img_centery);
2476 cpl_msg_info(fctid,
"Computed bias mean = %.4f, bias "
2478 results->mean, results->sigma);
2489 if (remove_bias == TRUE) {
2506 if (bad_pixel_map == NULL) {
2508 cpl_image_subtract(img, master_bias);
2510 if (biasdrift != 0.) {
2511 cpl_image_subtract_scalar(img, biasdrift);
2517 bad_pixel_dimx = cpl_image_get_size_x(bad_pixel_map);
2518 bad_pixel_dimy = cpl_image_get_size_y(bad_pixel_map);
2520 if ((bad_pixel_dimx != img_dimx) ||
2521 (bad_pixel_dimy != img_dimy)) {
2523 cpl_msg_error(fctid,
"Invalid bad pixel map size "
2524 "should be [%d, %d] but is [%d, %d].",
2526 bad_pixel_dimx, bad_pixel_dimy);
2531 if (option == GIBIAS_OPTION_CURVE) {
2533 const cxint* _bpix =
2534 cpl_image_get_data_int_const(bad_pixel_map);
2536 const cxdouble* _mbias =
2537 cpl_image_get_data_double_const(master_bias);
2539 cxdouble* _img = cpl_image_get_data_double(img);
2542 cpl_matrix_new(img_dimx * img_dimy, 1);
2544 cpl_matrix_new(img_dimx * img_dimy, 1);
2545 cpl_matrix* fit = NULL;
2548 for (j = 0; j < img_dimy; ++j) {
2550 register cxsize jj = j * img_dimx;
2552 for (i = 0; i < img_dimx; ++i) {
2554 cpl_matrix_set(x, jj + i, 0, i);
2555 cpl_matrix_set(y, jj + i, 0, j);
2560 fit = giraffe_chebyshev_fit2d(0., 0.,
2564 cpl_matrix_delete(y);
2567 cpl_matrix_delete(x);
2571 for (j = 0; j < img_dimy; ++j) {
2573 register cxsize jj = j * img_dimx;
2575 for (i = 0; i < img_dimx; ++i) {
2577 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2579 _img[jj + i] -= (_mbias[jj + i] +
2584 _img[jj + i] -= cpl_matrix_get(fit, i, j);
2591 cpl_matrix_delete(fit);
2595 else if (option == GIBIAS_OPTION_PLANE) {
2597 const cxint* _bpix =
2598 cpl_image_get_data_int_const(bad_pixel_map);
2600 const cxdouble* _mbias =
2601 cpl_image_get_data_double_const(master_bias);
2603 cxdouble* _img = cpl_image_get_data_double(img);
2606 for (j = 0; j < img_dimy; j++) {
2608 cxsize jj = j * img_dimx;
2610 for (i = 0; i < img_dimx; i++) {
2612 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2614 _img[jj + i] -= (_mbias[jj + i] +
2621 cpl_matrix_get(coeff, 0, 0) +
2622 cpl_matrix_get(coeff, 0, 1) * j +
2623 cpl_matrix_get(coeff, 0, 2) * i;
2633 const cxint* _bpix =
2634 cpl_image_get_data_int_const(bad_pixel_map);
2636 const cxdouble* _mbias =
2637 cpl_image_get_data_double_const(master_bias);
2639 cxdouble* _img = cpl_image_get_data_double(img);
2642 for (j = 0; j < img_dimy; j++) {
2644 cxsize jj = j * img_dimx;
2646 for (i = 0; i < img_dimx; i++) {
2648 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2650 _img[jj + i] -= (_mbias[jj + i] +
2656 _img[jj + i] -= results->mean;
2677 cpl_msg_error(fctid,
"Invalid method, aborting...");
2686 cpl_msg_info(fctid,
"Resulting biaslimits : %s",
2687 cx_string_get(results->limits));
2715 const cxchar *
const _id =
"giraffe_get_raw_areas";
2725 cxint32 tprescx = 0;
2726 cxint32 tprescy = 0;
2727 cxint32 tovrscx = 0;
2728 cxint32 tovrscy = 0;
2733 cpl_propertylist *properties;
2738 cpl_error_set(_id, CPL_ERROR_NULL_INPUT);
2742 winx = cpl_propertylist_get_int(properties, GIALIAS_WINX);
2743 winy = cpl_propertylist_get_int(properties, GIALIAS_WINY);
2745 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
2746 tprescx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
2751 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
2752 tprescy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
2757 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
2758 tovrscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
2763 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
2764 tovrscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
2771 m_tmp = cpl_matrix_new(1, 4);
2775 cpl_matrix_set(m_tmp, row, 0, 0.0);
2776 cpl_matrix_set(m_tmp, row, 1, (cxdouble) prescx - 1);
2777 cpl_matrix_set(m_tmp, row, 2, 0.0);
2778 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2780 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2786 cpl_matrix_set(m_tmp, row, 0, (cxdouble) winx - ovrscx);
2787 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2788 cpl_matrix_set(m_tmp, row, 2, 0.0);
2789 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2791 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2795 if (!(prescy == 0 || prescx == 0 || ovrscx == 0)) {
2797 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2798 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2799 cpl_matrix_set(m_tmp, row, 2, 0.0);
2800 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2802 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2806 else if (!(prescy == 0 || prescx == 0)) {
2808 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2809 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2810 cpl_matrix_set(m_tmp, row, 2, 0.0);
2811 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2813 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2817 else if (!(prescy == 0 || ovrscx == 0)) {
2819 cpl_matrix_set(m_tmp, row, 0, 0.0);
2820 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2821 cpl_matrix_set(m_tmp, row, 2, 0.0);
2822 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2824 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2828 else if (prescy > 0) {
2830 cpl_matrix_set(m_tmp, row, 0, 0.0);
2831 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2832 cpl_matrix_set(m_tmp, row, 2, 0.0);
2833 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2835 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2839 if (!(ovrscy == 0 || prescx == 0 || ovrscx == 0)) {
2841 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2842 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2843 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2844 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2846 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2850 else if (!(ovrscy == 0 || prescx == 0)) {
2852 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2853 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2854 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2855 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2857 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2861 else if (!(ovrscy == 0 || ovrscx == 0)) {
2863 cpl_matrix_set(m_tmp, row, 0, 0.0);
2864 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2865 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2866 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2868 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2872 else if (ovrscy > 0) {
2874 cpl_matrix_set(m_tmp, row, 0, 0.0);
2875 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2876 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2877 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2879 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2883 cpl_matrix_resize(m_tmp, 0, -1, 0, 0);
2887 cpl_matrix_delete(m_tmp);
2917 const cxchar *fctid =
"giraffe_trim_raw_areas";
2920 cxint newlx, newly, newux, newuy;
2934 cpl_msg_error(fctid,
"Image does not contain any properties!");
2938 nx = cpl_image_get_size_x(_image);
2939 ny = cpl_image_get_size_y(_image);
2941 if (cpl_propertylist_has(properties, GIALIAS_NAXIS1)) {
2942 cxint _nx = cpl_propertylist_get_int(properties, GIALIAS_NAXIS1);
2945 cpl_msg_warning(fctid,
"Image size (%d) and image property `%s' "
2946 "(%d) are not consistent! Using image size ...",
2947 nx, GIALIAS_NAXIS1, _nx);
2951 cpl_msg_warning(fctid,
"Image does not contain any property `%s'. "
2952 "Using image size (%d)", GIALIAS_NAXIS1, nx);
2956 if (cpl_propertylist_has(properties, GIALIAS_NAXIS2)) {
2957 cxint _ny = cpl_propertylist_get_int(properties, GIALIAS_NAXIS2);
2960 cpl_msg_warning(fctid,
"Image size (%d) and image property `%s' "
2961 "(%d) are not consistent! Using image size ...",
2962 ny, GIALIAS_NAXIS2, _ny);
2966 cpl_msg_warning(fctid,
"Image does not contain any property `%s'. "
2967 "Using image size (%d)", GIALIAS_NAXIS2, ny);
2970 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
2971 ovscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
2974 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
2975 ovscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
2978 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
2979 prscx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
2982 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
2983 prscy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
2995 tmp = cpl_image_extract(_image, newlx, newly, newux, newuy);
2998 cpl_image_delete(tmp);
3002 nx = cpl_image_get_size_x(_image);
3003 ny = cpl_image_get_size_y(_image);
3010 cpl_propertylist_set_int(properties, GIALIAS_NAXIS1, nx);
3011 cpl_propertylist_set_int(properties, GIALIAS_NAXIS2, ny);
3013 if (cpl_propertylist_has(properties, GIALIAS_CRPIX1)) {
3014 cxdouble crpix1 = cpl_propertylist_get_double(properties,
3022 crpix1 += (cxdouble)prscx;
3023 cpl_propertylist_set_double(properties, GIALIAS_CRPIX1, crpix1);
3026 if (cpl_propertylist_has(properties, GIALIAS_CRPIX2)) {
3027 cxdouble crpix2 = cpl_propertylist_get_double(properties,
3030 crpix2 -= (cxdouble) prscy;
3031 cpl_propertylist_set_double(properties, GIALIAS_CRPIX2, crpix2);
3034 cpl_propertylist_erase(properties, GIALIAS_OVSCX);
3035 cpl_propertylist_erase(properties, GIALIAS_OVSCY);
3036 cpl_propertylist_erase(properties, GIALIAS_PRSCX);
3037 cpl_propertylist_erase(properties, GIALIAS_PRSCY);
3107 const GiImage* master_bias,
const GiImage* bad_pixels,
3108 const cpl_matrix* biaslimits,
const GiBiasConfig* config)
3111 const cxchar*
const fctid =
"giraffe_bias_remove";
3115 cxdouble bias_value = 0.;
3119 cpl_matrix* areas = NULL;
3126 cpl_propertylist* properties;
3128 GiBiasResults bias_results = {0., 0., 0., NULL, NULL, NULL};
3131 cx_assert(raw != NULL);
3132 cx_assert(config != NULL);
3133 cx_assert(biaslimits == NULL);
3135 if (result == NULL) {
3144 if (strncmp(config->areas,
"None", 4) == 0) {
3146 cpl_msg_info(fctid,
"No bias areas specified, using pre/overscan"
3147 "regions of the raw frame.");
3151 if (areas == NULL) {
3152 if (cpl_error_get_code() == CPL_ERROR_NULL_INPUT) {
3153 cpl_msg_error(fctid,
"Invalid raw image! Image does not "
3154 "contain any properties!");
3157 cpl_msg_error(fctid,
"Invalid raw image! Image does not "
3158 "contain or has invalid pre- and overscan "
3168 areas = _giraffe_bias_get_areas(config->areas);
3170 if (areas == NULL) {
3172 cpl_msg_error(fctid,
"Invalid bias area specification '%s'!",
3178 cpl_msg_info(fctid,
"Using bias area(s) '%s' for bias computation",
3194 if (master_bias != NULL) {
3195 cpl_propertylist *_properties;
3196 cxbool is_same_size = FALSE;
3198 is_same_size = _giraffe_compare_overscans(raw, master_bias);
3205 if (is_same_size==FALSE) {
3206 cpl_msg_info(fctid,
"PRE/OVRSCAN Regions do not match between "
3207 "RAW Image and Masterbias Image.");
3214 if (cpl_propertylist_has(_properties, GIALIAS_BIASVALUE)) {
3215 bias_value = cpl_propertylist_get_double(_properties,
3226 if (bad_pixels != NULL) {
3227 cxbool is_same_size = FALSE;
3229 is_same_size = _giraffe_compare_overscans(raw, bad_pixels);
3236 if (is_same_size == FALSE) {
3237 cpl_msg_info(fctid,
"PRE/OVRSCAN Regions do not match between "
3238 "RAW Image and Bad Pixel Image.");
3253 timage = cpl_image_duplicate(_raw);
3261 error_code = _giraffe_bias_compute(&bias_results, timage,
3262 areas, _master_bias,
3263 _bad_pixels, config->method,
3264 config->option, config->model,
3265 config->remove, bias_value,
3266 config->xdeg, config->ydeg,
3267 config->xstep, config->ystep,
3268 config->sigma, config->iterations,
3271 cpl_matrix_delete(areas);
3281 _giraffe_biasresults_clear(&bias_results);
3283 cpl_msg_error(fctid,
"Bias computation failed with error %d!",
3296 if (config->remove == TRUE) {
3299 cpl_image_delete(timage);
3306 cpl_image_delete(timage);
3320 if (config->remove == TRUE) {
3322 cpl_propertylist_set_int(properties, GIALIAS_BITPIX, -32);
3323 cpl_propertylist_set_double(properties, GIALIAS_BZERO, 0.);
3324 cpl_propertylist_set_double(properties, GIALIAS_BSCALE, 1.);
3326 if (cpl_propertylist_has(properties, GIALIAS_GIRFTYPE)) {
3327 cpl_propertylist_set_string(properties,
3328 GIALIAS_GIRFTYPE,
"BRMIMG");
3331 cpl_propertylist_append_string(properties,
3332 GIALIAS_GIRFTYPE,
"BRMIMG");
3334 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3335 "GIRAFFE bias subtracted frame");
3352 if (cpl_propertylist_has(properties, GIALIAS_CONAD)) {
3353 cxdouble conad = cpl_propertylist_get_double(properties, GIALIAS_CONAD);
3359 cpl_msg_info(fctid,
"Converting bias subtracted frame to "
3360 "electrons; conversion factor is %.2f e-/ADU",
3362 cpl_image_multiply_scalar(_result, conad);
3365 cpl_msg_warning(fctid,
"Invalid conversion factor %.2f "
3366 "e-/ADU! Bias subtracted image will not be "
3367 "converted.", conad);
3371 cpl_msg_info(fctid,
"Conversion factor not found. Bias subtracted "
3372 "image will remain in ADU");
3377 s = cx_string_new();
3379 _giraffe_method_string(s, config->method, config->option);
3380 cpl_propertylist_append_string(properties, GIALIAS_BIASMETHOD, cx_string_get(s));
3382 cx_string_delete(s);
3384 cpl_propertylist_append_double(properties, GIALIAS_BIASVALUE,
3386 cpl_propertylist_append_double(properties, GIALIAS_BIASERROR,
3387 bias_results.sigma);
3388 cpl_propertylist_append_double(properties, GIALIAS_BIASSIGMA,
3391 if (bias_results.coeffs) {
3392 cpl_propertylist_append_double(properties, GIALIAS_BCLIPSIGMA,
3394 cpl_propertylist_append_int(properties, GIALIAS_BCLIPNITER,
3395 config->iterations);
3396 cpl_propertylist_append_double(properties, GIALIAS_BCLIPMFRAC,
3400 if (bias_results.limits) {
3401 cpl_propertylist_append_string(properties, GIALIAS_BIASAREAS,
3402 cx_string_get(bias_results.limits));
3405 if (bias_results.coeffs) {
3406 s = cx_string_new();
3408 _giraffe_stringify_coefficients(s, bias_results.coeffs);
3409 cpl_propertylist_append_string(properties, GIALIAS_BIASPLANE,
3412 cx_string_delete(s);
3420 _giraffe_biasresults_clear(&bias_results);
3444 GiBiasConfig* config = NULL;
3451 config = cx_calloc(1,
sizeof *config);
3458 config->method = GIBIAS_METHOD_UNDEFINED;
3459 config->option = GIBIAS_OPTION_UNDEFINED;
3460 config->model = GIBIAS_MODEL_UNDEFINED;
3466 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.remove");
3467 config->remove = cpl_parameter_get_bool(p);
3469 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.method");
3470 s = cpl_parameter_get_string(p);
3472 if (!strcmp(s,
"UNIFORM")) {
3473 config->method = GIBIAS_METHOD_UNIFORM;
3476 if (!strcmp(s,
"PLANE")) {
3477 config->method = GIBIAS_METHOD_PLANE;
3480 if (!strcmp(s,
"CURVE")) {
3481 config->method = GIBIAS_METHOD_CURVE;
3484 if (!strcmp(s,
"PROFILE")) {
3485 config->method = GIBIAS_METHOD_PROFILE;
3488 if (!strcmp(s,
"MASTER")) {
3489 config->method = GIBIAS_METHOD_MASTER;
3492 if (!strcmp(s,
"ZMASTER")) {
3493 config->method = GIBIAS_METHOD_ZMASTER;
3496 if (!strcmp(s,
"PROFILE+CURVE")) {
3497 config->method = GIBIAS_METHOD_PROFILE;
3498 config->option = GIBIAS_OPTION_CURVE;
3501 if (!strcmp(s,
"MASTER+PLANE")) {
3502 config->method = GIBIAS_METHOD_MASTER;
3503 config->option = GIBIAS_OPTION_PLANE;
3506 if (!strcmp(s,
"ZMASTER+PLANE")) {
3507 config->method = GIBIAS_METHOD_ZMASTER;
3508 config->option = GIBIAS_OPTION_PLANE;
3511 if (!strcmp(s,
"MASTER+CURVE")) {
3512 config->method = GIBIAS_METHOD_MASTER;
3513 config->option = GIBIAS_OPTION_CURVE;
3516 if (!strcmp(s,
"ZMASTER+CURVE")) {
3517 config->method = GIBIAS_METHOD_ZMASTER;
3518 config->option = GIBIAS_OPTION_CURVE;
3521 cx_assert(config->method != GIBIAS_METHOD_UNDEFINED);
3523 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.areas");
3524 s = cpl_parameter_get_string(p);
3525 config->areas = cx_strdup(s);
3527 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.sigma");
3528 config->sigma = cpl_parameter_get_double(p);
3530 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.iterations");
3531 config->iterations = cpl_parameter_get_int(p);
3533 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.fraction");
3534 config->fraction = cpl_parameter_get_double(p);
3536 if (config->method == GIBIAS_METHOD_CURVE ||
3537 config->option == GIBIAS_OPTION_CURVE) {
3538 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.xorder");
3539 config->xdeg = 1 + cpl_parameter_get_int(p);
3541 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.yorder");
3542 config->ydeg = 1 + cpl_parameter_get_int(p);
3545 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.xstep");
3546 config->xstep = cpl_parameter_get_int(p);
3548 p = cpl_parameterlist_find(list,
"giraffe.biasremoval.ystep");
3549 config->ystep = cpl_parameter_get_int(p);
3573 if (config->areas) {
3574 cx_free(config->areas);
3575 config->areas = NULL;
3607 p = cpl_parameter_new_value(
"giraffe.biasremoval.remove",
3609 "Enable bias removal",
3610 "giraffe.biasremoval",
3612 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"remove-bias");
3613 cpl_parameterlist_append(list, p);
3616 p = cpl_parameter_new_enum(
"giraffe.biasremoval.method",
3618 "Bias removal method",
3619 "giraffe.biasremoval",
3620 "PROFILE", 11,
"UNIFORM",
"PLANE",
"CURVE",
3621 "PROFILE",
"MASTER",
"ZMASTER",
"PROFILE+CURVE",
3622 "MASTER+PLANE",
"MASTER+CURVE",
"ZMASTER+PLANE",
3624 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-method");
3625 cpl_parameterlist_append(list, p);
3628 p = cpl_parameter_new_value(
"giraffe.biasremoval.areas",
3630 "Bias areas to use "
3631 "(Xl0:Xr0:Yl0:Yr0, ... ,Xln:Xrn:Yln:Yrn)",
3632 "giraffe.biasremoval",
3634 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-areas");
3635 cpl_parameterlist_append(list, p);
3638 p = cpl_parameter_new_value(
"giraffe.biasremoval.sigma",
3640 "Sigma Clipping: sigma threshold factor",
3641 "giraffe.biasremoval",
3643 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-sigma");
3644 cpl_parameterlist_append(list, p);
3647 p = cpl_parameter_new_value(
"giraffe.biasremoval.iterations",
3649 "Sigma Clipping: maximum number of "
3651 "giraffe.biasremoval",
3653 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-niter");
3654 cpl_parameterlist_append(list, p);
3657 p = cpl_parameter_new_value(
"giraffe.biasremoval.fraction",
3659 "Sigma Clipping: minimum fraction of points "
3660 "accepted/total [0.0..1.0]",
3661 "giraffe.biasremoval",
3663 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-mfrac");
3664 cpl_parameterlist_append(list, p);
3667 p = cpl_parameter_new_value(
"giraffe.biasremoval.xorder",
3669 "Order of X polynomial fit (CURVE method "
3671 "giraffe.biasremoval",
3673 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-xorder");
3674 cpl_parameterlist_append(list, p);
3676 p = cpl_parameter_new_value(
"giraffe.biasremoval.yorder",
3678 "Order of Y polynomial fit (CURVE method "
3680 "giraffe.biasremoval",
3682 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-yorder");
3683 cpl_parameterlist_append(list, p);
3686 p = cpl_parameter_new_value(
"giraffe.biasremoval.xstep",
3688 "Sampling step along X (CURVE method only)",
3689 "giraffe.biasremoval",
3691 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-xstep");
3692 cpl_parameterlist_append(list, p);
3695 p = cpl_parameter_new_value(
"giraffe.biasremoval.ystep",
3697 "Sampling step along Y (CURVE method only)",
3698 "giraffe.biasremoval",
3700 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI,
"bsremove-ystep");
3701 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.
cpl_matrix * giraffe_get_raw_areas(const GiImage *image)
Create bias areas from an image.
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.
void giraffe_bias_config_destroy(GiBiasConfig *config)
Destroys a bias removal setup structure.
cpl_propertylist * giraffe_image_get_properties(const GiImage *self)
Get the properties of an image.
cpl_image * giraffe_image_get(const GiImage *self)
Gets the image data.
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.
cxdouble giraffe_matrix_sigma_fit(const cpl_matrix *matrix, const cpl_matrix *matrix_fit)
Compute sigma of matrix fit.
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.
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.