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