lu {Matrix} | R Documentation |
Computes (generalized) triangular decompositions of square (sparse or dense) and non-square dense matrices.
lu(x, ...) ## S4 method for signature 'matrix' lu(x, warnSing = TRUE, ...) ## S4 method for signature 'dgeMatrix' lu(x, warnSing = TRUE, ...) ## S4 method for signature 'dgCMatrix' lu(x, errSing = TRUE, order = TRUE, tol = 1, keep.dimnames = TRUE, ...)
x |
a dense or sparse matrix, in the latter case of square dimension. No missing values or IEEE special values are allowed. |
warnSing |
(when |
errSing |
(when |
order |
logical or integer, used to choose which fill-reducing permutation technique will be used internally. Do not change unless you know what you are doing. |
tol |
positive number indicating the pivoting tolerance used in
|
keep.dimnames |
logical indicating that |
... |
further arguments passed to or from other methods. |
lu()
is a generic function with special methods for different types
of matrices. Use showMethods("lu")
to list all the methods
for the lu
generic.
The method for class dgeMatrix
(and all dense
matrices) is based on LAPACK's "dgetrf"
subroutine. It returns
a decomposition also for singular and non-square matrices.
The method for class dgCMatrix
(and all sparse
matrices) is based on functions from the CSparse library. It signals
an error (or returns NA
, when errSing = FALSE
, see
above) when the decomposition algorithm fails, as when x
is
(too close to) singular.
An object of class "LU"
, i.e., "denseLU"
(see its separate help page),
or "sparseLU"
, see sparseLU
; this is
a representation of a triangular decomposition of x
.
Because the underlying algorithm differ entirely,
in the dense case (class denseLU
), the
decomposition is
A = P L U,
where as in the sparse case (class
sparseLU
), it is
A = P' L U Q.
Golub, G., and Van Loan, C. F. (1989). Matrix Computations, 2nd edition, Johns Hopkins, Baltimore.
Timothy A. Davis (2006) Direct Methods for Sparse Linear Systems, SIAM Series “Fundamentals of Algorithms”.
Class definitions denseLU
and sparseLU
and function expand
;
qr
, chol
.
##--- Dense ------------------------- x <- Matrix(rnorm(9), 3, 3) lu(x) dim(x2 <- round(10 * x[,-3]))# non-square expand(lu2 <- lu(x2)) ##--- Sparse (see more in ?"sparseLU-class")----- % ./sparseLU-class.Rd pm <- as(readMM(system.file("external/pores_1.mtx", package = "Matrix")), "CsparseMatrix") str(pmLU <- lu(pm)) # p is a 0-based permutation of the rows # q is a 0-based permutation of the columns ## permute rows and columns of original matrix ppm <- pm[pmLU@p + 1L, pmLU@q + 1L] pLU <- drop0(pmLU@L %*% pmLU@U) # L %*% U -- dropping extra zeros ## equal up to "rounding" ppm[1:14, 1:5] pLU[1:14, 1:5]