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

◆ dgeqrf()

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

DGEQRF

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

Purpose:
 DGEQRF computes a QR factorization of a real M-by-N matrix A:

    A = Q * ( R ),
            ( 0 )

 where:

    Q is a M-by-M orthogonal matrix;
    R is an upper-triangular N-by-N matrix;
    0 is a (M-N)-by-N 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 above the diagonal of the array
          contain the min(M,N)-by-N upper trapezoidal matrix R (R is
          upper triangular if m >= n); the elements below the diagonal,
          with the array TAU, represent the orthogonal matrix Q as a
          product of min(m,n) 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 >= 1, if MIN(M,N) = 0, and LWORK >= N, otherwise.
          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(1) H(2) . . . H(k), 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:m) is stored on exit in A(i+1:m,i),
  and tau in TAU(i).

Definition at line 145 of file dgeqrf.f.

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