LAPACK  3.9.1
LAPACK: Linear Algebra PACKage

◆ sgemm()

subroutine sgemm ( character  TRANSA,
character  TRANSB,
integer  M,
integer  N,
integer  K,
real  ALPHA,
real, dimension(lda,*)  A,
integer  LDA,
real, dimension(ldb,*)  B,
integer  LDB,
real  BETA,
real, dimension(ldc,*)  C,
integer  LDC 
)

SGEMM

Purpose:
 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.
Parameters
[in]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.
[in]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.
[in]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.
[in]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.
[in]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.
[in]ALPHA
          ALPHA is REAL
           On entry, ALPHA specifies the scalar alpha.
[in]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.
[in]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 ).
[in]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.
[in]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 ).
[in]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.
[in,out]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 ).
[in]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 ).
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
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.

Definition at line 186 of file sgemm.f.

187 *
188 * -- Reference BLAS level3 routine --
189 * -- Reference BLAS is a software package provided by Univ. of Tennessee, --
190 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
191 *
192 * .. Scalar Arguments ..
193  REAL ALPHA,BETA
194  INTEGER K,LDA,LDB,LDC,M,N
195  CHARACTER TRANSA,TRANSB
196 * ..
197 * .. Array Arguments ..
198  REAL A(LDA,*),B(LDB,*),C(LDC,*)
199 * ..
200 *
201 * =====================================================================
202 *
203 * .. External Functions ..
204  LOGICAL LSAME
205  EXTERNAL lsame
206 * ..
207 * .. External Subroutines ..
208  EXTERNAL xerbla
209 * ..
210 * .. Intrinsic Functions ..
211  INTRINSIC max
212 * ..
213 * .. Local Scalars ..
214  REAL TEMP
215  INTEGER I,INFO,J,L,NROWA,NROWB
216  LOGICAL NOTA,NOTB
217 * ..
218 * .. Parameters ..
219  REAL ONE,ZERO
220  parameter(one=1.0e+0,zero=0.0e+0)
221 * ..
222 *
223 * Set NOTA and NOTB as true if A and B respectively are not
224 * transposed and set NROWA and NROWB as the number of rows of A
225 * and B respectively.
226 *
227  nota = lsame(transa,'N')
228  notb = lsame(transb,'N')
229  IF (nota) THEN
230  nrowa = m
231  ELSE
232  nrowa = k
233  END IF
234  IF (notb) THEN
235  nrowb = k
236  ELSE
237  nrowb = n
238  END IF
239 *
240 * Test the input parameters.
241 *
242  info = 0
243  IF ((.NOT.nota) .AND. (.NOT.lsame(transa,'C')) .AND.
244  + (.NOT.lsame(transa,'T'))) THEN
245  info = 1
246  ELSE IF ((.NOT.notb) .AND. (.NOT.lsame(transb,'C')) .AND.
247  + (.NOT.lsame(transb,'T'))) THEN
248  info = 2
249  ELSE IF (m.LT.0) THEN
250  info = 3
251  ELSE IF (n.LT.0) THEN
252  info = 4
253  ELSE IF (k.LT.0) THEN
254  info = 5
255  ELSE IF (lda.LT.max(1,nrowa)) THEN
256  info = 8
257  ELSE IF (ldb.LT.max(1,nrowb)) THEN
258  info = 10
259  ELSE IF (ldc.LT.max(1,m)) THEN
260  info = 13
261  END IF
262  IF (info.NE.0) THEN
263  CALL xerbla('SGEMM ',info)
264  RETURN
265  END IF
266 *
267 * Quick return if possible.
268 *
269  IF ((m.EQ.0) .OR. (n.EQ.0) .OR.
270  + (((alpha.EQ.zero).OR. (k.EQ.0)).AND. (beta.EQ.one))) RETURN
271 *
272 * And if alpha.eq.zero.
273 *
274  IF (alpha.EQ.zero) THEN
275  IF (beta.EQ.zero) THEN
276  DO 20 j = 1,n
277  DO 10 i = 1,m
278  c(i,j) = zero
279  10 CONTINUE
280  20 CONTINUE
281  ELSE
282  DO 40 j = 1,n
283  DO 30 i = 1,m
284  c(i,j) = beta*c(i,j)
285  30 CONTINUE
286  40 CONTINUE
287  END IF
288  RETURN
289  END IF
290 *
291 * Start the operations.
292 *
293  IF (notb) THEN
294  IF (nota) THEN
295 *
296 * Form C := alpha*A*B + beta*C.
297 *
298  DO 90 j = 1,n
299  IF (beta.EQ.zero) THEN
300  DO 50 i = 1,m
301  c(i,j) = zero
302  50 CONTINUE
303  ELSE IF (beta.NE.one) THEN
304  DO 60 i = 1,m
305  c(i,j) = beta*c(i,j)
306  60 CONTINUE
307  END IF
308  DO 80 l = 1,k
309  temp = alpha*b(l,j)
310  DO 70 i = 1,m
311  c(i,j) = c(i,j) + temp*a(i,l)
312  70 CONTINUE
313  80 CONTINUE
314  90 CONTINUE
315  ELSE
316 *
317 * Form C := alpha*A**T*B + beta*C
318 *
319  DO 120 j = 1,n
320  DO 110 i = 1,m
321  temp = zero
322  DO 100 l = 1,k
323  temp = temp + a(l,i)*b(l,j)
324  100 CONTINUE
325  IF (beta.EQ.zero) THEN
326  c(i,j) = alpha*temp
327  ELSE
328  c(i,j) = alpha*temp + beta*c(i,j)
329  END IF
330  110 CONTINUE
331  120 CONTINUE
332  END IF
333  ELSE
334  IF (nota) THEN
335 *
336 * Form C := alpha*A*B**T + beta*C
337 *
338  DO 170 j = 1,n
339  IF (beta.EQ.zero) THEN
340  DO 130 i = 1,m
341  c(i,j) = zero
342  130 CONTINUE
343  ELSE IF (beta.NE.one) THEN
344  DO 140 i = 1,m
345  c(i,j) = beta*c(i,j)
346  140 CONTINUE
347  END IF
348  DO 160 l = 1,k
349  temp = alpha*b(j,l)
350  DO 150 i = 1,m
351  c(i,j) = c(i,j) + temp*a(i,l)
352  150 CONTINUE
353  160 CONTINUE
354  170 CONTINUE
355  ELSE
356 *
357 * Form C := alpha*A**T*B**T + beta*C
358 *
359  DO 200 j = 1,n
360  DO 190 i = 1,m
361  temp = zero
362  DO 180 l = 1,k
363  temp = temp + a(l,i)*b(j,l)
364  180 CONTINUE
365  IF (beta.EQ.zero) THEN
366  c(i,j) = alpha*temp
367  ELSE
368  c(i,j) = alpha*temp + beta*c(i,j)
369  END IF
370  190 CONTINUE
371  200 CONTINUE
372  END IF
373  END IF
374 *
375  RETURN
376 *
377 * End of SGEMM .
378 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
Here is the call graph for this function: