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