00001 /* sgesvx.f -- translated by f2c (version 20061008). 00002 You must link the resulting object file with libf2c: 00003 on Microsoft Windows system, link with libf2c.lib; 00004 on Linux or Unix systems, link with .../path/to/libf2c.a -lm 00005 or, if you install libf2c.a in a standard place, with -lf2c -lm 00006 -- in that order, at the end of the command line, as in 00007 cc *.o -lf2c -lm 00008 Source for libf2c is in /netlib/f2c/libf2c.zip, e.g., 00009 00010 http://www.netlib.org/f2c/libf2c.zip 00011 */ 00012 00013 #include "f2c.h" 00014 #include "blaswrap.h" 00015 00016 /* Subroutine */ int sgesvx_(char *fact, char *trans, integer *n, integer * 00017 nrhs, real *a, integer *lda, real *af, integer *ldaf, integer *ipiv, 00018 char *equed, real *r__, real *c__, real *b, integer *ldb, real *x, 00019 integer *ldx, real *rcond, real *ferr, real *berr, real *work, 00020 integer *iwork, integer *info) 00021 { 00022 /* System generated locals */ 00023 integer a_dim1, a_offset, af_dim1, af_offset, b_dim1, b_offset, x_dim1, 00024 x_offset, i__1, i__2; 00025 real r__1, r__2; 00026 00027 /* Local variables */ 00028 integer i__, j; 00029 real amax; 00030 char norm[1]; 00031 extern logical lsame_(char *, char *); 00032 real rcmin, rcmax, anorm; 00033 logical equil; 00034 real colcnd; 00035 extern doublereal slamch_(char *), slange_(char *, integer *, 00036 integer *, real *, integer *, real *); 00037 logical nofact; 00038 extern /* Subroutine */ int slaqge_(integer *, integer *, real *, integer 00039 *, real *, real *, real *, real *, real *, char *), 00040 xerbla_(char *, integer *), sgecon_(char *, integer *, 00041 real *, integer *, real *, real *, real *, integer *, integer *); 00042 real bignum; 00043 integer infequ; 00044 logical colequ; 00045 extern /* Subroutine */ int sgeequ_(integer *, integer *, real *, integer 00046 *, real *, real *, real *, real *, real *, integer *), sgerfs_( 00047 char *, integer *, integer *, real *, integer *, real *, integer * 00048 , integer *, real *, integer *, real *, integer *, real *, real *, 00049 real *, integer *, integer *), sgetrf_(integer *, 00050 integer *, real *, integer *, integer *, integer *); 00051 real rowcnd; 00052 extern /* Subroutine */ int slacpy_(char *, integer *, integer *, real *, 00053 integer *, real *, integer *); 00054 logical notran; 00055 extern doublereal slantr_(char *, char *, char *, integer *, integer *, 00056 real *, integer *, real *); 00057 extern /* Subroutine */ int sgetrs_(char *, integer *, integer *, real *, 00058 integer *, integer *, real *, integer *, integer *); 00059 real smlnum; 00060 logical rowequ; 00061 real rpvgrw; 00062 00063 00064 /* -- LAPACK driver routine (version 3.2) -- */ 00065 /* Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd.. */ 00066 /* November 2006 */ 00067 00068 /* .. Scalar Arguments .. */ 00069 /* .. */ 00070 /* .. Array Arguments .. */ 00071 /* .. */ 00072 00073 /* Purpose */ 00074 /* ======= */ 00075 00076 /* SGESVX uses the LU factorization to compute the solution to a real */ 00077 /* system of linear equations */ 00078 /* A * X = B, */ 00079 /* where A is an N-by-N matrix and X and B are N-by-NRHS matrices. */ 00080 00081 /* Error bounds on the solution and a condition estimate are also */ 00082 /* provided. */ 00083 00084 /* Description */ 00085 /* =========== */ 00086 00087 /* The following steps are performed: */ 00088 00089 /* 1. If FACT = 'E', real scaling factors are computed to equilibrate */ 00090 /* the system: */ 00091 /* TRANS = 'N': diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B */ 00092 /* TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B */ 00093 /* TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B */ 00094 /* Whether or not the system will be equilibrated depends on the */ 00095 /* scaling of the matrix A, but if equilibration is used, A is */ 00096 /* overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N') */ 00097 /* or diag(C)*B (if TRANS = 'T' or 'C'). */ 00098 00099 /* 2. If FACT = 'N' or 'E', the LU decomposition is used to factor the */ 00100 /* matrix A (after equilibration if FACT = 'E') as */ 00101 /* A = P * L * U, */ 00102 /* where P is a permutation matrix, L is a unit lower triangular */ 00103 /* matrix, and U is upper triangular. */ 00104 00105 /* 3. If some U(i,i)=0, so that U is exactly singular, then the routine */ 00106 /* returns with INFO = i. Otherwise, the factored form of A is used */ 00107 /* to estimate the condition number of the matrix A. If the */ 00108 /* reciprocal of the condition number is less than machine precision, */ 00109 /* INFO = N+1 is returned as a warning, but the routine still goes on */ 00110 /* to solve for X and compute error bounds as described below. */ 00111 00112 /* 4. The system of equations is solved for X using the factored form */ 00113 /* of A. */ 00114 00115 /* 5. Iterative refinement is applied to improve the computed solution */ 00116 /* matrix and calculate error bounds and backward error estimates */ 00117 /* for it. */ 00118 00119 /* 6. If equilibration was used, the matrix X is premultiplied by */ 00120 /* diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so */ 00121 /* that it solves the original system before equilibration. */ 00122 00123 /* Arguments */ 00124 /* ========= */ 00125 00126 /* FACT (input) CHARACTER*1 */ 00127 /* Specifies whether or not the factored form of the matrix A is */ 00128 /* supplied on entry, and if not, whether the matrix A should be */ 00129 /* equilibrated before it is factored. */ 00130 /* = 'F': On entry, AF and IPIV contain the factored form of A. */ 00131 /* If EQUED is not 'N', the matrix A has been */ 00132 /* equilibrated with scaling factors given by R and C. */ 00133 /* A, AF, and IPIV are not modified. */ 00134 /* = 'N': The matrix A will be copied to AF and factored. */ 00135 /* = 'E': The matrix A will be equilibrated if necessary, then */ 00136 /* copied to AF and factored. */ 00137 00138 /* TRANS (input) CHARACTER*1 */ 00139 /* Specifies the form of the system of equations: */ 00140 /* = 'N': A * X = B (No transpose) */ 00141 /* = 'T': A**T * X = B (Transpose) */ 00142 /* = 'C': A**H * X = B (Transpose) */ 00143 00144 /* N (input) INTEGER */ 00145 /* The number of linear equations, i.e., the order of the */ 00146 /* matrix A. N >= 0. */ 00147 00148 /* NRHS (input) INTEGER */ 00149 /* The number of right hand sides, i.e., the number of columns */ 00150 /* of the matrices B and X. NRHS >= 0. */ 00151 00152 /* A (input/output) REAL array, dimension (LDA,N) */ 00153 /* On entry, the N-by-N matrix A. If FACT = 'F' and EQUED is */ 00154 /* not 'N', then A must have been equilibrated by the scaling */ 00155 /* factors in R and/or C. A is not modified if FACT = 'F' or */ 00156 /* 'N', or if FACT = 'E' and EQUED = 'N' on exit. */ 00157 00158 /* On exit, if EQUED .ne. 'N', A is scaled as follows: */ 00159 /* EQUED = 'R': A := diag(R) * A */ 00160 /* EQUED = 'C': A := A * diag(C) */ 00161 /* EQUED = 'B': A := diag(R) * A * diag(C). */ 00162 00163 /* LDA (input) INTEGER */ 00164 /* The leading dimension of the array A. LDA >= max(1,N). */ 00165 00166 /* AF (input or output) REAL array, dimension (LDAF,N) */ 00167 /* If FACT = 'F', then AF is an input argument and on entry */ 00168 /* contains the factors L and U from the factorization */ 00169 /* A = P*L*U as computed by SGETRF. If EQUED .ne. 'N', then */ 00170 /* AF is the factored form of the equilibrated matrix A. */ 00171 00172 /* If FACT = 'N', then AF is an output argument and on exit */ 00173 /* returns the factors L and U from the factorization A = P*L*U */ 00174 /* of the original matrix A. */ 00175 00176 /* If FACT = 'E', then AF is an output argument and on exit */ 00177 /* returns the factors L and U from the factorization A = P*L*U */ 00178 /* of the equilibrated matrix A (see the description of A for */ 00179 /* the form of the equilibrated matrix). */ 00180 00181 /* LDAF (input) INTEGER */ 00182 /* The leading dimension of the array AF. LDAF >= max(1,N). */ 00183 00184 /* IPIV (input or output) INTEGER array, dimension (N) */ 00185 /* If FACT = 'F', then IPIV is an input argument and on entry */ 00186 /* contains the pivot indices from the factorization A = P*L*U */ 00187 /* as computed by SGETRF; row i of the matrix was interchanged */ 00188 /* with row IPIV(i). */ 00189 00190 /* If FACT = 'N', then IPIV is an output argument and on exit */ 00191 /* contains the pivot indices from the factorization A = P*L*U */ 00192 /* of the original matrix A. */ 00193 00194 /* If FACT = 'E', then IPIV is an output argument and on exit */ 00195 /* contains the pivot indices from the factorization A = P*L*U */ 00196 /* of the equilibrated matrix A. */ 00197 00198 /* EQUED (input or output) 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 00210 /* R (input or output) REAL array, dimension (N) */ 00211 /* The row scale factors for A. If EQUED = 'R' or 'B', A is */ 00212 /* multiplied on the left by diag(R); if EQUED = 'N' or 'C', R */ 00213 /* is not accessed. R is an input argument if FACT = 'F'; */ 00214 /* otherwise, R is an output argument. If FACT = 'F' and */ 00215 /* EQUED = 'R' or 'B', each element of R must be positive. */ 00216 00217 /* C (input or output) REAL array, dimension (N) */ 00218 /* The column scale factors for A. If EQUED = 'C' or 'B', A is */ 00219 /* multiplied on the right by diag(C); if EQUED = 'N' or 'R', C */ 00220 /* is not accessed. C is an input argument if FACT = 'F'; */ 00221 /* otherwise, C is an output argument. If FACT = 'F' and */ 00222 /* EQUED = 'C' or 'B', each element of C must be positive. */ 00223 00224 /* B (input/output) REAL array, dimension (LDB,NRHS) */ 00225 /* On entry, the N-by-NRHS right hand side matrix B. */ 00226 /* On exit, */ 00227 /* if EQUED = 'N', B is not modified; */ 00228 /* if TRANS = 'N' and EQUED = 'R' or 'B', B is overwritten by */ 00229 /* diag(R)*B; */ 00230 /* if TRANS = 'T' or 'C' and EQUED = 'C' or 'B', B is */ 00231 /* overwritten by diag(C)*B. */ 00232 00233 /* LDB (input) INTEGER */ 00234 /* The leading dimension of the array B. LDB >= max(1,N). */ 00235 00236 /* X (output) REAL array, dimension (LDX,NRHS) */ 00237 /* If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X */ 00238 /* to the original system of equations. Note that A and B are */ 00239 /* modified on exit if EQUED .ne. 'N', and the solution to the */ 00240 /* equilibrated system is inv(diag(C))*X if TRANS = 'N' and */ 00241 /* EQUED = 'C' or 'B', or inv(diag(R))*X if TRANS = 'T' or 'C' */ 00242 /* and EQUED = 'R' or 'B'. */ 00243 00244 /* LDX (input) INTEGER */ 00245 /* The leading dimension of the array X. LDX >= max(1,N). */ 00246 00247 /* RCOND (output) REAL */ 00248 /* The estimate of the reciprocal condition number of the matrix */ 00249 /* A after equilibration (if done). If RCOND is less than the */ 00250 /* machine precision (in particular, if RCOND = 0), the matrix */ 00251 /* is singular to working precision. This condition is */ 00252 /* indicated by a return code of INFO > 0. */ 00253 00254 /* FERR (output) REAL array, dimension (NRHS) */ 00255 /* The estimated forward error bound for each solution vector */ 00256 /* X(j) (the j-th column of the solution matrix X). */ 00257 /* If XTRUE is the true solution corresponding to X(j), FERR(j) */ 00258 /* is an estimated upper bound for the magnitude of the largest */ 00259 /* element in (X(j) - XTRUE) divided by the magnitude of the */ 00260 /* largest element in X(j). The estimate is as reliable as */ 00261 /* the estimate for RCOND, and is almost always a slight */ 00262 /* overestimate of the true error. */ 00263 00264 /* BERR (output) REAL array, dimension (NRHS) */ 00265 /* The componentwise relative backward error of each solution */ 00266 /* vector X(j) (i.e., the smallest relative change in */ 00267 /* any element of A or B that makes X(j) an exact solution). */ 00268 00269 /* WORK (workspace/output) REAL array, dimension (4*N) */ 00270 /* On exit, WORK(1) contains the reciprocal pivot growth */ 00271 /* factor norm(A)/norm(U). The "max absolute element" norm is */ 00272 /* used. If WORK(1) is much less than 1, then the stability */ 00273 /* of the LU factorization of the (equilibrated) matrix A */ 00274 /* could be poor. This also means that the solution X, condition */ 00275 /* estimator RCOND, and forward error bound FERR could be */ 00276 /* unreliable. If factorization fails with 0<INFO<=N, then */ 00277 /* WORK(1) contains the reciprocal pivot growth factor for the */ 00278 /* leading INFO columns of A. */ 00279 00280 /* IWORK (workspace) INTEGER array, dimension (N) */ 00281 00282 /* INFO (output) INTEGER */ 00283 /* = 0: successful exit */ 00284 /* < 0: if INFO = -i, the i-th argument had an illegal value */ 00285 /* > 0: if INFO = i, and i is */ 00286 /* <= N: U(i,i) is exactly zero. The factorization has */ 00287 /* been completed, but the factor U is exactly */ 00288 /* singular, so the solution and error bounds */ 00289 /* could not be computed. RCOND = 0 is returned. */ 00290 /* = N+1: U is nonsingular, but RCOND is less than machine */ 00291 /* precision, meaning that the matrix is singular */ 00292 /* to working precision. Nevertheless, the */ 00293 /* solution and error bounds are computed because */ 00294 /* there are a number of situations where the */ 00295 /* computed solution can be more accurate than the */ 00296 /* value of RCOND would suggest. */ 00297 00298 /* ===================================================================== */ 00299 00300 /* .. Parameters .. */ 00301 /* .. */ 00302 /* .. Local Scalars .. */ 00303 /* .. */ 00304 /* .. External Functions .. */ 00305 /* .. */ 00306 /* .. External Subroutines .. */ 00307 /* .. */ 00308 /* .. Intrinsic Functions .. */ 00309 /* .. */ 00310 /* .. Executable Statements .. */ 00311 00312 /* Parameter adjustments */ 00313 a_dim1 = *lda; 00314 a_offset = 1 + a_dim1; 00315 a -= a_offset; 00316 af_dim1 = *ldaf; 00317 af_offset = 1 + af_dim1; 00318 af -= af_offset; 00319 --ipiv; 00320 --r__; 00321 --c__; 00322 b_dim1 = *ldb; 00323 b_offset = 1 + b_dim1; 00324 b -= b_offset; 00325 x_dim1 = *ldx; 00326 x_offset = 1 + x_dim1; 00327 x -= x_offset; 00328 --ferr; 00329 --berr; 00330 --work; 00331 --iwork; 00332 00333 /* Function Body */ 00334 *info = 0; 00335 nofact = lsame_(fact, "N"); 00336 equil = lsame_(fact, "E"); 00337 notran = lsame_(trans, "N"); 00338 if (nofact || equil) { 00339 *(unsigned char *)equed = 'N'; 00340 rowequ = FALSE_; 00341 colequ = FALSE_; 00342 } else { 00343 rowequ = lsame_(equed, "R") || lsame_(equed, 00344 "B"); 00345 colequ = lsame_(equed, "C") || lsame_(equed, 00346 "B"); 00347 smlnum = slamch_("Safe minimum"); 00348 bignum = 1.f / smlnum; 00349 } 00350 00351 /* Test the input parameters. */ 00352 00353 if (! nofact && ! equil && ! lsame_(fact, "F")) { 00354 *info = -1; 00355 } else if (! notran && ! lsame_(trans, "T") && ! 00356 lsame_(trans, "C")) { 00357 *info = -2; 00358 } else if (*n < 0) { 00359 *info = -3; 00360 } else if (*nrhs < 0) { 00361 *info = -4; 00362 } else if (*lda < max(1,*n)) { 00363 *info = -6; 00364 } else if (*ldaf < max(1,*n)) { 00365 *info = -8; 00366 } else if (lsame_(fact, "F") && ! (rowequ || colequ 00367 || lsame_(equed, "N"))) { 00368 *info = -10; 00369 } else { 00370 if (rowequ) { 00371 rcmin = bignum; 00372 rcmax = 0.f; 00373 i__1 = *n; 00374 for (j = 1; j <= i__1; ++j) { 00375 /* Computing MIN */ 00376 r__1 = rcmin, r__2 = r__[j]; 00377 rcmin = dmin(r__1,r__2); 00378 /* Computing MAX */ 00379 r__1 = rcmax, r__2 = r__[j]; 00380 rcmax = dmax(r__1,r__2); 00381 /* L10: */ 00382 } 00383 if (rcmin <= 0.f) { 00384 *info = -11; 00385 } else if (*n > 0) { 00386 rowcnd = dmax(rcmin,smlnum) / dmin(rcmax,bignum); 00387 } else { 00388 rowcnd = 1.f; 00389 } 00390 } 00391 if (colequ && *info == 0) { 00392 rcmin = bignum; 00393 rcmax = 0.f; 00394 i__1 = *n; 00395 for (j = 1; j <= i__1; ++j) { 00396 /* Computing MIN */ 00397 r__1 = rcmin, r__2 = c__[j]; 00398 rcmin = dmin(r__1,r__2); 00399 /* Computing MAX */ 00400 r__1 = rcmax, r__2 = c__[j]; 00401 rcmax = dmax(r__1,r__2); 00402 /* L20: */ 00403 } 00404 if (rcmin <= 0.f) { 00405 *info = -12; 00406 } else if (*n > 0) { 00407 colcnd = dmax(rcmin,smlnum) / dmin(rcmax,bignum); 00408 } else { 00409 colcnd = 1.f; 00410 } 00411 } 00412 if (*info == 0) { 00413 if (*ldb < max(1,*n)) { 00414 *info = -14; 00415 } else if (*ldx < max(1,*n)) { 00416 *info = -16; 00417 } 00418 } 00419 } 00420 00421 if (*info != 0) { 00422 i__1 = -(*info); 00423 xerbla_("SGESVX", &i__1); 00424 return 0; 00425 } 00426 00427 if (equil) { 00428 00429 /* Compute row and column scalings to equilibrate the matrix A. */ 00430 00431 sgeequ_(n, n, &a[a_offset], lda, &r__[1], &c__[1], &rowcnd, &colcnd, & 00432 amax, &infequ); 00433 if (infequ == 0) { 00434 00435 /* Equilibrate the matrix. */ 00436 00437 slaqge_(n, n, &a[a_offset], lda, &r__[1], &c__[1], &rowcnd, & 00438 colcnd, &amax, equed); 00439 rowequ = lsame_(equed, "R") || lsame_(equed, 00440 "B"); 00441 colequ = lsame_(equed, "C") || lsame_(equed, 00442 "B"); 00443 } 00444 } 00445 00446 /* Scale the right hand side. */ 00447 00448 if (notran) { 00449 if (rowequ) { 00450 i__1 = *nrhs; 00451 for (j = 1; j <= i__1; ++j) { 00452 i__2 = *n; 00453 for (i__ = 1; i__ <= i__2; ++i__) { 00454 b[i__ + j * b_dim1] = r__[i__] * b[i__ + j * b_dim1]; 00455 /* L30: */ 00456 } 00457 /* L40: */ 00458 } 00459 } 00460 } else if (colequ) { 00461 i__1 = *nrhs; 00462 for (j = 1; j <= i__1; ++j) { 00463 i__2 = *n; 00464 for (i__ = 1; i__ <= i__2; ++i__) { 00465 b[i__ + j * b_dim1] = c__[i__] * b[i__ + j * b_dim1]; 00466 /* L50: */ 00467 } 00468 /* L60: */ 00469 } 00470 } 00471 00472 if (nofact || equil) { 00473 00474 /* Compute the LU factorization of A. */ 00475 00476 slacpy_("Full", n, n, &a[a_offset], lda, &af[af_offset], ldaf); 00477 sgetrf_(n, n, &af[af_offset], ldaf, &ipiv[1], info); 00478 00479 /* Return if INFO is non-zero. */ 00480 00481 if (*info > 0) { 00482 00483 /* Compute the reciprocal pivot growth factor of the */ 00484 /* leading rank-deficient INFO columns of A. */ 00485 00486 rpvgrw = slantr_("M", "U", "N", info, info, &af[af_offset], ldaf, 00487 &work[1]); 00488 if (rpvgrw == 0.f) { 00489 rpvgrw = 1.f; 00490 } else { 00491 rpvgrw = slange_("M", n, info, &a[a_offset], lda, &work[1]) / rpvgrw; 00492 } 00493 work[1] = rpvgrw; 00494 *rcond = 0.f; 00495 return 0; 00496 } 00497 } 00498 00499 /* Compute the norm of the matrix A and the */ 00500 /* reciprocal pivot growth factor RPVGRW. */ 00501 00502 if (notran) { 00503 *(unsigned char *)norm = '1'; 00504 } else { 00505 *(unsigned char *)norm = 'I'; 00506 } 00507 anorm = slange_(norm, n, n, &a[a_offset], lda, &work[1]); 00508 rpvgrw = slantr_("M", "U", "N", n, n, &af[af_offset], ldaf, &work[1]); 00509 if (rpvgrw == 0.f) { 00510 rpvgrw = 1.f; 00511 } else { 00512 rpvgrw = slange_("M", n, n, &a[a_offset], lda, &work[1]) / 00513 rpvgrw; 00514 } 00515 00516 /* Compute the reciprocal of the condition number of A. */ 00517 00518 sgecon_(norm, n, &af[af_offset], ldaf, &anorm, rcond, &work[1], &iwork[1], 00519 info); 00520 00521 /* Compute the solution matrix X. */ 00522 00523 slacpy_("Full", n, nrhs, &b[b_offset], ldb, &x[x_offset], ldx); 00524 sgetrs_(trans, n, nrhs, &af[af_offset], ldaf, &ipiv[1], &x[x_offset], ldx, 00525 info); 00526 00527 /* Use iterative refinement to improve the computed solution and */ 00528 /* compute error bounds and backward error estimates for it. */ 00529 00530 sgerfs_(trans, n, nrhs, &a[a_offset], lda, &af[af_offset], ldaf, &ipiv[1], 00531 &b[b_offset], ldb, &x[x_offset], ldx, &ferr[1], &berr[1], &work[ 00532 1], &iwork[1], info); 00533 00534 /* Transform the solution matrix X to a solution of the original */ 00535 /* system. */ 00536 00537 if (notran) { 00538 if (colequ) { 00539 i__1 = *nrhs; 00540 for (j = 1; j <= i__1; ++j) { 00541 i__2 = *n; 00542 for (i__ = 1; i__ <= i__2; ++i__) { 00543 x[i__ + j * x_dim1] = c__[i__] * x[i__ + j * x_dim1]; 00544 /* L70: */ 00545 } 00546 /* L80: */ 00547 } 00548 i__1 = *nrhs; 00549 for (j = 1; j <= i__1; ++j) { 00550 ferr[j] /= colcnd; 00551 /* L90: */ 00552 } 00553 } 00554 } else if (rowequ) { 00555 i__1 = *nrhs; 00556 for (j = 1; j <= i__1; ++j) { 00557 i__2 = *n; 00558 for (i__ = 1; i__ <= i__2; ++i__) { 00559 x[i__ + j * x_dim1] = r__[i__] * x[i__ + j * x_dim1]; 00560 /* L100: */ 00561 } 00562 /* L110: */ 00563 } 00564 i__1 = *nrhs; 00565 for (j = 1; j <= i__1; ++j) { 00566 ferr[j] /= rowcnd; 00567 /* L120: */ 00568 } 00569 } 00570 00571 /* Set INFO = N+1 if the matrix is singular to working precision. */ 00572 00573 if (*rcond < slamch_("Epsilon")) { 00574 *info = *n + 1; 00575 } 00576 00577 work[1] = rpvgrw; 00578 return 0; 00579 00580 /* End of SGESVX */ 00581 00582 } /* sgesvx_ */