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