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