aboutsummaryrefslogtreecommitdiffhomepage
path: root/bench/BenchSparseUtil.h
blob: 13981f6b716b816c252c558957ae59a0b5e05a5f (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

#include <Eigen/Sparse>
#include <bench/BenchTimer.h>
#include <set>

using namespace std;
using namespace Eigen;
using namespace Eigen;

#ifndef SIZE
#define SIZE 1024
#endif

#ifndef DENSITY
#define DENSITY 0.01
#endif

#ifndef SCALAR
#define SCALAR double
#endif

typedef SCALAR Scalar;
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
typedef Matrix<Scalar,Dynamic,1> DenseVector;
typedef SparseMatrix<Scalar> EigenSparseMatrix;

void fillMatrix(float density, int rows, int cols,  EigenSparseMatrix& dst)
{
  dst.reserve(double(rows)*cols*density);
  for(int j = 0; j < cols; j++)
  {
    for(int i = 0; i < rows; i++)
    {
      Scalar v = (internal::random<float>(0,1) < density) ? internal::random<Scalar>() : 0;
      if (v!=0)
        dst.insert(i,j) = v;
    }
  }
  dst.finalize();
}

void fillMatrix2(int nnzPerCol, int rows, int cols,  EigenSparseMatrix& dst)
{
//   std::cout << "alloc " << nnzPerCol*cols << "\n";
  dst.reserve(nnzPerCol*cols);
  for(int j = 0; j < cols; j++)
  {
    std::set<int> aux;
    for(int i = 0; i < nnzPerCol; i++)
    {
      int k = internal::random<int>(0,rows-1);
      while (aux.find(k)!=aux.end())
        k = internal::random<int>(0,rows-1);
      aux.insert(k);

      dst.insert(k,j) = internal::random<Scalar>();
    }
  }
  dst.finalize();
}

void eiToDense(const EigenSparseMatrix& src, DenseMatrix& dst)
{
  dst.setZero();
  for (int j=0; j<src.cols(); ++j)
    for (EigenSparseMatrix::InnerIterator it(src.derived(), j); it; ++it)
      dst(it.index(),j) = it.value();
}

#ifndef NOGMM
#include "gmm/gmm.h"
typedef gmm::csc_matrix<Scalar> GmmSparse;
typedef gmm::col_matrix< gmm::wsvector<Scalar> > GmmDynSparse;
void eiToGmm(const EigenSparseMatrix& src, GmmSparse& dst)
{
  GmmDynSparse tmp(src.rows(), src.cols());
  for (int j=0; j<src.cols(); ++j)
    for (EigenSparseMatrix::InnerIterator it(src.derived(), j); it; ++it)
      tmp(it.index(),j) = it.value();
  gmm::copy(tmp, dst);
}
#endif

#ifndef NOMTL
#include <boost/numeric/mtl/mtl.hpp>
typedef mtl::compressed2D<Scalar, mtl::matrix::parameters<mtl::tag::col_major> > MtlSparse;
typedef mtl::compressed2D<Scalar, mtl::matrix::parameters<mtl::tag::row_major> > MtlSparseRowMajor;
void eiToMtl(const EigenSparseMatrix& src, MtlSparse& dst)
{
  mtl::matrix::inserter<MtlSparse> ins(dst);
  for (int j=0; j<src.cols(); ++j)
    for (EigenSparseMatrix::InnerIterator it(src.derived(), j); it; ++it)
      ins[it.index()][j] = it.value();
}
#endif

#ifdef CSPARSE
extern "C" {
#include "cs.h"
}
void eiToCSparse(const EigenSparseMatrix& src, cs* &dst)
{
  cs* aux = cs_spalloc (0, 0, 1, 1, 1);
  for (int j=0; j<src.cols(); ++j)
    for (EigenSparseMatrix::InnerIterator it(src.derived(), j); it; ++it)
      if (!cs_entry(aux, it.index(), j, it.value()))
      {
        std::cout << "cs_entry error\n";
        exit(2);
      }
   dst = cs_compress(aux);
//    cs_spfree(aux);
}
#endif // CSPARSE

#ifndef NOUBLAS
#include <boost/numeric/ublas/vector.hpp>
#include <boost/numeric/ublas/matrix.hpp>
#include <boost/numeric/ublas/io.hpp>
#include <boost/numeric/ublas/triangular.hpp>
#include <boost/numeric/ublas/vector_sparse.hpp>
#include <boost/numeric/ublas/matrix_sparse.hpp>
#include <boost/numeric/ublas/vector_of_vector.hpp>
#include <boost/numeric/ublas/operation.hpp>

typedef boost::numeric::ublas::compressed_matrix<Scalar,boost::numeric::ublas::column_major> UBlasSparse;

void eiToUblas(const EigenSparseMatrix& src, UBlasSparse& dst)
{
  dst.resize(src.rows(), src.cols(), false);
  for (int j=0; j<src.cols(); ++j)
    for (EigenSparseMatrix::InnerIterator it(src.derived(), j); it; ++it)
      dst(it.index(),j) = it.value();
}

template <typename EigenType, typename UblasType>
void eiToUblasVec(const EigenType& src, UblasType& dst)
{
  dst.resize(src.size());
  for (int j=0; j<src.size(); ++j)
      dst[j] = src.coeff(j);
}
#endif

#ifdef OSKI
extern "C" {
#include <oski/oski.h>
}
#endif