LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dgtt01 ( integer  N,
double precision, dimension( * )  DL,
double precision, dimension( * )  D,
double precision, dimension( * )  DU,
double precision, dimension( * )  DLF,
double precision, dimension( * )  DF,
double precision, dimension( * )  DUF,
double precision, dimension( * )  DU2,
integer, dimension( * )  IPIV,
double precision, dimension( ldwork, * )  WORK,
integer  LDWORK,
double precision, dimension( * )  RWORK,
double precision  RESID 
)

DGTT01

Purpose:
 DGTT01 reconstructs a tridiagonal matrix A from its LU factorization
 and computes the residual
    norm(L*U - A) / ( norm(A) * EPS ),
 where EPS is the machine epsilon.
Parameters
[in]N
          N is INTEGTER
          The order of the matrix A.  N >= 0.
[in]DL
          DL is DOUBLE PRECISION array, dimension (N-1)
          The (n-1) sub-diagonal elements of A.
[in]D
          D is DOUBLE PRECISION array, dimension (N)
          The diagonal elements of A.
[in]DU
          DU is DOUBLE PRECISION array, dimension (N-1)
          The (n-1) super-diagonal elements of A.
[in]DLF
          DLF is DOUBLE PRECISION array, dimension (N-1)
          The (n-1) multipliers that define the matrix L from the
          LU factorization of A.
[in]DF
          DF is DOUBLE PRECISION array, dimension (N)
          The n diagonal elements of the upper triangular matrix U from
          the LU factorization of A.
[in]DUF
          DUF is DOUBLE PRECISION array, dimension (N-1)
          The (n-1) elements of the first super-diagonal of U.
[in]DU2
          DU2 is DOUBLE PRECISION array, dimension (N-2)
          The (n-2) elements of the second super-diagonal of U.
[in]IPIV
          IPIV is INTEGER array, dimension (N)
          The pivot indices; for 1 <= i <= n, row i of the matrix was
          interchanged with row IPIV(i).  IPIV(i) will always be either
          i or i+1; IPIV(i) = i indicates a row interchange was not
          required.
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (LDWORK,N)
[in]LDWORK
          LDWORK is INTEGER
          The leading dimension of the array WORK.  LDWORK >= max(1,N).
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (N)
[out]RESID
          RESID is DOUBLE PRECISION
          The scaled residual:  norm(L*U - A) / (norm(A) * EPS)
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 136 of file dgtt01.f.

136 *
137 * -- LAPACK test routine (version 3.4.0) --
138 * -- LAPACK is a software package provided by Univ. of Tennessee, --
139 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
140 * November 2011
141 *
142 * .. Scalar Arguments ..
143  INTEGER ldwork, n
144  DOUBLE PRECISION resid
145 * ..
146 * .. Array Arguments ..
147  INTEGER ipiv( * )
148  DOUBLE PRECISION d( * ), df( * ), dl( * ), dlf( * ), du( * ),
149  $ du2( * ), duf( * ), rwork( * ),
150  $ work( ldwork, * )
151 * ..
152 *
153 * =====================================================================
154 *
155 * .. Parameters ..
156  DOUBLE PRECISION one, zero
157  parameter ( one = 1.0d+0, zero = 0.0d+0 )
158 * ..
159 * .. Local Scalars ..
160  INTEGER i, ip, j, lastj
161  DOUBLE PRECISION anorm, eps, li
162 * ..
163 * .. External Functions ..
164  DOUBLE PRECISION dlamch, dlangt, dlanhs
165  EXTERNAL dlamch, dlangt, dlanhs
166 * ..
167 * .. Intrinsic Functions ..
168  INTRINSIC min
169 * ..
170 * .. External Subroutines ..
171  EXTERNAL daxpy, dswap
172 * ..
173 * .. Executable Statements ..
174 *
175 * Quick return if possible
176 *
177  IF( n.LE.0 ) THEN
178  resid = zero
179  RETURN
180  END IF
181 *
182  eps = dlamch( 'Epsilon' )
183 *
184 * Copy the matrix U to WORK.
185 *
186  DO 20 j = 1, n
187  DO 10 i = 1, n
188  work( i, j ) = zero
189  10 CONTINUE
190  20 CONTINUE
191  DO 30 i = 1, n
192  IF( i.EQ.1 ) THEN
193  work( i, i ) = df( i )
194  IF( n.GE.2 )
195  $ work( i, i+1 ) = duf( i )
196  IF( n.GE.3 )
197  $ work( i, i+2 ) = du2( i )
198  ELSE IF( i.EQ.n ) THEN
199  work( i, i ) = df( i )
200  ELSE
201  work( i, i ) = df( i )
202  work( i, i+1 ) = duf( i )
203  IF( i.LT.n-1 )
204  $ work( i, i+2 ) = du2( i )
205  END IF
206  30 CONTINUE
207 *
208 * Multiply on the left by L.
209 *
210  lastj = n
211  DO 40 i = n - 1, 1, -1
212  li = dlf( i )
213  CALL daxpy( lastj-i+1, li, work( i, i ), ldwork,
214  $ work( i+1, i ), ldwork )
215  ip = ipiv( i )
216  IF( ip.EQ.i ) THEN
217  lastj = min( i+2, n )
218  ELSE
219  CALL dswap( lastj-i+1, work( i, i ), ldwork, work( i+1, i ),
220  $ ldwork )
221  END IF
222  40 CONTINUE
223 *
224 * Subtract the matrix A.
225 *
226  work( 1, 1 ) = work( 1, 1 ) - d( 1 )
227  IF( n.GT.1 ) THEN
228  work( 1, 2 ) = work( 1, 2 ) - du( 1 )
229  work( n, n-1 ) = work( n, n-1 ) - dl( n-1 )
230  work( n, n ) = work( n, n ) - d( n )
231  DO 50 i = 2, n - 1
232  work( i, i-1 ) = work( i, i-1 ) - dl( i-1 )
233  work( i, i ) = work( i, i ) - d( i )
234  work( i, i+1 ) = work( i, i+1 ) - du( i )
235  50 CONTINUE
236  END IF
237 *
238 * Compute the 1-norm of the tridiagonal matrix A.
239 *
240  anorm = dlangt( '1', n, dl, d, du )
241 *
242 * Compute the 1-norm of WORK, which is only guaranteed to be
243 * upper Hessenberg.
244 *
245  resid = dlanhs( '1', n, work, ldwork, rwork )
246 *
247 * Compute norm(L*U - A) / (norm(A) * EPS)
248 *
249  IF( anorm.LE.zero ) THEN
250  IF( resid.NE.zero )
251  $ resid = one / eps
252  ELSE
253  resid = ( resid / anorm ) / eps
254  END IF
255 *
256  RETURN
257 *
258 * End of DGTT01
259 *
double precision function dlanhs(NORM, N, A, LDA, WORK)
DLANHS returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: dlanhs.f:110
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine daxpy(N, DA, DX, INCX, DY, INCY)
DAXPY
Definition: daxpy.f:54
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
double precision function dlangt(NORM, N, DL, D, DU)
DLANGT returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition: dlangt.f:108

Here is the call graph for this function:

Here is the caller graph for this function: