nag_sparse_herm_sol (f11jsc) (PDF version)
f11 Chapter Contents
f11 Chapter Introduction
NAG Library Manual

NAG Library Function Document

nag_sparse_herm_sol (f11jsc)

 Contents

    1  Purpose
    7  Accuracy

1  Purpose

nag_sparse_herm_sol (f11jsc) solves a complex sparse Hermitian system of linear equations, represented in symmetric coordinate storage format, using a conjugate gradient or Lanczos method, without preconditioning, with Jacobi or with SSOR preconditioning.

2  Specification

#include <nag.h>
#include <nagf11.h>
void  nag_sparse_herm_sol (Nag_SparseSym_Method method, Nag_SparseSym_PrecType precon, Integer n, Integer nnz, const Complex a[], const Integer irow[], const Integer icol[], double omega, const Complex b[], double tol, Integer maxitn, Complex x[], double *rnorm, Integer *itn, double rdiag[], NagError *fail)

3  Description

nag_sparse_herm_sol (f11jsc) solves a complex sparse Hermitian linear system of equations
Ax=b,  
using a preconditioned conjugate gradient method (see Barrett et al. (1994)), or a preconditioned Lanczos method based on the algorithm SYMMLQ (see Paige and Saunders (1975)). The conjugate gradient method is more efficient if A is positive definite, but may fail to converge for indefinite matrices. In this case the Lanczos method should be used instead. For further details see Barrett et al. (1994).
nag_sparse_herm_sol (f11jsc) allows the following choices for the preconditioner:
For incomplete Cholesky (IC) preconditioning see nag_sparse_herm_chol_sol (f11jqc).
The matrix A is represented in symmetric coordinate storage (SCS) format (see Section 2.1.2 in the f11 Chapter Introduction) in the arrays a, irow and icol. The array a holds the nonzero entries in the lower triangular part of the matrix, while irow and icol hold the corresponding row and column indices.

4  References

Barrett R, Berry M, Chan T F, Demmel J, Donato J, Dongarra J, Eijkhout V, Pozo R, Romine C and Van der Vorst H (1994) Templates for the Solution of Linear Systems: Building Blocks for Iterative Methods SIAM, Philadelphia
Paige C C and Saunders M A (1975) Solution of sparse indefinite systems of linear equations SIAM J. Numer. Anal. 12 617–629
Young D (1971) Iterative Solution of Large Linear Systems Academic Press, New York

5  Arguments

1:     method Nag_SparseSym_MethodInput
On entry: specifies the iterative method to be used.
method=Nag_SparseSym_CG
Conjugate gradient method.
method=Nag_SparseSym_SYMMLQ
Lanczos method (SYMMLQ).
Constraint: method=Nag_SparseSym_CG or Nag_SparseSym_SYMMLQ.
2:     precon Nag_SparseSym_PrecTypeInput
On entry: specifies the type of preconditioning to be used.
precon=Nag_SparseSym_NoPrec
No preconditioning.
precon=Nag_SparseSym_JacPrec
Jacobi.
precon=Nag_SparseSym_SSORPrec
Symmetric successive-over-relaxation (SSOR).
Constraint: precon=Nag_SparseSym_NoPrec, Nag_SparseSym_JacPrec or Nag_SparseSym_SSORPrec.
3:     n IntegerInput
On entry: n, the order of the matrix A.
Constraint: n1.
4:     nnz IntegerInput
On entry: the number of nonzero elements in the lower triangular part of the matrix A.
Constraint: 1nnzn×n+1/2.
5:     a[nnz] const ComplexInput
On entry: the nonzero elements of 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_herm_sort (f11zpc) may be used to order the elements in this way.
6:     irow[nnz] const IntegerInput
7:     icol[nnz] const IntegerInput
On entry: 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_herm_sort (f11zpc)):
  • 1irow[i]n and 1icol[i]irow[i], for i=0,1,,nnz-1;
  • irow[i-1]<irow[i] or irow[i-1]=irow[i] and icol[i-1]<icol[i], for i=1,2,,nnz-1.
8:     omega doubleInput
On entry: if precon=Nag_SparseSym_SSORPrec, omega is the relaxation parameter ω to be used in the SSOR method. Otherwise omega need not be initialized.
Constraint: 0.0<omega<2.0.
9:     b[n] const ComplexInput
On entry: the right-hand side vector b.
10:   tol doubleInput
On entry: 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ε,10ε,nε is used, where ε is the machine precision. Otherwise τ=maxtol,10ε,nε is used.
Constraint: tol<1.0.
11:   maxitn IntegerInput
On entry: the maximum number of iterations allowed.
Constraint: maxitn1.
12:   x[n] ComplexInput/Output
On entry: an initial approximation to the solution vector x.
On exit: an improved approximation to the solution vector x.
13:   rnorm double *Output
On exit: the final value of the residual norm rk, where k is the output value of itn.
14:   itn Integer *Output
On exit: the number of iterations carried out.
15:   rdiag[n] doubleOutput
On exit: the elements of the diagonal matrix D-1, where D is the diagonal part of A. Note that since A is Hermitian the elements of D-1 are necessarily real.
16:   fail NagError *Input/Output
The NAG error argument (see Section 2.7 in How to Use the NAG Library and its Documentation).

6  Error Indicators and Warnings

NE_ACCURACY
The required accuracy could not be obtained. However, a reasonable accuracy has been achieved and further iterations could not improve the result.
NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 2.3.1.2 in How to Use the NAG Library and its Documentation for further information.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_COEFF_NOT_POS_DEF
The matrix of the coefficients a appears not to be positive definite. The computation cannot continue.
NE_CONVERGENCE
The solution has not converged after value iterations.
NE_INT
On entry, maxitn=value.
Constraint: maxitn1.
On entry, n=value.
Constraint: n1.
On entry, nnz=value.
Constraint: nnz1.
NE_INT_2
On entry, nnz=value and n=value.
Constraint: nnzn×n+1/2 
NE_INTERNAL_ERROR
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact NAG for assistance.
An unexpected error has been triggered by this function. Please contact NAG.
See Section 2.7.6 in How to Use the NAG Library and its Documentation for further information.
A serious error, code value, has occurred in an internal call to nag_sparse_herm_basic_solver (f11gsc). Check all function calls and array sizes. Seek expert help.
A serious error, code value, has occurred in an internal call to value. Check all function calls and array sizes. Seek expert help.
NE_INVALID_SCS
On entry, I=value, icol[I-1]=value and irow[I-1]=value.
Constraint: icol[I-1]1 and icol[I-1]irow[I-1].
On entry, i=value, irow[i-1]=value and n=value.
Constraint: irow[i-1]1 and irow[i-1]n.
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 2.7.5 in How to Use the NAG Library and its Documentation for further information.
NE_NOT_STRICTLY_INCREASING
On entry, a[i-1] is out of order: i=value.
On entry, the location (irow[I-1],icol[I-1]) is a duplicate: I=value. Consider calling nag_sparse_herm_sort (f11zpc) to reorder and sum or remove duplicates.
NE_PRECOND_NOT_POS_DEF
The preconditioner appears not to be positive definite. The computation cannot continue.
NE_REAL
On entry, omega=value.
Constraint: 0.0<omega<2.0.
On entry, tol=value.
Constraint: tol<1.0.
NE_ZERO_DIAG_ELEM
The matrix A has a non-real diagonal entry in row value.
The matrix A has a zero diagonal entry in row value.
The matrix A has no diagonal entry in row value.

7  Accuracy

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

8  Parallelism and Performance

nag_sparse_herm_sol (f11jsc) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
nag_sparse_herm_sol (f11jsc) makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the x06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this function. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9  Further Comments

The time taken by nag_sparse_herm_sol (f11jsc) for each iteration is roughly proportional to nnz. One iteration with the Lanczos method (SYMMLQ) requires a slightly larger number of operations than one iteration with the conjugate gradient method.
The number of iterations required to achieve a prescribed accuracy cannot easily be determined a priori, as it can depend dramatically on the conditioning and spectrum of the preconditioned matrix of the coefficients A-=M-1A.

10  Example

This example solves a complex sparse Hermitian positive definite system of equations using the conjugate gradient method, with SSOR preconditioning.

10.1  Program Text

Program Text (f11jsce.c)

10.2  Program Data

Program Data (f11jsce.d)

10.3  Program Results

Program Results (f11jsce.r)


nag_sparse_herm_sol (f11jsc) (PDF version)
f11 Chapter Contents
f11 Chapter Introduction
NAG Library Manual

© The Numerical Algorithms Group Ltd, Oxford, UK. 2016