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