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_CONF_H 00024 #define O2SCL_GSL_MMIN_CONF_H 00025 00026 #include <gsl/gsl_blas.h> 00027 #include <gsl/gsl_multimin.h> 00028 #include <o2scl/multi_min.h> 00029 00030 #ifndef DOXYGENP 00031 namespace o2scl { 00032 #endif 00033 00034 /** 00035 \brief Base minimization routines for gsl_mmin_conf and gsl_mmin_conp 00036 00037 Default template arguments 00038 - \c param_t - \ref multi_funct 00039 - \c vec_t - \ref ovector_view 00040 - \c alloc_vec_t - \ref ovector 00041 - \c alloc_t - \ref ovector_alloc 00042 - \c dfunc_t - \ref grad_funct 00043 - \c auto_grad_t - \ref gradient 00044 - \c def_auto_grad_t - \ref simple_grad 00045 00046 */ 00047 template<class param_t, class func_t = multi_funct<param_t>, 00048 class vec_t = ovector_view, class alloc_vec_t = ovector, 00049 class alloc_t = ovector_alloc, 00050 class dfunc_t = grad_funct<param_t,ovector_view>, 00051 class auto_grad_t = gradient<param_t,func_t,ovector_view>, 00052 class def_auto_grad_t = simple_grad<param_t,func_t,ovector_view> > 00053 class gsl_mmin_base : public multi_min<param_t,func_t,func_t,vec_t> 00054 { 00055 00056 #ifndef DOXYGEN_INTERNAL 00057 00058 protected: 00059 00060 /// User-specified function 00061 func_t *func; 00062 00063 /// User-specified gradient 00064 dfunc_t *grad; 00065 00066 /// Automatic gradient object 00067 auto_grad_t *agrad; 00068 00069 /// If true, a gradient has been specified 00070 bool grad_given; 00071 00072 /// User-specified parameter 00073 param_t *params; 00074 /// Memory size 00075 size_t dim; 00076 /// Memory allocation 00077 alloc_t ao; 00078 /// Temporary vector 00079 alloc_vec_t avt; 00080 /// Temporary vector 00081 alloc_vec_t avt2; 00082 00083 /// Take a step 00084 void take_step(const gsl_vector * x, const gsl_vector *px, 00085 double stepx, double lambda, gsl_vector *x1x, 00086 gsl_vector *dx) { 00087 00088 gsl_vector_set_zero(dx); 00089 gsl_blas_daxpy(-stepx * lambda, px, dx); 00090 00091 gsl_vector_memcpy (x1x, x); 00092 gsl_blas_daxpy(1.0, dx, x1x); 00093 } 00094 00095 /** 00096 \brief Line minimization 00097 00098 Do a line minimisation in the region (xa,fa) (xc,fc) to find 00099 an intermediate (xb,fb) satisifying fa > fb < fc. Choose an 00100 initial xb based on parabolic interpolation 00101 */ 00102 void intermediate_point(const gsl_vector *x, const gsl_vector *px, 00103 double lambda, double pg, double stepa, 00104 double stepc, double fa, double fc, 00105 gsl_vector *x1x, gsl_vector *dx, 00106 gsl_vector *gradient, double *stepx, double *f) { 00107 double stepb, fb; 00108 00109 bool trial_failed; 00110 00111 do { 00112 00113 double u = fabs (pg * lambda * stepc); 00114 stepb = 0.5 * stepc * u / ((fc - fa) + u); 00115 00116 take_step (x, px, stepb, lambda, x1x, dx); 00117 00118 for(size_t i=0;i<dim;i++) avt[i]=gsl_vector_get(x1x,i); 00119 (*func)(dim,avt,fb,*params); 00120 00121 trial_failed=false; 00122 if (fb >= fa && stepb > 0.0) { 00123 /* downhill step failed, reduce step-size and try again */ 00124 fc = fb; 00125 stepc = stepb; 00126 trial_failed=true; 00127 } 00128 00129 } while (trial_failed); 00130 00131 *stepx = stepb; 00132 *f = fb; 00133 00134 for(size_t i=0;i<dim;i++) avt[i]=gsl_vector_get(x1x,i); 00135 if (grad_given) { 00136 (*grad)(dim,avt,avt2,*params); 00137 } else { 00138 (*agrad)(dim,avt,avt2,*params); 00139 } 00140 for(size_t i=0;i<dim;i++) gsl_vector_set(gradient,i,avt2[i]); 00141 00142 return; 00143 } 00144 00145 /** 00146 \brief Perform the minimization 00147 00148 Starting at (x0, f0) move along the direction p to find a minimum 00149 f(x0 - lambda * p), returning the new point x1 = x0-lambda*p, 00150 f1=f(x1) and g1 = grad(f) at x1. 00151 00152 */ 00153 void minimize(const gsl_vector *x, const gsl_vector *xp, double lambda, 00154 double stepa, double stepb, double stepc, double fa, 00155 double fb, double fc, double xtol, gsl_vector *x1x, 00156 gsl_vector *dx1x, gsl_vector *x2x, gsl_vector *dx2x, 00157 gsl_vector *gradient, double *xstep, double *f, 00158 double *gnorm_u) { 00159 00160 double u = stepb; 00161 double v = stepa; 00162 double w = stepc; 00163 double fu = fb; 00164 double fv = fa; 00165 double fw = fc; 00166 00167 double old2 = fabs(w - v); 00168 double old1 = fabs(v - u); 00169 00170 double stepm, fm, pg, gnorm1; 00171 00172 int xiter = 0; 00173 00174 gsl_vector_memcpy (x2x, x1x); 00175 gsl_vector_memcpy (dx2x, dx1x); 00176 00177 *f = fb; 00178 *xstep = stepb; 00179 *gnorm_u = gsl_blas_dnrm2(gradient); 00180 00181 mid_trial: 00182 00183 xiter++; 00184 00185 if (xiter > nmaxiter) { 00186 // Exceeded maximum number of iterations 00187 return; 00188 } 00189 00190 { 00191 double dw = w - u; 00192 double dv = v - u; 00193 double du = 0.0; 00194 double e1 = ((fv - fu) * dw * dw + (fu - fw) * dv * dv); 00195 double e2 = 2.0 * ((fv - fu) * dw + (fu - fw) * dv); 00196 00197 if (e2 != 0.0) { 00198 du = e1 / e2; 00199 } 00200 00201 if (du > 0.0 && du < (stepc - stepb) && fabs(du) < 0.5 * old2) { 00202 stepm = u + du; 00203 } else if (du < 0.0 && du > (stepa - stepb) && 00204 fabs(du) < 0.5 * old2) { 00205 stepm = u + du; 00206 } else if ((stepc - stepb) > (stepb - stepa)) { 00207 stepm = 0.38 * (stepc - stepb) + stepb; 00208 } else { 00209 stepm = stepb - 0.38 * (stepb - stepa); 00210 } 00211 } 00212 00213 take_step (x, xp, stepm, lambda, x1x, dx1x); 00214 00215 for(size_t i=0;i<dim;i++) avt[i]=gsl_vector_get(x1x,i); 00216 (*func)(dim,avt,fm,*params); 00217 00218 if (fm > fb) { 00219 if (fm < fv) { 00220 w = v; 00221 v = stepm; 00222 fw = fv; 00223 fv = fm; 00224 } else if (fm < fw) { 00225 w = stepm; 00226 fw = fm; 00227 } 00228 if (stepm < stepb) { 00229 stepa = stepm; 00230 fa = fm; 00231 } else { 00232 stepc = stepm; 00233 fc = fm; 00234 } 00235 goto mid_trial; 00236 00237 } else if (fm <= fb) { 00238 00239 old2 = old1; 00240 old1 = fabs(u - stepm); 00241 w = v; 00242 v = u; 00243 u = stepm; 00244 fw = fv; 00245 fv = fu; 00246 fu = fm; 00247 00248 gsl_vector_memcpy (x2x, x1x); 00249 gsl_vector_memcpy (dx2x, dx1x); 00250 00251 for(size_t i=0;i<dim;i++) avt[i]=gsl_vector_get(x1x,i); 00252 if (grad_given) { 00253 (*grad)(dim,avt,avt2,*params); 00254 } else { 00255 (*agrad)(dim,avt,avt2,*params); 00256 } 00257 for(size_t i=0;i<dim;i++) gsl_vector_set(gradient,i,avt2[i]); 00258 00259 gsl_blas_ddot (xp, gradient, &pg); 00260 gnorm1 = gsl_blas_dnrm2 (gradient); 00261 00262 *f = fm; 00263 *xstep = stepm; 00264 *gnorm_u = gnorm1; 00265 00266 if (fabs (pg * lambda / gnorm1) < xtol) { 00267 /* SUCCESS */ 00268 return; 00269 } 00270 00271 if (stepm < stepb) { 00272 stepc = stepb; 00273 fc = fb; 00274 stepb = stepm; 00275 fb = fm; 00276 } else { 00277 stepa = stepb; 00278 fa = fb; 00279 stepb = stepm; 00280 fb = fm; 00281 } 00282 goto mid_trial; 00283 } 00284 } 00285 00286 #endif 00287 00288 public: 00289 00290 /// Stepsize for finite-differencing ( default \f$ 10^{-4} \f$ ) 00291 double deriv_h; 00292 00293 /// Maximum iterations for line minimization (default 10) 00294 int nmaxiter; 00295 00296 /// Default automatic Gradient object 00297 def_auto_grad_t def_grad; 00298 00299 gsl_mmin_base() { 00300 deriv_h=1.0e-4; 00301 nmaxiter=10; 00302 grad_given=false; 00303 agrad=&def_grad; 00304 } 00305 00306 /// Set the function 00307 int base_set(func_t &ufunc, param_t &pa) { 00308 func=&ufunc; 00309 params=&pa; 00310 grad_given=false; 00311 return 0; 00312 } 00313 00314 /// Set the function and the gradient 00315 int base_set_de(func_t &ufunc, dfunc_t &udfunc, param_t &pa) { 00316 func=&ufunc; 00317 params=&pa; 00318 grad=&udfunc; 00319 grad_given=true; 00320 return 0; 00321 } 00322 00323 /// Allocate memory 00324 int base_allocate(size_t nn) { 00325 ao.allocate(avt,nn); 00326 ao.allocate(avt2,nn); 00327 dim=nn; 00328 return 0; 00329 } 00330 00331 /// Clear allocated memory 00332 int base_free() { 00333 ao.free(avt); 00334 ao.free(avt2); 00335 dim=0; 00336 return 0; 00337 } 00338 00339 00340 }; 00341 00342 /** \brief Multidimensional minimization by the Fletcher-Reeves 00343 conjugate gradient algorithm (GSL) 00344 00345 The variable multi_min::tolf is used as the maximum value 00346 of the gradient and is \f$ 10^{-4} \f$ be default. 00347 00348 The gsl iterate() function for this minimizer chooses to return 00349 \c GSL_ENOPROG if the iteration fails to make progress without 00350 calling the error handler. This is presumably because the 00351 iterate() function can fail to make progress when the algorithm 00352 has succeeded in finding the minimum. I prefer to return a 00353 non-zero value from a function only in cases where the error 00354 handler will also be called, so the user is clear on what an 00355 "error" means in the local context. Thus if iterate() is failing 00356 to make progress, instead of returning a non-zero value, it sets 00357 the value of \ref it_info to a non-zero value. 00358 00359 \todo A bit of needless copying is required in the function 00360 wrapper to convert from \c gsl_vector to the templated vector 00361 type. This can be fixed, probably by rewriting take_step 00362 to produce a vec_t &x1 rather than a gsl_vector *x1; 00363 00364 Those who look in detail at the code will note that the state 00365 variable \c max_iter has not been included here, because 00366 it was not really used in the original GSL code for 00367 these minimizers. 00368 00369 Default template arguments 00370 - \c param_t - \ref multi_funct 00371 - \c vec_t - \ref ovector_view 00372 - \c alloc_vec_t - \ref ovector 00373 - \c alloc_t - \ref ovector_alloc 00374 - \c dfunc_t - \ref grad_funct 00375 - \c auto_grad_t - \ref gradient 00376 - \c def_auto_grad_t - \ref simple_grad 00377 */ 00378 template<class param_t, class func_t=multi_funct<param_t>, 00379 class vec_t=ovector_view, class alloc_vec_t=ovector, 00380 class alloc_t=ovector_alloc, 00381 class dfunc_t=grad_funct<param_t,ovector_view>, 00382 class auto_grad_t=gradient<param_t,func_t,ovector_view>, 00383 class def_auto_grad_t=simple_grad<param_t,func_t,ovector_view> > 00384 class gsl_mmin_conf : 00385 public gsl_mmin_base<param_t,func_t,vec_t,alloc_vec_t,alloc_t,dfunc_t, 00386 auto_grad_t,def_auto_grad_t> { 00387 00388 #ifndef DOXYGEN_INTERNAL 00389 00390 protected: 00391 00392 /// \name The original variables from the GSL state structure 00393 //@{ 00394 /// Desc 00395 int iter; 00396 /// Desc 00397 double step; 00398 /// Desc 00399 double tol; 00400 /// Desc 00401 gsl_vector *x1; 00402 /// Desc 00403 gsl_vector *dx1; 00404 /// Desc 00405 gsl_vector *x2; 00406 /// Desc 00407 double pnorm; 00408 /// Desc 00409 gsl_vector *p; 00410 /// Desc 00411 double g0norm; 00412 /// Desc 00413 gsl_vector *g0; 00414 //@} 00415 00416 /// \name Store the arguments to set() so we can use them for iterate() 00417 //@{ 00418 /// Desc 00419 gsl_vector *ugx; 00420 /// Desc 00421 gsl_vector *ugg; 00422 /// Desc 00423 gsl_vector *udx; 00424 /// Desc 00425 double it_min; 00426 //@} 00427 00428 /// Temporary vector 00429 alloc_vec_t avt5; 00430 /// Temporary vector 00431 alloc_vec_t avt6; 00432 /// Temporary vector 00433 alloc_vec_t avt7; 00434 /// Temporary vector 00435 alloc_vec_t avt8; 00436 00437 #endif 00438 00439 public: 00440 00441 gsl_mmin_conf() { 00442 lmin_tol=1.0e-4; 00443 this->tolf=1.0e-4; 00444 step_size=0.01; 00445 it_info=0; 00446 } 00447 00448 virtual ~gsl_mmin_conf() {} 00449 00450 /// Tolerance for the line minimization (default \f$ 10^{-4} \f$) 00451 double lmin_tol; 00452 00453 /// Size of the initial step (default 0.01) 00454 double step_size; 00455 00456 /** 00457 \brief Information from the last call to iterate() 00458 00459 This is 1 if \c pnorm or \c gnorm are 0 and 00460 2 if \c stepb is zero. 00461 00462 \todo Document this better 00463 */ 00464 int it_info; 00465 00466 /// Perform an iteration 00467 virtual int iterate() { 00468 00469 if (this->dim==0) { 00470 set_err_ret("Memory not allocated in iterate().",gsl_efailed); 00471 } 00472 00473 gsl_vector *x=ugx; 00474 gsl_vector *gradient=ugg; 00475 gsl_vector *dx=udx; 00476 00477 double fa = it_min, fb, fc; 00478 double dir; 00479 double stepa = 0.0, stepb, stepc = step; 00480 00481 double g1norm; 00482 double pg; 00483 00484 if (pnorm == 0.0 || g0norm == 0.0) { 00485 gsl_vector_set_zero(dx); 00486 it_info=1; 00487 return 0; 00488 } 00489 00490 /* Determine which direction is downhill, +p or -p */ 00491 00492 gsl_blas_ddot(p, gradient, &pg); 00493 00494 dir = (pg >= 0.0) ? +1.0 : -1.0; 00495 00496 /* Compute new trial point at x_c= x - step * p, where p is the 00497 current direction */ 00498 00499 this->take_step (x, p, stepc, dir / pnorm, x1, dx); 00500 00501 /* Evaluate function and gradient at new point xc */ 00502 00503 for(size_t i=0;i<this->dim;i++) avt5[i]=gsl_vector_get(x1,i); 00504 (*this->func)(this->dim,avt5,fc,*this->params); 00505 00506 if (fc < fa) { 00507 00508 /* Success, reduced the function value */ 00509 step = stepc * 2.0; 00510 it_min = fc; 00511 gsl_vector_memcpy(x,x1); 00512 00513 for(size_t i=0;i<this->dim;i++) avt6[i]=gsl_vector_get(x1,i); 00514 if (this->grad_given) { 00515 (*this->grad)(this->dim,avt6,avt7,*this->params); 00516 } else { 00517 (*this->agrad)(this->dim,avt6,avt7,*this->params); 00518 } 00519 for(size_t i=0;i<this->dim;i++) gsl_vector_set(gradient,i,avt7[i]); 00520 00521 it_info=0; 00522 return gsl_success; 00523 } 00524 00525 /* Do a line minimisation in the region (xa,fa) (xc,fc) to find an 00526 intermediate (xb,fb) satisifying fa > fb < fc. Choose an initial 00527 xb based on parabolic interpolation */ 00528 00529 this->intermediate_point(x,p,dir / pnorm,pg,stepa,stepc,fa,fc,x1, 00530 dx1,gradient,&stepb,&fb); 00531 00532 if (stepb == 0.0) { 00533 it_info=2; 00534 return 0; 00535 } 00536 00537 this->minimize(x,p,dir / pnorm,stepa,stepb,stepc,fa,fb,fc,tol, 00538 x1,dx1,x2,dx,gradient,&step,&it_min,&g1norm); 00539 00540 gsl_vector_memcpy (x, x2); 00541 00542 /* Choose a new conjugate direction for the next step */ 00543 iter=(iter+1) % x->size; 00544 00545 if (iter==0) { 00546 00547 gsl_vector_memcpy (p, gradient); 00548 pnorm=g1norm; 00549 00550 } else { 00551 00552 /* p' = g1 - beta * p */ 00553 double beta = -pow (g1norm / g0norm, 2.0); 00554 gsl_blas_dscal (-beta, p); 00555 gsl_blas_daxpy (1.0, gradient, p); 00556 pnorm=gsl_blas_dnrm2 (p); 00557 } 00558 00559 g0norm=g1norm; 00560 gsl_vector_memcpy (g0, gradient); 00561 00562 it_info=0; 00563 return gsl_success; 00564 00565 } 00566 00567 /// Return string denoting type("gsl_mmin_conf") 00568 virtual const char *type() { return "gsl_mmin_conf";} 00569 00570 /// Allocate the memory 00571 virtual int allocate(size_t n) { 00572 00573 x1=gsl_vector_calloc(n); 00574 if (x1==0) { 00575 set_err_ret("Failed to allocate x1 in gsl_mmin_conf::allocate().", 00576 gsl_enomem); 00577 } 00578 00579 dx1=gsl_vector_calloc(n); 00580 if (dx1==0) { 00581 gsl_vector_free(x1); 00582 set_err_ret("Failed to allocate dx1 in gsl_mmin_conf::allocate().", 00583 gsl_enomem); 00584 } 00585 00586 x2=gsl_vector_calloc(n); 00587 if (x2==0) { 00588 gsl_vector_free(dx1); 00589 gsl_vector_free(x1); 00590 set_err_ret("Failed to allocate x2 in gsl_mmin_conf::allocate().", 00591 gsl_enomem); 00592 } 00593 00594 p=gsl_vector_calloc(n); 00595 if (p==0) { 00596 gsl_vector_free(x2); 00597 gsl_vector_free(dx1); 00598 gsl_vector_free(x1); 00599 set_err_ret("Failed to allocate p in gsl_mmin_conf::allocate().", 00600 gsl_enomem); 00601 } 00602 00603 g0=gsl_vector_calloc(n); 00604 if (g0==0) { 00605 gsl_vector_free(p); 00606 gsl_vector_free(x2); 00607 gsl_vector_free(dx1); 00608 gsl_vector_free(x1); 00609 set_err_ret("Failed to allocate g0 in gsl_mmin_conf::allocate().", 00610 gsl_enomem); 00611 } 00612 00613 udx=gsl_vector_calloc(n); 00614 if (udx==0) { 00615 gsl_vector_free(g0); 00616 gsl_vector_free(p); 00617 gsl_vector_free(x2); 00618 gsl_vector_free(dx1); 00619 gsl_vector_free(x1); 00620 set_err_ret("Failed to allocate udx in gsl_mmin_conf::allocate().", 00621 gsl_enomem); 00622 } 00623 00624 ugg=gsl_vector_calloc(n); 00625 if (ugg==0) { 00626 gsl_vector_free(udx); 00627 gsl_vector_free(g0); 00628 gsl_vector_free(p); 00629 gsl_vector_free(x2); 00630 gsl_vector_free(dx1); 00631 gsl_vector_free(x1); 00632 set_err_ret("Failed to allocate ugg in gsl_mmin_conf::allocate().", 00633 gsl_enomem); 00634 } 00635 00636 ugx=gsl_vector_calloc(n); 00637 if (ugx==0) { 00638 gsl_vector_free(ugg); 00639 gsl_vector_free(udx); 00640 gsl_vector_free(g0); 00641 gsl_vector_free(p); 00642 gsl_vector_free(x2); 00643 gsl_vector_free(dx1); 00644 gsl_vector_free(x1); 00645 set_err_ret("Failed to allocate ugx in gsl_mmin_conf::allocate().", 00646 gsl_enomem); 00647 } 00648 00649 this->ao.allocate(avt5,n); 00650 this->ao.allocate(avt6,n); 00651 this->ao.allocate(avt7,n); 00652 this->ao.allocate(avt8,n); 00653 00654 this->base_allocate(n); 00655 00656 return gsl_success; 00657 00658 } 00659 00660 /// Free the allocated memory 00661 virtual int free() { 00662 gsl_vector_free(ugx); 00663 gsl_vector_free(ugg); 00664 gsl_vector_free(udx); 00665 gsl_vector_free(g0); 00666 gsl_vector_free(p); 00667 gsl_vector_free(x2); 00668 gsl_vector_free(dx1); 00669 gsl_vector_free(x1); 00670 00671 this->ao.free(avt5); 00672 this->ao.free(avt6); 00673 this->ao.free(avt7); 00674 this->ao.free(avt8); 00675 00676 this->base_free(); 00677 00678 return 0; 00679 } 00680 00681 /// Reset the minimizer to use the current point as a new starting point 00682 int restart() { 00683 iter=0; 00684 return gsl_success; 00685 } 00686 00687 /// Set the function and initial guess 00688 virtual int set(vec_t &x, double u_step_size, double tol_u, 00689 func_t &ufunc, param_t &pa) { 00690 00691 if (this->dim==0) { 00692 set_err_ret("Memory not allocated in set().",gsl_efailed); 00693 } 00694 00695 this->base_set(ufunc,pa); 00696 for(size_t i=0;i<this->dim;i++) gsl_vector_set(ugx,i,x[i]); 00697 00698 iter=0; 00699 step=u_step_size; 00700 tol=tol_u; 00701 00702 /// Evaluate the function and its gradient 00703 ufunc(this->dim,x,it_min,pa); 00704 this->agrad->set_function(ufunc); 00705 (*this->agrad)(this->dim,x,avt8,*this->params); 00706 for(size_t i=0;i<this->dim;i++) { 00707 gsl_vector_set(ugg,i,avt8[i]); 00708 } 00709 00710 /* Use the gradient as the initial direction */ 00711 00712 gsl_vector_memcpy(p,ugg); 00713 gsl_vector_memcpy(g0,ugg); 00714 00715 double gnorm=gsl_blas_dnrm2(ugg); 00716 00717 pnorm=gnorm; 00718 g0norm=gnorm; 00719 00720 return gsl_success; 00721 } 00722 00723 /// Set the function and initial guess 00724 virtual int set_de(vec_t &x, double u_step_size, double tol_u, 00725 func_t &ufunc, dfunc_t &udfunc, param_t &pa) { 00726 00727 if (this->dim==0) { 00728 set_err_ret("Memory not allocated in set().",gsl_efailed); 00729 } 00730 00731 this->base_set_de(ufunc,udfunc,pa); 00732 for(size_t i=0;i<this->dim;i++) gsl_vector_set(ugx,i,x[i]); 00733 00734 iter=0; 00735 step=u_step_size; 00736 tol =tol_u; 00737 00738 /// Evaluate the function and its gradient 00739 ufunc(this->dim,x,it_min,pa); 00740 (*this->agrad)(this->dim,x,avt8,*this->params); 00741 for(size_t i=0;i<this->dim;i++) { 00742 gsl_vector_set(ugg,i,avt8[i]); 00743 } 00744 00745 /* Use the gradient as the initial direction */ 00746 00747 gsl_vector_memcpy(p,ugg); 00748 gsl_vector_memcpy(g0,ugg); 00749 00750 double gnorm=gsl_blas_dnrm2(ugg); 00751 00752 pnorm=gnorm; 00753 g0norm=gnorm; 00754 00755 return gsl_success; 00756 } 00757 00758 /** \brief Calculate the minimum \c min of \c func w.r.t the 00759 array \c x of size \c nvar. 00760 */ 00761 virtual int mmin(size_t nn, vec_t &xx, double &fmin, param_t &pa, 00762 func_t &ufunc) { 00763 00764 int xiter=0, status; 00765 00766 allocate(nn); 00767 00768 set(xx,step_size,lmin_tol,ufunc,pa); 00769 00770 do { 00771 xiter++; 00772 00773 status=iterate(); 00774 00775 if (status || it_info!=0) { 00776 break; 00777 } 00778 00779 // Equivalent to gsl_multimin_test_gradient with 00780 // additional code to print out present iteration 00781 00782 double norm=gsl_blas_dnrm2(ugg); 00783 00784 if(this->verbose>0) { 00785 this->print_iter(nn,*((ovector *)ugx),it_min,xiter, 00786 norm,this->tolf,type()); 00787 } 00788 00789 if (norm<this->tolf) status=gsl_success; 00790 else status=gsl_continue; 00791 00792 } 00793 while (status==GSL_CONTINUE && xiter < this->ntrial); 00794 00795 for(size_t i=0;i<nn;i++) xx[i]=gsl_vector_get(ugx,i); 00796 fmin=it_min; 00797 00798 free(); 00799 this->last_ntrial=xiter; 00800 00801 if (xiter==this->ntrial) { 00802 set_err_ret("Too many iterations in gsl_mmin_conf::mmin().", 00803 gsl_emaxiter); 00804 } 00805 00806 return status; 00807 } 00808 00809 /** \brief Calculate the minimum \c min of \c func w.r.t the 00810 array \c x of size \c nvar. 00811 */ 00812 virtual int mmin_de(size_t nn, vec_t &xx, double &fmin, param_t &pa, 00813 func_t &ufunc, dfunc_t &udfunc) { 00814 00815 int xiter=0, status; 00816 00817 allocate(nn); 00818 00819 set_de(xx,step_size,lmin_tol,ufunc,udfunc,pa); 00820 00821 do { 00822 xiter++; 00823 00824 status=iterate(); 00825 00826 if (status || it_info!=0) { 00827 break; 00828 } 00829 00830 // Equivalent to gsl_multimin_test_gradient with 00831 // additional code to print out present iteration 00832 00833 double norm=gsl_blas_dnrm2(ugg); 00834 00835 if(this->verbose>0) { 00836 this->print_iter(nn,*((ovector *)ugx),it_min,xiter, 00837 norm,this->tolf,type()); 00838 } 00839 00840 if (norm<this->tolf) status=gsl_success; 00841 else status=gsl_continue; 00842 00843 } 00844 while (status==GSL_CONTINUE && xiter < this->ntrial); 00845 00846 for(size_t i=0;i<nn;i++) xx[i]=gsl_vector_get(ugx,i); 00847 fmin=it_min; 00848 00849 free(); 00850 this->last_ntrial=xiter; 00851 00852 if (xiter==this->ntrial) { 00853 set_err_ret("Too many iterations in gsl_mmin_conf::mmin_de().", 00854 gsl_emaxiter); 00855 } 00856 00857 return status; 00858 } 00859 00860 }; 00861 00862 /** 00863 \brief An array version of gsl_mmin_conf 00864 00865 \comment 00866 Doxygen has had trouble parsing this template, so I have 00867 simplified it slightly for the documentation below. 00868 \endcomment 00869 */ 00870 #ifdef DOXYGENP 00871 template<class param_t, size_t nv> class gsl_mmin_conf_array : 00872 public gsl_mmin_conf<param_t,func_t,vec_t,alloc_vec_t, 00873 alloc_t,dfunc_t,auto_grad_t,def_auto_grad_t> { 00874 }; 00875 #else 00876 template<class param_t, size_t nv> class gsl_mmin_conf_array : 00877 public gsl_mmin_conf<param_t,multi_vfunct<param_t,nv>, 00878 double[nv],double[nv],array_alloc<double[nv]>,grad_vfunct<param_t,nv>, 00879 gradient_array<param_t,multi_vfunct<param_t,nv>,nv>, 00880 simple_grad_array<param_t,multi_vfunct<param_t,nv>,nv> > { 00881 }; 00882 #endif 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