sparse_block.cpp 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008-2015 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #include "sparse.h"
  10. #include "AnnoyingScalar.h"
  11. template<typename T>
  12. typename Eigen::internal::enable_if<(T::Flags&RowMajorBit)==RowMajorBit, typename T::RowXpr>::type
  13. innervec(T& A, Index i)
  14. {
  15. return A.row(i);
  16. }
  17. template<typename T>
  18. typename Eigen::internal::enable_if<(T::Flags&RowMajorBit)==0, typename T::ColXpr>::type
  19. innervec(T& A, Index i)
  20. {
  21. return A.col(i);
  22. }
  23. template<typename SparseMatrixType> void sparse_block(const SparseMatrixType& ref)
  24. {
  25. const Index rows = ref.rows();
  26. const Index cols = ref.cols();
  27. const Index inner = ref.innerSize();
  28. const Index outer = ref.outerSize();
  29. typedef typename SparseMatrixType::Scalar Scalar;
  30. typedef typename SparseMatrixType::RealScalar RealScalar;
  31. typedef typename SparseMatrixType::StorageIndex StorageIndex;
  32. double density = (std::max)(8./(rows*cols), 0.01);
  33. typedef Matrix<Scalar,Dynamic,Dynamic,SparseMatrixType::IsRowMajor?RowMajor:ColMajor> DenseMatrix;
  34. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  35. typedef Matrix<Scalar,1,Dynamic> RowDenseVector;
  36. typedef SparseVector<Scalar> SparseVectorType;
  37. Scalar s1 = internal::random<Scalar>();
  38. {
  39. SparseMatrixType m(rows, cols);
  40. DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
  41. initSparse<Scalar>(density, refMat, m);
  42. VERIFY_IS_APPROX(m, refMat);
  43. // test InnerIterators and Block expressions
  44. for (int t=0; t<10; ++t)
  45. {
  46. Index j = internal::random<Index>(0,cols-2);
  47. Index i = internal::random<Index>(0,rows-2);
  48. Index w = internal::random<Index>(1,cols-j);
  49. Index h = internal::random<Index>(1,rows-i);
  50. VERIFY_IS_APPROX(m.block(i,j,h,w), refMat.block(i,j,h,w));
  51. for(Index 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(Index 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. VERIFY_IS_APPROX(m.block(i,j,h,w).coeff(r,c), refMat.block(i,j,h,w).coeff(r,c));
  58. }
  59. }
  60. for(Index r=0; r<h; r++)
  61. {
  62. VERIFY_IS_APPROX(m.block(i,j,h,w).row(r), refMat.block(i,j,h,w).row(r));
  63. for(Index c=0; c<w; c++)
  64. {
  65. VERIFY_IS_APPROX(m.block(i,j,h,w).row(r).coeff(c), refMat.block(i,j,h,w).row(r).coeff(c));
  66. VERIFY_IS_APPROX(m.block(i,j,h,w).coeff(r,c), refMat.block(i,j,h,w).coeff(r,c));
  67. }
  68. }
  69. VERIFY_IS_APPROX(m.middleCols(j,w), refMat.middleCols(j,w));
  70. VERIFY_IS_APPROX(m.middleRows(i,h), refMat.middleRows(i,h));
  71. for(Index r=0; r<h; r++)
  72. {
  73. VERIFY_IS_APPROX(m.middleCols(j,w).row(r), refMat.middleCols(j,w).row(r));
  74. VERIFY_IS_APPROX(m.middleRows(i,h).row(r), refMat.middleRows(i,h).row(r));
  75. for(Index c=0; c<w; c++)
  76. {
  77. VERIFY_IS_APPROX(m.col(c).coeff(r), refMat.col(c).coeff(r));
  78. VERIFY_IS_APPROX(m.row(r).coeff(c), refMat.row(r).coeff(c));
  79. VERIFY_IS_APPROX(m.middleCols(j,w).coeff(r,c), refMat.middleCols(j,w).coeff(r,c));
  80. VERIFY_IS_APPROX(m.middleRows(i,h).coeff(r,c), refMat.middleRows(i,h).coeff(r,c));
  81. if(m.middleCols(j,w).coeff(r,c) != Scalar(0))
  82. {
  83. VERIFY_IS_APPROX(m.middleCols(j,w).coeffRef(r,c), refMat.middleCols(j,w).coeff(r,c));
  84. }
  85. if(m.middleRows(i,h).coeff(r,c) != Scalar(0))
  86. {
  87. VERIFY_IS_APPROX(m.middleRows(i,h).coeff(r,c), refMat.middleRows(i,h).coeff(r,c));
  88. }
  89. }
  90. }
  91. for(Index c=0; c<w; c++)
  92. {
  93. VERIFY_IS_APPROX(m.middleCols(j,w).col(c), refMat.middleCols(j,w).col(c));
  94. VERIFY_IS_APPROX(m.middleRows(i,h).col(c), refMat.middleRows(i,h).col(c));
  95. }
  96. }
  97. for(Index c=0; c<cols; c++)
  98. {
  99. VERIFY_IS_APPROX(m.col(c) + m.col(c), (m + m).col(c));
  100. VERIFY_IS_APPROX(m.col(c) + m.col(c), refMat.col(c) + refMat.col(c));
  101. }
  102. for(Index r=0; r<rows; r++)
  103. {
  104. VERIFY_IS_APPROX(m.row(r) + m.row(r), (m + m).row(r));
  105. VERIFY_IS_APPROX(m.row(r) + m.row(r), refMat.row(r) + refMat.row(r));
  106. }
  107. }
  108. // test innerVector()
  109. {
  110. DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
  111. SparseMatrixType m2(rows, cols);
  112. initSparse<Scalar>(density, refMat2, m2);
  113. Index j0 = internal::random<Index>(0,outer-1);
  114. Index j1 = internal::random<Index>(0,outer-1);
  115. Index r0 = internal::random<Index>(0,rows-1);
  116. Index c0 = internal::random<Index>(0,cols-1);
  117. VERIFY_IS_APPROX(m2.innerVector(j0), innervec(refMat2,j0));
  118. VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), innervec(refMat2,j0)+innervec(refMat2,j1));
  119. m2.innerVector(j0) *= Scalar(2);
  120. innervec(refMat2,j0) *= Scalar(2);
  121. VERIFY_IS_APPROX(m2, refMat2);
  122. m2.row(r0) *= Scalar(3);
  123. refMat2.row(r0) *= Scalar(3);
  124. VERIFY_IS_APPROX(m2, refMat2);
  125. m2.col(c0) *= Scalar(4);
  126. refMat2.col(c0) *= Scalar(4);
  127. VERIFY_IS_APPROX(m2, refMat2);
  128. m2.row(r0) /= Scalar(3);
  129. refMat2.row(r0) /= Scalar(3);
  130. VERIFY_IS_APPROX(m2, refMat2);
  131. m2.col(c0) /= Scalar(4);
  132. refMat2.col(c0) /= Scalar(4);
  133. VERIFY_IS_APPROX(m2, refMat2);
  134. SparseVectorType v1;
  135. VERIFY_IS_APPROX(v1 = m2.col(c0) * 4, refMat2.col(c0)*4);
  136. VERIFY_IS_APPROX(v1 = m2.row(r0) * 4, refMat2.row(r0).transpose()*4);
  137. SparseMatrixType m3(rows,cols);
  138. m3.reserve(VectorXi::Constant(outer,int(inner/2)));
  139. for(Index j=0; j<outer; ++j)
  140. for(Index k=0; k<(std::min)(j,inner); ++k)
  141. m3.insertByOuterInner(j,k) = internal::convert_index<StorageIndex>(k+1);
  142. for(Index j=0; j<(std::min)(outer, inner); ++j)
  143. {
  144. VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
  145. if(j>0)
  146. VERIFY(RealScalar(j)==numext::real(m3.innerVector(j).lastCoeff()));
  147. }
  148. m3.makeCompressed();
  149. for(Index j=0; j<(std::min)(outer, inner); ++j)
  150. {
  151. VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
  152. if(j>0)
  153. VERIFY(RealScalar(j)==numext::real(m3.innerVector(j).lastCoeff()));
  154. }
  155. VERIFY(m3.innerVector(j0).nonZeros() == m3.transpose().innerVector(j0).nonZeros());
  156. // m2.innerVector(j0) = 2*m2.innerVector(j1);
  157. // refMat2.col(j0) = 2*refMat2.col(j1);
  158. // VERIFY_IS_APPROX(m2, refMat2);
  159. }
  160. // test innerVectors()
  161. {
  162. DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
  163. SparseMatrixType m2(rows, cols);
  164. initSparse<Scalar>(density, refMat2, m2);
  165. if(internal::random<float>(0,1)>0.5f) m2.makeCompressed();
  166. Index j0 = internal::random<Index>(0,outer-2);
  167. Index j1 = internal::random<Index>(0,outer-2);
  168. Index n0 = internal::random<Index>(1,outer-(std::max)(j0,j1));
  169. if(SparseMatrixType::IsRowMajor)
  170. VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(j0,0,n0,cols));
  171. else
  172. VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0));
  173. if(SparseMatrixType::IsRowMajor)
  174. VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
  175. refMat2.middleRows(j0,n0)+refMat2.middleRows(j1,n0));
  176. else
  177. VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
  178. refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
  179. VERIFY_IS_APPROX(m2, refMat2);
  180. VERIFY(m2.innerVectors(j0,n0).nonZeros() == m2.transpose().innerVectors(j0,n0).nonZeros());
  181. m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0);
  182. if(SparseMatrixType::IsRowMajor)
  183. refMat2.middleRows(j0,n0) = (refMat2.middleRows(j0,n0) + refMat2.middleRows(j1,n0)).eval();
  184. else
  185. refMat2.middleCols(j0,n0) = (refMat2.middleCols(j0,n0) + refMat2.middleCols(j1,n0)).eval();
  186. VERIFY_IS_APPROX(m2, refMat2);
  187. }
  188. // test generic blocks
  189. {
  190. DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
  191. SparseMatrixType m2(rows, cols);
  192. initSparse<Scalar>(density, refMat2, m2);
  193. Index j0 = internal::random<Index>(0,outer-2);
  194. Index j1 = internal::random<Index>(0,outer-2);
  195. Index n0 = internal::random<Index>(1,outer-(std::max)(j0,j1));
  196. if(SparseMatrixType::IsRowMajor)
  197. VERIFY_IS_APPROX(m2.block(j0,0,n0,cols), refMat2.block(j0,0,n0,cols));
  198. else
  199. VERIFY_IS_APPROX(m2.block(0,j0,rows,n0), refMat2.block(0,j0,rows,n0));
  200. if(SparseMatrixType::IsRowMajor)
  201. VERIFY_IS_APPROX(m2.block(j0,0,n0,cols)+m2.block(j1,0,n0,cols),
  202. refMat2.block(j0,0,n0,cols)+refMat2.block(j1,0,n0,cols));
  203. else
  204. VERIFY_IS_APPROX(m2.block(0,j0,rows,n0)+m2.block(0,j1,rows,n0),
  205. refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
  206. Index i = internal::random<Index>(0,m2.outerSize()-1);
  207. if(SparseMatrixType::IsRowMajor) {
  208. m2.innerVector(i) = m2.innerVector(i) * s1;
  209. refMat2.row(i) = refMat2.row(i) * s1;
  210. VERIFY_IS_APPROX(m2,refMat2);
  211. } else {
  212. m2.innerVector(i) = m2.innerVector(i) * s1;
  213. refMat2.col(i) = refMat2.col(i) * s1;
  214. VERIFY_IS_APPROX(m2,refMat2);
  215. }
  216. Index r0 = internal::random<Index>(0,rows-2);
  217. Index c0 = internal::random<Index>(0,cols-2);
  218. Index r1 = internal::random<Index>(1,rows-r0);
  219. Index c1 = internal::random<Index>(1,cols-c0);
  220. VERIFY_IS_APPROX(DenseVector(m2.col(c0)), refMat2.col(c0));
  221. VERIFY_IS_APPROX(m2.col(c0), refMat2.col(c0));
  222. VERIFY_IS_APPROX(RowDenseVector(m2.row(r0)), refMat2.row(r0));
  223. VERIFY_IS_APPROX(m2.row(r0), refMat2.row(r0));
  224. VERIFY_IS_APPROX(m2.block(r0,c0,r1,c1), refMat2.block(r0,c0,r1,c1));
  225. VERIFY_IS_APPROX((2*m2).block(r0,c0,r1,c1), (2*refMat2).block(r0,c0,r1,c1));
  226. if(m2.nonZeros()>0)
  227. {
  228. VERIFY_IS_APPROX(m2, refMat2);
  229. SparseMatrixType m3(rows, cols);
  230. DenseMatrix refMat3(rows, cols); refMat3.setZero();
  231. Index n = internal::random<Index>(1,10);
  232. for(Index k=0; k<n; ++k)
  233. {
  234. Index o1 = internal::random<Index>(0,outer-1);
  235. Index o2 = internal::random<Index>(0,outer-1);
  236. if(SparseMatrixType::IsRowMajor)
  237. {
  238. m3.innerVector(o1) = m2.row(o2);
  239. refMat3.row(o1) = refMat2.row(o2);
  240. }
  241. else
  242. {
  243. m3.innerVector(o1) = m2.col(o2);
  244. refMat3.col(o1) = refMat2.col(o2);
  245. }
  246. if(internal::random<bool>())
  247. m3.makeCompressed();
  248. }
  249. if(m3.nonZeros()>0)
  250. VERIFY_IS_APPROX(m3, refMat3);
  251. }
  252. }
  253. }
  254. EIGEN_DECLARE_TEST(sparse_block)
  255. {
  256. for(int i = 0; i < g_repeat; i++) {
  257. int r = Eigen::internal::random<int>(1,200), c = Eigen::internal::random<int>(1,200);
  258. if(Eigen::internal::random<int>(0,4) == 0) {
  259. r = c; // check square matrices in 25% of tries
  260. }
  261. EIGEN_UNUSED_VARIABLE(r+c);
  262. CALL_SUBTEST_1(( sparse_block(SparseMatrix<double>(1, 1)) ));
  263. CALL_SUBTEST_1(( sparse_block(SparseMatrix<double>(8, 8)) ));
  264. CALL_SUBTEST_1(( sparse_block(SparseMatrix<double>(r, c)) ));
  265. CALL_SUBTEST_2(( sparse_block(SparseMatrix<std::complex<double>, ColMajor>(r, c)) ));
  266. CALL_SUBTEST_2(( sparse_block(SparseMatrix<std::complex<double>, RowMajor>(r, c)) ));
  267. CALL_SUBTEST_3(( sparse_block(SparseMatrix<double,ColMajor,long int>(r, c)) ));
  268. CALL_SUBTEST_3(( sparse_block(SparseMatrix<double,RowMajor,long int>(r, c)) ));
  269. r = Eigen::internal::random<int>(1,100);
  270. c = Eigen::internal::random<int>(1,100);
  271. if(Eigen::internal::random<int>(0,4) == 0) {
  272. r = c; // check square matrices in 25% of tries
  273. }
  274. CALL_SUBTEST_4(( sparse_block(SparseMatrix<double,ColMajor,short int>(short(r), short(c))) ));
  275. CALL_SUBTEST_4(( sparse_block(SparseMatrix<double,RowMajor,short int>(short(r), short(c))) ));
  276. #ifndef EIGEN_TEST_ANNOYING_SCALAR_DONT_THROW
  277. AnnoyingScalar::dont_throw = true;
  278. #endif
  279. CALL_SUBTEST_5(( sparse_block(SparseMatrix<AnnoyingScalar>(r,c)) ));
  280. }
  281. }