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

◆ zhst01()

subroutine zhst01 ( integer n,
integer ilo,
integer ihi,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldh, * ) h,
integer ldh,
complex*16, dimension( ldq, * ) q,
integer ldq,
complex*16, dimension( lwork ) work,
integer lwork,
double precision, dimension( * ) rwork,
double precision, dimension( 2 ) result )

ZHST01

Purpose:
!>
!> ZHST01 tests the reduction of a general matrix A to upper Hessenberg
!> form:  A = Q*H*Q'.  Two test ratios are computed;
!>
!> RESULT(1) = norm( A - Q*H*Q' ) / ( norm(A) * N * EPS )
!> RESULT(2) = norm( I - Q'*Q ) / ( N * EPS )
!>
!> The matrix Q is assumed to be given explicitly as it would be
!> following ZGEHRD + ZUNGHR.
!>
!> In this version, ILO and IHI are not used, but they could be used
!> to save some work if this is desired.
!> 
Parameters
[in]N
!>          N is INTEGER
!>          The order of the matrix A.  N >= 0.
!> 
[in]ILO
!>          ILO is INTEGER
!> 
[in]IHI
!>          IHI is INTEGER
!>
!>          A is assumed to be upper triangular in rows and columns
!>          1:ILO-1 and IHI+1:N, so Q differs from the identity only in
!>          rows and columns ILO+1:IHI.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension (LDA,N)
!>          The original n by n matrix A.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in]H
!>          H is COMPLEX*16 array, dimension (LDH,N)
!>          The upper Hessenberg matrix H from the reduction A = Q*H*Q'
!>          as computed by ZGEHRD.  H is assumed to be zero below the
!>          first subdiagonal.
!> 
[in]LDH
!>          LDH is INTEGER
!>          The leading dimension of the array H.  LDH >= max(1,N).
!> 
[in]Q
!>          Q is COMPLEX*16 array, dimension (LDQ,N)
!>          The orthogonal matrix Q from the reduction A = Q*H*Q' as
!>          computed by ZGEHRD + ZUNGHR.
!> 
[in]LDQ
!>          LDQ is INTEGER
!>          The leading dimension of the array Q.  LDQ >= max(1,N).
!> 
[out]WORK
!>          WORK is COMPLEX*16 array, dimension (LWORK)
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The length of the array WORK.  LWORK >= 2*N*N.
!> 
[out]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension (N)
!> 
[out]RESULT
!>          RESULT is DOUBLE PRECISION array, dimension (2)
!>          RESULT(1) = norm( A - Q*H*Q' ) / ( norm(A) * N * EPS )
!>          RESULT(2) = norm( I - Q'*Q ) / ( N * EPS )
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 138 of file zhst01.f.

140*
141* -- LAPACK test routine --
142* -- LAPACK is a software package provided by Univ. of Tennessee, --
143* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
144*
145* .. Scalar Arguments ..
146 INTEGER IHI, ILO, LDA, LDH, LDQ, LWORK, N
147* ..
148* .. Array Arguments ..
149 DOUBLE PRECISION RESULT( 2 ), RWORK( * )
150 COMPLEX*16 A( LDA, * ), H( LDH, * ), Q( LDQ, * ),
151 $ WORK( LWORK )
152* ..
153*
154* =====================================================================
155*
156* .. Parameters ..
157 DOUBLE PRECISION ONE, ZERO
158 parameter( one = 1.0d+0, zero = 0.0d+0 )
159* ..
160* .. Local Scalars ..
161 INTEGER LDWORK
162 DOUBLE PRECISION ANORM, EPS, OVFL, SMLNUM, UNFL, WNORM
163* ..
164* .. External Functions ..
165 DOUBLE PRECISION DLAMCH, ZLANGE
166 EXTERNAL dlamch, zlange
167* ..
168* .. External Subroutines ..
169 EXTERNAL zgemm, zlacpy, zunt01
170* ..
171* .. Intrinsic Functions ..
172 INTRINSIC dcmplx, max, min
173* ..
174* .. Executable Statements ..
175*
176* Quick return if possible
177*
178 IF( n.LE.0 ) THEN
179 result( 1 ) = zero
180 result( 2 ) = zero
181 RETURN
182 END IF
183*
184 unfl = dlamch( 'Safe minimum' )
185 eps = dlamch( 'Precision' )
186 ovfl = one / unfl
187 smlnum = unfl*n / eps
188*
189* Test 1: Compute norm( A - Q*H*Q' ) / ( norm(A) * N * EPS )
190*
191* Copy A to WORK
192*
193 ldwork = max( 1, n )
194 CALL zlacpy( ' ', n, n, a, lda, work, ldwork )
195*
196* Compute Q*H
197*
198 CALL zgemm( 'No transpose', 'No transpose', n, n, n,
199 $ dcmplx( one ), q, ldq, h, ldh, dcmplx( zero ),
200 $ work( ldwork*n+1 ), ldwork )
201*
202* Compute A - Q*H*Q'
203*
204 CALL zgemm( 'No transpose', 'Conjugate transpose', n, n, n,
205 $ dcmplx( -one ), work( ldwork*n+1 ), ldwork, q, ldq,
206 $ dcmplx( one ), work, ldwork )
207*
208 anorm = max( zlange( '1', n, n, a, lda, rwork ), unfl )
209 wnorm = zlange( '1', n, n, work, ldwork, rwork )
210*
211* Note that RESULT(1) cannot overflow and is bounded by 1/(N*EPS)
212*
213 result( 1 ) = min( wnorm, anorm ) / max( smlnum, anorm*eps ) / n
214*
215* Test 2: Compute norm( I - Q'*Q ) / ( N * EPS )
216*
217 CALL zunt01( 'Columns', n, n, q, ldq, work, lwork, rwork,
218 $ result( 2 ) )
219*
220 RETURN
221*
222* End of ZHST01
223*
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:188
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlange.f:113
subroutine zunt01(rowcol, m, n, u, ldu, work, lwork, rwork, resid)
ZUNT01
Definition zunt01.f:126
Here is the call graph for this function:
Here is the caller graph for this function: