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

◆ sgeqrf()

subroutine sgeqrf ( integer  m,
integer  n,
real, dimension( lda, * )  a,
integer  lda,
real, dimension( * )  tau,
real, dimension( * )  work,
integer  lwork,
integer  info 
)

SGEQRF

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

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