irplib_strehl.c

00001 /* $Id: irplib_strehl.c,v 1.43 2009/11/18 21:37:48 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: 2009/11/18 21:37:48 $
00024  * $Revision: 1.43 $
00025  * $Name: HEAD $
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 static cpl_error_code update_bad_pixel_map(cpl_image* im);
00076 
00077 
00078 /*-----------------------------------------------------------------------------
00079                                    Functions code
00080  -----------------------------------------------------------------------------*/
00088 cpl_error_code update_bad_pixel_map(cpl_image* im)
00089 {
00090     int szx = cpl_image_get_size_x(im);
00091     int szy = cpl_image_get_size_y(im);
00092     int x = 0;
00093     cpl_mask* bpm = cpl_image_get_bpm(im);
00094 
00095     for (x = 1; x <=szx; x++)
00096     {
00097         int y = 0;
00098         for(y = 1; y <= szy; y++)
00099         {
00100             int isnull = 0;
00101             double value = cpl_image_get(im, x, y, &isnull);
00102             if (isnan(value))
00103             {
00104                 cpl_mask_set(bpm, x, y, CPL_BINARY_1);
00105             }
00106         }
00107     }
00108     return cpl_error_get_code();
00109 }
00140 cpl_error_code irplib_strehl_mark_bad_and_compute(cpl_image *   im,
00141                                      double              m1,
00142                                      double              m2,
00143                                      double              lam,
00144                                      double              dlam,
00145                                      double              pscale,
00146                                      int                 size,
00147                                      double              xpos,
00148                                      double              ypos,
00149                                      double              r1,
00150                                      double              r2,
00151                                      double              r3,
00152                                      int                 noise_box_sz,
00153                                      int                 noise_nsamples,
00154                                      double          *   strehl,
00155                                      double          *   strehl_err,
00156                                      double          *   star_bg,
00157                                      double          *   star_peak,
00158                                      double          *   star_flux,
00159                                      double          *   psf_peak,
00160                                      double          *   psf_flux,
00161                                      double          *   bg_noise)
00162 {
00163     cpl_ensure_code(!update_bad_pixel_map(im), cpl_error_get_code());
00164     return irplib_strehl_compute(im, m1, m2, lam, dlam, pscale, size, xpos, ypos,
00165             r1,
00166             r2,
00167             r3,
00168             noise_box_sz,
00169             noise_nsamples,
00170             strehl,
00171             strehl_err,
00172             star_bg,
00173             star_peak,
00174             star_flux,
00175             psf_peak,
00176             psf_flux,
00177             bg_noise);
00178 }
00179 
00180 /*----------------------------------------------------------------------------*/
00211 /*----------------------------------------------------------------------------*/
00212 cpl_error_code irplib_strehl_compute(const cpl_image *   im,
00213                                      double              m1,
00214                                      double              m2,
00215                                      double              lam,
00216                                      double              dlam,
00217                                      double              pscale,
00218                                      int                 size,
00219                                      double              xpos,
00220                                      double              ypos,
00221                                      double              r1,
00222                                      double              r2,
00223                                      double              r3,
00224                                      int                 noise_box_sz,
00225                                      int                 noise_nsamples,
00226                                      double          *   strehl,
00227                                      double          *   strehl_err,
00228                                      double          *   star_bg,
00229                                      double          *   star_peak,
00230                                      double          *   star_flux,
00231                                      double          *   psf_peak,
00232                                      double          *   psf_flux,
00233                                      double          *   bg_noise)
00234 {
00235     cpl_image  * psf;
00236     double       star_radius, max_radius;
00237 
00238     /* FIXME: Arbitrary choice of image border */
00239     const double window_size = (double)(IRPLIB_STREHL_RAD_CENTRAL);
00240 
00241     /* Determined empirically by C. Lidman for Strehl error computation */
00242     const double strehl_error_coefficient = CPL_MATH_PI * 0.007 / 0.0271;
00243     double       ring[4];
00244     /* cpl_flux_get_noise_ring() must succeed with this many tries */
00245     int          ring_tries = 3;
00246     cpl_errorstate prestate;
00247 
00248     /* Check compile-time constant */
00249     cpl_ensure_code(window_size > 0.0,  CPL_ERROR_ILLEGAL_INPUT);
00250 
00251     /* Test inputs */
00252     cpl_ensure_code(im != NULL,         CPL_ERROR_NULL_INPUT);
00253     cpl_ensure_code(strehl != NULL,     CPL_ERROR_NULL_INPUT);
00254     cpl_ensure_code(strehl_err != NULL, CPL_ERROR_NULL_INPUT);
00255     cpl_ensure_code(star_bg != NULL,    CPL_ERROR_NULL_INPUT);
00256     cpl_ensure_code(star_peak != NULL,  CPL_ERROR_NULL_INPUT);
00257     cpl_ensure_code(star_flux != NULL,  CPL_ERROR_NULL_INPUT);
00258     cpl_ensure_code(psf_peak != NULL,   CPL_ERROR_NULL_INPUT);
00259     cpl_ensure_code(psf_flux != NULL,   CPL_ERROR_NULL_INPUT);
00260 
00261     cpl_ensure_code(pscale > 0.0,      CPL_ERROR_ILLEGAL_INPUT);
00262 
00263     cpl_ensure_code(r1 > 0.0,      CPL_ERROR_ILLEGAL_INPUT);
00264     cpl_ensure_code(r2 > 0.0,      CPL_ERROR_ILLEGAL_INPUT);
00265     cpl_ensure_code(r3 > r2,       CPL_ERROR_ILLEGAL_INPUT);
00266 
00267     /* Computing a Strehl ratio is a story between an ideal PSF */
00268     /* and a candidate image supposed to approximate this ideal PSF. */
00269 
00270     /* Generate first appropriate PSF to find max peak */
00271     psf = irplib_strehl_generate_psf(m1, m2, lam, dlam, pscale, size);
00272     cpl_ensure_code(psf != NULL,      CPL_ERROR_ILLEGAL_OUTPUT);
00273 
00274     /* Compute flux in PSF and find max peak */
00275     *psf_peak = cpl_image_get_max(psf);
00276     cpl_image_delete(psf);
00277 
00278     assert( *psf_peak > 0.0); /* The ideal PSF has a positive maximum */
00279     *psf_flux = 1.0; /* The psf flux, cpl_image_get_flux(psf), is always 1 */
00280 
00281     /* Measure the background in the candidate image */
00282     *star_bg = irplib_strehl_ring_background(im, xpos, ypos, r2/pscale, r3/pscale,
00283                                              IRPLIB_BG_METHOD_AVER_REJ);
00284 
00285     /* Compute star_radius in pixels */
00286     star_radius = r1/pscale;
00287 
00288     /* Measure the flux on the candidate image */
00289     *star_flux = irplib_strehl_disk_flux(im, xpos, ypos, star_radius, *star_bg);
00290 
00291     cpl_ensure_code(*star_flux > 0.0,      CPL_ERROR_ILLEGAL_OUTPUT);
00292 
00293     /* Find the peak value on the central part of the candidate image */
00294     max_radius = window_size < star_radius ? window_size : star_radius;
00295     cpl_ensure_code(!irplib_strehl_disk_max(im, xpos, ypos, max_radius,
00296                                             star_peak), cpl_error_get_code());
00297     *star_peak -= *star_bg;
00298 
00299     cpl_ensure_code(*star_peak > 0.0,      CPL_ERROR_ILLEGAL_OUTPUT);
00300 
00301     /* Compute Strehl */
00302     /* (StarPeak / StarFlux) / (PsfPeak / PsfFlux) */
00303     *strehl = (*star_peak * *psf_flux ) / ( *star_flux * *psf_peak);
00304 
00305     if (*strehl > 1)
00306         cpl_msg_warning(cpl_func, "Extreme Strehl-ratio=%g, star_peak=%g, "
00307                         "star_flux=%g, psf_peak=%g, psf_flux=%g", *strehl,
00308                         *star_peak, *star_flux, *psf_peak, *psf_flux);
00309 
00310     /* Compute Strehl error */
00311     /* computation could fail if the image contains pixels with NaN value*/
00312     ring[0] = xpos;
00313     ring[1] = ypos;
00314     ring[2] = r2/pscale;
00315     ring[3] = r3/pscale;
00316 
00317     /* FIXME: With CPL 5.1 the recoverable error
00318        will be CPL_ERROR_DATA_NOT_FOUND */
00319     prestate = cpl_errorstate_get();
00320     while (cpl_flux_get_noise_ring(im, ring, noise_box_sz, noise_nsamples,
00321                                    bg_noise, NULL) && --ring_tries > 0);
00322     if (ring_tries > 0) {
00323         cpl_errorstate_set(prestate); /* Recover, if an error happened */
00324     } else {
00325         return cpl_error_set_where(cpl_func);
00326     }
00327 
00328     *strehl_err = strehl_error_coefficient * (*bg_noise) * pscale *
00329         star_radius * star_radius / *star_flux;
00330 
00331     /* This check should not be able to fail, but just to be sure */
00332     cpl_ensure_code(*strehl_err >= 0.0,       CPL_ERROR_ILLEGAL_OUTPUT);
00333 
00334     return CPL_ERROR_NONE;
00335 }
00336 
00337 /*----------------------------------------------------------------------------*/
00353 /*----------------------------------------------------------------------------*/
00354 double irplib_strehl_disk_flux(const cpl_image * im,
00355                                double            xpos,
00356                                double            ypos,
00357                                double            rad,
00358                                double            bg)
00359 {
00360     const double    sqr = rad * rad;
00361     double          sqrest;
00362     const float *   pim;
00363     double          flux = 0.0;
00364     double          yj, xi;
00365     int             nx, ny;
00366     int             lx, ly, ux, uy;
00367     int             i, j;
00368 
00369 
00370     /* Check entries */
00371     cpl_ensure(im != NULL, CPL_ERROR_NULL_INPUT, 0.0);
00372     cpl_ensure(cpl_image_get_type(im) == CPL_TYPE_FLOAT,
00373                CPL_ERROR_UNSUPPORTED_MODE, 0.0);
00374     cpl_ensure(rad > 0.0, CPL_ERROR_ILLEGAL_INPUT, 0.0);
00375 
00376     nx = cpl_image_get_size_x(im);
00377     ny = cpl_image_get_size_y(im);
00378 
00379     /* Round down */
00380     lx = (int)(xpos - rad);
00381     ly = (int)(ypos - rad);
00382     if (lx < 0) lx = 0;
00383     if (ly < 0) ly = 0;
00384 
00385     /* Round up */
00386     ux = (int)(xpos + rad) + 1;
00387     uy = (int)(ypos + rad) + 1;
00388     if (ux > (nx-1)) ux = nx-1;
00389     if (uy > (ny-1)) uy = ny-1;
00390 
00391     pim = cpl_image_get_data_float_const(im);
00392     for (j=ly ; j<uy ; j++) {
00393         yj = (double)j - ypos;
00394         sqrest = sqr - yj * yj;
00395         for (i=lx; i<ux ; i++) {
00396             xi = (double)i - xpos;
00397             if (sqrest >= xi * xi && irplib_isnan(pim[i+j*nx]) == 0) {
00398                 flux += (double)pim[i+j*nx] - bg;
00399             }
00400         }
00401     }
00402     return flux;
00403 }
00404 
00405 /*----------------------------------------------------------------------------*/
00419 /*----------------------------------------------------------------------------*/
00420 double irplib_strehl_ring_background(const cpl_image *       im,
00421                                      double                  xpos,
00422                                      double                  ypos,
00423                                      double                  rad_int,
00424                                      double                  rad_ext,
00425                                      irplib_strehl_bg_method mode)
00426 {
00427     int             npix;
00428     const double    sqr_int = rad_int * rad_int;
00429     const double    sqr_ext = rad_ext * rad_ext;
00430     double          dist;
00431     cpl_vector  *   pix_arr;
00432     const float *   pim;
00433     double          flux = 0.0;
00434     double          yj, xi;
00435     int             lx, ly, ux, uy;
00436     int             nx, ny;
00437     int             i, j;
00438 
00439     /* Check entries */
00440     cpl_ensure(im != NULL, CPL_ERROR_NULL_INPUT, 0.0);
00441     cpl_ensure(cpl_image_get_type(im) == CPL_TYPE_FLOAT,
00442                CPL_ERROR_UNSUPPORTED_MODE, 0.0);
00443     cpl_ensure(rad_int > 0.0, CPL_ERROR_ILLEGAL_INPUT, 0.0);
00444     cpl_ensure(rad_ext > rad_int, CPL_ERROR_ILLEGAL_INPUT, 0.0);
00445 
00446     cpl_ensure(mode == IRPLIB_BG_METHOD_AVER_REJ ||
00447                mode == IRPLIB_BG_METHOD_MEDIAN,
00448                CPL_ERROR_UNSUPPORTED_MODE, 0.0);
00449 
00450     nx = cpl_image_get_size_x(im);
00451     ny = cpl_image_get_size_y(im);
00452 
00453     /* Round down */
00454     lx = (int)(xpos - rad_ext);
00455     ly = (int)(ypos - rad_ext);
00456     if (lx < 0) lx = 0;
00457     if (ly < 0) ly = 0;
00458 
00459     /* Round up */
00460     ux = (int)(xpos + rad_ext) + 1;
00461     uy = (int)(ypos + rad_ext) + 1;
00462     if (ux > (nx-1)) ux = nx-1;
00463     if (uy > (ny-1)) uy = ny-1;
00464 
00465     npix = (ux - lx + 1) * (uy - ly + 1);
00466     cpl_ensure(npix >= IRPLIB_DISK_BG_MIN_PIX_NB, CPL_ERROR_DATA_NOT_FOUND, 0.0);
00467 
00468     /* Allocate pixel array to hold values in the ring */
00469     pix_arr = cpl_vector_new(npix);
00470 
00471     /* Count number of pixels in the ring */
00472     /* Retrieve all pixels which belong to the ring */
00473     pim = cpl_image_get_data_float_const(im);
00474     npix = 0;
00475     for (j=ly ; j<uy ; j++) {
00476         yj = (double)j - ypos;
00477         for (i=lx ; i<ux; i++) {
00478             xi = (double)i - xpos;
00479             dist = yj * yj + xi * xi;
00480             if (sqr_int <= dist && dist <= sqr_ext &&
00481                 irplib_isnan(pim[i+j*nx]) == 0) {
00482                 cpl_vector_set(pix_arr, npix, (double)pim[i+j*nx]);
00483                 npix++;
00484             }
00485         }
00486     }
00487 
00488     if (npix < IRPLIB_DISK_BG_MIN_PIX_NB) {
00489         cpl_vector_delete(pix_arr);
00490         cpl_ensure(0, CPL_ERROR_DATA_NOT_FOUND, 0.0);
00491     }
00492 
00493     /* Should not be able to fail now */
00494 
00495     /* Resize pixel array to actual number of values within the ring */
00496     cpl_vector_set_size(pix_arr, npix);
00497 
00498     if (mode == IRPLIB_BG_METHOD_AVER_REJ) {
00499         const int low_ind  = (int)((double)npix * IRPLIB_DISK_BG_REJ_LOW);
00500         const int high_ind = (int)((double)npix
00501                                    * (1.0 - IRPLIB_DISK_BG_REJ_HIGH));
00502 
00503         /* Sort the array */
00504         cpl_vector_sort(pix_arr, 1);
00505 
00506         for (i=low_ind ; i<high_ind ; i++) {
00507             flux += cpl_vector_get(pix_arr, i);
00508         }
00509         if (high_ind - low_ind > 1) flux /= (double)(high_ind - low_ind);
00510     } else /* if (mode == IRPLIB_BG_METHOD_MEDIAN) */ {
00511         flux = cpl_vector_get_median(pix_arr);
00512     }
00513 
00514     cpl_vector_delete(pix_arr);
00515 
00516     return flux;
00517 }
00518 
00519 /*----------------------------------------------------------------------------*/
00539 /*----------------------------------------------------------------------------*/
00540 cpl_image * irplib_strehl_generate_psf(
00541         double   m1,
00542         double   m2,
00543         double   lam,
00544         double   dlam,
00545         double   pscale,
00546         int      size)
00547 {
00548     cpl_image * otf_image = irplib_strehl_generate_otf(m1, m2, lam, dlam,
00549             size, pscale);
00550 
00551     if (otf_image == NULL) return NULL;
00552 
00553     /* Transform back to real space
00554        - Normalization is unnecessary, due to the subsequent normalisation.
00555        - An OTF is point symmetric about its center, i.e. it is even,
00556          i.e. the real space image is real.
00557        - Because of this a forward FFT works as well.
00558        - If the PSF ever needs to have its images halves swapped add
00559          CPL_FFT_SWAP_HALVES to the FFT call.
00560      */
00561 
00562     if (cpl_image_fft(otf_image, NULL, CPL_FFT_UNNORMALIZED) ||
00563 
00564         /* Compute absolute values of PSF */
00565         cpl_image_abs(otf_image) ||
00566 
00567         /* Normalize PSF to get flux=1  */
00568         cpl_image_normalise(otf_image, CPL_NORM_FLUX)) {
00569 
00570         cpl_image_delete(otf_image);
00571         return NULL;
00572     }
00573 
00574     return otf_image;
00575 }
00576 
00579 /*----------------------------------------------------------------------------*/
00595 /*----------------------------------------------------------------------------*/
00596 static cpl_image * irplib_strehl_generate_otf(
00597         double  m1,
00598         double  m2,
00599         double  lam,
00600         double  dlam,
00601         int     size,
00602         double  pscale)
00603 {
00604     cpl_image   *   otf_image;
00605     double      *   otf_data;
00606     double          obs_ratio ;  /* m1 / m2    */
00607     double          f_max ;      /* cut-off frequency        */
00608     int             pix0 ;      /* Pixel corresponding to the zero frequency */
00609     double          a, x, y;
00610     double          f, rsq, fc, invfc, lambda;
00611     double          sincy;
00612     double          invsize;
00613     register int    pos;
00614     int             i, j, k;
00615 
00616 
00617     cpl_ensure(m2   > 0.0,      CPL_ERROR_ILLEGAL_INPUT, NULL);
00618     cpl_ensure(m1   > m2,       CPL_ERROR_ILLEGAL_INPUT, NULL);
00619     cpl_ensure(lam  > 0.0,      CPL_ERROR_ILLEGAL_INPUT, NULL);
00620     cpl_ensure(dlam > 0.0,      CPL_ERROR_ILLEGAL_INPUT, NULL);
00621     cpl_ensure(size > 0,        CPL_ERROR_ILLEGAL_INPUT, NULL);
00622     cpl_ensure(pscale > 0.0,    CPL_ERROR_ILLEGAL_INPUT, NULL);
00623 
00624     /* Convert pixel scale from sec to radians, microns in meters    */
00625     pscale /= (double)206265;
00626     lam /= (double)1.0e6;
00627     dlam /= (double)1.0e6;
00628 
00629     /* Obscuration ratio    */
00630     obs_ratio = m2 / m1;
00631 
00632     /* Pixel corresponding to the zero frequency    */
00633     pix0 = size/2;
00634     invsize = (double)1.0 / (double)size;
00635 
00636     /* Cut-off frequency in pixels  */
00637     f_max = m1 * pscale * (double)size / lam;
00638 
00639     /* Allocate for output image    */
00640     otf_image = cpl_image_new(size, size, CPL_TYPE_DOUBLE);
00641     if (otf_image==NULL) return NULL;
00642     otf_data = cpl_image_get_data_double(otf_image);
00643 
00644     /* Now compute the OTF  */
00645     /* OPTIMIZED CODE !!! LIMITED READABILITY !!!   */
00646 
00647     for (k=1 ; k<=9 ; k++) {    /* iteration on the wavelength  */
00648         /* Compute intermediate cut-off frequency   */
00649         lambda = (double)(lam - dlam*(double)(k-5)/8.0);
00650         fc = (double)f_max * (double)lam / lambda;
00651         invfc = 1.0 / fc;
00652 
00653         /* Convolution with the detector pixels */
00654         pos = 0;
00655         for (j=0 ; j<size ; j++) {
00656             y = (double)(j-pix0);
00657             sincy = PSF_sinc(CPL_MATH_PI * y * invsize);
00658             for (i=0 ; i<size ; i++) {
00659                 x = (double)(i-pix0);
00660                 rsq = x*x + y*y;
00661                 if (rsq < fc*fc) {
00662                     if (rsq < 0.01)
00663                         a = 1.0;
00664                     else {
00665                         f = sqrt(rsq) * invfc;
00666                         a = PSF_TelOTF(f,obs_ratio) *
00667                             PSF_sinc(CPL_MATH_PI * x * invsize) * sincy;
00668                     }
00669                 } else {
00670                     a = 0.0;
00671                 }
00672                 otf_data[pos++] += a / 9.0;
00673             }
00674         }
00675     }
00676     return otf_image;
00677 }
00678 
00679 /*----------------------------------------------------------------------------*
00680  * H1 function
00681  *----------------------------------------------------------------------------*/
00682 static double PSF_H1(
00683         double  f,
00684         double  u,
00685         double  v)
00686 {
00687     const double e = fabs(1.0-v) > 0.0 ? -1.0 : 1.0; /* e = 1.0 iff v = 1.0 */
00688 
00689     return((v*v/CPL_MATH_PI)*acos((f/v)*(1.0+e*(1.0-u*u)/(4.0*f*f))));
00690 }
00691 
00692 /*----------------------------------------------------------------------------*
00693  * H2 function
00694  *----------------------------------------------------------------------------*/
00695 static double PSF_H2(double  f,
00696                      double  u)
00697 {
00698     const double tmp1 = (2.0 * f) / (1.0 + u);
00699     const double tmp2 = (1.0 - u) / (2.0 * f);
00700 
00701     return -1.0 * (f/CPL_MATH_PI) * (1.0+u)
00702         * sqrt((1.0-tmp1*tmp1)*(1.0-tmp2*tmp2));
00703 }
00704 
00705 /*----------------------------------------------------------------------------*
00706  * G function
00707  *----------------------------------------------------------------------------*/
00708 static double PSF_G(double  f,
00709                     double  u)
00710 {
00711     if (f <= (1.0-u)/2.0) return(u*u);
00712     if (f >= (1.0+u)/2.0) return(0.0);
00713     else return(PSF_H1(f,u,1.0) + PSF_H1(f,u,u) + PSF_H2(f,u));
00714 }
00715 
00716 /*----------------------------------------------------------------------------*
00717  * sinc function
00718  *----------------------------------------------------------------------------*/
00719 static double PSF_sinc(double x)
00720 {
00721   return fabs(x) > fabs(sin(x)) ? sin(x)/x : 1.0;
00722 }
00723 
00724 /*----------------------------------------------------------------------------*
00725  * Telescope OTF function
00726  *----------------------------------------------------------------------------*/
00727 static double PSF_TelOTF(double  f,
00728                          double  u)
00729 {
00730     return((PSF_G(f,1.0)+u*u*PSF_G(f/u,1.0)-2.0*PSF_G(f,u))/(1.0-u*u));
00731 }
00732 
00733 /*----------------------------------------------------------------------------*/
00745 /*----------------------------------------------------------------------------*/
00746 cpl_error_code irplib_strehl_disk_max(const cpl_image * self,
00747                                              double            xpos,
00748                                              double            ypos,
00749                                              double            radius,
00750                                              double          * ppeak)
00751 {
00752 
00753     const double    sqr = radius * radius;
00754     double          sqrest;
00755     const float *   pself;
00756     float           peak = FLT_MAX;  /* Avoid (false) uninit warning */
00757     double          yj, xi;
00758     int             nx, ny;
00759     int             lx, ly, ux, uy;
00760     int             i, j;
00761     cpl_boolean     first = CPL_TRUE;
00762 
00763 
00764     /* Check entries */
00765     cpl_ensure_code(ppeak != NULL, CPL_ERROR_NULL_INPUT);
00766     cpl_ensure_code(self  != NULL, CPL_ERROR_NULL_INPUT);
00767     cpl_ensure_code(cpl_image_get_type(self) == CPL_TYPE_FLOAT,
00768                     CPL_ERROR_UNSUPPORTED_MODE);
00769     cpl_ensure_code(radius > 0.0, CPL_ERROR_ILLEGAL_INPUT);
00770 
00771     nx = cpl_image_get_size_x(self);
00772     ny = cpl_image_get_size_y(self);
00773 
00774     /* Round down */
00775     lx = (int)(xpos - radius);
00776     ly = (int)(ypos - radius);
00777     if (lx < 0) lx = 0;
00778     if (ly < 0) ly = 0;
00779 
00780     /* Round up */
00781     ux = (int)(xpos + radius) + 1;
00782     uy = (int)(ypos + radius) + 1;
00783     if (ux > (nx-1)) ux = nx-1;
00784     if (uy > (ny-1)) uy = ny-1;
00785 
00786     pself = cpl_image_get_data_float_const(self);
00787     for (j=ly ; j<uy ; j++) {
00788         yj = (double)j - ypos;
00789         sqrest = sqr - yj * yj;
00790         for (i=lx; i<ux ; i++) {
00791             xi = (double)i - xpos;
00792             if (sqrest >= xi * xi && irplib_isnan(pself[i+j*nx]) == 0) {
00793                 if (first || pself[i+j*nx] > peak) {
00794                     first = CPL_FALSE;
00795                     peak = pself[i+j*nx];
00796                 }
00797             }
00798         }
00799     }
00800 
00801     cpl_ensure_code(!first, CPL_ERROR_DATA_NOT_FOUND);
00802 
00803     *ppeak = (double)peak;
00804 
00805     return CPL_ERROR_NONE;
00806 }

Generated on Wed Mar 9 15:46:17 2011 for NACO Pipeline Reference Manual by  doxygen 1.5.8