LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ sgbrfs()

subroutine sgbrfs ( character trans,
integer n,
integer kl,
integer ku,
integer nrhs,
real, dimension( ldab, * ) ab,
integer ldab,
real, dimension( ldafb, * ) afb,
integer ldafb,
integer, dimension( * ) ipiv,
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 )

SGBRFS

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

Purpose:
!>
!> SGBRFS improves the computed solution to a system of linear
!> equations when the coefficient matrix is banded, and provides
!> error bounds and backward error estimates for the solution.
!> 
Parameters
[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]N
!>          N is INTEGER
!>          The order of the matrix A.  N >= 0.
!> 
[in]KL
!>          KL is INTEGER
!>          The number of subdiagonals within the band of A.  KL >= 0.
!> 
[in]KU
!>          KU is INTEGER
!>          The number of superdiagonals within the band of A.  KU >= 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]AB
!>          AB is REAL array, dimension (LDAB,N)
!>          The original band matrix A, stored in rows 1 to KL+KU+1.
!>          The j-th column of A is stored in the j-th column of the
!>          array AB as follows:
!>          AB(ku+1+i-j,j) = A(i,j) for max(1,j-ku)<=i<=min(n,j+kl).
!> 
[in]LDAB
!>          LDAB is INTEGER
!>          The leading dimension of the array AB.  LDAB >= KL+KU+1.
!> 
[in]AFB
!>          AFB is REAL array, dimension (LDAFB,N)
!>          Details of the LU factorization of the band matrix A, as
!>          computed by SGBTRF.  U is stored as an upper triangular band
!>          matrix with KL+KU superdiagonals in rows 1 to KL+KU+1, and
!>          the multipliers used during the factorization are stored in
!>          rows KL+KU+2 to 2*KL+KU+1.
!> 
[in]LDAFB
!>          LDAFB is INTEGER
!>          The leading dimension of the array AFB.  LDAFB >= 2*KL*KU+1.
!> 
[in]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          The pivot indices from SGBTRF; for 1<=i<=N, row i of the
!>          matrix was interchanged with row IPIV(i).
!> 
[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,out]X
!>          X is REAL array, dimension (LDX,NRHS)
!>          On entry, the solution matrix X, as computed by SGBTRS.
!>          On exit, the improved 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
!> 
Internal Parameters:
!>  ITMAX is the maximum number of steps of iterative refinement.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 200 of file sgbrfs.f.

204*
205* -- LAPACK computational routine --
206* -- LAPACK is a software package provided by Univ. of Tennessee, --
207* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
208*
209* .. Scalar Arguments ..
210 CHARACTER TRANS
211 INTEGER INFO, KL, KU, LDAB, LDAFB, LDB, LDX, N, NRHS
212* ..
213* .. Array Arguments ..
214 INTEGER IPIV( * ), IWORK( * )
215 REAL AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ),
216 $ BERR( * ), FERR( * ), WORK( * ), X( LDX, * )
217* ..
218*
219* =====================================================================
220*
221* .. Parameters ..
222 INTEGER ITMAX
223 parameter( itmax = 5 )
224 REAL ZERO
225 parameter( zero = 0.0e+0 )
226 REAL ONE
227 parameter( one = 1.0e+0 )
228 REAL TWO
229 parameter( two = 2.0e+0 )
230 REAL THREE
231 parameter( three = 3.0e+0 )
232* ..
233* .. Local Scalars ..
234 LOGICAL NOTRAN
235 CHARACTER TRANST
236 INTEGER COUNT, I, J, K, KASE, KK, NZ
237 REAL EPS, LSTRES, S, SAFE1, SAFE2, SAFMIN, XK
238* ..
239* .. Local Arrays ..
240 INTEGER ISAVE( 3 )
241* ..
242* .. External Subroutines ..
243 EXTERNAL saxpy, scopy, sgbmv, sgbtrs, slacn2,
244 $ xerbla
245* ..
246* .. Intrinsic Functions ..
247 INTRINSIC abs, max, min
248* ..
249* .. External Functions ..
250 LOGICAL LSAME
251 REAL SLAMCH
252 EXTERNAL lsame, slamch
253* ..
254* .. Executable Statements ..
255*
256* Test the input parameters.
257*
258 info = 0
259 notran = lsame( trans, 'N' )
260 IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
261 $ lsame( trans, 'C' ) ) THEN
262 info = -1
263 ELSE IF( n.LT.0 ) THEN
264 info = -2
265 ELSE IF( kl.LT.0 ) THEN
266 info = -3
267 ELSE IF( ku.LT.0 ) THEN
268 info = -4
269 ELSE IF( nrhs.LT.0 ) THEN
270 info = -5
271 ELSE IF( ldab.LT.kl+ku+1 ) THEN
272 info = -7
273 ELSE IF( ldafb.LT.2*kl+ku+1 ) THEN
274 info = -9
275 ELSE IF( ldb.LT.max( 1, n ) ) THEN
276 info = -12
277 ELSE IF( ldx.LT.max( 1, n ) ) THEN
278 info = -14
279 END IF
280 IF( info.NE.0 ) THEN
281 CALL xerbla( 'SGBRFS', -info )
282 RETURN
283 END IF
284*
285* Quick return if possible
286*
287 IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
288 DO 10 j = 1, nrhs
289 ferr( j ) = zero
290 berr( j ) = zero
291 10 CONTINUE
292 RETURN
293 END IF
294*
295 IF( notran ) THEN
296 transt = 'T'
297 ELSE
298 transt = 'N'
299 END IF
300*
301* NZ = maximum number of nonzero elements in each row of A, plus 1
302*
303 nz = min( kl+ku+2, n+1 )
304 eps = slamch( 'Epsilon' )
305 safmin = slamch( 'Safe minimum' )
306 safe1 = real( nz )*safmin
307 safe2 = safe1 / eps
308*
309* Do for each right hand side
310*
311 DO 140 j = 1, nrhs
312*
313 count = 1
314 lstres = three
315 20 CONTINUE
316*
317* Loop until stopping criterion is satisfied.
318*
319* Compute residual R = B - op(A) * X,
320* where op(A) = A, A**T, or A**H, depending on TRANS.
321*
322 CALL scopy( n, b( 1, j ), 1, work( n+1 ), 1 )
323 CALL sgbmv( trans, n, n, kl, ku, -one, ab, ldab, x( 1, j ),
324 $ 1,
325 $ one, work( n+1 ), 1 )
326*
327* Compute componentwise relative backward error from formula
328*
329* max(i) ( abs(R(i)) / ( abs(op(A))*abs(X) + abs(B) )(i) )
330*
331* where abs(Z) is the componentwise absolute value of the matrix
332* or vector Z. If the i-th component of the denominator is less
333* than SAFE2, then SAFE1 is added to the i-th components of the
334* numerator and denominator before dividing.
335*
336 DO 30 i = 1, n
337 work( i ) = abs( b( i, j ) )
338 30 CONTINUE
339*
340* Compute abs(op(A))*abs(X) + abs(B).
341*
342 IF( notran ) THEN
343 DO 50 k = 1, n
344 kk = ku + 1 - k
345 xk = abs( x( k, j ) )
346 DO 40 i = max( 1, k-ku ), min( n, k+kl )
347 work( i ) = work( i ) + abs( ab( kk+i, k ) )*xk
348 40 CONTINUE
349 50 CONTINUE
350 ELSE
351 DO 70 k = 1, n
352 s = zero
353 kk = ku + 1 - k
354 DO 60 i = max( 1, k-ku ), min( n, k+kl )
355 s = s + abs( ab( kk+i, k ) )*abs( x( i, j ) )
356 60 CONTINUE
357 work( k ) = work( k ) + s
358 70 CONTINUE
359 END IF
360 s = zero
361 DO 80 i = 1, n
362 IF( work( i ).GT.safe2 ) THEN
363 s = max( s, abs( work( n+i ) ) / work( i ) )
364 ELSE
365 s = max( s, ( abs( work( n+i ) )+safe1 ) /
366 $ ( work( i )+safe1 ) )
367 END IF
368 80 CONTINUE
369 berr( j ) = s
370*
371* Test stopping criterion. Continue iterating if
372* 1) The residual BERR(J) is larger than machine epsilon, and
373* 2) BERR(J) decreased by at least a factor of 2 during the
374* last iteration, and
375* 3) At most ITMAX iterations tried.
376*
377 IF( berr( j ).GT.eps .AND. two*berr( j ).LE.lstres .AND.
378 $ count.LE.itmax ) THEN
379*
380* Update solution and try again.
381*
382 CALL sgbtrs( trans, n, kl, ku, 1, afb, ldafb, ipiv,
383 $ work( n+1 ), n, info )
384 CALL saxpy( n, one, work( n+1 ), 1, x( 1, j ), 1 )
385 lstres = berr( j )
386 count = count + 1
387 GO TO 20
388 END IF
389*
390* Bound error from formula
391*
392* norm(X - XTRUE) / norm(X) .le. FERR =
393* norm( abs(inv(op(A)))*
394* ( abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) ))) / norm(X)
395*
396* where
397* norm(Z) is the magnitude of the largest component of Z
398* inv(op(A)) is the inverse of op(A)
399* abs(Z) is the componentwise absolute value of the matrix or
400* vector Z
401* NZ is the maximum number of nonzeros in any row of A, plus 1
402* EPS is machine epsilon
403*
404* The i-th component of abs(R)+NZ*EPS*(abs(op(A))*abs(X)+abs(B))
405* is incremented by SAFE1 if the i-th component of
406* abs(op(A))*abs(X) + abs(B) is less than SAFE2.
407*
408* Use SLACN2 to estimate the infinity-norm of the matrix
409* inv(op(A)) * diag(W),
410* where W = abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) )))
411*
412 DO 90 i = 1, n
413 IF( work( i ).GT.safe2 ) THEN
414 work( i ) = abs( work( n+i ) ) + real( nz )*eps*work( i )
415 ELSE
416 work( i ) = abs( work( n+i ) ) + real( nz )*eps*work( i )
417 $ + safe1
418 END IF
419 90 CONTINUE
420*
421 kase = 0
422 100 CONTINUE
423 CALL slacn2( n, work( 2*n+1 ), work( n+1 ), iwork,
424 $ ferr( j ),
425 $ kase, isave )
426 IF( kase.NE.0 ) THEN
427 IF( kase.EQ.1 ) THEN
428*
429* Multiply by diag(W)*inv(op(A)**T).
430*
431 CALL sgbtrs( transt, n, kl, ku, 1, afb, ldafb, ipiv,
432 $ work( n+1 ), n, info )
433 DO 110 i = 1, n
434 work( n+i ) = work( n+i )*work( i )
435 110 CONTINUE
436 ELSE
437*
438* Multiply by inv(op(A))*diag(W).
439*
440 DO 120 i = 1, n
441 work( n+i ) = work( n+i )*work( i )
442 120 CONTINUE
443 CALL sgbtrs( trans, n, kl, ku, 1, afb, ldafb, ipiv,
444 $ work( n+1 ), n, info )
445 END IF
446 GO TO 100
447 END IF
448*
449* Normalize error.
450*
451 lstres = zero
452 DO 130 i = 1, n
453 lstres = max( lstres, abs( x( i, j ) ) )
454 130 CONTINUE
455 IF( lstres.NE.zero )
456 $ ferr( j ) = ferr( j ) / lstres
457*
458 140 CONTINUE
459*
460 RETURN
461*
462* End of SGBRFS
463*
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 sgbmv(trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy)
SGBMV
Definition sgbmv.f:188
subroutine sgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
SGBTRS
Definition sgbtrs.f:137
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
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: