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