zgbmv Subroutine

public subroutine zgbmv(trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy)

NAME

zgbmv(3f) - [BLAS:COMPLEX_16_BLAS_LEVEL2]

SYNOPSIS

 subroutine zgbmv(trans,m,n,kl,ku,alpha,a,lda,x,incx,beta,y,incy)

   .. Scalar Arguments ..
   complex(kind=real64),intent(in)    :: alpha,beta
   integer,intent(in)                 :: incx,incy,kl,ku,lda,m,n
   character,intent(in)               :: trans
   ..
   .. Array Arguments ..
   complex(kind=real64),intent(in)    :: a(lda,*),x(*)
   complex(kind=real64),intent(inout) :: y(*)
   ..

DEFINITION

ZGBMV performs one of the matrix-vector operations

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

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

where alpha and beta are scalars, x and y are vectors and A is an m by n band matrix, with kl sub-diagonals and ku super-diagonals.

OPTIONS

TRANS

       TRANS is CHARACTER*1
        On entry, TRANS specifies the operation to be performed as
        follows:

           TRANS = 'N' or 'n'   y := alpha*A*x + beta*y.

           TRANS = 'T' or 't'   y := alpha*A**T*x + beta*y.

           TRANS = 'C' or 'c'   y := alpha*A**H*x + beta*y.

M

       M is INTEGER
        On entry, M specifies the number of rows of the matrix A.
        M must be at least zero.

N

       N is INTEGER
        On entry, N specifies the number of columns of the matrix A.
        N must be at least zero.

KL

       KL is INTEGER
        On entry, KL specifies the number of sub-diagonals of the
        matrix A. KL must satisfy 0 .le. KL.

KU

       KU is INTEGER
        On entry, KU specifies the number of super-diagonals of the
        matrix A. KU must satisfy 0 .le. KU.

ALPHA

       ALPHA is complex(kind=real64)
        On entry, ALPHA specifies the scalar alpha.

A

       A is complex(kind=real64) array, dimension ( LDA, N )
        Before entry, the leading ( kl + ku + 1 ) by n part of the
        array A must contain the matrix of coefficients, supplied
        column by column, with the leading diagonal of the matrix in
        row ( ku + 1 ) of the array, the first super-diagonal
        starting at position 2 in row ku, the first sub-diagonal
        starting at position 1 in row ( ku + 2 ), and so on.
        Elements in the array A that do not correspond to elements
        in the band matrix (such as the top left ku by ku triangle)
        are not referenced.
        The following program segment will transfer a band matrix
        from conventional full matrix storage to band storage:

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

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
        ( kl + ku + 1 ).

X

       X is complex(kind=real64) array, dimension at least
        ( 1 + ( n - 1 )*abs( INCX ) ) when TRANS = 'N' or 'n'
        and at least
        ( 1 + ( m - 1 )*abs( INCX ) ) otherwise.
        Before entry, the incremented array X must contain the
        vector x.

INCX

       INCX is INTEGER
        On entry, INCX specifies the increment for the elements of
        X. INCX must not be zero.

BETA

       BETA is complex(kind=real64)
        On entry, BETA specifies the scalar beta. When BETA is
        supplied as zero then Y need not be set on input.

Y

       Y is complex(kind=real64) array, dimension at least
        ( 1 + ( m - 1 )*abs( INCY ) ) when TRANS = 'N' or 'n'
        and at least
        ( 1 + ( n - 1 )*abs( INCY ) ) otherwise.
        Before entry, the incremented array Y must contain the
        vector y. On exit, Y is overwritten by the updated vector y.

INCY

       INCY is INTEGER
        On entry, INCY specifies the increment for the elements of
        Y. INCY must not be zero.

AUTHORS

  • Univ. of Tennessee
  • Univ. of California Berkeley
  • Univ. of Colorado Denver
  • NAG Ltd.

date:December 2016

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.

SEE ALSO

Online html documentation available at
http://www.netlib.org/lapack/explore-html/

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: trans
integer, intent(in) :: m
integer, intent(in) :: n
integer, intent(in) :: kl
integer, intent(in) :: ku
complex(kind=real64), intent(in) :: alpha
complex(kind=real64), intent(in) :: a(lda,*)
integer, intent(in) :: lda
complex(kind=real64), intent(in) :: x(*)
integer, intent(in) :: incx
complex(kind=real64), intent(in) :: beta
complex(kind=real64), intent(inout) :: y(*)
integer, intent(in) :: incy

Contents

Source Code


Variables

Type Visibility Attributes Name Initial
integer, public :: i
integer, public :: info
integer, public :: ix
integer, public :: iy
integer, public :: j
integer, public :: jx
integer, public :: jy
integer, public :: k
integer, public :: kup1
integer, public :: kx
integer, public :: ky
integer, public :: lenx
integer, public :: leny
logical, public :: noconj
complex(kind=real64), public, parameter :: one = (1.0d+0,0.0d+0)
complex(kind=real64), public :: temp
complex(kind=real64), public, parameter :: zero = (0.0d+0,0.0d+0)

Source Code

       subroutine zgbmv(trans,m,n,kl,ku,alpha,a,lda,x,incx,beta,y,incy)
      implicit none
!
!  -- Reference BLAS level2 routine (version 3.7.0) --
!  -- Reference BLAS is a software package provided by Univ. of Tennessee,    --
!  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
!     December 2016
!
!     .. Scalar Arguments ..
      complex(kind=real64),intent(in)    :: alpha,beta
      integer,intent(in)                 :: incx,incy,kl,ku,lda,m,n
      character,intent(in)               :: trans
!     ..
!     .. Array Arguments ..
      complex(kind=real64),intent(in)    :: a(lda,*),x(*)
      complex(kind=real64),intent(inout) :: y(*)
!     ..
!
!  =====================================================================
!
!     .. Parameters ..
      complex(kind=real64) :: one
      parameter (one= (1.0d+0,0.0d+0))
      complex(kind=real64) :: zero
      parameter (zero= (0.0d+0,0.0d+0))
!     ..
!     .. Local Scalars ..
      complex(kind=real64) :: temp
      integer i,info,ix,iy,j,jx,jy,k,kup1,kx,ky,lenx,leny
      logical noconj
!     ..
!     .. External Functions ..
!      LOGICAL LSAME
!      EXTERNAL LSAME
!     ..
!     .. External Subroutines ..
!      EXTERNAL XERBLA
!     ..
!     .. Intrinsic Functions ..
      intrinsic dconjg,max,min
!     ..
!
!     Test the input parameters.
!
      info = 0
      if (.not.lsame(trans,'N') .and. .not.lsame(trans,'T') .and.  .not.lsame(trans,'C')) then
          info = 1
      elseif (m.lt.0) then
          info = 2
      elseif (n.lt.0) then
          info = 3
      elseif (kl.lt.0) then
          info = 4
      elseif (ku.lt.0) then
          info = 5
      elseif (lda.lt. (kl+ku+1)) then
          info = 8
      elseif (incx.eq.0) then
          info = 10
      elseif (incy.eq.0) then
          info = 13
      endif
      if (info.ne.0) then
          call xerbla('ZGBMV ',info)
          return
      endif
!
!     Quick return if possible.
!
      if ((m.eq.0) .or. (n.eq.0) .or.  ((alpha.eq.zero).and. (beta.eq.one))) return
!
      noconj = lsame(trans,'T')
!
!     Set  LENX  and  LENY, the lengths of the vectors x and y, and set
!     up the start points in  X  and  Y.
!
      if (lsame(trans,'N')) then
          lenx = n
          leny = m
      else
          lenx = m
          leny = n
      endif
      if (incx.gt.0) then
          kx = 1
      else
          kx = 1 - (lenx-1)*incx
      endif
      if (incy.gt.0) then
          ky = 1
      else
          ky = 1 - (leny-1)*incy
      endif
!
!     Start the operations. In this version the elements of A are
!     accessed sequentially with one pass through the band part of A.
!
!     First form  y := beta*y.
!
      if (beta.ne.one) then
          if (incy.eq.1) then
              if (beta.eq.zero) then
                  y(1:leny) = zero
              else
                  y(1:leny) = beta*y(1:leny)
              endif
          else
              iy = ky
              if (beta.eq.zero) then
                  do i = 1,leny
                      y(iy) = zero
                      iy = iy + incy
                  enddo
              else
                  do i = 1,leny
                      y(iy) = beta*y(iy)
                      iy = iy + incy
                  enddo
              endif
          endif
      endif
      if (alpha.eq.zero) return
      kup1 = ku + 1
      if (lsame(trans,'N')) then
!
!        Form  y := alpha*A*x + y.
!
          jx = kx
          if (incy.eq.1) then
              do j = 1,n
                  temp = alpha*x(jx)
                  k = kup1 - j
                  do i = max(1,j-ku),min(m,j+kl)
                      y(i) = y(i) + temp*a(k+i,j)
                  enddo
                  jx = jx + incx
              enddo
          else
              do j = 1,n
                  temp = alpha*x(jx)
                  iy = ky
                  k = kup1 - j
                  do i = max(1,j-ku),min(m,j+kl)
                      y(iy) = y(iy) + temp*a(k+i,j)
                      iy = iy + incy
                  enddo
                  jx = jx + incx
                  if (j.gt.ku) ky = ky + incy
              enddo
          endif
      else
!
!        Form  y := alpha*A**T*x + y  or  y := alpha*A**H*x + y.
!
          jy = ky
          if (incx.eq.1) then
              do j = 1,n
                  temp = zero
                  k = kup1 - j
                  if (noconj) then
                      do i = max(1,j-ku),min(m,j+kl)
                          temp = temp + a(k+i,j)*x(i)
                      enddo
                  else
                      do i = max(1,j-ku),min(m,j+kl)
                          temp = temp + dconjg(a(k+i,j))*x(i)
                      enddo
                  endif
                  y(jy) = y(jy) + alpha*temp
                  jy = jy + incy
              enddo
          else
              do j = 1,n
                  temp = zero
                  ix = kx
                  k = kup1 - j
                  if (noconj) then
                      do i = max(1,j-ku),min(m,j+kl)
                          temp = temp + a(k+i,j)*x(ix)
                          ix = ix + incx
                      enddo
                  else
                      do i = max(1,j-ku),min(m,j+kl)
                          temp = temp + dconjg(a(k+i,j))*x(ix)
                          ix = ix + incx
                      enddo
                  endif
                  y(jy) = y(jy) + alpha*temp
                  jy = jy + incy
                  if (j.gt.ku) kx = kx + incx
              enddo
          endif
      endif
!
      end subroutine zgbmv