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