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

◆ cungql()

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

CUNGQL

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

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

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

 as returned by CGEQLF.
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. M >= N >= 0.
[in]K
          K is INTEGER
          The number of elementary reflectors whose product defines the
          matrix Q. N >= K >= 0.
[in,out]A
          A is COMPLEX array, dimension (LDA,N)
          On entry, the (n-k+i)-th column must contain the vector which
          defines the elementary reflector H(i), for i = 1,2,...,k, as
          returned by CGEQLF in the last k columns 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 CGEQLF.
[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,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 has an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 127 of file cungql.f.

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