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.

dgetrff.f 4.6 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156
  1. SUBROUTINE DGETRFF( M, N, A, LDA, IPIV, INFO )
  2. *
  3. * -- LAPACK routine (version 3.0) --
  4. * Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
  5. * Courant Institute, Argonne National Lab, and Rice University
  6. * March 31, 1993
  7. *
  8. * .. Scalar Arguments ..
  9. INTEGER INFO, LDA, M, N
  10. * ..
  11. * .. Array Arguments ..
  12. INTEGER IPIV( * )
  13. DOUBLE PRECISION A( LDA, * )
  14. * ..
  15. *
  16. * Purpose
  17. * =======
  18. *
  19. * DGETRF computes an LU factorization of a general M-by-N matrix A
  20. * using partial pivoting with row interchanges.
  21. *
  22. * The factorization has the form
  23. * A = P * L * U
  24. * where P is a permutation matrix, L is lower triangular with unit
  25. * diagonal elements (lower trapezoidal if m > n), and U is upper
  26. * triangular (upper trapezoidal if m < n).
  27. *
  28. * This is the right-looking Level 3 BLAS version of the algorithm.
  29. *
  30. * Arguments
  31. * =========
  32. *
  33. * M (input) INTEGER
  34. * The number of rows of the matrix A. M >= 0.
  35. *
  36. * N (input) INTEGER
  37. * The number of columns of the matrix A. N >= 0.
  38. *
  39. * A (input/output) DOUBLE PRECISION array, dimension (LDA,N)
  40. * On entry, the M-by-N matrix to be factored.
  41. * On exit, the factors L and U from the factorization
  42. * A = P*L*U; the unit diagonal elements of L are not stored.
  43. *
  44. * LDA (input) INTEGER
  45. * The leading dimension of the array A. LDA >= max(1,M).
  46. *
  47. * IPIV (output) INTEGER array, dimension (min(M,N))
  48. * The pivot indices; for 1 <= i <= min(M,N), row i of the
  49. * matrix was interchanged with row IPIV(i).
  50. *
  51. * INFO (output) INTEGER
  52. * = 0: successful exit
  53. * < 0: if INFO = -i, the i-th argument had an illegal value
  54. * > 0: if INFO = i, U(i,i) is exactly zero. The factorization
  55. * has been completed, but the factor U is exactly
  56. * singular, and division by zero will occur if it is used
  57. * to solve a system of equations.
  58. *
  59. * =====================================================================
  60. *
  61. * .. Parameters ..
  62. DOUBLE PRECISION ONE
  63. PARAMETER ( ONE = 1.0D+0 )
  64. * ..
  65. * .. Local Scalars ..
  66. INTEGER I, IINFO, J, JB, NB
  67. * ..
  68. * .. External Subroutines ..
  69. EXTERNAL DGEMM, DGETF2, DLASWP, DTRSM, XERBLA
  70. * ..
  71. * .. Intrinsic Functions ..
  72. INTRINSIC MAX, MIN
  73. * ..
  74. * .. Executable Statements ..
  75. *
  76. * Test the input parameters.
  77. *
  78. INFO = 0
  79. IF( M.LT.0 ) THEN
  80. INFO = -1
  81. ELSE IF( N.LT.0 ) THEN
  82. INFO = -2
  83. ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
  84. INFO = -4
  85. END IF
  86. IF( INFO.NE.0 ) THEN
  87. CALL XERBLA( 'DGETRF', -INFO )
  88. RETURN
  89. END IF
  90. *
  91. * Quick return if possible
  92. *
  93. IF( M.EQ.0 .OR. N.EQ.0 )
  94. $ RETURN
  95. *
  96. * Determine the block size for this environment.
  97. *
  98. NB = 64
  99. IF( NB.LE.1 .OR. NB.GE.MIN( M, N ) ) THEN
  100. *
  101. * Use unblocked code.
  102. *
  103. CALL DGETF2( M, N, A, LDA, IPIV, INFO )
  104. ELSE
  105. *
  106. * Use blocked code.
  107. *
  108. DO 20 J = 1, MIN( M, N ), NB
  109. JB = MIN( MIN( M, N )-J+1, NB )
  110. *
  111. * Factor diagonal and subdiagonal blocks and test for exact
  112. * singularity.
  113. *
  114. CALL DGETF2( M-J+1, JB, A( J, J ), LDA, IPIV( J ), IINFO )
  115. *
  116. * Adjust INFO and the pivot indices.
  117. *
  118. IF( INFO.EQ.0 .AND. IINFO.GT.0 )
  119. $ INFO = IINFO + J - 1
  120. DO 10 I = J, MIN( M, J+JB-1 )
  121. IPIV( I ) = J - 1 + IPIV( I )
  122. 10 CONTINUE
  123. *
  124. * Apply interchanges to columns 1:J-1.
  125. *
  126. CALL DLASWP( J-1, A, LDA, J, J+JB-1, IPIV, 1 )
  127. *
  128. IF( J+JB.LE.N ) THEN
  129. *
  130. * Apply interchanges to columns J+JB:N.
  131. *
  132. CALL DLASWP( N-J-JB+1, A( 1, J+JB ), LDA, J, J+JB-1,
  133. $ IPIV, 1 )
  134. *
  135. * Compute block row of U.
  136. *
  137. CALL DTRSM( 'Left', 'Lower', 'No transpose', 'Unit', JB,
  138. $ N-J-JB+1, ONE, A( J, J ), LDA, A( J, J+JB ),
  139. $ LDA )
  140. IF( J+JB.LE.M ) THEN
  141. *
  142. * Update trailing submatrix.
  143. *
  144. CALL DGEMM( 'No transpose', 'No transpose', M-J-JB+1,
  145. $ N-J-JB+1, JB, -ONE, A( J+JB, J ), LDA,
  146. $ A( J, J+JB ), LDA, ONE, A( J+JB, J+JB ),
  147. $ LDA )
  148. END IF
  149. END IF
  150. 20 CONTINUE
  151. END IF
  152. RETURN
  153. *
  154. * End of DGETRF
  155. *
  156. END