LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
spprfs.f
Go to the documentation of this file.
1*> \brief \b SPPRFS
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download SPPRFS + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/spprfs.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/spprfs.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/spprfs.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE SPPRFS( UPLO, N, NRHS, AP, AFP, B, LDB, X, LDX, FERR,
22* BERR, WORK, IWORK, INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER UPLO
26* INTEGER INFO, LDB, LDX, N, NRHS
27* ..
28* .. Array Arguments ..
29* INTEGER IWORK( * )
30* REAL AFP( * ), AP( * ), B( LDB, * ), BERR( * ),
31* $ FERR( * ), WORK( * ), X( LDX, * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> SPPRFS improves the computed solution to a system of linear
41*> equations when the coefficient matrix is symmetric 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 REAL array, dimension (N*(N+1)/2)
72*> The upper or lower triangle of the symmetric 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 REAL array, dimension (N*(N+1)/2)
82*> The triangular factor U or L from the Cholesky factorization
83*> A = U**T*U or A = L*L**T, as computed by SPPTRF/CPPTRF,
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 REAL 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 REAL array, dimension (LDX,NRHS)
103*> On entry, the solution matrix X, as computed by SPPTRS.
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 REAL 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 REAL 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 REAL array, dimension (3*N)
137*> \endverbatim
138*>
139*> \param[out] IWORK
140*> \verbatim
141*> IWORK is INTEGER 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 pprfs
167*
168* =====================================================================
169 SUBROUTINE spprfs( UPLO, N, NRHS, AP, AFP, B, LDB, X, LDX, FERR,
170 $ BERR, WORK, IWORK, 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 INTEGER IWORK( * )
182 REAL AFP( * ), AP( * ), B( LDB, * ), BERR( * ),
183 $ ferr( * ), work( * ), x( ldx, * )
184* ..
185*
186* =====================================================================
187*
188* .. Parameters ..
189 INTEGER ITMAX
190 parameter( itmax = 5 )
191 REAL ZERO
192 parameter( zero = 0.0e+0 )
193 REAL ONE
194 parameter( one = 1.0e+0 )
195 REAL TWO
196 parameter( two = 2.0e+0 )
197 REAL THREE
198 parameter( three = 3.0e+0 )
199* ..
200* .. Local Scalars ..
201 LOGICAL UPPER
202 INTEGER COUNT, I, IK, J, K, KASE, KK, NZ
203 REAL EPS, LSTRES, S, SAFE1, SAFE2, SAFMIN, XK
204* ..
205* .. Local Arrays ..
206 INTEGER ISAVE( 3 )
207* ..
208* .. External Subroutines ..
209 EXTERNAL saxpy, scopy, slacn2, spptrs, sspmv, xerbla
210* ..
211* .. Intrinsic Functions ..
212 INTRINSIC abs, max
213* ..
214* .. External Functions ..
215 LOGICAL LSAME
216 REAL SLAMCH
217 EXTERNAL lsame, slamch
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( 'SPPRFS', -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 = slamch( 'Epsilon' )
255 safmin = slamch( '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 scopy( n, b( 1, j ), 1, work( n+1 ), 1 )
272 CALL sspmv( uplo, n, -one, ap, x( 1, j ), 1, one, work( n+1 ),
273 $ 1 )
274*
275* Compute componentwise relative backward error from formula
276*
277* max(i) ( abs(R(i)) / ( abs(A)*abs(X) + abs(B) )(i) )
278*
279* where abs(Z) is the componentwise absolute value of the matrix
280* or vector Z. If the i-th component of the denominator is less
281* than SAFE2, then SAFE1 is added to the i-th components of the
282* numerator and denominator before dividing.
283*
284 DO 30 i = 1, n
285 work( i ) = abs( b( i, j ) )
286 30 CONTINUE
287*
288* Compute abs(A)*abs(X) + abs(B).
289*
290 kk = 1
291 IF( upper ) THEN
292 DO 50 k = 1, n
293 s = zero
294 xk = abs( x( k, j ) )
295 ik = kk
296 DO 40 i = 1, k - 1
297 work( i ) = work( i ) + abs( ap( ik ) )*xk
298 s = s + abs( ap( ik ) )*abs( x( i, j ) )
299 ik = ik + 1
300 40 CONTINUE
301 work( k ) = work( k ) + abs( ap( kk+k-1 ) )*xk + s
302 kk = kk + k
303 50 CONTINUE
304 ELSE
305 DO 70 k = 1, n
306 s = zero
307 xk = abs( x( k, j ) )
308 work( k ) = work( k ) + abs( ap( kk ) )*xk
309 ik = kk + 1
310 DO 60 i = k + 1, n
311 work( i ) = work( i ) + abs( ap( ik ) )*xk
312 s = s + abs( ap( ik ) )*abs( x( i, j ) )
313 ik = ik + 1
314 60 CONTINUE
315 work( k ) = work( k ) + s
316 kk = kk + ( n-k+1 )
317 70 CONTINUE
318 END IF
319 s = zero
320 DO 80 i = 1, n
321 IF( work( i ).GT.safe2 ) THEN
322 s = max( s, abs( work( n+i ) ) / work( i ) )
323 ELSE
324 s = max( s, ( abs( work( n+i ) )+safe1 ) /
325 $ ( work( i )+safe1 ) )
326 END IF
327 80 CONTINUE
328 berr( j ) = s
329*
330* Test stopping criterion. Continue iterating if
331* 1) The residual BERR(J) is larger than machine epsilon, and
332* 2) BERR(J) decreased by at least a factor of 2 during the
333* last iteration, and
334* 3) At most ITMAX iterations tried.
335*
336 IF( berr( j ).GT.eps .AND. two*berr( j ).LE.lstres .AND.
337 $ count.LE.itmax ) THEN
338*
339* Update solution and try again.
340*
341 CALL spptrs( uplo, n, 1, afp, work( n+1 ), n, info )
342 CALL saxpy( n, one, work( n+1 ), 1, x( 1, j ), 1 )
343 lstres = berr( j )
344 count = count + 1
345 GO TO 20
346 END IF
347*
348* Bound error from formula
349*
350* norm(X - XTRUE) / norm(X) .le. FERR =
351* norm( abs(inv(A))*
352* ( abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) ))) / norm(X)
353*
354* where
355* norm(Z) is the magnitude of the largest component of Z
356* inv(A) is the inverse of A
357* abs(Z) is the componentwise absolute value of the matrix or
358* vector Z
359* NZ is the maximum number of nonzeros in any row of A, plus 1
360* EPS is machine epsilon
361*
362* The i-th component of abs(R)+NZ*EPS*(abs(A)*abs(X)+abs(B))
363* is incremented by SAFE1 if the i-th component of
364* abs(A)*abs(X) + abs(B) is less than SAFE2.
365*
366* Use SLACN2 to estimate the infinity-norm of the matrix
367* inv(A) * diag(W),
368* where W = abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) )))
369*
370 DO 90 i = 1, n
371 IF( work( i ).GT.safe2 ) THEN
372 work( i ) = abs( work( n+i ) ) + nz*eps*work( i )
373 ELSE
374 work( i ) = abs( work( n+i ) ) + nz*eps*work( i ) + safe1
375 END IF
376 90 CONTINUE
377*
378 kase = 0
379 100 CONTINUE
380 CALL slacn2( n, work( 2*n+1 ), work( n+1 ), iwork, ferr( j ),
381 $ kase, isave )
382 IF( kase.NE.0 ) THEN
383 IF( kase.EQ.1 ) THEN
384*
385* Multiply by diag(W)*inv(A**T).
386*
387 CALL spptrs( uplo, n, 1, afp, work( n+1 ), n, info )
388 DO 110 i = 1, n
389 work( n+i ) = work( i )*work( n+i )
390 110 CONTINUE
391 ELSE IF( kase.EQ.2 ) THEN
392*
393* Multiply by inv(A)*diag(W).
394*
395 DO 120 i = 1, n
396 work( n+i ) = work( i )*work( n+i )
397 120 CONTINUE
398 CALL spptrs( uplo, n, 1, afp, work( n+1 ), n, info )
399 END IF
400 GO TO 100
401 END IF
402*
403* Normalize error.
404*
405 lstres = zero
406 DO 130 i = 1, n
407 lstres = max( lstres, abs( x( i, j ) ) )
408 130 CONTINUE
409 IF( lstres.NE.zero )
410 $ ferr( j ) = ferr( j ) / lstres
411*
412 140 CONTINUE
413*
414 RETURN
415*
416* End of SPPRFS
417*
418 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine saxpy(n, sa, sx, incx, sy, incy)
SAXPY
Definition saxpy.f:89
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine sspmv(uplo, n, alpha, ap, x, incx, beta, y, incy)
SSPMV
Definition sspmv.f:147
subroutine slacn2(n, v, x, isgn, est, kase, isave)
SLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition slacn2.f:136
subroutine spprfs(uplo, n, nrhs, ap, afp, b, ldb, x, ldx, ferr, berr, work, iwork, info)
SPPRFS
Definition spprfs.f:171
subroutine spptrs(uplo, n, nrhs, ap, b, ldb, info)
SPPTRS
Definition spptrs.f:108