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