GIRAFFE Pipeline Reference Manual

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

This file is part of the GIRAFFE Pipeline Reference Manual 2.14.
Documentation copyright © 2002-2006 European Southern Observatory.
Generated on Wed Mar 11 2015 13:19:41 by doxygen 1.8.9.1 written by Dimitri van Heesch, © 1997-2004