Previous: , Up: Linear Algebra   [Contents][Index]


18.5 Specialized Solvers

: x = bicg (A, b)
: x = bicg (A, b, tol)
: x = bicg (A, b, tol, maxit)
: x = bicg (A, b, tol, maxit, M)
: x = bicg (A, b, tol, maxit, M1, M2)
: x = bicg (A, b, tol, maxit, M, [], x0)
: x = bicg (A, b, tol, maxit, M1, M2, x0)
: x = bicg (A, b, tol, maxit, M, [], x0, …)
: x = bicg (A, b, tol, maxit, M1, M2, x0, …)
: [x, flag, relres, iter, resvec] = bicg (A, b, …)

Solve the linear system of equations A * x = b by means of the Bi-Conjugate Gradient iterative method.

The input arguments are:

Any arguments which follow x0 are treated as parameters, and passed in an appropriate manner to any of the functions (Afun or Mfun) or that have been given to bicg.

The output parameters are:

Consider a trivial problem with a tridiagonal matrix

n = 20;
A = toeplitz (sparse ([1, 1], [1, 2], [2, 1] * n ^ 2, 1, n)) + ...
    toeplitz (sparse (1, 2, -1, 1, n) * n / 2, ...
              sparse (1, 2, 1, 1, n) * n / 2);
b = A * ones (n, 1);
restart = 5;
[M1, M2] = ilu (A);  # in this tridiag case, it corresponds to lu (A)
M = M1 * M2;
Afun = @(x, string) strcmp (string, "notransp") * (A * x) + ...
                     strcmp (string, "transp") * (A' * x);
Mfun = @(x, string) strcmp (string, "notransp") * (M \ x) + ...
                     strcmp (string, "transp") * (M' \ x);
M1fun = @(x, string) strcmp (string, "notransp") * (M1 \ x) + ...
                     strcmp (string, "transp") * (M1' \ x);
M2fun = @(x, string) strcmp (string, "notransp") * (M2 \ x) + ...
                     strcmp (string, "transp") * (M2' \ x);

EXAMPLE 1: simplest usage of bicg

x = bicg (A, b)

EXAMPLE 2: bicg with a function that computes A*x and A'*x

x = bicg (Afun, b, [], n)

EXAMPLE 3: bicg with a preconditioner matrix M

x = bicg (A, b, 1e-6, n, M)

EXAMPLE 4: bicg with a function as preconditioner

x = bicg (Afun, b, 1e-6, n, Mfun)

EXAMPLE 5: bicg with preconditioner matrices M1 and M2

x = bicg (A, b, 1e-6, n, M1, M2)

EXAMPLE 6: bicg with functions as preconditioners

x = bicg (Afun, b, 1e-6, n, M1fun, M2fun)

EXAMPLE 7: bicg with as input a function requiring an argument

function y = Ap (A, x, string, z)
  ## compute A^z * x or (A^z)' * x
  y = x;
  if (strcmp (string, "notransp"))
    for i = 1:z
      y = A * y;
    endfor
  elseif (strcmp (string, "transp"))
    for i = 1:z
      y = A' * y;
    endfor
  endif
endfunction

Apfun = @(x, string, p) Ap (A, x, string, p);
x = bicg (Apfun, b, [], [], [], [], [], 2);

Reference:

Y. Saad, Iterative Methods for Sparse Linear Systems, Second edition, 2003, SIAM.

See also: bicgstab, cgs, gmres, pcg, qmr, tfqmr.

: x = bicgstab (A, b, tol, maxit, M1, M2, x0, …)
: x = bicgstab (A, b, tol, maxit, M, [], x0, …)
: [x, flag, relres, iter, resvec] = bicgstab (A, b, …)

Solve A x = b using the stabilized Bi-conjugate gradient iterative method.

The input parameters are:

The arguments which follow x0 are treated as parameters, and passed in a proper way to any of the functions (A or M) which are passed to bicstab.

The output parameters are:

Let us consider a trivial problem with a tridiagonal matrix

n = 20;
A = toeplitz (sparse ([1, 1], [1, 2], [2, 1] * n ^ 2, 1, n))  + ...
    toeplitz (sparse (1, 2, -1, 1, n) * n / 2, ...
    sparse (1, 2, 1, 1, n) * n / 2);
b = A * ones (n, 1);
restart = 5;
[M1, M2] = ilu (A); # in this tridiag case, it corresponds to lu (A)
M = M1 * M2;
Afun = @(x) A * x;
Mfun = @(x) M \ x;
M1fun = @(x) M1 \ x;
M2fun = @(x) M2 \ x;

EXAMPLE 1: simplest usage of bicgstab

x = bicgstab (A, b, [], n)

EXAMPLE 2: bicgstab with a function which computes A * x

x = bicgstab (Afun, b, [], n)

EXAMPLE 3: bicgstab with a preconditioner matrix M

x = bicgstab (A, b, [], 1e-06, n, M)

EXAMPLE 4: bicgstab with a function as preconditioner

x = bicgstab (Afun, b, 1e-6, n, Mfun)

EXAMPLE 5: bicgstab with preconditioner matrices M1 and M2

x = bicgstab (A, b, [], 1e-6, n, M1, M2)

EXAMPLE 6: bicgstab with functions as preconditioners

x = bicgstab (Afun, b, 1e-6, n, M1fun, M2fun)

EXAMPLE 7: bicgstab with as input a function requiring an argument

function y = Ap (A, x, z) # compute A^z * x
   y = x;
   for i = 1:z
     y = A * y;
   endfor
 endfunction
Apfun = @(x, string, p) Ap (A, x, string, p);
x = bicgstab (Apfun, b, [], [], [], [], [], 2);

EXAMPLE 8: explicit example to show that bicgstab uses a right preconditioner

[M1, M2] = ilu (A + 0.1 * eye (n)); # factorization of A perturbed
M = M1 * M2;

## reference solution computed by bicgstab after one iteration
[x_ref, fl] = bicgstab (A, b, [], 1, M)

## right preconditioning
[y, fl] = bicgstab (A / M, b, [], 1)
x = M \ y # compare x and x_ref

Reference:

Y. Saad, Iterative Methods for Sparse Linear Systems, Second edition, 2003, SIAM

See also: bicg, cgs, gmres, pcg, qmr, tfqmr.

: x = cgs (A, b, tol, maxit, M1, M2, x0, …)
: x = cgs (A, b, tol, maxit, M, [], x0, …)
: [x, flag, relres, iter, resvec] = cgs (A, b, …)

Solve A x = b, where A is a square matrix, using the Conjugate Gradients Squared method.

The input arguments are:

The arguments which follow x0 are treated as parameters, and passed in a proper way to any of the functions (A or P) which are passed to cgs.

The output parameters are:

Let us consider a trivial problem with a tridiagonal matrix

n = 20;
A = toeplitz (sparse ([1, 1], [1, 2], [2, 1] * n ^ 2, 1, n))  + ...
    toeplitz (sparse (1, 2, -1, 1, n) * n / 2, ...
    sparse (1, 2, 1, 1, n) * n / 2);
b = A * ones (n, 1);
restart = 5;
[M1, M2] = ilu (A); # in this tridiag case it corresponds to chol (A)'
M = M1 * M2;
Afun = @(x) A * x;
Mfun = @(x) M \ x;
M1fun = @(x) M1 \ x;
M2fun = @(x) M2 \ x;

EXAMPLE 1: simplest usage of cgs

x = cgs (A, b, [], n)

EXAMPLE 2: cgs with a function which computes A * x

x = cgs (Afun, b, [], n)

EXAMPLE 3: cgs with a preconditioner matrix M

x = cgs (A, b, [], 1e-06, n, M)

EXAMPLE 4: cgs with a function as preconditioner

x = cgs (Afun, b, 1e-6, n, Mfun)

EXAMPLE 5: cgs with preconditioner matrices M1 and M2

x = cgs (A, b, [], 1e-6, n, M1, M2)

EXAMPLE 6: cgs with functions as preconditioners

x = cgs (Afun, b, 1e-6, n, M1fun, M2fun)

EXAMPLE 7: cgs with as input a function requiring an argument

function y = Ap (A, x, z) # compute A^z * x
   y = x;
   for i = 1:z
     y = A * y;
   endfor
 endfunction
Apfun = @(x, string, p) Ap (A, x, string, p);
x = cgs (Apfun, b, [], [], [], [], [], 2);

EXAMPLE 8: explicit example to show that cgs uses a right preconditioner

[M1, M2] = ilu (A + 0.3 * eye (n)); # factorization of A perturbed
M = M1 * M2;

## reference solution computed by cgs after one iteration
[x_ref, fl] = cgs (A, b, [], 1, M)

## right preconditioning
[y, fl] = cgs (A / M, b, [], 1)
x = M \ y # compare x and x_ref

References:

Y. Saad, Iterative Methods for Sparse Linear Systems, Second edition, 2003, SIAM

See also: pcg, bicgstab, bicg, gmres, qmr, tfqmr.

: x = gmres (A, b, restart, tol, maxit, M1, M2, x0, …)
: x = gmres (A, b, restart, tol, maxit, M, [], x0, …)
: [x, flag, relres, iter, resvec] = gmres (A, b, …)

Solve A x = b using the Preconditioned GMRES iterative method with restart, a.k.a. PGMRES(restart).

The input arguments are:

The arguments which follow x0 are treated as parameters, and passed in a proper way to any of the functions (A or M or M1 or M2) which are passed to gmres.

The outputs are:

Let us consider a trivial problem with a tridiagonal matrix

n = 20;
A = toeplitz (sparse ([1, 1], [1, 2], [2, 1] * n ^ 2, 1, n))  + ...
    toeplitz (sparse (1, 2, -1, 1, n) * n / 2, ...
    sparse (1, 2, 1, 1, n) * n / 2);
b = A * ones (n, 1);
restart = 5;
[M1, M2] = ilu (A); # in this tridiag case, it corresponds to lu (A)
M = M1 * M2;
Afun = @(x) A * x;
Mfun = @(x) M \ x;
M1fun = @(x) M1 \ x;
M2fun = @(x) M2 \ x;

EXAMPLE 1: simplest usage of gmres

x = gmres (A, b, [], [], n)

EXAMPLE 2: gmres with a function which computes A * x

x = gmres (Afun, b, [], [], n)

EXAMPLE 3: usage of gmres with the restart

x = gmres (A, b, restart);

EXAMPLE 4: gmres with a preconditioner matrix M with and without restart

x = gmres (A, b, [], 1e-06, n, M)
x = gmres (A, b, restart, 1e-06, n, M)

EXAMPLE 5: gmres with a function as preconditioner

x = gmres (Afun, b, [], 1e-6, n, Mfun)

EXAMPLE 6: gmres with preconditioner matrices M1 and M2

x = gmres (A, b, [], 1e-6, n, M1, M2)

EXAMPLE 7: gmres with functions as preconditioners

x = gmres (Afun, b, 1e-6, n, M1fun, M2fun)

EXAMPLE 8: gmres with as input a function requiring an argument

  function y = Ap (A, x, p) # compute A^p * x
     y = x;
     for i = 1:p
       y = A * y;
     endfor
  endfunction
Apfun = @(x, p) Ap (A, x, p);
x = gmres (Apfun, b, [], [], [], [], [], [], 2);

EXAMPLE 9: explicit example to show that gmres uses a left preconditioner

[M1, M2] = ilu (A + 0.1 * eye (n)); # factorization of A perturbed
M = M1 * M2;

## reference solution computed by gmres after two iterations
[x_ref, fl] = gmres (A, b, [], [], 1, M)

## left preconditioning
[x, fl] = gmres (M \ A, M \ b, [], [], 1)
x # compare x and x_ref

Reference:

Y. Saad, Iterative Methods for Sparse Linear Systems, Second edition, 2003, SIAM

See also: bicg, bicgstab, cgs, pcg, pcr, qmr, tfqmr.

: x = qmr (A, b, rtol, maxit, M1, M2, x0)
: x = qmr (A, b, rtol, maxit, P)
: [x, flag, relres, iter, resvec] = qmr (A, b, …)

Solve A x = b using the Quasi-Minimal Residual iterative method (without look-ahead).

A can be passed as a matrix or as a function handle or inline function f such that f(x, "notransp") = A*x and f(x, "transp") = A'*x.

The preconditioner P is given as P = M1 * M2. Both M1 and M2 can be passed as a matrix or as a function handle or inline function g such that g(x, "notransp") = M1 \ x or g(x, "notransp") = M2 \ x and g(x, "transp") = M1' \ x or g(x, "transp") = M2' \ x.

If called with more than one output parameter

References:

  1. R. Freund and N. Nachtigal, QMR: a quasi-minimal residual method for non-Hermitian linear systems, Numerische Mathematik, 1991, 60, pp. 315–339.
  2. R. Barrett, M. Berry, T. Chan, J. Demmel, J. Donato, J. Dongarra, V. Eijkhour, R. Pozo, C. Romine, and H. van der Vorst, Templates for the solution of linear systems: Building blocks for iterative methods, SIAM, 2nd ed., 1994.

See also: bicg, bicgstab, cgs, gmres, pcg.

: x = tfqmr (A, b, tol, maxit, M1, M2, x0, …)
: x = tfqmr (A, b, tol, maxit, M, [], x0, …)
: [x, flag, relres, iter, resvec] = tfqmr (A, b, …)

Solve A x = b using the Transpose-Tree qmr method, based on the cgs.

The input parameters are:

The arguments which follow x0 are treated as parameters, and passed in a proper way to any of the functions (A or M) which are passed to tfqmr.

The output parameters are:

Let us consider a trivial problem with a tridiagonal matrix

n = 20;
A = toeplitz (sparse ([1, 1], [1, 2], [2, 1] * n ^ 2, 1, n))  + ...
    toeplitz (sparse (1, 2, -1, 1, n) * n / 2, ...
    sparse (1, 2, 1, 1, n) * n / 2);
b = A * ones (n, 1);
restart = 5;
[M1, M2] = ilu (A); # in this tridiag case it corresponds to chol (A)'
M = M1 * M2;
Afun = @(x) A * x;
Mfun = @(x) M \ x;
M1fun = @(x) M1 \ x;
M2fun = @(x) M2 \ x;

EXAMPLE 1: simplest usage of tfqmr

x = tfqmr (A, b, [], n)

EXAMPLE 2: tfqmr with a function which computes A * x

x = tfqmr (Afun, b, [], n)

EXAMPLE 3: tfqmr with a preconditioner matrix M

x = tfqmr (A, b, [], 1e-06, n, M)

EXAMPLE 4: tfqmr with a function as preconditioner

x = tfqmr (Afun, b, 1e-6, n, Mfun)

EXAMPLE 5: tfqmr with preconditioner matrices M1 and M2

x = tfqmr (A, b, [], 1e-6, n, M1, M2)

EXAMPLE 6: tfmqr with functions as preconditioners

x = tfqmr (Afun, b, 1e-6, n, M1fun, M2fun)

EXAMPLE 7: tfqmr with as input a function requiring an argument

function y = Ap (A, x, z) # compute A^z * x
   y = x;
   for i = 1:z
     y = A * y;
   endfor
 endfunction
Apfun = @(x, string, p) Ap (A, x, string, p);
x = tfqmr (Apfun, b, [], [], [], [], [], 2);

EXAMPLE 8: explicit example to show that tfqmr uses a right preconditioner

[M1, M2] = ilu (A + 0.3 * eye (n)); # factorization of A perturbed
M = M1 * M2;

## reference solution computed by tfqmr after one iteration
[x_ref, fl] = tfqmr (A, b, [], 1, M)

## right preconditioning
[y, fl] = tfqmr (A / M, b, [], 1)
x = M \ y # compare x and x_ref

Reference:

Y. Saad, Iterative Methods for Sparse Linear Systems, Second edition, 2003, SIAM

See also: bicg, bicgstab, cgs, gmres, pcg, qmr, pcr.


Previous: , Up: Linear Algebra   [Contents][Index]