LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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 *> \date November 2011
167 *
168 *> \ingroup realOTHERcomputational
169 *
170 * =====================================================================
171  SUBROUTINE spprfs( UPLO, N, NRHS, AP, AFP, B, LDB, X, LDX, FERR,
172  $ berr, work, iwork, info )
173 *
174 * -- LAPACK computational routine (version 3.4.0) --
175 * -- LAPACK is a software package provided by Univ. of Tennessee, --
176 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
177 * November 2011
178 *
179 * .. Scalar Arguments ..
180  CHARACTER UPLO
181  INTEGER INFO, LDB, LDX, N, NRHS
182 * ..
183 * .. Array Arguments ..
184  INTEGER IWORK( * )
185  REAL AFP( * ), AP( * ), B( ldb, * ), BERR( * ),
186  $ ferr( * ), work( * ), x( ldx, * )
187 * ..
188 *
189 * =====================================================================
190 *
191 * .. Parameters ..
192  INTEGER ITMAX
193  parameter ( itmax = 5 )
194  REAL ZERO
195  parameter ( zero = 0.0e+0 )
196  REAL ONE
197  parameter ( one = 1.0e+0 )
198  REAL TWO
199  parameter ( two = 2.0e+0 )
200  REAL THREE
201  parameter ( three = 3.0e+0 )
202 * ..
203 * .. Local Scalars ..
204  LOGICAL UPPER
205  INTEGER COUNT, I, IK, J, K, KASE, KK, NZ
206  REAL EPS, LSTRES, S, SAFE1, SAFE2, SAFMIN, XK
207 * ..
208 * .. Local Arrays ..
209  INTEGER ISAVE( 3 )
210 * ..
211 * .. External Subroutines ..
212  EXTERNAL saxpy, scopy, slacn2, spptrs, sspmv, xerbla
213 * ..
214 * .. Intrinsic Functions ..
215  INTRINSIC abs, max
216 * ..
217 * .. External Functions ..
218  LOGICAL LSAME
219  REAL SLAMCH
220  EXTERNAL lsame, slamch
221 * ..
222 * .. Executable Statements ..
223 *
224 * Test the input parameters.
225 *
226  info = 0
227  upper = lsame( uplo, 'U' )
228  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
229  info = -1
230  ELSE IF( n.LT.0 ) THEN
231  info = -2
232  ELSE IF( nrhs.LT.0 ) THEN
233  info = -3
234  ELSE IF( ldb.LT.max( 1, n ) ) THEN
235  info = -7
236  ELSE IF( ldx.LT.max( 1, n ) ) THEN
237  info = -9
238  END IF
239  IF( info.NE.0 ) THEN
240  CALL xerbla( 'SPPRFS', -info )
241  RETURN
242  END IF
243 *
244 * Quick return if possible
245 *
246  IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
247  DO 10 j = 1, nrhs
248  ferr( j ) = zero
249  berr( j ) = zero
250  10 CONTINUE
251  RETURN
252  END IF
253 *
254 * NZ = maximum number of nonzero elements in each row of A, plus 1
255 *
256  nz = n + 1
257  eps = slamch( 'Epsilon' )
258  safmin = slamch( 'Safe minimum' )
259  safe1 = nz*safmin
260  safe2 = safe1 / eps
261 *
262 * Do for each right hand side
263 *
264  DO 140 j = 1, nrhs
265 *
266  count = 1
267  lstres = three
268  20 CONTINUE
269 *
270 * Loop until stopping criterion is satisfied.
271 *
272 * Compute residual R = B - A * X
273 *
274  CALL scopy( n, b( 1, j ), 1, work( n+1 ), 1 )
275  CALL sspmv( uplo, n, -one, ap, x( 1, j ), 1, one, work( n+1 ),
276  $ 1 )
277 *
278 * Compute componentwise relative backward error from formula
279 *
280 * max(i) ( abs(R(i)) / ( abs(A)*abs(X) + abs(B) )(i) )
281 *
282 * where abs(Z) is the componentwise absolute value of the matrix
283 * or vector Z. If the i-th component of the denominator is less
284 * than SAFE2, then SAFE1 is added to the i-th components of the
285 * numerator and denominator before dividing.
286 *
287  DO 30 i = 1, n
288  work( i ) = abs( b( i, j ) )
289  30 CONTINUE
290 *
291 * Compute abs(A)*abs(X) + abs(B).
292 *
293  kk = 1
294  IF( upper ) THEN
295  DO 50 k = 1, n
296  s = zero
297  xk = abs( x( k, j ) )
298  ik = kk
299  DO 40 i = 1, k - 1
300  work( i ) = work( i ) + abs( ap( ik ) )*xk
301  s = s + abs( ap( ik ) )*abs( x( i, j ) )
302  ik = ik + 1
303  40 CONTINUE
304  work( k ) = work( k ) + abs( ap( kk+k-1 ) )*xk + s
305  kk = kk + k
306  50 CONTINUE
307  ELSE
308  DO 70 k = 1, n
309  s = zero
310  xk = abs( x( k, j ) )
311  work( k ) = work( k ) + abs( ap( kk ) )*xk
312  ik = kk + 1
313  DO 60 i = k + 1, n
314  work( i ) = work( i ) + abs( ap( ik ) )*xk
315  s = s + abs( ap( ik ) )*abs( x( i, j ) )
316  ik = ik + 1
317  60 CONTINUE
318  work( k ) = work( k ) + s
319  kk = kk + ( n-k+1 )
320  70 CONTINUE
321  END IF
322  s = zero
323  DO 80 i = 1, n
324  IF( work( i ).GT.safe2 ) THEN
325  s = max( s, abs( work( n+i ) ) / work( i ) )
326  ELSE
327  s = max( s, ( abs( work( n+i ) )+safe1 ) /
328  $ ( work( i )+safe1 ) )
329  END IF
330  80 CONTINUE
331  berr( j ) = s
332 *
333 * Test stopping criterion. Continue iterating if
334 * 1) The residual BERR(J) is larger than machine epsilon, and
335 * 2) BERR(J) decreased by at least a factor of 2 during the
336 * last iteration, and
337 * 3) At most ITMAX iterations tried.
338 *
339  IF( berr( j ).GT.eps .AND. two*berr( j ).LE.lstres .AND.
340  $ count.LE.itmax ) THEN
341 *
342 * Update solution and try again.
343 *
344  CALL spptrs( uplo, n, 1, afp, work( n+1 ), n, info )
345  CALL saxpy( n, one, work( n+1 ), 1, x( 1, j ), 1 )
346  lstres = berr( j )
347  count = count + 1
348  GO TO 20
349  END IF
350 *
351 * Bound error from formula
352 *
353 * norm(X - XTRUE) / norm(X) .le. FERR =
354 * norm( abs(inv(A))*
355 * ( abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) ))) / norm(X)
356 *
357 * where
358 * norm(Z) is the magnitude of the largest component of Z
359 * inv(A) is the inverse of A
360 * abs(Z) is the componentwise absolute value of the matrix or
361 * vector Z
362 * NZ is the maximum number of nonzeros in any row of A, plus 1
363 * EPS is machine epsilon
364 *
365 * The i-th component of abs(R)+NZ*EPS*(abs(A)*abs(X)+abs(B))
366 * is incremented by SAFE1 if the i-th component of
367 * abs(A)*abs(X) + abs(B) is less than SAFE2.
368 *
369 * Use SLACN2 to estimate the infinity-norm of the matrix
370 * inv(A) * diag(W),
371 * where W = abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) )))
372 *
373  DO 90 i = 1, n
374  IF( work( i ).GT.safe2 ) THEN
375  work( i ) = abs( work( n+i ) ) + nz*eps*work( i )
376  ELSE
377  work( i ) = abs( work( n+i ) ) + nz*eps*work( i ) + safe1
378  END IF
379  90 CONTINUE
380 *
381  kase = 0
382  100 CONTINUE
383  CALL slacn2( n, work( 2*n+1 ), work( n+1 ), iwork, ferr( j ),
384  $ kase, isave )
385  IF( kase.NE.0 ) THEN
386  IF( kase.EQ.1 ) THEN
387 *
388 * Multiply by diag(W)*inv(A**T).
389 *
390  CALL spptrs( uplo, n, 1, afp, work( n+1 ), n, info )
391  DO 110 i = 1, n
392  work( n+i ) = work( i )*work( n+i )
393  110 CONTINUE
394  ELSE IF( kase.EQ.2 ) THEN
395 *
396 * Multiply by inv(A)*diag(W).
397 *
398  DO 120 i = 1, n
399  work( n+i ) = work( i )*work( n+i )
400  120 CONTINUE
401  CALL spptrs( uplo, n, 1, afp, work( n+1 ), n, info )
402  END IF
403  GO TO 100
404  END IF
405 *
406 * Normalize error.
407 *
408  lstres = zero
409  DO 130 i = 1, n
410  lstres = max( lstres, abs( x( i, j ) ) )
411  130 CONTINUE
412  IF( lstres.NE.zero )
413  $ ferr( j ) = ferr( j ) / lstres
414 *
415  140 CONTINUE
416 *
417  RETURN
418 *
419 * End of SPPRFS
420 *
421  END
subroutine sspmv(UPLO, N, ALPHA, AP, X, INCX, BETA, Y, INCY)
SSPMV
Definition: sspmv.f:149
subroutine spptrs(UPLO, N, NRHS, AP, B, LDB, INFO)
SPPTRS
Definition: spptrs.f:110
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
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:138
subroutine saxpy(N, SA, SX, INCX, SY, INCY)
SAXPY
Definition: saxpy.f:54
subroutine spprfs(UPLO, N, NRHS, AP, AFP, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, INFO)
SPPRFS
Definition: spprfs.f:173
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:53