LAPACK 3.12.0
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 "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 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 176 of file dtplqt2.f.

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