colSums {Matrix} | R Documentation |
Form row and column sums and means for
objects, for sparseMatrix
the result may
optionally be sparse (sparseVector
), too.
Row or column names are kept respectively as for base matrices
and colSums
methods, when the result is
numeric
vector.
colSums (x, na.rm = FALSE, dims = 1, ...) rowSums (x, na.rm = FALSE, dims = 1, ...) colMeans(x, na.rm = FALSE, dims = 1, ...) rowMeans(x, na.rm = FALSE, dims = 1, ...) ## S4 method for signature 'CsparseMatrix' colSums(x, na.rm = FALSE, dims = 1, sparseResult = FALSE) ## S4 method for signature 'CsparseMatrix' rowSums(x, na.rm = FALSE, dims = 1, sparseResult = FALSE) ## S4 method for signature 'CsparseMatrix' colMeans(x, na.rm = FALSE, dims = 1, sparseResult = FALSE) ## S4 method for signature 'CsparseMatrix' rowMeans(x, na.rm = FALSE, dims = 1, sparseResult = FALSE)
x |
a Matrix, i.e., inheriting from |
na.rm |
logical. Should missing values (including |
dims |
completely ignored by the |
... |
potentially further arguments, for method |
sparseResult |
logical indicating if the result should be sparse,
i.e., inheriting from class |
returns a numeric vector if sparseResult
is FALSE
as per
default. Otherwise, returns a sparseVector
.
dimnames(x)
are only kept (as names(v)
)
when the resulting v
is numeric
, since
sparseVector
s do not have names.
colSums
and the
sparseVector
classes.
(M <- bdiag(Diagonal(2), matrix(1:3, 3,4), diag(3:2))) # 7 x 8 colSums(M) d <- Diagonal(10, c(0,0,10,0,2,rep(0,5))) MM <- kronecker(d, M) dim(MM) # 70 80 length(MM@x) # 160, but many are '0' ; drop those: MM <- drop0(MM) length(MM@x) # 32 cm <- colSums(MM) (scm <- colSums(MM, sparseResult = TRUE)) stopifnot(is(scm, "sparseVector"), identical(cm, as.numeric(scm))) rowSums (MM, sparseResult = TRUE) # 14 of 70 are not zero colMeans(MM, sparseResult = TRUE) # 16 of 80 are not zero ## Since we have no 'NA's, these two are equivalent : stopifnot(identical(rowMeans(MM, sparseResult = TRUE), rowMeans(MM, sparseResult = TRUE, na.rm = TRUE)), rowMeans(Diagonal(16)) == 1/16, colSums(Diagonal(7)) == 1) ## dimnames(x) --> names( <value> ) : dimnames(M) <- list(paste0("r", 1:7), paste0("V",1:8)) M colSums(M) rowMeans(M) ## Assertions : stopifnot(all.equal(colSums(M), setNames(c(1,1,6,6,6,6,3,2), colnames(M))), all.equal(rowMeans(M), structure(c(1,1,4,8,12,3,2) / 8, .Names = paste0("r", 1:7))))