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