LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
dtplqt2.f
Go to the documentation of this file.
1 *> \brief \b 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.
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download DTPLQT2 + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dtplqt2.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dtplqt2.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dtplqt2.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE DTPLQT2( M, N, L, A, LDA, B, LDB, T, LDT, INFO )
22 *
23 * .. Scalar Arguments ..
24 * INTEGER INFO, LDA, LDB, LDT, N, M, L
25 * ..
26 * .. Array Arguments ..
27 * DOUBLE PRECISION A( LDA, * ), B( LDB, * ), T( LDT, * )
28 * ..
29 *
30 *
31 *> \par Purpose:
32 * =============
33 *>
34 *> \verbatim
35 *>
36 *> DTPLQT2 computes a LQ a factorization of a real "triangular-pentagonal"
37 *> matrix C, which is composed of a triangular block A and pentagonal block B,
38 *> using the compact WY representation for Q.
39 *> \endverbatim
40 *
41 * Arguments:
42 * ==========
43 *
44 *> \param[in] M
45 *> \verbatim
46 *> M is INTEGER
47 *> The total number of rows of the matrix B.
48 *> M >= 0.
49 *> \endverbatim
50 *>
51 *> \param[in] N
52 *> \verbatim
53 *> N is INTEGER
54 *> The number of columns of the matrix B, and the order of
55 *> the triangular matrix A.
56 *> N >= 0.
57 *> \endverbatim
58 *>
59 *> \param[in] L
60 *> \verbatim
61 *> L is INTEGER
62 *> The number of rows of the lower trapezoidal part of B.
63 *> MIN(M,N) >= L >= 0. See Further Details.
64 *> \endverbatim
65 *>
66 *> \param[in,out] A
67 *> \verbatim
68 *> A is DOUBLE PRECISION array, dimension (LDA,N)
69 *> On entry, the lower triangular M-by-M matrix A.
70 *> On exit, the elements on and below the diagonal of the array
71 *> contain the lower triangular matrix L.
72 *> \endverbatim
73 *>
74 *> \param[in] LDA
75 *> \verbatim
76 *> LDA is INTEGER
77 *> The leading dimension of the array A. LDA >= max(1,N).
78 *> \endverbatim
79 *>
80 *> \param[in,out] B
81 *> \verbatim
82 *> B is DOUBLE PRECISION array, dimension (LDB,N)
83 *> On entry, the pentagonal M-by-N matrix B. The first N-L columns
84 *> are rectangular, and the last L columns are lower trapezoidal.
85 *> On exit, B contains the pentagonal matrix V. See Further Details.
86 *> \endverbatim
87 *>
88 *> \param[in] LDB
89 *> \verbatim
90 *> LDB is INTEGER
91 *> The leading dimension of the array B. LDB >= max(1,M).
92 *> \endverbatim
93 *>
94 *> \param[out] T
95 *> \verbatim
96 *> T is DOUBLE PRECISION array, dimension (LDT,M)
97 *> The N-by-N upper triangular factor T of the block reflector.
98 *> See Further Details.
99 *> \endverbatim
100 *>
101 *> \param[in] LDT
102 *> \verbatim
103 *> LDT is INTEGER
104 *> The leading dimension of the array T. LDT >= max(1,M)
105 *> \endverbatim
106 *>
107 *> \param[out] INFO
108 *> \verbatim
109 *> INFO is INTEGER
110 *> = 0: successful exit
111 *> < 0: if INFO = -i, the i-th argument had an illegal value
112 *> \endverbatim
113 *
114 * Authors:
115 * ========
116 *
117 *> \author Univ. of Tennessee
118 *> \author Univ. of California Berkeley
119 *> \author Univ. of Colorado Denver
120 *> \author NAG Ltd.
121 *
122 *> \date December 2016
123 *
124 *> \ingroup doubleOTHERcomputational
125 *
126 *> \par Further Details:
127 * =====================
128 *>
129 *> \verbatim
130 *>
131 *> The input matrix C is a M-by-(M+N) matrix
132 *>
133 *> C = [ A ][ B ]
134 *>
135 *>
136 *> where A is an lower triangular N-by-N matrix, and B is M-by-N pentagonal
137 *> matrix consisting of a M-by-(N-L) rectangular matrix B1 left of a M-by-L
138 *> upper trapezoidal matrix B2:
139 *>
140 *> B = [ B1 ][ B2 ]
141 *> [ B1 ] <- M-by-(N-L) rectangular
142 *> [ B2 ] <- M-by-L lower trapezoidal.
143 *>
144 *> The lower trapezoidal matrix B2 consists of the first L columns of a
145 *> N-by-N lower triangular matrix, where 0 <= L <= MIN(M,N). If L=0,
146 *> B is rectangular M-by-N; if M=L=N, B is lower triangular.
147 *>
148 *> The matrix W stores the elementary reflectors H(i) in the i-th row
149 *> above the diagonal (of A) in the M-by-(M+N) input matrix C
150 *>
151 *> C = [ A ][ B ]
152 *> [ A ] <- lower triangular N-by-N
153 *> [ B ] <- M-by-N pentagonal
154 *>
155 *> so that W can be represented as
156 *>
157 *> W = [ I ][ V ]
158 *> [ I ] <- identity, N-by-N
159 *> [ V ] <- M-by-N, same form as B.
160 *>
161 *> Thus, all of information needed for W is contained on exit in B, which
162 *> we call V above. Note that V has the same form as B; that is,
163 *>
164 *> W = [ V1 ][ V2 ]
165 *> [ V1 ] <- M-by-(N-L) rectangular
166 *> [ V2 ] <- M-by-L lower trapezoidal.
167 *>
168 *> The rows of V represent the vectors which define the H(i)'s.
169 *> The (M+N)-by-(M+N) block reflector H is then given by
170 *>
171 *> H = I - W**T * T * W
172 *>
173 *> where W^H is the conjugate transpose of W and T is the upper triangular
174 *> factor of the block reflector.
175 *> \endverbatim
176 *>
177 * =====================================================================
178  SUBROUTINE dtplqt2( M, N, L, A, LDA, B, LDB, T, LDT, INFO )
179 *
180 * -- LAPACK computational routine (version 3.7.0) --
181 * -- LAPACK is a software package provided by Univ. of Tennessee, --
182 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
183 * December 2016
184 *
185 * .. Scalar Arguments ..
186  INTEGER INFO, LDA, LDB, LDT, N, M, L
187 * ..
188 * .. Array Arguments ..
189  DOUBLE PRECISION A( lda, * ), B( ldb, * ), T( ldt, * )
190 * ..
191 *
192 * =====================================================================
193 *
194 * .. Parameters ..
195  DOUBLE PRECISION ONE, ZERO
196  parameter( one = 1.0, zero = 0.0 )
197 * ..
198 * .. Local Scalars ..
199  INTEGER I, J, P, MP, NP
200  DOUBLE PRECISION ALPHA
201 * ..
202 * .. External Subroutines ..
203  EXTERNAL dlarfg, dgemv, dger, dtrmv, xerbla
204 * ..
205 * .. Intrinsic Functions ..
206  INTRINSIC max, min
207 * ..
208 * .. Executable Statements ..
209 *
210 * Test the input arguments
211 *
212  info = 0
213  IF( m.LT.0 ) THEN
214  info = -1
215  ELSE IF( n.LT.0 ) THEN
216  info = -2
217  ELSE IF( l.LT.0 .OR. l.GT.min(m,n) ) THEN
218  info = -3
219  ELSE IF( lda.LT.max( 1, m ) ) THEN
220  info = -5
221  ELSE IF( ldb.LT.max( 1, m ) ) THEN
222  info = -7
223  ELSE IF( ldt.LT.max( 1, m ) ) THEN
224  info = -9
225  END IF
226  IF( info.NE.0 ) THEN
227  CALL xerbla( 'DTPLQT2', -info )
228  RETURN
229  END IF
230 *
231 * Quick return if possible
232 *
233  IF( n.EQ.0 .OR. m.EQ.0 ) RETURN
234 *
235  DO i = 1, m
236 *
237 * Generate elementary reflector H(I) to annihilate B(I,:)
238 *
239  p = n-l+min( l, i )
240  CALL dlarfg( p+1, a( i, i ), b( i, 1 ), ldb, t( 1, i ) )
241  IF( i.LT.m ) THEN
242 *
243 * W(M-I:1) := C(I+1:M,I:N) * C(I,I:N) [use W = T(M,:)]
244 *
245  DO j = 1, m-i
246  t( m, j ) = (a( i+j, i ))
247  END DO
248  CALL dgemv( 'N', m-i, p, one, b( i+1, 1 ), ldb,
249  $ b( i, 1 ), ldb, one, t( m, 1 ), ldt )
250 *
251 * C(I+1:M,I:N) = C(I+1:M,I:N) + alpha * C(I,I:N)*W(M-1:1)^H
252 *
253  alpha = -(t( 1, i ))
254  DO j = 1, m-i
255  a( i+j, i ) = a( i+j, i ) + alpha*(t( m, j ))
256  END DO
257  CALL dger( m-i, p, alpha, t( m, 1 ), ldt,
258  $ b( i, 1 ), ldb, b( i+1, 1 ), ldb )
259  END IF
260  END DO
261 *
262  DO i = 2, m
263 *
264 * T(I,1:I-1) := C(I:I-1,1:N) * (alpha * C(I,I:N)^H)
265 *
266  alpha = -t( 1, i )
267 
268  DO j = 1, i-1
269  t( i, j ) = zero
270  END DO
271  p = min( i-1, l )
272  np = min( n-l+1, n )
273  mp = min( p+1, m )
274 *
275 * Triangular part of B2
276 *
277  DO j = 1, p
278  t( i, j ) = alpha*b( i, n-l+j )
279  END DO
280  CALL dtrmv( 'L', 'N', 'N', p, b( 1, np ), ldb,
281  $ t( i, 1 ), ldt )
282 *
283 * Rectangular part of B2
284 *
285  CALL dgemv( 'N', i-1-p, l, alpha, b( mp, np ), ldb,
286  $ b( i, np ), ldb, zero, t( i,mp ), ldt )
287 *
288 * B1
289 *
290  CALL dgemv( 'N', i-1, n-l, alpha, b, ldb, b( i, 1 ), ldb,
291  $ one, t( i, 1 ), ldt )
292 *
293 * T(1:I-1,I) := T(1:I-1,1:I-1) * T(I,1:I-1)
294 *
295  CALL dtrmv( 'L', 'T', 'N', i-1, t, ldt, t( i, 1 ), ldt )
296 *
297 * T(I,I) = tau(I)
298 *
299  t( i, i ) = t( 1, i )
300  t( 1, i ) = zero
301  END DO
302  DO i=1,m
303  DO j= i+1,m
304  t(i,j)=t(j,i)
305  t(j,i)= zero
306  END DO
307  END DO
308 
309 *
310 * End of DTPLQT2
311 *
312  END
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:132
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dlarfg(N, ALPHA, X, INCX, TAU)
DLARFG generates an elementary reflector (Householder matrix).
Definition: dlarfg.f:108
subroutine dtplqt2(M, N, L, A, LDA, B, LDB, T, LDT, 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.
Definition: dtplqt2.f:179
subroutine dtrmv(UPLO, TRANS, DIAG, N, A, LDA, X, INCX)
DTRMV
Definition: dtrmv.f:149