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_lapack_zppsvx (f07gp)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_lapack_zppsvx (f07gp) 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 by n Hermitian positive definite matrix stored in packed format and X and B are n by r matrices. Error bounds on the solution and a condition estimate are also provided.

Syntax

[ap, afp, equed, s, b, x, rcond, ferr, berr, info] = f07gp(fact, uplo, ap, afp, equed, s, b, 'n', n, 'nrhs_p', nrhs_p)
[ap, afp, equed, s, b, x, rcond, ferr, berr, info] = nag_lapack_zppsvx(fact, uplo, ap, afp, equed, s, b, 'n', n, 'nrhs_p', nrhs_p)

Description

nag_lapack_zppsvx (f07gp) performs the following steps:
1. If fact='E', real diagonal scaling factors, DS , are computed to equilibrate the system:
DS A DS DS-1 X = 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 by DS B.
2. If fact='N' or 'E', the Cholesky decomposition is used to factor the matrix A (after equilibration if fact='E') as A=UHU if uplo='U' or A=LLH if uplo='L', where U is an upper triangular matrix and L is a lower triangular matrix.
3. If the leading i by i principal minor of A is not positive definite, then the function returns with info=i. 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, infon+1 is returned as a warning, but the function still goes on to solve for X and compute error bounds as described below.
4. The system of equations is solved for X using the factored form of A.
5. Iterative refinement is applied to improve the computed solution matrix and to calculate error bounds and backward error estimates for it.
6. If equilibration was used, the matrix X is premultiplied by DS  so that it solves the original system before equilibration.

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 http://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

Parameters

Compulsory Input Parameters

1:     fact – string (length ≥ 1)
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='F'
afp contains the factorized form of A. If equed='Y', the matrix A has been equilibrated with scaling factors given by s. ap and afp will not be modified.
fact='N'
The matrix A will be copied to afp and factorized.
fact='E'
The matrix A will be equilibrated if necessary, then copied to afp and factorized.
Constraint: fact='F', 'N' or 'E'.
2:     uplo – string (length ≥ 1)
If uplo='U', the upper triangle of A is stored.
If uplo='L', the lower triangle of A is stored.
Constraint: uplo='U' or 'L'.
3:     ap: – complex array
The dimension of the array ap must be at least max1,n×n+1/2
If fact='F' and equed='Y', ap must contain the equilibrated matrix DSADS; otherwise, ap must contain the n by n Hermitian matrix A, packed by columns.
More precisely,
  • if uplo='U', the upper triangle of A must be stored with element Aij in api+jj-1/2 for ij;
  • if uplo='L', the lower triangle of A must be stored with element Aij in api+2n-jj-1/2 for ij.
4:     afp: – complex array
The dimension of the array afp must be at least max1,n×n+1/2
If fact='F', afp contains the triangular factor U or L from the Cholesky factorization A=UHU or A=LLH, in the same storage format as ap. If equed='Y', afp is the factorized form of the equilibrated matrix DSADS.
5:     equed – string (length ≥ 1)
If fact='N' or 'E', equed need not be set.
If fact='F', equed must specify the form of the equilibration that was performed as follows:
  • if equed='N', no equilibration;
  • if equed='Y', equilibration was performed, i.e., A has been replaced by DSADS.
Constraint: if fact='F', equed='N' or 'Y'.
6:     s: – double array
The dimension of the array s must be at least max1,n
If fact='N' or 'E', s need not be set.
If fact='F' and equed='Y', s must contain the scale factors, DS, for A; each element of s must be positive.
7:     bldb: – complex array
The first dimension of the array b must be at least max1,n.
The second dimension of the array b must be at least max1,nrhs_p.
The n by r right-hand side matrix B.

Optional Input Parameters

1:     n int64int32nag_int scalar
Default: the first dimension of the array b and the dimension of the array s.
n, the number of linear equations, i.e., the order of the matrix A.
Constraint: n0.
2:     nrhs_p int64int32nag_int scalar
Default: the second dimension of the array b.
r, the number of right-hand sides, i.e., the number of columns of the matrix B.
Constraint: nrhs_p0.

Output Parameters

1:     ap: – complex array
The dimension of the array ap will be max1,n×n+1/2
If fact='F' or 'N', or if fact='E' and equed='N', ap is not modified.
If fact='E' and equed='Y', ap stores DSADS.
2:     afp: – complex array
The dimension of the array afp will be max1,n×n+1/2
If fact='N' or if fact='E' and equed='N', afp returns the triangular factor U or L from the Cholesky factorization A=UHU or A=LLH of the original matrix A.
If fact='E' and equed='Y', afp 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 ap for the form of the equilibrated matrix).
3:     equed – string (length ≥ 1)
If fact='F', equed is unchanged from entry.
Otherwise, if no constraints are violated, equed specifies the form of the equilibration that was performed as specified above.
4:     s: – double array
The dimension of the array s will be max1,n
If fact='F', s is unchanged from entry.
Otherwise, if no constraints are violated and equed='Y', s contains the scale factors, DS, for A; each element of s is positive.
5:     bldb: – complex array
The first dimension of the array b will be max1,n.
The second dimension of the array b will be max1,nrhs_p.
If equed='N', b is not modified.
If equed='Y', b stores DSB.
6:     xldx: – complex array
The first dimension of the array x will be max1,n.
The second dimension of the array x will be max1,nrhs_p.
If info=0 or n+1, the n by r solution matrix X to the original system of equations. Note that the arrays A and B are modified on exit if equed='Y', and the solution to the equilibrated system is DS-1X.
7:     rcond – double scalar
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/A1 A-11 .
8:     ferrnrhs_p – double array
If info=0 or n+1, an estimate of the forward error bound for each computed solution vector, such that x^j-xj/xjferrj 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.
9:     berrnrhs_p – double array
If info=0 or n+1, 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).
10:   info int64int32nag_int scalar
info=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

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

   info<0
If info=-i, argument i had an illegal value. An explanatory message is output, and execution of the program is terminated.
   info>0andinfon
The leading minor of order _ 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.
W  info=n+1
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.

Accuracy

For each right-hand side vector b, the computed solution x is the exact solution of a perturbed system of equations A+Ex=b, where cn 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 condA,x^,b ,  
where condA,x^,b = A-1 A x^ + b / x^ condA = A-1 A κ A . If x^  is the j th column of X , then wc  is returned in berrj  and a bound on x - x^ / x^  is returned in ferrj . See Section 4.4 of Anderson et al. (1999) for further details.

Further Comments

The factorization of A  requires approximately 43 n3  floating-point operations.
For each right-hand side, computation of the backward error involves a minimum of 16n2  floating-point operations. Each step of iterative refinement involves an additional 24n2  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 8n2  operations.
The real analogue of this function is nag_lapack_dppsvx (f07gb).

Example

This example solves the equations
AX=B ,  
where A  is the Hermitian positive definite matrix
A = 3.23i+0.00 1.51-1.92i 1.90+0.84i 0.42+2.50i 1.51+1.92i 3.58i+0.00 -0.23+1.11i -1.18+1.37i 1.90-0.84i -0.23-1.11i 4.09i+0.00 2.33-0.14i 0.42-2.50i -1.18-1.37i 2.33+0.14i 4.29i+0.00  
and
B = 3.93-06.14i 1.48+06.58i 6.17+09.42i 4.65-04.75i -7.17-21.83i -4.91+02.29i 1.99-14.38i 7.64-10.79i .  
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.
function f07gp_example


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

% Upper triangular part of Hermitian matrix A
uplo = 'Upper';
n = int64(4);
ap = [ 3.23 + 0i,  ...
       1.51 - 1.92i,  3.58 + 0i,    ...
       1.90 + 0.84i, -0.23 + 1.11i, 4.09 + 0i,     ...
       0.42 + 2.50i, -1.18 + 1.37i, 2.33 - 0.14i,  4.29 + 0i];

% RHS
b = [ 3.93 -  6.14i,  1.48 +  6.58i;
      6.17 +  9.42i,  4.65 -  4.75i;
     -7.17 - 21.83i, -4.91 +  2.29i;
      1.99 - 14.38i,  7.64 - 10.79i];

% Input parameters
fact  = 'Equilibration';
afp   = complex(zeros(10,1));
equed = ' ';
s     = zeros(n,1);
[ap, afp, equed, s, b, x, rcond, ferr, berr, info] = ...
  f07gp( ...
         fact, uplo, ap, afp, equed, s, b);

disp('Solution(s)');
disp(x);
disp('Backward errors (machine-dependent)');
fprintf('%10.1e',berr);
fprintf('\n');
disp('Estimated forward error bounds (machine-dependent)');
fprintf('%10.1e',ferr);
fprintf('\n\n');
disp('Estimate of reciprocal condition number');
fprintf('%10.1e\n\n',rcond);
if equed=='N'
  fprintf('A has not been equilibrated\n');
else
  fprintf('A has been equilibrated\n');
end


f07gp example results

Solution(s)
   1.0000 - 1.0000i  -1.0000 + 2.0000i
  -0.0000 + 3.0000i   3.0000 - 4.0000i
  -4.0000 - 5.0000i  -2.0000 + 3.0000i
   2.0000 + 1.0000i   4.0000 - 5.0000i

Backward errors (machine-dependent)
   8.1e-17   9.9e-17
Estimated forward error bounds (machine-dependent)
   6.2e-14   7.5e-14

Estimate of reciprocal condition number
   6.6e-03

A has not been equilibrated

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