• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 20015 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 // This unit test cannot be easily written to work with EIGEN_DEFAULT_TO_ROW_MAJOR
11 #ifdef EIGEN_DEFAULT_TO_ROW_MAJOR
12 #undef EIGEN_DEFAULT_TO_ROW_MAJOR
13 #endif
14 
15 static long int nb_temporaries;
16 
on_temporary_creation()17 inline void on_temporary_creation() {
18   // here's a great place to set a breakpoint when debugging failures in this test!
19   nb_temporaries++;
20 }
21 
22 #define EIGEN_SPARSE_CREATE_TEMPORARY_PLUGIN { on_temporary_creation(); }
23 
24 #include "main.h"
25 #include <Eigen/SparseCore>
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 
check_const_correctness(const PlainObjectType &)34 template<typename PlainObjectType> void check_const_correctness(const PlainObjectType&)
35 {
36   // verify that ref-to-const don't have LvalueBit
37   typedef typename internal::add_const<PlainObjectType>::type ConstPlainObjectType;
38   VERIFY( !(internal::traits<Ref<ConstPlainObjectType> >::Flags & LvalueBit) );
39   VERIFY( !(internal::traits<Ref<ConstPlainObjectType, Aligned> >::Flags & LvalueBit) );
40   VERIFY( !(Ref<ConstPlainObjectType>::Flags & LvalueBit) );
41   VERIFY( !(Ref<ConstPlainObjectType, Aligned>::Flags & LvalueBit) );
42 }
43 
44 template<typename B>
call_ref_1(Ref<SparseMatrix<float>> a,const B & b)45 EIGEN_DONT_INLINE void call_ref_1(Ref<SparseMatrix<float> > a, const B &b) { VERIFY_IS_EQUAL(a.toDense(),b.toDense()); }
46 
47 template<typename B>
call_ref_2(const Ref<const SparseMatrix<float>> & a,const B & b)48 EIGEN_DONT_INLINE void call_ref_2(const Ref<const SparseMatrix<float> >& a, const B &b) { VERIFY_IS_EQUAL(a.toDense(),b.toDense()); }
49 
50 template<typename B>
call_ref_3(const Ref<const SparseMatrix<float>,StandardCompressedFormat> & a,const B & b)51 EIGEN_DONT_INLINE void call_ref_3(const Ref<const SparseMatrix<float>, StandardCompressedFormat>& a, const B &b) {
52   VERIFY(a.isCompressed());
53   VERIFY_IS_EQUAL(a.toDense(),b.toDense());
54 }
55 
56 template<typename B>
call_ref_4(Ref<SparseVector<float>> a,const B & b)57 EIGEN_DONT_INLINE void call_ref_4(Ref<SparseVector<float> > a, const B &b) { VERIFY_IS_EQUAL(a.toDense(),b.toDense()); }
58 
59 template<typename B>
call_ref_5(const Ref<const SparseVector<float>> & a,const B & b)60 EIGEN_DONT_INLINE void call_ref_5(const Ref<const SparseVector<float> >& a, const B &b) { VERIFY_IS_EQUAL(a.toDense(),b.toDense()); }
61 
call_ref()62 void call_ref()
63 {
64   SparseMatrix<float>               A = MatrixXf::Random(10,10).sparseView(0.5,1);
65   SparseMatrix<float,RowMajor>      B = MatrixXf::Random(10,10).sparseView(0.5,1);
66   SparseMatrix<float>               C = MatrixXf::Random(10,10).sparseView(0.5,1);
67   C.reserve(VectorXi::Constant(C.outerSize(), 2));
68   const SparseMatrix<float>&        Ac(A);
69   Block<SparseMatrix<float> >       Ab(A,0,1, 3,3);
70   const Block<SparseMatrix<float> > Abc(A,0,1,3,3);
71   SparseVector<float>               vc =  VectorXf::Random(10).sparseView(0.5,1);
72   SparseVector<float,RowMajor>      vr =  VectorXf::Random(10).sparseView(0.5,1);
73   SparseMatrix<float> AA = A*A;
74 
75 
76   VERIFY_EVALUATION_COUNT( call_ref_1(A, A),  0);
77 //   VERIFY_EVALUATION_COUNT( call_ref_1(Ac, Ac),  0); // does not compile on purpose
78   VERIFY_EVALUATION_COUNT( call_ref_2(A, A),  0);
79   VERIFY_EVALUATION_COUNT( call_ref_3(A, A),  0);
80   VERIFY_EVALUATION_COUNT( call_ref_2(A.transpose(), A.transpose()),  1);
81   VERIFY_EVALUATION_COUNT( call_ref_3(A.transpose(), A.transpose()),  1);
82   VERIFY_EVALUATION_COUNT( call_ref_2(Ac,Ac), 0);
83   VERIFY_EVALUATION_COUNT( call_ref_3(Ac,Ac), 0);
84   VERIFY_EVALUATION_COUNT( call_ref_2(A+A,2*Ac), 1);
85   VERIFY_EVALUATION_COUNT( call_ref_3(A+A,2*Ac), 1);
86   VERIFY_EVALUATION_COUNT( call_ref_2(B, B),  1);
87   VERIFY_EVALUATION_COUNT( call_ref_3(B, B),  1);
88   VERIFY_EVALUATION_COUNT( call_ref_2(B.transpose(), B.transpose()),  0);
89   VERIFY_EVALUATION_COUNT( call_ref_3(B.transpose(), B.transpose()),  0);
90   VERIFY_EVALUATION_COUNT( call_ref_2(A*A, AA),  3);
91   VERIFY_EVALUATION_COUNT( call_ref_3(A*A, AA),  3);
92 
93   VERIFY(!C.isCompressed());
94   VERIFY_EVALUATION_COUNT( call_ref_3(C, C),  1);
95 
96   Ref<SparseMatrix<float> > Ar(A);
97   VERIFY_IS_APPROX(Ar+Ar, A+A);
98   VERIFY_EVALUATION_COUNT( call_ref_1(Ar, A),  0);
99   VERIFY_EVALUATION_COUNT( call_ref_2(Ar, A),  0);
100 
101   Ref<SparseMatrix<float,RowMajor> > Br(B);
102   VERIFY_EVALUATION_COUNT( call_ref_1(Br.transpose(), Br.transpose()),  0);
103   VERIFY_EVALUATION_COUNT( call_ref_2(Br, Br),  1);
104   VERIFY_EVALUATION_COUNT( call_ref_2(Br.transpose(), Br.transpose()),  0);
105 
106   Ref<const SparseMatrix<float> > Arc(A);
107 //   VERIFY_EVALUATION_COUNT( call_ref_1(Arc, Arc),  0); // does not compile on purpose
108   VERIFY_EVALUATION_COUNT( call_ref_2(Arc, Arc),  0);
109 
110   VERIFY_EVALUATION_COUNT( call_ref_2(A.middleCols(1,3), A.middleCols(1,3)),  0);
111 
112   VERIFY_EVALUATION_COUNT( call_ref_2(A.col(2), A.col(2)),  0);
113   VERIFY_EVALUATION_COUNT( call_ref_2(vc, vc),  0);
114   VERIFY_EVALUATION_COUNT( call_ref_2(vr.transpose(), vr.transpose()),  0);
115   VERIFY_EVALUATION_COUNT( call_ref_2(vr, vr.transpose()),  0);
116 
117   VERIFY_EVALUATION_COUNT( call_ref_2(A.block(1,1,3,3), A.block(1,1,3,3)),  1); // should be 0 (allocate starts/nnz only)
118 
119   VERIFY_EVALUATION_COUNT( call_ref_4(vc, vc),  0);
120   VERIFY_EVALUATION_COUNT( call_ref_4(vr, vr.transpose()),  0);
121   VERIFY_EVALUATION_COUNT( call_ref_5(vc, vc),  0);
122   VERIFY_EVALUATION_COUNT( call_ref_5(vr, vr.transpose()),  0);
123   VERIFY_EVALUATION_COUNT( call_ref_4(A.col(2), A.col(2)),  0);
124   VERIFY_EVALUATION_COUNT( call_ref_5(A.col(2), A.col(2)),  0);
125   // VERIFY_EVALUATION_COUNT( call_ref_4(A.row(2), A.row(2).transpose()),  1); // does not compile on purpose
126   VERIFY_EVALUATION_COUNT( call_ref_5(A.row(2), A.row(2).transpose()),  1);
127 }
128 
test_sparse_ref()129 void test_sparse_ref()
130 {
131   for(int i = 0; i < g_repeat; i++) {
132     CALL_SUBTEST_1( check_const_correctness(SparseMatrix<float>()) );
133     CALL_SUBTEST_1( check_const_correctness(SparseMatrix<double,RowMajor>()) );
134     CALL_SUBTEST_2( call_ref() );
135 
136     CALL_SUBTEST_3( check_const_correctness(SparseVector<float>()) );
137     CALL_SUBTEST_3( check_const_correctness(SparseVector<double,RowMajor>()) );
138   }
139 }
140