1*> \brief \b ZLAVSY_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 ZLAVSY_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*16         A( LDA, * ), B( LDB, * )
21*       ..
22*
23*
24*> \par Purpose:
25*  =============
26*>
27*> \verbatim
28*>
29*> ZLAVSY_ROOK performs one of the matrix-vector operations
30*>    x := A*x  or  x := A'*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 ZSYTRF_ROOK.
33*>
34*> If TRANS = 'N', multiplies by U  or U * D  (or L  or L * D)
35*> If TRANS = 'T', 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*>          = 'T':  x := A'*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*16 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 ZSYTRF_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[in] 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 ZSYTRF_ROOK.
100*>
101*>          If UPLO = 'U':
102*>               If IPIV(k) > 0, then rows and columns k and IPIV(k)
103*>               were interchanged and D(k,k) is a 1-by-1 diagonal block.
104*>               (If IPIV( k ) = k, no interchange was done).
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*>               If IPIV(k) > 0, then rows and columns k and IPIV(k)
113*>               were interchanged and D(k,k) is a 1-by-1 diagonal block.
114*>               (If IPIV( k ) = k, no interchange was done).
115*>
116*>               If IPIV(k) < 0 and IPIV(k+1) < 0, then rows and
117*>               columns k and -IPIV(k) were interchanged and rows and
118*>               columns k+1 and -IPIV(k+1) were inerchaged,
119*>               D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
120*> \endverbatim
121*>
122*> \param[in,out] B
123*> \verbatim
124*>          B is COMPLEX*16 array, dimension (LDB,NRHS)
125*>          On entry, B contains NRHS vectors of length N.
126*>          On exit, B is overwritten with the product A * B.
127*> \endverbatim
128*>
129*> \param[in] LDB
130*> \verbatim
131*>          LDB is INTEGER
132*>          The leading dimension of the array B.  LDB >= max(1,N).
133*> \endverbatim
134*>
135*> \param[out] INFO
136*> \verbatim
137*>          INFO is INTEGER
138*>          = 0: successful exit
139*>          < 0: if INFO = -k, the k-th argument had an illegal value
140*> \endverbatim
141*
142*  Authors:
143*  ========
144*
145*> \author Univ. of Tennessee
146*> \author Univ. of California Berkeley
147*> \author Univ. of Colorado Denver
148*> \author NAG Ltd.
149*
150*> \ingroup complex16_lin
151*
152*  =====================================================================
153      SUBROUTINE ZLAVSY_ROOK( UPLO, TRANS, DIAG, N, NRHS, A, LDA, IPIV,
154     $                        B, LDB, INFO )
155*
156*  -- LAPACK test routine --
157*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
158*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
159*
160*     .. Scalar Arguments ..
161      CHARACTER          DIAG, TRANS, UPLO
162      INTEGER            INFO, LDA, LDB, N, NRHS
163*     ..
164*     .. Array Arguments ..
165      INTEGER            IPIV( * )
166      COMPLEX*16         A( LDA, * ), B( LDB, * )
167*     ..
168*
169*  =====================================================================
170*
171*     .. Parameters ..
172      COMPLEX*16         CONE
173      PARAMETER          ( CONE = ( 1.0D+0, 0.0D+0 ) )
174*     ..
175*     .. Local Scalars ..
176      LOGICAL            NOUNIT
177      INTEGER            J, K, KP
178      COMPLEX*16         D11, D12, D21, D22, T1, T2
179*     ..
180*     .. External Functions ..
181      LOGICAL            LSAME
182      EXTERNAL           LSAME
183*     ..
184*     .. External Subroutines ..
185      EXTERNAL           XERBLA, ZGEMV, ZGERU, ZSCAL, ZSWAP
186*     ..
187*     .. Intrinsic Functions ..
188      INTRINSIC          ABS, MAX
189*     ..
190*     .. Executable Statements ..
191*
192*     Test the input parameters.
193*
194      INFO = 0
195      IF( .NOT.LSAME( UPLO, 'U' ) .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
196         INFO = -1
197      ELSE IF( .NOT.LSAME( TRANS, 'N' ) .AND. .NOT.LSAME( TRANS, 'T' ) )
198     $          THEN
199         INFO = -2
200      ELSE IF( .NOT.LSAME( DIAG, 'U' ) .AND. .NOT.LSAME( DIAG, 'N' ) )
201     $          THEN
202         INFO = -3
203      ELSE IF( N.LT.0 ) THEN
204         INFO = -4
205      ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
206         INFO = -6
207      ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
208         INFO = -9
209      END IF
210      IF( INFO.NE.0 ) THEN
211         CALL XERBLA( 'ZLAVSY_ROOK ', -INFO )
212         RETURN
213      END IF
214*
215*     Quick return if possible.
216*
217      IF( N.EQ.0 )
218     $   RETURN
219*
220      NOUNIT = LSAME( DIAG, 'N' )
221*------------------------------------------
222*
223*     Compute  B := A * B  (No transpose)
224*
225*------------------------------------------
226      IF( LSAME( TRANS, 'N' ) ) THEN
227*
228*        Compute  B := U*B
229*        where U = P(m)*inv(U(m))* ... *P(1)*inv(U(1))
230*
231         IF( LSAME( UPLO, 'U' ) ) THEN
232*
233*        Loop forward applying the transformations.
234*
235            K = 1
236   10       CONTINUE
237            IF( K.GT.N )
238     $         GO TO 30
239            IF( IPIV( K ).GT.0 ) THEN
240*
241*              1 x 1 pivot block
242*
243*              Multiply by the diagonal element if forming U * D.
244*
245               IF( NOUNIT )
246     $            CALL ZSCAL( NRHS, A( K, K ), B( K, 1 ), LDB )
247*
248*              Multiply by  P(K) * inv(U(K))  if K > 1.
249*
250               IF( K.GT.1 ) THEN
251*
252*                 Apply the transformation.
253*
254                  CALL ZGERU( K-1, NRHS, CONE, A( 1, K ), 1, B( K, 1 ),
255     $                        LDB, B( 1, 1 ), LDB )
256*
257*                 Interchange if P(K) != I.
258*
259                  KP = IPIV( K )
260                  IF( KP.NE.K )
261     $               CALL ZSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
262               END IF
263               K = K + 1
264            ELSE
265*
266*              2 x 2 pivot block
267*
268*              Multiply by the diagonal block if forming U * D.
269*
270               IF( NOUNIT ) THEN
271                  D11 = A( K, K )
272                  D22 = A( K+1, K+1 )
273                  D12 = A( K, K+1 )
274                  D21 = D12
275                  DO 20 J = 1, NRHS
276                     T1 = B( K, J )
277                     T2 = B( K+1, J )
278                     B( K, J ) = D11*T1 + D12*T2
279                     B( K+1, J ) = D21*T1 + D22*T2
280   20             CONTINUE
281               END IF
282*
283*              Multiply by  P(K) * inv(U(K))  if K > 1.
284*
285               IF( K.GT.1 ) THEN
286*
287*                 Apply the transformations.
288*
289                  CALL ZGERU( K-1, NRHS, CONE, A( 1, K ), 1, B( K, 1 ),
290     $                        LDB, B( 1, 1 ), LDB )
291                  CALL ZGERU( K-1, NRHS, CONE, A( 1, K+1 ), 1,
292     $                        B( K+1, 1 ), LDB, B( 1, 1 ), LDB )
293*
294*                 Interchange if a permutation was applied at the
295*                 K-th step of the factorization.
296*
297*                 Swap the first of pair with IMAXth
298*
299                  KP = ABS( IPIV( K ) )
300                  IF( KP.NE.K )
301     $               CALL ZSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
302*
303*                 NOW swap the first of pair with Pth
304*
305                  KP = ABS( IPIV( K+1 ) )
306                  IF( KP.NE.K+1 )
307     $               CALL ZSWAP( NRHS, B( K+1, 1 ), LDB, B( KP, 1 ),
308     $                           LDB )
309               END IF
310               K = K + 2
311            END IF
312            GO TO 10
313   30       CONTINUE
314*
315*        Compute  B := L*B
316*        where L = P(1)*inv(L(1))* ... *P(m)*inv(L(m)) .
317*
318         ELSE
319*
320*           Loop backward applying the transformations to B.
321*
322            K = N
323   40       CONTINUE
324            IF( K.LT.1 )
325     $         GO TO 60
326*
327*           Test the pivot index.  If greater than zero, a 1 x 1
328*           pivot was used, otherwise a 2 x 2 pivot was used.
329*
330            IF( IPIV( K ).GT.0 ) THEN
331*
332*              1 x 1 pivot block:
333*
334*              Multiply by the diagonal element if forming L * D.
335*
336               IF( NOUNIT )
337     $            CALL ZSCAL( NRHS, A( K, K ), B( K, 1 ), LDB )
338*
339*              Multiply by  P(K) * inv(L(K))  if K < N.
340*
341               IF( K.NE.N ) THEN
342                  KP = IPIV( K )
343*
344*                 Apply the transformation.
345*
346                  CALL ZGERU( N-K, NRHS, CONE, A( K+1, K ), 1,
347     $                        B( K, 1 ), LDB, B( K+1, 1 ), LDB )
348*
349*                 Interchange if a permutation was applied at the
350*                 K-th step of the factorization.
351*
352                  IF( KP.NE.K )
353     $               CALL ZSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
354               END IF
355               K = K - 1
356*
357            ELSE
358*
359*              2 x 2 pivot block:
360*
361*              Multiply by the diagonal block if forming L * D.
362*
363               IF( NOUNIT ) THEN
364                  D11 = A( K-1, K-1 )
365                  D22 = A( K, K )
366                  D21 = A( K, K-1 )
367                  D12 = D21
368                  DO 50 J = 1, NRHS
369                     T1 = B( K-1, J )
370                     T2 = B( K, J )
371                     B( K-1, J ) = D11*T1 + D12*T2
372                     B( K, J ) = D21*T1 + D22*T2
373   50             CONTINUE
374               END IF
375*
376*              Multiply by  P(K) * inv(L(K))  if K < N.
377*
378               IF( K.NE.N ) THEN
379*
380*                 Apply the transformation.
381*
382                  CALL ZGERU( N-K, NRHS, CONE, A( K+1, K ), 1,
383     $                        B( K, 1 ), LDB, B( K+1, 1 ), LDB )
384                  CALL ZGERU( N-K, NRHS, CONE, A( K+1, K-1 ), 1,
385     $                        B( K-1, 1 ), LDB, B( K+1, 1 ), LDB )
386*
387*                 Interchange if a permutation was applied at the
388*                 K-th step of the factorization.
389*
390*                 Swap the second of pair with IMAXth
391*
392                  KP = ABS( IPIV( K ) )
393                  IF( KP.NE.K )
394     $               CALL ZSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
395*
396*                 NOW swap the first of pair with Pth
397*
398                  KP = ABS( IPIV( K-1 ) )
399                  IF( KP.NE.K-1 )
400     $               CALL ZSWAP( NRHS, B( K-1, 1 ), LDB, B( KP, 1 ),
401     $                           LDB )
402               END IF
403               K = K - 2
404            END IF
405            GO TO 40
406   60       CONTINUE
407         END IF
408*----------------------------------------
409*
410*     Compute  B := A' * B  (transpose)
411*
412*----------------------------------------
413      ELSE IF( LSAME( TRANS, 'T' ) ) THEN
414*
415*        Form  B := U'*B
416*        where U  = P(m)*inv(U(m))* ... *P(1)*inv(U(1))
417*        and   U' = inv(U'(1))*P(1)* ... *inv(U'(m))*P(m)
418*
419         IF( LSAME( UPLO, 'U' ) ) THEN
420*
421*           Loop backward applying the transformations.
422*
423            K = N
424   70       CONTINUE
425            IF( K.LT.1 )
426     $         GO TO 90
427*
428*           1 x 1 pivot block.
429*
430            IF( IPIV( K ).GT.0 ) THEN
431               IF( K.GT.1 ) THEN
432*
433*                 Interchange if P(K) != I.
434*
435                  KP = IPIV( K )
436                  IF( KP.NE.K )
437     $               CALL ZSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
438*
439*                 Apply the transformation
440*
441                  CALL ZGEMV( 'Transpose', K-1, NRHS, CONE, B, LDB,
442     $                        A( 1, K ), 1, CONE, B( K, 1 ), LDB )
443               END IF
444               IF( NOUNIT )
445     $            CALL ZSCAL( NRHS, A( K, K ), B( K, 1 ), LDB )
446               K = K - 1
447*
448*           2 x 2 pivot block.
449*
450            ELSE
451               IF( K.GT.2 ) THEN
452*
453*                 Swap the second of pair with Pth
454*
455                  KP = ABS( IPIV( K ) )
456                  IF( KP.NE.K )
457     $               CALL ZSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
458*
459*                 Now swap the first of pair with IMAX(r)th
460*
461                  KP = ABS( IPIV( K-1 ) )
462                  IF( KP.NE.K-1 )
463     $               CALL ZSWAP( NRHS, B( K-1, 1 ), LDB, B( KP, 1 ),
464     $                           LDB )
465*
466*                 Apply the transformations
467*
468                  CALL ZGEMV( 'Transpose', K-2, NRHS, CONE, B, LDB,
469     $                        A( 1, K ), 1, CONE, B( K, 1 ), LDB )
470                  CALL ZGEMV( 'Transpose', K-2, NRHS, CONE, B, LDB,
471     $                        A( 1, K-1 ), 1, CONE, B( K-1, 1 ), LDB )
472               END IF
473*
474*              Multiply by the diagonal block if non-unit.
475*
476               IF( NOUNIT ) THEN
477                  D11 = A( K-1, K-1 )
478                  D22 = A( K, K )
479                  D12 = A( K-1, K )
480                  D21 = D12
481                  DO 80 J = 1, NRHS
482                     T1 = B( K-1, J )
483                     T2 = B( K, J )
484                     B( K-1, J ) = D11*T1 + D12*T2
485                     B( K, J ) = D21*T1 + D22*T2
486   80             CONTINUE
487               END IF
488               K = K - 2
489            END IF
490            GO TO 70
491   90       CONTINUE
492*
493*        Form  B := L'*B
494*        where L  = P(1)*inv(L(1))* ... *P(m)*inv(L(m))
495*        and   L' = inv(L'(m))*P(m)* ... *inv(L'(1))*P(1)
496*
497         ELSE
498*
499*           Loop forward applying the L-transformations.
500*
501            K = 1
502  100       CONTINUE
503            IF( K.GT.N )
504     $         GO TO 120
505*
506*           1 x 1 pivot block
507*
508            IF( IPIV( K ).GT.0 ) THEN
509               IF( K.LT.N ) THEN
510*
511*                 Interchange if P(K) != I.
512*
513                  KP = IPIV( K )
514                  IF( KP.NE.K )
515     $               CALL ZSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
516*
517*                 Apply the transformation
518*
519                  CALL ZGEMV( 'Transpose', N-K, NRHS, CONE, B( K+1, 1 ),
520     $                       LDB, A( K+1, K ), 1, CONE, B( K, 1 ), LDB )
521               END IF
522               IF( NOUNIT )
523     $            CALL ZSCAL( NRHS, A( K, K ), B( K, 1 ), LDB )
524               K = K + 1
525*
526*           2 x 2 pivot block.
527*
528            ELSE
529               IF( K.LT.N-1 ) THEN
530*
531*                 Swap the first of pair with Pth
532*
533                  KP = ABS( IPIV( K ) )
534                  IF( KP.NE.K )
535     $               CALL ZSWAP( NRHS, B( K, 1 ), LDB, B( KP, 1 ), LDB )
536*
537*                 Now swap the second of pair with IMAX(r)th
538*
539                  KP = ABS( IPIV( K+1 ) )
540                  IF( KP.NE.K+1 )
541     $               CALL ZSWAP( NRHS, B( K+1, 1 ), LDB, B( KP, 1 ),
542     $                           LDB )
543*
544*                 Apply the transformation
545*
546                  CALL ZGEMV( 'Transpose', N-K-1, NRHS, CONE,
547     $                        B( K+2, 1 ), LDB, A( K+2, K+1 ), 1, CONE,
548     $                        B( K+1, 1 ), LDB )
549                  CALL ZGEMV( 'Transpose', N-K-1, NRHS, CONE,
550     $                        B( K+2, 1 ), LDB, A( K+2, K ), 1, CONE,
551     $                        B( K, 1 ), LDB )
552               END IF
553*
554*              Multiply by the diagonal block if non-unit.
555*
556               IF( NOUNIT ) THEN
557                  D11 = A( K, K )
558                  D22 = A( K+1, K+1 )
559                  D21 = A( K+1, K )
560                  D12 = D21
561                  DO 110 J = 1, NRHS
562                     T1 = B( K, J )
563                     T2 = B( K+1, J )
564                     B( K, J ) = D11*T1 + D12*T2
565                     B( K+1, J ) = D21*T1 + D22*T2
566  110             CONTINUE
567               END IF
568               K = K + 2
569            END IF
570            GO TO 100
571  120       CONTINUE
572         END IF
573      END IF
574      RETURN
575*
576*     End of ZLAVSY_ROOK
577*
578      END
579