sparse_product.cpp 25 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475
  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. //
  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. #if defined(_MSC_VER) && (_MSC_VER==1800)
  10. // This unit test takes forever to compile in Release mode with MSVC 2013,
  11. // multiple hours. So let's switch off optimization for this one.
  12. #pragma optimize("",off)
  13. #endif
  14. static long int nb_temporaries;
  15. inline void on_temporary_creation() {
  16. // here's a great place to set a breakpoint when debugging failures in this test!
  17. nb_temporaries++;
  18. }
  19. #define EIGEN_SPARSE_CREATE_TEMPORARY_PLUGIN { on_temporary_creation(); }
  20. #include "sparse.h"
  21. #define VERIFY_EVALUATION_COUNT(XPR,N) {\
  22. nb_temporaries = 0; \
  23. CALL_SUBTEST( XPR ); \
  24. if(nb_temporaries!=N) std::cerr << "nb_temporaries == " << nb_temporaries << "\n"; \
  25. VERIFY( (#XPR) && nb_temporaries==N ); \
  26. }
  27. template<typename SparseMatrixType> void sparse_product()
  28. {
  29. typedef typename SparseMatrixType::StorageIndex StorageIndex;
  30. Index n = 100;
  31. const Index rows = internal::random<Index>(1,n);
  32. const Index cols = internal::random<Index>(1,n);
  33. const Index depth = internal::random<Index>(1,n);
  34. typedef typename SparseMatrixType::Scalar Scalar;
  35. enum { Flags = SparseMatrixType::Flags };
  36. double density = (std::max)(8./(rows*cols), 0.2);
  37. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  38. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  39. typedef Matrix<Scalar,1,Dynamic> RowDenseVector;
  40. typedef SparseVector<Scalar,0,StorageIndex> ColSpVector;
  41. typedef SparseVector<Scalar,RowMajor,StorageIndex> RowSpVector;
  42. Scalar s1 = internal::random<Scalar>();
  43. Scalar s2 = internal::random<Scalar>();
  44. // test matrix-matrix product
  45. {
  46. DenseMatrix refMat2 = DenseMatrix::Zero(rows, depth);
  47. DenseMatrix refMat2t = DenseMatrix::Zero(depth, rows);
  48. DenseMatrix refMat3 = DenseMatrix::Zero(depth, cols);
  49. DenseMatrix refMat3t = DenseMatrix::Zero(cols, depth);
  50. DenseMatrix refMat4 = DenseMatrix::Zero(rows, cols);
  51. DenseMatrix refMat4t = DenseMatrix::Zero(cols, rows);
  52. DenseMatrix refMat5 = DenseMatrix::Random(depth, cols);
  53. DenseMatrix refMat6 = DenseMatrix::Random(rows, rows);
  54. DenseMatrix dm4 = DenseMatrix::Zero(rows, rows);
  55. // DenseVector dv1 = DenseVector::Random(rows);
  56. SparseMatrixType m2 (rows, depth);
  57. SparseMatrixType m2t(depth, rows);
  58. SparseMatrixType m3 (depth, cols);
  59. SparseMatrixType m3t(cols, depth);
  60. SparseMatrixType m4 (rows, cols);
  61. SparseMatrixType m4t(cols, rows);
  62. SparseMatrixType m6(rows, rows);
  63. initSparse(density, refMat2, m2);
  64. initSparse(density, refMat2t, m2t);
  65. initSparse(density, refMat3, m3);
  66. initSparse(density, refMat3t, m3t);
  67. initSparse(density, refMat4, m4);
  68. initSparse(density, refMat4t, m4t);
  69. initSparse(density, refMat6, m6);
  70. // int c = internal::random<int>(0,depth-1);
  71. // sparse * sparse
  72. VERIFY_IS_APPROX(m4=m2*m3, refMat4=refMat2*refMat3);
  73. VERIFY_IS_APPROX(m4=m2t.transpose()*m3, refMat4=refMat2t.transpose()*refMat3);
  74. VERIFY_IS_APPROX(m4=m2t.transpose()*m3t.transpose(), refMat4=refMat2t.transpose()*refMat3t.transpose());
  75. VERIFY_IS_APPROX(m4=m2*m3t.transpose(), refMat4=refMat2*refMat3t.transpose());
  76. VERIFY_IS_APPROX(m4 = m2*m3/s1, refMat4 = refMat2*refMat3/s1);
  77. VERIFY_IS_APPROX(m4 = m2*m3*s1, refMat4 = refMat2*refMat3*s1);
  78. VERIFY_IS_APPROX(m4 = s2*m2*m3*s1, refMat4 = s2*refMat2*refMat3*s1);
  79. VERIFY_IS_APPROX(m4 = (m2+m2)*m3, refMat4 = (refMat2+refMat2)*refMat3);
  80. VERIFY_IS_APPROX(m4 = m2*m3.leftCols(cols/2), refMat4 = refMat2*refMat3.leftCols(cols/2));
  81. VERIFY_IS_APPROX(m4 = m2*(m3+m3).leftCols(cols/2), refMat4 = refMat2*(refMat3+refMat3).leftCols(cols/2));
  82. VERIFY_IS_APPROX(m4=(m2*m3).pruned(0), refMat4=refMat2*refMat3);
  83. VERIFY_IS_APPROX(m4=(m2t.transpose()*m3).pruned(0), refMat4=refMat2t.transpose()*refMat3);
  84. VERIFY_IS_APPROX(m4=(m2t.transpose()*m3t.transpose()).pruned(0), refMat4=refMat2t.transpose()*refMat3t.transpose());
  85. VERIFY_IS_APPROX(m4=(m2*m3t.transpose()).pruned(0), refMat4=refMat2*refMat3t.transpose());
  86. // make sure the right product implementation is called:
  87. if((!SparseMatrixType::IsRowMajor) && m2.rows()<=m3.cols())
  88. {
  89. VERIFY_EVALUATION_COUNT(m4 = m2*m3, 3); // 1 temp for the result + 2 for transposing and get a sorted result.
  90. VERIFY_EVALUATION_COUNT(m4 = (m2*m3).pruned(0), 1);
  91. VERIFY_EVALUATION_COUNT(m4 = (m2*m3).eval().pruned(0), 4);
  92. }
  93. // and that pruning is effective:
  94. {
  95. DenseMatrix Ad(2,2);
  96. Ad << -1, 1, 1, 1;
  97. SparseMatrixType As(Ad.sparseView()), B(2,2);
  98. VERIFY_IS_EQUAL( (As*As.transpose()).eval().nonZeros(), 4);
  99. VERIFY_IS_EQUAL( (Ad*Ad.transpose()).eval().sparseView().eval().nonZeros(), 2);
  100. VERIFY_IS_EQUAL( (As*As.transpose()).pruned(1e-6).eval().nonZeros(), 2);
  101. }
  102. // dense ?= sparse * sparse
  103. VERIFY_IS_APPROX(dm4 =m2*m3, refMat4 =refMat2*refMat3);
  104. VERIFY_IS_APPROX(dm4+=m2*m3, refMat4+=refMat2*refMat3);
  105. VERIFY_IS_APPROX(dm4-=m2*m3, refMat4-=refMat2*refMat3);
  106. VERIFY_IS_APPROX(dm4 =m2t.transpose()*m3, refMat4 =refMat2t.transpose()*refMat3);
  107. VERIFY_IS_APPROX(dm4+=m2t.transpose()*m3, refMat4+=refMat2t.transpose()*refMat3);
  108. VERIFY_IS_APPROX(dm4-=m2t.transpose()*m3, refMat4-=refMat2t.transpose()*refMat3);
  109. VERIFY_IS_APPROX(dm4 =m2t.transpose()*m3t.transpose(), refMat4 =refMat2t.transpose()*refMat3t.transpose());
  110. VERIFY_IS_APPROX(dm4+=m2t.transpose()*m3t.transpose(), refMat4+=refMat2t.transpose()*refMat3t.transpose());
  111. VERIFY_IS_APPROX(dm4-=m2t.transpose()*m3t.transpose(), refMat4-=refMat2t.transpose()*refMat3t.transpose());
  112. VERIFY_IS_APPROX(dm4 =m2*m3t.transpose(), refMat4 =refMat2*refMat3t.transpose());
  113. VERIFY_IS_APPROX(dm4+=m2*m3t.transpose(), refMat4+=refMat2*refMat3t.transpose());
  114. VERIFY_IS_APPROX(dm4-=m2*m3t.transpose(), refMat4-=refMat2*refMat3t.transpose());
  115. VERIFY_IS_APPROX(dm4 = m2*m3*s1, refMat4 = refMat2*refMat3*s1);
  116. // test aliasing
  117. m4 = m2; refMat4 = refMat2;
  118. VERIFY_IS_APPROX(m4=m4*m3, refMat4=refMat4*refMat3);
  119. // sparse * dense matrix
  120. VERIFY_IS_APPROX(dm4=m2*refMat3, refMat4=refMat2*refMat3);
  121. VERIFY_IS_APPROX(dm4=m2*refMat3t.transpose(), refMat4=refMat2*refMat3t.transpose());
  122. VERIFY_IS_APPROX(dm4=m2t.transpose()*refMat3, refMat4=refMat2t.transpose()*refMat3);
  123. VERIFY_IS_APPROX(dm4=m2t.transpose()*refMat3t.transpose(), refMat4=refMat2t.transpose()*refMat3t.transpose());
  124. VERIFY_IS_APPROX(dm4=m2*refMat3, refMat4=refMat2*refMat3);
  125. VERIFY_IS_APPROX(dm4=dm4+m2*refMat3, refMat4=refMat4+refMat2*refMat3);
  126. VERIFY_IS_APPROX(dm4+=m2*refMat3, refMat4+=refMat2*refMat3);
  127. VERIFY_IS_APPROX(dm4-=m2*refMat3, refMat4-=refMat2*refMat3);
  128. VERIFY_IS_APPROX(dm4.noalias()+=m2*refMat3, refMat4+=refMat2*refMat3);
  129. VERIFY_IS_APPROX(dm4.noalias()-=m2*refMat3, refMat4-=refMat2*refMat3);
  130. VERIFY_IS_APPROX(dm4=m2*(refMat3+refMat3), refMat4=refMat2*(refMat3+refMat3));
  131. VERIFY_IS_APPROX(dm4=m2t.transpose()*(refMat3+refMat5)*0.5, refMat4=refMat2t.transpose()*(refMat3+refMat5)*0.5);
  132. // sparse * dense vector
  133. VERIFY_IS_APPROX(dm4.col(0)=m2*refMat3.col(0), refMat4.col(0)=refMat2*refMat3.col(0));
  134. VERIFY_IS_APPROX(dm4.col(0)=m2*refMat3t.transpose().col(0), refMat4.col(0)=refMat2*refMat3t.transpose().col(0));
  135. VERIFY_IS_APPROX(dm4.col(0)=m2t.transpose()*refMat3.col(0), refMat4.col(0)=refMat2t.transpose()*refMat3.col(0));
  136. VERIFY_IS_APPROX(dm4.col(0)=m2t.transpose()*refMat3t.transpose().col(0), refMat4.col(0)=refMat2t.transpose()*refMat3t.transpose().col(0));
  137. // dense * sparse
  138. VERIFY_IS_APPROX(dm4=refMat2*m3, refMat4=refMat2*refMat3);
  139. VERIFY_IS_APPROX(dm4=dm4+refMat2*m3, refMat4=refMat4+refMat2*refMat3);
  140. VERIFY_IS_APPROX(dm4+=refMat2*m3, refMat4+=refMat2*refMat3);
  141. VERIFY_IS_APPROX(dm4-=refMat2*m3, refMat4-=refMat2*refMat3);
  142. VERIFY_IS_APPROX(dm4.noalias()+=refMat2*m3, refMat4+=refMat2*refMat3);
  143. VERIFY_IS_APPROX(dm4.noalias()-=refMat2*m3, refMat4-=refMat2*refMat3);
  144. VERIFY_IS_APPROX(dm4=refMat2*m3t.transpose(), refMat4=refMat2*refMat3t.transpose());
  145. VERIFY_IS_APPROX(dm4=refMat2t.transpose()*m3, refMat4=refMat2t.transpose()*refMat3);
  146. VERIFY_IS_APPROX(dm4=refMat2t.transpose()*m3t.transpose(), refMat4=refMat2t.transpose()*refMat3t.transpose());
  147. // sparse * dense and dense * sparse outer product
  148. {
  149. Index c = internal::random<Index>(0,depth-1);
  150. Index r = internal::random<Index>(0,rows-1);
  151. Index c1 = internal::random<Index>(0,cols-1);
  152. Index r1 = internal::random<Index>(0,depth-1);
  153. DenseMatrix dm5 = DenseMatrix::Random(depth, cols);
  154. VERIFY_IS_APPROX( m4=m2.col(c)*dm5.col(c1).transpose(), refMat4=refMat2.col(c)*dm5.col(c1).transpose());
  155. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  156. VERIFY_IS_APPROX( m4=m2.middleCols(c,1)*dm5.col(c1).transpose(), refMat4=refMat2.col(c)*dm5.col(c1).transpose());
  157. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  158. VERIFY_IS_APPROX(dm4=m2.col(c)*dm5.col(c1).transpose(), refMat4=refMat2.col(c)*dm5.col(c1).transpose());
  159. VERIFY_IS_APPROX(m4=dm5.col(c1)*m2.col(c).transpose(), refMat4=dm5.col(c1)*refMat2.col(c).transpose());
  160. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  161. VERIFY_IS_APPROX(m4=dm5.col(c1)*m2.middleCols(c,1).transpose(), refMat4=dm5.col(c1)*refMat2.col(c).transpose());
  162. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  163. VERIFY_IS_APPROX(dm4=dm5.col(c1)*m2.col(c).transpose(), refMat4=dm5.col(c1)*refMat2.col(c).transpose());
  164. VERIFY_IS_APPROX( m4=dm5.row(r1).transpose()*m2.col(c).transpose(), refMat4=dm5.row(r1).transpose()*refMat2.col(c).transpose());
  165. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  166. VERIFY_IS_APPROX(dm4=dm5.row(r1).transpose()*m2.col(c).transpose(), refMat4=dm5.row(r1).transpose()*refMat2.col(c).transpose());
  167. VERIFY_IS_APPROX( m4=m2.row(r).transpose()*dm5.col(c1).transpose(), refMat4=refMat2.row(r).transpose()*dm5.col(c1).transpose());
  168. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  169. VERIFY_IS_APPROX( m4=m2.middleRows(r,1).transpose()*dm5.col(c1).transpose(), refMat4=refMat2.row(r).transpose()*dm5.col(c1).transpose());
  170. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  171. VERIFY_IS_APPROX(dm4=m2.row(r).transpose()*dm5.col(c1).transpose(), refMat4=refMat2.row(r).transpose()*dm5.col(c1).transpose());
  172. VERIFY_IS_APPROX( m4=dm5.col(c1)*m2.row(r), refMat4=dm5.col(c1)*refMat2.row(r));
  173. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  174. VERIFY_IS_APPROX( m4=dm5.col(c1)*m2.middleRows(r,1), refMat4=dm5.col(c1)*refMat2.row(r));
  175. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  176. VERIFY_IS_APPROX(dm4=dm5.col(c1)*m2.row(r), refMat4=dm5.col(c1)*refMat2.row(r));
  177. VERIFY_IS_APPROX( m4=dm5.row(r1).transpose()*m2.row(r), refMat4=dm5.row(r1).transpose()*refMat2.row(r));
  178. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  179. VERIFY_IS_APPROX(dm4=dm5.row(r1).transpose()*m2.row(r), refMat4=dm5.row(r1).transpose()*refMat2.row(r));
  180. }
  181. VERIFY_IS_APPROX(m6=m6*m6, refMat6=refMat6*refMat6);
  182. // sparse matrix * sparse vector
  183. ColSpVector cv0(cols), cv1;
  184. DenseVector dcv0(cols), dcv1;
  185. initSparse(2*density,dcv0, cv0);
  186. RowSpVector rv0(depth), rv1;
  187. RowDenseVector drv0(depth), drv1(rv1);
  188. initSparse(2*density,drv0, rv0);
  189. VERIFY_IS_APPROX(cv1=m3*cv0, dcv1=refMat3*dcv0);
  190. VERIFY_IS_APPROX(rv1=rv0*m3, drv1=drv0*refMat3);
  191. VERIFY_IS_APPROX(cv1=m3t.adjoint()*cv0, dcv1=refMat3t.adjoint()*dcv0);
  192. VERIFY_IS_APPROX(cv1=rv0*m3, dcv1=drv0*refMat3);
  193. VERIFY_IS_APPROX(rv1=m3*cv0, drv1=refMat3*dcv0);
  194. }
  195. // test matrix - diagonal product
  196. {
  197. DenseMatrix refM2 = DenseMatrix::Zero(rows, cols);
  198. DenseMatrix refM3 = DenseMatrix::Zero(rows, cols);
  199. DenseMatrix d3 = DenseMatrix::Zero(rows, cols);
  200. DiagonalMatrix<Scalar,Dynamic> d1(DenseVector::Random(cols));
  201. DiagonalMatrix<Scalar,Dynamic> d2(DenseVector::Random(rows));
  202. SparseMatrixType m2(rows, cols);
  203. SparseMatrixType m3(rows, cols);
  204. initSparse<Scalar>(density, refM2, m2);
  205. initSparse<Scalar>(density, refM3, m3);
  206. VERIFY_IS_APPROX(m3=m2*d1, refM3=refM2*d1);
  207. VERIFY_IS_APPROX(m3=m2.transpose()*d2, refM3=refM2.transpose()*d2);
  208. VERIFY_IS_APPROX(m3=d2*m2, refM3=d2*refM2);
  209. VERIFY_IS_APPROX(m3=d1*m2.transpose(), refM3=d1*refM2.transpose());
  210. // also check with a SparseWrapper:
  211. DenseVector v1 = DenseVector::Random(cols);
  212. DenseVector v2 = DenseVector::Random(rows);
  213. DenseVector v3 = DenseVector::Random(rows);
  214. VERIFY_IS_APPROX(m3=m2*v1.asDiagonal(), refM3=refM2*v1.asDiagonal());
  215. VERIFY_IS_APPROX(m3=m2.transpose()*v2.asDiagonal(), refM3=refM2.transpose()*v2.asDiagonal());
  216. VERIFY_IS_APPROX(m3=v2.asDiagonal()*m2, refM3=v2.asDiagonal()*refM2);
  217. VERIFY_IS_APPROX(m3=v1.asDiagonal()*m2.transpose(), refM3=v1.asDiagonal()*refM2.transpose());
  218. VERIFY_IS_APPROX(m3=v2.asDiagonal()*m2*v1.asDiagonal(), refM3=v2.asDiagonal()*refM2*v1.asDiagonal());
  219. VERIFY_IS_APPROX(v2=m2*v1.asDiagonal()*v1, refM2*v1.asDiagonal()*v1);
  220. VERIFY_IS_APPROX(v3=v2.asDiagonal()*m2*v1, v2.asDiagonal()*refM2*v1);
  221. // evaluate to a dense matrix to check the .row() and .col() iterator functions
  222. VERIFY_IS_APPROX(d3=m2*d1, refM3=refM2*d1);
  223. VERIFY_IS_APPROX(d3=m2.transpose()*d2, refM3=refM2.transpose()*d2);
  224. VERIFY_IS_APPROX(d3=d2*m2, refM3=d2*refM2);
  225. VERIFY_IS_APPROX(d3=d1*m2.transpose(), refM3=d1*refM2.transpose());
  226. }
  227. // test self-adjoint and triangular-view products
  228. {
  229. DenseMatrix b = DenseMatrix::Random(rows, rows);
  230. DenseMatrix x = DenseMatrix::Random(rows, rows);
  231. DenseMatrix refX = DenseMatrix::Random(rows, rows);
  232. DenseMatrix refUp = DenseMatrix::Zero(rows, rows);
  233. DenseMatrix refLo = DenseMatrix::Zero(rows, rows);
  234. DenseMatrix refS = DenseMatrix::Zero(rows, rows);
  235. DenseMatrix refA = DenseMatrix::Zero(rows, rows);
  236. SparseMatrixType mUp(rows, rows);
  237. SparseMatrixType mLo(rows, rows);
  238. SparseMatrixType mS(rows, rows);
  239. SparseMatrixType mA(rows, rows);
  240. initSparse<Scalar>(density, refA, mA);
  241. do {
  242. initSparse<Scalar>(density, refUp, mUp, ForceRealDiag|/*ForceNonZeroDiag|*/MakeUpperTriangular);
  243. } while (refUp.isZero());
  244. refLo = refUp.adjoint();
  245. mLo = mUp.adjoint();
  246. refS = refUp + refLo;
  247. refS.diagonal() *= 0.5;
  248. mS = mUp + mLo;
  249. // TODO be able to address the diagonal....
  250. for (int k=0; k<mS.outerSize(); ++k)
  251. for (typename SparseMatrixType::InnerIterator it(mS,k); it; ++it)
  252. if (it.index() == k)
  253. it.valueRef() *= Scalar(0.5);
  254. VERIFY_IS_APPROX(refS.adjoint(), refS);
  255. VERIFY_IS_APPROX(mS.adjoint(), mS);
  256. VERIFY_IS_APPROX(mS, refS);
  257. VERIFY_IS_APPROX(x=mS*b, refX=refS*b);
  258. // sparse selfadjointView with dense matrices
  259. VERIFY_IS_APPROX(x=mUp.template selfadjointView<Upper>()*b, refX=refS*b);
  260. VERIFY_IS_APPROX(x=mLo.template selfadjointView<Lower>()*b, refX=refS*b);
  261. VERIFY_IS_APPROX(x=mS.template selfadjointView<Upper|Lower>()*b, refX=refS*b);
  262. VERIFY_IS_APPROX(x=b * mUp.template selfadjointView<Upper>(), refX=b*refS);
  263. VERIFY_IS_APPROX(x=b * mLo.template selfadjointView<Lower>(), refX=b*refS);
  264. VERIFY_IS_APPROX(x=b * mS.template selfadjointView<Upper|Lower>(), refX=b*refS);
  265. VERIFY_IS_APPROX(x.noalias()+=mUp.template selfadjointView<Upper>()*b, refX+=refS*b);
  266. VERIFY_IS_APPROX(x.noalias()-=mLo.template selfadjointView<Lower>()*b, refX-=refS*b);
  267. VERIFY_IS_APPROX(x.noalias()+=mS.template selfadjointView<Upper|Lower>()*b, refX+=refS*b);
  268. // sparse selfadjointView with sparse matrices
  269. SparseMatrixType mSres(rows,rows);
  270. VERIFY_IS_APPROX(mSres = mLo.template selfadjointView<Lower>()*mS,
  271. refX = refLo.template selfadjointView<Lower>()*refS);
  272. VERIFY_IS_APPROX(mSres = mS * mLo.template selfadjointView<Lower>(),
  273. refX = refS * refLo.template selfadjointView<Lower>());
  274. // sparse triangularView with dense matrices
  275. VERIFY_IS_APPROX(x=mA.template triangularView<Upper>()*b, refX=refA.template triangularView<Upper>()*b);
  276. VERIFY_IS_APPROX(x=mA.template triangularView<Lower>()*b, refX=refA.template triangularView<Lower>()*b);
  277. VERIFY_IS_APPROX(x=b*mA.template triangularView<Upper>(), refX=b*refA.template triangularView<Upper>());
  278. VERIFY_IS_APPROX(x=b*mA.template triangularView<Lower>(), refX=b*refA.template triangularView<Lower>());
  279. // sparse triangularView with sparse matrices
  280. VERIFY_IS_APPROX(mSres = mA.template triangularView<Lower>()*mS, refX = refA.template triangularView<Lower>()*refS);
  281. VERIFY_IS_APPROX(mSres = mS * mA.template triangularView<Lower>(), refX = refS * refA.template triangularView<Lower>());
  282. VERIFY_IS_APPROX(mSres = mA.template triangularView<Upper>()*mS, refX = refA.template triangularView<Upper>()*refS);
  283. VERIFY_IS_APPROX(mSres = mS * mA.template triangularView<Upper>(), refX = refS * refA.template triangularView<Upper>());
  284. }
  285. }
  286. // New test for Bug in SparseTimeDenseProduct
  287. template<typename SparseMatrixType, typename DenseMatrixType> void sparse_product_regression_test()
  288. {
  289. // This code does not compile with afflicted versions of the bug
  290. SparseMatrixType sm1(3,2);
  291. DenseMatrixType m2(2,2);
  292. sm1.setZero();
  293. m2.setZero();
  294. DenseMatrixType m3 = sm1*m2;
  295. // This code produces a segfault with afflicted versions of another SparseTimeDenseProduct
  296. // bug
  297. SparseMatrixType sm2(20000,2);
  298. sm2.setZero();
  299. DenseMatrixType m4(sm2*m2);
  300. VERIFY_IS_APPROX( m4(0,0), 0.0 );
  301. }
  302. template<typename Scalar>
  303. void bug_942()
  304. {
  305. typedef Matrix<Scalar, Dynamic, 1> Vector;
  306. typedef SparseMatrix<Scalar, ColMajor> ColSpMat;
  307. typedef SparseMatrix<Scalar, RowMajor> RowSpMat;
  308. ColSpMat cmA(1,1);
  309. cmA.insert(0,0) = 1;
  310. RowSpMat rmA(1,1);
  311. rmA.insert(0,0) = 1;
  312. Vector d(1);
  313. d[0] = 2;
  314. double res = 2;
  315. VERIFY_IS_APPROX( ( cmA*d.asDiagonal() ).eval().coeff(0,0), res );
  316. VERIFY_IS_APPROX( ( d.asDiagonal()*rmA ).eval().coeff(0,0), res );
  317. VERIFY_IS_APPROX( ( rmA*d.asDiagonal() ).eval().coeff(0,0), res );
  318. VERIFY_IS_APPROX( ( d.asDiagonal()*cmA ).eval().coeff(0,0), res );
  319. }
  320. template<typename Real>
  321. void test_mixing_types()
  322. {
  323. typedef std::complex<Real> Cplx;
  324. typedef SparseMatrix<Real> SpMatReal;
  325. typedef SparseMatrix<Cplx> SpMatCplx;
  326. typedef SparseMatrix<Cplx,RowMajor> SpRowMatCplx;
  327. typedef Matrix<Real,Dynamic,Dynamic> DenseMatReal;
  328. typedef Matrix<Cplx,Dynamic,Dynamic> DenseMatCplx;
  329. Index n = internal::random<Index>(1,100);
  330. double density = (std::max)(8./(n*n), 0.2);
  331. SpMatReal sR1(n,n);
  332. SpMatCplx sC1(n,n), sC2(n,n), sC3(n,n);
  333. SpRowMatCplx sCR(n,n);
  334. DenseMatReal dR1(n,n);
  335. DenseMatCplx dC1(n,n), dC2(n,n), dC3(n,n);
  336. initSparse<Real>(density, dR1, sR1);
  337. initSparse<Cplx>(density, dC1, sC1);
  338. initSparse<Cplx>(density, dC2, sC2);
  339. VERIFY_IS_APPROX( sC2 = (sR1 * sC1), dC3 = dR1.template cast<Cplx>() * dC1 );
  340. VERIFY_IS_APPROX( sC2 = (sC1 * sR1), dC3 = dC1 * dR1.template cast<Cplx>() );
  341. VERIFY_IS_APPROX( sC2 = (sR1.transpose() * sC1), dC3 = dR1.template cast<Cplx>().transpose() * dC1 );
  342. VERIFY_IS_APPROX( sC2 = (sC1.transpose() * sR1), dC3 = dC1.transpose() * dR1.template cast<Cplx>() );
  343. VERIFY_IS_APPROX( sC2 = (sR1 * sC1.transpose()), dC3 = dR1.template cast<Cplx>() * dC1.transpose() );
  344. VERIFY_IS_APPROX( sC2 = (sC1 * sR1.transpose()), dC3 = dC1 * dR1.template cast<Cplx>().transpose() );
  345. VERIFY_IS_APPROX( sC2 = (sR1.transpose() * sC1.transpose()), dC3 = dR1.template cast<Cplx>().transpose() * dC1.transpose() );
  346. VERIFY_IS_APPROX( sC2 = (sC1.transpose() * sR1.transpose()), dC3 = dC1.transpose() * dR1.template cast<Cplx>().transpose() );
  347. VERIFY_IS_APPROX( sCR = (sR1 * sC1), dC3 = dR1.template cast<Cplx>() * dC1 );
  348. VERIFY_IS_APPROX( sCR = (sC1 * sR1), dC3 = dC1 * dR1.template cast<Cplx>() );
  349. VERIFY_IS_APPROX( sCR = (sR1.transpose() * sC1), dC3 = dR1.template cast<Cplx>().transpose() * dC1 );
  350. VERIFY_IS_APPROX( sCR = (sC1.transpose() * sR1), dC3 = dC1.transpose() * dR1.template cast<Cplx>() );
  351. VERIFY_IS_APPROX( sCR = (sR1 * sC1.transpose()), dC3 = dR1.template cast<Cplx>() * dC1.transpose() );
  352. VERIFY_IS_APPROX( sCR = (sC1 * sR1.transpose()), dC3 = dC1 * dR1.template cast<Cplx>().transpose() );
  353. VERIFY_IS_APPROX( sCR = (sR1.transpose() * sC1.transpose()), dC3 = dR1.template cast<Cplx>().transpose() * dC1.transpose() );
  354. VERIFY_IS_APPROX( sCR = (sC1.transpose() * sR1.transpose()), dC3 = dC1.transpose() * dR1.template cast<Cplx>().transpose() );
  355. VERIFY_IS_APPROX( sC2 = (sR1 * sC1).pruned(), dC3 = dR1.template cast<Cplx>() * dC1 );
  356. VERIFY_IS_APPROX( sC2 = (sC1 * sR1).pruned(), dC3 = dC1 * dR1.template cast<Cplx>() );
  357. VERIFY_IS_APPROX( sC2 = (sR1.transpose() * sC1).pruned(), dC3 = dR1.template cast<Cplx>().transpose() * dC1 );
  358. VERIFY_IS_APPROX( sC2 = (sC1.transpose() * sR1).pruned(), dC3 = dC1.transpose() * dR1.template cast<Cplx>() );
  359. VERIFY_IS_APPROX( sC2 = (sR1 * sC1.transpose()).pruned(), dC3 = dR1.template cast<Cplx>() * dC1.transpose() );
  360. VERIFY_IS_APPROX( sC2 = (sC1 * sR1.transpose()).pruned(), dC3 = dC1 * dR1.template cast<Cplx>().transpose() );
  361. VERIFY_IS_APPROX( sC2 = (sR1.transpose() * sC1.transpose()).pruned(), dC3 = dR1.template cast<Cplx>().transpose() * dC1.transpose() );
  362. VERIFY_IS_APPROX( sC2 = (sC1.transpose() * sR1.transpose()).pruned(), dC3 = dC1.transpose() * dR1.template cast<Cplx>().transpose() );
  363. VERIFY_IS_APPROX( sCR = (sR1 * sC1).pruned(), dC3 = dR1.template cast<Cplx>() * dC1 );
  364. VERIFY_IS_APPROX( sCR = (sC1 * sR1).pruned(), dC3 = dC1 * dR1.template cast<Cplx>() );
  365. VERIFY_IS_APPROX( sCR = (sR1.transpose() * sC1).pruned(), dC3 = dR1.template cast<Cplx>().transpose() * dC1 );
  366. VERIFY_IS_APPROX( sCR = (sC1.transpose() * sR1).pruned(), dC3 = dC1.transpose() * dR1.template cast<Cplx>() );
  367. VERIFY_IS_APPROX( sCR = (sR1 * sC1.transpose()).pruned(), dC3 = dR1.template cast<Cplx>() * dC1.transpose() );
  368. VERIFY_IS_APPROX( sCR = (sC1 * sR1.transpose()).pruned(), dC3 = dC1 * dR1.template cast<Cplx>().transpose() );
  369. VERIFY_IS_APPROX( sCR = (sR1.transpose() * sC1.transpose()).pruned(), dC3 = dR1.template cast<Cplx>().transpose() * dC1.transpose() );
  370. VERIFY_IS_APPROX( sCR = (sC1.transpose() * sR1.transpose()).pruned(), dC3 = dC1.transpose() * dR1.template cast<Cplx>().transpose() );
  371. VERIFY_IS_APPROX( dC2 = (sR1 * sC1), dC3 = dR1.template cast<Cplx>() * dC1 );
  372. VERIFY_IS_APPROX( dC2 = (sC1 * sR1), dC3 = dC1 * dR1.template cast<Cplx>() );
  373. VERIFY_IS_APPROX( dC2 = (sR1.transpose() * sC1), dC3 = dR1.template cast<Cplx>().transpose() * dC1 );
  374. VERIFY_IS_APPROX( dC2 = (sC1.transpose() * sR1), dC3 = dC1.transpose() * dR1.template cast<Cplx>() );
  375. VERIFY_IS_APPROX( dC2 = (sR1 * sC1.transpose()), dC3 = dR1.template cast<Cplx>() * dC1.transpose() );
  376. VERIFY_IS_APPROX( dC2 = (sC1 * sR1.transpose()), dC3 = dC1 * dR1.template cast<Cplx>().transpose() );
  377. VERIFY_IS_APPROX( dC2 = (sR1.transpose() * sC1.transpose()), dC3 = dR1.template cast<Cplx>().transpose() * dC1.transpose() );
  378. VERIFY_IS_APPROX( dC2 = (sC1.transpose() * sR1.transpose()), dC3 = dC1.transpose() * dR1.template cast<Cplx>().transpose() );
  379. VERIFY_IS_APPROX( dC2 = dR1 * sC1, dC3 = dR1.template cast<Cplx>() * sC1 );
  380. VERIFY_IS_APPROX( dC2 = sR1 * dC1, dC3 = sR1.template cast<Cplx>() * dC1 );
  381. VERIFY_IS_APPROX( dC2 = dC1 * sR1, dC3 = dC1 * sR1.template cast<Cplx>() );
  382. VERIFY_IS_APPROX( dC2 = sC1 * dR1, dC3 = sC1 * dR1.template cast<Cplx>() );
  383. VERIFY_IS_APPROX( dC2 = dR1.row(0) * sC1, dC3 = dR1.template cast<Cplx>().row(0) * sC1 );
  384. VERIFY_IS_APPROX( dC2 = sR1 * dC1.col(0), dC3 = sR1.template cast<Cplx>() * dC1.col(0) );
  385. VERIFY_IS_APPROX( dC2 = dC1.row(0) * sR1, dC3 = dC1.row(0) * sR1.template cast<Cplx>() );
  386. VERIFY_IS_APPROX( dC2 = sC1 * dR1.col(0), dC3 = sC1 * dR1.template cast<Cplx>().col(0) );
  387. }
  388. void test_sparse_product()
  389. {
  390. for(int i = 0; i < g_repeat; i++) {
  391. CALL_SUBTEST_1( (sparse_product<SparseMatrix<double,ColMajor> >()) );
  392. CALL_SUBTEST_1( (sparse_product<SparseMatrix<double,RowMajor> >()) );
  393. CALL_SUBTEST_1( (bug_942<double>()) );
  394. CALL_SUBTEST_2( (sparse_product<SparseMatrix<std::complex<double>, ColMajor > >()) );
  395. CALL_SUBTEST_2( (sparse_product<SparseMatrix<std::complex<double>, RowMajor > >()) );
  396. CALL_SUBTEST_3( (sparse_product<SparseMatrix<float,ColMajor,long int> >()) );
  397. CALL_SUBTEST_4( (sparse_product_regression_test<SparseMatrix<double,RowMajor>, Matrix<double, Dynamic, Dynamic, RowMajor> >()) );
  398. CALL_SUBTEST_5( (test_mixing_types<float>()) );
  399. }
  400. }