uves_wavecal_search.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: jmlarsen $
00022  * $Date: 2007/08/21 13:08:26 $
00023  * $Revision: 1.26 $
00024  * $Name: uves-3_9_0 $
00025  * $Log: uves_wavecal_search.c,v $
00026  * Revision 1.26  2007/08/21 13:08:26  jmlarsen
00027  * Removed irplib_access module, largely deprecated by CPL-4
00028  *
00029  * Revision 1.25  2007/06/06 08:17:34  amodigli
00030  * replace tab with 4 spaces
00031  *
00032  * Revision 1.24  2007/05/02 13:20:01  jmlarsen
00033  * Take error bars into account in line searching if arclamp was flat-fielded
00034  *
00035  * Revision 1.23  2007/04/24 12:50:29  jmlarsen
00036  * Replaced cpl_propertylist -> uves_propertylist which is much faster
00037  *
00038  * Revision 1.22  2007/04/20 14:46:45  jmlarsen
00039  * Added commented out code
00040  *
00041  * Revision 1.21  2007/03/05 10:25:08  jmlarsen
00042  * Include slope in Gaussian fit
00043  *
00044  * Revision 1.20  2007/02/23 13:33:38  jmlarsen
00045  * Added code to test unweighted fitting
00046  *
00047  * Revision 1.19  2007/02/22 15:38:26  jmlarsen
00048  * Use linear background term in Gaussian fit
00049  *
00050  * Revision 1.18  2006/11/15 15:02:15  jmlarsen
00051  * Implemented const safe workarounds for CPL functions
00052  *
00053  * Revision 1.16  2006/11/15 14:04:08  jmlarsen
00054  * Removed non-const version of parameterlist_get_first/last/next which is already
00055  * in CPL, added const-safe wrapper, unwrapper and deallocator functions
00056  *
00057  * Revision 1.15  2006/11/06 15:19:42  jmlarsen
00058  * Removed unused include directives
00059  *
00060  * Revision 1.14  2006/08/18 13:51:01  jmlarsen
00061  * Moved one message from info to debug level
00062  *
00063  * Revision 1.13  2006/08/17 13:56:53  jmlarsen
00064  * Reduced max line length
00065  *
00066  * Revision 1.12  2006/08/17 09:18:47  jmlarsen
00067  * Removed CPL2 code
00068  *
00069  * Revision 1.11  2006/08/11 14:38:24  jmlarsen
00070  * Minor text output change
00071  *
00072  * Revision 1.10  2006/08/11 09:01:17  jmlarsen
00073  * Set unextracted bins to zero flux rather than marking them as bad
00074  *
00075  * Revision 1.9  2006/07/14 12:45:58  jmlarsen
00076  * Removed a few messages
00077  *
00078  * Revision 1.8  2006/07/03 13:29:30  jmlarsen
00079  * Adapted to new 1d-fitting function interface
00080  *
00081  * Revision 1.7  2006/06/13 12:05:11  jmlarsen
00082  * Shortened max line length
00083  *
00084  * Revision 1.6  2006/05/12 15:06:30  jmlarsen
00085  * Killed code for method = gravity
00086  *
00087  * Revision 1.5  2006/04/24 09:34:26  jmlarsen
00088  * Adapted to new interface of gaussian fitting routine
00089  *
00090  * Revision 1.4  2006/03/03 13:54:11  jmlarsen
00091  * Changed syntax of check macro
00092  *
00093  * Revision 1.3  2006/02/15 13:19:15  jmlarsen
00094  * Reduced source code max. line length
00095  *
00096  * Revision 1.2  2006/02/08 07:52:16  jmlarsen
00097  * Added function returning library version
00098  *
00099  * Revision 1.34  2006/01/12 15:41:14  jmlarsen
00100  * Moved gauss. fitting to irplib
00101  *
00102  * Revision 1.33  2005/12/20 08:11:44  jmlarsen
00103  * Added CVS  entry
00104  *
00105  */
00106 
00107 /*----------------------------------------------------------------------------*/
00111 /*----------------------------------------------------------------------------*/
00114 #ifdef HAVE_CONFIG_H
00115 #  include <config.h>
00116 #endif
00117 
00118 #include <uves_wavecal_search.h>
00119 #include <uves_utils.h>
00120 #include <uves_utils_wrappers.h>
00121 #include <uves_utils_cpl.h>
00122 #include <uves_pfits.h>
00123 #include <uves_dump.h>
00124 #include <uves_error.h>
00125 #include <uves_msg.h>
00126 
00127 #include <cpl.h>
00128 #include <float.h>
00129 
00130 #define FIT_SLOPE 1
00131 #define WEIGHTED_FIT 1   /* Define to zero to get unweighted fit of emmision line
00132                             (like MIDAS) */
00133 
00134 static double
00135 xcenter(const cpl_image *image, const cpl_image *noise, int xlo, int xhi, int row,
00136     centering_method CENTERING_METHOD, int bin_disp,
00137     double *sigma, double *intensity, double *dx0, double *slope, double *background);
00138 
00139 static cpl_error_code
00140 detect_lines(const cpl_image *spectrum, const cpl_image *noise, 
00141          const uves_propertylist *spectrum_header, 
00142              bool flat_fielded,
00143          int RANGE, double THRESHOLD, centering_method CENTERING_METHOD, 
00144              int bin_disp,
00145          const polynomial *order_locations, cpl_image *arcframe, 
00146          cpl_table *linetable, 
00147          int *ndetected, int *nrows);
00148 
00149 /*----------------------------------------------------------------------------*/
00179 /*----------------------------------------------------------------------------*/
00180 cpl_table *
00181 uves_wavecal_search(const cpl_image *spectrum, const cpl_image *noise,
00182             const uves_propertylist *spectrum_header,
00183                     bool flat_fielded,
00184             const polynomial *order_locations, cpl_image *arcframe,
00185             int RANGE, int MINLINES, int MAXLINES,
00186             centering_method CENTERING_METHOD,
00187                     int bin_disp)
00188 {
00189     cpl_table *linetable = NULL;       /* Result */
00190 
00191     int nx, ny, norders;               /* Dimensions of raw image, number of orders */
00192     double threshold_low;              /* Threshold limits used for binary search */
00193     double threshold_high;
00194     double threshold = 0;             /* Current threshold */
00195     int lines_in_table;               /* Number of lines in line table */
00196     int lines_detected;               /* Number of lines actually found */
00197     bool max_thresh_found = false;    /* Is 'threshold_high' large enough? */
00198 
00199     passure( spectrum        != NULL, "Null input spectrum");
00200     passure( order_locations != NULL, "Null polynomial");
00201     passure( arcframe        != NULL, "Null raw image");
00202 
00203     if (flat_fielded) {
00204         assure( cpl_image_get_type(spectrum) == CPL_TYPE_DOUBLE,
00205                 CPL_ERROR_TYPE_MISMATCH,
00206                 "Spectrum image type is %s, must be double",
00207                 uves_tostring_cpl_type(cpl_image_get_type(spectrum)));
00208     }
00209     
00210     check(( nx      = cpl_image_get_size_x(spectrum),
00211         norders = cpl_image_get_size_y(spectrum)), "Error reading input spectrum");
00212     check( ny      = cpl_image_get_size_y(arcframe), "Error reading input image");
00213     assure(nx == cpl_image_get_size_x(arcframe), CPL_ERROR_INCOMPATIBLE_INPUT, 
00214        "Spectrum and image widths are different (%d and %d)",
00215        nx, cpl_image_get_size_x(arcframe));
00216     
00217     assure( MINLINES <= MAXLINES, CPL_ERROR_ILLEGAL_INPUT, 
00218         "minlines=%d maxlines=%d", MINLINES, MAXLINES );
00219     
00220     /* Initialize result line table */
00221     check(( linetable = cpl_table_new(MAXLINES),
00222         cpl_table_new_column(linetable, "X"     , CPL_TYPE_DOUBLE),
00223         cpl_table_new_column(linetable, "dX"    , CPL_TYPE_DOUBLE),
00224         cpl_table_new_column(linetable, "Xwidth", CPL_TYPE_DOUBLE),
00225         cpl_table_new_column(linetable, "Y"     , CPL_TYPE_INT),
00226         cpl_table_new_column(linetable, "Peak"  , CPL_TYPE_DOUBLE),
00227         cpl_table_new_column(linetable, "Background" , CPL_TYPE_DOUBLE),
00228         cpl_table_new_column(linetable, "Slope" , CPL_TYPE_DOUBLE)),
00229       "Could not create line table");
00230     
00231     uves_msg("Searching for emission lines");
00232 
00233     threshold_low  = 0.0;
00234 
00235     /* This start (guess) value is doubled until too few lines are detected */
00236     if (flat_fielded) {
00237         threshold_high = 10.0; /* dimensionless, number of stdevs above continuum */
00238     }
00239     else {
00240         threshold_high = cpl_image_get_mean(spectrum);
00241 
00242         assure( threshold_high > 0, CPL_ERROR_ILLEGAL_INPUT,
00243                 "Spectrum median flux is %e. Must be positive",
00244                 cpl_image_get_median(spectrum));
00245     }
00246     
00247     max_thresh_found = false;
00248 
00249     /* Detect lines and adjust threshold
00250        until MINLINES <= lines_detected <= MAXLINES */
00251     lines_detected = 0;
00252     while( (lines_detected < MINLINES || MAXLINES < lines_detected) && 
00253         fabs(threshold_low - threshold_high) > DBL_EPSILON )
00254     {
00255         threshold = (threshold_low + threshold_high)/2.0;
00256 
00257         check( detect_lines(spectrum, noise, spectrum_header,
00258                                 flat_fielded,
00259                 RANGE, threshold, CENTERING_METHOD,
00260                                 bin_disp,
00261                 order_locations,
00262                 NULL,
00263                 linetable,
00264                 &lines_detected,
00265                 &lines_in_table),
00266            "Could not search for emission lines");
00267         
00268         /* Update threshold */
00269         /* 'threshold_high' is doubled until the threshold is too high.
00270            Then a binary search is performed. */
00271         if (lines_detected < MINLINES)
00272         {
00273             max_thresh_found = true;
00274             threshold_high = threshold;
00275         }
00276         else if (MAXLINES < lines_detected) 
00277         {
00278             if (!max_thresh_found)
00279             {
00280                 threshold_high *= 2;
00281             }
00282             else
00283             {                
00284                 threshold_low = threshold;
00285             }
00286         }
00287     }
00288 
00289     assure( MINLINES <= lines_detected && lines_detected <= MAXLINES, 
00290         CPL_ERROR_CONTINUE,
00291         "Could not detect between %d and %d lines. Try to increase search range",
00292         MINLINES, MAXLINES);
00293     
00294     /* Draw detections on input image  */
00295     check( detect_lines(spectrum, noise, spectrum_header,
00296                         flat_fielded,
00297             RANGE, threshold, CENTERING_METHOD,
00298                         bin_disp,
00299             order_locations,
00300             arcframe,
00301             linetable,
00302             &lines_detected,
00303             &lines_in_table),
00304        "Could not search for emission lines");
00305     
00306     /* Remove the last part of the line table (garbage) */
00307     check( cpl_table_set_size(linetable, lines_in_table), 
00308        "Could not resize line table");
00309     
00310     uves_sort_table_1(linetable, "X", false);
00311     
00312   cleanup:
00313 #if 0 /* if flat-fielded */
00314     uves_free_image(&temp);
00315 #endif
00316     if (cpl_error_get_code() != CPL_ERROR_NONE)
00317     {
00318         uves_free_table(&linetable);
00319     }
00320     else
00321     {
00322         /* Returned is... */
00323         passure( cpl_table_get_ncol(linetable) == 7, "%d", 
00324              cpl_table_get_ncol(linetable));
00325         passure( cpl_table_has_column(linetable, "X"     ), " ");
00326         passure( cpl_table_has_column(linetable, "dX"    ), " ");
00327         passure( cpl_table_has_column(linetable, "Xwidth"), " ");
00328         passure( cpl_table_has_column(linetable, "Y"     ), " ");
00329         passure( cpl_table_has_column(linetable, "Peak"  ), " ");
00330         passure( cpl_table_has_column(linetable, "Background" ), " ");
00331         passure( cpl_table_has_column(linetable, "Slope" ), " ");
00332         
00333     }
00334     return linetable;
00335 }
00336 
00337 /*----------------------------------------------------------------------------*/
00388 /*----------------------------------------------------------------------------*/
00389 static cpl_error_code
00390 detect_lines(const cpl_image *spectrum, const cpl_image *noise, 
00391          const uves_propertylist *spectrum_header, 
00392              bool flat_fielded,
00393          int RANGE, double THRESHOLD, centering_method CENTERING_METHOD, 
00394              int bin_disp,
00395          const polynomial *order_locations, cpl_image *arcframe, 
00396          cpl_table *linetable, 
00397          int *ndetected, int *nrows)
00398 {
00399     int norders;      /* Number of orders */
00400     int minorder;     /* Relative order number of first row in spectrum image */
00401     int MAXLINES;     /* Number of rows in line table (max no. of 
00402              lines to search for) */
00403     int nx;           /* Width of spectrum (and raw image) */
00404     int x, order;     /* 'order' always counts from 1 */
00405     
00406     const double *spectrum_data;
00407     const double *noise_data;
00408 
00409     /* Check input */
00410     passure( spectrum        != NULL, " ");
00411     passure( noise           != NULL, " ");
00412     passure( spectrum_header != NULL, " ");
00413     nx      = cpl_image_get_size_x(spectrum);
00414     norders = cpl_image_get_size_y(spectrum);
00415     
00416     /* For efficiency reasons, get direct pointer to buffer,
00417        support only CPL_TYPE_DOUBLE */
00418     assure( cpl_image_get_type(spectrum) == CPL_TYPE_DOUBLE,
00419         CPL_ERROR_UNSUPPORTED_MODE,
00420         "Image type must be double. It is %s", 
00421         uves_tostring_cpl_type(cpl_image_get_type(spectrum)));
00422 
00423     spectrum_data = cpl_image_get_data_double_const(spectrum);
00424     noise_data    = cpl_image_get_data_double_const(noise);
00425 
00426     passure( RANGE > 0, "%d", RANGE);
00427 
00428     if (arcframe != NULL)
00429     {
00430         passure( order_locations != NULL, " ");
00431         passure( nx == cpl_image_get_size_x(arcframe), 
00432              "%d %d", nx, cpl_image_get_size_x(arcframe));
00433     }
00434     
00435     passure( linetable != NULL, " ");
00436     MAXLINES = cpl_table_get_nrow(linetable);
00437     passure( cpl_table_get_ncol(linetable) == 7, "%d", 
00438          cpl_table_get_ncol(linetable));
00439     passure( cpl_table_has_column(linetable, "X"     ), " ");
00440     passure( cpl_table_has_column(linetable, "dX"    ), " ");
00441     passure( cpl_table_has_column(linetable, "Xwidth"), " ");
00442     passure( cpl_table_has_column(linetable, "Y"     ), " ");
00443     passure( cpl_table_has_column(linetable, "Peak"  ), " ");
00444     passure( cpl_table_has_column(linetable, "Background" ), " ");
00445     passure( cpl_table_has_column(linetable, "Slope" ), " ");
00446     
00447     assure( THRESHOLD > 0, CPL_ERROR_ILLEGAL_INPUT, "Illegal threshold: %e",
00448         THRESHOLD);
00449 
00450     check( minorder = uves_pfits_get_crval2(spectrum_header), 
00451        "Error reading order number of first row");
00452 
00453     *ndetected = 0;    /* Counts the number of lines detected so far. */
00454     *nrows = 0;        /* A pointer to the first unused row in the
00455                   line table */
00456     
00457     for (order = minorder; order < minorder + norders; order++) {
00458         int spectrum_row = order - minorder + 1;
00459         int ndetected_order = 0;
00460         for (x = 1; x <= nx; x++) {
00461         double flux, dflux;
00462         int peak_width = 0;
00463         int xlo, xhi;
00464         double local_median;
00465         
00466         /* Check if there is a peak and determine its position and width */
00467         // flux = cpl_image_get(spectrum, x, spectrum_row, &pis_rejected);
00468         flux  = spectrum_data[(x-1) + (spectrum_row - 1) * nx];
00469                 dflux = noise_data   [(x-1) + (spectrum_row - 1) * nx];
00470         
00471         xlo = uves_max_int(x - RANGE, 1);
00472         xhi = uves_min_int(x + RANGE, nx);
00473         
00474         local_median = cpl_image_get_median_window(
00475             spectrum,
00476             uves_max_int(xlo, 1 ), spectrum_row,
00477             uves_min_int(xhi, nx), spectrum_row);
00478         
00479         while(x <= nx && 
00480                       (
00481                           (!flat_fielded && flux - local_median > THRESHOLD) 
00482                           ||
00483                           (flat_fielded && (flux - local_median) > THRESHOLD * dflux)
00484                           )
00485                     ) {
00486 #if WANT_BIG_LOGFILE
00487             uves_msg_debug("threshold = %f\tx = %d\tflux = %f\tmedian = %f", 
00488                    THRESHOLD, x, flux, local_median);
00489 #endif
00490             
00491             x += 1;
00492             peak_width += 1;
00493             
00494             if (x <= nx) {
00495             /* flux =
00496                cpl_image_get(spectrum, x,
00497                spectrum_row, &pis_rejected);
00498             */
00499             flux = spectrum_data[(x-1) + (spectrum_row - 1) * nx];
00500             xlo = uves_max_int(x - RANGE, 1);
00501             xhi = uves_min_int(x + RANGE, nx);
00502             local_median = cpl_image_get_median_window(
00503                 spectrum,
00504                 uves_max_int(xlo, 1 ), spectrum_row,
00505                 uves_min_int(xhi, nx), spectrum_row);
00506             }
00507         }
00508         /* x is now first position that is below (median + threshold) */
00509         
00510         if (peak_width > 0) {
00511             double x_peak, dx = 0, sigma, slope, back;
00512             check( x_peak = xcenter(spectrum, noise,
00513                         uves_max_int(1, x - peak_width), 
00514                         /* First position above threshold */ 
00515                         uves_max_int(1, x - 1),          
00516                         /* Last  position above threshold */ 
00517                         spectrum_row,
00518                         CENTERING_METHOD,
00519                                             bin_disp,
00520                         &sigma,
00521                         &flux,
00522                         &dx,
00523                                             &slope,
00524                                             &back),
00525                "Could not locate peak center");
00526             
00527 #if WANT_BIG_LOGFILE
00528             uves_msg_debug("(Order, x, flux) = (%d, %f, %f)", 
00529                    order, x_peak, flux);
00530 #endif        
00531             /* Add line to line table, but only if less
00532                lines that MAXLINES have been detected */
00533             if (*nrows < MAXLINES) {
00534             check(( cpl_table_set_int   (linetable, "Y"     , *nrows, order),
00535                 cpl_table_set_double(linetable, "X"     , *nrows, x_peak),
00536                 cpl_table_set_double(linetable, "dX"    , *nrows, dx),
00537                 cpl_table_set_double(linetable, "Xwidth", *nrows, sigma),
00538                 cpl_table_set_double(linetable, "Peak"  , *nrows, flux),
00539                 cpl_table_set_double(linetable, "Background" , *nrows, back),
00540                 cpl_table_set_double(linetable, "Slope" , *nrows, slope)),
00541                   "Could not update line table row %d", *nrows);
00542             (*nrows)++;
00543             }
00544             
00545             ndetected_order++;
00546             (*ndetected)++;
00547             
00548             if (arcframe != NULL) {
00549             int x1;
00550             int pen = 0;  /* Value to write */
00551             int ny = cpl_image_get_size_y(arcframe);
00552             /* We already know 'nx' from the width of the spectrum image */
00553             
00554             for (x1  = uves_max_int(
00555                  1 , uves_round_double(
00556                      x_peak - peak_width - 0*RANGE/2.0)); 
00557                  x1 <= uves_min_int(
00558                  nx, uves_round_double(
00559                      x_peak + peak_width + 0*RANGE/2.0)); 
00560                  x1++) {
00561                 check( cpl_image_set(
00562                        arcframe,
00563                        x1,
00564                        uves_min_int(
00565                        ny, 
00566                        uves_max_int(
00567                            1, 
00568                            (int) uves_polynomial_evaluate_2d(
00569                            order_locations, x1, order)     
00570                            )), 
00571                        pen),
00572                    "Error writing input image");
00573                 check( cpl_image_set(
00574                        arcframe,
00575                        uves_min_int(
00576                        nx,
00577                        uves_max_int((int) x_peak, 1)),
00578                        uves_min_int(
00579                        ny, 
00580                        uves_max_int(
00581                            1,
00582                            (int) uves_polynomial_evaluate_2d(
00583                            order_locations, x1, order)
00584                            - 10)), 
00585                        pen),
00586                    "Error writing input image");
00587             }
00588             }
00589         } /* line found */
00590         }/* for x */
00591         if (arcframe != NULL) uves_msg_debug("Order #%d: %d lines detected", 
00592                          order, ndetected_order);
00593     }/* for order */
00594     
00595     /* Remove doublets */
00596     {
00597     int i;
00598     int doublets_removed = 0;
00599     for (i = 0; i+1 < *nrows; i++) {
00600         if (fabs(cpl_table_get_double(linetable, "X", i  , NULL) - 
00601              cpl_table_get_double(linetable, "X", i+1, NULL))  < 2.0) 
00602         {
00603             /* If a doublet was found, delete it. 
00604                Make sure the table stays the same size
00605                by adding two rows at the end. */
00606 
00607             check( cpl_table_erase_window(linetable, i, 2),
00608                "Error removing rows");
00609             *nrows -= 2;
00610             *ndetected -= 2;
00611             
00612             check( cpl_table_set_size(linetable, 
00613                           cpl_table_get_nrow(linetable) + 2),
00614                "Could not resize line table");
00615             
00616             doublets_removed++;
00617         }
00618     }
00619     if (doublets_removed > 0)
00620         {
00621         uves_msg_debug("%d doublet%s removed", 
00622                    doublets_removed, doublets_removed > 1 ? "s" : "");
00623         }
00624     }
00625     
00626             uves_msg("Range = %d pixels; threshold = %.2f %s; %d lines detected", 
00627          RANGE, THRESHOLD, flat_fielded ? "stdev" : "ADU", *ndetected);
00628     
00629   cleanup:    
00630     return cpl_error_get_code();
00631 }
00632 
00633 /*----------------------------------------------------------------------------*/
00658 /*----------------------------------------------------------------------------*/
00659 static double
00660 xcenter(const cpl_image *image, const cpl_image *noise, int xlo, int xhi, int row,
00661     centering_method CENTERING_METHOD, int bin_disp,
00662     double *sigma, double *intensity, double *dx0, double *slope, double *background)
00663 {
00664     double x0;           /* Result */
00665     cpl_matrix *covariance = NULL;
00666     const double *image_data;
00667     bool converged;
00668     int lo_r, hi_r;
00669     
00670     int nx = cpl_image_get_size_x(image);
00671 
00672     passure(cpl_image_get_type(image) == CPL_TYPE_DOUBLE, " ");
00673 
00674     image_data = cpl_image_get_data_double_const(image);
00675 
00676     /* Make sure fit window is 13-19 pixels
00677        (7-9 pixels for binned CCD) */
00678     lo_r = 6;
00679     hi_r = 8;
00680     if (bin_disp >= 2) 
00681         {
00682             lo_r = 4;
00683             hi_r = 5;
00684         }
00685 
00686     {
00687         int xm = (xlo+xhi)/2;
00688 
00689         xlo = uves_max_int(1, xm - lo_r);
00690         xhi = uves_min_int(nx, xm + lo_r);
00691     }
00692 
00693     /* Increase fit window (up to 19 pixels) */
00694     do {
00695         converged = true;
00696         if (1 < xlo && 0 <
00697             //cpl_image_get(image, xlo - 1, row, &pis_rejected) &&
00698             //cpl_image_get(image, xlo - 1, row, &pis_rejected) <
00699             //cpl_image_get(image, xlo    , row, &pis_rejected) )
00700             image_data[(xlo-1-1) + (row - 1) * nx] &&
00701             image_data[(xlo-1-1) + (row - 1) * nx] <
00702             image_data[(xlo  -1) + (row - 1) * nx] )
00703             {
00704                 converged = false;
00705                 xlo -= 1;
00706             }
00707         
00708         if (xhi < nx && 0 <
00709             //cpl_image_get(image, xhi + 1, row, &pis_rejected) && 
00710             //cpl_image_get(image, xhi + 1, row, &pis_rejected) <
00711             //cpl_image_get(image, xhi    , row, &pis_rejected) )
00712             image_data[(xhi+1-1) + (row - 1) * nx] &&
00713             image_data[(xhi+1-1) + (row - 1) * nx] <
00714             image_data[(xhi  -1) + (row - 1) * nx] )
00715             {
00716                 converged = false;
00717                 xhi += 1;
00718             }
00719 
00720         if ((xhi-xlo+1) >= hi_r) 
00721             {
00722                 converged = true;
00723             }
00724 
00725     } while (!converged);
00726 
00727     /* Get precise location */
00728     if (CENTERING_METHOD == CENTERING_GAUSSIAN)
00729     {
00730 #if WEIGHTED_FIT
00731         uves_fit_1d_image(image, noise, NULL,
00732 #else /* Unweighted fit like MIDAS which gives larger sigma */
00733         uves_fit_1d_image(image, NULL, NULL,
00734 #endif
00735                   true, false, false,
00736                   xlo, xhi, row,
00737                   &x0, sigma, intensity, background, slope,
00738 #if WEIGHTED_FIT
00739                   NULL, NULL, &covariance,
00740 #else
00741                   NULL, NULL, NULL,
00742 #endif
00743 
00744 #if FIT_SLOPE
00745                               uves_gauss_linear, uves_gauss_linear_derivative, 5);
00746 #else
00747                               uves_gauss, uves_gauss_derivative, 4);
00748             *slope = 0;
00749 #endif
00750 
00751         /* The fitting routine sometimes (i.e. regularly) fails
00752          * because of low statistics.
00753          * Recover from specific fitting errors */
00754         if (cpl_error_get_code() == CPL_ERROR_NONE)
00755         {
00756             /* Variance is guaranteed to be positive */
00757 #if WEIGHTED_FIT
00758             *dx0 = sqrt(cpl_matrix_get(covariance, 0, 0));
00759 #else
00760             *dx0 = *sigma / sqrt(*intensity);
00761 #endif
00762 
00763 #if WANT_BIG_LOGFILE
00764              uves_msg_debug("Gaussian fit succeeded at (x, row, N) = (%f, %d, %d)",
00765                     x0, row, xhi-xlo+1);
00766 #endif
00767         }
00768         else if (cpl_error_get_code() == CPL_ERROR_CONTINUE)
00769         {
00770             /* Fitting failed */
00771             uves_error_reset();
00772 #if WANT_BIG_LOGFILE
00773             uves_msg_debug("Gaussian fit failed at (x, row, N) ="
00774                    " (%f, %d, %d), using centroid", 
00775                    x0, row, xhi-xlo+1);
00776 #endif
00777             *dx0 = *sigma / sqrt(*intensity);
00778         }
00779         else if (cpl_error_get_code() == CPL_ERROR_SINGULAR_MATRIX)
00780         {
00781             uves_error_reset();
00782             
00783             /* Fitting succeeded but covariance computation failed */
00784             uves_msg_debug("Covariance matrix computation failed");
00785             *dx0 = *sigma / sqrt(*intensity);
00786         }
00787         
00788         assure(cpl_error_get_code() == CPL_ERROR_NONE, cpl_error_get_code(),
00789            "Gaussian fitting failed");
00790 
00791 #if WANT_BIG_LOGFILE
00792         uves_msg_debug("Fit   = (x0=%f, sigma=%f, norm=%f, backg=%f, N=%d)", 
00793                x0,
00794                *sigma,
00795                *intensity,
00796                background,
00797                xhi - xlo + 1);
00798 #endif
00799 
00800         /* 'intensity' is the norm (area) of the gaussian fit.
00801            But we need to return the height above zero
00802            (for MIDAS compatibility). 
00803            height = f(x0) = background + norm/sqrt(2pi sigma^2)
00804         */
00805 
00806         *intensity = *background + (*intensity)/(sqrt(2*M_PI) * (*sigma));
00807 
00808     }
00809     else   /*  if (CENTERING_METHOD == CENTERING_GRAVITY) */
00810     {
00811         assure (false, CPL_ERROR_UNSUPPORTED_MODE,
00812             "Centering method (no. %d) is unsupported", 
00813             CENTERING_METHOD);
00814     }
00815 
00816   cleanup:
00817     uves_free_matrix(&covariance);
00818     return x0;
00819 }
00820 

Generated on Fri Apr 18 14:11:44 2008 for UVES Pipeline Reference Manual by  doxygen 1.5.1