Next: , Previous: Sparse Linear Algebra, Up: Sparse Matrices


22.3 Iterative Techniques applied to sparse matrices

The left division \ and right division / operators, discussed in the previous section, use direct solvers to resolve a linear equation of the form x = A \ b or x = b / A. Octave equally includes a number of functions to solve sparse linear equations using iterative techniques.

— Function File: x = pcg (A, b, tol, maxit, m1, m2, x0, ...)
— Function File: [x, flag, relres, iter, resvec, eigest] = pcg (...)

Solve the linear system of equations A * x = b by means of the Preconditioned Conjugate Gradient iterative 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 m) which are passed to pcg. See the examples below for further details. The output arguments are

Let us consider a trivial problem with a diagonal matrix (we exploit the sparsity of A)

               n = 10;
               A = diag (sparse (1:n));
               b = rand (n, 1);
               [l, u, p, q] = luinc (A, 1.e-3);

Example 1: Simplest use of pcg

            x = pcg(A,b)

Example 2: pcg with a function which computes A * x

            function y = apply_a (x)
              y = [1:N]'.*x;
            endfunction
          
            x = pcg ("apply_a", b)

Example 3: pcg with a preconditioner: l * u

          x = pcg (A, b, 1.e-6, 500, l*u);

Example 4: pcg with a preconditioner: l * u. Faster than Example 3 since lower and upper triangular matrices are easier to invert

          x = pcg (A, b, 1.e-6, 500, l, u);

Example 5: Preconditioned iteration, with full diagnostics. The preconditioner (quite strange, because even the original matrix A is trivial) is defined as a function

            function y = apply_m (x)
              k = floor (length (x) - 2);
              y = x;
              y(1:k) = x(1:k)./[1:k]';
            endfunction
          
            [x, flag, relres, iter, resvec, eigest] = ...
                               pcg (A, b, [], [], "apply_m");
            semilogy (1:iter+1, resvec);

Example 6: Finally, a preconditioner which depends on a parameter k.

            function y = apply_M (x, varargin)
            K = varargin{1};
            y = x;
            y(1:K) = x(1:K)./[1:K]';
            endfunction
          
            [x, flag, relres, iter, resvec, eigest] = ...
                 pcg (A, b, [], [], "apply_m", [], [], 3)

References:

  1. C.T. Kelley, Iterative Methods for Linear and Nonlinear Equations, SIAM, 1995. (the base PCG algorithm)
  2. Y. Saad, Iterative Methods for Sparse Linear Systems, PWS 1996. (condition number estimate from PCG) Revised version of this book is available online at http://www-users.cs.umn.edu/~saad/books.html

See also: sparse, pcr.

— Function File: x = pcr (A, b, tol, maxit, m, x0, ...)
— Function File: [x, flag, relres, iter, resvec] = pcr (...)

Solve the linear system of equations A * x = b by means of the Preconditioned Conjugate Residuals iterative 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 m) which are passed to pcr. See the examples below for further details. The output arguments are

Let us consider a trivial problem with a diagonal matrix (we exploit the sparsity of A)

               n = 10;
               A = sparse (diag (1:n));
               b = rand (N, 1);

Example 1: Simplest use of pcr

            x = pcr(A, b)

Example 2: pcr with a function which computes A * x.

            function y = apply_a (x)
              y = [1:10]'.*x;
            endfunction
          
            x = pcr ("apply_a", b)

Example 3: Preconditioned iteration, with full diagnostics. The preconditioner (quite strange, because even the original matrix A is trivial) is defined as a function

            function y = apply_m (x)
              k = floor (length(x)-2);
              y = x;
              y(1:k) = x(1:k)./[1:k]';
            endfunction
          
            [x, flag, relres, iter, resvec] = ...
                               pcr (A, b, [], [], "apply_m")
            semilogy([1:iter+1], resvec);

Example 4: Finally, a preconditioner which depends on a parameter k.

            function y = apply_m (x, varargin)
              k = varargin{1};
              y = x; y(1:k) = x(1:k)./[1:k]';
            endfunction
          
            [x, flag, relres, iter, resvec] = ...
                               pcr (A, b, [], [], "apply_m"', [], 3)

References:

[1] W. Hackbusch, Iterative Solution of Large Sparse Systems of Equations, section 9.5.4; Springer, 1994

See also: sparse, pcg.

The speed with which an iterative solver converges to a solution can be accelerated with the use of a pre-conditioning matrix M. In this case the linear equation M^-1 * x = M^-1 * A \ b is solved instead. Typical pre-conditioning matrices are partial factorizations of the original matrix.

— Loadable Function: [L, U, P, Q] = luinc (A, '0')
— Loadable Function: [L, U, P, Q] = luinc (A, droptol)
— Loadable Function: [L, U, P, Q] = luinc (A, opts)

Produce the incomplete LU factorization of the sparse matrix A. Two types of incomplete factorization are possible, and the type is determined by the second argument to luinc.

Called with a second argument of '0', the zero-level incomplete LU factorization is produced. This creates a factorization of A where the position of the non-zero arguments correspond to the same positions as in the matrix A.

Alternatively, the fill-in of the incomplete LU factorization can be controlled through the variable droptol or the structure opts. The umfpack multifrontal factorization code by Tim A. Davis is used for the incomplete LU factorization, (availability http://www.cise.ufl.edu/research/sparse/umfpack/)

droptol determines the values below which the values in the LU  factorization are dropped and replaced by zero. It must be a positive scalar, and any values in the factorization whose absolute value are less than this value are dropped, expect if leaving them increase the sparsity of the matrix. Setting droptol to zero results in a complete LU factorization which is the default.

opts is a structure containing one or more of the fields

droptol
The drop tolerance as above. If opts only contains droptol then this is equivalent to using the variable droptol.
milu
A logical variable flagging whether to use the modified incomplete LU  factorization. In the case that milu is true, the dropped values are subtracted from the diagonal of the matrix U of the factorization. The default is false.
udiag
A logical variable that flags whether zero elements on the diagonal of U should be replaced with droptol to attempt to avoid singular factors. The default is false.
thresh
Defines the pivot threshold in the interval [0,1]. Values outside that range are ignored.

All other fields in opts are ignored. The outputs from luinc are the same as for lu.

Given the string argument 'vector', luinc returns the values of p q as vector values.

See also: sparse, lu.