00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 #ifdef HAVE_CONFIG_H
00029 # include <config.h>
00030 #endif
00031
00032 #include <string.h>
00033 #include <math.h>
00034 #include <errno.h>
00035
00036 #include <cxmemory.h>
00037 #include <cxstrutils.h>
00038 #include <cxutils.h>
00039
00040 #include <cpl_msg.h>
00041 #include <cpl_parameter.h>
00042
00043 #include "gimacros.h"
00044 #include "gialias.h"
00045 #include "giarray.h"
00046 #include "gimatrix.h"
00047 #include "gichebyshev.h"
00048 #include "gimath.h"
00049 #include "gibias.h"
00050
00051
00060 struct GiBiasResults {
00061
00062 cxdouble mean;
00063 cxdouble sigma;
00064
00065 cx_string* limits;
00066
00067 cpl_matrix* coeffs;
00068
00069 cpl_image* model;
00070
00071 };
00072
00073 typedef struct GiBiasResults GiBiasResults;
00074
00075
00076
00077
00078
00079
00080 inline static void
00081 _giraffe_biasresults_clear(GiBiasResults *self)
00082 {
00083
00084 if (self != NULL) {
00085
00086 self->mean = 0.;
00087 self->sigma = 0.;
00088
00089 if (self->limits) {
00090 cx_string_delete(self->limits);
00091 self->limits = NULL;
00092 }
00093
00094 if (self->coeffs) {
00095 cpl_matrix_delete(self->coeffs);
00096 self->coeffs = NULL;
00097 }
00098
00099 if (self->model) {
00100 cpl_image_delete(self->model);
00101 self->model = NULL;
00102 }
00103
00104 }
00105
00106 return;
00107
00108 }
00109
00110
00111
00112
00113
00114
00115 inline static void
00116 _giraffe_method_string(cx_string *string, GiBiasMethod method,
00117 GiBiasOption option)
00118 {
00119
00120 switch (method) {
00121 case GIBIAS_METHOD_UNIFORM:
00122 cx_string_set(string, "UNIFORM");
00123 break;
00124
00125 case GIBIAS_METHOD_PLANE:
00126 cx_string_set(string, "PLANE");
00127 break;
00128
00129 case GIBIAS_METHOD_CURVE:
00130 cx_string_set(string, "CURVE");
00131 break;
00132
00133 case GIBIAS_METHOD_PROFILE:
00134 cx_string_set(string, "PROFILE");
00135 break;
00136
00137 case GIBIAS_METHOD_MASTER:
00138 cx_string_set(string, "MASTER");
00139 break;
00140
00141 case GIBIAS_METHOD_ZMASTER:
00142 cx_string_set(string, "ZMASTER");
00143 break;
00144
00145 default:
00146 break;
00147 }
00148
00149 if (option != GIBIAS_OPTION_UNDEFINED) {
00150 switch (option) {
00151 case GIBIAS_OPTION_PLANE:
00152 cx_string_append(string, "+PLANE");
00153 break;
00154
00155 case GIBIAS_OPTION_CURVE:
00156 cx_string_append(string, "+CURVE");
00157 break;
00158
00159 default:
00160 break;
00161 }
00162 }
00163
00164 return;
00165
00166 }
00167
00168
00169 inline static void
00170 _giraffe_stringify_coefficients(cx_string *string, cpl_matrix *matrix)
00171 {
00172
00173 register cxint i, j;
00174
00175 cxchar *tmp = cx_line_alloc();
00176
00177 cxint nr = cpl_matrix_get_nrow(matrix);
00178 cxint nc = cpl_matrix_get_ncol(matrix);
00179
00180 cxsize sz = cx_line_max();
00181
00182 cxdouble *data = cpl_matrix_get_data(matrix);
00183
00184
00185 for (i = 0; i < nr; i++) {
00186 for (j = 0; j < nc; j++) {
00187 snprintf(tmp, sz, "%g", data[i * nc + j]);
00188 cx_string_append(string, tmp);
00189
00190 if (i != nr - 1 || j < nc - 1) {
00191 cx_string_append(string, ":");
00192 }
00193 }
00194 }
00195
00196 cx_free(tmp);
00197
00198 return;
00199
00200 }
00201
00202
00218 inline static cxbool
00219 _giraffe_compare_overscans(const GiImage* image1, const GiImage* image2)
00220 {
00221
00222 cxint32 l1ovscx = -1;
00223 cxint32 l1ovscy = -1;
00224 cxint32 l1prscx = -1;
00225 cxint32 l1prscy = -1;
00226 cxint32 l2ovscx = -1;
00227 cxint32 l2ovscy = -1;
00228 cxint32 l2prscx = -1;
00229 cxint32 l2prscy = -1;
00230
00231 cpl_propertylist *l1, *l2;
00232
00233
00234 cx_assert(image1 != NULL && image2 != NULL);
00235
00236 l1 = giraffe_image_get_properties(image1);
00237 l2 = giraffe_image_get_properties(image2);
00238
00239 if (cpl_propertylist_has(l1, GIALIAS_OVSCX)) {
00240 l1ovscx = cpl_propertylist_get_int(l1, GIALIAS_OVSCX);
00241 }
00242 if (cpl_propertylist_has(l1, GIALIAS_OVSCY)) {
00243 l1ovscy = cpl_propertylist_get_int(l1, GIALIAS_OVSCY);
00244 }
00245 if (cpl_propertylist_has(l1, GIALIAS_PRSCX)) {
00246 l1prscx = cpl_propertylist_get_int(l1, GIALIAS_PRSCX);
00247 }
00248 if (cpl_propertylist_has(l1, GIALIAS_PRSCY)) {
00249 l1prscy = cpl_propertylist_get_int(l1, GIALIAS_PRSCY);
00250 }
00251
00252 if (cpl_propertylist_has(l2, GIALIAS_OVSCX)) {
00253 l2ovscx = cpl_propertylist_get_int(l2, GIALIAS_OVSCX);
00254 }
00255 if (cpl_propertylist_has(l2, GIALIAS_OVSCY)) {
00256 l2ovscy = cpl_propertylist_get_int(l2, GIALIAS_OVSCY);
00257 }
00258 if (cpl_propertylist_has(l2, GIALIAS_PRSCX)) {
00259 l2prscx = cpl_propertylist_get_int(l2, GIALIAS_PRSCX);
00260 }
00261 if (cpl_propertylist_has(l2, GIALIAS_PRSCY)) {
00262 l2prscy = cpl_propertylist_get_int(l2, GIALIAS_PRSCY);
00263 }
00264
00265 if (l1ovscx != l2ovscx || l1ovscy != l2ovscy) {
00266 return FALSE;
00267 }
00268
00269 if (l1prscx != l2prscx || l1prscy != l2prscy) {
00270 return FALSE;
00271 }
00272
00273 return TRUE;
00274
00275 }
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291 inline static cpl_matrix*
00292 _giraffe_bias_get_areas(const cxchar* string)
00293 {
00294
00295 cxchar** regions = NULL;
00296
00297 cpl_matrix* areas = NULL;
00298
00299
00300 cx_assert(string != NULL);
00301
00302 regions = cx_strsplit(string, ",", -1);
00303
00304 if (regions == NULL) {
00305
00306 return NULL;
00307
00308 }
00309 else {
00310
00311 const cxsize nvalues = 4;
00312
00313 cxsize i = 0;
00314 cxsize nregions = 0;
00315
00316 while (regions[i] != NULL) {
00317 ++i;
00318 }
00319
00320 nregions = i;
00321 areas = cpl_matrix_new(nregions, nvalues);
00322
00323 i = 0;
00324 while (regions[i] != NULL) {
00325
00326 register cxsize j = 0;
00327
00328 cxchar** limits = cx_strsplit(regions[i], ":", nvalues);
00329
00330
00331 if (limits == NULL) {
00332
00333 cpl_matrix_delete(areas);
00334 areas = NULL;
00335
00336 return NULL;
00337
00338 }
00339
00340 for (j = 0; j < nvalues; ++j) {
00341
00342 cxchar* last = NULL;
00343
00344 cxint status = 0;
00345 cxlong value = 0;
00346
00347
00348 if (limits[j] == NULL || *limits[j] == '\0') {
00349 break;
00350 }
00351
00352 status = errno;
00353 errno = 0;
00354
00355 value = strtol(limits[j], &last, 10);
00356
00357
00358
00359
00360
00361 if ((errno == ERANGE &&
00362 (value == LONG_MAX || value == LONG_MIN)) ||
00363 (errno != 0 && value == 0)) {
00364
00365 break;
00366
00367 }
00368
00369 errno = status;
00370
00371
00372
00373
00374
00375
00376 if (*last != '\0') {
00377 break;
00378 }
00379
00380 cpl_matrix_set(areas, i, j, value);
00381
00382 }
00383
00384 cx_strfreev(limits);
00385 limits = NULL;
00386
00387 if (j != nvalues) {
00388
00389 cpl_matrix_delete(areas);
00390 areas = NULL;
00391
00392 cx_strfreev(regions);
00393 regions = NULL;
00394
00395 return NULL;
00396
00397 }
00398
00399 ++i;
00400
00401 }
00402
00403 cx_strfreev(regions);
00404 regions = NULL;
00405
00406 }
00407
00408 return areas;
00409
00410 }
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426
00427
00428
00429
00430
00431
00432
00433
00434
00435
00436
00437
00438
00439
00440
00441
00442
00443
00444
00445
00446
00447
00448
00449
00450
00451
00452
00453
00454
00455
00456
00457
00458
00459 inline static cxint
00460 _giraffe_bias_compute_mean(GiBiasResults* results, const cpl_image* image,
00461 const cpl_matrix* areas, cxdouble sigma, cxint numiter,
00462 cxdouble maxfraction)
00463 {
00464
00465 const cxchar* const fctid = "giraffe_bias_compute_mean";
00466
00467
00468 const cxdouble* pdimg = NULL;
00469
00470 cxint32 j, l, k;
00471 cxint32 img_dimx, img_dimy;
00472 cxint32 ba_num_cols;
00473 cxint32 ba_num_rows;
00474 cxint32 curriter = 0;
00475 cxint32 x0, x1, x2, x3;
00476
00477 cxlong ntotal = 0L;
00478 cxlong naccepted = 0L;
00479 cxlong n = 0L;
00480 cxlong pixcount = 0L;
00481
00482 cxdouble currfraction = 2.0;
00483
00484 cx_string* tmp = NULL;
00485
00486 cpl_matrix* matrix_zz1;
00487 cpl_matrix* matrix_zz1diff;
00488
00489
00490
00491
00492
00493
00494 if (results->limits == NULL) {
00495 cpl_msg_info(fctid, "Unable to store biaslimits return parameter, "
00496 "aborting...");
00497 return -3;
00498 }
00499
00500 if (cpl_image_get_type(image) != CPL_TYPE_DOUBLE) {
00501 cpl_msg_info(fctid, "Only images allowed of type double, "
00502 "aborting ...");
00503 return -3;
00504 }
00505
00506
00507
00508
00509
00510
00511
00512
00513
00514
00515 if (areas == NULL) {
00516 cpl_msg_info(fctid, "Bias Areas: Missing bias areas, "
00517 "aborting ...");
00518 return -1;
00519 }
00520
00521 ba_num_cols = cpl_matrix_get_ncol(areas);
00522 ba_num_rows = cpl_matrix_get_nrow(areas);
00523
00524 for (j = 0; j < ba_num_rows; j++) {
00525 x3 = (cxint32) cpl_matrix_get(areas, j, 3);
00526 x2 = (cxint32) cpl_matrix_get(areas, j, 2);
00527 x1 = (cxint32) cpl_matrix_get(areas, j, 1);
00528 x0 = (cxint32) cpl_matrix_get(areas, j, 0);
00529
00530 ntotal += (cxulong) ((x3 - x2 + 1) * (x1 - x0 + 1));
00531 }
00532
00533 if (ntotal <= 0) {
00534 cpl_msg_info(fctid, "Bias Areas: Inconsistent specification, "
00535 "aborting ...");
00536 return -1;
00537 }
00538
00539 matrix_zz1 = cpl_matrix_new(ntotal, 1);
00540 matrix_zz1diff = cpl_matrix_new(ntotal, 1);
00541
00542
00543
00544
00545
00546
00547 img_dimx = cpl_image_get_size_x(image);
00548 img_dimy = cpl_image_get_size_y(image);
00549
00550 cx_string_set(results->limits, "");
00551
00552 for (j = 0; j < ba_num_rows; j++) {
00553
00554 x3 = (cxint32) cpl_matrix_get(areas, j, 3);
00555 x2 = (cxint32) cpl_matrix_get(areas, j, 2);
00556 x1 = (cxint32) cpl_matrix_get(areas, j, 1);
00557 x0 = (cxint32) cpl_matrix_get(areas, j, 0);
00558
00559 if ((x0 > img_dimx) || (x1 > img_dimx) ||
00560 (x2 > img_dimy) || (x3 > img_dimy)) {
00561 continue;
00562 }
00563
00564 tmp = cx_string_new();
00565
00566 cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
00567 cx_string_append(results->limits, cx_string_get(tmp));
00568
00569 cx_string_delete(tmp);
00570 tmp = NULL;
00571
00572 pdimg = cpl_image_get_data_double_const(image);
00573
00574 for (l = x2; l < x3 + 1; l++) {
00575 for (k = x0; k < x1 + 1; k++) {
00576 cpl_matrix_set(matrix_zz1, n, 1, pdimg[k + l * img_dimx]);
00577 n++;
00578 }
00579 }
00580 }
00581
00582 if (n != ntotal) {
00583 cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting ...");
00584
00585 cpl_matrix_delete(matrix_zz1);
00586 cpl_matrix_delete(matrix_zz1diff);
00587
00588 return -3;
00589 }
00590
00591 cpl_msg_info(fctid, "Bias Areas: Using %s",
00592 cx_string_get(results->limits));
00593
00594
00595
00596
00597
00598
00599
00600
00601
00602
00603 cpl_msg_info(fctid, "Sigma Clipping : Start");
00604
00605 naccepted = ntotal;
00606
00607 while ((naccepted > 0) && (curriter < numiter) &&
00608 (currfraction > maxfraction)) {
00609
00610 results->mean = cpl_matrix_get_mean(matrix_zz1);
00611 results->sigma = giraffe_matrix_sigma_mean(matrix_zz1, results->mean);
00612
00613 for (k = 0; k < cpl_matrix_get_nrow(matrix_zz1); k++) {
00614 cpl_matrix_set(matrix_zz1diff, k, 0,
00615 cpl_matrix_get(matrix_zz1diff, k, 0) -
00616 results->mean);
00617 }
00618
00619 cpl_msg_info(fctid, "bias[%d]: mean = %5g, sigma = %6.3g, "
00620 "ratio = %6.3g, accepted = %ld\n", curriter,
00621 results->mean, results->sigma, currfraction, naccepted);
00622
00623 pixcount = 0L;
00624 for (l = 0; l < cpl_matrix_get_nrow(matrix_zz1); l++) {
00625 if (fabs(cpl_matrix_get(matrix_zz1diff, l, 0)) >
00626 (results->sigma * sigma)) {
00627 continue;
00628 }
00629
00630 cpl_matrix_set(matrix_zz1, pixcount, 0,
00631 cpl_matrix_get(matrix_zz1, l, 0));
00632 pixcount++;
00633 }
00634
00635 cpl_matrix_resize(matrix_zz1, 0, 0, 0, pixcount -
00636 cpl_matrix_get_nrow(matrix_zz1));
00637
00638 cpl_matrix_resize(matrix_zz1diff, 0, 0, 0, pixcount -
00639 cpl_matrix_get_nrow(matrix_zz1diff));
00640
00641 if (pixcount == naccepted) {
00642 break;
00643 }
00644
00645 naccepted = pixcount;
00646
00647 currfraction = (cxdouble) naccepted / (cxdouble) ntotal;
00648 curriter++;
00649 }
00650
00651 cpl_msg_info(fctid, "Sigma Clipping : End");
00652
00653
00654
00655
00656
00657
00658 results->mean = cpl_matrix_get_mean(matrix_zz1);
00659 results->sigma = giraffe_matrix_sigma_mean(matrix_zz1, results->mean);
00660
00661 cpl_msg_info(fctid, "Sigma Clipping Results : bias[%d]: mean = %5g, "
00662 "sigma = %6.3g, ratio = %6.3g, accepted=%ld\n", curriter,
00663 results->mean, results->sigma, currfraction, naccepted);
00664
00665
00666
00667
00668
00669
00670 if (matrix_zz1 != NULL) {
00671 cpl_matrix_delete(matrix_zz1);
00672 }
00673
00674 if (matrix_zz1diff != NULL) {
00675 cpl_matrix_delete(matrix_zz1diff);
00676 }
00677
00678 return EXIT_SUCCESS;
00679
00680 }
00681
00682
00683
00684
00685
00686
00687
00688
00689
00690
00691
00692
00693
00694
00695
00696
00697
00698
00699
00700
00701
00702
00703
00704
00705
00706
00707
00708
00709
00710
00711
00712
00713
00714
00715
00716
00717
00718
00719
00720
00721
00722
00723
00724
00725
00726
00727
00728
00729
00730
00731
00732
00733
00734
00735 inline static cxint
00736 _giraffe_bias_compute_plane(GiBiasResults* results, const cpl_image* image,
00737 const cpl_matrix* areas, cxdouble kappa,
00738 cxint numiter, cxdouble maxfraction)
00739 {
00740
00741 const cxchar* const fctid = "giraffe_bias_compute_plane";
00742
00743
00744 cxint j = 0;
00745 cxint nx = 0;
00746 cxint ny = 0;
00747 cxint nareas = 0;
00748 cxint iteration = 0;
00749
00750 cxsize n = 0;
00751 cxsize ntotal = 0;
00752 cxsize naccepted = 0;
00753
00754 cxdouble fraction = 1.;
00755
00756 cpl_matrix* xbs = NULL;
00757 cpl_matrix* ybs = NULL;
00758 cpl_matrix* zbs = NULL;
00759 cpl_matrix* coeffs = NULL;
00760
00761
00762 cx_assert(results->limits != NULL);
00763 cx_assert(results->coeffs == NULL);
00764
00765 cx_assert(areas != NULL);
00766
00767 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
00768
00769
00770
00771
00772
00773
00774 nareas = cpl_matrix_get_nrow(areas);
00775
00776 for (j = 0; j < nareas; j++) {
00777
00778 cxint x3 = (cxint) cpl_matrix_get(areas, j, 3);
00779 cxint x2 = (cxint) cpl_matrix_get(areas, j, 2);
00780 cxint x1 = (cxint) cpl_matrix_get(areas, j, 1);
00781 cxint x0 = (cxint) cpl_matrix_get(areas, j, 0);
00782
00783 ntotal += (cxsize) ((x3 - x2 + 1) * (x1 - x0 + 1));
00784
00785 }
00786
00787 if (ntotal <= 0) {
00788
00789 cpl_msg_info(fctid, "Bias Areas: Inconsistent specification, "
00790 "aborting ...");
00791 return -1;
00792
00793 }
00794
00795 nx = cpl_image_get_size_x(image);
00796 ny = cpl_image_get_size_y(image);
00797
00798 cpl_msg_info(fctid, "Bias Areas: specified are %zd points in %dx%d "
00799 "image", ntotal, nx, ny);
00800
00801
00802
00803
00804
00805
00806 results->mean = 0.;
00807
00808 cx_string_set(results->limits, "");
00809
00810 xbs = cpl_matrix_new(ntotal, 1);
00811 ybs = cpl_matrix_new(ntotal, 1);
00812 zbs = cpl_matrix_new(1, ntotal);
00813
00814 for (j = 0; j < nareas; ++j) {
00815
00816 const cxdouble* _img = cpl_image_get_data_double_const(image);
00817
00818 cxint k = 0;
00819
00820 cx_string* tmp = NULL;
00821
00822
00823 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
00824 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
00825 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
00826 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
00827
00828 if ((x0 > nx) || (x1 > nx) || (x2 > ny) || (x3 > ny)) {
00829 continue;
00830 }
00831
00832 tmp = cx_string_new();
00833
00834 cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
00835 cx_string_append(results->limits, cx_string_get(tmp));
00836
00837 cx_string_delete(tmp);
00838 tmp = NULL;
00839
00840 for (k = x2; k < x3 + 1; ++k) {
00841
00842 register cxint l = 0;
00843
00844
00845 for (l = x0; l < x1 + 1; ++l) {
00846
00847 cpl_matrix_set(xbs, n, 0, l);
00848 cpl_matrix_set(ybs, n, 0, k);
00849 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
00850 ++n;
00851
00852 }
00853
00854 }
00855
00856 }
00857
00858 cpl_matrix_set_size(xbs, n, 1);
00859 cpl_matrix_set_size(ybs, n, 1);
00860 cpl_matrix_set_size(zbs, 1, n);
00861
00862 if (n != ntotal) {
00863
00864 cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting...");
00865
00866 cpl_matrix_delete(xbs);
00867 cpl_matrix_delete(ybs);
00868 cpl_matrix_delete(zbs);
00869
00870 return -1;
00871
00872 }
00873
00874 ntotal = n;
00875
00876 cpl_msg_info(fctid, "Bias Areas: Using %s [%zd pixels]",
00877 cx_string_get(results->limits), ntotal);
00878
00879
00880
00881
00882
00883
00884 cpl_msg_info(fctid, "Sigma Clipping : Start");
00885
00886 naccepted = ntotal;
00887
00888 while ((naccepted > 0) && (iteration < numiter) &&
00889 (fraction > maxfraction)) {
00890
00891 cxsize k = 0;
00892
00893 cxdouble sigma = 0.;
00894
00895 cpl_matrix* base = NULL;
00896 cpl_matrix* fit = NULL;
00897
00898
00899 base = cpl_matrix_new(3, naccepted);
00900
00901 if (base == NULL) {
00902
00903 cpl_msg_info(fctid, "Sigma Clipping: Error creating design "
00904 "matrix");
00905
00906 cpl_matrix_delete(zbs);
00907 cpl_matrix_delete(ybs);
00908 cpl_matrix_delete(xbs);
00909
00910 return -2;
00911 }
00912
00913 for (k = 0; k < naccepted; ++k) {
00914
00915 cpl_matrix_set(base, 0, k, 1.);
00916 cpl_matrix_set(base, 1, k, cpl_matrix_get(xbs, k, 0));
00917 cpl_matrix_set(base, 2, k, cpl_matrix_get(ybs, k, 0));
00918
00919 }
00920
00921 cpl_matrix_delete(coeffs);
00922 coeffs = NULL;
00923
00924 coeffs = giraffe_matrix_leastsq(base, zbs);
00925
00926 if (coeffs == NULL) {
00927
00928 cpl_msg_info(fctid, "Sigma Clipping : Error in least square "
00929 "solution, aborting...");
00930
00931 cpl_matrix_delete(base);
00932 base = NULL;
00933
00934 cpl_matrix_delete(xbs);
00935 cpl_matrix_delete(ybs);
00936 cpl_matrix_delete(zbs);
00937
00938 return -2;
00939
00940 }
00941
00942
00943
00944
00945
00946
00947 fit = cpl_matrix_product_create(coeffs, base);
00948
00949 cpl_matrix_delete(base);
00950 base = NULL;
00951
00952 results->mean = cpl_matrix_get_mean(fit);
00953
00954 sigma = giraffe_matrix_sigma_fit(zbs, fit);
00955
00956 cpl_msg_info(fctid, "Sigma Clipping : bias plane[%d]: %g + "
00957 "%g * x + %g * y, sigma = %.5g, ratio = %.4g, "
00958 "accepted = %zd\n", iteration,
00959 cpl_matrix_get(coeffs, 0, 0),
00960 cpl_matrix_get(coeffs, 0, 1),
00961 cpl_matrix_get(coeffs, 0, 2),
00962 sigma, fraction, naccepted);
00963
00964
00965
00966
00967
00968
00969 sigma *= kappa;
00970
00971 n = 0;
00972
00973 for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
00974
00975 register cxdouble z = cpl_matrix_get(zbs, 0, j);
00976
00977 if (fabs(cpl_matrix_get(fit, 0, j) - z) > sigma) {
00978 continue;
00979 }
00980
00981 cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
00982
00983 cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
00984
00985 cpl_matrix_set(zbs, 0, n, z);
00986 ++n;
00987
00988 }
00989
00990 cpl_matrix_set_size(xbs, n, 1);
00991 cpl_matrix_set_size(ybs, n, 1);
00992 cpl_matrix_set_size(zbs, 1, n);
00993
00994 cpl_matrix_delete(fit);
00995 fit = NULL;
00996
00997 if (n == naccepted) {
00998 break;
00999 }
01000
01001 naccepted = n;
01002
01003 fraction = (cxdouble)naccepted / (cxdouble)ntotal;
01004 ++iteration;
01005
01006 }
01007
01008 cpl_msg_info(fctid, "Sigma Clipping : End");
01009
01010
01011
01012
01013
01014
01015 results->coeffs = coeffs;
01016 results->sigma = cpl_matrix_get_stdev(zbs);
01017
01018 cpl_msg_info(fctid, "Sigma Clipping Results (%d/%zd, sigma = %g)",
01019 iteration, naccepted, results->sigma);
01020
01021
01022
01023
01024
01025
01026 cpl_matrix_delete(xbs);
01027 xbs = NULL;
01028
01029 cpl_matrix_delete(ybs);
01030 ybs = NULL;
01031
01032 cpl_matrix_delete(zbs);
01033 zbs = NULL;
01034
01035 return EXIT_SUCCESS;
01036
01037 }
01038
01039
01040
01041
01042
01043
01044
01045
01046
01047
01048
01049
01050
01051
01052
01053
01054
01055
01056
01057
01058
01059
01060
01061
01062
01063
01064
01065
01066
01067
01068
01069
01070
01071
01072
01073
01074
01075
01076
01077
01078
01079
01080
01081 inline static cxint
01082 _giraffe_bias_compute_curve(GiBiasResults* results, const cpl_image* image,
01083 const cpl_matrix *areas, cxdouble kappa,
01084 cxint numiter, cxdouble maxfraction,
01085 cxdouble xdeg, cxdouble ydeg,
01086 cxdouble xstep, cxdouble ystep)
01087 {
01088
01089 const cxchar* const fctid = "giraffe_bias_compute_curve";
01090
01091 cxint j = 0;
01092 cxint nx = 0;
01093 cxint ny = 0;
01094 cxint nareas = 0;
01095 cxint iteration = 0;
01096
01097 cxsize n = 0;
01098 cxsize ntotal = 0;
01099 cxsize naccepted = 0;
01100
01101 cxdouble fraction = 1.;
01102
01103 cpl_matrix* xbs = NULL;
01104 cpl_matrix* ybs = NULL;
01105 cpl_matrix* zbs = NULL;
01106
01107 GiChebyshev2D* fit = NULL;
01108
01109
01110 cx_assert(results != NULL);
01111 cx_assert(results->limits != NULL);
01112 cx_assert(results->coeffs == NULL);
01113
01114 cx_assert(areas != NULL);
01115
01116 cx_assert(image != NULL);
01117 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
01118
01119
01120
01121
01122
01123
01124 nareas = cpl_matrix_get_nrow(areas);
01125
01126 for (j = 0; j < nareas; ++j) {
01127
01128 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
01129 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
01130 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
01131 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
01132
01133 ntotal += (cxulong) (ceil((1. + x3 - x2) / ystep) *
01134 ceil((1. + x1 - x0) / xstep));
01135 }
01136
01137 nx = cpl_image_get_size_x(image);
01138 ny = cpl_image_get_size_y(image);
01139
01140 cpl_msg_info(fctid, "Bias Areas: Found %zd points in %dx%d image",
01141 ntotal, nx, ny);
01142
01143
01144
01145
01146
01147
01148 results->mean = 0.;
01149
01150 cx_string_set(results->limits, "");
01151
01152 xbs = cpl_matrix_new(ntotal, 1);
01153 ybs = cpl_matrix_new(ntotal, 1);
01154 zbs = cpl_matrix_new(1, ntotal);
01155
01156 for (j = 0; j < nareas; ++j) {
01157
01158 const cxdouble* _img = cpl_image_get_data_double_const(image);
01159
01160 cxint k = 0;
01161
01162 cx_string* tmp = NULL;
01163
01164
01165 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
01166 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
01167 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
01168 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
01169
01170 if ((x0 > nx) || (x1 > nx) ||
01171 (x2 > ny) || (x3 > ny)) {
01172 continue;
01173 }
01174
01175 tmp = cx_string_new();
01176
01177 cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
01178 cx_string_append(results->limits, cx_string_get(tmp));
01179
01180 cx_string_delete(tmp);
01181 tmp = NULL;
01182
01183 for (k = x2; k < x3 + 1; k += ystep) {
01184
01185 register cxint l = 0;
01186
01187
01188 for (l = x0; l < x1 + 1; l += xstep) {
01189
01190 cpl_matrix_set(xbs, n, 0, l);
01191 cpl_matrix_set(ybs, n, 0, k);
01192 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
01193 ++n;
01194
01195 }
01196
01197 }
01198
01199 }
01200
01201 cpl_matrix_set_size(xbs, n, 1);
01202 cpl_matrix_set_size(ybs, n, 1);
01203 cpl_matrix_set_size(zbs, 1, n);
01204
01205 if (n <= 0) {
01206
01207 cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting...");
01208
01209 cpl_matrix_delete(xbs);
01210 cpl_matrix_delete(ybs);
01211 cpl_matrix_delete(zbs);
01212
01213 return -1;
01214
01215 }
01216
01217 ntotal = n;
01218
01219 cpl_msg_info(fctid, "Bias Areas: Using %s [%zd pixels]",
01220 cx_string_get(results->limits), ntotal);
01221
01222
01223
01224
01225
01226
01227 cpl_msg_info(fctid, "Sigma Clipping : Start");
01228
01229 naccepted = ntotal;
01230
01231 while ((naccepted > 0) && (iteration < numiter) &&
01232 (fraction > maxfraction)) {
01233
01234 cxint status = 0;
01235
01236 cxdouble sigma = 0.;
01237
01238 cpl_matrix* base = NULL;
01239 cpl_matrix* coeffs = NULL;
01240 cpl_matrix* _coeffs = NULL;
01241 cpl_matrix* _fit = NULL;
01242
01243
01244 base = giraffe_chebyshev_base2d(0., 0., nx, ny,
01245 xdeg, ydeg, xbs, ybs);
01246
01247 if (base == NULL) {
01248
01249 cpl_msg_info(fctid, "Sigma Clipping: Error creating design "
01250 "matrix");
01251
01252 cpl_matrix_delete(zbs);
01253 cpl_matrix_delete(ybs);
01254 cpl_matrix_delete(xbs);
01255
01256 return -2;
01257 }
01258
01259 _coeffs = giraffe_matrix_leastsq(base, zbs);
01260
01261 cpl_matrix_delete(base);
01262 base = NULL;
01263
01264 if (_coeffs == NULL) {
01265
01266 cpl_msg_info(fctid, "Sigma Clipping : Error in least square "
01267 "solution, aborting...");
01268
01269 cpl_matrix_delete(xbs);
01270 cpl_matrix_delete(ybs);
01271 cpl_matrix_delete(zbs);
01272
01273 return -2;
01274
01275 }
01276
01277
01278
01279
01280
01281
01282 coeffs = cpl_matrix_wrap(xdeg, ydeg, cpl_matrix_get_data(_coeffs));
01283
01284 if (fit != NULL) {
01285 giraffe_chebyshev2d_delete(fit);
01286 fit = NULL;
01287 }
01288
01289 fit = giraffe_chebyshev2d_new(xdeg - 1, ydeg - 1);
01290 status = giraffe_chebyshev2d_set(fit, 0., nx, 0., ny,
01291 coeffs);
01292
01293 if (status != 0) {
01294
01295 giraffe_chebyshev2d_delete(fit);
01296 fit = NULL;
01297
01298 cpl_matrix_unwrap(coeffs);
01299 coeffs = NULL;
01300
01301 cpl_matrix_delete(_coeffs);
01302 _coeffs = NULL;
01303
01304 cpl_matrix_delete(xbs);
01305 cpl_matrix_delete(ybs);
01306 cpl_matrix_delete(zbs);
01307
01308 return -2;
01309
01310 }
01311
01312 cpl_matrix_unwrap(coeffs);
01313 coeffs = NULL;
01314
01315 cpl_matrix_delete(_coeffs);
01316 _coeffs = NULL;
01317
01318 _fit = cpl_matrix_new(1, cpl_matrix_get_ncol(zbs));
01319
01320 for (j = 0; j < cpl_matrix_get_ncol(_fit); ++j) {
01321
01322 cxdouble x = cpl_matrix_get(xbs, n, 0);
01323 cxdouble y = cpl_matrix_get(ybs, n, 0);
01324 cxdouble z = giraffe_chebyshev2d_eval(fit, x, y);
01325
01326 cpl_matrix_set(_fit, 0, j, z);
01327
01328 }
01329
01330 results->mean = cpl_matrix_get_mean(_fit);
01331
01332 sigma = giraffe_matrix_sigma_fit(zbs, _fit);
01333
01334 cpl_msg_info(fctid, "Sigma Clipping : bias surface[%d]: "
01335 "sigma = %8.5g, ratio = %7.4g, accepted = %zd\n",
01336 iteration, sigma, fraction, naccepted);
01337
01338
01339
01340
01341
01342
01343 sigma *= kappa;
01344
01345 n = 0;
01346
01347 for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
01348
01349 register cxdouble z = cpl_matrix_get(zbs, 0, j);
01350
01351 if (fabs(cpl_matrix_get(_fit, 0, j) - z) > sigma) {
01352 continue;
01353 }
01354
01355 cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
01356 cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
01357 cpl_matrix_set(zbs, 0, n, z);
01358 ++n;
01359
01360 }
01361
01362 cpl_matrix_set_size(xbs, n, 1);
01363 cpl_matrix_set_size(ybs, n, 1);
01364 cpl_matrix_set_size(zbs, 1, n);
01365
01366 cpl_matrix_delete(_fit);
01367 _fit = NULL;
01368
01369
01370 if (n == naccepted) {
01371 break;
01372 }
01373
01374 naccepted = n;
01375
01376 fraction = (cxdouble)naccepted / (cxdouble)ntotal;
01377 ++iteration;
01378
01379 }
01380
01381 cpl_msg_info(fctid, "Sigma Clipping : End");
01382
01383
01384
01385
01386
01387
01388 results->coeffs = cpl_matrix_duplicate(giraffe_chebyshev2d_coeffs(fit));
01389 results->sigma = cpl_matrix_get_stdev(zbs);
01390
01391 cpl_msg_info(fctid, "Sigma Clipping Results (%d/%zd, sigma = %g)",
01392 iteration, naccepted, results->sigma);
01393
01394
01395
01396
01397
01398
01399 giraffe_chebyshev2d_delete(fit);
01400 fit = NULL;
01401
01402 cpl_matrix_delete(xbs);
01403 xbs = NULL;
01404
01405 cpl_matrix_delete(ybs);
01406 ybs = NULL;
01407
01408 cpl_matrix_delete(zbs);
01409 zbs = NULL;
01410
01411 return EXIT_SUCCESS;
01412
01413 }
01414
01415
01416
01417
01418
01419
01420
01421
01422
01423
01424
01425
01426
01427
01428
01429
01430
01431
01432
01433
01434
01435
01436
01437
01438
01439 inline static cxint
01440 _giraffe_bias_compute_profile(GiBiasResults* results, const cpl_image* image,
01441 const cpl_matrix* areas, cxchar axis)
01442 {
01443
01444 const cxchar* const fctid = "_giraffe_bias_compute_profile";
01445
01446
01447 cxint nx = 0;
01448 cxint ny = 0;
01449 cxint sx = 0;
01450 cxint sy = 0;
01451
01452 cxsize j = 0;
01453 cxsize npixel = 0;
01454 cxsize ntotal = 0;
01455 cxsize nareas = 0;
01456
01457 cxdouble rms = 0;
01458
01459 cpl_matrix* _areas = NULL;
01460
01461 cpl_image* profile = NULL;
01462 cpl_image* model = NULL;
01463
01464
01465 cx_assert(results != NULL);
01466 cx_assert(results->limits != NULL);
01467 cx_assert(results->model == NULL);
01468
01469 cx_assert(areas != NULL);
01470
01471 cx_assert(image != NULL);
01472 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
01473
01474 cx_assert((axis == 'x') || (axis == 'y'));
01475
01476
01477 nx = cpl_image_get_size_x(image);
01478 ny = cpl_image_get_size_y(image);
01479
01480
01481
01482
01483
01484
01485 _areas = cpl_matrix_duplicate(areas);
01486 cx_assert(_areas != NULL);
01487
01488 nareas = cpl_matrix_get_nrow(_areas);
01489
01490 if (axis == 'y') {
01491
01492 for (j = 0; j < nareas; ++j) {
01493
01494 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
01495 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
01496 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
01497 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
01498
01499
01500 if (x_0 < 0) {
01501 x_0 = 0;
01502 cpl_matrix_set(_areas, j, 0, x_0);
01503 }
01504
01505 if (x_1 >= nx) {
01506 x_1 = nx - 1;
01507 cpl_matrix_set(_areas, j, 1, x_1);
01508 }
01509
01510 if (y_0 != 0) {
01511 y_0 = 0;
01512 cpl_matrix_set(_areas, j, 2, y_0);
01513 }
01514
01515 if (y_1 != ny - 1) {
01516 y_1 = ny - 1;
01517 cpl_matrix_set(_areas, j, 3, y_1);
01518 }
01519
01520 if ((x_1 >= x_0) && (y_1 >= y_0)) {
01521 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
01522 }
01523
01524 }
01525
01526 sx = nareas;
01527 sy = ny;
01528 }
01529 else {
01530
01531 for (j = 0; j < nareas; ++j) {
01532
01533 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
01534 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
01535 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
01536 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
01537
01538
01539 if (x_0 != 0) {
01540 x_0 = 0;
01541 cpl_matrix_set(_areas, j, 0, x_0);
01542 }
01543
01544 if (x_1 != nx) {
01545 x_1 = nx - 1;
01546 cpl_matrix_set(_areas, j, 1, x_1);
01547 }
01548
01549 if (y_0 < 0) {
01550 y_0 = 0;
01551 cpl_matrix_set(_areas, j, 2, y_0);
01552 }
01553
01554 if (y_1 >= ny) {
01555 y_1 = ny - 1;
01556 cpl_matrix_set(_areas, j, 3, y_1);
01557 }
01558
01559 if ((x_1 >= x_0) && (y_1 >= y_0)) {
01560 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
01561 }
01562
01563 }
01564
01565 sx = nareas;
01566 sy = nx;
01567
01568 }
01569
01570 cpl_msg_info(fctid, "Bias Areas: Found %zd points in %dx%d image",
01571 ntotal, nx, ny);
01572
01573
01574
01575
01576
01577
01578 results->mean = 0.;
01579 results->sigma = 0.;
01580
01581 cx_string_set(results->limits, "");
01582
01583 profile = cpl_image_new(sx, sy, CPL_TYPE_DOUBLE);
01584
01585 for (j = 0; j < nareas; ++j) {
01586
01587 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
01588 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
01589 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
01590 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
01591
01592 if ((x_1 >= x_0) && (y_1 >= y_0)) {
01593
01594 cx_string* tmp = cx_string_new();
01595
01596
01597 cx_string_sprintf(tmp, "%d:%d:%d:%d;", x_0, x_1, y_0, y_1);
01598 cx_string_append(results->limits, cx_string_get(tmp));
01599
01600 cx_string_delete(tmp);
01601 tmp = NULL;
01602
01603
01604
01605
01606
01607
01608 if (axis == 'y') {
01609
01610 cxint i = 0;
01611 cxint sz = x_1 - x_0 + 1;
01612
01613 const cxdouble* _img = cpl_image_get_data_double_const(image);
01614
01615 cxdouble* _profile = cpl_image_get_data_double(profile);
01616
01617
01618 for (i = y_0; i < y_1 + 1; ++i) {
01619
01620 register cxint k = i * nx;
01621 register cxint l = 0;
01622
01623 cxdouble m = giraffe_array_mean(&_img[k + x_0], sz);
01624
01625 _profile[i * sx + j] = m;
01626
01627 for (l = x_0; l < x_1 + 1; ++l) {
01628 rms += (_img[k + l] - m) * (_img[k + l] - m);
01629 ++npixel;
01630 }
01631
01632 }
01633
01634 }
01635 else {
01636
01637 cxint i = 0;
01638 cxint sz = y_1 - y_0 + 1;
01639
01640 const cxdouble* _img = cpl_image_get_data_double_const(image);
01641
01642 cxdouble* _profile = cpl_image_get_data_double(profile);
01643 cxdouble* buffer = cx_calloc(sz, sizeof(cxdouble));
01644
01645
01646 for (i = x_0; i < x_1 + 1; ++i) {
01647
01648 register cxint l = 0;
01649
01650 register cxdouble m = 0.;
01651
01652
01653 for (l = y_0; l < y_1 + 1; ++l) {
01654 buffer[l - y_0] = _img[l * nx + i];
01655 }
01656
01657 m = giraffe_array_mean(buffer, sz);
01658 _profile[i * sx + j] = m;
01659
01660 for (l = 0; l < sz; ++l) {
01661 rms += (buffer[l] - m) * (buffer[l] - m);
01662 ++npixel;
01663 }
01664
01665 }
01666
01667 cx_free(buffer);
01668 buffer = NULL;
01669
01670 }
01671
01672 }
01673
01674 }
01675
01676 cpl_matrix_delete(_areas);
01677 _areas = NULL;
01678
01679
01680
01681
01682
01683
01684 model = cpl_image_collapse_create(profile, 1);
01685
01686 cpl_image_delete(profile);
01687 profile = NULL;
01688
01689 cpl_image_divide_scalar(model, nareas);
01690 results->model = model;
01691
01692 model = NULL;
01693
01694 results->mean = cpl_image_get_mean(results->model);
01695
01696 if ( npixel > 1) {
01697 results->sigma = sqrt(rms / (npixel - 1));
01698 }
01699
01700 return EXIT_SUCCESS;
01701
01702 }
01703
01704
01705
01706
01707
01708
01709
01710
01711
01712
01713
01714
01715
01716
01717
01718
01719
01720
01721
01722
01723
01724
01725
01726
01727
01728
01729
01730
01731
01732
01733
01734 inline static cxint
01735 _giraffe_bias_compute(GiBiasResults* results, cpl_image* img,
01736 const cpl_matrix* biasareas,
01737 const cpl_image* master_bias,
01738 const cpl_image* bad_pixel_map,
01739 GiBiasMethod method, GiBiasOption option,
01740 GiBiasModel model, cxbool remove_bias, cxdouble mbias,
01741 cxdouble xdeg, cxdouble ydeg, cxdouble xstep,
01742 cxdouble ystep, cxdouble sigma, cxint numiter,
01743 cxdouble maxfraction)
01744 {
01745
01746 const cxchar* const fctid = "giraffe_bias_compute";
01747
01748 cxint i;
01749 cxint j;
01750 cxint img_dimx;
01751 cxint img_dimy;
01752 cxint img_centerx;
01753 cxint img_centery;
01754 cxint master_bias_dimx;
01755 cxint master_bias_dimy;
01756 cxint bad_pixel_dimx;
01757 cxint bad_pixel_dimy;
01758
01759 cxint error_code = 0;
01760
01761 const cpl_matrix* coeff = NULL;
01762
01763 cpl_matrix* tmp_x;
01764 cpl_matrix* tmp_y;
01765 cpl_matrix* tmp_fit;
01766
01767 cx_assert(results != NULL);
01768 cx_assert(results->limits == NULL);
01769 cx_assert(results->coeffs == NULL);
01770 cx_assert(results->model == NULL);
01771
01772 cx_assert(img != NULL);
01773 cx_assert(cpl_image_get_type(img) == CPL_TYPE_DOUBLE);
01774
01775 cx_assert(biasareas != NULL);
01776
01777
01778
01779
01780
01781
01782 img_dimx = cpl_image_get_size_x(img);
01783 img_dimy = cpl_image_get_size_y(img);
01784
01785 img_centerx = img_dimx / 2;
01786 img_centery = img_dimy / 2;
01787
01788 results->limits = cx_string_new();
01789
01790
01791
01792
01793
01794
01795 if (remove_bias) {
01796 cpl_msg_info(fctid, "removing bias value");
01797 }
01798 else {
01799 cpl_msg_info(fctid, "NOT removing bias value");
01800 }
01801
01802
01803 if (model == GIBIAS_MODEL_MEAN) {
01804
01805
01806
01807
01808
01809 cpl_msg_info(fctid, "using MEAN Model...");
01810
01811 if ((option == GIBIAS_OPTION_PLANE) ||
01812 (method == GIBIAS_METHOD_PLANE)) {
01813
01814 cpl_msg_error(fctid, "Can not use MEAN model with PLANE method.");
01815 return -3;
01816
01817 }
01818
01819 error_code = _giraffe_bias_compute_mean(results, img, biasareas,
01820 sigma, numiter, maxfraction);
01821
01822
01823 }
01824 else if ((method == GIBIAS_METHOD_CURVE) ||
01825 (option == GIBIAS_OPTION_CURVE)) {
01826
01827
01828
01829
01830
01831
01832 cpl_msg_info(fctid, "using CURVE Model...");
01833
01834 error_code = _giraffe_bias_compute_curve(results, img, biasareas,
01835 sigma, numiter, maxfraction, xdeg, ydeg, xstep, ystep);
01836
01837 if (error_code != EXIT_SUCCESS) {
01838
01839 cpl_msg_error(fctid, "Error during calculation of bias curve, "
01840 "aborting...");
01841 return error_code;
01842
01843 }
01844
01845 coeff = results->coeffs;
01846
01847 }
01848 else {
01849
01850
01851
01852
01853
01854 cpl_msg_info(fctid, "using PLANE (FITTED) model\n");
01855
01856 error_code = _giraffe_bias_compute_plane(results, img, biasareas,
01857 sigma, numiter, maxfraction);
01858
01859 if (error_code == EXIT_SUCCESS) {
01860
01861 coeff = results->coeffs;
01862
01863 results->mean = cpl_matrix_get(coeff, 0, 0) +
01864 cpl_matrix_get(coeff, 0, 1) * img_centerx +
01865 cpl_matrix_get(coeff, 0, 2) * img_centery;
01866
01867 }
01868 else {
01869
01870 cpl_msg_error(fctid, "Error during calculation of bias plane, "
01871 "aborting...");
01872 return error_code;
01873
01874 }
01875
01876 }
01877
01878
01879
01880
01881
01882
01883 switch (method) {
01884 case GIBIAS_METHOD_UNIFORM:
01885 {
01886
01887 cpl_msg_info(fctid, "using bias method == UNIFORM");
01888
01889
01890
01891
01892
01893
01894 if (model == GIBIAS_MODEL_MEAN) {
01895
01896 cpl_msg_info(fctid, "bias value(areas mean) = %f, "
01897 "bias sigma = %f", results->mean, results->sigma);
01898
01899 }
01900 else {
01901
01902 cpl_msg_info(fctid, "bias value at (%d, %d) = %f, "
01903 "bias sigma = %f", img_centerx, img_centery,
01904 results->mean, results->sigma);
01905
01906 }
01907
01908 if (remove_bias == TRUE) {
01909
01910 cpl_image_subtract_scalar(img, results->mean);
01911
01912 }
01913
01914 break;
01915 }
01916
01917 case GIBIAS_METHOD_PLANE:
01918 {
01919
01920 cpl_msg_info(fctid, "using bias method == PLANE");
01921
01922 if (coeff == NULL) {
01923
01924 error_code = _giraffe_bias_compute_plane(results, img,
01925 biasareas, sigma, numiter, maxfraction);
01926
01927 if (error_code == EXIT_SUCCESS) {
01928
01929 coeff = results->coeffs;
01930
01931 results->mean = cpl_matrix_get(coeff, 0, 0) +
01932 cpl_matrix_get(coeff, 0, 1) * img_centerx +
01933 cpl_matrix_get(coeff, 0, 2) * img_centery;
01934
01935 }
01936 else {
01937
01938 cpl_msg_error(fctid, "Error during calculation of bias "
01939 "plane, aborting...");
01940 return error_code;
01941
01942 }
01943
01944 }
01945
01946 cpl_msg_info(fctid, "bias value at (%d, %d) = %f, bias sigma = %f, "
01947 "bias plane = %g + %g * x + %g * y", img_centerx,
01948 img_centery, results->mean, results->sigma,
01949 cpl_matrix_get(coeff, 0, 0),
01950 cpl_matrix_get(coeff, 0, 1),
01951 cpl_matrix_get(coeff, 0, 2));
01952
01953 if (remove_bias == TRUE) {
01954
01955
01956
01957 cxdouble* _img = cpl_image_get_data_double(img);
01958
01959 #if 0
01960 cpl_image* bsmodel = cpl_image_new(img_dimx, img_dimy,
01961 CPL_TYPE_DOUBLE);
01962 cxdouble* _bsmodel = cpl_image_get_data_double(bsmodel);
01963 #endif
01964
01965
01966 for (j = 0; j < img_dimy; j++) {
01967
01968 cxsize jj = j * img_dimx;
01969
01970 for (i = 0; i < img_dimx; i++) {
01971
01972 #if 0
01973 _bsmodel[jj + i] = cpl_matrix_get(coeff, 0, 0)
01974 + cpl_matrix_get(coeff, 0, 1) * i
01975 + cpl_matrix_get(coeff, 0, 2) * j;
01976 #endif
01977
01978 _img[jj + i] -= cpl_matrix_get(coeff, 0, 0)
01979 + cpl_matrix_get(coeff, 0, 1) * i
01980 + cpl_matrix_get(coeff, 0, 2) * j;
01981
01982 }
01983
01984 }
01985
01986 #if 0
01987 cpl_image_save(bsmodel, "idiot.fits", CPL_BPP_IEEE_FLOAT,
01988 0, CPL_IO_DEFAULT);
01989 cpl_image_delete(bsmodel);
01990 #endif
01991
01992 }
01993
01994 break;
01995
01996 }
01997
01998 case GIBIAS_METHOD_CURVE:
01999 {
02000
02001 cpl_msg_info(fctid, "using bias method == CURVE");
02002
02003 if (coeff == NULL) {
02004
02005 error_code = _giraffe_bias_compute_curve(results, img,
02006 biasareas, sigma, numiter, maxfraction, xdeg, ydeg,
02007 xstep, ystep);
02008
02009 if (error_code != EXIT_SUCCESS) {
02010
02011 cpl_msg_error(fctid, "Error during calculation of bias "
02012 "surface curve, aborting...");
02013 return error_code;
02014
02015 }
02016
02017 coeff = results->coeffs;
02018
02019 }
02020
02021 cpl_msg_info(fctid, "bias value %f, bias sigma = %f\n",
02022 results->mean, results->sigma);
02023
02024
02025
02026
02027 if (remove_bias == TRUE) {
02028
02029 cxdouble* _img = cpl_image_get_data_double(img);
02030
02031 #if 0
02032 cpl_image* bsmodel = cpl_image_new(img_dimx, img_dimy,
02033 CPL_TYPE_DOUBLE);
02034 cxdouble* _bsmodel = cpl_image_get_data_double(bsmodel);
02035 #endif
02036
02037
02038 tmp_x = cpl_matrix_new(1, 1);
02039 tmp_y = cpl_matrix_new(1, 1);
02040
02041 for (j = 0; j < img_dimy; j++) {
02042
02043 cxsize jj = j * img_dimx;
02044
02045 cpl_matrix_set(tmp_y, 0, 0, j);
02046
02047 for (i = 0; i < img_dimx; i++) {
02048
02049 cpl_matrix_set(tmp_x, 0, 0, i);
02050 tmp_fit = giraffe_chebyshev_fit2d(0., 0., img_dimx,
02051 img_dimy, coeff, tmp_x, tmp_y);
02052
02053 #if 0
02054 _bsmodel[jj + i] = cpl_matrix_get(tmp_fit, 0, 0);
02055 #endif
02056 _img[jj + i] -= cpl_matrix_get(tmp_fit, 0, 0);
02057
02058 cpl_matrix_delete(tmp_fit);
02059
02060 }
02061
02062 }
02063
02064 cpl_matrix_delete(tmp_x);
02065 cpl_matrix_delete(tmp_y);
02066
02067 #if 0
02068 cpl_image_save(bsmodel, "idiot.fits", CPL_BPP_IEEE_FLOAT,
02069 0, CPL_IO_DEFAULT);
02070 cpl_image_delete(bsmodel);
02071 #endif
02072 }
02073
02074 break;
02075
02076 }
02077
02078 case GIBIAS_METHOD_PROFILE:
02079 {
02080
02081 cpl_msg_info(fctid, "using bias method 'PROFILE'");
02082
02083 error_code = _giraffe_bias_compute_profile(results, img,
02084 biasareas, 'y');
02085
02086 if (error_code != EXIT_SUCCESS) {
02087
02088 cpl_msg_error(fctid, "Error computing the bias area(s) "
02089 "profile");
02090 return error_code;
02091
02092 }
02093
02094 if (remove_bias == TRUE) {
02095
02096 const cxdouble* _bias =
02097 cpl_image_get_data_double(results->model);
02098
02099 cxdouble* _img = cpl_image_get_data_double(img);
02100
02101
02102 cx_assert(_bias != NULL);
02103 cx_assert(_img != NULL);
02104
02105 for (j = 0; j < img_dimy; ++j) {
02106
02107 cxsize jj = j * img_dimx;
02108
02109
02110 for (i = 0; i < img_dimx; ++i) {
02111 _img[jj + i] -= _bias[j];
02112 }
02113
02114 }
02115
02116 }
02117
02118 break;
02119
02120 }
02121
02122 case GIBIAS_METHOD_MASTER:
02123 case GIBIAS_METHOD_ZMASTER:
02124 {
02125
02126 cxdouble biasdrift = 0.;
02127
02128 if (master_bias == NULL) {
02129
02130 cpl_msg_error(fctid, "Master Bias Frame required with MASTER "
02131 "method.");
02132 return -5;
02133
02134 }
02135
02136 master_bias_dimx = cpl_image_get_size_x(master_bias);
02137 master_bias_dimy = cpl_image_get_size_y(master_bias);
02138
02139 if ((master_bias_dimx != img_dimx) ||
02140 (master_bias_dimy != img_dimy)) {
02141
02142 cpl_msg_error( fctid, "Invalid Masterbias size should "
02143 "be [%d, %d] but is [%d, %d].", img_dimx, img_dimy,
02144 master_bias_dimx, master_bias_dimy);
02145 return -7;
02146
02147 }
02148
02149 if (coeff == NULL) {
02150
02151 if (option == GIBIAS_OPTION_CURVE) {
02152
02153 error_code = _giraffe_bias_compute_curve(results, img,
02154 biasareas, sigma, numiter, maxfraction, xdeg,
02155 ydeg, xstep, ystep);
02156
02157 if (error_code != EXIT_SUCCESS) {
02158
02159 cpl_msg_error(fctid, "Error during calculation of "
02160 "bias surface curve, aborting...");
02161 return error_code;
02162
02163 }
02164
02165 coeff = results->coeffs;
02166
02167 }
02168 else {
02169
02170
02171
02172
02173
02174 error_code = _giraffe_bias_compute_plane(results, img,
02175 biasareas, sigma, numiter, maxfraction);
02176
02177 if (error_code == EXIT_SUCCESS) {
02178
02179 coeff = results->coeffs;
02180
02181 results->mean = cpl_matrix_get(coeff, 0, 0) +
02182 cpl_matrix_get(coeff, 0, 1) * img_centerx +
02183 cpl_matrix_get(coeff, 0, 2) * img_centery;
02184
02185 }
02186 else {
02187
02188 cpl_msg_error(fctid, "Error during calculation of "
02189 "bias surface plane, aborting...");
02190 return error_code;
02191
02192 }
02193
02194 }
02195
02196 }
02197
02198 if (method == GIBIAS_METHOD_ZMASTER) {
02199
02200
02201
02202
02203
02204
02205
02206 biasdrift = results->mean - mbias;
02207
02208 cpl_msg_info(fctid, "METHOD == ZMASTER, drift value %g",
02209 biasdrift);
02210
02211 }
02212 else {
02213
02214 cpl_msg_info(fctid, "METHOD == MASTER");
02215
02216 }
02217
02218
02219
02220
02221
02222
02223 if (option == GIBIAS_OPTION_CURVE) {
02224
02225
02226
02227 cpl_msg_info(fctid, "OPTION == CURVE, bias mean = %f, bias "
02228 "sigma = %f", results->mean, results->sigma);
02229
02230 }
02231 else {
02232
02233 if (option == GIBIAS_OPTION_PLANE) {
02234
02235 cpl_msg_info(fctid, "OPTION == PLANE, bias plane = %g + "
02236 "%g * x + %g * y, bias mean at (%d, %d) = %f, "
02237 "bias sigma = %f\n",
02238 cpl_matrix_get(coeff, 0, 0),
02239 cpl_matrix_get(coeff, 0, 1),
02240 cpl_matrix_get(coeff, 0, 2),
02241 img_centerx, img_centery, results->mean,
02242 results->sigma);
02243
02244 }
02245 else {
02246
02247 cpl_msg_info(fctid, "bias mean = %f, bias sigma = %f\n",
02248 results->mean, results->sigma);
02249
02250 }
02251
02252 }
02253
02254
02255
02256
02257
02258
02259 if (bad_pixel_map == NULL) {
02260
02261
02262
02263
02264
02265
02266
02267 if (remove_bias == TRUE) {
02268
02269 const cxdouble* _mbias =
02270 cpl_image_get_data_double_const(master_bias);
02271
02272 cxdouble* _img = cpl_image_get_data_double(img);
02273
02274
02275 for (j = 0; j < img_dimy; j++) {
02276
02277 cxsize jj = j * img_dimx;
02278
02279 for (i = 0; i < img_dimx; i++) {
02280 _img[jj + i] -= (_mbias[jj + i] + biasdrift);
02281 }
02282
02283 }
02284
02285 }
02286
02287 }
02288 else {
02289
02290
02291
02292
02293
02294
02295
02296
02297
02298 bad_pixel_dimx = cpl_image_get_size_x(bad_pixel_map);
02299 bad_pixel_dimy = cpl_image_get_size_y(bad_pixel_map);
02300
02301 if ((bad_pixel_dimx != img_dimx) ||
02302 (bad_pixel_dimy != img_dimy)) {
02303
02304 cpl_msg_error(fctid, "Invalid bad pixel map size should "
02305 "be [%d, %d] but is [%d, %d].", img_dimx,
02306 img_dimy, bad_pixel_dimx, bad_pixel_dimy);
02307 return -6;
02308
02309 }
02310
02311 if (remove_bias == TRUE) {
02312
02313 if (option == GIBIAS_OPTION_CURVE) {
02314
02315 const cxint* _bpix =
02316 cpl_image_get_data_int_const(bad_pixel_map);
02317
02318 const cxdouble* _mbias =
02319 cpl_image_get_data_double_const(master_bias);
02320
02321 cxdouble* _img = cpl_image_get_data_double(img);
02322
02323
02324 tmp_x = cpl_matrix_new(1, 1);
02325 tmp_y = cpl_matrix_new(1, 1);
02326
02327 for (j = 0; j < img_dimy; j++) {
02328
02329 cxsize jj = j * img_dimx;
02330
02331 cpl_matrix_set(tmp_y, 0, 0, j);
02332
02333 for (i = 0; i < img_dimx; i++) {
02334
02335 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
02336
02337 _img[jj + i] -= (_mbias[jj + i] +
02338 biasdrift);
02339
02340 }
02341 else {
02342
02343 cpl_matrix_set(tmp_x, 0, 0, i);
02344 tmp_fit = giraffe_chebyshev_fit2d(0., 0.,
02345 img_dimx, img_dimy, coeff,
02346 tmp_x, tmp_y);
02347
02348 _img[jj + i] -= cpl_matrix_get(tmp_fit,
02349 0, 0);
02350
02351 cpl_matrix_delete(tmp_fit);
02352
02353 }
02354
02355 }
02356
02357 }
02358
02359 cpl_matrix_delete(tmp_y);
02360 cpl_matrix_delete(tmp_x);
02361
02362 }
02363 else if (option == GIBIAS_OPTION_PLANE) {
02364
02365 const cxint* _bpix =
02366 cpl_image_get_data_int_const(bad_pixel_map);
02367
02368 const cxdouble* _mbias =
02369 cpl_image_get_data_double_const(master_bias);
02370
02371 cxdouble* _img = cpl_image_get_data_double(img);
02372
02373
02374 for (j = 0; j < img_dimy; j++) {
02375
02376 cxsize jj = j * img_dimx;
02377
02378 for (i = 0; i < img_dimx; i++) {
02379
02380 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
02381
02382 _img[jj + i] -= (_mbias[jj + i] +
02383 biasdrift);
02384
02385 }
02386 else {
02387
02388 _img[jj + i] -=
02389 cpl_matrix_get(coeff, 0, 0) +
02390 cpl_matrix_get(coeff, 0, 1) * j +
02391 cpl_matrix_get(coeff, 0, 2) * i;
02392 }
02393
02394 }
02395
02396 }
02397
02398 }
02399 else {
02400
02401 const cxint* _bpix =
02402 cpl_image_get_data_int_const(bad_pixel_map);
02403
02404 const cxdouble* _mbias =
02405 cpl_image_get_data_double_const(master_bias);
02406
02407 cxdouble* _img = cpl_image_get_data_double(img);
02408
02409
02410 cpl_msg_info(fctid, "OPTION==UNIFORM\n");
02411
02412 for (j = 0; j < img_dimy; j++) {
02413
02414 cxsize jj = j * img_dimx;
02415
02416 for (i = 0; i < img_dimx; i++) {
02417
02418 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
02419
02420 _img[jj + i] -= (_mbias[jj + i] +
02421 biasdrift);
02422
02423 }
02424 else {
02425
02426 _img[jj + i] -= results->mean;
02427
02428 }
02429
02430 }
02431
02432 }
02433
02434 }
02435
02436 }
02437
02438 }
02439
02440 break;
02441
02442 }
02443
02444 default:
02445 {
02446
02447 cpl_msg_error(fctid, "Invalid method, aborting...");
02448 return -4;
02449 break;
02450
02451 }
02452
02453 }
02454
02455
02456 cpl_msg_info(fctid, "Resulting biaslimits : %s",
02457 cx_string_get(results->limits));
02458
02459 return 0;
02460
02461 }
02462
02463
02481 cpl_matrix *
02482 giraffe_get_raw_areas(const GiImage *image)
02483 {
02484
02485 const cxchar *const _id = "giraffe_get_raw_areas";
02486
02487 cxint32 prescx = 0;
02488 cxint32 prescy = 0;
02489 cxint32 ovrscx = 0;
02490 cxint32 ovrscy = 0;
02491
02492 cxint32 winx = 0;
02493 cxint32 winy = 0;
02494
02495 cxint32 tprescx = 0;
02496 cxint32 tprescy = 0;
02497 cxint32 tovrscx = 0;
02498 cxint32 tovrscy = 0;
02499
02500 cxint32 row = 0;
02501
02502 cpl_matrix *m_tmp;
02503 cpl_propertylist *properties;
02504
02505
02506 properties = giraffe_image_get_properties(image);
02507 if (!properties) {
02508 cpl_error_set(_id, CPL_ERROR_NULL_INPUT);
02509 return NULL;
02510 }
02511
02512 winx = cpl_propertylist_get_int(properties, GIALIAS_WINX);
02513 winy = cpl_propertylist_get_int(properties, GIALIAS_WINY);
02514
02515 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
02516 tprescx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
02517 if (tprescx > 0) {
02518 prescx = tprescx;
02519 }
02520 }
02521 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
02522 tprescy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
02523 if (tprescy > 0) {
02524 prescy = tprescy;
02525 }
02526 }
02527 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
02528 tovrscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
02529 if (tovrscx > 0) {
02530 ovrscx = tovrscx;
02531 }
02532 }
02533 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
02534 tovrscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
02535 if (tovrscy > 0) {
02536 ovrscy = tovrscy;
02537 }
02538 }
02539
02540
02541 m_tmp = cpl_matrix_new(1, 4);
02542
02543 if (prescx > 0) {
02544
02545 cpl_matrix_set(m_tmp, row, 0, 0.0);
02546 cpl_matrix_set(m_tmp, row, 1, (cxdouble) prescx - 1);
02547 cpl_matrix_set(m_tmp, row, 2, 0.0);
02548 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02549
02550 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02551 row++;
02552 }
02553
02554 if (ovrscx > 0) {
02555
02556 cpl_matrix_set(m_tmp, row, 0, (cxdouble) winx - ovrscx);
02557 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
02558 cpl_matrix_set(m_tmp, row, 2, 0.0);
02559 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02560
02561 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02562 row++;
02563 }
02564
02565 if (!(prescy == 0 || prescx == 0 || ovrscx == 0)) {
02566
02567 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
02568 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
02569 cpl_matrix_set(m_tmp, row, 2, 0.0);
02570 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
02571
02572 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02573 row++;
02574
02575 }
02576 else if (!(prescy == 0 || prescx == 0)) {
02577
02578 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
02579 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
02580 cpl_matrix_set(m_tmp, row, 2, 0.0);
02581 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
02582
02583 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02584 row++;
02585
02586 }
02587 else if (!(prescy == 0 || ovrscx == 0)) {
02588
02589 cpl_matrix_set(m_tmp, row, 0, 0.0);
02590 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
02591 cpl_matrix_set(m_tmp, row, 2, 0.0);
02592 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
02593
02594 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02595 row++;
02596
02597 }
02598 else if (prescy > 0) {
02599
02600 cpl_matrix_set(m_tmp, row, 0, 0.0);
02601 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
02602 cpl_matrix_set(m_tmp, row, 2, 0.0);
02603 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
02604
02605 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02606 row++;
02607 }
02608
02609 if (!(ovrscy == 0 || prescx == 0 || ovrscx == 0)) {
02610
02611 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
02612 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
02613 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
02614 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02615
02616 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02617 row++;
02618
02619 }
02620 else if (!(ovrscy == 0 || prescx == 0)) {
02621
02622 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
02623 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
02624 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
02625 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02626
02627 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02628 row++;
02629
02630 }
02631 else if (!(ovrscy == 0 || ovrscx == 0)) {
02632
02633 cpl_matrix_set(m_tmp, row, 0, 0.0);
02634 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
02635 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
02636 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02637
02638 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02639 row++;
02640
02641 }
02642 else if (ovrscy > 0) {
02643
02644 cpl_matrix_set(m_tmp, row, 0, 0.0);
02645 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
02646 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
02647 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02648
02649 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02650 row++;
02651 }
02652
02653 cpl_matrix_resize(m_tmp, 0, -1, 0, 0);
02654 row--;
02655
02656 if (row == 0) {
02657 cpl_matrix_delete(m_tmp);
02658 return NULL;
02659 }
02660
02661 return m_tmp;
02662 }
02663
02664
02683 cxint
02684 giraffe_trim_raw_areas(GiImage *image)
02685 {
02686
02687 const cxchar *fctid = "giraffe_trim_raw_areas";
02688
02689 cxint nx, ny;
02690 cxint newlx, newly, newux, newuy;
02691
02692 cxint ovscx = 0;
02693 cxint ovscy = 0;
02694 cxint prscx = 0;
02695 cxint prscy = 0;
02696
02697 cpl_propertylist *properties = giraffe_image_get_properties(image);
02698
02699 cpl_image *_image = giraffe_image_get(image);
02700 cpl_image *tmp;
02701
02702
02703 if (!properties) {
02704 cpl_msg_error(fctid, "Image does not contain any properties!");
02705 return 1;
02706 }
02707
02708 nx = cpl_image_get_size_x(_image);
02709 ny = cpl_image_get_size_y(_image);
02710
02711 if (cpl_propertylist_has(properties, GIALIAS_NAXIS1)) {
02712 cxint _nx = cpl_propertylist_get_int(properties, GIALIAS_NAXIS1);
02713
02714 if (_nx != nx) {
02715 cpl_msg_warning(fctid, "Image size (%d) and image property `%s' "
02716 "(%d) are not consistent! Using image size ...",
02717 nx, GIALIAS_NAXIS1, _nx);
02718 }
02719 }
02720 else {
02721 cpl_msg_warning(fctid, "Image does not contain any property `%s'. "
02722 "Using image size (%d)", GIALIAS_NAXIS1, nx);
02723 }
02724
02725
02726 if (cpl_propertylist_has(properties, GIALIAS_NAXIS2)) {
02727 cxint _ny = cpl_propertylist_get_int(properties, GIALIAS_NAXIS2);
02728
02729 if (_ny != ny) {
02730 cpl_msg_warning(fctid, "Image size (%d) and image property `%s' "
02731 "(%d) are not consistent! Using image size ...",
02732 ny, GIALIAS_NAXIS2, _ny);
02733 }
02734 }
02735 else {
02736 cpl_msg_warning(fctid, "Image does not contain any property `%s'. "
02737 "Using image size (%d)", GIALIAS_NAXIS2, ny);
02738 }
02739
02740 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
02741 ovscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
02742 }
02743
02744 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
02745 ovscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
02746 }
02747
02748 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
02749 prscx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
02750 }
02751
02752 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
02753 prscy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
02754 }
02755
02756
02757
02758
02759
02760 newlx = prscx + 1;
02761 newly = prscy + 1;
02762 newux = nx - ovscx;
02763 newuy = ny - ovscy;
02764
02765 tmp = cpl_image_extract(_image, newlx, newly, newux, newuy);
02766
02767 giraffe_image_set(image, tmp);
02768 cpl_image_delete(tmp);
02769
02770 _image = giraffe_image_get(image);
02771
02772 nx = cpl_image_get_size_x(_image);
02773 ny = cpl_image_get_size_y(_image);
02774
02775
02776
02777
02778
02779
02780 cpl_propertylist_set_int(properties, GIALIAS_NAXIS1, nx);
02781 cpl_propertylist_set_int(properties, GIALIAS_NAXIS2, ny);
02782
02783 if (cpl_propertylist_has(properties, GIALIAS_CRPIX1)) {
02784 cxdouble crpix1 = cpl_propertylist_get_double(properties,
02785 GIALIAS_CRPIX1);
02786
02787
02788
02789
02790
02791
02792 crpix1 += (cxdouble)prscx;
02793 cpl_propertylist_set_double(properties, GIALIAS_CRPIX1, crpix1);
02794 }
02795
02796 if (cpl_propertylist_has(properties, GIALIAS_CRPIX2)) {
02797 cxdouble crpix2 = cpl_propertylist_get_double(properties,
02798 GIALIAS_CRPIX2);
02799
02800 crpix2 -= (cxdouble) prscy;
02801 cpl_propertylist_set_double(properties, GIALIAS_CRPIX2, crpix2);
02802 }
02803
02804 cpl_propertylist_erase(properties, GIALIAS_OVSCX);
02805 cpl_propertylist_erase(properties, GIALIAS_OVSCY);
02806 cpl_propertylist_erase(properties, GIALIAS_PRSCX);
02807 cpl_propertylist_erase(properties, GIALIAS_PRSCY);
02808
02809 return 0;
02810
02811 }
02812
02813
02875 cxint
02876 giraffe_bias_remove(GiImage* result, const GiImage* raw,
02877 const GiImage* master_bias, const GiImage* bad_pixels,
02878 const cpl_matrix* biaslimits, const GiBiasConfig* config)
02879 {
02880
02881 const cxchar* const fctid = "giraffe_bias_remove";
02882
02883 cxint error_code;
02884
02885 cxdouble bias_value = 0.;
02886
02887 cx_string* s;
02888
02889 cpl_matrix* areas = NULL;
02890
02891 cpl_image* _raw = giraffe_image_get(raw);
02892 cpl_image* _master_bias = giraffe_image_get(master_bias);
02893 cpl_image* _bad_pixels = giraffe_image_get(bad_pixels);
02894 cpl_image* timage;
02895
02896 cpl_propertylist* properties;
02897
02898 GiBiasResults bias_results = {0., 0., NULL, NULL, NULL};
02899
02900
02901 cx_assert(raw != NULL);
02902 cx_assert(config != NULL);
02903 cx_assert(biaslimits == NULL);
02904
02905 if (result == NULL) {
02906 return 1;
02907 }
02908
02909
02910
02911
02912
02913
02914 if (strncmp(config->areas, "None", 4) == 0) {
02915
02916 cpl_msg_info(fctid, "No bias areas specified, using pre/overscan"
02917 "regions of the raw frame.");
02918
02919 cpl_error_reset();
02920 areas = giraffe_get_raw_areas(raw);
02921 if (areas == NULL) {
02922 if (cpl_error_get_code() == CPL_ERROR_NULL_INPUT) {
02923 cpl_msg_error(fctid, "Invalid raw image! Image does not "
02924 "contain any properties!");
02925 }
02926 else {
02927 cpl_msg_error(fctid, "Invalid raw image! Image does not "
02928 "contain or has invalid pre- and overscan "
02929 "properties.");
02930 }
02931
02932 return 1;
02933 }
02934
02935 }
02936 else {
02937
02938 areas = _giraffe_bias_get_areas(config->areas);
02939
02940 if (areas == NULL) {
02941
02942 cpl_msg_error(fctid, "Invalid bias area specification '%s'!",
02943 config->areas);
02944 return 1;
02945
02946 }
02947
02948 cpl_msg_info(fctid, "Using bias area(s) '%s' for bias computation",
02949 config->areas);
02950
02951 }
02952
02953
02954
02955
02956
02957
02958
02959
02960
02961
02962
02963
02964 if (master_bias != NULL) {
02965 cpl_propertylist *_properties;
02966 cxbool is_same_size = FALSE;
02967
02968 is_same_size = _giraffe_compare_overscans(raw, master_bias);
02969
02970
02971
02972
02973
02974
02975 if (is_same_size==FALSE) {
02976 cpl_msg_info(fctid, "PRE/OVRSCAN Regions do not match between "
02977 "RAW Image and Masterbias Image.");
02978
02979 return 1;
02980 }
02981
02982 _properties = giraffe_image_get_properties(master_bias);
02983
02984 if (cpl_propertylist_has(_properties, GIALIAS_BIASVALUE)) {
02985 bias_value = cpl_propertylist_get_double(_properties,
02986 GIALIAS_BIASVALUE);
02987 }
02988 }
02989
02990
02991
02992
02993
02994
02995
02996 if (bad_pixels != NULL) {
02997 cxbool is_same_size = FALSE;
02998
02999 is_same_size = _giraffe_compare_overscans(raw, bad_pixels);
03000
03001
03002
03003
03004
03005
03006 if (is_same_size == FALSE) {
03007 cpl_msg_info(fctid, "PRE/OVRSCAN Regions do not match between "
03008 "RAW Image and Bad Pixel Image.");
03009
03010 return 1;
03011 }
03012 }
03013
03014
03015
03016
03017
03018
03019
03020
03021
03022
03023 timage = cpl_image_duplicate(_raw);
03024
03025
03026
03027
03028
03029
03030
03031 error_code = _giraffe_bias_compute(&bias_results, timage,
03032 areas, _master_bias,
03033 _bad_pixels, config->method,
03034 config->option, config->model,
03035 config->remove, bias_value,
03036 config->xdeg, config->ydeg,
03037 config->xstep, config->ystep,
03038 config->sigma, config->iterations,
03039 config->fraction);
03040
03041 cpl_matrix_delete(areas);
03042 areas = NULL;
03043
03044
03045
03046
03047
03048
03049 if (error_code) {
03050
03051 _giraffe_biasresults_clear(&bias_results);
03052
03053 cpl_msg_error(fctid, "Bias computation failed with error %d!",
03054 error_code);
03055 return 1;
03056
03057 }
03058
03059
03060
03061
03062
03063
03064 properties = giraffe_image_get_properties(raw);
03065
03066 if (config->remove == TRUE) {
03067
03068 giraffe_image_set(result, timage);
03069 cpl_image_delete(timage);
03070
03071 giraffe_image_set_properties(result, properties);
03072
03073 }
03074 else {
03075
03076 cpl_image_delete(timage);
03077
03078 giraffe_image_set(result, _raw);
03079 giraffe_image_set_properties(result, properties);
03080
03081 }
03082
03083
03084
03085
03086
03087
03088 properties = giraffe_image_get_properties(result);
03089
03090 if (config->remove == TRUE) {
03091
03092 cpl_propertylist_set_int(properties, GIALIAS_BITPIX, -32);
03093 cpl_propertylist_set_double(properties, GIALIAS_BZERO, 0.);
03094 cpl_propertylist_set_double(properties, GIALIAS_BSCALE, 1.);
03095
03096 if (cpl_propertylist_has(properties, GIALIAS_GIRFTYPE)) {
03097 cpl_propertylist_set_string(properties,
03098 GIALIAS_GIRFTYPE, "BRMIMG");
03099 }
03100 else {
03101 cpl_propertylist_append_string(properties,
03102 GIALIAS_GIRFTYPE, "BRMIMG");
03103 }
03104 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
03105 "GIRAFFE bias subtracted frame");
03106
03107
03108
03109
03110
03111 giraffe_trim_raw_areas(result);
03112
03113
03114
03115
03116
03117
03118
03119
03120 #if 0
03121 if (cpl_propertylist_has(properties, GIALIAS_CONAD)) {
03122 cxdouble conad = cpl_propertylist_get_double(properties, GIALIAS_CONAD);
03123
03124 if (conad > 0.) {
03125 cpl_image* _result = giraffe_image_get(result);
03126
03127 cpl_msg_info(fctid, "Converting bias subtracted frame to "
03128 "electrons; conversion factor is %.2f e-/ADU",
03129 conad);
03130 cpl_image_const_op_local(_result, conad, '*');
03131 }
03132 else {
03133 cpl_msg_warning(fctid, "Invalid conversion factor %.2f "
03134 "e-/ADU! Bias subtracted image will not be "
03135 "converted.", conad);
03136 }
03137 }
03138 else {
03139 cpl_msg_info(fctid, "Conversion factor not found. Bias subtracted "
03140 "image will remain in ADU");
03141 }
03142 #endif
03143 }
03144
03145 s = cx_string_new();
03146
03147 _giraffe_method_string(s, config->method, config->option);
03148 cpl_propertylist_append_string(properties, GIALIAS_BIASMETHOD, cx_string_get(s));
03149
03150 cx_string_delete(s);
03151
03152 cpl_propertylist_append_double(properties, GIALIAS_BIASVALUE,
03153 bias_results.mean);
03154 cpl_propertylist_append_double(properties, GIALIAS_BIASSIGMA,
03155 bias_results.sigma);
03156
03157 if (bias_results.coeffs) {
03158 cpl_propertylist_append_double(properties, GIALIAS_BCLIPSIGMA,
03159 config->sigma);
03160 cpl_propertylist_append_int(properties, GIALIAS_BCLIPNITER,
03161 config->iterations);
03162 cpl_propertylist_append_double(properties, GIALIAS_BCLIPMFRAC,
03163 config->fraction);
03164 }
03165
03166 if (bias_results.limits) {
03167 cpl_propertylist_append_string(properties, GIALIAS_BIASAREAS,
03168 cx_string_get(bias_results.limits));
03169 }
03170
03171 if (bias_results.coeffs) {
03172 s = cx_string_new();
03173
03174 _giraffe_stringify_coefficients(s, bias_results.coeffs);
03175 cpl_propertylist_append_string(properties, GIALIAS_BIASPLANE,
03176 cx_string_get(s));
03177
03178 cx_string_delete(s);
03179 }
03180
03181
03182
03183
03184
03185
03186 _giraffe_biasresults_clear(&bias_results);
03187
03188 return 0;
03189
03190 }
03191
03192
03203 GiBiasConfig*
03204 giraffe_bias_config_create(cpl_parameterlist* list)
03205 {
03206
03207 const cxchar* s;
03208 cpl_parameter* p;
03209
03210 GiBiasConfig* config = NULL;
03211
03212
03213 if (!list) {
03214 return NULL;
03215 }
03216
03217 config = cx_calloc(1, sizeof *config);
03218
03219
03220
03221
03222
03223
03224 config->method = GIBIAS_METHOD_UNDEFINED;
03225 config->option = GIBIAS_OPTION_UNDEFINED;
03226 config->model = GIBIAS_MODEL_UNDEFINED;
03227 config->mbias = 0.;
03228 config->xdeg = 1;
03229 config->ydeg = 1;
03230
03231
03232 p = cpl_parameterlist_find(list, "giraffe.biasremoval.remove");
03233 config->remove = cpl_parameter_get_bool(p);
03234
03235 p = cpl_parameterlist_find(list, "giraffe.biasremoval.method");
03236 s = cpl_parameter_get_string(p);
03237
03238 if (!strcmp(s, "UNIFORM")) {
03239 config->method = GIBIAS_METHOD_UNIFORM;
03240 }
03241
03242 if (!strcmp(s, "PLANE")) {
03243 config->method = GIBIAS_METHOD_PLANE;
03244 }
03245
03246 if (!strcmp(s, "PROFILE")) {
03247 config->method = GIBIAS_METHOD_PROFILE;
03248 }
03249
03250 if (!strcmp(s, "MASTER")) {
03251 config->method = GIBIAS_METHOD_MASTER;
03252 }
03253
03254 if (!strcmp(s, "ZMASTER")) {
03255 config->method = GIBIAS_METHOD_ZMASTER;
03256 }
03257
03258 if (!strcmp(s, "CURVE")) {
03259 config->method = GIBIAS_METHOD_CURVE;
03260 }
03261
03262 if (!strcmp(s, "MASTER+PLANE")) {
03263 config->method = GIBIAS_METHOD_MASTER;
03264 config->option = GIBIAS_OPTION_PLANE;
03265 }
03266
03267 if (!strcmp(s, "ZMASTER+PLANE")) {
03268 config->method = GIBIAS_METHOD_ZMASTER;
03269 config->option = GIBIAS_OPTION_PLANE;
03270 }
03271
03272 if (!strcmp(s, "MASTER+CURVE")) {
03273 config->method = GIBIAS_METHOD_MASTER;
03274 config->option = GIBIAS_OPTION_CURVE;
03275 }
03276
03277 if (!strcmp(s, "ZMASTER+CURVE")) {
03278 config->method = GIBIAS_METHOD_ZMASTER;
03279 config->option = GIBIAS_OPTION_CURVE;
03280 }
03281
03282 cx_assert(config->method != GIBIAS_METHOD_UNDEFINED);
03283
03284 p = cpl_parameterlist_find(list, "giraffe.biasremoval.areas");
03285 s = cpl_parameter_get_string(p);
03286 config->areas = cx_strdup(s);
03287
03288 p = cpl_parameterlist_find(list, "giraffe.biasremoval.sigma");
03289 config->sigma = cpl_parameter_get_double(p);
03290
03291 p = cpl_parameterlist_find(list, "giraffe.biasremoval.iterations");
03292 config->iterations = cpl_parameter_get_int(p);
03293
03294 p = cpl_parameterlist_find(list, "giraffe.biasremoval.fraction");
03295 config->fraction = cpl_parameter_get_double(p);
03296
03297 if (config->method == GIBIAS_METHOD_CURVE ||
03298 config->option == GIBIAS_OPTION_CURVE) {
03299 p = cpl_parameterlist_find(list, "giraffe.biasremoval.xorder");
03300 config->xdeg = 1 + cpl_parameter_get_int(p);
03301
03302 p = cpl_parameterlist_find(list, "giraffe.biasremoval.yorder");
03303 config->ydeg = 1 + cpl_parameter_get_int(p);
03304 }
03305
03306 p = cpl_parameterlist_find(list, "giraffe.biasremoval.xstep");
03307 config->xstep = cpl_parameter_get_int(p);
03308
03309 p = cpl_parameterlist_find(list, "giraffe.biasremoval.ystep");
03310 config->ystep = cpl_parameter_get_int(p);
03311
03312 return config;
03313
03314 }
03315
03316
03329 void
03330 giraffe_bias_config_destroy(GiBiasConfig* config)
03331 {
03332
03333 if (config) {
03334 if (config->areas) {
03335 cx_free(config->areas);
03336 config->areas = NULL;
03337 }
03338
03339 cx_free(config);
03340 }
03341
03342 return;
03343 }
03344
03345
03357 void
03358 giraffe_bias_config_add(cpl_parameterlist* list)
03359 {
03360
03361 cpl_parameter* p;
03362
03363
03364 if (!list) {
03365 return;
03366 }
03367
03368 p = cpl_parameter_new_value("giraffe.biasremoval.remove",
03369 CPL_TYPE_BOOL,
03370 "Enable bias removal",
03371 "giraffe.biasremoval",
03372 TRUE);
03373 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "remove-bias");
03374 cpl_parameterlist_append(list, p);
03375
03376
03377 p = cpl_parameter_new_enum("giraffe.biasremoval.method",
03378 CPL_TYPE_STRING,
03379 "Bias removal method",
03380 "giraffe.biasremoval",
03381 "PROFILE", 10, "UNIFORM", "PLANE", "PROFILE",
03382 "MASTER", "ZMASTER", "MASTER+PLANE",
03383 "ZMASTER+PLANE", "CURVE", "MASTER+CURVE",
03384 "ZMASTER+CURVE");
03385 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-method");
03386 cpl_parameterlist_append(list, p);
03387
03388
03389 p = cpl_parameter_new_value("giraffe.biasremoval.areas",
03390 CPL_TYPE_STRING,
03391 "Bias areas to use "
03392 "(Xl0:Xr0:Yl0:Yr0, ... ,Xln:Xrn:Yln:Yrn)",
03393 "giraffe.biasremoval",
03394 "5:40:0:4095");
03395 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-areas");
03396 cpl_parameterlist_append(list, p);
03397
03398
03399 p = cpl_parameter_new_value("giraffe.biasremoval.sigma",
03400 CPL_TYPE_DOUBLE,
03401 "Sigma Clipping: sigma threshold factor",
03402 "giraffe.biasremoval",
03403 2.5);
03404 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-sigma");
03405 cpl_parameterlist_append(list, p);
03406
03407
03408 p = cpl_parameter_new_value("giraffe.biasremoval.iterations",
03409 CPL_TYPE_INT,
03410 "Sigma Clipping: maximum number of "
03411 "iterations",
03412 "giraffe.biasremoval",
03413 5);
03414 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-niter");
03415 cpl_parameterlist_append(list, p);
03416
03417
03418 p = cpl_parameter_new_value("giraffe.biasremoval.fraction",
03419 CPL_TYPE_DOUBLE,
03420 "Sigma Clipping: minimum fraction of points "
03421 "accepted/total [0.0..1.0]",
03422 "giraffe.biasremoval",
03423 0.8);
03424 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-mfrac");
03425 cpl_parameterlist_append(list, p);
03426
03427
03428 p = cpl_parameter_new_value("giraffe.biasremoval.xorder",
03429 CPL_TYPE_INT,
03430 "Order of X polynomial fit (CURVE method "
03431 "only)",
03432 "giraffe.biasremoval",
03433 1);
03434 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-xorder");
03435 cpl_parameterlist_append(list, p);
03436
03437 p = cpl_parameter_new_value("giraffe.biasremoval.yorder",
03438 CPL_TYPE_INT,
03439 "Order of Y polynomial fit (CURVE method "
03440 "only)",
03441 "giraffe.biasremoval",
03442 1);
03443 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-yorder");
03444 cpl_parameterlist_append(list, p);
03445
03446
03447 p = cpl_parameter_new_value("giraffe.biasremoval.xstep",
03448 CPL_TYPE_INT,
03449 "Sampling step along X (CURVE method only)",
03450 "giraffe.biasremoval",
03451 1);
03452 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-xstep");
03453 cpl_parameterlist_append(list, p);
03454
03455
03456 p = cpl_parameter_new_value("giraffe.biasremoval.ystep",
03457 CPL_TYPE_INT,
03458 "Sampling step along Y (CURVE method only)",
03459 "giraffe.biasremoval",
03460 1);
03461 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-ystep");
03462 cpl_parameterlist_append(list, p);
03463
03464 return;
03465 }