You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

509 lines
17 KiB

  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. #include "sparse.h"
  12. template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& ref)
  13. {
  14. typedef typename SparseMatrixType::Index Index;
  15. typedef Matrix<Index,2,1> Vector2;
  16. const Index rows = ref.rows();
  17. const Index cols = ref.cols();
  18. typedef typename SparseMatrixType::Scalar Scalar;
  19. enum { Flags = SparseMatrixType::Flags };
  20. double density = (std::max)(8./(rows*cols), 0.01);
  21. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  22. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  23. Scalar eps = 1e-6;
  24. Scalar s1 = internal::random<Scalar>();
  25. {
  26. SparseMatrixType m(rows, cols);
  27. DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
  28. DenseVector vec1 = DenseVector::Random(rows);
  29. std::vector<Vector2> zeroCoords;
  30. std::vector<Vector2> nonzeroCoords;
  31. initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords);
  32. if (zeroCoords.size()==0 || nonzeroCoords.size()==0)
  33. return;
  34. // test coeff and coeffRef
  35. for (int i=0; i<(int)zeroCoords.size(); ++i)
  36. {
  37. VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps );
  38. if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value)
  39. VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[0].x(),zeroCoords[0].y()) = 5 );
  40. }
  41. VERIFY_IS_APPROX(m, refMat);
  42. m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
  43. refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
  44. VERIFY_IS_APPROX(m, refMat);
  45. /*
  46. // test InnerIterators and Block expressions
  47. for (int t=0; t<10; ++t)
  48. {
  49. int j = internal::random<int>(0,cols-1);
  50. int i = internal::random<int>(0,rows-1);
  51. int w = internal::random<int>(1,cols-j-1);
  52. int h = internal::random<int>(1,rows-i-1);
  53. // VERIFY_IS_APPROX(m.block(i,j,h,w), refMat.block(i,j,h,w));
  54. for(int c=0; c<w; c++)
  55. {
  56. VERIFY_IS_APPROX(m.block(i,j,h,w).col(c), refMat.block(i,j,h,w).col(c));
  57. for(int r=0; r<h; r++)
  58. {
  59. // VERIFY_IS_APPROX(m.block(i,j,h,w).col(c).coeff(r), refMat.block(i,j,h,w).col(c).coeff(r));
  60. }
  61. }
  62. // for(int r=0; r<h; r++)
  63. // {
  64. // VERIFY_IS_APPROX(m.block(i,j,h,w).row(r), refMat.block(i,j,h,w).row(r));
  65. // for(int c=0; c<w; c++)
  66. // {
  67. // VERIFY_IS_APPROX(m.block(i,j,h,w).row(r).coeff(c), refMat.block(i,j,h,w).row(r).coeff(c));
  68. // }
  69. // }
  70. }
  71. for(int c=0; c<cols; c++)
  72. {
  73. VERIFY_IS_APPROX(m.col(c) + m.col(c), (m + m).col(c));
  74. VERIFY_IS_APPROX(m.col(c) + m.col(c), refMat.col(c) + refMat.col(c));
  75. }
  76. for(int r=0; r<rows; r++)
  77. {
  78. VERIFY_IS_APPROX(m.row(r) + m.row(r), (m + m).row(r));
  79. VERIFY_IS_APPROX(m.row(r) + m.row(r), refMat.row(r) + refMat.row(r));
  80. }
  81. */
  82. // test assertion
  83. VERIFY_RAISES_ASSERT( m.coeffRef(-1,1) = 0 );
  84. VERIFY_RAISES_ASSERT( m.coeffRef(0,m.cols()) = 0 );
  85. }
  86. // test insert (inner random)
  87. {
  88. DenseMatrix m1(rows,cols);
  89. m1.setZero();
  90. SparseMatrixType m2(rows,cols);
  91. if(internal::random<int>()%2)
  92. m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
  93. for (Index j=0; j<cols; ++j)
  94. {
  95. for (Index k=0; k<rows/2; ++k)
  96. {
  97. Index i = internal::random<Index>(0,rows-1);
  98. if (m1.coeff(i,j)==Scalar(0))
  99. m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
  100. }
  101. }
  102. m2.finalize();
  103. VERIFY_IS_APPROX(m2,m1);
  104. }
  105. // test insert (fully random)
  106. {
  107. DenseMatrix m1(rows,cols);
  108. m1.setZero();
  109. SparseMatrixType m2(rows,cols);
  110. if(internal::random<int>()%2)
  111. m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
  112. for (int k=0; k<rows*cols; ++k)
  113. {
  114. Index i = internal::random<Index>(0,rows-1);
  115. Index j = internal::random<Index>(0,cols-1);
  116. if ((m1.coeff(i,j)==Scalar(0)) && (internal::random<int>()%2))
  117. m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
  118. else
  119. {
  120. Scalar v = internal::random<Scalar>();
  121. m2.coeffRef(i,j) += v;
  122. m1(i,j) += v;
  123. }
  124. }
  125. VERIFY_IS_APPROX(m2,m1);
  126. }
  127. // test insert (un-compressed)
  128. for(int mode=0;mode<4;++mode)
  129. {
  130. DenseMatrix m1(rows,cols);
  131. m1.setZero();
  132. SparseMatrixType m2(rows,cols);
  133. VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? m2.innerSize() : std::max<int>(1,m2.innerSize()/8)));
  134. m2.reserve(r);
  135. for (int k=0; k<rows*cols; ++k)
  136. {
  137. Index i = internal::random<Index>(0,rows-1);
  138. Index j = internal::random<Index>(0,cols-1);
  139. if (m1.coeff(i,j)==Scalar(0))
  140. m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
  141. if(mode==3)
  142. m2.reserve(r);
  143. }
  144. if(internal::random<int>()%2)
  145. m2.makeCompressed();
  146. VERIFY_IS_APPROX(m2,m1);
  147. }
  148. // test innerVector()
  149. {
  150. DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
  151. SparseMatrixType m2(rows, rows);
  152. initSparse<Scalar>(density, refMat2, m2);
  153. Index j0 = internal::random<Index>(0,rows-1);
  154. Index j1 = internal::random<Index>(0,rows-1);
  155. if(SparseMatrixType::IsRowMajor)
  156. VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.row(j0));
  157. else
  158. VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.col(j0));
  159. if(SparseMatrixType::IsRowMajor)
  160. VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.row(j0)+refMat2.row(j1));
  161. else
  162. VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1));
  163. SparseMatrixType m3(rows,rows);
  164. m3.reserve(VectorXi::Constant(rows,rows/2));
  165. for(Index j=0; j<rows; ++j)
  166. for(Index k=0; k<j; ++k)
  167. m3.insertByOuterInner(j,k) = k+1;
  168. for(Index j=0; j<rows; ++j)
  169. {
  170. VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
  171. if(j>0)
  172. VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
  173. }
  174. m3.makeCompressed();
  175. for(Index j=0; j<rows; ++j)
  176. {
  177. VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
  178. if(j>0)
  179. VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
  180. }
  181. //m2.innerVector(j0) = 2*m2.innerVector(j1);
  182. //refMat2.col(j0) = 2*refMat2.col(j1);
  183. //VERIFY_IS_APPROX(m2, refMat2);
  184. }
  185. // test innerVectors()
  186. {
  187. DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
  188. SparseMatrixType m2(rows, rows);
  189. initSparse<Scalar>(density, refMat2, m2);
  190. if(internal::random<float>(0,1)>0.5) m2.makeCompressed();
  191. Index j0 = internal::random<Index>(0,rows-2);
  192. Index j1 = internal::random<Index>(0,rows-2);
  193. Index n0 = internal::random<Index>(1,rows-(std::max)(j0,j1));
  194. if(SparseMatrixType::IsRowMajor)
  195. VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(j0,0,n0,cols));
  196. else
  197. VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0));
  198. if(SparseMatrixType::IsRowMajor)
  199. VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
  200. refMat2.middleRows(j0,n0)+refMat2.middleRows(j1,n0));
  201. else
  202. VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
  203. refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
  204. VERIFY_IS_APPROX(m2, refMat2);
  205. m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0);
  206. if(SparseMatrixType::IsRowMajor)
  207. refMat2.middleRows(j0,n0) = (refMat2.middleRows(j0,n0) + refMat2.middleRows(j1,n0)).eval();
  208. else
  209. refMat2.middleCols(j0,n0) = (refMat2.middleCols(j0,n0) + refMat2.middleCols(j1,n0)).eval();
  210. VERIFY_IS_APPROX(m2, refMat2);
  211. }
  212. // test basic computations
  213. {
  214. DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
  215. DenseMatrix refM2 = DenseMatrix::Zero(rows, rows);
  216. DenseMatrix refM3 = DenseMatrix::Zero(rows, rows);
  217. DenseMatrix refM4 = DenseMatrix::Zero(rows, rows);
  218. SparseMatrixType m1(rows, rows);
  219. SparseMatrixType m2(rows, rows);
  220. SparseMatrixType m3(rows, rows);
  221. SparseMatrixType m4(rows, rows);
  222. initSparse<Scalar>(density, refM1, m1);
  223. initSparse<Scalar>(density, refM2, m2);
  224. initSparse<Scalar>(density, refM3, m3);
  225. initSparse<Scalar>(density, refM4, m4);
  226. VERIFY_IS_APPROX(m1+m2, refM1+refM2);
  227. VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3);
  228. VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2));
  229. VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2);
  230. VERIFY_IS_APPROX(m1*=s1, refM1*=s1);
  231. VERIFY_IS_APPROX(m1/=s1, refM1/=s1);
  232. VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
  233. VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
  234. if(SparseMatrixType::IsRowMajor)
  235. VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0)));
  236. else
  237. VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.col(0).dot(refM2.row(0)));
  238. VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate());
  239. VERIFY_IS_APPROX(m1.real(), refM1.real());
  240. refM4.setRandom();
  241. // sparse cwise* dense
  242. VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4));
  243. // VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4);
  244. // test aliasing
  245. VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1));
  246. VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval()));
  247. VERIFY_IS_APPROX((m1 = -m1.transpose()), (refM1 = -refM1.transpose().eval()));
  248. VERIFY_IS_APPROX((m1 += -m1), (refM1 += -refM1));
  249. }
  250. // test transpose
  251. {
  252. DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
  253. SparseMatrixType m2(rows, rows);
  254. initSparse<Scalar>(density, refMat2, m2);
  255. VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval());
  256. VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
  257. VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
  258. }
  259. // test generic blocks
  260. {
  261. DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
  262. SparseMatrixType m2(rows, rows);
  263. initSparse<Scalar>(density, refMat2, m2);
  264. Index j0 = internal::random<Index>(0,rows-2);
  265. Index j1 = internal::random<Index>(0,rows-2);
  266. Index n0 = internal::random<Index>(1,rows-(std::max)(j0,j1));
  267. if(SparseMatrixType::IsRowMajor)
  268. VERIFY_IS_APPROX(m2.block(j0,0,n0,cols), refMat2.block(j0,0,n0,cols));
  269. else
  270. VERIFY_IS_APPROX(m2.block(0,j0,rows,n0), refMat2.block(0,j0,rows,n0));
  271. if(SparseMatrixType::IsRowMajor)
  272. VERIFY_IS_APPROX(m2.block(j0,0,n0,cols)+m2.block(j1,0,n0,cols),
  273. refMat2.block(j0,0,n0,cols)+refMat2.block(j1,0,n0,cols));
  274. else
  275. VERIFY_IS_APPROX(m2.block(0,j0,rows,n0)+m2.block(0,j1,rows,n0),
  276. refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
  277. Index i = internal::random<Index>(0,m2.outerSize()-1);
  278. if(SparseMatrixType::IsRowMajor) {
  279. m2.innerVector(i) = m2.innerVector(i) * s1;
  280. refMat2.row(i) = refMat2.row(i) * s1;
  281. VERIFY_IS_APPROX(m2,refMat2);
  282. } else {
  283. m2.innerVector(i) = m2.innerVector(i) * s1;
  284. refMat2.col(i) = refMat2.col(i) * s1;
  285. VERIFY_IS_APPROX(m2,refMat2);
  286. }
  287. }
  288. // test prune
  289. {
  290. SparseMatrixType m2(rows, rows);
  291. DenseMatrix refM2(rows, rows);
  292. refM2.setZero();
  293. int countFalseNonZero = 0;
  294. int countTrueNonZero = 0;
  295. for (Index j=0; j<m2.outerSize(); ++j)
  296. {
  297. m2.startVec(j);
  298. for (Index i=0; i<m2.innerSize(); ++i)
  299. {
  300. float x = internal::random<float>(0,1);
  301. if (x<0.1)
  302. {
  303. // do nothing
  304. }
  305. else if (x<0.5)
  306. {
  307. countFalseNonZero++;
  308. m2.insertBackByOuterInner(j,i) = Scalar(0);
  309. }
  310. else
  311. {
  312. countTrueNonZero++;
  313. m2.insertBackByOuterInner(j,i) = Scalar(1);
  314. if(SparseMatrixType::IsRowMajor)
  315. refM2(j,i) = Scalar(1);
  316. else
  317. refM2(i,j) = Scalar(1);
  318. }
  319. }
  320. }
  321. m2.finalize();
  322. VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros());
  323. VERIFY_IS_APPROX(m2, refM2);
  324. m2.prune(Scalar(1));
  325. VERIFY(countTrueNonZero==m2.nonZeros());
  326. VERIFY_IS_APPROX(m2, refM2);
  327. }
  328. // test setFromTriplets
  329. {
  330. typedef Triplet<Scalar,Index> TripletType;
  331. std::vector<TripletType> triplets;
  332. int ntriplets = rows*cols;
  333. triplets.reserve(ntriplets);
  334. DenseMatrix refMat(rows,cols);
  335. refMat.setZero();
  336. for(int i=0;i<ntriplets;++i)
  337. {
  338. Index r = internal::random<Index>(0,rows-1);
  339. Index c = internal::random<Index>(0,cols-1);
  340. Scalar v = internal::random<Scalar>();
  341. triplets.push_back(TripletType(r,c,v));
  342. refMat(r,c) += v;
  343. }
  344. SparseMatrixType m(rows,cols);
  345. m.setFromTriplets(triplets.begin(), triplets.end());
  346. VERIFY_IS_APPROX(m, refMat);
  347. }
  348. // test triangularView
  349. {
  350. DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
  351. SparseMatrixType m2(rows, rows), m3(rows, rows);
  352. initSparse<Scalar>(density, refMat2, m2);
  353. refMat3 = refMat2.template triangularView<Lower>();
  354. m3 = m2.template triangularView<Lower>();
  355. VERIFY_IS_APPROX(m3, refMat3);
  356. refMat3 = refMat2.template triangularView<Upper>();
  357. m3 = m2.template triangularView<Upper>();
  358. VERIFY_IS_APPROX(m3, refMat3);
  359. refMat3 = refMat2.template triangularView<UnitUpper>();
  360. m3 = m2.template triangularView<UnitUpper>();
  361. VERIFY_IS_APPROX(m3, refMat3);
  362. refMat3 = refMat2.template triangularView<UnitLower>();
  363. m3 = m2.template triangularView<UnitLower>();
  364. VERIFY_IS_APPROX(m3, refMat3);
  365. refMat3 = refMat2.template triangularView<StrictlyUpper>();
  366. m3 = m2.template triangularView<StrictlyUpper>();
  367. VERIFY_IS_APPROX(m3, refMat3);
  368. refMat3 = refMat2.template triangularView<StrictlyLower>();
  369. m3 = m2.template triangularView<StrictlyLower>();
  370. VERIFY_IS_APPROX(m3, refMat3);
  371. }
  372. // test selfadjointView
  373. if(!SparseMatrixType::IsRowMajor)
  374. {
  375. DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
  376. SparseMatrixType m2(rows, rows), m3(rows, rows);
  377. initSparse<Scalar>(density, refMat2, m2);
  378. refMat3 = refMat2.template selfadjointView<Lower>();
  379. m3 = m2.template selfadjointView<Lower>();
  380. VERIFY_IS_APPROX(m3, refMat3);
  381. }
  382. // test sparseView
  383. {
  384. DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
  385. SparseMatrixType m2(rows, rows);
  386. initSparse<Scalar>(density, refMat2, m2);
  387. VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval());
  388. }
  389. // test diagonal
  390. {
  391. DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
  392. SparseMatrixType m2(rows, rows);
  393. initSparse<Scalar>(density, refMat2, m2);
  394. VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval());
  395. }
  396. // test conservative resize
  397. {
  398. std::vector< std::pair<Index,Index> > inc;
  399. inc.push_back(std::pair<Index,Index>(-3,-2));
  400. inc.push_back(std::pair<Index,Index>(0,0));
  401. inc.push_back(std::pair<Index,Index>(3,2));
  402. inc.push_back(std::pair<Index,Index>(3,0));
  403. inc.push_back(std::pair<Index,Index>(0,3));
  404. for(size_t i = 0; i< inc.size(); i++) {
  405. Index incRows = inc[i].first;
  406. Index incCols = inc[i].second;
  407. SparseMatrixType m1(rows, cols);
  408. DenseMatrix refMat1 = DenseMatrix::Zero(rows, cols);
  409. initSparse<Scalar>(density, refMat1, m1);
  410. m1.conservativeResize(rows+incRows, cols+incCols);
  411. refMat1.conservativeResize(rows+incRows, cols+incCols);
  412. if (incRows > 0) refMat1.bottomRows(incRows).setZero();
  413. if (incCols > 0) refMat1.rightCols(incCols).setZero();
  414. VERIFY_IS_APPROX(m1, refMat1);
  415. // Insert new values
  416. if (incRows > 0)
  417. m1.insert(m1.rows()-1, 0) = refMat1(refMat1.rows()-1, 0) = 1;
  418. if (incCols > 0)
  419. m1.insert(0, m1.cols()-1) = refMat1(0, refMat1.cols()-1) = 1;
  420. VERIFY_IS_APPROX(m1, refMat1);
  421. }
  422. }
  423. // test Identity matrix
  424. {
  425. DenseMatrix refMat1 = DenseMatrix::Identity(rows, rows);
  426. SparseMatrixType m1(rows, rows);
  427. m1.setIdentity();
  428. VERIFY_IS_APPROX(m1, refMat1);
  429. }
  430. }
  431. void test_sparse_basic()
  432. {
  433. for(int i = 0; i < g_repeat; i++) {
  434. int s = Eigen::internal::random<int>(1,50);
  435. EIGEN_UNUSED_VARIABLE(s);
  436. CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) ));
  437. CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(s, s)) ));
  438. CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(s, s)) ));
  439. CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(s, s)) ));
  440. CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,long int>(s, s)) ));
  441. CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,long int>(s, s)) ));
  442. CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,short int>(short(s), short(s))) ));
  443. CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,short int>(short(s), short(s))) ));
  444. }
  445. }