sp_solve_sym {bigsparser}R Documentation

Solver for symmetric SFBM

Description

Solve Ax=b where A is a symmetric SFBM, and b is a vector.

Usage

sp_solve_sym(
  A,
  b,
  add_to_diag = rep(0, ncol(A)),
  tol = 1e-10,
  maxiter = 10 * ncol(A)
)

Arguments

A

A symmetric SFBM.

b

A vector.

add_to_diag

Vector (or single value) to virtually add to the diagonal of A. Default is 0s.

tol

Tolerance for convergence. Default is 1e-10.

maxiter

Maximum number of iterations for convergence.

Value

The vector x, solution of Ax=b.

Examples

N <- 100
spmat <- Matrix::rsparsematrix(N, N, 0.01, symmetric = TRUE)
X <- bigsparser::as_SFBM(as(spmat, "dgCMatrix"))
b <- runif(N)

test <- tryCatch(as.vector(Matrix::solve(spmat, b)), error = function(e) print(e))
test2 <- tryCatch(sp_solve_sym(X, b), error = function(e) print(e))

test3 <- as.vector(Matrix::solve(spmat + Matrix::Diagonal(N, 1:N), b))
test4 <- sp_solve_sym(X, b, add_to_diag = 1:N)
all.equal(test3, test4)


[Package bigsparser version 0.7.1 Index]