LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dget01 ( integer  M,
integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( ldafac, * )  AFAC,
integer  LDAFAC,
integer, dimension( * )  IPIV,
double precision, dimension( * )  RWORK,
double precision  RESID 
)

DGET01

Purpose:
 DGET01 reconstructs a matrix A from its L*U factorization and
 computes the residual
    norm(L*U - A) / ( N * norm(A) * EPS ),
 where EPS is the machine epsilon.
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix A.  N >= 0.
[in]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          The original M x N matrix A.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[in,out]AFAC
          AFAC is DOUBLE PRECISION array, dimension (LDAFAC,N)
          The factored form of the matrix A.  AFAC contains the factors
          L and U from the L*U factorization as computed by DGETRF.
          Overwritten with the reconstructed matrix, and then with the
          difference L*U - A.
[in]LDAFAC
          LDAFAC is INTEGER
          The leading dimension of the array AFAC.  LDAFAC >= max(1,M).
[in]IPIV
          IPIV is INTEGER array, dimension (N)
          The pivot indices from DGETRF.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (M)
[out]RESID
          RESID is DOUBLE PRECISION
          norm(L*U - A) / ( N * norm(A) * EPS )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 109 of file dget01.f.

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 *
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
double precision function ddot(N, DX, INCX, DY, INCY)
DDOT
Definition: ddot.f:53
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
double precision function dlange(NORM, M, N, A, LDA, WORK)
DLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: dlange.f:116
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

Here is the call graph for this function:

Here is the caller graph for this function: