sparse_solver.h 19 KB

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