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

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).

Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads.