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

◆ ztzrqf()

subroutine ztzrqf ( integer  m,
integer  n,
complex*16, dimension( lda, * )  a,
integer  lda,
complex*16, dimension( * )  tau,
integer  info 
)

ZTZRQF

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

Purpose:
 This routine is deprecated and has been replaced by routine ZTZRZF.

 ZTZRQF 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]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.
Further Details:
  The  factorization is obtained by Householder's method.  The kth
  transformation matrix, Z( k ), whose conjugate transpose is used to
  introduce zeros into the (m - k + 1)th row of A, is given in the form

     Z( k ) = ( I     0   ),
              ( 0  T( k ) )

  where

     T( k ) = I - tau*u( k )*u( k )**H,   u( k ) = (   1    ),
                                                   (   0    )
                                                   ( z( k ) )

  tau is a scalar and z( k ) is an ( n - m ) element vector.
  tau and z( k ) are chosen to annihilate the elements of the kth row
  of X.

  The scalar tau is returned in the kth element of TAU and the vector
  u( k ) in the kth row of A, such that the elements of z( k ) are
  in  a( k, m + 1 ), ..., a( k, n ). The elements of R are returned in
  the upper triangular part of A.

  Z is given by

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

Definition at line 137 of file ztzrqf.f.

138*
139* -- LAPACK computational routine --
140* -- LAPACK is a software package provided by Univ. of Tennessee, --
141* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
142*
143* .. Scalar Arguments ..
144 INTEGER INFO, LDA, M, N
145* ..
146* .. Array Arguments ..
147 COMPLEX*16 A( LDA, * ), TAU( * )
148* ..
149*
150* =====================================================================
151*
152* .. Parameters ..
153 COMPLEX*16 CONE, CZERO
154 parameter( cone = ( 1.0d+0, 0.0d+0 ),
155 $ czero = ( 0.0d+0, 0.0d+0 ) )
156* ..
157* .. Local Scalars ..
158 INTEGER I, K, M1
159 COMPLEX*16 ALPHA
160* ..
161* .. Intrinsic Functions ..
162 INTRINSIC dconjg, max, min
163* ..
164* .. External Subroutines ..
165 EXTERNAL xerbla, zaxpy, zcopy, zgemv, zgerc, zlacgv,
166 $ zlarfg
167* ..
168* .. Executable Statements ..
169*
170* Test the input parameters.
171*
172 info = 0
173 IF( m.LT.0 ) THEN
174 info = -1
175 ELSE IF( n.LT.m ) THEN
176 info = -2
177 ELSE IF( lda.LT.max( 1, m ) ) THEN
178 info = -4
179 END IF
180 IF( info.NE.0 ) THEN
181 CALL xerbla( 'ZTZRQF', -info )
182 RETURN
183 END IF
184*
185* Perform the factorization.
186*
187 IF( m.EQ.0 )
188 $ RETURN
189 IF( m.EQ.n ) THEN
190 DO 10 i = 1, n
191 tau( i ) = czero
192 10 CONTINUE
193 ELSE
194 m1 = min( m+1, n )
195 DO 20 k = m, 1, -1
196*
197* Use a Householder reflection to zero the kth row of A.
198* First set up the reflection.
199*
200 a( k, k ) = dconjg( a( k, k ) )
201 CALL zlacgv( n-m, a( k, m1 ), lda )
202 alpha = a( k, k )
203 CALL zlarfg( n-m+1, alpha, a( k, m1 ), lda, tau( k ) )
204 a( k, k ) = alpha
205 tau( k ) = dconjg( tau( k ) )
206*
207 IF( tau( k ).NE.czero .AND. k.GT.1 ) THEN
208*
209* We now perform the operation A := A*P( k )**H.
210*
211* Use the first ( k - 1 ) elements of TAU to store a( k ),
212* where a( k ) consists of the first ( k - 1 ) elements of
213* the kth column of A. Also let B denote the first
214* ( k - 1 ) rows of the last ( n - m ) columns of A.
215*
216 CALL zcopy( k-1, a( 1, k ), 1, tau, 1 )
217*
218* Form w = a( k ) + B*z( k ) in TAU.
219*
220 CALL zgemv( 'No transpose', k-1, n-m, cone, a( 1, m1 ),
221 $ lda, a( k, m1 ), lda, cone, tau, 1 )
222*
223* Now form a( k ) := a( k ) - conjg(tau)*w
224* and B := B - conjg(tau)*w*z( k )**H.
225*
226 CALL zaxpy( k-1, -dconjg( tau( k ) ), tau, 1, a( 1, k ),
227 $ 1 )
228 CALL zgerc( k-1, n-m, -dconjg( tau( k ) ), tau, 1,
229 $ a( k, m1 ), lda, a( 1, m1 ), lda )
230 END IF
231 20 CONTINUE
232 END IF
233*
234 RETURN
235*
236* End of ZTZRQF
237*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zaxpy(n, za, zx, incx, zy, incy)
ZAXPY
Definition zaxpy.f:88
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
ZGEMV
Definition zgemv.f:160
subroutine zgerc(m, n, alpha, x, incx, y, incy, a, lda)
ZGERC
Definition zgerc.f:130
subroutine zlacgv(n, x, incx)
ZLACGV conjugates a complex vector.
Definition zlacgv.f:74
subroutine zlarfg(n, alpha, x, incx, tau)
ZLARFG generates an elementary reflector (Householder matrix).
Definition zlarfg.f:106
Here is the call graph for this function:
Here is the caller graph for this function: