BenchSparseUtil.h File Reference
#include <Eigen/Sparse>
#include <bench/BenchTimer.h>
#include <set>
#include "gmm/gmm.h"
#include <boost/numeric/mtl/mtl.hpp>
#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>

Go to the source code of this file.

Macros

#define SIZE   1024
 
#define DENSITY   0.01
 
#define SCALAR   double
 

Typedefs

typedef SCALAR Scalar
 
typedef Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
 
typedef Matrix< Scalar, Dynamic, 1 > DenseVector
 
typedef SparseMatrix< ScalarEigenSparseMatrix
 
typedef gmm::csc_matrix< ScalarGmmSparse
 
typedef gmm::col_matrix< gmm::wsvector< Scalar > > GmmDynSparse
 
typedef mtl::compressed2D< Scalar, mtl::matrix::parameters< mtl::tag::col_major > > MtlSparse
 
typedef mtl::compressed2D< Scalar, mtl::matrix::parameters< mtl::tag::row_major > > MtlSparseRowMajor
 
typedef boost::numeric::ublas::compressed_matrix< Scalar, boost::numeric::ublas::column_major > UBlasSparse
 

Functions

void fillMatrix (float density, int rows, int cols, EigenSparseMatrix &dst)
 
void fillMatrix2 (int nnzPerCol, int rows, int cols, EigenSparseMatrix &dst)
 
void eiToDense (const EigenSparseMatrix &src, DenseMatrix &dst)
 
void eiToGmm (const EigenSparseMatrix &src, GmmSparse &dst)
 
void eiToMtl (const EigenSparseMatrix &src, MtlSparse &dst)
 
void eiToUblas (const EigenSparseMatrix &src, UBlasSparse &dst)
 
template<typename EigenType , typename UblasType >
void eiToUblasVec (const EigenType &src, UblasType &dst)
 

Macro Definition Documentation

◆ DENSITY

#define DENSITY   0.01

◆ SCALAR

#define SCALAR   double

◆ SIZE

#define SIZE   1024

Typedef Documentation

◆ DenseMatrix

typedef Matrix<Scalar, Dynamic, Dynamic> DenseMatrix

◆ DenseVector

typedef Matrix<Scalar, Dynamic, 1> DenseVector

◆ EigenSparseMatrix

◆ GmmDynSparse

typedef gmm::col_matrix<gmm::wsvector<Scalar> > GmmDynSparse

◆ GmmSparse

typedef gmm::csc_matrix<Scalar> GmmSparse

◆ MtlSparse

typedef mtl::compressed2D<Scalar, mtl::matrix::parameters<mtl::tag::col_major> > MtlSparse

◆ MtlSparseRowMajor

typedef mtl::compressed2D<Scalar, mtl::matrix::parameters<mtl::tag::row_major> > MtlSparseRowMajor

◆ Scalar

typedef SCALAR Scalar

◆ UBlasSparse

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

Function Documentation

◆ eiToDense()

void eiToDense ( const EigenSparseMatrix src,
DenseMatrix dst 
)
54  {
55  dst.setZero();
56  for (int j = 0; j < src.cols(); ++j)
57  for (EigenSparseMatrix::InnerIterator it(src.derived(), j); it; ++it) dst(it.index(), j) = it.value();
58 }
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Definition: CwiseNullaryOp.h:569
const Derived & derived() const
Definition: SparseMatrixBase.h:144
Index cols() const
Definition: SparseMatrix.h:161
Base::InnerIterator InnerIterator
Definition: SparseMatrix.h:138
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), Eigen::SparseMatrixBase< Derived >::derived(), j, and Eigen::PlainObjectBase< Derived >::setZero().

Referenced by main().

◆ eiToGmm()

void eiToGmm ( const EigenSparseMatrix src,
GmmSparse dst 
)
64  {
65  GmmDynSparse tmp(src.rows(), src.cols());
66  for (int j = 0; j < src.cols(); ++j)
67  for (EigenSparseMatrix::InnerIterator it(src.derived(), j); it; ++it) tmp(it.index(), j) = it.value();
68  gmm::copy(tmp, dst);
69 }
gmm::col_matrix< gmm::wsvector< Scalar > > GmmDynSparse
Definition: BenchSparseUtil.h:63
Index rows() const
Definition: SparseMatrix.h:159
EIGEN_BLAS_FUNC() copy(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
Definition: level1_impl.h:32
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365

References Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), copy(), Eigen::SparseMatrixBase< Derived >::derived(), j, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::rows(), and tmp.

Referenced by __attribute__(), and main().

◆ eiToMtl()

void eiToMtl ( const EigenSparseMatrix src,
MtlSparse dst 
)
76  {
77  mtl::matrix::inserter<MtlSparse> ins(dst);
78  for (int j = 0; j < src.cols(); ++j)
79  for (EigenSparseMatrix::InnerIterator it(src.derived(), j); it; ++it) ins[it.index()][j] = it.value();
80 }

References Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), Eigen::SparseMatrixBase< Derived >::derived(), and j.

Referenced by main().

◆ eiToUblas()

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

References Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), Eigen::SparseMatrixBase< Derived >::derived(), j, and Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::rows().

Referenced by main().

◆ eiToUblasVec()

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

References j.

Referenced by main().

◆ fillMatrix()

void fillMatrix ( float  density,
int  rows,
int  cols,
EigenSparseMatrix dst 
)
27  {
28  dst.reserve(double(rows) * cols * density);
29  for (int j = 0; j < cols; j++) {
30  for (int i = 0; i < rows; i++) {
31  Scalar v = (internal::random<float>(0, 1) < density) ? internal::random<Scalar>() : 0;
32  if (v != 0) dst.insert(i, j) = v;
33  }
34  }
35  dst.finalize();
36 }
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
SCALAR Scalar
Definition: BenchSparseUtil.h:22
int i
Definition: BiCGSTAB_step_by_step.cpp:9
int rows
Definition: Tutorial_commainit_02.cpp:1
int cols
Definition: Tutorial_commainit_02.cpp:1
void finalize()
Definition: SparseMatrix.h:461
void reserve(Index reserveSize)
Definition: SparseMatrix.h:315
Scalar & insert(Index row, Index col)
Definition: SparseMatrix.h:1586
density
Definition: UniformPSDSelfTest.py:19

References cols, UniformPSDSelfTest::density, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::finalize(), i, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::insert(), j, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::reserve(), rows, and v.

Referenced by main().

◆ fillMatrix2()

void fillMatrix2 ( int  nnzPerCol,
int  rows,
int  cols,
EigenSparseMatrix dst 
)
38  {
39  // std::cout << "alloc " << nnzPerCol*cols << "\n";
40  dst.reserve(nnzPerCol * cols);
41  for (int j = 0; j < cols; j++) {
42  std::set<int> aux;
43  for (int i = 0; i < nnzPerCol; i++) {
44  int k = internal::random<int>(0, rows - 1);
45  while (aux.find(k) != aux.end()) k = internal::random<int>(0, rows - 1);
46  aux.insert(k);
47 
48  dst.insert(k, j) = internal::random<Scalar>();
49  }
50  }
51  dst.finalize();
52 }
char char char int int * k
Definition: level2_impl.h:374

References cols, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::finalize(), i, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::insert(), j, k, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::reserve(), and rows.

Referenced by main().