00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 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 00031 #ifndef DOXYGENP 00032 namespace o2scl { 00033 #endif 00034 00035 /** 00036 \brief Virtual base for the gsl_mmin_bfgs2 wrapper 00037 00038 This is useful so that the gsl_mmin_linmin class doesn't need to 00039 depend on any template parameters, even though it will need a 00040 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, void *params)=0; 00048 /// Derivative 00049 virtual double wrap_df(double alpha, void *params)=0; 00050 /// Function and derivative 00051 virtual void wrap_fdf(double alpha, void *params, double *f, double *df)=0; 00052 }; 00053 00054 /** 00055 \brief Wrapper class for the gsl_mmin_bfgs2 minimizer 00056 00057 This is a reimplmentation of the internal GSL wrapper for 00058 function calls in the BFGS minimizer. 00059 00060 \future There's a bit of extra vector copying here which 00061 could potentially be avoided. 00062 */ 00063 template<class param_t, class func_t, 00064 class vec_t=ovector_view, class alloc_vec_t=ovector, 00065 class alloc_t=ovector_alloc, 00066 class dfunc_t=grad_funct<param_t,ovector_view>, 00067 class auto_grad_t=gradient<param_t,func_t,ovector_view>, 00068 class def_auto_grad_t=simple_grad<param_t,func_t,ovector_view> > 00069 class gsl_mmin_wrapper : public gsl_mmin_wrap_base { 00070 00071 #ifndef DOXYGEN_INTERNAL 00072 00073 protected: 00074 00075 /// Function 00076 func_t *func; 00077 00078 /// Derivative 00079 dfunc_t *dfunc; 00080 00081 /// Test 00082 auto_grad_t *agrad; 00083 00084 /// Parameters 00085 param_t *pa; 00086 00087 /** \name fixed values 00088 */ 00089 //@{ 00090 gsl_vector *x; 00091 gsl_vector *g; 00092 gsl_vector *p; 00093 //@} 00094 00095 /** \name cached values, for x(alpha) = x + alpha * p 00096 */ 00097 //@{ 00098 double f_alpha; 00099 double df_alpha; 00100 //@} 00101 00102 /** \name cache keys 00103 */ 00104 //@{ 00105 double f_cache_key; 00106 double df_cache_key; 00107 double x_cache_key; 00108 double g_cache_key; 00109 //@} 00110 00111 /// Move to a new point, using the cached value if possible 00112 void moveto(double alpha) { 00113 00114 if (alpha == x_cache_key) { 00115 /* using previously cached position */ 00116 return; 00117 } 00118 00119 /* set x_alpha = x + alpha * p */ 00120 00121 gsl_vector *tmp=gsl_vector_alloc(dim); 00122 for(size_t i=0;i<dim;i++) { 00123 av_x_alpha[i]=gsl_vector_get(x,i); 00124 gsl_vector_set(tmp,i,av_x_alpha[i]); 00125 } 00126 00127 gsl_blas_daxpy(alpha,p,tmp); 00128 for(size_t i=0;i<dim;i++) { 00129 av_x_alpha[i]=gsl_vector_get(tmp,i); 00130 } 00131 gsl_vector_free(tmp); 00132 00133 x_cache_key = alpha; 00134 } 00135 00136 /// Compute the slope 00137 double slope() { 00138 double df; 00139 gsl_vector *tmp=gsl_vector_alloc(dim); 00140 for(size_t i=0;i<dim;i++) { 00141 gsl_vector_set(tmp,i,av_g_alpha[i]); 00142 } 00143 gsl_blas_ddot(tmp,p,&df); 00144 for(size_t i=0;i<dim;i++) { 00145 av_g_alpha[i]=gsl_vector_get(tmp,i); 00146 } 00147 gsl_vector_free(tmp); 00148 return df; 00149 } 00150 00151 /// Evaluate the function 00152 virtual double wrap_f(double alpha, void *params) { 00153 if (alpha==f_cache_key) { 00154 return f_alpha; 00155 } 00156 moveto(alpha); 00157 (*func)(dim,av_x_alpha,f_alpha,*pa); 00158 00159 f_cache_key = alpha; 00160 return f_alpha; 00161 } 00162 00163 /// Evaluate the derivative 00164 virtual double wrap_df(double alpha, void *params) { 00165 00166 /* using previously cached df(alpha) */ 00167 if (alpha==df_cache_key) return df_alpha; 00168 00169 moveto(alpha); 00170 if (alpha!=g_cache_key) { 00171 simple_df(av_x_alpha,av_g_alpha); 00172 g_cache_key=alpha; 00173 } 00174 df_alpha=slope(); 00175 df_cache_key=alpha; 00176 00177 return df_alpha; 00178 } 00179 00180 /// A simple derivative 00181 int simple_df(vec_t &x2, vec_t &g2) { 00182 00183 double fv1, fv2, deriv_h=1.0e-4; 00184 00185 (*func)(dim,x2,fv1,*pa); 00186 00187 for(size_t i=0;i<dim;i++) { 00188 x2[i]+=deriv_h; 00189 (*func)(dim,x2,fv2,*pa); 00190 x2[i]-=deriv_h; 00191 g2[i]=(fv2-fv1)/deriv_h; 00192 } 00193 00194 return 0; 00195 } 00196 00197 /// Evaluate the function and the derivative 00198 virtual void wrap_fdf(double alpha, void *params, double *f, double *df) { 00199 00200 /* Check for previously cached values */ 00201 if (alpha == f_cache_key && alpha == df_cache_key) { 00202 *f = f_alpha; 00203 *df = df_alpha; 00204 return; 00205 } 00206 if (alpha == f_cache_key || alpha == df_cache_key) { 00207 *f = wrap_f (alpha, params); 00208 *df = wrap_df (alpha, params); 00209 return; 00210 } 00211 00212 moveto(alpha); 00213 (*func)(dim,av_x_alpha,f_alpha,*pa); 00214 simple_df(av_x_alpha,av_g_alpha); 00215 f_cache_key = alpha; 00216 g_cache_key = alpha; 00217 00218 df_alpha = slope(); 00219 df_cache_key = alpha; 00220 00221 *f = f_alpha; 00222 *df = df_alpha; 00223 00224 return; 00225 } 00226 00227 #endif 00228 00229 public: 00230 00231 /// Temporary storage 00232 alloc_vec_t av_x_alpha; 00233 00234 /// Temporary storage 00235 alloc_vec_t av_g_alpha; 00236 00237 /// Number of minimization dimensions 00238 size_t dim; 00239 00240 /// Initialize wrapper 00241 void prepare_wrapper(func_t &ufunc, param_t &upa, gsl_vector *t_x, 00242 double f, gsl_vector *t_g, gsl_vector *t_p) { 00243 00244 func=&ufunc; 00245 pa=&upa; 00246 00247 x=t_x; 00248 g=t_g; 00249 p=t_p; 00250 00251 x_cache_key=0.0; 00252 f_cache_key=0.0; 00253 g_cache_key=0.0; 00254 df_cache_key=0.0; 00255 00256 for(size_t i=0;i<dim;i++) { 00257 av_x_alpha[i]=gsl_vector_get(x,i); 00258 av_g_alpha[i]=gsl_vector_get(g,i); 00259 } 00260 00261 f_alpha=f; 00262 df_alpha=slope(); 00263 00264 return; 00265 } 00266 00267 /// Update position 00268 void update_position(double alpha, vec_t &t_x, double *t_f, 00269 gsl_vector *t_g) { 00270 00271 /* ensure that everything is fully cached */ 00272 { 00273 double t_f_alpha, t_df_alpha; 00274 wrap_fdf(alpha,0,&t_f_alpha,&t_df_alpha); 00275 } 00276 00277 *t_f = f_alpha; 00278 for(size_t i=0;i<dim;i++) { 00279 t_x[i]=av_x_alpha[i]; 00280 gsl_vector_set(t_g,i,av_g_alpha[i]); 00281 } 00282 } 00283 00284 /** 00285 \brief Convert cache values to the new minimizer direction 00286 00287 Convert the cache values from the end of the current minimisation 00288 to those needed for the start of the next minimisation, alpha=0 00289 */ 00290 void change_direction() { 00291 00292 /* The new x_alpha for alpha=0 is the current position and the 00293 new g_alpha for alpha=0 is the current gradient at the 00294 endpoint */ 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 Fletcher are 00367 rho = 0.01, sigma = 0.1, tau1 = 9, tau2 = 0.05, tau3 = 0.5 00368 00369 \todo Properly reference Fletcher here. 00370 */ 00371 int minimize(gsl_mmin_wrap_base &wrap, double rho, double sigma, 00372 double tau1, double tau2, double tau3, 00373 int order, double alpha1, double *alpha_new); 00374 }; 00375 00376 /** \brief Multidimensional minimization by the BFGS 00377 algorithm (GSL) 00378 00379 This class includes the optimizations from the GSL minimizer \c 00380 vector_bfgs2. 00381 00382 \todo Works with generic vector objects, but 00383 doesn't allow specification of jacobian yet. 00384 */ 00385 template<class param_t, class func_t=multi_funct<param_t>, 00386 class vec_t=ovector_view, class alloc_vec_t=ovector, 00387 class alloc_t=ovector_alloc, 00388 class dfunc_t=grad_funct<param_t,ovector_view>, 00389 class auto_grad_t=gradient<param_t,func_t,ovector_view>, 00390 class def_auto_grad_t=simple_grad<param_t,func_t,ovector_view> > 00391 class gsl_mmin_bfgs2 : public multi_min<param_t,func_t,func_t,vec_t> { 00392 00393 #ifndef DOXYGEN_INTERNAL 00394 00395 protected: 00396 00397 /// \name The original variables from the GSL state structure 00398 //@{ 00399 int iter; 00400 double step; 00401 double g0norm; 00402 double pnorm; 00403 double delta_f; 00404 /* f'(0) for f(x-alpha*p) */ 00405 double fp0; 00406 gsl_vector *x0; 00407 gsl_vector *g0; 00408 gsl_vector *p; 00409 /* work space */ 00410 gsl_vector *dx0; 00411 gsl_vector *dg0; 00412 /* wrapper function */ 00413 gsl_mmin_wrapper<param_t,func_t,vec_t,alloc_vec_t,alloc_t, 00414 dfunc_t,auto_grad_t,def_auto_grad_t> wrap; 00415 /* minimization parameters */ 00416 double rho; 00417 double sigma; 00418 double tau1; 00419 double tau2; 00420 double tau3; 00421 int order; 00422 //@} 00423 00424 /// The line minimizer 00425 gsl_mmin_linmin lm; 00426 00427 /// \name Store the arguments to set() so we can use them for iterate() 00428 //@{ 00429 vec_t *st_x; 00430 gsl_vector *st_dx; 00431 gsl_vector *st_grad; 00432 double st_f; 00433 //@} 00434 00435 /// Memory size 00436 size_t dim; 00437 00438 /// Memory allocation 00439 alloc_t ao; 00440 00441 #endif 00442 00443 public: 00444 00445 gsl_mmin_bfgs2() { 00446 lmin_tol=1.0e-4; 00447 this->tolf=1.0e-3; 00448 step_size=0.01; 00449 } 00450 00451 virtual ~gsl_mmin_bfgs2() {} 00452 00453 /// Perform an iteration 00454 virtual int iterate() { 00455 00456 double alpha = 0.0, alpha1; 00457 00458 double pg, dir; 00459 int status; 00460 00461 double f0 = st_f; 00462 00463 if (pnorm == 0.0 || g0norm == 0.0 || fp0 == 0) { 00464 gsl_vector_set_zero (st_dx); 00465 set_err_ret("Not making progress in gsl_mmin_bfgs2::iterate().", 00466 gsl_enoprog); 00467 } 00468 00469 if (delta_f < 0) { 00470 double del = GSL_MAX_DBL (-delta_f, 10 * GSL_DBL_EPSILON * fabs(f0)); 00471 alpha1 = GSL_MIN_DBL (1.0, 2.0 * del / (- fp0)); 00472 } else { 00473 alpha1 = fabs( step); 00474 } 00475 00476 /* line minimisation, with cubic interpolation (order = 3) */ 00477 00478 status=lm.minimize(wrap,rho,sigma,tau1,tau2,tau3,order, 00479 alpha1,&alpha); 00480 00481 if (status != GSL_SUCCESS) { 00482 return status; 00483 } 00484 00485 wrap.update_position(alpha,*st_x,&st_f,st_grad); 00486 00487 delta_f = st_f - f0; 00488 00489 /* Choose a new direction for the next step */ 00490 00491 { 00492 /* This is the BFGS update: */ 00493 /* p' = g1 - A dx - B dg */ 00494 /* A = - (1+ dg.dg/dx.dg) B + dg.g/dx.dg */ 00495 /* B = dx.g/dx.dg */ 00496 00497 double dxg, dgg, dxdg, dgnorm, A, B; 00498 00499 /* dx0 = x - x0 */ 00500 for(size_t i=0;i<dim;i++) gsl_vector_set(dx0,i,(*st_x)[i]); 00501 gsl_blas_daxpy (-1.0, x0, dx0); 00502 00503 gsl_vector_memcpy (st_dx, dx0); /* keep a copy */ 00504 00505 /* dg0 = g - g0 */ 00506 gsl_vector_memcpy (dg0, st_grad); 00507 gsl_blas_daxpy (-1.0, g0, dg0); 00508 00509 gsl_blas_ddot (dx0, st_grad, &dxg); 00510 gsl_blas_ddot (dg0, st_grad, &dgg); 00511 gsl_blas_ddot (dx0, dg0, &dxdg); 00512 00513 dgnorm = gsl_blas_dnrm2 (dg0); 00514 00515 if (dxdg != 0) { 00516 B = dxg / dxdg; 00517 A = -(1.0 + dgnorm * dgnorm / dxdg) * B + dgg / dxdg; 00518 } else { 00519 B = 0; 00520 A = 0; 00521 } 00522 00523 gsl_vector_memcpy (p, st_grad); 00524 gsl_blas_daxpy (-A, dx0, p); 00525 gsl_blas_daxpy (-B, dg0, p); 00526 } 00527 00528 gsl_vector_memcpy (g0, st_grad); 00529 00530 for(size_t i=0;i<dim;i++) gsl_vector_set(x0,i,(*st_x)[i]); 00531 00532 g0norm = gsl_blas_dnrm2 (g0); 00533 pnorm = gsl_blas_dnrm2 (p); 00534 /* update direction and fp0 */ 00535 00536 gsl_blas_ddot (p, st_grad, &pg); 00537 dir = (pg >= 0.0) ? -1.0 : +1.0; 00538 gsl_blas_dscal (dir / pnorm, p); 00539 pnorm = gsl_blas_dnrm2 (p); 00540 gsl_blas_ddot (p, g0, & fp0); 00541 00542 wrap.change_direction(); 00543 00544 return GSL_SUCCESS; 00545 00546 } 00547 00548 /// Return string denoting type("gsl_mmin_bfgs2") 00549 virtual const char *type() { return "gsl_mmin_bfgs2";} 00550 00551 /// Allocate the memory 00552 virtual int allocate(size_t n) { 00553 00554 p=gsl_vector_calloc(n); 00555 if (p == 0) { 00556 set_err_ret("Failed to allocate p in gsl_mmin_bfgs2::allocate().", 00557 gsl_enomem); 00558 } 00559 00560 x0=gsl_vector_calloc(n); 00561 if (x0 == 0) { 00562 gsl_vector_free(p); 00563 set_err_ret("Failed to allocate x0 in gsl_mmin_bfgs2::allocate().", 00564 gsl_enomem); 00565 } 00566 00567 g0=gsl_vector_calloc(n); 00568 if (g0 == 0) { 00569 gsl_vector_free(x0); 00570 gsl_vector_free(p); 00571 set_err_ret("Failed to allocate g0 in gsl_mmin_bfgs2::allocate().", 00572 gsl_enomem); 00573 } 00574 00575 dx0=gsl_vector_calloc(n); 00576 if (dx0 == 0) { 00577 gsl_vector_free(g0); 00578 gsl_vector_free(x0); 00579 gsl_vector_free(p); 00580 set_err_ret("Failed to allocate dx0 in gsl_mmin_bfgs2::allocate().", 00581 gsl_enomem); 00582 } 00583 00584 dg0=gsl_vector_calloc(n); 00585 if (dg0 == 0) { 00586 gsl_vector_free(dx0); 00587 gsl_vector_free(g0); 00588 gsl_vector_free(x0); 00589 gsl_vector_free(p); 00590 set_err_ret("Failed to allocate dg0 in gsl_mmin_bfgs2::allocate().", 00591 gsl_enomem); 00592 } 00593 00594 st_dx=gsl_vector_alloc(n); 00595 st_grad=gsl_vector_alloc(n); 00596 00597 ao.allocate(wrap.av_x_alpha,n); 00598 ao.allocate(wrap.av_g_alpha,n); 00599 wrap.dim=n; 00600 dim=n; 00601 00602 return GSL_SUCCESS; 00603 } 00604 00605 /// Free the allocated memory 00606 virtual int free() { 00607 ao.free(wrap.av_x_alpha); 00608 ao.free(wrap.av_g_alpha); 00609 gsl_vector_free(dg0); 00610 gsl_vector_free(dx0); 00611 gsl_vector_free(g0); 00612 gsl_vector_free(x0); 00613 gsl_vector_free(p); 00614 gsl_vector_free(st_dx); 00615 gsl_vector_free(st_grad); 00616 wrap.dim=0; 00617 dim=0; 00618 return 0; 00619 } 00620 00621 /// Reset the minimizer to use the current point as a new starting point 00622 int restart() { 00623 iter=0; 00624 return gsl_success; 00625 } 00626 00627 /// Set the function and initial guess 00628 virtual int set(vec_t &x, double u_step_size, double tol_u, 00629 func_t &ufunc, param_t &upa) { 00630 00631 iter=0; 00632 step=u_step_size; 00633 delta_f=0; 00634 00635 st_x=&x; 00636 00637 ufunc(dim,x,st_f,upa); 00638 { 00639 double fv2, deriv_h=1.0e-4; 00640 00641 for(size_t i=0;i<dim;i++) { 00642 x[i]+=deriv_h; 00643 ufunc(dim,x,fv2,upa); 00644 x[i]-=deriv_h; 00645 gsl_vector_set(st_grad,i,(fv2-st_f)/deriv_h); 00646 } 00647 } 00648 00649 /* Use the gradient as the initial direction */ 00650 00651 for(size_t i=0;i<dim;i++) { 00652 gsl_vector_set(x0,i,x[i]); 00653 } 00654 gsl_vector_memcpy(g0,st_grad); 00655 g0norm=gsl_blas_dnrm2(g0); 00656 00657 gsl_vector_memcpy(p,st_grad); 00658 gsl_blas_dscal(-1/g0norm,p); 00659 // pnorm should be 1 00660 pnorm=gsl_blas_dnrm2(p); 00661 fp0=- g0norm; 00662 00663 /* Prepare the wrapper */ 00664 00665 wrap.prepare_wrapper(ufunc,upa,x0,st_f,g0,p); 00666 00667 /* Prepare 1d minimisation parameters */ 00668 00669 rho=0.01; 00670 sigma=tol_u; 00671 tau1=9; 00672 tau2=0.05; 00673 tau3=0.5; 00674 // use cubic interpolation where possible 00675 order=3; 00676 00677 return GSL_SUCCESS; 00678 00679 } 00680 00681 /// Set the function, the gradient, and the initial guess 00682 virtual int set_de(vec_t &x, double u_step_size, double tol_u, 00683 func_t &ufunc, dfunc_t &udfunc, param_t &upa) { 00684 00685 iter=0; 00686 step=u_step_size; 00687 delta_f=0; 00688 00689 st_x=&x; 00690 00691 ufunc(dim,x,st_f,upa); 00692 //udfunc(dim,x,st_grad,upa); 00693 00694 /* Use the gradient as the initial direction */ 00695 00696 for(size_t i=0;i<dim;i++) { 00697 gsl_vector_set(x0,i,x[i]); 00698 } 00699 gsl_vector_memcpy(g0,st_grad); 00700 g0norm=gsl_blas_dnrm2(g0); 00701 00702 gsl_vector_memcpy(p,st_grad); 00703 gsl_blas_dscal(-1/g0norm,p); 00704 // pnorm should be 1 00705 pnorm=gsl_blas_dnrm2(p); 00706 fp0=- g0norm; 00707 00708 /* Prepare the wrapper */ 00709 00710 wrap.prepare_wrapper(ufunc,upa,x0,st_f,g0,p); 00711 00712 /* Prepare 1d minimisation parameters */ 00713 00714 rho=0.01; 00715 sigma=tol_u; 00716 tau1=9; 00717 tau2=0.05; 00718 tau3=0.5; 00719 // use cubic interpolation where possible 00720 order=3; 00721 00722 return GSL_SUCCESS; 00723 00724 } 00725 00726 /// The size of the first trial step 00727 double step_size; 00728 00729 /// The tolerance for the 1-dimensional minimizer 00730 double lmin_tol; 00731 00732 /** \brief Calculate the minimum \c min of \c func w.r.t the 00733 array \c x of size \c nn. 00734 */ 00735 virtual int mmin(size_t nn, vec_t &xx, double &fmin, param_t &pa, 00736 func_t &ufunc) { 00737 00738 int xiter=0, status; 00739 00740 allocate(nn); 00741 00742 set(xx,step_size,lmin_tol,ufunc,pa); 00743 00744 do { 00745 xiter++; 00746 00747 status=iterate(); 00748 if (status) { 00749 break; 00750 } 00751 00752 // Equivalent to gsl_multimin_test_gradient with 00753 // additional code to print out present iteration 00754 00755 double norm=gsl_blas_dnrm2(st_grad); 00756 00757 if(this->verbose>0) { 00758 this->print_iter(nn,*st_x,st_f,xiter, 00759 norm,this->tolf,"gsl_mmin_bfgs2"); 00760 } 00761 00762 if (norm<this->tolf) status=gsl_success; 00763 else status=gsl_continue; 00764 00765 } 00766 while (status == GSL_CONTINUE && xiter < this->ntrial); 00767 00768 for(size_t i=0;i<nn;i++) xx[i]=(*st_x)[i]; 00769 00770 fmin=st_f; 00771 00772 free(); 00773 00774 this->last_ntrial=xiter; 00775 00776 return status; 00777 } 00778 00779 /** \brief Calculate the minimum \c min of \c func w.r.t the 00780 array \c x of size \c nn. 00781 */ 00782 virtual int mmin_de(size_t nn, vec_t &xx, double &fmin, param_t &pa, 00783 func_t &ufunc, dfunc_t &udfunc) { 00784 00785 int xiter=0, status; 00786 00787 allocate(nn); 00788 00789 set_de(xx,step_size,lmin_tol,ufunc,udfunc,pa); 00790 00791 do { 00792 xiter++; 00793 00794 status=iterate(); 00795 00796 if (status) { 00797 break; 00798 } 00799 00800 // Equivalent to gsl_multimin_test_gradient with 00801 // additional code to print out present iteration 00802 00803 double norm=gsl_blas_dnrm2(st_grad); 00804 00805 if(this->verbose>0) { 00806 this->print_iter(nn,*st_x,st_f,xiter, 00807 norm,this->tolf,"gsl_mmin_bfgs2"); 00808 } 00809 00810 if (norm<this->tolf) status=gsl_success; 00811 else status=gsl_continue; 00812 00813 } 00814 while (status == GSL_CONTINUE && xiter < this->ntrial); 00815 00816 for(size_t i=0;i<nn;i++) xx[i]=(*st_x)[i]; 00817 00818 fmin=st_f; 00819 00820 free(); 00821 00822 this->last_ntrial=xiter; 00823 00824 return status; 00825 } 00826 00827 }; 00828 00829 #ifndef DOXYGENP 00830 } 00831 #endif 00832 00833 #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