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

◆ cunglq()

subroutine cunglq ( integer  m,
integer  n,
integer  k,
complex, dimension( lda, * )  a,
integer  lda,
complex, dimension( * )  tau,
complex, dimension( * )  work,
integer  lwork,
integer  info 
)

CUNGLQ

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

Purpose:
 CUNGLQ generates an M-by-N complex matrix Q with orthonormal rows,
 which is defined as the first M rows of a product of K elementary
 reflectors of order N

       Q  =  H(k)**H . . . H(2)**H H(1)**H

 as returned by CGELQF.
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrix Q. M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix Q. N >= M.
[in]K
          K is INTEGER
          The number of elementary reflectors whose product defines the
          matrix Q. M >= K >= 0.
[in,out]A
          A is COMPLEX array, dimension (LDA,N)
          On entry, the i-th row must contain the vector which defines
          the elementary reflector H(i), for i = 1,2,...,k, as returned
          by CGELQF in the first k rows of its array argument A.
          On exit, the M-by-N matrix Q.
[in]LDA
          LDA is INTEGER
          The first dimension of the array A. LDA >= max(1,M).
[in]TAU
          TAU is COMPLEX array, dimension (K)
          TAU(i) must contain the scalar factor of the elementary
          reflector H(i), as returned by CGELQF.
[out]WORK
          WORK is COMPLEX 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 has an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 126 of file cunglq.f.

127*
128* -- LAPACK computational routine --
129* -- LAPACK is a software package provided by Univ. of Tennessee, --
130* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
131*
132* .. Scalar Arguments ..
133 INTEGER INFO, K, LDA, LWORK, M, N
134* ..
135* .. Array Arguments ..
136 COMPLEX A( LDA, * ), TAU( * ), WORK( * )
137* ..
138*
139* =====================================================================
140*
141* .. Parameters ..
142 COMPLEX ZERO
143 parameter( zero = ( 0.0e+0, 0.0e+0 ) )
144* ..
145* .. Local Scalars ..
146 LOGICAL LQUERY
147 INTEGER I, IB, IINFO, IWS, J, KI, KK, L, LDWORK,
148 $ LWKOPT, NB, NBMIN, NX
149* ..
150* .. External Subroutines ..
151 EXTERNAL clarfb, clarft, cungl2, xerbla
152* ..
153* .. Intrinsic Functions ..
154 INTRINSIC max, min
155* ..
156* .. External Functions ..
157 INTEGER ILAENV
158 REAL SROUNDUP_LWORK
159 EXTERNAL ilaenv, sroundup_lwork
160* ..
161* .. Executable Statements ..
162*
163* Test the input arguments
164*
165 info = 0
166 nb = ilaenv( 1, 'CUNGLQ', ' ', m, n, k, -1 )
167 lwkopt = max( 1, m )*nb
168 work( 1 ) = sroundup_lwork(lwkopt)
169 lquery = ( lwork.EQ.-1 )
170 IF( m.LT.0 ) THEN
171 info = -1
172 ELSE IF( n.LT.m ) THEN
173 info = -2
174 ELSE IF( k.LT.0 .OR. k.GT.m ) THEN
175 info = -3
176 ELSE IF( lda.LT.max( 1, m ) ) THEN
177 info = -5
178 ELSE IF( lwork.LT.max( 1, m ) .AND. .NOT.lquery ) THEN
179 info = -8
180 END IF
181 IF( info.NE.0 ) THEN
182 CALL xerbla( 'CUNGLQ', -info )
183 RETURN
184 ELSE IF( lquery ) THEN
185 RETURN
186 END IF
187*
188* Quick return if possible
189*
190 IF( m.LE.0 ) THEN
191 work( 1 ) = 1
192 RETURN
193 END IF
194*
195 nbmin = 2
196 nx = 0
197 iws = m
198 IF( nb.GT.1 .AND. nb.LT.k ) THEN
199*
200* Determine when to cross over from blocked to unblocked code.
201*
202 nx = max( 0, ilaenv( 3, 'CUNGLQ', ' ', m, n, k, -1 ) )
203 IF( nx.LT.k ) THEN
204*
205* Determine if workspace is large enough for blocked code.
206*
207 ldwork = m
208 iws = ldwork*nb
209 IF( lwork.LT.iws ) THEN
210*
211* Not enough workspace to use optimal NB: reduce NB and
212* determine the minimum value of NB.
213*
214 nb = lwork / ldwork
215 nbmin = max( 2, ilaenv( 2, 'CUNGLQ', ' ', m, n, k, -1 ) )
216 END IF
217 END IF
218 END IF
219*
220 IF( nb.GE.nbmin .AND. nb.LT.k .AND. nx.LT.k ) THEN
221*
222* Use blocked code after the last block.
223* The first kk rows are handled by the block method.
224*
225 ki = ( ( k-nx-1 ) / nb )*nb
226 kk = min( k, ki+nb )
227*
228* Set A(kk+1:m,1:kk) to zero.
229*
230 DO 20 j = 1, kk
231 DO 10 i = kk + 1, m
232 a( i, j ) = zero
233 10 CONTINUE
234 20 CONTINUE
235 ELSE
236 kk = 0
237 END IF
238*
239* Use unblocked code for the last or only block.
240*
241 IF( kk.LT.m )
242 $ CALL cungl2( m-kk, n-kk, k-kk, a( kk+1, kk+1 ), lda,
243 $ tau( kk+1 ), work, iinfo )
244*
245 IF( kk.GT.0 ) THEN
246*
247* Use blocked code
248*
249 DO 50 i = ki + 1, 1, -nb
250 ib = min( nb, k-i+1 )
251 IF( i+ib.LE.m ) THEN
252*
253* Form the triangular factor of the block reflector
254* H = H(i) H(i+1) . . . H(i+ib-1)
255*
256 CALL clarft( 'Forward', 'Rowwise', n-i+1, ib, a( i, i ),
257 $ lda, tau( i ), work, ldwork )
258*
259* Apply H**H to A(i+ib:m,i:n) from the right
260*
261 CALL clarfb( 'Right', 'Conjugate transpose', 'Forward',
262 $ 'Rowwise', m-i-ib+1, n-i+1, ib, a( i, i ),
263 $ lda, work, ldwork, a( i+ib, i ), lda,
264 $ work( ib+1 ), ldwork )
265 END IF
266*
267* Apply H**H to columns i:n of current block
268*
269 CALL cungl2( ib, n-i+1, ib, a( i, i ), lda, tau( i ), work,
270 $ iinfo )
271*
272* Set columns 1:i-1 of current block to zero
273*
274 DO 40 j = 1, i - 1
275 DO 30 l = i, i + ib - 1
276 a( l, j ) = zero
277 30 CONTINUE
278 40 CONTINUE
279 50 CONTINUE
280 END IF
281*
282 work( 1 ) = sroundup_lwork(iws)
283 RETURN
284*
285* End of CUNGLQ
286*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine clarfb(side, trans, direct, storev, m, n, k, v, ldv, t, ldt, c, ldc, work, ldwork)
CLARFB applies a block reflector or its conjugate-transpose to a general rectangular matrix.
Definition clarfb.f:197
subroutine clarft(direct, storev, n, k, v, ldv, tau, t, ldt)
CLARFT forms the triangular factor T of a block reflector H = I - vtvH
Definition clarft.f:163
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine cungl2(m, n, k, a, lda, tau, work, info)
CUNGL2 generates all or part of the unitary matrix Q from an LQ factorization determined by cgelqf (u...
Definition cungl2.f:113
Here is the call graph for this function:
Here is the caller graph for this function: