blob: 6954ba7f97a43247a777f8ef2eac019d74b46830 [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>
Desire NUENTSA4cd82452013-06-11 14:42:29 +02006// Copyright (C) 2013 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
Gael Guennebaud86ccd992008-11-05 13:47:55 +00007//
Benoit Jacob69124cf2012-07-13 14:42:47 -04008// This Source Code Form is subject to the terms of the Mozilla
9// Public License v. 2.0. If a copy of the MPL was not distributed
10// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
Gael Guennebaud86ccd992008-11-05 13:47:55 +000011
12#include "sparse.h"
13
Gael Guennebaud178858f2009-01-19 15:20:45 +000014template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& ref)
15{
Hauke Heibelf1679c72010-06-20 17:37:56 +020016 typedef typename SparseMatrixType::Index Index;
17
18 const Index rows = ref.rows();
19 const Index cols = ref.cols();
Gael Guennebaud178858f2009-01-19 15:20:45 +000020 typedef typename SparseMatrixType::Scalar Scalar;
21 enum { Flags = SparseMatrixType::Flags };
Gael Guennebaud9f795582009-12-16 19:18:40 +010022
Gael Guennebaud42e25782011-08-19 14:18:05 +020023 double density = (std::max)(8./(rows*cols), 0.01);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000024 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
25 typedef Matrix<Scalar,Dynamic,1> DenseVector;
26 Scalar eps = 1e-6;
27
Benoit Jacob47160402010-10-25 10:15:22 -040028 Scalar s1 = internal::random<Scalar>();
Gael Guennebaud86ccd992008-11-05 13:47:55 +000029 {
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020030 SparseMatrixType m(rows, cols);
31 DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
32 DenseVector vec1 = DenseVector::Random(rows);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000033
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020034 std::vector<Vector2i> zeroCoords;
35 std::vector<Vector2i> nonzeroCoords;
36 initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000037
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020038 if (zeroCoords.size()==0 || nonzeroCoords.size()==0)
39 return;
Gael Guennebaud86ccd992008-11-05 13:47:55 +000040
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020041 // test coeff and coeffRef
42 for (int i=0; i<(int)zeroCoords.size(); ++i)
Gael Guennebaud86ccd992008-11-05 13:47:55 +000043 {
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020044 VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps );
45 if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value)
46 VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[0].x(),zeroCoords[0].y()) = 5 );
Gael Guennebaud86ccd992008-11-05 13:47:55 +000047 }
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020048 VERIFY_IS_APPROX(m, refMat);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000049
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020050 m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
51 refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000052
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020053 VERIFY_IS_APPROX(m, refMat);
54 /*
55 // test InnerIterators and Block expressions
56 for (int t=0; t<10; ++t)
57 {
58 int j = internal::random<int>(0,cols-1);
59 int i = internal::random<int>(0,rows-1);
60 int w = internal::random<int>(1,cols-j-1);
61 int h = internal::random<int>(1,rows-i-1);
62
63 // VERIFY_IS_APPROX(m.block(i,j,h,w), refMat.block(i,j,h,w));
64 for(int c=0; c<w; c++)
65 {
66 VERIFY_IS_APPROX(m.block(i,j,h,w).col(c), refMat.block(i,j,h,w).col(c));
67 for(int r=0; r<h; r++)
68 {
69 // VERIFY_IS_APPROX(m.block(i,j,h,w).col(c).coeff(r), refMat.block(i,j,h,w).col(c).coeff(r));
70 }
71 }
72 // for(int r=0; r<h; r++)
73 // {
74 // VERIFY_IS_APPROX(m.block(i,j,h,w).row(r), refMat.block(i,j,h,w).row(r));
75 // for(int c=0; c<w; c++)
76 // {
77 // VERIFY_IS_APPROX(m.block(i,j,h,w).row(r).coeff(c), refMat.block(i,j,h,w).row(r).coeff(c));
78 // }
79 // }
80 }
81
82 for(int c=0; c<cols; c++)
83 {
84 VERIFY_IS_APPROX(m.col(c) + m.col(c), (m + m).col(c));
85 VERIFY_IS_APPROX(m.col(c) + m.col(c), refMat.col(c) + refMat.col(c));
86 }
87
88 for(int r=0; r<rows; r++)
89 {
90 VERIFY_IS_APPROX(m.row(r) + m.row(r), (m + m).row(r));
91 VERIFY_IS_APPROX(m.row(r) + m.row(r), refMat.row(r) + refMat.row(r));
92 }
93 */
94 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +000095
Gael Guennebaud28293142009-05-04 14:25:12 +000096 // test insert (inner random)
Gael Guennebaud5015e482008-12-11 18:26:24 +000097 {
98 DenseMatrix m1(rows,cols);
99 m1.setZero();
Gael Guennebaud178858f2009-01-19 15:20:45 +0000100 SparseMatrixType m2(rows,cols);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100101 if(internal::random<int>()%2)
102 m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
Gael Guennebaud5015e482008-12-11 18:26:24 +0000103 for (int j=0; j<cols; ++j)
104 {
105 for (int k=0; k<rows/2; ++k)
106 {
Benoit Jacob47160402010-10-25 10:15:22 -0400107 int i = internal::random<int>(0,rows-1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000108 if (m1.coeff(i,j)==Scalar(0))
Benoit Jacob47160402010-10-25 10:15:22 -0400109 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud5015e482008-12-11 18:26:24 +0000110 }
111 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000112 m2.finalize();
113 VERIFY_IS_APPROX(m2,m1);
114 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100115
Gael Guennebaud28293142009-05-04 14:25:12 +0000116 // test insert (fully random)
117 {
118 DenseMatrix m1(rows,cols);
119 m1.setZero();
120 SparseMatrixType m2(rows,cols);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100121 if(internal::random<int>()%2)
122 m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
Gael Guennebaud28293142009-05-04 14:25:12 +0000123 for (int k=0; k<rows*cols; ++k)
124 {
Benoit Jacob47160402010-10-25 10:15:22 -0400125 int i = internal::random<int>(0,rows-1);
126 int j = internal::random<int>(0,cols-1);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100127 if ((m1.coeff(i,j)==Scalar(0)) && (internal::random<int>()%2))
Benoit Jacob47160402010-10-25 10:15:22 -0400128 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100129 else
130 {
131 Scalar v = internal::random<Scalar>();
132 m2.coeffRef(i,j) += v;
133 m1(i,j) += v;
134 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000135 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000136 VERIFY_IS_APPROX(m2,m1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000137 }
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200138
139 // test insert (un-compressed)
140 for(int mode=0;mode<4;++mode)
141 {
142 DenseMatrix m1(rows,cols);
143 m1.setZero();
144 SparseMatrixType m2(rows,cols);
145 VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? m2.innerSize() : std::max<int>(1,m2.innerSize()/8)));
146 m2.reserve(r);
147 for (int k=0; k<rows*cols; ++k)
148 {
149 int i = internal::random<int>(0,rows-1);
150 int j = internal::random<int>(0,cols-1);
151 if (m1.coeff(i,j)==Scalar(0))
152 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
153 if(mode==3)
154 m2.reserve(r);
155 }
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100156 if(internal::random<int>()%2)
157 m2.makeCompressed();
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200158 VERIFY_IS_APPROX(m2,m1);
159 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100160
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000161 // test innerVector()
162 {
163 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
Gael Guennebaud178858f2009-01-19 15:20:45 +0000164 SparseMatrixType m2(rows, rows);
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000165 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200166 int j0 = internal::random<int>(0,rows-1);
167 int j1 = internal::random<int>(0,rows-1);
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100168 if(SparseMatrixType::IsRowMajor)
169 VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.row(j0));
170 else
171 VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.col(j0));
172
173 if(SparseMatrixType::IsRowMajor)
174 VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.row(j0)+refMat2.row(j1));
175 else
176 VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100177
178 SparseMatrixType m3(rows,rows);
179 m3.reserve(VectorXi::Constant(rows,rows/2));
180 for(int j=0; j<rows; ++j)
181 for(int k=0; k<j; ++k)
182 m3.insertByOuterInner(j,k) = k+1;
183 for(int j=0; j<rows; ++j)
184 {
Gael Guennebaud65c59302013-06-12 10:37:15 +0200185 VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100186 if(j>0)
Gael Guennebaud65c59302013-06-12 10:37:15 +0200187 VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100188 }
189 m3.makeCompressed();
190 for(int j=0; j<rows; ++j)
191 {
Gael Guennebaud65c59302013-06-12 10:37:15 +0200192 VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100193 if(j>0)
Gael Guennebaud65c59302013-06-12 10:37:15 +0200194 VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100195 }
196
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000197 //m2.innerVector(j0) = 2*m2.innerVector(j1);
198 //refMat2.col(j0) = 2*refMat2.col(j1);
199 //VERIFY_IS_APPROX(m2, refMat2);
200 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100201
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000202 // test innerVectors()
203 {
204 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
205 SparseMatrixType m2(rows, rows);
206 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud7450b232013-04-12 13:20:13 +0200207 if(internal::random<float>(0,1)>0.5) m2.makeCompressed();
208
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200209 int j0 = internal::random<int>(0,rows-2);
210 int j1 = internal::random<int>(0,rows-2);
Gael Guennebaud42e25782011-08-19 14:18:05 +0200211 int n0 = internal::random<int>(1,rows-(std::max)(j0,j1));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100212 if(SparseMatrixType::IsRowMajor)
213 VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(j0,0,n0,cols));
214 else
215 VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0));
216 if(SparseMatrixType::IsRowMajor)
217 VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
Gael Guennebaud7450b232013-04-12 13:20:13 +0200218 refMat2.middleRows(j0,n0)+refMat2.middleRows(j1,n0));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100219 else
220 VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
221 refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
Gael Guennebaud7450b232013-04-12 13:20:13 +0200222
223 VERIFY_IS_APPROX(m2, refMat2);
224
225 m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0);
226 if(SparseMatrixType::IsRowMajor)
227 refMat2.middleRows(j0,n0) = (refMat2.middleRows(j0,n0) + refMat2.middleRows(j1,n0)).eval();
228 else
229 refMat2.middleCols(j0,n0) = (refMat2.middleCols(j0,n0) + refMat2.middleCols(j1,n0)).eval();
230
231 VERIFY_IS_APPROX(m2, refMat2);
232
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000233 }
Gael Guennebaud4e602832012-11-16 09:02:50 +0100234
235 // test basic computations
236 {
237 DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
238 DenseMatrix refM2 = DenseMatrix::Zero(rows, rows);
239 DenseMatrix refM3 = DenseMatrix::Zero(rows, rows);
240 DenseMatrix refM4 = DenseMatrix::Zero(rows, rows);
241 SparseMatrixType m1(rows, rows);
242 SparseMatrixType m2(rows, rows);
243 SparseMatrixType m3(rows, rows);
244 SparseMatrixType m4(rows, rows);
245 initSparse<Scalar>(density, refM1, m1);
246 initSparse<Scalar>(density, refM2, m2);
247 initSparse<Scalar>(density, refM3, m3);
248 initSparse<Scalar>(density, refM4, m4);
249
250 VERIFY_IS_APPROX(m1+m2, refM1+refM2);
251 VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3);
252 VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2));
253 VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2);
254
255 VERIFY_IS_APPROX(m1*=s1, refM1*=s1);
256 VERIFY_IS_APPROX(m1/=s1, refM1/=s1);
257
258 VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
259 VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
260
261 if(SparseMatrixType::IsRowMajor)
262 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0)));
263 else
264 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.col(0).dot(refM2.row(0)));
265
266 VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate());
267 VERIFY_IS_APPROX(m1.real(), refM1.real());
268
269 refM4.setRandom();
270 // sparse cwise* dense
271 VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4));
272// VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4);
273
274 // test aliasing
275 VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1));
276 VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval()));
277 VERIFY_IS_APPROX((m1 = -m1.transpose()), (refM1 = -refM1.transpose().eval()));
278 VERIFY_IS_APPROX((m1 += -m1), (refM1 += -refM1));
279 }
280
281 // test transpose
282 {
283 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
284 SparseMatrixType m2(rows, rows);
285 initSparse<Scalar>(density, refMat2, m2);
286 VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval());
287 VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
288
289 VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
290 }
291
292
293
294 // test generic blocks
295 {
296 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
297 SparseMatrixType m2(rows, rows);
298 initSparse<Scalar>(density, refMat2, m2);
299 int j0 = internal::random<int>(0,rows-2);
300 int j1 = internal::random<int>(0,rows-2);
301 int n0 = internal::random<int>(1,rows-(std::max)(j0,j1));
302 if(SparseMatrixType::IsRowMajor)
303 VERIFY_IS_APPROX(m2.block(j0,0,n0,cols), refMat2.block(j0,0,n0,cols));
304 else
305 VERIFY_IS_APPROX(m2.block(0,j0,rows,n0), refMat2.block(0,j0,rows,n0));
306
307 if(SparseMatrixType::IsRowMajor)
308 VERIFY_IS_APPROX(m2.block(j0,0,n0,cols)+m2.block(j1,0,n0,cols),
309 refMat2.block(j0,0,n0,cols)+refMat2.block(j1,0,n0,cols));
310 else
311 VERIFY_IS_APPROX(m2.block(0,j0,rows,n0)+m2.block(0,j1,rows,n0),
312 refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
Gael Guennebaud5f178e52013-06-14 10:52:19 +0200313
314 int i = internal::random<int>(0,m2.outerSize()-1);
315 if(SparseMatrixType::IsRowMajor) {
316 m2.innerVector(i) = m2.innerVector(i) * s1;
317 refMat2.row(i) = refMat2.row(i) * s1;
318 VERIFY_IS_APPROX(m2,refMat2);
319 } else {
320 m2.innerVector(i) = m2.innerVector(i) * s1;
321 refMat2.col(i) = refMat2.col(i) * s1;
322 VERIFY_IS_APPROX(m2,refMat2);
323 }
Gael Guennebaud4e602832012-11-16 09:02:50 +0100324 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000325
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000326 // test prune
327 {
328 SparseMatrixType m2(rows, rows);
329 DenseMatrix refM2(rows, rows);
330 refM2.setZero();
331 int countFalseNonZero = 0;
332 int countTrueNonZero = 0;
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000333 for (int j=0; j<m2.outerSize(); ++j)
Gael Guennebaud28293142009-05-04 14:25:12 +0000334 {
335 m2.startVec(j);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000336 for (int i=0; i<m2.innerSize(); ++i)
337 {
Benoit Jacob47160402010-10-25 10:15:22 -0400338 float x = internal::random<float>(0,1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000339 if (x<0.1)
340 {
341 // do nothing
342 }
343 else if (x<0.5)
344 {
345 countFalseNonZero++;
Gael Guennebaud8710bd22010-06-02 13:32:13 +0200346 m2.insertBackByOuterInner(j,i) = Scalar(0);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000347 }
348 else
349 {
350 countTrueNonZero++;
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100351 m2.insertBackByOuterInner(j,i) = Scalar(1);
352 if(SparseMatrixType::IsRowMajor)
353 refM2(j,i) = Scalar(1);
354 else
355 refM2(i,j) = Scalar(1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000356 }
357 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000358 }
359 m2.finalize();
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000360 VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros());
361 VERIFY_IS_APPROX(m2, refM2);
Hauke Heibeld204ec42010-11-02 14:33:33 +0100362 m2.prune(Scalar(1));
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000363 VERIFY(countTrueNonZero==m2.nonZeros());
364 VERIFY_IS_APPROX(m2, refM2);
365 }
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100366
Gael Guennebaud87138072012-01-28 11:13:59 +0100367 // test setFromTriplets
368 {
369 typedef Triplet<Scalar,Index> TripletType;
370 std::vector<TripletType> triplets;
371 int ntriplets = rows*cols;
372 triplets.reserve(ntriplets);
373 DenseMatrix refMat(rows,cols);
374 refMat.setZero();
375 for(int i=0;i<ntriplets;++i)
376 {
Gael Guennebaud18e3ac02012-01-31 09:14:01 +0100377 int r = internal::random<int>(0,rows-1);
378 int c = internal::random<int>(0,cols-1);
Gael Guennebaud87138072012-01-28 11:13:59 +0100379 Scalar v = internal::random<Scalar>();
Gael Guennebaud18e3ac02012-01-31 09:14:01 +0100380 triplets.push_back(TripletType(r,c,v));
381 refMat(r,c) += v;
Gael Guennebaud87138072012-01-28 11:13:59 +0100382 }
383 SparseMatrixType m(rows,cols);
384 m.setFromTriplets(triplets.begin(), triplets.end());
385 VERIFY_IS_APPROX(m, refMat);
386 }
387
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100388 // test triangularView
389 {
390 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
391 SparseMatrixType m2(rows, rows), m3(rows, rows);
392 initSparse<Scalar>(density, refMat2, m2);
393 refMat3 = refMat2.template triangularView<Lower>();
394 m3 = m2.template triangularView<Lower>();
395 VERIFY_IS_APPROX(m3, refMat3);
396
397 refMat3 = refMat2.template triangularView<Upper>();
398 m3 = m2.template triangularView<Upper>();
399 VERIFY_IS_APPROX(m3, refMat3);
400
401 refMat3 = refMat2.template triangularView<UnitUpper>();
402 m3 = m2.template triangularView<UnitUpper>();
403 VERIFY_IS_APPROX(m3, refMat3);
404
405 refMat3 = refMat2.template triangularView<UnitLower>();
406 m3 = m2.template triangularView<UnitLower>();
407 VERIFY_IS_APPROX(m3, refMat3);
Desire NUENTSA4cd82452013-06-11 14:42:29 +0200408
409 refMat3 = refMat2.template triangularView<StrictlyUpper>();
410 m3 = m2.template triangularView<StrictlyUpper>();
411 VERIFY_IS_APPROX(m3, refMat3);
412
413 refMat3 = refMat2.template triangularView<StrictlyLower>();
414 m3 = m2.template triangularView<StrictlyLower>();
415 VERIFY_IS_APPROX(m3, refMat3);
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100416 }
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200417
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100418 // test selfadjointView
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100419 if(!SparseMatrixType::IsRowMajor)
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100420 {
421 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
422 SparseMatrixType m2(rows, rows), m3(rows, rows);
423 initSparse<Scalar>(density, refMat2, m2);
424 refMat3 = refMat2.template selfadjointView<Lower>();
425 m3 = m2.template selfadjointView<Lower>();
426 VERIFY_IS_APPROX(m3, refMat3);
427 }
428
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200429 // test sparseView
430 {
431 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
432 SparseMatrixType m2(rows, rows);
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100433 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200434 VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval());
435 }
Gael Guennebaud82f9aa12011-12-04 21:49:21 +0100436
437 // test diagonal
438 {
439 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
440 SparseMatrixType m2(rows, rows);
441 initSparse<Scalar>(density, refMat2, m2);
442 VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval());
443 }
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200444
445 // test conservative resize
446 {
447 std::vector< std::pair<int,int> > inc;
448 inc.push_back(std::pair<int,int>(-3,-2));
449 inc.push_back(std::pair<int,int>(0,0));
450 inc.push_back(std::pair<int,int>(3,2));
451 inc.push_back(std::pair<int,int>(3,0));
452 inc.push_back(std::pair<int,int>(0,3));
453
454 for(size_t i = 0; i< inc.size(); i++) {
455 int incRows = inc[i].first;
456 int incCols = inc[i].second;
457 SparseMatrixType m1(rows, cols);
458 DenseMatrix refMat1 = DenseMatrix::Zero(rows, cols);
459 initSparse<Scalar>(density, refMat1, m1);
460
461 m1.conservativeResize(rows+incRows, cols+incCols);
462 refMat1.conservativeResize(rows+incRows, cols+incCols);
463 if (incRows > 0) refMat1.bottomRows(incRows).setZero();
464 if (incCols > 0) refMat1.rightCols(incCols).setZero();
465
466 VERIFY_IS_APPROX(m1, refMat1);
467
468 // Insert new values
469 if (incRows > 0)
470 m1.insert(refMat1.rows()-1, 0) = refMat1(refMat1.rows()-1, 0) = 1;
471 if (incCols > 0)
472 m1.insert(0, refMat1.cols()-1) = refMat1(0, refMat1.cols()-1) = 1;
473
474 VERIFY_IS_APPROX(m1, refMat1);
475
476
477 }
478 }
Desire NUENTSA4cd82452013-06-11 14:42:29 +0200479
480 // test Identity matrix
481 {
482 DenseMatrix refMat1 = DenseMatrix::Identity(rows, rows);
483 SparseMatrixType m1(rows, rows);
484 m1.setIdentity();
485 VERIFY_IS_APPROX(m1, refMat1);
486 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000487}
488
489void test_sparse_basic()
490{
491 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200492 int s = Eigen::internal::random<int>(1,50);
Gael Guennebauda1091ca2013-02-15 14:05:37 +0100493 EIGEN_UNUSED_VARIABLE(s);
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200494 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) ));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100495 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(s, s)) ));
496 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(s, s)) ));
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200497 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(s, s)) ));
498 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,long int>(s, s)) ));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100499 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,long int>(s, s)) ));
Gael Guennebaud47e89022013-06-10 10:34:03 +0200500
501 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,short int>(s, s)) ));
502 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,short int>(s, s)) ));
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000503 }
504}