All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
mmin_simp2.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 /* multimin/simplex2.c
24  *
25  * Copyright (C) 2007 Brian Gough
26  * Copyright (C) 2002 Tuomo Keskitalo, Ivo Alxneit
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 #ifndef O2SCL_MMIN_SIMP2_H
44 #define O2SCL_MMIN_SIMP2_H
45 
46 /** \file mmin_simp2.h
47  \brief File defining \ref o2scl::mmin_simp2
48 */
49 #include <string>
50 
51 #include <boost/numeric/ublas/vector.hpp>
52 
53 #include <gsl/gsl_multimin.h>
54 
55 #include <o2scl/mmin.h>
56 #include <o2scl/cblas.h>
57 #include <o2scl/vec_stats.h>
58 
59 #ifndef DOXYGEN_NO_O2NS
60 namespace o2scl {
61 #endif
62 
63  /** \brief Multidimensional minimization by the Simplex method (v2) (GSL)
64 
65  This class mins a function using Nelder and Mead's Simplex
66  algorithm. A simplex in a N-dimensional space is defined as a
67  set of N+1 points which describe an N-dimensional volume
68  surrounding the minimum. The algorithm proceeds by shifting the
69  simplex points until the simplex is sufficiently small and thus
70  the minimum is known with sufficient accuracy.
71 
72  This class has a high-level interface using mmin(),
73  mmin_twovec() or mmin_simplex() which automatically performs the
74  memory allocation and minimization, or a GSL-like interface
75  using allocate(), free(), interate() and set() or set_simplex().
76 
77  The simplex can be completely specified by the user (see
78  mmin_simplex() and set_simplex()). Alternatively, the simplex is
79  automatically specified given initial guess \f$ x_j \f$ and a
80  step size vector \f$ s_k \f$ for \f$ 0\leq k<n_s \f$. The
81  simplex \f$ p_{ij} \f$ with \f$ 0\leq i\leq n \f$ and \f$ 0\leq
82  j<n \f$ is chosen with \f$ p_{0j} = x_j \f$ and
83  \f{eqnarray*}
84  p_{i+1,j} &=& x_j\quad\mathrm{for}\quad i\neq j \\
85  p_{i+1,j} &=& x_j+s_{j~\mathrm{mod}~n_s}\quad\mathrm{for}\quad i=j
86  \f}
87  for \f$ 0<i<n \f$. The step size vector \f$ s \f$ is set by the
88  set_step() member function. The prsence of \f$ \mathrm{mod} \f$
89  in the recipe above just indicates that elements of the step
90  size vector are automatically re-used if there are less step
91  sizes than dimensions in the minimization problem.
92 
93  \note It is important that the initial simplex contains
94  sufficient variation in every direction of the parameter space
95  over which one is minimizing. For example, if all three points
96  in a simplex for minimizing over a two-dimensional space contain
97  nearly the same value for the second parameter, then the
98  minimizer may only min the function with respect to the
99  first parameter.
100 
101  See an example for the usage of this class
102  in \ref ex_mmin_sect .
103 
104  Default template arguments
105  - \c param_t - no default
106  - \c func_t - \ref multi_funct11
107  - \c vec_t - \ref boost::numeric::ublas::vector < double >
108 
109  Based on \ref Nelder65 .
110 
111  A variable <tt>count</tt> originally defined in the GSL simplex
112  state is not present here, because it was unused.
113 
114  \future Double check that the updates in gsl-1.13 are included
115  here, and also add support for the nmsimplex2rand algorithm
116  in GSL.
117  */
118  template<class func_t=multi_funct11,
120  public mmin_base<func_t,func_t,vec_t> {
121 
122  public:
123 
125 
126 #ifndef DOXYGEN_INTERNAL
127 
128  protected:
129 
130  /// An array of n+1 vectors containing the simplex
131  vec_t *x1;
132  /** \brief The n+1 function values at the simplex points
133 
134  \comment
135  This can't be the template type because it has to be
136  a different size (n+1) rather than n.
137  \endcomment
138  */
140  /// Workspace vector 1
141  vec_t ws1;
142  /// Workspace vector 2
143  vec_t ws2;
144  /// Workspace vector 3
145  vec_t ws3;
146  /// Center of simplex
147  vec_t center;
148  /// Desc
149  vec_t delta;
150  /// Distance of vector from center
151  vec_t xmc;
152  /// Squared simplex size
153  double S2;
154 
155  /// Compute the center of the simplex
157  size_t i, j;
158  double val;
159 
160  for (j = 0; j < dim; j++) {
161  val = 0.0;
162  for (i = 0; i < dim+1; i++) {
163  val += x1[i][j];
164  }
165  val /= ((double)dim)+1;
166  center[j]=val;
167  }
168 
169  return 0;
170  }
171 
172  /** \brief Compute the size of the simplex
173 
174  Calculates simplex size as average sum of length of vectors
175  from simplex center to corner points:
176 
177  \f$ (1/n) \sum || y - y_{\mathrm{center}} || \f$
178  */
179  double compute_size() {
180 
181  size_t P=dim+1;
182  double ss=0.0;
183 
184  for(size_t i=0;i<P;i++) {
185  for(size_t j=0;j<dim;j++) ws1[j]=x1[i][j];
186  o2scl_cblas::daxpy(-1.0,dim,center,ws1);
187  double t=o2scl_cblas::dnrm2(dim,ws1);
188  ss += t*t;
189  }
190 
191  /* Store squared size in the state */
192  S2=ss/((double)P);
193  return sqrt(ss/(double)(P));
194  }
195 
196  /** \brief Move a corner of a simplex
197 
198  Moves a simplex corner scaled by coeff (negative value
199  represents mirroring by the middle point of the "other" corner
200  points) and gives new corner in xc and function value at xc
201  in \c newval.
202  */
203  virtual int try_corner_move(const double coeff, size_t corner,
204  vec_t &xc, func_t &f, size_t nvar,
205  double &newval) {
206 
207  size_t P=dim+1;
208 
209  /*
210  xc = (1-coeff)*(P/(P-1)) * center(all) +
211  ((P*coeff-1)/(P-1))*x_corner
212  */
213  double alpha=(1.0-coeff)*((double)P)/((double)dim);
214  double beta=(((double)P)*coeff-1.0)/((double)dim);
215 
216  for(size_t j=0;j<dim;j++) {
217  xc[j]=center[j]*alpha;
218  xc[j]+=x1[corner][j]*beta;
219  }
220 
221  newval=f(nvar,xc);
222 
223  return 0;
224  }
225 
226  /// Update point \c i in the simplex with values \c xx
227  virtual int update_point(size_t i, vec_t &xx, double val) {
228 
229  const size_t P=dim+1;
230 
231  /* Compute delta = x - x_orig */
232  for(size_t j=0;j<dim;j++) {
233  delta[j]=xx[j];
234  delta[j]-=x1[i][j];
235  }
236 
237  /* Compute xmc = x_orig - c */
238  for(size_t j=0;j<dim;j++) {
239  xmc[j]=x1[i][j];
240  xmc[j]-=center[j];
241  }
242 
243  /* Update size: S2' = S2 + (2/P) *
244  (x_orig - c).delta + (P-1)*(delta/P)^2
245  */
246  double d=o2scl_cblas::dnrm2(dim,delta);
247  double xmcd=o2scl_cblas::ddot(dim,xmc,delta);
248  S2 += (2.0 / ((double)P)) * xmcd +
249  ((((double)P) - 1.0) / ((double)P)) * (d * d / ((double)P));
250 
251  /* Update center: c' = c + (x - x_orig) / P */
252  double alpha=1.0/((double)P);
253  for(size_t j=0;j<dim;j++) {
254  center[j]-=alpha*x1[i][j];
255  center[j]+=alpha*xx[j];
256  }
257 
258  for(size_t j=0;j<dim;j++) {
259  x1[i][j]=xx[j];
260  }
261  y1[i]=val;
262 
263  return 0;
264  }
265 
266  /** \brief Contract the simplex towards the best point
267 
268  Function contracts the simplex in respect to best valued
269  corner. All corners besides the best corner are moved.
270 
271  \comment
272  The GSL version requires a vector for workspace
273  named 'xc' which is not required here.
274  \endcomment
275  */
276  virtual int contract_by_best(size_t best, func_t &f,
277  size_t nvar) {
278 
279  /* Function contracts the simplex in respect to best valued
280  corner. That is, all corners besides the best corner are
281  moved. (This function is rarely called in practice, since
282  it is the last choice, hence not optimised - BJG)
283  */
284  size_t i,j, it;
285  double newval;
286  bool failed;
287 
288  int status=success;
289 
290  for (i=0;i<dim+1;i++) {
291 
292  if (i!=best) {
293 
294  for (j=0;j<dim;j++) {
295  x1[i][j]=0.5*(x1[i][j]+x1[best][j]);
296  }
297  y1[i]=f(nvar,x1[i]);
298  if (!o2scl::is_finite(y1[i])) {
299  std::string err=((std::string)"Function not finite (returned ")+
300  dtos(y1[i])+" in mmin_simp2::contract_by_best().";
301  O2SCL_ERR(err.c_str(),exc_ebadfunc);
302  }
303 
304  /*
305  if (avoid_nonzero==true && ret!=0) {
306  std::cout << "Found problem in contract." << std::endl;
307 
308  // copy the old point to ws3
309  for (j=0;j<dim;j++) {
310  ws3[j]=2.0*x1[i][j]-x1[best][j];
311  }
312 
313  // Try (21*best+20*xold)/41, (22*best+19*xold)/41, ...,
314  // (40*best+xold)/41 to see if we can find a contraction
315  // that works
316  for(it=0;ret!=0 && it<20;it++) {
317  for (j=0;j<dim;j++) {
318  x1[i][j]=((20-it)*ws3[j]+(it+21)*x1[best][j])/41.0;
319  }
320  y1[i]=f(nvar,x1[i]);
321  std::cout << "it, ret, x: " << it << " " << ret << " "
322  << x << std::endl;
323  }
324  if (ret!=0) {
325  O2SCL_CONV2_RET("Failed to find suitable contraction ",
326  "in mmin_simp2::contract_by_best().",
327  exc_efailed,this->err_nonconv);
328  }
329  }
330  */
331  }
332  }
333 
334  compute_center();
335  size=compute_size();
336 
337  return success;
338  }
339 
340  /// Number of variables to be mind over
341  size_t dim;
342 
343  /// Function
344  func_t *func;
345 
346  /// True if set() has been called
348 
349  /// Vector of step sizes
351 
352  /** \brief If true, try to automatically avoid regions where the
353  function returns a non-zero value (default false)
354 
355  \note This option doesn't work yet, so I've made the variable
356  protected to prevent the casual user from changing it.
357  */
359 
360 #endif
361 
362  public:
363 
364  mmin_simp2() {
365  dim=0;
366  print_simplex=0;
367  step_vec.resize(1);
368  step_vec[0]=1.0;
369  avoid_nonzero=false;
370  }
371 
372  virtual ~mmin_simp2() {
373  }
374 
375  /// Set the step sizes for each independent variable
376  template<class vec2_t> int set_step(size_t nv, vec2_t &step) {
377  if (nv>0) {
378  step_vec.resize(nv);
379  for(size_t i=0;i<nv;i++) step_vec[i]=step[i];
380  }
381  return 0;
382  }
383 
384  /// Size of current simplex computed by iterate()
385  double size;
386 
387  /// Present minimum vector computed by iterate()
388  vec_t x;
389 
390  /// Function value at minimum computed by iterate()
391  double fval;
392 
393  /** \brief Print simplex information in print_iter() (default 0)
394 
395  If this is 1 and \ref verbose is greater than 0,
396  then print_iter() will print the function values at all
397  the simplex points. If this is 2 and \ref verbose is
398  greater than 0, then print_iter() will print the
399  simplex coordinates in addition to the function values.
400  */
402 
403  /** \brief Calculate the minimum \c min of \c func w.r.t the
404  array \c x of size \c nvar.
405  */
406  virtual int mmin(size_t nn, vec_t &xx, double &fmin,
407  func_t &ufunc) {
408 
409  if (nn==0) {
410  O2SCL_ERR2("Tried to min over zero variables ",
411  " in mmin_simp2::mmin().",exc_einval);
412  }
413 
414  int ret=0,status,iter=0;
415 
416  allocate(nn);
417 
418  vec_t ss(nn);
419  for (size_t is=0;is<nn;is++) ss[is]=step_vec[is % step_vec.size()];
420  ret=set(ufunc,nn,xx,ss);
421 
422  if(ret!=0) {
423  return ret;
424  }
425 
426  do {
427  iter++;
428 
429  status=iterate();
430  if(status) break;
431 
432  if(this->verbose>0) {
433  print_iter(nn,x,x1,fval,iter,size,this->tol_abs,
434  "mmin_simp2");
435  }
436 
437  status=gsl_multimin_test_size(size,this->tol_abs);
438 
439  } while(status == GSL_CONTINUE && iter<this->ntrial);
440 
441  for (size_t i=0;i<nn;i++) xx[i]=x[i];
442  fmin=fval;
443 
444  this->last_ntrial=iter;
445 
446  if(iter>=this->ntrial) {
447  std::string str="Exceeded maximum number of iterations ("+
448  itos(this->ntrial)+") in mmin_simp2::mmin().";
449  O2SCL_CONV_RET(str.c_str(),exc_emaxiter,this->err_nonconv);
450  }
451 
452  return status;
453  }
454 
455  /** \brief Calculate the minimum \c min of \c func w.r.t the
456  array \c x of size \c nvar, using \c xx and \c xx2 to specify
457  the simplex
458  */
459  virtual int mmin_twovec(size_t nn, vec_t &xx, vec_t &xx2, double &fmin,
460  func_t &ufunc) {
461 
462  int ret=0,i,status,iter=0;
463 
464  allocate(nn);
465 
466  vec_t ss(nn);
467  for (size_t is=0;is<nn;is++) ss[is]=xx2[is]-xx[is];
468  ret=set(ufunc,nn,xx,ss);
469 
470  if(ret!=0) {
471  return ret;
472  }
473 
474  do {
475  iter++;
476 
477  status=iterate();
478  if(status) break;
479 
480  if(this->verbose>0) {
481  print_iter(nn,x,x1,fval,iter,size,this->tol_abs,
482  "mmin_simp2");
483  }
484 
485  status=gsl_multimin_test_size(size,this->tol_abs);
486 
487  } while(status == GSL_CONTINUE && iter<this->ntrial);
488 
489  for (i=0;i<((int)nn);i++) xx[i]=x[i];
490  fmin=fval;
491 
492  this->last_ntrial=iter;
493 
494  if(iter>=this->ntrial) {
495  std::string str="Exceeded maximum number of iterations ("+
496  itos(this->ntrial)+") in mmin_simp2::mmin_twovec().";
497  O2SCL_CONV_RET(str.c_str(),exc_emaxiter,this->err_nonconv);
498  }
499 
500  return status;
501  }
502 
503  /** \brief Calculate the minimum \c min of \c func w.r.t the
504  array \c x of size \c nvar, given an initial simplex
505  */
506  template<class mat_t>
507  int mmin_simplex(size_t nn, mat_t &sx, double &fmin,
508  func_t &ufunc) {
509 
510  int ret=0,i,status,iter=0;
511 
512  allocate(nn);
513 
514  ret=set_simplex(ufunc,sx);
515  if(ret!=0) {
516  return ret;
517  }
518 
519  do {
520  iter++;
521 
522  status=iterate();
523  if(status) break;
524 
525  if(this->verbose>0) {
526  print_iter(nn,x,x1,fval,iter,size,this->tol_abs,
527  "mmin_simp2");
528  }
529 
530  status=gsl_multimin_test_size(size,this->tol_abs);
531 
532  } while(status == GSL_CONTINUE && iter<this->ntrial);
533 
534  for (i=0;i<((int)nn);i++) sx(0,i)=x[i];
535  fmin=fval;
536 
537  this->last_ntrial=iter;
538 
539  if (iter>=this->ntrial) {
540  std::string str="Exceeded maximum number of iterations ("+
541  itos(this->ntrial)+") in mmin_simp2::mmin_simplex().";
542  O2SCL_CONV_RET(str.c_str(),exc_emaxiter,this->err_nonconv);
543  }
544 
545  return status;
546  }
547 
548  /// Allocate the memory
549  virtual int allocate(size_t n) {
550 
551  set_called=false;
552  dim=n;
553 
554  x.resize(n);
555  x1=new vec_t[n+1];
556  for(size_t i=0;i<n+1;i++) x1[i].resize(n);
557  y1.resize(n+1);
558  ws1.resize(n);
559  ws2.resize(n);
560  ws3.resize(n);
561  center.resize(n);
562  delta.resize(n);
563  xmc.resize(n);
564 
565  return success;
566  }
567 
568  /// Set the function and initial guess
569  virtual int set(func_t &ufunc, size_t n, vec_t &ax,
570  vec_t &step_size) {
571  size_t i;
572 
573  if(dim!=n) allocate(n);
574  func=&ufunc;
575 
576  // Copy initial guess to x
577  for (i=0;i<dim;i++) x[i]=ax[i];
578 
579  // first point is the original x0
580 
581  y1[0]=ufunc(dim,ax);
582  if (!o2scl::is_finite(y1[0])) {
583  std::string err=((std::string)"Function not finite (returned ")+
584  dtos(y1[0])+" in mmin_simp2::set().";
585  O2SCL_ERR(err.c_str(),exc_ebadfunc);
586  }
587  for(i=0;i<dim;i++) x1[0][i]=ax[i];
588 
589  /* following points are initialized to x0+step_size */
590 
591  for (i=1;i<dim+1;i++) {
592  for(size_t j=0;j<dim;j++) x1[i][j]=x[j];
593  x1[i][i-1]=x1[i][i-1]+step_size[i-1];
594  y1[i]=ufunc(dim,x1[i]);
595  }
596 
597  /* Initialize simplex size */
598 
599  compute_center();
600  size=compute_size();
601 
602  set_called=true;
603 
604  return success;
605  }
606 
607  /// Set the function and initial simplex
608  template<class mat_t>
609  int set_simplex(func_t &ufunc, mat_t &sx) {
610 
611  if(dim==0) {
612  O2SCL_ERR2("Memory not allocated in ",
613  "mmin_simp2::set_simplex().",exc_ebadlen);
614  }
615 
616  func=&ufunc;
617 
618  for(size_t i=0;i<dim+1;i++) {
619  for(size_t j=0;j<dim;j++) {
620  x1[i][j]=sx(i,j);
621  }
622  y1[i]=ufunc(dim,x1[i]);
623  if (!o2scl::is_finite(y1[i])) {
624  std::string err=((std::string)"Function not finite (returned ")+
625  dtos(y1[i])+" in mmin_simp2::set_simplex().";
626  O2SCL_ERR(err.c_str(),exc_ebadfunc);
627  }
628  }
629 
630  /* Initialize simplex size */
631 
632  compute_center();
633  size=compute_size();
634 
635  set_called=true;
636 
637  return success;
638  }
639 
640  /// Perform an iteration
641  virtual int iterate() {
642 
643  size_t n=dim+1;
644  size_t i;
645  size_t lo, hi, s_hi;
646  double dhi,ds_hi,dlo;
647  int status;
648  double val,val2;
649 
650  /* get index of highest, second highest and lowest point */
651 
652  lo=hi=0;
653  dlo=dhi=y1[0];
654  s_hi=1;
655  ds_hi=y1[1];
656 
657  for (i=1;i<n;i++) {
658  val=y1[i];
659  if(val<dlo) {
660  dlo=val;
661  lo=i;
662  } else if (val > dhi) {
663  ds_hi=dhi;
664  s_hi=hi;
665  dhi=val;
666  hi=i;
667  } else if(val > ds_hi) {
668  ds_hi=val;
669  s_hi=i;
670  }
671  }
672 
673  /* reflect the highest value */
674 
675  int ret1=try_corner_move(-1.0,hi,ws1,*func,dim,val);
676 
677  if (avoid_nonzero && ret1!=0) {
678  std::cout << "Found problem move1: " << std::endl;
679  for (size_t it=0;it<20 && ret1!=0;it++) {
680  ret1=try_corner_move(-1.0+((double)it)/10.0,hi,ws1,
681  *func,dim,val);
682  std::cout << "it,ret: " << it << " " << ret1 << std::endl;
683  }
684  if (ret1!=0) {
685  O2SCL_ERR2("Failed to move corner (1) in ",
686  "mmin_simp2::iterate().",exc_efailed);
687  }
688  }
689 
690  if (o2scl::is_finite(val) && val<y1[lo]) {
691 
692  /* reflected point becomes lowest point,try expansion */
693 
694  int ret2=try_corner_move(-2.0,hi,ws2,*func,dim,val2);
695 
696  if (avoid_nonzero && ret2!=0) {
697  std::cout << "Found problem move2: " << std::endl;
698  for (size_t it=0;it<20 && ret2!=0;it++) {
699  ret2=try_corner_move(-2.0+((double)it)/6.7,hi,ws2,
700  *func,dim,val2);
701  std::cout << "it,ret: " << it << " " << ret2 << std::endl;
702  }
703  if (ret2!=0) {
704  O2SCL_ERR2("Failed to move corner (2) in ",
705  "mmin_simp2::iterate().",exc_efailed);
706  }
707  }
708 
709  if (o2scl::is_finite(val2) && val2<y1[lo]) {
710  update_point(hi,ws2,val2);
711  } else {
712  update_point(hi,ws1,val);
713  }
714 
715  } else if (!o2scl::is_finite(val) || val > y1[s_hi]) {
716 
717  /* reflection does not improve things enough */
718 
719  if (o2scl::is_finite(val) && val <= y1[hi]) {
720 
721  /* if trial point is better than highest point,replace
722  highest point */
723 
724  update_point(hi,ws1,val);
725  }
726 
727  /* try one dimensional contraction */
728 
729  int ret3=try_corner_move(0.5,hi,ws2,*func,dim,val2);
730 
731  if (avoid_nonzero && ret3!=0) {
732  std::cout << "Found problem move3: " << std::endl;
733  for (size_t it=0;it<20 && ret3!=0;it++) {
734  ret3=try_corner_move(0.025*((double)(19-it)),hi,ws2,
735  *func,dim,val2);
736  std::cout << "it,ret: " << it << " " << ret3 << std::endl;
737  }
738  if (ret3!=0) {
739  O2SCL_ERR2("Failed to move corner (2) in ",
740  "mmin_simp2::iterate().",exc_efailed);
741  }
742  }
743 
744  if (o2scl::is_finite(val2) && val2 <= y1[hi]) {
745 
746  update_point(hi,ws2,val2);
747 
748  } else {
749 
750  /* contract the whole simplex in respect to the best point */
751  status=contract_by_best(lo,*func,dim);
752  if(status != 0) {
753  O2SCL_ERR("Function contract_by_best() failed in iterate().",
754  exc_efailed);
755  }
756 
757  }
758 
759  } else {
760 
761  /* trial point is better than second highest point.
762  Replace highest point by it */
763 
764  update_point(hi,ws1,val);
765  }
766 
767  /* return lowest point of simplex as x */
768 
769  vector_min(dim+1,y1,lo,val);
770  for(i=0;i<dim;i++) x[i]=x1[lo][i];
771  fval=y1[lo];
772 
773  /* Update simplex size */
774 
775  if (S2 > 0) {
776  size=sqrt(S2);
777  } else {
778  /* recompute if accumulated error has made size invalid */
779  size=compute_size();
780  }
781 
782  return success;
783  }
784 
785  /** \brief Print out iteration information.
786 
787  Depending on the value of the variable verbose, this prints out
788  the iteration information. If verbose=0, then no information is
789  printed, while if verbose>1, then after each iteration, the
790  present values of x and y are output to std::cout along with the
791  iteration number. If verbose>=2 then each iteration waits for a
792  character.
793  */
794  virtual int print_iter(size_t nv, vec_t &xx, vec_t *simp,
795  double y, int iter,
796  double value, double limit,
797  std::string comment) {
798 
799  if (this->verbose<=0) return 0;
800 
801  size_t i;
802  char ch;
803 
804  (*this->outs) << comment << " Iteration: " << iter << std::endl;
805  (*this->outs) << "x: ";
806  for(i=0;i<nv;i++) (*this->outs) << xx[i] << " ";
807  (*this->outs) << std::endl;
808  if (print_simplex>0) {
809  if (print_simplex==2) {
810  (*this->outs) << "Simplex Values:" << std::endl;
811  for(i=0;i<nv+1;i++) {
812  (*this->outs) << i << ": ";
813  for(size_t j=0;j<nv;j++) {
814  (*this->outs) << simp[i][j] << " ";
815  }
816  (*this->outs) << ": " << y1[i] << std::endl;
817  }
818  } else {
819  (*this->outs) << "Simplex Values:" << std::endl;
820  for(i=0;i<nv+1;i++) {
821  (*this->outs) << y1[i] << " ";
822  }
823  (*this->outs) << std::endl;
824  }
825  }
826  (*this->outs) << "y: " << y << " Val: " << value << " Lim: "
827  << limit << std::endl;
828  if (this->verbose>1) {
829  (*this->outs) << "Press a key and type enter to continue. ";
830  (*this->ins) >> ch;
831  }
832 
833  return 0;
834  }
835 
836  /// Return string denoting type("mmin_simp2")
837  virtual const char *type() { return "mmin_simp2";}
838 
839 #ifndef DOXYGEN_INTERNAL
840 
841  private:
842 
844  (const mmin_simp2<func_t,vec_t> &);
845  mmin_simp2<func_t,vec_t>& operator=
846  (const mmin_simp2<func_t,vec_t>&);
847 
848 #endif
849 
850  };
851 
852 #ifndef DOXYGEN_NO_O2NS
853 }
854 #endif
855 
856 #endif
vec_t ws1
Workspace vector 1.
Definition: mmin_simp2.h:141
double S2
Squared simplex size.
Definition: mmin_simp2.h:153
virtual int update_point(size_t i, vec_t &xx, double val)
Update point i in the simplex with values xx.
Definition: mmin_simp2.h:227
std::istream * ins
Stream for verbose input.
Definition: mmin.h:173
ubvector step_vec
Vector of step sizes.
Definition: mmin_simp2.h:350
func_t * func
Function.
Definition: mmin_simp2.h:344
#define O2SCL_CONV_RET(d, n, b)
Set a "convergence" error and return the error value.
Definition: err_hnd.h:292
std::ostream * outs
Stream for verbose output.
Definition: mmin.h:170
virtual int allocate(size_t n)
Allocate the memory.
Definition: mmin_simp2.h:549
virtual int set(func_t &ufunc, size_t n, vec_t &ax, vec_t &step_size)
Set the function and initial guess.
Definition: mmin_simp2.h:569
bool is_finite(double x)
Return false if x is infinite or not a number.
vec_t delta
Desc.
Definition: mmin_simp2.h:149
invalid argument supplied by user
Definition: err_hnd.h:59
virtual int mmin_twovec(size_t nn, vec_t &xx, vec_t &xx2, double &fmin, func_t &ufunc)
Calculate the minimum min of func w.r.t the array x of size nvar, using xx and xx2 to specify the sim...
Definition: mmin_simp2.h:459
exceeded max number of iterations
Definition: err_hnd.h:73
Multidimensional minimization by the Simplex method (v2) (GSL)
Definition: mmin_simp2.h:119
int mmin_simplex(size_t nn, mat_t &sx, double &fmin, func_t &ufunc)
Calculate the minimum min of func w.r.t the array x of size nvar, given an initial simplex...
Definition: mmin_simp2.h:507
void vector_min(size_t n, const vec_t &data, size_t &index, data_t &val)
Compute the minimum of the first n elements of a vector.
Definition: vector.h:740
generic failure
Definition: err_hnd.h:61
double size
Size of current simplex computed by iterate()
Definition: mmin_simp2.h:385
double dnrm2(const size_t N, const vec_t &X)
Compute the norm of the vector X.
Definition: cblas_base.h:156
virtual int iterate()
Perform an iteration.
Definition: mmin_simp2.h:641
bool err_nonconv
If true, call the error handler if the routine does not "converge".
Definition: mmin.h:208
matrix, vector lengths are not conformant
Definition: err_hnd.h:89
Multidimensional minimization [abstract base].
Definition: mmin.h:163
vec_t xmc
Distance of vector from center.
Definition: mmin_simp2.h:151
bool avoid_nonzero
If true, try to automatically avoid regions where the function returns a non-zero value (default fals...
Definition: mmin_simp2.h:358
vec_t x
Present minimum vector computed by iterate()
Definition: mmin_simp2.h:388
#define O2SCL_ERR2(d, d2, n)
Set an error, two-string version.
Definition: err_hnd.h:281
std::string dtos(double x, int prec=6, bool auto_prec=false)
Convert a double to a string.
size_t dim
Number of variables to be mind over.
Definition: mmin_simp2.h:341
virtual int try_corner_move(const double coeff, size_t corner, vec_t &xc, func_t &f, size_t nvar, double &newval)
Move a corner of a simplex.
Definition: mmin_simp2.h:203
double ddot(const size_t N, const vec_t &X, const vec2_t &Y)
Compute .
Definition: cblas_base.h:123
vec_t * x1
An array of n+1 vectors containing the simplex.
Definition: mmin_simp2.h:131
double compute_size()
Compute the size of the simplex.
Definition: mmin_simp2.h:179
virtual int mmin(size_t nn, vec_t &xx, double &fmin, func_t &ufunc)
Calculate the minimum min of func w.r.t the array x of size nvar.
Definition: mmin_simp2.h:406
#define O2SCL_ERR(d, n)
Set an error with message d and code n.
Definition: err_hnd.h:273
int last_ntrial
The number of iterations for in the most recent minimization.
Definition: mmin.h:205
double tol_abs
The independent variable tolerance.
Definition: mmin.h:202
int set_simplex(func_t &ufunc, mat_t &sx)
Set the function and initial simplex.
Definition: mmin_simp2.h:609
void daxpy(const double alpha, const size_t N, const vec_t &X, vec2_t &Y)
Compute .
Definition: cblas_base.h:98
int set_step(size_t nv, vec2_t &step)
Set the step sizes for each independent variable.
Definition: mmin_simp2.h:376
vec_t center
Center of simplex.
Definition: mmin_simp2.h:147
vec_t ws3
Workspace vector 3.
Definition: mmin_simp2.h:145
problem with user-supplied function
Definition: err_hnd.h:69
double fval
Function value at minimum computed by iterate()
Definition: mmin_simp2.h:391
ubvector y1
The n+1 function values at the simplex points.
Definition: mmin_simp2.h:139
virtual int contract_by_best(size_t best, func_t &f, size_t nvar)
Contract the simplex towards the best point.
Definition: mmin_simp2.h:276
std::function< double(size_t, const boost::numeric::ublas::vector< double > &)> multi_funct11
Multi-dimensional function typedef.
Definition: multi_funct.h:45
vec_t ws2
Workspace vector 2.
Definition: mmin_simp2.h:143
std::string itos(int x)
Convert an integer to a string.
int print_simplex
Print simplex information in print_iter() (default 0)
Definition: mmin_simp2.h:401
int compute_center()
Compute the center of the simplex.
Definition: mmin_simp2.h:156
bool set_called
True if set() has been called.
Definition: mmin_simp2.h:347
Success.
Definition: err_hnd.h:47
virtual int print_iter(size_t nv, vec_t &xx, vec_t *simp, double y, int iter, double value, double limit, std::string comment)
Print out iteration information.
Definition: mmin_simp2.h:794
virtual const char * type()
Return string denoting type("mmin_simp2")
Definition: mmin_simp2.h:837
int ntrial
Maximum number of iterations.
Definition: mmin.h:196

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..