| lsolve.bicg {Rlinsolve} | R 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 |
B |
a vector of length |
xinit |
a length- |
reltol |
tolerance level for stopping iterations. |
maxiter |
maximum number of iterations allowed. |
preconditioner |
an |
verbose |
a logical; |
Value
a named list containing
- x
solution; a vector of length
nor a matrix of size(n\times k).- iter
the number of iterations required.
- errors
a vector of errors for stopping criterion.
References
Fletcher R (1976). “Conjugate gradient methods for indefinite systems.” In Watson GA (ed.), Numerical Analysis, volume 506, 73–89. Springer Berlin Heidelberg, Berlin, Heidelberg. ISBN 978-3-540-07610-0 978-3-540-38129-7.
Voevodin VV (1983). “The question of non-self-adjoint extension of the conjugate gradients method is closed.” USSR Computational Mathematics and Mathematical Physics, 23(2), 143–144. ISSN 00415553.
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)