irplib_strehl.c

00001 /* $Id: irplib_strehl.c,v 1.37 2008/03/25 13:31:33 llundin Exp $
00002  *
00003  * This file is part of the irplib package
00004  * Copyright (C) 2002,2003 European Southern Observatory
00005  *
00006  * This program is free software; you can redistribute it and/or modify
00007  * it under the terms of the GNU General Public License as published by
00008  * the Free Software Foundation; either version 2 of the License, or
00009  * (at your option) any later version.
00010  *
00011  * This program is distributed in the hope that it will be useful,
00012  * but WITHOUT ANY WARRANTY; without even the implied warranty of
00013  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014  * GNU General Public License for more details.
00015  *
00016  * You should have received a copy of the GNU General Public License
00017  * along with this program; if not, write to the Free Software
00018  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02111-1307  USA
00019  */
00020 
00021 /*
00022  * $Author: llundin $
00023  * $Date: 2008/03/25 13:31:33 $
00024  * $Revision: 1.37 $
00025  * $Name: uves-4_2_2 $
00026  */
00027 
00028 #ifdef HAVE_CONFIG_H
00029 #include <config.h>
00030 #endif
00031 
00032 /*-----------------------------------------------------------------------------
00033                                    Includes
00034  -----------------------------------------------------------------------------*/
00035 
00036 #include <string.h>
00037 #include <assert.h>
00038 #include <math.h>
00039 #include <float.h>
00040 
00041 #include <cpl.h>
00042 
00043 #include "irplib_utils.h"
00044 #include "irplib_strehl.h"
00045 
00046 /*----------------------------------------------------------------------------*/
00050 /*----------------------------------------------------------------------------*/
00051 
00052 /*-----------------------------------------------------------------------------
00053                                    Define
00054  -----------------------------------------------------------------------------*/
00055 
00056 #ifndef IRPLIB_STREHL_RAD_CENTRAL
00057 #define IRPLIB_STREHL_RAD_CENTRAL 5
00058 #endif
00059 
00060 #define IRPLIB_DISK_BG_MIN_PIX_NB    30
00061 #define IRPLIB_DISK_BG_REJ_LOW       0.1
00062 #define IRPLIB_DISK_BG_REJ_HIGH      0.1
00063 
00064 /*-----------------------------------------------------------------------------
00065                                    Functions prototypes
00066  -----------------------------------------------------------------------------*/
00067 
00068 static cpl_image * irplib_strehl_generate_otf(double, double, double, double, 
00069                                               int, double);
00070 static double PSF_H1(double, double, double);
00071 static double PSF_H2(double, double);
00072 static double PSF_G(double, double);
00073 static double PSF_sinc(double);
00074 static double PSF_TelOTF(double, double);
00075 
00076 
00077 /*-----------------------------------------------------------------------------
00078                                    Functions code
00079  -----------------------------------------------------------------------------*/
00082 /*----------------------------------------------------------------------------*/
00111 /*----------------------------------------------------------------------------*/
00112 cpl_error_code irplib_strehl_compute(const cpl_image *   im,
00113                                      double              m1,
00114                                      double              m2,
00115                                      double              lam,
00116                                      double              dlam,
00117                                      double              pscale,
00118                                      int                 size,
00119                                      double              xpos,
00120                                      double              ypos,
00121                                      double              r1,
00122                                      double              r2,
00123                                      double              r3,
00124                                      int                 noise_box_sz,
00125                                      int                 noise_nsamples,
00126                                      double          *   strehl,
00127                                      double          *   strehl_err,
00128                                      double          *   star_bg,
00129                                      double          *   star_peak,
00130                                      double          *   star_flux,
00131                                      double          *   psf_peak,
00132                                      double          *   psf_flux,
00133                                      double          *   bg_noise)
00134 {
00135     cpl_image  * psf;
00136     double       star_radius, max_radius;
00137 
00138     /* FIXME: Arbitrary choice of image border */
00139     const double window_size = (double)(IRPLIB_STREHL_RAD_CENTRAL);
00140 
00141     /* Determined empirically by C. Lidman for Strehl error computation */
00142     const double strehl_error_coefficient = CPL_MATH_PI * 0.007 / 0.0271;
00143     double       ring[4];
00144 
00145     /* Check compile-time constant */
00146     cpl_ensure_code(window_size > 0.0,  CPL_ERROR_ILLEGAL_INPUT);
00147 
00148     /* Test inputs */
00149     cpl_ensure_code(im != NULL,         CPL_ERROR_NULL_INPUT);
00150     cpl_ensure_code(strehl != NULL,     CPL_ERROR_NULL_INPUT);
00151     cpl_ensure_code(strehl_err != NULL, CPL_ERROR_NULL_INPUT);
00152     cpl_ensure_code(star_bg != NULL,    CPL_ERROR_NULL_INPUT);
00153     cpl_ensure_code(star_peak != NULL,  CPL_ERROR_NULL_INPUT);
00154     cpl_ensure_code(star_flux != NULL,  CPL_ERROR_NULL_INPUT);
00155     cpl_ensure_code(psf_peak != NULL,   CPL_ERROR_NULL_INPUT);
00156     cpl_ensure_code(psf_flux != NULL,   CPL_ERROR_NULL_INPUT);
00157 
00158     cpl_ensure_code(pscale > 0.0,      CPL_ERROR_ILLEGAL_INPUT);
00159 
00160     cpl_ensure_code(r1 > 0.0,      CPL_ERROR_ILLEGAL_INPUT);
00161     cpl_ensure_code(r2 > 0.0,      CPL_ERROR_ILLEGAL_INPUT);
00162     cpl_ensure_code(r3 > r2,       CPL_ERROR_ILLEGAL_INPUT);
00163 
00164     /* Computing a Strehl ratio is a story between an ideal PSF */
00165     /* and a candidate image supposed to approximate this ideal PSF. */
00166 
00167     /* Generate first appropriate PSF to find max peak */
00168     psf = irplib_strehl_generate_psf(m1, m2, lam, dlam, pscale, size);
00169     cpl_ensure_code(psf != NULL,      CPL_ERROR_ILLEGAL_OUTPUT);
00170 
00171     /* Compute flux in PSF and find max peak */
00172     *psf_peak = cpl_image_get_max(psf);
00173     cpl_image_delete(psf);
00174 
00175     assert( *psf_peak > 0.0); /* The ideal PSF has a positive maximum */
00176     *psf_flux = 1.0; /* The psf flux, cpl_image_get_flux(psf), is always 1 */
00177 
00178     /* Measure the background in the candidate image */
00179     *star_bg = irplib_strehl_ring_background(im, xpos, ypos, r2/pscale, r3/pscale,
00180                                              IRPLIB_BG_METHOD_AVER_REJ);
00181 
00182     /* Compute star_radius in pixels */
00183     star_radius = r1/pscale;
00184 
00185     /* Measure the flux on the candidate image */
00186     *star_flux = irplib_strehl_disk_flux(im, xpos, ypos, star_radius, *star_bg);
00187 
00188     cpl_ensure_code(*star_flux > 0.0,      CPL_ERROR_ILLEGAL_OUTPUT);
00189 
00190     /* Find the peak value on the central part of the candidate image */
00191     max_radius = window_size < star_radius ? window_size : star_radius;
00192     cpl_ensure_code(!irplib_strehl_disk_max(im, xpos, ypos, max_radius,
00193                                             star_peak), cpl_error_get_code());
00194     *star_peak -= *star_bg;
00195 
00196     cpl_ensure_code(*star_peak > 0.0,      CPL_ERROR_ILLEGAL_OUTPUT);
00197 
00198     /* Compute Strehl */
00199     /* (StarPeak / StarFlux) / (PsfPeak / PsfFlux) */
00200     *strehl = (*star_peak * *psf_flux ) / ( *star_flux * *psf_peak);
00201 
00202     if (*strehl > 1)
00203         cpl_msg_warning(cpl_func, "Extreme Strehl-ratio=%g, star_peak=%g, "
00204                         "star_flux=%g, psf_peak=%g, psf_flux=%g", *strehl,
00205                         *star_peak, *star_flux, *psf_peak, *psf_flux);
00206 
00207     /* Compute Strehl error */
00208     ring[0] = xpos;
00209     ring[1] = ypos;
00210     ring[2] = r2/pscale;
00211     ring[3] = r3/pscale;
00212     cpl_ensure_code(!cpl_flux_get_noise_ring(im, ring, noise_box_sz,
00213                                              noise_nsamples, bg_noise, NULL),
00214                     cpl_error_get_code());
00215 
00216     *strehl_err = strehl_error_coefficient * (*bg_noise) * pscale * 
00217         star_radius * star_radius / *star_flux;
00218 
00219     /* This check should not be able to fail, but just to be sure */
00220     cpl_ensure_code(*strehl_err >= 0.0,       CPL_ERROR_ILLEGAL_OUTPUT);
00221 
00222     return CPL_ERROR_NONE;
00223 }
00224 
00225 /*----------------------------------------------------------------------------*/
00241 /*----------------------------------------------------------------------------*/
00242 double irplib_strehl_disk_flux(const cpl_image * im,
00243                                double            xpos,
00244                                double            ypos,
00245                                double            rad,
00246                                double            bg)
00247 {
00248     const double    sqr = rad * rad;
00249     double          sqrest;
00250     const float *   pim;
00251     double          flux = 0.0;
00252     double          yj, xi; 
00253     int             nx, ny;
00254     int             lx, ly, ux, uy;
00255     int             i, j;
00256 
00257 
00258     /* Check entries */
00259     cpl_ensure(im != NULL, CPL_ERROR_NULL_INPUT, 0.0);
00260     cpl_ensure(cpl_image_get_type(im) == CPL_TYPE_FLOAT,
00261                CPL_ERROR_UNSUPPORTED_MODE, 0.0);
00262     cpl_ensure(rad > 0.0, CPL_ERROR_ILLEGAL_INPUT, 0.0);
00263  
00264     nx = cpl_image_get_size_x(im);
00265     ny = cpl_image_get_size_y(im);
00266 
00267     /* Round down */
00268     lx = (int)(xpos - rad);
00269     ly = (int)(ypos - rad);
00270     if (lx < 0) lx = 0;
00271     if (ly < 0) ly = 0;
00272 
00273     /* Round up */
00274     ux = (int)(xpos + rad) + 1;
00275     uy = (int)(ypos + rad) + 1;
00276     if (ux > (nx-1)) ux = nx-1;
00277     if (uy > (ny-1)) uy = ny-1;
00278 
00279     /* im is actually not modified */
00280     pim = cpl_image_get_data_float((cpl_image*)im);
00281     for (j=ly ; j<uy ; j++) {
00282         yj = (double)j - ypos;
00283         sqrest = sqr - yj * yj;
00284         for (i=lx; i<ux ; i++) {
00285             xi = (double)i - xpos;
00286             if (sqrest >= xi * xi && irplib_isnan(pim[i+j*nx]) == 0) {
00287                 flux += (double)pim[i+j*nx] - bg;
00288             }
00289         }
00290     }
00291     return flux;
00292 }
00293 
00294 /*----------------------------------------------------------------------------*/
00308 /*----------------------------------------------------------------------------*/
00309 double irplib_strehl_ring_background(const cpl_image *       im,
00310                                      double                  xpos,
00311                                      double                  ypos,
00312                                      double                  rad_int,
00313                                      double                  rad_ext,
00314                                      irplib_strehl_bg_method mode)
00315 {
00316     int             npix;
00317     const double    sqr_int = rad_int * rad_int;
00318     const double    sqr_ext = rad_ext * rad_ext;
00319     double          dist;
00320     cpl_vector  *   pix_arr;
00321     const float *   pim;
00322     double          flux = 0.0;
00323     double          yj, xi; 
00324     int             lx, ly, ux, uy;
00325     int             nx, ny;
00326     int             i, j;
00327 
00328     /* Check entries */
00329     cpl_ensure(im != NULL, CPL_ERROR_NULL_INPUT, 0.0);
00330     cpl_ensure(cpl_image_get_type(im) == CPL_TYPE_FLOAT,
00331                CPL_ERROR_UNSUPPORTED_MODE, 0.0);
00332     cpl_ensure(rad_int > 0.0, CPL_ERROR_ILLEGAL_INPUT, 0.0);
00333     cpl_ensure(rad_ext > rad_int, CPL_ERROR_ILLEGAL_INPUT, 0.0);
00334 
00335     cpl_ensure(mode == IRPLIB_BG_METHOD_AVER_REJ ||
00336                mode == IRPLIB_BG_METHOD_MEDIAN,
00337                CPL_ERROR_UNSUPPORTED_MODE, 0.0);
00338 
00339     nx = cpl_image_get_size_x(im);
00340     ny = cpl_image_get_size_y(im);
00341 
00342     /* Round down */
00343     lx = (int)(xpos - rad_ext);
00344     ly = (int)(ypos - rad_ext);
00345     if (lx < 0) lx = 0;
00346     if (ly < 0) ly = 0;
00347 
00348     /* Round up */
00349     ux = (int)(xpos + rad_ext) + 1;
00350     uy = (int)(ypos + rad_ext) + 1;
00351     if (ux > (nx-1)) ux = nx-1;
00352     if (uy > (ny-1)) uy = ny-1;
00353 
00354     npix = (ux - lx + 1) * (uy - ly + 1);
00355     cpl_ensure(npix >= IRPLIB_DISK_BG_MIN_PIX_NB, CPL_ERROR_DATA_NOT_FOUND, 0.0);
00356 
00357     /* Allocate pixel array to hold values in the ring */
00358     pix_arr = cpl_vector_new(npix);
00359 
00360     /* Count number of pixels in the ring */
00361     /* Retrieve all pixels which belong to the ring */
00362     /* im is actually _not_ modified */
00363     pim = cpl_image_get_data_float((cpl_image*)im);
00364     npix = 0;
00365     for (j=ly ; j<uy ; j++) {
00366         yj = (double)j - ypos;
00367         for (i=lx ; i<ux; i++) {
00368             xi = (double)i - xpos;
00369             dist = yj * yj + xi * xi;
00370             if (sqr_int <= dist && dist <= sqr_ext && 
00371                 irplib_isnan(pim[i+j*nx]) == 0) {
00372                 cpl_vector_set(pix_arr, npix, (double)pim[i+j*nx]);
00373                 npix++;
00374             }
00375         }
00376     }
00377 
00378     if (npix < IRPLIB_DISK_BG_MIN_PIX_NB) {
00379         cpl_vector_delete(pix_arr);
00380         cpl_ensure(0, CPL_ERROR_DATA_NOT_FOUND, 0.0);
00381     }
00382 
00383     /* Should not be able to fail now */
00384 
00385     /* Resize pixel array to actual number of values within the ring */
00386     cpl_vector_set_size(pix_arr, npix);
00387 
00388     if (mode == IRPLIB_BG_METHOD_AVER_REJ) {
00389         const int low_ind  = (int)((double)npix * IRPLIB_DISK_BG_REJ_LOW);
00390         const int high_ind = (int)((double)npix
00391                                    * (1.0 - IRPLIB_DISK_BG_REJ_HIGH));
00392 
00393         /* Sort the array */
00394         cpl_vector_sort(pix_arr, 1);
00395 
00396         for (i=low_ind ; i<high_ind ; i++) {
00397             flux += cpl_vector_get(pix_arr, i);
00398         }
00399         if (high_ind - low_ind > 1) flux /= (double)(high_ind - low_ind);
00400     } else /* if (mode == IRPLIB_BG_METHOD_MEDIAN) */ {
00401         flux = cpl_vector_get_median(pix_arr);
00402     }
00403 
00404     cpl_vector_delete(pix_arr);
00405 
00406     return flux;
00407 }
00408 
00409 /*----------------------------------------------------------------------------*/
00429 /*----------------------------------------------------------------------------*/
00430 cpl_image * irplib_strehl_generate_psf(
00431         double   m1,
00432         double   m2,
00433         double   lam,
00434         double   dlam,
00435         double   pscale,
00436         int      size)
00437 {
00438     cpl_image * otf_image = irplib_strehl_generate_otf(m1, m2, lam, dlam, 
00439             size, pscale);
00440 
00441     if (otf_image == NULL) return NULL;
00442 
00443     /* Transform back to real space
00444        - Normalization is unnecessary, due to the subsequent normalisation.
00445        - An OTF is point symmetric about its center, i.e. it is even,
00446          i.e. the real space image is real.
00447        - Because of this a forward FFT works as well.
00448        - If the PSF ever needs to have its images halves swapped add
00449          CPL_FFT_SWAP_HALVES to the FFT call.
00450      */
00451 
00452     if (cpl_image_fft(otf_image, NULL, CPL_FFT_UNNORMALIZED) ||
00453 
00454         /* Compute absolute values of PSF */
00455         cpl_image_abs(otf_image) ||
00456 
00457         /* Normalize PSF to get flux=1  */
00458         cpl_image_normalise(otf_image, CPL_NORM_FLUX)) {
00459 
00460         cpl_image_delete(otf_image);
00461         return NULL;
00462     }
00463     
00464     return otf_image;
00465 }
00466 
00469 /*----------------------------------------------------------------------------*/
00485 /*----------------------------------------------------------------------------*/
00486 static cpl_image * irplib_strehl_generate_otf(
00487         double  m1,
00488         double  m2,
00489         double  lam,
00490         double  dlam,
00491         int     size,
00492         double  pscale)
00493 {
00494     cpl_image   *   otf_image;
00495     double      *   otf_data;
00496     double          obs_ratio ;  /* m1 / m2    */
00497     double          f_max ;      /* cut-off frequency        */
00498     int             pix0 ;      /* Pixel corresponding to the zero frequency */
00499     double          a, x, y;
00500     double          f, rsq, fc, invfc, lambda;
00501     double          sincy;
00502     double          invsize;
00503     register int    pos;
00504     int             i, j, k;
00505 
00506 
00507     cpl_ensure(m2   > 0.0,      CPL_ERROR_ILLEGAL_INPUT, NULL);
00508     cpl_ensure(m1   > m2,       CPL_ERROR_ILLEGAL_INPUT, NULL);
00509     cpl_ensure(lam  > 0.0,      CPL_ERROR_ILLEGAL_INPUT, NULL);
00510     cpl_ensure(dlam > 0.0,      CPL_ERROR_ILLEGAL_INPUT, NULL);
00511     cpl_ensure(size > 0,        CPL_ERROR_ILLEGAL_INPUT, NULL);
00512     cpl_ensure(pscale > 0.0,    CPL_ERROR_ILLEGAL_INPUT, NULL);
00513 
00514     /* Convert pixel scale from sec to radians, microns in meters    */
00515     pscale /= (double)206265;
00516     lam /= (double)1.0e6;
00517     dlam /= (double)1.0e6;
00518 
00519     /* Obscuration ratio    */
00520     obs_ratio = m2 / m1;
00521 
00522     /* Pixel corresponding to the zero frequency    */
00523     pix0 = size/2;
00524     invsize = (double)1.0 / (double)size;
00525 
00526     /* Cut-off frequency in pixels  */
00527     f_max = m1 * pscale * (double)size / lam;
00528 
00529     /* Allocate for output image    */
00530     otf_image = cpl_image_new(size, size, CPL_TYPE_DOUBLE);
00531     if (otf_image==NULL) return NULL;
00532     otf_data = cpl_image_get_data_double(otf_image);
00533 
00534     /* Now compute the OTF  */
00535     /* OPTIMIZED CODE !!! LIMITED READABILITY !!!   */
00536 
00537     for (k=1 ; k<=9 ; k++) {    /* iteration on the wavelength  */
00538         /* Compute intermediate cut-off frequency   */
00539         lambda = (double)(lam - dlam*(double)(k-5)/8.0);
00540         fc = (double)f_max * (double)lam / lambda;
00541         invfc = 1.0 / fc;
00542 
00543         /* Convolution with the detector pixels */
00544         pos = 0;
00545         for (j=0 ; j<size ; j++) {
00546             y = (double)(j-pix0);
00547             sincy = PSF_sinc(CPL_MATH_PI * y * invsize);
00548             for (i=0 ; i<size ; i++) {
00549                 x = (double)(i-pix0);
00550                 rsq = x*x + y*y;
00551                 if (rsq < fc*fc) {
00552                     if (rsq < 0.01)
00553                         a = 1.0;
00554                     else {
00555                         f = sqrt(rsq) * invfc;
00556                         a = PSF_TelOTF(f,obs_ratio) *
00557                             PSF_sinc(CPL_MATH_PI * x * invsize) * sincy;
00558                     }
00559                 } else {
00560                     a = 0.0;
00561                 }
00562                 otf_data[pos++] += a / 9.0;
00563             }
00564         }
00565     }
00566     return otf_image;
00567 }
00568 
00569 /*----------------------------------------------------------------------------*
00570  * H1 function
00571  *----------------------------------------------------------------------------*/
00572 static double PSF_H1(
00573         double  f,
00574         double  u,
00575         double  v)
00576 {
00577     const double e = fabs(1.0-v) > 0.0 ? -1.0 : 1.0; /* e = 1.0 iff v = 1.0 */
00578 
00579     return((v*v/CPL_MATH_PI)*acos((f/v)*(1.0+e*(1.0-u*u)/(4.0*f*f))));
00580 }
00581 
00582 /*----------------------------------------------------------------------------*
00583  * H2 function
00584  *----------------------------------------------------------------------------*/
00585 static double PSF_H2(double  f,
00586                      double  u)
00587 {
00588     const double tmp1 = (2.0 * f) / (1.0 + u);
00589     const double tmp2 = (1.0 - u) / (2.0 * f);
00590 
00591     return -1.0 * (f/CPL_MATH_PI) * (1.0+u)
00592         * sqrt((1.0-tmp1*tmp1)*(1.0-tmp2*tmp2));
00593 }   
00594 
00595 /*----------------------------------------------------------------------------*
00596  * G function
00597  *----------------------------------------------------------------------------*/
00598 static double PSF_G(double  f, 
00599                     double  u)
00600 {
00601     if (f <= (1.0-u)/2.0) return(u*u);
00602     if (f >= (1.0+u)/2.0) return(0.0);
00603     else return(PSF_H1(f,u,1.0) + PSF_H1(f,u,u) + PSF_H2(f,u));
00604 }   
00605 
00606 /*----------------------------------------------------------------------------*
00607  * sinc function
00608  *----------------------------------------------------------------------------*/
00609 static double PSF_sinc(double x)
00610 {
00611   return fabs(x) > fabs(sin(x)) ? sin(x)/x : 1.0;
00612 }
00613 
00614 /*----------------------------------------------------------------------------*
00615  * Telescope OTF function 
00616  *----------------------------------------------------------------------------*/
00617 static double PSF_TelOTF(double  f,  
00618                          double  u)
00619 {
00620     return((PSF_G(f,1.0)+u*u*PSF_G(f/u,1.0)-2.0*PSF_G(f,u))/(1.0-u*u));
00621 }
00622 
00623 /*----------------------------------------------------------------------------*/
00635 /*----------------------------------------------------------------------------*/
00636 cpl_error_code irplib_strehl_disk_max(const cpl_image * self,
00637                                              double            xpos,
00638                                              double            ypos,
00639                                              double            radius,
00640                                              double          * ppeak)
00641 {
00642 
00643     const double    sqr = radius * radius;
00644     double          sqrest;
00645     const float *   pself;
00646     float           peak = FLT_MAX;  /* Avoid (false) uninit warning */
00647     double          yj, xi; 
00648     int             nx, ny;
00649     int             lx, ly, ux, uy;
00650     int             i, j;
00651     cpl_boolean     first = CPL_TRUE;
00652 
00653 
00654     /* Check entries */
00655     cpl_ensure_code(ppeak != NULL, CPL_ERROR_NULL_INPUT);
00656     cpl_ensure_code(self  != NULL, CPL_ERROR_NULL_INPUT);
00657     cpl_ensure_code(cpl_image_get_type(self) == CPL_TYPE_FLOAT,
00658                     CPL_ERROR_UNSUPPORTED_MODE);
00659     cpl_ensure_code(radius > 0.0, CPL_ERROR_ILLEGAL_INPUT);
00660  
00661     nx = cpl_image_get_size_x(self);
00662     ny = cpl_image_get_size_y(self);
00663 
00664     /* Round down */
00665     lx = (int)(xpos - radius);
00666     ly = (int)(ypos - radius);
00667     if (lx < 0) lx = 0;
00668     if (ly < 0) ly = 0;
00669 
00670     /* Round up */
00671     ux = (int)(xpos + radius) + 1;
00672     uy = (int)(ypos + radius) + 1;
00673     if (ux > (nx-1)) ux = nx-1;
00674     if (uy > (ny-1)) uy = ny-1;
00675 
00676     /* im is actually not modified */
00677     pself = cpl_image_get_data_float((cpl_image*)self);
00678     for (j=ly ; j<uy ; j++) {
00679         yj = (double)j - ypos;
00680         sqrest = sqr - yj * yj;
00681         for (i=lx; i<ux ; i++) {
00682             xi = (double)i - xpos;
00683             if (sqrest >= xi * xi && irplib_isnan(pself[i+j*nx]) == 0) {
00684                 if (first || pself[i+j*nx] > peak) {
00685                     first = CPL_FALSE;
00686                     peak = pself[i+j*nx];
00687                 }
00688             }
00689         }
00690     }
00691 
00692     cpl_ensure_code(!first, CPL_ERROR_DATA_NOT_FOUND);
00693 
00694     *ppeak = (double)peak;
00695 
00696     return CPL_ERROR_NONE;
00697 }

Generated on Mon Apr 21 10:56:54 2008 for UVES Pipeline Reference Manual by  doxygen 1.5.1