LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.
Date
September 2012
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 175 of file dtpqrt2.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: