libLBFGS: a library of Limited-memory Broyden-Fletcher-Goldfarb-Shanno (L-BFGS)

Introduction

This library is a C port of the implementation of Limited-memory Broyden-Fletcher-Goldfarb-Shanno (L-BFGS) method written by Jorge Nocedal. The original FORTRAN source code is available at: http://www.ece.northwestern.edu/~nocedal/lbfgs.html

The L-BFGS method solves the unconstrainted minimization problem,

    minimize F(x), x = (x1, x2, ..., xN),

only if the objective function F(x) and its gradient G(x) are computable. The well-known Newton's method requires computation of the inverse of the hessian matrix of the objective function. However, the computational cost for the inverse hessian matrix is expensive especially when the objective function takes a large number of variables. The L-BFGS method iteratively finds a minimizer by approximating the inverse hessian matrix by information from last m iterations. This innovation saves the memory storage and computational time drastically for large-scaled problems.

Among the various ports of L-BFGS, this library provides several features:

This library is used by:

Download

libLBFGS is distributed under the term of the MIT license.

Third-party modules

History

Documentation

Sample code

#include <stdio.h>
#include <lbfgs.h>
static lbfgsfloatval_t evaluate(
void *instance,
const lbfgsfloatval_t *x,
lbfgsfloatval_t *g,
const int n,
const lbfgsfloatval_t step
)
{
int i;
lbfgsfloatval_t fx = 0.0;
for (i = 0;i < n;i += 2) {
lbfgsfloatval_t t1 = 1.0 - x[i];
lbfgsfloatval_t t2 = 10.0 * (x[i+1] - x[i] * x[i]);
g[i+1] = 20.0 * t2;
g[i] = -2.0 * (x[i] * g[i+1] + t1);
fx += t1 * t1 + t2 * t2;
}
return fx;
}
static int progress(
void *instance,
const lbfgsfloatval_t *x,
const lbfgsfloatval_t *g,
const lbfgsfloatval_t fx,
const lbfgsfloatval_t xnorm,
const lbfgsfloatval_t gnorm,
const lbfgsfloatval_t step,
int n,
int k,
int ls
)
{
printf("Iteration %d:\n", k);
printf(" fx = %f, x[0] = %f, x[1] = %f\n", fx, x[0], x[1]);
printf(" xnorm = %f, gnorm = %f, step = %f\n", xnorm, gnorm, step);
printf("\n");
return 0;
}
#define N 100
int main(int argc, char *argv[])
{
int i, ret = 0;
lbfgsfloatval_t fx;
lbfgsfloatval_t *x = lbfgs_malloc(N);
if (x == NULL) {
printf("ERROR: Failed to allocate a memory block for variables.\n");
return 1;
}
/* Initialize the variables. */
for (i = 0;i < N;i += 2) {
x[i] = -1.2;
x[i+1] = 1.0;
}
/* Initialize the parameters for the L-BFGS optimization. */
/*param.linesearch = LBFGS_LINESEARCH_BACKTRACKING;*/
/*
Start the L-BFGS optimization; this will invoke the callback functions
evaluate() and progress() when necessary.
*/
ret = lbfgs(N, x, &fx, evaluate, progress, NULL, &param);
/* Report the result. */
printf("L-BFGS optimization terminated with status code = %d\n", ret);
printf(" fx = %f, x[0] = %f, x[1] = %f\n", fx, x[0], x[1]);
return 0;
}

Acknowledgements

The L-BFGS algorithm is described in:

The line search algorithms used in this implementation are described in:

This library also implements Orthant-Wise Limited-memory Quasi-Newton (OWL-QN) method presented in:

Special thanks go to:

Finally I would like to thank the original author, Jorge Nocedal, who has been distributing the effieicnt and explanatory implementation in an open source licence.

Reference


Copyright (c) 2002-2014 by Naoaki Okazaki
Tue Sep 2 2014 11:35:37