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