1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154
|
#### Methods for the sparse QR decomposition
## TODO: qr.R() generic that allows optional args ['backPermute']
## --- so we can add it to our qr.R() method, *instead* of this :
qrR <- function(qr, complete = FALSE, backPermute = TRUE, row.names = TRUE) {
ir <- seq_len(qr@Dim[if(complete) 1L else 2L])
r <- if(backPermute <- backPermute && (n <- length(qr@q)) && !isSeq(qr@q, n-1L))
qr@R[ir, order(qr@q), drop = FALSE] else
qr@R[ir, , drop = FALSE]
if(row.names && !is.null(rn <- qr@V@Dimnames[[1]])) # qr.R() in 'base' gives rownames
r@Dimnames[[1]] <- rn[seq_len(r@Dim[1L])]
if(complete || backPermute) r else as(r, "triangularMatrix")
}
setMethod("qr.R", signature(qr = "sparseQR"),
function(qr, complete = FALSE) {
if(nonTRUEoption("Matrix.quiet.qr.R") && nonTRUEoption("Matrix.quiet"))
warning("qr.R(<sparse>) may differ from qr.R(<dense>) because of permutations. Possibly use our qrR() instead")
qrR(qr, complete=complete, backPermute=FALSE)
})
## if(identical("", as.character(formals(qr.Q)$Dvec))) { # "new"
setMethod("qr.Q", "sparseQR",
function(qr, complete=FALSE, Dvec)
{
d <- qr@Dim
## ir <- seq_len(d[k <- if(complete) 1L else 2L])
k <- if(complete) 1L else 2L
if(missing(Dvec)) Dvec <- rep.int(1, if(complete) d[1] else min(d))
D <- .sparseDiagonal(d[1], x=Dvec, cols=0L:(d[k] -1L))
qr.qy(qr, D)
})
## } else {
## setMethod("qr.Q", "sparseQR",
## function(qr, complete=FALSE, Dvec = rep.int(1, if(complete) d[1] else min(d)))
## {
## d <- qr@Dim
## ir <- seq_len(d[k <- if(complete) 1L else 2L])
## D <- .sparseDiagonal(d[1], x=Dvec, cols=0L:(d[k] -1L))
## qr.qy(qr, D)
## })
## }
## NB: Here, the .Call()s to sparseQR_qty all set keep_names = TRUE
## --- instead of allowing it to become an argument,
## mainly because the base functions qr.qy() / qr.qty() have no '...' formal argument
## To change, would make these *implicit* generics in 'methods' - as qr.R
## Also, qr() itself has keep.names = TRUE/FALSE -- should be enough
setMethod("qr.qy", signature(qr = "sparseQR", y = "ddenseMatrix"),
function(qr, y) .Call(sparseQR_qty, qr, y, FALSE, TRUE),
valueClass = "dgeMatrix")
setMethod("qr.qy", signature(qr = "sparseQR", y = "matrix"),
function(qr, y) .Call(sparseQR_qty, qr, y, FALSE, TRUE),
valueClass = "dgeMatrix")
setMethod("qr.qy", signature(qr = "sparseQR", y = "numeric"),
## drop to vector {to be 100% standard-R-matrix compatible} :
function(qr, y) .Call(sparseQR_qty, qr, y, FALSE, TRUE)@x)
setMethod("qr.qy", signature(qr = "sparseQR", y = "Matrix"),
function(qr, y) .Call(sparseQR_qty, qr,
as(as(y, "denseMatrix"),"dgeMatrix"), FALSE, TRUE),
valueClass = "dgeMatrix")
setMethod("qr.qty", signature(qr = "sparseQR", y = "ddenseMatrix"),
function(qr, y) .Call(sparseQR_qty, qr, y, TRUE, TRUE),
valueClass = "dgeMatrix")
setMethod("qr.qty", signature(qr = "sparseQR", y = "matrix"),
function(qr, y) .Call(sparseQR_qty, qr, y, TRUE, TRUE),
valueClass = "dgeMatrix")
setMethod("qr.qty", signature(qr = "sparseQR", y = "numeric"),
function(qr, y) .Call(sparseQR_qty, qr, y, TRUE, TRUE)@x)
setMethod("qr.qty", signature(qr = "sparseQR", y = "Matrix"),
function(qr, y) .Call(sparseQR_qty, qr,
as(as(y, "denseMatrix"),"dgeMatrix"), TRUE, TRUE),
valueClass = "dgeMatrix")
## FIXME: really should happen in C, i.e sparseQR_coef() in ../src/sparseQR.c :
.coef.trunc <- function(qr, res, drop=FALSE) {
if(!all((d <- lengths(res@Dimnames)) == 0L) && !identical(d, D <- res@Dim)) {
## Fix dimnames from dim (when not NULL !) :
if(d[[1]]) length(res@Dimnames[[1]]) <- D[[1]]
if(d[[2]]) length(res@Dimnames[[2]]) <- D[[2]]
}
res[seq_len(ncol(qr@R)),,drop=drop]
}
setMethod("qr.coef", signature(qr = "sparseQR", y = "ddenseMatrix"),
function(qr, y)
.coef.trunc(qr, .Call(sparseQR_coef, qr, y)),
valueClass = "dgeMatrix")
setMethod("qr.coef", signature(qr = "sparseQR", y = "matrix"),
function(qr, y)
.coef.trunc(qr, .Call(sparseQR_coef, qr, y)),
valueClass = "dgeMatrix")
setMethod("qr.coef", signature(qr = "sparseQR", y = "numeric"),
function(qr, y)
.coef.trunc(qr, .Call(sparseQR_coef, qr, y), drop=TRUE))
setMethod("qr.coef", signature(qr = "sparseQR", y = "Matrix"),
function(qr, y)
.coef.trunc(qr, .Call(sparseQR_coef, qr,
as(as(y, "denseMatrix"),"dgeMatrix"))),
valueClass = "dgeMatrix")
## qr.resid() & qr.fitted() : ---------------------------
setMethod("qr.resid", signature(qr = "sparseQR", y = "ddenseMatrix"),
function(qr, y)
.Call(sparseQR_resid_fitted, qr, y, TRUE),
valueClass = "dgeMatrix")
setMethod("qr.resid", signature(qr = "sparseQR", y = "matrix"),
function(qr, y)
.Call(sparseQR_resid_fitted, qr, y, TRUE),
valueClass = "dgeMatrix")
setMethod("qr.resid", signature(qr = "sparseQR", y = "numeric"),
function(qr, y) drop(.Call(sparseQR_resid_fitted, qr, y, TRUE)))
setMethod("qr.resid", signature(qr = "sparseQR", y = "Matrix"),
function(qr, y)
.Call(sparseQR_resid_fitted, qr,
as(as(y, "denseMatrix"),"dgeMatrix"), TRUE),
valueClass = "dgeMatrix")
setMethod("qr.fitted", signature(qr = "sparseQR", y = "ddenseMatrix"),
function(qr, y, k)
.Call(sparseQR_resid_fitted, qr, y, FALSE),
valueClass = "dgeMatrix")
setMethod("qr.fitted", signature(qr = "sparseQR", y = "matrix"),
function(qr, y, k)
.Call(sparseQR_resid_fitted, qr, y, FALSE),
valueClass = "dgeMatrix")
setMethod("qr.fitted", signature(qr = "sparseQR", y = "numeric"),
function(qr, y, k) drop(.Call(sparseQR_resid_fitted, qr, y, FALSE)))
setMethod("qr.fitted", signature(qr = "sparseQR", y = "Matrix"),
function(qr, y, k)
.Call(sparseQR_resid_fitted, qr,
as(as(y, "denseMatrix"),"dgeMatrix"), FALSE),
valueClass = "dgeMatrix")
##
setMethod("solve", signature(a = "sparseQR", b = "ANY"),
function(a, b, ...) qr.coef(a, b))
|