zgemm Subroutine

public subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)

NAME

zgemm(3f) - [BLAS:COMPLEX16_BLAS_LEVEL3]

SYNOPSIS

 subroutine zgemm(transa,transb,m,n,k,alpha,a,lda,b,ldb,beta,c,ldc)

   .. Scalar Arguments ..
   complex(kind=real64),intent(in)    :: alpha,beta
   integer,intent(in)                 :: k,lda,ldb,ldc,m,n
   character,intent(in)               :: transa,transb
   ..
   .. Array Arguments ..
   complex(kind=real64),intent(in)    :: a(lda,*),b(ldb,*)
   complex(kind=real64),intent(inout) :: c(ldc,*)
   ..

DEFINITION

ZGEMM performs one of the matrix-matrix operations

 C := alpha*op( A )*op( B ) + beta*C,

where op( X ) is one of

 op( X ) = X   or   op( X ) = X**T   or   op( X ) = X**H,

alpha and beta are scalars, and A, B and C are matrices, with op( A ) an m by k matrix, op( B ) a k by n matrix and C an m by n matrix.

OPTIONS

TRANSA

       TRANSA is CHARACTER*1
        On entry, TRANSA specifies the form of op( A ) to be used in
        the matrix multiplication as follows:

           TRANSA = 'N' or 'n',  op( A ) = A.

           TRANSA = 'T' or 't',  op( A ) = A**T.

           TRANSA = 'C' or 'c',  op( A ) = A**H.

TRANSB

       TRANSB is CHARACTER*1
        On entry, TRANSB specifies the form of op( B ) to be used in
        the matrix multiplication as follows:

           TRANSB = 'N' or 'n',  op( B ) = B.

           TRANSB = 'T' or 't',  op( B ) = B**T.

           TRANSB = 'C' or 'c',  op( B ) = B**H.

M

       M is INTEGER
        On entry, M specifies the number of rows of the matrix
        op( A ) and 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
        op( B ) and the number of columns of the matrix C. N must be
        at least zero.

K

       K is INTEGER
        On entry, K specifies the number of columns of the matrix
        op( A ) and the number of rows of the matrix op( B ). K 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, ka ), where ka is
        k when TRANSA = 'N' or 'n', and is m otherwise.
        Before entry with TRANSA = 'N' or 'n', the leading m by k
        part of the array A must contain the matrix A, otherwise
        the leading k by m part of the array A must contain the
        matrix A.

LDA

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

B

       B is complex(kind=real64) array, dimension ( LDB, kb ), where kb is
        n when TRANSB = 'N' or 'n', and is k otherwise.
        Before entry with TRANSB = 'N' or 'n', the leading k by n
        part of the array B must contain the matrix B, otherwise
        the leading n by k 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. When TRANSB = 'N' or 'n' then
        LDB must be at least max( 1, k ), otherwise LDB must be at
        least max( 1, n ).

BETA

       BETA is complex(kind=real64)
        On entry, BETA specifies the scalar beta. When BETA is
        supplied as zero then C need not be set on input.

C

       C is complex(kind=real64) 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 matrix
        ( alpha*op( A )*op( B ) + beta*C ).

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) :: transa
character(len=1), intent(in) :: transb
integer, intent(in) :: m
integer, intent(in) :: n
integer, intent(in) :: k
complex(kind=real64), intent(in) :: alpha
complex(kind=real64), intent(in) :: a(lda,*)
integer, intent(in) :: lda
complex(kind=real64), intent(in) :: b(ldb,*)
integer, intent(in) :: ldb
complex(kind=real64), intent(in) :: beta
complex(kind=real64), intent(inout) :: c(ldc,*)
integer, intent(in) :: ldc

Contents

Source Code


Variables

Type Visibility Attributes Name Initial
logical, public :: conja
logical, public :: conjb
integer, public :: i
integer, public :: info
integer, public :: j
integer, public :: l
logical, public :: nota
logical, public :: notb
integer, public :: nrowa
integer, public :: nrowb
complex(kind=real64), public, parameter :: one = (1.0d+0,0.0d+0)
complex(kind=real64), public :: temp
complex(kind=real64), public, parameter :: zero = (0.0d+0,0.0d+0)

Source Code

       subroutine zgemm(transa,transb,m,n,k,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 ..
      complex(kind=real64),intent(in)    :: alpha,beta
      integer,intent(in)                 :: k,lda,ldb,ldc,m,n
      character,intent(in)               :: transa,transb
!     ..
!     .. Array Arguments ..
      complex(kind=real64),intent(in)    :: a(lda,*),b(ldb,*)
      complex(kind=real64),intent(inout) :: c(ldc,*)
!     ..
!
!  =====================================================================
!
!     .. External Functions ..
!      LOGICAL LSAME
!      EXTERNAL LSAME
!     ..
!     .. External Subroutines ..
!      EXTERNAL XERBLA
!     ..
!     .. Intrinsic Functions ..
      intrinsic dconjg,max
!     ..
!     .. Local Scalars ..
      complex(kind=real64) :: temp
      integer i,info,j,l,nrowa,nrowb
      logical conja,conjb,nota,notb
!     ..
!     .. 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))
!     ..
!
!     Set  NOTA  and  NOTB  as  true if  A  and  B  respectively are not
!     conjugated or transposed, set  CONJA and CONJB  as true if  A  and
!     B  respectively are to be  transposed but  not conjugated  and set
!     NROWA and NROWB  as the number of rows  of  A  and  B  respectively.
!
      nota = lsame(transa,'N')
      notb = lsame(transb,'N')
      conja = lsame(transa,'C')
      conjb = lsame(transb,'C')
      if (nota) then
          nrowa = m
      else
          nrowa = k
      endif
      if (notb) then
          nrowb = k
      else
          nrowb = n
      endif
!
!     Test the input parameters.
!
      info = 0
      if ((.not.nota) .and. (.not.conja) .and.  (.not.lsame(transa,'T'))) then
          info = 1
      elseif ((.not.notb) .and. (.not.conjb) .and.  (.not.lsame(transb,'T'))) then
          info = 2
      elseif (m.lt.0) then
          info = 3
      elseif (n.lt.0) then
          info = 4
      elseif (k.lt.0) then
          info = 5
      elseif (lda.lt.max(1,nrowa)) then
          info = 8
      elseif (ldb.lt.max(1,nrowb)) then
          info = 10
      elseif (ldc.lt.max(1,m)) then
          info = 13
      endif
      if (info.ne.0) then
          call xerbla('ZGEMM ',info)
          return
      endif
!
!     Quick return if possible.
!
      if ((m.eq.0) .or. (n.eq.0) .or.  (((alpha.eq.zero).or. (k.eq.0)).and. (beta.eq.one))) return
!
!     And when  alpha.eq.zero.
!
      if (alpha.eq.zero) then
          if (beta.eq.zero) then
              do j = 1,n
                  c(1:m,j) = zero
              enddo
          else
              do j = 1,n
                  c(1:m,j) = beta*c(1:m,j)
              enddo
          endif
          return
      endif
!
!     Start the operations.
!
      if (notb) then
          if (nota) then
!
!           Form  C := alpha*A*B + beta*C.
!
              do j = 1,n
                  if (beta.eq.zero) then
                      c(1:m,j) = zero
                  elseif (beta.ne.one) then
                      c(1:m,j) = beta*c(1:m,j)
                  endif
                  do l = 1,k
                      temp = alpha*b(l,j)
                      c(1:m,j) = c(1:m,j) + temp*a(1:m,l)
                  enddo
              enddo
          elseif (conja) then
!
!           Form  C := alpha*A**H*B + beta*C.
!
              do j = 1,n
                  do i = 1,m
                      temp = zero
                      do l = 1,k
                          temp = temp + dconjg(a(l,i))*b(l,j)
                      enddo
                      if (beta.eq.zero) then
                          c(i,j) = alpha*temp
                      else
                          c(i,j) = alpha*temp + beta*c(i,j)
                      endif
                  enddo
              enddo
          else
!
!           Form  C := alpha*A**T*B + beta*C
!
              do j = 1,n
                  do i = 1,m
                      temp = zero
                      do l = 1,k
                          temp = temp + a(l,i)*b(l,j)
                      enddo
                      if (beta.eq.zero) then
                          c(i,j) = alpha*temp
                      else
                          c(i,j) = alpha*temp + beta*c(i,j)
                      endif
                  enddo
              enddo
          endif
      elseif (nota) then
          if (conjb) then
!
!           Form  C := alpha*A*B**H + beta*C.
!
              do j = 1,n
                  if (beta.eq.zero) then
                      c(1:m,j) = zero
                  elseif (beta.ne.one) then
                      c(1:m,j) = beta*c(1:m,j)
                  endif
                  do l = 1,k
                      temp = alpha*dconjg(b(j,l))
                      do i = 1,m
                          c(i,j) = c(i,j) + temp*a(i,l)
                      enddo
                  enddo
              enddo
          else
!
!           Form  C := alpha*A*B**T + beta*C
!
              do j = 1,n
                  if (beta.eq.zero) then
                      c(1:m,j) = zero
                  elseif (beta.ne.one) then
                      c(1:m,j) = beta*c(1:m,j)
                  endif
                  do l = 1,k
                      temp = alpha*b(j,l)
                      c(1:m,j) = c(1:m,j) + temp*a(1:m,l)
                  enddo
              enddo
          endif
      elseif (conja) then
          if (conjb) then
!
!           Form  C := alpha*A**H*B**H + beta*C.
!
              do j = 1,n
                  do i = 1,m
                      temp = zero
                      do l = 1,k
                          temp = temp + dconjg(a(l,i))*dconjg(b(j,l))
                      enddo
                      if (beta.eq.zero) then
                          c(i,j) = alpha*temp
                      else
                          c(i,j) = alpha*temp + beta*c(i,j)
                      endif
                  enddo
              enddo
          else
!
!           Form  C := alpha*A**H*B**T + beta*C
!
              do j = 1,n
                  do i = 1,m
                      temp = zero
                      do l = 1,k
                          temp = temp + dconjg(a(l,i))*b(j,l)
                      enddo
                      if (beta.eq.zero) then
                          c(i,j) = alpha*temp
                      else
                          c(i,j) = alpha*temp + beta*c(i,j)
                      endif
                  enddo
              enddo
          endif
      else
          if (conjb) then
!
!           Form  C := alpha*A**T*B**H + beta*C
!
              do j = 1,n
                  do i = 1,m
                      temp = zero
                      do l = 1,k
                          temp = temp + a(l,i)*dconjg(b(j,l))
                      enddo
                      if (beta.eq.zero) then
                          c(i,j) = alpha*temp
                      else
                          c(i,j) = alpha*temp + beta*c(i,j)
                      endif
                  enddo
              enddo
          else
!
!           Form  C := alpha*A**T*B**T + beta*C
!
              do j = 1,n
                  do i = 1,m
                      temp = zero
                      do l = 1,k
                          temp = temp + a(l,i)*b(j,l)
                      enddo
                      if (beta.eq.zero) then
                          c(i,j) = alpha*temp
                      else
                          c(i,j) = alpha*temp + beta*c(i,j)
                      endif
                  enddo
              enddo
          endif
      endif
!
      end subroutine zgemm