• 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) 2008-2010 Gael Guennebaud <g.gael@free.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 
11 // import basic and product tests for deprectaed DynamicSparseMatrix
12 #define EIGEN_NO_DEPRECATED_WARNING
13 #include "sparse_basic.cpp"
14 #include "sparse_product.cpp"
15 #include <Eigen/SparseExtra>
16 
17 template<typename SetterType,typename DenseType, typename Scalar, int Options>
test_random_setter(SparseMatrix<Scalar,Options> & sm,const DenseType & ref,const std::vector<Vector2i> & nonzeroCoords)18 bool test_random_setter(SparseMatrix<Scalar,Options>& sm, const DenseType& ref, const std::vector<Vector2i>& nonzeroCoords)
19 {
20   {
21     sm.setZero();
22     SetterType w(sm);
23     std::vector<Vector2i> remaining = nonzeroCoords;
24     while(!remaining.empty())
25     {
26       int i = internal::random<int>(0,static_cast<int>(remaining.size())-1);
27       w(remaining[i].x(),remaining[i].y()) = ref.coeff(remaining[i].x(),remaining[i].y());
28       remaining[i] = remaining.back();
29       remaining.pop_back();
30     }
31   }
32   return sm.isApprox(ref);
33 }
34 
35 template<typename SetterType,typename DenseType, typename T>
test_random_setter(DynamicSparseMatrix<T> & sm,const DenseType & ref,const std::vector<Vector2i> & nonzeroCoords)36 bool test_random_setter(DynamicSparseMatrix<T>& sm, const DenseType& ref, const std::vector<Vector2i>& nonzeroCoords)
37 {
38   sm.setZero();
39   std::vector<Vector2i> remaining = nonzeroCoords;
40   while(!remaining.empty())
41   {
42     int i = internal::random<int>(0,static_cast<int>(remaining.size())-1);
43     sm.coeffRef(remaining[i].x(),remaining[i].y()) = ref.coeff(remaining[i].x(),remaining[i].y());
44     remaining[i] = remaining.back();
45     remaining.pop_back();
46   }
47   return sm.isApprox(ref);
48 }
49 
sparse_extra(const SparseMatrixType & ref)50 template<typename SparseMatrixType> void sparse_extra(const SparseMatrixType& ref)
51 {
52   typedef typename SparseMatrixType::Index Index;
53   const Index rows = ref.rows();
54   const Index cols = ref.cols();
55   typedef typename SparseMatrixType::Scalar Scalar;
56   enum { Flags = SparseMatrixType::Flags };
57 
58   double density = (std::max)(8./(rows*cols), 0.01);
59   typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
60   typedef Matrix<Scalar,Dynamic,1> DenseVector;
61   Scalar eps = 1e-6;
62 
63   SparseMatrixType m(rows, cols);
64   DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
65   DenseVector vec1 = DenseVector::Random(rows);
66 
67   std::vector<Vector2i> zeroCoords;
68   std::vector<Vector2i> nonzeroCoords;
69   initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords);
70 
71   if (zeroCoords.size()==0 || nonzeroCoords.size()==0)
72     return;
73 
74   // test coeff and coeffRef
75   for (int i=0; i<(int)zeroCoords.size(); ++i)
76   {
77     VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps );
78     if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value)
79       VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[0].x(),zeroCoords[0].y()) = 5 );
80   }
81   VERIFY_IS_APPROX(m, refMat);
82 
83   m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
84   refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
85 
86   VERIFY_IS_APPROX(m, refMat);
87 
88   // random setter
89 //   {
90 //     m.setZero();
91 //     VERIFY_IS_NOT_APPROX(m, refMat);
92 //     SparseSetter<SparseMatrixType, RandomAccessPattern> w(m);
93 //     std::vector<Vector2i> remaining = nonzeroCoords;
94 //     while(!remaining.empty())
95 //     {
96 //       int i = internal::random<int>(0,remaining.size()-1);
97 //       w->coeffRef(remaining[i].x(),remaining[i].y()) = refMat.coeff(remaining[i].x(),remaining[i].y());
98 //       remaining[i] = remaining.back();
99 //       remaining.pop_back();
100 //     }
101 //   }
102 //   VERIFY_IS_APPROX(m, refMat);
103 
104     VERIFY(( test_random_setter<RandomSetter<SparseMatrixType, StdMapTraits> >(m,refMat,nonzeroCoords) ));
105     #ifdef EIGEN_UNORDERED_MAP_SUPPORT
106     VERIFY(( test_random_setter<RandomSetter<SparseMatrixType, StdUnorderedMapTraits> >(m,refMat,nonzeroCoords) ));
107     #endif
108     #ifdef _DENSE_HASH_MAP_H_
109     VERIFY(( test_random_setter<RandomSetter<SparseMatrixType, GoogleDenseHashMapTraits> >(m,refMat,nonzeroCoords) ));
110     #endif
111     #ifdef _SPARSE_HASH_MAP_H_
112     VERIFY(( test_random_setter<RandomSetter<SparseMatrixType, GoogleSparseHashMapTraits> >(m,refMat,nonzeroCoords) ));
113     #endif
114 
115 
116   // test RandomSetter
117   /*{
118     SparseMatrixType m1(rows,cols), m2(rows,cols);
119     DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
120     initSparse<Scalar>(density, refM1, m1);
121     {
122       Eigen::RandomSetter<SparseMatrixType > setter(m2);
123       for (int j=0; j<m1.outerSize(); ++j)
124         for (typename SparseMatrixType::InnerIterator i(m1,j); i; ++i)
125           setter(i.index(), j) = i.value();
126     }
127     VERIFY_IS_APPROX(m1, m2);
128   }*/
129 
130 
131 }
132 
test_sparse_extra()133 void test_sparse_extra()
134 {
135   for(int i = 0; i < g_repeat; i++) {
136     int s = Eigen::internal::random<int>(1,50);
137     CALL_SUBTEST_1( sparse_extra(SparseMatrix<double>(8, 8)) );
138     CALL_SUBTEST_2( sparse_extra(SparseMatrix<std::complex<double> >(s, s)) );
139     CALL_SUBTEST_1( sparse_extra(SparseMatrix<double>(s, s)) );
140 
141     CALL_SUBTEST_3( sparse_extra(DynamicSparseMatrix<double>(s, s)) );
142 //    CALL_SUBTEST_3(( sparse_basic(DynamicSparseMatrix<double>(s, s)) ));
143 //    CALL_SUBTEST_3(( sparse_basic(DynamicSparseMatrix<double,ColMajor,long int>(s, s)) ));
144 
145     CALL_SUBTEST_3( (sparse_product<DynamicSparseMatrix<float, ColMajor> >()) );
146     CALL_SUBTEST_3( (sparse_product<DynamicSparseMatrix<float, RowMajor> >()) );
147   }
148 }
149