00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 Some of the source code in this file has been adapted from OOL 00009 library written by Sergio Drumond Ventura, Luis Alberto D'Afonseca, 00010 and Ricardo Biloti available from http://ool.sourceforge.net. 00011 00012 O2scl is free software; you can redistribute it and/or modify 00013 it under the terms of the GNU General Public License as published by 00014 the Free Software Foundation; either version 3 of the License, or 00015 (at your option) any later version. 00016 00017 O2scl is distributed in the hope that it will be useful, 00018 but WITHOUT ANY WARRANTY; without even the implied warranty of 00019 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00020 GNU General Public License for more details. 00021 00022 You should have received a copy of the GNU General Public License 00023 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00024 00025 ------------------------------------------------------------------- 00026 */ 00027 #ifndef O2SCL_OOL_MMIN_SPG_H 00028 #define O2SCL_OOL_MMIN_SPG_H 00029 00030 #include <o2scl/text_file.h> 00031 #include <o2scl/ovector_tlate.h> 00032 #include <o2scl/multi_funct.h> 00033 #include <o2scl/ool_constr_mmin.h> 00034 #include <gsl/gsl_math.h> 00035 00036 #ifndef DOXYGENP 00037 namespace o2scl { 00038 #endif 00039 00040 /** 00041 \brief Constrained minimization by the 00042 spectral projected gradient method (OOL) 00043 */ 00044 template<class param_t, class func_t, class dfunc_t, 00045 class vec_t=ovector_view, class alloc_vec_t=ovector, 00046 class alloc_t=ovector_alloc> class ool_mmin_spg : 00047 public ool_constr_mmin<param_t,func_t,dfunc_t,ool_hfunct<void *>,vec_t, 00048 alloc_vec_t,alloc_t> { 00049 00050 #ifndef DOXYGEN_INTERNAL 00051 00052 protected: 00053 00054 /// A convenient typedef for the unused Hessian product type 00055 typedef ool_hfunct<void *> hfunc_t; 00056 00057 /// Armijo parameter 00058 double alpha; 00059 00060 /// Temporary vector 00061 alloc_vec_t xx; 00062 /// Temporary vector 00063 alloc_vec_t d; 00064 /// Temporary vector 00065 alloc_vec_t s; 00066 /// Temporary vector 00067 alloc_vec_t y; 00068 /// Temporary vector 00069 alloc_vec_t fvec; 00070 /// Non-monotone parameter 00071 size_t m; 00072 /// Desc 00073 int tail; 00074 00075 /// Line search 00076 int line_search() { 00077 00078 double lambda, lambdanew; 00079 double fmax, faux, fxx, dTg; 00080 short armijo_failure; 00081 size_t ii; 00082 00083 // Saving the previous gradient and 00084 // compute d = P(x - alpha * g) - x */ 00085 for(size_t i=0;i<this->dim;i++) { 00086 y[i]=this->gradient[i]; 00087 d[i]=this->x[i]; 00088 } 00089 for(size_t i=0;i<this->dim;i++) { 00090 d[i]-=alpha*this->gradient[i]; 00091 } 00092 proj(d); 00093 for(size_t i=0;i<this->dim;i++) { 00094 d[i]-=this->x[i]; 00095 } 00096 00097 dTg=0.0; 00098 for(size_t i=0;i<this->dim;i++) dTg+=d[i]*this->gradient[i]; 00099 00100 lambda=1; 00101 armijo_failure=1; 00102 while (armijo_failure) { 00103 00104 /* x trial */ 00105 for(size_t i=0;i<this->dim;i++) xx[i]=this->x[i]; 00106 for(size_t i=0;i<this->dim;i++) xx[i]+=lambda*d[i]; 00107 (*this->func)(this->dim,xx,fxx,*this->param); 00108 00109 fmax=0; 00110 for(ii=0;ii<m;ii++) { 00111 faux=fvec[ii]+gamma*lambda*dTg; 00112 fmax=GSL_MAX(fmax,faux); 00113 } 00114 00115 if (fxx>fmax) { 00116 armijo_failure=1; 00117 lambdanew=-lambda*lambda*dTg/(2*(fxx-this->f-lambda*dTg)); 00118 lambda=GSL_MAX(sigma1*lambda,GSL_MIN(sigma2*lambda,lambdanew)); 00119 } else { 00120 armijo_failure=0; 00121 } 00122 } 00123 00124 /* st->s = x_{k+1} - x_k */ 00125 for(size_t i=0;i<this->dim;i++) s[i]=xx[i]; 00126 for(size_t i=0;i<this->dim;i++) s[i]-=this->x[i]; 00127 00128 /* M->x = x_{k+1} */ 00129 for(size_t i=0;i<this->dim;i++) this->x[i]=xx[i]; 00130 (*this->dfunc)(this->dim,this->x,this->gradient,*this->param); 00131 this->f=fxx; 00132 00133 /* Infinite norm of g1 = d <- [P(x - g) - x] */ 00134 for(size_t i=0;i<this->dim;i++) d[i]=this->x[i]; 00135 for(size_t i=0;i<this->dim;i++) d[i]-=this->gradient[i]; 00136 proj(d); 00137 for(size_t i=0;i<this->dim;i++) d[i]-=this->x[i]; 00138 00139 double maxval=fabs(d[0]); 00140 for(size_t i=1;i<this->dim;i++) { 00141 if (fabs(d[i])>maxval) { 00142 maxval=fabs(d[i]); 00143 } 00144 } 00145 this->size=fabs(maxval); 00146 00147 /* st->y = - (g(x_{k+1}) - g(x_k)) */ 00148 for(size_t i=0;i<this->dim;i++) y[i]-=this->gradient[i]; 00149 00150 m++; 00151 m=GSL_MIN(M,m); 00152 tail++; 00153 tail=tail % M; 00154 fvec[tail]=this->f; 00155 00156 return 0; 00157 } 00158 00159 /// Project into feasible region 00160 int proj(vec_t &xt) { 00161 size_t ii, n=this->dim; 00162 00163 for(ii=0;ii<n;ii++) { 00164 xt[ii]=GSL_MAX(this->L[ii],GSL_MIN(xt[ii],this->U[ii])); 00165 } 00166 return 0; 00167 } 00168 00169 #endif 00170 00171 public: 00172 00173 ool_mmin_spg() { 00174 fmin=-1.0e99; 00175 tol=1.0e-4; 00176 alphamin=1.0e-30; 00177 alphamax=1.0e30; 00178 gamma=1.0e-4; 00179 sigma1=0.1; 00180 sigma2=0.9; 00181 M=10; 00182 } 00183 00184 ~ool_mmin_spg() { 00185 this->ao.free(fvec); 00186 } 00187 00188 /** 00189 \brief Minimum function value (default \f$ 10^{-99} \f$) 00190 00191 If the function value is below this value, then the algorithm 00192 assumes that the function is not bounded and exits. 00193 */ 00194 double fmin; 00195 /// Tolerance on infinite norm (default \f$ 10^{-4} \f$) 00196 double tol; 00197 /// Lower bound to spectral step size (default \f$ 10^{-30} \f$) 00198 double alphamin; 00199 /// Upper bound to spectral step size (default \f$ 10^{30} \f$) 00200 double alphamax; 00201 /// Sufficient decrease parameter (default \f$ 10^{-4} \f$) 00202 double gamma; 00203 /// Lower bound to the step length reduction (default 0.1) 00204 double sigma1; 00205 /// Upper bound to the step length reduction (default 0.9) 00206 double sigma2; 00207 /// Monotonicity parameter (M=1 forces monotonicity) (default 10) 00208 size_t M; 00209 00210 /// Allocate memory 00211 virtual int allocate(const size_t n) { 00212 if (this->dim>0) free(); 00213 this->ao.allocate(xx,n); 00214 this->ao.allocate(d,n); 00215 this->ao.allocate(s,n); 00216 this->ao.allocate(y,n); 00217 return ool_constr_mmin<param_t,func_t,dfunc_t,hfunc_t,vec_t,alloc_vec_t, 00218 alloc_t>::allocate(n); 00219 } 00220 00221 /// Free previously allocated memory 00222 virtual int free() { 00223 if (this->dim>0) this->ao.free(xx); 00224 if (this->dim>0) this->ao.free(d); 00225 if (this->dim>0) this->ao.free(s); 00226 if (this->dim>0) this->ao.free(y); 00227 return ool_constr_mmin<param_t,func_t,dfunc_t,hfunc_t,vec_t,alloc_vec_t, 00228 alloc_t>::free(); 00229 } 00230 00231 /// Set the function, the initial guess, and the parameters 00232 virtual int set(func_t &fn, dfunc_t &dfn, vec_t &init, param_t &par) { 00233 00234 int ret=ool_constr_mmin<param_t,func_t,dfunc_t,hfunc_t,vec_t,alloc_vec_t, 00235 alloc_t>::set(fn,dfn,init,par); 00236 00237 // Turn initial guess into a feasible point 00238 proj(this->x); 00239 00240 // Evaluate function and gradient 00241 (*this->func)(this->dim,this->x,this->f,*this->param); 00242 (*this->dfunc)(this->dim,this->x,this->gradient,*this->param); 00243 00244 /* Infinite norm of d <- g1 = [P(x - g) - x] */ 00245 o2scl::vector_copy(this->dim,this->x,d); 00246 for(size_t i=0;i<this->dim;i++) { 00247 d[i]-=this->gradient[i]; 00248 } 00249 proj(d); 00250 for(size_t i=0;i<this->dim;i++) d[i]-=this->x[i]; 00251 00252 double maxval=fabs(d[0]); 00253 for(size_t i=1;i<this->dim;i++) { 00254 if (fabs(d[i])>maxval) { 00255 maxval=fabs(d[i]); 00256 } 00257 } 00258 this->size=fabs(maxval); 00259 00260 /* alpha_0 */ 00261 maxval=fabs(this->gradient[0]); 00262 for(size_t i=1;i<this->dim;i++) { 00263 if (fabs(this->gradient[i])>maxval) { 00264 maxval=fabs(this->gradient[i]); 00265 } 00266 } 00267 alpha=1.0/fabs(maxval); 00268 00269 // Fixme: at the moment, this memory allocation is 00270 // left hanging. It's taken care of the destructor, 00271 // but this isn't so elegant. 00272 this->ao.allocate(fvec,M); 00273 m=1; 00274 tail=0; 00275 fvec[tail]=this->f; 00276 00277 return 0; 00278 } 00279 00280 /// Restart the minimizer 00281 virtual int restart() { 00282 00283 proj(this->x); 00284 00285 // Evaluate function and gradient 00286 (*this->func)(this->dim,this->x,this->f,*this->param); 00287 (*this->dfunc)(this->dim,this->x,this->gradient,*this->param); 00288 00289 /* alpha_0 */ 00290 double maxval=fabs(this->gradient[0]); 00291 for(size_t i=1;i<this->dim;i++) { 00292 if (fabs(this->gradient[i])>maxval) { 00293 maxval=fabs(this->gradient[i]); 00294 } 00295 } 00296 alpha=1.0/fabs(maxval); 00297 00298 /* Infinite norm of d <- g1 = [P(x - g) - x] */ 00299 o2scl::vector_copy(this->dim,this->x,d); 00300 for(size_t i=0;i<this->dim;i++) { 00301 d[i]-=this->gradient[i]; 00302 } 00303 proj(d); 00304 for(size_t i=0;i<this->dim;i++) d[i]-=this->x[i]; 00305 00306 maxval=fabs(d[0]); 00307 for(size_t i=1;i<this->dim;i++) { 00308 if (fabs(d[i])>maxval) { 00309 maxval=fabs(d[i]); 00310 } 00311 } 00312 this->size=fabs(maxval); 00313 00314 m=1; 00315 tail=0; 00316 fvec[tail]=this->f; 00317 00318 return 0; 00319 } 00320 00321 /// Perform an iteration 00322 virtual int iterate() { 00323 double t; 00324 00325 line_search(); 00326 00327 double bk=0.0; 00328 for(size_t i=0;i<this->dim;i++) bk+=s[i]*y[i]; 00329 00330 if (bk>=0.0) { 00331 alpha=alphamax; 00332 } else { 00333 double ak=0.0; 00334 for(size_t i=0;i<this->dim;i++) ak+=s[i]*s[i]; 00335 ak=-ak/bk; 00336 alpha=GSL_MIN(alphamax,GSL_MAX(alphamin,ak)); 00337 } 00338 00339 return 0; 00340 } 00341 00342 /// See if we're finished 00343 virtual int is_optimal() { 00344 if (this->size>tol && this->f>fmin) { 00345 return GSL_CONTINUE; 00346 } 00347 return GSL_SUCCESS; 00348 } 00349 00350 /// Return string denoting type ("ool_mmin_spg") 00351 const char *type() { return "ool_mmin_spg"; } 00352 00353 }; 00354 00355 #ifndef DOXYGENP 00356 } 00357 #endif 00358 00359 #endif 00360
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