LAPACK  3.11.0
LAPACK: Linear Algebra PACKage
sgemm.f
1 *> \brief \b SGEMM
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * SUBROUTINE SGEMM(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
12 *
13 * .. Scalar Arguments ..
14 * REAL ALPHA,BETA
15 * INTEGER K,LDA,LDB,LDC,M,N
16 * CHARACTER TRANSA,TRANSB
17 * ..
18 * .. Array Arguments ..
19 * REAL A(LDA,*),B(LDB,*),C(LDC,*)
20 * ..
21 *
22 *
23 *> \par Purpose:
24 * =============
25 *>
26 *> \verbatim
27 *>
28 *> SGEMM performs one of the matrix-matrix operations
29 *>
30 *> C := alpha*op( A )*op( B ) + beta*C,
31 *>
32 *> where op( X ) is one of
33 *>
34 *> op( X ) = X or op( X ) = X**T,
35 *>
36 *> alpha and beta are scalars, and A, B and C are matrices, with op( A )
37 *> an m by k matrix, op( B ) a k by n matrix and C an m by n matrix.
38 *> \endverbatim
39 *
40 * Arguments:
41 * ==========
42 *
43 *> \param[in] TRANSA
44 *> \verbatim
45 *> TRANSA is CHARACTER*1
46 *> On entry, TRANSA specifies the form of op( A ) to be used in
47 *> the matrix multiplication as follows:
48 *>
49 *> TRANSA = 'N' or 'n', op( A ) = A.
50 *>
51 *> TRANSA = 'T' or 't', op( A ) = A**T.
52 *>
53 *> TRANSA = 'C' or 'c', op( A ) = A**T.
54 *> \endverbatim
55 *>
56 *> \param[in] TRANSB
57 *> \verbatim
58 *> TRANSB is CHARACTER*1
59 *> On entry, TRANSB specifies the form of op( B ) to be used in
60 *> the matrix multiplication as follows:
61 *>
62 *> TRANSB = 'N' or 'n', op( B ) = B.
63 *>
64 *> TRANSB = 'T' or 't', op( B ) = B**T.
65 *>
66 *> TRANSB = 'C' or 'c', op( B ) = B**T.
67 *> \endverbatim
68 *>
69 *> \param[in] M
70 *> \verbatim
71 *> M is INTEGER
72 *> On entry, M specifies the number of rows of the matrix
73 *> op( A ) and of the matrix C. M must be at least zero.
74 *> \endverbatim
75 *>
76 *> \param[in] N
77 *> \verbatim
78 *> N is INTEGER
79 *> On entry, N specifies the number of columns of the matrix
80 *> op( B ) and the number of columns of the matrix C. N must be
81 *> at least zero.
82 *> \endverbatim
83 *>
84 *> \param[in] K
85 *> \verbatim
86 *> K is INTEGER
87 *> On entry, K specifies the number of columns of the matrix
88 *> op( A ) and the number of rows of the matrix op( B ). K must
89 *> be at least zero.
90 *> \endverbatim
91 *>
92 *> \param[in] ALPHA
93 *> \verbatim
94 *> ALPHA is REAL
95 *> On entry, ALPHA specifies the scalar alpha.
96 *> \endverbatim
97 *>
98 *> \param[in] A
99 *> \verbatim
100 *> A is REAL array, dimension ( LDA, ka ), where ka is
101 *> k when TRANSA = 'N' or 'n', and is m otherwise.
102 *> Before entry with TRANSA = 'N' or 'n', the leading m by k
103 *> part of the array A must contain the matrix A, otherwise
104 *> the leading k by m part of the array A must contain the
105 *> matrix A.
106 *> \endverbatim
107 *>
108 *> \param[in] LDA
109 *> \verbatim
110 *> LDA is INTEGER
111 *> On entry, LDA specifies the first dimension of A as declared
112 *> in the calling (sub) program. When TRANSA = 'N' or 'n' then
113 *> LDA must be at least max( 1, m ), otherwise LDA must be at
114 *> least max( 1, k ).
115 *> \endverbatim
116 *>
117 *> \param[in] B
118 *> \verbatim
119 *> B is REAL array, dimension ( LDB, kb ), where kb is
120 *> n when TRANSB = 'N' or 'n', and is k otherwise.
121 *> Before entry with TRANSB = 'N' or 'n', the leading k by n
122 *> part of the array B must contain the matrix B, otherwise
123 *> the leading n by k part of the array B must contain the
124 *> matrix B.
125 *> \endverbatim
126 *>
127 *> \param[in] LDB
128 *> \verbatim
129 *> LDB is INTEGER
130 *> On entry, LDB specifies the first dimension of B as declared
131 *> in the calling (sub) program. When TRANSB = 'N' or 'n' then
132 *> LDB must be at least max( 1, k ), otherwise LDB must be at
133 *> least max( 1, n ).
134 *> \endverbatim
135 *>
136 *> \param[in] BETA
137 *> \verbatim
138 *> BETA is REAL
139 *> On entry, BETA specifies the scalar beta. When BETA is
140 *> supplied as zero then C need not be set on input.
141 *> \endverbatim
142 *>
143 *> \param[in,out] C
144 *> \verbatim
145 *> C is REAL array, dimension ( LDC, N )
146 *> Before entry, the leading m by n part of the array C must
147 *> contain the matrix C, except when beta is zero, in which
148 *> case C need not be set on entry.
149 *> On exit, the array C is overwritten by the m by n matrix
150 *> ( alpha*op( A )*op( B ) + beta*C ).
151 *> \endverbatim
152 *>
153 *> \param[in] LDC
154 *> \verbatim
155 *> LDC is INTEGER
156 *> On entry, LDC specifies the first dimension of C as declared
157 *> in the calling (sub) program. LDC must be at least
158 *> max( 1, m ).
159 *> \endverbatim
160 *
161 * Authors:
162 * ========
163 *
164 *> \author Univ. of Tennessee
165 *> \author Univ. of California Berkeley
166 *> \author Univ. of Colorado Denver
167 *> \author NAG Ltd.
168 *
169 *> \ingroup gemm
170 *
171 *> \par Further Details:
172 * =====================
173 *>
174 *> \verbatim
175 *>
176 *> Level 3 Blas routine.
177 *>
178 *> -- Written on 8-February-1989.
179 *> Jack Dongarra, Argonne National Laboratory.
180 *> Iain Duff, AERE Harwell.
181 *> Jeremy Du Croz, Numerical Algorithms Group Ltd.
182 *> Sven Hammarling, Numerical Algorithms Group Ltd.
183 *> \endverbatim
184 *>
185 * =====================================================================
186  SUBROUTINE sgemm(TRANSA,TRANSB,M,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
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 *
379  END
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine sgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
SGEMM
Definition: sgemm.f:187