LAPACK
3.4.2
LAPACK: Linear Algebra PACKage
|
Functions/Subroutines | |
subroutine | cposv (UPLO, N, NRHS, A, LDA, B, LDB, INFO) |
CPOSV computes the solution to system of linear equations A * X = B for PO matrices | |
subroutine | cposvx (FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, RWORK, INFO) |
CPOSVX computes the solution to system of linear equations A * X = B for PO matrices | |
subroutine | cposvxx (FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, S, B, LDB, X, LDX, RCOND, RPVGRW, BERR, N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, RWORK, INFO) |
CPOSVXX computes the solution to system of linear equations A * X = B for PO matrices |
This is the group of complex solve driver functions for PO matrices
subroutine cposv | ( | character | UPLO, |
integer | N, | ||
integer | NRHS, | ||
complex, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
integer | INFO | ||
) |
CPOSV computes the solution to system of linear equations A * X = B for PO matrices
Download CPOSV + dependencies [TGZ] [ZIP] [TXT]CPOSV computes the solution to a complex system of linear equations A * X = B, where A is an N-by-N Hermitian positive definite matrix and X and B are N-by-NRHS matrices. The Cholesky decomposition is used to factor A as A = U**H* U, if UPLO = 'U', or A = L * L**H, 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 A * X = B.
[in] | UPLO | UPLO is CHARACTER*1 = 'U': Upper triangle of A is stored; = 'L': Lower triangle of A is stored. |
[in] | N | N is INTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0. |
[in] | NRHS | NRHS is INTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0. |
[in,out] | A | A is COMPLEX array, dimension (LDA,N) On entry, the Hermitian matrix A. If UPLO = 'U', the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = 'L', the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. On exit, if INFO = 0, the factor U or L from the Cholesky factorization A = U**H*U or A = L*L**H. |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | B | B is COMPLEX array, dimension (LDB,NRHS) On entry, the N-by-NRHS right hand side matrix B. On exit, if INFO = 0, the N-by-NRHS solution matrix X. |
[in] | LDB | LDB is INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value > 0: if INFO = i, 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. |
Definition at line 131 of file cposv.f.
subroutine cposvx | ( | character | FACT, |
character | UPLO, | ||
integer | N, | ||
integer | NRHS, | ||
complex, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex, dimension( ldaf, * ) | AF, | ||
integer | LDAF, | ||
character | EQUED, | ||
real, dimension( * ) | S, | ||
complex, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
complex, dimension( ldx, * ) | X, | ||
integer | LDX, | ||
real | RCOND, | ||
real, dimension( * ) | FERR, | ||
real, dimension( * ) | BERR, | ||
complex, dimension( * ) | WORK, | ||
real, dimension( * ) | RWORK, | ||
integer | INFO | ||
) |
CPOSVX computes the solution to system of linear equations A * X = B for PO matrices
Download CPOSVX + dependencies [TGZ] [ZIP] [TXT]CPOSVX uses the Cholesky factorization A = U**H*U or A = L*L**H to compute the solution to a complex system of linear equations A * X = B, where A is an N-by-N Hermitian positive definite matrix and X and B are N-by-NRHS matrices. Error bounds on the solution and a condition estimate are also provided.
The following steps are performed: 1. 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. 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to factor the matrix A (after equilibration if FACT = 'E') as A = U**H* U, if UPLO = 'U', or A = L * L**H, 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 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. 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 calculate error bounds and backward error estimates for it. 6. If equilibration was used, the matrix X is premultiplied by diag(S) so that it solves the original system before equilibration.
[in] | FACT | FACT is CHARACTER*1 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': On entry, AF contains the factored 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. = 'N': The matrix A will be copied to AF and factored. = 'E': The matrix A will be equilibrated if necessary, then copied to AF and factored. |
[in] | UPLO | UPLO is CHARACTER*1 = 'U': Upper triangle of A is stored; = 'L': Lower triangle of A is stored. |
[in] | N | N is INTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0. |
[in] | NRHS | NRHS is INTEGER The number of right hand sides, i.e., the number of columns of the matrices B and X. NRHS >= 0. |
[in,out] | A | A is COMPLEX array, dimension (LDA,N) On entry, the Hermitian matrix A, except if FACT = 'F' and EQUED = 'Y', then A must contain the equilibrated matrix diag(S)*A*diag(S). If UPLO = 'U', the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = 'L', the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. A is not modified if FACT = 'F' or 'N', or if FACT = 'E' and EQUED = 'N' on exit. On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by diag(S)*A*diag(S). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | AF | AF is COMPLEX array, dimension (LDAF,N) If FACT = 'F', then AF is an input argument and on entry contains the triangular factor U or L from the Cholesky factorization A = U**H*U or A = L*L**H, in the same storage format as A. If EQUED .ne. 'N', then AF is the factored form of the equilibrated matrix diag(S)*A*diag(S). If FACT = 'N', then AF is an output argument and on exit returns the triangular factor U or L from the Cholesky factorization A = U**H*U or A = L*L**H of the original matrix A. If FACT = 'E', then AF is an output argument and on exit returns the triangular factor U or L from the Cholesky factorization A = U**H*U or A = L*L**H of the equilibrated matrix A (see the description of A for the form of the equilibrated matrix). |
[in] | LDAF | LDAF is INTEGER The leading dimension of the array AF. LDAF >= max(1,N). |
[in,out] | EQUED | EQUED is CHARACTER*1 Specifies the form of equilibration that was done. = 'N': No equilibration (always true if FACT = 'N'). = 'Y': Equilibration was done, i.e., A has been replaced by diag(S) * A * diag(S). EQUED is an input argument if FACT = 'F'; otherwise, it is an output argument. |
[in,out] | S | S is REAL array, dimension (N) The scale factors for A; not accessed if EQUED = 'N'. S is an input argument if FACT = 'F'; otherwise, S is an output argument. If FACT = 'F' and EQUED = 'Y', each element of S must be positive. |
[in,out] | B | B is COMPLEX array, dimension (LDB,NRHS) On entry, the N-by-NRHS righthand side matrix B. On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y', B is overwritten by diag(S) * B. |
[in] | LDB | LDB is INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[out] | X | X is COMPLEX array, dimension (LDX,NRHS) If INFO = 0 or INFO = N+1, the N-by-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. |
[in] | LDX | LDX is INTEGER The leading dimension of the array X. LDX >= max(1,N). |
[out] | RCOND | RCOND is REAL 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 | FERR is REAL array, dimension (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) 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 | BERR is REAL array, dimension (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 | WORK is COMPLEX array, dimension (2*N) |
[out] | RWORK | RWORK is REAL array, dimension (N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value > 0: if INFO = i, and i is <= 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. |
Definition at line 305 of file cposvx.f.
subroutine cposvxx | ( | character | FACT, |
character | UPLO, | ||
integer | N, | ||
integer | NRHS, | ||
complex, dimension( lda, * ) | A, | ||
integer | LDA, | ||
complex, dimension( ldaf, * ) | AF, | ||
integer | LDAF, | ||
character | EQUED, | ||
real, dimension( * ) | S, | ||
complex, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
complex, dimension( ldx, * ) | X, | ||
integer | LDX, | ||
real | RCOND, | ||
real | RPVGRW, | ||
real, dimension( * ) | BERR, | ||
integer | N_ERR_BNDS, | ||
real, dimension( nrhs, * ) | ERR_BNDS_NORM, | ||
real, dimension( nrhs, * ) | ERR_BNDS_COMP, | ||
integer | NPARAMS, | ||
real, dimension( * ) | PARAMS, | ||
complex, dimension( * ) | WORK, | ||
real, dimension( * ) | RWORK, | ||
integer | INFO | ||
) |
CPOSVXX computes the solution to system of linear equations A * X = B for PO matrices
Download CPOSVXX + dependencies [TGZ] [ZIP] [TXT]CPOSVXX uses the Cholesky factorization A = U**T*U or A = L*L**T to compute the solution to a complex system of linear equations A * X = B, where A is an N-by-N symmetric positive definite matrix and X and B are N-by-NRHS matrices. If requested, both normwise and maximum componentwise error bounds are returned. CPOSVXX will return a solution with a tiny guaranteed error (O(eps) where eps is the working machine precision) unless the matrix is very ill-conditioned, in which case a warning is returned. Relevant condition numbers also are calculated and returned. CPOSVXX accepts user-provided factorizations and equilibration factors; see the definitions of the FACT and EQUED options. Solving with refinement and using a factorization from a previous CPOSVXX call will also produce a solution with either O(eps) errors or warnings, but we cannot make that claim for general user-provided factorizations and equilibration factors if they differ from what CPOSVXX would itself produce.
The following steps are performed: 1. 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. 2. 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. 3. If the leading i-by-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 (see argument RCOND). If the reciprocal of the condition number is less than machine precision, the routine 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. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero), the routine will use iterative refinement to try to get a small error and error bounds. Refinement calculates the residual to at least twice the working precision. 6. If equilibration was used, the matrix X is premultiplied by diag(S) so that it solves the original system before equilibration.
Some optional parameters are bundled in the PARAMS array. These settings determine how refinement is performed, but often the defaults are acceptable. If the defaults are acceptable, users can pass NPARAMS = 0 which prevents the source code from accessing the PARAMS argument.
[in] | FACT | FACT is CHARACTER*1 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': On entry, AF contains the factored form of A. If EQUED is not 'N', the matrix A has been equilibrated with scaling factors given by S. A and AF are not modified. = 'N': The matrix A will be copied to AF and factored. = 'E': The matrix A will be equilibrated if necessary, then copied to AF and factored. |
[in] | UPLO | UPLO is CHARACTER*1 = 'U': Upper triangle of A is stored; = 'L': Lower triangle of A is stored. |
[in] | N | N is INTEGER The number of linear equations, i.e., the order of the matrix A. N >= 0. |
[in] | NRHS | NRHS is INTEGER The number of right hand sides, i.e., the number of columns of the matrices B and X. NRHS >= 0. |
[in,out] | A | A is COMPLEX array, dimension (LDA,N) On entry, the symmetric matrix A, except if FACT = 'F' and EQUED = 'Y', then A must contain the equilibrated matrix diag(S)*A*diag(S). If UPLO = 'U', the leading N-by-N upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If UPLO = 'L', the leading N-by-N lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced. A is not modified if FACT = 'F' or 'N', or if FACT = 'E' and EQUED = 'N' on exit. On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by diag(S)*A*diag(S). |
[in] | LDA | LDA is INTEGER The leading dimension of the array A. LDA >= max(1,N). |
[in,out] | AF | AF is COMPLEX array, dimension (LDAF,N) If FACT = 'F', then AF is an input argument and on entry 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 A. If EQUED .ne. 'N', then AF is the factored form of the equilibrated matrix diag(S)*A*diag(S). If FACT = 'N', then AF is an output argument and on exit 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', then AF is an output argument and on exit 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 A for the form of the equilibrated matrix). |
[in] | LDAF | LDAF is INTEGER The leading dimension of the array AF. LDAF >= max(1,N). |
[in,out] | EQUED | EQUED is CHARACTER*1 Specifies the form of equilibration that was done. = 'N': No equilibration (always true if FACT = 'N'). = 'Y': Both row and column equilibration, i.e., A has been replaced by diag(S) * A * diag(S). EQUED is an input argument if FACT = 'F'; otherwise, it is an output argument. |
[in,out] | S | S is REAL array, dimension (N) The row scale factors for A. If EQUED = 'Y', A is multiplied on the left and right by diag(S). S is an input argument if FACT = 'F'; otherwise, S is an output argument. If FACT = 'F' and EQUED = 'Y', each element of S must be positive. If S is output, each element of S is a power of the radix. If S is input, each element of S should be a power of the radix to ensure a reliable solution and error estimates. Scaling by powers of the radix does not cause rounding errors unless the result underflows or overflows. Rounding errors during scaling lead to refining with a matrix that is not equivalent to the input matrix, producing error estimates that may not be reliable. |
[in,out] | B | B is COMPLEX array, dimension (LDB,NRHS) On entry, the N-by-NRHS right hand side matrix B. On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y', B is overwritten by diag(S)*B; |
[in] | LDB | LDB is INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[out] | X | X is COMPLEX array, dimension (LDX,NRHS) If INFO = 0, the N-by-NRHS solution matrix X to the original system of equations. Note that A and B are modified on exit if EQUED .ne. 'N', and the solution to the equilibrated system is inv(diag(S))*X. |
[in] | LDX | LDX is INTEGER The leading dimension of the array X. LDX >= max(1,N). |
[out] | RCOND | RCOND is REAL Reciprocal scaled condition number. This is an estimate of the reciprocal Skeel condition number of the matrix A after equilibration (if done). If this is less than the machine precision (in particular, if it is zero), the matrix is singular to working precision. Note that the error may still be small even if this number is very small and the matrix appears ill- conditioned. |
[out] | RPVGRW | RPVGRW is REAL Reciprocal pivot growth. On exit, this contains the reciprocal pivot growth factor norm(A)/norm(U). The "max absolute element" norm is used. If this is much less than 1, then the stability of the LU factorization of the (equilibrated) matrix A could be poor. This also means that the solution X, estimated condition numbers, and error bounds could be unreliable. If factorization fails with 0<INFO<=N, then this contains the reciprocal pivot growth factor for the leading INFO columns of A. |
[out] | BERR | BERR is REAL array, dimension (NRHS) Componentwise relative backward error. This is 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). |
[in] | N_ERR_BNDS | N_ERR_BNDS is INTEGER Number of error bounds to return for each right hand side and each type (normwise or componentwise). See ERR_BNDS_NORM and ERR_BNDS_COMP below. |
[out] | ERR_BNDS_NORM | ERR_BNDS_NORM is REAL array, dimension (NRHS, N_ERR_BNDS) For each right-hand side, this array contains information about various error bounds and condition numbers corresponding to the normwise relative error, which is defined as follows: Normwise relative error in the ith solution vector: max_j (abs(XTRUE(j,i) - X(j,i))) ------------------------------ max_j abs(X(j,i)) The array is indexed by the type of error information as described below. There currently are up to three pieces of information returned. The first index in ERR_BNDS_NORM(i,:) corresponds to the ith right-hand side. The second index in ERR_BNDS_NORM(:,err) contains the following three fields: err = 1 "Trust/don't trust" boolean. Trust the answer if the reciprocal condition number is less than the threshold sqrt(n) * slamch('Epsilon'). err = 2 "Guaranteed" error bound: The estimated forward error, almost certainly within a factor of 10 of the true error so long as the next entry is greater than the threshold sqrt(n) * slamch('Epsilon'). This error bound should only be trusted if the previous boolean is true. err = 3 Reciprocal condition number: Estimated normwise reciprocal condition number. Compared with the threshold sqrt(n) * slamch('Epsilon') to determine if the error estimate is "guaranteed". These reciprocal condition numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some appropriately scaled matrix Z. Let Z = S*A, where S scales each row by a power of the radix so all absolute row sums of Z are approximately 1. See Lapack Working Note 165 for further details and extra cautions. |
[out] | ERR_BNDS_COMP | ERR_BNDS_COMP is REAL array, dimension (NRHS, N_ERR_BNDS) For each right-hand side, this array contains information about various error bounds and condition numbers corresponding to the componentwise relative error, which is defined as follows: Componentwise relative error in the ith solution vector: abs(XTRUE(j,i) - X(j,i)) max_j ---------------------- abs(X(j,i)) The array is indexed by the right-hand side i (on which the componentwise relative error depends), and the type of error information as described below. There currently are up to three pieces of information returned for each right-hand side. If componentwise accuracy is not requested (PARAMS(3) = 0.0), then ERR_BNDS_COMP is not accessed. If N_ERR_BNDS .LT. 3, then at most the first (:,N_ERR_BNDS) entries are returned. The first index in ERR_BNDS_COMP(i,:) corresponds to the ith right-hand side. The second index in ERR_BNDS_COMP(:,err) contains the following three fields: err = 1 "Trust/don't trust" boolean. Trust the answer if the reciprocal condition number is less than the threshold sqrt(n) * slamch('Epsilon'). err = 2 "Guaranteed" error bound: The estimated forward error, almost certainly within a factor of 10 of the true error so long as the next entry is greater than the threshold sqrt(n) * slamch('Epsilon'). This error bound should only be trusted if the previous boolean is true. err = 3 Reciprocal condition number: Estimated componentwise reciprocal condition number. Compared with the threshold sqrt(n) * slamch('Epsilon') to determine if the error estimate is "guaranteed". These reciprocal condition numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some appropriately scaled matrix Z. Let Z = S*(A*diag(x)), where x is the solution for the current right-hand side and S scales each row of A*diag(x) by a power of the radix so all absolute row sums of Z are approximately 1. See Lapack Working Note 165 for further details and extra cautions. |
[in] | NPARAMS | NPARAMS is INTEGER Specifies the number of parameters set in PARAMS. If .LE. 0, the PARAMS array is never referenced and default values are used. |
[in,out] | PARAMS | PARAMS is REAL array, dimension NPARAMS Specifies algorithm parameters. If an entry is .LT. 0.0, then that entry will be filled with default value used for that parameter. Only positions up to NPARAMS are accessed; defaults are used for higher-numbered parameters. PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative refinement or not. Default: 1.0 = 0.0 : No refinement is performed, and no error bounds are computed. = 1.0 : Use the double-precision refinement algorithm, possibly with doubled-single computations if the compilation environment does not support DOUBLE PRECISION. (other values are reserved for future use) PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual computations allowed for refinement. Default: 10 Aggressive: Set to 100 to permit convergence using approximate factorizations or factorizations other than LU. If the factorization uses a technique other than Gaussian elimination, the guarantees in err_bnds_norm and err_bnds_comp may no longer be trustworthy. PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code will attempt to find a solution with small componentwise relative error in the double-precision algorithm. Positive is true, 0.0 is false. Default: 1.0 (attempt componentwise convergence) |
[out] | WORK | WORK is COMPLEX array, dimension (2*N) |
[out] | RWORK | RWORK is REAL array, dimension (2*N) |
[out] | INFO | INFO is INTEGER = 0: Successful exit. The solution to every right-hand side is guaranteed. < 0: If INFO = -i, the i-th argument had an illegal value > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization has been completed, but the factor U is exactly singular, so the solution and error bounds could not be computed. RCOND = 0 is returned. = N+J: The solution corresponding to the Jth right-hand side is not guaranteed. The solutions corresponding to other right- hand sides K with K > J may not be guaranteed as well, but only the first such right-hand side is reported. If a small componentwise error is not requested (PARAMS(3) = 0.0) then the Jth right-hand side is the first with a normwise error bound that is not guaranteed (the smallest J such that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0) the Jth right-hand side is the first with either a normwise or componentwise error bound that is not guaranteed (the smallest J such that either ERR_BNDS_NORM(J,1) = 0.0 or ERR_BNDS_COMP(J,1) = 0.0). See the definition of ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information about all of the right-hand sides check ERR_BNDS_NORM or ERR_BNDS_COMP. |
Definition at line 494 of file cposvxx.f.