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

◆ zhetrd()

subroutine zhetrd ( character  uplo,
integer  n,
complex*16, dimension( lda, * )  a,
integer  lda,
double precision, dimension( * )  d,
double precision, dimension( * )  e,
complex*16, dimension( * )  tau,
complex*16, dimension( * )  work,
integer  lwork,
integer  info 
)

ZHETRD

Download ZHETRD + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 ZHETRD reduces a complex Hermitian matrix A to real symmetric
 tridiagonal form T by a unitary similarity transformation:
 Q**H * A * Q = T.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]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.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]D
          D is DOUBLE PRECISION array, dimension (N)
          The diagonal elements of the tridiagonal matrix T:
          D(i) = A(i,i).
[out]E
          E is DOUBLE PRECISION array, dimension (N-1)
          The off-diagonal elements of the tridiagonal matrix T:
          E(i) = A(i,i+1) if UPLO = 'U', E(i) = A(i+1,i) if UPLO = 'L'.
[out]TAU
          TAU is COMPLEX*16 array, dimension (N-1)
          The scalar factors of the elementary reflectors (see Further
          Details).
[out]WORK
          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.  LWORK >= 1.
          For optimum performance LWORK >= N*NB, where NB is the
          optimal blocksize.

          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.
[out]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:
  If UPLO = 'U', the matrix Q is represented as a product of elementary
  reflectors

     Q = H(n-1) . . . H(2) H(1).

  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(i+1:n) = 0 and v(i) = 1; v(1:i-1) is stored on exit in
  A(1:i-1,i+1), 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(n-1).

  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) = 0 and v(i+1) = 1; v(i+2:n) is stored on exit in A(i+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':

    (  d   e   v2  v3  v4 )              (  d                  )
    (      d   e   v3  v4 )              (  e   d              )
    (          d   e   v4 )              (  v1  e   d          )
    (              d   e  )              (  v1  v2  e   d      )
    (                  d  )              (  v1  v2  v3  e   d  )

  where d and e denote diagonal and off-diagonal elements of T, and vi
  denotes an element of the vector defining H(i).

Definition at line 191 of file zhetrd.f.

192*
193* -- LAPACK computational routine --
194* -- LAPACK is a software package provided by Univ. of Tennessee, --
195* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
196*
197* .. Scalar Arguments ..
198 CHARACTER UPLO
199 INTEGER INFO, LDA, LWORK, N
200* ..
201* .. Array Arguments ..
202 DOUBLE PRECISION D( * ), E( * )
203 COMPLEX*16 A( LDA, * ), TAU( * ), WORK( * )
204* ..
205*
206* =====================================================================
207*
208* .. Parameters ..
209 DOUBLE PRECISION ONE
210 parameter( one = 1.0d+0 )
211 COMPLEX*16 CONE
212 parameter( cone = ( 1.0d+0, 0.0d+0 ) )
213* ..
214* .. Local Scalars ..
215 LOGICAL LQUERY, UPPER
216 INTEGER I, IINFO, IWS, J, KK, LDWORK, LWKOPT, NB,
217 $ NBMIN, NX
218* ..
219* .. External Subroutines ..
220 EXTERNAL xerbla, zher2k, zhetd2, zlatrd
221* ..
222* .. Intrinsic Functions ..
223 INTRINSIC max
224* ..
225* .. External Functions ..
226 LOGICAL LSAME
227 INTEGER ILAENV
228 EXTERNAL lsame, ilaenv
229* ..
230* .. Executable Statements ..
231*
232* Test the input parameters
233*
234 info = 0
235 upper = lsame( uplo, 'U' )
236 lquery = ( lwork.EQ.-1 )
237 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
238 info = -1
239 ELSE IF( n.LT.0 ) THEN
240 info = -2
241 ELSE IF( lda.LT.max( 1, n ) ) THEN
242 info = -4
243 ELSE IF( lwork.LT.1 .AND. .NOT.lquery ) THEN
244 info = -9
245 END IF
246*
247 IF( info.EQ.0 ) THEN
248*
249* Determine the block size.
250*
251 nb = ilaenv( 1, 'ZHETRD', uplo, n, -1, -1, -1 )
252 lwkopt = n*nb
253 work( 1 ) = lwkopt
254 END IF
255*
256 IF( info.NE.0 ) THEN
257 CALL xerbla( 'ZHETRD', -info )
258 RETURN
259 ELSE IF( lquery ) THEN
260 RETURN
261 END IF
262*
263* Quick return if possible
264*
265 IF( n.EQ.0 ) THEN
266 work( 1 ) = 1
267 RETURN
268 END IF
269*
270 nx = n
271 iws = 1
272 IF( nb.GT.1 .AND. nb.LT.n ) THEN
273*
274* Determine when to cross over from blocked to unblocked code
275* (last block is always handled by unblocked code).
276*
277 nx = max( nb, ilaenv( 3, 'ZHETRD', uplo, n, -1, -1, -1 ) )
278 IF( nx.LT.n ) THEN
279*
280* Determine if workspace is large enough for blocked code.
281*
282 ldwork = n
283 iws = ldwork*nb
284 IF( lwork.LT.iws ) THEN
285*
286* Not enough workspace to use optimal NB: determine the
287* minimum value of NB, and reduce NB or force use of
288* unblocked code by setting NX = N.
289*
290 nb = max( lwork / ldwork, 1 )
291 nbmin = ilaenv( 2, 'ZHETRD', uplo, n, -1, -1, -1 )
292 IF( nb.LT.nbmin )
293 $ nx = n
294 END IF
295 ELSE
296 nx = n
297 END IF
298 ELSE
299 nb = 1
300 END IF
301*
302 IF( upper ) THEN
303*
304* Reduce the upper triangle of A.
305* Columns 1:kk are handled by the unblocked method.
306*
307 kk = n - ( ( n-nx+nb-1 ) / nb )*nb
308 DO 20 i = n - nb + 1, kk + 1, -nb
309*
310* Reduce columns i:i+nb-1 to tridiagonal form and form the
311* matrix W which is needed to update the unreduced part of
312* the matrix
313*
314 CALL zlatrd( uplo, i+nb-1, nb, a, lda, e, tau, work,
315 $ ldwork )
316*
317* Update the unreduced submatrix A(1:i-1,1:i-1), using an
318* update of the form: A := A - V*W**H - W*V**H
319*
320 CALL zher2k( uplo, 'No transpose', i-1, nb, -cone,
321 $ a( 1, i ), lda, work, ldwork, one, a, lda )
322*
323* Copy superdiagonal elements back into A, and diagonal
324* elements into D
325*
326 DO 10 j = i, i + nb - 1
327 a( j-1, j ) = e( j-1 )
328 d( j ) = dble( a( j, j ) )
329 10 CONTINUE
330 20 CONTINUE
331*
332* Use unblocked code to reduce the last or only block
333*
334 CALL zhetd2( uplo, kk, a, lda, d, e, tau, iinfo )
335 ELSE
336*
337* Reduce the lower triangle of A
338*
339 DO 40 i = 1, n - nx, nb
340*
341* Reduce columns i:i+nb-1 to tridiagonal form and form the
342* matrix W which is needed to update the unreduced part of
343* the matrix
344*
345 CALL zlatrd( uplo, n-i+1, nb, a( i, i ), lda, e( i ),
346 $ tau( i ), work, ldwork )
347*
348* Update the unreduced submatrix A(i+nb:n,i+nb:n), using
349* an update of the form: A := A - V*W**H - W*V**H
350*
351 CALL zher2k( uplo, 'No transpose', n-i-nb+1, nb, -cone,
352 $ a( i+nb, i ), lda, work( nb+1 ), ldwork, one,
353 $ a( i+nb, i+nb ), lda )
354*
355* Copy subdiagonal elements back into A, and diagonal
356* elements into D
357*
358 DO 30 j = i, i + nb - 1
359 a( j+1, j ) = e( j )
360 d( j ) = dble( a( j, j ) )
361 30 CONTINUE
362 40 CONTINUE
363*
364* Use unblocked code to reduce the last or only block
365*
366 CALL zhetd2( uplo, n-i+1, a( i, i ), lda, d( i ), e( i ),
367 $ tau( i ), iinfo )
368 END IF
369*
370 work( 1 ) = lwkopt
371 RETURN
372*
373* End of ZHETRD
374*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zher2k(uplo, trans, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZHER2K
Definition zher2k.f:198
subroutine zhetd2(uplo, n, a, lda, d, e, tau, info)
ZHETD2 reduces a Hermitian matrix to real symmetric tridiagonal form by an unitary similarity transfo...
Definition zhetd2.f:175
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine zlatrd(uplo, n, nb, a, lda, e, tau, w, ldw)
ZLATRD reduces the first nb rows and columns of a symmetric/Hermitian matrix A to real tridiagonal fo...
Definition zlatrd.f:199
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: