sparse.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_TESTSPARSE_H
11 #define EIGEN_TESTSPARSE_H
12 
13 #define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET
14 
15 #include "main.h"
16 
17 #ifdef min
18 #undef min
19 #endif
20 
21 #ifdef max
22 #undef max
23 #endif
24 
25 #include <unordered_map>
26 #define EIGEN_UNORDERED_MAP_SUPPORT
27 
28 #include <Eigen/Cholesky>
29 #include <Eigen/LU>
30 #include <Eigen/Sparse>
31 
33 
34 /* Initializes both a sparse and dense matrix with same random values,
35  * and a ratio of \a density non zero entries.
36  * \param flags is a union of ForceNonZeroDiag, MakeLowerTriangular and MakeUpperTriangular
37  * allowing to control the shape of the matrix.
38  * \param zeroCoords and nonzeroCoords allows to get the coordinate lists of the non zero,
39  * and zero coefficients respectively.
40  */
41 template <typename Scalar, int Opt1, int Opt2, typename StorageIndex>
43  SparseMatrix<Scalar, Opt2, StorageIndex>& sparseMat, int flags = 0,
44  std::vector<Matrix<StorageIndex, 2, 1> >* zeroCoords = 0,
45  std::vector<Matrix<StorageIndex, 2, 1> >* nonzeroCoords = 0) {
47  sparseMat.setZero();
48  // sparseMat.reserve(int(refMat.rows()*refMat.cols()*density));
49  int nnz = static_cast<int>((1.5 * density) * static_cast<double>(IsRowMajor ? refMat.cols() : refMat.rows()));
50  sparseMat.reserve(VectorXi::Constant(IsRowMajor ? refMat.rows() : refMat.cols(), nnz));
51 
52  Index insert_count = 0;
53  for (Index j = 0; j < sparseMat.outerSize(); j++) {
54  // sparseMat.startVec(j);
55  for (Index i = 0; i < sparseMat.innerSize(); i++) {
56  Index ai(i), aj(j);
57  if (IsRowMajor) std::swap(ai, aj);
58  Scalar v = (internal::random<double>(0, 1) < density) ? internal::random<Scalar>() : Scalar(0);
59  if ((flags & ForceNonZeroDiag) && (i == j)) {
60  // FIXME: the following is too conservative
61  v = internal::random<Scalar>() * Scalar(3.);
62  v = v * v;
63  if (numext::real(v) > 0)
64  v += Scalar(5);
65  else
66  v -= Scalar(5);
67  }
68  if ((flags & MakeLowerTriangular) && aj > ai)
69  v = Scalar(0);
70  else if ((flags & MakeUpperTriangular) && aj < ai)
71  v = Scalar(0);
72 
73  if ((flags & ForceRealDiag) && (i == j)) v = numext::real(v);
74 
75  if (!numext::is_exactly_zero(v)) {
76  // sparseMat.insertBackByOuterInner(j,i) = v;
77  sparseMat.insertByOuterInner(j, i) = v;
78  ++insert_count;
79  if (nonzeroCoords) nonzeroCoords->push_back(Matrix<StorageIndex, 2, 1>(ai, aj));
80  } else if (zeroCoords) {
81  zeroCoords->push_back(Matrix<StorageIndex, 2, 1>(ai, aj));
82  }
83  refMat(ai, aj) = v;
84 
85  // make sure we only insert as many as the sparse matrix supports
86  if (insert_count == NumTraits<StorageIndex>::highest()) return;
87  }
88  }
89  // sparseMat.finalize();
90 }
91 
92 template <typename Scalar, int Options, typename Index>
94  std::vector<int>* zeroCoords = 0, std::vector<int>* nonzeroCoords = 0) {
95  sparseVec.reserve(int(refVec.size() * density));
96  sparseVec.setZero();
97  for (int i = 0; i < refVec.size(); i++) {
98  Scalar v = (internal::random<double>(0, 1) < density) ? internal::random<Scalar>() : Scalar(0);
99  if (!numext::is_exactly_zero(v)) {
100  sparseVec.insertBack(i) = v;
101  if (nonzeroCoords) nonzeroCoords->push_back(i);
102  } else if (zeroCoords)
103  zeroCoords->push_back(i);
104  refVec[i] = v;
105  }
106 }
107 
108 template <typename Scalar, int Options, typename Index>
110  std::vector<int>* zeroCoords = 0, std::vector<int>* nonzeroCoords = 0) {
111  sparseVec.reserve(int(refVec.size() * density));
112  sparseVec.setZero();
113  for (int i = 0; i < refVec.size(); i++) {
114  Scalar v = (internal::random<double>(0, 1) < density) ? internal::random<Scalar>() : Scalar(0);
115  if (v != Scalar(0)) {
116  sparseVec.insertBack(i) = v;
117  if (nonzeroCoords) nonzeroCoords->push_back(i);
118  } else if (zeroCoords)
119  zeroCoords->push_back(i);
120  refVec[i] = v;
121  }
122 }
123 
124 #endif // EIGEN_TESTSPARSE_H
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
int i
Definition: BiCGSTAB_step_by_step.cpp:9
SCALAR Scalar
Definition: bench_gemm.cpp:45
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:192
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:191
A versatible sparse matrix representation.
Definition: SparseMatrix.h:121
void setZero()
Definition: SparseMatrix.h:303
Index outerSize() const
Definition: SparseMatrix.h:166
Scalar & insertByOuterInner(Index j, Index i)
Definition: SparseMatrix.h:566
Index innerSize() const
Definition: SparseMatrix.h:164
void reserve(Index reserveSize)
Definition: SparseMatrix.h:315
a sparse vector class
Definition: SparseVector.h:62
Scalar & insertBack(Index i)
Definition: SparseVector.h:142
void setZero()
Definition: SparseVector.h:127
void reserve(Index reserveSize)
Definition: SparseVector.h:186
float real
Definition: datatypes.h:10
EIGEN_BLAS_FUNC() swap(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
Definition: level1_impl.h:117
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bool is_exactly_zero(const X &x)
Definition: Meta.h:592
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
density
Definition: UniformPSDSelfTest.py:19
@ ForceRealDiag
Definition: sparse.h:32
@ ForceNonZeroDiag
Definition: sparse.h:32
@ MakeUpperTriangular
Definition: sparse.h:32
@ MakeLowerTriangular
Definition: sparse.h:32
void initSparse(double density, Matrix< Scalar, Dynamic, Dynamic, Opt1 > &refMat, SparseMatrix< Scalar, Opt2, StorageIndex > &sparseMat, int flags=0, std::vector< Matrix< StorageIndex, 2, 1 > > *zeroCoords=0, std::vector< Matrix< StorageIndex, 2, 1 > > *nonzeroCoords=0)
Definition: sparse.h:42
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2