sparse_basic.cpp 25 KB

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