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

◆ sgelqf()

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

SGELQF

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

Purpose:
 SGELQF computes an LQ factorization of a real M-by-N matrix A:

    A = ( L 0 ) *  Q

 where:

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

Definition at line 142 of file sgelqf.f.

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