00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 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_PGRAD_H 00028 #define O2SCL_OOL_MMIN_PGRAD_H 00029 00030 #include <gsl/gsl_math.h> 00031 00032 #include <o2scl/text_file.h> 00033 #include <o2scl/ovector_tlate.h> 00034 #include <o2scl/multi_funct.h> 00035 #include <o2scl/ool_constr_mmin.h> 00036 00037 #ifndef DOXYGENP 00038 namespace o2scl { 00039 #endif 00040 00041 /** 00042 \brief Constrained minimization by the projected gradient method (OOL) 00043 00044 \todo Complete the mmin() interface with automatic gradient 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_pgrad : 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 /// Temporary vector 00060 alloc_vec_t xx; 00061 00062 /// Project into feasible region 00063 int proj(vec_t &xt) { 00064 size_t ii, n=this->dim; 00065 00066 for(ii=0;ii<n;ii++) { 00067 xt[ii]=GSL_MAX(this->L[ii],GSL_MIN(xt[ii],this->U[ii])); 00068 } 00069 return 0; 00070 } 00071 00072 /// Line search 00073 int line_search() { 00074 00075 double t, tnew, fx, fxx, dif2, gtd; 00076 00077 fx=this->f; 00078 00079 tnew=1; 00080 00081 do { 00082 00083 t = tnew; 00084 /* xx = x - t df */ 00085 o2scl_arith::vector_copy(this->dim,this->x,xx); 00086 for(size_t i=0;i<this->dim;i++) xx[i]-=t*this->gradient[i]; 00087 proj(xx); 00088 00089 (*this->func)(this->dim,xx,fxx,*this->param); 00090 00091 o2scl_arith::vector_copy(this->dim,xx,this->dx); 00092 for(size_t i=0;i<this->dim;i++) this->dx[i]-=this->x[i]; 00093 00094 dif2=0.0; 00095 for(size_t i=0;i<this->dim;i++) dif2+=this->dx[i]*this->dx[i]; 00096 00097 gtd=0.0; 00098 for(size_t i=0;i<this->dim;i++) { 00099 gtd+=this->gradient[i]*this->dx[i]; 00100 } 00101 00102 tnew = - t * t * gtd / (2*(fxx - this->f - gtd)); 00103 00104 double arg1=sigma1*t; 00105 double arg2=GSL_MIN(sigma2*t,tnew); 00106 tnew=GSL_MAX(arg1,arg2); 00107 00108 // sufficient decrease condition (Armijo) 00109 } while(fxx > fx - (alpha/t) * dif2 ); 00110 00111 o2scl_arith::vector_copy(this->dim,xx,this->x); 00112 this->f = fxx; 00113 00114 (*this->dfunc)(this->dim,this->x,this->gradient,*this->param); 00115 00116 return 0; 00117 } 00118 00119 #endif 00120 00121 public: 00122 00123 ool_mmin_pgrad() { 00124 fmin=-1.0e99; 00125 tol=1.0e-4; 00126 alpha=1.0e-4; 00127 sigma1=0.1; 00128 sigma2=0.9; 00129 this->ntrial=1000; 00130 } 00131 00132 /** 00133 \brief Minimum function value (default \f$ 10^{-99} \f$) 00134 00135 If the function value is below this value, then the algorithm 00136 assumes that the function is not bounded and exits. 00137 */ 00138 double fmin; 00139 /// Tolerance on infinite norm 00140 double tol; 00141 /** 00142 \brief Constant for the sufficient decrease condition 00143 (default \f$ 10^{-4} \f$) 00144 */ 00145 double alpha; 00146 /// Lower bound to the step length reduction 00147 double sigma1; 00148 /// Upper bound to the step length reduction 00149 double sigma2; 00150 00151 /// Allocate memory 00152 virtual int allocate(const size_t n) { 00153 if (this->dim>0) free(); 00154 this->ao.allocate(xx,n); 00155 return ool_constr_mmin<param_t,func_t,dfunc_t,hfunc_t,vec_t,alloc_vec_t, 00156 alloc_t>::allocate(n); 00157 } 00158 00159 /// Free previously allocated memory 00160 virtual int free() { 00161 if (this->dim>0) this->ao.free(xx); 00162 return ool_constr_mmin<param_t,func_t,dfunc_t,hfunc_t,vec_t,alloc_vec_t, 00163 alloc_t>::free(); 00164 } 00165 00166 /// Set the function, the initial guess, and the parameters 00167 virtual int set(func_t &fn, dfunc_t &dfn, vec_t &init, param_t &par) { 00168 00169 int ret=ool_constr_mmin<param_t,func_t,dfunc_t,hfunc_t,vec_t, 00170 alloc_vec_t,alloc_t>::set(fn,dfn,init,par); 00171 00172 // Turn initial guess into a feasible point 00173 proj(this->x); 00174 00175 // Evaluate function and gradient 00176 (*this->func)(this->dim,this->x,this->f,*this->param); 00177 (*this->dfunc)(this->dim,this->x,this->gradient,*this->param); 00178 00179 return 0; 00180 } 00181 00182 /// Restart the minimizer 00183 virtual int restart() { 00184 // Turn x into a feasible point 00185 proj(this->x); 00186 00187 // Evaluate function and gradient 00188 (*this->func)(this->dim,this->x,this->f,*this->param); 00189 (*this->dfunc)(this->dim,this->x,this->gradient,*this->param); 00190 00191 return 0; 00192 } 00193 00194 /// Perform an iteration 00195 virtual int iterate() { 00196 double t; 00197 00198 line_search(); 00199 00200 /* Infinite norm of g1 = d <- [P(x - g) - x] */ 00201 o2scl_arith::vector_copy(this->dim,this->x,xx); 00202 for(size_t i=0;i<this->dim;i++) { 00203 xx[i]-=this->gradient[i]; 00204 } 00205 proj(xx); 00206 for(size_t i=0;i<this->dim;i++) xx[i]-=this->x[i]; 00207 00208 double maxval=fabs(xx[0]); 00209 for(size_t i=1;i<this->dim;i++) { 00210 if (fabs(xx[i])>maxval) { 00211 maxval=fabs(xx[i]); 00212 } 00213 } 00214 this->size=fabs(maxval); 00215 00216 return GSL_SUCCESS; 00217 } 00218 00219 /// See if we're finished 00220 virtual int is_optimal() { 00221 if (this->size>tol && this->f>fmin) { 00222 return GSL_CONTINUE; 00223 } 00224 return GSL_SUCCESS; 00225 } 00226 00227 /// Return string denoting type ("ool_mmin_pgrad") 00228 const char *type() { return "ool_mmin_pgrad"; } 00229 00230 }; 00231 00232 #ifndef DOXYGENP 00233 } 00234 #endif 00235 00236 #endif 00237
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