Annotation of rpl/lapack/lapack/dlaqps.f, revision 1.10

1.10    ! bertrand    1: *> \brief \b DLAQPS
        !             2: *
        !             3: *  =========== DOCUMENTATION ===========
        !             4: *
        !             5: * Online html documentation available at 
        !             6: *            http://www.netlib.org/lapack/explore-html/ 
        !             7: *
        !             8: *> \htmlonly
        !             9: *> Download DLAQPS + dependencies 
        !            10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlaqps.f"> 
        !            11: *> [TGZ]</a> 
        !            12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlaqps.f"> 
        !            13: *> [ZIP]</a> 
        !            14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlaqps.f"> 
        !            15: *> [TXT]</a>
        !            16: *> \endhtmlonly 
        !            17: *
        !            18: *  Definition:
        !            19: *  ===========
        !            20: *
        !            21: *       SUBROUTINE DLAQPS( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
        !            22: *                          VN2, AUXV, F, LDF )
        !            23: * 
        !            24: *       .. Scalar Arguments ..
        !            25: *       INTEGER            KB, LDA, LDF, M, N, NB, OFFSET
        !            26: *       ..
        !            27: *       .. Array Arguments ..
        !            28: *       INTEGER            JPVT( * )
        !            29: *       DOUBLE PRECISION   A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * ),
        !            30: *      $                   VN1( * ), VN2( * )
        !            31: *       ..
        !            32: *  
        !            33: *
        !            34: *> \par Purpose:
        !            35: *  =============
        !            36: *>
        !            37: *> \verbatim
        !            38: *>
        !            39: *> DLAQPS computes a step of QR factorization with column pivoting
        !            40: *> of a real M-by-N matrix A by using Blas-3.  It tries to factorize
        !            41: *> NB columns from A starting from the row OFFSET+1, and updates all
        !            42: *> of the matrix with Blas-3 xGEMM.
        !            43: *>
        !            44: *> In some cases, due to catastrophic cancellations, it cannot
        !            45: *> factorize NB columns.  Hence, the actual number of factorized
        !            46: *> columns is returned in KB.
        !            47: *>
        !            48: *> Block A(1:OFFSET,1:N) is accordingly pivoted, but not factorized.
        !            49: *> \endverbatim
        !            50: *
        !            51: *  Arguments:
        !            52: *  ==========
        !            53: *
        !            54: *> \param[in] M
        !            55: *> \verbatim
        !            56: *>          M is INTEGER
        !            57: *>          The number of rows of the matrix A. M >= 0.
        !            58: *> \endverbatim
        !            59: *>
        !            60: *> \param[in] N
        !            61: *> \verbatim
        !            62: *>          N is INTEGER
        !            63: *>          The number of columns of the matrix A. N >= 0
        !            64: *> \endverbatim
        !            65: *>
        !            66: *> \param[in] OFFSET
        !            67: *> \verbatim
        !            68: *>          OFFSET is INTEGER
        !            69: *>          The number of rows of A that have been factorized in
        !            70: *>          previous steps.
        !            71: *> \endverbatim
        !            72: *>
        !            73: *> \param[in] NB
        !            74: *> \verbatim
        !            75: *>          NB is INTEGER
        !            76: *>          The number of columns to factorize.
        !            77: *> \endverbatim
        !            78: *>
        !            79: *> \param[out] KB
        !            80: *> \verbatim
        !            81: *>          KB is INTEGER
        !            82: *>          The number of columns actually factorized.
        !            83: *> \endverbatim
        !            84: *>
        !            85: *> \param[in,out] A
        !            86: *> \verbatim
        !            87: *>          A is DOUBLE PRECISION array, dimension (LDA,N)
        !            88: *>          On entry, the M-by-N matrix A.
        !            89: *>          On exit, block A(OFFSET+1:M,1:KB) is the triangular
        !            90: *>          factor obtained and block A(1:OFFSET,1:N) has been
        !            91: *>          accordingly pivoted, but no factorized.
        !            92: *>          The rest of the matrix, block A(OFFSET+1:M,KB+1:N) has
        !            93: *>          been updated.
        !            94: *> \endverbatim
        !            95: *>
        !            96: *> \param[in] LDA
        !            97: *> \verbatim
        !            98: *>          LDA is INTEGER
        !            99: *>          The leading dimension of the array A. LDA >= max(1,M).
        !           100: *> \endverbatim
        !           101: *>
        !           102: *> \param[in,out] JPVT
        !           103: *> \verbatim
        !           104: *>          JPVT is INTEGER array, dimension (N)
        !           105: *>          JPVT(I) = K <==> Column K of the full matrix A has been
        !           106: *>          permuted into position I in AP.
        !           107: *> \endverbatim
        !           108: *>
        !           109: *> \param[out] TAU
        !           110: *> \verbatim
        !           111: *>          TAU is DOUBLE PRECISION array, dimension (KB)
        !           112: *>          The scalar factors of the elementary reflectors.
        !           113: *> \endverbatim
        !           114: *>
        !           115: *> \param[in,out] VN1
        !           116: *> \verbatim
        !           117: *>          VN1 is DOUBLE PRECISION array, dimension (N)
        !           118: *>          The vector with the partial column norms.
        !           119: *> \endverbatim
        !           120: *>
        !           121: *> \param[in,out] VN2
        !           122: *> \verbatim
        !           123: *>          VN2 is DOUBLE PRECISION array, dimension (N)
        !           124: *>          The vector with the exact column norms.
        !           125: *> \endverbatim
        !           126: *>
        !           127: *> \param[in,out] AUXV
        !           128: *> \verbatim
        !           129: *>          AUXV is DOUBLE PRECISION array, dimension (NB)
        !           130: *>          Auxiliar vector.
        !           131: *> \endverbatim
        !           132: *>
        !           133: *> \param[in,out] F
        !           134: *> \verbatim
        !           135: *>          F is DOUBLE PRECISION array, dimension (LDF,NB)
        !           136: *>          Matrix F**T = L*Y**T*A.
        !           137: *> \endverbatim
        !           138: *>
        !           139: *> \param[in] LDF
        !           140: *> \verbatim
        !           141: *>          LDF is INTEGER
        !           142: *>          The leading dimension of the array F. LDF >= max(1,N).
        !           143: *> \endverbatim
        !           144: *
        !           145: *  Authors:
        !           146: *  ========
        !           147: *
        !           148: *> \author Univ. of Tennessee 
        !           149: *> \author Univ. of California Berkeley 
        !           150: *> \author Univ. of Colorado Denver 
        !           151: *> \author NAG Ltd. 
        !           152: *
        !           153: *> \date November 2011
        !           154: *
        !           155: *> \ingroup doubleOTHERauxiliary
        !           156: *
        !           157: *> \par Contributors:
        !           158: *  ==================
        !           159: *>
        !           160: *>    G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain
        !           161: *>    X. Sun, Computer Science Dept., Duke University, USA
        !           162: *> \n
        !           163: *>  Partial column norm updating strategy modified on April 2011
        !           164: *>    Z. Drmac and Z. Bujanovic, Dept. of Mathematics,
        !           165: *>    University of Zagreb, Croatia.
        !           166: *
        !           167: *> \par References:
        !           168: *  ================
        !           169: *>
        !           170: *> LAPACK Working Note 176
        !           171: *
        !           172: *> \htmlonly
        !           173: *> <a href="http://www.netlib.org/lapack/lawnspdf/lawn176.pdf">[PDF]</a> 
        !           174: *> \endhtmlonly 
        !           175: *
        !           176: *  =====================================================================
1.1       bertrand  177:       SUBROUTINE DLAQPS( M, N, OFFSET, NB, KB, A, LDA, JPVT, TAU, VN1,
                    178:      $                   VN2, AUXV, F, LDF )
                    179: *
1.10    ! bertrand  180: *  -- LAPACK auxiliary routine (version 3.4.0) --
1.1       bertrand  181: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
                    182: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
1.10    ! bertrand  183: *     November 2011
1.1       bertrand  184: *
                    185: *     .. Scalar Arguments ..
                    186:       INTEGER            KB, LDA, LDF, M, N, NB, OFFSET
                    187: *     ..
                    188: *     .. Array Arguments ..
                    189:       INTEGER            JPVT( * )
                    190:       DOUBLE PRECISION   A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * ),
                    191:      $                   VN1( * ), VN2( * )
                    192: *     ..
                    193: *
                    194: *  =====================================================================
                    195: *
                    196: *     .. Parameters ..
                    197:       DOUBLE PRECISION   ZERO, ONE
                    198:       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
                    199: *     ..
                    200: *     .. Local Scalars ..
                    201:       INTEGER            ITEMP, J, K, LASTRK, LSTICC, PVT, RK
                    202:       DOUBLE PRECISION   AKK, TEMP, TEMP2, TOL3Z
                    203: *     ..
                    204: *     .. External Subroutines ..
1.5       bertrand  205:       EXTERNAL           DGEMM, DGEMV, DLARFG, DSWAP
1.1       bertrand  206: *     ..
                    207: *     .. Intrinsic Functions ..
                    208:       INTRINSIC          ABS, DBLE, MAX, MIN, NINT, SQRT
                    209: *     ..
                    210: *     .. External Functions ..
                    211:       INTEGER            IDAMAX
                    212:       DOUBLE PRECISION   DLAMCH, DNRM2
                    213:       EXTERNAL           IDAMAX, DLAMCH, DNRM2
                    214: *     ..
                    215: *     .. Executable Statements ..
                    216: *
                    217:       LASTRK = MIN( M, N+OFFSET )
                    218:       LSTICC = 0
                    219:       K = 0
                    220:       TOL3Z = SQRT(DLAMCH('Epsilon'))
                    221: *
                    222: *     Beginning of while loop.
                    223: *
                    224:    10 CONTINUE
                    225:       IF( ( K.LT.NB ) .AND. ( LSTICC.EQ.0 ) ) THEN
                    226:          K = K + 1
                    227:          RK = OFFSET + K
                    228: *
                    229: *        Determine ith pivot column and swap if necessary
                    230: *
                    231:          PVT = ( K-1 ) + IDAMAX( N-K+1, VN1( K ), 1 )
                    232:          IF( PVT.NE.K ) THEN
                    233:             CALL DSWAP( M, A( 1, PVT ), 1, A( 1, K ), 1 )
                    234:             CALL DSWAP( K-1, F( PVT, 1 ), LDF, F( K, 1 ), LDF )
                    235:             ITEMP = JPVT( PVT )
                    236:             JPVT( PVT ) = JPVT( K )
                    237:             JPVT( K ) = ITEMP
                    238:             VN1( PVT ) = VN1( K )
                    239:             VN2( PVT ) = VN2( K )
                    240:          END IF
                    241: *
                    242: *        Apply previous Householder reflectors to column K:
1.9       bertrand  243: *        A(RK:M,K) := A(RK:M,K) - A(RK:M,1:K-1)*F(K,1:K-1)**T.
1.1       bertrand  244: *
                    245:          IF( K.GT.1 ) THEN
                    246:             CALL DGEMV( 'No transpose', M-RK+1, K-1, -ONE, A( RK, 1 ),
                    247:      $                  LDA, F( K, 1 ), LDF, ONE, A( RK, K ), 1 )
                    248:          END IF
                    249: *
                    250: *        Generate elementary reflector H(k).
                    251: *
                    252:          IF( RK.LT.M ) THEN
1.5       bertrand  253:             CALL DLARFG( M-RK+1, A( RK, K ), A( RK+1, K ), 1, TAU( K ) )
1.1       bertrand  254:          ELSE
1.5       bertrand  255:             CALL DLARFG( 1, A( RK, K ), A( RK, K ), 1, TAU( K ) )
1.1       bertrand  256:          END IF
                    257: *
                    258:          AKK = A( RK, K )
                    259:          A( RK, K ) = ONE
                    260: *
                    261: *        Compute Kth column of F:
                    262: *
1.9       bertrand  263: *        Compute  F(K+1:N,K) := tau(K)*A(RK:M,K+1:N)**T*A(RK:M,K).
1.1       bertrand  264: *
                    265:          IF( K.LT.N ) THEN
                    266:             CALL DGEMV( 'Transpose', M-RK+1, N-K, TAU( K ),
                    267:      $                  A( RK, K+1 ), LDA, A( RK, K ), 1, ZERO,
                    268:      $                  F( K+1, K ), 1 )
                    269:          END IF
                    270: *
                    271: *        Padding F(1:K,K) with zeros.
                    272: *
                    273:          DO 20 J = 1, K
                    274:             F( J, K ) = ZERO
                    275:    20    CONTINUE
                    276: *
                    277: *        Incremental updating of F:
1.9       bertrand  278: *        F(1:N,K) := F(1:N,K) - tau(K)*F(1:N,1:K-1)*A(RK:M,1:K-1)**T
1.1       bertrand  279: *                    *A(RK:M,K).
                    280: *
                    281:          IF( K.GT.1 ) THEN
                    282:             CALL DGEMV( 'Transpose', M-RK+1, K-1, -TAU( K ), A( RK, 1 ),
                    283:      $                  LDA, A( RK, K ), 1, ZERO, AUXV( 1 ), 1 )
                    284: *
                    285:             CALL DGEMV( 'No transpose', N, K-1, ONE, F( 1, 1 ), LDF,
                    286:      $                  AUXV( 1 ), 1, ONE, F( 1, K ), 1 )
                    287:          END IF
                    288: *
                    289: *        Update the current row of A:
1.9       bertrand  290: *        A(RK,K+1:N) := A(RK,K+1:N) - A(RK,1:K)*F(K+1:N,1:K)**T.
1.1       bertrand  291: *
                    292:          IF( K.LT.N ) THEN
                    293:             CALL DGEMV( 'No transpose', N-K, K, -ONE, F( K+1, 1 ), LDF,
                    294:      $                  A( RK, 1 ), LDA, ONE, A( RK, K+1 ), LDA )
                    295:          END IF
                    296: *
                    297: *        Update partial column norms.
                    298: *
                    299:          IF( RK.LT.LASTRK ) THEN
                    300:             DO 30 J = K + 1, N
                    301:                IF( VN1( J ).NE.ZERO ) THEN
                    302: *
                    303: *                 NOTE: The following 4 lines follow from the analysis in
                    304: *                 Lapack Working Note 176.
                    305: *
                    306:                   TEMP = ABS( A( RK, J ) ) / VN1( J )
                    307:                   TEMP = MAX( ZERO, ( ONE+TEMP )*( ONE-TEMP ) )
                    308:                   TEMP2 = TEMP*( VN1( J ) / VN2( J ) )**2
                    309:                   IF( TEMP2 .LE. TOL3Z ) THEN
                    310:                      VN2( J ) = DBLE( LSTICC )
                    311:                      LSTICC = J
                    312:                   ELSE
                    313:                      VN1( J ) = VN1( J )*SQRT( TEMP )
                    314:                   END IF
                    315:                END IF
                    316:    30       CONTINUE
                    317:          END IF
                    318: *
                    319:          A( RK, K ) = AKK
                    320: *
                    321: *        End of while loop.
                    322: *
                    323:          GO TO 10
                    324:       END IF
                    325:       KB = K
                    326:       RK = OFFSET + KB
                    327: *
                    328: *     Apply the block reflector to the rest of the matrix:
                    329: *     A(OFFSET+KB+1:M,KB+1:N) := A(OFFSET+KB+1:M,KB+1:N) -
1.9       bertrand  330: *                         A(OFFSET+KB+1:M,1:KB)*F(KB+1:N,1:KB)**T.
1.1       bertrand  331: *
                    332:       IF( KB.LT.MIN( N, M-OFFSET ) ) THEN
                    333:          CALL DGEMM( 'No transpose', 'Transpose', M-RK, N-KB, KB, -ONE,
                    334:      $               A( RK+1, 1 ), LDA, F( KB+1, 1 ), LDF, ONE,
                    335:      $               A( RK+1, KB+1 ), LDA )
                    336:       END IF
                    337: *
                    338: *     Recomputation of difficult columns.
                    339: *
                    340:    40 CONTINUE
                    341:       IF( LSTICC.GT.0 ) THEN
                    342:          ITEMP = NINT( VN2( LSTICC ) )
                    343:          VN1( LSTICC ) = DNRM2( M-RK, A( RK+1, LSTICC ), 1 )
                    344: *
                    345: *        NOTE: The computation of VN1( LSTICC ) relies on the fact that 
                    346: *        SNRM2 does not fail on vectors with norm below the value of
                    347: *        SQRT(DLAMCH('S')) 
                    348: *
                    349:          VN2( LSTICC ) = VN1( LSTICC )
                    350:          LSTICC = ITEMP
                    351:          GO TO 40
                    352:       END IF
                    353: *
                    354:       RETURN
                    355: *
                    356: *     End of DLAQPS
                    357: *
                    358:       END

CVSweb interface <joel.bertrand@systella.fr>