Title: | Limited-memory BFGS Optimization |
---|---|
Description: | A wrapper built around the libLBFGS optimization library by Naoaki Okazaki. The lbfgs package implements both the Limited-memory Broyden-Fletcher-Goldfarb-Shanno (L-BFGS) and the Orthant-Wise Quasi-Newton Limited-Memory (OWL-QN) optimization algorithms. The L-BFGS algorithm solves the problem of minimizing an objective, given its gradient, by iteratively computing approximations of the inverse Hessian matrix. The OWL-QN algorithm finds the optimum of an objective plus the L1-norm of the problem's parameters. The package offers a fast and memory-efficient implementation of these optimization routines, which is particularly suited for high-dimensional problems. |
Authors: | Antonio Coppola [aut, cre, cph], Brandon Stewart [aut, cph], Naoaki Okazaki [aut, cph], David Ardia [ctb, cph], Dirk Eddelbuettel [ctb, cph], Katharine Mullen [ctb, cph], Jorge Nocedal [ctb, cph] |
Maintainer: | Antonio Coppola <[email protected]> |
License: | GPL (>= 2) |
Version: | 1.2.1.2 |
Built: | 2024-12-07 06:29:54 UTC |
Source: | CRAN |
Performs function optimization using the Limited-memory Broyden-Fletcher-Goldfarb-Shanno (L-BFGS) and Orthant-Wise Limited-memory Quasi-Newton optimization (OWL-QN) algorithms. A wrapper to the libLBFGS library by Naoaki Okazaki, based on an implementation of the L-BFGS method written by Jorge Nocedal. Please note that significant portions of this help file are taken from Okazaki's original documentation. For further information, please refer to the libLBFGS page.
lbfgs(call_eval, call_grad, vars, environment=NULL, ..., invisible = 0, m = 6, epsilon = 1e-5, past = 0, delta = 0, max_iterations = 0, linesearch_algorithm = "LBFGS_LINESEARCH_DEFAULT", max_linesearch = 20, min_step = 1e-20, max_step = 1e+20, ftol = 1e-4, wolfe = 0.9, gtol = 0.9, orthantwise_c = 0, orthantwise_start = 0, orthantwise_end = length(vars))
lbfgs(call_eval, call_grad, vars, environment=NULL, ..., invisible = 0, m = 6, epsilon = 1e-5, past = 0, delta = 0, max_iterations = 0, linesearch_algorithm = "LBFGS_LINESEARCH_DEFAULT", max_linesearch = 20, min_step = 1e-20, max_step = 1e+20, ftol = 1e-4, wolfe = 0.9, gtol = 0.9, orthantwise_c = 0, orthantwise_start = 0, orthantwise_end = length(vars))
call_eval |
The function to be optimized. This should be either an R object taking in a numeric vector as its first parameter, and returning a scalar output, or an external pointer to a C++ function compiled using the |
call_grad |
A function returning the gradient vector of the objective. This should be either an R object taking in a numeric vector as its first parameter, and returning a scalar output, or an external pointer to a C++ function compiled using the |
vars |
A numeric vector containing the initial values for all variables. |
environment |
An R environment containing all extra arguments to be passed to the objective function and to the gradient, which must be matched exactly. If the objective function and the gradient are implemented in C++, extra arguments must be passed using this option, rather than the |
... |
Other arguments to be passed to |
invisible |
Defaults to |
m |
The number of corrections to approximate the inverse Hessian matrix. The L-BFGS routine stores the computation results of previous |
epsilon |
Epsilon for convergence test. This parameter determines the accuracy with which the solution is to be found. A minimization terminates when |
past |
Distance for delta-based convergence test. This parameter determines the distance, in iterations, to compute the rate of decrease of the objective function. If the value of this parameter is zero, the library does not perform the delta-based convergence test. The default value is |
delta |
Delta for convergence test. This parameter determines the minimum rate of decrease of the objective function. The library stops iterations when the following condition is met: |
max_iterations |
The maximum number of iterations. The |
linesearch_algorithm |
The line search algorithm. This parameter specifies a line search algorithm to be used by the L-BFGS routine. Valid arguments are the following:
If OWL-QN is invoked ( |
max_linesearch |
The maximum number of trials for the line search.This parameter controls the number of function and gradients evaluations per iteration for the line search routine. The default value is |
min_step |
The minimum step of the line search routine. The default value is |
max_step |
The maximum step of the line search. The default value is |
ftol |
A parameter to control the accuracy of the line search routine. The default value is |
wolfe |
A coefficient for the Wolfe condition. This parameter is valid only when the backtracking line-search algorithm is used with the Wolfe condition. The default value is |
gtol |
A parameter to control the accuracy of the line search routine. The default value is |
orthantwise_c |
Coefficient for the |
orthantwise_start |
Start index for computing L1 norm of the variables. This parameter is valid only for OWL-QN method (i.e., |
orthantwise_end |
End index for computing |
A list with the following components:
value |
The minimized value of the objective function. |
par |
A numerical array. The best set of parameters found. |
convergence |
An integer code. Zero indicates that convergence was reached without issues. Negative values indicate errors in the execution of the L-BFGS routine. |
message |
A character object detailing execution errors. This component is only returned if the convergence code is different form zero. |
# Rosenbrock Banana function objective <- function(x) { x1 <- x[1] x2 <- x[2] 100 * (x2 - x1 * x1)^2 + (1 - x1)^2 } gradient <- function(x) { ## Gradient of 'fr' x1 <- x[1] x2 <- x[2] c(-400 * x1 * (x2 - x1 * x1) - 2 * (1 - x1), 200 * (x2 - x1 * x1)) } output <- lbfgs(objective, gradient, c(-1.2,1)) # An example using OWL-QN to perform a Poisson regression using data from # Golub, Todd R., et al. "Molecular classification of cancer: class discovery # and class prediction by gene expression monitoring." Science 286.5439 (1999): # 531-537. A workspace with the dataset ("Leukemia.RData") is included # in the package distribution. # data(Leukemia) # X <- Leukemia$x # y <- Leukemia$y # X1 <- cbind(1, X) # pois.likelihood <- function(par, X, y, prec=0) { # Xbeta <- X%*%par # -(sum(y*Xbeta - exp(Xbeta)) -.5*sum(par^2*prec)) # } # pois.gradient <- function(par, X, y, prec=0) { # Xbeta <- X%*%par # expXbeta <- exp(Xbeta) # -(crossprod(X,(y-expXbeta)) -par*prec) # } # output <- lbfgs(pois.likelihood,pois.gradient, X=X1, y=y, prec=0, # rep(0, ncol(X1)), invisible=1, orthantwise_c=10, # linesearch_algorithm="LBFGS_LINESEARCH_BACKTRACKING", # orthantwise_start = 1, orthantwise_end = ncol(X1)) # Trivial Example objective <- function(x){ a <- x[1] b <- x[2] return(a^2 + b^2) } gradient <- function(x){ return(2*x) } output <- lbfgs(objective, gradient, c(100,13))
# Rosenbrock Banana function objective <- function(x) { x1 <- x[1] x2 <- x[2] 100 * (x2 - x1 * x1)^2 + (1 - x1)^2 } gradient <- function(x) { ## Gradient of 'fr' x1 <- x[1] x2 <- x[2] c(-400 * x1 * (x2 - x1 * x1) - 2 * (1 - x1), 200 * (x2 - x1 * x1)) } output <- lbfgs(objective, gradient, c(-1.2,1)) # An example using OWL-QN to perform a Poisson regression using data from # Golub, Todd R., et al. "Molecular classification of cancer: class discovery # and class prediction by gene expression monitoring." Science 286.5439 (1999): # 531-537. A workspace with the dataset ("Leukemia.RData") is included # in the package distribution. # data(Leukemia) # X <- Leukemia$x # y <- Leukemia$y # X1 <- cbind(1, X) # pois.likelihood <- function(par, X, y, prec=0) { # Xbeta <- X%*%par # -(sum(y*Xbeta - exp(Xbeta)) -.5*sum(par^2*prec)) # } # pois.gradient <- function(par, X, y, prec=0) { # Xbeta <- X%*%par # expXbeta <- exp(Xbeta) # -(crossprod(X,(y-expXbeta)) -par*prec) # } # output <- lbfgs(pois.likelihood,pois.gradient, X=X1, y=y, prec=0, # rep(0, ncol(X1)), invisible=1, orthantwise_c=10, # linesearch_algorithm="LBFGS_LINESEARCH_BACKTRACKING", # orthantwise_start = 1, orthantwise_end = ncol(X1)) # Trivial Example objective <- function(x){ a <- x[1] b <- x[2] return(a^2 + b^2) } gradient <- function(x){ return(2*x) } output <- lbfgs(objective, gradient, c(100,13))
Data from Golub, Todd R., et al. "Molecular classification of cancer: class discovery and class prediction by gene expression monitoring." Science 286.5439 (1999): 531-537. The study uses microarray data to perform cancer classification based on gene expression monitoring.
data(Leukemia)
data(Leukemia)
y |
A vector of binary values specifying the cancer class for 72 leukemia patients. A value of 0 corresponds to patients with acute lymphoblastic leukemia (ALL), and 1 corresponds to patients with acute myeloid leukemia (AML). |
x |
A 72-by-3571 matrix specifying the levels of expressions of 3571 genes for the 72 different patients. |