LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zgtt01()

subroutine zgtt01 ( integer  n,
complex*16, dimension( * )  dl,
complex*16, dimension( * )  d,
complex*16, dimension( * )  du,
complex*16, dimension( * )  dlf,
complex*16, dimension( * )  df,
complex*16, dimension( * )  duf,
complex*16, dimension( * )  du2,
integer, dimension( * )  ipiv,
complex*16, dimension( ldwork, * )  work,
integer  ldwork,
double precision, dimension( * )  rwork,
double precision  resid 
)

ZGTT01

Purpose:
 ZGTT01 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 INTEGER
          The order of the matrix A.  N >= 0.
[in]DL
          DL is COMPLEX*16 array, dimension (N-1)
          The (n-1) sub-diagonal elements of A.
[in]D
          D is COMPLEX*16 array, dimension (N)
          The diagonal elements of A.
[in]DU
          DU is COMPLEX*16 array, dimension (N-1)
          The (n-1) super-diagonal elements of A.
[in]DLF
          DLF is COMPLEX*16 array, dimension (N-1)
          The (n-1) multipliers that define the matrix L from the
          LU factorization of A.
[in]DF
          DF is COMPLEX*16 array, dimension (N)
          The n diagonal elements of the upper triangular matrix U from
          the LU factorization of A.
[in]DUF
          DUF is COMPLEX*16 array, dimension (N-1)
          The (n-1) elements of the first super-diagonal of U.
[in]DU2
          DU2 is COMPLEX*16 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 COMPLEX*16 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.

Definition at line 132 of file zgtt01.f.

134*
135* -- LAPACK test routine --
136* -- LAPACK is a software package provided by Univ. of Tennessee, --
137* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
138*
139* .. Scalar Arguments ..
140 INTEGER LDWORK, N
141 DOUBLE PRECISION RESID
142* ..
143* .. Array Arguments ..
144 INTEGER IPIV( * )
145 DOUBLE PRECISION RWORK( * )
146 COMPLEX*16 D( * ), DF( * ), DL( * ), DLF( * ), DU( * ),
147 $ DU2( * ), DUF( * ), WORK( LDWORK, * )
148* ..
149*
150* =====================================================================
151*
152* .. Parameters ..
153 DOUBLE PRECISION ONE, ZERO
154 parameter( one = 1.0d+0, zero = 0.0d+0 )
155* ..
156* .. Local Scalars ..
157 INTEGER I, IP, J, LASTJ
158 DOUBLE PRECISION ANORM, EPS
159 COMPLEX*16 LI
160* ..
161* .. External Functions ..
162 DOUBLE PRECISION DLAMCH, ZLANGT, ZLANHS
163 EXTERNAL dlamch, zlangt, zlanhs
164* ..
165* .. Intrinsic Functions ..
166 INTRINSIC min
167* ..
168* .. External Subroutines ..
169 EXTERNAL zaxpy, zswap
170* ..
171* .. Executable Statements ..
172*
173* Quick return if possible
174*
175 IF( n.LE.0 ) THEN
176 resid = zero
177 RETURN
178 END IF
179*
180 eps = dlamch( 'Epsilon' )
181*
182* Copy the matrix U to WORK.
183*
184 DO 20 j = 1, n
185 DO 10 i = 1, n
186 work( i, j ) = zero
187 10 CONTINUE
188 20 CONTINUE
189 DO 30 i = 1, n
190 IF( i.EQ.1 ) THEN
191 work( i, i ) = df( i )
192 IF( n.GE.2 )
193 $ work( i, i+1 ) = duf( i )
194 IF( n.GE.3 )
195 $ work( i, i+2 ) = du2( i )
196 ELSE IF( i.EQ.n ) THEN
197 work( i, i ) = df( i )
198 ELSE
199 work( i, i ) = df( i )
200 work( i, i+1 ) = duf( i )
201 IF( i.LT.n-1 )
202 $ work( i, i+2 ) = du2( i )
203 END IF
204 30 CONTINUE
205*
206* Multiply on the left by L.
207*
208 lastj = n
209 DO 40 i = n - 1, 1, -1
210 li = dlf( i )
211 CALL zaxpy( lastj-i+1, li, work( i, i ), ldwork,
212 $ work( i+1, i ), ldwork )
213 ip = ipiv( i )
214 IF( ip.EQ.i ) THEN
215 lastj = min( i+2, n )
216 ELSE
217 CALL zswap( lastj-i+1, work( i, i ), ldwork, work( i+1, i ),
218 $ ldwork )
219 END IF
220 40 CONTINUE
221*
222* Subtract the matrix A.
223*
224 work( 1, 1 ) = work( 1, 1 ) - d( 1 )
225 IF( n.GT.1 ) THEN
226 work( 1, 2 ) = work( 1, 2 ) - du( 1 )
227 work( n, n-1 ) = work( n, n-1 ) - dl( n-1 )
228 work( n, n ) = work( n, n ) - d( n )
229 DO 50 i = 2, n - 1
230 work( i, i-1 ) = work( i, i-1 ) - dl( i-1 )
231 work( i, i ) = work( i, i ) - d( i )
232 work( i, i+1 ) = work( i, i+1 ) - du( i )
233 50 CONTINUE
234 END IF
235*
236* Compute the 1-norm of the tridiagonal matrix A.
237*
238 anorm = zlangt( '1', n, dl, d, du )
239*
240* Compute the 1-norm of WORK, which is only guaranteed to be
241* upper Hessenberg.
242*
243 resid = zlanhs( '1', n, work, ldwork, rwork )
244*
245* Compute norm(L*U - A) / (norm(A) * EPS)
246*
247 IF( anorm.LE.zero ) THEN
248 IF( resid.NE.zero )
249 $ resid = one / eps
250 ELSE
251 resid = ( resid / anorm ) / eps
252 END IF
253*
254 RETURN
255*
256* End of ZGTT01
257*
subroutine zaxpy(n, za, zx, incx, zy, incy)
ZAXPY
Definition zaxpy.f:88
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlangt(norm, n, dl, d, du)
ZLANGT returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlangt.f:106
double precision function zlanhs(norm, n, a, lda, work)
ZLANHS returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlanhs.f:109
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81
Here is the call graph for this function:
Here is the caller graph for this function: