• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 
2 #ifndef EIGEN_BENCH_UTIL_H
3 #define EIGEN_BENCH_UTIL_H
4 
5 #include <Eigen/Core>
6 #include "BenchTimer.h"
7 
8 using namespace std;
9 using namespace Eigen;
10 
11 #include <boost/preprocessor/repetition/enum_params.hpp>
12 #include <boost/preprocessor/repetition.hpp>
13 #include <boost/preprocessor/seq.hpp>
14 #include <boost/preprocessor/array.hpp>
15 #include <boost/preprocessor/arithmetic.hpp>
16 #include <boost/preprocessor/comparison.hpp>
17 #include <boost/preprocessor/punctuation.hpp>
18 #include <boost/preprocessor/punctuation/comma.hpp>
19 #include <boost/preprocessor/stringize.hpp>
20 
21 template<typename MatrixType> void initMatrix_random(MatrixType& mat) __attribute__((noinline));
initMatrix_random(MatrixType & mat)22 template<typename MatrixType> void initMatrix_random(MatrixType& mat)
23 {
24   mat.setRandom();// = MatrixType::random(mat.rows(), mat.cols());
25 }
26 
27 template<typename MatrixType> void initMatrix_identity(MatrixType& mat) __attribute__((noinline));
initMatrix_identity(MatrixType & mat)28 template<typename MatrixType> void initMatrix_identity(MatrixType& mat)
29 {
30   mat.setIdentity();
31 }
32 
33 #ifndef __INTEL_COMPILER
34 #define DISABLE_SSE_EXCEPTIONS()  { \
35   int aux; \
36   asm( \
37   "stmxcsr   %[aux]           \n\t" \
38   "orl       $32832, %[aux]   \n\t" \
39   "ldmxcsr   %[aux]           \n\t" \
40   : : [aux] "m" (aux)); \
41 }
42 #else
43 #define DISABLE_SSE_EXCEPTIONS()
44 #endif
45 
46 #ifdef BENCH_GMM
47 #include <gmm/gmm.h>
48 template <typename EigenMatrixType, typename GmmMatrixType>
eiToGmm(const EigenMatrixType & src,GmmMatrixType & dst)49 void eiToGmm(const EigenMatrixType& src, GmmMatrixType& dst)
50 {
51   dst.resize(src.rows(),src.cols());
52   for (int j=0; j<src.cols(); ++j)
53     for (int i=0; i<src.rows(); ++i)
54       dst(i,j) = src.coeff(i,j);
55 }
56 #endif
57 
58 
59 #ifdef BENCH_GSL
60 #include <gsl/gsl_matrix.h>
61 #include <gsl/gsl_linalg.h>
62 #include <gsl/gsl_eigen.h>
63 template <typename EigenMatrixType>
eiToGsl(const EigenMatrixType & src,gsl_matrix ** dst)64 void eiToGsl(const EigenMatrixType& src, gsl_matrix** dst)
65 {
66   for (int j=0; j<src.cols(); ++j)
67     for (int i=0; i<src.rows(); ++i)
68       gsl_matrix_set(*dst, i, j, src.coeff(i,j));
69 }
70 #endif
71 
72 #ifdef BENCH_UBLAS
73 #include <boost/numeric/ublas/matrix.hpp>
74 #include <boost/numeric/ublas/vector.hpp>
75 template <typename EigenMatrixType, typename UblasMatrixType>
eiToUblas(const EigenMatrixType & src,UblasMatrixType & dst)76 void eiToUblas(const EigenMatrixType& src, UblasMatrixType& dst)
77 {
78   dst.resize(src.rows(),src.cols());
79   for (int j=0; j<src.cols(); ++j)
80     for (int i=0; i<src.rows(); ++i)
81       dst(i,j) = src.coeff(i,j);
82 }
83 template <typename EigenType, typename UblasType>
eiToUblasVec(const EigenType & src,UblasType & dst)84 void eiToUblasVec(const EigenType& src, UblasType& dst)
85 {
86   dst.resize(src.size());
87   for (int j=0; j<src.size(); ++j)
88       dst[j] = src.coeff(j);
89 }
90 #endif
91 
92 #endif // EIGEN_BENCH_UTIL_H
93