NAG CL Interface
f07hpc (zpbsvx)

Settings help

CL Name Style:


1 Purpose

f07hpc uses the Cholesky factorization
A=UHU   or   A=LLH  
to compute the solution to a complex system of linear equations
AX=B ,  
where A is an n×n Hermitian positive definite band matrix of bandwidth (2kd+1) and X and B are n×r matrices. Error bounds on the solution and a condition estimate are also provided.

2 Specification

#include <nag.h>
void  f07hpc (Nag_OrderType order, Nag_FactoredFormType fact, Nag_UploType uplo, Integer n, Integer kd, Integer nrhs, Complex ab[], Integer pdab, Complex afb[], Integer pdafb, Nag_EquilibrationType *equed, double s[], Complex b[], Integer pdb, Complex x[], Integer pdx, double *rcond, double ferr[], double berr[], NagError *fail)
The function may be called by the names: f07hpc, nag_lapacklin_zpbsvx or nag_zpbsvx.

3 Description

f07hpc performs the following steps:
  1. 1.If fact=Nag_EquilibrateAndFactor, real diagonal scaling factors, DS , are computed to equilibrate the system:
    (DSADS) (DS-1X) = DS B .  
    Whether or not the system will be equilibrated depends on the scaling of the matrix A, but if equilibration is used, A is overwritten by DS A DS and B× DS B.
  2. 2.If fact=Nag_NotFactored or Nag_EquilibrateAndFactor, the Cholesky decomposition is used to factor the matrix A (after equilibration if fact=Nag_EquilibrateAndFactor) as A=UHU if uplo=Nag_Upper or A=LLH if uplo=Nag_Lower, where U is an upper triangular matrix and L is a lower triangular matrix.
  3. 3.If the leading i×i principal minor of A is not positive definite, then the function returns with fail.errnum=i and fail.code= NE_MAT_NOT_POS_DEF. Otherwise, the factored form of A is used to estimate the condition number of the matrix A. If the reciprocal of the condition number is less than machine precision, fail.code= NE_SINGULAR_WP is returned as a warning, but the function still goes on to solve for X and compute error bounds as described below.
  4. 4.The system of equations is solved for X using the factored form of A.
  5. 5.Iterative refinement is applied to improve the computed solution matrix and to calculate error bounds and backward error estimates for it.
  6. 6.If equilibration was used, the matrix X is premultiplied by DS so that it solves the original system before equilibration.

4 References

Anderson E, Bai Z, Bischof C, Blackford S, Demmel J, Dongarra J J, Du Croz J J, Greenbaum A, Hammarling S, McKenney A and Sorensen D (1999) LAPACK Users' Guide (3rd Edition) SIAM, Philadelphia https://www.netlib.org/lapack/lug
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
Higham N J (2002) Accuracy and Stability of Numerical Algorithms (2nd Edition) SIAM, Philadelphia

5 Arguments

