blob: e300f2537081997cb1ed37959e293800abc8ff58 [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
Gael Guennebaud178858f2009-01-19 15:20:45 +000028 SparseMatrixType m(rows, cols);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000029 DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
30 DenseVector vec1 = DenseVector::Random(rows);
Benoit Jacob47160402010-10-25 10:15:22 -040031 Scalar s1 = internal::random<Scalar>();
Gael Guennebaud86ccd992008-11-05 13:47:55 +000032
33 std::vector<Vector2i> zeroCoords;
34 std::vector<Vector2i> nonzeroCoords;
35 initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords);
Gael Guennebaud9f795582009-12-16 19:18:40 +010036
Gael Guennebaud86ccd992008-11-05 13:47:55 +000037 if (zeroCoords.size()==0 || nonzeroCoords.size()==0)
38 return;
39
40 // test coeff and coeffRef
41 for (int i=0; i<(int)zeroCoords.size(); ++i)
42 {
43 VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps );
Hauke Heibel7bc8e3a2010-10-25 22:13:49 +020044 if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value)
Gael Guennebaud178858f2009-01-19 15:20:45 +000045 VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[0].x(),zeroCoords[0].y()) = 5 );
Gael Guennebaud86ccd992008-11-05 13:47:55 +000046 }
47 VERIFY_IS_APPROX(m, refMat);
48
49 m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
50 refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
51
52 VERIFY_IS_APPROX(m, refMat);
Gael Guennebaudc4c70662009-01-14 14:24:10 +000053 /*
Gael Guennebaud86ccd992008-11-05 13:47:55 +000054 // test InnerIterators and Block expressions
55 for (int t=0; t<10; ++t)
56 {
Benoit Jacob47160402010-10-25 10:15:22 -040057 int j = internal::random<int>(0,cols-1);
58 int i = internal::random<int>(0,rows-1);
59 int w = internal::random<int>(1,cols-j-1);
60 int h = internal::random<int>(1,rows-i-1);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000061
Gael Guennebaudc4c70662009-01-14 14:24:10 +000062// VERIFY_IS_APPROX(m.block(i,j,h,w), refMat.block(i,j,h,w));
Gael Guennebaud86ccd992008-11-05 13:47:55 +000063 for(int c=0; c<w; c++)
64 {
65 VERIFY_IS_APPROX(m.block(i,j,h,w).col(c), refMat.block(i,j,h,w).col(c));
66 for(int r=0; r<h; r++)
67 {
Gael Guennebaudc4c70662009-01-14 14:24:10 +000068// 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 +000069 }
70 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +000071// for(int r=0; r<h; r++)
72// {
73// VERIFY_IS_APPROX(m.block(i,j,h,w).row(r), refMat.block(i,j,h,w).row(r));
74// for(int c=0; c<w; c++)
75// {
76// VERIFY_IS_APPROX(m.block(i,j,h,w).row(r).coeff(c), refMat.block(i,j,h,w).row(r).coeff(c));
77// }
78// }
Gael Guennebaud86ccd992008-11-05 13:47:55 +000079 }
80
81 for(int c=0; c<cols; c++)
82 {
83 VERIFY_IS_APPROX(m.col(c) + m.col(c), (m + m).col(c));
84 VERIFY_IS_APPROX(m.col(c) + m.col(c), refMat.col(c) + refMat.col(c));
85 }
86
87 for(int r=0; r<rows; r++)
88 {
89 VERIFY_IS_APPROX(m.row(r) + m.row(r), (m + m).row(r));
90 VERIFY_IS_APPROX(m.row(r) + m.row(r), refMat.row(r) + refMat.row(r));
91 }
92 */
93
Gael Guennebaud28293142009-05-04 14:25:12 +000094 // test insert (inner random)
Gael Guennebaud5015e482008-12-11 18:26:24 +000095 {
96 DenseMatrix m1(rows,cols);
97 m1.setZero();
Gael Guennebaud178858f2009-01-19 15:20:45 +000098 SparseMatrixType m2(rows,cols);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +010099 if(internal::random<int>()%2)
100 m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
Gael Guennebaud5015e482008-12-11 18:26:24 +0000101 for (int j=0; j<cols; ++j)
102 {
103 for (int k=0; k<rows/2; ++k)
104 {
Benoit Jacob47160402010-10-25 10:15:22 -0400105 int i = internal::random<int>(0,rows-1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000106 if (m1.coeff(i,j)==Scalar(0))
Benoit Jacob47160402010-10-25 10:15:22 -0400107 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud5015e482008-12-11 18:26:24 +0000108 }
109 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000110 m2.finalize();
111 VERIFY_IS_APPROX(m2,m1);
112 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100113
Gael Guennebaud28293142009-05-04 14:25:12 +0000114 // test insert (fully random)
115 {
116 DenseMatrix m1(rows,cols);
117 m1.setZero();
118 SparseMatrixType m2(rows,cols);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100119 if(internal::random<int>()%2)
120 m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
Gael Guennebaud28293142009-05-04 14:25:12 +0000121 for (int k=0; k<rows*cols; ++k)
122 {
Benoit Jacob47160402010-10-25 10:15:22 -0400123 int i = internal::random<int>(0,rows-1);
124 int j = internal::random<int>(0,cols-1);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100125 if ((m1.coeff(i,j)==Scalar(0)) && (internal::random<int>()%2))
Benoit Jacob47160402010-10-25 10:15:22 -0400126 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100127 else
128 {
129 Scalar v = internal::random<Scalar>();
130 m2.coeffRef(i,j) += v;
131 m1(i,j) += v;
132 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000133 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000134 VERIFY_IS_APPROX(m2,m1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000135 }
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200136
137 // test insert (un-compressed)
138 for(int mode=0;mode<4;++mode)
139 {
140 DenseMatrix m1(rows,cols);
141 m1.setZero();
142 SparseMatrixType m2(rows,cols);
143 VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? m2.innerSize() : std::max<int>(1,m2.innerSize()/8)));
144 m2.reserve(r);
145 for (int k=0; k<rows*cols; ++k)
146 {
147 int i = internal::random<int>(0,rows-1);
148 int j = internal::random<int>(0,cols-1);
149 if (m1.coeff(i,j)==Scalar(0))
150 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
151 if(mode==3)
152 m2.reserve(r);
153 }
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100154 if(internal::random<int>()%2)
155 m2.makeCompressed();
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200156 VERIFY_IS_APPROX(m2,m1);
157 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100158
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000159 // test innerVector()
160 {
161 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
Gael Guennebaud178858f2009-01-19 15:20:45 +0000162 SparseMatrixType m2(rows, rows);
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000163 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200164 int j0 = internal::random<int>(0,rows-1);
165 int j1 = internal::random<int>(0,rows-1);
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100166 if(SparseMatrixType::IsRowMajor)
167 VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.row(j0));
168 else
169 VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.col(j0));
170
171 if(SparseMatrixType::IsRowMajor)
172 VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.row(j0)+refMat2.row(j1));
173 else
174 VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100175
176 SparseMatrixType m3(rows,rows);
177 m3.reserve(VectorXi::Constant(rows,rows/2));
178 for(int j=0; j<rows; ++j)
179 for(int k=0; k<j; ++k)
180 m3.insertByOuterInner(j,k) = k+1;
181 for(int j=0; j<rows; ++j)
182 {
Gael Guennebaud65c59302013-06-12 10:37:15 +0200183 VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100184 if(j>0)
Gael Guennebaud65c59302013-06-12 10:37:15 +0200185 VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100186 }
187 m3.makeCompressed();
188 for(int j=0; j<rows; ++j)
189 {
Gael Guennebaud65c59302013-06-12 10:37:15 +0200190 VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100191 if(j>0)
Gael Guennebaud65c59302013-06-12 10:37:15 +0200192 VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100193 }
194
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000195 //m2.innerVector(j0) = 2*m2.innerVector(j1);
196 //refMat2.col(j0) = 2*refMat2.col(j1);
197 //VERIFY_IS_APPROX(m2, refMat2);
198 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100199
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000200 // test innerVectors()
201 {
202 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
203 SparseMatrixType m2(rows, rows);
204 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud7450b232013-04-12 13:20:13 +0200205 if(internal::random<float>(0,1)>0.5) m2.makeCompressed();
206
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200207 int j0 = internal::random<int>(0,rows-2);
208 int j1 = internal::random<int>(0,rows-2);
Gael Guennebaud42e25782011-08-19 14:18:05 +0200209 int n0 = internal::random<int>(1,rows-(std::max)(j0,j1));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100210 if(SparseMatrixType::IsRowMajor)
211 VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(j0,0,n0,cols));
212 else
213 VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0));
214 if(SparseMatrixType::IsRowMajor)
215 VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
Gael Guennebaud7450b232013-04-12 13:20:13 +0200216 refMat2.middleRows(j0,n0)+refMat2.middleRows(j1,n0));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100217 else
218 VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
219 refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
Gael Guennebaud7450b232013-04-12 13:20:13 +0200220
221 VERIFY_IS_APPROX(m2, refMat2);
222
223 m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0);
224 if(SparseMatrixType::IsRowMajor)
225 refMat2.middleRows(j0,n0) = (refMat2.middleRows(j0,n0) + refMat2.middleRows(j1,n0)).eval();
226 else
227 refMat2.middleCols(j0,n0) = (refMat2.middleCols(j0,n0) + refMat2.middleCols(j1,n0)).eval();
228
229 VERIFY_IS_APPROX(m2, refMat2);
230
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000231 }
Gael Guennebaud4e602832012-11-16 09:02:50 +0100232
233 // test basic computations
234 {
235 DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
236 DenseMatrix refM2 = DenseMatrix::Zero(rows, rows);
237 DenseMatrix refM3 = DenseMatrix::Zero(rows, rows);
238 DenseMatrix refM4 = DenseMatrix::Zero(rows, rows);
239 SparseMatrixType m1(rows, rows);
240 SparseMatrixType m2(rows, rows);
241 SparseMatrixType m3(rows, rows);
242 SparseMatrixType m4(rows, rows);
243 initSparse<Scalar>(density, refM1, m1);
244 initSparse<Scalar>(density, refM2, m2);
245 initSparse<Scalar>(density, refM3, m3);
246 initSparse<Scalar>(density, refM4, m4);
247
248 VERIFY_IS_APPROX(m1+m2, refM1+refM2);
249 VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3);
250 VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2));
251 VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2);
252
253 VERIFY_IS_APPROX(m1*=s1, refM1*=s1);
254 VERIFY_IS_APPROX(m1/=s1, refM1/=s1);
255
256 VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
257 VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
258
259 if(SparseMatrixType::IsRowMajor)
260 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0)));
261 else
262 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.col(0).dot(refM2.row(0)));
263
264 VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate());
265 VERIFY_IS_APPROX(m1.real(), refM1.real());
266
267 refM4.setRandom();
268 // sparse cwise* dense
269 VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4));
270// VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4);
271
272 // test aliasing
273 VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1));
274 VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval()));
275 VERIFY_IS_APPROX((m1 = -m1.transpose()), (refM1 = -refM1.transpose().eval()));
276 VERIFY_IS_APPROX((m1 += -m1), (refM1 += -refM1));
277 }
278
279 // test transpose
280 {
281 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
282 SparseMatrixType m2(rows, rows);
283 initSparse<Scalar>(density, refMat2, m2);
284 VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval());
285 VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
286
287 VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
288 }
289
290
291
292 // test generic blocks
293 {
294 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
295 SparseMatrixType m2(rows, rows);
296 initSparse<Scalar>(density, refMat2, m2);
297 int j0 = internal::random<int>(0,rows-2);
298 int j1 = internal::random<int>(0,rows-2);
299 int n0 = internal::random<int>(1,rows-(std::max)(j0,j1));
300 if(SparseMatrixType::IsRowMajor)
301 VERIFY_IS_APPROX(m2.block(j0,0,n0,cols), refMat2.block(j0,0,n0,cols));
302 else
303 VERIFY_IS_APPROX(m2.block(0,j0,rows,n0), refMat2.block(0,j0,rows,n0));
304
305 if(SparseMatrixType::IsRowMajor)
306 VERIFY_IS_APPROX(m2.block(j0,0,n0,cols)+m2.block(j1,0,n0,cols),
307 refMat2.block(j0,0,n0,cols)+refMat2.block(j1,0,n0,cols));
308 else
309 VERIFY_IS_APPROX(m2.block(0,j0,rows,n0)+m2.block(0,j1,rows,n0),
310 refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
Gael Guennebaud5f178e52013-06-14 10:52:19 +0200311
312 int i = internal::random<int>(0,m2.outerSize()-1);
313 if(SparseMatrixType::IsRowMajor) {
314 m2.innerVector(i) = m2.innerVector(i) * s1;
315 refMat2.row(i) = refMat2.row(i) * s1;
316 VERIFY_IS_APPROX(m2,refMat2);
317 } else {
318 m2.innerVector(i) = m2.innerVector(i) * s1;
319 refMat2.col(i) = refMat2.col(i) * s1;
320 VERIFY_IS_APPROX(m2,refMat2);
321 }
Gael Guennebaud4e602832012-11-16 09:02:50 +0100322 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000323
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000324 // test prune
325 {
326 SparseMatrixType m2(rows, rows);
327 DenseMatrix refM2(rows, rows);
328 refM2.setZero();
329 int countFalseNonZero = 0;
330 int countTrueNonZero = 0;
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000331 for (int j=0; j<m2.outerSize(); ++j)
Gael Guennebaud28293142009-05-04 14:25:12 +0000332 {
333 m2.startVec(j);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000334 for (int i=0; i<m2.innerSize(); ++i)
335 {
Benoit Jacob47160402010-10-25 10:15:22 -0400336 float x = internal::random<float>(0,1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000337 if (x<0.1)
338 {
339 // do nothing
340 }
341 else if (x<0.5)
342 {
343 countFalseNonZero++;
Gael Guennebaud8710bd22010-06-02 13:32:13 +0200344 m2.insertBackByOuterInner(j,i) = Scalar(0);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000345 }
346 else
347 {
348 countTrueNonZero++;
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100349 m2.insertBackByOuterInner(j,i) = Scalar(1);
350 if(SparseMatrixType::IsRowMajor)
351 refM2(j,i) = Scalar(1);
352 else
353 refM2(i,j) = Scalar(1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000354 }
355 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000356 }
357 m2.finalize();
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000358 VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros());
359 VERIFY_IS_APPROX(m2, refM2);
Hauke Heibeld204ec42010-11-02 14:33:33 +0100360 m2.prune(Scalar(1));
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000361 VERIFY(countTrueNonZero==m2.nonZeros());
362 VERIFY_IS_APPROX(m2, refM2);
363 }
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100364
Gael Guennebaud87138072012-01-28 11:13:59 +0100365 // test setFromTriplets
366 {
367 typedef Triplet<Scalar,Index> TripletType;
368 std::vector<TripletType> triplets;
369 int ntriplets = rows*cols;
370 triplets.reserve(ntriplets);
371 DenseMatrix refMat(rows,cols);
372 refMat.setZero();
373 for(int i=0;i<ntriplets;++i)
374 {
Gael Guennebaud18e3ac02012-01-31 09:14:01 +0100375 int r = internal::random<int>(0,rows-1);
376 int c = internal::random<int>(0,cols-1);
Gael Guennebaud87138072012-01-28 11:13:59 +0100377 Scalar v = internal::random<Scalar>();
Gael Guennebaud18e3ac02012-01-31 09:14:01 +0100378 triplets.push_back(TripletType(r,c,v));
379 refMat(r,c) += v;
Gael Guennebaud87138072012-01-28 11:13:59 +0100380 }
381 SparseMatrixType m(rows,cols);
382 m.setFromTriplets(triplets.begin(), triplets.end());
383 VERIFY_IS_APPROX(m, refMat);
384 }
385
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100386 // test triangularView
387 {
388 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
389 SparseMatrixType m2(rows, rows), m3(rows, rows);
390 initSparse<Scalar>(density, refMat2, m2);
391 refMat3 = refMat2.template triangularView<Lower>();
392 m3 = m2.template triangularView<Lower>();
393 VERIFY_IS_APPROX(m3, refMat3);
394
395 refMat3 = refMat2.template triangularView<Upper>();
396 m3 = m2.template triangularView<Upper>();
397 VERIFY_IS_APPROX(m3, refMat3);
398
399 refMat3 = refMat2.template triangularView<UnitUpper>();
400 m3 = m2.template triangularView<UnitUpper>();
401 VERIFY_IS_APPROX(m3, refMat3);
402
403 refMat3 = refMat2.template triangularView<UnitLower>();
404 m3 = m2.template triangularView<UnitLower>();
405 VERIFY_IS_APPROX(m3, refMat3);
Desire NUENTSA4cd82452013-06-11 14:42:29 +0200406
407 refMat3 = refMat2.template triangularView<StrictlyUpper>();
408 m3 = m2.template triangularView<StrictlyUpper>();
409 VERIFY_IS_APPROX(m3, refMat3);
410
411 refMat3 = refMat2.template triangularView<StrictlyLower>();
412 m3 = m2.template triangularView<StrictlyLower>();
413 VERIFY_IS_APPROX(m3, refMat3);
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100414 }
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200415
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100416 // test selfadjointView
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100417 if(!SparseMatrixType::IsRowMajor)
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100418 {
419 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
420 SparseMatrixType m2(rows, rows), m3(rows, rows);
421 initSparse<Scalar>(density, refMat2, m2);
422 refMat3 = refMat2.template selfadjointView<Lower>();
423 m3 = m2.template selfadjointView<Lower>();
424 VERIFY_IS_APPROX(m3, refMat3);
425 }
426
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200427 // test sparseView
428 {
429 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
430 SparseMatrixType m2(rows, rows);
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100431 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200432 VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval());
433 }
Gael Guennebaud82f9aa12011-12-04 21:49:21 +0100434
435 // test diagonal
436 {
437 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
438 SparseMatrixType m2(rows, rows);
439 initSparse<Scalar>(density, refMat2, m2);
440 VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval());
441 }
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200442
443 // test conservative resize
444 {
445 std::vector< std::pair<int,int> > inc;
446 inc.push_back(std::pair<int,int>(-3,-2));
447 inc.push_back(std::pair<int,int>(0,0));
448 inc.push_back(std::pair<int,int>(3,2));
449 inc.push_back(std::pair<int,int>(3,0));
450 inc.push_back(std::pair<int,int>(0,3));
451
452 for(size_t i = 0; i< inc.size(); i++) {
453 int incRows = inc[i].first;
454 int incCols = inc[i].second;
455 SparseMatrixType m1(rows, cols);
456 DenseMatrix refMat1 = DenseMatrix::Zero(rows, cols);
457 initSparse<Scalar>(density, refMat1, m1);
458
459 m1.conservativeResize(rows+incRows, cols+incCols);
460 refMat1.conservativeResize(rows+incRows, cols+incCols);
461 if (incRows > 0) refMat1.bottomRows(incRows).setZero();
462 if (incCols > 0) refMat1.rightCols(incCols).setZero();
463
464 VERIFY_IS_APPROX(m1, refMat1);
465
466 // Insert new values
467 if (incRows > 0)
468 m1.insert(refMat1.rows()-1, 0) = refMat1(refMat1.rows()-1, 0) = 1;
469 if (incCols > 0)
470 m1.insert(0, refMat1.cols()-1) = refMat1(0, refMat1.cols()-1) = 1;
471
472 VERIFY_IS_APPROX(m1, refMat1);
473
474
475 }
476 }
Desire NUENTSA4cd82452013-06-11 14:42:29 +0200477
478 // test Identity matrix
479 {
480 DenseMatrix refMat1 = DenseMatrix::Identity(rows, rows);
481 SparseMatrixType m1(rows, rows);
482 m1.setIdentity();
483 VERIFY_IS_APPROX(m1, refMat1);
484 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000485}
486
487void test_sparse_basic()
488{
489 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200490 int s = Eigen::internal::random<int>(1,50);
Gael Guennebauda1091ca2013-02-15 14:05:37 +0100491 EIGEN_UNUSED_VARIABLE(s);
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200492 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) ));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100493 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(s, s)) ));
494 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(s, s)) ));
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200495 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(s, s)) ));
496 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,long int>(s, s)) ));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100497 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,long int>(s, s)) ));
Gael Guennebaud47e89022013-06-10 10:34:03 +0200498
499 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,short int>(s, s)) ));
500 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,short int>(s, s)) ));
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000501 }
502}