All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
mroot_hybrids.h
Go to the documentation of this file.
1 /*
2  -------------------------------------------------------------------
3 
4  Copyright (C) 2006-2014, Andrew W. Steiner
5 
6  This file is part of O2scl.
7 
8  O2scl is free software; you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version.
12 
13  O2scl is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with O2scl. If not, see <http://www.gnu.org/licenses/>.
20 
21  -------------------------------------------------------------------
22 */
23 /* multiroots/hybrid.c
24  * multiroots/dogleg.c
25  *
26  * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Brian Gough
27  *
28  * This program is free software; you can redistribute it and/or modify
29  * it under the terms of the GNU General Public License as published by
30  * the Free Software Foundation; either version 3 of the License, or (at
31  * your option) any later version.
32  *
33  * This program is distributed in the hope that it will be useful, but
34  * WITHOUT ANY WARRANTY; without even the implied warranty of
35  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
36  * General Public License for more details.
37  *
38  * You should have received a copy of the GNU General Public License
39  * along with this program; if not, write to the Free Software
40  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
41  * 02110-1301, USA.
42  */
43 /** \file mroot_hybrids.h
44  \brief File defining \ref o2scl::mroot_hybrids and
45  specializations
46 */
47 #ifndef O2SCL_MROOT_HYBRIDS_H
48 #define O2SCL_MROOT_HYBRIDS_H
49 
50 #include <string>
51 
52 #include <gsl/gsl_multiroots.h>
53 #include <gsl/gsl_linalg.h>
54 
55 #include <boost/numeric/ublas/vector.hpp>
56 #include <boost/numeric/ublas/vector_proxy.hpp>
57 #include <boost/numeric/ublas/matrix.hpp>
58 #include <boost/numeric/ublas/matrix_proxy.hpp>
59 
60 #include <o2scl/vector.h>
61 #include <o2scl/mroot.h>
62 #include <o2scl/jacobian.h>
63 #include <o2scl/qr.h>
64 // For matrix_out() below
65 #include <o2scl/columnify.h>
66 
67 #ifndef DOXYGEN_NO_O2NS
68 namespace o2scl {
69 #endif
70 
71  /** \brief Base functions for \ref mroot_hybrids
72  */
74 
75  public:
76 
79 
80  /// Compute the actual reduction
81  double compute_actual_reduction(double fnorm0, double fnorm1) {
82  double actred;
83  if (fnorm1<fnorm0) {
84  double u=fnorm1/fnorm0;
85  actred=1-u*u;
86  } else {
87  actred=-1;
88  }
89  return actred;
90  }
91 
92  /// Compute the predicted reduction phi1p=|Q^T f + R dx|
93  double compute_predicted_reduction(double fnorm0, double fnorm1) {
94  double prered;
95  if (fnorm1<fnorm0) {
96  double u=fnorm1/fnorm0;
97  prered=1-u*u;
98  } else {
99  prered=0;
100  }
101  return prered;
102  }
103 
104  /// Compute \f$ R \cdot g \f$ where \c g is the \gradient
105  template<class vec2_t, class mat_t>
106  void compute_Rg(size_t N, const mat_t &r2,
107  const ubvector &gradient2, vec2_t &Rg) {
108 
109  for (size_t i=0;i<N;i++) {
110  double sum=0.0;
111  for (size_t j=i;j<N;j++) {
112  sum+=r2(i,j)*gradient2[j];
113  }
114  Rg[i]=sum;
115  }
116 
117  return;
118  }
119 
120  /// Compute \c w and \c v
121  template<class vec2_t>
122  void compute_wv(size_t n, const ubvector &qtdf2,
123  const ubvector &rdx2, const vec2_t &dx2,
124  const ubvector &diag2, double pnorm,
125  ubvector &w2, ubvector &v2) {
126 
127  for (size_t i=0;i<n;i++) {
128  double diagi=diag2[i];
129  w2[i]=(qtdf2[i]-rdx2[i])/pnorm;
130  v2[i]=diagi*diagi*dx2[i]/pnorm;
131  }
132 
133  return;
134  }
135 
136  /// Compute \f$ R \cdot \mathrm{dx} \f$
137  template<class vec2_t, class mat_t>
138  void compute_rdx(size_t N, const mat_t &r2,
139  const vec2_t &dx2, ubvector &rdx2) {
140 
141  for (size_t i=0;i<N;i++) {
142  double sum=0.0;
143  for (size_t j=i;j<N;j++) {
144  sum+=r2(i,j)*dx2[j];
145  }
146  rdx2[i]=sum;
147  }
148 
149  return;
150  }
151 
152  /** \brief Compute the norm of the vector \f$ \vec{v} \f$ defined
153  by \f$ v_i = d_i ff_i \f$
154  */
155  template<class vec2_t, class vec3_t>
156  double scaled_enorm(size_t n, const vec2_t &d, const vec3_t &ff) {
157  double e2=0.0;
158  for (size_t i=0;i<n;i++) {
159  double u=d[i]*ff[i];
160  e2+=u*u;
161  }
162  return sqrt(e2);
163  }
164 
165  /// Compute delta
166  template<class vec2_t>
167  double compute_delta(size_t n, ubvector &diag2, vec2_t &x2) {
168  double Dx=scaled_enorm(n,diag2,x2);
169  double factor=100.0;
170  return (Dx > 0) ? factor*Dx : factor;
171  }
172 
173  /** \brief Compute the Euclidean norm of \c ff
174 
175  \future Replace this with \c dnrm2 from \ref cblas_base.h
176  */
177  template<class vec2_t> double enorm(size_t N, const vec2_t &ff) {
178  double e2=0.0;
179  for (size_t i=0;i<N;i++) {
180  double fi=ff[i];
181  e2+=fi*fi;
182  }
183  return sqrt(e2);
184  }
185 
186  /// Compute the Euclidean norm of the sum of \c a and \c b
187  double enorm_sum(size_t n, const ubvector &a,
188  const ubvector &b) {
189  double e2=0.0;
190  for (size_t i=0;i<n;i++) {
191  double u=a[i]+b[i];
192  e2+=u*u;
193  }
194  return sqrt(e2);
195  }
196 
197  /** \brief Compute a trial step and put the result in \c xx_trial
198 
199  \future Replace with daxpy.
200  */
201  template<class vec2_t>
202  int compute_trial_step(size_t N, vec2_t &xl, vec2_t &dxl,
203  vec2_t &xx_trial) {
204  for (size_t i=0;i<N;i++) {
205  xx_trial[i]=xl[i]+dxl[i];
206  }
207  return 0;
208  }
209 
210  /** \brief Compute the change in the function value
211 
212  \future Replace with daxpy?
213  */
214  template<class vec2_t>
215  int compute_df(size_t n, const vec2_t &ff_trial,
216  const vec2_t &fl, ubvector &dfl) {
217  for (size_t i=0;i<n;i++) {
218  dfl[i]=ff_trial[i]-fl[i];
219  }
220 
221  return 0;
222  }
223 
224  /// Compute \c diag, the norm of the columns of the Jacobian
225  template<class mat2_t>
226  void compute_diag(size_t n, const mat2_t &J2, ubvector &diag2) {
227  for (size_t j=0;j<n;j++) {
228  double sum=0;
229  for (size_t i=0;i<n;i++) {
230  const double Jij=J2(i,j);
231  sum+=Jij*Jij;
232  }
233  if (sum == 0) {
234  sum=1.0;
235  }
236  diag2[j]=sqrt(sum);
237  }
238 
239  return;
240  }
241 
242  /** \brief Compute \f$ Q^{T} f \f$
243 
244  \future This is just right-multiplication, so we could
245  use the \o2 cblas routines instead.
246  */
247  template<class vec2_t, class vec3_t, class vec4_t>
248  void compute_qtf(size_t N, const vec2_t &q2,
249  const vec3_t &f2, vec4_t &qtf2) {
250  for (size_t j=0;j<N;j++) {
251  double sum=0.0;
252  for (size_t i=0;i<N;i++) {
253  sum+=q2(i,j)*f2[i];
254  }
255  qtf2[j]=sum;
256  }
257 
258  return;
259  }
260 
261  /// Update \c diag
262  template<class mat2_t>
263  void update_diag(size_t n, const mat2_t &J2, ubvector &diag2) {
264  for (size_t j=0;j<n;j++) {
265  double cnorm, diagj, sum=0.0;
266  for (size_t i=0;i<n;i++) {
267  double Jij=J2(i,j);
268  sum+=Jij*Jij;
269  }
270  if (sum == 0) {
271  sum=1.0;
272  }
273  cnorm=sqrt(sum);
274  diagj=diag2[j];
275  if (cnorm > diagj) {
276  diag2[j]=cnorm;
277  }
278  }
279 
280  return;
281  }
282 
283  /** \brief Form appropriate convex combination of the Gauss-Newton
284  direction and the scaled gradient direction.
285 
286  Using the Gauss-Newton direction given in \c newton (a vector of
287  size N), and the gradient direction in \c gradient (also a
288  vector of size N), this computes
289  \f[
290  \mathrm{pp}=\alpha \mathrm{newton}+\beta \mathrm{gradient}
291  \f]
292  */
293  template<class vec2_t>
294  void scaled_addition(size_t N, double alpha, ubvector &newton2,
295  double beta, ubvector &gradient2, vec2_t &pp) {
296  for (size_t i=0;i<N;i++) {
297  pp[i]=alpha*newton2[i]+beta*gradient2[i];
298  }
299 
300  return;
301  }
302 
303  /// Compute the Gauss-Newton direction
304  template<class mat_t>
305  int newton_direction(const size_t N, const mat_t &r2,
306  const ubvector &qtf2, ubvector &p) {
307  size_t i;
308  int status=0;
309 
310  vector_copy(N,qtf2,p);
311  o2scl_cblas::dtrsv(o2scl_cblas::o2cblas_RowMajor,
312  o2scl_cblas::o2cblas_Upper,
313  o2scl_cblas::o2cblas_NoTrans,
314  o2scl_cblas::o2cblas_NonUnit,N,N,r2,p);
315 
316  for (i=0;i<N;i++) {
317  p[i]*=-1.0;
318  }
319 
320  return status;
321  }
322 
323  /// Compute the gradient direction
324  template<class mat_t>
325  void gradient_direction(const size_t M, const size_t N,
326  const mat_t &r2, const ubvector &qtf2,
327  const ubvector &diag2, ubvector &g) {
328  for (size_t j=0;j<M;j++) {
329  double sum=0.0;
330  for (size_t i=0;i<N;i++) {
331  sum+=r2(i,j)*qtf2[i];
332  }
333  g[j]=-sum/diag2[j];
334  }
335 
336  return;
337  }
338 
339  /// Compute the point at which the gradient is minimized
340  void minimum_step(const size_t N, double gnorm,
341  const ubvector &diag2, ubvector &g) {
342  for (size_t i=0;i<N;i++) {
343  g[i]=(g[i]/gnorm)/diag2[i];
344  }
345 
346  return;
347  }
348 
349  /** \brief Take a dogleg step
350 
351  Given the QR decomposition of an \c n by \c n matrix "A",
352  a diagonal matrix \c diag, a vector \c b, and a positive
353  number \c delta, this function determines the convex
354  combination \c x of the Gauss-Newton and scaled gradient
355  directions which minimizes \f$ A x-b \f$ in the least
356  squares sense, subject to the restriction that the
357  Euclidean norm of \f$ d x \f$ is smaller than the
358  value of \c delta.
359  */
360  template<class vec2_t, class mat_t>
361  int dogleg(size_t n, const mat_t &r2, const ubvector &qtf2,
362  const ubvector &diag2, double delta2,
363  ubvector &newton2, ubvector &gradient2,
364  vec2_t &p) {
365 
366  double qnorm, gnorm, sgnorm, bnorm, temp;
367 
368  // Compute the Gauss-Newton direction
369  newton_direction(n,r2,qtf2,newton2);
370  qnorm=scaled_enorm(n,diag2,newton2);
371 
372  // Test whether the gauss-newton direction is acceptable
373  if (qnorm <= delta2) {
374  for(size_t i=0;i<n;i++) p[i]=newton2[i];
375  return success;
376  }
377 
378  // Compute the gradient direction
379  gradient_direction(n,n,r2,qtf2,diag2,gradient2);
380  gnorm=enorm(n,gradient2);
381 
382  // The special case that the scaled gradient is zero
383  if (gnorm == 0) {
384  double alpha=delta2/qnorm;
385  double beta=0;
386  scaled_addition(n,alpha,newton2,beta,gradient2,p);
387  return success;
388  }
389 
390  // -----------------------------------------------------
391  // The next set of statements calculates the point along the
392  // scaled gradient at which the quadratic is minimized.
393 
394  minimum_step(n,gnorm,diag2,gradient2);
395 
396  // Use p as temporary space to compute Rg
397  compute_Rg(n,r2,gradient2,p);
398 
399  temp=enorm(n,p);
400  sgnorm=(gnorm/temp)/temp;
401 
402  // -----------------------------------------------------
403  // Test whether the scaled gradient direction is acceptable
404 
405  if (sgnorm>delta2) {
406 
407  double alpha=0;
408  double beta=delta2;
409 
410  scaled_addition(n,alpha,newton2,beta,gradient2,p);
411 
412  return success;
413  }
414 
415  // If not, calculate the point along the dogleg step
416  // at which the quadratic is minimized
417  bnorm=enorm(n,qtf2);
418 
419  double bg=bnorm/gnorm;
420  double bq=bnorm/qnorm;
421  double dq=delta2/qnorm;
422  double dq2=dq*dq;
423  double sd=sgnorm/delta2;
424  double sd2=sd*sd;
425 
426  double t1=bg*bq*sd;
427  double u=t1-dq;
428  double t2=t1-dq*sd2+sqrt(u*u+(1-dq2)*(1-sd2));
429 
430  double alpha=dq*(1-sd2)/t2;
431  double beta=(1-alpha)*sgnorm;
432 
433  // Form the appropriate convex combination of the gauss-newton
434  // direction and the scaled gradient direction.
435 
436  scaled_addition(n,alpha,newton2,beta,gradient2,p);
437 
438  return success;
439  }
440 
441  };
442 
443  /** \brief Multidimensional root-finding algorithm using
444  Powell's Hybrid method (GSL)
445 
446  This is a recasted version of the GSL routines which use a
447  modified version of Powell's Hybrid method as implemented in the
448  HYBRJ algorithm in MINPACK (\ref Garbow80). Both the scaled and
449  unscaled options are available by setting \ref int_scaling (the
450  scaled version is the default). If derivatives are not provided,
451  they will be computed automatically. This class provides the
452  GSL-like interface using allocate(), set() (or set_de() in case
453  where derivatives are available), iterate(), and the
454  higher-level interfaces, msolve() and msolve_de(), which perform
455  the solution and the memory allocation automatically. Some
456  additional checking is performed in case the user calls the
457  functions out of order (i.e. set() without allocate()).
458 
459  The functions msolve() and msolve_de() use the condition \f$
460  \sum_i |f_i|<\f$ \ref mroot::tol_rel to determine if the solver has
461  succeeded.
462 
463  See the \ref multisolve_subsect section of the User's guide for
464  general information about \o2 solvers.
465  There is an example for the usage of the multidimensional solver
466  classes given in <tt>examples/ex_mroot.cpp</tt>, see the \ref
467  ex_mroot_sect .
468 
469  The original GSL algorithm has been modified to shrink the
470  stepsize if a proposed step causes the function to return a
471  non-zero value. This allows the routine to automatically try to
472  avoid regions where the function is not defined. The algorithm
473  is also modified to check that it is not sending non-finite
474  values to the user-specified function. To return to the default
475  GSL behavior, set \ref shrink_step and \ref extra_finite_check
476  to false.
477 
478  The default method for numerically computing the Jacobian is
479  from \ref jacobian_gsl. This default is identical to the GSL
480  approach, except that the default value of \ref
481  jacobian_gsl::epsmin is non-zero. See \ref jacobian_gsl
482  for more details.
483 
484  By default convergence failures result in calling the exception
485  handler, but this can be turned off by setting \ref
486  mroot::err_nonconv to false. If \ref mroot::err_nonconv is
487  false, then the functions \ref iterate(), \ref msolve() and \ref
488  msolve_de() will return a non-zero value if convergence fails.
489  Note that the default Jacobian object, \ref def_jac also has a
490  data member \ref jacobian_gsl::err_nonconv which separately
491  handles the case where the one row of the Jacobian is all zero.
492 
493  \future Is all the setting of vectors and matrices to zero really
494  necessary? Do they need to be executed even if memory hasn't
495  been recently allocated?
496 
497  \future Convert more ubvectors to vec_t.
498 
499  \future Some more of the element-wise vector manipulation could be
500  converted to BLAS routines.
501 
502  \future It's kind of strange that set() sets jac_given to false
503  and set_de() has to reset it to true. Can this be simplified?
504 
505  \future Many of these minpack functions could be put in their
506  own "minpack_tools" class, or possibly moved to be
507  linear algebra routines instead.
508 
509  \future There are still some numbers in here which the user
510  could have control over, for example, the <tt>nslow2</tt>
511  threshold which indicates failure.
512  */
513  template<
514  class func_t=mm_funct11,
517  class jfunc_t=jac_funct11 > class mroot_hybrids :
518  public mroot<func_t,vec_t,jfunc_t>, mroot_hybrids_base {
519 
520 #ifndef DOXYGEN_INTERNAL
521 
522  protected:
523 
524  /// Number of iterations
525  int iter;
526  /// Compute the number of failures
527  size_t ncfail;
528  /// Compute the number of successes
529  size_t ncsuc;
530  /// The number of times the actual reduction is less than 0.001
531  size_t nslow1;
532  /// The number of times the actual reduction is less than 0.1
533  size_t nslow2;
534  /// The norm of the current function value
535  double fnorm;
536  /// The limit of the Nuclidean norm
537  double delta;
538  /// Jacobian
539  mat_t J;
540  /// Q matrix from QR decomposition
541  mat_t q;
542  /// R matrix from QR decomposition
543  mat_t r;
544  /// The diagonal elements
546  /// The value of \f$ Q^T f \f$
548  /// The Newton direction
550  /// The gradient direction
552  /// The change in the function value
554  /// The value of \f$ Q^T \cdot \mathrm{df} \f$
556  /// The value of \f$ R \cdot \mathrm{dx} \f$
558  /// The value of \f$ w=(Q^T df-R dx)/|dx| \f$
560  /// The value of \f$ v=D^2 dx/|dx| \f$
562 
563  /// The user-specified Jacobian
564  jfunc_t *jac;
565 
566  /// The automatic Jacobian
568 
569  /// The value of the derivative
570  vec_t dx;
571 
572  /// Trial root
573  vec_t x_trial;
574 
575  /// Trial function value
576  vec_t f_trial;
577 
578  /// The number of equations and unknowns
579  size_t dim;
580 
581  /// True if the jacobian has been given
582  bool jac_given;
583 
584  /// The user-specified function
585  func_t *fnewp;
586 
587  /// True if "set" has been called
589 
590  /// Finish the solution after set() or set_de() has been called
591  int solve_set(size_t nn, vec_t &xx, func_t &ufunc) {
592 
593  int status;
594  iter=0;
595 
596  do {
597  iter++;
598 
599  if (iterate()!=0) {
600  O2SCL_CONV2_RET("Function iterate() failed in mroot_hybrids::",
601  "solve_set().",exc_efailed,this->err_nonconv);
602  }
603 
604  // ------------------------------------------------------
605  // The equivalent of the statement:
606  //
607  // status=gsl_multiroot_test_residual(f,this->tol_rel);
608 
609  double resid=0.0;
610  for(size_t i=0;i<nn;i++) {
611  resid+=fabs(f[i]);
612  }
613  if (resid<this->tol_rel) status=success;
614  else status=gsl_continue;
615 
616  // ------------------------------------------------------
617 
618  if (this->verbose>0) {
619  this->print_iter(nn,x,f,iter,resid,this->tol_rel,
620  "mroot_hybrids");
621  }
622 
623  } while (status==gsl_continue && iter<this->ntrial);
624 
625  for(size_t i=0;i<nn;i++) {
626  xx[i]=x[i];
627  }
628 
629  this->last_ntrial=iter;
630 
631  if (((int)iter)>=this->ntrial) {
632  O2SCL_CONV2_RET("Function mroot_hybrids::msolve() ",
633  "exceeded max. number of iterations.",
634  exc_emaxiter,this->err_nonconv);
635  }
636 
637  return success;
638  }
639 
640 #endif
641 
642  public:
643 
644  mroot_hybrids() {
645  shrink_step=true;
646  dim=0;
647  ajac=&def_jac;
648  def_jac.epsrel=sqrt(std::numeric_limits<double>::epsilon());
649  int_scaling=true;
650  jac_given=false;
651  set_called=false;
652  extra_finite_check=true;
653  }
654 
655  virtual ~mroot_hybrids() {
656  }
657 
658  /** \brief If true, iterate() will shrink the step-size automatically if
659  the function returns a non-zero value (default true)
660 
661  The original GSL behavior can be obtained by setting
662  this to \c false.
663  */
665 
666  /** \brief If true, double check that the input function values are
667  finite (default true)
668  */
670 
671  /// If true, use the internal scaling method (default true)
673 
674  /// Default automatic Jacobian object
676 
677  /// Set the automatic Jacobian object
679  ajac=&j;
680  return 0;
681  }
682 
683  /** \brief The value of the function at the present iteration
684 
685  \comment
686  We need this to be public so that the user can see if
687  iterate() is working
688  \endcomment
689  */
690  vec_t f;
691 
692  /// The present solution
693  vec_t x;
694 
695  /** \brief Perform an iteration
696 
697  At the end of the iteration, the current value of the solution
698  is stored in \ref x.
699 
700  \todo This now just returns zero during convergence errors,
701  because it caused problems with tov_eos_ts, but this should be
702  fixed so that it returns a non-zero value, and functions which
703  call iterate() need to handle the return value appropriately.
704  */
705  int iterate() {
706 
707  if (!set_called) {
708  O2SCL_ERR2("Function set() not called or most recent call ",
709  "failed in mroot_hybrids::iterate().",
710  exc_efailed);
711  }
712 
713  double prered, actred;
714  double pnorm, fnorm1, fnorm1p;
715  double ratio;
716  double p1=0.1, p5=0.5, p001=0.001, p0001=0.0001;
717 
718  /* Compute qtf=Q^T f */
719 
720  compute_qtf(dim,q,f,qtf);
721 
722  /* Compute dogleg step */
723 
725 
726  /* Take a trial step */
727 
729  pnorm=scaled_enorm(dim,diag,dx);
730 
731  if (iter==1) {
732  if (pnorm<delta) {
733  delta=pnorm;
734  }
735  }
736 
737  if (extra_finite_check) {
738  // AWS, 11/14/13: This is not strictly necessary, because if
739  // x_trial is not finite the iteration will fail to converge
740  // anyway, but it is disconcerting to the user to have
741  // non-finite values sent to the user-specified function for
742  // no apparent reason. In reality what appears to be happening
743  // is that pnorm was previously zero (because of a vanishingly
744  // small step size), and the call to compute_wv() below then
745  // leads to non-finite values. On the other hand, checking all
746  // the vector values is time consuming, so I perform this
747  // check only if extra_finite_check is true.
748  for(size_t ik=0;ik<dim;ik++) {
749  if (!o2scl::is_finite(x_trial[ik])) {
750  O2SCL_CONV2_RET("Iteration lead to non-finite values in ",
751  "mroot_hybrids::iterate().",exc_efailed,
752  this->err_nonconv);
753  }
754  }
755  }
756 
757  /* Evaluate function at x+p */
758 
759  int status;
760 
761  if (shrink_step==false) {
762 
763  // Evaluate the function at the new point, exit if it fails
764 
765  status=(*fnewp)(dim,x_trial,f_trial);
766 
767  if (status != success) {
768  std::string str="Function returned non-zero value ("+itos(status)+
769  ") in mroot_hybrids::iterate().";
771  }
772 
773  } else {
774 
775  // Evaluate the function at the new point, try to recover
776  // if it fails
777 
778  status=(*fnewp)(dim,x_trial,f_trial);
779 
780  int bit=0;
781  while(status!=0 && bit<20) {
782  for(size_t ib=0;ib<dim;ib++) {
783  x_trial[ib]=(x_trial[ib]+x[ib])/2.0;
784  }
785  status=(*fnewp)(dim,x_trial,f_trial);
786  bit++;
787  }
788 
789  // Exit if we weren't able to find a new good point
790 
791  if (status != success) {
792  std::string str="No suitable step found, function returned ("+
793  itos(status)+") in mroot_hybrids::iterate().";
795  }
796 
797  }
798 
799  /* Set df=f_trial-f */
800 
802 
803  /* Compute the scaled actual reduction */
804 
805  fnorm1=enorm(dim,f_trial);
806  actred=compute_actual_reduction(fnorm,fnorm1);
807 
808  /* Compute rdx=R dx */
809 
811 
812  /* Compute the scaled predicted reduction phi1p=|Q^T f+R dx| */
813 
814  fnorm1p=enorm_sum(dim,qtf,rdx);
815  prered=compute_predicted_reduction(fnorm,fnorm1p);
816 
817  /* Compute the ratio of the actual to predicted reduction */
818 
819  if (prered > 0) {
820  ratio=actred/prered;
821  } else {
822  ratio=0;
823  }
824 
825  /* Update the step bound */
826 
827  if (ratio<p1) {
828  ncsuc=0;
829  ncfail++;
830  delta*=p5;
831  } else {
832  ncfail=0;
833  ncsuc++;
834 
835  if (ratio >= p5 || ncsuc > 1) {
836  delta=GSL_MAX(delta,pnorm/p5);
837  }
838  if (fabs (ratio-1) <= p1) {
839  delta=pnorm/p5;
840  }
841  }
842 
843  /* Test for successful iteration */
844 
845  if (ratio >= p0001) {
846  for(size_t i=0;i<dim;i++) {
847  x[i]=x_trial[i];
848  f[i]=f_trial[i];
849  }
850  fnorm=fnorm1;
851  iter++;
852  }
853 
854  /* Determine the progress of the iteration */
855 
856  nslow1++;
857  if (actred >= p001) {
858  nslow1=0;
859  }
860 
861  if (actred >= p1) {
862  nslow2=0;
863  }
864 
865  if (ncfail==2) {
866 
867  int jac_ret;
868 
869  if (jac_given) jac_ret=(*jac)(dim,x,dim,f,J);
870  else jac_ret=(*ajac)(dim,x,dim,f,J);
871 
872  if (jac_ret!=0) {
873  std::string str="Jacobian failed and returned ("+
874  itos(jac_ret)+") in mroot_hybrids::iterate().";
875  O2SCL_CONV_RET(str.c_str(),exc_efailed,this->err_nonconv);
876  }
877 
878  nslow2++;
879 
880  if (iter==1) {
881  if (int_scaling) {
882  compute_diag(dim,J,diag);
883  }
884  delta=compute_delta(dim,diag,x);
885  } else {
886  if (int_scaling) {
887  update_diag(dim,J,diag);
888  }
889  }
890 
891  o2scl_linalg::QR_decomp_unpack(dim,dim,this->J,this->q,this->r);
892 
893  return success;
894  }
895 
896  /* Compute qtdf=Q^T df, w=(Q^T df-R dx)/|dx|, v=D^2 dx/|dx| */
897 
899  compute_wv(dim,qtdf,rdx,dx,diag,pnorm,w,v);
900 
901  /* Rank-1 update of the jacobian Q'R'=Q(R+w v^T) */
902 
904 
905  /* No progress as measured by jacobian evaluations */
906 
907  if (nslow2==5) {
908  O2SCL_CONV2_RET("No progress in Jacobian in mroot_hybrids::",
909  "iterate().",exc_enoprogj,this->err_nonconv);
910  }
911 
912  /* No progress as measured by function evaluations */
913 
914  if (nslow1==10) {
915  O2SCL_CONV2_RET("No progress in function in mroot_hybrids::",
916  "iterate().",exc_enoprog,this->err_nonconv);
917  }
918 
919  return success;
920  }
921 
922  /// Allocate the memory
923  void allocate(size_t n) {
924 
925  q.resize(n,n);
926  r.resize(n,n);
927  diag.resize(n);
928  qtf.resize(n);
929  newton.resize(n);
930  gradient.resize(n);
931  df.resize(n);
932  qtdf.resize(n);
933  rdx.resize(n);
934  w.resize(n);
935  v.resize(n);
936 
937  for(size_t ii=0;ii<n;ii++) {
938  for(size_t jj=0;jj<n;jj++) {
939  q(ii,jj)=0.0;
940  }
941  }
942  for(size_t ii=0;ii<n;ii++) {
943  for(size_t jj=0;jj<n;jj++) {
944  r(ii,jj)=0.0;
945  }
946  }
947 
948  for(size_t ii=0;ii<diag.size();ii++) diag[ii]=0.0;
949  for(size_t ii=0;ii<qtf.size();ii++) qtf[ii]=0.0;
950  for(size_t ii=0;ii<newton.size();ii++) newton[ii]=0.0;
951  for(size_t ii=0;ii<gradient.size();ii++) gradient[ii]=0.0;
952  for(size_t ii=0;ii<df.size();ii++) df[ii]=0.0;
953  for(size_t ii=0;ii<qtdf.size();ii++) qtdf[ii]=0.0;
954  for(size_t ii=0;ii<rdx.size();ii++) rdx[ii]=0.0;
955  for(size_t ii=0;ii<w.size();ii++) w[ii]=0.0;
956  for(size_t ii=0;ii<v.size();ii++) v[ii]=0.0;
957 
958  J.resize(n,n);
959 
960  x.resize(n);
961  dx.resize(n);
962  f.resize(n);
963  x_trial.resize(n);
964  f_trial.resize(n);
965 
966  for(size_t i=0;i<n;i++) {
967  x[i]=0.0;
968  dx[i]=0.0;
969  f[i]=0.0;
970  }
971 
972  dim=n;
973 
974  return;
975  }
976 
977  /// Return the type,\c "mroot_hybrids".
978  virtual const char *type() { return "mroot_hybrids"; }
979 
980  /** \brief Solve \c func with derivatives \c dfunc using \c x as
981  an initial guess, returning \c x.
982 
983  */
984  virtual int msolve_de(size_t nn, vec_t &xx, func_t &ufunc,
985  jfunc_t &dfunc) {
986 
987  int ret=set_de(nn,xx,ufunc,dfunc);
988  if (ret!=success) {
989  O2SCL_CONV2_RET("Function set_de() failed in mroot_hybrids::",
990  "msolve_de().",exc_efailed,this->err_nonconv);
991  }
992 
993  return solve_set(nn,xx,ufunc);
994  }
995 
996  /// Solve \c ufunc using \c xx as an initial guess, returning \c xx.
997  virtual int msolve(size_t nn, vec_t &xx, func_t &ufunc) {
998 
999  int ret=set(nn,xx,ufunc);
1000  if (ret!=success) {
1001  O2SCL_CONV2_RET("Function set() failed in mroot_hybrids::",
1002  "msolve().",exc_efailed,this->err_nonconv);
1003  }
1004  int status=solve_set(nn,xx,ufunc);
1005  for(size_t i=0;i<nn;i++) xx[i]=x[i];
1006 
1007  return status;
1008  }
1009 
1010  /** \brief Set the function, the parameters, and the initial guess
1011  */
1012  int set(size_t nn, vec_t &ax, func_t &ufunc) {
1013 
1014  int status;
1015 
1016  if (nn!=dim) {
1017  allocate(nn);
1018  }
1019 
1020  fnewp=&ufunc;
1021 
1022  // Specify function for automatic Jacobian
1023  ajac->set_function(ufunc);
1024 
1025  if (dim==0) {
1026  O2SCL_ERR2("No memory allocated in ",
1027  "mroot_hybrids::set().",o2scl::exc_ebadlen);
1028  }
1029 
1030  // Copy the user-specified solution
1031  for(size_t i=0;i<dim;i++) x[i]=ax[i];
1032 
1033  status=ufunc(dim,ax,f);
1034  if (status!=0) {
1035  O2SCL_ERR2("Function returned non-zero in ",
1036  "mroot_hybrids::set().",exc_ebadfunc);
1037  }
1038 
1039  if (jac_given) status=(*jac)(dim,ax,dim,f,J);
1040  else status=(*ajac)(dim,ax,dim,f,J);
1041 
1042  if (status!=0) {
1043  O2SCL_CONV2_RET("Jacobian failed in ",
1044  "mroot_hybrids::set().",exc_efailed,this->err_nonconv);
1045  }
1046 
1047  iter=1;
1048  fnorm=enorm(dim,f);
1049  ncfail=0;
1050  ncsuc=0;
1051  nslow1=0;
1052  nslow2=0;
1053 
1054  for(size_t i=0;i<dim;i++) dx[i]=0.0;
1055 
1056  /* Store column norms in diag */
1057 
1058  if (int_scaling) {
1059  compute_diag(dim,J,diag);
1060  } else {
1061  for(size_t ii=0;ii<diag.size();ii++) diag[ii]=0.0;
1062  }
1063 
1064  /* Set delta to factor |D x| or to factor if |D x| is zero */
1065 
1066  delta=compute_delta(dim,diag,x);
1067 
1068  /* Factorize J into QR decomposition */
1069  o2scl_linalg::QR_decomp_unpack(dim,dim,this->J,this->q,this->r);
1070  set_called=true;
1071  jac_given=false;
1072 
1073  return 0;
1074  }
1075 
1076  /** \brief Set the function, the Jacobian, the parameters,
1077  and the initial guess
1078  */
1079  int set_de(size_t nn, vec_t &ax, func_t &ufunc, jfunc_t &dfunc) {
1080 
1081  fnewp=&ufunc;
1082  jac=&dfunc;
1083 
1084  // Make sure set() uses the right Jacobian
1085  jac_given=true;
1086 
1087  int ret=set(nn,ax,ufunc);
1088 
1089  // Reset jac_given since set() will set it back to false
1090  jac_given=true;
1091  set_called=true;
1092 
1093  return ret;
1094  }
1095 
1096 #ifndef DOXYGEN_INTERNAL
1097 
1098  private:
1099 
1104 
1105 #endif
1106 
1107  };
1108 
1109 #ifndef DOXYGEN_NO_O2NS
1110 }
1111 #endif
1112 
1113 #if defined (O2SCL_COND_FLAG) || defined (DOXYGEN)
1114 #if defined (O2SCL_ARMA) || defined (DOXYGEN)
1115 #include <armadillo>
1116 namespace o2scl {
1117  /** \brief A version of \ref mroot_hybrids which uses
1118  Armadillo for the QR decomposition
1119 
1120  \note This class template is only defined if Armadillo
1121  was enabled when \o2 was installed.
1122  */
1123  template<class func_t, class vec_t, class mat_t, class jfunc_t>
1125  public mroot_hybrids<func_t,vec_t,mat_t,jfunc_t> {
1126 
1127  virtual void qr_decomp_unpack() {
1128  arma::qr_econ(this->q,this->r,this->J);
1129  return;
1130  }
1131 
1132  };
1133 }
1134 #endif
1135 #if defined (O2SCL_HAVE_EIGEN) || defined (DOXYGEN)
1136 #include <Eigen/Dense>
1137 namespace o2scl {
1138  /** \brief A version of \ref mroot_hybrids
1139  which uses Eigen for the QR decomposition
1140 
1141  \note This class template is only defined if Eigen
1142  was enabled when \o2 was installed.
1143  */
1144  template<class func_t, class vec_t, class mat_t, class jfunc_t>
1146  public mroot_hybrids<func_t,vec_t,mat_t,jfunc_t> {
1147 
1148  virtual void qr_decomp_unpack() {
1149  Eigen::HouseholderQR<Eigen::MatrixXd> hqr(this->J);
1150  this->q=hqr.householderQ();
1151  this->r=hqr.matrixQR().triangularView<Eigen::Upper>();
1152  return;
1153  }
1154 
1155  };
1156 }
1157 #endif
1158 #else
1159 #include <o2scl/mroot_special.h>
1160 #endif
1161 
1162 #endif
Class for automatically computing gradients [abstract base].
Definition: mmin.h:53
int iter
Number of iterations.
double enorm(size_t N, const vec2_t &ff)
Compute the Euclidean norm of ff.
bool shrink_step
If true, iterate() will shrink the step-size automatically if the function returns a non-zero value (...
double compute_delta(size_t n, ubvector &diag2, vec2_t &x2)
Compute delta.
vec_t dx
The value of the derivative.
void dtrsv(const enum o2cblas_order order, const enum o2cblas_uplo Uplo, const enum o2cblas_transpose TransA, const enum o2cblas_diag Diag, const size_t M, const size_t N, const mat_t &A, vec_t &X)
Compute .
Definition: cblas_base.h:308
ubvector qtdf
The value of .
int verbose
Output control (default 0)
Definition: mroot.h:88
size_t ncsuc
Compute the number of successes.
vec_t x_trial
Trial root.
ubvector v
The value of .
ubvector w
The value of .
#define O2SCL_CONV_RET(d, n, b)
Set a "convergence" error and return the error value.
Definition: err_hnd.h:292
jacobian_gsl< func_t, vec_t, mat_t > def_jac
Default automatic Jacobian object.
size_t ncfail
Compute the number of failures.
void allocate(size_t n)
Allocate the memory.
bool is_finite(double x)
Return false if x is infinite or not a number.
double fnorm
The norm of the current function value.
mat_t J
Jacobian.
void compute_Rg(size_t N, const mat_t &r2, const ubvector &gradient2, vec2_t &Rg)
Compute where g is the gradient gradient .
std::function< int(size_t, const boost::numeric::ublas::vector< double > &, boost::numeric::ublas::vector< double > &) > mm_funct11
Array of multi-dimensional functions typedef.
Definition: mm_funct.h:43
int compute_df(size_t n, const vec2_t &ff_trial, const vec2_t &fl, ubvector &dfl)
Compute the change in the function value.
ubvector diag
The diagonal elements.
double compute_predicted_reduction(double fnorm0, double fnorm1)
Compute the predicted reduction phi1p=|Q^T f + R dx|.
Definition: mroot_hybrids.h:93
exceeded max number of iterations
Definition: err_hnd.h:73
void compute_wv(size_t n, const ubvector &qtdf2, const ubvector &rdx2, const vec2_t &dx2, const ubvector &diag2, double pnorm, ubvector &w2, ubvector &v2)
Compute w and v.
void update_diag(size_t n, const mat2_t &J2, ubvector &diag2)
Update diag.
#define O2SCL_CONV2_RET(d, d2, n, b)
Set an error and return the error value, two-string version.
Definition: err_hnd.h:298
Multidimensional root-finding [abstract base].
Definition: mroot.h:66
A version of mroot_hybrids which uses Eigen for the QR decomposition.
int set_de(size_t nn, vec_t &ax, func_t &ufunc, jfunc_t &dfunc)
Set the function, the Jacobian, the parameters, and the initial guess.
jacobian< func_t, vec_t, mat_t > * ajac
The automatic Jacobian.
iteration has not converged
Definition: err_hnd.h:51
size_t nslow1
The number of times the actual reduction is less than 0.001.
double enorm_sum(size_t n, const ubvector &a, const ubvector &b)
Compute the Euclidean norm of the sum of a and b.
vec_t f_trial
Trial function value.
int newton_direction(const size_t N, const mat_t &r2, const ubvector &qtf2, ubvector &p)
Compute the Gauss-Newton direction.
int solve_set(size_t nn, vec_t &xx, func_t &ufunc)
Finish the solution after set() or set_de() has been called.
void compute_rdx(size_t N, const mat_t &r2, const vec2_t &dx2, ubvector &rdx2)
Compute .
generic failure
Definition: err_hnd.h:61
int print_iter(size_t n, const vec2_t &x, const vec3_t &y, int iter, double value=0.0, double limit=0.0, std::string comment="")
Print out iteration information.
Definition: mroot.h:133
size_t nslow2
The number of times the actual reduction is less than 0.1.
int ntrial
Maximum number of iterations (default 100)
Definition: mroot.h:91
void QR_decomp_unpack(const size_t M, const size_t N, mat_t &A, mat2_t &Q, mat3_t &R)
Compute the unpacked QR decomposition of matrix A.
Definition: qr_base.h:231
double delta
The limit of the Nuclidean norm.
vec_t f
The value of the function at the present iteration.
void minimum_step(const size_t N, double gnorm, const ubvector &diag2, ubvector &g)
Compute the point at which the gradient is minimized.
mat_t q
Q matrix from QR decomposition.
virtual int set_jacobian(jacobian< func_t, vec_t, mat_t > &j)
Set the automatic Jacobian object.
bool int_scaling
If true, use the internal scaling method (default true)
Multidimensional root-finding algorithm using Powell's Hybrid method (GSL)
void compute_diag(size_t n, const mat2_t &J2, ubvector &diag2)
Compute diag, the norm of the columns of the Jacobian.
iteration is not making progress toward solution
Definition: err_hnd.h:105
matrix, vector lengths are not conformant
Definition: err_hnd.h:89
vec_t x
The present solution.
bool err_nonconv
If true, call the error handler if msolve() or msolve_de() does not converge (default true) ...
Definition: mroot.h:99
A version of mroot_hybrids which uses Armadillo for the QR decomposition.
#define O2SCL_ERR2(d, d2, n)
Set an error, two-string version.
Definition: err_hnd.h:281
void scaled_addition(size_t N, double alpha, ubvector &newton2, double beta, ubvector &gradient2, vec2_t &pp)
Form appropriate convex combination of the Gauss-Newton direction and the scaled gradient direction...
int last_ntrial
The number of iterations for in the most recent minimization.
Definition: mroot.h:94
int iterate()
Perform an iteration.
jacobian jacobian evaluations are not improving the solution
Definition: err_hnd.h:107
std::function< int(size_t, boost::numeric::ublas::vector< double > &, size_t, boost::numeric::ublas::vector< double > &, boost::numeric::ublas::matrix< double > &) > jac_funct11
Jacobian function (not necessarily square)
Definition: jacobian.h:44
int set(size_t nn, vec_t &ax, func_t &ufunc)
Set the function, the parameters, and the initial guess.
Base functions for mroot_hybrids.
Definition: mroot_hybrids.h:73
jfunc_t * jac
The user-specified Jacobian.
func_t * fnewp
The user-specified function.
ubvector df
The change in the function value.
ubvector gradient
The gradient direction.
ubvector qtf
The value of .
bool extra_finite_check
If true, double check that the input function values are finite (default true)
void QR_update(size_t M, size_t N, mat1_t &Q, mat2_t &R, vec1_t &w, vec2_t &v)
Update a QR factorisation for A= Q R, A' = A + u v^T,.
Definition: qr_base.h:165
void gradient_direction(const size_t M, const size_t N, const mat_t &r2, const ubvector &qtf2, const ubvector &diag2, ubvector &g)
Compute the gradient direction.
ubvector newton
The Newton direction.
virtual const char * type()
Return the type,"mroot_hybrids".
virtual int msolve_de(size_t nn, vec_t &xx, func_t &ufunc, jfunc_t &dfunc)
Solve func with derivatives dfunc using x as an initial guess, returning x.
problem with user-supplied function
Definition: err_hnd.h:69
virtual int msolve(size_t nn, vec_t &xx, func_t &ufunc)
Solve ufunc using xx as an initial guess, returning xx.
int dogleg(size_t n, const mat_t &r2, const ubvector &qtf2, const ubvector &diag2, double delta2, ubvector &newton2, ubvector &gradient2, vec2_t &p)
Take a dogleg step.
mat_t r
R matrix from QR decomposition.
void vector_copy(vec_t &src, vec2_t &dest)
Simple generic vector copy.
Definition: vector.h:73
int compute_trial_step(size_t N, vec2_t &xl, vec2_t &dxl, vec2_t &xx_trial)
Compute a trial step and put the result in xx_trial.
bool set_called
True if "set" has been called.
double tol_rel
The maximum value of the functions for success (default 1.0e-8)
Definition: mroot.h:82
bool jac_given
True if the jacobian has been given.
size_t dim
The number of equations and unknowns.
std::string itos(int x)
Convert an integer to a string.
static const double x2[5]
Definition: inte_qng_gsl.h:66
double scaled_enorm(size_t n, const vec2_t &d, const vec3_t &ff)
Compute the norm of the vector defined by .
double compute_actual_reduction(double fnorm0, double fnorm1)
Compute the actual reduction.
Definition: mroot_hybrids.h:81
Success.
Definition: err_hnd.h:47
void compute_qtf(size_t N, const vec2_t &q2, const vec3_t &f2, vec4_t &qtf2)
Compute .
Simple automatic Jacobian.
Definition: jacobian.h:126
Base for providing a numerical jacobian [abstract base].
Definition: jacobian.h:64
ubvector rdx
The value of .

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).
Hosted at Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads..