LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine strrfs ( character  UPLO,
character  TRANS,
character  DIAG,
integer  N,
integer  NRHS,
real, dimension( lda, * )  A,
integer  LDA,
real, dimension( ldb, * )  B,
integer  LDB,
real, dimension( ldx, * )  X,
integer  LDX,
real, dimension( * )  FERR,
real, dimension( * )  BERR,
real, dimension( * )  WORK,
integer, dimension( * )  IWORK,
integer  INFO 
)

STRRFS

Download STRRFS + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 STRRFS provides error bounds and backward error estimates for the
 solution to a system of linear equations with a triangular
 coefficient matrix.

 The solution matrix X must be computed by STRTRS or some other
 means before entering this routine.  STRRFS does not do iterative
 refinement because doing so cannot improve the backward error.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  A is upper triangular;
          = 'L':  A is lower triangular.
[in]TRANS
          TRANS is CHARACTER*1
          Specifies the form of the system of equations:
          = 'N':  A * X = B  (No transpose)
          = 'T':  A**T * X = B  (Transpose)
          = 'C':  A**H * X = B  (Conjugate transpose = Transpose)
[in]DIAG
          DIAG is CHARACTER*1
          = 'N':  A is non-unit triangular;
          = 'U':  A is unit triangular.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]NRHS
          NRHS is INTEGER
          The number of right hand sides, i.e., the number of columns
          of the matrices B and X.  NRHS >= 0.
[in]A
          A is REAL array, dimension (LDA,N)
          The triangular matrix A.  If UPLO = 'U', the leading N-by-N
          upper triangular part of the array A contains the upper
          triangular matrix, and the strictly lower triangular part of
          A is not referenced.  If UPLO = 'L', the leading N-by-N lower
          triangular part of the array A contains the lower triangular
          matrix, and the strictly upper triangular part of A is not
          referenced.  If DIAG = 'U', the diagonal elements of A are
          also not referenced and are assumed to be 1.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in]B
          B is REAL array, dimension (LDB,NRHS)
          The right hand side matrix B.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[in]X
          X is REAL array, dimension (LDX,NRHS)
          The solution matrix X.
[in]LDX
          LDX is INTEGER
          The leading dimension of the array X.  LDX >= max(1,N).
[out]FERR
          FERR is REAL array, dimension (NRHS)
          The estimated forward error bound for each solution vector
          X(j) (the j-th column of the solution matrix X).
          If XTRUE is the true solution corresponding to X(j), FERR(j)
          is an estimated upper bound for the magnitude of the largest
          element in (X(j) - XTRUE) divided by the magnitude of the
          largest element in X(j).  The estimate is as reliable as
          the estimate for RCOND, and is almost always a slight
          overestimate of the true error.
[out]BERR
          BERR is REAL array, dimension (NRHS)
          The componentwise relative backward error of each solution
          vector X(j) (i.e., the smallest relative change in
          any element of A or B that makes X(j) an exact solution).
[out]WORK
          WORK is REAL array, dimension (3*N)
[out]IWORK
          IWORK is INTEGER array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 184 of file strrfs.f.

184 *
185 * -- LAPACK computational routine (version 3.4.0) --
186 * -- LAPACK is a software package provided by Univ. of Tennessee, --
187 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
188 * November 2011
189 *
190 * .. Scalar Arguments ..
191  CHARACTER diag, trans, uplo
192  INTEGER info, lda, ldb, ldx, n, nrhs
193 * ..
194 * .. Array Arguments ..
195  INTEGER iwork( * )
196  REAL a( lda, * ), b( ldb, * ), berr( * ), ferr( * ),
197  $ work( * ), x( ldx, * )
198 * ..
199 *
200 * =====================================================================
201 *
202 * .. Parameters ..
203  REAL zero
204  parameter ( zero = 0.0e+0 )
205  REAL one
206  parameter ( one = 1.0e+0 )
207 * ..
208 * .. Local Scalars ..
209  LOGICAL notran, nounit, upper
210  CHARACTER transt
211  INTEGER i, j, k, kase, nz
212  REAL eps, lstres, s, safe1, safe2, safmin, xk
213 * ..
214 * .. Local Arrays ..
215  INTEGER isave( 3 )
216 * ..
217 * .. External Subroutines ..
218  EXTERNAL saxpy, scopy, slacn2, strmv, strsv, xerbla
219 * ..
220 * .. Intrinsic Functions ..
221  INTRINSIC abs, max
222 * ..
223 * .. External Functions ..
224  LOGICAL lsame
225  REAL slamch
226  EXTERNAL lsame, slamch
227 * ..
228 * .. Executable Statements ..
229 *
230 * Test the input parameters.
231 *
232  info = 0
233  upper = lsame( uplo, 'U' )
234  notran = lsame( trans, 'N' )
235  nounit = lsame( diag, 'N' )
236 *
237  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
238  info = -1
239  ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
240  $ lsame( trans, 'C' ) ) THEN
241  info = -2
242  ELSE IF( .NOT.nounit .AND. .NOT.lsame( diag, 'U' ) ) THEN
243  info = -3
244  ELSE IF( n.LT.0 ) THEN
245  info = -4
246  ELSE IF( nrhs.LT.0 ) THEN
247  info = -5
248  ELSE IF( lda.LT.max( 1, n ) ) THEN
249  info = -7
250  ELSE IF( ldb.LT.max( 1, n ) ) THEN
251  info = -9
252  ELSE IF( ldx.LT.max( 1, n ) ) THEN
253  info = -11
254  END IF
255  IF( info.NE.0 ) THEN
256  CALL xerbla( 'STRRFS', -info )
257  RETURN
258  END IF
259 *
260 * Quick return if possible
261 *
262  IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
263  DO 10 j = 1, nrhs
264  ferr( j ) = zero
265  berr( j ) = zero
266  10 CONTINUE
267  RETURN
268  END IF
269 *
270  IF( notran ) THEN
271  transt = 'T'
272  ELSE
273  transt = 'N'
274  END IF
275 *
276 * NZ = maximum number of nonzero elements in each row of A, plus 1
277 *
278  nz = n + 1
279  eps = slamch( 'Epsilon' )
280  safmin = slamch( 'Safe minimum' )
281  safe1 = nz*safmin
282  safe2 = safe1 / eps
283 *
284 * Do for each right hand side
285 *
286  DO 250 j = 1, nrhs
287 *
288 * Compute residual R = B - op(A) * X,
289 * where op(A) = A or A**T, depending on TRANS.
290 *
291  CALL scopy( n, x( 1, j ), 1, work( n+1 ), 1 )
292  CALL strmv( uplo, trans, diag, n, a, lda, work( n+1 ), 1 )
293  CALL saxpy( n, -one, b( 1, j ), 1, work( n+1 ), 1 )
294 *
295 * Compute componentwise relative backward error from formula
296 *
297 * max(i) ( abs(R(i)) / ( abs(op(A))*abs(X) + abs(B) )(i) )
298 *
299 * where abs(Z) is the componentwise absolute value of the matrix
300 * or vector Z. If the i-th component of the denominator is less
301 * than SAFE2, then SAFE1 is added to the i-th components of the
302 * numerator and denominator before dividing.
303 *
304  DO 20 i = 1, n
305  work( i ) = abs( b( i, j ) )
306  20 CONTINUE
307 *
308  IF( notran ) THEN
309 *
310 * Compute abs(A)*abs(X) + abs(B).
311 *
312  IF( upper ) THEN
313  IF( nounit ) THEN
314  DO 40 k = 1, n
315  xk = abs( x( k, j ) )
316  DO 30 i = 1, k
317  work( i ) = work( i ) + abs( a( i, k ) )*xk
318  30 CONTINUE
319  40 CONTINUE
320  ELSE
321  DO 60 k = 1, n
322  xk = abs( x( k, j ) )
323  DO 50 i = 1, k - 1
324  work( i ) = work( i ) + abs( a( i, k ) )*xk
325  50 CONTINUE
326  work( k ) = work( k ) + xk
327  60 CONTINUE
328  END IF
329  ELSE
330  IF( nounit ) THEN
331  DO 80 k = 1, n
332  xk = abs( x( k, j ) )
333  DO 70 i = k, n
334  work( i ) = work( i ) + abs( a( i, k ) )*xk
335  70 CONTINUE
336  80 CONTINUE
337  ELSE
338  DO 100 k = 1, n
339  xk = abs( x( k, j ) )
340  DO 90 i = k + 1, n
341  work( i ) = work( i ) + abs( a( i, k ) )*xk
342  90 CONTINUE
343  work( k ) = work( k ) + xk
344  100 CONTINUE
345  END IF
346  END IF
347  ELSE
348 *
349 * Compute abs(A**T)*abs(X) + abs(B).
350 *
351  IF( upper ) THEN
352  IF( nounit ) THEN
353  DO 120 k = 1, n
354  s = zero
355  DO 110 i = 1, k
356  s = s + abs( a( i, k ) )*abs( x( i, j ) )
357  110 CONTINUE
358  work( k ) = work( k ) + s
359  120 CONTINUE
360  ELSE
361  DO 140 k = 1, n
362  s = abs( x( k, j ) )
363  DO 130 i = 1, k - 1
364  s = s + abs( a( i, k ) )*abs( x( i, j ) )
365  130 CONTINUE
366  work( k ) = work( k ) + s
367  140 CONTINUE
368  END IF
369  ELSE
370  IF( nounit ) THEN
371  DO 160 k = 1, n
372  s = zero
373  DO 150 i = k, n
374  s = s + abs( a( i, k ) )*abs( x( i, j ) )
375  150 CONTINUE
376  work( k ) = work( k ) + s
377  160 CONTINUE
378  ELSE
379  DO 180 k = 1, n
380  s = abs( x( k, j ) )
381  DO 170 i = k + 1, n
382  s = s + abs( a( i, k ) )*abs( x( i, j ) )
383  170 CONTINUE
384  work( k ) = work( k ) + s
385  180 CONTINUE
386  END IF
387  END IF
388  END IF
389  s = zero
390  DO 190 i = 1, n
391  IF( work( i ).GT.safe2 ) THEN
392  s = max( s, abs( work( n+i ) ) / work( i ) )
393  ELSE
394  s = max( s, ( abs( work( n+i ) )+safe1 ) /
395  $ ( work( i )+safe1 ) )
396  END IF
397  190 CONTINUE
398  berr( j ) = s
399 *
400 * Bound error from formula
401 *
402 * norm(X - XTRUE) / norm(X) .le. FERR =
403 * norm( abs(inv(op(A)))*
404 * ( abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) ))) / norm(X)
405 *
406 * where
407 * norm(Z) is the magnitude of the largest component of Z
408 * inv(op(A)) is the inverse of op(A)
409 * abs(Z) is the componentwise absolute value of the matrix or
410 * vector Z
411 * NZ is the maximum number of nonzeros in any row of A, plus 1
412 * EPS is machine epsilon
413 *
414 * The i-th component of abs(R)+NZ*EPS*(abs(op(A))*abs(X)+abs(B))
415 * is incremented by SAFE1 if the i-th component of
416 * abs(op(A))*abs(X) + abs(B) is less than SAFE2.
417 *
418 * Use SLACN2 to estimate the infinity-norm of the matrix
419 * inv(op(A)) * diag(W),
420 * where W = abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) )))
421 *
422  DO 200 i = 1, n
423  IF( work( i ).GT.safe2 ) THEN
424  work( i ) = abs( work( n+i ) ) + nz*eps*work( i )
425  ELSE
426  work( i ) = abs( work( n+i ) ) + nz*eps*work( i ) + safe1
427  END IF
428  200 CONTINUE
429 *
430  kase = 0
431  210 CONTINUE
432  CALL slacn2( n, work( 2*n+1 ), work( n+1 ), iwork, ferr( j ),
433  $ kase, isave )
434  IF( kase.NE.0 ) THEN
435  IF( kase.EQ.1 ) THEN
436 *
437 * Multiply by diag(W)*inv(op(A)**T).
438 *
439  CALL strsv( uplo, transt, diag, n, a, lda, 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 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine strmv(UPLO, TRANS, DIAG, N, A, LDA, X, INCX)
STRMV
Definition: strmv.f:149
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:138
subroutine strsv(UPLO, TRANS, DIAG, N, A, LDA, X, INCX)
STRSV
Definition: strsv.f:151
subroutine saxpy(N, SA, SX, INCX, SY, INCY)
SAXPY
Definition: saxpy.f:54
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:53

Here is the call graph for this function:

Here is the caller graph for this function: