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

Generated on Tue Jun 19 14:39:20 2007 for UVES Pipeline Reference Manual by  doxygen 1.4.6