You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

309 lines
9.2 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2011 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. #include "sparse.h"
  10. #include <Eigen/SparseCore>
  11. template<typename Solver, typename Rhs, typename DenseMat, typename DenseRhs>
  12. void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const DenseMat& dA, const DenseRhs& db)
  13. {
  14. typedef typename Solver::MatrixType Mat;
  15. typedef typename Mat::Scalar Scalar;
  16. DenseRhs refX = dA.lu().solve(db);
  17. Rhs x(b.rows(), b.cols());
  18. Rhs oldb = b;
  19. solver.compute(A);
  20. if (solver.info() != Success)
  21. {
  22. std::cerr << "sparse solver testing: factorization failed (check_sparse_solving)\n";
  23. exit(0);
  24. return;
  25. }
  26. x = solver.solve(b);
  27. if (solver.info() != Success)
  28. {
  29. std::cerr << "sparse solver testing: solving failed\n";
  30. return;
  31. }
  32. VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
  33. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  34. x.setZero();
  35. // test the analyze/factorize API
  36. solver.analyzePattern(A);
  37. solver.factorize(A);
  38. if (solver.info() != Success)
  39. {
  40. std::cerr << "sparse solver testing: factorization failed (check_sparse_solving)\n";
  41. exit(0);
  42. return;
  43. }
  44. x = solver.solve(b);
  45. if (solver.info() != Success)
  46. {
  47. std::cerr << "sparse solver testing: solving failed\n";
  48. return;
  49. }
  50. VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
  51. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  52. // test Block as the result and rhs:
  53. {
  54. DenseRhs x(db.rows(), db.cols());
  55. DenseRhs b(db), oldb(db);
  56. x.setZero();
  57. x.block(0,0,x.rows(),x.cols()) = solver.solve(b.block(0,0,b.rows(),b.cols()));
  58. VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
  59. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  60. }
  61. }
  62. template<typename Solver, typename Rhs>
  63. void check_sparse_solving_real_cases(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const Rhs& refX)
  64. {
  65. typedef typename Solver::MatrixType Mat;
  66. typedef typename Mat::Scalar Scalar;
  67. typedef typename Mat::RealScalar RealScalar;
  68. Rhs x(b.rows(), b.cols());
  69. solver.compute(A);
  70. if (solver.info() != Success)
  71. {
  72. std::cerr << "sparse solver testing: factorization failed (check_sparse_solving_real_cases)\n";
  73. exit(0);
  74. return;
  75. }
  76. x = solver.solve(b);
  77. if (solver.info() != Success)
  78. {
  79. std::cerr << "sparse solver testing: solving failed\n";
  80. return;
  81. }
  82. RealScalar res_error;
  83. // Compute the norm of the relative error
  84. if(refX.size() != 0)
  85. res_error = (refX - x).norm()/refX.norm();
  86. else
  87. {
  88. // Compute the relative residual norm
  89. res_error = (b - A * x).norm()/b.norm();
  90. }
  91. if (res_error > test_precision<Scalar>() ){
  92. std::cerr << "Test " << g_test_stack.back() << " failed in "EI_PP_MAKE_STRING(__FILE__)
  93. << " (" << EI_PP_MAKE_STRING(__LINE__) << ")" << std::endl << std::endl;
  94. abort();
  95. }
  96. }
  97. template<typename Solver, typename DenseMat>
  98. void check_sparse_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
  99. {
  100. typedef typename Solver::MatrixType Mat;
  101. typedef typename Mat::Scalar Scalar;
  102. typedef typename Mat::RealScalar RealScalar;
  103. solver.compute(A);
  104. if (solver.info() != Success)
  105. {
  106. std::cerr << "sparse solver testing: factorization failed (check_sparse_determinant)\n";
  107. return;
  108. }
  109. Scalar refDet = dA.determinant();
  110. VERIFY_IS_APPROX(refDet,solver.determinant());
  111. }
  112. template<typename Solver, typename DenseMat>
  113. int generate_sparse_spd_problem(Solver& , typename Solver::MatrixType& A, typename Solver::MatrixType& halfA, DenseMat& dA, int maxSize = 300)
  114. {
  115. typedef typename Solver::MatrixType Mat;
  116. typedef typename Mat::Scalar Scalar;
  117. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  118. int size = internal::random<int>(1,maxSize);
  119. double density = (std::max)(8./(size*size), 0.01);
  120. Mat M(size, size);
  121. DenseMatrix dM(size, size);
  122. initSparse<Scalar>(density, dM, M, ForceNonZeroDiag);
  123. A = M * M.adjoint();
  124. dA = dM * dM.adjoint();
  125. halfA.resize(size,size);
  126. halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
  127. return size;
  128. }
  129. #ifdef TEST_REAL_CASES
  130. template<typename Scalar>
  131. inline std::string get_matrixfolder()
  132. {
  133. std::string mat_folder = TEST_REAL_CASES;
  134. if( internal::is_same<Scalar, std::complex<float> >::value || internal::is_same<Scalar, std::complex<double> >::value )
  135. mat_folder = mat_folder + static_cast<string>("/complex/");
  136. else
  137. mat_folder = mat_folder + static_cast<string>("/real/");
  138. return mat_folder;
  139. }
  140. #endif
  141. template<typename Solver> void check_sparse_spd_solving(Solver& solver)
  142. {
  143. typedef typename Solver::MatrixType Mat;
  144. typedef typename Mat::Scalar Scalar;
  145. typedef typename Mat::Index Index;
  146. typedef SparseMatrix<Scalar,ColMajor> SpMat;
  147. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  148. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  149. // generate the problem
  150. Mat A, halfA;
  151. DenseMatrix dA;
  152. int size = generate_sparse_spd_problem(solver, A, halfA, dA);
  153. // generate the right hand sides
  154. int rhsCols = internal::random<int>(1,16);
  155. double density = (std::max)(8./(size*rhsCols), 0.1);
  156. SpMat B(size,rhsCols);
  157. DenseVector b = DenseVector::Random(size);
  158. DenseMatrix dB(size,rhsCols);
  159. initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
  160. for (int i = 0; i < g_repeat; i++) {
  161. check_sparse_solving(solver, A, b, dA, b);
  162. check_sparse_solving(solver, halfA, b, dA, b);
  163. check_sparse_solving(solver, A, dB, dA, dB);
  164. check_sparse_solving(solver, halfA, dB, dA, dB);
  165. check_sparse_solving(solver, A, B, dA, dB);
  166. check_sparse_solving(solver, halfA, B, dA, dB);
  167. }
  168. // First, get the folder
  169. #ifdef TEST_REAL_CASES
  170. if (internal::is_same<Scalar, float>::value
  171. || internal::is_same<Scalar, std::complex<float> >::value)
  172. return ;
  173. std::string mat_folder = get_matrixfolder<Scalar>();
  174. MatrixMarketIterator<Scalar> it(mat_folder);
  175. for (; it; ++it)
  176. {
  177. if (it.sym() == SPD){
  178. Mat halfA;
  179. PermutationMatrix<Dynamic, Dynamic, Index> pnull;
  180. halfA.template selfadjointView<Solver::UpLo>() = it.matrix().template triangularView<Eigen::Lower>().twistedBy(pnull);
  181. std::cout<< " ==== SOLVING WITH MATRIX " << it.matname() << " ==== \n";
  182. check_sparse_solving_real_cases(solver, it.matrix(), it.rhs(), it.refX());
  183. check_sparse_solving_real_cases(solver, halfA, it.rhs(), it.refX());
  184. }
  185. }
  186. #endif
  187. }
  188. template<typename Solver> void check_sparse_spd_determinant(Solver& solver)
  189. {
  190. typedef typename Solver::MatrixType Mat;
  191. typedef typename Mat::Scalar Scalar;
  192. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  193. // generate the problem
  194. Mat A, halfA;
  195. DenseMatrix dA;
  196. generate_sparse_spd_problem(solver, A, halfA, dA, 30);
  197. for (int i = 0; i < g_repeat; i++) {
  198. check_sparse_determinant(solver, A, dA);
  199. check_sparse_determinant(solver, halfA, dA );
  200. }
  201. }
  202. template<typename Solver, typename DenseMat>
  203. int generate_sparse_square_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300)
  204. {
  205. typedef typename Solver::MatrixType Mat;
  206. typedef typename Mat::Scalar Scalar;
  207. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  208. int size = internal::random<int>(1,maxSize);
  209. double density = (std::max)(8./(size*size), 0.01);
  210. A.resize(size,size);
  211. dA.resize(size,size);
  212. initSparse<Scalar>(density, dA, A, ForceNonZeroDiag);
  213. return size;
  214. }
  215. template<typename Solver> void check_sparse_square_solving(Solver& solver)
  216. {
  217. typedef typename Solver::MatrixType Mat;
  218. typedef typename Mat::Scalar Scalar;
  219. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  220. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  221. int rhsCols = internal::random<int>(1,16);
  222. Mat A;
  223. DenseMatrix dA;
  224. int size = generate_sparse_square_problem(solver, A, dA);
  225. DenseVector b = DenseVector::Random(size);
  226. DenseMatrix dB = DenseMatrix::Random(size,rhsCols);
  227. A.makeCompressed();
  228. for (int i = 0; i < g_repeat; i++) {
  229. check_sparse_solving(solver, A, b, dA, b);
  230. check_sparse_solving(solver, A, dB, dA, dB);
  231. }
  232. // First, get the folder
  233. #ifdef TEST_REAL_CASES
  234. if (internal::is_same<Scalar, float>::value
  235. || internal::is_same<Scalar, std::complex<float> >::value)
  236. return ;
  237. std::string mat_folder = get_matrixfolder<Scalar>();
  238. MatrixMarketIterator<Scalar> it(mat_folder);
  239. for (; it; ++it)
  240. {
  241. std::cout<< " ==== SOLVING WITH MATRIX " << it.matname() << " ==== \n";
  242. check_sparse_solving_real_cases(solver, it.matrix(), it.rhs(), it.refX());
  243. }
  244. #endif
  245. }
  246. template<typename Solver> void check_sparse_square_determinant(Solver& solver)
  247. {
  248. typedef typename Solver::MatrixType Mat;
  249. typedef typename Mat::Scalar Scalar;
  250. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  251. // generate the problem
  252. Mat A;
  253. DenseMatrix dA;
  254. generate_sparse_square_problem(solver, A, dA, 30);
  255. A.makeCompressed();
  256. for (int i = 0; i < g_repeat; i++) {
  257. check_sparse_determinant(solver, A, dA);
  258. }
  259. }