|
123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384 |
- *> \brief \b SGEMM
- *
- * =========== DOCUMENTATION ===========
- *
- * Online html documentation available at
- * http://www.netlib.org/lapack/explore-html/
- *
- * Definition:
- * ===========
- *
- * SUBROUTINE SGEMM(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
- *
- * .. Scalar Arguments ..
- * REAL ALPHA,BETA
- * INTEGER K,LDA,LDB,LDC,M,N
- * CHARACTER TRANSA,TRANSB
- * ..
- * .. Array Arguments ..
- * REAL A(LDA,*),B(LDB,*),C(LDC,*)
- * ..
- *
- *
- *> \par Purpose:
- * =============
- *>
- *> \verbatim
- *>
- *> 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.
- *> \endverbatim
- *
- * Arguments:
- * ==========
- *
- *> \param[in] TRANSA
- *> \verbatim
- *> 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.
- *> \endverbatim
- *>
- *> \param[in] TRANSB
- *> \verbatim
- *> 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.
- *> \endverbatim
- *>
- *> \param[in] M
- *> \verbatim
- *> 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.
- *> \endverbatim
- *>
- *> \param[in] N
- *> \verbatim
- *> 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.
- *> \endverbatim
- *>
- *> \param[in] K
- *> \verbatim
- *> 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.
- *> \endverbatim
- *>
- *> \param[in] ALPHA
- *> \verbatim
- *> ALPHA is REAL
- *> On entry, ALPHA specifies the scalar alpha.
- *> \endverbatim
- *>
- *> \param[in] A
- *> \verbatim
- *> 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.
- *> \endverbatim
- *>
- *> \param[in] LDA
- *> \verbatim
- *> 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 ).
- *> \endverbatim
- *>
- *> \param[in] B
- *> \verbatim
- *> 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.
- *> \endverbatim
- *>
- *> \param[in] LDB
- *> \verbatim
- *> 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 ).
- *> \endverbatim
- *>
- *> \param[in] BETA
- *> \verbatim
- *> BETA is REAL
- *> On entry, BETA specifies the scalar beta. When BETA is
- *> supplied as zero then C need not be set on input.
- *> \endverbatim
- *>
- *> \param[in,out] C
- *> \verbatim
- *> 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 ).
- *> \endverbatim
- *>
- *> \param[in] LDC
- *> \verbatim
- *> 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 ).
- *> \endverbatim
- *
- * Authors:
- * ========
- *
- *> \author Univ. of Tennessee
- *> \author Univ. of California Berkeley
- *> \author Univ. of Colorado Denver
- *> \author NAG Ltd.
- *
- *> \date December 2016
- *
- *> \ingroup single_blas_level3
- *
- *> \par Further Details:
- * =====================
- *>
- *> \verbatim
- *>
- *> 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.
- *> \endverbatim
- *>
- * =====================================================================
- SUBROUTINE SGEMM(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
- *
- * -- 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 ALPHA,BETA
- INTEGER K,LDA,LDB,LDC,M,N
- CHARACTER TRANSA,TRANSB
- * ..
- * .. Array Arguments ..
- REAL A(LDA,*),B(LDB,*),C(LDC,*)
- * ..
- *
- * =====================================================================
- *
- * .. External Functions ..
- LOGICAL LSAME
- EXTERNAL LSAME
- * ..
- * .. External Subroutines ..
- EXTERNAL XERBLA
- * ..
- * .. Intrinsic Functions ..
- INTRINSIC MAX
- * ..
- * .. Local Scalars ..
- REAL TEMP
- INTEGER I,INFO,J,L,NCOLA,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, NCOLA and NROWB as the number of rows
- * and columns of A and the number of rows of B respectively.
- *
- NOTA = LSAME(TRANSA,'N')
- NOTB = LSAME(TRANSB,'N')
- IF (NOTA) THEN
- NROWA = M
- NCOLA = K
- ELSE
- NROWA = K
- NCOLA = M
- END IF
- IF (NOTB) THEN
- NROWB = K
- ELSE
- NROWB = N
- END IF
- *
- * Test the input parameters.
- *
- INFO = 0
- IF ((.NOT.NOTA) .AND. (.NOT.LSAME(TRANSA,'C')) .AND.
- + (.NOT.LSAME(TRANSA,'T'))) THEN
- INFO = 1
- ELSE IF ((.NOT.NOTB) .AND. (.NOT.LSAME(TRANSB,'C')) .AND.
- + (.NOT.LSAME(TRANSB,'T'))) THEN
- INFO = 2
- ELSE IF (M.LT.0) THEN
- INFO = 3
- ELSE IF (N.LT.0) THEN
- INFO = 4
- ELSE IF (K.LT.0) THEN
- INFO = 5
- ELSE IF (LDA.LT.MAX(1,NROWA)) THEN
- INFO = 8
- ELSE IF (LDB.LT.MAX(1,NROWB)) THEN
- INFO = 10
- ELSE IF (LDC.LT.MAX(1,M)) THEN
- INFO = 13
- END IF
- IF (INFO.NE.0) THEN
- CALL XERBLA('SGEMM ',INFO)
- RETURN
- END IF
- *
- * 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 20 J = 1,N
- DO 10 I = 1,M
- C(I,J) = ZERO
- 10 CONTINUE
- 20 CONTINUE
- ELSE
- DO 40 J = 1,N
- DO 30 I = 1,M
- C(I,J) = BETA*C(I,J)
- 30 CONTINUE
- 40 CONTINUE
- END IF
- RETURN
- END IF
- *
- * Start the operations.
- *
- IF (NOTB) THEN
- IF (NOTA) THEN
- *
- * Form C := alpha*A*B + beta*C.
- *
- DO 90 J = 1,N
- IF (BETA.EQ.ZERO) THEN
- DO 50 I = 1,M
- C(I,J) = ZERO
- 50 CONTINUE
- ELSE IF (BETA.NE.ONE) THEN
- DO 60 I = 1,M
- C(I,J) = BETA*C(I,J)
- 60 CONTINUE
- END IF
- DO 80 L = 1,K
- TEMP = ALPHA*B(L,J)
- DO 70 I = 1,M
- C(I,J) = C(I,J) + TEMP*A(I,L)
- 70 CONTINUE
- 80 CONTINUE
- 90 CONTINUE
- ELSE
- *
- * Form C := alpha*A**T*B + beta*C
- *
- DO 120 J = 1,N
- DO 110 I = 1,M
- TEMP = ZERO
- DO 100 L = 1,K
- TEMP = TEMP + A(L,I)*B(L,J)
- 100 CONTINUE
- IF (BETA.EQ.ZERO) THEN
- C(I,J) = ALPHA*TEMP
- ELSE
- C(I,J) = ALPHA*TEMP + BETA*C(I,J)
- END IF
- 110 CONTINUE
- 120 CONTINUE
- END IF
- ELSE
- IF (NOTA) THEN
- *
- * Form C := alpha*A*B**T + beta*C
- *
- DO 170 J = 1,N
- IF (BETA.EQ.ZERO) THEN
- DO 130 I = 1,M
- C(I,J) = ZERO
- 130 CONTINUE
- ELSE IF (BETA.NE.ONE) THEN
- DO 140 I = 1,M
- C(I,J) = BETA*C(I,J)
- 140 CONTINUE
- END IF
- DO 160 L = 1,K
- TEMP = ALPHA*B(J,L)
- DO 150 I = 1,M
- C(I,J) = C(I,J) + TEMP*A(I,L)
- 150 CONTINUE
- 160 CONTINUE
- 170 CONTINUE
- ELSE
- *
- * Form C := alpha*A**T*B**T + beta*C
- *
- DO 200 J = 1,N
- DO 190 I = 1,M
- TEMP = ZERO
- DO 180 L = 1,K
- TEMP = TEMP + A(L,I)*B(J,L)
- 180 CONTINUE
- IF (BETA.EQ.ZERO) THEN
- C(I,J) = ALPHA*TEMP
- ELSE
- C(I,J) = ALPHA*TEMP + BETA*C(I,J)
- END IF
- 190 CONTINUE
- 200 CONTINUE
- END IF
- END IF
- *
- RETURN
- *
- * End of SGEMM .
- *
- END
|