File: __bfgsmin.cc

package info (click to toggle)
octave-optim 1.5.3-2
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 2,368 kB
  • sloc: cpp: 1,437; makefile: 185; perl: 169; xml: 29; sh: 3
file content (519 lines) | stat: -rw-r--r-- 19,462 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
// Copyright (C) 2004,2005,2006,2007,2010 Michael Creel <michael.creel@uab.es>
//
// This program is free software; you can redistribute it and/or modify it under
// the terms of the GNU General Public License as published by the Free Software
// Foundation; either version 3 of the License, or (at your option) any later
// version.
//
// This program is distributed in the hope that it will be useful, but WITHOUT
// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
// FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
// details.
//
// You should have received a copy of the GNU General Public License along with
// this program; if not, see <http://www.gnu.org/licenses/>.

// the functions defined in this file are:
// __bfgsmin_obj: bulletproofed objective function that allows checking for availability of analytic gradient
// __numgradient: numeric gradient, used only if analytic not supplied
// __bisectionstep: fallback stepsize algorithm
// __newtonstep: default stepsize algorithm
// __bfgsmin: the DLD function that does the minimization, to be called from bfgsmin.m

#include <oct.h>
#include <octave/parse.h>
#include <octave/Cell.h>
#include <octave/lo-mappers.h>
#include <float.h>
#include <climits>

#include "error-helpers.h"

int __bfgsmin_obj(double &obj, const std::string f, const octave_value_list f_args, const ColumnVector theta, const int minarg)
{
  octave_value_list f_return, f_args_new;
  int success = 1;
  f_args_new = f_args;
  f_args_new(minarg - 1) = theta;
  bool err;
  SET_ERR (f_return = OCTAVE__FEVAL (f, f_args_new), err);
  if (! err)
    {
      SET_ERR (obj = f_return(0).double_value(), err);
    }
  // bullet-proof the objective function
  if (err)
    {
      warning("__bfgsmin_obj: objective function could not be evaluated - setting to DBL_MAX");
      obj = DBL_MAX;
      success = 0;
    }
  return success;
}


// __numgradient: numeric central difference gradient for bfgs.
// This is the same as numgradient, except the derivative is known to be a vector, it's defined as a column,
// and the finite difference delta is incorporated directly rather than called from a function
int __numgradient(ColumnVector &derivative, const std::string f, const octave_value_list f_args, const int minarg)
{
  double SQRT_EPS, diff, delta, obj_left, obj_right, p;
  int j, test, success;
  ColumnVector parameter = f_args(minarg - 1).column_vector_value();
  int k = parameter.rows();
  ColumnVector g(k);
  SQRT_EPS = sqrt(DBL_EPSILON);
  diff = exp(log(DBL_EPSILON)/3.0);
  // get 1st derivative by central difference
  for (j=0; j<k; j++)
    {
      p = parameter(j);
      // determine delta for finite differencing
      test = (fabs(p) + SQRT_EPS) * SQRT_EPS > diff;
      if (test) delta = (fabs(p) + SQRT_EPS) * SQRT_EPS;
      else delta = diff;
      // right side
      parameter(j) = p + delta;
      success = __bfgsmin_obj(obj_right, f, f_args, parameter, minarg);
      if (! success)
        _p_error ("__numgradient: objective function failed, can't compute numeric gradient");
      // left size
      parameter(j) = p - delta;
      success = __bfgsmin_obj(obj_left, f, f_args, parameter, minarg);
      if (! success)
        _p_error ("__numgradient: objective function failed, can't compute numeric gradient");
      parameter(j) = p;  // restore original parameter for next round
      g(j) = (obj_right - obj_left) / (2.0*delta);
    }
  derivative = g;
  return success;
}

int __bfgsmin_gradient(ColumnVector &derivative, const std::string f, octave_value_list f_args, const ColumnVector theta, const int minarg, int try_analytic_gradient, int &have_analytic_gradient) {
  octave_value_list f_return;
  int k = theta.rows();
  int success;
  ColumnVector g(k);
  Matrix check_gradient(k,1);
  if (have_analytic_gradient)
    {
      f_args(minarg - 1) = theta;
      f_return = OCTAVE__FEVAL (f, f_args);
      g = f_return(1).column_vector_value();
    }
  else if (try_analytic_gradient)
    {
      f_args(minarg - 1) = theta;
      f_return = OCTAVE__FEVAL (f, f_args);
      if (f_return.length() > 1)
        {
          if (f_return(1).is_real_matrix())
            {
              if ((f_return(1).rows() == k) & (f_return(1).columns() == 1))
                {
                  g = f_return(1).column_vector_value();
                  have_analytic_gradient = 1;
                }
              else
                have_analytic_gradient = 0;
            }
          else
            have_analytic_gradient = 0;
        }
      else
        have_analytic_gradient = 0;
      if (!have_analytic_gradient) __numgradient(g, f, f_args, minarg);
    }
  else
    __numgradient(g, f, f_args, minarg);
  // check that gradient is ok
  check_gradient.column(0) = g;
  if (check_gradient.any_element_is_inf_or_nan ())
    {
      _p_error("__bfgsmin_gradient: gradient contains NaNs or Inf");
      success = 0;
    }
  else success = 1;
  derivative = g;
  return success;
}


// this is the lbfgs direction, used if control has 5 elements
ColumnVector lbfgs_recursion(const int memory, const Matrix sigmas, const Matrix gammas, ColumnVector d)
{
  if (memory == 0)
    {
      const int n = sigmas.columns();
      ColumnVector sig = sigmas.column(n-1);
      ColumnVector gam = gammas.column(n-1);
      // do conditioning if there is any memory
      double cond = gam.transpose()*gam;
      if (cond > 0)
        {
          cond = (sig.transpose()*gam) / cond;
          d = cond*d;
        }
      return d;
    }
  else
    {
      ColumnVector sig = sigmas.column(memory-1);
      ColumnVector gam = gammas.column(memory-1);
      double rho;
      rho = 1.0 / (gam.transpose() * sig);
      double alpha;
      alpha = rho * (sig.transpose() * d);
      d = d - alpha*gam;
      d = lbfgs_recursion(memory - 1, sigmas, gammas, d);
      d = d + (alpha - rho * gam.transpose() * d) * sig;
    }
  return d;
}

// __bisectionstep: fallback stepsize method if __newtonstep fails
int __bisectionstep(double &step, double &obj, const std::string f, const octave_value_list f_args, const ColumnVector x, const ColumnVector dx, const int minarg, const int verbose)
{
  double best_obj, improvement, improvement_0;
  ColumnVector trial;
  // initial values
  best_obj = obj;
  improvement_0 = 0.0;
  step = 1.0;
  trial = x + step*dx;
  __bfgsmin_obj(obj, f, f_args, trial, minarg);
  if (verbose) printf("bisectionstep: trial step: %g  obj value: %g\n", step, obj);
  // this first loop goes until an improvement is found
  while (obj >= best_obj)
    {
      if (step < 2.0*DBL_EPSILON)
        {
          if (verbose) warning("bisectionstep: unable to find improvement, setting step to zero");
          step = 0.0;
          return 0;
        }
      step = 0.5*step;
      trial = x + step*dx;
      __bfgsmin_obj(obj, f, f_args, trial, minarg);
      if (verbose)
        printf("bisectionstep: trial step: %g  obj value: %g  best_value: %g\n", step, obj, best_obj);
    }
  // now keep going until rate of improvement is too low, or reach max trials
  best_obj = obj;
  while (step > 2.0*DBL_EPSILON)
    {
      step = 0.5*step;
      trial = x + step*dx;
      __bfgsmin_obj(obj, f, f_args, trial, minarg);
      if (verbose) printf("bisectionstep: trial step: %g  obj value: %g\n", step, obj);
      // if improved, record new best and try another step
      if (obj < best_obj)
        {
          improvement = best_obj - obj;
          best_obj = obj;
          if (improvement > 0.5*improvement_0) {
            improvement_0 = improvement;
          }
          else break;
        }
      else
        {
          step = 2.0*step; // put it back to best found
          obj = best_obj;
          break;
        }
    }
  return 1;
}

