sparse_vector.cpp
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 #include "sparse.h"
11 
12 template<typename Scalar,typename StorageIndex> void sparse_vector(int rows, int cols)
13 {
14  double densityMat = (std::max)(8./(rows*cols), 0.01);
15  double densityVec = (std::max)(8./(rows), 0.1);
18  typedef SparseVector<Scalar,0,StorageIndex> SparseVectorType;
19  typedef SparseMatrix<Scalar,0,StorageIndex> SparseMatrixType;
20  Scalar eps = 1e-6;
21 
22  SparseMatrixType m1(rows,rows);
23  SparseVectorType v1(rows), v2(rows), v3(rows);
24  DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
25  DenseVector refV1 = DenseVector::Random(rows),
26  refV2 = DenseVector::Random(rows),
27  refV3 = DenseVector::Random(rows);
28 
29  std::vector<int> zerocoords, nonzerocoords;
30  initSparse<Scalar>(densityVec, refV1, v1, &zerocoords, &nonzerocoords);
31  initSparse<Scalar>(densityMat, refM1, m1);
32 
33  initSparse<Scalar>(densityVec, refV2, v2);
34  initSparse<Scalar>(densityVec, refV3, v3);
35 
36  Scalar s1 = internal::random<Scalar>();
37 
38  // test coeff and coeffRef
39  for (unsigned int i=0; i<zerocoords.size(); ++i)
40  {
41  VERIFY_IS_MUCH_SMALLER_THAN( v1.coeff(zerocoords[i]), eps );
42  //VERIFY_RAISES_ASSERT( v1.coeffRef(zerocoords[i]) = 5 );
43  }
44  {
45  VERIFY(int(nonzerocoords.size()) == v1.nonZeros());
46  int j=0;
47  for (typename SparseVectorType::InnerIterator it(v1); it; ++it,++j)
48  {
49  VERIFY(nonzerocoords[j]==it.index());
50  VERIFY(it.value()==v1.coeff(it.index()));
51  VERIFY(it.value()==refV1.coeff(it.index()));
52  }
53  }
54  VERIFY_IS_APPROX(v1, refV1);
55 
56  // test coeffRef with reallocation
57  {
58  SparseVectorType v4(rows);
59  DenseVector v5 = DenseVector::Zero(rows);
60  for(int k=0; k<rows; ++k)
61  {
62  int i = internal::random<int>(0,rows-1);
63  Scalar v = internal::random<Scalar>();
64  v4.coeffRef(i) += v;
65  v5.coeffRef(i) += v;
66  }
67  VERIFY_IS_APPROX(v4,v5);
68  }
69 
70  v1.coeffRef(nonzerocoords[0]) = Scalar(5);
71  refV1.coeffRef(nonzerocoords[0]) = Scalar(5);
72  VERIFY_IS_APPROX(v1, refV1);
73 
74  VERIFY_IS_APPROX(v1+v2, refV1+refV2);
75  VERIFY_IS_APPROX(v1+v2+v3, refV1+refV2+refV3);
76 
77  VERIFY_IS_APPROX(v1*s1-v2, refV1*s1-refV2);
78 
79  VERIFY_IS_APPROX(v1*=s1, refV1*=s1);
80  VERIFY_IS_APPROX(v1/=s1, refV1/=s1);
81 
82  VERIFY_IS_APPROX(v1+=v2, refV1+=refV2);
83  VERIFY_IS_APPROX(v1-=v2, refV1-=refV2);
84 
85  VERIFY_IS_APPROX(v1.dot(v2), refV1.dot(refV2));
86  VERIFY_IS_APPROX(v1.dot(refV2), refV1.dot(refV2));
87 
88  VERIFY_IS_APPROX(m1*v2, refM1*refV2);
89  VERIFY_IS_APPROX(v1.dot(m1*v2), refV1.dot(refM1*refV2));
90  {
91  int i = internal::random<int>(0,rows-1);
92  VERIFY_IS_APPROX(v1.dot(m1.col(i)), refV1.dot(refM1.col(i)));
93  }
94 
95 
96  VERIFY_IS_APPROX(v1.squaredNorm(), refV1.squaredNorm());
97 
98  VERIFY_IS_APPROX(v1.blueNorm(), refV1.blueNorm());
99 
100  // test aliasing
101  VERIFY_IS_APPROX((v1 = -v1), (refV1 = -refV1));
102  VERIFY_IS_APPROX((v1 = v1.transpose()), (refV1 = refV1.transpose().eval()));
103  VERIFY_IS_APPROX((v1 += -v1), (refV1 += -refV1));
104 
105  // sparse matrix to sparse vector
106  SparseMatrixType mv1;
107  VERIFY_IS_APPROX((mv1=v1),v1);
108  VERIFY_IS_APPROX(mv1,(v1=mv1));
109  VERIFY_IS_APPROX(mv1,(v1=mv1.transpose()));
110 
111  // check copy to dense vector with transpose
112  refV3.resize(0);
113  VERIFY_IS_APPROX(refV3 = v1.transpose(),v1.toDense());
114  VERIFY_IS_APPROX(DenseVector(v1),v1.toDense());
115 
116  // test conservative resize
117  {
118  std::vector<StorageIndex> inc;
119  if(rows > 3)
120  inc.push_back(-3);
121  inc.push_back(0);
122  inc.push_back(3);
123  inc.push_back(1);
124  inc.push_back(10);
125 
126  for(std::size_t i = 0; i< inc.size(); i++) {
127  StorageIndex incRows = inc[i];
128  SparseVectorType vec1(rows);
129  DenseVector refVec1 = DenseVector::Zero(rows);
130  initSparse<Scalar>(densityVec, refVec1, vec1);
131 
132  vec1.conservativeResize(rows+incRows);
133  refVec1.conservativeResize(rows+incRows);
134  if (incRows > 0) refVec1.tail(incRows).setZero();
135 
136  VERIFY_IS_APPROX(vec1, refVec1);
137 
138  // Insert new values
139  if (incRows > 0)
140  vec1.insert(vec1.rows()-1) = refVec1(refVec1.rows()-1) = 1;
141 
142  VERIFY_IS_APPROX(vec1, refVec1);
143  }
144  }
145 
146 }
147 
149 {
150  for(int i = 0; i < g_repeat; i++) {
151  int r = Eigen::internal::random<int>(1,500), c = Eigen::internal::random<int>(1,500);
152  if(Eigen::internal::random<int>(0,4) == 0) {
153  r = c; // check square matrices in 25% of tries
154  }
156 
157  CALL_SUBTEST_1(( sparse_vector<double,int>(8, 8) ));
158  CALL_SUBTEST_2(( sparse_vector<std::complex<double>, int>(r, c) ));
159  CALL_SUBTEST_1(( sparse_vector<double,long int>(r, c) ));
160  CALL_SUBTEST_1(( sparse_vector<double,short>(r, c) ));
161  }
162 }
163 
VERIFY_IS_MUCH_SMALLER_THAN
#define VERIFY_IS_MUCH_SMALLER_THAN(a, b)
Definition: main.h:390
Eigen::SparseMatrix
A versatible sparse matrix representation.
Definition: SparseMatrix.h:96
e
Array< double, 1, 3 > e(1./3., 0.5, 2.)
c
Scalar Scalar * c
Definition: benchVecAdd.cpp:17
m1
Matrix3d m1
Definition: IOFormat.cpp:2
EIGEN_DECLARE_TEST
EIGEN_DECLARE_TEST(sparse_vector)
Definition: sparse_vector.cpp:148
rows
int rows
Definition: Tutorial_commainit_02.cpp:1
CALL_SUBTEST_1
#define CALL_SUBTEST_1(FUNC)
Definition: split_test_helper.h:4
EIGEN_UNUSED_VARIABLE
#define EIGEN_UNUSED_VARIABLE(var)
Definition: Macros.h:1076
j
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2
Eigen::g_repeat
static int g_repeat
Definition: main.h:169
CALL_SUBTEST_2
#define CALL_SUBTEST_2(FUNC)
Definition: split_test_helper.h:10
size_t
std::size_t size_t
Definition: wrap/pybind11/include/pybind11/detail/common.h:490
DenseVector
Matrix< Scalar, Dynamic, 1 > DenseVector
Definition: BenchSparseUtil.h:24
VERIFY_IS_APPROX
#define VERIFY_IS_APPROX(a, b)
Definition: integer_types.cpp:15
Eigen::SparseVector
a sparse vector class
Definition: SparseUtil.h:54
v2
Vector v2
Definition: testSerializationBase.cpp:39
v
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
DenseMatrix
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
Definition: BenchSparseUtil.h:23
Eigen::Matrix< Scalar, Dynamic, Dynamic >
v3
Vector v3
Definition: testSerializationBase.cpp:40
cols
int cols
Definition: Tutorial_commainit_02.cpp:1
max
#define max(a, b)
Definition: datatypes.h:20
sparse_vector
void sparse_vector(int rows, int cols)
Definition: sparse_vector.cpp:12
vec1
RowVectorXd vec1(3)
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
sparse.h
v1
Vector v1
Definition: testSerializationBase.cpp:38
Scalar
SCALAR Scalar
Definition: bench_gemm.cpp:46
VERIFY
#define VERIFY(a)
Definition: main.h:380


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:04:24