Object-oriented Scientific Computing Library: Version 0.910
gsl_mmin_simp2.h
00001 /*
00002   -------------------------------------------------------------------
00003 
00004   Copyright (C) 2006-2012, Andrew W. Steiner
00005 
00006   This file is part of O2scl.
00007 
00008   O2scl is free software; you can redistribute it and/or modify
00009   it under the terms of the GNU General Public License as published by
00010   the Free Software Foundation; either version 3 of the License, or
00011   (at your option) any later version.
00012 
00013   O2scl is distributed in the hope that it will be useful,
00014   but WITHOUT ANY WARRANTY; without even the implied warranty of
00015   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016   GNU General Public License for more details.
00017 
00018   You should have received a copy of the GNU General Public License
00019   along with O2scl. If not, see <http://www.gnu.org/licenses/>.
00020 
00021   -------------------------------------------------------------------
00022 */
00023 /* multimin/simplex2.c
00024  * 
00025  * Copyright (C) 2007 Brian Gough
00026  * Copyright (C) 2002 Tuomo Keskitalo, Ivo Alxneit
00027  * 
00028  * This program is free software; you can redistribute it and/or modify
00029  * it under the terms of the GNU General Public License as published by
00030  * the Free Software Foundation; either version 3 of the License, or (at
00031  * your option) any later version.
00032  * 
00033  * This program is distributed in the hope that it will be useful, but
00034  * WITHOUT ANY WARRANTY; without even the implied warranty of
00035  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00036  * General Public License for more details.
00037  * 
00038  * You should have received a copy of the GNU General Public License
00039  * along with this program; if not, write to the Free Software
00040  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 
00041  * 02110-1301, USA.
00042  */
00043 #ifndef O2SCL_GSL_MMIN_SIMP2_H
00044 #define O2SCL_GSL_MMIN_SIMP2_H
00045 
00046 #include <string>
00047 
00048 #include <gsl/gsl_multimin.h>
00049 
00050 #include <o2scl/multi_min.h>
00051 #include <o2scl/cblas.h>
00052 #include <o2scl/vec_stats.h>
00053 
00054 #ifndef DOXYGENP
00055 namespace o2scl {
00056 #endif
00057 
00058   /** \brief Multidimensional minimization by the Simplex method (v2) (GSL)
00059 
00060       This class minimizes a function using Nelder and Mead's Simplex
00061       algorithm. A simplex in a N-dimensional space is defined as a
00062       set of N+1 points which describe an N-dimensional volume
00063       surrounding the minimum. The algorithm proceeds by shifting the
00064       simplex points until the simplex is sufficiently small and thus
00065       the minimum is known with sufficient accuracy.
00066 
00067       For the earlier method used in GSL, see \ref gsl_mmin_simp .
00068 
00069       This class has a high-level interface using mmin(),
00070       mmin_twovec() or mmin_simplex() which automatically performs the
00071       memory allocation and minimization, or a GSL-like interface
00072       using allocate(), free(), interate() and set() or set_simplex().
00073 
00074       The simplex can be completely specified by the user (see
00075       mmin_simplex() and set_simplex()). Alternatively, the simplex is
00076       automatically specified given initial guess \f$ x_j \f$ and a
00077       step size vector \f$ s_k \f$ for \f$ 0\leq k<n_s \f$. The
00078       simplex \f$ p_{ij} \f$ with \f$ 0\leq i\leq n \f$ and \f$ 0\leq
00079       j<n \f$ is chosen with \f$ p_{0j} = x_j \f$ and
00080       \f{eqnarray*}
00081       p_{i+1,j} &=& x_j\quad\mathrm{for}\quad i\neq j \\
00082       p_{i+1,j} &=& x_j+s_{j~\mathrm{mod}~n_s}\quad\mathrm{for}\quad i=j
00083       \f}
00084       for \f$ 0<i<n \f$. The step size vector \f$ s \f$ is set by the
00085       set_step() member function. The prsence of \f$ \mathrm{mod} \f$
00086       in the recipe above just indicates that elements of the step
00087       size vector are automatically re-used if there are less step
00088       sizes than dimensions in the minimization problem.
00089 
00090       \note It is important that the initial simplex contains
00091       sufficient variation in every direction of the parameter space
00092       over which one is minimizing. For example, if all three points
00093       in a simplex for minimizing over a two-dimensional space contain
00094       nearly the same value for the second parameter, then the
00095       minimizer may only minimize the function with respect to the
00096       first parameter. 
00097 
00098       See an example for the usage of this class
00099       in \ref ex_mmin_sect .
00100 
00101       Default template arguments
00102       - \c param_t - no default
00103       - \c func_t - \ref multi_funct<param_t>
00104       - \c vec_t - \ref ovector_base
00105       - \c alloc_vec_t - \ref ovector
00106       - \c alloc_t - \ref ovector_alloc
00107 
00108       Based on \ref Nelder65 .
00109 
00110       A variable <tt>count</tt> originally defined in the GSL simplex
00111       state is not present here, because it was unused.
00112 
00113       \future Double check that the updates in gsl-1.13 are included
00114       here, and also add support for the nmsimplex2rand algorithm
00115       in GSL.
00116   */
00117 #ifdef DOXYGENP
00118   template<class func_t=multi_funct<>,
00119     class vec_t=ovector_base, class alloc_vec_t=ovector,
00120     class alloc_t=ovector_alloc> class gsl_mmin_simp2 :
00121   public multi_min
00122 #else
00123     template<class func_t=multi_funct<>,
00124     class vec_t=ovector_base, class alloc_vec_t=ovector,
00125     class alloc_t=ovector_alloc> class gsl_mmin_simp2 :
00126   public multi_min<func_t,func_t,vec_t>
00127 #endif
00128     {
00129     
00130 #ifndef DOXYGEN_INTERNAL
00131     
00132       protected:
00133 
00134       /// An array of n+1 vectors containing the simplex
00135       alloc_vec_t *x1;
00136       /** \brief The n+1 function values at the simplex points
00137           
00138           \comment
00139           This can't be the template type because it has to be 
00140           a different size (n+1) rather than n.
00141           \endcomment
00142       */
00143       ovector y1;
00144       /// Workspace vector 1
00145       alloc_vec_t ws1;
00146       /// Workspace vector 2
00147       alloc_vec_t ws2;
00148       /// Workspace vector 3
00149       alloc_vec_t ws3;
00150       /// Center of simplex
00151       alloc_vec_t center;
00152       /// Desc
00153       alloc_vec_t delta;
00154       /// Distance of vector from center
00155       alloc_vec_t xmc;
00156       /// Squared simplex size
00157       double S2;
00158     
00159       /// Compute the center of the simplex 
00160       int compute_center() {
00161         size_t i, j;
00162         double val;
00163       
00164         for (j = 0; j < dim; j++) {
00165           val = 0.0;
00166           for (i = 0; i < dim+1; i++) {
00167             val += x1[i][j];
00168           }
00169           val /= dim+1;
00170           center[j]=val;
00171         }
00172       
00173         return 0;
00174       }
00175     
00176       /** \brief Compute the size of the simplex
00177         
00178           Calculates simplex size as average sum of length of vectors
00179           from simplex center to corner points:
00180           
00181           \f$ (1/n) \sum || y - y_{\mathrm{center}} || \f$
00182       */
00183       double compute_size() {
00184         
00185         size_t P=dim+1;
00186         double ss=0.0;
00187         
00188         for(size_t i=0;i<P;i++) {
00189           for(size_t j=0;j<dim;j++) ws1[j]=x1[i][j];
00190           o2scl_cblas::daxpy(-1.0,dim,center,ws1);
00191           double t=o2scl_cblas::dnrm2(dim,ws1);
00192           ss += t*t;
00193         }
00194         
00195         /* Store squared size in the state */
00196         S2=ss/P;
00197         return sqrt(ss/(double)(P));
00198       }
00199 
00200       /** \brief Move a corner of a simplex
00201 
00202           Moves a simplex corner scaled by coeff (negative value
00203           represents mirroring by the middle point of the "other" corner
00204           points) and gives new corner in xc and function value at xc 
00205           in \c newval.
00206       */
00207       virtual int try_corner_move(const double coeff, size_t corner, 
00208                                   vec_t &xc, func_t &f, size_t nvar, 
00209                                   double &newval) {
00210         
00211         size_t P=dim+1;
00212         
00213         /* 
00214            xc = (1-coeff)*(P/(P-1)) * center(all) + 
00215            ((P*coeff-1)/(P-1))*x_corner 
00216         */
00217         double alpha=(1.0-coeff)*P/dim;
00218         double beta=(P*coeff-1.0)/dim;
00219         
00220         for(size_t j=0;j<dim;j++) {
00221           xc[j]=center[j]*alpha;
00222           xc[j]+=x1[corner][j]*beta;
00223         }
00224         
00225         newval=f(nvar,xc);
00226 
00227         return 0;
00228       }
00229 
00230       /// Update point \c i in the simplex with values \c xx
00231       virtual int update_point(size_t i, vec_t &xx, double val) {
00232 
00233         const size_t P=dim+1;
00234 
00235         /* Compute delta = x - x_orig */
00236         for(size_t j=0;j<dim;j++) {
00237           delta[j]=xx[j];
00238           delta[j]-=x1[i][j];
00239         }
00240 
00241         /* Compute xmc = x_orig - c */
00242         for(size_t j=0;j<dim;j++) {
00243           xmc[j]=x1[i][j];
00244           xmc[j]-=center[j];
00245         }
00246         
00247         /* Update size: S2' = S2 + (2/P) * 
00248            (x_orig - c).delta + (P-1)*(delta/P)^2 
00249         */
00250         double d=o2scl_cblas::dnrm2(dim,delta);
00251         double xmcd=o2scl_cblas::ddot(dim,xmc,delta);
00252         S2 += (2.0 / P) * xmcd + ((P - 1.0) / P) * (d * d / P);
00253 
00254         /* Update center:  c' = c + (x - x_orig) / P */
00255         double alpha=1.0/P;
00256         for(size_t j=0;j<dim;j++) {
00257           center[j]-=alpha*x1[i][j];
00258           center[j]+=alpha*xx[j];
00259         }
00260         
00261         for(size_t j=0;j<dim;j++) {
00262           x1[i][j]=xx[j];
00263         }
00264         y1[i]=val;
00265 
00266         return 0;
00267       }
00268       
00269       /** \brief Contract the simplex towards the best point
00270           
00271           Function contracts the simplex in respect to best valued
00272           corner. All corners besides the best corner are moved.
00273         
00274           \comment 
00275           The GSL version requires a vector for workspace
00276           named 'xc' which is not required here.
00277           \endcomment
00278       */
00279       virtual int contract_by_best(size_t best, func_t &f, 
00280                                    size_t nvar) {
00281         
00282         /* Function contracts the simplex in respect to best valued
00283            corner. That is, all corners besides the best corner are
00284            moved. (This function is rarely called in practice, since
00285            it is the last choice, hence not optimised - BJG) 
00286         */
00287         size_t i,j, it;
00288         double newval;
00289         bool failed;
00290 
00291         int status=gsl_success;
00292         
00293         for (i=0;i<dim+1;i++) {
00294         
00295           if (i!=best) {
00296           
00297             for (j=0;j<dim;j++) {
00298               x1[i][j]=0.5*(x1[i][j]+x1[best][j]);
00299             }
00300             y1[i]=f(nvar,x1[i]);
00301             if (!gsl_finite(y1[i])) {
00302               std::string err=((std::string)"Function not finite (returned ")+
00303                 dtos(y1[i])+" in gsl_mmin_simp2::contract_by_best().";
00304               O2SCL_ERR_RET(err.c_str(),gsl_ebadfunc);
00305             }
00306 
00307             /*
00308             if (avoid_nonzero==true && ret!=0) {
00309               std::cout << "Found problem in contract." << std::endl;
00310 
00311               // copy the old point to ws3
00312               for (j=0;j<dim;j++) {
00313                 ws3[j]=2.0*x1[i][j]-x1[best][j];
00314               }
00315             
00316               // Try (21*best+20*xold)/41, (22*best+19*xold)/41, ...,
00317               // (40*best+xold)/41 to see if we can find a contraction
00318               // that works
00319               for(it=0;ret!=0 && it<20;it++) {
00320                 for (j=0;j<dim;j++) {
00321                   x1[i][j]=((20-it)*ws3[j]+(it+21)*x1[best][j])/41.0;
00322                 }
00323                 y1[i]=f(nvar,x1[i]);
00324                 std::cout << "it, ret, x: " << it << " " << ret << " "
00325                           << x << std::endl;
00326               }
00327               if (ret!=0) {
00328                 this->last_conv=gsl_efailed;
00329                 O2SCL_CONV2_RET("Failed to find suitable contraction ",
00330                                 "in gsl_mmin_simp2::contract_by_best().",
00331                                 gsl_efailed,this->err_nonconv);
00332               }
00333             }
00334             */
00335           }
00336         }
00337         
00338         compute_center();
00339         size=compute_size();
00340         
00341         return gsl_success;
00342       }
00343 
00344       /// Number of variables to be minimized over
00345       size_t dim;
00346 
00347       /// Function
00348       func_t *func;
00349 
00350       /// True if set() has been called
00351       bool set_called;
00352 
00353       /// Vector of step sizes
00354       ovector step_vec;
00355 
00356       /// Vector allocator
00357       alloc_t ao;
00358 
00359       /** \brief If true, try to automatically avoid regions where the 
00360           function returns a non-zero value (default false)
00361           
00362           \note This option doesn't work yet, so I've made the variable
00363           protected to prevent the casual user from changing it.
00364       */
00365       bool avoid_nonzero;
00366 
00367 #endif
00368 
00369       public:
00370 
00371       gsl_mmin_simp2() {
00372         dim=0;
00373         print_simplex=0;
00374         step_vec.allocate(1);
00375         step_vec[0]=1.0;
00376         avoid_nonzero=false;
00377       }
00378     
00379       virtual ~gsl_mmin_simp2() {
00380         free();
00381         step_vec.free();
00382       }
00383 
00384       /// Set the step sizes for each independent variable
00385       template<class vec2_t> int set_step(size_t nv, vec2_t &step) {
00386         if (nv>0) {
00387           step_vec.free();
00388           step_vec.allocate(nv);
00389           for(size_t i=0;i<nv;i++) step_vec[i]=step[i];
00390         }
00391         return 0;
00392       }
00393 
00394       /// Size of current simplex computed by iterate()
00395       double size;
00396 
00397       /// Present minimum vector computed by iterate()
00398       alloc_vec_t x;
00399 
00400       /// Function value at minimum computed by iterate()
00401       double fval;
00402 
00403       /** \brief Print simplex information in print_iter() (default 0)
00404         
00405           If this is 1 and \ref verbose is greater than 0, 
00406           then print_iter() will print the function values at all 
00407           the simplex points. If this is 2 and \ref verbose is
00408           greater than 0, then print_iter() will print the 
00409           simplex coordinates in addition to the function values.
00410       */
00411       int print_simplex;
00412 
00413       /** \brief Calculate the minimum \c min of \c func w.r.t the
00414           array \c x of size \c nvar.
00415       */
00416       virtual int mmin(size_t nn, vec_t &xx, double &fmin, 
00417                        func_t &ufunc) {
00418 
00419         if (nn==0) {
00420           O2SCL_ERR2_RET("Tried to minimize over zero variables ",
00421                          " in gsl_mmin_simp2::mmin().",gsl_einval);
00422         }
00423 
00424         this->last_conv=0;
00425         int ret=0,status,iter=0;
00426       
00427         allocate(nn);
00428 
00429         alloc_vec_t ss;
00430         ao.allocate(ss,nn);
00431         for (size_t is=0;is<nn;is++) ss[is]=step_vec[is % step_vec.size()];
00432         ret=set(ufunc,nn,xx,ss);
00433         ao.free(ss);
00434 
00435         if(ret!=0) {
00436           free();
00437           return ret;
00438         }
00439   
00440         do {
00441           iter++;
00442           
00443           status=iterate();
00444           if(status) break;
00445           
00446           if(this->verbose>0) {
00447             print_iter(nn,x,x1,fval,iter,size,this->tol_abs,
00448                        "gsl_mmin_simp2");
00449           }
00450     
00451           status=gsl_multimin_test_size(size,this->tol_abs);
00452           
00453         } while(status == GSL_CONTINUE && iter<this->ntrial);
00454         
00455         for (size_t i=0;i<nn;i++) xx[i]=x[i];
00456         fmin=fval;
00457   
00458         free();
00459         this->last_ntrial=iter;
00460         
00461         if(iter>=this->ntrial) {
00462           this->last_conv=gsl_emaxiter;
00463           std::string str="Exceeded maximum number of iterations ("+
00464             itos(this->ntrial)+") in mmin_simp2::mmin().";
00465           O2SCL_CONV_RET(str.c_str(),gsl_emaxiter,this->err_nonconv);
00466         }
00467 
00468         return status;
00469       }
00470     
00471       /** \brief Calculate the minimum \c min of \c func w.r.t the
00472           array \c x of size \c nvar, using \c xx and \c xx2 to specify
00473           the simplex
00474       */
00475       virtual int mmin_twovec(size_t nn, vec_t &xx, vec_t &xx2, double &fmin, 
00476                               func_t &ufunc) {
00477       
00478         this->last_conv=0;
00479         int ret=0,i,status,iter=0;
00480       
00481         allocate(nn);
00482 
00483         alloc_vec_t ss;
00484         ao.allocate(ss,nn);
00485 
00486         for (size_t is=0;is<nn;is++) ss[is]=xx2[is]-xx[is];
00487         ret=set(ufunc,nn,xx,ss);
00488         ao.free(ss);
00489       
00490         if(ret!=0) {
00491           free();
00492           return ret;
00493         }
00494   
00495         do {
00496           iter++;
00497           
00498           status=iterate();
00499           if(status) break;
00500 
00501           if(this->verbose>0) {
00502             print_iter(nn,x,x1,fval,iter,size,this->tol_abs,
00503                        "gsl_mmin_simp2");
00504           }
00505         
00506           status=gsl_multimin_test_size(size,this->tol_abs);
00507           
00508         } while(status == GSL_CONTINUE && iter<this->ntrial);
00509         
00510         for (i=0;i<((int)nn);i++) xx[i]=x[i];
00511         fmin=fval;
00512   
00513         free();
00514         this->last_ntrial=iter;
00515 
00516         if(iter>=this->ntrial) {
00517           this->last_conv=gsl_emaxiter;
00518           std::string str="Exceeded maximum number of iterations ("+
00519             itos(this->ntrial)+") in mmin_simp2::mmin_twovec().";
00520           O2SCL_CONV_RET(str.c_str(),gsl_emaxiter,this->err_nonconv);
00521         }
00522 
00523         return status;
00524       }
00525     
00526       /** \brief Calculate the minimum \c min of \c func w.r.t the
00527           array \c x of size \c nvar, given an initial simplex
00528       */
00529       template<class mat_t> 
00530       int mmin_simplex(size_t nn, mat_t &sx, double &fmin, 
00531                        func_t &ufunc) {
00532       
00533         int ret=0,i,status,iter=0;
00534       
00535         allocate(nn);
00536         
00537         ret=set_simplex(ufunc,sx);
00538         if(ret!=0) {
00539           free();
00540           return ret;
00541         }
00542   
00543         do {
00544           iter++;
00545           
00546           status=iterate();
00547           if(status) break;
00548 
00549           if(this->verbose>0) {
00550             print_iter(nn,x,x1,fval,iter,size,this->tol_abs,
00551                        "gsl_mmin_simp2");
00552           }
00553     
00554           status=gsl_multimin_test_size(size,this->tol_abs);
00555           
00556         } while(status == GSL_CONTINUE && iter<this->ntrial);
00557         
00558         for (i=0;i<((int)nn);i++) sx[0][i]=x[i];
00559         fmin=fval;
00560   
00561         free();
00562         this->last_ntrial=iter;
00563 
00564         if (iter>=this->ntrial) {
00565           this->last_conv=gsl_emaxiter;
00566           std::string str="Exceeded maximum number of iterations ("+
00567             itos(this->ntrial)+") in gsl_mmin_simp2::mmin_simplex().";
00568           O2SCL_CONV_RET(str.c_str(),gsl_emaxiter,this->err_nonconv);
00569         }
00570 
00571         return status;
00572       }
00573     
00574       /// Allocate the memory
00575       virtual int allocate(size_t n) {
00576         
00577         if(dim!=0) free();
00578         set_called=false;
00579         dim=n;
00580 
00581         ao.allocate(x,n);
00582         x1=new alloc_vec_t[n+1];
00583         for(size_t i=0;i<n+1;i++) ao.allocate(x1[i],n);
00584         y1.allocate(n+1);
00585         ao.allocate(ws1,n);
00586         ao.allocate(ws2,n);
00587         ao.allocate(ws3,n);
00588         ao.allocate(center,n);
00589         ao.allocate(delta,n);
00590         ao.allocate(xmc,n);
00591 
00592         return gsl_success;
00593       }
00594     
00595       /// Free the allocated memory
00596       virtual int free() {
00597 
00598         if (dim>0) {
00599           for(size_t i=0;i<dim+1;i++) ao.free(x1[i]);
00600           delete[] x1;
00601           y1.free();
00602           ao.free(ws1);
00603           ao.free(ws2);
00604           ao.free(ws3);
00605           ao.free(center);
00606           ao.free(delta);
00607           ao.free(xmc);
00608           ao.free(x);
00609         }
00610 
00611         dim=0;
00612 
00613         return 0;
00614       }
00615 
00616       /// Set the function and initial guess
00617       virtual int set(func_t &ufunc, size_t n, vec_t &ax, 
00618                       vec_t &step_size) {
00619         size_t i;
00620       
00621         if(dim!=n) allocate(n);
00622         func=&ufunc;
00623 
00624         // Copy initial guess to x
00625         for (i=0;i<dim;i++) x[i]=ax[i];
00626       
00627         // first point is the original x0 
00628       
00629         y1[0]=ufunc(dim,ax);
00630         if (!gsl_finite(y1[0])) {
00631           std::string err=((std::string)"Function not finite (returned ")+
00632             dtos(y1[0])+" in gsl_mmin_simp2::set().";
00633           O2SCL_ERR_RET(err.c_str(),gsl_ebadfunc);
00634         }
00635         for(i=0;i<dim;i++) x1[0][i]=ax[i];
00636   
00637         /* following points are initialized to x0+step_size */
00638       
00639         for (i=1;i<dim+1;i++) {
00640           for(size_t j=0;j<dim;j++) x1[i][j]=x[j];
00641           x1[i][i-1]=x1[i][i-1]+step_size[i-1];
00642           y1[i]=ufunc(dim,x1[i]);
00643         }
00644  
00645         /* Initialize simplex size */
00646       
00647         compute_center();
00648         size=compute_size();
00649         
00650         set_called=true;
00651 
00652         return gsl_success;
00653       }
00654 
00655       /// Set the function and initial simplex
00656       template<class mat_t> 
00657       int set_simplex(func_t &ufunc, mat_t &sx) {
00658 
00659         if(dim==0) {
00660           O2SCL_ERR2_RET("Memory not allocated in ",
00661                          "gsl_mmin_simp2::set_simplex().",gsl_ebadlen);
00662         }
00663         
00664         func=&ufunc;
00665 
00666         for(size_t i=0;i<dim+1;i++) {
00667           for(size_t j=0;j<dim;j++) {
00668             x1[i][j]=sx[i][j];
00669           }
00670           y1[i]=ufunc(dim,x1[i]);
00671           if (!gsl_finite(y1[i])) {
00672             std::string err=((std::string)"Function not finite (returned ")+
00673               dtos(y1[i])+" in gsl_mmin_simp2::set_simplex().";
00674             O2SCL_ERR_RET(err.c_str(),gsl_ebadfunc);
00675           }
00676         }
00677         
00678         /* Initialize simplex size */
00679         
00680         compute_center();
00681         size=compute_size();
00682         
00683         set_called=true;
00684  
00685         return gsl_success;
00686       }
00687 
00688       /// Perform an iteration
00689       virtual int iterate() {
00690 
00691         size_t n=dim+1;
00692         size_t i;
00693         size_t lo, hi, s_hi;
00694         double dhi,ds_hi,dlo;
00695         int status;
00696         double val,val2;
00697  
00698         /* get index of highest, second highest and lowest point */
00699 
00700         lo=hi=0;
00701         dlo=dhi=y1[0];
00702         s_hi=1;
00703         ds_hi=y1[1];
00704       
00705         for (i=1;i<n;i++) {
00706           val=y1[i];
00707           if(val<dlo) {
00708             dlo=val;
00709             lo=i;
00710           } else if (val > dhi) {
00711             ds_hi=dhi;
00712             s_hi=hi;
00713             dhi=val;
00714             hi=i;
00715           } else if(val > ds_hi) {
00716             ds_hi=val;
00717             s_hi=i;
00718           }
00719         }
00720 
00721         /* reflect the highest value */
00722         
00723         int ret1=try_corner_move(-1.0,hi,ws1,*func,dim,val);
00724 
00725         if (avoid_nonzero && ret1!=0) {
00726           std::cout << "Found problem move1: " << std::endl;
00727           for (size_t it=0;it<20 && ret1!=0;it++) {
00728             ret1=try_corner_move(-1.0+((double)it)/10.0,hi,ws1,
00729                                  *func,dim,val);
00730             std::cout << "it,ret: " << it << " " << ret1 << std::endl;
00731           }
00732           if (ret1!=0) {
00733             O2SCL_ERR2_RET("Failed to move corner (1) in ",
00734                            "gsl_mmin_simp2::iterate().",gsl_efailed);
00735           }
00736         }
00737       
00738         if (gsl_finite(val) && val<y1[lo]) {
00739 
00740           /* reflected point becomes lowest point,try expansion */
00741         
00742           int ret2=try_corner_move(-2.0,hi,ws2,*func,dim,val2);
00743 
00744           if (avoid_nonzero && ret2!=0) {
00745             std::cout << "Found problem move2: " << std::endl;
00746             for (size_t it=0;it<20 && ret2!=0;it++) {
00747               ret2=try_corner_move(-2.0+((double)it)/6.7,hi,ws2,
00748                                    *func,dim,val2);
00749               std::cout << "it,ret: " << it << " " << ret2 << std::endl;
00750             }
00751             if (ret2!=0) {
00752               O2SCL_ERR2_RET("Failed to move corner (2) in ",
00753                              "gsl_mmin_simp2::iterate().",gsl_efailed);
00754             }
00755           }
00756 
00757           if (gsl_finite(val2) && val2<y1[lo]) {
00758             update_point(hi,ws2,val2);
00759           } else {
00760             update_point(hi,ws1,val);
00761           }
00762           
00763         } else if (!gsl_finite(val) || val > y1[s_hi]) {
00764         
00765           /* reflection does not improve things enough */
00766         
00767           if (gsl_finite(val) && val <= y1[hi]) {
00768             
00769             /* if trial point is better than highest point,replace
00770                highest point */
00771             
00772             update_point(hi,ws1,val);
00773           }
00774       
00775           /* try one dimensional contraction */
00776         
00777           int ret3=try_corner_move(0.5,hi,ws2,*func,dim,val2);
00778 
00779           if (avoid_nonzero && ret3!=0) {
00780             std::cout << "Found problem move3: " << std::endl;
00781             for (size_t it=0;it<20 && ret3!=0;it++) {
00782               ret3=try_corner_move(0.025*((double)(19-it)),hi,ws2,
00783                                    *func,dim,val2);
00784               std::cout << "it,ret: " << it << " " << ret3 << std::endl;
00785             }
00786             if (ret3!=0) {
00787               O2SCL_ERR2_RET("Failed to move corner (2) in ",
00788                              "gsl_mmin_simp2::iterate().",gsl_efailed);
00789             }
00790           }
00791           
00792           if (gsl_finite(val2) && val2 <= y1[hi]) {
00793 
00794             update_point(hi,ws2,val2);
00795 
00796           } else {
00797 
00798             /* contract the whole simplex in respect to the best point */
00799             status=contract_by_best(lo,*func,dim);
00800             if(status != 0) {
00801               O2SCL_ERR("Function contract_by_best() failed in iterate().",
00802                         gsl_efailed);
00803             }
00804             
00805           }
00806 
00807         } else {
00808 
00809           /* trial point is better than second highest point.
00810              Replace highest point by it */
00811           
00812           update_point(hi,ws1,val);
00813         }
00814   
00815         /* return lowest point of simplex as x */
00816       
00817         vector_min(dim+1,y1,lo,val);
00818         for(i=0;i<dim;i++) x[i]=x1[lo][i];
00819         fval=y1[lo];
00820       
00821         /* Update simplex size */
00822         
00823         if (S2 > 0) {
00824           size=sqrt(S2);
00825         } else {
00826           /* recompute if accumulated error has made size invalid */
00827           size=compute_size();
00828         }
00829 
00830         return gsl_success;
00831       }
00832       
00833       /** \brief Print out iteration information.
00834           
00835           Depending on the value of the variable verbose, this prints out
00836           the iteration information. If verbose=0, then no information is
00837           printed, while if verbose>1, then after each iteration, the
00838           present values of x and y are output to std::cout along with the
00839           iteration number. If verbose>=2 then each iteration waits for a
00840           character.
00841       */
00842       virtual int print_iter(size_t nv, vec_t &xx, alloc_vec_t *simp,
00843                              double y, int iter,
00844                              double value, double limit,
00845                              std::string comment) {
00846         
00847         if (this->verbose<=0) return 0;
00848         
00849         size_t i;
00850         char ch;
00851       
00852         (*this->outs) << comment << " Iteration: " << iter << std::endl;
00853         (*this->outs) << "x: ";
00854         for(i=0;i<nv;i++) (*this->outs) << xx[i] << " ";
00855         (*this->outs) << std::endl;
00856         if (print_simplex>0) {
00857           if (print_simplex==2) {
00858             (*this->outs) << "Simplex Values:" << std::endl;
00859             for(i=0;i<nv+1;i++) {
00860               (*this->outs) << i << ": ";
00861               for(size_t j=0;j<nv;j++) {
00862                 (*this->outs) << simp[i][j] << " ";
00863               }
00864               (*this->outs) << ": " << y1[i] << std::endl;
00865             }
00866           } else {
00867             (*this->outs) << "Simplex Values:" << std::endl;
00868             for(i=0;i<nv+1;i++) {
00869               (*this->outs) << y1[i] << " ";
00870             }
00871             (*this->outs) << std::endl;
00872           }
00873         }
00874         (*this->outs) << "y: " << y << " Val: " << value << " Lim: " 
00875         << limit << std::endl;
00876         if (this->verbose>1) {
00877           (*this->outs) << "Press a key and type enter to continue. ";
00878           (*this->ins) >> ch;
00879         }
00880         
00881         return 0;
00882       }
00883 
00884       /// Return string denoting type("gsl_mmin_simp2")
00885       virtual const char *type() { return "gsl_mmin_simp2";}
00886 
00887     };
00888   
00889 #ifndef DOXYGENP
00890 }
00891 #endif
00892 
00893 #endif
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).

Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads.