![]() |
LAPACK
3.4.1
LAPACK: Linear Algebra PACKage
|
00001 *> \brief <b> DGESVX computes the solution to system of linear equations A * X = B for GE matrices</b> 00002 * 00003 * =========== DOCUMENTATION =========== 00004 * 00005 * Online html documentation available at 00006 * http://www.netlib.org/lapack/explore-html/ 00007 * 00008 *> \htmlonly 00009 *> Download DGESVX + dependencies 00010 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgesvx.f"> 00011 *> [TGZ]</a> 00012 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgesvx.f"> 00013 *> [ZIP]</a> 00014 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgesvx.f"> 00015 *> [TXT]</a> 00016 *> \endhtmlonly 00017 * 00018 * Definition: 00019 * =========== 00020 * 00021 * SUBROUTINE DGESVX( FACT, TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV, 00022 * EQUED, R, C, B, LDB, X, LDX, RCOND, FERR, BERR, 00023 * WORK, IWORK, INFO ) 00024 * 00025 * .. Scalar Arguments .. 00026 * CHARACTER EQUED, FACT, TRANS 00027 * INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS 00028 * DOUBLE PRECISION RCOND 00029 * .. 00030 * .. Array Arguments .. 00031 * INTEGER IPIV( * ), IWORK( * ) 00032 * DOUBLE PRECISION A( LDA, * ), AF( LDAF, * ), B( LDB, * ), 00033 * $ BERR( * ), C( * ), FERR( * ), R( * ), 00034 * $ WORK( * ), X( LDX, * ) 00035 * .. 00036 * 00037 * 00038 *> \par Purpose: 00039 * ============= 00040 *> 00041 *> \verbatim 00042 *> 00043 *> DGESVX uses the LU factorization to compute the solution to a real 00044 *> system of linear equations 00045 *> A * X = B, 00046 *> where A is an N-by-N matrix and X and B are N-by-NRHS matrices. 00047 *> 00048 *> Error bounds on the solution and a condition estimate are also 00049 *> provided. 00050 *> \endverbatim 00051 * 00052 *> \par Description: 00053 * ================= 00054 *> 00055 *> \verbatim 00056 *> 00057 *> The following steps are performed: 00058 *> 00059 *> 1. If FACT = 'E', real scaling factors are computed to equilibrate 00060 *> the system: 00061 *> TRANS = 'N': diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B 00062 *> TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B 00063 *> TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B 00064 *> Whether or not the system will be equilibrated depends on the 00065 *> scaling of the matrix A, but if equilibration is used, A is 00066 *> overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N') 00067 *> or diag(C)*B (if TRANS = 'T' or 'C'). 00068 *> 00069 *> 2. If FACT = 'N' or 'E', the LU decomposition is used to factor the 00070 *> matrix A (after equilibration if FACT = 'E') as 00071 *> A = P * L * U, 00072 *> where P is a permutation matrix, L is a unit lower triangular 00073 *> matrix, and U is upper triangular. 00074 *> 00075 *> 3. If some U(i,i)=0, so that U is exactly singular, then the routine 00076 *> returns with INFO = i. Otherwise, the factored form of A is used 00077 *> to estimate the condition number of the matrix A. If the 00078 *> reciprocal of the condition number is less than machine precision, 00079 *> INFO = N+1 is returned as a warning, but the routine still goes on 00080 *> to solve for X and compute error bounds as described below. 00081 *> 00082 *> 4. The system of equations is solved for X using the factored form 00083 *> of A. 00084 *> 00085 *> 5. Iterative refinement is applied to improve the computed solution 00086 *> matrix and calculate error bounds and backward error estimates 00087 *> for it. 00088 *> 00089 *> 6. If equilibration was used, the matrix X is premultiplied by 00090 *> diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so 00091 *> that it solves the original system before equilibration. 00092 *> \endverbatim 00093 * 00094 * Arguments: 00095 * ========== 00096 * 00097 *> \param[in] FACT 00098 *> \verbatim 00099 *> FACT is CHARACTER*1 00100 *> Specifies whether or not the factored form of the matrix A is 00101 *> supplied on entry, and if not, whether the matrix A should be 00102 *> equilibrated before it is factored. 00103 *> = 'F': On entry, AF and IPIV contain the factored form of A. 00104 *> If EQUED is not 'N', the matrix A has been 00105 *> equilibrated with scaling factors given by R and C. 00106 *> A, AF, and IPIV are not modified. 00107 *> = 'N': The matrix A will be copied to AF and factored. 00108 *> = 'E': The matrix A will be equilibrated if necessary, then 00109 *> copied to AF and factored. 00110 *> \endverbatim 00111 *> 00112 *> \param[in] TRANS 00113 *> \verbatim 00114 *> TRANS is CHARACTER*1 00115 *> Specifies the form of the system of equations: 00116 *> = 'N': A * X = B (No transpose) 00117 *> = 'T': A**T * X = B (Transpose) 00118 *> = 'C': A**H * X = B (Transpose) 00119 *> \endverbatim 00120 *> 00121 *> \param[in] N 00122 *> \verbatim 00123 *> N is INTEGER 00124 *> The number of linear equations, i.e., the order of the 00125 *> matrix A. N >= 0. 00126 *> \endverbatim 00127 *> 00128 *> \param[in] NRHS 00129 *> \verbatim 00130 *> NRHS is INTEGER 00131 *> The number of right hand sides, i.e., the number of columns 00132 *> of the matrices B and X. NRHS >= 0. 00133 *> \endverbatim 00134 *> 00135 *> \param[in,out] A 00136 *> \verbatim 00137 *> A is DOUBLE PRECISION array, dimension (LDA,N) 00138 *> On entry, the N-by-N matrix A. If FACT = 'F' and EQUED is 00139 *> not 'N', then A must have been equilibrated by the scaling 00140 *> factors in R and/or C. A is not modified if FACT = 'F' or 00141 *> 'N', or if FACT = 'E' and EQUED = 'N' on exit. 00142 *> 00143 *> On exit, if EQUED .ne. 'N', A is scaled as follows: 00144 *> EQUED = 'R': A := diag(R) * A 00145 *> EQUED = 'C': A := A * diag(C) 00146 *> EQUED = 'B': A := diag(R) * A * diag(C). 00147 *> \endverbatim 00148 *> 00149 *> \param[in] LDA 00150 *> \verbatim 00151 *> LDA is INTEGER 00152 *> The leading dimension of the array A. LDA >= max(1,N). 00153 *> \endverbatim 00154 *> 00155 *> \param[in,out] AF 00156 *> \verbatim 00157 *> AF is DOUBLE PRECISION array, dimension (LDAF,N) 00158 *> If FACT = 'F', then AF is an input argument and on entry 00159 *> contains the factors L and U from the factorization 00160 *> A = P*L*U as computed by DGETRF. If EQUED .ne. 'N', then 00161 *> AF is the factored form of the equilibrated matrix A. 00162 *> 00163 *> If FACT = 'N', then AF is an output argument and on exit 00164 *> returns the factors L and U from the factorization A = P*L*U 00165 *> of the original matrix A. 00166 *> 00167 *> If FACT = 'E', then AF is an output argument and on exit 00168 *> returns the factors L and U from the factorization A = P*L*U 00169 *> of the equilibrated matrix A (see the description of A for 00170 *> the form of the equilibrated matrix). 00171 *> \endverbatim 00172 *> 00173 *> \param[in] LDAF 00174 *> \verbatim 00175 *> LDAF is INTEGER 00176 *> The leading dimension of the array AF. LDAF >= max(1,N). 00177 *> \endverbatim 00178 *> 00179 *> \param[in,out] IPIV 00180 *> \verbatim 00181 *> IPIV is INTEGER array, dimension (N) 00182 *> If FACT = 'F', then IPIV is an input argument and on entry 00183 *> contains the pivot indices from the factorization A = P*L*U 00184 *> as computed by DGETRF; row i of the matrix was interchanged 00185 *> with row IPIV(i). 00186 *> 00187 *> If FACT = 'N', then IPIV is an output argument and on exit 00188 *> contains the pivot indices from the factorization A = P*L*U 00189 *> of the original matrix A. 00190 *> 00191 *> If FACT = 'E', then IPIV is an output argument and on exit 00192 *> contains the pivot indices from the factorization A = P*L*U 00193 *> of the equilibrated matrix A. 00194 *> \endverbatim 00195 *> 00196 *> \param[in,out] EQUED 00197 *> \verbatim 00198 *> EQUED is CHARACTER*1 00199 *> Specifies the form of equilibration that was done. 00200 *> = 'N': No equilibration (always true if FACT = 'N'). 00201 *> = 'R': Row equilibration, i.e., A has been premultiplied by 00202 *> diag(R). 00203 *> = 'C': Column equilibration, i.e., A has been postmultiplied 00204 *> by diag(C). 00205 *> = 'B': Both row and column equilibration, i.e., A has been 00206 *> replaced by diag(R) * A * diag(C). 00207 *> EQUED is an input argument if FACT = 'F'; otherwise, it is an 00208 *> output argument. 00209 *> \endverbatim 00210 *> 00211 *> \param[in,out] R 00212 *> \verbatim 00213 *> R is DOUBLE PRECISION array, dimension (N) 00214 *> The row scale factors for A. If EQUED = 'R' or 'B', A is 00215 *> multiplied on the left by diag(R); if EQUED = 'N' or 'C', R 00216 *> is not accessed. R is an input argument if FACT = 'F'; 00217 *> otherwise, R is an output argument. If FACT = 'F' and 00218 *> EQUED = 'R' or 'B', each element of R must be positive. 00219 *> \endverbatim 00220 *> 00221 *> \param[in,out] C 00222 *> \verbatim 00223 *> C is DOUBLE PRECISION array, dimension (N) 00224 *> The column scale factors for A. If EQUED = 'C' or 'B', A is 00225 *> multiplied on the right by diag(C); if EQUED = 'N' or 'R', C 00226 *> is not accessed. C is an input argument if FACT = 'F'; 00227 *> otherwise, C is an output argument. If FACT = 'F' and 00228 *> EQUED = 'C' or 'B', each element of C must be positive. 00229 *> \endverbatim 00230 *> 00231 *> \param[in,out] B 00232 *> \verbatim 00233 *> B is DOUBLE PRECISION array, dimension (LDB,NRHS) 00234 *> On entry, the N-by-NRHS right hand side matrix B. 00235 *> On exit, 00236 *> if EQUED = 'N', B is not modified; 00237 *> if TRANS = 'N' and EQUED = 'R' or 'B', B is overwritten by 00238 *> diag(R)*B; 00239 *> if TRANS = 'T' or 'C' and EQUED = 'C' or 'B', B is 00240 *> overwritten by diag(C)*B. 00241 *> \endverbatim 00242 *> 00243 *> \param[in] LDB 00244 *> \verbatim 00245 *> LDB is INTEGER 00246 *> The leading dimension of the array B. LDB >= max(1,N). 00247 *> \endverbatim 00248 *> 00249 *> \param[out] X 00250 *> \verbatim 00251 *> X is DOUBLE PRECISION array, dimension (LDX,NRHS) 00252 *> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X 00253 *> to the original system of equations. Note that A and B are 00254 *> modified on exit if EQUED .ne. 'N', and the solution to the 00255 *> equilibrated system is inv(diag(C))*X if TRANS = 'N' and 00256 *> EQUED = 'C' or 'B', or inv(diag(R))*X if TRANS = 'T' or 'C' 00257 *> and EQUED = 'R' or 'B'. 00258 *> \endverbatim 00259 *> 00260 *> \param[in] LDX 00261 *> \verbatim 00262 *> LDX is INTEGER 00263 *> The leading dimension of the array X. LDX >= max(1,N). 00264 *> \endverbatim 00265 *> 00266 *> \param[out] RCOND 00267 *> \verbatim 00268 *> RCOND is DOUBLE PRECISION 00269 *> The estimate of the reciprocal condition number of the matrix 00270 *> A after equilibration (if done). If RCOND is less than the 00271 *> machine precision (in particular, if RCOND = 0), the matrix 00272 *> is singular to working precision. This condition is 00273 *> indicated by a return code of INFO > 0. 00274 *> \endverbatim 00275 *> 00276 *> \param[out] FERR 00277 *> \verbatim 00278 *> FERR is DOUBLE PRECISION array, dimension (NRHS) 00279 *> The estimated forward error bound for each solution vector 00280 *> X(j) (the j-th column of the solution matrix X). 00281 *> If XTRUE is the true solution corresponding to X(j), FERR(j) 00282 *> is an estimated upper bound for the magnitude of the largest 00283 *> element in (X(j) - XTRUE) divided by the magnitude of the 00284 *> largest element in X(j). The estimate is as reliable as 00285 *> the estimate for RCOND, and is almost always a slight 00286 *> overestimate of the true error. 00287 *> \endverbatim 00288 *> 00289 *> \param[out] BERR 00290 *> \verbatim 00291 *> BERR is DOUBLE PRECISION array, dimension (NRHS) 00292 *> The componentwise relative backward error of each solution 00293 *> vector X(j) (i.e., the smallest relative change in 00294 *> any element of A or B that makes X(j) an exact solution). 00295 *> \endverbatim 00296 *> 00297 *> \param[out] WORK 00298 *> \verbatim 00299 *> WORK is DOUBLE PRECISION array, dimension (4*N) 00300 *> On exit, WORK(1) contains the reciprocal pivot growth 00301 *> factor norm(A)/norm(U). The "max absolute element" norm is 00302 *> used. If WORK(1) is much less than 1, then the stability 00303 *> of the LU factorization of the (equilibrated) matrix A 00304 *> could be poor. This also means that the solution X, condition 00305 *> estimator RCOND, and forward error bound FERR could be 00306 *> unreliable. If factorization fails with 0<INFO<=N, then 00307 *> WORK(1) contains the reciprocal pivot growth factor for the 00308 *> leading INFO columns of A. 00309 *> \endverbatim 00310 *> 00311 *> \param[out] IWORK 00312 *> \verbatim 00313 *> IWORK is INTEGER array, dimension (N) 00314 *> \endverbatim 00315 *> 00316 *> \param[out] INFO 00317 *> \verbatim 00318 *> INFO is INTEGER 00319 *> = 0: successful exit 00320 *> < 0: if INFO = -i, the i-th argument had an illegal value 00321 *> > 0: if INFO = i, and i is 00322 *> <= N: U(i,i) is exactly zero. The factorization has 00323 *> been completed, but the factor U is exactly 00324 *> singular, so the solution and error bounds 00325 *> could not be computed. RCOND = 0 is returned. 00326 *> = N+1: U is nonsingular, but RCOND is less than machine 00327 *> precision, meaning that the matrix is singular 00328 *> to working precision. Nevertheless, the 00329 *> solution and error bounds are computed because 00330 *> there are a number of situations where the 00331 *> computed solution can be more accurate than the 00332 *> value of RCOND would suggest. 00333 *> \endverbatim 00334 * 00335 * Authors: 00336 * ======== 00337 * 00338 *> \author Univ. of Tennessee 00339 *> \author Univ. of California Berkeley 00340 *> \author Univ. of Colorado Denver 00341 *> \author NAG Ltd. 00342 * 00343 *> \date April 2012 00344 * 00345 *> \ingroup doubleGEsolve 00346 * 00347 * ===================================================================== 00348 SUBROUTINE DGESVX( FACT, TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV, 00349 $ EQUED, R, C, B, LDB, X, LDX, RCOND, FERR, BERR, 00350 $ WORK, IWORK, INFO ) 00351 * 00352 * -- LAPACK driver routine (version 3.4.1) -- 00353 * -- LAPACK is a software package provided by Univ. of Tennessee, -- 00354 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 00355 * April 2012 00356 * 00357 * .. Scalar Arguments .. 00358 CHARACTER EQUED, FACT, TRANS 00359 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS 00360 DOUBLE PRECISION RCOND 00361 * .. 00362 * .. Array Arguments .. 00363 INTEGER IPIV( * ), IWORK( * ) 00364 DOUBLE PRECISION A( LDA, * ), AF( LDAF, * ), B( LDB, * ), 00365 $ BERR( * ), C( * ), FERR( * ), R( * ), 00366 $ WORK( * ), X( LDX, * ) 00367 * .. 00368 * 00369 * ===================================================================== 00370 * 00371 * .. Parameters .. 00372 DOUBLE PRECISION ZERO, ONE 00373 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 ) 00374 * .. 00375 * .. Local Scalars .. 00376 LOGICAL COLEQU, EQUIL, NOFACT, NOTRAN, ROWEQU 00377 CHARACTER NORM 00378 INTEGER I, INFEQU, J 00379 DOUBLE PRECISION AMAX, ANORM, BIGNUM, COLCND, RCMAX, RCMIN, 00380 $ ROWCND, RPVGRW, SMLNUM 00381 * .. 00382 * .. External Functions .. 00383 LOGICAL LSAME 00384 DOUBLE PRECISION DLAMCH, DLANGE, DLANTR 00385 EXTERNAL LSAME, DLAMCH, DLANGE, DLANTR 00386 * .. 00387 * .. External Subroutines .. 00388 EXTERNAL DGECON, DGEEQU, DGERFS, DGETRF, DGETRS, DLACPY, 00389 $ DLAQGE, XERBLA 00390 * .. 00391 * .. Intrinsic Functions .. 00392 INTRINSIC MAX, MIN 00393 * .. 00394 * .. Executable Statements .. 00395 * 00396 INFO = 0 00397 NOFACT = LSAME( FACT, 'N' ) 00398 EQUIL = LSAME( FACT, 'E' ) 00399 NOTRAN = LSAME( TRANS, 'N' ) 00400 IF( NOFACT .OR. EQUIL ) THEN 00401 EQUED = 'N' 00402 ROWEQU = .FALSE. 00403 COLEQU = .FALSE. 00404 ELSE 00405 ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' ) 00406 COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' ) 00407 SMLNUM = DLAMCH( 'Safe minimum' ) 00408 BIGNUM = ONE / SMLNUM 00409 END IF 00410 * 00411 * Test the input parameters. 00412 * 00413 IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.LSAME( FACT, 'F' ) ) 00414 $ THEN 00415 INFO = -1 00416 ELSE IF( .NOT.NOTRAN .AND. .NOT.LSAME( TRANS, 'T' ) .AND. .NOT. 00417 $ LSAME( TRANS, 'C' ) ) THEN 00418 INFO = -2 00419 ELSE IF( N.LT.0 ) THEN 00420 INFO = -3 00421 ELSE IF( NRHS.LT.0 ) THEN 00422 INFO = -4 00423 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN 00424 INFO = -6 00425 ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN 00426 INFO = -8 00427 ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT. 00428 $ ( ROWEQU .OR. COLEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN 00429 INFO = -10 00430 ELSE 00431 IF( ROWEQU ) THEN 00432 RCMIN = BIGNUM 00433 RCMAX = ZERO 00434 DO 10 J = 1, N 00435 RCMIN = MIN( RCMIN, R( J ) ) 00436 RCMAX = MAX( RCMAX, R( J ) ) 00437 10 CONTINUE 00438 IF( RCMIN.LE.ZERO ) THEN 00439 INFO = -11 00440 ELSE IF( N.GT.0 ) THEN 00441 ROWCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM ) 00442 ELSE 00443 ROWCND = ONE 00444 END IF 00445 END IF 00446 IF( COLEQU .AND. INFO.EQ.0 ) THEN 00447 RCMIN = BIGNUM 00448 RCMAX = ZERO 00449 DO 20 J = 1, N 00450 RCMIN = MIN( RCMIN, C( J ) ) 00451 RCMAX = MAX( RCMAX, C( J ) ) 00452 20 CONTINUE 00453 IF( RCMIN.LE.ZERO ) THEN 00454 INFO = -12 00455 ELSE IF( N.GT.0 ) THEN 00456 COLCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM ) 00457 ELSE 00458 COLCND = ONE 00459 END IF 00460 END IF 00461 IF( INFO.EQ.0 ) THEN 00462 IF( LDB.LT.MAX( 1, N ) ) THEN 00463 INFO = -14 00464 ELSE IF( LDX.LT.MAX( 1, N ) ) THEN 00465 INFO = -16 00466 END IF 00467 END IF 00468 END IF 00469 * 00470 IF( INFO.NE.0 ) THEN 00471 CALL XERBLA( 'DGESVX', -INFO ) 00472 RETURN 00473 END IF 00474 * 00475 IF( EQUIL ) THEN 00476 * 00477 * Compute row and column scalings to equilibrate the matrix A. 00478 * 00479 CALL DGEEQU( N, N, A, LDA, R, C, ROWCND, COLCND, AMAX, INFEQU ) 00480 IF( INFEQU.EQ.0 ) THEN 00481 * 00482 * Equilibrate the matrix. 00483 * 00484 CALL DLAQGE( N, N, A, LDA, R, C, ROWCND, COLCND, AMAX, 00485 $ EQUED ) 00486 ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' ) 00487 COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' ) 00488 END IF 00489 END IF 00490 * 00491 * Scale the right hand side. 00492 * 00493 IF( NOTRAN ) THEN 00494 IF( ROWEQU ) THEN 00495 DO 40 J = 1, NRHS 00496 DO 30 I = 1, N 00497 B( I, J ) = R( I )*B( I, J ) 00498 30 CONTINUE 00499 40 CONTINUE 00500 END IF 00501 ELSE IF( COLEQU ) THEN 00502 DO 60 J = 1, NRHS 00503 DO 50 I = 1, N 00504 B( I, J ) = C( I )*B( I, J ) 00505 50 CONTINUE 00506 60 CONTINUE 00507 END IF 00508 * 00509 IF( NOFACT .OR. EQUIL ) THEN 00510 * 00511 * Compute the LU factorization of A. 00512 * 00513 CALL DLACPY( 'Full', N, N, A, LDA, AF, LDAF ) 00514 CALL DGETRF( N, N, AF, LDAF, IPIV, INFO ) 00515 * 00516 * Return if INFO is non-zero. 00517 * 00518 IF( INFO.GT.0 ) THEN 00519 * 00520 * Compute the reciprocal pivot growth factor of the 00521 * leading rank-deficient INFO columns of A. 00522 * 00523 RPVGRW = DLANTR( 'M', 'U', 'N', INFO, INFO, AF, LDAF, 00524 $ WORK ) 00525 IF( RPVGRW.EQ.ZERO ) THEN 00526 RPVGRW = ONE 00527 ELSE 00528 RPVGRW = DLANGE( 'M', N, INFO, A, LDA, WORK ) / RPVGRW 00529 END IF 00530 WORK( 1 ) = RPVGRW 00531 RCOND = ZERO 00532 RETURN 00533 END IF 00534 END IF 00535 * 00536 * Compute the norm of the matrix A and the 00537 * reciprocal pivot growth factor RPVGRW. 00538 * 00539 IF( NOTRAN ) THEN 00540 NORM = '1' 00541 ELSE 00542 NORM = 'I' 00543 END IF 00544 ANORM = DLANGE( NORM, N, N, A, LDA, WORK ) 00545 RPVGRW = DLANTR( 'M', 'U', 'N', N, N, AF, LDAF, WORK ) 00546 IF( RPVGRW.EQ.ZERO ) THEN 00547 RPVGRW = ONE 00548 ELSE 00549 RPVGRW = DLANGE( 'M', N, N, A, LDA, WORK ) / RPVGRW 00550 END IF 00551 * 00552 * Compute the reciprocal of the condition number of A. 00553 * 00554 CALL DGECON( NORM, N, AF, LDAF, ANORM, RCOND, WORK, IWORK, INFO ) 00555 * 00556 * Compute the solution matrix X. 00557 * 00558 CALL DLACPY( 'Full', N, NRHS, B, LDB, X, LDX ) 00559 CALL DGETRS( TRANS, N, NRHS, AF, LDAF, IPIV, X, LDX, INFO ) 00560 * 00561 * Use iterative refinement to improve the computed solution and 00562 * compute error bounds and backward error estimates for it. 00563 * 00564 CALL DGERFS( TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV, B, LDB, X, 00565 $ LDX, FERR, BERR, WORK, IWORK, INFO ) 00566 * 00567 * Transform the solution matrix X to a solution of the original 00568 * system. 00569 * 00570 IF( NOTRAN ) THEN 00571 IF( COLEQU ) THEN 00572 DO 80 J = 1, NRHS 00573 DO 70 I = 1, N 00574 X( I, J ) = C( I )*X( I, J ) 00575 70 CONTINUE 00576 80 CONTINUE 00577 DO 90 J = 1, NRHS 00578 FERR( J ) = FERR( J ) / COLCND 00579 90 CONTINUE 00580 END IF 00581 ELSE IF( ROWEQU ) THEN 00582 DO 110 J = 1, NRHS 00583 DO 100 I = 1, N 00584 X( I, J ) = R( I )*X( I, J ) 00585 100 CONTINUE 00586 110 CONTINUE 00587 DO 120 J = 1, NRHS 00588 FERR( J ) = FERR( J ) / ROWCND 00589 120 CONTINUE 00590 END IF 00591 * 00592 WORK( 1 ) = RPVGRW 00593 * 00594 * Set INFO = N+1 if the matrix is singular to working precision. 00595 * 00596 IF( RCOND.LT.DLAMCH( 'Epsilon' ) ) 00597 $ INFO = N + 1 00598 RETURN 00599 * 00600 * End of DGESVX 00601 * 00602 END