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