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.

337 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. //
  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_product()
  11. {
  12. typedef typename SparseMatrixType::StorageIndex StorageIndex;
  13. Index n = 100;
  14. const Index rows = internal::random<Index>(1,n);
  15. const Index cols = internal::random<Index>(1,n);
  16. const Index depth = internal::random<Index>(1,n);
  17. typedef typename SparseMatrixType::Scalar Scalar;
  18. enum { Flags = SparseMatrixType::Flags };
  19. double density = (std::max)(8./(rows*cols), 0.2);
  20. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  21. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  22. typedef Matrix<Scalar,1,Dynamic> RowDenseVector;
  23. typedef SparseVector<Scalar,0,StorageIndex> ColSpVector;
  24. typedef SparseVector<Scalar,RowMajor,StorageIndex> RowSpVector;
  25. Scalar s1 = internal::random<Scalar>();
  26. Scalar s2 = internal::random<Scalar>();
  27. // test matrix-matrix product
  28. {
  29. DenseMatrix refMat2 = DenseMatrix::Zero(rows, depth);
  30. DenseMatrix refMat2t = DenseMatrix::Zero(depth, rows);
  31. DenseMatrix refMat3 = DenseMatrix::Zero(depth, cols);
  32. DenseMatrix refMat3t = DenseMatrix::Zero(cols, depth);
  33. DenseMatrix refMat4 = DenseMatrix::Zero(rows, cols);
  34. DenseMatrix refMat4t = DenseMatrix::Zero(cols, rows);
  35. DenseMatrix refMat5 = DenseMatrix::Random(depth, cols);
  36. DenseMatrix refMat6 = DenseMatrix::Random(rows, rows);
  37. DenseMatrix dm4 = DenseMatrix::Zero(rows, rows);
  38. // DenseVector dv1 = DenseVector::Random(rows);
  39. SparseMatrixType m2 (rows, depth);
  40. SparseMatrixType m2t(depth, rows);
  41. SparseMatrixType m3 (depth, cols);
  42. SparseMatrixType m3t(cols, depth);
  43. SparseMatrixType m4 (rows, cols);
  44. SparseMatrixType m4t(cols, rows);
  45. SparseMatrixType m6(rows, rows);
  46. initSparse(density, refMat2, m2);
  47. initSparse(density, refMat2t, m2t);
  48. initSparse(density, refMat3, m3);
  49. initSparse(density, refMat3t, m3t);
  50. initSparse(density, refMat4, m4);
  51. initSparse(density, refMat4t, m4t);
  52. initSparse(density, refMat6, m6);
  53. // int c = internal::random<int>(0,depth-1);
  54. // sparse * sparse
  55. VERIFY_IS_APPROX(m4=m2*m3, refMat4=refMat2*refMat3);
  56. VERIFY_IS_APPROX(m4=m2t.transpose()*m3, refMat4=refMat2t.transpose()*refMat3);
  57. VERIFY_IS_APPROX(m4=m2t.transpose()*m3t.transpose(), refMat4=refMat2t.transpose()*refMat3t.transpose());
  58. VERIFY_IS_APPROX(m4=m2*m3t.transpose(), refMat4=refMat2*refMat3t.transpose());
  59. VERIFY_IS_APPROX(m4 = m2*m3/s1, refMat4 = refMat2*refMat3/s1);
  60. VERIFY_IS_APPROX(m4 = m2*m3*s1, refMat4 = refMat2*refMat3*s1);
  61. VERIFY_IS_APPROX(m4 = s2*m2*m3*s1, refMat4 = s2*refMat2*refMat3*s1);
  62. VERIFY_IS_APPROX(m4 = (m2+m2)*m3, refMat4 = (refMat2+refMat2)*refMat3);
  63. VERIFY_IS_APPROX(m4 = m2*m3.leftCols(cols/2), refMat4 = refMat2*refMat3.leftCols(cols/2));
  64. VERIFY_IS_APPROX(m4 = m2*(m3+m3).leftCols(cols/2), refMat4 = refMat2*(refMat3+refMat3).leftCols(cols/2));
  65. VERIFY_IS_APPROX(m4=(m2*m3).pruned(0), refMat4=refMat2*refMat3);
  66. VERIFY_IS_APPROX(m4=(m2t.transpose()*m3).pruned(0), refMat4=refMat2t.transpose()*refMat3);
  67. VERIFY_IS_APPROX(m4=(m2t.transpose()*m3t.transpose()).pruned(0), refMat4=refMat2t.transpose()*refMat3t.transpose());
  68. VERIFY_IS_APPROX(m4=(m2*m3t.transpose()).pruned(0), refMat4=refMat2*refMat3t.transpose());
  69. // dense ?= sparse * sparse
  70. VERIFY_IS_APPROX(dm4 =m2*m3, refMat4 =refMat2*refMat3);
  71. VERIFY_IS_APPROX(dm4+=m2*m3, refMat4+=refMat2*refMat3);
  72. VERIFY_IS_APPROX(dm4-=m2*m3, refMat4-=refMat2*refMat3);
  73. VERIFY_IS_APPROX(dm4 =m2t.transpose()*m3, refMat4 =refMat2t.transpose()*refMat3);
  74. VERIFY_IS_APPROX(dm4+=m2t.transpose()*m3, refMat4+=refMat2t.transpose()*refMat3);
  75. VERIFY_IS_APPROX(dm4-=m2t.transpose()*m3, refMat4-=refMat2t.transpose()*refMat3);
  76. VERIFY_IS_APPROX(dm4 =m2t.transpose()*m3t.transpose(), refMat4 =refMat2t.transpose()*refMat3t.transpose());
  77. VERIFY_IS_APPROX(dm4+=m2t.transpose()*m3t.transpose(), refMat4+=refMat2t.transpose()*refMat3t.transpose());
  78. VERIFY_IS_APPROX(dm4-=m2t.transpose()*m3t.transpose(), refMat4-=refMat2t.transpose()*refMat3t.transpose());
  79. VERIFY_IS_APPROX(dm4 =m2*m3t.transpose(), refMat4 =refMat2*refMat3t.transpose());
  80. VERIFY_IS_APPROX(dm4+=m2*m3t.transpose(), refMat4+=refMat2*refMat3t.transpose());
  81. VERIFY_IS_APPROX(dm4-=m2*m3t.transpose(), refMat4-=refMat2*refMat3t.transpose());
  82. VERIFY_IS_APPROX(dm4 = m2*m3*s1, refMat4 = refMat2*refMat3*s1);
  83. // test aliasing
  84. m4 = m2; refMat4 = refMat2;
  85. VERIFY_IS_APPROX(m4=m4*m3, refMat4=refMat4*refMat3);
  86. // sparse * dense matrix
  87. VERIFY_IS_APPROX(dm4=m2*refMat3, refMat4=refMat2*refMat3);
  88. VERIFY_IS_APPROX(dm4=m2*refMat3t.transpose(), refMat4=refMat2*refMat3t.transpose());
  89. VERIFY_IS_APPROX(dm4=m2t.transpose()*refMat3, refMat4=refMat2t.transpose()*refMat3);
  90. VERIFY_IS_APPROX(dm4=m2t.transpose()*refMat3t.transpose(), refMat4=refMat2t.transpose()*refMat3t.transpose());
  91. VERIFY_IS_APPROX(dm4=m2*refMat3, refMat4=refMat2*refMat3);
  92. VERIFY_IS_APPROX(dm4=dm4+m2*refMat3, refMat4=refMat4+refMat2*refMat3);
  93. VERIFY_IS_APPROX(dm4+=m2*refMat3, refMat4+=refMat2*refMat3);
  94. VERIFY_IS_APPROX(dm4-=m2*refMat3, refMat4-=refMat2*refMat3);
  95. VERIFY_IS_APPROX(dm4.noalias()+=m2*refMat3, refMat4+=refMat2*refMat3);
  96. VERIFY_IS_APPROX(dm4.noalias()-=m2*refMat3, refMat4-=refMat2*refMat3);
  97. VERIFY_IS_APPROX(dm4=m2*(refMat3+refMat3), refMat4=refMat2*(refMat3+refMat3));
  98. VERIFY_IS_APPROX(dm4=m2t.transpose()*(refMat3+refMat5)*0.5, refMat4=refMat2t.transpose()*(refMat3+refMat5)*0.5);
  99. // sparse * dense vector
  100. VERIFY_IS_APPROX(dm4.col(0)=m2*refMat3.col(0), refMat4.col(0)=refMat2*refMat3.col(0));
  101. VERIFY_IS_APPROX(dm4.col(0)=m2*refMat3t.transpose().col(0), refMat4.col(0)=refMat2*refMat3t.transpose().col(0));
  102. VERIFY_IS_APPROX(dm4.col(0)=m2t.transpose()*refMat3.col(0), refMat4.col(0)=refMat2t.transpose()*refMat3.col(0));
  103. VERIFY_IS_APPROX(dm4.col(0)=m2t.transpose()*refMat3t.transpose().col(0), refMat4.col(0)=refMat2t.transpose()*refMat3t.transpose().col(0));
  104. // dense * sparse
  105. VERIFY_IS_APPROX(dm4=refMat2*m3, refMat4=refMat2*refMat3);
  106. VERIFY_IS_APPROX(dm4=dm4+refMat2*m3, refMat4=refMat4+refMat2*refMat3);
  107. VERIFY_IS_APPROX(dm4+=refMat2*m3, refMat4+=refMat2*refMat3);
  108. VERIFY_IS_APPROX(dm4-=refMat2*m3, refMat4-=refMat2*refMat3);
  109. VERIFY_IS_APPROX(dm4.noalias()+=refMat2*m3, refMat4+=refMat2*refMat3);
  110. VERIFY_IS_APPROX(dm4.noalias()-=refMat2*m3, refMat4-=refMat2*refMat3);
  111. VERIFY_IS_APPROX(dm4=refMat2*m3t.transpose(), refMat4=refMat2*refMat3t.transpose());
  112. VERIFY_IS_APPROX(dm4=refMat2t.transpose()*m3, refMat4=refMat2t.transpose()*refMat3);
  113. VERIFY_IS_APPROX(dm4=refMat2t.transpose()*m3t.transpose(), refMat4=refMat2t.transpose()*refMat3t.transpose());
  114. // sparse * dense and dense * sparse outer product
  115. {
  116. Index c = internal::random<Index>(0,depth-1);
  117. Index r = internal::random<Index>(0,rows-1);
  118. Index c1 = internal::random<Index>(0,cols-1);
  119. Index r1 = internal::random<Index>(0,depth-1);
  120. DenseMatrix dm5 = DenseMatrix::Random(depth, cols);
  121. VERIFY_IS_APPROX( m4=m2.col(c)*dm5.col(c1).transpose(), refMat4=refMat2.col(c)*dm5.col(c1).transpose());
  122. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  123. VERIFY_IS_APPROX( m4=m2.middleCols(c,1)*dm5.col(c1).transpose(), refMat4=refMat2.col(c)*dm5.col(c1).transpose());
  124. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  125. VERIFY_IS_APPROX(dm4=m2.col(c)*dm5.col(c1).transpose(), refMat4=refMat2.col(c)*dm5.col(c1).transpose());
  126. VERIFY_IS_APPROX(m4=dm5.col(c1)*m2.col(c).transpose(), refMat4=dm5.col(c1)*refMat2.col(c).transpose());
  127. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  128. VERIFY_IS_APPROX(m4=dm5.col(c1)*m2.middleCols(c,1).transpose(), refMat4=dm5.col(c1)*refMat2.col(c).transpose());
  129. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  130. VERIFY_IS_APPROX(dm4=dm5.col(c1)*m2.col(c).transpose(), refMat4=dm5.col(c1)*refMat2.col(c).transpose());
  131. VERIFY_IS_APPROX( m4=dm5.row(r1).transpose()*m2.col(c).transpose(), refMat4=dm5.row(r1).transpose()*refMat2.col(c).transpose());
  132. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  133. VERIFY_IS_APPROX(dm4=dm5.row(r1).transpose()*m2.col(c).transpose(), refMat4=dm5.row(r1).transpose()*refMat2.col(c).transpose());
  134. VERIFY_IS_APPROX( m4=m2.row(r).transpose()*dm5.col(c1).transpose(), refMat4=refMat2.row(r).transpose()*dm5.col(c1).transpose());
  135. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  136. VERIFY_IS_APPROX( m4=m2.middleRows(r,1).transpose()*dm5.col(c1).transpose(), refMat4=refMat2.row(r).transpose()*dm5.col(c1).transpose());
  137. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  138. VERIFY_IS_APPROX(dm4=m2.row(r).transpose()*dm5.col(c1).transpose(), refMat4=refMat2.row(r).transpose()*dm5.col(c1).transpose());
  139. VERIFY_IS_APPROX( m4=dm5.col(c1)*m2.row(r), refMat4=dm5.col(c1)*refMat2.row(r));
  140. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  141. VERIFY_IS_APPROX( m4=dm5.col(c1)*m2.middleRows(r,1), refMat4=dm5.col(c1)*refMat2.row(r));
  142. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  143. VERIFY_IS_APPROX(dm4=dm5.col(c1)*m2.row(r), refMat4=dm5.col(c1)*refMat2.row(r));
  144. VERIFY_IS_APPROX( m4=dm5.row(r1).transpose()*m2.row(r), refMat4=dm5.row(r1).transpose()*refMat2.row(r));
  145. VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
  146. VERIFY_IS_APPROX(dm4=dm5.row(r1).transpose()*m2.row(r), refMat4=dm5.row(r1).transpose()*refMat2.row(r));
  147. }
  148. VERIFY_IS_APPROX(m6=m6*m6, refMat6=refMat6*refMat6);
  149. // sparse matrix * sparse vector
  150. ColSpVector cv0(cols), cv1;
  151. DenseVector dcv0(cols), dcv1;
  152. initSparse(2*density,dcv0, cv0);
  153. RowSpVector rv0(depth), rv1;
  154. RowDenseVector drv0(depth), drv1(rv1);
  155. initSparse(2*density,drv0, rv0);
  156. VERIFY_IS_APPROX(cv1=m3*cv0, dcv1=refMat3*dcv0);
  157. VERIFY_IS_APPROX(rv1=rv0*m3, drv1=drv0*refMat3);
  158. VERIFY_IS_APPROX(cv1=m3t.adjoint()*cv0, dcv1=refMat3t.adjoint()*dcv0);
  159. VERIFY_IS_APPROX(cv1=rv0*m3, dcv1=drv0*refMat3);
  160. VERIFY_IS_APPROX(rv1=m3*cv0, drv1=refMat3*dcv0);
  161. }
  162. // test matrix - diagonal product
  163. {
  164. DenseMatrix refM2 = DenseMatrix::Zero(rows, cols);
  165. DenseMatrix refM3 = DenseMatrix::Zero(rows, cols);
  166. DenseMatrix d3 = DenseMatrix::Zero(rows, cols);
  167. DiagonalMatrix<Scalar,Dynamic> d1(DenseVector::Random(cols));
  168. DiagonalMatrix<Scalar,Dynamic> d2(DenseVector::Random(rows));
  169. SparseMatrixType m2(rows, cols);
  170. SparseMatrixType m3(rows, cols);
  171. initSparse<Scalar>(density, refM2, m2);
  172. initSparse<Scalar>(density, refM3, m3);
  173. VERIFY_IS_APPROX(m3=m2*d1, refM3=refM2*d1);
  174. VERIFY_IS_APPROX(m3=m2.transpose()*d2, refM3=refM2.transpose()*d2);
  175. VERIFY_IS_APPROX(m3=d2*m2, refM3=d2*refM2);
  176. VERIFY_IS_APPROX(m3=d1*m2.transpose(), refM3=d1*refM2.transpose());
  177. // also check with a SparseWrapper:
  178. DenseVector v1 = DenseVector::Random(cols);
  179. DenseVector v2 = DenseVector::Random(rows);
  180. VERIFY_IS_APPROX(m3=m2*v1.asDiagonal(), refM3=refM2*v1.asDiagonal());
  181. VERIFY_IS_APPROX(m3=m2.transpose()*v2.asDiagonal(), refM3=refM2.transpose()*v2.asDiagonal());
  182. VERIFY_IS_APPROX(m3=v2.asDiagonal()*m2, refM3=v2.asDiagonal()*refM2);
  183. VERIFY_IS_APPROX(m3=v1.asDiagonal()*m2.transpose(), refM3=v1.asDiagonal()*refM2.transpose());
  184. VERIFY_IS_APPROX(m3=v2.asDiagonal()*m2*v1.asDiagonal(), refM3=v2.asDiagonal()*refM2*v1.asDiagonal());
  185. // evaluate to a dense matrix to check the .row() and .col() iterator functions
  186. VERIFY_IS_APPROX(d3=m2*d1, refM3=refM2*d1);
  187. VERIFY_IS_APPROX(d3=m2.transpose()*d2, refM3=refM2.transpose()*d2);
  188. VERIFY_IS_APPROX(d3=d2*m2, refM3=d2*refM2);
  189. VERIFY_IS_APPROX(d3=d1*m2.transpose(), refM3=d1*refM2.transpose());
  190. }
  191. // test self-adjoint and triangular-view products
  192. {
  193. DenseMatrix b = DenseMatrix::Random(rows, rows);
  194. DenseMatrix x = DenseMatrix::Random(rows, rows);
  195. DenseMatrix refX = DenseMatrix::Random(rows, rows);
  196. DenseMatrix refUp = DenseMatrix::Zero(rows, rows);
  197. DenseMatrix refLo = DenseMatrix::Zero(rows, rows);
  198. DenseMatrix refS = DenseMatrix::Zero(rows, rows);
  199. DenseMatrix refA = DenseMatrix::Zero(rows, rows);
  200. SparseMatrixType mUp(rows, rows);
  201. SparseMatrixType mLo(rows, rows);
  202. SparseMatrixType mS(rows, rows);
  203. SparseMatrixType mA(rows, rows);
  204. initSparse<Scalar>(density, refA, mA);
  205. do {
  206. initSparse<Scalar>(density, refUp, mUp, ForceRealDiag|/*ForceNonZeroDiag|*/MakeUpperTriangular);
  207. } while (refUp.isZero());
  208. refLo = refUp.adjoint();
  209. mLo = mUp.adjoint();
  210. refS = refUp + refLo;
  211. refS.diagonal() *= 0.5;
  212. mS = mUp + mLo;
  213. // TODO be able to address the diagonal....
  214. for (int k=0; k<mS.outerSize(); ++k)
  215. for (typename SparseMatrixType::InnerIterator it(mS,k); it; ++it)
  216. if (it.index() == k)
  217. it.valueRef() *= 0.5;
  218. VERIFY_IS_APPROX(refS.adjoint(), refS);
  219. VERIFY_IS_APPROX(mS.adjoint(), mS);
  220. VERIFY_IS_APPROX(mS, refS);
  221. VERIFY_IS_APPROX(x=mS*b, refX=refS*b);
  222. // sparse selfadjointView with dense matrices
  223. VERIFY_IS_APPROX(x=mUp.template selfadjointView<Upper>()*b, refX=refS*b);
  224. VERIFY_IS_APPROX(x=mLo.template selfadjointView<Lower>()*b, refX=refS*b);
  225. VERIFY_IS_APPROX(x=mS.template selfadjointView<Upper|Lower>()*b, refX=refS*b);
  226. // sparse selfadjointView with sparse matrices
  227. SparseMatrixType mSres(rows,rows);
  228. VERIFY_IS_APPROX(mSres = mLo.template selfadjointView<Lower>()*mS,
  229. refX = refLo.template selfadjointView<Lower>()*refS);
  230. VERIFY_IS_APPROX(mSres = mS * mLo.template selfadjointView<Lower>(),
  231. refX = refS * refLo.template selfadjointView<Lower>());
  232. // sparse triangularView with dense matrices
  233. VERIFY_IS_APPROX(x=mA.template triangularView<Upper>()*b, refX=refA.template triangularView<Upper>()*b);
  234. VERIFY_IS_APPROX(x=mA.template triangularView<Lower>()*b, refX=refA.template triangularView<Lower>()*b);
  235. VERIFY_IS_APPROX(x=b*mA.template triangularView<Upper>(), refX=b*refA.template triangularView<Upper>());
  236. VERIFY_IS_APPROX(x=b*mA.template triangularView<Lower>(), refX=b*refA.template triangularView<Lower>());
  237. // sparse triangularView with sparse matrices
  238. VERIFY_IS_APPROX(mSres = mA.template triangularView<Lower>()*mS, refX = refA.template triangularView<Lower>()*refS);
  239. VERIFY_IS_APPROX(mSres = mS * mA.template triangularView<Lower>(), refX = refS * refA.template triangularView<Lower>());
  240. VERIFY_IS_APPROX(mSres = mA.template triangularView<Upper>()*mS, refX = refA.template triangularView<Upper>()*refS);
  241. VERIFY_IS_APPROX(mSres = mS * mA.template triangularView<Upper>(), refX = refS * refA.template triangularView<Upper>());
  242. }
  243. }
  244. // New test for Bug in SparseTimeDenseProduct
  245. template<typename SparseMatrixType, typename DenseMatrixType> void sparse_product_regression_test()
  246. {
  247. // This code does not compile with afflicted versions of the bug
  248. SparseMatrixType sm1(3,2);
  249. DenseMatrixType m2(2,2);
  250. sm1.setZero();
  251. m2.setZero();
  252. DenseMatrixType m3 = sm1*m2;
  253. // This code produces a segfault with afflicted versions of another SparseTimeDenseProduct
  254. // bug
  255. SparseMatrixType sm2(20000,2);
  256. sm2.setZero();
  257. DenseMatrixType m4(sm2*m2);
  258. VERIFY_IS_APPROX( m4(0,0), 0.0 );
  259. }
  260. template<typename Scalar>
  261. void bug_942()
  262. {
  263. typedef Matrix<Scalar, Dynamic, 1> Vector;
  264. typedef SparseMatrix<Scalar, ColMajor> ColSpMat;
  265. typedef SparseMatrix<Scalar, RowMajor> RowSpMat;
  266. ColSpMat cmA(1,1);
  267. cmA.insert(0,0) = 1;
  268. RowSpMat rmA(1,1);
  269. rmA.insert(0,0) = 1;
  270. Vector d(1);
  271. d[0] = 2;
  272. double res = 2;
  273. VERIFY_IS_APPROX( ( cmA*d.asDiagonal() ).eval().coeff(0,0), res );
  274. VERIFY_IS_APPROX( ( d.asDiagonal()*rmA ).eval().coeff(0,0), res );
  275. VERIFY_IS_APPROX( ( rmA*d.asDiagonal() ).eval().coeff(0,0), res );
  276. VERIFY_IS_APPROX( ( d.asDiagonal()*cmA ).eval().coeff(0,0), res );
  277. }
  278. void test_sparse_product()
  279. {
  280. for(int i = 0; i < g_repeat; i++) {
  281. CALL_SUBTEST_1( (sparse_product<SparseMatrix<double,ColMajor> >()) );
  282. CALL_SUBTEST_1( (sparse_product<SparseMatrix<double,RowMajor> >()) );
  283. CALL_SUBTEST_1( (bug_942<double>()) );
  284. CALL_SUBTEST_2( (sparse_product<SparseMatrix<std::complex<double>, ColMajor > >()) );
  285. CALL_SUBTEST_2( (sparse_product<SparseMatrix<std::complex<double>, RowMajor > >()) );
  286. CALL_SUBTEST_3( (sparse_product<SparseMatrix<float,ColMajor,long int> >()) );
  287. CALL_SUBTEST_4( (sparse_product_regression_test<SparseMatrix<double,RowMajor>, Matrix<double, Dynamic, Dynamic, RowMajor> >()) );
  288. }
  289. }