GIRAFFE Pipeline Reference Manual

gibias.c

00001 /* $Id: gibias.c,v 1.27 2007/08/23 08:39:06 rpalsa Exp $
00002  *
00003  * This file is part of the GIRAFFE Pipeline
00004  * Copyright (C) 2002-2006 European Southern Observatory
00005  *
00006  * This program is free software; you can redistribute it and/or modify
00007  * it under the terms of the GNU General Public License as published by
00008  * the Free Software Foundation; either version 2 of the License, or
00009  * (at your option) any later version.
00010  *
00011  * This program is distributed in the hope that it will be useful,
00012  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00013  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014  * GNU General Public License for more details.
00015  *
00016  * You should have received a copy of the GNU General Public License
00017  * along with this program; if not, write to the Free Software
00018  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00019  */
00020 
00021 /*
00022  * $Author: rpalsa $
00023  * $Date: 2007/08/23 08:39:06 $
00024  * $Revision: 1.27 $
00025  * $Name: giraffe-2_5_2 $
00026  */
00027 
00028 #ifdef HAVE_CONFIG_H
00029 #  include <config.h>
00030 #endif
00031 
00032 #include <string.h>
00033 #include <math.h>
00034 
00035 #include <cxmemory.h>
00036 #include <cxstrutils.h>
00037 #include <cxutils.h>
00038 
00039 #include <cpl_msg.h>
00040 #include <cpl_parameter.h>
00041 
00042 #include "gimacros.h"
00043 #include "gialias.h"
00044 #include "gimatrix.h"
00045 #include "gimath.h"
00046 #include "gibias.h"
00047 
00048 
00057 struct GiBiasResults {
00058     cxdouble   mean;
00059     cxdouble   sigma;
00060     cx_string  *limits;
00061     cpl_matrix *coeffs;
00062 };
00063 
00064 typedef struct GiBiasResults GiBiasResults;
00065 
00066 
00067 /*
00068  * Destroys a bias results structure
00069  */
00070 
00071 inline static void
00072 _giraffe_bias_results_destroy(GiBiasResults *results)
00073 {
00074 
00075     if (results) {
00076         if (results->limits) {
00077             cx_string_delete(results->limits);
00078             results->limits = NULL;
00079         }
00080 
00081         if (results->coeffs) {
00082             cpl_matrix_delete(results->coeffs);
00083             results->coeffs = NULL;
00084         }
00085 
00086         cx_free(results);
00087     }
00088 
00089     return;
00090 
00091 }
00092 
00093 
00094 /*
00095  * Transfor bias method and option value into a string
00096  */
00097 
00098 inline static void
00099 _giraffe_method_string(cx_string *string, GiBiasMethod method,
00100                        GiBiasOption option)
00101 {
00102 
00103     switch (method) {
00104     case GIBIAS_METHOD_UNIFORM:
00105         cx_string_set(string, "UNIFORM");
00106         break;
00107 
00108     case GIBIAS_METHOD_PLANE:
00109         cx_string_set(string, "PLANE");
00110         break;
00111 
00112     case GIBIAS_METHOD_CURVE:
00113         cx_string_set(string, "CURVE");
00114         break;
00115 
00116     case GIBIAS_METHOD_MASTER:
00117         cx_string_set(string, "MASTER");
00118         break;
00119 
00120     case GIBIAS_METHOD_ZMASTER:
00121         cx_string_set(string, "ZMASTER");
00122         break;
00123 
00124     default:
00125         break;
00126     }
00127 
00128     if (option != GIBIAS_OPTION_UNDEFINED) {
00129         switch (option) {
00130         case GIBIAS_OPTION_PLANE:
00131             cx_string_append(string, "+PLANE");
00132             break;
00133 
00134         case GIBIAS_OPTION_CURVE:
00135             cx_string_append(string, "+CURVE");
00136             break;
00137 
00138         default:
00139             break;
00140         }
00141     }
00142 
00143     return;
00144 
00145 }
00146 
00147 
00148 inline static void
00149 _giraffe_stringify_coefficients(cx_string *string, cpl_matrix *matrix)
00150 {
00151 
00152     register cxint i, j;
00153 
00154     cxchar *tmp = cx_line_alloc();
00155 
00156     cxint nr = cpl_matrix_get_nrow(matrix);
00157     cxint nc = cpl_matrix_get_ncol(matrix);
00158 
00159     cxsize sz = cx_line_max();
00160 
00161     cxdouble *data = cpl_matrix_get_data(matrix);
00162 
00163 
00164     for (i = 0; i < nr; i++) {
00165         for (j = 0; j < nc; j++) {
00166             snprintf(tmp, sz, "%g", data[i * nc  + j]);
00167             cx_string_append(string, tmp);
00168 
00169             if (i != nr - 1 || j < nc - 1) {
00170                 cx_string_append(string, ":");
00171             }
00172         }
00173     }
00174 
00175     cx_free(tmp);
00176 
00177     return;
00178 
00179 }
00180 
00181 
00197 inline static cxbool
00198 _giraffe_compare_overscans(GiImage *image1, GiImage *image2)
00199 {
00200 
00201     cxint32 l1ovscx = -1;
00202     cxint32 l1ovscy = -1;
00203     cxint32 l1prscx = -1;
00204     cxint32 l1prscy = -1;
00205     cxint32 l2ovscx = -1;
00206     cxint32 l2ovscy = -1;
00207     cxint32 l2prscx = -1;
00208     cxint32 l2prscy = -1;
00209 
00210     cpl_propertylist *l1, *l2;
00211 
00212 
00213     cx_assert(image1 != NULL && image2 != NULL);
00214 
00215     l1 = giraffe_image_get_properties(image1);
00216     l2 = giraffe_image_get_properties(image2);
00217 
00218     if (cpl_propertylist_has(l1, GIALIAS_OVSCX)) {
00219         l1ovscx = cpl_propertylist_get_int(l1, GIALIAS_OVSCX);
00220     }
00221     if (cpl_propertylist_has(l1, GIALIAS_OVSCY)) {
00222         l1ovscy = cpl_propertylist_get_int(l1, GIALIAS_OVSCY);
00223     }
00224     if (cpl_propertylist_has(l1, GIALIAS_PRSCX)) {
00225         l1prscx = cpl_propertylist_get_int(l1, GIALIAS_PRSCX);
00226     }
00227     if (cpl_propertylist_has(l1, GIALIAS_PRSCY)) {
00228         l1prscy = cpl_propertylist_get_int(l1, GIALIAS_PRSCY);
00229     }
00230 
00231     if (cpl_propertylist_has(l2, GIALIAS_OVSCX)) {
00232         l2ovscx = cpl_propertylist_get_int(l2, GIALIAS_OVSCX);
00233     }
00234     if (cpl_propertylist_has(l2, GIALIAS_OVSCY)) {
00235         l2ovscy = cpl_propertylist_get_int(l2, GIALIAS_OVSCY);
00236     }
00237     if (cpl_propertylist_has(l2, GIALIAS_PRSCX)) {
00238         l2prscx = cpl_propertylist_get_int(l2, GIALIAS_PRSCX);
00239     }
00240     if (cpl_propertylist_has(l2, GIALIAS_PRSCY)) {
00241         l2prscy = cpl_propertylist_get_int(l2, GIALIAS_PRSCY);
00242     }
00243 
00244     if (l1ovscx != l2ovscx || l1ovscy != l2ovscy) {
00245         return FALSE;
00246     }
00247 
00248     if (l1prscx != l2prscx || l1prscy != l2prscy) {
00249         return FALSE;
00250     }
00251 
00252     return TRUE;
00253 
00254 }
00255 
00256 
00257 /*
00258  * @brief
00259  *   Compute bias mean and sigma values over bias areas.
00260  *
00261  * @param  image        Image over which to calculate mean and sigma
00262  *                      inside bias areas
00263  * @param  areas        Bias areas specifications [N,4] as a matrix
00264  *                      where N is the number of areas
00265  * @param  sigma        Sigma clipping algorithm: sigma value
00266  * @param  numiter      Sigma clipping algorithm: max number of iterations
00267  * @param  maxfraction  Sigma clipping algorithm: max. fraction of pixels
00268  *
00269  * @param  results      A @em GiBiasResults structure which contains
00270  *                      the computed mean bias value, its standard deviation
00271  *                      and the coordinates specifications of the areas used
00272  *                      on return.
00273  *
00274  * @returns EXIT_SUCCESS on sucess, negative value otherwise
00275  *
00276  * Computes mean average of pixels value over areas defined by @em areas.
00277  * @em areas is a matrix [N, 4] specifying the lower left and upper right
00278  * corners of each of the @em N areas of the CCD used for bias determination.
00279  *
00280  * @code
00281  *   areas = [[xmin1, xmax1, ymin1, ymax1],
00282  *            [xmin2, xmax2, ymin2, ymax2],
00283  *             ....., ....., ....., .....
00284  *            [xminN, xmaxN, yminN, ymaxN]]
00285  * @endcode
00286  *
00287  * Areas pixels whose value is too far from mean average are rejected by
00288  * kappa-sigma clipping defined by parameters sigma, numiter and maxfraction.
00289  * This rejection loops until good/bad points ratio is enough or the
00290  * maximum number of iteration has been reached.
00291  *
00292  * The function returns the bias mean @em results->mean, the bias sigma
00293  * @em results->sigma computed over the kappa-sigma clipped areas pixels
00294  * value and @em results->limits string specifying the areas used as
00295  * follows:
00296  * @code
00297  *   "xmin1:xmax1:ymin1:ymax1;...;xminN:xmaxN:yminN:ymaxN"
00298  * @endcode
00299  *
00300  * @see giraffe_bias_compute_plane, giraffe_compute_remove_bias
00301  */
00302 
00303 inline static cxint
00304 _giraffe_bias_compute_mean(cpl_image *image, cpl_matrix *areas, cxdouble sigma,
00305                            cxint numiter, cxdouble maxfraction,
00306                            GiBiasResults *results)
00307 {
00308 
00309     const cxchar *fctid = "giraffe_bias_compute_mean";
00310 
00311     cxbool not_ok = FALSE;
00312 
00313     cxint32 j, l, k;
00314     cxint32 img_dimx, img_dimy;
00315     cxint32 ba_num_cols;
00316     cxint32 ba_num_rows;
00317     cxint32 curriter = 0;
00318     cxint32 x0, x1, x2, x3;
00319 
00320     cxlong ntotal    = 0L;
00321     cxlong naccepted = 0L;
00322     cxlong n         = 0L;
00323     cxlong pixcount  = 0L;
00324 
00325     cxdouble *pdimg;
00326     cxdouble currfraction = 2.0;
00327 
00328     cx_string *tmp = NULL;
00329 
00330     cpl_matrix *matrix_zz1;
00331     cpl_matrix *matrix_zz1diff;
00332 
00333 
00334     /*
00335      * Initialization
00336      */
00337 
00338     if (results->limits == NULL) {
00339         cpl_msg_info(fctid, "Unable to store biaslimits return parameter, "
00340                      "aborting...");
00341         return -3;
00342     }
00343 
00344     if (cpl_image_get_type(image) != CPL_TYPE_DOUBLE) {
00345         cpl_msg_info(fctid, "Only images allowed of type double, "
00346                      "aborting ...");
00347         return -3;
00348     }
00349 
00350 
00351     /*
00352      * Preprocessing
00353      */
00354 
00355     /*
00356      *  Validate Bias Areas and calculate number of points
00357      */
00358 
00359     if (areas == NULL) {
00360         cpl_msg_info(fctid, "Bias Areas: Missing bias areas, "
00361                      "aborting ...");
00362         return -1;
00363     }
00364 
00365     ba_num_cols = cpl_matrix_get_ncol(areas);
00366     ba_num_rows = cpl_matrix_get_nrow(areas);
00367 
00368     for (j = 0; j < ba_num_rows; j++) {
00369         x3 = (cxint32) cpl_matrix_get(areas, j, 3);
00370         x2 = (cxint32) cpl_matrix_get(areas, j, 2);
00371         x1 = (cxint32) cpl_matrix_get(areas, j, 1);
00372         x0 = (cxint32) cpl_matrix_get(areas, j, 0);
00373 
00374         ntotal += (cxulong) ((x3 - x2) * (x1 - x0));
00375     }
00376 
00377     if (ntotal <= 0) {
00378         cpl_msg_info(fctid, "Bias Areas: Inconsistent specification, "
00379                      "aborting ...");
00380         return -1;
00381     }
00382 
00383     matrix_zz1     = cpl_matrix_new(ntotal, 1);
00384     matrix_zz1diff = cpl_matrix_new(ntotal, 1);
00385 
00386 
00387     /*
00388      *  Compute Bias Areas Values
00389      */
00390 
00391     img_dimx = cpl_image_get_size_x(image);
00392     img_dimy = cpl_image_get_size_y(image);
00393 
00394     cx_string_set(results->limits, "");
00395 
00396     for (j = 0; j < ba_num_rows; j++) {
00397 
00398         x3 = (cxint32) cpl_matrix_get(areas, j, 3);
00399         x2 = (cxint32) cpl_matrix_get(areas, j, 2);
00400         x1 = (cxint32) cpl_matrix_get(areas, j, 1);
00401         x0 = (cxint32) cpl_matrix_get(areas, j, 0);
00402 
00403         not_ok = ((x0 > img_dimx) || (x1 > img_dimx) ||
00404                   (x2 > img_dimy) || (x3 > img_dimy));
00405 
00406         if (not_ok) {
00407             continue;
00408         }
00409 
00410         tmp = cx_string_new();
00411 
00412         cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
00413         cx_string_append(results->limits, cx_string_get(tmp));
00414 
00415         cx_string_delete(tmp);
00416         tmp = NULL;
00417 
00418         pdimg = cpl_image_get_data_double(image);
00419 
00420         for (l = x2; l < x3; l++) {
00421             for (k = x0; k < x1; k++) {
00422                 cpl_matrix_set(matrix_zz1, n, 1, pdimg[k + l * img_dimx]);
00423                 n++;
00424             }
00425         }
00426     }
00427 
00428     if (n != ntotal) {
00429         cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting ...");
00430 
00431         cpl_matrix_delete(matrix_zz1);
00432         cpl_matrix_delete(matrix_zz1diff);
00433 
00434         return -3;
00435     }
00436 
00437     cpl_msg_info(fctid, "Bias Areas: Using %s",
00438                  cx_string_get(results->limits));
00439 
00440 
00441     /*
00442      * Processing
00443      */
00444 
00445     /*
00446      *  Perform Sigma Clipping...
00447      */
00448 
00449     cpl_msg_info(fctid, "Sigma Clipping : Start");
00450 
00451     naccepted = ntotal;
00452 
00453     while ((naccepted > 0) && (curriter < numiter) &&
00454            (currfraction > maxfraction)) {
00455 
00456         results->mean  = cpl_matrix_get_mean(matrix_zz1);
00457         results->sigma = giraffe_matrix_sigma_mean(matrix_zz1, results->mean);
00458 
00459         for (k = 0; k < cpl_matrix_get_nrow(matrix_zz1); k++) {
00460             cpl_matrix_set(matrix_zz1diff, k, 0,
00461                            cpl_matrix_get(matrix_zz1diff, k, 0) -
00462                            results->mean);
00463         }
00464 
00465         cpl_msg_info(fctid, "bias[%d]: mean = %5g, sigma = %6.3g, "
00466                      "ratio = %6.3g, accepted = %ld\n", curriter,
00467                      results->mean, results->sigma, currfraction, naccepted);
00468 
00469         pixcount = 0L;
00470         for (l = 0; l < cpl_matrix_get_nrow(matrix_zz1); l++) {
00471             if (fabs(cpl_matrix_get(matrix_zz1diff, l, 0)) >
00472                 (results->sigma * sigma)) {
00473                 continue;
00474             }
00475 
00476             cpl_matrix_set(matrix_zz1, pixcount, 0,
00477                            cpl_matrix_get(matrix_zz1, l, 0));
00478             pixcount++;
00479         }
00480 
00481         cpl_matrix_resize(matrix_zz1, 0, 0, 0, pixcount -
00482                           cpl_matrix_get_nrow(matrix_zz1));
00483 
00484         cpl_matrix_resize(matrix_zz1diff, 0, 0, 0, pixcount -
00485                           cpl_matrix_get_nrow(matrix_zz1diff));
00486 
00487         if (pixcount == naccepted) {
00488             break;
00489         }
00490 
00491         naccepted = pixcount;
00492 
00493         currfraction = (cxdouble) naccepted / (cxdouble) ntotal;
00494         curriter++;
00495     }
00496 
00497     cpl_msg_info(fctid, "Sigma Clipping : End");
00498 
00499 
00500     /*
00501      * Postprocessing
00502      */
00503 
00504     results->mean = cpl_matrix_get_mean(matrix_zz1);
00505     results->sigma = giraffe_matrix_sigma_mean(matrix_zz1, results->mean);
00506 
00507     cpl_msg_info(fctid, "Sigma Clipping Results : bias[%d]: mean = %5g, "
00508                  "sigma = %6.3g, ratio = %6.3g, accepted=%ld\n", curriter,
00509                  results->mean, results->sigma, currfraction, naccepted);
00510 
00511 
00512     /*
00513      * Cleanup
00514      */
00515 
00516     if (matrix_zz1 != NULL) {
00517         cpl_matrix_delete(matrix_zz1);
00518     }
00519 
00520     if (matrix_zz1diff != NULL) {
00521         cpl_matrix_delete(matrix_zz1diff);
00522     }
00523 
00524     return EXIT_SUCCESS;
00525 
00526 }
00527 
00528 
00529 /*
00530  * @brief
00531  *   Compute bias plane over bias areas.
00532  *
00533  * @param  image        Image over which to calculate mean and sigma
00534  *                      inside bias areas
00535  * @param  areas        Bias areas specifications [N,4] as a matrix
00536  *                      where N is the number of areas
00537  * @param  sigma        Sigma clipping algorithm: sigma value
00538  * @param  numiter      Sigma clipping algorithm: max number of iterations
00539  * @param  maxfraction  Sigma clipping algorithm: max. fraction of pixels
00540  *
00541  * @param  results      A @em GiBiasResults structure which contains
00542  *                      the computed coefficients of the fitted bias plane,
00543  *                      their standard deviation and the coordinates
00544  *                      specifications of the areas used on return.
00545  *
00546  * @returns EXIT_SUCCESS on sucess, negative value otherwise
00547  *
00548  * Fit a plane through the pixels value over areas defined by @em areas.
00549  * @em areas is a matrix [N, 4] specifying the lower left and upper right
00550  * corners of each of the @em N areas of the CCD used for bias determination.
00551  *
00552  * @code
00553  *   areas = [[xmin1, xmax1, ymin1, ymax1],
00554  *            [xmin2, xmax2, ymin2, ymax2],
00555  *             ....., ....., ....., .....
00556  *            [xminN, xmaxN, yminN, ymaxN]]
00557  * @endcode
00558  *
00559  * Areas pixels whose value is too far from mean average are rejected by
00560  * kappa-sigma clipping defined by parameters sigma, numiter and maxfraction.
00561  * This rejection loops until good/bad points ratio is enough or the
00562  * maximum number of iteration has been reached.
00563  *
00564  * The function returns the bias fitted plane coefficients @em BP1, BP2, BP3:
00565  * @code
00566  *   [BiasPlane] = [img] * inv([1XY]) over the bias areas where:
00567  *   [BiasPlane] = [BP1,BP2,BP3] coefficients of plane z = BP1+BP2*x+BP3*y
00568  *   [Z] = 1-column matrix of pixel values
00569  *   [1XY] = 3-columns matrix of all pixels coordinates [1,X,Y]
00570  * @endcode
00571  *
00572  * the bias sigma computed over the kappa sigma clipped areas pixels value
00573  * and a string specifying the areas used as follows
00574  * @code
00575  *   "xmin1:xmax1:ymin1:ymax1;...;xminN:xmaxN:yminN:ymaxN"
00576  * @endcode
00577  *
00578  * @see giraffe_bias_compute_mean, giraffe_compute_remove_bias
00579  */
00580 
00581 inline static cxint
00582 _giraffe_bias_compute_plane(cpl_image *image, cpl_matrix *areas,
00583                             cxdouble sigma, cxint numiter,
00584                             cxdouble maxfraction, GiBiasResults *results)
00585 {
00586 
00587     const cxchar *fctid = "giraffe_bias_compute_plane";
00588 
00589     cxbool not_ok = FALSE;
00590 
00591     cxint32 j, l, k;
00592     cxint32 img_dimx, img_dimy;
00593     cxint32 ba_num_cols, ba_num_rows;
00594     cxint32 curriter = 0;
00595     cxint32 x0, x1, x2, x3;
00596 
00597     cxlong ntotal   = 0L;
00598     cxlong naccepted  = 0L;
00599     cxlong n                  = 0L;
00600     cxlong pixcount           = 0L;
00601     cxlong kindex;
00602 
00603     cxdouble *pdimg;
00604     cxdouble currfraction = 2.0;
00605     cxdouble d_fit_sigma;
00606 
00607     cx_string *tmp = NULL;
00608 
00609     cpl_matrix *matrix_xbs = NULL;
00610     cpl_matrix *matrix_ybs = NULL;
00611     cpl_matrix *matrix_zbs = NULL;
00612     cpl_matrix *matrix_coeff = NULL;
00613     cpl_matrix *matrix_base = NULL;
00614     cpl_matrix *matrix_fit = NULL;
00615 
00616 
00617     /*
00618      * Initialization
00619      */
00620 
00621     if (results->limits == NULL) {
00622         cpl_msg_info(fctid, "Unable to store biaslimits return parameter, "
00623                      "aborting...");
00624         return -3;
00625     }
00626 
00627     if (cpl_image_get_type(image) != CPL_TYPE_DOUBLE) {
00628         cpl_msg_info(fctid, "Only image allowed of type double, aborting ...");
00629         return -3;
00630     }
00631 
00632 
00633     /*
00634      * Preprocessing
00635      */
00636 
00637     /*
00638      *  Validate Bias Areas and calculate number of points
00639      */
00640 
00641     if (areas == NULL) {
00642         cpl_msg_info(fctid, "Bias Areas: Missing bias areas, "
00643                      "aborting ...");
00644         return -1;
00645     }
00646 
00647     ba_num_cols = cpl_matrix_get_ncol(areas);
00648     ba_num_rows = cpl_matrix_get_nrow(areas);
00649 
00650     for (j = 0; j < ba_num_rows; j++) {
00651         x3 = (cxint32) cpl_matrix_get(areas, j, 3);
00652         x2 = (cxint32) cpl_matrix_get(areas, j, 2);
00653         x1 = (cxint32) cpl_matrix_get(areas, j, 1);
00654         x0 = (cxint32) cpl_matrix_get(areas, j, 0);
00655 
00656         ntotal += (cxulong) ((x3 - x2) * (x1 - x0));
00657     }
00658 
00659     if (ntotal <= 0) {
00660         cpl_msg_info(fctid, "Bias Areas: Inconsistent specification, "
00661                      "aborting ...");
00662         return -1;
00663     }
00664 
00665     matrix_xbs = cpl_matrix_new(ntotal, 1);
00666     matrix_ybs = cpl_matrix_new(ntotal, 1);
00667     matrix_zbs = cpl_matrix_new(1, ntotal);
00668 
00669     img_dimx = cpl_image_get_size_x(image);
00670     img_dimy = cpl_image_get_size_y(image);
00671 
00672     cpl_msg_info(fctid, "Bias Areas: specified are %ld points in %dx%d "
00673                  "image", ntotal, img_dimx, img_dimy);
00674 
00675 
00676     /*
00677      *  Compute Bias Areas Values
00678      */
00679 
00680     cx_string_set(results->limits, "");
00681 
00682     for (j = 0; j < ba_num_rows; j++) {
00683 
00684         x3 = (cxint32) cpl_matrix_get(areas, j, 3);
00685         x2 = (cxint32) cpl_matrix_get(areas, j, 2);
00686         x1 = (cxint32) cpl_matrix_get(areas, j, 1);
00687         x0 = (cxint32) cpl_matrix_get(areas, j, 0);
00688 
00689         not_ok = ((x0 > img_dimx) || (x1 > img_dimx) ||
00690                   (x2 > img_dimy) || (x3 > img_dimy));
00691 
00692         if (not_ok) {
00693             continue;
00694         }
00695 
00696         tmp = cx_string_new();
00697 
00698         cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
00699         cx_string_append(results->limits, cx_string_get(tmp));
00700 
00701         cx_string_delete(tmp);
00702         tmp = NULL;
00703 
00704         pdimg = cpl_image_get_data_double(image);
00705 
00706         for (k = x0; k < x1; k++) {
00707             for (l = x2; l < x3; l++) {
00708                 cpl_matrix_set(matrix_xbs, n, 0, (cxdouble)k);
00709                 cpl_matrix_set(matrix_ybs, n, 0, (cxdouble)l);
00710                 cpl_matrix_set(matrix_zbs, 0, n, pdimg[l * img_dimx + k]);
00711                 n++;
00712             }
00713         }
00714     }
00715 
00716     if (n != ntotal) {
00717         cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting ...");
00718 
00719         cpl_matrix_delete(matrix_xbs);
00720         cpl_matrix_delete(matrix_ybs);
00721         cpl_matrix_delete(matrix_zbs);
00722 
00723         return -3;
00724     }
00725 
00726     cpl_msg_info(fctid, "Bias Areas: Using %s [%ld pixels]",
00727                  cx_string_get(results->limits), ntotal);
00728 
00729 
00730     /*
00731      * Processing
00732      */
00733 
00734     cpl_msg_info(fctid, "Sigma Clipping : Start");
00735 
00736 
00737     /*
00738      *  Perform Sigma Clipping
00739      */
00740 
00741     naccepted = ntotal;
00742 
00743     while ((naccepted > 0) && (curriter < numiter) &&
00744            (currfraction > maxfraction)) {
00745 
00746         matrix_base = cpl_matrix_new(3, naccepted);
00747 
00748         for (kindex = 0; kindex < naccepted; kindex++) {
00749 
00750             cpl_matrix_set(matrix_base, 0, kindex, 1.0);
00751             cpl_matrix_set(matrix_base, 1, kindex,
00752                            cpl_matrix_get(matrix_xbs, kindex, 0));
00753             cpl_matrix_set(matrix_base, 2, kindex,
00754                            cpl_matrix_get(matrix_ybs, kindex, 0));
00755 
00756         }
00757 
00758         if (matrix_coeff != NULL) {
00759             cpl_matrix_delete(matrix_coeff);
00760             matrix_coeff = NULL;
00761         }
00762 
00763         matrix_coeff = giraffe_matrix_leastsq(matrix_base, matrix_zbs);
00764 
00765         if (matrix_coeff == NULL) {
00766 
00767             cpl_msg_info(fctid, "Sigma Clipping : Error in least square "
00768                          "solution, aborting...");
00769 
00770             cpl_matrix_delete(matrix_base);
00771             cpl_matrix_delete(matrix_xbs);
00772             cpl_matrix_delete(matrix_ybs);
00773             cpl_matrix_delete(matrix_zbs);
00774 
00775             return -2;
00776 
00777         }
00778 
00779         matrix_fit = cpl_matrix_product_create(matrix_coeff, matrix_base);
00780 
00781         if (matrix_base != NULL) {
00782             cpl_matrix_delete(matrix_base);
00783             matrix_base = NULL;
00784         }
00785 
00786         d_fit_sigma = giraffe_matrix_sigma_fit(matrix_zbs, matrix_fit);
00787 
00788         cpl_msg_info(fctid, "Sigma Clipping : bias plane[%d]: %g + "
00789                      "%g * x + %g * y, sigma = %8.5g, ratio = %7.4g, "
00790                      "accepted = %ld\n", curriter,
00791                      cpl_matrix_get(matrix_coeff, 0, 0),
00792                      cpl_matrix_get(matrix_coeff, 0, 1),
00793                      cpl_matrix_get(matrix_coeff, 0, 2),
00794                      d_fit_sigma, currfraction, naccepted);
00795 
00796 
00797         /*
00798          *  Perform Clipping
00799          */
00800 
00801         pixcount = 0;
00802         for (l = 0;l < cpl_matrix_get_ncol(matrix_zbs); l++) {
00803             if (fabs(cpl_matrix_get(matrix_fit, 0, l) -
00804                      cpl_matrix_get(matrix_zbs, 0, l)) >
00805                 (d_fit_sigma * sigma)) {
00806                 continue;
00807             }
00808 
00809             cpl_matrix_set(matrix_xbs, pixcount, 0,
00810                            cpl_matrix_get(matrix_xbs, l, 0));
00811 
00812             cpl_matrix_set(matrix_ybs, pixcount, 0,
00813                            cpl_matrix_get(matrix_ybs, l, 0));
00814 
00815             cpl_matrix_set(matrix_zbs, 0, pixcount,
00816                            cpl_matrix_get(matrix_zbs, 0, l));
00817 
00818             pixcount++;
00819         }
00820 
00821         cpl_matrix_resize(matrix_xbs, 0, pixcount -
00822                           cpl_matrix_get_ncol(matrix_xbs), 0, 0);
00823 
00824         cpl_matrix_resize(matrix_ybs, 0, pixcount -
00825                           cpl_matrix_get_ncol(matrix_ybs), 0, 0);
00826 
00827         cpl_matrix_resize(matrix_zbs, 0, 0, 0, pixcount -
00828                           cpl_matrix_get_ncol(matrix_zbs));
00829 
00830         if (matrix_fit != NULL)  {
00831             cpl_matrix_delete(matrix_fit);
00832             matrix_fit = NULL;
00833         }
00834 
00835         if (pixcount == naccepted) {
00836             break;
00837         }
00838 
00839         naccepted = pixcount;
00840 
00841         currfraction = (cxdouble)naccepted / (cxdouble)ntotal;
00842         curriter++;
00843     }
00844 
00845     cpl_msg_info(fctid, "Sigma Clipping : End");
00846 
00847 
00848     /*
00849      * Postprocessing
00850      */
00851 
00852     if (results->coeffs != NULL) {
00853         cpl_matrix_delete(results->coeffs);
00854     }
00855 
00856     results->coeffs = matrix_coeff;
00857     results->sigma = cpl_matrix_get_stdev(matrix_zbs);
00858 
00859     cpl_msg_info(fctid, "Sigma Clipping Results (%d/%ld, sigma = %g)",
00860                  curriter, naccepted, results->sigma);
00861 
00862 
00863     /*
00864      * Cleanup
00865      */
00866 
00867     if (matrix_xbs != NULL) {
00868         cpl_matrix_delete(matrix_xbs);
00869     }
00870 
00871     if (matrix_ybs != NULL) {
00872         cpl_matrix_delete(matrix_ybs  );
00873     }
00874 
00875     if (matrix_zbs != NULL) {
00876         cpl_matrix_delete(matrix_zbs  );
00877     }
00878 
00879     return EXIT_SUCCESS;
00880 
00881 }
00882 
00883 
00884 /*
00885  * @brief
00886  *    Compute bias curve over bias areas.
00887  *
00888  * @param image        Image over which to calculate bias curve
00889  *                     inside bias areas
00890  * @param areas        bias areas specifications [N,4] as a matrix
00891  *                     where N is the number of areas
00892  * @param sigma        sigma clipping algorithm: sigma value
00893  * @param numiter      sigma clipping algorithm: max number of iterations
00894  * @param maxfraction  sigma clipping algorithm: max. fraction of pixels
00895  * @param xdeg         Polynom order for fit in x direction
00896  * @param ydeg         Polynom order for fit in y direction
00897  * @param xstep        Step in x direction for polynomial fit
00898  * @param ystep        Step in y direction for polynomial fit
00899  *
00900  * @param results      computed bias curve over bias areas coordinates
00901  *                     specifications
00902  *
00903  * @returns EXIT_SUCCESS on sucess, negative value otherwise
00904  *
00905  * Fit a surface using a rectangular mesh defined by @em xstep and
00906  * @em ystep in @em areas.
00907  * Mesh pixels whose value is too far from fitted surface are rejected by
00908  * kappa-sigma clipping defined by parameters of @em sigma.
00909  * This rejection loops until good/bad points ratio is enough
00910  * (maxfraction) or the maximum number of iteration (numiter) is reached.
00911  *
00912  * The function returns the bias fitted surface coefficients
00913  * @a results->biascoeffs computed using a 2D Chebyshev polynomial fit
00914  * whose X and Y order are given respectively by @a xdeg and  @a ydeg,
00915  * The bias sigma @a results->biassigma computed over the
00916  * kappa-sigma clipped pixels value and a string specifying the areas
00917  * used as follows:
00918  * @code
00919  *   "xstart:xend:xstep;ymin,ystart,ystep"
00920  * @endcode
00921  *
00922  * @see giraffe_bias_compute_plane, giraffe_bias_compute_remove_bias
00923  */
00924 
00925 inline static cxint
00926 _giraffe_bias_compute_curve(cpl_image *image, cpl_matrix *areas,
00927                             cxdouble sigma, cxint numiter,
00928                             cxdouble maxfraction, cxdouble xdeg, cxdouble ydeg,
00929                             cxdouble xstep, cxdouble ystep,
00930                             GiBiasResults *results)
00931 {
00932 
00933     const cxchar *fctid = "giraffe_bias_compute_curve";
00934 
00935     cxbool not_ok = FALSE;
00936 
00937     cxint32 j, l, k;
00938     cxint32 img_dimx, img_dimy;
00939     cxint32 ba_num_cols, ba_num_rows;
00940     cxint32 curriter = 0;
00941     cxint32 x0, x1, x2, x3;
00942 
00943     cxlong ntotal    = 0L;
00944     cxlong naccepted = 0L;
00945     cxlong n         = 0L;
00946     cxlong pixcount  = 0L;
00947 
00948     cxdouble *pdimg;
00949     cxdouble currfraction = 2.0;
00950     cxdouble d_fit_sigma;
00951 
00952     cx_string *tmp = NULL;
00953 
00954     cpl_matrix *matrix_xbs = NULL;
00955     cpl_matrix *matrix_ybs = NULL;
00956     cpl_matrix *matrix_zbs = NULL;
00957     cpl_matrix *matrix_base = NULL;
00958     cpl_matrix *matrix_fit = NULL;
00959     cpl_matrix *matrix_chebyshev = NULL;
00960 
00961 
00962     /*
00963      * Initialization
00964      */
00965 
00966     if (results->limits == NULL) {
00967         cpl_msg_info(fctid, "Unable to store biaslimits return parameter, "
00968                      "aborting...");
00969         return -3;
00970     }
00971 
00972     if (cpl_image_get_type(image) != CPL_TYPE_DOUBLE) {
00973         cpl_msg_info(fctid, "Only images allowed of type double, aborting...");
00974         return -3;
00975     }
00976 
00977 
00978     /*
00979      * Preprocessing
00980      */
00981 
00982     /*
00983      *  Validate Bias Areas and calculate number of points
00984      */
00985 
00986     if (areas == NULL) {
00987         cpl_msg_info(fctid, "Bias Areas: Missing bias areas, aborting...");
00988         return -1;
00989     }
00990 
00991     ba_num_cols = cpl_matrix_get_ncol(areas);
00992     ba_num_rows = cpl_matrix_get_nrow(areas);
00993 
00994     for (j = 0; j < ba_num_rows; j++) {
00995 
00996         x3 = (cxint32) cpl_matrix_get(areas, j, 3);
00997         x2 = (cxint32) cpl_matrix_get(areas, j, 2);
00998         x1 = (cxint32) cpl_matrix_get(areas, j, 1);
00999         x0 = (cxint32) cpl_matrix_get(areas, j, 0);
01000 
01001         ntotal += (cxulong) (ceil((1.0 + x3 - x2) / ystep) *
01002                              ceil((1.0 + x1 - x0) / xstep));
01003     }
01004 
01005     matrix_xbs = cpl_matrix_new(ntotal, 1);
01006     matrix_ybs = cpl_matrix_new(ntotal, 1);
01007     matrix_zbs = cpl_matrix_new(1, ntotal);
01008 
01009     img_dimx = cpl_image_get_size_x(image);
01010     img_dimy = cpl_image_get_size_y(image);
01011 
01012     cpl_msg_info(fctid, "Bias Areas: Found %ld points in %dx%d image",
01013                  ntotal, img_dimx, img_dimy);
01014 
01015 
01016     /*
01017      *  Compute Bias Areas Values
01018      */
01019 
01020     results->mean = 0.0;
01021 
01022     cx_string_set(results->limits, "");
01023 
01024     for (j = 0; j < ba_num_rows; j++) {
01025 
01026         x3 = (cxint32) cpl_matrix_get(areas, j, 3);
01027         x2 = (cxint32) cpl_matrix_get(areas, j, 2);
01028         x1 = (cxint32) cpl_matrix_get(areas, j, 1);
01029         x0 = (cxint32) cpl_matrix_get(areas, j, 0);
01030 
01031         not_ok = ((x0 > img_dimx) || (x1 > img_dimx) ||
01032                   (x2 > img_dimy) || (x3 > img_dimy));
01033 
01034         if (not_ok) {
01035             continue;
01036         }
01037 
01038         tmp = cx_string_new();
01039 
01040         cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
01041         cx_string_append(results->limits, cx_string_get(tmp));
01042 
01043         cx_string_delete(tmp);
01044         tmp = NULL;
01045 
01046         pdimg = cpl_image_get_data_double(image);
01047 
01048         for (k = x2; k < x3; k += xstep) {
01049             for (l = x0; l < x1; l += ystep) {
01050                 cpl_matrix_set(matrix_xbs, n, 0, (cxdouble)k);
01051                 cpl_matrix_set(matrix_ybs, n, 0, (cxdouble)l);
01052                 cpl_matrix_set(matrix_zbs, 0, n, pdimg[k * img_dimx + l]);
01053                 n++;
01054             }
01055         }
01056     }
01057 
01058     cpl_matrix_resize(matrix_xbs, 0, 0, 0, n - ntotal);
01059     cpl_matrix_resize(matrix_ybs, 0, 0, 0, n - ntotal);
01060     cpl_matrix_resize(matrix_zbs, n - ntotal, 0, 0, 0);
01061 
01062     ntotal = n;
01063 
01064     if (ntotal <= 0) {
01065         cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting...");
01066 
01067         cpl_matrix_delete(matrix_xbs);
01068         cpl_matrix_delete(matrix_ybs);
01069         cpl_matrix_delete(matrix_zbs);
01070 
01071         return -1;
01072     }
01073 
01074     cpl_msg_info(fctid, "Bias Areas: Using %s [%ld pixels]",
01075                  cx_string_get(results->limits), ntotal);
01076 
01077 
01078     /*
01079      * Processing
01080      */
01081 
01082     /*
01083      *  Perform Sigma Clipping
01084      */
01085 
01086     cpl_msg_info(fctid, "Sigma Clipping : Start");
01087 
01088     naccepted = ntotal;
01089 
01090     while ((naccepted > 0) && (curriter < numiter) &&
01091            (currfraction > maxfraction)) {
01092 
01093         matrix_base = giraffe_chebyshev_base2d(0.0, 0.0, img_dimy, img_dimx,
01094                                                xdeg, ydeg, matrix_xbs,
01095                                                matrix_ybs);
01096 
01097         if (matrix_chebyshev != NULL) {
01098             cpl_matrix_delete(matrix_chebyshev);
01099             matrix_chebyshev = NULL;
01100         }
01101 
01102         matrix_chebyshev = giraffe_matrix_leastsq(matrix_base, matrix_zbs);
01103 
01104         if (matrix_base != NULL) {
01105             cpl_matrix_delete(matrix_base);
01106             matrix_base = NULL;
01107         }
01108 
01109         if (matrix_chebyshev == NULL) {
01110 
01111             cpl_msg_info(fctid, "Sigma Clipping : Error in least square "
01112                          "solution, aborting...");
01113 
01114             cpl_matrix_delete(matrix_base);
01115             cpl_matrix_delete(matrix_xbs);
01116             cpl_matrix_delete(matrix_ybs);
01117             cpl_matrix_delete(matrix_zbs);
01118 
01119             return -2;
01120 
01121         }
01122 
01123         cpl_matrix_resize(matrix_chebyshev,
01124                           xdeg - cpl_matrix_get_nrow(matrix_chebyshev), 0, 0,
01125                           ydeg - cpl_matrix_get_ncol(matrix_chebyshev));
01126 
01127         matrix_fit = giraffe_chebyshev_fit2d(0.0, 0.0, img_dimy, img_dimx,
01128                                              matrix_chebyshev, matrix_xbs,
01129                                              matrix_ybs);
01130 
01131         results->mean = cpl_matrix_get_mean(matrix_fit);
01132 
01133         d_fit_sigma = giraffe_matrix_sigma_fit(matrix_zbs, matrix_fit);
01134 
01135         cpl_msg_info(fctid, "Sigma Clipping : bias surface[%d]: "
01136                      "sigma = %8.5g, ratio = %7.4g, accepted = %ld\n",
01137                      curriter, d_fit_sigma, currfraction, naccepted);
01138 
01139 
01140         /*
01141          *  Perform Clipping
01142          */
01143 
01144         pixcount = 0;
01145         for (l = 0; l < cpl_matrix_get_ncol(matrix_zbs); l++) {
01146 
01147             if (fabs(cpl_matrix_get(matrix_fit, l, 0) -
01148                      cpl_matrix_get(matrix_zbs, 0, l)) >
01149                 (d_fit_sigma * sigma)) {
01150                 continue;
01151             }
01152 
01153             cpl_matrix_set(matrix_xbs, pixcount, 0,
01154                            cpl_matrix_get(matrix_xbs, l, 0));
01155 
01156             cpl_matrix_set(matrix_ybs, pixcount, 0,
01157                            cpl_matrix_get(matrix_ybs, l, 0));
01158 
01159             cpl_matrix_set(matrix_zbs, 0, pixcount,
01160                            cpl_matrix_get(matrix_zbs, 0, l));
01161 
01162             pixcount++;
01163         }
01164 
01165         cpl_matrix_resize(matrix_xbs, pixcount -
01166                           cpl_matrix_get_nrow(matrix_xbs), 0, 0, 0);
01167 
01168         cpl_matrix_resize(matrix_ybs, pixcount -
01169                           cpl_matrix_get_nrow(matrix_ybs), 0, 0, 0);
01170 
01171         cpl_matrix_resize(matrix_zbs, 0, 0, 0, pixcount -
01172                           cpl_matrix_get_ncol(matrix_zbs));
01173 
01174         if (matrix_fit != NULL) {
01175             cpl_matrix_delete(matrix_fit);
01176             matrix_fit = NULL;
01177         }
01178 
01179         if (pixcount == naccepted) {
01180             break;
01181         }
01182 
01183         naccepted = pixcount;
01184 
01185         currfraction = (cxdouble)naccepted / (cxdouble)ntotal;
01186         curriter++;
01187     }
01188 
01189     cpl_msg_info(fctid, "Sigma Clipping : End");
01190 
01191 
01192     /*
01193      * Postprocessing
01194      */
01195 
01196     if (results->coeffs != NULL) {
01197         cpl_matrix_delete(results->coeffs);
01198     }
01199 
01200     results->coeffs = matrix_chebyshev;
01201     results->sigma  = cpl_matrix_get_stdev(matrix_zbs);
01202 
01203     cpl_msg_info(fctid, "Sigma Clipping Results (%d/%ld, sigma = %g)",
01204                  curriter, naccepted, results->sigma);
01205 
01206 
01207     /*
01208      * Cleanup
01209      */
01210 
01211     if (matrix_xbs != NULL) {
01212         cpl_matrix_delete(matrix_xbs);
01213     }
01214 
01215     if (matrix_ybs != NULL) {
01216         cpl_matrix_delete(matrix_ybs);
01217     }
01218 
01219     if (matrix_zbs != NULL) {
01220         cpl_matrix_delete(matrix_zbs);
01221     }
01222 
01223     return EXIT_SUCCESS;
01224 
01225 }
01226 
01227 
01228 /*
01229  * @brief
01230  *   Actual Calculation/Logic Module of Bias Removal
01231  *
01232  * @param img            Input image
01233  * @param biasareas      Bias area to use
01234  * @param master_bias    Master bias Image
01235  * @param bad_pixel_map  Bad pixel map Image
01236  * @param method         Method to use for bias removal
01237  * @param option         Option to use for bias removal
01238  * @param model          Model to use for bias removal
01239  * @param remove_bias    Remove bias from result image
01240  * @param mbias          Master Bias value read from Master Bias
01241  *                       FITS Keywords
01242  * @param xdeg           X Degree of Polynomial of Chebyshev fit
01243  * @param ydeg           Y Degree of Polynomial of Chebyshev fit
01244  * @param xstep          X Step to use during Chebyshev fit
01245  * @param ystep          Y Step to use during Chebyshev fit
01246  * @param sigma          Sigma to use during Kappa Sigma Clipping
01247  * @param numiter        Max Num Iterations to use during Kappa Sigma
01248  *                       Clipping
01249  * @param maxfraction    Maximum Fraction to use during Kappa Sigma Clipping
01250  * @param results        Values determined
01251  *
01252  * @return The function returns 0 on success, or an error code otherwise.
01253  *
01254  * TBD
01255  */
01256 
01257 inline static cxint
01258 _giraffe_bias_compute(cpl_image *img, cpl_matrix *biasareas,
01259                       cpl_image *master_bias, cpl_image *bad_pixel_map,
01260                       GiBiasMethod method, GiBiasOption option,
01261                       GiBiasModel model, cxbool remove_bias, cxdouble mbias,
01262                       cxdouble xdeg, cxdouble ydeg, cxdouble xstep,
01263                       cxdouble ystep, cxdouble sigma, cxint numiter,
01264                       cxdouble maxfraction, GiBiasResults *results)
01265 {
01266 
01267     const cxchar *fctid = "giraffe_bias_compute";
01268 
01269     cxint32 i, j;
01270     cxint32 img_dimx, img_dimy;
01271     cxint32 img_centerx, img_centery;
01272     cxint32 master_bias_dimx, master_bias_dimy;
01273     cxint32 bad_pixel_dimx, bad_pixel_dimy;
01274     cxint32 *pibadpixelmap;
01275 
01276     cxint32 error_code = 0;
01277 
01278     cxulong jj;
01279 
01280     cxdouble *pdimg, *pdmbias, biasdrift;
01281 
01282     cpl_matrix *tmp_x, *tmp_y, *tmp_fit;
01283     cpl_matrix *matrix_coeff = NULL;              //TODO do we need this????
01284 
01285 
01286     /*
01287      * Initialization
01288      */
01289 
01290     if (results->limits == NULL) {
01291         cpl_msg_info(fctid, "Unable to store bias limits, aborting...");
01292         return -3;
01293     }
01294 
01295     if (cpl_image_get_type(img) != CPL_TYPE_DOUBLE) {
01296         cpl_msg_info(fctid, "Only images of type double allowed, aborting...");
01297         return -3;
01298     }
01299 
01300 
01301     /*
01302      * Preprocessing
01303      */
01304 
01305     img_dimx = cpl_image_get_size_x(img);
01306     img_dimy = cpl_image_get_size_y(img);
01307 
01308     img_centerx = img_dimx >> 1;
01309     img_centery = img_dimy >> 1;
01310 
01311 
01312     /*
01313      *  Validate Bias Areas and calculate number of points
01314      */
01315 
01316     if (biasareas == NULL) {
01317         cpl_msg_info(fctid, "Bias Areas: Missing bias areas, aborting...");
01318         return -1;
01319     }
01320 
01321     cx_string_set(results->limits, "");
01322 
01323     if (remove_bias) {
01324         cpl_msg_info(fctid, "removing bias value\n");
01325     }
01326     else {
01327         cpl_msg_info(fctid, "NOT removing bias value");
01328     }
01329 
01330 
01331     /*
01332      * Processing
01333      */
01334 
01335     if (model == GIBIAS_MODEL_MEAN) {
01336 
01337         /*
01338          *  Compute bias average value over biaslimits areas
01339          */
01340 
01341         cpl_msg_info(fctid, "using MEAN Model...");
01342 
01343         if ((option == GIBIAS_OPTION_PLANE) ||
01344             (method == GIBIAS_METHOD_PLANE)) {
01345             cpl_msg_error(fctid, "Can not use MEAN model with PLANE method.");
01346             return -3;
01347         }
01348 
01349         error_code = _giraffe_bias_compute_mean(img, biasareas, sigma, numiter,
01350                                                 maxfraction, results);
01351 
01352 
01353     }
01354     else if ((method == GIBIAS_METHOD_CURVE) ||
01355              (option == GIBIAS_OPTION_CURVE)) {
01356 
01357         /*
01358          * A 2D Bias Surface is fitted on a grid of the pixels in the
01359          * bias areas
01360          */
01361 
01362         cpl_msg_info(fctid, "using CURVE Model...");
01363 
01364         error_code = _giraffe_bias_compute_curve(img, biasareas, sigma,
01365                                                  numiter, maxfraction, xdeg,
01366                                                  ydeg, xstep, ystep, results);
01367 
01368         if (error_code != EXIT_SUCCESS) {
01369             cpl_msg_error(fctid, "Error during calculation of bias curve, "
01370                           "aborting...");
01371             return error_code;
01372         }
01373 
01374         matrix_coeff = cpl_matrix_duplicate(results->coeffs);
01375 
01376     }
01377     else {
01378 
01379         /*
01380          *  A 2D Bias Plane is fitted on the pixels of the biaslimits area
01381          */
01382 
01383         cpl_msg_info(fctid, "using PLANE (FITTED) model\n");
01384 
01385         error_code = _giraffe_bias_compute_plane(img, biasareas, sigma,
01386                                                  numiter, maxfraction,
01387                                                  results);
01388 
01389         if (error_code == EXIT_SUCCESS) {
01390             matrix_coeff = cpl_matrix_duplicate(results->coeffs);
01391 
01392             results->mean = cpl_matrix_get(matrix_coeff, 0, 0) +
01393                 cpl_matrix_get(matrix_coeff, 0, 1) * img_centerx +
01394                 cpl_matrix_get(matrix_coeff, 0, 2) * img_centery;
01395 
01396         }
01397         else {
01398             cpl_msg_error(fctid, "Error during calculation of bias plane, "
01399                           "aborting...");
01400             return error_code;
01401         }
01402 
01403     }
01404 
01405 
01406     /*
01407      *  Perform Computation based on method....
01408      */
01409 
01410     if (method == GIBIAS_METHOD_UNIFORM) {
01411 
01412         /*
01413          *  METHOD == UNIFORM
01414          */
01415 
01416         cpl_msg_info(fctid, "using bias method == UNIFORM");
01417 
01418         /* we subtract the bias average value */
01419 
01420         if (model == GIBIAS_MODEL_MEAN) {
01421             cpl_msg_info(fctid, "bias value(areas mean) = %f, bias sigma = %f",
01422                          results->mean, results->sigma);
01423         }
01424         else {
01425             cpl_msg_info(fctid, "bias value at (%d, %d) = %f, bias sigma = %f",
01426                          img_centerx, img_centery, results->mean,
01427                          results->sigma);
01428         }
01429 
01430         if (remove_bias == TRUE) {
01431             cpl_image_subtract_scalar(img, results->mean);
01432         }
01433 
01434     }
01435     else if (method == GIBIAS_METHOD_PLANE) {
01436 
01437         /*
01438          *  METHOD == PLANE
01439          */
01440 
01441         cpl_msg_info(fctid, "using bias method == PLANE");
01442 
01443         if (matrix_coeff == NULL) {
01444 
01445             error_code = _giraffe_bias_compute_plane(img, biasareas, sigma,
01446                                                      numiter, maxfraction,
01447                                                      results);
01448 
01449             if (error_code == EXIT_SUCCESS) {
01450                 matrix_coeff = cpl_matrix_duplicate(results->coeffs);
01451 
01452                 results->mean = cpl_matrix_get(matrix_coeff, 0, 0) +
01453                     cpl_matrix_get(matrix_coeff, 0, 1) * img_centerx +
01454                     cpl_matrix_get(matrix_coeff, 0, 2) * img_centery;
01455 
01456             }
01457             else {
01458                 cpl_msg_error(fctid, "Error during calculation of bias "
01459                               "plane, aborting...");
01460                 return error_code;
01461             }
01462         }
01463 
01464         cpl_msg_info(fctid, "bias value at (%d, %d) = %f, bias sigma = %f, "
01465                      "bias plane = %g + %g * x + %g * y", img_centerx,
01466                      img_centery, results->mean, results->sigma,
01467                      cpl_matrix_get(matrix_coeff, 0, 0),
01468                      cpl_matrix_get(matrix_coeff, 0, 1),
01469                      cpl_matrix_get(matrix_coeff, 0, 2));
01470 
01471         if (remove_bias == TRUE) {
01472 
01473             //TODO Move plane removal to gir_image...
01474 
01475             pdimg = cpl_image_get_data_double(img);
01476 
01477             for (j = 0; j < img_dimy; j++) {
01478                 jj = j * img_dimx;
01479                 for (i = 0;i < img_dimx; i++) {
01480                     pdimg[jj + i] = (cpl_matrix_get(matrix_coeff, 0, 0) +
01481                                      cpl_matrix_get(matrix_coeff, 0, 1) * j +
01482                                      cpl_matrix_get(matrix_coeff, 0, 2) * i);
01483                 }
01484             }
01485         }
01486 
01487     }
01488     else if (method == GIBIAS_METHOD_CURVE) {
01489 
01490         /*
01491          *  METHOD == CURVE
01492          */
01493 
01494         cpl_msg_info(fctid, "using bias method == CURVE");
01495 
01496         if (matrix_coeff == NULL) {
01497 
01498             error_code = _giraffe_bias_compute_curve(img, biasareas, sigma,
01499                                                      numiter, maxfraction,
01500                                                      xdeg, ydeg, xstep,ystep,
01501                                                      results);
01502 
01503             if (error_code != EXIT_SUCCESS) {
01504                 cpl_msg_error(fctid, "Error during calculation of bias "
01505                               "surface curve, aborting...");
01506                 return error_code;
01507             }
01508 
01509             matrix_coeff = cpl_matrix_duplicate(results->coeffs);
01510 
01511         }
01512 
01513         cpl_msg_info(fctid, "bias value %f, bias sigma = %f\n",
01514                      results->mean, results->sigma);
01515 
01516         //TODO Should we print the matrix as well????
01517 
01518         if (remove_bias == TRUE) {
01519 
01520             pdimg = cpl_image_get_data_double(img);
01521 
01522             tmp_x = cpl_matrix_new(1, 1);
01523             tmp_y = cpl_matrix_new(1, 1);
01524 
01525             for (j = 0; j < img_dimy; j++) {
01526                 jj = j * img_dimx;
01527                 cpl_matrix_set(tmp_x, 0, 0, j);
01528                 for (i = 0; i < img_dimx; i++) {
01529                     cpl_matrix_set(tmp_y, 0, 0, i);
01530                     tmp_fit = giraffe_chebyshev_fit2d(0.0, 0.0, img_dimx,
01531                                                       img_dimy, matrix_coeff,
01532                                                       tmp_x, tmp_y);
01533 
01534                     pdimg[jj + i] -= cpl_matrix_get(tmp_fit, 0, 0);
01535 
01536                     cpl_matrix_delete(tmp_fit);
01537                 }
01538             }
01539 
01540             cpl_matrix_delete(tmp_x);
01541             cpl_matrix_delete(tmp_y);
01542 
01543         }
01544 
01545     }
01546     else if ((method == GIBIAS_METHOD_MASTER) ||
01547              (method == GIBIAS_METHOD_ZMASTER)) {
01548 
01549         /*
01550          *  METHOD == MASTER || ZMASTER
01551          */
01552 
01553         biasdrift = 0.0;
01554 
01555         if (master_bias == NULL) {
01556             cpl_msg_error(fctid, "Master Bias Frame required with MASTER "
01557                           "method.");
01558             return -5;
01559         }
01560 
01561         master_bias_dimx = cpl_image_get_size_x(master_bias);
01562         master_bias_dimy = cpl_image_get_size_y(master_bias);
01563 
01564         if ((master_bias_dimx != img_dimx) || (master_bias_dimy != img_dimy)) {
01565             cpl_msg_error( fctid, "Invalid Masterbias size should be [%d, %d] "
01566                            "but is [%d, %d].", img_dimx, img_dimy,
01567                            master_bias_dimx, master_bias_dimy);
01568             return -7;
01569         }
01570 
01571         if (matrix_coeff == NULL) {
01572 
01573             if (option == GIBIAS_OPTION_CURVE) {
01574 
01575                 error_code = _giraffe_bias_compute_curve(img, biasareas, sigma,
01576                                                          numiter, maxfraction,
01577                                                          xdeg, ydeg, xstep,
01578                                                          ystep, results);
01579 
01580                 if (error_code != EXIT_SUCCESS) {
01581                     cpl_msg_error(fctid, "Error during calculation of bias "
01582                                   "surface curve, aborting...");
01583                     return error_code;
01584                 }
01585 
01586                 matrix_coeff = cpl_matrix_duplicate(results->coeffs);
01587 
01588             }
01589             else {
01590 
01591                 /*
01592                  *  Default to plane...
01593                  */
01594 
01595                 error_code = _giraffe_bias_compute_plane(img, biasareas, sigma,
01596                                                          numiter,maxfraction,
01597                                                          results);
01598 
01599                 if (error_code == EXIT_SUCCESS) {
01600                     matrix_coeff = cpl_matrix_duplicate(results->coeffs);
01601 
01602                     results->mean = cpl_matrix_get(matrix_coeff, 0, 0) +
01603                         cpl_matrix_get(matrix_coeff, 0, 1) * img_centerx +
01604                         cpl_matrix_get(matrix_coeff, 0, 2) * img_centery;
01605                 }
01606                 else {
01607                     cpl_msg_error(fctid, "Error during calculation of bias "
01608                                   "surface plane, aborting...");
01609                     return error_code;
01610                 }
01611             }
01612         }
01613 
01614         if (method == GIBIAS_METHOD_ZMASTER) {
01615 
01616             /*
01617              * A possible drift of the average bias is compensated using
01618              * the pre/overscans reference value as an additive correction
01619              *
01620              */
01621 
01622             biasdrift = results->mean - mbias;
01623 
01624             cpl_msg_info(fctid, "METHOD == ZMASTER, drift value %g",
01625                          biasdrift);
01626 
01627         }
01628         else {
01629             cpl_msg_info(fctid, "METHOD == MASTER");
01630         }
01631 
01632         /*
01633          *  Write some info into log file...
01634          */
01635 
01636         if (option == GIBIAS_OPTION_CURVE) {
01637 
01638             //TODO dump matrix_coeff to log
01639 
01640             cpl_msg_info(fctid, "OPTION == CURVE, bias mean = %f, bias "
01641                          "sigma = %f", results->mean, results->sigma);
01642 
01643         }
01644         else if (option == GIBIAS_OPTION_PLANE) {
01645             cpl_msg_info(fctid, "OPTION == PLANE, bias plane = %g + "
01646                          "%g * x + %g * y, bias mean at (%d, %d) = %f, "
01647                          "bias sigma = %f\n",
01648                          cpl_matrix_get(matrix_coeff, 0, 0),
01649                          cpl_matrix_get(matrix_coeff, 0, 1),
01650                          cpl_matrix_get(matrix_coeff, 0, 2),
01651                          img_centerx, img_centery, results->mean,
01652                          results->sigma);
01653         }
01654         else {
01655             cpl_msg_info(fctid, "bias mean = %f, bias sigma = %f\n",
01656                          results->mean, results->sigma);
01657         }
01658 
01659 
01660         /*
01661          *  Handle bad pixel bitmap...
01662          */
01663 
01664         if (bad_pixel_map == NULL) {
01665 
01666             /*
01667              *  No bad pixel map was given...
01668              *  Subtract masterbias pixels without modification except
01669              *  pre/ovrscan trimming.
01670              */
01671 
01672             if (remove_bias == TRUE) {
01673                 pdimg = cpl_image_get_data_double(img);
01674                 pdmbias = cpl_image_get_data_double(master_bias);
01675 
01676                 for (j = 0; j < img_dimy; j++) {
01677                     jj = j * img_dimx;
01678                     for (i = 0; i < img_dimx; i++) {
01679                         pdimg[jj + i] -= (pdmbias[jj + i] + biasdrift);
01680                     }
01681                 }
01682             }
01683 
01684         }
01685         else {
01686 
01687             /*
01688              *  Bad pixel map was given...
01689              *  For all pixels not flagged in bad_pixel_map the master_bias
01690              *  pixel values are subtracted, otherwise results->mean or
01691              *  biasmodelimage pixels are subtracted depending on model.
01692              */
01693 
01694             bad_pixel_dimx = cpl_image_get_size_x(bad_pixel_map);
01695             bad_pixel_dimy = cpl_image_get_size_y(bad_pixel_map);
01696 
01697             if ((bad_pixel_dimx != img_dimx) || (bad_pixel_dimy != img_dimy)) {
01698                 cpl_msg_error(fctid, "Invalid Badpixel map size should be "
01699                               "[%d, %d] but is [%d, %d].", img_dimx, img_dimy,
01700                               bad_pixel_dimx, bad_pixel_dimy);
01701                 return -6;
01702             }
01703 
01704             if (remove_bias == TRUE) {
01705 
01706                 if (option == GIBIAS_OPTION_CURVE) {
01707 
01708                     pdimg = cpl_image_get_data_double(img);
01709                     pdmbias = cpl_image_get_data_double(master_bias);
01710                     pibadpixelmap = cpl_image_get_data_int(bad_pixel_map);
01711 
01712                     tmp_x = cpl_matrix_new(1, 1);
01713                     tmp_y = cpl_matrix_new(1, 1);
01714 
01715                     for (j = 0; j < img_dimy; j++) {
01716                         jj = j * img_dimx;
01717                         cpl_matrix_set(tmp_x, 0, 0, j);
01718                         for (i = 0; i < img_dimx; i++) {
01719                             if (!(pibadpixelmap[jj + i] & GIR_M_PIX_SET)) {
01720                                 pdimg[jj + i] -= (pdmbias[jj + i] + biasdrift);
01721                             }
01722                             else {
01723                                 cpl_matrix_set(tmp_y, 0, 0, i);
01724                                 tmp_fit =
01725                                     giraffe_chebyshev_fit2d(0.0, 0.0,
01726                                                             img_dimx, img_dimy,
01727                                                             matrix_coeff,
01728                                                             tmp_x, tmp_y);
01729 
01730                                 pdimg[jj + i] -= cpl_matrix_get(tmp_fit, 0, 0);
01731 
01732                                 cpl_matrix_delete(tmp_fit);
01733 
01734                             }
01735                         }
01736                     }
01737 
01738                     cpl_matrix_delete(tmp_y);
01739                     cpl_matrix_delete(tmp_x);
01740 
01741                 }
01742                 else if (option == GIBIAS_OPTION_PLANE) {
01743 
01744                     pdimg = cpl_image_get_data_double(img);
01745                     pdmbias = cpl_image_get_data_double(master_bias);
01746                     pibadpixelmap = cpl_image_get_data_int(bad_pixel_map);
01747 
01748                     for (j = 0; j < img_dimy; j++) {
01749                         jj = j * img_dimx;
01750                         for (i = 0; i < img_dimx; i++) {
01751                             if (!(pibadpixelmap[jj + i] & GIR_M_PIX_SET)) {
01752                                 pdimg[jj + i] -= (pdmbias[jj + i] + biasdrift);
01753                             }
01754                             else {
01755                                 pdimg[jj + i] -=
01756                                     (cpl_matrix_get(matrix_coeff, 0, 0) +
01757                                      cpl_matrix_get(matrix_coeff, 0, 1) * j +
01758                                      cpl_matrix_get(matrix_coeff, 0, 2) * i);
01759                             }
01760                         }
01761                     }
01762 
01763                 }
01764                 else {
01765 
01766                     cpl_msg_info(fctid, "OPTION==UNIFORM\n");
01767 
01768                     pdimg = cpl_image_get_data_double(img);
01769                     pdmbias = cpl_image_get_data_double(master_bias);
01770                     pibadpixelmap = cpl_image_get_data_int(bad_pixel_map);
01771 
01772                     for (j = 0; j < img_dimy; j++) {
01773                         jj = j * img_dimx;
01774                         for (i = 0; i < img_dimx; i++) {
01775                             if (!(pibadpixelmap[jj + i] & GIR_M_PIX_SET)) {
01776                                 pdimg[jj + i] -= (pdmbias[jj + i] + biasdrift);
01777                             }
01778                             else {
01779                                 pdimg[jj + i] -= (results->mean);
01780                             }
01781                         }
01782                     }
01783 
01784                 }
01785 
01786             }
01787 
01788         }
01789 
01790     }
01791     else {
01792 
01793         cpl_msg_error(fctid, "Invalid method, aborting...");
01794         return -4;
01795 
01796     }
01797 
01798 
01799     /*
01800      * Postprocessing
01801      */
01802 
01803     cpl_msg_info(fctid, "Resulting biaslimits : %s",
01804                  cx_string_get(results->limits));
01805 
01806     /*
01807      * Cleanup
01808      */
01809 
01810     if (matrix_coeff) {
01811         cpl_matrix_delete(matrix_coeff);
01812     }
01813 
01814     return 0;
01815 
01816 }
01817 
01818 
01836 cpl_matrix *
01837 giraffe_get_raw_areas(GiImage *image)
01838 {
01839 
01840     const cxchar *const _id = "giraffe_get_raw_areas";
01841 
01842     cxint32 prescx = 0;
01843     cxint32 prescy = 0;
01844     cxint32 ovrscx = 0;
01845     cxint32 ovrscy = 0;
01846 
01847     cxint32 winx = 0;
01848     cxint32 winy = 0;
01849 
01850     cxint32 tprescx = 0;
01851     cxint32 tprescy = 0;
01852     cxint32 tovrscx = 0;
01853     cxint32 tovrscy = 0;
01854 
01855     cxint32 row = 0;
01856 
01857     cpl_matrix *m_tmp;
01858     cpl_propertylist *properties;
01859 
01860 
01861     properties = giraffe_image_get_properties(image);
01862     if (!properties) {
01863         cpl_error_set(_id, CPL_ERROR_NULL_INPUT);
01864         return NULL;
01865     }
01866 
01867     winx = cpl_propertylist_get_int(properties, GIALIAS_WINX);
01868     winy = cpl_propertylist_get_int(properties, GIALIAS_WINY);
01869 
01870     if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
01871         tprescx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
01872         if (tprescx > 0) {
01873             prescx = tprescx;
01874         }
01875     }
01876     if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
01877         tprescy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
01878         if (tprescy > 0) {
01879             prescy = tprescy;
01880         }
01881     }
01882     if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
01883         tovrscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
01884         if (tovrscx > 0) {
01885             ovrscx = tovrscx;
01886         }
01887     }
01888     if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
01889         tovrscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
01890         if (tovrscy > 0) {
01891             ovrscy = tovrscy;
01892         }
01893     }
01894 
01895 
01896     m_tmp = cpl_matrix_new(1, 4);
01897 
01898     if (prescx > 0) {
01899 
01900         cpl_matrix_set(m_tmp, row, 0, 0.0);
01901         cpl_matrix_set(m_tmp, row, 1, (cxdouble) prescx);
01902         cpl_matrix_set(m_tmp, row, 2, 0.0);
01903         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy);
01904 
01905         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
01906         row++;
01907     }
01908 
01909     if (ovrscx > 0) {
01910 
01911         cpl_matrix_set(m_tmp, row, 0, (cxdouble) winx - ovrscx);
01912         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx);
01913         cpl_matrix_set(m_tmp, row, 2, 0.0);
01914         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy);
01915 
01916         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
01917         row++;
01918     }
01919 
01920     if (!(prescy == 0 || prescx == 0 || ovrscx == 0)) {
01921 
01922         cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx + 1);
01923         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
01924         cpl_matrix_set(m_tmp, row, 2, 0.0);
01925         cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy);
01926 
01927         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
01928         row++;
01929 
01930     }
01931     else if (!(prescy == 0 || prescx == 0)) {
01932 
01933         cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx + 1);
01934         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx);
01935         cpl_matrix_set(m_tmp, row, 2, 0.0);
01936         cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy);
01937 
01938         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
01939         row++;
01940 
01941     }
01942     else if (!(prescy == 0 || ovrscx == 0)) {
01943 
01944         cpl_matrix_set(m_tmp, row, 0, 0.0);
01945         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
01946         cpl_matrix_set(m_tmp, row, 2, 0.0);
01947         cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy);
01948 
01949         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
01950         row++;
01951 
01952     }
01953     else if (prescy > 0) {
01954 
01955         cpl_matrix_set(m_tmp, row, 0, 0.0);
01956         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx);
01957         cpl_matrix_set(m_tmp, row, 2, 0.0);
01958         cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy);
01959 
01960         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
01961         row++;
01962     }
01963 
01964     if (!(ovrscy == 0 || prescx == 0 || ovrscx == 0)) {
01965 
01966         cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx + 1);
01967         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
01968         cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
01969         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy);
01970 
01971         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
01972         row++;
01973 
01974     }
01975     else if (!(ovrscy == 0 || prescx == 0)) {
01976 
01977         cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx + 1);
01978         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx);
01979         cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
01980         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy);
01981 
01982         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
01983         row++;
01984 
01985     }
01986     else if (!(ovrscy == 0 || ovrscx == 0)) {
01987 
01988         cpl_matrix_set(m_tmp, row, 0, 0.0);
01989         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
01990         cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
01991         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy);
01992 
01993         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
01994         row++;
01995 
01996     }
01997     else if (ovrscy > 0) {
01998 
01999         cpl_matrix_set(m_tmp, row, 0, 0.0);
02000         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx);
02001         cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
02002         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy);
02003 
02004         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02005         row++;
02006     }
02007 
02008     cpl_matrix_resize(m_tmp, 0, -1, 0, 0);
02009     row--;
02010 
02011     if (row == 0) {
02012         cpl_matrix_delete(m_tmp);
02013         return NULL;
02014     }
02015 
02016     return m_tmp;
02017 }
02018 
02019 
02038 cxint
02039 giraffe_trim_raw_areas(GiImage *image)
02040 {
02041 
02042     const cxchar *fctid = "giraffe_trim_raw_areas";
02043 
02044     cxint nx, ny;
02045     cxint newlx, newly, newux, newuy;
02046 
02047     cxint ovscx = 0;
02048     cxint ovscy = 0;
02049     cxint prscx = 0;
02050     cxint prscy = 0;
02051 
02052     cpl_propertylist *properties = giraffe_image_get_properties(image);
02053 
02054     cpl_image *_image = giraffe_image_get(image);
02055     cpl_image *tmp;
02056 
02057 
02058     if (!properties) {
02059         cpl_msg_error(fctid, "Image does not contain any properties!");
02060         return 1;
02061     }
02062 
02063     nx = cpl_image_get_size_x(_image);
02064     ny = cpl_image_get_size_y(_image);
02065 
02066     if (cpl_propertylist_has(properties, GIALIAS_NAXIS1)) {
02067         cxint _nx = cpl_propertylist_get_int(properties, GIALIAS_NAXIS1);
02068 
02069         if (_nx != nx) {
02070             cpl_msg_warning(fctid, "Image size (%d) and image property `%s' "
02071                             "(%d) are not consistent! Using image size ...",
02072                             nx, GIALIAS_NAXIS1, _nx);
02073         }
02074     }
02075     else {
02076         cpl_msg_warning(fctid, "Image does not contain any property `%s'. "
02077                         "Using image size (%d)", GIALIAS_NAXIS1, nx);
02078     }
02079 
02080 
02081     if (cpl_propertylist_has(properties, GIALIAS_NAXIS2)) {
02082         cxint _ny = cpl_propertylist_get_int(properties, GIALIAS_NAXIS2);
02083 
02084         if (_ny != ny) {
02085             cpl_msg_warning(fctid, "Image size (%d) and image property `%s' "
02086                             "(%d) are not consistent! Using image size ...",
02087                             ny, GIALIAS_NAXIS2, _ny);
02088         }
02089     }
02090     else {
02091         cpl_msg_warning(fctid, "Image does not contain any property `%s'. "
02092                         "Using image size (%d)", GIALIAS_NAXIS2, ny);
02093     }
02094 
02095     if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
02096         ovscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
02097     }
02098 
02099     if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
02100         ovscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
02101     }
02102 
02103     if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
02104         prscx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
02105     }
02106 
02107     if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
02108         prscy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
02109     }
02110 
02111     // FIXME:
02112     //  Check that the one pixel offset is ok. cpl_image_extract starts
02113     //  counting from 1,1 with right and to boundaries inclusive.
02114 
02115     newlx = prscx + 1;
02116     newly = prscy + 1;
02117     newux = nx - ovscx;
02118     newuy = ny - ovscy;
02119 
02120     tmp = cpl_image_extract(_image, newlx, newly, newux, newuy);
02121 
02122     giraffe_image_set(image, tmp);
02123     cpl_image_delete(tmp);
02124 
02125     _image = giraffe_image_get(image);
02126 
02127     nx = cpl_image_get_size_x(_image);
02128     ny = cpl_image_get_size_y(_image);
02129 
02130 
02131     /*
02132      * Update image properties
02133      */
02134 
02135     cpl_propertylist_set_int(properties, GIALIAS_NAXIS1, nx);
02136     cpl_propertylist_set_int(properties, GIALIAS_NAXIS2, ny);
02137 
02138     if (cpl_propertylist_has(properties, GIALIAS_CRPIX1)) {
02139         cxdouble crpix1 = cpl_propertylist_get_double(properties,
02140             GIALIAS_CRPIX1);
02141 
02142         /*
02143          * For GIRAFFE the reference pixel is defined as
02144          * 1 - width(prescan)
02145          */
02146 
02147         crpix1 += (cxdouble)prscx;
02148         cpl_propertylist_set_double(properties, GIALIAS_CRPIX1, crpix1);
02149     }
02150 
02151     if (cpl_propertylist_has(properties, GIALIAS_CRPIX2)) {
02152         cxdouble crpix2 = cpl_propertylist_get_double(properties,
02153             GIALIAS_CRPIX2);
02154 
02155         crpix2 -= (cxdouble) prscy;
02156         cpl_propertylist_set_double(properties, GIALIAS_CRPIX2, crpix2);
02157     }
02158 
02159     cpl_propertylist_erase(properties, GIALIAS_OVSCX);
02160     cpl_propertylist_erase(properties, GIALIAS_OVSCY);
02161     cpl_propertylist_erase(properties, GIALIAS_PRSCX);
02162     cpl_propertylist_erase(properties, GIALIAS_PRSCY);
02163 
02164     return 0;
02165 
02166 }
02167 
02168 
02230 cxint
02231 giraffe_bias_remove(GiImage *result, GiImage *raw,
02232                     GiImage *master_bias, GiImage *bad_pixels,
02233                     cpl_matrix *biaslimits, GiBiasConfig *config)
02234 {
02235 
02236     const cxchar *fctid = "giraffe_bias_remove";
02237 
02238     cxint32 error_code;
02239 
02240     cxdouble bias_value = 0.;
02241 
02242     cx_string *s;
02243 
02244     cpl_matrix *areas = biaslimits;
02245 
02246     cpl_image *_raw = giraffe_image_get(raw);
02247     cpl_image *_master_bias = giraffe_image_get(master_bias);
02248     cpl_image *_bad_pixels = giraffe_image_get(bad_pixels);
02249     cpl_image *timage;
02250 
02251     cpl_propertylist *properties;
02252 
02253     GiBiasResults  *bias_results;
02254 
02255 
02256     cx_assert(raw != NULL);
02257     cx_assert(config != NULL);
02258 
02259     if (result == NULL) {
02260         return 1;
02261     }
02262 
02263 
02264     /*
02265      * Preprocessing
02266      */
02267 
02268 
02269     if (areas == NULL) {
02270 
02271         cpl_msg_info(fctid, "No Bias Areas specified, using Frame "
02272                      "Pre/Overscan regions.");
02273 
02274         cpl_error_reset();
02275         areas = giraffe_get_raw_areas(raw);
02276         if (areas == NULL) {
02277             if (cpl_error_get_code() == CPL_ERROR_NULL_INPUT) {
02278                 cpl_msg_error(fctid, "Invalid raw image! Image does not "
02279                               "contain any properties!");
02280             }
02281             else {
02282                 cpl_msg_error(fctid, "Invalid raw image! Image does not "
02283                               "contain or has invalid pre- and overscan "
02284                               "properties.");
02285             }
02286 
02287             return 1;
02288         }
02289 
02290     }
02291 
02292 
02293     /*
02294      * Processing
02295      */
02296 
02297 
02298     /*
02299      *  Check whether a masterbias image was specified
02300      *  If so, compare pre/ovrscan keywords/areas
02301      */
02302 
02303     if (master_bias != NULL) {
02304         cpl_propertylist *_properties;
02305         cxbool is_same_size = FALSE;
02306 
02307         is_same_size = _giraffe_compare_overscans(raw, master_bias);
02308 
02309         // FIXME:
02310         //  Fully processed calibrations usually do not have overscans
02311         //  any longer. Instead of failing at this point some dummies
02312         //  could be created.
02313 
02314         if (is_same_size==FALSE) {
02315             cpl_msg_info(fctid, "PRE/OVRSCAN Regions do not match between "
02316                          "RAW Image and Masterbias Image.");
02317 
02318             return 1;
02319         }
02320 
02321         _properties = giraffe_image_get_properties(master_bias);
02322 
02323         if (cpl_propertylist_has(_properties, GIALIAS_BIASVALUE)) {
02324             bias_value = cpl_propertylist_get_double(_properties, GIALIAS_BIASVALUE);
02325         }
02326     }
02327 
02328 
02329     /*
02330      *  Check whether a Bad Pixel Map image was specified
02331      *  If so, compare pre/ovrscan keywords/areas
02332      */
02333 
02334     if (bad_pixels != NULL) {
02335         cxbool is_same_size = FALSE;
02336 
02337         is_same_size = _giraffe_compare_overscans(raw, bad_pixels);
02338 
02339         // FIXME:
02340         //  Fully processed calibrations usually do not have overscans
02341         //  any longer. Instead of failing at this point some dummies
02342         //  could be created.
02343 
02344         if (is_same_size == FALSE) {
02345             cpl_msg_info(fctid, "PRE/OVRSCAN Regions do not match between "
02346                          "RAW Image and Bad Pixel Image.");
02347 
02348             return 1;
02349         }
02350     }
02351 
02352 
02353     /*
02354      *  Create a copy of the raw frame. We are not working on the raw
02355      *  frame directly.
02356      */
02357 
02358     // FIXME:
02359     //   Can we use the raw frame directly? (RP)
02360 
02361     timage = cpl_image_duplicate(_raw);
02362 
02363 
02364     /*
02365      * Initialize the results structure
02366      */
02367 
02368     bias_results = cx_malloc(sizeof(GiBiasResults));
02369 
02370     bias_results->mean   = 0;
02371     bias_results->sigma  = 0;
02372     bias_results->limits = cx_string_new();
02373     bias_results->coeffs = NULL;
02374 
02375 
02376     /*
02377      *  Remove Bias Calculation. The bias corrected image will be
02378      *  stored in timage.
02379      */
02380 
02381     error_code = _giraffe_bias_compute(timage, areas, _master_bias,
02382                                        _bad_pixels, config->method,
02383                                        config->option, config->model,
02384                                        config->remove, bias_value,
02385                                        config->xdeg, config->ydeg,
02386                                        config->xstep, config->ystep,
02387                                        config->sigma, config->iterations,
02388                                        config->fraction, bias_results);
02389 
02390     // FIXME:
02391     //  Posslibly check returned error code and do a dedicated
02392     //  exception handling (RP)
02393 
02394     if (error_code) {
02395         cpl_msg_error(fctid, "Bias computation failed with error %d!",
02396                       error_code);
02397         _giraffe_bias_results_destroy(bias_results);
02398         bias_results = NULL;
02399 
02400         if (!biaslimits) {
02401             cpl_matrix_delete(areas);
02402         }
02403 
02404         return 1;
02405     }
02406 
02407     if (!biaslimits) {
02408         cpl_matrix_delete(areas);
02409     }
02410 
02411 
02412     /*
02413      *  Should we remove Bias or not?
02414      */
02415 
02416     properties = giraffe_image_get_properties(raw);
02417 
02418     if (config->remove == TRUE) {
02419 
02420         giraffe_image_set(result, timage);
02421         cpl_image_delete(timage);
02422 
02423         giraffe_image_set_properties(result, properties);
02424 
02425     }
02426     else {
02427 
02428         cpl_image_delete(timage);
02429 
02430         giraffe_image_set(result, _raw);
02431         giraffe_image_set_properties(result, properties);
02432 
02433     }
02434 
02435 
02436     /*
02437      * Postprocessing
02438      */
02439 
02440     properties = giraffe_image_get_properties(result);
02441 
02442     if (config->remove == TRUE) {
02443 
02444         cpl_propertylist_set_int(properties, GIALIAS_BITPIX, -32);
02445         cpl_propertylist_set_double(properties, GIALIAS_BZERO, 0.);
02446         cpl_propertylist_set_double(properties, GIALIAS_BSCALE, 1.);
02447 
02448         if (cpl_propertylist_has(properties, GIALIAS_GIRFTYPE)) {
02449             cpl_propertylist_set_string(properties,
02450                                         GIALIAS_GIRFTYPE, "BRMIMG");
02451         }
02452         else {
02453             cpl_propertylist_append_string(properties,
02454                                            GIALIAS_GIRFTYPE, "BRMIMG");
02455         }
02456         cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
02457                               "GIRAFFE bias subtracted frame");
02458 
02459         /*
02460          * Remove pre- and overscans
02461          */
02462 
02463         giraffe_trim_raw_areas(result);
02464 
02465 
02466         /*
02467          * Convert bias corrected image to electrons
02468          */
02469 
02470         // FIXME:
02471         //  Is this really needed? Disabled on request of DFO (RP)
02472 #if 0
02473         if (cpl_propertylist_has(properties, GIALIAS_CONAD)) {
02474             cxdouble conad = cpl_propertylist_get_double(properties, GIALIAS_CONAD);
02475 
02476             if (conad > 0.) {
02477                 cpl_image *_result = giraffe_image_get(result);
02478 
02479                 cpl_msg_info(fctid, "Converting bias subtracted frame to "
02480                              "electrons; conversion factor is %.2f e-/ADU",
02481                              conad);
02482                 cpl_image_const_op_local(_result, conad, '*');
02483             }
02484             else {
02485                 cpl_msg_warning(fctid, "Invalid conversion factor %.2f "
02486                                 "e-/ADU! Bias subtracted image will not be "
02487                                 "converted.", conad);
02488             }
02489         }
02490         else {
02491             cpl_msg_info(fctid, "Conversion factor not found. Bias subtracted "
02492                          "image will remain in ADU");
02493         }
02494 #endif
02495     }
02496 
02497     s = cx_string_new();
02498 
02499     _giraffe_method_string(s, config->method, config->option);
02500     cpl_propertylist_append_string(properties, GIALIAS_BIASMETHOD, cx_string_get(s));
02501 
02502     cx_string_delete(s);
02503 
02504     cpl_propertylist_append_double(properties, GIALIAS_BIASVALUE,
02505                             bias_results->mean);
02506     cpl_propertylist_append_double(properties, GIALIAS_BIASSIGMA,
02507                             bias_results->sigma);
02508 
02509     if (bias_results->coeffs) {
02510         cpl_propertylist_append_double(properties, GIALIAS_BCLIPSIGMA,
02511                                 config->sigma);
02512         cpl_propertylist_append_int(properties, GIALIAS_BCLIPNITER,
02513                                 config->iterations);
02514         cpl_propertylist_append_double(properties, GIALIAS_BCLIPMFRAC,
02515                                 config->fraction);
02516     }
02517 
02518     if (bias_results->limits) {
02519         cpl_propertylist_append_string(properties, GIALIAS_BIASAREAS,
02520                                 cx_string_get(bias_results->limits));
02521     }
02522 
02523     if (bias_results->coeffs) {
02524         s = cx_string_new();
02525 
02526         _giraffe_stringify_coefficients(s, bias_results->coeffs);
02527         cpl_propertylist_append_string(properties, GIALIAS_BIASPLANE,
02528                                 cx_string_get(s));
02529 
02530         cx_string_delete(s);
02531     }
02532 
02533 
02534     /*
02535      * Cleanup
02536      */
02537 
02538     _giraffe_bias_results_destroy(bias_results);
02539 
02540     return 0;
02541 
02542 }
02543 
02544 
02555 GiBiasConfig *
02556 giraffe_bias_config_create(cpl_parameterlist *list)
02557 {
02558 
02559     const cxchar *s;
02560     cpl_parameter *p;
02561 
02562     GiBiasConfig *config = NULL;
02563 
02564 
02565     if (!list) {
02566         return NULL;
02567     }
02568 
02569     config = cx_calloc(1, sizeof *config);
02570 
02571 
02572     /*
02573      * Some defaults
02574      */
02575 
02576     config->method = GIBIAS_METHOD_UNDEFINED;
02577     config->option = GIBIAS_OPTION_UNDEFINED;
02578     config->model = GIBIAS_MODEL_UNDEFINED;
02579     config->mbias = 0.;
02580     config->xdeg = 1;
02581     config->ydeg = 1;
02582 
02583 
02584     p = cpl_parameterlist_find(list, "giraffe.biasremoval.remove");
02585     config->remove = cpl_parameter_get_bool(p);
02586 
02587     p = cpl_parameterlist_find(list, "giraffe.biasremoval.method");
02588     s = cpl_parameter_get_string(p);
02589 
02590     if (!strcmp(s, "UNIFORM")) {
02591         config->method = GIBIAS_METHOD_UNIFORM;
02592     }
02593     if (!strcmp(s, "PLANE")) {
02594         config->method = GIBIAS_METHOD_PLANE;
02595     }
02596     if (!strcmp(s, "MASTER")) {
02597         config->method = GIBIAS_METHOD_MASTER;
02598     }
02599     if (!strcmp(s, "ZMASTER")) {
02600         config->method = GIBIAS_METHOD_ZMASTER;
02601     }
02602     if (!strcmp(s, "CURVE")) {
02603         config->method = GIBIAS_METHOD_CURVE;
02604     }
02605     if (!strcmp(s, "MASTER+PLANE")) {
02606         config->method = GIBIAS_METHOD_MASTER;
02607         config->option = GIBIAS_OPTION_PLANE;
02608     }
02609     if (!strcmp(s, "ZMASTER+PLANE")) {
02610         config->method = GIBIAS_METHOD_ZMASTER;
02611         config->option = GIBIAS_OPTION_PLANE;
02612     }
02613     if (!strcmp(s, "MASTER+CURVE")) {
02614         config->method = GIBIAS_METHOD_MASTER;
02615         config->option = GIBIAS_OPTION_CURVE;
02616     }
02617     if (!strcmp(s, "ZMASTER+CURVE")) {
02618         config->method = GIBIAS_METHOD_ZMASTER;
02619         config->option = GIBIAS_OPTION_CURVE;
02620     }
02621 
02622     cx_assert(config->method != GIBIAS_METHOD_UNDEFINED);
02623 
02624     p = cpl_parameterlist_find(list, "giraffe.biasremoval.areas");
02625     s = cpl_parameter_get_string(p);
02626     config->areas = cx_strdup(s);
02627 
02628     p = cpl_parameterlist_find(list, "giraffe.biasremoval.sigma");
02629     config->sigma = cpl_parameter_get_double(p);
02630 
02631     p = cpl_parameterlist_find(list, "giraffe.biasremoval.iterations");
02632     config->iterations = cpl_parameter_get_int(p);
02633 
02634     p = cpl_parameterlist_find(list, "giraffe.biasremoval.fraction");
02635     config->fraction = cpl_parameter_get_double(p);
02636 
02637     if (config->method == GIBIAS_METHOD_CURVE ||
02638         config->option == GIBIAS_OPTION_CURVE) {
02639         p = cpl_parameterlist_find(list, "giraffe.biasremoval.xorder");
02640         config->xdeg = 1 + cpl_parameter_get_int(p);
02641 
02642         p = cpl_parameterlist_find(list, "giraffe.biasremoval.yorder");
02643         config->ydeg = 1 + cpl_parameter_get_int(p);
02644     }
02645 
02646     p = cpl_parameterlist_find(list, "giraffe.biasremoval.xstep");
02647     config->xstep = cpl_parameter_get_int(p);
02648 
02649     p = cpl_parameterlist_find(list, "giraffe.biasremoval.ystep");
02650     config->ystep = cpl_parameter_get_int(p);
02651 
02652     return config;
02653 
02654 }
02655 
02656 
02669 void
02670 giraffe_bias_config_destroy(GiBiasConfig *config)
02671 {
02672 
02673     if (config) {
02674         if (config->areas) {
02675             cx_free(config->areas);
02676             config->areas = NULL;
02677         }
02678 
02679         cx_free(config);
02680     }
02681 
02682     return;
02683 }
02684 
02685 
02697 void
02698 giraffe_bias_config_add(cpl_parameterlist *list)
02699 {
02700 
02701     cpl_parameter *p;
02702 
02703 
02704     if (!list) {
02705         return;
02706     }
02707 
02708     p = cpl_parameter_new_value("giraffe.biasremoval.remove",
02709                                 CPL_TYPE_BOOL,
02710                                 "Enable bias removal",
02711                                 "giraffe.biasremoval",
02712                                 TRUE);
02713     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "remove-bias");
02714     cpl_parameterlist_append(list, p);
02715 
02716 
02717     p = cpl_parameter_new_enum("giraffe.biasremoval.method",
02718                                CPL_TYPE_STRING,
02719                                "Bias removal method",
02720                                "giraffe.biasremoval",
02721                                "MASTER", 9, "UNIFORM", "PLANE", "MASTER",
02722                                "ZMASTER", "MASTER+PLANE", "ZMASTER+PLANE",
02723                                "CURVE", "MASTER+CURVE", "ZMASTER+CURVE");
02724     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-method");
02725     cpl_parameterlist_append(list, p);
02726 
02727 
02728     p = cpl_parameter_new_value("giraffe.biasremoval.areas",
02729                                 CPL_TYPE_STRING,
02730                                 "Bias areas to use "
02731                                 "(Xl0:Xr0:Yl0:Yr0, ... ,Xln:Xrn:Yln:Yrn)",
02732                                 "giraffe.biasremoval",
02733                                 "None");
02734     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-areas");
02735     cpl_parameterlist_append(list, p);
02736 
02737 
02738     p = cpl_parameter_new_value("giraffe.biasremoval.sigma",
02739                                 CPL_TYPE_DOUBLE,
02740                                 "Sigma Clipping: sigma threshold factor",
02741                                 "giraffe.biasremoval",
02742                                 2.5);
02743     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-sigma");
02744     cpl_parameterlist_append(list, p);
02745 
02746 
02747     p = cpl_parameter_new_value("giraffe.biasremoval.iterations",
02748                                 CPL_TYPE_INT,
02749                                 "Sigma Clipping: maximum number of "
02750                                 "iterations",
02751                                 "giraffe.biasremoval",
02752                                 5);
02753     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-niter");
02754     cpl_parameterlist_append(list, p);
02755 
02756 
02757     p = cpl_parameter_new_value("giraffe.biasremoval.fraction",
02758                                 CPL_TYPE_DOUBLE,
02759                                 "Sigma Clipping: minimum fraction of points "
02760                                 "accepted/total [0.0..1.0]",
02761                                 "giraffe.biasremoval",
02762                                 0.8);
02763     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-mfrac");
02764     cpl_parameterlist_append(list, p);
02765 
02766 
02767     p = cpl_parameter_new_value("giraffe.biasremoval.xorder",
02768                                 CPL_TYPE_INT,
02769                                 "Order of X polynomial fit (CURVE method "
02770                                 "only)",
02771                                 "giraffe.biasremoval",
02772                                 1);
02773     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-xorder");
02774     cpl_parameterlist_append(list, p);
02775 
02776     p = cpl_parameter_new_value("giraffe.biasremoval.yorder",
02777                                 CPL_TYPE_INT,
02778                                 "Order of Y polynomial fit (CURVE method "
02779                                 "only)",
02780                                 "giraffe.biasremoval",
02781                                 1);
02782     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-yorder");
02783     cpl_parameterlist_append(list, p);
02784 
02785 
02786     p = cpl_parameter_new_value("giraffe.biasremoval.xstep",
02787                                 CPL_TYPE_INT,
02788                                 "Sampling step along X (CURVE method only)",
02789                                 "giraffe.biasremoval",
02790                                 1);
02791     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-xstep");
02792     cpl_parameterlist_append(list, p);
02793 
02794 
02795     p = cpl_parameter_new_value("giraffe.biasremoval.ystep",
02796                                 CPL_TYPE_INT,
02797                                 "Sampling step along Y (CURVE method only)",
02798                                 "giraffe.biasremoval",
02799                                 1);
02800     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-ystep");
02801     cpl_parameterlist_append(list, p);
02802 
02803     return;
02804 }

This file is part of the GIRAFFE Pipeline Reference Manual 2.5.2.
Documentation copyright © 2002-2006 European Southern Observatory.
Generated on Fri Jun 13 14:36:21 2008 by doxygen 1.4.6 written by Dimitri van Heesch, © 1997-2004