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

◆ slaorhr_col_getrfnp2()

recursive subroutine slaorhr_col_getrfnp2 ( integer m,
integer n,
real, dimension( lda, * ) a,
integer lda,
real, dimension( * ) d,
integer info )

SLAORHR_COL_GETRFNP2

Download DLAORHR_GETRF2NP + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> SLAORHR_COL_GETRFNP2 computes the modified LU factorization without
!> pivoting of a real general M-by-N matrix A. The factorization has
!> the form:
!>
!>     A - S = L * U,
!>
!> where:
!>    S is a m-by-n diagonal sign matrix with the diagonal D, so that
!>    D(i) = S(i,i), 1 <= i <= min(M,N). The diagonal D is constructed
!>    as D(i)=-SIGN(A(i,i)), where A(i,i) is the value after performing
!>    i-1 steps of Gaussian elimination. This means that the diagonal
!>    element at each step of  Gaussian elimination is at
!>    least one in absolute value (so that division-by-zero not
!>    possible during the division by the diagonal element);
!>
!>    L is a M-by-N lower triangular matrix with unit diagonal elements
!>    (lower trapezoidal if M > N);
!>
!>    and U is a M-by-N upper triangular matrix
!>    (upper trapezoidal if M < N).
!>
!> This routine is an auxiliary routine used in the Householder
!> reconstruction routine SORHR_COL. In SORHR_COL, this routine is
!> applied to an M-by-N matrix A with orthonormal columns, where each
!> element is bounded by one in absolute value. With the choice of
!> the matrix S above, one can show that the diagonal element at each
!> step of Gaussian elimination is the largest (in absolute value) in
!> the column on or below the diagonal, so that no pivoting is required
!> for numerical stability [1].
!>
!> For more details on the Householder reconstruction algorithm,
!> including the modified LU factorization, see [1].
!>
!> This is the recursive version of the LU factorization algorithm.
!> Denote A - S by B. The algorithm divides the matrix B into four
!> submatrices:
!>
!>        [  B11 | B12  ]  where B11 is n1 by n1,
!>    B = [ -----|----- ]        B21 is (m-n1) by n1,
!>        [  B21 | B22  ]        B12 is n1 by n2,
!>                               B22 is (m-n1) by n2,
!>                               with n1 = min(m,n)/2, n2 = n-n1.
!>
!>
!> The subroutine calls itself to factor B11, solves for B21,
!> solves for B12, updates B22, then calls itself to factor B22.
!>
!> For more details on the recursive LU algorithm, see [2].
!>
!> SLAORHR_COL_GETRFNP2 is called to factorize a block by the blocked
!> routine SLAORHR_COL_GETRFNP, which uses blocked code calling
!> Level 3 BLAS to update the submatrix. However, SLAORHR_COL_GETRFNP2
!> is self-sufficient and can be used without SLAORHR_COL_GETRFNP.
!>
!> [1] ,
!>     G. Ballard, J. Demmel, L. Grigori, M. Jacquelin, H.D. Nguyen,
!>     E. Solomonik, J. Parallel Distrib. Comput.,
!>     vol. 85, pp. 3-31, 2015.
!>
!> [2] , F. Gustavson, IBM J. of Res. and Dev.,
!>     vol. 41, no. 6, pp. 737-755, 1997.
!> 
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-S=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]D
!>          D is REAL array, dimension min(M,N)
!>          The diagonal elements of the diagonal M-by-N sign matrix S,
!>          D(i) = S(i,i), where 1 <= i <= min(M,N). The elements can
!>          be only plus or minus one.
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
!>
!> November 2019, Igor Kozachenko,
!>                Computer Science Division,
!>                University of California, Berkeley
!>
!> 

Definition at line 164 of file slaorhr_col_getrfnp2.f.

166 IMPLICIT NONE
167*
168* -- LAPACK computational routine --
169* -- LAPACK is a software package provided by Univ. of Tennessee, --
170* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
171*
172* .. Scalar Arguments ..
173 INTEGER INFO, LDA, M, N
174* ..
175* .. Array Arguments ..
176 REAL A( LDA, * ), D( * )
177* ..
178*
179* =====================================================================
180*
181* .. Parameters ..
182 REAL ONE
183 parameter( one = 1.0e+0 )
184* ..
185* .. Local Scalars ..
186 REAL SFMIN
187 INTEGER I, IINFO, N1, N2
188* ..
189* .. External Functions ..
190 REAL SLAMCH
191 EXTERNAL slamch
192* ..
193* .. External Subroutines ..
194 EXTERNAL sgemm, sscal, strsm, xerbla
195* ..
196* .. Intrinsic Functions ..
197 INTRINSIC abs, sign, max, min
198* ..
199* .. Executable Statements ..
200*
201* Test the input parameters
202*
203 info = 0
204 IF( m.LT.0 ) THEN
205 info = -1
206 ELSE IF( n.LT.0 ) THEN
207 info = -2
208 ELSE IF( lda.LT.max( 1, m ) ) THEN
209 info = -4
210 END IF
211 IF( info.NE.0 ) THEN
212 CALL xerbla( 'SLAORHR_COL_GETRFNP2', -info )
213 RETURN
214 END IF
215*
216* Quick return if possible
217*
218 IF( min( m, n ).EQ.0 )
219 $ RETURN
220
221 IF ( m.EQ.1 ) THEN
222*
223* One row case, (also recursion termination case),
224* use unblocked code
225*
226* Transfer the sign
227*
228 d( 1 ) = -sign( one, a( 1, 1 ) )
229*
230* Construct the row of U
231*
232 a( 1, 1 ) = a( 1, 1 ) - d( 1 )
233*
234 ELSE IF( n.EQ.1 ) THEN
235*
236* One column case, (also recursion termination case),
237* use unblocked code
238*
239* Transfer the sign
240*
241 d( 1 ) = -sign( one, a( 1, 1 ) )
242*
243* Construct the row of U
244*
245 a( 1, 1 ) = a( 1, 1 ) - d( 1 )
246*
247* Scale the elements 2:M of the column
248*
249* Determine machine safe minimum
250*
251 sfmin = slamch('S')
252*
253* Construct the subdiagonal elements of L
254*
255 IF( abs( a( 1, 1 ) ) .GE. sfmin ) THEN
256 CALL sscal( m-1, one / a( 1, 1 ), a( 2, 1 ), 1 )
257 ELSE
258 DO i = 2, m
259 a( i, 1 ) = a( i, 1 ) / a( 1, 1 )
260 END DO
261 END IF
262*
263 ELSE
264*
265* Divide the matrix B into four submatrices
266*
267 n1 = min( m, n ) / 2
268 n2 = n-n1
269
270*
271* Factor B11, recursive call
272*
273 CALL slaorhr_col_getrfnp2( n1, n1, a, lda, d, iinfo )
274*
275* Solve for B21
276*
277 CALL strsm( 'R', 'U', 'N', 'N', m-n1, n1, one, a, lda,
278 $ a( n1+1, 1 ), lda )
279*
280* Solve for B12
281*
282 CALL strsm( 'L', 'L', 'N', 'U', n1, n2, one, a, lda,
283 $ a( 1, n1+1 ), lda )
284*
285* Update B22, i.e. compute the Schur complement
286* B22 := B22 - B21*B12
287*
288 CALL sgemm( 'N', 'N', m-n1, n2, n1, -one, a( n1+1, 1 ), lda,
289 $ a( 1, n1+1 ), lda, one, a( n1+1, n1+1 ), lda )
290*
291* Factor B22, recursive call
292*
293 CALL slaorhr_col_getrfnp2( m-n1, n2, a( n1+1, n1+1 ), lda,
294 $ d( n1+1 ), iinfo )
295*
296 END IF
297 RETURN
298*
299* End of SLAORHR_COL_GETRFNP2
300*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:188
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
recursive subroutine slaorhr_col_getrfnp2(m, n, a, lda, d, info)
SLAORHR_COL_GETRFNP2
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
Here is the call graph for this function:
Here is the caller graph for this function: