00001 /* dgesvxx.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 dgesvxx_(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 *rpvgrw, doublereal *berr, integer *n_err_bnds__, 00021 doublereal *err_bnds_norm__, doublereal *err_bnds_comp__, integer * 00022 nparams, doublereal *params, doublereal *work, integer *iwork, 00023 integer *info) 00024 { 00025 /* System generated locals */ 00026 integer a_dim1, a_offset, af_dim1, af_offset, b_dim1, b_offset, x_dim1, 00027 x_offset, err_bnds_norm_dim1, err_bnds_norm_offset, 00028 err_bnds_comp_dim1, err_bnds_comp_offset, i__1; 00029 doublereal d__1, d__2; 00030 00031 /* Local variables */ 00032 integer j; 00033 extern doublereal dla_rpvgrw__(integer *, integer *, doublereal *, 00034 integer *, doublereal *, integer *); 00035 doublereal amax; 00036 extern logical lsame_(char *, char *); 00037 doublereal rcmin, rcmax; 00038 logical equil; 00039 extern doublereal dlamch_(char *); 00040 extern /* Subroutine */ int dlaqge_(integer *, integer *, doublereal *, 00041 integer *, doublereal *, doublereal *, doublereal *, doublereal *, 00042 doublereal *, char *); 00043 doublereal colcnd; 00044 logical nofact; 00045 extern /* Subroutine */ int dgetrf_(integer *, integer *, doublereal *, 00046 integer *, integer *, integer *), dlacpy_(char *, integer *, 00047 integer *, doublereal *, integer *, doublereal *, integer *), xerbla_(char *, integer *); 00048 doublereal bignum; 00049 integer infequ; 00050 logical colequ; 00051 extern /* Subroutine */ int dgetrs_(char *, integer *, integer *, 00052 doublereal *, integer *, integer *, doublereal *, integer *, 00053 integer *); 00054 doublereal rowcnd; 00055 logical notran; 00056 doublereal smlnum; 00057 logical rowequ; 00058 extern /* Subroutine */ int dlascl2_(integer *, integer *, doublereal *, 00059 doublereal *, integer *), dgeequb_(integer *, integer *, 00060 doublereal *, integer *, doublereal *, doublereal *, doublereal *, 00061 doublereal *, doublereal *, integer *), dgerfsx_(char *, char *, 00062 integer *, integer *, doublereal *, integer *, doublereal *, 00063 integer *, integer *, doublereal *, doublereal *, doublereal *, 00064 integer *, doublereal *, integer *, doublereal *, doublereal *, 00065 integer *, doublereal *, doublereal *, integer *, doublereal *, 00066 doublereal *, integer *, integer *); 00067 00068 00069 /* -- LAPACK driver routine (version 3.2) -- */ 00070 /* -- Contributed by James Demmel, Deaglan Halligan, Yozo Hida and -- */ 00071 /* -- Jason Riedy of Univ. of California Berkeley. -- */ 00072 /* -- November 2008 -- */ 00073 00074 /* -- LAPACK is a software package provided by Univ. of Tennessee, -- */ 00075 /* -- Univ. of California Berkeley and NAG Ltd. -- */ 00076 00077 /* .. */ 00078 /* .. Scalar Arguments .. */ 00079 /* .. */ 00080 /* .. Array Arguments .. */ 00081 /* .. */ 00082 00083 /* Purpose */ 00084 /* ======= */ 00085 00086 /* DGESVXX uses the LU factorization to compute the solution to a */ 00087 /* double precision system of linear equations A * X = B, where A is an */ 00088 /* N-by-N matrix and X and B are N-by-NRHS matrices. */ 00089 00090 /* If requested, both normwise and maximum componentwise error bounds */ 00091 /* are returned. DGESVXX will return a solution with a tiny */ 00092 /* guaranteed error (O(eps) where eps is the working machine */ 00093 /* precision) unless the matrix is very ill-conditioned, in which */ 00094 /* case a warning is returned. Relevant condition numbers also are */ 00095 /* calculated and returned. */ 00096 00097 /* DGESVXX accepts user-provided factorizations and equilibration */ 00098 /* factors; see the definitions of the FACT and EQUED options. */ 00099 /* Solving with refinement and using a factorization from a previous */ 00100 /* DGESVXX call will also produce a solution with either O(eps) */ 00101 /* errors or warnings, but we cannot make that claim for general */ 00102 /* user-provided factorizations and equilibration factors if they */ 00103 /* differ from what DGESVXX would itself produce. */ 00104 00105 /* Description */ 00106 /* =========== */ 00107 00108 /* The following steps are performed: */ 00109 00110 /* 1. If FACT = 'E', double precision scaling factors are computed to equilibrate */ 00111 /* the system: */ 00112 00113 /* TRANS = 'N': diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B */ 00114 /* TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B */ 00115 /* TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B */ 00116 00117 /* Whether or not the system will be equilibrated depends on the */ 00118 /* scaling of the matrix A, but if equilibration is used, A is */ 00119 /* overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N') */ 00120 /* or diag(C)*B (if TRANS = 'T' or 'C'). */ 00121 00122 /* 2. If FACT = 'N' or 'E', the LU decomposition is used to factor */ 00123 /* the matrix A (after equilibration if FACT = 'E') as */ 00124 00125 /* A = P * L * U, */ 00126 00127 /* where P is a permutation matrix, L is a unit lower triangular */ 00128 /* matrix, and U is upper triangular. */ 00129 00130 /* 3. If some U(i,i)=0, so that U is exactly singular, then the */ 00131 /* routine returns with INFO = i. Otherwise, the factored form of A */ 00132 /* is used to estimate the condition number of the matrix A (see */ 00133 /* argument RCOND). If the reciprocal of the condition number is less */ 00134 /* than machine precision, the routine still goes on to solve for X */ 00135 /* and compute error bounds as described below. */ 00136 00137 /* 4. The system of equations is solved for X using the factored form */ 00138 /* of A. */ 00139 00140 /* 5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero), */ 00141 /* the routine will use iterative refinement to try to get a small */ 00142 /* error and error bounds. Refinement calculates the residual to at */ 00143 /* least twice the working precision. */ 00144 00145 /* 6. If equilibration was used, the matrix X is premultiplied by */ 00146 /* diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so */ 00147 /* that it solves the original system before equilibration. */ 00148 00149 /* Arguments */ 00150 /* ========= */ 00151 00152 /* Some optional parameters are bundled in the PARAMS array. These */ 00153 /* settings determine how refinement is performed, but often the */ 00154 /* defaults are acceptable. If the defaults are acceptable, users */ 00155 /* can pass NPARAMS = 0 which prevents the source code from accessing */ 00156 /* the PARAMS argument. */ 00157 00158 /* FACT (input) CHARACTER*1 */ 00159 /* Specifies whether or not the factored form of the matrix A is */ 00160 /* supplied on entry, and if not, whether the matrix A should be */ 00161 /* equilibrated before it is factored. */ 00162 /* = 'F': On entry, AF and IPIV contain the factored form of A. */ 00163 /* If EQUED is not 'N', the matrix A has been */ 00164 /* equilibrated with scaling factors given by R and C. */ 00165 /* A, AF, and IPIV are not modified. */ 00166 /* = 'N': The matrix A will be copied to AF and factored. */ 00167 /* = 'E': The matrix A will be equilibrated if necessary, then */ 00168 /* copied to AF and factored. */ 00169 00170 /* TRANS (input) CHARACTER*1 */ 00171 /* Specifies the form of the system of equations: */ 00172 /* = 'N': A * X = B (No transpose) */ 00173 /* = 'T': A**T * X = B (Transpose) */ 00174 /* = 'C': A**H * X = B (Conjugate Transpose = Transpose) */ 00175 00176 /* N (input) INTEGER */ 00177 /* The number of linear equations, i.e., the order of the */ 00178 /* matrix A. N >= 0. */ 00179 00180 /* NRHS (input) INTEGER */ 00181 /* The number of right hand sides, i.e., the number of columns */ 00182 /* of the matrices B and X. NRHS >= 0. */ 00183 00184 /* A (input/output) DOUBLE PRECISION array, dimension (LDA,N) */ 00185 /* On entry, the N-by-N matrix A. If FACT = 'F' and EQUED is */ 00186 /* not 'N', then A must have been equilibrated by the scaling */ 00187 /* factors in R and/or C. A is not modified if FACT = 'F' or */ 00188 /* 'N', or if FACT = 'E' and EQUED = 'N' on exit. */ 00189 00190 /* On exit, if EQUED .ne. 'N', A is scaled as follows: */ 00191 /* EQUED = 'R': A := diag(R) * A */ 00192 /* EQUED = 'C': A := A * diag(C) */ 00193 /* EQUED = 'B': A := diag(R) * A * diag(C). */ 00194 00195 /* LDA (input) INTEGER */ 00196 /* The leading dimension of the array A. LDA >= max(1,N). */ 00197 00198 /* AF (input or output) DOUBLE PRECISION array, dimension (LDAF,N) */ 00199 /* If FACT = 'F', then AF is an input argument and on entry */ 00200 /* contains the factors L and U from the factorization */ 00201 /* A = P*L*U as computed by DGETRF. If EQUED .ne. 'N', then */ 00202 /* AF is the factored form of the equilibrated matrix A. */ 00203 00204 /* If FACT = 'N', then AF is an output argument and on exit */ 00205 /* returns the factors L and U from the factorization A = P*L*U */ 00206 /* of the original matrix A. */ 00207 00208 /* If FACT = 'E', then AF is an output argument and on exit */ 00209 /* returns the factors L and U from the factorization A = P*L*U */ 00210 /* of the equilibrated matrix A (see the description of A for */ 00211 /* the form of the equilibrated matrix). */ 00212 00213 /* LDAF (input) INTEGER */ 00214 /* The leading dimension of the array AF. LDAF >= max(1,N). */ 00215 00216 /* IPIV (input or output) INTEGER array, dimension (N) */ 00217 /* If FACT = 'F', then IPIV is an input argument and on entry */ 00218 /* contains the pivot indices from the factorization A = P*L*U */ 00219 /* as computed by DGETRF; row i of the matrix was interchanged */ 00220 /* with row IPIV(i). */ 00221 00222 /* If FACT = 'N', then IPIV is an output argument and on exit */ 00223 /* contains the pivot indices from the factorization A = P*L*U */ 00224 /* of the original matrix A. */ 00225 00226 /* If FACT = 'E', then IPIV is an output argument and on exit */ 00227 /* contains the pivot indices from the factorization A = P*L*U */ 00228 /* of the equilibrated matrix A. */ 00229 00230 /* EQUED (input or output) CHARACTER*1 */ 00231 /* Specifies the form of equilibration that was done. */ 00232 /* = 'N': No equilibration (always true if FACT = 'N'). */ 00233 /* = 'R': Row equilibration, i.e., A has been premultiplied by */ 00234 /* diag(R). */ 00235 /* = 'C': Column equilibration, i.e., A has been postmultiplied */ 00236 /* by diag(C). */ 00237 /* = 'B': Both row and column equilibration, i.e., A has been */ 00238 /* replaced by diag(R) * A * diag(C). */ 00239 /* EQUED is an input argument if FACT = 'F'; otherwise, it is an */ 00240 /* output argument. */ 00241 00242 /* R (input or output) DOUBLE PRECISION array, dimension (N) */ 00243 /* The row scale factors for A. If EQUED = 'R' or 'B', A is */ 00244 /* multiplied on the left by diag(R); if EQUED = 'N' or 'C', R */ 00245 /* is not accessed. R is an input argument if FACT = 'F'; */ 00246 /* otherwise, R is an output argument. If FACT = 'F' and */ 00247 /* EQUED = 'R' or 'B', each element of R must be positive. */ 00248 /* If R is output, each element of R is a power of the radix. */ 00249 /* If R is input, each element of R should be a power of the radix */ 00250 /* to ensure a reliable solution and error estimates. Scaling by */ 00251 /* powers of the radix does not cause rounding errors unless the */ 00252 /* result underflows or overflows. Rounding errors during scaling */ 00253 /* lead to refining with a matrix that is not equivalent to the */ 00254 /* input matrix, producing error estimates that may not be */ 00255 /* reliable. */ 00256 00257 /* C (input or output) DOUBLE PRECISION array, dimension (N) */ 00258 /* The column scale factors for A. If EQUED = 'C' or 'B', A is */ 00259 /* multiplied on the right by diag(C); if EQUED = 'N' or 'R', C */ 00260 /* is not accessed. C is an input argument if FACT = 'F'; */ 00261 /* otherwise, C is an output argument. If FACT = 'F' and */ 00262 /* EQUED = 'C' or 'B', each element of C must be positive. */ 00263 /* If C is output, each element of C is a power of the radix. */ 00264 /* If C is input, each element of C should be a power of the radix */ 00265 /* to ensure a reliable solution and error estimates. Scaling by */ 00266 /* powers of the radix does not cause rounding errors unless the */ 00267 /* result underflows or overflows. Rounding errors during scaling */ 00268 /* lead to refining with a matrix that is not equivalent to the */ 00269 /* input matrix, producing error estimates that may not be */ 00270 /* reliable. */ 00271 00272 /* B (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS) */ 00273 /* On entry, the N-by-NRHS right hand side matrix B. */ 00274 /* On exit, */ 00275 /* if EQUED = 'N', B is not modified; */ 00276 /* if TRANS = 'N' and EQUED = 'R' or 'B', B is overwritten by */ 00277 /* diag(R)*B; */ 00278 /* if TRANS = 'T' or 'C' and EQUED = 'C' or 'B', B is */ 00279 /* overwritten by diag(C)*B. */ 00280 00281 /* LDB (input) INTEGER */ 00282 /* The leading dimension of the array B. LDB >= max(1,N). */ 00283 00284 /* X (output) DOUBLE PRECISION array, dimension (LDX,NRHS) */ 00285 /* If INFO = 0, the N-by-NRHS solution matrix X to the original */ 00286 /* system of equations. Note that A and B are modified on exit */ 00287 /* if EQUED .ne. 'N', and the solution to the equilibrated system is */ 00288 /* inv(diag(C))*X if TRANS = 'N' and EQUED = 'C' or 'B', or */ 00289 /* inv(diag(R))*X if TRANS = 'T' or 'C' and EQUED = 'R' or 'B'. */ 00290 00291 /* LDX (input) INTEGER */ 00292 /* The leading dimension of the array X. LDX >= max(1,N). */ 00293 00294 /* RCOND (output) DOUBLE PRECISION */ 00295 /* Reciprocal scaled condition number. This is an estimate of the */ 00296 /* reciprocal Skeel condition number of the matrix A after */ 00297 /* equilibration (if done). If this is less than the machine */ 00298 /* precision (in particular, if it is zero), the matrix is singular */ 00299 /* to working precision. Note that the error may still be small even */ 00300 /* if this number is very small and the matrix appears ill- */ 00301 /* conditioned. */ 00302 00303 /* RPVGRW (output) DOUBLE PRECISION */ 00304 /* Reciprocal pivot growth. On exit, this contains the reciprocal */ 00305 /* pivot growth factor norm(A)/norm(U). The "max absolute element" */ 00306 /* norm is used. If this is much less than 1, then the stability of */ 00307 /* the LU factorization of the (equilibrated) matrix A could be poor. */ 00308 /* This also means that the solution X, estimated condition numbers, */ 00309 /* and error bounds could be unreliable. If factorization fails with */ 00310 /* 0<INFO<=N, then this contains the reciprocal pivot growth factor */ 00311 /* for the leading INFO columns of A. In DGESVX, this quantity is */ 00312 /* returned in WORK(1). */ 00313 00314 /* BERR (output) DOUBLE PRECISION array, dimension (NRHS) */ 00315 /* Componentwise relative backward error. This is the */ 00316 /* componentwise relative backward error of each solution vector X(j) */ 00317 /* (i.e., the smallest relative change in any element of A or B that */ 00318 /* makes X(j) an exact solution). */ 00319 00320 /* N_ERR_BNDS (input) INTEGER */ 00321 /* Number of error bounds to return for each right hand side */ 00322 /* and each type (normwise or componentwise). See ERR_BNDS_NORM and */ 00323 /* ERR_BNDS_COMP below. */ 00324 00325 /* ERR_BNDS_NORM (output) DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS) */ 00326 /* For each right-hand side, this array contains information about */ 00327 /* various error bounds and condition numbers corresponding to the */ 00328 /* normwise relative error, which is defined as follows: */ 00329 00330 /* Normwise relative error in the ith solution vector: */ 00331 /* max_j (abs(XTRUE(j,i) - X(j,i))) */ 00332 /* ------------------------------ */ 00333 /* max_j abs(X(j,i)) */ 00334 00335 /* The array is indexed by the type of error information as described */ 00336 /* below. There currently are up to three pieces of information */ 00337 /* returned. */ 00338 00339 /* The first index in ERR_BNDS_NORM(i,:) corresponds to the ith */ 00340 /* right-hand side. */ 00341 00342 /* The second index in ERR_BNDS_NORM(:,err) contains the following */ 00343 /* three fields: */ 00344 /* err = 1 "Trust/don't trust" boolean. Trust the answer if the */ 00345 /* reciprocal condition number is less than the threshold */ 00346 /* sqrt(n) * dlamch('Epsilon'). */ 00347 00348 /* err = 2 "Guaranteed" error bound: The estimated forward error, */ 00349 /* almost certainly within a factor of 10 of the true error */ 00350 /* so long as the next entry is greater than the threshold */ 00351 /* sqrt(n) * dlamch('Epsilon'). This error bound should only */ 00352 /* be trusted if the previous boolean is true. */ 00353 00354 /* err = 3 Reciprocal condition number: Estimated normwise */ 00355 /* reciprocal condition number. Compared with the threshold */ 00356 /* sqrt(n) * dlamch('Epsilon') to determine if the error */ 00357 /* estimate is "guaranteed". These reciprocal condition */ 00358 /* numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some */ 00359 /* appropriately scaled matrix Z. */ 00360 /* Let Z = S*A, where S scales each row by a power of the */ 00361 /* radix so all absolute row sums of Z are approximately 1. */ 00362 00363 /* See Lapack Working Note 165 for further details and extra */ 00364 /* cautions. */ 00365 00366 /* ERR_BNDS_COMP (output) DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS) */ 00367 /* For each right-hand side, this array contains information about */ 00368 /* various error bounds and condition numbers corresponding to the */ 00369 /* componentwise relative error, which is defined as follows: */ 00370 00371 /* Componentwise relative error in the ith solution vector: */ 00372 /* abs(XTRUE(j,i) - X(j,i)) */ 00373 /* max_j ---------------------- */ 00374 /* abs(X(j,i)) */ 00375 00376 /* The array is indexed by the right-hand side i (on which the */ 00377 /* componentwise relative error depends), and the type of error */ 00378 /* information as described below. There currently are up to three */ 00379 /* pieces of information returned for each right-hand side. If */ 00380 /* componentwise accuracy is not requested (PARAMS(3) = 0.0), then */ 00381 /* ERR_BNDS_COMP is not accessed. If N_ERR_BNDS .LT. 3, then at most */ 00382 /* the first (:,N_ERR_BNDS) entries are returned. */ 00383 00384 /* The first index in ERR_BNDS_COMP(i,:) corresponds to the ith */ 00385 /* right-hand side. */ 00386 00387 /* The second index in ERR_BNDS_COMP(:,err) contains the following */ 00388 /* three fields: */ 00389 /* err = 1 "Trust/don't trust" boolean. Trust the answer if the */ 00390 /* reciprocal condition number is less than the threshold */ 00391 /* sqrt(n) * dlamch('Epsilon'). */ 00392 00393 /* err = 2 "Guaranteed" error bound: The estimated forward error, */ 00394 /* almost certainly within a factor of 10 of the true error */ 00395 /* so long as the next entry is greater than the threshold */ 00396 /* sqrt(n) * dlamch('Epsilon'). This error bound should only */ 00397 /* be trusted if the previous boolean is true. */ 00398 00399 /* err = 3 Reciprocal condition number: Estimated componentwise */ 00400 /* reciprocal condition number. Compared with the threshold */ 00401 /* sqrt(n) * dlamch('Epsilon') to determine if the error */ 00402 /* estimate is "guaranteed". These reciprocal condition */ 00403 /* numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some */ 00404 /* appropriately scaled matrix Z. */ 00405 /* Let Z = S*(A*diag(x)), where x is the solution for the */ 00406 /* current right-hand side and S scales each row of */ 00407 /* A*diag(x) by a power of the radix so all absolute row */ 00408 /* sums of Z are approximately 1. */ 00409 00410 /* See Lapack Working Note 165 for further details and extra */ 00411 /* cautions. */ 00412 00413 /* NPARAMS (input) INTEGER */ 00414 /* Specifies the number of parameters set in PARAMS. If .LE. 0, the */ 00415 /* PARAMS array is never referenced and default values are used. */ 00416 00417 /* PARAMS (input / output) DOUBLE PRECISION array, dimension NPARAMS */ 00418 /* Specifies algorithm parameters. If an entry is .LT. 0.0, then */ 00419 /* that entry will be filled with default value used for that */ 00420 /* parameter. Only positions up to NPARAMS are accessed; defaults */ 00421 /* are used for higher-numbered parameters. */ 00422 00423 /* PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative */ 00424 /* refinement or not. */ 00425 /* Default: 1.0D+0 */ 00426 /* = 0.0 : No refinement is performed, and no error bounds are */ 00427 /* computed. */ 00428 /* = 1.0 : Use the extra-precise refinement algorithm. */ 00429 /* (other values are reserved for future use) */ 00430 00431 /* PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual */ 00432 /* computations allowed for refinement. */ 00433 /* Default: 10 */ 00434 /* Aggressive: Set to 100 to permit convergence using approximate */ 00435 /* factorizations or factorizations other than LU. If */ 00436 /* the factorization uses a technique other than */ 00437 /* Gaussian elimination, the guarantees in */ 00438 /* err_bnds_norm and err_bnds_comp may no longer be */ 00439 /* trustworthy. */ 00440 00441 /* PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code */ 00442 /* will attempt to find a solution with small componentwise */ 00443 /* relative error in the double-precision algorithm. Positive */ 00444 /* is true, 0.0 is false. */ 00445 /* Default: 1.0 (attempt componentwise convergence) */ 00446 00447 /* WORK (workspace) DOUBLE PRECISION array, dimension (4*N) */ 00448 00449 /* IWORK (workspace) INTEGER array, dimension (N) */ 00450 00451 /* INFO (output) INTEGER */ 00452 /* = 0: Successful exit. The solution to every right-hand side is */ 00453 /* guaranteed. */ 00454 /* < 0: If INFO = -i, the i-th argument had an illegal value */ 00455 /* > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization */ 00456 /* has been completed, but the factor U is exactly singular, so */ 00457 /* the solution and error bounds could not be computed. RCOND = 0 */ 00458 /* is returned. */ 00459 /* = N+J: The solution corresponding to the Jth right-hand side is */ 00460 /* not guaranteed. The solutions corresponding to other right- */ 00461 /* hand sides K with K > J may not be guaranteed as well, but */ 00462 /* only the first such right-hand side is reported. If a small */ 00463 /* componentwise error is not requested (PARAMS(3) = 0.0) then */ 00464 /* the Jth right-hand side is the first with a normwise error */ 00465 /* bound that is not guaranteed (the smallest J such */ 00466 /* that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0) */ 00467 /* the Jth right-hand side is the first with either a normwise or */ 00468 /* componentwise error bound that is not guaranteed (the smallest */ 00469 /* J such that either ERR_BNDS_NORM(J,1) = 0.0 or */ 00470 /* ERR_BNDS_COMP(J,1) = 0.0). See the definition of */ 00471 /* ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information */ 00472 /* about all of the right-hand sides check ERR_BNDS_NORM or */ 00473 /* ERR_BNDS_COMP. */ 00474 00475 /* ================================================================== */ 00476 00477 /* .. Parameters .. */ 00478 /* .. */ 00479 /* .. Local Scalars .. */ 00480 /* .. */ 00481 /* .. External Functions .. */ 00482 /* .. */ 00483 /* .. External Subroutines .. */ 00484 /* .. */ 00485 /* .. Intrinsic Functions .. */ 00486 /* .. */ 00487 /* .. Executable Statements .. */ 00488 00489 /* Parameter adjustments */ 00490 err_bnds_comp_dim1 = *nrhs; 00491 err_bnds_comp_offset = 1 + err_bnds_comp_dim1; 00492 err_bnds_comp__ -= err_bnds_comp_offset; 00493 err_bnds_norm_dim1 = *nrhs; 00494 err_bnds_norm_offset = 1 + err_bnds_norm_dim1; 00495 err_bnds_norm__ -= err_bnds_norm_offset; 00496 a_dim1 = *lda; 00497 a_offset = 1 + a_dim1; 00498 a -= a_offset; 00499 af_dim1 = *ldaf; 00500 af_offset = 1 + af_dim1; 00501 af -= af_offset; 00502 --ipiv; 00503 --r__; 00504 --c__; 00505 b_dim1 = *ldb; 00506 b_offset = 1 + b_dim1; 00507 b -= b_offset; 00508 x_dim1 = *ldx; 00509 x_offset = 1 + x_dim1; 00510 x -= x_offset; 00511 --berr; 00512 --params; 00513 --work; 00514 --iwork; 00515 00516 /* Function Body */ 00517 *info = 0; 00518 nofact = lsame_(fact, "N"); 00519 equil = lsame_(fact, "E"); 00520 notran = lsame_(trans, "N"); 00521 smlnum = dlamch_("Safe minimum"); 00522 bignum = 1. / smlnum; 00523 if (nofact || equil) { 00524 *(unsigned char *)equed = 'N'; 00525 rowequ = FALSE_; 00526 colequ = FALSE_; 00527 } else { 00528 rowequ = lsame_(equed, "R") || lsame_(equed, 00529 "B"); 00530 colequ = lsame_(equed, "C") || lsame_(equed, 00531 "B"); 00532 } 00533 00534 /* Default is failure. If an input parameter is wrong or */ 00535 /* factorization fails, make everything look horrible. Only the */ 00536 /* pivot growth is set here, the rest is initialized in DGERFSX. */ 00537 00538 *rpvgrw = 0.; 00539 00540 /* Test the input parameters. PARAMS is not tested until DGERFSX. */ 00541 00542 if (! nofact && ! equil && ! lsame_(fact, "F")) { 00543 *info = -1; 00544 } else if (! notran && ! lsame_(trans, "T") && ! 00545 lsame_(trans, "C")) { 00546 *info = -2; 00547 } else if (*n < 0) { 00548 *info = -3; 00549 } else if (*nrhs < 0) { 00550 *info = -4; 00551 } else if (*lda < max(1,*n)) { 00552 *info = -6; 00553 } else if (*ldaf < max(1,*n)) { 00554 *info = -8; 00555 } else if (lsame_(fact, "F") && ! (rowequ || colequ 00556 || lsame_(equed, "N"))) { 00557 *info = -10; 00558 } else { 00559 if (rowequ) { 00560 rcmin = bignum; 00561 rcmax = 0.; 00562 i__1 = *n; 00563 for (j = 1; j <= i__1; ++j) { 00564 /* Computing MIN */ 00565 d__1 = rcmin, d__2 = r__[j]; 00566 rcmin = min(d__1,d__2); 00567 /* Computing MAX */ 00568 d__1 = rcmax, d__2 = r__[j]; 00569 rcmax = max(d__1,d__2); 00570 /* L10: */ 00571 } 00572 if (rcmin <= 0.) { 00573 *info = -11; 00574 } else if (*n > 0) { 00575 rowcnd = max(rcmin,smlnum) / min(rcmax,bignum); 00576 } else { 00577 rowcnd = 1.; 00578 } 00579 } 00580 if (colequ && *info == 0) { 00581 rcmin = bignum; 00582 rcmax = 0.; 00583 i__1 = *n; 00584 for (j = 1; j <= i__1; ++j) { 00585 /* Computing MIN */ 00586 d__1 = rcmin, d__2 = c__[j]; 00587 rcmin = min(d__1,d__2); 00588 /* Computing MAX */ 00589 d__1 = rcmax, d__2 = c__[j]; 00590 rcmax = max(d__1,d__2); 00591 /* L20: */ 00592 } 00593 if (rcmin <= 0.) { 00594 *info = -12; 00595 } else if (*n > 0) { 00596 colcnd = max(rcmin,smlnum) / min(rcmax,bignum); 00597 } else { 00598 colcnd = 1.; 00599 } 00600 } 00601 if (*info == 0) { 00602 if (*ldb < max(1,*n)) { 00603 *info = -14; 00604 } else if (*ldx < max(1,*n)) { 00605 *info = -16; 00606 } 00607 } 00608 } 00609 00610 if (*info != 0) { 00611 i__1 = -(*info); 00612 xerbla_("DGESVXX", &i__1); 00613 return 0; 00614 } 00615 00616 if (equil) { 00617 00618 /* Compute row and column scalings to equilibrate the matrix A. */ 00619 00620 dgeequb_(n, n, &a[a_offset], lda, &r__[1], &c__[1], &rowcnd, &colcnd, 00621 &amax, &infequ); 00622 if (infequ == 0) { 00623 00624 /* Equilibrate the matrix. */ 00625 00626 dlaqge_(n, n, &a[a_offset], lda, &r__[1], &c__[1], &rowcnd, & 00627 colcnd, &amax, equed); 00628 rowequ = lsame_(equed, "R") || lsame_(equed, 00629 "B"); 00630 colequ = lsame_(equed, "C") || lsame_(equed, 00631 "B"); 00632 } 00633 00634 /* If the scaling factors are not applied, set them to 1.0. */ 00635 00636 if (! rowequ) { 00637 i__1 = *n; 00638 for (j = 1; j <= i__1; ++j) { 00639 r__[j] = 1.; 00640 } 00641 } 00642 if (! colequ) { 00643 i__1 = *n; 00644 for (j = 1; j <= i__1; ++j) { 00645 c__[j] = 1.; 00646 } 00647 } 00648 } 00649 00650 /* Scale the right-hand side. */ 00651 00652 if (notran) { 00653 if (rowequ) { 00654 dlascl2_(n, nrhs, &r__[1], &b[b_offset], ldb); 00655 } 00656 } else { 00657 if (colequ) { 00658 dlascl2_(n, nrhs, &c__[1], &b[b_offset], ldb); 00659 } 00660 } 00661 00662 if (nofact || equil) { 00663 00664 /* Compute the LU factorization of A. */ 00665 00666 dlacpy_("Full", n, n, &a[a_offset], lda, &af[af_offset], ldaf); 00667 dgetrf_(n, n, &af[af_offset], ldaf, &ipiv[1], info); 00668 00669 /* Return if INFO is non-zero. */ 00670 00671 if (*info > 0) { 00672 00673 /* Pivot in column INFO is exactly 0 */ 00674 /* Compute the reciprocal pivot growth factor of the */ 00675 /* leading rank-deficient INFO columns of A. */ 00676 00677 *rpvgrw = dla_rpvgrw__(n, info, &a[a_offset], lda, &af[af_offset], 00678 ldaf); 00679 return 0; 00680 } 00681 } 00682 00683 /* Compute the reciprocal pivot growth factor RPVGRW. */ 00684 00685 *rpvgrw = dla_rpvgrw__(n, n, &a[a_offset], lda, &af[af_offset], ldaf); 00686 00687 /* Compute the solution matrix X. */ 00688 00689 dlacpy_("Full", n, nrhs, &b[b_offset], ldb, &x[x_offset], ldx); 00690 dgetrs_(trans, n, nrhs, &af[af_offset], ldaf, &ipiv[1], &x[x_offset], ldx, 00691 info); 00692 00693 /* Use iterative refinement to improve the computed solution and */ 00694 /* compute error bounds and backward error estimates for it. */ 00695 00696 dgerfsx_(trans, equed, n, nrhs, &a[a_offset], lda, &af[af_offset], ldaf, & 00697 ipiv[1], &r__[1], &c__[1], &b[b_offset], ldb, &x[x_offset], ldx, 00698 rcond, &berr[1], n_err_bnds__, &err_bnds_norm__[ 00699 err_bnds_norm_offset], &err_bnds_comp__[err_bnds_comp_offset], 00700 nparams, ¶ms[1], &work[1], &iwork[1], info); 00701 00702 /* Scale solutions. */ 00703 00704 if (colequ && notran) { 00705 dlascl2_(n, nrhs, &c__[1], &x[x_offset], ldx); 00706 } else if (rowequ && ! notran) { 00707 dlascl2_(n, nrhs, &r__[1], &x[x_offset], ldx); 00708 } 00709 00710 return 0; 00711 00712 /* End of DGESVXX */ 00713 } /* dgesvxx_ */