1 SUBROUTINE pzhentrd( UPLO, N, A, IA, JA, DESCA, D, E, TAU, WORK,
2 $ LWORK, RWORK, LRWORK, INFO )
11 INTEGER IA, INFO, JA, LRWORK, LWORK, N
15 DOUBLE PRECISION D( * ), E( * ), RWORK( * )
16 COMPLEX*16 A( * ), TAU( * ), WORK( * )
266 INTEGER BLOCK_CYCLIC_2D, DLEN_, DTYPE_, CTXT_, M_, N_,
267 $ mb_, nb_, rsrc_, csrc_, lld_
268 parameter( block_cyclic_2d = 1, dlen_ = 9, dtype_ = 1,
269 $ ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
270 $ rsrc_ = 7, csrc_ = 8, lld_ = 9 )
272 parameter( one = 1.0d+0 )
274 parameter( cone = ( 1.0d+0, 0.0d+0 ) )
277 LOGICAL LQUERY, UPPER
278 CHARACTER COLCTOP, ROWCTOP
279 INTEGER ANB, CTXTB, I, IACOL, IAROW, ICOFFA, ICTXT,
280 $ iinfo, indb, indrd, indre, indtau, indw, ipw,
281 $ iroffa, j, jb, jx, k, kk, llrwork, llwork,
282 $ lrwmin, lwmin, minsz, mycol, mycolb, myrow,
283 $ myrowb, nb, np, npcol, npcolb, nprow, nprowb,
284 $ nps, nq, onepmin, oneprmin, sqnpc, ttlrwmin,
288 INTEGER DESCB( DLEN_ ), DESCW( DLEN_ ), IDUM1( 3 ),
292 EXTERNAL blacs_get, blacs_gridexit, blacs_gridinfo,
300 INTEGER INDXG2L, INDXG2P, NUMROC, PJLAENV
301 EXTERNAL lsame, indxg2l, indxg2p, numroc, pjlaenv
304 INTRINSIC dble, dcmplx, ichar, int,
max,
min, mod, sqrt
309 IF( block_cyclic_2d*csrc_*ctxt_*dlen_*dtype_*lld_*mb_*m_*nb_*n_*
313 ictxt = desca( ctxt_ )
314 CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
319 IF( nprow.EQ.-1 )
THEN
320 info = -( 600+ctxt_ )
322 CALL chk1mat( n, 2, n, 2, ia, ja, desca, 6, info )
323 upper = lsame( uplo,
'U' )
326 iroffa = mod( ia-1, desca( mb_ ) )
327 icoffa = mod( ja-1, desca( nb_ ) )
328 iarow = indxg2p( ia, nb, myrow, desca( rsrc_ ), nprow )
329 iacol = indxg2p( ja, nb, mycol, desca( csrc_ ), npcol )
330 np = numroc( n, nb, myrow, iarow, nprow )
331 nq =
max( 1, numroc( n+ja-1, nb, mycol, desca( csrc_ ),
333 lwmin =
max( ( np+1 )*nb, 3*nb )
334 anb = pjlaenv( ictxt, 3,
'PZHETTRD',
'L', 0, 0, 0, 0 )
335 minsz = pjlaenv( ictxt, 5,
'PZHETTRD',
'L', 0, 0, 0, 0 )
336 sqnpc = int( sqrt( dble( nprow*npcol ) ) )
337 nps =
max( numroc( n, 1, 0, 0, sqnpc ), 2*anb )
338 ttlwmin = 2*( anb+1 )*( 4*nps+2 ) + ( nps+2 )*nps
342 work( 1 ) = dcmplx( dble( ttlwmin ) )
343 rwork( 1 ) = dble( ttlrwmin )
344 lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 )
345 IF( .NOT.upper .AND. .NOT.lsame( uplo,
'L' ) )
THEN
351 ELSE IF( iroffa.NE.icoffa .OR. icoffa.NE.0 )
THEN
353 ELSE IF( desca( mb_ ).NE.desca( nb_ ) )
THEN
355 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery )
THEN
357 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery )
THEN
362 idum1( 1 ) = ichar(
'U' )
364 idum1( 1 ) = ichar(
'L' )
367 IF( lwork.EQ.-1 )
THEN
373 IF( lrwork.EQ.-1 )
THEN
379 CALL pchk1mat( n, 2, n, 2, ia, ja, desca, 6, 3, idum1, idum2,
384 CALL pxerbla( ictxt,
'PZHENTRD', -info )
386 ELSE IF( lquery )
THEN
396 onepmin = n*n + 3*n + 1
398 CALL igamn2d( ictxt,
'A',
' ', 1, 1, llwork, 1, 1, -1, -1, -1,
403 CALL igamn2d( ictxt,
'A',
' ', 1, 1, llrwork, 1, 1, -1, -1, -1,
411 IF( ( n.LT.minsz .OR. sqnpc.EQ.1 ) .AND. llwork.GE.onepmin .AND.
412 $ llrwork.GE.oneprmin .AND. .NOT.upper )
THEN
416 IF( llwork.GE.ttlwmin .AND. llrwork.GE.ttlrwmin .AND. .NOT.
422 IF( nprowb.GE.1 )
THEN
428 indtau = indb + nps*nps
430 llwork = llwork - indw + 1
432 CALL blacs_get( ictxt, 10, ctxtb )
433 CALL blacs_gridinit( ctxtb,
'Row major', sqnpc, sqnpc )
434 CALL blacs_gridinfo( ctxtb, nprowb, npcolb, myrowb, mycolb )
435 CALL descset( descb, n, n, 1, 1, 0, 0, ctxtb, nps )
437 CALL pztrmr2d( uplo,
'N', n, n, a, ia, ja, desca, work( indb ),
438 $ 1, 1, descb, ictxt )
443 IF( nprowb.GT.0 )
THEN
445 IF( nprowb.EQ.1 )
THEN
446 CALL zhetrd( uplo, n, work( indb ), nps, rwork( indrd ),
447 $ rwork( indre ), work( indtau ),
448 $ work( indw ), llwork, info )
451 CALL pzhettrd(
'L', n, work( indb ), 1, 1, descb,
452 $ rwork( indrd ), rwork( indre ),
453 $ work( indtau ), work( indw ), llwork,
462 CALL pdlamr1d( n-1, rwork( indre ), 1, 1, descb, e, 1, ja,
465 CALL pdlamr1d( n, rwork( indrd ), 1, 1, descb, d, 1, ja,
468 CALL pzlamr1d( n, work( indtau ), 1, 1, descb, tau, 1, ja,
471 CALL pztrmr2d( uplo,
'N', n, n, work( indb ), 1, 1, descb, a,
472 $ ia, ja, desca, ictxt )
475 $
CALL blacs_gridexit( ctxtb )
479 CALL pb_topget( ictxt,
'Combine',
'Columnwise', colctop )
480 CALL pb_topget( ictxt,
'Combine',
'Rowwise', rowctop )
481 CALL pb_topset( ictxt,
'Combine',
'Columnwise',
'1-tree' )
482 CALL pb_topset( ictxt,
'Combine',
'Rowwise',
'1-tree' )
490 kk = mod( ja+n-1, nb )
493 CALL descset( descw, n, nb, nb, nb, iarow,
494 $ indxg2p( ja+n-kk, nb, mycol, desca( csrc_ ),
495 $ npcol ), ictxt,
max( 1, np ) )
497 DO 10 k = n - kk + 1, nb + 1, -nb
498 jb =
min( n-k+1, nb )
506 CALL pzlatrd( uplo, k+jb-1, jb, a, ia, ja, desca, d, e,
507 $ tau, work, 1, 1, descw, work( ipw ) )
513 CALL pzher2k( uplo,
'No transpose', k-1, jb, -cone, a,
514 $ ia, j, desca, work, 1, 1, descw, one, a,
519 jx =
min( indxg2l( j, nb, 0, iacol, npcol ), nq )
520 CALL pzelset( a, i-1, j, desca, dcmplx( e( jx ) ) )
522 descw( csrc_ ) = mod( descw( csrc_ )+npcol-1, npcol )
528 CALL pzhetd2( uplo,
min( n, nb ), a, ia, ja, desca, d, e,
529 $ tau, work, lwork, iinfo )
535 kk = mod( ja+n-1, nb )
538 CALL descset( descw, n, nb, nb, nb, iarow, iacol, ictxt,
541 DO 20 k = 1, n - nb, nb
549 CALL pzlatrd( uplo, n-k+1, nb, a, i, j, desca, d, e, tau,
550 $ work, k, 1, descw, work( ipw ) )
556 CALL pzher2k( uplo,
'No transpose', n-k-nb+1, nb, -cone,
557 $ a, i+nb, j, desca, work, k+nb, 1, descw,
558 $ one, a, i+nb, j+nb, desca )
562 jx =
min( indxg2l( j+nb-1, nb, 0, iacol, npcol ), nq )
563 CALL pzelset( a, i+nb, j+nb-1, desca, dcmplx( e( jx ) ) )
565 descw( csrc_ ) = mod( descw( csrc_ )+1, npcol )
571 CALL pzhetd2( uplo, kk, a, ia+k-1, ja+k-1, desca, d, e, tau,
572 $ work, lwork, iinfo )
575 CALL pb_topset( ictxt,
'Combine',
'Columnwise', colctop )
576 CALL pb_topset( ictxt,
'Combine',
'Rowwise', rowctop )
580 work( 1 ) = dcmplx( dble( ttlwmin ) )
581 rwork( 1 ) = dble( ttlrwmin )