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