1*bf2c3715SXin Li
2*bf2c3715SXin Li // g++ -DNDEBUG -O3 -I.. benchEigenSolver.cpp -o benchEigenSolver && ./benchEigenSolver
3*bf2c3715SXin Li // options:
4*bf2c3715SXin Li // -DBENCH_GMM
5*bf2c3715SXin Li // -DBENCH_GSL -lgsl /usr/lib/libcblas.so.3
6*bf2c3715SXin Li // -DEIGEN_DONT_VECTORIZE
7*bf2c3715SXin Li // -msse2
8*bf2c3715SXin Li // -DREPEAT=100
9*bf2c3715SXin Li // -DTRIES=10
10*bf2c3715SXin Li // -DSCALAR=double
11*bf2c3715SXin Li
12*bf2c3715SXin Li #include <iostream>
13*bf2c3715SXin Li
14*bf2c3715SXin Li #include <Eigen/Core>
15*bf2c3715SXin Li #include <Eigen/QR>
16*bf2c3715SXin Li #include <bench/BenchUtil.h>
17*bf2c3715SXin Li using namespace Eigen;
18*bf2c3715SXin Li
19*bf2c3715SXin Li #ifndef REPEAT
20*bf2c3715SXin Li #define REPEAT 1000
21*bf2c3715SXin Li #endif
22*bf2c3715SXin Li
23*bf2c3715SXin Li #ifndef TRIES
24*bf2c3715SXin Li #define TRIES 4
25*bf2c3715SXin Li #endif
26*bf2c3715SXin Li
27*bf2c3715SXin Li #ifndef SCALAR
28*bf2c3715SXin Li #define SCALAR float
29*bf2c3715SXin Li #endif
30*bf2c3715SXin Li
31*bf2c3715SXin Li typedef SCALAR Scalar;
32*bf2c3715SXin Li
33*bf2c3715SXin Li template <typename MatrixType>
benchEigenSolver(const MatrixType & m)34*bf2c3715SXin Li __attribute__ ((noinline)) void benchEigenSolver(const MatrixType& m)
35*bf2c3715SXin Li {
36*bf2c3715SXin Li int rows = m.rows();
37*bf2c3715SXin Li int cols = m.cols();
38*bf2c3715SXin Li
39*bf2c3715SXin Li int stdRepeats = std::max(1,int((REPEAT*1000)/(rows*rows*sqrt(rows))));
40*bf2c3715SXin Li int saRepeats = stdRepeats * 4;
41*bf2c3715SXin Li
42*bf2c3715SXin Li typedef typename MatrixType::Scalar Scalar;
43*bf2c3715SXin Li typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
44*bf2c3715SXin Li
45*bf2c3715SXin Li MatrixType a = MatrixType::Random(rows,cols);
46*bf2c3715SXin Li SquareMatrixType covMat = a * a.adjoint();
47*bf2c3715SXin Li
48*bf2c3715SXin Li BenchTimer timerSa, timerStd;
49*bf2c3715SXin Li
50*bf2c3715SXin Li Scalar acc = 0;
51*bf2c3715SXin Li int r = internal::random<int>(0,covMat.rows()-1);
52*bf2c3715SXin Li int c = internal::random<int>(0,covMat.cols()-1);
53*bf2c3715SXin Li {
54*bf2c3715SXin Li SelfAdjointEigenSolver<SquareMatrixType> ei(covMat);
55*bf2c3715SXin Li for (int t=0; t<TRIES; ++t)
56*bf2c3715SXin Li {
57*bf2c3715SXin Li timerSa.start();
58*bf2c3715SXin Li for (int k=0; k<saRepeats; ++k)
59*bf2c3715SXin Li {
60*bf2c3715SXin Li ei.compute(covMat);
61*bf2c3715SXin Li acc += ei.eigenvectors().coeff(r,c);
62*bf2c3715SXin Li }
63*bf2c3715SXin Li timerSa.stop();
64*bf2c3715SXin Li }
65*bf2c3715SXin Li }
66*bf2c3715SXin Li
67*bf2c3715SXin Li {
68*bf2c3715SXin Li EigenSolver<SquareMatrixType> ei(covMat);
69*bf2c3715SXin Li for (int t=0; t<TRIES; ++t)
70*bf2c3715SXin Li {
71*bf2c3715SXin Li timerStd.start();
72*bf2c3715SXin Li for (int k=0; k<stdRepeats; ++k)
73*bf2c3715SXin Li {
74*bf2c3715SXin Li ei.compute(covMat);
75*bf2c3715SXin Li acc += ei.eigenvectors().coeff(r,c);
76*bf2c3715SXin Li }
77*bf2c3715SXin Li timerStd.stop();
78*bf2c3715SXin Li }
79*bf2c3715SXin Li }
80*bf2c3715SXin Li
81*bf2c3715SXin Li if (MatrixType::RowsAtCompileTime==Dynamic)
82*bf2c3715SXin Li std::cout << "dyn ";
83*bf2c3715SXin Li else
84*bf2c3715SXin Li std::cout << "fixed ";
85*bf2c3715SXin Li std::cout << covMat.rows() << " \t"
86*bf2c3715SXin Li << timerSa.value() * REPEAT / saRepeats << "s \t"
87*bf2c3715SXin Li << timerStd.value() * REPEAT / stdRepeats << "s";
88*bf2c3715SXin Li
89*bf2c3715SXin Li #ifdef BENCH_GMM
90*bf2c3715SXin Li if (MatrixType::RowsAtCompileTime==Dynamic)
91*bf2c3715SXin Li {
92*bf2c3715SXin Li timerSa.reset();
93*bf2c3715SXin Li timerStd.reset();
94*bf2c3715SXin Li
95*bf2c3715SXin Li gmm::dense_matrix<Scalar> gmmCovMat(covMat.rows(),covMat.cols());
96*bf2c3715SXin Li gmm::dense_matrix<Scalar> eigvect(covMat.rows(),covMat.cols());
97*bf2c3715SXin Li std::vector<Scalar> eigval(covMat.rows());
98*bf2c3715SXin Li eiToGmm(covMat, gmmCovMat);
99*bf2c3715SXin Li for (int t=0; t<TRIES; ++t)
100*bf2c3715SXin Li {
101*bf2c3715SXin Li timerSa.start();
102*bf2c3715SXin Li for (int k=0; k<saRepeats; ++k)
103*bf2c3715SXin Li {
104*bf2c3715SXin Li gmm::symmetric_qr_algorithm(gmmCovMat, eigval, eigvect);
105*bf2c3715SXin Li acc += eigvect(r,c);
106*bf2c3715SXin Li }
107*bf2c3715SXin Li timerSa.stop();
108*bf2c3715SXin Li }
109*bf2c3715SXin Li // the non-selfadjoint solver does not compute the eigen vectors
110*bf2c3715SXin Li // for (int t=0; t<TRIES; ++t)
111*bf2c3715SXin Li // {
112*bf2c3715SXin Li // timerStd.start();
113*bf2c3715SXin Li // for (int k=0; k<stdRepeats; ++k)
114*bf2c3715SXin Li // {
115*bf2c3715SXin Li // gmm::implicit_qr_algorithm(gmmCovMat, eigval, eigvect);
116*bf2c3715SXin Li // acc += eigvect(r,c);
117*bf2c3715SXin Li // }
118*bf2c3715SXin Li // timerStd.stop();
119*bf2c3715SXin Li // }
120*bf2c3715SXin Li
121*bf2c3715SXin Li std::cout << " | \t"
122*bf2c3715SXin Li << timerSa.value() * REPEAT / saRepeats << "s"
123*bf2c3715SXin Li << /*timerStd.value() * REPEAT / stdRepeats << "s"*/ " na ";
124*bf2c3715SXin Li }
125*bf2c3715SXin Li #endif
126*bf2c3715SXin Li
127*bf2c3715SXin Li #ifdef BENCH_GSL
128*bf2c3715SXin Li if (MatrixType::RowsAtCompileTime==Dynamic)
129*bf2c3715SXin Li {
130*bf2c3715SXin Li timerSa.reset();
131*bf2c3715SXin Li timerStd.reset();
132*bf2c3715SXin Li
133*bf2c3715SXin Li gsl_matrix* gslCovMat = gsl_matrix_alloc(covMat.rows(),covMat.cols());
134*bf2c3715SXin Li gsl_matrix* gslCopy = gsl_matrix_alloc(covMat.rows(),covMat.cols());
135*bf2c3715SXin Li gsl_matrix* eigvect = gsl_matrix_alloc(covMat.rows(),covMat.cols());
136*bf2c3715SXin Li gsl_vector* eigval = gsl_vector_alloc(covMat.rows());
137*bf2c3715SXin Li gsl_eigen_symmv_workspace* eisymm = gsl_eigen_symmv_alloc(covMat.rows());
138*bf2c3715SXin Li
139*bf2c3715SXin Li gsl_matrix_complex* eigvectz = gsl_matrix_complex_alloc(covMat.rows(),covMat.cols());
140*bf2c3715SXin Li gsl_vector_complex* eigvalz = gsl_vector_complex_alloc(covMat.rows());
141*bf2c3715SXin Li gsl_eigen_nonsymmv_workspace* einonsymm = gsl_eigen_nonsymmv_alloc(covMat.rows());
142*bf2c3715SXin Li
143*bf2c3715SXin Li eiToGsl(covMat, &gslCovMat);
144*bf2c3715SXin Li for (int t=0; t<TRIES; ++t)
145*bf2c3715SXin Li {
146*bf2c3715SXin Li timerSa.start();
147*bf2c3715SXin Li for (int k=0; k<saRepeats; ++k)
148*bf2c3715SXin Li {
149*bf2c3715SXin Li gsl_matrix_memcpy(gslCopy,gslCovMat);
150*bf2c3715SXin Li gsl_eigen_symmv(gslCopy, eigval, eigvect, eisymm);
151*bf2c3715SXin Li acc += gsl_matrix_get(eigvect,r,c);
152*bf2c3715SXin Li }
153*bf2c3715SXin Li timerSa.stop();
154*bf2c3715SXin Li }
155*bf2c3715SXin Li for (int t=0; t<TRIES; ++t)
156*bf2c3715SXin Li {
157*bf2c3715SXin Li timerStd.start();
158*bf2c3715SXin Li for (int k=0; k<stdRepeats; ++k)
159*bf2c3715SXin Li {
160*bf2c3715SXin Li gsl_matrix_memcpy(gslCopy,gslCovMat);
161*bf2c3715SXin Li gsl_eigen_nonsymmv(gslCopy, eigvalz, eigvectz, einonsymm);
162*bf2c3715SXin Li acc += GSL_REAL(gsl_matrix_complex_get(eigvectz,r,c));
163*bf2c3715SXin Li }
164*bf2c3715SXin Li timerStd.stop();
165*bf2c3715SXin Li }
166*bf2c3715SXin Li
167*bf2c3715SXin Li std::cout << " | \t"
168*bf2c3715SXin Li << timerSa.value() * REPEAT / saRepeats << "s \t"
169*bf2c3715SXin Li << timerStd.value() * REPEAT / stdRepeats << "s";
170*bf2c3715SXin Li
171*bf2c3715SXin Li gsl_matrix_free(gslCovMat);
172*bf2c3715SXin Li gsl_vector_free(gslCopy);
173*bf2c3715SXin Li gsl_matrix_free(eigvect);
174*bf2c3715SXin Li gsl_vector_free(eigval);
175*bf2c3715SXin Li gsl_matrix_complex_free(eigvectz);
176*bf2c3715SXin Li gsl_vector_complex_free(eigvalz);
177*bf2c3715SXin Li gsl_eigen_symmv_free(eisymm);
178*bf2c3715SXin Li gsl_eigen_nonsymmv_free(einonsymm);
179*bf2c3715SXin Li }
180*bf2c3715SXin Li #endif
181*bf2c3715SXin Li
182*bf2c3715SXin Li std::cout << "\n";
183*bf2c3715SXin Li
184*bf2c3715SXin Li // make sure the compiler does not optimize too much
185*bf2c3715SXin Li if (acc==123)
186*bf2c3715SXin Li std::cout << acc;
187*bf2c3715SXin Li }
188*bf2c3715SXin Li
main(int argc,char * argv[])189*bf2c3715SXin Li int main(int argc, char* argv[])
190*bf2c3715SXin Li {
191*bf2c3715SXin Li const int dynsizes[] = {4,6,8,12,16,24,32,64,128,256,512,0};
192*bf2c3715SXin Li std::cout << "size selfadjoint generic";
193*bf2c3715SXin Li #ifdef BENCH_GMM
194*bf2c3715SXin Li std::cout << " GMM++ ";
195*bf2c3715SXin Li #endif
196*bf2c3715SXin Li #ifdef BENCH_GSL
197*bf2c3715SXin Li std::cout << " GSL (double + ATLAS) ";
198*bf2c3715SXin Li #endif
199*bf2c3715SXin Li std::cout << "\n";
200*bf2c3715SXin Li for (uint i=0; dynsizes[i]>0; ++i)
201*bf2c3715SXin Li benchEigenSolver(Matrix<Scalar,Dynamic,Dynamic>(dynsizes[i],dynsizes[i]));
202*bf2c3715SXin Li
203*bf2c3715SXin Li benchEigenSolver(Matrix<Scalar,2,2>());
204*bf2c3715SXin Li benchEigenSolver(Matrix<Scalar,3,3>());
205*bf2c3715SXin Li benchEigenSolver(Matrix<Scalar,4,4>());
206*bf2c3715SXin Li benchEigenSolver(Matrix<Scalar,6,6>());
207*bf2c3715SXin Li benchEigenSolver(Matrix<Scalar,8,8>());
208*bf2c3715SXin Li benchEigenSolver(Matrix<Scalar,12,12>());
209*bf2c3715SXin Li benchEigenSolver(Matrix<Scalar,16,16>());
210*bf2c3715SXin Li return 0;
211*bf2c3715SXin Li }
212*bf2c3715SXin Li
213