![]() |
Object-oriented Scientific Computing Library: Version 0.910
|
00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006-2012, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 O2scl is free software; you can redistribute it and/or modify 00009 it under the terms of the GNU General Public License as published by 00010 the Free Software Foundation; either version 3 of the License, or 00011 (at your option) any later version. 00012 00013 O2scl is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 00018 You should have received a copy of the GNU General Public License 00019 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00020 00021 ------------------------------------------------------------------- 00022 */ 00023 #ifndef O2SCL_GSL_MMIN_BFGS2_H 00024 #define O2SCL_GSL_MMIN_BFGS2_H 00025 00026 #include <gsl/gsl_blas.h> 00027 #include <gsl/gsl_poly.h> 00028 #include <gsl/gsl_multimin.h> 00029 #include <o2scl/multi_min.h> 00030 #include <o2scl/cblas.h> 00031 00032 #ifndef DOXYGENP 00033 namespace o2scl { 00034 #endif 00035 00036 /** \brief Virtual base for the gsl_mmin_bfgs2 wrapper 00037 00038 This class is useful so that the gsl_mmin_linmin class doesn't 00039 need to depend on any template parameters, even though it will 00040 need a wrapping object as an argument for the 00041 gsl_mmin_linmin::minimize() function. 00042 */ 00043 class gsl_mmin_wrap_base { 00044 public: 00045 virtual ~gsl_mmin_wrap_base() {} 00046 /// Function 00047 virtual double wrap_f(double alpha)=0; 00048 /// Derivative 00049 virtual double wrap_df(double alpha)=0; 00050 /// Function and derivative 00051 virtual void wrap_fdf(double alpha, double *f, double *df)=0; 00052 }; 00053 00054 /** \brief Wrapper class for the gsl_mmin_bfgs2 minimizer 00055 00056 This is a reimplmentation of the internal GSL wrapper for 00057 function calls in the BFGS minimizer. 00058 */ 00059 template<class func_t, class vec_t=ovector_base, 00060 class alloc_vec_t=ovector, class alloc_t=ovector_alloc, 00061 class dfunc_t=grad_funct<ovector_base>, 00062 class auto_grad_t=gradient<func_t,ovector_base> > 00063 class gsl_mmin_wrapper : public gsl_mmin_wrap_base { 00064 00065 #ifndef DOXYGEN_INTERNAL 00066 00067 protected: 00068 00069 /// Function 00070 func_t *func; 00071 00072 /// Derivative 00073 dfunc_t *dfunc; 00074 00075 /// The automatic gradient object 00076 auto_grad_t *agrad; 00077 00078 /// True if the gradient was given by the user 00079 bool grad_given; 00080 00081 /** \name fixed values 00082 */ 00083 //@{ 00084 vec_t *x; 00085 vec_t *g; 00086 vec_t *p; 00087 //@} 00088 00089 /** \name cached values, for x(alpha)=x + alpha * p 00090 */ 00091 //@{ 00092 double f_alpha; 00093 double df_alpha; 00094 //@} 00095 00096 /** \name cache keys 00097 */ 00098 //@{ 00099 double f_cache_key; 00100 double df_cache_key; 00101 double x_cache_key; 00102 double g_cache_key; 00103 //@} 00104 00105 /// Move to a new point, using the cached value if possible 00106 void moveto(double alpha) { 00107 00108 if (alpha == x_cache_key) { 00109 /* using previously cached position */ 00110 return; 00111 } 00112 00113 /* set x_alpha=x + alpha * p */ 00114 00115 for(size_t i=0;i<dim;i++) { 00116 av_x_alpha[i]=(*x)[i]; 00117 } 00118 o2scl_cblas::daxpy(alpha,dim,*p,av_x_alpha); 00119 00120 x_cache_key=alpha; 00121 } 00122 00123 /// Compute the slope 00124 double slope() { 00125 return o2scl_cblas::ddot(dim,av_g_alpha,*p); 00126 } 00127 00128 /// Evaluate the function 00129 virtual double wrap_f(double alpha) { 00130 if (alpha==f_cache_key) { 00131 return f_alpha; 00132 } 00133 moveto(alpha); 00134 f_alpha=(*func)(dim,av_x_alpha); 00135 00136 f_cache_key=alpha; 00137 return f_alpha; 00138 } 00139 00140 /// Evaluate the derivative 00141 virtual double wrap_df(double alpha) { 00142 00143 /* using previously cached df(alpha) */ 00144 if (alpha==df_cache_key) return df_alpha; 00145 00146 moveto(alpha); 00147 if (alpha!=g_cache_key) { 00148 if (grad_given) { 00149 (*dfunc)(dim,av_x_alpha,av_g_alpha); 00150 } else { 00151 (*agrad)(dim,av_x_alpha,av_g_alpha); 00152 } 00153 g_cache_key=alpha; 00154 } 00155 df_alpha=slope(); 00156 df_cache_key=alpha; 00157 00158 return df_alpha; 00159 } 00160 00161 /// Evaluate the function and the derivative 00162 virtual void wrap_fdf(double alpha, double *f, double *df) { 00163 00164 /* Check for previously cached values */ 00165 if (alpha == f_cache_key && alpha == df_cache_key) { 00166 *f=f_alpha; 00167 *df=df_alpha; 00168 return; 00169 } 00170 if (alpha == f_cache_key || alpha == df_cache_key) { 00171 *f=wrap_f(alpha); 00172 *df=wrap_df(alpha); 00173 return; 00174 } 00175 00176 moveto(alpha); 00177 f_alpha=(*func)(dim,av_x_alpha); 00178 if (grad_given) { 00179 (*dfunc)(dim,av_x_alpha,av_g_alpha); 00180 } else { 00181 (*agrad)(dim,av_x_alpha,av_g_alpha); 00182 } 00183 f_cache_key=alpha; 00184 g_cache_key=alpha; 00185 00186 df_alpha=slope(); 00187 df_cache_key=alpha; 00188 00189 *f=f_alpha; 00190 *df=df_alpha; 00191 00192 return; 00193 } 00194 00195 #endif 00196 00197 public: 00198 00199 /// Temporary storage 00200 alloc_vec_t av_x_alpha; 00201 00202 /// Temporary storage 00203 alloc_vec_t av_g_alpha; 00204 00205 /// Number of minimization dimensions 00206 size_t dim; 00207 00208 /// Initialize wrapper 00209 void prepare_wrapper(func_t &ufunc, dfunc_t *udfunc, vec_t &t_x, 00210 double f, vec_t &t_g, vec_t &t_p, auto_grad_t *ag) { 00211 00212 func=&ufunc; 00213 dfunc=udfunc; 00214 if (dfunc!=0) grad_given=true; 00215 else grad_given=false; 00216 agrad=ag; 00217 00218 x=&t_x; 00219 g=&t_g; 00220 p=&t_p; 00221 00222 x_cache_key=0.0; 00223 f_cache_key=0.0; 00224 g_cache_key=0.0; 00225 df_cache_key=0.0; 00226 00227 for(size_t i=0;i<dim;i++) { 00228 av_x_alpha[i]=(*x)[i]; 00229 av_g_alpha[i]=(*g)[i]; 00230 } 00231 00232 f_alpha=f; 00233 df_alpha=slope(); 00234 00235 return; 00236 } 00237 00238 /// Update position 00239 void update_position(double alpha, vec_t &t_x, double *t_f, vec_t &t_g) { 00240 00241 /* ensure that everything is fully cached */ 00242 { 00243 double t_f_alpha, t_df_alpha; 00244 wrap_fdf(alpha,&t_f_alpha,&t_df_alpha); 00245 } 00246 00247 *t_f=f_alpha; 00248 for(size_t i=0;i<dim;i++) { 00249 t_x[i]=av_x_alpha[i]; 00250 t_g[i]=av_g_alpha[i]; 00251 } 00252 } 00253 00254 /** \brief Convert cache values to the new minimizer direction 00255 00256 Convert the cache values from the end of the current minimisation 00257 to those needed for the start of the next minimisation, alpha=0 00258 */ 00259 void change_direction() { 00260 00261 /* The new x_alpha for alpha=0 is the current position and the 00262 new g_alpha for alpha=0 is the current gradient at the 00263 endpoint 00264 */ 00265 for(size_t i=0;i<dim;i++) { 00266 av_x_alpha[i]=(*x)[i]; 00267 av_g_alpha[i]=(*g)[i]; 00268 } 00269 x_cache_key=0.0; 00270 g_cache_key=0.0; 00271 00272 /* The function value does not change */ 00273 f_cache_key=0.0; 00274 00275 /* Calculate the slope along the new direction vector, p */ 00276 df_alpha=slope(); 00277 df_cache_key=0.0; 00278 00279 return; 00280 } 00281 00282 }; 00283 00284 /** \brief The line minimizer for gsl_mmin_bfgs2 00285 */ 00286 class gsl_mmin_linmin2 { 00287 00288 #ifndef DOXYGEN_INTERNAL 00289 00290 protected: 00291 00292 /** \brief Minimize the interpolating quadratic 00293 00294 Find a minimum in x=[0,1] of the interpolating quadratic through 00295 (0,f0) (1,f1) with derivative fp0 at x=0. The interpolating 00296 polynomial is q(x)=f0 + fp0 * z + (f1-f0-fp0) * z^2 00297 */ 00298 double interp_quad(double f0, double fp0, double f1, double zl, 00299 double zh); 00300 00301 /** \brief Minimize the interpolating cubic 00302 00303 Find a minimum in x=[0,1] of the interpolating cubic through 00304 (0,f0) (1,f1) with derivatives fp0 at x=0 and fp1 at x=1. 00305 00306 The interpolating polynomial is: 00307 00308 c(x)=f0 + fp0 * z + eta * z^2 + xi * z^3 00309 00310 where eta=3*(f1-f0)-2*fp0-fp1, xi=fp0+fp1-2*(f1-f0). 00311 */ 00312 double cubic(double c0, double c1, double c2, double c3, double z); 00313 00314 /// Test to see curvature is positive 00315 void check_extremum(double c0, double c1, double c2, double c3, double z, 00316 double *zmin, double *fmin); 00317 00318 /// Interpolate using a cubic 00319 double interp_cubic(double f0, double fp0, double f1, 00320 double fp1, double zl, double zh); 00321 00322 /// Perform the interpolation 00323 double interpolate(double a, double fa, double fpa, double b, 00324 double fb, double fpb, double xmin, double xmax, 00325 int order); 00326 #endif 00327 00328 public: 00329 00330 /** \brief The line minimization 00331 00332 Recommended values from \ref Fletcher87 are 00333 <tt>rho=0.01, sigma=0.1, tau1=9, tau2=0.05, 00334 tau3=0.5 </tt> 00335 */ 00336 int minimize(gsl_mmin_wrap_base &wrap, double rho, double sigma, 00337 double tau1, double tau2, double tau3, 00338 int order, double alpha1, double *alpha_new); 00339 }; 00340 00341 /** \brief Multidimensional minimization by the BFGS 00342 algorithm (GSL) 00343 00344 The functions mmin() and mmin_de() minimize a given function 00345 until the gradient is smaller than the value of multi_min::tol_rel 00346 (which defaults to \f$ 10^{-4} \f$ ). 00347 00348 See an example for the usage of this class in 00349 \ref ex_mmin_sect . 00350 00351 This class includes the optimizations from the GSL minimizer \c 00352 vector_bfgs2. 00353 00354 Default template arguments 00355 - \c func_t - \ref multi_funct< \ref ovector_base > 00356 - \c vec_t - \ref ovector_base 00357 - \c alloc_vec_t - \ref ovector 00358 - \c alloc_t - \ref ovector_alloc 00359 - \c dfunc_t - \ref mm_funct< \ref ovector_base > 00360 - \c auto_grad_t - \ref gradient<func_t, \ref ovector_base > 00361 - \c def_auto_grad_t - \ref simple_grad<func_t,\ref ovector_base > 00362 00363 \todo While BFGS does well in the \c ex_mmin example with the 00364 initial guess of \f$ (1,0,7\pi) \f$ it seems to converge more 00365 poorly for the spring function than the other minimizers with 00366 other initial guesses, and I think this will happen in the GSL 00367 versions too. I need to examine this more closely with some code 00368 designed to clearly show this. 00369 */ 00370 template<class func_t=multi_funct<>, 00371 class vec_t=ovector_base, class alloc_vec_t=ovector, 00372 class alloc_t=ovector_alloc, class dfunc_t=grad_funct<ovector_base>, 00373 class auto_grad_t=gradient<func_t,ovector_base>, 00374 class def_auto_grad_t=simple_grad<func_t,ovector_base> > 00375 class gsl_mmin_bfgs2 : public multi_min<func_t,func_t,vec_t> { 00376 00377 #ifndef DOXYGEN_INTERNAL 00378 00379 protected: 00380 00381 /// \name The original variables from the GSL state structure 00382 //@{ 00383 int iter; 00384 double step; 00385 double g0norm; 00386 double pnorm; 00387 double delta_f; 00388 /* f'(0) for f(x-alpha*p) */ 00389 double fp0; 00390 alloc_vec_t x0; 00391 alloc_vec_t g0; 00392 alloc_vec_t p; 00393 /* work space */ 00394 alloc_vec_t dx0; 00395 alloc_vec_t dg0; 00396 /* wrapper function */ 00397 gsl_mmin_wrapper<func_t,vec_t,alloc_vec_t,alloc_t, 00398 dfunc_t,auto_grad_t> wrap; 00399 /* minimization parameters */ 00400 double rho; 00401 double sigma; 00402 double tau1; 00403 double tau2; 00404 double tau3; 00405 int order; 00406 //@} 00407 00408 /// The line minimizer 00409 gsl_mmin_linmin2 lm; 00410 00411 /// \name Store the arguments to set() so we can use them for iterate() 00412 //@{ 00413 vec_t *st_x; 00414 alloc_vec_t st_dx; 00415 alloc_vec_t st_grad; 00416 double st_f; 00417 //@} 00418 00419 /// Memory size 00420 size_t dim; 00421 00422 /// Memory allocation 00423 alloc_t ao; 00424 00425 /// Automatic gradient object 00426 auto_grad_t *agrad; 00427 00428 #endif 00429 00430 public: 00431 00432 gsl_mmin_bfgs2() { 00433 // The default values of lmin_tol and step_size from the 00434 // example in the GSL documentation 00435 lmin_tol=1.0e-4; 00436 step_size=0.01; 00437 agrad=&def_grad; 00438 } 00439 00440 virtual ~gsl_mmin_bfgs2() {} 00441 00442 /// Perform an iteration 00443 virtual int iterate() { 00444 00445 this->last_conv=0; 00446 00447 double alpha=0.0, alpha1; 00448 00449 double pg, dir; 00450 int status; 00451 00452 double f0=st_f; 00453 00454 if (pnorm == 0.0 || g0norm == 0.0 || fp0 == 0) { 00455 for(size_t i=0;i<this->dim;i++) st_dx[i]=0.0; 00456 this->last_conv=gsl_enoprog; 00457 O2SCL_CONV2("Either pnorm, g0norm, or fp0 vanished in ", 00458 "gsl_mmin_bfgs2::iterate().", 00459 gsl_enoprog,this->err_nonconv); 00460 // The functions mmin() and mmin_de() use this return value, so 00461 // when err_nonconv is false, we still want to return a non-zero 00462 // value. 00463 return gsl_enoprog; 00464 } 00465 00466 if (delta_f < 0) { 00467 double del=GSL_MAX_DBL(-delta_f, 10 * GSL_DBL_EPSILON * fabs(f0)); 00468 alpha1=GSL_MIN_DBL(1.0, 2.0 * del / (- fp0)); 00469 } else { 00470 alpha1=fabs(step); 00471 } 00472 00473 /* line minimisation, with cubic interpolation (order=3) */ 00474 00475 status=lm.minimize(wrap,rho,sigma,tau1,tau2,tau3,order, 00476 alpha1,&alpha); 00477 00478 if (status != gsl_success) { 00479 this->last_conv=gsl_enoprog; 00480 O2SCL_CONV2("Variable stepb vanished in ", 00481 "gsl_mmin_bfgs2::iterate().", 00482 gsl_enoprog,this->err_nonconv); 00483 // The functions mmin() and mmin_de() use this return value, so 00484 // when err_nonconv is false, we still want to return a non-zero 00485 // value. 00486 return gsl_enoprog; 00487 } 00488 00489 wrap.update_position(alpha,*st_x,&st_f,st_grad); 00490 00491 delta_f=st_f - f0; 00492 00493 /* Choose a new direction for the next step */ 00494 00495 { 00496 /* This is the BFGS update: */ 00497 /* p'=g1 - A dx - B dg */ 00498 /* A=- (1+ dg.dg/dx.dg) B + dg.g/dx.dg */ 00499 /* B=dx.g/dx.dg */ 00500 00501 double dxg, dgg, dxdg, dgnorm, A, B; 00502 00503 /* dx0=x - x0 */ 00504 for(size_t i=0;i<dim;i++) dx0[i]=(*st_x)[i]; 00505 o2scl_cblas::daxpy(-1.0,dim,x0,dx0); 00506 00507 /* keep a copy */ 00508 for(size_t i=0;i<dim;i++) st_dx[i]=dx0[i]; 00509 00510 /* dg0=g - g0 */ 00511 for(size_t i=0;i<dim;i++) dg0[i]=st_grad[i]; 00512 o2scl_cblas::daxpy(-1.0,dim,g0,dg0); 00513 00514 dxg=o2scl_cblas::ddot(dim,dx0,st_grad); 00515 dgg=o2scl_cblas::ddot(dim,dg0,st_grad); 00516 dxdg=o2scl_cblas::ddot(dim,dx0,dg0); 00517 00518 dgnorm=o2scl_cblas::dnrm2(dim,dg0); 00519 00520 if (dxdg != 0) { 00521 B=dxg / dxdg; 00522 A=-(1.0 + dgnorm * dgnorm / dxdg) * B + dgg / dxdg; 00523 } else { 00524 B=0; 00525 A=0; 00526 } 00527 00528 for(size_t i=0;i<dim;i++) p[i]=st_grad[i]; 00529 o2scl_cblas::daxpy(-A,dim,dx0,p); 00530 o2scl_cblas::daxpy(-B,dim,dg0,p); 00531 } 00532 00533 for(size_t i=0;i<dim;i++) { 00534 g0[i]=st_grad[i]; 00535 x0[i]=(*st_x)[i]; 00536 } 00537 00538 g0norm=o2scl_cblas::dnrm2(dim,g0); 00539 pnorm=o2scl_cblas::dnrm2(dim,p); 00540 00541 /* update direction and fp0 */ 00542 00543 pg=o2scl_cblas::ddot(dim,st_grad,p); 00544 dir=(pg >= 0.0) ? -1.0 : +1.0; 00545 o2scl_cblas::dscal(dir/pnorm,dim,p); 00546 pnorm=o2scl_cblas::dnrm2(dim,p); 00547 fp0=o2scl_cblas::ddot(dim,p,g0); 00548 00549 wrap.change_direction(); 00550 00551 return gsl_success; 00552 } 00553 00554 /// Return string denoting type("gsl_mmin_bfgs2") 00555 virtual const char *type() { return "gsl_mmin_bfgs2";} 00556 00557 /// Allocate the memory 00558 virtual int allocate(size_t n) { 00559 00560 ao.allocate(p,n); 00561 ao.allocate(x0,n); 00562 ao.allocate(g0,n); 00563 ao.allocate(dx0,n); 00564 ao.allocate(dg0,n); 00565 00566 for(size_t i=0;i<n;i++) { 00567 p[i]=0.0; 00568 x0[i]=0.0; 00569 g0[i]=0.0; 00570 dx0[i]=0.0; 00571 dg0[i]=0.0; 00572 } 00573 00574 ao.allocate(st_dx,n); 00575 ao.allocate(st_grad,n); 00576 00577 ao.allocate(wrap.av_x_alpha,n); 00578 ao.allocate(wrap.av_g_alpha,n); 00579 wrap.dim=n; 00580 dim=n; 00581 00582 return gsl_success; 00583 } 00584 00585 /// Free the allocated memory 00586 virtual int free() { 00587 ao.free(wrap.av_x_alpha); 00588 ao.free(wrap.av_g_alpha); 00589 ao.free(st_grad); 00590 ao.free(dg0); 00591 ao.free(dx0); 00592 ao.free(g0); 00593 ao.free(x0); 00594 ao.free(p); 00595 ao.free(st_dx); 00596 wrap.dim=0; 00597 dim=0; 00598 return 0; 00599 } 00600 00601 /// Reset the minimizer to use the current point as a new starting point 00602 int restart() { 00603 iter=0; 00604 return gsl_success; 00605 } 00606 00607 /// Set the function and initial guess 00608 virtual int set(vec_t &x, double u_step_size, double tol_u, 00609 func_t &ufunc) { 00610 00611 iter=0; 00612 step=u_step_size; 00613 delta_f=0; 00614 00615 st_x=&x; 00616 00617 st_f=ufunc(dim,x); 00618 agrad->set_function(ufunc); 00619 (*agrad)(dim,x,st_grad); 00620 00621 /* Use the gradient as the initial direction */ 00622 00623 for(size_t i=0;i<dim;i++) { 00624 x0[i]=x[i]; 00625 g0[i]=st_grad[i]; 00626 } 00627 g0norm=o2scl_cblas::dnrm2(dim,g0); 00628 for(size_t i=0;i<dim;i++) { 00629 p[i]=st_grad[i]; 00630 } 00631 o2scl_cblas::dscal(-1.0/g0norm,dim,p); 00632 00633 // pnorm should be 1 00634 pnorm=o2scl_cblas::dnrm2(dim,p); 00635 fp0=-g0norm; 00636 00637 /* Prepare the wrapper */ 00638 00639 wrap.prepare_wrapper(ufunc,0,x0,st_f,g0,p,agrad); 00640 00641 /* Prepare 1d minimisation parameters */ 00642 00643 rho=0.01; 00644 sigma=tol_u; 00645 tau1=9; 00646 tau2=0.05; 00647 tau3=0.5; 00648 // use cubic interpolation where possible 00649 order=3; 00650 00651 return gsl_success; 00652 00653 } 00654 00655 /// Set the function, the gradient, and the initial guess 00656 virtual int set_de(vec_t &x, double u_step_size, double tol_u, 00657 func_t &ufunc, dfunc_t &udfunc) { 00658 00659 iter=0; 00660 step=u_step_size; 00661 delta_f=0; 00662 00663 st_x=&x; 00664 00665 st_f=ufunc(dim,x); 00666 udfunc(dim,x,st_grad); 00667 00668 /* Use the gradient as the initial direction */ 00669 00670 for(size_t i=0;i<dim;i++) { 00671 x0[i]=x[i]; 00672 g0[i]=st_grad[i]; 00673 } 00674 g0norm=o2scl_cblas::dnrm2(dim,g0); 00675 for(size_t i=0;i<dim;i++) { 00676 p[i]=st_grad[i]; 00677 } 00678 o2scl_cblas::dscal(-1.0/g0norm,dim,p); 00679 00680 // pnorm should be 1 00681 pnorm=o2scl_cblas::dnrm2(dim,p); 00682 fp0=-g0norm; 00683 00684 /* Prepare the wrapper */ 00685 00686 wrap.prepare_wrapper(ufunc,&udfunc,x0,st_f,g0,p,agrad); 00687 00688 /* Prepare 1d minimisation parameters */ 00689 00690 rho=0.01; 00691 sigma=tol_u; 00692 tau1=9; 00693 tau2=0.05; 00694 tau3=0.5; 00695 // use cubic interpolation where possible 00696 order=3; 00697 00698 return gsl_success; 00699 00700 } 00701 00702 /// The size of the first trial step (default 0.01) 00703 double step_size; 00704 00705 /// The tolerance for the 1-dimensional minimizer 00706 double lmin_tol; 00707 00708 /// Default automatic gradient object 00709 def_auto_grad_t def_grad; 00710 00711 /** \brief Calculate the minimum \c min of \c func w.r.t the 00712 array \c x of size \c nn. 00713 */ 00714 virtual int mmin(size_t nn, vec_t &xx, double &fmin, 00715 func_t &ufunc) { 00716 00717 if (nn==0) { 00718 O2SCL_ERR2_RET("Tried to minimize over zero variables ", 00719 " in gsl_mmin_bfgs2::mmin().",gsl_einval); 00720 } 00721 00722 int xiter=0, status; 00723 00724 allocate(nn); 00725 00726 set(xx,step_size,lmin_tol,ufunc); 00727 00728 do { 00729 00730 xiter++; 00731 00732 status=iterate(); 00733 00734 if (status) { 00735 break; 00736 } 00737 00738 // Equivalent to gsl_multimin_test_gradient with 00739 // additional code to print out present iteration 00740 00741 double norm=o2scl_cblas::dnrm2(nn,st_grad); 00742 00743 if(this->verbose>0) { 00744 this->print_iter(nn,*st_x,st_f,xiter, 00745 norm,this->tol_rel,"gsl_mmin_bfgs2"); 00746 } 00747 00748 if (norm<this->tol_rel) { 00749 status=gsl_success; 00750 } else { 00751 status=gsl_continue; 00752 } 00753 00754 } while (status == gsl_continue && xiter < this->ntrial); 00755 00756 for(size_t i=0;i<nn;i++) xx[i]=(*st_x)[i]; 00757 fmin=st_f; 00758 00759 free(); 00760 this->last_ntrial=xiter; 00761 00762 if (status==gsl_continue && xiter==this->ntrial) { 00763 this->last_conv=gsl_emaxiter; 00764 O2SCL_CONV_RET("Too many iterations in gsl_mmin_bfgs2::mmin().", 00765 gsl_emaxiter,this->err_nonconv); 00766 } 00767 return status; 00768 } 00769 00770 /** \brief Calculate the minimum \c min of \c func w.r.t the 00771 array \c x of size \c nn. 00772 */ 00773 virtual int mmin_de(size_t nn, vec_t &xx, double &fmin, 00774 func_t &ufunc, dfunc_t &udfunc) { 00775 00776 if (nn==0) { 00777 O2SCL_ERR2_RET("Tried to minimize over zero variables ", 00778 "in gsl_mmin_bfgs2::mmin().",gsl_einval); 00779 } 00780 00781 int xiter=0, status; 00782 00783 allocate(nn); 00784 00785 set_de(xx,step_size,lmin_tol,ufunc,udfunc); 00786 00787 do { 00788 xiter++; 00789 00790 status=iterate(); 00791 00792 if (status) { 00793 break; 00794 } 00795 00796 // Equivalent to gsl_multimin_test_gradient with 00797 // additional code to print out present iteration 00798 00799 double norm=o2scl_cblas::dnrm2(nn,st_grad); 00800 00801 if(this->verbose>0) { 00802 this->print_iter(nn,*st_x,st_f,xiter, 00803 norm,this->tol_rel,"gsl_mmin_bfgs2"); 00804 } 00805 00806 if (norm<this->tol_rel) status=gsl_success; 00807 else status=gsl_continue; 00808 00809 } while (status == gsl_continue && xiter < this->ntrial); 00810 00811 for(size_t i=0;i<nn;i++) xx[i]=(*st_x)[i]; 00812 fmin=st_f; 00813 00814 free(); 00815 this->last_ntrial=xiter; 00816 00817 if (status==gsl_continue && xiter==this->ntrial) { 00818 this->last_conv=gsl_emaxiter; 00819 O2SCL_CONV_RET("Too many iterations in gsl_mmin_bfgs2::mmin_de().", 00820 gsl_emaxiter,this->err_nonconv); 00821 } 00822 00823 return status; 00824 } 00825 00826 }; 00827 00828 #ifndef DOXYGENP 00829 } 00830 #endif 00831 00832 #endif
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).