LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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 of 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 of 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 of 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.
Date
November 2015
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 189 of file sgemm.f.

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

Here is the call graph for this function: