LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
stprfs.f
Go to the documentation of this file.
1*> \brief \b STPRFS
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download STPRFS + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/stprfs.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/stprfs.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/stprfs.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE STPRFS( UPLO, TRANS, DIAG, N, NRHS, AP, B, LDB, X, LDX,
20* FERR, BERR, WORK, IWORK, INFO )
21*
22* .. Scalar Arguments ..
23* CHARACTER DIAG, TRANS, UPLO
24* INTEGER INFO, LDB, LDX, N, NRHS
25* ..
26* .. Array Arguments ..
27* INTEGER IWORK( * )
28* REAL AP( * ), B( LDB, * ), BERR( * ), FERR( * ),
29* $ WORK( * ), X( LDX, * )
30* ..
31*
32*
33*> \par Purpose:
34* =============
35*>
36*> \verbatim
37*>
38*> STPRFS provides error bounds and backward error estimates for the
39*> solution to a system of linear equations with a triangular packed
40*> coefficient matrix.
41*>
42*> The solution matrix X must be computed by STPTRS or some other
43*> means before entering this routine. STPRFS does not do iterative
44*> refinement because doing so cannot improve the backward error.
45*> \endverbatim
46*
47* Arguments:
48* ==========
49*
50*> \param[in] UPLO
51*> \verbatim
52*> UPLO is CHARACTER*1
53*> = 'U': A is upper triangular;
54*> = 'L': A is lower triangular.
55*> \endverbatim
56*>
57*> \param[in] TRANS
58*> \verbatim
59*> TRANS is CHARACTER*1
60*> Specifies the form of the system of equations:
61*> = 'N': A * X = B (No transpose)
62*> = 'T': A**T * X = B (Transpose)
63*> = 'C': A**H * X = B (Conjugate transpose = Transpose)
64*> \endverbatim
65*>
66*> \param[in] DIAG
67*> \verbatim
68*> DIAG is CHARACTER*1
69*> = 'N': A is non-unit triangular;
70*> = 'U': A is unit triangular.
71*> \endverbatim
72*>
73*> \param[in] N
74*> \verbatim
75*> N is INTEGER
76*> The order of the matrix A. N >= 0.
77*> \endverbatim
78*>
79*> \param[in] NRHS
80*> \verbatim
81*> NRHS is INTEGER
82*> The number of right hand sides, i.e., the number of columns
83*> of the matrices B and X. NRHS >= 0.
84*> \endverbatim
85*>
86*> \param[in] AP
87*> \verbatim
88*> AP is REAL array, dimension (N*(N+1)/2)
89*> The upper or lower triangular matrix A, packed columnwise in
90*> a linear array. The j-th column of A is stored in the array
91*> AP as follows:
92*> if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
93*> if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.
94*> If DIAG = 'U', the diagonal elements of A are not referenced
95*> and are assumed to be 1.
96*> \endverbatim
97*>
98*> \param[in] B
99*> \verbatim
100*> B is REAL array, dimension (LDB,NRHS)
101*> The right hand side matrix B.
102*> \endverbatim
103*>
104*> \param[in] LDB
105*> \verbatim
106*> LDB is INTEGER
107*> The leading dimension of the array B. LDB >= max(1,N).
108*> \endverbatim
109*>
110*> \param[in] X
111*> \verbatim
112*> X is REAL array, dimension (LDX,NRHS)
113*> The solution matrix X.
114*> \endverbatim
115*>
116*> \param[in] LDX
117*> \verbatim
118*> LDX is INTEGER
119*> The leading dimension of the array X. LDX >= max(1,N).
120*> \endverbatim
121*>
122*> \param[out] FERR
123*> \verbatim
124*> FERR is REAL array, dimension (NRHS)
125*> The estimated forward error bound for each solution vector
126*> X(j) (the j-th column of the solution matrix X).
127*> If XTRUE is the true solution corresponding to X(j), FERR(j)
128*> is an estimated upper bound for the magnitude of the largest
129*> element in (X(j) - XTRUE) divided by the magnitude of the
130*> largest element in X(j). The estimate is as reliable as
131*> the estimate for RCOND, and is almost always a slight
132*> overestimate of the true error.
133*> \endverbatim
134*>
135*> \param[out] BERR
136*> \verbatim
137*> BERR is REAL array, dimension (NRHS)
138*> The componentwise relative backward error of each solution
139*> vector X(j) (i.e., the smallest relative change in
140*> any element of A or B that makes X(j) an exact solution).
141*> \endverbatim
142*>
143*> \param[out] WORK
144*> \verbatim
145*> WORK is REAL array, dimension (3*N)
146*> \endverbatim
147*>
148*> \param[out] IWORK
149*> \verbatim
150*> IWORK is INTEGER array, dimension (N)
151*> \endverbatim
152*>
153*> \param[out] INFO
154*> \verbatim
155*> INFO is INTEGER
156*> = 0: successful exit
157*> < 0: if INFO = -i, the i-th argument had an illegal value
158*> \endverbatim
159*
160* Authors:
161* ========
162*
163*> \author Univ. of Tennessee
164*> \author Univ. of California Berkeley
165*> \author Univ. of Colorado Denver
166*> \author NAG Ltd.
167*
168*> \ingroup tprfs
169*
170* =====================================================================
171 SUBROUTINE stprfs( UPLO, TRANS, DIAG, N, NRHS, AP, B, LDB, X,
172 $ LDX,
173 $ FERR, BERR, WORK, IWORK, INFO )
174*
175* -- LAPACK computational routine --
176* -- LAPACK is a software package provided by Univ. of Tennessee, --
177* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
178*
179* .. Scalar Arguments ..
180 CHARACTER DIAG, TRANS, UPLO
181 INTEGER INFO, LDB, LDX, N, NRHS
182* ..
183* .. Array Arguments ..
184 INTEGER IWORK( * )
185 REAL AP( * ), B( LDB, * ), BERR( * ), FERR( * ),
186 $ work( * ), x( ldx, * )
187* ..
188*
189* =====================================================================
190*
191* .. Parameters ..
192 REAL ZERO
193 PARAMETER ( ZERO = 0.0e+0 )
194 REAL ONE
195 parameter( one = 1.0e+0 )
196* ..
197* .. Local Scalars ..
198 LOGICAL NOTRAN, NOUNIT, UPPER
199 CHARACTER TRANST
200 INTEGER I, J, K, KASE, KC, NZ
201 REAL EPS, LSTRES, S, SAFE1, SAFE2, SAFMIN, XK
202* ..
203* .. Local Arrays ..
204 INTEGER ISAVE( 3 )
205* ..
206* .. External Subroutines ..
207 EXTERNAL saxpy, scopy, slacn2, stpmv, stpsv,
208 $ xerbla
209* ..
210* .. Intrinsic Functions ..
211 INTRINSIC abs, max
212* ..
213* .. External Functions ..
214 LOGICAL LSAME
215 REAL SLAMCH
216 EXTERNAL lsame, slamch
217* ..
218* .. Executable Statements ..
219*
220* Test the input parameters.
221*
222 info = 0
223 upper = lsame( uplo, 'U' )
224 notran = lsame( trans, 'N' )
225 nounit = lsame( diag, 'N' )
226*
227 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
228 info = -1
229 ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
230 $ lsame( trans, 'C' ) ) THEN
231 info = -2
232 ELSE IF( .NOT.nounit .AND. .NOT.lsame( diag, 'U' ) ) THEN
233 info = -3
234 ELSE IF( n.LT.0 ) THEN
235 info = -4
236 ELSE IF( nrhs.LT.0 ) THEN
237 info = -5
238 ELSE IF( ldb.LT.max( 1, n ) ) THEN
239 info = -8
240 ELSE IF( ldx.LT.max( 1, n ) ) THEN
241 info = -10
242 END IF
243 IF( info.NE.0 ) THEN
244 CALL xerbla( 'STPRFS', -info )
245 RETURN
246 END IF
247*
248* Quick return if possible
249*
250 IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
251 DO 10 j = 1, nrhs
252 ferr( j ) = zero
253 berr( j ) = zero
254 10 CONTINUE
255 RETURN
256 END IF
257*
258 IF( notran ) THEN
259 transt = 'T'
260 ELSE
261 transt = 'N'
262 END IF
263*
264* NZ = maximum number of nonzero elements in each row of A, plus 1
265*
266 nz = n + 1
267 eps = slamch( 'Epsilon' )
268 safmin = slamch( 'Safe minimum' )
269 safe1 = real( nz )*safmin
270 safe2 = safe1 / eps
271*
272* Do for each right hand side
273*
274 DO 250 j = 1, nrhs
275*
276* Compute residual R = B - op(A) * X,
277* where op(A) = A or A**T, depending on TRANS.
278*
279 CALL scopy( n, x( 1, j ), 1, work( n+1 ), 1 )
280 CALL stpmv( uplo, trans, diag, n, ap, work( n+1 ), 1 )
281 CALL saxpy( n, -one, b( 1, j ), 1, work( n+1 ), 1 )
282*
283* Compute componentwise relative backward error from formula
284*
285* max(i) ( abs(R(i)) / ( abs(op(A))*abs(X) + abs(B) )(i) )
286*
287* where abs(Z) is the componentwise absolute value of the matrix
288* or vector Z. If the i-th component of the denominator is less
289* than SAFE2, then SAFE1 is added to the i-th components of the
290* numerator and denominator before dividing.
291*
292 DO 20 i = 1, n
293 work( i ) = abs( b( i, j ) )
294 20 CONTINUE
295*
296 IF( notran ) THEN
297*
298* Compute abs(A)*abs(X) + abs(B).
299*
300 IF( upper ) THEN
301 kc = 1
302 IF( nounit ) THEN
303 DO 40 k = 1, n
304 xk = abs( x( k, j ) )
305 DO 30 i = 1, k
306 work( i ) = work( i ) + abs( ap( kc+i-1 ) )*xk
307 30 CONTINUE
308 kc = kc + k
309 40 CONTINUE
310 ELSE
311 DO 60 k = 1, n
312 xk = abs( x( k, j ) )
313 DO 50 i = 1, k - 1
314 work( i ) = work( i ) + abs( ap( kc+i-1 ) )*xk
315 50 CONTINUE
316 work( k ) = work( k ) + xk
317 kc = kc + k
318 60 CONTINUE
319 END IF
320 ELSE
321 kc = 1
322 IF( nounit ) THEN
323 DO 80 k = 1, n
324 xk = abs( x( k, j ) )
325 DO 70 i = k, n
326 work( i ) = work( i ) + abs( ap( kc+i-k ) )*xk
327 70 CONTINUE
328 kc = kc + n - k + 1
329 80 CONTINUE
330 ELSE
331 DO 100 k = 1, n
332 xk = abs( x( k, j ) )
333 DO 90 i = k + 1, n
334 work( i ) = work( i ) + abs( ap( kc+i-k ) )*xk
335 90 CONTINUE
336 work( k ) = work( k ) + xk
337 kc = kc + n - k + 1
338 100 CONTINUE
339 END IF
340 END IF
341 ELSE
342*
343* Compute abs(A**T)*abs(X) + abs(B).
344*
345 IF( upper ) THEN
346 kc = 1
347 IF( nounit ) THEN
348 DO 120 k = 1, n
349 s = zero
350 DO 110 i = 1, k
351 s = s + abs( ap( kc+i-1 ) )*abs( x( i, j ) )
352 110 CONTINUE
353 work( k ) = work( k ) + s
354 kc = kc + k
355 120 CONTINUE
356 ELSE
357 DO 140 k = 1, n
358 s = abs( x( k, j ) )
359 DO 130 i = 1, k - 1
360 s = s + abs( ap( kc+i-1 ) )*abs( x( i, j ) )
361 130 CONTINUE
362 work( k ) = work( k ) + s
363 kc = kc + k
364 140 CONTINUE
365 END IF
366 ELSE
367 kc = 1
368 IF( nounit ) THEN
369 DO 160 k = 1, n
370 s = zero
371 DO 150 i = k, n
372 s = s + abs( ap( kc+i-k ) )*abs( x( i, j ) )
373 150 CONTINUE
374 work( k ) = work( k ) + s
375 kc = kc + n - k + 1
376 160 CONTINUE
377 ELSE
378 DO 180 k = 1, n
379 s = abs( x( k, j ) )
380 DO 170 i = k + 1, n
381 s = s + abs( ap( kc+i-k ) )*abs( x( i, j ) )
382 170 CONTINUE
383 work( k ) = work( k ) + s
384 kc = kc + n - k + 1
385 180 CONTINUE
386 END IF
387 END IF
388 END IF
389 s = zero
390 DO 190 i = 1, n
391 IF( work( i ).GT.safe2 ) THEN
392 s = max( s, abs( work( n+i ) ) / work( i ) )
393 ELSE
394 s = max( s, ( abs( work( n+i ) )+safe1 ) /
395 $ ( work( i )+safe1 ) )
396 END IF
397 190 CONTINUE
398 berr( j ) = s
399*
400* Bound error from formula
401*
402* norm(X - XTRUE) / norm(X) .le. FERR =
403* norm( abs(inv(op(A)))*
404* ( abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) ))) / norm(X)
405*
406* where
407* norm(Z) is the magnitude of the largest component of Z
408* inv(op(A)) is the inverse of op(A)
409* abs(Z) is the componentwise absolute value of the matrix or
410* vector Z
411* NZ is the maximum number of nonzeros in any row of A, plus 1
412* EPS is machine epsilon
413*
414* The i-th component of abs(R)+NZ*EPS*(abs(op(A))*abs(X)+abs(B))
415* is incremented by SAFE1 if the i-th component of
416* abs(op(A))*abs(X) + abs(B) is less than SAFE2.
417*
418* Use SLACN2 to estimate the infinity-norm of the matrix
419* inv(op(A)) * diag(W),
420* where W = abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) )))
421*
422 DO 200 i = 1, n
423 IF( work( i ).GT.safe2 ) THEN
424 work( i ) = abs( work( n+i ) ) + real( nz )*eps*work( i )
425 ELSE
426 work( i ) = abs( work( n+i ) ) + real( nz )*eps*work( i )
427 $ + safe1
428 END IF
429 200 CONTINUE
430*
431 kase = 0
432 210 CONTINUE
433 CALL slacn2( n, work( 2*n+1 ), work( n+1 ), iwork,
434 $ ferr( j ),
435 $ kase, isave )
436 IF( kase.NE.0 ) THEN
437 IF( kase.EQ.1 ) THEN
438*
439* Multiply by diag(W)*inv(op(A)**T).
440*
441 CALL stpsv( uplo, transt, diag, n, ap, work( n+1 ),
442 $ 1 )
443 DO 220 i = 1, n
444 work( n+i ) = work( i )*work( n+i )
445 220 CONTINUE
446 ELSE
447*
448* Multiply by inv(op(A))*diag(W).
449*
450 DO 230 i = 1, n
451 work( n+i ) = work( i )*work( n+i )
452 230 CONTINUE
453 CALL stpsv( uplo, trans, diag, n, ap, work( n+1 ), 1 )
454 END IF
455 GO TO 210
456 END IF
457*
458* Normalize error.
459*
460 lstres = zero
461 DO 240 i = 1, n
462 lstres = max( lstres, abs( x( i, j ) ) )
463 240 CONTINUE
464 IF( lstres.NE.zero )
465 $ ferr( j ) = ferr( j ) / lstres
466*
467 250 CONTINUE
468*
469 RETURN
470*
471* End of STPRFS
472*
473 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 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
subroutine stpmv(uplo, trans, diag, n, ap, x, incx)
STPMV
Definition stpmv.f:142
subroutine stprfs(uplo, trans, diag, n, nrhs, ap, b, ldb, x, ldx, ferr, berr, work, iwork, info)
STPRFS
Definition stprfs.f:174
subroutine stpsv(uplo, trans, diag, n, ap, x, incx)
STPSV
Definition stpsv.f:144