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.

254 lines
9.8 KiB

  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. template<typename SparseMatrixType> void sparse_block(const SparseMatrixType& ref)
  11. {
  12. const Index rows = ref.rows();
  13. const Index cols = ref.cols();
  14. const Index inner = ref.innerSize();
  15. const Index outer = ref.outerSize();
  16. typedef typename SparseMatrixType::Scalar Scalar;
  17. double density = (std::max)(8./(rows*cols), 0.01);
  18. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  19. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  20. typedef Matrix<Scalar,1,Dynamic> RowDenseVector;
  21. Scalar s1 = internal::random<Scalar>();
  22. {
  23. SparseMatrixType m(rows, cols);
  24. DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
  25. initSparse<Scalar>(density, refMat, m);
  26. VERIFY_IS_APPROX(m, refMat);
  27. // test InnerIterators and Block expressions
  28. for (int t=0; t<10; ++t)
  29. {
  30. Index j = internal::random<Index>(0,cols-2);
  31. Index i = internal::random<Index>(0,rows-2);
  32. Index w = internal::random<Index>(1,cols-j);
  33. Index h = internal::random<Index>(1,rows-i);
  34. VERIFY_IS_APPROX(m.block(i,j,h,w), refMat.block(i,j,h,w));
  35. for(Index c=0; c<w; c++)
  36. {
  37. VERIFY_IS_APPROX(m.block(i,j,h,w).col(c), refMat.block(i,j,h,w).col(c));
  38. for(Index r=0; r<h; r++)
  39. {
  40. VERIFY_IS_APPROX(m.block(i,j,h,w).col(c).coeff(r), refMat.block(i,j,h,w).col(c).coeff(r));
  41. VERIFY_IS_APPROX(m.block(i,j,h,w).coeff(r,c), refMat.block(i,j,h,w).coeff(r,c));
  42. }
  43. }
  44. for(Index r=0; r<h; r++)
  45. {
  46. VERIFY_IS_APPROX(m.block(i,j,h,w).row(r), refMat.block(i,j,h,w).row(r));
  47. for(Index c=0; c<w; c++)
  48. {
  49. VERIFY_IS_APPROX(m.block(i,j,h,w).row(r).coeff(c), refMat.block(i,j,h,w).row(r).coeff(c));
  50. VERIFY_IS_APPROX(m.block(i,j,h,w).coeff(r,c), refMat.block(i,j,h,w).coeff(r,c));
  51. }
  52. }
  53. VERIFY_IS_APPROX(m.middleCols(j,w), refMat.middleCols(j,w));
  54. VERIFY_IS_APPROX(m.middleRows(i,h), refMat.middleRows(i,h));
  55. for(Index r=0; r<h; r++)
  56. {
  57. VERIFY_IS_APPROX(m.middleCols(j,w).row(r), refMat.middleCols(j,w).row(r));
  58. VERIFY_IS_APPROX(m.middleRows(i,h).row(r), refMat.middleRows(i,h).row(r));
  59. for(Index c=0; c<w; c++)
  60. {
  61. VERIFY_IS_APPROX(m.col(c).coeff(r), refMat.col(c).coeff(r));
  62. VERIFY_IS_APPROX(m.row(r).coeff(c), refMat.row(r).coeff(c));
  63. VERIFY_IS_APPROX(m.middleCols(j,w).coeff(r,c), refMat.middleCols(j,w).coeff(r,c));
  64. VERIFY_IS_APPROX(m.middleRows(i,h).coeff(r,c), refMat.middleRows(i,h).coeff(r,c));
  65. if(m.middleCols(j,w).coeff(r,c) != Scalar(0))
  66. {
  67. VERIFY_IS_APPROX(m.middleCols(j,w).coeffRef(r,c), refMat.middleCols(j,w).coeff(r,c));
  68. }
  69. if(m.middleRows(i,h).coeff(r,c) != Scalar(0))
  70. {
  71. VERIFY_IS_APPROX(m.middleRows(i,h).coeff(r,c), refMat.middleRows(i,h).coeff(r,c));
  72. }
  73. }
  74. }
  75. for(Index c=0; c<w; c++)
  76. {
  77. VERIFY_IS_APPROX(m.middleCols(j,w).col(c), refMat.middleCols(j,w).col(c));
  78. VERIFY_IS_APPROX(m.middleRows(i,h).col(c), refMat.middleRows(i,h).col(c));
  79. }
  80. }
  81. for(Index c=0; c<cols; c++)
  82. {
  83. VERIFY_IS_APPROX(m.col(c) + m.col(c), (m + m).col(c));
  84. VERIFY_IS_APPROX(m.col(c) + m.col(c), refMat.col(c) + refMat.col(c));
  85. }
  86. for(Index r=0; r<rows; r++)
  87. {
  88. VERIFY_IS_APPROX(m.row(r) + m.row(r), (m + m).row(r));
  89. VERIFY_IS_APPROX(m.row(r) + m.row(r), refMat.row(r) + refMat.row(r));
  90. }
  91. }
  92. // test innerVector()
  93. {
  94. DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
  95. SparseMatrixType m2(rows, cols);
  96. initSparse<Scalar>(density, refMat2, m2);
  97. Index j0 = internal::random<Index>(0,outer-1);
  98. Index j1 = internal::random<Index>(0,outer-1);
  99. if(SparseMatrixType::IsRowMajor)
  100. VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.row(j0));
  101. else
  102. VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.col(j0));
  103. if(SparseMatrixType::IsRowMajor)
  104. VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.row(j0)+refMat2.row(j1));
  105. else
  106. VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1));
  107. SparseMatrixType m3(rows,cols);
  108. m3.reserve(VectorXi::Constant(outer,int(inner/2)));
  109. for(Index j=0; j<outer; ++j)
  110. for(Index k=0; k<(std::min)(j,inner); ++k)
  111. m3.insertByOuterInner(j,k) = k+1;
  112. for(Index j=0; j<(std::min)(outer, inner); ++j)
  113. {
  114. VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
  115. if(j>0)
  116. VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
  117. }
  118. m3.makeCompressed();
  119. for(Index j=0; j<(std::min)(outer, inner); ++j)
  120. {
  121. VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
  122. if(j>0)
  123. VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
  124. }
  125. VERIFY(m3.innerVector(j0).nonZeros() == m3.transpose().innerVector(j0).nonZeros());
  126. // m2.innerVector(j0) = 2*m2.innerVector(j1);
  127. // refMat2.col(j0) = 2*refMat2.col(j1);
  128. // VERIFY_IS_APPROX(m2, refMat2);
  129. }
  130. // test innerVectors()
  131. {
  132. DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
  133. SparseMatrixType m2(rows, cols);
  134. initSparse<Scalar>(density, refMat2, m2);
  135. if(internal::random<float>(0,1)>0.5) m2.makeCompressed();
  136. Index j0 = internal::random<Index>(0,outer-2);
  137. Index j1 = internal::random<Index>(0,outer-2);
  138. Index n0 = internal::random<Index>(1,outer-(std::max)(j0,j1));
  139. if(SparseMatrixType::IsRowMajor)
  140. VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(j0,0,n0,cols));
  141. else
  142. VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0));
  143. if(SparseMatrixType::IsRowMajor)
  144. VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
  145. refMat2.middleRows(j0,n0)+refMat2.middleRows(j1,n0));
  146. else
  147. VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
  148. refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
  149. VERIFY_IS_APPROX(m2, refMat2);
  150. VERIFY(m2.innerVectors(j0,n0).nonZeros() == m2.transpose().innerVectors(j0,n0).nonZeros());
  151. m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0);
  152. if(SparseMatrixType::IsRowMajor)
  153. refMat2.middleRows(j0,n0) = (refMat2.middleRows(j0,n0) + refMat2.middleRows(j1,n0)).eval();
  154. else
  155. refMat2.middleCols(j0,n0) = (refMat2.middleCols(j0,n0) + refMat2.middleCols(j1,n0)).eval();
  156. VERIFY_IS_APPROX(m2, refMat2);
  157. }
  158. // test generic blocks
  159. {
  160. DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
  161. SparseMatrixType m2(rows, cols);
  162. initSparse<Scalar>(density, refMat2, m2);
  163. Index j0 = internal::random<Index>(0,outer-2);
  164. Index j1 = internal::random<Index>(0,outer-2);
  165. Index n0 = internal::random<Index>(1,outer-(std::max)(j0,j1));
  166. if(SparseMatrixType::IsRowMajor)
  167. VERIFY_IS_APPROX(m2.block(j0,0,n0,cols), refMat2.block(j0,0,n0,cols));
  168. else
  169. VERIFY_IS_APPROX(m2.block(0,j0,rows,n0), refMat2.block(0,j0,rows,n0));
  170. if(SparseMatrixType::IsRowMajor)
  171. VERIFY_IS_APPROX(m2.block(j0,0,n0,cols)+m2.block(j1,0,n0,cols),
  172. refMat2.block(j0,0,n0,cols)+refMat2.block(j1,0,n0,cols));
  173. else
  174. VERIFY_IS_APPROX(m2.block(0,j0,rows,n0)+m2.block(0,j1,rows,n0),
  175. refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
  176. Index i = internal::random<Index>(0,m2.outerSize()-1);
  177. if(SparseMatrixType::IsRowMajor) {
  178. m2.innerVector(i) = m2.innerVector(i) * s1;
  179. refMat2.row(i) = refMat2.row(i) * s1;
  180. VERIFY_IS_APPROX(m2,refMat2);
  181. } else {
  182. m2.innerVector(i) = m2.innerVector(i) * s1;
  183. refMat2.col(i) = refMat2.col(i) * s1;
  184. VERIFY_IS_APPROX(m2,refMat2);
  185. }
  186. Index r0 = internal::random<Index>(0,rows-2);
  187. Index c0 = internal::random<Index>(0,cols-2);
  188. Index r1 = internal::random<Index>(1,rows-r0);
  189. Index c1 = internal::random<Index>(1,cols-c0);
  190. VERIFY_IS_APPROX(DenseVector(m2.col(c0)), refMat2.col(c0));
  191. VERIFY_IS_APPROX(m2.col(c0), refMat2.col(c0));
  192. VERIFY_IS_APPROX(RowDenseVector(m2.row(r0)), refMat2.row(r0));
  193. VERIFY_IS_APPROX(m2.row(r0), refMat2.row(r0));
  194. VERIFY_IS_APPROX(m2.block(r0,c0,r1,c1), refMat2.block(r0,c0,r1,c1));
  195. VERIFY_IS_APPROX((2*m2).block(r0,c0,r1,c1), (2*refMat2).block(r0,c0,r1,c1));
  196. }
  197. }
  198. void test_sparse_block()
  199. {
  200. for(int i = 0; i < g_repeat; i++) {
  201. int r = StormEigen::internal::random<int>(1,200), c = StormEigen::internal::random<int>(1,200);
  202. if(StormEigen::internal::random<int>(0,4) == 0) {
  203. r = c; // check square matrices in 25% of tries
  204. }
  205. EIGEN_UNUSED_VARIABLE(r+c);
  206. CALL_SUBTEST_1(( sparse_block(SparseMatrix<double>(1, 1)) ));
  207. CALL_SUBTEST_1(( sparse_block(SparseMatrix<double>(8, 8)) ));
  208. CALL_SUBTEST_1(( sparse_block(SparseMatrix<double>(r, c)) ));
  209. CALL_SUBTEST_2(( sparse_block(SparseMatrix<std::complex<double>, ColMajor>(r, c)) ));
  210. CALL_SUBTEST_2(( sparse_block(SparseMatrix<std::complex<double>, RowMajor>(r, c)) ));
  211. CALL_SUBTEST_3(( sparse_block(SparseMatrix<double,ColMajor,long int>(r, c)) ));
  212. CALL_SUBTEST_3(( sparse_block(SparseMatrix<double,RowMajor,long int>(r, c)) ));
  213. r = StormEigen::internal::random<int>(1,100);
  214. c = StormEigen::internal::random<int>(1,100);
  215. if(StormEigen::internal::random<int>(0,4) == 0) {
  216. r = c; // check square matrices in 25% of tries
  217. }
  218. CALL_SUBTEST_4(( sparse_block(SparseMatrix<double,ColMajor,short int>(short(r), short(c))) ));
  219. CALL_SUBTEST_4(( sparse_block(SparseMatrix<double,RowMajor,short int>(short(r), short(c))) ));
  220. }
  221. }