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