LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dtbrfs.f
Go to the documentation of this file.
1*> \brief \b DTBRFS
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DTBRFS + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dtbrfs.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dtbrfs.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dtbrfs.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DTBRFS( UPLO, TRANS, DIAG, N, KD, NRHS, AB, LDAB, B,
20* LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO )
21*
22* .. Scalar Arguments ..
23* CHARACTER DIAG, TRANS, UPLO
24* INTEGER INFO, KD, LDAB, LDB, LDX, N, NRHS
25* ..
26* .. Array Arguments ..
27* INTEGER IWORK( * )
28* DOUBLE PRECISION AB( LDAB, * ), B( LDB, * ), BERR( * ),
29* $ FERR( * ), WORK( * ), X( LDX, * )
30* ..
31*
32*
33*> \par Purpose:
34* =============
35*>
36*> \verbatim
37*>
38*> DTBRFS provides error bounds and backward error estimates for the
39*> solution to a system of linear equations with a triangular band
40*> coefficient matrix.
41*>
42*> The solution matrix X must be computed by DTBTRS or some other
43*> means before entering this routine. DTBRFS does not do iterative
44*> refinement because doing so cannot improve the backward error.
45*> \endverbatim
46*
47* Arguments:
48* ==========
49*
50*> \param[in] UPLO
51*> \verbatim
52*> UPLO is CHARACTER*1
53*> = 'U': A is upper triangular;
54*> = 'L': A is lower triangular.
55*> \endverbatim
56*>
57*> \param[in] TRANS
58*> \verbatim
59*> TRANS is CHARACTER*1
60*> Specifies the form of the system of equations:
61*> = 'N': A * X = B (No transpose)
62*> = 'T': A**T * X = B (Transpose)
63*> = 'C': A**H * X = B (Conjugate transpose = Transpose)
64*> \endverbatim
65*>
66*> \param[in] DIAG
67*> \verbatim
68*> DIAG is CHARACTER*1
69*> = 'N': A is non-unit triangular;
70*> = 'U': A is unit triangular.
71*> \endverbatim
72*>
73*> \param[in] N
74*> \verbatim
75*> N is INTEGER
76*> The order of the matrix A. N >= 0.
77*> \endverbatim
78*>
79*> \param[in] KD
80*> \verbatim
81*> KD is INTEGER
82*> The number of superdiagonals or subdiagonals of the
83*> triangular band matrix A. KD >= 0.
84*> \endverbatim
85*>
86*> \param[in] NRHS
87*> \verbatim
88*> NRHS is INTEGER
89*> The number of right hand sides, i.e., the number of columns
90*> of the matrices B and X. NRHS >= 0.
91*> \endverbatim
92*>
93*> \param[in] AB
94*> \verbatim
95*> AB is DOUBLE PRECISION array, dimension (LDAB,N)
96*> The upper or lower triangular band matrix A, stored in the
97*> first kd+1 rows of the array. The j-th column of A is stored
98*> in the j-th column of the array AB as follows:
99*> if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
100*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+kd).
101*> If DIAG = 'U', the diagonal elements of A are not referenced
102*> and are assumed to be 1.
103*> \endverbatim
104*>
105*> \param[in] LDAB
106*> \verbatim
107*> LDAB is INTEGER
108*> The leading dimension of the array AB. LDAB >= KD+1.
109*> \endverbatim
110*>
111*> \param[in] B
112*> \verbatim
113*> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
114*> The right hand side matrix B.
115*> \endverbatim
116*>
117*> \param[in] LDB
118*> \verbatim
119*> LDB is INTEGER
120*> The leading dimension of the array B. LDB >= max(1,N).
121*> \endverbatim
122*>
123*> \param[in] X
124*> \verbatim
125*> X is DOUBLE PRECISION array, dimension (LDX,NRHS)
126*> The solution matrix X.
127*> \endverbatim
128*>
129*> \param[in] LDX
130*> \verbatim
131*> LDX is INTEGER
132*> The leading dimension of the array X. LDX >= max(1,N).
133*> \endverbatim
134*>
135*> \param[out] FERR
136*> \verbatim
137*> FERR is DOUBLE PRECISION array, dimension (NRHS)
138*> The estimated forward error bound for each solution vector
139*> X(j) (the j-th column of the solution matrix X).
140*> If XTRUE is the true solution corresponding to X(j), FERR(j)
141*> is an estimated upper bound for the magnitude of the largest
142*> element in (X(j) - XTRUE) divided by the magnitude of the
143*> largest element in X(j). The estimate is as reliable as
144*> the estimate for RCOND, and is almost always a slight
145*> overestimate of the true error.
146*> \endverbatim
147*>
148*> \param[out] BERR
149*> \verbatim
150*> BERR is DOUBLE PRECISION array, dimension (NRHS)
151*> The componentwise relative backward error of each solution
152*> vector X(j) (i.e., the smallest relative change in
153*> any element of A or B that makes X(j) an exact solution).
154*> \endverbatim
155*>
156*> \param[out] WORK
157*> \verbatim
158*> WORK is DOUBLE PRECISION array, dimension (3*N)
159*> \endverbatim
160*>
161*> \param[out] IWORK
162*> \verbatim
163*> IWORK is INTEGER array, dimension (N)
164*> \endverbatim
165*>
166*> \param[out] INFO
167*> \verbatim
168*> INFO is INTEGER
169*> = 0: successful exit
170*> < 0: if INFO = -i, the i-th argument had an illegal value
171*> \endverbatim
172*
173* Authors:
174* ========
175*
176*> \author Univ. of Tennessee
177*> \author Univ. of California Berkeley
178*> \author Univ. of Colorado Denver
179*> \author NAG Ltd.
180*
181*> \ingroup tbrfs
182*
183* =====================================================================
184 SUBROUTINE dtbrfs( UPLO, TRANS, DIAG, N, KD, NRHS, AB, LDAB, B,
185 $ LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO )
186*
187* -- LAPACK computational routine --
188* -- LAPACK is a software package provided by Univ. of Tennessee, --
189* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
190*
191* .. Scalar Arguments ..
192 CHARACTER DIAG, TRANS, UPLO
193 INTEGER INFO, KD, LDAB, LDB, LDX, N, NRHS
194* ..
195* .. Array Arguments ..
196 INTEGER IWORK( * )
197 DOUBLE PRECISION AB( LDAB, * ), B( LDB, * ), BERR( * ),
198 $ ferr( * ), work( * ), x( ldx, * )
199* ..
200*
201* =====================================================================
202*
203* .. Parameters ..
204 DOUBLE PRECISION ZERO
205 parameter( zero = 0.0d+0 )
206 DOUBLE PRECISION ONE
207 parameter( one = 1.0d+0 )
208* ..
209* .. Local Scalars ..
210 LOGICAL NOTRAN, NOUNIT, UPPER
211 CHARACTER TRANST
212 INTEGER I, J, K, KASE, NZ
213 DOUBLE PRECISION EPS, LSTRES, S, SAFE1, SAFE2, SAFMIN, XK
214* ..
215* .. Local Arrays ..
216 INTEGER ISAVE( 3 )
217* ..
218* .. External Subroutines ..
219 EXTERNAL daxpy, dcopy, dlacn2, dtbmv, dtbsv,
220 $ xerbla
221* ..
222* .. Intrinsic Functions ..
223 INTRINSIC abs, max, min
224* ..
225* .. External Functions ..
226 LOGICAL LSAME
227 DOUBLE PRECISION DLAMCH
228 EXTERNAL lsame, dlamch
229* ..
230* .. Executable Statements ..
231*
232* Test the input parameters.
233*
234 info = 0
235 upper = lsame( uplo, 'U' )
236 notran = lsame( trans, 'N' )
237 nounit = lsame( diag, 'N' )
238*
239 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
240 info = -1
241 ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
242 $ lsame( trans, 'C' ) ) THEN
243 info = -2
244 ELSE IF( .NOT.nounit .AND. .NOT.lsame( diag, 'U' ) ) THEN
245 info = -3
246 ELSE IF( n.LT.0 ) THEN
247 info = -4
248 ELSE IF( kd.LT.0 ) THEN
249 info = -5
250 ELSE IF( nrhs.LT.0 ) THEN
251 info = -6
252 ELSE IF( ldab.LT.kd+1 ) THEN
253 info = -8
254 ELSE IF( ldb.LT.max( 1, n ) ) THEN
255 info = -10
256 ELSE IF( ldx.LT.max( 1, n ) ) THEN
257 info = -12
258 END IF
259 IF( info.NE.0 ) THEN
260 CALL xerbla( 'DTBRFS', -info )
261 RETURN
262 END IF
263*
264* Quick return if possible
265*
266 IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
267 DO 10 j = 1, nrhs
268 ferr( j ) = zero
269 berr( j ) = zero
270 10 CONTINUE
271 RETURN
272 END IF
273*
274 IF( notran ) THEN
275 transt = 'T'
276 ELSE
277 transt = 'N'
278 END IF
279*
280* NZ = maximum number of nonzero elements in each row of A, plus 1
281*
282 nz = kd + 2
283 eps = dlamch( 'Epsilon' )
284 safmin = dlamch( 'Safe minimum' )
285 safe1 = nz*safmin
286 safe2 = safe1 / eps
287*
288* Do for each right hand side
289*
290 DO 250 j = 1, nrhs
291*
292* Compute residual R = B - op(A) * X,
293* where op(A) = A or A**T, depending on TRANS.
294*
295 CALL dcopy( n, x( 1, j ), 1, work( n+1 ), 1 )
296 CALL dtbmv( uplo, trans, diag, n, kd, ab, ldab, work( n+1 ),
297 $ 1 )
298 CALL daxpy( n, -one, b( 1, j ), 1, work( n+1 ), 1 )
299*
300* Compute componentwise relative backward error from formula
301*
302* max(i) ( abs(R(i)) / ( abs(op(A))*abs(X) + abs(B) )(i) )
303*
304* where abs(Z) is the componentwise absolute value of the matrix
305* or vector Z. If the i-th component of the denominator is less
306* than SAFE2, then SAFE1 is added to the i-th components of the
307* numerator and denominator before dividing.
308*
309 DO 20 i = 1, n
310 work( i ) = abs( b( i, j ) )
311 20 CONTINUE
312*
313 IF( notran ) THEN
314*
315* Compute abs(A)*abs(X) + abs(B).
316*
317 IF( upper ) THEN
318 IF( nounit ) THEN
319 DO 40 k = 1, n
320 xk = abs( x( k, j ) )
321 DO 30 i = max( 1, k-kd ), k
322 work( i ) = work( i ) +
323 $ abs( ab( kd+1+i-k, k ) )*xk
324 30 CONTINUE
325 40 CONTINUE
326 ELSE
327 DO 60 k = 1, n
328 xk = abs( x( k, j ) )
329 DO 50 i = max( 1, k-kd ), k - 1
330 work( i ) = work( i ) +
331 $ abs( ab( kd+1+i-k, k ) )*xk
332 50 CONTINUE
333 work( k ) = work( k ) + xk
334 60 CONTINUE
335 END IF
336 ELSE
337 IF( nounit ) THEN
338 DO 80 k = 1, n
339 xk = abs( x( k, j ) )
340 DO 70 i = k, min( n, k+kd )
341 work( i ) = work( i ) + abs( ab( 1+i-k, k ) )*xk
342 70 CONTINUE
343 80 CONTINUE
344 ELSE
345 DO 100 k = 1, n
346 xk = abs( x( k, j ) )
347 DO 90 i = k + 1, min( n, k+kd )
348 work( i ) = work( i ) + abs( ab( 1+i-k, k ) )*xk
349 90 CONTINUE
350 work( k ) = work( k ) + xk
351 100 CONTINUE
352 END IF
353 END IF
354 ELSE
355*
356* Compute abs(A**T)*abs(X) + abs(B).
357*
358 IF( upper ) THEN
359 IF( nounit ) THEN
360 DO 120 k = 1, n
361 s = zero
362 DO 110 i = max( 1, k-kd ), k
363 s = s + abs( ab( kd+1+i-k, k ) )*
364 $ abs( x( i, j ) )
365 110 CONTINUE
366 work( k ) = work( k ) + s
367 120 CONTINUE
368 ELSE
369 DO 140 k = 1, n
370 s = abs( x( k, j ) )
371 DO 130 i = max( 1, k-kd ), k - 1
372 s = s + abs( ab( kd+1+i-k, k ) )*
373 $ abs( x( i, j ) )
374 130 CONTINUE
375 work( k ) = work( k ) + s
376 140 CONTINUE
377 END IF
378 ELSE
379 IF( nounit ) THEN
380 DO 160 k = 1, n
381 s = zero
382 DO 150 i = k, min( n, k+kd )
383 s = s + abs( ab( 1+i-k, k ) )*abs( x( i, j ) )
384 150 CONTINUE
385 work( k ) = work( k ) + s
386 160 CONTINUE
387 ELSE
388 DO 180 k = 1, n
389 s = abs( x( k, j ) )
390 DO 170 i = k + 1, min( n, k+kd )
391 s = s + abs( ab( 1+i-k, k ) )*abs( x( i, j ) )
392 170 CONTINUE
393 work( k ) = work( k ) + s
394 180 CONTINUE
395 END IF
396 END IF
397 END IF
398 s = zero
399 DO 190 i = 1, n
400 IF( work( i ).GT.safe2 ) THEN
401 s = max( s, abs( work( n+i ) ) / work( i ) )
402 ELSE
403 s = max( s, ( abs( work( n+i ) )+safe1 ) /
404 $ ( work( i )+safe1 ) )
405 END IF
406 190 CONTINUE
407 berr( j ) = s
408*
409* Bound error from formula
410*
411* norm(X - XTRUE) / norm(X) .le. FERR =
412* norm( abs(inv(op(A)))*
413* ( abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) ))) / norm(X)
414*
415* where
416* norm(Z) is the magnitude of the largest component of Z
417* inv(op(A)) is the inverse of op(A)
418* abs(Z) is the componentwise absolute value of the matrix or
419* vector Z
420* NZ is the maximum number of nonzeros in any row of A, plus 1
421* EPS is machine epsilon
422*
423* The i-th component of abs(R)+NZ*EPS*(abs(op(A))*abs(X)+abs(B))
424* is incremented by SAFE1 if the i-th component of
425* abs(op(A))*abs(X) + abs(B) is less than SAFE2.
426*
427* Use DLACN2 to estimate the infinity-norm of the matrix
428* inv(op(A)) * diag(W),
429* where W = abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) )))
430*
431 DO 200 i = 1, n
432 IF( work( i ).GT.safe2 ) THEN
433 work( i ) = abs( work( n+i ) ) + nz*eps*work( i )
434 ELSE
435 work( i ) = abs( work( n+i ) ) + nz*eps*work( i ) + safe1
436 END IF
437 200 CONTINUE
438*
439 kase = 0
440 210 CONTINUE
441 CALL dlacn2( n, work( 2*n+1 ), work( n+1 ), iwork,
442 $ ferr( j ),
443 $ kase, isave )
444 IF( kase.NE.0 ) THEN
445 IF( kase.EQ.1 ) THEN
446*
447* Multiply by diag(W)*inv(op(A)**T).
448*
449 CALL dtbsv( uplo, transt, diag, n, kd, ab, ldab,
450 $ work( n+1 ), 1 )
451 DO 220 i = 1, n
452 work( n+i ) = work( i )*work( n+i )
453 220 CONTINUE
454 ELSE
455*
456* Multiply by inv(op(A))*diag(W).
457*
458 DO 230 i = 1, n
459 work( n+i ) = work( i )*work( n+i )
460 230 CONTINUE
461 CALL dtbsv( uplo, trans, diag, n, kd, ab, ldab,
462 $ work( n+1 ), 1 )
463 END IF
464 GO TO 210
465 END IF
466*
467* Normalize error.
468*
469 lstres = zero
470 DO 240 i = 1, n
471 lstres = max( lstres, abs( x( i, j ) ) )
472 240 CONTINUE
473 IF( lstres.NE.zero )
474 $ ferr( j ) = ferr( j ) / lstres
475*
476 250 CONTINUE
477*
478 RETURN
479*
480* End of DTBRFS
481*
482 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine daxpy(n, da, dx, incx, dy, incy)
DAXPY
Definition daxpy.f:89
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dlacn2(n, v, x, isgn, est, kase, isave)
DLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition dlacn2.f:134
subroutine dtbmv(uplo, trans, diag, n, k, a, lda, x, incx)
DTBMV
Definition dtbmv.f:186
subroutine dtbrfs(uplo, trans, diag, n, kd, nrhs, ab, ldab, b, ldb, x, ldx, ferr, berr, work, iwork, info)
DTBRFS
Definition dtbrfs.f:186
subroutine dtbsv(uplo, trans, diag, n, k, a, lda, x, incx)
DTBSV
Definition dtbsv.f:189