cher2 Subroutine

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

NAME

cher2(3f) - [BLAS:COMPLEX_BLAS_LEVEL2] A := A + alpha*CX*CONJUGATE-TRANSPOSE(CY)n + CONJUGATE(alpha)*CY*CONJUGATE-TRANSPOSE(CX);
==> n A a (square) hermitian matrix.
(performs the hermitian rank 2 operation)

SYNOPSIS

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

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

DEFINITION

CHER2 performs the hermitian rank 2 operation

 A := alpha*x*y**H + conjg( alpha )*y*x**H + A,

where alpha is a scalar, 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
        On entry, ALPHA specifies the scalar alpha.

X

       X is COMPLEX 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.

Y

       Y is COMPLEX array, dimension at least
        ( 1 + ( n - 1 )*abs( INCY ) ).
        Before entry, the incremented array Y must contain the n
        element vector y.

INCY

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

A

       A is COMPLEX 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. On exit, the
        upper triangular part of the array A is overwritten by the
        upper triangular part of the updated matrix.
        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. On exit, the
        lower triangular part of the array A is overwritten by the
        lower triangular part of the updated matrix.
        Note that the imaginary parts of the diagonal elements need
        not be set, they are assumed to be zero, and on exit they
        are set to 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 ).

AUTHORS

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

date:December 2016

FURTHER DETAILS

Level 2 Blas routine.

– 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, intent(in) :: alpha
complex, intent(in) :: x(*)
integer, intent(in) :: incx
complex, intent(in) :: y(*)
integer, intent(in) :: incy
complex, intent(inout) :: a(lda,*)
integer, intent(in) :: lda

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, public :: temp1
complex, public :: temp2
complex, public, parameter :: zero = (0.0e+0,0.0e+0)

Source Code

       subroutine cher2(uplo,n,alpha,x,incx,y,incy,a,lda)
      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,intent(in)     :: alpha
      integer,intent(in)     :: incx,incy,lda,n
      character,intent(in)   :: uplo
!     ..
!     .. Array Arguments ..
      complex,intent(inout)  :: a(lda,*)
      complex,intent(in)     :: x(*),y(*)
!     ..
!
!  =====================================================================
!
!     .. Parameters ..
      complex zero
      parameter (zero= (0.0e+0,0.0e+0))
!     ..
!     .. Local Scalars ..
      complex temp1,temp2
      integer i,info,ix,iy,j,jx,jy,kx,ky
!     ..
!     .. External Functions ..
!     ..
!     .. External Subroutines ..
!     ..
!     .. Intrinsic Functions ..
      intrinsic conjg,max,real
!     ..
!
!     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 = 5
      elseif (incy.eq.0) then
          info = 7
      elseif (lda.lt.max(1,n)) then
          info = 9
      endif
      if (info.ne.0) then
          call xerbla('CHER2 ',info)
          return
      endif
!
!     Quick return if possible.
!
      if ((n.eq.0) .or. (alpha.eq.zero)) return
!
!     Set up the start points in X and Y if the increments are not both
!     unity.
!
      if ((incx.ne.1) .or. (incy.ne.1)) then
          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
          jx = kx
          jy = ky
      endif
!
!     Start the operations. In this version the elements of A are
!     accessed sequentially with one pass through the triangular part
!     of A.
!
      if (lsame(uplo,'U')) then
!
!        Form  A  when A is stored in the upper triangle.
!
          if ((incx.eq.1) .and. (incy.eq.1)) then
              do j = 1,n
                  if ((x(j).ne.zero) .or. (y(j).ne.zero)) then
                      temp1 = alpha*conjg(y(j))
                      temp2 = conjg(alpha*x(j))
                      do i = 1,j - 1
                          a(i,j) = a(i,j) + x(i)*temp1 + y(i)*temp2
                      enddo
                      a(j,j) = real(a(j,j)) + real(x(j)*temp1+y(j)*temp2)
                  else
                      a(j,j) = real(a(j,j))
                  endif
              enddo
          else
              do j = 1,n
                  if ((x(jx).ne.zero) .or. (y(jy).ne.zero)) then
                      temp1 = alpha*conjg(y(jy))
                      temp2 = conjg(alpha*x(jx))
                      ix = kx
                      iy = ky
                      do i = 1,j - 1
                          a(i,j) = a(i,j) + x(ix)*temp1 + y(iy)*temp2
                          ix = ix + incx
                          iy = iy + incy
                      enddo
                      a(j,j) = real(a(j,j)) + real(x(jx)*temp1+y(jy)*temp2)
                  else
                      a(j,j) = real(a(j,j))
                  endif
                  jx = jx + incx
                  jy = jy + incy
              enddo
          endif
      else
!
!        Form  A  when A is stored in the lower triangle.
!
          if ((incx.eq.1) .and. (incy.eq.1)) then
              do j = 1,n
                  if ((x(j).ne.zero) .or. (y(j).ne.zero)) then
                      temp1 = alpha*conjg(y(j))
                      temp2 = conjg(alpha*x(j))
                      a(j,j) = real(a(j,j)) + real(x(j)*temp1+y(j)*temp2)
                      do i = j + 1,n
                          a(i,j) = a(i,j) + x(i)*temp1 + y(i)*temp2
                      enddo
                  else
                      a(j,j) = real(a(j,j))
                  endif
              enddo
          else
              do j = 1,n
                  if ((x(jx).ne.zero) .or. (y(jy).ne.zero)) then
                      temp1 = alpha*conjg(y(jy))
                      temp2 = conjg(alpha*x(jx))
                      a(j,j) = real(a(j,j)) + real(x(jx)*temp1+y(jy)*temp2)
                      ix = jx
                      iy = jy
                      do i = j + 1,n
                          ix = ix + incx
                          iy = iy + incy
                          a(i,j) = a(i,j) + x(ix)*temp1 + y(iy)*temp2
                      enddo
                  else
                      a(j,j) = real(a(j,j))
                  endif
                  jx = jx + incx
                  jy = jy + incy
              enddo
          endif
      endif
!
      end subroutine cher2