hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_sparse_complex_gen_solve_bdilu (f11du)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_sparse_complex_gen_solve_bdilu (f11du) solves a complex sparse non-Hermitian system of linear equations, represented in coordinate storage format, using a restarted generalized minimal residual (RGMRES), conjugate gradient squared (CGS), stabilized bi-conjugate gradient (Bi-CGSTAB), or transpose-free quasi-minimal residual (TFQMR) method, with block Jacobi or additive Schwarz preconditioning.

Syntax

[x, rnorm, itn, ifail] = f11du(method, nz, a, irow, icol, nb, istb, indb, ipivp, ipivq, istr, idiag, b, m, tol, maxitn, x, 'n', n, 'la', la, 'lindb', lindb)
[x, rnorm, itn, ifail] = nag_sparse_complex_gen_solve_bdilu(method, nz, a, irow, icol, nb, istb, indb, ipivp, ipivq, istr, idiag, b, m, tol, maxitn, x, 'n', n, 'la', la, 'lindb', lindb)

Description

nag_sparse_complex_gen_solve_bdilu (f11du) solves a complex sparse non-Hermitian linear system of equations
Ax=b,  
using a preconditioned RGMRES (see Saad and Schultz (1986)), CGS (see Sonneveld (1989)), Bi-CGSTAB() (see Van der Vorst (1989) and Sleijpen and Fokkema (1993)), or TFQMR (see Freund and Nachtigal (1991) and Freund (1993)) method.
nag_sparse_complex_gen_solve_bdilu (f11du) uses the incomplete (possibly overlapping) block LU factorization determined by nag_sparse_complex_gen_precon_bdilu (f11dt) as the preconditioning matrix. A call to nag_sparse_complex_gen_solve_bdilu (f11du) must always be preceded by a call to nag_sparse_complex_gen_precon_bdilu (f11dt). Alternative preconditioners for the same storage scheme are available by calling nag_sparse_complex_gen_solve_ilu (f11dq) or nag_sparse_complex_gen_solve_jacssor (f11ds).
The matrix A, and the preconditioning matrix M, are represented in coordinate storage (CS) format (see Coordinate storage (CS) format in the F11 Chapter Introduction) in the arrays a, irow and icol, as returned from nag_sparse_complex_gen_precon_bdilu (f11dt). The array a holds the nonzero entries in these matrices, while irow and icol hold the corresponding row and column indices.
nag_sparse_complex_gen_solve_bdilu (f11du) is a Black Box function which calls nag_sparse_complex_gen_basic_setup (f11br), nag_sparse_complex_gen_basic_solver (f11bs) and nag_sparse_complex_gen_basic_diag (f11bt). If you wish to use an alternative storage scheme, preconditioner, or termination criterion, or require additional diagnostic information, you should call these underlying functions directly.

References

Freund R W (1993) A transpose-free quasi-minimal residual algorithm for non-Hermitian linear systems SIAM J. Sci. Comput. 14 470–482
Freund R W and Nachtigal N (1991) QMR: a Quasi-Minimal Residual Method for Non-Hermitian Linear Systems Numer. Math. 60 315–339
Saad Y and Schultz M (1986) GMRES: a generalized minimal residual algorithm for solving nonsymmetric linear systems SIAM J. Sci. Statist. Comput. 7 856–869
Sleijpen G L G and Fokkema D R (1993) BiCGSTAB for linear equations involving matrices with complex spectrum ETNA 1 11–32
Sonneveld P (1989) CGS, a fast Lanczos-type solver for nonsymmetric linear systems SIAM J. Sci. Statist. Comput. 10 36–52
Van der Vorst H (1989) Bi-CGSTAB, a fast and smoothly converging variant of Bi-CG for the solution of nonsymmetric linear systems SIAM J. Sci. Statist. Comput. 13 631–644

Parameters

Compulsory Input Parameters

1:     method – string
Specifies the iterative method to be used.
method='RGMRES'
Restarted generalized minimum residual method.
method='CGS'
Conjugate gradient squared method.
method='BICGSTAB'
Bi-conjugate gradient stabilized () method.
method='TFQMR'
Transpose-free quasi-minimal residual method.
Constraint: method='RGMRES', 'CGS', 'BICGSTAB' or 'TFQMR'.
2:     nz int64int32nag_int scalar
3:     ala – complex array
4:     irowla int64int32nag_int array
5:     icolla int64int32nag_int array
6:     nb int64int32nag_int scalar
7:     istbnb+1 int64int32nag_int array
8:     indblindb int64int32nag_int array
9:     ipivplindb int64int32nag_int array
10:   ipivqlindb int64int32nag_int array
11:   istrlindb+1 int64int32nag_int array
12:   idiaglindb int64int32nag_int array
The values returned in arrays irow, icol, ipivp, ipivq, istr and idiag by a previous call to nag_sparse_complex_gen_precon_bdilu (f11dt).
The arrays istb, indb and a together with the scalars n, nz, la, nb and lindb must be the same values that were supplied in the preceding call to nag_sparse_complex_gen_precon_bdilu (f11dt).
The values returned in arrays irow, icol, ipivp, ipivq, istr and idiag by a previous call to nag_sparse_complex_gen_precon_bdilu (f11dt).
The arrays istb, indb and the scalars nb and lindb must be the same values that were supplied in the preceding call to nag_sparse_complex_gen_precon_bdilu (f11dt).
13:   bn – complex array
The right-hand side vector b.
14:   m int64int32nag_int scalar
If method='RGMRES', m is the dimension of the restart subspace.
If method='BICGSTAB', m is the order  of the polynomial Bi-CGSTAB method.
Otherwise, m is not referenced.
Constraints:
  • if method='RGMRES', 0<mminn,50;
  • if method='BICGSTAB', 0<mminn,10.
15:   tol – double scalar
The required tolerance. Let xk denote the approximate solution at iteration k, and rk the corresponding residual. The algorithm is considered to have converged at iteration k if
rkτ×b+Axk.  
If tol0.0, τ=maxε,nε is used, where ε is the machine precision. Otherwise τ=maxtol,10ε,nε is used.
Constraint: tol<1.0.
16:   maxitn int64int32nag_int scalar
The maximum number of iterations allowed.
Constraint: maxitn1.
17:   xn – complex array
An initial approximation to the solution vector x.

Optional Input Parameters

1:     n int64int32nag_int scalar
2:     la int64int32nag_int scalar
3:     lindb int64int32nag_int scalar
Default: For la, the dimension of the arrays a, icol, irow. (An error is raised if these dimensions are not equal.)For lindb, the dimension of the arrays idiag, indb, ipivp, ipivq. (An error is raised if these dimensions are not equal.)For n, the dimension of the arrays b, x. (An error is raised if these dimensions are not equal.)
The values returned in arrays irow, icol, ipivp, ipivq, istr and idiag by a previous call to nag_sparse_complex_gen_precon_bdilu (f11dt).
The arrays istb, indb and a together with the scalars n, nz, la, nb and lindb must be the same values that were supplied in the preceding call to nag_sparse_complex_gen_precon_bdilu (f11dt).
Default: For la, the dimension of the arrays a, icol, irow. (An error is raised if these dimensions are not equal.)For lindb, the dimension of the arrays idiag, indb, ipivp, ipivq. (An error is raised if these dimensions are not equal.)For n, the dimension of the arrays b, x. (An error is raised if these dimensions are not equal.)
The values returned in arrays irow, icol, ipivp, ipivq, istr and idiag by a previous call to nag_sparse_complex_gen_precon_bdilu (f11dt).
The arrays istb, indb and the scalars nb and lindb must be the same values that were supplied in the preceding call to nag_sparse_complex_gen_precon_bdilu (f11dt).

Output Parameters

1:     xn – complex array
An improved approximation to the solution vector x.
2:     rnorm – double scalar
The final value of the residual norm rk, where k is the output value of itn.
3:     itn int64int32nag_int scalar
The number of iterations carried out.
4:     ifail int64int32nag_int scalar
ifail=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

Errors or warnings detected by the function:
   ifail=1
Constraint: 1nbn.
Constraint: istbb+1>istbb, for b=1,2,,nb.
Constraint: if method='RGMRES', 1mminn,_.
Constraint: istb11.
Constraint: la2×nz.
Constraint: lindbistbnb+1-1.
Constraint: 1indbmn, for m=1,2,,istbnb+1-1 
Constraint: maxitn1.
Constraint: method='RGMRES', 'CGS' or 'BICGSTAB'.
Constraint: nzn2.
Constraint: nz1.
Constraint: n1.
Constraint: tol<1.0.
lwork is too small.
   ifail=2
Constraint: 1icolin, for i=1,2,,nz.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_complex_gen_precon_bdilu (f11dt) and nag_sparse_complex_gen_solve_bdilu (f11du).
Constraint: 1irowin, for i=1,2,,nz.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_complex_gen_precon_bdilu (f11dt) and nag_sparse_complex_gen_solve_bdilu (f11du).
On entry, element _ of a was out of order.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_complex_gen_precon_bdilu (f11dt) and nag_sparse_complex_gen_solve_bdilu (f11du).
On entry, location _ of irow,icol was a duplicate.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_complex_gen_precon_bdilu (f11dt) and nag_sparse_complex_gen_solve_bdilu (f11du).
   ifail=3
The CS representation of the preconditioner is invalid.
Check that a, irow, icol, ipivp, ipivq, istr and idiag have not been corrupted between calls to nag_sparse_complex_gen_precon_bdilu (f11dt) and nag_sparse_complex_gen_solve_bdilu (f11du).
   ifail=4
The required accuracy could not be obtained. However, a reasonable accuracy may have been achieved. You should check the output value of rnorm for acceptability. This error code usually implies that your problem has been fully and satisfactorily solved to within or close to the accuracy available on your system. Further iterations are unlikely to improve on this situation.
   ifail=5
The solution has not converged after _ iterations.
   ifail=6
Algorithmic breakdown. A solution is returned, although it is possible that it is completely inaccurate.
   ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
   ifail=-399
Your licence key may have expired or may not have been installed correctly.
   ifail=-999
Dynamic memory allocation failed.

Accuracy

On successful termination, the final residual rk=b-Axk, where k=itn, satisfies the termination criterion
rkτ×b+Axk.  
The value of the final residual norm is returned in rnorm.

Further Comments

The time taken by nag_sparse_complex_gen_solve_bdilu (f11du) for each iteration is roughly proportional to the value of nnzc returned from the preceding call to nag_sparse_complex_gen_precon_bdilu (f11dt).
The number of iterations required to achieve a prescribed accuracy cannot be easily determined a priori, as it can depend dramatically on the conditioning and spectrum of the preconditioned coefficient matrix A-=M-1A.

Example

This example program reads in a sparse matrix A and a vector b. It calls nag_sparse_complex_gen_precon_bdilu (f11dt), with the array lfill=0 and the array dtol=0.0, to compute an overlapping incomplete LU factorization. This is then used as an additive Schwarz preconditioner on a call to nag_sparse_complex_gen_solve_bdilu (f11du) which uses the RGMRES method to solve Ax=b.
function f11du_example


fprintf('f11du example results\n\n');

n   = int64(9);
nz  = int64(33);
a    = zeros(20*nz, 1);
irow = zeros(20*nz, 1, 'int64');
icol = zeros(20*nz, 1, 'int64');
a(1:nz) = [ 96 - 64i; -20 + 22i; -36 + 14i; 
           -12 + 10i;  96 - 64i; -20 + 22i; -36 + 14i;
           -12 + 10i;  96 - 64i; -36 + 14i;
           -28 + 18i;  96 - 64i; -20 + 22i; -36 + 14i;
           -28 + 18i; -12 + 10i;  96 - 64i; -20 + 22i; -36 + 14i;
           -28 + 18i; -12 + 10i;  96 - 64i; -36 + 14i;
           -28 + 18i;  96 - 64i; -20 + 22i;
           -28 + 18i; -12 + 10i;  96 - 64i; -20 + 22i; 
           -28 + 18i; -12 + 10i;  96 - 64i];
irow(1:nz) = [1; 1; 1;  2; 2; 2; 2;  3; 3; 3;  4; 4; 4; 4;  5; 5; 5; 5; 5; 
              6; 6; 6; 6;  7; 7; 7;  8; 8; 8; 8;  9; 9; 9];
icol(1:nz) = [1; 2; 4;  1; 2; 3; 5;  2; 3; 6;  1; 4; 5; 7;  2; 4; 5; 6; 8; 
              3; 5; 6; 9;  4; 7; 8;  5; 7; 8; 9;  6; 8; 9];

% 3 blocks
nb    = int64(3);
nover = 1;
lfill = [int64(0); 0; 0];
dtol  = [0; 0; 0];
pstrat = {'n'; 'n'; 'n'};
milu   = {'n'; 'n'; 'n'};

% Define diagonal block indices.
% In this example use blocks of MB consecutive rows and initialise
% assuming no overlap.
mb    = idivide(n+nb-1, nb);
istb  = zeros(nb+1, 1, 'int64');
indb  = zeros(3*n, 1, 'int64');
ipivp = zeros(3*n, 1, 'int64');
ipivq = zeros(3*n, 1, 'int64');
istb(1:nb) = [1:mb:nb*mb];
istb(nb+1) = n+1;
indb(1:n)  = [1:n];

% Modify indb and istb to account for overlap.
[istb, indb, ifail] = f11du_overlap(n, nz, irow, icol, nb, ...
                                    istb, indb, 3*n, nover);
if (ifail == -999)
  error('indb is too small, size of indb = %d', numel(indb));
end

% Calculate Factorisation
[a, irow, icol, ipivp, ipivq, istr, idiag, nnzc, npivm, ifail] = ...
     f11dt(n, nz, a, irow, icol, istb, indb, ...
           lfill, dtol, milu, ipivp, ipivq, 'pstrat', pstrat);

% RHS b and initial guess x
x = complex(zeros(n,1));
b = x;
b(1:n) = 100+4i;

% Solver algorithmic parameters
method = 'rgmres';
m      = int64(2);
tol    = 1e-6;
maxitn = int64(100);

% Solve Ax = b
[x, rnorm, itn, ifail] = f11du( ...
                                method, nz, a, irow, icol, nb, istb, indb, ...
                                ipivp, ipivq, istr, idiag, b, m, tol, ...
                                maxitn, x);

fprintf('Converged in %d iterations\n', itn);
fprintf('Final redidual norm = %16.3d\n\n', rnorm);
disp('Solution');
disp(x);



function [istb, indb, ifail] = f11du_overlap(n, nz, irow, icol, nb, ...
                                             istb, indb, lindb, nover)

  ifail = 0;

  % This function takes a set of row indices indb defining the diagonal
  % blocks to be used in f11dt to define a block Jacobi or additive Schwarz
  % preconditioner, and expands them to allow for nover levels of overlap.
  % The pointer array istb is also updated accordingly, so that the returned
  % values of istb and indb can be passed to f11dt to define overlapping
  % diagonal blocks.
  iwork = zeros(3*n+1, 1, 'int64');

  % Find the number of non-zero elements in each row of the matrix A, and
  % the start address of each row. Store the start addresses in
  % iwork(n+1,...,2*n+1).
  for k=1:nz
    iwork(irow(k)) = iwork(irow(k)) + 1;
  end
  iwork(n+1) = 1;
  for i = 1:n
    iwork(n+i+1) = iwork(n+i) + iwork(i);
  end

  % Loop over blocks
  for k=1:nb
    % Initialize marker array.
    iwork(1:n) = 0;

    % Mark the rows already in block k in the workspace array.
    for l = istb(k):istb(k+1)-1
      iwork(indb(l)) = 1;
    end

    % Loop over levels of overlap.
    for iover=1:nover
      % Initialize counter of new row indices to be added.
      ind = 0;

      % Loop over the rows currently in the diagonal block.
      for l = istb(k):istb(k+1)-1
        row = indb(l);

        % Loop over non-zero elements in row
        for i = iwork(n+row):iwork(n+row+1)-1

          % If the column index of the non-zero element is not in the
          % existing set for this block, store it to be added later, and
          % mark it in the marker array.
          if (iwork(icol(i))==0)
            iwork(icol(i)) = 1;
            ind = ind + 1;
            iwork(2*n+1+ind) = icol(i);
          end
        end
      end

      % Shift the indices in indb and add the new entries for block k.
      % Change istb accordingly.
      nadd = ind;
      if (istb(nb+1)+nadd-1>lindb) Then
        ifail = -999;
        return;
      end

      for i = istb(nb+1) - 1:-1:istb(k+1)
        indb(i+nadd) = indb(i);
      end
      n21 = 2*n + 1;
      ik = istb(k+1) - 1;
      indb(ik+1:ik+nadd) = iwork(n21+1:n21+nadd);
      istb(k+1:nb+1) = istb(k+1:nb+1) + nadd;
    end
  end
f11du example results

Converged in 8 iterations
Final redidual norm =        6.492e-04

Solution
   2.2040 + 1.6972i
   2.3511 + 1.9275i
   1.5931 + 1.4368i
   2.8641 + 1.9762i
   3.0687 + 2.2645i
   2.0467 + 1.6948i
   2.2065 + 1.3244i
   2.3724 + 1.5170i
   1.6254 + 1.1783i


PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015