1*> \brief \b CLAVHE_ROOK
2*
3*  =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6*            http://www.netlib.org/lapack/explore-html/
7*
8*  Definition:
9*  ===========
10*
11*       SUBROUTINE CLAVHE_ROOK( UPLO, TRANS, DIAG, N, NRHS, A, LDA, IPIV, B,
12*                               LDB, INFO )
13*
14*       .. Scalar Arguments ..
15*       CHARACTER          DIAG, TRANS, UPLO
16*       INTEGER            INFO, LDA, LDB, N, NRHS
17*       ..
18*       .. Array Arguments ..
19*       INTEGER            IPIV( * )
20*       COMPLEX            A( LDA, * ), B( LDB, * )
21*       ..
22*
23*
24*> \par Purpose:
25*  =============
26*>
27*> \verbatim
28*>
29*> CLAVHE_ROOK performs one of the matrix-vector operations
30*>    x := A*x  or  x := A^H*x,
31*> where x is an N element vector and  A is one of the factors
32*> from the block U*D*U' or L*D*L' factorization computed by CHETRF_ROOK.
33*>
34*> If TRANS = 'N', multiplies by U  or U * D  (or L  or L * D)
35*> If TRANS = 'C', multiplies by U' or D * U' (or L' or D * L')
36*> \endverbatim
37*
38*  Arguments:
39*  ==========
40*
41*> \param[in] UPLO
42*> \verbatim
43*>          UPLO is CHARACTER*1
44*>          Specifies whether the factor stored in A is upper or lower
45*>          triangular.
46*>          = 'U':  Upper triangular
47*>          = 'L':  Lower triangular
48*> \endverbatim
49*>
50*> \param[in] TRANS
51*> \verbatim
52*>          TRANS is CHARACTER*1
53*>          Specifies the operation to be performed:
54*>          = 'N':  x := A*x
55*>          = 'C':   x := A^H*x
56*> \endverbatim
57*>
58*> \param[in] DIAG
59*> \verbatim
60*>          DIAG is CHARACTER*1
61*>          Specifies whether or not the diagonal blocks are unit
62*>          matrices.  If the diagonal blocks are assumed to be unit,
63*>          then A = U or A = L, otherwise A = U*D or A = L*D.
64*>          = 'U':  Diagonal blocks are assumed to be unit matrices.
65*>          = 'N':  Diagonal blocks are assumed to be non-unit matrices.
66*> \endverbatim
67*>
68*> \param[in] N
69*> \verbatim
70*>          N is INTEGER
71*>          The number of rows and columns of the matrix A.  N >= 0.
72*> \endverbatim
73*>
74*> \param[in] NRHS
75*> \verbatim
76*>          NRHS is INTEGER
77*>          The number of right hand sides, i.e., the number of vectors
78*>          x to be multiplied by A.  NRHS >= 0.
79*> \endverbatim
80*>
81*> \param[in] A
82*> \verbatim
83*>          A is COMPLEX array, dimension (LDA,N)
84*>          The block diagonal matrix D and the multipliers used to
85*>          obtain the factor U or L as computed by CHETRF_ROOK.
86*>          Stored as a 2-D triangular matrix.
87*> \endverbatim
88*>
89*> \param[in] LDA
90*> \verbatim
91*>          LDA is INTEGER
92*>          The leading dimension of the array A.  LDA >= max(1,N).
93*> \endverbatim
94*>
95*> \param[out] IPIV
96*> \verbatim
97*>          IPIV is INTEGER array, dimension (N)
98*>          Details of the interchanges and the block structure of D,
99*>          as determined by CHETRF_ROOK.
100*>          If UPLO = 'U':
101*>             Only the last KB elements of IPIV are set.
102*>
103*>             If IPIV(k) > 0, then rows and columns k and IPIV(k) were
104*>             interchanged and D(k,k) is a 1-by-1 diagonal block.
105*>
106*>             If IPIV(k) < 0 and IPIV(k-1) < 0, then rows and
107*>             columns k and -IPIV(k) were interchanged and rows and
108*>             columns k-1 and -IPIV(k-1) were inerchaged,
109*>             D(k-1:k,k-1:k) is a 2-by-2 diagonal block.
110*>
111*>          If UPLO = 'L':
112*>             Only the first KB elements of IPIV are set.
113*>
114*>             If IPIV(k) > 0, then rows and columns k and IPIV(k)
115*>             were interchanged and D(k,k) is a 1-by-1 diagonal block.
116*>
117*>             If IPIV(k) < 0 and IPIV(k+1) < 0, then rows and
118*>             columns k and -IPIV(k) were interchanged and rows and
119*>             columns k+1 and -IPIV(k+1) were inerchaged,
120*>             D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
121*> \endverbatim
122*>
123*> \param[in,out] B
124*> \verbatim
125*>          B is COMPLEX array, dimension (LDB,NRHS)
126*>          On entry, B contains NRHS vectors of length N.
127*>          On exit, B is overwritten with the product A * B.
128*> \endverbatim
129*>
130*> \param[in] LDB
131*> \verbatim
132*>          LDB is INTEGER
133*>          The leading dimension of the array B.  LDB >= max(1,N).
134*> \endverbatim
135*>
136*> \param[out] INFO
137*> \verbatim
138*>          INFO is INTEGER
139*>          = 0: successful exit
140*>          < 0: if INFO = -k, the k-th argument had an illegal value
141*> \endverbatim
142*
143*  Authors:
144*  ========
145*
146*> \author Univ. of Tennessee
147*> \author Univ. of California Berkeley
148*> \author Univ. of Colorado Denver
149*> \author NAG Ltd.
150*
151*> \ingroup complex_lin
152*
153*  =====================================================================
154      SUBROUTINE CLAVHE_ROOK( UPLO, TRANS, DIAG, N, NRHS, A, LDA, IPIV,
155     $                        B, LDB, INFO )
156*
157*  -- LAPACK test routine --
158*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
159*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
160*
161*     .. Scalar Arguments ..
162      CHARACTER          DIAG, TRANS, UPLO
163      INTEGER            INFO, LDA, LDB, N, NRHS
164*     ..
165*     .. Array Arguments ..
166      INTEGER            IPIV( * )
167      COMPLEX            A( LDA, * ), B( LDB, * )
168*     ..
169*
170*  =====================================================================
171*
172*     .. Parameters ..
173      COMPLEX            CONE
174      PARAMETER          ( CONE = ( 1.0E+0, 0.0E+0 ) )
175*     ..
176*     .. Local Scalars ..
177      LOGICAL            NOUNIT
178      INTEGER            J, K, KP
179      COMPLEX            D11, D12, D21, D22, T1, T2
180*     ..
181*     .. External Functions ..
182      LOGICAL            LSAME
183      EXTERNAL           LSAME
184*     ..
185*     .. External Subroutines ..
186      EXTERNAL           CGEMV, CGERU, CLACGV, CSCAL, CSWAP, XERBLA
187*     ..
188*     .. Intrinsic Functions ..
189      INTRINSIC          ABS, CONJG, MAX
190*     ..
191*     .. Executable Statements ..
192*
193*     Test the input parameters.
194*
195      INFO = 0
196      IF( .NOT.LSAME( UPLO, 'U' ) .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
197         INFO = -1
198      ELSE IF( .NOT.LSAME( TRANS, 'N' ) .AND. .NOT.LSAME( TRANS, 'C' ) )
199     $          THEN
200         INFO = -2
201      ELSE IF( .NOT.LSAME( DIAG, 'U' ) .AND. .NOT.LSAME( DIAG, 'N' ) )
202     $          THEN
203         INFO = -3
204      ELSE IF( N.LT.0 ) THEN
205         INFO = -4
206      ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
207         INFO = -6
208      ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
209         INFO = -9
210      END IF
211      IF( INFO.NE.0 ) THEN
212         CALL XERBLA( 'CLAVHE_ROOK ', -INFO )
213         RETURN
214      END IF
215*
216*     Quick return if possible.
217*
218      IF( N.EQ.0 )
219     $   RETURN
220*
221      NOUNIT = LSAME( DIAG, 'N' )
222*------------------------------------------
223*
224*     Compute  B := A * B  (No transpose)
225*
226*------------------------------------------
227      IF( LSAME( TRANS, 'N' ) ) THEN
228*
229*        Compute  B := U*B
230*        where U = P(m)*inv(U(m))* ... *P(1)*inv(U(1))
231*
232         IF( LSAME( UPLO, 'U' ) ) THEN
233*
234*        Loop forward applying the transformations.
235*
236            K = 1
237   10       CONTINUE
238            IF( K.GT.N )
239     $         GO TO 30
240            IF( IPIV( K ).GT.0 ) THEN
241*
242*              1 x 1 pivot block
243*
244*              Multiply by the diagonal element if forming U * D.
245*
246               IF( NOUNIT )
247     $            CALL CSCAL( NRHS, A( K, K ), B( K, 1 ), LDB )
248*
249*              Multiply by  P(K) * inv(U(K))  if K > 1.
250*
251               IF( K.GT.1 ) THEN
252*
253*                 Apply the transformation.
254*
255                  CALL CGERU( K-1, NRHS, CONE, A( 1, K ), 1, B( K, 1 ),
256     $                        LDB, B( 1, 1 ), LDB )
257*
258*                 Interchange if P(K) != I.
259*
260                  KP = IPIV( K )
261                  IF( KP.NE.K )
262     $               CALL CSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
263               END IF
264               K = K + 1
265            ELSE
266*
267*              2 x 2 pivot block
268*
269*              Multiply by the diagonal block if forming U * D.
270*
271               IF( NOUNIT ) THEN
272                  D11 = A( K, K )
273                  D22 = A( K+1, K+1 )
274                  D12 = A( K, K+1 )
275                  D21 = CONJG( D12 )
276                  DO 20 J = 1, NRHS
277                     T1 = B( K, J )
278                     T2 = B( K+1, J )
279                     B( K, J ) = D11*T1 + D12*T2
280                     B( K+1, J ) = D21*T1 + D22*T2
281   20             CONTINUE
282               END IF
283*
284*              Multiply by  P(K) * inv(U(K))  if K > 1.
285*
286               IF( K.GT.1 ) THEN
287*
288*                 Apply the transformations.
289*
290                  CALL CGERU( K-1, NRHS, CONE, A( 1, K ), 1, B( K, 1 ),
291     $                        LDB, B( 1, 1 ), LDB )
292                  CALL CGERU( K-1, NRHS, CONE, A( 1, K+1 ), 1,
293     $                        B( K+1, 1 ), LDB, B( 1, 1 ), LDB )
294*
295*                 Interchange if a permutation was applied at the
296*                 K-th step of the factorization.
297*
298*                 Swap the first of pair with IMAXth
299*
300                  KP = ABS( IPIV( K ) )
301                  IF( KP.NE.K )
302     $               CALL CSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
303*
304*                 NOW swap the first of pair with Pth
305*
306                  KP = ABS( IPIV( K+1 ) )
307                  IF( KP.NE.K+1 )
308     $               CALL CSWAP( NRHS, B( K+1, 1 ), LDB, B( KP, 1 ),
309     $                           LDB )
310               END IF
311               K = K + 2
312            END IF
313            GO TO 10
314   30       CONTINUE
315*
316*        Compute  B := L*B
317*        where L = P(1)*inv(L(1))* ... *P(m)*inv(L(m)) .
318*
319         ELSE
320*
321*           Loop backward applying the transformations to B.
322*
323            K = N
324   40       CONTINUE
325            IF( K.LT.1 )
326     $         GO TO 60
327*
328*           Test the pivot index.  If greater than zero, a 1 x 1
329*           pivot was used, otherwise a 2 x 2 pivot was used.
330*
331            IF( IPIV( K ).GT.0 ) THEN
332*
333*              1 x 1 pivot block:
334*
335*              Multiply by the diagonal element if forming L * D.
336*
337               IF( NOUNIT )
338     $            CALL CSCAL( NRHS, A( K, K ), B( K, 1 ), LDB )
339*
340*              Multiply by  P(K) * inv(L(K))  if K < N.
341*
342               IF( K.NE.N ) THEN
343                  KP = IPIV( K )
344*
345*                 Apply the transformation.
346*
347                  CALL CGERU( N-K, NRHS, CONE, A( K+1, K ), 1,
348     $                        B( K, 1 ), LDB, B( K+1, 1 ), LDB )
349*
350*                 Interchange if a permutation was applied at the
351*                 K-th step of the factorization.
352*
353                  IF( KP.NE.K )
354     $               CALL CSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
355               END IF
356               K = K - 1
357*
358            ELSE
359*
360*              2 x 2 pivot block:
361*
362*              Multiply by the diagonal block if forming L * D.
363*
364               IF( NOUNIT ) THEN
365                  D11 = A( K-1, K-1 )
366                  D22 = A( K, K )
367                  D21 = A( K, K-1 )
368                  D12 = CONJG( D21 )
369                  DO 50 J = 1, NRHS
370                     T1 = B( K-1, J )
371                     T2 = B( K, J )
372                     B( K-1, J ) = D11*T1 + D12*T2
373                     B( K, J ) = D21*T1 + D22*T2
374   50             CONTINUE
375               END IF
376*
377*              Multiply by  P(K) * inv(L(K))  if K < N.
378*
379               IF( K.NE.N ) THEN
380*
381*                 Apply the transformation.
382*
383                  CALL CGERU( N-K, NRHS, CONE, A( K+1, K ), 1,
384     $                        B( K, 1 ), LDB, B( K+1, 1 ), LDB )
385                  CALL CGERU( N-K, NRHS, CONE, A( K+1, K-1 ), 1,
386     $                        B( K-1, 1 ), LDB, B( K+1, 1 ), LDB )
387*
388*                 Interchange if a permutation was applied at the
389*                 K-th step of the factorization.
390*
391*
392*                 Swap the second of pair with IMAXth
393*
394                  KP = ABS( IPIV( K ) )
395                  IF( KP.NE.K )
396     $               CALL CSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
397*
398*                 NOW swap the first of pair with Pth
399*
400                  KP = ABS( IPIV( K-1 ) )
401                  IF( KP.NE.K-1 )
402     $               CALL CSWAP( NRHS, B( K-1, 1 ), LDB, B( KP, 1 ),
403     $                           LDB )
404*
405               END IF
406               K = K - 2
407            END IF
408            GO TO 40
409   60       CONTINUE
410         END IF
411*--------------------------------------------------
412*
413*     Compute  B := A^H * B  (conjugate transpose)
414*
415*--------------------------------------------------
416      ELSE
417*
418*        Form  B := U^H*B
419*        where U  = P(m)*inv(U(m))* ... *P(1)*inv(U(1))
420*        and   U^H = inv(U^H(1))*P(1)* ... *inv(U^H(m))*P(m)
421*
422         IF( LSAME( UPLO, 'U' ) ) THEN
423*
424*           Loop backward applying the transformations.
425*
426            K = N
427   70       IF( K.LT.1 )
428     $         GO TO 90
429*
430*           1 x 1 pivot block.
431*
432            IF( IPIV( K ).GT.0 ) THEN
433               IF( K.GT.1 ) THEN
434*
435*                 Interchange if P(K) != I.
436*
437                  KP = IPIV( K )
438                  IF( KP.NE.K )
439     $               CALL CSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
440*
441*                 Apply the transformation
442*                    y = y - B' conjg(x),
443*                 where x is a column of A and y is a row of B.
444*
445                  CALL CLACGV( NRHS, B( K, 1 ), LDB )
446                  CALL CGEMV( 'Conjugate', K-1, NRHS, CONE, B, LDB,
447     $                        A( 1, K ), 1, CONE, B( K, 1 ), LDB )
448                  CALL CLACGV( NRHS, B( K, 1 ), LDB )
449               END IF
450               IF( NOUNIT )
451     $            CALL CSCAL( NRHS, A( K, K ), B( K, 1 ), LDB )
452               K = K - 1
453*
454*           2 x 2 pivot block.
455*
456            ELSE
457               IF( K.GT.2 ) THEN
458*
459*                 Swap the second of pair with Pth
460*
461                  KP = ABS( IPIV( K ) )
462                  IF( KP.NE.K )
463     $               CALL CSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
464*
465*                 Now swap the first of pair with IMAX(r)th
466*
467                  KP = ABS( IPIV( K-1 ) )
468                  IF( KP.NE.K-1 )
469     $               CALL CSWAP( NRHS, B( K-1, 1 ), LDB, B( KP, 1 ),
470     $                           LDB )
471*
472*                 Apply the transformations
473*                    y = y - B' conjg(x),
474*                 where x is a block column of A and y is a block
475*                 row of B.
476*
477                  CALL CLACGV( NRHS, B( K, 1 ), LDB )
478                  CALL CGEMV( 'Conjugate', K-2, NRHS, CONE, B, LDB,
479     $                        A( 1, K ), 1, CONE, B( K, 1 ), LDB )
480                  CALL CLACGV( NRHS, B( K, 1 ), LDB )
481*
482                  CALL CLACGV( NRHS, B( K-1, 1 ), LDB )
483                  CALL CGEMV( 'Conjugate', K-2, NRHS, CONE, B, LDB,
484     $                        A( 1, K-1 ), 1, CONE, B( K-1, 1 ), LDB )
485                  CALL CLACGV( NRHS, B( K-1, 1 ), LDB )
486               END IF
487*
488*              Multiply by the diagonal block if non-unit.
489*
490               IF( NOUNIT ) THEN
491                  D11 = A( K-1, K-1 )
492                  D22 = A( K, K )
493                  D12 = A( K-1, K )
494                  D21 = CONJG( D12 )
495                  DO 80 J = 1, NRHS
496                     T1 = B( K-1, J )
497                     T2 = B( K, J )
498                     B( K-1, J ) = D11*T1 + D12*T2
499                     B( K, J ) = D21*T1 + D22*T2
500   80             CONTINUE
501               END IF
502               K = K - 2
503            END IF
504            GO TO 70
505   90       CONTINUE
506*
507*        Form  B := L^H*B
508*        where L  = P(1)*inv(L(1))* ... *P(m)*inv(L(m))
509*        and   L^H = inv(L^H(m))*P(m)* ... *inv(L^H(1))*P(1)
510*
511         ELSE
512*
513*           Loop forward applying the L-transformations.
514*
515            K = 1
516  100       CONTINUE
517            IF( K.GT.N )
518     $         GO TO 120
519*
520*           1 x 1 pivot block
521*
522            IF( IPIV( K ).GT.0 ) THEN
523               IF( K.LT.N ) THEN
524*
525*                 Interchange if P(K) != I.
526*
527                  KP = IPIV( K )
528                  IF( KP.NE.K )
529     $               CALL CSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
530*
531*                 Apply the transformation
532*
533                  CALL CLACGV( NRHS, B( K, 1 ), LDB )
534                  CALL CGEMV( 'Conjugate', N-K, NRHS, CONE, B( K+1, 1 ),
535     $                       LDB, A( K+1, K ), 1, CONE, B( K, 1 ), LDB )
536                  CALL CLACGV( NRHS, B( K, 1 ), LDB )
537               END IF
538               IF( NOUNIT )
539     $            CALL CSCAL( NRHS, A( K, K ), B( K, 1 ), LDB )
540               K = K + 1
541*
542*           2 x 2 pivot block.
543*
544            ELSE
545               IF( K.LT.N-1 ) THEN
546*
547*                 Swap the first of pair with Pth
548*
549                  KP = ABS( IPIV( K ) )
550                  IF( KP.NE.K )
551     $               CALL CSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
552*
553*                 Now swap the second of pair with IMAX(r)th
554*
555                  KP = ABS( IPIV( K+1 ) )
556                  IF( KP.NE.K+1 )
557     $               CALL CSWAP( NRHS, B( K+1, 1 ), LDB, B( KP, 1 ),
558     $                           LDB )
559*
560*                 Apply the transformation
561*
562                  CALL CLACGV( NRHS, B( K+1, 1 ), LDB )
563                  CALL CGEMV( 'Conjugate', N-K-1, NRHS, CONE,
564     $                        B( K+2, 1 ), LDB, A( K+2, K+1 ), 1, CONE,
565     $                        B( K+1, 1 ), LDB )
566                  CALL CLACGV( NRHS, B( K+1, 1 ), LDB )
567*
568                  CALL CLACGV( NRHS, B( K, 1 ), LDB )
569                  CALL CGEMV( 'Conjugate', N-K-1, NRHS, CONE,
570     $                        B( K+2, 1 ), LDB, A( K+2, K ), 1, CONE,
571     $                        B( K, 1 ), LDB )
572                  CALL CLACGV( NRHS, B( K, 1 ), LDB )
573               END IF
574*
575*              Multiply by the diagonal block if non-unit.
576*
577               IF( NOUNIT ) THEN
578                  D11 = A( K, K )
579                  D22 = A( K+1, K+1 )
580                  D21 = A( K+1, K )
581                  D12 = CONJG( D21 )
582                  DO 110 J = 1, NRHS
583                     T1 = B( K, J )
584                     T2 = B( K+1, J )
585                     B( K, J ) = D11*T1 + D12*T2
586                     B( K+1, J ) = D21*T1 + D22*T2
587  110             CONTINUE
588               END IF
589               K = K + 2
590            END IF
591            GO TO 100
592  120       CONTINUE
593         END IF
594*
595      END IF
596      RETURN
597*
598*     End of CLAVHE_ROOK
599*
600      END
601