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

◆ dgeqlf()

subroutine dgeqlf ( integer  m,
integer  n,
double precision, dimension( lda, * )  a,
integer  lda,
double precision, dimension( * )  tau,
double precision, dimension( * )  work,
integer  lwork,
integer  info 
)

DGEQLF

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

Purpose:
 DGEQLF computes a QL factorization of a real 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 DOUBLE PRECISION 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
          orthogonal 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 DOUBLE PRECISION array, dimension (min(M,N))
          The scalar factors of the elementary reflectors (see Further
          Details).
[out]WORK
          WORK is DOUBLE PRECISION 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**T

  where tau is a real scalar, and v is a real 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 dgeqlf.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 DOUBLE PRECISION 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 dgeql2, dlarfb, dlarft, xerbla
159* ..
160* .. Intrinsic Functions ..
161 INTRINSIC max, min
162* ..
163* .. External Functions ..
164 INTEGER ILAENV
165 EXTERNAL ilaenv
166* ..
167* .. Executable Statements ..
168*
169* Test the input arguments
170*
171 info = 0
172 lquery = ( lwork.EQ.-1 )
173 IF( m.LT.0 ) THEN
174 info = -1
175 ELSE IF( n.LT.0 ) THEN
176 info = -2
177 ELSE IF( lda.LT.max( 1, m ) ) THEN
178 info = -4
179 END IF
180*
181 IF( info.EQ.0 ) THEN
182 k = min( m, n )
183 IF( k.EQ.0 ) THEN
184 lwkopt = 1
185 ELSE
186 nb = ilaenv( 1, 'DGEQLF', ' ', m, n, -1, -1 )
187 lwkopt = n*nb
188 END IF
189 work( 1 ) = lwkopt
190*
191 IF( lwork.LT.max( 1, n ) .AND. .NOT.lquery ) THEN
192 info = -7
193 END IF
194 END IF
195*
196 IF( info.NE.0 ) THEN
197 CALL xerbla( 'DGEQLF', -info )
198 RETURN
199 ELSE IF( lquery ) THEN
200 RETURN
201 END IF
202*
203* Quick return if possible
204*
205 IF( k.EQ.0 ) THEN
206 RETURN
207 END IF
208*
209 nbmin = 2
210 nx = 1
211 iws = n
212 IF( nb.GT.1 .AND. nb.LT.k ) THEN
213*
214* Determine when to cross over from blocked to unblocked code.
215*
216 nx = max( 0, ilaenv( 3, 'DGEQLF', ' ', m, n, -1, -1 ) )
217 IF( nx.LT.k ) THEN
218*
219* Determine if workspace is large enough for blocked code.
220*
221 ldwork = n
222 iws = ldwork*nb
223 IF( lwork.LT.iws ) THEN
224*
225* Not enough workspace to use optimal NB: reduce NB and
226* determine the minimum value of NB.
227*
228 nb = lwork / ldwork
229 nbmin = max( 2, ilaenv( 2, 'DGEQLF', ' ', m, n, -1,
230 $ -1 ) )
231 END IF
232 END IF
233 END IF
234*
235 IF( nb.GE.nbmin .AND. nb.LT.k .AND. nx.LT.k ) THEN
236*
237* Use blocked code initially.
238* The last kk columns are handled by the block method.
239*
240 ki = ( ( k-nx-1 ) / nb )*nb
241 kk = min( k, ki+nb )
242*
243 DO 10 i = k - kk + ki + 1, k - kk + 1, -nb
244 ib = min( k-i+1, nb )
245*
246* Compute the QL factorization of the current block
247* A(1:m-k+i+ib-1,n-k+i:n-k+i+ib-1)
248*
249 CALL dgeql2( m-k+i+ib-1, ib, a( 1, n-k+i ), lda, tau( i ),
250 $ work, iinfo )
251 IF( n-k+i.GT.1 ) THEN
252*
253* Form the triangular factor of the block reflector
254* H = H(i+ib-1) . . . H(i+1) H(i)
255*
256 CALL dlarft( 'Backward', 'Columnwise', m-k+i+ib-1, ib,
257 $ a( 1, n-k+i ), lda, tau( i ), work, ldwork )
258*
259* Apply H**T to A(1:m-k+i+ib-1,1:n-k+i-1) from the left
260*
261 CALL dlarfb( 'Left', 'Transpose', 'Backward',
262 $ 'Columnwise', m-k+i+ib-1, n-k+i-1, ib,
263 $ a( 1, n-k+i ), lda, work, ldwork, a, lda,
264 $ work( ib+1 ), ldwork )
265 END IF
266 10 CONTINUE
267 mu = m - k + i + nb - 1
268 nu = n - k + i + nb - 1
269 ELSE
270 mu = m
271 nu = n
272 END IF
273*
274* Use unblocked code to factor the last or only block
275*
276 IF( mu.GT.0 .AND. nu.GT.0 )
277 $ CALL dgeql2( mu, nu, a, lda, tau, work, iinfo )
278*
279 work( 1 ) = iws
280 RETURN
281*
282* End of DGEQLF
283*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgeql2(m, n, a, lda, tau, work, info)
DGEQL2 computes the QL factorization of a general rectangular matrix using an unblocked algorithm.
Definition dgeql2.f:123
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine dlarfb(side, trans, direct, storev, m, n, k, v, ldv, t, ldt, c, ldc, work, ldwork)
DLARFB applies a block reflector or its transpose to a general rectangular matrix.
Definition dlarfb.f:197
subroutine dlarft(direct, storev, n, k, v, ldv, tau, t, ldt)
DLARFT forms the triangular factor T of a block reflector H = I - vtvH
Definition dlarft.f:163
Here is the call graph for this function:
Here is the caller graph for this function: