You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

dgetrf.f 6.7 kB

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