PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_sparse_complex_gen_precon_ssor_solve (f11dr)
Purpose
nag_sparse_complex_gen_precon_ssor_solve (f11dr) solves a system of linear equations involving the preconditioning matrix corresponding to SSOR applied to a complex sparse non-Hermitian matrix, represented in coordinate storage format.
Syntax
[
x,
ifail] = f11dr(
trans,
a,
irow,
icol,
rdiag,
omega,
check,
y, 'n',
n, 'nz',
nz)
[
x,
ifail] = nag_sparse_complex_gen_precon_ssor_solve(
trans,
a,
irow,
icol,
rdiag,
omega,
check,
y, 'n',
n, 'nz',
nz)
Description
nag_sparse_complex_gen_precon_ssor_solve (f11dr) solves a system of linear equations
according to the value of the argument
trans, where the matrix
corresponds to symmetric successive-over-relaxation (SSOR)
Young (1971) applied to a linear system
, where
is a complex sparse non-Hermitian matrix stored in coordinate storage (CS) format (see
Coordinate storage (CS) format in the F11 Chapter Introduction).
In the definition of given above is the diagonal part of , is the strictly lower triangular part of , is the strictly upper triangular part of , and is a user-defined relaxation parameter.
It is envisaged that a common use of
nag_sparse_complex_gen_precon_ssor_solve (f11dr) will be to carry out the preconditioning step required in the application of
nag_sparse_complex_gen_basic_solver (f11bs) to sparse linear systems. For an illustration of this use of
nag_sparse_complex_gen_precon_ssor_solve (f11dr) see the example program given in
Example.
nag_sparse_complex_gen_precon_ssor_solve (f11dr) is also used for this purpose by the Black Box function
nag_sparse_complex_gen_solve_jacssor (f11ds).
References
Young D (1971) Iterative Solution of Large Linear Systems Academic Press, New York
Parameters
Compulsory Input Parameters
- 1:
– string (length ≥ 1)
-
Specifies whether or not the matrix
is transposed.
- is solved.
- is solved.
Constraint:
or .
- 2:
– complex array
-
The nonzero elements in the matrix
, 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_complex_gen_sort (f11zn) may be used to order the elements in this way.
- 3:
– int64int32nag_int array
- 4:
– int64int32nag_int array
-
The row and column indices of the nonzero elements supplied in
a.
Constraints:
irow and
icol must satisfy the following constraints (which may be imposed by a call to
nag_sparse_complex_gen_sort (f11zn)):
- and , for ;
- either or both and , for .
- 5:
– complex array
-
The elements of the diagonal matrix , where is the diagonal part of .
- 6:
– double scalar
-
The relaxation parameter .
Constraint:
.
- 7:
– string (length ≥ 1)
-
Specifies whether or not the CS representation of the matrix
should be checked.
- Checks are carried on the values of n, nz, irow, icol and omega.
- None of these checks are carried out.
Constraint:
or .
- 8:
– complex array
-
The right-hand side vector .
Optional Input Parameters
- 1:
– int64int32nag_int scalar
-
Default:
the dimension of the arrays
rdiag,
y. (An error is raised if these dimensions are not equal.)
, the order of the matrix .
Constraint:
.
- 2:
– 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 matrix .
Constraint:
.
Output Parameters
- 1:
– complex array
-
The solution vector .
- 2:
– int64int32nag_int scalar
unless the function detects an error (see
Error Indicators and Warnings).
Error Indicators and Warnings
Errors or warnings detected by the function:
-
-
On entry, | or , |
or | or . |
-
-
On entry, | , |
or | , |
or | , |
or | omega lies outside the interval , |
-
-
On entry, the arrays
irow and
icol fail to satisfy the following constraints:
- and , for ;
- or and , for .
Therefore a nonzero element has been supplied which does not lie in the matrix
, is out of order, or has duplicate row and column indices. Call
nag_sparse_complex_gen_sort (f11zn) to reorder and sum or remove duplicates.
-
-
On entry, the matrix has a zero diagonal element. The SSOR preconditioner is not appropriate for this problem.
-
An unexpected error has been triggered by this routine. Please
contact
NAG.
-
Your licence key may have expired or may not have been installed correctly.
-
Dynamic memory allocation failed.
Accuracy
If
the computed solution
is the exact solution of a perturbed system of equations
, where
is a modest linear function of
, and
is the
machine precision. An equivalent result holds when
.
Further Comments
Timing
The time taken for a call to
nag_sparse_complex_gen_precon_ssor_solve (f11dr) is proportional to
nz.
Use of check
It is expected that a common use of
nag_sparse_complex_gen_precon_ssor_solve (f11dr) will be to carry out the preconditioning step required in the application of
nag_sparse_complex_gen_basic_solver (f11bs) to sparse linear systems. In this situation
nag_sparse_complex_gen_precon_ssor_solve (f11dr) is likely to be called many times with the same matrix
. In the interests of both reliability and efficiency, you are recommended to set
for the first of such calls, and
for all subsequent calls.
Example
This example solves a complex sparse linear system of equations
using RGMRES with SSOR preconditioning.
The RGMRES algorithm itself is implemented by the reverse communication function
nag_sparse_complex_gen_basic_solver (f11bs), 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.
- If , a matrix-vector product is required. This is implemented by a call to nag_sparse_complex_gen_matvec (f11xn).
- If , a conjugate transposed matrix-vector product is required in the estimation of the norm of . This is implemented by a call to nag_sparse_complex_gen_matvec (f11xn).
- If , a solution of the preconditioning equation is required. This is achieved by a call to nag_sparse_complex_gen_precon_ssor_solve (f11dr).
- If , nag_sparse_complex_gen_basic_solver (f11bs) has completed its tasks. Either the iteration has terminated, or an error condition has arisen.
For further details see the function document for
nag_sparse_complex_gen_basic_solver (f11bs).
Open in the MATLAB editor:
f11dr_example
function f11dr_example
fprintf('f11dr example results\n\n');
n = int64(5);
m = int64(2);
nz = int64(16);
sparseA = [ 2 3 1 1;
1 -1 1 2;
-1 0 1 4;
0 2 2 2;
-2 1 2 3;
1 0 2 5;
0 -1 3 1;
5 4 3 3;
3 -1 3 4;
1 0 3 5;
-2 2 4 1;
-3 1 4 4;
0 3 4 5;
4 -2 5 2;
-2 0 5 3;
-6 1 5 5];
a(1:nz) = sparseA(:,1) + i*sparseA(:,2);
irow(1:nz) = int64(sparseA(:,3));
icol(1:nz) = int64(sparseA(:,4));
b = [ -3 + 3i;
-11 + 5i;
23 + 48i;
-41 + 2i;
-28 - 31i];
method = 'CGS';
precon = 'P';
norm_p = 'I';
tol = 1e-10;
maxitn = int64(1000);
anorm = 0;
sigmax = 0;
monit = int64(0);
lwork = max([121+n*(3+m)+m*(m+5),120+7*n,120+(2*n+m)*(m+2)+2*n,120+10*n]);
[lwreq, work, ifail] = f11br( ...
method, precon, n, m, tol, maxitn, anorm, ...
sigmax, monit, lwork, 'norm_p', norm_p);
rdiag = zeros(n, 1);
if precon == 'P' || precon == 'p'
dcount = zeros(n, 1, 'int64');
for j = 1:nz
if irow(j) == icol(j)
dcount(irow(j)) = dcount(irow(j)) + 1;
if a(j) ~= 0
rdiag(irow(j)) = 1/a(j);
else
error('Matrix has a zero diagonal element');
end
end
end
for j = 1:n
if dcount(j)==0
error('Matrix has a missing diagonal element');
elseif dcount(j) > 1
error('Matrix has a multiple diagonal element');
end
end
end
omega = 1.4;
ckdrf = 'C';
irevcm = int64(0);
wgt = zeros(n, 1);
x = complex(zeros(n, 1));
ckxnf = 'C';
while (irevcm ~= 4)
[irevcm, x, b, work, ifail] = f11bs( ...
irevcm, x, b, wgt, work);
if (irevcm == -1)
[b, ifail] = f11xn( ...
'T', a, irow, icol, ckxnf, x);
ckxnf = 'N';
elseif (irevcm == 1)
[b, ifail] = f11xn( ...
'N', a, irow, icol, ckxnf, x);
ckxnf = 'N';
elseif (irevcm == 2)
[b, ifail] = f11dr( ...
'N', a, irow, icol, rdiag, omega, ckdrf, x);
ckdrf = 'N';
end
end
if ifail == 0
[itn, stplhs, stprhs, anorm, sigmax, work, ifail] = ...
f11bt(work);
fprintf('Converged in %d iterations\n', itn);
fprintf('Matrix norm = %16.3e\n', anorm);
fprintf('Final residual norm = %16.4e\n\n', stplhs);
disp('Solution');
disp(x);
end
f11dr example results
Converged in 5 iterations
Matrix norm = 1.500e+01
Final residual norm = 4.6185e-14
Solution
1.0000 + 2.0000i
2.0000 + 3.0000i
3.0000 + 4.0000i
4.0000 + 5.0000i
5.0000 + 6.0000i
PDF version (NAG web site
, 64-bit version, 64-bit version)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015