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.

523 lines
17 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. #include <sstream>
  12. template<typename Solver, typename Rhs, typename DenseMat, typename DenseRhs>
  13. void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const DenseMat& dA, const DenseRhs& db)
  14. {
  15. typedef typename Solver::MatrixType Mat;
  16. typedef typename Mat::Scalar Scalar;
  17. typedef typename Mat::StorageIndex StorageIndex;
  18. DenseRhs refX = dA.householderQr().solve(db);
  19. {
  20. Rhs x(A.cols(), b.cols());
  21. Rhs oldb = b;
  22. solver.compute(A);
  23. if (solver.info() != Success)
  24. {
  25. std::cerr << "ERROR | sparse solver testing, factorization failed (" << typeid(Solver).name() << ")\n";
  26. VERIFY(solver.info() == Success);
  27. }
  28. x = solver.solve(b);
  29. if (solver.info() != Success)
  30. {
  31. std::cerr << "WARNING | sparse solver testing: solving failed (" << typeid(Solver).name() << ")\n";
  32. return;
  33. }
  34. VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
  35. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  36. x.setZero();
  37. // test the analyze/factorize API
  38. solver.analyzePattern(A);
  39. solver.factorize(A);
  40. VERIFY(solver.info() == Success && "factorization failed when using analyzePattern/factorize API");
  41. x = solver.solve(b);
  42. VERIFY(solver.info() == Success && "solving failed when using analyzePattern/factorize API");
  43. VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
  44. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  45. x.setZero();
  46. // test with Map
  47. MappedSparseMatrix<Scalar,Mat::Options,StorageIndex> Am(A.rows(), A.cols(), A.nonZeros(), const_cast<StorageIndex*>(A.outerIndexPtr()), const_cast<StorageIndex*>(A.innerIndexPtr()), const_cast<Scalar*>(A.valuePtr()));
  48. solver.compute(Am);
  49. VERIFY(solver.info() == Success && "factorization failed when using Map");
  50. DenseRhs dx(refX);
  51. dx.setZero();
  52. Map<DenseRhs> xm(dx.data(), dx.rows(), dx.cols());
  53. Map<const DenseRhs> bm(db.data(), db.rows(), db.cols());
  54. xm = solver.solve(bm);
  55. VERIFY(solver.info() == Success && "solving failed when using Map");
  56. VERIFY(oldb.isApprox(bm) && "sparse solver testing: the rhs should not be modified!");
  57. VERIFY(xm.isApprox(refX,test_precision<Scalar>()));
  58. }
  59. // if not too large, do some extra check:
  60. if(A.rows()<2000)
  61. {
  62. // test initialization ctor
  63. {
  64. Rhs x(b.rows(), b.cols());
  65. Solver solver2(A);
  66. VERIFY(solver2.info() == Success);
  67. x = solver2.solve(b);
  68. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  69. }
  70. // test dense Block as the result and rhs:
  71. {
  72. DenseRhs x(refX.rows(), refX.cols());
  73. DenseRhs oldb(db);
  74. x.setZero();
  75. x.block(0,0,x.rows(),x.cols()) = solver.solve(db.block(0,0,db.rows(),db.cols()));
  76. VERIFY(oldb.isApprox(db) && "sparse solver testing: the rhs should not be modified!");
  77. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  78. }
  79. // test uncompressed inputs
  80. {
  81. Mat A2 = A;
  82. A2.reserve((ArrayXf::Random(A.outerSize())+2).template cast<typename Mat::StorageIndex>().eval());
  83. solver.compute(A2);
  84. Rhs x = solver.solve(b);
  85. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  86. }
  87. // test expression as input
  88. {
  89. solver.compute(0.5*(A+A));
  90. Rhs x = solver.solve(b);
  91. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  92. Solver solver2(0.5*(A+A));
  93. Rhs x2 = solver2.solve(b);
  94. VERIFY(x2.isApprox(refX,test_precision<Scalar>()));
  95. }
  96. }
  97. }
  98. template<typename Solver, typename Rhs>
  99. void check_sparse_solving_real_cases(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const typename Solver::MatrixType& fullA, const Rhs& refX)
  100. {
  101. typedef typename Solver::MatrixType Mat;
  102. typedef typename Mat::Scalar Scalar;
  103. typedef typename Mat::RealScalar RealScalar;
  104. Rhs x(A.cols(), b.cols());
  105. solver.compute(A);
  106. if (solver.info() != Success)
  107. {
  108. std::cerr << "ERROR | sparse solver testing, factorization failed (" << typeid(Solver).name() << ")\n";
  109. VERIFY(solver.info() == Success);
  110. }
  111. x = solver.solve(b);
  112. if (solver.info() != Success)
  113. {
  114. std::cerr << "WARNING | sparse solver testing, solving failed (" << typeid(Solver).name() << ")\n";
  115. return;
  116. }
  117. RealScalar res_error = (fullA*x-b).norm()/b.norm();
  118. VERIFY( (res_error <= test_precision<Scalar>() ) && "sparse solver failed without noticing it");
  119. if(refX.size() != 0 && (refX - x).norm()/refX.norm() > test_precision<Scalar>())
  120. {
  121. std::cerr << "WARNING | found solution is different from the provided reference one\n";
  122. }
  123. }
  124. template<typename Solver, typename DenseMat>
  125. void check_sparse_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
  126. {
  127. typedef typename Solver::MatrixType Mat;
  128. typedef typename Mat::Scalar Scalar;
  129. solver.compute(A);
  130. if (solver.info() != Success)
  131. {
  132. std::cerr << "WARNING | sparse solver testing: factorization failed (check_sparse_determinant)\n";
  133. return;
  134. }
  135. Scalar refDet = dA.determinant();
  136. VERIFY_IS_APPROX(refDet,solver.determinant());
  137. }
  138. template<typename Solver, typename DenseMat>
  139. void check_sparse_abs_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
  140. {
  141. using std::abs;
  142. typedef typename Solver::MatrixType Mat;
  143. typedef typename Mat::Scalar Scalar;
  144. solver.compute(A);
  145. if (solver.info() != Success)
  146. {
  147. std::cerr << "WARNING | sparse solver testing: factorization failed (check_sparse_abs_determinant)\n";
  148. return;
  149. }
  150. Scalar refDet = abs(dA.determinant());
  151. VERIFY_IS_APPROX(refDet,solver.absDeterminant());
  152. }
  153. template<typename Solver, typename DenseMat>
  154. int generate_sparse_spd_problem(Solver& , typename Solver::MatrixType& A, typename Solver::MatrixType& halfA, DenseMat& dA, int maxSize = 300)
  155. {
  156. typedef typename Solver::MatrixType Mat;
  157. typedef typename Mat::Scalar Scalar;
  158. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  159. int size = internal::random<int>(1,maxSize);
  160. double density = (std::max)(8./(size*size), 0.01);
  161. Mat M(size, size);
  162. DenseMatrix dM(size, size);
  163. initSparse<Scalar>(density, dM, M, ForceNonZeroDiag);
  164. A = M * M.adjoint();
  165. dA = dM * dM.adjoint();
  166. halfA.resize(size,size);
  167. if(Solver::UpLo==(Lower|Upper))
  168. halfA = A;
  169. else
  170. halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
  171. return size;
  172. }
  173. #ifdef TEST_REAL_CASES
  174. template<typename Scalar>
  175. inline std::string get_matrixfolder()
  176. {
  177. std::string mat_folder = TEST_REAL_CASES;
  178. if( internal::is_same<Scalar, std::complex<float> >::value || internal::is_same<Scalar, std::complex<double> >::value )
  179. mat_folder = mat_folder + static_cast<std::string>("/complex/");
  180. else
  181. mat_folder = mat_folder + static_cast<std::string>("/real/");
  182. return mat_folder;
  183. }
  184. std::string sym_to_string(int sym)
  185. {
  186. if(sym==Symmetric) return "Symmetric ";
  187. if(sym==SPD) return "SPD ";
  188. return "";
  189. }
  190. template<typename Derived>
  191. std::string solver_stats(const IterativeSolverBase<Derived> &solver)
  192. {
  193. std::stringstream ss;
  194. ss << solver.iterations() << " iters, error: " << solver.error();
  195. return ss.str();
  196. }
  197. template<typename Derived>
  198. std::string solver_stats(const SparseSolverBase<Derived> &/*solver*/)
  199. {
  200. return "";
  201. }
  202. #endif
  203. template<typename Solver> void check_sparse_spd_solving(Solver& solver, int maxSize = 300, int maxRealWorldSize = 100000)
  204. {
  205. typedef typename Solver::MatrixType Mat;
  206. typedef typename Mat::Scalar Scalar;
  207. typedef typename Mat::StorageIndex StorageIndex;
  208. typedef SparseMatrix<Scalar,ColMajor, StorageIndex> SpMat;
  209. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  210. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  211. // generate the problem
  212. Mat A, halfA;
  213. DenseMatrix dA;
  214. for (int i = 0; i < g_repeat; i++) {
  215. int size = generate_sparse_spd_problem(solver, A, halfA, dA, maxSize);
  216. // generate the right hand sides
  217. int rhsCols = internal::random<int>(1,16);
  218. double density = (std::max)(8./(size*rhsCols), 0.1);
  219. SpMat B(size,rhsCols);
  220. DenseVector b = DenseVector::Random(size);
  221. DenseMatrix dB(size,rhsCols);
  222. initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
  223. CALL_SUBTEST( check_sparse_solving(solver, A, b, dA, b) );
  224. CALL_SUBTEST( check_sparse_solving(solver, halfA, b, dA, b) );
  225. CALL_SUBTEST( check_sparse_solving(solver, A, dB, dA, dB) );
  226. CALL_SUBTEST( check_sparse_solving(solver, halfA, dB, dA, dB) );
  227. CALL_SUBTEST( check_sparse_solving(solver, A, B, dA, dB) );
  228. CALL_SUBTEST( check_sparse_solving(solver, halfA, B, dA, dB) );
  229. // check only once
  230. if(i==0)
  231. {
  232. b = DenseVector::Zero(size);
  233. check_sparse_solving(solver, A, b, dA, b);
  234. }
  235. }
  236. // First, get the folder
  237. #ifdef TEST_REAL_CASES
  238. // Test real problems with double precision only
  239. if (internal::is_same<typename NumTraits<Scalar>::Real, double>::value)
  240. {
  241. std::string mat_folder = get_matrixfolder<Scalar>();
  242. MatrixMarketIterator<Scalar> it(mat_folder);
  243. for (; it; ++it)
  244. {
  245. if (it.sym() == SPD){
  246. A = it.matrix();
  247. if(A.diagonal().size() <= maxRealWorldSize)
  248. {
  249. DenseVector b = it.rhs();
  250. DenseVector refX = it.refX();
  251. PermutationMatrix<Dynamic, Dynamic, StorageIndex> pnull;
  252. halfA.resize(A.rows(), A.cols());
  253. if(Solver::UpLo == (Lower|Upper))
  254. halfA = A;
  255. else
  256. halfA.template selfadjointView<Solver::UpLo>() = A.template triangularView<StormEigen::Lower>().twistedBy(pnull);
  257. std::cout << "INFO | Testing " << sym_to_string(it.sym()) << "sparse problem " << it.matname()
  258. << " (" << A.rows() << "x" << A.cols() << ") using " << typeid(Solver).name() << "..." << std::endl;
  259. CALL_SUBTEST( check_sparse_solving_real_cases(solver, A, b, A, refX) );
  260. std::string stats = solver_stats(solver);
  261. if(stats.size()>0)
  262. std::cout << "INFO | " << stats << std::endl;
  263. CALL_SUBTEST( check_sparse_solving_real_cases(solver, halfA, b, A, refX) );
  264. }
  265. else
  266. {
  267. std::cout << "INFO | Skip sparse problem \"" << it.matname() << "\" (too large)" << std::endl;
  268. }
  269. }
  270. }
  271. }
  272. #else
  273. EIGEN_UNUSED_VARIABLE(maxRealWorldSize);
  274. #endif
  275. }
  276. template<typename Solver> void check_sparse_spd_determinant(Solver& solver)
  277. {
  278. typedef typename Solver::MatrixType Mat;
  279. typedef typename Mat::Scalar Scalar;
  280. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  281. // generate the problem
  282. Mat A, halfA;
  283. DenseMatrix dA;
  284. generate_sparse_spd_problem(solver, A, halfA, dA, 30);
  285. for (int i = 0; i < g_repeat; i++) {
  286. check_sparse_determinant(solver, A, dA);
  287. check_sparse_determinant(solver, halfA, dA );
  288. }
  289. }
  290. template<typename Solver, typename DenseMat>
  291. Index generate_sparse_square_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300, int options = ForceNonZeroDiag)
  292. {
  293. typedef typename Solver::MatrixType Mat;
  294. typedef typename Mat::Scalar Scalar;
  295. Index size = internal::random<int>(1,maxSize);
  296. double density = (std::max)(8./(size*size), 0.01);
  297. A.resize(size,size);
  298. dA.resize(size,size);
  299. initSparse<Scalar>(density, dA, A, options);
  300. return size;
  301. }
  302. struct prune_column {
  303. Index m_col;
  304. prune_column(Index col) : m_col(col) {}
  305. template<class Scalar>
  306. bool operator()(Index, Index col, const Scalar&) const {
  307. return col != m_col;
  308. }
  309. };
  310. template<typename Solver> void check_sparse_square_solving(Solver& solver, int maxSize = 300, int maxRealWorldSize = 100000, bool checkDeficient = false)
  311. {
  312. typedef typename Solver::MatrixType Mat;
  313. typedef typename Mat::Scalar Scalar;
  314. typedef SparseMatrix<Scalar,ColMajor, typename Mat::StorageIndex> SpMat;
  315. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  316. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  317. int rhsCols = internal::random<int>(1,16);
  318. Mat A;
  319. DenseMatrix dA;
  320. for (int i = 0; i < g_repeat; i++) {
  321. Index size = generate_sparse_square_problem(solver, A, dA, maxSize);
  322. A.makeCompressed();
  323. DenseVector b = DenseVector::Random(size);
  324. DenseMatrix dB(size,rhsCols);
  325. SpMat B(size,rhsCols);
  326. double density = (std::max)(8./(size*rhsCols), 0.1);
  327. initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
  328. B.makeCompressed();
  329. CALL_SUBTEST(check_sparse_solving(solver, A, b, dA, b));
  330. CALL_SUBTEST(check_sparse_solving(solver, A, dB, dA, dB));
  331. CALL_SUBTEST(check_sparse_solving(solver, A, B, dA, dB));
  332. // check only once
  333. if(i==0)
  334. {
  335. b = DenseVector::Zero(size);
  336. check_sparse_solving(solver, A, b, dA, b);
  337. }
  338. // regression test for Bug 792 (structurally rank deficient matrices):
  339. if(checkDeficient && size>1) {
  340. Index col = internal::random<int>(0,int(size-1));
  341. A.prune(prune_column(col));
  342. solver.compute(A);
  343. VERIFY_IS_EQUAL(solver.info(), NumericalIssue);
  344. }
  345. }
  346. // First, get the folder
  347. #ifdef TEST_REAL_CASES
  348. // Test real problems with double precision only
  349. if (internal::is_same<typename NumTraits<Scalar>::Real, double>::value)
  350. {
  351. std::string mat_folder = get_matrixfolder<Scalar>();
  352. MatrixMarketIterator<Scalar> it(mat_folder);
  353. for (; it; ++it)
  354. {
  355. A = it.matrix();
  356. if(A.diagonal().size() <= maxRealWorldSize)
  357. {
  358. DenseVector b = it.rhs();
  359. DenseVector refX = it.refX();
  360. std::cout << "INFO | Testing " << sym_to_string(it.sym()) << "sparse problem " << it.matname()
  361. << " (" << A.rows() << "x" << A.cols() << ") using " << typeid(Solver).name() << "..." << std::endl;
  362. CALL_SUBTEST(check_sparse_solving_real_cases(solver, A, b, A, refX));
  363. std::string stats = solver_stats(solver);
  364. if(stats.size()>0)
  365. std::cout << "INFO | " << stats << std::endl;
  366. }
  367. else
  368. {
  369. std::cout << "INFO | SKIP sparse problem \"" << it.matname() << "\" (too large)" << std::endl;
  370. }
  371. }
  372. }
  373. #else
  374. EIGEN_UNUSED_VARIABLE(maxRealWorldSize);
  375. #endif
  376. }
  377. template<typename Solver> void check_sparse_square_determinant(Solver& solver)
  378. {
  379. typedef typename Solver::MatrixType Mat;
  380. typedef typename Mat::Scalar Scalar;
  381. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  382. for (int i = 0; i < g_repeat; i++) {
  383. // generate the problem
  384. Mat A;
  385. DenseMatrix dA;
  386. int size = internal::random<int>(1,30);
  387. dA.setRandom(size,size);
  388. dA = (dA.array().abs()<0.3).select(0,dA);
  389. dA.diagonal() = (dA.diagonal().array()==0).select(1,dA.diagonal());
  390. A = dA.sparseView();
  391. A.makeCompressed();
  392. check_sparse_determinant(solver, A, dA);
  393. }
  394. }
  395. template<typename Solver> void check_sparse_square_abs_determinant(Solver& solver)
  396. {
  397. typedef typename Solver::MatrixType Mat;
  398. typedef typename Mat::Scalar Scalar;
  399. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  400. for (int i = 0; i < g_repeat; i++) {
  401. // generate the problem
  402. Mat A;
  403. DenseMatrix dA;
  404. generate_sparse_square_problem(solver, A, dA, 30);
  405. A.makeCompressed();
  406. check_sparse_abs_determinant(solver, A, dA);
  407. }
  408. }
  409. template<typename Solver, typename DenseMat>
  410. void generate_sparse_leastsquare_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300, int options = ForceNonZeroDiag)
  411. {
  412. typedef typename Solver::MatrixType Mat;
  413. typedef typename Mat::Scalar Scalar;
  414. int rows = internal::random<int>(1,maxSize);
  415. int cols = internal::random<int>(1,rows);
  416. double density = (std::max)(8./(rows*cols), 0.01);
  417. A.resize(rows,cols);
  418. dA.resize(rows,cols);
  419. initSparse<Scalar>(density, dA, A, options);
  420. }
  421. template<typename Solver> void check_sparse_leastsquare_solving(Solver& solver)
  422. {
  423. typedef typename Solver::MatrixType Mat;
  424. typedef typename Mat::Scalar Scalar;
  425. typedef SparseMatrix<Scalar,ColMajor, typename Mat::StorageIndex> SpMat;
  426. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  427. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  428. int rhsCols = internal::random<int>(1,16);
  429. Mat A;
  430. DenseMatrix dA;
  431. for (int i = 0; i < g_repeat; i++) {
  432. generate_sparse_leastsquare_problem(solver, A, dA);
  433. A.makeCompressed();
  434. DenseVector b = DenseVector::Random(A.rows());
  435. DenseMatrix dB(A.rows(),rhsCols);
  436. SpMat B(A.rows(),rhsCols);
  437. double density = (std::max)(8./(A.rows()*rhsCols), 0.1);
  438. initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
  439. B.makeCompressed();
  440. check_sparse_solving(solver, A, b, dA, b);
  441. check_sparse_solving(solver, A, dB, dA, dB);
  442. check_sparse_solving(solver, A, B, dA, dB);
  443. // check only once
  444. if(i==0)
  445. {
  446. b = DenseVector::Zero(A.rows());
  447. check_sparse_solving(solver, A, b, dA, b);
  448. }
  449. }
  450. }