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