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