LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
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*> \ingroup getrf
105*
106* =====================================================================
107 SUBROUTINE dgetrf( M, N, A, LDA, IPIV, INFO )
108*
109* -- LAPACK computational routine --
110* -- LAPACK is a software package provided by Univ. of Tennessee, --
111* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
112*
113* .. Scalar Arguments ..
114 INTEGER INFO, LDA, M, N
115* ..
116* .. Array Arguments ..
117 INTEGER IPIV( * )
118 DOUBLE PRECISION A( LDA, * )
119* ..
120*
121* =====================================================================
122*
123* .. Parameters ..
124 DOUBLE PRECISION ONE
125 parameter( one = 1.0d+0 )
126* ..
127* .. Local Scalars ..
128 INTEGER I, IINFO, J, JB, NB
129* ..
130* .. External Subroutines ..
131 EXTERNAL dgemm, dgetrf2, dlaswp, dtrsm, xerbla
132* ..
133* .. External Functions ..
134 INTEGER ILAENV
135 EXTERNAL ilaenv
136* ..
137* .. Intrinsic Functions ..
138 INTRINSIC max, min
139* ..
140* .. Executable Statements ..
141*
142* Test the input parameters.
143*
144 info = 0
145 IF( m.LT.0 ) THEN
146 info = -1
147 ELSE IF( n.LT.0 ) THEN
148 info = -2
149 ELSE IF( lda.LT.max( 1, m ) ) THEN
150 info = -4
151 END IF
152 IF( info.NE.0 ) THEN
153 CALL xerbla( 'DGETRF', -info )
154 RETURN
155 END IF
156*
157* Quick return if possible
158*
159 IF( m.EQ.0 .OR. n.EQ.0 )
160 $ RETURN
161*
162* Determine the block size for this environment.
163*
164 nb = ilaenv( 1, 'DGETRF', ' ', m, n, -1, -1 )
165 IF( nb.LE.1 .OR. nb.GE.min( m, n ) ) THEN
166*
167* Use unblocked code.
168*
169 CALL dgetrf2( m, n, a, lda, ipiv, info )
170 ELSE
171*
172* Use blocked code.
173*
174 DO 20 j = 1, min( m, n ), nb
175 jb = min( min( m, n )-j+1, nb )
176*
177* Factor diagonal and subdiagonal blocks and test for exact
178* singularity.
179*
180 CALL dgetrf2( m-j+1, jb, a( j, j ), lda, ipiv( j ), iinfo )
181*
182* Adjust INFO and the pivot indices.
183*
184 IF( info.EQ.0 .AND. iinfo.GT.0 )
185 $ info = iinfo + j - 1
186 DO 10 i = j, min( m, j+jb-1 )
187 ipiv( i ) = j - 1 + ipiv( i )
188 10 CONTINUE
189*
190* Apply interchanges to columns 1:J-1.
191*
192 CALL dlaswp( j-1, a, lda, j, j+jb-1, ipiv, 1 )
193*
194 IF( j+jb.LE.n ) THEN
195*
196* Apply interchanges to columns J+JB:N.
197*
198 CALL dlaswp( n-j-jb+1, a( 1, j+jb ), lda, j, j+jb-1,
199 $ ipiv, 1 )
200*
201* Compute block row of U.
202*
203 CALL dtrsm( 'Left', 'Lower', 'No transpose', 'Unit', jb,
204 $ n-j-jb+1, one, a( j, j ), lda, a( j, j+jb ),
205 $ lda )
206 IF( j+jb.LE.m ) THEN
207*
208* Update trailing submatrix.
209*
210 CALL dgemm( 'No transpose', 'No transpose', m-j-jb+1,
211 $ n-j-jb+1, jb, -one, a( j+jb, j ), lda,
212 $ a( j, j+jb ), lda, one, a( j+jb, j+jb ),
213 $ lda )
214 END IF
215 END IF
216 20 CONTINUE
217 END IF
218 RETURN
219*
220* End of DGETRF
221*
222 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
recursive subroutine dgetrf2(m, n, a, lda, ipiv, info)
DGETRF2
Definition dgetrf2.f:113
subroutine dgetrf(m, n, a, lda, ipiv, info)
DGETRF
Definition dgetrf.f:108
subroutine dlaswp(n, a, lda, k1, k2, ipiv, incx)
DLASWP performs a series of row interchanges on a general rectangular matrix.
Definition dlaswp.f:115
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181