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 ntotal = 0;
1478 cxsize nareas = 0;
1479 cxsize nvalid = 0;
1480
1481 cxdouble rms = 0.;
1482 cxdouble sigma = 0.;
1483
1484 cpl_matrix* _areas = NULL;
1485
1486 cpl_image* profile = NULL;
1487 cpl_image* model = NULL;
1488
1489
1490 cx_assert(results != NULL);
1491 cx_assert(results->limits != NULL);
1492 cx_assert(results->model == NULL);
1493
1494 cx_assert(areas != NULL);
1495
1496 cx_assert(image != NULL);
1497 cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
1498
1499 cx_assert((axis == 'x') || (axis == 'y'));
1500
1501
1502 nx = cpl_image_get_size_x(image);
1503 ny = cpl_image_get_size_y(image);
1504
1505
1506 /*
1507 * Validate Bias Areas and calculate number of points
1508 */
1509
1510 _areas = cpl_matrix_duplicate(areas);
1511 cx_assert(_areas != NULL);
1512
1513 nareas = cpl_matrix_get_nrow(_areas);
1514
1515 if (axis == 'y') {
1516
1517 for (j = 0; j < nareas; ++j) {
1518
1519 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1520 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1521 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1522 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1523
1524
1525 if (x_0 < 0) {
1526 x_0 = 0;
1527 cpl_matrix_set(_areas, j, 0, x_0);
1528 }
1529
1530 if (x_1 >= nx) {
1531 x_1 = nx - 1;
1532 cpl_matrix_set(_areas, j, 1, x_1);
1533 }
1534
1535 if (y_0 != 0) {
1536 y_0 = 0;
1537 cpl_matrix_set(_areas, j, 2, y_0);
1538 }
1539
1540 if (y_1 != ny - 1) {
1541 y_1 = ny - 1;
1542 cpl_matrix_set(_areas, j, 3, y_1);
1543 }
1544
1545 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1546 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
1547 }
1548
1549 }
1550
1551 sx = nareas;
1552 sy = ny;
1553 }
1554 else {
1555
1556 for (j = 0; j < nareas; ++j) {
1557
1558 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1559 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1560 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1561 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1562
1563
1564 if (x_0 != 0) {
1565 x_0 = 0;
1566 cpl_matrix_set(_areas, j, 0, x_0);
1567 }
1568
1569 if (x_1 != nx) {
1570 x_1 = nx - 1;
1571 cpl_matrix_set(_areas, j, 1, x_1);
1572 }
1573
1574 if (y_0 < 0) {
1575 y_0 = 0;
1576 cpl_matrix_set(_areas, j, 2, y_0);
1577 }
1578
1579 if (y_1 >= ny) {
1580 y_1 = ny - 1;
1581 cpl_matrix_set(_areas, j, 3, y_1);
1582 }
1583
1584 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1585 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
1586 }
1587
1588 }
1589
1590 sx = nareas;
1591 sy = nx;
1592
1593 }
1594
1595 cpl_msg_info(fctid, "Bias Areas: Found %zu points in %dx%d image",
1596 ntotal, nx, ny);
1597
1598
1599 /*
1600 * Compute Bias Areas Values
1601 */
1602
1603 results->mean = 0.;
1604 results->sigma = 0.;
1605 results->rms = 0.;
1606
1607 cx_string_set(results->limits, "");
1608
1609 profile = cpl_image_new(sx, sy, CPL_TYPE_DOUBLE);
1610
1611 for (j = 0; j < nareas; ++j) {
1612
1613 cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
1614 cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
1615 cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
1616 cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
1617
1618 if ((x_1 >= x_0) && (y_1 >= y_0)) {
1619
1620 cx_string* tmp = cx_string_new();
1621
1622
1623 cx_string_sprintf(tmp, "%d:%d:%d:%d;", x_0, x_1, y_0, y_1);
1624 cx_string_append(results->limits, cx_string_get(tmp));
1625
1626 cx_string_delete(tmp);
1627 tmp = NULL;
1628
1629
1630 /*
1631 * Update the number of actually used (valid) areas.
1632 */
1633
1634 ++nvalid;
1635
1636
1637 /*
1638 * Compute the profile along the given axis for each bias region.
1639 */
1640
1641 if (axis == 'y') {
1642
1643 cxint i = 0;
1644 cxint sz = x_1 - x_0 + 1;
1645
1646 cxdouble residual = 0.;
1647
1648 const cxdouble* _img = cpl_image_get_data_double_const(image);
1649
1650 cxdouble* _profile = cpl_image_get_data_double(profile);
1651
1652
1653 for (i = y_0; i < y_1 + 1; ++i) {
1654
1655 register cxint k = i * nx;
1656 register cxint l = 0;
1657
1658 cxdouble m = giraffe_array_mean(&_img[k + x_0], sz);
1659
1660 _profile[i * sx + j] = m;
1661
1662 for (l = x_0; l < x_1 + 1; ++l) {
1663 residual += (_img[k + l] - m) * (_img[k + l] - m);
1664 }
1665
1666 }
1667
1668
1669 /*
1670 * Compute the mean RMS and the mean bias error
1671 */
1672
1673 rms += residual / (sz - 1);
1674 sigma += residual / (sz * (sz - 1));
1675
1676 }
1677 else {
1678
1679 cxint i = 0;
1680 cxint sz = y_1 - y_0 + 1;
1681
1682 cxdouble residual = 0.;
1683
1684 const cxdouble* _img = cpl_image_get_data_double_const(image);
1685
1686 cxdouble* _profile = cpl_image_get_data_double(profile);
1687 cxdouble* buffer = cx_calloc(sz, sizeof(cxdouble));
1688
1689
1690 for (i = x_0; i < x_1 + 1; ++i) {
1691
1692 register cxint l = 0;
1693
1694 register cxdouble m = 0.;
1695
1696
1697 for (l = y_0; l < y_1 + 1; ++l) {
1698 buffer[l - y_0] = _img[l * nx + i];
1699 }
1700
1701 m = giraffe_array_mean(buffer, sz);
1702 _profile[i * sx + j] = m;
1703
1704 for (l = 0; l < sz; ++l) {
1705 residual += (buffer[l] - m) * (buffer[l] - m);
1706 }
1707
1708 }
1709
1710
1711 /*
1712 * Compute the mean RMS and the mean bias error
1713 */
1714
1715 rms += residual / (sz - 1);
1716 sigma += residual / (sz * (sz - 1));
1717
1718 cx_free(buffer);
1719 buffer = NULL;
1720
1721 }
1722
1723 }
1724
1725 }
1726
1727 cpl_matrix_delete(_areas);
1728 _areas = NULL;
1729
1730
1731 /*
1732 * Compute average profile of all bias areas
1733 */
1734
1735 model = cpl_image_collapse_create(profile, 1);
1736
1737 cpl_image_delete(profile);
1738 profile = NULL;
1739
1740 cpl_image_divide_scalar(model, nvalid);
1741 results->model = model;
1742
1743 model = NULL;
1744
1745 results->mean = cpl_image_get_mean(results->model);
1746
1747 if (nvalid > 0) {
1748 results->rms = sqrt(rms / (sy * nvalid));
1749 results->sigma = sqrt(sigma / sy) / nvalid;
1750 }
1751
1752 return EXIT_SUCCESS;
1753
1754}
1755
1756
1757/*
1758 * @brief
1759 * Fit a Chebyshev polynomial to a bias profile model.
1760 *
1761 * @param results Container for the computed mean bias, rms and the profile
1762 * model.
1763 * @param image The image of the raw profile model.
1764 * @param order The order of the Chebyshev polynomial.
1765 *
1766 * @return
1767 * The function returns @c EXIT_SUCCESS on success, and a non-zero
1768 * value otherwise.
1769 *
1770 * The function fits a one-dimensional Chebyshev polynomial of the order
1771 * @em order to the one-dimensional, raw profile image @em profile.
1772 *
1773 * The raw profile image is expected to be an image with a single column!
1774 * It may be an alias for the profile model which is present in the results
1775 * container.
1776 *
1777 * If the fit is successful, the result is stored as the new profile model
1778 * in the results container. Any previous model is overwritten!
1779 */
1780
1781inline static cxint
1782_giraffe_bias_fit_profile(GiBiasResults* results, const cpl_image* profile,
1783 cxint order)
1784{
1785
1786 cxint i = 0;
1787
1788 cxint sy = cpl_image_get_size_y(profile);
1789 cxint nc = order + 1;
1790
1791 cxdouble* _profile = (cxdouble*)cpl_image_get_data_double_const(profile);
1792
1793 cpl_matrix* base = NULL;
1794 cpl_matrix* baseT = NULL;
1795 cpl_matrix* coeff = NULL;
1796 cpl_matrix* fit = NULL;
1797 cpl_matrix* x = cpl_matrix_new(sy, 1);
1798 cpl_matrix* y = cpl_matrix_wrap(sy, 1, _profile);
1799 cpl_matrix* covy = cpl_matrix_new(sy, sy);
1800 cpl_matrix* covc = cpl_matrix_new(nc, nc);
1801
1802
1803 if ((x == NULL) || (y == NULL) || (covy == NULL) || (covc == NULL)) {
1804
1805 cpl_matrix_delete(covc);
1806 cpl_matrix_delete(covy);
1807 cpl_matrix_unwrap(y);
1808 cpl_matrix_delete(x);
1809
1810 return -1;
1811 }
1812
1813
1814 /*
1815 * Set abscissa values for the least square fit, and create the
1816 * design matrix of the fit (the basis polynomials at the abscissa values).
1817 */
1818
1819 for (i = 0; i < sy; ++i)
1820 {
1821 cpl_matrix_set(x, i, 0, i);
1822 }
1823
1824 baseT = giraffe_chebyshev_base1d(0., sy, nc, x);
1825 base = cpl_matrix_transpose_create(baseT);
1826
1827 cpl_matrix_delete(baseT);
1828 baseT = NULL;
1829
1830 cpl_matrix_delete(x);
1831 x = NULL;
1832
1833 if (base == NULL) {
1834 cpl_matrix_unwrap(y);
1835 return -2;
1836 }
1837
1838 cpl_matrix_fill_diagonal(covy, results->rms * results->rms, 0);
1839
1840
1841 /*
1842 * Compute the coefficients of the fitting polynomial, and evaluate the
1843 * fit at the given abscissa values.
1844 */
1845
1846 coeff = giraffe_matrix_solve_cholesky(base, y, covy, covc);
1847
1848 cpl_matrix_delete(covy);
1849 covy = NULL;
1850
1851 if (coeff == NULL) {
1852
1853 cpl_matrix_delete(base);
1854 base = NULL;
1855
1856 cpl_matrix_unwrap(y);
1857 y = NULL;
1858
1859 return -3;
1860
1861 }
1862
1863
1864 fit = cpl_matrix_product_create(base, coeff);
1865
1866 cpl_matrix_delete(coeff);
1867 coeff = NULL;
1868
1869 cpl_matrix_delete(base);
1870 base = NULL;
1871
1872 if (fit == NULL) {
1873
1874 cpl_matrix_unwrap(y);
1875 y = NULL;
1876
1877 return -3;
1878
1879 }
1880
1881
1882 /*
1883 * Update the model in the results container
1884 */
1885
1886 for (i = 0; i < sy; ++i)
1887 {
1888 cpl_matrix_set(y, i, 0, cpl_matrix_get(fit, i, 0));
1889 }
1890
1891 cpl_matrix_delete(fit);
1892 fit = NULL;
1893
1894 cpl_matrix_unwrap(y);
1895 y = NULL;
1896
1897
1898 /*
1899 * Update mean bias and the bias error estimates. The error of the bias
1900 * is adopted to be the error of the constant term of the bias polynomial
1901 * model.
1902 */
1903
1904 results->mean = cpl_image_get_mean(results->model);
1905 results->sigma = sqrt(cpl_matrix_get(covc, 0, 0));
1906
1907#if 0
1908 cpl_image_save(results->model, "idiot.fits", CPL_BPP_IEEE_FLOAT,
1909 0, CPL_IO_DEFAULT);
1910#endif
1911
1912 /*
1913 * Cleanup
1914 */
1915
1916 cpl_matrix_delete(covc);
1917 covc = NULL;
1918
1919 return EXIT_SUCCESS;
1920}
1921
1922
1923/*
1924 * @brief
1925 * Actual Calculation/Logic Module of Bias Removal
1926 *
1927 * @param img Input image
1928 * @param biasareas Bias area to use
1929 * @param master_bias Master bias Image
1930 * @param bad_pixel_map Bad pixel map Image
1931 * @param method Method to use for bias removal
1932 * @param option Option to use for bias removal
1933 * @param model Model to use for bias removal
1934 * @param remove_bias Remove bias from result image
1935 * @param mbias Master Bias value read from Master Bias
1936 * FITS Keywords
1937 * @param xdeg X Degree of Polynomial of Chebyshev fit
1938 * @param ydeg Y Degree of Polynomial of Chebyshev fit
1939 * @param xstep X Step to use during Chebyshev fit
1940 * @param ystep Y Step to use during Chebyshev fit
1941 * @param sigma Sigma to use during Kappa Sigma Clipping
1942 * @param numiter Max Num Iterations to use during Kappa Sigma
1943 * Clipping
1944 * @param maxfraction Maximum Fraction to use during Kappa Sigma Clipping
1945 * @param results Values determined
1946 *
1947 * @return The function returns 0 on success, or an error code otherwise.
1948 *
1949 * TBD
1950 */
1951
1952inline static cxint
1953_giraffe_bias_compute(GiBiasResults* results, cpl_image* img,
1954 const cpl_matrix* biasareas,
1955 const cpl_image* master_bias,
1956 const cpl_image* bad_pixel_map,
1957 GiBiasMethod method, GiBiasOption option,
1958 GiBiasModel model, cxbool remove_bias, cxdouble mbias,
1959 cxdouble xdeg, cxdouble ydeg, cxdouble xstep,
1960 cxdouble ystep, cxdouble sigma, cxint numiter,
1961 cxdouble maxfraction)
1962{
1963
1964 const cxchar* const fctid = "giraffe_bias_compute";
1965
1966 cxint i;
1967 cxint j;
1968 cxint img_dimx;
1969 cxint img_dimy;
1970 cxint img_centerx;
1971 cxint img_centery;
1972 cxint master_bias_dimx;
1973 cxint master_bias_dimy;
1974 cxint bad_pixel_dimx;
1975 cxint bad_pixel_dimy;
1976
1977 cxint error_code = 0;
1978
1979 cx_string* s = NULL;
1980
1981 const cpl_matrix* coeff = NULL;
1982
1983
1984 cx_assert(results != NULL);
1985 cx_assert(results->limits == NULL);
1986 cx_assert(results->coeffs == NULL);
1987 cx_assert(results->model == NULL);
1988
1989 cx_assert(img != NULL);
1990 cx_assert(cpl_image_get_type(img) == CPL_TYPE_DOUBLE);
1991
1992 cx_assert(biasareas != NULL);
1993
1994
1995 /*
1996 * Preprocessing
1997 */
1998
1999 img_dimx = cpl_image_get_size_x(img);
2000 img_dimy = cpl_image_get_size_y(img);
2001
2002 img_centerx = img_dimx / 2;
2003 img_centery = img_dimy / 2;
2004
2005 results->limits = cx_string_new();
2006
2007
2008 /*
2009 * Processing
2010 */
2011
2012 if (remove_bias) {
2013 cpl_msg_info(fctid, "Bias correction will be done.");
2014 }
2015 else {
2016 cpl_msg_warning(fctid, "No Bias correction will be done!");
2017 }
2018
2019
2020 if (model == GIBIAS_MODEL_MEAN) {
2021
2022 /*
2023 * Compute bias average value over biaslimits areas
2024 */
2025
2026 cpl_msg_info(fctid, "Using bias model 'MEAN' ...");
2027
2028 if ((option == GIBIAS_OPTION_PLANE) ||
2029 (method == GIBIAS_METHOD_PLANE)) {
2030
2031 cpl_msg_error(fctid, "Can not use MEAN model with PLANE method.");
2032 return -3;
2033
2034 }
2035
2036 error_code = _giraffe_bias_compute_mean(results, img, biasareas,
2037 sigma, numiter, maxfraction);
2038
2039
2040 }
2041 else if ((method == GIBIAS_METHOD_CURVE) ||
2042 ((method != GIBIAS_METHOD_PROFILE) &&
2043 (option == GIBIAS_OPTION_CURVE))) {
2044
2045 /*
2046 * A 2D Bias Surface is fitted on a grid of the pixels in the
2047 * bias areas
2048 */
2049
2050 cpl_msg_info(fctid, "Using bias model 'CURVE' ...");
2051
2052 error_code = _giraffe_bias_compute_curve(results, img, biasareas,
2053 sigma, numiter, maxfraction, xdeg, ydeg, xstep, ystep);
2054
2055 if (error_code != EXIT_SUCCESS) {
2056
2057 cpl_msg_error(fctid, "Error during calculation of bias curve, "
2058 "aborting...");
2059 return error_code;
2060
2061 }
2062
2063 coeff = results->coeffs;
2064
2065 }
2066 else {
2067
2068 /*
2069 * A 2D Bias Plane is fitted on the pixels of the biaslimits area
2070 */
2071
2072 cpl_msg_info(fctid, "Using bias model 'PLANE (FITTED)' ...");
2073
2074 error_code = _giraffe_bias_compute_plane(results, img, biasareas,
2075 sigma, numiter, maxfraction);
2076
2077 if (error_code == EXIT_SUCCESS) {
2078
2079 coeff = results->coeffs;
2080
2081 results->mean = cpl_matrix_get(coeff, 0, 0) +
2082 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2083 cpl_matrix_get(coeff, 0, 2) * img_centery;
2084
2085 }
2086 else {
2087
2088 cpl_msg_error(fctid, "Error during calculation of bias plane, "
2089 "aborting...");
2090 return error_code;
2091
2092 }
2093
2094 }
2095
2096
2097 /*
2098 * Perform computation based on method....
2099 */
2100
2101 s = cx_string_new();
2102 _giraffe_method_string(s, method, option);
2103
2104 cpl_msg_info(fctid, "Using bias method '%s'", cx_string_get(s));
2105
2106 cx_string_delete(s);
2107 s = NULL;
2108
2109
2110 switch (method) {
2111 case GIBIAS_METHOD_UNIFORM:
2112 {
2113
2114 /*
2115 * Subtract the average bias value
2116 */
2117
2118 if (model == GIBIAS_MODEL_MEAN) {
2119
2120 cpl_msg_info(fctid, "bias value (areas mean) = %f, "
2121 "bias sigma = %f", results->mean, results->sigma);
2122
2123 }
2124 else {
2125
2126 cpl_msg_info(fctid, "bias value at (%d, %d) = %f, "
2127 "bias sigma = %f", img_centerx, img_centery,
2128 results->mean, results->sigma);
2129
2130 }
2131
2132 if (remove_bias == TRUE) {
2133
2134 cpl_image_subtract_scalar(img, results->mean);
2135
2136 }
2137
2138 break;
2139 }
2140
2141 case GIBIAS_METHOD_PLANE:
2142 {
2143
2144 if (coeff == NULL) {
2145
2146 error_code = _giraffe_bias_compute_plane(results, img,
2147 biasareas, sigma, numiter, maxfraction);
2148
2149 if (error_code == EXIT_SUCCESS) {
2150
2151 coeff = results->coeffs;
2152
2153 results->mean = cpl_matrix_get(coeff, 0, 0) +
2154 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2155 cpl_matrix_get(coeff, 0, 2) * img_centery;
2156
2157 }
2158 else {
2159
2160 cpl_msg_error(fctid, "Error during calculation of bias "
2161 "plane, aborting...");
2162 return error_code;
2163
2164 }
2165
2166 }
2167
2168 cpl_msg_info(fctid, "bias value at (%d, %d) = %f, bias sigma = %f, "
2169 "bias plane = %g + %g * x + %g * y", img_centerx,
2170 img_centery, results->mean, results->sigma,
2171 cpl_matrix_get(coeff, 0, 0),
2172 cpl_matrix_get(coeff, 0, 1),
2173 cpl_matrix_get(coeff, 0, 2));
2174
2175 if (remove_bias == TRUE) {
2176
2177 //TODO Move plane removal to gir_image...
2178
2179 cxdouble* _img = cpl_image_get_data_double(img);
2180
2181#if 0
2182 cpl_image* bsmodel = cpl_image_new(img_dimx, img_dimy,
2183 CPL_TYPE_DOUBLE);
2184 cxdouble* _bsmodel = cpl_image_get_data_double(bsmodel);
2185#endif
2186
2187
2188 for (j = 0; j < img_dimy; j++) {
2189
2190 cxsize jj = j * img_dimx;
2191
2192 for (i = 0; i < img_dimx; i++) {
2193
2194#if 0
2195 _bsmodel[jj + i] = cpl_matrix_get(coeff, 0, 0)
2196 + cpl_matrix_get(coeff, 0, 1) * i
2197 + cpl_matrix_get(coeff, 0, 2) * j;
2198#endif
2199
2200 _img[jj + i] -= cpl_matrix_get(coeff, 0, 0)
2201 + cpl_matrix_get(coeff, 0, 1) * i
2202 + cpl_matrix_get(coeff, 0, 2) * j;
2203
2204 }
2205
2206 }
2207
2208#if 0
2209 cpl_image_save(bsmodel, "idiot.fits", CPL_BPP_IEEE_FLOAT,
2210 0, CPL_IO_DEFAULT);
2211 cpl_image_delete(bsmodel);
2212#endif
2213
2214 }
2215
2216 break;
2217
2218 }
2219
2220 case GIBIAS_METHOD_CURVE:
2221 {
2222
2223 if (coeff == NULL) {
2224
2225 error_code =
2226 _giraffe_bias_compute_curve(results, img, biasareas,
2227 sigma, numiter, maxfraction,
2228 xdeg, ydeg, xstep, ystep);
2229
2230 if (error_code != EXIT_SUCCESS) {
2231
2232 cpl_msg_error(fctid, "Error during calculation of bias "
2233 "surface curve, aborting...");
2234 return error_code;
2235
2236 }
2237
2238 coeff = results->coeffs;
2239
2240 }
2241
2242 cpl_msg_info(fctid, "bias value %f, bias sigma = %f\n",
2243 results->mean, results->sigma);
2244
2245
2246 if (remove_bias == TRUE) {
2247
2248 cpl_image* bsmodel = NULL;
2249
2250 cpl_matrix* x = cpl_matrix_new(img_dimx * img_dimy, 1);
2251 cpl_matrix* y = cpl_matrix_new(img_dimx * img_dimy, 1);
2252 cpl_matrix* fit = NULL;
2253
2254
2255 for (j = 0; j < img_dimy; ++j) {
2256
2257 register cxsize jj = j * img_dimx;
2258
2259 for (i = 0; i < img_dimx; ++i) {
2260
2261 cpl_matrix_set(x, jj + i, 0, i);
2262 cpl_matrix_set(y, jj + i, 0, j);
2263
2264 }
2265
2266 }
2267
2268 fit = giraffe_chebyshev_fit2d(0., 0., img_dimx, img_dimy,
2269 coeff, x, y);
2270
2271 cpl_matrix_delete(y);
2272 y = NULL;
2273
2274 cpl_matrix_delete(x);
2275 x = NULL;
2276
2277 bsmodel = cpl_image_wrap_double(img_dimx, img_dimy,
2278 cpl_matrix_get_data(fit));
2279
2280#if 0
2281 cpl_image_save(bsmodel, "idiot.fits", CPL_BPP_IEEE_FLOAT,
2282 0, CPL_IO_DEFAULT);
2283#endif
2284
2285 cpl_image_subtract(img, bsmodel);
2286
2287 cpl_image_unwrap(bsmodel);
2288 bsmodel = NULL;
2289
2290 cpl_matrix_delete(fit);
2291 fit = NULL;
2292
2293 }
2294
2295 break;
2296
2297 }
2298
2299 case GIBIAS_METHOD_PROFILE:
2300 {
2301
2302 error_code = _giraffe_bias_compute_profile(results, img,
2303 biasareas, 'y');
2304
2305 if (error_code != EXIT_SUCCESS) {
2306
2307 cpl_msg_error(fctid, "Error computing the bias area(s) "
2308 "profile");
2309 return error_code;
2310
2311 }
2312
2313 if (option == GIBIAS_OPTION_CURVE) {
2314
2315 error_code = _giraffe_bias_fit_profile(results, results->model,
2316 ydeg - 1);
2317 if (error_code != EXIT_SUCCESS) {
2318
2319 cpl_msg_error(fctid, "Error fitting the bias profile");
2320 return error_code;
2321
2322 }
2323
2324 }
2325
2326 if (remove_bias == TRUE) {
2327
2328 const cxdouble* _bias =
2329 cpl_image_get_data_double(results->model);
2330
2331 cxdouble* _img = cpl_image_get_data_double(img);
2332
2333
2334 cx_assert(_bias != NULL);
2335 cx_assert(_img != NULL);
2336
2337 for (j = 0; j < img_dimy; ++j) {
2338
2339 cxsize jj = j * img_dimx;
2340
2341
2342 for (i = 0; i < img_dimx; ++i) {
2343 _img[jj + i] -= _bias[j];
2344 }
2345
2346 }
2347
2348 }
2349
2350 break;
2351
2352 }
2353
2354 case GIBIAS_METHOD_MASTER:
2355 case GIBIAS_METHOD_ZMASTER:
2356 {
2357
2358 cxdouble biasdrift = 0.;
2359
2360 if (master_bias == NULL) {
2361
2362 cpl_msg_error(fctid, "Master Bias Frame required with MASTER "
2363 "method.");
2364 return -5;
2365
2366 }
2367
2368 master_bias_dimx = cpl_image_get_size_x(master_bias);
2369 master_bias_dimy = cpl_image_get_size_y(master_bias);
2370
2371 if ((master_bias_dimx != img_dimx) ||
2372 (master_bias_dimy != img_dimy)) {
2373
2374 cpl_msg_error(fctid, "Invalid master bias! Size should "
2375 "be [%d, %d] but is [%d, %d].", img_dimx, img_dimy,
2376 master_bias_dimx, master_bias_dimy);
2377 return -7;
2378
2379 }
2380
2381 if (coeff == NULL) {
2382
2383 if (option == GIBIAS_OPTION_CURVE) {
2384
2385 error_code = _giraffe_bias_compute_curve(results, img,
2386 biasareas, sigma, numiter, maxfraction, xdeg,
2387 ydeg, xstep, ystep);
2388
2389 if (error_code != EXIT_SUCCESS) {
2390
2391 cpl_msg_error(fctid, "Error during calculation of "
2392 "bias surface curve, aborting...");
2393 return error_code;
2394
2395 }
2396
2397 coeff = results->coeffs;
2398
2399 }
2400 else {
2401
2402 /*
2403 * Default to plane...
2404 */
2405
2406 error_code = _giraffe_bias_compute_plane(results, img,
2407 biasareas, sigma, numiter, maxfraction);
2408
2409 if (error_code == EXIT_SUCCESS) {
2410
2411 coeff = results->coeffs;
2412
2413 results->mean = cpl_matrix_get(coeff, 0, 0) +
2414 cpl_matrix_get(coeff, 0, 1) * img_centerx +
2415 cpl_matrix_get(coeff, 0, 2) * img_centery;
2416
2417 }
2418 else {
2419
2420 cpl_msg_error(fctid, "Error during calculation of "
2421 "bias surface plane, aborting...");
2422 return error_code;
2423
2424 }
2425
2426 }
2427
2428 }
2429
2430 if (method == GIBIAS_METHOD_ZMASTER) {
2431
2432 /*
2433 * A possible drift of the average bias is compensated using
2434 * the pre/overscans reference value as an additive
2435 * correction
2436 */
2437
2438 biasdrift = results->mean - mbias;
2439
2440 cpl_msg_info(fctid, "Using bias drift value %.4g", biasdrift);
2441
2442 }
2443
2444
2445 /*
2446 * Write some info into log file...
2447 */
2448
2449 if (option == GIBIAS_OPTION_CURVE) {
2450
2451 //TODO dump coeff to log
2452
2453 cpl_msg_info(fctid, "Computed bias mean = %.4f, bias "
2454 "sigma = %.4e", results->mean, results->sigma);
2455
2456 }
2457 else {
2458
2459 if (option == GIBIAS_OPTION_PLANE) {
2460
2461 cpl_msg_info(fctid, "Bias plane = %.4e + "
2462 "%.4e * x + %.4e * y",
2463 cpl_matrix_get(coeff, 0, 0),
2464 cpl_matrix_get(coeff, 0, 1),
2465 cpl_matrix_get(coeff, 0, 2));
2466 cpl_msg_info(fctid, "Computed bias mean = %.4f, "
2467 "bias sigma = %.4e at (%d, %d)",
2468 results->mean, results->sigma,
2469 img_centerx, img_centery);
2470
2471 }
2472 else {
2473
2474 cpl_msg_info(fctid, "Computed bias mean = %.4f, bias "
2475 "sigma = %.4e",
2476 results->mean, results->sigma);
2477
2478 }
2479
2480 }
2481
2482
2483 /*
2484 * Remove bias if requested
2485 */
2486
2487 if (remove_bias == TRUE) {
2488
2489
2490 /*
2491 * No bad pixel map was given...
2492 * Subtract master bias pixels without modification except
2493 * pre/overscan trimming.
2494 */
2495
2496 /*
2497 * Bad pixel map was given...
2498 * For all pixels not flagged in bad_pixel_map the
2499 * master_bias pixel values are subtracted, otherwise
2500 * results->mean or bias model image pixels are subtracted
2501 * depending on model.
2502 */
2503
2504 if (bad_pixel_map == NULL) {
2505
2506 cpl_image_subtract(img, master_bias);
2507
2508 if (biasdrift != 0.) {
2509 cpl_image_subtract_scalar(img, biasdrift);
2510 }
2511
2512 }
2513 else {
2514
2515 bad_pixel_dimx = cpl_image_get_size_x(bad_pixel_map);
2516 bad_pixel_dimy = cpl_image_get_size_y(bad_pixel_map);
2517
2518 if ((bad_pixel_dimx != img_dimx) ||
2519 (bad_pixel_dimy != img_dimy)) {
2520
2521 cpl_msg_error(fctid, "Invalid bad pixel map size "
2522 "should be [%d, %d] but is [%d, %d].",
2523 img_dimx, img_dimy,
2524 bad_pixel_dimx, bad_pixel_dimy);
2525 return -6;
2526
2527 }
2528
2529 if (option == GIBIAS_OPTION_CURVE) {
2530
2531 const cxint* _bpix =
2532 cpl_image_get_data_int_const(bad_pixel_map);
2533
2534 const cxdouble* _mbias =
2535 cpl_image_get_data_double_const(master_bias);
2536
2537 cxdouble* _img = cpl_image_get_data_double(img);
2538
2539 cpl_matrix* x =
2540 cpl_matrix_new(img_dimx * img_dimy, 1);
2541 cpl_matrix* y =
2542 cpl_matrix_new(img_dimx * img_dimy, 1);
2543 cpl_matrix* fit = NULL;
2544
2545
2546 for (j = 0; j < img_dimy; ++j) {
2547
2548 register cxsize jj = j * img_dimx;
2549
2550 for (i = 0; i < img_dimx; ++i) {
2551
2552 cpl_matrix_set(x, jj + i, 0, i);
2553 cpl_matrix_set(y, jj + i, 0, j);
2554
2555 }
2556 }
2557
2558 fit = giraffe_chebyshev_fit2d(0., 0.,
2559 img_dimx, img_dimy,
2560 coeff, x, y);
2561
2562 cpl_matrix_delete(y);
2563 y = NULL;
2564
2565 cpl_matrix_delete(x);
2566 x = NULL;
2567
2568
2569 for (j = 0; j < img_dimy; ++j) {
2570
2571 register cxsize jj = j * img_dimx;
2572
2573 for (i = 0; i < img_dimx; ++i) {
2574
2575 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2576
2577 _img[jj + i] -= (_mbias[jj + i] +
2578 biasdrift);
2579
2580 }
2581 else {
2582 _img[jj + i] -= cpl_matrix_get(fit, i, j);
2583 }
2584
2585 }
2586
2587 }
2588
2589 cpl_matrix_delete(fit);
2590 fit = NULL;
2591
2592 }
2593 else if (option == GIBIAS_OPTION_PLANE) {
2594
2595 const cxint* _bpix =
2596 cpl_image_get_data_int_const(bad_pixel_map);
2597
2598 const cxdouble* _mbias =
2599 cpl_image_get_data_double_const(master_bias);
2600
2601 cxdouble* _img = cpl_image_get_data_double(img);
2602
2603
2604 for (j = 0; j < img_dimy; j++) {
2605
2606 cxsize jj = j * img_dimx;
2607
2608 for (i = 0; i < img_dimx; i++) {
2609
2610 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2611
2612 _img[jj + i] -= (_mbias[jj + i] +
2613 biasdrift);
2614
2615 }
2616 else {
2617
2618 _img[jj + i] -=
2619 cpl_matrix_get(coeff, 0, 0) +
2620 cpl_matrix_get(coeff, 0, 1) * j +
2621 cpl_matrix_get(coeff, 0, 2) * i;
2622 }
2623
2624 }
2625
2626 }
2627
2628 }
2629 else {
2630
2631 const cxint* _bpix =
2632 cpl_image_get_data_int_const(bad_pixel_map);
2633
2634 const cxdouble* _mbias =
2635 cpl_image_get_data_double_const(master_bias);
2636
2637 cxdouble* _img = cpl_image_get_data_double(img);
2638
2639
2640 for (j = 0; j < img_dimy; j++) {
2641
2642 cxsize jj = j * img_dimx;
2643
2644 for (i = 0; i < img_dimx; i++) {
2645
2646 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
2647
2648 _img[jj + i] -= (_mbias[jj + i] +
2649 biasdrift);
2650
2651 }
2652 else {
2653
2654 _img[jj + i] -= results->mean;
2655
2656 }
2657
2658 }
2659
2660 }
2661
2662 }
2663
2664 }
2665
2666 }
2667
2668 break;
2669
2670 }
2671
2672 default:
2673 {
2674
2675 cpl_msg_error(fctid, "Invalid method, aborting...");
2676 return -4;
2677 break;
2678
2679 }
2680
2681 }
2682
2683
2684 cpl_msg_info(fctid, "Resulting biaslimits : %s",
2685 cx_string_get(results->limits));
2686
2687 return 0;
2688
2689}
2690
2691
2709cpl_matrix *
2710giraffe_get_raw_areas(const GiImage *image)
2711{
2712
2713 const cxchar *const _id = "giraffe_get_raw_areas";
2714
2715 cxint32 prescx = 0;
2716 cxint32 prescy = 0;
2717 cxint32 ovrscx = 0;
2718 cxint32 ovrscy = 0;
2719
2720 cxint32 winx = 0;
2721 cxint32 winy = 0;
2722
2723 cxint32 tprescx = 0;
2724 cxint32 tprescy = 0;
2725 cxint32 tovrscx = 0;
2726 cxint32 tovrscy = 0;
2727
2728 cxint32 row = 0;
2729
2730 cpl_matrix *m_tmp;
2731 cpl_propertylist *properties;
2732
2733
2734 properties = giraffe_image_get_properties(image);
2735 if (!properties) {
2736 cpl_error_set(_id, CPL_ERROR_NULL_INPUT);
2737 return NULL;
2738 }
2739
2740 winx = cpl_propertylist_get_int(properties, GIALIAS_WINX);
2741 winy = cpl_propertylist_get_int(properties, GIALIAS_WINY);
2742
2743 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
2744 tprescx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
2745 if (tprescx > 0) {
2746 prescx = tprescx;
2747 }
2748 }
2749 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
2750 tprescy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
2751 if (tprescy > 0) {
2752 prescy = tprescy;
2753 }
2754 }
2755 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
2756 tovrscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
2757 if (tovrscx > 0) {
2758 ovrscx = tovrscx;
2759 }
2760 }
2761 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
2762 tovrscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
2763 if (tovrscy > 0) {
2764 ovrscy = tovrscy;
2765 }
2766 }
2767
2768
2769 m_tmp = cpl_matrix_new(1, 4);
2770
2771 if (prescx > 0) {
2772
2773 cpl_matrix_set(m_tmp, row, 0, 0.0);
2774 cpl_matrix_set(m_tmp, row, 1, (cxdouble) prescx - 1);
2775 cpl_matrix_set(m_tmp, row, 2, 0.0);
2776 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2777
2778 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2779 row++;
2780 }
2781
2782 if (ovrscx > 0) {
2783
2784 cpl_matrix_set(m_tmp, row, 0, (cxdouble) winx - ovrscx);
2785 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2786 cpl_matrix_set(m_tmp, row, 2, 0.0);
2787 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2788
2789 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2790 row++;
2791 }
2792
2793 if (!(prescy == 0 || prescx == 0 || ovrscx == 0)) {
2794
2795 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2796 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2797 cpl_matrix_set(m_tmp, row, 2, 0.0);
2798 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2799
2800 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2801 row++;
2802
2803 }
2804 else if (!(prescy == 0 || prescx == 0)) {
2805
2806 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2807 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2808 cpl_matrix_set(m_tmp, row, 2, 0.0);
2809 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2810
2811 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2812 row++;
2813
2814 }
2815 else if (!(prescy == 0 || ovrscx == 0)) {
2816
2817 cpl_matrix_set(m_tmp, row, 0, 0.0);
2818 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2819 cpl_matrix_set(m_tmp, row, 2, 0.0);
2820 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2821
2822 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2823 row++;
2824
2825 }
2826 else if (prescy > 0) {
2827
2828 cpl_matrix_set(m_tmp, row, 0, 0.0);
2829 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2830 cpl_matrix_set(m_tmp, row, 2, 0.0);
2831 cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
2832
2833 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2834 row++;
2835 }
2836
2837 if (!(ovrscy == 0 || prescx == 0 || ovrscx == 0)) {
2838
2839 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2840 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2841 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2842 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2843
2844 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2845 row++;
2846
2847 }
2848 else if (!(ovrscy == 0 || prescx == 0)) {
2849
2850 cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
2851 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2852 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2853 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2854
2855 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2856 row++;
2857
2858 }
2859 else if (!(ovrscy == 0 || ovrscx == 0)) {
2860
2861 cpl_matrix_set(m_tmp, row, 0, 0.0);
2862 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
2863 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2864 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2865
2866 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2867 row++;
2868
2869 }
2870 else if (ovrscy > 0) {
2871
2872 cpl_matrix_set(m_tmp, row, 0, 0.0);
2873 cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
2874 cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
2875 cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
2876
2877 cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
2878 row++;
2879 }
2880
2881 cpl_matrix_resize(m_tmp, 0, -1, 0, 0);
2882 row--;
2883
2884 if (row == 0) {
2885 cpl_matrix_delete(m_tmp);
2886 return NULL;
2887 }
2888
2889 return m_tmp;
2890}
2891
2892
2911cxint
2913{
2914
2915 const cxchar *fctid = "giraffe_trim_raw_areas";
2916
2917 cxint nx, ny;
2918 cxint newlx, newly, newux, newuy;
2919
2920 cxint ovscx = 0;
2921 cxint ovscy = 0;
2922 cxint prscx = 0;
2923 cxint prscy = 0;
2924
2925 cpl_propertylist *properties = giraffe_image_get_properties(image);
2926
2927 cpl_image *_image = giraffe_image_get(image);
2928 cpl_image *tmp;
2929
2930
2931 if (!properties) {
2932 cpl_msg_error(fctid, "Image does not contain any properties!");
2933 return 1;
2934 }
2935
2936 nx = cpl_image_get_size_x(_image);
2937 ny = cpl_image_get_size_y(_image);
2938
2939 if (cpl_propertylist_has(properties, GIALIAS_NAXIS1)) {
2940 cxint _nx = cpl_propertylist_get_int(properties, GIALIAS_NAXIS1);
2941
2942 if (_nx != nx) {
2943 cpl_msg_warning(fctid, "Image size (%d) and image property `%s' "
2944 "(%d) are not consistent! Using image size ...",
2945 nx, GIALIAS_NAXIS1, _nx);
2946 }
2947 }
2948 else {
2949 cpl_msg_warning(fctid, "Image does not contain any property `%s'. "
2950 "Using image size (%d)", GIALIAS_NAXIS1, nx);
2951 }
2952
2953
2954 if (cpl_propertylist_has(properties, GIALIAS_NAXIS2)) {
2955 cxint _ny = cpl_propertylist_get_int(properties, GIALIAS_NAXIS2);
2956
2957 if (_ny != ny) {
2958 cpl_msg_warning(fctid, "Image size (%d) and image property `%s' "
2959 "(%d) are not consistent! Using image size ...",
2960 ny, GIALIAS_NAXIS2, _ny);
2961 }
2962 }
2963 else {
2964 cpl_msg_warning(fctid, "Image does not contain any property `%s'. "
2965 "Using image size (%d)", GIALIAS_NAXIS2, ny);
2966 }
2967
2968 if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
2969 ovscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
2970 }
2971
2972 if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
2973 ovscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
2974 }
2975
2976 if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
2977 prscx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
2978 }
2979
2980 if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
2981 prscy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
2982 }
2983
2984 // FIXME:
2985 // Check that the one pixel offset is ok. cpl_image_extract starts
2986 // counting from 1,1 with right and to boundaries inclusive.
2987
2988 newlx = prscx + 1;
2989 newly = prscy + 1;
2990 newux = nx - ovscx;
2991 newuy = ny - ovscy;
2992
2993 tmp = cpl_image_extract(_image, newlx, newly, newux, newuy);
2994
2995 giraffe_image_set(image, tmp);
2996 cpl_image_delete(tmp);
2997
2998 _image = giraffe_image_get(image);
2999
3000 nx = cpl_image_get_size_x(_image);
3001 ny = cpl_image_get_size_y(_image);
3002
3003
3004 /*
3005 * Update image properties
3006 */
3007
3008 cpl_propertylist_set_int(properties, GIALIAS_NAXIS1, nx);
3009 cpl_propertylist_set_int(properties, GIALIAS_NAXIS2, ny);
3010
3011 if (cpl_propertylist_has(properties, GIALIAS_CRPIX1)) {
3012 cxdouble crpix1 = cpl_propertylist_get_double(properties,
3013 GIALIAS_CRPIX1);
3014
3015 /*
3016 * For GIRAFFE the reference pixel is defined as
3017 * 1 - width(prescan)
3018 */
3019
3020 crpix1 += (cxdouble)prscx;
3021 cpl_propertylist_set_double(properties, GIALIAS_CRPIX1, crpix1);
3022 }
3023
3024 if (cpl_propertylist_has(properties, GIALIAS_CRPIX2)) {
3025 cxdouble crpix2 = cpl_propertylist_get_double(properties,
3026 GIALIAS_CRPIX2);
3027
3028 crpix2 -= (cxdouble) prscy;
3029 cpl_propertylist_set_double(properties, GIALIAS_CRPIX2, crpix2);
3030 }
3031
3032 cpl_propertylist_erase(properties, GIALIAS_OVSCX);
3033 cpl_propertylist_erase(properties, GIALIAS_OVSCY);
3034 cpl_propertylist_erase(properties, GIALIAS_PRSCX);
3035 cpl_propertylist_erase(properties, GIALIAS_PRSCY);
3036
3037 return 0;
3038
3039}
3040
3041
3103cxint
3104giraffe_bias_remove(GiImage* result, const GiImage* raw,
3105 const GiImage* master_bias, const GiImage* bad_pixels,
3106 const cpl_matrix* biaslimits, const GiBiasConfig* config)
3107{
3108
3109 const cxchar* const fctid = "giraffe_bias_remove";
3110
3111 cxint error_code;
3112
3113 cxdouble bias_value = 0.;
3114
3115 cx_string* s;
3116
3117 cpl_matrix* areas = NULL;
3118
3119 cpl_image* _raw = giraffe_image_get(raw);
3120 cpl_image* _master_bias = giraffe_image_get(master_bias);
3121 cpl_image* _bad_pixels = giraffe_image_get(bad_pixels);
3122 cpl_image* timage;
3123
3124 cpl_propertylist* properties;
3125
3126 GiBiasResults bias_results = {0., 0., 0., NULL, NULL, NULL};
3127
3128
3129 cx_assert(raw != NULL);
3130 cx_assert(config != NULL);
3131 cx_assert(biaslimits == NULL);
3132
3133 if (result == NULL) {
3134 return 1;
3135 }
3136
3137
3138 /*
3139 * Setup user defined areas to use for the bias computation
3140 */
3141
3142 if (strncmp(config->areas, "None", 4) == 0) {
3143
3144 cpl_msg_info(fctid, "No bias areas specified, using pre/overscan"
3145 "regions of the raw frame.");
3146
3147 cpl_error_reset();
3148 areas = giraffe_get_raw_areas(raw);
3149 if (areas == NULL) {
3150 if (cpl_error_get_code() == CPL_ERROR_NULL_INPUT) {
3151 cpl_msg_error(fctid, "Invalid raw image! Image does not "
3152 "contain any properties!");
3153 }
3154 else {
3155 cpl_msg_error(fctid, "Invalid raw image! Image does not "
3156 "contain or has invalid pre- and overscan "
3157 "properties.");
3158 }
3159
3160 return 1;
3161 }
3162
3163 }
3164 else {
3165
3166 areas = _giraffe_bias_get_areas(config->areas);
3167
3168 if (areas == NULL) {
3169
3170 cpl_msg_error(fctid, "Invalid bias area specification '%s'!",
3171 config->areas);
3172 return 1;
3173
3174 }
3175
3176 cpl_msg_info(fctid, "Using bias area(s) '%s' for bias computation",
3177 config->areas);
3178
3179 }
3180
3181
3182 /*
3183 * Processing
3184 */
3185
3186
3187 /*
3188 * Check whether a masterbias image was specified
3189 * If so, compare pre/ovrscan keywords/areas
3190 */
3191
3192 if (master_bias != NULL) {
3193 cpl_propertylist *_properties;
3194 cxbool is_same_size = FALSE;
3195
3196 is_same_size = _giraffe_compare_overscans(raw, master_bias);
3197
3198 // FIXME:
3199 // Fully processed calibrations usually do not have overscans
3200 // any longer. Instead of failing at this point some dummies
3201 // could be created.
3202
3203 if (is_same_size==FALSE) {
3204 cpl_msg_info(fctid, "PRE/OVRSCAN Regions do not match between "
3205 "RAW Image and Masterbias Image.");
3206
3207 return 1;
3208 }
3209
3210 _properties = giraffe_image_get_properties(master_bias);
3211
3212 if (cpl_propertylist_has(_properties, GIALIAS_BIASVALUE)) {
3213 bias_value = cpl_propertylist_get_double(_properties,
3214 GIALIAS_BIASVALUE);
3215 }
3216 }
3217
3218
3219 /*
3220 * Check whether a Bad Pixel Map image was specified
3221 * If so, compare pre/ovrscan keywords/areas
3222 */
3223
3224 if (bad_pixels != NULL) {
3225 cxbool is_same_size = FALSE;
3226
3227 is_same_size = _giraffe_compare_overscans(raw, bad_pixels);
3228
3229 // FIXME:
3230 // Fully processed calibrations usually do not have overscans
3231 // any longer. Instead of failing at this point some dummies
3232 // could be created.
3233
3234 if (is_same_size == FALSE) {
3235 cpl_msg_info(fctid, "PRE/OVRSCAN Regions do not match between "
3236 "RAW Image and Bad Pixel Image.");
3237
3238 return 1;
3239 }
3240 }
3241
3242
3243 /*
3244 * Create a copy of the raw frame. We are not working on the raw
3245 * frame directly.
3246 */
3247
3248 // FIXME:
3249 // Can we use the raw frame directly? (RP)
3250
3251 timage = cpl_image_duplicate(_raw);
3252
3253
3254 /*
3255 * Remove Bias Calculation. The bias corrected image will be
3256 * stored in timage.
3257 */
3258
3259 error_code = _giraffe_bias_compute(&bias_results, timage,
3260 areas, _master_bias,
3261 _bad_pixels, config->method,
3262 config->option, config->model,
3263 config->remove, bias_value,
3264 config->xdeg, config->ydeg,
3265 config->xstep, config->ystep,
3266 config->sigma, config->iterations,
3267 config->fraction);
3268
3269 cpl_matrix_delete(areas);
3270 areas = NULL;
3271
3272
3273 // FIXME:
3274 // Possibly check returned error code and do a dedicated
3275 // exception handling (RP)
3276
3277 if (error_code) {
3278
3279 _giraffe_biasresults_clear(&bias_results);
3280
3281 cpl_msg_error(fctid, "Bias computation failed with error %d!",
3282 error_code);
3283 return 1;
3284
3285 }
3286
3287
3288 /*
3289 * Should we remove Bias or not?
3290 */
3291
3292 properties = giraffe_image_get_properties(raw);
3293
3294 if (config->remove == TRUE) {
3295
3296 giraffe_image_set(result, timage);
3297 cpl_image_delete(timage);
3298
3299 giraffe_image_set_properties(result, properties);
3300
3301 }
3302 else {
3303
3304 cpl_image_delete(timage);
3305
3306 giraffe_image_set(result, _raw);
3307 giraffe_image_set_properties(result, properties);
3308
3309 }
3310
3311
3312 /*
3313 * Postprocessing
3314 */
3315
3316 properties = giraffe_image_get_properties(result);
3317
3318 if (config->remove == TRUE) {
3319
3320 cpl_propertylist_set_int(properties, GIALIAS_BITPIX, -32);
3321 cpl_propertylist_set_double(properties, GIALIAS_BZERO, 0.);
3322 cpl_propertylist_set_double(properties, GIALIAS_BSCALE, 1.);
3323
3324 if (cpl_propertylist_has(properties, GIALIAS_GIRFTYPE)) {
3325 cpl_propertylist_set_string(properties,
3326 GIALIAS_GIRFTYPE, "BRMIMG");
3327 }
3328 else {
3329 cpl_propertylist_append_string(properties,
3330 GIALIAS_GIRFTYPE, "BRMIMG");
3331 }
3332 cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
3333 "GIRAFFE bias subtracted frame");
3334
3335 /*
3336 * Remove pre- and overscans
3337 */
3338
3339 giraffe_trim_raw_areas(result);
3340
3341
3342 /*
3343 * Convert bias corrected image to electrons
3344 */
3345
3346 // FIXME:
3347 // Is this really needed? Disabled on request of DFO (RP)
3348 // If it has to be put back use giraffe_propertylist_get_conad()!
3349#if 0
3350 if (cpl_propertylist_has(properties, GIALIAS_CONAD)) {
3351 cxdouble conad = cpl_propertylist_get_double(properties, GIALIAS_CONAD);
3352
3353 if (conad > 0.) {
3354
3355 cpl_image* _result = giraffe_image_get(result);
3356
3357 cpl_msg_info(fctid, "Converting bias subtracted frame to "
3358 "electrons; conversion factor is %.2f e-/ADU",
3359 conad);
3360 cpl_image_multiply_scalar(_result, conad);
3361 }
3362 else {
3363 cpl_msg_warning(fctid, "Invalid conversion factor %.2f "
3364 "e-/ADU! Bias subtracted image will not be "
3365 "converted.", conad);
3366 }
3367 }
3368 else {
3369 cpl_msg_info(fctid, "Conversion factor not found. Bias subtracted "
3370 "image will remain in ADU");
3371 }
3372#endif
3373 }
3374
3375 s = cx_string_new();
3376
3377 _giraffe_method_string(s, config->method, config->option);
3378 cpl_propertylist_append_string(properties, GIALIAS_BIASMETHOD, cx_string_get(s));
3379
3380 cx_string_delete(s);
3381
3382 cpl_propertylist_append_double(properties, GIALIAS_BIASVALUE,
3383 bias_results.mean);
3384 cpl_propertylist_append_double(properties, GIALIAS_BIASERROR,
3385 bias_results.sigma);
3386 cpl_propertylist_append_double(properties, GIALIAS_BIASSIGMA,
3387 bias_results.rms);
3388
3389 if (bias_results.coeffs) {
3390 cpl_propertylist_append_double(properties, GIALIAS_BCLIPSIGMA,
3391 config->sigma);
3392 cpl_propertylist_append_int(properties, GIALIAS_BCLIPNITER,
3393 config->iterations);
3394 cpl_propertylist_append_double(properties, GIALIAS_BCLIPMFRAC,
3395 config->fraction);
3396 }
3397
3398 if (bias_results.limits) {
3399 cpl_propertylist_append_string(properties, GIALIAS_BIASAREAS,
3400 cx_string_get(bias_results.limits));
3401 }
3402
3403 if (bias_results.coeffs) {
3404 s = cx_string_new();
3405
3406 _giraffe_stringify_coefficients(s, bias_results.coeffs);
3407 cpl_propertylist_append_string(properties, GIALIAS_BIASPLANE,
3408 cx_string_get(s));
3409
3410 cx_string_delete(s);
3411 }
3412
3413
3414 /*
3415 * Cleanup
3416 */
3417
3418 _giraffe_biasresults_clear(&bias_results);
3419
3420 return 0;
3421
3422}
3423
3424
3435GiBiasConfig*
3436giraffe_bias_config_create(cpl_parameterlist* list)
3437{
3438
3439 const cxchar* s;
3440 cpl_parameter* p;
3441
3442 GiBiasConfig* config = NULL;
3443
3444
3445 if (!list) {
3446 return NULL;
3447 }
3448
3449 config = cx_calloc(1, sizeof *config);
3450
3451
3452 /*
3453 * Some defaults
3454 */
3455
3456 config->method = GIBIAS_METHOD_UNDEFINED;
3457 config->option = GIBIAS_OPTION_UNDEFINED;
3458 config->model = GIBIAS_MODEL_UNDEFINED;
3459 config->mbias = 0.;
3460 config->xdeg = 1;
3461 config->ydeg = 1;
3462
3463
3464 p = cpl_parameterlist_find(list, "giraffe.biasremoval.remove");
3465 config->remove = cpl_parameter_get_bool(p);
3466
3467 p = cpl_parameterlist_find(list, "giraffe.biasremoval.method");
3468 s = cpl_parameter_get_string(p);
3469
3470 if (!strcmp(s, "UNIFORM")) {
3471 config->method = GIBIAS_METHOD_UNIFORM;
3472 }
3473
3474 if (!strcmp(s, "PLANE")) {
3475 config->method = GIBIAS_METHOD_PLANE;
3476 }
3477
3478 if (!strcmp(s, "CURVE")) {
3479 config->method = GIBIAS_METHOD_CURVE;
3480 }
3481
3482 if (!strcmp(s, "PROFILE")) {
3483 config->method = GIBIAS_METHOD_PROFILE;
3484 }
3485
3486 if (!strcmp(s, "MASTER")) {
3487 config->method = GIBIAS_METHOD_MASTER;
3488 }
3489
3490 if (!strcmp(s, "ZMASTER")) {
3491 config->method = GIBIAS_METHOD_ZMASTER;
3492 }
3493
3494 if (!strcmp(s, "PROFILE+CURVE")) {
3495 config->method = GIBIAS_METHOD_PROFILE;
3496 config->option = GIBIAS_OPTION_CURVE;
3497 }
3498
3499 if (!strcmp(s, "MASTER+PLANE")) {
3500 config->method = GIBIAS_METHOD_MASTER;
3501 config->option = GIBIAS_OPTION_PLANE;
3502 }
3503
3504 if (!strcmp(s, "ZMASTER+PLANE")) {
3505 config->method = GIBIAS_METHOD_ZMASTER;
3506 config->option = GIBIAS_OPTION_PLANE;
3507 }
3508
3509 if (!strcmp(s, "MASTER+CURVE")) {
3510 config->method = GIBIAS_METHOD_MASTER;
3511 config->option = GIBIAS_OPTION_CURVE;
3512 }
3513
3514 if (!strcmp(s, "ZMASTER+CURVE")) {
3515 config->method = GIBIAS_METHOD_ZMASTER;
3516 config->option = GIBIAS_OPTION_CURVE;
3517 }
3518
3519 cx_assert(config->method != GIBIAS_METHOD_UNDEFINED);
3520
3521 p = cpl_parameterlist_find(list, "giraffe.biasremoval.areas");
3522 s = cpl_parameter_get_string(p);
3523 config->areas = cx_strdup(s);
3524
3525 p = cpl_parameterlist_find(list, "giraffe.biasremoval.sigma");
3526 config->sigma = cpl_parameter_get_double(p);
3527
3528 p = cpl_parameterlist_find(list, "giraffe.biasremoval.iterations");
3529 config->iterations = cpl_parameter_get_int(p);
3530
3531 p = cpl_parameterlist_find(list, "giraffe.biasremoval.fraction");
3532 config->fraction = cpl_parameter_get_double(p);
3533
3534 if (config->method == GIBIAS_METHOD_CURVE ||
3535 config->option == GIBIAS_OPTION_CURVE) {
3536 p = cpl_parameterlist_find(list, "giraffe.biasremoval.xorder");
3537 config->xdeg = 1 + cpl_parameter_get_int(p);
3538
3539 p = cpl_parameterlist_find(list, "giraffe.biasremoval.yorder");
3540 config->ydeg = 1 + cpl_parameter_get_int(p);
3541 }
3542
3543 p = cpl_parameterlist_find(list, "giraffe.biasremoval.xstep");
3544 config->xstep = cpl_parameter_get_int(p);
3545
3546 p = cpl_parameterlist_find(list, "giraffe.biasremoval.ystep");
3547 config->ystep = cpl_parameter_get_int(p);
3548
3549 return config;
3550
3551}
3552
3553
3566void
3567giraffe_bias_config_destroy(GiBiasConfig* config)
3568{
3569
3570 if (config) {
3571 if (config->areas) {
3572 cx_free(config->areas);
3573 config->areas = NULL;
3574 }
3575
3576 cx_free(config);
3577 }
3578
3579 return;
3580}
3581
3582
3594void
3595giraffe_bias_config_add(cpl_parameterlist* list)
3596{
3597
3598 cpl_parameter* p;
3599
3600
3601 if (!list) {
3602 return;
3603 }
3604
3605 p = cpl_parameter_new_value("giraffe.biasremoval.remove",
3606 CPL_TYPE_BOOL,
3607 "Enable bias removal",
3608 "giraffe.biasremoval",
3609 TRUE);
3610 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "remove-bias");
3611 cpl_parameterlist_append(list, p);
3612
3613
3614 p = cpl_parameter_new_enum("giraffe.biasremoval.method",
3615 CPL_TYPE_STRING,
3616 "Bias removal method",
3617 "giraffe.biasremoval",
3618 "PROFILE", 11, "UNIFORM", "PLANE", "CURVE",
3619 "PROFILE", "MASTER", "ZMASTER", "PROFILE+CURVE",
3620 "MASTER+PLANE", "MASTER+CURVE", "ZMASTER+PLANE",
3621 "ZMASTER+CURVE");
3622 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-method");
3623 cpl_parameterlist_append(list, p);
3624
3625
3626 p = cpl_parameter_new_value("giraffe.biasremoval.areas",
3627 CPL_TYPE_STRING,
3628 "Bias areas to use "
3629 "(Xl0:Xr0:Yl0:Yr0, ... ,Xln:Xrn:Yln:Yrn)",
3630 "giraffe.biasremoval",
3631 "5:40:0:4095");
3632 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-areas");
3633 cpl_parameterlist_append(list, p);
3634
3635
3636 p = cpl_parameter_new_value("giraffe.biasremoval.sigma",
3637 CPL_TYPE_DOUBLE,
3638 "Sigma Clipping: sigma threshold factor",
3639 "giraffe.biasremoval",
3640 2.5);
3641 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-sigma");
3642 cpl_parameterlist_append(list, p);
3643
3644
3645 p = cpl_parameter_new_value("giraffe.biasremoval.iterations",
3646 CPL_TYPE_INT,
3647 "Sigma Clipping: maximum number of "
3648 "iterations",
3649 "giraffe.biasremoval",
3650 5);
3651 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-niter");
3652 cpl_parameterlist_append(list, p);
3653
3654
3655 p = cpl_parameter_new_value("giraffe.biasremoval.fraction",
3656 CPL_TYPE_DOUBLE,
3657 "Sigma Clipping: minimum fraction of points "
3658 "accepted/total [0.0..1.0]",
3659 "giraffe.biasremoval",
3660 0.8);
3661 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-mfrac");
3662 cpl_parameterlist_append(list, p);
3663
3664
3665 p = cpl_parameter_new_value("giraffe.biasremoval.xorder",
3666 CPL_TYPE_INT,
3667 "Order of X polynomial fit (CURVE method "
3668 "only)",
3669 "giraffe.biasremoval",
3670 1);
3671 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-xorder");
3672 cpl_parameterlist_append(list, p);
3673
3674 p = cpl_parameter_new_value("giraffe.biasremoval.yorder",
3675 CPL_TYPE_INT,
3676 "Order of Y polynomial fit (CURVE method "
3677 "only)",
3678 "giraffe.biasremoval",
3679 1);
3680 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-yorder");
3681 cpl_parameterlist_append(list, p);
3682
3683
3684 p = cpl_parameter_new_value("giraffe.biasremoval.xstep",
3685 CPL_TYPE_INT,
3686 "Sampling step along X (CURVE method only)",
3687 "giraffe.biasremoval",
3688 1);
3689 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-xstep");
3690 cpl_parameterlist_append(list, p);
3691
3692
3693 p = cpl_parameter_new_value("giraffe.biasremoval.ystep",
3694 CPL_TYPE_INT,
3695 "Sampling step along Y (CURVE method only)",
3696 "giraffe.biasremoval",
3697 1);
3698 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-ystep");
3699 cpl_parameterlist_append(list, p);
3700
3701 return;
3702}
cxint giraffe_trim_raw_areas(GiImage *image)
Remove pre- and overscan ares from an image.
Definition: gibias.c:2912
GiBiasConfig * giraffe_bias_config_create(cpl_parameterlist *list)
Creates a setup structure for a bias removal task.
Definition: gibias.c:3436
void giraffe_bias_config_add(cpl_parameterlist *list)
Adds parameters for the bias removal.
Definition: gibias.c:3595
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:3104
cpl_matrix * giraffe_get_raw_areas(const GiImage *image)
Create bias areas from an image.
Definition: gibias.c:2710
void giraffe_bias_config_destroy(GiBiasConfig *config)
Destroys a bias removal setup structure.
Definition: gibias.c:3567
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 13:47:22 by doxygen 1.9.6 written by Dimitri van Heesch, © 1997-2004