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

◆ dgeqrfp()

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

DGEQRFP

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

Purpose:
 DGEQR2P 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 with nonnegative diagonal
    entries;
    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 diagonal entries of R
          are nonnegative; 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).

 See Lapack Working Note 203 for details

Definition at line 148 of file dgeqrfp.f.

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