sparse_solver.h 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699
  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 <Eigen/SparseLU>
  12. #include <sstream>
  13. template<typename Solver, typename Rhs, typename Guess,typename Result>
  14. void solve_with_guess(IterativeSolverBase<Solver>& solver, const MatrixBase<Rhs>& b, const Guess& g, Result &x) {
  15. if(internal::random<bool>())
  16. {
  17. // With a temporary through evaluator<SolveWithGuess>
  18. x = solver.derived().solveWithGuess(b,g) + Result::Zero(x.rows(), x.cols());
  19. }
  20. else
  21. {
  22. // direct evaluation within x through Assignment<Result,SolveWithGuess>
  23. x = solver.derived().solveWithGuess(b.derived(),g);
  24. }
  25. }
  26. template<typename Solver, typename Rhs, typename Guess,typename Result>
  27. void solve_with_guess(SparseSolverBase<Solver>& solver, const MatrixBase<Rhs>& b, const Guess& , Result& x) {
  28. if(internal::random<bool>())
  29. x = solver.derived().solve(b) + Result::Zero(x.rows(), x.cols());
  30. else
  31. x = solver.derived().solve(b);
  32. }
  33. template<typename Solver, typename Rhs, typename Guess,typename Result>
  34. void solve_with_guess(SparseSolverBase<Solver>& solver, const SparseMatrixBase<Rhs>& b, const Guess& , Result& x) {
  35. x = solver.derived().solve(b);
  36. }
  37. template<typename Solver, typename Rhs, typename DenseMat, typename DenseRhs>
  38. void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const DenseMat& dA, const DenseRhs& db)
  39. {
  40. typedef typename Solver::MatrixType Mat;
  41. typedef typename Mat::Scalar Scalar;
  42. typedef typename Mat::StorageIndex StorageIndex;
  43. DenseRhs refX = dA.householderQr().solve(db);
  44. {
  45. Rhs x(A.cols(), b.cols());
  46. Rhs oldb = b;
  47. solver.compute(A);
  48. if (solver.info() != Success)
  49. {
  50. std::cerr << "ERROR | sparse solver testing, factorization failed (" << typeid(Solver).name() << ")\n";
  51. VERIFY(solver.info() == Success);
  52. }
  53. x = solver.solve(b);
  54. if (solver.info() != Success)
  55. {
  56. std::cerr << "WARNING: sparse solver testing: solving failed (" << typeid(Solver).name() << ")\n";
  57. // dump call stack:
  58. g_test_level++;
  59. VERIFY(solver.info() == Success);
  60. g_test_level--;
  61. return;
  62. }
  63. VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
  64. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  65. x.setZero();
  66. solve_with_guess(solver, b, x, x);
  67. VERIFY(solver.info() == Success && "solving failed when using solve_with_guess API");
  68. VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
  69. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  70. x.setZero();
  71. // test the analyze/factorize API
  72. solver.analyzePattern(A);
  73. solver.factorize(A);
  74. VERIFY(solver.info() == Success && "factorization failed when using analyzePattern/factorize API");
  75. x = solver.solve(b);
  76. VERIFY(solver.info() == Success && "solving failed when using analyzePattern/factorize API");
  77. VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
  78. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  79. x.setZero();
  80. // test with Map
  81. 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()));
  82. solver.compute(Am);
  83. VERIFY(solver.info() == Success && "factorization failed when using Map");
  84. DenseRhs dx(refX);
  85. dx.setZero();
  86. Map<DenseRhs> xm(dx.data(), dx.rows(), dx.cols());
  87. Map<const DenseRhs> bm(db.data(), db.rows(), db.cols());
  88. xm = solver.solve(bm);
  89. VERIFY(solver.info() == Success && "solving failed when using Map");
  90. VERIFY(oldb.isApprox(bm) && "sparse solver testing: the rhs should not be modified!");
  91. VERIFY(xm.isApprox(refX,test_precision<Scalar>()));
  92. }
  93. // if not too large, do some extra check:
  94. if(A.rows()<2000)
  95. {
  96. // test initialization ctor
  97. {
  98. Rhs x(b.rows(), b.cols());
  99. Solver solver2(A);
  100. VERIFY(solver2.info() == Success);
  101. x = solver2.solve(b);
  102. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  103. }
  104. // test dense Block as the result and rhs:
  105. {
  106. DenseRhs x(refX.rows(), refX.cols());
  107. DenseRhs oldb(db);
  108. x.setZero();
  109. x.block(0,0,x.rows(),x.cols()) = solver.solve(db.block(0,0,db.rows(),db.cols()));
  110. VERIFY(oldb.isApprox(db) && "sparse solver testing: the rhs should not be modified!");
  111. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  112. }
  113. // test uncompressed inputs
  114. {
  115. Mat A2 = A;
  116. A2.reserve((ArrayXf::Random(A.outerSize())+2).template cast<typename Mat::StorageIndex>().eval());
  117. solver.compute(A2);
  118. Rhs x = solver.solve(b);
  119. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  120. }
  121. // test expression as input
  122. {
  123. solver.compute(0.5*(A+A));
  124. Rhs x = solver.solve(b);
  125. VERIFY(x.isApprox(refX,test_precision<Scalar>()));
  126. Solver solver2(0.5*(A+A));
  127. Rhs x2 = solver2.solve(b);
  128. VERIFY(x2.isApprox(refX,test_precision<Scalar>()));
  129. }
  130. }
  131. }
  132. // specialization of generic check_sparse_solving for SuperLU in order to also test adjoint and transpose solves
  133. template<typename Scalar, typename Rhs, typename DenseMat, typename DenseRhs>
  134. void check_sparse_solving(Eigen::SparseLU<Eigen::SparseMatrix<Scalar> >& solver, const typename Eigen::SparseMatrix<Scalar>& A, const Rhs& b, const DenseMat& dA, const DenseRhs& db)
  135. {
  136. typedef typename Eigen::SparseMatrix<Scalar> Mat;
  137. typedef typename Mat::StorageIndex StorageIndex;
  138. typedef typename Eigen::SparseLU<Eigen::SparseMatrix<Scalar> > Solver;
  139. // reference solutions computed by dense QR solver
  140. DenseRhs refX1 = dA.householderQr().solve(db); // solution of A x = db
  141. DenseRhs refX2 = dA.transpose().householderQr().solve(db); // solution of A^T * x = db (use transposed matrix A^T)
  142. DenseRhs refX3 = dA.adjoint().householderQr().solve(db); // solution of A^* * x = db (use adjoint matrix A^*)
  143. {
  144. Rhs x1(A.cols(), b.cols());
  145. Rhs x2(A.cols(), b.cols());
  146. Rhs x3(A.cols(), b.cols());
  147. Rhs oldb = b;
  148. solver.compute(A);
  149. if (solver.info() != Success)
  150. {
  151. std::cerr << "ERROR | sparse solver testing, factorization failed (" << typeid(Solver).name() << ")\n";
  152. VERIFY(solver.info() == Success);
  153. }
  154. x1 = solver.solve(b);
  155. if (solver.info() != Success)
  156. {
  157. std::cerr << "WARNING | sparse solver testing: solving failed (" << typeid(Solver).name() << ")\n";
  158. return;
  159. }
  160. VERIFY(oldb.isApprox(b,0.0) && "sparse solver testing: the rhs should not be modified!");
  161. VERIFY(x1.isApprox(refX1,test_precision<Scalar>()));
  162. // test solve with transposed
  163. x2 = solver.transpose().solve(b);
  164. VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
  165. VERIFY(x2.isApprox(refX2,test_precision<Scalar>()));
  166. // test solve with adjoint
  167. //solver.template _solve_impl_transposed<true>(b, x3);
  168. x3 = solver.adjoint().solve(b);
  169. VERIFY(oldb.isApprox(b,0.0) && "sparse solver testing: the rhs should not be modified!");
  170. VERIFY(x3.isApprox(refX3,test_precision<Scalar>()));
  171. x1.setZero();
  172. solve_with_guess(solver, b, x1, x1);
  173. VERIFY(solver.info() == Success && "solving failed when using analyzePattern/factorize API");
  174. VERIFY(oldb.isApprox(b,0.0) && "sparse solver testing: the rhs should not be modified!");
  175. VERIFY(x1.isApprox(refX1,test_precision<Scalar>()));
  176. x1.setZero();
  177. x2.setZero();
  178. x3.setZero();
  179. // test the analyze/factorize API
  180. solver.analyzePattern(A);
  181. solver.factorize(A);
  182. VERIFY(solver.info() == Success && "factorization failed when using analyzePattern/factorize API");
  183. x1 = solver.solve(b);
  184. x2 = solver.transpose().solve(b);
  185. x3 = solver.adjoint().solve(b);
  186. VERIFY(solver.info() == Success && "solving failed when using analyzePattern/factorize API");
  187. VERIFY(oldb.isApprox(b,0.0) && "sparse solver testing: the rhs should not be modified!");
  188. VERIFY(x1.isApprox(refX1,test_precision<Scalar>()));
  189. VERIFY(x2.isApprox(refX2,test_precision<Scalar>()));
  190. VERIFY(x3.isApprox(refX3,test_precision<Scalar>()));
  191. x1.setZero();
  192. // test with Map
  193. 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()));
  194. solver.compute(Am);
  195. VERIFY(solver.info() == Success && "factorization failed when using Map");
  196. DenseRhs dx(refX1);
  197. dx.setZero();
  198. Map<DenseRhs> xm(dx.data(), dx.rows(), dx.cols());
  199. Map<const DenseRhs> bm(db.data(), db.rows(), db.cols());
  200. xm = solver.solve(bm);
  201. VERIFY(solver.info() == Success && "solving failed when using Map");
  202. VERIFY(oldb.isApprox(bm,0.0) && "sparse solver testing: the rhs should not be modified!");
  203. VERIFY(xm.isApprox(refX1,test_precision<Scalar>()));
  204. }
  205. // if not too large, do some extra check:
  206. if(A.rows()<2000)
  207. {
  208. // test initialization ctor
  209. {
  210. Rhs x(b.rows(), b.cols());
  211. Solver solver2(A);
  212. VERIFY(solver2.info() == Success);
  213. x = solver2.solve(b);
  214. VERIFY(x.isApprox(refX1,test_precision<Scalar>()));
  215. }
  216. // test dense Block as the result and rhs:
  217. {
  218. DenseRhs x(refX1.rows(), refX1.cols());
  219. DenseRhs oldb(db);
  220. x.setZero();
  221. x.block(0,0,x.rows(),x.cols()) = solver.solve(db.block(0,0,db.rows(),db.cols()));
  222. VERIFY(oldb.isApprox(db,0.0) && "sparse solver testing: the rhs should not be modified!");
  223. VERIFY(x.isApprox(refX1,test_precision<Scalar>()));
  224. }
  225. // test uncompressed inputs
  226. {
  227. Mat A2 = A;
  228. A2.reserve((ArrayXf::Random(A.outerSize())+2).template cast<typename Mat::StorageIndex>().eval());
  229. solver.compute(A2);
  230. Rhs x = solver.solve(b);
  231. VERIFY(x.isApprox(refX1,test_precision<Scalar>()));
  232. }
  233. // test expression as input
  234. {
  235. solver.compute(0.5*(A+A));
  236. Rhs x = solver.solve(b);
  237. VERIFY(x.isApprox(refX1,test_precision<Scalar>()));
  238. Solver solver2(0.5*(A+A));
  239. Rhs x2 = solver2.solve(b);
  240. VERIFY(x2.isApprox(refX1,test_precision<Scalar>()));
  241. }
  242. }
  243. }
  244. template<typename Solver, typename Rhs>
  245. void check_sparse_solving_real_cases(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const typename Solver::MatrixType& fullA, const Rhs& refX)
  246. {
  247. typedef typename Solver::MatrixType Mat;
  248. typedef typename Mat::Scalar Scalar;
  249. typedef typename Mat::RealScalar RealScalar;
  250. Rhs x(A.cols(), b.cols());
  251. solver.compute(A);
  252. if (solver.info() != Success)
  253. {
  254. std::cerr << "ERROR | sparse solver testing, factorization failed (" << typeid(Solver).name() << ")\n";
  255. VERIFY(solver.info() == Success);
  256. }
  257. x = solver.solve(b);
  258. if (solver.info() != Success)
  259. {
  260. std::cerr << "WARNING | sparse solver testing, solving failed (" << typeid(Solver).name() << ")\n";
  261. return;
  262. }
  263. RealScalar res_error = (fullA*x-b).norm()/b.norm();
  264. VERIFY( (res_error <= test_precision<Scalar>() ) && "sparse solver failed without noticing it");
  265. if(refX.size() != 0 && (refX - x).norm()/refX.norm() > test_precision<Scalar>())
  266. {
  267. std::cerr << "WARNING | found solution is different from the provided reference one\n";
  268. }
  269. }
  270. template<typename Solver, typename DenseMat>
  271. void check_sparse_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
  272. {
  273. typedef typename Solver::MatrixType Mat;
  274. typedef typename Mat::Scalar Scalar;
  275. solver.compute(A);
  276. if (solver.info() != Success)
  277. {
  278. std::cerr << "WARNING | sparse solver testing: factorization failed (check_sparse_determinant)\n";
  279. return;
  280. }
  281. Scalar refDet = dA.determinant();
  282. VERIFY_IS_APPROX(refDet,solver.determinant());
  283. }
  284. template<typename Solver, typename DenseMat>
  285. void check_sparse_abs_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
  286. {
  287. using std::abs;
  288. typedef typename Solver::MatrixType Mat;
  289. typedef typename Mat::Scalar Scalar;
  290. solver.compute(A);
  291. if (solver.info() != Success)
  292. {
  293. std::cerr << "WARNING | sparse solver testing: factorization failed (check_sparse_abs_determinant)\n";
  294. return;
  295. }
  296. Scalar refDet = abs(dA.determinant());
  297. VERIFY_IS_APPROX(refDet,solver.absDeterminant());
  298. }
  299. template<typename Solver, typename DenseMat>
  300. int generate_sparse_spd_problem(Solver& , typename Solver::MatrixType& A, typename Solver::MatrixType& halfA, DenseMat& dA, int maxSize = 300)
  301. {
  302. typedef typename Solver::MatrixType Mat;
  303. typedef typename Mat::Scalar Scalar;
  304. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  305. int size = internal::random<int>(1,maxSize);
  306. double density = (std::max)(8./(size*size), 0.01);
  307. Mat M(size, size);
  308. DenseMatrix dM(size, size);
  309. initSparse<Scalar>(density, dM, M, ForceNonZeroDiag);
  310. A = M * M.adjoint();
  311. dA = dM * dM.adjoint();
  312. halfA.resize(size,size);
  313. if(Solver::UpLo==(Lower|Upper))
  314. halfA = A;
  315. else
  316. halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
  317. return size;
  318. }
  319. #ifdef TEST_REAL_CASES
  320. template<typename Scalar>
  321. inline std::string get_matrixfolder()
  322. {
  323. std::string mat_folder = TEST_REAL_CASES;
  324. if( internal::is_same<Scalar, std::complex<float> >::value || internal::is_same<Scalar, std::complex<double> >::value )
  325. mat_folder = mat_folder + static_cast<std::string>("/complex/");
  326. else
  327. mat_folder = mat_folder + static_cast<std::string>("/real/");
  328. return mat_folder;
  329. }
  330. std::string sym_to_string(int sym)
  331. {
  332. if(sym==Symmetric) return "Symmetric ";
  333. if(sym==SPD) return "SPD ";
  334. return "";
  335. }
  336. template<typename Derived>
  337. std::string solver_stats(const IterativeSolverBase<Derived> &solver)
  338. {
  339. std::stringstream ss;
  340. ss << solver.iterations() << " iters, error: " << solver.error();
  341. return ss.str();
  342. }
  343. template<typename Derived>
  344. std::string solver_stats(const SparseSolverBase<Derived> &/*solver*/)
  345. {
  346. return "";
  347. }
  348. #endif
  349. template<typename Solver> void check_sparse_spd_solving(Solver& solver, int maxSize = (std::min)(300,EIGEN_TEST_MAX_SIZE), int maxRealWorldSize = 100000)
  350. {
  351. typedef typename Solver::MatrixType Mat;
  352. typedef typename Mat::Scalar Scalar;
  353. typedef typename Mat::StorageIndex StorageIndex;
  354. typedef SparseMatrix<Scalar,ColMajor, StorageIndex> SpMat;
  355. typedef SparseVector<Scalar, 0, StorageIndex> SpVec;
  356. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  357. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  358. // generate the problem
  359. Mat A, halfA;
  360. DenseMatrix dA;
  361. for (int i = 0; i < g_repeat; i++) {
  362. int size = generate_sparse_spd_problem(solver, A, halfA, dA, maxSize);
  363. // generate the right hand sides
  364. int rhsCols = internal::random<int>(1,16);
  365. double density = (std::max)(8./(size*rhsCols), 0.1);
  366. SpMat B(size,rhsCols);
  367. DenseVector b = DenseVector::Random(size);
  368. DenseMatrix dB(size,rhsCols);
  369. initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
  370. SpVec c = B.col(0);
  371. DenseVector dc = dB.col(0);
  372. CALL_SUBTEST( check_sparse_solving(solver, A, b, dA, b) );
  373. CALL_SUBTEST( check_sparse_solving(solver, halfA, b, dA, b) );
  374. CALL_SUBTEST( check_sparse_solving(solver, A, dB, dA, dB) );
  375. CALL_SUBTEST( check_sparse_solving(solver, halfA, dB, dA, dB) );
  376. CALL_SUBTEST( check_sparse_solving(solver, A, B, dA, dB) );
  377. CALL_SUBTEST( check_sparse_solving(solver, halfA, B, dA, dB) );
  378. CALL_SUBTEST( check_sparse_solving(solver, A, c, dA, dc) );
  379. CALL_SUBTEST( check_sparse_solving(solver, halfA, c, dA, dc) );
  380. // check only once
  381. if(i==0)
  382. {
  383. b = DenseVector::Zero(size);
  384. check_sparse_solving(solver, A, b, dA, b);
  385. }
  386. }
  387. // First, get the folder
  388. #ifdef TEST_REAL_CASES
  389. // Test real problems with double precision only
  390. if (internal::is_same<typename NumTraits<Scalar>::Real, double>::value)
  391. {
  392. std::string mat_folder = get_matrixfolder<Scalar>();
  393. MatrixMarketIterator<Scalar> it(mat_folder);
  394. for (; it; ++it)
  395. {
  396. if (it.sym() == SPD){
  397. A = it.matrix();
  398. if(A.diagonal().size() <= maxRealWorldSize)
  399. {
  400. DenseVector b = it.rhs();
  401. DenseVector refX = it.refX();
  402. PermutationMatrix<Dynamic, Dynamic, StorageIndex> pnull;
  403. halfA.resize(A.rows(), A.cols());
  404. if(Solver::UpLo == (Lower|Upper))
  405. halfA = A;
  406. else
  407. halfA.template selfadjointView<Solver::UpLo>() = A.template triangularView<Eigen::Lower>().twistedBy(pnull);
  408. std::cout << "INFO | Testing " << sym_to_string(it.sym()) << "sparse problem " << it.matname()
  409. << " (" << A.rows() << "x" << A.cols() << ") using " << typeid(Solver).name() << "..." << std::endl;
  410. CALL_SUBTEST( check_sparse_solving_real_cases(solver, A, b, A, refX) );
  411. std::string stats = solver_stats(solver);
  412. if(stats.size()>0)
  413. std::cout << "INFO | " << stats << std::endl;
  414. CALL_SUBTEST( check_sparse_solving_real_cases(solver, halfA, b, A, refX) );
  415. }
  416. else
  417. {
  418. std::cout << "INFO | Skip sparse problem \"" << it.matname() << "\" (too large)" << std::endl;
  419. }
  420. }
  421. }
  422. }
  423. #else
  424. EIGEN_UNUSED_VARIABLE(maxRealWorldSize);
  425. #endif
  426. }
  427. template<typename Solver> void check_sparse_spd_determinant(Solver& solver)
  428. {
  429. typedef typename Solver::MatrixType Mat;
  430. typedef typename Mat::Scalar Scalar;
  431. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  432. // generate the problem
  433. Mat A, halfA;
  434. DenseMatrix dA;
  435. generate_sparse_spd_problem(solver, A, halfA, dA, 30);
  436. for (int i = 0; i < g_repeat; i++) {
  437. check_sparse_determinant(solver, A, dA);
  438. check_sparse_determinant(solver, halfA, dA );
  439. }
  440. }
  441. template<typename Solver, typename DenseMat>
  442. Index generate_sparse_square_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300, int options = ForceNonZeroDiag)
  443. {
  444. typedef typename Solver::MatrixType Mat;
  445. typedef typename Mat::Scalar Scalar;
  446. Index size = internal::random<int>(1,maxSize);
  447. double density = (std::max)(8./(size*size), 0.01);
  448. A.resize(size,size);
  449. dA.resize(size,size);
  450. initSparse<Scalar>(density, dA, A, options);
  451. return size;
  452. }
  453. struct prune_column {
  454. Index m_col;
  455. prune_column(Index col) : m_col(col) {}
  456. template<class Scalar>
  457. bool operator()(Index, Index col, const Scalar&) const {
  458. return col != m_col;
  459. }
  460. };
  461. template<typename Solver> void check_sparse_square_solving(Solver& solver, int maxSize = 300, int maxRealWorldSize = 100000, bool checkDeficient = false)
  462. {
  463. typedef typename Solver::MatrixType Mat;
  464. typedef typename Mat::Scalar Scalar;
  465. typedef SparseMatrix<Scalar,ColMajor, typename Mat::StorageIndex> SpMat;
  466. typedef SparseVector<Scalar, 0, typename Mat::StorageIndex> SpVec;
  467. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  468. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  469. int rhsCols = internal::random<int>(1,16);
  470. Mat A;
  471. DenseMatrix dA;
  472. for (int i = 0; i < g_repeat; i++) {
  473. Index size = generate_sparse_square_problem(solver, A, dA, maxSize);
  474. A.makeCompressed();
  475. DenseVector b = DenseVector::Random(size);
  476. DenseMatrix dB(size,rhsCols);
  477. SpMat B(size,rhsCols);
  478. double density = (std::max)(8./(size*rhsCols), 0.1);
  479. initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
  480. B.makeCompressed();
  481. SpVec c = B.col(0);
  482. DenseVector dc = dB.col(0);
  483. CALL_SUBTEST(check_sparse_solving(solver, A, b, dA, b));
  484. CALL_SUBTEST(check_sparse_solving(solver, A, dB, dA, dB));
  485. CALL_SUBTEST(check_sparse_solving(solver, A, B, dA, dB));
  486. CALL_SUBTEST(check_sparse_solving(solver, A, c, dA, dc));
  487. // check only once
  488. if(i==0)
  489. {
  490. CALL_SUBTEST(b = DenseVector::Zero(size); check_sparse_solving(solver, A, b, dA, b));
  491. }
  492. // regression test for Bug 792 (structurally rank deficient matrices):
  493. if(checkDeficient && size>1) {
  494. Index col = internal::random<int>(0,int(size-1));
  495. A.prune(prune_column(col));
  496. solver.compute(A);
  497. VERIFY_IS_EQUAL(solver.info(), NumericalIssue);
  498. }
  499. }
  500. // First, get the folder
  501. #ifdef TEST_REAL_CASES
  502. // Test real problems with double precision only
  503. if (internal::is_same<typename NumTraits<Scalar>::Real, double>::value)
  504. {
  505. std::string mat_folder = get_matrixfolder<Scalar>();
  506. MatrixMarketIterator<Scalar> it(mat_folder);
  507. for (; it; ++it)
  508. {
  509. A = it.matrix();
  510. if(A.diagonal().size() <= maxRealWorldSize)
  511. {
  512. DenseVector b = it.rhs();
  513. DenseVector refX = it.refX();
  514. std::cout << "INFO | Testing " << sym_to_string(it.sym()) << "sparse problem " << it.matname()
  515. << " (" << A.rows() << "x" << A.cols() << ") using " << typeid(Solver).name() << "..." << std::endl;
  516. CALL_SUBTEST(check_sparse_solving_real_cases(solver, A, b, A, refX));
  517. std::string stats = solver_stats(solver);
  518. if(stats.size()>0)
  519. std::cout << "INFO | " << stats << std::endl;
  520. }
  521. else
  522. {
  523. std::cout << "INFO | SKIP sparse problem \"" << it.matname() << "\" (too large)" << std::endl;
  524. }
  525. }
  526. }
  527. #else
  528. EIGEN_UNUSED_VARIABLE(maxRealWorldSize);
  529. #endif
  530. }
  531. template<typename Solver> void check_sparse_square_determinant(Solver& solver)
  532. {
  533. typedef typename Solver::MatrixType Mat;
  534. typedef typename Mat::Scalar Scalar;
  535. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  536. for (int i = 0; i < g_repeat; i++) {
  537. // generate the problem
  538. Mat A;
  539. DenseMatrix dA;
  540. int size = internal::random<int>(1,30);
  541. dA.setRandom(size,size);
  542. dA = (dA.array().abs()<0.3).select(0,dA);
  543. dA.diagonal() = (dA.diagonal().array()==0).select(1,dA.diagonal());
  544. A = dA.sparseView();
  545. A.makeCompressed();
  546. check_sparse_determinant(solver, A, dA);
  547. }
  548. }
  549. template<typename Solver> void check_sparse_square_abs_determinant(Solver& solver)
  550. {
  551. typedef typename Solver::MatrixType Mat;
  552. typedef typename Mat::Scalar Scalar;
  553. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  554. for (int i = 0; i < g_repeat; i++) {
  555. // generate the problem
  556. Mat A;
  557. DenseMatrix dA;
  558. generate_sparse_square_problem(solver, A, dA, 30);
  559. A.makeCompressed();
  560. check_sparse_abs_determinant(solver, A, dA);
  561. }
  562. }
  563. template<typename Solver, typename DenseMat>
  564. void generate_sparse_leastsquare_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300, int options = ForceNonZeroDiag)
  565. {
  566. typedef typename Solver::MatrixType Mat;
  567. typedef typename Mat::Scalar Scalar;
  568. int rows = internal::random<int>(1,maxSize);
  569. int cols = internal::random<int>(1,rows);
  570. double density = (std::max)(8./(rows*cols), 0.01);
  571. A.resize(rows,cols);
  572. dA.resize(rows,cols);
  573. initSparse<Scalar>(density, dA, A, options);
  574. }
  575. template<typename Solver> void check_sparse_leastsquare_solving(Solver& solver)
  576. {
  577. typedef typename Solver::MatrixType Mat;
  578. typedef typename Mat::Scalar Scalar;
  579. typedef SparseMatrix<Scalar,ColMajor, typename Mat::StorageIndex> SpMat;
  580. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  581. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  582. int rhsCols = internal::random<int>(1,16);
  583. Mat A;
  584. DenseMatrix dA;
  585. for (int i = 0; i < g_repeat; i++) {
  586. generate_sparse_leastsquare_problem(solver, A, dA);
  587. A.makeCompressed();
  588. DenseVector b = DenseVector::Random(A.rows());
  589. DenseMatrix dB(A.rows(),rhsCols);
  590. SpMat B(A.rows(),rhsCols);
  591. double density = (std::max)(8./(A.rows()*rhsCols), 0.1);
  592. initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
  593. B.makeCompressed();
  594. check_sparse_solving(solver, A, b, dA, b);
  595. check_sparse_solving(solver, A, dB, dA, dB);
  596. check_sparse_solving(solver, A, B, dA, dB);
  597. // check only once
  598. if(i==0)
  599. {
  600. b = DenseVector::Zero(A.rows());
  601. check_sparse_solving(solver, A, b, dA, b);
  602. }
  603. }
  604. }