c@427: /* dgetrf.f -- translated by f2c (version 20061008). c@427: You must link the resulting object file with libf2c: c@427: on Microsoft Windows system, link with libf2c.lib; c@427: on Linux or Unix systems, link with .../path/to/libf2c.a -lm c@427: or, if you install libf2c.a in a standard place, with -lf2c -lm c@427: -- in that order, at the end of the command line, as in c@427: cc *.o -lf2c -lm c@427: Source for libf2c is in /netlib/f2c/libf2c.zip, e.g., c@427: c@427: http://www.netlib.org/f2c/libf2c.zip c@427: */ c@427: c@427: #include "f2c.h" c@427: #include "blaswrap.h" c@427: c@427: /* Table of constant values */ c@427: c@427: static integer c__1 = 1; c@427: static integer c_n1 = -1; c@427: static doublereal c_b16 = 1.; c@427: static doublereal c_b19 = -1.; c@427: c@427: /* Subroutine */ int dgetrf_(integer *m, integer *n, doublereal *a, integer * c@427: lda, integer *ipiv, integer *info) c@427: { c@427: /* System generated locals */ c@427: integer a_dim1, a_offset, i__1, i__2, i__3, i__4, i__5; c@427: c@427: /* Local variables */ c@427: integer i__, j, jb, nb; c@427: extern /* Subroutine */ int dgemm_(char *, char *, integer *, integer *, c@427: integer *, doublereal *, doublereal *, integer *, doublereal *, c@427: integer *, doublereal *, doublereal *, integer *); c@427: integer iinfo; c@427: extern /* Subroutine */ int dtrsm_(char *, char *, char *, char *, c@427: integer *, integer *, doublereal *, doublereal *, integer *, c@427: doublereal *, integer *), dgetf2_( c@427: integer *, integer *, doublereal *, integer *, integer *, integer c@427: *), xerbla_(char *, integer *); c@427: extern integer ilaenv_(integer *, char *, char *, integer *, integer *, c@427: integer *, integer *); c@427: extern /* Subroutine */ int dlaswp_(integer *, doublereal *, integer *, c@427: integer *, integer *, integer *, integer *); c@427: c@427: c@427: /* -- LAPACK routine (version 3.2) -- */ c@427: /* Univ. of Tennessee, Univ. of California Berkeley and NAG Ltd.. */ c@427: /* November 2006 */ c@427: c@427: /* .. Scalar Arguments .. */ c@427: /* .. */ c@427: /* .. Array Arguments .. */ c@427: /* .. */ c@427: c@427: /* Purpose */ c@427: /* ======= */ c@427: c@427: /* DGETRF computes an LU factorization of a general M-by-N matrix A */ c@427: /* using partial pivoting with row interchanges. */ c@427: c@427: /* The factorization has the form */ c@427: /* A = P * L * U */ c@427: /* where P is a permutation matrix, L is lower triangular with unit */ c@427: /* diagonal elements (lower trapezoidal if m > n), and U is upper */ c@427: /* triangular (upper trapezoidal if m < n). */ c@427: c@427: /* This is the right-looking Level 3 BLAS version of the algorithm. */ c@427: c@427: /* Arguments */ c@427: /* ========= */ c@427: c@427: /* M (input) INTEGER */ c@427: /* The number of rows of the matrix A. M >= 0. */ c@427: c@427: /* N (input) INTEGER */ c@427: /* The number of columns of the matrix A. N >= 0. */ c@427: c@427: /* A (input/output) DOUBLE PRECISION array, dimension (LDA,N) */ c@427: /* On entry, the M-by-N matrix to be factored. */ c@427: /* On exit, the factors L and U from the factorization */ c@427: /* A = P*L*U; the unit diagonal elements of L are not stored. */ c@427: c@427: /* LDA (input) INTEGER */ c@427: /* The leading dimension of the array A. LDA >= max(1,M). */ c@427: c@427: /* IPIV (output) INTEGER array, dimension (min(M,N)) */ c@427: /* The pivot indices; for 1 <= i <= min(M,N), row i of the */ c@427: /* matrix was interchanged with row IPIV(i). */ c@427: c@427: /* INFO (output) INTEGER */ c@427: /* = 0: successful exit */ c@427: /* < 0: if INFO = -i, the i-th argument had an illegal value */ c@427: /* > 0: if INFO = i, U(i,i) is exactly zero. The factorization */ c@427: /* has been completed, but the factor U is exactly */ c@427: /* singular, and division by zero will occur if it is used */ c@427: /* to solve a system of equations. */ c@427: c@427: /* ===================================================================== */ c@427: c@427: /* .. Parameters .. */ c@427: /* .. */ c@427: /* .. Local Scalars .. */ c@427: /* .. */ c@427: /* .. External Subroutines .. */ c@427: /* .. */ c@427: /* .. External Functions .. */ c@427: /* .. */ c@427: /* .. Intrinsic Functions .. */ c@427: /* .. */ c@427: /* .. Executable Statements .. */ c@427: c@427: /* Test the input parameters. */ c@427: c@427: /* Parameter adjustments */ c@427: a_dim1 = *lda; c@427: a_offset = 1 + a_dim1; c@427: a -= a_offset; c@427: --ipiv; c@427: c@427: /* Function Body */ c@427: *info = 0; c@427: if (*m < 0) { c@427: *info = -1; c@427: } else if (*n < 0) { c@427: *info = -2; c@427: } else if (*lda < max(1,*m)) { c@427: *info = -4; c@427: } c@427: if (*info != 0) { c@427: i__1 = -(*info); c@427: xerbla_("DGETRF", &i__1); c@427: return 0; c@427: } c@427: c@427: /* Quick return if possible */ c@427: c@427: if (*m == 0 || *n == 0) { c@427: return 0; c@427: } c@427: c@427: /* Determine the block size for this environment. */ c@427: c@427: nb = ilaenv_(&c__1, "DGETRF", " ", m, n, &c_n1, &c_n1); c@427: if (nb <= 1 || nb >= min(*m,*n)) { c@427: c@427: /* Use unblocked code. */ c@427: c@427: dgetf2_(m, n, &a[a_offset], lda, &ipiv[1], info); c@427: } else { c@427: c@427: /* Use blocked code. */ c@427: c@427: i__1 = min(*m,*n); c@427: i__2 = nb; c@427: for (j = 1; i__2 < 0 ? j >= i__1 : j <= i__1; j += i__2) { c@427: /* Computing MIN */ c@427: i__3 = min(*m,*n) - j + 1; c@427: jb = min(i__3,nb); c@427: c@427: /* Factor diagonal and subdiagonal blocks and test for exact */ c@427: /* singularity. */ c@427: c@427: i__3 = *m - j + 1; c@427: dgetf2_(&i__3, &jb, &a[j + j * a_dim1], lda, &ipiv[j], &iinfo); c@427: c@427: /* Adjust INFO and the pivot indices. */ c@427: c@427: if (*info == 0 && iinfo > 0) { c@427: *info = iinfo + j - 1; c@427: } c@427: /* Computing MIN */ c@427: i__4 = *m, i__5 = j + jb - 1; c@427: i__3 = min(i__4,i__5); c@427: for (i__ = j; i__ <= i__3; ++i__) { c@427: ipiv[i__] = j - 1 + ipiv[i__]; c@427: /* L10: */ c@427: } c@427: c@427: /* Apply interchanges to columns 1:J-1. */ c@427: c@427: i__3 = j - 1; c@427: i__4 = j + jb - 1; c@427: dlaswp_(&i__3, &a[a_offset], lda, &j, &i__4, &ipiv[1], &c__1); c@427: c@427: if (j + jb <= *n) { c@427: c@427: /* Apply interchanges to columns J+JB:N. */ c@427: c@427: i__3 = *n - j - jb + 1; c@427: i__4 = j + jb - 1; c@427: dlaswp_(&i__3, &a[(j + jb) * a_dim1 + 1], lda, &j, &i__4, & c@427: ipiv[1], &c__1); c@427: c@427: /* Compute block row of U. */ c@427: c@427: i__3 = *n - j - jb + 1; c@427: dtrsm_("Left", "Lower", "No transpose", "Unit", &jb, &i__3, & c@427: c_b16, &a[j + j * a_dim1], lda, &a[j + (j + jb) * c@427: a_dim1], lda); c@427: if (j + jb <= *m) { c@427: c@427: /* Update trailing submatrix. */ c@427: c@427: i__3 = *m - j - jb + 1; c@427: i__4 = *n - j - jb + 1; c@427: dgemm_("No transpose", "No transpose", &i__3, &i__4, &jb, c@427: &c_b19, &a[j + jb + j * a_dim1], lda, &a[j + (j + c@427: jb) * a_dim1], lda, &c_b16, &a[j + jb + (j + jb) * c@427: a_dim1], lda); c@427: } c@427: } c@427: /* L20: */ c@427: } c@427: } c@427: return 0; c@427: c@427: /* End of DGETRF */ c@427: c@427: } /* dgetrf_ */