00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 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 #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 functions mmin() and mmin_de() minimize a given function 00038 until the gradientis smaller than the value of multi_min::tolf 00039 (which defaults to \f$ 10^{-4} \f$ ). 00040 00041 See an example for the usage of this class 00042 in \ref ex_mmin_sect . 00043 00044 \future A bit of needless copying is required in the function 00045 wrapper to convert from \c gsl_vector to the templated vector 00046 type. This can be fixed. 00047 00048 */ 00049 template<class param_t, class func_t=multi_funct<param_t>, 00050 class vec_t=ovector_base, class alloc_vec_t=ovector, 00051 class alloc_t=ovector_alloc, 00052 class dfunc_t=grad_funct<param_t,ovector_base>, 00053 class auto_grad_t=gradient<param_t,func_t,ovector_base>, 00054 class def_auto_grad_t=simple_grad<param_t,func_t,ovector_base> > 00055 class gsl_mmin_conp : 00056 public gsl_mmin_conf<param_t,func_t,vec_t,alloc_vec_t,alloc_t,dfunc_t, 00057 auto_grad_t,def_auto_grad_t> { 00058 00059 public: 00060 00061 /// Perform an iteration 00062 virtual int iterate() { 00063 00064 if (this->dim==0) { 00065 O2SCL_ERR_RET("Memory not allocated in iterate().",gsl_efailed); 00066 } 00067 00068 gsl_vector *x=this->ugx; 00069 gsl_vector *gradient=this->ugg; 00070 gsl_vector *dx=this->udx; 00071 00072 double fa = this->it_min, fb, fc; 00073 double dir; 00074 double stepa = 0.0, stepb, stepc=this->step; 00075 00076 double g1norm; 00077 double pg; 00078 00079 if (this->pnorm == 0.0 || this->g0norm == 0.0) { 00080 gsl_vector_set_zero(dx); 00081 this->last_conv=gsl_enoprog; 00082 O2SCL_CONV2_RET("Either pnorm or g0norm vanished in ", 00083 "in gsl_mmin_conp::iterate().", 00084 gsl_enoprog,this->err_nonconv); 00085 return 0; 00086 } 00087 00088 /* Determine which direction is downhill, +p or -p */ 00089 00090 gsl_blas_ddot (this->p, gradient, &pg); 00091 00092 dir = (pg >= 0.0) ? +1.0 : -1.0; 00093 00094 /* Compute new trial point at x_c= x - step * p, where p is the 00095 current direction */ 00096 00097 take_step (x, this->p, stepc, dir / this->pnorm, this->x1, dx); 00098 00099 /* Evaluate function and gradient at new point xc */ 00100 00101 for(size_t i=0;i<this->dim;i++) { 00102 this->avt5[i]=gsl_vector_get(this->x1,i); 00103 } 00104 (*this->func)(this->dim,this->avt5,fc,*this->params); 00105 00106 if (fc < fa) { 00107 00108 /* Success, reduced the function value */ 00109 this->step = stepc * 2.0; 00110 this->it_min = fc; 00111 gsl_vector_memcpy (x, this->x1); 00112 00113 for(size_t i=0;i<this->dim;i++) { 00114 this->avt6[i]=gsl_vector_get(this->x1,i); 00115 } 00116 if (this->grad_given) { 00117 (*this->grad)(this->dim,this->avt6,this->avt7,*this->params); 00118 } else { 00119 (*this->agrad)(this->dim,this->avt6,this->avt7,*this->params); 00120 } 00121 for(size_t i=0;i<this->dim;i++) { 00122 gsl_vector_set(gradient,i,this->avt7[i]); 00123 } 00124 00125 return gsl_success; 00126 } 00127 00128 /* Do a line minimisation in the region (xa,fa) (xc,fc) to find an 00129 intermediate (xb,fb) satisifying fa > fb < fc. Choose an initial 00130 xb based on parabolic interpolation */ 00131 00132 this->intermediate_point(x, this->p, dir / this->pnorm, pg, 00133 stepa, stepc, fa, fc, this->x1, this->dx1, 00134 gradient, &stepb, &fb); 00135 00136 if (stepb == 0.0) { 00137 this->last_conv=gsl_enoprog; 00138 O2SCL_CONV2_RET("Variable stepb vanished in ", 00139 "in gsl_mmin_conp::iterate().", 00140 gsl_enoprog,this->err_nonconv); 00141 } 00142 00143 this->minimize(x,this->p,dir / this->pnorm,stepa,stepb,stepc, fa, 00144 fb, fc, this->tol, this->x1, this->dx1, this->x2, 00145 dx, gradient, &(this->step), &(this->it_min), &g1norm); 00146 00147 gsl_vector_memcpy (x,this->x2); 00148 00149 /* Choose a new conjugate direction for the next step */ 00150 00151 this->iter = (this->iter + 1) % x->size; 00152 00153 if (this->iter == 0) { 00154 gsl_vector_memcpy (this->p, gradient); 00155 this->pnorm = g1norm; 00156 } else { 00157 /* p' = g1 - beta * p */ 00158 double g0g1, beta; 00159 00160 /* g0' = g0 - g1 */ 00161 gsl_blas_daxpy (-1.0, gradient, this->g0); 00162 /* g1g0 = (g0-g1).g1 */ 00163 gsl_blas_ddot(this->g0, gradient, &g0g1); 00164 /* beta = -((g1 - g0).g1)/(g0.g0) */ 00165 beta = g0g1 / (this->g0norm*this->g0norm); 00166 00167 gsl_blas_dscal (-beta, this->p); 00168 gsl_blas_daxpy (1.0, gradient, this->p); 00169 this->pnorm = gsl_blas_dnrm2 (this->p); 00170 } 00171 00172 this->g0norm = g1norm; 00173 gsl_vector_memcpy (this->g0, gradient); 00174 00175 return gsl_success; 00176 } 00177 00178 /// Return string denoting type("gsl_mmin_conp") 00179 virtual const char *type() { return "gsl_mmin_conp";} 00180 00181 }; 00182 00183 /** 00184 \brief An array version of gsl_mmin_conp 00185 00186 \comment 00187 Doxygen has had trouble parsing this template, so I have 00188 simplified it slightly for the documentation below. 00189 \endcomment 00190 */ 00191 #ifdef DOXYGENP 00192 template<class param_t, size_t nv> class gsl_mmin_conp_array : 00193 public gsl_mmin_conp<param_t,func_t,vec_t,alloc_vec_t, 00194 alloc_t,dfunc_t,auto_grad_t,def_auto_grad_t> { 00195 }; 00196 #else 00197 template<class param_t, size_t nv> class gsl_mmin_conp_array : 00198 public gsl_mmin_conp<param_t,multi_vfunct<param_t,nv>, 00199 double[nv],double[nv],array_alloc<double[nv]>,grad_vfunct<param_t,nv>, 00200 gradient_array<param_t,multi_vfunct<param_t,nv>,nv>, 00201 simple_grad_array<param_t,multi_vfunct<param_t,nv>,nv> > { 00202 }; 00203 #endif 00204 00205 00206 #ifndef DOXYGENP 00207 } 00208 #endif 00209 00210 #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