LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zlarfb_gett.f
Go to the documentation of this file.
1*> \brief \b ZLARFB_GETT
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZLARFB_GETT + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zlarfb_gett.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zlarfb_gett.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zlarfb_gett.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZLARFB_GETT( IDENT, M, N, K, T, LDT, A, LDA, B, LDB,
20* $ WORK, LDWORK )
21* IMPLICIT NONE
22*
23* .. Scalar Arguments ..
24* CHARACTER IDENT
25* INTEGER K, LDA, LDB, LDT, LDWORK, M, N
26* ..
27* .. Array Arguments ..
28* COMPLEX*16 A( LDA, * ), B( LDB, * ), T( LDT, * ),
29* $ WORK( LDWORK, * )
30* ..
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> ZLARFB_GETT applies a complex Householder block reflector H from the
38*> left to a complex (K+M)-by-N "triangular-pentagonal" matrix
39*> composed of two block matrices: an upper trapezoidal K-by-N matrix A
40*> stored in the array A, and a rectangular M-by-(N-K) matrix B, stored
41*> in the array B. The block reflector H is stored in a compact
42*> WY-representation, where the elementary reflectors are in the
43*> arrays A, B and T. See Further Details section.
44*> \endverbatim
45*
46* Arguments:
47* ==========
48*
49*> \param[in] IDENT
50*> \verbatim
51*> IDENT is CHARACTER*1
52*> If IDENT = not 'I', or not 'i', then V1 is unit
53*> lower-triangular and stored in the left K-by-K block of
54*> the input matrix A,
55*> If IDENT = 'I' or 'i', then V1 is an identity matrix and
56*> not stored.
57*> See Further Details section.
58*> \endverbatim
59*>
60*> \param[in] M
61*> \verbatim
62*> M is INTEGER
63*> The number of rows of the matrix B.
64*> M >= 0.
65*> \endverbatim
66*>
67*> \param[in] N
68*> \verbatim
69*> N is INTEGER
70*> The number of columns of the matrices A and B.
71*> N >= 0.
72*> \endverbatim
73*>
74*> \param[in] K
75*> \verbatim
76*> K is INTEGER
77*> The number or rows of the matrix A.
78*> K is also order of the matrix T, i.e. the number of
79*> elementary reflectors whose product defines the block
80*> reflector. 0 <= K <= N.
81*> \endverbatim
82*>
83*> \param[in] T
84*> \verbatim
85*> T is COMPLEX*16 array, dimension (LDT,K)
86*> The upper-triangular K-by-K matrix T in the representation
87*> of the block reflector.
88*> \endverbatim
89*>
90*> \param[in] LDT
91*> \verbatim
92*> LDT is INTEGER
93*> The leading dimension of the array T. LDT >= K.
94*> \endverbatim
95*>
96*> \param[in,out] A
97*> \verbatim
98*> A is COMPLEX*16 array, dimension (LDA,N)
99*>
100*> On entry:
101*> a) In the K-by-N upper-trapezoidal part A: input matrix A.
102*> b) In the columns below the diagonal: columns of V1
103*> (ones are not stored on the diagonal).
104*>
105*> On exit:
106*> A is overwritten by rectangular K-by-N product H*A.
107*>
108*> See Further Details section.
109*> \endverbatim
110*>
111*> \param[in] LDA
112*> \verbatim
113*> LDB is INTEGER
114*> The leading dimension of the array A. LDA >= max(1,K).
115*> \endverbatim
116*>
117*> \param[in,out] B
118*> \verbatim
119*> B is COMPLEX*16 array, dimension (LDB,N)
120*>
121*> On entry:
122*> a) In the M-by-(N-K) right block: input matrix B.
123*> b) In the M-by-N left block: columns of V2.
124*>
125*> On exit:
126*> B is overwritten by rectangular M-by-N product H*B.
127*>
128*> See Further Details section.
129*> \endverbatim
130*>
131*> \param[in] LDB
132*> \verbatim
133*> LDB is INTEGER
134*> The leading dimension of the array B. LDB >= max(1,M).
135*> \endverbatim
136*>
137*> \param[out] WORK
138*> \verbatim
139*> WORK is COMPLEX*16 array,
140*> dimension (LDWORK,max(K,N-K))
141*> \endverbatim
142*>
143*> \param[in] LDWORK
144*> \verbatim
145*> LDWORK is INTEGER
146*> The leading dimension of the array WORK. LDWORK>=max(1,K).
147*>
148*> \endverbatim
149*
150* Authors:
151* ========
152*
153*> \author Univ. of Tennessee
154*> \author Univ. of California Berkeley
155*> \author Univ. of Colorado Denver
156*> \author NAG Ltd.
157*
158*> \ingroup larfb_gett
159*
160*> \par Contributors:
161* ==================
162*>
163*> \verbatim
164*>
165*> November 2020, Igor Kozachenko,
166*> Computer Science Division,
167*> University of California, Berkeley
168*>
169*> \endverbatim
170*
171*> \par Further Details:
172* =====================
173*>
174*> \verbatim
175*>
176*> (1) Description of the Algebraic Operation.
177*>
178*> The matrix A is a K-by-N matrix composed of two column block
179*> matrices, A1, which is K-by-K, and A2, which is K-by-(N-K):
180*> A = ( A1, A2 ).
181*> The matrix B is an M-by-N matrix composed of two column block
182*> matrices, B1, which is M-by-K, and B2, which is M-by-(N-K):
183*> B = ( B1, B2 ).
184*>
185*> Perform the operation:
186*>
187*> ( A_out ) := H * ( A_in ) = ( I - V * T * V**H ) * ( A_in ) =
188*> ( B_out ) ( B_in ) ( B_in )
189*> = ( I - ( V1 ) * T * ( V1**H, V2**H ) ) * ( A_in )
190*> ( V2 ) ( B_in )
191*> On input:
192*>
193*> a) ( A_in ) consists of two block columns:
194*> ( B_in )
195*>
196*> ( A_in ) = (( A1_in ) ( A2_in )) = (( A1_in ) ( A2_in ))
197*> ( B_in ) (( B1_in ) ( B2_in )) (( 0 ) ( B2_in )),
198*>
199*> where the column blocks are:
200*>
201*> ( A1_in ) is a K-by-K upper-triangular matrix stored in the
202*> upper triangular part of the array A(1:K,1:K).
203*> ( B1_in ) is an M-by-K rectangular ZERO matrix and not stored.
204*>
205*> ( A2_in ) is a K-by-(N-K) rectangular matrix stored
206*> in the array A(1:K,K+1:N).
207*> ( B2_in ) is an M-by-(N-K) rectangular matrix stored
208*> in the array B(1:M,K+1:N).
209*>
210*> b) V = ( V1 )
211*> ( V2 )
212*>
213*> where:
214*> 1) if IDENT == 'I',V1 is a K-by-K identity matrix, not stored;
215*> 2) if IDENT != 'I',V1 is a K-by-K unit lower-triangular matrix,
216*> stored in the lower-triangular part of the array
217*> A(1:K,1:K) (ones are not stored),
218*> and V2 is an M-by-K rectangular stored the array B(1:M,1:K),
219*> (because on input B1_in is a rectangular zero
220*> matrix that is not stored and the space is
221*> used to store V2).
222*>
223*> c) T is a K-by-K upper-triangular matrix stored
224*> in the array T(1:K,1:K).
225*>
226*> On output:
227*>
228*> a) ( A_out ) consists of two block columns:
229*> ( B_out )
230*>
231*> ( A_out ) = (( A1_out ) ( A2_out ))
232*> ( B_out ) (( B1_out ) ( B2_out )),
233*>
234*> where the column blocks are:
235*>
236*> ( A1_out ) is a K-by-K square matrix, or a K-by-K
237*> upper-triangular matrix, if V1 is an
238*> identity matrix. AiOut is stored in
239*> the array A(1:K,1:K).
240*> ( B1_out ) is an M-by-K rectangular matrix stored
241*> in the array B(1:M,K:N).
242*>
243*> ( A2_out ) is a K-by-(N-K) rectangular matrix stored
244*> in the array A(1:K,K+1:N).
245*> ( B2_out ) is an M-by-(N-K) rectangular matrix stored
246*> in the array B(1:M,K+1:N).
247*>
248*>
249*> The operation above can be represented as the same operation
250*> on each block column:
251*>
252*> ( A1_out ) := H * ( A1_in ) = ( I - V * T * V**H ) * ( A1_in )
253*> ( B1_out ) ( 0 ) ( 0 )
254*>
255*> ( A2_out ) := H * ( A2_in ) = ( I - V * T * V**H ) * ( A2_in )
256*> ( B2_out ) ( B2_in ) ( B2_in )
257*>
258*> If IDENT != 'I':
259*>
260*> The computation for column block 1:
261*>
262*> A1_out: = A1_in - V1*T*(V1**H)*A1_in
263*>
264*> B1_out: = - V2*T*(V1**H)*A1_in
265*>
266*> The computation for column block 2, which exists if N > K:
267*>
268*> A2_out: = A2_in - V1*T*( (V1**H)*A2_in + (V2**H)*B2_in )
269*>
270*> B2_out: = B2_in - V2*T*( (V1**H)*A2_in + (V2**H)*B2_in )
271*>
272*> If IDENT == 'I':
273*>
274*> The operation for column block 1:
275*>
276*> A1_out: = A1_in - V1*T*A1_in
277*>
278*> B1_out: = - V2*T*A1_in
279*>
280*> The computation for column block 2, which exists if N > K:
281*>
282*> A2_out: = A2_in - T*( A2_in + (V2**H)*B2_in )
283*>
284*> B2_out: = B2_in - V2*T*( A2_in + (V2**H)*B2_in )
285*>
286*> (2) Description of the Algorithmic Computation.
287*>
288*> In the first step, we compute column block 2, i.e. A2 and B2.
289*> Here, we need to use the K-by-(N-K) rectangular workspace
290*> matrix W2 that is of the same size as the matrix A2.
291*> W2 is stored in the array WORK(1:K,1:(N-K)).
292*>
293*> In the second step, we compute column block 1, i.e. A1 and B1.
294*> Here, we need to use the K-by-K square workspace matrix W1
295*> that is of the same size as the as the matrix A1.
296*> W1 is stored in the array WORK(1:K,1:K).
297*>
298*> NOTE: Hence, in this routine, we need the workspace array WORK
299*> only of size WORK(1:K,1:max(K,N-K)) so it can hold both W2 from
300*> the first step and W1 from the second step.
301*>
302*> Case (A), when V1 is unit lower-triangular, i.e. IDENT != 'I',
303*> more computations than in the Case (B).
304*>
305*> if( IDENT != 'I' ) then
306*> if ( N > K ) then
307*> (First Step - column block 2)
308*> col2_(1) W2: = A2
309*> col2_(2) W2: = (V1**H) * W2 = (unit_lower_tr_of_(A1)**H) * W2
310*> col2_(3) W2: = W2 + (V2**H) * B2 = W2 + (B1**H) * B2
311*> col2_(4) W2: = T * W2
312*> col2_(5) B2: = B2 - V2 * W2 = B2 - B1 * W2
313*> col2_(6) W2: = V1 * W2 = unit_lower_tr_of_(A1) * W2
314*> col2_(7) A2: = A2 - W2
315*> else
316*> (Second Step - column block 1)
317*> col1_(1) W1: = A1
318*> col1_(2) W1: = (V1**H) * W1 = (unit_lower_tr_of_(A1)**H) * W1
319*> col1_(3) W1: = T * W1
320*> col1_(4) B1: = - V2 * W1 = - B1 * W1
321*> col1_(5) square W1: = V1 * W1 = unit_lower_tr_of_(A1) * W1
322*> col1_(6) square A1: = A1 - W1
323*> end if
324*> end if
325*>
326*> Case (B), when V1 is an identity matrix, i.e. IDENT == 'I',
327*> less computations than in the Case (A)
328*>
329*> if( IDENT == 'I' ) then
330*> if ( N > K ) then
331*> (First Step - column block 2)
332*> col2_(1) W2: = A2
333*> col2_(3) W2: = W2 + (V2**H) * B2 = W2 + (B1**H) * B2
334*> col2_(4) W2: = T * W2
335*> col2_(5) B2: = B2 - V2 * W2 = B2 - B1 * W2
336*> col2_(7) A2: = A2 - W2
337*> else
338*> (Second Step - column block 1)
339*> col1_(1) W1: = A1
340*> col1_(3) W1: = T * W1
341*> col1_(4) B1: = - V2 * W1 = - B1 * W1
342*> col1_(6) upper-triangular_of_(A1): = A1 - W1
343*> end if
344*> end if
345*>
346*> Combine these cases (A) and (B) together, this is the resulting
347*> algorithm:
348*>
349*> if ( N > K ) then
350*>
351*> (First Step - column block 2)
352*>
353*> col2_(1) W2: = A2
354*> if( IDENT != 'I' ) then
355*> col2_(2) W2: = (V1**H) * W2
356*> = (unit_lower_tr_of_(A1)**H) * W2
357*> end if
358*> col2_(3) W2: = W2 + (V2**H) * B2 = W2 + (B1**H) * B2]
359*> col2_(4) W2: = T * W2
360*> col2_(5) B2: = B2 - V2 * W2 = B2 - B1 * W2
361*> if( IDENT != 'I' ) then
362*> col2_(6) W2: = V1 * W2 = unit_lower_tr_of_(A1) * W2
363*> end if
364*> col2_(7) A2: = A2 - W2
365*>
366*> else
367*>
368*> (Second Step - column block 1)
369*>
370*> col1_(1) W1: = A1
371*> if( IDENT != 'I' ) then
372*> col1_(2) W1: = (V1**H) * W1
373*> = (unit_lower_tr_of_(A1)**H) * W1
374*> end if
375*> col1_(3) W1: = T * W1
376*> col1_(4) B1: = - V2 * W1 = - B1 * W1
377*> if( IDENT != 'I' ) then
378*> col1_(5) square W1: = V1 * W1 = unit_lower_tr_of_(A1) * W1
379*> col1_(6_a) below_diag_of_(A1): = - below_diag_of_(W1)
380*> end if
381*> col1_(6_b) up_tr_of_(A1): = up_tr_of_(A1) - up_tr_of_(W1)
382*>
383*> end if
384*>
385*> \endverbatim
386*>
387* =====================================================================
388 SUBROUTINE zlarfb_gett( IDENT, M, N, K, T, LDT, A, LDA, B, LDB,
389 $ WORK, LDWORK )
390 IMPLICIT NONE
391*
392* -- LAPACK auxiliary routine --
393* -- LAPACK is a software package provided by Univ. of Tennessee, --
394* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
395*
396* .. Scalar Arguments ..
397 CHARACTER IDENT
398 INTEGER K, LDA, LDB, LDT, LDWORK, M, N
399* ..
400* .. Array Arguments ..
401 COMPLEX*16 A( LDA, * ), B( LDB, * ), T( LDT, * ),
402 $ work( ldwork, * )
403* ..
404*
405* =====================================================================
406*
407* .. Parameters ..
408 COMPLEX*16 CONE, CZERO
409 parameter( cone = ( 1.0d+0, 0.0d+0 ),
410 $ czero = ( 0.0d+0, 0.0d+0 ) )
411* ..
412* .. Local Scalars ..
413 LOGICAL LNOTIDENT
414 INTEGER I, J
415* ..
416* .. EXTERNAL FUNCTIONS ..
417 LOGICAL LSAME
418 EXTERNAL lsame
419* ..
420* .. External Subroutines ..
421 EXTERNAL zcopy, zgemm, ztrmm
422* ..
423* .. Executable Statements ..
424*
425* Quick return if possible
426*
427 IF( m.LT.0 .OR. n.LE.0 .OR. k.EQ.0 .OR. k.GT.n )
428 $ RETURN
429*
430 lnotident = .NOT.lsame( ident, 'I' )
431*
432* ------------------------------------------------------------------
433*
434* First Step. Computation of the Column Block 2:
435*
436* ( A2 ) := H * ( A2 )
437* ( B2 ) ( B2 )
438*
439* ------------------------------------------------------------------
440*
441 IF( n.GT.k ) THEN
442*
443* col2_(1) Compute W2: = A2. Therefore, copy A2 = A(1:K, K+1:N)
444* into W2=WORK(1:K, 1:N-K) column-by-column.
445*
446 DO j = 1, n-k
447 CALL zcopy( k, a( 1, k+j ), 1, work( 1, j ), 1 )
448 END DO
449
450 IF( lnotident ) THEN
451*
452* col2_(2) Compute W2: = (V1**H) * W2 = (A1**H) * W2,
453* V1 is not an identity matrix, but unit lower-triangular
454* V1 stored in A1 (diagonal ones are not stored).
455*
456*
457 CALL ztrmm( 'L', 'L', 'C', 'U', k, n-k, cone, a, lda,
458 $ work, ldwork )
459 END IF
460*
461* col2_(3) Compute W2: = W2 + (V2**H) * B2 = W2 + (B1**H) * B2
462* V2 stored in B1.
463*
464 IF( m.GT.0 ) THEN
465 CALL zgemm( 'C', 'N', k, n-k, m, cone, b, ldb,
466 $ b( 1, k+1 ), ldb, cone, work, ldwork )
467 END IF
468*
469* col2_(4) Compute W2: = T * W2,
470* T is upper-triangular.
471*
472 CALL ztrmm( 'L', 'U', 'N', 'N', k, n-k, cone, t, ldt,
473 $ work, ldwork )
474*
475* col2_(5) Compute B2: = B2 - V2 * W2 = B2 - B1 * W2,
476* V2 stored in B1.
477*
478 IF( m.GT.0 ) THEN
479 CALL zgemm( 'N', 'N', m, n-k, k, -cone, b, ldb,
480 $ work, ldwork, cone, b( 1, k+1 ), ldb )
481 END IF
482*
483 IF( lnotident ) THEN
484*
485* col2_(6) Compute W2: = V1 * W2 = A1 * W2,
486* V1 is not an identity matrix, but unit lower-triangular,
487* V1 stored in A1 (diagonal ones are not stored).
488*
489 CALL ztrmm( 'L', 'L', 'N', 'U', k, n-k, cone, a, lda,
490 $ work, ldwork )
491 END IF
492*
493* col2_(7) Compute A2: = A2 - W2 =
494* = A(1:K, K+1:N-K) - WORK(1:K, 1:N-K),
495* column-by-column.
496*
497 DO j = 1, n-k
498 DO i = 1, k
499 a( i, k+j ) = a( i, k+j ) - work( i, j )
500 END DO
501 END DO
502*
503 END IF
504*
505* ------------------------------------------------------------------
506*
507* Second Step. Computation of the Column Block 1:
508*
509* ( A1 ) := H * ( A1 )
510* ( B1 ) ( 0 )
511*
512* ------------------------------------------------------------------
513*
514* col1_(1) Compute W1: = A1. Copy the upper-triangular
515* A1 = A(1:K, 1:K) into the upper-triangular
516* W1 = WORK(1:K, 1:K) column-by-column.
517*
518 DO j = 1, k
519 CALL zcopy( j, a( 1, j ), 1, work( 1, j ), 1 )
520 END DO
521*
522* Set the subdiagonal elements of W1 to zero column-by-column.
523*
524 DO j = 1, k - 1
525 DO i = j + 1, k
526 work( i, j ) = czero
527 END DO
528 END DO
529*
530 IF( lnotident ) THEN
531*
532* col1_(2) Compute W1: = (V1**H) * W1 = (A1**H) * W1,
533* V1 is not an identity matrix, but unit lower-triangular
534* V1 stored in A1 (diagonal ones are not stored),
535* W1 is upper-triangular with zeroes below the diagonal.
536*
537 CALL ztrmm( 'L', 'L', 'C', 'U', k, k, cone, a, lda,
538 $ work, ldwork )
539 END IF
540*
541* col1_(3) Compute W1: = T * W1,
542* T is upper-triangular,
543* W1 is upper-triangular with zeroes below the diagonal.
544*
545 CALL ztrmm( 'L', 'U', 'N', 'N', k, k, cone, t, ldt,
546 $ work, ldwork )
547*
548* col1_(4) Compute B1: = - V2 * W1 = - B1 * W1,
549* V2 = B1, W1 is upper-triangular with zeroes below the diagonal.
550*
551 IF( m.GT.0 ) THEN
552 CALL ztrmm( 'R', 'U', 'N', 'N', m, k, -cone, work, ldwork,
553 $ b, ldb )
554 END IF
555*
556 IF( lnotident ) THEN
557*
558* col1_(5) Compute W1: = V1 * W1 = A1 * W1,
559* V1 is not an identity matrix, but unit lower-triangular
560* V1 stored in A1 (diagonal ones are not stored),
561* W1 is upper-triangular on input with zeroes below the diagonal,
562* and square on output.
563*
564 CALL ztrmm( 'L', 'L', 'N', 'U', k, k, cone, a, lda,
565 $ work, ldwork )
566*
567* col1_(6) Compute A1: = A1 - W1 = A(1:K, 1:K) - WORK(1:K, 1:K)
568* column-by-column. A1 is upper-triangular on input.
569* If IDENT, A1 is square on output, and W1 is square,
570* if NOT IDENT, A1 is upper-triangular on output,
571* W1 is upper-triangular.
572*
573* col1_(6)_a Compute elements of A1 below the diagonal.
574*
575 DO j = 1, k - 1
576 DO i = j + 1, k
577 a( i, j ) = - work( i, j )
578 END DO
579 END DO
580*
581 END IF
582*
583* col1_(6)_b Compute elements of A1 on and above the diagonal.
584*
585 DO j = 1, k
586 DO i = 1, j
587 a( i, j ) = a( i, j ) - work( i, j )
588 END DO
589 END DO
590*
591 RETURN
592*
593* End of ZLARFB_GETT
594*
595 END
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:188
subroutine zlarfb_gett(ident, m, n, k, t, ldt, a, lda, b, ldb, work, ldwork)
ZLARFB_GETT
subroutine ztrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRMM
Definition ztrmm.f:177