ADMB Documentation  -a65f1c97
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
lmnewton.cpp
Go to the documentation of this file.
1 
8 #include <admodel.h>
9 
10 #ifndef OPT_LIB
11  #include <cassert>
12  #include <climits>
13 #endif
14 
15 #ifdef ISZERO
16  #undef ISZERO
17 #endif
18 #define ISZERO(d) ((d)==0.0)
19 
20 #ifdef __cplusplus
21 extern "C" {
22 #endif
23 
24 int lbfgs_(long int *, long int *, double *, double *, double *, long int *,
25  double *, long int *, double *, double *, double *, long int *,long int *,
26  long int *);
27 
28 #ifdef __cplusplus
29 }
30 #endif
31 
33  const independent_variables& _x, int m)
34 {
36  if (m<=0)
37  {
38  cerr << "illegal value for number of steps in limited memory"
39  " quasi-newton method -- set equal to 10" << endl;
40  }
41  int noprintx=0;
42  int on;
43  int maxfn_option=0;
44  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-maxfn"))>-1)
45  {
46  maxfn_option=1;
47  }
48  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-nox"))>-1)
49  {
50  noprintx=1;
51  }
52 
53  // get the number of active
54  int nvar = initial_params::nvarcalc();
55 
56  double _crit=0;
57  int itn=0;
58  int ifn=0;
59  int nopt = 0;
60  // set the convergence criterion by command line
61  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-crit",nopt))>-1)
62  {
63  if (!nopt)
64  {
65  cerr << "Usage -crit option needs number -- ignored" << endl;
66  }
67  else
68  {
69  char * end;
70  _crit=strtod(ad_comm::argv[on+1],&end);
71  if (_crit<=0)
72  {
73  cerr << "Usage -crit option needs positive number -- ignored" << endl;
74  _crit=0.0;
75  }
76  }
77  }
78 
80  // set convergence criterion for this phase
81  if (!(!convergence_criteria))
82  {
86  }
87  if (!ISZERO(_crit))
88  {
89  crit = _crit;
90  }
91  if (!(!maximum_function_evaluations) && !maxfn_option)
92  {
96  }
97  dvariable vf=0.0;
98 
99  double xtol, f;
100  dvector diag(1,nvar);
101  //int j, n, iflag, icall;
102  int iflag, icall;
103  double fbest=1.e+100;
104  dvector g(1,nvar);
105  g.initialize();
106  dvector xbest(1,nvar);
107  dvector gbest(1,nvar);
108  //double t1, t2;
109  long int diagco=0;
110  int iprintx[2];
111  //double epsx;
112  //m = 35;
113  dvector w(1,nvar+2*m+2*nvar*m);
114  nopt=0;
115  int on1;
116  if ( (on1=option_match(ad_comm::argc,ad_comm::argv,"-iprint",nopt))>-1)
117  {
118  if (!nopt)
119  {
120  cerr << "Usage -iprint option needs integer -- ignored" << endl;
121  }
122  else
123  {
124  int jj=atoi(ad_comm::argv[on1+1]);
125  iprint=jj;
126  }
127  }
128  iprintx[0] = iprint;
129  iprintx[1] = 0;
130  crit = 1e-5;
131  xtol = 1e-16;
132  icall = 0;
133  iflag = 0;
134  long int linfo=0;
135 
136 L20:
137  f = 0.;
140  userfunction();
142  f=value(vf);
143  ifn++;
144  if (f<fbest)
145  {
146  fbest=f;
147  xbest=x;
148  gbest=g;
149  }
150 
151  gradcalc(nvar,g);
152  if(fmod(double(itn),double(iprint)) == 0)
153  {
154  if (iprint>0)
155  {
156  if (!itn)
157  {
158  ad_printf("\nInitial statistics: ");
159  }
160  else
161  {
162  ad_printf("\nIntermediate statistics: ");
163  }
164 
165  ad_printf(
166  "%d variables; iteration %ld; function evaluation %ld\n",
167  nvar, itn, ifn);
168 
169  if (!itn)
170  {
171  ad_printf(
172  "Function value %12.4le; maximum gradient component mag %12.4le\n",
173  f, max(fabs(g)));
174  }
175  else
176  {
177  ad_printf(
178  "Function value %12.4le; maximum gradient component mag %12.4le\n",
179  fbest, max(gbest)
180  );
181  }
182  if (!noprintx)
183  {
184  if (!itn)
185  fmmdisp(x, g, nvar, 0,noprintx);
186  else
187  fmmdisp(xbest, gbest, nvar, 0,noprintx);
188  }
189  }
190  }
191 
192  long int lnvar=nvar;
193  long int litn=itn;
194  long int liprintx= *iprintx;
195  long int liflag=iflag;
196  long int lm=m;
197  lbfgs_(&lnvar, &lm, &(x[1]) , &f, &(g[1]), &diagco, &(diag[1]),
198  &liprintx, &crit, &xtol, &(w[1]), &liflag,&litn,&linfo);
199  itn=int(litn);
200  iflag=int(liflag);
201 
202  if (iflag <= 0)
203  {
204  goto L50;
205  }
206  ++icall;
207  if (icall > maxfn)
208  {
209  goto L50;
210  }
211  goto L20;
212 L50:
213  if (iprint>0)
214  {
215  ad_printf("\nfinal statistics: ");
216 
217  ad_printf(
218  "%d variables; iteration %ld; function evaluation %ld\n",
219  nvar, itn, ifn);
220  ad_printf(
221  "Function value %12.4le; maximum gradient component mag %12.4le\n",
222  f, max(g));
223  fmmdisp(x, g, nvar, 0,noprintx);
224  }
226  x=xbest;
228 }
229 
231  (double& f, const independent_variables& _x, int m, int noprintx,
232  int maxfn, double crit)
233 {
235  if (m<=0)
236  {
237  cerr << "illegal value for number of steps in limited memory"
238  " quasi-newton method -- set equal to 10" << endl;
239  }
240  int on;
241 
242  // get the number of active
243  int nvar = initial_params::nvarcalc();
244 
245  double _crit=0;
246  int itn=0;
247  int ifn=0;
248  int nopt = 0;
249  // set the convergence criterion by command line
250  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-crit",nopt))>-1)
251  {
252  if (!nopt)
253  {
254  cerr << "Usage -crit option needs number -- ignored" << endl;
255  }
256  else
257  {
258  char* end;
259  _crit=strtod(ad_comm::argv[on+1],&end);
260  if (_crit<=0)
261  {
262  cerr << "Usage -crit option needs positive number -- ignored" << endl;
263  }
264  }
265  }
267  dvariable vf=0.0;
268 
269  double xtol;
270  dvector diag(1,nvar);
271  //int j, n, iflag, icall;
272  int iflag, icall;
273  double fbest=1.e+100;
274  dvector g(1,nvar);
275  g.initialize();
276  dvector xbest(1,nvar);
277  dvector gbest(1,nvar);
278  //double t1, t2;
279  long int diagco=0;
280  int iprintx[2];
281  //double epsx;
282  //m = 35;
283  dvector w(1,nvar+2*m+2*nvar*m);
284  iprintx[0] = iprint;
285  iprintx[1] = 0;
286  crit = 1e-5;
287  xtol = 1e-16;
288  icall = 0;
289  iflag = 0;
290  long int linfo=0;
291 
292 L20:
293  f = 0.;
296  userfunction();
298  f=value(vf);
299  ifn++;
300  if (f<fbest)
301  {
302  fbest=f;
303  xbest=x;
304  gbest=g;
305  }
306 
307  gradcalc(nvar,g);
308  if(fmod(double(itn),double(iprint)) == 0)
309  {
310  if (iprint>0)
311  {
312  if (!itn)
313  {
314  ad_printf("\nInitial statistics: ");
315  }
316  else
317  {
318  ad_printf("\nIntermediate statistics: ");
319  }
320 
321  ad_printf(
322  "%d variables; iteration %ld; function evaluation %ld\n",
323  nvar, itn, ifn);
324 
325  if (!itn)
326  {
327  ad_printf(
328  "Function value %12.4le; maximum gradient component mag %12.4le\n",
329  f, max(g));
330  }
331  else
332  {
333  ad_printf(
334  "Function value %12.4le; maximum gradient component mag %12.4le\n",
335  fbest, max(gbest));
336  }
337  if (!noprintx)
338  {
339  if (!itn)
340  fmmdisp(x, g, nvar, 0,noprintx);
341  else
342  fmmdisp(xbest, gbest, nvar, 0,noprintx);
343  }
344  }
345  }
346  long int lnvar=nvar;
347  long int litn=itn;
348  long int liprintx= *iprintx;
349  long int liflag=iflag;
350  long int lm=m;
351  lbfgs_(&lnvar, &lm, &(x[1]) , &f, &(g[1]), &diagco, &(diag[1]),
352  &liprintx, &crit, &xtol, &(w[1]), &liflag,&litn,&linfo);
353  itn=int(litn);
354  iflag=int(liflag);
355 
356  if (iflag <= 0)
357  {
358  goto L50;
359  }
360  ++icall;
361  if (icall > maxfn)
362  {
363  goto L50;
364  }
365  goto L20;
366 L50:
367  if (iprint>0)
368  {
369  ad_printf("\nfinal statistics: ");
370  ad_printf("%d variables; iteration %ld; function evaluation %ld\n",
371  nvar, itn, ifn);
372  ad_printf(
373  "Function value %12.4le; maximum gradient component mag %12.4le\n",
374  f, max(g));
375  fmmdisp(x, g, nvar, 0,noprintx);
376  }
377  //gradient_structure::set_NO_DERIVATIVES();
378  x=xbest;
379  f=fbest;
381 }
void fmmdisp(const dvector &x, const dvector &g, const int &nvar, int scroll_flag, int noprintx)
Description not yet available.
Definition: fmm_disp.cpp:102
static void set_NO_DERIVATIVES(void)
Disable accumulation of derivative information.
Definition: gradstrc.cpp:641
#define x
Vector of double precision numbers.
Definition: dvector.h:50
int lbfgs_(integer *n, integer *m, doublereal *x, doublereal *f, doublereal *g, logical *diagco, doublereal *diag, integer *iprint, doublereal *eps, doublereal *xtol, doublereal *w, integer *iflag, integer *iter, integer *info)
Definition: lbfgs.cpp:261
static char ** argv
Definition: fvar.hpp:8866
#define ISZERO(d)
Definition: lmnewton.cpp:18
df1_two_variable fabs(const df1_two_variable &x)
Definition: df12fun.cpp:891
static dvariable reset(const dvar_vector &x)
Definition: model.cpp:345
const int iprint
Definition: fvar.hpp:9504
ADMB variable vector.
Definition: fvar.hpp:2172
void gradcalc(int nvar, const dvector &g)
Definition: sgradclc.cpp:77
int atoi(adstring &s)
Returns a integer converted from input s.
Definition: atoi.cpp:20
static int nvarcalc()
Definition: model.cpp:152
void limited_memory_quasi_newton(const independent_variables &, int)
Definition: lmnewton.cpp:32
prnstream & endl(prnstream &)
Description not yet available.
Definition: fvar.hpp:1937
static int current_phase
Definition: admodel.h:842
static dvector maximum_function_evaluations
Definition: admodel.h:1917
#define min(a, b)
Definition: cbivnorm.cpp:188
static int argc
Definition: fvar.hpp:8863
int indexmax() const
Get maximum valid index.
Definition: dvector.h:204
static double gmax
Definition: admodel.h:2396
static objective_function_value * pobjfun
Definition: admodel.h:2394
Description not yet available.
void initialize(void)
Initialze all elements of dvector to zero.
Definition: dvect5.cpp:10
int option_match(int argc, char *argv[], const char *string)
Checks if the program has been invoked with a particular command line argument (&quot;string&quot;).
Definition: optmatch.cpp:25
#define w
static void set_YES_DERIVATIVES(void)
Enable accumulation of derivative information.
Definition: gradstrc.cpp:650
dvector value(const df1_one_vector &v)
Definition: df11fun.cpp:69
virtual void userfunction(void)=0
#define max(a, b)
Definition: cbivnorm.cpp:189
static dvector convergence_criteria
Definition: admodel.h:1915
Fundamental data type for reverse mode automatic differentiation.
Definition: fvar.hpp:1518
int ad_printf(FILE *stream, const char *format, Args...args)
Definition: fvar.hpp:9487