Go to main content
Oracle Developer Studio 12.5 Man Pages

Exit Print View

Updated: June 2017
 
 

sp_sgemm (3p)

Name

sp_sgemm - a SuperLU routine that 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' or op(X) = conjg(X'), alpha and beta are scalars, A is a sparse matrix of type SuperMatrix, and B and C are dense matrices, with op( A ) an m by k matrix, op( B ) a k by n matrix and C an m by n matrix.

Synopsis

#include <sunperf.h>


int sp_sgemm(char *transa, char *transb, int m, int  n,  int  k,  float
alpha,  SuperMatrix  *A, float *B, int ldb, float beta, float
*C, int ldc)


int sp_dgemm(char *transa, char *transb, int m, int n,  int  k,  double
alpha,  SuperMatrix *A, double *B, int ldb, double beta, dou-
ble *C, int ldc)


int sp_cgemm(char *transa, char *transb, int m, int n, int  k,  complex
alpha,  SuperMatrix  *A,  complex  *B, int ldb, complex beta,
complex *C, int ldc)


int sp_zgemm(char *transa, char *transb, int m, int n, int  k,  double-
complex  alpha,  SuperMatrix  *A,  doublecomplex *B, int ldb,
doublecomplex beta, doublecomplex *C, int ldc)


long sp_sgemm_64(char *transa, char *transb, long m, long  n,  long  k,
float  alpha,  SuperMatrix_64  *A,  float *B, long ldb, float
beta, float *C, long ldc)


long sp_dgemm_64(char *transa, char *transb, long m, long  n,  long  k,
double  alpha, SuperMatrix_64 *A, double *B, long ldb, double
beta, double *C, long ldc)


long sp_cgemm_64(char *transa, char *transb, long m, long  n,  long  k,
complex  alpha, SuperMatrix_64 *A, complex *B, long ldb, com-
plex beta, complex *C, long ldc)


long sp_zgemm_64(char *transa, char *transb, long m, long  n,  long  k,
doublecomplex  alpha,  SuperMatrix_64  *A,  doublecomplex *B,
long ldb, doublecomplex beta, doublecomplex *C, long ldc)

Description

Oracle Solaris Studio Performance Library                          sp_gemm(3P)



NAME
       sp_gemm:  sp_cgemm,  sp_dgemm,  sp_sgemm,  sp_zgemm - a SuperLU routine
       that 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' or op(X) = conjg(X'),

       alpha and beta are scalars, A is a sparse matrix of  type  SuperMatrix,
       and B and C are dense matrices, with op( A ) an m by k matrix,  op( B )
       a  k by n matrix and  C an m by n matrix.


SYNOPSIS
       #include <sunperf.h>


       int sp_sgemm(char *transa, char *transb, int m, int  n,  int  k,  float
                 alpha,  SuperMatrix  *A, float *B, int ldb, float beta, float
                 *C, int ldc)


       int sp_dgemm(char *transa, char *transb, int m, int n,  int  k,  double
                 alpha,  SuperMatrix *A, double *B, int ldb, double beta, dou-
                 ble *C, int ldc)


       int sp_cgemm(char *transa, char *transb, int m, int n, int  k,  complex
                 alpha,  SuperMatrix  *A,  complex  *B, int ldb, complex beta,
                 complex *C, int ldc)


       int sp_zgemm(char *transa, char *transb, int m, int n, int  k,  double-
                 complex  alpha,  SuperMatrix  *A,  doublecomplex *B, int ldb,
                 doublecomplex beta, doublecomplex *C, int ldc)


       long sp_sgemm_64(char *transa, char *transb, long m, long  n,  long  k,
                 float  alpha,  SuperMatrix_64  *A,  float *B, long ldb, float
                 beta, float *C, long ldc)


       long sp_dgemm_64(char *transa, char *transb, long m, long  n,  long  k,
                 double  alpha, SuperMatrix_64 *A, double *B, long ldb, double
                 beta, double *C, long ldc)


       long sp_cgemm_64(char *transa, char *transb, long m, long  n,  long  k,
                 complex  alpha, SuperMatrix_64 *A, complex *B, long ldb, com-
                 plex beta, complex *C, long ldc)


       long sp_zgemm_64(char *transa, char *transb, long m, long  n,  long  k,
                 doublecomplex  alpha,  SuperMatrix_64  *A,  doublecomplex *B,
                 long ldb, doublecomplex beta, doublecomplex *C, long ldc)


PURPOSE
       sp_gemm 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',

       alpha and beta are scalars; A is sparse and of type SuperMatrix; B  and
       C are dense matrices, with op( A ) an m by k matrix, op( B )  a  k by n
       matrix and  C an m by n matrix.

       sp_gemm returns 0 on exit.


ARGUMENTS
       char *transa (input)
              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'.
              transa = 'C' or 'c',  op( A ) = A'.


       char *transb (input)
              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'.
              transb = 'C' or 'c',  op( B ) = B'.


       int m (input)
              On entry,  m  specifies  the number  of rows  of the  matrix op(
              A )  and of the  matrix  C.  m  must  be at least  zero.


       int n (input)
              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.


       int k (input)
              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.


       float alpha (input)
              On entry, alpha specifies the scaling value of matrix A.


       SuperMatrix *A (input)
              General  matrix  A  in  sparse  format with dimensions (A->nrow,
              A->ncol).  Currently, the type of A can be:
              Stype = NC or NCP; Dtype = SLU_C; Mtype = GE.
              In the future, more general A can be handled.


       int lda (input)
              On entry, lda specifies the first dimension of A as declared  in
              the  calling  routine.  When   transa  =  'N' or 'n' then lda >=
              max(1, m), otherwise  lda >= max(1, k).


       float *B (input)
              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 con-
              tain  the matrix  B,  otherwise the leading  n by k  part of the
              array  B  must contain  the matrix B.


       int ldb (input)
              On entry, ldb specifies the first dimension of B as declared  in
              the  calling  routine.  When   transb  =  'N' or 'n' then ldb >=
              max(1, k); otherwise ldb >= max(1, n).


       float beta (input)
              On entry,  beta  specifies the scaling value of matrix C.   When
              beta   is  supplied  as  zero  then C need not be initialized on
              input.


       float *C (input/output)
              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).


       int ldc (input)
              On  entry, ldc specifies the first dimension of C as declared in
              the  calling routine.  ldc >= max(1, m).


SEE ALSO
       SuperMatrix

       http://crd.lbl.gov/~xiaoye/SuperLU/

       James W. Demmel, Stanley C. Eisenstat, John R. Gilbert,  Xiaoye  S.  Li
       and  Joseph  W. H. Liu, "A supernodal approach to sparse partial pivot-
       ing", SIAM J. Matrix Analysis and Applications, Vol. 20, Num. 3,  1999,
       pp. 720-755.



                                  7 Nov 2015                       sp_gemm(3P)