LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
zgtrfs.f
Go to the documentation of this file.
1 *> \brief \b ZGTRFS
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download ZGTRFS + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgtrfs.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgtrfs.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgtrfs.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZGTRFS( TRANS, N, NRHS, DL, D, DU, DLF, DF, DUF, DU2,
22 * IPIV, B, LDB, X, LDX, FERR, BERR, WORK, RWORK,
23 * INFO )
24 *
25 * .. Scalar Arguments ..
26 * CHARACTER TRANS
27 * INTEGER INFO, LDB, LDX, N, NRHS
28 * ..
29 * .. Array Arguments ..
30 * INTEGER IPIV( * )
31 * DOUBLE PRECISION BERR( * ), FERR( * ), RWORK( * )
32 * COMPLEX*16 B( LDB, * ), D( * ), DF( * ), DL( * ),
33 * $ DLF( * ), DU( * ), DU2( * ), DUF( * ),
34 * $ WORK( * ), X( LDX, * )
35 * ..
36 *
37 *
38 *> \par Purpose:
39 * =============
40 *>
41 *> \verbatim
42 *>
43 *> ZGTRFS improves the computed solution to a system of linear
44 *> equations when the coefficient matrix is tridiagonal, and provides
45 *> error bounds and backward error estimates for the solution.
46 *> \endverbatim
47 *
48 * Arguments:
49 * ==========
50 *
51 *> \param[in] TRANS
52 *> \verbatim
53 *> TRANS is CHARACTER*1
54 *> Specifies the form of the system of equations:
55 *> = 'N': A * X = B (No transpose)
56 *> = 'T': A**T * X = B (Transpose)
57 *> = 'C': A**H * X = B (Conjugate transpose)
58 *> \endverbatim
59 *>
60 *> \param[in] N
61 *> \verbatim
62 *> N is INTEGER
63 *> The order of the matrix A. N >= 0.
64 *> \endverbatim
65 *>
66 *> \param[in] NRHS
67 *> \verbatim
68 *> NRHS is INTEGER
69 *> The number of right hand sides, i.e., the number of columns
70 *> of the matrix B. NRHS >= 0.
71 *> \endverbatim
72 *>
73 *> \param[in] DL
74 *> \verbatim
75 *> DL is COMPLEX*16 array, dimension (N-1)
76 *> The (n-1) subdiagonal elements of A.
77 *> \endverbatim
78 *>
79 *> \param[in] D
80 *> \verbatim
81 *> D is COMPLEX*16 array, dimension (N)
82 *> The diagonal elements of A.
83 *> \endverbatim
84 *>
85 *> \param[in] DU
86 *> \verbatim
87 *> DU is COMPLEX*16 array, dimension (N-1)
88 *> The (n-1) superdiagonal elements of A.
89 *> \endverbatim
90 *>
91 *> \param[in] DLF
92 *> \verbatim
93 *> DLF is COMPLEX*16 array, dimension (N-1)
94 *> The (n-1) multipliers that define the matrix L from the
95 *> LU factorization of A as computed by ZGTTRF.
96 *> \endverbatim
97 *>
98 *> \param[in] DF
99 *> \verbatim
100 *> DF is COMPLEX*16 array, dimension (N)
101 *> The n diagonal elements of the upper triangular matrix U from
102 *> the LU factorization of A.
103 *> \endverbatim
104 *>
105 *> \param[in] DUF
106 *> \verbatim
107 *> DUF is COMPLEX*16 array, dimension (N-1)
108 *> The (n-1) elements of the first superdiagonal of U.
109 *> \endverbatim
110 *>
111 *> \param[in] DU2
112 *> \verbatim
113 *> DU2 is COMPLEX*16 array, dimension (N-2)
114 *> The (n-2) elements of the second superdiagonal of U.
115 *> \endverbatim
116 *>
117 *> \param[in] IPIV
118 *> \verbatim
119 *> IPIV is INTEGER array, dimension (N)
120 *> The pivot indices; for 1 <= i <= n, row i of the matrix was
121 *> interchanged with row IPIV(i). IPIV(i) will always be either
122 *> i or i+1; IPIV(i) = i indicates a row interchange was not
123 *> required.
124 *> \endverbatim
125 *>
126 *> \param[in] B
127 *> \verbatim
128 *> B is COMPLEX*16 array, dimension (LDB,NRHS)
129 *> The right hand side matrix B.
130 *> \endverbatim
131 *>
132 *> \param[in] LDB
133 *> \verbatim
134 *> LDB is INTEGER
135 *> The leading dimension of the array B. LDB >= max(1,N).
136 *> \endverbatim
137 *>
138 *> \param[in,out] X
139 *> \verbatim
140 *> X is COMPLEX*16 array, dimension (LDX,NRHS)
141 *> On entry, the solution matrix X, as computed by ZGTTRS.
142 *> On exit, the improved solution matrix X.
143 *> \endverbatim
144 *>
145 *> \param[in] LDX
146 *> \verbatim
147 *> LDX is INTEGER
148 *> The leading dimension of the array X. LDX >= max(1,N).
149 *> \endverbatim
150 *>
151 *> \param[out] FERR
152 *> \verbatim
153 *> FERR is DOUBLE PRECISION array, dimension (NRHS)
154 *> The estimated forward error bound for each solution vector
155 *> X(j) (the j-th column of the solution matrix X).
156 *> If XTRUE is the true solution corresponding to X(j), FERR(j)
157 *> is an estimated upper bound for the magnitude of the largest
158 *> element in (X(j) - XTRUE) divided by the magnitude of the
159 *> largest element in X(j). The estimate is as reliable as
160 *> the estimate for RCOND, and is almost always a slight
161 *> overestimate of the true error.
162 *> \endverbatim
163 *>
164 *> \param[out] BERR
165 *> \verbatim
166 *> BERR is DOUBLE PRECISION array, dimension (NRHS)
167 *> The componentwise relative backward error of each solution
168 *> vector X(j) (i.e., the smallest relative change in
169 *> any element of A or B that makes X(j) an exact solution).
170 *> \endverbatim
171 *>
172 *> \param[out] WORK
173 *> \verbatim
174 *> WORK is COMPLEX*16 array, dimension (2*N)
175 *> \endverbatim
176 *>
177 *> \param[out] RWORK
178 *> \verbatim
179 *> RWORK is DOUBLE PRECISION array, dimension (N)
180 *> \endverbatim
181 *>
182 *> \param[out] INFO
183 *> \verbatim
184 *> INFO is INTEGER
185 *> = 0: successful exit
186 *> < 0: if INFO = -i, the i-th argument had an illegal value
187 *> \endverbatim
188 *
189 *> \par Internal Parameters:
190 * =========================
191 *>
192 *> \verbatim
193 *> ITMAX is the maximum number of steps of iterative refinement.
194 *> \endverbatim
195 *
196 * Authors:
197 * ========
198 *
199 *> \author Univ. of Tennessee
200 *> \author Univ. of California Berkeley
201 *> \author Univ. of Colorado Denver
202 *> \author NAG Ltd.
203 *
204 *> \date September 2012
205 *
206 *> \ingroup complex16GTcomputational
207 *
208 * =====================================================================
209  SUBROUTINE zgtrfs( TRANS, N, NRHS, DL, D, DU, DLF, DF, DUF, DU2,
210  $ ipiv, b, ldb, x, ldx, ferr, berr, work, rwork,
211  $ info )
212 *
213 * -- LAPACK computational routine (version 3.4.2) --
214 * -- LAPACK is a software package provided by Univ. of Tennessee, --
215 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
216 * September 2012
217 *
218 * .. Scalar Arguments ..
219  CHARACTER TRANS
220  INTEGER INFO, LDB, LDX, N, NRHS
221 * ..
222 * .. Array Arguments ..
223  INTEGER IPIV( * )
224  DOUBLE PRECISION BERR( * ), FERR( * ), RWORK( * )
225  COMPLEX*16 B( ldb, * ), D( * ), DF( * ), DL( * ),
226  $ dlf( * ), du( * ), du2( * ), duf( * ),
227  $ work( * ), x( ldx, * )
228 * ..
229 *
230 * =====================================================================
231 *
232 * .. Parameters ..
233  INTEGER ITMAX
234  parameter ( itmax = 5 )
235  DOUBLE PRECISION ZERO, ONE
236  parameter ( zero = 0.0d+0, one = 1.0d+0 )
237  DOUBLE PRECISION TWO
238  parameter ( two = 2.0d+0 )
239  DOUBLE PRECISION THREE
240  parameter ( three = 3.0d+0 )
241 * ..
242 * .. Local Scalars ..
243  LOGICAL NOTRAN
244  CHARACTER TRANSN, TRANST
245  INTEGER COUNT, I, J, KASE, NZ
246  DOUBLE PRECISION EPS, LSTRES, S, SAFE1, SAFE2, SAFMIN
247  COMPLEX*16 ZDUM
248 * ..
249 * .. Local Arrays ..
250  INTEGER ISAVE( 3 )
251 * ..
252 * .. External Subroutines ..
253  EXTERNAL xerbla, zaxpy, zcopy, zgttrs, zlacn2, zlagtm
254 * ..
255 * .. Intrinsic Functions ..
256  INTRINSIC abs, dble, dcmplx, dimag, max
257 * ..
258 * .. External Functions ..
259  LOGICAL LSAME
260  DOUBLE PRECISION DLAMCH
261  EXTERNAL lsame, dlamch
262 * ..
263 * .. Statement Functions ..
264  DOUBLE PRECISION CABS1
265 * ..
266 * .. Statement Function definitions ..
267  cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
268 * ..
269 * .. Executable Statements ..
270 *
271 * Test the input parameters.
272 *
273  info = 0
274  notran = lsame( trans, 'N' )
275  IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
276  $ lsame( trans, 'C' ) ) THEN
277  info = -1
278  ELSE IF( n.LT.0 ) THEN
279  info = -2
280  ELSE IF( nrhs.LT.0 ) THEN
281  info = -3
282  ELSE IF( ldb.LT.max( 1, n ) ) THEN
283  info = -13
284  ELSE IF( ldx.LT.max( 1, n ) ) THEN
285  info = -15
286  END IF
287  IF( info.NE.0 ) THEN
288  CALL xerbla( 'ZGTRFS', -info )
289  RETURN
290  END IF
291 *
292 * Quick return if possible
293 *
294  IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
295  DO 10 j = 1, nrhs
296  ferr( j ) = zero
297  berr( j ) = zero
298  10 CONTINUE
299  RETURN
300  END IF
301 *
302  IF( notran ) THEN
303  transn = 'N'
304  transt = 'C'
305  ELSE
306  transn = 'C'
307  transt = 'N'
308  END IF
309 *
310 * NZ = maximum number of nonzero elements in each row of A, plus 1
311 *
312  nz = 4
313  eps = dlamch( 'Epsilon' )
314  safmin = dlamch( 'Safe minimum' )
315  safe1 = nz*safmin
316  safe2 = safe1 / eps
317 *
318 * Do for each right hand side
319 *
320  DO 110 j = 1, nrhs
321 *
322  count = 1
323  lstres = three
324  20 CONTINUE
325 *
326 * Loop until stopping criterion is satisfied.
327 *
328 * Compute residual R = B - op(A) * X,
329 * where op(A) = A, A**T, or A**H, depending on TRANS.
330 *
331  CALL zcopy( n, b( 1, j ), 1, work, 1 )
332  CALL zlagtm( trans, n, 1, -one, dl, d, du, x( 1, j ), ldx, one,
333  $ work, n )
334 *
335 * Compute abs(op(A))*abs(x) + abs(b) for use in the backward
336 * error bound.
337 *
338  IF( notran ) THEN
339  IF( n.EQ.1 ) THEN
340  rwork( 1 ) = cabs1( b( 1, j ) ) +
341  $ cabs1( d( 1 ) )*cabs1( x( 1, j ) )
342  ELSE
343  rwork( 1 ) = cabs1( b( 1, j ) ) +
344  $ cabs1( d( 1 ) )*cabs1( x( 1, j ) ) +
345  $ cabs1( du( 1 ) )*cabs1( x( 2, j ) )
346  DO 30 i = 2, n - 1
347  rwork( i ) = cabs1( b( i, j ) ) +
348  $ cabs1( dl( i-1 ) )*cabs1( x( i-1, j ) ) +
349  $ cabs1( d( i ) )*cabs1( x( i, j ) ) +
350  $ cabs1( du( i ) )*cabs1( x( i+1, j ) )
351  30 CONTINUE
352  rwork( n ) = cabs1( b( n, j ) ) +
353  $ cabs1( dl( n-1 ) )*cabs1( x( n-1, j ) ) +
354  $ cabs1( d( n ) )*cabs1( x( n, j ) )
355  END IF
356  ELSE
357  IF( n.EQ.1 ) THEN
358  rwork( 1 ) = cabs1( b( 1, j ) ) +
359  $ cabs1( d( 1 ) )*cabs1( x( 1, j ) )
360  ELSE
361  rwork( 1 ) = cabs1( b( 1, j ) ) +
362  $ cabs1( d( 1 ) )*cabs1( x( 1, j ) ) +
363  $ cabs1( dl( 1 ) )*cabs1( x( 2, j ) )
364  DO 40 i = 2, n - 1
365  rwork( i ) = cabs1( b( i, j ) ) +
366  $ cabs1( du( i-1 ) )*cabs1( x( i-1, j ) ) +
367  $ cabs1( d( i ) )*cabs1( x( i, j ) ) +
368  $ cabs1( dl( i ) )*cabs1( x( i+1, j ) )
369  40 CONTINUE
370  rwork( n ) = cabs1( b( n, j ) ) +
371  $ cabs1( du( n-1 ) )*cabs1( x( n-1, j ) ) +
372  $ cabs1( d( n ) )*cabs1( x( n, j ) )
373  END IF
374  END IF
375 *
376 * Compute componentwise relative backward error from formula
377 *
378 * max(i) ( abs(R(i)) / ( abs(op(A))*abs(X) + abs(B) )(i) )
379 *
380 * where abs(Z) is the componentwise absolute value of the matrix
381 * or vector Z. If the i-th component of the denominator is less
382 * than SAFE2, then SAFE1 is added to the i-th components of the
383 * numerator and denominator before dividing.
384 *
385  s = zero
386  DO 50 i = 1, n
387  IF( rwork( i ).GT.safe2 ) THEN
388  s = max( s, cabs1( work( i ) ) / rwork( i ) )
389  ELSE
390  s = max( s, ( cabs1( work( i ) )+safe1 ) /
391  $ ( rwork( i )+safe1 ) )
392  END IF
393  50 CONTINUE
394  berr( j ) = s
395 *
396 * Test stopping criterion. Continue iterating if
397 * 1) The residual BERR(J) is larger than machine epsilon, and
398 * 2) BERR(J) decreased by at least a factor of 2 during the
399 * last iteration, and
400 * 3) At most ITMAX iterations tried.
401 *
402  IF( berr( j ).GT.eps .AND. two*berr( j ).LE.lstres .AND.
403  $ count.LE.itmax ) THEN
404 *
405 * Update solution and try again.
406 *
407  CALL zgttrs( trans, n, 1, dlf, df, duf, du2, ipiv, work, n,
408  $ info )
409  CALL zaxpy( n, dcmplx( one ), work, 1, x( 1, j ), 1 )
410  lstres = berr( j )
411  count = count + 1
412  GO TO 20
413  END IF
414 *
415 * Bound error from formula
416 *
417 * norm(X - XTRUE) / norm(X) .le. FERR =
418 * norm( abs(inv(op(A)))*
419 * ( abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) ))) / norm(X)
420 *
421 * where
422 * norm(Z) is the magnitude of the largest component of Z
423 * inv(op(A)) is the inverse of op(A)
424 * abs(Z) is the componentwise absolute value of the matrix or
425 * vector Z
426 * NZ is the maximum number of nonzeros in any row of A, plus 1
427 * EPS is machine epsilon
428 *
429 * The i-th component of abs(R)+NZ*EPS*(abs(op(A))*abs(X)+abs(B))
430 * is incremented by SAFE1 if the i-th component of
431 * abs(op(A))*abs(X) + abs(B) is less than SAFE2.
432 *
433 * Use ZLACN2 to estimate the infinity-norm of the matrix
434 * inv(op(A)) * diag(W),
435 * where W = abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) )))
436 *
437  DO 60 i = 1, n
438  IF( rwork( i ).GT.safe2 ) THEN
439  rwork( i ) = cabs1( work( i ) ) + nz*eps*rwork( i )
440  ELSE
441  rwork( i ) = cabs1( work( i ) ) + nz*eps*rwork( i ) +
442  $ safe1
443  END IF
444  60 CONTINUE
445 *
446  kase = 0
447  70 CONTINUE
448  CALL zlacn2( n, work( n+1 ), work, ferr( j ), kase, isave )
449  IF( kase.NE.0 ) THEN
450  IF( kase.EQ.1 ) THEN
451 *
452 * Multiply by diag(W)*inv(op(A)**H).
453 *
454  CALL zgttrs( transt, n, 1, dlf, df, duf, du2, ipiv, work,
455  $ n, info )
456  DO 80 i = 1, n
457  work( i ) = rwork( i )*work( i )
458  80 CONTINUE
459  ELSE
460 *
461 * Multiply by inv(op(A))*diag(W).
462 *
463  DO 90 i = 1, n
464  work( i ) = rwork( i )*work( i )
465  90 CONTINUE
466  CALL zgttrs( transn, n, 1, dlf, df, duf, du2, ipiv, work,
467  $ n, info )
468  END IF
469  GO TO 70
470  END IF
471 *
472 * Normalize error.
473 *
474  lstres = zero
475  DO 100 i = 1, n
476  lstres = max( lstres, cabs1( x( i, j ) ) )
477  100 CONTINUE
478  IF( lstres.NE.zero )
479  $ ferr( j ) = ferr( j ) / lstres
480 *
481  110 CONTINUE
482 *
483  RETURN
484 *
485 * End of ZGTRFS
486 *
487  END
subroutine zlagtm(TRANS, N, NRHS, ALPHA, DL, D, DU, X, LDX, BETA, B, LDB)
ZLAGTM performs a matrix-matrix product of the form C = αAB+βC, where A is a tridiagonal matrix...
Definition: zlagtm.f:147
subroutine zgttrs(TRANS, N, NRHS, DL, D, DU, DU2, IPIV, B, LDB, INFO)
ZGTTRS
Definition: zgttrs.f:140
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zlacn2(N, V, X, EST, KASE, ISAVE)
ZLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition: zlacn2.f:135
subroutine zaxpy(N, ZA, ZX, INCX, ZY, INCY)
ZAXPY
Definition: zaxpy.f:53
subroutine zgtrfs(TRANS, N, NRHS, DL, D, DU, DLF, DF, DUF, DU2, IPIV, B, LDB, X, LDX, FERR, BERR, WORK, RWORK, INFO)
ZGTRFS
Definition: zgtrfs.f:212