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