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