LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ sgemmtr()

subroutine sgemmtr ( character uplo,
character transa,
character transb,
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 )

SGEMMTR

Purpose:
!>
!> SGEMMTR  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 n by k matrix,  op( B )  a  k by n matrix and  C an n by n matrix.
!> Thereby, the routine only accesses and updates the upper or lower
!> triangular part of the result matrix C. This behaviour can be used if
!> the resulting matrix C is known to be symmetric.
!> 
Parameters
[in]UPLO
!>          UPLO is CHARACTER*1
!>           On entry, UPLO specifies whether the lower or the upper
!>           triangular part of C is access and updated.
!>
!>              UPLO = 'L' or 'l', the lower triangular part of C is used.
!>
!>              UPLO = 'U' or 'u', the upper triangular part of C is used.
!> 
[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]N
!>          N is INTEGER
!>           On entry,  N specifies the number of rows and columns of
!>           the matrix C, the number of columns of op(B) and the number
!>           of rows of op(A).  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  n  otherwise.
!>           Before entry with  TRANSA = 'N' or 'n',  the leading  n 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, n ), 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  n 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 upper or lower triangular part of the matrix
!>           C  is overwritten by the n 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, n ).
!> 
Author
Martin Koehler
Further Details:
!>
!>  Level 3 Blas routine.
!>
!>  -- Written on 19-July-2023.
!>     Martin Koehler, MPI Magdeburg
!> 

Definition at line 189 of file sgemmtr.f.

191 IMPLICIT NONE
192*
193* -- Reference BLAS level3 routine --
194* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
195* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
196*
197* .. Scalar Arguments ..
198 REAL ALPHA,BETA
199 INTEGER K,LDA,LDB,LDC,N
200 CHARACTER TRANSA,TRANSB,UPLO
201* ..
202* .. Array Arguments ..
203 REAL A(LDA,*),B(LDB,*),C(LDC,*)
204* ..
205*
206* =====================================================================
207*
208* .. External Functions ..
209 LOGICAL LSAME
210 EXTERNAL lsame
211* ..
212* .. External Subroutines ..
213 EXTERNAL xerbla
214* ..
215* .. Intrinsic Functions ..
216 INTRINSIC max
217* ..
218* .. Local Scalars ..
219 REAL TEMP
220 INTEGER I,INFO,J,L,NROWA,NROWB, ISTART, ISTOP
221 LOGICAL NOTA,NOTB, UPPER
222* ..
223* .. Parameters ..
224 REAL ONE,ZERO
225 parameter(one=1.0d+0,zero=0.0d+0)
226* ..
227*
228* Set NOTA and NOTB as true if A and B respectively are not
229* transposed and set NROWA and NROWB as the number of rows of A
230* and B respectively.
231*
232 nota = lsame(transa,'N')
233 notb = lsame(transb,'N')
234 IF (nota) THEN
235 nrowa = n
236 ELSE
237 nrowa = k
238 END IF
239 IF (notb) THEN
240 nrowb = k
241 ELSE
242 nrowb = n
243 END IF
244 upper = lsame(uplo, 'U')
245*
246* Test the input parameters.
247*
248 info = 0
249 IF ((.NOT. upper) .AND. (.NOT. lsame(uplo, 'L'))) THEN
250 info = 1
251 ELSE IF ((.NOT.nota) .AND. (.NOT.lsame(transa,'C')) .AND.
252 + (.NOT.lsame(transa,'T'))) THEN
253 info = 2
254 ELSE IF ((.NOT.notb) .AND. (.NOT.lsame(transb,'C')) .AND.
255 + (.NOT.lsame(transb,'T'))) THEN
256 info = 3
257 ELSE IF (n.LT.0) THEN
258 info = 4
259 ELSE IF (k.LT.0) THEN
260 info = 5
261 ELSE IF (lda.LT.max(1,nrowa)) THEN
262 info = 8
263 ELSE IF (ldb.LT.max(1,nrowb)) THEN
264 info = 10
265 ELSE IF (ldc.LT.max(1,n)) THEN
266 info = 13
267 END IF
268 IF (info.NE.0) THEN
269 CALL xerbla('SGEMMTR',info)
270 RETURN
271 END IF
272*
273* Quick return if possible.
274*
275 IF (n.EQ.0) 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 IF (upper) THEN
283 istart = 1
284 istop = j
285 ELSE
286 istart = j
287 istop = n
288 END IF
289
290 DO 10 i = istart, istop
291 c(i,j) = zero
292 10 CONTINUE
293 20 CONTINUE
294 ELSE
295 DO 40 j = 1,n
296 IF (upper) THEN
297 istart = 1
298 istop = j
299 ELSE
300 istart = j
301 istop = n
302 END IF
303
304 DO 30 i = istart, istop
305 c(i,j) = beta*c(i,j)
306 30 CONTINUE
307 40 CONTINUE
308 END IF
309 RETURN
310 END IF
311*
312* Start the operations.
313*
314 IF (notb) THEN
315 IF (nota) THEN
316*
317* Form C := alpha*A*B + beta*C.
318*
319 DO 90 j = 1,n
320 IF (upper) THEN
321 istart = 1
322 istop = j
323 ELSE
324 istart = j
325 istop = n
326 END IF
327 IF (beta.EQ.zero) THEN
328 DO 50 i = istart, istop
329 c(i,j) = zero
330 50 CONTINUE
331 ELSE IF (beta.NE.one) THEN
332 DO 60 i = istart, istop
333 c(i,j) = beta*c(i,j)
334 60 CONTINUE
335 END IF
336 DO 80 l = 1,k
337 temp = alpha*b(l,j)
338 DO 70 i = istart, istop
339 c(i,j) = c(i,j) + temp*a(i,l)
340 70 CONTINUE
341 80 CONTINUE
342 90 CONTINUE
343 ELSE
344*
345* Form C := alpha*A**T*B + beta*C
346*
347 DO 120 j = 1,n
348 IF (upper) THEN
349 istart = 1
350 istop = j
351 ELSE
352 istart = j
353 istop = n
354 END IF
355
356 DO 110 i = istart, istop
357 temp = zero
358 DO 100 l = 1,k
359 temp = temp + a(l,i)*b(l,j)
360 100 CONTINUE
361 IF (beta.EQ.zero) THEN
362 c(i,j) = alpha*temp
363 ELSE
364 c(i,j) = alpha*temp + beta*c(i,j)
365 END IF
366 110 CONTINUE
367 120 CONTINUE
368 END IF
369 ELSE
370 IF (nota) THEN
371*
372* Form C := alpha*A*B**T + beta*C
373*
374 DO 170 j = 1,n
375 IF (upper) THEN
376 istart = 1
377 istop = j
378 ELSE
379 istart = j
380 istop = n
381 END IF
382
383 IF (beta.EQ.zero) THEN
384 DO 130 i = istart,istop
385 c(i,j) = zero
386 130 CONTINUE
387 ELSE IF (beta.NE.one) THEN
388 DO 140 i = istart,istop
389 c(i,j) = beta*c(i,j)
390 140 CONTINUE
391 END IF
392 DO 160 l = 1,k
393 temp = alpha*b(j,l)
394 DO 150 i = istart,istop
395 c(i,j) = c(i,j) + temp*a(i,l)
396 150 CONTINUE
397 160 CONTINUE
398 170 CONTINUE
399 ELSE
400*
401* Form C := alpha*A**T*B**T + beta*C
402*
403 DO 200 j = 1,n
404 IF (upper) THEN
405 istart = 1
406 istop = j
407 ELSE
408 istart = j
409 istop = n
410 END IF
411
412 DO 190 i = istart, istop
413 temp = zero
414 DO 180 l = 1,k
415 temp = temp + a(l,i)*b(j,l)
416 180 CONTINUE
417 IF (beta.EQ.zero) THEN
418 c(i,j) = alpha*temp
419 ELSE
420 c(i,j) = alpha*temp + beta*c(i,j)
421 END IF
422 190 CONTINUE
423 200 CONTINUE
424 END IF
425 END IF
426*
427 RETURN
428*
429* End of SGEMMTR
430*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: