sparse_basic.cpp 29 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. // Copyright (C) 2008 Daniel Gomez Ferro <dgomezferro@gmail.com>
  6. // Copyright (C) 2013 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
  7. //
  8. // 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/.
  11. #ifndef EIGEN_SPARSE_TEST_INCLUDED_FROM_SPARSE_EXTRA
  12. static long g_realloc_count = 0;
  13. #define EIGEN_SPARSE_COMPRESSED_STORAGE_REALLOCATE_PLUGIN g_realloc_count++;
  14. static long g_dense_op_sparse_count = 0;
  15. #define EIGEN_SPARSE_ASSIGNMENT_FROM_DENSE_OP_SPARSE_PLUGIN g_dense_op_sparse_count++;
  16. #define EIGEN_SPARSE_ASSIGNMENT_FROM_SPARSE_ADD_DENSE_PLUGIN g_dense_op_sparse_count+=10;
  17. #define EIGEN_SPARSE_ASSIGNMENT_FROM_SPARSE_SUB_DENSE_PLUGIN g_dense_op_sparse_count+=20;
  18. #endif
  19. #include "sparse.h"
  20. template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& ref)
  21. {
  22. typedef typename SparseMatrixType::StorageIndex StorageIndex;
  23. typedef Matrix<StorageIndex,2,1> Vector2;
  24. const Index rows = ref.rows();
  25. const Index cols = ref.cols();
  26. //const Index inner = ref.innerSize();
  27. //const Index outer = ref.outerSize();
  28. typedef typename SparseMatrixType::Scalar Scalar;
  29. typedef typename SparseMatrixType::RealScalar RealScalar;
  30. enum { Flags = SparseMatrixType::Flags };
  31. double density = (std::max)(8./(rows*cols), 0.01);
  32. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  33. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  34. Scalar eps = 1e-6;
  35. Scalar s1 = internal::random<Scalar>();
  36. {
  37. SparseMatrixType m(rows, cols);
  38. DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
  39. DenseVector vec1 = DenseVector::Random(rows);
  40. std::vector<Vector2> zeroCoords;
  41. std::vector<Vector2> nonzeroCoords;
  42. initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords);
  43. // test coeff and coeffRef
  44. for (std::size_t i=0; i<zeroCoords.size(); ++i)
  45. {
  46. VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps );
  47. if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value)
  48. VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[i].x(),zeroCoords[i].y()) = 5 );
  49. }
  50. VERIFY_IS_APPROX(m, refMat);
  51. if(!nonzeroCoords.empty()) {
  52. m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
  53. refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
  54. }
  55. VERIFY_IS_APPROX(m, refMat);
  56. // test assertion
  57. VERIFY_RAISES_ASSERT( m.coeffRef(-1,1) = 0 );
  58. VERIFY_RAISES_ASSERT( m.coeffRef(0,m.cols()) = 0 );
  59. }
  60. // test insert (inner random)
  61. {
  62. DenseMatrix m1(rows,cols);
  63. m1.setZero();
  64. SparseMatrixType m2(rows,cols);
  65. bool call_reserve = internal::random<int>()%2;
  66. Index nnz = internal::random<int>(1,int(rows)/2);
  67. if(call_reserve)
  68. {
  69. if(internal::random<int>()%2)
  70. m2.reserve(VectorXi::Constant(m2.outerSize(), int(nnz)));
  71. else
  72. m2.reserve(m2.outerSize() * nnz);
  73. }
  74. g_realloc_count = 0;
  75. for (Index j=0; j<cols; ++j)
  76. {
  77. for (Index k=0; k<nnz; ++k)
  78. {
  79. Index i = internal::random<Index>(0,rows-1);
  80. if (m1.coeff(i,j)==Scalar(0))
  81. m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
  82. }
  83. }
  84. if(call_reserve && !SparseMatrixType::IsRowMajor)
  85. {
  86. VERIFY(g_realloc_count==0);
  87. }
  88. m2.finalize();
  89. VERIFY_IS_APPROX(m2,m1);
  90. }
  91. // test insert (fully random)
  92. {
  93. DenseMatrix m1(rows,cols);
  94. m1.setZero();
  95. SparseMatrixType m2(rows,cols);
  96. if(internal::random<int>()%2)
  97. m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
  98. for (int k=0; k<rows*cols; ++k)
  99. {
  100. Index i = internal::random<Index>(0,rows-1);
  101. Index j = internal::random<Index>(0,cols-1);
  102. if ((m1.coeff(i,j)==Scalar(0)) && (internal::random<int>()%2))
  103. m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
  104. else
  105. {
  106. Scalar v = internal::random<Scalar>();
  107. m2.coeffRef(i,j) += v;
  108. m1(i,j) += v;
  109. }
  110. }
  111. VERIFY_IS_APPROX(m2,m1);
  112. }
  113. // test insert (un-compressed)
  114. for(int mode=0;mode<4;++mode)
  115. {
  116. DenseMatrix m1(rows,cols);
  117. m1.setZero();
  118. SparseMatrixType m2(rows,cols);
  119. VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? int(m2.innerSize()) : std::max<int>(1,int(m2.innerSize())/8)));
  120. m2.reserve(r);
  121. for (Index k=0; k<rows*cols; ++k)
  122. {
  123. Index i = internal::random<Index>(0,rows-1);
  124. Index j = internal::random<Index>(0,cols-1);
  125. if (m1.coeff(i,j)==Scalar(0))
  126. m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
  127. if(mode==3)
  128. m2.reserve(r);
  129. }
  130. if(internal::random<int>()%2)
  131. m2.makeCompressed();
  132. VERIFY_IS_APPROX(m2,m1);
  133. }
  134. // test basic computations
  135. {
  136. DenseMatrix refM1 = DenseMatrix::Zero(rows, cols);
  137. DenseMatrix refM2 = DenseMatrix::Zero(rows, cols);
  138. DenseMatrix refM3 = DenseMatrix::Zero(rows, cols);
  139. DenseMatrix refM4 = DenseMatrix::Zero(rows, cols);
  140. SparseMatrixType m1(rows, cols);
  141. SparseMatrixType m2(rows, cols);
  142. SparseMatrixType m3(rows, cols);
  143. SparseMatrixType m4(rows, cols);
  144. initSparse<Scalar>(density, refM1, m1);
  145. initSparse<Scalar>(density, refM2, m2);
  146. initSparse<Scalar>(density, refM3, m3);
  147. initSparse<Scalar>(density, refM4, m4);
  148. if(internal::random<bool>())
  149. m1.makeCompressed();
  150. Index m1_nnz = m1.nonZeros();
  151. VERIFY_IS_APPROX(m1*s1, refM1*s1);
  152. VERIFY_IS_APPROX(m1+m2, refM1+refM2);
  153. VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3);
  154. VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2));
  155. VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2);
  156. VERIFY_IS_APPROX(m4=m1/s1, refM1/s1);
  157. VERIFY_IS_EQUAL(m4.nonZeros(), m1_nnz);
  158. if(SparseMatrixType::IsRowMajor)
  159. VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0)));
  160. else
  161. VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.col(0)), refM1.col(0).dot(refM2.col(0)));
  162. DenseVector rv = DenseVector::Random(m1.cols());
  163. DenseVector cv = DenseVector::Random(m1.rows());
  164. Index r = internal::random<Index>(0,m1.rows()-2);
  165. Index c = internal::random<Index>(0,m1.cols()-1);
  166. VERIFY_IS_APPROX(( m1.template block<1,Dynamic>(r,0,1,m1.cols()).dot(rv)) , refM1.row(r).dot(rv));
  167. VERIFY_IS_APPROX(m1.row(r).dot(rv), refM1.row(r).dot(rv));
  168. VERIFY_IS_APPROX(m1.col(c).dot(cv), refM1.col(c).dot(cv));
  169. VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate());
  170. VERIFY_IS_APPROX(m1.real(), refM1.real());
  171. refM4.setRandom();
  172. // sparse cwise* dense
  173. VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4));
  174. // dense cwise* sparse
  175. VERIFY_IS_APPROX(refM4.cwiseProduct(m3), refM4.cwiseProduct(refM3));
  176. // VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4);
  177. // mixed sparse-dense
  178. VERIFY_IS_APPROX(refM4 + m3, refM4 + refM3);
  179. VERIFY_IS_APPROX(m3 + refM4, refM3 + refM4);
  180. VERIFY_IS_APPROX(refM4 - m3, refM4 - refM3);
  181. VERIFY_IS_APPROX(m3 - refM4, refM3 - refM4);
  182. VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + RealScalar(0.5)*m3).eval(), RealScalar(0.5)*refM4 + RealScalar(0.5)*refM3);
  183. VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + m3*RealScalar(0.5)).eval(), RealScalar(0.5)*refM4 + RealScalar(0.5)*refM3);
  184. VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + m3.cwiseProduct(m3)).eval(), RealScalar(0.5)*refM4 + refM3.cwiseProduct(refM3));
  185. VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + RealScalar(0.5)*m3).eval(), RealScalar(0.5)*refM4 + RealScalar(0.5)*refM3);
  186. VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + m3*RealScalar(0.5)).eval(), RealScalar(0.5)*refM4 + RealScalar(0.5)*refM3);
  187. VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + (m3+m3)).eval(), RealScalar(0.5)*refM4 + (refM3+refM3));
  188. VERIFY_IS_APPROX(((refM3+m3)+RealScalar(0.5)*m3).eval(), RealScalar(0.5)*refM3 + (refM3+refM3));
  189. VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + (refM3+m3)).eval(), RealScalar(0.5)*refM4 + (refM3+refM3));
  190. VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + (m3+refM3)).eval(), RealScalar(0.5)*refM4 + (refM3+refM3));
  191. VERIFY_IS_APPROX(m1.sum(), refM1.sum());
  192. m4 = m1; refM4 = m4;
  193. VERIFY_IS_APPROX(m1*=s1, refM1*=s1);
  194. VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
  195. VERIFY_IS_APPROX(m1/=s1, refM1/=s1);
  196. VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
  197. VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
  198. VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
  199. refM3 = refM1;
  200. VERIFY_IS_APPROX(refM1+=m2, refM3+=refM2);
  201. VERIFY_IS_APPROX(refM1-=m2, refM3-=refM2);
  202. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1 =m2+refM4, refM3 =refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,10);
  203. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1+=m2+refM4, refM3+=refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
  204. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1-=m2+refM4, refM3-=refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
  205. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1 =refM4+m2, refM3 =refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
  206. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1+=refM4+m2, refM3+=refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
  207. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1-=refM4+m2, refM3-=refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
  208. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1 =m2-refM4, refM3 =refM2-refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,20);
  209. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1+=m2-refM4, refM3+=refM2-refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
  210. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1-=m2-refM4, refM3-=refM2-refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
  211. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1 =refM4-m2, refM3 =refM4-refM2); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
  212. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1+=refM4-m2, refM3+=refM4-refM2); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
  213. g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1-=refM4-m2, refM3-=refM4-refM2); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
  214. refM3 = m3;
  215. if (rows>=2 && cols>=2)
  216. {
  217. VERIFY_RAISES_ASSERT( m1 += m1.innerVector(0) );
  218. VERIFY_RAISES_ASSERT( m1 -= m1.innerVector(0) );
  219. VERIFY_RAISES_ASSERT( refM1 -= m1.innerVector(0) );
  220. VERIFY_RAISES_ASSERT( refM1 += m1.innerVector(0) );
  221. }
  222. m1 = m4; refM1 = refM4;
  223. // test aliasing
  224. VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1));
  225. VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
  226. m1 = m4; refM1 = refM4;
  227. VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval()));
  228. VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
  229. m1 = m4; refM1 = refM4;
  230. VERIFY_IS_APPROX((m1 = -m1.transpose()), (refM1 = -refM1.transpose().eval()));
  231. VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
  232. m1 = m4; refM1 = refM4;
  233. VERIFY_IS_APPROX((m1 += -m1), (refM1 += -refM1));
  234. VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
  235. m1 = m4; refM1 = refM4;
  236. if(m1.isCompressed())
  237. {
  238. VERIFY_IS_APPROX(m1.coeffs().sum(), m1.sum());
  239. m1.coeffs() += s1;
  240. for(Index j = 0; j<m1.outerSize(); ++j)
  241. for(typename SparseMatrixType::InnerIterator it(m1,j); it; ++it)
  242. refM1(it.row(), it.col()) += s1;
  243. VERIFY_IS_APPROX(m1, refM1);
  244. }
  245. // and/or
  246. {
  247. typedef SparseMatrix<bool, SparseMatrixType::Options, typename SparseMatrixType::StorageIndex> SpBool;
  248. SpBool mb1 = m1.real().template cast<bool>();
  249. SpBool mb2 = m2.real().template cast<bool>();
  250. VERIFY_IS_EQUAL(mb1.template cast<int>().sum(), refM1.real().template cast<bool>().count());
  251. VERIFY_IS_EQUAL((mb1 && mb2).template cast<int>().sum(), (refM1.real().template cast<bool>() && refM2.real().template cast<bool>()).count());
  252. VERIFY_IS_EQUAL((mb1 || mb2).template cast<int>().sum(), (refM1.real().template cast<bool>() || refM2.real().template cast<bool>()).count());
  253. SpBool mb3 = mb1 && mb2;
  254. if(mb1.coeffs().all() && mb2.coeffs().all())
  255. {
  256. VERIFY_IS_EQUAL(mb3.nonZeros(), (refM1.real().template cast<bool>() && refM2.real().template cast<bool>()).count());
  257. }
  258. }
  259. }
  260. // test reverse iterators
  261. {
  262. DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
  263. SparseMatrixType m2(rows, cols);
  264. initSparse<Scalar>(density, refMat2, m2);
  265. std::vector<Scalar> ref_value(m2.innerSize());
  266. std::vector<Index> ref_index(m2.innerSize());
  267. if(internal::random<bool>())
  268. m2.makeCompressed();
  269. for(Index j = 0; j<m2.outerSize(); ++j)
  270. {
  271. Index count_forward = 0;
  272. for(typename SparseMatrixType::InnerIterator it(m2,j); it; ++it)
  273. {
  274. ref_value[ref_value.size()-1-count_forward] = it.value();
  275. ref_index[ref_index.size()-1-count_forward] = it.index();
  276. count_forward++;
  277. }
  278. Index count_reverse = 0;
  279. for(typename SparseMatrixType::ReverseInnerIterator it(m2,j); it; --it)
  280. {
  281. VERIFY_IS_APPROX( std::abs(ref_value[ref_value.size()-count_forward+count_reverse])+1, std::abs(it.value())+1);
  282. VERIFY_IS_EQUAL( ref_index[ref_index.size()-count_forward+count_reverse] , it.index());
  283. count_reverse++;
  284. }
  285. VERIFY_IS_EQUAL(count_forward, count_reverse);
  286. }
  287. }
  288. // test transpose
  289. {
  290. DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
  291. SparseMatrixType m2(rows, cols);
  292. initSparse<Scalar>(density, refMat2, m2);
  293. VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval());
  294. VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
  295. VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
  296. // check isApprox handles opposite storage order
  297. typename Transpose<SparseMatrixType>::PlainObject m3(m2);
  298. VERIFY(m2.isApprox(m3));
  299. }
  300. // test prune
  301. {
  302. SparseMatrixType m2(rows, cols);
  303. DenseMatrix refM2(rows, cols);
  304. refM2.setZero();
  305. int countFalseNonZero = 0;
  306. int countTrueNonZero = 0;
  307. m2.reserve(VectorXi::Constant(m2.outerSize(), int(m2.innerSize())));
  308. for (Index j=0; j<m2.cols(); ++j)
  309. {
  310. for (Index i=0; i<m2.rows(); ++i)
  311. {
  312. float x = internal::random<float>(0,1);
  313. if (x<0.1f)
  314. {
  315. // do nothing
  316. }
  317. else if (x<0.5f)
  318. {
  319. countFalseNonZero++;
  320. m2.insert(i,j) = Scalar(0);
  321. }
  322. else
  323. {
  324. countTrueNonZero++;
  325. m2.insert(i,j) = Scalar(1);
  326. refM2(i,j) = Scalar(1);
  327. }
  328. }
  329. }
  330. if(internal::random<bool>())
  331. m2.makeCompressed();
  332. VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros());
  333. if(countTrueNonZero>0)
  334. VERIFY_IS_APPROX(m2, refM2);
  335. m2.prune(Scalar(1));
  336. VERIFY(countTrueNonZero==m2.nonZeros());
  337. VERIFY_IS_APPROX(m2, refM2);
  338. }
  339. // test setFromTriplets
  340. {
  341. typedef Triplet<Scalar,StorageIndex> TripletType;
  342. std::vector<TripletType> triplets;
  343. Index ntriplets = rows*cols;
  344. triplets.reserve(ntriplets);
  345. DenseMatrix refMat_sum = DenseMatrix::Zero(rows,cols);
  346. DenseMatrix refMat_prod = DenseMatrix::Zero(rows,cols);
  347. DenseMatrix refMat_last = DenseMatrix::Zero(rows,cols);
  348. for(Index i=0;i<ntriplets;++i)
  349. {
  350. StorageIndex r = internal::random<StorageIndex>(0,StorageIndex(rows-1));
  351. StorageIndex c = internal::random<StorageIndex>(0,StorageIndex(cols-1));
  352. Scalar v = internal::random<Scalar>();
  353. triplets.push_back(TripletType(r,c,v));
  354. refMat_sum(r,c) += v;
  355. if(std::abs(refMat_prod(r,c))==0)
  356. refMat_prod(r,c) = v;
  357. else
  358. refMat_prod(r,c) *= v;
  359. refMat_last(r,c) = v;
  360. }
  361. SparseMatrixType m(rows,cols);
  362. m.setFromTriplets(triplets.begin(), triplets.end());
  363. VERIFY_IS_APPROX(m, refMat_sum);
  364. m.setFromTriplets(triplets.begin(), triplets.end(), std::multiplies<Scalar>());
  365. VERIFY_IS_APPROX(m, refMat_prod);
  366. #if (EIGEN_COMP_CXXVER >= 11)
  367. m.setFromTriplets(triplets.begin(), triplets.end(), [] (Scalar,Scalar b) { return b; });
  368. VERIFY_IS_APPROX(m, refMat_last);
  369. #endif
  370. }
  371. // test Map
  372. {
  373. DenseMatrix refMat2(rows, cols), refMat3(rows, cols);
  374. SparseMatrixType m2(rows, cols), m3(rows, cols);
  375. initSparse<Scalar>(density, refMat2, m2);
  376. initSparse<Scalar>(density, refMat3, m3);
  377. {
  378. Map<SparseMatrixType> mapMat2(m2.rows(), m2.cols(), m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr());
  379. Map<SparseMatrixType> mapMat3(m3.rows(), m3.cols(), m3.nonZeros(), m3.outerIndexPtr(), m3.innerIndexPtr(), m3.valuePtr(), m3.innerNonZeroPtr());
  380. VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
  381. VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
  382. }
  383. {
  384. MappedSparseMatrix<Scalar,SparseMatrixType::Options,StorageIndex> mapMat2(m2.rows(), m2.cols(), m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr());
  385. MappedSparseMatrix<Scalar,SparseMatrixType::Options,StorageIndex> mapMat3(m3.rows(), m3.cols(), m3.nonZeros(), m3.outerIndexPtr(), m3.innerIndexPtr(), m3.valuePtr(), m3.innerNonZeroPtr());
  386. VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
  387. VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
  388. }
  389. Index i = internal::random<Index>(0,rows-1);
  390. Index j = internal::random<Index>(0,cols-1);
  391. m2.coeffRef(i,j) = 123;
  392. if(internal::random<bool>())
  393. m2.makeCompressed();
  394. Map<SparseMatrixType> mapMat2(rows, cols, m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr());
  395. VERIFY_IS_EQUAL(m2.coeff(i,j),Scalar(123));
  396. VERIFY_IS_EQUAL(mapMat2.coeff(i,j),Scalar(123));
  397. mapMat2.coeffRef(i,j) = -123;
  398. VERIFY_IS_EQUAL(m2.coeff(i,j),Scalar(-123));
  399. }
  400. // test triangularView
  401. {
  402. DenseMatrix refMat2(rows, cols), refMat3(rows, cols);
  403. SparseMatrixType m2(rows, cols), m3(rows, cols);
  404. initSparse<Scalar>(density, refMat2, m2);
  405. refMat3 = refMat2.template triangularView<Lower>();
  406. m3 = m2.template triangularView<Lower>();
  407. VERIFY_IS_APPROX(m3, refMat3);
  408. refMat3 = refMat2.template triangularView<Upper>();
  409. m3 = m2.template triangularView<Upper>();
  410. VERIFY_IS_APPROX(m3, refMat3);
  411. {
  412. refMat3 = refMat2.template triangularView<UnitUpper>();
  413. m3 = m2.template triangularView<UnitUpper>();
  414. VERIFY_IS_APPROX(m3, refMat3);
  415. refMat3 = refMat2.template triangularView<UnitLower>();
  416. m3 = m2.template triangularView<UnitLower>();
  417. VERIFY_IS_APPROX(m3, refMat3);
  418. }
  419. refMat3 = refMat2.template triangularView<StrictlyUpper>();
  420. m3 = m2.template triangularView<StrictlyUpper>();
  421. VERIFY_IS_APPROX(m3, refMat3);
  422. refMat3 = refMat2.template triangularView<StrictlyLower>();
  423. m3 = m2.template triangularView<StrictlyLower>();
  424. VERIFY_IS_APPROX(m3, refMat3);
  425. // check sparse-triangular to dense
  426. refMat3 = m2.template triangularView<StrictlyUpper>();
  427. VERIFY_IS_APPROX(refMat3, DenseMatrix(refMat2.template triangularView<StrictlyUpper>()));
  428. }
  429. // test selfadjointView
  430. if(!SparseMatrixType::IsRowMajor)
  431. {
  432. DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
  433. SparseMatrixType m2(rows, rows), m3(rows, rows);
  434. initSparse<Scalar>(density, refMat2, m2);
  435. refMat3 = refMat2.template selfadjointView<Lower>();
  436. m3 = m2.template selfadjointView<Lower>();
  437. VERIFY_IS_APPROX(m3, refMat3);
  438. refMat3 += refMat2.template selfadjointView<Lower>();
  439. m3 += m2.template selfadjointView<Lower>();
  440. VERIFY_IS_APPROX(m3, refMat3);
  441. refMat3 -= refMat2.template selfadjointView<Lower>();
  442. m3 -= m2.template selfadjointView<Lower>();
  443. VERIFY_IS_APPROX(m3, refMat3);
  444. // selfadjointView only works for square matrices:
  445. SparseMatrixType m4(rows, rows+1);
  446. VERIFY_RAISES_ASSERT(m4.template selfadjointView<Lower>());
  447. VERIFY_RAISES_ASSERT(m4.template selfadjointView<Upper>());
  448. }
  449. // test sparseView
  450. {
  451. DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
  452. SparseMatrixType m2(rows, rows);
  453. initSparse<Scalar>(density, refMat2, m2);
  454. VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval());
  455. // sparse view on expressions:
  456. VERIFY_IS_APPROX((s1*m2).eval(), (s1*refMat2).sparseView().eval());
  457. VERIFY_IS_APPROX((m2+m2).eval(), (refMat2+refMat2).sparseView().eval());
  458. VERIFY_IS_APPROX((m2*m2).eval(), (refMat2.lazyProduct(refMat2)).sparseView().eval());
  459. VERIFY_IS_APPROX((m2*m2).eval(), (refMat2*refMat2).sparseView().eval());
  460. }
  461. // test diagonal
  462. {
  463. DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
  464. SparseMatrixType m2(rows, cols);
  465. initSparse<Scalar>(density, refMat2, m2);
  466. VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval());
  467. DenseVector d = m2.diagonal();
  468. VERIFY_IS_APPROX(d, refMat2.diagonal().eval());
  469. d = m2.diagonal().array();
  470. VERIFY_IS_APPROX(d, refMat2.diagonal().eval());
  471. VERIFY_IS_APPROX(const_cast<const SparseMatrixType&>(m2).diagonal(), refMat2.diagonal().eval());
  472. initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag);
  473. m2.diagonal() += refMat2.diagonal();
  474. refMat2.diagonal() += refMat2.diagonal();
  475. VERIFY_IS_APPROX(m2, refMat2);
  476. }
  477. // test diagonal to sparse
  478. {
  479. DenseVector d = DenseVector::Random(rows);
  480. DenseMatrix refMat2 = d.asDiagonal();
  481. SparseMatrixType m2;
  482. m2 = d.asDiagonal();
  483. VERIFY_IS_APPROX(m2, refMat2);
  484. SparseMatrixType m3(d.asDiagonal());
  485. VERIFY_IS_APPROX(m3, refMat2);
  486. refMat2 += d.asDiagonal();
  487. m2 += d.asDiagonal();
  488. VERIFY_IS_APPROX(m2, refMat2);
  489. m2.setZero(); m2 += d.asDiagonal();
  490. refMat2.setZero(); refMat2 += d.asDiagonal();
  491. VERIFY_IS_APPROX(m2, refMat2);
  492. m2.setZero(); m2 -= d.asDiagonal();
  493. refMat2.setZero(); refMat2 -= d.asDiagonal();
  494. VERIFY_IS_APPROX(m2, refMat2);
  495. initSparse<Scalar>(density, refMat2, m2);
  496. m2.makeCompressed();
  497. m2 += d.asDiagonal();
  498. refMat2 += d.asDiagonal();
  499. VERIFY_IS_APPROX(m2, refMat2);
  500. initSparse<Scalar>(density, refMat2, m2);
  501. m2.makeCompressed();
  502. VectorXi res(rows);
  503. for(Index i=0; i<rows; ++i)
  504. res(i) = internal::random<int>(0,3);
  505. m2.reserve(res);
  506. m2 -= d.asDiagonal();
  507. refMat2 -= d.asDiagonal();
  508. VERIFY_IS_APPROX(m2, refMat2);
  509. }
  510. // test conservative resize
  511. {
  512. std::vector< std::pair<StorageIndex,StorageIndex> > inc;
  513. if(rows > 3 && cols > 2)
  514. inc.push_back(std::pair<StorageIndex,StorageIndex>(-3,-2));
  515. inc.push_back(std::pair<StorageIndex,StorageIndex>(0,0));
  516. inc.push_back(std::pair<StorageIndex,StorageIndex>(3,2));
  517. inc.push_back(std::pair<StorageIndex,StorageIndex>(3,0));
  518. inc.push_back(std::pair<StorageIndex,StorageIndex>(0,3));
  519. inc.push_back(std::pair<StorageIndex,StorageIndex>(0,-1));
  520. inc.push_back(std::pair<StorageIndex,StorageIndex>(-1,0));
  521. inc.push_back(std::pair<StorageIndex,StorageIndex>(-1,-1));
  522. for(size_t i = 0; i< inc.size(); i++) {
  523. StorageIndex incRows = inc[i].first;
  524. StorageIndex incCols = inc[i].second;
  525. SparseMatrixType m1(rows, cols);
  526. DenseMatrix refMat1 = DenseMatrix::Zero(rows, cols);
  527. initSparse<Scalar>(density, refMat1, m1);
  528. SparseMatrixType m2 = m1;
  529. m2.makeCompressed();
  530. m1.conservativeResize(rows+incRows, cols+incCols);
  531. m2.conservativeResize(rows+incRows, cols+incCols);
  532. refMat1.conservativeResize(rows+incRows, cols+incCols);
  533. if (incRows > 0) refMat1.bottomRows(incRows).setZero();
  534. if (incCols > 0) refMat1.rightCols(incCols).setZero();
  535. VERIFY_IS_APPROX(m1, refMat1);
  536. VERIFY_IS_APPROX(m2, refMat1);
  537. // Insert new values
  538. if (incRows > 0)
  539. m1.insert(m1.rows()-1, 0) = refMat1(refMat1.rows()-1, 0) = 1;
  540. if (incCols > 0)
  541. m1.insert(0, m1.cols()-1) = refMat1(0, refMat1.cols()-1) = 1;
  542. VERIFY_IS_APPROX(m1, refMat1);
  543. }
  544. }
  545. // test Identity matrix
  546. {
  547. DenseMatrix refMat1 = DenseMatrix::Identity(rows, rows);
  548. SparseMatrixType m1(rows, rows);
  549. m1.setIdentity();
  550. VERIFY_IS_APPROX(m1, refMat1);
  551. for(int k=0; k<rows*rows/4; ++k)
  552. {
  553. Index i = internal::random<Index>(0,rows-1);
  554. Index j = internal::random<Index>(0,rows-1);
  555. Scalar v = internal::random<Scalar>();
  556. m1.coeffRef(i,j) = v;
  557. refMat1.coeffRef(i,j) = v;
  558. VERIFY_IS_APPROX(m1, refMat1);
  559. if(internal::random<Index>(0,10)<2)
  560. m1.makeCompressed();
  561. }
  562. m1.setIdentity();
  563. refMat1.setIdentity();
  564. VERIFY_IS_APPROX(m1, refMat1);
  565. }
  566. // test array/vector of InnerIterator
  567. {
  568. typedef typename SparseMatrixType::InnerIterator IteratorType;
  569. DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
  570. SparseMatrixType m2(rows, cols);
  571. initSparse<Scalar>(density, refMat2, m2);
  572. IteratorType static_array[2];
  573. static_array[0] = IteratorType(m2,0);
  574. static_array[1] = IteratorType(m2,m2.outerSize()-1);
  575. VERIFY( static_array[0] || m2.innerVector(static_array[0].outer()).nonZeros() == 0 );
  576. VERIFY( static_array[1] || m2.innerVector(static_array[1].outer()).nonZeros() == 0 );
  577. if(static_array[0] && static_array[1])
  578. {
  579. ++(static_array[1]);
  580. static_array[1] = IteratorType(m2,0);
  581. VERIFY( static_array[1] );
  582. VERIFY( static_array[1].index() == static_array[0].index() );
  583. VERIFY( static_array[1].outer() == static_array[0].outer() );
  584. VERIFY( static_array[1].value() == static_array[0].value() );
  585. }
  586. std::vector<IteratorType> iters(2);
  587. iters[0] = IteratorType(m2,0);
  588. iters[1] = IteratorType(m2,m2.outerSize()-1);
  589. }
  590. // test reserve with empty rows/columns
  591. {
  592. SparseMatrixType m1(0,cols);
  593. m1.reserve(ArrayXi::Constant(m1.outerSize(),1));
  594. SparseMatrixType m2(rows,0);
  595. m2.reserve(ArrayXi::Constant(m2.outerSize(),1));
  596. }
  597. }
  598. template<typename SparseMatrixType>
  599. void big_sparse_triplet(Index rows, Index cols, double density) {
  600. typedef typename SparseMatrixType::StorageIndex StorageIndex;
  601. typedef typename SparseMatrixType::Scalar Scalar;
  602. typedef Triplet<Scalar,Index> TripletType;
  603. std::vector<TripletType> triplets;
  604. double nelements = density * rows*cols;
  605. VERIFY(nelements>=0 && nelements < static_cast<double>(NumTraits<StorageIndex>::highest()));
  606. Index ntriplets = Index(nelements);
  607. triplets.reserve(ntriplets);
  608. Scalar sum = Scalar(0);
  609. for(Index i=0;i<ntriplets;++i)
  610. {
  611. Index r = internal::random<Index>(0,rows-1);
  612. Index c = internal::random<Index>(0,cols-1);
  613. // use positive values to prevent numerical cancellation errors in sum
  614. Scalar v = numext::abs(internal::random<Scalar>());
  615. triplets.push_back(TripletType(r,c,v));
  616. sum += v;
  617. }
  618. SparseMatrixType m(rows,cols);
  619. m.setFromTriplets(triplets.begin(), triplets.end());
  620. VERIFY(m.nonZeros() <= ntriplets);
  621. VERIFY_IS_APPROX(sum, m.sum());
  622. }
  623. template<int>
  624. void bug1105()
  625. {
  626. // Regression test for bug 1105
  627. int n = Eigen::internal::random<int>(200,600);
  628. SparseMatrix<std::complex<double>,0, long> mat(n, n);
  629. std::complex<double> val;
  630. for(int i=0; i<n; ++i)
  631. {
  632. mat.coeffRef(i, i%(n/10)) = val;
  633. VERIFY(mat.data().allocatedSize()<20*n);
  634. }
  635. }
  636. #ifndef EIGEN_SPARSE_TEST_INCLUDED_FROM_SPARSE_EXTRA
  637. EIGEN_DECLARE_TEST(sparse_basic)
  638. {
  639. g_dense_op_sparse_count = 0; // Suppresses compiler warning.
  640. for(int i = 0; i < g_repeat; i++) {
  641. int r = Eigen::internal::random<int>(1,200), c = Eigen::internal::random<int>(1,200);
  642. if(Eigen::internal::random<int>(0,4) == 0) {
  643. r = c; // check square matrices in 25% of tries
  644. }
  645. EIGEN_UNUSED_VARIABLE(r+c);
  646. CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(1, 1)) ));
  647. CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) ));
  648. CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(r, c)) ));
  649. CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(r, c)) ));
  650. CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(r, c)) ));
  651. CALL_SUBTEST_5(( sparse_basic(SparseMatrix<double,ColMajor,long int>(r, c)) ));
  652. CALL_SUBTEST_5(( sparse_basic(SparseMatrix<double,RowMajor,long int>(r, c)) ));
  653. r = Eigen::internal::random<int>(1,100);
  654. c = Eigen::internal::random<int>(1,100);
  655. if(Eigen::internal::random<int>(0,4) == 0) {
  656. r = c; // check square matrices in 25% of tries
  657. }
  658. CALL_SUBTEST_6(( sparse_basic(SparseMatrix<double,ColMajor,short int>(short(r), short(c))) ));
  659. CALL_SUBTEST_6(( sparse_basic(SparseMatrix<double,RowMajor,short int>(short(r), short(c))) ));
  660. }
  661. // Regression test for bug 900: (manually insert higher values here, if you have enough RAM):
  662. CALL_SUBTEST_3((big_sparse_triplet<SparseMatrix<float, RowMajor, int> >(10000, 10000, 0.125)));
  663. CALL_SUBTEST_4((big_sparse_triplet<SparseMatrix<double, ColMajor, long int> >(10000, 10000, 0.125)));
  664. CALL_SUBTEST_7( bug1105<0>() );
  665. }
  666. #endif