001:       SUBROUTINE SGESVXX( FACT, TRANS, N, NRHS, A, LDA, AF, LDAF, IPIV,
002:      $                    EQUED, R, C, B, LDB, X, LDX, RCOND, RPVGRW,
003:      $                    BERR, N_ERR_BNDS, ERR_BNDS_NORM,
004:      $                    ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, IWORK,
005:      $                    INFO )
006: *
007: *     -- LAPACK driver routine (version 3.2.1)                          --
008: *     -- Contributed by James Demmel, Deaglan Halligan, Yozo Hida and --
009: *     -- Jason Riedy of Univ. of California Berkeley.                 --
010: *     -- April 2009                                                   --
011: *
012: *     -- LAPACK is a software package provided by Univ. of Tennessee, --
013: *     -- Univ. of California Berkeley and NAG Ltd.                    --
014: *
015:       IMPLICIT NONE
016: *     ..
017: *     .. Scalar Arguments ..
018:       CHARACTER          EQUED, FACT, TRANS
019:       INTEGER            INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
020:      $                   N_ERR_BNDS
021:       REAL               RCOND, RPVGRW
022: *     ..
023: *     .. Array Arguments ..
024:       INTEGER            IPIV( * ), IWORK( * )
025:       REAL               A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
026:      $                   X( LDX , * ),WORK( * )
027:       REAL               R( * ), C( * ), PARAMS( * ), BERR( * ),
028:      $                   ERR_BNDS_NORM( NRHS, * ),
029:      $                   ERR_BNDS_COMP( NRHS, * )
030: *     ..
031: *
032: *     Purpose
033: *     =======
034: *
035: *     SGESVXX uses the LU factorization to compute the solution to a
036: *     real system of linear equations  A * X = B,  where A is an
037: *     N-by-N matrix and X and B are N-by-NRHS matrices.
038: *
039: *     If requested, both normwise and maximum componentwise error bounds
040: *     are returned. SGESVXX will return a solution with a tiny
041: *     guaranteed error (O(eps) where eps is the working machine
042: *     precision) unless the matrix is very ill-conditioned, in which
043: *     case a warning is returned. Relevant condition numbers also are
044: *     calculated and returned.
045: *
046: *     SGESVXX accepts user-provided factorizations and equilibration
047: *     factors; see the definitions of the FACT and EQUED options.
048: *     Solving with refinement and using a factorization from a previous
049: *     SGESVXX call will also produce a solution with either O(eps)
050: *     errors or warnings, but we cannot make that claim for general
051: *     user-provided factorizations and equilibration factors if they
052: *     differ from what SGESVXX would itself produce.
053: *
054: *     Description
055: *     ===========
056: *
057: *     The following steps are performed:
058: *
059: *     1. If FACT = 'E', real scaling factors are computed to equilibrate
060: *     the system:
061: *
062: *       TRANS = 'N':  diag(R)*A*diag(C)     *inv(diag(C))*X = diag(R)*B
063: *       TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
064: *       TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B
065: *
066: *     Whether or not the system will be equilibrated depends on the
067: *     scaling of the matrix A, but if equilibration is used, A is
068: *     overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N')
069: *     or diag(C)*B (if TRANS = 'T' or 'C').
070: *
071: *     2. If FACT = 'N' or 'E', the LU decomposition is used to factor
072: *     the matrix A (after equilibration if FACT = 'E') as
073: *
074: *       A = P * L * U,
075: *
076: *     where P is a permutation matrix, L is a unit lower triangular
077: *     matrix, and U is upper triangular.
078: *
079: *     3. If some U(i,i)=0, so that U is exactly singular, then the
080: *     routine returns with INFO = i. Otherwise, the factored form of A
081: *     is used to estimate the condition number of the matrix A (see
082: *     argument RCOND). If the reciprocal of the condition number is less
083: *     than machine precision, the routine still goes on to solve for X
084: *     and compute error bounds as described below.
085: *
086: *     4. The system of equations is solved for X using the factored form
087: *     of A.
088: *
089: *     5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
090: *     the routine will use iterative refinement to try to get a small
091: *     error and error bounds.  Refinement calculates the residual to at
092: *     least twice the working precision.
093: *
094: *     6. If equilibration was used, the matrix X is premultiplied by
095: *     diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so
096: *     that it solves the original system before equilibration.
097: *
098: *     Arguments
099: *     =========
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) REAL 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) REAL 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 SGETRF.  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 SGETRF; 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) REAL 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) REAL 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) REAL 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) REAL 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) REAL
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) REAL
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 SGESVX, this quantity is
261: *     returned in WORK(1).
262: *
263: *     BERR    (output) REAL 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) REAL 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) * slamch('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) * slamch('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) * slamch('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) REAL 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) * slamch('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) * slamch('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) * slamch('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) REAL 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.0
375: *            = 0.0 : No refinement is performed, and no error bounds are
376: *                    computed.
377: *            = 1.0 : Use the double-precision refinement algorithm,
378: *                    possibly with doubled-single computations if the
379: *                    compilation environment does not support DOUBLE
380: *                    PRECISION.
381: *              (other values are reserved for future use)
382: *
383: *       PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
384: *            computations allowed for refinement.
385: *         Default: 10
386: *         Aggressive: Set to 100 to permit convergence using approximate
387: *                     factorizations or factorizations other than LU. If
388: *                     the factorization uses a technique other than
389: *                     Gaussian elimination, the guarantees in
390: *                     err_bnds_norm and err_bnds_comp may no longer be
391: *                     trustworthy.
392: *
393: *       PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
394: *            will attempt to find a solution with small componentwise
395: *            relative error in the double-precision algorithm.  Positive
396: *            is true, 0.0 is false.
397: *         Default: 1.0 (attempt componentwise convergence)
398: *
399: *     WORK    (workspace) REAL array, dimension (4*N)
400: *
401: *     IWORK   (workspace) INTEGER array, dimension (N)
402: *
403: *     INFO    (output) INTEGER
404: *       = 0:  Successful exit. The solution to every right-hand side is
405: *         guaranteed.
406: *       < 0:  If INFO = -i, the i-th argument had an illegal value
407: *       > 0 and <= N:  U(INFO,INFO) is exactly zero.  The factorization
408: *         has been completed, but the factor U is exactly singular, so
409: *         the solution and error bounds could not be computed. RCOND = 0
410: *         is returned.
411: *       = N+J: The solution corresponding to the Jth right-hand side is
412: *         not guaranteed. The solutions corresponding to other right-
413: *         hand sides K with K > J may not be guaranteed as well, but
414: *         only the first such right-hand side is reported. If a small
415: *         componentwise error is not requested (PARAMS(3) = 0.0) then
416: *         the Jth right-hand side is the first with a normwise error
417: *         bound that is not guaranteed (the smallest J such
418: *         that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
419: *         the Jth right-hand side is the first with either a normwise or
420: *         componentwise error bound that is not guaranteed (the smallest
421: *         J such that either ERR_BNDS_NORM(J,1) = 0.0 or
422: *         ERR_BNDS_COMP(J,1) = 0.0). See the definition of
423: *         ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
424: *         about all of the right-hand sides check ERR_BNDS_NORM or
425: *         ERR_BNDS_COMP.
426: *
427: *     ==================================================================
428: *
429: *     .. Parameters ..
430:       REAL               ZERO, ONE
431:       PARAMETER          ( ZERO = 0.0E+0, ONE = 1.0E+0 )
432:       INTEGER            FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
433:       INTEGER            RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
434:       INTEGER            CMP_ERR_I, PIV_GROWTH_I
435:       PARAMETER          ( FINAL_NRM_ERR_I = 1, FINAL_CMP_ERR_I = 2,
436:      $                   BERR_I = 3 )
437:       PARAMETER          ( RCOND_I = 4, NRM_RCOND_I = 5, NRM_ERR_I = 6 )
438:       PARAMETER          ( CMP_RCOND_I = 7, CMP_ERR_I = 8,
439:      $                   PIV_GROWTH_I = 9 )
440: *     ..
441: *     .. Local Scalars ..
442:       LOGICAL            COLEQU, EQUIL, NOFACT, NOTRAN, ROWEQU
443:       INTEGER            INFEQU, J
444:       REAL               AMAX, BIGNUM, COLCND, RCMAX, RCMIN, ROWCND,
445:      $                   SMLNUM
446: *     ..
447: *     .. External Functions ..
448:       EXTERNAL           LSAME, SLAMCH, SLA_RPVGRW
449:       LOGICAL            LSAME
450:       REAL               SLAMCH, SLA_RPVGRW
451: *     ..
452: *     .. External Subroutines ..
453:       EXTERNAL           SGEEQUB, SGETRF, SGETRS, SLACPY, SLAQGE,
454:      $                   XERBLA, SLASCL2, SGERFSX
455: *     ..
456: *     .. Intrinsic Functions ..
457:       INTRINSIC          MAX, MIN
458: *     ..
459: *     .. Executable Statements ..
460: *
461:       INFO = 0
462:       NOFACT = LSAME( FACT, 'N' )
463:       EQUIL = LSAME( FACT, 'E' )
464:       NOTRAN = LSAME( TRANS, 'N' )
465:       SMLNUM = SLAMCH( 'Safe minimum' )
466:       BIGNUM = ONE / SMLNUM
467:       IF( NOFACT .OR. EQUIL ) THEN
468:          EQUED = 'N'
469:          ROWEQU = .FALSE.
470:          COLEQU = .FALSE.
471:       ELSE
472:          ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
473:          COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
474:       END IF
475: *
476: *     Default is failure.  If an input parameter is wrong or
477: *     factorization fails, make everything look horrible.  Only the
478: *     pivot growth is set here, the rest is initialized in SGERFSX.
479: *
480:       RPVGRW = ZERO
481: *
482: *     Test the input parameters.  PARAMS is not tested until SGERFSX.
483: *
484:       IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.
485:      $     LSAME( FACT, 'F' ) ) THEN
486:          INFO = -1
487:       ELSE IF( .NOT.NOTRAN .AND. .NOT.LSAME( TRANS, 'T' ) .AND. .NOT.
488:      $        LSAME( TRANS, 'C' ) ) THEN
489:          INFO = -2
490:       ELSE IF( N.LT.0 ) THEN
491:          INFO = -3
492:       ELSE IF( NRHS.LT.0 ) THEN
493:          INFO = -4
494:       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
495:          INFO = -6
496:       ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN
497:          INFO = -8
498:       ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
499:      $        ( ROWEQU .OR. COLEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
500:          INFO = -10
501:       ELSE
502:          IF( ROWEQU ) THEN
503:             RCMIN = BIGNUM
504:             RCMAX = ZERO
505:             DO 10 J = 1, N
506:                RCMIN = MIN( RCMIN, R( J ) )
507:                RCMAX = MAX( RCMAX, R( J ) )
508:  10         CONTINUE
509:             IF( RCMIN.LE.ZERO ) THEN
510:                INFO = -11
511:             ELSE IF( N.GT.0 ) THEN
512:                ROWCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
513:             ELSE
514:                ROWCND = ONE
515:             END IF
516:          END IF
517:          IF( COLEQU .AND. INFO.EQ.0 ) THEN
518:             RCMIN = BIGNUM
519:             RCMAX = ZERO
520:             DO 20 J = 1, N
521:                RCMIN = MIN( RCMIN, C( J ) )
522:                RCMAX = MAX( RCMAX, C( J ) )
523:  20         CONTINUE
524:             IF( RCMIN.LE.ZERO ) THEN
525:                INFO = -12
526:             ELSE IF( N.GT.0 ) THEN
527:                COLCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
528:             ELSE
529:                COLCND = ONE
530:             END IF
531:          END IF
532:          IF( INFO.EQ.0 ) THEN
533:             IF( LDB.LT.MAX( 1, N ) ) THEN
534:                INFO = -14
535:             ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
536:                INFO = -16
537:             END IF
538:          END IF
539:       END IF
540: *
541:       IF( INFO.NE.0 ) THEN
542:          CALL XERBLA( 'SGESVXX', -INFO )
543:          RETURN
544:       END IF
545: *
546:       IF( EQUIL ) THEN
547: *
548: *     Compute row and column scalings to equilibrate the matrix A.
549: *
550:          CALL SGEEQUB( N, N, A, LDA, R, C, ROWCND, COLCND, AMAX,
551:      $        INFEQU )
552:          IF( INFEQU.EQ.0 ) THEN
553: *
554: *     Equilibrate the matrix.
555: *
556:             CALL SLAQGE( N, N, A, LDA, R, C, ROWCND, COLCND, AMAX,
557:      $           EQUED )
558:             ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
559:             COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
560:          END IF
561: *
562: *     If the scaling factors are not applied, set them to 1.0.
563: *
564:          IF ( .NOT.ROWEQU ) THEN
565:             DO J = 1, N
566:                R( J ) = 1.0
567:             END DO
568:          END IF
569:          IF ( .NOT.COLEQU ) THEN
570:             DO J = 1, N
571:                C( J ) = 1.0
572:             END DO
573:          END IF
574:       END IF
575: *
576: *     Scale the right-hand side.
577: *
578:       IF( NOTRAN ) THEN
579:          IF( ROWEQU ) CALL SLASCL2( N, NRHS, R, B, LDB )
580:       ELSE
581:          IF( COLEQU ) CALL SLASCL2( N, NRHS, C, B, LDB )
582:       END IF
583: *
584:       IF( NOFACT .OR. EQUIL ) THEN
585: *
586: *        Compute the LU factorization of A.
587: *
588:          CALL SLACPY( 'Full', N, N, A, LDA, AF, LDAF )
589:          CALL SGETRF( N, N, AF, LDAF, IPIV, INFO )
590: *
591: *        Return if INFO is non-zero.
592: *
593:          IF( INFO.GT.0 ) THEN
594: *
595: *           Pivot in column INFO is exactly 0
596: *           Compute the reciprocal pivot growth factor of the
597: *           leading rank-deficient INFO columns of A.
598: *
599:             RPVGRW = SLA_RPVGRW( N, INFO, A, LDA, AF, LDAF )
600:             RETURN
601:          END IF
602:       END IF
603: *
604: *     Compute the reciprocal pivot growth factor RPVGRW.
605: *
606:       RPVGRW = SLA_RPVGRW( N, N, A, LDA, AF, LDAF )
607: *
608: *     Compute the solution matrix X.
609: *
610:       CALL SLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
611:       CALL SGETRS( TRANS, N, NRHS, AF, LDAF, IPIV, X, LDX, INFO )
612: *
613: *     Use iterative refinement to improve the computed solution and
614: *     compute error bounds and backward error estimates for it.
615: *
616:       CALL SGERFSX( TRANS, EQUED, N, NRHS, A, LDA, AF, LDAF,
617:      $     IPIV, R, C, B, LDB, X, LDX, RCOND, BERR,
618:      $     N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS,
619:      $     WORK, IWORK, INFO )
620: *
621: *     Scale solutions.
622: *
623:       IF ( COLEQU .AND. NOTRAN ) THEN
624:          CALL SLASCL2 ( N, NRHS, C, X, LDX )
625:       ELSE IF ( ROWEQU .AND. .NOT.NOTRAN ) THEN
626:          CALL SLASCL2 ( N, NRHS, R, X, LDX )
627:       END IF
628: *
629:       RETURN
630: *
631: *     End of SGESVXX
632: 
633:       END
634: