f07fp 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 and X and B are n by r matrices. Error bounds on the solution and a condition estimate are also provided.

Syntax

C#
public static void f07fp(
	string fact,
	string uplo,
	int n,
	int nrhs,
	Complex[,] a,
	Complex[,] af,
	ref string equed,
	double[] s,
	Complex[,] b,
	Complex[,] x,
	out double rcond,
	double[] ferr,
	double[] berr,
	out int info
)
Visual Basic
Public Shared Sub f07fp ( _
	fact As String, _
	uplo As String, _
	n As Integer, _
	nrhs As Integer, _
	a As Complex(,), _
	af As Complex(,), _
	ByRef equed As String, _
	s As Double(), _
	b As Complex(,), _
	x As Complex(,), _
	<OutAttribute> ByRef rcond As Double, _
	ferr As Double(), _
	berr As Double(), _
	<OutAttribute> ByRef info As Integer _
)
Visual C++
public:
static void f07fp(
	String^ fact, 
	String^ uplo, 
	int n, 
	int nrhs, 
	array<Complex,2>^ a, 
	array<Complex,2>^ af, 
	String^% equed, 
	array<double>^ s, 
	array<Complex,2>^ b, 
	array<Complex,2>^ x, 
	[OutAttribute] double% rcond, 
	array<double>^ ferr, 
	array<double>^ berr, 
	[OutAttribute] int% info
)
F#
static member f07fp : 
        fact : string * 
        uplo : string * 
        n : int * 
        nrhs : int * 
        a : Complex[,] * 
        af : Complex[,] * 
        equed : string byref * 
        s : float[] * 
        b : Complex[,] * 
        x : Complex[,] * 
        rcond : float byref * 
        ferr : float[] * 
        berr : float[] * 
        info : int byref -> unit 

Parameters

fact
Type: System..::..String
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="F"
af contains the factorized form of A. If equed="Y", the matrix A has been equilibrated with scaling factors given by s. a and af will not be modified.
fact="N"
The matrix A will be copied to af and factorized.
fact="E"
The matrix A will be equilibrated if necessary, then copied to af and factorized.
Constraint: fact="F", "N" or "E".
uplo
Type: System..::..String
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".
n
Type: System..::..Int32
On entry: n, the number of linear equations, i.e., the order of the matrix A.
Constraint: n0.
nrhs
Type: System..::..Int32
On entry: r, the number of right-hand sides, i.e., the number of columns of the matrix B.
Constraint: nrhs0.
a
Type: array<NagLibrary..::..Complex,2>[,](,)[,][,]
An array of size [dim1, dim2]
Note: dim1 must satisfy the constraint: dim1max1,n
Note: the second dimension of the array a must be at least max1,n.
On entry: the n by n Hermitian matrix A.
If fact="F" and equed="Y", a must have been equilibrated by the scaling factor in s as DSADS.
  • If uplo="U", the upper triangular part of A must be stored and the elements of the array below the diagonal are not referenced.
  • If uplo="L", the lower triangular part of A must be stored and the elements of the array above the diagonal are not referenced.
On exit: if fact="F" or "N", or if fact="E" and equed="N", a is not modified.
If fact="E" and equed="Y", a is overwritten by DSADS.
af
Type: array<NagLibrary..::..Complex,2>[,](,)[,][,]
An array of size [dim1, dim2]
Note: dim1 must satisfy the constraint: dim1max1,n
Note: the second dimension of the array af must be at least max1,n.
On entry: if fact="F", af contains the triangular factor U or L from the Cholesky factorization A=UHU or A=LLH, in the same storage format as a. If equed"N", af is the factorized form of the equilibrated matrix DSADS.
On exit: if fact="N", af returns the triangular factor U or L from the Cholesky factorization A=UHU or A=LLH of the original matrix A.
If fact="E", af 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 a for the form of the equilibrated matrix).
equed
Type: System..::..String%
On entry: 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.
On exit: 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.
Constraint: if fact="F", equed="N" or "Y".
s
Type: array<System..::..Double>[]()[][]
An array of size [dim1]
Note: the dimension of the array s must be at least max1,n.
On entry: 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.
On exit: 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.
b
Type: array<NagLibrary..::..Complex,2>[,](,)[,][,]
An array of size [dim1, dim2]
Note: dim1 must satisfy the constraint: dim1max1,n
Note: the second dimension of the array b must be at least max1,nrhs.
On entry: the n by r right-hand side matrix B.
On exit: if equed="N", b is not modified.
If equed="Y", b is overwritten by DSB.
x
Type: array<NagLibrary..::..Complex,2>[,](,)[,][,]
An array of size [dim1, dim2]
Note: dim1 must satisfy the constraint: dim1max1,n
Note: the second dimension of the array x must be at least max1,nrhs.
On exit: 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.
rcond
Type: System..::..Double%
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.
ferr
Type: array<System..::..Double>[]()[][]
An array of size [nrhs]
On exit: if info=0 or n+1, 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.
berr
Type: array<System..::..Double>[]()[][]
An array of size [nrhs]
On exit: 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).
info
Type: System..::..Int32%
On exit: info=0 unless the method detects an error (see [Error Indicators and Warnings]).

Description

f07fp performs the following steps:
1. If fact="E", real diagonal scaling factors, DS, are computed to equilibrate the system:
DSADSDS-1X=DSB.
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 DSADS and B by DSB.
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 method 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, info=n+1 is returned as a warning, but the method 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

Error Indicators and Warnings

Some error messages may refer to parameters that are dropped from this interface (LDA, LDAF, LDB, LDX) In these cases, an error in another parameter has usually caused an incorrect value to be inferred.
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 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.
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.
ifail=-9000
An error occured, see message report.
ifail=-6000
Invalid Parameters value
ifail=-4000
Invalid dimension for array value
ifail=-8000
Negative dimension for array value
ifail=-6000
Invalid Parameters value
ifail=-6000
Invalid Parameters value

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
  • if uplo="U", EcnεUHU;
  • if uplo="L", EcnεLLH,
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^wccondA,x^,b
where condA,x^,b=A-1Ax^+b/x^condA=A-1AκA. If x^ is the jth 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.

Parallelism and Performance

None.

Further Comments

The factorization of A requires approximately 43n3 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 method is f07fb.

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.64-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.

Example program (C#): f07fpe.cs

Example program data: f07fpe.d

Example program results: f07fpe.r

See Also