Ice Sheet System Model  4.18
Code documentation
BrentSearch.cpp
Go to the documentation of this file.
1 
5 #ifdef HAVE_CONFIG_H
6  #include <config.h>
7 #else
8 #error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!"
9 #endif
10 
11 #include <float.h>
12 #include <iomanip>
13 #include <cmath>
14 
15 #include "../Exceptions/exceptions.h"
16 #include "../io/io.h"
17 #include "../MemOps/MemOps.h"
18 #include "./Verbosity.h"
19 #include "./OptPars.h"
20 #include "./types.h"
21 #include "./isnan.h"
22 
23 void BrentSearch(IssmDouble** pJ,OptPars optpars,IssmDouble* X0,IssmDouble (*f)(IssmDouble*,void*),IssmDouble (*g)(IssmDouble**,IssmDouble*,void*),void* usr){
24 
25  /* This routine is optimizing a given function using Brent's method
26  * (Golden or parabolic procedure)*/
27 
28  /*Intermediary*/
29  int iter;
30  IssmDouble si,gold,intervalgold,oldintervalgold;
31  IssmDouble parab_num,parab_den,distance;
32  IssmDouble fxmax,fxmin,fxbest;
33  IssmDouble fx,fx1,fx2;
34  IssmDouble x,x1,x2,xm,xbest;
35  IssmDouble tol1,tol2,seps;
36  IssmDouble tolerance = 1.e-4;
37 
38  /*Recover parameters:*/
39  int nsteps = optpars.nsteps;
40  int nsize = optpars.nsize;
41  IssmDouble xmin = optpars.xmin;
42  IssmDouble xmax = optpars.xmax;
43  int *maxiter = optpars.maxiter;
44  IssmDouble *cm_jump = optpars.cm_jump;
45 
46  /*Initialize gradient and controls*/
47  IssmDouble* G = NULL;
48  IssmDouble* J = xNew<IssmDouble>(nsteps);
49  IssmDouble* X = xNew<IssmDouble>(nsize);
50 
51  /*Header of printf*/
52  _printf0_("\n");
53  _printf0_(" x | Cost function f(x) | List of contributions\n");
54 
55  /*start iterations*/
56  for(int n=0;n<nsteps;n++){
57 
58  /*Print iteration number*/
59  _printf0_("====================== step "<< n+1 << "/" << nsteps <<" ===============================\n");
60 
61  /*Reset some variables*/
62  iter = 0;
63  xmin = 0.;
64  xmax = 1.;
65  bool loop = true;
66  cout<<setprecision(5);
67 
68  /*Get current Gradient at xmin=0*/
69  _printf0_(" x = "<<setw(9)<<xmin<<" | ");
70  fxmin = (*g)(&G,X0,usr); if(xIsNan<IssmDouble>(fxmin)) _error_("Function evaluation returned NaN");
71 
72  /*Get f(xmax)*/
73  _printf0_(" x = "<<setw(9)<<xmax<<" | ");
74  for(int i=0;i<nsize;i++) X[i]=X0[i]+xmax*G[i];
75  fxmax = (*f)(X,usr); if (xIsNan<IssmDouble>(fxmax)) _error_("Function evaluation returned NaN");
76  //if(VerboseControl()) _printf0_(" N/A "<<setw(12)<<xmax<<" "<<setw(12)<<fxmax<<" N/A boundary\n");
77 
78  /*test if jump option activated and xmin==0*/
79  if(!xIsNan<IssmDouble>(cm_jump[n]) && (xmin==0) && (fxmax/fxmin)<cm_jump[n]){
80  for(int i=0;i<nsize;i++) X0[i]=X0[i]+xmax*G[i];
81  xDelete<IssmDouble>(G);
82  J[n]=fxmax;
83  continue;
84  }
85 
86  /*initialize optimization variables*/
87  seps=sqrt(DBL_EPSILON); //precision of a IssmDouble
88  distance=0.0; //new_x=old_x + distance
89  gold=0.5*(3.0-sqrt(5.0)); //gold = 1 - golden ratio
90  intervalgold=0.0; //distance used by Golden procedure
91 
92  /*1: initialize the values of the 4 x needed (x1,x2,x,xbest)*/
93  x1=xmin+gold*(xmax-xmin);
94  x2=x1;
95  xbest=x1;
96  x=xbest;
97 
98  /*2: call the function to be evaluated*/
99  _printf0_(" x = "<<setw(9)<<x<<" | ");
100  for(int i=0;i<nsize;i++) X[i]=X0[i]+x*G[i];
101  fxbest = (*f)(X,usr); if(xIsNan<IssmDouble>(fxbest)) _error_("Function evaluation returned NaN");
102  iter++;
103 
104  /*3: update the other variables*/
105  fx1=fxbest;
106  fx2=fxbest;
107  /*xm is always in the middle of a and b*/
108  xm=0.5*(xmin+xmax);
109  /*update tolerances*/
110  tol1=seps*sqrt(pow(xbest,2))+tolerance/3.0;
111  tol2=2.0*tol1;
112 
113  /*4: print result*/
114  //if(VerboseControl())
115  //_printf0_(" "<<setw(5)<<iter<<" "<<setw(12)<<xbest<<" "<<setw(12)<<fxbest<<" "<<setw(12)<<pow(pow(xbest-xm,2),0.5)<<" initial\n");
116  if (!xIsNan<IssmDouble>(cm_jump[n]) && (xmin==0) && ((fxbest/fxmin)<cm_jump[n])){
117  //if(VerboseControl()) _printf0_(" optimization terminated: current x satisfies criteria 'cm_jump'=" << cm_jump[n] << "\n");
118  loop=false;
119  }
120 
121  while(loop){
122 
123  bool goldenflag=true;
124 
125  // Is a parabolic fit possible ?
126  if (sqrt(pow(intervalgold,2))>tol1){
127 
128  // Yes, so fit parabola
129  goldenflag=false;
130  parab_num=(xbest-x1)*(xbest-x1)*(fxbest-fx2)-(xbest-x2)*(xbest-x2)*(fxbest-fx1);;
131  parab_den=2.0*(xbest-x1)*(fxbest-fx2)-2.0*(xbest-x2)*(fxbest-fx1);
132 
133  //reverse p if necessary
134  if(parab_den>0.0){
135  parab_num=-parab_num;
136  }
137  parab_den=sqrt(pow(parab_den,2));
138  oldintervalgold=intervalgold;
139  intervalgold=distance;
140 
141  // Is the parabola acceptable (we use seps here because in some configuration parab_num==parab_den*(xmax-xbest)
142  // and the result is not repeatable anymore
143  if (( sqrt(pow(parab_num,2)) < sqrt(pow(0.5*parab_den*oldintervalgold,2))) &&
144  (parab_num>parab_den*(xmin-xbest)+seps) &&
145  (parab_num<parab_den*(xmax-xbest)-seps)){
146 
147  // Yes, parabolic interpolation step
148  distance=parab_num/parab_den;
149  x=xbest+distance;
150 
151  // f must not be evaluated too close to min_x or max_x
152  if (((x-xmin)<tol2) || ((xmax-x)<tol2)){
153  if ((xm-xbest)<0.0) si=-1;
154  else si=1;
155  //compute new distance
156  distance=tol1*si;
157  }
158  }
159  else{
160  // Not acceptable, must do a golden section step
161  goldenflag=true;
162  }
163  }
164 
165  //Golden procedure
166  if(goldenflag){
167  // compute the new distance d
168  if(xbest>=xm){
169  intervalgold=xmin-xbest;
170  }
171  else{
172  intervalgold=xmax-xbest;
173  }
174  distance=gold*intervalgold;
175  }
176 
177  // The function must not be evaluated too close to xbest
178  if(distance<0) si=-1;
179  else si=1;
180  if(sqrt(pow(distance,2))>tol1) x=xbest+si*sqrt(pow(distance,2));
181  else x=xbest+si*tol1;
182 
183  //evaluate function on x
184  _printf0_(" x = "<<setw(9)<<x<<" | ");
185  for(int i=0;i<nsize;i++) X[i]=X0[i]+x*G[i];
186  fx = (*f)(X,usr); if(xIsNan<IssmDouble>(fx)) _error_("Function evaluation returned NaN");
187  iter++;
188 
189  // Update a, b, xm, x1, x2, tol1, tol2
190  if (fx<=fxbest){
191  if (x>=xbest) xmin=xbest;
192  else xmax=xbest;
193  x1=x2; fx1=fx2;
194  x2=xbest; fx2=fxbest;
195  xbest=x; fxbest=fx;
196  }
197  else{ // fx > fxbest
198  if (x<xbest) xmin=x;
199  else xmax=x;
200  if ((fx<=fx2) || (x2==xbest)){
201  x1=x2; fx1=fx2;
202  x2=x; fx2=fx;
203  }
204  else if ( (fx <= fx1) || (x1 == xbest) || (x1 == x2) ){
205  x1=x; fx1=fx;
206  }
207  }
208  xm = 0.5*(xmin+xmax);
209  tol1=seps*pow(pow(xbest,2),0.5)+tolerance/3.0;
210  tol2=2.0*tol1;
211  //if(VerboseControl())
212  // _printf0_(" "<<setw(5)<<iter<<" "<<setw(12)<<x<<" "<<setw(12)<<fx<<" "<<setw(12)<<pow(pow(xbest-xm,2),0.5)<<
213  // " "<<(goldenflag?"golden":"parabolic")<<"\n");
214 
215  /*Stop the optimization?*/
216  if (sqrt(pow(xbest-xm,2)) < (tol2-0.5*(xmax-xmin))){
217  //if(VerboseControl()) _printf0_(" optimization terminated: current x satisfies criteria 'tolx'=" << tolerance << "\n");
218  loop=false;
219  }
220  else if (iter>=maxiter[n]){
221  //if(VerboseControl()) _printf0_(" exiting: Maximum number of iterations has been exceeded ('maxiter'=" << maxiter[n] << ")\n");
222  loop=false;
223  }
224  else if (!xIsNan<IssmDouble>(cm_jump[n]) && (xmin==0) && ((fxbest/fxmin)<cm_jump[n])){
225  //if(VerboseControl()) _printf0_(" optimization terminated: current x satisfies criteria 'cm_jump'=" << cm_jump[n] << "\n");
226  loop=false;
227  }
228  else{
229  //continue
230  loop=true;
231  }
232  }//end while
233 
234  //Now, check that the value on the boundaries are not better than current fxbest
235  if (fxbest>fxmin){
236  xbest=optpars.xmin; fxbest=fxmin;
237  }
238  if (fxbest>fxmax){
239  xbest=optpars.xmax; fxbest=fxmax;
240  }
241 
242  /*Assign output pointers: */
243  for(int i=0;i<nsize;i++) X0[i]=X0[i]+xbest*G[i];
244  xDelete<IssmDouble>(G);
245  J[n]=fxbest;
246  }
247 
248  /*return*/
249  xDelete<IssmDouble>(X);
250  *pJ=J;
251 }
OptPars::nsize
int nsize
Definition: OptPars.h:17
IssmDouble
double IssmDouble
Definition: types.h:37
_printf0_
#define _printf0_(StreamArgs)
Definition: Print.h:29
types.h
prototypes for types.h
OptPars::maxiter
int * maxiter
Definition: OptPars.h:15
BrentSearch
void BrentSearch(IssmDouble **pJ, OptPars optpars, IssmDouble *X0, IssmDouble(*f)(IssmDouble *, void *), IssmDouble(*g)(IssmDouble **, IssmDouble *, void *), void *usr)
Definition: BrentSearch.cpp:23
OptPars::xmax
IssmDouble xmax
Definition: OptPars.h:13
OptPars.h
place holder for optimization parameters
fx
IssmDouble fx(IssmDouble x, IssmDouble y, IssmDouble z, int testid)
Definition: fsanalyticals.cpp:218
OptPars
Definition: OptPars.h:10
OptPars::nsteps
int nsteps
Definition: OptPars.h:16
OptPars::xmin
IssmDouble xmin
Definition: OptPars.h:12
isnan.h
_error_
#define _error_(StreamArgs)
Definition: exceptions.h:49
Verbosity.h
OptPars::cm_jump
IssmDouble * cm_jump
Definition: OptPars.h:14