uves_utils_wrappers.c

00001 /*                                                                              *
00002  *   This file is part of the ESO UVES Pipeline                                 *
00003  *   Copyright (C) 2004,2005 European Southern Observatory                      *
00004  *                                                                              *
00005  *   This library is free software; you can redistribute it and/or modify       *
00006  *   it under the terms of the GNU General Public License as published by       *
00007  *   the Free Software Foundation; either version 2 of the License, or          *
00008  *   (at your option) any later version.                                        *
00009  *                                                                              *
00010  *   This program is distributed in the hope that it will be useful,            *
00011  *   but WITHOUT ANY WARRANTY; without even the implied warranty of             *
00012  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the              *
00013  *   GNU General Public License for more details.                               *
00014  *                                                                              *
00015  *   You should have received a copy of the GNU General Public License          *
00016  *   along with this program; if not, write to the Free Software                *
00017  *   Foundation, 51 Franklin St, Fifth Floor, Boston, MA  02111-1307  USA       *
00018  *                                                                              */
00019 
00020 /*
00021  * $Author: amodigli $
00022  * $Date: 2010/09/24 09:32:09 $
00023  * $Revision: 1.120 $
00024  * $Name: uves-4_9_1 $
00025  */
00026 
00027 #ifdef HAVE_CONFIG_H
00028 #  include <config.h>
00029 #endif
00030 
00031 /*----------------------------------------------------------------------------*/
00039 /*----------------------------------------------------------------------------*/
00040 
00041 #include <uves_utils_wrappers.h>
00042 
00043 #include <uves_utils.h>
00044 #include <uves_utils_cpl.h>
00045 #include <uves_error.h>
00046 #include <uves_dump.h>
00047 #include <cpl.h>
00048 
00049 #include <irplib_utils.h>
00050 #include <stdbool.h>
00051 #include <assert.h>
00052 
00053 /*-----------------------------------------------------------------------------
00054                                    Local functions
00055  -----------------------------------------------------------------------------*/
00056 
00057 static int get_candidate(const double *a, const int ia[],
00058              int M, int N, int D,
00059              double lambda,
00060              int    (*f)(const double x[], const double a[], 
00061                      double *result),
00062              int (*dfda)(const double x[], const double a[], 
00063                      double result[]),
00064              const double *x,
00065              const double *y,
00066              const double *sigma,
00067              double *partials,
00068              cpl_matrix *alpha,
00069              cpl_matrix *beta,
00070              double *a_da);
00071 
00072 static double get_chisq(int N, int D,
00073             int (*f)(const double x[], const double a[], 
00074                  double *result),
00075             const double *a,
00076             const double *x,
00077             const double *y,
00078             const double *sigma);
00079 
00080 
00081 static cpl_image * 
00082 uves_image_filter_wrapper(const cpl_image *b, 
00083                           const cpl_matrix *k, 
00084                           cpl_filter_mode mode);
00085 cpl_image * 
00086 uves_image_filter_median(const cpl_image * img, const cpl_matrix * mx)
00087 {
00088     return uves_image_filter_wrapper(img, mx, CPL_FILTER_MEDIAN);
00089 }
00090 
00091 cpl_image * 
00092 uves_image_filter_linear(const cpl_image *img, const cpl_matrix * mx)
00093 {
00094     return uves_image_filter_wrapper(img, mx, CPL_FILTER_LINEAR);
00095 
00096 }
00097 
00098 /*----------------------------------------------------------------------------
00099                                    Defines
00100  ---------------------------------------------------------------------------*/
00101 
00103 #define image_is_rejected(badmap, x, y) \
00104   ((badmap) != NULL && (badmap)[((x)-1) + ((y)-1)*nx] == CPL_BINARY_1)
00105 
00106 #ifndef UVES_FIT_MAXITER
00107 #define UVES_FIT_MAXITER 1000
00108 #endif
00109 
00110 /*-----------------------------------------------------------------------------
00111                                    Types
00112  -----------------------------------------------------------------------------*/
00113 /* @cond */
00114 typedef struct {
00115     double x;
00116     double y;
00117 } uves_fit_1d_input;
00118 /* @endcond */
00119 
00120 /*-----------------------------------------------------------------------------
00121                                    Implementation
00122  -----------------------------------------------------------------------------*/
00127 static cpl_image * 
00128 uves_image_filter_wrapper(const cpl_image *b, 
00129                           const cpl_matrix *k, 
00130                           cpl_filter_mode mode)
00131 {
00132     const double EPSILON = 1E-5;
00133     int nx   = cpl_image_get_size_x(b);
00134     int ny   = cpl_image_get_size_y(b);
00135     int nrow = cpl_matrix_get_nrow(k);
00136     int ncol = cpl_matrix_get_ncol(k);
00137     int i, j;
00138     cpl_type type = cpl_image_get_type(b);
00139     cpl_image * a = cpl_image_new(nx, ny, type);
00140     // where m is a cpl_mask with a CPL_BINARY_1 whereever k has a 1.0.
00141     cpl_mask* m = cpl_mask_new(ncol, nrow);
00142     cpl_msg_warning(cpl_func, "nx[%d], ny[%d], ncol[%d], nrow[%d]", nx, ny, ncol, nrow);
00143     for (i = 0; i < ncol ; i++)
00144     {
00145         for (j = 0; j < nrow ; j++)
00146         {
00147             double value = cpl_matrix_get(k, j, i);
00148             if (fabs(value - 1.0) < EPSILON)
00149             {
00150                 cpl_mask_set(m, i + 1, j + 1, CPL_BINARY_1);
00151             }
00152         }
00153     }
00154 
00155     cpl_image_filter_mask(a, b, m, mode, CPL_BORDER_FILTER);
00156     cpl_mask_delete(m);
00157     return a;
00158 }
00159 
00160 cpl_image*
00161 uves_image_filter_mode(const cpl_image* b,
00162                       const cpl_matrix * ker,
00163                       cpl_filter_mode filter)
00164 {
00165   int nx   = cpl_image_get_size_x(b);
00166   int ny   = cpl_image_get_size_y(b);
00167   int type = cpl_image_get_type(b);
00168   cpl_image * a = cpl_image_new(nx, ny, type);
00169 
00170   switch(filter) {
00171   case CPL_FILTER_MEDIAN:
00172     check_nomsg(cpl_image_filter(a, b, ker, CPL_FILTER_MEDIAN, CPL_BORDER_FILTER));
00173     break;
00174   case CPL_FILTER_LINEAR:
00175     check_nomsg(cpl_image_filter(a, b, ker, CPL_FILTER_LINEAR, CPL_BORDER_FILTER));
00176     break;
00177   case CPL_FILTER_STDEV:
00178     cpl_image_filter(a, b, ker, CPL_FILTER_STDEV, CPL_BORDER_FILTER);
00179     break;
00180   case CPL_FILTER_MORPHO:
00181     cpl_image_filter(a, b, ker, CPL_FILTER_MORPHO, CPL_BORDER_FILTER);
00182     break;
00183   default:
00184     uves_msg_error("Filter type not supported");
00185     return NULL;
00186   }
00187  cleanup:
00188 
00189   return a;
00190 
00191 }
00192 
00193 
00194 /*----------------------------------------------------------------------------*/
00201 /*----------------------------------------------------------------------------*/
00202 void uves_image_reject_all(cpl_image *image)
00203 {
00204     int nx, ny;
00205     int x, y;
00206 
00207     assure_nomsg( image != NULL, CPL_ERROR_NULL_INPUT );
00208     nx = cpl_image_get_size_x(image);
00209     ny = cpl_image_get_size_y(image);
00210 
00211     for (y = 1; y <= ny; y++) {
00212         for (x = 1; x <= nx; x++) {
00213             cpl_image_reject(image, x, y);
00214         }
00215     }
00216     
00217   cleanup:
00218     return;    
00219 }
00220 
00221 
00222 /*----------------------------------------------------------------------------*/
00281 /*----------------------------------------------------------------------------*/
00282 
00283 cpl_error_code
00284 uves_fit_1d_image(const cpl_image *image, const cpl_image *noise,
00285           const cpl_binary *image_badmap,
00286           bool horizontal, bool fix_back, bool fit_back,
00287           int xlo, int xhi, int y_0,
00288           double *x0, double *sigma, double *norm, double *background,
00289                   double *slope,
00290           double *mse, double *red_chisq,
00291           cpl_matrix **covariance,
00292           int (*f)   (const double x[], const double a[], double *result),
00293           int (*dfda)(const double x[], const double a[], double result[]),
00294           int M)
00295 {
00296     cpl_vector *x = NULL;
00297     cpl_vector *y = NULL;
00298     cpl_vector *sigma_y = NULL;     /* Noise vector          */
00299     
00300     cpl_fit_mode fit_pattern;
00301     int nx, ny;                     /* Image size            */
00302     int N;                          /* Number of good pixels */
00303     int i, j;
00304     cpl_type image_type;
00305 
00306     const double *image_data       = NULL;        /* Pointer to data */
00307     const double *noise_data       = NULL;        /* Pointer to data */
00308 
00309     assure( x0 != NULL        , CPL_ERROR_NULL_INPUT, "Null fit parameter");
00310     assure( sigma != NULL     , CPL_ERROR_NULL_INPUT, "Null fit parameter");
00311     assure( norm != NULL      , CPL_ERROR_NULL_INPUT, "Null fit parameter");
00312     assure( background != NULL, CPL_ERROR_NULL_INPUT, "Null fit parameter");
00313     /* mse, red_chisq, covariance may be NULL */
00314 
00315     assure( image != NULL, CPL_ERROR_NULL_INPUT, "Null image");
00316    
00317     image_type = cpl_image_get_type(image);
00318 
00319     /* To support the following types, use cpl_image_get() or
00320        more multiple type pointers to data buffer.
00321        cpl_ensure_code( image_type == CPL_TYPE_INT ||
00322        image_type == CPL_TYPE_FLOAT ||
00323        image_type == CPL_TYPE_DOUBLE,
00324        CPL_ERROR_TYPE_MISMATCH);
00325     */
00326 
00327     assure( image_type == CPL_TYPE_DOUBLE, CPL_ERROR_UNSUPPORTED_MODE,
00328         "Unsupported type: %s", uves_tostring_cpl_type(image_type));
00329 
00330     image_data = cpl_image_get_data_double_const(image);
00331 
00332     if (noise != NULL)
00333     {
00334         image_type = cpl_image_get_type(noise);
00335         /*  See comment above.
00336           cpl_ensure_code( image_type == CPL_TYPE_INT ||
00337           image_type == CPL_TYPE_FLOAT ||
00338           image_type == CPL_TYPE_DOUBLE,
00339           CPL_ERROR_TYPE_MISMATCH);
00340         */
00341         assure( image_type == CPL_TYPE_DOUBLE, CPL_ERROR_UNSUPPORTED_MODE, 
00342             "Unsupported type: %s", uves_tostring_cpl_type(image_type));
00343 
00344         noise_data = cpl_image_get_data_double_const(noise);
00345     }   
00346 
00347     nx = cpl_image_get_size_x(image);
00348     ny = cpl_image_get_size_y(image);
00349 
00350     if (horizontal)
00351     {
00352         assure( 1 <= xlo && xlo <= xhi && xhi <= nx &&
00353             1 <= y_0  && y_0  <= ny, CPL_ERROR_ACCESS_OUT_OF_RANGE,
00354             "Illegal window (%d, %d)-(%d, %d), image: %dx%d",
00355             xlo, y_0, xhi, y_0,
00356             nx, ny);
00357     }
00358     else
00359     {
00360         assure( 1 <= xlo && xlo <= xhi && xhi <= ny &&
00361             1 <= y_0  && y_0  <= nx,
00362             CPL_ERROR_ACCESS_OUT_OF_RANGE,
00363             "Illegal window (%d, %d)-(%d, %d), image: %dx%d",
00364             y_0, xlo, y_0, xhi,
00365             nx, ny);
00366     }
00367     
00368     /* Noise image must have same size
00369      * as the input image
00370      */
00371     if (noise != NULL)
00372     {
00373         assure( cpl_image_get_size_x(noise) == nx &&
00374             cpl_image_get_size_y(noise) == ny,
00375             CPL_ERROR_INCOMPATIBLE_INPUT, "Noise image: %dx%d, image: %dx%d:",
00376             cpl_image_get_size_x(noise),
00377             cpl_image_get_size_y(noise),
00378             nx, ny);
00379     }
00380    
00381     /* Count good pixels in sub-window, check that noise image is positive */
00382     N = 0;
00383     for (i = xlo; i <= xhi; i++)
00384     {
00385         if ( !image_is_rejected(image_badmap,
00386                     (horizontal) ? i : y_0,
00387                     (horizontal) ? y_0 : i) )
00388         {
00389             if ( noise != NULL)
00390             {
00391                 if( !image_is_rejected(image_badmap,
00392                            (horizontal) ? i : y_0,
00393                            (horizontal) ? y_0 : i) )
00394                 {
00395                     /* Noise image must be positive (only check
00396                        pixels that are actually used) */
00397                     assure( /*cpl_image_get(noise,
00398                           (horizontal) ? i : y_0,
00399                           (horizontal) ? y_0 : i,
00400                           &pis_rejected)*/
00401                     noise_data[(((horizontal) ? i : y_0) - 1) +
00402                            (((horizontal) ? y_0 : i) - 1) * nx]
00403                     > 0,
00404                     CPL_ERROR_ILLEGAL_INPUT,
00405                     "Non-positive noise at (%d, %d): %e",
00406                     (horizontal) ? i : y_0,
00407                     (horizontal) ? y_0 : i,
00408                     noise_data[(((horizontal) ? i : y_0) - 1) +
00409                            (((horizontal) ? y_0 : i) - 1) * nx]
00410                     );
00411                     
00412                     N += 1;
00413                 }
00414                 else
00415                 {
00416                     /* Pixel value is good, but noise pixel is
00417                        bad. Don't use. */
00418                 }
00419             }
00420             else
00421             {
00422                 /* Pixel is good. No noise image */
00423                 N += 1;
00424             }
00425         }
00426     }
00427    
00428     /* Check that there is at least one good pixel */
00429     assure( N >= 1, CPL_ERROR_ILLEGAL_INPUT, "Only %d good pixel(s)", N);
00430 
00431     /* Allocate space */
00432     x = cpl_vector_new(N);
00433     y = cpl_vector_new(N);
00434     if (noise != NULL)
00435     {
00436         sigma_y = cpl_vector_new(N);
00437         assure_mem( sigma_y );
00438     }
00439 
00440     if (fix_back)
00441     {
00442         fit_pattern = CPL_FIT_CENTROID | CPL_FIT_STDEV | CPL_FIT_AREA;
00443     }
00444     else if (fit_back)
00445     {
00446         fit_pattern = CPL_FIT_AREA | CPL_FIT_OFFSET;
00447     }
00448     else
00449     {
00450         fit_pattern = CPL_FIT_ALL;
00451     }
00452    
00453     assure_mem( x );
00454     assure_mem( y );
00455     
00456     /* Copy the N good pixels from the input image to vectors,
00457        j count good pixels */
00458     for (i = xlo, j = 0;
00459      i <= xhi;
00460      i++)
00461     {
00462         double flux;
00463         
00464         /*
00465           flux = cpl_image_get(image,
00466           (horizontal) ? xlo+i : y_0,
00467           (horizontal) ? y_0 : xlo+i,
00468           &pis_rejected);
00469         */
00470         
00471         flux = image_data[(((horizontal) ? i : y_0) - 1) +
00472                   (((horizontal) ? y_0 : i) - 1) * nx];
00473        
00474         //if (!pis_rejected)
00475         if ( !image_is_rejected(image_badmap,
00476                     (horizontal) ? i : y_0,
00477                     (horizontal) ? y_0 : i) )
00478         {
00479             if (noise != NULL)
00480             {
00481                 double flux_noise;
00482                
00483                 /* flux_noise = cpl_image_get(noise,
00484                    (horizontal) ? xlo+i : y_0,
00485                    (horizontal) ? y_0 : xlo+i,
00486                    &pis_rejected);
00487                 */
00488                
00489                 flux_noise = noise_data
00490                 [(((horizontal) ? i : y_0) - 1) +
00491                  (((horizontal) ? y_0 : i) - 1) * nx];
00492                
00493                 //if (!pis_rejected)
00494                 if ( !image_is_rejected(image_badmap,
00495                             (horizontal) ?
00496                             i : y_0,
00497                             (horizontal)
00498                             ? y_0 : i) )
00499                 {
00500                     cpl_vector_set(x,       j, i);
00501                     cpl_vector_set(y,       j, flux);
00502                     cpl_vector_set(sigma_y, j, flux_noise);
00503                     j++;
00504                 }
00505             }
00506             else
00507             {
00508                 cpl_vector_set(x, j, i);
00509                 cpl_vector_set(y, j, flux);
00510                 j++;
00511             }
00512         }
00513     }
00514     passure( j == N, "%d %d", j, N);
00515     
00516     check( uves_fit_1d(x, NULL,      /* x, sigma_x */
00517                y, sigma_y,
00518                fit_pattern, fit_back,
00519                x0, sigma, norm, background,
00520                        slope,
00521                mse, red_chisq,
00522                covariance,
00523                f, dfda, M),
00524        "Fit failed");
00525 
00526     
00527   cleanup:
00528     uves_free_vector(&x);
00529     uves_free_vector(&y);
00530     uves_free_vector(&sigma_y);
00531 
00532     return cpl_error_get_code();
00533 }
00534 
00535 /*----------------------------------------------------------------------------*/
00543 /*----------------------------------------------------------------------------*/
00544 static int uves_fit_1d_compare(const void *left,
00545                    const void *right)
00546 {
00547     return 
00548     (((uves_fit_1d_input *)left )->x <  
00549      ((uves_fit_1d_input *)right)->x) ? -1 :
00550     (((uves_fit_1d_input *)left )->x == 
00551      ((uves_fit_1d_input *)right)->x) ? 0  : 1;
00552 }
00553 
00554 /*----------------------------------------------------------------------------*/
00577 /*----------------------------------------------------------------------------*/
00578 cpl_error_code
00579 uves_fit_1d(cpl_vector *x, const cpl_vector *sigma_x,
00580             cpl_vector *y, const cpl_vector *sigma_y,
00581             cpl_fit_mode fit_pars, bool fit_back,
00582             double *x0, double *sigma, double *area, double *offset, double *slope,
00583             double *mse, double *red_chisq,
00584             cpl_matrix **covariance,
00585             int (*f)   (const double x[], const double a[], double *result),
00586             int (*dfda)(const double x[], const double a[], double result[]),
00587             int M)
00588 {
00589     cpl_matrix *x_matrix = NULL; /* LM algorithm needs a matrix,
00590                       not a vector                 */
00591 
00592     int N;                          /* Number of data points        */
00593     double xlo, xhi;                /* Min/max x                    */
00594 
00595     /* Initial parameter values */
00596     double x0_guess    = 0;  /* Avoid warnings about uninitialized variables */
00597     double sigma_guess = 0;
00598     double area_guess;
00599     double offset_guess;
00600 
00601     cpl_ensure_code( M == 4 || M == 5, CPL_ERROR_UNSUPPORTED_MODE);
00602 
00603     /* Validate input */
00604     cpl_ensure_code( x       != NULL, CPL_ERROR_NULL_INPUT);
00605     cpl_ensure_code( sigma_x == NULL, CPL_ERROR_UNSUPPORTED_MODE);
00606     cpl_ensure_code( y       != NULL, CPL_ERROR_NULL_INPUT);
00607     /* sigma_y may be NULL or non-NULL */
00608     
00609     N = cpl_vector_get_size(x);
00610 
00611     cpl_ensure_code( N == cpl_vector_get_size(y),
00612              CPL_ERROR_INCOMPATIBLE_INPUT);
00613 
00614     if (sigma_x != NULL)
00615     {
00616         cpl_ensure_code( N == cpl_vector_get_size(sigma_x),
00617                  CPL_ERROR_INCOMPATIBLE_INPUT);
00618     }
00619     if (sigma_y != NULL)
00620     {
00621         cpl_ensure_code( N == cpl_vector_get_size(sigma_y),
00622                  CPL_ERROR_INCOMPATIBLE_INPUT);
00623     }
00624 
00625     cpl_ensure_code( x0     != NULL &&
00626              sigma  != NULL &&
00627              area   != NULL &&
00628              offset != NULL &&
00629                      (M != 5 || slope != NULL), CPL_ERROR_NULL_INPUT);
00630 
00631     if (! (fit_pars & CPL_FIT_STDEV))
00632     {
00633         cpl_ensure_code( *sigma > 0, CPL_ERROR_ILLEGAL_INPUT);
00634     }
00635 
00636     cpl_ensure_code( !fit_back || fit_pars == (CPL_FIT_OFFSET | CPL_FIT_AREA),
00637              CPL_ERROR_INCOMPATIBLE_INPUT);
00638 
00639     /* Input area must be positive if fit_back is false */
00640     if (! (fit_pars & CPL_FIT_AREA) && !fit_back)
00641     {
00642         cpl_ensure_code( *area > 0, CPL_ERROR_ILLEGAL_INPUT);
00643     }
00644 
00645     /* mse, red_chisq may be NULL */
00646 
00647     /* Need more than number_of_parameters points to calculate chi^2.
00648      * There are less than 5 parameters. */
00649     cpl_ensure_code( red_chisq == NULL || N >= 5, CPL_ERROR_ILLEGAL_INPUT);
00650     
00651     if (covariance != NULL) *covariance = NULL;
00652     /* If covariance computation is requested, then
00653      * return either the covariance matrix or NULL
00654      * (don't return undefined pointer).
00655      */
00656     
00657     /* Cannot compute chi square & covariance without sigma_y */
00658     cpl_ensure_code( (red_chisq == NULL && covariance == NULL) || 
00659              sigma_y != NULL,
00660              CPL_ERROR_INCOMPATIBLE_INPUT);
00661     
00662     /* Create matrix from x-data */
00663     x_matrix = cpl_matrix_wrap(N, 1, cpl_vector_get_data(x));
00664     if (x_matrix == NULL)
00665     {
00666         cpl_ensure_code(
00667                 CPL_FALSE,
00668         CPL_ERROR_ILLEGAL_OUTPUT);
00669     }
00670     
00671     /* Check that any provided sigmas are positive. */
00672     if (sigma_x != NULL && cpl_vector_get_min(sigma_x) <= 0)
00673     {
00674         cpl_matrix_unwrap(x_matrix);
00675         cpl_ensure_code(
00676         CPL_FALSE,
00677         CPL_ERROR_ILLEGAL_INPUT);
00678     }
00679     if (sigma_y != NULL && cpl_vector_get_min(sigma_y) <= 0)
00680     {
00681         cpl_matrix_unwrap(x_matrix);
00682         cpl_ensure_code(
00683         CPL_FALSE,
00684         CPL_ERROR_ILLEGAL_INPUT);
00685     }
00686 
00687     /* Compute guess parameters using robust estimation
00688      * (This step might be improved by taking into account the 
00689      * uncertainties but for simplicity's sake do not)
00690      */
00691     if (fit_back)
00692     {
00693         /* We need to estimate only these two parameters */
00694         assert( fit_pars == CPL_FIT_OFFSET || CPL_FIT_AREA);
00695 
00696 #if defined CPL_VERSION_CODE && CPL_VERSION_CODE >= CPL_VERSION(4, 0, 0)
00697         offset_guess = cpl_vector_get_median_const(y);
00698 #else
00699         offset_guess = cpl_vector_get_median(y);
00700 #endif
00701         area_guess = N*(cpl_vector_get_mean(y) - offset_guess);
00702         /* Sum of (flux-offset) */
00703 
00704         x0_guess = *x0;
00705         sigma_guess = *sigma;
00706     }
00707     else {
00708     double sum  = 0.0;
00709     double quartile[3];
00710     double fraction[3] = {0.25 , 0.50 , 0.75};
00711     const double *y_data = cpl_vector_get_data_const(y);
00712 
00713     if (fit_pars & CPL_FIT_OFFSET)
00714         {
00715         /* Estimate offset as 25% percentile of y-values.
00716            (The minimum value may be too low for noisy input,
00717            the median is too high if there is not much
00718            background in the supplied data, so use
00719            something inbetween).
00720         */
00721 
00722         cpl_vector *y_dup = cpl_vector_duplicate(y);
00723         
00724         if (y_dup == NULL)
00725             {
00726             cpl_matrix_unwrap(x_matrix);
00727             cpl_ensure_code(
00728                 CPL_FALSE,
00729                 CPL_ERROR_ILLEGAL_OUTPUT);
00730             }
00731         
00732         offset_guess = uves_utils_get_kth_double(
00733             cpl_vector_get_data(y_dup), N, N/4);
00734         
00735         cpl_vector_delete(y_dup);
00736         }
00737     else
00738         {
00739         offset_guess = *offset;
00740         }
00741     
00742     /* Get quartiles of distribution
00743        (only bother if it's needed for estimation of x0 or sigma) */
00744     if ( (fit_pars & CPL_FIT_CENTROID) ||
00745          (fit_pars & CPL_FIT_STDEV   )
00746         )
00747         {
00748         /* The algorithm requires the input to be sorted
00749            as function of x, so do that (using qsort), and
00750            work on the sorted copy. Of course, the y-vector
00751            must be re-ordered along with x.
00752            sigma_x and sigma_y are not used, so don't copy those.
00753         */
00754         
00755         uves_fit_1d_input
00756             *sorted_input = cpl_malloc(N * sizeof(*sorted_input));
00757         const double *x_data = cpl_matrix_get_data_const(x_matrix);
00758         cpl_boolean is_sorted = CPL_TRUE;
00759         int i;
00760         
00761         if (sorted_input == NULL)
00762             {
00763             cpl_matrix_unwrap(x_matrix);
00764             cpl_ensure_code(
00765                 CPL_FALSE,
00766                 CPL_ERROR_ILLEGAL_INPUT);
00767             }
00768         
00769         for (i = 0; i < N; i++)
00770             {
00771             sorted_input[i].x = x_data[i];
00772             sorted_input[i].y = y_data[i];
00773             
00774             is_sorted = is_sorted && 
00775                 (i==0 || (x_data[i-1] < x_data[i]));
00776             }            
00777         
00778         if (!is_sorted)
00779             {
00780             qsort(sorted_input, N, sizeof(*sorted_input),
00781                   &uves_fit_1d_compare);
00782             }
00783 
00784         for(i = 0; i < N; i++)
00785             {
00786             double flux = sorted_input[i].y;
00787             
00788             sum += (flux - offset_guess);
00789             }
00790         /* Note that 'sum' must be calculated the same way as
00791            'running_sum' below, Otherwise (due to round-off error)
00792            'running_sum' might end up being different from 'sum'(!).
00793            Specifically, it will not work to calculate 'sum' as
00794            
00795            (flux1 + ... + fluxN)  -  N*offset_guess
00796         */
00797         
00798         if (sum > 0.0)
00799             {
00800             double flux, x1, x2;
00801             double running_sum = 0.0;
00802             int j;
00803             
00804             i = 0;
00805             flux = sorted_input[i].y - offset_guess;
00806             
00807             for (j = 0; j < 3; j++)
00808                 {
00809                 double limit = fraction[j] * sum;
00810                 double k;
00811                 
00812                 while (running_sum + flux < limit && i < N-1)
00813                     {
00814                     running_sum += flux;
00815                     i++;
00816                     flux = sorted_input[i].y - offset_guess;
00817                     }
00818 
00819                 /* Fraction [0;1] of current flux needed
00820                    to reach the quartile */
00821                 k = (limit - running_sum)/flux;
00822                 
00823                 if (k <= 0.5)
00824                     {
00825                     /* Interpolate linearly between
00826                        current and previous position
00827                     */
00828                     if (0 < i)
00829                         {
00830                         x1 = sorted_input[i-1].x;
00831                         x2 = sorted_input[i].x;
00832                         
00833                         quartile[j] = 
00834                             x1*(0.5-k) +
00835                             x2*(0.5+k);
00836                         /*
00837                           k=0   => quartile = midpoint,
00838                           k=0.5 => quartile = x2
00839                         */
00840                         }
00841                     else
00842                         {
00843                         quartile[j] = sorted_input[i].x;
00844                         }
00845                     }
00846                 else
00847                     {
00848                     /* Interpolate linearly between
00849                        current and next position */
00850                     if (i < N-1)
00851                         {
00852                         x1 = sorted_input[i].x;
00853                         x2 = sorted_input[i+1].x;
00854                         
00855                         quartile[j] = 
00856                             x1*( 1.5-k) +
00857                             x2*(-0.5+k);
00858                         /*
00859                           k=0.5 => quartile = x1,
00860                           k=1.0 => quartile = midpoint
00861                         */
00862                         }
00863                     else
00864                         {
00865                         quartile[j] = sorted_input[i].x;
00866                         }
00867                     }
00868                 }
00869             }
00870         else
00871             {
00872             /* If there's no flux (sum = 0) then
00873                set quartiles to something that's not 
00874                completely insensible.
00875             */
00876             quartile[1] = cpl_matrix_get_mean(x_matrix);
00877             
00878             quartile[2] = quartile[1];
00879             quartile[0] = quartile[2] - 1.0;
00880             /* Then sigma_guess = 1.0 */
00881             }
00882 
00883         cpl_free(sorted_input);
00884         } /* If need to compute quartiles */
00885         
00886     /* x0_guess = median of distribution */
00887     x0_guess = (fit_pars & CPL_FIT_CENTROID) ? quartile[1] : *x0;
00888     
00889     /* sigma_guess = median of absolute residuals
00890      *
00891      *  (68% is within 1 sigma, and 50% is within 0.6744
00892      *  sigma,  => quartile3-quartile1 = 2*0.6744 sigma)
00893      */
00894     sigma_guess = (fit_pars & CPL_FIT_STDEV) ? 
00895         (quartile[2] - quartile[0]) / (2*0.6744) : *sigma;
00896     
00897     area_guess  = (fit_pars & CPL_FIT_AREA) ?
00898         (cpl_vector_get_max(y) - offset_guess) * sqrt(2*M_PI) * sigma_guess : *area;
00899     /* This formula makes sense only if the area is positive */
00900     
00901     /* Make sure that the area is a positive number */
00902     if ( area_guess <= 0)  area_guess = 1.0;
00903     if (sigma_guess <= 0) sigma_guess = 1.0;
00904     }
00905     
00906     /* Wrap parameters, fit, unwrap */
00907     {
00908     cpl_vector *a;
00909     int ia[5];            /* The last element
00910                  is ignored if
00911                  M = 4 */
00912     cpl_error_code ec;
00913 
00914     assert(M == 4 || M == 5);
00915     a = cpl_vector_new(M);
00916 
00917     if (a == NULL)
00918         {
00919         cpl_matrix_unwrap(x_matrix);
00920         cpl_ensure_code(
00921             CPL_FALSE,
00922             CPL_ERROR_ILLEGAL_OUTPUT);
00923         }
00924 
00925     cpl_vector_set(a, 0, x0_guess);
00926     cpl_vector_set(a, 1, sigma_guess);
00927     cpl_vector_set(a, 2, area_guess);
00928     cpl_vector_set(a, 3, offset_guess);
00929     
00930     ia[0] = fit_pars & CPL_FIT_CENTROID;
00931     ia[1] = fit_pars & CPL_FIT_STDEV;
00932     ia[2] = fit_pars & CPL_FIT_AREA;
00933     ia[3] = fit_pars & CPL_FIT_OFFSET;
00934 
00935     if (M == 5)
00936         {
00937         /* linear sky-term, first hold it constant,
00938          * then call LM-fitting a second time where
00939          * it is non-constant */
00940         if (fit_pars & CPL_FIT_OFFSET)
00941                     {
00942                         cpl_vector_set(a, 4, 0);
00943                     }
00944                 else
00945                     {
00946                         cpl_vector_set(a, 4, *slope);
00947                     }
00948 
00949         ia[4] = 0;
00950         }
00951     
00952     ec = uves_fit(x_matrix, NULL,
00953             y, sigma_y, 
00954             a, ia, f, dfda,
00955             mse, red_chisq,
00956             covariance);
00957     
00958 //    printf("Sky: e='%s'\n", cpl_error_get_message()); cpl_vector_dump(a, stdout);
00959     if (M == 5)
00960         {
00961         ia[4] = fit_pars & CPL_FIT_OFFSET;
00962 
00963         if (covariance != NULL)
00964             {
00965             uves_free_matrix(covariance);
00966             }
00967 
00968         ec = uves_fit(x_matrix, NULL,
00969                 y, sigma_y, 
00970                 a, ia, f, dfda,
00971                 mse, red_chisq,
00972                 covariance);
00973 //    printf("Sky: e='%s'\n", cpl_error_get_message()); cpl_vector_dump(a, stdout);
00974         }
00975 
00976     cpl_matrix_unwrap(x_matrix);
00977     
00978     /* Check return status of fitting routine */
00979     if (ec == CPL_ERROR_NONE      ||
00980         ec == CPL_ERROR_SINGULAR_MATRIX)
00981         {
00982         /* The LM algorithm converged. The computation
00983          *  of the covariance matrix might have failed.
00984          */
00985         
00986         /* In principle, the LM algorithm might have converged
00987          * to a negative sigma (even if the guess value was
00988          * positive). Make sure that the returned sigma is positive
00989          * (by convention).
00990          */
00991 
00992         if (CPL_FIT_CENTROID) *x0     = cpl_vector_get(a, 0);
00993         if (CPL_FIT_STDEV   ) *sigma  = fabs(cpl_vector_get(a, 1));
00994         if (CPL_FIT_AREA    ) *area   = cpl_vector_get(a, 2);
00995         if (CPL_FIT_OFFSET  ) {
00996                     *offset = cpl_vector_get(a, 3);
00997                     if (M == 5) *slope = cpl_vector_get(a, 4);
00998                 }
00999         }
01000     
01001     cpl_vector_delete(a);
01002     
01003     xlo = cpl_vector_get_min(x);
01004     xhi = cpl_vector_get_max(x);
01005 
01006     if (ec == CPL_ERROR_CONTINUE ||
01007         !(
01008         !irplib_isnan(*x0    ) && !irplib_isinf(*x0    ) &&
01009         !irplib_isnan(*sigma ) && !irplib_isinf(*sigma ) &&
01010         !irplib_isnan(*area  ) && !irplib_isinf(*area  ) &&
01011         !irplib_isnan(*offset) && !irplib_isinf(*offset) &&
01012         ((M != 5) || (!irplib_isnan(*slope ) && !irplib_isinf(*slope ))) &&
01013         xlo <= *x0 && *x0 <= xhi &&
01014         0 < *sigma && *sigma < (xhi - xlo + 1) &&
01015         (fit_back || 0 < *area)
01016         /* This extra check on the background level makes sense
01017            iff the input flux is assumed to be positive
01018            && *offset > - *area  */
01019         )
01020         )
01021         {
01022         /* The LM algorithm did not converge, or it converged to
01023          * a non-sensical result. Return the guess parameter values
01024          * in order to enable the caller to recover. */
01025 
01026         *x0         = x0_guess;
01027         *sigma      = sigma_guess;
01028         *area       = area_guess;
01029         *offset     = offset_guess;
01030                 if (M == 5) *slope = 0;
01031         
01032         /* In this case the covariance matrix will not make sense
01033            (because the LM algorithm failed), so delete it */
01034         if (covariance != NULL && *covariance != NULL)
01035             {
01036             cpl_matrix_delete(*covariance);
01037             *covariance = NULL;
01038             }
01039 
01040         /* Return CPL_ERROR_CONTINUE if the fitting routine failed */
01041         cpl_ensure_code(
01042             CPL_FALSE,
01043             CPL_ERROR_CONTINUE);
01044         }
01045     }
01046     
01047     return CPL_ERROR_NONE;
01048 }
01049 /* @endcond */
01050 
01051 #define DEBUG_LM 0   /* Set to non-zero to print info on the error msg level */
01052 /*----------------------------------------------------------------------------*/
01119 /*----------------------------------------------------------------------------*/
01120 cpl_error_code
01121 uves_fit(const cpl_matrix *x, const cpl_matrix *sigma_x,
01122      const cpl_vector *y, const cpl_vector *sigma_y,
01123      cpl_vector *a, const int ia[],
01124      int    (*f)(const double x[], const double a[], double *result),
01125      int (*dfda)(const double x[], const double a[], double result[]),
01126      double *mse,
01127      double *red_chisq,
01128      cpl_matrix **covariance)
01129 {
01130     const double *x_data     = NULL; /* Pointer to input data                  */
01131     const double *y_data     = NULL; /* Pointer to input data                  */
01132     const double *sigma_data = NULL; /* Pointer to input data                  */
01133     int N    = 0;                    /* Number of data points                  */
01134     int D    = 0;                    /* Dimension of x-points                  */
01135     int M    = 0;                    /* Number of fit parameters               */
01136     int Mfit = 0;                    /* Number of non-constant fit
01137                         parameters                             */
01138 
01139     double lambda    = 0.0;          /* Lambda in L-M algorithm                */
01140     double MAXLAMBDA = 10e40;        /* Parameter to control the graceful exit
01141                     if steepest descent unexpectedly fails */
01142     double chi_sq    = 0.0;          /* Current  chi^2                         */
01143     int count        = 0;            /* Number of successive small improvements
01144                     in chi^2 */
01145     int iterations   = 0;
01146    
01147     cpl_matrix *alpha  = NULL;       /* The MxM ~curvature matrix used in L-M  */
01148     cpl_matrix *beta   = NULL;       /* Mx1 matrix = -.5 grad(chi^2)           */
01149     double *a_data     = NULL;       /* Parameters, a                          */
01150     double *a_da       = NULL;       /* Candidate position a+da                */
01151     double *part       = NULL;       /* The partial derivatives df/da          */
01152     int *ia_local      = NULL;       /* non-NULL version of ia                 */
01153    
01154     /* If covariance computation is requested, then either
01155      * return the covariance matrix or return NULL.
01156      */
01157     if (covariance != NULL) *covariance = NULL;
01158 
01159     /* Validate input */
01160     cpl_ensure_code(x       != NULL, CPL_ERROR_NULL_INPUT);
01161     cpl_ensure_code(sigma_x == NULL, CPL_ERROR_UNSUPPORTED_MODE);
01162     cpl_ensure_code(y       != NULL, CPL_ERROR_NULL_INPUT);
01163     cpl_ensure_code(a       != NULL, CPL_ERROR_NULL_INPUT);
01164     /* ia may be NULL */
01165     cpl_ensure_code(f       != NULL, CPL_ERROR_NULL_INPUT);
01166     cpl_ensure_code(dfda    != NULL, CPL_ERROR_NULL_INPUT);
01167 
01168     /* Chi^2 and covariance computations require sigmas to be known */
01169     cpl_ensure_code( sigma_y != NULL || (red_chisq == NULL && covariance == NULL),
01170              CPL_ERROR_INCOMPATIBLE_INPUT);
01171 
01172     D = cpl_matrix_get_ncol(x);
01173     N = cpl_matrix_get_nrow(x);
01174     M = cpl_vector_get_size(a);
01175     cpl_ensure_code(N > 0 && D > 0 && M > 0, CPL_ERROR_ILLEGAL_INPUT);
01176 
01177     cpl_ensure_code( cpl_vector_get_size(y) == N,
01178              CPL_ERROR_INCOMPATIBLE_INPUT);
01179 
01180     x_data = cpl_matrix_get_data_const(x);
01181     y_data = cpl_vector_get_data_const(y);
01182     a_data = cpl_vector_get_data(a);
01183 
01184     if (sigma_y != NULL)
01185     {
01186         cpl_ensure_code( cpl_vector_get_size(sigma_y) == N,
01187                  CPL_ERROR_INCOMPATIBLE_INPUT);
01188         /* Sigmas must be positive */
01189         cpl_ensure_code( cpl_vector_get_min (sigma_y) > 0,
01190                  CPL_ERROR_ILLEGAL_INPUT);
01191         sigma_data = cpl_vector_get_data_const(sigma_y);
01192     }
01193 
01194     ia_local = cpl_malloc(M * sizeof(int));
01195     cpl_ensure_code(ia_local != NULL, CPL_ERROR_ILLEGAL_OUTPUT);
01196 
01197     /* Count non-constant fit parameters, copy ia */
01198     if (ia != NULL)
01199     {
01200         int i;
01201 
01202         Mfit = 0;
01203         for (i = 0; i < M; i++)
01204         {
01205             ia_local[i] = ia[i];
01206 
01207             if (ia[i] != 0) 
01208             {
01209                 Mfit += 1;
01210             }
01211         }
01212         
01213         if (! (Mfit > 0))
01214         {
01215             cpl_free(ia_local);
01216             cpl_ensure_code( CPL_FALSE,
01217                      CPL_ERROR_ILLEGAL_INPUT);
01218         }
01219     }
01220     else
01221     {
01222         /* All parameters participate */
01223         int i;
01224         
01225         Mfit = M;
01226         
01227         for (i = 0; i < M; i++)
01228         {
01229             ia_local[i] = 1;
01230         }
01231     }
01232 
01233     /* To compute reduced chi^2, we need N > Mfit */
01234     if (! ( red_chisq == NULL || N > Mfit ) )
01235     {
01236         cpl_free(ia_local);
01237         cpl_ensure_code(
01238         CPL_FALSE,
01239         CPL_ERROR_ILLEGAL_INPUT);
01240     }
01241 
01242     /* Create alpha, beta, a_da, part  work space */
01243     alpha = cpl_matrix_new(Mfit, Mfit);
01244     if (alpha == NULL)
01245     {
01246         cpl_free(ia_local);
01247         cpl_ensure_code(
01248         CPL_FALSE,
01249         CPL_ERROR_ILLEGAL_OUTPUT);
01250     }
01251    
01252     beta = cpl_matrix_new(Mfit, 1);
01253     if (beta == NULL)
01254     {
01255         cpl_free(ia_local);
01256         cpl_matrix_delete(alpha);
01257         cpl_ensure_code(
01258         CPL_FALSE,
01259         CPL_ERROR_ILLEGAL_OUTPUT);
01260     }
01261 
01262     a_da = cpl_malloc(M * sizeof(double));
01263     if (a_da == NULL)
01264     {
01265         cpl_free(ia_local);
01266         cpl_matrix_delete(alpha);
01267         cpl_matrix_delete(beta);
01268         cpl_ensure_code(
01269         CPL_FALSE,
01270         CPL_ERROR_ILLEGAL_OUTPUT);
01271     }
01272 
01273     part = cpl_malloc(M * sizeof(double));
01274     if (part == NULL)
01275     {
01276         cpl_free(ia_local);
01277         cpl_matrix_delete(alpha);
01278         cpl_matrix_delete(beta);
01279         cpl_free(a_da);
01280         cpl_ensure_code(
01281         CPL_FALSE,
01282         CPL_ERROR_ILLEGAL_OUTPUT);
01283     }
01284 
01285     /* Initialize loop variables */
01286     lambda = 0.001;
01287     count = 0;
01288     iterations = 0;
01289     if( (chi_sq = get_chisq(N, D, f, a_data, x_data, y_data, sigma_data)) < 0)
01290     {
01291         cpl_free(ia_local);
01292         cpl_matrix_delete(alpha);
01293         cpl_matrix_delete(beta);
01294         cpl_free(a_da);
01295         cpl_free(part);
01296         cpl_ensure_code(
01297         CPL_FALSE,
01298         cpl_error_get_code());
01299     }
01300 
01301 #if DEBUG_LM    
01302      uves_msg_error("Initial chi^2 = %f", chi_sq); 
01303      {int i;
01304      for (i = 0; i < M; i++) 
01305      {
01306          uves_msg_error("Initial a[%d] = %e", i, a_data[i]); 
01307      }
01308      }
01309 #endif
01310 
01311     /* Iterate until chi^2 didn't improve substantially many (say, 5)
01312        times in a row */
01313     while (count < 5 && lambda < MAXLAMBDA && iterations < UVES_FIT_MAXITER)
01314     {
01315         /* In each iteration lambda increases, or chi^2 decreases or
01316            count increases. Because chi^2 is bounded from below
01317            (and lambda and count from above), the loop will terminate */
01318 
01319         double chi_sq_candidate = 0.0;
01320         int returncode = 0;
01321 
01322         /* Get candidate position in parameter space = a+da,
01323          * where  alpha * da = beta .
01324          * Increase lambda until alpha is non-singular
01325          */
01326        
01327         while( (returncode = get_candidate(a_data, ia_local,
01328                            M, N, D,
01329                            lambda, f, dfda,
01330                            x_data, y_data, sigma_data,
01331                            part, alpha, beta, a_da)
01332                ) != 0
01333            && cpl_error_get_code() == CPL_ERROR_SINGULAR_MATRIX
01334            && lambda < MAXLAMBDA)
01335         {
01336 #if DEBUG_LM    
01337             uves_msg_error("Singular matrix. lambda = %e", lambda);
01338 #endif
01339             cpl_error_reset();
01340             lambda *= 9.0;
01341         }
01342        
01343         /* Set error if lambda diverged */
01344         if ( !( lambda < MAXLAMBDA ) )
01345         {
01346             cpl_free(ia_local);
01347             cpl_matrix_delete(alpha);
01348             cpl_matrix_delete(beta);
01349             cpl_free(a_da);
01350             cpl_free(part);
01351             cpl_ensure_code(
01352             CPL_FALSE,
01353             CPL_ERROR_CONTINUE);
01354         }
01355        
01356         if (returncode != 0)
01357         {
01358             cpl_free(ia_local);
01359             cpl_matrix_delete(alpha);
01360             cpl_matrix_delete(beta);
01361             cpl_free(a_da);
01362             cpl_free(part);
01363             cpl_ensure_code(
01364             CPL_FALSE,
01365             cpl_error_get_code());
01366         }
01367 
01368         /* Get chi^2(a+da) */
01369         if ((chi_sq_candidate = get_chisq(N, D, f, a_da,
01370                           x_data, y_data, sigma_data)) < 0)
01371         {
01372             cpl_free(ia_local);
01373             cpl_matrix_delete(alpha);
01374             cpl_matrix_delete(beta);
01375             cpl_free(a_da);
01376             cpl_free(part);
01377             cpl_ensure_code(
01378             CPL_FALSE,
01379             cpl_error_get_code());
01380         }
01381 
01382         if (chi_sq_candidate > 1.000001 * chi_sq)
01383         {
01384             /* Move towards steepest descent */
01385 #if DEBUG_LM
01386             uves_msg_error("Chi^2 = %f  Candidate = %f  Lambda = %e",
01387                chi_sq, chi_sq_candidate, lambda); 
01388 #endif            
01389             lambda *= 9.0;
01390         }
01391         else
01392         {
01393 #if DEBUG_LM
01394             uves_msg_error("Chi^2 = %f  Candidate = %f* Lambda = %e count = %d",
01395                chi_sq, chi_sq_candidate, lambda, count);
01396 #endif
01397        
01398             /* Move towards Newton's algorithm */
01399             lambda /= 10.0;
01400 
01401             /* Count the number of successive improvements in chi^2 of
01402                less than 0.01 relative */
01403             if ( chi_sq == 0 ||
01404              (chi_sq - chi_sq_candidate)/chi_sq < .01)
01405             {
01406                 count += 1;
01407             }
01408             else
01409             {
01410                 /* Chi^2 improved by a significant amount,
01411                    reset counter */
01412                 count = 0;
01413             }
01414 
01415             if (chi_sq_candidate < chi_sq)
01416             {
01417                 /* chi^2 improved, update a and chi^2 */
01418                 int i;
01419                 for (i = 0; i < M; i++) 
01420                 {
01421                     a_data[i] = a_da[i];
01422 #if DEBUG_LM
01423                     uves_msg_error("-> a[%d] = %e", i, a_da[i]); 
01424 #endif
01425                 }
01426                 chi_sq = chi_sq_candidate;
01427             }
01428         }
01429         iterations++;
01430     }
01431 
01432     /* Set error if we didn't converge */
01433     if ( !( lambda < MAXLAMBDA && iterations < UVES_FIT_MAXITER ) )
01434     {
01435 #if DEBUG_LM
01436         uves_msg_error("Failed to converge, lambda=%f iterations=%d",
01437                lambda, iterations);
01438 #endif
01439         cpl_free(ia_local);
01440         cpl_matrix_delete(alpha);
01441         cpl_matrix_delete(beta);
01442         cpl_free(a_da);
01443         cpl_free(part);
01444         cpl_ensure_code(
01445         CPL_FALSE,
01446         CPL_ERROR_CONTINUE);
01447     }
01448 
01449     /* Compute mse if requested */
01450     if (mse != NULL)
01451     {
01452         int i;
01453 
01454         *mse = 0.0;
01455        
01456         for(i = 0; i < N; i++)
01457         {
01458             double fx_i = 0.0;
01459             double residual = 0.0;
01460            
01461             /* Evaluate f(x_i) at the best fit parameters */
01462             if( f(&(x_data[i*D]),
01463               a_data,
01464               &fx_i) != 0)
01465             {
01466                 cpl_free(ia_local);
01467                 cpl_matrix_delete(alpha);
01468                 cpl_matrix_delete(beta);
01469                 cpl_free(a_da);
01470                 cpl_free(part);
01471                 cpl_ensure_code(
01472                 CPL_FALSE,
01473                 CPL_ERROR_ILLEGAL_INPUT);
01474             }
01475 
01476             residual = y_data[i] - fx_i;
01477             *mse += residual * residual;
01478         }
01479         *mse /= N;
01480     }
01481 
01482     /* Compute reduced chi^2 if requested */
01483     if (red_chisq != NULL)
01484     {
01485         /* We already know the optimal chi^2 (and that N > Mfit)*/
01486         *red_chisq = chi_sq / (N-Mfit);
01487     }
01488 
01489     /* Compute covariance matrix if requested
01490      * cov = alpha(lambda=0)^-1              
01491      */
01492     if (covariance != NULL)
01493     {
01494         cpl_matrix *cov;
01495 
01496         if( get_candidate(a_data, ia_local, 
01497                   M, N, D, 0.0, f, dfda, 
01498                   x_data, y_data, sigma_data,
01499                   part, alpha, beta, a_da)
01500         != 0)
01501         {
01502             cpl_free(ia_local);
01503             cpl_matrix_delete(alpha);
01504             cpl_matrix_delete(beta);
01505             cpl_free(a_da);
01506             cpl_free(part);
01507             cpl_ensure_code(
01508             CPL_FALSE,
01509             cpl_error_get_code());
01510         }
01511        
01512         cov = cpl_matrix_invert_create(alpha);
01513         if (cov == NULL)
01514         {
01515             cpl_free(ia_local);
01516             cpl_matrix_delete(alpha);
01517             cpl_matrix_delete(beta);
01518             cpl_free(a_da);
01519             cpl_free(part);
01520             cpl_ensure_code(
01521             CPL_FALSE,
01522             cpl_error_get_code());
01523         }
01524        
01525         /* Make sure that variances are positive */
01526         {
01527         int i;
01528         for (i = 0; i < Mfit; i++)
01529             {
01530             if ( !(cpl_matrix_get(cov, i, i) > 0) )
01531                 {
01532                 cpl_free(ia_local);
01533                 cpl_matrix_delete(alpha);
01534                 cpl_matrix_delete(beta);
01535                 cpl_free(a_da);
01536                 cpl_free(part);
01537                 cpl_matrix_delete(cov);
01538                 *covariance = NULL;
01539                 cpl_ensure_code(
01540                     CPL_FALSE,
01541                     CPL_ERROR_SINGULAR_MATRIX);
01542                 }
01543             }
01544         }
01545 
01546         /* Expand covariance matrix from Mfit x Mfit
01547            to M x M. Set rows/columns corresponding to fixed
01548            parameters to zero */
01549 
01550         *covariance = cpl_matrix_new(M, M);
01551         if (*covariance == NULL)
01552         {
01553             cpl_free(ia_local);
01554             cpl_matrix_delete(alpha);
01555             cpl_matrix_delete(beta);
01556             cpl_free(a_da);
01557             cpl_free(part);
01558             cpl_matrix_delete(cov);
01559             cpl_ensure_code(
01560             CPL_FALSE,
01561             CPL_ERROR_ILLEGAL_OUTPUT);
01562         }
01563 
01564         {
01565         int j, jmfit;
01566 
01567         for (j = 0, jmfit = 0; j < M; j++)
01568             if (ia_local[j] != 0)
01569             {
01570                 int i, imfit;
01571 
01572                 for (i = 0, imfit = 0; i < M; i++)
01573                 if (ia_local[i] != 0)
01574                     {
01575                     cpl_matrix_set(*covariance, i, j,
01576                                cpl_matrix_get(
01577                                cov, imfit, jmfit));
01578                     imfit += 1;
01579                     }
01580                 
01581                 assert( imfit == Mfit );
01582 
01583                 jmfit += 1;
01584             }
01585         
01586         assert( jmfit == Mfit );
01587         }
01588 
01589         cpl_matrix_delete(cov);
01590     }
01591     
01592     cpl_free(ia_local);
01593     cpl_matrix_delete(alpha);
01594     cpl_matrix_delete(beta);
01595     cpl_free(a_da);
01596     cpl_free(part);
01597    
01598     return CPL_ERROR_NONE;
01599 }
01600 
01601 /*----------------------------------------------------------------------------*/
01611 /*----------------------------------------------------------------------------*/
01612 cpl_error_code
01613 uves_cast_image(cpl_image **image, cpl_type to_type)
01614 {
01615     cpl_image *temp = NULL;
01616     
01617     assure( image != NULL, CPL_ERROR_NULL_INPUT, "Null image");
01618 
01619     temp = *image;
01620     
01621     check( *image = cpl_image_cast(temp, to_type), "Couldn't convert image to %s", 
01622        uves_tostring_cpl_type(to_type));
01623     
01624   cleanup:
01625     uves_free_image(&temp);
01626     return cpl_error_get_code();    
01627 }
01628 
01629 
01630 /*-----------------------------------------------------------------------------*/
01645 /*-----------------------------------------------------------------------------*/
01646 cpl_error_code
01647 uves_crop_image(cpl_image **image, int x1, int y_1, int x2, int y2)
01648 {
01649     cpl_image *temp = NULL;
01650     
01651     assure( image != NULL, CPL_ERROR_NULL_INPUT, "Null image");
01652 
01653     temp = *image;
01654     
01655     check( *image = cpl_image_extract(temp, x1, y_1, x2, y2), 
01656        "Could not extract image");
01657     
01658   cleanup:
01659     uves_free_image(&temp);
01660     return cpl_error_get_code();    
01661 }
01662 
01663 /*----------------------------------------------------------------------------*/
01687 /*----------------------------------------------------------------------------*/
01688 cpl_error_code
01689 uves_get_property_value(const uves_propertylist *plist, const char *keyword, 
01690             cpl_type keywordtype, void *result)
01691 {
01692     cpl_type t;
01693     
01694     /* Check input */
01695     assure( plist != NULL  , CPL_ERROR_NULL_INPUT, "Null property list");
01696     assure( keyword != NULL, CPL_ERROR_NULL_INPUT, "Null keyword");
01697     
01698     /* Check for existence... */
01699     assure( uves_propertylist_contains(plist, keyword), CPL_ERROR_DATA_NOT_FOUND,
01700         "Keyword %s does not exist", keyword );
01701     /* ...and type of keyword */
01702     check( t = uves_propertylist_get_type(plist, keyword) ,
01703        "Could not read type of keyword '%s'", keyword );
01704     assure(t == keywordtype , CPL_ERROR_TYPE_MISMATCH, 
01705        "Keyword '%s' has wrong type (%s). %s expected",
01706        keyword, uves_tostring_cpl_type(t), uves_tostring_cpl_type(keywordtype));
01707     
01708     /* Read the keyword */
01709     switch (keywordtype)
01710     {
01711     case CPL_TYPE_INT   : 
01712         check( *((      int    *)result) =
01713         uves_propertylist_get_int(plist, keyword),
01714         "Could not get (integer) value of %s", keyword );
01715         break;
01716     case CPL_TYPE_BOOL  : 
01717         check( *((      bool   *)result) =
01718            uves_propertylist_get_bool(plist, keyword), 
01719            "Could not get (boolean) value of %s", keyword ); 
01720         break;
01721     case CPL_TYPE_DOUBLE:
01722         check( *((      double *)result) = 
01723            uves_propertylist_get_double(plist, keyword), 
01724            "Could not get (double) value of %s" , keyword );
01725         break;
01726     case CPL_TYPE_STRING: 
01727         check( *((const char * *)result) = 
01728            uves_propertylist_get_string(plist, keyword), 
01729            "Could not get (string) value of %s" , keyword ); 
01730         break;
01731     default:
01732         assure( false, CPL_ERROR_INVALID_TYPE, "Unknown type");
01733         break;
01734     }
01735     
01736   cleanup:
01737     return cpl_error_get_code();
01738 }
01739 
01740 
01741 /*----------------------------------------------------------------------------*/
01766 /*----------------------------------------------------------------------------*/
01767 
01768 const char *
01769 uves_find_frame(const cpl_frameset *frames, const char **wanted, int N, int *found,
01770         const cpl_frame **frame)
01771 {
01772     const char *filename = NULL;  /* Return NULL in case of error */
01773     int i;
01774     const cpl_frame *f = NULL;
01775     
01776     /* Return well-defined pointers in case of error */
01777     *found = 0;
01778     if (frame != NULL)
01779     {
01780         *frame = NULL;
01781     }
01782 
01783     for (i = 0; i < N; i++)
01784     {
01785         check( f = cpl_frameset_find_const(frames, wanted[i]), 
01786            "Could not search through frame set");
01787         if (f != NULL) 
01788         {
01789             check( filename = cpl_frame_get_filename(f), 
01790                "Could not read frame filename");
01791             *found = i;
01792             if (frame != NULL)
01793             {
01794                 *frame = f;
01795             }
01796             /* break */
01797             i = N;
01798         }
01799     }
01800     
01801     /* Set an error if a matching frame wasn't found */
01802     assure(filename != NULL, CPL_ERROR_DATA_NOT_FOUND, "No such frame in frame set");
01803 
01804   cleanup:
01805     return filename;
01806 }
01807 
01808 
01809 /*----------------------------------------------------------------------------*/
01816 /*----------------------------------------------------------------------------*/
01817 int
01818 uves_get_nextensions(const char *filename)
01819 {
01820     int result = 0;
01821     cpl_frame *f = NULL;
01822     
01823     /* CPL only supports reading the number of extensions of a FITS
01824        file if this is in a frame, so create a frame for the specified file */
01825       
01826     check(( f = cpl_frame_new(),
01827         cpl_frame_set_filename(f, filename)),
01828       "Could not create frame");
01829 
01830     check( result = cpl_frame_get_nextensions(f),
01831        "Error reading number of extensions of file '%s'", filename);
01832   cleanup:
01833     cpl_frame_delete(f);
01834     return result;
01835 }
01836 
01837 /*----------------------------------------------------------------------------*/
01858 /*----------------------------------------------------------------------------*/
01859 cpl_error_code
01860 uves_get_parameter(const cpl_parameterlist *parameters, const char *context,
01861            const char *recipe_id, const char *name, cpl_type type, 
01862            void *value)
01863 {
01864     char *fullname = NULL;
01865     const cpl_parameter *p = NULL;
01866     
01867     passure( parameters != NULL, " ");
01868     /* 'context' may be NULL */
01869     passure( recipe_id != NULL, " ");
01870     passure( name != NULL, " ");
01871     passure( value != NULL, " ");
01872 
01873     if (context != NULL)
01874     {
01875         check( fullname = uves_sprintf("%s.%s.%s", context, recipe_id, name),
01876            "Error getting full parameter name");
01877     }
01878     else
01879     {
01880         check( fullname = uves_sprintf("%s.%s", recipe_id, name),
01881            "Error getting full parameter name");
01882     }
01883 
01884     /* Const cast */
01885     check( p = cpl_parameterlist_find_const(parameters, fullname), 
01886        "Error searching for parameter '%s'", fullname);
01887     assure( p != NULL, CPL_ERROR_DATA_NOT_FOUND, 
01888         "No parameter '%s' in parameter list", fullname);
01889     
01890     /* Check that parameter has the correct type */
01891     {
01892     cpl_type partype;
01893     check(  partype = cpl_parameter_get_type(p), 
01894         "Could not read type of parameter '%s'", fullname);
01895     assure( partype == type, CPL_ERROR_TYPE_MISMATCH,
01896         "Parameter '%s' has type %s. Expected type was %s", 
01897         fullname,
01898         uves_tostring_cpl_type(partype), uves_tostring_cpl_type(type));    
01899     }
01900 
01901     /* Read the parameter */
01902     switch(type)
01903     {
01904     case CPL_TYPE_INT:
01905         check( *(int *         )value = 
01906            cpl_parameter_get_int   (p),
01907            "Could not read integer parameter '%s'", fullname);  
01908         break;
01909     case CPL_TYPE_BOOL:
01910         check( *(bool *        )value = 
01911            cpl_parameter_get_bool  (p),
01912            "Could not read boolean parameter '%s'", fullname);  
01913         break;
01914     case CPL_TYPE_DOUBLE:    
01915         check( *(double *      )value = 
01916            cpl_parameter_get_double(p),
01917            "Could not read double parameter '%s'" , fullname );  
01918         break;
01919     case CPL_TYPE_STRING:
01920         check( *(const char ** )value = 
01921            cpl_parameter_get_string(p),
01922            "Could not read string parameter '%s'" , fullname );  
01923         break;
01924     default:
01925         assure(false, CPL_ERROR_UNSUPPORTED_MODE,
01926            "Don't know how to read parameter '%s' of type %s",
01927            fullname, uves_tostring_cpl_type(type));
01928         break;
01929     }
01930     
01931   cleanup:
01932     cpl_free(fullname); fullname = NULL;
01933     return cpl_error_get_code();
01934 }
01935 
01936 /*----------------------------------------------------------------------------*/
01952 /*----------------------------------------------------------------------------*/
01953 cpl_error_code
01954 uves_set_parameter(cpl_parameterlist *parameters, 
01955            const char *context,
01956            const char *name, cpl_type type, void *value)
01957 {
01958     char *fullname = NULL;
01959     cpl_parameter *p = NULL;
01960     
01961     check( fullname = uves_sprintf("%s.%s", context, name),
01962        "Error getting full parameter name");
01963 
01964     /* Const cast */
01965     check( p = cpl_parameterlist_find(parameters, fullname), 
01966        "Error searching for parameter '%s'", fullname);
01967     assure( p != NULL, CPL_ERROR_DATA_NOT_FOUND, 
01968         "No parameter '%s' in parameter list", fullname);
01969 
01970     /* Check that parameter has the correct type */
01971     {
01972     cpl_type partype;
01973     check(  partype = cpl_parameter_get_type(p), 
01974         "Could not read type of parameter '%s'", fullname);
01975     assure( partype == type, CPL_ERROR_TYPE_MISMATCH,
01976         "Parameter '%s' has type %s. Expected type was %s", 
01977         fullname, uves_tostring_cpl_type(partype),
01978         uves_tostring_cpl_type(type));    
01979     }
01980 
01981     /* Set the parameter */
01982     switch(type)
01983     {
01984     case CPL_TYPE_INT:
01985         check( cpl_parameter_set_int   (p, *((int *)    value)), 
01986            "Could not set integer parameter '%s'", fullname);  
01987         break;
01988     case CPL_TYPE_BOOL:    
01989         check( cpl_parameter_set_bool  (p, *((bool *)    value)), 
01990            "Could not set boolean parameter '%s'", fullname);
01991         break;
01992     case CPL_TYPE_DOUBLE:
01993         check( cpl_parameter_set_double(p, *((double *) value)), 
01994            "Could not set double parameter '%s'" , fullname);
01995         break;
01996     case CPL_TYPE_STRING:
01997         check( cpl_parameter_set_string(p, *((char **)  value)),
01998            "Could not set string parameter '%s'" , fullname);  
01999         break;
02000     default:
02001         assure(false, CPL_ERROR_UNSUPPORTED_MODE,
02002            "Don't know how to set parameter of type %s", 
02003            uves_tostring_cpl_type(type));
02004         break;
02005     }
02006     
02007   cleanup:
02008     cpl_free(fullname); fullname = NULL;
02009     return cpl_error_get_code();
02010 }
02011 
02012 /*----------------------------------------------------------------------------*/
02026 /*----------------------------------------------------------------------------*/
02027 
02028 cpl_error_code
02029 uves_set_parameter_default(cpl_parameterlist *parameters, const char *context,
02030                const char *parname, 
02031                cpl_type type, void *value)
02032 {
02033     const char *full_name = NULL;
02034     cpl_parameter *p = NULL;
02035     cpl_type partype;
02036 
02037     if (context != NULL)
02038     {
02039         full_name = uves_sprintf("%s.%s", context, parname);
02040     }
02041     else
02042     {
02043         full_name = uves_sprintf("%s", parname);
02044     }
02045 
02046     if (full_name == NULL)
02047     {
02048         return CPL_ERROR_ILLEGAL_OUTPUT;
02049     }
02050 
02051 
02052     if ( (p = cpl_parameterlist_find(parameters, full_name)) == NULL)
02053     {
02054         uves_msg_error("Missing parameter: '%s'", full_name);
02055 
02056         uves_free_string_const(&full_name);
02057         if (cpl_error_get_code() != CPL_ERROR_NONE) return cpl_error_get_code();
02058         else return CPL_ERROR_DATA_NOT_FOUND;
02059     }
02060     
02061     partype = cpl_parameter_get_type(p);
02062     
02063     if (partype != type)
02064     {
02065         uves_msg_error("Parameter '%s' has type %s. Expected type was %s", 
02066               full_name, uves_tostring_cpl_type(partype), 
02067               uves_tostring_cpl_type(type));
02068 
02069         uves_free_string_const(&full_name);
02070         return CPL_ERROR_TYPE_MISMATCH;
02071     }
02072 
02073     switch(type)
02074     {
02075     case CPL_TYPE_INT:
02076         cpl_parameter_set_default_int   (p, *((int *)    value)); 
02077         break;
02078     case CPL_TYPE_BOOL:
02079         cpl_parameter_set_default_bool  (p, *((bool *)   value)); 
02080         break;
02081     case CPL_TYPE_DOUBLE:
02082         cpl_parameter_set_default_double(p, *((double *) value)); 
02083         break;
02084     case CPL_TYPE_STRING:
02085         cpl_parameter_set_default_string(p, *((char **)  value)); 
02086         break;
02087     default:
02088         uves_msg_error("Unknown type: %s", uves_tostring_cpl_type(type));
02089 
02090         uves_free_string_const(&full_name);
02091         return CPL_ERROR_INVALID_TYPE;
02092     }
02093 
02094     if (cpl_error_get_code() != CPL_ERROR_NONE)
02095     {
02096         uves_msg_error("Error changing value of parameter '%s'", 
02097                full_name);
02098 
02099         uves_free_string_const(&full_name);
02100         return cpl_error_get_code();
02101     }
02102 
02103 
02104     uves_free_string_const(&full_name);
02105     return CPL_ERROR_NONE;
02106 }
02107 
02108 
02109 /*----------------------------------------------------------------------------*/
02121 /*----------------------------------------------------------------------------*/
02122 void 
02123 uves_raise_to_median_frac(cpl_table *t, const char *column, double fraction)
02124 {
02125     int i = 0;
02126     double threshold;
02127 
02128     assure_nomsg(t != NULL, CPL_ERROR_NULL_INPUT);
02129     assure(cpl_table_has_column(t, column), CPL_ERROR_DATA_NOT_FOUND,
02130        "No such column: %s", column);
02131 
02132     assure(cpl_table_get_column_type(t, column) == CPL_TYPE_DOUBLE,
02133        CPL_ERROR_UNSUPPORTED_MODE, "Column %s has type %s. %s expected",
02134        column,
02135        uves_tostring_cpl_type(cpl_table_get_column_type(t, column)),
02136        uves_tostring_cpl_type(CPL_TYPE_DOUBLE));
02137 
02138 
02139     threshold = fraction * cpl_table_get_column_median(t, column);
02140     for (i = 0; i < cpl_table_get_nrow(t); i++)
02141     {
02142         if (cpl_table_get_double(t, column, i, NULL) < threshold)
02143         {
02144             cpl_table_set_double(t, column, i, threshold);
02145         }
02146     }
02147 
02148   cleanup:
02149     return;
02150 }
02151 
02152 
02153 /*----------------------------------------------------------------------------*/
02170 /*----------------------------------------------------------------------------*/
02171 
02172 int
02173 uves_select_table_rows(cpl_table *t,  const char *column, 
02174                cpl_table_select_operator operator, double value)
02175 {
02176     int result = 0;
02177     cpl_type type;
02178     
02179     assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
02180     assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT, 
02181         "No such column: %s", column);
02182 
02183     type = cpl_table_get_column_type(t, column);
02184 
02185     assure( type == CPL_TYPE_DOUBLE || type == CPL_TYPE_FLOAT ||
02186         type == CPL_TYPE_INT, CPL_ERROR_INVALID_TYPE,
02187         "Column '%s' must be double or int. %s found", column, 
02188         uves_tostring_cpl_type(type));
02189 
02190     check( cpl_table_select_all(t), "Error selecting rows");
02191     
02192     if      (type == CPL_TYPE_DOUBLE)
02193     {
02194         result = cpl_table_and_selected_double(t, column, operator, value);
02195     }
02196     else if (type == CPL_TYPE_FLOAT)
02197     {
02198         result = cpl_table_and_selected_float(t, column, operator, value);
02199     }
02200     else if (type == CPL_TYPE_INT)
02201     {
02202         result = cpl_table_and_selected_int(t, column, operator, 
02203                                                 uves_round_double(value));
02204     }
02205     else { /*impossible*/ passure(false, ""); }
02206     
02207   cleanup:
02208     return result;
02209 
02210 }
02211 
02212 /*----------------------------------------------------------------------------*/
02226 /*----------------------------------------------------------------------------*/
02227 int
02228 uves_extract_table_rows_local(cpl_table *t, const char *column,
02229                   cpl_table_select_operator operator, double value)
02230 {
02231     int result = 0;
02232     
02233     assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
02234     assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT, 
02235         "No such column: %s", column);
02236 
02237     check( result = uves_select_table_rows(t, column, operator, value),
02238        "Error selecting rows");
02239 
02240     check( cpl_table_not_selected(t), "Error selecting rows");
02241 
02242     check( cpl_table_erase_selected(t), "Error deleting rows");
02243 
02244   cleanup:
02245     return result;
02246 }
02247 
02248 /*----------------------------------------------------------------------------*/
02267 /*----------------------------------------------------------------------------*/
02268 int
02269 uves_erase_table_rows(cpl_table *t, const char *column,
02270               cpl_table_select_operator operator, double value)
02271 {
02272     int result = 0;
02273     
02274     assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
02275     assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT, 
02276         "No such column: %s", column);
02277 
02278     check( result = uves_select_table_rows(t, column, operator, value),
02279        "Error selecting rows");
02280 
02281     check( cpl_table_erase_selected(t), "Error deleting rows");
02282 
02283   cleanup:
02284     return result;
02285 }
02286 
02287 
02288 
02289 
02290 
02291 /*----------------------------------------------------------------------------*/
02300 /*----------------------------------------------------------------------------*/
02301 
02302 void uves_propertylist_append_property(uves_propertylist *plist, const cpl_property *p)
02303 {
02304     switch(cpl_property_get_type(p)) {
02305     case CPL_TYPE_CHAR:
02306         uves_propertylist_append_char(plist, cpl_property_get_name(p), cpl_property_get_char(p));
02307         break;
02308     case CPL_TYPE_BOOL:
02309         uves_propertylist_append_bool(plist, cpl_property_get_name(p), cpl_property_get_bool(p));
02310         break;
02311     case CPL_TYPE_INT:
02312         uves_propertylist_append_int(plist, cpl_property_get_name(p), cpl_property_get_int(p));
02313         break;
02314     case CPL_TYPE_LONG:
02315         uves_propertylist_append_long(plist, cpl_property_get_name(p), cpl_property_get_long(p));
02316         break;
02317     case CPL_TYPE_FLOAT:
02318         uves_propertylist_append_float(plist, cpl_property_get_name(p), cpl_property_get_float(p));
02319         break;
02320     case CPL_TYPE_DOUBLE:
02321         uves_propertylist_append_double(plist, cpl_property_get_name(p), cpl_property_get_double(p));
02322         break;
02323     case CPL_TYPE_STRING:
02324         uves_propertylist_append_string(plist, cpl_property_get_name(p), cpl_property_get_string(p));
02325         break;
02326     default:
02327         assure( false, CPL_ERROR_UNSUPPORTED_MODE,
02328                 "Type is %s", uves_tostring_cpl_type(cpl_property_get_type(p)));
02329         break;
02330     }
02331   cleanup:
02332     return;
02333 }
02334 
02335 
02336 
02337 /*----------------------------------------------------------------------------*/
02343 /*----------------------------------------------------------------------------*/
02344 int
02345 uves_table_and_selected_invalid(cpl_table *t, const char *column)
02346 {
02347     if (cpl_table_get_column_type(t, column) != CPL_TYPE_STRING)
02348         {
02349             return cpl_table_and_selected_invalid(t, column);
02350         }
02351     else
02352         {
02353             int i = 0;
02354             for (i = 0; i < cpl_table_get_nrow(t); i++)
02355                 {
02356                     if (cpl_table_is_selected(t, i))
02357                         {
02358                             if (cpl_table_is_valid(t, column, i))
02359                                 {
02360                                     cpl_table_unselect_row(t, i);
02361                                 }
02362                             /* else keep it selected */
02363                         }
02364                     /* else unselected, don't change */
02365                 }
02366 
02367             return cpl_table_count_selected(t);
02368         }
02369 }
02370 
02371     
02372 /*----------------------------------------------------------------------------*/
02389 /*----------------------------------------------------------------------------*/
02390 int
02391 uves_erase_invalid_table_rows(cpl_table *t, const char *column)
02392 {
02393     int result = 0;
02394 
02395     assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
02396 
02397     if (column == NULL)
02398     /* Loop through all columns */
02399     {
02400         const char *name;
02401         cpl_table *argument = t;
02402         result = 0;
02403         while ( (name = cpl_table_get_column_name(argument)) != NULL)
02404         {
02405             int n_deleted_rows;
02406 
02407             argument = NULL;
02408             n_deleted_rows = uves_erase_invalid_table_rows(t, name);
02409 
02410             if (n_deleted_rows > 0) 
02411             {
02412                 uves_msg_low("%d rows with invalid '%s' removed", 
02413                      n_deleted_rows, name);
02414             }
02415             result += n_deleted_rows;
02416         }
02417     }
02418     else
02419     {
02420         assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT,
02421             "No such column: %s", column);
02422         
02423         check(( cpl_table_select_all(t),
02424             result = uves_table_and_selected_invalid(t, column), /* workaround here */
02425             cpl_table_erase_selected(t)),              /* and here */
02426                   "Error deleting rows");
02427     }
02428     
02429   cleanup:
02430     return result;
02431 }
02432 
02433 /*----------------------------------------------------------------------------*/
02450 /*----------------------------------------------------------------------------*/
02451 cpl_table *
02452 uves_extract_table_rows(const cpl_table *t, const char *column,
02453             cpl_table_select_operator operator, double value)
02454 {
02455     cpl_table *result = NULL;
02456     
02457     assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
02458     assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT,
02459         "No such column: %s", column);
02460     
02461     /* 1. Extract (duplicate) the entire table
02462        2. remove rows *not* satisfying the criterion */
02463     check(( result = cpl_table_duplicate(t),
02464 
02465         uves_select_table_rows(result, column, operator, value),
02466         cpl_table_not_selected(result),  /* Inverses selection */
02467         
02468         cpl_table_erase_selected(result)),
02469        
02470        "Error extracting rows");
02471 
02472     passure( cpl_table_count_selected(result) == cpl_table_get_nrow(result),
02473              "%d %d",
02474              cpl_table_count_selected(result), cpl_table_get_nrow(result) );
02475     
02476   cleanup:
02477     if (cpl_error_get_code() != CPL_ERROR_NONE)
02478     {
02479         uves_free_table(&result);
02480     }
02481     return result;
02482 }
02483 
02484 /*----------------------------------------------------------------------------*/
02496 /*----------------------------------------------------------------------------*/
02497 void
02498 uves_sort_table_1(cpl_table *t, const char *column, bool reverse)
02499 {
02500     uves_propertylist *plist = NULL;
02501     
02502     assure(t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
02503     assure(cpl_table_has_column(t, column), CPL_ERROR_ILLEGAL_INPUT, 
02504        "No column '%s'", column);
02505 
02506     check(( plist = uves_propertylist_new(),
02507         uves_propertylist_append_bool(plist, column, reverse)),
02508        "Could not create property list for sorting");
02509 
02510     check( uves_table_sort(t, plist), "Could not sort table");
02511 
02512   cleanup:
02513     uves_free_propertylist(&plist);
02514     return;
02515 }
02516 
02517 /*----------------------------------------------------------------------------*/
02533 /*----------------------------------------------------------------------------*/
02534 void
02535 uves_sort_table_2(cpl_table *t, const char *column1, const char *column2, 
02536           bool reverse1, bool reverse2)
02537 {
02538     uves_propertylist *plist = NULL;
02539     
02540     assure(t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
02541     assure(cpl_table_has_column(t, column1), CPL_ERROR_ILLEGAL_INPUT, 
02542        "No column '%s'", column1);
02543     assure(cpl_table_has_column(t, column2), CPL_ERROR_ILLEGAL_INPUT,
02544        "No column '%s'", column2);
02545 
02546     check(( plist = uves_propertylist_new(),
02547         uves_propertylist_append_bool(plist, column1, reverse1),
02548         uves_propertylist_append_bool(plist, column2, reverse2)),
02549        "Could not create property list for sorting");
02550     check( uves_table_sort(t, plist), "Could not sort table");
02551     
02552   cleanup:
02553     uves_free_propertylist(&plist);
02554     return;
02555 }
02556 /*----------------------------------------------------------------------------*/
02571 /*----------------------------------------------------------------------------*/
02572 void
02573 uves_sort_table_3(cpl_table *t, const char *column1, const char *column2, 
02574                   const char *column3,
02575           bool reverse1, bool reverse2, bool reverse3)
02576 {
02577     uves_propertylist *plist = NULL;
02578     
02579     assure(t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
02580     assure(cpl_table_has_column(t, column1), CPL_ERROR_ILLEGAL_INPUT, 
02581        "No column '%s'", column1);
02582     assure(cpl_table_has_column(t, column2), CPL_ERROR_ILLEGAL_INPUT,
02583        "No column '%s'", column2);
02584     assure(cpl_table_has_column(t, column3), CPL_ERROR_ILLEGAL_INPUT,
02585        "No column '%s'", column3);
02586 
02587     check(( plist = uves_propertylist_new(),
02588         uves_propertylist_append_bool(plist, column1, reverse1),
02589         uves_propertylist_append_bool(plist, column2, reverse2),
02590             uves_propertylist_append_bool(plist, column3, reverse3)),
02591         "Could not create property list for sorting");
02592     check( uves_table_sort(t, plist), "Could not sort table");
02593     
02594   cleanup:
02595     uves_free_propertylist(&plist);
02596     return;
02597 }
02598 /*----------------------------------------------------------------------------*/
02603 /*----------------------------------------------------------------------------*/
02604 void uves_free(const void *mem)
02605 {
02606     cpl_free((void *)mem); /* No, it is not a bug. The cast is safe */
02607     return;
02608 }
02609 
02610 /*----------------------------------------------------------------------------*/
02615 /*----------------------------------------------------------------------------*/
02616 void uves_free_image(cpl_image **i) 
02617 {if(i){cpl_image_delete(*i); *i = NULL;}}
02618 /*----------------------------------------------------------------------------*/
02623 /*----------------------------------------------------------------------------*/
02624 void uves_free_mask(cpl_mask **m)
02625 {if(m){cpl_mask_delete(*m); *m = NULL;}}
02626 /*----------------------------------------------------------------------------*/
02631 /*----------------------------------------------------------------------------*/
02632 void uves_free_imagelist(cpl_imagelist **i)
02633 {if(i){cpl_imagelist_delete(*i);        *i = NULL;}}
02634 
02635 /*----------------------------------------------------------------------------*/
02640 /*----------------------------------------------------------------------------*/
02641 void uves_free_table(cpl_table **t)
02642 {if(t){cpl_table_delete(*t);            *t = NULL;}}
02643 
02644 /*----------------------------------------------------------------------------*/
02649 /*----------------------------------------------------------------------------*/
02650 void uves_free_table_const(const cpl_table **t)
02651 {if(t){cpl_table_delete((cpl_table*) (*t));            *t = NULL;}}
02652 
02653 /*----------------------------------------------------------------------------*/
02658 /*----------------------------------------------------------------------------*/
02659 void uves_free_propertylist(uves_propertylist **p)
02660 {if(p){uves_propertylist_delete(*p);     *p = NULL;}}
02661 
02662 /*----------------------------------------------------------------------------*/
02667 /*----------------------------------------------------------------------------*/
02668 void uves_free_propertylist_const(const uves_propertylist **p)
02669 {if(p){uves_propertylist_delete(*p);     *p = NULL;}}
02670 
02671 /*----------------------------------------------------------------------------*/
02676 /*----------------------------------------------------------------------------*/
02677 void uves_free_property(cpl_property **p)
02678 {if(p){cpl_property_delete(*p);     *p = NULL;}}
02679 /*----------------------------------------------------------------------------*/
02684 /*----------------------------------------------------------------------------*/
02685 void uves_free_polynomial(cpl_polynomial **p)
02686 {if(p){cpl_polynomial_delete(*p);       *p = NULL;}}
02687 /*----------------------------------------------------------------------------*/
02692 /*----------------------------------------------------------------------------*/
02693 void uves_free_matrix(cpl_matrix **m)
02694 {if(m){cpl_matrix_delete(*m);           *m = NULL;}}
02695 /*----------------------------------------------------------------------------*/
02700 /*----------------------------------------------------------------------------*/
02701 void uves_free_parameterlist(cpl_parameterlist **p)
02702 {if(p){cpl_parameterlist_delete(*p);    *p = NULL;}}
02703 /*----------------------------------------------------------------------------*/
02708 /*----------------------------------------------------------------------------*/
02709 void uves_free_frameset(cpl_frameset **f)
02710 {if(f){cpl_frameset_delete(*f);    *f = NULL;}}
02711 /*----------------------------------------------------------------------------*/
02716 /*----------------------------------------------------------------------------*/
02717 void uves_free_frame(cpl_frame **f)
02718 {if(f){cpl_frame_delete(*f);    *f = NULL;}}
02719 /*----------------------------------------------------------------------------*/
02724 /*----------------------------------------------------------------------------*/
02725 void uves_free_bivector(cpl_bivector **b)
02726 {if(b){cpl_bivector_delete(*b);           *b = NULL;}}
02727 /*----------------------------------------------------------------------------*/
02732 /*----------------------------------------------------------------------------*/
02733 void uves_free_vector(cpl_vector **v)
02734 {if(v){cpl_vector_delete(*v);           *v = NULL;}}
02735 /*----------------------------------------------------------------------------*/
02740 /*----------------------------------------------------------------------------*/
02741 void uves_free_stats(cpl_stats **s)
02742 {if(s){cpl_stats_delete(*s);            *s = NULL;}}
02743 /*----------------------------------------------------------------------------*/
02748 /*----------------------------------------------------------------------------*/
02749 void uves_unwrap_vector(cpl_vector **v)
02750 {if(v){cpl_vector_unwrap(*v);           *v = NULL;}}
02751 /*----------------------------------------------------------------------------*/
02756 /*----------------------------------------------------------------------------*/
02757 void uves_unwrap_vector_const(const cpl_vector **v)
02758 {if(v){cpl_vector_unwrap((cpl_vector*) (*v));           *v = NULL;}}
02759 /*----------------------------------------------------------------------------*/
02764 /*----------------------------------------------------------------------------*/
02765 void uves_unwrap_bivector_vectors(cpl_bivector **b)
02766 {if(b){cpl_bivector_unwrap_vectors(*b); *b = NULL;}}
02767 /*----------------------------------------------------------------------------*/
02772 /*----------------------------------------------------------------------------*/
02773 void uves_free_array(cpl_array **a)
02774 {if(a){cpl_array_delete(*a);           *a = NULL;}}
02775 
02776 /*----------------------------------------------------------------------------*/
02781 /*----------------------------------------------------------------------------*/
02782 void uves_free_int(int **i)
02783 {if(i){cpl_free(*i);           *i = NULL;}}
02784 
02785 /*----------------------------------------------------------------------------*/
02790 /*----------------------------------------------------------------------------*/
02791 void uves_free_int_const(const int **i)
02792 {if(i){uves_free(*i);           *i = NULL;}}
02793 /*----------------------------------------------------------------------------*/
02798 /*----------------------------------------------------------------------------*/
02799 void uves_free_float(float **f)
02800 {if(f){cpl_free(*f);           *f = NULL;}}
02801 
02802 /*----------------------------------------------------------------------------*/
02807 /*----------------------------------------------------------------------------*/
02808 void uves_free_double(double **d)
02809 {if(d){cpl_free(*d);           *d = NULL;}}
02810 
02811 /*----------------------------------------------------------------------------*/
02816 /*----------------------------------------------------------------------------*/
02817 void uves_free_string(char **s)
02818 {if(s){cpl_free(*s);           *s = NULL;}}
02819 /*----------------------------------------------------------------------------*/
02824 /*----------------------------------------------------------------------------*/
02825 void uves_free_string_const(const char **s)
02826 {if(s){uves_free(*s);           *s = NULL;}}
02827 
02828 /*----------------------------------------------------------------------------*/
02849 /*----------------------------------------------------------------------------*/
02850 
02851 static double
02852 get_chisq(int N, int D,
02853       int (*f)(const double x[], const double a[], double *result),
02854       const double *a,
02855       const double *x,
02856       const double *y,
02857       const double *sigma)
02858 {
02859     double chi_sq;     /* Result */
02860     int i = 0;
02861 
02862     /* For efficiency, don't check input in this static function */
02863 
02864     chi_sq = 0.0;
02865     for (i = 0; i < N; i++)
02866     {
02867         double fx_i;
02868         double residual;                 /* Residual in units of uncertainty */
02869         const double *x_i = &(x[0+i*D]);
02870 
02871         /* Evaluate */
02872         cpl_ensure( f(x_i,
02873               a,
02874               &fx_i) == 0, CPL_ERROR_ILLEGAL_INPUT, -1.0);
02875 
02876         /* Accumulate */
02877         if (sigma == NULL)
02878         {
02879             residual = (fx_i - y[i]);
02880         }
02881         else
02882         {
02883             residual = (fx_i - y[i]) / sigma[i];
02884         }
02885 
02886         chi_sq += residual*residual;
02887        
02888     }
02889 
02890     return chi_sq;
02891 }
02892 /*----------------------------------------------------------------------------*/
02926 /*----------------------------------------------------------------------------*/
02927 static int
02928 get_candidate(const double *a, const int ia[],
02929           int M, int N, int D,
02930           double lambda,
02931           int    (*f)(const double x[], const double a[], double *result),
02932           int (*dfda)(const double x[], const double a[], double result[]),
02933           const double *x,
02934           const double *y,
02935           const double *sigma,
02936           double *partials,
02937           cpl_matrix *alpha,
02938           cpl_matrix *beta,
02939           double *a_da)
02940 {
02941     int Mfit = 0;         /* Number of non-constant fit parameters */
02942     cpl_matrix *da;       /* Solution of   alpha * da = beta */
02943     double *alpha_data;
02944     double *beta_data;
02945     double *da_data;
02946     int i, imfit = 0;
02947     int j, jmfit = 0;
02948     int k = 0;
02949 
02950     /* For efficiency, don't check input in this static function */
02951 
02952     Mfit = cpl_matrix_get_nrow(alpha);
02953 
02954     alpha_data    = cpl_matrix_get_data(alpha);
02955     beta_data     = cpl_matrix_get_data(beta);
02956    
02957     /* Build alpha, beta:
02958      *
02959      *  alpha[i,j] = sum_{k=1,N} (sigma_k)^-2 * df/da_i * df/da_j  *
02960      *                           (1 + delta_ij lambda) ,
02961      *
02962      *   beta[i]   = sum_{k=1,N} (sigma_k)^-2 * ( y_k - f(x_k) ) * df/da_i
02963      *
02964      * where (i,j) loop over the non-constant parameters (0 to Mfit-1),
02965      * delta is Kronecker's delta, and all df/da are evaluated in x_k
02966      */
02967 
02968     cpl_matrix_fill(alpha, 0.0);
02969     cpl_matrix_fill(beta , 0.0);
02970 
02971     for (k = 0; k < N; k++)
02972     {
02973         double sm2 = 0.0;                /* (sigma_k)^-2 */
02974         double fx_k = 0.0;               /* f(x_k)       */
02975         const double *x_k = &(x[0+k*D]); /* x_k          */
02976 
02977         if (sigma == NULL)
02978         {
02979             sm2 = 1.0;
02980         }
02981         else
02982         {
02983             sm2 = 1.0 / (sigma[k] * sigma[k]);
02984         }
02985         
02986         /* Evaluate f(x_k) */
02987         cpl_ensure( f(x_k, a, &fx_k) == 0, CPL_ERROR_ILLEGAL_INPUT, -1);
02988 
02989         /* Evaluate (all) df/da (x_k) */
02990         cpl_ensure( dfda(x_k, a, partials) == 0, 
02991             CPL_ERROR_ILLEGAL_INPUT, -1);
02992 
02993         for (i = 0, imfit = 0; i < M; i++)
02994         {
02995             if (ia[i] != 0)
02996             {
02997                 /* Beta */
02998                 beta_data[imfit] +=
02999                 sm2 * (y[k] - fx_k) * partials[i];
03000                 
03001                 /* Alpha is symmetrical, so compute
03002                    only lower-left part */
03003                 for (j = 0, jmfit = 0; j < i; j++)
03004                 {
03005                     if (ia[j] != 0)
03006                     {
03007                         alpha_data[jmfit + imfit*Mfit] +=
03008                         sm2 * partials[i] * 
03009                         partials[j];
03010                         
03011                         jmfit += 1;
03012                     }
03013                 }
03014                 
03015                 /* Alpha, diagonal terms */
03016                 j = i;
03017                 jmfit = imfit;
03018                 
03019                 alpha_data[jmfit + imfit*Mfit] += 
03020                 sm2 * partials[i] *
03021                 partials[j] * (1 + lambda);
03022                 
03023                 imfit += 1;
03024             }
03025         }
03026 
03027         assert( imfit == Mfit );
03028     }
03029 
03030     /* Create upper-right part of alpha */
03031     for (i = 0, imfit = 0; i < M; i++)
03032     {
03033         if (ia[i] != 0)
03034         {
03035             for (j = i+1, jmfit = imfit+1; j < M; j++)
03036             {
03037                 if (ia[j] != 0)
03038                 {
03039                     alpha_data[jmfit + imfit*Mfit] = 
03040                     alpha_data[imfit + jmfit*Mfit];
03041                     
03042                     jmfit += 1;
03043                 }
03044             }
03045             assert( jmfit == Mfit );
03046 
03047             imfit += 1;
03048         }
03049     }
03050     assert( imfit == Mfit );
03051     
03052     da = cpl_matrix_solve(alpha, beta);
03053 
03054     cpl_ensure(da != NULL, cpl_error_get_code(), -1);
03055 
03056     /* Create a+da vector by adding a and da */
03057     da_data   = cpl_matrix_get_data(da);
03058 
03059     for (i = 0, imfit = 0; i < M; i++)
03060     {
03061         if (ia[i] != 0)
03062         {
03063             a_da[i] = a[i] + da_data[0 + imfit*1];
03064             
03065             imfit += 1;
03066         }
03067         else
03068         {
03069             a_da[i] = a[i];
03070         }
03071     }
03072     
03073     assert( imfit == Mfit );
03074 
03075     cpl_matrix_delete(da);
03076 
03077     return 0;
03078 }
03079 

Generated on 8 Mar 2011 for UVES Pipeline Reference Manual by  doxygen 1.6.1