blob: 8fc1904b1dc320d3446d1e687a2e3521ed56c87c [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 */
Gael Guennebauda915f022013-06-28 16:16:02 +020094
95 // test assertion
96 VERIFY_RAISES_ASSERT( m.coeffRef(-1,1) = 0 );
97 VERIFY_RAISES_ASSERT( m.coeffRef(0,m.cols()) = 0 );
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020098 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +000099
Gael Guennebaud28293142009-05-04 14:25:12 +0000100 // test insert (inner random)
Gael Guennebaud5015e482008-12-11 18:26:24 +0000101 {
102 DenseMatrix m1(rows,cols);
103 m1.setZero();
Gael Guennebaud178858f2009-01-19 15:20:45 +0000104 SparseMatrixType m2(rows,cols);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100105 if(internal::random<int>()%2)
106 m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
Gael Guennebaud5015e482008-12-11 18:26:24 +0000107 for (int j=0; j<cols; ++j)
108 {
109 for (int k=0; k<rows/2; ++k)
110 {
Benoit Jacob47160402010-10-25 10:15:22 -0400111 int i = internal::random<int>(0,rows-1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000112 if (m1.coeff(i,j)==Scalar(0))
Benoit Jacob47160402010-10-25 10:15:22 -0400113 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud5015e482008-12-11 18:26:24 +0000114 }
115 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000116 m2.finalize();
117 VERIFY_IS_APPROX(m2,m1);
118 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100119
Gael Guennebaud28293142009-05-04 14:25:12 +0000120 // test insert (fully random)
121 {
122 DenseMatrix m1(rows,cols);
123 m1.setZero();
124 SparseMatrixType m2(rows,cols);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100125 if(internal::random<int>()%2)
126 m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
Gael Guennebaud28293142009-05-04 14:25:12 +0000127 for (int k=0; k<rows*cols; ++k)
128 {
Benoit Jacob47160402010-10-25 10:15:22 -0400129 int i = internal::random<int>(0,rows-1);
130 int j = internal::random<int>(0,cols-1);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100131 if ((m1.coeff(i,j)==Scalar(0)) && (internal::random<int>()%2))
Benoit Jacob47160402010-10-25 10:15:22 -0400132 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100133 else
134 {
135 Scalar v = internal::random<Scalar>();
136 m2.coeffRef(i,j) += v;
137 m1(i,j) += v;
138 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000139 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000140 VERIFY_IS_APPROX(m2,m1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000141 }
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200142
143 // test insert (un-compressed)
144 for(int mode=0;mode<4;++mode)
145 {
146 DenseMatrix m1(rows,cols);
147 m1.setZero();
148 SparseMatrixType m2(rows,cols);
149 VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? m2.innerSize() : std::max<int>(1,m2.innerSize()/8)));
150 m2.reserve(r);
151 for (int k=0; k<rows*cols; ++k)
152 {
153 int i = internal::random<int>(0,rows-1);
154 int j = internal::random<int>(0,cols-1);
155 if (m1.coeff(i,j)==Scalar(0))
156 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
157 if(mode==3)
158 m2.reserve(r);
159 }
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100160 if(internal::random<int>()%2)
161 m2.makeCompressed();
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200162 VERIFY_IS_APPROX(m2,m1);
163 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100164
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000165 // test innerVector()
166 {
167 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
Gael Guennebaud178858f2009-01-19 15:20:45 +0000168 SparseMatrixType m2(rows, rows);
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000169 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200170 int j0 = internal::random<int>(0,rows-1);
171 int j1 = internal::random<int>(0,rows-1);
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100172 if(SparseMatrixType::IsRowMajor)
173 VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.row(j0));
174 else
175 VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.col(j0));
176
177 if(SparseMatrixType::IsRowMajor)
178 VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.row(j0)+refMat2.row(j1));
179 else
180 VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100181
182 SparseMatrixType m3(rows,rows);
183 m3.reserve(VectorXi::Constant(rows,rows/2));
184 for(int j=0; j<rows; ++j)
185 for(int k=0; k<j; ++k)
186 m3.insertByOuterInner(j,k) = k+1;
187 for(int j=0; j<rows; ++j)
188 {
Gael Guennebaud65c59302013-06-12 10:37:15 +0200189 VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100190 if(j>0)
Gael Guennebaud65c59302013-06-12 10:37:15 +0200191 VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100192 }
193 m3.makeCompressed();
194 for(int j=0; j<rows; ++j)
195 {
Gael Guennebaud65c59302013-06-12 10:37:15 +0200196 VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100197 if(j>0)
Gael Guennebaud65c59302013-06-12 10:37:15 +0200198 VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100199 }
200
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000201 //m2.innerVector(j0) = 2*m2.innerVector(j1);
202 //refMat2.col(j0) = 2*refMat2.col(j1);
203 //VERIFY_IS_APPROX(m2, refMat2);
204 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100205
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000206 // test innerVectors()
207 {
208 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
209 SparseMatrixType m2(rows, rows);
210 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud7450b232013-04-12 13:20:13 +0200211 if(internal::random<float>(0,1)>0.5) m2.makeCompressed();
212
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200213 int j0 = internal::random<int>(0,rows-2);
214 int j1 = internal::random<int>(0,rows-2);
Gael Guennebaud42e25782011-08-19 14:18:05 +0200215 int n0 = internal::random<int>(1,rows-(std::max)(j0,j1));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100216 if(SparseMatrixType::IsRowMajor)
217 VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(j0,0,n0,cols));
218 else
219 VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0));
220 if(SparseMatrixType::IsRowMajor)
221 VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
Gael Guennebaud7450b232013-04-12 13:20:13 +0200222 refMat2.middleRows(j0,n0)+refMat2.middleRows(j1,n0));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100223 else
224 VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
225 refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
Gael Guennebaud7450b232013-04-12 13:20:13 +0200226
227 VERIFY_IS_APPROX(m2, refMat2);
228
229 m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0);
230 if(SparseMatrixType::IsRowMajor)
231 refMat2.middleRows(j0,n0) = (refMat2.middleRows(j0,n0) + refMat2.middleRows(j1,n0)).eval();
232 else
233 refMat2.middleCols(j0,n0) = (refMat2.middleCols(j0,n0) + refMat2.middleCols(j1,n0)).eval();
234
235 VERIFY_IS_APPROX(m2, refMat2);
236
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000237 }
Gael Guennebaud4e602832012-11-16 09:02:50 +0100238
239 // test basic computations
240 {
241 DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
242 DenseMatrix refM2 = DenseMatrix::Zero(rows, rows);
243 DenseMatrix refM3 = DenseMatrix::Zero(rows, rows);
244 DenseMatrix refM4 = DenseMatrix::Zero(rows, rows);
245 SparseMatrixType m1(rows, rows);
246 SparseMatrixType m2(rows, rows);
247 SparseMatrixType m3(rows, rows);
248 SparseMatrixType m4(rows, rows);
249 initSparse<Scalar>(density, refM1, m1);
250 initSparse<Scalar>(density, refM2, m2);
251 initSparse<Scalar>(density, refM3, m3);
252 initSparse<Scalar>(density, refM4, m4);
253
254 VERIFY_IS_APPROX(m1+m2, refM1+refM2);
255 VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3);
256 VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2));
257 VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2);
258
259 VERIFY_IS_APPROX(m1*=s1, refM1*=s1);
260 VERIFY_IS_APPROX(m1/=s1, refM1/=s1);
261
262 VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
263 VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
264
265 if(SparseMatrixType::IsRowMajor)
266 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0)));
267 else
268 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.col(0).dot(refM2.row(0)));
269
270 VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate());
271 VERIFY_IS_APPROX(m1.real(), refM1.real());
272
273 refM4.setRandom();
274 // sparse cwise* dense
275 VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4));
276// VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4);
277
278 // test aliasing
279 VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1));
280 VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval()));
281 VERIFY_IS_APPROX((m1 = -m1.transpose()), (refM1 = -refM1.transpose().eval()));
282 VERIFY_IS_APPROX((m1 += -m1), (refM1 += -refM1));
283 }
284
285 // test transpose
286 {
287 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
288 SparseMatrixType m2(rows, rows);
289 initSparse<Scalar>(density, refMat2, m2);
290 VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval());
291 VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
292
293 VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
294 }
295
296
297
298 // test generic blocks
299 {
300 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
301 SparseMatrixType m2(rows, rows);
302 initSparse<Scalar>(density, refMat2, m2);
303 int j0 = internal::random<int>(0,rows-2);
304 int j1 = internal::random<int>(0,rows-2);
305 int n0 = internal::random<int>(1,rows-(std::max)(j0,j1));
306 if(SparseMatrixType::IsRowMajor)
307 VERIFY_IS_APPROX(m2.block(j0,0,n0,cols), refMat2.block(j0,0,n0,cols));
308 else
309 VERIFY_IS_APPROX(m2.block(0,j0,rows,n0), refMat2.block(0,j0,rows,n0));
310
311 if(SparseMatrixType::IsRowMajor)
312 VERIFY_IS_APPROX(m2.block(j0,0,n0,cols)+m2.block(j1,0,n0,cols),
313 refMat2.block(j0,0,n0,cols)+refMat2.block(j1,0,n0,cols));
314 else
315 VERIFY_IS_APPROX(m2.block(0,j0,rows,n0)+m2.block(0,j1,rows,n0),
316 refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
Gael Guennebaud5f178e52013-06-14 10:52:19 +0200317
318 int i = internal::random<int>(0,m2.outerSize()-1);
319 if(SparseMatrixType::IsRowMajor) {
320 m2.innerVector(i) = m2.innerVector(i) * s1;
321 refMat2.row(i) = refMat2.row(i) * s1;
322 VERIFY_IS_APPROX(m2,refMat2);
323 } else {
324 m2.innerVector(i) = m2.innerVector(i) * s1;
325 refMat2.col(i) = refMat2.col(i) * s1;
326 VERIFY_IS_APPROX(m2,refMat2);
327 }
Gael Guennebaud4e602832012-11-16 09:02:50 +0100328 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000329
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000330 // test prune
331 {
332 SparseMatrixType m2(rows, rows);
333 DenseMatrix refM2(rows, rows);
334 refM2.setZero();
335 int countFalseNonZero = 0;
336 int countTrueNonZero = 0;
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000337 for (int j=0; j<m2.outerSize(); ++j)
Gael Guennebaud28293142009-05-04 14:25:12 +0000338 {
339 m2.startVec(j);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000340 for (int i=0; i<m2.innerSize(); ++i)
341 {
Benoit Jacob47160402010-10-25 10:15:22 -0400342 float x = internal::random<float>(0,1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000343 if (x<0.1)
344 {
345 // do nothing
346 }
347 else if (x<0.5)
348 {
349 countFalseNonZero++;
Gael Guennebaud8710bd22010-06-02 13:32:13 +0200350 m2.insertBackByOuterInner(j,i) = Scalar(0);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000351 }
352 else
353 {
354 countTrueNonZero++;
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100355 m2.insertBackByOuterInner(j,i) = Scalar(1);
356 if(SparseMatrixType::IsRowMajor)
357 refM2(j,i) = Scalar(1);
358 else
359 refM2(i,j) = Scalar(1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000360 }
361 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000362 }
363 m2.finalize();
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000364 VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros());
365 VERIFY_IS_APPROX(m2, refM2);
Hauke Heibeld204ec42010-11-02 14:33:33 +0100366 m2.prune(Scalar(1));
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000367 VERIFY(countTrueNonZero==m2.nonZeros());
368 VERIFY_IS_APPROX(m2, refM2);
369 }
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100370
Gael Guennebaud87138072012-01-28 11:13:59 +0100371 // test setFromTriplets
372 {
373 typedef Triplet<Scalar,Index> TripletType;
374 std::vector<TripletType> triplets;
375 int ntriplets = rows*cols;
376 triplets.reserve(ntriplets);
377 DenseMatrix refMat(rows,cols);
378 refMat.setZero();
379 for(int i=0;i<ntriplets;++i)
380 {
Gael Guennebaud18e3ac02012-01-31 09:14:01 +0100381 int r = internal::random<int>(0,rows-1);
382 int c = internal::random<int>(0,cols-1);
Gael Guennebaud87138072012-01-28 11:13:59 +0100383 Scalar v = internal::random<Scalar>();
Gael Guennebaud18e3ac02012-01-31 09:14:01 +0100384 triplets.push_back(TripletType(r,c,v));
385 refMat(r,c) += v;
Gael Guennebaud87138072012-01-28 11:13:59 +0100386 }
387 SparseMatrixType m(rows,cols);
388 m.setFromTriplets(triplets.begin(), triplets.end());
389 VERIFY_IS_APPROX(m, refMat);
390 }
391
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100392 // test triangularView
393 {
394 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
395 SparseMatrixType m2(rows, rows), m3(rows, rows);
396 initSparse<Scalar>(density, refMat2, m2);
397 refMat3 = refMat2.template triangularView<Lower>();
398 m3 = m2.template triangularView<Lower>();
399 VERIFY_IS_APPROX(m3, refMat3);
400
401 refMat3 = refMat2.template triangularView<Upper>();
402 m3 = m2.template triangularView<Upper>();
403 VERIFY_IS_APPROX(m3, refMat3);
404
405 refMat3 = refMat2.template triangularView<UnitUpper>();
406 m3 = m2.template triangularView<UnitUpper>();
407 VERIFY_IS_APPROX(m3, refMat3);
408
409 refMat3 = refMat2.template triangularView<UnitLower>();
410 m3 = m2.template triangularView<UnitLower>();
411 VERIFY_IS_APPROX(m3, refMat3);
Desire NUENTSA4cd82452013-06-11 14:42:29 +0200412
413 refMat3 = refMat2.template triangularView<StrictlyUpper>();
414 m3 = m2.template triangularView<StrictlyUpper>();
415 VERIFY_IS_APPROX(m3, refMat3);
416
417 refMat3 = refMat2.template triangularView<StrictlyLower>();
418 m3 = m2.template triangularView<StrictlyLower>();
419 VERIFY_IS_APPROX(m3, refMat3);
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100420 }
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200421
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100422 // test selfadjointView
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100423 if(!SparseMatrixType::IsRowMajor)
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100424 {
425 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
426 SparseMatrixType m2(rows, rows), m3(rows, rows);
427 initSparse<Scalar>(density, refMat2, m2);
428 refMat3 = refMat2.template selfadjointView<Lower>();
429 m3 = m2.template selfadjointView<Lower>();
430 VERIFY_IS_APPROX(m3, refMat3);
431 }
432
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200433 // test sparseView
434 {
435 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
436 SparseMatrixType m2(rows, rows);
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100437 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200438 VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval());
439 }
Gael Guennebaud82f9aa12011-12-04 21:49:21 +0100440
441 // test diagonal
442 {
443 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
444 SparseMatrixType m2(rows, rows);
445 initSparse<Scalar>(density, refMat2, m2);
446 VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval());
447 }
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200448
449 // test conservative resize
450 {
451 std::vector< std::pair<int,int> > inc;
452 inc.push_back(std::pair<int,int>(-3,-2));
453 inc.push_back(std::pair<int,int>(0,0));
454 inc.push_back(std::pair<int,int>(3,2));
455 inc.push_back(std::pair<int,int>(3,0));
456 inc.push_back(std::pair<int,int>(0,3));
457
458 for(size_t i = 0; i< inc.size(); i++) {
459 int incRows = inc[i].first;
460 int incCols = inc[i].second;
461 SparseMatrixType m1(rows, cols);
462 DenseMatrix refMat1 = DenseMatrix::Zero(rows, cols);
463 initSparse<Scalar>(density, refMat1, m1);
464
465 m1.conservativeResize(rows+incRows, cols+incCols);
466 refMat1.conservativeResize(rows+incRows, cols+incCols);
467 if (incRows > 0) refMat1.bottomRows(incRows).setZero();
468 if (incCols > 0) refMat1.rightCols(incCols).setZero();
469
470 VERIFY_IS_APPROX(m1, refMat1);
471
472 // Insert new values
473 if (incRows > 0)
474 m1.insert(refMat1.rows()-1, 0) = refMat1(refMat1.rows()-1, 0) = 1;
475 if (incCols > 0)
476 m1.insert(0, refMat1.cols()-1) = refMat1(0, refMat1.cols()-1) = 1;
477
478 VERIFY_IS_APPROX(m1, refMat1);
479
480
481 }
482 }
Desire NUENTSA4cd82452013-06-11 14:42:29 +0200483
484 // test Identity matrix
485 {
486 DenseMatrix refMat1 = DenseMatrix::Identity(rows, rows);
487 SparseMatrixType m1(rows, rows);
488 m1.setIdentity();
489 VERIFY_IS_APPROX(m1, refMat1);
490 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000491}
492
493void test_sparse_basic()
494{
495 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200496 int s = Eigen::internal::random<int>(1,50);
Gael Guennebauda1091ca2013-02-15 14:05:37 +0100497 EIGEN_UNUSED_VARIABLE(s);
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200498 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) ));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100499 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(s, s)) ));
500 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(s, s)) ));
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200501 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(s, s)) ));
502 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,long int>(s, s)) ));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100503 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,long int>(s, s)) ));
Gael Guennebaud47e89022013-06-10 10:34:03 +0200504
505 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,short int>(s, s)) ));
506 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,short int>(s, s)) ));
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000507 }
508}