2010-02-22 18:09:57 +08:00
|
|
|
|
|
|
|
// g++-4.4 bench_gemm.cpp -I .. -O2 -DNDEBUG -lrt -fopenmp && OMP_NUM_THREADS=2 ./a.out
|
|
|
|
// icpc bench_gemm.cpp -I .. -O3 -DNDEBUG -lrt -openmp && OMP_NUM_THREADS=2 ./a.out
|
|
|
|
|
2010-03-05 16:57:04 +08:00
|
|
|
#include <iostream>
|
2010-06-22 17:10:38 +08:00
|
|
|
#include <Eigen/Core>
|
2010-02-23 00:57:15 +08:00
|
|
|
#include <bench/BenchTimer.h>
|
2010-02-22 18:09:57 +08:00
|
|
|
|
|
|
|
using namespace std;
|
|
|
|
using namespace Eigen;
|
|
|
|
|
|
|
|
#ifndef SCALAR
|
2010-07-22 22:29:35 +08:00
|
|
|
// #define SCALAR std::complex<float>
|
2010-07-16 02:45:45 +08:00
|
|
|
#define SCALAR float
|
2010-02-22 18:09:57 +08:00
|
|
|
#endif
|
|
|
|
|
|
|
|
typedef SCALAR Scalar;
|
2010-07-08 01:49:48 +08:00
|
|
|
typedef NumTraits<Scalar>::Real RealScalar;
|
|
|
|
typedef Matrix<RealScalar,Dynamic,Dynamic> A;
|
2010-07-22 19:19:09 +08:00
|
|
|
typedef Matrix</*Real*/Scalar,Dynamic,Dynamic> B;
|
2010-07-08 01:49:48 +08:00
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> C;
|
2010-07-22 19:19:09 +08:00
|
|
|
typedef Matrix<RealScalar,Dynamic,Dynamic> M;
|
2010-02-22 18:09:57 +08:00
|
|
|
|
2010-02-26 19:32:00 +08:00
|
|
|
#ifdef HAVE_BLAS
|
|
|
|
|
|
|
|
extern "C" {
|
2012-07-08 22:18:51 +08:00
|
|
|
#include <Eigen/src/misc/blas.h>
|
2010-02-26 19:32:00 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
static float fone = 1;
|
|
|
|
static float fzero = 0;
|
|
|
|
static double done = 1;
|
|
|
|
static double szero = 0;
|
2010-07-05 22:18:09 +08:00
|
|
|
static std::complex<float> cfone = 1;
|
|
|
|
static std::complex<float> cfzero = 0;
|
2010-07-07 22:41:29 +08:00
|
|
|
static std::complex<double> cdone = 1;
|
|
|
|
static std::complex<double> cdzero = 0;
|
2010-02-26 19:32:00 +08:00
|
|
|
static char notrans = 'N';
|
2010-07-22 19:19:09 +08:00
|
|
|
static char trans = 'T';
|
2010-02-26 19:32:00 +08:00
|
|
|
static char nonunit = 'N';
|
|
|
|
static char lower = 'L';
|
|
|
|
static char right = 'R';
|
|
|
|
static int intone = 1;
|
|
|
|
|
|
|
|
void blas_gemm(const MatrixXf& a, const MatrixXf& b, MatrixXf& c)
|
|
|
|
{
|
|
|
|
int M = c.rows(); int N = c.cols(); int K = a.cols();
|
|
|
|
int lda = a.rows(); int ldb = b.rows(); int ldc = c.rows();
|
|
|
|
|
|
|
|
sgemm_(¬rans,¬rans,&M,&N,&K,&fone,
|
|
|
|
const_cast<float*>(a.data()),&lda,
|
|
|
|
const_cast<float*>(b.data()),&ldb,&fone,
|
|
|
|
c.data(),&ldc);
|
|
|
|
}
|
|
|
|
|
2010-12-10 02:46:26 +08:00
|
|
|
EIGEN_DONT_INLINE void blas_gemm(const MatrixXd& a, const MatrixXd& b, MatrixXd& c)
|
|
|
|
{
|
|
|
|
int M = c.rows(); int N = c.cols(); int K = a.cols();
|
|
|
|
int lda = a.rows(); int ldb = b.rows(); int ldc = c.rows();
|
|
|
|
|
|
|
|
dgemm_(¬rans,¬rans,&M,&N,&K,&done,
|
|
|
|
const_cast<double*>(a.data()),&lda,
|
|
|
|
const_cast<double*>(b.data()),&ldb,&done,
|
|
|
|
c.data(),&ldc);
|
|
|
|
}
|
|
|
|
|
2010-07-05 22:18:09 +08:00
|
|
|
void blas_gemm(const MatrixXcf& a, const MatrixXcf& b, MatrixXcf& c)
|
|
|
|
{
|
|
|
|
int M = c.rows(); int N = c.cols(); int K = a.cols();
|
|
|
|
int lda = a.rows(); int ldb = b.rows(); int ldc = c.rows();
|
|
|
|
|
|
|
|
cgemm_(¬rans,¬rans,&M,&N,&K,(float*)&cfone,
|
|
|
|
const_cast<float*>((const float*)a.data()),&lda,
|
|
|
|
const_cast<float*>((const float*)b.data()),&ldb,(float*)&cfone,
|
|
|
|
(float*)c.data(),&ldc);
|
|
|
|
}
|
|
|
|
|
2010-07-07 22:41:29 +08:00
|
|
|
void blas_gemm(const MatrixXcd& a, const MatrixXcd& b, MatrixXcd& c)
|
|
|
|
{
|
|
|
|
int M = c.rows(); int N = c.cols(); int K = a.cols();
|
|
|
|
int lda = a.rows(); int ldb = b.rows(); int ldc = c.rows();
|
|
|
|
|
|
|
|
zgemm_(¬rans,¬rans,&M,&N,&K,(double*)&cdone,
|
|
|
|
const_cast<double*>((const double*)a.data()),&lda,
|
|
|
|
const_cast<double*>((const double*)b.data()),&ldb,(double*)&cdone,
|
|
|
|
(double*)c.data(),&ldc);
|
|
|
|
}
|
|
|
|
|
2010-02-26 19:32:00 +08:00
|
|
|
|
|
|
|
|
|
|
|
#endif
|
|
|
|
|
2010-07-22 19:19:09 +08:00
|
|
|
void matlab_cplx_cplx(const M& ar, const M& ai, const M& br, const M& bi, M& cr, M& ci)
|
|
|
|
{
|
|
|
|
cr.noalias() += ar * br;
|
|
|
|
cr.noalias() -= ai * bi;
|
|
|
|
ci.noalias() += ar * bi;
|
|
|
|
ci.noalias() += ai * br;
|
|
|
|
}
|
|
|
|
|
|
|
|
void matlab_real_cplx(const M& a, const M& br, const M& bi, M& cr, M& ci)
|
|
|
|
{
|
|
|
|
cr.noalias() += a * br;
|
|
|
|
ci.noalias() += a * bi;
|
|
|
|
}
|
|
|
|
|
|
|
|
void matlab_cplx_real(const M& ar, const M& ai, const M& b, M& cr, M& ci)
|
|
|
|
{
|
|
|
|
cr.noalias() += ar * b;
|
|
|
|
ci.noalias() += ai * b;
|
|
|
|
}
|
|
|
|
|
2010-07-08 01:49:48 +08:00
|
|
|
template<typename A, typename B, typename C>
|
|
|
|
EIGEN_DONT_INLINE void gemm(const A& a, const B& b, C& c)
|
2010-02-22 18:09:57 +08:00
|
|
|
{
|
2010-07-22 19:19:09 +08:00
|
|
|
c.noalias() += a * b;
|
2010-02-22 18:09:57 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
int main(int argc, char ** argv)
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
std::ptrdiff_t l1 = internal::queryL1CacheSize();
|
|
|
|
std::ptrdiff_t l2 = internal::queryTopLevelCacheSize();
|
2010-07-07 22:41:29 +08:00
|
|
|
std::cout << "L1 cache size = " << (l1>0 ? l1/1024 : -1) << " KB\n";
|
|
|
|
std::cout << "L2/L3 cache size = " << (l2>0 ? l2/1024 : -1) << " KB\n";
|
2010-10-25 22:15:22 +08:00
|
|
|
typedef internal::gebp_traits<Scalar,Scalar> Traits;
|
2010-07-22 19:19:09 +08:00
|
|
|
std::cout << "Register blocking = " << Traits::mr << " x " << Traits::nr << "\n";
|
2010-06-22 05:28:50 +08:00
|
|
|
|
2010-03-01 17:57:32 +08:00
|
|
|
int rep = 1; // number of repetitions per try
|
2010-06-22 05:28:50 +08:00
|
|
|
int tries = 2; // number of tries, we keep the best
|
2010-02-26 19:32:00 +08:00
|
|
|
|
2010-06-19 05:25:57 +08:00
|
|
|
int s = 2048;
|
|
|
|
int cache_size = -1;
|
|
|
|
|
|
|
|
bool need_help = false;
|
|
|
|
for (int i=1; i<argc; ++i)
|
|
|
|
{
|
|
|
|
if(argv[i][0]=='s')
|
|
|
|
s = atoi(argv[i]+1);
|
|
|
|
else if(argv[i][0]=='c')
|
|
|
|
cache_size = atoi(argv[i]+1);
|
2010-06-22 17:10:38 +08:00
|
|
|
else if(argv[i][0]=='t')
|
|
|
|
tries = atoi(argv[i]+1);
|
|
|
|
else if(argv[i][0]=='p')
|
|
|
|
rep = atoi(argv[i]+1);
|
2010-06-19 05:25:57 +08:00
|
|
|
else
|
|
|
|
need_help = true;
|
|
|
|
}
|
|
|
|
|
|
|
|
if(need_help)
|
|
|
|
{
|
2010-06-22 17:10:38 +08:00
|
|
|
std::cout << argv[0] << " s<matrix size> c<cache size> t<nb tries> p<nb repeats>\n";
|
2010-06-19 05:25:57 +08:00
|
|
|
return 1;
|
|
|
|
}
|
|
|
|
|
|
|
|
if(cache_size>0)
|
2010-07-05 22:18:09 +08:00
|
|
|
setCpuCacheSizes(cache_size,96*cache_size);
|
2010-06-19 05:25:57 +08:00
|
|
|
|
2010-02-22 18:09:57 +08:00
|
|
|
int m = s;
|
|
|
|
int n = s;
|
|
|
|
int p = s;
|
2010-07-22 19:19:09 +08:00
|
|
|
A a(m,p); a.setRandom();
|
|
|
|
B b(p,n); b.setRandom();
|
|
|
|
C c(m,n); c.setOnes();
|
2011-09-09 16:36:20 +08:00
|
|
|
C rc = c;
|
2010-02-22 18:09:57 +08:00
|
|
|
|
2010-06-21 18:07:05 +08:00
|
|
|
std::cout << "Matrix sizes = " << m << "x" << p << " * " << p << "x" << n << "\n";
|
2010-07-22 19:19:09 +08:00
|
|
|
std::ptrdiff_t mc(m), nc(n), kc(p);
|
2011-09-09 16:36:20 +08:00
|
|
|
internal::computeProductBlockingSizes<Scalar,Scalar>(kc, mc, nc);
|
2010-07-22 19:19:09 +08:00
|
|
|
std::cout << "blocking size (mc x kc) = " << mc << " x " << kc << "\n";
|
2010-06-21 18:07:05 +08:00
|
|
|
|
2010-07-08 01:49:48 +08:00
|
|
|
C r = c;
|
2010-02-26 19:32:00 +08:00
|
|
|
|
|
|
|
// check the parallel product is correct
|
2010-07-22 19:19:09 +08:00
|
|
|
#if defined EIGEN_HAS_OPENMP
|
2010-02-26 19:32:00 +08:00
|
|
|
int procs = omp_get_max_threads();
|
2010-03-01 17:57:32 +08:00
|
|
|
if(procs>1)
|
|
|
|
{
|
|
|
|
#ifdef HAVE_BLAS
|
|
|
|
blas_gemm(a,b,r);
|
|
|
|
#else
|
|
|
|
omp_set_num_threads(1);
|
|
|
|
r.noalias() += a * b;
|
|
|
|
omp_set_num_threads(procs);
|
|
|
|
#endif
|
|
|
|
c.noalias() += a * b;
|
|
|
|
if(!r.isApprox(c)) std::cerr << "Warning, your parallel product is crap!\n\n";
|
|
|
|
}
|
2010-07-22 19:19:09 +08:00
|
|
|
#elif defined HAVE_BLAS
|
|
|
|
blas_gemm(a,b,r);
|
|
|
|
c.noalias() += a * b;
|
|
|
|
if(!r.isApprox(c)) std::cerr << "Warning, your product is crap!\n\n";
|
|
|
|
#else
|
|
|
|
gemm(a,b,c);
|
|
|
|
r.noalias() += a.cast<Scalar>() * b.cast<Scalar>();
|
|
|
|
if(!r.isApprox(c)) std::cerr << "Warning, your product is crap!\n\n";
|
2010-02-26 19:32:00 +08:00
|
|
|
#endif
|
|
|
|
|
|
|
|
#ifdef HAVE_BLAS
|
2010-03-01 17:57:32 +08:00
|
|
|
BenchTimer tblas;
|
2011-09-09 16:36:20 +08:00
|
|
|
c = rc;
|
2010-03-01 17:57:32 +08:00
|
|
|
BENCH(tblas, tries, rep, blas_gemm(a,b,c));
|
|
|
|
std::cout << "blas cpu " << tblas.best(CPU_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/tblas.best(CPU_TIMER))*1e-9 << " GFLOPS \t(" << tblas.total(CPU_TIMER) << "s)\n";
|
|
|
|
std::cout << "blas real " << tblas.best(REAL_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/tblas.best(REAL_TIMER))*1e-9 << " GFLOPS \t(" << tblas.total(REAL_TIMER) << "s)\n";
|
2010-02-26 19:32:00 +08:00
|
|
|
#endif
|
2010-02-22 18:09:57 +08:00
|
|
|
|
2010-03-01 17:57:32 +08:00
|
|
|
BenchTimer tmt;
|
2011-09-09 16:36:20 +08:00
|
|
|
c = rc;
|
2010-03-01 17:57:32 +08:00
|
|
|
BENCH(tmt, tries, rep, gemm(a,b,c));
|
|
|
|
std::cout << "eigen cpu " << tmt.best(CPU_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/tmt.best(CPU_TIMER))*1e-9 << " GFLOPS \t(" << tmt.total(CPU_TIMER) << "s)\n";
|
|
|
|
std::cout << "eigen real " << tmt.best(REAL_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/tmt.best(REAL_TIMER))*1e-9 << " GFLOPS \t(" << tmt.total(REAL_TIMER) << "s)\n";
|
|
|
|
|
|
|
|
#ifdef EIGEN_HAS_OPENMP
|
|
|
|
if(procs>1)
|
|
|
|
{
|
|
|
|
BenchTimer tmono;
|
2011-09-09 16:36:20 +08:00
|
|
|
omp_set_num_threads(1);
|
|
|
|
Eigen::internal::setNbThreads(1);
|
|
|
|
c = rc;
|
2010-03-01 17:57:32 +08:00
|
|
|
BENCH(tmono, tries, rep, gemm(a,b,c));
|
|
|
|
std::cout << "eigen mono cpu " << tmono.best(CPU_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/tmono.best(CPU_TIMER))*1e-9 << " GFLOPS \t(" << tmono.total(CPU_TIMER) << "s)\n";
|
|
|
|
std::cout << "eigen mono real " << tmono.best(REAL_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/tmono.best(REAL_TIMER))*1e-9 << " GFLOPS \t(" << tmono.total(REAL_TIMER) << "s)\n";
|
|
|
|
std::cout << "mt speed up x" << tmono.best(CPU_TIMER) / tmt.best(REAL_TIMER) << " => " << (100.0*tmono.best(CPU_TIMER) / tmt.best(REAL_TIMER))/procs << "%\n";
|
|
|
|
}
|
|
|
|
#endif
|
2010-07-22 19:19:09 +08:00
|
|
|
|
|
|
|
#ifdef DECOUPLED
|
|
|
|
if((NumTraits<A::Scalar>::IsComplex) && (NumTraits<B::Scalar>::IsComplex))
|
|
|
|
{
|
|
|
|
M ar(m,p); ar.setRandom();
|
|
|
|
M ai(m,p); ai.setRandom();
|
|
|
|
M br(p,n); br.setRandom();
|
|
|
|
M bi(p,n); bi.setRandom();
|
|
|
|
M cr(m,n); cr.setRandom();
|
|
|
|
M ci(m,n); ci.setRandom();
|
|
|
|
|
|
|
|
BenchTimer t;
|
|
|
|
BENCH(t, tries, rep, matlab_cplx_cplx(ar,ai,br,bi,cr,ci));
|
|
|
|
std::cout << "\"matlab\" cpu " << t.best(CPU_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/t.best(CPU_TIMER))*1e-9 << " GFLOPS \t(" << t.total(CPU_TIMER) << "s)\n";
|
|
|
|
std::cout << "\"matlab\" real " << t.best(REAL_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/t.best(REAL_TIMER))*1e-9 << " GFLOPS \t(" << t.total(REAL_TIMER) << "s)\n";
|
|
|
|
}
|
|
|
|
if((!NumTraits<A::Scalar>::IsComplex) && (NumTraits<B::Scalar>::IsComplex))
|
|
|
|
{
|
|
|
|
M a(m,p); a.setRandom();
|
|
|
|
M br(p,n); br.setRandom();
|
|
|
|
M bi(p,n); bi.setRandom();
|
|
|
|
M cr(m,n); cr.setRandom();
|
|
|
|
M ci(m,n); ci.setRandom();
|
|
|
|
|
|
|
|
BenchTimer t;
|
|
|
|
BENCH(t, tries, rep, matlab_real_cplx(a,br,bi,cr,ci));
|
|
|
|
std::cout << "\"matlab\" cpu " << t.best(CPU_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/t.best(CPU_TIMER))*1e-9 << " GFLOPS \t(" << t.total(CPU_TIMER) << "s)\n";
|
|
|
|
std::cout << "\"matlab\" real " << t.best(REAL_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/t.best(REAL_TIMER))*1e-9 << " GFLOPS \t(" << t.total(REAL_TIMER) << "s)\n";
|
|
|
|
}
|
|
|
|
if((NumTraits<A::Scalar>::IsComplex) && (!NumTraits<B::Scalar>::IsComplex))
|
|
|
|
{
|
|
|
|
M ar(m,p); ar.setRandom();
|
|
|
|
M ai(m,p); ai.setRandom();
|
|
|
|
M b(p,n); b.setRandom();
|
|
|
|
M cr(m,n); cr.setRandom();
|
|
|
|
M ci(m,n); ci.setRandom();
|
|
|
|
|
|
|
|
BenchTimer t;
|
|
|
|
BENCH(t, tries, rep, matlab_cplx_real(ar,ai,b,cr,ci));
|
|
|
|
std::cout << "\"matlab\" cpu " << t.best(CPU_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/t.best(CPU_TIMER))*1e-9 << " GFLOPS \t(" << t.total(CPU_TIMER) << "s)\n";
|
|
|
|
std::cout << "\"matlab\" real " << t.best(REAL_TIMER)/rep << "s \t" << (double(m)*n*p*rep*2/t.best(REAL_TIMER))*1e-9 << " GFLOPS \t(" << t.total(REAL_TIMER) << "s)\n";
|
|
|
|
}
|
|
|
|
#endif
|
2010-02-22 18:09:57 +08:00
|
|
|
|
|
|
|
return 0;
|
|
|
|
}
|
|
|
|
|