LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cptrfs ( character  UPLO,
integer  N,
integer  NRHS,
real, dimension( * )  D,
complex, dimension( * )  E,
real, dimension( * )  DF,
complex, dimension( * )  EF,
complex, dimension( ldb, * )  B,
integer  LDB,
complex, dimension( ldx, * )  X,
integer  LDX,
real, dimension( * )  FERR,
real, dimension( * )  BERR,
complex, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer  INFO 
)

CPTRFS

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

Purpose:
 CPTRFS improves the computed solution to a system of linear
 equations when the coefficient matrix is Hermitian positive definite
 and tridiagonal, and provides error bounds and backward error
 estimates for the solution.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the superdiagonal or the subdiagonal of the
          tridiagonal matrix A is stored and the form of the
          factorization:
          = 'U':  E is the superdiagonal of A, and A = U**H*D*U;
          = 'L':  E is the subdiagonal of A, and A = L*D*L**H.
          (The two forms are equivalent if A is real.)
[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 matrix B.  NRHS >= 0.
[in]D
          D is REAL array, dimension (N)
          The n real diagonal elements of the tridiagonal matrix A.
[in]E
          E is COMPLEX array, dimension (N-1)
          The (n-1) off-diagonal elements of the tridiagonal matrix A
          (see UPLO).
[in]DF
          DF is REAL array, dimension (N)
          The n diagonal elements of the diagonal matrix D from
          the factorization computed by CPTTRF.
[in]EF
          EF is COMPLEX array, dimension (N-1)
          The (n-1) off-diagonal elements of the unit bidiagonal
          factor U or L from the factorization computed by CPTTRF
          (see UPLO).
[in]B
          B is COMPLEX 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 COMPLEX array, dimension (LDX,NRHS)
          On entry, the solution matrix X, as computed by CPTTRS.
          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 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).
[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 COMPLEX array, dimension (N)
[out]RWORK
          RWORK is REAL 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.
Date
September 2012

Definition at line 185 of file cptrfs.f.

185 *
186 * -- LAPACK computational routine (version 3.4.2) --
187 * -- LAPACK is a software package provided by Univ. of Tennessee, --
188 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
189 * September 2012
190 *
191 * .. Scalar Arguments ..
192  CHARACTER uplo
193  INTEGER info, ldb, ldx, n, nrhs
194 * ..
195 * .. Array Arguments ..
196  REAL berr( * ), d( * ), df( * ), ferr( * ),
197  $ rwork( * )
198  COMPLEX b( ldb, * ), e( * ), ef( * ), work( * ),
199  $ x( ldx, * )
200 * ..
201 *
202 * =====================================================================
203 *
204 * .. Parameters ..
205  INTEGER itmax
206  parameter ( itmax = 5 )
207  REAL zero
208  parameter ( zero = 0.0e+0 )
209  REAL one
210  parameter ( one = 1.0e+0 )
211  REAL two
212  parameter ( two = 2.0e+0 )
213  REAL three
214  parameter ( three = 3.0e+0 )
215 * ..
216 * .. Local Scalars ..
217  LOGICAL upper
218  INTEGER count, i, ix, j, nz
219  REAL eps, lstres, s, safe1, safe2, safmin
220  COMPLEX bi, cx, dx, ex, zdum
221 * ..
222 * .. External Functions ..
223  LOGICAL lsame
224  INTEGER isamax
225  REAL slamch
226  EXTERNAL lsame, isamax, slamch
227 * ..
228 * .. External Subroutines ..
229  EXTERNAL caxpy, cpttrs, xerbla
230 * ..
231 * .. Intrinsic Functions ..
232  INTRINSIC abs, aimag, cmplx, conjg, max, real
233 * ..
234 * .. Statement Functions ..
235  REAL cabs1
236 * ..
237 * .. Statement Function definitions ..
238  cabs1( zdum ) = abs( REAL( ZDUM ) ) + abs( aimag( zdum ) )
239 * ..
240 * .. Executable Statements ..
241 *
242 * Test the input parameters.
243 *
244  info = 0
245  upper = lsame( uplo, 'U' )
246  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
247  info = -1
248  ELSE IF( n.LT.0 ) THEN
249  info = -2
250  ELSE IF( nrhs.LT.0 ) THEN
251  info = -3
252  ELSE IF( ldb.LT.max( 1, n ) ) THEN
253  info = -9
254  ELSE IF( ldx.LT.max( 1, n ) ) THEN
255  info = -11
256  END IF
257  IF( info.NE.0 ) THEN
258  CALL xerbla( 'CPTRFS', -info )
259  RETURN
260  END IF
261 *
262 * Quick return if possible
263 *
264  IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
265  DO 10 j = 1, nrhs
266  ferr( j ) = zero
267  berr( j ) = zero
268  10 CONTINUE
269  RETURN
270  END IF
271 *
272 * NZ = maximum number of nonzero elements in each row of A, plus 1
273 *
274  nz = 4
275  eps = slamch( 'Epsilon' )
276  safmin = slamch( 'Safe minimum' )
277  safe1 = nz*safmin
278  safe2 = safe1 / eps
279 *
280 * Do for each right hand side
281 *
282  DO 100 j = 1, nrhs
283 *
284  count = 1
285  lstres = three
286  20 CONTINUE
287 *
288 * Loop until stopping criterion is satisfied.
289 *
290 * Compute residual R = B - A * X. Also compute
291 * abs(A)*abs(x) + abs(b) for use in the backward error bound.
292 *
293  IF( upper ) THEN
294  IF( n.EQ.1 ) THEN
295  bi = b( 1, j )
296  dx = d( 1 )*x( 1, j )
297  work( 1 ) = bi - dx
298  rwork( 1 ) = cabs1( bi ) + cabs1( dx )
299  ELSE
300  bi = b( 1, j )
301  dx = d( 1 )*x( 1, j )
302  ex = e( 1 )*x( 2, j )
303  work( 1 ) = bi - dx - ex
304  rwork( 1 ) = cabs1( bi ) + cabs1( dx ) +
305  $ cabs1( e( 1 ) )*cabs1( x( 2, j ) )
306  DO 30 i = 2, n - 1
307  bi = b( i, j )
308  cx = conjg( e( i-1 ) )*x( i-1, j )
309  dx = d( i )*x( i, j )
310  ex = e( i )*x( i+1, j )
311  work( i ) = bi - cx - dx - ex
312  rwork( i ) = cabs1( bi ) +
313  $ cabs1( e( i-1 ) )*cabs1( x( i-1, j ) ) +
314  $ cabs1( dx ) + cabs1( e( i ) )*
315  $ cabs1( x( i+1, j ) )
316  30 CONTINUE
317  bi = b( n, j )
318  cx = conjg( e( n-1 ) )*x( n-1, j )
319  dx = d( n )*x( n, j )
320  work( n ) = bi - cx - dx
321  rwork( n ) = cabs1( bi ) + cabs1( e( n-1 ) )*
322  $ cabs1( x( n-1, j ) ) + cabs1( dx )
323  END IF
324  ELSE
325  IF( n.EQ.1 ) THEN
326  bi = b( 1, j )
327  dx = d( 1 )*x( 1, j )
328  work( 1 ) = bi - dx
329  rwork( 1 ) = cabs1( bi ) + cabs1( dx )
330  ELSE
331  bi = b( 1, j )
332  dx = d( 1 )*x( 1, j )
333  ex = conjg( e( 1 ) )*x( 2, j )
334  work( 1 ) = bi - dx - ex
335  rwork( 1 ) = cabs1( bi ) + cabs1( dx ) +
336  $ cabs1( e( 1 ) )*cabs1( x( 2, j ) )
337  DO 40 i = 2, n - 1
338  bi = b( i, j )
339  cx = e( i-1 )*x( i-1, j )
340  dx = d( i )*x( i, j )
341  ex = conjg( e( i ) )*x( i+1, j )
342  work( i ) = bi - cx - dx - ex
343  rwork( i ) = cabs1( bi ) +
344  $ cabs1( e( i-1 ) )*cabs1( x( i-1, j ) ) +
345  $ cabs1( dx ) + cabs1( e( i ) )*
346  $ cabs1( x( i+1, j ) )
347  40 CONTINUE
348  bi = b( n, j )
349  cx = e( n-1 )*x( n-1, j )
350  dx = d( n )*x( n, j )
351  work( n ) = bi - cx - dx
352  rwork( n ) = cabs1( bi ) + cabs1( e( n-1 ) )*
353  $ cabs1( x( n-1, j ) ) + cabs1( dx )
354  END IF
355  END IF
356 *
357 * Compute componentwise relative backward error from formula
358 *
359 * max(i) ( abs(R(i)) / ( abs(A)*abs(X) + abs(B) )(i) )
360 *
361 * where abs(Z) is the componentwise absolute value of the matrix
362 * or vector Z. If the i-th component of the denominator is less
363 * than SAFE2, then SAFE1 is added to the i-th components of the
364 * numerator and denominator before dividing.
365 *
366  s = zero
367  DO 50 i = 1, n
368  IF( rwork( i ).GT.safe2 ) THEN
369  s = max( s, cabs1( work( i ) ) / rwork( i ) )
370  ELSE
371  s = max( s, ( cabs1( work( i ) )+safe1 ) /
372  $ ( rwork( i )+safe1 ) )
373  END IF
374  50 CONTINUE
375  berr( j ) = s
376 *
377 * Test stopping criterion. Continue iterating if
378 * 1) The residual BERR(J) is larger than machine epsilon, and
379 * 2) BERR(J) decreased by at least a factor of 2 during the
380 * last iteration, and
381 * 3) At most ITMAX iterations tried.
382 *
383  IF( berr( j ).GT.eps .AND. two*berr( j ).LE.lstres .AND.
384  $ count.LE.itmax ) THEN
385 *
386 * Update solution and try again.
387 *
388  CALL cpttrs( uplo, n, 1, df, ef, work, n, info )
389  CALL caxpy( n, cmplx( one ), work, 1, x( 1, j ), 1 )
390  lstres = berr( j )
391  count = count + 1
392  GO TO 20
393  END IF
394 *
395 * Bound error from formula
396 *
397 * norm(X - XTRUE) / norm(X) .le. FERR =
398 * norm( abs(inv(A))*
399 * ( abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) ))) / norm(X)
400 *
401 * where
402 * norm(Z) is the magnitude of the largest component of Z
403 * inv(A) is the inverse of A
404 * abs(Z) is the componentwise absolute value of the matrix or
405 * vector Z
406 * NZ is the maximum number of nonzeros in any row of A, plus 1
407 * EPS is machine epsilon
408 *
409 * The i-th component of abs(R)+NZ*EPS*(abs(A)*abs(X)+abs(B))
410 * is incremented by SAFE1 if the i-th component of
411 * abs(A)*abs(X) + abs(B) is less than SAFE2.
412 *
413  DO 60 i = 1, n
414  IF( rwork( i ).GT.safe2 ) THEN
415  rwork( i ) = cabs1( work( i ) ) + nz*eps*rwork( i )
416  ELSE
417  rwork( i ) = cabs1( work( i ) ) + nz*eps*rwork( i ) +
418  $ safe1
419  END IF
420  60 CONTINUE
421  ix = isamax( n, rwork, 1 )
422  ferr( j ) = rwork( ix )
423 *
424 * Estimate the norm of inv(A).
425 *
426 * Solve M(A) * x = e, where M(A) = (m(i,j)) is given by
427 *
428 * m(i,j) = abs(A(i,j)), i = j,
429 * m(i,j) = -abs(A(i,j)), i .ne. j,
430 *
431 * and e = [ 1, 1, ..., 1 ]**T. Note M(A) = M(L)*D*M(L)**H.
432 *
433 * Solve M(L) * x = e.
434 *
435  rwork( 1 ) = one
436  DO 70 i = 2, n
437  rwork( i ) = one + rwork( i-1 )*abs( ef( i-1 ) )
438  70 CONTINUE
439 *
440 * Solve D * M(L)**H * x = b.
441 *
442  rwork( n ) = rwork( n ) / df( n )
443  DO 80 i = n - 1, 1, -1
444  rwork( i ) = rwork( i ) / df( i ) +
445  $ rwork( i+1 )*abs( ef( i ) )
446  80 CONTINUE
447 *
448 * Compute norm(inv(A)) = max(x(i)), 1<=i<=n.
449 *
450  ix = isamax( n, rwork, 1 )
451  ferr( j ) = ferr( j )*abs( rwork( ix ) )
452 *
453 * Normalize error.
454 *
455  lstres = zero
456  DO 90 i = 1, n
457  lstres = max( lstres, abs( x( i, j ) ) )
458  90 CONTINUE
459  IF( lstres.NE.zero )
460  $ ferr( j ) = ferr( j ) / lstres
461 *
462  100 CONTINUE
463 *
464  RETURN
465 *
466 * End of CPTRFS
467 *
integer function isamax(N, SX, INCX)
ISAMAX
Definition: isamax.f:53
subroutine cpttrs(UPLO, N, NRHS, D, E, B, LDB, INFO)
CPTTRS
Definition: cpttrs.f:123
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine caxpy(N, CA, CX, INCX, CY, INCY)
CAXPY
Definition: caxpy.f:53
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: