GIRAFFE Pipeline Reference Manual

gibias.c

00001 /* $Id: gibias.c,v 1.30 2010/06/21 09:07:21 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: 2010/06/21 09:07:21 $
00024  * $Revision: 1.30 $
00025  * $Name: giraffe-2_8_6 $
00026  */
00027 
00028 #ifdef HAVE_CONFIG_H
00029 #  include <config.h>
00030 #endif
00031 
00032 #include <string.h>
00033 #include <math.h>
00034 #include <errno.h>
00035 
00036 #include <cxmemory.h>
00037 #include <cxstrutils.h>
00038 #include <cxutils.h>
00039 
00040 #include <cpl_msg.h>
00041 #include <cpl_parameter.h>
00042 
00043 #include "gimacros.h"
00044 #include "gialias.h"
00045 #include "giarray.h"
00046 #include "gimatrix.h"
00047 #include "gichebyshev.h"
00048 #include "gimath.h"
00049 #include "gibias.h"
00050 
00051 
00060 struct GiBiasResults {
00061 
00062     cxdouble mean;
00063     cxdouble sigma;
00064     cxdouble rms;
00065 
00066     cx_string* limits;
00067 
00068     cpl_matrix* coeffs;
00069 
00070     cpl_image* model;
00071 
00072 };
00073 
00074 typedef struct GiBiasResults GiBiasResults;
00075 
00076 
00077 /*
00078  * Clear a bias results structure
00079  */
00080 
00081 inline static void
00082 _giraffe_biasresults_clear(GiBiasResults *self)
00083 {
00084 
00085     if (self != NULL) {
00086 
00087         self->mean  = 0.;
00088         self->sigma = 0.;
00089         self->rms   = 0.;
00090 
00091         if (self->limits) {
00092             cx_string_delete(self->limits);
00093             self->limits = NULL;
00094         }
00095 
00096         if (self->coeffs) {
00097             cpl_matrix_delete(self->coeffs);
00098             self->coeffs = NULL;
00099         }
00100 
00101         if (self->model) {
00102             cpl_image_delete(self->model);
00103             self->model = NULL;
00104         }
00105 
00106     }
00107 
00108     return;
00109 
00110 }
00111 
00112 
00113 /*
00114  * Transfor bias method and option value into a string
00115  */
00116 
00117 inline static void
00118 _giraffe_method_string(cx_string *string, GiBiasMethod method,
00119                        GiBiasOption option)
00120 {
00121 
00122     switch (method) {
00123     case GIBIAS_METHOD_UNIFORM:
00124         cx_string_set(string, "UNIFORM");
00125         break;
00126 
00127     case GIBIAS_METHOD_PLANE:
00128         cx_string_set(string, "PLANE");
00129         break;
00130 
00131     case GIBIAS_METHOD_CURVE:
00132         cx_string_set(string, "CURVE");
00133         break;
00134 
00135     case GIBIAS_METHOD_PROFILE:
00136         cx_string_set(string, "PROFILE");
00137         break;
00138 
00139     case GIBIAS_METHOD_MASTER:
00140         cx_string_set(string, "MASTER");
00141         break;
00142 
00143     case GIBIAS_METHOD_ZMASTER:
00144         cx_string_set(string, "ZMASTER");
00145         break;
00146 
00147     default:
00148         break;
00149     }
00150 
00151     if (option != GIBIAS_OPTION_UNDEFINED) {
00152         switch (option) {
00153         case GIBIAS_OPTION_PLANE:
00154             cx_string_append(string, "+PLANE");
00155             break;
00156 
00157         case GIBIAS_OPTION_CURVE:
00158             cx_string_append(string, "+CURVE");
00159             break;
00160 
00161         default:
00162             break;
00163         }
00164     }
00165 
00166     return;
00167 
00168 }
00169 
00170 
00171 inline static void
00172 _giraffe_stringify_coefficients(cx_string *string, cpl_matrix *matrix)
00173 {
00174 
00175     register cxint i, j;
00176 
00177     cxchar *tmp = cx_line_alloc();
00178 
00179     cxint nr = cpl_matrix_get_nrow(matrix);
00180     cxint nc = cpl_matrix_get_ncol(matrix);
00181 
00182     cxsize sz = cx_line_max();
00183 
00184     cxdouble *data = cpl_matrix_get_data(matrix);
00185 
00186 
00187     for (i = 0; i < nr; i++) {
00188         for (j = 0; j < nc; j++) {
00189             snprintf(tmp, sz, "%g", data[i * nc  + j]);
00190             cx_string_append(string, tmp);
00191 
00192             if (i != nr - 1 || j < nc - 1) {
00193                 cx_string_append(string, ":");
00194             }
00195         }
00196     }
00197 
00198     cx_free(tmp);
00199 
00200     return;
00201 
00202 }
00203 
00204 
00220 inline static cxbool
00221 _giraffe_compare_overscans(const GiImage* image1, const GiImage* image2)
00222 {
00223 
00224     cxint32 l1ovscx = -1;
00225     cxint32 l1ovscy = -1;
00226     cxint32 l1prscx = -1;
00227     cxint32 l1prscy = -1;
00228     cxint32 l2ovscx = -1;
00229     cxint32 l2ovscy = -1;
00230     cxint32 l2prscx = -1;
00231     cxint32 l2prscy = -1;
00232 
00233     cpl_propertylist *l1, *l2;
00234 
00235 
00236     cx_assert(image1 != NULL && image2 != NULL);
00237 
00238     l1 = giraffe_image_get_properties(image1);
00239     l2 = giraffe_image_get_properties(image2);
00240 
00241     if (cpl_propertylist_has(l1, GIALIAS_OVSCX)) {
00242         l1ovscx = cpl_propertylist_get_int(l1, GIALIAS_OVSCX);
00243     }
00244     if (cpl_propertylist_has(l1, GIALIAS_OVSCY)) {
00245         l1ovscy = cpl_propertylist_get_int(l1, GIALIAS_OVSCY);
00246     }
00247     if (cpl_propertylist_has(l1, GIALIAS_PRSCX)) {
00248         l1prscx = cpl_propertylist_get_int(l1, GIALIAS_PRSCX);
00249     }
00250     if (cpl_propertylist_has(l1, GIALIAS_PRSCY)) {
00251         l1prscy = cpl_propertylist_get_int(l1, GIALIAS_PRSCY);
00252     }
00253 
00254     if (cpl_propertylist_has(l2, GIALIAS_OVSCX)) {
00255         l2ovscx = cpl_propertylist_get_int(l2, GIALIAS_OVSCX);
00256     }
00257     if (cpl_propertylist_has(l2, GIALIAS_OVSCY)) {
00258         l2ovscy = cpl_propertylist_get_int(l2, GIALIAS_OVSCY);
00259     }
00260     if (cpl_propertylist_has(l2, GIALIAS_PRSCX)) {
00261         l2prscx = cpl_propertylist_get_int(l2, GIALIAS_PRSCX);
00262     }
00263     if (cpl_propertylist_has(l2, GIALIAS_PRSCY)) {
00264         l2prscy = cpl_propertylist_get_int(l2, GIALIAS_PRSCY);
00265     }
00266 
00267     if (l1ovscx != l2ovscx || l1ovscy != l2ovscy) {
00268         return FALSE;
00269     }
00270 
00271     if (l1prscx != l2prscx || l1prscy != l2prscy) {
00272         return FALSE;
00273     }
00274 
00275     return TRUE;
00276 
00277 }
00278 
00279 
00280 /*
00281  * @brief
00282  *   Parse a bias area specification string.
00283  *
00284  * @param string  The string to parse.
00285  *
00286  * @return
00287  *   The function returns a newly allocated matrix containing the limits
00288  *   of the specified areas.
00289  *
00290  * TBD
00291  */
00292 
00293 inline static cpl_matrix*
00294 _giraffe_bias_get_areas(const cxchar* string)
00295 {
00296 
00297     cxchar** regions = NULL;
00298 
00299     cpl_matrix* areas = NULL;
00300 
00301 
00302     cx_assert(string != NULL);
00303 
00304     regions = cx_strsplit(string, ",", -1);
00305 
00306     if (regions == NULL) {
00307 
00308         return NULL;
00309 
00310     }
00311     else {
00312 
00313         const cxsize nvalues = 4;
00314 
00315         cxsize i = 0;
00316         cxsize nregions = 0;
00317 
00318         while (regions[i] != NULL) {
00319             ++i;
00320         }
00321 
00322         nregions = i;
00323         areas = cpl_matrix_new(nregions, nvalues);
00324 
00325         i = 0;
00326         while (regions[i] != NULL) {
00327 
00328             register cxsize j = 0;
00329 
00330             cxchar** limits = cx_strsplit(regions[i], ":", nvalues);
00331 
00332 
00333             if (limits == NULL) {
00334 
00335                 cpl_matrix_delete(areas);
00336                 areas = NULL;
00337 
00338                 return NULL;
00339 
00340             }
00341 
00342             for (j = 0; j < nvalues; ++j) {
00343 
00344                 cxchar* last = NULL;
00345 
00346                 cxint status = 0;
00347                 cxlong value = 0;
00348 
00349 
00350                 if (limits[j] == NULL || *limits[j] == '\0') {
00351                     break;
00352                 }
00353 
00354                 status = errno;
00355                 errno = 0;
00356 
00357                 value = strtol(limits[j], &last, 10);
00358 
00359                 /*
00360                  * Check for various errors
00361                  */
00362 
00363                 if ((errno == ERANGE &&
00364                     (value == LONG_MAX || value == LONG_MIN)) ||
00365                     (errno != 0 && value == 0)) {
00366 
00367                     break;
00368 
00369                 }
00370 
00371                 errno = status;
00372 
00373 
00374                 /*
00375                  * Check that the input string was parsed completely.
00376                  */
00377 
00378                 if (*last != '\0') {
00379                     break;
00380                 }
00381 
00382                 cpl_matrix_set(areas, i, j, value);
00383 
00384             }
00385 
00386             cx_strfreev(limits);
00387             limits = NULL;
00388 
00389             if (j != nvalues) {
00390 
00391                 cpl_matrix_delete(areas);
00392                 areas = NULL;
00393 
00394                 cx_strfreev(regions);
00395                 regions = NULL;
00396 
00397                 return NULL;
00398 
00399             }
00400 
00401             ++i;
00402 
00403         }
00404 
00405         cx_strfreev(regions);
00406         regions = NULL;
00407 
00408     }
00409 
00410     return areas;
00411 
00412 }
00413 
00414 
00415 /*
00416  * @brief
00417  *   Compute bias mean and sigma values over bias areas.
00418  *
00419  * @param  image        Image over which to calculate mean and sigma
00420  *                      inside bias areas
00421  * @param  areas        Bias areas specifications [N,4] as a matrix
00422  *                      where N is the number of areas
00423  * @param  kappa        Sigma clipping algorithm: kappa value
00424  * @param  numiter      Sigma clipping algorithm: max number of iterations
00425  * @param  maxfraction  Sigma clipping algorithm: max. fraction of pixels
00426  *
00427  * @param  results      A @em GiBiasResults structure which contains
00428  *                      the computed mean bias value, its standard deviation
00429  *                      and the coordinates specifications of the areas used
00430  *                      on return.
00431  *
00432  * @returns EXIT_SUCCESS on sucess, negative value otherwise
00433  *
00434  * Computes mean average of pixels value over areas defined by @em areas.
00435  * @em areas is a matrix [N, 4] specifying the lower left and upper right
00436  * corners of each of the @em N areas of the CCD used for bias determination.
00437  *
00438  * @code
00439  *   areas = [[xmin1, xmax1, ymin1, ymax1],
00440  *            [xmin2, xmax2, ymin2, ymax2],
00441  *             ....., ....., ....., .....
00442  *            [xminN, xmaxN, yminN, ymaxN]]
00443  * @endcode
00444  *
00445  * Areas pixels whose value is too far from mean average are rejected by
00446  * kappa-sigma clipping defined by parameters kappa, numiter and maxfraction.
00447  * This rejection loops until good/bad points ratio is enough or the
00448  * maximum number of iteration has been reached.
00449  *
00450  * The function returns the bias mean @em results->mean, the bias sigma
00451  * @em results->sigma computed over the kappa-sigma clipped areas pixels
00452  * value and @em results->limits string specifying the areas used as
00453  * follows:
00454  * @code
00455  *   "xmin1:xmax1:ymin1:ymax1;...;xminN:xmaxN:yminN:ymaxN"
00456  * @endcode
00457  *
00458  * @see giraffe_bias_compute_plane, giraffe_compute_remove_bias
00459  */
00460 
00461 inline static cxint
00462 _giraffe_bias_compute_mean(GiBiasResults* results, const cpl_image* image,
00463         const cpl_matrix* areas, cxdouble kappa, cxint numiter,
00464         cxdouble maxfraction)
00465 {
00466 
00467     const cxchar* const fctid = "giraffe_bias_compute_mean";
00468 
00469 
00470     const cxdouble* pdimg = NULL;
00471 
00472     cxint j, l, k;
00473     cxint img_dimx, img_dimy;
00474     cxint ba_num_cols;
00475     cxint ba_num_rows;
00476     cxint curriter = 0;
00477     cxint x0, x1, x2, x3;
00478 
00479     cxdouble sigma = 0.;
00480 
00481     cxlong ntotal    = 0L;
00482     cxlong naccepted = 0L;
00483     cxlong n         = 0L;
00484     cxlong pixcount  = 0L;
00485 
00486     cxdouble currfraction = 2.;
00487 
00488     cx_string* tmp = NULL;
00489 
00490     cpl_matrix* matrix_zz1;
00491     cpl_matrix* matrix_zz1diff;
00492 
00493 
00494     /*
00495      * Initialization
00496      */
00497 
00498     if (results->limits == NULL) {
00499         cpl_msg_info(fctid, "Unable to store biaslimits return parameter, "
00500                      "aborting...");
00501         return -3;
00502     }
00503 
00504     if (cpl_image_get_type(image) != CPL_TYPE_DOUBLE) {
00505         cpl_msg_info(fctid, "Only images allowed of type double, "
00506                      "aborting ...");
00507         return -3;
00508     }
00509 
00510 
00511     /*
00512      * Preprocessing
00513      */
00514 
00515     /*
00516      *  Validate Bias Areas and calculate number of points
00517      */
00518 
00519     if (areas == NULL) {
00520         cpl_msg_info(fctid, "Bias Areas: Missing bias areas, "
00521                      "aborting ...");
00522         return -1;
00523     }
00524 
00525     ba_num_cols = cpl_matrix_get_ncol(areas);
00526     ba_num_rows = cpl_matrix_get_nrow(areas);
00527 
00528     for (j = 0; j < ba_num_rows; j++) {
00529         x3 = (cxint) cpl_matrix_get(areas, j, 3);
00530         x2 = (cxint) cpl_matrix_get(areas, j, 2);
00531         x1 = (cxint) cpl_matrix_get(areas, j, 1);
00532         x0 = (cxint) cpl_matrix_get(areas, j, 0);
00533 
00534         ntotal += (cxulong) ((x3 - x2 + 1) * (x1 - x0 + 1));
00535     }
00536 
00537     if (ntotal <= 0) {
00538         cpl_msg_info(fctid, "Bias Areas: Inconsistent specification, "
00539                      "aborting ...");
00540         return -1;
00541     }
00542 
00543     matrix_zz1     = cpl_matrix_new(ntotal, 1);
00544     matrix_zz1diff = cpl_matrix_new(ntotal, 1);
00545 
00546 
00547     /*
00548      *  Compute Bias Areas Values
00549      */
00550 
00551     img_dimx = cpl_image_get_size_x(image);
00552     img_dimy = cpl_image_get_size_y(image);
00553 
00554     cx_string_set(results->limits, "");
00555 
00556     for (j = 0; j < ba_num_rows; j++) {
00557 
00558         x3 = (cxint) cpl_matrix_get(areas, j, 3);
00559         x2 = (cxint) cpl_matrix_get(areas, j, 2);
00560         x1 = (cxint) cpl_matrix_get(areas, j, 1);
00561         x0 = (cxint) cpl_matrix_get(areas, j, 0);
00562 
00563         if ((x0 > img_dimx) || (x1 > img_dimx) ||
00564             (x2 > img_dimy) || (x3 > img_dimy)) {
00565             continue;
00566         }
00567 
00568         tmp = cx_string_new();
00569 
00570         cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
00571         cx_string_append(results->limits, cx_string_get(tmp));
00572 
00573         cx_string_delete(tmp);
00574         tmp = NULL;
00575 
00576         pdimg = cpl_image_get_data_double_const(image);
00577 
00578         for (l = x2; l < x3 + 1; l++) {
00579             for (k = x0; k < x1 + 1; k++) {
00580                 cpl_matrix_set(matrix_zz1, n, 1, pdimg[k + l * img_dimx]);
00581                 n++;
00582             }
00583         }
00584     }
00585 
00586     if (n != ntotal) {
00587         cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting ...");
00588 
00589         cpl_matrix_delete(matrix_zz1);
00590         cpl_matrix_delete(matrix_zz1diff);
00591 
00592         return -3;
00593     }
00594 
00595     cpl_msg_info(fctid, "Bias Areas: Using %s",
00596                  cx_string_get(results->limits));
00597 
00598 
00599     /*
00600      * Processing
00601      */
00602 
00603     /*
00604      *  Perform Sigma Clipping...
00605      */
00606 
00607     cpl_msg_info(fctid, "Sigma Clipping : Start");
00608 
00609     naccepted = ntotal;
00610 
00611     while ((naccepted > 0) && (curriter < numiter) &&
00612            (currfraction > maxfraction)) {
00613 
00614         cxdouble ksigma = 0.;
00615 
00616         results->mean  = cpl_matrix_get_mean(matrix_zz1);
00617         sigma = giraffe_matrix_sigma_mean(matrix_zz1, results->mean);
00618 
00619         for (k = 0; k < cpl_matrix_get_nrow(matrix_zz1); k++) {
00620             cpl_matrix_set(matrix_zz1diff, k, 0,
00621                            cpl_matrix_get(matrix_zz1, k, 0) - results->mean);
00622         }
00623 
00624         cpl_msg_info(fctid, "bias[%d]: mean = %5g, sigma = %6.3g, "
00625                      "ratio = %6.3g, accepted = %ld\n", curriter,
00626                      results->mean, sigma, currfraction, naccepted);
00627 
00628         ksigma = sigma * kappa;
00629 
00630         pixcount = 0L;
00631         for (l = 0; l < cpl_matrix_get_nrow(matrix_zz1); l++) {
00632             if (fabs(cpl_matrix_get(matrix_zz1diff, l, 0)) > ksigma) {
00633                 continue;
00634             }
00635 
00636             cpl_matrix_set(matrix_zz1, pixcount, 0,
00637                            cpl_matrix_get(matrix_zz1, l, 0));
00638             ++pixcount;
00639         }
00640 
00641         cpl_matrix_resize(matrix_zz1, 0, 0, 0, pixcount -
00642                           cpl_matrix_get_nrow(matrix_zz1));
00643 
00644         cpl_matrix_resize(matrix_zz1diff, 0, 0, 0, pixcount -
00645                           cpl_matrix_get_nrow(matrix_zz1diff));
00646 
00647         if (pixcount == naccepted) {
00648             break;
00649         }
00650 
00651         naccepted = pixcount;
00652 
00653         currfraction = (cxdouble) naccepted / (cxdouble) ntotal;
00654         ++curriter;
00655     }
00656 
00657     cpl_msg_info(fctid, "Sigma Clipping : End");
00658 
00659 
00660     /*
00661      * Postprocessing
00662      */
00663 
00664     results->mean = cpl_matrix_get_mean(matrix_zz1);
00665     results->rms = giraffe_matrix_sigma_mean(matrix_zz1, results->mean);
00666 
00667     results->sigma = results->rms / sqrt(cpl_matrix_get_nrow(matrix_zz1));
00668 
00669 
00670     cpl_msg_info(fctid, "Sigma Clipping Results : bias[%d]: mean = %5g, "
00671                  "sigma = %6.3g, ratio = %6.3g, accepted=%ld\n", curriter,
00672                  results->mean, results->rms, currfraction, naccepted);
00673 
00674 
00675     /*
00676      * Cleanup
00677      */
00678 
00679     if (matrix_zz1 != NULL) {
00680         cpl_matrix_delete(matrix_zz1);
00681     }
00682 
00683     if (matrix_zz1diff != NULL) {
00684         cpl_matrix_delete(matrix_zz1diff);
00685     }
00686 
00687     return EXIT_SUCCESS;
00688 
00689 }
00690 
00691 
00692 /*
00693  * @brief
00694  *   Compute bias plane over bias areas.
00695  *
00696  * @param  image        Image over which to calculate mean and sigma
00697  *                      inside bias areas
00698  * @param  areas        Bias areas specifications [N,4] as a matrix
00699  *                      where N is the number of areas
00700  * @param  kappa        Sigma clipping algorithm: kappa value
00701  * @param  numiter      Sigma clipping algorithm: max number of iterations
00702  * @param  maxfraction  Sigma clipping algorithm: max. fraction of pixels
00703  *
00704  * @param  results      A @em GiBiasResults structure which contains
00705  *                      the computed coefficients of the fitted bias plane,
00706  *                      their standard deviation and the coordinates
00707  *                      specifications of the areas used on return.
00708  *
00709  * @returns EXIT_SUCCESS on sucess, negative value otherwise
00710  *
00711  * Fit a plane through the pixels value over areas defined by @em areas.
00712  * @em areas is a matrix [N, 4] specifying the lower left and upper right
00713  * corners of each of the @em N areas of the CCD used for bias determination.
00714  *
00715  * @code
00716  *   areas = [[xmin1, xmax1, ymin1, ymax1],
00717  *            [xmin2, xmax2, ymin2, ymax2],
00718  *             ....., ....., ....., .....
00719  *            [xminN, xmaxN, yminN, ymaxN]]
00720  * @endcode
00721  *
00722  * Areas pixels whose value is too far from mean average are rejected by
00723  * kappa-sigma clipping defined by parameters kappa, numiter and maxfraction.
00724  * This rejection loops until good/bad points ratio is enough or the
00725  * maximum number of iteration has been reached.
00726  *
00727  * The function returns the bias fitted plane coefficients @em BP1, BP2, BP3:
00728  * @code
00729  *   [BiasPlane] = [img] * inv([1XY]) over the bias areas where:
00730  *   [BiasPlane] = [BP1,BP2,BP3] coefficients of plane z = BP1+BP2*x+BP3*y
00731  *   [Z] = 1-column matrix of pixel values
00732  *   [1XY] = 3-columns matrix of all pixels coordinates [1,X,Y]
00733  * @endcode
00734  *
00735  * the bias sigma computed over the kappa sigma clipped areas pixels value
00736  * and a string specifying the areas used as follows
00737  * @code
00738  *   "xmin1:xmax1:ymin1:ymax1;...;xminN:xmaxN:yminN:ymaxN"
00739  * @endcode
00740  *
00741  * @see giraffe_bias_compute_mean, giraffe_compute_remove_bias
00742  */
00743 
00744 inline static cxint
00745 _giraffe_bias_compute_plane(GiBiasResults* results, const cpl_image* image,
00746                             const cpl_matrix* areas, cxdouble kappa,
00747                             cxint numiter, cxdouble maxfraction)
00748 {
00749 
00750     const cxchar* const fctid = "giraffe_bias_compute_plane";
00751 
00752 
00753     cxint j = 0;
00754     cxint nx = 0;
00755     cxint ny = 0;
00756     cxint nareas = 0;
00757     cxint iteration = 0;
00758 
00759     cxsize n = 0;
00760     cxsize ntotal = 0;
00761     cxsize naccepted = 0;
00762 
00763     cxdouble fraction = 1.;
00764     cxdouble sigma    = 0.;
00765 
00766     cpl_matrix* xbs = NULL;
00767     cpl_matrix* ybs = NULL;
00768     cpl_matrix* zbs = NULL;
00769     cpl_matrix* coeffs = NULL;
00770 
00771 
00772     cx_assert(results->limits != NULL);
00773     cx_assert(results->coeffs == NULL);
00774 
00775     cx_assert(areas != NULL);
00776 
00777     cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
00778 
00779 
00780     /*
00781      *  Validate Bias Areas and calculate number of points
00782      */
00783 
00784     nareas = cpl_matrix_get_nrow(areas);
00785 
00786     for (j = 0; j < nareas; j++) {
00787 
00788         cxint x3 = (cxint) cpl_matrix_get(areas, j, 3);
00789         cxint x2 = (cxint) cpl_matrix_get(areas, j, 2);
00790         cxint x1 = (cxint) cpl_matrix_get(areas, j, 1);
00791         cxint x0 = (cxint) cpl_matrix_get(areas, j, 0);
00792 
00793         ntotal += (cxsize) ((x3 - x2 + 1) * (x1 - x0 + 1));
00794 
00795     }
00796 
00797     if (ntotal <= 0) {
00798 
00799         cpl_msg_info(fctid, "Bias Areas: Inconsistent specification, "
00800                      "aborting ...");
00801         return -1;
00802 
00803     }
00804 
00805     nx = cpl_image_get_size_x(image);
00806     ny = cpl_image_get_size_y(image);
00807 
00808     cpl_msg_info(fctid, "Bias Areas: specified are %zu points in %dx%d "
00809                  "image", ntotal, nx, ny);
00810 
00811 
00812     /*
00813      *  Compute Bias Areas Values
00814      */
00815 
00816     results->mean  = 0.;
00817     results->sigma = 0.;
00818     results->rms   = 0.;
00819 
00820     cx_string_set(results->limits, "");
00821 
00822     xbs = cpl_matrix_new(ntotal, 1);
00823     ybs = cpl_matrix_new(ntotal, 1);
00824     zbs = cpl_matrix_new(1, ntotal);
00825 
00826     for (j = 0; j < nareas; ++j) {
00827 
00828         const cxdouble* _img =  cpl_image_get_data_double_const(image);
00829 
00830         cxint k = 0;
00831 
00832         cx_string* tmp = NULL;
00833 
00834 
00835         cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
00836         cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
00837         cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
00838         cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
00839 
00840         if ((x0 > nx) || (x1 > nx) || (x2 > ny) || (x3 > ny)) {
00841             continue;
00842         }
00843 
00844         tmp = cx_string_new();
00845 
00846         cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
00847         cx_string_append(results->limits, cx_string_get(tmp));
00848 
00849         cx_string_delete(tmp);
00850         tmp = NULL;
00851 
00852         for (k = x2; k < x3 + 1; ++k) {
00853 
00854             register cxint l = 0;
00855 
00856 
00857             for (l = x0; l < x1 + 1; ++l) {
00858 
00859                 cpl_matrix_set(xbs, n, 0, l);
00860                 cpl_matrix_set(ybs, n, 0, k);
00861                 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
00862                 ++n;
00863 
00864             }
00865 
00866         }
00867 
00868     }
00869 
00870     cpl_matrix_set_size(xbs, n, 1);
00871     cpl_matrix_set_size(ybs, n, 1);
00872     cpl_matrix_set_size(zbs, 1, n);
00873 
00874     if (n != ntotal) {
00875 
00876         cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting...");
00877 
00878         cpl_matrix_delete(xbs);
00879         cpl_matrix_delete(ybs);
00880         cpl_matrix_delete(zbs);
00881 
00882         return -1;
00883 
00884     }
00885 
00886     ntotal = n;
00887 
00888     cpl_msg_info(fctid, "Bias Areas: Using %s [%zu pixels]",
00889                  cx_string_get(results->limits), ntotal);
00890 
00891 
00892     /*
00893      *  Perform Sigma Clipping
00894      */
00895 
00896     cpl_msg_info(fctid, "Sigma Clipping : Start");
00897 
00898     naccepted = ntotal;
00899 
00900     while ((naccepted > 0) && (iteration < numiter) &&
00901            (fraction > maxfraction)) {
00902 
00903         cxsize k = 0;
00904 
00905         cxdouble ksigma = 0.;
00906 
00907         cpl_matrix* base = NULL;
00908         cpl_matrix* fit = NULL;
00909 
00910 
00911         base = cpl_matrix_new(3, naccepted);
00912 
00913         if (base == NULL) {
00914 
00915             cpl_msg_info(fctid, "Sigma Clipping: Error creating design "
00916                          "matrix");
00917 
00918             cpl_matrix_delete(zbs);
00919             cpl_matrix_delete(ybs);
00920             cpl_matrix_delete(xbs);
00921 
00922             return -2;
00923         }
00924 
00925         for (k = 0; k < naccepted; ++k) {
00926 
00927             cpl_matrix_set(base, 0, k, 1.);
00928             cpl_matrix_set(base, 1, k, cpl_matrix_get(xbs, k, 0));
00929             cpl_matrix_set(base, 2, k, cpl_matrix_get(ybs, k, 0));
00930 
00931         }
00932 
00933         cpl_matrix_delete(coeffs);
00934         coeffs = NULL;
00935 
00936         coeffs = giraffe_matrix_leastsq(base, zbs);
00937 
00938         if (coeffs == NULL) {
00939 
00940             cpl_msg_info(fctid, "Sigma Clipping : Error in least square "
00941                          "solution, aborting...");
00942 
00943             cpl_matrix_delete(base);
00944             base = NULL;
00945 
00946             cpl_matrix_delete(xbs);
00947             cpl_matrix_delete(ybs);
00948             cpl_matrix_delete(zbs);
00949 
00950             return -2;
00951 
00952         }
00953 
00954 
00955         /*
00956          * Compute bias model fit and reject deviating pixels
00957          */
00958 
00959         fit = cpl_matrix_product_create(coeffs, base);
00960 
00961         cpl_matrix_delete(base);
00962         base = NULL;
00963 
00964         results->mean = cpl_matrix_get_mean(fit);
00965 
00966         sigma = giraffe_matrix_sigma_fit(zbs, fit);
00967 
00968         cpl_msg_info(fctid, "Sigma Clipping : bias plane[%d]: %g + "
00969                      "%g * x + %g * y, sigma = %.5g, ratio = %.4g, "
00970                      "accepted = %zu\n", iteration,
00971                      cpl_matrix_get(coeffs, 0, 0),
00972                      cpl_matrix_get(coeffs, 0, 1),
00973                      cpl_matrix_get(coeffs, 0, 2),
00974                      sigma, fraction, naccepted);
00975 
00976 
00977         /*
00978          *  Perform Clipping
00979          */
00980 
00981         ksigma = sigma * kappa;
00982 
00983         n = 0;
00984 
00985         for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
00986 
00987             register cxdouble z = cpl_matrix_get(zbs, 0, j);
00988 
00989             if (fabs(cpl_matrix_get(fit, 0, j) - z) > ksigma) {
00990                 continue;
00991             }
00992 
00993             cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
00994 
00995             cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
00996 
00997             cpl_matrix_set(zbs, 0, n, z);
00998             ++n;
00999 
01000         }
01001 
01002         cpl_matrix_set_size(xbs, n, 1);
01003         cpl_matrix_set_size(ybs, n, 1);
01004         cpl_matrix_set_size(zbs, 1, n);
01005 
01006         cpl_matrix_delete(fit);
01007         fit = NULL;
01008 
01009         if (n == naccepted) {
01010             break;
01011         }
01012 
01013         naccepted = n;
01014 
01015         fraction = (cxdouble)naccepted / (cxdouble)ntotal;
01016         ++iteration;
01017 
01018     }
01019 
01020     cpl_msg_info(fctid, "Sigma Clipping : End");
01021 
01022 
01023     /*
01024      * Store fit results
01025      */
01026 
01027     results->coeffs = coeffs;
01028     results->rms    = sigma;
01029 
01030     // FIXME: The following is not correct. It should be the error of
01031     //        results->mean which has to be obtained from error propagation.
01032     //        At least the following is connected to the fitted model,
01033     //        instead of the original code which computed only the standard
01034     //        deviation of the raw data.
01035 
01036     results->sigma  = sigma / sqrt(cpl_matrix_get_ncol(zbs));
01037 
01038     cpl_msg_info(fctid, "Sigma Clipping Results (%d/%zu, sigma = %g)",
01039                  iteration, naccepted, results->rms);
01040 
01041 
01042     /*
01043      * Cleanup
01044      */
01045 
01046     cpl_matrix_delete(xbs);
01047     xbs = NULL;
01048 
01049     cpl_matrix_delete(ybs);
01050     ybs = NULL;
01051 
01052     cpl_matrix_delete(zbs);
01053     zbs = NULL;
01054 
01055     return EXIT_SUCCESS;
01056 
01057 }
01058 
01059 
01060 /*
01061  * @brief
01062  *   Compute bias curve over bias areas.
01063  *
01064  * @param results      computed bias curve over bias areas coordinates
01065  *                     specifications
01066  * @param image        Image over which to calculate bias curve
01067  *                     inside bias areas
01068  * @param areas        bias areas specifications [N,4] as a matrix
01069  *                     where N is the number of areas
01070  * @param kappa        sigma clipping algorithm: kappa value
01071  * @param numiter      sigma clipping algorithm: max number of iterations
01072  * @param maxfraction  sigma clipping algorithm: max. fraction of pixels
01073  * @param xdeg         Polynom order for fit in x direction
01074  * @param ydeg         Polynom order for fit in y direction
01075  * @param xstep        Step in x direction for polynomial fit
01076  * @param ystep        Step in y direction for polynomial fit
01077  *
01078  * @return
01079  *   EXIT_SUCCESS on sucess, negative value otherwise
01080  *
01081  * Fit a surface using a rectangular mesh defined by @em xstep and
01082  * @em ystep in @em areas.
01083  * Mesh pixels whose value is too far from fitted surface are rejected by
01084  * kappa-sigma clipping defined by parameters of @em kappa.
01085  * This rejection loops until good/bad points ratio is enough
01086  * (maxfraction) or the maximum number of iteration (numiter) is reached.
01087  *
01088  * The function returns the bias fitted surface coefficients
01089  * @a results->biascoeffs computed using a 2D Chebyshev polynomial fit
01090  * whose X and Y order are given respectively by @a xdeg and  @a ydeg,
01091  * The bias sigma @a results->biassigma computed over the
01092  * kappa-sigma clipped pixels value and a string specifying the areas
01093  * used as follows:
01094  * @code
01095  *   "xstart:xend:xstep;ymin,ystart,ystep"
01096  * @endcode
01097  *
01098  * @see giraffe_bias_compute_plane, giraffe_bias_compute_remove_bias
01099  */
01100 
01101 inline static cxint
01102 _giraffe_bias_compute_curve(GiBiasResults* results, const cpl_image* image,
01103                             const cpl_matrix *areas, cxdouble kappa,
01104                             cxint numiter, cxdouble maxfraction,
01105                             cxdouble xdeg, cxdouble ydeg,
01106                             cxdouble xstep, cxdouble ystep)
01107 {
01108 
01109     const cxchar* const fctid = "giraffe_bias_compute_curve";
01110 
01111     cxint j = 0;
01112     cxint nx = 0;
01113     cxint ny = 0;
01114     cxint nareas = 0;
01115     cxint iteration = 0;
01116 
01117     cxsize n = 0;
01118     cxsize ntotal = 0;
01119     cxsize naccepted = 0;
01120 
01121     cxdouble fraction = 1.;
01122     cxdouble sigma    = 0.;
01123 
01124     cpl_matrix* xbs = NULL;
01125     cpl_matrix* ybs = NULL;
01126     cpl_matrix* zbs = NULL;
01127 
01128     GiChebyshev2D* fit = NULL;
01129 
01130 
01131     cx_assert(results != NULL);
01132     cx_assert(results->limits != NULL);
01133     cx_assert(results->coeffs == NULL);
01134 
01135     cx_assert(areas != NULL);
01136 
01137     cx_assert(image != NULL);
01138     cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
01139 
01140 
01141     /*
01142      *  Validate Bias Areas and calculate number of points
01143      */
01144 
01145     nareas = cpl_matrix_get_nrow(areas);
01146 
01147     for (j = 0; j < nareas; ++j) {
01148 
01149         cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
01150         cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
01151         cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
01152         cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
01153 
01154         ntotal += (cxulong) (ceil((1. + x3 - x2) / ystep) *
01155                              ceil((1. + x1 - x0) / xstep));
01156     }
01157 
01158     nx = cpl_image_get_size_x(image);
01159     ny = cpl_image_get_size_y(image);
01160 
01161     cpl_msg_info(fctid, "Bias Areas: Found %zu points in %dx%d image",
01162                  ntotal, nx, ny);
01163 
01164 
01165     /*
01166      *  Compute Bias Areas Values
01167      */
01168 
01169     results->mean  = 0.;
01170     results->sigma = 0.;
01171     results->rms   = 0.;
01172 
01173     cx_string_set(results->limits, "");
01174 
01175     xbs = cpl_matrix_new(ntotal, 1);
01176     ybs = cpl_matrix_new(ntotal, 1);
01177     zbs = cpl_matrix_new(1, ntotal);
01178 
01179     for (j = 0; j < nareas; ++j) {
01180 
01181         const cxdouble* _img =  cpl_image_get_data_double_const(image);
01182 
01183         cxint k = 0;
01184 
01185         cx_string* tmp = NULL;
01186 
01187 
01188         cxint x3 = (cxint)cpl_matrix_get(areas, j, 3);
01189         cxint x2 = (cxint)cpl_matrix_get(areas, j, 2);
01190         cxint x1 = (cxint)cpl_matrix_get(areas, j, 1);
01191         cxint x0 = (cxint)cpl_matrix_get(areas, j, 0);
01192 
01193         if ((x0 > nx) || (x1 > nx) ||
01194             (x2 > ny) || (x3 > ny)) {
01195             continue;
01196         }
01197 
01198         tmp = cx_string_new();
01199 
01200         cx_string_sprintf(tmp, "%d:%d:%d:%d;", x0, x1, x2, x3);
01201         cx_string_append(results->limits, cx_string_get(tmp));
01202 
01203         cx_string_delete(tmp);
01204         tmp = NULL;
01205 
01206         for (k = x2; k < x3 + 1; k += ystep) {
01207 
01208             register cxint l = 0;
01209 
01210 
01211             for (l = x0; l < x1 + 1; l += xstep) {
01212 
01213                 cpl_matrix_set(xbs, n, 0, l);
01214                 cpl_matrix_set(ybs, n, 0, k);
01215                 cpl_matrix_set(zbs, 0, n, _img[k * nx + l]);
01216                 ++n;
01217 
01218             }
01219 
01220         }
01221 
01222     }
01223 
01224     cpl_matrix_set_size(xbs, n, 1);
01225     cpl_matrix_set_size(ybs, n, 1);
01226     cpl_matrix_set_size(zbs, 1, n);
01227 
01228     if (n <= 0) {
01229 
01230         cpl_msg_info(fctid, "Bias Areas: Validation failed, aborting...");
01231 
01232         cpl_matrix_delete(xbs);
01233         cpl_matrix_delete(ybs);
01234         cpl_matrix_delete(zbs);
01235 
01236         return -1;
01237 
01238     }
01239 
01240     ntotal = n;
01241 
01242     cpl_msg_info(fctid, "Bias Areas: Using %s [%zu pixels]",
01243                  cx_string_get(results->limits), ntotal);
01244 
01245 
01246     /*
01247      *  Perform Sigma Clipping
01248      */
01249 
01250     cpl_msg_info(fctid, "Sigma Clipping : Start");
01251 
01252     naccepted = ntotal;
01253 
01254     while ((naccepted > 0) && (iteration < numiter) &&
01255            (fraction > maxfraction)) {
01256 
01257         cxint status = 0;
01258 
01259         cxdouble ksigma = 0.;
01260 
01261         cpl_matrix* base = NULL;
01262         cpl_matrix* coeffs = NULL;
01263         cpl_matrix* _coeffs = NULL;
01264         cpl_matrix* _fit = NULL;
01265 
01266 
01267         base = giraffe_chebyshev_base2d(0., 0., nx, ny,
01268                                         xdeg, ydeg, xbs, ybs);
01269 
01270         if (base == NULL) {
01271 
01272             cpl_msg_info(fctid, "Sigma Clipping: Error creating design "
01273                          "matrix");
01274 
01275             cpl_matrix_delete(zbs);
01276             cpl_matrix_delete(ybs);
01277             cpl_matrix_delete(xbs);
01278 
01279             return -2;
01280         }
01281 
01282         _coeffs = giraffe_matrix_leastsq(base, zbs);
01283 
01284         cpl_matrix_delete(base);
01285         base = NULL;
01286 
01287         if (_coeffs == NULL) {
01288 
01289             cpl_msg_info(fctid, "Sigma Clipping : Error in least square "
01290                          "solution, aborting...");
01291 
01292             cpl_matrix_delete(xbs);
01293             cpl_matrix_delete(ybs);
01294             cpl_matrix_delete(zbs);
01295 
01296             return -2;
01297 
01298         }
01299 
01300 
01301         /*
01302          * Compute bias model fit and reject deviating pixels
01303          */
01304 
01305         coeffs = cpl_matrix_wrap(xdeg, ydeg, cpl_matrix_get_data(_coeffs));
01306 
01307         if (fit != NULL) {
01308             giraffe_chebyshev2d_delete(fit);
01309             fit = NULL;
01310         }
01311 
01312         fit = giraffe_chebyshev2d_new(xdeg - 1, ydeg - 1);
01313         status = giraffe_chebyshev2d_set(fit, 0., nx, 0., ny,
01314                                          coeffs);
01315 
01316         if (status != 0) {
01317 
01318             giraffe_chebyshev2d_delete(fit);
01319             fit = NULL;
01320 
01321             cpl_matrix_unwrap(coeffs);
01322             coeffs = NULL;
01323 
01324             cpl_matrix_delete(_coeffs);
01325             _coeffs = NULL;
01326 
01327             cpl_matrix_delete(xbs);
01328             cpl_matrix_delete(ybs);
01329             cpl_matrix_delete(zbs);
01330 
01331             return -2;
01332 
01333         }
01334 
01335         cpl_matrix_unwrap(coeffs);
01336         coeffs = NULL;
01337 
01338         cpl_matrix_delete(_coeffs);
01339         _coeffs = NULL;
01340 
01341         _fit = cpl_matrix_new(1, cpl_matrix_get_ncol(zbs));
01342 
01343         for (j = 0; j < cpl_matrix_get_ncol(_fit); ++j) {
01344 
01345             cxdouble x = cpl_matrix_get(xbs, n, 0);
01346             cxdouble y = cpl_matrix_get(ybs, n, 0);
01347             cxdouble z = giraffe_chebyshev2d_eval(fit, x, y);
01348 
01349             cpl_matrix_set(_fit, 0, j, z);
01350 
01351         }
01352 
01353         results->mean = cpl_matrix_get_mean(_fit);
01354 
01355         sigma = giraffe_matrix_sigma_fit(zbs, _fit);
01356 
01357         cpl_msg_info(fctid, "Sigma Clipping : bias surface[%d]: "
01358                      "sigma = %8.5g, ratio = %7.4g, accepted = %zu\n",
01359                      iteration, sigma, fraction, naccepted);
01360 
01361 
01362         /*
01363          *  Perform Clipping
01364          */
01365 
01366         ksigma = sigma * kappa;
01367 
01368         n = 0;
01369 
01370         for (j = 0; j < cpl_matrix_get_ncol(zbs); ++j) {
01371 
01372             register cxdouble z = cpl_matrix_get(zbs, 0, j);
01373 
01374             if (fabs(cpl_matrix_get(_fit, 0, j) - z) > ksigma) {
01375                 continue;
01376             }
01377 
01378             cpl_matrix_set(xbs, n, 0, cpl_matrix_get(xbs, j, 0));
01379             cpl_matrix_set(ybs, n, 0, cpl_matrix_get(ybs, j, 0));
01380             cpl_matrix_set(zbs, 0, n, z);
01381             ++n;
01382 
01383         }
01384 
01385         cpl_matrix_set_size(xbs, n, 1);
01386         cpl_matrix_set_size(ybs, n, 1);
01387         cpl_matrix_set_size(zbs, 1, n);
01388 
01389         cpl_matrix_delete(_fit);
01390         _fit = NULL;
01391 
01392 
01393         if (n == naccepted) {
01394             break;
01395         }
01396 
01397         naccepted = n;
01398 
01399         fraction = (cxdouble)naccepted / (cxdouble)ntotal;
01400         ++iteration;
01401 
01402     }
01403 
01404     cpl_msg_info(fctid, "Sigma Clipping : End");
01405 
01406 
01407     /*
01408      * Store fit results
01409      */
01410 
01411     results->coeffs = cpl_matrix_duplicate(giraffe_chebyshev2d_coeffs(fit));
01412     results->rms    = sigma;
01413 
01414     // FIXME: The following is not correct. It should be the error of
01415     //        results->mean which has to be obtained from error propagation.
01416     //        At least the following is connected to the fitted model,
01417     //        instead of the original code which computed only the standard
01418     //        deviation of the raw data.
01419 
01420     results->sigma  = sigma / sqrt(cpl_matrix_get_ncol(zbs));
01421 
01422     cpl_msg_info(fctid, "Sigma Clipping Results (%d/%zu, sigma = %g)",
01423                  iteration, naccepted, results->rms);
01424 
01425 
01426     /*
01427      * Cleanup
01428      */
01429 
01430     giraffe_chebyshev2d_delete(fit);
01431     fit = NULL;
01432 
01433     cpl_matrix_delete(xbs);
01434     xbs = NULL;
01435 
01436     cpl_matrix_delete(ybs);
01437     ybs = NULL;
01438 
01439     cpl_matrix_delete(zbs);
01440     zbs = NULL;
01441 
01442     return EXIT_SUCCESS;
01443 
01444 }
01445 
01446 
01447 /*
01448  * @brief
01449  *   Compute a bias profile model from the collapsed bias areas.
01450  *
01451  * @param results  Container for the computed mean bias, rms and the profile
01452  *                 model.
01453  * @param image    The image from which the profile model is computed.
01454  * @param areas    Image regions used to compute the bias model.
01455  *
01456  * @return
01457  *   The function returns @c EXIT_SUCCESS on success, and a non-zero
01458  *   value otherwise.
01459  *
01460  * The function computes a bias model profile along the axis @em axis from
01461  * the image areas @em areas of the input image @em image. The profile model
01462  * is computed by computing the mean pixel value along the axis
01463  * perpendicular to @em axis for each image area listed in @em areas.
01464  * The bias model value is then the average mean value of all given bias
01465  * areas for each pixel along @em axis.
01466  *
01467  * The given bias areas must have the same size along the axis @em axis!
01468  */
01469 
01470 inline static cxint
01471 _giraffe_bias_compute_profile(GiBiasResults* results, const cpl_image* image,
01472                               const cpl_matrix* areas, cxchar axis)
01473 {
01474 
01475     const cxchar* const fctid = "_giraffe_bias_compute_profile";
01476 
01477 
01478     cxint nx = 0;
01479     cxint ny = 0;
01480     cxint sx = 0;
01481     cxint sy = 0;
01482 
01483     cxsize j = 0;
01484     cxsize npixel = 0;
01485     cxsize ntotal = 0;
01486     cxsize nareas = 0;
01487     cxsize nvalid = 0;
01488 
01489     cxdouble rms   = 0.;
01490     cxdouble sigma = 0.;
01491 
01492     cpl_matrix* _areas = NULL;
01493 
01494     cpl_image* profile = NULL;
01495     cpl_image* model = NULL;
01496 
01497 
01498     cx_assert(results != NULL);
01499     cx_assert(results->limits != NULL);
01500     cx_assert(results->model == NULL);
01501 
01502     cx_assert(areas != NULL);
01503 
01504     cx_assert(image != NULL);
01505     cx_assert(cpl_image_get_type(image) == CPL_TYPE_DOUBLE);
01506 
01507     cx_assert((axis == 'x') || (axis == 'y'));
01508 
01509 
01510     nx = cpl_image_get_size_x(image);
01511     ny = cpl_image_get_size_y(image);
01512 
01513 
01514     /*
01515      *  Validate Bias Areas and calculate number of points
01516      */
01517 
01518     _areas = cpl_matrix_duplicate(areas);
01519     cx_assert(_areas != NULL);
01520 
01521     nareas = cpl_matrix_get_nrow(_areas);
01522 
01523     if (axis == 'y') {
01524 
01525         for (j = 0; j < nareas; ++j) {
01526 
01527             cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
01528             cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
01529             cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
01530             cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
01531 
01532 
01533             if (x_0 < 0) {
01534                 x_0 = 0;
01535                 cpl_matrix_set(_areas, j, 0, x_0);
01536             }
01537 
01538             if (x_1 >= nx) {
01539                 x_1 = nx - 1;
01540                 cpl_matrix_set(_areas, j, 1, x_1);
01541             }
01542 
01543             if (y_0 != 0) {
01544                 y_0 = 0;
01545                 cpl_matrix_set(_areas, j, 2, y_0);
01546             }
01547 
01548             if (y_1 != ny - 1) {
01549                 y_1 = ny - 1;
01550                 cpl_matrix_set(_areas, j, 3, y_1);
01551             }
01552 
01553             if ((x_1 >= x_0) && (y_1 >= y_0)) {
01554                 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
01555             }
01556 
01557         }
01558 
01559         sx = nareas;
01560         sy = ny;
01561     }
01562     else {
01563 
01564         for (j = 0; j < nareas; ++j) {
01565 
01566             cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
01567             cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
01568             cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
01569             cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
01570 
01571 
01572             if (x_0 != 0) {
01573                 x_0 = 0;
01574                 cpl_matrix_set(_areas, j, 0, x_0);
01575             }
01576 
01577             if (x_1 != nx) {
01578                 x_1 = nx - 1;
01579                 cpl_matrix_set(_areas, j, 1, x_1);
01580             }
01581 
01582             if (y_0 < 0) {
01583                 y_0 = 0;
01584                 cpl_matrix_set(_areas, j, 2, y_0);
01585             }
01586 
01587             if (y_1 >= ny) {
01588                 y_1 = ny - 1;
01589                 cpl_matrix_set(_areas, j, 3, y_1);
01590             }
01591 
01592             if ((x_1 >= x_0) && (y_1 >= y_0)) {
01593                 ntotal += (y_1 - y_0 + 1) * (x_1 - x_0 + 1);
01594             }
01595 
01596         }
01597 
01598         sx = nareas;
01599         sy = nx;
01600 
01601     }
01602 
01603     cpl_msg_info(fctid, "Bias Areas: Found %zu points in %dx%d image",
01604                  ntotal, nx, ny);
01605 
01606 
01607     /*
01608      *  Compute Bias Areas Values
01609      */
01610 
01611     results->mean  = 0.;
01612     results->sigma = 0.;
01613     results->rms   = 0.;
01614 
01615     cx_string_set(results->limits, "");
01616 
01617     profile = cpl_image_new(sx, sy, CPL_TYPE_DOUBLE);
01618 
01619     for (j = 0; j < nareas; ++j) {
01620 
01621         cxint y_1 = (cxint)cpl_matrix_get(_areas, j, 3);
01622         cxint y_0 = (cxint)cpl_matrix_get(_areas, j, 2);
01623         cxint x_1 = (cxint)cpl_matrix_get(_areas, j, 1);
01624         cxint x_0 = (cxint)cpl_matrix_get(_areas, j, 0);
01625 
01626         if ((x_1 >= x_0) && (y_1 >= y_0)) {
01627 
01628             cx_string* tmp = cx_string_new();
01629 
01630 
01631             cx_string_sprintf(tmp, "%d:%d:%d:%d;", x_0, x_1, y_0, y_1);
01632             cx_string_append(results->limits, cx_string_get(tmp));
01633 
01634             cx_string_delete(tmp);
01635             tmp = NULL;
01636 
01637 
01638             /*
01639              * Update the number of actually used (valid) areas.
01640              */
01641 
01642             ++nvalid;
01643 
01644 
01645             /*
01646              * Compute the profile along the given axis for each bias region.
01647              */
01648 
01649             if (axis == 'y') {
01650 
01651                 cxint i = 0;
01652                 cxint sz = x_1 - x_0 + 1;
01653 
01654                 cxdouble residual = 0.;
01655 
01656                 const cxdouble* _img = cpl_image_get_data_double_const(image);
01657 
01658                 cxdouble* _profile = cpl_image_get_data_double(profile);
01659 
01660 
01661                 for (i = y_0; i < y_1 + 1; ++i) {
01662 
01663                     register cxint k = i * nx;
01664                     register cxint l = 0;
01665 
01666                     cxdouble m = giraffe_array_mean(&_img[k + x_0], sz);
01667 
01668                     _profile[i * sx + j] = m;
01669 
01670                     for (l = x_0; l < x_1 + 1; ++l) {
01671                         residual += (_img[k + l] - m) * (_img[k + l] - m);
01672                     }
01673                     ++npixel;
01674 
01675                 }
01676 
01677 
01678                 /*
01679                  * Compute the mean RMS and the mean bias error
01680                  */
01681 
01682                 rms   += residual / (sz - 1);
01683                 sigma += residual / (sz * (sz - 1));
01684 
01685             }
01686             else {
01687 
01688                 cxint i = 0;
01689                 cxint sz = y_1 - y_0 + 1;
01690 
01691                 cxdouble residual = 0.;
01692 
01693                 const cxdouble* _img = cpl_image_get_data_double_const(image);
01694 
01695                 cxdouble* _profile = cpl_image_get_data_double(profile);
01696                 cxdouble* buffer = cx_calloc(sz, sizeof(cxdouble));
01697 
01698 
01699                 for (i = x_0; i < x_1 + 1; ++i) {
01700 
01701                     register cxint l = 0;
01702 
01703                     register cxdouble m = 0.;
01704 
01705 
01706                     for (l = y_0; l < y_1 + 1; ++l) {
01707                         buffer[l - y_0] = _img[l * nx + i];
01708                     }
01709 
01710                     m = giraffe_array_mean(buffer, sz);
01711                     _profile[i * sx + j] = m;
01712 
01713                     for (l = 0; l < sz; ++l) {
01714                         residual += (buffer[l] - m) * (buffer[l] - m);
01715                     }
01716                     ++npixel;
01717 
01718                 }
01719 
01720 
01721                 /*
01722                  * Compute the mean RMS and the mean bias error
01723                  */
01724 
01725                 rms   += residual / (sz - 1);
01726                 sigma += residual / (sz * (sz - 1));
01727 
01728                 cx_free(buffer);
01729                 buffer = NULL;
01730 
01731             }
01732 
01733         }
01734 
01735     }
01736 
01737     cpl_matrix_delete(_areas);
01738     _areas = NULL;
01739 
01740 
01741     /*
01742      * Compute average profile of all bias areas
01743      */
01744 
01745     model = cpl_image_collapse_create(profile, 1);
01746 
01747     cpl_image_delete(profile);
01748     profile = NULL;
01749 
01750     cpl_image_divide_scalar(model, nvalid);
01751     results->model = model;
01752 
01753     model = NULL;
01754 
01755     results->mean = cpl_image_get_mean(results->model);
01756 
01757     if (nvalid > 0) {
01758         results->rms   = sqrt(rms / (sy * nvalid));
01759         results->sigma = sqrt(sigma / sy) / nvalid;
01760     }
01761 
01762     return EXIT_SUCCESS;
01763 
01764 }
01765 
01766 
01767 /*
01768  * @brief
01769  *   Actual Calculation/Logic Module of Bias Removal
01770  *
01771  * @param img            Input image
01772  * @param biasareas      Bias area to use
01773  * @param master_bias    Master bias Image
01774  * @param bad_pixel_map  Bad pixel map Image
01775  * @param method         Method to use for bias removal
01776  * @param option         Option to use for bias removal
01777  * @param model          Model to use for bias removal
01778  * @param remove_bias    Remove bias from result image
01779  * @param mbias          Master Bias value read from Master Bias
01780  *                       FITS Keywords
01781  * @param xdeg           X Degree of Polynomial of Chebyshev fit
01782  * @param ydeg           Y Degree of Polynomial of Chebyshev fit
01783  * @param xstep          X Step to use during Chebyshev fit
01784  * @param ystep          Y Step to use during Chebyshev fit
01785  * @param sigma          Sigma to use during Kappa Sigma Clipping
01786  * @param numiter        Max Num Iterations to use during Kappa Sigma
01787  *                       Clipping
01788  * @param maxfraction    Maximum Fraction to use during Kappa Sigma Clipping
01789  * @param results        Values determined
01790  *
01791  * @return The function returns 0 on success, or an error code otherwise.
01792  *
01793  * TBD
01794  */
01795 
01796 inline static cxint
01797 _giraffe_bias_compute(GiBiasResults* results, cpl_image* img,
01798                       const cpl_matrix* biasareas,
01799                       const cpl_image* master_bias,
01800                       const cpl_image* bad_pixel_map,
01801                       GiBiasMethod method, GiBiasOption option,
01802                       GiBiasModel model, cxbool remove_bias, cxdouble mbias,
01803                       cxdouble xdeg, cxdouble ydeg, cxdouble xstep,
01804                       cxdouble ystep, cxdouble sigma, cxint numiter,
01805                       cxdouble maxfraction)
01806 {
01807 
01808     const cxchar* const fctid = "giraffe_bias_compute";
01809 
01810     cxint i;
01811     cxint j;
01812     cxint img_dimx;
01813     cxint img_dimy;
01814     cxint img_centerx;
01815     cxint img_centery;
01816     cxint master_bias_dimx;
01817     cxint master_bias_dimy;
01818     cxint bad_pixel_dimx;
01819     cxint bad_pixel_dimy;
01820 
01821     cxint error_code = 0;
01822 
01823     const cpl_matrix* coeff = NULL;              //TODO do we need this????
01824 
01825     cpl_matrix* tmp_x;
01826     cpl_matrix* tmp_y;
01827     cpl_matrix* tmp_fit;
01828 
01829     cx_assert(results != NULL);
01830     cx_assert(results->limits == NULL);
01831     cx_assert(results->coeffs == NULL);
01832     cx_assert(results->model == NULL);
01833 
01834     cx_assert(img != NULL);
01835     cx_assert(cpl_image_get_type(img) == CPL_TYPE_DOUBLE);
01836 
01837     cx_assert(biasareas != NULL);
01838 
01839 
01840     /*
01841      * Preprocessing
01842      */
01843 
01844     img_dimx = cpl_image_get_size_x(img);
01845     img_dimy = cpl_image_get_size_y(img);
01846 
01847     img_centerx = img_dimx / 2;
01848     img_centery = img_dimy / 2;
01849 
01850     results->limits = cx_string_new();
01851 
01852 
01853     /*
01854      * Processing
01855      */
01856 
01857     if (remove_bias) {
01858         cpl_msg_info(fctid, "removing bias value");
01859     }
01860     else {
01861         cpl_msg_info(fctid, "NOT removing bias value");
01862     }
01863 
01864 
01865     if (model == GIBIAS_MODEL_MEAN) {
01866 
01867         /*
01868          *  Compute bias average value over biaslimits areas
01869          */
01870 
01871         cpl_msg_info(fctid, "using MEAN Model...");
01872 
01873         if ((option == GIBIAS_OPTION_PLANE) ||
01874             (method == GIBIAS_METHOD_PLANE)) {
01875 
01876             cpl_msg_error(fctid, "Can not use MEAN model with PLANE method.");
01877             return -3;
01878 
01879         }
01880 
01881         error_code = _giraffe_bias_compute_mean(results, img, biasareas,
01882                 sigma, numiter, maxfraction);
01883 
01884 
01885     }
01886     else if ((method == GIBIAS_METHOD_CURVE) ||
01887              (option == GIBIAS_OPTION_CURVE)) {
01888 
01889         /*
01890          * A 2D Bias Surface is fitted on a grid of the pixels in the
01891          * bias areas
01892          */
01893 
01894         cpl_msg_info(fctid, "using CURVE Model...");
01895 
01896         error_code = _giraffe_bias_compute_curve(results, img, biasareas,
01897                 sigma, numiter, maxfraction, xdeg, ydeg, xstep, ystep);
01898 
01899         if (error_code != EXIT_SUCCESS) {
01900 
01901             cpl_msg_error(fctid, "Error during calculation of bias curve, "
01902                           "aborting...");
01903             return error_code;
01904 
01905         }
01906 
01907         coeff = results->coeffs;
01908 
01909     }
01910     else {
01911 
01912         /*
01913          *  A 2D Bias Plane is fitted on the pixels of the biaslimits area
01914          */
01915 
01916         cpl_msg_info(fctid, "using PLANE (FITTED) model\n");
01917 
01918         error_code = _giraffe_bias_compute_plane(results, img, biasareas,
01919                 sigma, numiter, maxfraction);
01920 
01921         if (error_code == EXIT_SUCCESS) {
01922 
01923             coeff = results->coeffs;
01924 
01925             results->mean = cpl_matrix_get(coeff, 0, 0) +
01926                 cpl_matrix_get(coeff, 0, 1) * img_centerx +
01927                 cpl_matrix_get(coeff, 0, 2) * img_centery;
01928 
01929         }
01930         else {
01931 
01932             cpl_msg_error(fctid, "Error during calculation of bias plane, "
01933                           "aborting...");
01934             return error_code;
01935 
01936         }
01937 
01938     }
01939 
01940 
01941     /*
01942      *  Perform computation based on method....
01943      */
01944 
01945     switch (method) {
01946         case GIBIAS_METHOD_UNIFORM:
01947         {
01948 
01949             cpl_msg_info(fctid, "using bias method == UNIFORM");
01950 
01951 
01952             /*
01953              * Subtract the average bias value
01954              */
01955 
01956             if (model == GIBIAS_MODEL_MEAN) {
01957 
01958                 cpl_msg_info(fctid, "bias value(areas mean) = %f, "
01959                         "bias sigma = %f", results->mean, results->sigma);
01960 
01961             }
01962             else {
01963 
01964                 cpl_msg_info(fctid, "bias value at (%d, %d) = %f, "
01965                         "bias sigma = %f", img_centerx, img_centery,
01966                         results->mean, results->sigma);
01967 
01968             }
01969 
01970             if (remove_bias == TRUE) {
01971 
01972                 cpl_image_subtract_scalar(img, results->mean);
01973 
01974             }
01975 
01976             break;
01977         }
01978 
01979         case GIBIAS_METHOD_PLANE:
01980         {
01981 
01982             cpl_msg_info(fctid, "using bias method == PLANE");
01983 
01984             if (coeff == NULL) {
01985 
01986                 error_code = _giraffe_bias_compute_plane(results, img,
01987                         biasareas, sigma, numiter, maxfraction);
01988 
01989                 if (error_code == EXIT_SUCCESS) {
01990 
01991                     coeff = results->coeffs;
01992 
01993                     results->mean = cpl_matrix_get(coeff, 0, 0) +
01994                         cpl_matrix_get(coeff, 0, 1) * img_centerx +
01995                         cpl_matrix_get(coeff, 0, 2) * img_centery;
01996 
01997                 }
01998                 else {
01999 
02000                     cpl_msg_error(fctid, "Error during calculation of bias "
02001                                   "plane, aborting...");
02002                     return error_code;
02003 
02004                 }
02005 
02006             }
02007 
02008             cpl_msg_info(fctid, "bias value at (%d, %d) = %f, bias sigma = %f, "
02009                          "bias plane = %g + %g * x + %g * y", img_centerx,
02010                          img_centery, results->mean, results->sigma,
02011                          cpl_matrix_get(coeff, 0, 0),
02012                          cpl_matrix_get(coeff, 0, 1),
02013                          cpl_matrix_get(coeff, 0, 2));
02014 
02015             if (remove_bias == TRUE) {
02016 
02017                 //TODO Move plane removal to gir_image...
02018 
02019                 cxdouble* _img = cpl_image_get_data_double(img);
02020 
02021 #if 0
02022                 cpl_image* bsmodel = cpl_image_new(img_dimx, img_dimy,
02023                                                    CPL_TYPE_DOUBLE);
02024                 cxdouble* _bsmodel = cpl_image_get_data_double(bsmodel);
02025 #endif
02026 
02027 
02028                 for (j = 0; j < img_dimy; j++) {
02029 
02030                     cxsize jj = j * img_dimx;
02031 
02032                     for (i = 0; i < img_dimx; i++) {
02033 
02034 #if 0
02035                         _bsmodel[jj + i] = cpl_matrix_get(coeff, 0, 0)
02036                                 + cpl_matrix_get(coeff, 0, 1) * i
02037                                 + cpl_matrix_get(coeff, 0, 2) * j;
02038 #endif
02039 
02040                         _img[jj + i] -= cpl_matrix_get(coeff, 0, 0)
02041                                 + cpl_matrix_get(coeff, 0, 1) * i
02042                                 + cpl_matrix_get(coeff, 0, 2) * j;
02043 
02044                     }
02045 
02046                 }
02047 
02048 #if 0
02049                 cpl_image_save(bsmodel, "idiot.fits", CPL_BPP_IEEE_FLOAT,
02050                                0, CPL_IO_DEFAULT);
02051                 cpl_image_delete(bsmodel);
02052 #endif
02053 
02054             }
02055 
02056             break;
02057 
02058         }
02059 
02060         case GIBIAS_METHOD_CURVE:
02061         {
02062 
02063             cpl_msg_info(fctid, "using bias method == CURVE");
02064 
02065             if (coeff == NULL) {
02066 
02067                 error_code = _giraffe_bias_compute_curve(results, img,
02068                         biasareas, sigma, numiter, maxfraction, xdeg, ydeg,
02069                         xstep, ystep);
02070 
02071                 if (error_code != EXIT_SUCCESS) {
02072 
02073                     cpl_msg_error(fctid, "Error during calculation of bias "
02074                                   "surface curve, aborting...");
02075                     return error_code;
02076 
02077                 }
02078 
02079                 coeff = results->coeffs;
02080 
02081             }
02082 
02083             cpl_msg_info(fctid, "bias value %f, bias sigma = %f\n",
02084                     results->mean, results->sigma);
02085 
02086             //TODO Should we print the matrix as well????
02087 
02088 
02089             if (remove_bias == TRUE) {
02090 
02091                 cxdouble* _img = cpl_image_get_data_double(img);
02092 
02093 #if 0
02094                 cpl_image* bsmodel = cpl_image_new(img_dimx, img_dimy,
02095                                                    CPL_TYPE_DOUBLE);
02096                 cxdouble* _bsmodel = cpl_image_get_data_double(bsmodel);
02097 #endif
02098 
02099 
02100                 tmp_x = cpl_matrix_new(1, 1);
02101                 tmp_y = cpl_matrix_new(1, 1);
02102 
02103                 for (j = 0; j < img_dimy; j++) {
02104 
02105                     cxsize jj = j * img_dimx;
02106 
02107                     cpl_matrix_set(tmp_y, 0, 0, j);
02108 
02109                     for (i = 0; i < img_dimx; i++) {
02110 
02111                         cpl_matrix_set(tmp_x, 0, 0, i);
02112                         tmp_fit = giraffe_chebyshev_fit2d(0., 0., img_dimx,
02113                                 img_dimy, coeff, tmp_x, tmp_y);
02114 
02115 #if 0
02116                         _bsmodel[jj + i] = cpl_matrix_get(tmp_fit, 0, 0);
02117 #endif
02118                         _img[jj + i] -= cpl_matrix_get(tmp_fit, 0, 0);
02119 
02120                         cpl_matrix_delete(tmp_fit);
02121 
02122                     }
02123 
02124                 }
02125 
02126                 cpl_matrix_delete(tmp_x);
02127                 cpl_matrix_delete(tmp_y);
02128 
02129 #if 0
02130                 cpl_image_save(bsmodel, "idiot.fits", CPL_BPP_IEEE_FLOAT,
02131                                0, CPL_IO_DEFAULT);
02132                 cpl_image_delete(bsmodel);
02133 #endif
02134             }
02135 
02136             break;
02137 
02138         }
02139 
02140         case GIBIAS_METHOD_PROFILE:
02141         {
02142 
02143             cpl_msg_info(fctid, "using bias method 'PROFILE'");
02144 
02145             error_code = _giraffe_bias_compute_profile(results, img,
02146                                                        biasareas, 'y');
02147 
02148             if (error_code != EXIT_SUCCESS) {
02149 
02150                 cpl_msg_error(fctid, "Error computing the bias area(s) "
02151                               "profile");
02152                 return error_code;
02153 
02154             }
02155 
02156             if (remove_bias == TRUE) {
02157 
02158                 const cxdouble* _bias =
02159                     cpl_image_get_data_double(results->model);
02160 
02161                 cxdouble* _img = cpl_image_get_data_double(img);
02162 
02163 
02164                 cx_assert(_bias != NULL);
02165                 cx_assert(_img != NULL);
02166 
02167                 for (j = 0; j < img_dimy; ++j) {
02168 
02169                     cxsize jj = j * img_dimx;
02170 
02171 
02172                     for (i = 0; i < img_dimx; ++i) {
02173                         _img[jj + i] -= _bias[j];
02174                     }
02175 
02176                 }
02177 
02178             }
02179 
02180             break;
02181 
02182         }
02183 
02184         case GIBIAS_METHOD_MASTER:
02185         case GIBIAS_METHOD_ZMASTER:
02186         {
02187 
02188             cxdouble biasdrift = 0.;
02189 
02190             if (master_bias == NULL) {
02191 
02192                 cpl_msg_error(fctid, "Master Bias Frame required with MASTER "
02193                               "method.");
02194                 return -5;
02195 
02196             }
02197 
02198             master_bias_dimx = cpl_image_get_size_x(master_bias);
02199             master_bias_dimy = cpl_image_get_size_y(master_bias);
02200 
02201             if ((master_bias_dimx != img_dimx) ||
02202                 (master_bias_dimy != img_dimy)) {
02203 
02204                 cpl_msg_error(fctid, "Invalid master bias! Size should "
02205                         "be [%d, %d] but is [%d, %d].", img_dimx, img_dimy,
02206                         master_bias_dimx, master_bias_dimy);
02207                 return -7;
02208 
02209             }
02210 
02211             if (coeff == NULL) {
02212 
02213                 if (option == GIBIAS_OPTION_CURVE) {
02214 
02215                     error_code = _giraffe_bias_compute_curve(results, img,
02216                             biasareas, sigma, numiter, maxfraction, xdeg,
02217                             ydeg, xstep, ystep);
02218 
02219                     if (error_code != EXIT_SUCCESS) {
02220 
02221                         cpl_msg_error(fctid, "Error during calculation of "
02222                                 "bias surface curve, aborting...");
02223                         return error_code;
02224 
02225                     }
02226 
02227                     coeff = results->coeffs;
02228 
02229                 }
02230                 else {
02231 
02232                     /*
02233                      *  Default to plane...
02234                      */
02235 
02236                     error_code = _giraffe_bias_compute_plane(results, img,
02237                             biasareas, sigma, numiter, maxfraction);
02238 
02239                     if (error_code == EXIT_SUCCESS) {
02240 
02241                         coeff = results->coeffs;
02242 
02243                         results->mean = cpl_matrix_get(coeff, 0, 0) +
02244                             cpl_matrix_get(coeff, 0, 1) * img_centerx +
02245                             cpl_matrix_get(coeff, 0, 2) * img_centery;
02246 
02247                     }
02248                     else {
02249 
02250                         cpl_msg_error(fctid, "Error during calculation of "
02251                                 "bias surface plane, aborting...");
02252                         return error_code;
02253 
02254                     }
02255 
02256                 }
02257 
02258             }
02259 
02260             if (method == GIBIAS_METHOD_ZMASTER) {
02261 
02262                 /*
02263                  * A possible drift of the average bias is compensated using
02264                  * the pre/overscans reference value as an additive
02265                  * correction
02266                  */
02267 
02268                 biasdrift = results->mean - mbias;
02269 
02270                 cpl_msg_info(fctid, "METHOD == ZMASTER, drift value %g",
02271                         biasdrift);
02272 
02273             }
02274             else {
02275 
02276                 cpl_msg_info(fctid, "METHOD == MASTER");
02277 
02278             }
02279 
02280 
02281             /*
02282              *  Write some info into log file...
02283              */
02284 
02285             if (option == GIBIAS_OPTION_CURVE) {
02286 
02287                 //TODO dump coeff to log
02288 
02289                 cpl_msg_info(fctid, "OPTION == CURVE, bias mean = %f, bias "
02290                         "sigma = %f", results->mean, results->sigma);
02291 
02292             }
02293             else {
02294 
02295                 if (option == GIBIAS_OPTION_PLANE) {
02296 
02297                     cpl_msg_info(fctid, "OPTION == PLANE, bias plane = %g + "
02298                             "%g * x + %g * y, bias mean at (%d, %d) = %f, "
02299                             "bias sigma = %f\n",
02300                             cpl_matrix_get(coeff, 0, 0),
02301                             cpl_matrix_get(coeff, 0, 1),
02302                             cpl_matrix_get(coeff, 0, 2),
02303                             img_centerx, img_centery, results->mean,
02304                             results->sigma);
02305 
02306                 }
02307                 else {
02308 
02309                     cpl_msg_info(fctid, "bias mean = %f, bias sigma = %f\n",
02310                             results->mean, results->sigma);
02311 
02312                 }
02313 
02314             }
02315 
02316 
02317             /*
02318              *  Handle bad pixel bitmap...
02319              */
02320 
02321             if (bad_pixel_map == NULL) {
02322 
02323                 /*
02324                  *  No bad pixel map was given...
02325                  *  Subtract master bias pixels without modification except
02326                  *  pre/overscan trimming.
02327                  */
02328 
02329                 if (remove_bias == TRUE) {
02330 
02331                     const cxdouble* _mbias =
02332                         cpl_image_get_data_double_const(master_bias);
02333 
02334                     cxdouble* _img = cpl_image_get_data_double(img);
02335 
02336 
02337                     for (j = 0; j < img_dimy; j++) {
02338 
02339                         cxsize jj = j * img_dimx;
02340 
02341                         for (i = 0; i < img_dimx; i++) {
02342                             _img[jj + i] -= (_mbias[jj + i] + biasdrift);
02343                         }
02344 
02345                     }
02346 
02347                 }
02348 
02349             }
02350             else {
02351 
02352                 /*
02353                  *  Bad pixel map was given...
02354                  *  For all pixels not flagged in bad_pixel_map the
02355                  *  master_bias pixel values are subtracted, otherwise
02356                  *  results->mean or bias model image pixels are subtracted
02357                  *  depending on model.
02358                  */
02359 
02360                 bad_pixel_dimx = cpl_image_get_size_x(bad_pixel_map);
02361                 bad_pixel_dimy = cpl_image_get_size_y(bad_pixel_map);
02362 
02363                 if ((bad_pixel_dimx != img_dimx) ||
02364                     (bad_pixel_dimy != img_dimy)) {
02365 
02366                     cpl_msg_error(fctid, "Invalid bad pixel map size should "
02367                             "be [%d, %d] but is [%d, %d].", img_dimx,
02368                             img_dimy, bad_pixel_dimx, bad_pixel_dimy);
02369                     return -6;
02370 
02371                 }
02372 
02373                 if (remove_bias == TRUE) {
02374 
02375                     if (option == GIBIAS_OPTION_CURVE) {
02376 
02377                         const cxint* _bpix =
02378                             cpl_image_get_data_int_const(bad_pixel_map);
02379 
02380                         const cxdouble* _mbias =
02381                             cpl_image_get_data_double_const(master_bias);
02382 
02383                         cxdouble* _img = cpl_image_get_data_double(img);
02384 
02385 
02386                         tmp_x = cpl_matrix_new(1, 1);
02387                         tmp_y = cpl_matrix_new(1, 1);
02388 
02389                         for (j = 0; j < img_dimy; j++) {
02390 
02391                             cxsize jj = j * img_dimx;
02392 
02393                             cpl_matrix_set(tmp_y, 0, 0, j);
02394 
02395                             for (i = 0; i < img_dimx; i++) {
02396 
02397                                 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
02398 
02399                                     _img[jj + i] -= (_mbias[jj + i] +
02400                                             biasdrift);
02401 
02402                                 }
02403                                 else {
02404 
02405                                     cpl_matrix_set(tmp_x, 0, 0, i);
02406                                     tmp_fit = giraffe_chebyshev_fit2d(0., 0.,
02407                                             img_dimx, img_dimy, coeff,
02408                                             tmp_x, tmp_y);
02409 
02410                                     _img[jj + i] -= cpl_matrix_get(tmp_fit,
02411                                             0, 0);
02412 
02413                                     cpl_matrix_delete(tmp_fit);
02414 
02415                                 }
02416 
02417                             }
02418 
02419                         }
02420 
02421                         cpl_matrix_delete(tmp_y);
02422                         cpl_matrix_delete(tmp_x);
02423 
02424                     }
02425                     else if (option == GIBIAS_OPTION_PLANE) {
02426 
02427                         const cxint* _bpix =
02428                             cpl_image_get_data_int_const(bad_pixel_map);
02429 
02430                         const cxdouble* _mbias =
02431                             cpl_image_get_data_double_const(master_bias);
02432 
02433                         cxdouble* _img = cpl_image_get_data_double(img);
02434 
02435 
02436                         for (j = 0; j < img_dimy; j++) {
02437 
02438                             cxsize jj = j * img_dimx;
02439 
02440                             for (i = 0; i < img_dimx; i++) {
02441 
02442                                 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
02443 
02444                                     _img[jj + i] -= (_mbias[jj + i] +
02445                                             biasdrift);
02446 
02447                                 }
02448                                 else {
02449 
02450                                     _img[jj + i] -=
02451                                         cpl_matrix_get(coeff, 0, 0) +
02452                                         cpl_matrix_get(coeff, 0, 1) * j +
02453                                         cpl_matrix_get(coeff, 0, 2) * i;
02454                                 }
02455 
02456                             }
02457 
02458                         }
02459 
02460                     }
02461                     else {
02462 
02463                         const cxint* _bpix =
02464                             cpl_image_get_data_int_const(bad_pixel_map);
02465 
02466                         const cxdouble* _mbias =
02467                             cpl_image_get_data_double_const(master_bias);
02468 
02469                         cxdouble* _img = cpl_image_get_data_double(img);
02470 
02471 
02472                         cpl_msg_info(fctid, "OPTION==UNIFORM\n");
02473 
02474                         for (j = 0; j < img_dimy; j++) {
02475 
02476                             cxsize jj = j * img_dimx;
02477 
02478                             for (i = 0; i < img_dimx; i++) {
02479 
02480                                 if (!(_bpix[jj + i] & GIR_M_PIX_SET)) {
02481 
02482                                     _img[jj + i] -= (_mbias[jj + i] +
02483                                             biasdrift);
02484 
02485                                 }
02486                                 else {
02487 
02488                                     _img[jj + i] -= results->mean;
02489 
02490                                 }
02491 
02492                             }
02493 
02494                         }
02495 
02496                     }
02497 
02498                 }
02499 
02500             }
02501 
02502             break;
02503 
02504         }
02505 
02506         default:
02507         {
02508 
02509             cpl_msg_error(fctid, "Invalid method, aborting...");
02510             return -4;
02511             break;
02512 
02513         }
02514 
02515     }
02516 
02517 
02518     cpl_msg_info(fctid, "Resulting biaslimits : %s",
02519             cx_string_get(results->limits));
02520 
02521     return 0;
02522 
02523 }
02524 
02525 
02543 cpl_matrix *
02544 giraffe_get_raw_areas(const GiImage *image)
02545 {
02546 
02547     const cxchar *const _id = "giraffe_get_raw_areas";
02548 
02549     cxint32 prescx = 0;
02550     cxint32 prescy = 0;
02551     cxint32 ovrscx = 0;
02552     cxint32 ovrscy = 0;
02553 
02554     cxint32 winx = 0;
02555     cxint32 winy = 0;
02556 
02557     cxint32 tprescx = 0;
02558     cxint32 tprescy = 0;
02559     cxint32 tovrscx = 0;
02560     cxint32 tovrscy = 0;
02561 
02562     cxint32 row = 0;
02563 
02564     cpl_matrix *m_tmp;
02565     cpl_propertylist *properties;
02566 
02567 
02568     properties = giraffe_image_get_properties(image);
02569     if (!properties) {
02570         cpl_error_set(_id, CPL_ERROR_NULL_INPUT);
02571         return NULL;
02572     }
02573 
02574     winx = cpl_propertylist_get_int(properties, GIALIAS_WINX);
02575     winy = cpl_propertylist_get_int(properties, GIALIAS_WINY);
02576 
02577     if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
02578         tprescx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
02579         if (tprescx > 0) {
02580             prescx = tprescx;
02581         }
02582     }
02583     if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
02584         tprescy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
02585         if (tprescy > 0) {
02586             prescy = tprescy;
02587         }
02588     }
02589     if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
02590         tovrscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
02591         if (tovrscx > 0) {
02592             ovrscx = tovrscx;
02593         }
02594     }
02595     if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
02596         tovrscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
02597         if (tovrscy > 0) {
02598             ovrscy = tovrscy;
02599         }
02600     }
02601 
02602 
02603     m_tmp = cpl_matrix_new(1, 4);
02604 
02605     if (prescx > 0) {
02606 
02607         cpl_matrix_set(m_tmp, row, 0, 0.0);
02608         cpl_matrix_set(m_tmp, row, 1, (cxdouble) prescx - 1);
02609         cpl_matrix_set(m_tmp, row, 2, 0.0);
02610         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02611 
02612         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02613         row++;
02614     }
02615 
02616     if (ovrscx > 0) {
02617 
02618         cpl_matrix_set(m_tmp, row, 0, (cxdouble) winx - ovrscx);
02619         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
02620         cpl_matrix_set(m_tmp, row, 2, 0.0);
02621         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02622 
02623         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02624         row++;
02625     }
02626 
02627     if (!(prescy == 0 || prescx == 0 || ovrscx == 0)) {
02628 
02629         cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
02630         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
02631         cpl_matrix_set(m_tmp, row, 2, 0.0);
02632         cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
02633 
02634         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02635         row++;
02636 
02637     }
02638     else if (!(prescy == 0 || prescx == 0)) {
02639 
02640         cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
02641         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
02642         cpl_matrix_set(m_tmp, row, 2, 0.0);
02643         cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
02644 
02645         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02646         row++;
02647 
02648     }
02649     else if (!(prescy == 0 || ovrscx == 0)) {
02650 
02651         cpl_matrix_set(m_tmp, row, 0, 0.0);
02652         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
02653         cpl_matrix_set(m_tmp, row, 2, 0.0);
02654         cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
02655 
02656         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02657         row++;
02658 
02659     }
02660     else if (prescy > 0) {
02661 
02662         cpl_matrix_set(m_tmp, row, 0, 0.0);
02663         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
02664         cpl_matrix_set(m_tmp, row, 2, 0.0);
02665         cpl_matrix_set(m_tmp, row, 3, (cxdouble) prescy - 1);
02666 
02667         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02668         row++;
02669     }
02670 
02671     if (!(ovrscy == 0 || prescx == 0 || ovrscx == 0)) {
02672 
02673         cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
02674         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
02675         cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
02676         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02677 
02678         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02679         row++;
02680 
02681     }
02682     else if (!(ovrscy == 0 || prescx == 0)) {
02683 
02684         cpl_matrix_set(m_tmp, row, 0, (cxdouble) prescx);
02685         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
02686         cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
02687         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02688 
02689         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02690         row++;
02691 
02692     }
02693     else if (!(ovrscy == 0 || ovrscx == 0)) {
02694 
02695         cpl_matrix_set(m_tmp, row, 0, 0.0);
02696         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - ovrscx - 1);
02697         cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
02698         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02699 
02700         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02701         row++;
02702 
02703     }
02704     else if (ovrscy > 0) {
02705 
02706         cpl_matrix_set(m_tmp, row, 0, 0.0);
02707         cpl_matrix_set(m_tmp, row, 1, (cxdouble) winx - 1);
02708         cpl_matrix_set(m_tmp, row, 2, (cxdouble) winy - ovrscy);
02709         cpl_matrix_set(m_tmp, row, 3, (cxdouble) winy - 1);
02710 
02711         cpl_matrix_resize(m_tmp, 0, 1, 0, 0);
02712         row++;
02713     }
02714 
02715     cpl_matrix_resize(m_tmp, 0, -1, 0, 0);
02716     row--;
02717 
02718     if (row == 0) {
02719         cpl_matrix_delete(m_tmp);
02720         return NULL;
02721     }
02722 
02723     return m_tmp;
02724 }
02725 
02726 
02745 cxint
02746 giraffe_trim_raw_areas(GiImage *image)
02747 {
02748 
02749     const cxchar *fctid = "giraffe_trim_raw_areas";
02750 
02751     cxint nx, ny;
02752     cxint newlx, newly, newux, newuy;
02753 
02754     cxint ovscx = 0;
02755     cxint ovscy = 0;
02756     cxint prscx = 0;
02757     cxint prscy = 0;
02758 
02759     cpl_propertylist *properties = giraffe_image_get_properties(image);
02760 
02761     cpl_image *_image = giraffe_image_get(image);
02762     cpl_image *tmp;
02763 
02764 
02765     if (!properties) {
02766         cpl_msg_error(fctid, "Image does not contain any properties!");
02767         return 1;
02768     }
02769 
02770     nx = cpl_image_get_size_x(_image);
02771     ny = cpl_image_get_size_y(_image);
02772 
02773     if (cpl_propertylist_has(properties, GIALIAS_NAXIS1)) {
02774         cxint _nx = cpl_propertylist_get_int(properties, GIALIAS_NAXIS1);
02775 
02776         if (_nx != nx) {
02777             cpl_msg_warning(fctid, "Image size (%d) and image property `%s' "
02778                             "(%d) are not consistent! Using image size ...",
02779                             nx, GIALIAS_NAXIS1, _nx);
02780         }
02781     }
02782     else {
02783         cpl_msg_warning(fctid, "Image does not contain any property `%s'. "
02784                         "Using image size (%d)", GIALIAS_NAXIS1, nx);
02785     }
02786 
02787 
02788     if (cpl_propertylist_has(properties, GIALIAS_NAXIS2)) {
02789         cxint _ny = cpl_propertylist_get_int(properties, GIALIAS_NAXIS2);
02790 
02791         if (_ny != ny) {
02792             cpl_msg_warning(fctid, "Image size (%d) and image property `%s' "
02793                             "(%d) are not consistent! Using image size ...",
02794                             ny, GIALIAS_NAXIS2, _ny);
02795         }
02796     }
02797     else {
02798         cpl_msg_warning(fctid, "Image does not contain any property `%s'. "
02799                         "Using image size (%d)", GIALIAS_NAXIS2, ny);
02800     }
02801 
02802     if (cpl_propertylist_has(properties, GIALIAS_OVSCX)) {
02803         ovscx = cpl_propertylist_get_int(properties, GIALIAS_OVSCX);
02804     }
02805 
02806     if (cpl_propertylist_has(properties, GIALIAS_OVSCY)) {
02807         ovscy = cpl_propertylist_get_int(properties, GIALIAS_OVSCY);
02808     }
02809 
02810     if (cpl_propertylist_has(properties, GIALIAS_PRSCX)) {
02811         prscx = cpl_propertylist_get_int(properties, GIALIAS_PRSCX);
02812     }
02813 
02814     if (cpl_propertylist_has(properties, GIALIAS_PRSCY)) {
02815         prscy = cpl_propertylist_get_int(properties, GIALIAS_PRSCY);
02816     }
02817 
02818     // FIXME:
02819     //  Check that the one pixel offset is ok. cpl_image_extract starts
02820     //  counting from 1,1 with right and to boundaries inclusive.
02821 
02822     newlx = prscx + 1;
02823     newly = prscy + 1;
02824     newux = nx - ovscx;
02825     newuy = ny - ovscy;
02826 
02827     tmp = cpl_image_extract(_image, newlx, newly, newux, newuy);
02828 
02829     giraffe_image_set(image, tmp);
02830     cpl_image_delete(tmp);
02831 
02832     _image = giraffe_image_get(image);
02833 
02834     nx = cpl_image_get_size_x(_image);
02835     ny = cpl_image_get_size_y(_image);
02836 
02837 
02838     /*
02839      * Update image properties
02840      */
02841 
02842     cpl_propertylist_set_int(properties, GIALIAS_NAXIS1, nx);
02843     cpl_propertylist_set_int(properties, GIALIAS_NAXIS2, ny);
02844 
02845     if (cpl_propertylist_has(properties, GIALIAS_CRPIX1)) {
02846         cxdouble crpix1 = cpl_propertylist_get_double(properties,
02847             GIALIAS_CRPIX1);
02848 
02849         /*
02850          * For GIRAFFE the reference pixel is defined as
02851          * 1 - width(prescan)
02852          */
02853 
02854         crpix1 += (cxdouble)prscx;
02855         cpl_propertylist_set_double(properties, GIALIAS_CRPIX1, crpix1);
02856     }
02857 
02858     if (cpl_propertylist_has(properties, GIALIAS_CRPIX2)) {
02859         cxdouble crpix2 = cpl_propertylist_get_double(properties,
02860             GIALIAS_CRPIX2);
02861 
02862         crpix2 -= (cxdouble) prscy;
02863         cpl_propertylist_set_double(properties, GIALIAS_CRPIX2, crpix2);
02864     }
02865 
02866     cpl_propertylist_erase(properties, GIALIAS_OVSCX);
02867     cpl_propertylist_erase(properties, GIALIAS_OVSCY);
02868     cpl_propertylist_erase(properties, GIALIAS_PRSCX);
02869     cpl_propertylist_erase(properties, GIALIAS_PRSCY);
02870 
02871     return 0;
02872 
02873 }
02874 
02875 
02937 cxint
02938 giraffe_bias_remove(GiImage* result, const GiImage* raw,
02939                     const GiImage* master_bias, const GiImage* bad_pixels,
02940                     const cpl_matrix* biaslimits, const GiBiasConfig* config)
02941 {
02942 
02943     const cxchar* const fctid = "giraffe_bias_remove";
02944 
02945     cxint error_code;
02946 
02947     cxdouble bias_value = 0.;
02948 
02949     cx_string* s;
02950 
02951     cpl_matrix* areas = NULL;
02952 
02953     cpl_image* _raw = giraffe_image_get(raw);
02954     cpl_image* _master_bias = giraffe_image_get(master_bias);
02955     cpl_image* _bad_pixels = giraffe_image_get(bad_pixels);
02956     cpl_image* timage;
02957 
02958     cpl_propertylist* properties;
02959 
02960     GiBiasResults  bias_results = {0., 0., 0., NULL, NULL, NULL};
02961 
02962 
02963     cx_assert(raw != NULL);
02964     cx_assert(config != NULL);
02965     cx_assert(biaslimits == NULL);
02966 
02967     if (result == NULL) {
02968         return 1;
02969     }
02970 
02971 
02972     /*
02973      * Setup user defined areas to use for the bias computation
02974      */
02975 
02976     if (strncmp(config->areas, "None", 4) == 0) {
02977 
02978         cpl_msg_info(fctid, "No bias areas specified, using pre/overscan"
02979                 "regions of the raw frame.");
02980 
02981         cpl_error_reset();
02982         areas = giraffe_get_raw_areas(raw);
02983         if (areas == NULL) {
02984             if (cpl_error_get_code() == CPL_ERROR_NULL_INPUT) {
02985                 cpl_msg_error(fctid, "Invalid raw image! Image does not "
02986                               "contain any properties!");
02987             }
02988             else {
02989                 cpl_msg_error(fctid, "Invalid raw image! Image does not "
02990                               "contain or has invalid pre- and overscan "
02991                               "properties.");
02992             }
02993 
02994             return 1;
02995         }
02996 
02997     }
02998     else {
02999 
03000         areas = _giraffe_bias_get_areas(config->areas);
03001 
03002         if (areas == NULL) {
03003 
03004             cpl_msg_error(fctid, "Invalid bias area specification '%s'!",
03005                     config->areas);
03006             return 1;
03007 
03008         }
03009 
03010         cpl_msg_info(fctid, "Using bias area(s) '%s' for bias computation",
03011                 config->areas);
03012 
03013     }
03014 
03015 
03016     /*
03017      * Processing
03018      */
03019 
03020 
03021     /*
03022      *  Check whether a masterbias image was specified
03023      *  If so, compare pre/ovrscan keywords/areas
03024      */
03025 
03026     if (master_bias != NULL) {
03027         cpl_propertylist *_properties;
03028         cxbool is_same_size = FALSE;
03029 
03030         is_same_size = _giraffe_compare_overscans(raw, master_bias);
03031 
03032         // FIXME:
03033         //  Fully processed calibrations usually do not have overscans
03034         //  any longer. Instead of failing at this point some dummies
03035         //  could be created.
03036 
03037         if (is_same_size==FALSE) {
03038             cpl_msg_info(fctid, "PRE/OVRSCAN Regions do not match between "
03039                          "RAW Image and Masterbias Image.");
03040 
03041             return 1;
03042         }
03043 
03044         _properties = giraffe_image_get_properties(master_bias);
03045 
03046         if (cpl_propertylist_has(_properties, GIALIAS_BIASVALUE)) {
03047             bias_value = cpl_propertylist_get_double(_properties,
03048                                                      GIALIAS_BIASVALUE);
03049         }
03050     }
03051 
03052 
03053     /*
03054      *  Check whether a Bad Pixel Map image was specified
03055      *  If so, compare pre/ovrscan keywords/areas
03056      */
03057 
03058     if (bad_pixels != NULL) {
03059         cxbool is_same_size = FALSE;
03060 
03061         is_same_size = _giraffe_compare_overscans(raw, bad_pixels);
03062 
03063         // FIXME:
03064         //  Fully processed calibrations usually do not have overscans
03065         //  any longer. Instead of failing at this point some dummies
03066         //  could be created.
03067 
03068         if (is_same_size == FALSE) {
03069             cpl_msg_info(fctid, "PRE/OVRSCAN Regions do not match between "
03070                          "RAW Image and Bad Pixel Image.");
03071 
03072             return 1;
03073         }
03074     }
03075 
03076 
03077     /*
03078      *  Create a copy of the raw frame. We are not working on the raw
03079      *  frame directly.
03080      */
03081 
03082     // FIXME:
03083     //   Can we use the raw frame directly? (RP)
03084 
03085     timage = cpl_image_duplicate(_raw);
03086 
03087 
03088     /*
03089      *  Remove Bias Calculation. The bias corrected image will be
03090      *  stored in timage.
03091      */
03092 
03093     error_code = _giraffe_bias_compute(&bias_results, timage,
03094                                        areas, _master_bias,
03095                                        _bad_pixels, config->method,
03096                                        config->option, config->model,
03097                                        config->remove, bias_value,
03098                                        config->xdeg, config->ydeg,
03099                                        config->xstep, config->ystep,
03100                                        config->sigma, config->iterations,
03101                                        config->fraction);
03102 
03103     cpl_matrix_delete(areas);
03104     areas = NULL;
03105 
03106 
03107     // FIXME:
03108     //  Posslibly check returned error code and do a dedicated
03109     //  exception handling (RP)
03110 
03111     if (error_code) {
03112 
03113         _giraffe_biasresults_clear(&bias_results);
03114 
03115         cpl_msg_error(fctid, "Bias computation failed with error %d!",
03116                       error_code);
03117         return 1;
03118 
03119     }
03120 
03121 
03122     /*
03123      *  Should we remove Bias or not?
03124      */
03125 
03126     properties = giraffe_image_get_properties(raw);
03127 
03128     if (config->remove == TRUE) {
03129 
03130         giraffe_image_set(result, timage);
03131         cpl_image_delete(timage);
03132 
03133         giraffe_image_set_properties(result, properties);
03134 
03135     }
03136     else {
03137 
03138         cpl_image_delete(timage);
03139 
03140         giraffe_image_set(result, _raw);
03141         giraffe_image_set_properties(result, properties);
03142 
03143     }
03144 
03145 
03146     /*
03147      * Postprocessing
03148      */
03149 
03150     properties = giraffe_image_get_properties(result);
03151 
03152     if (config->remove == TRUE) {
03153 
03154         cpl_propertylist_set_int(properties, GIALIAS_BITPIX, -32);
03155         cpl_propertylist_set_double(properties, GIALIAS_BZERO, 0.);
03156         cpl_propertylist_set_double(properties, GIALIAS_BSCALE, 1.);
03157 
03158         if (cpl_propertylist_has(properties, GIALIAS_GIRFTYPE)) {
03159             cpl_propertylist_set_string(properties,
03160                                         GIALIAS_GIRFTYPE, "BRMIMG");
03161         }
03162         else {
03163             cpl_propertylist_append_string(properties,
03164                                            GIALIAS_GIRFTYPE, "BRMIMG");
03165         }
03166         cpl_propertylist_set_comment(properties, GIALIAS_GIRFTYPE,
03167                               "GIRAFFE bias subtracted frame");
03168 
03169         /*
03170          * Remove pre- and overscans
03171          */
03172 
03173         giraffe_trim_raw_areas(result);
03174 
03175 
03176         /*
03177          * Convert bias corrected image to electrons
03178          */
03179 
03180         // FIXME:
03181         //  Is this really needed? Disabled on request of DFO (RP)
03182 #if 0
03183         if (cpl_propertylist_has(properties, GIALIAS_CONAD)) {
03184             cxdouble conad = cpl_propertylist_get_double(properties, GIALIAS_CONAD);
03185 
03186             if (conad > 0.) {
03187 
03188                 cpl_image* _result = giraffe_image_get(result);
03189 
03190                 cpl_msg_info(fctid, "Converting bias subtracted frame to "
03191                              "electrons; conversion factor is %.2f e-/ADU",
03192                              conad);
03193                 cpl_image_multiply_scalar(_result, conad);
03194             }
03195             else {
03196                 cpl_msg_warning(fctid, "Invalid conversion factor %.2f "
03197                                 "e-/ADU! Bias subtracted image will not be "
03198                                 "converted.", conad);
03199             }
03200         }
03201         else {
03202             cpl_msg_info(fctid, "Conversion factor not found. Bias subtracted "
03203                          "image will remain in ADU");
03204         }
03205 #endif
03206     }
03207 
03208     s = cx_string_new();
03209 
03210     _giraffe_method_string(s, config->method, config->option);
03211     cpl_propertylist_append_string(properties, GIALIAS_BIASMETHOD, cx_string_get(s));
03212 
03213     cx_string_delete(s);
03214 
03215     cpl_propertylist_append_double(properties, GIALIAS_BIASVALUE,
03216                             bias_results.mean);
03217     cpl_propertylist_append_double(properties, GIALIAS_BIASERROR,
03218                             bias_results.sigma);
03219     cpl_propertylist_append_double(properties, GIALIAS_BIASSIGMA,
03220                             bias_results.rms);
03221 
03222     if (bias_results.coeffs) {
03223         cpl_propertylist_append_double(properties, GIALIAS_BCLIPSIGMA,
03224                                 config->sigma);
03225         cpl_propertylist_append_int(properties, GIALIAS_BCLIPNITER,
03226                                 config->iterations);
03227         cpl_propertylist_append_double(properties, GIALIAS_BCLIPMFRAC,
03228                                 config->fraction);
03229     }
03230 
03231     if (bias_results.limits) {
03232         cpl_propertylist_append_string(properties, GIALIAS_BIASAREAS,
03233                                 cx_string_get(bias_results.limits));
03234     }
03235 
03236     if (bias_results.coeffs) {
03237         s = cx_string_new();
03238 
03239         _giraffe_stringify_coefficients(s, bias_results.coeffs);
03240         cpl_propertylist_append_string(properties, GIALIAS_BIASPLANE,
03241                                 cx_string_get(s));
03242 
03243         cx_string_delete(s);
03244     }
03245 
03246 
03247     /*
03248      * Cleanup
03249      */
03250 
03251     _giraffe_biasresults_clear(&bias_results);
03252 
03253     return 0;
03254 
03255 }
03256 
03257 
03268 GiBiasConfig*
03269 giraffe_bias_config_create(cpl_parameterlist* list)
03270 {
03271 
03272     const cxchar* s;
03273     cpl_parameter* p;
03274 
03275     GiBiasConfig* config = NULL;
03276 
03277 
03278     if (!list) {
03279         return NULL;
03280     }
03281 
03282     config = cx_calloc(1, sizeof *config);
03283 
03284 
03285     /*
03286      * Some defaults
03287      */
03288 
03289     config->method = GIBIAS_METHOD_UNDEFINED;
03290     config->option = GIBIAS_OPTION_UNDEFINED;
03291     config->model = GIBIAS_MODEL_UNDEFINED;
03292     config->mbias = 0.;
03293     config->xdeg = 1;
03294     config->ydeg = 1;
03295 
03296 
03297     p = cpl_parameterlist_find(list, "giraffe.biasremoval.remove");
03298     config->remove = cpl_parameter_get_bool(p);
03299 
03300     p = cpl_parameterlist_find(list, "giraffe.biasremoval.method");
03301     s = cpl_parameter_get_string(p);
03302 
03303     if (!strcmp(s, "UNIFORM")) {
03304         config->method = GIBIAS_METHOD_UNIFORM;
03305     }
03306 
03307     if (!strcmp(s, "PLANE")) {
03308         config->method = GIBIAS_METHOD_PLANE;
03309     }
03310 
03311     if (!strcmp(s, "PROFILE")) {
03312         config->method = GIBIAS_METHOD_PROFILE;
03313     }
03314 
03315     if (!strcmp(s, "MASTER")) {
03316         config->method = GIBIAS_METHOD_MASTER;
03317     }
03318 
03319     if (!strcmp(s, "ZMASTER")) {
03320         config->method = GIBIAS_METHOD_ZMASTER;
03321     }
03322 
03323     if (!strcmp(s, "CURVE")) {
03324         config->method = GIBIAS_METHOD_CURVE;
03325     }
03326 
03327     if (!strcmp(s, "MASTER+PLANE")) {
03328         config->method = GIBIAS_METHOD_MASTER;
03329         config->option = GIBIAS_OPTION_PLANE;
03330     }
03331 
03332     if (!strcmp(s, "ZMASTER+PLANE")) {
03333         config->method = GIBIAS_METHOD_ZMASTER;
03334         config->option = GIBIAS_OPTION_PLANE;
03335     }
03336 
03337     if (!strcmp(s, "MASTER+CURVE")) {
03338         config->method = GIBIAS_METHOD_MASTER;
03339         config->option = GIBIAS_OPTION_CURVE;
03340     }
03341 
03342     if (!strcmp(s, "ZMASTER+CURVE")) {
03343         config->method = GIBIAS_METHOD_ZMASTER;
03344         config->option = GIBIAS_OPTION_CURVE;
03345     }
03346 
03347     cx_assert(config->method != GIBIAS_METHOD_UNDEFINED);
03348 
03349     p = cpl_parameterlist_find(list, "giraffe.biasremoval.areas");
03350     s = cpl_parameter_get_string(p);
03351     config->areas = cx_strdup(s);
03352 
03353     p = cpl_parameterlist_find(list, "giraffe.biasremoval.sigma");
03354     config->sigma = cpl_parameter_get_double(p);
03355 
03356     p = cpl_parameterlist_find(list, "giraffe.biasremoval.iterations");
03357     config->iterations = cpl_parameter_get_int(p);
03358 
03359     p = cpl_parameterlist_find(list, "giraffe.biasremoval.fraction");
03360     config->fraction = cpl_parameter_get_double(p);
03361 
03362     if (config->method == GIBIAS_METHOD_CURVE ||
03363         config->option == GIBIAS_OPTION_CURVE) {
03364         p = cpl_parameterlist_find(list, "giraffe.biasremoval.xorder");
03365         config->xdeg = 1 + cpl_parameter_get_int(p);
03366 
03367         p = cpl_parameterlist_find(list, "giraffe.biasremoval.yorder");
03368         config->ydeg = 1 + cpl_parameter_get_int(p);
03369     }
03370 
03371     p = cpl_parameterlist_find(list, "giraffe.biasremoval.xstep");
03372     config->xstep = cpl_parameter_get_int(p);
03373 
03374     p = cpl_parameterlist_find(list, "giraffe.biasremoval.ystep");
03375     config->ystep = cpl_parameter_get_int(p);
03376 
03377     return config;
03378 
03379 }
03380 
03381 
03394 void
03395 giraffe_bias_config_destroy(GiBiasConfig* config)
03396 {
03397 
03398     if (config) {
03399         if (config->areas) {
03400             cx_free(config->areas);
03401             config->areas = NULL;
03402         }
03403 
03404         cx_free(config);
03405     }
03406 
03407     return;
03408 }
03409 
03410 
03422 void
03423 giraffe_bias_config_add(cpl_parameterlist* list)
03424 {
03425 
03426     cpl_parameter* p;
03427 
03428 
03429     if (!list) {
03430         return;
03431     }
03432 
03433     p = cpl_parameter_new_value("giraffe.biasremoval.remove",
03434                                 CPL_TYPE_BOOL,
03435                                 "Enable bias removal",
03436                                 "giraffe.biasremoval",
03437                                 TRUE);
03438     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "remove-bias");
03439     cpl_parameterlist_append(list, p);
03440 
03441 
03442     p = cpl_parameter_new_enum("giraffe.biasremoval.method",
03443                                CPL_TYPE_STRING,
03444                                "Bias removal method",
03445                                "giraffe.biasremoval",
03446                                "PROFILE", 10, "UNIFORM", "PLANE", "PROFILE",
03447                                "MASTER", "ZMASTER", "MASTER+PLANE",
03448                                "ZMASTER+PLANE", "CURVE", "MASTER+CURVE",
03449                                "ZMASTER+CURVE");
03450     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-method");
03451     cpl_parameterlist_append(list, p);
03452 
03453 
03454     p = cpl_parameter_new_value("giraffe.biasremoval.areas",
03455                                 CPL_TYPE_STRING,
03456                                 "Bias areas to use "
03457                                 "(Xl0:Xr0:Yl0:Yr0, ... ,Xln:Xrn:Yln:Yrn)",
03458                                 "giraffe.biasremoval",
03459                                 "5:40:0:4095");
03460     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-areas");
03461     cpl_parameterlist_append(list, p);
03462 
03463 
03464     p = cpl_parameter_new_value("giraffe.biasremoval.sigma",
03465                                 CPL_TYPE_DOUBLE,
03466                                 "Sigma Clipping: sigma threshold factor",
03467                                 "giraffe.biasremoval",
03468                                 2.5);
03469     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-sigma");
03470     cpl_parameterlist_append(list, p);
03471 
03472 
03473     p = cpl_parameter_new_value("giraffe.biasremoval.iterations",
03474                                 CPL_TYPE_INT,
03475                                 "Sigma Clipping: maximum number of "
03476                                 "iterations",
03477                                 "giraffe.biasremoval",
03478                                 5);
03479     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-niter");
03480     cpl_parameterlist_append(list, p);
03481 
03482 
03483     p = cpl_parameter_new_value("giraffe.biasremoval.fraction",
03484                                 CPL_TYPE_DOUBLE,
03485                                 "Sigma Clipping: minimum fraction of points "
03486                                 "accepted/total [0.0..1.0]",
03487                                 "giraffe.biasremoval",
03488                                 0.8);
03489     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-mfrac");
03490     cpl_parameterlist_append(list, p);
03491 
03492 
03493     p = cpl_parameter_new_value("giraffe.biasremoval.xorder",
03494                                 CPL_TYPE_INT,
03495                                 "Order of X polynomial fit (CURVE method "
03496                                 "only)",
03497                                 "giraffe.biasremoval",
03498                                 1);
03499     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-xorder");
03500     cpl_parameterlist_append(list, p);
03501 
03502     p = cpl_parameter_new_value("giraffe.biasremoval.yorder",
03503                                 CPL_TYPE_INT,
03504                                 "Order of Y polynomial fit (CURVE method "
03505                                 "only)",
03506                                 "giraffe.biasremoval",
03507                                 1);
03508     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-yorder");
03509     cpl_parameterlist_append(list, p);
03510 
03511 
03512     p = cpl_parameter_new_value("giraffe.biasremoval.xstep",
03513                                 CPL_TYPE_INT,
03514                                 "Sampling step along X (CURVE method only)",
03515                                 "giraffe.biasremoval",
03516                                 1);
03517     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-xstep");
03518     cpl_parameterlist_append(list, p);
03519 
03520 
03521     p = cpl_parameter_new_value("giraffe.biasremoval.ystep",
03522                                 CPL_TYPE_INT,
03523                                 "Sampling step along Y (CURVE method only)",
03524                                 "giraffe.biasremoval",
03525                                 1);
03526     cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, "bsremove-ystep");
03527     cpl_parameterlist_append(list, p);
03528 
03529     return;
03530 }

This file is part of the GIRAFFE Pipeline Reference Manual 2.8.6.
Documentation copyright © 2002-2006 European Southern Observatory.
Generated on Mon Jun 28 12:56:56 2010 by doxygen 1.5.1 written by Dimitri van Heesch, © 1997-2004