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 \note Possibly not working quite yet 00045 */ 00046 template<class param_t, class func_t, class dfunc_t, 00047 class vec_t=ovector_view, class alloc_vec_t=ovector, 00048 class alloc_t=ovector_alloc> class ool_mmin_spg : 00049 public ool_constr_mmin<param_t,func_t,dfunc_t,ool_hfunct<void *>,vec_t, 00050 alloc_vec_t,alloc_t> { 00051 00052 #ifndef DOXYGEN_INTERNAL 00053 00054 protected: 00055 00056 /// A convenient typedef for the unused Hessian product type 00057 typedef ool_hfunct<void *> hfunc_t; 00058 00059 /// Armijo parameter 00060 double alpha; 00061 00062 /// Temporary vector 00063 alloc_vec_t xx; 00064 /// Temporary vector 00065 alloc_vec_t d; 00066 /// Temporary vector 00067 alloc_vec_t s; 00068 /// Temporary vector 00069 alloc_vec_t y; 00070 /// Temporary vector 00071 alloc_vec_t fvec; 00072 /// Desc 00073 size_t m; 00074 /// Desc 00075 int tail; 00076 00077 /// Line search 00078 int line_search() { 00079 00080 double lambda, lambdanew; 00081 double fmax, faux, fxx, dTg; 00082 short armijo_failure; 00083 size_t ii; 00084 00085 // Saving the previous gradient and 00086 // compute d = P(x - alpha * g) - x */ 00087 for(size_t i=0;i<this->dim;i++) { 00088 y[i]=this->gradient[i]; 00089 d[i]=this->x[i]; 00090 } 00091 for(size_t i=0;i<this->dim;i++) { 00092 d[i]-=alpha*this->gradient[i]; 00093 } 00094 proj(d); 00095 for(size_t i=0;i<this->dim;i++) { 00096 d[i]-=this->x[i]; 00097 } 00098 00099 dTg=0.0; 00100 for(size_t i=0;i<this->dim;i++) dTg+=d[i]*this->gradient[i]; 00101 00102 lambda=1; 00103 armijo_failure=1; 00104 while (armijo_failure) { 00105 00106 /* x trial */ 00107 for(size_t i=0;i<this->dim;i++) xx[i]=this->x[i]; 00108 for(size_t i=0;i<this->dim;i++) xx[i]+=lambda*d[i]; 00109 (*this->func)(this->dim,xx,fxx,*this->param); 00110 00111 fmax=0; 00112 for(ii=0;ii<m;ii++) { 00113 faux=fvec[ii]+gamma*lambda*dTg; 00114 fmax=GSL_MAX(fmax,faux); 00115 } 00116 00117 if (fxx>fmax) { 00118 armijo_failure=1; 00119 lambdanew=-lambda*lambda*dTg/(2*(fxx-this->f-lambda*dTg)); 00120 lambda=GSL_MAX(sigma1*lambda,GSL_MIN(sigma2*lambda,lambdanew)); 00121 } else { 00122 armijo_failure=0; 00123 } 00124 } 00125 00126 /* st->s = x_{k+1} - x_k */ 00127 for(size_t i=0;i<this->dim;i++) s[i]=xx[i]; 00128 for(size_t i=0;i<this->dim;i++) s[i]-=this->x[i]; 00129 00130 /* M->x = x_{k+1} */ 00131 for(size_t i=0;i<this->dim;i++) this->x[i]=xx[i]; 00132 (*this->dfunc)(this->dim,this->x,this->gradient,*this->param); 00133 this->f=fxx; 00134 00135 /* Infinite norm of g1 = d <- [P(x - g) - x] */ 00136 for(size_t i=0;i<this->dim;i++) d[i]=this->x[i]; 00137 for(size_t i=0;i<this->dim;i++) d[i]-=this->gradient[i]; 00138 proj(d); 00139 for(size_t i=0;i<this->dim;i++) d[i]-=this->x[i]; 00140 00141 double maxval=fabs(d[0]); 00142 for(size_t i=1;i<this->dim;i++) { 00143 if (fabs(d[i])>maxval) { 00144 maxval=fabs(d[i]); 00145 } 00146 } 00147 this->size=fabs(maxval); 00148 00149 /* st->y = - (g(x_{k+1}) - g(x_k)) */ 00150 for(size_t i=0;i<this->dim;i++) y[i]-=this->gradient[i]; 00151 00152 m++; 00153 m=GSL_MIN(M,m); 00154 tail++; 00155 tail=tail % M; 00156 fvec[tail]=this->f; 00157 00158 return 0; 00159 } 00160 00161 /// Project into feasible region 00162 int proj(vec_t &xt) { 00163 size_t ii, n=this->dim; 00164 00165 for(ii=0;ii<n;ii++) { 00166 xt[ii]=GSL_MAX(this->L[ii],GSL_MIN(xt[ii],this->U[ii])); 00167 } 00168 return 0; 00169 } 00170 00171 #endif 00172 00173 public: 00174 00175 ool_mmin_spg() { 00176 fmin=-1.0e99; 00177 tol=1.0e-4; 00178 alphamin=1.0e-30; 00179 alphamax=1.0e30; 00180 gamma=1.0e-4; 00181 sigma1=0.1; 00182 sigma2=0.9; 00183 M=10; 00184 } 00185 00186 ~ool_mmin_spg() { 00187 this->ao.free(fvec); 00188 } 00189 00190 /** 00191 \brief Minimum function value (default \f$ 10^{-99} \f$) 00192 00193 If the function value is below this value, then the algorithm 00194 assumes that the function is not bounded and exits. 00195 */ 00196 double fmin; 00197 /// Tolerance on infinite norm (default \f$ 10^{-4} \f$) 00198 double tol; 00199 /// Lower bound to spectral step size (default \f$ 10^{-30} \f$) 00200 double alphamin; 00201 /// Upper bound to spectral step size (default \f$ 10^{30} \f$) 00202 double alphamax; 00203 /// Sufficient decrease parameter (default \f$ 10^{-4} \f$) 00204 double gamma; 00205 /// Lower bound to the step length reduction (default 0.1) 00206 double sigma1; 00207 /// Upper bound to the step length reduction (default 0.9) 00208 double sigma2; 00209 /// Monotonicity parameter (M=1 forces monotonicity) (default 10) 00210 size_t M; 00211 00212 /// Allocate memory 00213 virtual int allocate(const size_t n) { 00214 if (this->dim>0) free(); 00215 this->ao.allocate(xx,n); 00216 this->ao.allocate(d,n); 00217 this->ao.allocate(s,n); 00218 this->ao.allocate(y,n); 00219 return ool_constr_mmin<param_t,func_t,dfunc_t,hfunc_t,vec_t,alloc_vec_t, 00220 alloc_t>::allocate(n); 00221 } 00222 00223 /// Free previously allocated memory 00224 virtual int free() { 00225 if (this->dim>0) this->ao.free(xx); 00226 if (this->dim>0) this->ao.free(d); 00227 if (this->dim>0) this->ao.free(s); 00228 if (this->dim>0) this->ao.free(y); 00229 return ool_constr_mmin<param_t,func_t,dfunc_t,hfunc_t,vec_t,alloc_vec_t, 00230 alloc_t>::free(); 00231 } 00232 00233 /// Set the function, the initial guess, and the parameters 00234 virtual int set(func_t &fn, dfunc_t &dfn, vec_t &init, param_t &par) { 00235 00236 int ret=ool_constr_mmin<param_t,func_t,dfunc_t,hfunc_t,vec_t,alloc_vec_t, 00237 alloc_t>::set(fn,dfn,init,par); 00238 00239 // Turn initial guess into a feasible point 00240 proj(this->x); 00241 00242 // Evaluate function and gradient 00243 (*this->func)(this->dim,this->x,this->f,*this->param); 00244 (*this->dfunc)(this->dim,this->x,this->gradient,*this->param); 00245 00246 /* Infinite norm of d <- g1 = [P(x - g) - x] */ 00247 o2scl_arith::vector_copy(this->dim,this->x,d); 00248 for(size_t i=0;i<this->dim;i++) { 00249 d[i]-=this->gradient[i]; 00250 } 00251 proj(d); 00252 for(size_t i=0;i<this->dim;i++) d[i]-=this->x[i]; 00253 00254 double maxval=fabs(d[0]); 00255 for(size_t i=1;i<this->dim;i++) { 00256 if (fabs(d[i])>maxval) { 00257 maxval=fabs(d[i]); 00258 } 00259 } 00260 this->size=fabs(maxval); 00261 00262 /* alpha_0 */ 00263 maxval=fabs(this->gradient[0]); 00264 for(size_t i=1;i<this->dim;i++) { 00265 if (fabs(this->gradient[i])>maxval) { 00266 maxval=fabs(this->gradient[i]); 00267 } 00268 } 00269 alpha=1.0/fabs(maxval); 00270 00271 // Fixme: at the moment, this memory allocation is 00272 // left hanging. It's taken care of the destructor, 00273 // but this isn't so elegant. 00274 this->ao.allocate(fvec,M); 00275 m=1; 00276 tail=0; 00277 fvec[tail]=this->f; 00278 00279 return 0; 00280 } 00281 00282 /// Restart the minimizer 00283 virtual int restart() { 00284 00285 proj(this->x); 00286 00287 // Evaluate function and gradient 00288 (*this->func)(this->dim,this->x,this->f,*this->param); 00289 (*this->dfunc)(this->dim,this->x,this->gradient,*this->param); 00290 00291 /* alpha_0 */ 00292 double maxval=fabs(this->gradient[0]); 00293 for(size_t i=1;i<this->dim;i++) { 00294 if (fabs(this->gradient[i])>maxval) { 00295 maxval=fabs(this->gradient[i]); 00296 } 00297 } 00298 alpha=1.0/fabs(maxval); 00299 00300 /* Infinite norm of d <- g1 = [P(x - g) - x] */ 00301 o2scl_arith::vector_copy(this->dim,this->x,d); 00302 for(size_t i=0;i<this->dim;i++) { 00303 d[i]-=this->gradient[i]; 00304 } 00305 proj(d); 00306 for(size_t i=0;i<this->dim;i++) d[i]-=this->x[i]; 00307 00308 maxval=fabs(d[0]); 00309 for(size_t i=1;i<this->dim;i++) { 00310 if (fabs(d[i])>maxval) { 00311 maxval=fabs(d[i]); 00312 } 00313 } 00314 this->size=fabs(maxval); 00315 00316 m=1; 00317 tail=0; 00318 fvec[tail]=this->f; 00319 00320 return 0; 00321 } 00322 00323 /// Perform an iteration 00324 virtual int iterate() { 00325 double t; 00326 00327 line_search(); 00328 00329 double bk=0.0; 00330 for(size_t i=0;i<this->dim;i++) bk+=s[i]*y[i]; 00331 00332 if (bk>=0.0) { 00333 alpha=alphamax; 00334 } else { 00335 double ak=0.0; 00336 for(size_t i=0;i<this->dim;i++) ak+=s[i]*s[i]; 00337 ak=-ak/bk; 00338 alpha=GSL_MIN(alphamax,GSL_MAX(alphamin,ak)); 00339 } 00340 00341 return 0; 00342 } 00343 00344 /// See if we're finished 00345 virtual int is_optimal() { 00346 if (this->size>tol && this->f>fmin) { 00347 return GSL_CONTINUE; 00348 } 00349 return GSL_SUCCESS; 00350 } 00351 00352 /// Return string denoting type ("ool_mmin_spg") 00353 const char *type() { return "ool_mmin_spg"; } 00354 00355 }; 00356 00357 #ifndef DOXYGENP 00358 } 00359 #endif 00360 00361 #endif 00362
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