zhemv Subroutine

public subroutine zhemv(uplo, n, alpha, a, lda, x, incx, beta, y, incy)

NAME

zhemv(3f) - [BLAS:COMPLEX_16_BLAS_LEVEL2]

SYNOPSIS

 subroutine zhemv(uplo,n,alpha,a,lda,x,incx,beta,y,incy)

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

DEFINITION

ZHEMV 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 matrix.

OPTIONS

UPLO

       UPLO is CHARACTER*1
        On entry, UPLO specifies whether the upper or lower
        triangular part of the array A is to be referenced as
        follows:

           UPLO = 'U' or 'u'   Only the upper triangular part of A
                               is to be referenced.

           UPLO = 'L' or 'l'   Only the lower triangular part of A
                               is to be referenced.

N

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

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 with UPLO = 'U' or 'u', the leading n by n
        upper triangular part of the array A must contain the upper
        triangular part of the hermitian matrix and the strictly
        lower triangular part of A is not referenced.
        Before entry with UPLO = 'L' or 'l', the leading n by n
        lower triangular part of the array A must contain the lower
        triangular part of the hermitian matrix and the strictly
        upper triangular part of A is not referenced.
        Note that the imaginary parts of the diagonal elements need
        not be set and are assumed to be zero.

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
        max( 1, n ).

X

       X is complex(kind=real64) array, dimension at least
        ( 1 + ( n - 1 )*abs( INCX ) ).
        Before entry, the incremented array X must contain the n
        element 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 + ( n - 1 )*abs( INCY ) ).
        Before entry, the incremented array Y must contain the n
        element 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) :: uplo
integer, intent(in) :: n
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 :: kx
integer, public :: ky
complex(kind=real64), public, parameter :: one = (1.0d+0,0.0d+0)
complex(kind=real64), public :: temp1
complex(kind=real64), public :: temp2
complex(kind=real64), public, parameter :: zero = (0.0d+0,0.0d+0)

Source Code

subroutine zhemv(uplo,n,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,lda,n
      character,intent(in)               :: uplo
!     ..
!     .. 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) :: temp1,temp2
      integer i,info,ix,iy,j,jx,jy,kx,ky
!     ..
!     .. External Functions ..
!      logical lsame
!      external lsame
!     ..
!     .. External Subroutines ..
!      external xerbla
!     ..
!     .. Intrinsic Functions ..
      intrinsic dble,dconjg,max
!     ..
!
!     Test the input parameters.
!
      info = 0
      if (.not.lsame(uplo,'U') .and. .not.lsame(uplo,'L')) then
          info = 1
      elseif (n.lt.0) then
          info = 2
      elseif (lda.lt.max(1,n)) then
          info = 5
      elseif (incx.eq.0) then
          info = 7
      elseif (incy.eq.0) then
          info = 10
      endif
      if (info.ne.0) then
          call xerbla('ZHEMV ',info)
          return
      endif
!
!     Quick return if possible.
!
      if ((n.eq.0) .or. ((alpha.eq.zero).and. (beta.eq.one))) return
!
!     Set up the start points in  X  and  Y.
!
      if (incx.gt.0) then
          kx = 1
      else
          kx = 1 - (n-1)*incx
      endif
      if (incy.gt.0) then
          ky = 1
      else
          ky = 1 - (n-1)*incy
      endif
!
!     Start the operations. In this version the elements of A are
!     accessed sequentially with one pass through the triangular 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:n) = zero
              else
                  y(1:n) = beta*y(1:n)
              endif
          else
              iy = ky
              if (beta.eq.zero) then
                  do i = 1,n
                      y(iy) = zero
                      iy = iy + incy
                  enddo
              else
                  do i = 1,n
                      y(iy) = beta*y(iy)
                      iy = iy + incy
                  enddo
              endif
          endif
      endif
      if (alpha.eq.zero) return
      if (lsame(uplo,'U')) then
!
!        Form  y  when A is stored in upper triangle.
!
          if ((incx.eq.1) .and. (incy.eq.1)) then
              do j = 1,n
                  temp1 = alpha*x(j)
                  temp2 = zero
                  do i = 1,j - 1
                      y(i) = y(i) + temp1*a(i,j)
                      temp2 = temp2 + dconjg(a(i,j))*x(i)
                  enddo
                  y(j) = y(j) + temp1*dble(a(j,j)) + alpha*temp2
              enddo
          else
              jx = kx
              jy = ky
              do j = 1,n
                  temp1 = alpha*x(jx)
                  temp2 = zero
                  ix = kx
                  iy = ky
                  do i = 1,j - 1
                      y(iy) = y(iy) + temp1*a(i,j)
                      temp2 = temp2 + dconjg(a(i,j))*x(ix)
                      ix = ix + incx
                      iy = iy + incy
                  enddo
                  y(jy) = y(jy) + temp1*dble(a(j,j)) + alpha*temp2
                  jx = jx + incx
                  jy = jy + incy
              enddo
          endif
      else
!
!        Form  y  when A is stored in lower triangle.
!
          if ((incx.eq.1) .and. (incy.eq.1)) then
              do j = 1,n
                  temp1 = alpha*x(j)
                  temp2 = zero
                  y(j) = y(j) + temp1*dble(a(j,j))
                  do i = j + 1,n
                      y(i) = y(i) + temp1*a(i,j)
                      temp2 = temp2 + dconjg(a(i,j))*x(i)
                  enddo
                  y(j) = y(j) + alpha*temp2
              enddo
          else
              jx = kx
              jy = ky
              do j = 1,n
                  temp1 = alpha*x(jx)
                  temp2 = zero
                  y(jy) = y(jy) + temp1*dble(a(j,j))
                  ix = jx
                  iy = jy
                  do i = j + 1,n
                      ix = ix + incx
                      iy = iy + incy
                      y(iy) = y(iy) + temp1*a(i,j)
                      temp2 = temp2 + dconjg(a(i,j))*x(ix)
                  enddo
                  y(jy) = y(jy) + alpha*temp2
                  jx = jx + incx
                  jy = jy + incy
              enddo
          endif
      endif

end subroutine zhemv