ADMB Documentation  -a65f1c97
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
amoeba.cpp
Go to the documentation of this file.
1 
5 #include <admodel.h>
6 #include <math.h>
7 #define NRANSI
8 //#include "nrutil.h"
9 #define NMAX 5000
10 
14 #define GET_PSUM \
15  for (j=1;j<=ndim;j++) {\
16  for (sum=0.0,i=1;i<=mpts;i++) sum += p[i][j];\
17  psum[j]=sum;}
18 
22 #define SWAP(a,b) {swap=(a);(a)=(b);(b)=swap;}
23 
24 
38 void function_minimizer::adamoeba(const dmatrix& _p, const dvector& _y,
39  int ndim, double ftol,int nfunk)
40 {
41  dmatrix& p=(dmatrix&) _p;
42  dvector& y=(dvector&) _y;
43  int i,ihi,ilo,inhi,j,mpts=ndim+1;
44  double rtol,sum,swap,ysave,ytry;
45 
46  dvector psum(1,ndim);
47  nfunk=0;
48  GET_PSUM
49  for (;;) {
50  ilo=1;
51  ihi = y[1]>y[2] ? (inhi=2,1) : (inhi=1,2);
52  for (i=1;i<=mpts;i++) {
53  if (y[i] <= y[ilo]) ilo=i;
54  if (y[i] > y[ihi]) {
55  inhi=ihi;
56  ihi=i;
57  } else if (y[i] > y[inhi] && i != ihi) inhi=i;
58  }
59  rtol=2.0*fabs(y[ihi]-y[ilo])/(fabs(y[ihi])+fabs(y[ilo]));
60  if (rtol < ftol) {
61  SWAP(y[1],y[ilo])
62  for (i=1;i<=ndim;i++) SWAP(p[1][i],p[ilo][i])
63  break;
64  }
65  if (nfunk >= NMAX)
66  {
67  cerr << "NMAX exceeded" << endl;
68  }
69  nfunk += 2;
70  ytry=amxxx(p,y,psum,ndim,ihi,-1.0);
71  if (ytry <= y[ilo])
72  ytry=amxxx(p,y,psum,ndim,ihi,2.0);
73  else if (ytry >= y[inhi]) {
74  ysave=y[ihi];
75  ytry=amxxx(p,y,psum,ndim,ihi,0.5);
76  if (ytry >= ysave) {
77  for (i=1;i<=mpts;i++) {
78  if (i != ilo) {
79  for (j=1;j<=ndim;j++)
80  p[i][j]=psum[j]=0.5*(p[i][j]+p[ilo][j]);
81 
82  dvariable vf=0.0;
85  userfunction();
87 
88  y[i]=value(vf);
89  }
90  }
91  nfunk += ndim;
92  GET_PSUM
93  }
94  } else --(nfunk);
95  }
96 }
97 #undef SWAP
98 #undef GET_PSUM
99 #undef NMAX
100 #undef NRANSI
101 #define NRANSI
102 
107 double function_minimizer::amxxx(const dmatrix& _p, const dvector& _y,
108  const dvector& _psum, int ndim, int ihi, double fac)
109 {
110  dmatrix& p=(dmatrix&) _p;
111  dvector& y=(dvector&) _y;
112  dvector& psum=(dvector&) _psum;
113  int j;
114  double fac1,fac2,ytry;
115 
116  dvector ptry(1,ndim);
117  fac1=(1.0-fac)/ndim;
118  fac2=fac1-fac;
119  for (j=1;j<=ndim;j++) ptry[j]=psum[j]*fac1-p[ihi][j]*fac2;
120 
121  dvariable vf=0.0;
124  userfunction();
126  ytry=value(vf);
127 
128  if (ytry < y[ihi]) {
129  y[ihi]=ytry;
130  for (j=1;j<=ndim;j++) {
131  psum[j] += ptry[j]-p[ihi][j];
132  p[ihi][j]=ptry[j];
133  }
134  }
135  return ytry;
136 }
137 #undef NRANSI
138 
150 void function_minimizer::neldmead(int n, dvector& _start, dvector& _xmin,
151  double *ynewlo, double reqmin, double delta,int *icount, int *numres,
152  int *ifault)
153 //
154 // Purpose:
155 //
156 // NELMIN minimizes a function using the Nelder-Mead algorithm.
157 //
158 // Discussion:
159 //
160 // This routine seeks the minimum value of a user-specified function.
161 //
162 // Simplex function minimisation procedure due to Nelder+Mead(1965),
163 // as implemented by O'Neill(1971, Appl.Statist. 20, 338-45), with
164 // subsequent comments by Chambers+Ertel(1974, 23, 250-1), Benyon(1976,
165 // 25, 97) and Hill(1978, 27, 380-2)
166 //
167 // This routine does not include a termination test using the
168 // fitting of a quadratic surface.
169 //
170 // Modified:
171 //
172 // October 2010 by Derek Seiple
173 //
174 // Author:
175 //
176 // FORTRAN77 version by R ONeill
177 // C++ version by John Burkardt
178 //
179 // Reference:
180 //
181 // John Nelder, Roger Mead,
182 // A simplex method for function minimization,
183 // Computer Journal,
184 // Volume 7, 1965, pages 308-313.
185 //
186 // R ONeill,
187 // Algorithm AS 47:
188 // Function Minimization Using a Simplex Procedure,
189 // Applied Statistics,
190 // Volume 20, Number 3, 1971, pages 338-345.
191 //
192 // Parameters:
193 //
194 // Input, int N, the number of variables.
195 //
196 // Input/output, double START[N]. On input, a starting point
197 // for the iteration. On output, this data may have been overwritten.
198 //
199 // Output, double XMIN[N], the coordinates of the point which
200 // is estimated to minimize the function.
201 //
202 // Output, double YNEWLO, the minimum value of the function.
203 //
204 // Input, double REQMIN, the terminating limit for the variance
205 // of function values.
206 //
207 // Input, int KONVGE, the convergence check is carried out
208 // every KONVGE iterations.
209 //
210 // Input, int KCOUNT, the maximum number of function
211 // evaluations.
212 //
213 // Output, int *ICOUNT, the number of function evaluations
214 // used.
215 //
216 // Output, int *NUMRES, the number of restarts.
217 //
218 // Output, int *IFAULT, error indicator.
219 // 0, no errors detected.
220 // 1, REQMIN, N, or KONVGE has an illegal value.
221 // 2, iteration terminated because KCOUNT was exceeded without convergence.
222 //
223 {
224  dvector& start=(dvector&) _start;
225  dvector& xmin=(dvector&) _xmin;
226  int slb = start.indexmin();
227  int xlb = xmin.indexmin();
228  start.shift(0);
229  xmin.shift(0);
230 
231  int i;
232  dvector step(0,n-1);
233  for(i=0;i<n;i++)
234  {
235  step(i)=delta;
236  }
237 
238  int konvge = 10;
239  int kcount = 5000;
240  double ccoeff = 0.5;
241  double del;
242  double dn;
243  double dnn;
244  double ecoeff = 2.0;
245  double eps = 0.001;
246 
247  int ihi;
248  int ilo;
249  int j;
250  int jcount;
251  int l;
252  int nn;
253  double rcoeff = 1.0;
254  double rq;
255  double x;
256  double y2star;
257  double ylo;
258  double ystar;
259  double z;
260 
261  // Check the input parameters.
262  if ( reqmin <= 0.0 )
263  {
264  *ifault = 1;
265  return;
266  }
267 
268  if ( n < 1 )
269  {
270  *ifault = 1;
271  return;
272  }
273 
274  if ( konvge < 1 )
275  {
276  *ifault = 1;
277  return;
278  }
279 
280  dvector p(0,n*(n+1)-1);
281  dvector pstar(0,n-1);
282  dvector p2star(0,n-1);
283  dvector pbar(0,n-1);
284  dvector y(0,n);
285 
286  *icount = 0;
287  *numres = 0;
288 
289  jcount = konvge;
290  dn = ( double ) ( n );
291  nn = n + 1;
292  dnn = ( double ) ( nn );
293  del = 1.0;
294  rq = reqmin * dn;
295 
296  // Initial or restarted loop.
297  dvariable vf=0.0;
298 
299  for ( ; ; )
300  {
301  for ( i = 0; i < n; i++ )
302  {
303  p[i+n*n] = start[i];
304  }
305  start.shift(1);
306  vf=0.0;
309  userfunction();
311  y[n] = value(vf);
312  start.shift(0);
313 
314  *icount = *icount + 1;
315 
316  for ( j = 0; j < n; j++ )
317  {
318  x = start[j];
319  start[j] = start[j] + step[j] * del;
320  for ( i = 0; i < n; i++ )
321  {
322  p[i+j*n] = start[i];
323  }
324  start.shift(1);
325  vf=0.0;
328  userfunction();
330  y[j]=value(vf);
331  start.shift(0);
332 
333  *icount = *icount + 1;
334  start[j] = x;
335  }
336  // The simplex construction is complete.
337 
338  // Find highest and lowest Y values. YNEWLO = Y(IHI) indicates
339  // the vertex of the simplex to be replaced.
340  ylo = y[0];
341  ilo = 0;
342 
343  for ( i = 1; i < nn; i++ )
344  {
345  if ( y[i] < ylo )
346  {
347  ylo = y[i];
348  ilo = i;
349  }
350  }
351 
352  // Inner loop.
353  for ( ; ; )
354  {
355  if ( kcount <= *icount )
356  {
357  break;
358  }
359  *ynewlo = y[0];
360  ihi = 0;
361 
362  for ( i = 1; i < nn; i++ )
363  {
364  if ( *ynewlo < y[i] )
365  {
366  *ynewlo = y[i];
367  ihi = i;
368  }
369  }
370 
371  // Calculate PBAR, the centroid of the simplex vertices
372  // excepting the vertex with Y value YNEWLO.
373  for ( i = 0; i < n; i++ )
374  {
375  z = 0.0;
376  for ( j = 0; j < nn; j++ )
377  {
378  z = z + p[i+j*n];
379  }
380  z = z - p[i+ihi*n];
381  pbar[i] = z / dn;
382  }
383 
384  // Reflection through the centroid.
385  for ( i = 0; i < n; i++ )
386  {
387  pstar[i] = pbar[i] + rcoeff * ( pbar[i] - p[i+ihi*n] );
388  }
389  pstar.shift(1);
390  vf=0.0;
393  userfunction();
395  ystar = value(vf);
396  pstar.shift(0);
397 
398  *icount = *icount + 1;
399 
400  // Successful reflection, so extension.
401  if ( ystar < ylo )
402  {
403  for ( i = 0; i < n; i++ )
404  {
405  p2star[i] = pbar[i] + ecoeff * ( pstar[i] - pbar[i] );
406  }
407  p2star.shift(1);
408  vf=0.0;
411  userfunction();
413  y2star = value(vf);
414  p2star.shift(0);
415 
416  *icount = *icount + 1;
417 
418  // Check extension.
419  if ( ystar < y2star )
420  {
421  for ( i = 0; i < n; i++ )
422  {
423  p[i+ihi*n] = pstar[i];
424  }
425  y[ihi] = ystar;
426  }
427 
428  // Retain extension or contraction.
429  else
430  {
431  for ( i = 0; i < n; i++ )
432  {
433  p[i+ihi*n] = p2star[i];
434  }
435  y[ihi] = y2star;
436  }
437  }
438 
439  // No extension.
440  else
441  {
442  l = 0;
443  for ( i = 0; i < nn; i++ )
444  {
445  if ( ystar < y[i] )
446  {
447  l = l + 1;
448  }
449  }
450 
451  if ( 1 < l )
452  {
453  for ( i = 0; i < n; i++ )
454  {
455  p[i+ihi*n] = pstar[i];
456  }
457  y[ihi] = ystar;
458  }
459 
460  // Contraction on the Y(IHI) side of the centroid.
461  else if ( l == 0 )
462  {
463  for ( i = 0; i < n; i++ )
464  {
465  p2star[i] = pbar[i] + ccoeff * ( p[i+ihi*n] - pbar[i] );
466  }
467  p2star.shift(1);
468  vf=0.0;
471  userfunction();
473  y2star = value(vf);
474  p2star.shift(0);
475 
476  *icount = *icount + 1;
477 
478  // Contract the whole simplex.
479  if ( y[ihi] < y2star )
480  {
481  for ( j = 0; j < nn; j++ )
482  {
483  for ( i = 0; i < n; i++ )
484  {
485  p[i+j*n] = ( p[i+j*n] + p[i+ilo*n] ) * 0.5;
486  xmin[i] = p[i+j*n];
487  }
488  xmin.shift(1);
489  vf=0.0;
492  userfunction();
494  y[j] = value(vf);
495  xmin.shift(0);
496 
497  *icount = *icount + 1;
498  }
499  ylo = y[0];
500  ilo = 0;
501 
502  for ( i = 1; i < nn; i++ )
503  {
504  if ( y[i] < ylo )
505  {
506  ylo = y[i];
507  ilo = i;
508  }
509  }
510  continue;
511  }
512 
513  // Retain contraction.
514  else
515  {
516  for ( i = 0; i < n; i++ )
517  {
518  p[i+ihi*n] = p2star[i];
519  }
520  y[ihi] = y2star;
521  }
522  }
523 
524  // Contraction on the reflection side of the centroid.
525  else if ( l == 1 )
526  {
527  for ( i = 0; i < n; i++ )
528  {
529  p2star[i] = pbar[i] + ccoeff * ( pstar[i] - pbar[i] );
530  }
531  p2star.shift(1);
532  vf=0.0;
535  userfunction();
537  y2star = value(vf);
538  p2star.shift(0);
539 
540  *icount = *icount + 1;
541 
542  // Retain reflection?
543  if ( y2star <= ystar )
544  {
545  for ( i = 0; i < n; i++ )
546  {
547  p[i+ihi*n] = p2star[i];
548  }
549  y[ihi] = y2star;
550  }
551  else
552  {
553  for ( i = 0; i < n; i++ )
554  {
555  p[i+ihi*n] = pstar[i];
556  }
557  y[ihi] = ystar;
558  }
559  }
560  }
561 
562  // Check if YLO improved.
563  if ( y[ihi] < ylo )
564  {
565  ylo = y[ihi];
566  ilo = ihi;
567  }
568  jcount = jcount - 1;
569 
570  if ( 0 < jcount )
571  {
572  continue;
573  }
574 
575  // Check to see if minimum reached.
576  if ( *icount <= kcount )
577  {
578  jcount = konvge;
579 
580  z = 0.0;
581  for ( i = 0; i < nn; i++ )
582  {
583  z = z + y[i];
584  }
585  x = z / dnn;
586 
587  z = 0.0;
588  for ( i = 0; i < nn; i++ )
589  {
590  z = z + pow ( y[i] - x, 2 );
591  }
592 
593  if ( z <= rq )
594  {
595  break;
596  }
597  }
598  }
599 
600  // Factorial tests to check that YNEWLO is a local minimum.
601  for ( i = 0; i < n; i++ )
602  {
603  xmin[i] = p[i+ilo*n];
604  }
605  *ynewlo = y[ilo];
606 
607  if ( kcount < *icount )
608  {
609  *ifault = 2;
610  break;
611  }
612 
613  *ifault = 0;
614 
615  for ( i = 0; i < n; i++ )
616  {
617  del = step[i] * eps;
618  xmin[i] = xmin[i] + del;
619  xmin.shift(1);
620  vf=0.0;
623  userfunction();
625  z = value(vf);
626  xmin.shift(0);
627 
628  *icount = *icount + 1;
629  if ( z < *ynewlo )
630  {
631  *ifault = 2;
632  break;
633  }
634  xmin[i] = xmin[i] - del - del;
635  xmin.shift(1);
636  vf=0.0;
639  userfunction();
641  z = value(vf);
642  xmin.shift(0);
643 
644  *icount = *icount + 1;
645  if ( z < *ynewlo )
646  {
647  *ifault = 2;
648  break;
649  }
650  xmin[i] = xmin[i] + del;
651  }
652 
653  if ( *ifault == 0 )
654  {
655  break;
656  }
657 
658  // Restart the procedure.
659  for ( i = 0; i < n; i++ )
660  {
661  start[i] = xmin[i];
662  }
663  del = eps;
664  *numres = *numres + 1;
665  }
666 
667  start.shift(slb);
668  xmin.shift(xlb);
669 
670  return;
671 }
#define x
Vector of double precision numbers.
Definition: dvector.h:50
int indexmin() const
Get minimum valid index.
Definition: dvector.h:199
void neldmead(int n, dvector &_start, dvector &_xmin, double *ynewlo, double reqmin, double delta, int *icount, int *numres, int *ifault)
Nelder-Mead simplex alogrithm.
Definition: amoeba.cpp:150
double sum(const d3_array &darray)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: d3arr.cpp:21
df1_two_variable fabs(const df1_two_variable &x)
Definition: df12fun.cpp:891
static dvariable reset(const dvar_vector &x)
Definition: model.cpp:345
ADMB variable vector.
Definition: fvar.hpp:2172
void adamoeba(const dmatrix &p, const dvector &y, int ndim, double ftol, int maxfn)
Nelder-Mead simplex alogrithm.
Definition: amoeba.cpp:38
#define NMAX
Definition: amoeba.cpp:9
prnstream & endl(prnstream &)
static objective_function_value * pobjfun
Definition: admodel.h:2394
Description not yet available.
Description not yet available.
Definition: fvar.hpp:2819
dvector & shift(int min)
Shift valid range of subscripts.
Definition: dvector.cpp:52
double amxxx(const dmatrix &p, const dvector &y, const dvector &psum, int ndim, int ihi, double fac)
Description not yet available.
Definition: amoeba.cpp:107
#define SWAP(a, b)
Reverses a and b.
Definition: amoeba.cpp:22
double eps
Definition: ftweak.cpp:13
dvector value(const df1_one_vector &v)
Definition: df11fun.cpp:69
virtual void userfunction(void)=0
#define GET_PSUM
Computes psum.
Definition: amoeba.cpp:14
Fundamental data type for reverse mode automatic differentiation.
Definition: fvar.hpp:1518
d3_array pow(const d3_array &m, int e)
Description not yet available.
Definition: d3arr6.cpp:17