LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ stbmv()

subroutine stbmv ( character  uplo,
character  trans,
character  diag,
integer  n,
integer  k,
real, dimension(lda,*)  a,
integer  lda,
real, dimension(*)  x,
integer  incx 
)

STBMV

Purpose:
 STBMV  performs one of the matrix-vector operations

    x := A*x,   or   x := A**T*x,

 where x is an n element vector and  A is an n by n unit, or non-unit,
 upper or lower triangular band matrix, with ( k + 1 ) diagonals.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
           On entry, UPLO specifies whether the matrix is an upper or
           lower triangular matrix as follows:

              UPLO = 'U' or 'u'   A is an upper triangular matrix.

              UPLO = 'L' or 'l'   A is a lower triangular matrix.
[in]TRANS
          TRANS is CHARACTER*1
           On entry, TRANS specifies the operation to be performed as
           follows:

              TRANS = 'N' or 'n'   x := A*x.

              TRANS = 'T' or 't'   x := A**T*x.

              TRANS = 'C' or 'c'   x := A**T*x.
[in]DIAG
          DIAG is CHARACTER*1
           On entry, DIAG specifies whether or not A is unit
           triangular as follows:

              DIAG = 'U' or 'u'   A is assumed to be unit triangular.

              DIAG = 'N' or 'n'   A is not assumed to be unit
                                  triangular.
[in]N
          N is INTEGER
           On entry, N specifies the order of the matrix A.
           N must be at least zero.
[in]K
          K is INTEGER
           On entry with UPLO = 'U' or 'u', K specifies the number of
           super-diagonals of the matrix A.
           On entry with UPLO = 'L' or 'l', K specifies the number of
           sub-diagonals of the matrix A.
           K must satisfy  0 .le. K.
[in]A
          A is REAL array, dimension ( LDA, N )
           Before entry with UPLO = 'U' or 'u', the leading ( k + 1 )
           by n part of the array A must contain the upper triangular
           band part of the matrix of coefficients, supplied column by
           column, with the leading diagonal of the matrix in row
           ( k + 1 ) of the array, the first super-diagonal starting at
           position 2 in row k, and so on. The top left k by k triangle
           of the array A is not referenced.
           The following program segment will transfer an upper
           triangular band matrix from conventional full matrix storage
           to band storage:

                 DO 20, J = 1, N
                    M = K + 1 - J
                    DO 10, I = MAX( 1, J - K ), J
                       A( M + I, J ) = matrix( I, J )
              10    CONTINUE
              20 CONTINUE

           Before entry with UPLO = 'L' or 'l', the leading ( k + 1 )
           by n part of the array A must contain the lower triangular
           band part of the matrix of coefficients, supplied column by
           column, with the leading diagonal of the matrix in row 1 of
           the array, the first sub-diagonal starting at position 1 in
           row 2, and so on. The bottom right k by k triangle of the
           array A is not referenced.
           The following program segment will transfer a lower
           triangular band matrix from conventional full matrix storage
           to band storage:

                 DO 20, J = 1, N
                    M = 1 - J
                    DO 10, I = J, MIN( N, J + K )
                       A( M + I, J ) = matrix( I, J )
              10    CONTINUE
              20 CONTINUE

           Note that when DIAG = 'U' or 'u' the elements of the array A
           corresponding to the diagonal elements of the matrix are not
           referenced, but are assumed to be unity.
[in]LDA
          LDA is INTEGER
           On entry, LDA specifies the first dimension of A as declared
           in the calling (sub) program. LDA must be at least
           ( k + 1 ).
[in,out]X
          X is REAL array, dimension at least
           ( 1 + ( n - 1 )*abs( INCX ) ).
           Before entry, the incremented array X must contain the n
           element vector x. On exit, X is overwritten with the
           transformed vector x.
[in]INCX
          INCX is INTEGER
           On entry, INCX specifies the increment for the elements of
           X. INCX must not be zero.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
  Level 2 Blas routine.
  The vector and matrix arguments are not referenced when N = 0, or M = 0

  -- Written on 22-October-1986.
     Jack Dongarra, Argonne National Lab.
     Jeremy Du Croz, Nag Central Office.
     Sven Hammarling, Nag Central Office.
     Richard Hanson, Sandia National Labs.

Definition at line 185 of file stbmv.f.

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.
227 + .NOT.lsame(trans,'T') .AND.
228 + .NOT.lsame(trans,'C')) THEN
229 info = 2
230 ELSE IF (.NOT.lsame(diag,'U') .AND.
231 + .NOT.lsame(diag,'N')) THEN
232 info = 3
233 ELSE IF (n.LT.0) THEN
234 info = 4
235 ELSE IF (k.LT.0) THEN
236 info = 5
237 ELSE IF (lda.LT. (k+1)) THEN
238 info = 7
239 ELSE IF (incx.EQ.0) THEN
240 info = 9
241 END IF
242 IF (info.NE.0) THEN
243 CALL xerbla('STBMV ',info)
244 RETURN
245 END IF
246*
247* Quick return if possible.
248*
249 IF (n.EQ.0) RETURN
250*
251 nounit = lsame(diag,'N')
252*
253* Set up the start point in X if the increment is not unity. This
254* will be ( N - 1 )*INCX too small for descending loops.
255*
256 IF (incx.LE.0) THEN
257 kx = 1 - (n-1)*incx
258 ELSE IF (incx.NE.1) THEN
259 kx = 1
260 END IF
261*
262* Start the operations. In this version the elements of A are
263* accessed sequentially with one pass through A.
264*
265 IF (lsame(trans,'N')) THEN
266*
267* Form x := A*x.
268*
269 IF (lsame(uplo,'U')) THEN
270 kplus1 = k + 1
271 IF (incx.EQ.1) THEN
272 DO 20 j = 1,n
273 IF (x(j).NE.zero) THEN
274 temp = x(j)
275 l = kplus1 - j
276 DO 10 i = max(1,j-k),j - 1
277 x(i) = x(i) + temp*a(l+i,j)
278 10 CONTINUE
279 IF (nounit) x(j) = x(j)*a(kplus1,j)
280 END IF
281 20 CONTINUE
282 ELSE
283 jx = kx
284 DO 40 j = 1,n
285 IF (x(jx).NE.zero) THEN
286 temp = x(jx)
287 ix = kx
288 l = kplus1 - j
289 DO 30 i = max(1,j-k),j - 1
290 x(ix) = x(ix) + temp*a(l+i,j)
291 ix = ix + incx
292 30 CONTINUE
293 IF (nounit) x(jx) = x(jx)*a(kplus1,j)
294 END IF
295 jx = jx + incx
296 IF (j.GT.k) kx = kx + incx
297 40 CONTINUE
298 END IF
299 ELSE
300 IF (incx.EQ.1) THEN
301 DO 60 j = n,1,-1
302 IF (x(j).NE.zero) THEN
303 temp = x(j)
304 l = 1 - j
305 DO 50 i = min(n,j+k),j + 1,-1
306 x(i) = x(i) + temp*a(l+i,j)
307 50 CONTINUE
308 IF (nounit) x(j) = x(j)*a(1,j)
309 END IF
310 60 CONTINUE
311 ELSE
312 kx = kx + (n-1)*incx
313 jx = kx
314 DO 80 j = n,1,-1
315 IF (x(jx).NE.zero) THEN
316 temp = x(jx)
317 ix = kx
318 l = 1 - j
319 DO 70 i = min(n,j+k),j + 1,-1
320 x(ix) = x(ix) + temp*a(l+i,j)
321 ix = ix - incx
322 70 CONTINUE
323 IF (nounit) x(jx) = x(jx)*a(1,j)
324 END IF
325 jx = jx - incx
326 IF ((n-j).GE.k) kx = kx - incx
327 80 CONTINUE
328 END IF
329 END IF
330 ELSE
331*
332* Form x := A**T*x.
333*
334 IF (lsame(uplo,'U')) THEN
335 kplus1 = k + 1
336 IF (incx.EQ.1) THEN
337 DO 100 j = n,1,-1
338 temp = x(j)
339 l = kplus1 - j
340 IF (nounit) temp = temp*a(kplus1,j)
341 DO 90 i = j - 1,max(1,j-k),-1
342 temp = temp + a(l+i,j)*x(i)
343 90 CONTINUE
344 x(j) = temp
345 100 CONTINUE
346 ELSE
347 kx = kx + (n-1)*incx
348 jx = kx
349 DO 120 j = n,1,-1
350 temp = x(jx)
351 kx = kx - incx
352 ix = kx
353 l = kplus1 - j
354 IF (nounit) temp = temp*a(kplus1,j)
355 DO 110 i = j - 1,max(1,j-k),-1
356 temp = temp + a(l+i,j)*x(ix)
357 ix = ix - incx
358 110 CONTINUE
359 x(jx) = temp
360 jx = jx - incx
361 120 CONTINUE
362 END IF
363 ELSE
364 IF (incx.EQ.1) THEN
365 DO 140 j = 1,n
366 temp = x(j)
367 l = 1 - j
368 IF (nounit) temp = temp*a(1,j)
369 DO 130 i = j + 1,min(n,j+k)
370 temp = temp + a(l+i,j)*x(i)
371 130 CONTINUE
372 x(j) = temp
373 140 CONTINUE
374 ELSE
375 jx = kx
376 DO 160 j = 1,n
377 temp = x(jx)
378 kx = kx + incx
379 ix = kx
380 l = 1 - j
381 IF (nounit) temp = temp*a(1,j)
382 DO 150 i = j + 1,min(n,j+k)
383 temp = temp + a(l+i,j)*x(ix)
384 ix = ix + incx
385 150 CONTINUE
386 x(jx) = temp
387 jx = jx + incx
388 160 CONTINUE
389 END IF
390 END IF
391 END IF
392*
393 RETURN
394*
395* End of STBMV
396*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: