1*> \brief \b STBMV
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 STBMV(UPLO,TRANS,DIAG,N,K,A,LDA,X,INCX)
12*
13*       .. Scalar Arguments ..
14*       INTEGER INCX,K,LDA,N
15*       CHARACTER DIAG,TRANS,UPLO
16*       ..
17*       .. Array Arguments ..
18*       REAL A(LDA,*),X(*)
19*       ..
20*
21*
22*> \par Purpose:
23*  =============
24*>
25*> \verbatim
26*>
27*> STBMV  performs one of the matrix-vector operations
28*>
29*>    x := A*x,   or   x := A**T*x,
30*>
31*> where x is an n element vector and  A is an n by n unit, or non-unit,
32*> upper or lower triangular band matrix, with ( k + 1 ) diagonals.
33*> \endverbatim
34*
35*  Arguments:
36*  ==========
37*
38*> \param[in] UPLO
39*> \verbatim
40*>          UPLO is CHARACTER*1
41*>           On entry, UPLO specifies whether the matrix is an upper or
42*>           lower triangular matrix as follows:
43*>
44*>              UPLO = 'U' or 'u'   A is an upper triangular matrix.
45*>
46*>              UPLO = 'L' or 'l'   A is a lower triangular matrix.
47*> \endverbatim
48*>
49*> \param[in] TRANS
50*> \verbatim
51*>          TRANS is CHARACTER*1
52*>           On entry, TRANS specifies the operation to be performed as
53*>           follows:
54*>
55*>              TRANS = 'N' or 'n'   x := A*x.
56*>
57*>              TRANS = 'T' or 't'   x := A**T*x.
58*>
59*>              TRANS = 'C' or 'c'   x := A**T*x.
60*> \endverbatim
61*>
62*> \param[in] DIAG
63*> \verbatim
64*>          DIAG is CHARACTER*1
65*>           On entry, DIAG specifies whether or not A is unit
66*>           triangular as follows:
67*>
68*>              DIAG = 'U' or 'u'   A is assumed to be unit triangular.
69*>
70*>              DIAG = 'N' or 'n'   A is not assumed to be unit
71*>                                  triangular.
72*> \endverbatim
73*>
74*> \param[in] N
75*> \verbatim
76*>          N is INTEGER
77*>           On entry, N specifies the order of the matrix A.
78*>           N must be at least zero.
79*> \endverbatim
80*>
81*> \param[in] K
82*> \verbatim
83*>          K is INTEGER
84*>           On entry with UPLO = 'U' or 'u', K specifies the number of
85*>           super-diagonals of the matrix A.
86*>           On entry with UPLO = 'L' or 'l', K specifies the number of
87*>           sub-diagonals of the matrix A.
88*>           K must satisfy  0 .le. K.
89*> \endverbatim
90*>
91*> \param[in] A
92*> \verbatim
93*>          A is REAL array, dimension ( LDA, N )
94*>           Before entry with UPLO = 'U' or 'u', the leading ( k + 1 )
95*>           by n part of the array A must contain the upper triangular
96*>           band part of the matrix of coefficients, supplied column by
97*>           column, with the leading diagonal of the matrix in row
98*>           ( k + 1 ) of the array, the first super-diagonal starting at
99*>           position 2 in row k, and so on. The top left k by k triangle
100*>           of the array A is not referenced.
101*>           The following program segment will transfer an upper
102*>           triangular band matrix from conventional full matrix storage
103*>           to band storage:
104*>
105*>                 DO 20, J = 1, N
106*>                    M = K + 1 - J
107*>                    DO 10, I = MAX( 1, J - K ), J
108*>                       A( M + I, J ) = matrix( I, J )
109*>              10    CONTINUE
110*>              20 CONTINUE
111*>
112*>           Before entry with UPLO = 'L' or 'l', the leading ( k + 1 )
113*>           by n part of the array A must contain the lower triangular
114*>           band part of the matrix of coefficients, supplied column by
115*>           column, with the leading diagonal of the matrix in row 1 of
116*>           the array, the first sub-diagonal starting at position 1 in
117*>           row 2, and so on. The bottom right k by k triangle of the
118*>           array A is not referenced.
119*>           The following program segment will transfer a lower
120*>           triangular band matrix from conventional full matrix storage
121*>           to band storage:
122*>
123*>                 DO 20, J = 1, N
124*>                    M = 1 - J
125*>                    DO 10, I = J, MIN( N, J + K )
126*>                       A( M + I, J ) = matrix( I, J )
127*>              10    CONTINUE
128*>              20 CONTINUE
129*>
130*>           Note that when DIAG = 'U' or 'u' the elements of the array A
131*>           corresponding to the diagonal elements of the matrix are not
132*>           referenced, but are assumed to be unity.
133*> \endverbatim
134*>
135*> \param[in] LDA
136*> \verbatim
137*>          LDA is INTEGER
138*>           On entry, LDA specifies the first dimension of A as declared
139*>           in the calling (sub) program. LDA must be at least
140*>           ( k + 1 ).
141*> \endverbatim
142*>
143*> \param[in,out] X
144*> \verbatim
145*>          X is REAL array, dimension at least
146*>           ( 1 + ( n - 1 )*abs( INCX ) ).
147*>           Before entry, the incremented array X must contain the n
148*>           element vector x. On exit, X is overwritten with the
149*>           transformed vector x.
150*> \endverbatim
151*>
152*> \param[in] INCX
153*> \verbatim
154*>          INCX is INTEGER
155*>           On entry, INCX specifies the increment for the elements of
156*>           X. INCX must not be zero.
157*> \endverbatim
158*
159*  Authors:
160*  ========
161*
162*> \author Univ. of Tennessee
163*> \author Univ. of California Berkeley
164*> \author Univ. of Colorado Denver
165*> \author NAG Ltd.
166*
167*> \ingroup single_blas_level2
168*
169*> \par Further Details:
170*  =====================
171*>
172*> \verbatim
173*>
174*>  Level 2 Blas routine.
175*>  The vector and matrix arguments are not referenced when N = 0, or M = 0
176*>
177*>  -- Written on 22-October-1986.
178*>     Jack Dongarra, Argonne National Lab.
179*>     Jeremy Du Croz, Nag Central Office.
180*>     Sven Hammarling, Nag Central Office.
181*>     Richard Hanson, Sandia National Labs.
182*> \endverbatim
183*>
184*  =====================================================================
185      SUBROUTINE STBMV(UPLO,TRANS,DIAG,N,K,A,LDA,X,INCX)
186*
187*  -- Reference BLAS level2 routine --
188*  -- Reference BLAS is a software package provided by Univ. of Tennessee,    --
189*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
190*
191*     .. Scalar Arguments ..
192      INTEGER INCX,K,LDA,N
193      CHARACTER DIAG,TRANS,UPLO
194*     ..
195*     .. Array Arguments ..
196      REAL A(LDA,*),X(*)
197*     ..
198*
199*  =====================================================================
200*
201*     .. Parameters ..
202      REAL ZERO
203      PARAMETER (ZERO=0.0E+0)
204*     ..
205*     .. Local Scalars ..
206      REAL TEMP
207      INTEGER I,INFO,IX,J,JX,KPLUS1,KX,L
208      LOGICAL NOUNIT
209*     ..
210*     .. External Functions ..
211      LOGICAL LSAME
212      EXTERNAL LSAME
213*     ..
214*     .. External Subroutines ..
215      EXTERNAL XERBLA
216*     ..
217*     .. Intrinsic Functions ..
218      INTRINSIC MAX,MIN
219*     ..
220*
221*     Test the input parameters.
222*
223      INFO = 0
224      IF (.NOT.LSAME(UPLO,'U') .AND. .NOT.LSAME(UPLO,'L')) THEN
225          INFO = 1
226      ELSE IF (.NOT.LSAME(TRANS,'N') .AND. .NOT.LSAME(TRANS,'T') .AND.
227     +         .NOT.LSAME(TRANS,'C')) THEN
228          INFO = 2
229      ELSE IF (.NOT.LSAME(DIAG,'U') .AND. .NOT.LSAME(DIAG,'N')) THEN
230          INFO = 3
231      ELSE IF (N.LT.0) THEN
232          INFO = 4
233      ELSE IF (K.LT.0) THEN
234          INFO = 5
235      ELSE IF (LDA.LT. (K+1)) THEN
236          INFO = 7
237      ELSE IF (INCX.EQ.0) THEN
238          INFO = 9
239      END IF
240      IF (INFO.NE.0) THEN
241          CALL XERBLA('STBMV ',INFO)
242          RETURN
243      END IF
244*
245*     Quick return if possible.
246*
247      IF (N.EQ.0) RETURN
248*
249      NOUNIT = LSAME(DIAG,'N')
250*
251*     Set up the start point in X if the increment is not unity. This
252*     will be  ( N - 1 )*INCX   too small for descending loops.
253*
254      IF (INCX.LE.0) THEN
255          KX = 1 - (N-1)*INCX
256      ELSE IF (INCX.NE.1) THEN
257          KX = 1
258      END IF
259*
260*     Start the operations. In this version the elements of A are
261*     accessed sequentially with one pass through A.
262*
263      IF (LSAME(TRANS,'N')) THEN
264*
265*         Form  x := A*x.
266*
267          IF (LSAME(UPLO,'U')) THEN
268              KPLUS1 = K + 1
269              IF (INCX.EQ.1) THEN
270                  DO 20 J = 1,N
271                      IF (X(J).NE.ZERO) THEN
272                          TEMP = X(J)
273                          L = KPLUS1 - J
274                          DO 10 I = MAX(1,J-K),J - 1
275                              X(I) = X(I) + TEMP*A(L+I,J)
276   10                     CONTINUE
277                          IF (NOUNIT) X(J) = X(J)*A(KPLUS1,J)
278                      END IF
279   20             CONTINUE
280              ELSE
281                  JX = KX
282                  DO 40 J = 1,N
283                      IF (X(JX).NE.ZERO) THEN
284                          TEMP = X(JX)
285                          IX = KX
286                          L = KPLUS1 - J
287                          DO 30 I = MAX(1,J-K),J - 1
288                              X(IX) = X(IX) + TEMP*A(L+I,J)
289                              IX = IX + INCX
290   30                     CONTINUE
291                          IF (NOUNIT) X(JX) = X(JX)*A(KPLUS1,J)
292                      END IF
293                      JX = JX + INCX
294                      IF (J.GT.K) KX = KX + INCX
295   40             CONTINUE
296              END IF
297          ELSE
298              IF (INCX.EQ.1) THEN
299                  DO 60 J = N,1,-1
300                      IF (X(J).NE.ZERO) THEN
301                          TEMP = X(J)
302                          L = 1 - J
303                          DO 50 I = MIN(N,J+K),J + 1,-1
304                              X(I) = X(I) + TEMP*A(L+I,J)
305   50                     CONTINUE
306                          IF (NOUNIT) X(J) = X(J)*A(1,J)
307                      END IF
308   60             CONTINUE
309              ELSE
310                  KX = KX + (N-1)*INCX
311                  JX = KX
312                  DO 80 J = N,1,-1
313                      IF (X(JX).NE.ZERO) THEN
314                          TEMP = X(JX)
315                          IX = KX
316                          L = 1 - J
317                          DO 70 I = MIN(N,J+K),J + 1,-1
318                              X(IX) = X(IX) + TEMP*A(L+I,J)
319                              IX = IX - INCX
320   70                     CONTINUE
321                          IF (NOUNIT) X(JX) = X(JX)*A(1,J)
322                      END IF
323                      JX = JX - INCX
324                      IF ((N-J).GE.K) KX = KX - INCX
325   80             CONTINUE
326              END IF
327          END IF
328      ELSE
329*
330*        Form  x := A**T*x.
331*
332          IF (LSAME(UPLO,'U')) THEN
333              KPLUS1 = K + 1
334              IF (INCX.EQ.1) THEN
335                  DO 100 J = N,1,-1
336                      TEMP = X(J)
337                      L = KPLUS1 - J
338                      IF (NOUNIT) TEMP = TEMP*A(KPLUS1,J)
339                      DO 90 I = J - 1,MAX(1,J-K),-1
340                          TEMP = TEMP + A(L+I,J)*X(I)
341   90                 CONTINUE
342                      X(J) = TEMP
343  100             CONTINUE
344              ELSE
345                  KX = KX + (N-1)*INCX
346                  JX = KX
347                  DO 120 J = N,1,-1
348                      TEMP = X(JX)
349                      KX = KX - INCX
350                      IX = KX
351                      L = KPLUS1 - J
352                      IF (NOUNIT) TEMP = TEMP*A(KPLUS1,J)
353                      DO 110 I = J - 1,MAX(1,J-K),-1
354                          TEMP = TEMP + A(L+I,J)*X(IX)
355                          IX = IX - INCX
356  110                 CONTINUE
357                      X(JX) = TEMP
358                      JX = JX - INCX
359  120             CONTINUE
360              END IF
361          ELSE
362              IF (INCX.EQ.1) THEN
363                  DO 140 J = 1,N
364                      TEMP = X(J)
365                      L = 1 - J
366                      IF (NOUNIT) TEMP = TEMP*A(1,J)
367                      DO 130 I = J + 1,MIN(N,J+K)
368                          TEMP = TEMP + A(L+I,J)*X(I)
369  130                 CONTINUE
370                      X(J) = TEMP
371  140             CONTINUE
372              ELSE
373                  JX = KX
374                  DO 160 J = 1,N
375                      TEMP = X(JX)
376                      KX = KX + INCX
377                      IX = KX
378                      L = 1 - J
379                      IF (NOUNIT) TEMP = TEMP*A(1,J)
380                      DO 150 I = J + 1,MIN(N,J+K)
381                          TEMP = TEMP + A(L+I,J)*X(IX)
382                          IX = IX + INCX
383  150                 CONTINUE
384                      X(JX) = TEMP
385                      JX = JX + INCX
386  160             CONTINUE
387              END IF
388          END IF
389      END IF
390*
391      RETURN
392*
393*     End of STBMV
394*
395      END
396