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

◆ ztzrzf()

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

ZTZRZF

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

Purpose:
 ZTZRZF 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*16 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*16 array, dimension (M)
          The scalar factors of the elementary reflectors.
[out]WORK
          WORK is COMPLEX*16 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 ZTZRZF,
  and tau(k) is the kth element of the array TAU.

Definition at line 150 of file ztzrzf.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*16 A( LDA, * ), TAU( * ), WORK( * )
161* ..
162*
163* =====================================================================
164*
165* .. Parameters ..
166 COMPLEX*16 ZERO
167 parameter( zero = ( 0.0d+0, 0.0d+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, zlarzb, zlarzt, zlatrz
176* ..
177* .. Intrinsic Functions ..
178 INTRINSIC max, min
179* ..
180* .. External Functions ..
181 INTEGER ILAENV
182 EXTERNAL ilaenv
183* ..
184* .. Executable Statements ..
185*
186* Test the input arguments
187*
188 info = 0
189 lquery = ( lwork.EQ.-1 )
190 IF( m.LT.0 ) THEN
191 info = -1
192 ELSE IF( n.LT.m ) THEN
193 info = -2
194 ELSE IF( lda.LT.max( 1, m ) ) THEN
195 info = -4
196 END IF
197*
198 IF( info.EQ.0 ) THEN
199 IF( m.EQ.0 .OR. m.EQ.n ) THEN
200 lwkopt = 1
201 lwkmin = 1
202 ELSE
203*
204* Determine the block size.
205*
206 nb = ilaenv( 1, 'ZGERQF', ' ', m, n, -1, -1 )
207 lwkopt = m*nb
208 lwkmin = max( 1, m )
209 END IF
210 work( 1 ) = lwkopt
211*
212 IF( lwork.LT.lwkmin .AND. .NOT.lquery ) THEN
213 info = -7
214 END IF
215 END IF
216*
217 IF( info.NE.0 ) THEN
218 CALL xerbla( 'ZTZRZF', -info )
219 RETURN
220 ELSE IF( lquery ) THEN
221 RETURN
222 END IF
223*
224* Quick return if possible
225*
226 IF( m.EQ.0 ) THEN
227 RETURN
228 ELSE IF( m.EQ.n ) THEN
229 DO 10 i = 1, n
230 tau( i ) = zero
231 10 CONTINUE
232 RETURN
233 END IF
234*
235 nbmin = 2
236 nx = 1
237 iws = m
238 IF( nb.GT.1 .AND. nb.LT.m ) THEN
239*
240* Determine when to cross over from blocked to unblocked code.
241*
242 nx = max( 0, ilaenv( 3, 'ZGERQF', ' ', m, n, -1, -1 ) )
243 IF( nx.LT.m ) THEN
244*
245* Determine if workspace is large enough for blocked code.
246*
247 ldwork = m
248 iws = ldwork*nb
249 IF( lwork.LT.iws ) THEN
250*
251* Not enough workspace to use optimal NB: reduce NB and
252* determine the minimum value of NB.
253*
254 nb = lwork / ldwork
255 nbmin = max( 2, ilaenv( 2, 'ZGERQF', ' ', m, n, -1,
256 $ -1 ) )
257 END IF
258 END IF
259 END IF
260*
261 IF( nb.GE.nbmin .AND. nb.LT.m .AND. nx.LT.m ) THEN
262*
263* Use blocked code initially.
264* The last kk rows are handled by the block method.
265*
266 m1 = min( m+1, n )
267 ki = ( ( m-nx-1 ) / nb )*nb
268 kk = min( m, ki+nb )
269*
270 DO 20 i = m - kk + ki + 1, m - kk + 1, -nb
271 ib = min( m-i+1, nb )
272*
273* Compute the TZ factorization of the current block
274* A(i:i+ib-1,i:n)
275*
276 CALL zlatrz( ib, n-i+1, n-m, a( i, i ), lda, tau( i ),
277 $ work )
278 IF( i.GT.1 ) THEN
279*
280* Form the triangular factor of the block reflector
281* H = H(i+ib-1) . . . H(i+1) H(i)
282*
283 CALL zlarzt( 'Backward', 'Rowwise', n-m, ib, a( i, m1 ),
284 $ lda, tau( i ), work, ldwork )
285*
286* Apply H to A(1:i-1,i:n) from the right
287*
288 CALL zlarzb( 'Right', 'No transpose', 'Backward',
289 $ 'Rowwise', i-1, n-i+1, ib, n-m, a( i, m1 ),
290 $ lda, work, ldwork, a( 1, i ), lda,
291 $ work( ib+1 ), ldwork )
292 END IF
293 20 CONTINUE
294 mu = i + nb - 1
295 ELSE
296 mu = m
297 END IF
298*
299* Use unblocked code to factor the last or only block
300*
301 IF( mu.GT.0 )
302 $ CALL zlatrz( mu, n, n-m, a, lda, tau, work )
303*
304 work( 1 ) = lwkopt
305*
306 RETURN
307*
308* End of ZTZRZF
309*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine zlarzb(side, trans, direct, storev, m, n, k, l, v, ldv, t, ldt, c, ldc, work, ldwork)
ZLARZB applies a block reflector or its conjugate-transpose to a general matrix.
Definition zlarzb.f:183
subroutine zlarzt(direct, storev, n, k, v, ldv, tau, t, ldt)
ZLARZT forms the triangular factor T of a block reflector H = I - vtvH.
Definition zlarzt.f:185
subroutine zlatrz(m, n, l, a, lda, tau, work)
ZLATRZ factors an upper trapezoidal matrix by means of unitary transformations.
Definition zlatrz.f:140
Here is the call graph for this function:
Here is the caller graph for this function: