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