LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zlahef ( character  UPLO,
integer  N,
integer  NB,
integer  KB,
complex*16, dimension( lda, * )  A,
integer  LDA,
integer, dimension( * )  IPIV,
complex*16, dimension( ldw, * )  W,
integer  LDW,
integer  INFO 
)

ZLAHEF computes a partial factorization of a complex Hermitian indefinite matrix using the Bunch-Kaufman diagonal pivoting method (blocked algorithm, calling Level 3 BLAS).

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

Purpose:
 ZLAHEF computes a partial factorization of a complex Hermitian
 matrix A using the Bunch-Kaufman diagonal pivoting method. The
 partial factorization has the form:

 A  =  ( I  U12 ) ( A11  0  ) (  I      0     )  if UPLO = 'U', or:
       ( 0  U22 ) (  0   D  ) ( U12**H U22**H )

 A  =  ( L11  0 ) (  D   0  ) ( L11**H L21**H )  if UPLO = 'L'
       ( L21  I ) (  0  A22 ) (  0      I     )

 where the order of D is at most NB. The actual order is returned in
 the argument KB, and is either NB or NB-1, or N if N <= NB.
 Note that U**H denotes the conjugate transpose of U.

 ZLAHEF is an auxiliary routine called by ZHETRF. It uses blocked code
 (calling Level 3 BLAS) to update the submatrix A11 (if UPLO = 'U') or
 A22 (if UPLO = 'L').
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the upper or lower triangular part of the
          Hermitian matrix A is stored:
          = 'U':  Upper triangular
          = 'L':  Lower triangular
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]NB
          NB is INTEGER
          The maximum number of columns of the matrix A that should be
          factored.  NB should be at least 2 to allow for 2-by-2 pivot
          blocks.
[out]KB
          KB is INTEGER
          The number of columns of A that were actually factored.
          KB is either NB-1 or NB, or N if N <= NB.
[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, A contains details of the partial factorization.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]IPIV
          IPIV is INTEGER array, dimension (N)
          Details of the interchanges and the block structure of D.

          If UPLO = 'U':
             Only the last KB elements of IPIV are set.

             If IPIV(k) > 0, then rows and columns k and IPIV(k) were
             interchanged and D(k,k) is a 1-by-1 diagonal block.

             If IPIV(k) = IPIV(k-1) < 0, then rows and columns
             k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
             is a 2-by-2 diagonal block.

          If UPLO = 'L':
             Only the first KB elements of IPIV are set.

             If IPIV(k) > 0, then rows and columns k and IPIV(k) were
             interchanged and D(k,k) is a 1-by-1 diagonal block.

             If IPIV(k) = IPIV(k+1) < 0, then rows and columns
             k+1 and -IPIV(k) were interchanged and D(k:k+1,k:k+1)
             is a 2-by-2 diagonal block.
[out]W
          W is COMPLEX*16 array, dimension (LDW,NB)
[in]LDW
          LDW is INTEGER
          The leading dimension of the array W.  LDW >= max(1,N).
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          > 0: if INFO = k, D(k,k) is exactly zero.  The factorization
               has been completed, but the block diagonal matrix D is
               exactly singular.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2013
Contributors:
  November 2013,  Igor Kozachenko,
                  Computer Science Division,
                  University of California, Berkeley

Definition at line 179 of file zlahef.f.

179 *
180 * -- LAPACK computational routine (version 3.5.0) --
181 * -- LAPACK is a software package provided by Univ. of Tennessee, --
182 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
183 * November 2013
184 *
185 * .. Scalar Arguments ..
186  CHARACTER uplo
187  INTEGER info, kb, lda, ldw, n, nb
188 * ..
189 * .. Array Arguments ..
190  INTEGER ipiv( * )
191  COMPLEX*16 a( lda, * ), w( ldw, * )
192 * ..
193 *
194 * =====================================================================
195 *
196 * .. Parameters ..
197  DOUBLE PRECISION zero, one
198  parameter ( zero = 0.0d+0, one = 1.0d+0 )
199  COMPLEX*16 cone
200  parameter ( cone = ( 1.0d+0, 0.0d+0 ) )
201  DOUBLE PRECISION eight, sevten
202  parameter ( eight = 8.0d+0, sevten = 17.0d+0 )
203 * ..
204 * .. Local Scalars ..
205  INTEGER imax, j, jb, jj, jmax, jp, k, kk, kkw, kp,
206  $ kstep, kw
207  DOUBLE PRECISION absakk, alpha, colmax, r1, rowmax, t
208  COMPLEX*16 d11, d21, d22, z
209 * ..
210 * .. External Functions ..
211  LOGICAL lsame
212  INTEGER izamax
213  EXTERNAL lsame, izamax
214 * ..
215 * .. External Subroutines ..
216  EXTERNAL zcopy, zdscal, zgemm, zgemv, zlacgv, zswap
217 * ..
218 * .. Intrinsic Functions ..
219  INTRINSIC abs, dble, dconjg, dimag, max, min, sqrt
220 * ..
221 * .. Statement Functions ..
222  DOUBLE PRECISION cabs1
223 * ..
224 * .. Statement Function definitions ..
225  cabs1( z ) = abs( dble( z ) ) + abs( dimag( z ) )
226 * ..
227 * .. Executable Statements ..
228 *
229  info = 0
230 *
231 * Initialize ALPHA for use in choosing pivot block size.
232 *
233  alpha = ( one+sqrt( sevten ) ) / eight
234 *
235  IF( lsame( uplo, 'U' ) ) THEN
236 *
237 * Factorize the trailing columns of A using the upper triangle
238 * of A and working backwards, and compute the matrix W = U12*D
239 * for use in updating A11 (note that conjg(W) is actually stored)
240 *
241 * K is the main loop index, decreasing from N in steps of 1 or 2
242 *
243 * KW is the column of W which corresponds to column K of A
244 *
245  k = n
246  10 CONTINUE
247  kw = nb + k - n
248 *
249 * Exit from loop
250 *
251  IF( ( k.LE.n-nb+1 .AND. nb.LT.n ) .OR. k.LT.1 )
252  $ GO TO 30
253 *
254  kstep = 1
255 *
256 * Copy column K of A to column KW of W and update it
257 *
258  CALL zcopy( k-1, a( 1, k ), 1, w( 1, kw ), 1 )
259  w( k, kw ) = dble( a( k, k ) )
260  IF( k.LT.n ) THEN
261  CALL zgemv( 'No transpose', k, n-k, -cone, a( 1, k+1 ), lda,
262  $ w( k, kw+1 ), ldw, cone, w( 1, kw ), 1 )
263  w( k, kw ) = dble( w( k, kw ) )
264  END IF
265 *
266 * Determine rows and columns to be interchanged and whether
267 * a 1-by-1 or 2-by-2 pivot block will be used
268 *
269  absakk = abs( dble( w( k, kw ) ) )
270 *
271 * IMAX is the row-index of the largest off-diagonal element in
272 * column K, and COLMAX is its absolute value.
273 * Determine both COLMAX and IMAX.
274 *
275  IF( k.GT.1 ) THEN
276  imax = izamax( k-1, w( 1, kw ), 1 )
277  colmax = cabs1( w( imax, kw ) )
278  ELSE
279  colmax = zero
280  END IF
281 *
282  IF( max( absakk, colmax ).EQ.zero ) THEN
283 *
284 * Column K is zero or underflow: set INFO and continue
285 *
286  IF( info.EQ.0 )
287  $ info = k
288  kp = k
289  a( k, k ) = dble( a( k, k ) )
290  ELSE
291 *
292 * ============================================================
293 *
294 * BEGIN pivot search
295 *
296 * Case(1)
297  IF( absakk.GE.alpha*colmax ) THEN
298 *
299 * no interchange, use 1-by-1 pivot block
300 *
301  kp = k
302  ELSE
303 *
304 * BEGIN pivot search along IMAX row
305 *
306 *
307 * Copy column IMAX to column KW-1 of W and update it
308 *
309  CALL zcopy( imax-1, a( 1, imax ), 1, w( 1, kw-1 ), 1 )
310  w( imax, kw-1 ) = dble( a( imax, imax ) )
311  CALL zcopy( k-imax, a( imax, imax+1 ), lda,
312  $ w( imax+1, kw-1 ), 1 )
313  CALL zlacgv( k-imax, w( imax+1, kw-1 ), 1 )
314  IF( k.LT.n ) THEN
315  CALL zgemv( 'No transpose', k, n-k, -cone,
316  $ a( 1, k+1 ), lda, w( imax, kw+1 ), ldw,
317  $ cone, w( 1, kw-1 ), 1 )
318  w( imax, kw-1 ) = dble( w( imax, kw-1 ) )
319  END IF
320 *
321 * JMAX is the column-index of the largest off-diagonal
322 * element in row IMAX, and ROWMAX is its absolute value.
323 * Determine only ROWMAX.
324 *
325  jmax = imax + izamax( k-imax, w( imax+1, kw-1 ), 1 )
326  rowmax = cabs1( w( jmax, kw-1 ) )
327  IF( imax.GT.1 ) THEN
328  jmax = izamax( imax-1, w( 1, kw-1 ), 1 )
329  rowmax = max( rowmax, cabs1( w( jmax, kw-1 ) ) )
330  END IF
331 *
332 * Case(2)
333  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
334 *
335 * no interchange, use 1-by-1 pivot block
336 *
337  kp = k
338 *
339 * Case(3)
340  ELSE IF( abs( dble( w( imax, kw-1 ) ) ).GE.alpha*rowmax )
341  $ THEN
342 *
343 * interchange rows and columns K and IMAX, use 1-by-1
344 * pivot block
345 *
346  kp = imax
347 *
348 * copy column KW-1 of W to column KW of W
349 *
350  CALL zcopy( k, w( 1, kw-1 ), 1, w( 1, kw ), 1 )
351 *
352 * Case(4)
353  ELSE
354 *
355 * interchange rows and columns K-1 and IMAX, use 2-by-2
356 * pivot block
357 *
358  kp = imax
359  kstep = 2
360  END IF
361 *
362 *
363 * END pivot search along IMAX row
364 *
365  END IF
366 *
367 * END pivot search
368 *
369 * ============================================================
370 *
371 * KK is the column of A where pivoting step stopped
372 *
373  kk = k - kstep + 1
374 *
375 * KKW is the column of W which corresponds to column KK of A
376 *
377  kkw = nb + kk - n
378 *
379 * Interchange rows and columns KP and KK.
380 * Updated column KP is already stored in column KKW of W.
381 *
382  IF( kp.NE.kk ) THEN
383 *
384 * Copy non-updated column KK to column KP of submatrix A
385 * at step K. No need to copy element into column K
386 * (or K and K-1 for 2-by-2 pivot) of A, since these columns
387 * will be later overwritten.
388 *
389  a( kp, kp ) = dble( a( kk, kk ) )
390  CALL zcopy( kk-1-kp, a( kp+1, kk ), 1, a( kp, kp+1 ),
391  $ lda )
392  CALL zlacgv( kk-1-kp, a( kp, kp+1 ), lda )
393  IF( kp.GT.1 )
394  $ CALL zcopy( kp-1, a( 1, kk ), 1, a( 1, kp ), 1 )
395 *
396 * Interchange rows KK and KP in last K+1 to N columns of A
397 * (columns K (or K and K-1 for 2-by-2 pivot) of A will be
398 * later overwritten). Interchange rows KK and KP
399 * in last KKW to NB columns of W.
400 *
401  IF( k.LT.n )
402  $ CALL zswap( n-k, a( kk, k+1 ), lda, a( kp, k+1 ),
403  $ lda )
404  CALL zswap( n-kk+1, w( kk, kkw ), ldw, w( kp, kkw ),
405  $ ldw )
406  END IF
407 *
408  IF( kstep.EQ.1 ) THEN
409 *
410 * 1-by-1 pivot block D(k): column kw of W now holds
411 *
412 * W(kw) = U(k)*D(k),
413 *
414 * where U(k) is the k-th column of U
415 *
416 * (1) Store subdiag. elements of column U(k)
417 * and 1-by-1 block D(k) in column k of A.
418 * (NOTE: Diagonal element U(k,k) is a UNIT element
419 * and not stored)
420 * A(k,k) := D(k,k) = W(k,kw)
421 * A(1:k-1,k) := U(1:k-1,k) = W(1:k-1,kw)/D(k,k)
422 *
423 * (NOTE: No need to use for Hermitian matrix
424 * A( K, K ) = DBLE( W( K, K) ) to separately copy diagonal
425 * element D(k,k) from W (potentially saves only one load))
426  CALL zcopy( k, w( 1, kw ), 1, a( 1, k ), 1 )
427  IF( k.GT.1 ) THEN
428 *
429 * (NOTE: No need to check if A(k,k) is NOT ZERO,
430 * since that was ensured earlier in pivot search:
431 * case A(k,k) = 0 falls into 2x2 pivot case(4))
432 *
433  r1 = one / dble( a( k, k ) )
434  CALL zdscal( k-1, r1, a( 1, k ), 1 )
435 *
436 * (2) Conjugate column W(kw)
437 *
438  CALL zlacgv( k-1, w( 1, kw ), 1 )
439  END IF
440 *
441  ELSE
442 *
443 * 2-by-2 pivot block D(k): columns kw and kw-1 of W now hold
444 *
445 * ( W(kw-1) W(kw) ) = ( U(k-1) U(k) )*D(k)
446 *
447 * where U(k) and U(k-1) are the k-th and (k-1)-th columns
448 * of U
449 *
450 * (1) Store U(1:k-2,k-1) and U(1:k-2,k) and 2-by-2
451 * block D(k-1:k,k-1:k) in columns k-1 and k of A.
452 * (NOTE: 2-by-2 diagonal block U(k-1:k,k-1:k) is a UNIT
453 * block and not stored)
454 * A(k-1:k,k-1:k) := D(k-1:k,k-1:k) = W(k-1:k,kw-1:kw)
455 * A(1:k-2,k-1:k) := U(1:k-2,k:k-1:k) =
456 * = W(1:k-2,kw-1:kw) * ( D(k-1:k,k-1:k)**(-1) )
457 *
458  IF( k.GT.2 ) THEN
459 *
460 * Factor out the columns of the inverse of 2-by-2 pivot
461 * block D, so that each column contains 1, to reduce the
462 * number of FLOPS when we multiply panel
463 * ( W(kw-1) W(kw) ) by this inverse, i.e. by D**(-1).
464 *
465 * D**(-1) = ( d11 cj(d21) )**(-1) =
466 * ( d21 d22 )
467 *
468 * = 1/(d11*d22-|d21|**2) * ( ( d22) (-cj(d21) ) ) =
469 * ( (-d21) ( d11 ) )
470 *
471 * = 1/(|d21|**2) * 1/((d11/cj(d21))*(d22/d21)-1) *
472 *
473 * * ( d21*( d22/d21 ) conj(d21)*( - 1 ) ) =
474 * ( ( -1 ) ( d11/conj(d21) ) )
475 *
476 * = 1/(|d21|**2) * 1/(D22*D11-1) *
477 *
478 * * ( d21*( D11 ) conj(d21)*( -1 ) ) =
479 * ( ( -1 ) ( D22 ) )
480 *
481 * = (1/|d21|**2) * T * ( d21*( D11 ) conj(d21)*( -1 ) ) =
482 * ( ( -1 ) ( D22 ) )
483 *
484 * = ( (T/conj(d21))*( D11 ) (T/d21)*( -1 ) ) =
485 * ( ( -1 ) ( D22 ) )
486 *
487 * = ( conj(D21)*( D11 ) D21*( -1 ) )
488 * ( ( -1 ) ( D22 ) ),
489 *
490 * where D11 = d22/d21,
491 * D22 = d11/conj(d21),
492 * D21 = T/d21,
493 * T = 1/(D22*D11-1).
494 *
495 * (NOTE: No need to check for division by ZERO,
496 * since that was ensured earlier in pivot search:
497 * (a) d21 != 0, since in 2x2 pivot case(4)
498 * |d21| should be larger than |d11| and |d22|;
499 * (b) (D22*D11 - 1) != 0, since from (a),
500 * both |D11| < 1, |D22| < 1, hence |D22*D11| << 1.)
501 *
502  d21 = w( k-1, kw )
503  d11 = w( k, kw ) / dconjg( d21 )
504  d22 = w( k-1, kw-1 ) / d21
505  t = one / ( dble( d11*d22 )-one )
506  d21 = t / d21
507 *
508 * Update elements in columns A(k-1) and A(k) as
509 * dot products of rows of ( W(kw-1) W(kw) ) and columns
510 * of D**(-1)
511 *
512  DO 20 j = 1, k - 2
513  a( j, k-1 ) = d21*( d11*w( j, kw-1 )-w( j, kw ) )
514  a( j, k ) = dconjg( d21 )*
515  $ ( d22*w( j, kw )-w( j, kw-1 ) )
516  20 CONTINUE
517  END IF
518 *
519 * Copy D(k) to A
520 *
521  a( k-1, k-1 ) = w( k-1, kw-1 )
522  a( k-1, k ) = w( k-1, kw )
523  a( k, k ) = w( k, kw )
524 *
525 * (2) Conjugate columns W(kw) and W(kw-1)
526 *
527  CALL zlacgv( k-1, w( 1, kw ), 1 )
528  CALL zlacgv( k-2, w( 1, kw-1 ), 1 )
529 *
530  END IF
531 *
532  END IF
533 *
534 * Store details of the interchanges in IPIV
535 *
536  IF( kstep.EQ.1 ) THEN
537  ipiv( k ) = kp
538  ELSE
539  ipiv( k ) = -kp
540  ipiv( k-1 ) = -kp
541  END IF
542 *
543 * Decrease K and return to the start of the main loop
544 *
545  k = k - kstep
546  GO TO 10
547 *
548  30 CONTINUE
549 *
550 * Update the upper triangle of A11 (= A(1:k,1:k)) as
551 *
552 * A11 := A11 - U12*D*U12**H = A11 - U12*W**H
553 *
554 * computing blocks of NB columns at a time (note that conjg(W) is
555 * actually stored)
556 *
557  DO 50 j = ( ( k-1 ) / nb )*nb + 1, 1, -nb
558  jb = min( nb, k-j+1 )
559 *
560 * Update the upper triangle of the diagonal block
561 *
562  DO 40 jj = j, j + jb - 1
563  a( jj, jj ) = dble( a( jj, jj ) )
564  CALL zgemv( 'No transpose', jj-j+1, n-k, -cone,
565  $ a( j, k+1 ), lda, w( jj, kw+1 ), ldw, cone,
566  $ a( j, jj ), 1 )
567  a( jj, jj ) = dble( a( jj, jj ) )
568  40 CONTINUE
569 *
570 * Update the rectangular superdiagonal block
571 *
572  CALL zgemm( 'No transpose', 'Transpose', j-1, jb, n-k,
573  $ -cone, a( 1, k+1 ), lda, w( j, kw+1 ), ldw,
574  $ cone, a( 1, j ), lda )
575  50 CONTINUE
576 *
577 * Put U12 in standard form by partially undoing the interchanges
578 * in columns k+1:n looping backwards from k+1 to n
579 *
580  j = k + 1
581  60 CONTINUE
582 *
583 * Undo the interchanges (if any) of rows JJ and JP at each
584 * step J
585 *
586 * (Here, J is a diagonal index)
587  jj = j
588  jp = ipiv( j )
589  IF( jp.LT.0 ) THEN
590  jp = -jp
591 * (Here, J is a diagonal index)
592  j = j + 1
593  END IF
594 * (NOTE: Here, J is used to determine row length. Length N-J+1
595 * of the rows to swap back doesn't include diagonal element)
596  j = j + 1
597  IF( jp.NE.jj .AND. j.LE.n )
598  $ CALL zswap( n-j+1, a( jp, j ), lda, a( jj, j ), lda )
599  IF( j.LT.n )
600  $ GO TO 60
601 *
602 * Set KB to the number of columns factorized
603 *
604  kb = n - k
605 *
606  ELSE
607 *
608 * Factorize the leading columns of A using the lower triangle
609 * of A and working forwards, and compute the matrix W = L21*D
610 * for use in updating A22 (note that conjg(W) is actually stored)
611 *
612 * K is the main loop index, increasing from 1 in steps of 1 or 2
613 *
614  k = 1
615  70 CONTINUE
616 *
617 * Exit from loop
618 *
619  IF( ( k.GE.nb .AND. nb.LT.n ) .OR. k.GT.n )
620  $ GO TO 90
621 *
622  kstep = 1
623 *
624 * Copy column K of A to column K of W and update it
625 *
626  w( k, k ) = dble( a( k, k ) )
627  IF( k.LT.n )
628  $ CALL zcopy( n-k, a( k+1, k ), 1, w( k+1, k ), 1 )
629  CALL zgemv( 'No transpose', n-k+1, k-1, -cone, a( k, 1 ), lda,
630  $ w( k, 1 ), ldw, cone, w( k, k ), 1 )
631  w( k, k ) = dble( w( k, k ) )
632 *
633 * Determine rows and columns to be interchanged and whether
634 * a 1-by-1 or 2-by-2 pivot block will be used
635 *
636  absakk = abs( dble( w( k, k ) ) )
637 *
638 * IMAX is the row-index of the largest off-diagonal element in
639 * column K, and COLMAX is its absolute value.
640 * Determine both COLMAX and IMAX.
641 *
642  IF( k.LT.n ) THEN
643  imax = k + izamax( n-k, w( k+1, k ), 1 )
644  colmax = cabs1( w( imax, k ) )
645  ELSE
646  colmax = zero
647  END IF
648 *
649  IF( max( absakk, colmax ).EQ.zero ) THEN
650 *
651 * Column K is zero or underflow: set INFO and continue
652 *
653  IF( info.EQ.0 )
654  $ info = k
655  kp = k
656  a( k, k ) = dble( a( k, k ) )
657  ELSE
658 *
659 * ============================================================
660 *
661 * BEGIN pivot search
662 *
663 * Case(1)
664  IF( absakk.GE.alpha*colmax ) THEN
665 *
666 * no interchange, use 1-by-1 pivot block
667 *
668  kp = k
669  ELSE
670 *
671 * BEGIN pivot search along IMAX row
672 *
673 *
674 * Copy column IMAX to column K+1 of W and update it
675 *
676  CALL zcopy( imax-k, a( imax, k ), lda, w( k, k+1 ), 1 )
677  CALL zlacgv( imax-k, w( k, k+1 ), 1 )
678  w( imax, k+1 ) = dble( a( imax, imax ) )
679  IF( imax.LT.n )
680  $ CALL zcopy( n-imax, a( imax+1, imax ), 1,
681  $ w( imax+1, k+1 ), 1 )
682  CALL zgemv( 'No transpose', n-k+1, k-1, -cone, a( k, 1 ),
683  $ lda, w( imax, 1 ), ldw, cone, w( k, k+1 ),
684  $ 1 )
685  w( imax, k+1 ) = dble( w( imax, k+1 ) )
686 *
687 * JMAX is the column-index of the largest off-diagonal
688 * element in row IMAX, and ROWMAX is its absolute value.
689 * Determine only ROWMAX.
690 *
691  jmax = k - 1 + izamax( imax-k, w( k, k+1 ), 1 )
692  rowmax = cabs1( w( jmax, k+1 ) )
693  IF( imax.LT.n ) THEN
694  jmax = imax + izamax( n-imax, w( imax+1, k+1 ), 1 )
695  rowmax = max( rowmax, cabs1( w( jmax, k+1 ) ) )
696  END IF
697 *
698 * Case(2)
699  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
700 *
701 * no interchange, use 1-by-1 pivot block
702 *
703  kp = k
704 *
705 * Case(3)
706  ELSE IF( abs( dble( w( imax, k+1 ) ) ).GE.alpha*rowmax )
707  $ THEN
708 *
709 * interchange rows and columns K and IMAX, use 1-by-1
710 * pivot block
711 *
712  kp = imax
713 *
714 * copy column K+1 of W to column K of W
715 *
716  CALL zcopy( n-k+1, w( k, k+1 ), 1, w( k, k ), 1 )
717 *
718 * Case(4)
719  ELSE
720 *
721 * interchange rows and columns K+1 and IMAX, use 2-by-2
722 * pivot block
723 *
724  kp = imax
725  kstep = 2
726  END IF
727 *
728 *
729 * END pivot search along IMAX row
730 *
731  END IF
732 *
733 * END pivot search
734 *
735 * ============================================================
736 *
737 * KK is the column of A where pivoting step stopped
738 *
739  kk = k + kstep - 1
740 *
741 * Interchange rows and columns KP and KK.
742 * Updated column KP is already stored in column KK of W.
743 *
744  IF( kp.NE.kk ) THEN
745 *
746 * Copy non-updated column KK to column KP of submatrix A
747 * at step K. No need to copy element into column K
748 * (or K and K+1 for 2-by-2 pivot) of A, since these columns
749 * will be later overwritten.
750 *
751  a( kp, kp ) = dble( a( kk, kk ) )
752  CALL zcopy( kp-kk-1, a( kk+1, kk ), 1, a( kp, kk+1 ),
753  $ lda )
754  CALL zlacgv( kp-kk-1, a( kp, kk+1 ), lda )
755  IF( kp.LT.n )
756  $ CALL zcopy( n-kp, a( kp+1, kk ), 1, a( kp+1, kp ), 1 )
757 *
758 * Interchange rows KK and KP in first K-1 columns of A
759 * (columns K (or K and K+1 for 2-by-2 pivot) of A will be
760 * later overwritten). Interchange rows KK and KP
761 * in first KK columns of W.
762 *
763  IF( k.GT.1 )
764  $ CALL zswap( k-1, a( kk, 1 ), lda, a( kp, 1 ), lda )
765  CALL zswap( kk, w( kk, 1 ), ldw, w( kp, 1 ), ldw )
766  END IF
767 *
768  IF( kstep.EQ.1 ) THEN
769 *
770 * 1-by-1 pivot block D(k): column k of W now holds
771 *
772 * W(k) = L(k)*D(k),
773 *
774 * where L(k) is the k-th column of L
775 *
776 * (1) Store subdiag. elements of column L(k)
777 * and 1-by-1 block D(k) in column k of A.
778 * (NOTE: Diagonal element L(k,k) is a UNIT element
779 * and not stored)
780 * A(k,k) := D(k,k) = W(k,k)
781 * A(k+1:N,k) := L(k+1:N,k) = W(k+1:N,k)/D(k,k)
782 *
783 * (NOTE: No need to use for Hermitian matrix
784 * A( K, K ) = DBLE( W( K, K) ) to separately copy diagonal
785 * element D(k,k) from W (potentially saves only one load))
786  CALL zcopy( n-k+1, w( k, k ), 1, a( k, k ), 1 )
787  IF( k.LT.n ) THEN
788 *
789 * (NOTE: No need to check if A(k,k) is NOT ZERO,
790 * since that was ensured earlier in pivot search:
791 * case A(k,k) = 0 falls into 2x2 pivot case(4))
792 *
793  r1 = one / dble( a( k, k ) )
794  CALL zdscal( n-k, r1, a( k+1, k ), 1 )
795 *
796 * (2) Conjugate column W(k)
797 *
798  CALL zlacgv( n-k, w( k+1, k ), 1 )
799  END IF
800 *
801  ELSE
802 *
803 * 2-by-2 pivot block D(k): columns k and k+1 of W now hold
804 *
805 * ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k)
806 *
807 * where L(k) and L(k+1) are the k-th and (k+1)-th columns
808 * of L
809 *
810 * (1) Store L(k+2:N,k) and L(k+2:N,k+1) and 2-by-2
811 * block D(k:k+1,k:k+1) in columns k and k+1 of A.
812 * (NOTE: 2-by-2 diagonal block L(k:k+1,k:k+1) is a UNIT
813 * block and not stored)
814 * A(k:k+1,k:k+1) := D(k:k+1,k:k+1) = W(k:k+1,k:k+1)
815 * A(k+2:N,k:k+1) := L(k+2:N,k:k+1) =
816 * = W(k+2:N,k:k+1) * ( D(k:k+1,k:k+1)**(-1) )
817 *
818  IF( k.LT.n-1 ) THEN
819 *
820 * Factor out the columns of the inverse of 2-by-2 pivot
821 * block D, so that each column contains 1, to reduce the
822 * number of FLOPS when we multiply panel
823 * ( W(kw-1) W(kw) ) by this inverse, i.e. by D**(-1).
824 *
825 * D**(-1) = ( d11 cj(d21) )**(-1) =
826 * ( d21 d22 )
827 *
828 * = 1/(d11*d22-|d21|**2) * ( ( d22) (-cj(d21) ) ) =
829 * ( (-d21) ( d11 ) )
830 *
831 * = 1/(|d21|**2) * 1/((d11/cj(d21))*(d22/d21)-1) *
832 *
833 * * ( d21*( d22/d21 ) conj(d21)*( - 1 ) ) =
834 * ( ( -1 ) ( d11/conj(d21) ) )
835 *
836 * = 1/(|d21|**2) * 1/(D22*D11-1) *
837 *
838 * * ( d21*( D11 ) conj(d21)*( -1 ) ) =
839 * ( ( -1 ) ( D22 ) )
840 *
841 * = (1/|d21|**2) * T * ( d21*( D11 ) conj(d21)*( -1 ) ) =
842 * ( ( -1 ) ( D22 ) )
843 *
844 * = ( (T/conj(d21))*( D11 ) (T/d21)*( -1 ) ) =
845 * ( ( -1 ) ( D22 ) )
846 *
847 * = ( conj(D21)*( D11 ) D21*( -1 ) )
848 * ( ( -1 ) ( D22 ) ),
849 *
850 * where D11 = d22/d21,
851 * D22 = d11/conj(d21),
852 * D21 = T/d21,
853 * T = 1/(D22*D11-1).
854 *
855 * (NOTE: No need to check for division by ZERO,
856 * since that was ensured earlier in pivot search:
857 * (a) d21 != 0, since in 2x2 pivot case(4)
858 * |d21| should be larger than |d11| and |d22|;
859 * (b) (D22*D11 - 1) != 0, since from (a),
860 * both |D11| < 1, |D22| < 1, hence |D22*D11| << 1.)
861 *
862  d21 = w( k+1, k )
863  d11 = w( k+1, k+1 ) / d21
864  d22 = w( k, k ) / dconjg( d21 )
865  t = one / ( dble( d11*d22 )-one )
866  d21 = t / d21
867 *
868 * Update elements in columns A(k) and A(k+1) as
869 * dot products of rows of ( W(k) W(k+1) ) and columns
870 * of D**(-1)
871 *
872  DO 80 j = k + 2, n
873  a( j, k ) = dconjg( d21 )*
874  $ ( d11*w( j, k )-w( j, k+1 ) )
875  a( j, k+1 ) = d21*( d22*w( j, k+1 )-w( j, k ) )
876  80 CONTINUE
877  END IF
878 *
879 * Copy D(k) to A
880 *
881  a( k, k ) = w( k, k )
882  a( k+1, k ) = w( k+1, k )
883  a( k+1, k+1 ) = w( k+1, k+1 )
884 *
885 * (2) Conjugate columns W(k) and W(k+1)
886 *
887  CALL zlacgv( n-k, w( k+1, k ), 1 )
888  CALL zlacgv( n-k-1, w( k+2, k+1 ), 1 )
889 *
890  END IF
891 *
892  END IF
893 *
894 * Store details of the interchanges in IPIV
895 *
896  IF( kstep.EQ.1 ) THEN
897  ipiv( k ) = kp
898  ELSE
899  ipiv( k ) = -kp
900  ipiv( k+1 ) = -kp
901  END IF
902 *
903 * Increase K and return to the start of the main loop
904 *
905  k = k + kstep
906  GO TO 70
907 *
908  90 CONTINUE
909 *
910 * Update the lower triangle of A22 (= A(k:n,k:n)) as
911 *
912 * A22 := A22 - L21*D*L21**H = A22 - L21*W**H
913 *
914 * computing blocks of NB columns at a time (note that conjg(W) is
915 * actually stored)
916 *
917  DO 110 j = k, n, nb
918  jb = min( nb, n-j+1 )
919 *
920 * Update the lower triangle of the diagonal block
921 *
922  DO 100 jj = j, j + jb - 1
923  a( jj, jj ) = dble( a( jj, jj ) )
924  CALL zgemv( 'No transpose', j+jb-jj, k-1, -cone,
925  $ a( jj, 1 ), lda, w( jj, 1 ), ldw, cone,
926  $ a( jj, jj ), 1 )
927  a( jj, jj ) = dble( a( jj, jj ) )
928  100 CONTINUE
929 *
930 * Update the rectangular subdiagonal block
931 *
932  IF( j+jb.LE.n )
933  $ CALL zgemm( 'No transpose', 'Transpose', n-j-jb+1, jb,
934  $ k-1, -cone, a( j+jb, 1 ), lda, w( j, 1 ),
935  $ ldw, cone, a( j+jb, j ), lda )
936  110 CONTINUE
937 *
938 * Put L21 in standard form by partially undoing the interchanges
939 * of rows in columns 1:k-1 looping backwards from k-1 to 1
940 *
941  j = k - 1
942  120 CONTINUE
943 *
944 * Undo the interchanges (if any) of rows JJ and JP at each
945 * step J
946 *
947 * (Here, J is a diagonal index)
948  jj = j
949  jp = ipiv( j )
950  IF( jp.LT.0 ) THEN
951  jp = -jp
952 * (Here, J is a diagonal index)
953  j = j - 1
954  END IF
955 * (NOTE: Here, J is used to determine row length. Length J
956 * of the rows to swap back doesn't include diagonal element)
957  j = j - 1
958  IF( jp.NE.jj .AND. j.GE.1 )
959  $ CALL zswap( j, a( jp, 1 ), lda, a( jj, 1 ), lda )
960  IF( j.GT.1 )
961  $ GO TO 120
962 *
963 * Set KB to the number of columns factorized
964 *
965  kb = k - 1
966 *
967  END IF
968  RETURN
969 *
970 * End of ZLAHEF
971 *
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:52
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:189
integer function izamax(N, ZX, INCX)
IZAMAX
Definition: izamax.f:53
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zlacgv(N, X, INCX)
ZLACGV conjugates a complex vector.
Definition: zlacgv.f:76

Here is the call graph for this function:

Here is the caller graph for this function: