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

◆ dgelqf()

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

DGELQF

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

Purpose:
 DGELQF computes an LQ factorization of a real M-by-N matrix A:

    A = ( L 0 ) *  Q

 where:

    Q is a N-by-N orthogonal matrix;
    L is a lower-triangular M-by-M matrix;
    0 is a M-by-(N-M) zero matrix, if M < N.
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, the elements on and below the diagonal of the array
          contain the m-by-min(m,n) lower trapezoidal matrix L (L is
          lower triangular if m <= n); the elements above the diagonal,
          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,M).
          For optimum performance LWORK >= M*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(1:i-1) = 0 and v(i) = 1; v(i+1:n) is stored on exit in A(i,i+1:n),
  and tau in TAU(i).

Definition at line 142 of file dgelqf.f.

143*
144* -- LAPACK computational routine --
145* -- LAPACK is a software package provided by Univ. of Tennessee, --
146* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
147*
148* .. Scalar Arguments ..
149 INTEGER INFO, LDA, LWORK, M, N
150* ..
151* .. Array Arguments ..
152 DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * )
153* ..
154*
155* =====================================================================
156*
157* .. Local Scalars ..
158 LOGICAL LQUERY
159 INTEGER I, IB, IINFO, IWS, K, LDWORK, LWKOPT, NB,
160 $ NBMIN, NX
161* ..
162* .. External Subroutines ..
163 EXTERNAL dgelq2, dlarfb, dlarft, xerbla
164* ..
165* .. Intrinsic Functions ..
166 INTRINSIC max, min
167* ..
168* .. External Functions ..
169 INTEGER ILAENV
170 EXTERNAL ilaenv
171* ..
172* .. Executable Statements ..
173*
174* Test the input arguments
175*
176 info = 0
177 nb = ilaenv( 1, 'DGELQF', ' ', m, n, -1, -1 )
178 lwkopt = m*nb
179 work( 1 ) = lwkopt
180 lquery = ( lwork.EQ.-1 )
181 IF( m.LT.0 ) THEN
182 info = -1
183 ELSE IF( n.LT.0 ) THEN
184 info = -2
185 ELSE IF( lda.LT.max( 1, m ) ) THEN
186 info = -4
187 ELSE IF( lwork.LT.max( 1, m ) .AND. .NOT.lquery ) THEN
188 info = -7
189 END IF
190 IF( info.NE.0 ) THEN
191 CALL xerbla( 'DGELQF', -info )
192 RETURN
193 ELSE IF( lquery ) THEN
194 RETURN
195 END IF
196*
197* Quick return if possible
198*
199 k = min( m, n )
200 IF( k.EQ.0 ) THEN
201 work( 1 ) = 1
202 RETURN
203 END IF
204*
205 nbmin = 2
206 nx = 0
207 iws = m
208 IF( nb.GT.1 .AND. nb.LT.k ) THEN
209*
210* Determine when to cross over from blocked to unblocked code.
211*
212 nx = max( 0, ilaenv( 3, 'DGELQF', ' ', m, n, -1, -1 ) )
213 IF( nx.LT.k ) THEN
214*
215* Determine if workspace is large enough for blocked code.
216*
217 ldwork = m
218 iws = ldwork*nb
219 IF( lwork.LT.iws ) THEN
220*
221* Not enough workspace to use optimal NB: reduce NB and
222* determine the minimum value of NB.
223*
224 nb = lwork / ldwork
225 nbmin = max( 2, ilaenv( 2, 'DGELQF', ' ', m, n, -1,
226 $ -1 ) )
227 END IF
228 END IF
229 END IF
230*
231 IF( nb.GE.nbmin .AND. nb.LT.k .AND. nx.LT.k ) THEN
232*
233* Use blocked code initially
234*
235 DO 10 i = 1, k - nx, nb
236 ib = min( k-i+1, nb )
237*
238* Compute the LQ factorization of the current block
239* A(i:i+ib-1,i:n)
240*
241 CALL dgelq2( ib, n-i+1, a( i, i ), lda, tau( i ), work,
242 $ iinfo )
243 IF( i+ib.LE.m ) THEN
244*
245* Form the triangular factor of the block reflector
246* H = H(i) H(i+1) . . . H(i+ib-1)
247*
248 CALL dlarft( 'Forward', 'Rowwise', n-i+1, ib, a( i, i ),
249 $ lda, tau( i ), work, ldwork )
250*
251* Apply H to A(i+ib:m,i:n) from the right
252*
253 CALL dlarfb( 'Right', 'No transpose', 'Forward',
254 $ 'Rowwise', m-i-ib+1, n-i+1, ib, a( i, i ),
255 $ lda, work, ldwork, a( i+ib, i ), lda,
256 $ work( ib+1 ), ldwork )
257 END IF
258 10 CONTINUE
259 ELSE
260 i = 1
261 END IF
262*
263* Use unblocked code to factor the last or only block.
264*
265 IF( i.LE.k )
266 $ CALL dgelq2( m-i+1, n-i+1, a( i, i ), lda, tau( i ), work,
267 $ iinfo )
268*
269 work( 1 ) = iws
270 RETURN
271*
272* End of DGELQF
273*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgelq2(m, n, a, lda, tau, work, info)
DGELQ2 computes the LQ factorization of a general rectangular matrix using an unblocked algorithm.
Definition dgelq2.f:129
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: