LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dtpmqrt.f
Go to the documentation of this file.
1*> \brief \b DTPMQRT
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DTPMQRT + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dtpmqrt.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dtpmqrt.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dtpmqrt.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DTPMQRT( SIDE, TRANS, M, N, K, L, NB, V, LDV, T, LDT,
20* A, LDA, B, LDB, WORK, INFO )
21*
22* .. Scalar Arguments ..
23* CHARACTER SIDE, TRANS
24* INTEGER INFO, K, LDV, LDA, LDB, M, N, L, NB, LDT
25* ..
26* .. Array Arguments ..
27* DOUBLE PRECISION V( LDV, * ), A( LDA, * ), B( LDB, * ),
28* $ T( LDT, * ), WORK( * )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> DTPMQRT applies a real orthogonal matrix Q obtained from a
38*> "triangular-pentagonal" real block reflector H to a general
39*> real matrix C, which consists of two blocks A and B.
40*> \endverbatim
41*
42* Arguments:
43* ==========
44*
45*> \param[in] SIDE
46*> \verbatim
47*> SIDE is CHARACTER*1
48*> = 'L': apply Q or Q**T from the Left;
49*> = 'R': apply Q or Q**T from the Right.
50*> \endverbatim
51*>
52*> \param[in] TRANS
53*> \verbatim
54*> TRANS is CHARACTER*1
55*> = 'N': No transpose, apply Q;
56*> = 'T': Transpose, apply Q**T.
57*> \endverbatim
58*>
59*> \param[in] M
60*> \verbatim
61*> M is INTEGER
62*> The number of rows of the matrix B. M >= 0.
63*> \endverbatim
64*>
65*> \param[in] N
66*> \verbatim
67*> N is INTEGER
68*> The number of columns of the matrix B. N >= 0.
69*> \endverbatim
70*>
71*> \param[in] K
72*> \verbatim
73*> K is INTEGER
74*> The number of elementary reflectors whose product defines
75*> the matrix Q.
76*> \endverbatim
77*>
78*> \param[in] L
79*> \verbatim
80*> L is INTEGER
81*> The order of the trapezoidal part of V.
82*> K >= L >= 0. See Further Details.
83*> \endverbatim
84*>
85*> \param[in] NB
86*> \verbatim
87*> NB is INTEGER
88*> The block size used for the storage of T. K >= NB >= 1.
89*> This must be the same value of NB used to generate T
90*> in CTPQRT.
91*> \endverbatim
92*>
93*> \param[in] V
94*> \verbatim
95*> V is DOUBLE PRECISION array, dimension (LDV,K)
96*> The i-th column must contain the vector which defines the
97*> elementary reflector H(i), for i = 1,2,...,k, as returned by
98*> CTPQRT in B. See Further Details.
99*> \endverbatim
100*>
101*> \param[in] LDV
102*> \verbatim
103*> LDV is INTEGER
104*> The leading dimension of the array V.
105*> If SIDE = 'L', LDV >= max(1,M);
106*> if SIDE = 'R', LDV >= max(1,N).
107*> \endverbatim
108*>
109*> \param[in] T
110*> \verbatim
111*> T is DOUBLE PRECISION array, dimension (LDT,K)
112*> The upper triangular factors of the block reflectors
113*> as returned by CTPQRT, stored as a NB-by-K matrix.
114*> \endverbatim
115*>
116*> \param[in] LDT
117*> \verbatim
118*> LDT is INTEGER
119*> The leading dimension of the array T. LDT >= NB.
120*> \endverbatim
121*>
122*> \param[in,out] A
123*> \verbatim
124*> A is DOUBLE PRECISION array, dimension
125*> (LDA,N) if SIDE = 'L' or
126*> (LDA,K) if SIDE = 'R'
127*> On entry, the K-by-N or M-by-K matrix A.
128*> On exit, A is overwritten by the corresponding block of
129*> Q*C or Q**T*C or C*Q or C*Q**T. See Further Details.
130*> \endverbatim
131*>
132*> \param[in] LDA
133*> \verbatim
134*> LDA is INTEGER
135*> The leading dimension of the array A.
136*> If SIDE = 'L', LDC >= max(1,K);
137*> If SIDE = 'R', LDC >= max(1,M).
138*> \endverbatim
139*>
140*> \param[in,out] B
141*> \verbatim
142*> B is DOUBLE PRECISION array, dimension (LDB,N)
143*> On entry, the M-by-N matrix B.
144*> On exit, B is overwritten by the corresponding block of
145*> Q*C or Q**T*C or C*Q or C*Q**T. See Further Details.
146*> \endverbatim
147*>
148*> \param[in] LDB
149*> \verbatim
150*> LDB is INTEGER
151*> The leading dimension of the array B.
152*> LDB >= max(1,M).
153*> \endverbatim
154*>
155*> \param[out] WORK
156*> \verbatim
157*> WORK is DOUBLE PRECISION array. The dimension of WORK is
158*> N*NB if SIDE = 'L', or M*NB if SIDE = 'R'.
159*> \endverbatim
160*>
161*> \param[out] INFO
162*> \verbatim
163*> INFO is INTEGER
164*> = 0: successful exit
165*> < 0: if INFO = -i, the i-th argument had an illegal value
166*> \endverbatim
167*
168* Authors:
169* ========
170*
171*> \author Univ. of Tennessee
172*> \author Univ. of California Berkeley
173*> \author Univ. of Colorado Denver
174*> \author NAG Ltd.
175*
176*> \ingroup tpmqrt
177*
178*> \par Further Details:
179* =====================
180*>
181*> \verbatim
182*>
183*> The columns of the pentagonal matrix V contain the elementary reflectors
184*> H(1), H(2), ..., H(K); V is composed of a rectangular block V1 and a
185*> trapezoidal block V2:
186*>
187*> V = [V1]
188*> [V2].
189*>
190*> The size of the trapezoidal block V2 is determined by the parameter L,
191*> where 0 <= L <= K; V2 is upper trapezoidal, consisting of the first L
192*> rows of a K-by-K upper triangular matrix. If L=K, V2 is upper triangular;
193*> if L=0, there is no trapezoidal block, hence V = V1 is rectangular.
194*>
195*> If SIDE = 'L': C = [A] where A is K-by-N, B is M-by-N and V is M-by-K.
196*> [B]
197*>
198*> If SIDE = 'R': C = [A B] where A is M-by-K, B is M-by-N and V is N-by-K.
199*>
200*> The real orthogonal matrix Q is formed from V and T.
201*>
202*> If TRANS='N' and SIDE='L', C is on exit replaced with Q * C.
203*>
204*> If TRANS='T' and SIDE='L', C is on exit replaced with Q**T * C.
205*>
206*> If TRANS='N' and SIDE='R', C is on exit replaced with C * Q.
207*>
208*> If TRANS='T' and SIDE='R', C is on exit replaced with C * Q**T.
209*> \endverbatim
210*>
211* =====================================================================
212 SUBROUTINE dtpmqrt( SIDE, TRANS, M, N, K, L, NB, V, LDV, T,
213 $ LDT,
214 $ A, LDA, B, LDB, WORK, INFO )
215*
216* -- LAPACK computational routine --
217* -- LAPACK is a software package provided by Univ. of Tennessee, --
218* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
219*
220* .. Scalar Arguments ..
221 CHARACTER SIDE, TRANS
222 INTEGER INFO, K, LDV, LDA, LDB, M, N, L, NB, LDT
223* ..
224* .. Array Arguments ..
225 DOUBLE PRECISION V( LDV, * ), A( LDA, * ), B( LDB, * ),
226 $ T( LDT, * ), WORK( * )
227* ..
228*
229* =====================================================================
230*
231* ..
232* .. Local Scalars ..
233 LOGICAL LEFT, RIGHT, TRAN, NOTRAN
234 INTEGER I, IB, MB, LB, KF, LDAQ, LDVQ
235* ..
236* .. External Functions ..
237 LOGICAL LSAME
238 EXTERNAL LSAME
239* ..
240* .. External Subroutines ..
241 EXTERNAL dtprfb, xerbla
242* ..
243* .. Intrinsic Functions ..
244 INTRINSIC max, min
245* ..
246* .. Executable Statements ..
247*
248* .. Test the input arguments ..
249*
250 info = 0
251 left = lsame( side, 'L' )
252 right = lsame( side, 'R' )
253 tran = lsame( trans, 'T' )
254 notran = lsame( trans, 'N' )
255*
256 IF ( left ) THEN
257 ldvq = max( 1, m )
258 ldaq = max( 1, k )
259 ELSE IF ( right ) THEN
260 ldvq = max( 1, n )
261 ldaq = max( 1, m )
262 END IF
263 IF( .NOT.left .AND. .NOT.right ) THEN
264 info = -1
265 ELSE IF( .NOT.tran .AND. .NOT.notran ) THEN
266 info = -2
267 ELSE IF( m.LT.0 ) THEN
268 info = -3
269 ELSE IF( n.LT.0 ) THEN
270 info = -4
271 ELSE IF( k.LT.0 ) THEN
272 info = -5
273 ELSE IF( l.LT.0 .OR. l.GT.k ) THEN
274 info = -6
275 ELSE IF( nb.LT.1 .OR. (nb.GT.k .AND. k.GT.0) ) THEN
276 info = -7
277 ELSE IF( ldv.LT.ldvq ) THEN
278 info = -9
279 ELSE IF( ldt.LT.nb ) THEN
280 info = -11
281 ELSE IF( lda.LT.ldaq ) THEN
282 info = -13
283 ELSE IF( ldb.LT.max( 1, m ) ) THEN
284 info = -15
285 END IF
286*
287 IF( info.NE.0 ) THEN
288 CALL xerbla( 'DTPMQRT', -info )
289 RETURN
290 END IF
291*
292* .. Quick return if possible ..
293*
294 IF( m.EQ.0 .OR. n.EQ.0 .OR. k.EQ.0 ) RETURN
295*
296 IF( left .AND. tran ) THEN
297*
298 DO i = 1, k, nb
299 ib = min( nb, k-i+1 )
300 mb = min( m-l+i+ib-1, m )
301 IF( i.GE.l ) THEN
302 lb = 0
303 ELSE
304 lb = mb-m+l-i+1
305 END IF
306 CALL dtprfb( 'L', 'T', 'F', 'C', mb, n, ib, lb,
307 $ v( 1, i ), ldv, t( 1, i ), ldt,
308 $ a( i, 1 ), lda, b, ldb, work, ib )
309 END DO
310*
311 ELSE IF( right .AND. notran ) THEN
312*
313 DO i = 1, k, nb
314 ib = min( nb, k-i+1 )
315 mb = min( n-l+i+ib-1, n )
316 IF( i.GE.l ) THEN
317 lb = 0
318 ELSE
319 lb = mb-n+l-i+1
320 END IF
321 CALL dtprfb( 'R', 'N', 'F', 'C', m, mb, ib, lb,
322 $ v( 1, i ), ldv, t( 1, i ), ldt,
323 $ a( 1, i ), lda, b, ldb, work, m )
324 END DO
325*
326 ELSE IF( left .AND. notran ) THEN
327*
328 kf = ((k-1)/nb)*nb+1
329 DO i = kf, 1, -nb
330 ib = min( nb, k-i+1 )
331 mb = min( m-l+i+ib-1, m )
332 IF( i.GE.l ) THEN
333 lb = 0
334 ELSE
335 lb = mb-m+l-i+1
336 END IF
337 CALL dtprfb( 'L', 'N', 'F', 'C', mb, n, ib, lb,
338 $ v( 1, i ), ldv, t( 1, i ), ldt,
339 $ a( i, 1 ), lda, b, ldb, work, ib )
340 END DO
341*
342 ELSE IF( right .AND. tran ) THEN
343*
344 kf = ((k-1)/nb)*nb+1
345 DO i = kf, 1, -nb
346 ib = min( nb, k-i+1 )
347 mb = min( n-l+i+ib-1, n )
348 IF( i.GE.l ) THEN
349 lb = 0
350 ELSE
351 lb = mb-n+l-i+1
352 END IF
353 CALL dtprfb( 'R', 'T', 'F', 'C', m, mb, ib, lb,
354 $ v( 1, i ), ldv, t( 1, i ), ldt,
355 $ a( 1, i ), lda, b, ldb, work, m )
356 END DO
357*
358 END IF
359*
360 RETURN
361*
362* End of DTPMQRT
363*
364 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dtpmqrt(side, trans, m, n, k, l, nb, v, ldv, t, ldt, a, lda, b, ldb, work, info)
DTPMQRT
Definition dtpmqrt.f:215
subroutine dtprfb(side, trans, direct, storev, m, n, k, l, v, ldv, t, ldt, a, lda, b, ldb, work, ldwork)
DTPRFB applies a real "triangular-pentagonal" block reflector to a real matrix, which is composed of ...
Definition dtprfb.f:249