00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 O2scl is free software; you can redistribute it and/or modify 00009 it under the terms of the GNU General Public License as published by 00010 the Free Software Foundation; either version 3 of the License, or 00011 (at your option) any later version. 00012 00013 O2scl is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 00018 You should have received a copy of the GNU General Public License 00019 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00020 00021 ------------------------------------------------------------------- 00022 */ 00023 #ifndef O2SCL_GSL_MROOT_HYBRIDS_H 00024 #define O2SCL_GSL_MROOT_HYBRIDS_H 00025 00026 #include <gsl/gsl_multiroots.h> 00027 #include <gsl/gsl_linalg.h> 00028 00029 #include <o2scl/ovector_tlate.h> 00030 #include <o2scl/gsl_mroot_hybrids_b.h> 00031 #include <o2scl/mroot.h> 00032 #include <o2scl/jacobian.h> 00033 00034 #ifndef DOXYGENP 00035 namespace o2scl { 00036 #endif 00037 00038 /** 00039 \brief Multidimensional root-finding algorithm using 00040 Powell's Hybrid method (GSL) 00041 00042 This is a recasted version of the GSL routines which use a 00043 modified version of Powell's Hybrid method as implemented in the 00044 HYBRJ algorithm in MINPACK. Both the scaled and unscaled options 00045 are available by setting \ref int_scaling (the scaled version is 00046 the default). If derivatives are not provided, they will be 00047 computed automatically. This class provides the GSL-like 00048 interface using allocate(), set() (or set_de() in case where 00049 derivatives are available), iterate(), and free() and 00050 higher-level interfaces, msolve() and msolve_de(), which perform 00051 the solution automatically. Some additional checking is 00052 performed in case the user calls the functions out of order 00053 (i.e. set() without allocate()). 00054 00055 The functions msolve() and msolve_de() use the condition \f$ 00056 \sum_i |f_i| < \f$ \ref mroot::tolf to determine if the solver has 00057 succeeded. 00058 00059 The original GSL algorithm has been modified to shrink the 00060 stepsize if a proposed step causes the function to return a 00061 non-zero value. This allows the routine to automatically try to 00062 avoid regions where the function is not defined. To return to 00063 the default GSL behavior, set \ref shrink_step to false. 00064 00065 \todo Need to complete generalization to generic matrix 00066 types 00067 */ 00068 template<class param_t, class func_t=mm_funct<param_t>, 00069 class vec_t=ovector_view, class alloc_vec_t=ovector, 00070 class alloc_t=ovector_alloc, 00071 class jfunc_t=jac_funct<param_t,vec_t,omatrix_view> > 00072 class gsl_mroot_hybrids : 00073 public mroot<param_t,func_t,vec_t,jfunc_t>, public hybrids_base { 00074 00075 #ifndef DOXYGEN_INTERNAL 00076 00077 protected: 00078 00079 /// Desc 00080 double enorm_new(const vec_t &ff) { 00081 double e2 = 0 ; 00082 size_t i; 00083 for (i = 0; i < dim ; i++) { 00084 double fi= ff[i]; 00085 e2 += fi * fi ; 00086 } 00087 return sqrt(e2); 00088 } 00089 00090 /// Desc 00091 int compute_trial_step_new(gsl_vector *xl, gsl_vector *dxl, 00092 vec_t &x_trial) { 00093 size_t i, N = xl->size; 00094 00095 for (i = 0; i < N; i++) { 00096 double pi = gsl_vector_get (dxl, i); 00097 double xi = gsl_vector_get (xl, i); 00098 x_trial[i]=xi+pi; 00099 } 00100 return 0; 00101 } 00102 00103 /// Desc 00104 int compute_df_new(const vec_t &f_trial, 00105 const gsl_vector * fl, gsl_vector * dfl) { 00106 size_t i, n = fl->size; 00107 00108 for (i = 0; i < n; i++) { 00109 double dfi =f_trial[i] - gsl_vector_get (fl, i); 00110 gsl_vector_set (dfl, i, dfi); 00111 } 00112 00113 return 0; 00114 } 00115 00116 /// Pointer to the user-specified Jacobian object 00117 jfunc_t *jac; 00118 00119 /// Pointer to the automatic Jacobian object 00120 jacobian<param_t,func_t,vec_t,omatrix_view> *ajac; 00121 00122 /// Memory allocator for objects of type \c alloc_vec_t 00123 alloc_t ao; 00124 00125 /// A structure for \ref gsl_mroot_hybrids 00126 typedef struct 00127 { 00128 size_t iter; 00129 size_t ncfail; 00130 size_t ncsuc; 00131 size_t nslow1; 00132 size_t nslow2; 00133 double fnorm; 00134 double delta; 00135 gsl_matrix *J; 00136 gsl_matrix *q; 00137 gsl_matrix *r; 00138 gsl_vector *tau; 00139 gsl_vector *diag; 00140 gsl_vector *qtf; 00141 gsl_vector *newton; 00142 gsl_vector *gradient; 00143 gsl_vector *df; 00144 gsl_vector *qtdf; 00145 gsl_vector *rdx; 00146 gsl_vector *w; 00147 gsl_vector *v; 00148 } 00149 o2scl_hybrid_state_t; 00150 00151 /// The present solution 00152 gsl_vector *x; 00153 00154 /// The value of the derivative 00155 gsl_vector *dx; 00156 00157 /// The solver state 00158 o2scl_hybrid_state_t *state; 00159 00160 /// The function parameters 00161 param_t *params; 00162 00163 /// The number of equations and unknowns 00164 size_t dim; 00165 00166 /// True if the jacobian has been given 00167 bool jac_given; 00168 00169 /// Pointer to the user-specified function 00170 func_t *fnewp; 00171 00172 /// True if "set" has been called 00173 bool set_called; 00174 00175 /// Finish the solution after set() or set_de() has been called 00176 int solve_set(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc) { 00177 int iter=0, status; 00178 00179 do { 00180 iter++; 00181 00182 status=iterate(xx); 00183 if (status) break; 00184 00185 // ------------------------------------------------------ 00186 // The effective equivalent of the statement: 00187 // 00188 // status=gsl_multiroot_test_residual(f,this->tolf); 00189 00190 double resid=0.0; 00191 for(size_t i=0;i<nn;i++) { 00192 resid+=fabs(gsl_vector_get(f,i)); 00193 } 00194 if (resid<this->tolf) status=gsl_success; 00195 else status=gsl_continue; 00196 00197 // ------------------------------------------------------ 00198 00199 if (this->verbose>0) { 00200 ovector *rx=(ovector *)x; 00201 ovector *rf=(ovector *)f; 00202 this->print_iter(nn,*rx,*rf,iter,resid,this->tolf, 00203 "gsl_mroot_hybrids"); 00204 } 00205 00206 } while (status == GSL_CONTINUE && iter < this->ntrial); 00207 00208 for(size_t i=0;i<nn;i++) { 00209 xx[i]=gsl_vector_get(x,i); 00210 } 00211 00212 this->last_ntrial=iter; 00213 00214 if (iter>=this->ntrial) { 00215 set_err_ret("Function msolve() exceeded max number of iterations.", 00216 o2scl::gsl_emaxiter); 00217 } 00218 00219 if (status!=0) { 00220 add_err_ret("Function iterate() failed in msolve().", 00221 gsl_efailed); 00222 } 00223 00224 return gsl_success; 00225 } 00226 00227 #endif 00228 00229 public: 00230 00231 gsl_mroot_hybrids() { 00232 shrink_step=true; 00233 dim=0; 00234 ajac=&def_jac; 00235 def_jac.epsrel=GSL_SQRT_DBL_EPSILON; 00236 int_scaling=true; 00237 jac_given=false; 00238 set_called=false; 00239 } 00240 00241 virtual ~gsl_mroot_hybrids() {} 00242 00243 /** \brief If true, iterate() will shrink the step-size automatically if 00244 the function returns a non-zero value (default true) 00245 00246 The original GSL behavior can be obtained by setting 00247 this to \c false. 00248 */ 00249 bool shrink_step; 00250 00251 /// If true, use the internal scaling method (default true) 00252 bool int_scaling; 00253 00254 /// Default automatic Jacobian object 00255 simple_jacobian<param_t,func_t,vec_t,omatrix_view, 00256 alloc_vec_t,alloc_t> def_jac; 00257 00258 /// Set the automatic Jacobian object 00259 virtual int set_jacobian(jacobian<param_t,func_t,vec_t,omatrix_view> &j) { 00260 ajac=&j; 00261 return 0; 00262 } 00263 00264 /** 00265 \brief The value of the function at the present iteration 00266 00267 \comment 00268 We need this to be public so that the user can see if 00269 iterate() is working 00270 \endcomment 00271 00272 */ 00273 gsl_vector *f; 00274 00275 /** 00276 \brief Perform an iteration 00277 00278 At the end of the iteration, the current value of the solution 00279 is stored in \c ux. 00280 */ 00281 int iterate(vec_t &ux) { 00282 00283 if (!set_called) { 00284 std::string s2="Function set() not called or most recent call failed"; 00285 s2+="in gsl_mroot_hybrids::iterate()."; 00286 set_err_ret(s2.c_str(),gsl_efailed); 00287 } 00288 00289 const double fnorm=state->fnorm; 00290 00291 gsl_matrix *J=state->J; 00292 gsl_matrix *q=state->q; 00293 gsl_matrix *r=state->r; 00294 gsl_vector *tau=state->tau; 00295 gsl_vector *diag=state->diag; 00296 gsl_vector *qtf=state->qtf; 00297 gsl_vector *df=state->df; 00298 gsl_vector *qtdf=state->qtdf; 00299 gsl_vector *rdx=state->rdx; 00300 gsl_vector *w=state->w; 00301 gsl_vector *v=state->v; 00302 00303 double prered,actred; 00304 double pnorm,fnorm1,fnorm1p; 00305 double ratio; 00306 double p1=0.1,p5=0.5,p001=0.001,p0001=0.0001; 00307 00308 /* Compute qtf=Q^T f */ 00309 00310 compute_qtf(q,f,qtf); 00311 00312 /* Compute dogleg step */ 00313 00314 dogleg(r,qtf,diag,state->delta,state->newton,state->gradient,dx); 00315 00316 /* Take a trial step */ 00317 00318 alloc_vec_t x_trial,f_trial; 00319 ao.allocate(x_trial,dim); 00320 ao.allocate(f_trial,dim); 00321 00322 compute_trial_step_new(x,dx,x_trial); 00323 00324 pnorm=scaled_enorm(diag,dx); 00325 00326 if (state->iter == 1) { 00327 if (pnorm < state->delta) { 00328 state->delta=pnorm; 00329 } 00330 } 00331 00332 /* Evaluate function at x + p */ 00333 00334 int status; 00335 00336 if (shrink_step==false) { 00337 00338 // Evaluate the function at the new point, exit if it fails 00339 00340 status=(*fnewp)(dim,x_trial,f_trial,*params); 00341 00342 if (status != gsl_success) { 00343 ao.free(x_trial); 00344 ao.free(f_trial); 00345 add_err_ret 00346 ("Returned non-zero value in gsl_mroot_hybrids::iterate().", 00347 o2scl::gsl_ebadfunc); 00348 } 00349 00350 } else { 00351 00352 // Evaluate the function at the new point, try to recover 00353 // if it fails 00354 00355 status=(*fnewp)(dim,x_trial,f_trial,*params); 00356 00357 int bit=0; 00358 while(status!=0 && bit<20) { 00359 for(size_t ib=0;ib<dim;ib++) { 00360 x_trial[ib]=(x_trial[ib]+gsl_vector_get(x,ib))/2.0; 00361 } 00362 status=(*fnewp)(dim,x_trial,f_trial,*params); 00363 bit++; 00364 } 00365 00366 // Exit if we weren't able to find a new good point 00367 00368 if (status != gsl_success) { 00369 ao.free(x_trial); 00370 ao.free(f_trial); 00371 add_err_ret 00372 ("No suitable step found in gsl_mroot_hybrids::iterate().", 00373 o2scl::gsl_ebadfunc); 00374 } 00375 00376 } 00377 00378 /* Set df=f_trial-f */ 00379 00380 compute_df_new(f_trial,f,df); 00381 00382 /* Compute the scaled actual reduction */ 00383 00384 fnorm1=enorm_new(f_trial); 00385 00386 actred=compute_actual_reduction(fnorm,fnorm1); 00387 00388 /* Compute rdx=R dx */ 00389 00390 compute_rdx(r,dx,rdx); 00391 00392 /* Compute the scaled predicted reduction phi1p=|Q^T f + R dx| */ 00393 00394 fnorm1p=enorm_sum(qtf,rdx); 00395 00396 prered=compute_predicted_reduction(fnorm,fnorm1p); 00397 00398 /* Compute the ratio of the actual to predicted reduction */ 00399 00400 if (prered > 0) { 00401 ratio=actred / prered; 00402 } else { 00403 ratio=0; 00404 } 00405 00406 /* Update the step bound */ 00407 00408 if (ratio < p1) { 00409 state->ncsuc=0; 00410 state->ncfail++; 00411 state->delta *= p5; 00412 } else { 00413 state->ncfail=0; 00414 state->ncsuc++; 00415 00416 if (ratio >= p5 || state->ncsuc > 1) { 00417 state->delta=GSL_MAX (state->delta,pnorm / p5); 00418 } 00419 if (fabs (ratio - 1) <= p1) { 00420 state->delta=pnorm / p5; 00421 } 00422 } 00423 00424 /* Test for successful iteration */ 00425 00426 if (ratio >= p0001) { 00427 for(size_t i=0;i<dim;i++) { 00428 gsl_vector_set(x,i,x_trial[i]); 00429 gsl_vector_set(f,i,f_trial[i]); 00430 } 00431 state->fnorm=fnorm1; 00432 state->iter++; 00433 } 00434 00435 ao.free(x_trial); 00436 ao.free(f_trial); 00437 00438 /* Determine the progress of the iteration */ 00439 00440 state->nslow1++; 00441 if (actred >= p001) { 00442 state->nslow1=0; 00443 } 00444 00445 if (actred >= p1) { 00446 state->nslow2=0; 00447 } 00448 00449 if (state->ncfail == 2) { 00450 00451 int rx; 00452 00453 alloc_vec_t bx; 00454 alloc_vec_t bf; 00455 ao.allocate(bx,dim); 00456 ao.allocate(bf,dim); 00457 for(size_t ikx=0;ikx<dim;ikx++) { 00458 bx[ikx]=gsl_vector_get(x,ikx); 00459 bf[ikx]=gsl_vector_get(f,ikx); 00460 } 00461 omatrix_view *oJ2=(omatrix_view *)J; 00462 00463 if (jac_given) rx=(*jac)(dim,bx,bf,*oJ2,*params); 00464 else rx=(*ajac)(dim,bx,bf,*oJ2,*params); 00465 00466 ao.free(bx); 00467 ao.free(bf); 00468 00469 if (rx!=0) { 00470 add_err_ret("Jacobian failed in iterate().",gsl_efailed); 00471 } 00472 00473 state->nslow2++; 00474 00475 if (state->iter == 1) { 00476 if (int_scaling) { 00477 compute_diag (J,diag); 00478 } 00479 state->delta=compute_delta (diag, x); 00480 } else { 00481 if (int_scaling) { 00482 update_diag (J,diag); 00483 } 00484 } 00485 00486 /* Factorize J into QR decomposition */ 00487 00488 gsl_linalg_QR_decomp (J,tau); 00489 gsl_linalg_QR_unpack (J,tau,q,r); 00490 00491 return gsl_success; 00492 } 00493 00494 /* Compute qtdf=Q^T df,w=(Q^T df - R dx)/|dx|, v=D^2 dx/|dx| */ 00495 00496 compute_qtf(q,df,qtdf); 00497 00498 compute_wv(qtdf,rdx,dx,diag,pnorm,w,v); 00499 00500 /* Rank-1 update of the jacobian Q'R'=Q(R + w v^T) */ 00501 00502 gsl_linalg_QR_update(q,r,w,v); 00503 00504 /* No progress as measured by jacobian evaluations */ 00505 00506 if (state->nslow2 == 5) { 00507 set_err_ret("No progress in Jacobian in iterate().", 00508 gsl_enoprogj); 00509 } 00510 00511 /* No progress as measured by function evaluations */ 00512 00513 if (state->nslow1 == 10) { 00514 set_err_ret("No progress in function in iterate().", 00515 gsl_enoprog); 00516 } 00517 00518 for(size_t i=0;i<dim;i++) ux[i]=gsl_vector_get(x,i); 00519 00520 return gsl_success; 00521 } 00522 00523 /// Allocate the memory 00524 int allocate(size_t n) { 00525 int status; 00526 00527 if (dim>0) free(); 00528 00529 x=gsl_vector_calloc(n); 00530 if (x == 0) { 00531 set_err_ret("failed to allocate space for x",gsl_enomem); 00532 } 00533 00534 f=gsl_vector_calloc(n); 00535 if (f == 0) { 00536 gsl_vector_free(x); 00537 set_err_ret("failed to allocate space for f",gsl_enomem); 00538 } 00539 dx=gsl_vector_calloc(n); 00540 00541 if (dx == 0) { 00542 gsl_vector_free(x); 00543 gsl_vector_free(f); 00544 set_err_ret("failed to allocate space for dx",gsl_enomem); 00545 } 00546 00547 state=(o2scl_hybrid_state_t *) 00548 malloc(gsl_multiroot_fsolver_hybrids->size); 00549 00550 if (state == 0) { 00551 gsl_vector_free(dx); 00552 gsl_vector_free(x); 00553 gsl_vector_free(f); 00554 set_err_ret("failed to allocate space for multiroot solver state", 00555 o2scl::gsl_enomem); 00556 } 00557 00558 status=(gsl_multiroot_fsolver_hybrids->alloc)(state,n); 00559 00560 if (status != gsl_success) { 00561 (gsl_multiroot_fsolver_hybrids->free)(state); 00562 std::free(state); 00563 gsl_vector_free(dx); 00564 gsl_vector_free(x); 00565 gsl_vector_free(f); 00566 set_err_ret("Solver allocation failed.",status); 00567 } 00568 00569 dim=n; 00570 00571 return o2scl::gsl_success; 00572 } 00573 00574 /// Free the allocated memory 00575 int free() { 00576 if (dim>0) { 00577 (gsl_multiroot_fsolver_hybrids->free)(state); 00578 std::free(state); 00579 gsl_vector_free(dx); 00580 gsl_vector_free(x); 00581 gsl_vector_free(f); 00582 dim=0; 00583 } 00584 return o2scl::gsl_success; 00585 } 00586 00587 /// Return the type,\c "gsl_mroot_hybrids". 00588 virtual const char *type() { return "gsl_mroot_hybrids"; } 00589 00590 /** \brief Solve \c func with derivatives \c dfunc using \c x as 00591 an initial guess,returning \c x. 00592 00593 */ 00594 virtual int msolve_de(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc, 00595 jfunc_t &dfunc) { 00596 int status; 00597 00598 status=set_de(nn,xx,ufunc,dfunc,pa); 00599 if (status!=0) { 00600 add_err_ret("Function set() failed in msolve().",gsl_efailed); 00601 } 00602 00603 return solve_set(nn,xx,pa,ufunc); 00604 } 00605 00606 /// Solve \c ufunc using \c xx as an initial guess,returning \c xx. 00607 virtual int msolve(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc) { 00608 00609 int status; 00610 00611 status=set(nn,xx,ufunc,pa); 00612 if (status!=0) { 00613 add_err_ret("Function set() failed in msolve().",gsl_efailed); 00614 } 00615 00616 return solve_set(nn,xx,pa,ufunc); 00617 } 00618 00619 /** \brief Set the function, the parameters, and the initial guess 00620 */ 00621 int set(size_t nn, vec_t &ax, func_t &ufunc, param_t &pa) { 00622 00623 int status; 00624 00625 if (nn!=dim) { 00626 status=allocate(nn); 00627 if (status!=0) { 00628 add_err_ret("Function allocate() failed in msolve().",gsl_efailed); 00629 } 00630 } 00631 00632 fnewp=&ufunc; 00633 00634 // Specify function for automatic Jacobian 00635 ajac->set_function(ufunc); 00636 00637 if (dim==0) { 00638 set_err_ret("Memory has not been allocated. Use allocate().", 00639 o2scl::gsl_ebadlen); 00640 } 00641 00642 params=&pa; 00643 00644 // We can't use gsl_vector_memcpy here,since ax may have a 00645 // different size than dim. 00646 for(size_t i=0;i<dim;i++) gsl_vector_set(x,i,ax[i]); 00647 00648 o2scl_hybrid_state_t *hstate=(o2scl_hybrid_state_t *)state; 00649 00650 gsl_matrix *J=state->J; 00651 gsl_matrix *q=state->q; 00652 gsl_matrix *r=state->r; 00653 gsl_vector *tau=state->tau; 00654 gsl_vector *diag=state->diag; 00655 00656 alloc_vec_t aaf; 00657 ao.allocate(aaf,dim); 00658 for(size_t ikx=0;ikx<dim;ikx++) { 00659 aaf[ikx]=gsl_vector_get(f,ikx); 00660 } 00661 omatrix_view *oJ=(omatrix_view *)J; 00662 00663 status=ufunc(dim,ax,aaf,*params); 00664 if (status!=0) { 00665 set_err_ret("Function returned non-zero in set().",gsl_ebadfunc); 00666 } 00667 00668 if (jac_given) status=(*jac)(dim,ax,aaf,*oJ,*params); 00669 else status=(*ajac)(dim,ax,aaf,*oJ,*params); 00670 if (status!=0) { 00671 add_err_ret("Jacobian failed in set().",gsl_efailed); 00672 } 00673 00674 for(size_t ikx=0;ikx<dim;ikx++) { 00675 gsl_vector_set(f,ikx,aaf[ikx]); 00676 } 00677 ao.free(aaf); 00678 00679 state->iter=1; 00680 state->fnorm=enorm (f); 00681 state->ncfail=0; 00682 state->ncsuc=0; 00683 state->nslow1=0; 00684 state->nslow2=0; 00685 00686 gsl_vector_set_all(dx,0.0); 00687 00688 /* Store column norms in diag */ 00689 00690 if (int_scaling) { 00691 compute_diag(J,diag); 00692 } else { 00693 gsl_vector_set_all(diag,1.0); 00694 } 00695 00696 /* Set delta to factor |D x| or to factor if |D x| is zero */ 00697 00698 state->delta=compute_delta(diag,x); 00699 00700 /* Factorize J into QR decomposition */ 00701 00702 status=gsl_linalg_QR_decomp(J,tau); 00703 if (status!=0) { 00704 set_err_ret("QR decomposition failed in set().",gsl_efailed); 00705 } 00706 00707 status=gsl_linalg_QR_unpack(J,tau,q,r); 00708 if (status!=0) { 00709 set_err_ret("QR unpacking failed in set().",gsl_efailed); 00710 } 00711 00712 set_called=true; 00713 jac_given=false; 00714 00715 return gsl_success; 00716 } 00717 00718 /** \brief Set the function, the Jacobian, the parameters, 00719 and the initial guess 00720 */ 00721 int set_de(size_t nn, vec_t &ax, func_t &ufunc, jfunc_t &dfunc, 00722 param_t &pa) { 00723 00724 fnewp=&ufunc; 00725 jac=&dfunc; 00726 00727 /// Make sure set() uses the right Jacobian 00728 jac_given=true; 00729 00730 int ret=set(nn,ax,ufunc,pa); 00731 00732 /// Reset jac_given since set() will set it back to false 00733 jac_given=true; 00734 set_called=true; 00735 00736 return ret; 00737 } 00738 00739 }; 00740 00741 #ifndef DOXYGENP 00742 } 00743 #endif 00744 00745 #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