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