LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sporfs.f
Go to the documentation of this file.
1*> \brief \b SPORFS
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download SPORFS + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sporfs.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sporfs.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sporfs.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE SPORFS( UPLO, N, NRHS, A, LDA, AF, LDAF, B, LDB, X,
22* LDX, FERR, BERR, WORK, IWORK, INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER UPLO
26* INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
27* ..
28* .. Array Arguments ..
29* INTEGER IWORK( * )
30* REAL A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
31* $ BERR( * ), FERR( * ), WORK( * ), X( LDX, * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> SPORFS improves the computed solution to a system of linear
41*> equations when the coefficient matrix is symmetric 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 REAL array, dimension (LDA,N)
72*> The symmetric 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 REAL array, dimension (LDAF,N)
90*> The triangular factor U or L from the Cholesky factorization
91*> A = U**T*U or A = L*L**T, as computed by SPOTRF.
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 REAL 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 REAL array, dimension (LDX,NRHS)
115*> On entry, the solution matrix X, as computed by SPOTRS.
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 REAL 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 REAL 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 REAL array, dimension (3*N)
149*> \endverbatim
150*>
151*> \param[out] IWORK
152*> \verbatim
153*> IWORK is INTEGER 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*> \ingroup porfs
179*
180* =====================================================================
181 SUBROUTINE sporfs( UPLO, N, NRHS, A, LDA, AF, LDAF, B, LDB, X,
182 $ LDX, FERR, BERR, WORK, IWORK, INFO )
183*
184* -- LAPACK computational routine --
185* -- LAPACK is a software package provided by Univ. of Tennessee, --
186* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
187*
188* .. Scalar Arguments ..
189 CHARACTER UPLO
190 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
191* ..
192* .. Array Arguments ..
193 INTEGER IWORK( * )
194 REAL A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
195 $ berr( * ), ferr( * ), work( * ), x( ldx, * )
196* ..
197*
198* =====================================================================
199*
200* .. Parameters ..
201 INTEGER ITMAX
202 parameter( itmax = 5 )
203 REAL ZERO
204 parameter( zero = 0.0e+0 )
205 REAL ONE
206 parameter( one = 1.0e+0 )
207 REAL TWO
208 parameter( two = 2.0e+0 )
209 REAL THREE
210 parameter( three = 3.0e+0 )
211* ..
212* .. Local Scalars ..
213 LOGICAL UPPER
214 INTEGER COUNT, I, J, K, KASE, NZ
215 REAL EPS, LSTRES, S, SAFE1, SAFE2, SAFMIN, XK
216* ..
217* .. Local Arrays ..
218 INTEGER ISAVE( 3 )
219* ..
220* .. External Subroutines ..
221 EXTERNAL saxpy, scopy, slacn2, spotrs, ssymv, xerbla
222* ..
223* .. Intrinsic Functions ..
224 INTRINSIC abs, max
225* ..
226* .. External Functions ..
227 LOGICAL LSAME
228 REAL SLAMCH
229 EXTERNAL lsame, slamch
230* ..
231* .. Executable Statements ..
232*
233* Test the input parameters.
234*
235 info = 0
236 upper = lsame( uplo, 'U' )
237 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
238 info = -1
239 ELSE IF( n.LT.0 ) THEN
240 info = -2
241 ELSE IF( nrhs.LT.0 ) THEN
242 info = -3
243 ELSE IF( lda.LT.max( 1, n ) ) THEN
244 info = -5
245 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
246 info = -7
247 ELSE IF( ldb.LT.max( 1, n ) ) THEN
248 info = -9
249 ELSE IF( ldx.LT.max( 1, n ) ) THEN
250 info = -11
251 END IF
252 IF( info.NE.0 ) THEN
253 CALL xerbla( 'SPORFS', -info )
254 RETURN
255 END IF
256*
257* Quick return if possible
258*
259 IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
260 DO 10 j = 1, nrhs
261 ferr( j ) = zero
262 berr( j ) = zero
263 10 CONTINUE
264 RETURN
265 END IF
266*
267* NZ = maximum number of nonzero elements in each row of A, plus 1
268*
269 nz = n + 1
270 eps = slamch( 'Epsilon' )
271 safmin = slamch( 'Safe minimum' )
272 safe1 = nz*safmin
273 safe2 = safe1 / eps
274*
275* Do for each right hand side
276*
277 DO 140 j = 1, nrhs
278*
279 count = 1
280 lstres = three
281 20 CONTINUE
282*
283* Loop until stopping criterion is satisfied.
284*
285* Compute residual R = B - A * X
286*
287 CALL scopy( n, b( 1, j ), 1, work( n+1 ), 1 )
288 CALL ssymv( uplo, n, -one, a, lda, x( 1, j ), 1, one,
289 $ work( n+1 ), 1 )
290*
291* Compute componentwise relative backward error from formula
292*
293* max(i) ( abs(R(i)) / ( abs(A)*abs(X) + abs(B) )(i) )
294*
295* where abs(Z) is the componentwise absolute value of the matrix
296* or vector Z. If the i-th component of the denominator is less
297* than SAFE2, then SAFE1 is added to the i-th components of the
298* numerator and denominator before dividing.
299*
300 DO 30 i = 1, n
301 work( i ) = abs( b( i, j ) )
302 30 CONTINUE
303*
304* Compute abs(A)*abs(X) + abs(B).
305*
306 IF( upper ) THEN
307 DO 50 k = 1, n
308 s = zero
309 xk = abs( x( k, j ) )
310 DO 40 i = 1, k - 1
311 work( i ) = work( i ) + abs( a( i, k ) )*xk
312 s = s + abs( a( i, k ) )*abs( x( i, j ) )
313 40 CONTINUE
314 work( k ) = work( k ) + abs( a( k, k ) )*xk + s
315 50 CONTINUE
316 ELSE
317 DO 70 k = 1, n
318 s = zero
319 xk = abs( x( k, j ) )
320 work( k ) = work( k ) + abs( a( k, k ) )*xk
321 DO 60 i = k + 1, n
322 work( i ) = work( i ) + abs( a( i, k ) )*xk
323 s = s + abs( a( i, k ) )*abs( x( i, j ) )
324 60 CONTINUE
325 work( k ) = work( k ) + s
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 spotrs( uplo, n, 1, af, ldaf, work( n+1 ), n, info )
351 CALL saxpy( n, one, work( n+1 ), 1, x( 1, j ), 1 )
352 lstres = berr( j )
353 count = count + 1
354 GO TO 20
355 END IF
356*
357* Bound error from formula
358*
359* norm(X - XTRUE) / norm(X) .le. FERR =
360* norm( abs(inv(A))*
361* ( abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) ))) / norm(X)
362*
363* where
364* norm(Z) is the magnitude of the largest component of Z
365* inv(A) is the inverse of A
366* abs(Z) is the componentwise absolute value of the matrix or
367* vector Z
368* NZ is the maximum number of nonzeros in any row of A, plus 1
369* EPS is machine epsilon
370*
371* The i-th component of abs(R)+NZ*EPS*(abs(A)*abs(X)+abs(B))
372* is incremented by SAFE1 if the i-th component of
373* abs(A)*abs(X) + abs(B) is less than SAFE2.
374*
375* Use SLACN2 to estimate the infinity-norm of the matrix
376* inv(A) * diag(W),
377* where W = abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) )))
378*
379 DO 90 i = 1, n
380 IF( work( i ).GT.safe2 ) THEN
381 work( i ) = abs( work( n+i ) ) + nz*eps*work( i )
382 ELSE
383 work( i ) = abs( work( n+i ) ) + nz*eps*work( i ) + safe1
384 END IF
385 90 CONTINUE
386*
387 kase = 0
388 100 CONTINUE
389 CALL slacn2( n, work( 2*n+1 ), work( n+1 ), iwork, ferr( j ),
390 $ kase, isave )
391 IF( kase.NE.0 ) THEN
392 IF( kase.EQ.1 ) THEN
393*
394* Multiply by diag(W)*inv(A**T).
395*
396 CALL spotrs( uplo, n, 1, af, ldaf, work( n+1 ), n, info )
397 DO 110 i = 1, n
398 work( n+i ) = work( i )*work( n+i )
399 110 CONTINUE
400 ELSE IF( kase.EQ.2 ) THEN
401*
402* Multiply by inv(A)*diag(W).
403*
404 DO 120 i = 1, n
405 work( n+i ) = work( i )*work( n+i )
406 120 CONTINUE
407 CALL spotrs( uplo, n, 1, af, ldaf, work( n+1 ), n, info )
408 END IF
409 GO TO 100
410 END IF
411*
412* Normalize error.
413*
414 lstres = zero
415 DO 130 i = 1, n
416 lstres = max( lstres, abs( x( i, j ) ) )
417 130 CONTINUE
418 IF( lstres.NE.zero )
419 $ ferr( j ) = ferr( j ) / lstres
420*
421 140 CONTINUE
422*
423 RETURN
424*
425* End of SPORFS
426*
427 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 ssymv(uplo, n, alpha, a, lda, x, incx, beta, y, incy)
SSYMV
Definition ssymv.f:152
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 sporfs(uplo, n, nrhs, a, lda, af, ldaf, b, ldb, x, ldx, ferr, berr, work, iwork, info)
SPORFS
Definition sporfs.f:183
subroutine spotrs(uplo, n, nrhs, a, lda, b, ldb, info)
SPOTRS
Definition spotrs.f:110