aboutsummaryrefslogtreecommitdiffhomepage
path: root/bench/benchEigenSolver.cpp
blob: 395bff3f94a323d1fea8c0466dedcbd1b1617c4d (plain)
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
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210

// g++ -DNDEBUG -O3 -I.. benchEigenSolver.cpp  -o benchEigenSolver && ./benchEigenSolver
// options:
//  -DBENCH_GMM
//  -DBENCH_GSL -lgsl /usr/lib/libcblas.so.3
//  -DEIGEN_DONT_VECTORIZE
//  -msse2
//  -DREPEAT=100
//  -DTRIES=10
//  -DSCALAR=double

#include <Eigen/Array>
#include <Eigen/QR>
#include <bench/BenchUtil.h>
using namespace Eigen;

#ifndef REPEAT
#define REPEAT 1000
#endif

#ifndef TRIES
#define TRIES 4
#endif

#ifndef SCALAR
#define SCALAR float
#endif

typedef SCALAR Scalar;

template <typename MatrixType>
__attribute__ ((noinline)) void benchEigenSolver(const MatrixType& m)
{
  int rows = m.rows();
  int cols = m.cols();

  int stdRepeats = std::max(1,int((REPEAT*1000)/(rows*rows*sqrt(rows))));
  int saRepeats = stdRepeats * 4;

  typedef typename MatrixType::Scalar Scalar;
  typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;

  MatrixType a = MatrixType::Random(rows,cols);
  SquareMatrixType covMat =  a * a.adjoint();

  BenchTimer timerSa, timerStd;

  Scalar acc = 0;
  int r = ei_random<int>(0,covMat.rows()-1);
  int c = ei_random<int>(0,covMat.cols()-1);
  {
    SelfAdjointEigenSolver<SquareMatrixType> ei(covMat);
    for (int t=0; t<TRIES; ++t)
    {
      timerSa.start();
      for (int k=0; k<saRepeats; ++k)
      {
        ei.compute(covMat);
        acc += ei.eigenvectors().coeff(r,c);
      }
      timerSa.stop();
    }
  }

  {
    EigenSolver<SquareMatrixType> ei(covMat);
    for (int t=0; t<TRIES; ++t)
    {
      timerStd.start();
      for (int k=0; k<stdRepeats; ++k)
      {
        ei.compute(covMat);
        acc += ei.eigenvectors().coeff(r,c);
      }
      timerStd.stop();
    }
  }

  if (MatrixType::RowsAtCompileTime==Dynamic)
    std::cout << "dyn   ";
  else
    std::cout << "fixed ";
  std::cout << covMat.rows() << " \t"
            << timerSa.value() * REPEAT / saRepeats << "s \t"
            << timerStd.value() * REPEAT / stdRepeats << "s";

  #ifdef BENCH_GMM
  if (MatrixType::RowsAtCompileTime==Dynamic)
  {
    timerSa.reset();
    timerStd.reset();

    gmm::dense_matrix<Scalar> gmmCovMat(covMat.rows(),covMat.cols());
    gmm::dense_matrix<Scalar> eigvect(covMat.rows(),covMat.cols());
    std::vector<Scalar> eigval(covMat.rows());
    eiToGmm(covMat, gmmCovMat);
    for (int t=0; t<TRIES; ++t)
    {
      timerSa.start();
      for (int k=0; k<saRepeats; ++k)
      {
        gmm::symmetric_qr_algorithm(gmmCovMat, eigval, eigvect);
        acc += eigvect(r,c);
      }
      timerSa.stop();
    }
    // the non-selfadjoint solver does not compute the eigen vectors
//     for (int t=0; t<TRIES; ++t)
//     {
//       timerStd.start();
//       for (int k=0; k<stdRepeats; ++k)
//       {
//         gmm::implicit_qr_algorithm(gmmCovMat, eigval, eigvect);
//         acc += eigvect(r,c);
//       }
//       timerStd.stop();
//     }

    std::cout << " | \t"
              << timerSa.value() * REPEAT / saRepeats << "s"
              << /*timerStd.value() * REPEAT / stdRepeats << "s"*/ "   na   ";
  }
  #endif

  #ifdef BENCH_GSL
  if (MatrixType::RowsAtCompileTime==Dynamic)
  {
    timerSa.reset();
    timerStd.reset();

    gsl_matrix* gslCovMat = gsl_matrix_alloc(covMat.rows(),covMat.cols());
    gsl_matrix* gslCopy = gsl_matrix_alloc(covMat.rows(),covMat.cols());
    gsl_matrix* eigvect = gsl_matrix_alloc(covMat.rows(),covMat.cols());
    gsl_vector* eigval  = gsl_vector_alloc(covMat.rows());
    gsl_eigen_symmv_workspace* eisymm = gsl_eigen_symmv_alloc(covMat.rows());
    
    gsl_matrix_complex* eigvectz = gsl_matrix_complex_alloc(covMat.rows(),covMat.cols());
    gsl_vector_complex* eigvalz  = gsl_vector_complex_alloc(covMat.rows());
    gsl_eigen_nonsymmv_workspace* einonsymm = gsl_eigen_nonsymmv_alloc(covMat.rows());
    
    eiToGsl(covMat, &gslCovMat);
    for (int t=0; t<TRIES; ++t)
    {
      timerSa.start();
      for (int k=0; k<saRepeats; ++k)
      {
        gsl_matrix_memcpy(gslCopy,gslCovMat);
        gsl_eigen_symmv(gslCopy, eigval, eigvect, eisymm);
        acc += gsl_matrix_get(eigvect,r,c);
      }
      timerSa.stop();
    }
    for (int t=0; t<TRIES; ++t)
    {
      timerStd.start();
      for (int k=0; k<stdRepeats; ++k)
      {
        gsl_matrix_memcpy(gslCopy,gslCovMat);
        gsl_eigen_nonsymmv(gslCopy, eigvalz, eigvectz, einonsymm);
        acc += GSL_REAL(gsl_matrix_complex_get(eigvectz,r,c));
      }
      timerStd.stop();
    }

    std::cout << " | \t"
              << timerSa.value() * REPEAT / saRepeats << "s \t"
              << timerStd.value() * REPEAT / stdRepeats << "s";

    gsl_matrix_free(gslCovMat);
    gsl_vector_free(gslCopy);
    gsl_matrix_free(eigvect);
    gsl_vector_free(eigval);
    gsl_matrix_complex_free(eigvectz);
    gsl_vector_complex_free(eigvalz);
    gsl_eigen_symmv_free(eisymm);
    gsl_eigen_nonsymmv_free(einonsymm);
  }
  #endif

  std::cout << "\n";
  
  // make sure the compiler does not optimize too much
  if (acc==123)
    std::cout << acc;
}

int main(int argc, char* argv[])
{
  const int dynsizes[] = {4,6,8,12,16,24,32,64,128,256,512,0};
  std::cout << "size            selfadjoint       generic";
  #ifdef BENCH_GMM
  std::cout << "        GMM++          ";
  #endif
  #ifdef BENCH_GSL
  std::cout << "       GSL (double + ATLAS)  ";
  #endif
  std::cout << "\n";
  for (uint i=0; dynsizes[i]>0; ++i)
    benchEigenSolver(Matrix<Scalar,Dynamic,Dynamic>(dynsizes[i],dynsizes[i]));

  benchEigenSolver(Matrix<Scalar,2,2>());
  benchEigenSolver(Matrix<Scalar,3,3>());
  benchEigenSolver(Matrix<Scalar,4,4>());
  benchEigenSolver(Matrix<Scalar,6,6>());
  benchEigenSolver(Matrix<Scalar,8,8>());
  benchEigenSolver(Matrix<Scalar,12,12>());
  benchEigenSolver(Matrix<Scalar,16,16>());
  return 0;
}