LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ sgetrf()

subroutine sgetrf ( integer  M,
integer  N,
real, dimension( lda, * )  A,
integer  LDA,
integer, dimension( * )  IPIV,
integer  INFO 
)

SGETRF VARIANT: Crout Level 3 BLAS version of the algorithm.

SGETRF VARIANT: iterative version of Sivan Toledo's recursive LU algorithm

SGETRF VARIANT: left-looking Level 3 BLAS version of the algorithm.

Purpose:

 SGETRF computes an LU factorization of a general M-by-N matrix A
 using partial pivoting with row interchanges.

 The factorization has the form
    A = P * L * U
 where P is a permutation matrix, L is lower triangular with unit
 diagonal elements (lower trapezoidal if m > n), and U is upper
 triangular (upper trapezoidal if m < n).

 This is the Crout Level 3 BLAS version of the algorithm.
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix A.  N >= 0.
[in,out]A
          A is REAL array, dimension (LDA,N)
          On entry, the M-by-N matrix to be factored.
          On exit, the factors L and U from the factorization
          A = P*L*U; the unit diagonal elements of L are not stored.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[out]IPIV
          IPIV is INTEGER array, dimension (min(M,N))
          The pivot indices; for 1 <= i <= min(M,N), row i of the
          matrix was interchanged with row IPIV(i).
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, U(i,i) is exactly zero. The factorization
                has been completed, but the factor U is exactly
                singular, and division by zero will occur if it is used
                to solve a system of equations.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
December 2016

Purpose:

 SGETRF computes an LU factorization of a general M-by-N matrix A
 using partial pivoting with row interchanges.

 The factorization has the form
    A = P * L * U
 where P is a permutation matrix, L is lower triangular with unit
 diagonal elements (lower trapezoidal if m > n), and U is upper
 triangular (upper trapezoidal if m < n).

 This is the left-looking Level 3 BLAS version of the algorithm.
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix A.  N >= 0.
[in,out]A
          A is REAL array, dimension (LDA,N)
          On entry, the M-by-N matrix to be factored.
          On exit, the factors L and U from the factorization
          A = P*L*U; the unit diagonal elements of L are not stored.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[out]IPIV
          IPIV is INTEGER array, dimension (min(M,N))
          The pivot indices; for 1 <= i <= min(M,N), row i of the
          matrix was interchanged with row IPIV(i).
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, U(i,i) is exactly zero. The factorization
                has been completed, but the factor U is exactly
                singular, and division by zero will occur if it is used
                to solve a system of equations.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
December 2016

Purpose:

 SGETRF computes an LU factorization of a general M-by-N matrix A
 using partial pivoting with row interchanges.

 The factorization has the form
    A = P * L * U
 where P is a permutation matrix, L is lower triangular with unit
 diagonal elements (lower trapezoidal if m > n), and U is upper
 triangular (upper trapezoidal if m < n).

 This code implements an iterative version of Sivan Toledo's recursive
 LU algorithm[1].  For square matrices, this iterative versions should
 be within a factor of two of the optimum number of memory transfers.

 The pattern is as follows, with the large blocks of U being updated
 in one call to STRSM, and the dotted lines denoting sections that
 have had all pending permutations applied:

  1 2 3 4 5 6 7 8
 +-+-+---+-------+------
 | |1|   |       |
 |.+-+ 2 |       |
 | | |   |       |
 |.|.+-+-+   4   |
 | | | |1|       |
 | | |.+-+       |
 | | | | |       |
 |.|.|.|.+-+-+---+  8
 | | | | | |1|   |
 | | | | |.+-+ 2 |
 | | | | | | |   |
 | | | | |.|.+-+-+
 | | | | | | | |1|
 | | | | | | |.+-+
 | | | | | | | | |
 |.|.|.|.|.|.|.|.+-----
 | | | | | | | | |

 The 1-2-1-4-1-2-1-8-... pattern is the position of the last 1 bit in
 the binary expansion of the current column.  Each Schur update is
 applied as soon as the necessary portion of U is available.

 [1] Toledo, S. 1997. Locality of Reference in LU Decomposition with
 Partial Pivoting. SIAM J. Matrix Anal. Appl. 18, 4 (Oct. 1997),
 1065-1081. http://dx.doi.org/10.1137/S0895479896297744
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix A.  N >= 0.
[in,out]A
          A is REAL array, dimension (LDA,N)
          On entry, the M-by-N matrix to be factored.
          On exit, the factors L and U from the factorization
          A = P*L*U; the unit diagonal elements of L are not stored.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[out]IPIV
          IPIV is INTEGER array, dimension (min(M,N))
          The pivot indices; for 1 <= i <= min(M,N), row i of the
          matrix was interchanged with row IPIV(i).
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, U(i,i) is exactly zero. The factorization
                has been completed, but the factor U is exactly
                singular, and division by zero will occur if it is used
                to solve a system of equations.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
December 2016

Definition at line 101 of file sgetrf.f.

102*
103* -- LAPACK computational routine --
104* -- LAPACK is a software package provided by Univ. of Tennessee, --
105* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
106*
107* .. Scalar Arguments ..
108 INTEGER INFO, LDA, M, N
109* ..
110* .. Array Arguments ..
111 INTEGER IPIV( * )
112 REAL A( LDA, * )
113* ..
114*
115* =====================================================================
116*
117* .. Parameters ..
118 REAL ONE
119 parameter( one = 1.0e+0 )
120* ..
121* .. Local Scalars ..
122 INTEGER I, IINFO, J, JB, NB
123* ..
124* .. External Subroutines ..
125 EXTERNAL sgemm, sgetf2, slaswp, strsm, xerbla
126* ..
127* .. External Functions ..
128 INTEGER ILAENV
129 EXTERNAL ilaenv
130* ..
131* .. Intrinsic Functions ..
132 INTRINSIC max, min
133* ..
134* .. Executable Statements ..
135*
136* Test the input parameters.
137*
138 info = 0
139 IF( m.LT.0 ) THEN
140 info = -1
141 ELSE IF( n.LT.0 ) THEN
142 info = -2
143 ELSE IF( lda.LT.max( 1, m ) ) THEN
144 info = -4
145 END IF
146 IF( info.NE.0 ) THEN
147 CALL xerbla( 'SGETRF', -info )
148 RETURN
149 END IF
150*
151* Quick return if possible
152*
153 IF( m.EQ.0 .OR. n.EQ.0 )
154 $ RETURN
155*
156* Determine the block size for this environment.
157*
158 nb = ilaenv( 1, 'SGETRF', ' ', m, n, -1, -1 )
159 IF( nb.LE.1 .OR. nb.GE.min( m, n ) ) THEN
160*
161* Use unblocked code.
162*
163 CALL sgetf2( m, n, a, lda, ipiv, info )
164 ELSE
165*
166* Use blocked code.
167*
168 DO 20 j = 1, min( m, n ), nb
169 jb = min( min( m, n )-j+1, nb )
170*
171* Update current block.
172*
173 CALL sgemm( 'No transpose', 'No transpose',
174 $ m-j+1, jb, j-1, -one,
175 $ a( j, 1 ), lda, a( 1, j ), lda, one,
176 $ a( j, j ), lda )
177
178*
179* Factor diagonal and subdiagonal blocks and test for exact
180* singularity.
181*
182 CALL sgetf2( m-j+1, jb, a( j, j ), lda, ipiv( j ), iinfo )
183*
184* Adjust INFO and the pivot indices.
185*
186 IF( info.EQ.0 .AND. iinfo.GT.0 )
187 $ info = iinfo + j - 1
188 DO 10 i = j, min( m, j+jb-1 )
189 ipiv( i ) = j - 1 + ipiv( i )
190 10 CONTINUE
191*
192* Apply interchanges to column 1:J-1
193*
194 CALL slaswp( j-1, a, lda, j, j+jb-1, ipiv, 1 )
195*
196 IF ( j+jb.LE.n ) THEN
197*
198* Apply interchanges to column J+JB:N
199*
200 CALL slaswp( n-j-jb+1, a( 1, j+jb ), lda, j, j+jb-1,
201 $ ipiv, 1 )
202*
203 CALL sgemm( 'No transpose', 'No transpose',
204 $ jb, n-j-jb+1, j-1, -one,
205 $ a( j, 1 ), lda, a( 1, j+jb ), lda, one,
206 $ a( j, j+jb ), lda )
207*
208* Compute block row of U.
209*
210 CALL strsm( 'Left', 'Lower', 'No transpose', 'Unit',
211 $ jb, n-j-jb+1, one, a( j, j ), lda,
212 $ a( j, j+jb ), lda )
213 END IF
214
215 20 CONTINUE
216
217 END IF
218 RETURN
219*
220* End of SGETRF
221*
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: ilaenv.f:162
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine sgetf2(M, N, A, LDA, IPIV, INFO)
SGETF2 computes the LU factorization of a general m-by-n matrix using partial pivoting with row inter...
Definition: sgetf2.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 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
Here is the call graph for this function: