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

◆ cgeqlf()

subroutine cgeqlf ( integer  m,
integer  n,
complex, dimension( lda, * )  a,
integer  lda,
complex, dimension( * )  tau,
complex, dimension( * )  work,
integer  lwork,
integer  info 
)

CGEQLF

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

Purpose:
 CGEQLF computes a QL factorization of a complex M-by-N matrix A:
 A = Q * L.
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 COMPLEX array, dimension (LDA,N)
          On entry, the M-by-N matrix A.
          On exit,
          if m >= n, the lower triangle of the subarray
          A(m-n+1:m,1:n) contains the N-by-N lower triangular matrix L;
          if m <= n, the elements on and below the (n-m)-th
          superdiagonal contain the M-by-N lower trapezoidal matrix L;
          the remaining elements, with the array TAU, represent the
          unitary matrix Q as a product of elementary reflectors
          (see Further Details).
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[out]TAU
          TAU is COMPLEX array, dimension (min(M,N))
          The scalar factors of the elementary reflectors (see Further
          Details).
[out]WORK
          WORK is COMPLEX array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.  LWORK >= max(1,N).
          For optimum performance LWORK >= N*NB, where NB is
          the optimal blocksize.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[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.
Further Details:
  The matrix Q is represented as a product of elementary reflectors

     Q = H(k) . . . H(2) H(1), where k = min(m,n).

  Each H(i) has the form

     H(i) = I - tau * v * v**H

  where tau is a complex scalar, and v is a complex vector with
  v(m-k+i+1:m) = 0 and v(m-k+i) = 1; v(1:m-k+i-1) is stored on exit in
  A(1:m-k+i-1,n-k+i), and tau in TAU(i).

Definition at line 137 of file cgeqlf.f.

138*
139* -- LAPACK computational routine --
140* -- LAPACK is a software package provided by Univ. of Tennessee, --
141* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
142*
143* .. Scalar Arguments ..
144 INTEGER INFO, LDA, LWORK, M, N
145* ..
146* .. Array Arguments ..
147 COMPLEX A( LDA, * ), TAU( * ), WORK( * )
148* ..
149*
150* =====================================================================
151*
152* .. Local Scalars ..
153 LOGICAL LQUERY
154 INTEGER I, IB, IINFO, IWS, K, KI, KK, LDWORK, LWKOPT,
155 $ MU, NB, NBMIN, NU, NX
156* ..
157* .. External Subroutines ..
158 EXTERNAL cgeql2, clarfb, clarft, xerbla
159* ..
160* .. Intrinsic Functions ..
161 INTRINSIC max, min
162* ..
163* .. External Functions ..
164 INTEGER ILAENV
165 REAL SROUNDUP_LWORK
166 EXTERNAL ilaenv, sroundup_lwork
167* ..
168* .. Executable Statements ..
169*
170* Test the input arguments
171*
172 info = 0
173 lquery = ( lwork.EQ.-1 )
174 IF( m.LT.0 ) THEN
175 info = -1
176 ELSE IF( n.LT.0 ) THEN
177 info = -2
178 ELSE IF( lda.LT.max( 1, m ) ) THEN
179 info = -4
180 END IF
181*
182 IF( info.EQ.0 ) THEN
183 k = min( m, n )
184 IF( k.EQ.0 ) THEN
185 lwkopt = 1
186 ELSE
187 nb = ilaenv( 1, 'CGEQLF', ' ', m, n, -1, -1 )
188 lwkopt = n*nb
189 END IF
190 work( 1 ) = sroundup_lwork(lwkopt)
191*
192 IF( lwork.LT.max( 1, n ) .AND. .NOT.lquery ) THEN
193 info = -7
194 END IF
195 END IF
196*
197 IF( info.NE.0 ) THEN
198 CALL xerbla( 'CGEQLF', -info )
199 RETURN
200 ELSE IF( lquery ) THEN
201 RETURN
202 END IF
203*
204* Quick return if possible
205*
206 IF( k.EQ.0 ) THEN
207 RETURN
208 END IF
209*
210 nbmin = 2
211 nx = 1
212 iws = n
213 IF( nb.GT.1 .AND. nb.LT.k ) THEN
214*
215* Determine when to cross over from blocked to unblocked code.
216*
217 nx = max( 0, ilaenv( 3, 'CGEQLF', ' ', m, n, -1, -1 ) )
218 IF( nx.LT.k ) THEN
219*
220* Determine if workspace is large enough for blocked code.
221*
222 ldwork = n
223 iws = ldwork*nb
224 IF( lwork.LT.iws ) THEN
225*
226* Not enough workspace to use optimal NB: reduce NB and
227* determine the minimum value of NB.
228*
229 nb = lwork / ldwork
230 nbmin = max( 2, ilaenv( 2, 'CGEQLF', ' ', m, n, -1,
231 $ -1 ) )
232 END IF
233 END IF
234 END IF
235*
236 IF( nb.GE.nbmin .AND. nb.LT.k .AND. nx.LT.k ) THEN
237*
238* Use blocked code initially.
239* The last kk columns are handled by the block method.
240*
241 ki = ( ( k-nx-1 ) / nb )*nb
242 kk = min( k, ki+nb )
243*
244 DO 10 i = k - kk + ki + 1, k - kk + 1, -nb
245 ib = min( k-i+1, nb )
246*
247* Compute the QL factorization of the current block
248* A(1:m-k+i+ib-1,n-k+i:n-k+i+ib-1)
249*
250 CALL cgeql2( m-k+i+ib-1, ib, a( 1, n-k+i ), lda, tau( i ),
251 $ work, iinfo )
252 IF( n-k+i.GT.1 ) THEN
253*
254* Form the triangular factor of the block reflector
255* H = H(i+ib-1) . . . H(i+1) H(i)
256*
257 CALL clarft( 'Backward', 'Columnwise', m-k+i+ib-1, ib,
258 $ a( 1, n-k+i ), lda, tau( i ), work, ldwork )
259*
260* Apply H**H to A(1:m-k+i+ib-1,1:n-k+i-1) from the left
261*
262 CALL clarfb( 'Left', 'Conjugate transpose', 'Backward',
263 $ 'Columnwise', m-k+i+ib-1, n-k+i-1, ib,
264 $ a( 1, n-k+i ), lda, work, ldwork, a, lda,
265 $ work( ib+1 ), ldwork )
266 END IF
267 10 CONTINUE
268 mu = m - k + i + nb - 1
269 nu = n - k + i + nb - 1
270 ELSE
271 mu = m
272 nu = n
273 END IF
274*
275* Use unblocked code to factor the last or only block
276*
277 IF( mu.GT.0 .AND. nu.GT.0 )
278 $ CALL cgeql2( mu, nu, a, lda, tau, work, iinfo )
279*
280 work( 1 ) = sroundup_lwork(iws)
281 RETURN
282*
283* End of CGEQLF
284*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgeql2(m, n, a, lda, tau, work, info)
CGEQL2 computes the QL factorization of a general rectangular matrix using an unblocked algorithm.
Definition cgeql2.f:123
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
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.
Definition clarfb.f:197
subroutine clarft(direct, storev, n, k, v, ldv, tau, t, ldt)
CLARFT forms the triangular factor T of a block reflector H = I - vtvH
Definition clarft.f:163
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
Here is the call graph for this function:
Here is the caller graph for this function: