blob: 7bd57cdbea6ebeb58f4ae95bf7e8045dddedefc6 [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;
Gael Guennebaud869b4442016-01-25 11:55:39 +010018 typedef SparseVector<Scalar,0,StorageIndex> SparseVectorType;
19 typedef SparseMatrix<Scalar,0,StorageIndex> SparseMatrixType;
Gael Guennebaud709e9032009-01-07 17:01:57 +000020 Scalar eps = 1e-6;
21
Gael Guennebaud2c03e6f2011-12-22 14:01:06 +010022 SparseMatrixType m1(rows,rows);
Gael Guennebaud709e9032009-01-07 17:01:57 +000023 SparseVectorType v1(rows), v2(rows), v3(rows);
Gael Guennebaud2c03e6f2011-12-22 14:01:06 +010024 DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
Gael Guennebaud709e9032009-01-07 17:01:57 +000025 DenseVector refV1 = DenseVector::Random(rows),
Gael Guennebaud349c2c92014-10-09 23:35:49 +020026 refV2 = DenseVector::Random(rows),
27 refV3 = DenseVector::Random(rows);
Gael Guennebaud709e9032009-01-07 17:01:57 +000028
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
Benoit Jacob47160402010-10-25 10:15:22 -040036 Scalar s1 = internal::random<Scalar>();
Gael Guennebaud709e9032009-01-07 17:01:57 +000037
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 );
Gael Guennebaud36c478c2009-01-19 22:29:28 +000042 //VERIFY_RAISES_ASSERT( v1.coeffRef(zerocoords[i]) = 5 );
Gael Guennebaud709e9032009-01-07 17:01:57 +000043 }
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());
Gael Guennebaudc4c70662009-01-14 14:24:10 +000050 VERIFY(it.value()==v1.coeff(it.index()));
51 VERIFY(it.value()==refV1.coeff(it.index()));
Gael Guennebaud709e9032009-01-07 17:01:57 +000052 }
53 }
54 VERIFY_IS_APPROX(v1, refV1);
Gael Guennebaud349c2c92014-10-09 23:35:49 +020055
56 // test coeffRef with reallocation
57 {
Gael Guennebaud9a2447b2015-06-09 09:11:12 +020058 SparseVectorType v4(rows);
59 DenseVector v5 = DenseVector::Zero(rows);
Gael Guennebaud349c2c92014-10-09 23:35:49 +020060 for(int k=0; k<rows; ++k)
61 {
62 int i = internal::random<int>(0,rows-1);
63 Scalar v = internal::random<Scalar>();
Gael Guennebaud9a2447b2015-06-09 09:11:12 +020064 v4.coeffRef(i) += v;
65 v5.coeffRef(i) += v;
Gael Guennebaud349c2c92014-10-09 23:35:49 +020066 }
Gael Guennebaud9a2447b2015-06-09 09:11:12 +020067 VERIFY_IS_APPROX(v4,v5);
Gael Guennebaud349c2c92014-10-09 23:35:49 +020068 }
Gael Guennebaud709e9032009-01-07 17:01:57 +000069
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
Gael Guennebaud2d534662009-01-14 21:27:54 +000079 VERIFY_IS_APPROX(v1*=s1, refV1*=s1);
80 VERIFY_IS_APPROX(v1/=s1, refV1/=s1);
Gael Guennebauda2324d62010-04-13 10:40:55 +020081
Gael Guennebaude7c48fa2009-01-23 13:59:32 +000082 VERIFY_IS_APPROX(v1+=v2, refV1+=refV2);
83 VERIFY_IS_APPROX(v1-=v2, refV1-=refV2);
Gael Guennebaud2d534662009-01-14 21:27:54 +000084
Gael Guennebaud709e9032009-01-07 17:01:57 +000085 VERIFY_IS_APPROX(v1.dot(v2), refV1.dot(refV2));
Gael Guennebauda9688f02009-02-09 09:59:30 +000086 VERIFY_IS_APPROX(v1.dot(refV2), refV1.dot(refV2));
Gael Guennebaud709e9032009-01-07 17:01:57 +000087
Gael Guennebaudfbb53b62014-09-01 16:53:52 +020088 VERIFY_IS_APPROX(m1*v2, refM1*refV2);
Gael Guennebaud2c03e6f2011-12-22 14:01:06 +010089 VERIFY_IS_APPROX(v1.dot(m1*v2), refV1.dot(refM1*refV2));
Gael Guennebaud869b4442016-01-25 11:55:39 +010090 {
91 int i = internal::random<int>(0,rows-1);
92 VERIFY_IS_APPROX(v1.dot(m1.col(i)), refV1.dot(refM1.col(i)));
93 }
Gael Guennebaud2c03e6f2011-12-22 14:01:06 +010094
95
Gael Guennebauda2324d62010-04-13 10:40:55 +020096 VERIFY_IS_APPROX(v1.squaredNorm(), refV1.squaredNorm());
Desire NUENTSA5b9bb002013-01-21 15:37:06 +010097
98 VERIFY_IS_APPROX(v1.blueNorm(), refV1.blueNorm());
Gael Guennebauda2324d62010-04-13 10:40:55 +020099
Gael Guennebaude75b1eb2012-07-25 09:33:50 +0200100 // 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));
Gael Guennebaud98ce4452013-03-06 11:58:22 +0100104
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()));
Gael Guennebaude3929482013-06-10 00:06:40 +0200110
111 // check copy to dense vector with transpose
112 refV3.resize(0);
113 VERIFY_IS_APPROX(refV3 = v1.transpose(),v1.toDense());
Erik Schultheis13954c42021-11-15 22:16:01 +0000114 VERIFY_IS_APPROX(DenseVector(v1),v1.toDense());
Gael Guennebaude75b1eb2012-07-25 09:33:50 +0200115
Gael Guennebaud869b4442016-01-25 11:55:39 +0100116 // 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
Gael Guennebaud709e9032009-01-07 17:01:57 +0000146}
Erik Schultheis13954c42021-11-15 22:16:01 +0000147void test_pruning() {
148 using SparseVectorType = SparseVector<double, 0, int>;
149
150 SparseVectorType vec;
151 auto init_vec = [&](){;
152 vec.resize(10);
153 vec.insert(3) = 0.1;
154 vec.insert(5) = 1.0;
155 vec.insert(8) = -0.1;
156 vec.insert(9) = -0.2;
157 };
158 init_vec();
159
160 VERIFY_IS_EQUAL(vec.nonZeros(), 4);
161 VERIFY_IS_EQUAL(vec.prune(0.1, 1.0), 2);
162 VERIFY_IS_EQUAL(vec.nonZeros(), 2);
163 VERIFY_IS_EQUAL(vec.coeff(5), 1.0);
164 VERIFY_IS_EQUAL(vec.coeff(9), -0.2);
165
166 init_vec();
167 VERIFY_IS_EQUAL(vec.prune([](double v) { return v >= 0; }), 2);
168 VERIFY_IS_EQUAL(vec.nonZeros(), 2);
169 VERIFY_IS_EQUAL(vec.coeff(3), 0.1);
170 VERIFY_IS_EQUAL(vec.coeff(5), 1.0);
171}
Gael Guennebaud709e9032009-01-07 17:01:57 +0000172
Gael Guennebaud82f0ce22018-07-17 14:46:15 +0200173EIGEN_DECLARE_TEST(sparse_vector)
Gael Guennebaud709e9032009-01-07 17:01:57 +0000174{
175 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaud869b4442016-01-25 11:55:39 +0100176 int r = Eigen::internal::random<int>(1,500), c = Eigen::internal::random<int>(1,500);
177 if(Eigen::internal::random<int>(0,4) == 0) {
178 r = c; // check square matrices in 25% of tries
179 }
180 EIGEN_UNUSED_VARIABLE(r+c);
181
Gael Guennebaud4b6b3f32014-02-15 09:35:23 +0100182 CALL_SUBTEST_1(( sparse_vector<double,int>(8, 8) ));
Gael Guennebaud869b4442016-01-25 11:55:39 +0100183 CALL_SUBTEST_2(( sparse_vector<std::complex<double>, int>(r, c) ));
184 CALL_SUBTEST_1(( sparse_vector<double,long int>(r, c) ));
185 CALL_SUBTEST_1(( sparse_vector<double,short>(r, c) ));
Gael Guennebaud709e9032009-01-07 17:01:57 +0000186 }
Erik Schultheis13954c42021-11-15 22:16:01 +0000187
188 CALL_SUBTEST_1(test_pruning());
Gael Guennebaud709e9032009-01-07 17:01:57 +0000189}
190