OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
sgetrf.f
Go to the documentation of this file.
1C> \brief \b SGETRF VARIANT: iterative version of Sivan Toledo's recursive LU algorithm
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE SGETRF( M, N, A, LDA, IPIV, INFO )
12*
13* .. Scalar Arguments ..
14* INTEGER INFO, LDA, M, N
15* ..
16* .. Array Arguments ..
17* INTEGER IPIV( * )
18* REAL A( LDA, * )
19* ..
20*
21* Purpose
22* =======
23*
24C>\details \b Purpose:
25C>\verbatim
26C>
27C> SGETRF computes an LU factorization of a general M-by-N matrix A
28C> using partial pivoting with row interchanges.
29C>
30C> The factorization has the form
31C> A = P * L * U
32C> where P is a permutation matrix, L is lower triangular with unit
33C> diagonal elements (lower trapezoidal if m > n), and U is upper
34C> triangular (upper trapezoidal if m < n).
35C>
36C> This code implements an iterative version of Sivan Toledo's recursive
37C> LU algorithm[1]. For square matrices, this iterative versions should
38C> be within a factor of two of the optimum number of memory transfers.
39C>
40C> The pattern is as follows, with the large blocks of U being updated
41C> in one call to STRSM, and the dotted lines denoting sections that
42C> have had all pending permutations applied:
43C>
44C> 1 2 3 4 5 6 7 8
45C> +-+-+---+-------+------
46C> | |1| | |
47C> |.+-+ 2 | |
48C> | | | | |
49C> |.|.+-+-+ 4 |
50C> | | | |1| |
51C> | | |.+-+ |
52C> | | | | | |
53C> |.|.|.|.+-+-+---+ 8
54C> | | | | | |1| |
55C> | | | | |.+-+ 2 |
56C> | | | | | | | |
57C> | | | | |.|.+-+-+
58C> | | | | | | | |1|
59C> | | | | | | |.+-+
60C> | | | | | | | | |
61C> |.|.|.|.|.|.|.|.+-----
62C> | | | | | | | | |
63C>
64C> The 1-2-1-4-1-2-1-8-... pattern is the position of the last 1 bit in
65C> the binary expansion of the current column. Each Schur update is
66C> applied as soon as the necessary portion of U is available.
67C>
68C> [1] Toledo, S. 1997. Locality of Reference in LU Decomposition with
69C> Partial Pivoting. SIAM J. Matrix Anal. Appl. 18, 4 (Oct. 1997),
70C> 1065-1081. http://dx.doi.org/10.1137/S0895479896297744
71C>
72C>\endverbatim
73*
74* Arguments:
75* ==========
76*
77C> \param[in] M
78C> \verbatim
79C> M is INTEGER
80C> The number of rows of the matrix A. M >= 0.
81C> \endverbatim
82C>
83C> \param[in] N
84C> \verbatim
85C> N is INTEGER
86C> The number of columns of the matrix A. N >= 0.
87C> \endverbatim
88C>
89C> \param[in,out] A
90C> \verbatim
91C> A is REAL array, dimension (LDA,N)
92C> On entry, the M-by-N matrix to be factored.
93C> On exit, the factors L and U from the factorization
94C> A = P*L*U; the unit diagonal elements of L are not stored.
95C> \endverbatim
96C>
97C> \param[in] LDA
98C> \verbatim
99C> LDA is INTEGER
100C> The leading dimension of the array A. LDA >= max(1,M).
101C> \endverbatim
102C>
103C> \param[out] IPIV
104C> \verbatim
105C> IPIV is INTEGER array, dimension (min(M,N))
106C> The pivot indices; for 1 <= i <= min(M,N), row i of the
107C> matrix was interchanged with row IPIV(i).
108C> \endverbatim
109C>
110C> \param[out] INFO
111C> \verbatim
112C> INFO is INTEGER
113C> = 0: successful exit
114C> < 0: if INFO = -i, the i-th argument had an illegal value
115C> > 0: if INFO = i, U(i,i) is exactly zero. The factorization
116C> has been completed, but the factor U is exactly
117C> singular, and division by zero will occur if it is used
118C> to solve a system of equations.
119C> \endverbatim
120C>
121*
122* Authors:
123* ========
124*
125C> \author Univ. of Tennessee
126C> \author Univ. of California Berkeley
127C> \author Univ. of Colorado Denver
128C> \author NAG Ltd.
129*
130C> \date December 2016
131*
132C> \ingroup variantsGEcomputational
133*
134* =====================================================================
135 SUBROUTINE sgetrf( M, N, A, LDA, IPIV, INFO )
136*
137* -- LAPACK computational routine (version 3.X) --
138* -- LAPACK is a software package provided by Univ. of Tennessee, --
139* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
140*
141* .. Scalar Arguments ..
142 INTEGER INFO, LDA, M, N
143* ..
144* .. Array Arguments ..
145 INTEGER IPIV( * )
146 REAL A( LDA, * )
147* ..
148*
149* =====================================================================
150*
151* .. Parameters ..
152 REAL ONE, ZERO, NEGONE
153 parameter( one = 1.0e+0, zero = 0.0e+0 )
154 parameter( negone = -1.0e+0 )
155* ..
156* .. Local Scalars ..
157 REAL SFMIN, TMP
158 INTEGER I, J, JP, NSTEP, NTOPIV, NPIVED, KAHEAD
159 INTEGER KSTART, IPIVSTART, JPIVSTART, KCOLS
160* ..
161* .. External Functions ..
162 REAL SLAMCH
163 INTEGER ISAMAX
164 LOGICAL SISNAN
165 EXTERNAL slamch, isamax, sisnan
166* ..
167* .. External Subroutines ..
168 EXTERNAL strsm, sscal, xerbla, slaswp
169* ..
170* .. Intrinsic Functions ..
171 INTRINSIC max, min, iand
172* ..
173* .. Executable Statements ..
174*
175* Test the input parameters.
176*
177 info = 0
178 IF( m.LT.0 ) THEN
179 info = -1
180 ELSE IF( n.LT.0 ) THEN
181 info = -2
182 ELSE IF( lda.LT.max( 1, m ) ) THEN
183 info = -4
184 END IF
185 IF( info.NE.0 ) THEN
186 CALL xerbla( 'SGETRF', -info )
187 RETURN
188 END IF
189*
190* Quick return if possible
191*
192 IF( m.EQ.0 .OR. n.EQ.0 )
193 $ RETURN
194*
195* Compute machine safe minimum
196*
197 sfmin = slamch( 'S' )
198*
199 nstep = min( m, n )
200 DO j = 1, nstep
201 kahead = iand( j, -j )
202 kstart = j + 1 - kahead
203 kcols = min( kahead, m-j )
204*
205* Find pivot.
206*
207 jp = j - 1 + isamax( m-j+1, a( j, j ), 1 )
208 ipiv( j ) = jp
209
210! Permute just this column.
211 IF (jp .NE. j) THEN
212 tmp = a( j, j )
213 a( j, j ) = a( jp, j )
214 a( jp, j ) = tmp
215 END IF
216
217! Apply pending permutations to L
218 ntopiv = 1
219 ipivstart = j
220 jpivstart = j - ntopiv
221 DO WHILE ( ntopiv .LT. kahead )
222 CALL slaswp( ntopiv, a( 1, jpivstart ), lda, ipivstart, j,
223 $ ipiv, 1 )
224 ipivstart = ipivstart - ntopiv;
225 ntopiv = ntopiv * 2;
226 jpivstart = jpivstart - ntopiv;
227 END DO
228
229! Permute U block to match L
230 CALL slaswp( kcols, a( 1,j+1 ), lda, kstart, j, ipiv, 1 )
231
232! Factor the current column
233 IF( a( j, j ).NE.zero .AND. .NOT.sisnan( a( j, j ) ) ) THEN
234 IF( abs(a( j, j )) .GE. sfmin ) THEN
235 CALL sscal( m-j, one / a( j, j ), a( j+1, j ), 1 )
236 ELSE
237 DO i = 1, m-j
238 a( j+i, j ) = a( j+i, j ) / a( j, j )
239 END DO
240 END IF
241 ELSE IF( a( j,j ) .EQ. zero .AND. info .EQ. 0 ) THEN
242 info = j
243 END IF
244
245! Solve for U block.
246 CALL strsm( 'Left', 'Lower', 'No transpose', 'Unit', kahead,
247 $ kcols, one, a( kstart, kstart ), lda,
248 $ a( kstart, j+1 ), lda )
249! Schur complement.
250 CALL sgemm( 'No transpose', 'No transpose', m-j,
251 $ kcols, kahead, negone, a( j+1, kstart ), lda,
252 $ a( kstart, j+1 ), lda, one, a( j+1, j+1 ), lda )
253 END DO
254
255! Handle pivot permutations on the way out of the recursion
256 npived = iand( nstep, -nstep )
257 j = nstep - npived
258 DO WHILE ( j .GT. 0 )
259 ntopiv = iand( j, -j )
260 CALL slaswp( ntopiv, a( 1, j-ntopiv+1 ), lda, j+1, nstep,
261 $ ipiv, 1 )
262 j = j - ntopiv
263 END DO
264
265! If short and wide, handle the rest of the columns.
266 IF ( m .LT. n ) THEN
267 CALL slaswp( n-m, a( 1, m+kcols+1 ), lda, 1, m, ipiv, 1 )
268 CALL strsm( 'Left', 'Lower', 'No transpose', 'Unit', m,
269 $ n-m, one, a, lda, a( 1,m+kcols+1 ), lda )
270 END IF
271
272 RETURN
273*
274* End of SGETRF
275*
276 END
subroutine xerbla(srname, info)
XERBLA
Definition xerbla.f:60
subroutine sgetrf(m, n, a, lda, ipiv, info)
SGETRF
Definition sgetrf.f:108
subroutine slaswp(n, a, lda, k1, k2, ipiv, incx)
SLASWP performs a series of row interchanges on a general rectangular matrix.
Definition slaswp.f:115
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine strsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRSM
Definition strsm.f:181
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:187
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21