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
#include <stddef.h>
#include <stdio.h>
#include <stdlib.h>

void
initGeMatrix(size_t m, size_t n,
             double *A,
             ptrdiff_t incRowA, ptrdiff_t incColA)
{
    for (size_t i=0; i<m; ++i) {
        for (size_t j=0; j<n; ++j) {
            A[i*incRowA+j*incColA] = i*n + j + 1;
        }
    }
}

void
printGeMatrix(size_t m, size_t n,
              const double *A,
              ptrdiff_t incRowA, ptrdiff_t incColA)
{
    for (size_t i=0; i<m; ++i) {
        for (size_t j=0; j<n; ++j) {
            printf("%10.4lf ", A[i*incRowA+j*incColA]);
        }
        printf("\n");
    }
    printf("\n");
}

void
dgemm_ref(size_t m, size_t n, size_t k,
          double alpha,
          const double *A, ptrdiff_t incRowA, ptrdiff_t incColA,
          const double *B, ptrdiff_t incRowB, ptrdiff_t incColB,
          double beta,
          double *C, ptrdiff_t incRowC, ptrdiff_t incColC)
{
    size_t i, j, l;

    if (beta!=1) {
        if (beta!=0) {
            for (i=0; i<m; ++i) {
                for (j=0; j<n; ++j) {
                    C[i*incRowC+j*incColC] *= beta;
                }
            }
        } else {
            for (i=0; i<m; ++i) {
                for (j=0; j<n; ++j) {
                    C[i*incRowC+j*incColC] = 0;
                }
            }
        }
    }
    if (alpha!=0) {
        for (i=0; i<m; ++i) {
            for (j=0; j<n; ++j) {
                for (l=0; l<k; ++l) {
                    C[i*incRowC+j*incColC] += alpha*A[i*incRowA+l*incColA]
                                                   *B[l*incRowB+j*incColB];
                }
            }
        }
    }
}


#ifndef DGEMM_MR
#define DGEMM_MR    4
#endif

#ifndef DGEMM_NR
#define DGEMM_NR    5
#endif

void
dgemm_micro(size_t k, double alpha,
            const double *A, const double *B,
            double beta,
            double *C, ptrdiff_t incRowC, ptrdiff_t incColC)
{
    double AB[DGEMM_MR*DGEMM_NR];

    // AB <- A*B
    for (size_t i=0; i<DGEMM_MR*DGEMM_NR; ++i) {
        AB[i] = 0;
    }
    for (size_t l=0; l<k; ++l) {
        for (size_t i=0; i<DGEMM_MR; ++i) {
            for (size_t j=0; j<DGEMM_NR; ++j) {
                AB[i+j*DGEMM_MR] += A[i+l*DGEMM_MR]*B[l*DGEMM_NR+j];
            }
        }
    }
    // C <- beta*C
    if (beta!=1) {
        if (beta!=0) {
            for (size_t i=0; i<DGEMM_MR; ++i) {
                for (size_t j=0; j<DGEMM_NR; ++j) {
                    C[i*incRowC+j*incColC] *= beta;
                }
            }
        } else {
            for (size_t i=0; i<DGEMM_MR; ++i) {
                for (size_t j=0; j<DGEMM_NR; ++j) {
                    C[i*incRowC+j*incColC] = 0;
                }
            }
        }
    }
    // C <- C + alpha*AB
    for (size_t i=0; i<DGEMM_MR; ++i) {
        for (size_t j=0; j<DGEMM_NR; ++j) {
            C[i*incRowC+j*incColC] += alpha*AB[i+j*DGEMM_MR];
        }
    }
}


int
main()
{
    size_t      k       = 7;

    double      *A      = malloc(DGEMM_MR*k*sizeof(*A));
    double      *B      = malloc(k*DGEMM_NR*sizeof(*B));

    double      C0[DGEMM_MR*DGEMM_NR];
    double      C1[DGEMM_MR*DGEMM_NR];

    initGeMatrix(DGEMM_MR, k, A, 1, DGEMM_MR);
    initGeMatrix(k, DGEMM_NR, B, DGEMM_NR, 1);
    initGeMatrix(DGEMM_MR, DGEMM_NR, C0, 1, DGEMM_MR);
    initGeMatrix(DGEMM_MR, DGEMM_NR, C1, 1, DGEMM_MR);

    printf("A=\n");
    printGeMatrix(DGEMM_MR, k, A, 1, DGEMM_MR);
    printf("B=\n");
    printGeMatrix(k, DGEMM_NR, B, DGEMM_NR, 1);
    printf("C=\n");
    printGeMatrix(DGEMM_MR, DGEMM_NR, C0, 1, DGEMM_MR);

    double alpha = 1;
    double beta  = 1;
    dgemm_ref(DGEMM_MR, DGEMM_NR, k,
              alpha,
              A, 1, DGEMM_MR,
              B, DGEMM_NR, 1,
              beta,
              C0, 1, DGEMM_MR);

    dgemm_micro(k, alpha, A, B, beta, C1, 1, DGEMM_MR);


    printf("gemm_ref computed C=\n");
    printGeMatrix(DGEMM_MR, DGEMM_NR, C0, 1, DGEMM_MR);

    printf("gemm_micro computed C=\n");
    printGeMatrix(DGEMM_MR, DGEMM_NR, C1, 1, DGEMM_MR);

    free(A);
    free(B);
}