1: order Nag_OrderType Input
On entry: the order argument specifies the two-dimensional storage scheme being used, i.e., row-major ordering or column-major ordering. C language defined storage is specified by order=Nag_RowMajor. See Section 3.1.3 in the Introduction to the NAG Library CL Interface for a more detailed explanation of the use of this argument.
Constraint: order=Nag_RowMajor or Nag_ColMajor.
2: fact Nag_FactoredFormType Input
On entry: specifies whether or not the factorized form of the matrix A is supplied on entry, and if not, whether the matrix A should be equilibrated before it is factorized.
fact=Nag_Factored
afb contains the factorized form of A. If equed=Nag_Equilibrated, the matrix A has been equilibrated with scaling factors given by s. ab and afb will not be modified.
fact=Nag_NotFactored
The matrix A will be copied to afb and factorized.
fact=Nag_EquilibrateAndFactor
The matrix A will be equilibrated if necessary, then copied to afb and factorized.
Constraint: fact=Nag_Factored, Nag_NotFactored or Nag_EquilibrateAndFactor.
3: uplo Nag_UploType Input
On entry: if uplo=Nag_Upper, the upper triangle of A is stored.
If uplo=Nag_Lower, the lower triangle of A is stored.
Constraint: uplo=Nag_Upper or Nag_Lower.
4: n Integer Input
On entry: n, the number of linear equations, i.e., the order of the matrix A.
Constraint: n0.
5: kd Integer Input
On entry: kd, the number of superdiagonals of the matrix A if uplo=Nag_Upper, or the number of subdiagonals if uplo=Nag_Lower.
Constraint: kd0.
6: nrhs Integer Input
On entry: r, the number of right-hand sides, i.e., the number of columns of the matrix B.
Constraint: nrhs0.
7: ab[dim] Complex Input/Output
Note: the dimension, dim, of the array ab must be at least max(1,pdab×n).
On entry: the upper or lower triangle of the Hermitian band matrix A, except if fact=Nag_Factored and equed=Nag_Equilibrated, in which case ab must contain the equilibrated matrix DSADS.
This is stored as a notional two-dimensional array with row elements or column elements stored contiguously. The storage of elements of Aij, depends on the order and uplo arguments as follows:
if order=Nag_ColMajor and uplo=Nag_Upper,
Aij is stored in ab[kd+i-j+(j-1)×pdab], for j=1,,n and i=max(1,j-kd),,j;
if order=Nag_ColMajor and uplo=Nag_Lower,
Aij is stored in ab[i-j+(j-1)×pdab], for j=1,,n and i=j,,min(n,j+kd);
if order=Nag_RowMajor and uplo=Nag_Upper,
Aij is stored in ab[j-i+(i-1)×pdab], for i=1,,n and j=i,,min(n,i+kd);
if order=Nag_RowMajor and uplo=Nag_Lower,
Aij is stored in ab[kd+j-i+(i-1)×pdab], for i=1,,n and j=max(1,i-kd),,i.
On exit: if fact=Nag_EquilibrateAndFactor and equed=Nag_Equilibrated, ab is overwritten by DSADS.
8: pdab Integer Input
On entry: the stride separating row or column elements (depending on the value of order) of the matrix A in the array ab.
Constraint: pdabkd+1.
9: afb[dim] Complex Input/Output
Note: the dimension, dim, of the array afb must be at least max(1,pdafb×n).
On entry: if fact=Nag_Factored, afb contains the triangular factor U or L from the Cholesky factorization A=UHU or A=LLH of the band matrix A, in the same storage format as A. If equed=Nag_Equilibrated, afb is the factorized form of the equilibrated matrix A.
On exit: if fact=Nag_NotFactored, afb returns the triangular factor U or L from the Cholesky factorization A=UHU or A=LLH.
If fact=Nag_EquilibrateAndFactor, afb returns the triangular factor U or L from the Cholesky factorization A=UHU or A=LLH of the equilibrated matrix A (see the description of ab for the form of the equilibrated matrix).
10: pdafb Integer Input
On entry: the stride separating row or column elements (depending on the value of order) of the matrix A in the array afb.
Constraint: pdafbkd+1.
11: equed Nag_EquilibrationType * Input/Output
On entry: if fact=Nag_NotFactored or Nag_EquilibrateAndFactor, equed need not be set.
If fact=Nag_Factored, equed must specify the form of the equilibration that was performed as follows:
  • if equed=Nag_NoEquilibration, no equilibration;
  • if equed=Nag_Equilibrated, equilibration was performed, i.e., A has been replaced by DSADS.
On exit: if fact=Nag_Factored, equed is unchanged from entry.
Otherwise, if no constraints are violated, equed specifies the form of the equilibration that was performed as specified above.
Constraint: if fact=Nag_Factored, equed=Nag_NoEquilibration or Nag_Equilibrated.
12: s[dim] double Input/Output
Note: the dimension, dim, of the array s must be at least max(1,n).
On entry: if fact=Nag_NotFactored or Nag_EquilibrateAndFactor, s need not be set.
If fact=Nag_Factored and equed=Nag_Equilibrated, s must contain the scale factors, DS, for A; each element of s must be positive.
On exit: if fact=Nag_Factored, s is unchanged from entry.
Otherwise, if no constraints are violated and equed=Nag_Equilibrated, s contains the scale factors, DS, for A; each element of s is positive.
13: b[dim] Complex Input/Output
Note: the dimension, dim, of the array b must be at least
  • max(1,pdb×nrhs) when order=Nag_ColMajor;
  • max(1,n×pdb) when order=Nag_RowMajor.
The (i,j)th element of the matrix B is stored in
  • b[(j-1)×pdb+i-1] when order=Nag_ColMajor;
  • b[(i-1)×pdb+j-1] when order=Nag_RowMajor.
On entry: the n×r right-hand side matrix B.
On exit: if equed=Nag_NoEquilibration, b is not modified.
If equed=Nag_Equilibrated, b is overwritten by DSB.
14: pdb Integer Input
On entry: the stride separating row or column elements (depending on the value of order) in the array b.
Constraints:
  • if order=Nag_ColMajor, pdbmax(1,n);
  • if order=Nag_RowMajor, pdbmax(1,nrhs).
15: x[dim] Complex Output
Note: the dimension, dim, of the array x must be at least
  • max(1,pdx×nrhs) when order=Nag_ColMajor;
  • max(1,n×pdx) when order=Nag_RowMajor.
The (i,j)th element of the matrix X is stored in
  • x[(j-1)×pdx+i-1] when order=Nag_ColMajor;
  • x[(i-1)×pdx+j-1] when order=Nag_RowMajor.
On exit: if fail.code= NE_NOERROR or NE_SINGULAR_WP, the n×r solution matrix X to the original system of equations. Note that the arrays A and B are modified on exit if equed=Nag_Equilibrated, and the solution to the equilibrated system is DS-1X.
16: pdx Integer Input
On entry: the stride separating row or column elements (depending on the value of order) in the array x.
Constraints:
  • if order=Nag_ColMajor, pdxmax(1,n);
  • if order=Nag_RowMajor, pdxmax(1,nrhs).
