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