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.

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