17: rcond double * Output
On exit: if no constraints are violated, an estimate of the reciprocal condition number of the matrix A (after equilibration if that is performed), computed as rcond=1.0/(A1A-11).
18: ferr[nrhs] double Output
On exit: if fail.code= NE_NOERROR or NE_SINGULAR_WP, an estimate of the forward error bound for each computed solution vector, such that x^j-xj/xjferr[j-1] where x^j is the jth column of the computed solution returned in the array x and xj is the corresponding column of the exact solution X. The estimate is as reliable as the estimate for rcond, and is almost always a slight overestimate of the true error.
19: berr[nrhs] double Output
On exit: if fail.code= NE_NOERROR or NE_SINGULAR_WP, an estimate of the component-wise relative backward error of each computed solution vector x^j (i.e., the smallest relative change in any element of A or B that makes x^j an exact solution).
20: fail NagError * Input/Output
The NAG error argument (see Section 7 in the Introduction to the NAG Library CL Interface).

6 Error Indicators and Warnings

NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 3.1.2 in the Introduction to the NAG Library CL Interface for further information.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_INT
On entry, kd=value.
Constraint: kd0.
On entry, n=value.
Constraint: n0.
On entry, nrhs=value.
Constraint: nrhs0.
On entry, pdab=value.
Constraint: pdab>0.
On entry, pdafb=value.
Constraint: pdafb>0.
On entry, pdb=value.
Constraint: pdb>0.
On entry, pdx=value.
Constraint: pdx>0.
NE_INT_2
On entry, pdab=value and kd=value.
Constraint: pdabkd+1.
On entry, pdafb=value and kd=value.
Constraint: pdafbkd+1.
On entry, pdb=value and n=value.
Constraint: pdbmax(1,n).
On entry, pdb=value and nrhs=value.
Constraint: pdbmax(1,nrhs).
On entry, pdx=value and n=value.
Constraint: pdxmax(1,n).
On entry, pdx=value and nrhs=value.
Constraint: pdxmax(1,nrhs).
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.
See Section 7.5 in the Introduction to the NAG Library CL Interface for further information.
NE_MAT_NOT_POS_DEF
The leading minor of order value of A is not positive definite, so the factorization could not be completed, and the solution has not been computed. rcond=0.0 is returned.
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library CL Interface for further information.
NE_SINGULAR_WP
U (or L) is nonsingular, but rcond is less than machine precision, meaning that the matrix is singular to working precision. Nevertheless, the solution and error bounds are computed because there are a number of situations where the computed solution can be more accurate than the value of rcond would suggest.

7 Accuracy

For each right-hand side vector b, the computed solution x is the exact solution of a perturbed system of equations (A+E) x = b, where c(n) is a modest linear function of n, and ε is the machine precision. See Section 10.1 of Higham (2002) for further details.
If x^ is the true solution, then the computed solution x satisfies a forward error bound of the form
x-x^ x^ wc cond(A,x^,b)  
where cond(A,x^,b) = |A-1|(|A||x^|+|b|)/ x^ cond(A) = |A-1||A|κ (A). If x^ is the j th column of X , then wc is returned in berr[j-1] and a bound on x-x^ / x^ is returned in ferr[j-1] . See Section 4.4 of Anderson et al. (1999) for further details.

8 Parallelism and Performance

f07hpc is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f07hpc 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

When nk , the factorization of A requires approximately 4 n (k+1) 2 floating-point operations, where k is the number of superdiagonals.
For each right-hand side, computation of the backward error involves a minimum of 32nk floating-point operations. Each step of iterative refinement involves an additional 48nk operations. At most five steps of iterative refinement are performed, but usually only one or two steps are required. Estimating the forward error involves solving a number of systems of equations of the form Ax=b ; the number is usually 4 or 5 and never more than 11. Each solution involves approximately 16nk operations.
The real analogue of this function is f07hbc.

10 Example

This example solves the equations
AX=B ,  
where A is the Hermitian positive definite band matrix
A = ( 9.39i+0.00 1.08-1.73i 0.00i+0.00 0.00i+0.00 1.08+1.73i 1.69i+0.00 -0.04+0.29i 0.00i+0.00 0.00i+0.00 -0.04-0.29i 2.65i+0.00 -0.33+2.24i 0.00i+0.00 0.00i+0.00 -0.33-2.24i 2.17i+0.00 )  
and
B = ( -12.42+68.42i 54.30-56.56i -9.93+00.88i 18.32+04.76i -27.30-00.01i -4.40+09.97i 5.31+23.63i 9.43+01.41i ) .  
Error estimates for the solutions, information on equilibration and an estimate of the reciprocal of the condition number of the scaled matrix A are also output.

10.1 Program Text

Program Text (f07hpce.c)

10.2 Program Data

Program Data (f07hpce.d)

10.3 Program Results

Program Results (f07hpce.r)