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

◆ dsbmv()

subroutine dsbmv ( character  uplo,
integer  n,
integer  k,
double precision  alpha,
double precision, dimension(lda,*)  a,
integer  lda,
double precision, dimension(*)  x,
integer  incx,
double precision  beta,
double precision, dimension(*)  y,
integer  incy 
)

DSBMV

Purpose:
 DSBMV  performs the matrix-vector  operation

    y := alpha*A*x + beta*y,

 where alpha and beta are scalars, x and y are n element vectors and
 A is an n by n symmetric band matrix, with k super-diagonals.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
           On entry, UPLO specifies whether the upper or lower
           triangular part of the band matrix A is being supplied as
           follows:

              UPLO = 'U' or 'u'   The upper triangular part of A is
                                  being supplied.

              UPLO = 'L' or 'l'   The lower triangular part of A is
                                  being supplied.
[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, K specifies the number of super-diagonals of the
           matrix A. K must satisfy  0 .le. K.
[in]ALPHA
          ALPHA is DOUBLE PRECISION.
           On entry, ALPHA specifies the scalar alpha.
[in]A
          A is DOUBLE PRECISION 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 symmetric matrix, 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 the upper
           triangular part of a symmetric 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 symmetric matrix, 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 the lower
           triangular part of a symmetric 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
[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]X
          X is DOUBLE PRECISION array, dimension at least
           ( 1 + ( n - 1 )*abs( INCX ) ).
           Before entry, the incremented array X must contain the
           vector x.
[in]INCX
          INCX is INTEGER
           On entry, INCX specifies the increment for the elements of
           X. INCX must not be zero.
[in]BETA
          BETA is DOUBLE PRECISION.
           On entry, BETA specifies the scalar beta.
[in,out]Y
          Y is DOUBLE PRECISION array, dimension at least
           ( 1 + ( n - 1 )*abs( INCY ) ).
           Before entry, the incremented array Y must contain the
           vector y. On exit, Y is overwritten by the updated vector y.
[in]INCY
          INCY is INTEGER
           On entry, INCY specifies the increment for the elements of
           Y. INCY 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 183 of file dsbmv.f.

184*
185* -- Reference BLAS level2 routine --
186* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
187* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
188*
189* .. Scalar Arguments ..
190 DOUBLE PRECISION ALPHA,BETA
191 INTEGER INCX,INCY,K,LDA,N
192 CHARACTER UPLO
193* ..
194* .. Array Arguments ..
195 DOUBLE PRECISION A(LDA,*),X(*),Y(*)
196* ..
197*
198* =====================================================================
199*
200* .. Parameters ..
201 DOUBLE PRECISION ONE,ZERO
202 parameter(one=1.0d+0,zero=0.0d+0)
203* ..
204* .. Local Scalars ..
205 DOUBLE PRECISION TEMP1,TEMP2
206 INTEGER I,INFO,IX,IY,J,JX,JY,KPLUS1,KX,KY,L
207* ..
208* .. External Functions ..
209 LOGICAL LSAME
210 EXTERNAL lsame
211* ..
212* .. External Subroutines ..
213 EXTERNAL xerbla
214* ..
215* .. Intrinsic Functions ..
216 INTRINSIC max,min
217* ..
218*
219* Test the input parameters.
220*
221 info = 0
222 IF (.NOT.lsame(uplo,'U') .AND. .NOT.lsame(uplo,'L')) THEN
223 info = 1
224 ELSE IF (n.LT.0) THEN
225 info = 2
226 ELSE IF (k.LT.0) THEN
227 info = 3
228 ELSE IF (lda.LT. (k+1)) THEN
229 info = 6
230 ELSE IF (incx.EQ.0) THEN
231 info = 8
232 ELSE IF (incy.EQ.0) THEN
233 info = 11
234 END IF
235 IF (info.NE.0) THEN
236 CALL xerbla('DSBMV ',info)
237 RETURN
238 END IF
239*
240* Quick return if possible.
241*
242 IF ((n.EQ.0) .OR. ((alpha.EQ.zero).AND. (beta.EQ.one))) RETURN
243*
244* Set up the start points in X and Y.
245*
246 IF (incx.GT.0) THEN
247 kx = 1
248 ELSE
249 kx = 1 - (n-1)*incx
250 END IF
251 IF (incy.GT.0) THEN
252 ky = 1
253 ELSE
254 ky = 1 - (n-1)*incy
255 END IF
256*
257* Start the operations. In this version the elements of the array A
258* are accessed sequentially with one pass through A.
259*
260* First form y := beta*y.
261*
262 IF (beta.NE.one) THEN
263 IF (incy.EQ.1) THEN
264 IF (beta.EQ.zero) THEN
265 DO 10 i = 1,n
266 y(i) = zero
267 10 CONTINUE
268 ELSE
269 DO 20 i = 1,n
270 y(i) = beta*y(i)
271 20 CONTINUE
272 END IF
273 ELSE
274 iy = ky
275 IF (beta.EQ.zero) THEN
276 DO 30 i = 1,n
277 y(iy) = zero
278 iy = iy + incy
279 30 CONTINUE
280 ELSE
281 DO 40 i = 1,n
282 y(iy) = beta*y(iy)
283 iy = iy + incy
284 40 CONTINUE
285 END IF
286 END IF
287 END IF
288 IF (alpha.EQ.zero) RETURN
289 IF (lsame(uplo,'U')) THEN
290*
291* Form y when upper triangle of A is stored.
292*
293 kplus1 = k + 1
294 IF ((incx.EQ.1) .AND. (incy.EQ.1)) THEN
295 DO 60 j = 1,n
296 temp1 = alpha*x(j)
297 temp2 = zero
298 l = kplus1 - j
299 DO 50 i = max(1,j-k),j - 1
300 y(i) = y(i) + temp1*a(l+i,j)
301 temp2 = temp2 + a(l+i,j)*x(i)
302 50 CONTINUE
303 y(j) = y(j) + temp1*a(kplus1,j) + alpha*temp2
304 60 CONTINUE
305 ELSE
306 jx = kx
307 jy = ky
308 DO 80 j = 1,n
309 temp1 = alpha*x(jx)
310 temp2 = zero
311 ix = kx
312 iy = ky
313 l = kplus1 - j
314 DO 70 i = max(1,j-k),j - 1
315 y(iy) = y(iy) + temp1*a(l+i,j)
316 temp2 = temp2 + a(l+i,j)*x(ix)
317 ix = ix + incx
318 iy = iy + incy
319 70 CONTINUE
320 y(jy) = y(jy) + temp1*a(kplus1,j) + alpha*temp2
321 jx = jx + incx
322 jy = jy + incy
323 IF (j.GT.k) THEN
324 kx = kx + incx
325 ky = ky + incy
326 END IF
327 80 CONTINUE
328 END IF
329 ELSE
330*
331* Form y when lower triangle of A is stored.
332*
333 IF ((incx.EQ.1) .AND. (incy.EQ.1)) THEN
334 DO 100 j = 1,n
335 temp1 = alpha*x(j)
336 temp2 = zero
337 y(j) = y(j) + temp1*a(1,j)
338 l = 1 - j
339 DO 90 i = j + 1,min(n,j+k)
340 y(i) = y(i) + temp1*a(l+i,j)
341 temp2 = temp2 + a(l+i,j)*x(i)
342 90 CONTINUE
343 y(j) = y(j) + alpha*temp2
344 100 CONTINUE
345 ELSE
346 jx = kx
347 jy = ky
348 DO 120 j = 1,n
349 temp1 = alpha*x(jx)
350 temp2 = zero
351 y(jy) = y(jy) + temp1*a(1,j)
352 l = 1 - j
353 ix = jx
354 iy = jy
355 DO 110 i = j + 1,min(n,j+k)
356 ix = ix + incx
357 iy = iy + incy
358 y(iy) = y(iy) + temp1*a(l+i,j)
359 temp2 = temp2 + a(l+i,j)*x(ix)
360 110 CONTINUE
361 y(jy) = y(jy) + alpha*temp2
362 jx = jx + incx
363 jy = jy + incy
364 120 CONTINUE
365 END IF
366 END IF
367*
368 RETURN
369*
370* End of DSBMV
371*
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: