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

◆ ctzrzf()

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

CTZRZF

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

Purpose:
 CTZRZF reduces the M-by-N ( M<=N ) complex upper trapezoidal matrix A
 to upper triangular form by means of unitary transformations.

 The upper trapezoidal matrix A is factored as

    A = ( R  0 ) * Z,

 where Z is an N-by-N unitary matrix and R is an M-by-M upper
 triangular matrix.
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 >= M.
[in,out]A
          A is COMPLEX array, dimension (LDA,N)
          On entry, the leading M-by-N upper trapezoidal part of the
          array A must contain the matrix to be factorized.
          On exit, the leading M-by-M upper triangular part of A
          contains the upper triangular matrix R, and elements M+1 to
          N of the first M rows of A, with the array TAU, represent the
          unitary matrix Z as a product of M elementary reflectors.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[out]TAU
          TAU is COMPLEX array, dimension (M)
          The scalar factors of the elementary reflectors.
[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 had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
A. Petitet, Computer Science Dept., Univ. of Tenn., Knoxville, USA
Further Details:
  The N-by-N matrix Z can be computed by

     Z =  Z(1)*Z(2)* ... *Z(M)

  where each N-by-N Z(k) is given by

     Z(k) = I - tau(k)*v(k)*v(k)**H

  with v(k) is the kth row vector of the M-by-N matrix

     V = ( I   A(:,M+1:N) )

  I is the M-by-M identity matrix, A(:,M+1:N)
  is the output stored in A on exit from CTZRZF,
  and tau(k) is the kth element of the array TAU.

Definition at line 150 of file ctzrzf.f.

151*
152* -- LAPACK computational routine --
153* -- LAPACK is a software package provided by Univ. of Tennessee, --
154* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
155*
156* .. Scalar Arguments ..
157 INTEGER INFO, LDA, LWORK, M, N
158* ..
159* .. Array Arguments ..
160 COMPLEX A( LDA, * ), TAU( * ), WORK( * )
161* ..
162*
163* =====================================================================
164*
165* .. Parameters ..
166 COMPLEX ZERO
167 parameter( zero = ( 0.0e+0, 0.0e+0 ) )
168* ..
169* .. Local Scalars ..
170 LOGICAL LQUERY
171 INTEGER I, IB, IWS, KI, KK, LDWORK, LWKMIN, LWKOPT,
172 $ M1, MU, NB, NBMIN, NX
173* ..
174* .. External Subroutines ..
175 EXTERNAL xerbla, clarzb, clarzt, clatrz
176* ..
177* .. Intrinsic Functions ..
178 INTRINSIC max, min
179* ..
180* .. External Functions ..
181 INTEGER ILAENV
182 REAL SROUNDUP_LWORK
183 EXTERNAL ilaenv, sroundup_lwork
184* ..
185* .. Executable Statements ..
186*
187* Test the input arguments
188*
189 info = 0
190 lquery = ( lwork.EQ.-1 )
191 IF( m.LT.0 ) THEN
192 info = -1
193 ELSE IF( n.LT.m ) THEN
194 info = -2
195 ELSE IF( lda.LT.max( 1, m ) ) THEN
196 info = -4
197 END IF
198*
199 IF( info.EQ.0 ) THEN
200 IF( m.EQ.0 .OR. m.EQ.n ) THEN
201 lwkopt = 1
202 lwkmin = 1
203 ELSE
204*
205* Determine the block size.
206*
207 nb = ilaenv( 1, 'CGERQF', ' ', m, n, -1, -1 )
208 lwkopt = m*nb
209 lwkmin = max( 1, m )
210 END IF
211 work( 1 ) = sroundup_lwork(lwkopt)
212*
213 IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
214 info = -7
215 END IF
216 END IF
217*
218 IF( info.NE.0 ) THEN
219 CALL xerbla( 'CTZRZF', -info )
220 RETURN
221 ELSE IF( lquery ) THEN
222 RETURN
223 END IF
224*
225* Quick return if possible
226*
227 IF( m.EQ.0 ) THEN
228 RETURN
229 ELSE IF( m.EQ.n ) THEN
230 DO 10 i = 1, n
231 tau( i ) = zero
232 10 CONTINUE
233 RETURN
234 END IF
235*
236 nbmin = 2
237 nx = 1
238 iws = m
239 IF( nb.GT.1 .AND. nb.LT.m ) THEN
240*
241* Determine when to cross over from blocked to unblocked code.
242*
243 nx = max( 0, ilaenv( 3, 'CGERQF', ' ', m, n, -1, -1 ) )
244 IF( nx.LT.m ) THEN
245*
246* Determine if workspace is large enough for blocked code.
247*
248 ldwork = m
249 iws = ldwork*nb
250 IF( lwork.LT.iws ) THEN
251*
252* Not enough workspace to use optimal NB: reduce NB and
253* determine the minimum value of NB.
254*
255 nb = lwork / ldwork
256 nbmin = max( 2, ilaenv( 2, 'CGERQF', ' ', m, n, -1,
257 $ -1 ) )
258 END IF
259 END IF
260 END IF
261*
262 IF( nb.GE.nbmin .AND. nb.LT.m .AND. nx.LT.m ) THEN
263*
264* Use blocked code initially.
265* The last kk rows are handled by the block method.
266*
267 m1 = min( m+1, n )
268 ki = ( ( m-nx-1 ) / nb )*nb
269 kk = min( m, ki+nb )
270*
271 DO 20 i = m - kk + ki + 1, m - kk + 1, -nb
272 ib = min( m-i+1, nb )
273*
274* Compute the TZ factorization of the current block
275* A(i:i+ib-1,i:n)
276*
277 CALL clatrz( ib, n-i+1, n-m, a( i, i ), lda, tau( i ),
278 $ work )
279 IF( i.GT.1 ) THEN
280*
281* Form the triangular factor of the block reflector
282* H = H(i+ib-1) . . . H(i+1) H(i)
283*
284 CALL clarzt( 'Backward', 'Rowwise', n-m, ib, a( i, m1 ),
285 $ lda, tau( i ), work, ldwork )
286*
287* Apply H to A(1:i-1,i:n) from the right
288*
289 CALL clarzb( 'Right', 'No transpose', 'Backward',
290 $ 'Rowwise', i-1, n-i+1, ib, n-m, a( i, m1 ),
291 $ lda, work, ldwork, a( 1, i ), lda,
292 $ work( ib+1 ), ldwork )
293 END IF
294 20 CONTINUE
295 mu = i + nb - 1
296 ELSE
297 mu = m
298 END IF
299*
300* Use unblocked code to factor the last or only block
301*
302 IF( mu.GT.0 )
303 $ CALL clatrz( mu, n, n-m, a, lda, tau, work )
304*
305 work( 1 ) = sroundup_lwork(lwkopt)
306*
307 RETURN
308*
309* End of CTZRZF
310*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine clarzb(side, trans, direct, storev, m, n, k, l, v, ldv, t, ldt, c, ldc, work, ldwork)
CLARZB applies a block reflector or its conjugate-transpose to a general matrix.
Definition clarzb.f:183
subroutine clarzt(direct, storev, n, k, v, ldv, tau, t, ldt)
CLARZT forms the triangular factor T of a block reflector H = I - vtvH.
Definition clarzt.f:185
subroutine clatrz(m, n, l, a, lda, tau, work)
CLATRZ factors an upper trapezoidal matrix by means of unitary transformations.
Definition clatrz.f:140
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
Here is the call graph for this function:
Here is the caller graph for this function: