260 SUBROUTINE zlasyf_rk( UPLO, N, NB, KB, A, LDA, E, IPIV, W, LDW,
269 INTEGER INFO, KB, LDA, LDW, N, NB
273 COMPLEX*16 A( LDA, * ), E( * ), W( LDW, * )
279 DOUBLE PRECISION ZERO, ONE
280 parameter( zero = 0.0d+0, one = 1.0d+0 )
281 DOUBLE PRECISION EIGHT, SEVTEN
282 parameter( eight = 8.0d+0, sevten = 17.0d+0 )
283 COMPLEX*16 CONE, CZERO
284 parameter( cone = ( 1.0d+0, 0.0d+0 ),
285 $ czero = ( 0.0d+0, 0.0d+0 ) )
289 INTEGER IMAX, ITEMP, J, JB, JJ, JMAX, K, KK, KW, KKW,
291 DOUBLE PRECISION ABSAKK, ALPHA, COLMAX, ROWMAX, SFMIN, DTEMP
292 COMPLEX*16 D11, D12, D21, D22, R1, T, Z
297 DOUBLE PRECISION DLAMCH
298 EXTERNAL lsame, izamax, dlamch
304 INTRINSIC abs, dble, dimag, max, min, sqrt
307 DOUBLE PRECISION CABS1
310 cabs1( z ) = abs( dble( z ) ) + abs( dimag( z ) )
318 alpha = ( one+sqrt( sevten ) ) / eight
322 sfmin = dlamch(
'S' )
324 IF( lsame( uplo,
'U' ) )
THEN
346 IF( ( k.LE.n-nb+1 .AND. nb.LT.n ) .OR. k.LT.1 )
354 CALL zcopy( k, a( 1, k ), 1, w( 1, kw ), 1 )
356 $
CALL zgemv(
'No transpose', k, n-k, -cone, a( 1, k+1 ),
357 $ lda, w( k, kw+1 ), ldw, cone, w( 1, kw ), 1 )
362 absakk = cabs1( w( k, kw ) )
369 imax = izamax( k-1, w( 1, kw ), 1 )
370 colmax = cabs1( w( imax, kw ) )
375 IF( max( absakk, colmax ).EQ.zero )
THEN
382 CALL zcopy( k, w( 1, kw ), 1, a( 1, k ), 1 )
398 IF( .NOT.( absakk.LT.alpha*colmax ) )
THEN
417 CALL zcopy( imax, a( 1, imax ), 1, w( 1, kw-1 ), 1 )
418 CALL zcopy( k-imax, a( imax, imax+1 ), lda,
419 $ w( imax+1, kw-1 ), 1 )
422 $
CALL zgemv(
'No transpose', k, n-k, -cone,
423 $ a( 1, k+1 ), lda, w( imax, kw+1 ), ldw,
424 $ cone, w( 1, kw-1 ), 1 )
431 jmax = imax + izamax( k-imax, w( imax+1, kw-1 ),
433 rowmax = cabs1( w( jmax, kw-1 ) )
439 itemp = izamax( imax-1, w( 1, kw-1 ), 1 )
440 dtemp = cabs1( w( itemp, kw-1 ) )
441 IF( dtemp.GT.rowmax )
THEN
451 IF( .NOT.(cabs1( w( imax, kw-1 ) ).LT.alpha*rowmax ) )
461 CALL zcopy( k, w( 1, kw-1 ), 1, w( 1, kw ), 1 )
468 ELSE IF( ( p.EQ.jmax ) .OR. ( rowmax.LE.colmax ) )
487 CALL zcopy( k, w( 1, kw-1 ), 1, w( 1, kw ), 1 )
493 IF( .NOT. done )
GOTO 12
505 IF( ( kstep.EQ.2 ) .AND. ( p.NE.k ) )
THEN
509 CALL zcopy( k-p, a( p+1, k ), 1, a( p, p+1 ), lda )
510 CALL zcopy( p, a( 1, k ), 1, a( 1, p ), 1 )
515 CALL zswap( n-k+1, a( k, k ), lda, a( p, k ), lda )
516 CALL zswap( n-kk+1, w( k, kkw ), ldw, w( p, kkw ), ldw )
525 a( kp, k ) = a( kk, k )
526 CALL zcopy( k-1-kp, a( kp+1, kk ), 1, a( kp, kp+1 ),
528 CALL zcopy( kp, a( 1, kk ), 1, a( 1, kp ), 1 )
533 CALL zswap( n-kk+1, a( kk, kk ), lda, a( kp, kk ), lda )
534 CALL zswap( n-kk+1, w( kk, kkw ), ldw, w( kp, kkw ),
538 IF( kstep.EQ.1 )
THEN
548 CALL zcopy( k, w( 1, kw ), 1, a( 1, k ), 1 )
550 IF( cabs1( a( k, k ) ).GE.sfmin )
THEN
551 r1 = cone / a( k, k )
552 CALL zscal( k-1, r1, a( 1, k ), 1 )
553 ELSE IF( a( k, k ).NE.czero )
THEN
555 a( ii, k ) = a( ii, k ) / a( k, k )
580 d11 = w( k, kw ) / d12
581 d22 = w( k-1, kw-1 ) / d12
582 t = cone / ( d11*d22-cone )
584 a( j, k-1 ) = t*( (d11*w( j, kw-1 )-w( j, kw ) ) /
586 a( j, k ) = t*( ( d22*w( j, kw )-w( j, kw-1 ) ) /
595 a( k-1, k-1 ) = w( k-1, kw-1 )
597 a( k, k ) = w( k, kw )
598 e( k ) = w( k-1, kw )
609 IF( kstep.EQ.1 )
THEN
629 DO 50 j = ( ( k-1 ) / nb )*nb + 1, 1, -nb
630 jb = min( nb, k-j+1 )
634 DO 40 jj = j, j + jb - 1
635 CALL zgemv(
'No transpose', jj-j+1, n-k, -cone,
636 $ a( j, k+1 ), lda, w( jj, kw+1 ), ldw, cone,
643 $
CALL zgemm(
'No transpose',
'Transpose', j-1, jb,
644 $ n-k, -cone, a( 1, k+1 ), lda, w( j, kw+1 ),
645 $ ldw, cone, a( 1, j ), lda )
669 IF( ( k.GE.nb .AND. nb.LT.n ) .OR. k.GT.n )
677 CALL zcopy( n-k+1, a( k, k ), 1, w( k, k ), 1 )
679 $
CALL zgemv(
'No transpose', n-k+1, k-1, -cone, a( k, 1 ),
680 $ lda, w( k, 1 ), ldw, cone, w( k, k ), 1 )
685 absakk = cabs1( w( k, k ) )
692 imax = k + izamax( n-k, w( k+1, k ), 1 )
693 colmax = cabs1( w( imax, k ) )
698 IF( max( absakk, colmax ).EQ.zero )
THEN
705 CALL zcopy( n-k+1, w( k, k ), 1, a( k, k ), 1 )
721 IF( .NOT.( absakk.LT.alpha*colmax ) )
THEN
740 CALL zcopy( imax-k, a( imax, k ), lda, w( k, k+1 ), 1)
741 CALL zcopy( n-imax+1, a( imax, imax ), 1,
742 $ w( imax, k+1 ), 1 )
744 $
CALL zgemv(
'No transpose', n-k+1, k-1, -cone,
745 $ a( k, 1 ), lda, w( imax, 1 ), ldw,
746 $ cone, w( k, k+1 ), 1 )
753 jmax = k - 1 + izamax( imax-k, w( k, k+1 ), 1 )
754 rowmax = cabs1( w( jmax, k+1 ) )
760 itemp = imax + izamax( n-imax, w( imax+1, k+1 ), 1)
761 dtemp = cabs1( w( itemp, k+1 ) )
762 IF( dtemp.GT.rowmax )
THEN
772 IF( .NOT.( cabs1( w( imax, k+1 ) ).LT.alpha*rowmax ) )
782 CALL zcopy( n-k+1, w( k, k+1 ), 1, w( k, k ), 1 )
789 ELSE IF( ( p.EQ.jmax ) .OR. ( rowmax.LE.colmax ) )
808 CALL zcopy( n-k+1, w( k, k+1 ), 1, w( k, k ), 1 )
814 IF( .NOT. done )
GOTO 72
822 IF( ( kstep.EQ.2 ) .AND. ( p.NE.k ) )
THEN
826 CALL zcopy( p-k, a( k, k ), 1, a( p, k ), lda )
827 CALL zcopy( n-p+1, a( p, k ), 1, a( p, p ), 1 )
832 CALL zswap( k, a( k, 1 ), lda, a( p, 1 ), lda )
833 CALL zswap( kk, w( k, 1 ), ldw, w( p, 1 ), ldw )
842 a( kp, k ) = a( kk, k )
843 CALL zcopy( kp-k-1, a( k+1, kk ), 1, a( kp, k+1 ), lda )
844 CALL zcopy( n-kp+1, a( kp, kk ), 1, a( kp, kp ), 1 )
848 CALL zswap( kk, a( kk, 1 ), lda, a( kp, 1 ), lda )
849 CALL zswap( kk, w( kk, 1 ), ldw, w( kp, 1 ), ldw )
852 IF( kstep.EQ.1 )
THEN
862 CALL zcopy( n-k+1, w( k, k ), 1, a( k, k ), 1 )
864 IF( cabs1( a( k, k ) ).GE.sfmin )
THEN
865 r1 = cone / a( k, k )
866 CALL zscal( n-k, r1, a( k+1, k ), 1 )
867 ELSE IF( a( k, k ).NE.czero )
THEN
869 a( ii, k ) = a( ii, k ) / a( k, k )
893 d11 = w( k+1, k+1 ) / d21
894 d22 = w( k, k ) / d21
895 t = cone / ( d11*d22-cone )
897 a( j, k ) = t*( ( d11*w( j, k )-w( j, k+1 ) ) /
899 a( j, k+1 ) = t*( ( d22*w( j, k+1 )-w( j, k ) ) /
908 a( k, k ) = w( k, k )
910 a( k+1, k+1 ) = w( k+1, k+1 )
922 IF( kstep.EQ.1 )
THEN
943 jb = min( nb, n-j+1 )
947 DO 100 jj = j, j + jb - 1
948 CALL zgemv(
'No transpose', j+jb-jj, k-1, -cone,
949 $ a( jj, 1 ), lda, w( jj, 1 ), ldw, cone,
956 $
CALL zgemm(
'No transpose',
'Transpose', n-j-jb+1, jb,
957 $ k-1, -cone, a( j+jb, 1 ), lda, w( j, 1 ),
958 $ ldw, cone, a( j+jb, j ), lda )
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
subroutine zgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
ZGEMV
subroutine zlasyf_rk(uplo, n, nb, kb, a, lda, e, ipiv, w, ldw, info)
ZLASYF_RK computes a partial factorization of a complex symmetric indefinite matrix using bounded Bun...
subroutine zscal(n, za, zx, incx)
ZSCAL
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP