Possible Solution

#include <float.h>
#include <math.h>
#include <stddef.h>
#include <stdio.h>
#include <stdlib.h>
#include <sys/times.h>
#include <unistd.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("%11.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    4
#endif

#ifndef DGEMM_MC
#define DGEMM_MC    256
#endif

#ifndef DGEMM_NC
#define DGEMM_NC    256
#endif

#ifndef DGEMM_KC
#define DGEMM_KC    512
#endif

void
dgepack_A(size_t m, size_t k,
          const double *A, ptrdiff_t incRowA, ptrdiff_t incColA,
          double *p)
{
    size_t mb = (m+DGEMM_MR-1)/DGEMM_MR;

    for (size_t l=0; l<k; ++l) {
        for (size_t i1=0; i1<mb; ++i1) {
            for (size_t i0=0; i0<DGEMM_MR; ++i0) {
                size_t i  = i1*DGEMM_MR + i0;
                size_t nu = i1*DGEMM_MR*k + l*DGEMM_MR + i0;
                p[nu] = (i<m) ? A[i*incRowA + l*incColA]
                              : 0;
            }
        }
    }
}

void
dgepack_B(size_t k, size_t n,
          const double *B, ptrdiff_t incRowB, ptrdiff_t incColB,
          double *p)
{
    size_t nb = (n+DGEMM_NR-1)/DGEMM_NR;

    for (size_t j1=0; j1<nb; ++j1) {
        for (size_t j0=0; j0<DGEMM_NR; ++j0) {
            for (size_t l=0; l<k; ++l) {
                size_t j  = j1*DGEMM_NR + j0;
                size_t nu = j1*DGEMM_NR*k + l*DGEMM_NR + j0;
                p[nu] = (j<n) ? B[l*incRowB + j*incColB]
                              : 0;
            }
        }
    }
}

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];
        }
    }
}

void
dgescal(size_t m, size_t n,
        double alpha,
        double *X, size_t incRowX, size_t incColX)
{
    if (alpha==1) {
        return;
    }
    if (alpha!=0) {
        for (size_t i=0; i<m; ++i) {
            for (size_t j=0; j<n; ++j) {
                X[i*incRowX+j*incColX] *= alpha;
            }
        }
    } else {
        for (size_t i=0; i<m; ++i) {
            for (size_t j=0; j<n; ++j) {
                X[i*incRowX+j*incColX] = 0;
            }
        }
    }
}

void
dgeaxpy(size_t m, size_t n,
        double alpha,
        const double *X, size_t incRowX, size_t incColX,
        double *Y, size_t incRowY, size_t incColY)
{
    if (alpha==0) {
        return;
    }
    for (size_t i=0; i<m; ++i) {
        for (size_t j=0; j<n; ++j) {
            Y[i*incRowY+j*incColY] += alpha*X[i*incRowX+j*incColX];
        }
    }
}

void
dgemm_macro(size_t m, size_t n, 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];

    size_t mb = (m+DGEMM_MR-1) / DGEMM_MR;
    size_t nb = (n+DGEMM_NR-1) / DGEMM_NR;

    size_t mr = m % DGEMM_MR;
    size_t nr = n % DGEMM_NR;

    for (size_t i=0; i<mb; ++i) {
        size_t m_ = (i<mb-1 || mr==0) ? DGEMM_MR
                                      : mr;
        for (size_t j=0; j<nb; ++j) {
            size_t n_ = (j<nb-1 || nr==0) ? DGEMM_NR
                                          : nr;
            if (m_==DGEMM_MR && n_==DGEMM_NR) {
                dgemm_micro(k, alpha,
                            &A[i*DGEMM_MR*k], &B[j*k*DGEMM_NR],
                            beta,
                            &C[i*DGEMM_MR*incRowC+j*DGEMM_NR*incColC],
                            incRowC, incColC);
            } else {
                dgemm_micro(k, alpha,
                            &A[i*DGEMM_MR*k], &B[j*k*DGEMM_NR],
                            0,
                            AB, 1, DGEMM_MR);
                dgescal(DGEMM_MR, DGEMM_NR,
                        beta,
                        &C[i*DGEMM_MR*incRowC+j*DGEMM_NR*incColC],
                        incRowC, incColC);
                dgeaxpy(DGEMM_MR, DGEMM_NR,
                        1,
                        AB, 1, DGEMM_MR,
                        &C[i*DGEMM_MR*incRowC+j*DGEMM_NR*incColC],
                        incRowC, incColC);
            }
        }
    }
}

