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

◆ cungrq()

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

CUNGRQ

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

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

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

 as returned by CGERQF.
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 (m-k+i)-th row must contain the vector which
          defines the elementary reflector H(i), for i = 1,2,...,k, as
          returned by CGERQF in the last 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 CGERQF.
[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 127 of file cungrq.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, II, IINFO, IWS, J, KK, L, LDWORK,
149 $ LWKOPT, NB, NBMIN, NX
150* ..
151* .. External Subroutines ..
152 EXTERNAL clarfb, clarft, cungr2, 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.m ) THEN
171 info = -2
172 ELSE IF( k.LT.0 .OR. k.GT.m ) 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( m.LE.0 ) THEN
180 lwkopt = 1
181 ELSE
182 nb = ilaenv( 1, 'CUNGRQ', ' ', m, n, k, -1 )
183 lwkopt = m*nb
184 END IF
185 work( 1 ) = sroundup_lwork(lwkopt)
186*
187 IF( lwork.LT.max( 1, m ) .AND. .NOT.lquery ) THEN
188 info = -8
189 END IF
190 END IF
191*
192 IF( info.NE.0 ) THEN
193 CALL xerbla( 'CUNGRQ', -info )
194 RETURN
195 ELSE IF( lquery ) THEN
196 RETURN
197 END IF
198*
199* Quick return if possible
200*
201 IF( m.LE.0 ) THEN
202 RETURN
203 END IF
204*
205 nbmin = 2
206 nx = 0
207 iws = m
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, 'CUNGRQ', ' ', m, n, k, -1 ) )
213 IF( nx.LT.k ) THEN
214*
215* Determine if workspace is large enough for blocked code.
216*
217 ldwork = m
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, 'CUNGRQ', ' ', 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 rows are handled by the block method.
234*
235 kk = min( k, ( ( k-nx+nb-1 ) / nb )*nb )
236*
237* Set A(1:m-kk,n-kk+1:n) to zero.
238*
239 DO 20 j = n - kk + 1, n
240 DO 10 i = 1, m - kk
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 cungr2( 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 ii = m - k + i
259 IF( ii.GT.1 ) THEN
260*
261* Form the triangular factor of the block reflector
262* H = H(i+ib-1) . . . H(i+1) H(i)
263*
264 CALL clarft( 'Backward', 'Rowwise', n-k+i+ib-1, ib,
265 $ a( ii, 1 ), lda, tau( i ), work, ldwork )
266*
267* Apply H**H to A(1:m-k+i-1,1:n-k+i+ib-1) from the right
268*
269 CALL clarfb( 'Right', 'Conjugate transpose', 'Backward',
270 $ 'Rowwise', ii-1, n-k+i+ib-1, ib, a( ii, 1 ),
271 $ lda, work, ldwork, a, lda, work( ib+1 ),
272 $ ldwork )
273 END IF
274*
275* Apply H**H to columns 1:n-k+i+ib-1 of current block
276*
277 CALL cungr2( ib, n-k+i+ib-1, ib, a( ii, 1 ), lda, tau( i ),
278 $ work, iinfo )
279*
280* Set columns n-k+i+ib:n of current block to zero
281*
282 DO 40 l = n - k + i + ib, n
283 DO 30 j = ii, ii + ib - 1
284 a( j, l ) = zero
285 30 CONTINUE
286 40 CONTINUE
287 50 CONTINUE
288 END IF
289*
290 work( 1 ) = sroundup_lwork(iws)
291 RETURN
292*
293* End of CUNGRQ
294*
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 cungr2(m, n, k, a, lda, tau, work, info)
CUNGR2 generates all or part of the unitary matrix Q from an RQ factorization determined by cgerqf (u...
Definition cungr2.f:114
Here is the call graph for this function:
Here is the caller graph for this function: