LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
subroutine ztplqt2 ( integer  M,
integer  N,
integer  L,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( ldb, * )  B,
integer  LDB,
complex*16, dimension( ldt, * )  T,
integer  LDT,
integer  INFO 
)

ZTPLQT2 computes a LQ factorization of a real or complex "triangular-pentagonal" matrix, which is composed of a triangular block and a pentagonal block, using the compact WY representation for Q.

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

Purpose:
 ZTPLQT2 computes a LQ a factorization of a complex "triangular-pentagonal"
 matrix C, which is composed of a triangular block A and pentagonal block B,
 using the compact WY representation for Q.
Parameters
[in]M
          M is INTEGER
          The total number of rows of the matrix B.
          M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix B, and the order of
          the triangular matrix A.
          N >= 0.
[in]L
          L is INTEGER
          The number of rows of the lower trapezoidal part of B.
          MIN(M,N) >= L >= 0.  See Further Details.
[in,out]A
          A is COMPLEX*16 array, dimension (LDA,N)
          On entry, the lower triangular M-by-M matrix A.
          On exit, the elements on and below the diagonal of the array
          contain the lower triangular matrix L.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in,out]B
          B is COMPLEX*16 array, dimension (LDB,N)
          On entry, the pentagonal M-by-N matrix B.  The first N-L columns
          are rectangular, and the last L columns are lower trapezoidal.
          On exit, B contains the pentagonal matrix V.  See Further Details.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,M).
[out]T
          T is COMPLEX*16 array, dimension (LDT,M)
          The N-by-N upper triangular factor T of the block reflector.
          See Further Details.
[in]LDT
          LDT is INTEGER
          The leading dimension of the array T.  LDT >= max(1,M)
