dsyr Subroutine

public subroutine dsyr(uplo, n, alpha, x, incx, a, lda)

NAME

dsyr(3f) - [BLAS:DOUBLE_BLAS_LEVEL3]

SYNOPSIS

 subroutine dsyr(uplo,n,alpha,x,incx,a,lda)

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

DEFINITION

DSYR performs the symmetric rank 1 operation

 A := alpha*x*x**T + A,

where alpha is a real scalar, x is an n element vector and A is an n by n symmetric 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 DOUBLE PRECISION.
        On entry, ALPHA specifies the scalar alpha.

X

       X is DOUBLE PRECISION 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.

A

       A is DOUBLE PRECISION 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 symmetric 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 symmetric 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.

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

\ingroup double_blas_level2

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
double precision, intent(in) :: alpha
double precision, intent(in) :: x(*)
integer, intent(in) :: incx
double precision, 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 :: j
integer, public :: jx
integer, public :: kx
double precision, public :: temp
double precision, public, parameter :: zero = 0.0d+0

Source Code

       subroutine dsyr(uplo,n,alpha,x,incx,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 ..
      double precision,intent(in)       :: alpha
      integer,intent(in)                :: incx,lda,n
      character,intent(in)              :: uplo
!     ..
!     .. Array Arguments ..
      double precision,intent(inout)    :: a(lda,*)
      double precision,intent(in)       :: x(*)
!     ..
!
!  =====================================================================
!
!     .. Parameters ..
      double precision zero
      parameter (zero=0.0d+0)
!     ..
!     .. Local Scalars ..
      double precision temp
      integer i,info,ix,j,jx,kx
!     ..
!     .. External Functions ..
!     ..
!     .. External Subroutines ..
!     ..
!     .. Intrinsic Functions ..
      intrinsic 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 (incx.eq.0) then
          info = 5
      elseif (lda.lt.max(1,n)) then
          info = 7
      endif
      if (info.ne.0) then
          call xerbla('DSYR  ',info)
          return
      endif
!
!     Quick return if possible.
!
      if ((n.eq.0) .or. (alpha.eq.zero)) return
!
!     Set the start point in X if the increment is not unity.
!
      if (incx.le.0) then
          kx = 1 - (n-1)*incx
      elseif (incx.ne.1) then
          kx = 1
      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 upper triangle.
!
          if (incx.eq.1) then
              do j = 1,n
                  if (x(j).ne.zero) then
                      temp = alpha*x(j)
                      a(1:j,j) = a(1:j,j) + x(1:j)*temp
                  endif
              enddo
          else
              jx = kx
              do j = 1,n
                  if (x(jx).ne.zero) then
                      temp = alpha*x(jx)
                      ix = kx
                      do i = 1,j
                          a(i,j) = a(i,j) + x(ix)*temp
                          ix = ix + incx
                      enddo
                  endif
                  jx = jx + incx
              enddo
          endif
      else
!
!        Form  A  when A is stored in lower triangle.
!
          if (incx.eq.1) then
              do j = 1,n
                  if (x(j).ne.zero) then
                      temp = alpha*x(j)
                      a(j:n,j) = a(j:n,j) + x(j:n)*temp
                  endif
              enddo
          else
              jx = kx
              do j = 1,n
                  if (x(jx).ne.zero) then
                      temp = alpha*x(jx)
                      ix = jx
                      do i = j,n
                          a(i,j) = a(i,j) + x(ix)*temp
                          ix = ix + incx
                      enddo
                  endif
                  jx = jx + incx
              enddo
          endif
      endif
!
      end subroutine dsyr