ssymm Subroutine

public subroutine ssymm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)

NAME

ssymm(3f) - [BLAS:SINGLE_BLAS_LEVEL3] C:=alpha*A*B+beta*C, A symmetric, B, C rectangular.

SYNOPSIS

 subroutine ssymm(side,uplo,m,n,alpha,a,lda,b,ldb,beta,c,ldc)

   .. Scalar Arguments ..
   real,intent(in)      :: alpha,beta
   integer,intent(in)   :: lda,ldb,ldc,m,n
   character,intent(in) :: side,uplo
   ..
   .. Array Arguments ..
   real,intent(in)      :: a(lda,*),b(ldb,*)
   real,intent(inout)   :: c(ldc,*)
   ..

DEFINITION

SSYMM performs one of the matrix-matrix operations

 C := alpha*A*B + beta*C,

or

 C := alpha*B*A + beta*C,

where alpha and beta are scalars, A is a symmetric matrix and B and C are m by n matrices.

OPTIONS

SIDE

       SIDE is CHARACTER*1
        On entry, SIDE specifies whether the symmetric matrix A
        appears on the left or right in the operation as follows:

           SIDE = 'L' or 'l'   C := alpha*A*B + beta*C,

           SIDE = 'R' or 'r'   C := alpha*B*A + beta*C,

UPLO

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

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

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

M

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

N

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

ALPHA

       ALPHA is REAL
        On entry, ALPHA specifies the scalar alpha.

A

       A is REAL array, dimension ( LDA, ka ), where ka is
        m when SIDE = 'L' or 'l' and is n otherwise.
        Before entry with SIDE = 'L' or 'l', the m by m part of
        the array A must contain the symmetric matrix, such that
        when UPLO = 'U' or 'u', the leading m by m 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, and when UPLO = 'L' or 'l',
        the leading m by m 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.
        Before entry with SIDE = 'R' or 'r', the n by n part of
        the array A must contain the symmetric matrix, such that
        when 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, and when 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.

LDA

       LDA is INTEGER
        On entry, LDA specifies the first dimension of A as declared
        in the calling (sub) program. When SIDE = 'L' or 'l' then
        LDA must be at least max( 1, m ), otherwise LDA must be at
        least max( 1, n ).

B

       B is REAL array, dimension ( LDB, N )
        Before entry, the leading m by n part of the array B must
        contain the matrix B.

LDB

       LDB is INTEGER
        On entry, LDB specifies the first dimension of B as declared
        in the calling (sub) program. LDB must be at least
        max( 1, m ).

BETA

       BETA is REAL
        On entry, BETA specifies the scalar beta. When BETA is
        supplied as zero then C need not be set on input.

C

       C is REAL array, dimension ( LDC, N )
        Before entry, the leading m by n part of the array C must
        contain the matrix C, except when beta is zero, in which
        case C need not be set on entry.
        On exit, the array C is overwritten by the m by n updated
        matrix.

LDC

       LDC is INTEGER
        On entry, LDC specifies the first dimension of C as declared
        in the calling (sub) program. LDC must be at least
        max( 1, m ).

AUTHORS

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

date:December 2016

FURTHER DETAILS

Level 3 Blas routine.

– Written on 8-February-1989. Jack Dongarra, Argonne National Laboratory. Iain Duff, AERE Harwell. Jeremy Du Croz, Numerical Algorithms Group Ltd. Sven Hammarling, Numerical Algorithms Group Ltd.

SEE ALSO

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

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: side
character(len=1), intent(in) :: uplo
integer, intent(in) :: m
integer, intent(in) :: n
real, intent(in) :: alpha
real, intent(in) :: a(lda,*)
integer, intent(in) :: lda
real, intent(in) :: b(ldb,*)
integer, intent(in) :: ldb
real, intent(in) :: beta
real, intent(inout) :: c(ldc,*)
integer, intent(in) :: ldc

Contents

Source Code


Variables

Type Visibility Attributes Name Initial
integer, public :: i
integer, public :: info
integer, public :: j
integer, public :: k
integer, public :: nrowa
real, public, parameter :: one = 1.0e+0
real, public :: temp1
real, public :: temp2
logical, public :: upper
real, public, parameter :: zero = 0.0e+0

Source Code

       subroutine ssymm(side,uplo,m,n,alpha,a,lda,b,ldb,beta,c,ldc)
      implicit none
