1
       2
       3
       4
       5
       6
       7
       8
       9
      10
      11
      12
      13
      14
      15
      16
      17
      18
      19
      20
      21
      22
      23
      24
      25
      26
      27
      28
      29
      30
      31
      32
      33
      34
      35
      36
      37
      38
      39
      40
      41
      42
      43
      44
      45
      46
      47
      48
      49
      50
      51
      52
      53
      54
      55
      56
      57
      58
      59
      60
      61
      62
      63
      64
      65
      66
      67
      68
      69
      70
      71
      72
      73
      74
      75
      76
      77
      78
      79
      80
      81
      82
      83
      84
      85
      86
      87
      88
      89
      90
      91
      92
      93
      94
      95
      96
      97
      98
      99
     100
     101
     102
     103
     104
     105
     106
     107
     108
     109
     110
     111
     112
     113
     114
     115
     116
     117
     118
     119
     120
     121
     122
     123
     124
     125
     126
     127
     128
     129
     130
     131
     132
     133
     134
     135
     136
     137
     138
     139
     140
     141
     142
     143
     144
     145
     146
     147
     148
     149
     150
     151
     152
     153
     154
     155
     156
     157
     158
     159
     160
     161
     162
     163
     164
     165
     166
     167
     168
     169
     170
#include <ulmblas.h>
#include <auxiliary/xerbla.h>
#include <level3/dgemm_nn.h>

void
ULMBLAS(dgemm)(const enum Trans  transA,
               const enum Trans  transB,
               const int         m,
               const int         n,
               const int         k,
               const double      alpha,
               const double      *A,
               const int         ldA,
               const double      *B,
               const int         ldB,
               const double      beta,
               double            *C,
               const int         ldC)
{
//
//  Local scalars
//
    int     i, j;

//
//  Quick return if possible.
//
    if (m==0 || n==0 || ((alpha==0.0 || k==0) && (beta==1.0))) {
        return;
    }

//
//  And if alpha is exactly zero
//
    if (alpha==0.0) {
        if (beta==0.0) {
            for (j=0; j<n; ++j) {
                for (i=0; i<m; ++i) {
                    C[i+j*ldC] = 0.0;
                }
            }
        } else {
            for (j=0; j<n; ++j) {
                for (i=0; i<m; ++i) {
                    C[i+j*ldC] *= beta;
                }
            }
        }
        return;
    }

//
//  Start the operations.
//
    if (transB==NoTrans || transB==Conj) {
        if (transA==NoTrans || transA==Conj) {
//
//          Form  C := alpha*A*B + beta*C.
//
            dgemm_nn(m, n, k,
                     alpha,
                     A, 1, ldA,
                     B, 1, ldB,
                     beta,
                     C, 1, ldC);
        } else {
//
//          Form  C := alpha*A**T*B + beta*C
//
            dgemm_nn(m, n, k,
                     alpha,
                     A, ldA, 1,
                     B, 1, ldB,
                     beta,
                     C, 1, ldC);
        }
    } else {
        if (transA==NoTrans || transA==Conj) {
//
//          Form  C := alpha*A*B**T + beta*C
//
            dgemm_nn(m, n, k,
                     alpha,
                     A, 1, ldA,
                     B, ldB, 1,
                     beta,
                     C, 1, ldC);
        } else {
//
//          Form  C := alpha*A**T*B**T + beta*C
//
            dgemm_nn(m, n, k,
                     alpha,
                     A, ldA, 1,
                     B, ldB, 1,
                     beta,
                     C, 1, ldC);
        }
    }
}

void
F77BLAS(dgemm)(const char     *_transA,
               const char     *_transB,
               const int      *_m,
               const int      *_n,
               const int      *_k,
               const double   *_alpha,
               const double   *A,
               const int      *_ldA,
               const double   *B,
               const int      *_ldB,
               const double   *_beta,
               double         *C,
               const int      *_ldC)
{
//
//  Dereference scalar parameters
//
    enum Trans transA = charToTranspose(*_transA);
    enum Trans transB = charToTranspose(*_transB);
    int m             = *_m;
    int n             = *_n;
    int k             = *_k;
    double alpha      = *_alpha;
    int ldA           = *_ldA;
    int ldB           = *_ldB;
    double beta       = *_beta;
    int ldC           = *_ldC;

//
//  Set  numRowsA and numRowsB as the number of rows of A and B
//
    int numRowsA = (transA==NoTrans || transA==Conj) ? m : k;
    int numRowsB = (transB==NoTrans || transB==Conj) ? k : n;

//
//  Test the input parameters
//
    int info = 0;
    if (transA==0) {
        info = 1;
    } else if (transB==0) {
        info = 2;
    } else if (m<0) {
        info = 3;
    } else if (n<0) {
        info = 4;
    } else if (k<0) {
        info = 5;
    } else if (ldA<max(1,numRowsA)) {
        info = 8;
    } else if (ldB<max(1,numRowsB)) {
        info = 10;
    } else if (ldC<max(1,m)) {
        info = 13;
    }

    if (info!=0) {
        F77BLAS(xerbla)("DGEMM ", &info);
    }

    ULMBLAS(dgemm)(transA, transB,
                   m, n, k,
                   alpha,
                   A, ldA,
                   B, ldB,
                   beta,
                   C, ldC);
}