![]() |
Object-oriented Scientific Computing Library: Version 0.910
|
00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2011-2012, 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 #ifndef GSL_MROOT_BROYDEN_H 00025 #define GSL_MROOT_BROYDEN_H 00026 00027 #include <o2scl/mroot.h> 00028 #include <o2scl/permutation.h> 00029 #include <o2scl/lu.h> 00030 00031 #ifndef DOXYGENP 00032 namespace o2scl { 00033 #endif 00034 00035 /** \brief Multidimensional root-finding using Broyden's method (GSL) 00036 00037 Experimental. 00038 00039 See \ref Broyden65. 00040 */ 00041 template<class func_t=mm_funct<>, class vec_t=ovector_base, 00042 class alloc_vec_t=ovector, class alloc_t=ovector_alloc, 00043 class mat_t=omatrix_base, class alloc_mat_t=omatrix, 00044 class mat_alloc_t=omatrix_alloc, class jfunc_t=jac_funct<vec_t,mat_t> > 00045 class gsl_mroot_broyden : public mroot<func_t,vec_t,jfunc_t> { 00046 00047 protected: 00048 00049 /// Desc 00050 omatrix H; 00051 00052 /// LU decomposition 00053 omatrix lu; 00054 00055 /// Permutation object for the LU decomposition 00056 permutation perm; 00057 00058 /// Desc 00059 ovector v; 00060 00061 /// Desc 00062 ovector w; 00063 00064 /// Desc 00065 ovector y; 00066 00067 /// Desc 00068 ovector p; 00069 00070 /// Desc 00071 ovector fnew; 00072 00073 /// Desc 00074 ovector x_trial; 00075 00076 /// Desc 00077 double phi; 00078 00079 /// Memory allocation object 00080 alloc_t ao; 00081 00082 /// Stepsize vector 00083 alloc_vec_t dx_int; 00084 00085 /// Function value vector 00086 alloc_vec_t f_int; 00087 00088 /// A pointer to the user-specified function 00089 func_t *user_func; 00090 00091 /// Function values 00092 vec_t *user_f; 00093 00094 /// Initial guess and current solution 00095 vec_t *user_x; 00096 00097 /// Initial and current step 00098 vec_t *user_dx; 00099 00100 /// Number of variables 00101 size_t user_nvar; 00102 00103 /// Size of memory allocated 00104 size_t mem_size; 00105 00106 /// Jacobian 00107 jacobian<func_t,vec_t,mat_t> *ajac; 00108 00109 /** \brief Clear allocated vectors and matrices 00110 00111 This function is called by set() before each solve. 00112 */ 00113 void clear() { 00114 if (mem_size>0) { 00115 lu.set_all(0.0); 00116 for(size_t i=0;i<mem_size;i++) { 00117 perm[i]=0; 00118 } 00119 H.set_all(0.0); 00120 v.set_all(0.0); 00121 w.set_all(0.0); 00122 y.set_all(0.0); 00123 fnew.set_all(0.0); 00124 x_trial.set_all(0.0); 00125 p.set_all(0.0); 00126 } 00127 return; 00128 } 00129 00130 public: 00131 00132 gsl_mroot_broyden() { 00133 mem_size=0; 00134 ajac=&def_jac; 00135 def_jac.epsrel=GSL_SQRT_DBL_EPSILON; 00136 } 00137 00138 /// Default Jacobian object 00139 simple_jacobian<func_t,vec_t,mat_t,alloc_vec_t,alloc_t> def_jac; 00140 00141 virtual ~gsl_mroot_broyden() { 00142 free(); 00143 } 00144 00145 /// Allocate memory 00146 void allocate(size_t n) { 00147 00148 if (n!=mem_size) { 00149 00150 free(); 00151 00152 lu.allocate(n,n); 00153 perm.allocate(n); 00154 H.allocate(n,n); 00155 v.allocate(n); 00156 w.allocate(n); 00157 y.allocate(n); 00158 fnew.allocate(n); 00159 x_trial.allocate(n); 00160 p.allocate(n); 00161 ao.allocate(dx_int,n); 00162 ao.allocate(f_int,n); 00163 00164 mem_size=n; 00165 00166 } 00167 00168 return; 00169 } 00170 00171 /// Euclidean norm 00172 double enorm(size_t nvar, const vec_t &ff) { 00173 double e2=0; 00174 size_t i; 00175 for (i=0;i<nvar;i++) { 00176 double fi=ff[i]; 00177 e2+=fi*fi; 00178 } 00179 return sqrt(e2); 00180 } 00181 00182 /** \brief Set the function, initial guess, and provide vectors to 00183 store function values and stepsize 00184 00185 The initial values of \c f and \c dx are ignored. 00186 */ 00187 void set(func_t &func, size_t nvar, vec_t &x, vec_t &f, vec_t &dx) { 00188 00189 if (nvar!=mem_size) allocate(nvar); 00190 clear(); 00191 00192 user_func=&func; 00193 user_nvar=nvar; 00194 user_x=&x; 00195 user_f=&f; 00196 user_dx=&dx; 00197 00198 size_t i, j; 00199 int signum=0; 00200 00201 func(nvar,x,f); 00202 00203 ajac->set_function(func); 00204 (*ajac)(nvar,x,f,lu); 00205 00206 o2scl_linalg::LU_decomp<>(nvar,lu,perm,signum); 00207 o2scl_linalg::LU_invert<omatrix,omatrix,omatrix_col> 00208 (nvar,lu,perm,H); 00209 00210 for (i=0;i<nvar;i++) { 00211 for (j=0;j<nvar;j++) { 00212 H[i][j]=-H[i][j]; 00213 } 00214 } 00215 00216 dx.set_all(0.0); 00217 00218 phi=enorm(nvar,f); 00219 00220 return; 00221 } 00222 00223 /// Perform an iteration 00224 int iterate() { 00225 00226 func_t &func=*user_func; 00227 size_t nvar=user_nvar; 00228 vec_t &x=*user_x; 00229 vec_t &f=*user_f; 00230 vec_t &dx=*user_dx; 00231 00232 double phi0, phi1, t, lambda; 00233 00234 size_t i, j, iter; 00235 00236 /* p=H f */ 00237 00238 for (i=0;i<nvar;i++) { 00239 double sum=0; 00240 for (j=0;j<nvar;j++) { 00241 sum+=H[i][j]*f[j]; 00242 } 00243 p[i]=sum; 00244 } 00245 00246 t=1.0; 00247 iter=0; 00248 00249 phi0=phi; 00250 00251 bool new_step=true; 00252 while (new_step) { 00253 new_step=false; 00254 00255 for (i=0;i<nvar;i++) { 00256 x_trial[i]=x[i]+t*p[i]; 00257 } 00258 00259 { 00260 int status=func(nvar,x_trial,fnew); 00261 if (status != gsl_success) { 00262 return gsl_ebadfunc; 00263 } 00264 } 00265 00266 phi1=enorm(nvar,fnew); 00267 00268 iter++; 00269 00270 if (phi1>phi0 && iter<10 && t>0.1) { 00271 00272 // [GSL] Full step goes uphill, take a reduced step instead. 00273 double theta=phi1/phi0; 00274 t*=(sqrt(1.0+6.0*theta)-1.0)/(3.0*theta); 00275 new_step=true; 00276 } 00277 00278 if (new_step==false && phi1>phi0) { 00279 00280 // [GSL] Need to recompute Jacobian 00281 int signum=0; 00282 00283 (*ajac)(nvar,x,f,lu); 00284 00285 for (i=0;i<nvar;i++) { 00286 for (j=0;j<nvar;j++) { 00287 lu[i][j]=-lu[i][j]; 00288 } 00289 } 00290 00291 o2scl_linalg::LU_decomp(nvar,lu,perm,signum); 00292 o2scl_linalg::LU_invert<omatrix,omatrix,omatrix_col> 00293 (nvar,lu,perm,H); 00294 o2scl_linalg::LU_solve(nvar,lu,perm,f,p); 00295 00296 t=1.0; 00297 00298 for (i=0;i<nvar;i++) { 00299 x_trial[i]=x[i]+t*p[i]; 00300 } 00301 00302 { 00303 int status=func(nvar,x_trial,fnew); 00304 if (status != gsl_success) { 00305 return gsl_ebadfunc; 00306 } 00307 } 00308 00309 phi1=enorm(nvar,fnew); 00310 } 00311 00312 } 00313 00314 /* y=f'-f */ 00315 for (i=0;i<nvar;i++) { 00316 y[i]=fnew[i]-f[i]; 00317 } 00318 00319 /* v=H y */ 00320 for (i=0;i<nvar;i++) { 00321 double sum=0.0; 00322 for (j=0;j<nvar;j++) { 00323 sum+=H[i][j]*y[j]; 00324 } 00325 v[i]=sum; 00326 } 00327 00328 /* lambda=p dot v */ 00329 lambda=0.0; 00330 for (i=0;i<nvar;i++) { 00331 lambda+=p[i]*v[i]; 00332 } 00333 if (lambda==0) { 00334 O2SCL_ERR2("Approximation to Jacobian has collapsed in ", 00335 "gsl_mroot_broyden::iterate().",gsl_ezerodiv); 00336 } 00337 00338 /* v'=v+t*p */ 00339 for (i=0;i<nvar;i++) { 00340 v[i]=v[i]+t*p[i]; 00341 } 00342 00343 /* w^T=p^T H */ 00344 for (i=0;i<nvar;i++) { 00345 double sum=0; 00346 for (j=0;j<nvar;j++) { 00347 sum+=H[j][i]*p[j]; 00348 } 00349 w[i]=sum; 00350 } 00351 00352 /* Hij -> Hij-(vi wj/lambda) */ 00353 for (i=0;i<nvar;i++) { 00354 double vi=v[i]; 00355 for (j=0;j<nvar;j++) { 00356 double wj=w[j]; 00357 H[i][j]-=vi*wj/lambda; 00358 } 00359 } 00360 00361 /* copy fnew into f */ 00362 vector_copy(nvar,fnew,f); 00363 00364 /* copy x_trial into x */ 00365 vector_copy(nvar,x_trial,x); 00366 for (i=0;i<nvar;i++) { 00367 dx[i]=t*p[i]; 00368 } 00369 00370 phi=phi1; 00371 00372 return gsl_success; 00373 } 00374 00375 /// Desc 00376 virtual int msolve(size_t n, vec_t &x, func_t &func) { 00377 00378 int status; 00379 00380 set(func,n,x,f_int,dx_int); 00381 00382 int iter=0; 00383 00384 do { 00385 iter++; 00386 00387 status=iterate(); 00388 if (status) break; 00389 00390 // ------------------------------------------------------ 00391 // The equivalent of the statement: 00392 // 00393 // status=gsl_multiroot_test_residual(f,this->tol_rel); 00394 00395 double resid=0.0; 00396 for(size_t i=0;i<n;i++) { 00397 resid+=fabs(f_int[i]); 00398 } 00399 if (resid<this->tol_rel) status=gsl_success; 00400 else status=gsl_continue; 00401 00402 // ------------------------------------------------------ 00403 00404 if (this->verbose>0) { 00405 this->print_iter(n,x,f_int,iter,resid,this->tol_rel, 00406 "gsl_mroot_broyden"); 00407 } 00408 00409 } while (status==gsl_continue && iter<this->ntrial); 00410 00411 this->last_ntrial=iter; 00412 00413 if (iter>=this->ntrial) { 00414 if (this->last_conv==0) this->last_conv=gsl_emaxiter; 00415 O2SCL_CONV2_RET("Function gsl_mroot_broyden::msolve() ", 00416 "exceeded max. number of iterations.", 00417 gsl_emaxiter,this->err_nonconv); 00418 } 00419 00420 if (status!=0) { 00421 O2SCL_ERR2_RET("Function iterate() failed in ", 00422 "gsl_mroot_broyden::solve_set().",gsl_efailed); 00423 } 00424 00425 return gsl_success; 00426 } 00427 00428 /// Desc 00429 void free() { 00430 if (mem_size>0) { 00431 H.free(); 00432 lu.free(); 00433 perm.free(); 00434 v.free(); 00435 w.free(); 00436 y.free(); 00437 p.free(); 00438 fnew.free(); 00439 x_trial.free(); 00440 ao.free(dx_int); 00441 ao.free(f_int); 00442 } 00443 return; 00444 } 00445 00446 }; 00447 00448 #ifndef DOXYGENP 00449 } 00450 #endif 00451 00452 #endif
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).