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

◆ ctpqrt2()

subroutine ctpqrt2 ( integer m,
integer n,
integer l,
complex, dimension( lda, * ) a,
integer lda,
complex, dimension( ldb, * ) b,
integer ldb,
complex, dimension( ldt, * ) t,
integer ldt,
integer info )

CTPQRT2 computes a QR 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 CTPQRT2 + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> CTPQRT2 computes a QR factorization of a complex 
!> 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 upper trapezoidal part of B.
!>          MIN(M,N) >= L >= 0.  See Further Details.
!> 
[in,out]A
!>          A is COMPLEX array, dimension (LDA,N)
!>          On entry, the upper triangular N-by-N matrix A.
!>          On exit, the elements on and above the diagonal of the array
!>          contain the upper triangular matrix R.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[in,out]B
!>          B is COMPLEX array, dimension (LDB,N)
!>          On entry, the pentagonal M-by-N matrix B.  The first M-L rows
!>          are rectangular, and the last L rows are upper 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 array, dimension (LDT,N)
!>          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,N)
!> 
[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 input matrix C is a (N+M)-by-N matrix
!>
!>               C = [ A ]
!>                   [ B ]
!>
!>  where A is an upper triangular N-by-N matrix, and B is M-by-N pentagonal
!>  matrix consisting of a (M-L)-by-N rectangular matrix B1 on top of a L-by-N
!>  upper trapezoidal matrix B2:
!>
!>               B = [ B1 ]  <- (M-L)-by-N rectangular
!>                   [ B2 ]  <-     L-by-N upper trapezoidal.
!>
!>  The upper trapezoidal matrix B2 consists of the first L rows of a
!>  N-by-N upper triangular matrix, where 0 <= L <= MIN(M,N).  If L=0,
!>  B is rectangular M-by-N; if M=L=N, B is upper triangular.
!>
!>  The matrix W stores the elementary reflectors H(i) in the i-th column
!>  below the diagonal (of A) in the (N+M)-by-N input matrix C
!>
!>               C = [ A ]  <- upper triangular N-by-N
!>                   [ B ]  <- M-by-N pentagonal
!>
!>  so that W can be represented as
!>
!>               W = [ 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,
!>
!>               V = [ V1 ] <- (M-L)-by-N rectangular
!>                   [ V2 ] <-     L-by-N upper trapezoidal.
!>
!>  The columns 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 * W**H
!>
!>  where W**H is the conjugate transpose of W and T is the upper triangular
!>  factor of the block reflector.
!> 

Definition at line 170 of file ctpqrt2.f.

171*
172* -- LAPACK computational routine --
173* -- LAPACK is a software package provided by Univ. of Tennessee, --
174* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
175*
176* .. Scalar Arguments ..
177 INTEGER INFO, LDA, LDB, LDT, N, M, L
178* ..
179* .. Array Arguments ..
180 COMPLEX A( LDA, * ), B( LDB, * ), T( LDT, * )
181* ..
182*
183* =====================================================================
184*
185* .. Parameters ..
186 COMPLEX ONE, ZERO
187 parameter( one = (1.0,0.0), zero = (0.0,0.0) )
188* ..
189* .. Local Scalars ..
190 INTEGER I, J, P, MP, NP
191 COMPLEX ALPHA
192* ..
193* .. External Subroutines ..
194 EXTERNAL clarfg, cgemv, cgerc, ctrmv, xerbla
195* ..
196* .. Intrinsic Functions ..
197 INTRINSIC max, min
198* ..
199* .. Executable Statements ..
200*
201* Test the input arguments
202*
203 info = 0
204 IF( m.LT.0 ) THEN
205 info = -1
206 ELSE IF( n.LT.0 ) THEN
207 info = -2
208 ELSE IF( l.LT.0 .OR. l.GT.min(m,n) ) THEN
209 info = -3
210 ELSE IF( lda.LT.max( 1, n ) ) THEN
211 info = -5
212 ELSE IF( ldb.LT.max( 1, m ) ) THEN
213 info = -7
214 ELSE IF( ldt.LT.max( 1, n ) ) THEN
215 info = -9
216 END IF
217 IF( info.NE.0 ) THEN
218 CALL xerbla( 'CTPQRT2', -info )
219 RETURN
220 END IF
221*
222* Quick return if possible
223*
224 IF( n.EQ.0 .OR. m.EQ.0 ) RETURN
225*
226 DO i = 1, n
227*
228* Generate elementary reflector H(I) to annihilate B(:,I)
229*
230 p = m-l+min( l, i )
231 CALL clarfg( p+1, a( i, i ), b( 1, i ), 1, t( i, 1 ) )
232 IF( i.LT.n ) THEN
233*
234* W(1:N-I) := C(I:M,I+1:N)**H * C(I:M,I) [use W = T(:,N)]
235*
236 DO j = 1, n-i
237 t( j, n ) = conjg(a( i, i+j ))
238 END DO
239 CALL cgemv( 'C', p, n-i, one, b( 1, i+1 ), ldb,
240 $ b( 1, i ), 1, one, t( 1, n ), 1 )
241*
242* C(I:M,I+1:N) = C(I:m,I+1:N) + alpha*C(I:M,I)*W(1:N-1)**H
243*
244 alpha = -conjg(t( i, 1 ))
245 DO j = 1, n-i
246 a( i, i+j ) = a( i, i+j ) + alpha*conjg(t( j, n ))
247 END DO
248 CALL cgerc( p, n-i, alpha, b( 1, i ), 1,
249 $ t( 1, n ), 1, b( 1, i+1 ), ldb )
250 END IF
251 END DO
252*
253 DO i = 2, n
254*
255* T(1:I-1,I) := C(I:M,1:I-1)**H * (alpha * C(I:M,I))
256*
257 alpha = -t( i, 1 )
258
259 DO j = 1, i-1
260 t( j, i ) = zero
261 END DO
262 p = min( i-1, l )
263 mp = min( m-l+1, m )
264 np = min( p+1, n )
265*
266* Triangular part of B2
267*
268 DO j = 1, p
269 t( j, i ) = alpha*b( m-l+j, i )
270 END DO
271 CALL ctrmv( 'U', 'C', 'N', p, b( mp, 1 ), ldb,
272 $ t( 1, i ), 1 )
273*
274* Rectangular part of B2
275*
276 CALL cgemv( 'C', l, i-1-p, alpha, b( mp, np ), ldb,
277 $ b( mp, i ), 1, zero, t( np, i ), 1 )
278*
279* B1
280*
281 CALL cgemv( 'C', m-l, i-1, alpha, b, ldb, b( 1, i ), 1,
282 $ one, t( 1, i ), 1 )
283*
284* T(1:I-1,I) := T(1:I-1,1:I-1) * T(1:I-1,I)
285*
286 CALL ctrmv( 'U', 'N', 'N', i-1, t, ldt, t( 1, i ), 1 )
287*
288* T(I,I) = tau(I)
289*
290 t( i, i ) = t( i, 1 )
291 t( i, 1 ) = zero
292 END DO
293
294*
295* End of CTPQRT2
296*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
CGEMV
Definition cgemv.f:160
subroutine cgerc(m, n, alpha, x, incx, y, incy, a, lda)
CGERC
Definition cgerc.f:130
subroutine clarfg(n, alpha, x, incx, tau)
CLARFG generates an elementary reflector (Householder matrix).
Definition clarfg.f:104
subroutine ctrmv(uplo, trans, diag, n, a, lda, x, incx)
CTRMV
Definition ctrmv.f:147
Here is the call graph for this function:
Here is the caller graph for this function: