sparseQR-class {Matrix}R Documentation

Sparse QR decomposition of a sparse matrix

Description

Objects class "sparseQR" represent a QR decomposition of a sparse m x n (“long”: m >= n) rectangular matrix A, typically resulting from qr(), see ‘Details’ notably about row and column permutations for pivoting.

Details

For a sparse m x n (“long”: m >= n) rectangular matrix A, the sparse QR decomposition is either
of the form P A = Q R with a (row) permutation matrix P, (encoded in the p slot of the result) if the q slot is of length 0,
or of the form P A P* = Q R with an extra (column) permutation matrix P* (encoded in the q slot). Note that the row permutation P A in R is simply A[p+1, ] where p is the p-slot, a 0-based permutation of 1:m applied to the rows of the original matrix.

If the q slot has length n it is a 0-based permutation of 1:n applied to the columns of the original matrix to reduce the amount of “fill-in” in the matrix R, and A P* in R is simply A[ , q+1].

R is an m by n matrix that is zero below the main diagonal, i.e., upper triangular (m by n) with m-n extra zero rows.

The matrix Q is a "virtual matrix". It is the product of n Householder transformations. The information to generate these Householder transformations is stored in the V and beta slots.
Note however that qr.Q() returns the row permuted matrix Q* := P^(-1) Q = P'Q as permutation matrices are orthogonal; and Q* is orthogonal itself because Q and P are. This is useful because then, as in the dense matrix and base R matrix qr case, we have the mathematical identity

P A = Q* R,

in R as

            A[p+1,] == qr.Q(*) %*% R .

The "sparseQR" methods for the qr.* functions return objects of class "dgeMatrix" (see dgeMatrix). Results from qr.coef, qr.resid and qr.fitted (when k == ncol(R)) are well-defined and should match those from the corresponding dense matrix calculations. However, because the matrix Q is not uniquely defined, the results of qr.qy and qr.qty do not necessarily match those from the corresponding dense matrix calculations.

Also, the results of qr.qy and qr.qty apply to the permuted column order when the q slot has length n.

Objects from the Class

Objects can be created by calls of the form new("sparseQR", ...) but are more commonly created by function qr applied to a sparse matrix such as a matrix of class dgCMatrix.

Slots

V:

Object of class "dgCMatrix". The columns of V are the vectors that generate the Householder transformations of which the matrix Q is composed.

beta:

Object of class "numeric", the normalizing factors for the Householder transformations.

p:

Object of class "integer": Permutation (of 0:(n-1)) applied to the rows of the original matrix.

R:

Object of class "dgCMatrix": An upper triangular matrix of the same dimension as X.

q:

Object of class "integer": Permutation applied from the right, i.e., to the columns of the original matrix. Can be of length 0 which implies no permutation.

Methods

qr.R

signature(qr = "sparseQR"): compute the upper triangular R matrix of the QR decomposition. Note that this currently warns because of possible permutation mismatch with the classical qr.R() result, and you can suppress these warnings by setting options() either "Matrix.quiet.qr.R" or (the more general) either "Matrix.quiet" to TRUE.

qr.Q

signature(qr = "sparseQR"): compute the orthogonal Q matrix of the QR decomposition.

qr.coef

signature(qr = "sparseQR", y = "ddenseMatrix"): ...

qr.coef

signature(qr = "sparseQR", y = "matrix"): ...

qr.coef

signature(qr = "sparseQR", y = "numeric"): ...

qr.fitted

signature(qr = "sparseQR", y = "ddenseMatrix"): ...

qr.fitted

signature(qr = "sparseQR", y = "matrix"): ...

qr.fitted

signature(qr = "sparseQR", y = "numeric"): ...

qr.qty

signature(qr = "sparseQR", y = "ddenseMatrix"): ...

qr.qty

signature(qr = "sparseQR", y = "matrix"): ...

qr.qty

signature(qr = "sparseQR", y = "numeric"): ...

qr.qy

signature(qr = "sparseQR", y = "ddenseMatrix"): ...

qr.qy

signature(qr = "sparseQR", y = "matrix"): ...

qr.qy

signature(qr = "sparseQR", y = "numeric"): ...

qr.resid

signature(qr = "sparseQR", y = "ddenseMatrix"): ...

qr.resid

signature(qr = "sparseQR", y = "matrix"): ...

qr.resid

signature(qr = "sparseQR", y = "numeric"): ...

solve

signature(a = "sparseQR", b = "ANY"): For solve(a,b), simply uses qr.coef(a,b).

See Also

qr, qr.Q, qr.R, qr.fitted, qr.resid, qr.coef, qr.qty, qr.qy,

Permutation matrices in the Matrix package: pMatrix; dgCMatrix, dgeMatrix.

Examples

data(KNex)
mm <- KNex $ mm
 y <- KNex $  y
 y. <- as(as.matrix(y), "dgCMatrix")
str(qrm <- qr(mm))
 qc  <- qr.coef  (qrm, y); qc. <- qr.coef  (qrm, y.) # 2nd failed in Matrix <= 1.1-0
 qf  <- qr.fitted(qrm, y); qf. <- qr.fitted(qrm, y.)
 qs  <- qr.resid (qrm, y); qs. <- qr.resid (qrm, y.)
stopifnot(all.equal(qc, as.numeric(qc.),  tolerance=1e-12),
          all.equal(qf, as.numeric(qf.),  tolerance=1e-12),
          all.equal(qs, as.numeric(qs.),  tolerance=1e-12),
          all.equal(qf+qs, y, tolerance=1e-12))


[Package Matrix version 1.2-17 Index]