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