krylov_BICG: Biconjugate Gradient method

lsolve.bicgR Documentation

Biconjugate Gradient method

Description

Biconjugate Gradient(BiCG) method is a modification of Conjugate Gradient for nonsymmetric systems using evaluations with respect to A^T as well as A in matrix-vector multiplications. For an overdetermined system where nrow(A)>ncol(A), it is automatically transformed to the normal equation. Underdetermined system - nrow(A)<ncol(A) - is not supported. Preconditioning matrix M, in theory, should be symmetric and positive definite with fast computability for inverse, though it is not limited until the solver level.

Usage

lsolve.bicg(
  A,
  B,
  xinit = NA,
  reltol = 1e-05,
  maxiter = 10000,
  preconditioner = diag(ncol(A)),
  verbose = TRUE
)

Arguments

A

an (m\times n) dense or sparse matrix. See also sparseMatrix.

B

a vector of length m or an (m\times k) matrix (dense or sparse) for solving k systems simultaneously.

xinit

a length-n vector for initial starting point. NA to start from a random initial point near 0.

reltol

tolerance level for stopping iterations.

maxiter

maximum number of iterations allowed.

preconditioner

an (n\times n) preconditioning matrix; default is an identity matrix.

verbose

a logical; TRUE to show progress of computation.

Value

a named list containing

x

solution; a vector of length n or a matrix of size (n\times k).

iter

the number of iterations required.

errors

a vector of errors for stopping criterion.

References

\insertRef

watson_conjugate_1976Rlinsolve

\insertRef

voevodin_question_1983Rlinsolve

Examples

## Overdetermined System
set.seed(100)
A = matrix(rnorm(10*5),nrow=10)
x = rnorm(5)
b = A%*%x

out1 = lsolve.cg(A,b)
out2 = lsolve.bicg(A,b)
matout = cbind(matrix(x),out1$x, out2$x);
colnames(matout) = c("true x","CG result", "BiCG result")
print(matout)


kisungyou/Rlinsolve documentation built on Aug. 21, 2023, 4:52 a.m.