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

◆ dtplqt2()

subroutine dtplqt2 ( integer m,
integer n,
integer l,
double precision, dimension( lda, * ) a,
integer lda,
double precision, dimension( ldb, * ) b,
integer ldb,
double precision, dimension( ldt, * ) t,
integer ldt,
integer info )

DTPLQT2 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 DTPLQT2 + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> DTPLQT2 computes a LQ a factorization of a real 
!> 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 DOUBLE PRECISION array, dimension (LDA,M)
!>          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,M).
!> 
[in,out]B
!>          B is DOUBLE PRECISION 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 DOUBLE PRECISION 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.
Further Details:
!>
!>  The input matrix C is a M-by-(M+N) matrix
!>
!>               C = [ A ][ B ]
!>
!>
!>  where A is an lower triangular M-by-M 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 M-by-M
!>                   [ B ]  <- M-by-N pentagonal
!>
!>  so that W can be represented as
!>
!>               W = [ I ][ V ]
!>                   [ I ]  <- identity, M-by-M
!>                   [ 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 174 of file dtplqt2.f.

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