bandSparse.RdConstruct a sparse banded matrix by specifying its non-zero sup- and super-diagonals.
bandSparse(n, m = n, k, diagonals, symmetric = FALSE,
repr = "C", giveCsparse = (repr == "C"))the matrix dimension \((n,m) = (nrow, ncol)\).
integer vector of “diagonal numbers”, with identical
meaning as in band(*, k), i.e., relative to the main diagonal,
which is k=0.
optional list of sub-/super- diagonals; if missing,
the result will be a pattern matrix, i.e., inheriting from
class nMatrix.
diagonals can also be \(n' \times d\) matrix, where
d <- length(k) and \(n' >= min(n,m)\). In that case,
the sub-/super- diagonals are taken from the columns of
diagonals, where only the first several rows will be used
(typically) for off-diagonals.
logical; if true the result will be symmetric
(inheriting from class symmetricMatrix) and
only the upper or lower triangle must be specified (via k and
diagonals).
character string, one of "C",
"T", or "R", specifying the sparse representation to
be used for the result, i.e., one from the super classes
CsparseMatrix, TsparseMatrix, or
RsparseMatrix.
(deprecated, replaced with repr):
logical indicating if the result should be a
CsparseMatrix or a
TsparseMatrix, where the default was TRUE,
and now is determined from repr; very often Csparse matrices are
more efficient subsequently, but not always.
a sparse matrix (of class
CsparseMatrix) of dimension \(n \times m\)
with diagonal “bands” as specified.
band, for extraction of matrix bands;
bdiag, diag,
sparseMatrix,
Matrix.
diags <- list(1:30, 10*(1:20), 100*(1:20))
s1 <- bandSparse(13, k = -c(0:2, 6), diag = c(diags, diags[2]), symm=TRUE)
s1
#> 13 x 13 sparse Matrix of class "dsCMatrix"
#>
#> [1,] 1 10 100 . . . 10 . . . . . .
#> [2,] 10 2 20 200 . . . 20 . . . . .
#> [3,] 100 20 3 30 300 . . . 30 . . . .
#> [4,] . 200 30 4 40 400 . . . 40 . . .
#> [5,] . . 300 40 5 50 500 . . . 50 . .
#> [6,] . . . 400 50 6 60 600 . . . 60 .
#> [7,] 10 . . . 500 60 7 70 700 . . . 70
#> [8,] . 20 . . . 600 70 8 80 800 . . .
#> [9,] . . 30 . . . 700 80 9 90 900 . .
#> [10,] . . . 40 . . . 800 90 10 100 1000 .
#> [11,] . . . . 50 . . . 900 100 11 110 1100
#> [12,] . . . . . 60 . . . 1000 110 12 120
#> [13,] . . . . . . 70 . . . 1100 120 13
s2 <- bandSparse(13, k = c(0:2, 6), diag = c(diags, diags[2]), symm=TRUE)
stopifnot(identical(s1, t(s2)), is(s1,"dsCMatrix"))
## a pattern Matrix of *full* (sub-)diagonals:
bk <- c(0:4, 7,9)
(s3 <- bandSparse(30, k = bk, symm = TRUE))
#> 30 x 30 sparse Matrix of class "nsCMatrix"
#>
#> [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,] . . . . . . . . . . . . . . . . . . . . | . | . . | | | | |
## If you want a pattern matrix, but with "sparse"-diagonals,
## you currently need to go via logical sparse:
lLis <- lapply(list(rpois(20, 2), rpois(20, 1), rpois(20, 3))[c(1:3, 2:3, 3:2)],
as.logical)
(s4 <- bandSparse(20, k = bk, symm = TRUE, diag = lLis))
#> 20 x 20 sparse Matrix of class "lsCMatrix"
#>
#> [1,] | : | : | . . | . : . . . . . . . . . .
#> [2,] : | : | : | . . | . : . . . . . . . . .
#> [3,] | : | | | | | . . | . | . . . . . . . .
#> [4,] : | | | | : | : . . : . | . . . . . . .
#> [5,] | : | | | | | | | . . | . | . . . . . .
#> [6,] . | | : | | | | | | . . | . | . . . . .
#> [7,] . . | | | | | | | | | . . | . | . . . .
#> [8,] | . . : | | | | | | | | . . | . | . . .
#> [9,] . | . . | | | | | | | | | . . | . | . .
#> [10,] : . | . . | | | | | | | | | . . | . | .
#> [11,] . : . : . . | | | | | | | | | . . | . |
#> [12,] . . | . | . . | | | | | | | | | . . | .
#> [13,] . . . | . | . . | | | | | | | | | . . |
#> [14,] . . . . | . | . . | | | | | : | : | . .
#> [15,] . . . . . | . | . . | | | : | | | | | .
#> [16,] . . . . . . | . | . . | | | | : : | : |
#> [17,] . . . . . . . | . | . . | : | : | | | |
#> [18,] . . . . . . . . | . | . . | | | | : | |
#> [19,] . . . . . . . . . | . | . . | : | | | |
#> [20,] . . . . . . . . . . | . | . . | | | | |
(s4. <- as(drop0(s4), "nsparseMatrix"))
#> 20 x 20 sparse Matrix of class "nsCMatrix"
#>
#> [1,] | . | . | . . | . . . . . . . . . . . .
#> [2,] . | . | . | . . | . . . . . . . . . . .
#> [3,] | . | | | | | . . | . | . . . . . . . .
#> [4,] . | | | | . | . . . . . | . . . . . . .
#> [5,] | . | | | | | | | . . | . | . . . . . .
#> [6,] . | | . | | | | | | . . | . | . . . . .
#> [7,] . . | | | | | | | | | . . | . | . . . .
#> [8,] | . . . | | | | | | | | . . | . | . . .
#> [9,] . | . . | | | | | | | | | . . | . | . .
#> [10,] . . | . . | | | | | | | | | . . | . | .
#> [11,] . . . . . . | | | | | | | | | . . | . |
#> [12,] . . | . | . . | | | | | | | | | . . | .
#> [13,] . . . | . | . . | | | | | | | | | . . |
#> [14,] . . . . | . | . . | | | | | . | . | . .
#> [15,] . . . . . | . | . . | | | . | | | | | .
#> [16,] . . . . . . | . | . . | | | | . . | . |
#> [17,] . . . . . . . | . | . . | . | . | | | |
#> [18,] . . . . . . . . | . | . . | | | | . | |
#> [19,] . . . . . . . . . | . | . . | . | | | |
#> [20,] . . . . . . . . . . | . | . . | | | | |
n <- 1e4
bk <- c(0:5, 7,11)
bMat <- matrix(1:8, n, 8, byrow=TRUE)
bLis <- as.data.frame(bMat)
B <- bandSparse(n, k = bk, diag = bLis)
Bs <- bandSparse(n, k = bk, diag = bLis, symmetric=TRUE)
B [1:15, 1:30]
#> 15 x 30 sparse Matrix of class "dgCMatrix"
#>
#> [1,] 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . . . . . . .
#> [2,] . 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . . . . . .
#> [3,] . . 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . . . . .
#> [4,] . . . 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . . . .
#> [5,] . . . . 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . . .
#> [6,] . . . . . 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . .
#> [7,] . . . . . . 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . .
#> [8,] . . . . . . . 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . .
#> [9,] . . . . . . . . 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . .
#> [10,] . . . . . . . . . 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . .
#> [11,] . . . . . . . . . . 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . .
#> [12,] . . . . . . . . . . . 1 2 3 4 5 6 . 7 . . . 8 . . . . . . .
#> [13,] . . . . . . . . . . . . 1 2 3 4 5 6 . 7 . . . 8 . . . . . .
#> [14,] . . . . . . . . . . . . . 1 2 3 4 5 6 . 7 . . . 8 . . . . .
#> [15,] . . . . . . . . . . . . . . 1 2 3 4 5 6 . 7 . . . 8 . . . .
Bs[1:15, 1:30]
#> 15 x 30 sparse Matrix of class "dgCMatrix"
#>
#> [1,] 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . . . . . . .
#> [2,] 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . . . . . .
#> [3,] 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . . . . .
#> [4,] 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . . . .
#> [5,] 5 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . . .
#> [6,] 6 5 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . . .
#> [7,] . 6 5 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . . .
#> [8,] 7 . 6 5 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . . .
#> [9,] . 7 . 6 5 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . . .
#> [10,] . . 7 . 6 5 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . . .
#> [11,] . . . 7 . 6 5 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . . . .
#> [12,] 8 . . . 7 . 6 5 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . . .
#> [13,] . 8 . . . 7 . 6 5 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . . .
#> [14,] . . 8 . . . 7 . 6 5 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . . .
#> [15,] . . . 8 . . . 7 . 6 5 4 3 2 1 2 3 4 5 6 . 7 . . . 8 . . . .
## can use a list *or* a matrix for specifying the diagonals:
stopifnot(identical(B, bandSparse(n, k = bk, diag = bMat)),
identical(Bs, bandSparse(n, k = bk, diag = bMat, symmetric=TRUE))
, inherits(B, "dtCMatrix") # triangular!
)