GIRAFFE Pipeline Reference Manual

gibias.c
1/*
2 * This file is part of the GIRAFFE Pipeline
3 * Copyright (C) 2002-2019 European Southern Observatory
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18 */
19
20#ifdef HAVE_CONFIG_H
21# include <config.h>
22#endif
23
24#include <string.h>
25#include <math.h>
26#include <errno.h>
27
28#include <cxmemory.h>
29#include <cxmessages.h>
30#include <cxstrutils.h>
31#include <cxutils.h>
32
33#include <cpl_msg.h>
34#include <cpl_parameter.h>
35
36#include "gimacros.h"
37#include "gialias.h"
38#include "giarray.h"
39#include "gimatrix.h"
40#include "gichebyshev.h"
41#include "gimath.h"
42#include "gibias.h"
43
44
53struct GiBiasResults {
54
55 cxdouble mean;
56 cxdouble sigma;
57 cxdouble rms;
58
59 cx_string* limits;
60
61 cpl_matrix* coeffs;
62
63 cpl_image* model;
64
65};
66
67typedef struct GiBiasResults GiBiasResults;
68
69
70/*
71 * Clear a bias results structure
72 */
73
74inline static void
75_giraffe_biasresults_clear(GiBiasResults *self)
76{
77
78 if (self != NULL) {
79
80 self->mean = 0.;
81 self->sigma = 0.;
82 self->rms = 0.;
83
84 if (self->limits) {
85 cx_string_delete(self->limits);
86 self->limits = NULL;
87 }
88
89 if (self->coeffs) {
90 cpl_matrix_delete(self->coeffs);
91 self->coeffs = NULL;
92 }
93
94 if (self->model) {
95 cpl_image_delete(self->model);
96 self->model = NULL;
97 }
98
99 }
100
101 return;
102
103}
104
105
106/*
107 * Transform bias method and option value into a string
108 */
109
110inline static void
111_giraffe_method_string(cx_string *string, GiBiasMethod method,
112 GiBiasOption option)
113{
114
115 switch (method) {
116 case GIBIAS_METHOD_UNIFORM:
117 cx_string_set(string, "UNIFORM");
118 break;
119
120 case GIBIAS_METHOD_PLANE:
121 cx_string_set(string, "PLANE");
122 break;
123
124 case GIBIAS_METHOD_CURVE:
125 cx_string_set(string, "CURVE");
126 break;
127
128 case GIBIAS_METHOD_PROFILE:
129 cx_string_set(string, "PROFILE");
130 break;
131
132 case GIBIAS_METHOD_MASTER:
133 cx_string_set(string, "MASTER");
134 break;
135
136 case GIBIAS_METHOD_ZMASTER:
137 cx_string_set(string, "ZMASTER");
138 break;
139
140 default:
141 break;
142 }
143
144 if (option != GIBIAS_OPTION_UNDEFINED) {
145 switch (option) {
146 case GIBIAS_OPTION_PLANE:
147 cx_string_append(string, "+PLANE");
148 break;
149
150 case GIBIAS_OPTION_CURVE:
151 cx_string_append(string, "+CURVE");
152 break;
153
154 default:
155 break;
156 }
157 }
158
159 return;
160
161}
162
163
164inline static void
165_giraffe_stringify_coefficients(cx_string *string, cpl_matrix *matrix)
166{
167
168 register cxint i, j;
169
170 cxchar *tmp = cx_line_alloc();
171
172 cxint nr = cpl_matrix_get_nrow(matrix);
173 cxint nc = cpl_matrix_get_ncol(matrix);
174
175 cxsize sz = cx_line_max();
176
177 cxdouble *data = cpl_matrix_get_data(matrix);
178
179
180 for (i = 0; i < nr; i++) {
181 for (j = 0; j < nc; j++) {
182 snprintf(tmp, sz, "%g", data[i * nc + j]);
183 cx_string_append(string, tmp);
184
185 if (i != nr - 1 || j < nc - 1) {
186 cx_string_append(string, ":");
187 }
188 }
189 }
190
191 cx_free(tmp);
192
193 return;
194
195}
196
197
213inline static cxbool
214_giraffe_compare_overscans(const GiImage* image1, const GiImage* image2)
215{
216
217 cxint32 l1ovscx = -1;
218 cxint32 l1ovscy = -1;
219 cxint32 l1prscx = -1;
220 cxint32 l1prscy = -1;
221 cxint32 l2ovscx = -1;
222 cxint32 l2ovscy = -1;
223 cxint32 l2prscx = -1;
224 cxint32 l2prscy = -1;
225
226 cpl_propertylist *l1, *l2;
227
228
229 cx_assert(image1 != NULL && image2 != NULL);
230
231 l1 = giraffe_image_get_properties(image1);
232 l2 = giraffe_image_get_properties(image2);
233
234 if (cpl_propertylist_has(l1, GIALIAS_OVSCX)) {
235 l1ovscx = cpl_propertylist_get_int(l1, GIALIAS_OVSCX);
236 }
237 if (cpl_propertylist_has(l1, GIALIAS_OVSCY)) {
238 l1ovscy = cpl_propertylist_get_int(l1, GIALIAS_OVSCY);
239 }
240 if (cpl_propertylist_has(l1, GIALIAS_PRSCX)) {
241 l1prscx = cpl_propertylist_get_int(l1, GIALIAS_PRSCX);
242 }
243 if (cpl_propertylist_has(l1, GIALIAS_PRSCY)) {
244 l1prscy = cpl_propertylist_get_int(l1, GIALIAS_PRSCY);
245 }
246
247 if (cpl_propertylist_has(l2, GIALIAS_OVSCX)) {
248 l2ovscx = cpl_propertylist_get_int(l2, GIALIAS_OVSCX);
249 }
250 if (cpl_propertylist_has(l2, GIALIAS_OVSCY)) {
251 l2ovscy = cpl_propertylist_get_int(l2, GIALIAS_OVSCY);
252 }
253 if (cpl_propertylist_has(l2, GIALIAS_PRSCX)) {
254 l2prscx = cpl_propertylist_get_int(l2, GIALIAS_PRSCX);
255 }
256 if (cpl_propertylist_has(l2, GIALIAS_PRSCY)) {
257 l2prscy = cpl_propertylist_get_int(l2, GIALIAS_PRSCY);
258 }
259
260 if (l1ovscx != l2ovscx || l1ovscy != l2ovscy) {
261 return FALSE;
262 }
263
264 if (l1prscx != l2prscx || l1prscy != l2prscy) {
265 return FALSE;
266 }
267
268 return TRUE;
269
270}
271
272
273/*
274 * @brief
275 * Parse a bias area specification string.
276 *
277 * @param string The string to parse.
278 *
279 * @return
280 * The function returns a newly allocated matrix containing the limits
281 * of the specified areas.
282 *
283 * TBD
284 */
285
286inline static cpl_matrix*
287_giraffe_bias_get_areas(const cxchar* string)
288{
289
290 cxchar** regions = NULL;
291
292 cpl_matrix* areas = NULL;
293
294
295 cx_assert(string != NULL);
296
297 regions = cx_strsplit(string, ",", -1);
298
299 if (regions == NULL) {
300
301 return NULL;
302
303 }
304 else {
305
306 const cxsize nvalues = 4;
307
308 cxsize i = 0;
309 cxsize nregions = 0;
310
311 while (regions[i] != NULL) {
312 ++i;
313 }
314
315 nregions = i;
316 areas = cpl_matrix_new(nregions, nvalues);
317
318 i = 0;
319 while (regions[i] != NULL) {
320
321 register cxsize j = 0;
322
323 cxchar** limits = cx_strsplit(regions[i], ":", nvalues);
324
325
326 if (limits == NULL) {
327
328 cpl_matrix_delete(areas);
329 areas = NULL;
330
331 return NULL;
332
333 }
334
335 for (j = 0; j < nvalues; ++j) {
336
337 cxchar* last = NULL;
338
339 cxint status = 0;
340 cxlong value = 0;
341
342
343 if (limits[j] == NULL || *limits[j] == '\0') {
344 break;
345 }
346
347 status = errno;
348 errno = 0;
349
350 value = strtol(limits[j], &last, 10);
351
352 /*
353 * Check for various errors
354 */
355
356 if ((errno == ERANGE &&
357 (value == LONG_MAX || value == LONG_MIN)) ||
358 (errno != 0 && value == 0)) {
359
360 break;
361
362 }
363
364 errno = status;
365
366
367 /*
368 * Check that the input string was parsed completely.
369 */
370
371 if (*last != '\0') {
372 break;
373 }
374
375 cpl_matrix_set(areas, i, j, value);
376
377 }
378
379 cx_strfreev(limits);
380 limits = NULL;
381
382 if (j != nvalues) {
383
384 cpl_matrix_delete(areas);
385 areas = NULL;
386
387 cx_strfreev(regions);
388 regions = NULL;
389
390 return NULL;
391
392 }
393
394 ++i;
395
396 }
397
398 cx_strfreev(regions);
399 regions = NULL;
400
401 }
402
403 return areas;
404
405}
406
407
408/*
409 * @brief
410 * Compute bias mean and sigma values over bias areas.
411 *
412 * @param image Image over which to calculate mean and sigma
413 * inside bias areas
414 * @param areas Bias areas specifications [N,4] as a matrix
415 * where N is the number of areas
416 * @param kappa Sigma clipping algorithm: kappa value
417 * @param numiter Sigma clipping algorithm: max number of iterations
418 * @param maxfraction Sigma clipping algorithm: max. fraction of pixels
419 *
420 * @param results A @em GiBiasResults structure which contains
421 * the computed mean bias value, its standard deviation
422 * and the coordinates specifications of the areas used
423 * on return.
424 *
425 * @returns EXIT_SUCCESS on sucess, negative value otherwise
426 *
427 * Computes mean average of pixels value over areas defined by @em areas.
428 * @em areas is a matrix [N, 4] specifying the lower left and upper right
429 * corners of each of the @em N areas of the CCD used for bias determination.
430 *
431 * @code
432 * areas = [[xmin1, xmax1, ymin1, ymax1],
433 * [xmin2, xmax2, ymin2, ymax2],
434 * ....., ....., ....., .....
435 * [xminN, xmaxN, yminN, ymaxN]]
436 * @endcode
437 *
438 * Areas pixels whose value is too far from mean average are rejected by
439 * kappa-sigma clipping defined by parameters kappa, numiter and maxfraction.
440 * This rejection loops until good/bad points ratio is enough or the
441 * maximum number of iteration has been reached.
442 *
443 * The function returns the bias mean @em results->mean, the bias sigma
444 * @em results->sigma computed over the kappa-sigma clipped areas pixels
445 * value and @em results->limits string specifying the areas used as
446 * follows:
447 * @code
448 * "xmin1:xmax1:ymin1:ymax1;...;xminN:xmaxN:yminN:ymaxN"
449 * @endcode
450 *
451 * @see giraffe_bias_compute_plane, giraffe_compute_remove_bias
452 */
453
454inline static cxint
455_giraffe_bias_compute_mean(GiBiasResults* results, const cpl_image* image,
456 const cpl_matrix* areas, cxdouble kappa, cxint numiter,
457 cxdouble maxfraction)
458{
459
460 const cxchar* const fctid = "giraffe_bias_compute_mean";
461
462
463 const cxdouble* pdimg = NULL;
464
465 cxint j, l, k;
466 cxint img_dimx, img_dimy;
467 /*cxint ba_num_cols;*/
468 cxint ba_num_rows;
469 cxint curriter = 0;
470 cxint x0, x1, x2, x3;
471
472 cxdouble sigma = 0.;
473
474 cxlong ntotal = 0L;
475 cxlong naccepted = 0L;
476 cxlong n = 0L;
477 cxlong pixcount = 0L;
478
479 cxdouble currfraction = 2.;
480
481 cx_string* tmp = NULL;
482
483 cpl_matrix* matrix_zz1;
484 cpl_matrix* matrix_zz1diff;
485
486
487 /*
488 * Initialization
489 */
490
491 if (results->limits == NULL) {
492 cpl_msg_info(fctid, "Unable to store biaslimits return parameter, "
493 "aborting...");
494 return -3;
495 }
496
497 if (cpl_image_get_type(image) != CPL_TYPE_DOUBLE) {
498 cpl_msg_info(fctid, "Only images allowed of type double, "
499 "aborting ...");
500 return -3;
501 }
502
503
504 /*
505 * Preprocessing
506 */
507
508 /*
509 * Validate Bias Areas and calculate number of points
510 */
511
512 if (areas == NULL) {
513 cpl_msg_info(fctid, "Bias Areas: Missing bias areas, "
514 "aborting ...");
515 return -1;
516 }
517
518 /*ba_num_cols = cpl_matrix_get_ncol(areas);*/
519 ba_num_rows = cpl_matrix_get_nrow(areas);
520
521 for (j = 0; j < ba_num_rows; j++) {
522 x3 = (cxint) cpl_matrix_get(areas, j, 3);
523 x2 = (cxint) cpl_matrix_get(areas, j, 2);
524 x1 = (cxint) cpl_matrix_get(areas, j, 1);
525 x0 = (cxint) cpl_matrix_get(areas, j, 0);
526
527 ntotal += (cxulong) ((x3 - x2 + 1) * (x1 - x0 + 1));
528 }
529
530 if (ntotal <= 0) {
531 cpl_msg_info(fctid, "Bias Areas: Inconsistent specification, "
532 "aborting ...");
533 return -1;
534 }
535
536 matrix_zz1 = cpl_matrix_new(ntotal, 1);
537 matrix_zz1diff = cpl_matrix_new(ntotal, 1);
538
539
540 /*
541 * Compute Bias Areas Values
542 */
543
544 img_dimx = cpl_image_get_size_x(image);
545 img_dimy = cpl_image_get_size_y(image);
546
547 cx_string_set(results->limits, "");
548
549 for (j = 0; j < ba_num_rows; j++) {
550
551 x3 = (cxint) cpl_matrix_get(areas, j, 3);
552 x2 = (cxint) cpl_matrix_get(areas, j, 2);
553 x1 = (cxint) cpl_matrix_get(areas, j, 1);
554 x0 = (cxint) cpl_matrix_get(areas, j, 0);
555
556 if ((x0 > img_dimx) || (x1 > img_dimx) ||
557 (x2 > img_dimy) || (x3 > img_dimy)) {
558 continue;
559 }
560
561 tmp = cx_string_new();
562
563 cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
564 cx_string_append(results->limits, cx_string_get(tmp));
565
566 cx_string_delete(tmp);
567 tmp = NULL;
568
569 pdimg = cpl_image_get_data_double_const(image);
570
571 for (l = x2; l < x3 + 1; l++) {
572 for (k = x0; k < x1 + 1; k++) {
573 cpl_matrix_set(matrix_zz1, n, 1, pdimg[k + l * img_dimx]);
574 n++;
575 }
576 }
577 }
578
579 if (n != ntotal) {
580 cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting ...");
581
582 cpl_matrix_delete(matrix_zz1);
583 cpl_matrix_delete(matrix_zz1diff);
584
585 return -3;
586 }
587
588 cpl_msg_info(fctid, "Bias Areas: Using %s",
589 cx_string_get(results->limits));
590
591
592 /*
593 * Processing
594 */
595
596 /*
597 * Perform Sigma Clipping...
598 */
599
600 cpl_msg_info(fctid, "Sigma Clipping : Start");
601
602 naccepted = ntotal;
603
604 while ((naccepted > 0) && (curriter < numiter) &&
605 (currfraction > maxfraction)) {
606
607 cxdouble ksigma = 0.;
608
609 results->mean = cpl_matrix_get_mean(matrix_zz1);
610 sigma = giraffe_matrix_sigma_mean(matrix_zz1, results->mean);
611
612 for (k = 0; k < cpl_matrix_get_nrow(matrix_zz1); k++) {
613 cpl_matrix_set(matrix_zz1diff, k, 0,
614 cpl_matrix_get(matrix_zz1, k, 0) - results->mean);
615 }
616
617 cpl_msg_info(fctid, "bias[%d]: mean = %5g, sigma = %6.3g, "
618 "ratio = %6.3g, accepted = %ld\n", curriter,
619 results->mean, sigma, currfraction, naccepted);
620
621 ksigma = sigma * kappa;
622
623 pixcount = 0L;
624 for (l = 0; l < cpl_matrix_get_nrow(matrix_zz1); l++) {
625 if (fabs(cpl_matrix_get(matrix_zz1diff, l, 0)) > ksigma) {
626 continue;
627 }
628
629 cpl_matrix_set(matrix_zz1, pixcount, 0,
630 cpl_matrix_get(matrix_zz1, l, 0));
631 ++pixcount;
632 }
633
634 cpl_matrix_resize(matrix_zz1, 0, 0, 0, pixcount -
635 cpl_matrix_get_nrow(matrix_zz1));
636
637 cpl_matrix_resize(matrix_zz1diff, 0, 0, 0, pixcount -
638 cpl_matrix_get_nrow(matrix_zz1diff));
639
640 if (pixcount == naccepted) {
641 break;
642 }
643
644 naccepted = pixcount;
645
646 currfraction = (cxdouble) naccepted / (cxdouble) ntotal;
647 ++curriter;
648 }
649
650 cpl_msg_info(fctid, "Sigma Clipping : End");
651
652
653 /*
654 * Postprocessing
655 */
656
657 results->mean = cpl_matrix_get_mean(matrix_zz1);
658 results->rms = giraffe_matrix_sigma_mean(matrix_zz1, results->mean);
659
660 results->sigma = results->rms / sqrt(cpl_matrix_get_nrow(matrix_zz1));
661
662
663 cpl_msg_info(fctid, "Sigma Clipping Results : bias[%d]: mean = %5g, "
664 "sigma = %6.3g, ratio = %6.3g, accepted=%ld\n", curriter,
665 results->mean, results->rms, currfraction, naccepted);
666
667
668 /*
669 * Cleanup
670 */
671
672 if (matrix_zz1 != NULL) {
673 cpl_matrix_delete(matrix_zz1);
674 }
675
676 if (matrix_zz1diff != NULL) {
677 cpl_matrix_delete(matrix_zz1diff);
678 }
679
680 return EXIT_SUCCESS;
681
682}
683
684
685/*
686 * @brief
687 * Compute bias plane over bias areas.
688 *
689 * @param image Image over which to calculate mean and sigma
690 * inside bias areas
691 * @param areas Bias areas specifications [N,4] as a matrix
692 * where N is the number of areas
693 * @param kappa Sigma clipping algorithm: kappa value
694 * @param numiter Sigma clipping algorithm: max number of iterations
695 * @param maxfraction Sigma clipping algorithm: max. fraction of pixels
696 *
697 * @param results A @em GiBiasResults structure which contains
698 * the computed coefficients of the fitted bias plane,
699 * their standard deviation and the coordinates
700 * specifications of the areas used on return.
701 *
702 * @returns EXIT_SUCCESS on sucess, negative value otherwise
703 *
704 * Fit a plane through the pixels value over areas defined by @em areas.
705 * @em areas is a matrix [N, 4] specifying the lower left and upper right
706 * corners of each of the @em N areas of the CCD used for bias determination.
707 *
708 * @code
709 * areas = [[xmin1, xmax1, ymin1, ymax1],
710 * [xmin2, xmax2, ymin2, ymax2],
711 * ....., ....., ....., .....
712 * [xminN, xmaxN, yminN, ymaxN]]
713 * @endcode
714 *
715 * Areas pixels whose value is too far from mean average are rejected by
716 * kappa-sigma clipping defined by parameters kappa, numiter and maxfraction.
717 * This rejection loops until good/bad points ratio is enough or the
718 * maximum number of iteration has been reached.
719 *
720 * The function returns the bias fitted plane coefficients @em BP1, BP2, BP3:
721 * @code
722 * [BiasPlane] = [img] * inv([1XY]) over the bias areas where:
723 * [BiasPlane] = [BP1,BP2,BP3] coefficients of plane z = BP1+BP2*x+BP3*y
724 * [Z] = 1-column matrix of pixel values
725 * [1XY] = 3-columns matrix of all pixels coordinates [1,X,Y]
726 * @endcode
727 *
728 * the bias sigma computed over the kappa sigma clipped areas pixels value
729 * and a string specifying the areas used as follows
730 * @code
731 * "xmin1:xmax1:ymin1:ymax1;...;xminN:xmaxN:yminN:ymaxN"
732 * @endcode
733 *
734 * @see giraffe_bias_compute_mean, giraffe_compute_remove_bias
735 */
736
737inline static cxint
738_giraffe_bias_compute_plane(GiBiasResults* results, const cpl_image* image,
739 const cpl_matrix* areas, cxdouble kappa,
740 cxint numiter, cxdouble maxfraction)
741{
742
743 const cxchar* const fctid = "giraffe_bias_compute_plane";
744
745
746 cxint j = 0;
747 cxint nx = 0;
748 cxint ny = 0;
749 cxint nareas = 0;
750 cxint iteration = 0;
751
752 cxsize n = 0;
753 cxsize ntotal = 0;
754 cxsize naccepted = 0;
755
756 cxdouble fraction = 1.;
757 cxdouble sigma = 0.;
758
759 cpl_matrix* xbs = NULL;
760 cpl_matrix* ybs = NULL;
761 cpl_matrix* zbs = NULL;
762 cpl_matrix* coeffs = NULL;
763
764
765 cx_assert(results->limits != NULL);
766 cx_assert(results->coeffs == NULL);
767
768 cx_assert(areas != NULL);
769
770 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
771
772
773 /*
774 * Validate Bias Areas and calculate number of points
775 */
776
777 nareas = cpl_matrix_get_nrow(areas);
778
779 for (j = 0; j < nareas; j++) {
780
781 cxint x3 = (cxint) cpl_matrix_get(areas, j, 3);
782 cxint x2 = (cxint) cpl_matrix_get(areas, j, 2);
783 cxint x1 = (cxint) cpl_matrix_get(areas, j, 1);
784 cxint x0 = (cxint) cpl_matrix_get(areas, j, 0);
785
786 ntotal += (cxsize) ((x3 - x2 + 1) * (x1 - x0 + 1));
787
788 }
789
790 if (ntotal <= 0) {
791
792 cpl_msg_info(fctid, "Bias Areas: Inconsistent specification, "
793 "aborting ...");
794 return -1;
795
796 }
797
798 nx = cpl_image_get_size_x(image);
799 ny = cpl_image_get_size_y(image);
800
801 cpl_msg_info(fctid, "Bias Areas: specified are %zu points in %dx%d "
802 "image", ntotal, nx, ny);
803
804
805 /*
806 * Compute Bias Areas Values
807 */
808
809 results->mean = 0.;
810 results->sigma = 0.;
811 results->rms = 0.;
812
813 cx_string_set(results->limits, "");
814
815 xbs = cpl_matrix_new(ntotal, 1);
816 ybs = cpl_matrix_new(ntotal, 1);
817 zbs = cpl_matrix_new(1, ntotal);
818
819 for (j = 0; j < nareas; ++j) {
820
821 const cxdouble* _img = cpl_image_get_data_double_const(image);
822
823 cxint k = 0;
824
825 cx_string* tmp = NULL;
826
827
828 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
829 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
830 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
831 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
832
833 if ((x0 > nx) || (x1 > nx) || (x2 > ny) || (x3 > ny)) {
834 continue;
835 }
836
837 tmp = cx_string_new();
838
839 cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
840 cx_string_append(results->limits, cx_string_get(tmp));
841
842 cx_string_delete(tmp);
843 tmp = NULL;
844
845 for (k = x2; k < x3 + 1; ++k) {
846
847 register cxint l = 0;
848
849
850 for (l = x0; l < x1 + 1; ++l) {
851
852 cpl_matrix_set(xbs, n, 0, l);
853 cpl_matrix_set(ybs, n, 0, k);
854 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
855 ++n;
856
857 }
858
859 }
860
861 }
862
863 cpl_matrix_set_size(xbs, n, 1);
864 cpl_matrix_set_size(ybs, n, 1);
865 cpl_matrix_set_size(zbs, 1, n);
866
867 if (n != ntotal) {
868
869 cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting...");
870
871 cpl_matrix_delete(xbs);
872 cpl_matrix_delete(ybs);
873 cpl_matrix_delete(zbs);
874
875 return -1;
876
877 }
878
879 ntotal = n;
880
881 cpl_msg_info(fctid, "Bias Areas: Using %s [%zu pixels]",
882 cx_string_get(results->limits), ntotal);
883
884
885 /*
886 * Perform Sigma Clipping
887 */
888
889 cpl_msg_info(fctid, "Sigma Clipping : Start");
890
891 naccepted = ntotal;
892
893 while ((naccepted > 0) && (iteration < numiter) &&
894 (fraction > maxfraction)) {
895
896 cxsize k = 0;
897
898 cxdouble ksigma = 0.;
899
900 cpl_matrix* base = NULL;
901 cpl_matrix* fit = NULL;
902
903
904 base = cpl_matrix_new(3, naccepted);
905
906 if (base == NULL) {
907
908 cpl_msg_info(fctid, "Sigma Clipping: Error creating design "
909 "matrix");
910
911 cpl_matrix_delete(zbs);
912 cpl_matrix_delete(ybs);
913 cpl_matrix_delete(xbs);
914
915 return -2;
916 }
917
918 for (k = 0; k < naccepted; ++k) {
919
920 cpl_matrix_set(base, 0, k, 1.);
921 cpl_matrix_set(base, 1, k, cpl_matrix_get(xbs, k, 0));
922 cpl_matrix_set(base, 2, k, cpl_matrix_get(ybs, k, 0));
923
924 }
925
926 cpl_matrix_delete(coeffs);
927 coeffs = NULL;
928
929 coeffs = giraffe_matrix_leastsq(base, zbs);
930
931 if (coeffs == NULL) {
932
933 cpl_msg_info(fctid, "Sigma Clipping : Error in least square "
934 "solution, aborting...");
935
936 cpl_matrix_delete(base);
937 base = NULL;
938
939 cpl_matrix_delete(xbs);
940 cpl_matrix_delete(ybs);
941 cpl_matrix_delete(zbs);
942
943 return -2;
944
945 }
946
947
948 /*
949 * Compute bias model fit and reject deviating pixels
950 */
951
952 fit = cpl_matrix_product_create(coeffs, base);
953
954 cpl_matrix_delete(base);
955 base = NULL;
956
957 results->mean = cpl_matrix_get_mean(fit);
958
959 sigma = giraffe_matrix_sigma_fit(zbs, fit);
960
961 cpl_msg_info(fctid, "Sigma Clipping : bias plane[%d]: %g + "
962 "%g * x + %g * y, sigma = %.5g, ratio = %.4g, "
963 "accepted = %zu\n", iteration,
964 cpl_matrix_get(coeffs, 0, 0),
965 cpl_matrix_get(coeffs, 0, 1),
966 cpl_matrix_get(coeffs, 0, 2),
967 sigma, fraction, naccepted);
968
969
970 /*
971 * Perform Clipping
972 */
973
974 ksigma = sigma * kappa;
975
976 n = 0;
977
978 for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
979
980 register cxdouble z = cpl_matrix_get(zbs, 0, j);
981
982 if (fabs(cpl_matrix_get(fit, 0, j) - z) > ksigma) {
983 continue;
984 }
985
986 cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
987
988 cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
989
990 cpl_matrix_set(zbs, 0, n, z);
991 ++n;
992
993 }
994
995 cpl_matrix_set_size(xbs, n, 1);
996 cpl_matrix_set_size(ybs, n, 1);
997 cpl_matrix_set_size(zbs, 1, n);
998
999 cpl_matrix_delete(fit);
1000 fit = NULL;
1001
1002 if (n == naccepted) {
1003 break;
1004 }
1005
1006 naccepted = n;
1007
1008 fraction = (cxdouble)naccepted / (cxdouble)ntotal;
1009 ++iteration;
1010
1011 }
1012
1013 cpl_msg_info(fctid, "Sigma Clipping : End");
1014
1015
1016 /*
1017 * Store fit results
1018 */
1019
1020 results->coeffs = coeffs;
1021 results->rms = sigma;
1022
1023 // FIXME: The following is not correct. It should be the error of
1024 // results->mean which has to be obtained from error propagation.
1025 // At least the following is connected to the fitted model,
1026 // instead of the original code which computed only the standard
1027 // deviation of the raw data.
1028
1029 results->sigma = sigma / sqrt(cpl_matrix_get_ncol(zbs));
1030
1031 cpl_msg_info(fctid, "Sigma Clipping Results (%d/%zu, sigma = %g)",
1032 iteration, naccepted, results->rms);
1033
1034
1035 /*
1036 * Cleanup
1037 */
1038
1039 cpl_matrix_delete(xbs);
1040 xbs = NULL;
1041
1042 cpl_matrix_delete(ybs);
1043 ybs = NULL;
1044
1045 cpl_matrix_delete(zbs);
1046 zbs = NULL;
1047
1048 return EXIT_SUCCESS;
1049
1050}
1051
1052
1053/*
1054 * @brief
1055 * Compute bias curve over bias areas.
1056 *
1057 * @param results computed bias curve over bias areas coordinates
1058 * specifications
1059 * @param image Image over which to calculate bias curve
1060 * inside bias areas
1061 * @param areas bias areas specifications [N,4] as a matrix
1062 * where N is the number of areas
1063 * @param kappa sigma clipping algorithm: kappa value
1064 * @param numiter sigma clipping algorithm: max number of iterations
1065 * @param maxfraction sigma clipping algorithm: max. fraction of pixels
1066 * @param xdeg Polynom order for fit in x direction
1067 * @param ydeg Polynom order for fit in y direction
1068 * @param xstep Step in x direction for polynomial fit
1069 * @param ystep Step in y direction for polynomial fit
1070 *
1071 * @return
1072 * EXIT_SUCCESS on sucess, negative value otherwise
1073 *
1074 * Fit a surface using a rectangular mesh defined by @em xstep and
1075 * @em ystep in @em areas.
1076 * Mesh pixels whose value is too far from fitted surface are rejected by
1077 * kappa-sigma clipping defined by parameters of @em kappa.
1078 * This rejection loops until good/bad points ratio is enough
1079 * (maxfraction) or the maximum number of iteration (numiter) is reached.
1080 *
1081 * The function returns the bias fitted surface coefficients
1082 * @a results->biascoeffs computed using a 2D Chebyshev polynomial fit
1083 * whose X and Y order are given respectively by @a xdeg and @a ydeg,
1084 * The bias sigma @a results->biassigma computed over the
1085 * kappa-sigma clipped pixels value and a string specifying the areas
1086 * used as follows:
1087 * @code
1088 * "xstart:xend:xstep;ymin,ystart,ystep"
1089 * @endcode
1090 *
1091 * @see giraffe_bias_compute_plane, giraffe_bias_compute_remove_bias
1092 */
1093
1094inline static cxint
1095_giraffe_bias_compute_curve(GiBiasResults* results, const cpl_image* image,
1096 const cpl_matrix *areas, cxdouble kappa,
1097 cxint numiter, cxdouble maxfraction,
1098 cxdouble xdeg, cxdouble ydeg,
1099 cxdouble xstep, cxdouble ystep)
1100{
1101
1102 const cxchar* const fctid = "giraffe_bias_compute_curve";
1103
1104 cxint j = 0;
1105 cxint nx = 0;
1106 cxint ny = 0;
1107 cxint nareas = 0;
1108 cxint iteration = 0;
1109
1110 cxsize n = 0;
1111 cxsize ntotal = 0;
1112 cxsize naccepted = 0;
1113
1114 cxdouble fraction = 1.;
1115 cxdouble sigma = 0.;
1116
1117 cpl_matrix* xbs = NULL;
1118 cpl_matrix* ybs = NULL;
1119 cpl_matrix* zbs = NULL;
1120
1121 GiChebyshev2D* fit = NULL;
1122
1123
1124 cx_assert(results != NULL);
1125 cx_assert(results->limits != NULL);
1126 cx_assert(results->coeffs == NULL);
1127
1128 cx_assert(areas != NULL);
1129
1130 cx_assert(image != NULL);
1131 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
1132
1133
1134 /*
1135 * Validate Bias Areas and calculate number of points
1136 */
1137
1138 nareas = cpl_matrix_get_nrow(areas);
1139
1140 for (j = 0; j < nareas; ++j) {
1141
1142 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
1143 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
1144 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
1145 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
1146
1147 ntotal += (cxulong) (ceil((1. + x3 - x2) / ystep) *
1148 ceil((1. + x1 - x0) / xstep));
1149 }
1150
1151 nx = cpl_image_get_size_x(image);
1152 ny = cpl_image_get_size_y(image);
1153
1154 cpl_msg_info(fctid, "Bias Areas: Found %zu points in %dx%d image",
1155 ntotal, nx, ny);
1156
1157
1158 /*
1159 * Compute Bias Areas Values
1160 */
1161
1162 results->mean = 0.;
1163 results->sigma = 0.;
1164 results->rms = 0.;
1165
1166 cx_string_set(results->limits, "");
1167
1168 xbs = cpl_matrix_new(ntotal, 1);
1169 ybs = cpl_matrix_new(ntotal, 1);
1170 zbs = cpl_matrix_new(1, ntotal);
1171
1172 for (j = 0; j < nareas; ++j) {
1173
1174 const cxdouble* _img = cpl_image_get_data_double_const(image);
1175
1176 cxint k = 0;
1177
1178 cx_string* tmp = NULL;
1179
1180
1181 cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
1182 cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
1183 cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
1184 cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
1185
1186 if ((x0 > nx) || (x1 > nx) ||
1187 (x2 > ny) || (x3 > ny)) {
1188 continue;
1189 }
1190
1191 tmp = cx_string_new();
1192
1193 cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
1194 cx_string_append(results->limits, cx_string_get(tmp));
1195
1196 cx_string_delete(tmp);
1197 tmp = NULL;
1198
1199 for (k = x2; k < x3 + 1; k += ystep) {
1200
1201 register cxint l = 0;
1202
1203
1204 for (l = x0; l < x1 + 1; l += xstep) {
1205
1206 cpl_matrix_set(xbs, n, 0, l);
1207 cpl_matrix_set(ybs, n, 0, k);
1208 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
1209 ++n;
1210
1211 }
1212
1213 }
1214
1215 }
1216
1217 cpl_matrix_set_size(xbs, n, 1);
1218 cpl_matrix_set_size(ybs, n, 1);
1219 cpl_matrix_set_size(zbs, 1, n);
1220
1221 if (n <= 0) {
1222
1223 cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting...");
1224
1225 cpl_matrix_delete(xbs);
1226 cpl_matrix_delete(ybs);
1227 cpl_matrix_delete(zbs);
1228
1229 return -1;
1230
1231 }
1232
1233 ntotal = n;
1234
1235 cpl_msg_info(fctid, "Bias Areas: Using %s [%zu pixels]",
1236 cx_string_get(results->limits), ntotal);
1237
1238
1239 /*
1240 * Perform Sigma Clipping
1241 */
1242
1243 cpl_msg_info(fctid, "Sigma Clipping : Start");
1244
1245 naccepted = ntotal;
1246
1247 while ((naccepted > 0) && (iteration < numiter) &&
1248 (fraction > maxfraction)) {
1249
1250 cxint status = 0;
1251
1252 cxdouble ksigma = 0.;
1253
1254 cpl_matrix* base = NULL;
1255 cpl_matrix* coeffs = NULL;
1256 cpl_matrix* _coeffs = NULL;
1257 cpl_matrix* _fit = NULL;
1258
1259
1260 base = giraffe_chebyshev_base2d(0., 0., nx, ny,
1261 xdeg, ydeg, xbs, ybs);
1262
1263 if (base == NULL) {
1264
1265 cpl_msg_info(fctid, "Sigma Clipping: Error creating design "
1266 "matrix");
1267
1268 cpl_matrix_delete(zbs);
1269 cpl_matrix_delete(ybs);
1270 cpl_matrix_delete(xbs);
1271
1272 return -2;
1273 }
1274
1275 _coeffs = giraffe_matrix_leastsq(base, zbs);
1276
1277 cpl_matrix_delete(base);
1278 base = NULL;
1279
1280 if (_coeffs == NULL) {
1281
1282 cpl_msg_info(fctid, "Sigma Clipping : Error in least square "
1283 "solution, aborting...");
1284
1285 cpl_matrix_delete(xbs);
1286 cpl_matrix_delete(ybs);
1287 cpl_matrix_delete(zbs);
1288
1289 return -2;
1290
1291 }
1292
1293
1294 /*
1295 * Compute bias model fit and reject deviating pixels
1296 */
1297
1298 coeffs = cpl_matrix_wrap(xdeg, ydeg, cpl_matrix_get_data(_coeffs));
1299
1300 if (fit != NULL) {
1301 giraffe_chebyshev2d_delete(fit);
1302 fit = NULL;
1303 }
1304
1305 fit = giraffe_chebyshev2d_new(xdeg - 1, ydeg - 1);
1306 status = giraffe_chebyshev2d_set(fit, 0., nx, 0., ny,
1307 coeffs);
1308
1309 if (status != 0) {
1310
1311 giraffe_chebyshev2d_delete(fit);
1312 fit = NULL;
1313
1314 cpl_matrix_unwrap(coeffs);
1315 coeffs = NULL;
1316
1317 cpl_matrix_delete(_coeffs);
1318 _coeffs = NULL;
1319
1320 cpl_matrix_delete(xbs);
1321 cpl_matrix_delete(ybs);
1322 cpl_matrix_delete(zbs);
1323
1324 return -2;
1325
1326 }
1327
1328 cpl_matrix_unwrap(coeffs);
1329 coeffs = NULL;
1330
1331 cpl_matrix_delete(_coeffs);
1332 _coeffs = NULL;
1333
1334 _fit = cpl_matrix_new(1, cpl_matrix_get_ncol(zbs));
1335
1336 for (j = 0; j < cpl_matrix_get_ncol(_fit); ++j) {
1337
1338 cxdouble x = cpl_matrix_get(xbs, n, 0);
1339 cxdouble y = cpl_matrix_get(ybs, n, 0);
1340 cxdouble z = giraffe_chebyshev2d_eval(fit, x, y);
1341
1342 cpl_matrix_set(_fit, 0, j, z);
1343
1344 }
1345
1346 results->mean = cpl_matrix_get_mean(_fit);
1347
1348 sigma = giraffe_matrix_sigma_fit(zbs, _fit);
1349
1350 cpl_msg_info(fctid, "Sigma Clipping : bias surface[%d]: "
1351 "sigma = %8.5g, ratio = %7.4g, accepted = %zu\n",
1352 iteration, sigma, fraction, naccepted);
1353
1354
1355 /*
1356 * Perform Clipping
1357 */
1358
1359 ksigma = sigma * kappa;
1360
1361 n = 0;
1362
1363 for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
1364
1365 register cxdouble z = cpl_matrix_get(zbs, 0, j);
1366
1367 if (fabs(cpl_matrix_get(_fit, 0, j) - z) > ksigma) {
1368 continue;
1369 }
1370
1371 cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
1372 cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
1373 cpl_matrix_set(zbs, 0, n, z);
1374 ++n;
1375
1376 }
1377
1378 cpl_matrix_set_size(xbs, n, 1);
1379 cpl_matrix_set_size(ybs, n, 1);
1380 cpl_matrix_set_size(zbs, 1, n);
1381
1382 cpl_matrix_delete(_fit);
1383 _fit = NULL;
1384
1385
1386 if (n == naccepted) {
1387 break;
1388 }
1389
1390 naccepted = n;
1391
1392 fraction = (cxdouble)naccepted / (cxdouble)ntotal;
1393 ++iteration;
1394
1395 }
1396
1397 cpl_msg_info(fctid, "Sigma Clipping : End");
1398
1399
1400 /*
1401 * Store fit results
1402 */
1403
1404 results->coeffs = cpl_matrix_duplicate(giraffe_chebyshev2d_coeffs(fit));
1405 results->rms = sigma;
1406
1407 // FIXME: The following is not correct. It should be the error of
1408 // results->mean which has to be obtained from error propagation.
1409 // At least the following is connected to the fitted model,
1410 // instead of the original code which computed only the standard
1411 // deviation of the raw data.
1412
1413 results->sigma = sigma / sqrt(cpl_matrix_get_ncol(zbs));
1414
1415 cpl_msg_info(fctid, "Sigma Clipping Results (%d/%zu, sigma = %g)",
1416 iteration, naccepted, results->rms);
1417
1418
1419 /*
1420 * Cleanup
1421 */
1422
1423 giraffe_chebyshev2d_delete(fit);
1424 fit = NULL;
1425
1426 cpl_matrix_delete(xbs);
1427 xbs = NULL;
1428
1429 cpl_matrix_delete(ybs);
1430 ybs = NULL;
1431
1432 cpl_matrix_delete(zbs);
1433 zbs = NULL;
1434
1435 return EXIT_SUCCESS;
1436
1437}
1438
1439
1440/*
1441 * @brief
1442 * Compute a bias profile model from the collapsed bias areas.
1443 *
1444 * @param results Container for the computed mean bias, rms and the profile
1445 * model.
1446 * @param image The image from which the profile model is computed.
1447 * @param areas Image regions used to compute the bias model.
1448 *
1449 * @return
1450 * The function returns @c EXIT_SUCCESS on success, and a non-zero
1451 * value otherwise.
1452 *
1453 * The function computes a bias model profile along the axis @em axis from
1454 * the image areas @em areas of the input image @em image. The profile model
1455 * is computed by computing the mean pixel value along the axis
1456 * perpendicular to @em axis for each image area listed in @em areas.
1457 * The bias model value is then the average mean value of all given bias
1458 * areas for each pixel along @em axis.
1459 *
1460 * The given bias areas must have the same size along the axis @em axis!
1461 */
1462
1463inline static cxint
1464_giraffe_bias_compute_profile(GiBiasResults* results, const cpl_image* image,
1465 const cpl_matrix* areas, cxchar axis)
1466{
1467
1468 const cxchar* const fctid = "_giraffe_bias_compute_profile";
1469
1470
1471 cxint nx = 0;
1472 cxint ny = 0;
1473 cxint sx = 0;
1474 cxint sy = 0;
1475
1476 cxsize j = 0;
1477 cxsize npixel = 0;
1478 cxsize ntotal = 0;
1479 cxsize nareas = 0;
1480 cxsize nvalid = 0;
1481
1482 cxdouble rms = 0.;
1483 cxdouble sigma = 0.;
1484
1485 cpl_matrix* _areas = NULL;
1486
1487 cpl_image* profile = NULL;
1488 cpl_image* model = NULL;
1489
1490
1491 cx_assert(results != NULL);
1492 cx_assert(results->limits != NULL);
1493 cx_assert(results->model == NULL);
1494
1495 cx_assert(areas != NULL);
1496
1497 cx_assert(image != NULL);
1498 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
1499
1500 cx_assert((axis == 'x') || (axis == 'y'));
1501
1502
1503 nx = cpl_image_get_size_x(image);
1504 ny = cpl_image_get_size_y(image);
1505
1506
1507 /*
1508 * Validate Bias Areas and calculate number of points
1509 */
1510
1511 _areas = cpl_matrix_duplicate(areas);
1512 cx_assert(_areas != NULL);
1513
1514 nareas = cpl_matrix_get_nrow(_areas);
1515
1516 if (axis == 'y') {
1517
1518 for (j = 0; j < nareas; ++j) {
1519
1520 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1521 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1522 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1523 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1524
1525
1526 if (x_0 < 0) {
1527 x_0 = 0;
1528 cpl_matrix_set(_areas, j, 0, x_0);
1529 }
1530
1531 if (x_1 >= nx) {
1532 x_1 = nx - 1;
1533 cpl_matrix_set(_areas, j, 1, x_1);
1534 }
1535
1536 if (y_0 != 0) {
1537 y_0 = 0;
1538 cpl_matrix_set(_areas, j, 2, y_0);
1539 }
1540
1541 if (y_1 != ny - 1) {
1542 y_1 = ny - 1;
1543 cpl_matrix_set(_areas, j, 3, y_1);
1544 }
1545
1546 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1547 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
1548 }
1549
1550 }
1551
1552 sx = nareas;
1553 sy = ny;
1554 }
1555 else {
1556
1557 for (j = 0; j < nareas; ++j) {
1558
1559 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1560 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1561 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1562 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1563
1564
1565 if (x_0 != 0) {
1566 x_0 = 0;
1567 cpl_matrix_set(_areas, j, 0, x_0);
1568 }
1569
1570 if (x_1 != nx) {
1571 x_1 = nx - 1;
1572 cpl_matrix_set(_areas, j, 1, x_1);
1573 }
1574
1575 if (y_0 < 0) {
1576 y_0 = 0;
1577 cpl_matrix_set(_areas, j, 2, y_0);
1578 }
1579
1580 if (y_1 >= ny) {
1581 y_1 = ny - 1;
1582 cpl_matrix_set(_areas, j, 3, y_1);
1583 }
1584
1585 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1586 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
1587 }
1588
1589 }
1590
1591 sx = nareas;
1592 sy = nx;
1593
1594 }
1595
1596 cpl_msg_info(fctid, "Bias Areas: Found %zu points in %dx%d image",
1597 ntotal, nx, ny);
1598
1599
1600 /*
1601 * Compute Bias Areas Values
1602 */
1603
1604 results->mean = 0.;
1605 results->sigma = 0.;
1606 results->rms = 0.;
1607
1608 cx_string_set(results->limits, "");
1609
1610 profile = cpl_image_new(sx, sy, CPL_TYPE_DOUBLE);
1611
1612 for (j = 0; j < nareas; ++j) {
1613
1614 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1615 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1616 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1617 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1618
1619 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1620
1621 cx_string* tmp = cx_string_new();
1622
1623
1624 cx_string_sprintf(tmp, "%d:%d:%d:%d;", x_0, x_1, y_0, y_1);
1625 cx_string_append(results->limits, cx_string_get(tmp));
1626
1627 cx_string_delete(tmp);
1628 tmp = NULL;
1629
1630
1631 /*
1632 * Update the number of actually used (valid) areas.
1633 */
1634
1635 ++nvalid;
1636
1637
1638 /*
1639 * Compute the profile along the given axis for each bias region.
1640 */
1641
1642 if (axis == 'y') {
1643
1644 cxint i = 0;
1645 cxint sz = x_1 - x_0 + 1;
1646
1647 cxdouble residual = 0.;
1648
1649 const cxdouble* _img = cpl_image_get_data_double_const(image);
1650
1651 cxdouble* _profile = cpl_image_get_data_double(profile);
1652
1653
1654 for (i = y_0; i < y_1 + 1; ++i) {
1655
1656 register cxint k = i * nx;
1657 register cxint l = 0;
1658
1659 cxdouble m = giraffe_array_mean(&_img[k + x_0], sz);
1660
1661 _profile[i * sx + j] = m;
1662
1663 for (l = x_0; l < x_1 + 1; ++l) {
1664 residual += (_img[k + l] - m) * (_img[k + l] - m);
1665 }
1666 ++npixel;
1667
1668 }
1669
1670
1671 /*
1672 * Compute the mean RMS and the mean bias error
1673 */
1674
1675 rms += residual / (sz - 1);
1676 sigma += residual / (sz * (sz - 1));
1677
1678 }
1679 else {
1680
1681 cxint i = 0;
1682 cxint sz = y_1 - y_0 + 1;
1683
1684 cxdouble residual = 0.;
1685
1686 const cxdouble* _img = cpl_image_get_data_double_const(image);
1687
1688 cxdouble* _profile = cpl_image_get_data_double(profile);
1689 cxdouble* buffer = cx_calloc(sz, sizeof(cxdouble));
1690
1691
1692 for (i = x_0; i < x_1 + 1; ++i) {
1693
1694 register cxint l = 0;
1695
1696 register cxdouble m = 0.;
1697
1698
1699 for (l = y_0; l < y_1 + 1; ++l) {
1700 buffer[l - y_0] = _img[l * nx + i];
1701 }
1702
1703 m = giraffe_array_mean(buffer, sz);
1704 _profile[i * sx + j] = m;
1705
1706 for (l = 0; l < sz; ++l) {
1707 residual += (buffer[l] - m) * (buffer[l] - m);
1708 }
1709 ++npixel;
1710
1711 }
1712
1713
1714 /*
1715 * Compute the mean RMS and the mean bias error
1716 */
1717
1718 rms += residual / (sz - 1);
1719 sigma += residual / (sz * (sz - 1));
1720
1721 cx_free(buffer);
1722 buffer = NULL;
1723
1724 }
1725
1726 }
1727
1728 }
1729
1730 cpl_matrix_delete(_areas);
1731 _areas = NULL;
1732
1733
1734 /*
1735 * Compute average profile of all bias areas
1736 */
1737
1738 model = cpl_image_collapse_create(profile, 1);
1739
1740 cpl_image_delete(profile);
1741 profile = NULL;
1742
1743 cpl_image_divide_scalar(model, nvalid);
1744 results->model = model;
1745
1746 model = NULL;
1747
1748 results->mean = cpl_image_get_mean(results->model);
1749
1750 if (nvalid > 0) {
1751 results->rms = sqrt(rms / (sy * nvalid));
1752 results->sigma = sqrt(sigma / sy) / nvalid;
1753 }
1754
1755 return EXIT_SUCCESS;
1756
1757}
1758
1759
1760/*
1761 * @brief
1762 * Fit a Chebyshev polynomial to a bias profile model.
1763 *
1764 * @param results Container for the computed mean bias, rms and the profile
1765 * model.
1766 * @param image The image of the raw profile model.
1767 * @param order The order of the Chebyshev polynomial.
1768 *
1769 * @return
1770 * The function returns @c EXIT_SUCCESS on success, and a non-zero
1771 * value otherwise.
1772 *
1773 * The function fits a one-dimensional Chebyshev polynomial of the order
1774 * @em order to the one-dimensional, raw profile image @em profile.
1775 *
1776 * The raw profile image is expected to be an image with a single column!
1777 * It may be an alias for the profile model which is present in the results
1778 * container.
1779 *
1780 * If the fit is successful, the result is stored as the new profile model
1781 * in the results container. Any previous model is overwritten!
1782 */
1783
1784inline static cxint
1785_giraffe_bias_fit_profile(GiBiasResults* results, const cpl_image* profile,
1786 cxint order)
1787{
1788
1789 cxint i = 0;
1790
1791 cxint sy = cpl_image_get_size_y(profile);
1792 cxint nc = order + 1;
1793
1794 cxdouble* _profile = (cxdouble*)cpl_image_get_data_double_const(profile);
1795
1796 cpl_matrix* base = NULL;
1797 cpl_matrix* baseT = NULL;
1798 cpl_matrix* coeff = NULL;
1799 cpl_matrix* fit = NULL;
1800 cpl_matrix* x = cpl_matrix_new(sy, 1);
1801 cpl_matrix* y = cpl_matrix_wrap(sy, 1, _profile);
1802 cpl_matrix* covy = cpl_matrix_new(sy, sy);
1803 cpl_matrix* covc = cpl_matrix_new(nc, nc);
1804
1805
1806 if ((x == NULL) || (y == NULL) || (covy == NULL) || (covc == NULL)) {
1807
1808 cpl_matrix_delete(covc);
1809 cpl_matrix_delete(covy);
1810 cpl_matrix_unwrap(y);
1811 cpl_matrix_delete(x);
1812
1813 return -1;
1814 }
1815
1816
1817 /*
1818 * Set abscissa values for the least square fit, and create the
1819 * design matrix of the fit (the basis polynomials at the abscissa values).
1820 */
1821
1822 for (i = 0; i < sy; ++i)
1823 {
1824 cpl_matrix_set(x, i, 0, i);
1825 }
1826
1827 baseT = giraffe_chebyshev_base1d(0., sy, nc, x);
1828 base = cpl_matrix_transpose_create(baseT);
1829
1830 cpl_matrix_delete(baseT);
1831 baseT = NULL;
1832
1833 cpl_matrix_delete(x);
1834 x = NULL;
1835
1836 if (base == NULL) {
1837 cpl_matrix_unwrap(y);
1838 return -2;
1839 }
1840
1841 cpl_matrix_fill_diagonal(covy, results->rms * results->rms, 0);
1842
1843
1844 /*
1845 * Compute the coefficients of the fitting polynomial, and evaluate the
1846 * fit at the given abscissa values.
1847 */
1848
1849 coeff = giraffe_matrix_solve_cholesky(base, y, covy, covc);
1850
1851 cpl_matrix_delete(covy);
1852 covy = NULL;
1853
1854 if (coeff == NULL) {
1855
1856 cpl_matrix_delete(base);
1857 base = NULL;
1858
1859 cpl_matrix_unwrap(y);
1860 y = NULL;
1861
1862 return -3;
1863
1864 }
1865
1866
1867 fit = cpl_matrix_product_create(base, coeff);
1868
1869 cpl_matrix_delete(coeff);
1870 coeff = NULL;
1871
1872 cpl_matrix_delete(base);
1873 base = NULL;
1874
1875 if (fit == NULL) {
1876
1877 cpl_matrix_unwrap(y);
1878 y = NULL;
1879
1880 return -3;
1881
1882 }
1883
1884
1885 /*
1886 * Update the model in the results container
1887 */
1888
1889 for (i = 0; i < sy; ++i)
1890 {
1891 cpl_matrix_set(y, i, 0, cpl_matrix_get(fit, i, 0));
1892 }
1893
1894 cpl_matrix_delete(fit);
1895 fit = NULL;
1896
1897 cpl_matrix_unwrap(y);
1898 y = NULL;
1899
1900
1901 /*
1902 * Update mean bias and the bias error estimates. The error of the bias
1903 * is adopted to be the error of the constant term of the bias polynomial
1904 * model.
1905 */
1906
1907 results->mean = cpl_image_get_mean(results->model);
1908 results->sigma = sqrt(cpl_matrix_get(covc, 0, 0));
1909
1910#if 0
1911 cpl_image_save(results->model, "idiot.fits", CPL_BPP_IEEE_FLOAT,
1912 0, CPL_IO_DEFAULT);
1913#endif
1914
1915 /*
1916 * Cleanup
1917 */
1918
1919 cpl_matrix_delete(covc);
1920 covc = NULL;
1921
1922 return EXIT_SUCCESS;
1923}
1924
1925
1926/*
1927 * @brief
1928 * Actual Calculation/Logic Module of Bias Removal
1929 *
1930 * @param img Input image
1931 * @param biasareas Bias area to use
1932 * @param master_bias Master bias Image
1933 * @param bad_pixel_map Bad pixel map Image
1934 * @param method Method to use for bias removal
1935 * @param option Option to use for bias removal
1936 * @param model Model to use for bias removal
1937 * @param remove_bias Remove bias from result image
1938 * @param mbias Master Bias value read from Master Bias
1939 * FITS Keywords
1940 * @param xdeg X Degree of Polynomial of Chebyshev fit
1941 * @param ydeg Y Degree of Polynomial of Chebyshev fit
1942 * @param xstep X Step to use during Chebyshev fit
1943 * @param ystep Y Step to use during Chebyshev fit
1944 * @param sigma Sigma to use during Kappa Sigma Clipping
1945 * @param numiter Max Num Iterations to use during Kappa Sigma
1946 * Clipping
1947 * @param maxfraction Maximum Fraction to use during Kappa Sigma Clipping
1948 * @param results Values determined
1949 *
1950 * @return The function returns 0 on success, or an error code otherwise.
1951 *
1952 * TBD
1953 */
1954
1955inline static cxint
1956_giraffe_bias_compute(GiBiasResults* results, cpl_image* img,
1957 const cpl_matrix* biasareas,
1958 const cpl_image* master_bias,
1959 const cpl_image* bad_pixel_map,
1960 GiBiasMethod method, GiBiasOption option,
1961 GiBiasModel model, cxbool remove_bias, cxdouble mbias,
1962 cxdouble xdeg, cxdouble ydeg, cxdouble xstep,
1963 cxdouble ystep, cxdouble sigma, cxint numiter,
1964 cxdouble maxfraction)
1965{
1966
1967 const cxchar* const fctid = "giraffe_bias_compute";
1968
1969 cxint i;
1970 cxint j;
1971 cxint img_dimx;
1972 cxint img_dimy;
1973 cxint img_centerx;
1974 cxint img_centery;
1975 cxint master_bias_dimx;
1976 cxint master_bias_dimy;
1977 cxint bad_pixel_dimx;
1978 cxint bad_pixel_dimy;
1979
1980 cxint error_code = 0;
1981
1982 cx_string* s = NULL;
1983
1984 const cpl_matrix* coeff = NULL;
1985
1986
1987 cx_assert(results != NULL);
1988 cx_assert(results->limits == NULL);
1989 cx_assert(results->coeffs == NULL);
1990 cx_assert(results->model == NULL);
1991
1992 cx_assert(img != NULL);
1993 cx_assert(cpl_image_get_type(img) == CPL_TYPE_DOUBLE);
1994
1995 cx_assert(biasareas != NULL);
1996
1997
1998 /*
1999 * Preprocessing
2000 */
2001
2002 img_dimx = cpl_image_get_size_x(img);
2003 img_dimy = cpl_image_get_size_y(img);
2004
2005 img_centerx = img_dimx / 2;
2006 img_centery = img_dimy / 2;
2007
2008 results->limits = cx_string_new();
2009
2010
2011 /*
2012 * Processing
2013 */
2014
2015 if (remove_bias) {
2016 cpl_msg_info(fctid, "Bias correction will be done.");
2017 }
2018 else {
2019 cpl_msg_warning(fctid, "No Bias correction will be done!");
2020 }
2021
2022
2023 if (model == GIBIAS_MODEL_MEAN) {
2024
2025 /*
2026 * Compute bias average value over biaslimits areas
2027 */
2028
2029 cpl_msg_info(fctid, "Using bias model 'MEAN' ...");
2030
2031 if ((option == GIBIAS_OPTION_PLANE) ||
2032 (method == GIBIAS_METHOD_PLANE)) {
2033
2034 cpl_msg_error(fctid, "Can not use MEAN model with PLANE method.");
2035 return -3;
2036
2037 }
2038
2039 error_code = _giraffe_bias_compute_mean(results, img, biasareas,
2040 sigma, numiter, maxfraction);
2041
2042
2043 }
2044 else if ((method == GIBIAS_METHOD_CURVE) ||
2045 ((method != GIBIAS_METHOD_PROFILE) &&
2046 (option == GIBIAS_OPTION_CURVE))) {
2047
2048 /*
2049 * A 2D Bias Surface is fitted on a grid of the pixels in the
2050 * bias areas
2051 */
2052
2053 cpl_msg_info(fctid, "Using bias model 'CURVE' ...");
2054
2055 error_code = _giraffe_bias_compute_curve(results, img, biasareas,
2056 sigma, numiter, maxfraction, xdeg, ydeg, xstep, ystep);
2057
2058 if (error_code != EXIT_SUCCESS) {
2059
2060 cpl_msg_error(fctid, "Error during calculation of bias curve, "
2061 "aborting...");
2062 return error_code;
2063
2064 }
2065
2066 coeff = results->coeffs;
2067
2068 }
2069 else {
2070
2071 /*
2072 * A 2D Bias Plane is fitted on the pixels of the biaslimits area
2073 */
2074
2075 cpl_msg_info(fctid, "Using bias model 'PLANE (FITTED)' ...");
2076
2077 error_code = _giraffe_bias_compute_plane(results, img, biasareas,
2078 sigma, numiter, maxfraction);
2079
2080 if (error_code == EXIT_SUCCESS) {
2081
2082 coeff = results->coeffs;
2083
2084 results->mean = cpl_matrix_get(coeff, 0, 0) +
2085 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2086 cpl_matrix_get(coeff, 0, 2) * img_centery;
2087
2088 }
2089 else {
2090
2091 cpl_msg_error(fctid, "Error during calculation of bias plane, "
2092 "aborting...");
2093 return error_code;
2094
2095 }
2096
2097 }
2098
2099
2100 /*
2101 * Perform computation based on method....
2102 */
2103
2104 s = cx_string_new();
2105 _giraffe_method_string(s, method, option);
2106
2107 cpl_msg_info(fctid, "Using bias method '%s'", cx_string_get(s));
2108
2109 cx_string_delete(s);
2110 s = NULL;
2111
2112
2113 switch (method) {
2114 case GIBIAS_METHOD_UNIFORM:
2115 {
2116
2117 /*
2118 * Subtract the average bias value
2119 */
2120
2121 if (model == GIBIAS_MODEL_MEAN) {
2122
2123 cpl_msg_info(fctid, "bias value (areas mean) = %f, "
2124 "bias sigma = %f", results->mean, results->sigma);
2125
2126 }
2127 else {
2128
2129 cpl_msg_info(fctid, "bias value at (%d, %d) = %f, "
2130 "bias sigma = %f", img_centerx, img_centery,
2131 results->mean, results->sigma);
2132
2133 }
2134
2135 if (remove_bias == TRUE) {
2136
2137 cpl_image_subtract_scalar(img, results->mean);
2138
2139 }
2140
2141 break;
2142 }
2143
2144 case GIBIAS_METHOD_PLANE:
2145 {
2146
2147 if (coeff == NULL) {
2148
2149 error_code = _giraffe_bias_compute_plane(results, img,
2150 biasareas, sigma, numiter, maxfraction);
2151
2152 if (error_code == EXIT_SUCCESS) {
2153
2154 coeff = results->coeffs;
2155
2156 results->mean = cpl_matrix_get(coeff, 0, 0) +
2157 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2158 cpl_matrix_get(coeff, 0, 2) * img_centery;
2159
2160 }
2161 else {
2162
2163 cpl_msg_error(fctid, "Error during calculation of bias "
2164 "plane, aborting...");
2165 return error_code;
2166
2167 }
2168
2169 }
2170
2171 cpl_msg_info(fctid, "bias value at (%d, %d) = %f, bias sigma = %f, "
2172 "bias plane = %g + %g * x + %g * y", img_centerx,
2173 img_centery, results->mean, results->sigma,
2174 cpl_matrix_get(coeff, 0, 0),
2175 cpl_matrix_get(coeff, 0, 1),
2176 cpl_matrix_get(coeff, 0, 2));
2177
2178 if (remove_bias == TRUE) {
2179
2180 //TODO Move plane removal to gir_image...
2181
2182 cxdouble* _img = cpl_image_get_data_double(img);
2183
2184#if 0
2185 cpl_image* bsmodel = cpl_image_new(img_dimx, img_dimy,
2186 CPL_TYPE_DOUBLE);
2187 cxdouble* _bsmodel = cpl_image_get_data_double(bsmodel);
2188#endif
2189
2190
2191 for (j = 0; j < img_dimy; j++) {
2192
2193 cxsize jj = j * img_dimx;
2194
2195 for (i = 0; i < img_dimx; i++) {
2196
2197#if 0
2198 _bsmodel[jj + i] = cpl_matrix_get(coeff, 0, 0)
2199 + cpl_matrix_get(coeff, 0, 1) * i
2200 + cpl_matrix_get(coeff, 0, 2) * j;
2201#endif
2202
2203 _img[jj + i] -= cpl_matrix_get(coeff, 0, 0)
2204 + cpl_matrix_get(coeff, 0, 1) * i
2205 + cpl_matrix_get(coeff, 0, 2) * j;
2206
2207 }
2208
2209 }
2210
2211#if 0
2212 cpl_image_save(bsmodel, "idiot.fits", CPL_BPP_IEEE_FLOAT,
2213 0, CPL_IO_DEFAULT);
2214 cpl_image_delete(bsmodel);
2215#endif
2216
2217 }
2218
2219 break;
2220
2221 }
2222
2223 case GIBIAS_METHOD_CURVE:
2224 {
2225
2226 if (coeff == NULL) {
2227
2228 error_code =
2229 _giraffe_bias_compute_curve(results, img, biasareas,
2230 sigma, numiter, maxfraction,
2231 xdeg, ydeg, xstep, ystep);
2232
2233 if (error_code != EXIT_SUCCESS) {
2234
2235 cpl_msg_error(fctid, "Error during calculation of bias "
2236 "surface curve, aborting...");
2237 return error_code;
2238
2239 }
2240
2241 coeff = results->coeffs;
2242
2243 }
2244
2245 cpl_msg_info(fctid, "bias value %f, bias sigma = %f\n",
2246 results->mean, results->sigma);
2247
2248
2249 if (remove_bias == TRUE) {
2250
2251 cpl_image* bsmodel = NULL;
2252
2253 cpl_matrix* x = cpl_matrix_new(img_dimx * img_dimy, 1);
2254 cpl_matrix* y = cpl_matrix_new(img_dimx * img_dimy, 1);
2255 cpl_matrix* fit = NULL;
2256
2257
2258 for (j = 0; j < img_dimy; ++j) {
2259
2260 register cxsize jj = j * img_dimx;
2261
2262 for (i = 0; i < img_dimx; ++i) {
2263
2264 cpl_matrix_set(x, jj + i, 0, i);
2265 cpl_matrix_set(y, jj + i, 0, j);
2266
2267 }
2268
2269 }
2270
2271 fit = giraffe_chebyshev_fit2d(0., 0., img_dimx, img_dimy,
2272 coeff, x, y);
2273
2274 cpl_matrix_delete(y);
2275 y = NULL;
2276
2277 cpl_matrix_delete(x);
2278 x = NULL;
2279
2280 bsmodel = cpl_image_wrap_double(img_dimx, img_dimy,
2281 cpl_matrix_get_data(fit));
2282
2283#if 0
2284 cpl_image_save(bsmodel, "idiot.fits", CPL_BPP_IEEE_FLOAT,
2285 0, CPL_IO_DEFAULT);
2286#endif
2287
2288 cpl_image_subtract(img, bsmodel);
2289
2290 cpl_image_unwrap(bsmodel);
2291 bsmodel = NULL;
2292
2293 cpl_matrix_delete(fit);
2294 fit = NULL;
2295
2296 }
2297
2298 break;
2299
2300 }
2301
2302 case GIBIAS_METHOD_PROFILE:
2303 {
2304
2305 error_code = _giraffe_bias_compute_profile(results, img,
2306 biasareas, 'y');
2307
2308 if (error_code != EXIT_SUCCESS) {
2309
2310 cpl_msg_error(fctid, "Error computing the bias area(s) "
2311 "profile");
2312 return error_code;
2313
2314 }
2315
2316 if (option == GIBIAS_OPTION_CURVE) {
2317
2318 error_code = _giraffe_bias_fit_profile(results, results->model,
2319 ydeg - 1);
2320 if (error_code != EXIT_SUCCESS) {
2321
2322 cpl_msg_error(fctid, "Error fitting the bias profile");
2323 return error_code;
2324
2325 }
2326
2327 }
2328
2329 if (remove_bias == TRUE) {
2330
2331 const cxdouble* _bias =
2332 cpl_image_get_data_double(results->model);
2333
2334 cxdouble* _img = cpl_image_get_data_double(img);
2335
2336
2337 cx_assert(_bias != NULL);
2338 cx_assert(_img != NULL);
2339
2340 for (j = 0; j < img_dimy; ++j) {
2341
2342 cxsize jj = j * img_dimx;
2343
2344
2345 for (i = 0; i < img_dimx; ++i) {
2346 _img[jj + i] -= _bias[j];
2347 }
2348
2349 }
2350
2351 }
2352
2353 break;
2354
2355 }
2356
2357 case GIBIAS_METHOD_MASTER:
2358 case GIBIAS_METHOD_ZMASTER:
2359 {
2360
2361 cxdouble biasdrift = 0.;
2362
2363 if (master_bias == NULL) {
2364
2365 cpl_msg_error(fctid, "Master Bias Frame required with MASTER "
2366 "method.");
2367 return -5;
2368
2369 }
2370
2371 master_bias_dimx = cpl_image_get_size_x(master_bias);
2372 master_bias_dimy = cpl_image_get_size_y(master_bias);
2373
2374 if ((master_bias_dimx != img_dimx) ||
2375 (master_bias_dimy != img_dimy)) {
2376
2377 cpl_msg_error(fctid, "Invalid master bias! Size should "
2378 "be [%d, %d] but is [%d, %d].", img_dimx, img_dimy,
2379 master_bias_dimx, master_bias_dimy);
2380 return -7;
2381
2382 }
2383
2384 if (coeff == NULL) {
2385
2386 if (option == GIBIAS_OPTION_CURVE) {
2387
2388 error_code = _giraffe_bias_compute_curve(results, img,
2389 biasareas, sigma, numiter, maxfraction, xdeg,
2390 ydeg, xstep, ystep);
2391
2392 if (error_code != EXIT_SUCCESS) {
2393
2394 cpl_msg_error(fctid, "Error during calculation of "
2395 "bias surface curve, aborting...");
2396 return error_code;
2397
2398 }
2399
2400 coeff = results->coeffs;
2401
2402 }
2403 else {
2404
2405 /*
2406 * Default to plane...
2407 */
2408
2409 error_code = _giraffe_bias_compute_plane(results, img,
2410 biasareas, sigma, numiter, maxfraction);
2411
2412 if (error_code == EXIT_SUCCESS) {
2413
2414 coeff = results->coeffs;
2415
2416 results->mean = cpl_matrix_get(coeff, 0, 0) +
2417 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2418 cpl_matrix_get(coeff, 0, 2) * img_centery;
2419
2420 }
2421 else {
2422
2423 cpl_msg_error(fctid, "Error during calculation of "
2424 "bias surface plane, aborting...");
2425 return error_code;
2426
2427 }
2428
2429 }
2430
2431 }
2432
2433 if (method == GIBIAS_METHOD_ZMASTER) {
2434
2435 /*
2436 * A possible drift of the average bias is compensated using
2437 * the pre/overscans reference value as an additive
2438 * correction
2439 */
2440
2441 biasdrift = results->mean - mbias;
2442
2443 cpl_msg_info(fctid, "Using bias drift value %.4g", biasdrift);
2444
2445 }
2446
2447
2448 /*
2449 * Write some info into log file...
2450 */
2451
2452 if (option == GIBIAS_OPTION_CURVE) {
2453
2454 //TODO dump coeff to log
2455
2456 cpl_msg_info(fctid, "Computed bias mean = %.4f, bias "
2457 "sigma = %.4e", results->mean, results->sigma);
2458
2459 }
2460 else {
2461
2462 if (option == GIBIAS_OPTION_PLANE) {
2463
2464 cpl_msg_info(fctid, "Bias plane = %.4e + "
2465 "%.4e * x + %.4e * y",
2466 cpl_matrix_get(coeff, 0, 0),
2467 cpl_matrix_get(coeff, 0, 1),
2468 cpl_matrix_get(coeff, 0, 2));
2469 cpl_msg_info(fctid, "Computed bias mean = %.4f, "
2470 "bias sigma = %.4e at (%d, %d)",
2471 results->mean, results->sigma,
2472 img_centerx, img_centery);
2473
2474 }
2475 else {
2476
2477 cpl_msg_info(fctid, "Computed bias mean = %.4f, bias "
2478 "sigma = %.4e",
2479 results->mean, results->sigma);
2480
2481 }
2482
2483 }
2484
2485
2486 /*
2487 * Remove bias if requested
2488 */
2489
2490 if (remove_bias == TRUE) {
2491
2492
2493 /*
2494 * No bad pixel map was given...
2495 * Subtract master bias pixels without modification except
2496 * pre/overscan trimming.
2497 */
2498
2499 /*
2500 * Bad pixel map was given...
2501 * For all pixels not flagged in bad_pixel_map the
2502 * master_bias pixel values are subtracted, otherwise
2503 * results->mean or bias model image pixels are subtracted
2504 * depending on model.
2505 */
2506
2507 if (bad_pixel_map == NULL) {
2508
2509 cpl_image_subtract(img, master_bias);
2510
2511 if (biasdrift != 0.) {
2512 cpl_image_subtract_scalar(img, biasdrift);
2513 }
2514
2515 }
2516 else {
2517
2518 bad_pixel_dimx = cpl_image_get_size_x(bad_pixel_map);
2519 bad_pixel_dimy = cpl_image_get_size_y(bad_pixel_map);
2520
2521 if ((bad_pixel_dimx != img_dimx) ||
2522 (bad_pixel_dimy != img_dimy)) {
2523
2524 cpl_msg_error(fctid, "Invalid bad pixel map size "
2525 "should be [%d, %d] but is [%d, %d].",
2526 img_dimx, img_dimy,
2527 bad_pixel_dimx, bad_pixel_dimy);
2528 return -6;
2529
2530 }
2531
2532 if (option == GIBIAS_OPTION_CURVE) {
2533
2534 const cxint* _bpix =
2535 cpl_image_get_data_int_const(bad_pixel_map);
2536
2537 const cxdouble* _mbias =
2538 cpl_image_get_data_double_const(master_bias);
2539
2540 cxdouble* _img = cpl_image_get_data_double(img);
2541
2542 cpl_matrix* x =
2543 cpl_matrix_new(img_dimx * img_dimy, 1);
2544 cpl_matrix* y =
2545 cpl_matrix_new(img_dimx * img_dimy, 1);
2546 cpl_matrix* fit = NULL;
2547
2548
2549 for (j = 0; j < img_dimy; ++j) {
2550
2551 register cxsize jj = j * img_dimx;
2552
2553 for (i = 0; i < img_dimx; ++i) {
2554
2555 cpl_matrix_set(x, jj + i, 0, i);
2556 cpl_matrix_set(y, jj + i, 0, j);
2557
2558 }
2559 }
2560
2561 fit = giraffe_chebyshev_fit2d(0., 0.,
2562 img_dimx, img_dimy,
2563 coeff, x, y);
2564
2565 cpl_matrix_delete(y);
2566 y = NULL;
2567
2568 cpl_matrix_delete(x);
2569 x = NULL;
2570
2571
2572 for (j = 0; j < img_dimy; ++j) {
2573
2574 register cxsize jj = j * img_dimx;
2575
2576 for (i = 0; i < img_dimx; ++i) {
2577
2578 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2579
2580 _img[jj + i] -= (_mbias[jj + i] +
2581 biasdrift);
2582
2583 }
2584 else {
2585 _img[jj + i] -= cpl_matrix_get(fit, i, j);
2586 }
2587
2588 }
2589
2590 }
2591
2592 cpl_matrix_delete(fit);
2593 fit = NULL;
2594
2595 }
2596 else if (option == GIBIAS_OPTION_PLANE) {
2597
2598 const cxint* _bpix =
2599 cpl_image_get_data_int_const(bad_pixel_map);
2600
2601 const cxdouble* _mbias =
2602 cpl_image_get_data_double_const(master_bias);
2603
2604 cxdouble* _img = cpl_image_get_data_double(img);
2605
2606
2607 for (j = 0; j < img_dimy; j++) {
2608
2609 cxsize jj = j * img_dimx;
2610
2611 for (i = 0; i < img_dimx; i++) {
2612
2613 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2614
2615 _img[jj + i] -= (_mbias[jj + i] +
2616 biasdrift);
2617
2618 }
2619 else {
2620
2621 _img[jj + i] -=
2622 cpl_matrix_get(coeff, 0, 0) +
2623 cpl_matrix_get(coeff, 0, 1) * j +
2624 cpl_matrix_get(coeff, 0, 2) * i;
2625 }
2626
2627 }
2628
2629 }
2630
2631 }
2632 else {
2633
2634 const cxint* _bpix =
2635 cpl_image_get_data_int_const(bad_pixel_map);
2636
2637 const cxdouble* _mbias =
2638 cpl_image_get_data_double_const(master_bias);
2639
2640 cxdouble* _img = cpl_image_get_data_double(img);
2641
2642
2643 for (j = 0; j < img_dimy; j++) {
2644
2645 cxsize jj = j * img_dimx;
2646
2647 for (i = 0; i < img_dimx; i++) {
2648
2649 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2650
2651 _img[jj + i] -= (_mbias[jj + i] +
2652 biasdrift);
2653
2654 }
2655 else {
2656
2657 _img[jj + i] -= results->mean;
2658
2659 }
2660
2661 }
2662
2663 }
2664
2665 }
2666
2667 }
2668
2669 }
2670
2671 break;
2672
2673 }
2674
2675 default:
2676 {
2677
2678 cpl_msg_error(fctid, "Invalid method, aborting...");
2679 return -4;
2680 break;
2681
2682 }
2683
2684 }
2685
2686
2687 cpl_msg_info(fctid, "Resulting biaslimits : %s",
2688 cx_string_get(results->limits));
2689
2690 return 0;
2691
2692}
2693
2694
2712cpl_matrix *
2713giraffe_get_raw_areas(const GiImage *image)
2714{
2715
2716 const cxchar *const _id = "giraffe_get_raw_areas";
2717
2718 cxint32 prescx = 0;
2719 cxint32 prescy = 0;
2720 cxint32 ovrscx = 0;
2721 cxint32 ovrscy = 0;
2722
2723 cxint32 winx = 0;
2724 cxint32 winy = 0;
2725
2726 cxint32 tprescx = 0;
2727 cxint32 tprescy = 0;
2728 cxint32 tovrscx = 0;
2729 cxint32 tovrscy = 0;
2730
2731 cxint32 row = 0;
2732
2733 cpl_matrix *m_tmp;
2734 cpl_propertylist *properties;
2735
2736
2737 properties = giraffe_image_get_properties(image);
2738 if (!properties) {
2739 cpl_error_set(_id, CPL_ERROR_NULL_INPUT);
2740 return NULL;
2741 }
2742
2743 winx = cpl_propertylist_get_int(properties, GIALIAS_WINX);
2744 winy = cpl_propertylist_get_int(properties, GIALIAS_WINY);
2745
2746 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
2747 tprescx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
2748 if (tprescx > 0) {
2749 prescx = tprescx;
2750 }
2751 }
2752 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
2753 tprescy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
2754 if (tprescy > 0) {
2755 prescy = tprescy;
2756 }
2757 }
2758 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
2759 tovrscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
2760 if (tovrscx > 0) {
2761 ovrscx = tovrscx;
2762 }
2763 }
2764 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
2765 tovrscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
2766 if (tovrscy > 0) {
2767 ovrscy = tovrscy;
2768 }
2769 }
2770
2771
2772 m_tmp = cpl_matrix_new(1, 4);
2773
2774 if (prescx > 0) {
2775
2776 cpl_matrix_set(m_tmp, row, 0, 0.0);
2777 cpl_matrix_set(m_tmp, row, 1, (cxdouble) prescx - 1);
2778 cpl_matrix_set(m_tmp, row, 2, 0.0);
2779 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2780
2781 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2782 row++;
2783 }
2784
2785 if (ovrscx > 0) {
2786
2787 cpl_matrix_set(m_tmp, row, 0, (cxdouble) winx - ovrscx);
2788 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2789 cpl_matrix_set(m_tmp, row, 2, 0.0);
2790 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2791
2792 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2793 row++;
2794 }
2795
2796 if (!(prescy == 0 || prescx == 0 || ovrscx == 0)) {
2797
2798 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2799 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2800 cpl_matrix_set(m_tmp, row, 2, 0.0);
2801 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2802
2803 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2804 row++;
2805
2806 }
2807 else if (!(prescy == 0 || prescx == 0)) {
2808
2809 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2810 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2811 cpl_matrix_set(m_tmp, row, 2, 0.0);
2812 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2813
2814 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2815 row++;
2816
2817 }
2818 else if (!(prescy == 0 || ovrscx == 0)) {
2819
2820 cpl_matrix_set(m_tmp, row, 0, 0.0);
2821 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2822 cpl_matrix_set(m_tmp, row, 2, 0.0);
2823 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2824
2825 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2826 row++;
2827
2828 }
2829 else if (prescy > 0) {
2830
2831 cpl_matrix_set(m_tmp, row, 0, 0.0);
2832 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2833 cpl_matrix_set(m_tmp, row, 2, 0.0);
2834 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2835
2836 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2837 row++;
2838 }
2839
2840 if (!(ovrscy == 0 || prescx == 0 || ovrscx == 0)) {
2841
2842 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2843 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2844 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2845 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2846
2847 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2848 row++;
2849
2850 }
2851 else if (!(ovrscy == 0 || prescx == 0)) {
2852
2853 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2854 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2855 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2856 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2857
2858 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2859 row++;
2860
2861 }
2862 else if (!(ovrscy == 0 || ovrscx == 0)) {
2863
2864 cpl_matrix_set(m_tmp, row, 0, 0.0);
2865 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2866 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2867 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2868
2869 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2870 row++;
2871
2872 }
2873 else if (ovrscy > 0) {
2874
2875 cpl_matrix_set(m_tmp, row, 0, 0.0);
2876 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2877 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2878 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2879
2880 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2881 row++;
2882 }
2883
2884 cpl_matrix_resize(m_tmp, 0, -1, 0, 0);
2885 row--;
2886
2887 if (row == 0) {
2888 cpl_matrix_delete(m_tmp);
2889 return NULL;
2890 }
2891
2892 return m_tmp;
2893}
2894
2895
2914cxint
2916{
2917
2918 const cxchar *fctid = "giraffe_trim_raw_areas";
2919
2920 cxint nx, ny;
2921 cxint newlx, newly, newux, newuy;
2922
2923 cxint ovscx = 0;
2924 cxint ovscy = 0;
2925 cxint prscx = 0;
2926 cxint prscy = 0;
2927
2928 cpl_propertylist *properties = giraffe_image_get_properties(image);
2929
2930 cpl_image *_image = giraffe_image_get(image);
2931 cpl_image *tmp;
2932
2933
2934 if (!properties) {
2935 cpl_msg_error(fctid, "Image does not contain any properties!");
2936 return 1;
2937 }
2938
2939 nx = cpl_image_get_size_x(_image);
2940 ny = cpl_image_get_size_y(_image);
2941
2942 if (cpl_propertylist_has(properties, GIALIAS_NAXIS1)) {
2943 cxint _nx = cpl_propertylist_get_int(properties, GIALIAS_NAXIS1);
2944
2945 if (_nx != nx) {
2946 cpl_msg_warning(fctid, "Image size (%d) and image property `%s' "
2947 "(%d) are not consistent! Using image size ...",
2948 nx, GIALIAS_NAXIS1, _nx);
2949 }
2950 }
2951 else {
2952 cpl_msg_warning(fctid, "Image does not contain any property `%s'. "
2953 "Using image size (%d)", GIALIAS_NAXIS1, nx);
2954 }
2955
2956
2957 if (cpl_propertylist_has(properties, GIALIAS_NAXIS2)) {
2958 cxint _ny = cpl_propertylist_get_int(properties, GIALIAS_NAXIS2);
2959
2960 if (_ny != ny) {
2961 cpl_msg_warning(fctid, "Image size (%d) and image property `%s' "
2962 "(%d) are not consistent! Using image size ...",
2963 ny, GIALIAS_NAXIS2, _ny);
2964 }
2965 }
2966 else {
2967 cpl_msg_warning(fctid, "Image does not contain any property `%s'. "
2968 "Using image size (%d)", GIALIAS_NAXIS2, ny);
2969 }
2970
2971 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
2972 ovscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
2973 }
2974
2975 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
2976 ovscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
2977 }
2978
2979 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
2980 prscx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
2981 }
2982
2983 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
2984 prscy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
2985 }
2986
2987 // FIXME:
2988 // Check that the one pixel offset is ok. cpl_image_extract starts
2989 // counting from 1,1 with right and to boundaries inclusive.
2990
2991 newlx = prscx + 1;
2992 newly = prscy + 1;
2993 newux = nx - ovscx;
2994 newuy = ny - ovscy;
2995
2996 tmp = cpl_image_extract(_image, newlx, newly, newux, newuy);
2997
2998 giraffe_image_set(image, tmp);
2999 cpl_image_delete(tmp);
3000
3001 _image = giraffe_image_get(image);
3002
3003 nx = cpl_image_get_size_x(_image);
3004 ny = cpl_image_get_size_y(_image);
3005
3006
3007 /*
3008 * Update image properties
3009 */
3010
3011 cpl_propertylist_set_int(properties, GIALIAS_NAXIS1, nx);
3012 cpl_propertylist_set_int(properties, GIALIAS_NAXIS2, ny);
3013
3014 if (cpl_propertylist_has(properties, GIALIAS_CRPIX1)) {
3015 cxdouble crpix1 = cpl_propertylist_get_double(properties,
3016 GIALIAS_CRPIX1);
3017
3018 /*
3019 * For GIRAFFE the reference pixel is defined as
3020 * 1 - width(prescan)
3021 */
3022
3023 crpix1 += (cxdouble)prscx;
3024 cpl_propertylist_set_double(properties, GIALIAS_CRPIX1, crpix1);
3025 }
3026
3027 if (cpl_propertylist_has(properties, GIALIAS_CRPIX2)) {
3028 cxdouble crpix2 = cpl_propertylist_get_double(properties,
3029 GIALIAS_CRPIX2);
3030
3031 crpix2 -= (cxdouble) prscy;
3032 cpl_propertylist_set_double(properties, GIALIAS_CRPIX2, crpix2);
3033 }
3034
3035 cpl_propertylist_erase(properties, GIALIAS_OVSCX);
3036 cpl_propertylist_erase(properties, GIALIAS_OVSCY);
3037 cpl_propertylist_erase(properties, GIALIAS_PRSCX);
3038 cpl_propertylist_erase(properties, GIALIAS_PRSCY);
3039
3040 return 0;
3041
3042}
3043
3044
3106cxint
3107giraffe_bias_remove(GiImage* result, const GiImage* raw,
3108 const GiImage* master_bias, const GiImage* bad_pixels,
3109 const cpl_matrix* biaslimits, const GiBiasConfig* config)
3110{
3111
3112 const cxchar* const fctid = "giraffe_bias_remove";
3113
3114 cxint error_code;
3115
3116 cxdouble bias_value = 0.;
3117
3118 cx_string* s;
3119
3120 cpl_matrix* areas = NULL;
3121
3122 cpl_image* _raw = giraffe_image_get(raw);
3123 cpl_image* _master_bias = giraffe_image_get(master_bias);
3124 cpl_image* _bad_pixels = giraffe_image_get(bad_pixels);
3125 cpl_image* timage;
3126
3127 cpl_propertylist* properties;
3128
3129 GiBiasResults bias_results = {0., 0., 0., NULL, NULL, NULL};
3130
3131
3132 cx_assert(raw != NULL);
3133 cx_assert(config != NULL);
3134 cx_assert(biaslimits == NULL);
3135
3136 if (result == NULL) {
3137 return 1;
3138 }
3139
3140
3141 /*
3142 * Setup user defined areas to use for the bias computation
3143 */
3144
3145 if (strncmp(config->areas, "None", 4) == 0) {
3146
3147 cpl_msg_info(fctid, "No bias areas specified, using pre/overscan"
3148 "regions of the raw frame.");
3149
3150 cpl_error_reset();
3151 areas = giraffe_get_raw_areas(raw);
3152 if (areas == NULL) {
3153 if (cpl_error_get_code() == CPL_ERROR_NULL_INPUT) {
3154 cpl_msg_error(fctid, "Invalid raw image! Image does not "
3155 "contain any properties!");
3156 }
3157 else {
3158 cpl_msg_error(fctid, "Invalid raw image! Image does not "
3159 "contain or has invalid pre- and overscan "
3160 "properties.");
3161 }
3162
3163 return 1;
3164 }
3165
3166 }
3167 else {
3168
3169 areas = _giraffe_bias_get_areas(config->areas);
3170
3171 if (areas == NULL) {
3172
3173 cpl_msg_error(fctid, "Invalid bias area specification '%s'!",
3174 config->areas);
3175 return 1;
3176
3177 }
3178
3179 cpl_msg_info(fctid, "Using bias area(s) '%s' for bias computation",
3180 config->areas);
3181
3182 }
3183
3184
3185 /*
3186 * Processing
3187 */
3188
3189
3190 /*
3191 * Check whether a masterbias image was specified
3192 * If so, compare pre/ovrscan keywords/areas
3193 */
3194
3195 if (master_bias != NULL) {
3196 cpl_propertylist *_properties;
3197 cxbool is_same_size = FALSE;
3198
3199 is_same_size = _giraffe_compare_overscans(raw, master_bias);
3200
3201 // FIXME:
3202 // Fully processed calibrations usually do not have overscans
3203 // any longer. Instead of failing at this point some dummies
3204 // could be created.
3205
3206 if (is_same_size==FALSE) {
3207 cpl_msg_info(fctid, "PRE/OVRSCAN Regions do not match between "
3208 "RAW Image and Masterbias Image.");
3209
3210 return 1;
3211 }
3212
3213 _properties = giraffe_image_get_properties(master_bias);
3214
3215 if (cpl_propertylist_has(_properties, GIALIAS_BIASVALUE)) {
3216 bias_value = cpl_propertylist_get_double(_properties,
3217 GIALIAS_BIASVALUE);
3218 }
3219 }
3220
3221
3222 /*
3223 * Check whether a Bad Pixel Map image was specified
3224 * If so, compare pre/ovrscan keywords/areas
3225 */
3226
3227 if (bad_pixels != NULL) {
3228 cxbool is_same_size = FALSE;
3229
3230 is_same_size = _giraffe_compare_overscans(raw, bad_pixels);
3231
3232 // FIXME:
3233 // Fully processed calibrations usually do not have overscans
3234 // any longer. Instead of failing at this point some dummies
3235 // could be created.
3236
3237 if (is_same_size == FALSE) {
3238 cpl_msg_info(fctid, "PRE/OVRSCAN Regions do not match between "
3239 "RAW Image and Bad Pixel Image.");
3240
3241 return 1;
3242 }
3243 }
3244
3245
3246 /*
3247 * Create a copy of the raw frame. We are not working on the raw
3248 * frame directly.
3249 */
3250
3251 // FIXME:
3252 // Can we use the raw frame directly? (RP)
3253
3254 timage = cpl_image_duplicate(_raw);
3255
3256
3257 /*
3258 * Remove Bias Calculation. The bias corrected image will be
3259 * stored in timage.
3260 */
3261
3262 error_code = _giraffe_bias_compute(&bias_results, timage,
3263 areas, _master_bias,
3264 _bad_pixels, config->method,
3265 config->option, config->model,
3266 config->remove, bias_value,
3267 config->xdeg, config->ydeg,
3268 config->xstep, config->ystep,
3269 config->sigma, config->iterations,
3270 config->fraction);
3271
3272 cpl_matrix_delete(areas);
3273 areas = NULL;
3274
3275
3276 // FIXME:
3277 // Possibly check returned error code and do a dedicated
3278 // exception handling (RP)
3279
3280 if (error_code) {
3281
3282 _giraffe_biasresults_clear(&bias_results);
3283
3284 cpl_msg_error(fctid, "Bias computation failed with error %d!",
3285 error_code);
3286 return 1;
3287
3288 }
3289
3290
3291 /*
3292 * Should we remove Bias or not?
3293 */
3294
3295 properties = giraffe_image_get_properties(raw);
3296
3297 if (config->remove == TRUE) {
3298
3299 giraffe_image_set(result, timage);
3300 cpl_image_delete(timage);
3301
3302 giraffe_image_set_properties(result, properties);
3303
3304 }
3305 else {
3306
3307 cpl_image_delete(timage);
3308
3309 giraffe_image_set(result, _raw);
3310 giraffe_image_set_properties(result, properties);
3311
3312 }
3313
3314
3315 /*
3316 * Postprocessing
3317 */
3318
3319 properties = giraffe_image_get_properties(result);
3320
3321 if (config->remove == TRUE) {
3322
3323 cpl_propertylist_set_int(properties, GIALIAS_BITPIX, -32);
3324 cpl_propertylist_set_double(properties, GIALIAS_BZERO, 0.);
3325 cpl_propertylist_set_double(properties, GIALIAS_BSCALE, 1.);
3326
3327 if (cpl_propertylist_has(properties, GIALIAS_GIRFTYPE)) {
3328 cpl_propertylist_set_string(properties,
3329 GIALIAS_GIRFTYPE, "BRMIMG");
3330 }
3331 else {
3332 cpl_propertylist_append_string(properties,
3333 GIALIAS_GIRFTYPE, "BRMIMG");
3334 }
3335 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3336 "GIRAFFE bias subtracted frame");
3337
3338 /*
3339 * Remove pre- and overscans
3340 */
3341
3342 giraffe_trim_raw_areas(result);
3343
3344
3345 /*
3346 * Convert bias corrected image to electrons
3347 */
3348
3349 // FIXME:
3350 // Is this really needed? Disabled on request of DFO (RP)
3351 // If it has to be put back use giraffe_propertylist_get_conad()!
3352#if 0
3353 if (cpl_propertylist_has(properties, GIALIAS_CONAD)) {
3354 cxdouble conad = cpl_propertylist_get_double(properties, GIALIAS_CONAD);
3355
3356 if (conad > 0.) {
3357
3358 cpl_image* _result = giraffe_image_get(result);
3359
3360 cpl_msg_info(fctid, "Converting bias subtracted frame to "
3361 "electrons; conversion factor is %.2f e-/ADU",
3362 conad);
3363 cpl_image_multiply_scalar(_result, conad);
3364 }
3365 else {
3366 cpl_msg_warning(fctid, "Invalid conversion factor %.2f "
3367 "e-/ADU! Bias subtracted image will not be "
3368 "converted.", conad);
3369 }
3370 }
3371 else {
3372 cpl_msg_info(fctid, "Conversion factor not found. Bias subtracted "
3373 "image will remain in ADU");
3374 }
3375#endif
3376 }
3377
3378 s = cx_string_new();
3379
3380 _giraffe_method_string(s, config->method, config->option);
3381 cpl_propertylist_append_string(properties, GIALIAS_BIASMETHOD, cx_string_get(s));
3382
3383 cx_string_delete(s);
3384
3385 cpl_propertylist_append_double(properties, GIALIAS_BIASVALUE,
3386 bias_results.mean);
3387 cpl_propertylist_append_double(properties, GIALIAS_BIASERROR,
3388 bias_results.sigma);
3389 cpl_propertylist_append_double(properties, GIALIAS_BIASSIGMA,
3390 bias_results.rms);
3391
3392 if (bias_results.coeffs) {
3393 cpl_propertylist_append_double(properties, GIALIAS_BCLIPSIGMA,
3394 config->sigma);
3395 cpl_propertylist_append_int(properties, GIALIAS_BCLIPNITER,
3396 config->iterations);
3397 cpl_propertylist_append_double(properties, GIALIAS_BCLIPMFRAC,
3398 config->fraction);
3399 }
3400
3401 if (bias_results.limits) {
3402 cpl_propertylist_append_string(properties, GIALIAS_BIASAREAS,
3403 cx_string_get(bias_results.limits));
3404 }
3405
3406 if (bias_results.coeffs) {
3407 s = cx_string_new();
3408
3409 _giraffe_stringify_coefficients(s, bias_results.coeffs);
3410 cpl_propertylist_append_string(properties, GIALIAS_BIASPLANE,
3411 cx_string_get(s));
3412
3413 cx_string_delete(s);
3414 }
3415
3416
3417 /*
3418 * Cleanup
3419 */
3420
3421 _giraffe_biasresults_clear(&bias_results);
3422
3423 return 0;
3424
3425}
3426
3427
3438GiBiasConfig*
3439giraffe_bias_config_create(cpl_parameterlist* list)
3440{
3441
3442 const cxchar* s;
3443 cpl_parameter* p;
3444
3445 GiBiasConfig* config = NULL;
3446
3447
3448 if (!list) {
3449 return NULL;
3450 }
3451
3452 config = cx_calloc(1, sizeof *config);
3453
3454
3455 /*
3456 * Some defaults
3457 */
3458
3459 config->method = GIBIAS_METHOD_UNDEFINED;
3460 config->option = GIBIAS_OPTION_UNDEFINED;
3461 config->model = GIBIAS_MODEL_UNDEFINED;
3462 config->mbias = 0.;
3463 config->xdeg = 1;
3464 config->ydeg = 1;
3465
3466
3467 p = cpl_parameterlist_find(list, "giraffe.biasremoval.remove");
3468 config->remove = cpl_parameter_get_bool(p);
3469
3470 p = cpl_parameterlist_find(list, "giraffe.biasremoval.method");
3471 s = cpl_parameter_get_string(p);
3472
3473 if (!strcmp(s, "UNIFORM")) {
3474 config->method = GIBIAS_METHOD_UNIFORM;
3475 }
3476
3477 if (!strcmp(s, "PLANE")) {
3478 config->method = GIBIAS_METHOD_PLANE;
3479 }
3480
3481 if (!strcmp(s, "CURVE")) {
3482 config->method = GIBIAS_METHOD_CURVE;
3483 }
3484
3485 if (!strcmp(s, "PROFILE")) {
3486 config->method = GIBIAS_METHOD_PROFILE;
3487 }
3488
3489 if (!strcmp(s, "MASTER")) {
3490 config->method = GIBIAS_METHOD_MASTER;
3491 }
3492
3493 if (!strcmp(s, "ZMASTER")) {
3494 config->method = GIBIAS_METHOD_ZMASTER;
3495 }
3496
3497 if (!strcmp(s, "PROFILE+CURVE")) {
3498 config->method = GIBIAS_METHOD_PROFILE;
3499 config->option = GIBIAS_OPTION_CURVE;
3500 }
3501
3502 if (!strcmp(s, "MASTER+PLANE")) {
3503 config->method = GIBIAS_METHOD_MASTER;
3504 config->option = GIBIAS_OPTION_PLANE;
3505 }
3506
3507 if (!strcmp(s, "ZMASTER+PLANE")) {
3508 config->method = GIBIAS_METHOD_ZMASTER;
3509 config->option = GIBIAS_OPTION_PLANE;
3510 }
3511
3512 if (!strcmp(s, "MASTER+CURVE")) {
3513 config->method = GIBIAS_METHOD_MASTER;
3514 config->option = GIBIAS_OPTION_CURVE;
3515 }
3516
3517 if (!strcmp(s, "ZMASTER+CURVE")) {
3518 config->method = GIBIAS_METHOD_ZMASTER;
3519 config->option = GIBIAS_OPTION_CURVE;
3520 }
3521
3522 cx_assert(config->method != GIBIAS_METHOD_UNDEFINED);
3523
3524 p = cpl_parameterlist_find(list, "giraffe.biasremoval.areas");
3525 s = cpl_parameter_get_string(p);
3526 config->areas = cx_strdup(s);
3527
3528 p = cpl_parameterlist_find(list, "giraffe.biasremoval.sigma");
3529 config->sigma = cpl_parameter_get_double(p);
3530
3531 p = cpl_parameterlist_find(list, "giraffe.biasremoval.iterations");
3532 config->iterations = cpl_parameter_get_int(p);
3533
3534 p = cpl_parameterlist_find(list, "giraffe.biasremoval.fraction");
3535 config->fraction = cpl_parameter_get_double(p);
3536
3537 if (config->method == GIBIAS_METHOD_CURVE ||
3538 config->option == GIBIAS_OPTION_CURVE) {
3539 p = cpl_parameterlist_find(list, "giraffe.biasremoval.xorder");
3540 config->xdeg = 1 + cpl_parameter_get_int(p);
3541
3542 p = cpl_parameterlist_find(list, "giraffe.biasremoval.yorder");
3543 config->ydeg = 1 + cpl_parameter_get_int(p);
3544 }
3545
3546 p = cpl_parameterlist_find(list, "giraffe.biasremoval.xstep");
3547 config->xstep = cpl_parameter_get_int(p);
3548
3549 p = cpl_parameterlist_find(list, "giraffe.biasremoval.ystep");
3550 config->ystep = cpl_parameter_get_int(p);
3551
3552 return config;
3553
3554}
3555
3556
3569void
3570giraffe_bias_config_destroy(GiBiasConfig* config)
3571{
3572
3573 if (config) {
3574 if (config->areas) {
3575 cx_free(config->areas);
3576 config->areas = NULL;
3577 }
3578
3579 cx_free(config);
3580 }
3581
3582 return;
3583}
3584
3585
3597void
3598giraffe_bias_config_add(cpl_parameterlist* list)
3599{
3600
3601 cpl_parameter* p;
3602
3603
3604 if (!list) {
3605 return;
3606 }
3607
3608 p = cpl_parameter_new_value("giraffe.biasremoval.remove",
3609 CPL_TYPE_BOOL,
3610 "Enable bias removal",
3611 "giraffe.biasremoval",
3612 TRUE);
3613 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "remove-bias");
3614 cpl_parameterlist_append(list, p);
3615
3616
3617 p = cpl_parameter_new_enum("giraffe.biasremoval.method",
3618 CPL_TYPE_STRING,
3619 "Bias removal method",
3620 "giraffe.biasremoval",
3621 "PROFILE", 11, "UNIFORM", "PLANE", "CURVE",
3622 "PROFILE", "MASTER", "ZMASTER", "PROFILE+CURVE",
3623 "MASTER+PLANE", "MASTER+CURVE", "ZMASTER+PLANE",
3624 "ZMASTER+CURVE");
3625 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-method");
3626 cpl_parameterlist_append(list, p);
3627
3628
3629 p = cpl_parameter_new_value("giraffe.biasremoval.areas",
3630 CPL_TYPE_STRING,
3631 "Bias areas to use "
3632 "(Xl0:Xr0:Yl0:Yr0, ... ,Xln:Xrn:Yln:Yrn)",
3633 "giraffe.biasremoval",
3634 "5:40:0:4095");
3635 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-areas");
3636 cpl_parameterlist_append(list, p);
3637
3638
3639 p = cpl_parameter_new_value("giraffe.biasremoval.sigma",
3640 CPL_TYPE_DOUBLE,
3641 "Sigma Clipping: sigma threshold factor",
3642 "giraffe.biasremoval",
3643 2.5);
3644 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-sigma");
3645 cpl_parameterlist_append(list, p);
3646
3647
3648 p = cpl_parameter_new_value("giraffe.biasremoval.iterations",
3649 CPL_TYPE_INT,
3650 "Sigma Clipping: maximum number of "
3651 "iterations",
3652 "giraffe.biasremoval",
3653 5);
3654 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-niter");
3655 cpl_parameterlist_append(list, p);
3656
3657
3658 p = cpl_parameter_new_value("giraffe.biasremoval.fraction",
3659 CPL_TYPE_DOUBLE,
3660 "Sigma Clipping: minimum fraction of points "
3661 "accepted/total [0.0..1.0]",
3662 "giraffe.biasremoval",
3663 0.8);
3664 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-mfrac");
3665 cpl_parameterlist_append(list, p);
3666
3667
3668 p = cpl_parameter_new_value("giraffe.biasremoval.xorder",
3669 CPL_TYPE_INT,
3670 "Order of X polynomial fit (CURVE method "
3671 "only)",
3672 "giraffe.biasremoval",
3673 1);
3674 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-xorder");
3675 cpl_parameterlist_append(list, p);
3676
3677 p = cpl_parameter_new_value("giraffe.biasremoval.yorder",
3678 CPL_TYPE_INT,
3679 "Order of Y polynomial fit (CURVE method "
3680 "only)",
3681 "giraffe.biasremoval",
3682 1);
3683 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-yorder");
3684 cpl_parameterlist_append(list, p);
3685
3686
3687 p = cpl_parameter_new_value("giraffe.biasremoval.xstep",
3688 CPL_TYPE_INT,
3689 "Sampling step along X (CURVE method only)",
3690 "giraffe.biasremoval",
3691 1);
3692 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-xstep");
3693 cpl_parameterlist_append(list, p);
3694
3695
3696 p = cpl_parameter_new_value("giraffe.biasremoval.ystep",
3697 CPL_TYPE_INT,
3698 "Sampling step along Y (CURVE method only)",
3699 "giraffe.biasremoval",
3700 1);
3701 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-ystep");
3702 cpl_parameterlist_append(list, p);
3703
3704 return;
3705}
cxint giraffe_trim_raw_areas(GiImage *image)
Remove pre- and overscan ares from an image.
Definition: gibias.c:2915
GiBiasConfig * giraffe_bias_config_create(cpl_parameterlist *list)
Creates a setup structure for a bias removal task.
Definition: gibias.c:3439
void giraffe_bias_config_add(cpl_parameterlist *list)
Adds parameters for the bias removal.
Definition: gibias.c:3598
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.
Definition: gibias.c:3107
cpl_matrix * giraffe_get_raw_areas(const GiImage *image)
Create bias areas from an image.
Definition: gibias.c:2713
void giraffe_bias_config_destroy(GiBiasConfig *config)
Destroys a bias removal setup structure.
Definition: gibias.c:3570
cpl_image * giraffe_image_get(const GiImage *self)
Gets the image data.
Definition: giimage.c:218
cpl_propertylist * giraffe_image_get_properties(const GiImage *self)
Get the properties of an image.
Definition: giimage.c:282
cxint giraffe_image_set(GiImage *self, cpl_image *image)
Sets the image data.
Definition: giimage.c:244
cxint giraffe_image_set_properties(GiImage *self, cpl_propertylist *properties)
Attaches a property list to an image.
Definition: giimage.c:312
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.
Definition: gimatrix.c:579
cxdouble giraffe_matrix_sigma_fit(const cpl_matrix *matrix, const cpl_matrix *matrix_fit)
Compute sigma of matrix fit.
Definition: gimatrix.c:275
cpl_matrix * giraffe_matrix_leastsq(const cpl_matrix *mA, const cpl_matrix *mB)
Computes the solution of an equation using a pseudo-inverse.
Definition: gimatrix.c:503
cxdouble giraffe_matrix_sigma_mean(const cpl_matrix *matrix, cxdouble mean)
Compute sigma of matrix elements, with a given mean value.
Definition: gimatrix.c:229

This file is part of the GIRAFFE Pipeline Reference Manual 2.19.4.
Documentation copyright © 2002-2006 European Southern Observatory.
Generated on Fri Feb 6 2026 11:30:08 by doxygen 1.9.6 written by Dimitri van Heesch, © 1997-2004