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