LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zhpt21 ( integer  ITYPE,
character  UPLO,
integer  N,
integer  KBAND,
complex*16, dimension( * )  AP,
double precision, dimension( * )  D,
double precision, dimension( * )  E,
complex*16, dimension( ldu, * )  U,
integer  LDU,
complex*16, dimension( * )  VP,
complex*16, dimension( * )  TAU,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
double precision, dimension( 2 )  RESULT 
)

ZHPT21

Purpose:
 ZHPT21  generally checks a decomposition of the form

         A = U S UC>
 where * means conjugate transpose, A is hermitian, U is
 unitary, and S is diagonal (if KBAND=0) or (real) symmetric
 tridiagonal (if KBAND=1).  If ITYPE=1, then U is represented as
 a dense matrix, otherwise the U is expressed as a product of
 Householder transformations, whose vectors are stored in the
 array "V" and whose scaling constants are in "TAU"; we shall
 use the letter "V" to refer to the product of Householder
 transformations (which should be equal to U).

 Specifically, if ITYPE=1, then:

         RESULT(1) = | A - U S U* | / ( |A| n ulp ) *andC>         RESULT(2) = | I - UU* | / ( n ulp )

 If ITYPE=2, then:

         RESULT(1) = | A - V S V* | / ( |A| n ulp )

 If ITYPE=3, then:

         RESULT(1) = | I - UV* | / ( n ulp )

 Packed storage means that, for example, if UPLO='U', then the columns
 of the upper triangle of A are stored one after another, so that
 A(1,j+1) immediately follows A(j,j) in the array AP.  Similarly, if
 UPLO='L', then the columns of the lower triangle of A are stored one
 after another in AP, so that A(j+1,j+1) immediately follows A(n,j)
 in the array AP.  This means that A(i,j) is stored in:

    AP( i + j*(j-1)/2 )                 if UPLO='U'

    AP( i + (2*n-j)*(j-1)/2 )           if UPLO='L'

 The array VP bears the same relation to the matrix V that A does to
 AP.

 For ITYPE > 1, the transformation U is expressed as a product
 of Householder transformations:

    If UPLO='U', then  V = H(n-1)...H(1),  where

        H(j) = I  -  tau(j) v(j) v(j)C>
    and the first j-1 elements of v(j) are stored in V(1:j-1,j+1),
    (i.e., VP( j*(j+1)/2 + 1 : j*(j+1)/2 + j-1 ) ),
    the j-th element is 1, and the last n-j elements are 0.

    If UPLO='L', then  V = H(1)...H(n-1),  where

        H(j) = I  -  tau(j) v(j) v(j)C>
    and the first j elements of v(j) are 0, the (j+1)-st is 1, and the
    (j+2)-nd through n-th elements are stored in V(j+2:n,j) (i.e.,
    in VP( (2*n-j)*(j-1)/2 + j+2 : (2*n-j)*(j-1)/2 + n ) .)
Parameters
[in]ITYPE
          ITYPE is INTEGER
          Specifies the type of tests to be performed.
          1: U expressed as a dense unitary matrix:
             RESULT(1) = | A - U S U* | / ( |A| n ulp )   *andC>             RESULT(2) = | I - UU* | / ( n ulp )

          2: U expressed as a product V of Housholder transformations:
             RESULT(1) = | A - V S V* | / ( |A| n ulp )

          3: U expressed both as a dense unitary matrix and
             as a product of Housholder transformations:
             RESULT(1) = | I - UV* | / ( n ulp )
[in]UPLO
          UPLO is CHARACTER
          If UPLO='U', the upper triangle of A and V will be used and
          the (strictly) lower triangle will not be referenced.
          If UPLO='L', the lower triangle of A and V will be used and
          the (strictly) upper triangle will not be referenced.
[in]N
          N is INTEGER
          The size of the matrix.  If it is zero, ZHPT21 does nothing.
          It must be at least zero.
[in]KBAND
          KBAND is INTEGER
          The bandwidth of the matrix.  It may only be zero or one.
          If zero, then S is diagonal, and E is not referenced.  If
          one, then S is symmetric tri-diagonal.
[in]AP
          AP is COMPLEX*16 array, dimension (N*(N+1)/2)
          The original (unfactored) matrix.  It is assumed to be
          hermitian, and contains the columns of just the upper
          triangle (UPLO='U') or only the lower triangle (UPLO='L'),
          packed one after another.
[in]D
          D is DOUBLE PRECISION array, dimension (N)
          The diagonal of the (symmetric tri-) diagonal matrix.
[in]E
          E is DOUBLE PRECISION array, dimension (N)
          The off-diagonal of the (symmetric tri-) diagonal matrix.
          E(1) is the (1,2) and (2,1) element, E(2) is the (2,3) and
          (3,2) element, etc.
          Not referenced if KBAND=0.
[in]U
          U is COMPLEX*16 array, dimension (LDU, N)
          If ITYPE=1 or 3, this contains the unitary matrix in
          the decomposition, expressed as a dense matrix.  If ITYPE=2,
          then it is not referenced.
[in]LDU
          LDU is INTEGER
          The leading dimension of U.  LDU must be at least N and
          at least 1.
[in]VP
          VP is DOUBLE PRECISION array, dimension (N*(N+1)/2)
          If ITYPE=2 or 3, the columns of this array contain the
          Householder vectors used to describe the unitary matrix
          in the decomposition, as described in purpose.
          *NOTE* If ITYPE=2 or 3, V is modified and restored.  The
          subdiagonal (if UPLO='L') or the superdiagonal (if UPLO='U')
          is set to one, and later reset to its original value, during
          the course of the calculation.
          If ITYPE=1, then it is neither referenced nor modified.
[in]TAU
          TAU is COMPLEX*16 array, dimension (N)
          If ITYPE >= 2, then TAU(j) is the scalar factor of
          v(j) v(j)* in the Householder transformation H(j) of
          the product  U = H(1)...H(n-2)
          If ITYPE < 2, then TAU is not referenced.
[out]WORK
          WORK is COMPLEX*16 array, dimension (N**2)
          Workspace.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (N)
          Workspace.
[out]RESULT
          RESULT is DOUBLE PRECISION array, dimension (2)
          The values computed by the two tests described above.  The
          values are currently limited to 1/ulp, to avoid overflow.
          RESULT(1) is always modified.  RESULT(2) is modified only
          if ITYPE=1.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 225 of file zhpt21.f.

225 *
226 * -- LAPACK test routine (version 3.4.0) --
227 * -- LAPACK is a software package provided by Univ. of Tennessee, --
228 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
229 * November 2011
230 *
231 * .. Scalar Arguments ..
232  CHARACTER uplo
233  INTEGER itype, kband, ldu, n
234 * ..
235 * .. Array Arguments ..
236  DOUBLE PRECISION d( * ), e( * ), result( 2 ), rwork( * )
237  COMPLEX*16 ap( * ), tau( * ), u( ldu, * ), vp( * ),
238  $ work( * )
239 * ..
240 *
241 * =====================================================================
242 *
243 * .. Parameters ..
244  DOUBLE PRECISION zero, one, ten
245  parameter ( zero = 0.0d+0, one = 1.0d+0, ten = 10.0d+0 )
246  DOUBLE PRECISION half
247  parameter ( half = 1.0d+0 / 2.0d+0 )
248  COMPLEX*16 czero, cone
249  parameter ( czero = ( 0.0d+0, 0.0d+0 ),
250  $ cone = ( 1.0d+0, 0.0d+0 ) )
251 * ..
252 * .. Local Scalars ..
253  LOGICAL lower
254  CHARACTER cuplo
255  INTEGER iinfo, j, jp, jp1, jr, lap
256  DOUBLE PRECISION anorm, ulp, unfl, wnorm
257  COMPLEX*16 temp, vsave
258 * ..
259 * .. External Functions ..
260  LOGICAL lsame
261  DOUBLE PRECISION dlamch, zlange, zlanhp
262  COMPLEX*16 zdotc
263  EXTERNAL lsame, dlamch, zlange, zlanhp, zdotc
264 * ..
265 * .. External Subroutines ..
266  EXTERNAL zaxpy, zcopy, zgemm, zhpmv, zhpr, zhpr2,
267  $ zlacpy, zlaset, zupmtr
268 * ..
269 * .. Intrinsic Functions ..
270  INTRINSIC dble, dcmplx, max, min
271 * ..
272 * .. Executable Statements ..
273 *
274 * Constants
275 *
276  result( 1 ) = zero
277  IF( itype.EQ.1 )
278  $ result( 2 ) = zero
279  IF( n.LE.0 )
280  $ RETURN
281 *
282  lap = ( n*( n+1 ) ) / 2
283 *
284  IF( lsame( uplo, 'U' ) ) THEN
285  lower = .false.
286  cuplo = 'U'
287  ELSE
288  lower = .true.
289  cuplo = 'L'
290  END IF
291 *
292  unfl = dlamch( 'Safe minimum' )
293  ulp = dlamch( 'Epsilon' )*dlamch( 'Base' )
294 *
295 * Some Error Checks
296 *
297  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
298  result( 1 ) = ten / ulp
299  RETURN
300  END IF
301 *
302 * Do Test 1
303 *
304 * Norm of A:
305 *
306  IF( itype.EQ.3 ) THEN
307  anorm = one
308  ELSE
309  anorm = max( zlanhp( '1', cuplo, n, ap, rwork ), unfl )
310  END IF
311 *
312 * Compute error matrix:
313 *
314  IF( itype.EQ.1 ) THEN
315 *
316 * ITYPE=1: error = A - U S U*
317 *
318  CALL zlaset( 'Full', n, n, czero, czero, work, n )
319  CALL zcopy( lap, ap, 1, work, 1 )
320 *
321  DO 10 j = 1, n
322  CALL zhpr( cuplo, n, -d( j ), u( 1, j ), 1, work )
323  10 CONTINUE
324 *
325  IF( n.GT.1 .AND. kband.EQ.1 ) THEN
326  DO 20 j = 1, n - 1
327  CALL zhpr2( cuplo, n, -dcmplx( e( j ) ), u( 1, j ), 1,
328  $ u( 1, j-1 ), 1, work )
329  20 CONTINUE
330  END IF
331  wnorm = zlanhp( '1', cuplo, n, work, rwork )
332 *
333  ELSE IF( itype.EQ.2 ) THEN
334 *
335 * ITYPE=2: error = V S V* - A
336 *
337  CALL zlaset( 'Full', n, n, czero, czero, work, n )
338 *
339  IF( lower ) THEN
340  work( lap ) = d( n )
341  DO 40 j = n - 1, 1, -1
342  jp = ( ( 2*n-j )*( j-1 ) ) / 2
343  jp1 = jp + n - j
344  IF( kband.EQ.1 ) THEN
345  work( jp+j+1 ) = ( cone-tau( j ) )*e( j )
346  DO 30 jr = j + 2, n
347  work( jp+jr ) = -tau( j )*e( j )*vp( jp+jr )
348  30 CONTINUE
349  END IF
350 *
351  IF( tau( j ).NE.czero ) THEN
352  vsave = vp( jp+j+1 )
353  vp( jp+j+1 ) = cone
354  CALL zhpmv( 'L', n-j, cone, work( jp1+j+1 ),
355  $ vp( jp+j+1 ), 1, czero, work( lap+1 ), 1 )
356  temp = -half*tau( j )*zdotc( n-j, work( lap+1 ), 1,
357  $ vp( jp+j+1 ), 1 )
358  CALL zaxpy( n-j, temp, vp( jp+j+1 ), 1, work( lap+1 ),
359  $ 1 )
360  CALL zhpr2( 'L', n-j, -tau( j ), vp( jp+j+1 ), 1,
361  $ work( lap+1 ), 1, work( jp1+j+1 ) )
362 *
363  vp( jp+j+1 ) = vsave
364  END IF
365  work( jp+j ) = d( j )
366  40 CONTINUE
367  ELSE
368  work( 1 ) = d( 1 )
369  DO 60 j = 1, n - 1
370  jp = ( j*( j-1 ) ) / 2
371  jp1 = jp + j
372  IF( kband.EQ.1 ) THEN
373  work( jp1+j ) = ( cone-tau( j ) )*e( j )
374  DO 50 jr = 1, j - 1
375  work( jp1+jr ) = -tau( j )*e( j )*vp( jp1+jr )
376  50 CONTINUE
377  END IF
378 *
379  IF( tau( j ).NE.czero ) THEN
380  vsave = vp( jp1+j )
381  vp( jp1+j ) = cone
382  CALL zhpmv( 'U', j, cone, work, vp( jp1+1 ), 1, czero,
383  $ work( lap+1 ), 1 )
384  temp = -half*tau( j )*zdotc( j, work( lap+1 ), 1,
385  $ vp( jp1+1 ), 1 )
386  CALL zaxpy( j, temp, vp( jp1+1 ), 1, work( lap+1 ),
387  $ 1 )
388  CALL zhpr2( 'U', j, -tau( j ), vp( jp1+1 ), 1,
389  $ work( lap+1 ), 1, work )
390  vp( jp1+j ) = vsave
391  END IF
392  work( jp1+j+1 ) = d( j+1 )
393  60 CONTINUE
394  END IF
395 *
396  DO 70 j = 1, lap
397  work( j ) = work( j ) - ap( j )
398  70 CONTINUE
399  wnorm = zlanhp( '1', cuplo, n, work, rwork )
400 *
401  ELSE IF( itype.EQ.3 ) THEN
402 *
403 * ITYPE=3: error = U V* - I
404 *
405  IF( n.LT.2 )
406  $ RETURN
407  CALL zlacpy( ' ', n, n, u, ldu, work, n )
408  CALL zupmtr( 'R', cuplo, 'C', n, n, vp, tau, work, n,
409  $ work( n**2+1 ), iinfo )
410  IF( iinfo.NE.0 ) THEN
411  result( 1 ) = ten / ulp
412  RETURN
413  END IF
414 *
415  DO 80 j = 1, n
416  work( ( n+1 )*( j-1 )+1 ) = work( ( n+1 )*( j-1 )+1 ) - cone
417  80 CONTINUE
418 *
419  wnorm = zlange( '1', n, n, work, n, rwork )
420  END IF
421 *
422  IF( anorm.GT.wnorm ) THEN
423  result( 1 ) = ( wnorm / anorm ) / ( n*ulp )
424  ELSE
425  IF( anorm.LT.one ) THEN
426  result( 1 ) = ( min( wnorm, n*anorm ) / anorm ) / ( n*ulp )
427  ELSE
428  result( 1 ) = min( wnorm / anorm, dble( n ) ) / ( n*ulp )
429  END IF
430  END IF
431 *
432 * Do Test 2
433 *
434 * Compute UU* - I
435 *
436  IF( itype.EQ.1 ) THEN
437  CALL zgemm( 'N', 'C', n, n, n, cone, u, ldu, u, ldu, czero,
438  $ work, n )
439 *
440  DO 90 j = 1, n
441  work( ( n+1 )*( j-1 )+1 ) = work( ( n+1 )*( j-1 )+1 ) - cone
442  90 CONTINUE
443 *
444  result( 2 ) = min( zlange( '1', n, n, work, n, rwork ),
445  $ dble( n ) ) / ( n*ulp )
446  END IF
447 *
448  RETURN
449 *
450 * End of ZHPT21
451 *
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
subroutine zhpr2(UPLO, N, ALPHA, X, INCX, Y, INCY, AP)
ZHPR2
Definition: zhpr2.f:147
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
double precision function zlanhp(NORM, UPLO, N, AP, WORK)
ZLANHP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a complex Hermitian matrix supplied in packed form.
Definition: zlanhp.f:119
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:189
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: zlaset.f:108
complex *16 function zdotc(N, ZX, INCX, ZY, INCY)
ZDOTC
Definition: zdotc.f:54
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
subroutine zhpr(UPLO, N, ALPHA, X, INCX, AP)
ZHPR
Definition: zhpr.f:132
subroutine zupmtr(SIDE, UPLO, TRANS, M, N, AP, TAU, C, LDC, WORK, INFO)
ZUPMTR
Definition: zupmtr.f:152
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zhpmv(UPLO, N, ALPHA, AP, X, INCX, BETA, Y, INCY)
ZHPMV
Definition: zhpmv.f:151
subroutine zaxpy(N, ZA, ZX, INCX, ZY, INCY)
ZAXPY
Definition: zaxpy.f:53

Here is the call graph for this function:

Here is the caller graph for this function: