bandSparse function

Construct Sparse Banded Matrix from (Sup-/Super-) Diagonals

Construct Sparse Banded Matrix from (Sup-/Super-) Diagonals

Construct 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"))

Arguments

  • n,m: the matrix dimension (n,m)=(nrow,ncol)(n,m) = (nrow, ncol).

  • k: integer vector of diagonal numbers , with identical meaning as in band(*, k), i.e., relative to the main diagonal, which is k=0.

  • diagonals: 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 nxdn' x d matrix, where d <- length(k) and n>=min(n,m)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.

  • symmetric: 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).

  • repr: 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.

  • giveCsparse: (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.

Returns

a sparse matrix (of class

CsparseMatrix) of dimension nxmn x m

with diagonal bands as specified.

See Also

band, for extraction of matrix bands; bdiag, diag, sparseMatrix, Matrix.

Examples

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 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)) ## 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)) (s4. <- as(drop0(s4), "nsparseMatrix")) 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] Bs[1:15, 1:30] ## 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! )
  • Maintainer: Martin Maechler
  • License: GPL (>= 2) | file LICENCE
  • Last published: 2025-03-11