LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
dget01.f
Go to the documentation of this file.
1 *> \brief \b DGET01
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * SUBROUTINE DGET01( M, N, A, LDA, AFAC, LDAFAC, IPIV, RWORK,
12 * RESID )
13 *
14 * .. Scalar Arguments ..
15 * INTEGER LDA, LDAFAC, M, N
16 * DOUBLE PRECISION RESID
17 * ..
18 * .. Array Arguments ..
19 * INTEGER IPIV( * )
20 * DOUBLE PRECISION A( LDA, * ), AFAC( LDAFAC, * ), RWORK( * )
21 * ..
22 *
23 *
24 *> \par Purpose:
25 * =============
26 *>
27 *> \verbatim
28 *>
29 *> DGET01 reconstructs a matrix A from its L*U factorization and
30 *> computes the residual
31 *> norm(L*U - A) / ( N * norm(A) * EPS ),
32 *> where EPS is the machine epsilon.
33 *> \endverbatim
34 *
35 * Arguments:
36 * ==========
37 *
38 *> \param[in] M
39 *> \verbatim
40 *> M is INTEGER
41 *> The number of rows of the matrix A. M >= 0.
42 *> \endverbatim
43 *>
44 *> \param[in] N
45 *> \verbatim
46 *> N is INTEGER
47 *> The number of columns of the matrix A. N >= 0.
48 *> \endverbatim
49 *>
50 *> \param[in] A
51 *> \verbatim
52 *> A is DOUBLE PRECISION array, dimension (LDA,N)
53 *> The original M x N matrix A.
54 *> \endverbatim
55 *>
56 *> \param[in] LDA
57 *> \verbatim
58 *> LDA is INTEGER
59 *> The leading dimension of the array A. LDA >= max(1,M).
60 *> \endverbatim
61 *>
62 *> \param[in,out] AFAC
63 *> \verbatim
64 *> AFAC is DOUBLE PRECISION array, dimension (LDAFAC,N)
65 *> The factored form of the matrix A. AFAC contains the factors
66 *> L and U from the L*U factorization as computed by DGETRF.
67 *> Overwritten with the reconstructed matrix, and then with the
68 *> difference L*U - A.
69 *> \endverbatim
70 *>
71 *> \param[in] LDAFAC
72 *> \verbatim
73 *> LDAFAC is INTEGER
74 *> The leading dimension of the array AFAC. LDAFAC >= max(1,M).
75 *> \endverbatim
76 *>
77 *> \param[in] IPIV
78 *> \verbatim
79 *> IPIV is INTEGER array, dimension (N)
80 *> The pivot indices from DGETRF.
81 *> \endverbatim
82 *>
83 *> \param[out] RWORK
84 *> \verbatim
85 *> RWORK is DOUBLE PRECISION array, dimension (M)
86 *> \endverbatim
87 *>
88 *> \param[out] RESID
89 *> \verbatim
90 *> RESID is DOUBLE PRECISION
91 *> norm(L*U - A) / ( N * norm(A) * EPS )
92 *> \endverbatim
93 *
94 * Authors:
95 * ========
96 *
97 *> \author Univ. of Tennessee
98 *> \author Univ. of California Berkeley
99 *> \author Univ. of Colorado Denver
100 *> \author NAG Ltd.
101 *
102 *> \date November 2011
103 *
104 *> \ingroup double_lin
105 *
106 * =====================================================================
107  SUBROUTINE dget01( M, N, A, LDA, AFAC, LDAFAC, IPIV, RWORK,
108  $ resid )
109 *
110 * -- LAPACK test routine (version 3.4.0) --
111 * -- LAPACK is a software package provided by Univ. of Tennessee, --
112 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
113 * November 2011
114 *
115 * .. Scalar Arguments ..
116  INTEGER LDA, LDAFAC, M, N
117  DOUBLE PRECISION RESID
118 * ..
119 * .. Array Arguments ..
120  INTEGER IPIV( * )
121  DOUBLE PRECISION A( lda, * ), AFAC( ldafac, * ), RWORK( * )
122 * ..
123 *
124 * =====================================================================
125 *
126 *
127 * .. Parameters ..
128  DOUBLE PRECISION ZERO, ONE
129  parameter ( zero = 0.0d+0, one = 1.0d+0 )
130 * ..
131 * .. Local Scalars ..
132  INTEGER I, J, K
133  DOUBLE PRECISION ANORM, EPS, T
134 * ..
135 * .. External Functions ..
136  DOUBLE PRECISION DDOT, DLAMCH, DLANGE
137  EXTERNAL ddot, dlamch, dlange
138 * ..
139 * .. External Subroutines ..
140  EXTERNAL dgemv, dlaswp, dscal, dtrmv
141 * ..
142 * .. Intrinsic Functions ..
143  INTRINSIC dble, min
144 * ..
145 * .. Executable Statements ..
146 *
147 * Quick exit if M = 0 or N = 0.
148 *
149  IF( m.LE.0 .OR. n.LE.0 ) THEN
150  resid = zero
151  RETURN
152  END IF
153 *
154 * Determine EPS and the norm of A.
155 *
156  eps = dlamch( 'Epsilon' )
157  anorm = dlange( '1', m, n, a, lda, rwork )
158 *
159 * Compute the product L*U and overwrite AFAC with the result.
160 * A column at a time of the product is obtained, starting with
161 * column N.
162 *
163  DO 10 k = n, 1, -1
164  IF( k.GT.m ) THEN
165  CALL dtrmv( 'Lower', 'No transpose', 'Unit', m, afac,
166  $ ldafac, afac( 1, k ), 1 )
167  ELSE
168 *
169 * Compute elements (K+1:M,K)
170 *
171  t = afac( k, k )
172  IF( k+1.LE.m ) THEN
173  CALL dscal( m-k, t, afac( k+1, k ), 1 )
174  CALL dgemv( 'No transpose', m-k, k-1, one,
175  $ afac( k+1, 1 ), ldafac, afac( 1, k ), 1, one,
176  $ afac( k+1, k ), 1 )
177  END IF
178 *
179 * Compute the (K,K) element
180 *
181  afac( k, k ) = t + ddot( k-1, afac( k, 1 ), ldafac,
182  $ afac( 1, k ), 1 )
183 *
184 * Compute elements (1:K-1,K)
185 *
186  CALL dtrmv( 'Lower', 'No transpose', 'Unit', k-1, afac,
187  $ ldafac, afac( 1, k ), 1 )
188  END IF
189  10 CONTINUE
190  CALL dlaswp( n, afac, ldafac, 1, min( m, n ), ipiv, -1 )
191 *
192 * Compute the difference L*U - A and store in AFAC.
193 *
194  DO 30 j = 1, n
195  DO 20 i = 1, m
196  afac( i, j ) = afac( i, j ) - a( i, j )
197  20 CONTINUE
198  30 CONTINUE
199 *
200 * Compute norm( L*U - A ) / ( N * norm(A) * EPS )
201 *
202  resid = dlange( '1', m, n, afac, ldafac, rwork )
203 *
204  IF( anorm.LE.zero ) THEN
205  IF( resid.NE.zero )
206  $ resid = one / eps
207  ELSE
208  resid = ( ( resid / dble( n ) ) / anorm ) / eps
209  END IF
210 *
211  RETURN
212 *
213 * End of DGET01
214 *
215  END
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
subroutine dget01(M, N, A, LDA, AFAC, LDAFAC, IPIV, RWORK, RESID)
DGET01
Definition: dget01.f:109
subroutine dlaswp(N, A, LDA, K1, K2, IPIV, INCX)
DLASWP performs a series of row interchanges on a general rectangular matrix.
Definition: dlaswp.f:116
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
subroutine dtrmv(UPLO, TRANS, DIAG, N, A, LDA, X, INCX)
DTRMV
Definition: dtrmv.f:149