LAPACK  3.10.1
LAPACK: Linear Algebra PACKage
dgbsvx.f
Go to the documentation of this file.
1 *> \brief <b> DGBSVX computes the solution to system of linear equations A * X = B for GB matrices</b>
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download DGBSVX + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgbsvx.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgbsvx.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgbsvx.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE DGBSVX( FACT, TRANS, N, KL, KU, NRHS, AB, LDAB, AFB,
22 * LDAFB, IPIV, EQUED, R, C, B, LDB, X, LDX,
23 * RCOND, FERR, BERR, WORK, IWORK, INFO )
24 *
25 * .. Scalar Arguments ..
26 * CHARACTER EQUED, FACT, TRANS
27 * INTEGER INFO, KL, KU, LDAB, LDAFB, LDB, LDX, N, NRHS
28 * DOUBLE PRECISION RCOND
29 * ..
30 * .. Array Arguments ..
31 * INTEGER IPIV( * ), IWORK( * )
32 * DOUBLE PRECISION AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ),
33 * $ BERR( * ), C( * ), FERR( * ), R( * ),
34 * $ WORK( * ), X( LDX, * )
35 * ..
36 *
37 *
38 *> \par Purpose:
39 * =============
40 *>
41 *> \verbatim
42 *>
43 *> DGBSVX uses the LU factorization to compute the solution to a real
44 *> system of linear equations A * X = B, A**T * X = B, or A**H * X = B,
45 *> where A is a band matrix of order N with KL subdiagonals and KU
46 *> superdiagonals, and X and B are N-by-NRHS matrices.
47 *>
48 *> Error bounds on the solution and a condition estimate are also
49 *> provided.
50 *> \endverbatim
51 *
52 *> \par Description:
53 * =================
54 *>
55 *> \verbatim
56 *>
57 *> The following steps are performed by this subroutine:
58 *>
59 *> 1. If FACT = 'E', real scaling factors are computed to equilibrate
60 *> the system:
61 *> TRANS = 'N': diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B
62 *> TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
63 *> TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B
64 *> Whether or not the system will be equilibrated depends on the
65 *> scaling of the matrix A, but if equilibration is used, A is
66 *> overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N')
67 *> or diag(C)*B (if TRANS = 'T' or 'C').
68 *>
69 *> 2. If FACT = 'N' or 'E', the LU decomposition is used to factor the
70 *> matrix A (after equilibration if FACT = 'E') as
71 *> A = L * U,
72 *> where L is a product of permutation and unit lower triangular
73 *> matrices with KL subdiagonals, and U is upper triangular with
74 *> KL+KU superdiagonals.
75 *>
76 *> 3. If some U(i,i)=0, so that U is exactly singular, then the routine
77 *> returns with INFO = i. Otherwise, the factored form of A is used
78 *> to estimate the condition number of the matrix A. If the
79 *> reciprocal of the condition number is less than machine precision,
80 *> INFO = N+1 is returned as a warning, but the routine still goes on
81 *> to solve for X and compute error bounds as described below.
82 *>
83 *> 4. The system of equations is solved for X using the factored form
84 *> of A.
85 *>
86 *> 5. Iterative refinement is applied to improve the computed solution
87 *> matrix and calculate error bounds and backward error estimates
88 *> for it.
89 *>
90 *> 6. If equilibration was used, the matrix X is premultiplied by
91 *> diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so
92 *> that it solves the original system before equilibration.
93 *> \endverbatim
94 *
95 * Arguments:
96 * ==========
97 *
98 *> \param[in] FACT
99 *> \verbatim
100 *> FACT is CHARACTER*1
101 *> Specifies whether or not the factored form of the matrix A is
102 *> supplied on entry, and if not, whether the matrix A should be
103 *> equilibrated before it is factored.
104 *> = 'F': On entry, AFB and IPIV contain the factored form of
105 *> A. If EQUED is not 'N', the matrix A has been
106 *> equilibrated with scaling factors given by R and C.
107 *> AB, AFB, and IPIV are not modified.
108 *> = 'N': The matrix A will be copied to AFB and factored.
109 *> = 'E': The matrix A will be equilibrated if necessary, then
110 *> copied to AFB and factored.
111 *> \endverbatim
112 *>
113 *> \param[in] TRANS
114 *> \verbatim
115 *> TRANS is CHARACTER*1
116 *> Specifies the form of the system of equations.
117 *> = 'N': A * X = B (No transpose)
118 *> = 'T': A**T * X = B (Transpose)
119 *> = 'C': A**H * X = B (Transpose)
120 *> \endverbatim
121 *>
122 *> \param[in] N
123 *> \verbatim
124 *> N is INTEGER
125 *> The number of linear equations, i.e., the order of the
126 *> matrix A. N >= 0.
127 *> \endverbatim
128 *>
129 *> \param[in] KL
130 *> \verbatim
131 *> KL is INTEGER
132 *> The number of subdiagonals within the band of A. KL >= 0.
133 *> \endverbatim
134 *>
135 *> \param[in] KU
136 *> \verbatim
137 *> KU is INTEGER
138 *> The number of superdiagonals within the band of A. KU >= 0.
139 *> \endverbatim
140 *>
141 *> \param[in] NRHS
142 *> \verbatim
143 *> NRHS is INTEGER
144 *> The number of right hand sides, i.e., the number of columns
145 *> of the matrices B and X. NRHS >= 0.
146 *> \endverbatim
147 *>
148 *> \param[in,out] AB
149 *> \verbatim
150 *> AB is DOUBLE PRECISION array, dimension (LDAB,N)
151 *> On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
152 *> The j-th column of A is stored in the j-th column of the
153 *> array AB as follows:
154 *> AB(KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+kl)
155 *>
156 *> If FACT = 'F' and EQUED is not 'N', then A must have been
157 *> equilibrated by the scaling factors in R and/or C. AB is not
158 *> modified if FACT = 'F' or 'N', or if FACT = 'E' and
159 *> EQUED = 'N' on exit.
160 *>
161 *> On exit, if EQUED .ne. 'N', A is scaled as follows:
162 *> EQUED = 'R': A := diag(R) * A
163 *> EQUED = 'C': A := A * diag(C)
164 *> EQUED = 'B': A := diag(R) * A * diag(C).
165 *> \endverbatim
166 *>
167 *> \param[in] LDAB
168 *> \verbatim
169 *> LDAB is INTEGER
170 *> The leading dimension of the array AB. LDAB >= KL+KU+1.
171 *> \endverbatim
172 *>
173 *> \param[in,out] AFB
174 *> \verbatim
175 *> AFB is DOUBLE PRECISION array, dimension (LDAFB,N)
176 *> If FACT = 'F', then AFB is an input argument and on entry
177 *> contains details of the LU factorization of the band matrix
178 *> A, as computed by DGBTRF. U is stored as an upper triangular
179 *> band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1,
180 *> and the multipliers used during the factorization are stored
181 *> in rows KL+KU+2 to 2*KL+KU+1. If EQUED .ne. 'N', then AFB is
182 *> the factored form of the equilibrated matrix A.
183 *>
184 *> If FACT = 'N', then AFB is an output argument and on exit
185 *> returns details of the LU factorization of A.
186 *>
187 *> If FACT = 'E', then AFB is an output argument and on exit
188 *> returns details of the LU factorization of the equilibrated
189 *> matrix A (see the description of AB for the form of the
190 *> equilibrated matrix).
191 *> \endverbatim
192 *>
193 *> \param[in] LDAFB
194 *> \verbatim
195 *> LDAFB is INTEGER
196 *> The leading dimension of the array AFB. LDAFB >= 2*KL+KU+1.
197 *> \endverbatim
198 *>
199 *> \param[in,out] IPIV
200 *> \verbatim
201 *> IPIV is INTEGER array, dimension (N)
202 *> If FACT = 'F', then IPIV is an input argument and on entry
203 *> contains the pivot indices from the factorization A = L*U
204 *> as computed by DGBTRF; row i of the matrix was interchanged
205 *> with row IPIV(i).
206 *>
207 *> If FACT = 'N', then IPIV is an output argument and on exit
208 *> contains the pivot indices from the factorization A = L*U
209 *> of the original matrix A.
210 *>
211 *> If FACT = 'E', then IPIV is an output argument and on exit
212 *> contains the pivot indices from the factorization A = L*U
213 *> of the equilibrated matrix A.
214 *> \endverbatim
215 *>
216 *> \param[in,out] EQUED
217 *> \verbatim
218 *> EQUED is CHARACTER*1
219 *> Specifies the form of equilibration that was done.
220 *> = 'N': No equilibration (always true if FACT = 'N').
221 *> = 'R': Row equilibration, i.e., A has been premultiplied by
222 *> diag(R).
223 *> = 'C': Column equilibration, i.e., A has been postmultiplied
224 *> by diag(C).
225 *> = 'B': Both row and column equilibration, i.e., A has been
226 *> replaced by diag(R) * A * diag(C).
227 *> EQUED is an input argument if FACT = 'F'; otherwise, it is an
228 *> output argument.
229 *> \endverbatim
230 *>
231 *> \param[in,out] R
232 *> \verbatim
233 *> R is DOUBLE PRECISION array, dimension (N)
234 *> The row scale factors for A. If EQUED = 'R' or 'B', A is
235 *> multiplied on the left by diag(R); if EQUED = 'N' or 'C', R
236 *> is not accessed. R is an input argument if FACT = 'F';
237 *> otherwise, R is an output argument. If FACT = 'F' and
238 *> EQUED = 'R' or 'B', each element of R must be positive.
239 *> \endverbatim
240 *>
241 *> \param[in,out] C
242 *> \verbatim
243 *> C is DOUBLE PRECISION array, dimension (N)
244 *> The column scale factors for A. If EQUED = 'C' or 'B', A is
245 *> multiplied on the right by diag(C); if EQUED = 'N' or 'R', C
246 *> is not accessed. C is an input argument if FACT = 'F';
247 *> otherwise, C is an output argument. If FACT = 'F' and
248 *> EQUED = 'C' or 'B', each element of C must be positive.
249 *> \endverbatim
250 *>
251 *> \param[in,out] B
252 *> \verbatim
253 *> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
254 *> On entry, the right hand side matrix B.
255 *> On exit,
256 *> if EQUED = 'N', B is not modified;
257 *> if TRANS = 'N' and EQUED = 'R' or 'B', B is overwritten by
258 *> diag(R)*B;
259 *> if TRANS = 'T' or 'C' and EQUED = 'C' or 'B', B is
260 *> overwritten by diag(C)*B.
261 *> \endverbatim
262 *>
263 *> \param[in] LDB
264 *> \verbatim
265 *> LDB is INTEGER
266 *> The leading dimension of the array B. LDB >= max(1,N).
267 *> \endverbatim
268 *>
269 *> \param[out] X
270 *> \verbatim
271 *> X is DOUBLE PRECISION array, dimension (LDX,NRHS)
272 *> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X
273 *> to the original system of equations. Note that A and B are
274 *> modified on exit if EQUED .ne. 'N', and the solution to the
275 *> equilibrated system is inv(diag(C))*X if TRANS = 'N' and
276 *> EQUED = 'C' or 'B', or inv(diag(R))*X if TRANS = 'T' or 'C'
277 *> and EQUED = 'R' or 'B'.
278 *> \endverbatim
279 *>
280 *> \param[in] LDX
281 *> \verbatim
282 *> LDX is INTEGER
283 *> The leading dimension of the array X. LDX >= max(1,N).
284 *> \endverbatim
285 *>
286 *> \param[out] RCOND
287 *> \verbatim
288 *> RCOND is DOUBLE PRECISION
289 *> The estimate of the reciprocal condition number of the matrix
290 *> A after equilibration (if done). If RCOND is less than the
291 *> machine precision (in particular, if RCOND = 0), the matrix
292 *> is singular to working precision. This condition is
293 *> indicated by a return code of INFO > 0.
294 *> \endverbatim
295 *>
296 *> \param[out] FERR
297 *> \verbatim
298 *> FERR is DOUBLE PRECISION array, dimension (NRHS)
299 *> The estimated forward error bound for each solution vector
300 *> X(j) (the j-th column of the solution matrix X).
301 *> If XTRUE is the true solution corresponding to X(j), FERR(j)
302 *> is an estimated upper bound for the magnitude of the largest
303 *> element in (X(j) - XTRUE) divided by the magnitude of the
304 *> largest element in X(j). The estimate is as reliable as
305 *> the estimate for RCOND, and is almost always a slight
306 *> overestimate of the true error.
307 *> \endverbatim
308 *>
309 *> \param[out] BERR
310 *> \verbatim
311 *> BERR is DOUBLE PRECISION array, dimension (NRHS)
312 *> The componentwise relative backward error of each solution
313 *> vector X(j) (i.e., the smallest relative change in
314 *> any element of A or B that makes X(j) an exact solution).
315 *> \endverbatim
316 *>
317 *> \param[out] WORK
318 *> \verbatim
319 *> WORK is DOUBLE PRECISION array, dimension (3*N)
320 *> On exit, WORK(1) contains the reciprocal pivot growth
321 *> factor norm(A)/norm(U). The "max absolute element" norm is
322 *> used. If WORK(1) is much less than 1, then the stability
323 *> of the LU factorization of the (equilibrated) matrix A
324 *> could be poor. This also means that the solution X, condition
325 *> estimator RCOND, and forward error bound FERR could be
326 *> unreliable. If factorization fails with 0<INFO<=N, then
327 *> WORK(1) contains the reciprocal pivot growth factor for the
328 *> leading INFO columns of A.
329 *> \endverbatim
330 *>
331 *> \param[out] IWORK
332 *> \verbatim
333 *> IWORK is INTEGER array, dimension (N)
334 *> \endverbatim
335 *>
336 *> \param[out] INFO
337 *> \verbatim
338 *> INFO is INTEGER
339 *> = 0: successful exit
340 *> < 0: if INFO = -i, the i-th argument had an illegal value
341 *> > 0: if INFO = i, and i is
342 *> <= N: U(i,i) is exactly zero. The factorization
343 *> has been completed, but the factor U is exactly
344 *> singular, so the solution and error bounds
345 *> could not be computed. RCOND = 0 is returned.
346 *> = N+1: U is nonsingular, but RCOND is less than machine
347 *> precision, meaning that the matrix is singular
348 *> to working precision. Nevertheless, the
349 *> solution and error bounds are computed because
350 *> there are a number of situations where the
351 *> computed solution can be more accurate than the
352 *> value of RCOND would suggest.
353 *> \endverbatim
354 *
355 * Authors:
356 * ========
357 *
358 *> \author Univ. of Tennessee
359 *> \author Univ. of California Berkeley
360 *> \author Univ. of Colorado Denver
361 *> \author NAG Ltd.
362 *
363 *> \ingroup doubleGBsolve
364 *
365 * =====================================================================
366  SUBROUTINE dgbsvx( FACT, TRANS, N, KL, KU, NRHS, AB, LDAB, AFB,
367  $ LDAFB, IPIV, EQUED, R, C, B, LDB, X, LDX,
368  $ RCOND, FERR, BERR, WORK, IWORK, INFO )
369 *
370 * -- LAPACK driver routine --
371 * -- LAPACK is a software package provided by Univ. of Tennessee, --
372 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
373 *
374 * .. Scalar Arguments ..
375  CHARACTER EQUED, FACT, TRANS
376  INTEGER INFO, KL, KU, LDAB, LDAFB, LDB, LDX, N, NRHS
377  DOUBLE PRECISION RCOND
378 * ..
379 * .. Array Arguments ..
380  INTEGER IPIV( * ), IWORK( * )
381  DOUBLE PRECISION AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ),
382  $ berr( * ), c( * ), ferr( * ), r( * ),
383  $ work( * ), x( ldx, * )
384 * ..
385 *
386 * =====================================================================
387 *
388 * .. Parameters ..
389  DOUBLE PRECISION ZERO, ONE
390  PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
391 * ..
392 * .. Local Scalars ..
393  LOGICAL COLEQU, EQUIL, NOFACT, NOTRAN, ROWEQU
394  CHARACTER NORM
395  INTEGER I, INFEQU, J, J1, J2
396  DOUBLE PRECISION AMAX, ANORM, BIGNUM, COLCND, RCMAX, RCMIN,
397  $ rowcnd, rpvgrw, smlnum
398 * ..
399 * .. External Functions ..
400  LOGICAL LSAME
401  DOUBLE PRECISION DLAMCH, DLANGB, DLANTB
402  EXTERNAL lsame, dlamch, dlangb, dlantb
403 * ..
404 * .. External Subroutines ..
405  EXTERNAL dcopy, dgbcon, dgbequ, dgbrfs, dgbtrf, dgbtrs,
406  $ dlacpy, dlaqgb, xerbla
407 * ..
408 * .. Intrinsic Functions ..
409  INTRINSIC abs, max, min
410 * ..
411 * .. Executable Statements ..
412 *
413  info = 0
414  nofact = lsame( fact, 'N' )
415  equil = lsame( fact, 'E' )
416  notran = lsame( trans, 'N' )
417  IF( nofact .OR. equil ) THEN
418  equed = 'N'
419  rowequ = .false.
420  colequ = .false.
421  ELSE
422  rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
423  colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
424  smlnum = dlamch( 'Safe minimum' )
425  bignum = one / smlnum
426  END IF
427 *
428 * Test the input parameters.
429 *
430  IF( .NOT.nofact .AND. .NOT.equil .AND. .NOT.lsame( fact, 'F' ) )
431  $ THEN
432  info = -1
433  ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
434  $ lsame( trans, 'C' ) ) THEN
435  info = -2
436  ELSE IF( n.LT.0 ) THEN
437  info = -3
438  ELSE IF( kl.LT.0 ) THEN
439  info = -4
440  ELSE IF( ku.LT.0 ) THEN
441  info = -5
442  ELSE IF( nrhs.LT.0 ) THEN
443  info = -6
444  ELSE IF( ldab.LT.kl+ku+1 ) THEN
445  info = -8
446  ELSE IF( ldafb.LT.2*kl+ku+1 ) THEN
447  info = -10
448  ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
449  $ ( rowequ .OR. colequ .OR. lsame( equed, 'N' ) ) ) THEN
450  info = -12
451  ELSE
452  IF( rowequ ) THEN
453  rcmin = bignum
454  rcmax = zero
455  DO 10 j = 1, n
456  rcmin = min( rcmin, r( j ) )
457  rcmax = max( rcmax, r( j ) )
458  10 CONTINUE
459  IF( rcmin.LE.zero ) THEN
460  info = -13
461  ELSE IF( n.GT.0 ) THEN
462  rowcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
463  ELSE
464  rowcnd = one
465  END IF
466  END IF
467  IF( colequ .AND. info.EQ.0 ) THEN
468  rcmin = bignum
469  rcmax = zero
470  DO 20 j = 1, n
471  rcmin = min( rcmin, c( j ) )
472  rcmax = max( rcmax, c( j ) )
473  20 CONTINUE
474  IF( rcmin.LE.zero ) THEN
475  info = -14
476  ELSE IF( n.GT.0 ) THEN
477  colcnd = max( rcmin, smlnum ) / min( rcmax, bignum )
478  ELSE
479  colcnd = one
480  END IF
481  END IF
482  IF( info.EQ.0 ) THEN
483  IF( ldb.LT.max( 1, n ) ) THEN
484  info = -16
485  ELSE IF( ldx.LT.max( 1, n ) ) THEN
486  info = -18
487  END IF
488  END IF
489  END IF
490 *
491  IF( info.NE.0 ) THEN
492  CALL xerbla( 'DGBSVX', -info )
493  RETURN
494  END IF
495 *
496  IF( equil ) THEN
497 *
498 * Compute row and column scalings to equilibrate the matrix A.
499 *
500  CALL dgbequ( n, n, kl, ku, ab, ldab, r, c, rowcnd, colcnd,
501  $ amax, infequ )
502  IF( infequ.EQ.0 ) THEN
503 *
504 * Equilibrate the matrix.
505 *
506  CALL dlaqgb( n, n, kl, ku, ab, ldab, r, c, rowcnd, colcnd,
507  $ amax, equed )
508  rowequ = lsame( equed, 'R' ) .OR. lsame( equed, 'B' )
509  colequ = lsame( equed, 'C' ) .OR. lsame( equed, 'B' )
510  END IF
511  END IF
512 *
513 * Scale the right hand side.
514 *
515  IF( notran ) THEN
516  IF( rowequ ) THEN
517  DO 40 j = 1, nrhs
518  DO 30 i = 1, n
519  b( i, j ) = r( i )*b( i, j )
520  30 CONTINUE
521  40 CONTINUE
522  END IF
523  ELSE IF( colequ ) THEN
524  DO 60 j = 1, nrhs
525  DO 50 i = 1, n
526  b( i, j ) = c( i )*b( i, j )
527  50 CONTINUE
528  60 CONTINUE
529  END IF
530 *
531  IF( nofact .OR. equil ) THEN
532 *
533 * Compute the LU factorization of the band matrix A.
534 *
535  DO 70 j = 1, n
536  j1 = max( j-ku, 1 )
537  j2 = min( j+kl, n )
538  CALL dcopy( j2-j1+1, ab( ku+1-j+j1, j ), 1,
539  $ afb( kl+ku+1-j+j1, j ), 1 )
540  70 CONTINUE
541 *
542  CALL dgbtrf( n, n, kl, ku, afb, ldafb, ipiv, info )
543 *
544 * Return if INFO is non-zero.
545 *
546  IF( info.GT.0 ) THEN
547 *
548 * Compute the reciprocal pivot growth factor of the
549 * leading rank-deficient INFO columns of A.
550 *
551  anorm = zero
552  DO 90 j = 1, info
553  DO 80 i = max( ku+2-j, 1 ), min( n+ku+1-j, kl+ku+1 )
554  anorm = max( anorm, abs( ab( i, j ) ) )
555  80 CONTINUE
556  90 CONTINUE
557  rpvgrw = dlantb( 'M', 'U', 'N', info, min( info-1, kl+ku ),
558  $ afb( max( 1, kl+ku+2-info ), 1 ), ldafb,
559  $ work )
560  IF( rpvgrw.EQ.zero ) THEN
561  rpvgrw = one
562  ELSE
563  rpvgrw = anorm / rpvgrw
564  END IF
565  work( 1 ) = rpvgrw
566  rcond = zero
567  RETURN
568  END IF
569  END IF
570 *
571 * Compute the norm of the matrix A and the
572 * reciprocal pivot growth factor RPVGRW.
573 *
574  IF( notran ) THEN
575  norm = '1'
576  ELSE
577  norm = 'I'
578  END IF
579  anorm = dlangb( norm, n, kl, ku, ab, ldab, work )
580  rpvgrw = dlantb( 'M', 'U', 'N', n, kl+ku, afb, ldafb, work )
581  IF( rpvgrw.EQ.zero ) THEN
582  rpvgrw = one
583  ELSE
584  rpvgrw = dlangb( 'M', n, kl, ku, ab, ldab, work ) / rpvgrw
585  END IF
586 *
587 * Compute the reciprocal of the condition number of A.
588 *
589  CALL dgbcon( norm, n, kl, ku, afb, ldafb, ipiv, anorm, rcond,
590  $ work, iwork, info )
591 *
592 * Compute the solution matrix X.
593 *
594  CALL dlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
595  CALL dgbtrs( trans, n, kl, ku, nrhs, afb, ldafb, ipiv, x, ldx,
596  $ info )
597 *
598 * Use iterative refinement to improve the computed solution and
599 * compute error bounds and backward error estimates for it.
600 *
601  CALL dgbrfs( trans, n, kl, ku, nrhs, ab, ldab, afb, ldafb, ipiv,
602  $ b, ldb, x, ldx, ferr, berr, work, iwork, info )
603 *
604 * Transform the solution matrix X to a solution of the original
605 * system.
606 *
607  IF( notran ) THEN
608  IF( colequ ) THEN
609  DO 110 j = 1, nrhs
610  DO 100 i = 1, n
611  x( i, j ) = c( i )*x( i, j )
612  100 CONTINUE
613  110 CONTINUE
614  DO 120 j = 1, nrhs
615  ferr( j ) = ferr( j ) / colcnd
616  120 CONTINUE
617  END IF
618  ELSE IF( rowequ ) THEN
619  DO 140 j = 1, nrhs
620  DO 130 i = 1, n
621  x( i, j ) = r( i )*x( i, j )
622  130 CONTINUE
623  140 CONTINUE
624  DO 150 j = 1, nrhs
625  ferr( j ) = ferr( j ) / rowcnd
626  150 CONTINUE
627  END IF
628 *
629 * Set INFO = N+1 if the matrix is singular to working precision.
630 *
631  IF( rcond.LT.dlamch( 'Epsilon' ) )
632  $ info = n + 1
633 *
634  work( 1 ) = rpvgrw
635  RETURN
636 *
637 * End of DGBSVX
638 *
639  END
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:103
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:82
subroutine dlaqgb(M, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND, AMAX, EQUED)
DLAQGB scales a general band matrix, using row and column scaling factors computed by sgbequ.
Definition: dlaqgb.f:159
subroutine dgbtrs(TRANS, N, KL, KU, NRHS, AB, LDAB, IPIV, B, LDB, INFO)
DGBTRS
Definition: dgbtrs.f:138
subroutine dgbtrf(M, N, KL, KU, AB, LDAB, IPIV, INFO)
DGBTRF
Definition: dgbtrf.f:144
subroutine dgbequ(M, N, KL, KU, AB, LDAB, R, C, ROWCND, COLCND, AMAX, INFO)
DGBEQU
Definition: dgbequ.f:153
subroutine dgbcon(NORM, N, KL, KU, AB, LDAB, IPIV, ANORM, RCOND, WORK, IWORK, INFO)
DGBCON
Definition: dgbcon.f:146
subroutine dgbrfs(TRANS, N, KL, KU, NRHS, AB, LDAB, AFB, LDAFB, IPIV, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
DGBRFS
Definition: dgbrfs.f:205
subroutine dgbsvx(FACT, TRANS, N, KL, KU, NRHS, AB, LDAB, AFB, LDAFB, IPIV, EQUED, R, C, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO)
DGBSVX computes the solution to system of linear equations A * X = B for GB matrices
Definition: dgbsvx.f:369