LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine chbmv ( character  UPLO,
integer  N,
integer  K,
complex  ALPHA,
complex, dimension(lda,*)  A,
integer  LDA,
complex, dimension(*)  X,
integer  INCX,
complex  BETA,
complex, dimension(*)  Y,
integer  INCY 
)

CHBMV

Purpose:
 CHBMV  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 hermitian 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 COMPLEX
           On entry, ALPHA specifies the scalar alpha.
[in]A
          A is COMPLEX array of 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 hermitian 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 hermitian 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 hermitian 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 hermitian 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 the imaginary parts of the diagonal elements need
           not be set and are assumed to be zero.
[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 COMPLEX array of 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 COMPLEX
           On entry, BETA specifies the scalar beta.
[in,out]Y
          Y is COMPLEX array of 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.
Date
November 2011
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 189 of file chbmv.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: