A reference class for computing sparse Hessians

Details

The sparseHessianFD function calls the initializer for the sparseHessianFD class, and returns a sparseHessianFD object.

sparseHessianFD(x, fn, gr, rows, cols, delta, index1, complex, ...)

The function, gradient and sparsity pattern are declared as part of the initialization.

Once initialized, the $hessian method will evaluate the Hessian at x.


obj <- sparseHessian(x, fn, gr, rows, cols, ...)
obj$hessian(x)

For convenience, the class provides wrapper methods to the fn and gr functions that were specified in the initializer.


obj$fn(P) ## wrapper to objective function
obj$gr(P) ## wrapper to gradient
obj$fngr(P) ## list of obj function and gradient
obj$fngrhs(P) ## list of obj function, gradient and Hessian.

Arguments to initializer

x

an vector at which the function, gradient and Hessian are initialized and tested.

fn, gr

R functions that return the function value and gradient, evaluated at x.

rows, cols

Numeric vectors with row and column indices of the non-zero elements in the lower triangle (including diagonal) of the Hessian.

delta

The perturbation amount for finite difference (or complex step) of the gradient to compute the Hessian. Defaults to 1e-07.

index1

TRUE if rows and cols use 1-based (R format) indexing (FALSE for 0-based (C format) indexing.

complex

TRUE if Hessian will be computed using the complex step method, and FALSE (default) if using finite differences. If TRUE, both fn and gr must accept complex arguments and return complex values.

...

other arguments to be passed to fn and gr.

Other methods are described below. Do not access any of the fields directly. The internal structure is subject to change in future versions.

Fields

fn1

A closure for calling fn(x, ...).

gr1

A closure for calling gr(x, ...).

iRow,jCol

Numeric vectors with row and column indices of the non-zero elements in the lower triangle (including diagonal) of the Hessian.

delta

The perturbation amount for finite differencing of the gradient to compute the Hessian. Defaults to 1e-07.

index1

TRUE if rows and cols use 1-based (R format) indexing (FALSE for 0-based (C format) indexing.

complex

TRUE if Hessian will be computed using the complex step method, and FALSE (default) if using finite differences.

D

raw finite differences (internal use only)

nvars

Number of variables (length of x)

nnz

Number of non-zero elements in the lower triangle of the Hessian.

ready

TRUE if object has been initialized, and Hessian has been partitioned.

idx,pntr

Column indices and row pointers for non-zero elements in lower triangle of the permuted Hessian. Row-oriented compressed storage.

colors

A vector representation of the partitioning of the columns. There are nvars elements, one for each column of the permuted Hessian. The value corresponds to the "color" for that column.

perm,invperm

Permutation vector and its inverse

Methods

fn(x)

Return function value, evaluated at x: fn(x, ...)

fngr(x)

Return list of function value and gradient, evaluated at x

fngrhs(x)

Return list of function value, gradient, and Hessian, evaluated at x

get_invperm()

Return integer vector of inverse of permutation used for computing Hessian

get_nnz()

Return number of non-zero elements in lower triangle of Hessian

get_nvars()

Return dimension (number of rows or columns) of Hessian

get_pattern()

Return pattern matrix of lower triangle of Hessian

get_perm()

Return integer vector of permutation used for computing Hessian

get_perm_pattern()

Return pattern matrix of lower triangle of *permuted* Hessian

gr(x)

Return gradient, evaluated at x: gr(x,...)

hessian(x)

Return sparse Hessian, evaluated at x, as a dgCMatrix object.

initialize( x, fn, gr, rows, cols, delta = 1e-07, index1 = TRUE, complex = FALSE, ... )

Initialize object with functions to compute the objective function and gradient (fn and gr), row and column indices of non-zero elements (rows and cols), an initial variable vector x at which fn and gr can be evaluated, a finite differencing parameter delta, flags for 0 or 1-based indexing (index1), whether the complex step method will be used, and other arguments (...) to be passed to fn and gr.

partition()

Return the partitioning used to compute finite differences

pointers(out.index1 = index1)

Return list with indices (idx) and pointers (pntr) for sparsity pattern of the compressed sparse Hessian. Since the Hessian is symmetric, the indices and pointers for row-oriented and column-oriented storage patterns are the same.

Examples

## Log posterior density of hierarchical binary choice model. See vignette.
set.seed(123)
data("binary_small")
N <- length(binary[["Y"]])
k <- NROW(binary[["X"]])
T <- binary[["T"]]
P <- rnorm((N+1)*k)
priors <- list(inv.Sigma = rWishart(1,k+5,diag(k))[,,1],
               inv.Omega = diag(k))
true.hess <- binary.hess(P, binary, priors)
pattern <- Matrix.to.Coord(Matrix::tril(true.hess))
str(pattern)
#> List of 2
#>  $ rows: int [1:1310] 1 2 3 4 201 202 203 204 2 3 ...
#>  $ cols: int [1:1310] 1 1 1 1 1 1 1 1 2 2 ...
obj <- sparseHessianFD(P, fn=binary.f, gr=binary.grad,
       rows=pattern[["rows"]], cols=pattern[["cols"]],
                      data=binary, priors=priors)
hs <- obj$hessian(P)
all.equal(hs, true.hess)
#> [1] TRUE

f <- obj$fn(P) ## obj function
df <- obj$gr(P) ## gradient
fdf <- obj$fngr(P) ## list of obj function and gradient
fdfhs <- obj$fngrhs(P) ## list of obj function, gradient and Hessian.