![]() |
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_CONP_H 00024 #define O2SCL_GSL_MMIN_CONP_H 00025 00026 #include <gsl/gsl_blas.h> 00027 #include <gsl/gsl_multimin.h> 00028 #include <o2scl/gsl_mmin_conf.h> 00029 00030 #ifndef DOXYGENP 00031 namespace o2scl { 00032 #endif 00033 00034 /** \brief Multidimensional minimization by the Polak-Ribiere 00035 conjugate gradient algorithm (GSL) 00036 00037 The functions mmin() and mmin_de() minimize a given function 00038 until the gradient is smaller than the value of multi_min::tol_rel 00039 (which defaults to \f$ 10^{-4} \f$ ). 00040 00041 See an example for the usage of this class in 00042 \ref ex_mmin_sect . 00043 */ 00044 template<class func_t=multi_funct<>, 00045 class vec_t=ovector_base, class alloc_vec_t=ovector, 00046 class alloc_t=ovector_alloc, class dfunc_t=grad_funct<ovector_base>, 00047 class auto_grad_t=gradient<func_t,ovector_base>, 00048 class def_auto_grad_t=simple_grad<func_t,ovector_base> > 00049 class gsl_mmin_conp : 00050 public gsl_mmin_conf<func_t,vec_t,alloc_vec_t,alloc_t,dfunc_t, 00051 auto_grad_t,def_auto_grad_t> { 00052 00053 public: 00054 00055 /// Perform an iteration 00056 virtual int iterate() { 00057 00058 if (this->dim==0) { 00059 O2SCL_ERR_RET("Memory not allocated in iterate().",gsl_efailed); 00060 } 00061 00062 vec_t &x=this->ugx; 00063 vec_t &gradient=this->ugg; 00064 vec_t &dx=this->udx; 00065 00066 double fa = this->it_min, fb, fc; 00067 double dir; 00068 double stepa = 0.0, stepb, stepc=this->step; 00069 00070 double g1norm; 00071 double pg; 00072 00073 if (this->pnorm == 0.0 || this->g0norm == 0.0) { 00074 for(size_t i=0;i<this->dim;i++) dx[i]=0.0; 00075 this->last_conv=gsl_enoprog; 00076 O2SCL_CONV2_RET("Either pnorm or g0norm vanished in ", 00077 "in gsl_mmin_conp::iterate().", 00078 gsl_enoprog,this->err_nonconv); 00079 return 0; 00080 } 00081 00082 /* Determine which direction is downhill, +p or -p */ 00083 00084 pg=o2scl_cblas::ddot(this->dim,this->p,gradient); 00085 00086 dir = (pg >= 0.0) ? +1.0 : -1.0; 00087 00088 /* Compute new trial point at x_c= x - step * p, where p is the 00089 current direction */ 00090 00091 take_step (x, this->p, stepc, dir / this->pnorm, this->x1, dx); 00092 00093 /* Evaluate function and gradient at new point xc */ 00094 00095 fc=(*this->func)(this->dim,this->x1); 00096 00097 if (fc < fa) { 00098 00099 /* Success, reduced the function value */ 00100 this->step = stepc * 2.0; 00101 this->it_min = fc; 00102 for(size_t i=0;i<this->dim;i++) { 00103 x[i]=this->x1[i]; 00104 } 00105 00106 if (this->grad_given) { 00107 (*this->grad)(this->dim,this->x1,gradient); 00108 } else { 00109 (*this->agrad)(this->dim,this->x1,gradient); 00110 } 00111 00112 return gsl_success; 00113 } 00114 00115 /* Do a line minimisation in the region (xa,fa) (xc,fc) to find an 00116 intermediate (xb,fb) satisifying fa > fb < fc. Choose an initial 00117 xb based on parabolic interpolation */ 00118 00119 this->intermediate_point(x, this->p, dir / this->pnorm, pg, 00120 stepa, stepc, fa, fc, this->x1, this->dx1, 00121 gradient, &stepb, &fb); 00122 00123 if (stepb == 0.0) { 00124 this->last_conv=gsl_enoprog; 00125 O2SCL_CONV2_RET("Variable stepb vanished in gsl_mmin_conp::", 00126 "iterate().",gsl_enoprog,this->err_nonconv); 00127 } 00128 00129 this->minimize(x,this->p,dir / this->pnorm,stepa,stepb,stepc, fa, 00130 fb, fc, this->tol, this->x1, this->dx1, this->x2, 00131 dx, gradient, &(this->step), &(this->it_min), &g1norm); 00132 00133 for(size_t i=0;i<this->dim;i++) x[i]=this->x2[i]; 00134 00135 /* Choose a new conjugate direction for the next step */ 00136 00137 this->iter = (this->iter + 1) % this->dim; 00138 00139 if (this->iter == 0) { 00140 for(size_t i=0;i<this->dim;i++) this->p[i]=gradient[i]; 00141 this->pnorm = g1norm; 00142 } else { 00143 /* p' = g1 - beta * p */ 00144 double g0g1, beta; 00145 00146 /* g0' = g0 - g1 */ 00147 o2scl_cblas::daxpy(-1.0,this->dim,gradient,this->g0); 00148 /* g1g0 = (g0-g1).g1 */ 00149 g0g1=o2scl_cblas::ddot(this->dim,this->g0,gradient); 00150 /* beta = -((g1 - g0).g1)/(g0.g0) */ 00151 beta = g0g1 / (this->g0norm*this->g0norm); 00152 00153 o2scl_cblas::dscal(-beta,this->dim,this->p); 00154 o2scl_cblas::daxpy(1.0,this->dim,gradient,this->p); 00155 this->pnorm = o2scl_cblas::dnrm2(this->dim,this->p); 00156 } 00157 00158 this->g0norm = g1norm; 00159 for(size_t i=0;i<this->dim;i++) { 00160 this->g0[i]=gradient[i]; 00161 } 00162 00163 return gsl_success; 00164 } 00165 00166 /// Return string denoting type("gsl_mmin_conp") 00167 virtual const char *type() { return "gsl_mmin_conp";} 00168 00169 }; 00170 00171 #ifndef DOXYGENP 00172 } 00173 #endif 00174 00175 #endif
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).