LAPACK 3.3.0
|
00001 SUBROUTINE ZPORFSX( UPLO, EQUED, N, NRHS, A, LDA, AF, LDAF, S, B, 00002 $ LDB, X, LDX, RCOND, BERR, N_ERR_BNDS, 00003 $ ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS, 00004 $ WORK, RWORK, INFO ) 00005 * 00006 * -- LAPACK routine (version 3.2.2) -- 00007 * -- Contributed by James Demmel, Deaglan Halligan, Yozo Hida and -- 00008 * -- Jason Riedy of Univ. of California Berkeley. -- 00009 * -- June 2010 -- 00010 * 00011 * -- LAPACK is a software package provided by Univ. of Tennessee, -- 00012 * -- Univ. of California Berkeley and NAG Ltd. -- 00013 * 00014 IMPLICIT NONE 00015 * .. 00016 * .. Scalar Arguments .. 00017 CHARACTER UPLO, EQUED 00018 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS, 00019 $ N_ERR_BNDS 00020 DOUBLE PRECISION RCOND 00021 * .. 00022 * .. Array Arguments .. 00023 COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ), 00024 $ X( LDX, * ), WORK( * ) 00025 DOUBLE PRECISION RWORK( * ), S( * ), PARAMS(*), BERR( * ), 00026 $ ERR_BNDS_NORM( NRHS, * ), 00027 $ ERR_BNDS_COMP( NRHS, * ) 00028 * .. 00029 * 00030 * Purpose 00031 * ======= 00032 * 00033 * ZPORFSX improves the computed solution to a system of linear 00034 * equations when the coefficient matrix is symmetric positive 00035 * definite, and provides error bounds and backward error estimates 00036 * for the solution. In addition to normwise error bound, the code 00037 * provides maximum componentwise error bound if possible. See 00038 * comments for ERR_BNDS_NORM and ERR_BNDS_COMP for details of the 00039 * error bounds. 00040 * 00041 * The original system of linear equations may have been equilibrated 00042 * before calling this routine, as described by arguments EQUED and S 00043 * below. In this case, the solution and error bounds returned are 00044 * for the original unequilibrated system. 00045 * 00046 * Arguments 00047 * ========= 00048 * 00049 * Some optional parameters are bundled in the PARAMS array. These 00050 * settings determine how refinement is performed, but often the 00051 * defaults are acceptable. If the defaults are acceptable, users 00052 * can pass NPARAMS = 0 which prevents the source code from accessing 00053 * the PARAMS argument. 00054 * 00055 * UPLO (input) CHARACTER*1 00056 * = 'U': Upper triangle of A is stored; 00057 * = 'L': Lower triangle of A is stored. 00058 * 00059 * EQUED (input) CHARACTER*1 00060 * Specifies the form of equilibration that was done to A 00061 * before calling this routine. This is needed to compute 00062 * the solution and error bounds correctly. 00063 * = 'N': No equilibration 00064 * = 'Y': Both row and column equilibration, i.e., A has been 00065 * replaced by diag(S) * A * diag(S). 00066 * The right hand side B has been changed accordingly. 00067 * 00068 * N (input) INTEGER 00069 * The order of the matrix A. N >= 0. 00070 * 00071 * NRHS (input) INTEGER 00072 * The number of right hand sides, i.e., the number of columns 00073 * of the matrices B and X. NRHS >= 0. 00074 * 00075 * A (input) COMPLEX*16 array, dimension (LDA,N) 00076 * The symmetric matrix A. If UPLO = 'U', the leading N-by-N 00077 * upper triangular part of A contains the upper triangular part 00078 * of the matrix A, and the strictly lower triangular part of A 00079 * is not referenced. If UPLO = 'L', the leading N-by-N lower 00080 * triangular part of A contains the lower triangular part of 00081 * the matrix A, and the strictly upper triangular part of A is 00082 * not referenced. 00083 * 00084 * LDA (input) INTEGER 00085 * The leading dimension of the array A. LDA >= max(1,N). 00086 * 00087 * AF (input) COMPLEX*16 array, dimension (LDAF,N) 00088 * The triangular factor U or L from the Cholesky factorization 00089 * A = U**T*U or A = L*L**T, as computed by DPOTRF. 00090 * 00091 * LDAF (input) INTEGER 00092 * The leading dimension of the array AF. LDAF >= max(1,N). 00093 * 00094 * S (input or output) DOUBLE PRECISION array, dimension (N) 00095 * The row scale factors for A. If EQUED = 'Y', A is multiplied on 00096 * the left and right by diag(S). S is an input argument if FACT = 00097 * 'F'; otherwise, S is an output argument. If FACT = 'F' and EQUED 00098 * = 'Y', each element of S must be positive. If S is output, each 00099 * element of S is a power of the radix. If S is input, each element 00100 * of S should be a power of the radix to ensure a reliable solution 00101 * and error estimates. Scaling by powers of the radix does not cause 00102 * rounding errors unless the result underflows or overflows. 00103 * Rounding errors during scaling lead to refining with a matrix that 00104 * is not equivalent to the input matrix, producing error estimates 00105 * that may not be reliable. 00106 * 00107 * B (input) COMPLEX*16 array, dimension (LDB,NRHS) 00108 * The right hand side matrix B. 00109 * 00110 * LDB (input) INTEGER 00111 * The leading dimension of the array B. LDB >= max(1,N). 00112 * 00113 * X (input/output) COMPLEX*16 array, dimension (LDX,NRHS) 00114 * On entry, the solution matrix X, as computed by DGETRS. 00115 * On exit, the improved solution matrix X. 00116 * 00117 * LDX (input) INTEGER 00118 * The leading dimension of the array X. LDX >= max(1,N). 00119 * 00120 * RCOND (output) DOUBLE PRECISION 00121 * Reciprocal scaled condition number. This is an estimate of the 00122 * reciprocal Skeel condition number of the matrix A after 00123 * equilibration (if done). If this is less than the machine 00124 * precision (in particular, if it is zero), the matrix is singular 00125 * to working precision. Note that the error may still be small even 00126 * if this number is very small and the matrix appears ill- 00127 * conditioned. 00128 * 00129 * BERR (output) DOUBLE PRECISION array, dimension (NRHS) 00130 * Componentwise relative backward error. This is the 00131 * componentwise relative backward error of each solution vector X(j) 00132 * (i.e., the smallest relative change in any element of A or B that 00133 * makes X(j) an exact solution). 00134 * 00135 * N_ERR_BNDS (input) INTEGER 00136 * Number of error bounds to return for each right hand side 00137 * and each type (normwise or componentwise). See ERR_BNDS_NORM and 00138 * ERR_BNDS_COMP below. 00139 * 00140 * ERR_BNDS_NORM (output) DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS) 00141 * For each right-hand side, this array contains information about 00142 * various error bounds and condition numbers corresponding to the 00143 * normwise relative error, which is defined as follows: 00144 * 00145 * Normwise relative error in the ith solution vector: 00146 * max_j (abs(XTRUE(j,i) - X(j,i))) 00147 * ------------------------------ 00148 * max_j abs(X(j,i)) 00149 * 00150 * The array is indexed by the type of error information as described 00151 * below. There currently are up to three pieces of information 00152 * returned. 00153 * 00154 * The first index in ERR_BNDS_NORM(i,:) corresponds to the ith 00155 * right-hand side. 00156 * 00157 * The second index in ERR_BNDS_NORM(:,err) contains the following 00158 * three fields: 00159 * err = 1 "Trust/don't trust" boolean. Trust the answer if the 00160 * reciprocal condition number is less than the threshold 00161 * sqrt(n) * dlamch('Epsilon'). 00162 * 00163 * err = 2 "Guaranteed" error bound: The estimated forward error, 00164 * almost certainly within a factor of 10 of the true error 00165 * so long as the next entry is greater than the threshold 00166 * sqrt(n) * dlamch('Epsilon'). This error bound should only 00167 * be trusted if the previous boolean is true. 00168 * 00169 * err = 3 Reciprocal condition number: Estimated normwise 00170 * reciprocal condition number. Compared with the threshold 00171 * sqrt(n) * dlamch('Epsilon') to determine if the error 00172 * estimate is "guaranteed". These reciprocal condition 00173 * numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some 00174 * appropriately scaled matrix Z. 00175 * Let Z = S*A, where S scales each row by a power of the 00176 * radix so all absolute row sums of Z are approximately 1. 00177 * 00178 * See Lapack Working Note 165 for further details and extra 00179 * cautions. 00180 * 00181 * ERR_BNDS_COMP (output) DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS) 00182 * For each right-hand side, this array contains information about 00183 * various error bounds and condition numbers corresponding to the 00184 * componentwise relative error, which is defined as follows: 00185 * 00186 * Componentwise relative error in the ith solution vector: 00187 * abs(XTRUE(j,i) - X(j,i)) 00188 * max_j ---------------------- 00189 * abs(X(j,i)) 00190 * 00191 * The array is indexed by the right-hand side i (on which the 00192 * componentwise relative error depends), and the type of error 00193 * information as described below. There currently are up to three 00194 * pieces of information returned for each right-hand side. If 00195 * componentwise accuracy is not requested (PARAMS(3) = 0.0), then 00196 * ERR_BNDS_COMP is not accessed. If N_ERR_BNDS .LT. 3, then at most 00197 * the first (:,N_ERR_BNDS) entries are returned. 00198 * 00199 * The first index in ERR_BNDS_COMP(i,:) corresponds to the ith 00200 * right-hand side. 00201 * 00202 * The second index in ERR_BNDS_COMP(:,err) contains the following 00203 * three fields: 00204 * err = 1 "Trust/don't trust" boolean. Trust the answer if the 00205 * reciprocal condition number is less than the threshold 00206 * sqrt(n) * dlamch('Epsilon'). 00207 * 00208 * err = 2 "Guaranteed" error bound: The estimated forward error, 00209 * almost certainly within a factor of 10 of the true error 00210 * so long as the next entry is greater than the threshold 00211 * sqrt(n) * dlamch('Epsilon'). This error bound should only 00212 * be trusted if the previous boolean is true. 00213 * 00214 * err = 3 Reciprocal condition number: Estimated componentwise 00215 * reciprocal condition number. Compared with the threshold 00216 * sqrt(n) * dlamch('Epsilon') to determine if the error 00217 * estimate is "guaranteed". These reciprocal condition 00218 * numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some 00219 * appropriately scaled matrix Z. 00220 * Let Z = S*(A*diag(x)), where x is the solution for the 00221 * current right-hand side and S scales each row of 00222 * A*diag(x) by a power of the radix so all absolute row 00223 * sums of Z are approximately 1. 00224 * 00225 * See Lapack Working Note 165 for further details and extra 00226 * cautions. 00227 * 00228 * NPARAMS (input) INTEGER 00229 * Specifies the number of parameters set in PARAMS. If .LE. 0, the 00230 * PARAMS array is never referenced and default values are used. 00231 * 00232 * PARAMS (input / output) DOUBLE PRECISION array, dimension NPARAMS 00233 * Specifies algorithm parameters. If an entry is .LT. 0.0, then 00234 * that entry will be filled with default value used for that 00235 * parameter. Only positions up to NPARAMS are accessed; defaults 00236 * are used for higher-numbered parameters. 00237 * 00238 * PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative 00239 * refinement or not. 00240 * Default: 1.0D+0 00241 * = 0.0 : No refinement is performed, and no error bounds are 00242 * computed. 00243 * = 1.0 : Use the double-precision refinement algorithm, 00244 * possibly with doubled-single computations if the 00245 * compilation environment does not support DOUBLE 00246 * PRECISION. 00247 * (other values are reserved for future use) 00248 * 00249 * PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual 00250 * computations allowed for refinement. 00251 * Default: 10 00252 * Aggressive: Set to 100 to permit convergence using approximate 00253 * factorizations or factorizations other than LU. If 00254 * the factorization uses a technique other than 00255 * Gaussian elimination, the guarantees in 00256 * err_bnds_norm and err_bnds_comp may no longer be 00257 * trustworthy. 00258 * 00259 * PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code 00260 * will attempt to find a solution with small componentwise 00261 * relative error in the double-precision algorithm. Positive 00262 * is true, 0.0 is false. 00263 * Default: 1.0 (attempt componentwise convergence) 00264 * 00265 * WORK (workspace) COMPLEX*16 array, dimension (2*N) 00266 * 00267 * RWORK (workspace) DOUBLE PRECISION array, dimension (2*N) 00268 * 00269 * INFO (output) INTEGER 00270 * = 0: Successful exit. The solution to every right-hand side is 00271 * guaranteed. 00272 * < 0: If INFO = -i, the i-th argument had an illegal value 00273 * > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization 00274 * has been completed, but the factor U is exactly singular, so 00275 * the solution and error bounds could not be computed. RCOND = 0 00276 * is returned. 00277 * = N+J: The solution corresponding to the Jth right-hand side is 00278 * not guaranteed. The solutions corresponding to other right- 00279 * hand sides K with K > J may not be guaranteed as well, but 00280 * only the first such right-hand side is reported. If a small 00281 * componentwise error is not requested (PARAMS(3) = 0.0) then 00282 * the Jth right-hand side is the first with a normwise error 00283 * bound that is not guaranteed (the smallest J such 00284 * that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0) 00285 * the Jth right-hand side is the first with either a normwise or 00286 * componentwise error bound that is not guaranteed (the smallest 00287 * J such that either ERR_BNDS_NORM(J,1) = 0.0 or 00288 * ERR_BNDS_COMP(J,1) = 0.0). See the definition of 00289 * ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information 00290 * about all of the right-hand sides check ERR_BNDS_NORM or 00291 * ERR_BNDS_COMP. 00292 * 00293 * ================================================================== 00294 * 00295 * .. Parameters .. 00296 DOUBLE PRECISION ZERO, ONE 00297 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 ) 00298 DOUBLE PRECISION ITREF_DEFAULT, ITHRESH_DEFAULT 00299 DOUBLE PRECISION COMPONENTWISE_DEFAULT, RTHRESH_DEFAULT 00300 DOUBLE PRECISION DZTHRESH_DEFAULT 00301 PARAMETER ( ITREF_DEFAULT = 1.0D+0 ) 00302 PARAMETER ( ITHRESH_DEFAULT = 10.0D+0 ) 00303 PARAMETER ( COMPONENTWISE_DEFAULT = 1.0D+0 ) 00304 PARAMETER ( RTHRESH_DEFAULT = 0.5D+0 ) 00305 PARAMETER ( DZTHRESH_DEFAULT = 0.25D+0 ) 00306 INTEGER LA_LINRX_ITREF_I, LA_LINRX_ITHRESH_I, 00307 $ LA_LINRX_CWISE_I 00308 PARAMETER ( LA_LINRX_ITREF_I = 1, 00309 $ LA_LINRX_ITHRESH_I = 2 ) 00310 PARAMETER ( LA_LINRX_CWISE_I = 3 ) 00311 INTEGER LA_LINRX_TRUST_I, LA_LINRX_ERR_I, 00312 $ LA_LINRX_RCOND_I 00313 PARAMETER ( LA_LINRX_TRUST_I = 1, LA_LINRX_ERR_I = 2 ) 00314 PARAMETER ( LA_LINRX_RCOND_I = 3 ) 00315 * .. 00316 * .. Local Scalars .. 00317 CHARACTER(1) NORM 00318 LOGICAL RCEQU 00319 INTEGER J, PREC_TYPE, REF_TYPE 00320 INTEGER N_NORMS 00321 DOUBLE PRECISION ANORM, RCOND_TMP 00322 DOUBLE PRECISION ILLRCOND_THRESH, ERR_LBND, CWISE_WRONG 00323 LOGICAL IGNORE_CWISE 00324 INTEGER ITHRESH 00325 DOUBLE PRECISION RTHRESH, UNSTABLE_THRESH 00326 * .. 00327 * .. External Subroutines .. 00328 EXTERNAL XERBLA, ZPOCON, ZLA_PORFSX_EXTENDED 00329 * .. 00330 * .. Intrinsic Functions .. 00331 INTRINSIC MAX, SQRT, TRANSFER 00332 * .. 00333 * .. External Functions .. 00334 EXTERNAL LSAME, BLAS_FPINFO_X, ILATRANS, ILAPREC 00335 EXTERNAL DLAMCH, ZLANHE, ZLA_PORCOND_X, ZLA_PORCOND_C 00336 DOUBLE PRECISION DLAMCH, ZLANHE, ZLA_PORCOND_X, ZLA_PORCOND_C 00337 LOGICAL LSAME 00338 INTEGER BLAS_FPINFO_X 00339 INTEGER ILATRANS, ILAPREC 00340 * .. 00341 * .. Executable Statements .. 00342 * 00343 * Check the input parameters. 00344 * 00345 INFO = 0 00346 REF_TYPE = INT( ITREF_DEFAULT ) 00347 IF ( NPARAMS .GE. LA_LINRX_ITREF_I ) THEN 00348 IF ( PARAMS( LA_LINRX_ITREF_I ) .LT. 0.0D+0 ) THEN 00349 PARAMS( LA_LINRX_ITREF_I ) = ITREF_DEFAULT 00350 ELSE 00351 REF_TYPE = PARAMS( LA_LINRX_ITREF_I ) 00352 END IF 00353 END IF 00354 * 00355 * Set default parameters. 00356 * 00357 ILLRCOND_THRESH = DBLE( N ) * DLAMCH( 'Epsilon' ) 00358 ITHRESH = INT( ITHRESH_DEFAULT ) 00359 RTHRESH = RTHRESH_DEFAULT 00360 UNSTABLE_THRESH = DZTHRESH_DEFAULT 00361 IGNORE_CWISE = COMPONENTWISE_DEFAULT .EQ. 0.0D+0 00362 * 00363 IF ( NPARAMS.GE.LA_LINRX_ITHRESH_I ) THEN 00364 IF ( PARAMS(LA_LINRX_ITHRESH_I ).LT.0.0D+0 ) THEN 00365 PARAMS( LA_LINRX_ITHRESH_I ) = ITHRESH 00366 ELSE 00367 ITHRESH = INT( PARAMS( LA_LINRX_ITHRESH_I ) ) 00368 END IF 00369 END IF 00370 IF ( NPARAMS.GE.LA_LINRX_CWISE_I ) THEN 00371 IF ( PARAMS(LA_LINRX_CWISE_I ).LT.0.0D+0 ) THEN 00372 IF ( IGNORE_CWISE ) THEN 00373 PARAMS( LA_LINRX_CWISE_I ) = 0.0D+0 00374 ELSE 00375 PARAMS( LA_LINRX_CWISE_I ) = 1.0D+0 00376 END IF 00377 ELSE 00378 IGNORE_CWISE = PARAMS( LA_LINRX_CWISE_I ) .EQ. 0.0D+0 00379 END IF 00380 END IF 00381 IF ( REF_TYPE .EQ. 0 .OR. N_ERR_BNDS .EQ. 0 ) THEN 00382 N_NORMS = 0 00383 ELSE IF ( IGNORE_CWISE ) THEN 00384 N_NORMS = 1 00385 ELSE 00386 N_NORMS = 2 00387 END IF 00388 * 00389 RCEQU = LSAME( EQUED, 'Y' ) 00390 * 00391 * Test input parameters. 00392 * 00393 IF (.NOT.LSAME( UPLO, 'U' ) .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN 00394 INFO = -1 00395 ELSE IF( .NOT.RCEQU .AND. .NOT.LSAME( EQUED, 'N' ) ) THEN 00396 INFO = -2 00397 ELSE IF( N.LT.0 ) THEN 00398 INFO = -3 00399 ELSE IF( NRHS.LT.0 ) THEN 00400 INFO = -4 00401 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN 00402 INFO = -6 00403 ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN 00404 INFO = -8 00405 ELSE IF( LDB.LT.MAX( 1, N ) ) THEN 00406 INFO = -11 00407 ELSE IF( LDX.LT.MAX( 1, N ) ) THEN 00408 INFO = -13 00409 END IF 00410 IF( INFO.NE.0 ) THEN 00411 CALL XERBLA( 'ZPORFSX', -INFO ) 00412 RETURN 00413 END IF 00414 * 00415 * Quick return if possible. 00416 * 00417 IF( N.EQ.0 .OR. NRHS.EQ.0 ) THEN 00418 RCOND = 1.0D+0 00419 DO J = 1, NRHS 00420 BERR( J ) = 0.0D+0 00421 IF ( N_ERR_BNDS .GE. 1 ) THEN 00422 ERR_BNDS_NORM( J, LA_LINRX_TRUST_I ) = 1.0D+0 00423 ERR_BNDS_COMP( J, LA_LINRX_TRUST_I ) = 1.0D+0 00424 END IF 00425 IF ( N_ERR_BNDS .GE. 2 ) THEN 00426 ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) = 0.0D+0 00427 ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) = 0.0D+0 00428 END IF 00429 IF ( N_ERR_BNDS .GE. 3 ) THEN 00430 ERR_BNDS_NORM( J, LA_LINRX_RCOND_I ) = 1.0D+0 00431 ERR_BNDS_COMP( J, LA_LINRX_RCOND_I ) = 1.0D+0 00432 END IF 00433 END DO 00434 RETURN 00435 END IF 00436 * 00437 * Default to failure. 00438 * 00439 RCOND = 0.0D+0 00440 DO J = 1, NRHS 00441 BERR( J ) = 1.0D+0 00442 IF ( N_ERR_BNDS .GE. 1 ) THEN 00443 ERR_BNDS_NORM( J, LA_LINRX_TRUST_I ) = 1.0D+0 00444 ERR_BNDS_COMP( J, LA_LINRX_TRUST_I ) = 1.0D+0 00445 END IF 00446 IF ( N_ERR_BNDS .GE. 2 ) THEN 00447 ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) = 1.0D+0 00448 ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) = 1.0D+0 00449 END IF 00450 IF ( N_ERR_BNDS .GE. 3 ) THEN 00451 ERR_BNDS_NORM( J, LA_LINRX_RCOND_I ) = 0.0D+0 00452 ERR_BNDS_COMP( J, LA_LINRX_RCOND_I ) = 0.0D+0 00453 END IF 00454 END DO 00455 * 00456 * Compute the norm of A and the reciprocal of the condition 00457 * number of A. 00458 * 00459 NORM = 'I' 00460 ANORM = ZLANHE( NORM, UPLO, N, A, LDA, RWORK ) 00461 CALL ZPOCON( UPLO, N, AF, LDAF, ANORM, RCOND, WORK, RWORK, 00462 $ INFO ) 00463 * 00464 * Perform refinement on each right-hand side 00465 * 00466 IF ( REF_TYPE .NE. 0 ) THEN 00467 00468 PREC_TYPE = ILAPREC( 'E' ) 00469 00470 CALL ZLA_PORFSX_EXTENDED( PREC_TYPE, UPLO, N, 00471 $ NRHS, A, LDA, AF, LDAF, RCEQU, S, B, 00472 $ LDB, X, LDX, BERR, N_NORMS, ERR_BNDS_NORM, ERR_BNDS_COMP, 00473 $ WORK, RWORK, WORK(N+1), 00474 $ TRANSFER (RWORK(1:2*N), (/ (ZERO, ZERO) /), N), RCOND, 00475 $ ITHRESH, RTHRESH, UNSTABLE_THRESH, IGNORE_CWISE, 00476 $ INFO ) 00477 END IF 00478 00479 ERR_LBND = MAX( 10.0D+0, SQRT( DBLE( N ) ) ) * DLAMCH( 'Epsilon' ) 00480 IF ( N_ERR_BNDS .GE. 1 .AND. N_NORMS .GE. 1 ) THEN 00481 * 00482 * Compute scaled normwise condition number cond(A*C). 00483 * 00484 IF ( RCEQU ) THEN 00485 RCOND_TMP = ZLA_PORCOND_C( UPLO, N, A, LDA, AF, LDAF, 00486 $ S, .TRUE., INFO, WORK, RWORK ) 00487 ELSE 00488 RCOND_TMP = ZLA_PORCOND_C( UPLO, N, A, LDA, AF, LDAF, 00489 $ S, .FALSE., INFO, WORK, RWORK ) 00490 END IF 00491 DO J = 1, NRHS 00492 * 00493 * Cap the error at 1.0. 00494 * 00495 IF ( N_ERR_BNDS .GE. LA_LINRX_ERR_I 00496 $ .AND. ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) .GT. 1.0D+0 ) 00497 $ ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) = 1.0D+0 00498 * 00499 * Threshold the error (see LAWN). 00500 * 00501 IF ( RCOND_TMP .LT. ILLRCOND_THRESH ) THEN 00502 ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) = 1.0D+0 00503 ERR_BNDS_NORM( J, LA_LINRX_TRUST_I ) = 0.0D+0 00504 IF ( INFO .LE. N ) INFO = N + J 00505 ELSE IF ( ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) .LT. ERR_LBND ) 00506 $ THEN 00507 ERR_BNDS_NORM( J, LA_LINRX_ERR_I ) = ERR_LBND 00508 ERR_BNDS_NORM( J, LA_LINRX_TRUST_I ) = 1.0D+0 00509 END IF 00510 * 00511 * Save the condition number. 00512 * 00513 IF ( N_ERR_BNDS .GE. LA_LINRX_RCOND_I ) THEN 00514 ERR_BNDS_NORM( J, LA_LINRX_RCOND_I ) = RCOND_TMP 00515 END IF 00516 00517 END DO 00518 END IF 00519 00520 IF (N_ERR_BNDS .GE. 1 .AND. N_NORMS .GE. 2) THEN 00521 * 00522 * Compute componentwise condition number cond(A*diag(Y(:,J))) for 00523 * each right-hand side using the current solution as an estimate of 00524 * the true solution. If the componentwise error estimate is too 00525 * large, then the solution is a lousy estimate of truth and the 00526 * estimated RCOND may be too optimistic. To avoid misleading users, 00527 * the inverse condition number is set to 0.0 when the estimated 00528 * cwise error is at least CWISE_WRONG. 00529 * 00530 CWISE_WRONG = SQRT( DLAMCH( 'Epsilon' ) ) 00531 DO J = 1, NRHS 00532 IF (ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) .LT. CWISE_WRONG ) 00533 $ THEN 00534 RCOND_TMP = ZLA_PORCOND_X( UPLO, N, A, LDA, AF, LDAF, 00535 $ X(1,J), INFO, WORK, RWORK ) 00536 ELSE 00537 RCOND_TMP = 0.0D+0 00538 END IF 00539 * 00540 * Cap the error at 1.0. 00541 * 00542 IF ( N_ERR_BNDS .GE. LA_LINRX_ERR_I 00543 $ .AND. ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) .GT. 1.0D+0 ) 00544 $ ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) = 1.0D+0 00545 * 00546 * Threshold the error (see LAWN). 00547 * 00548 IF (RCOND_TMP .LT. ILLRCOND_THRESH) THEN 00549 ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) = 1.0D+0 00550 ERR_BNDS_COMP( J, LA_LINRX_TRUST_I ) = 0.0D+0 00551 IF ( PARAMS( LA_LINRX_CWISE_I ) .EQ. 1.0D+0 00552 $ .AND. INFO.LT.N + J ) INFO = N + J 00553 ELSE IF ( ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) 00554 $ .LT. ERR_LBND ) THEN 00555 ERR_BNDS_COMP( J, LA_LINRX_ERR_I ) = ERR_LBND 00556 ERR_BNDS_COMP( J, LA_LINRX_TRUST_I ) = 1.0D+0 00557 END IF 00558 * 00559 * Save the condition number. 00560 * 00561 IF ( N_ERR_BNDS .GE. LA_LINRX_RCOND_I ) THEN 00562 ERR_BNDS_COMP( J, LA_LINRX_RCOND_I ) = RCOND_TMP 00563 END IF 00564 00565 END DO 00566 END IF 00567 * 00568 RETURN 00569 * 00570 * End of ZPORFSX 00571 * 00572 END