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

◆ chetrd_he2hb()

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

CHETRD_HE2HB

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

Purpose:
!>
!> CHETRD_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
[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]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.
!> 
[in,out]A
!>          A is COMPLEX 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]AB
!>          AB is COMPLEX 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).
!> 
[in]LDAB
!>          LDAB is INTEGER
!>          The leading dimension of the array AB.  LDAB >= KD+1.
!> 
[out]TAU
!>          TAU is COMPLEX array, dimension (N-KD)
!>          The scalar factors of the elementary reflectors (see Further
!>          Details).
!> 
[out]WORK
!>          WORK is COMPLEX array, dimension (MAX(1,LWORK))
!>          On exit, if INFO = 0, or if LWORK = -1,
!>          WORK(1) returns the size of LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The dimension of the array WORK which should be calculated
!>          by a workspace query.
!>          If N <= KD+1, LWORK >= 1, else 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.
!> 
[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:
!>
!>  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).
!> 

Definition at line 241 of file chetrd_he2hb.f.

243*
244 IMPLICIT NONE
245*
246* -- LAPACK computational routine --
247* -- LAPACK is a software package provided by Univ. of Tennessee, --
248* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
249*
250* .. Scalar Arguments ..
251 CHARACTER UPLO
252 INTEGER INFO, LDA, LDAB, LWORK, N, KD
253* ..
254* .. Array Arguments ..
255 COMPLEX A( LDA, * ), AB( LDAB, * ),
256 $ TAU( * ), WORK( * )
257* ..
258*
259* =====================================================================
260*
261* .. Parameters ..
262 REAL RONE
263 COMPLEX ZERO, ONE, HALF
264 parameter( rone = 1.0e+0,
265 $ zero = ( 0.0e+0, 0.0e+0 ),
266 $ one = ( 1.0e+0, 0.0e+0 ),
267 $ half = ( 0.5e+0, 0.0e+0 ) )
268* ..
269* .. Local Scalars ..
270 LOGICAL LQUERY, UPPER
271 INTEGER I, J, IINFO, LWMIN, PN, PK, LK,
272 $ LDT, LDW, LDS2, LDS1,
273 $ LS2, LS1, LW, LT,
274 $ TPOS, WPOS, S2POS, S1POS
275* ..
276* .. External Subroutines ..
277 EXTERNAL xerbla, cher2k, chemm, cgemm,
278 $ ccopy,
280* ..
281* .. Intrinsic Functions ..
282 INTRINSIC min, max
283* ..
284* .. External Functions ..
285 LOGICAL LSAME
286 INTEGER ILAENV2STAGE
287 REAL SROUNDUP_LWORK
289* ..
290* .. Executable Statements ..
291*
292* Determine the minimal workspace size required
293* and test the input parameters
294*
295 info = 0
296 upper = lsame( uplo, 'U' )
297 lquery = ( lwork.EQ.-1 )
298 IF( n.LE.kd+1 ) THEN
299 lwmin = 1
300 ELSE
301 lwmin = ilaenv2stage( 4, 'CHETRD_HE2HB', '', n, kd, -1, -1 )
302 END IF
303*
304 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
305 info = -1
306 ELSE IF( n.LT.0 ) THEN
307 info = -2
308 ELSE IF( kd.LT.0 ) THEN
309 info = -3
310 ELSE IF( lda.LT.max( 1, n ) ) THEN
311 info = -5
312 ELSE IF( ldab.LT.max( 1, kd+1 ) ) THEN
313 info = -7
314 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
315 info = -10
316 END IF
317*
318 IF( info.NE.0 ) THEN
319 CALL xerbla( 'CHETRD_HE2HB', -info )
320 RETURN
321 ELSE IF( lquery ) THEN
322 work( 1 ) = sroundup_lwork( lwmin )
323 RETURN
324 END IF
325*
326* Quick return if possible
327* Copy the upper/lower portion of A into AB
328*
329 IF( n.LE.kd+1 ) THEN
330 IF( upper ) THEN
331 DO 100 i = 1, n
332 lk = min( kd+1, i )
333 CALL ccopy( lk, a( i-lk+1, i ), 1,
334 $ ab( kd+1-lk+1, i ), 1 )
335 100 CONTINUE
336 ELSE
337 DO 110 i = 1, n
338 lk = min( kd+1, n-i+1 )
339 CALL ccopy( lk, a( i, i ), 1, ab( 1, i ), 1 )
340 110 CONTINUE
341 ENDIF
342 work( 1 ) = 1
343 RETURN
344 END IF
345*
346* Determine the pointer position for the workspace
347*
348 ldt = kd
349 lds1 = kd
350 lt = ldt*kd
351 lw = n*kd
352 ls1 = lds1*kd
353 ls2 = lwmin - lt - lw - ls1
354* LS2 = N*MAX(KD,FACTOPTNB)
355 tpos = 1
356 wpos = tpos + lt
357 s1pos = wpos + lw
358 s2pos = s1pos + ls1
359 IF( upper ) THEN
360 ldw = kd
361 lds2 = kd
362 ELSE
363 ldw = n
364 lds2 = n
365 ENDIF
366*
367*
368* Set the workspace of the triangular matrix T to zero once such a
369* way every time T is generated the upper/lower portion will be always zero
370*
371 CALL claset( "A", ldt, kd, zero, zero, work( tpos ), ldt )
372*
373 IF( upper ) THEN
374 DO 10 i = 1, n - kd, kd
375 pn = n-i-kd+1
376 pk = min( n-i-kd+1, kd )
377*
378* Compute the LQ factorization of the current block
379*
380 CALL cgelqf( kd, pn, a( i, i+kd ), lda,
381 $ tau( i ), work( s2pos ), ls2, iinfo )
382*
383* Copy the upper portion of A into AB
384*
385 DO 20 j = i, i+pk-1
386 lk = min( kd, n-j ) + 1
387 CALL ccopy( lk, a( j, j ), lda, ab( kd+1, j ),
388 $ ldab-1 )
389 20 CONTINUE
390*
391 CALL claset( 'Lower', pk, pk, zero, one,
392 $ a( i, i+kd ), lda )
393*
394* Form the matrix T
395*
396 CALL clarft( 'Forward', 'Rowwise', pn, pk,
397 $ a( i, i+kd ), lda, tau( i ),
398 $ work( tpos ), ldt )
399*
400* Compute W:
401*
402 CALL cgemm( 'Conjugate', 'No transpose', pk, pn, pk,
403 $ one, work( tpos ), ldt,
404 $ a( i, i+kd ), lda,
405 $ zero, work( s2pos ), lds2 )
406*
407 CALL chemm( 'Right', uplo, pk, pn,
408 $ one, a( i+kd, i+kd ), lda,
409 $ work( s2pos ), lds2,
410 $ zero, work( wpos ), ldw )
411*
412 CALL cgemm( 'No transpose', 'Conjugate', pk, pk, pn,
413 $ one, work( wpos ), ldw,
414 $ work( s2pos ), lds2,
415 $ zero, work( s1pos ), lds1 )
416*
417 CALL cgemm( 'No transpose', 'No transpose', pk, pn, pk,
418 $ -half, work( s1pos ), lds1,
419 $ a( i, i+kd ), lda,
420 $ one, work( wpos ), ldw )
421*
422*
423* Update the unreduced submatrix A(i+kd:n,i+kd:n), using
424* an update of the form: A := A - V'*W - W'*V
425*
426 CALL cher2k( uplo, 'Conjugate', pn, pk,
427 $ -one, a( i, i+kd ), lda,
428 $ work( wpos ), ldw,
429 $ rone, a( i+kd, i+kd ), lda )
430 10 CONTINUE
431*
432* Copy the upper band to AB which is the band storage matrix
433*
434 DO 30 j = n-kd+1, n
435 lk = min(kd, n-j) + 1
436 CALL ccopy( lk, a( j, j ), lda, ab( kd+1, j ), ldab-1 )
437 30 CONTINUE
438*
439 ELSE
440*
441* Reduce the lower triangle of A to lower band matrix
442*
443 DO 40 i = 1, n - kd, kd
444 pn = n-i-kd+1
445 pk = min( n-i-kd+1, kd )
446*
447* Compute the QR factorization of the current block
448*
449 CALL cgeqrf( pn, kd, a( i+kd, i ), lda,
450 $ tau( i ), work( s2pos ), ls2, iinfo )
451*
452* Copy the upper portion of A into AB
453*
454 DO 50 j = i, i+pk-1
455 lk = min( kd, n-j ) + 1
456 CALL ccopy( lk, a( j, j ), 1, ab( 1, j ), 1 )
457 50 CONTINUE
458*
459 CALL claset( 'Upper', pk, pk, zero, one,
460 $ a( i+kd, i ), lda )
461*
462* Form the matrix T
463*
464 CALL clarft( 'Forward', 'Columnwise', pn, pk,
465 $ a( i+kd, i ), lda, tau( i ),
466 $ work( tpos ), ldt )
467*
468* Compute W:
469*
470 CALL cgemm( 'No transpose', 'No transpose', pn, pk, pk,
471 $ one, a( i+kd, i ), lda,
472 $ work( tpos ), ldt,
473 $ zero, work( s2pos ), lds2 )
474*
475 CALL chemm( 'Left', uplo, pn, pk,
476 $ one, a( i+kd, i+kd ), lda,
477 $ work( s2pos ), lds2,
478 $ zero, work( wpos ), ldw )
479*
480 CALL cgemm( 'Conjugate', 'No transpose', pk, pk, pn,
481 $ one, work( s2pos ), lds2,
482 $ work( wpos ), ldw,
483 $ zero, work( s1pos ), lds1 )
484*
485 CALL cgemm( 'No transpose', 'No transpose', pn, pk, pk,
486 $ -half, a( i+kd, i ), lda,
487 $ work( s1pos ), lds1,
488 $ one, work( wpos ), ldw )
489*
490*
491* Update the unreduced submatrix A(i+kd:n,i+kd:n), using
492* an update of the form: A := A - V*W' - W*V'
493*
494 CALL cher2k( uplo, 'No transpose', pn, pk,
495 $ -one, a( i+kd, i ), lda,
496 $ work( wpos ), ldw,
497 $ rone, a( i+kd, i+kd ), lda )
498* ==================================================================
499* RESTORE A FOR COMPARISON AND CHECKING TO BE REMOVED
500* DO 45 J = I, I+PK-1
501* LK = MIN( KD, N-J ) + 1
502* CALL CCOPY( LK, AB( 1, J ), 1, A( J, J ), 1 )
503* 45 CONTINUE
504* ==================================================================
505 40 CONTINUE
506*
507* Copy the lower band to AB which is the band storage matrix
508*
509 DO 60 j = n-kd+1, n
510 lk = min(kd, n-j) + 1
511 CALL ccopy( lk, a( j, j ), 1, ab( 1, j ), 1 )
512 60 CONTINUE
513
514 END IF
515*
516 work( 1 ) = sroundup_lwork( lwmin )
517 RETURN
518*
519* End of CHETRD_HE2HB
520*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
subroutine cgelqf(m, n, a, lda, tau, work, lwork, info)
CGELQF
Definition cgelqf.f:142
subroutine cgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CGEMM
Definition cgemm.f:188
subroutine cgeqrf(m, n, a, lda, tau, work, lwork, info)
CGEQRF
Definition cgeqrf.f:144
subroutine chemm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
CHEMM
Definition chemm.f:191
subroutine cher2k(uplo, trans, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CHER2K
Definition cher2k.f:197
integer function ilaenv2stage(ispec, name, opts, n1, n2, n3, n4)
ILAENV2STAGE
recursive subroutine clarft(direct, storev, n, k, v, ldv, tau, t, ldt)
CLARFT forms the triangular factor T of a block reflector H = I - vtvH
Definition clarft.f:162
subroutine claset(uplo, m, n, alpha, beta, a, lda)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition claset.f:104
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
Here is the call graph for this function:
Here is the caller graph for this function: