benchCholesky.cpp File Reference
#include <iostream>
#include <Eigen/Core>
#include <Eigen/Cholesky>
#include <bench/BenchUtil.h>

Macros

#define REPEAT   10000
 
#define TRIES   10
 

Typedefs

typedef float Scalar
 

Functions

template<typename MatrixType >
 __attribute__ ((noinline)) void benchLLT(const MatrixType &m)
 
int main (int argc, char *argv[])
 

Macro Definition Documentation

◆ REPEAT

#define REPEAT   10000

◆ TRIES

#define TRIES   10

Typedef Documentation

◆ Scalar

typedef float Scalar

Function Documentation

◆ __attribute__()

template<typename MatrixType >
__attribute__ ( (noinline)  ) const &
28  {
29  int rows = m.rows();
30  int cols = m.cols();
31 
32  double cost = 0;
33  for (int j = 0; j < rows; ++j) {
34  int r = std::max(rows - j - 1, 0);
35  cost += 2 * (r * j + r + j);
36  }
37 
38  int repeats = (REPEAT * 1000) / (rows * rows);
39 
40  typedef typename MatrixType::Scalar Scalar;
42 
43  MatrixType a = MatrixType::Random(rows, cols);
44  SquareMatrixType covMat = a * a.adjoint();
45 
46  BenchTimer timerNoSqrt, timerSqrt;
47 
48  Scalar acc = 0;
49  int r = internal::random<int>(0, covMat.rows() - 1);
50  int c = internal::random<int>(0, covMat.cols() - 1);
51  for (int t = 0; t < TRIES; ++t) {
52  timerNoSqrt.start();
53  for (int k = 0; k < repeats; ++k) {
54  LDLT<SquareMatrixType> cholnosqrt(covMat);
55  acc += cholnosqrt.matrixL().coeff(r, c);
56  }
57  timerNoSqrt.stop();
58  }
59 
60  for (int t = 0; t < TRIES; ++t) {
61  timerSqrt.start();
62  for (int k = 0; k < repeats; ++k) {
63  LLT<SquareMatrixType> chol(covMat);
64  acc += chol.matrixL().coeff(r, c);
65  }
66  timerSqrt.stop();
67  }
68 
69  if (MatrixType::RowsAtCompileTime == Dynamic)
70  std::cout << "dyn ";
71  else
72  std::cout << "fixed ";
73  std::cout << covMat.rows() << " \t" << (timerNoSqrt.best()) / repeats << "s "
74  << "(" << 1e-9 * cost * repeats / timerNoSqrt.best() << " GFLOPS)\t" << (timerSqrt.best()) / repeats << "s "
75  << "(" << 1e-9 * cost * repeats / timerSqrt.best() << " GFLOPS)\n";
76 
77 #ifdef BENCH_GSL
78  if (MatrixType::RowsAtCompileTime == Dynamic) {
79  timerSqrt.reset();
80 
81  gsl_matrix* gslCovMat = gsl_matrix_alloc(covMat.rows(), covMat.cols());
82  gsl_matrix* gslCopy = gsl_matrix_alloc(covMat.rows(), covMat.cols());
83 
84  eiToGsl(covMat, &gslCovMat);
85  for (int t = 0; t < TRIES; ++t) {
86  timerSqrt.start();
87  for (int k = 0; k < repeats; ++k) {
88  gsl_matrix_memcpy(gslCopy, gslCovMat);
89  gsl_linalg_cholesky_decomp(gslCopy);
90  acc += gsl_matrix_get(gslCopy, r, c);
91  }
92  timerSqrt.stop();
93  }
94 
95  std::cout << " | \t" << timerSqrt.value() * REPEAT / repeats << "s";
96 
97  gsl_matrix_free(gslCovMat);
98  }
99 #endif
100  std::cout << "\n";
101  // make sure the compiler does not optimize too much
102  if (acc == 123) std::cout << acc;
103 }
Array< double, 1, 3 > e(1./3., 0.5, 2.)
int rows
Definition: Tutorial_commainit_02.cpp:1
int cols
Definition: Tutorial_commainit_02.cpp:1
#define REPEAT
Definition: benchCholesky.cpp:18
#define TRIES
Definition: benchCholesky.cpp:22
float Scalar
Definition: benchCholesky.cpp:25
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
Definition: BenchTimer.h:55
double value(int TIMER=CPU_TIMER) const
Definition: BenchTimer.h:94
void reset()
Definition: BenchTimer.h:68
void stop()
Definition: BenchTimer.h:77
void start()
Definition: BenchTimer.h:73
double best(int TIMER=CPU_TIMER) const
Definition: BenchTimer.h:98
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:63
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:70
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
#define max(a, b)
Definition: datatypes.h:23
const Scalar * a
Definition: level2_cplx_impl.h:32
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
const int Dynamic
Definition: Constants.h:25
r
Definition: UniformPSDSelfTest.py:20
int c
Definition: calibrate.py:100
t
Definition: plotPSD.py:36
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References a, Eigen::BenchTimer::best(), calibrate::c, cols, Eigen::Dynamic, e(), j, k, m, Eigen::LDLT< MatrixType_, UpLo_ >::matrixL(), Eigen::LLT< MatrixType_, UpLo_ >::matrixL(), max, UniformPSDSelfTest::r, REPEAT, Eigen::BenchTimer::reset(), rows, Eigen::BenchTimer::start(), Eigen::BenchTimer::stop(), plotPSD::t, TRIES, and Eigen::BenchTimer::value().

◆ main()

int main ( int argc  ,
char argv[] 
)
105  {
106  const int dynsizes[] = {4, 6, 8, 16, 24, 32, 49, 64, 128, 256, 512, 900, 1500, 0};
107  std::cout << "size LDLT LLT";
108  // #ifdef BENCH_GSL
109  // std::cout << " GSL (standard + double + ATLAS) ";
110  // #endif
111  std::cout << "\n";
112  for (int i = 0; dynsizes[i] > 0; ++i) benchLLT(Matrix<Scalar, Dynamic, Dynamic>(dynsizes[i], dynsizes[i]));
113 
114  benchLLT(Matrix<Scalar, 2, 2>());
115  benchLLT(Matrix<Scalar, 3, 3>());
116  benchLLT(Matrix<Scalar, 4, 4>());
117  benchLLT(Matrix<Scalar, 5, 5>());
118  benchLLT(Matrix<Scalar, 6, 6>());
119  benchLLT(Matrix<Scalar, 7, 7>());
120  benchLLT(Matrix<Scalar, 8, 8>());
121  benchLLT(Matrix<Scalar, 12, 12>());
122  benchLLT(Matrix<Scalar, 16, 16>());
123  return 0;
124 }
int i
Definition: BiCGSTAB_step_by_step.cpp:9

References i.