sinfo_focus.c

00001 /*
00002  * This file is part of the ESO SINFONI Pipeline
00003  * Copyright (C) 2004,2005 European Southern Observatory
00004  *
00005  * This program 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 * E.S.O. - VLT project
00021 *
00022 *
00023 *
00024 * who       when      what
00025 * --------  --------  ----------------------------------------------
00026 * schreib  16/01/02  created
00027 */
00028 
00029 /************************************************************************
00030 *   NAME
00031 *        sinfo_focus.c -
00032 *        routines to determine the focus position of the detector
00033 *
00034 *   SYNOPSIS
00035 *   #include "sinfo_focus.h"
00036 *
00037 *   1) double sinfo_new_gaussian_ellipse ( double * xdat, double * parlist )
00038 *   2) void sinfo_new_gaussian_ellipse_deriv( double * xdat, 
00039                                               double * parlist, 
00040                                               double * dervs )
00041 *   3) static int new_inv_mat (void)
00042 *   4) static void new_get_mat ( double * xdat,
00043 *                            int    * xdim,
00044 *                            double * ydat,
00045 *                            double * wdat,
00046 *                            int    * ndat,
00047 *                            double * fpar,
00048 *                            double * epar,
00049 *                            int    * npar )
00050 *   5) static int new_get_vec ( double * xdat,
00051 *                           int    * xdim,
00052 *                           double * ydat,
00053 *                           double * wdat,
00054 *                           int    * ndat,
00055 *                           double * fpar,
00056 *                           double * epar,
00057 *                           int    * npar )
00058 *   6) int new_lsqfit ( double * xdat,
00059 *                   int    * xdim,
00060 *                   double * ydat,
00061 *                   double * wdat,
00062 *                   int    * ndat,
00063 *                   double * fpar,
00064 *                   double * epar,
00065 *                   int    * mpar,
00066 *                   int    * npar,
00067 *                   double * tol ,
00068 *                   int    * its ,
00069 *                   double * lab  )
00070 *   7) int sinfo_new_fit_2d_gaussian( cpl_image   * lineImage, 
00071 *                         double     * fit_par, 
00072 *                         double     * derv_par   
00073 *                         int        * mpar,
00074 *                         int          lleftx,
00075 *                         int          llefty,
00076 *                         int          halfbox_x,
00077 *                         int          halfbox_y, int* check )
00078 *   8) cpl_image * sinfo_new_plot_gaussian ( cpl_image   * image, 
00079 *                                double     * parlist )
00080 *   9) static int new_gauss2ellipse ( double     * parlist ,
00081 *  10) float sinfo_new_determine_conversion_factor ( cpl_imagelist * cube, 
00082 *                                        float     mag,
00083 *                                        float     exptime,
00084 *                                        int       lleftx,
00085 *                                        int       llefty,
00086 *                                        int       halfbox_x,
00087 *                                        int       halfbox_y, 
00088 *                                        int* check )
00089 *
00090 *   DESCRIPTION
00091 *   1) Compute the value of a 2d Gaussian function at a given point.
00092 *      The ellptical 2D Gaussian is:
00093 *      F(x,y) = par(2) * EXP( -4.0*log(2.0)*[(xr/par(4))^2+(yr/par(5))^2]) + 
00094                 par(3),
00095 *      where: xr = xo * cos(par(6)) + yo * sin(par(6))
00096 *      yr = -xo * sin(par(6)) + yo * cos(par(6))
00097 *      and:   x0 = x - par(0)
00098 *             y0 = y - par(1)
00099 *   2) calculates the partial derivatives for a 2d Gaussian function with
00100 *      parameters parlist at position xdat 
00101 *   3) calculates the inverse of matrix2. The algorithm used 
00102 *      is the Gauss-Jordan algorithm described in Stoer,
00103 *      Numerische Mathematik, 1. Teil.
00104 *   4) builds the sinfo_matrix 
00105 *   5) calculates the correction sinfo_vector. The sinfo_matrix has been
00106 *      built by get_mat(), we only have to rescale it for the 
00107 *      current value of labda. The sinfo_matrix is rescaled so that
00108 *      the diagonal gets the value 1 + labda.
00109 *      Next we calculate the inverse of the sinfo_matrix and then
00110 *      the correction sinfo_vector.
00111 *   6) this is a routine for making a least-squares fit of a
00112 *      function to a set of data points. The method used is
00113 *      described in: Marquardt, J.Soc.Ind.Appl.Math. 11. 431 (1963).
00114 *      This method is a mixture of the steepest descent method 
00115 *      and the Taylor method.
00116 *   7) fits the image of a point source by using a 2-D Gaussian
00117 *      fit.
00118 *   8) plots an image of a given 2D-Gaussian 
00119 *   9) converts gauss parameters to ellipse parameters. 
00120 *  10) determines an intensity conversion factor for the instrument
00121 *      by fitting a 2D-Gaussian to an collapsed image of a standard star
00122 *      with known brightness (only for non-AO observations).
00123 *      Then the resulting Gaussian is integrated and the counts
00124 *      are divided by the exposure time (Fits header information) 
00125 *
00126 *   FILES
00127 *
00128 *   ENVIRONMENT
00129 *
00130 *   RETURN VALUES
00131 *
00132 *   CAUTIONS
00133 *
00134 *   EXAMPLES
00135 *
00136 *   SEE ALSO
00137 *
00138 *   BUGS
00139 *
00140 *------------------------------------------------------------------------
00141 */
00142 
00143 #ifdef HAVE_CONFIG_H
00144 #  include <config.h>
00145 #endif
00146 #include "sinfo_vltPort.h"
00147 
00148 /*
00149  * System Headers
00150  */
00151 
00152 /*
00153  * Local Headers
00154  */
00155 
00156 #include "sinfo_focus.h"
00157 #include "sinfo_recipes.h"
00158 #include <float.h>
00159 
00160 /*----------------------------------------------------------------------------
00161  *                                 Defines
00162  *--------------------------------------------------------------------------*/
00163 
00164 #define XDIMG          2         /* dimension of the x values */
00165 #define TOLG           0.001     /* fitting tolerance */
00166 #define LABG           0.1       /* labda parameter */
00167 #define ITSG           200       /* maximum number of iterations */
00168 #define LABFACG        10.0      /* labda step factor */
00169 #define LABMAXG        1.0e+10   /* maximum value for labda */
00170 #define LABMING        1.0e-10   /* minimum value for labda */
00171 #define NPAR           7         /* number of fit parameters */
00172 #define PI_NUMB        (3.1415926535897932384626433832795) /* pi */
00173 
00174 
00175 /*----------------------------------------------------------------------------
00176  *                                    Local variables
00177  *--------------------------------------------------------------------------*/
00178 
00179 static double chi1 ;                    /* old reduced chi-squared */
00180 static double chi2 ;                    /* new reduced chi-squared */
00181 static double labda ;                   /* mixing parameter */
00182 static double vec[NPAR] ;               /* correction sinfo_vector */
00183 static double matrix1[NPAR][NPAR] ;     /* original sinfo_matrix */
00184 static double matrix2[NPAR][NPAR] ;     /* inverse of matrix1 */
00185 static int    nfree ;                   /* number of free parameters */
00186 static int    parptr[NPAR] ;            /* parameter pointer */
00187 
00188 /*----------------------------------------------------------------------------
00189  *                    Functions private to this module
00190  *--------------------------------------------------------------------------*/
00191 
00192 static int new_inv_mat (void) ;
00193 
00194 static void new_get_mat ( double * xdat,
00195                       int    * xdim,
00196                       double * ydat,
00197                       double * wdat,
00198                       int    * ndat,
00199                       double * fpar,
00200                       double * epar/*,
00201                       int    * npar */) ;
00202 
00203 static int new_get_vec ( double * xdat,
00204                      int    * xdim,
00205                      double * ydat,
00206                      double * wdat,
00207                      int    * ndat,
00208                      double * fpar,
00209                      double * epar,
00210                      int    * npar ) ;
00211 
00212 static int new_gauss2Ellipse ( double  * parlist ) ;
00221 /*----------------------------------------------------------------------------
00222  *                            Function codes
00223  *--------------------------------------------------------------------------*/
00224 
00225 /*-------------------------------------------------------------------------*/
00249 /*--------------------------------------------------------------------------*/
00250 
00251 double sinfo_new_gaussian_ellipse(double * xdat, double * parlist)
00252 {
00253     double  result ;
00254     double x ;
00255     double y ;
00256     double fwhmx ;
00257     double fwhmy ;
00258     double costheta ;
00259     double sintheta ;
00260     double argX ;           /* arguments in the exponent */
00261     double argY ;
00262 
00263     /* some abbreviations */
00264     x =  xdat[0] -  parlist[0] ;
00265     y =  xdat[1] -  parlist[1] ;
00266 
00267     fwhmx = fabs(parlist[4]) ;
00268     fwhmy = fabs(parlist[5]) ;
00269     
00270     costheta = cos ( parlist[6] ) ;
00271     sintheta = sin ( parlist[6] ) ;
00272 
00273     argX = x * costheta + y * sintheta ;
00274     argY = -x * sintheta + y * costheta ;
00275     
00276     /* function */
00277     result =  parlist[2] * exp(-4.*log(2.0)*((argX/fwhmx)*(argX/fwhmx)+
00278                                              (argY/fwhmy)*(argY/fwhmy))) +
00279               parlist[3] ; 
00280 
00281     return result ;
00282 }
00283 
00308 void 
00309 sinfo_new_gaussian_ellipse_deriv(double * xdat, 
00310                                  double * parlist, 
00311                                  double * dervs )
00312 {
00313     double x ;
00314     double y ;
00315     double fwhmx ;
00316     double fwhmy ;
00317     double argX ;
00318     double argY ;
00319     double expon ;
00320     double e8log2 ;
00321     double fwx2 ;
00322     double fwy2 ;
00323     double costheta ;
00324     double sintheta ;
00325 
00326     /* some abbreviations */
00327     x = xdat[0] - parlist[0] ;
00328     y = xdat[1] - parlist[1] ;
00329 
00330     fwhmx = fabs(parlist[4]) ;
00331     fwhmy = fabs(parlist[5]) ;
00332     fwx2 = fwhmx * fwhmx ;
00333     fwy2 = fwhmy * fwhmy ;
00334 
00335     costheta = cos ( parlist[6] ) ;
00336     sintheta = sin ( parlist[6] ) ;
00337 
00338     argX = x * costheta + y * sintheta ;
00339     argY = -x * sintheta + y * costheta ;
00340 
00341     expon = exp ( -4.0 * log(2.0) * ((argX/fwhmx)*(argX/fwhmx) + 
00342                                      (argY/fwhmy)*(argY/fwhmy)) ) ;
00343     e8log2 = expon * 8.0 * log(2.0) ;
00344 
00345     /* determine the derivatives */
00346     /* partial derivative x-position */
00347     dervs[0] = -parlist[2]*e8log2 * (-argX*costheta/fwx2 + argY*sintheta/fwy2);
00348     /* partial derivative y-position */
00349     dervs[1] = -parlist[2]*e8log2 * (-argX*sintheta/fwx2 - argY*costheta/fwy2);
00350     /* partial derivative amplitude */
00351     dervs[2] = expon ;
00352     /* partial derivative background */
00353     dervs[3] = 1. ;
00354     /* partial derivative fwhmx */
00355     dervs[4] = parlist[2]*e8log2 * argX*argX/(fwx2*fwhmx) ;
00356     /* partial derivative fwhmy */
00357     dervs[5] = parlist[2]*e8log2 * argY*argY/(fwy2*fwhmy) ;
00358     /* partial derivative theta */
00359     dervs[6] = -parlist[2]*e8log2 * argY * argX * (1.0/fwx2 - 1.0/fwy2) ;
00360     
00361 }
00362 
00373 static int new_inv_mat (void)
00374 {
00375     double even ;
00376     double hv[NPAR] ;
00377     double mjk ;
00378     double rowmax ;
00379     int evin ;
00380     int i, j, k, row ;
00381     int per[NPAR] ;
00382    
00383     /* set permutation array */
00384     for ( i = 0 ; i < nfree ; i++ )
00385     {
00386         per[i] = i ;
00387     }
00388     
00389     for ( j = 0 ; j < nfree ; j++ ) /* in j-th column */
00390     {
00391         /* determine largest element of a row */                                
00392         rowmax = fabs ( matrix2[j][j] ) ;
00393         row = j ;                         
00394 
00395         for ( i = j + 1 ; i < nfree ; i++ )
00396         {
00397             if ( fabs ( matrix2[i][j] ) > rowmax )
00398             {
00399                 rowmax = fabs( matrix2[i][j] ) ;
00400                 row = i ;
00401             }
00402         }
00403 
00404         /* determinant is zero! */
00405         if ( matrix2[row][j] == 0.0 )
00406         {
00407             return -6 ;
00408         }
00409         
00410         /* if the largest element is not on the diagonal, 
00411            then permutate rows */
00412         if ( row > j )
00413         {
00414             for ( k = 0 ; k < nfree ; k++ )
00415             {
00416                 even = matrix2[j][k] ;
00417                 matrix2[j][k] = matrix2[row][k] ;
00418                 matrix2[row][k] = even ;
00419             }
00420             /* keep track of permutation */
00421             evin = per[j] ;
00422             per[j] = per[row] ;
00423             per[row] = evin ;
00424         }
00425         
00426         /* modify column */
00427         even = 1.0 / matrix2[j][j] ;
00428         for ( i = 0 ; i < nfree ; i++ )
00429         {
00430             matrix2[i][j] *= even ;
00431         }
00432         matrix2[j][j] = even ;
00433         
00434         for ( k = 0 ; k < j ; k++ )
00435         {
00436             mjk = matrix2[j][k] ;
00437             for ( i = 0 ; i < j ; i++ )
00438             {
00439                 matrix2[i][k] -= matrix2[i][j] * mjk ;
00440             }
00441             for ( i = j + 1 ; i < nfree ; i++ )
00442             {
00443                 matrix2[i][k] -= matrix2[i][j] * mjk ;
00444             }
00445             matrix2[j][k] = -even * mjk ;
00446         }
00447     
00448         for ( k = j + 1 ; k < nfree ; k++ )
00449         {
00450             mjk = matrix2[j][k] ;
00451             for ( i = 0 ; i < j ; i++ )
00452             {
00453                 matrix2[i][k]  -= matrix2[i][j] * mjk ;
00454             }
00455             for ( i = j + 1 ; i < nfree ; i++ )
00456             {
00457                 matrix2[i][k]  -= matrix2[i][j] * mjk ;
00458             }
00459             matrix2[j][k] = -even * mjk ;
00460         }
00461     }
00462     
00463     /* finally, repermute the columns */
00464     for ( i = 0 ; i < nfree ; i++ )
00465     {
00466         for ( k = 0 ; k < nfree ; k++ )
00467         {
00468             hv[per[k]] = matrix2[i][k] ;
00469         }
00470         for ( k = 0 ; k < nfree ; k++ )
00471         {
00472             matrix2[i][k] = hv[k] ;
00473         }
00474     }
00475         
00476     /* all is well */
00477     return 0 ;
00478 }
00479     
00495 static void new_get_mat ( double * xdat,
00496                       int    * xdim,
00497                       double * ydat,
00498                       double * wdat,
00499                       int    * ndat,
00500                       double * fpar,
00501                       double * epar/*,
00502                       int    * npar */)
00503 {
00504     double wd ;
00505     double wn ;
00506     double yd ;
00507     int i, j, n ;
00508 
00509     for ( j = 0 ; j < nfree ; j++ )
00510     {
00511         vec[j] = 0.0 ; /* zero sinfo_vector */
00512         for ( i = 0 ; i<= j ; i++ )   /* zero sinfo_matrix only on 
00513                                          and below diagonal */
00514         {
00515             matrix1[j][i] = 0.0 ;
00516         }
00517     }
00518     chi2 = 0.0 ;  /* reset reduced chi-squared */
00519     
00520     /* loop through data points */
00521     for ( n = 0 ; n < (*ndat) ; n++ )
00522     {
00523         wn = wdat[n] ;
00524         if ( wn > 0.0 )  /* legal weight ? */
00525         {
00526             yd=ydat[n] - sinfo_new_gaussian_ellipse(&xdat[(*xdim) * n],fpar) ;
00527             sinfo_new_gaussian_ellipse_deriv( &xdat[(*xdim) * n], fpar, epar ) ;
00528             chi2 += yd * yd * wn ; /* add to chi-squared */
00529             for ( j = 0 ; j < nfree ; j++ )
00530             {
00531                 wd = epar[parptr[j]] * wn ;  /* weighted derivative */
00532                 vec[j] += yd * wd ;       /* fill sinfo_vector */
00533                 for ( i = 0 ; i <= j ; i++ ) /* fill sinfo_matrix */
00534                 {
00535                     matrix1[j][i] += epar[parptr[i]] * wd ;
00536                 }
00537             }
00538         }
00539     }                   
00540 }  
00541                 
00542             
00568 static int new_get_vec ( double * xdat,
00569                      int    * xdim,
00570                      double * ydat,
00571                      double * wdat,
00572                      int    * ndat,
00573                      double * fpar,
00574                      double * epar,
00575                      int    * npar )
00576 {
00577     double dj ;
00578     double dy ;
00579     double mii ;
00580     double mji ;
00581     double mjj ;
00582     double wn ;
00583     int i, j, n, r ;
00584 
00585     /* loop to modify and scale the sinfo_matrix */
00586     for ( j = 0 ; j < nfree ; j++ )
00587     {
00588         mjj = matrix1[j][j] ;
00589         if ( mjj <= 0.0 )             /* diagonal element wrong */
00590         {
00591             return -5 ;
00592         }
00593         mjj = sqrt( mjj ) ;
00594         for ( i = 0 ; i < j ; i++ )
00595         {
00596             mji = matrix1[j][i] / mjj / sqrt( matrix1[i][i] ) ;
00597             matrix2[i][j] = matrix2[j][i] = mji ;
00598         }
00599         matrix2[j][j] = 1.0 + labda ;  /* scaled value on diagonal */
00600     }    
00601     
00602     if ( (r = new_inv_mat()) ) /* sinfo_invert sinfo_matrix inlace */
00603     {
00604         return r ;
00605     }
00606     
00607     for ( i = 0 ; i < (*npar) ; i ++ )
00608     {
00609         epar[i] = fpar[i] ;
00610     }
00611     
00612     /* loop to calculate correction sinfo_vector */
00613     for ( j = 0 ; j < nfree ; j++ )
00614     {
00615         dj = 0.0 ;
00616         mjj = matrix1[j][j] ;
00617         if ( mjj <= 0.0)               /* not allowed */
00618         {
00619             return -7 ;
00620         }
00621         mjj = sqrt ( mjj ) ;
00622         for ( i = 0 ; i < nfree ; i++ )
00623         {
00624             mii = matrix1[i][i] ;
00625             if ( mii <= 0.0 )
00626             {
00627                 return -7 ;
00628             }
00629             mii = sqrt( mii ) ;
00630             dj += vec[i] * matrix2[j][i] / mjj / mii ;
00631         }
00632         epar[parptr[j]] += dj ;       /* new parameters */
00633     }    
00634     chi1 = 0.0 ;                      /* reset reduced chi-squared */
00635  
00636     /* loop through the data points */
00637     for ( n = 0 ; n < (*ndat) ; n++ )
00638     {
00639         wn = wdat[n] ;               /* get weight */
00640         if ( wn > 0.0 )              /* legal weight */
00641         {
00642             dy=ydat[n] - sinfo_new_gaussian_ellipse(&xdat[(*xdim) * n],epar);
00643             chi1 += wdat[n] * dy * dy ;
00644         }
00645     }
00646     return 0 ;
00647 }   
00648     
00649         
00650 
00698 int sinfo_new_lsqfitd ( double * xdat,
00699               int    * xdim,
00700               double * ydat,
00701               double * wdat,
00702               int    * ndat,
00703               double * fpar,
00704               double * epar,
00705               int    * mpar,
00706               int    * npar,
00707               double * tol ,
00708               int    * its ,
00709               double * lab  )
00710 {
00711     int i, n, r ;
00712     int itc ;                      /* fate of fit */
00713     int found ;                    /* fit converged: 1, not yet converged: 0 */
00714     int  nuse ;                    /* number of useable data points */
00715     double tolerance ;             /* accuracy */
00716 
00717     itc   = 0 ;                    /* fate of fit */
00718     found = 0 ;                    /* reset */
00719     nfree = 0 ;                    /* number of free parameters */
00720     nuse  = 0 ;                    /* number of legal data points */
00721 
00722     if ( *tol < (DBL_EPSILON * 10.0 ) )
00723     {
00724         tolerance = DBL_EPSILON * 10.0 ;  /* default tolerance */
00725     }
00726     else
00727     {
00728         tolerance = *tol ;                /* tolerance */
00729     }
00730     
00731     labda = fabs( *lab ) * LABFACG ; /* start value for mixing parameter */
00732     for ( i = 0 ; i < (*npar) ; i++ )
00733     {
00734         if ( mpar[i] )
00735         {
00736             if ( nfree > NPAR )         /* too many free parameters */
00737             {
00738                 return -1 ;
00739             }
00740             parptr[nfree++] = i ;         /* a free parameter */
00741         }
00742     }
00743     
00744     if (nfree == 0)                       /* no free parameters */     
00745     {
00746         return -2 ;
00747     }
00748     
00749     for ( n = 0 ; n < (*ndat) ; n++ )
00750     {
00751         if ( wdat[n] > 0.0 )              /* legal weight */
00752         {
00753             nuse ++ ;
00754         }
00755     }
00756     
00757     if ( nfree >= nuse )
00758     {
00759         return -3 ;                       /* no degrees of freedom */
00760     }
00761     if ( labda == 0.0 )                   /* linear fit */
00762     {
00763         /* initialize fpar array */
00764         for ( i = 0 ; i < nfree ; fpar[parptr[i++]] = 0.0 ) ;  
00765         new_get_mat ( xdat, xdim, ydat, wdat, ndat, fpar, epar/*, npar */) ;
00766         r =  new_get_vec ( xdat, xdim, ydat, wdat, ndat, fpar, epar, npar ) ;
00767         if ( r )                         /* error */
00768         {
00769             return r ;
00770         }
00771         for ( i = 0 ; i < (*npar) ; i++ )
00772         {
00773             fpar[i] = epar[i] ;           /* save new parameters */
00774             epar[i] = 0.0 ;               /* and set errors to zero */
00775         }
00776         chi1 = sqrt( chi1 / (double) (nuse - nfree) ) ;
00777         for ( i = 0 ; i < nfree ; i++ )
00778         {
00779             if ( (matrix1[i][i] <= 0.0 ) || (matrix2[i][i] <= 0.0) )
00780             {
00781                 return -7 ;
00782             }
00783             epar[parptr[i]] = chi1 * sqrt( matrix2[i][i] ) / 
00784                                      sqrt( matrix1[i][i] ) ;
00785         }
00786     }
00787     else                                  /* non-linear fit */
00788     {
00789         /*----------------------------------------------------------------
00790          * the non-linear fit uses the steepest descent method in combination
00791          * with the Taylor method. The mixing of these methods is controlled
00792          * by labda. In the outer loop ( called the iteration loop ) we build
00793          * the matrix and calculate the correction sinfo_vector. In the 
00794            inner loop
00795          * (called the interpolation loop) we check whether we have obtained a
00796          * better solution than the previous one. If so, we leave the inner loop
00797          * else we increase lambda ( give more weight to the steepest descent 
00798          * method) calculate the correction vector and check again. After the 
00799          * inner loop
00800          * we do a final check on the goodness of the fit and if this satisfies
00801          * the tolerance we calculate the errors of the fitted parameters.
00802          */
00803         while ( !found )                  /* iteration loop */
00804         {      
00805             if ( itc++ == (*its) )        /* increase iteration counter */
00806             {
00807                 return -4 ;               
00808             }
00809             new_get_mat( xdat, xdim, ydat, wdat, ndat, fpar, epar/*, npar*/ ) ;
00810             
00811             /*-------------------------------------------------------------
00812              * here we decrease labda since we may assume that each iteration
00813              * brings us closer to the answer.
00814              */
00815             if ( labda > LABMING )
00816             {
00817                 labda = labda / LABFACG ;         /* decrease labda */
00818             }
00819             r = new_get_vec ( xdat, xdim, ydat, wdat, ndat, fpar, epar, npar ) ;
00820 
00821             if ( r )                      /* error */
00822             {
00823                 return r ;
00824             }
00825 
00826             while ( chi1 >= chi2 )        /* interpolation loop */
00827             {
00828                 /*-----------------------------------------------------------
00829                  * The next statement is based on experience, not on the 
00830                  * mathematics of the problem. It is assumed that we have 
00831                  * reached convergence when the pure steepest descent method 
00832                  * does not produce a better solution.
00833                  */
00834                 if ( labda > LABMAXG )    /* assume solution found */
00835                 {
00836                     break ;
00837                 }
00838                 labda = labda * LABFACG ; /* increase mixing parameter */
00839                 r = new_get_vec(xdat,xdim,ydat,wdat,ndat,fpar,epar,npar) ;
00840 
00841                 if ( r )                  /* error */
00842                 {
00843                     return r ;
00844                 }
00845             }
00846 
00847             if ( labda <= LABMAXG )        /* save old parameters */
00848             {
00849                 for ( i = 0 ; i < *npar ; i++ )
00850                 {
00851                     fpar[i] = epar[i] ;
00852                 }
00853             }
00854             if ( (fabs( chi2 - chi1 ) <= (tolerance * chi1)) || 
00855                       (labda > LABMAXG) )
00856             {
00857                 /*-----------------------------------------------------------
00858                  * we have a satisfying solution, so now we need to calculate 
00859                  * the correct errors of the fitted parameters. This we do by 
00860                  * using the pure Taylor
00861                  * method because we are very close to the real solution.
00862                  */
00863                 labda = LABMING ;              /* for Taylor solution */
00864                 new_get_mat(xdat,xdim,ydat,wdat,ndat,fpar,epar/*, npar */) ;
00865                 r=new_get_vec(xdat,xdim,ydat,wdat,ndat,fpar,epar,npar ) ;
00866 
00867                 if ( r )                    /* error */
00868                 {
00869                     return r ;
00870                 }
00871                 for ( i = 0 ; i < (*npar) ; i++ )
00872                 {
00873                     epar[i] = 0.0 ;          /* set error to zero */
00874                 }
00875                 chi2 = sqrt ( chi2 / (double) (nuse - nfree) ) ;
00876 
00877                 for ( i = 0 ; i < nfree ; i++ )
00878                 {
00879                     if ( (matrix1[i][i] <= 0.0) || (matrix2[i][i] <= 0.0) )
00880                     {
00881                         return -7 ;
00882                     }
00883                     epar[parptr[i]] = chi2 * sqrt( matrix2[i][i] ) / 
00884                                              sqrt( matrix1[i][i] ) ;
00885                 }
00886                 found = 1 ;                  /* we found a solution */
00887             }
00888         }
00889     }
00890     return itc ;                             /* return number of iterations */
00891 }
00892 
00922 int 
00923 sinfo_new_fit_2d_gaussian ( cpl_image   * image, 
00924                     double     * fit_par, 
00925                     double     * derv_par,   
00926                     int        * mpar,
00927                     int          lleftx,
00928                     int          llefty,
00929                     int          halfbox_x,
00930                     int          halfbox_y, 
00931                     int* check )
00932 {
00933     int i, j, n ;
00934     int col, row ;
00935     int boxi, boxj ;
00936     int iters ;
00937     int ndata ;
00938     int xdim ;
00939     int npar ;
00940     int its ;
00941     double lab ;
00942     double tol ;
00943     double maxval ;
00944     double background ;
00945     double amplitude ;
00946     float * backarray=NULL ;
00947     double M, Mx, My ;
00948     double Mxx, Mxy, Myy ; 
00949     double X0, Y0 ;
00950     double xydat[4 *halfbox_x*halfbox_y][XDIMG] ;
00951     double zdat[4*halfbox_x*halfbox_y] ;
00952     double wdat[4*halfbox_x*halfbox_y] ;
00953     double xco, yco ;
00954     double value ;
00955     double denom ;
00956     double temp ;
00957     int llx, lly ;
00958     int foundrow ;
00959     int foundcol ;
00960     int k ;
00961     int ilx=0;
00962     int ily=0;
00963     float* pidata=NULL;
00964 
00965 
00966     if ( NULL == image )
00967     {
00968         sinfo_msg_error("no image given") ;
00969         return -1 ;
00970     }
00971     ilx=cpl_image_get_size_x(image);
00972     ily=cpl_image_get_size_y(image);
00973 
00974     if ( NULL == fit_par )
00975     {
00976         sinfo_msg_error("no fit parameters given") ;
00977         return -1 ;
00978     }
00979     if ( NULL == derv_par )
00980     {
00981         sinfo_msg_error("no derivatives of fit parameters given") ;
00982         return -1 ;
00983     }
00984     if ( lleftx < 0 || lleftx + 2*halfbox_x >= ilx ||
00985          llefty < 0 || llefty + 2*halfbox_y >= ily )
00986     {
00987         sinfo_msg_error("wrong lower left point of fitting box given!") ;
00988         return -1 ;
00989     }
00990     if ( halfbox_x <= 1 || halfbox_y <= 1 )
00991     {
00992         sinfo_msg_error("wrong box dimensions given") ;
00993         return -1 ;
00994     }
00995     /* allocate memory */
00996     if ( NULL == (backarray = (float*) cpl_calloc(4*halfbox_x+4*halfbox_y, 
00997                   sizeof(float) ) ) ) 
00998     {
00999         sinfo_msg_error("could not allocate memory") ;
01000         return -1 ;
01001     }
01002 
01003     /* -------------------------------------------------------------------
01004      * find the initial estimates for the free parameters
01005      */
01006 
01007     /* first search for the position of the maximum intensity */
01008     foundrow = 0 ;
01009     foundcol = 0 ;
01010     maxval   = -SINFO_DBL_MAX ;
01011     pidata=cpl_image_get_data_float(image);
01012     for ( col = lleftx ; col < lleftx + 2*halfbox_x ; col++ )
01013     {
01014         for ( row = llefty ; row < llefty + 2*halfbox_y ; row++ )
01015         {
01016             if ( isnan(pidata[col+row*ilx]) )
01017             {
01018                 continue ;
01019             }
01020             if ( maxval < pidata[col+row*ilx] ) 
01021             {
01022                 maxval = pidata[col+row*ilx] ;
01023                 foundrow = row ;
01024                 foundcol = col ;
01025             }
01026         }
01027     }
01028 
01029     if ( foundrow == 0 || foundcol == 0 || maxval <= 0. ||
01030          foundrow == ilx-1 || foundcol == ily-1 )
01031     {
01032         sinfo_msg_warning("no maximum found") ;
01033         cpl_free(backarray) ;
01034         return -1 ;
01035     }
01036 
01037     /* determine the lower left sinfo_edge of the fitting box, center it 
01038        on the maximum value */
01039     llx = foundcol - halfbox_x ;
01040     lly = foundrow - halfbox_y ;
01041     if ((foundcol - halfbox_x) > 0) {
01042       llx = (foundcol - halfbox_x);
01043     } else {
01044       llx=1;
01045       check++;
01046     }
01047 
01048     if ((foundrow - halfbox_y) > 0) {
01049       lly = (foundrow - halfbox_y);
01050     } else {
01051       lly=1;
01052       check++;
01053     } 
01054 
01055     if ( ( llx + 2*halfbox_x) <  ilx-1 ) {
01056       halfbox_x=halfbox_x;
01057     } else {
01058        halfbox_x=(int) (ilx-2-llx)/2;
01059       check++;
01060     }
01061 
01062     if ( ( lly + 2*halfbox_y) <  ily-1 ) {
01063       halfbox_y= halfbox_y;
01064     } else {
01065       halfbox_y=(int) (ily-2-lly)/2;
01066       check++;
01067     }
01068 
01069     if ( llx <= 0 || lly < 0 || llx + 2*halfbox_x >= ilx-1 ||
01070          lly + 2*halfbox_y >= ily )
01071     {
01072         sinfo_msg_error("box does not fit into image") ;
01073         cpl_free(backarray) ;
01074         return -1 ;
01075     }
01076     
01077     /* determine the zeroth and first order moments of the image  
01078        within the fitting box */ 
01079     M = Mx = My = 0. ;
01080     n = 0 ;
01081     boxi = boxj = 0 ;
01082     for ( j = lly ; j < lly + 2*halfbox_y ; j++ )
01083     {
01084         boxj = j - lly ;
01085         for ( i = llx ; i < llx + 2*halfbox_x ; i++ )
01086         {
01087             boxi = i - llx ;
01088             if ( !isnan(pidata[i+j*ilx]) )
01089             {
01090                 M  += pidata[i+j*ilx] ;  
01091                 Mx += (double)boxi * pidata[i+j*ilx] ;
01092                 My += (double)boxj * pidata[i+j*ilx] ;
01093                 /*-----------------------------------------------------------
01094                  * estimate the amplitude and the background 
01095                  * go through the margins of the fitting box 
01096                  * and calculate the clean mean to
01097                  * determine the background 
01098                  */
01099                 if ( i == llx || i == llx + 2*halfbox_x -1 ||
01100                      j == lly || j == lly + 2*halfbox_y -1 )
01101                 {
01102                     backarray[n] = pidata[i+j*ilx] ;
01103                     n++ ;
01104                 }
01105             }
01106         }
01107     }
01108     if ( M <= 0. )
01109     {
01110         sinfo_msg_warning("only negative or zero values") ;
01111         cpl_free(backarray) ;
01112         return -1 ;
01113     }
01114     if ( n < 3 )
01115     {
01116         sinfo_msg_error("not enough data points to calculate background") ;
01117         cpl_free(backarray) ;
01118         return -1 ;
01119     }
01120     /* determine the background as sinfo_median of the surrounding pixels */
01121     if (FLT_MAX==(background=sinfo_new_clean_mean(backarray,n,10.,10.))) 
01122     {
01123         sinfo_msg_error("it was not possible to compute the "
01124                         "clean mean of the background values") ;
01125         cpl_free(backarray) ;
01126         return -1 ;
01127     }
01128     cpl_free (backarray) ;
01129     /* now calculate the amplitude estimation */
01130     amplitude = maxval - background ;
01131     if ( amplitude < 1e-12 )
01132     {
01133         sinfo_msg_warning("amplitude is too small") ;
01134         return -1 ;
01135     }
01136 
01137     /* determine the center of gravity = centroid */
01138     X0 = Mx / M ;
01139     Y0 = My / M ;
01140     /* if one of the values is outside the fitting box return with error */
01141     if ( X0 <= 0. || Y0 <= 0. || X0 >= 2.*(double)halfbox_x || 
01142          Y0 >= 2.*(double)halfbox_y )
01143     {
01144         sinfo_msg_warning("center of gravity is outside the fitting box!") ;
01145         return -1 ;
01146     }
01147     
01148     /*------------------------------------------------------------------------ 
01149      * put the data in the 2-d array xydat[][] (pixel position) and zdat[] 
01150      * (data values) additionally, determine the second order momentum
01151      */
01152     n = 0 ;
01153     M = Mx = Mxx = My = Myy = Mxy = 0. ;
01154     boxi = boxj = 0 ;
01155     for ( j = lly ; j < lly + 2*halfbox_y ; j++ )
01156     {
01157         boxj = j - lly ;
01158         for ( i = llx ; i < llx + 2*halfbox_x ; i++ )
01159         {
01160             boxi = i - llx ;
01161             value = pidata[i+j*ilx] ;
01162             if ( !isnan(value) )
01163             {
01164                 xydat[n][0] = (double) boxi ;
01165                 xydat[n][1] = (double) boxj ;
01166                 zdat[n]     = value ;
01167                 wdat[n]     = 1. ;
01168                 n++ ; 
01169                
01170                 /* now calculate the moments without background in the 
01171                    centroid coordinate system */
01172                 value -= background ;
01173                 xco = (double) boxi - X0 ;
01174                 yco = (double) boxj - Y0 ;
01175                 M   += value ;
01176                 Mx  += xco * value ;
01177                 My  += yco * value ;
01178                 Mxx += xco * xco * value ;
01179                 Myy += yco * yco * value ;
01180                 Mxy += xco * yco * value ;
01181             }
01182         }
01183     }
01184     if ( M <= 0. )
01185     {
01186         sinfo_msg_warning("only negative or zero values") ;
01187         return -1 ;
01188     }
01189   
01190     /* ----------------------------------------------------------------
01191      * estimate the fwhm_x and fwhm_y and theta 
01192      */ 
01193   
01194     /* first scale the moments */
01195     Mx  /= M ;
01196     My  /= M ;
01197     Mxx /= M ;
01198     Myy /= M ;
01199     Mxy /= M ;
01200    
01201     denom = 2. * (Mxx*Myy - Mxy*Mxy) ;
01202     if ( denom == 0. )
01203     {
01204         sinfo_msg_error("denominator is zero!") ;
01205         return -1 ;
01206     }
01207     
01208     /* now associate the parameter list with the found estimates */
01209     fit_par[0] = X0 ;
01210     fit_par[1] = Y0 ;
01211     fit_par[2] = amplitude ;
01212     fit_par[3] = background ;
01213     fit_par[4] = Myy/denom ;
01214     fit_par[5] = Mxx/denom ;
01215     fit_par[6] = -Mxy/denom ;
01216 
01217     /* convert the moments to ellipse paramters */
01218     if ( 0 > new_gauss2Ellipse (fit_par) )
01219     {
01220         sinfo_msg_warning("gauss2Ellipse does not run!") ;
01221         return -1 ;
01222     }
01223 
01224     /* total number of data points */
01225     ndata = 4 * halfbox_x * halfbox_y ;
01226     xdim = XDIMG ; /* dimension of xydat array */
01227     npar = NPAR ; /* number of parameters in the fit */
01228     its = ITSG ;
01229     lab = LABG ;
01230     tol = TOLG ;
01231     for ( i = 0 ; i < NPAR ; i++ )
01232     {
01233         derv_par[i] = 0. ;
01234     }
01235     
01236     if ( 0 > ( iters = sinfo_new_lsqfitd ( &xydat[0][0],
01237                                  &xdim,
01238                                  zdat,
01239                                  wdat, 
01240                                  &ndata,
01241                                  fit_par,
01242                                  derv_par,
01243                                  mpar,
01244                                  &npar,
01245                                  &tol,
01246                                  &its,
01247                                  &lab )) ) 
01248     {
01249         sinfo_msg_warning(" least squares fit failed, error no: %d!", iters) ;
01250         return -1 ;
01251     }
01252 
01253     /* exclude impossible fit results */
01254     if ( fit_par[2] <= 0. || fit_par[4] < 0. || fit_par[5] < 0. )
01255     {
01256         sinfo_msg_error("sorry, some impossible negative fit results!") ;
01257         return -1 ;
01258     }
01259     fit_par[0] += llx ;
01260     fit_par[1] += lly ;
01261     if ( fit_par[0] < llx || fit_par[0] >= llx + 2*halfbox_x ||
01262          fit_par[1] < lly || fit_par[1] >= lly + 2*halfbox_y )
01263     {
01264         sinfo_msg_error("sorry, centroid after the fit "
01265                         "outside the fitting box") ;
01266         return -1 ;
01267     }
01268 
01269     /* exchange fwhmx and fwhmy if |theta| is bigger than 
01270        pi/4 and subtract pi/2 from theta */
01271     if ( fabs ( fit_par[6] ) > PI_NUMB / 4. )
01272     {
01273         /* first convert angle to smaller than 2 pi */
01274         if ( fabs (fit_par[6]) >= 2. * PI_NUMB )
01275         { 
01276             k = (int) (fit_par[6] / (2.*PI_NUMB)) ;
01277             if ( k > 0 ) 
01278             {
01279                 fit_par[6] -= k*2.*PI_NUMB ;
01280             }
01281             else
01282             {
01283                 fit_par[6] += k*2.*PI_NUMB ;
01284             }
01285         }
01286         /* first convert angle to smaller than pi/2 */
01287         if ( fabs (fit_par[6]) > PI_NUMB / 2. )
01288         {
01289             if ( fit_par[6] > 0. )
01290             {
01291                 fit_par[6] -= PI_NUMB ;
01292             }
01293             else
01294             {
01295                 fit_par[6] += PI_NUMB ;
01296             }
01297         }
01298 
01299         if ( fabs (fit_par[6]) > PI_NUMB / 4. )
01300         {
01301             temp       = fit_par[4] ;
01302             fit_par[4] = fit_par[5] ;
01303             fit_par[5] = temp ;
01304             if ( fit_par[6] < 0. )
01305             { 
01306                 fit_par[6] += PI_NUMB / 2. ;
01307             }
01308             else
01309             {
01310                 fit_par[6] -= PI_NUMB / 2. ;
01311             }  
01312         }
01313     }
01314     
01315     return iters ;
01316 }
01317 
01327 cpl_image * 
01328 sinfo_new_plot_gaussian (cpl_image   * image, 
01329                          double     * parlist )
01330 {
01331     int col, row ;
01332     cpl_image * retImage ;
01333     double xdat[2] ;
01334     int ilx=0;
01335     int ily=0;
01336     float* podata=NULL;
01337 
01338     if ( image == NULL )
01339     {
01340         sinfo_msg_error("no input image given!") ;
01341         return NULL ;
01342     }
01343     ilx=cpl_image_get_size_x(image);
01344     ily=cpl_image_get_size_y(image);
01345 
01346     if ( parlist == NULL ) 
01347     {
01348         sinfo_msg_error("no Gaussian parameters given!") ;
01349         return NULL ;
01350     }
01351  
01352     retImage = cpl_image_new (ilx, ily, CPL_TYPE_FLOAT) ;
01353     podata=cpl_image_get_data_float(retImage);
01354     for ( row = 0 ; row < ily ; row++ )
01355     {
01356         for ( col = 0 ; col < ilx ; col++ )
01357         {
01358             xdat[0] = (double) col ;
01359             xdat[1] = (double) row ;
01360             podata[col+row*ilx] = sinfo_new_gaussian_ellipse( xdat , parlist) ; 
01361         }
01362     }
01363 
01364     return retImage ;
01365 }
01366 
01374 static int new_gauss2Ellipse ( double     * parlist )
01375 {
01376     double a, b, c ;
01377     double ellipseconst ;
01378     double axisX, axisY, phi ;
01379     double p ;
01380     
01381     if ( parlist == NULL )
01382     {
01383         sinfo_msg_error(" no parameters given!\n") ;
01384         return -1 ;
01385     }
01386 
01387     a = parlist[4] ; /* fwhmx */
01388     b = parlist[5] ; /* fwhmy */
01389     c = parlist[6] ; /* theta */
01390 
01391     ellipseconst = 2. * log(2.) ;
01392 
01393     if ( a*b - c*c <= 0. )
01394     {
01395         sinfo_msg_warning("estimates of moments are unusable, "
01396                           "they do not make an ellipse!") ;
01397         return -1 ;
01398     }
01399     
01400     if ( a == b )
01401     { 
01402         phi = 0. ;
01403     }
01404     else
01405     {
01406         phi = 0.5 * atan( 2. * c / (a-b) ) ;
01407     }
01408 
01409     p = sqrt ( (a-b) * (a-b) + 4. * c*c ) ;
01410 
01411     if ( a > b )
01412     {
01413         axisX = 2. * sqrt ( ellipseconst / (a+b+p) ) ;
01414         axisY = 2. * sqrt ( ellipseconst / (a+b-p) ) ;
01415     }
01416     else
01417     {
01418         axisX = 2. * sqrt ( ellipseconst / (a+b-p) ) ;
01419         axisY = 2. * sqrt ( ellipseconst / (a+b+p) ) ;
01420     }
01421 
01422     parlist[4] = axisX ;
01423     parlist[5] = axisY ;
01424     parlist[6] = phi ;
01425 
01426     return 0 ;
01427 }
01428 
01452 float sinfo_new_determine_conversion_factor ( cpl_imagelist * cube, 
01453                                   float     mag,
01454                                   float     exptime,
01455                                   int       llx,
01456                                   int       lly,
01457                                   int       halfbox_x,
01458                                   int       halfbox_y, 
01459                                   int* check )
01460 {
01461     int row, col, i ;
01462     int first_row, first_col ;
01463     int last_row, last_col ;
01464     float factor ;
01465     int mpar[7] ;
01466     double fit_par[7] ;
01467     double derv_par[7] ;
01468     int fitInd ;
01469     double sum ;
01470     double xdat[2] ;
01471     cpl_image * summedIm ;
01472 
01473     int ilx=0;
01474     int ily=0;
01475     int inp=0;
01476 
01477     if ( NULL == cube )
01478     {
01479         sinfo_msg_error(" no cube given!\n") ;
01480         return -FLT_MAX ;
01481     }
01482 
01483     ilx=cpl_image_get_size_x(cpl_imagelist_get(cube,0));
01484     ily=cpl_image_get_size_y(cpl_imagelist_get(cube,0));
01485     inp=cpl_imagelist_get_size(cube);
01486 
01487     if ( halfbox_x <= 0 || halfbox_y <= 0 || 
01488        2*halfbox_x > ilx || 2*halfbox_y > ily)
01489     {
01490         sinfo_msg_error(" wrong width of halfbox given!") ;
01491         return -FLT_MAX ;
01492     }
01493     if ( exptime <= 0. )
01494     {
01495         sinfo_msg_error(" impossible exposure time given !") ;
01496         return -FLT_MAX ;
01497     }
01498 
01499     /* collapse the cube to be able to do 2D-Gaussian fitting */
01500     if ( NULL == (summedIm = sinfo_new_sum_cube_to_image(cube)) )
01501     {
01502         sinfo_msg_error(" sinfo_averageCubeToImage failed!") ;
01503         return -FLT_MAX ;
01504     }
01505 
01506     /* call the 2D-Gaussian fit routine */
01507     for ( i = 0 ; i < 7 ; i++ )
01508     {
01509         mpar[i] = 1 ;
01510     }
01511     if ( -1 == (fitInd = sinfo_new_fit_2d_gaussian(summedIm, fit_par, derv_par,
01512                                                    mpar, llx, lly, halfbox_x, 
01513                                                    halfbox_y, check)) )
01514     {
01515         sinfo_msg_warning("sinfo_fit2dGaussian failed!") ;
01516         cpl_image_delete( summedIm) ;
01517         return -FLT_MAX ;
01518     }
01519     cpl_image_delete(summedIm) ;
01520 
01521     /* now integrate the found 2D Gaussian by first 
01522        subtracting the background */
01523     if  ((fit_par[0] - halfbox_x) < 0) {
01524       first_col=0;
01525       check++;
01526     } else {
01527       first_col=(fit_par[0] - halfbox_x);
01528     }
01529 
01530     if ((fit_par[0] + halfbox_x) < ilx) {
01531       last_col  = (fit_par[0] + halfbox_x);
01532     } else {
01533       last_col = (ilx-1) ;
01534       check++;
01535     }
01536 
01537     if ((fit_par[1] - halfbox_y) < 0) {
01538       first_row=0;
01539       check++;
01540     } else {
01541       first_row=(fit_par[1] - halfbox_y) ;
01542     }
01543 
01544     if ((fit_par[1] + halfbox_y) < ily) {
01545       last_row=(fit_par[1] + halfbox_y);
01546     } else {
01547       last_row= (ily-1);
01548       check++;
01549     }
01550 
01551 
01552     if ( first_col < 0 || first_row < 0 || last_col >= ilx || last_row >= ily )
01553     {
01554         sinfo_msg_error("star badly centered in FOV or fitting box too big!") ;
01555         return -FLT_MAX ;
01556     }
01557     sum = 0. ;
01558     for ( row = first_row ; row < last_row ; row++ )  
01559     {
01560         for( col = first_col ; col < last_col ; col++ )
01561         {
01562             xdat[0] = (double) col ;
01563             xdat[1] = (double) row ;
01564             sum += (sinfo_new_gaussian_ellipse( xdat, fit_par ) - fit_par[3]) ;
01565         }
01566     }
01567     if ( sum <= 0. )
01568     {
01569         sinfo_msg_error("zero or negative sum of counts!") ;
01570         return -FLT_MAX ;
01571     }
01572     factor = mag / (float)sum * exptime ;
01573     return factor ;
01574 }
01575                               
01576 /*--------------------------------------------------------------------------*/

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