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