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