Scroll to navigation

/home/abuild/rpmbuild/BUILD/lapack-3.12.0/SRC/zhetrd_he2hb.f(3) Library Functions Manual /home/abuild/rpmbuild/BUILD/lapack-3.12.0/SRC/zhetrd_he2hb.f(3)

NAME

/home/abuild/rpmbuild/BUILD/lapack-3.12.0/SRC/zhetrd_he2hb.f

SYNOPSIS

Functions/Subroutines


subroutine ZHETRD_HE2HB (uplo, n, kd, a, lda, ab, ldab, tau, work, lwork, info)
ZHETRD_HE2HB

Function/Subroutine Documentation

subroutine ZHETRD_HE2HB (character uplo, integer n, integer kd, complex*16, dimension( lda, * ) a, integer lda, complex*16, dimension( ldab, * ) ab, integer ldab, complex*16, dimension( * ) tau, complex*16, dimension( * ) work, integer lwork, integer info)

ZHETRD_HE2HB

Purpose:

!>
!> ZHETRD_HE2HB reduces a complex Hermitian matrix A to complex Hermitian
!> band-diagonal form AB by a unitary similarity transformation:
!> Q**H * A * Q = AB.
!> 

Parameters

UPLO

!>          UPLO is CHARACTER*1
!>          = 'U':  Upper triangle of A is stored;
!>          = 'L':  Lower triangle of A is stored.
!> 

N

!>          N is INTEGER
!>          The order of the matrix A.  N >= 0.
!> 

KD

!>          KD is INTEGER
!>          The number of superdiagonals of the reduced matrix if UPLO = 'U',
!>          or the number of subdiagonals if UPLO = 'L'.  KD >= 0.
!>          The reduced matrix is stored in the array AB.
!> 

A

!>          A is COMPLEX*16 array, dimension (LDA,N)
!>          On entry, the Hermitian matrix A.  If UPLO = 'U', the leading
!>          N-by-N upper triangular part of A contains the upper
!>          triangular part of the matrix A, and the strictly lower
!>          triangular part of A is not referenced.  If UPLO = 'L', the
!>          leading N-by-N lower triangular part of A contains the lower
!>          triangular part of the matrix A, and the strictly upper
!>          triangular part of A is not referenced.
!>          On exit, if UPLO = 'U', the diagonal and first superdiagonal
!>          of A are overwritten by the corresponding elements of the
!>          tridiagonal matrix T, and the elements above the first
!>          superdiagonal, with the array TAU, represent the unitary
!>          matrix Q as a product of elementary reflectors; if UPLO
!>          = 'L', the diagonal and first subdiagonal of A are over-
!>          written by the corresponding elements of the tridiagonal
!>          matrix T, and the elements below the first subdiagonal, with
!>          the array TAU, represent the unitary matrix Q as a product
!>          of elementary reflectors. See Further Details.
!> 

LDA

!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 

AB

!>          AB is COMPLEX*16 array, dimension (LDAB,N)
!>          On exit, the upper or lower triangle of the Hermitian band
!>          matrix A, stored in the first KD+1 rows of the array.  The
!>          j-th column of A is stored in the j-th column of the array AB
!>          as follows:
!>          if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
!>          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+kd).
!> 

LDAB

!>          LDAB is INTEGER
!>          The leading dimension of the array AB.  LDAB >= KD+1.
!> 

TAU

!>          TAU is COMPLEX*16 array, dimension (N-KD)
!>          The scalar factors of the elementary reflectors (see Further
!>          Details).
!> 

WORK

!>          WORK is COMPLEX*16 array, dimension (LWORK)
!>          On exit, if INFO = 0, or if LWORK=-1, 
!>          WORK(1) returns the size of LWORK.
!> 

LWORK

!>          LWORK is INTEGER
!>          The dimension of the array WORK which should be calculated
!>          by a workspace query. LWORK = MAX(1, LWORK_QUERY)
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!>          LWORK_QUERY = N*KD + N*max(KD,FACTOPTNB) + 2*KD*KD
!>          where FACTOPTNB is the blocking used by the QR or LQ
!>          algorithm, usually FACTOPTNB=128 is a good choice otherwise
!>          putting LWORK=-1 will provide the size of WORK.
!> 

INFO

!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value
!> 

Author

Univ. of Tennessee

Univ. of California Berkeley

Univ. of Colorado Denver

NAG Ltd.

Further Details:

!>
!>  Implemented by Azzam Haidar.
!>
!>  All details are available on technical report, SC11, SC13 papers.
!>
!>  Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
!>  Parallel reduction to condensed forms for symmetric eigenvalue problems
!>  using aggregated fine-grained and memory-aware kernels. In Proceedings
!>  of 2011 International Conference for High Performance Computing,
!>  Networking, Storage and Analysis (SC '11), New York, NY, USA,
!>  Article 8 , 11 pages.
!>  http://doi.acm.org/10.1145/2063384.2063394
!>
!>  A. Haidar, J. Kurzak, P. Luszczek, 2013.
!>  An improved parallel singular value algorithm and its implementation 
!>  for multicore hardware, In Proceedings of 2013 International Conference
!>  for High Performance Computing, Networking, Storage and Analysis (SC '13).
!>  Denver, Colorado, USA, 2013.
!>  Article 90, 12 pages.
!>  http://doi.acm.org/10.1145/2503210.2503292
!>
!>  A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
!>  A novel hybrid CPU-GPU generalized eigensolver for electronic structure 
!>  calculations based on fine-grained memory aware tasks.
!>  International Journal of High Performance Computing Applications.
!>  Volume 28 Issue 2, Pages 196-209, May 2014.
!>  http://hpc.sagepub.com/content/28/2/196 
!>
!> 

!>
!>  If UPLO = 'U', the matrix Q is represented as a product of elementary
!>  reflectors
!>
!>     Q = H(k)**H . . . H(2)**H H(1)**H, where k = n-kd.
!>
!>  Each H(i) has the form
!>
!>     H(i) = I - tau * v * v**H
!>
!>  where tau is a complex scalar, and v is a complex vector with
!>  v(1:i+kd-1) = 0 and v(i+kd) = 1; conjg(v(i+kd+1:n)) is stored on exit in
!>  A(i,i+kd+1:n), and tau in TAU(i).
!>
!>  If UPLO = 'L', the matrix Q is represented as a product of elementary
!>  reflectors
!>
!>     Q = H(1) H(2) . . . H(k), where k = n-kd.
!>
!>  Each H(i) has the form
!>
!>     H(i) = I - tau * v * v**H
!>
!>  where tau is a complex scalar, and v is a complex vector with
!>  v(kd+1:i) = 0 and v(i+kd+1) = 1; v(i+kd+2:n) is stored on exit in
!>  A(i+kd+2:n,i), and tau in TAU(i).
!>
!>  The contents of A on exit are illustrated by the following examples
!>  with n = 5:
!>
!>  if UPLO = 'U':                       if UPLO = 'L':
!>
!>    (  ab  ab/v1  v1      v1     v1    )              (  ab                            )
!>    (      ab     ab/v2   v2     v2    )              (  ab/v1  ab                     )
!>    (             ab      ab/v3  v3    )              (  v1     ab/v2  ab              )
!>    (                     ab     ab/v4 )              (  v1     v2     ab/v3  ab       )
!>    (                            ab    )              (  v1     v2     v3     ab/v4 ab )
!>
!>  where d and e denote diagonal and off-diagonal elements of T, and vi
!>  denotes an element of the vector defining H(i).
!> .fi

Definition at line 241 of file zhetrd_he2hb.f.

Author

Generated automatically by Doxygen for LAPACK from the source code.

Version 3.12.0 LAPACK