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 
on_temporary_creation()18 inline void on_temporary_creation() {
19   // here's a great place to set a breakpoint when debugging failures in this test!
20   nb_temporaries++;
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 
sparse_product()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);
47   typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
48   typedef Matrix<Scalar,Dynamic,1> DenseVector;
49   typedef Matrix<Scalar,1,Dynamic> RowDenseVector;
50   typedef SparseVector<Scalar,0,StorageIndex> ColSpVector;
51   typedef SparseVector<Scalar,RowMajor,StorageIndex> RowSpVector;
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     // make sure the right product implementation is called:
104     if((!SparseMatrixType::IsRowMajor) && m2.rows()<=m3.cols())
105     {
106       VERIFY_EVALUATION_COUNT(m4 = m2*m3, 3); // 1 temp for the result + 2 for transposing and get a sorted result.
107       VERIFY_EVALUATION_COUNT(m4 = (m2*m3).pruned(0), 1);
108       VERIFY_EVALUATION_COUNT(m4 = (m2*m3).eval().pruned(0), 4);
109     }
110 
111     // and that pruning is effective:
112     {
113       DenseMatrix Ad(2,2);
114       Ad << -1, 1, 1, 1;
115       SparseMatrixType As(Ad.sparseView()), B(2,2);
116       VERIFY_IS_EQUAL( (As*As.transpose()).eval().nonZeros(), 4);
117       VERIFY_IS_EQUAL( (Ad*Ad.transpose()).eval().sparseView().eval().nonZeros(), 2);
118       VERIFY_IS_EQUAL( (As*As.transpose()).pruned(1e-6).eval().nonZeros(), 2);
119     }
120 
121     // dense ?= sparse * sparse
122     VERIFY_IS_APPROX(dm4 =m2*m3, refMat4 =refMat2*refMat3);
123     VERIFY_IS_APPROX(dm4+=m2*m3, refMat4+=refMat2*refMat3);
124     VERIFY_IS_APPROX(dm4-=m2*m3, refMat4-=refMat2*refMat3);
125     VERIFY_IS_APPROX(dm4 =m2t.transpose()*m3, refMat4 =refMat2t.transpose()*refMat3);
126     VERIFY_IS_APPROX(dm4+=m2t.transpose()*m3, refMat4+=refMat2t.transpose()*refMat3);
127     VERIFY_IS_APPROX(dm4-=m2t.transpose()*m3, refMat4-=refMat2t.transpose()*refMat3);
128     VERIFY_IS_APPROX(dm4 =m2t.transpose()*m3t.transpose(), refMat4 =refMat2t.transpose()*refMat3t.transpose());
129     VERIFY_IS_APPROX(dm4+=m2t.transpose()*m3t.transpose(), refMat4+=refMat2t.transpose()*refMat3t.transpose());
130     VERIFY_IS_APPROX(dm4-=m2t.transpose()*m3t.transpose(), refMat4-=refMat2t.transpose()*refMat3t.transpose());
131     VERIFY_IS_APPROX(dm4 =m2*m3t.transpose(), refMat4 =refMat2*refMat3t.transpose());
132     VERIFY_IS_APPROX(dm4+=m2*m3t.transpose(), refMat4+=refMat2*refMat3t.transpose());
133     VERIFY_IS_APPROX(dm4-=m2*m3t.transpose(), refMat4-=refMat2*refMat3t.transpose());
134     VERIFY_IS_APPROX(dm4 = m2*m3*s1, refMat4 = refMat2*refMat3*s1);
135 
136     // test aliasing
137     m4 = m2; refMat4 = refMat2;
138     VERIFY_IS_APPROX(m4=m4*m3, refMat4=refMat4*refMat3);
139 
140     // sparse * dense matrix
141     VERIFY_IS_APPROX(dm4=m2*refMat3, refMat4=refMat2*refMat3);
142     VERIFY_IS_APPROX(dm4=m2*refMat3t.transpose(), refMat4=refMat2*refMat3t.transpose());
143     VERIFY_IS_APPROX(dm4=m2t.transpose()*refMat3, refMat4=refMat2t.transpose()*refMat3);
144     VERIFY_IS_APPROX(dm4=m2t.transpose()*refMat3t.transpose(), refMat4=refMat2t.transpose()*refMat3t.transpose());
145 
146     VERIFY_IS_APPROX(dm4=m2*refMat3, refMat4=refMat2*refMat3);
147     VERIFY_IS_APPROX(dm4=dm4+m2*refMat3, refMat4=refMat4+refMat2*refMat3);
148     VERIFY_IS_APPROX(dm4+=m2*refMat3, refMat4+=refMat2*refMat3);
149     VERIFY_IS_APPROX(dm4-=m2*refMat3, refMat4-=refMat2*refMat3);
150     VERIFY_IS_APPROX(dm4.noalias()+=m2*refMat3, refMat4+=refMat2*refMat3);
151     VERIFY_IS_APPROX(dm4.noalias()-=m2*refMat3, refMat4-=refMat2*refMat3);
152     VERIFY_IS_APPROX(dm4=m2*(refMat3+refMat3), refMat4=refMat2*(refMat3+refMat3));
153     VERIFY_IS_APPROX(dm4=m2t.transpose()*(refMat3+refMat5)*0.5, refMat4=refMat2t.transpose()*(refMat3+refMat5)*0.5);
154 
155     // sparse * dense vector
156     VERIFY_IS_APPROX(dm4.col(0)=m2*refMat3.col(0), refMat4.col(0)=refMat2*refMat3.col(0));
157     VERIFY_IS_APPROX(dm4.col(0)=m2*refMat3t.transpose().col(0), refMat4.col(0)=refMat2*refMat3t.transpose().col(0));
158     VERIFY_IS_APPROX(dm4.col(0)=m2t.transpose()*refMat3.col(0), refMat4.col(0)=refMat2t.transpose()*refMat3.col(0));
159     VERIFY_IS_APPROX(dm4.col(0)=m2t.transpose()*refMat3t.transpose().col(0), refMat4.col(0)=refMat2t.transpose()*refMat3t.transpose().col(0));
160 
161     // dense * sparse
162     VERIFY_IS_APPROX(dm4=refMat2*m3, refMat4=refMat2*refMat3);
163     VERIFY_IS_APPROX(dm4=dm4+refMat2*m3, refMat4=refMat4+refMat2*refMat3);
164     VERIFY_IS_APPROX(dm4+=refMat2*m3, refMat4+=refMat2*refMat3);
165     VERIFY_IS_APPROX(dm4-=refMat2*m3, refMat4-=refMat2*refMat3);
166     VERIFY_IS_APPROX(dm4.noalias()+=refMat2*m3, refMat4+=refMat2*refMat3);
167     VERIFY_IS_APPROX(dm4.noalias()-=refMat2*m3, refMat4-=refMat2*refMat3);
168     VERIFY_IS_APPROX(dm4=refMat2*m3t.transpose(), refMat4=refMat2*refMat3t.transpose());
169     VERIFY_IS_APPROX(dm4=refMat2t.transpose()*m3, refMat4=refMat2t.transpose()*refMat3);
170     VERIFY_IS_APPROX(dm4=refMat2t.transpose()*m3t.transpose(), refMat4=refMat2t.transpose()*refMat3t.transpose());
171 
172     // sparse * dense and dense * sparse outer product
173     {
174       Index c  = internal::random<Index>(0,depth-1);
175       Index r  = internal::random<Index>(0,rows-1);
176       Index c1 = internal::random<Index>(0,cols-1);
177       Index r1 = internal::random<Index>(0,depth-1);
178       DenseMatrix dm5  = DenseMatrix::Random(depth, cols);
179 
180       VERIFY_IS_APPROX( m4=m2.col(c)*dm5.col(c1).transpose(), refMat4=refMat2.col(c)*dm5.col(c1).transpose());
181       VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
182       VERIFY_IS_APPROX( m4=m2.middleCols(c,1)*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(dm4=m2.col(c)*dm5.col(c1).transpose(), refMat4=refMat2.col(c)*dm5.col(c1).transpose());
185 
186       VERIFY_IS_APPROX(m4=dm5.col(c1)*m2.col(c).transpose(), refMat4=dm5.col(c1)*refMat2.col(c).transpose());
187       VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
188       VERIFY_IS_APPROX(m4=dm5.col(c1)*m2.middleCols(c,1).transpose(), refMat4=dm5.col(c1)*refMat2.col(c).transpose());
189       VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
190       VERIFY_IS_APPROX(dm4=dm5.col(c1)*m2.col(c).transpose(), refMat4=dm5.col(c1)*refMat2.col(c).transpose());
191 
192       VERIFY_IS_APPROX( m4=dm5.row(r1).transpose()*m2.col(c).transpose(), refMat4=dm5.row(r1).transpose()*refMat2.col(c).transpose());
193       VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
194       VERIFY_IS_APPROX(dm4=dm5.row(r1).transpose()*m2.col(c).transpose(), refMat4=dm5.row(r1).transpose()*refMat2.col(c).transpose());
195 
196       VERIFY_IS_APPROX( m4=m2.row(r).transpose()*dm5.col(c1).transpose(), refMat4=refMat2.row(r).transpose()*dm5.col(c1).transpose());
197       VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
198       VERIFY_IS_APPROX( m4=m2.middleRows(r,1).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(dm4=m2.row(r).transpose()*dm5.col(c1).transpose(), refMat4=refMat2.row(r).transpose()*dm5.col(c1).transpose());
201 
202       VERIFY_IS_APPROX( m4=dm5.col(c1)*m2.row(r), refMat4=dm5.col(c1)*refMat2.row(r));
203       VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
204       VERIFY_IS_APPROX( m4=dm5.col(c1)*m2.middleRows(r,1), refMat4=dm5.col(c1)*refMat2.row(r));
205       VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
206       VERIFY_IS_APPROX(dm4=dm5.col(c1)*m2.row(r), refMat4=dm5.col(c1)*refMat2.row(r));
207 
208       VERIFY_IS_APPROX( m4=dm5.row(r1).transpose()*m2.row(r), refMat4=dm5.row(r1).transpose()*refMat2.row(r));
209       VERIFY_IS_EQUAL(m4.nonZeros(), (refMat4.array()!=0).count());
210       VERIFY_IS_APPROX(dm4=dm5.row(r1).transpose()*m2.row(r), refMat4=dm5.row(r1).transpose()*refMat2.row(r));
211     }
212 
213     VERIFY_IS_APPROX(m6=m6*m6, refMat6=refMat6*refMat6);
214 
215     // sparse matrix * sparse vector
216     ColSpVector cv0(cols), cv1;
217     DenseVector dcv0(cols), dcv1;
218     initSparse(2*density,dcv0, cv0);
219 
220     RowSpVector rv0(depth), rv1;
221     RowDenseVector drv0(depth), drv1(rv1);
222     initSparse(2*density,drv0, rv0);
223 
224     VERIFY_IS_APPROX(cv1=m3*cv0, dcv1=refMat3*dcv0);
225     VERIFY_IS_APPROX(rv1=rv0*m3, drv1=drv0*refMat3);
226     VERIFY_IS_APPROX(cv1=m3t.adjoint()*cv0, dcv1=refMat3t.adjoint()*dcv0);
227     VERIFY_IS_APPROX(cv1=rv0*m3, dcv1=drv0*refMat3);
228     VERIFY_IS_APPROX(rv1=m3*cv0, drv1=refMat3*dcv0);
229   }
230 
231   // test matrix - diagonal product
232   {
233     DenseMatrix refM2 = DenseMatrix::Zero(rows, cols);
234     DenseMatrix refM3 = DenseMatrix::Zero(rows, cols);
235     DenseMatrix d3 = DenseMatrix::Zero(rows, cols);
236     DiagonalMatrix<Scalar,Dynamic> d1(DenseVector::Random(cols));
237     DiagonalMatrix<Scalar,Dynamic> d2(DenseVector::Random(rows));
238     SparseMatrixType m2(rows, cols);
239     SparseMatrixType m3(rows, cols);
240     initSparse<Scalar>(density, refM2, m2);
241     initSparse<Scalar>(density, refM3, m3);
242     VERIFY_IS_APPROX(m3=m2*d1, refM3=refM2*d1);
243     VERIFY_IS_APPROX(m3=m2.transpose()*d2, refM3=refM2.transpose()*d2);
244     VERIFY_IS_APPROX(m3=d2*m2, refM3=d2*refM2);
245     VERIFY_IS_APPROX(m3=d1*m2.transpose(), refM3=d1*refM2.transpose());
246 
247     // also check with a SparseWrapper:
248     DenseVector v1 = DenseVector::Random(cols);
249     DenseVector v2 = DenseVector::Random(rows);
250     DenseVector v3 = DenseVector::Random(rows);
251     VERIFY_IS_APPROX(m3=m2*v1.asDiagonal(), refM3=refM2*v1.asDiagonal());
252     VERIFY_IS_APPROX(m3=m2.transpose()*v2.asDiagonal(), refM3=refM2.transpose()*v2.asDiagonal());
253     VERIFY_IS_APPROX(m3=v2.asDiagonal()*m2, refM3=v2.asDiagonal()*refM2);
254     VERIFY_IS_APPROX(m3=v1.asDiagonal()*m2.transpose(), refM3=v1.asDiagonal()*refM2.transpose());
255 
256     VERIFY_IS_APPROX(m3=v2.asDiagonal()*m2*v1.asDiagonal(), refM3=v2.asDiagonal()*refM2*v1.asDiagonal());
257 
258     VERIFY_IS_APPROX(v2=m2*v1.asDiagonal()*v1, refM2*v1.asDiagonal()*v1);
259     VERIFY_IS_APPROX(v3=v2.asDiagonal()*m2*v1, v2.asDiagonal()*refM2*v1);
260 
261     // evaluate to a dense matrix to check the .row() and .col() iterator functions
262     VERIFY_IS_APPROX(d3=m2*d1, refM3=refM2*d1);
263     VERIFY_IS_APPROX(d3=m2.transpose()*d2, refM3=refM2.transpose()*d2);
264     VERIFY_IS_APPROX(d3=d2*m2, refM3=d2*refM2);
265     VERIFY_IS_APPROX(d3=d1*m2.transpose(), refM3=d1*refM2.transpose());
266   }
267 
268   // test self-adjoint and triangular-view products
269   {
270     DenseMatrix b = DenseMatrix::Random(rows, rows);
271     DenseMatrix x = DenseMatrix::Random(rows, rows);
272     DenseMatrix refX = DenseMatrix::Random(rows, rows);
273     DenseMatrix refUp = DenseMatrix::Zero(rows, rows);
274     DenseMatrix refLo = DenseMatrix::Zero(rows, rows);
275     DenseMatrix refS = DenseMatrix::Zero(rows, rows);
276     DenseMatrix refA = DenseMatrix::Zero(rows, rows);
277     SparseMatrixType mUp(rows, rows);
278     SparseMatrixType mLo(rows, rows);
279     SparseMatrixType mS(rows, rows);
280     SparseMatrixType mA(rows, rows);
281     initSparse<Scalar>(density, refA, mA);
282     do {
283       initSparse<Scalar>(density, refUp, mUp, ForceRealDiag|/*ForceNonZeroDiag|*/MakeUpperTriangular);
284     } while (refUp.isZero());
285     refLo = refUp.adjoint();
286     mLo = mUp.adjoint();
287     refS = refUp + refLo;
288     refS.diagonal() *= 0.5;
289     mS = mUp + mLo;
290     // TODO be able to address the diagonal....
291     for (int k=0; k<mS.outerSize(); ++k)
292       for (typename SparseMatrixType::InnerIterator it(mS,k); it; ++it)
293         if (it.index() == k)
294           it.valueRef() *= Scalar(0.5);
295 
296     VERIFY_IS_APPROX(refS.adjoint(), refS);
297     VERIFY_IS_APPROX(mS.adjoint(), mS);
298     VERIFY_IS_APPROX(mS, refS);
299     VERIFY_IS_APPROX(x=mS*b, refX=refS*b);
300 
301     // sparse selfadjointView with dense matrices
302     VERIFY_IS_APPROX(x=mUp.template selfadjointView<Upper>()*b, refX=refS*b);
303     VERIFY_IS_APPROX(x=mLo.template selfadjointView<Lower>()*b, refX=refS*b);
304     VERIFY_IS_APPROX(x=mS.template selfadjointView<Upper|Lower>()*b, refX=refS*b);
305 
306     VERIFY_IS_APPROX(x=b * mUp.template selfadjointView<Upper>(),       refX=b*refS);
307     VERIFY_IS_APPROX(x=b * mLo.template selfadjointView<Lower>(),       refX=b*refS);
308     VERIFY_IS_APPROX(x=b * mS.template selfadjointView<Upper|Lower>(),  refX=b*refS);
309 
310     VERIFY_IS_APPROX(x.noalias()+=mUp.template selfadjointView<Upper>()*b, refX+=refS*b);
311     VERIFY_IS_APPROX(x.noalias()-=mLo.template selfadjointView<Lower>()*b, refX-=refS*b);
312     VERIFY_IS_APPROX(x.noalias()+=mS.template selfadjointView<Upper|Lower>()*b, refX+=refS*b);
313 
314     // sparse selfadjointView with sparse matrices
315     SparseMatrixType mSres(rows,rows);
316     VERIFY_IS_APPROX(mSres = mLo.template selfadjointView<Lower>()*mS,
317                      refX = refLo.template selfadjointView<Lower>()*refS);
318     VERIFY_IS_APPROX(mSres = mS * mLo.template selfadjointView<Lower>(),
319                      refX = refS * refLo.template selfadjointView<Lower>());
320 
321     // sparse triangularView with dense matrices
322     VERIFY_IS_APPROX(x=mA.template triangularView<Upper>()*b, refX=refA.template triangularView<Upper>()*b);
323     VERIFY_IS_APPROX(x=mA.template triangularView<Lower>()*b, refX=refA.template triangularView<Lower>()*b);
324     VERIFY_IS_APPROX(x=b*mA.template triangularView<Upper>(), refX=b*refA.template triangularView<Upper>());
325     VERIFY_IS_APPROX(x=b*mA.template triangularView<Lower>(), refX=b*refA.template triangularView<Lower>());
326 
327     // sparse triangularView with sparse matrices
328     VERIFY_IS_APPROX(mSres = mA.template triangularView<Lower>()*mS,   refX = refA.template triangularView<Lower>()*refS);
329     VERIFY_IS_APPROX(mSres = mS * mA.template triangularView<Lower>(), refX = refS * refA.template triangularView<Lower>());
330     VERIFY_IS_APPROX(mSres = mA.template triangularView<Upper>()*mS,   refX = refA.template triangularView<Upper>()*refS);
331     VERIFY_IS_APPROX(mSres = mS * mA.template triangularView<Upper>(), refX = refS * refA.template triangularView<Upper>());
332   }
333 }
334 
335 // New test for Bug in SparseTimeDenseProduct
sparse_product_regression_test()336 template<typename SparseMatrixType, typename DenseMatrixType> void sparse_product_regression_test()
337 {
338   // This code does not compile with afflicted versions of the bug
339   SparseMatrixType sm1(3,2);
340   DenseMatrixType m2(2,2);
341   sm1.setZero();
342   m2.setZero();
343 
344   DenseMatrixType m3 = sm1*m2;
345 
346 
347   // This code produces a segfault with afflicted versions of another SparseTimeDenseProduct
348   // bug
349 
350   SparseMatrixType sm2(20000,2);
351   sm2.setZero();
352   DenseMatrixType m4(sm2*m2);
353 
354   VERIFY_IS_APPROX( m4(0,0), 0.0 );
355 }
356 
357 template<typename Scalar>
bug_942()358 void bug_942()
359 {
360   typedef Matrix<Scalar, Dynamic, 1>     Vector;
361   typedef SparseMatrix<Scalar, ColMajor> ColSpMat;
362   typedef SparseMatrix<Scalar, RowMajor> RowSpMat;
363   ColSpMat cmA(1,1);
364   cmA.insert(0,0) = 1;
365 
366   RowSpMat rmA(1,1);
367   rmA.insert(0,0) = 1;
368 
369   Vector d(1);
370   d[0] = 2;
371 
372   double res = 2;
373 
374   VERIFY_IS_APPROX( ( cmA*d.asDiagonal() ).eval().coeff(0,0), res );
375   VERIFY_IS_APPROX( ( d.asDiagonal()*rmA ).eval().coeff(0,0), res );
376   VERIFY_IS_APPROX( ( rmA*d.asDiagonal() ).eval().coeff(0,0), res );
377   VERIFY_IS_APPROX( ( d.asDiagonal()*cmA ).eval().coeff(0,0), res );
378 }
379 
380 template<typename Real>
test_mixing_types()381 void test_mixing_types()
382 {
383   typedef std::complex<Real> Cplx;
384   typedef SparseMatrix<Real> SpMatReal;
385   typedef SparseMatrix<Cplx> SpMatCplx;
386   typedef SparseMatrix<Cplx,RowMajor> SpRowMatCplx;
387   typedef Matrix<Real,Dynamic,Dynamic> DenseMatReal;
388   typedef Matrix<Cplx,Dynamic,Dynamic> DenseMatCplx;
389 
390   Index n = internal::random<Index>(1,100);
391   double density = (std::max)(8./(n*n), 0.2);
392 
393   SpMatReal sR1(n,n);
394   SpMatCplx sC1(n,n), sC2(n,n), sC3(n,n);
395   SpRowMatCplx sCR(n,n);
396   DenseMatReal dR1(n,n);
397   DenseMatCplx dC1(n,n), dC2(n,n), dC3(n,n);
398 
399   initSparse<Real>(density, dR1, sR1);
400   initSparse<Cplx>(density, dC1, sC1);
401   initSparse<Cplx>(density, dC2, sC2);
402 
403   VERIFY_IS_APPROX( sC2 = (sR1 * sC1),                         dC3 = dR1.template cast<Cplx>() * dC1 );
404   VERIFY_IS_APPROX( sC2 = (sC1 * sR1),                         dC3 = dC1 * dR1.template cast<Cplx>() );
405   VERIFY_IS_APPROX( sC2 = (sR1.transpose() * sC1),             dC3 = dR1.template cast<Cplx>().transpose() * dC1 );
406   VERIFY_IS_APPROX( sC2 = (sC1.transpose() * sR1),             dC3 = dC1.transpose() * dR1.template cast<Cplx>() );
407   VERIFY_IS_APPROX( sC2 = (sR1 * sC1.transpose()),             dC3 = dR1.template cast<Cplx>() * dC1.transpose() );
408   VERIFY_IS_APPROX( sC2 = (sC1 * sR1.transpose()),             dC3 = dC1 * dR1.template cast<Cplx>().transpose() );
409   VERIFY_IS_APPROX( sC2 = (sR1.transpose() * sC1.transpose()), dC3 = dR1.template cast<Cplx>().transpose() * dC1.transpose() );
410   VERIFY_IS_APPROX( sC2 = (sC1.transpose() * sR1.transpose()), dC3 = dC1.transpose() * dR1.template cast<Cplx>().transpose() );
411 
412   VERIFY_IS_APPROX( sCR = (sR1 * sC1),                         dC3 = dR1.template cast<Cplx>() * dC1 );
413   VERIFY_IS_APPROX( sCR = (sC1 * sR1),                         dC3 = dC1 * dR1.template cast<Cplx>() );
414   VERIFY_IS_APPROX( sCR = (sR1.transpose() * sC1),             dC3 = dR1.template cast<Cplx>().transpose() * dC1 );
415   VERIFY_IS_APPROX( sCR = (sC1.transpose() * sR1),             dC3 = dC1.transpose() * dR1.template cast<Cplx>() );
416   VERIFY_IS_APPROX( sCR = (sR1 * sC1.transpose()),             dC3 = dR1.template cast<Cplx>() * dC1.transpose() );
417   VERIFY_IS_APPROX( sCR = (sC1 * sR1.transpose()),             dC3 = dC1 * dR1.template cast<Cplx>().transpose() );
418   VERIFY_IS_APPROX( sCR = (sR1.transpose() * sC1.transpose()), dC3 = dR1.template cast<Cplx>().transpose() * dC1.transpose() );
419   VERIFY_IS_APPROX( sCR = (sC1.transpose() * sR1.transpose()), dC3 = dC1.transpose() * dR1.template cast<Cplx>().transpose() );
420 
421 
422   VERIFY_IS_APPROX( sC2 = (sR1 * sC1).pruned(),                         dC3 = dR1.template cast<Cplx>() * dC1 );
423   VERIFY_IS_APPROX( sC2 = (sC1 * sR1).pruned(),                         dC3 = dC1 * dR1.template cast<Cplx>() );
424   VERIFY_IS_APPROX( sC2 = (sR1.transpose() * sC1).pruned(),             dC3 = dR1.template cast<Cplx>().transpose() * dC1 );
425   VERIFY_IS_APPROX( sC2 = (sC1.transpose() * sR1).pruned(),             dC3 = dC1.transpose() * dR1.template cast<Cplx>() );
426   VERIFY_IS_APPROX( sC2 = (sR1 * sC1.transpose()).pruned(),             dC3 = dR1.template cast<Cplx>() * dC1.transpose() );
427   VERIFY_IS_APPROX( sC2 = (sC1 * sR1.transpose()).pruned(),             dC3 = dC1 * dR1.template cast<Cplx>().transpose() );
428   VERIFY_IS_APPROX( sC2 = (sR1.transpose() * sC1.transpose()).pruned(), dC3 = dR1.template cast<Cplx>().transpose() * dC1.transpose() );
429   VERIFY_IS_APPROX( sC2 = (sC1.transpose() * sR1.transpose()).pruned(), dC3 = dC1.transpose() * dR1.template cast<Cplx>().transpose() );
430 
431   VERIFY_IS_APPROX( sCR = (sR1 * sC1).pruned(),                         dC3 = dR1.template cast<Cplx>() * dC1 );
432   VERIFY_IS_APPROX( sCR = (sC1 * sR1).pruned(),                         dC3 = dC1 * dR1.template cast<Cplx>() );
433   VERIFY_IS_APPROX( sCR = (sR1.transpose() * sC1).pruned(),             dC3 = dR1.template cast<Cplx>().transpose() * dC1 );
434   VERIFY_IS_APPROX( sCR = (sC1.transpose() * sR1).pruned(),             dC3 = dC1.transpose() * dR1.template cast<Cplx>() );
435   VERIFY_IS_APPROX( sCR = (sR1 * sC1.transpose()).pruned(),             dC3 = dR1.template cast<Cplx>() * dC1.transpose() );
436   VERIFY_IS_APPROX( sCR = (sC1 * sR1.transpose()).pruned(),             dC3 = dC1 * dR1.template cast<Cplx>().transpose() );
437   VERIFY_IS_APPROX( sCR = (sR1.transpose() * sC1.transpose()).pruned(), dC3 = dR1.template cast<Cplx>().transpose() * dC1.transpose() );
438   VERIFY_IS_APPROX( sCR = (sC1.transpose() * sR1.transpose()).pruned(), dC3 = dC1.transpose() * dR1.template cast<Cplx>().transpose() );
439 
440 
441   VERIFY_IS_APPROX( dC2 = (sR1 * sC1),                         dC3 = dR1.template cast<Cplx>() * dC1 );
442   VERIFY_IS_APPROX( dC2 = (sC1 * sR1),                         dC3 = dC1 * dR1.template cast<Cplx>() );
443   VERIFY_IS_APPROX( dC2 = (sR1.transpose() * sC1),             dC3 = dR1.template cast<Cplx>().transpose() * dC1 );
444   VERIFY_IS_APPROX( dC2 = (sC1.transpose() * sR1),             dC3 = dC1.transpose() * dR1.template cast<Cplx>() );
445   VERIFY_IS_APPROX( dC2 = (sR1 * sC1.transpose()),             dC3 = dR1.template cast<Cplx>() * dC1.transpose() );
446   VERIFY_IS_APPROX( dC2 = (sC1 * sR1.transpose()),             dC3 = dC1 * dR1.template cast<Cplx>().transpose() );
447   VERIFY_IS_APPROX( dC2 = (sR1.transpose() * sC1.transpose()), dC3 = dR1.template cast<Cplx>().transpose() * dC1.transpose() );
448   VERIFY_IS_APPROX( dC2 = (sC1.transpose() * sR1.transpose()), dC3 = dC1.transpose() * dR1.template cast<Cplx>().transpose() );
449 
450 
451   VERIFY_IS_APPROX( dC2 = dR1 * sC1, dC3 = dR1.template cast<Cplx>() * sC1 );
452   VERIFY_IS_APPROX( dC2 = sR1 * dC1, dC3 = sR1.template cast<Cplx>() * dC1 );
453   VERIFY_IS_APPROX( dC2 = dC1 * sR1, dC3 = dC1 * sR1.template cast<Cplx>() );
454   VERIFY_IS_APPROX( dC2 = sC1 * dR1, dC3 = sC1 * dR1.template cast<Cplx>() );
455 
456   VERIFY_IS_APPROX( dC2 = dR1.row(0) * sC1, dC3 = dR1.template cast<Cplx>().row(0) * sC1 );
457   VERIFY_IS_APPROX( dC2 = sR1 * dC1.col(0), dC3 = sR1.template cast<Cplx>() * dC1.col(0) );
458   VERIFY_IS_APPROX( dC2 = dC1.row(0) * sR1, dC3 = dC1.row(0) * sR1.template cast<Cplx>() );
459   VERIFY_IS_APPROX( dC2 = sC1 * dR1.col(0), dC3 = sC1 * dR1.template cast<Cplx>().col(0) );
460 }
461 
test_sparse_product()462 void test_sparse_product()
463 {
464   for(int i = 0; i < g_repeat; i++) {
465     CALL_SUBTEST_1( (sparse_product<SparseMatrix<double,ColMajor> >()) );
466     CALL_SUBTEST_1( (sparse_product<SparseMatrix<double,RowMajor> >()) );
467     CALL_SUBTEST_1( (bug_942<double>()) );
468     CALL_SUBTEST_2( (sparse_product<SparseMatrix<std::complex<double>, ColMajor > >()) );
469     CALL_SUBTEST_2( (sparse_product<SparseMatrix<std::complex<double>, RowMajor > >()) );
470     CALL_SUBTEST_3( (sparse_product<SparseMatrix<float,ColMajor,long int> >()) );
471     CALL_SUBTEST_4( (sparse_product_regression_test<SparseMatrix<double,RowMajor>, Matrix<double, Dynamic, Dynamic, RowMajor> >()) );
472 
473     CALL_SUBTEST_5( (test_mixing_types<float>()) );
474   }
475 }
476