blob: c86534bad275bf5e44451aebb732d6a8c948a241 [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;
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +020017 typedef Matrix<Index,2,1> Vector2;
18
Hauke Heibelf1679c72010-06-20 17:37:56 +020019 const Index rows = ref.rows();
20 const Index cols = ref.cols();
Gael Guennebaud178858f2009-01-19 15:20:45 +000021 typedef typename SparseMatrixType::Scalar Scalar;
22 enum { Flags = SparseMatrixType::Flags };
Gael Guennebaud9f795582009-12-16 19:18:40 +010023
Gael Guennebaud42e25782011-08-19 14:18:05 +020024 double density = (std::max)(8./(rows*cols), 0.01);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000025 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
26 typedef Matrix<Scalar,Dynamic,1> DenseVector;
27 Scalar eps = 1e-6;
28
Benoit Jacob47160402010-10-25 10:15:22 -040029 Scalar s1 = internal::random<Scalar>();
Gael Guennebaud86ccd992008-11-05 13:47:55 +000030 {
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020031 SparseMatrixType m(rows, cols);
32 DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
33 DenseVector vec1 = DenseVector::Random(rows);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000034
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +020035 std::vector<Vector2> zeroCoords;
36 std::vector<Vector2> nonzeroCoords;
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020037 initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000038
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020039 if (zeroCoords.size()==0 || nonzeroCoords.size()==0)
40 return;
Gael Guennebaud86ccd992008-11-05 13:47:55 +000041
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020042 // test coeff and coeffRef
43 for (int i=0; i<(int)zeroCoords.size(); ++i)
Gael Guennebaud86ccd992008-11-05 13:47:55 +000044 {
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020045 VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps );
46 if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value)
47 VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[0].x(),zeroCoords[0].y()) = 5 );
Gael Guennebaud86ccd992008-11-05 13:47:55 +000048 }
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020049 VERIFY_IS_APPROX(m, refMat);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000050
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020051 m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
52 refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000053
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020054 VERIFY_IS_APPROX(m, refMat);
55 /*
56 // test InnerIterators and Block expressions
57 for (int t=0; t<10; ++t)
58 {
59 int j = internal::random<int>(0,cols-1);
60 int i = internal::random<int>(0,rows-1);
61 int w = internal::random<int>(1,cols-j-1);
62 int h = internal::random<int>(1,rows-i-1);
63
64 // VERIFY_IS_APPROX(m.block(i,j,h,w), refMat.block(i,j,h,w));
65 for(int c=0; c<w; c++)
66 {
67 VERIFY_IS_APPROX(m.block(i,j,h,w).col(c), refMat.block(i,j,h,w).col(c));
68 for(int r=0; r<h; r++)
69 {
70 // VERIFY_IS_APPROX(m.block(i,j,h,w).col(c).coeff(r), refMat.block(i,j,h,w).col(c).coeff(r));
71 }
72 }
73 // for(int r=0; r<h; r++)
74 // {
75 // VERIFY_IS_APPROX(m.block(i,j,h,w).row(r), refMat.block(i,j,h,w).row(r));
76 // for(int c=0; c<w; c++)
77 // {
78 // VERIFY_IS_APPROX(m.block(i,j,h,w).row(r).coeff(c), refMat.block(i,j,h,w).row(r).coeff(c));
79 // }
80 // }
81 }
82
83 for(int c=0; c<cols; c++)
84 {
85 VERIFY_IS_APPROX(m.col(c) + m.col(c), (m + m).col(c));
86 VERIFY_IS_APPROX(m.col(c) + m.col(c), refMat.col(c) + refMat.col(c));
87 }
88
89 for(int r=0; r<rows; r++)
90 {
91 VERIFY_IS_APPROX(m.row(r) + m.row(r), (m + m).row(r));
92 VERIFY_IS_APPROX(m.row(r) + m.row(r), refMat.row(r) + refMat.row(r));
93 }
94 */
Gael Guennebauda915f022013-06-28 16:16:02 +020095
96 // test assertion
97 VERIFY_RAISES_ASSERT( m.coeffRef(-1,1) = 0 );
98 VERIFY_RAISES_ASSERT( m.coeffRef(0,m.cols()) = 0 );
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020099 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000100
Gael Guennebaud28293142009-05-04 14:25:12 +0000101 // test insert (inner random)
Gael Guennebaud5015e482008-12-11 18:26:24 +0000102 {
103 DenseMatrix m1(rows,cols);
104 m1.setZero();
Gael Guennebaud178858f2009-01-19 15:20:45 +0000105 SparseMatrixType m2(rows,cols);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100106 if(internal::random<int>()%2)
107 m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200108 for (Index j=0; j<cols; ++j)
Gael Guennebaud5015e482008-12-11 18:26:24 +0000109 {
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200110 for (Index k=0; k<rows/2; ++k)
Gael Guennebaud5015e482008-12-11 18:26:24 +0000111 {
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200112 Index i = internal::random<Index>(0,rows-1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000113 if (m1.coeff(i,j)==Scalar(0))
Benoit Jacob47160402010-10-25 10:15:22 -0400114 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud5015e482008-12-11 18:26:24 +0000115 }
116 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000117 m2.finalize();
118 VERIFY_IS_APPROX(m2,m1);
119 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100120
Gael Guennebaud28293142009-05-04 14:25:12 +0000121 // test insert (fully random)
122 {
123 DenseMatrix m1(rows,cols);
124 m1.setZero();
125 SparseMatrixType m2(rows,cols);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100126 if(internal::random<int>()%2)
127 m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
Gael Guennebaud28293142009-05-04 14:25:12 +0000128 for (int k=0; k<rows*cols; ++k)
129 {
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200130 Index i = internal::random<Index>(0,rows-1);
131 Index j = internal::random<Index>(0,cols-1);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100132 if ((m1.coeff(i,j)==Scalar(0)) && (internal::random<int>()%2))
Benoit Jacob47160402010-10-25 10:15:22 -0400133 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100134 else
135 {
136 Scalar v = internal::random<Scalar>();
137 m2.coeffRef(i,j) += v;
138 m1(i,j) += v;
139 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000140 }
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);
Gael Guennebaudb47ef142014-07-08 16:47:11 +0200150 VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? int(m2.innerSize()) : std::max<int>(1,int(m2.innerSize())/8)));
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200151 m2.reserve(r);
152 for (int k=0; k<rows*cols; ++k)
153 {
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200154 Index i = internal::random<Index>(0,rows-1);
155 Index j = internal::random<Index>(0,cols-1);
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200156 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 }
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100161 if(internal::random<int>()%2)
162 m2.makeCompressed();
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200163 VERIFY_IS_APPROX(m2,m1);
164 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100165
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000166 // test innerVector()
167 {
168 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
Gael Guennebaud178858f2009-01-19 15:20:45 +0000169 SparseMatrixType m2(rows, rows);
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000170 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200171 Index j0 = internal::random<Index>(0,rows-1);
172 Index j1 = internal::random<Index>(0,rows-1);
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100173 if(SparseMatrixType::IsRowMajor)
174 VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.row(j0));
175 else
176 VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.col(j0));
177
178 if(SparseMatrixType::IsRowMajor)
179 VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.row(j0)+refMat2.row(j1));
180 else
181 VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100182
183 SparseMatrixType m3(rows,rows);
Gael Guennebaudb47ef142014-07-08 16:47:11 +0200184 m3.reserve(VectorXi::Constant(rows,int(rows/2)));
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200185 for(Index j=0; j<rows; ++j)
186 for(Index k=0; k<j; ++k)
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100187 m3.insertByOuterInner(j,k) = k+1;
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200188 for(Index j=0; j<rows; ++j)
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100189 {
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 m3.makeCompressed();
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200195 for(Index j=0; j<rows; ++j)
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100196 {
Gael Guennebaud65c59302013-06-12 10:37:15 +0200197 VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100198 if(j>0)
Gael Guennebaud65c59302013-06-12 10:37:15 +0200199 VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100200 }
Gael Guennebaud5c2d1b42013-11-10 15:26:07 +0100201
202 VERIFY(m3.innerVector(j0).nonZeros() == m3.transpose().innerVector(j0).nonZeros());
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100203
Gael Guennebaud4aac8722014-07-22 12:54:03 +0200204// m2.innerVector(j0) = 2*m2.innerVector(j1);
205// refMat2.col(j0) = 2*refMat2.col(j1);
206// VERIFY_IS_APPROX(m2, refMat2);
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000207 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100208
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000209 // test innerVectors()
210 {
211 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
212 SparseMatrixType m2(rows, rows);
213 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud7450b232013-04-12 13:20:13 +0200214 if(internal::random<float>(0,1)>0.5) m2.makeCompressed();
215
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200216 Index j0 = internal::random<Index>(0,rows-2);
217 Index j1 = internal::random<Index>(0,rows-2);
218 Index n0 = internal::random<Index>(1,rows-(std::max)(j0,j1));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100219 if(SparseMatrixType::IsRowMajor)
220 VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(j0,0,n0,cols));
221 else
222 VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0));
223 if(SparseMatrixType::IsRowMajor)
224 VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
Gael Guennebaud7450b232013-04-12 13:20:13 +0200225 refMat2.middleRows(j0,n0)+refMat2.middleRows(j1,n0));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100226 else
227 VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
228 refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
Gael Guennebaud7450b232013-04-12 13:20:13 +0200229
230 VERIFY_IS_APPROX(m2, refMat2);
231
Gael Guennebaud5c2d1b42013-11-10 15:26:07 +0100232 VERIFY(m2.innerVectors(j0,n0).nonZeros() == m2.transpose().innerVectors(j0,n0).nonZeros());
233
Gael Guennebaud7450b232013-04-12 13:20:13 +0200234 m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0);
235 if(SparseMatrixType::IsRowMajor)
236 refMat2.middleRows(j0,n0) = (refMat2.middleRows(j0,n0) + refMat2.middleRows(j1,n0)).eval();
237 else
238 refMat2.middleCols(j0,n0) = (refMat2.middleCols(j0,n0) + refMat2.middleCols(j1,n0)).eval();
239
240 VERIFY_IS_APPROX(m2, refMat2);
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000241 }
Gael Guennebaud4aac8722014-07-22 12:54:03 +0200242
Gael Guennebaud4e602832012-11-16 09:02:50 +0100243 // test basic computations
244 {
245 DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
246 DenseMatrix refM2 = DenseMatrix::Zero(rows, rows);
247 DenseMatrix refM3 = DenseMatrix::Zero(rows, rows);
248 DenseMatrix refM4 = DenseMatrix::Zero(rows, rows);
249 SparseMatrixType m1(rows, rows);
250 SparseMatrixType m2(rows, rows);
251 SparseMatrixType m3(rows, rows);
252 SparseMatrixType m4(rows, rows);
253 initSparse<Scalar>(density, refM1, m1);
254 initSparse<Scalar>(density, refM2, m2);
255 initSparse<Scalar>(density, refM3, m3);
256 initSparse<Scalar>(density, refM4, m4);
257
Gael Guennebaud4aac8722014-07-22 12:54:03 +0200258 VERIFY_IS_APPROX(m1*s1, refM1*s1);
Gael Guennebaud4e602832012-11-16 09:02:50 +0100259 VERIFY_IS_APPROX(m1+m2, refM1+refM2);
260 VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3);
261 VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2));
262 VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2);
263
264 VERIFY_IS_APPROX(m1*=s1, refM1*=s1);
265 VERIFY_IS_APPROX(m1/=s1, refM1/=s1);
266
267 VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
268 VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
269
270 if(SparseMatrixType::IsRowMajor)
271 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0)));
272 else
273 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.col(0).dot(refM2.row(0)));
Gael Guennebaud3573a102014-02-17 13:46:17 +0100274
275 DenseVector rv = DenseVector::Random(m1.cols());
276 DenseVector cv = DenseVector::Random(m1.rows());
277 Index r = internal::random<Index>(0,m1.rows()-2);
278 Index c = internal::random<Index>(0,m1.cols()-1);
279 VERIFY_IS_APPROX(( m1.template block<1,Dynamic>(r,0,1,m1.cols()).dot(rv)) , refM1.row(r).dot(rv));
280 VERIFY_IS_APPROX(m1.row(r).dot(rv), refM1.row(r).dot(rv));
281 VERIFY_IS_APPROX(m1.col(c).dot(cv), refM1.col(c).dot(cv));
Gael Guennebaud4e602832012-11-16 09:02:50 +0100282
283 VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate());
284 VERIFY_IS_APPROX(m1.real(), refM1.real());
285
286 refM4.setRandom();
287 // sparse cwise* dense
288 VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4));
289// VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4);
290
291 // test aliasing
292 VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1));
293 VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval()));
294 VERIFY_IS_APPROX((m1 = -m1.transpose()), (refM1 = -refM1.transpose().eval()));
295 VERIFY_IS_APPROX((m1 += -m1), (refM1 += -refM1));
296 }
297
298 // test transpose
299 {
300 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
301 SparseMatrixType m2(rows, rows);
302 initSparse<Scalar>(density, refMat2, m2);
303 VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval());
304 VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
305
306 VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
307 }
308
309
310
311 // test generic blocks
312 {
313 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
314 SparseMatrixType m2(rows, rows);
315 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200316 Index j0 = internal::random<Index>(0,rows-2);
317 Index j1 = internal::random<Index>(0,rows-2);
318 Index n0 = internal::random<Index>(1,rows-(std::max)(j0,j1));
Gael Guennebaud4e602832012-11-16 09:02:50 +0100319 if(SparseMatrixType::IsRowMajor)
320 VERIFY_IS_APPROX(m2.block(j0,0,n0,cols), refMat2.block(j0,0,n0,cols));
321 else
322 VERIFY_IS_APPROX(m2.block(0,j0,rows,n0), refMat2.block(0,j0,rows,n0));
323
324 if(SparseMatrixType::IsRowMajor)
325 VERIFY_IS_APPROX(m2.block(j0,0,n0,cols)+m2.block(j1,0,n0,cols),
326 refMat2.block(j0,0,n0,cols)+refMat2.block(j1,0,n0,cols));
327 else
328 VERIFY_IS_APPROX(m2.block(0,j0,rows,n0)+m2.block(0,j1,rows,n0),
329 refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
Gael Guennebaud5f178e52013-06-14 10:52:19 +0200330
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200331 Index i = internal::random<Index>(0,m2.outerSize()-1);
Gael Guennebaud5f178e52013-06-14 10:52:19 +0200332 if(SparseMatrixType::IsRowMajor) {
333 m2.innerVector(i) = m2.innerVector(i) * s1;
334 refMat2.row(i) = refMat2.row(i) * s1;
335 VERIFY_IS_APPROX(m2,refMat2);
336 } else {
337 m2.innerVector(i) = m2.innerVector(i) * s1;
338 refMat2.col(i) = refMat2.col(i) * s1;
339 VERIFY_IS_APPROX(m2,refMat2);
340 }
Gael Guennebaud4e602832012-11-16 09:02:50 +0100341 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000342
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000343 // test prune
344 {
345 SparseMatrixType m2(rows, rows);
346 DenseMatrix refM2(rows, rows);
347 refM2.setZero();
348 int countFalseNonZero = 0;
349 int countTrueNonZero = 0;
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200350 for (Index j=0; j<m2.outerSize(); ++j)
Gael Guennebaud28293142009-05-04 14:25:12 +0000351 {
352 m2.startVec(j);
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200353 for (Index i=0; i<m2.innerSize(); ++i)
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000354 {
Benoit Jacob47160402010-10-25 10:15:22 -0400355 float x = internal::random<float>(0,1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000356 if (x<0.1)
357 {
358 // do nothing
359 }
360 else if (x<0.5)
361 {
362 countFalseNonZero++;
Gael Guennebaud8710bd22010-06-02 13:32:13 +0200363 m2.insertBackByOuterInner(j,i) = Scalar(0);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000364 }
365 else
366 {
367 countTrueNonZero++;
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100368 m2.insertBackByOuterInner(j,i) = Scalar(1);
369 if(SparseMatrixType::IsRowMajor)
370 refM2(j,i) = Scalar(1);
371 else
372 refM2(i,j) = Scalar(1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000373 }
374 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000375 }
376 m2.finalize();
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000377 VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros());
378 VERIFY_IS_APPROX(m2, refM2);
Hauke Heibeld204ec42010-11-02 14:33:33 +0100379 m2.prune(Scalar(1));
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000380 VERIFY(countTrueNonZero==m2.nonZeros());
381 VERIFY_IS_APPROX(m2, refM2);
382 }
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100383
Gael Guennebaud87138072012-01-28 11:13:59 +0100384 // test setFromTriplets
385 {
386 typedef Triplet<Scalar,Index> TripletType;
387 std::vector<TripletType> triplets;
Gael Guennebaudb47ef142014-07-08 16:47:11 +0200388 Index ntriplets = rows*cols;
Gael Guennebaud87138072012-01-28 11:13:59 +0100389 triplets.reserve(ntriplets);
390 DenseMatrix refMat(rows,cols);
391 refMat.setZero();
Gael Guennebaudb47ef142014-07-08 16:47:11 +0200392 for(Index i=0;i<ntriplets;++i)
Gael Guennebaud87138072012-01-28 11:13:59 +0100393 {
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200394 Index r = internal::random<Index>(0,rows-1);
395 Index c = internal::random<Index>(0,cols-1);
Gael Guennebaud87138072012-01-28 11:13:59 +0100396 Scalar v = internal::random<Scalar>();
Gael Guennebaud18e3ac02012-01-31 09:14:01 +0100397 triplets.push_back(TripletType(r,c,v));
398 refMat(r,c) += v;
Gael Guennebaud87138072012-01-28 11:13:59 +0100399 }
400 SparseMatrixType m(rows,cols);
401 m.setFromTriplets(triplets.begin(), triplets.end());
402 VERIFY_IS_APPROX(m, refMat);
403 }
404
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100405 // test triangularView
406 {
407 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
408 SparseMatrixType m2(rows, rows), m3(rows, rows);
409 initSparse<Scalar>(density, refMat2, m2);
410 refMat3 = refMat2.template triangularView<Lower>();
411 m3 = m2.template triangularView<Lower>();
412 VERIFY_IS_APPROX(m3, refMat3);
413
414 refMat3 = refMat2.template triangularView<Upper>();
415 m3 = m2.template triangularView<Upper>();
416 VERIFY_IS_APPROX(m3, refMat3);
417
418 refMat3 = refMat2.template triangularView<UnitUpper>();
419 m3 = m2.template triangularView<UnitUpper>();
420 VERIFY_IS_APPROX(m3, refMat3);
421
422 refMat3 = refMat2.template triangularView<UnitLower>();
423 m3 = m2.template triangularView<UnitLower>();
424 VERIFY_IS_APPROX(m3, refMat3);
Desire NUENTSA4cd82452013-06-11 14:42:29 +0200425
426 refMat3 = refMat2.template triangularView<StrictlyUpper>();
427 m3 = m2.template triangularView<StrictlyUpper>();
428 VERIFY_IS_APPROX(m3, refMat3);
429
430 refMat3 = refMat2.template triangularView<StrictlyLower>();
431 m3 = m2.template triangularView<StrictlyLower>();
432 VERIFY_IS_APPROX(m3, refMat3);
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100433 }
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200434
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100435 // test selfadjointView
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100436 if(!SparseMatrixType::IsRowMajor)
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100437 {
438 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
439 SparseMatrixType m2(rows, rows), m3(rows, rows);
440 initSparse<Scalar>(density, refMat2, m2);
441 refMat3 = refMat2.template selfadjointView<Lower>();
442 m3 = m2.template selfadjointView<Lower>();
443 VERIFY_IS_APPROX(m3, refMat3);
444 }
445
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200446 // test sparseView
447 {
448 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
449 SparseMatrixType m2(rows, rows);
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100450 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200451 VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval());
452 }
Gael Guennebaud82f9aa12011-12-04 21:49:21 +0100453
454 // test diagonal
455 {
456 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
457 SparseMatrixType m2(rows, rows);
458 initSparse<Scalar>(density, refMat2, m2);
459 VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval());
460 }
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200461
462 // test conservative resize
463 {
Gael Guennebaud7ee378d2013-07-12 16:40:02 +0200464 std::vector< std::pair<Index,Index> > inc;
465 inc.push_back(std::pair<Index,Index>(-3,-2));
466 inc.push_back(std::pair<Index,Index>(0,0));
467 inc.push_back(std::pair<Index,Index>(3,2));
468 inc.push_back(std::pair<Index,Index>(3,0));
469 inc.push_back(std::pair<Index,Index>(0,3));
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200470
471 for(size_t i = 0; i< inc.size(); i++) {
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200472 Index incRows = inc[i].first;
473 Index incCols = inc[i].second;
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200474 SparseMatrixType m1(rows, cols);
475 DenseMatrix refMat1 = DenseMatrix::Zero(rows, cols);
476 initSparse<Scalar>(density, refMat1, m1);
477
478 m1.conservativeResize(rows+incRows, cols+incCols);
479 refMat1.conservativeResize(rows+incRows, cols+incCols);
480 if (incRows > 0) refMat1.bottomRows(incRows).setZero();
481 if (incCols > 0) refMat1.rightCols(incCols).setZero();
482
483 VERIFY_IS_APPROX(m1, refMat1);
484
485 // Insert new values
486 if (incRows > 0)
Gael Guennebaud7ee378d2013-07-12 16:40:02 +0200487 m1.insert(m1.rows()-1, 0) = refMat1(refMat1.rows()-1, 0) = 1;
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200488 if (incCols > 0)
Gael Guennebaud7ee378d2013-07-12 16:40:02 +0200489 m1.insert(0, m1.cols()-1) = refMat1(0, refMat1.cols()-1) = 1;
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200490
491 VERIFY_IS_APPROX(m1, refMat1);
492
493
494 }
495 }
Desire NUENTSA4cd82452013-06-11 14:42:29 +0200496
497 // test Identity matrix
498 {
499 DenseMatrix refMat1 = DenseMatrix::Identity(rows, rows);
500 SparseMatrixType m1(rows, rows);
501 m1.setIdentity();
502 VERIFY_IS_APPROX(m1, refMat1);
503 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000504}
505
506void test_sparse_basic()
507{
508 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200509 int s = Eigen::internal::random<int>(1,50);
Gael Guennebauda1091ca2013-02-15 14:05:37 +0100510 EIGEN_UNUSED_VARIABLE(s);
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200511 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) ));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100512 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(s, s)) ));
513 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(s, s)) ));
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200514 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(s, s)) ));
515 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,long int>(s, s)) ));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100516 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,long int>(s, s)) ));
Gael Guennebaud47e89022013-06-10 10:34:03 +0200517
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200518 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,short int>(short(s), short(s))) ));
519 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,short int>(short(s), short(s))) ));
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000520 }
521}