00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 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 #ifndef O2SCL_GSL_MMIN_CONP_H 00024 #define O2SCL_GSL_MMIN_CONP_H 00025 00026 #include <gsl/gsl_blas.h> 00027 #include <gsl/gsl_multimin.h> 00028 #include <o2scl/gsl_mmin_conf.h> 00029 00030 #ifndef DOXYGENP 00031 namespace o2scl { 00032 #endif 00033 00034 /** \brief Multidimensional minimization by the Polak-Ribiere 00035 conjugate gradient algorithm (GSL) 00036 00037 The variable multi_min::tolf is used as the maximum value 00038 of the gradient and is \f$ 10^{-4} \f$ be default. 00039 00040 The gsl iterate() function for this minimizer chooses to return 00041 \c GSL_ENOPROG if the iteration fails to make progress without 00042 calling the error handler. This is presumably because the 00043 iterate() function can fail to make progress when the algorithm 00044 has succeeded in finding the minimum. I prefer to return a 00045 non-zero value from a function only in cases where the error 00046 handler will also be called, so the user is clear on what an 00047 "error" means in the local context. Thus if iterate() is failing 00048 to make progress, instead of returning a non-zero value, it sets 00049 the value of \ref it_info to a non-zero value. 00050 00051 \todo A bit of needless copying is required in the function 00052 wrapper to convert from \c gsl_vector to the templated vector 00053 type. This can be fixed. 00054 \todo Document stopping conditions 00055 00056 */ 00057 template<class param_t, class func_t=multi_funct<param_t>, 00058 class vec_t=ovector_view, class alloc_vec_t=ovector, 00059 class alloc_t=ovector_alloc, 00060 class dfunc_t=grad_funct<param_t,ovector_view>, 00061 class auto_grad_t=gradient<param_t,func_t,ovector_view>, 00062 class def_auto_grad_t=simple_grad<param_t,func_t,ovector_view> > 00063 class gsl_mmin_conp : 00064 public gsl_mmin_conf<param_t,func_t,vec_t,alloc_vec_t,alloc_t,dfunc_t, 00065 auto_grad_t,def_auto_grad_t> { 00066 00067 public: 00068 00069 /// Perform an iteration 00070 virtual int iterate() { 00071 00072 if (this->dim==0) { 00073 set_err_ret("Memory not allocated in iterate().",gsl_efailed); 00074 } 00075 00076 gsl_vector *x=this->ugx; 00077 gsl_vector *gradient=this->ugg; 00078 gsl_vector *dx=this->udx; 00079 00080 double fa = this->it_min, fb, fc; 00081 double dir; 00082 double stepa = 0.0, stepb, stepc=this->step; 00083 00084 double g1norm; 00085 double pg; 00086 00087 if (this->pnorm == 0.0 || this->g0norm == 0.0) { 00088 gsl_vector_set_zero(dx); 00089 this->it_info=1; 00090 return 0; 00091 } 00092 00093 /* Determine which direction is downhill, +p or -p */ 00094 00095 gsl_blas_ddot (this->p, gradient, &pg); 00096 00097 dir = (pg >= 0.0) ? +1.0 : -1.0; 00098 00099 /* Compute new trial point at x_c= x - step * p, where p is the 00100 current direction */ 00101 00102 take_step (x, this->p, stepc, dir / this->pnorm, this->x1, dx); 00103 00104 /* Evaluate function and gradient at new point xc */ 00105 00106 for(size_t i=0;i<this->dim;i++) { 00107 this->avt5[i]=gsl_vector_get(this->x1,i); 00108 } 00109 (*this->func)(this->dim,this->avt5,fc,*this->params); 00110 00111 if (fc < fa) { 00112 /* Success, reduced the function value */ 00113 this->step = stepc * 2.0; 00114 this->it_min = fc; 00115 gsl_vector_memcpy (x, this->x1); 00116 00117 for(size_t i=0;i<this->dim;i++) { 00118 this->avt6[i]=gsl_vector_get(this->x1,i); 00119 } 00120 //simple_df(this->avt6,gradient); 00121 00122 this->it_info=0; 00123 return gsl_success; 00124 } 00125 00126 /* Do a line minimisation in the region (xa,fa) (xc,fc) to find an 00127 intermediate (xb,fb) satisifying fa > fb < fc. Choose an initial 00128 xb based on parabolic interpolation */ 00129 00130 this->intermediate_point(x, this->p, dir / this->pnorm, pg, 00131 stepa, stepc, fa, fc, this->x1, this->dx1, 00132 gradient, &stepb, &fb); 00133 00134 if (stepb == 0.0) { 00135 this->it_info=2; 00136 return 0; 00137 } 00138 00139 this->minimize(x,this->p,dir / this->pnorm,stepa,stepb,stepc, fa, 00140 fb, fc, this->tol, this->x1, this->dx1, this->x2, 00141 dx, gradient, &(this->step), &(this->it_min), &g1norm); 00142 00143 gsl_vector_memcpy (x,this->x2); 00144 00145 /* Choose a new conjugate direction for the next step */ 00146 00147 this->iter = (this->iter + 1) % x->size; 00148 00149 if (this->iter == 0) { 00150 gsl_vector_memcpy (this->p, gradient); 00151 this->pnorm = g1norm; 00152 } else { 00153 /* p' = g1 - beta * p */ 00154 double g0g1, beta; 00155 00156 /* g0' = g0 - g1 */ 00157 gsl_blas_daxpy (-1.0, gradient, this->g0); 00158 /* g1g0 = (g0-g1).g1 */ 00159 gsl_blas_ddot(this->g0, gradient, &g0g1); 00160 /* beta = -((g1 - g0).g1)/(g0.g0) */ 00161 beta = g0g1 / (this->g0norm*this->g0norm); 00162 00163 gsl_blas_dscal (-beta, this->p); 00164 gsl_blas_daxpy (1.0, gradient, this->p); 00165 this->pnorm = gsl_blas_dnrm2 (this->p); 00166 } 00167 00168 this->g0norm = g1norm; 00169 gsl_vector_memcpy (this->g0, gradient); 00170 00171 this->it_info=0; 00172 return gsl_success; 00173 } 00174 00175 /// Return string denoting type("gsl_mmin_conp") 00176 virtual const char *type() { return "gsl_mmin_conp";} 00177 00178 }; 00179 00180 /** 00181 \brief An array version of gsl_mmin_conp 00182 00183 \comment 00184 Doxygen has had trouble parsing this template, so I have 00185 simplified it slightly for the documentation below. 00186 \endcomment 00187 */ 00188 #ifdef DOXYGENP 00189 template<class param_t, size_t nv> class gsl_mmin_conp_array : 00190 public gsl_mmin_conp<param_t,func_t,vec_t,alloc_vec_t, 00191 alloc_t,dfunc_t,auto_grad_t,def_auto_grad_t> { 00192 }; 00193 #else 00194 template<class param_t, size_t nv> class gsl_mmin_conp_array : 00195 public gsl_mmin_conp<param_t,multi_vfunct<param_t,nv>, 00196 double[nv],double[nv],array_alloc<double[nv]>,grad_vfunct<param_t,nv>, 00197 gradient_array<param_t,multi_vfunct<param_t,nv>,nv>, 00198 simple_grad_array<param_t,multi_vfunct<param_t,nv>,nv> > { 00199 }; 00200 #endif 00201 00202 00203 #ifndef DOXYGENP 00204 } 00205 #endif 00206 00207 #endif
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