sparse_permutations.cpp 9.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2011-2015 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. static long int nb_transposed_copies;
  10. #define EIGEN_SPARSE_TRANSPOSED_COPY_PLUGIN {nb_transposed_copies++;}
  11. #define VERIFY_TRANSPOSITION_COUNT(XPR,N) {\
  12. nb_transposed_copies = 0; \
  13. XPR; \
  14. if(nb_transposed_copies!=N) std::cerr << "nb_transposed_copies == " << nb_transposed_copies << "\n"; \
  15. VERIFY( (#XPR) && nb_transposed_copies==N ); \
  16. }
  17. #include "sparse.h"
  18. template<typename T>
  19. bool is_sorted(const T& mat) {
  20. for(Index k = 0; k<mat.outerSize(); ++k)
  21. {
  22. Index prev = -1;
  23. for(typename T::InnerIterator it(mat,k); it; ++it)
  24. {
  25. if(prev>=it.index())
  26. return false;
  27. prev = it.index();
  28. }
  29. }
  30. return true;
  31. }
  32. template<typename T>
  33. typename internal::nested_eval<T,1>::type eval(const T &xpr)
  34. {
  35. VERIFY( int(internal::nested_eval<T,1>::type::Flags&RowMajorBit) == int(internal::evaluator<T>::Flags&RowMajorBit) );
  36. return xpr;
  37. }
  38. template<int OtherStorage, typename SparseMatrixType> void sparse_permutations(const SparseMatrixType& ref)
  39. {
  40. const Index rows = ref.rows();
  41. const Index cols = ref.cols();
  42. typedef typename SparseMatrixType::Scalar Scalar;
  43. typedef typename SparseMatrixType::StorageIndex StorageIndex;
  44. typedef SparseMatrix<Scalar, OtherStorage, StorageIndex> OtherSparseMatrixType;
  45. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  46. typedef Matrix<StorageIndex,Dynamic,1> VectorI;
  47. // bool IsRowMajor1 = SparseMatrixType::IsRowMajor;
  48. // bool IsRowMajor2 = OtherSparseMatrixType::IsRowMajor;
  49. double density = (std::max)(8./(rows*cols), 0.01);
  50. SparseMatrixType mat(rows, cols), up(rows,cols), lo(rows,cols);
  51. OtherSparseMatrixType res;
  52. DenseMatrix mat_d = DenseMatrix::Zero(rows, cols), up_sym_d, lo_sym_d, res_d;
  53. initSparse<Scalar>(density, mat_d, mat, 0);
  54. up = mat.template triangularView<Upper>();
  55. lo = mat.template triangularView<Lower>();
  56. up_sym_d = mat_d.template selfadjointView<Upper>();
  57. lo_sym_d = mat_d.template selfadjointView<Lower>();
  58. VERIFY_IS_APPROX(mat, mat_d);
  59. VERIFY_IS_APPROX(up, DenseMatrix(mat_d.template triangularView<Upper>()));
  60. VERIFY_IS_APPROX(lo, DenseMatrix(mat_d.template triangularView<Lower>()));
  61. PermutationMatrix<Dynamic> p, p_null;
  62. VectorI pi;
  63. randomPermutationVector(pi, cols);
  64. p.indices() = pi;
  65. VERIFY( is_sorted( ::eval(mat*p) ));
  66. VERIFY( is_sorted( res = mat*p ));
  67. VERIFY_TRANSPOSITION_COUNT( ::eval(mat*p), 0);
  68. //VERIFY_TRANSPOSITION_COUNT( res = mat*p, IsRowMajor ? 1 : 0 );
  69. res_d = mat_d*p;
  70. VERIFY(res.isApprox(res_d) && "mat*p");
  71. VERIFY( is_sorted( ::eval(p*mat) ));
  72. VERIFY( is_sorted( res = p*mat ));
  73. VERIFY_TRANSPOSITION_COUNT( ::eval(p*mat), 0);
  74. res_d = p*mat_d;
  75. VERIFY(res.isApprox(res_d) && "p*mat");
  76. VERIFY( is_sorted( (mat*p).eval() ));
  77. VERIFY( is_sorted( res = mat*p.inverse() ));
  78. VERIFY_TRANSPOSITION_COUNT( ::eval(mat*p.inverse()), 0);
  79. res_d = mat*p.inverse();
  80. VERIFY(res.isApprox(res_d) && "mat*inv(p)");
  81. VERIFY( is_sorted( (p*mat+p*mat).eval() ));
  82. VERIFY( is_sorted( res = p.inverse()*mat ));
  83. VERIFY_TRANSPOSITION_COUNT( ::eval(p.inverse()*mat), 0);
  84. res_d = p.inverse()*mat_d;
  85. VERIFY(res.isApprox(res_d) && "inv(p)*mat");
  86. VERIFY( is_sorted( (p * mat * p.inverse()).eval() ));
  87. VERIFY( is_sorted( res = mat.twistedBy(p) ));
  88. VERIFY_TRANSPOSITION_COUNT( ::eval(p * mat * p.inverse()), 0);
  89. res_d = (p * mat_d) * p.inverse();
  90. VERIFY(res.isApprox(res_d) && "p*mat*inv(p)");
  91. VERIFY( is_sorted( res = mat.template selfadjointView<Upper>().twistedBy(p_null) ));
  92. res_d = up_sym_d;
  93. VERIFY(res.isApprox(res_d) && "full selfadjoint upper to full");
  94. VERIFY( is_sorted( res = mat.template selfadjointView<Lower>().twistedBy(p_null) ));
  95. res_d = lo_sym_d;
  96. VERIFY(res.isApprox(res_d) && "full selfadjoint lower to full");
  97. VERIFY( is_sorted( res = up.template selfadjointView<Upper>().twistedBy(p_null) ));
  98. res_d = up_sym_d;
  99. VERIFY(res.isApprox(res_d) && "upper selfadjoint to full");
  100. VERIFY( is_sorted( res = lo.template selfadjointView<Lower>().twistedBy(p_null) ));
  101. res_d = lo_sym_d;
  102. VERIFY(res.isApprox(res_d) && "lower selfadjoint full");
  103. VERIFY( is_sorted( res = mat.template selfadjointView<Upper>() ));
  104. res_d = up_sym_d;
  105. VERIFY(res.isApprox(res_d) && "full selfadjoint upper to full");
  106. VERIFY( is_sorted( res = mat.template selfadjointView<Lower>() ));
  107. res_d = lo_sym_d;
  108. VERIFY(res.isApprox(res_d) && "full selfadjoint lower to full");
  109. VERIFY( is_sorted( res = up.template selfadjointView<Upper>() ));
  110. res_d = up_sym_d;
  111. VERIFY(res.isApprox(res_d) && "upper selfadjoint to full");
  112. VERIFY( is_sorted( res = lo.template selfadjointView<Lower>() ));
  113. res_d = lo_sym_d;
  114. VERIFY(res.isApprox(res_d) && "lower selfadjoint full");
  115. res.template selfadjointView<Upper>() = mat.template selfadjointView<Upper>();
  116. res_d = up_sym_d.template triangularView<Upper>();
  117. VERIFY(res.isApprox(res_d) && "full selfadjoint upper to upper");
  118. res.template selfadjointView<Lower>() = mat.template selfadjointView<Upper>();
  119. res_d = up_sym_d.template triangularView<Lower>();
  120. VERIFY(res.isApprox(res_d) && "full selfadjoint upper to lower");
  121. res.template selfadjointView<Upper>() = mat.template selfadjointView<Lower>();
  122. res_d = lo_sym_d.template triangularView<Upper>();
  123. VERIFY(res.isApprox(res_d) && "full selfadjoint lower to upper");
  124. res.template selfadjointView<Lower>() = mat.template selfadjointView<Lower>();
  125. res_d = lo_sym_d.template triangularView<Lower>();
  126. VERIFY(res.isApprox(res_d) && "full selfadjoint lower to lower");
  127. res.template selfadjointView<Upper>() = mat.template selfadjointView<Upper>().twistedBy(p);
  128. res_d = ((p * up_sym_d) * p.inverse()).eval().template triangularView<Upper>();
  129. VERIFY(res.isApprox(res_d) && "full selfadjoint upper twisted to upper");
  130. res.template selfadjointView<Upper>() = mat.template selfadjointView<Lower>().twistedBy(p);
  131. res_d = ((p * lo_sym_d) * p.inverse()).eval().template triangularView<Upper>();
  132. VERIFY(res.isApprox(res_d) && "full selfadjoint lower twisted to upper");
  133. res.template selfadjointView<Lower>() = mat.template selfadjointView<Lower>().twistedBy(p);
  134. res_d = ((p * lo_sym_d) * p.inverse()).eval().template triangularView<Lower>();
  135. VERIFY(res.isApprox(res_d) && "full selfadjoint lower twisted to lower");
  136. res.template selfadjointView<Lower>() = mat.template selfadjointView<Upper>().twistedBy(p);
  137. res_d = ((p * up_sym_d) * p.inverse()).eval().template triangularView<Lower>();
  138. VERIFY(res.isApprox(res_d) && "full selfadjoint upper twisted to lower");
  139. res.template selfadjointView<Upper>() = up.template selfadjointView<Upper>().twistedBy(p);
  140. res_d = ((p * up_sym_d) * p.inverse()).eval().template triangularView<Upper>();
  141. VERIFY(res.isApprox(res_d) && "upper selfadjoint twisted to upper");
  142. res.template selfadjointView<Upper>() = lo.template selfadjointView<Lower>().twistedBy(p);
  143. res_d = ((p * lo_sym_d) * p.inverse()).eval().template triangularView<Upper>();
  144. VERIFY(res.isApprox(res_d) && "lower selfadjoint twisted to upper");
  145. res.template selfadjointView<Lower>() = lo.template selfadjointView<Lower>().twistedBy(p);
  146. res_d = ((p * lo_sym_d) * p.inverse()).eval().template triangularView<Lower>();
  147. VERIFY(res.isApprox(res_d) && "lower selfadjoint twisted to lower");
  148. res.template selfadjointView<Lower>() = up.template selfadjointView<Upper>().twistedBy(p);
  149. res_d = ((p * up_sym_d) * p.inverse()).eval().template triangularView<Lower>();
  150. VERIFY(res.isApprox(res_d) && "upper selfadjoint twisted to lower");
  151. VERIFY( is_sorted( res = mat.template selfadjointView<Upper>().twistedBy(p) ));
  152. res_d = (p * up_sym_d) * p.inverse();
  153. VERIFY(res.isApprox(res_d) && "full selfadjoint upper twisted to full");
  154. VERIFY( is_sorted( res = mat.template selfadjointView<Lower>().twistedBy(p) ));
  155. res_d = (p * lo_sym_d) * p.inverse();
  156. VERIFY(res.isApprox(res_d) && "full selfadjoint lower twisted to full");
  157. VERIFY( is_sorted( res = up.template selfadjointView<Upper>().twistedBy(p) ));
  158. res_d = (p * up_sym_d) * p.inverse();
  159. VERIFY(res.isApprox(res_d) && "upper selfadjoint twisted to full");
  160. VERIFY( is_sorted( res = lo.template selfadjointView<Lower>().twistedBy(p) ));
  161. res_d = (p * lo_sym_d) * p.inverse();
  162. VERIFY(res.isApprox(res_d) && "lower selfadjoint twisted to full");
  163. }
  164. template<typename Scalar> void sparse_permutations_all(int size)
  165. {
  166. CALL_SUBTEST(( sparse_permutations<ColMajor>(SparseMatrix<Scalar, ColMajor>(size,size)) ));
  167. CALL_SUBTEST(( sparse_permutations<ColMajor>(SparseMatrix<Scalar, RowMajor>(size,size)) ));
  168. CALL_SUBTEST(( sparse_permutations<RowMajor>(SparseMatrix<Scalar, ColMajor>(size,size)) ));
  169. CALL_SUBTEST(( sparse_permutations<RowMajor>(SparseMatrix<Scalar, RowMajor>(size,size)) ));
  170. }
  171. void test_sparse_permutations()
  172. {
  173. for(int i = 0; i < g_repeat; i++) {
  174. int s = Eigen::internal::random<int>(1,50);
  175. CALL_SUBTEST_1(( sparse_permutations_all<double>(s) ));
  176. CALL_SUBTEST_2(( sparse_permutations_all<std::complex<double> >(s) ));
  177. }
  178. VERIFY((internal::is_same<internal::permutation_matrix_product<SparseMatrix<double>,OnTheRight,false,SparseShape>::ReturnType,
  179. internal::nested_eval<Product<SparseMatrix<double>,PermutationMatrix<Dynamic,Dynamic>,AliasFreeProduct>,1>::type>::value));
  180. VERIFY((internal::is_same<internal::permutation_matrix_product<SparseMatrix<double>,OnTheLeft,false,SparseShape>::ReturnType,
  181. internal::nested_eval<Product<PermutationMatrix<Dynamic,Dynamic>,SparseMatrix<double>,AliasFreeProduct>,1>::type>::value));
  182. }