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