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


gtsam
Author(s):
autogenerated on Tue Jul 4 2023 02:35:54