LAPACK  3.10.0
LAPACK: Linear Algebra PACKage

◆ 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 >= 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(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 144 of file dgeqrf.f.

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