3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download ZPTRFS + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zptrfs.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zptrfs.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zptrfs.f">
21 * SUBROUTINE ZPTRFS( UPLO, N, NRHS, D, E, DF, EF, B, LDB, X, LDX,
22 * FERR, BERR, WORK, RWORK, INFO )
24 * .. Scalar Arguments ..
26 * INTEGER INFO, LDB, LDX, N, NRHS
28 * .. Array Arguments ..
29 * DOUBLE PRECISION BERR( * ), D( * ), DF( * ), FERR( * ),
31 * COMPLEX*16 B( LDB, * ), E( * ), EF( * ), WORK( * ),
41 *> ZPTRFS improves the computed solution to a system of linear
42 *> equations when the coefficient matrix is Hermitian positive definite
43 *> and tridiagonal, and provides error bounds and backward error
44 *> estimates for the solution.
52 *> UPLO is CHARACTER*1
53 *> Specifies whether the superdiagonal or the subdiagonal of the
54 *> tridiagonal matrix A is stored and the form of the
56 *> = 'U': E is the superdiagonal of A, and A = U**H*D*U;
57 *> = 'L': E is the subdiagonal of A, and A = L*D*L**H.
58 *> (The two forms are equivalent if A is real.)
64 *> The order of the matrix A. N >= 0.
70 *> The number of right hand sides, i.e., the number of columns
71 *> of the matrix B. NRHS >= 0.
76 *> D is DOUBLE PRECISION array, dimension (N)
77 *> The n real diagonal elements of the tridiagonal matrix A.
82 *> E is COMPLEX*16 array, dimension (N-1)
83 *> The (n-1) off-diagonal elements of the tridiagonal matrix A
89 *> DF is DOUBLE PRECISION array, dimension (N)
90 *> The n diagonal elements of the diagonal matrix D from
91 *> the factorization computed by ZPTTRF.
96 *> EF is COMPLEX*16 array, dimension (N-1)
97 *> The (n-1) off-diagonal elements of the unit bidiagonal
98 *> factor U or L from the factorization computed by ZPTTRF
104 *> B is COMPLEX*16 array, dimension (LDB,NRHS)
105 *> The right hand side matrix B.
111 *> The leading dimension of the array B. LDB >= max(1,N).
116 *> X is COMPLEX*16 array, dimension (LDX,NRHS)
117 *> On entry, the solution matrix X, as computed by ZPTTRS.
118 *> On exit, the improved solution matrix X.
124 *> The leading dimension of the array X. LDX >= max(1,N).
129 *> FERR is DOUBLE PRECISION array, dimension (NRHS)
130 *> The forward error bound for each solution vector
131 *> X(j) (the j-th column of the solution matrix X).
132 *> If XTRUE is the true solution corresponding to X(j), FERR(j)
133 *> is an estimated upper bound for the magnitude of the largest
134 *> element in (X(j) - XTRUE) divided by the magnitude of the
135 *> largest element in X(j).
140 *> BERR is DOUBLE PRECISION array, dimension (NRHS)
141 *> The componentwise relative backward error of each solution
142 *> vector X(j) (i.e., the smallest relative change in
143 *> any element of A or B that makes X(j) an exact solution).
148 *> WORK is COMPLEX*16 array, dimension (N)
153 *> RWORK is DOUBLE PRECISION array, dimension (N)
159 *> = 0: successful exit
160 *> < 0: if INFO = -i, the i-th argument had an illegal value
163 *> \par Internal Parameters:
164 * =========================
167 *> ITMAX is the maximum number of steps of iterative refinement.
173 *> \author Univ. of Tennessee
174 *> \author Univ. of California Berkeley
175 *> \author Univ. of Colorado Denver
178 *> \date September 2012
180 *> \ingroup complex16PTcomputational
182 * =====================================================================
183 SUBROUTINE ZPTRFS( UPLO, N, NRHS, D, E, DF, EF, B, LDB, X, LDX,
184 $ FERR, BERR, WORK, RWORK, INFO )
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..--
191 * .. Scalar Arguments ..
193 INTEGER INFO, LDB, LDX, N, NRHS
195 * .. Array Arguments ..
196 DOUBLE PRECISION BERR( * ), D( * ), DF( * ), FERR( * ),
198 COMPLEX*16 B( LDB, * ), E( * ), EF( * ), WORK( * ),
202 * =====================================================================
206 PARAMETER ( ITMAX = 5 )
207 DOUBLE PRECISION ZERO
208 PARAMETER ( ZERO = 0.0D+0 )
210 PARAMETER ( ONE = 1.0D+0 )
212 PARAMETER ( TWO = 2.0D+0 )
213 DOUBLE PRECISION THREE
214 PARAMETER ( THREE = 3.0D+0 )
216 * .. Local Scalars ..
218 INTEGER COUNT, I, IX, J, NZ
219 DOUBLE PRECISION EPS, LSTRES, S, SAFE1, SAFE2, SAFMIN
220 COMPLEX*16 BI, CX, DX, EX, ZDUM
222 * .. External Functions ..
225 DOUBLE PRECISION DLAMCH
226 EXTERNAL LSAME, IDAMAX, DLAMCH
228 * .. External Subroutines ..
229 EXTERNAL XERBLA, ZAXPY, ZPTTRS
231 * .. Intrinsic Functions ..
232 INTRINSIC ABS, DBLE, DCMPLX, DCONJG, DIMAG, MAX
234 * .. Statement Functions ..
235 DOUBLE PRECISION CABS1
237 * .. Statement Function definitions ..
238 CABS1( ZDUM ) = ABS( DBLE( ZDUM ) ) + ABS( DIMAG( ZDUM ) )
240 * .. Executable Statements ..
242 * Test the input parameters.
245 UPPER = LSAME( UPLO, 'U' )
246 IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
248 ELSE IF( N.LT.0 ) THEN
250 ELSE IF( NRHS.LT.0 ) THEN
252 ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
254 ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
258 CALL XERBLA( 'ZPTRFS', -INFO )
262 * Quick return if possible
264 IF( N.EQ.0 .OR. NRHS.EQ.0 ) THEN
272 * NZ = maximum number of nonzero elements in each row of A, plus 1
275 EPS = DLAMCH( 'Epsilon' )
276 SAFMIN = DLAMCH( 'Safe minimum' )
280 * Do for each right hand side
288 * Loop until stopping criterion is satisfied.
290 * Compute residual R = B - A * X. Also compute
291 * abs(A)*abs(x) + abs(b) for use in the backward error bound.
296 DX = D( 1 )*X( 1, J )
298 RWORK( 1 ) = CABS1( BI ) + CABS1( DX )
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 ) )
308 CX = DCONJG( 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 ) )
318 CX = DCONJG( 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 )
327 DX = D( 1 )*X( 1, J )
329 RWORK( 1 ) = CABS1( BI ) + CABS1( DX )
332 DX = D( 1 )*X( 1, J )
333 EX = DCONJG( 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 ) )
339 CX = E( I-1 )*X( I-1, J )
340 DX = D( I )*X( I, J )
341 EX = DCONJG( 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 ) )
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 )
357 * Compute componentwise relative backward error from formula
359 * max(i) ( abs(R(i)) / ( abs(A)*abs(X) + abs(B) )(i) )
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.
368 IF( RWORK( I ).GT.SAFE2 ) THEN
369 S = MAX( S, CABS1( WORK( I ) ) / RWORK( I ) )
371 S = MAX( S, ( CABS1( WORK( I ) )+SAFE1 ) /
372 $ ( RWORK( I )+SAFE1 ) )
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.
383 IF( BERR( J ).GT.EPS .AND. TWO*BERR( J ).LE.LSTRES .AND.
384 $ COUNT.LE.ITMAX ) THEN
386 * Update solution and try again.
388 CALL ZPTTRS( UPLO, N, 1, DF, EF, WORK, N, INFO )
389 CALL ZAXPY( N, DCMPLX( ONE ), WORK, 1, X( 1, J ), 1 )
395 * Bound error from formula
397 * norm(X - XTRUE) / norm(X) .le. FERR =
399 * ( abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) ))) / norm(X)
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
406 * NZ is the maximum number of nonzeros in any row of A, plus 1
407 * EPS is machine epsilon
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.
414 IF( RWORK( I ).GT.SAFE2 ) THEN
415 RWORK( I ) = CABS1( WORK( I ) ) + NZ*EPS*RWORK( I )
417 RWORK( I ) = CABS1( WORK( I ) ) + NZ*EPS*RWORK( I ) +
421 IX = IDAMAX( N, RWORK, 1 )
422 FERR( J ) = RWORK( IX )
424 * Estimate the norm of inv(A).
426 * Solve M(A) * x = e, where M(A) = (m(i,j)) is given by
428 * m(i,j) = abs(A(i,j)), i = j,
429 * m(i,j) = -abs(A(i,j)), i .ne. j,
431 * and e = [ 1, 1, ..., 1 ]**T. Note M(A) = M(L)*D*M(L)**H.
433 * Solve M(L) * x = e.
437 RWORK( I ) = ONE + RWORK( I-1 )*ABS( EF( I-1 ) )
440 * Solve D * M(L)**H * x = b.
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 ) )
448 * Compute norm(inv(A)) = max(x(i)), 1<=i<=n.
450 IX = IDAMAX( N, RWORK, 1 )
451 FERR( J ) = FERR( J )*ABS( RWORK( IX ) )
457 LSTRES = MAX( LSTRES, ABS( X( I, J ) ) )
460 $ FERR( J ) = FERR( J ) / LSTRES