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