krylov_CGS: Conjugate Gradient Squared method

lsolve.cgsR Documentation

Conjugate Gradient Squared method

Description

Conjugate Gradient Squared(CGS) method is an extension of Conjugate Gradient method where the system is symmetric and positive definite. It aims at achieving faster convergence using an idea of contraction operator twice. For a square matrix A,it is required to be symmetric and positive definite. 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.cgs(
  A,
  B,
  xinit = NA,
  reltol = 1e-05,
  maxiter = 10000,
  preconditioner = diag(ncol(A)),
  adjsym = TRUE,
  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.

adjsym

a logical; TRUE to symmetrize the system by transforming the system into normal equation, FALSE otherwise.

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

sonneveld_cgs_1989Rlinsolve

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.cgs(A,b)
matout = cbind(matrix(x),out1$x, out2$x);
colnames(matout) = c("true x","CG result", "CGS result")
print(matout)


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