zhpmv Subroutine

public subroutine zhpmv(uplo, n, alpha, ap, x, incx, beta, y, incy)

NAME

zhpmv(3f) - [BLAS:COMPLEX_16_BLAS_LEVEL2]

SYNOPSIS

 subroutine zhpmv(uplo,n,alpha,ap,x,incx,beta,y,incy)

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

DEFINITION

ZHPMV 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, supplied in packed form.

OPTIONS

UPLO

       UPLO is CHARACTER*1
        On entry, UPLO specifies whether the upper or lower
        triangular part of the matrix A is supplied in the packed
        array AP as follows:

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

           UPLO = 'L' or 'l'   The lower triangular part of A is
                               supplied in AP.

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.

AP

       AP is complex(kind=real64) array, dimension at least
        ( ( n*( n + 1 ) )/2 ).
        Before entry with UPLO = 'U' or 'u', the array AP must
        contain the upper triangular part of the hermitian matrix
        packed sequentially, column by column, so that AP( 1 )
        contains a( 1, 1 ), AP( 2 ) and AP( 3 ) contain a( 1, 2 )
        and a( 2, 2 ) respectively, and so on.
        Before entry with UPLO = 'L' or 'l', the array AP must
        contain the lower triangular part of the hermitian matrix
        packed sequentially, column by column, so that AP( 1 )
        contains a( 1, 1 ), AP( 2 ) and AP( 3 ) contain a( 2, 1 )
        and a( 3, 1 ) respectively, and so on.
        Note that the imaginary parts of the diagonal elements need
        not be set and are assumed to be zero.

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) :: ap(*)
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 :: kk
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 zhpmv(uplo,n,alpha,ap,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,n
      character,intent(in)               :: uplo
!     ..
!     .. Array Arguments ..
      complex(kind=real64),intent(in)    :: ap(*),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,k,kk,kx,ky
!     ..
!     .. External Functions ..
!      LOGICAL LSAME
!      EXTERNAL LSAME
!     ..
!     .. External Subroutines ..
!      EXTERNAL XERBLA
!     ..
!     .. Intrinsic Functions ..
      intrinsic dble,dconjg
!     ..
!
!     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 (incx.eq.0) then
          info = 6
      elseif (incy.eq.0) then
          info = 9
      endif
      if (info.ne.0) then
          call xerbla('ZHPMV ',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 the array AP
!     are accessed sequentially with one pass through AP.
!
!     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
      kk = 1
      if (lsame(uplo,'U')) then
!
!        Form  y  when AP contains the upper triangle.
!
          if ((incx.eq.1) .and. (incy.eq.1)) then
              do j = 1,n
                  temp1 = alpha*x(j)
                  temp2 = zero
                  k = kk
                  do i = 1,j - 1
                      y(i) = y(i) + temp1*ap(k)
                      temp2 = temp2 + dconjg(ap(k))*x(i)
                      k = k + 1
                  enddo
                  y(j) = y(j) + temp1*dble(ap(kk+j-1)) + alpha*temp2
                  kk = kk + j
              enddo
          else
              jx = kx
              jy = ky
              do j = 1,n
                  temp1 = alpha*x(jx)
                  temp2 = zero
                  ix = kx
                  iy = ky
                  do k = kk,kk + j - 2
                      y(iy) = y(iy) + temp1*ap(k)
                      temp2 = temp2 + dconjg(ap(k))*x(ix)
                      ix = ix + incx
                      iy = iy + incy
                  enddo
                  y(jy) = y(jy) + temp1*dble(ap(kk+j-1)) + alpha*temp2
                  jx = jx + incx
                  jy = jy + incy
                  kk = kk + j
              enddo
          endif
      else
!
!        Form  y  when AP contains the 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(ap(kk))
                  k = kk + 1
                  do i = j + 1,n
                      y(i) = y(i) + temp1*ap(k)
                      temp2 = temp2 + dconjg(ap(k))*x(i)
                      k = k + 1
                  enddo
                  y(j) = y(j) + alpha*temp2
                  kk = kk + (n-j+1)
              enddo
          else
              jx = kx
              jy = ky
              do j = 1,n
                  temp1 = alpha*x(jx)
                  temp2 = zero
                  y(jy) = y(jy) + temp1*dble(ap(kk))
                  ix = jx
                  iy = jy
                  do k = kk + 1,kk + n - j
                      ix = ix + incx
                      iy = iy + incy
                      y(iy) = y(iy) + temp1*ap(k)
                      temp2 = temp2 + dconjg(ap(k))*x(ix)
                  enddo
                  y(jy) = y(jy) + alpha*temp2
                  jx = jx + incx
                  jy = jy + incy
                  kk = kk + (n-j+1)
              enddo
          endif
      endif

      end subroutine zhpmv