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