void
dgemm_frame(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)
{
    if (k==0 || alpha==0) {
        dgescal(m, n, beta, C, incRowC, incColC);
        return;
    }

    size_t mb = (m+DGEMM_MC-1) / DGEMM_MC;
    size_t nb = (n+DGEMM_NC-1) / DGEMM_NC;
    size_t kb = (k+DGEMM_KC-1) / DGEMM_KC;

    size_t mr = m % DGEMM_MC;
    size_t nr = n % DGEMM_NC;
    size_t kr = k % DGEMM_KC;

    double      *Ap      = malloc(DGEMM_MC*DGEMM_KC*sizeof(*Ap));
    double      *Bp      = malloc(DGEMM_KC*DGEMM_NC*sizeof(*Bp));

    for (size_t j=0; j<nb; ++j) {
        size_t n_ = (j<nb-1 || nr==0) ? DGEMM_NC
                                      : nr;
        for (size_t l=0; l<kb; ++l) {
            size_t k_ = (l<kb-1 || kr==0) ? DGEMM_KC
                                          : kr;
            double beta_ = (l==0) ? beta
                                  : 1;

            dgepack_B(k_, n_,
                      &B[l*DGEMM_KC*incRowB+j*DGEMM_NC*incColB],
                      incRowB, incColB,
                      Bp);

            for (size_t i=0; i<mb; ++i) {
                size_t m_ = (i<mb-1 || mr==0) ? DGEMM_MC
                                              : mr;
                dgepack_A(m_, k_,
                          &A[i*DGEMM_MC*incRowA+l*DGEMM_KC*incColA],
                          incRowA, incColA,
                          Ap);
                dgemm_macro(m_, n_, k_,
                            alpha,
                            Ap, Bp,
                            beta_,
                            &C[i*DGEMM_MC*incRowC+j*DGEMM_NC*incColC],
                            incRowC, incColC);
            }
        }
    }

    free(Bp);
    free(Ap);
}

//-- Function for benchmarking and testing -------------------------------------

double
walltime()
{
   struct tms    ts;
   static double ClockTick=0.0;

   if (ClockTick==0.0) {
        ClockTick = 1.0 / ((double) sysconf(_SC_CLK_TCK));
   }
   return ((double) times(&ts)) * ClockTick;
}

void
randGeMatrix(size_t m, size_t n, double *A, ptrdiff_t incRowA, ptrdiff_t incColA)
{
    for (size_t j=0; j<n; ++j) {
        for (size_t i=0; i<m; ++i) {
            A[i*incRowA+j*incColA] = ((double)rand()-RAND_MAX/2)*200/RAND_MAX;
        }
    }
}

#define MIN(X,Y)   ((X)<(Y) ? (X) : (Y))
#define MAX(X,Y)   ((X)>(Y) ? (X) : (Y))

double
dgenrm1(size_t m, size_t n, const double *A, ptrdiff_t incRowA, ptrdiff_t incColA)
{
    double  result = 0;

    for (size_t j=0; j<n; ++j) {
        double sum = 0;
        for (size_t i=0; i<m; ++i) {
            sum += fabs(A[i*incRowA+j*incColA]);
        }
        if (sum>result) {
            result = sum;
        }
    }
    return result;
}

double
err_dgemm(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,
          const double *C0, ptrdiff_t incRowC0, ptrdiff_t incColC0,
          double *C, ptrdiff_t incRowC, ptrdiff_t incColC)
{
    double normA = dgenrm1(m, k, A, incRowA, incColA);
    double normB = dgenrm1(k, n, B, incRowB, incColB);
    double normC = dgenrm1(m, n, C, incRowC0, incColC0);
    double normD;
    size_t mn    = (m>n)  ? m  : n;
    size_t mnk   = (mn>k) ? mn : k;

    normA = MAX(normA, fabs(alpha)*normA);
    normC = MAX(normC, fabs(beta)*normC);

    dgeaxpy(m, n, -1.0, C0, incRowC0, incColC0, C, incRowC, incColC);
    normD = dgenrm1(m, n, C, incRowC, incColC);

    return normD/(mnk*normA*normB*normC);
}

void
dcopy(size_t n,
      const double *x, ptrdiff_t incX,
      double *y, ptrdiff_t incY)
{
    for (size_t i=0; i<n; ++i) {
        y[i*incY] = x[i*incX];
    }
}

void
dgecopy(size_t m, size_t n,
        const double *X, ptrdiff_t incRowX, ptrdiff_t incColX,
        double *Y, ptrdiff_t incRowY, ptrdiff_t incColY)
{
    if (incRowX<incColX) {
        for (size_t j=0; j<n; ++j) {
            dcopy(m, &X[j*incColX], incRowX, &Y[j*incColY], incRowY);
        }
    } else {
        for (size_t i=0; i<m; ++i) {
            dcopy(n, &X[i*incRowX], incColX, &Y[i*incRowY], incColY);
        }
    }
}


//------------------------------------------------------------------------------

#ifndef MIN_N
#define MIN_N 100
#endif

#ifndef MAX_N
#define MAX_N 4000
#endif

#ifndef INC_N
#define INC_N 100
#endif

#ifndef MIN_M
#define MIN_M 100
#endif

#ifndef MAX_M
#define MAX_M 4000
#endif

#ifndef INC_M
#define INC_M 100
#endif

#ifndef MIN_K
#define MIN_K 100
#endif

#ifndef MAX_K
#define MAX_K 4000
#endif

#ifndef INC_K
#define INC_K 100
#endif

#ifndef ALPHA
#define ALPHA 1
#endif

#ifndef BETA
#define BETA 1
#endif

#ifndef ROWMAJOR_A
#define ROWMAJOR_A 0
#endif

#ifndef ROWMAJOR_B
#define ROWMAJOR_B 0
#endif

#ifndef ROWMAJOR_C
#define ROWMAJOR_C 0
#endif

double A_[MAX_M*MAX_K];
double B_[MAX_K*MAX_N];
double C_[MAX_M*MAX_N];

double C0[MAX_M*MAX_N];     // reference solution
double C1[MAX_M*MAX_N];     // tested solution


int
main()
{
    randGeMatrix(MAX_M, MAX_K, A_, 1, MAX_M);
    randGeMatrix(MAX_K, MAX_N, B_, 1, MAX_K);
    randGeMatrix(MAX_N, MAX_M, C_, 1, MAX_M);

    printf("#%9s %9s %9s", "m", "n", "k");
    printf(" %12s %12s %17s", "t", "MFLOPS", "Residual Error");
    printf("\n");

    for (size_t m=MIN_M, n=MIN_N, k=MIN_K; n<=MAX_N && m<=MAX_M && k<=MAX_K;
         m+=INC_M, n+=INC_N, k+=INC_K)
    {
        double t, dt, err;
        size_t runs  = 1;
        double ops   = 2.0*m/1000*n/1000*k;

        ptrdiff_t incRowA = (ROWMAJOR_A==1) ? k : 1;
        ptrdiff_t incColA = (ROWMAJOR_A==1) ? 1 : m;

        ptrdiff_t incRowB = (ROWMAJOR_B==1) ? n : 1;
        ptrdiff_t incColB = (ROWMAJOR_B==1) ? 1 : k;

        ptrdiff_t incRowC = (ROWMAJOR_C==1) ? n : 1;
        ptrdiff_t incColC = (ROWMAJOR_C==1) ? 1 : m;

        printf(" %9zu %9zu %9td", m, n, k);

        // compute reference solution
        dgecopy(m, n, C_, 1, MAX_M, C0, incRowC, incColC);
        dgemm_ref(m, n, k,
                  ALPHA,
                  A_, incRowA, incRowA,
                  B_, incRowB, incColB,
                  BETA,
                  C0, incRowC, incColC);

        // benchmark dgemm_frame
        t    = 0;
        runs = 0;
        do {
            dgecopy(m, n, C_, 1, MAX_M, C1, incRowC, incColC);
            dt = walltime();
            dgemm_frame(m, n, k,
                        ALPHA,
                        A_, incRowA, incRowA,
                        B_, incRowB, incColB,
                        BETA,
                        C1, incRowC, incColC);
            dt = walltime() - dt;
            t += dt;
            ++runs;
        } while (t<0.3);
        t /= runs;

        err = err_dgemm(m, n, k,
                        ALPHA,
                        A_, incRowA, incColA,
                        B_, incRowB, incColB,
                        BETA,
                        C0, incRowC, incColC,
                        C1, incRowC, incColC);

        printf(" %12.2e %12.2lf %12.2e %4s", t, ops/t,
               err, (err<DBL_EPSILON) ? "PASS" : "FAIL");

        printf("\n");
    }

    return 0;
}