sgemm(3f) - [BLAS:SINGLE_BLAS_LEVEL3]
C:=alphaAB+beta*C, A, B, C rectangular.
subroutine sgemm(transa,transb,m,n,k,alpha,a,lda,b,ldb,beta,c,ldc)
.. Scalar Arguments ..
real,intent(in) :: alpha,beta
integer,intent(in) :: k,lda,ldb,ldc,m,n
character,intent(in) :: transa,transb
..
.. Array Arguments ..
real,intent(in) :: a(lda,*),b(ldb,*)
real,intent(inout) :: c(ldc,*)
..
SGEMM 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,
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.
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**T.
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**T.
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 REAL
On entry, ALPHA specifies the scalar alpha.
A
A is REAL 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 REAL 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 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 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 ).
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.
Online html documentation available at
http://www.netlib.org/lapack/explore-html/
Type | Intent | Optional | 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 | |||
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 |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public | :: | i | ||||
integer, | public | :: | info | ||||
integer, | public | :: | j | ||||
integer, | public | :: | l | ||||
logical, | public | :: | nota | ||||
logical, | public | :: | notb | ||||
integer, | public | :: | nrowa | ||||
integer, | public | :: | nrowb | ||||
real, | public, | parameter | :: | one | = | 1.0e+0 | |
real, | public | :: | temp | ||||
real, | public, | parameter | :: | zero | = | 0.0e+0 |
subroutine sgemm(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 ..
real,intent(in) :: alpha,beta
integer,intent(in) :: k,lda,ldb,ldc,m,n
character,intent(in) :: transa,transb
! ..
! .. 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 temp
integer i,info,j,l,nrowa,nrowb
logical nota,notb
! ..
! .. Parameters ..
real one,zero
parameter (one=1.0e+0,zero=0.0e+0)
! ..
!
! Set NOTA and NOTB as true if A and B respectively are not
! transposed and set NROWA and NROWB as the number of rows of A
! and B respectively.
!
nota = lsame(transa,'N')
notb = lsame(transb,'N')
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.lsame(transa,'C')) .and. (.not.lsame(transa,'T'))) then
info = 1
elseif ((.not.notb) .and. (.not.lsame(transb,'C')) .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('SGEMM ',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 if 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
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
else
if (nota) then
!
! 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
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 sgemm