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