// __newtonstep: default stepsize algorithm
int __newtonstep(double &step, double &obj, const std::string f, const octave_value_list f_args, const ColumnVector x, const ColumnVector dx, const int minarg, const int verbose)
{
  double obj_0, obj_left, obj_right, delta, inv_delta_sq, gradient, hessian;
  int found_improvement = 0;
  obj_0 = obj;
  delta = 0.001; // experimentation shows that this is a good choice
  inv_delta_sq = 1.0 / (delta*delta);
  ColumnVector x_right = x + delta*dx;
  ColumnVector x_left = x  - delta*dx;
  // right
  __bfgsmin_obj(obj_right, f, f_args, x_right, minarg);
  // left
  __bfgsmin_obj(obj_left, f, f_args, x_left, minarg);
  gradient = (obj_right - obj_left) / (2.0*delta);  // take central difference
  hessian =  inv_delta_sq*(obj_right - 2.0*obj_0 + obj_left);
  hessian = fabs(hessian); // ensures we're going in a decreasing direction
  if (hessian < 2.0*DBL_EPSILON) hessian = 1.0; // avoid div by zero
  step = - gradient / hessian;  // hessian inverse gradient: the Newton step
  //      step = (step < 1.0)*step + 1.0*(step >= 1.0); // maximum stepsize is 1.0 - conservative
  // ensure that this is improvement, and if not, fall back to bisection
  __bfgsmin_obj(obj, f, f_args, x + step*dx, minarg);
  if (verbose) printf("newtonstep: trial step: %g  obj value: %g\n", step, obj);
  if (obj > obj_0)
    {
      obj = obj_0;
      if (verbose) warning("__stepsize: no improvement with Newton step, falling back to bisection");
      found_improvement = __bisectionstep(step, obj, f, f_args, x, dx, minarg, verbose);
    }
  else found_improvement = 1;
  if (OCTAVE__MATH__ISNAN (obj)) {
    obj = obj_0;
    if (verbose) warning("__stepsize: objective function crash in Newton step, falling back to bisection");
    found_improvement = __bisectionstep(step, obj, f, f_args, x, dx, minarg, verbose);
  }
  else found_improvement = 1;
  return found_improvement;
}

