136 SUBROUTINE cgelqf( M, N, A, LDA, TAU, WORK, LWORK, INFO )
144 INTEGER INFO, LDA, LWORK, M, N
147 COMPLEX A( lda, * ), TAU( * ), WORK( * )
154 INTEGER I, IB, IINFO, IWS, K, LDWORK, LWKOPT, NB,
172 nb = ilaenv( 1,
'CGELQF',
' ', m, n, -1, -1 )
175 lquery = ( lwork.EQ.-1 )
178 ELSE IF( n.LT.0 )
THEN
180 ELSE IF( lda.LT.max( 1, m ) )
THEN
182 ELSE IF( lwork.LT.max( 1, m ) .AND. .NOT.lquery )
THEN
186 CALL xerbla(
'CGELQF', -info )
188 ELSE IF( lquery )
THEN
203 IF( nb.GT.1 .AND. nb.LT.k )
THEN
207 nx = max( 0, ilaenv( 3,
'CGELQF',
' ', m, n, -1, -1 ) )
214 IF( lwork.LT.iws )
THEN
220 nbmin = max( 2, ilaenv( 2,
'CGELQF',
' ', m, n, -1,
226 IF( nb.GE.nbmin .AND. nb.LT.k .AND. nx.LT.k )
THEN
230 DO 10 i = 1, k - nx, nb
231 ib = min( k-i+1, nb )
236 CALL cgelq2( ib, n-i+1, a( i, i ), lda, tau( i ), work,
243 CALL clarft(
'Forward',
'Rowwise', n-i+1, ib, a( i, i ),
244 $ lda, tau( i ), work, ldwork )
248 CALL clarfb(
'Right',
'No transpose',
'Forward',
249 $
'Rowwise', m-i-ib+1, n-i+1, ib, a( i, i ),
250 $ lda, work, ldwork, a( i+ib, i ), lda,
251 $ work( ib+1 ), ldwork )
261 $
CALL cgelq2( m-i+1, n-i+1, a( i, i ), lda, tau( i ), work,
subroutine clarft(DIRECT, STOREV, N, K, V, LDV, TAU, T, LDT)
CLARFT forms the triangular factor T of a block reflector H = I - vtvH
subroutine xerbla(SRNAME, INFO)
XERBLA
subroutine cgelqf(M, N, A, LDA, TAU, WORK, LWORK, INFO)
CGELQF
subroutine cgelq2(M, N, A, LDA, TAU, WORK, INFO)
CGELQ2 computes the LQ factorization of a general rectangular matrix using an unblocked algorithm...
subroutine clarfb(SIDE, TRANS, DIRECT, STOREV, M, N, K, V, LDV, T, LDT, C, LDC, WORK, LDWORK)
CLARFB applies a block reflector or its conjugate-transpose to a general rectangular matrix...