!
!  -- Reference BLAS level3 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 ..
      real,intent(in)      :: alpha,beta
      integer,intent(in)   :: lda,ldb,ldc,m,n
      character,intent(in) :: side,uplo
!     ..
!     .. Array Arguments ..
      real,intent(in)      :: a(lda,*),b(ldb,*)
      real,intent(inout)   :: c(ldc,*)
!     ..
!
!  =====================================================================
!
!     .. External Functions ..  LOGICAL LSAME
!     ..
!     .. External Subroutines ..  EXTERNAL XERBLA
!     ..
!     .. Intrinsic Functions ..
      intrinsic max
!     ..
!     .. Local Scalars ..
      real temp1,temp2
      integer i,info,j,k,nrowa
      logical upper
!     ..
!     .. Parameters ..
      real,parameter :: one=1.0e+0,zero=0.0e+0
!     ..
!
!     Set NROWA as the number of rows of A.
!
      if (lsame(side,'L')) then
          nrowa = m
      else
          nrowa = n
      endif
      upper = lsame(uplo,'U')
!
!     Test the input parameters.
!
      info = 0
      if ((.not.lsame(side,'L')) .and. (.not.lsame(side,'R'))) then
          info = 1
      elseif ((.not.upper) .and. (.not.lsame(uplo,'L'))) then
          info = 2
      elseif (m.lt.0) then
          info = 3
      elseif (n.lt.0) then
          info = 4
      elseif (lda.lt.max(1,nrowa)) then
          info = 7
      elseif (ldb.lt.max(1,m)) then
          info = 9
      elseif (ldc.lt.max(1,m)) then
          info = 12
      endif
      if (info.ne.0) then
          call xerbla('SSYMM ',info)
          return
      endif
!
!     Quick return if possible.
!
      if ((m.eq.0) .or. (n.eq.0) .or.  ((alpha.eq.zero).and. (beta.eq.one))) return
!
!     And when  alpha.eq.zero.
!
      if (alpha.eq.zero) then
          if (beta.eq.zero) then
              c(1:m,1:n) = zero
          else
              c(1:m,1:n) = beta*c(1:m,1:n)
          endif
          return
      endif
!
!     Start the operations.
!
      if (lsame(side,'L')) then
!
!        Form  C := alpha*A*B + beta*C.
!
          if (upper) then
              do j = 1,n
                  do i = 1,m
                      temp1 = alpha*b(i,j)
                      temp2 = zero
                      do k = 1,i - 1
                          c(k,j) = c(k,j) + temp1*a(k,i)
                          temp2 = temp2 + b(k,j)*a(k,i)
                      enddo
                      if (beta.eq.zero) then
                          c(i,j) = temp1*a(i,i) + alpha*temp2
                      else
                          c(i,j) = beta*c(i,j) + temp1*a(i,i) + alpha*temp2
                      endif
                  enddo
              enddo
          else
              do j = 1,n
                  do i = m,1,-1
                      temp1 = alpha*b(i,j)
                      temp2 = zero
                      do k = i + 1,m
                          c(k,j) = c(k,j) + temp1*a(k,i)
                          temp2 = temp2 + b(k,j)*a(k,i)
                      enddo
                      if (beta.eq.zero) then
                          c(i,j) = temp1*a(i,i) + alpha*temp2
                      else
                          c(i,j) = beta*c(i,j) + temp1*a(i,i) + alpha*temp2
                      endif
                  enddo
              enddo
          endif
      else
!
!        Form  C := alpha*B*A + beta*C.
!
          do j = 1,n
              temp1 = alpha*a(j,j)
              if (beta.eq.zero) then
                  do i = 1,m
                      c(i,j) = temp1*b(i,j)
                  enddo
              else
                  do i = 1,m
                      c(i,j) = beta*c(i,j) + temp1*b(i,j)
                  enddo
              endif
              do k = 1,j - 1
                  if (upper) then
                      temp1 = alpha*a(k,j)
                  else
                      temp1 = alpha*a(j,k)
                  endif
                  do i = 1,m
                      c(i,j) = c(i,j) + temp1*b(i,k)
                  enddo
              enddo
              do k = j + 1,n
                  if (upper) then
                      temp1 = alpha*a(j,k)
                  else
                      temp1 = alpha*a(k,j)
                  endif
                  do i = 1,m
                      c(i,j) = c(i,j) + temp1*b(i,k)
                  enddo
              enddo
          enddo
      endif

      end subroutine ssymm