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_real_symm_precon_ssor_solve (f11jd)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_sparse_real_symm_precon_ssor_solve (f11jd) solves a system of linear equations involving the preconditioning matrix corresponding to SSOR applied to a real sparse symmetric matrix, represented in symmetric coordinate storage format.

Syntax

[x, ifail] = f11jd(a, irow, icol, rdiag, omega, check, y, 'n', n, 'nz', nz)
[x, ifail] = nag_sparse_real_symm_precon_ssor_solve(a, irow, icol, rdiag, omega, check, y, 'n', n, 'nz', nz)

Description

nag_sparse_real_symm_precon_ssor_solve (f11jd) solves a system of equations
Mx=y  
involving the preconditioning matrix
M=1ω2-ω D+ω L D-1 D+ω LT  
corresponding to symmetric successive-over-relaxation (SSOR) (see Young (1971)) on a linear system Ax=b, where A is a sparse symmetric matrix stored in symmetric coordinate storage (SCS) format (see Symmetric coordinate storage (SCS) format in the F11 Chapter Introduction).
In the definition of M given above D is the diagonal part of A, L is the strictly lower triangular part of A, and ω is a user-defined relaxation parameter.
It is envisaged that a common use of nag_sparse_real_symm_precon_ssor_solve (f11jd) will be to carry out the preconditioning step required in the application of nag_sparse_real_symm_basic_solver (f11ge) to sparse linear systems. For an illustration of this use of nag_sparse_real_symm_precon_ssor_solve (f11jd) see the example program given in Example. nag_sparse_real_symm_precon_ssor_solve (f11jd) is also used for this purpose by the Black Box function nag_sparse_real_symm_solve_jacssor (f11je).

References

Young D (1971) Iterative Solution of Large Linear Systems Academic Press, New York

Parameters

Compulsory Input Parameters

1:     anz – double array
The nonzero elements in the lower triangular part of the matrix A, ordered by increasing row index, and by increasing column index within each row. Multiple entries for the same row and column indices are not permitted. The function nag_sparse_real_symm_sort (f11zb) may be used to order the elements in this way.
2:     irownz int64int32nag_int array
3:     icolnz int64int32nag_int array
The row and column indices of the nonzero elements supplied in array a.
Constraints:
irow and icol must satisfy these constraints (which may be imposed by a call to nag_sparse_real_symm_sort (f11zb)):
  • 1irowin and 1icoliirowi, for i=1,2,,nz;
  • irowi-1<irowi or irowi-1=irowi and icoli-1<icoli, for i=2,3,,nz.
4:     rdiagn – double array
The elements of the diagonal matrix D-1, where D is the diagonal part of A.
5:     omega – double scalar
The relaxation parameter ω.
Constraint: 0.0<omega<2.0.
6:     check – string (length ≥ 1)
Specifies whether or not the input data should be checked.
check='C'
Checks are carried out on the values of n, nz, irow, icol and omega.
check='N'
None of these checks are carried out.
See also Use of check.
Constraint: check='C' or 'N'.
7:     yn – double array
The right-hand side vector y.

Optional Input Parameters

1:     n int64int32nag_int scalar
Default: the dimension of the arrays rdiag, y. (An error is raised if these dimensions are not equal.)
n, the order of the matrix A.
Constraint: n1.
2:     nz int64int32nag_int scalar
Default: the dimension of the arrays a, irow, icol. (An error is raised if these dimensions are not equal.)
The number of nonzero elements in the lower triangular part of A.
Constraint: 1nzn×n+1/2.

Output Parameters

1:     xn – double array
The solution vector x.
2:     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
On entry,check'C' or 'N'.
   ifail=2
On entry,n<1,
ornz<1,
ornz>n×n+1/2,
oromega lies outside the interval 0.0,2.0,
   ifail=3
On entry, the arrays irow and icol fail to satisfy the following constraints:
  • 1irowin and 1icoliirowi, for i=1,2,,nz;
  • irowi-1<irowi or irowi-1=irowi and icoli-1<icoli, for i=2,3,,nz.
Therefore a nonzero element has been supplied which does not lie in the lower triangular part of A, is out of order, or has duplicate row and column indices. Call nag_sparse_real_symm_sort (f11zb) to reorder and sum or remove duplicates.
   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

The computed solution x is the exact solution of a perturbed system of equations M+δMx=y, where
δMcnεD+ωLD-1D+ωLT,  
cn is a modest linear function of n, and ε is the machine precision.

Further Comments

Timing

The time taken for a call to nag_sparse_real_symm_precon_ssor_solve (f11jd) is proportional to nz.

Use of check

It is expected that a common use of nag_sparse_real_symm_precon_ssor_solve (f11jd) will be to carry out the preconditioning step required in the application of nag_sparse_real_symm_basic_solver (f11ge) to sparse symmetric linear systems. In this situation nag_sparse_real_symm_precon_ssor_solve (f11jd) is likely to be called many times with the same matrix M. In the interests of both reliability and efficiency, you are recommended to set check='C' for the first of such calls, and to set check='N' for all subsequent calls.

Example

This example solves a sparse symmetric linear system of equations
Ax=b,  
using the conjugate-gradient (CG) method with SSOR preconditioning.
The CG algorithm itself is implemented by the reverse communication function nag_sparse_real_symm_basic_solver (f11ge), which returns repeatedly to the calling program with various values of the argument irevcm. This argument indicates the action to be taken by the calling program.
For further details see the function document for nag_sparse_real_symm_basic_solver (f11ge).
function f11jd_example


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

% Solve Ax = b by CG with SSOR preconditioning

% Matrix A and RHS b

n    = int64(7);
nz   = int64(16);
a = zeros(nz, 1);
irow = zeros(nz, 1, 'int64');
icol = irow;
a          = [4; 1; 5; 2; 2; 3;-1; 1; 4; 1;-2; 3; 2;-1;-2; 5];
irow(1:nz) = [1; 2; 2; 3; 4; 4; 5; 5; 5; 6; 6; 6; 7; 7; 7; 7];
icol(1:nz) = [1; 1; 2; 3; 2; 4; 1; 4; 5; 2; 5; 6; 1; 2; 3; 7];

b = [15; 18; -8; 21; 11; 10; 29];

% Initialise solver
method = 'CG';
precon = 'P';
tol    = 1e-6;
maxitn = int64(100);
anorm  = 0;
sigmax = 0;
maxits = int64(10);
monit = int64(0);
[lwreq, work, ifail] = ...
  f11gd( ...
         method, precon, n, tol, maxitn, anorm, sigmax, ...
         maxits, monit, 'sigcmp', 'N', 'norm_p', 'I');

% SSOR preconditioner
% Calculate reciprocal diagonal matrix elements.
iwork = zeros(n+1, 1, 'int64');
rdiag = zeros(n, 1);
for i=1:nz
  if irow(i) == icol(i)
    iwork(irow(i)) = iwork(irow(i)) + 1;
    if a(i) ~= 0
      rdiag(irow(i)) = 1/a(i);
    else
      error('Matrix has a zero diagonal element');
    end
  end
end
for i=1:n
  if iwork(i) == 0
    error('Matrix has a missing diagonal elemen');
  elseif iwork(i) > 2
    error('Matrix has a multiple diagonal element');
  end
end

% Solver inputs
irevcm = int64(0);
x      = zeros(n, 1);
wgt    = zeros(n, 1);

% Preconditioner inputs
omega  = 1;
check  = 'C';

% Solve the linear system
while irevcm ~= 4
  [irevcm, x, b, work, ifail] = f11ge( ...
                                       irevcm, x, b, wgt, work);
  if (irevcm == 1)
    % Compute matrix vector product
    [b, ifail] = f11xe( ...
                        a, irow, icol, 'N', x);
  elseif (irevcm == 2)
    % SSOR preconditioning
    [b, ifail] = f11jd( ...
                        a, irow, icol, rdiag, omega, check, x);
  end
end

% Get information about the computation
[itn, stplhs, stprhs, anorm, sigmax, its, sigerr, ifail] = ...
f11gf(work);
fprintf('Converged in %d iterations\n', itn);
fprintf('Final residual norm = %14.4e\n\n', stplhs);
disp('Solution');
disp(x);


f11jd example results

Converged in 6 iterations
Final residual norm =     7.1054e-15

Solution
    1.0000
    2.0000
    3.0000
    4.0000
    5.0000
    6.0000
    7.0000


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