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

◆ dlasyf_rk()

subroutine dlasyf_rk ( character uplo,
integer n,
integer nb,
integer kb,
double precision, dimension( lda, * ) a,
integer lda,
double precision, dimension( * ) e,
integer, dimension( * ) ipiv,
double precision, dimension( ldw, * ) w,
integer ldw,
integer info )

DLASYF_RK computes a partial factorization of a real symmetric indefinite matrix using bounded Bunch-Kaufman (rook) diagonal pivoting method.

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

Purpose:
!> DLASYF_RK computes a partial factorization of a real symmetric
!> matrix A using the bounded Bunch-Kaufman (rook) 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.
!>
!> DLASYF_RK is an auxiliary routine called by DSYTRF_RK. 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 DOUBLE PRECISION 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, contains:
!>            a) ONLY diagonal elements of the symmetric block diagonal
!>               matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
!>               (superdiagonal (or subdiagonal) elements of D
!>                are stored on exit in array E), and
!>            b) If UPLO = 'U': factor U in the superdiagonal part of A.
!>               If UPLO = 'L': factor L in the subdiagonal part of A.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]E
!>          E is DOUBLE PRECISION array, dimension (N)
!>          On exit, contains the superdiagonal (or subdiagonal)
!>          elements of the symmetric block diagonal matrix D
!>          with 1-by-1 or 2-by-2 diagonal blocks, where
!>          If UPLO = 'U': E(i) = D(i-1,i), i=2:N, E(1) is set to 0;
!>          If UPLO = 'L': E(i) = D(i+1,i), i=1:N-1, E(N) is set to 0.
!>
!>          NOTE: For 1-by-1 diagonal block D(k), where
!>          1 <= k <= N, the element E(k) is set to 0 in both
!>          UPLO = 'U' or UPLO = 'L' cases.
!> 
[out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          IPIV describes the permutation matrix P in the factorization
!>          of matrix A as follows. The absolute value of IPIV(k)
!>          represents the index of row and column that were
!>          interchanged with the k-th row and column. The value of UPLO
!>          describes the order in which the interchanges were applied.
!>          Also, the sign of IPIV represents the block structure of
!>          the symmetric block diagonal matrix D with 1-by-1 or 2-by-2
!>          diagonal blocks which correspond to 1 or 2 interchanges
!>          at each factorization step.
!>
!>          If UPLO = 'U',
!>          ( in factorization order, k decreases from N to 1 ):
!>            a) A single positive entry IPIV(k) > 0 means:
!>               D(k,k) is a 1-by-1 diagonal block.
!>               If IPIV(k) != k, rows and columns k and IPIV(k) were
!>               interchanged in the submatrix A(1:N,N-KB+1:N);
!>               If IPIV(k) = k, no interchange occurred.
!>
!>
!>            b) A pair of consecutive negative entries
!>               IPIV(k) < 0 and IPIV(k-1) < 0 means:
!>               D(k-1:k,k-1:k) is a 2-by-2 diagonal block.
!>               (NOTE: negative entries in IPIV appear ONLY in pairs).
!>               1) If -IPIV(k) != k, rows and columns
!>                  k and -IPIV(k) were interchanged
!>                  in the matrix A(1:N,N-KB+1:N).
!>                  If -IPIV(k) = k, no interchange occurred.
!>               2) If -IPIV(k-1) != k-1, rows and columns
!>                  k-1 and -IPIV(k-1) were interchanged
!>                  in the submatrix A(1:N,N-KB+1:N).
!>                  If -IPIV(k-1) = k-1, no interchange occurred.
!>
!>            c) In both cases a) and b) is always ABS( IPIV(k) ) <= k.
!>
!>            d) NOTE: Any entry IPIV(k) is always NONZERO on output.
!>
!>          If UPLO = 'L',
!>          ( in factorization order, k increases from 1 to N ):
!>            a) A single positive entry IPIV(k) > 0 means:
!>               D(k,k) is a 1-by-1 diagonal block.
!>               If IPIV(k) != k, rows and columns k and IPIV(k) were
!>               interchanged in the submatrix A(1:N,1:KB).
!>               If IPIV(k) = k, no interchange occurred.
!>
!>            b) A pair of consecutive negative entries
!>               IPIV(k) < 0 and IPIV(k+1) < 0 means:
!>               D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
!>               (NOTE: negative entries in IPIV appear ONLY in pairs).
!>               1) If -IPIV(k) != k, rows and columns
!>                  k and -IPIV(k) were interchanged
!>                  in the submatrix A(1:N,1:KB).
!>                  If -IPIV(k) = k, no interchange occurred.
!>               2) If -IPIV(k+1) != k+1, rows and columns
!>                  k-1 and -IPIV(k-1) were interchanged
!>                  in the submatrix A(1:N,1:KB).
!>                  If -IPIV(k+1) = k+1, no interchange occurred.
!>
!>            c) In both cases a) and b) is always ABS( IPIV(k) ) >= k.
!>
!>            d) NOTE: Any entry IPIV(k) is always NONZERO on output.
!> 
[out]W
!>          W is DOUBLE PRECISION 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, the k-th argument had an illegal value
!>
!>          > 0: If INFO = k, the matrix A is singular, because:
!>                 If UPLO = 'U': column k in the upper
!>                 triangular part of A contains all zeros.
!>                 If UPLO = 'L': column k in the lower
!>                 triangular part of A contains all zeros.
!>
!>               Therefore D(k,k) is exactly zero, and superdiagonal
!>               elements of column k of U (or subdiagonal elements of
!>               column k of L ) are all zeros. The factorization has
!>               been completed, but the block diagonal matrix D is
!>               exactly singular, and division by zero will occur if
!>               it is used to solve a system of equations.
!>
!>               NOTE: INFO only stores the first occurrence of
!>               a singularity, any subsequent occurrence of singularity
!>               is not stored in INFO even though the factorization
!>               always completes.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
!>
!>  December 2016,  Igor Kozachenko,
!>                  Computer Science Division,
!>                  University of California, Berkeley
!>
!>  September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
!>                  School of Mathematics,
!>                  University of Manchester
!>
!> 

Definition at line 258 of file dlasyf_rk.f.

260*
261* -- LAPACK computational routine --
262* -- LAPACK is a software package provided by Univ. of Tennessee, --
263* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
264*
265* .. Scalar Arguments ..
266 CHARACTER UPLO
267 INTEGER INFO, KB, LDA, LDW, N, NB
268* ..
269* .. Array Arguments ..
270 INTEGER IPIV( * )
271 DOUBLE PRECISION A( LDA, * ), E( * ), W( LDW, * )
272* ..
273*
274* =====================================================================
275*
276* .. Parameters ..
277 DOUBLE PRECISION ZERO, ONE
278 parameter( zero = 0.0d+0, one = 1.0d+0 )
279 DOUBLE PRECISION EIGHT, SEVTEN
280 parameter( eight = 8.0d+0, sevten = 17.0d+0 )
281* ..
282* .. Local Scalars ..
283 LOGICAL DONE
284 INTEGER IMAX, ITEMP, J, JMAX, K, KK, KW, KKW,
285 $ KP, KSTEP, P, II
286 DOUBLE PRECISION ABSAKK, ALPHA, COLMAX, D11, D12, D21, D22,
287 $ DTEMP, R1, ROWMAX, T, SFMIN
288* ..
289* .. External Functions ..
290 LOGICAL LSAME
291 INTEGER IDAMAX
292 DOUBLE PRECISION DLAMCH
293 EXTERNAL lsame, idamax, dlamch
294* ..
295* .. External Subroutines ..
296 EXTERNAL dcopy, dgemmtr, dgemv, dscal, dswap
297* ..
298* .. Intrinsic Functions ..
299 INTRINSIC abs, max, min, sqrt
300* ..
301* .. Executable Statements ..
302*
303 info = 0
304*
305* Initialize ALPHA for use in choosing pivot block size.
306*
307 alpha = ( one+sqrt( sevten ) ) / eight
308*
309* Compute machine safe minimum
310*
311 sfmin = dlamch( 'S' )
312*
313 IF( lsame( uplo, 'U' ) ) THEN
314*
315* Factorize the trailing columns of A using the upper triangle
316* of A and working backwards, and compute the matrix W = U12*D
317* for use in updating A11
318*
319* Initialize the first entry of array E, where superdiagonal
320* elements of D are stored
321*
322 e( 1 ) = zero
323*
324* K is the main loop index, decreasing from N in steps of 1 or 2
325*
326 k = n
327 10 CONTINUE
328*
329* KW is the column of W which corresponds to column K of A
330*
331 kw = nb + k - n
332*
333* Exit from loop
334*
335 IF( ( k.LE.n-nb+1 .AND. nb.LT.n ) .OR. k.LT.1 )
336 $ GO TO 30
337*
338 kstep = 1
339 p = k
340*
341* Copy column K of A to column KW of W and update it
342*
343 CALL dcopy( k, a( 1, k ), 1, w( 1, kw ), 1 )
344 IF( k.LT.n )
345 $ CALL dgemv( 'No transpose', k, n-k, -one, a( 1, k+1 ),
346 $ lda, w( k, kw+1 ), ldw, one, w( 1, kw ), 1 )
347*
348* Determine rows and columns to be interchanged and whether
349* a 1-by-1 or 2-by-2 pivot block will be used
350*
351 absakk = abs( w( k, kw ) )
352*
353* IMAX is the row-index of the largest off-diagonal element in
354* column K, and COLMAX is its absolute value.
355* Determine both COLMAX and IMAX.
356*
357 IF( k.GT.1 ) THEN
358 imax = idamax( k-1, w( 1, kw ), 1 )
359 colmax = abs( w( imax, kw ) )
360 ELSE
361 colmax = zero
362 END IF
363*
364 IF( max( absakk, colmax ).EQ.zero ) THEN
365*
366* Column K is zero or underflow: set INFO and continue
367*
368 IF( info.EQ.0 )
369 $ info = k
370 kp = k
371 CALL dcopy( k, w( 1, kw ), 1, a( 1, k ), 1 )
372*
373* Set E( K ) to zero
374*
375 IF( k.GT.1 )
376 $ e( k ) = zero
377*
378 ELSE
379*
380* ============================================================
381*
382* Test for interchange
383*
384* Equivalent to testing for ABSAKK.GE.ALPHA*COLMAX
385* (used to handle NaN and Inf)
386*
387 IF( .NOT.( absakk.LT.alpha*colmax ) ) THEN
388*
389* no interchange, use 1-by-1 pivot block
390*
391 kp = k
392*
393 ELSE
394*
395 done = .false.
396*
397* Loop until pivot found
398*
399 12 CONTINUE
400*
401* Begin pivot search loop body
402*
403*
404* Copy column IMAX to column KW-1 of W and update it
405*
406 CALL dcopy( imax, a( 1, imax ), 1, w( 1, kw-1 ),
407 $ 1 )
408 CALL dcopy( k-imax, a( imax, imax+1 ), lda,
409 $ w( imax+1, kw-1 ), 1 )
410*
411 IF( k.LT.n )
412 $ CALL dgemv( 'No transpose', k, n-k, -one,
413 $ a( 1, k+1 ), lda, w( imax, kw+1 ), ldw,
414 $ one, w( 1, kw-1 ), 1 )
415*
416* JMAX is the column-index of the largest off-diagonal
417* element in row IMAX, and ROWMAX is its absolute value.
418* Determine both ROWMAX and JMAX.
419*
420 IF( imax.NE.k ) THEN
421 jmax = imax + idamax( k-imax, w( imax+1, kw-1 ),
422 $ 1 )
423 rowmax = abs( w( jmax, kw-1 ) )
424 ELSE
425 rowmax = zero
426 END IF
427*
428 IF( imax.GT.1 ) THEN
429 itemp = idamax( imax-1, w( 1, kw-1 ), 1 )
430 dtemp = abs( w( itemp, kw-1 ) )
431 IF( dtemp.GT.rowmax ) THEN
432 rowmax = dtemp
433 jmax = itemp
434 END IF
435 END IF
436*
437* Equivalent to testing for
438* ABS( W( IMAX, KW-1 ) ).GE.ALPHA*ROWMAX
439* (used to handle NaN and Inf)
440*
441 IF( .NOT.(abs( w( imax, kw-1 ) ).LT.alpha*rowmax ) )
442 $ THEN
443*
444* interchange rows and columns K and IMAX,
445* use 1-by-1 pivot block
446*
447 kp = imax
448*
449* copy column KW-1 of W to column KW of W
450*
451 CALL dcopy( k, w( 1, kw-1 ), 1, w( 1, kw ), 1 )
452*
453 done = .true.
454*
455* Equivalent to testing for ROWMAX.EQ.COLMAX,
456* (used to handle NaN and Inf)
457*
458 ELSE IF( ( p.EQ.jmax ) .OR. ( rowmax.LE.colmax ) )
459 $ THEN
460*
461* interchange rows and columns K-1 and IMAX,
462* use 2-by-2 pivot block
463*
464 kp = imax
465 kstep = 2
466 done = .true.
467 ELSE
468*
469* Pivot not found: set params and repeat
470*
471 p = imax
472 colmax = rowmax
473 imax = jmax
474*
475* Copy updated JMAXth (next IMAXth) column to Kth of W
476*
477 CALL dcopy( k, w( 1, kw-1 ), 1, w( 1, kw ), 1 )
478*
479 END IF
480*
481* End pivot search loop body
482*
483 IF( .NOT. done ) GOTO 12
484*
485 END IF
486*
487* ============================================================
488*
489 kk = k - kstep + 1
490*
491* KKW is the column of W which corresponds to column KK of A
492*
493 kkw = nb + kk - n
494*
495 IF( ( kstep.EQ.2 ) .AND. ( p.NE.k ) ) THEN
496*
497* Copy non-updated column K to column P
498*
499 CALL dcopy( k-p, a( p+1, k ), 1, a( p, p+1 ), lda )
500 CALL dcopy( p, a( 1, k ), 1, a( 1, p ), 1 )
501*
502* Interchange rows K and P in last N-K+1 columns of A
503* and last N-K+2 columns of W
504*
505 CALL dswap( n-k+1, a( k, k ), lda, a( p, k ), lda )
506 CALL dswap( n-kk+1, w( k, kkw ), ldw, w( p, kkw ),
507 $ ldw )
508 END IF
509*
510* Updated column KP is already stored in column KKW of W
511*
512 IF( kp.NE.kk ) THEN
513*
514* Copy non-updated column KK to column KP
515*
516 a( kp, k ) = a( kk, k )
517 CALL dcopy( k-1-kp, a( kp+1, kk ), 1, a( kp, kp+1 ),
518 $ lda )
519 CALL dcopy( kp, a( 1, kk ), 1, a( 1, kp ), 1 )
520*
521* Interchange rows KK and KP in last N-KK+1 columns
522* of A and W
523*
524 CALL dswap( n-kk+1, a( kk, kk ), lda, a( kp, kk ),
525 $ lda )
526 CALL dswap( n-kk+1, w( kk, kkw ), ldw, w( kp, kkw ),
527 $ ldw )
528 END IF
529*
530 IF( kstep.EQ.1 ) THEN
531*
532* 1-by-1 pivot block D(k): column KW of W now holds
533*
534* W(k) = U(k)*D(k)
535*
536* where U(k) is the k-th column of U
537*
538* Store U(k) in column k of A
539*
540 CALL dcopy( k, w( 1, kw ), 1, a( 1, k ), 1 )
541 IF( k.GT.1 ) THEN
542 IF( abs( a( k, k ) ).GE.sfmin ) THEN
543 r1 = one / a( k, k )
544 CALL dscal( k-1, r1, a( 1, k ), 1 )
545 ELSE IF( a( k, k ).NE.zero ) THEN
546 DO 14 ii = 1, k - 1
547 a( ii, k ) = a( ii, k ) / a( k, k )
548 14 CONTINUE
549 END IF
550*
551* Store the superdiagonal element of D in array E
552*
553 e( k ) = zero
554*
555 END IF
556*
557 ELSE
558*
559* 2-by-2 pivot block D(k): columns KW and KW-1 of W now
560* hold
561*
562* ( W(k-1) W(k) ) = ( U(k-1) U(k) )*D(k)
563*
564* where U(k) and U(k-1) are the k-th and (k-1)-th columns
565* of U
566*
567 IF( k.GT.2 ) THEN
568*
569* Store U(k) and U(k-1) in columns k and k-1 of A
570*
571 d12 = w( k-1, kw )
572 d11 = w( k, kw ) / d12
573 d22 = w( k-1, kw-1 ) / d12
574 t = one / ( d11*d22-one )
575 DO 20 j = 1, k - 2
576 a( j, k-1 ) = t*( (d11*w( j, kw-1 )-w( j, kw ) ) /
577 $ d12 )
578 a( j, k ) = t*( ( d22*w( j, kw )-w( j, kw-1 ) ) /
579 $ d12 )
580 20 CONTINUE
581 END IF
582*
583* Copy diagonal elements of D(K) to A,
584* copy superdiagonal element of D(K) to E(K) and
585* ZERO out superdiagonal entry of A
586*
587 a( k-1, k-1 ) = w( k-1, kw-1 )
588 a( k-1, k ) = zero
589 a( k, k ) = w( k, kw )
590 e( k ) = w( k-1, kw )
591 e( k-1 ) = zero
592*
593 END IF
594*
595* End column K is nonsingular
596*
597 END IF
598*
599* Store details of the interchanges in IPIV
600*
601 IF( kstep.EQ.1 ) THEN
602 ipiv( k ) = kp
603 ELSE
604 ipiv( k ) = -p
605 ipiv( k-1 ) = -kp
606 END IF
607*
608* Decrease K and return to the start of the main loop
609*
610 k = k - kstep
611 GO TO 10
612*
613 30 CONTINUE
614*
615* Update the upper triangle of A11 (= A(1:k,1:k)) as
616*
617* A11 := A11 - U12*D*U12**T = A11 - U12*W**T
618*
619 CALL dgemmtr( 'Upper', 'No transpose', 'Transpose', k, n-k,
620 $ -one, a( 1, k+1 ), lda, w( 1, kw+1 ), ldw,
621 $ one, a( 1, 1 ), lda )
622*
623* Set KB to the number of columns factorized
624*
625 kb = n - k
626*
627 ELSE
628*
629* Factorize the leading columns of A using the lower triangle
630* of A and working forwards, and compute the matrix W = L21*D
631* for use in updating A22
632*
633* Initialize the unused last entry of the subdiagonal array E.
634*
635 e( n ) = zero
636*
637* K is the main loop index, increasing from 1 in steps of 1 or 2
638*
639 k = 1
640 70 CONTINUE
641*
642* Exit from loop
643*
644 IF( ( k.GE.nb .AND. nb.LT.n ) .OR. k.GT.n )
645 $ GO TO 90
646*
647 kstep = 1
648 p = k
649*
650* Copy column K of A to column K of W and update it
651*
652 CALL dcopy( n-k+1, a( k, k ), 1, w( k, k ), 1 )
653 IF( k.GT.1 )
654 $ CALL dgemv( 'No transpose', n-k+1, k-1, -one, a( k, 1 ),
655 $ lda, w( k, 1 ), ldw, one, w( k, k ), 1 )
656*
657* Determine rows and columns to be interchanged and whether
658* a 1-by-1 or 2-by-2 pivot block will be used
659*
660 absakk = abs( w( k, k ) )
661*
662* IMAX is the row-index of the largest off-diagonal element in
663* column K, and COLMAX is its absolute value.
664* Determine both COLMAX and IMAX.
665*
666 IF( k.LT.n ) THEN
667 imax = k + idamax( n-k, w( k+1, k ), 1 )
668 colmax = abs( w( imax, k ) )
669 ELSE
670 colmax = zero
671 END IF
672*
673 IF( max( absakk, colmax ).EQ.zero ) THEN
674*
675* Column K is zero or underflow: set INFO and continue
676*
677 IF( info.EQ.0 )
678 $ info = k
679 kp = k
680 CALL dcopy( n-k+1, w( k, k ), 1, a( k, k ), 1 )
681*
682* Set E( K ) to zero
683*
684 IF( k.LT.n )
685 $ e( k ) = zero
686*
687 ELSE
688*
689* ============================================================
690*
691* Test for interchange
692*
693* Equivalent to testing for ABSAKK.GE.ALPHA*COLMAX
694* (used to handle NaN and Inf)
695*
696 IF( .NOT.( absakk.LT.alpha*colmax ) ) THEN
697*
698* no interchange, use 1-by-1 pivot block
699*
700 kp = k
701*
702 ELSE
703*
704 done = .false.
705*
706* Loop until pivot found
707*
708 72 CONTINUE
709*
710* Begin pivot search loop body
711*
712*
713* Copy column IMAX to column K+1 of W and update it
714*
715 CALL dcopy( imax-k, a( imax, k ), lda, w( k, k+1 ),
716 $ 1)
717 CALL dcopy( n-imax+1, a( imax, imax ), 1,
718 $ w( imax, k+1 ), 1 )
719 IF( k.GT.1 )
720 $ CALL dgemv( 'No transpose', n-k+1, k-1, -one,
721 $ a( k, 1 ), lda, w( imax, 1 ), ldw,
722 $ one, w( k, k+1 ), 1 )
723*
724* JMAX is the column-index of the largest off-diagonal
725* element in row IMAX, and ROWMAX is its absolute value.
726* Determine both ROWMAX and JMAX.
727*
728 IF( imax.NE.k ) THEN
729 jmax = k - 1 + idamax( imax-k, w( k, k+1 ), 1 )
730 rowmax = abs( w( jmax, k+1 ) )
731 ELSE
732 rowmax = zero
733 END IF
734*
735 IF( imax.LT.n ) THEN
736 itemp = imax + idamax( n-imax, w( imax+1, k+1 ),
737 $ 1)
738 dtemp = abs( w( itemp, k+1 ) )
739 IF( dtemp.GT.rowmax ) THEN
740 rowmax = dtemp
741 jmax = itemp
742 END IF
743 END IF
744*
745* Equivalent to testing for
746* ABS( W( IMAX, K+1 ) ).GE.ALPHA*ROWMAX
747* (used to handle NaN and Inf)
748*
749 IF( .NOT.( abs( w( imax, k+1 ) ).LT.alpha*rowmax ) )
750 $ THEN
751*
752* interchange rows and columns K and IMAX,
753* use 1-by-1 pivot block
754*
755 kp = imax
756*
757* copy column K+1 of W to column K of W
758*
759 CALL dcopy( n-k+1, w( k, k+1 ), 1, w( k, k ),
760 $ 1 )
761*
762 done = .true.
763*
764* Equivalent to testing for ROWMAX.EQ.COLMAX,
765* (used to handle NaN and Inf)
766*
767 ELSE IF( ( p.EQ.jmax ) .OR. ( rowmax.LE.colmax ) )
768 $ THEN
769*
770* interchange rows and columns K+1 and IMAX,
771* use 2-by-2 pivot block
772*
773 kp = imax
774 kstep = 2
775 done = .true.
776 ELSE
777*
778* Pivot not found: set params and repeat
779*
780 p = imax
781 colmax = rowmax
782 imax = jmax
783*
784* Copy updated JMAXth (next IMAXth) column to Kth of W
785*
786 CALL dcopy( n-k+1, w( k, k+1 ), 1, w( k, k ),
787 $ 1 )
788*
789 END IF
790*
791* End pivot search loop body
792*
793 IF( .NOT. done ) GOTO 72
794*
795 END IF
796*
797* ============================================================
798*
799 kk = k + kstep - 1
800*
801 IF( ( kstep.EQ.2 ) .AND. ( p.NE.k ) ) THEN
802*
803* Copy non-updated column K to column P
804*
805 CALL dcopy( p-k, a( k, k ), 1, a( p, k ), lda )
806 CALL dcopy( n-p+1, a( p, k ), 1, a( p, p ), 1 )
807*
808* Interchange rows K and P in first K columns of A
809* and first K+1 columns of W
810*
811 CALL dswap( k, a( k, 1 ), lda, a( p, 1 ), lda )
812 CALL dswap( kk, w( k, 1 ), ldw, w( p, 1 ), ldw )
813 END IF
814*
815* Updated column KP is already stored in column KK of W
816*
817 IF( kp.NE.kk ) THEN
818*
819* Copy non-updated column KK to column KP
820*
821 a( kp, k ) = a( kk, k )
822 CALL dcopy( kp-k-1, a( k+1, kk ), 1, a( kp, k+1 ),
823 $ lda )
824 CALL dcopy( n-kp+1, a( kp, kk ), 1, a( kp, kp ), 1 )
825*
826* Interchange rows KK and KP in first KK columns of A and W
827*
828 CALL dswap( kk, a( kk, 1 ), lda, a( kp, 1 ), lda )
829 CALL dswap( kk, w( kk, 1 ), ldw, w( kp, 1 ), ldw )
830 END IF
831*
832 IF( kstep.EQ.1 ) THEN
833*
834* 1-by-1 pivot block D(k): column k of W now holds
835*
836* W(k) = L(k)*D(k)
837*
838* where L(k) is the k-th column of L
839*
840* Store L(k) in column k of A
841*
842 CALL dcopy( n-k+1, w( k, k ), 1, a( k, k ), 1 )
843 IF( k.LT.n ) THEN
844 IF( abs( a( k, k ) ).GE.sfmin ) THEN
845 r1 = one / a( k, k )
846 CALL dscal( n-k, r1, a( k+1, k ), 1 )
847 ELSE IF( a( k, k ).NE.zero ) THEN
848 DO 74 ii = k + 1, n
849 a( ii, k ) = a( ii, k ) / a( k, k )
850 74 CONTINUE
851 END IF
852*
853* Store the subdiagonal element of D in array E
854*
855 e( k ) = zero
856*
857 END IF
858*
859 ELSE
860*
861* 2-by-2 pivot block D(k): columns k and k+1 of W now hold
862*
863* ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k)
864*
865* where L(k) and L(k+1) are the k-th and (k+1)-th columns
866* of L
867*
868 IF( k.LT.n-1 ) THEN
869*
870* Store L(k) and L(k+1) in columns k and k+1 of A
871*
872 d21 = w( k+1, k )
873 d11 = w( k+1, k+1 ) / d21
874 d22 = w( k, k ) / d21
875 t = one / ( d11*d22-one )
876 DO 80 j = k + 2, n
877 a( j, k ) = t*( ( d11*w( j, k )-w( j, k+1 ) ) /
878 $ d21 )
879 a( j, k+1 ) = t*( ( d22*w( j, k+1 )-w( j, k ) ) /
880 $ d21 )
881 80 CONTINUE
882 END IF
883*
884* Copy diagonal elements of D(K) to A,
885* copy subdiagonal element of D(K) to E(K) and
886* ZERO out subdiagonal entry of A
887*
888 a( k, k ) = w( k, k )
889 a( k+1, k ) = zero
890 a( k+1, k+1 ) = w( k+1, k+1 )
891 e( k ) = w( k+1, k )
892 e( k+1 ) = zero
893*
894 END IF
895*
896* End column K is nonsingular
897*
898 END IF
899*
900* Store details of the interchanges in IPIV
901*
902 IF( kstep.EQ.1 ) THEN
903 ipiv( k ) = kp
904 ELSE
905 ipiv( k ) = -p
906 ipiv( k+1 ) = -kp
907 END IF
908*
909* Increase K and return to the start of the main loop
910*
911 k = k + kstep
912 GO TO 70
913*
914 90 CONTINUE
915*
916* Update the lower triangle of A22 (= A(k:n,k:n)) as
917*
918* A22 := A22 - L21*D*L21**T = A22 - L21*W**T
919*
920 CALL dgemmtr( 'Lower', 'No transpose', 'Transpose', n-k+1,
921 $ k-1, -one, a( k, 1 ), lda, w( k, 1 ), ldw,
922 $ one, a( k, k ), lda )
923*
924* Set KB to the number of columns factorized
925*
926 kb = k - 1
927*
928 END IF
929*
930 RETURN
931*
932* End of DLASYF_RK
933*
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dgemmtr(uplo, transa, transb, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMMTR
Definition dgemmtr.f:191
subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
DGEMV
Definition dgemv.f:158
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
Here is the call graph for this function:
Here is the caller graph for this function: