00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 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, class dfunc_t=grad_funct<param_t,vec_t> > 00060 class gsl_mmin_conp : 00061 public gsl_mmin_conf<param_t,func_t,vec_t,alloc_vec_t,alloc_t,dfunc_t> { 00062 00063 public: 00064 00065 /// Perform an iteration 00066 virtual int iterate() { 00067 00068 if (this->dim==0) { 00069 set_err_ret("Memory not allocated in iterate().",gsl_efailed); 00070 } 00071 00072 gsl_vector *x=this->ugx; 00073 gsl_vector *gradient=this->ugg; 00074 gsl_vector *dx=this->udx; 00075 00076 double fa = this->it_min, fb, fc; 00077 double dir; 00078 double stepa = 0.0, stepb, stepc=this->step; 00079 00080 double g1norm; 00081 double pg; 00082 00083 if (this->pnorm == 0.0 || this->g0norm == 0.0) { 00084 gsl_vector_set_zero(dx); 00085 this->it_info=1; 00086 return 0; 00087 } 00088 00089 /* Determine which direction is downhill, +p or -p */ 00090 00091 gsl_blas_ddot (this->p, gradient, &pg); 00092 00093 dir = (pg >= 0.0) ? +1.0 : -1.0; 00094 00095 /* Compute new trial point at x_c= x - step * p, where p is the 00096 current direction */ 00097 00098 take_step (x, this->p, stepc, dir / this->pnorm, this->x1, dx); 00099 00100 /* Evaluate function and gradient at new point xc */ 00101 00102 for(size_t i=0;i<this->dim;i++) { 00103 this->avt5[i]=gsl_vector_get(this->x1,i); 00104 } 00105 (*this->func)(this->dim,this->avt5,fc,*this->params); 00106 00107 if (fc < fa) { 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 //simple_df(this->avt6,gradient); 00117 00118 this->it_info=0; 00119 return gsl_success; 00120 } 00121 00122 /* Do a line minimisation in the region (xa,fa) (xc,fc) to find an 00123 intermediate (xb,fb) satisifying fa > fb < fc. Choose an initial 00124 xb based on parabolic interpolation */ 00125 00126 this->intermediate_point(x, this->p, dir / this->pnorm, pg, 00127 stepa, stepc, fa, fc, this->x1, this->dx1, 00128 gradient, &stepb, &fb); 00129 00130 if (stepb == 0.0) { 00131 this->it_info=2; 00132 return 0; 00133 } 00134 00135 minimize(x,this->p,dir / this->pnorm,stepa,stepb,stepc, fa, 00136 fb, fc, this->tol, this->x1, this->dx1, this->x2, 00137 dx, gradient, &(this->step), &(this->it_min), &g1norm); 00138 00139 gsl_vector_memcpy (x,this->x2); 00140 00141 /* Choose a new conjugate direction for the next step */ 00142 00143 this->iter = (this->iter + 1) % x->size; 00144 00145 if (this->iter == 0) { 00146 gsl_vector_memcpy (this->p, gradient); 00147 this->pnorm = g1norm; 00148 } else { 00149 /* p' = g1 - beta * p */ 00150 double g0g1, beta; 00151 00152 /* g0' = g0 - g1 */ 00153 gsl_blas_daxpy (-1.0, gradient, this->g0); 00154 /* g1g0 = (g0-g1).g1 */ 00155 gsl_blas_ddot(this->g0, gradient, &g0g1); 00156 /* beta = -((g1 - g0).g1)/(g0.g0) */ 00157 beta = g0g1 / (this->g0norm*this->g0norm); 00158 00159 gsl_blas_dscal (-beta, this->p); 00160 gsl_blas_daxpy (1.0, gradient, this->p); 00161 this->pnorm = gsl_blas_dnrm2 (this->p); 00162 } 00163 00164 this->g0norm = g1norm; 00165 gsl_vector_memcpy (this->g0, gradient); 00166 00167 this->it_info=0; 00168 return gsl_success; 00169 } 00170 00171 /// Return string denoting type("gsl_mmin_conp") 00172 virtual const char *type() { return "gsl_mmin_conp";} 00173 00174 }; 00175 00176 #ifndef DOXYGENP 00177 } 00178 #endif 00179 00180 #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