SQUIC: Sparse QUadratic approximate Inverse Covariance matrix...

Description Usage Arguments Details Value References Examples

View source: R/SQUIC.R

Description

SQUIC is a second-order, L1-regularized maximum likelihood method for performant large-scale sparse precision matrix estimation.

Usage

1
2
 SQUIC(Y, lambda, max_iter = 100, tol = 1e-3, verbose = 1, M = NULL,  X0 = NULL, W0 = NULL)
 

Arguments

Y

Input data in the form p (dimensions) by n (samples).

lambda

Scalar tuning parameter controlling the sparsity of the precision matrix (a non-zero positive number).

max_iter

Maximum number of Newton iterations of the outer loop. Default: 100.

tol

Tolerance for convergence and approximate inversion. Default: 1e-3.

verbose

Level of printing output (0 or 1). Default: 1.

M

The matrix encoding the sparsity pattern of the matrix tuning parameter, i.e., Lambda (p by p). Default: NULL.

X0

Initial guess of the precision matrix (p by p). Default: Identity.

W0

Initial guess of the inverse of the precision matrix (p by p). Default: Identity.

Details

SQUIC is a high-performance precision matrix estimation package intended for large-scale applications. The algorithm utilizes second-order information in conjunction with a highly optimized, supernodal Cholesky factorization routine (CHOLMOD), a block-oriented computational approach for both the coordinate descent and the approximate matrix inversion, and finally, a sparse representation of the sample covariance matrix. The library is implemented in C++ with shared-memory parallelization using OpenMP. For further details, please see the listed references.

NOTE: (1) If max_iter=0, the returned value for the inverse of the precision matrix (W) is the sparse sample covariance matrix (S). (2) The number of threads used by SQUIC can be defined by setting the enviroment variable OMP_NUM_THREADS (e.g., base> export OMP_NUM_THREADS=12). This may require a restart of the session).

Value

A list with components

X

Estimated precision matrix (p by p).

W

Estimated inverse of the precision matrix (p by p).

info_time_total

Total runtime.

info_time_sample_cov

Runtime of the sample covariance matrix.

info_time_optimize

Runtime of the Newton steps.

info_time_factor

Runtime of the Cholesky factorization.

info_time_approximate_inv

Runtime of the approximate matrix inversion.

info_time_coordinate_upd

Runtime of the coordinate descent update.

info_objective

Value of the objective at each Newton iteration.

info_logdetX

Value of the log determinant of the precision matrix.

info_trSX

Value of the trace of the sample covariance matrix times the precision matrix.

References

Bollhöfer, M., Eftekhari, A., Scheidegger, S. and Schenk, O., 2019. Large-Scale Sparse Inverse Covariance Matrix Estimation. SIAM Journal on Scientific Computing, 41(1), pp.A380-A401.

Eftekhari, A., Bollhöfer, M. and Schenk, O., 2018, November. Distributed Memory Sparse Inverse Covariance Matrix Estimation on High-performance Computing Architectures. In SC18: International Conference for High Performance Computing, Networking, Storage and Analysis (pp. 253-264). IEEE.

Eftekhari, A., Pasadakis, D., Bollhöfer, M., Scheidegger, S. and Schenk, O., 2021. Block-Enhanced Precision Matrix Estimation for Large-Scale Datasets. Journal of Computational Science, p. 101389.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
 library(SQUIC)

 p = 1024
 n = 100
 lambda = .4
 max_iter = 100
 tol = 1e-3

 # generate a tridiagonal matrix
 iC_star = Matrix::bandSparse(p, p, (-1):1, list(rep(-.5, p-1), rep(1.25, p), rep(-.5, p-1)));

 # generate the data
 z    = replicate(n,rnorm(p));
 iC_L = chol(iC_star);
 data = matrix(solve(iC_L,z),p,n);

 # Run SQUIC
 out <- SQUIC(data,lambda)
 

aefty/SQUIC_R documentation built on July 29, 2021, 5:17 p.m.