Annotation of rpl/lapack/lapack/dgesvxx.f, revision 1.4
1.1 bertrand 1: SUBROUTINE DGESVXX( FACT, TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV,
2: $ EQUED, R, C, B, LDB, X, LDX, RCOND, RPVGRW,
3: $ BERR, N_ERR_BNDS, ERR_BNDS_NORM,
4: $ ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, IWORK,
5: $ INFO )
6: *
7: * -- LAPACK driver routine (version 3.2.2) --
8: * -- Contributed by James Demmel, Deaglan Halligan, Yozo Hida and --
9: * -- Jason Riedy of Univ. of California Berkeley. --
10: * -- June 2010 --
11: *
12: * -- LAPACK is a software package provided by Univ. of Tennessee, --
13: * -- Univ. of California Berkeley and NAG Ltd. --
14: *
15: IMPLICIT NONE
16: * ..
17: * .. Scalar Arguments ..
18: CHARACTER EQUED, FACT, TRANS
19: INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
20: $ N_ERR_BNDS
21: DOUBLE PRECISION RCOND, RPVGRW
22: * ..
23: * .. Array Arguments ..
24: INTEGER IPIV( * ), IWORK( * )
25: DOUBLE PRECISION A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
26: $ X( LDX , * ),WORK( * )
27: DOUBLE PRECISION R( * ), C( * ), PARAMS( * ), BERR( * ),
28: $ ERR_BNDS_NORM( NRHS, * ),
29: $ ERR_BNDS_COMP( NRHS, * )
30: * ..
31: *
32: * Purpose
33: * =======
34: *
35: * DGESVXX uses the LU factorization to compute the solution to a
36: * double precision system of linear equations A * X = B, where A is an
37: * N-by-N matrix and X and B are N-by-NRHS matrices.
38: *
39: * If requested, both normwise and maximum componentwise error bounds
40: * are returned. DGESVXX will return a solution with a tiny
41: * guaranteed error (O(eps) where eps is the working machine
42: * precision) unless the matrix is very ill-conditioned, in which
43: * case a warning is returned. Relevant condition numbers also are
44: * calculated and returned.
45: *
46: * DGESVXX accepts user-provided factorizations and equilibration
47: * factors; see the definitions of the FACT and EQUED options.
48: * Solving with refinement and using a factorization from a previous
49: * DGESVXX call will also produce a solution with either O(eps)
50: * errors or warnings, but we cannot make that claim for general
51: * user-provided factorizations and equilibration factors if they
52: * differ from what DGESVXX would itself produce.
53: *
54: * Description
55: * ===========
56: *
57: * The following steps are performed:
58: *
59: * 1. If FACT = 'E', double precision scaling factors are computed to equilibrate
60: * the system:
61: *
62: * TRANS = 'N': diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B
63: * TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
64: * TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B
65: *
66: * Whether or not the system will be equilibrated depends on the
67: * scaling of the matrix A, but if equilibration is used, A is
68: * overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N')
69: * or diag(C)*B (if TRANS = 'T' or 'C').
70: *
71: * 2. If FACT = 'N' or 'E', the LU decomposition is used to factor
72: * the matrix A (after equilibration if FACT = 'E') as
73: *
74: * A = P * L * U,
75: *
76: * where P is a permutation matrix, L is a unit lower triangular
77: * matrix, and U is upper triangular.
78: *
79: * 3. If some U(i,i)=0, so that U is exactly singular, then the
80: * routine returns with INFO = i. Otherwise, the factored form of A
81: * is used to estimate the condition number of the matrix A (see
82: * argument RCOND). If the reciprocal of the condition number is less
83: * than machine precision, the routine still goes on to solve for X
84: * and compute error bounds as described below.
85: *
86: * 4. The system of equations is solved for X using the factored form
87: * of A.
88: *
89: * 5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
90: * the routine will use iterative refinement to try to get a small
91: * error and error bounds. Refinement calculates the residual to at
92: * least twice the working precision.
93: *
94: * 6. If equilibration was used, the matrix X is premultiplied by
95: * diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so
96: * that it solves the original system before equilibration.
97: *
98: * Arguments
99: * =========
100: *
101: * Some optional parameters are bundled in the PARAMS array. These
102: * settings determine how refinement is performed, but often the
103: * defaults are acceptable. If the defaults are acceptable, users
104: * can pass NPARAMS = 0 which prevents the source code from accessing
105: * the PARAMS argument.
106: *
107: * FACT (input) CHARACTER*1
108: * Specifies whether or not the factored form of the matrix A is
109: * supplied on entry, and if not, whether the matrix A should be
110: * equilibrated before it is factored.
111: * = 'F': On entry, AF and IPIV contain the factored form of A.
112: * If EQUED is not 'N', the matrix A has been
113: * equilibrated with scaling factors given by R and C.
114: * A, AF, and IPIV are not modified.
115: * = 'N': The matrix A will be copied to AF and factored.
116: * = 'E': The matrix A will be equilibrated if necessary, then
117: * copied to AF and factored.
118: *
119: * TRANS (input) CHARACTER*1
120: * Specifies the form of the system of equations:
121: * = 'N': A * X = B (No transpose)
122: * = 'T': A**T * X = B (Transpose)
123: * = 'C': A**H * X = B (Conjugate Transpose = Transpose)
124: *
125: * N (input) INTEGER
126: * The number of linear equations, i.e., the order of the
127: * matrix A. N >= 0.
128: *
129: * NRHS (input) INTEGER
130: * The number of right hand sides, i.e., the number of columns
131: * of the matrices B and X. NRHS >= 0.
132: *
133: * A (input/output) DOUBLE PRECISION array, dimension (LDA,N)
134: * On entry, the N-by-N matrix A. If FACT = 'F' and EQUED is
135: * not 'N', then A must have been equilibrated by the scaling
136: * factors in R and/or C. A is not modified if FACT = 'F' or
137: * 'N', or if FACT = 'E' and EQUED = 'N' on exit.
138: *
139: * On exit, if EQUED .ne. 'N', A is scaled as follows:
140: * EQUED = 'R': A := diag(R) * A
141: * EQUED = 'C': A := A * diag(C)
142: * EQUED = 'B': A := diag(R) * A * diag(C).
143: *
144: * LDA (input) INTEGER
145: * The leading dimension of the array A. LDA >= max(1,N).
146: *
147: * AF (input or output) DOUBLE PRECISION array, dimension (LDAF,N)
148: * If FACT = 'F', then AF is an input argument and on entry
149: * contains the factors L and U from the factorization
150: * A = P*L*U as computed by DGETRF. If EQUED .ne. 'N', then
151: * AF is the factored form of the equilibrated matrix A.
152: *
153: * If FACT = 'N', then AF is an output argument and on exit
154: * returns the factors L and U from the factorization A = P*L*U
155: * of the original matrix A.
156: *
157: * If FACT = 'E', then AF is an output argument and on exit
158: * returns the factors L and U from the factorization A = P*L*U
159: * of the equilibrated matrix A (see the description of A for
160: * the form of the equilibrated matrix).
161: *
162: * LDAF (input) INTEGER
163: * The leading dimension of the array AF. LDAF >= max(1,N).
164: *
165: * IPIV (input or output) INTEGER array, dimension (N)
166: * If FACT = 'F', then IPIV is an input argument and on entry
167: * contains the pivot indices from the factorization A = P*L*U
168: * as computed by DGETRF; row i of the matrix was interchanged
169: * with row IPIV(i).
170: *
171: * If FACT = 'N', then IPIV is an output argument and on exit
172: * contains the pivot indices from the factorization A = P*L*U
173: * of the original matrix A.
174: *
175: * If FACT = 'E', then IPIV is an output argument and on exit
176: * contains the pivot indices from the factorization A = P*L*U
177: * of the equilibrated matrix A.
178: *
179: * EQUED (input or output) CHARACTER*1
180: * Specifies the form of equilibration that was done.
181: * = 'N': No equilibration (always true if FACT = 'N').
182: * = 'R': Row equilibration, i.e., A has been premultiplied by
183: * diag(R).
184: * = 'C': Column equilibration, i.e., A has been postmultiplied
185: * by diag(C).
186: * = 'B': Both row and column equilibration, i.e., A has been
187: * replaced by diag(R) * A * diag(C).
188: * EQUED is an input argument if FACT = 'F'; otherwise, it is an
189: * output argument.
190: *
191: * R (input or output) DOUBLE PRECISION array, dimension (N)
192: * The row scale factors for A. If EQUED = 'R' or 'B', A is
193: * multiplied on the left by diag(R); if EQUED = 'N' or 'C', R
194: * is not accessed. R is an input argument if FACT = 'F';
195: * otherwise, R is an output argument. If FACT = 'F' and
196: * EQUED = 'R' or 'B', each element of R must be positive.
197: * If R is output, each element of R is a power of the radix.
198: * If R is input, each element of R should be a power of the radix
199: * to ensure a reliable solution and error estimates. Scaling by
200: * powers of the radix does not cause rounding errors unless the
201: * result underflows or overflows. Rounding errors during scaling
202: * lead to refining with a matrix that is not equivalent to the
203: * input matrix, producing error estimates that may not be
204: * reliable.
205: *
206: * C (input or output) DOUBLE PRECISION array, dimension (N)
207: * The column scale factors for A. If EQUED = 'C' or 'B', A is
208: * multiplied on the right by diag(C); if EQUED = 'N' or 'R', C
209: * is not accessed. C is an input argument if FACT = 'F';
210: * otherwise, C is an output argument. If FACT = 'F' and
211: * EQUED = 'C' or 'B', each element of C must be positive.
212: * If C is output, each element of C is a power of the radix.
213: * If C is input, each element of C should be a power of the radix
214: * to ensure a reliable solution and error estimates. Scaling by
215: * powers of the radix does not cause rounding errors unless the
216: * result underflows or overflows. Rounding errors during scaling
217: * lead to refining with a matrix that is not equivalent to the
218: * input matrix, producing error estimates that may not be
219: * reliable.
220: *
221: * B (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
222: * On entry, the N-by-NRHS right hand side matrix B.
223: * On exit,
224: * if EQUED = 'N', B is not modified;
225: * if TRANS = 'N' and EQUED = 'R' or 'B', B is overwritten by
226: * diag(R)*B;
227: * if TRANS = 'T' or 'C' and EQUED = 'C' or 'B', B is
228: * overwritten by diag(C)*B.
229: *
230: * LDB (input) INTEGER
231: * The leading dimension of the array B. LDB >= max(1,N).
232: *
233: * X (output) DOUBLE PRECISION array, dimension (LDX,NRHS)
234: * If INFO = 0, the N-by-NRHS solution matrix X to the original
235: * system of equations. Note that A and B are modified on exit
236: * if EQUED .ne. 'N', and the solution to the equilibrated system is
237: * inv(diag(C))*X if TRANS = 'N' and EQUED = 'C' or 'B', or
238: * inv(diag(R))*X if TRANS = 'T' or 'C' and EQUED = 'R' or 'B'.
239: *
240: * LDX (input) INTEGER
241: * The leading dimension of the array X. LDX >= max(1,N).
242: *
243: * RCOND (output) DOUBLE PRECISION
244: * Reciprocal scaled condition number. This is an estimate of the
245: * reciprocal Skeel condition number of the matrix A after
246: * equilibration (if done). If this is less than the machine
247: * precision (in particular, if it is zero), the matrix is singular
248: * to working precision. Note that the error may still be small even
249: * if this number is very small and the matrix appears ill-
250: * conditioned.
251: *
252: * RPVGRW (output) DOUBLE PRECISION
253: * Reciprocal pivot growth. On exit, this contains the reciprocal
254: * pivot growth factor norm(A)/norm(U). The "max absolute element"
255: * norm is used. If this is much less than 1, then the stability of
256: * the LU factorization of the (equilibrated) matrix A could be poor.
257: * This also means that the solution X, estimated condition numbers,
258: * and error bounds could be unreliable. If factorization fails with
259: * 0<INFO<=N, then this contains the reciprocal pivot growth factor
260: * for the leading INFO columns of A. In DGESVX, this quantity is
261: * returned in WORK(1).
262: *
263: * BERR (output) DOUBLE PRECISION array, dimension (NRHS)
264: * Componentwise relative backward error. This is the
265: * componentwise relative backward error of each solution vector X(j)
266: * (i.e., the smallest relative change in any element of A or B that
267: * makes X(j) an exact solution).
268: *
269: * N_ERR_BNDS (input) INTEGER
270: * Number of error bounds to return for each right hand side
271: * and each type (normwise or componentwise). See ERR_BNDS_NORM and
272: * ERR_BNDS_COMP below.
273: *
274: * ERR_BNDS_NORM (output) DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
275: * For each right-hand side, this array contains information about
276: * various error bounds and condition numbers corresponding to the
277: * normwise relative error, which is defined as follows:
278: *
279: * Normwise relative error in the ith solution vector:
280: * max_j (abs(XTRUE(j,i) - X(j,i)))
281: * ------------------------------
282: * max_j abs(X(j,i))
283: *
284: * The array is indexed by the type of error information as described
285: * below. There currently are up to three pieces of information
286: * returned.
287: *
288: * The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
289: * right-hand side.
290: *
291: * The second index in ERR_BNDS_NORM(:,err) contains the following
292: * three fields:
293: * err = 1 "Trust/don't trust" boolean. Trust the answer if the
294: * reciprocal condition number is less than the threshold
295: * sqrt(n) * dlamch('Epsilon').
296: *
297: * err = 2 "Guaranteed" error bound: The estimated forward error,
298: * almost certainly within a factor of 10 of the true error
299: * so long as the next entry is greater than the threshold
300: * sqrt(n) * dlamch('Epsilon'). This error bound should only
301: * be trusted if the previous boolean is true.
302: *
303: * err = 3 Reciprocal condition number: Estimated normwise
304: * reciprocal condition number. Compared with the threshold
305: * sqrt(n) * dlamch('Epsilon') to determine if the error
306: * estimate is "guaranteed". These reciprocal condition
307: * numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
308: * appropriately scaled matrix Z.
309: * Let Z = S*A, where S scales each row by a power of the
310: * radix so all absolute row sums of Z are approximately 1.
311: *
312: * See Lapack Working Note 165 for further details and extra
313: * cautions.
314: *
315: * ERR_BNDS_COMP (output) DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
316: * For each right-hand side, this array contains information about
317: * various error bounds and condition numbers corresponding to the
318: * componentwise relative error, which is defined as follows:
319: *
320: * Componentwise relative error in the ith solution vector:
321: * abs(XTRUE(j,i) - X(j,i))
322: * max_j ----------------------
323: * abs(X(j,i))
324: *
325: * The array is indexed by the right-hand side i (on which the
326: * componentwise relative error depends), and the type of error
327: * information as described below. There currently are up to three
328: * pieces of information returned for each right-hand side. If
329: * componentwise accuracy is not requested (PARAMS(3) = 0.0), then
330: * ERR_BNDS_COMP is not accessed. If N_ERR_BNDS .LT. 3, then at most
331: * the first (:,N_ERR_BNDS) entries are returned.
332: *
333: * The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
334: * right-hand side.
335: *
336: * The second index in ERR_BNDS_COMP(:,err) contains the following
337: * three fields:
338: * err = 1 "Trust/don't trust" boolean. Trust the answer if the
339: * reciprocal condition number is less than the threshold
340: * sqrt(n) * dlamch('Epsilon').
341: *
342: * err = 2 "Guaranteed" error bound: The estimated forward error,
343: * almost certainly within a factor of 10 of the true error
344: * so long as the next entry is greater than the threshold
345: * sqrt(n) * dlamch('Epsilon'). This error bound should only
346: * be trusted if the previous boolean is true.
347: *
348: * err = 3 Reciprocal condition number: Estimated componentwise
349: * reciprocal condition number. Compared with the threshold
350: * sqrt(n) * dlamch('Epsilon') to determine if the error
351: * estimate is "guaranteed". These reciprocal condition
352: * numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
353: * appropriately scaled matrix Z.
354: * Let Z = S*(A*diag(x)), where x is the solution for the
355: * current right-hand side and S scales each row of
356: * A*diag(x) by a power of the radix so all absolute row
357: * sums of Z are approximately 1.
358: *
359: * See Lapack Working Note 165 for further details and extra
360: * cautions.
361: *
362: * NPARAMS (input) INTEGER
363: * Specifies the number of parameters set in PARAMS. If .LE. 0, the
364: * PARAMS array is never referenced and default values are used.
365: *
366: * PARAMS (input / output) DOUBLE PRECISION array, dimension (NPARAMS)
367: * Specifies algorithm parameters. If an entry is .LT. 0.0, then
368: * that entry will be filled with default value used for that
369: * parameter. Only positions up to NPARAMS are accessed; defaults
370: * are used for higher-numbered parameters.
371: *
372: * PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
373: * refinement or not.
374: * Default: 1.0D+0
375: * = 0.0 : No refinement is performed, and no error bounds are
376: * computed.
377: * = 1.0 : Use the extra-precise refinement algorithm.
378: * (other values are reserved for future use)
379: *
380: * PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
381: * computations allowed for refinement.
382: * Default: 10
383: * Aggressive: Set to 100 to permit convergence using approximate
384: * factorizations or factorizations other than LU. If
385: * the factorization uses a technique other than
386: * Gaussian elimination, the guarantees in
387: * err_bnds_norm and err_bnds_comp may no longer be
388: * trustworthy.
389: *
390: * PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
391: * will attempt to find a solution with small componentwise
392: * relative error in the double-precision algorithm. Positive
393: * is true, 0.0 is false.
394: * Default: 1.0 (attempt componentwise convergence)
395: *
396: * WORK (workspace) DOUBLE PRECISION array, dimension (4*N)
397: *
398: * IWORK (workspace) INTEGER array, dimension (N)
399: *
400: * INFO (output) INTEGER
401: * = 0: Successful exit. The solution to every right-hand side is
402: * guaranteed.
403: * < 0: If INFO = -i, the i-th argument had an illegal value
404: * > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization
405: * has been completed, but the factor U is exactly singular, so
406: * the solution and error bounds could not be computed. RCOND = 0
407: * is returned.
408: * = N+J: The solution corresponding to the Jth right-hand side is
409: * not guaranteed. The solutions corresponding to other right-
410: * hand sides K with K > J may not be guaranteed as well, but
411: * only the first such right-hand side is reported. If a small
412: * componentwise error is not requested (PARAMS(3) = 0.0) then
413: * the Jth right-hand side is the first with a normwise error
414: * bound that is not guaranteed (the smallest J such
415: * that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
416: * the Jth right-hand side is the first with either a normwise or
417: * componentwise error bound that is not guaranteed (the smallest
418: * J such that either ERR_BNDS_NORM(J,1) = 0.0 or
419: * ERR_BNDS_COMP(J,1) = 0.0). See the definition of
420: * ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
421: * about all of the right-hand sides check ERR_BNDS_NORM or
422: * ERR_BNDS_COMP.
423: *
424: * ==================================================================
425: *
426: * .. Parameters ..
427: DOUBLE PRECISION ZERO, ONE
428: PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
429: INTEGER FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
430: INTEGER RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
431: INTEGER CMP_ERR_I, PIV_GROWTH_I
432: PARAMETER ( FINAL_NRM_ERR_I = 1, FINAL_CMP_ERR_I = 2,
433: $ BERR_I = 3 )
434: PARAMETER ( RCOND_I = 4, NRM_RCOND_I = 5, NRM_ERR_I = 6 )
435: PARAMETER ( CMP_RCOND_I = 7, CMP_ERR_I = 8,
436: $ PIV_GROWTH_I = 9 )
437: * ..
438: * .. Local Scalars ..
439: LOGICAL COLEQU, EQUIL, NOFACT, NOTRAN, ROWEQU
440: INTEGER INFEQU, J
441: DOUBLE PRECISION AMAX, BIGNUM, COLCND, RCMAX, RCMIN, ROWCND,
442: $ SMLNUM
443: * ..
444: * .. External Functions ..
445: EXTERNAL LSAME, DLAMCH, DLA_RPVGRW
446: LOGICAL LSAME
447: DOUBLE PRECISION DLAMCH, DLA_RPVGRW
448: * ..
449: * .. External Subroutines ..
450: EXTERNAL DGEEQUB, DGETRF, DGETRS, DLACPY, DLAQGE,
451: $ XERBLA, DLASCL2, DGERFSX
452: * ..
453: * .. Intrinsic Functions ..
454: INTRINSIC MAX, MIN
455: * ..
456: * .. Executable Statements ..
457: *
458: INFO = 0
459: NOFACT = LSAME( FACT, 'N' )
460: EQUIL = LSAME( FACT, 'E' )
461: NOTRAN = LSAME( TRANS, 'N' )
462: SMLNUM = DLAMCH( 'Safe minimum' )
463: BIGNUM = ONE / SMLNUM
464: IF( NOFACT .OR. EQUIL ) THEN
465: EQUED = 'N'
466: ROWEQU = .FALSE.
467: COLEQU = .FALSE.
468: ELSE
469: ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
470: COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
471: END IF
472: *
473: * Default is failure. If an input parameter is wrong or
474: * factorization fails, make everything look horrible. Only the
475: * pivot growth is set here, the rest is initialized in DGERFSX.
476: *
477: RPVGRW = ZERO
478: *
479: * Test the input parameters. PARAMS is not tested until DGERFSX.
480: *
481: IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.
482: $ LSAME( FACT, 'F' ) ) THEN
483: INFO = -1
484: ELSE IF( .NOT.NOTRAN .AND. .NOT.LSAME( TRANS, 'T' ) .AND. .NOT.
485: $ LSAME( TRANS, 'C' ) ) THEN
486: INFO = -2
487: ELSE IF( N.LT.0 ) THEN
488: INFO = -3
489: ELSE IF( NRHS.LT.0 ) THEN
490: INFO = -4
491: ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
492: INFO = -6
493: ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN
494: INFO = -8
495: ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
496: $ ( ROWEQU .OR. COLEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
497: INFO = -10
498: ELSE
499: IF( ROWEQU ) THEN
500: RCMIN = BIGNUM
501: RCMAX = ZERO
502: DO 10 J = 1, N
503: RCMIN = MIN( RCMIN, R( J ) )
504: RCMAX = MAX( RCMAX, R( J ) )
505: 10 CONTINUE
506: IF( RCMIN.LE.ZERO ) THEN
507: INFO = -11
508: ELSE IF( N.GT.0 ) THEN
509: ROWCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
510: ELSE
511: ROWCND = ONE
512: END IF
513: END IF
514: IF( COLEQU .AND. INFO.EQ.0 ) THEN
515: RCMIN = BIGNUM
516: RCMAX = ZERO
517: DO 20 J = 1, N
518: RCMIN = MIN( RCMIN, C( J ) )
519: RCMAX = MAX( RCMAX, C( J ) )
520: 20 CONTINUE
521: IF( RCMIN.LE.ZERO ) THEN
522: INFO = -12
523: ELSE IF( N.GT.0 ) THEN
524: COLCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
525: ELSE
526: COLCND = ONE
527: END IF
528: END IF
529: IF( INFO.EQ.0 ) THEN
530: IF( LDB.LT.MAX( 1, N ) ) THEN
531: INFO = -14
532: ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
533: INFO = -16
534: END IF
535: END IF
536: END IF
537: *
538: IF( INFO.NE.0 ) THEN
539: CALL XERBLA( 'DGESVXX', -INFO )
540: RETURN
541: END IF
542: *
543: IF( EQUIL ) THEN
544: *
545: * Compute row and column scalings to equilibrate the matrix A.
546: *
547: CALL DGEEQUB( N, N, A, LDA, R, C, ROWCND, COLCND, AMAX,
548: $ INFEQU )
549: IF( INFEQU.EQ.0 ) THEN
550: *
551: * Equilibrate the matrix.
552: *
553: CALL DLAQGE( N, N, A, LDA, R, C, ROWCND, COLCND, AMAX,
554: $ EQUED )
555: ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
556: COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
557: END IF
558: *
559: * If the scaling factors are not applied, set them to 1.0.
560: *
561: IF ( .NOT.ROWEQU ) THEN
562: DO J = 1, N
563: R( J ) = 1.0D+0
564: END DO
565: END IF
566: IF ( .NOT.COLEQU ) THEN
567: DO J = 1, N
568: C( J ) = 1.0D+0
569: END DO
570: END IF
571: END IF
572: *
573: * Scale the right-hand side.
574: *
575: IF( NOTRAN ) THEN
576: IF( ROWEQU ) CALL DLASCL2( N, NRHS, R, B, LDB )
577: ELSE
578: IF( COLEQU ) CALL DLASCL2( N, NRHS, C, B, LDB )
579: END IF
580: *
581: IF( NOFACT .OR. EQUIL ) THEN
582: *
583: * Compute the LU factorization of A.
584: *
585: CALL DLACPY( 'Full', N, N, A, LDA, AF, LDAF )
586: CALL DGETRF( N, N, AF, LDAF, IPIV, INFO )
587: *
588: * Return if INFO is non-zero.
589: *
590: IF( INFO.GT.0 ) THEN
591: *
592: * Pivot in column INFO is exactly 0
593: * Compute the reciprocal pivot growth factor of the
594: * leading rank-deficient INFO columns of A.
595: *
596: RPVGRW = DLA_RPVGRW( N, INFO, A, LDA, AF, LDAF )
597: RETURN
598: END IF
599: END IF
600: *
601: * Compute the reciprocal pivot growth factor RPVGRW.
602: *
603: RPVGRW = DLA_RPVGRW( N, N, A, LDA, AF, LDAF )
604: *
605: * Compute the solution matrix X.
606: *
607: CALL DLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
608: CALL DGETRS( TRANS, N, NRHS, AF, LDAF, IPIV, X, LDX, INFO )
609: *
610: * Use iterative refinement to improve the computed solution and
611: * compute error bounds and backward error estimates for it.
612: *
613: CALL DGERFSX( TRANS, EQUED, N, NRHS, A, LDA, AF, LDAF,
614: $ IPIV, R, C, B, LDB, X, LDX, RCOND, BERR,
615: $ N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS,
616: $ WORK, IWORK, INFO )
617: *
618: * Scale solutions.
619: *
620: IF ( COLEQU .AND. NOTRAN ) THEN
621: CALL DLASCL2 ( N, NRHS, C, X, LDX )
622: ELSE IF ( ROWEQU .AND. .NOT.NOTRAN ) THEN
623: CALL DLASCL2 ( N, NRHS, R, X, LDX )
624: END IF
625: *
626: RETURN
627: *
628: * End of DGESVXX
629:
630: END
CVSweb interface <joel.bertrand@systella.fr>