sparse_chol_crs {SparseChol} | R Documentation |
Sparse Cholesky decomposition with sparse representation
Description
Sparse Cholesky decomposition with sparse representation
Usage
sparse_chol_crs(n, Ap, Ai, Ax)
Arguments
n |
Integer specifying the dimension of the matrix |
Ap |
numeric (integer valued) vector of pointers, one for each column (or row), to the initial (zero-based) index of elements in the column (or row). |
Ai |
Integer vector specifying the row positions of the non-zero values of the matrix |
Ax |
values of the non-zero matrix entries |
Details
Generates the LDL decomposition of a symmetric, sparse matrix using the method described by Timothy Davis (see references). Required input is a matrix in sparse format from the matrix package, see sparseMatrix, or the package function dense_to_sparse. To instead use a matrix directly, see sparse_chol.
Value
A list with elements n, Ai, Ap, Ax (corresponding to above arguments) for matrix L, and element D, which contains the diagonal values of matrix D.
Examples
n <- 10
Ap <- c(0, 1, 2, 3, 4, 6, 7, 9, 11, 15, 19)
Ai <- c(1, 2, 3, 4, 2,5, 6, 5,7, 5,8, 1,5,8,9, 2,5,7,10)
Ax <- c(1.7, 1., 1.5, 1.1, .02,2.6, 1.2, .16,1.3, .09,1.6,
.13,.52,.11,1.4, .01,.53,.56,3.1)
out <-sparse_chol_crs(n,Ap,Ai,Ax)
sparse_L(out)
sparse_D(out)