1       SUBROUTINE ZGBSVXX( FACT, TRANS, N, KL, KU, NRHS, AB, LDAB, AFB,
  2      $                    LDAFB, IPIV, EQUED, R, C, B, LDB, X, LDX,
  3      $                    RCOND, RPVGRW, BERR, N_ERR_BNDS,
  4      $                    ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS,
  5      $                    WORK, RWORK, 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, LDAB, LDAFB, LDB, LDX, N, NRHS, NPARAMS,
 20      $                   N_ERR_BNDS
 21       DOUBLE PRECISION   RCOND, RPVGRW
 22 *     ..
 23 *     .. Array Arguments ..
 24       INTEGER            IPIV( * )
 25       COMPLEX*16         AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ),
 26      $                   X( LDX , * ),WORK( * )
 27       DOUBLE PRECISION   R( * ), C( * ), PARAMS( * ), BERR( * ),
 28      $                   ERR_BNDS_NORM( NRHS, * ),
 29      $                   ERR_BNDS_COMP( NRHS, * ), RWORK( * )
 30 *     ..
 31 *
 32 *     Purpose
 33 *     =======
 34 *
 35 *     ZGBSVXX uses the LU factorization to compute the solution to a
 36 *     complex*16 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. ZGBSVXX 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 *     ZGBSVXX 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 *     ZGBSVXX 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 ZGBSVXX 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 *     KL      (input) INTEGER
130 *     The number of subdiagonals within the band of A.  KL >= 0.
131 *
132 *     KU      (input) INTEGER
133 *     The number of superdiagonals within the band of A.  KU >= 0.
134 *
135 *     NRHS    (input) INTEGER
136 *     The number of right hand sides, i.e., the number of columns
137 *     of the matrices B and X.  NRHS >= 0.
138 *
139 *     AB      (input/output) COMPLEX*16 array, dimension (LDAB,N)
140 *     On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
141 *     The j-th column of A is stored in the j-th column of the
142 *     array AB as follows:
143 *     AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)
144 *
145 *     If FACT = 'F' and EQUED is not 'N', then AB must have been
146 *     equilibrated by the scaling factors in R and/or C.  AB is not
147 *     modified if FACT = 'F' or 'N', or if FACT = 'E' and
148 *     EQUED = 'N' on exit.
149 *
150 *     On exit, if EQUED .ne. 'N', A is scaled as follows:
151 *     EQUED = 'R':  A := diag(R) * A
152 *     EQUED = 'C':  A := A * diag(C)
153 *     EQUED = 'B':  A := diag(R) * A * diag(C).
154 *
155 *     LDAB    (input) INTEGER
156 *     The leading dimension of the array AB.  LDAB >= KL+KU+1.
157 *
158 *     AFB     (input or output) COMPLEX*16 array, dimension (LDAFB,N)
159 *     If FACT = 'F', then AFB is an input argument and on entry
160 *     contains details of the LU factorization of the band matrix
161 *     A, as computed by ZGBTRF.  U is stored as an upper triangular
162 *     band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
163 *     and the multipliers used during the factorization are stored
164 *     in rows KL+KU+2 to 2*KL+KU+1.  If EQUED .ne. 'N', then AFB is
165 *     the factored form of the equilibrated matrix A.
166 *
167 *     If FACT = 'N', then AF is an output argument and on exit
168 *     returns the factors L and U from the factorization A = P*L*U
169 *     of the original matrix A.
170 *
171 *     If FACT = 'E', then AF is an output argument and on exit
172 *     returns the factors L and U from the factorization A = P*L*U
173 *     of the equilibrated matrix A (see the description of A for
174 *     the form of the equilibrated matrix).
175 *
176 *     LDAFB   (input) INTEGER
177 *     The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1.
178 *
179 *     IPIV    (input or output) INTEGER array, dimension (N)
180 *     If FACT = 'F', then IPIV is an input argument and on entry
181 *     contains the pivot indices from the factorization A = P*L*U
182 *     as computed by DGETRF; row i of the matrix was interchanged
183 *     with row IPIV(i).
184 *
185 *     If FACT = 'N', then IPIV is an output argument and on exit
186 *     contains the pivot indices from the factorization A = P*L*U
187 *     of the original matrix A.
188 *
189 *     If FACT = 'E', then IPIV is an output argument and on exit
190 *     contains the pivot indices from the factorization A = P*L*U
191 *     of the equilibrated matrix A.
192 *
193 *     EQUED   (input or output) CHARACTER*1
194 *     Specifies the form of equilibration that was done.
195 *       = 'N':  No equilibration (always true if FACT = 'N').
196 *       = 'R':  Row equilibration, i.e., A has been premultiplied by
197 *               diag(R).
198 *       = 'C':  Column equilibration, i.e., A has been postmultiplied
199 *               by diag(C).
200 *       = 'B':  Both row and column equilibration, i.e., A has been
201 *               replaced by diag(R) * A * diag(C).
202 *     EQUED is an input argument if FACT = 'F'; otherwise, it is an
203 *     output argument.
204 *
205 *     R       (input or output) DOUBLE PRECISION array, dimension (N)
206 *     The row scale factors for A.  If EQUED = 'R' or 'B', A is
207 *     multiplied on the left by diag(R); if EQUED = 'N' or 'C', R
208 *     is not accessed.  R is an input argument if FACT = 'F';
209 *     otherwise, R is an output argument.  If FACT = 'F' and
210 *     EQUED = 'R' or 'B', each element of R must be positive.
211 *     If R is output, each element of R is a power of the radix.
212 *     If R is input, each element of R should be a power of the radix
213 *     to ensure a reliable solution and error estimates. Scaling by
214 *     powers of the radix does not cause rounding errors unless the
215 *     result underflows or overflows. Rounding errors during scaling
216 *     lead to refining with a matrix that is not equivalent to the
217 *     input matrix, producing error estimates that may not be
218 *     reliable.
219 *
220 *     C       (input or output) DOUBLE PRECISION array, dimension (N)
221 *     The column scale factors for A.  If EQUED = 'C' or 'B', A is
222 *     multiplied on the right by diag(C); if EQUED = 'N' or 'R', C
223 *     is not accessed.  C is an input argument if FACT = 'F';
224 *     otherwise, C is an output argument.  If FACT = 'F' and
225 *     EQUED = 'C' or 'B', each element of C must be positive.
226 *     If C is output, each element of C is a power of the radix.
227 *     If C is input, each element of C should be a power of the radix
228 *     to ensure a reliable solution and error estimates. Scaling by
229 *     powers of the radix does not cause rounding errors unless the
230 *     result underflows or overflows. Rounding errors during scaling
231 *     lead to refining with a matrix that is not equivalent to the
232 *     input matrix, producing error estimates that may not be
233 *     reliable.
234 *
235 *     B       (input/output) COMPLEX*16 array, dimension (LDB,NRHS)
236 *     On entry, the N-by-NRHS right hand side matrix B.
237 *     On exit,
238 *     if EQUED = 'N', B is not modified;
239 *     if TRANS = 'N' and EQUED = 'R' or 'B', B is overwritten by
240 *        diag(R)*B;
241 *     if TRANS = 'T' or 'C' and EQUED = 'C' or 'B', B is
242 *        overwritten by diag(C)*B.
243 *
244 *     LDB     (input) INTEGER
245 *     The leading dimension of the array B.  LDB >= max(1,N).
246 *
247 *     X       (output) COMPLEX*16 array, dimension (LDX,NRHS)
248 *     If INFO = 0, the N-by-NRHS solution matrix X to the original
249 *     system of equations.  Note that A and B are modified on exit
250 *     if EQUED .ne. 'N', and the solution to the equilibrated system is
251 *     inv(diag(C))*X if TRANS = 'N' and EQUED = 'C' or 'B', or
252 *     inv(diag(R))*X if TRANS = 'T' or 'C' and EQUED = 'R' or 'B'.
253 *
254 *     LDX     (input) INTEGER
255 *     The leading dimension of the array X.  LDX >= max(1,N).
256 *
257 *     RCOND   (output) DOUBLE PRECISION
258 *     Reciprocal scaled condition number.  This is an estimate of the
259 *     reciprocal Skeel condition number of the matrix A after
260 *     equilibration (if done).  If this is less than the machine
261 *     precision (in particular, if it is zero), the matrix is singular
262 *     to working precision.  Note that the error may still be small even
263 *     if this number is very small and the matrix appears ill-
264 *     conditioned.
265 *
266 *     RPVGRW  (output) DOUBLE PRECISION
267 *     Reciprocal pivot growth.  On exit, this contains the reciprocal
268 *     pivot growth factor norm(A)/norm(U). The "max absolute element"
269 *     norm is used.  If this is much less than 1, then the stability of
270 *     the LU factorization of the (equilibrated) matrix A could be poor.
271 *     This also means that the solution X, estimated condition numbers,
272 *     and error bounds could be unreliable. If factorization fails with
273 *     0<INFO<=N, then this contains the reciprocal pivot growth factor
274 *     for the leading INFO columns of A.  In DGESVX, this quantity is
275 *     returned in WORK(1).
276 *
277 *     BERR    (output) DOUBLE PRECISION array, dimension (NRHS)
278 *     Componentwise relative backward error.  This is the
279 *     componentwise relative backward error of each solution vector X(j)
280 *     (i.e., the smallest relative change in any element of A or B that
281 *     makes X(j) an exact solution).
282 *
283 *     N_ERR_BNDS (input) INTEGER
284 *     Number of error bounds to return for each right hand side
285 *     and each type (normwise or componentwise).  See ERR_BNDS_NORM and
286 *     ERR_BNDS_COMP below.
287 *
288 *     ERR_BNDS_NORM  (output) DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
289 *     For each right-hand side, this array contains information about
290 *     various error bounds and condition numbers corresponding to the
291 *     normwise relative error, which is defined as follows:
292 *
293 *     Normwise relative error in the ith solution vector:
294 *             max_j (abs(XTRUE(j,i) - X(j,i)))
295 *            ------------------------------
296 *                  max_j abs(X(j,i))
297 *
298 *     The array is indexed by the type of error information as described
299 *     below. There currently are up to three pieces of information
300 *     returned.
301 *
302 *     The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
303 *     right-hand side.
304 *
305 *     The second index in ERR_BNDS_NORM(:,err) contains the following
306 *     three fields:
307 *     err = 1 "Trust/don't trust" boolean. Trust the answer if the
308 *              reciprocal condition number is less than the threshold
309 *              sqrt(n) * dlamch('Epsilon').
310 *
311 *     err = 2 "Guaranteed" error bound: The estimated forward error,
312 *              almost certainly within a factor of 10 of the true error
313 *              so long as the next entry is greater than the threshold
314 *              sqrt(n) * dlamch('Epsilon'). This error bound should only
315 *              be trusted if the previous boolean is true.
316 *
317 *     err = 3  Reciprocal condition number: Estimated normwise
318 *              reciprocal condition number.  Compared with the threshold
319 *              sqrt(n) * dlamch('Epsilon') to determine if the error
320 *              estimate is "guaranteed". These reciprocal condition
321 *              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
322 *              appropriately scaled matrix Z.
323 *              Let Z = S*A, where S scales each row by a power of the
324 *              radix so all absolute row sums of Z are approximately 1.
325 *
326 *     See Lapack Working Note 165 for further details and extra
327 *     cautions.
328 *
329 *     ERR_BNDS_COMP  (output) DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
330 *     For each right-hand side, this array contains information about
331 *     various error bounds and condition numbers corresponding to the
332 *     componentwise relative error, which is defined as follows:
333 *
334 *     Componentwise relative error in the ith solution vector:
335 *                    abs(XTRUE(j,i) - X(j,i))
336 *             max_j ----------------------
337 *                         abs(X(j,i))
338 *
339 *     The array is indexed by the right-hand side i (on which the
340 *     componentwise relative error depends), and the type of error
341 *     information as described below. There currently are up to three
342 *     pieces of information returned for each right-hand side. If
343 *     componentwise accuracy is not requested (PARAMS(3) = 0.0), then
344 *     ERR_BNDS_COMP is not accessed.  If N_ERR_BNDS .LT. 3, then at most
345 *     the first (:,N_ERR_BNDS) entries are returned.
346 *
347 *     The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
348 *     right-hand side.
349 *
350 *     The second index in ERR_BNDS_COMP(:,err) contains the following
351 *     three fields:
352 *     err = 1 "Trust/don't trust" boolean. Trust the answer if the
353 *              reciprocal condition number is less than the threshold
354 *              sqrt(n) * dlamch('Epsilon').
355 *
356 *     err = 2 "Guaranteed" error bound: The estimated forward error,
357 *              almost certainly within a factor of 10 of the true error
358 *              so long as the next entry is greater than the threshold
359 *              sqrt(n) * dlamch('Epsilon'). This error bound should only
360 *              be trusted if the previous boolean is true.
361 *
362 *     err = 3  Reciprocal condition number: Estimated componentwise
363 *              reciprocal condition number.  Compared with the threshold
364 *              sqrt(n) * dlamch('Epsilon') to determine if the error
365 *              estimate is "guaranteed". These reciprocal condition
366 *              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
367 *              appropriately scaled matrix Z.
368 *              Let Z = S*(A*diag(x)), where x is the solution for the
369 *              current right-hand side and S scales each row of
370 *              A*diag(x) by a power of the radix so all absolute row
371 *              sums of Z are approximately 1.
372 *
373 *     See Lapack Working Note 165 for further details and extra
374 *     cautions.
375 *
376 *     NPARAMS (input) INTEGER
377 *     Specifies the number of parameters set in PARAMS.  If .LE. 0, the
378 *     PARAMS array is never referenced and default values are used.
379 *
380 *     PARAMS  (input / output) DOUBLE PRECISION array, dimension NPARAMS
381 *     Specifies algorithm parameters.  If an entry is .LT. 0.0, then
382 *     that entry will be filled with default value used for that
383 *     parameter.  Only positions up to NPARAMS are accessed; defaults
384 *     are used for higher-numbered parameters.
385 *
386 *       PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
387 *            refinement or not.
388 *         Default: 1.0D+0
389 *            = 0.0 : No refinement is performed, and no error bounds are
390 *                    computed.
391 *            = 1.0 : Use the extra-precise refinement algorithm.
392 *              (other values are reserved for future use)
393 *
394 *       PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
395 *            computations allowed for refinement.
396 *         Default: 10
397 *         Aggressive: Set to 100 to permit convergence using approximate
398 *                     factorizations or factorizations other than LU. If
399 *                     the factorization uses a technique other than
400 *                     Gaussian elimination, the guarantees in
401 *                     err_bnds_norm and err_bnds_comp may no longer be
402 *                     trustworthy.
403 *
404 *       PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
405 *            will attempt to find a solution with small componentwise
406 *            relative error in the double-precision algorithm.  Positive
407 *            is true, 0.0 is false.
408 *         Default: 1.0 (attempt componentwise convergence)
409 *
410 *     WORK    (workspace) COMPLEX*16 array, dimension (2*N)
411 *
412 *     RWORK   (workspace) DOUBLE PRECISION array, dimension (2*N)
413 *
414 *     INFO    (output) INTEGER
415 *       = 0:  Successful exit. The solution to every right-hand side is
416 *         guaranteed.
417 *       < 0:  If INFO = -i, the i-th argument had an illegal value
418 *       > 0 and <= N:  U(INFO,INFO) is exactly zero.  The factorization
419 *         has been completed, but the factor U is exactly singular, so
420 *         the solution and error bounds could not be computed. RCOND = 0
421 *         is returned.
422 *       = N+J: The solution corresponding to the Jth right-hand side is
423 *         not guaranteed. The solutions corresponding to other right-
424 *         hand sides K with K > J may not be guaranteed as well, but
425 *         only the first such right-hand side is reported. If a small
426 *         componentwise error is not requested (PARAMS(3) = 0.0) then
427 *         the Jth right-hand side is the first with a normwise error
428 *         bound that is not guaranteed (the smallest J such
429 *         that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
430 *         the Jth right-hand side is the first with either a normwise or
431 *         componentwise error bound that is not guaranteed (the smallest
432 *         J such that either ERR_BNDS_NORM(J,1) = 0.0 or
433 *         ERR_BNDS_COMP(J,1) = 0.0). See the definition of
434 *         ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
435 *         about all of the right-hand sides check ERR_BNDS_NORM or
436 *         ERR_BNDS_COMP.
437 *
438 *     ==================================================================
439 *
440 *     .. Parameters ..
441       DOUBLE PRECISION   ZERO, ONE
442       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
443       INTEGER            FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
444       INTEGER            RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
445       INTEGER            CMP_ERR_I, PIV_GROWTH_I
446       PARAMETER          ( FINAL_NRM_ERR_I = 1, FINAL_CMP_ERR_I = 2,
447      $                   BERR_I = 3 )
448       PARAMETER          ( RCOND_I = 4, NRM_RCOND_I = 5, NRM_ERR_I = 6 )
449       PARAMETER          ( CMP_RCOND_I = 7, CMP_ERR_I = 8,
450      $                   PIV_GROWTH_I = 9 )
451 *     ..
452 *     .. Local Scalars ..
453       LOGICAL            COLEQU, EQUIL, NOFACT, NOTRAN, ROWEQU
454       INTEGER            INFEQU, I, J, KL, KU
455       DOUBLE PRECISION   AMAX, BIGNUM, COLCND, RCMAX, RCMIN,
456      $                   ROWCND, SMLNUM
457 *     ..
458 *     .. External Functions ..
459       EXTERNAL           LSAME, DLAMCH, ZLA_GBRPVGRW
460       LOGICAL            LSAME
461       DOUBLE PRECISION   DLAMCH, ZLA_GBRPVGRW
462 *     ..
463 *     .. External Subroutines ..
464       EXTERNAL           ZGBEQUB, ZGBTRF, ZGBTRS, ZLACPY, ZLAQGB,
465      $                   XERBLA, ZLASCL2, ZGBRFSX
466 *     ..
467 *     .. Intrinsic Functions ..
468       INTRINSIC          MAXMIN
469 *     ..
470 *     .. Executable Statements ..
471 *
472       INFO = 0
473       NOFACT = LSAME( FACT, 'N' )
474       EQUIL = LSAME( FACT, 'E' )
475       NOTRAN = LSAME( TRANS, 'N' )
476       SMLNUM = DLAMCH( 'Safe minimum' )
477       BIGNUM = ONE / SMLNUM
478       IF( NOFACT .OR. EQUIL ) THEN
479          EQUED = 'N'
480          ROWEQU = .FALSE.
481          COLEQU = .FALSE.
482       ELSE
483          ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
484          COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
485       END IF
486 *
487 *     Default is failure.  If an input parameter is wrong or
488 *     factorization fails, make everything look horrible.  Only the
489 *     pivot growth is set here, the rest is initialized in ZGBRFSX.
490 *
491       RPVGRW = ZERO
492 *
493 *     Test the input parameters.  PARAMS is not tested until DGERFSX.
494 *
495       IF.NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.
496      $     LSAME( FACT, 'F' ) ) THEN
497          INFO = -1
498       ELSE IF.NOT.NOTRAN .AND. .NOT.LSAME( TRANS, 'T' ) .AND. .NOT.
499      $        LSAME( TRANS, 'C' ) ) THEN
500          INFO = -2
501       ELSE IF( N.LT.0 ) THEN
502          INFO = -3
503       ELSE IF( KL.LT.0 ) THEN
504          INFO = -4
505       ELSE IF( KU.LT.0 ) THEN
506          INFO = -5
507       ELSE IF( NRHS.LT.0 ) THEN
508          INFO = -6
509       ELSE IF( LDAB.LT.KL+KU+1 ) THEN
510          INFO = -8
511       ELSE IF( LDAFB.LT.2*KL+KU+1 ) THEN
512          INFO = -10
513       ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
514      $        ( ROWEQU .OR. COLEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
515          INFO = -12
516       ELSE
517          IF( ROWEQU ) THEN
518             RCMIN = BIGNUM
519             RCMAX = ZERO
520             DO 10 J = 1, N
521                RCMIN = MIN( RCMIN, R( J ) )
522                RCMAX = MAX( RCMAX, R( J ) )
523  10         CONTINUE
524             IF( RCMIN.LE.ZERO ) THEN
525                INFO = -13
526             ELSE IF( N.GT.0 ) THEN
527                ROWCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
528             ELSE
529                ROWCND = ONE
530             END IF
531          END IF
532          IF( COLEQU .AND. INFO.EQ.0 ) THEN
533             RCMIN = BIGNUM
534             RCMAX = ZERO
535             DO 20 J = 1, N
536                RCMIN = MIN( RCMIN, C( J ) )
537                RCMAX = MAX( RCMAX, C( J ) )
538  20         CONTINUE
539             IF( RCMIN.LE.ZERO ) THEN
540                INFO = -14
541             ELSE IF( N.GT.0 ) THEN
542                COLCND = MAX( RCMIN, SMLNUM ) / MIN( RCMAX, BIGNUM )
543             ELSE
544                COLCND = ONE
545             END IF
546          END IF
547          IF( INFO.EQ.0 ) THEN
548             IF( LDB.LT.MAX1, N ) ) THEN
549                INFO = -15
550             ELSE IF( LDX.LT.MAX1, N ) ) THEN
551                INFO = -16
552             END IF
553          END IF
554       END IF
555 *
556       IF( INFO.NE.0 ) THEN
557          CALL XERBLA( 'ZGBSVXX'-INFO )
558          RETURN
559       END IF
560 *
561       IF( EQUIL ) THEN
562 *
563 *     Compute row and column scalings to equilibrate the matrix A.
564 *
565          CALL ZGBEQUB( N, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND,
566      $        AMAX, INFEQU )
567          IF( INFEQU.EQ.0 ) THEN
568 *
569 *     Equilibrate the matrix.
570 *
571             CALL ZLAQGB( N, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND,
572      $           AMAX, EQUED )
573             ROWEQU = LSAME( EQUED, 'R' ) .OR. LSAME( EQUED, 'B' )
574             COLEQU = LSAME( EQUED, 'C' ) .OR. LSAME( EQUED, 'B' )
575          END IF
576 *
577 *     If the scaling factors are not applied, set them to 1.0.
578 *
579          IF ( .NOT.ROWEQU ) THEN
580             DO J = 1, N
581                R( J ) = 1.0D+0
582             END DO
583          END IF
584          IF ( .NOT.COLEQU ) THEN
585             DO J = 1, N
586                C( J ) = 1.0D+0
587             END DO
588          END IF
589       END IF
590 *
591 *     Scale the right-hand side.
592 *
593       IF( NOTRAN ) THEN
594          IF( ROWEQU ) CALL ZLASCL2( N, NRHS, R, B, LDB )
595       ELSE
596          IF( COLEQU ) CALL ZLASCL2( N, NRHS, C, B, LDB )
597       END IF
598 *
599       IF( NOFACT .OR. EQUIL ) THEN
600 *
601 *        Compute the LU factorization of A.
602 *
603          DO 40, J = 1, N
604             DO 30, I = KL+12*KL+KU+1
605                AFB( I, J ) = AB( I-KL, J )
606  30         CONTINUE
607  40      CONTINUE
608          CALL ZGBTRF( N, N, KL, KU, AFB, LDAFB, IPIV, INFO )
609 *
610 *        Return if INFO is non-zero.
611 *
612          IF( INFO.GT.0 ) THEN
613 *
614 *           Pivot in column INFO is exactly 0
615 *           Compute the reciprocal pivot growth factor of the
616 *           leading rank-deficient INFO columns of A.
617 *
618             RPVGRW = ZLA_GBRPVGRW( N, KL, KU, INFO, AB, LDAB, AFB,
619      $           LDAFB )
620             RETURN
621          END IF
622       END IF
623 *
624 *     Compute the reciprocal pivot growth factor RPVGRW.
625 *
626       RPVGRW = ZLA_GBRPVGRW( N, KL, KU, N, AB, LDAB, AFB, LDAFB )
627 *
628 *     Compute the solution matrix X.
629 *
630       CALL ZLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
631       CALL ZGBTRS( TRANS, N, KL, KU, NRHS, AFB, LDAFB, IPIV, X, LDX,
632      $     INFO )
633 *
634 *     Use iterative refinement to improve the computed solution and
635 *     compute error bounds and backward error estimates for it.
636 *
637       CALL ZGBRFSX( TRANS, EQUED, N, KL, KU, NRHS, AB, LDAB, AFB, LDAFB,
638      $     IPIV, R, C, B, LDB, X, LDX, RCOND, BERR,
639      $     N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS,
640      $     WORK, RWORK, INFO )
641 
642 *
643 *     Scale solutions.
644 *
645       IF ( COLEQU .AND. NOTRAN ) THEN
646          CALL ZLASCL2( N, NRHS, C, X, LDX )
647       ELSE IF ( ROWEQU .AND. .NOT.NOTRAN ) THEN
648          CALL ZLASCL2( N, NRHS, R, X, LDX )
649       END IF
650 *
651       RETURN
652 *
653 *     End of ZGBSVXX
654 *
655       END