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/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 <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 (v2) (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 the earlier method used in GSL, see \ref gsl_mmin_simp . 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 See an example for the usage of this class 00089 in \ref ex_mmin_sect . 00090 00091 Based on \ref Nelder65 . 00092 */ 00093 #ifdef DOXYGENP 00094 template<class param_t,class func_t=multi_funct<param_t>, 00095 class vec_t=ovector_base, class alloc_vec_t=ovector, 00096 class alloc_t=ovector_alloc> class gsl_mmin_simp2 : 00097 public multi_min 00098 #else 00099 template<class param_t,class func_t=multi_funct<param_t>, 00100 class vec_t=ovector_base, class alloc_vec_t=ovector, 00101 class alloc_t=ovector_alloc> class gsl_mmin_simp2 : 00102 public multi_min<param_t,func_t,func_t,vec_t> 00103 #endif 00104 { 00105 00106 #ifndef DOXYGEN_INTERNAL 00107 00108 protected: 00109 00110 /// An array of n+1 vectors containing the simplex 00111 alloc_vec_t *x1; 00112 /** 00113 \brief The n+1 function values at the simplex points 00114 00115 \comment 00116 This can't be the template type because it has to be 00117 a different size (n+1) rather than n. 00118 \endcomment 00119 */ 00120 ovector y1; 00121 /// Workspace vector 1 00122 alloc_vec_t ws1; 00123 /// Workspace vector 2 00124 alloc_vec_t ws2; 00125 /// Workspace vector 3 00126 alloc_vec_t ws3; 00127 /// Center of simplex 00128 alloc_vec_t center; 00129 /// Desc 00130 alloc_vec_t delta; 00131 /// Distance of vector from center 00132 alloc_vec_t xmc; 00133 /// Squared simplex size 00134 double S2; 00135 00136 /// Compute the center of the simplex 00137 int compute_center() { 00138 size_t i, j; 00139 double val; 00140 00141 for (j = 0; j < dim; j++) { 00142 val = 0.0; 00143 for (i = 0; i < dim+1; i++) { 00144 val += x1[i][j]; 00145 } 00146 val /= dim+1; 00147 center[j]=val; 00148 } 00149 00150 return 0; 00151 } 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 compute_size() { 00162 00163 size_t P=dim+1; 00164 double ss=0.0; 00165 00166 for(size_t i=0;i<P;i++) { 00167 for(size_t j=0;j<dim;j++) ws1[j]=x1[i][j]; 00168 o2scl_cblas::daxpy(dim,-1.0,center,ws1); 00169 double t=o2scl_cblas::dnrm2(dim,ws1); 00170 ss += t*t; 00171 } 00172 00173 /* Store squared size in the state */ 00174 S2=ss/P; 00175 return sqrt(ss/(double)(P)); 00176 } 00177 00178 /** 00179 \brief Move a corner of a simplex 00180 00181 Moves a simplex corner scaled by coeff (negative value 00182 represents mirroring by the middle point of the "other" corner 00183 points) and gives new corner in xc and function value at xc 00184 in \c newval. 00185 */ 00186 virtual int try_corner_move(const double coeff, size_t corner, 00187 vec_t &xc, func_t &f, size_t nvar, 00188 param_t &pa, double &newval) { 00189 00190 size_t P=dim+1; 00191 00192 /* 00193 xc = (1-coeff)*(P/(P-1)) * center(all) + 00194 ((P*coeff-1)/(P-1))*x_corner 00195 */ 00196 double alpha=(1.0-coeff)*P/dim; 00197 double beta=(P*coeff-1.0)/dim; 00198 00199 for(size_t j=0;j<dim;j++) { 00200 xc[j]=center[j]*alpha; 00201 xc[j]+=x1[corner][j]*beta; 00202 } 00203 00204 return f(nvar,xc,newval,pa); 00205 } 00206 00207 /// Desc 00208 virtual int update_point(size_t i, vec_t &xx, double val) { 00209 00210 const size_t P=dim+1; 00211 00212 /* Compute delta = x - x_orig */ 00213 for(size_t j=0;j<dim;j++) { 00214 delta[j]=xx[j]; 00215 delta[j]-=x1[i][j]; 00216 } 00217 00218 /* Compute xmc = x_orig - c */ 00219 for(size_t j=0;j<dim;j++) { 00220 xmc[j]=x1[i][j]; 00221 xmc[j]-=center[j]; 00222 } 00223 00224 /* Update size: S2' = S2 + (2/P) * 00225 (x_orig - c).delta + (P-1)*(delta/P)^2 00226 */ 00227 double d=o2scl_cblas::dnrm2(dim,delta); 00228 double xmcd=o2scl_cblas::ddot(dim,xmc,delta); 00229 S2 += (2.0 / P) * xmcd + ((P - 1.0) / P) * (d * d / P); 00230 00231 /* Update center: c' = c + (x - x_orig) / P */ 00232 double alpha=1.0/P; 00233 for(size_t j=0;j<dim;j++) { 00234 center[j]-=alpha*x1[i][j]; 00235 center[j]+=alpha*xx[j]; 00236 } 00237 00238 for(size_t j=0;j<dim;j++) { 00239 x1[i][j]=xx[j]; 00240 } 00241 y1[i]=val; 00242 00243 return 0; 00244 } 00245 00246 /** 00247 \brief Contract the simplex towards the best point 00248 00249 Function contracts the simplex in respect to best valued 00250 corner. All corners besides the best corner are moved. 00251 00252 \comment 00253 The GSL version requires a vector for workspace 00254 named 'xc' which is not required here. 00255 \endcomment 00256 */ 00257 virtual int contract_by_best(size_t best, func_t &f, 00258 size_t nvar, param_t &pa) { 00259 00260 /* Function contracts the simplex in respect to best valued 00261 corner. That is, all corners besides the best corner are 00262 moved. (This function is rarely called in practice, since 00263 it is the last choice, hence not optimised - BJG) 00264 */ 00265 size_t i,j, it; 00266 double newval; 00267 bool failed; 00268 00269 int status=gsl_success; 00270 00271 for (i=0;i<dim+1;i++) { 00272 00273 if (i!=best) { 00274 00275 for (j=0;j<dim;j++) { 00276 x1[i][j]=0.5*(x1[i][j]+x1[best][j]); 00277 } 00278 int ret=f(nvar,x1[i],y1[i],pa); 00279 if (!finite(y1[i])) { 00280 O2SCL_ERR2_RET("Function not finite in ", 00281 "gsl_mmin_simp2::contract_by_best().", 00282 gsl_ebadfunc); 00283 } 00284 00285 if (avoid_nonzero==true && ret!=0) { 00286 std::cout << "Found problem in contract." << std::endl; 00287 00288 // copy the old point to ws3 00289 for (j=0;j<dim;j++) { 00290 ws3[j]=2.0*x1[i][j]-x1[best][j]; 00291 } 00292 00293 // Try (21*best+20*xold)/41, (22*best+19*xold)/41, ..., 00294 // (40*best+xold)/41 to see if we can find a contraction 00295 // that works 00296 for(it=0;ret!=0 && it<20;it++) { 00297 for (j=0;j<dim;j++) { 00298 x1[i][j]=((20-it)*ws3[j]+(it+21)*x1[best][j])/41.0; 00299 } 00300 ret=f(nvar,x1[i],y1[i],pa); 00301 std::cout << "it, ret, x: " << it << " " << ret << " " 00302 << x << std::endl; 00303 } 00304 if (ret!=0) { 00305 this->last_conv=gsl_efailed; 00306 O2SCL_CONV2_RET("Failed to find suitable contraction ", 00307 "in gsl_mmin_simp2::contract_by_best().", 00308 gsl_efailed,this->err_nonconv); 00309 } 00310 } 00311 } 00312 } 00313 00314 return gsl_success; 00315 } 00316 00317 /// Number of variables to be minimized over 00318 size_t dim; 00319 00320 /// Function 00321 func_t *func; 00322 00323 /// Parameters 00324 param_t *params; 00325 00326 /// True if set() has been called 00327 bool set_called; 00328 00329 /// Vector of step sizes 00330 ovector step_vec; 00331 00332 /// Vector allocator 00333 alloc_t ao; 00334 00335 /** 00336 \brief If true, try to automatically avoid regions where the 00337 function returns a non-zero value (default false) 00338 00339 \note This option doesn't work yet, so I've made the variable 00340 protected to prevent the user from changing it. 00341 */ 00342 bool avoid_nonzero; 00343 00344 #endif 00345 00346 public: 00347 00348 gsl_mmin_simp2() { 00349 dim=0; 00350 print_simplex=0; 00351 step_vec.allocate(1); 00352 step_vec[0]=1.0; 00353 avoid_nonzero=false; 00354 } 00355 00356 virtual ~gsl_mmin_simp2() { 00357 free(); 00358 step_vec.free(); 00359 } 00360 00361 /// Set the step sizes for each independent variable 00362 template<class vec2_t> int set_step(size_t nv, vec2_t &step) { 00363 if (nv>0) { 00364 step_vec.free(); 00365 step_vec.allocate(nv); 00366 for(size_t i=0;i<nv;i++) step_vec[i]=step[i]; 00367 } 00368 return 0; 00369 } 00370 00371 /// Size of current simplex computed by iterate() 00372 double size; 00373 00374 /// Present minimum vector computed by iterate() 00375 alloc_vec_t x; 00376 00377 /// Function value at minimum computed by iterate() 00378 double fval; 00379 00380 /** 00381 \brief Print simplex information in print_iter() (default 0) 00382 00383 If this is 1 and \ref verbose is greater than 0, 00384 then print_iter() will print the function values at all 00385 the simplex points. If this is 2 and \ref verbose is 00386 greater than 0, then print_iter() will print the 00387 simplex coordinates in addition to the function values. 00388 */ 00389 int print_simplex; 00390 00391 /** \brief Calculate the minimum \c min of \c func w.r.t the 00392 array \c x of size \c nvar. 00393 */ 00394 virtual int mmin(size_t nn, vec_t &xx, double &fmin, param_t &pa, 00395 func_t &ufunc) { 00396 00397 if (nn==0) { 00398 O2SCL_ERR2_RET("Tried to minimize over zero variables ", 00399 " in gsl_mmin_simp2::mmin().",gsl_einval); 00400 } 00401 00402 this->last_conv=0; 00403 int ret=0,status,iter=0; 00404 00405 allocate(nn); 00406 00407 alloc_vec_t ss; 00408 ao.allocate(ss,nn); 00409 for (size_t is=0;is<nn;is++) ss[is]=step_vec[is % step_vec.size()]; 00410 ret=set(ufunc,pa,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->tolx, 00426 "gsl_mmin_simp2"); 00427 } 00428 00429 status=gsl_multimin_test_size(size,this->tolx); 00430 00431 } while(status == GSL_CONTINUE && iter<this->ntrial); 00432 00433 for (size_t i=0;i<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 O2SCL_CONV_RET("Exceeded max # of iterations in mmin_simp2::mmin().", 00442 gsl_emaxiter,this->err_nonconv); 00443 } 00444 00445 return status; 00446 } 00447 00448 /** \brief Calculate the minimum \c min of \c func w.r.t the 00449 array \c x of size \c nvar, using \c xx and \c xx2 to specify 00450 the simplex 00451 */ 00452 virtual int mmin_twovec(size_t nn, vec_t &xx, vec_t &xx2, double &fmin, 00453 param_t &pa, func_t &ufunc) { 00454 00455 this->last_conv=0; 00456 int ret=0,i,status,iter=0; 00457 00458 allocate(nn); 00459 00460 alloc_vec_t ss; 00461 ao.allocate(ss,nn); 00462 00463 for (size_t is=0;is<nn;is++) ss[is]=xx2[is]-xx[is]; 00464 ret=set(ufunc,pa,nn,xx,ss); 00465 ao.free(ss); 00466 00467 if(ret!=0) { 00468 free(); 00469 return ret; 00470 } 00471 00472 do { 00473 iter++; 00474 00475 status=iterate(); 00476 if(status) break; 00477 00478 if(this->verbose>0) { 00479 print_iter(nn,x,x1,fval,iter,size,this->tolx, 00480 "gsl_mmin_simp2"); 00481 } 00482 00483 status=gsl_multimin_test_size(size,this->tolx); 00484 00485 } while(status == GSL_CONTINUE && iter<this->ntrial); 00486 00487 for (i=0;i<((int)nn);i++) xx[i]=x[i]; 00488 fmin=fval; 00489 00490 free(); 00491 this->last_ntrial=iter; 00492 00493 if(iter>=this->ntrial) { 00494 this->last_conv=gsl_emaxiter; 00495 O2SCL_CONV2_RET("Exceeded max number of iterations in ", 00496 "gsl_mmin_simp2::mmin_twovec().", 00497 gsl_emaxiter,this->err_nonconv); 00498 } 00499 00500 return status; 00501 } 00502 00503 /** \brief Calculate the minimum \c min of \c func w.r.t the 00504 array \c x of size \c nvar, given an initial simplex 00505 */ 00506 template<class mat_t> 00507 int mmin_simplex(size_t nn, mat_t &sx, double &fmin, 00508 param_t &pa, func_t &ufunc) { 00509 00510 int ret=0,i,status,iter=0; 00511 00512 allocate(nn); 00513 00514 ret=set_simplex(ufunc,pa,sx); 00515 if(ret!=0) { 00516 free(); 00517 return ret; 00518 } 00519 00520 do { 00521 iter++; 00522 00523 status=iterate(); 00524 if(status) break; 00525 00526 if(this->verbose>0) { 00527 print_iter(nn,x,x1,fval,iter,size,this->tolx, 00528 "gsl_mmin_simp2"); 00529 } 00530 00531 status=gsl_multimin_test_size(size,this->tolx); 00532 00533 } while(status == GSL_CONTINUE && iter<this->ntrial); 00534 00535 for (i=0;i<((int)nn);i++) sx[0][i]=x[i]; 00536 fmin=fval; 00537 00538 free(); 00539 this->last_ntrial=iter; 00540 00541 if(iter>=this->ntrial) { 00542 O2SCL_ERR_RET 00543 ("Exceeded max # of iterations in gsl_mmin_simp2::mmin_twovec().", 00544 gsl_emaxiter); 00545 } 00546 00547 return status; 00548 } 00549 00550 /// Allocate the memory 00551 virtual int allocate(size_t n) { 00552 int status; 00553 if(dim!=0) free(); 00554 set_called=false; 00555 dim=n; 00556 00557 ao.allocate(x,n); 00558 x1=new alloc_vec_t[n+1]; 00559 for(size_t i=0;i<n+1;i++) ao.allocate(x1[i],n); 00560 y1.allocate(n+1); 00561 ao.allocate(ws1,n); 00562 ao.allocate(ws2,n); 00563 ao.allocate(ws3,n); 00564 ao.allocate(center,n); 00565 ao.allocate(delta,n); 00566 ao.allocate(xmc,n); 00567 00568 return gsl_success; 00569 } 00570 00571 /// Free the allocated memory 00572 virtual int free() { 00573 00574 if (dim>0) { 00575 for(size_t i=0;i<dim+1;i++) ao.free(x1[i]); 00576 delete[] x1; 00577 y1.free(); 00578 ao.free(ws1); 00579 ao.free(ws2); 00580 ao.free(ws3); 00581 ao.free(center); 00582 ao.free(delta); 00583 ao.free(xmc); 00584 ao.free(x); 00585 } 00586 00587 dim=0; 00588 00589 return 0; 00590 } 00591 00592 /// Set the function and initial guess 00593 virtual int set(func_t &ufunc, param_t &pa, size_t n, vec_t &ax, 00594 vec_t &step_size) { 00595 size_t i; 00596 00597 if(dim!=n) allocate(n); 00598 params=&pa; 00599 func=&ufunc; 00600 00601 // Copy initial guess to x 00602 for (i=0;i<dim;i++) x[i]=ax[i]; 00603 00604 // first point is the original x0 00605 00606 int ret=ufunc(dim,ax,y1[0],pa); 00607 if (!finite(y1[0])) { 00608 O2SCL_ERR2_RET("Function not finite in ", 00609 "gsl_mmin_simp2::set().",gsl_ebadfunc); 00610 } 00611 for(i=0;i<dim;i++) x1[0][i]=ax[i]; 00612 00613 /* following points are initialized to x0+step_size */ 00614 00615 for (i=1;i<dim+1;i++) { 00616 for(size_t j=0;j<dim;j++) x1[i][j]=x[j]; 00617 x1[i][i-1]=x1[i][i-1]+step_size[i-1]; 00618 ufunc(dim,x1[i],y1[i],pa); 00619 } 00620 00621 /* Initialize simplex size */ 00622 00623 compute_center(); 00624 size=compute_size(); 00625 00626 set_called=true; 00627 00628 return gsl_success; 00629 } 00630 00631 /// Set the function and initial simplex 00632 template<class mat_t> 00633 int set_simplex(func_t &ufunc, param_t &pa, mat_t &sx) { 00634 00635 if(dim==0) { 00636 O2SCL_ERR_RET("Memory not allocated in gsl_mmin_simp2::set().", 00637 gsl_ebadlen); 00638 } 00639 00640 params=&pa; 00641 func=&ufunc; 00642 00643 for(size_t i=0;i<dim+1;i++) { 00644 for(size_t j=0;j<dim;j++) { 00645 x1[i][j]=sx[i][j]; 00646 } 00647 ufunc(dim,x1[i],y1[i],pa); 00648 if (!finite(y1[i])) { 00649 O2SCL_ERR2_RET("Function not finite in ", 00650 "gsl_mmin_simp2::set_simplex().",gsl_ebadfunc); 00651 } 00652 } 00653 00654 /* Initialize simplex size */ 00655 00656 compute_center(); 00657 size=compute_size(); 00658 00659 set_called=true; 00660 00661 return gsl_success; 00662 } 00663 00664 /// Perform an iteration 00665 virtual int iterate() { 00666 00667 size_t n=dim+1; 00668 size_t i; 00669 size_t lo, hi, s_hi; 00670 double dhi,ds_hi,dlo; 00671 int status; 00672 double val,val2; 00673 00674 /* get index of highest, second highest and lowest point */ 00675 00676 lo=hi=0; 00677 dlo=dhi=y1[0]; 00678 s_hi=1; 00679 ds_hi=y1[1]; 00680 00681 for (i=1;i<n;i++) { 00682 val=y1[i]; 00683 if(val<dlo) { 00684 dlo=val; 00685 lo=i; 00686 } else if (val > dhi) { 00687 ds_hi=dhi; 00688 s_hi=hi; 00689 dhi=val; 00690 hi=i; 00691 } else if(val > ds_hi) { 00692 ds_hi=val; 00693 s_hi=i; 00694 } 00695 } 00696 00697 /* reflect the highest value */ 00698 00699 int ret1=try_corner_move(-1.0,hi,ws1,*func,dim,*params,val); 00700 00701 if (avoid_nonzero && ret1!=0) { 00702 std::cout << "Found problem move1: " << std::endl; 00703 for (size_t it=0;it<20 && ret1!=0;it++) { 00704 ret1=try_corner_move(-1.0+((double)it)/10.0,hi,ws1, 00705 *func,dim,*params,val); 00706 std::cout << "it,ret: " << it << " " << ret1 << std::endl; 00707 } 00708 if (ret1!=0) { 00709 O2SCL_ERR2_RET("Failed to move corner (1) in ", 00710 "gsl_mmin_simp2::iterate().",gsl_efailed); 00711 } 00712 } 00713 00714 if (gsl_finite(val) && val<y1[lo]) { 00715 00716 /* reflected point becomes lowest point,try expansion */ 00717 00718 int ret2=try_corner_move(-2.0,hi,ws2,*func,dim,*params,val2); 00719 00720 if (avoid_nonzero && ret2!=0) { 00721 std::cout << "Found problem move2: " << std::endl; 00722 for (size_t it=0;it<20 && ret2!=0;it++) { 00723 ret2=try_corner_move(-2.0+((double)it)/6.7,hi,ws2, 00724 *func,dim,*params,val2); 00725 std::cout << "it,ret: " << it << " " << ret2 << std::endl; 00726 } 00727 if (ret2!=0) { 00728 O2SCL_ERR2_RET("Failed to move corner (2) in ", 00729 "gsl_mmin_simp2::iterate().",gsl_efailed); 00730 } 00731 } 00732 00733 if (gsl_finite(val2) && val2<y1[lo]) { 00734 update_point(hi,ws2,val2); 00735 } else { 00736 update_point(hi,ws1,val); 00737 } 00738 00739 } else if (!gsl_finite(val) || val > y1[s_hi]) { 00740 00741 /* reflection does not improve things enough */ 00742 00743 if (gsl_finite(val) && val <= y1[hi]) { 00744 00745 /* if trial point is better than highest point,replace 00746 highest point */ 00747 00748 update_point(hi,ws1,val); 00749 } 00750 00751 /* try one dimensional contraction */ 00752 00753 int ret3=try_corner_move(0.5,hi,ws2,*func,dim,*params,val2); 00754 00755 if (avoid_nonzero && ret3!=0) { 00756 std::cout << "Found problem move3: " << std::endl; 00757 for (size_t it=0;it<20 && ret3!=0;it++) { 00758 ret3=try_corner_move(0.025*((double)(19-it)),hi,ws2, 00759 *func,dim,*params,val2); 00760 std::cout << "it,ret: " << it << " " << ret3 << std::endl; 00761 } 00762 if (ret3!=0) { 00763 O2SCL_ERR2_RET("Failed to move corner (2) in ", 00764 "gsl_mmin_simp2::iterate().",gsl_efailed); 00765 } 00766 } 00767 00768 if (gsl_finite(val2) && val2 <= y1[hi]) { 00769 00770 update_point(hi,ws2,val2); 00771 00772 } else { 00773 00774 /* contract the whole simplex in respect to the best point */ 00775 status=contract_by_best(lo,*func,dim,*params); 00776 if(status != 0) { 00777 O2SCL_ERR("Function contract_by_best() failed in iterate().", 00778 gsl_efailed); 00779 } 00780 00781 } 00782 00783 } else { 00784 00785 /* trial point is better than second highest point. 00786 Replace highest point by it */ 00787 00788 update_point(hi,ws1,val); 00789 } 00790 00791 /* return lowest point of simplex as x */ 00792 00793 lo=vector_min_index(dim+1,y1,val); 00794 for(i=0;i<dim;i++) x[i]=x1[lo][i]; 00795 fval=y1[lo]; 00796 00797 /* Update simplex size */ 00798 00799 if (S2 > 0) { 00800 size=sqrt(S2); 00801 } else { 00802 /* recompute if accumulated error has made size invalid */ 00803 size=compute_size(); 00804 } 00805 00806 return gsl_success; 00807 } 00808 00809 /** 00810 \brief Print out iteration information. 00811 00812 Depending on the value of the variable verbose, this prints out 00813 the iteration information. If verbose=0, then no information is 00814 printed, while if verbose>1, then after each iteration, the 00815 present values of x and y are output to std::cout along with the 00816 iteration number. If verbose>=2 then each iteration waits for a 00817 character. 00818 */ 00819 virtual int print_iter(size_t nv, vec_t &xx, alloc_vec_t *simp, 00820 double y, int iter, 00821 double value, double limit, 00822 std::string comment) { 00823 00824 if (this->verbose<=0) return 0; 00825 00826 size_t i; 00827 char ch; 00828 00829 (*this->outs) << comment << " Iteration: " << iter << std::endl; 00830 text_out_file outsf(this->outs,79); 00831 outsf.word_out("x:"); 00832 for(i=0;i<nv;i++) outsf.double_out(xx[i]); 00833 outsf.end_line(); 00834 if (print_simplex>0) { 00835 if (print_simplex==2) { 00836 (*this->outs) << "Simplex Values:" << std::endl; 00837 for(i=0;i<nv+1;i++) { 00838 outsf.int_out(i); 00839 outsf.char_out(':'); 00840 for(size_t j=0;j<nv;j++) { 00841 outsf.double_out(simp[i][j]); 00842 } 00843 outsf.char_out(':'); 00844 outsf.double_out(y1[i]); 00845 outsf.end_line(); 00846 } 00847 } else { 00848 outsf.word_out("Simplex Values:"); 00849 for(i=0;i<nv+1;i++) { 00850 outsf.double_out(y1[i]); 00851 } 00852 outsf.end_line(); 00853 } 00854 } 00855 (*this->outs) << "y: " << y << " Val: " << value << " Lim: " 00856 << limit << std::endl; 00857 if (this->verbose>1) { 00858 (*this->outs) << "Press a key and type enter to continue. "; 00859 (*this->ins) >> ch; 00860 } 00861 00862 return 0; 00863 } 00864 00865 /// Return string denoting type("gsl_mmin_simp2") 00866 virtual const char *type() { return "gsl_mmin_simp2";} 00867 00868 }; 00869 00870 #ifndef DOXYGENP 00871 } 00872 #endif 00873 00874 #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