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

◆ sgeqrfp()

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

SGEQRFP

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

Purpose:
 SGEQR2P 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 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 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 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 >= 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 sgeqrfp.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 REAL 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 sgeqr2p, slarfb, slarft, xerbla
170* ..
171* .. Intrinsic Functions ..
172 INTRINSIC max, min
173* ..
174* .. External Functions ..
175 INTEGER ILAENV
176 REAL SROUNDUP_LWORK
177 EXTERNAL ilaenv, sroundup_lwork
178* ..
179* .. Executable Statements ..
180*
181* Test the input arguments
182*
183 info = 0
184 nb = ilaenv( 1, 'SGEQRF', ' ', m, n, -1, -1 )
185 lwkopt = n*nb
186 work( 1 ) = sroundup_lwork(lwkopt)
187 lquery = ( lwork.EQ.-1 )
188 IF( m.LT.0 ) THEN
189 info = -1
190 ELSE IF( n.LT.0 ) THEN
191 info = -2
192 ELSE IF( lda.LT.max( 1, m ) ) THEN
193 info = -4
194 ELSE IF( lwork.LT.max( 1, n ) .AND. .NOT.lquery ) THEN
195 info = -7
196 END IF
197 IF( info.NE.0 ) THEN
198 CALL xerbla( 'SGEQRFP', -info )
199 RETURN
200 ELSE IF( lquery ) THEN
201 RETURN
202 END IF
203*
204* Quick return if possible
205*
206 k = min( m, n )
207 IF( k.EQ.0 ) THEN
208 work( 1 ) = 1
209 RETURN
210 END IF
211*
212 nbmin = 2
213 nx = 0
214 iws = n
215 IF( nb.GT.1 .AND. nb.LT.k ) THEN
216*
217* Determine when to cross over from blocked to unblocked code.
218*
219 nx = max( 0, ilaenv( 3, 'SGEQRF', ' ', m, n, -1, -1 ) )
220 IF( nx.LT.k ) THEN
221*
222* Determine if workspace is large enough for blocked code.
223*
224 ldwork = n
225 iws = ldwork*nb
226 IF( lwork.LT.iws ) THEN
227*
228* Not enough workspace to use optimal NB: reduce NB and
229* determine the minimum value of NB.
230*
231 nb = lwork / ldwork
232 nbmin = max( 2, ilaenv( 2, 'SGEQRF', ' ', m, n, -1,
233 $ -1 ) )
234 END IF
235 END IF
236 END IF
237*
238 IF( nb.GE.nbmin .AND. nb.LT.k .AND. nx.LT.k ) THEN
239*
240* Use blocked code initially
241*
242 DO 10 i = 1, k - nx, nb
243 ib = min( k-i+1, nb )
244*
245* Compute the QR factorization of the current block
246* A(i:m,i:i+ib-1)
247*
248 CALL sgeqr2p( m-i+1, ib, a( i, i ), lda, tau( i ), work,
249 $ iinfo )
250 IF( i+ib.LE.n ) THEN
251*
252* Form the triangular factor of the block reflector
253* H = H(i) H(i+1) . . . H(i+ib-1)
254*
255 CALL slarft( 'Forward', 'Columnwise', m-i+1, ib,
256 $ a( i, i ), lda, tau( i ), work, ldwork )
257*
258* Apply H**T to A(i:m,i+ib:n) from the left
259*
260 CALL slarfb( 'Left', 'Transpose', 'Forward',
261 $ 'Columnwise', m-i+1, n-i-ib+1, ib,
262 $ a( i, i ), lda, work, ldwork, a( i, i+ib ),
263 $ lda, work( ib+1 ), ldwork )
264 END IF
265 10 CONTINUE
266 ELSE
267 i = 1
268 END IF
269*
270* Use unblocked code to factor the last or only block.
271*
272 IF( i.LE.k )
273 $ CALL sgeqr2p( m-i+1, n-i+1, a( i, i ), lda, tau( i ), work,
274 $ iinfo )
275*
276 work( 1 ) = sroundup_lwork(iws)
277 RETURN
278*
279* End of SGEQRFP
280*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sgeqr2p(m, n, a, lda, tau, work, info)
SGEQR2P computes the QR factorization of a general rectangular matrix with non-negative diagonal elem...
Definition sgeqr2p.f:134
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: