NAG FL Interface
f07gnf (zppsv)

Settings help

FL Name Style:


FL Specification Language:


1 Purpose

f07gnf computes the solution to a complex system of linear equations
AX=B ,  
where A is an n×n Hermitian positive definite matrix stored in packed format and X and B are n×r matrices.

2 Specification

Fortran Interface
Subroutine f07gnf ( uplo, n, nrhs, ap, b, ldb, info)
Integer, Intent (In) :: n, nrhs, ldb
Integer, Intent (Out) :: info
Complex (Kind=nag_wp), Intent (Inout) :: ap(*), b(ldb,*)
Character (1), Intent (In) :: uplo
C Header Interface
#include <nag.h>
void  f07gnf_ (const char *uplo, const Integer *n, const Integer *nrhs, Complex ap[], Complex b[], const Integer *ldb, Integer *info, const Charlen length_uplo)
The routine may be called by the names f07gnf, nagf_lapacklin_zppsv or its LAPACK name zppsv.

3 Description

f07gnf uses the Cholesky decomposition to factor A 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. The factored form of A is then used to solve the system of equations AX=B.

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

5 Arguments

1: uplo Character(1) Input
On entry: 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'.
2: n Integer Input
On entry: n, the number of linear equations, i.e., the order of the matrix A.
Constraint: n0.
3: nrhs Integer Input
On entry: r, the number of right-hand sides, i.e., the number of columns of the matrix B.
Constraint: nrhs0.
4: ap(*) Complex (Kind=nag_wp) array Input/Output
Note: the dimension of the array ap must be at least max(1,n×(n+1)/2).
On entry: the n×n Hermitian matrix A, packed by columns.
More precisely,
  • if uplo='U', the upper triangle of A must be stored with element Aij in ap(i+j(j-1)/2) for ij;
  • if uplo='L', the lower triangle of A must be stored with element Aij in ap(i+(2n-j)(j-1)/2) for ij.
On exit: if info=0, the factor U or L from the Cholesky factorization A=UHU or A=LLH, in the same storage format as A.
5: b(ldb,*) Complex (Kind=nag_wp) array Input/Output
Note: the second dimension of the array b must be at least max(1,nrhs).
to solve the equations Ax=b, where b is a single right-hand side, b may be supplied as a one-dimensional array with length ldb=max(1,n).
On entry: the n×r right-hand side matrix B.
On exit: if info=0, the n×r solution matrix X.
6: ldb Integer Input
On entry: the first dimension of the array b as declared in the (sub)program from which f07gnf is called.
Constraint: ldbmax(1,n).
7: info Integer Output
On exit: info=0 unless the routine detects an error (see Section 6).

6 Error Indicators and 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>0
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.

7 Accuracy

The computed solution for a single right-hand side, x^ , satisfies an equation of the form
(A+E) x^=b ,  
where
E1 = O(ε) A1  
and ε is the machine precision. An approximate error bound for the computed solution is given by
x^-x1 x1 κ(A) E1 A1 ,  
where κ(A) = A-11 A1 , the condition number of A with respect to the solution of the linear equations. See Section 4.4 of Anderson et al. (1999) for further details.
f07gpf is a comprehensive LAPACK driver that returns forward and backward error bounds and an estimate of the condition number. Alternatively, f04cef solves Ax=b and returns a forward error bound and condition estimate. f04cef calls f07gnf to solve the equations.

8 Parallelism and Performance

f07gnf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f07gnf 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 routine. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9 Further Comments

The total number of floating-point operations is approximately 43 n3 + 8n2r , where r is the number of right-hand sides.
The real analogue of this routine is f07gaf.

10 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 6.17+09.42i -7.17-21.83i 1.99-14.38i ) .  
Details of the Cholesky factorization of A are also output.

10.1 Program Text

Program Text (f07gnfe.f90)

10.2 Program Data

Program Data (f07gnfe.d)

10.3 Program Results

Program Results (f07gnfe.r)