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 template<class vec_t=ovector_base, class alloc_vec_t=ovector, 00062 class alloc_t=ovector_alloc, class mat_t=omatrix_base, 00063 class alloc_mat_t=omatrix, class mat_alloc_t=omatrix_alloc> 00064 class mroot_hybrid_state_t { 00065 00066 public: 00067 00068 mroot_hybrid_state_t() { 00069 dim2=0; 00070 } 00071 00072 /** \brief Allocate memory for a solver with \c n variables 00073 00074 \future Convert to using gsl_alloc_arrays() 00075 */ 00076 int allocate(size_t n) { 00077 00078 ma.allocate(J,n,n); 00079 00080 q=gsl_matrix_calloc(n, n); 00081 if (q==0) { 00082 ma.free(J,n); 00083 O2SCL_ERR2_RET("Failed to allocate space for 'q'", 00084 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00085 } 00086 00087 r=gsl_matrix_calloc(n, n); 00088 if (r==0) { 00089 ma.free(J,n); 00090 gsl_matrix_free(q); 00091 O2SCL_ERR2_RET("Failed to allocate space for 'r'", 00092 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00093 } 00094 00095 tau=gsl_vector_calloc(n); 00096 if (tau==0) { 00097 ma.free(J,n); 00098 gsl_matrix_free(q); 00099 gsl_matrix_free(r); 00100 O2SCL_ERR2_RET("Failed to allocate space for tau", 00101 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00102 } 00103 00104 diag=gsl_vector_calloc(n); 00105 if (diag==0) { 00106 ma.free(J,n); 00107 gsl_matrix_free(q); 00108 gsl_matrix_free(r); 00109 gsl_vector_free(tau); 00110 O2SCL_ERR2_RET("Failed to allocate space for diag", 00111 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00112 } 00113 00114 qtf=gsl_vector_calloc(n); 00115 if (qtf==0) { 00116 ma.free(J,n); 00117 gsl_matrix_free(q); 00118 gsl_matrix_free(r); 00119 gsl_vector_free(tau); 00120 gsl_vector_free(diag); 00121 O2SCL_ERR2_RET("Failed to allocate space for qtf", 00122 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00123 } 00124 00125 newton=gsl_vector_calloc(n); 00126 if (newton==0) { 00127 ma.free(J,n); 00128 gsl_matrix_free(q); 00129 gsl_matrix_free(r); 00130 gsl_vector_free(tau); 00131 gsl_vector_free(diag); 00132 gsl_vector_free(qtf); 00133 O2SCL_ERR2_RET("Failed to allocate space for newton", 00134 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00135 } 00136 00137 gradient=gsl_vector_calloc(n); 00138 if (gradient==0) { 00139 ma.free(J,n); 00140 gsl_matrix_free(q); 00141 gsl_matrix_free(r); 00142 gsl_vector_free(tau); 00143 gsl_vector_free(diag); 00144 gsl_vector_free(qtf); 00145 gsl_vector_free(newton); 00146 O2SCL_ERR2_RET("Failed to allocate space for gradient", 00147 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00148 } 00149 00150 df=gsl_vector_calloc(n); 00151 if (df==0) { 00152 ma.free(J,n); 00153 gsl_matrix_free(q); 00154 gsl_matrix_free(r); 00155 gsl_vector_free(tau); 00156 gsl_vector_free(diag); 00157 gsl_vector_free(qtf); 00158 gsl_vector_free(newton); 00159 gsl_vector_free(gradient); 00160 O2SCL_ERR2_RET("Failed to allocate space for df", 00161 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00162 } 00163 00164 qtdf=gsl_vector_calloc(n); 00165 if (qtdf==0) { 00166 ma.free(J,n); 00167 gsl_matrix_free(q); 00168 gsl_matrix_free(r); 00169 gsl_vector_free(tau); 00170 gsl_vector_free(diag); 00171 gsl_vector_free(qtf); 00172 gsl_vector_free(newton); 00173 gsl_vector_free(gradient); 00174 gsl_vector_free(df); 00175 O2SCL_ERR2_RET("Failed to allocate space for qtdf", 00176 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00177 } 00178 00179 rdx=gsl_vector_calloc(n); 00180 if (rdx==0) { 00181 ma.free(J,n); 00182 gsl_matrix_free(q); 00183 gsl_matrix_free(r); 00184 gsl_vector_free(tau); 00185 gsl_vector_free(diag); 00186 gsl_vector_free(qtf); 00187 gsl_vector_free(newton); 00188 gsl_vector_free(gradient); 00189 gsl_vector_free(df); 00190 gsl_vector_free(qtdf); 00191 O2SCL_ERR2_RET("Failed to allocate space for rdx", 00192 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00193 } 00194 00195 w=gsl_vector_calloc(n); 00196 if (w==0) { 00197 ma.free(J,n); 00198 gsl_matrix_free(q); 00199 gsl_matrix_free(r); 00200 gsl_vector_free(tau); 00201 gsl_vector_free(diag); 00202 gsl_vector_free(qtf); 00203 gsl_vector_free(newton); 00204 gsl_vector_free(gradient); 00205 gsl_vector_free(df); 00206 gsl_vector_free(qtdf); 00207 gsl_vector_free(rdx); 00208 O2SCL_ERR2_RET("Failed to allocate space for w", 00209 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00210 } 00211 00212 v=gsl_vector_calloc(n); 00213 if (v==0) { 00214 ma.free(J,n); 00215 gsl_matrix_free(q); 00216 gsl_matrix_free(r); 00217 gsl_vector_free(tau); 00218 gsl_vector_free(diag); 00219 gsl_vector_free(qtf); 00220 gsl_vector_free(newton); 00221 gsl_vector_free(gradient); 00222 gsl_vector_free(df); 00223 gsl_vector_free(qtdf); 00224 gsl_vector_free(rdx); 00225 gsl_vector_free(w); 00226 O2SCL_ERR2_RET("Failed to allocate space for v", 00227 " in mroot_hybrid_state_t::allocate().",gsl_enomem); 00228 } 00229 00230 dim2=n; 00231 return gsl_success; 00232 } 00233 00234 /// Free allocated memory 00235 int free() { 00236 if (dim2>0) { 00237 gsl_vector_free(v); 00238 gsl_vector_free(w); 00239 gsl_vector_free(rdx); 00240 gsl_vector_free(qtdf); 00241 gsl_vector_free(df); 00242 gsl_vector_free(gradient); 00243 gsl_vector_free(newton); 00244 gsl_vector_free(qtf); 00245 gsl_vector_free(diag); 00246 gsl_vector_free(tau); 00247 gsl_matrix_free(r); 00248 gsl_matrix_free(q); 00249 ma.free(J,dim2); 00250 } 00251 return 0; 00252 } 00253 00254 /// Vector allocator 00255 alloc_t va; 00256 /// Matrix allocator 00257 mat_alloc_t ma; 00258 /// Number of iterations 00259 size_t iter; 00260 /// Compute the number of failures 00261 size_t ncfail; 00262 /// Compute the number of successes 00263 size_t ncsuc; 00264 /// The number of times the actual reduction is less than 0.001 00265 size_t nslow1; 00266 /// The number of times the actual reduction is less than 0.1 00267 size_t nslow2; 00268 /// The norm of the current function value 00269 double fnorm; 00270 /// The limit of the Nuclidean norm 00271 double delta; 00272 /// Jacobian 00273 alloc_mat_t J; 00274 /// Q matrix from QR decomposition 00275 gsl_matrix *q; 00276 /// R matrix from QR decomposition 00277 gsl_matrix *r; 00278 /// The tau vector from QR decomposition 00279 gsl_vector *tau; 00280 /// The diagonal elements 00281 gsl_vector *diag; 00282 /// The value of \f$ Q^T f \f$ 00283 gsl_vector *qtf; 00284 /// The Newton direction 00285 gsl_vector *newton; 00286 /// The gradient direction 00287 gsl_vector *gradient; 00288 /// The change in the function value 00289 gsl_vector *df; 00290 /// The value of \f$ Q^T \cdot \mathrm{df} \f$ 00291 gsl_vector *qtdf; 00292 /// The value of \f$ R \cdot \mathrm{dx} \f$ 00293 gsl_vector *rdx; 00294 /// The value of \f$ w=(Q^T df - R dx)/|dx| \f$ 00295 gsl_vector *w; 00296 /// The value of \f$ v=D^2 dx/|dx| \f$ 00297 gsl_vector *v; 00298 /// Number of variables 00299 size_t dim2; 00300 //@} 00301 }; 00302 00303 /** 00304 \brief Multidimensional root-finding algorithm using 00305 Powell's Hybrid method (GSL) 00306 00307 This is a recasted version of the GSL routines which use a 00308 modified version of Powell's Hybrid method as implemented in the 00309 HYBRJ algorithm in MINPACK (\ref Garbow80). Both the scaled and 00310 unscaled options are available by setting \ref int_scaling (the 00311 scaled version is the default). If derivatives are not provided, 00312 they will be computed automatically. This class provides the 00313 GSL-like interface using allocate(), set() (or set_de() in case 00314 where derivatives are available), iterate(), and free() and 00315 higher-level interfaces, msolve() and msolve_de(), which perform 00316 the solution and the memory allocation automatically. Some 00317 additional checking is performed in case the user calls the 00318 functions out of order (i.e. set() without allocate()). 00319 00320 The functions msolve() and msolve_de() use the condition \f$ 00321 \sum_i |f_i| < \f$ \ref mroot::tolf to determine if the solver has 00322 succeeded. 00323 00324 The original GSL algorithm has been modified to shrink the 00325 stepsize if a proposed step causes the function to return a 00326 non-zero value. This allows the routine to automatically try to 00327 avoid regions where the function is not defined. To return to 00328 the default GSL behavior, set \ref shrink_step to false. 00329 00330 The default method for numerically computing the Jacobian is 00331 from \ref simple_jacobian. This default is identical to the GSL 00332 approach, except that the default value of \ref 00333 simple_jacobian::epsmin is non-zero. See \ref simple_jacobian 00334 for more details. 00335 00336 There is an example for the usage of the multidimensional solver 00337 classes given in <tt>examples/ex_mroot.cpp</tt>, see the \ref 00338 ex_mroot_sect . 00339 00340 \future It's kind of strange that set() sets jac_given to false 00341 and set_de() has to reset it to true. Can this be simplified? 00342 00343 \future Replace scaled_enorm() with scaled_enorm_tvec() 00344 and enorm() with enorm_tvec() everywhere, and then 00345 remove the tvec suffix. 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 mroot_hybrids_base { 00354 00355 #ifndef DOXYGEN_INTERNAL 00356 00357 protected: 00358 00359 /// Compute \f$ R \cdot g \f$ where \c g is the \gradient 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 /// Compute \c w and \c v 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 /// Compute \f$ R \cdot \mathrm{dx} \f$ 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 /// Compute the scaled Euclidean norm 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 /// Compute delta 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 /// Compute the Euclidiean norm 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 /// Compute a trial step and put the result in \c xx_trial 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 /// Compute the change in the function value 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 /// Compute \c diag, the norm of the columns of the Jacobian 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 /// Compute \f$ Q^T \cdot f \f$ 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 /// Update \c diag 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 /** \brief Form appropriate convex combination of the Gauss-Newton 00524 direction and the scaled gradient direction. 00525 00526 Using the Gauss-Newton direction given in \c newton (a vector of 00527 size N), and the gradient direction in \c gradient (also a 00528 vector of size N), this computes 00529 \f[ 00530 \mathrm{pp} = \alpha \mathrm{newton} + \beta \mathrm{gradient} 00531 \f] 00532 */ 00533 void scaled_addition_tvec(size_t N, double alpha, gsl_vector *newton, 00534 double beta, gsl_vector *gradient, vec_t &pp) { 00535 size_t i; 00536 00537 for (i = 0; i < N; i++) { 00538 double ni = gsl_vector_get (newton, i); 00539 double gi = gsl_vector_get (gradient, i); 00540 pp[i]=alpha*ni + beta*gi; 00541 } 00542 } 00543 00544 /** \brief Take a dogleg step 00545 00546 Given the QR decomposition of an \c n by \c n matrix "A", 00547 a diagonal matrix \c diag, a vector \c b, and a positive 00548 number \c delta, this function determines the convex 00549 combination \c x of the Gauss-Newton and scaled gradient 00550 directions which minimizes \f$ A x - b \f$ in the least 00551 squares sense, subject to the restriction that the 00552 Euclidean norm of \f$ d x \f$ is smaller than the 00553 value of \c delta. 00554 */ 00555 int dogleg(size_t n, const gsl_matrix *r, const gsl_vector *qtf, 00556 const gsl_vector *diag, double delta, 00557 gsl_vector *newton, gsl_vector *gradient, 00558 vec_t &p) { 00559 00560 double qnorm, gnorm, sgnorm, bnorm, temp; 00561 00562 // Compute the Gauss-Newton direction 00563 newton_direction(r, qtf, newton); 00564 00565 qnorm=scaled_enorm(diag,newton); 00566 00567 // Test whether the gauss-newton direction is acceptable 00568 if (qnorm <= delta) { 00569 for(size_t i=0;i<n;i++) p[i]=gsl_vector_get(newton,i); 00570 return gsl_success; 00571 } 00572 00573 // Compute the gradient direction 00574 gradient_direction(r, qtf, diag, gradient); 00575 00576 gnorm = enorm(gradient); 00577 00578 // The special case that the scaled gradient is zero 00579 if (gnorm == 0) { 00580 double alpha = delta / qnorm; 00581 double beta = 0; 00582 scaled_addition_tvec(n,alpha, newton, beta, gradient, p); 00583 return gsl_success; 00584 } 00585 00586 // ----------------------------------------------------- 00587 // The next set of statements calculates the point along the 00588 // scaled gradient at which the quadratic is minimized. 00589 00590 minimum_step (gnorm, diag, gradient); 00591 00592 // Use p as temporary space to compute Rg 00593 compute_Rg(n,r,gradient,p); 00594 00595 temp = enorm_tvec(p); 00596 sgnorm = (gnorm / temp) / temp; 00597 00598 // ----------------------------------------------------- 00599 // Test whether the scaled gradient direction is acceptable 00600 00601 if (sgnorm > delta) { 00602 00603 double alpha = 0; 00604 double beta = delta; 00605 00606 scaled_addition_tvec(n,alpha, newton, beta, gradient, p); 00607 00608 return gsl_success; 00609 } 00610 00611 // If not, calculate the point along the dogleg step 00612 // at which the quadratic is minimized 00613 bnorm = enorm (qtf); 00614 00615 double bg = bnorm / gnorm; 00616 double bq = bnorm / qnorm; 00617 double dq = delta / qnorm; 00618 double dq2 = dq * dq; 00619 double sd = sgnorm / delta; 00620 double sd2 = sd * sd; 00621 00622 double t1 = bg * bq * sd; 00623 double u = t1 - dq; 00624 double t2 = t1 - dq * sd2 + sqrt (u * u + (1-dq2) * (1 - sd2)); 00625 00626 double alpha = dq * (1 - sd2) / t2; 00627 double beta = (1 - alpha) * sgnorm; 00628 00629 // Form the appropriate convex combination of the gauss-newton 00630 // direction and the scaled gradient direction. 00631 00632 scaled_addition_tvec(n,alpha, newton, beta, gradient, p); 00633 00634 return gsl_success; 00635 } 00636 00637 /// The user-specified Jacobian 00638 jfunc_t *jac; 00639 00640 /// The automatic Jacobian 00641 jacobian<param_t,func_t,vec_t,mat_t> *ajac; 00642 00643 /// Memory allocator for objects of type \c alloc_vec_t 00644 alloc_t ao; 00645 00646 /// The value of the derivative 00647 alloc_vec_t dx; 00648 00649 /// Trial root 00650 alloc_vec_t x_trial; 00651 00652 /// Trial function value 00653 alloc_vec_t f_trial; 00654 00655 /// The solver state 00656 mroot_hybrid_state_t<vec_t,alloc_vec_t,alloc_t, 00657 mat_t,alloc_mat_t,mat_alloc_t> state; 00658 00659 /// The function parameters 00660 param_t *params; 00661 00662 /// The number of equations and unknowns 00663 size_t dim; 00664 00665 /// True if the jacobian has been given 00666 bool jac_given; 00667 00668 /// The user-specified function 00669 func_t *fnewp; 00670 00671 /// True if "set" has been called 00672 bool set_called; 00673 00674 /// Finish the solution after set() or set_de() has been called 00675 int solve_set(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc) { 00676 00677 int iter=0, status; 00678 00679 do { 00680 iter++; 00681 00682 status=iterate(); 00683 if (status) break; 00684 00685 // ------------------------------------------------------ 00686 // The equivalent of the statement: 00687 // 00688 // status=gsl_multiroot_test_residual(f,this->tolf); 00689 00690 double resid=0.0; 00691 for(size_t i=0;i<nn;i++) { 00692 resid+=fabs(f[i]); 00693 } 00694 if (resid<this->tolf) status=gsl_success; 00695 else status=gsl_continue; 00696 00697 // ------------------------------------------------------ 00698 00699 if (this->verbose>0) { 00700 this->print_iter(nn,x,f,iter,resid,this->tolf, 00701 "gsl_mroot_hybrids"); 00702 } 00703 00704 } while (status==GSL_CONTINUE && iter < this->ntrial); 00705 00706 for(size_t i=0;i<nn;i++) { 00707 xx[i]=x[i]; 00708 } 00709 00710 this->last_ntrial=iter; 00711 00712 if (iter>=this->ntrial) { 00713 if (this->last_conv==0) this->last_conv=gsl_emaxiter; 00714 O2SCL_CONV2_RET("Function gsl_mroot_hybrids::msolve() ", 00715 "exceeded max. number of iterations.", 00716 gsl_emaxiter,this->err_nonconv); 00717 } 00718 00719 if (status!=0) { 00720 O2SCL_ERR_RET("Function iterate() failed in msolve().", 00721 gsl_efailed); 00722 } 00723 00724 return gsl_success; 00725 } 00726 00727 #endif 00728 00729 public: 00730 00731 gsl_mroot_hybrids() { 00732 shrink_step=true; 00733 dim=0; 00734 ajac=&def_jac; 00735 def_jac.epsrel=GSL_SQRT_DBL_EPSILON; 00736 int_scaling=true; 00737 jac_given=false; 00738 set_called=false; 00739 } 00740 00741 virtual ~gsl_mroot_hybrids() {} 00742 00743 /** \brief If true, iterate() will shrink the step-size automatically if 00744 the function returns a non-zero value (default true) 00745 00746 The original GSL behavior can be obtained by setting 00747 this to \c false. 00748 */ 00749 bool shrink_step; 00750 00751 /// If true, use the internal scaling method (default true) 00752 bool int_scaling; 00753 00754 /// Default automatic Jacobian object 00755 simple_jacobian<param_t,func_t,vec_t,mat_t, 00756 alloc_vec_t,alloc_t> def_jac; 00757 00758 /// Set the automatic Jacobian object 00759 virtual int set_jacobian(jacobian<param_t,func_t,vec_t,mat_t> &j) { 00760 ajac=&j; 00761 return 0; 00762 } 00763 00764 /** 00765 \brief The value of the function at the present iteration 00766 00767 \comment 00768 We need this to be public so that the user can see if 00769 iterate() is working 00770 \endcomment 00771 */ 00772 alloc_vec_t f; 00773 00774 /// The present solution 00775 alloc_vec_t x; 00776 00777 /** 00778 \brief Perform an iteration 00779 00780 At the end of the iteration, the current value of the solution 00781 is stored in \ref x. 00782 */ 00783 int iterate() { 00784 00785 if (!set_called) { 00786 O2SCL_ERR2_RET("Function set() not called or most recent call ", 00787 "failed in gsl_mroot_hybrids::iterate().", 00788 gsl_efailed); 00789 } 00790 00791 const double fnorm=state.fnorm; 00792 00793 this->last_conv=0; 00794 00795 gsl_matrix *q=state.q; 00796 gsl_matrix *r=state.r; 00797 gsl_vector *tau=state.tau; 00798 gsl_vector *diag=state.diag; 00799 gsl_vector *qtf=state.qtf; 00800 gsl_vector *df=state.df; 00801 gsl_vector *qtdf=state.qtdf; 00802 gsl_vector *rdx=state.rdx; 00803 gsl_vector *w=state.w; 00804 gsl_vector *v=state.v; 00805 00806 double prered,actred; 00807 double pnorm,fnorm1,fnorm1p; 00808 double ratio; 00809 double p1=0.1, p5=0.5, p001=0.001, p0001=0.0001; 00810 00811 /* Compute qtf=Q^T f */ 00812 00813 compute_qtf_tvec(dim,q,f,qtf); 00814 00815 /* Compute dogleg step */ 00816 00817 dogleg(dim,r,qtf,diag,state.delta,state.newton,state.gradient,dx); 00818 00819 /* Take a trial step */ 00820 00821 compute_trial_step_tvec(dim,x,dx,x_trial); 00822 pnorm=scaled_enorm_tvec(dim,diag,dx); 00823 if (state.iter==1) { 00824 if (pnorm < state.delta) { 00825 state.delta=pnorm; 00826 } 00827 } 00828 00829 /* Evaluate function at x + p */ 00830 00831 int status; 00832 00833 if (shrink_step==false) { 00834 00835 // Evaluate the function at the new point, exit if it fails 00836 00837 status=(*fnewp)(dim,x_trial,f_trial,*params); 00838 00839 if (status != gsl_success) { 00840 this->last_conv=gsl_ebadfunc; 00841 O2SCL_CONV_RET 00842 ("Returned non-zero value in gsl_mroot_hybrids::iterate().", 00843 o2scl::gsl_ebadfunc,this->err_nonconv); 00844 } 00845 00846 } else { 00847 00848 // Evaluate the function at the new point, try to recover 00849 // if it fails 00850 00851 status=(*fnewp)(dim,x_trial,f_trial,*params); 00852 00853 int bit=0; 00854 while(status!=0 && bit<20) { 00855 for(size_t ib=0;ib<dim;ib++) { 00856 x_trial[ib]=(x_trial[ib]+x[ib])/2.0; 00857 } 00858 status=(*fnewp)(dim,x_trial,f_trial,*params); 00859 bit++; 00860 } 00861 00862 // Exit if we weren't able to find a new good point 00863 00864 if (status != gsl_success) { 00865 this->last_conv=gsl_ebadfunc; 00866 O2SCL_CONV_RET 00867 ("No suitable step found in gsl_mroot_hybrids::iterate().", 00868 o2scl::gsl_ebadfunc,this->err_nonconv); 00869 } 00870 00871 } 00872 00873 /* Set df=f_trial-f */ 00874 00875 compute_df_tvec(dim,f_trial,f,df); 00876 00877 /* Compute the scaled actual reduction */ 00878 00879 fnorm1=enorm_tvec(f_trial); 00880 actred=compute_actual_reduction(fnorm,fnorm1); 00881 00882 /* Compute rdx=R dx */ 00883 00884 compute_rdx(dim,r,dx,rdx); 00885 00886 /* Compute the scaled predicted reduction phi1p=|Q^T f + R dx| */ 00887 00888 fnorm1p=enorm_sum(qtf,rdx); 00889 prered=compute_predicted_reduction(fnorm,fnorm1p); 00890 00891 /* Compute the ratio of the actual to predicted reduction */ 00892 00893 if (prered > 0) { 00894 ratio=actred / prered; 00895 } else { 00896 ratio=0; 00897 } 00898 00899 /* Update the step bound */ 00900 00901 if (ratio < p1) { 00902 state.ncsuc=0; 00903 state.ncfail++; 00904 state.delta *= p5; 00905 } else { 00906 state.ncfail=0; 00907 state.ncsuc++; 00908 00909 if (ratio >= p5 || state.ncsuc > 1) { 00910 state.delta=GSL_MAX (state.delta,pnorm / p5); 00911 } 00912 if (fabs (ratio - 1) <= p1) { 00913 state.delta=pnorm / p5; 00914 } 00915 } 00916 00917 /* Test for successful iteration */ 00918 00919 if (ratio >= p0001) { 00920 for(size_t i=0;i<dim;i++) { 00921 x[i]=x_trial[i]; 00922 f[i]=f_trial[i]; 00923 } 00924 state.fnorm=fnorm1; 00925 state.iter++; 00926 } 00927 00928 /* Determine the progress of the iteration */ 00929 00930 state.nslow1++; 00931 if (actred >= p001) { 00932 state.nslow1=0; 00933 } 00934 00935 if (actred >= p1) { 00936 state.nslow2=0; 00937 } 00938 00939 if (state.ncfail==2) { 00940 00941 int rx; 00942 00943 if (jac_given) rx=(*jac)(dim,x,f,state.J,*params); 00944 else rx=(*ajac)(dim,x,f,state.J,*params); 00945 00946 if (rx!=0) { 00947 this->last_conv=gsl_efailed; 00948 O2SCL_CONV_RET 00949 ("Jacobian failed in gsl_mroot_hybrids::iterate().", 00950 gsl_efailed,this->err_nonconv); 00951 } 00952 00953 state.nslow2++; 00954 00955 if (state.iter==1) { 00956 if (int_scaling) { 00957 compute_diag_tvec(dim,state.J,diag); 00958 } 00959 state.delta=compute_delta(dim,diag,x); 00960 } else { 00961 if (int_scaling) { 00962 update_diag_tvec(dim,state.J,diag); 00963 } 00964 } 00965 00966 /* Factorize J into QR decomposition */ 00967 ovector *otau=(ovector *)tau; 00968 omatrix *oq=(omatrix *)q; 00969 omatrix *oor=(omatrix *)r; 00970 o2scl_linalg::QR_decomp(dim,dim,state.J,*otau); 00971 o2scl_linalg::QR_unpack(dim,dim,state.J,*otau,*oq,*oor); 00972 00973 return gsl_success; 00974 } 00975 00976 /* Compute qtdf=Q^T df,w=(Q^T df - R dx)/|dx|, v=D^2 dx/|dx| */ 00977 00978 compute_qtf(q,df,qtdf); 00979 compute_wv(dim,qtdf,rdx,dx,diag,pnorm,w,v); 00980 00981 /* Rank-1 update of the jacobian Q'R'=Q(R + w v^T) */ 00982 00983 gsl_linalg_QR_update(q,r,w,v); 00984 00985 /* No progress as measured by jacobian evaluations */ 00986 00987 if (state.nslow2==5) { 00988 this->last_conv=gsl_enoprogj; 00989 O2SCL_CONV_RET 00990 ("No progress in Jacobian in gsl_mroot_hybrids::iterate().", 00991 gsl_enoprogj,this->err_nonconv); 00992 } 00993 00994 /* No progress as measured by function evaluations */ 00995 00996 if (state.nslow1==10) { 00997 this->last_conv=gsl_enoprog; 00998 O2SCL_CONV_RET 00999 ("No progress in function in gsl_mroot_hybrids::iterate().", 01000 gsl_enoprog,this->err_nonconv); 01001 } 01002 01003 return gsl_success; 01004 } 01005 01006 /// Allocate the memory 01007 int allocate(size_t n) { 01008 int status; 01009 01010 if (dim>0) free(); 01011 01012 ao.allocate(x,n); 01013 ao.allocate(dx,n); 01014 ao.allocate(f,n); 01015 ao.allocate(x_trial,n); 01016 ao.allocate(f_trial,n); 01017 status=state.allocate(n); 01018 01019 for(size_t i=0;i<n;i++) { 01020 x[i]=0.0; 01021 dx[i]=0.0; 01022 f[i]=0.0; 01023 } 01024 01025 dim=n; 01026 01027 return o2scl::gsl_success; 01028 } 01029 01030 /// Free the allocated memory 01031 int free() { 01032 if (dim>0) { 01033 state.free(); 01034 ao.free(x_trial); 01035 ao.free(f_trial); 01036 ao.free(dx); 01037 ao.free(x); 01038 ao.free(f); 01039 dim=0; 01040 } 01041 return o2scl::gsl_success; 01042 } 01043 01044 /// Return the type,\c "gsl_mroot_hybrids". 01045 virtual const char *type() { return "gsl_mroot_hybrids"; } 01046 01047 /** \brief Solve \c func with derivatives \c dfunc using \c x as 01048 an initial guess, returning \c x. 01049 01050 */ 01051 virtual int msolve_de(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc, 01052 jfunc_t &dfunc) { 01053 int status; 01054 01055 status=set_de(nn,xx,ufunc,dfunc,pa); 01056 if (status!=0) { 01057 O2SCL_ERR_RET("Function set() failed in msolve().",gsl_efailed); 01058 } 01059 01060 return solve_set(nn,xx,pa,ufunc); 01061 } 01062 01063 /// Solve \c ufunc using \c xx as an initial guess, returning \c xx. 01064 virtual int msolve(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc) { 01065 01066 int status; 01067 01068 status=set(nn,xx,ufunc,pa); 01069 if (status!=0) { 01070 O2SCL_ERR_RET("Function set() failed in msolve().",gsl_efailed); 01071 } 01072 01073 status=solve_set(nn,xx,pa,ufunc); 01074 for(size_t i=0;i<nn;i++) xx[i]=x[i]; 01075 01076 return status; 01077 } 01078 01079 /** \brief Set the function, the parameters, and the initial guess 01080 */ 01081 int set(size_t nn, vec_t &ax, func_t &ufunc, param_t &pa) { 01082 01083 int status; 01084 01085 if (nn!=dim) { 01086 status=allocate(nn); 01087 if (status!=0) { 01088 O2SCL_ERR_RET("Function allocate() failed in msolve().", 01089 gsl_efailed); 01090 } 01091 } 01092 01093 fnewp=&ufunc; 01094 01095 // Specify function for automatic Jacobian 01096 ajac->set_function(ufunc); 01097 01098 if (dim==0) { 01099 O2SCL_ERR_RET("Memory has not been allocated. Use allocate().", 01100 o2scl::gsl_ebadlen); 01101 } 01102 01103 params=&pa; 01104 01105 // Copy the user-specified solution 01106 for(size_t i=0;i<dim;i++) x[i]=ax[i]; 01107 01108 gsl_matrix *q=state.q; 01109 gsl_matrix *r=state.r; 01110 gsl_vector *tau=state.tau; 01111 gsl_vector *diag=state.diag; 01112 01113 status=ufunc(dim,ax,f,*params); 01114 if (status!=0) { 01115 O2SCL_ERR_RET("Function returned non-zero in set().", 01116 gsl_ebadfunc); 01117 } 01118 01119 if (jac_given) status=(*jac)(dim,ax,f,state.J,*params); 01120 else status=(*ajac)(dim,ax,f,state.J,*params); 01121 01122 if (status!=0) { 01123 O2SCL_ERR_RET("Jacobian failed in set().",gsl_efailed); 01124 } 01125 01126 state.iter=1; 01127 state.fnorm=enorm_tvec(f); 01128 state.ncfail=0; 01129 state.ncsuc=0; 01130 state.nslow1=0; 01131 state.nslow2=0; 01132 01133 for(size_t i=0;i<dim;i++) dx[i]=0.0; 01134 01135 /* Store column norms in diag */ 01136 01137 if (int_scaling) { 01138 compute_diag_tvec(dim,state.J,diag); 01139 } else { 01140 gsl_vector_set_all(diag,1.0); 01141 } 01142 01143 /* Set delta to factor |D x| or to factor if |D x| is zero */ 01144 01145 state.delta=compute_delta(dim,diag,x); 01146 01147 /* Factorize J into QR decomposition */ 01148 01149 ovector *otau=(ovector *)tau; 01150 omatrix *oq=(omatrix *)q; 01151 omatrix *oor=(omatrix *)r; 01152 01153 status=o2scl_linalg::QR_decomp(dim,dim,state.J,*otau); 01154 if (status!=0) { 01155 O2SCL_ERR_RET("QR decomposition failed in set().",gsl_efailed); 01156 } 01157 01158 status=o2scl_linalg::QR_unpack(dim,dim,state.J,*otau,*oq,*oor); 01159 if (status!=0) { 01160 O2SCL_ERR_RET("QR unpacking failed in set().",gsl_efailed); 01161 } 01162 01163 set_called=true; 01164 jac_given=false; 01165 01166 return gsl_success; 01167 } 01168 01169 /** \brief Set the function, the Jacobian, the parameters, 01170 and the initial guess 01171 */ 01172 int set_de(size_t nn, vec_t &ax, func_t &ufunc, jfunc_t &dfunc, 01173 param_t &pa) { 01174 01175 fnewp=&ufunc; 01176 jac=&dfunc; 01177 01178 // Make sure set() uses the right Jacobian 01179 jac_given=true; 01180 01181 int ret=set(nn,ax,ufunc,pa); 01182 01183 // Reset jac_given since set() will set it back to false 01184 jac_given=true; 01185 set_called=true; 01186 01187 return ret; 01188 } 01189 01190 }; 01191 01192 #ifndef DOXYGENP 01193 } 01194 #endif 01195 01196 #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