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