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_MROOT_HYBRIDS_H 00024 #define O2SCL_GSL_MROOT_HYBRIDS_H 00025 00026 #include <gsl/gsl_multiroots.h> 00027 #include <gsl/gsl_linalg.h> 00028 00029 #include <o2scl/ovector_tlate.h> 00030 #include <o2scl/gsl_mroot_hybrids_b.h> 00031 #include <o2scl/mroot.h> 00032 #include <o2scl/jacobian.h> 00033 #include <o2scl/qr.h> 00034 00035 #ifndef DOXYGENP 00036 namespace o2scl { 00037 #endif 00038 00039 /// State class for \ref gsl_mroot_hybrids 00040 template<class vec_t=ovector_view, class alloc_vec_t=ovector, 00041 class alloc_t=ovector_alloc, class mat_t=omatrix_view, 00042 class alloc_mat_t=omatrix, class mat_alloc_t=omatrix_alloc> 00043 class o2scl_hybrid_state_t { 00044 00045 public: 00046 00047 o2scl_hybrid_state_t() { 00048 dim2=0; 00049 } 00050 00051 /// Allocate memory for a solver with \c n variables 00052 int allocate(size_t n) { 00053 00054 ma.allocate(J,n,n); 00055 00056 q=gsl_matrix_calloc(n, n); 00057 if (q==0) { 00058 ma.free(J,n); 00059 set_err2_ret("Failed to allocate space for q", 00060 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00061 } 00062 00063 r=gsl_matrix_calloc(n, n); 00064 if (r==0) { 00065 ma.free(J,n); 00066 gsl_matrix_free(q); 00067 set_err2_ret("Failed to allocate space for r", 00068 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00069 } 00070 00071 tau=gsl_vector_calloc(n); 00072 if (tau==0) { 00073 ma.free(J,n); 00074 gsl_matrix_free(q); 00075 gsl_matrix_free(r); 00076 set_err2_ret("Failed to allocate space for tau", 00077 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00078 } 00079 00080 diag=gsl_vector_calloc(n); 00081 if (diag==0) { 00082 ma.free(J,n); 00083 gsl_matrix_free(q); 00084 gsl_matrix_free(r); 00085 gsl_vector_free(tau); 00086 set_err2_ret("Failed to allocate space for diag", 00087 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00088 } 00089 00090 qtf=gsl_vector_calloc(n); 00091 if (qtf==0) { 00092 ma.free(J,n); 00093 gsl_matrix_free(q); 00094 gsl_matrix_free(r); 00095 gsl_vector_free(tau); 00096 gsl_vector_free(diag); 00097 set_err2_ret("Failed to allocate space for qtf", 00098 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00099 } 00100 00101 newton=gsl_vector_calloc(n); 00102 if (newton==0) { 00103 ma.free(J,n); 00104 gsl_matrix_free(q); 00105 gsl_matrix_free(r); 00106 gsl_vector_free(tau); 00107 gsl_vector_free(diag); 00108 gsl_vector_free(qtf); 00109 set_err2_ret("Failed to allocate space for newton", 00110 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00111 } 00112 00113 gradient=gsl_vector_calloc(n); 00114 if (gradient==0) { 00115 ma.free(J,n); 00116 gsl_matrix_free(q); 00117 gsl_matrix_free(r); 00118 gsl_vector_free(tau); 00119 gsl_vector_free(diag); 00120 gsl_vector_free(qtf); 00121 gsl_vector_free(newton); 00122 set_err2_ret("Failed to allocate space for gradient", 00123 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00124 } 00125 00126 df=gsl_vector_calloc(n); 00127 if (df==0) { 00128 ma.free(J,n); 00129 gsl_matrix_free(q); 00130 gsl_matrix_free(r); 00131 gsl_vector_free(tau); 00132 gsl_vector_free(diag); 00133 gsl_vector_free(qtf); 00134 gsl_vector_free(newton); 00135 gsl_vector_free(gradient); 00136 set_err2_ret("Failed to allocate space for df", 00137 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00138 } 00139 00140 qtdf=gsl_vector_calloc(n); 00141 if (qtdf==0) { 00142 ma.free(J,n); 00143 gsl_matrix_free(q); 00144 gsl_matrix_free(r); 00145 gsl_vector_free(tau); 00146 gsl_vector_free(diag); 00147 gsl_vector_free(qtf); 00148 gsl_vector_free(newton); 00149 gsl_vector_free(gradient); 00150 gsl_vector_free(df); 00151 set_err2_ret("Failed to allocate space for qtdf", 00152 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00153 } 00154 00155 rdx=gsl_vector_calloc(n); 00156 if (rdx==0) { 00157 ma.free(J,n); 00158 gsl_matrix_free(q); 00159 gsl_matrix_free(r); 00160 gsl_vector_free(tau); 00161 gsl_vector_free(diag); 00162 gsl_vector_free(qtf); 00163 gsl_vector_free(newton); 00164 gsl_vector_free(gradient); 00165 gsl_vector_free(df); 00166 gsl_vector_free(qtdf); 00167 set_err2_ret("Failed to allocate space for rdx", 00168 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00169 } 00170 00171 w=gsl_vector_calloc(n); 00172 if (w==0) { 00173 ma.free(J,n); 00174 gsl_matrix_free(q); 00175 gsl_matrix_free(r); 00176 gsl_vector_free(tau); 00177 gsl_vector_free(diag); 00178 gsl_vector_free(qtf); 00179 gsl_vector_free(newton); 00180 gsl_vector_free(gradient); 00181 gsl_vector_free(df); 00182 gsl_vector_free(qtdf); 00183 gsl_vector_free(rdx); 00184 set_err2_ret("Failed to allocate space for w", 00185 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00186 } 00187 00188 v=gsl_vector_calloc(n); 00189 if (v==0) { 00190 ma.free(J,n); 00191 gsl_matrix_free(q); 00192 gsl_matrix_free(r); 00193 gsl_vector_free(tau); 00194 gsl_vector_free(diag); 00195 gsl_vector_free(qtf); 00196 gsl_vector_free(newton); 00197 gsl_vector_free(gradient); 00198 gsl_vector_free(df); 00199 gsl_vector_free(qtdf); 00200 gsl_vector_free(rdx); 00201 gsl_vector_free(w); 00202 set_err2_ret("Failed to allocate space for v", 00203 " in o2scl_hybrid_state_t::allocate().",gsl_enomem); 00204 } 00205 00206 dim2=n; 00207 return GSL_SUCCESS; 00208 } 00209 00210 /// Free allocated memory 00211 int free() { 00212 if (dim2>0) { 00213 gsl_vector_free(v); 00214 gsl_vector_free(w); 00215 gsl_vector_free(rdx); 00216 gsl_vector_free(qtdf); 00217 gsl_vector_free(df); 00218 gsl_vector_free(gradient); 00219 gsl_vector_free(newton); 00220 gsl_vector_free(qtf); 00221 gsl_vector_free(diag); 00222 gsl_vector_free(tau); 00223 gsl_matrix_free(r); 00224 gsl_matrix_free(q); 00225 ma.free(J,dim2); 00226 } 00227 return 0; 00228 } 00229 00230 /// Vector allocator 00231 alloc_t va; 00232 /// Matrix allocator 00233 mat_alloc_t ma; 00234 /// Number of iterations 00235 size_t iter; 00236 /// Desc 00237 size_t ncfail; 00238 /// Desc 00239 size_t ncsuc; 00240 /// Desc 00241 size_t nslow1; 00242 /// Desc 00243 size_t nslow2; 00244 /// Desc 00245 double fnorm; 00246 /// Desc 00247 double delta; 00248 /// Jacobian 00249 alloc_mat_t J; 00250 /// Q matrix from QR decomposition 00251 gsl_matrix *q; 00252 /// R matrix from QR decomposition 00253 gsl_matrix *r; 00254 /// tau vector from QR decomposition 00255 gsl_vector *tau; 00256 /// Desc 00257 gsl_vector *diag; 00258 /// Desc 00259 gsl_vector *qtf; 00260 /// Desc 00261 gsl_vector *newton; 00262 /// Desc 00263 gsl_vector *gradient; 00264 /// Desc 00265 gsl_vector *df; 00266 /// Desc 00267 gsl_vector *qtdf; 00268 /// Desc 00269 gsl_vector *rdx; 00270 /// Desc 00271 gsl_vector *w; 00272 /// Desc 00273 gsl_vector *v; 00274 /// Number of variables 00275 size_t dim2; 00276 //@} 00277 }; 00278 00279 /** 00280 \brief Multidimensional root-finding algorithm using 00281 Powell's Hybrid method (GSL) 00282 00283 This is a recasted version of the GSL routines which use a 00284 modified version of Powell's Hybrid method as implemented in the 00285 HYBRJ algorithm in MINPACK. Both the scaled and unscaled options 00286 are available by setting \ref int_scaling (the scaled version is 00287 the default). If derivatives are not provided, they will be 00288 computed automatically. This class provides the GSL-like 00289 interface using allocate(), set() (or set_de() in case where 00290 derivatives are available), iterate(), and free() and 00291 higher-level interfaces, msolve() and msolve_de(), which perform 00292 the solution and the memory allocation automatically. Some 00293 additional checking is performed in case the user calls the 00294 functions out of order (i.e. set() without allocate()). 00295 00296 The functions msolve() and msolve_de() use the condition \f$ 00297 \sum_i |f_i| < \f$ \ref mroot::tolf to determine if the solver has 00298 succeeded. 00299 00300 The original GSL algorithm has been modified to shrink the 00301 stepsize if a proposed step causes the function to return a 00302 non-zero value. This allows the routine to automatically try to 00303 avoid regions where the function is not defined. To return to 00304 the default GSL behavior, set \ref shrink_step to false. 00305 00306 The default method for numerically computing the Jacobian 00307 is from \ref simple_jacobian. This default is identical 00308 to the GSL approach, except that the default 00309 value of \ref simple_jacobian::epsmin is non-zero. See 00310 \ref simple_jacobian for more details. 00311 00312 There is an example for the usage of the multidimensional solver 00313 classes given in <tt>examples/ex_mroot.cpp</tt>, see \ref 00314 ex_mroot_sect . 00315 00316 */ 00317 template<class param_t, class func_t=mm_funct<param_t>, 00318 class vec_t=ovector_view, class alloc_vec_t=ovector, 00319 class alloc_t=ovector_alloc, class mat_t=omatrix_view, 00320 class alloc_mat_t=omatrix, class mat_alloc_t=omatrix_alloc, 00321 class jfunc_t=jac_funct<param_t,vec_t,mat_t> > 00322 class gsl_mroot_hybrids : 00323 public mroot<param_t,func_t,vec_t,jfunc_t>, public hybrids_base { 00324 00325 #ifndef DOXYGEN_INTERNAL 00326 00327 protected: 00328 00329 /// Desc 00330 void compute_Rg(size_t N, const gsl_matrix * r, 00331 const gsl_vector * gradient, vec_t &Rg) { 00332 size_t i, j; 00333 00334 for (i = 0; i < N; i++) { 00335 double sum = 0; 00336 00337 for (j = i; j < N; j++) { 00338 double gj = gsl_vector_get (gradient, j); 00339 double rij = gsl_matrix_get (r, i, j); 00340 sum += rij * gj; 00341 } 00342 00343 Rg[i]=sum; 00344 } 00345 } 00346 00347 /// Desc 00348 void compute_wv(size_t n, const gsl_vector * qtdf, 00349 const gsl_vector *rdx, const vec_t &dxx, 00350 const gsl_vector *diag, double pnorm, 00351 gsl_vector * w, gsl_vector * v) { 00352 00353 size_t i; 00354 00355 for (i = 0; i < n; i++) { 00356 00357 double qtdfi = gsl_vector_get (qtdf, i); 00358 double rdxi = gsl_vector_get (rdx, i); 00359 double dxi = dxx[i]; 00360 double diagi = gsl_vector_get (diag, i); 00361 00362 gsl_vector_set (w, i, (qtdfi - rdxi) / pnorm); 00363 gsl_vector_set (v, i, diagi * diagi * dxi / pnorm); 00364 } 00365 } 00366 00367 /// Desc 00368 void compute_rdx(size_t N, const gsl_matrix * r, 00369 const vec_t &dxx, gsl_vector * rdx) { 00370 size_t i, j; 00371 00372 for (i = 0; i < N; i++) { 00373 double sum = 0; 00374 for (j = i; j < N; j++) { 00375 sum += gsl_matrix_get (r, i, j) * dxx[j]; 00376 } 00377 gsl_vector_set (rdx, i, sum); 00378 } 00379 } 00380 00381 /// Desc 00382 double scaled_enorm_tvec(size_t n, const gsl_vector * d, const vec_t &ff) { 00383 double e2 = 0; 00384 size_t i; 00385 for (i = 0; i < n ; i++) { 00386 double fi=ff[i]; 00387 double di= gsl_vector_get(d, i); 00388 double u = di * fi; 00389 e2 += u * u ; 00390 } 00391 return sqrt(e2); 00392 } 00393 00394 /// Desc 00395 double compute_delta(size_t n, gsl_vector * diag, vec_t &xx) { 00396 double Dx = scaled_enorm_tvec(n,diag,xx); 00397 double factor = 100; 00398 00399 return (Dx > 0) ? factor * Dx : factor; 00400 } 00401 00402 /// Desc 00403 double enorm_tvec(const vec_t &ff) { 00404 double e2 = 0; 00405 size_t i; 00406 for (i = 0; i < dim ; i++) { 00407 double fi= ff[i]; 00408 e2 += fi * fi; 00409 } 00410 return sqrt(e2); 00411 } 00412 00413 /// Desc 00414 int compute_trial_step_tvec(size_t N, vec_t &xl, vec_t &dxl, 00415 vec_t &xx_trial) { 00416 size_t i; 00417 00418 for (i = 0; i < N; i++) { 00419 xx_trial[i]=xl[i]+dxl[i]; 00420 } 00421 return 0; 00422 } 00423 00424 /// Desc 00425 int compute_df_tvec(size_t n, const vec_t &ff_trial, 00426 const vec_t &fl, gsl_vector * dfl) { 00427 size_t i; 00428 00429 for (i = 0; i < n; i++) { 00430 double dfi=ff_trial[i]-fl[i]; 00431 gsl_vector_set(dfl,i,dfi); 00432 } 00433 00434 return 0; 00435 } 00436 00437 /// Desc 00438 void compute_diag_tvec(size_t n, const mat_t &J, gsl_vector *diag) { 00439 size_t i, j; 00440 00441 for (j = 0; j < n; j++) { 00442 double sum = 0; 00443 for (i = 0; i < n; i++) { 00444 double Jij = J[i][j]; 00445 sum += Jij * Jij; 00446 } 00447 if (sum == 0) { 00448 sum = 1.0; 00449 } 00450 00451 gsl_vector_set(diag,j,sqrt(sum)); 00452 } 00453 } 00454 00455 /// Desc 00456 void compute_qtf_tvec(size_t N, const gsl_matrix * q, const vec_t &ff, 00457 gsl_vector * qtf) { 00458 size_t i, j; 00459 00460 for (j = 0; j < N; j++) { 00461 double sum = 0; 00462 for (i = 0; i < N; i++) { 00463 sum += gsl_matrix_get (q, i, j)*ff[i]; 00464 } 00465 00466 gsl_vector_set (qtf, j, sum); 00467 } 00468 } 00469 00470 /// Desc 00471 void update_diag_tvec(size_t n, const mat_t &J, gsl_vector * diag) { 00472 size_t i, j; 00473 00474 for (j = 0; j < n; j++) { 00475 double cnorm, diagj, sum = 0; 00476 for (i = 0; i < n; i++) { 00477 double Jij = J[i][j]; 00478 sum += Jij * Jij; 00479 } 00480 if (sum == 0) { 00481 sum = 1.0; 00482 } 00483 00484 cnorm = sqrt (sum); 00485 diagj = gsl_vector_get (diag, j); 00486 00487 if (cnorm > diagj) { 00488 gsl_vector_set (diag, j, cnorm); 00489 } 00490 } 00491 } 00492 00493 /// Desc 00494 void scaled_addition_tvec(size_t N, double alpha, gsl_vector *newton, 00495 double beta, gsl_vector *gradient, vec_t &pp) { 00496 size_t i; 00497 00498 for (i = 0; i < N; i++) { 00499 double ni = gsl_vector_get (newton, i); 00500 double gi = gsl_vector_get (gradient, i); 00501 pp[i]=alpha*ni + beta*gi; 00502 } 00503 } 00504 00505 /// Take a dogleg step 00506 int dogleg(size_t n, const gsl_matrix *r, const gsl_vector *qtf, 00507 const gsl_vector *diag, double delta, 00508 gsl_vector *newton, gsl_vector *gradient, 00509 vec_t &p) { 00510 00511 double qnorm, gnorm, sgnorm, bnorm, temp; 00512 00513 newton_direction(r, qtf, newton); 00514 00515 qnorm = scaled_enorm (diag, newton); 00516 00517 if (qnorm <= delta) { 00518 for(size_t i=0;i<n;i++) p[i]=gsl_vector_get(newton,i); 00519 return GSL_SUCCESS; 00520 } 00521 00522 gradient_direction (r, qtf, diag, gradient); 00523 00524 gnorm = enorm (gradient); 00525 00526 if (gnorm == 0) { 00527 double alpha = delta / qnorm; 00528 double beta = 0; 00529 scaled_addition_tvec(n,alpha, newton, beta, gradient, p); 00530 return GSL_SUCCESS; 00531 } 00532 00533 minimum_step (gnorm, diag, gradient); 00534 00535 /* Use p as temporary space to compute Rg */ 00536 compute_Rg(n,r, gradient, p); 00537 00538 temp = enorm_tvec(p); 00539 sgnorm = (gnorm / temp) / temp; 00540 00541 if (sgnorm > delta) { 00542 double alpha = 0; 00543 double beta = delta; 00544 scaled_addition_tvec(n,alpha, newton, beta, gradient, p); 00545 return GSL_SUCCESS; 00546 } 00547 00548 bnorm = enorm (qtf); 00549 00550 double bg = bnorm / gnorm; 00551 double bq = bnorm / qnorm; 00552 double dq = delta / qnorm; 00553 double dq2 = dq * dq; 00554 double sd = sgnorm / delta; 00555 double sd2 = sd * sd; 00556 00557 double t1 = bg * bq * sd; 00558 double u = t1 - dq; 00559 double t2 = t1 - dq * sd2 + sqrt (u * u + (1-dq2) * (1 - sd2)); 00560 00561 double alpha = dq * (1 - sd2) / t2; 00562 double beta = (1 - alpha) * sgnorm; 00563 00564 scaled_addition_tvec(n,alpha, newton, beta, gradient, p); 00565 00566 return GSL_SUCCESS; 00567 } 00568 00569 /// Pointer to the user-specified Jacobian object 00570 jfunc_t *jac; 00571 00572 /// Pointer to the automatic Jacobian object 00573 jacobian<param_t,func_t,vec_t,mat_t> *ajac; 00574 00575 /// Memory allocator for objects of type \c alloc_vec_t 00576 alloc_t ao; 00577 00578 /// The value of the derivative 00579 alloc_vec_t dx; 00580 00581 /// Trial root 00582 alloc_vec_t x_trial; 00583 00584 /// Trial function value 00585 alloc_vec_t f_trial; 00586 00587 /// The solver state 00588 o2scl_hybrid_state_t<vec_t,alloc_vec_t,alloc_t, 00589 mat_t,alloc_mat_t,mat_alloc_t> state; 00590 00591 /// The function parameters 00592 param_t *params; 00593 00594 /// The number of equations and unknowns 00595 size_t dim; 00596 00597 /// True if the jacobian has been given 00598 bool jac_given; 00599 00600 /// Pointer to the user-specified function 00601 func_t *fnewp; 00602 00603 /// True if "set" has been called 00604 bool set_called; 00605 00606 /// Finish the solution after set() or set_de() has been called 00607 int solve_set(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc) { 00608 int iter=0, status; 00609 00610 do { 00611 iter++; 00612 00613 status=iterate(); 00614 if (status) break; 00615 00616 // ------------------------------------------------------ 00617 // The effective equivalent of the statement: 00618 // 00619 // status=gsl_multiroot_test_residual(f,this->tolf); 00620 00621 double resid=0.0; 00622 for(size_t i=0;i<nn;i++) { 00623 resid+=fabs(f[i]); 00624 } 00625 if (resid<this->tolf) status=gsl_success; 00626 else status=gsl_continue; 00627 00628 // ------------------------------------------------------ 00629 00630 if (this->verbose>0) { 00631 this->print_iter(nn,x,f,iter,resid,this->tolf, 00632 "gsl_mroot_hybrids"); 00633 } 00634 00635 } while (status==GSL_CONTINUE && iter < this->ntrial); 00636 00637 for(size_t i=0;i<nn;i++) { 00638 xx[i]=x[i]; 00639 } 00640 00641 this->last_ntrial=iter; 00642 00643 if (iter>=this->ntrial) { 00644 set_err2_ret("Function gsl_mroot_hybrids::msolve() ", 00645 "exceeded max. number of iterations.", 00646 o2scl::gsl_emaxiter); 00647 } 00648 00649 if (status!=0) { 00650 add_err_ret("Function iterate() failed in msolve().", 00651 gsl_efailed); 00652 } 00653 00654 return gsl_success; 00655 } 00656 00657 #endif 00658 00659 public: 00660 00661 gsl_mroot_hybrids() { 00662 shrink_step=true; 00663 dim=0; 00664 ajac=&def_jac; 00665 def_jac.epsrel=GSL_SQRT_DBL_EPSILON; 00666 int_scaling=true; 00667 jac_given=false; 00668 set_called=false; 00669 } 00670 00671 virtual ~gsl_mroot_hybrids() {} 00672 00673 /** \brief If true, iterate() will shrink the step-size automatically if 00674 the function returns a non-zero value (default true) 00675 00676 The original GSL behavior can be obtained by setting 00677 this to \c false. 00678 */ 00679 bool shrink_step; 00680 00681 /// If true, use the internal scaling method (default true) 00682 bool int_scaling; 00683 00684 /// Default automatic Jacobian object 00685 simple_jacobian<param_t,func_t,vec_t,mat_t, 00686 alloc_vec_t,alloc_t> def_jac; 00687 00688 /// Set the automatic Jacobian object 00689 virtual int set_jacobian(jacobian<param_t,func_t,vec_t,mat_t> &j) { 00690 ajac=&j; 00691 return 0; 00692 } 00693 00694 /** 00695 \brief The value of the function at the present iteration 00696 00697 \comment 00698 We need this to be public so that the user can see if 00699 iterate() is working 00700 \endcomment 00701 00702 */ 00703 alloc_vec_t f; 00704 00705 /// The present solution 00706 alloc_vec_t x; 00707 00708 /** 00709 \brief Perform an iteration 00710 00711 At the end of the iteration, the current value of the solution 00712 is stored in \ref x. 00713 */ 00714 int iterate() { 00715 00716 if (!set_called) { 00717 set_err2_ret("Function set() not called or most recent call ", 00718 "failed in gsl_mroot_hybrids::iterate().", 00719 gsl_efailed); 00720 } 00721 00722 const double fnorm=state.fnorm; 00723 00724 gsl_matrix *q=state.q; 00725 gsl_matrix *r=state.r; 00726 gsl_vector *tau=state.tau; 00727 gsl_vector *diag=state.diag; 00728 gsl_vector *qtf=state.qtf; 00729 gsl_vector *df=state.df; 00730 gsl_vector *qtdf=state.qtdf; 00731 gsl_vector *rdx=state.rdx; 00732 gsl_vector *w=state.w; 00733 gsl_vector *v=state.v; 00734 00735 double prered,actred; 00736 double pnorm,fnorm1,fnorm1p; 00737 double ratio; 00738 double p1=0.1,p5=0.5,p001=0.001,p0001=0.0001; 00739 00740 /* Compute qtf=Q^T f */ 00741 00742 compute_qtf_tvec(dim,q,f,qtf); 00743 00744 /* Compute dogleg step */ 00745 00746 dogleg(dim,r,qtf,diag,state.delta,state.newton,state.gradient,dx); 00747 00748 /* Take a trial step */ 00749 00750 compute_trial_step_tvec(dim,x,dx,x_trial); 00751 pnorm=scaled_enorm_tvec(dim,diag,dx); 00752 if (state.iter==1) { 00753 if (pnorm < state.delta) { 00754 state.delta=pnorm; 00755 } 00756 } 00757 00758 /* Evaluate function at x + p */ 00759 00760 int status; 00761 00762 if (shrink_step==false) { 00763 00764 // Evaluate the function at the new point, exit if it fails 00765 00766 status=(*fnewp)(dim,x_trial,f_trial,*params); 00767 00768 if (status != gsl_success) { 00769 add_err_ret 00770 ("Returned non-zero value in gsl_mroot_hybrids::iterate().", 00771 o2scl::gsl_ebadfunc); 00772 } 00773 00774 } else { 00775 00776 // Evaluate the function at the new point, try to recover 00777 // if it fails 00778 00779 status=(*fnewp)(dim,x_trial,f_trial,*params); 00780 00781 int bit=0; 00782 while(status!=0 && bit<20) { 00783 for(size_t ib=0;ib<dim;ib++) { 00784 x_trial[ib]=(x_trial[ib]+x[ib])/2.0; 00785 } 00786 status=(*fnewp)(dim,x_trial,f_trial,*params); 00787 bit++; 00788 } 00789 00790 // Exit if we weren't able to find a new good point 00791 00792 if (status != gsl_success) { 00793 add_err_ret 00794 ("No suitable step found in gsl_mroot_hybrids::iterate().", 00795 o2scl::gsl_ebadfunc); 00796 } 00797 00798 } 00799 00800 /* Set df=f_trial-f */ 00801 00802 compute_df_tvec(dim,f_trial,f,df); 00803 00804 /* Compute the scaled actual reduction */ 00805 00806 fnorm1=enorm_tvec(f_trial); 00807 actred=compute_actual_reduction(fnorm,fnorm1); 00808 00809 /* Compute rdx=R dx */ 00810 00811 compute_rdx(dim,r,dx,rdx); 00812 00813 /* Compute the scaled predicted reduction phi1p=|Q^T f + R dx| */ 00814 00815 fnorm1p=enorm_sum(qtf,rdx); 00816 prered=compute_predicted_reduction(fnorm,fnorm1p); 00817 00818 /* Compute the ratio of the actual to predicted reduction */ 00819 00820 if (prered > 0) { 00821 ratio=actred / prered; 00822 } else { 00823 ratio=0; 00824 } 00825 00826 /* Update the step bound */ 00827 00828 if (ratio < p1) { 00829 state.ncsuc=0; 00830 state.ncfail++; 00831 state.delta *= p5; 00832 } else { 00833 state.ncfail=0; 00834 state.ncsuc++; 00835 00836 if (ratio >= p5 || state.ncsuc > 1) { 00837 state.delta=GSL_MAX (state.delta,pnorm / p5); 00838 } 00839 if (fabs (ratio - 1) <= p1) { 00840 state.delta=pnorm / p5; 00841 } 00842 } 00843 00844 /* Test for successful iteration */ 00845 00846 if (ratio >= p0001) { 00847 for(size_t i=0;i<dim;i++) { 00848 x[i]=x_trial[i]; 00849 f[i]=f_trial[i]; 00850 } 00851 state.fnorm=fnorm1; 00852 state.iter++; 00853 } 00854 00855 /* Determine the progress of the iteration */ 00856 00857 state.nslow1++; 00858 if (actred >= p001) { 00859 state.nslow1=0; 00860 } 00861 00862 if (actred >= p1) { 00863 state.nslow2=0; 00864 } 00865 00866 if (state.ncfail==2) { 00867 00868 int rx; 00869 00870 if (jac_given) rx=(*jac)(dim,x,f,state.J,*params); 00871 else rx=(*ajac)(dim,x,f,state.J,*params); 00872 00873 if (rx!=0) { 00874 add_err_ret 00875 ("Jacobian failed in gsl_mroot_hybrids::iterate().", 00876 gsl_efailed); 00877 } 00878 00879 state.nslow2++; 00880 00881 if (state.iter==1) { 00882 if (int_scaling) { 00883 compute_diag_tvec(dim,state.J,diag); 00884 } 00885 state.delta=compute_delta(dim,diag,x); 00886 } else { 00887 if (int_scaling) { 00888 update_diag_tvec(dim,state.J,diag); 00889 } 00890 } 00891 00892 /* Factorize J into QR decomposition */ 00893 ovector *otau=(ovector *)tau; 00894 omatrix *oq=(omatrix *)q; 00895 omatrix *oor=(omatrix *)r; 00896 o2scl_linalg::QR_decomp(dim,dim,state.J,*otau); 00897 o2scl_linalg::QR_unpack(dim,dim,state.J,*otau,*oq,*oor); 00898 00899 return gsl_success; 00900 } 00901 00902 /* Compute qtdf=Q^T df,w=(Q^T df - R dx)/|dx|, v=D^2 dx/|dx| */ 00903 00904 compute_qtf(q,df,qtdf); 00905 compute_wv(dim,qtdf,rdx,dx,diag,pnorm,w,v); 00906 00907 /* Rank-1 update of the jacobian Q'R'=Q(R + w v^T) */ 00908 00909 gsl_linalg_QR_update(q,r,w,v); 00910 00911 /* No progress as measured by jacobian evaluations */ 00912 00913 if (state.nslow2==5) { 00914 set_err_ret 00915 ("No progress in Jacobian in gsl_mroot_hybrids::iterate().", 00916 gsl_enoprogj); 00917 } 00918 00919 /* No progress as measured by function evaluations */ 00920 00921 if (state.nslow1==10) { 00922 set_err_ret 00923 ("No progress in function in gsl_mroot_hybrids::iterate().", 00924 gsl_enoprog); 00925 } 00926 00927 return gsl_success; 00928 } 00929 00930 /// Allocate the memory 00931 int allocate(size_t n) { 00932 int status; 00933 00934 if (dim>0) free(); 00935 00936 ao.allocate(x,n); 00937 ao.allocate(dx,n); 00938 ao.allocate(f,n); 00939 ao.allocate(x_trial,n); 00940 ao.allocate(f_trial,n); 00941 status=state.allocate(n); 00942 00943 for(size_t i=0;i<n;i++) { 00944 x[i]=0.0; 00945 dx[i]=0.0; 00946 f[i]=0.0; 00947 } 00948 00949 dim=n; 00950 00951 return o2scl::gsl_success; 00952 } 00953 00954 /// Free the allocated memory 00955 int free() { 00956 if (dim>0) { 00957 state.free(); 00958 ao.free(x_trial); 00959 ao.free(f_trial); 00960 ao.free(dx); 00961 ao.free(x); 00962 ao.free(f); 00963 dim=0; 00964 } 00965 return o2scl::gsl_success; 00966 } 00967 00968 /// Return the type,\c "gsl_mroot_hybrids". 00969 virtual const char *type() { return "gsl_mroot_hybrids"; } 00970 00971 /** \brief Solve \c func with derivatives \c dfunc using \c x as 00972 an initial guess, returning \c x. 00973 00974 */ 00975 virtual int msolve_de(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc, 00976 jfunc_t &dfunc) { 00977 int status; 00978 00979 status=set_de(nn,xx,ufunc,dfunc,pa); 00980 if (status!=0) { 00981 add_err_ret("Function set() failed in msolve().",gsl_efailed); 00982 } 00983 00984 return solve_set(nn,xx,pa,ufunc); 00985 } 00986 00987 /// Solve \c ufunc using \c xx as an initial guess, returning \c xx. 00988 virtual int msolve(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc) { 00989 00990 int status; 00991 00992 status=set(nn,xx,ufunc,pa); 00993 if (status!=0) { 00994 add_err_ret("Function set() failed in msolve().",gsl_efailed); 00995 } 00996 00997 status=solve_set(nn,xx,pa,ufunc); 00998 for(size_t i=0;i<nn;i++) xx[i]=x[i]; 00999 01000 return status; 01001 } 01002 01003 /** \brief Set the function, the parameters, and the initial guess 01004 */ 01005 int set(size_t nn, vec_t &ax, func_t &ufunc, param_t &pa) { 01006 01007 int status; 01008 01009 if (nn!=dim) { 01010 status=allocate(nn); 01011 if (status!=0) { 01012 add_err_ret("Function allocate() failed in msolve().",gsl_efailed); 01013 } 01014 } 01015 01016 fnewp=&ufunc; 01017 01018 // Specify function for automatic Jacobian 01019 ajac->set_function(ufunc); 01020 01021 if (dim==0) { 01022 set_err_ret("Memory has not been allocated. Use allocate().", 01023 o2scl::gsl_ebadlen); 01024 } 01025 01026 params=&pa; 01027 01028 // Copy the user-specified solution 01029 for(size_t i=0;i<dim;i++) x[i]=ax[i]; 01030 01031 gsl_matrix *q=state.q; 01032 gsl_matrix *r=state.r; 01033 gsl_vector *tau=state.tau; 01034 gsl_vector *diag=state.diag; 01035 01036 status=ufunc(dim,ax,f,*params); 01037 if (status!=0) { 01038 set_err_ret("Function returned non-zero in set().",gsl_ebadfunc); 01039 } 01040 01041 if (jac_given) status=(*jac)(dim,ax,f,state.J,*params); 01042 else status=(*ajac)(dim,ax,f,state.J,*params); 01043 01044 if (status!=0) { 01045 add_err_ret("Jacobian failed in set().",gsl_efailed); 01046 } 01047 01048 state.iter=1; 01049 state.fnorm=enorm_tvec(f); 01050 state.ncfail=0; 01051 state.ncsuc=0; 01052 state.nslow1=0; 01053 state.nslow2=0; 01054 01055 for(size_t i=0;i<dim;i++) dx[i]=0.0; 01056 01057 /* Store column norms in diag */ 01058 01059 if (int_scaling) { 01060 compute_diag_tvec(dim,state.J,diag); 01061 } else { 01062 gsl_vector_set_all(diag,1.0); 01063 } 01064 01065 /* Set delta to factor |D x| or to factor if |D x| is zero */ 01066 01067 state.delta=compute_delta(dim,diag,x); 01068 01069 /* Factorize J into QR decomposition */ 01070 01071 ovector *otau=(ovector *)tau; 01072 omatrix *oq=(omatrix *)q; 01073 omatrix *oor=(omatrix *)r; 01074 01075 status=o2scl_linalg::QR_decomp(dim,dim,state.J,*otau); 01076 if (status!=0) { 01077 set_err_ret("QR decomposition failed in set().",gsl_efailed); 01078 } 01079 01080 status=o2scl_linalg::QR_unpack(dim,dim,state.J,*otau,*oq,*oor); 01081 if (status!=0) { 01082 set_err_ret("QR unpacking failed in set().",gsl_efailed); 01083 } 01084 01085 set_called=true; 01086 jac_given=false; 01087 01088 return gsl_success; 01089 } 01090 01091 /** \brief Set the function, the Jacobian, the parameters, 01092 and the initial guess 01093 */ 01094 int set_de(size_t nn, vec_t &ax, func_t &ufunc, jfunc_t &dfunc, 01095 param_t &pa) { 01096 01097 fnewp=&ufunc; 01098 jac=&dfunc; 01099 01100 // Make sure set() uses the right Jacobian 01101 jac_given=true; 01102 01103 int ret=set(nn,ax,ufunc,pa); 01104 01105 // Reset jac_given since set() will set it back to false 01106 jac_given=true; 01107 set_called=true; 01108 01109 return ret; 01110 } 01111 01112 }; 01113 01114 #ifndef DOXYGENP 01115 } 01116 #endif 01117 01118 #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