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