Takahashi_Davis {sparseinv} | R Documentation |
Takahashi equations
Description
Computes the sparse inverse subset of a sparse matrix Q
using the Takahashi equations.
Usage
Takahashi_Davis(Q = NULL, cholQp = NULL, return_perm_chol = 0, P = 0,
gc = 0)
Arguments
Q |
precision matrix of class |
cholQp |
the Cholesky factor of class |
return_perm_chol |
if 1, the Cholesky factor of the permuted |
P |
the permutation matrix of class |
gc |
do garbage collection throughout (may increase computational time but useful for small memory machines) |
Details
This function first computes the Cholesky factor of Q
. The fill-in reduction permutation is the approximate minimum degree permutation (amd) of Timothy Davis' SuiteSparse package configured to be slightly more aggressive than that in the Matrix
package. The function then uses the Takahashi equations to compute the variances at the non-zero locations in the Cholesky factor from the factor itself. The equations themselves are implemented in C using the SparseSuite package of Timothy Davis.
Value
if return_perm_chol == 0, the sparse inverse subset of Q is returned, where the non-zero elements correspond to those in the Cholesky factor of its permutation.
If !(return_perm_chol == 0), a list with three elements is returned: S
(the sparse inverse subset), Lp (the Cholesky factor of the permuted matrix) and P (the
permutation matrix)
Note
This package is a wrapper for C functions implemented by Timothy Davis in SuiteSparse. The author of this package has done no work on the sparse inverse routines themselves and any acknowledgment should include one to SuiteSparse (see below for reference). The author of this package was made aware of this methodology by Botond Cseke.
References
Takahashi, K., Fagan, J., Chin, M.-S., 1973. Formation of a sparse bus impedance matrix and its application to short circuit study. 8th PICA Conf. Proc. June 4–6, Minneapolis, Minn.
Davis, T. A., 2014. sparseinv: Sparse Inverse Subset. URL https://au.mathworks.com/matlabcentral/fileexchange/33966-sparseinv–sparse-inverse-subset Davis, T. A., 2006. Direct Methods for Sparse Linear Systems. SIAM, Philadelphia, PA.
Examples
require(Matrix)
Q = sparseMatrix(i = c(1, 1, 2, 2),
j = c(1, 2, 1, 2),
x = c(0.1, 0.2, 0.2, 1))
X <- cholPermute(Q)
S_partial = Takahashi_Davis(Q, cholQp = X$Qpermchol, P = X$P)