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.

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