DEFUN_DLD(__bfgsmin, args, ,"__bfgsmin: backend for bfgs minimization\n\
Users should not use this directly. Use bfgsmin.m instead")
{
  std::string f (args(0).string_value());
  Cell f_args_cell (args(1).cell_value());
  octave_value_list f_args, f_return; // holder for return items

  int max_iters, verbosity, criterion, minarg, convergence, iter, memory,
    gradient_ok, i, j, k, conv_fun, conv_param, conv_grad, have_gradient,
    try_gradient, warnings;
  double func_tol, param_tol, gradient_tol, stepsize, obj_value, obj_in,
    last_obj_value, denominator, test;
  Matrix H, H1, H2;
  ColumnVector thetain, d, g, g_new, p, q, sig, gam;

  // controls
  Cell control (args(2).cell_value());
  max_iters = control(0).int_value();
  if (max_iters == -1) max_iters = INT_MAX;
  verbosity = control(1).int_value();
  criterion = control(2).int_value();
  minarg = control(3).int_value();
  memory = control(4).int_value();
  func_tol = control(5).double_value();
  param_tol = control(6).double_value();
  gradient_tol = control(7).double_value();

  // want to see warnings?
  warnings = 0;
  if (verbosity == 3) warnings = 1;

  // copy cell contents over to octave_value_list to use feval()
  k = f_args_cell.numel ();
  f_args(k); // resize only once
  for (i = 0; i<k; i++) f_args(i) = f_args_cell(i);

  // get the minimization argument
  ColumnVector theta  = f_args(minarg - 1).column_vector_value();
  k = theta.rows();

  // containers for items in limited memory version
  Matrix sigmas(k, memory);
  Matrix gammas(k, memory);
  sigmas.fill(0.0);
  gammas.fill(0.0);

  // initialize things
  have_gradient = 0; // have analytic gradient
  try_gradient = 1;  // try to get analytic gradient
  convergence = -1; // if this doesn't change, it means that maxiters were exceeded
  thetain = theta;
  H = identity_matrix(k,k);

  // Initial obj_value
  __bfgsmin_obj(obj_in, f, f_args, theta, minarg);
  if (warnings) printf("initial obj_value %g\n", obj_in);

  // Initial gradient (try analytic, and use it if it's close enough to numeric)
  __bfgsmin_gradient(g, f, f_args, theta, minarg, 1, have_gradient);	// try analytic
  if (have_gradient)
    {
      // check equality if analytic available
      if (warnings) printf("function claims to provide analytic gradient\n");
      have_gradient = 0;				// force numeric
      __bfgsmin_gradient(g_new, f, f_args, theta, minarg, 0, have_gradient);
      p = g - g_new;
      have_gradient = sqrt(p.transpose() * p) < gradient_tol;
      if (have_gradient && warnings) printf("function claims to provide analytic gradient, and it agrees with numeric - using analytic\n");
      if (!have_gradient && warnings) printf("function claims to provide analytic gradient, but it does not agree with numeric - using numeric\n");
    }

  last_obj_value = obj_in; // initialize, is updated after each iteration
  // MAIN LOOP STARTS HERE
  for (iter = 0; iter < max_iters; iter++)
    {
      if(memory > 0)
        {  // lbfgs
          if (iter < memory) d = lbfgs_recursion(iter, sigmas, gammas, g);
          else d = lbfgs_recursion(memory, sigmas, gammas, g);
          d = -d;
        }
      else
        d = -H*g; // ordinary bfgs
      // convergence tests
      conv_fun = 0;
      conv_param = 0;
      conv_grad = 0;
      // function convergence
      p = theta+d;
      __bfgsmin_obj(obj_value, f, f_args, p, minarg);
      if (fabs(last_obj_value) > 1.0)	conv_fun=(fabs((obj_value/last_obj_value-1)))<func_tol;
      else conv_fun = fabs(obj_value - last_obj_value) < func_tol;
      // parameter change convergence
      test = sqrt(theta.transpose() * theta);
      if (test > 1.0) conv_param = sqrt(d.transpose() * d) / test < param_tol;
      else conv_param = sqrt(d.transpose() * d) < param_tol;		// Want intermediate results?
      // gradient convergence
      conv_grad = sqrt(g.transpose() * g) < gradient_tol;
      // Are we done?
      if (criterion == 1)
        {
          if (conv_fun && conv_param && conv_grad)
            {
              convergence = 1;
              break;
            }
        }
      else if (conv_fun)
        {
          convergence = 1;
          break;
        }
      // if not done, then take a step
      // stepsize: try (l)bfgs direction, then steepest descent if it fails
      f_args(minarg - 1) = theta;
      obj_value = last_obj_value;
      __newtonstep(stepsize, obj_value, f, f_args, theta, d, minarg, warnings);
      if (stepsize == 0.0)
        {
          // fall back to steepest descent
          if (warnings) warning("bfgsmin: BFGS direction fails, switch to steepest descent");
          d = -g; // try steepest descent
          H = identity_matrix(k,k); // accompany with Hessian reset, for good measure
          obj_value = last_obj_value;
          __newtonstep(stepsize, obj_value, f, f_args, theta, d, minarg, warnings);
          if (stepsize == 0.0)
            {  // if true, exit, we can't find a direction of descent
              warning("bfgsmin: failure, exiting. Try different start values?");
              f_return(0) = theta;
              f_return(1) = obj_value;
              f_return(2) = -1;
              f_return(3) = iter;
              return octave_value_list(f_return);
            }
        }
      p = stepsize*d;
      // Want intermediate results?
      if (verbosity > 1)
        {
          printf("------------------------------------------------\n");
          printf("bfgsmin iteration %d  convergence (f g p): %d %d %d\n", iter, conv_fun, conv_grad, conv_param);
          if (warnings)
            {
              if (memory > 0) printf("Using LBFGS, memory is last %d iterations\n",memory);
            }
          printf("\nfunction value: %g  stepsize: %g  \n\n", last_obj_value, stepsize);
          if (have_gradient) printf("used analytic gradient\n");
          else printf("used numeric gradient\n");
          for (j = 0; j<k; j++) printf("%15.5f %15.5f %15.5f\n",theta(j), g(j), p(j));
        }
      //--------------------------------------------------
      // // Are we done?
      // if (criterion == 1) {
      // 	if (conv_fun && conv_param && conv_grad) {
      // 		convergence = 1;
      // 		break;
      // 	}
      // }
      // else if (conv_fun) {
      // 	convergence = 1;
      // 	break;
      // }
      //-------------------------------------------------- 
      last_obj_value = obj_value;
      theta = theta + p;
      // new gradient
      gradient_ok = __bfgsmin_gradient(g_new, f, f_args, theta, minarg, try_gradient, have_gradient);
      if (memory == 0)
        {  //bfgs?
          // Hessian update if gradient ok
          if (gradient_ok)
            {
              q = g_new-g;
              g = g_new;
              denominator = q.transpose()*p;
              if ((fabs(denominator) < DBL_EPSILON))
                {  // reset Hessian if necessary
                  if (verbosity == 1) printf("bfgsmin: Hessian reset\n");
                  H = identity_matrix(k,k);
                }
              else
                {
                  H1 = (1.0+(q.transpose() * H * q) / denominator) / denominator
                    * (p * p.transpose());
                  H2 = (p * q.transpose() * H + H*q*p.transpose());
                  H2 = H2 / denominator;
                  H = H + H1 - H2;
                }
            }
          else
            H = identity_matrix(k,k); // reset hessian if gradient fails
          // then try to start again with steepest descent
        }
      else
        {  // otherwise lbfgs
          // save components for Hessian if gradient ok
          if (gradient_ok)
            {
              sig = p; // change in parameter
              gam = g_new - g; // change in gradient
              g = g_new;
              // shift remembered vectors to the right (forget last)
              for(j = memory - 1; j > 0; j--)
                {
                  for(i = 0; i < k; i++)
                    {
                      sigmas(i,j) = sigmas(i,j-1);
                      gammas(i,j) = gammas(i,j-1);
                    }
                }
              // insert new vectors in left-most column
              for(i = 0; i < k; i++)
                {
                  sigmas(i, 0) = sig(i);
                  gammas(i, 0) = gam(i);
                }
            }
          else
            {
              // failed gradient - loose memory and use previous theta
              sigmas.fill(0.0);
              gammas.fill(0.0);
              theta = theta - p;
            }
        }
    }
  // Want last iteration results?
  if (verbosity > 0)
    {
      printf("------------------------------------------------\n");
      printf("bfgsmin final results: %d iterations\n", iter);
      if (warnings)
        {
          if (memory > 0) printf("Used LBFGS, memory is last %d iterations\n",memory);
        }
      printf("\nfunction value: %g\n\n", obj_value);
      if (convergence == -1)                      printf("NO CONVERGENCE: max iters exceeded\n");
      if ((convergence == 1) & (criterion == 1))  printf("STRONG CONVERGENCE\n");
      if ((convergence == 1) & !(criterion == 1)) printf("WEAK CONVERGENCE\n");
      if (convergence == 2)                       printf("NO CONVERGENCE: algorithm failed\n");
      printf("Function conv %d  Param conv %d  Gradient conv %d\n\n", conv_fun, conv_param, conv_grad);
      if (have_gradient) printf("used analytic gradient\n");
      else printf("used numeric gradient\n");
      printf("          param    gradient (n)          change\n");
      for (j = 0; j<k; j++) printf("%15.5f %15.5f %15.5f\n",theta(j), g(j), d(j));
    }
  f_return(3) = iter;
  f_return(2) = convergence;
  f_return(1) = obj_value;
  f_return(0) = theta;
  return f_return;
}