LAPACK 3.3.0
|
00001 SUBROUTINE DGETRF ( M, N, A, LDA, IPIV, INFO) 00002 * 00003 * -- LAPACK routine (version 3.1) -- 00004 * Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd.. 00005 * March 2008 00006 * 00007 * .. Scalar Arguments .. 00008 INTEGER INFO, LDA, M, N 00009 * .. 00010 * .. Array Arguments .. 00011 INTEGER IPIV( * ) 00012 DOUBLE PRECISION A( LDA, * ) 00013 * .. 00014 * 00015 * Purpose 00016 * ======= 00017 * 00018 * DGETRF computes an LU factorization of a general M-by-N matrix A 00019 * using partial pivoting with row interchanges. 00020 * 00021 * The factorization has the form 00022 * A = P * L * U 00023 * where P is a permutation matrix, L is lower triangular with unit 00024 * diagonal elements (lower trapezoidal if m > n), and U is upper 00025 * triangular (upper trapezoidal if m < n). 00026 * 00027 * This is the left-looking Level 3 BLAS version of the algorithm. 00028 * 00029 * Arguments 00030 * ========= 00031 * 00032 * M (input) INTEGER 00033 * The number of rows of the matrix A. M >= 0. 00034 * 00035 * N (input) INTEGER 00036 * The number of columns of the matrix A. N >= 0. 00037 * 00038 * A (input/output) DOUBLE PRECISION array, dimension (LDA,N) 00039 * On entry, the M-by-N matrix to be factored. 00040 * On exit, the factors L and U from the factorization 00041 * A = P*L*U; the unit diagonal elements of L are not stored. 00042 * 00043 * LDA (input) INTEGER 00044 * The leading dimension of the array A. LDA >= max(1,M). 00045 * 00046 * IPIV (output) INTEGER array, dimension (min(M,N)) 00047 * The pivot indices; for 1 <= i <= min(M,N), row i of the 00048 * matrix was interchanged with row IPIV(i). 00049 * 00050 * INFO (output) INTEGER 00051 * = 0: successful exit 00052 * < 0: if INFO = -i, the i-th argument had an illegal value 00053 * > 0: if INFO = i, U(i,i) is exactly zero. The factorization 00054 * has been completed, but the factor U is exactly 00055 * singular, and division by zero will occur if it is used 00056 * to solve a system of equations. 00057 * 00058 * ===================================================================== 00059 * 00060 * .. Parameters .. 00061 DOUBLE PRECISION ONE 00062 PARAMETER ( ONE = 1.0D+0 ) 00063 * .. 00064 * .. Local Scalars .. 00065 INTEGER I, IINFO, J, JB, K, NB 00066 * .. 00067 * .. External Subroutines .. 00068 EXTERNAL DGEMM, DGETF2, DLASWP, DTRSM, XERBLA 00069 * .. 00070 * .. External Functions .. 00071 INTEGER ILAENV 00072 EXTERNAL ILAENV 00073 * .. 00074 * .. Intrinsic Functions .. 00075 INTRINSIC MAX, MIN 00076 * .. 00077 * .. Executable Statements .. 00078 * 00079 * Test the input parameters. 00080 * 00081 INFO = 0 00082 IF( M.LT.0 ) THEN 00083 INFO = -1 00084 ELSE IF( N.LT.0 ) THEN 00085 INFO = -2 00086 ELSE IF( LDA.LT.MAX( 1, M ) ) THEN 00087 INFO = -4 00088 END IF 00089 IF( INFO.NE.0 ) THEN 00090 CALL XERBLA( 'DGETRF', -INFO ) 00091 RETURN 00092 END IF 00093 * 00094 * Quick return if possible 00095 * 00096 IF( M.EQ.0 .OR. N.EQ.0 ) 00097 $ RETURN 00098 * 00099 * Determine the block size for this environment. 00100 * 00101 NB = ILAENV( 1, 'DGETRF', ' ', M, N, -1, -1 ) 00102 IF( NB.LE.1 .OR. NB.GE.MIN( M, N ) ) THEN 00103 * 00104 * Use unblocked code. 00105 * 00106 CALL DGETF2( M, N, A, LDA, IPIV, INFO ) 00107 00108 ELSE 00109 * 00110 * Use blocked code. 00111 * 00112 DO 20 J = 1, MIN( M, N ), NB 00113 JB = MIN( MIN( M, N )-J+1, NB ) 00114 * 00115 * Update before factoring the current panel 00116 * 00117 DO 30 K = 1, J-NB, NB 00118 * 00119 * Apply interchanges to rows K:K+NB-1. 00120 * 00121 CALL DLASWP( JB, A(1, J), LDA, K, K+NB-1, IPIV, 1 ) 00122 * 00123 * Compute block row of U. 00124 * 00125 CALL DTRSM( 'Left', 'Lower', 'No transpose', 'Unit', 00126 $ NB, JB, ONE, A( K, K ), LDA, 00127 $ A( K, J ), LDA ) 00128 * 00129 * Update trailing submatrix. 00130 * 00131 CALL DGEMM( 'No transpose', 'No transpose', 00132 $ M-K-NB+1, JB, NB, -ONE, 00133 $ A( K+NB, K ), LDA, A( K, J ), LDA, ONE, 00134 $ A( K+NB, J ), LDA ) 00135 30 CONTINUE 00136 * 00137 * Factor diagonal and subdiagonal blocks and test for exact 00138 * singularity. 00139 * 00140 CALL DGETF2( M-J+1, JB, A( J, J ), LDA, IPIV( J ), IINFO ) 00141 * 00142 * Adjust INFO and the pivot indices. 00143 * 00144 IF( INFO.EQ.0 .AND. IINFO.GT.0 ) 00145 $ INFO = IINFO + J - 1 00146 DO 10 I = J, MIN( M, J+JB-1 ) 00147 IPIV( I ) = J - 1 + IPIV( I ) 00148 10 CONTINUE 00149 * 00150 20 CONTINUE 00151 00152 * 00153 * Apply interchanges to the left-overs 00154 * 00155 DO 40 K = 1, MIN( M, N ), NB 00156 CALL DLASWP( K-1, A( 1, 1 ), LDA, K, 00157 $ MIN (K+NB-1, MIN ( M, N )), IPIV, 1 ) 00158 40 CONTINUE 00159 * 00160 * Apply update to the M+1:N columns when N > M 00161 * 00162 IF ( N.GT.M ) THEN 00163 00164 CALL DLASWP( N-M, A(1, M+1), LDA, 1, M, IPIV, 1 ) 00165 00166 DO 50 K = 1, M, NB 00167 00168 JB = MIN( M-K+1, NB ) 00169 * 00170 CALL DTRSM( 'Left', 'Lower', 'No transpose', 'Unit', 00171 $ JB, N-M, ONE, A( K, K ), LDA, 00172 $ A( K, M+1 ), LDA ) 00173 00174 * 00175 IF ( K+NB.LE.M ) THEN 00176 CALL DGEMM( 'No transpose', 'No transpose', 00177 $ M-K-NB+1, N-M, NB, -ONE, 00178 $ A( K+NB, K ), LDA, A( K, M+1 ), LDA, ONE, 00179 $ A( K+NB, M+1 ), LDA ) 00180 END IF 00181 50 CONTINUE 00182 END IF 00183 * 00184 END IF 00185 RETURN 00186 * 00187 * End of DGETRF 00188 * 00189 END