LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.
Date
November 2011

Definition at line 142 of file zhst01.f.

142 *
143 * -- LAPACK test routine (version 3.4.0) --
144 * -- LAPACK is a software package provided by Univ. of Tennessee, --
145 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
146 * November 2011
147 *
148 * .. Scalar Arguments ..
149  INTEGER ihi, ilo, lda, ldh, ldq, lwork, n
150 * ..
151 * .. Array Arguments ..
152  DOUBLE PRECISION result( 2 ), rwork( * )
153  COMPLEX*16 a( lda, * ), h( ldh, * ), q( ldq, * ),
154  $ work( lwork )
155 * ..
156 *
157 * =====================================================================
158 *
159 * .. Parameters ..
160  DOUBLE PRECISION one, zero
161  parameter ( one = 1.0d+0, zero = 0.0d+0 )
162 * ..
163 * .. Local Scalars ..
164  INTEGER ldwork
165  DOUBLE PRECISION anorm, eps, ovfl, smlnum, unfl, wnorm
166 * ..
167 * .. External Functions ..
168  DOUBLE PRECISION dlamch, zlange
169  EXTERNAL dlamch, zlange
170 * ..
171 * .. External Subroutines ..
172  EXTERNAL dlabad, zgemm, zlacpy, zunt01
173 * ..
174 * .. Intrinsic Functions ..
175  INTRINSIC dcmplx, max, min
176 * ..
177 * .. Executable Statements ..
178 *
179 * Quick return if possible
180 *
181  IF( n.LE.0 ) THEN
182  result( 1 ) = zero
183  result( 2 ) = zero
184  RETURN
185  END IF
186 *
187  unfl = dlamch( 'Safe minimum' )
188  eps = dlamch( 'Precision' )
189  ovfl = one / unfl
190  CALL dlabad( unfl, ovfl )
191  smlnum = unfl*n / eps
192 *
193 * Test 1: Compute norm( A - Q*H*Q' ) / ( norm(A) * N * EPS )
194 *
195 * Copy A to WORK
196 *
197  ldwork = max( 1, n )
198  CALL zlacpy( ' ', n, n, a, lda, work, ldwork )
199 *
200 * Compute Q*H
201 *
202  CALL zgemm( 'No transpose', 'No transpose', n, n, n,
203  $ dcmplx( one ), q, ldq, h, ldh, dcmplx( zero ),
204  $ work( ldwork*n+1 ), ldwork )
205 *
206 * Compute A - Q*H*Q'
207 *
208  CALL zgemm( 'No transpose', 'Conjugate transpose', n, n, n,
209  $ dcmplx( -one ), work( ldwork*n+1 ), ldwork, q, ldq,
210  $ dcmplx( one ), work, ldwork )
211 *
212  anorm = max( zlange( '1', n, n, a, lda, rwork ), unfl )
213  wnorm = zlange( '1', n, n, work, ldwork, rwork )
214 *
215 * Note that RESULT(1) cannot overflow and is bounded by 1/(N*EPS)
216 *
217  result( 1 ) = min( wnorm, anorm ) / max( smlnum, anorm*eps ) / n
218 *
219 * Test 2: Compute norm( I - Q'*Q ) / ( N * EPS )
220 *
221  CALL zunt01( 'Columns', n, n, q, ldq, work, lwork, rwork,
222  $ result( 2 ) )
223 *
224  RETURN
225 *
226 * End of ZHST01
227 *
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:189
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:76
subroutine zunt01(ROWCOL, M, N, U, LDU, WORK, LWORK, RWORK, RESID)
ZUNT01
Definition: zunt01.f:128
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:117

Here is the call graph for this function:

Here is the caller graph for this function: