LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zlasyf ( 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 
)

ZLASYF computes a partial factorization of a complex symmetric matrix using the Bunch-Kaufman diagonal pivoting method.

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

Purpose:
 ZLASYF computes a partial factorization of a complex symmetric 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**T U22**T )

 A  =  ( L11  0 ) ( D    0  ) ( L11**T L21**T )  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**T denotes the transpose of U.

 ZLASYF is an auxiliary routine called by ZSYTRF. 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
          symmetric 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 symmetric 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 zlasyf.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  DOUBLE PRECISION eight, sevten
200  parameter ( eight = 8.0d+0, sevten = 17.0d+0 )
201  COMPLEX*16 cone
202  parameter ( cone = ( 1.0d+0, 0.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, rowmax
208  COMPLEX*16 d11, d21, d22, r1, t, z
209 * ..
210 * .. External Functions ..
211  LOGICAL lsame
212  INTEGER izamax
213  EXTERNAL lsame, izamax
214 * ..
215 * .. External Subroutines ..
216  EXTERNAL zcopy, zgemm, zgemv, zscal, zswap
217 * ..
218 * .. Intrinsic Functions ..
219  INTRINSIC abs, dble, 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
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 * Copy column K of A to column KW of W and update it
255 *
256  CALL zcopy( k, a( 1, k ), 1, w( 1, kw ), 1 )
257  IF( k.LT.n )
258  $ CALL zgemv( 'No transpose', k, n-k, -cone, a( 1, k+1 ), lda,
259  $ w( k, kw+1 ), ldw, cone, w( 1, kw ), 1 )
260 *
261  kstep = 1
262 *
263 * Determine rows and columns to be interchanged and whether
264 * a 1-by-1 or 2-by-2 pivot block will be used
265 *
266  absakk = cabs1( w( k, kw ) )
267 *
268 * IMAX is the row-index of the largest off-diagonal element in
269 
270 *
271  IF( k.GT.1 ) THEN
272  imax = izamax( k-1, w( 1, kw ), 1 )
273  colmax = cabs1( w( imax, kw ) )
274  ELSE
275  colmax = zero
276  END IF
277 *
278  IF( max( absakk, colmax ).EQ.zero ) THEN
279 *
280 * Column K is zero or underflow: set INFO and continue
281 *
282  IF( info.EQ.0 )
283  $ info = k
284  kp = k
285  ELSE
286  IF( absakk.GE.alpha*colmax ) THEN
287 *
288 * no interchange, use 1-by-1 pivot block
289 *
290  kp = k
291  ELSE
292 *
293 * Copy column IMAX to column KW-1 of W and update it
294 *
295  CALL zcopy( imax, a( 1, imax ), 1, w( 1, kw-1 ), 1 )
296  CALL zcopy( k-imax, a( imax, imax+1 ), lda,
297  $ w( imax+1, kw-1 ), 1 )
298  IF( k.LT.n )
299  $ CALL zgemv( 'No transpose', k, n-k, -cone,
300  $ a( 1, k+1 ), lda, w( imax, kw+1 ), ldw,
301  $ cone, w( 1, kw-1 ), 1 )
302 *
303 * JMAX is the column-index of the largest off-diagonal
304 * element in row IMAX, and ROWMAX is its absolute value
305 *
306  jmax = imax + izamax( k-imax, w( imax+1, kw-1 ), 1 )
307  rowmax = cabs1( w( jmax, kw-1 ) )
308  IF( imax.GT.1 ) THEN
309  jmax = izamax( imax-1, w( 1, kw-1 ), 1 )
310  rowmax = max( rowmax, cabs1( w( jmax, kw-1 ) ) )
311  END IF
312 *
313  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
314 *
315 * no interchange, use 1-by-1 pivot block
316 *
317  kp = k
318  ELSE IF( cabs1( w( imax, kw-1 ) ).GE.alpha*rowmax ) THEN
319 *
320 * interchange rows and columns K and IMAX, use 1-by-1
321 * pivot block
322 *
323  kp = imax
324 *
325 * copy column KW-1 of W to column KW of W
326 *
327  CALL zcopy( k, w( 1, kw-1 ), 1, w( 1, kw ), 1 )
328  ELSE
329 *
330 * interchange rows and columns K-1 and IMAX, use 2-by-2
331 * pivot block
332 *
333  kp = imax
334  kstep = 2
335  END IF
336  END IF
337 *
338 * ============================================================
339 *
340 * KK is the column of A where pivoting step stopped
341 *
342  kk = k - kstep + 1
343 *
344 * KKW is the column of W which corresponds to column KK of A
345 *
346  kkw = nb + kk - n
347 *
348 * Interchange rows and columns KP and KK.
349 * Updated column KP is already stored in column KKW of W.
350 *
351  IF( kp.NE.kk ) THEN
352 *
353 * Copy non-updated column KK to column KP of submatrix A
354 * at step K. No need to copy element into column K
355 * (or K and K-1 for 2-by-2 pivot) of A, since these columns
356 * will be later overwritten.
357 *
358  a( kp, kp ) = a( kk, kk )
359  CALL zcopy( kk-1-kp, a( kp+1, kk ), 1, a( kp, kp+1 ),
360  $ lda )
361  IF( kp.GT.1 )
362  $ CALL zcopy( kp-1, a( 1, kk ), 1, a( 1, kp ), 1 )
363 *
364 * Interchange rows KK and KP in last K+1 to N columns of A
365 * (columns K (or K and K-1 for 2-by-2 pivot) of A will be
366 * later overwritten). Interchange rows KK and KP
367 * in last KKW to NB columns of W.
368 *
369  IF( k.LT.n )
370  $ CALL zswap( n-k, a( kk, k+1 ), lda, a( kp, k+1 ),
371  $ lda )
372  CALL zswap( n-kk+1, w( kk, kkw ), ldw, w( kp, kkw ),
373  $ ldw )
374  END IF
375 *
376  IF( kstep.EQ.1 ) THEN
377 *
378 * 1-by-1 pivot block D(k): column kw of W now holds
379 *
380 * W(kw) = U(k)*D(k),
381 *
382 * where U(k) is the k-th column of U
383 *
384 * Store subdiag. elements of column U(k)
385 * and 1-by-1 block D(k) in column k of A.
386 * NOTE: Diagonal element U(k,k) is a UNIT element
387 * and not stored.
388 * A(k,k) := D(k,k) = W(k,kw)
389 * A(1:k-1,k) := U(1:k-1,k) = W(1:k-1,kw)/D(k,k)
390 *
391  CALL zcopy( k, w( 1, kw ), 1, a( 1, k ), 1 )
392  r1 = cone / a( k, k )
393  CALL zscal( k-1, r1, a( 1, k ), 1 )
394 *
395  ELSE
396 *
397 * 2-by-2 pivot block D(k): columns kw and kw-1 of W now hold
398 *
399 * ( W(kw-1) W(kw) ) = ( U(k-1) U(k) )*D(k)
400 *
401 * where U(k) and U(k-1) are the k-th and (k-1)-th columns
402 * of U
403 *
404 * Store U(1:k-2,k-1) and U(1:k-2,k) and 2-by-2
405 * block D(k-1:k,k-1:k) in columns k-1 and k of A.
406 * NOTE: 2-by-2 diagonal block U(k-1:k,k-1:k) is a UNIT
407 * block and not stored.
408 * A(k-1:k,k-1:k) := D(k-1:k,k-1:k) = W(k-1:k,kw-1:kw)
409 * A(1:k-2,k-1:k) := U(1:k-2,k:k-1:k) =
410 * = W(1:k-2,kw-1:kw) * ( D(k-1:k,k-1:k)**(-1) )
411 *
412  IF( k.GT.2 ) THEN
413 *
414 * Compose the columns of the inverse of 2-by-2 pivot
415 * block D in the following way to reduce the number
416 * of FLOPS when we myltiply panel ( W(kw-1) W(kw) ) by
417 * this inverse
418 *
419 * D**(-1) = ( d11 d21 )**(-1) =
420 * ( d21 d22 )
421 *
422 * = 1/(d11*d22-d21**2) * ( ( d22 ) (-d21 ) ) =
423 * ( (-d21 ) ( d11 ) )
424 *
425 * = 1/d21 * 1/((d11/d21)*(d22/d21)-1) *
426 *
427 * * ( ( d22/d21 ) ( -1 ) ) =
428 * ( ( -1 ) ( d11/d21 ) )
429 *
430 * = 1/d21 * 1/(D22*D11-1) * ( ( D11 ) ( -1 ) ) =
431 * ( ( -1 ) ( D22 ) )
432 *
433 * = 1/d21 * T * ( ( D11 ) ( -1 ) )
434 * ( ( -1 ) ( D22 ) )
435 *
436 * = D21 * ( ( D11 ) ( -1 ) )
437 * ( ( -1 ) ( D22 ) )
438 *
439  d21 = w( k-1, kw )
440  d11 = w( k, kw ) / d21
441  d22 = w( k-1, kw-1 ) / d21
442  t = cone / ( d11*d22-cone )
443  d21 = t / d21
444 *
445 * Update elements in columns A(k-1) and A(k) as
446 * dot products of rows of ( W(kw-1) W(kw) ) and columns
447 * of D**(-1)
448 *
449  DO 20 j = 1, k - 2
450  a( j, k-1 ) = d21*( d11*w( j, kw-1 )-w( j, kw ) )
451  a( j, k ) = d21*( d22*w( j, kw )-w( j, kw-1 ) )
452  20 CONTINUE
453  END IF
454 *
455 * Copy D(k) to A
456 *
457  a( k-1, k-1 ) = w( k-1, kw-1 )
458  a( k-1, k ) = w( k-1, kw )
459  a( k, k ) = w( k, kw )
460 *
461  END IF
462 *
463  END IF
464 *
465 * Store details of the interchanges in IPIV
466 *
467  IF( kstep.EQ.1 ) THEN
468  ipiv( k ) = kp
469  ELSE
470  ipiv( k ) = -kp
471  ipiv( k-1 ) = -kp
472  END IF
473 *
474 * Decrease K and return to the start of the main loop
475 *
476  k = k - kstep
477  GO TO 10
478 *
479  30 CONTINUE
480 *
481 * Update the upper triangle of A11 (= A(1:k,1:k)) as
482 *
483 * A11 := A11 - U12*D*U12**T = A11 - U12*W**T
484 *
485 * computing blocks of NB columns at a time
486 *
487  DO 50 j = ( ( k-1 ) / nb )*nb + 1, 1, -nb
488  jb = min( nb, k-j+1 )
489 *
490 * Update the upper triangle of the diagonal block
491 *
492  DO 40 jj = j, j + jb - 1
493  CALL zgemv( 'No transpose', jj-j+1, n-k, -cone,
494  $ a( j, k+1 ), lda, w( jj, kw+1 ), ldw, cone,
495  $ a( j, jj ), 1 )
496  40 CONTINUE
497 *
498 * Update the rectangular superdiagonal block
499 *
500  CALL zgemm( 'No transpose', 'Transpose', j-1, jb, n-k,
501  $ -cone, a( 1, k+1 ), lda, w( j, kw+1 ), ldw,
502  $ cone, a( 1, j ), lda )
503  50 CONTINUE
504 *
505 * Put U12 in standard form by partially undoing the interchanges
506 * in columns k+1:n looping backwards from k+1 to n
507 *
508  j = k + 1
509  60 CONTINUE
510 *
511 * Undo the interchanges (if any) of rows JJ and JP at each
512 * step J
513 *
514 * (Here, J is a diagonal index)
515  jj = j
516  jp = ipiv( j )
517  IF( jp.LT.0 ) THEN
518  jp = -jp
519 * (Here, J is a diagonal index)
520  j = j + 1
521  END IF
522 * (NOTE: Here, J is used to determine row length. Length N-J+1
523 * of the rows to swap back doesn't include diagonal element)
524  j = j + 1
525  IF( jp.NE.jj .AND. j.LE.n )
526  $ CALL zswap( n-j+1, a( jp, j ), lda, a( jj, j ), lda )
527  IF( j.LT.n )
528  $ GO TO 60
529 *
530 * Set KB to the number of columns factorized
531 *
532  kb = n - k
533 *
534  ELSE
535 *
536 * Factorize the leading columns of A using the lower triangle
537 * of A and working forwards, and compute the matrix W = L21*D
538 * for use in updating A22
539 *
540 * K is the main loop index, increasing from 1 in steps of 1 or 2
541 *
542  k = 1
543  70 CONTINUE
544 *
545 * Exit from loop
546 *
547  IF( ( k.GE.nb .AND. nb.LT.n ) .OR. k.GT.n )
548  $ GO TO 90
549 *
550 * Copy column K of A to column K of W and update it
551 *
552  CALL zcopy( n-k+1, a( k, k ), 1, w( k, k ), 1 )
553  CALL zgemv( 'No transpose', n-k+1, k-1, -cone, a( k, 1 ), lda,
554  $ w( k, 1 ), ldw, cone, w( k, k ), 1 )
555 *
556  kstep = 1
557 *
558 * Determine rows and columns to be interchanged and whether
559 * a 1-by-1 or 2-by-2 pivot block will be used
560 *
561  absakk = cabs1( w( k, k ) )
562 *
563 * IMAX is the row-index of the largest off-diagonal element in
564 
565 *
566  IF( k.LT.n ) THEN
567  imax = k + izamax( n-k, w( k+1, k ), 1 )
568  colmax = cabs1( w( imax, k ) )
569  ELSE
570  colmax = zero
571  END IF
572 *
573  IF( max( absakk, colmax ).EQ.zero ) THEN
574 *
575 * Column K is zero or underflow: set INFO and continue
576 *
577  IF( info.EQ.0 )
578  $ info = k
579  kp = k
580  ELSE
581  IF( absakk.GE.alpha*colmax ) THEN
582 *
583 * no interchange, use 1-by-1 pivot block
584 *
585  kp = k
586  ELSE
587 *
588 * Copy column IMAX to column K+1 of W and update it
589 *
590  CALL zcopy( imax-k, a( imax, k ), lda, w( k, k+1 ), 1 )
591  CALL zcopy( n-imax+1, a( imax, imax ), 1, w( imax, k+1 ),
592  $ 1 )
593  CALL zgemv( 'No transpose', n-k+1, k-1, -cone, a( k, 1 ),
594  $ lda, w( imax, 1 ), ldw, cone, w( k, k+1 ),
595  $ 1 )
596 *
597 * JMAX is the column-index of the largest off-diagonal
598 * element in row IMAX, and ROWMAX is its absolute value
599 *
600  jmax = k - 1 + izamax( imax-k, w( k, k+1 ), 1 )
601  rowmax = cabs1( w( jmax, k+1 ) )
602  IF( imax.LT.n ) THEN
603  jmax = imax + izamax( n-imax, w( imax+1, k+1 ), 1 )
604  rowmax = max( rowmax, cabs1( w( jmax, k+1 ) ) )
605  END IF
606 *
607  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
608 *
609 * no interchange, use 1-by-1 pivot block
610 *
611  kp = k
612  ELSE IF( cabs1( w( imax, k+1 ) ).GE.alpha*rowmax ) THEN
613 *
614 * interchange rows and columns K and IMAX, use 1-by-1
615 * pivot block
616 *
617  kp = imax
618 *
619 * copy column K+1 of W to column K of W
620 *
621  CALL zcopy( n-k+1, w( k, k+1 ), 1, w( k, k ), 1 )
622  ELSE
623 *
624 * interchange rows and columns K+1 and IMAX, use 2-by-2
625 * pivot block
626 *
627  kp = imax
628  kstep = 2
629  END IF
630  END IF
631 *
632 * ============================================================
633 *
634 * KK is the column of A where pivoting step stopped
635 *
636  kk = k + kstep - 1
637 *
638 * Interchange rows and columns KP and KK.
639 * Updated column KP is already stored in column KK of W.
640 *
641  IF( kp.NE.kk ) THEN
642 *
643 * Copy non-updated column KK to column KP of submatrix A
644 * at step K. No need to copy element into column K
645 * (or K and K+1 for 2-by-2 pivot) of A, since these columns
646 * will be later overwritten.
647 *
648  a( kp, kp ) = a( kk, kk )
649  CALL zcopy( kp-kk-1, a( kk+1, kk ), 1, a( kp, kk+1 ),
650  $ lda )
651  IF( kp.LT.n )
652  $ CALL zcopy( n-kp, a( kp+1, kk ), 1, a( kp+1, kp ), 1 )
653 *
654 * Interchange rows KK and KP in first K-1 columns of A
655 * (columns K (or K and K+1 for 2-by-2 pivot) of A will be
656 * later overwritten). Interchange rows KK and KP
657 * in first KK columns of W.
658 *
659  IF( k.GT.1 )
660  $ CALL zswap( k-1, a( kk, 1 ), lda, a( kp, 1 ), lda )
661  CALL zswap( kk, w( kk, 1 ), ldw, w( kp, 1 ), ldw )
662  END IF
663 *
664  IF( kstep.EQ.1 ) THEN
665 *
666 * 1-by-1 pivot block D(k): column k of W now holds
667 *
668 * W(k) = L(k)*D(k),
669 *
670 * where L(k) is the k-th column of L
671 *
672 * Store subdiag. elements of column L(k)
673 * and 1-by-1 block D(k) in column k of A.
674 * (NOTE: Diagonal element L(k,k) is a UNIT element
675 * and not stored)
676 * A(k,k) := D(k,k) = W(k,k)
677 * A(k+1:N,k) := L(k+1:N,k) = W(k+1:N,k)/D(k,k)
678 *
679  CALL zcopy( n-k+1, w( k, k ), 1, a( k, k ), 1 )
680  IF( k.LT.n ) THEN
681  r1 = cone / a( k, k )
682  CALL zscal( n-k, r1, a( k+1, k ), 1 )
683  END IF
684 *
685  ELSE
686 *
687 * 2-by-2 pivot block D(k): columns k and k+1 of W now hold
688 *
689 * ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k)
690 *
691 * where L(k) and L(k+1) are the k-th and (k+1)-th columns
692 * of L
693 *
694 * Store L(k+2:N,k) and L(k+2:N,k+1) and 2-by-2
695 * block D(k:k+1,k:k+1) in columns k and k+1 of A.
696 * (NOTE: 2-by-2 diagonal block L(k:k+1,k:k+1) is a UNIT
697 * block and not stored)
698 * A(k:k+1,k:k+1) := D(k:k+1,k:k+1) = W(k:k+1,k:k+1)
699 * A(k+2:N,k:k+1) := L(k+2:N,k:k+1) =
700 * = W(k+2:N,k:k+1) * ( D(k:k+1,k:k+1)**(-1) )
701 *
702  IF( k.LT.n-1 ) THEN
703 *
704 * Compose the columns of the inverse of 2-by-2 pivot
705 * block D in the following way to reduce the number
706 * of FLOPS when we myltiply panel ( W(k) W(k+1) ) by
707 * this inverse
708 *
709 * D**(-1) = ( d11 d21 )**(-1) =
710 * ( d21 d22 )
711 *
712 * = 1/(d11*d22-d21**2) * ( ( d22 ) (-d21 ) ) =
713 * ( (-d21 ) ( d11 ) )
714 *
715 * = 1/d21 * 1/((d11/d21)*(d22/d21)-1) *
716 *
717 * * ( ( d22/d21 ) ( -1 ) ) =
718 * ( ( -1 ) ( d11/d21 ) )
719 *
720 * = 1/d21 * 1/(D22*D11-1) * ( ( D11 ) ( -1 ) ) =
721 * ( ( -1 ) ( D22 ) )
722 *
723 * = 1/d21 * T * ( ( D11 ) ( -1 ) )
724 * ( ( -1 ) ( D22 ) )
725 *
726 * = D21 * ( ( D11 ) ( -1 ) )
727 * ( ( -1 ) ( D22 ) )
728 *
729  d21 = w( k+1, k )
730  d11 = w( k+1, k+1 ) / d21
731  d22 = w( k, k ) / d21
732  t = cone / ( d11*d22-cone )
733  d21 = t / d21
734 *
735 * Update elements in columns A(k) and A(k+1) as
736 * dot products of rows of ( W(k) W(k+1) ) and columns
737 * of D**(-1)
738 *
739  DO 80 j = k + 2, n
740  a( j, k ) = d21*( d11*w( j, k )-w( j, k+1 ) )
741  a( j, k+1 ) = d21*( d22*w( j, k+1 )-w( j, k ) )
742  80 CONTINUE
743  END IF
744 *
745 * Copy D(k) to A
746 *
747  a( k, k ) = w( k, k )
748  a( k+1, k ) = w( k+1, k )
749  a( k+1, k+1 ) = w( k+1, k+1 )
750 *
751  END IF
752 *
753  END IF
754 *
755 * Store details of the interchanges in IPIV
756 *
757  IF( kstep.EQ.1 ) THEN
758  ipiv( k ) = kp
759  ELSE
760  ipiv( k ) = -kp
761  ipiv( k+1 ) = -kp
762  END IF
763 *
764 * Increase K and return to the start of the main loop
765 *
766  k = k + kstep
767  GO TO 70
768 *
769  90 CONTINUE
770 *
771 * Update the lower triangle of A22 (= A(k:n,k:n)) as
772 *
773 * A22 := A22 - L21*D*L21**T = A22 - L21*W**T
774 *
775 * computing blocks of NB columns at a time
776 *
777  DO 110 j = k, n, nb
778  jb = min( nb, n-j+1 )
779 *
780 * Update the lower triangle of the diagonal block
781 *
782  DO 100 jj = j, j + jb - 1
783  CALL zgemv( 'No transpose', j+jb-jj, k-1, -cone,
784  $ a( jj, 1 ), lda, w( jj, 1 ), ldw, cone,
785  $ a( jj, jj ), 1 )
786  100 CONTINUE
787 *
788 * Update the rectangular subdiagonal block
789 *
790  IF( j+jb.LE.n )
791  $ CALL zgemm( 'No transpose', 'Transpose', n-j-jb+1, jb,
792  $ k-1, -cone, a( j+jb, 1 ), lda, w( j, 1 ),
793  $ ldw, cone, a( j+jb, j ), lda )
794  110 CONTINUE
795 *
796 * Put L21 in standard form by partially undoing the interchanges
797 * of rows in columns 1:k-1 looping backwards from k-1 to 1
798 *
799  j = k - 1
800  120 CONTINUE
801 *
802 * Undo the interchanges (if any) of rows JJ and JP at each
803 * step J
804 *
805 * (Here, J is a diagonal index)
806  jj = j
807  jp = ipiv( j )
808  IF( jp.LT.0 ) THEN
809  jp = -jp
810 * (Here, J is a diagonal index)
811  j = j - 1
812  END IF
813 * (NOTE: Here, J is used to determine row length. Length J
814 * of the rows to swap back doesn't include diagonal element)
815  j = j - 1
816  IF( jp.NE.jj .AND. j.GE.1 )
817  $ CALL zswap( j, a( jp, 1 ), lda, a( jj, 1 ), lda )
818  IF( j.GT.1 )
819  $ GO TO 120
820 *
821 * Set KB to the number of columns factorized
822 *
823  kb = k - 1
824 *
825  END IF
826  RETURN
827 *
828 * End of ZLASYF
829 *
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
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zscal(N, ZA, ZX, INCX)
ZSCAL
Definition: zscal.f:54

Here is the call graph for this function:

Here is the caller graph for this function: