File:  [local] / rpl / lapack / lapack / dlasyf_rk.f
Revision 1.5: download - view: text, annotated - select for diffs - revision graph
Mon Aug 7 08:39:00 2023 UTC (9 months ago) by bertrand
Branches: MAIN
CVS tags: rpl-4_1_35, rpl-4_1_34, HEAD
Première mise à jour de lapack et blas.

    1: *> \brief \b DLASYF_RK computes a partial factorization of a real symmetric indefinite matrix using bounded Bunch-Kaufman (rook) diagonal pivoting method.
    2: *
    3: *  =========== DOCUMENTATION ===========
    4: *
    5: * Online html documentation available at
    6: *            http://www.netlib.org/lapack/explore-html/
    7: *
    8: *> \htmlonly
    9: *> Download DLASYF_RK + dependencies
   10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlasyf_rk.f">
   11: *> [TGZ]</a>
   12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlasyf_rk.f">
   13: *> [ZIP]</a>
   14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlasyf_rk.f">
   15: *> [TXT]</a>
   16: *> \endhtmlonly
   17: *
   18: *  Definition:
   19: *  ===========
   20: *
   21: *       SUBROUTINE DLASYF_RK( UPLO, N, NB, KB, A, LDA, E, IPIV, W, LDW,
   22: *                             INFO )
   23: *
   24: *       .. Scalar Arguments ..
   25: *       CHARACTER          UPLO
   26: *       INTEGER            INFO, KB, LDA, LDW, N, NB
   27: *       ..
   28: *       .. Array Arguments ..
   29: *       INTEGER            IPIV( * )
   30: *       DOUBLE PRECISION   A( LDA, * ), E( * ), W( LDW, * )
   31: *       ..
   32: *
   33: *
   34: *> \par Purpose:
   35: *  =============
   36: *>
   37: *> \verbatim
   38: *> DLASYF_RK computes a partial factorization of a real symmetric
   39: *> matrix A using the bounded Bunch-Kaufman (rook) diagonal
   40: *> pivoting method. The partial factorization has the form:
   41: *>
   42: *> A  =  ( I  U12 ) ( A11  0  ) (  I       0    )  if UPLO = 'U', or:
   43: *>       ( 0  U22 ) (  0   D  ) ( U12**T U22**T )
   44: *>
   45: *> A  =  ( L11  0 ) (  D   0  ) ( L11**T L21**T )  if UPLO = 'L',
   46: *>       ( L21  I ) (  0  A22 ) (  0       I    )
   47: *>
   48: *> where the order of D is at most NB. The actual order is returned in
   49: *> the argument KB, and is either NB or NB-1, or N if N <= NB.
   50: *>
   51: *> DLASYF_RK is an auxiliary routine called by DSYTRF_RK. It uses
   52: *> blocked code (calling Level 3 BLAS) to update the submatrix
   53: *> A11 (if UPLO = 'U') or A22 (if UPLO = 'L').
   54: *> \endverbatim
   55: *
   56: *  Arguments:
   57: *  ==========
   58: *
   59: *> \param[in] UPLO
   60: *> \verbatim
   61: *>          UPLO is CHARACTER*1
   62: *>          Specifies whether the upper or lower triangular part of the
   63: *>          symmetric matrix A is stored:
   64: *>          = 'U':  Upper triangular
   65: *>          = 'L':  Lower triangular
   66: *> \endverbatim
   67: *>
   68: *> \param[in] N
   69: *> \verbatim
   70: *>          N is INTEGER
   71: *>          The order of the matrix A.  N >= 0.
   72: *> \endverbatim
   73: *>
   74: *> \param[in] NB
   75: *> \verbatim
   76: *>          NB is INTEGER
   77: *>          The maximum number of columns of the matrix A that should be
   78: *>          factored.  NB should be at least 2 to allow for 2-by-2 pivot
   79: *>          blocks.
   80: *> \endverbatim
   81: *>
   82: *> \param[out] KB
   83: *> \verbatim
   84: *>          KB is INTEGER
   85: *>          The number of columns of A that were actually factored.
   86: *>          KB is either NB-1 or NB, or N if N <= NB.
   87: *> \endverbatim
   88: *>
   89: *> \param[in,out] A
   90: *> \verbatim
   91: *>          A is DOUBLE PRECISION array, dimension (LDA,N)
   92: *>          On entry, the symmetric matrix A.
   93: *>            If UPLO = 'U': the leading N-by-N upper triangular part
   94: *>            of A contains the upper triangular part of the matrix A,
   95: *>            and the strictly lower triangular part of A is not
   96: *>            referenced.
   97: *>
   98: *>            If UPLO = 'L': the leading N-by-N lower triangular part
   99: *>            of A contains the lower triangular part of the matrix A,
  100: *>            and the strictly upper triangular part of A is not
  101: *>            referenced.
  102: *>
  103: *>          On exit, contains:
  104: *>            a) ONLY diagonal elements of the symmetric block diagonal
  105: *>               matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
  106: *>               (superdiagonal (or subdiagonal) elements of D
  107: *>                are stored on exit in array E), and
  108: *>            b) If UPLO = 'U': factor U in the superdiagonal part of A.
  109: *>               If UPLO = 'L': factor L in the subdiagonal part of A.
  110: *> \endverbatim
  111: *>
  112: *> \param[in] LDA
  113: *> \verbatim
  114: *>          LDA is INTEGER
  115: *>          The leading dimension of the array A.  LDA >= max(1,N).
  116: *> \endverbatim
  117: *>
  118: *> \param[out] E
  119: *> \verbatim
  120: *>          E is DOUBLE PRECISION array, dimension (N)
  121: *>          On exit, contains the superdiagonal (or subdiagonal)
  122: *>          elements of the symmetric block diagonal matrix D
  123: *>          with 1-by-1 or 2-by-2 diagonal blocks, where
  124: *>          If UPLO = 'U': E(i) = D(i-1,i), i=2:N, E(1) is set to 0;
  125: *>          If UPLO = 'L': E(i) = D(i+1,i), i=1:N-1, E(N) is set to 0.
  126: *>
  127: *>          NOTE: For 1-by-1 diagonal block D(k), where
  128: *>          1 <= k <= N, the element E(k) is set to 0 in both
  129: *>          UPLO = 'U' or UPLO = 'L' cases.
  130: *> \endverbatim
  131: *>
  132: *> \param[out] IPIV
  133: *> \verbatim
  134: *>          IPIV is INTEGER array, dimension (N)
  135: *>          IPIV describes the permutation matrix P in the factorization
  136: *>          of matrix A as follows. The absolute value of IPIV(k)
  137: *>          represents the index of row and column that were
  138: *>          interchanged with the k-th row and column. The value of UPLO
  139: *>          describes the order in which the interchanges were applied.
  140: *>          Also, the sign of IPIV represents the block structure of
  141: *>          the symmetric block diagonal matrix D with 1-by-1 or 2-by-2
  142: *>          diagonal blocks which correspond to 1 or 2 interchanges
  143: *>          at each factorization step.
  144: *>
  145: *>          If UPLO = 'U',
  146: *>          ( in factorization order, k decreases from N to 1 ):
  147: *>            a) A single positive entry IPIV(k) > 0 means:
  148: *>               D(k,k) is a 1-by-1 diagonal block.
  149: *>               If IPIV(k) != k, rows and columns k and IPIV(k) were
  150: *>               interchanged in the submatrix A(1:N,N-KB+1:N);
  151: *>               If IPIV(k) = k, no interchange occurred.
  152: *>
  153: *>
  154: *>            b) A pair of consecutive negative entries
  155: *>               IPIV(k) < 0 and IPIV(k-1) < 0 means:
  156: *>               D(k-1:k,k-1:k) is a 2-by-2 diagonal block.
  157: *>               (NOTE: negative entries in IPIV appear ONLY in pairs).
  158: *>               1) If -IPIV(k) != k, rows and columns
  159: *>                  k and -IPIV(k) were interchanged
  160: *>                  in the matrix A(1:N,N-KB+1:N).
  161: *>                  If -IPIV(k) = k, no interchange occurred.
  162: *>               2) If -IPIV(k-1) != k-1, rows and columns
  163: *>                  k-1 and -IPIV(k-1) were interchanged
  164: *>                  in the submatrix A(1:N,N-KB+1:N).
  165: *>                  If -IPIV(k-1) = k-1, no interchange occurred.
  166: *>
  167: *>            c) In both cases a) and b) is always ABS( IPIV(k) ) <= k.
  168: *>
  169: *>            d) NOTE: Any entry IPIV(k) is always NONZERO on output.
  170: *>
  171: *>          If UPLO = 'L',
  172: *>          ( in factorization order, k increases from 1 to N ):
  173: *>            a) A single positive entry IPIV(k) > 0 means:
  174: *>               D(k,k) is a 1-by-1 diagonal block.
  175: *>               If IPIV(k) != k, rows and columns k and IPIV(k) were
  176: *>               interchanged in the submatrix A(1:N,1:KB).
  177: *>               If IPIV(k) = k, no interchange occurred.
  178: *>
  179: *>            b) A pair of consecutive negative entries
  180: *>               IPIV(k) < 0 and IPIV(k+1) < 0 means:
  181: *>               D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
  182: *>               (NOTE: negative entries in IPIV appear ONLY in pairs).
  183: *>               1) If -IPIV(k) != k, rows and columns
  184: *>                  k and -IPIV(k) were interchanged
  185: *>                  in the submatrix A(1:N,1:KB).
  186: *>                  If -IPIV(k) = k, no interchange occurred.
  187: *>               2) If -IPIV(k+1) != k+1, rows and columns
  188: *>                  k-1 and -IPIV(k-1) were interchanged
  189: *>                  in the submatrix A(1:N,1:KB).
  190: *>                  If -IPIV(k+1) = k+1, no interchange occurred.
  191: *>
  192: *>            c) In both cases a) and b) is always ABS( IPIV(k) ) >= k.
  193: *>
  194: *>            d) NOTE: Any entry IPIV(k) is always NONZERO on output.
  195: *> \endverbatim
  196: *>
  197: *> \param[out] W
  198: *> \verbatim
  199: *>          W is DOUBLE PRECISION array, dimension (LDW,NB)
  200: *> \endverbatim
  201: *>
  202: *> \param[in] LDW
  203: *> \verbatim
  204: *>          LDW is INTEGER
  205: *>          The leading dimension of the array W.  LDW >= max(1,N).
  206: *> \endverbatim
  207: *>
  208: *> \param[out] INFO
  209: *> \verbatim
  210: *>          INFO is INTEGER
  211: *>          = 0: successful exit
  212: *>
  213: *>          < 0: If INFO = -k, the k-th argument had an illegal value
  214: *>
  215: *>          > 0: If INFO = k, the matrix A is singular, because:
  216: *>                 If UPLO = 'U': column k in the upper
  217: *>                 triangular part of A contains all zeros.
  218: *>                 If UPLO = 'L': column k in the lower
  219: *>                 triangular part of A contains all zeros.
  220: *>
  221: *>               Therefore D(k,k) is exactly zero, and superdiagonal
  222: *>               elements of column k of U (or subdiagonal elements of
  223: *>               column k of L ) are all zeros. The factorization has
  224: *>               been completed, but the block diagonal matrix D is
  225: *>               exactly singular, and division by zero will occur if
  226: *>               it is used to solve a system of equations.
  227: *>
  228: *>               NOTE: INFO only stores the first occurrence of
  229: *>               a singularity, any subsequent occurrence of singularity
  230: *>               is not stored in INFO even though the factorization
  231: *>               always completes.
  232: *> \endverbatim
  233: *
  234: *  Authors:
  235: *  ========
  236: *
  237: *> \author Univ. of Tennessee
  238: *> \author Univ. of California Berkeley
  239: *> \author Univ. of Colorado Denver
  240: *> \author NAG Ltd.
  241: *
  242: *> \ingroup doubleSYcomputational
  243: *
  244: *> \par Contributors:
  245: *  ==================
  246: *>
  247: *> \verbatim
  248: *>
  249: *>  December 2016,  Igor Kozachenko,
  250: *>                  Computer Science Division,
  251: *>                  University of California, Berkeley
  252: *>
  253: *>  September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
  254: *>                  School of Mathematics,
  255: *>                  University of Manchester
  256: *>
  257: *> \endverbatim
  258: *
  259: *  =====================================================================
  260:       SUBROUTINE DLASYF_RK( UPLO, N, NB, KB, A, LDA, E, IPIV, W, LDW,
  261:      $                      INFO )
  262: *
  263: *  -- LAPACK computational routine --
  264: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
  265: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  266: *
  267: *     .. Scalar Arguments ..
  268:       CHARACTER          UPLO
  269:       INTEGER            INFO, KB, LDA, LDW, N, NB
  270: *     ..
  271: *     .. Array Arguments ..
  272:       INTEGER            IPIV( * )
  273:       DOUBLE PRECISION   A( LDA, * ), E( * ), W( LDW, * )
  274: *     ..
  275: *
  276: *  =====================================================================
  277: *
  278: *     .. Parameters ..
  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: *     ..
  284: *     .. Local Scalars ..
  285:       LOGICAL            DONE
  286:       INTEGER            IMAX, ITEMP, J, JB, JJ, JMAX, K, KK, KW, KKW,
  287:      $                   KP, KSTEP, P, II
  288:       DOUBLE PRECISION   ABSAKK, ALPHA, COLMAX, D11, D12, D21, D22,
  289:      $                   DTEMP, R1, ROWMAX, T, SFMIN
  290: *     ..
  291: *     .. External Functions ..
  292:       LOGICAL            LSAME
  293:       INTEGER            IDAMAX
  294:       DOUBLE PRECISION   DLAMCH
  295:       EXTERNAL           LSAME, IDAMAX, DLAMCH
  296: *     ..
  297: *     .. External Subroutines ..
  298:       EXTERNAL           DCOPY, DGEMM, DGEMV, DSCAL, DSWAP
  299: *     ..
  300: *     .. Intrinsic Functions ..
  301:       INTRINSIC          ABS, MAX, MIN, SQRT
  302: *     ..
  303: *     .. Executable Statements ..
  304: *
  305:       INFO = 0
  306: *
  307: *     Initialize ALPHA for use in choosing pivot block size.
  308: *
  309:       ALPHA = ( ONE+SQRT( SEVTEN ) ) / EIGHT
  310: *
  311: *     Compute machine safe minimum
  312: *
  313:       SFMIN = DLAMCH( 'S' )
  314: *
  315:       IF( LSAME( UPLO, 'U' ) ) THEN
  316: *
  317: *        Factorize the trailing columns of A using the upper triangle
  318: *        of A and working backwards, and compute the matrix W = U12*D
  319: *        for use in updating A11
  320: *
  321: *        Initialize the first entry of array E, where superdiagonal
  322: *        elements of D are stored
  323: *
  324:          E( 1 ) = ZERO
  325: *
  326: *        K is the main loop index, decreasing from N in steps of 1 or 2
  327: *
  328:          K = N
  329:    10    CONTINUE
  330: *
  331: *        KW is the column of W which corresponds to column K of A
  332: *
  333:          KW = NB + K - N
  334: *
  335: *        Exit from loop
  336: *
  337:          IF( ( K.LE.N-NB+1 .AND. NB.LT.N ) .OR. K.LT.1 )
  338:      $      GO TO 30
  339: *
  340:          KSTEP = 1
  341:          P = K
  342: *
  343: *        Copy column K of A to column KW of W and update it
  344: *
  345:          CALL DCOPY( K, A( 1, K ), 1, W( 1, KW ), 1 )
  346:          IF( K.LT.N )
  347:      $      CALL DGEMV( 'No transpose', K, N-K, -ONE, A( 1, K+1 ),
  348:      $                  LDA, W( K, KW+1 ), LDW, ONE, W( 1, KW ), 1 )
  349: *
  350: *        Determine rows and columns to be interchanged and whether
  351: *        a 1-by-1 or 2-by-2 pivot block will be used
  352: *
  353:          ABSAKK = ABS( W( K, KW ) )
  354: *
  355: *        IMAX is the row-index of the largest off-diagonal element in
  356: *        column K, and COLMAX is its absolute value.
  357: *        Determine both COLMAX and IMAX.
  358: *
  359:          IF( K.GT.1 ) THEN
  360:             IMAX = IDAMAX( K-1, W( 1, KW ), 1 )
  361:             COLMAX = ABS( W( IMAX, KW ) )
  362:          ELSE
  363:             COLMAX = ZERO
  364:          END IF
  365: *
  366:          IF( MAX( ABSAKK, COLMAX ).EQ.ZERO ) THEN
  367: *
  368: *           Column K is zero or underflow: set INFO and continue
  369: *
  370:             IF( INFO.EQ.0 )
  371:      $         INFO = K
  372:             KP = K
  373:             CALL DCOPY( K, W( 1, KW ), 1, A( 1, K ), 1 )
  374: *
  375: *           Set E( K ) to zero
  376: *
  377:             IF( K.GT.1 )
  378:      $         E( K ) = ZERO
  379: *
  380:          ELSE
  381: *
  382: *           ============================================================
  383: *
  384: *           Test for interchange
  385: *
  386: *           Equivalent to testing for ABSAKK.GE.ALPHA*COLMAX
  387: *           (used to handle NaN and Inf)
  388: *
  389:             IF( .NOT.( ABSAKK.LT.ALPHA*COLMAX ) ) THEN
  390: *
  391: *              no interchange, use 1-by-1 pivot block
  392: *
  393:                KP = K
  394: *
  395:             ELSE
  396: *
  397:                DONE = .FALSE.
  398: *
  399: *              Loop until pivot found
  400: *
  401:    12          CONTINUE
  402: *
  403: *                 Begin pivot search loop body
  404: *
  405: *
  406: *                 Copy column IMAX to column KW-1 of W and update it
  407: *
  408:                   CALL DCOPY( IMAX, A( 1, IMAX ), 1, W( 1, KW-1 ), 1 )
  409:                   CALL DCOPY( K-IMAX, A( IMAX, IMAX+1 ), LDA,
  410:      $                        W( IMAX+1, KW-1 ), 1 )
  411: *
  412:                   IF( K.LT.N )
  413:      $               CALL DGEMV( 'No transpose', K, N-K, -ONE,
  414:      $                           A( 1, K+1 ), LDA, W( IMAX, KW+1 ), LDW,
  415:      $                           ONE, W( 1, KW-1 ), 1 )
  416: *
  417: *                 JMAX is the column-index of the largest off-diagonal
  418: *                 element in row IMAX, and ROWMAX is its absolute value.
  419: *                 Determine both ROWMAX and JMAX.
  420: *
  421:                   IF( IMAX.NE.K ) THEN
  422:                      JMAX = IMAX + IDAMAX( K-IMAX, W( IMAX+1, KW-1 ),
  423:      $                                     1 )
  424:                      ROWMAX = ABS( W( JMAX, KW-1 ) )
  425:                   ELSE
  426:                      ROWMAX = ZERO
  427:                   END IF
  428: *
  429:                   IF( IMAX.GT.1 ) THEN
  430:                      ITEMP = IDAMAX( IMAX-1, W( 1, KW-1 ), 1 )
  431:                      DTEMP = ABS( W( ITEMP, KW-1 ) )
  432:                      IF( DTEMP.GT.ROWMAX ) THEN
  433:                         ROWMAX = DTEMP
  434:                         JMAX = ITEMP
  435:                      END IF
  436:                   END IF
  437: *
  438: *                 Equivalent to testing for
  439: *                 ABS( W( IMAX, KW-1 ) ).GE.ALPHA*ROWMAX
  440: *                 (used to handle NaN and Inf)
  441: *
  442:                   IF( .NOT.(ABS( W( IMAX, KW-1 ) ).LT.ALPHA*ROWMAX ) )
  443:      $            THEN
  444: *
  445: *                    interchange rows and columns K and IMAX,
  446: *                    use 1-by-1 pivot block
  447: *
  448:                      KP = IMAX
  449: *
  450: *                    copy column KW-1 of W to column KW of W
  451: *
  452:                      CALL DCOPY( K, W( 1, KW-1 ), 1, W( 1, KW ), 1 )
  453: *
  454:                      DONE = .TRUE.
  455: *
  456: *                 Equivalent to testing for ROWMAX.EQ.COLMAX,
  457: *                 (used to handle NaN and Inf)
  458: *
  459:                   ELSE IF( ( P.EQ.JMAX ) .OR. ( ROWMAX.LE.COLMAX ) )
  460:      $            THEN
  461: *
  462: *                    interchange rows and columns K-1 and IMAX,
  463: *                    use 2-by-2 pivot block
  464: *
  465:                      KP = IMAX
  466:                      KSTEP = 2
  467:                      DONE = .TRUE.
  468:                   ELSE
  469: *
  470: *                    Pivot not found: set params and repeat
  471: *
  472:                      P = IMAX
  473:                      COLMAX = ROWMAX
  474:                      IMAX = JMAX
  475: *
  476: *                    Copy updated JMAXth (next IMAXth) column to Kth of W
  477: *
  478:                      CALL DCOPY( K, W( 1, KW-1 ), 1, W( 1, KW ), 1 )
  479: *
  480:                   END IF
  481: *
  482: *                 End pivot search loop body
  483: *
  484:                IF( .NOT. DONE ) GOTO 12
  485: *
  486:             END IF
  487: *
  488: *           ============================================================
  489: *
  490:             KK = K - KSTEP + 1
  491: *
  492: *           KKW is the column of W which corresponds to column KK of A
  493: *
  494:             KKW = NB + KK - N
  495: *
  496:             IF( ( KSTEP.EQ.2 ) .AND. ( P.NE.K ) ) THEN
  497: *
  498: *              Copy non-updated column K to column P
  499: *
  500:                CALL DCOPY( K-P, A( P+1, K ), 1, A( P, P+1 ), LDA )
  501:                CALL DCOPY( P, A( 1, K ), 1, A( 1, P ), 1 )
  502: *
  503: *              Interchange rows K and P in last N-K+1 columns of A
  504: *              and last N-K+2 columns of W
  505: *
  506:                CALL DSWAP( N-K+1, A( K, K ), LDA, A( P, K ), LDA )
  507:                CALL DSWAP( N-KK+1, W( K, KKW ), LDW, W( P, KKW ), 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 ), LDA )
  525:                CALL DSWAP( N-KK+1, W( KK, KKW ), LDW, W( KP, KKW ),
  526:      $                     LDW )
  527:             END IF
  528: *
  529:             IF( KSTEP.EQ.1 ) THEN
  530: *
  531: *              1-by-1 pivot block D(k): column KW of W now holds
  532: *
  533: *              W(k) = U(k)*D(k)
  534: *
  535: *              where U(k) is the k-th column of U
  536: *
  537: *              Store U(k) in column k of A
  538: *
  539:                CALL DCOPY( K, W( 1, KW ), 1, A( 1, K ), 1 )
  540:                IF( K.GT.1 ) THEN
  541:                   IF( ABS( A( K, K ) ).GE.SFMIN ) THEN
  542:                      R1 = ONE / A( K, K )
  543:                      CALL DSCAL( K-1, R1, A( 1, K ), 1 )
  544:                   ELSE IF( A( K, K ).NE.ZERO ) THEN
  545:                      DO 14 II = 1, K - 1
  546:                         A( II, K ) = A( II, K ) / A( K, K )
  547:    14                CONTINUE
  548:                   END IF
  549: *
  550: *                 Store the superdiagonal element of D in array E
  551: *
  552:                   E( K ) = ZERO
  553: *
  554:                END IF
  555: *
  556:             ELSE
  557: *
  558: *              2-by-2 pivot block D(k): columns KW and KW-1 of W now
  559: *              hold
  560: *
  561: *              ( W(k-1) W(k) ) = ( U(k-1) U(k) )*D(k)
  562: *
  563: *              where U(k) and U(k-1) are the k-th and (k-1)-th columns
  564: *              of U
  565: *
  566:                IF( K.GT.2 ) THEN
  567: *
  568: *                 Store U(k) and U(k-1) in columns k and k-1 of A
  569: *
  570:                   D12 = W( K-1, KW )
  571:                   D11 = W( K, KW ) / D12
  572:                   D22 = W( K-1, KW-1 ) / D12
  573:                   T = ONE / ( D11*D22-ONE )
  574:                   DO 20 J = 1, K - 2
  575:                      A( J, K-1 ) = T*( (D11*W( J, KW-1 )-W( J, KW ) ) /
  576:      $                             D12 )
  577:                      A( J, K ) = T*( ( D22*W( J, KW )-W( J, KW-1 ) ) /
  578:      $                           D12 )
  579:    20             CONTINUE
  580:                END IF
  581: *
  582: *              Copy diagonal elements of D(K) to A,
  583: *              copy superdiagonal element of D(K) to E(K) and
  584: *              ZERO out superdiagonal entry of A
  585: *
  586:                A( K-1, K-1 ) = W( K-1, KW-1 )
  587:                A( K-1, K ) = ZERO
  588:                A( K, K ) = W( K, KW )
  589:                E( K ) = W( K-1, KW )
  590:                E( K-1 ) = ZERO
  591: *
  592:             END IF
  593: *
  594: *           End column K is nonsingular
  595: *
  596:          END IF
  597: *
  598: *        Store details of the interchanges in IPIV
  599: *
  600:          IF( KSTEP.EQ.1 ) THEN
  601:             IPIV( K ) = KP
  602:          ELSE
  603:             IPIV( K ) = -P
  604:             IPIV( K-1 ) = -KP
  605:          END IF
  606: *
  607: *        Decrease K and return to the start of the main loop
  608: *
  609:          K = K - KSTEP
  610:          GO TO 10
  611: *
  612:    30    CONTINUE
  613: *
  614: *        Update the upper triangle of A11 (= A(1:k,1:k)) as
  615: *
  616: *        A11 := A11 - U12*D*U12**T = A11 - U12*W**T
  617: *
  618: *        computing blocks of NB columns at a time
  619: *
  620:          DO 50 J = ( ( K-1 ) / NB )*NB + 1, 1, -NB
  621:             JB = MIN( NB, K-J+1 )
  622: *
  623: *           Update the upper triangle of the diagonal block
  624: *
  625:             DO 40 JJ = J, J + JB - 1
  626:                CALL DGEMV( 'No transpose', JJ-J+1, N-K, -ONE,
  627:      $                     A( J, K+1 ), LDA, W( JJ, KW+1 ), LDW, ONE,
  628:      $                     A( J, JJ ), 1 )
  629:    40       CONTINUE
  630: *
  631: *           Update the rectangular superdiagonal block
  632: *
  633:             IF( J.GE.2 )
  634:      $         CALL DGEMM( 'No transpose', 'Transpose', J-1, JB,
  635:      $                  N-K, -ONE, A( 1, K+1 ), LDA, W( J, KW+1 ),
  636:      $                  LDW, ONE, A( 1, J ), LDA )
  637:    50    CONTINUE
  638: *
  639: *        Set KB to the number of columns factorized
  640: *
  641:          KB = N - K
  642: *
  643:       ELSE
  644: *
  645: *        Factorize the leading columns of A using the lower triangle
  646: *        of A and working forwards, and compute the matrix W = L21*D
  647: *        for use in updating A22
  648: *
  649: *        Initialize the unused last entry of the subdiagonal array E.
  650: *
  651:          E( N ) = ZERO
  652: *
  653: *        K is the main loop index, increasing from 1 in steps of 1 or 2
  654: *
  655:          K = 1
  656:    70   CONTINUE
  657: *
  658: *        Exit from loop
  659: *
  660:          IF( ( K.GE.NB .AND. NB.LT.N ) .OR. K.GT.N )
  661:      $      GO TO 90
  662: *
  663:          KSTEP = 1
  664:          P = K
  665: *
  666: *        Copy column K of A to column K of W and update it
  667: *
  668:          CALL DCOPY( N-K+1, A( K, K ), 1, W( K, K ), 1 )
  669:          IF( K.GT.1 )
  670:      $      CALL DGEMV( 'No transpose', N-K+1, K-1, -ONE, A( K, 1 ),
  671:      $                  LDA, W( K, 1 ), LDW, ONE, W( K, K ), 1 )
  672: *
  673: *        Determine rows and columns to be interchanged and whether
  674: *        a 1-by-1 or 2-by-2 pivot block will be used
  675: *
  676:          ABSAKK = ABS( W( K, K ) )
  677: *
  678: *        IMAX is the row-index of the largest off-diagonal element in
  679: *        column K, and COLMAX is its absolute value.
  680: *        Determine both COLMAX and IMAX.
  681: *
  682:          IF( K.LT.N ) THEN
  683:             IMAX = K + IDAMAX( N-K, W( K+1, K ), 1 )
  684:             COLMAX = ABS( W( IMAX, K ) )
  685:          ELSE
  686:             COLMAX = ZERO
  687:          END IF
  688: *
  689:          IF( MAX( ABSAKK, COLMAX ).EQ.ZERO ) THEN
  690: *
  691: *           Column K is zero or underflow: set INFO and continue
  692: *
  693:             IF( INFO.EQ.0 )
  694:      $         INFO = K
  695:             KP = K
  696:             CALL DCOPY( N-K+1, W( K, K ), 1, A( K, K ), 1 )
  697: *
  698: *           Set E( K ) to zero
  699: *
  700:             IF( K.LT.N )
  701:      $         E( K ) = ZERO
  702: *
  703:          ELSE
  704: *
  705: *           ============================================================
  706: *
  707: *           Test for interchange
  708: *
  709: *           Equivalent to testing for ABSAKK.GE.ALPHA*COLMAX
  710: *           (used to handle NaN and Inf)
  711: *
  712:             IF( .NOT.( ABSAKK.LT.ALPHA*COLMAX ) ) THEN
  713: *
  714: *              no interchange, use 1-by-1 pivot block
  715: *
  716:                KP = K
  717: *
  718:             ELSE
  719: *
  720:                DONE = .FALSE.
  721: *
  722: *              Loop until pivot found
  723: *
  724:    72          CONTINUE
  725: *
  726: *                 Begin pivot search loop body
  727: *
  728: *
  729: *                 Copy column IMAX to column K+1 of W and update it
  730: *
  731:                   CALL DCOPY( IMAX-K, A( IMAX, K ), LDA, W( K, K+1 ), 1)
  732:                   CALL DCOPY( N-IMAX+1, A( IMAX, IMAX ), 1,
  733:      $                        W( IMAX, K+1 ), 1 )
  734:                   IF( K.GT.1 )
  735:      $               CALL DGEMV( 'No transpose', N-K+1, K-1, -ONE,
  736:      $                           A( K, 1 ), LDA, W( IMAX, 1 ), LDW,
  737:      $                           ONE, W( K, K+1 ), 1 )
  738: *
  739: *                 JMAX is the column-index of the largest off-diagonal
  740: *                 element in row IMAX, and ROWMAX is its absolute value.
  741: *                 Determine both ROWMAX and JMAX.
  742: *
  743:                   IF( IMAX.NE.K ) THEN
  744:                      JMAX = K - 1 + IDAMAX( IMAX-K, W( K, K+1 ), 1 )
  745:                      ROWMAX = ABS( W( JMAX, K+1 ) )
  746:                   ELSE
  747:                      ROWMAX = ZERO
  748:                   END IF
  749: *
  750:                   IF( IMAX.LT.N ) THEN
  751:                      ITEMP = IMAX + IDAMAX( N-IMAX, W( IMAX+1, K+1 ), 1)
  752:                      DTEMP = ABS( W( ITEMP, K+1 ) )
  753:                      IF( DTEMP.GT.ROWMAX ) THEN
  754:                         ROWMAX = DTEMP
  755:                         JMAX = ITEMP
  756:                      END IF
  757:                   END IF
  758: *
  759: *                 Equivalent to testing for
  760: *                 ABS( W( IMAX, K+1 ) ).GE.ALPHA*ROWMAX
  761: *                 (used to handle NaN and Inf)
  762: *
  763:                   IF( .NOT.( ABS( W( IMAX, K+1 ) ).LT.ALPHA*ROWMAX ) )
  764:      $            THEN
  765: *
  766: *                    interchange rows and columns K and IMAX,
  767: *                    use 1-by-1 pivot block
  768: *
  769:                      KP = IMAX
  770: *
  771: *                    copy column K+1 of W to column K of W
  772: *
  773:                      CALL DCOPY( N-K+1, W( K, K+1 ), 1, W( K, K ), 1 )
  774: *
  775:                      DONE = .TRUE.
  776: *
  777: *                 Equivalent to testing for ROWMAX.EQ.COLMAX,
  778: *                 (used to handle NaN and Inf)
  779: *
  780:                   ELSE IF( ( P.EQ.JMAX ) .OR. ( ROWMAX.LE.COLMAX ) )
  781:      $            THEN
  782: *
  783: *                    interchange rows and columns K+1 and IMAX,
  784: *                    use 2-by-2 pivot block
  785: *
  786:                      KP = IMAX
  787:                      KSTEP = 2
  788:                      DONE = .TRUE.
  789:                   ELSE
  790: *
  791: *                    Pivot not found: set params and repeat
  792: *
  793:                      P = IMAX
  794:                      COLMAX = ROWMAX
  795:                      IMAX = JMAX
  796: *
  797: *                    Copy updated JMAXth (next IMAXth) column to Kth of W
  798: *
  799:                      CALL DCOPY( N-K+1, W( K, K+1 ), 1, W( K, K ), 1 )
  800: *
  801:                   END IF
  802: *
  803: *                 End pivot search loop body
  804: *
  805:                IF( .NOT. DONE ) GOTO 72
  806: *
  807:             END IF
  808: *
  809: *           ============================================================
  810: *
  811:             KK = K + KSTEP - 1
  812: *
  813:             IF( ( KSTEP.EQ.2 ) .AND. ( P.NE.K ) ) THEN
  814: *
  815: *              Copy non-updated column K to column P
  816: *
  817:                CALL DCOPY( P-K, A( K, K ), 1, A( P, K ), LDA )
  818:                CALL DCOPY( N-P+1, A( P, K ), 1, A( P, P ), 1 )
  819: *
  820: *              Interchange rows K and P in first K columns of A
  821: *              and first K+1 columns of W
  822: *
  823:                CALL DSWAP( K, A( K, 1 ), LDA, A( P, 1 ), LDA )
  824:                CALL DSWAP( KK, W( K, 1 ), LDW, W( P, 1 ), LDW )
  825:             END IF
  826: *
  827: *           Updated column KP is already stored in column KK of W
  828: *
  829:             IF( KP.NE.KK ) THEN
  830: *
  831: *              Copy non-updated column KK to column KP
  832: *
  833:                A( KP, K ) = A( KK, K )
  834:                CALL DCOPY( KP-K-1, A( K+1, KK ), 1, A( KP, K+1 ), LDA )
  835:                CALL DCOPY( N-KP+1, A( KP, KK ), 1, A( KP, KP ), 1 )
  836: *
  837: *              Interchange rows KK and KP in first KK columns of A and W
  838: *
  839:                CALL DSWAP( KK, A( KK, 1 ), LDA, A( KP, 1 ), LDA )
  840:                CALL DSWAP( KK, W( KK, 1 ), LDW, W( KP, 1 ), LDW )
  841:             END IF
  842: *
  843:             IF( KSTEP.EQ.1 ) THEN
  844: *
  845: *              1-by-1 pivot block D(k): column k of W now holds
  846: *
  847: *              W(k) = L(k)*D(k)
  848: *
  849: *              where L(k) is the k-th column of L
  850: *
  851: *              Store L(k) in column k of A
  852: *
  853:                CALL DCOPY( N-K+1, W( K, K ), 1, A( K, K ), 1 )
  854:                IF( K.LT.N ) THEN
  855:                   IF( ABS( A( K, K ) ).GE.SFMIN ) THEN
  856:                      R1 = ONE / A( K, K )
  857:                      CALL DSCAL( N-K, R1, A( K+1, K ), 1 )
  858:                   ELSE IF( A( K, K ).NE.ZERO ) THEN
  859:                      DO 74 II = K + 1, N
  860:                         A( II, K ) = A( II, K ) / A( K, K )
  861:    74                CONTINUE
  862:                   END IF
  863: *
  864: *                 Store the subdiagonal element of D in array E
  865: *
  866:                   E( K ) = ZERO
  867: *
  868:                END IF
  869: *
  870:             ELSE
  871: *
  872: *              2-by-2 pivot block D(k): columns k and k+1 of W now hold
  873: *
  874: *              ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k)
  875: *
  876: *              where L(k) and L(k+1) are the k-th and (k+1)-th columns
  877: *              of L
  878: *
  879:                IF( K.LT.N-1 ) THEN
  880: *
  881: *                 Store L(k) and L(k+1) in columns k and k+1 of A
  882: *
  883:                   D21 = W( K+1, K )
  884:                   D11 = W( K+1, K+1 ) / D21
  885:                   D22 = W( K, K ) / D21
  886:                   T = ONE / ( D11*D22-ONE )
  887:                   DO 80 J = K + 2, N
  888:                      A( J, K ) = T*( ( D11*W( J, K )-W( J, K+1 ) ) /
  889:      $                           D21 )
  890:                      A( J, K+1 ) = T*( ( D22*W( J, K+1 )-W( J, K ) ) /
  891:      $                             D21 )
  892:    80             CONTINUE
  893:                END IF
  894: *
  895: *              Copy diagonal elements of D(K) to A,
  896: *              copy subdiagonal element of D(K) to E(K) and
  897: *              ZERO out subdiagonal entry of A
  898: *
  899:                A( K, K ) = W( K, K )
  900:                A( K+1, K ) = ZERO
  901:                A( K+1, K+1 ) = W( K+1, K+1 )
  902:                E( K ) = W( K+1, K )
  903:                E( K+1 ) = ZERO
  904: *
  905:             END IF
  906: *
  907: *           End column K is nonsingular
  908: *
  909:          END IF
  910: *
  911: *        Store details of the interchanges in IPIV
  912: *
  913:          IF( KSTEP.EQ.1 ) THEN
  914:             IPIV( K ) = KP
  915:          ELSE
  916:             IPIV( K ) = -P
  917:             IPIV( K+1 ) = -KP
  918:          END IF
  919: *
  920: *        Increase K and return to the start of the main loop
  921: *
  922:          K = K + KSTEP
  923:          GO TO 70
  924: *
  925:    90    CONTINUE
  926: *
  927: *        Update the lower triangle of A22 (= A(k:n,k:n)) as
  928: *
  929: *        A22 := A22 - L21*D*L21**T = A22 - L21*W**T
  930: *
  931: *        computing blocks of NB columns at a time
  932: *
  933:          DO 110 J = K, N, NB
  934:             JB = MIN( NB, N-J+1 )
  935: *
  936: *           Update the lower triangle of the diagonal block
  937: *
  938:             DO 100 JJ = J, J + JB - 1
  939:                CALL DGEMV( 'No transpose', J+JB-JJ, K-1, -ONE,
  940:      $                     A( JJ, 1 ), LDA, W( JJ, 1 ), LDW, ONE,
  941:      $                     A( JJ, JJ ), 1 )
  942:   100       CONTINUE
  943: *
  944: *           Update the rectangular subdiagonal block
  945: *
  946:             IF( J+JB.LE.N )
  947:      $         CALL DGEMM( 'No transpose', 'Transpose', N-J-JB+1, JB,
  948:      $                     K-1, -ONE, A( J+JB, 1 ), LDA, W( J, 1 ),
  949:      $                     LDW, ONE, A( J+JB, J ), LDA )
  950:   110    CONTINUE
  951: *
  952: *        Set KB to the number of columns factorized
  953: *
  954:          KB = K - 1
  955: *
  956:       END IF
  957: *
  958:       RETURN
  959: *
  960: *     End of DLASYF_RK
  961: *
  962:       END

CVSweb interface <joel.bertrand@systella.fr>