Annotation of rpl/lapack/lapack/dgetrf.f, revision 1.8

1.8     ! bertrand    1: *> \brief \b DGETRF
        !             2: *
        !             3: *  =========== DOCUMENTATION ===========
        !             4: *
        !             5: * Online html documentation available at 
        !             6: *            http://www.netlib.org/lapack/explore-html/ 
        !             7: *
        !             8: *> \htmlonly
        !             9: *> Download DGETRF + dependencies 
        !            10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgetrf.f"> 
        !            11: *> [TGZ]</a> 
        !            12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgetrf.f"> 
        !            13: *> [ZIP]</a> 
        !            14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgetrf.f"> 
        !            15: *> [TXT]</a>
        !            16: *> \endhtmlonly 
        !            17: *
        !            18: *  Definition:
        !            19: *  ===========
        !            20: *
        !            21: *       SUBROUTINE DGETRF( M, N, A, LDA, IPIV, INFO )
        !            22: * 
        !            23: *       .. Scalar Arguments ..
        !            24: *       INTEGER            INFO, LDA, M, N
        !            25: *       ..
        !            26: *       .. Array Arguments ..
        !            27: *       INTEGER            IPIV( * )
        !            28: *       DOUBLE PRECISION   A( LDA, * )
        !            29: *       ..
        !            30: *  
        !            31: *
        !            32: *> \par Purpose:
        !            33: *  =============
        !            34: *>
        !            35: *> \verbatim
        !            36: *>
        !            37: *> DGETRF computes an LU factorization of a general M-by-N matrix A
        !            38: *> using partial pivoting with row interchanges.
        !            39: *>
        !            40: *> The factorization has the form
        !            41: *>    A = P * L * U
        !            42: *> where P is a permutation matrix, L is lower triangular with unit
        !            43: *> diagonal elements (lower trapezoidal if m > n), and U is upper
        !            44: *> triangular (upper trapezoidal if m < n).
        !            45: *>
        !            46: *> This is the right-looking Level 3 BLAS version of the algorithm.
        !            47: *> \endverbatim
        !            48: *
        !            49: *  Arguments:
        !            50: *  ==========
        !            51: *
        !            52: *> \param[in] M
        !            53: *> \verbatim
        !            54: *>          M is INTEGER
        !            55: *>          The number of rows of the matrix A.  M >= 0.
        !            56: *> \endverbatim
        !            57: *>
        !            58: *> \param[in] N
        !            59: *> \verbatim
        !            60: *>          N is INTEGER
        !            61: *>          The number of columns of the matrix A.  N >= 0.
        !            62: *> \endverbatim
        !            63: *>
        !            64: *> \param[in,out] A
        !            65: *> \verbatim
        !            66: *>          A is DOUBLE PRECISION array, dimension (LDA,N)
        !            67: *>          On entry, the M-by-N matrix to be factored.
        !            68: *>          On exit, the factors L and U from the factorization
        !            69: *>          A = P*L*U; the unit diagonal elements of L are not stored.
        !            70: *> \endverbatim
        !            71: *>
        !            72: *> \param[in] LDA
        !            73: *> \verbatim
        !            74: *>          LDA is INTEGER
        !            75: *>          The leading dimension of the array A.  LDA >= max(1,M).
        !            76: *> \endverbatim
        !            77: *>
        !            78: *> \param[out] IPIV
        !            79: *> \verbatim
        !            80: *>          IPIV is INTEGER array, dimension (min(M,N))
        !            81: *>          The pivot indices; for 1 <= i <= min(M,N), row i of the
        !            82: *>          matrix was interchanged with row IPIV(i).
        !            83: *> \endverbatim
        !            84: *>
        !            85: *> \param[out] INFO
        !            86: *> \verbatim
        !            87: *>          INFO is INTEGER
        !            88: *>          = 0:  successful exit
        !            89: *>          < 0:  if INFO = -i, the i-th argument had an illegal value
        !            90: *>          > 0:  if INFO = i, U(i,i) is exactly zero. The factorization
        !            91: *>                has been completed, but the factor U is exactly
        !            92: *>                singular, and division by zero will occur if it is used
        !            93: *>                to solve a system of equations.
        !            94: *> \endverbatim
        !            95: *
        !            96: *  Authors:
        !            97: *  ========
        !            98: *
        !            99: *> \author Univ. of Tennessee 
        !           100: *> \author Univ. of California Berkeley 
        !           101: *> \author Univ. of Colorado Denver 
        !           102: *> \author NAG Ltd. 
        !           103: *
        !           104: *> \date November 2011
        !           105: *
        !           106: *> \ingroup doubleGEcomputational
        !           107: *
        !           108: *  =====================================================================
1.1       bertrand  109:       SUBROUTINE DGETRF( M, N, A, LDA, IPIV, INFO )
                    110: *
1.8     ! bertrand  111: *  -- LAPACK computational routine (version 3.4.0) --
1.1       bertrand  112: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
                    113: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
1.8     ! bertrand  114: *     November 2011
1.1       bertrand  115: *
                    116: *     .. Scalar Arguments ..
                    117:       INTEGER            INFO, LDA, M, N
                    118: *     ..
                    119: *     .. Array Arguments ..
                    120:       INTEGER            IPIV( * )
                    121:       DOUBLE PRECISION   A( LDA, * )
                    122: *     ..
                    123: *
                    124: *  =====================================================================
                    125: *
                    126: *     .. Parameters ..
                    127:       DOUBLE PRECISION   ONE
                    128:       PARAMETER          ( ONE = 1.0D+0 )
                    129: *     ..
                    130: *     .. Local Scalars ..
                    131:       INTEGER            I, IINFO, J, JB, NB
                    132: *     ..
                    133: *     .. External Subroutines ..
                    134:       EXTERNAL           DGEMM, DGETF2, DLASWP, DTRSM, XERBLA
                    135: *     ..
                    136: *     .. External Functions ..
                    137:       INTEGER            ILAENV
                    138:       EXTERNAL           ILAENV
                    139: *     ..
                    140: *     .. Intrinsic Functions ..
                    141:       INTRINSIC          MAX, MIN
                    142: *     ..
                    143: *     .. Executable Statements ..
                    144: *
                    145: *     Test the input parameters.
                    146: *
                    147:       INFO = 0
                    148:       IF( M.LT.0 ) THEN
                    149:          INFO = -1
                    150:       ELSE IF( N.LT.0 ) THEN
                    151:          INFO = -2
                    152:       ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
                    153:          INFO = -4
                    154:       END IF
                    155:       IF( INFO.NE.0 ) THEN
                    156:          CALL XERBLA( 'DGETRF', -INFO )
                    157:          RETURN
                    158:       END IF
                    159: *
                    160: *     Quick return if possible
                    161: *
                    162:       IF( M.EQ.0 .OR. N.EQ.0 )
                    163:      $   RETURN
                    164: *
                    165: *     Determine the block size for this environment.
                    166: *
                    167:       NB = ILAENV( 1, 'DGETRF', ' ', M, N, -1, -1 )
                    168:       IF( NB.LE.1 .OR. NB.GE.MIN( M, N ) ) THEN
                    169: *
                    170: *        Use unblocked code.
                    171: *
                    172:          CALL DGETF2( M, N, A, LDA, IPIV, INFO )
                    173:       ELSE
                    174: *
                    175: *        Use blocked code.
                    176: *
                    177:          DO 20 J = 1, MIN( M, N ), NB
                    178:             JB = MIN( MIN( M, N )-J+1, NB )
                    179: *
                    180: *           Factor diagonal and subdiagonal blocks and test for exact
                    181: *           singularity.
                    182: *
                    183:             CALL DGETF2( M-J+1, JB, A( J, J ), LDA, IPIV( J ), IINFO )
                    184: *
                    185: *           Adjust INFO and the pivot indices.
                    186: *
                    187:             IF( INFO.EQ.0 .AND. IINFO.GT.0 )
                    188:      $         INFO = IINFO + J - 1
                    189:             DO 10 I = J, MIN( M, J+JB-1 )
                    190:                IPIV( I ) = J - 1 + IPIV( I )
                    191:    10       CONTINUE
                    192: *
                    193: *           Apply interchanges to columns 1:J-1.
                    194: *
                    195:             CALL DLASWP( J-1, A, LDA, J, J+JB-1, IPIV, 1 )
                    196: *
                    197:             IF( J+JB.LE.N ) THEN
                    198: *
                    199: *              Apply interchanges to columns J+JB:N.
                    200: *
                    201:                CALL DLASWP( N-J-JB+1, A( 1, J+JB ), LDA, J, J+JB-1,
                    202:      $                      IPIV, 1 )
                    203: *
                    204: *              Compute block row of U.
                    205: *
                    206:                CALL DTRSM( 'Left', 'Lower', 'No transpose', 'Unit', JB,
                    207:      $                     N-J-JB+1, ONE, A( J, J ), LDA, A( J, J+JB ),
                    208:      $                     LDA )
                    209:                IF( J+JB.LE.M ) THEN
                    210: *
                    211: *                 Update trailing submatrix.
                    212: *
                    213:                   CALL DGEMM( 'No transpose', 'No transpose', M-J-JB+1,
                    214:      $                        N-J-JB+1, JB, -ONE, A( J+JB, J ), LDA,
                    215:      $                        A( J, J+JB ), LDA, ONE, A( J+JB, J+JB ),
                    216:      $                        LDA )
                    217:                END IF
                    218:             END IF
                    219:    20    CONTINUE
                    220:       END IF
                    221:       RETURN
                    222: *
                    223: *     End of DGETRF
                    224: *
                    225:       END

CVSweb interface <joel.bertrand@systella.fr>