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 Internally, there is a little unnecessary copying 00066 back and forth of vectors 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 /// Pointer to the user-specified Jacobian object 00080 jfunc_t *jac; 00081 00082 /// Pointer to the automatic Jacobian object 00083 jacobian<param_t,func_t,vec_t,omatrix_view> *ajac; 00084 00085 /// Memory allocator for objects of type \c alloc_vec_t 00086 alloc_t ao; 00087 00088 /// A structure for \ref gsl_mroot_hybrids 00089 typedef struct 00090 { 00091 size_t iter; 00092 size_t ncfail; 00093 size_t ncsuc; 00094 size_t nslow1; 00095 size_t nslow2; 00096 double fnorm; 00097 double delta; 00098 gsl_matrix *J; 00099 gsl_matrix *q; 00100 gsl_matrix *r; 00101 gsl_vector *tau; 00102 gsl_vector *diag; 00103 gsl_vector *qtf; 00104 gsl_vector *newton; 00105 gsl_vector *gradient; 00106 gsl_vector *x_trial; 00107 gsl_vector *f_trial; 00108 gsl_vector *df; 00109 gsl_vector *qtdf; 00110 gsl_vector *rdx; 00111 gsl_vector *w; 00112 gsl_vector *v; 00113 } 00114 o2scl_hybrid_state_t; 00115 00116 /// The present solution 00117 gsl_vector *x; 00118 00119 /// The value of the derivative 00120 gsl_vector *dx; 00121 00122 /// The solver state 00123 o2scl_hybrid_state_t *state; 00124 00125 /// The function parameters 00126 param_t *params; 00127 00128 /// The number of equations and unknowns 00129 size_t dim; 00130 00131 /// True if the jacobian has been given 00132 bool jac_given; 00133 00134 /// Pointer to the user-specified function 00135 func_t *fnewp; 00136 00137 /// True if "set" has been called 00138 bool set_called; 00139 00140 /// Finish the solution after set() or set_de() has been called 00141 int solve_set(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc) { 00142 int iter=0, status; 00143 00144 do { 00145 iter++; 00146 00147 status=iterate(xx); 00148 if (status) break; 00149 00150 // ------------------------------------------------------ 00151 // The effective equivalent of the statement: 00152 // 00153 // status=gsl_multiroot_test_residual(f,this->tolf); 00154 00155 double resid=0.0; 00156 for(size_t i=0;i<nn;i++) { 00157 resid+=fabs(gsl_vector_get(f,i)); 00158 } 00159 if (resid<this->tolf) status=gsl_success; 00160 else status=gsl_continue; 00161 00162 // ------------------------------------------------------ 00163 00164 if (this->verbose>0) { 00165 ovector *rx=(ovector *)x; 00166 ovector *rf=(ovector *)f; 00167 this->print_iter(nn,*rx,*rf,iter,resid,this->tolf, 00168 "gsl_mroot_hybrids"); 00169 } 00170 00171 } while (status == GSL_CONTINUE && iter < this->ntrial); 00172 00173 for(size_t i=0;i<nn;i++) { 00174 xx[i]=gsl_vector_get(x,i); 00175 } 00176 00177 this->last_ntrial=iter; 00178 00179 if (iter>=this->ntrial) { 00180 set_err_ret("Function msolve() exceeded max number of iterations.", 00181 o2scl::gsl_emaxiter); 00182 } 00183 00184 if (status!=0) { 00185 add_err_ret("Function iterate() failed in msolve().", 00186 gsl_efailed); 00187 } 00188 00189 return gsl_success; 00190 } 00191 00192 #endif 00193 00194 public: 00195 00196 gsl_mroot_hybrids() { 00197 shrink_step=true; 00198 dim=0; 00199 ajac=&def_jac; 00200 def_jac.epsrel=GSL_SQRT_DBL_EPSILON; 00201 int_scaling=true; 00202 jac_given=false; 00203 set_called=false; 00204 } 00205 00206 virtual ~gsl_mroot_hybrids() {} 00207 00208 /** \brief If true, iterate() will shrink the step-size automatically if 00209 the function returns a non-zero value (default true) 00210 00211 The original GSL behavior can be obtained by setting 00212 this to \c false. 00213 */ 00214 bool shrink_step; 00215 00216 /// If true, use the internal scaling method (default true) 00217 bool int_scaling; 00218 00219 /// Default automatic Jacobian object 00220 simple_jacobian<param_t,func_t,vec_t,omatrix_view, 00221 alloc_vec_t,alloc_t> def_jac; 00222 00223 /// Set the automatic Jacobian object 00224 virtual int set_jacobian(jacobian<param_t,func_t,vec_t,omatrix_view> &j) { 00225 ajac=&j; 00226 return 0; 00227 } 00228 00229 /** 00230 \brief The value of the function at the present iteration 00231 00232 \comment 00233 We need this to be public so that the user can see if 00234 iterate() is working 00235 \endcomment 00236 00237 */ 00238 gsl_vector *f; 00239 00240 /** 00241 \brief Perform an iteration 00242 00243 At the end of the iteration, the current value of the solution 00244 is stored in \c ux. 00245 */ 00246 int iterate(vec_t &ux) { 00247 00248 if (!set_called) { 00249 std::string s="Function set() not called or most recent call failed"; 00250 s+="in gsl_mroot_hybrids::iterate()."; 00251 set_err_ret(s.c_str(),gsl_efailed); 00252 } 00253 00254 const double fnorm=state->fnorm; 00255 00256 gsl_matrix *J=state->J; 00257 gsl_matrix *q=state->q; 00258 gsl_matrix *r=state->r; 00259 gsl_vector *tau=state->tau; 00260 gsl_vector *diag=state->diag; 00261 gsl_vector *qtf=state->qtf; 00262 gsl_vector *x_trial=state->x_trial; 00263 gsl_vector *f_trial=state->f_trial; 00264 gsl_vector *df=state->df; 00265 gsl_vector *qtdf=state->qtdf; 00266 gsl_vector *rdx=state->rdx; 00267 gsl_vector *w=state->w; 00268 gsl_vector *v=state->v; 00269 00270 double prered,actred; 00271 double pnorm,fnorm1,fnorm1p; 00272 double ratio; 00273 double p1=0.1,p5=0.5,p001=0.001,p0001=0.0001; 00274 00275 /* Compute qtf=Q^T f */ 00276 00277 compute_qtf(q,f,qtf); 00278 00279 /* Compute dogleg step */ 00280 00281 dogleg(r,qtf,diag,state->delta,state->newton,state->gradient,dx); 00282 00283 /* Take a trial step */ 00284 00285 compute_trial_step(x,dx,state->x_trial); 00286 00287 pnorm=scaled_enorm(diag,dx); 00288 00289 if (state->iter == 1) { 00290 if (pnorm < state->delta) { 00291 state->delta=pnorm; 00292 } 00293 } 00294 00295 /* Evaluate function at x + p */ 00296 00297 alloc_vec_t avx,avf; 00298 ao.allocate(avx,dim); 00299 ao.allocate(avf,dim); 00300 for(size_t il=0;il<dim;il++) { 00301 avx[il]=gsl_vector_get(x_trial,il); 00302 avf[il]=gsl_vector_get(f_trial,il); 00303 } 00304 00305 int status; 00306 00307 if (shrink_step==false) { 00308 00309 // Evaluate the function at the new point, exit if it fails 00310 00311 status=(*fnewp)(dim,avx,avf,*params); 00312 00313 if (status != gsl_success) { 00314 ao.free(avx); 00315 ao.free(avf); 00316 add_err_ret 00317 ("Returned non-zero value in gsl_mroot_hybrids::iterate().", 00318 o2scl::gsl_ebadfunc); 00319 } 00320 00321 } else { 00322 00323 // Evaluate the function at the new point, try to recover 00324 // if it fails 00325 00326 status=(*fnewp)(dim,avx,avf,*params); 00327 00328 int bit=0; 00329 while(status!=0 && bit<20) { 00330 for(size_t ib=0;ib<dim;ib++) { 00331 avx[ib]=(avx[ib]+gsl_vector_get(x,ib))/2.0; 00332 } 00333 status=(*fnewp)(dim,avx,avf,*params); 00334 bit++; 00335 } 00336 00337 // Exit if we weren't able to find a new good point 00338 00339 if (status != gsl_success) { 00340 ao.free(avx); 00341 ao.free(avf); 00342 add_err_ret 00343 ("No suitable step found in gsl_mroot_hybrids::iterate().", 00344 o2scl::gsl_ebadfunc); 00345 } 00346 00347 } 00348 00349 for(size_t il=0;il<dim;il++) { 00350 gsl_vector_set(x_trial,il,avx[il]); 00351 gsl_vector_set(f_trial,il,avf[il]); 00352 } 00353 ao.free(avx); 00354 ao.free(avf); 00355 00356 /* Set df=f_trial-f */ 00357 00358 compute_df(f_trial,f,df); 00359 00360 /* Compute the scaled actual reduction */ 00361 00362 fnorm1=enorm(f_trial); 00363 00364 actred=compute_actual_reduction(fnorm,fnorm1); 00365 00366 /* Compute rdx=R dx */ 00367 00368 compute_rdx(r,dx,rdx); 00369 00370 /* Compute the scaled predicted reduction phi1p=|Q^T f + R dx| */ 00371 00372 fnorm1p=enorm_sum(qtf,rdx); 00373 00374 prered=compute_predicted_reduction(fnorm,fnorm1p); 00375 00376 /* Compute the ratio of the actual to predicted reduction */ 00377 00378 if (prered > 0) { 00379 ratio=actred / prered; 00380 } else { 00381 ratio=0; 00382 } 00383 00384 /* Update the step bound */ 00385 00386 if (ratio < p1) { 00387 state->ncsuc=0; 00388 state->ncfail++; 00389 state->delta *= p5; 00390 } else { 00391 state->ncfail=0; 00392 state->ncsuc++; 00393 00394 if (ratio >= p5 || state->ncsuc > 1) { 00395 state->delta=GSL_MAX (state->delta,pnorm / p5); 00396 } 00397 if (fabs (ratio - 1) <= p1) { 00398 state->delta=pnorm / p5; 00399 } 00400 } 00401 00402 /* Test for successful iteration */ 00403 00404 if (ratio >= p0001) { 00405 gsl_vector_memcpy ( x,x_trial); 00406 gsl_vector_memcpy ( f,f_trial); 00407 state->fnorm=fnorm1; 00408 state->iter++; 00409 } 00410 00411 /* Determine the progress of the iteration */ 00412 00413 state->nslow1++; 00414 if (actred >= p001) { 00415 state->nslow1=0; 00416 } 00417 00418 if (actred >= p1) { 00419 state->nslow2=0; 00420 } 00421 00422 if (state->ncfail == 2) { 00423 00424 int rx; 00425 00426 alloc_vec_t bx; 00427 alloc_vec_t bf; 00428 ao.allocate(bx,dim); 00429 ao.allocate(bf,dim); 00430 for(size_t ikx=0;ikx<dim;ikx++) { 00431 bx[ikx]=gsl_vector_get(x,ikx); 00432 bf[ikx]=gsl_vector_get(f,ikx); 00433 } 00434 omatrix_view *oJ2=(omatrix_view *)J; 00435 00436 if (jac_given) rx=(*jac)(dim,bx,bf,*oJ2,*params); 00437 else rx=(*ajac)(dim,bx,bf,*oJ2,*params); 00438 00439 ao.free(bx); 00440 ao.free(bf); 00441 00442 if (rx!=0) { 00443 add_err_ret("Jacobian failed in iterate().",gsl_efailed); 00444 } 00445 00446 state->nslow2++; 00447 00448 if (state->iter == 1) { 00449 if (int_scaling) { 00450 compute_diag (J,diag); 00451 } 00452 state->delta=compute_delta (diag, x); 00453 } else { 00454 if (int_scaling) { 00455 update_diag (J,diag); 00456 } 00457 } 00458 00459 /* Factorize J into QR decomposition */ 00460 00461 gsl_linalg_QR_decomp (J,tau); 00462 gsl_linalg_QR_unpack (J,tau,q,r); 00463 00464 return gsl_success; 00465 } 00466 00467 /* Compute qtdf=Q^T df,w=(Q^T df - R dx)/|dx|, v=D^2 dx/|dx| */ 00468 00469 compute_qtf(q,df,qtdf); 00470 00471 compute_wv(qtdf,rdx,dx,diag,pnorm,w,v); 00472 00473 /* Rank-1 update of the jacobian Q'R'=Q(R + w v^T) */ 00474 00475 gsl_linalg_QR_update(q,r,w,v); 00476 00477 /* No progress as measured by jacobian evaluations */ 00478 00479 if (state->nslow2 == 5) { 00480 set_err_ret("No progress in Jacobian in iterate().", 00481 gsl_enoprogj); 00482 } 00483 00484 /* No progress as measured by function evaluations */ 00485 00486 if (state->nslow1 == 10) { 00487 set_err_ret("No progress in function in iterate().", 00488 gsl_enoprog); 00489 } 00490 00491 for(size_t i=0;i<dim;i++) ux[i]=gsl_vector_get(x,i); 00492 00493 return gsl_success; 00494 } 00495 00496 /// Allocate the memory 00497 int allocate(size_t n) { 00498 int status; 00499 00500 if (dim>0) free(); 00501 00502 x=gsl_vector_calloc(n); 00503 if (x == 0) { 00504 set_err_ret("failed to allocate space for x",gsl_enomem); 00505 } 00506 00507 f=gsl_vector_calloc(n); 00508 if (f == 0) { 00509 gsl_vector_free(x); 00510 set_err_ret("failed to allocate space for f",gsl_enomem); 00511 } 00512 dx=gsl_vector_calloc(n); 00513 00514 if (dx == 0) { 00515 gsl_vector_free(x); 00516 gsl_vector_free(f); 00517 set_err_ret("failed to allocate space for dx",gsl_enomem); 00518 } 00519 00520 state=(o2scl_hybrid_state_t *) 00521 malloc(gsl_multiroot_fsolver_hybrids->size); 00522 00523 if (state == 0) { 00524 gsl_vector_free(dx); 00525 gsl_vector_free(x); 00526 gsl_vector_free(f); 00527 set_err_ret("failed to allocate space for multiroot solver state", 00528 o2scl::gsl_enomem); 00529 } 00530 00531 status=(gsl_multiroot_fsolver_hybrids->alloc)(state,n); 00532 00533 if (status != gsl_success) { 00534 (gsl_multiroot_fsolver_hybrids->free)(state); 00535 std::free(state); 00536 gsl_vector_free(dx); 00537 gsl_vector_free(x); 00538 gsl_vector_free(f); 00539 set_err_ret("Solver allocation failed.",status); 00540 } 00541 00542 dim=n; 00543 00544 return o2scl::gsl_success; 00545 } 00546 00547 /// Free the allocated memory 00548 int free() { 00549 if (dim>0) { 00550 (gsl_multiroot_fsolver_hybrids->free)(state); 00551 std::free(state); 00552 gsl_vector_free(dx); 00553 gsl_vector_free(x); 00554 gsl_vector_free(f); 00555 dim=0; 00556 } 00557 return o2scl::gsl_success; 00558 } 00559 00560 /// Return the type,\c "gsl_mroot_hybrids". 00561 virtual const char *type() { return "gsl_mroot_hybrids"; } 00562 00563 /** \brief Solve \c func with derivatives \c dfunc using \c x as 00564 an initial guess,returning \c x. 00565 00566 */ 00567 virtual int msolve_de(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc, 00568 jfunc_t &dfunc) { 00569 int status; 00570 00571 status=set_de(nn,xx,ufunc,dfunc,pa); 00572 if (status!=0) { 00573 add_err_ret("Function set() failed in msolve().",gsl_efailed); 00574 } 00575 00576 return solve_set(nn,xx,pa,ufunc); 00577 } 00578 00579 /// Solve \c ufunc using \c xx as an initial guess,returning \c xx. 00580 virtual int msolve(size_t nn, vec_t &xx, param_t &pa, func_t &ufunc) { 00581 00582 int status; 00583 00584 status=set(nn,xx,ufunc,pa); 00585 if (status!=0) { 00586 add_err_ret("Function set() failed in msolve().",gsl_efailed); 00587 } 00588 00589 return solve_set(nn,xx,pa,ufunc); 00590 } 00591 00592 /** \brief Set the function, the parameters, and the initial guess 00593 */ 00594 int set(size_t nn, vec_t &ax, func_t &ufunc, param_t &pa) { 00595 00596 int status; 00597 00598 if (nn!=dim) { 00599 status=allocate(nn); 00600 if (status!=0) { 00601 add_err_ret("Function allocate() failed in msolve().",gsl_efailed); 00602 } 00603 } 00604 00605 fnewp=&ufunc; 00606 00607 // Specify function for automatic Jacobian 00608 ajac->set_function(ufunc); 00609 00610 if (dim==0) { 00611 set_err_ret("Memory has not been allocated. Use allocate().", 00612 o2scl::gsl_ebadlen); 00613 } 00614 00615 params=&pa; 00616 00617 // We can't use gsl_vector_memcpy here,since ax may have a 00618 // different size than dim. 00619 for(size_t i=0;i<dim;i++) gsl_vector_set(x,i,ax[i]); 00620 00621 o2scl_hybrid_state_t *hstate=(o2scl_hybrid_state_t *)state; 00622 00623 gsl_matrix *J=state->J; 00624 gsl_matrix *q=state->q; 00625 gsl_matrix *r=state->r; 00626 gsl_vector *tau=state->tau; 00627 gsl_vector *diag=state->diag; 00628 00629 alloc_vec_t aaf; 00630 ao.allocate(aaf,dim); 00631 for(size_t ikx=0;ikx<dim;ikx++) { 00632 aaf[ikx]=gsl_vector_get(f,ikx); 00633 } 00634 omatrix_view *oJ=(omatrix_view *)J; 00635 00636 status=ufunc(dim,ax,aaf,*params); 00637 if (status!=0) { 00638 set_err_ret("Function returned non-zero in set().",gsl_ebadfunc); 00639 } 00640 00641 if (jac_given) status=(*jac)(dim,ax,aaf,*oJ,*params); 00642 else status=(*ajac)(dim,ax,aaf,*oJ,*params); 00643 if (status!=0) { 00644 add_err_ret("Jacobian failed in set().",gsl_efailed); 00645 } 00646 00647 for(size_t ikx=0;ikx<dim;ikx++) { 00648 gsl_vector_set(f,ikx,aaf[ikx]); 00649 } 00650 ao.free(aaf); 00651 00652 state->iter=1; 00653 state->fnorm=enorm (f); 00654 state->ncfail=0; 00655 state->ncsuc=0; 00656 state->nslow1=0; 00657 state->nslow2=0; 00658 00659 gsl_vector_set_all(dx,0.0); 00660 00661 /* Store column norms in diag */ 00662 00663 if (int_scaling) { 00664 compute_diag(J,diag); 00665 } else { 00666 gsl_vector_set_all(diag,1.0); 00667 } 00668 00669 /* Set delta to factor |D x| or to factor if |D x| is zero */ 00670 00671 state->delta=compute_delta(diag,x); 00672 00673 /* Factorize J into QR decomposition */ 00674 00675 status=gsl_linalg_QR_decomp(J,tau); 00676 if (status!=0) { 00677 set_err_ret("QR decomposition failed in set().",gsl_efailed); 00678 } 00679 00680 status=gsl_linalg_QR_unpack(J,tau,q,r); 00681 if (status!=0) { 00682 set_err_ret("QR unpacking failed in set().",gsl_efailed); 00683 } 00684 00685 set_called=true; 00686 jac_given=false; 00687 00688 return gsl_success; 00689 } 00690 00691 /** \brief Set the function, the Jacobian, the parameters, 00692 and the initial guess 00693 */ 00694 int set_de(size_t nn, vec_t &ax, func_t &ufunc, jfunc_t &dfunc, 00695 param_t &pa) { 00696 00697 fnewp=&ufunc; 00698 jac=&dfunc; 00699 00700 /// Make sure set() uses the right Jacobian 00701 jac_given=true; 00702 00703 int ret=set(nn,ax,ufunc,pa); 00704 00705 /// Reset jac_given since set() will set it back to false 00706 jac_given=true; 00707 set_called=true; 00708 00709 return ret; 00710 } 00711 00712 }; 00713 00714 #ifndef DOXYGENP 00715 } 00716 #endif 00717 00718 #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