blob: 88cea09485401078b3b72bbf3030973c29a05211 [file] [log] [blame]
Gael Guennebaud709e9032009-01-07 17:01:57 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Gael Guennebaud709e9032009-01-07 17:01:57 +00003//
Gael Guennebaud22c76092011-03-22 11:58:22 +01004// Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
Gael Guennebaud709e9032009-01-07 17:01:57 +00005//
Benoit Jacob69124cf2012-07-13 14:42:47 -04006// 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/.
Gael Guennebaud709e9032009-01-07 17:01:57 +00009
10#include "sparse.h"
11
Gael Guennebaud869b4442016-01-25 11:55:39 +010012template<typename Scalar,typename StorageIndex> void sparse_vector(int rows, int cols)
Gael Guennebaud709e9032009-01-07 17:01:57 +000013{
Gael Guennebaud42e25782011-08-19 14:18:05 +020014 double densityMat = (std::max)(8./(rows*cols), 0.01);
Christoph Hertzbergdacb4692016-05-05 13:35:45 +020015 double densityVec = (std::max)(8./(rows), 0.1);
Gael Guennebaud709e9032009-01-07 17:01:57 +000016 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
17 typedef Matrix<Scalar,Dynamic,1> DenseVector;
Antonio Sánchezc614b2b2022-12-06 00:02:31 +000018 typedef Matrix<DenseIndex,Dynamic,1> DenseIndexVector;
Gael Guennebaud869b4442016-01-25 11:55:39 +010019 typedef SparseVector<Scalar,0,StorageIndex> SparseVectorType;
20 typedef SparseMatrix<Scalar,0,StorageIndex> SparseMatrixType;
Gael Guennebaud709e9032009-01-07 17:01:57 +000021 Scalar eps = 1e-6;
22
Gael Guennebaud2c03e6f2011-12-22 14:01:06 +010023 SparseMatrixType m1(rows,rows);
Gael Guennebaud709e9032009-01-07 17:01:57 +000024 SparseVectorType v1(rows), v2(rows), v3(rows);
Gael Guennebaud2c03e6f2011-12-22 14:01:06 +010025 DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
Gael Guennebaud709e9032009-01-07 17:01:57 +000026 DenseVector refV1 = DenseVector::Random(rows),
Gael Guennebaud349c2c92014-10-09 23:35:49 +020027 refV2 = DenseVector::Random(rows),
28 refV3 = DenseVector::Random(rows);
Gael Guennebaud709e9032009-01-07 17:01:57 +000029
30 std::vector<int> zerocoords, nonzerocoords;
31 initSparse<Scalar>(densityVec, refV1, v1, &zerocoords, &nonzerocoords);
32 initSparse<Scalar>(densityMat, refM1, m1);
33
34 initSparse<Scalar>(densityVec, refV2, v2);
35 initSparse<Scalar>(densityVec, refV3, v3);
36
Benoit Jacob47160402010-10-25 10:15:22 -040037 Scalar s1 = internal::random<Scalar>();
Gael Guennebaud709e9032009-01-07 17:01:57 +000038
39 // test coeff and coeffRef
40 for (unsigned int i=0; i<zerocoords.size(); ++i)
41 {
42 VERIFY_IS_MUCH_SMALLER_THAN( v1.coeff(zerocoords[i]), eps );
Gael Guennebaud36c478c2009-01-19 22:29:28 +000043 //VERIFY_RAISES_ASSERT( v1.coeffRef(zerocoords[i]) = 5 );
Gael Guennebaud709e9032009-01-07 17:01:57 +000044 }
45 {
46 VERIFY(int(nonzerocoords.size()) == v1.nonZeros());
47 int j=0;
48 for (typename SparseVectorType::InnerIterator it(v1); it; ++it,++j)
49 {
50 VERIFY(nonzerocoords[j]==it.index());
Erik Schultheisd271a7d2022-01-26 18:16:19 +000051 VERIFY_IS_EQUAL(it.value(), v1.coeff(it.index()));
52 VERIFY_IS_EQUAL(it.value(), refV1.coeff(it.index()));
Gael Guennebaud709e9032009-01-07 17:01:57 +000053 }
54 }
55 VERIFY_IS_APPROX(v1, refV1);
Gael Guennebaud349c2c92014-10-09 23:35:49 +020056
57 // test coeffRef with reallocation
58 {
Gael Guennebaud9a2447b2015-06-09 09:11:12 +020059 SparseVectorType v4(rows);
60 DenseVector v5 = DenseVector::Zero(rows);
Gael Guennebaud349c2c92014-10-09 23:35:49 +020061 for(int k=0; k<rows; ++k)
62 {
63 int i = internal::random<int>(0,rows-1);
64 Scalar v = internal::random<Scalar>();
Gael Guennebaud9a2447b2015-06-09 09:11:12 +020065 v4.coeffRef(i) += v;
66 v5.coeffRef(i) += v;
Gael Guennebaud349c2c92014-10-09 23:35:49 +020067 }
Gael Guennebaud9a2447b2015-06-09 09:11:12 +020068 VERIFY_IS_APPROX(v4,v5);
Gael Guennebaud349c2c92014-10-09 23:35:49 +020069 }
Gael Guennebaud709e9032009-01-07 17:01:57 +000070
71 v1.coeffRef(nonzerocoords[0]) = Scalar(5);
72 refV1.coeffRef(nonzerocoords[0]) = Scalar(5);
73 VERIFY_IS_APPROX(v1, refV1);
74
75 VERIFY_IS_APPROX(v1+v2, refV1+refV2);
76 VERIFY_IS_APPROX(v1+v2+v3, refV1+refV2+refV3);
77
78 VERIFY_IS_APPROX(v1*s1-v2, refV1*s1-refV2);
79
Gael Guennebaud2d534662009-01-14 21:27:54 +000080 VERIFY_IS_APPROX(v1*=s1, refV1*=s1);
81 VERIFY_IS_APPROX(v1/=s1, refV1/=s1);
Gael Guennebauda2324d62010-04-13 10:40:55 +020082
Gael Guennebaude7c48fa2009-01-23 13:59:32 +000083 VERIFY_IS_APPROX(v1+=v2, refV1+=refV2);
84 VERIFY_IS_APPROX(v1-=v2, refV1-=refV2);
Gael Guennebaud2d534662009-01-14 21:27:54 +000085
Gael Guennebaud709e9032009-01-07 17:01:57 +000086 VERIFY_IS_APPROX(v1.dot(v2), refV1.dot(refV2));
Gael Guennebauda9688f02009-02-09 09:59:30 +000087 VERIFY_IS_APPROX(v1.dot(refV2), refV1.dot(refV2));
Gael Guennebaud709e9032009-01-07 17:01:57 +000088
Gael Guennebaudfbb53b62014-09-01 16:53:52 +020089 VERIFY_IS_APPROX(m1*v2, refM1*refV2);
Gael Guennebaud2c03e6f2011-12-22 14:01:06 +010090 VERIFY_IS_APPROX(v1.dot(m1*v2), refV1.dot(refM1*refV2));
Gael Guennebaud869b4442016-01-25 11:55:39 +010091 {
92 int i = internal::random<int>(0,rows-1);
93 VERIFY_IS_APPROX(v1.dot(m1.col(i)), refV1.dot(refM1.col(i)));
94 }
Gael Guennebaud2c03e6f2011-12-22 14:01:06 +010095
96
Gael Guennebauda2324d62010-04-13 10:40:55 +020097 VERIFY_IS_APPROX(v1.squaredNorm(), refV1.squaredNorm());
Desire NUENTSA5b9bb002013-01-21 15:37:06 +010098
99 VERIFY_IS_APPROX(v1.blueNorm(), refV1.blueNorm());
Gael Guennebauda2324d62010-04-13 10:40:55 +0200100
Gael Guennebaude75b1eb2012-07-25 09:33:50 +0200101 // test aliasing
102 VERIFY_IS_APPROX((v1 = -v1), (refV1 = -refV1));
103 VERIFY_IS_APPROX((v1 = v1.transpose()), (refV1 = refV1.transpose().eval()));
104 VERIFY_IS_APPROX((v1 += -v1), (refV1 += -refV1));
Gael Guennebaud98ce4452013-03-06 11:58:22 +0100105
106 // sparse matrix to sparse vector
107 SparseMatrixType mv1;
108 VERIFY_IS_APPROX((mv1=v1),v1);
109 VERIFY_IS_APPROX(mv1,(v1=mv1));
110 VERIFY_IS_APPROX(mv1,(v1=mv1.transpose()));
Gael Guennebaude3929482013-06-10 00:06:40 +0200111
112 // check copy to dense vector with transpose
113 refV3.resize(0);
114 VERIFY_IS_APPROX(refV3 = v1.transpose(),v1.toDense());
Erik Schultheis13954c42021-11-15 22:16:01 +0000115 VERIFY_IS_APPROX(DenseVector(v1),v1.toDense());
Gael Guennebaude75b1eb2012-07-25 09:33:50 +0200116
Gael Guennebaud869b4442016-01-25 11:55:39 +0100117 // test conservative resize
118 {
119 std::vector<StorageIndex> inc;
120 if(rows > 3)
121 inc.push_back(-3);
122 inc.push_back(0);
123 inc.push_back(3);
124 inc.push_back(1);
125 inc.push_back(10);
126
127 for(std::size_t i = 0; i< inc.size(); i++) {
128 StorageIndex incRows = inc[i];
129 SparseVectorType vec1(rows);
130 DenseVector refVec1 = DenseVector::Zero(rows);
131 initSparse<Scalar>(densityVec, refVec1, vec1);
132
133 vec1.conservativeResize(rows+incRows);
134 refVec1.conservativeResize(rows+incRows);
135 if (incRows > 0) refVec1.tail(incRows).setZero();
136
137 VERIFY_IS_APPROX(vec1, refVec1);
138
139 // Insert new values
140 if (incRows > 0)
141 vec1.insert(vec1.rows()-1) = refVec1(refVec1.rows()-1) = 1;
142
143 VERIFY_IS_APPROX(vec1, refVec1);
144 }
145 }
146
Charles Schlosser44fe5392022-12-01 19:28:56 +0000147 // test sort
148 if(rows > 1)
149 {
150 SparseVectorType vec1(rows);
151 DenseVector refVec1 = DenseVector::Zero(rows);
Antonio Sánchezc614b2b2022-12-06 00:02:31 +0000152 DenseIndexVector innerIndices(rows);
Charles Schlosser44fe5392022-12-01 19:28:56 +0000153 innerIndices.setLinSpaced(0, rows - 1);
154 std::random_shuffle(innerIndices.begin(), innerIndices.end());
155 Index nz = internal::random<Index>(2, rows / 2);
156 for (Index k = 0; k < nz; k++)
157 {
158 Index i = innerIndices[k];
159 Scalar val = internal::random<Scalar>();
160 refVec1.coeffRef(i) = val;
161 vec1.insert(i) = val;
162 }
163
164 vec1.template sortInnerIndices<std::greater<>>();
165 VERIFY_IS_APPROX(vec1, refVec1);
166 VERIFY_IS_EQUAL(vec1.template innerIndicesAreSorted<std::greater<>>(), 1);
167 VERIFY_IS_EQUAL(vec1.template innerIndicesAreSorted<std::less<>>(), 0);
168 vec1.template sortInnerIndices<std::less<>>();
169 VERIFY_IS_APPROX(vec1, refVec1);
170 VERIFY_IS_EQUAL(vec1.template innerIndicesAreSorted<std::greater<>>(), 0);
171 VERIFY_IS_EQUAL(vec1.template innerIndicesAreSorted<std::less<>>(), 1);
172 }
173
Gael Guennebaud709e9032009-01-07 17:01:57 +0000174}
Erik Schultheis13954c42021-11-15 22:16:01 +0000175void test_pruning() {
176 using SparseVectorType = SparseVector<double, 0, int>;
177
178 SparseVectorType vec;
179 auto init_vec = [&](){;
180 vec.resize(10);
181 vec.insert(3) = 0.1;
182 vec.insert(5) = 1.0;
183 vec.insert(8) = -0.1;
184 vec.insert(9) = -0.2;
185 };
186 init_vec();
187
188 VERIFY_IS_EQUAL(vec.nonZeros(), 4);
189 VERIFY_IS_EQUAL(vec.prune(0.1, 1.0), 2);
190 VERIFY_IS_EQUAL(vec.nonZeros(), 2);
191 VERIFY_IS_EQUAL(vec.coeff(5), 1.0);
192 VERIFY_IS_EQUAL(vec.coeff(9), -0.2);
193
194 init_vec();
195 VERIFY_IS_EQUAL(vec.prune([](double v) { return v >= 0; }), 2);
196 VERIFY_IS_EQUAL(vec.nonZeros(), 2);
197 VERIFY_IS_EQUAL(vec.coeff(3), 0.1);
198 VERIFY_IS_EQUAL(vec.coeff(5), 1.0);
199}
Gael Guennebaud709e9032009-01-07 17:01:57 +0000200
Gael Guennebaud82f0ce22018-07-17 14:46:15 +0200201EIGEN_DECLARE_TEST(sparse_vector)
Gael Guennebaud709e9032009-01-07 17:01:57 +0000202{
203 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaud869b4442016-01-25 11:55:39 +0100204 int r = Eigen::internal::random<int>(1,500), c = Eigen::internal::random<int>(1,500);
205 if(Eigen::internal::random<int>(0,4) == 0) {
206 r = c; // check square matrices in 25% of tries
207 }
208 EIGEN_UNUSED_VARIABLE(r+c);
209
Gael Guennebaud4b6b3f32014-02-15 09:35:23 +0100210 CALL_SUBTEST_1(( sparse_vector<double,int>(8, 8) ));
Gael Guennebaud869b4442016-01-25 11:55:39 +0100211 CALL_SUBTEST_2(( sparse_vector<std::complex<double>, int>(r, c) ));
212 CALL_SUBTEST_1(( sparse_vector<double,long int>(r, c) ));
213 CALL_SUBTEST_1(( sparse_vector<double,short>(r, c) ));
Gael Guennebaud709e9032009-01-07 17:01:57 +0000214 }
Erik Schultheis13954c42021-11-15 22:16:01 +0000215
216 CALL_SUBTEST_1(test_pruning());
Gael Guennebaud709e9032009-01-07 17:01:57 +0000217}
218