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

◆ dtpqrt2()

subroutine dtpqrt2 ( 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 
)

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

Purpose:
 DTPQRT2 computes a QR 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 upper trapezoidal part of B.
          MIN(M,N) >= L >= 0.  See Further Details.
[in,out]A
          A is DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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**T

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

Definition at line 172 of file dtpqrt2.f.

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