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_basic_solver (f11bs)

Purpose

nag_sparse_complex_gen_basic_solver (f11bs) is an iterative solver for a complex general (non-Hermitian) system of simultaneous linear equations; nag_sparse_complex_gen_basic_solver (f11bs) is the second in a suite of three functions, where the first function, nag_sparse_complex_gen_basic_setup (f11br), must be called prior to nag_sparse_complex_gen_basic_solver (f11bs) to set up the suite, and the third function in the suite, nag_sparse_complex_gen_basic_diag (f11bt), can be used to return additional information about the computation.
These three functions are suitable for the solution of large sparse general (non-Hermitian) systems of equations.

Syntax

[irevcm, u, v, work, ifail] = f11bs(irevcm, u, v, wgt, work, 'lwork', lwork)
[irevcm, u, v, work, ifail] = nag_sparse_complex_gen_basic_solver(irevcm, u, v, wgt, work, 'lwork', lwork)

Description

nag_sparse_complex_gen_basic_solver (f11bs) solves the general (non-Hermitian) system of linear simultaneous equations Ax = bAx=b of order nn, where nn is large and the coefficient matrix AA is sparse, using one of four available methods: RGMRES, the preconditioned restarted generalized minimum residual method (see Saad and Schultz (1986)); CGS, the preconditioned conjugate gradient squared method (see Sonneveld (1989)); Bi-CGSTAB(), the bi-conjugate gradient stabilized method of order  (see Van der Vorst (1989) and Sleijpen and Fokkema (1993)); or TFQMR, the transpose-free quasi-minimal residual method (see Freund and Nachtigal (1991) and Freund (1993)).
For a general description of the methods employed you are referred to Section [Description] in (f11br).
nag_sparse_complex_gen_basic_solver (f11bs) can solve the system after the first function in the suite, nag_sparse_complex_gen_basic_setup (f11br), has been called to initialize the computation and specify the method of solution. The third function in the suite, nag_sparse_complex_gen_basic_diag (f11bt), can be used to return additional information generated by the computation during monitoring steps and after nag_sparse_complex_gen_basic_solver (f11bs) has completed its tasks.
nag_sparse_complex_gen_basic_solver (f11bs) uses reverse communication, i.e., it returns repeatedly to the calling program with the parameter irevcm (see Section [Parameters]) set to specified values which require the calling program to carry out one of the following tasks:
Through the parameter irevcm the calling program can cause immediate or tidy termination of the execution. On final exit, the last iterates of the solution and of the residual vectors of the original system of equations are returned.
Reverse communication has the following advantages.
  1. Maximum flexibility in the representation and storage of sparse matrices: all matrix operations are performed outside the solver function, thereby avoiding the need for a complicated interface with enough flexibility to cope with all types of storage schemes and sparsity patterns. This applies also to preconditioners.
  2. Enhanced user interaction: you can closely monitor the progress of the solution and tidy or immediate termination can be requested. This is useful, for example, when alternative termination criteria are to be employed or in case of failure of the external functions used to perform matrix operations.

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
Higham N J (1988) FORTRAN codes for estimating the one-norm of a real or complex matrix, with applications to condition estimation ACM Trans. Math. Software 14 381–396
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

Note: this function uses reverse communication. Its use involves an initial entry, intermediate exits and re-entries, and a final exit, as indicated by the parameter irevcm. Between intermediate exits and re-entries, all parameters other than irevcm and v must remain unchanged.

Compulsory Input Parameters

1:     irevcm – int64int32nag_int scalar
On initial entry: irevcm = 0irevcm=0, otherwise an error condition will be raised.
On intermediate re-entry: must either be unchanged from its previous exit value, or can have one of the following values.
irevcm = 5irevcm=5
Tidy termination: the computation will terminate at the end of the current iteration. Further reverse communication exits may occur depending on when the termination request is issued. nag_sparse_complex_gen_basic_solver (f11bs) will then return with the termination code irevcm = 4irevcm=4. Note that before calling nag_sparse_complex_gen_basic_solver (f11bs) with irevcm = 5irevcm=5 the calling program must have performed the tasks required by the value of irevcm returned by the previous call to nag_sparse_complex_gen_basic_solver (f11bs), otherwise subsequently returned values may be invalid.
irevcm = 6irevcm=6
Immediate termination: nag_sparse_complex_gen_basic_solver (f11bs) will return immediately with termination code irevcm = 4irevcm=4 and with any useful information available. This includes the last iterate of the solution.
Immediate termination may be useful, for example, when errors are detected during matrix-vector multiplication or during the solution of the preconditioning equation.
Changing irevcm to any other value between calls will result in an error.
Constraint: on initial entry, irevcm = 0irevcm=0; on re-entry, either irevcm must remain unchanged or be reset to irevcm = 5irevcm=5 or 66.
2:     u( : :) – complex array
Note: the dimension of the array u must be at least nn.
On initial entry: an initial estimate, x0x0, of the solution of the system of equations Ax = bAx=b.
On intermediate re-entry: must remain unchanged.
3:     v( : :) – complex array
Note: the dimension of the array v must be at least nn.
On initial entry: the right-hand side bb of the system of equations Ax = bAx=b.
On intermediate re-entry: the returned value of irevcm determines the contents of v as follows.
If irevcm = 1irevcm=-1, 11 or 22, v must store the vector vv, the result of the operation specified by the value of irevcm returned by the previous call to nag_sparse_complex_gen_basic_solver (f11bs).
If irevcm = 3irevcm=3, v must remain unchanged.
4:     wgt( : :) – double array
Note: the dimension of the array wgt must be at least max (1,n)max(1,n).
The user-supplied weights, if these are to be used in the computation of the vector norms in the termination criterion (see Sections [Description] and [Parameters] in (f11br)).
Constraint: if weights are to be used, at least one element of wgt must be nonzero.
5:     work(lwork) – complex array
lwork, the dimension of the array, must satisfy the constraint lworklwreqlworklwreq, where lwreq is returned by nag_sparse_complex_gen_basic_setup (f11br).
On initial entry: the array work as returned by nag_sparse_complex_gen_basic_setup (f11br) (see also Section [Parameters] in (f11br)).
lwork, the dimension of the array, must satisfy the constraint lworklwreqlworklwreq, where lwreq is returned by nag_sparse_complex_gen_basic_setup (f11br).
On intermediate re-entry: must remain unchanged.

Optional Input Parameters

1:     lwork – int64int32nag_int scalar
Default: The dimension of the array work.
On initial entry: the dimension of the array work as declared in the (sub)program from which nag_sparse_complex_gen_basic_solver (f11bs) is called (see also Sections [Description] and [Parameters] in (f11br)). The required amount of workspace is as follows:
Method Requirements
RGMRES lwork = 120 + n(m + 3) + m(m + 5) + 1lwork=120+n(m+3)+m(m+5)+1, where mm is the dimension of the basis
CGS lwork = 120 + 7nlwork=120+7n
Bi-CGSTAB() lwork = 120 + (2n + )( + 2) + plwork=120+(2n+)(+2)+p, where  is the order of the method
TFQMR lwork = 120 + 10nlwork=120+10n,
where
  • p = 2np=2n, if > 1>1 and iterm = 2iterm=2 was supplied;
  • p = np=n, if > 1>1 and a preconditioner is used or iterm = 2iterm=2 was supplied;
  • p = 0p=0, otherwise.
Constraint: lworklwreqlworklwreq, where lwreq is returned by nag_sparse_complex_gen_basic_setup (f11br).

Input Parameters Omitted from the MATLAB Interface

None.

Output Parameters

1:     irevcm – int64int32nag_int scalar
On intermediate exit: has the following meanings.
irevcm = -1irevcm=-1
The calling program must compute the matrix-vector product v = AHuv=AHu, where uu and vv are stored in u and v, respectively; RGMRES, CGS and Bi-CGSTAB() methods return irevcm = -1irevcm=-1 only if the matrix norm A1A1 or AA is estimated internally using Higham's method. This can only happen if iterm = 1iterm=1 in nag_sparse_complex_gen_basic_setup (f11br).
irevcm = 1irevcm=1
The calling program must compute the matrix-vector product v = Auv=Au, where uu and vv are stored in u and v, respectively.
irevcm = 2irevcm=2
The calling program must solve the preconditioning equation Mv = uMv=u, where uu and vv are stored in u and v, respectively.
irevcm = 3irevcm=3
Monitoring step: the solution and residual at the current iteration are returned in the arrays u and v, respectively. No action by the calling program is required. nag_sparse_complex_gen_basic_diag (f11bt) can be called at this step to return additional information.
On final exit: irevcm = 4irevcm=4: nag_sparse_complex_gen_basic_solver (f11bs) has completed its tasks. The value of ifail determines whether the iteration has been successfully completed, errors have been detected or the calling program has requested termination.
2:     u( : :) – complex array
Note: the dimension of the array u must be at least nn.
On intermediate exit: the returned value of irevcm determines the contents of u as follows.
If irevcm = -1irevcm=-1, 11 or 22, u holds the vector uu on which the operation specified by irevcm is to be carried out.
If irevcm = 3irevcm=3, u holds the current iterate of the solution vector.
On final exit: if ifail = 3ifail=3 or ifail < 0ifail<0, the array u is unchanged from the initial entry to nag_sparse_complex_gen_basic_solver (f11bs).
If ifail = 1ifail=1, the array u is unchanged from the last entry to nag_sparse_complex_gen_basic_solver (f11bs).
Otherwise, u holds the last available iterate of the solution of the system of equations, for all returned values of ifail.
3:     v( : :) – complex array
Note: the dimension of the array v must be at least nn.
On intermediate exit: if irevcm = 3irevcm=3, v holds the current iterate of the residual vector. Note that this is an approximation to the true residual vector. Otherwise, it does not contain any useful information.
On final exit: if ifail = 3ifail=3 or ifail < 0ifail<0, the array v is unchanged from the initial entry to nag_sparse_complex_gen_basic_solver (f11bs).
If ifail = 1ifail=1, the array v is unchanged from the last entry to nag_sparse_complex_gen_basic_solver (f11bs).
If ifail = 0ifail=0 or 22, the array v contains the true residual vector of the system of equations (see also Section [Error Indicators and Warnings]).
Otherwise, v stores the last available iterate of the residual vector unless ifail = 8ifail=8 is returned on last entry, in which case v is set to 0.00.0.
4:     work(lwork) – complex array
5:     ifail – int64int32nag_int scalar
ifail = 0ifail=0 unless the function detects an error (see [Error Indicators and Warnings]).

Error Indicators and Warnings

Errors or warnings detected by the function:

Cases prefixed with W are classified as warnings and do not generate an error of type NAG:error_n. See nag_issue_warnings.

  ifail = iifail=-i
If ifail = iifail=-i, parameter ii had an illegal value on entry. The parameters are numbered as follows:
1: irevcm, 2: u, 3: v, 4: wgt, 5: work, 6: lwork, 7: ifail.
W ifail = 1ifail=1
nag_sparse_complex_gen_basic_solver (f11bs) has been called again after returning the termination code irevcm = 4irevcm=4. No further computation has been carried out and all input data and data stored for access by nag_sparse_complex_gen_basic_diag (f11bt) have remained unchanged.
W ifail = 2ifail=2
The required accuracy could not be obtained. However, nag_sparse_complex_gen_basic_solver (f11bs) has terminated with reasonable accuracy: the last iterate of the residual satisfied the termination criterion but the exact residual r = bAxr=b-Ax, did not. After the first occurrence of this situation, the iteration was restarted once, but nag_sparse_complex_gen_basic_solver (f11bs) could not improve on the accuracy. 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. You should call nag_sparse_complex_gen_basic_diag (f11bt) to check the values of the left- and right-hand sides of the termination condition.
  ifail = 3ifail=3
nag_sparse_complex_gen_basic_setup (f11br) was either not called before calling nag_sparse_complex_gen_basic_solver (f11bs) or it returned an error. The arguments u and v remain unchanged.
W ifail = 4ifail=4
The calling program requested a tidy termination before the solution had converged. The arrays u and v return the last iterates available of the solution and of the residual vector, respectively.
W ifail = 5ifail=5
The solution did not converge within the maximum number of iterations allowed. The arrays u and v return the last iterates available of the solution and of the residual vector, respectively.
  ifail = 6ifail=6
Algorithmic breakdown. The last available iterates of the solution and residuals are returned, although it is possible that they are completely inaccurate.
W ifail = 8ifail=8
The calling program requested an immediate termination. However, the array u returns the last iterate of the solution, the array v returns the last iterate of the residual vector, for the CGS and TFQMR methods only.
  ifail = 10ifail=10
User-supplied weights are to be used, but all elements of the array wgt are zero.

Accuracy

On completion, i.e., irevcm = 4irevcm=4 on exit, the arrays u and v will return the solution and residual vectors, xkxk and rk = bAxkrk=b-Axk, respectively, at the kkth iteration, the last iteration performed, unless an immediate termination was requested.
On successful completion, the termination criterion is satisfied to within the user-specified tolerance, as described in nag_sparse_complex_gen_basic_setup (f11br). The computed values of the left- and right-hand sides of the termination criterion selected can be obtained by a call to nag_sparse_complex_gen_basic_diag (f11bt).

Further Comments

The number of operations carried out by nag_sparse_complex_gen_basic_solver (f11bs) for each iteration is likely to be principally determined by the computation of the matrix-vector products v = Auv=Au and by the solution of the preconditioning equation Mv = uMv=u in the calling program. Each of these operations is carried out once every iteration.
The number of the remaining operations in nag_sparse_complex_gen_basic_solver (f11bs) for each iteration is approximately proportional to nn.
The number of iterations required to achieve a prescribed accuracy cannot be easily determined at the onset, as it can depend dramatically on the conditioning and spectrum of the preconditioned matrix of the coefficients A = M1AA-=M-1A (RGMRES, CGS and TFQMR methods) or A = AM1A-=AM-1 (Bi-CGSTAB() method).
Additional matrix-vector products are required for the computation of A1A1 or AA, when this has not been supplied to nag_sparse_complex_gen_basic_setup (f11br) and is required by the termination criterion employed.
If the termination criterion rkp τ (bp + Ap × xkp) rkp τ ( bp+Ap × xkp )  is used (see Section [Description] in (f11br)) and x0xkx0xk, then the required accuracy cannot be obtained due to loss of significant digits. The iteration is restarted automatically at some suitable point: nag_sparse_complex_gen_basic_solver (f11bs) sets x0 = xkx0=xk and the computation begins again. For particularly badly scaled problems, more than one restart may be necessary. This does not apply to the RGMRES method which, by its own nature, self-restarts every super-iteration. Naturally, restarting adds to computational costs: it is recommended that the iteration should start from a value x0x0 which is as close to the true solution x~ as can be estimated. Otherwise, the iteration should start from x0 = 0x0=0.

Example

function nag_sparse_complex_gen_basic_solver_example
method = 'TFQMR   ';
precon = 'P';
n = 8;
m = int64(1);
tol = 1e-08;
maxitn = int64(20);
anorm = 0;
sigmax = 0;
monit = int64(2);
lwork = int64(6000);
nz=24;
a=complex(zeros(10000,1));
a(1:24)=[ 2 + 1.i,
          -1 + 1.i,
          1 - 3.i,
          4 + 7.i,
          -3 + 0.i,
          2 + 4.i,
          -7 - 5.i,
          2 + 1.i,
          3 + 2.i,
          -4 + 2.i,
          0 + 1.i,
          5 - 3.i,
          -1 + 2.i,
          8 + 6.i,
          -3 - 4.i,
          -6 - 2.i,
          5 - 2.i,
          2 + 0.i,
          0 - 5.i,
          -1 + 5.i,
          6 + 2.i,
          -1 + 4.i,
          2 + 0.i,
          3 + 3.i];
irow=zeros(10000,1,'int64');
irow(1:24) = [int64(1),1,1,2,2,2,3,3,4,4,4,4,5,5,5,6,6,6,7,7,7,8,8,8];
icol=zeros(10000,1,'int64');
icol(1:24) = [int64(1),4,8,1,2,5,3,6,1,3,4,7,2,5,7,1,3,6,3,5,7,2,6,8];
ipivp = zeros(n,1,'int64');
ipivq = zeros(n,1,'int64');
irevcm = int64(0);
lfill = int64(0);
dtol = 0;
milu = 'N';
wgt = zeros(8,1);
u = complex(zeros(8,1));
v = [ 7 + 11.i;
      1 + 24.i;
      -13 - 18.i;
      -10 + 3.i;
      23 + 14.i;
      17 - 7.i;
      15 - 3.i;
      -3 + 20.i];
[a, irow, icol, ipivp, ipivq, istr, idiag, nnzc, npivm, ifail] = ...
    nag_sparse_complex_gen_precon_ilu(int64(nz), a, irow, icol, lfill, ...
    dtol, milu, ipivp, ipivq);

[lwreq, work, ifail] = ...
     nag_sparse_complex_gen_basic_setup(method, precon, int64(n), m, tol, ...
     maxitn, anorm, sigmax, monit, lwork, 'norm_p', '1');

while (irevcm ~= 4)
  [irevcm, u, v, work, ifail] = ...
    nag_sparse_complex_gen_basic_solver(irevcm, u, v, wgt, work);

  if (irevcm == -1)
    [v, ifail] = ...
    nag_sparse_complex_gen_matvec('T', a(1:nz), irow(1:nz), icol(1:nz), ...
    'N', u);
  elseif (irevcm == 1)
    [v, ifail] = ...
    nag_sparse_complex_gen_matvec('N', a(1:nz), irow(1:nz), icol(1:nz), ...
    'N', u);
  elseif (irevcm == 2)
    [v, ifail] = ...
    nag_sparse_complex_gen_precon_ilu_solve('N', a, irow, icol, ipivp, ...
    ipivq, istr, idiag, 'N', u);
  elseif (irevcm == 3)
    [itn, stplhs, stprhs, anorm, sigmax, work, ifail] = ...
    nag_sparse_complex_gen_basic_diag(work);
    fprintf('\nMonitoring at iteration number %d\nresidual norm: %14.4e\n', ...
            itn, stplhs);
    fprintf('\n   Solution Vector\n');
    for i = 1:n
      fprintf('%+16.4e + %+16.4eI\n', real(u(i)), imag(u(i)));
    end
    fprintf('\n   Residual Vector\n');
    for i = 1:n
      fprintf('%+16.4e + %+16.4eI\n', real(v(i)), imag(v(i)));
    end
  end
end
% Get information about the computation
[itn, stplhs, stprhs, anorm, sigmax, work, ifail] = ...
    nag_sparse_complex_gen_basic_diag(work);
fprintf('\nNumber of iterations for convergence: %d\n', itn);
fprintf('Residual norm:                           %14.4e\n', stplhs);
fprintf('Right-hand side of termination criteria: %14.4e\n', stprhs);
fprintf('i-norm of matrix a:                      %14.4e\n', anorm);
fprintf('\n   Solution Vector\n');
for i = 1:n
  fprintf('%+16.4e + %+16.4eI\n', real(u(i)), imag(u(i)));
end
fprintf('\n   Residual Vector\n');
for i = 1:n
  fprintf('%+16.4e + %+16.4eI\n', real(v(i)), imag(v(i)));
end
 

Monitoring at iteration number 2
residual norm:     8.2345e+01

   Solution Vector
     +6.9055e-01 +      +1.4236e+00I
     +7.3931e-02 +      -1.1880e+00I
     +1.4778e+00 +      +4.7846e-01I
     +5.6572e+00 +      -3.0786e+00I
     +1.4243e+00 +      -1.1246e+00I
     +1.0374e-01 +      +1.9740e+00I
     +4.4985e-01 +      -1.2715e+00I
     +2.5704e+00 +      +1.7578e+00I

   Residual Vector
     +1.7772e+00 +      +4.6797e+00I
     +1.0774e+00 +      +6.4600e+00I
     -3.2812e+00 +      -1.1314e+01I
     -3.8698e+00 +      -1.6438e+00I
     +8.9912e+00 +      +1.1100e+01I
     +9.7428e+00 +      -4.6218e-01I
     +3.1668e+00 +      +2.8721e+00I
     -1.0323e+01 +      +1.5837e+00I

Number of iterations for convergence: 4
Residual norm:                               2.7448e-11
Right-hand side of termination criteria:     8.9100e-06
i-norm of matrix a:                          2.7000e+01

   Solution Vector
     +1.0000e+00 +      +1.0000e+00I
     +2.0000e+00 +      -1.0000e+00I
     +3.0000e+00 +      +1.0000e+00I
     +4.0000e+00 +      -1.0000e+00I
     +3.0000e+00 +      -1.0000e+00I
     +2.0000e+00 +      +1.0000e+00I
     +1.0000e+00 +      -1.0000e+00I
     -8.1379e-13 +      +3.0000e+00I

   Residual Vector
     -6.1817e-13 +      -1.9504e-12I
     +1.7746e-12 +      -4.7180e-12I
     +1.3003e-12 +      +1.9043e-12I
     +1.6733e-12 +      -8.1535e-13I
     -2.3448e-12 +      -1.4975e-12I
     -1.7870e-12 +      +6.5459e-13I
     -7.8337e-13 +      +1.1635e-12I
     +1.9051e-12 +      +2.5580e-12I

function f11bs_example
method = 'TFQMR   ';
precon = 'P';
n = 8;
m = int64(1);
tol = 1e-08;
maxitn = int64(20);
anorm = 0;
sigmax = 0;
monit = int64(2);
lwork = int64(6000);
nz=24;
a=complex(zeros(10000,1));
a(1:24)=[ 2 + 1.i,
          -1 + 1.i,
          1 - 3.i,
          4 + 7.i,
          -3 + 0.i,
          2 + 4.i,
          -7 - 5.i,
          2 + 1.i,
          3 + 2.i,
          -4 + 2.i,
          0 + 1.i,
          5 - 3.i,
          -1 + 2.i,
          8 + 6.i,
          -3 - 4.i,
          -6 - 2.i,
          5 - 2.i,
          2 + 0.i,
          0 - 5.i,
          -1 + 5.i,
          6 + 2.i,
          -1 + 4.i,
          2 + 0.i,
          3 + 3.i];
irow=zeros(10000,1,'int64');
irow(1:24) = [int64(1),1,1,2,2,2,3,3,4,4,4,4,5,5,5,6,6,6,7,7,7,8,8,8];
icol=zeros(10000,1,'int64');
icol(1:24) = [int64(1),4,8,1,2,5,3,6,1,3,4,7,2,5,7,1,3,6,3,5,7,2,6,8];
ipivp = zeros(n,1,'int64');
ipivq = zeros(n,1,'int64');
irevcm = int64(0);
lfill = int64(0);
dtol = 0;
milu = 'N';
wgt = zeros(8,1);
u = complex(zeros(8,1));
v = [ 7 + 11.i;
      1 + 24.i;
      -13 - 18.i;
      -10 + 3.i;
      23 + 14.i;
      17 - 7.i;
      15 - 3.i;
      -3 + 20.i];
[a, irow, icol, ipivp, ipivq, istr, idiag, nnzc, npivm, ifail] = ...
    f11dn(int64(nz), a, irow, icol, lfill, dtol, milu, ipivp, ipivq);

[lwreq, work, ifail] = ...
     f11br(method, precon, int64(n), m, tol, maxitn, anorm, sigmax, monit, lwork, 'norm_p', '1');

while (irevcm ~= 4)
  [irevcm, u, v, work, ifail] = f11bs(irevcm, u, v, wgt, work);

  if (irevcm == -1)
    [v, ifail] = f11xn('T', a(1:nz), irow(1:nz), icol(1:nz), 'N', u);
  elseif (irevcm == 1)
    [v, ifail] = f11xn('N', a(1:nz), irow(1:nz), icol(1:nz), 'N', u);
  elseif (irevcm == 2)
    [v, ifail] = f11dp('N', a, irow, icol, ipivp, ipivq, istr, idiag, 'N', u);
  elseif (irevcm == 3)
    [itn, stplhs, stprhs, anorm, sigmax, work, ifail] = f11bt(work);
    fprintf('\nMonitoring at iteration number %d\nresidual norm: %14.4e\n', itn, stplhs);
    fprintf('\n   Solution Vector\n');
    for i = 1:n
      fprintf('%+16.4e + %+16.4eI\n', real(u(i)), imag(u(i)));
    end
    fprintf('\n   Residual Vector\n');
    for i = 1:n
      fprintf('%+16.4e + %+16.4eI\n', real(v(i)), imag(v(i)));
    end
  end
end
% Get information about the computation
[itn, stplhs, stprhs, anorm, sigmax, work, ifail] = f11bt(work);
fprintf('\nNumber of iterations for convergence: %d\n', itn);
fprintf('Residual norm:                           %14.4e\n', stplhs);
fprintf('Right-hand side of termination criteria: %14.4e\n', stprhs);
fprintf('i-norm of matrix a:                      %14.4e\n', anorm);
fprintf('\n   Solution Vector\n');
for i = 1:n
  fprintf('%+16.4e + %+16.4eI\n', real(u(i)), imag(u(i)));
end
fprintf('\n   Residual Vector\n');
for i = 1:n
  fprintf('%+16.4e + %+16.4eI\n', real(v(i)), imag(v(i)));
end
 

Monitoring at iteration number 2
residual norm:     8.2345e+01

   Solution Vector
     +6.9055e-01 +      +1.4236e+00I
     +7.3931e-02 +      -1.1880e+00I
     +1.4778e+00 +      +4.7846e-01I
     +5.6572e+00 +      -3.0786e+00I
     +1.4243e+00 +      -1.1246e+00I
     +1.0374e-01 +      +1.9740e+00I
     +4.4985e-01 +      -1.2715e+00I
     +2.5704e+00 +      +1.7578e+00I

   Residual Vector
     +1.7772e+00 +      +4.6797e+00I
     +1.0774e+00 +      +6.4600e+00I
     -3.2812e+00 +      -1.1314e+01I
     -3.8698e+00 +      -1.6438e+00I
     +8.9912e+00 +      +1.1100e+01I
     +9.7428e+00 +      -4.6218e-01I
     +3.1668e+00 +      +2.8721e+00I
     -1.0323e+01 +      +1.5837e+00I

Number of iterations for convergence: 4
Residual norm:                               2.5802e-11
Right-hand side of termination criteria:     8.9100e-06
i-norm of matrix a:                          2.7000e+01

   Solution Vector
     +1.0000e+00 +      +1.0000e+00I
     +2.0000e+00 +      -1.0000e+00I
     +3.0000e+00 +      +1.0000e+00I
     +4.0000e+00 +      -1.0000e+00I
     +3.0000e+00 +      -1.0000e+00I
     +2.0000e+00 +      +1.0000e+00I
     +1.0000e+00 +      -1.0000e+00I
     -5.4061e-13 +      +3.0000e+00I

   Residual Vector
     -4.3165e-13 +      -1.9966e-12I
     +1.7355e-12 +      -4.6327e-12I
     +6.8212e-13 +      +2.2240e-12I
     +1.4957e-12 +      -5.4090e-13I
     -1.8190e-12 +      -2.1512e-12I
     -1.9149e-12 +      +1.3678e-13I
     -1.0107e-12 +      +1.1386e-12I
     +1.2381e-12 +      +2.6539e-12I


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–2013