|
|
◆ dppsvx()
| void dppsvx |
( |
char |
fact, |
|
|
char |
uplo, |
|
|
int |
n, |
|
|
int |
nrhs, |
|
|
double |
ap[], |
|
|
double |
afp[], |
|
|
char * |
equed, |
|
|
double |
s[], |
|
|
int |
ldb, |
|
|
double |
b[], |
|
|
int |
ldx, |
|
|
double |
x[], |
|
|
double * |
rcond, |
|
|
double |
ferr[], |
|
|
double |
berr[], |
|
|
double |
work[], |
|
|
int |
iwork[], |
|
|
int * |
info |
|
) |
| |
(Expert driver) Solution to system of linear equations AX = B for a symmetric positive definite matrix in packed form
- Purpose
- This routine uses the Cholesky factorization A = U^T*U or A = L*L^T to compute the solution to a real system of
linear equations
A * X = B,
where A is an n x n symmetric positive definite matrix stored in packed form and X and B are n x nrhs matrices.
Error bounds on the solution and a condition estimate are also provided.
- Description
- The following steps are performed:
- If fact = 'E', real scaling factors are computed to equilibrate the system:
diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * 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 diag(S)*A*diag(S) and B by diag(S)*B.
- If fact = 'N' or 'E', the Cholesky decomposition is used to factor the matrix A (after equilibration if fact = 'E') as
A = U^T * U, if uplo = 'U', or
A = L * L^T, if uplo = 'L',
where U is an upper triangular matrix and L is a lower triangular matrix.
- If the leading i x i principal minor 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.
- The system of equations is solved for X using the factored form of A.
- Iterative refinement is applied to improve the computed solution matrix and calculate error bounds and backward error estimates for it.
- If equilibration was used, the matrix X is premultiplied by diag(S) so that it solves the original system before equilibration.
- Parameters
-
| [in] | fact | Specifies whether or not the factored form of the matrix A is supplied on entry, and if not, whether the matrix A should be equilibrated before it is factored.
= 'F': afp[] contains the factored form of ap[]. If equed = 'Y', the matrix A has been equilibrated with scaling factors given by s[]. ap[] and afp[] will not be modified.
= 'N': The matrix A will be copied to afp[] and factored.
= 'E': The matrix A will be equilibrated if necessary, then copied to afp[] and factored. |
| [in] | uplo | = 'U': Upper triangle of A is stored.
= 'L': Lower triangle of A is stored. |
| [in] | n | Number of linear equations, i.e., order of the matrix A. (n >= 0) (If n = 0, returns without computation) |
| [in] | nrhs | Number of right hand sides, i.e., number of columns of the matrix B. (nrhs >= 0) (if nrhs = 0, returns without computation) |
| [in] | ap[] | Array ap[lap] (lap >= n(n + 1)/2)
[in] n x n symmetric positove definite matrix A in packed form, except if fact = 'F' and equed = 'Y', then ap[] must contain the equilibrated matrix diag(S)*A*diag(S). The upper or lower triangular part is to be stored in accordance with uplo.
[out] If fact = 'E' and equed = 'Y', ap[] is overwritten by diag(S)*A*diag(S). Otherwise, ap[] is not modified. |
| [in,out] | afp[] | Array afp[lafp] (lafp >= n(n + 1)/2)
[in] If fact = 'F', afp[] contains the triangular factor U or L from the Cholesky factorization A = U^T*U or A = L*L^T in the same storage format as ap[]. If equed = 'Y', then afp[] is the factored form of the equilibrated matrix diag(S)*A*diag(S).
[out] If fact = 'N', afp[] returns the triangular factor U or L from the Cholesky factorization A = U^T*U or A = L*L^T of the original matrix A.
If fact='E', afp[] returns the triangular factor U or L from the Cholesky factorization A = U^T*U or A = L*L^T of the equilibrated matrix A (see the description of ap[] for the form of the equilibrated matrix). |
| [in,out] | equed | Specifies the form of equilibration that was done.
= 'N': No equilibration.
= 'Y': Equilibration was done, i.e., A has been replaced by diag(S)*A*diag(S).
[in] If fact = 'F', the form of equilibration of given matrix A.
[out] If fact = 'E', returns the result of equilibration. If fact = 'N', always returns 'N'. |
| [in,out] | s[] | Array s[ls] (ls >= n)
The scale factors for A; not accessed if equed = 'N'.
[in] If fact = 'F', the scale factors for given matrix A. (Each element > 0)
[out] If fact = 'E', the resulted scale factors. |
| [in] | ldb | Leading dimension of the two dimensional array b[]. (ldb >= max(1, n)) |
| [in,out] | b[][] | Array b[lb][ldb] (lb >= nrhs)
[in] n x nrhs right hand side matrix B.
[out] If equed = 'N', b[][] is not modified; If equed = 'Y', b[][] is overwritten by diag(S)*B. |
| [in] | ldx | Leading dimension of the two dimensional array x[]. (ldx >= max(1, n)) |
| [out] | x[][] | Array x[lx][ldx] (lx >= nrhs)
If info = 0 or info = n+1, the n x nrhs solution matrix X to the original system of equations. Note that if equed = 'Y', a[][] and b[][] are modified on exit, and the solution to the equilibrated system is inv(diag(S))*X. |
| [out] | rcond | The estimate of the reciprocal condition number of the matrix A after equilibration (if done). If rcond is less than the machine precision (in particular, if rcond=0), the matrix is singular to working precision. This condition is indicated by a return code of info > 0. |
| [out] | ferr[] | Array ferr[lferr] (lferr >= nrhs)
The estimated forward error bound for each solution vector X(j) (the j-th column of the solution matrix X). If Xtrue is the true solution corresponding to X(j), ferr[j-1] is an estimated upper bound for the magnitude of the largest element in (X(j) - Xtrue) divided by the magnitude of the largest element in X(j). The estimate is as reliable as the estimate for rcond, and is almost always a slight overestimate of the true error. |
| [out] | berr[] | Array berr[lberr] (lberr >= nrhs)
The componentwise relative backward error of each solution vector X(j) (i.e., the smallest relative change in any element of A or B that makes X(j) an exact solution). |
| [out] | work[] | Array work[lwork] (lwork >= 3*n)
Work array. |
| [out] | iwork[] | Array iwork[liwork] (liwork >= n)
Work array. |
| [out] | info | = 0: Successful exit
= -1: The argument fact had an illegal value (fact != 'F', 'N' nor 'E')
= -2: The argument uplo had an illegal value (uplo != 'U' nor 'L')
= -3: The argument n had an illegal value (n < 0)
= -4: The argument nrhs had an illegal value (nrhs < 0)
= -7: The argument equed had an illegal value (fact = 'F' and equed != 'N' nor 'Y')
= -8: The argument s had an illegal value (s[i] <= 0 when fact = 'F' and equed = 'Y')
= -9: The argument ldb had an illegal value (ldb < max(1, n))
= -11: The argument ldx had an illegal value (ldx < max(1, n))
= i (0 < i <= n): The leading minor of order i of A is not positive definite, so the factorization could not be completed, and the solution has not been computed. rcond = 0 is returned.
= n+1: U 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. |
- Reference
- LAPACK
|