[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.
Date
December 2016
Further Details:
  The input matrix C is a M-by-(M+N) matrix

               C = [ A ][ B ]


  where A is an lower triangular N-by-N matrix, and B is M-by-N pentagonal
  matrix consisting of a M-by-(N-L) rectangular matrix B1 left of a M-by-L
  upper trapezoidal matrix B2:

               B = [ B1 ][ B2 ]
                   [ B1 ]  <-     M-by-(N-L) rectangular
                   [ B2 ]  <-     M-by-L lower trapezoidal.

  The lower trapezoidal matrix B2 consists of the first L columns of a
  N-by-N lower triangular matrix, where 0 <= L <= MIN(M,N).  If L=0,
  B is rectangular M-by-N; if M=L=N, B is lower triangular.

  The matrix W stores the elementary reflectors H(i) in the i-th row
  above the diagonal (of A) in the M-by-(M+N) input matrix C

               C = [ A ][ B ]
                   [ A ]  <- lower triangular N-by-N
                   [ B ]  <- M-by-N pentagonal

  so that W can be represented as

               W = [ I ][ V ]
                   [ I ]  <- identity, N-by-N
                   [ V ]  <- M-by-N, same form as B.

  Thus, all of information needed for W is contained on exit in B, which
  we call V above.  Note that V has the same form as B; that is,

               W = [ V1 ][ V2 ]
                   [ V1 ] <-     M-by-(N-L) rectangular
                   [ V2 ] <-     M-by-L lower trapezoidal.

  The rows of V represent the vectors which define the H(i)'s.
  The (M+N)-by-(M+N) block reflector H is then given by

               H = I - W**T * T * W

  where W^H is the conjugate transpose of W and T is the upper triangular
  factor of the block reflector.

Definition at line 179 of file ztplqt2.f.

179 *
180 * -- LAPACK computational routine (version 3.7.0) --
181 * -- LAPACK is a software package provided by Univ. of Tennessee, --
182 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
183 * December 2016
184 *
185 * .. Scalar Arguments ..
186  INTEGER info, lda, ldb, ldt, n, m, l
187 * ..
188 * .. Array Arguments ..
189  COMPLEX*16 a( lda, * ), b( ldb, * ), t( ldt, * )
190 * ..
191 *
192 * =====================================================================
193 *
194 * .. Parameters ..
195  COMPLEX*16 one, zero
196  parameter( zero = ( 0.0d+0, 0.0d+0 ),one = ( 1.0d+0, 0.0d+0 ) )
197 * ..
198 * .. Local Scalars ..
199  INTEGER i, j, p, mp, np
200  COMPLEX*16 alpha
201 * ..
202 * .. External Subroutines ..
203  EXTERNAL zlarfg, zgemv, zgerc, ztrmv, xerbla
204 * ..
205 * .. Intrinsic Functions ..
206  INTRINSIC max, min
207 * ..
208 * .. Executable Statements ..
209 *
210 * Test the input arguments
211 *
212  info = 0
213  IF( m.LT.0 ) THEN
214  info = -1
215  ELSE IF( n.LT.0 ) THEN
216  info = -2
217  ELSE IF( l.LT.0 .OR. l.GT.min(m,n) ) THEN
218  info = -3
219  ELSE IF( lda.LT.max( 1, m ) ) THEN
220  info = -5
221  ELSE IF( ldb.LT.max( 1, m ) ) THEN
222  info = -7
223  ELSE IF( ldt.LT.max( 1, m ) ) THEN
224  info = -9
225  END IF
226  IF( info.NE.0 ) THEN
227  CALL xerbla( 'ZTPLQT2', -info )
228  RETURN
229  END IF
230 *
231 * Quick return if possible
232 *
233  IF( n.EQ.0 .OR. m.EQ.0 ) RETURN
234 *
235  DO i = 1, m
236 *
237 * Generate elementary reflector H(I) to annihilate B(I,:)
238 *
239  p = n-l+min( l, i )
240  CALL zlarfg( p+1, a( i, i ), b( i, 1 ), ldb, t( 1, i ) )
241  t(1,i)=conjg(t(1,i))
242  IF( i.LT.m ) THEN
243  DO j = 1, p
244  b( i, j ) = conjg(b(i,j))
245  END DO
246 *
247 * W(M-I:1) := C(I+1:M,I:N) * C(I,I:N) [use W = T(M,:)]
248 *
249  DO j = 1, m-i
250  t( m, j ) = (a( i+j, i ))
251  END DO
252  CALL zgemv( 'N', m-i, p, one, b( i+1, 1 ), ldb,
253  $ b( i, 1 ), ldb, one, t( m, 1 ), ldt )
254 *
255 * C(I+1:M,I:N) = C(I+1:M,I:N) + alpha * C(I,I:N)*W(M-1:1)^H
256 *
257  alpha = -(t( 1, i ))
258  DO j = 1, m-i
259  a( i+j, i ) = a( i+j, i ) + alpha*(t( m, j ))
260  END DO
261  CALL zgerc( m-i, p, (alpha), t( m, 1 ), ldt,
262  $ b( i, 1 ), ldb, b( i+1, 1 ), ldb )
263  DO j = 1, p
264  b( i, j ) = conjg(b(i,j))
265  END DO
266  END IF
267  END DO
268 *
269  DO i = 2, m
270 *
271 * T(I,1:I-1) := C(I:I-1,1:N)**H * (alpha * C(I,I:N))
272 *
273  alpha = -(t( 1, i ))
274  DO j = 1, i-1
275  t( i, j ) = zero
276  END DO
277  p = min( i-1, l )
278  np = min( n-l+1, n )
279  mp = min( p+1, m )
280  DO j = 1, n-l+p
281  b(i,j)=conjg(b(i,j))
282  END DO
283 *
284 * Triangular part of B2
285 *
286  DO j = 1, p
287  t( i, j ) = (alpha*b( i, n-l+j ))
288  END DO
289  CALL ztrmv( 'L', 'N', 'N', p, b( 1, np ), ldb,
290  $ t( i, 1 ), ldt )
291 *
292 * Rectangular part of B2
293 *
294  CALL zgemv( 'N', i-1-p, l, alpha, b( mp, np ), ldb,
295  $ b( i, np ), ldb, zero, t( i,mp ), ldt )
296 *
297 * B1
298 
299 *
300  CALL zgemv( 'N', i-1, n-l, alpha, b, ldb, b( i, 1 ), ldb,
301  $ one, t( i, 1 ), ldt )
302 *
303 
304 *
305 * T(1:I-1,I) := T(1:I-1,1:I-1) * T(I,1:I-1)
306 *
307  DO j = 1, i-1
308  t(i,j)=conjg(t(i,j))
309  END DO
310  CALL ztrmv( 'L', 'C', 'N', i-1, t, ldt, t( i, 1 ), ldt )
311  DO j = 1, i-1
312  t(i,j)=conjg(t(i,j))
313  END DO
314  DO j = 1, n-l+p
315  b(i,j)=conjg(b(i,j))
316  END DO
317 *
318 * T(I,I) = tau(I)
319 *
320  t( i, i ) = t( 1, i )
321  t( 1, i ) = zero
322  END DO
323  DO i=1,m
324  DO j= i+1,m
325  t(i,j)=(t(j,i))
326  t(j,i)=zero
327  END DO
328  END DO
329 
330 *
331 * End of ZTPLQT2
332 *
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:132
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ztrmv(UPLO, TRANS, DIAG, N, A, LDA, X, INCX)
ZTRMV
Definition: ztrmv.f:149
subroutine zlarfg(N, ALPHA, X, INCX, TAU)
ZLARFG generates an elementary reflector (Householder matrix).
Definition: zlarfg.f:108

Here is the call graph for this function:

Here is the caller graph for this function: