NAG FL Interface
f07fbf (dposvx)

Settings help

FL Name Style:


FL Specification Language:


1 Purpose

f07fbf uses the Cholesky factorization
A=UTU   or   A=LLT  
to compute the solution to a real system of linear equations
AX=B ,  
where A is an n×n symmetric positive definite matrix and X and B are n×r matrices. Error bounds on the solution and a condition estimate are also provided.

2 Specification

Fortran Interface
Subroutine f07fbf ( fact, uplo, n, nrhs, a, lda, af, ldaf, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, iwork, info)
Integer, Intent (In) :: n, nrhs, lda, ldaf, ldb, ldx
Integer, Intent (Out) :: iwork(n), info
Real (Kind=nag_wp), Intent (Inout) :: a(lda,*), af(ldaf,*), s(*), b(ldb,*), x(ldx,*)
Real (Kind=nag_wp), Intent (Out) :: rcond, ferr(nrhs), berr(nrhs), work(3*n)
Character (1), Intent (In) :: fact, uplo
Character (1), Intent (InOut) :: equed
C Header Interface
#include <nag.h>
void  f07fbf_ (const char *fact, const char *uplo, const Integer *n, const Integer *nrhs, double a[], const Integer *lda, double af[], const Integer *ldaf, char *equed, double s[], double b[], const Integer *ldb, double x[], const Integer *ldx, double *rcond, double ferr[], double berr[], double work[], Integer iwork[], Integer *info, const Charlen length_fact, const Charlen length_uplo, const Charlen length_equed)
The routine may be called by the names f07fbf, nagf_lapacklin_dposvx or its LAPACK name dposvx.

3 Description

f07fbf performs the following steps:
  1. 1.If fact='E', 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='N' or 'E', the Cholesky decomposition is used to factor the matrix A (after equilibration if fact='E') as A=UTU if uplo='U' or A=LLT if uplo='L', 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 routine 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 routine 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: fact Character(1) 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='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'.
2: 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'.
3: n Integer Input
On entry: n, the number of linear equations, i.e., the order of the matrix A.
Constraint: n0.
4: nrhs Integer Input
On entry: r, the number of right-hand sides, i.e., the number of columns of the matrix B.
Constraint: nrhs0.
5: a(lda,*) Real (Kind=nag_wp) array Input/Output
Note: the second dimension of the array a must be at least max(1,n).
On entry: the n×n symmetric 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.
6: lda Integer Input
On entry: the first dimension of the array a as declared in the (sub)program from which f07fbf is called.
Constraint: ldamax(1,n).
7: af(ldaf,*) Real (Kind=nag_wp) array Input/Output
Note: the second dimension of the array af must be at least max(1,n).
On entry: if fact='F', af contains the triangular factor U or L from the Cholesky factorization A=UTU or A=LLT, 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=UTU or A=LLT of the original matrix A.
If fact='E', af returns the triangular factor U or L from the Cholesky factorization A=UTU or A=LLT of the equilibrated matrix A (see the description of a for the form of the equilibrated matrix).
8: ldaf Integer Input
On entry: the first dimension of the array af as declared in the (sub)program from which f07fbf is called.
Constraint: ldafmax(1,n).
9: equed Character(1) Input/Output
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'.
10: s(*) Real (Kind=nag_wp) array Input/Output
Note: the dimension of the array s must be at least max(1,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.
11: b(ldb,*) Real (Kind=nag_wp) array Input/Output
Note: the second dimension of the array b must be at least max(1,nrhs).
On entry: the n×r right-hand side matrix B.
On exit: if equed='N', b is not modified.
If equed='Y', b is overwritten by DSB.
12: ldb Integer Input
On entry: the first dimension of the array b as declared in the (sub)program from which f07fbf is called.
Constraint: ldbmax(1,n).
13: x(ldx,*) Real (Kind=nag_wp) array Output
Note: the second dimension of the array x must be at least max(1,nrhs).
On exit: if info=0 or n+1, 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='Y', and the solution to the equilibrated system is DS-1X.
14: ldx Integer Input
On entry: the first dimension of the array x as declared in the (sub)program from which f07fbf is called.
Constraint: ldxmax(1,n).
15: rcond Real (Kind=nag_wp) 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).
16: ferr(nrhs) Real (Kind=nag_wp) array Output
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) 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.
17: berr(nrhs) Real (Kind=nag_wp) array Output
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).
18: work(3×n) Real (Kind=nag_wp) array Workspace
19: iwork(n) Integer array Workspace
20: 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>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.

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) and a bound on x-x^ / x^ is returned in ferr(j) . See Section 4.4 of Anderson et al. (1999) for further details.

8 Parallelism and Performance

Background information to multithreading can be found in the Multithreading documentation.
f07fbf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f07fbf 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 factorization of A requires approximately 13 n3 floating-point operations.
For each right-hand side, computation of the backward error involves a minimum of 4n2 floating-point operations. Each step of iterative refinement involves an additional 6n2 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 2n2 operations.
The complex analogue of this routine is f07fpf.

10 Example

This example solves the equations
AX=B ,  
where A is the symmetric positive definite matrix
A = ( 4.16 -3.12 0.56 -0.10 -3.12 5.03 -0.83 1.18 0.56 -0.83 0.76 0.34 -0.10 1.18 0.34 1.18 )  
and
B = ( 8.70 8.30 -13.35 2.13 1.89 1.61 -4.14 5.00 ) .  
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 (f07fbfe.f90)

10.2 Program Data

Program Data (f07fbfe.d)

10.3 Program Results

Program Results (f07fbfe.r)