blob: d4579e4c9c9ab2bd9692286921fd45f09e9afc4b [file] [log] [blame]
Gael Guennebaud86ccd992008-11-05 13:47:55 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Gael Guennebaud86ccd992008-11-05 13:47:55 +00003//
Gael Guennebaud22c76092011-03-22 11:58:22 +01004// Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
Gael Guennebaud86ccd992008-11-05 13:47:55 +00005// Copyright (C) 2008 Daniel Gomez Ferro <dgomezferro@gmail.com>
6//
7// Eigen is free software; you can redistribute it and/or
8// modify it under the terms of the GNU Lesser General Public
9// License as published by the Free Software Foundation; either
10// version 3 of the License, or (at your option) any later version.
11//
12// Alternatively, you can redistribute it and/or
13// modify it under the terms of the GNU General Public License as
14// published by the Free Software Foundation; either version 2 of
15// the License, or (at your option) any later version.
16//
17// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
18// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
20// GNU General Public License for more details.
21//
22// You should have received a copy of the GNU Lesser General Public
23// License and a copy of the GNU General Public License along with
24// Eigen. If not, see <http://www.gnu.org/licenses/>.
25
26#include "sparse.h"
27
Gael Guennebaud178858f2009-01-19 15:20:45 +000028template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& ref)
29{
Hauke Heibelf1679c72010-06-20 17:37:56 +020030 typedef typename SparseMatrixType::Index Index;
31
32 const Index rows = ref.rows();
33 const Index cols = ref.cols();
Gael Guennebaud178858f2009-01-19 15:20:45 +000034 typedef typename SparseMatrixType::Scalar Scalar;
35 enum { Flags = SparseMatrixType::Flags };
Gael Guennebaud9f795582009-12-16 19:18:40 +010036
Gael Guennebaud42e25782011-08-19 14:18:05 +020037 double density = (std::max)(8./(rows*cols), 0.01);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000038 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
39 typedef Matrix<Scalar,Dynamic,1> DenseVector;
40 Scalar eps = 1e-6;
41
Gael Guennebaud178858f2009-01-19 15:20:45 +000042 SparseMatrixType m(rows, cols);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000043 DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
44 DenseVector vec1 = DenseVector::Random(rows);
Benoit Jacob47160402010-10-25 10:15:22 -040045 Scalar s1 = internal::random<Scalar>();
Gael Guennebaud86ccd992008-11-05 13:47:55 +000046
47 std::vector<Vector2i> zeroCoords;
48 std::vector<Vector2i> nonzeroCoords;
49 initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords);
Gael Guennebaud9f795582009-12-16 19:18:40 +010050
Gael Guennebaud86ccd992008-11-05 13:47:55 +000051 if (zeroCoords.size()==0 || nonzeroCoords.size()==0)
52 return;
53
54 // test coeff and coeffRef
55 for (int i=0; i<(int)zeroCoords.size(); ++i)
56 {
57 VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps );
Hauke Heibel7bc8e3a2010-10-25 22:13:49 +020058 if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value)
Gael Guennebaud178858f2009-01-19 15:20:45 +000059 VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[0].x(),zeroCoords[0].y()) = 5 );
Gael Guennebaud86ccd992008-11-05 13:47:55 +000060 }
61 VERIFY_IS_APPROX(m, refMat);
62
63 m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
64 refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
65
66 VERIFY_IS_APPROX(m, refMat);
Gael Guennebaudc4c70662009-01-14 14:24:10 +000067 /*
Gael Guennebaud86ccd992008-11-05 13:47:55 +000068 // test InnerIterators and Block expressions
69 for (int t=0; t<10; ++t)
70 {
Benoit Jacob47160402010-10-25 10:15:22 -040071 int j = internal::random<int>(0,cols-1);
72 int i = internal::random<int>(0,rows-1);
73 int w = internal::random<int>(1,cols-j-1);
74 int h = internal::random<int>(1,rows-i-1);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000075
Gael Guennebaudc4c70662009-01-14 14:24:10 +000076// VERIFY_IS_APPROX(m.block(i,j,h,w), refMat.block(i,j,h,w));
Gael Guennebaud86ccd992008-11-05 13:47:55 +000077 for(int c=0; c<w; c++)
78 {
79 VERIFY_IS_APPROX(m.block(i,j,h,w).col(c), refMat.block(i,j,h,w).col(c));
80 for(int r=0; r<h; r++)
81 {
Gael Guennebaudc4c70662009-01-14 14:24:10 +000082// VERIFY_IS_APPROX(m.block(i,j,h,w).col(c).coeff(r), refMat.block(i,j,h,w).col(c).coeff(r));
Gael Guennebaud86ccd992008-11-05 13:47:55 +000083 }
84 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +000085// for(int r=0; r<h; r++)
86// {
87// VERIFY_IS_APPROX(m.block(i,j,h,w).row(r), refMat.block(i,j,h,w).row(r));
88// for(int c=0; c<w; c++)
89// {
90// VERIFY_IS_APPROX(m.block(i,j,h,w).row(r).coeff(c), refMat.block(i,j,h,w).row(r).coeff(c));
91// }
92// }
Gael Guennebaud86ccd992008-11-05 13:47:55 +000093 }
94
95 for(int c=0; c<cols; c++)
96 {
97 VERIFY_IS_APPROX(m.col(c) + m.col(c), (m + m).col(c));
98 VERIFY_IS_APPROX(m.col(c) + m.col(c), refMat.col(c) + refMat.col(c));
99 }
100
101 for(int r=0; r<rows; r++)
102 {
103 VERIFY_IS_APPROX(m.row(r) + m.row(r), (m + m).row(r));
104 VERIFY_IS_APPROX(m.row(r) + m.row(r), refMat.row(r) + refMat.row(r));
105 }
106 */
107
Gael Guennebaud28293142009-05-04 14:25:12 +0000108 // test insert (inner random)
Gael Guennebaud5015e482008-12-11 18:26:24 +0000109 {
110 DenseMatrix m1(rows,cols);
111 m1.setZero();
Gael Guennebaud178858f2009-01-19 15:20:45 +0000112 SparseMatrixType m2(rows,cols);
Gael Guennebaud28293142009-05-04 14:25:12 +0000113 m2.reserve(10);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000114 for (int j=0; j<cols; ++j)
115 {
116 for (int k=0; k<rows/2; ++k)
117 {
Benoit Jacob47160402010-10-25 10:15:22 -0400118 int i = internal::random<int>(0,rows-1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000119 if (m1.coeff(i,j)==Scalar(0))
Benoit Jacob47160402010-10-25 10:15:22 -0400120 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud5015e482008-12-11 18:26:24 +0000121 }
122 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000123 m2.finalize();
124 VERIFY_IS_APPROX(m2,m1);
125 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100126
Gael Guennebaud28293142009-05-04 14:25:12 +0000127 // test insert (fully random)
128 {
129 DenseMatrix m1(rows,cols);
130 m1.setZero();
131 SparseMatrixType m2(rows,cols);
132 m2.reserve(10);
133 for (int k=0; k<rows*cols; ++k)
134 {
Benoit Jacob47160402010-10-25 10:15:22 -0400135 int i = internal::random<int>(0,rows-1);
136 int j = internal::random<int>(0,cols-1);
Gael Guennebaud28293142009-05-04 14:25:12 +0000137 if (m1.coeff(i,j)==Scalar(0))
Benoit Jacob47160402010-10-25 10:15:22 -0400138 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud28293142009-05-04 14:25:12 +0000139 }
140 m2.finalize();
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000141 VERIFY_IS_APPROX(m2,m1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000142 }
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200143
144 // test insert (un-compressed)
145 for(int mode=0;mode<4;++mode)
146 {
147 DenseMatrix m1(rows,cols);
148 m1.setZero();
149 SparseMatrixType m2(rows,cols);
150 VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? m2.innerSize() : std::max<int>(1,m2.innerSize()/8)));
151 m2.reserve(r);
152 for (int k=0; k<rows*cols; ++k)
153 {
154 int i = internal::random<int>(0,rows-1);
155 int j = internal::random<int>(0,cols-1);
156 if (m1.coeff(i,j)==Scalar(0))
157 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
158 if(mode==3)
159 m2.reserve(r);
160 }
161 m2.finalize();
162 m2.makeCompressed();
163 VERIFY_IS_APPROX(m2,m1);
164 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100165
Gael Guennebaud2d534662009-01-14 21:27:54 +0000166 // test basic computations
167 {
168 DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
169 DenseMatrix refM2 = DenseMatrix::Zero(rows, rows);
170 DenseMatrix refM3 = DenseMatrix::Zero(rows, rows);
171 DenseMatrix refM4 = DenseMatrix::Zero(rows, rows);
Gael Guennebaud178858f2009-01-19 15:20:45 +0000172 SparseMatrixType m1(rows, rows);
173 SparseMatrixType m2(rows, rows);
174 SparseMatrixType m3(rows, rows);
175 SparseMatrixType m4(rows, rows);
Gael Guennebaud2d534662009-01-14 21:27:54 +0000176 initSparse<Scalar>(density, refM1, m1);
177 initSparse<Scalar>(density, refM2, m2);
178 initSparse<Scalar>(density, refM3, m3);
179 initSparse<Scalar>(density, refM4, m4);
180
181 VERIFY_IS_APPROX(m1+m2, refM1+refM2);
182 VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3);
Gael Guennebaud9f795582009-12-16 19:18:40 +0100183 VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2));
Gael Guennebaud2d534662009-01-14 21:27:54 +0000184 VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2);
185
186 VERIFY_IS_APPROX(m1*=s1, refM1*=s1);
187 VERIFY_IS_APPROX(m1/=s1, refM1/=s1);
Gael Guennebaud9f795582009-12-16 19:18:40 +0100188
Gael Guennebaude7c48fa2009-01-23 13:59:32 +0000189 VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
190 VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
Gael Guennebaud9f795582009-12-16 19:18:40 +0100191
Gael Guennebauda9688f02009-02-09 09:59:30 +0000192 VERIFY_IS_APPROX(m1.col(0).dot(refM2.row(0)), refM1.col(0).dot(refM2.row(0)));
Gael Guennebaud9f795582009-12-16 19:18:40 +0100193
Gael Guennebaud2d534662009-01-14 21:27:54 +0000194 refM4.setRandom();
195 // sparse cwise* dense
Gael Guennebaud9f795582009-12-16 19:18:40 +0100196 VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4));
Gael Guennebaud2d534662009-01-14 21:27:54 +0000197// VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4);
198 }
199
Gael Guennebaudece48a62010-06-18 11:28:30 +0200200 // test transpose
201 {
202 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
203 SparseMatrixType m2(rows, rows);
204 initSparse<Scalar>(density, refMat2, m2);
205 VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval());
206 VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
207
208 VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
209 }
210
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000211 // test innerVector()
212 {
213 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
Gael Guennebaud178858f2009-01-19 15:20:45 +0000214 SparseMatrixType m2(rows, rows);
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000215 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200216 int j0 = internal::random<int>(0,rows-1);
217 int j1 = internal::random<int>(0,rows-1);
Gael Guennebaud2d534662009-01-14 21:27:54 +0000218 VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.col(j0));
219 VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1));
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000220 //m2.innerVector(j0) = 2*m2.innerVector(j1);
221 //refMat2.col(j0) = 2*refMat2.col(j1);
222 //VERIFY_IS_APPROX(m2, refMat2);
223 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100224
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000225 // test innerVectors()
226 {
227 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
228 SparseMatrixType m2(rows, rows);
229 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200230 int j0 = internal::random<int>(0,rows-2);
231 int j1 = internal::random<int>(0,rows-2);
Gael Guennebaud42e25782011-08-19 14:18:05 +0200232 int n0 = internal::random<int>(1,rows-(std::max)(j0,j1));
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000233 VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0));
234 VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
235 refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
236 //m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0);
237 //refMat2.block(0,j0,rows,n0) = refMat2.block(0,j0,rows,n0) + refMat2.block(0,j1,rows,n0);
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000238 }
239
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000240 // test prune
241 {
242 SparseMatrixType m2(rows, rows);
243 DenseMatrix refM2(rows, rows);
244 refM2.setZero();
245 int countFalseNonZero = 0;
246 int countTrueNonZero = 0;
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000247 for (int j=0; j<m2.outerSize(); ++j)
Gael Guennebaud28293142009-05-04 14:25:12 +0000248 {
249 m2.startVec(j);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000250 for (int i=0; i<m2.innerSize(); ++i)
251 {
Benoit Jacob47160402010-10-25 10:15:22 -0400252 float x = internal::random<float>(0,1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000253 if (x<0.1)
254 {
255 // do nothing
256 }
257 else if (x<0.5)
258 {
259 countFalseNonZero++;
Gael Guennebaud8710bd22010-06-02 13:32:13 +0200260 m2.insertBackByOuterInner(j,i) = Scalar(0);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000261 }
262 else
263 {
264 countTrueNonZero++;
Gael Guennebaud8710bd22010-06-02 13:32:13 +0200265 m2.insertBackByOuterInner(j,i) = refM2(i,j) = Scalar(1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000266 }
267 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000268 }
269 m2.finalize();
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000270 VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros());
271 VERIFY_IS_APPROX(m2, refM2);
Hauke Heibeld204ec42010-11-02 14:33:33 +0100272 m2.prune(Scalar(1));
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000273 VERIFY(countTrueNonZero==m2.nonZeros());
274 VERIFY_IS_APPROX(m2, refM2);
275 }
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200276
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100277 // test selfadjointView
278 {
279 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
280 SparseMatrixType m2(rows, rows), m3(rows, rows);
281 initSparse<Scalar>(density, refMat2, m2);
282 refMat3 = refMat2.template selfadjointView<Lower>();
283 m3 = m2.template selfadjointView<Lower>();
284 VERIFY_IS_APPROX(m3, refMat3);
285 }
286
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200287 // test sparseView
288 {
289 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
290 SparseMatrixType m2(rows, rows);
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100291 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200292 VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval());
293 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000294}
295
296void test_sparse_basic()
297{
298 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200299 int s = Eigen::internal::random<int>(1,50);
300 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) ));
301 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double> >(s, s)) ));
302 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(s, s)) ));
303 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,long int>(s, s)) ));
Gael Guennebaud9f795582009-12-16 19:18:40 +0100304
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200305 CALL_SUBTEST_3(( sparse_basic(DynamicSparseMatrix<double>(s, s)) ));
306 CALL_SUBTEST_3(( sparse_basic(DynamicSparseMatrix<double,ColMajor,long int>(s, s)) ));
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000307 }
308}