dgetrf.c 8.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269
  1. /* dgetrf.f -- translated by f2c (version 20061008).
  2. You must link the resulting object file with libf2c:
  3. on Microsoft Windows system, link with libf2c.lib;
  4. on Linux or Unix systems, link with .../path/to/libf2c.a -lm
  5. or, if you install libf2c.a in a standard place, with -lf2c -lm
  6. -- in that order, at the end of the command line, as in
  7. cc *.o -lf2c -lm
  8. Source for libf2c is in /netlib/f2c/libf2c.zip, e.g.,
  9. http://www.netlib.org/f2c/libf2c.zip
  10. */
  11. #include "f2c.h"
  12. #include "blaswrap.h"
  13. /* Table of constant values */
  14. static integer c__1 = 1;
  15. static doublereal c_b12 = 1.;
  16. static doublereal c_b15 = -1.;
  17. /* Subroutine */ int _starpu_dgetrf_(integer *m, integer *n, doublereal *a, integer *
  18. lda, integer *ipiv, integer *info)
  19. {
  20. /* System generated locals */
  21. integer a_dim1, a_offset, i__1, i__2, i__3;
  22. doublereal d__1;
  23. /* Local variables */
  24. integer i__, j, ipivstart, jpivstart, jp;
  25. doublereal tmp;
  26. extern /* Subroutine */ int _starpu_dscal_(integer *, doublereal *, doublereal *,
  27. integer *), _starpu_dgemm_(char *, char *, integer *, integer *, integer *
  28. , doublereal *, doublereal *, integer *, doublereal *, integer *,
  29. doublereal *, doublereal *, integer *);
  30. integer kcols;
  31. doublereal sfmin;
  32. integer nstep;
  33. extern /* Subroutine */ int _starpu_dtrsm_(char *, char *, char *, char *,
  34. integer *, integer *, doublereal *, doublereal *, integer *,
  35. doublereal *, integer *);
  36. integer kahead;
  37. extern doublereal _starpu_dlamch_(char *);
  38. extern integer _starpu_idamax_(integer *, doublereal *, integer *);
  39. extern logical _starpu_disnan_(doublereal *);
  40. extern /* Subroutine */ int _starpu_xerbla_(char *, integer *);
  41. integer npived;
  42. extern /* Subroutine */ int _starpu_dlaswp_(integer *, doublereal *, integer *,
  43. integer *, integer *, integer *, integer *);
  44. integer kstart, ntopiv;
  45. /* -- LAPACK routine (version 3.X) -- */
  46. /* Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd.. */
  47. /* May 2008 */
  48. /* .. Scalar Arguments .. */
  49. /* .. */
  50. /* .. Array Arguments .. */
  51. /* .. */
  52. /* Purpose */
  53. /* ======= */
  54. /* DGETRF computes an LU factorization of a general M-by-N matrix A */
  55. /* using partial pivoting with row interchanges. */
  56. /* The factorization has the form */
  57. /* A = P * L * U */
  58. /* where P is a permutation matrix, L is lower triangular with unit */
  59. /* diagonal elements (lower trapezoidal if m > n), and U is upper */
  60. /* triangular (upper trapezoidal if m < n). */
  61. /* This code implements an iterative version of Sivan Toledo's recursive */
  62. /* LU algorithm[1]. For square matrices, this iterative versions should */
  63. /* be within a factor of two of the optimum number of memory transfers. */
  64. /* The pattern is as follows, with the large blocks of U being updated */
  65. /* in one call to DTRSM, and the dotted lines denoting sections that */
  66. /* have had all pending permutations applied: */
  67. /* 1 2 3 4 5 6 7 8 */
  68. /* +-+-+---+-------+------ */
  69. /* | |1| | | */
  70. /* |.+-+ 2 | | */
  71. /* | | | | | */
  72. /* |.|.+-+-+ 4 | */
  73. /* | | | |1| | */
  74. /* | | |.+-+ | */
  75. /* | | | | | | */
  76. /* |.|.|.|.+-+-+---+ 8 */
  77. /* | | | | | |1| | */
  78. /* | | | | |.+-+ 2 | */
  79. /* | | | | | | | | */
  80. /* | | | | |.|.+-+-+ */
  81. /* | | | | | | | |1| */
  82. /* | | | | | | |.+-+ */
  83. /* | | | | | | | | | */
  84. /* |.|.|.|.|.|.|.|.+----- */
  85. /* | | | | | | | | | */
  86. /* The 1-2-1-4-1-2-1-8-... pattern is the position of the last 1 bit in */
  87. /* the binary expansion of the current column. Each Schur update is */
  88. /* applied as soon as the necessary portion of U is available. */
  89. /* [1] Toledo, S. 1997. Locality of Reference in LU Decomposition with */
  90. /* Partial Pivoting. SIAM J. Matrix Anal. Appl. 18, 4 (Oct. 1997), */
  91. /* 1065-1081. http://dx.doi.org/10.1137/S0895479896297744 */
  92. /* Arguments */
  93. /* ========= */
  94. /* M (input) INTEGER */
  95. /* The number of rows of the matrix A. M >= 0. */
  96. /* N (input) INTEGER */
  97. /* The number of columns of the matrix A. N >= 0. */
  98. /* A (input/output) DOUBLE PRECISION array, dimension (LDA,N) */
  99. /* On entry, the M-by-N matrix to be factored. */
  100. /* On exit, the factors L and U from the factorization */
  101. /* A = P*L*U; the unit diagonal elements of L are not stored. */
  102. /* LDA (input) INTEGER */
  103. /* The leading dimension of the array A. LDA >= max(1,M). */
  104. /* IPIV (output) INTEGER array, dimension (min(M,N)) */
  105. /* The pivot indices; for 1 <= i <= min(M,N), row i of the */
  106. /* matrix was interchanged with row IPIV(i). */
  107. /* INFO (output) INTEGER */
  108. /* = 0: successful exit */
  109. /* < 0: if INFO = -i, the i-th argument had an illegal value */
  110. /* > 0: if INFO = i, U(i,i) is exactly zero. The factorization */
  111. /* has been completed, but the factor U is exactly */
  112. /* singular, and division by zero will occur if it is used */
  113. /* to solve a system of equations. */
  114. /* ===================================================================== */
  115. /* .. Parameters .. */
  116. /* .. */
  117. /* .. Local Scalars .. */
  118. /* .. */
  119. /* .. External Functions .. */
  120. /* .. */
  121. /* .. External Subroutines .. */
  122. /* .. */
  123. /* .. Intrinsic Functions .. */
  124. /* .. */
  125. /* .. Executable Statements .. */
  126. /* Test the input parameters. */
  127. /* Parameter adjustments */
  128. a_dim1 = *lda;
  129. a_offset = 1 + a_dim1;
  130. a -= a_offset;
  131. --ipiv;
  132. /* Function Body */
  133. *info = 0;
  134. if (*m < 0) {
  135. *info = -1;
  136. } else if (*n < 0) {
  137. *info = -2;
  138. } else if (*lda < max(1,*m)) {
  139. *info = -4;
  140. }
  141. if (*info != 0) {
  142. i__1 = -(*info);
  143. _starpu_xerbla_("DGETRF", &i__1);
  144. return 0;
  145. }
  146. /* Quick return if possible */
  147. if (*m == 0 || *n == 0) {
  148. return 0;
  149. }
  150. /* Compute machine safe minimum */
  151. sfmin = _starpu_dlamch_("S");
  152. nstep = min(*m,*n);
  153. i__1 = nstep;
  154. for (j = 1; j <= i__1; ++j) {
  155. kahead = j & -j;
  156. kstart = j + 1 - kahead;
  157. /* Computing MIN */
  158. i__2 = kahead, i__3 = *m - j;
  159. kcols = min(i__2,i__3);
  160. /* Find pivot. */
  161. i__2 = *m - j + 1;
  162. jp = j - 1 + _starpu_idamax_(&i__2, &a[j + j * a_dim1], &c__1);
  163. ipiv[j] = jp;
  164. /* Permute just this column. */
  165. if (jp != j) {
  166. tmp = a[j + j * a_dim1];
  167. a[j + j * a_dim1] = a[jp + j * a_dim1];
  168. a[jp + j * a_dim1] = tmp;
  169. }
  170. /* Apply pending permutations to L */
  171. ntopiv = 1;
  172. ipivstart = j;
  173. jpivstart = j - ntopiv;
  174. while(ntopiv < kahead) {
  175. _starpu_dlaswp_(&ntopiv, &a[jpivstart * a_dim1 + 1], lda, &ipivstart, &j,
  176. &ipiv[1], &c__1);
  177. ipivstart -= ntopiv;
  178. ntopiv <<= 1;
  179. jpivstart -= ntopiv;
  180. }
  181. /* Permute U block to match L */
  182. _starpu_dlaswp_(&kcols, &a[(j + 1) * a_dim1 + 1], lda, &kstart, &j, &ipiv[1],
  183. &c__1);
  184. /* Factor the current column */
  185. if (a[j + j * a_dim1] != 0. && ! _starpu_disnan_(&a[j + j * a_dim1])) {
  186. if ((d__1 = a[j + j * a_dim1], abs(d__1)) >= sfmin) {
  187. i__2 = *m - j;
  188. d__1 = 1. / a[j + j * a_dim1];
  189. _starpu_dscal_(&i__2, &d__1, &a[j + 1 + j * a_dim1], &c__1);
  190. } else {
  191. i__2 = *m - j;
  192. for (i__ = 1; i__ <= i__2; ++i__) {
  193. a[j + i__ + j * a_dim1] /= a[j + j * a_dim1];
  194. }
  195. }
  196. } else if (a[j + j * a_dim1] == 0. && *info == 0) {
  197. *info = j;
  198. }
  199. /* Solve for U block. */
  200. _starpu_dtrsm_("Left", "Lower", "No transpose", "Unit", &kahead, &kcols, &
  201. c_b12, &a[kstart + kstart * a_dim1], lda, &a[kstart + (j + 1)
  202. * a_dim1], lda);
  203. /* Schur complement. */
  204. i__2 = *m - j;
  205. _starpu_dgemm_("No transpose", "No transpose", &i__2, &kcols, &kahead, &c_b15,
  206. &a[j + 1 + kstart * a_dim1], lda, &a[kstart + (j + 1) *
  207. a_dim1], lda, &c_b12, &a[j + 1 + (j + 1) * a_dim1], lda);
  208. }
  209. /* Handle pivot permutations on the way out of the recursion */
  210. npived = nstep & -nstep;
  211. j = nstep - npived;
  212. while(j > 0) {
  213. ntopiv = j & -j;
  214. i__1 = j + 1;
  215. _starpu_dlaswp_(&ntopiv, &a[(j - ntopiv + 1) * a_dim1 + 1], lda, &i__1, &
  216. nstep, &ipiv[1], &c__1);
  217. j -= ntopiv;
  218. }
  219. /* If short and wide, handle the rest of the columns. */
  220. if (*m < *n) {
  221. i__1 = *n - *m;
  222. _starpu_dlaswp_(&i__1, &a[(*m + kcols + 1) * a_dim1 + 1], lda, &c__1, m, &
  223. ipiv[1], &c__1);
  224. i__1 = *n - *m;
  225. _starpu_dtrsm_("Left", "Lower", "No transpose", "Unit", m, &i__1, &c_b12, &a[
  226. a_offset], lda, &a[(*m + kcols + 1) * a_dim1 + 1], lda);
  227. }
  228. return 0;
  229. /* End of DGETRF */
  230. } /* _starpu_dgetrf_ */