Annotation of rpl/lapack/lapack/dgesvx.f, revision 1.3
1.1 bertrand 1: SUBROUTINE DGESVX( FACT, TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV,
2: $ EQUED, R, C, B, LDB, X, LDX, RCOND, FERR, BERR,
3: $ WORK, IWORK, INFO )
4: *
5: * -- LAPACK driver routine (version 3.2) --
6: * -- LAPACK is a software package provided by Univ. of Tennessee, --
7: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
8: * November 2006
9: *
10: * .. Scalar Arguments ..
11: CHARACTER EQUED, FACT, TRANS
12: INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
13: DOUBLE PRECISION RCOND
14: * ..
15: * .. Array Arguments ..
16: INTEGER IPIV( * ), IWORK( * )
17: DOUBLE PRECISION A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
18: $ BERR( * ), C( * ), FERR( * ), R( * ),
19: $ WORK( * ), X( LDX, * )
20: * ..
21: *
22: * Purpose
23: * =======
24: *
25: * DGESVX uses the LU factorization to compute the solution to a real
26: * system of linear equations
27: * A * X = B,
28: * where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
29: *
30: * Error bounds on the solution and a condition estimate are also
31: * provided.
32: *
33: * Description
34: * ===========
35: *
36: * The following steps are performed:
37: *
38: * 1. If FACT = 'E', real scaling factors are computed to equilibrate
39: * the system:
40: * TRANS = 'N': diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B
41: * TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
42: * TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B
43: * Whether or not the system will be equilibrated depends on the
44: * scaling of the matrix A, but if equilibration is used, A is
45: * overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N')
46: * or diag(C)*B (if TRANS = 'T' or 'C').
47: *
48: * 2. If FACT = 'N' or 'E', the LU decomposition is used to factor the
49: * matrix A (after equilibration if FACT = 'E') as
50: * A = P * L * U,
51: * where P is a permutation matrix, L is a unit lower triangular
52: * matrix, and U is upper triangular.
53: *
54: * 3. If some U(i,i)=0, so that U is exactly singular, then the routine
55: * returns with INFO = i. Otherwise, the factored form of A is used
56: * to estimate the condition number of the matrix A. If the
57: * reciprocal of the condition number is less than machine precision,
58: * INFO = N+1 is returned as a warning, but the routine still goes on
59: * to solve for X and compute error bounds as described below.
60: *
61: * 4. The system of equations is solved for X using the factored form
62: * of A.
63: *
64: * 5. Iterative refinement is applied to improve the computed solution
65: * matrix and calculate error bounds and backward error estimates
66: * for it.
67: *
68: * 6. If equilibration was used, the matrix X is premultiplied by
69: * diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so
70: * that it solves the original system before equilibration.
71: *
72: * Arguments
73: * =========
74: *
75: * FACT (input) CHARACTER*1
76: * Specifies whether or not the factored form of the matrix A is
77: * supplied on entry, and if not, whether the matrix A should be
78: * equilibrated before it is factored.
79: * = 'F': On entry, AF and IPIV contain the factored form of A.
80: * If EQUED is not 'N', the matrix A has been
81: * equilibrated with scaling factors given by R and C.
82: * A, AF, and IPIV are not modified.
83: * = 'N': The matrix A will be copied to AF and factored.
84: * = 'E': The matrix A will be equilibrated if necessary, then
85: * copied to AF and factored.
86: *
87: * TRANS (input) CHARACTER*1
88: * Specifies the form of the system of equations:
89: * = 'N': A * X = B (No transpose)
90: * = 'T': A**T * X = B (Transpose)
91: * = 'C': A**H * X = B (Transpose)
92: *
93: * N (input) INTEGER
94: * The number of linear equations, i.e., the order of the
95: * matrix A. N >= 0.
96: *
97: * NRHS (input) INTEGER
98: * The number of right hand sides, i.e., the number of columns
99: * of the matrices B and X. NRHS >= 0.
100: *
101: * A (input/output) DOUBLE PRECISION array, dimension (LDA,N)
102: * On entry, the N-by-N matrix A. If FACT = 'F' and EQUED is
103: * not 'N', then A must have been equilibrated by the scaling
104: * factors in R and/or C. A is not modified if FACT = 'F' or
105: * 'N', or if FACT = 'E' and EQUED = 'N' on exit.
106: *
107: * On exit, if EQUED .ne. 'N', A is scaled as follows:
108: * EQUED = 'R': A := diag(R) * A
109: * EQUED = 'C': A := A * diag(C)
110: * EQUED = 'B': A := diag(R) * A * diag(C).
111: *
112: * LDA (input) INTEGER
113: * The leading dimension of the array A. LDA >= max(1,N).
114: *
115: * AF (input or output) DOUBLE PRECISION array, dimension (LDAF,N)
116: * If FACT = 'F', then AF is an input argument and on entry
117: * contains the factors L and U from the factorization
118: * A = P*L*U as computed by DGETRF. If EQUED .ne. 'N', then
119: * AF is the factored form of the equilibrated matrix A.
120: *
121: * If FACT = 'N', then AF is an output argument and on exit
122: * returns the factors L and U from the factorization A = P*L*U
123: * of the original matrix A.
124: *
125: * If FACT = 'E', then AF is an output argument and on exit
126: * returns the factors L and U from the factorization A = P*L*U
127: * of the equilibrated matrix A (see the description of A for
128: * the form of the equilibrated matrix).
129: *
130: * LDAF (input) INTEGER
131: * The leading dimension of the array AF. LDAF >= max(1,N).
132: *
133: * IPIV (input or output) INTEGER array, dimension (N)
134: * If FACT = 'F', then IPIV is an input argument and on entry
135: * contains the pivot indices from the factorization A = P*L*U
136: * as computed by DGETRF; row i of the matrix was interchanged
137: * with row IPIV(i).
138: *
139: * If FACT = 'N', then IPIV is an output argument and on exit
140: * contains the pivot indices from the factorization A = P*L*U
141: * of the original matrix A.
142: *
143: * If FACT = 'E', then IPIV is an output argument and on exit
144: * contains the pivot indices from the factorization A = P*L*U
145: * of the equilibrated matrix A.
146: *
147: * EQUED (input or output) CHARACTER*1
148: * Specifies the form of equilibration that was done.
149: * = 'N': No equilibration (always true if FACT = 'N').
150: * = 'R': Row equilibration, i.e., A has been premultiplied by
151: * diag(R).
152: * = 'C': Column equilibration, i.e., A has been postmultiplied
153: * by diag(C).
154: * = 'B': Both row and column equilibration, i.e., A has been
155: * replaced by diag(R) * A * diag(C).
156: * EQUED is an input argument if FACT = 'F'; otherwise, it is an
157: * output argument.
158: *
159: * R (input or output) DOUBLE PRECISION array, dimension (N)
160: * The row scale factors for A. If EQUED = 'R' or 'B', A is
161: * multiplied on the left by diag(R); if EQUED = 'N' or 'C', R
162: * is not accessed. R is an input argument if FACT = 'F';
163: * otherwise, R is an output argument. If FACT = 'F' and
164: * EQUED = 'R' or 'B', each element of R must be positive.
165: *
166: * C (input or output) DOUBLE PRECISION array, dimension (N)
167: * The column scale factors for A. If EQUED = 'C' or 'B', A is
168: * multiplied on the right by diag(C); if EQUED = 'N' or 'R', C
169: * is not accessed. C is an input argument if FACT = 'F';
170: * otherwise, C is an output argument. If FACT = 'F' and
171: * EQUED = 'C' or 'B', each element of C must be positive.
172: *
173: * B (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS)
174: * On entry, the N-by-NRHS right hand side matrix B.
175: * On exit,
176: * if EQUED = 'N', B is not modified;
177: * if TRANS = 'N' and EQUED = 'R' or 'B', B is overwritten by
178: * diag(R)*B;
179: * if TRANS = 'T' or 'C' and EQUED = 'C' or 'B', B is
180: * overwritten by diag(C)*B.
181: *
182: * LDB (input) INTEGER
183: * The leading dimension of the array B. LDB >= max(1,N).
184: *
185: * X (output) DOUBLE PRECISION array, dimension (LDX,NRHS)
186: * If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X
187: * to the original system of equations. Note that A and B are
188: * modified on exit if EQUED .ne. 'N', and the solution to the
189: * equilibrated system is inv(diag(C))*X if TRANS = 'N' and
190: * EQUED = 'C' or 'B', or inv(diag(R))*X if TRANS = 'T' or 'C'
191: * and EQUED = 'R' or 'B'.
192: *
193: * LDX (input) INTEGER
194: * The leading dimension of the array X. LDX >= max(1,N).
195: *
196: * RCOND (output) DOUBLE PRECISION
197: * The estimate of the reciprocal condition number of the matrix
198: * A after equilibration (if done). If RCOND is less than the
199: * machine precision (in particular, if RCOND = 0), the matrix
200: * is singular to working precision. This condition is
201: * indicated by a return code of INFO > 0.
202: *
203: * FERR (output) DOUBLE PRECISION array, dimension (NRHS)
204: * The estimated forward error bound for each solution vector
205: * X(j) (the j-th column of the solution matrix X).
206: * If XTRUE is the true solution corresponding to X(j), FERR(j)
207: * is an estimated upper bound for the magnitude of the largest
208: * element in (X(j) - XTRUE) divided by the magnitude of the
209: * largest element in X(j). The estimate is as reliable as
210: * the estimate for RCOND, and is almost always a slight
211: * overestimate of the true error.
212: *
213: * BERR (output) DOUBLE PRECISION array, dimension (NRHS)
214: * The componentwise relative backward error of each solution
215: * vector X(j) (i.e., the smallest relative change in
216: * any element of A or B that makes X(j) an exact solution).
217: *
218: * WORK (workspace/output) DOUBLE PRECISION array, dimension (4*N)
219: * On exit, WORK(1) contains the reciprocal pivot growth
220: * factor norm(A)/norm(U). The "max absolute element" norm is
221: * used. If WORK(1) is much less than 1, then the stability
222: * of the LU factorization of the (equilibrated) matrix A
223: * could be poor. This also means that the solution X, condition
224: * estimator RCOND, and forward error bound FERR could be
225: * unreliable. If factorization fails with 0<INFO<=N, then
226: * WORK(1) contains the reciprocal pivot growth factor for the
227: * leading INFO columns of A.
228: *
229: * IWORK (workspace) INTEGER array, dimension (N)
230: *
231: * INFO (output) INTEGER
232: * = 0: successful exit
233: * < 0: if INFO = -i, the i-th argument had an illegal value
234: * > 0: if INFO = i, and i is
235: * <= N: U(i,i) is exactly zero. The factorization has
236: * been completed, but the factor U is exactly
237: * singular, so the solution and error bounds
238: * could not be computed. RCOND = 0 is returned.
239: * = N+1: U is nonsingular, but RCOND is less than machine
240: * precision, meaning that the matrix is singular
241: * to working precision. Nevertheless, the
242: * solution and error bounds are computed because
243: * there are a number of situations where the
244: * computed solution can be more accurate than the
245: * value of RCOND would suggest.
246: *
247: * =====================================================================
248: *
249: * .. Parameters ..
250: DOUBLE PRECISION ZERO, ONE
251: PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
252: * ..
253: * .. Local Scalars ..
254: LOGICAL COLEQU, EQUIL, NOFACT, NOTRAN, ROWEQU
255: CHARACTER NORM
256: INTEGER I, INFEQU, J
257: DOUBLE PRECISION AMAX, ANORM, BIGNUM, COLCND, RCMAX, RCMIN,
258: $ ROWCND, RPVGRW, SMLNUM
259: * ..
260: * .. External Functions ..
261: LOGICAL LSAME
262: DOUBLE PRECISION DLAMCH, DLANGE, DLANTR
263: EXTERNAL LSAME, DLAMCH, DLANGE, DLANTR
264: * ..
265: * .. External Subroutines ..
266: EXTERNAL DGECON, DGEEQU, DGERFS, DGETRF, DGETRS, DLACPY,
267: $ DLAQGE, XERBLA
268: * ..
269: * .. Intrinsic Functions ..
270: INTRINSIC MAX, MIN
271: * ..
272: * .. Executable Statements ..
273: *
274: INFO = 0
275: NOFACT = LSAME( FACT, 'N' )
276: EQUIL = LSAME( FACT, 'E' )
277: NOTRAN = LSAME( TRANS, 'N' )
278: IF( NOFACT .OR. EQUIL ) THEN
279: EQUED = 'N'
280: ROWEQU = .FALSE.
281: COLEQU = .FALSE.
282: ELSE
283: ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
284: COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
285: SMLNUM = DLAMCH( 'Safe minimum' )
286: BIGNUM = ONE / SMLNUM
287: END IF
288: *
289: * Test the input parameters.
290: *
291: IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.LSAME( FACT, 'F' ) )
292: $ THEN
293: INFO = -1
294: ELSE IF( .NOT.NOTRAN .AND. .NOT.LSAME( TRANS, 'T' ) .AND. .NOT.
295: $ LSAME( TRANS, 'C' ) ) THEN
296: INFO = -2
297: ELSE IF( N.LT.0 ) THEN
298: INFO = -3
299: ELSE IF( NRHS.LT.0 ) THEN
300: INFO = -4
301: ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
302: INFO = -6
303: ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN
304: INFO = -8
305: ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
306: $ ( ROWEQU .OR. COLEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
307: INFO = -10
308: ELSE
309: IF( ROWEQU ) THEN
310: RCMIN = BIGNUM
311: RCMAX = ZERO
312: DO 10 J = 1, N
313: RCMIN = MIN( RCMIN, R( J ) )
314: RCMAX = MAX( RCMAX, R( J ) )
315: 10 CONTINUE
316: IF( RCMIN.LE.ZERO ) THEN
317: INFO = -11
318: ELSE IF( N.GT.0 ) THEN
319: ROWCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
320: ELSE
321: ROWCND = ONE
322: END IF
323: END IF
324: IF( COLEQU .AND. INFO.EQ.0 ) THEN
325: RCMIN = BIGNUM
326: RCMAX = ZERO
327: DO 20 J = 1, N
328: RCMIN = MIN( RCMIN, C( J ) )
329: RCMAX = MAX( RCMAX, C( J ) )
330: 20 CONTINUE
331: IF( RCMIN.LE.ZERO ) THEN
332: INFO = -12
333: ELSE IF( N.GT.0 ) THEN
334: COLCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
335: ELSE
336: COLCND = ONE
337: END IF
338: END IF
339: IF( INFO.EQ.0 ) THEN
340: IF( LDB.LT.MAX( 1, N ) ) THEN
341: INFO = -14
342: ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
343: INFO = -16
344: END IF
345: END IF
346: END IF
347: *
348: IF( INFO.NE.0 ) THEN
349: CALL XERBLA( 'DGESVX', -INFO )
350: RETURN
351: END IF
352: *
353: IF( EQUIL ) THEN
354: *
355: * Compute row and column scalings to equilibrate the matrix A.
356: *
357: CALL DGEEQU( N, N, A, LDA, R, C, ROWCND, COLCND, AMAX, INFEQU )
358: IF( INFEQU.EQ.0 ) THEN
359: *
360: * Equilibrate the matrix.
361: *
362: CALL DLAQGE( N, N, A, LDA, R, C, ROWCND, COLCND, AMAX,
363: $ EQUED )
364: ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
365: COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
366: END IF
367: END IF
368: *
369: * Scale the right hand side.
370: *
371: IF( NOTRAN ) THEN
372: IF( ROWEQU ) THEN
373: DO 40 J = 1, NRHS
374: DO 30 I = 1, N
375: B( I, J ) = R( I )*B( I, J )
376: 30 CONTINUE
377: 40 CONTINUE
378: END IF
379: ELSE IF( COLEQU ) THEN
380: DO 60 J = 1, NRHS
381: DO 50 I = 1, N
382: B( I, J ) = C( I )*B( I, J )
383: 50 CONTINUE
384: 60 CONTINUE
385: END IF
386: *
387: IF( NOFACT .OR. EQUIL ) THEN
388: *
389: * Compute the LU factorization of A.
390: *
391: CALL DLACPY( 'Full', N, N, A, LDA, AF, LDAF )
392: CALL DGETRF( N, N, AF, LDAF, IPIV, INFO )
393: *
394: * Return if INFO is non-zero.
395: *
396: IF( INFO.GT.0 ) THEN
397: *
398: * Compute the reciprocal pivot growth factor of the
399: * leading rank-deficient INFO columns of A.
400: *
401: RPVGRW = DLANTR( 'M', 'U', 'N', INFO, INFO, AF, LDAF,
402: $ WORK )
403: IF( RPVGRW.EQ.ZERO ) THEN
404: RPVGRW = ONE
405: ELSE
406: RPVGRW = DLANGE( 'M', N, INFO, A, LDA, WORK ) / RPVGRW
407: END IF
408: WORK( 1 ) = RPVGRW
409: RCOND = ZERO
410: RETURN
411: END IF
412: END IF
413: *
414: * Compute the norm of the matrix A and the
415: * reciprocal pivot growth factor RPVGRW.
416: *
417: IF( NOTRAN ) THEN
418: NORM = '1'
419: ELSE
420: NORM = 'I'
421: END IF
422: ANORM = DLANGE( NORM, N, N, A, LDA, WORK )
423: RPVGRW = DLANTR( 'M', 'U', 'N', N, N, AF, LDAF, WORK )
424: IF( RPVGRW.EQ.ZERO ) THEN
425: RPVGRW = ONE
426: ELSE
427: RPVGRW = DLANGE( 'M', N, N, A, LDA, WORK ) / RPVGRW
428: END IF
429: *
430: * Compute the reciprocal of the condition number of A.
431: *
432: CALL DGECON( NORM, N, AF, LDAF, ANORM, RCOND, WORK, IWORK, INFO )
433: *
434: * Compute the solution matrix X.
435: *
436: CALL DLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
437: CALL DGETRS( TRANS, N, NRHS, AF, LDAF, IPIV, X, LDX, INFO )
438: *
439: * Use iterative refinement to improve the computed solution and
440: * compute error bounds and backward error estimates for it.
441: *
442: CALL DGERFS( TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV, B, LDB, X,
443: $ LDX, FERR, BERR, WORK, IWORK, INFO )
444: *
445: * Transform the solution matrix X to a solution of the original
446: * system.
447: *
448: IF( NOTRAN ) THEN
449: IF( COLEQU ) THEN
450: DO 80 J = 1, NRHS
451: DO 70 I = 1, N
452: X( I, J ) = C( I )*X( I, J )
453: 70 CONTINUE
454: 80 CONTINUE
455: DO 90 J = 1, NRHS
456: FERR( J ) = FERR( J ) / COLCND
457: 90 CONTINUE
458: END IF
459: ELSE IF( ROWEQU ) THEN
460: DO 110 J = 1, NRHS
461: DO 100 I = 1, N
462: X( I, J ) = R( I )*X( I, J )
463: 100 CONTINUE
464: 110 CONTINUE
465: DO 120 J = 1, NRHS
466: FERR( J ) = FERR( J ) / ROWCND
467: 120 CONTINUE
468: END IF
469: *
470: WORK( 1 ) = RPVGRW
471: *
472: * Set INFO = N+1 if the matrix is singular to working precision.
473: *
474: IF( RCOND.LT.DLAMCH( 'Epsilon' ) )
475: $ INFO = N + 1
476: RETURN
477: *
478: * End of DGESVX
479: *
480: END
CVSweb interface <joel.bertrand@systella.fr>