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 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261
|
###################################
# Matrix multiplication.
# Note that "%*%" methods should NOT call other matrix multiplication methods
# directly on their ResidualMatrix arguments. Rather, any ResidualMatrix should
# be broken down into the seed or the underlying matrix before further multiplication.
# This reduces the risk of infinite S4 recursion when 'y' is also an S4 matrix class.
#
# Specifically, the .*_ResidualMatrix functions take a seed object that IGNORES
# any non-FALSE setting of @transposed. It will then break down the seed into
# its constituents and perform the multiplication, such that any further S4
# dispatch occurs on the lower components of the seed.
#
# We also coerce each matrix product to a full matrix - which it usually is, anyway
# - to avoid the unnecessary overhead of multiplying DelayedArray instances.
#' @export
#' @importFrom Matrix t
#' @importFrom DelayedArray DelayedArray seed
setMethod("%*%", c("ResidualMatrix", "ANY"), function(x, y) {
x_seed <- seed(x)
if (is_transposed(x_seed)) {
out <- t(.leftmult_ResidualMatrix(t(y), x_seed))
} else {
out <- .rightmult_ResidualMatrix(x_seed, y)
}
DelayedArray(out)
})
#' @importFrom Matrix crossprod
.rightmult_ResidualMatrix <- function(x_seed, y) {
# Order of operations chosen to minimize size of intermediates,
# under the assumption that ncol(y) is very small.
as.matrix(get_matrix2(x_seed) %*% y) - get_Q(x_seed) %*% as.matrix(get_Qty(x_seed) %*% y)
}
#' @export
#' @importFrom Matrix t
#' @importFrom DelayedArray DelayedArray seed
setMethod("%*%", c("ANY", "ResidualMatrix"), function(x, y) {
y_seed <- seed(y)
if (is_transposed(y_seed)) {
if (!is.null(dim(x))) {
# Vectors don't quite behave as 1-column matrices here.
# so we need to be a bit more careful.
x <- t(x)
}
out <- t(.rightmult_ResidualMatrix(y_seed, x))
} else {
out <- .leftmult_ResidualMatrix(x, y_seed)
}
DelayedArray(out)
})
#' @importFrom Matrix tcrossprod
.leftmult_ResidualMatrix <- function(x, y_seed) {
# Order of operations chosen to minimize size of intermediates,
# under the assumption that nrow(x) is very small.
as.matrix(x %*% get_matrix2(y_seed)) - as.matrix(x %*% get_Q(y_seed)) %*% get_Qty(y_seed)
}
#' @export
#' @importFrom Matrix t
#' @importFrom DelayedArray DelayedArray seed
setMethod("%*%", c("ResidualMatrix", "ResidualMatrix"), function(x, y) {
x_seed <- seed(x)
if (is_transposed(x_seed)) {
out <- t(.leftmult_ResidualMatrix(t(y), x_seed))
} else {
out <- .rightmult_ResidualMatrix(x_seed, y)
}
DelayedArray(out)
})
# Unlike DeferredMatrix, we are happy to delegate to the left/right %*%
# for dual ResidualMatrix multiplication. This is because there are no
# operations that will collapse a ResidualMatrix instance to a DelayedMatrix
# prior to multiplication. Thus, we do not have to worry about writing
# specialized methods to avoid the overhead of DelayedMatrix multiplication.
###################################
# Crossproduct.
# Technically, this could be defined in terms of '%*%'.
# However, we re-define it in terms of 'crossprod' for efficiency,
# to exploit the potential availability of optimized crossprod for .matrix.
#' @export
#' @importFrom Matrix crossprod
#' @importFrom DelayedArray DelayedArray seed
setMethod("crossprod", c("ResidualMatrix", "missing"), function(x, y) {
x_seed <- seed(x)
if (is_transposed(x_seed)) {
# No need for t(), it's symmetric.
out <- .tcp_ResidualMatrix(x_seed)
} else {
out <- .crossprod_ResidualMatrix(x_seed)
}
DelayedArray(out)
})
#' @importFrom Matrix crossprod
.crossprod_ResidualMatrix <- function(x_seed) {
mat <- get_matrix2(x_seed)
Qty <- get_Qty(x_seed)
Q <- get_Q(x_seed)
# We assume that nrow(Q) >> ncol(Q) and nrow(mat) >> ncol(mat)
# in order for this to be efficient.
ytQQty <- crossprod(Qty)
QtQ <- crossprod(Q)
# Using this addition order to minimize numeric instability.
(crossprod(mat) - ytQQty) + (crossprod(Qty, QtQ %*% Qty) - ytQQty)
}
#' @export
#' @importFrom Matrix crossprod
#' @importFrom DelayedArray DelayedArray seed
setMethod("crossprod", c("ResidualMatrix", "ANY"), function(x, y) {
x_seed <- seed(x)
if (is_transposed(x_seed)) {
out <- .rightmult_ResidualMatrix(x_seed, y)
} else {
out <- .rightcross_ResidualMatrix(x_seed, y)
}
DelayedArray(out)
})
#' @importFrom Matrix crossprod
.rightcross_ResidualMatrix <- function(x_seed, y) {
# Order of operations chosen to minimize size of intermediates,
# under the assumption that ncol(y) is very small.
as.matrix(crossprod(get_matrix2(x_seed), y)) - crossprod(get_Qty(x_seed), as.matrix(crossprod(get_Q(x_seed), y)))
}
#' @export
#' @importFrom Matrix crossprod t
#' @importFrom DelayedArray DelayedArray seed
setMethod("crossprod", c("ANY", "ResidualMatrix"), function(x, y) {
y_seed <- seed(y)
if (is_transposed(y_seed)) {
out <- t(.rightmult_ResidualMatrix(y_seed, x))
} else {
out <- .leftcross_ResidualMatrix(x, y_seed)
}
DelayedArray(out)
})
#' @importFrom Matrix crossprod
.leftcross_ResidualMatrix <- function(x, y_seed) {
# Order of operations chosen to minimize size of intermediates,
# under the assumption that ncol(x) is very small.
as.matrix(crossprod(x, get_matrix2(y_seed))) - as.matrix(crossprod(x, get_Q(y_seed))) %*% get_Qty(y_seed)
}
#' @export
#' @importFrom Matrix crossprod
#' @importFrom DelayedArray DelayedArray seed
setMethod("crossprod", c("ResidualMatrix", "ResidualMatrix"), function(x, y) {
x_seed <- seed(x)
if (is_transposed(x_seed)) {
out <- .rightmult_ResidualMatrix(x_seed, y)
} else {
out <- .rightcross_ResidualMatrix(x_seed, y)
}
DelayedArray(out)
})
###################################
# Transposed crossproduct.
# Technically, this could be defined in terms of '%*%'.
# However, we re-define it in terms of 'tcrossprod' for efficiency,
# to exploit the potential availability of optimized crossprod for .matrix.
#' @export
#' @importFrom Matrix tcrossprod
#' @importFrom DelayedArray DelayedArray seed
setMethod("tcrossprod", c("ResidualMatrix", "missing"), function(x, y) {
x_seed <- seed(x)
if (is_transposed(x_seed)) {
out <- .crossprod_ResidualMatrix(x_seed)
} else {
out <- .tcp_ResidualMatrix(x_seed)
}
DelayedArray(out)
})
#' @importFrom Matrix tcrossprod
.tcp_ResidualMatrix <- function(x_seed) {
mat <- get_matrix2(x_seed)
Qty <- get_Qty(x_seed)
Q <- get_Q(x_seed)
# We assume that ncol(mat) >> nrow(mat) for this to be efficient.
# We also try to avoid constructing QQt under the assumption that
# nrow(Q) >> ncol(Q) for Q derived from a full-rank design matrix.
QQtyyt<- Q %*% tcrossprod(Qty, mat)
QQtyytQQt <- tcrossprod(QQtyyt %*% Q, Q)
# Using this addition order to minimize numeric instability.
(tcrossprod(mat) - QQtyyt) + (QQtyytQQt - t(QQtyyt))
}
#' @export
#' @importFrom Matrix tcrossprod t
#' @importFrom DelayedArray DelayedArray seed
setMethod("tcrossprod", c("ResidualMatrix", "ANY"), function(x, y) {
if (is.null(dim(y))) { # for consistency with base::tcrossprod.
stop("non-conformable arguments")
}
x_seed <- seed(x)
if (is_transposed(x_seed)) {
out <- t(.leftmult_ResidualMatrix(y, x_seed))
} else {
out <- .righttcp_ResidualMatrix(x_seed, y)
}
DelayedArray(out)
})
#' @importFrom Matrix tcrossprod
.righttcp_ResidualMatrix <- function(x_seed, y) {
# Order of operations chosen to minimize size of intermediates,
# assuming that nrow(y) is very small.
as.matrix(tcrossprod(get_matrix2(x_seed), y)) - get_Q(x_seed) %*% as.matrix(tcrossprod(get_Qty(x_seed), y))
}
#' @export
#' @importFrom Matrix tcrossprod
#' @importFrom DelayedArray DelayedArray seed
setMethod("tcrossprod", c("ANY", "ResidualMatrix"), function(x, y) {
y_seed <- seed(y)
if (is_transposed(y_seed)) {
out <- .leftmult_ResidualMatrix(x, y_seed)
} else {
out <- .lefttcp_ResidualMatrix(x, y_seed)
}
DelayedArray(out)
})
#' @importFrom Matrix tcrossprod
.lefttcp_ResidualMatrix <- function(x, y_seed) {
# Order of operations chosen to minimize size of intermediates.
# assuming that nrow(x) is very small.
as.matrix(tcrossprod(x, get_matrix2(y_seed))) - tcrossprod(as.matrix(tcrossprod(x, get_Qty(y_seed))), get_Q(y_seed))
}
#' @export
#' @importFrom Matrix tcrossprod t
#' @importFrom DelayedArray DelayedArray seed
setMethod("tcrossprod", c("ResidualMatrix", "ResidualMatrix"), function(x, y) {
x_seed <- seed(x)
if (is_transposed(x_seed)) {
out <- t(.leftmult_ResidualMatrix(y, x_seed))
} else {
out <- .righttcp_ResidualMatrix(x_seed, y)
}
DelayedArray(out)
})
|