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.

533 lines
16 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.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 <iostream>
  10. #include <fstream>
  11. #include "Eigen/SparseCore"
  12. #include <bench/BenchTimer.h>
  13. #include <cstdlib>
  14. #include <string>
  15. #include <Eigen/Cholesky>
  16. #include <Eigen/Jacobi>
  17. #include <Eigen/Householder>
  18. #include <Eigen/IterativeLinearSolvers>
  19. #include <unsupported/Eigen/IterativeSolvers>
  20. #include <Eigen/LU>
  21. #include <unsupported/Eigen/SparseExtra>
  22. #ifdef EIGEN_CHOLMOD_SUPPORT
  23. #include <Eigen/CholmodSupport>
  24. #endif
  25. #ifdef EIGEN_UMFPACK_SUPPORT
  26. #include <Eigen/UmfPackSupport>
  27. #endif
  28. #ifdef EIGEN_PARDISO_SUPPORT
  29. #include <Eigen/PardisoSupport>
  30. #endif
  31. #ifdef EIGEN_SUPERLU_SUPPORT
  32. #include <Eigen/SuperLUSupport>
  33. #endif
  34. #ifdef EIGEN_PASTIX_SUPPORT
  35. #include <Eigen/PaStiXSupport>
  36. #endif
  37. // CONSTANTS
  38. #define EIGEN_UMFPACK 0
  39. #define EIGEN_SUPERLU 1
  40. #define EIGEN_PASTIX 2
  41. #define EIGEN_PARDISO 3
  42. #define EIGEN_BICGSTAB 4
  43. #define EIGEN_BICGSTAB_ILUT 5
  44. #define EIGEN_GMRES 6
  45. #define EIGEN_GMRES_ILUT 7
  46. #define EIGEN_SIMPLICIAL_LDLT 8
  47. #define EIGEN_CHOLMOD_LDLT 9
  48. #define EIGEN_PASTIX_LDLT 10
  49. #define EIGEN_PARDISO_LDLT 11
  50. #define EIGEN_SIMPLICIAL_LLT 12
  51. #define EIGEN_CHOLMOD_SUPERNODAL_LLT 13
  52. #define EIGEN_CHOLMOD_SIMPLICIAL_LLT 14
  53. #define EIGEN_PASTIX_LLT 15
  54. #define EIGEN_PARDISO_LLT 16
  55. #define EIGEN_CG 17
  56. #define EIGEN_CG_PRECOND 18
  57. #define EIGEN_ALL_SOLVERS 19
  58. using namespace Eigen;
  59. using namespace std;
  60. struct Stats{
  61. ComputationInfo info;
  62. double total_time;
  63. double compute_time;
  64. double solve_time;
  65. double rel_error;
  66. int memory_used;
  67. int iterations;
  68. int isavail;
  69. int isIterative;
  70. };
  71. // Global variables for input parameters
  72. int MaximumIters; // Maximum number of iterations
  73. double RelErr; // Relative error of the computed solution
  74. template<typename T> inline typename NumTraits<T>::Real test_precision() { return NumTraits<T>::dummy_precision(); }
  75. template<> inline float test_precision<float>() { return 1e-3f; }
  76. template<> inline double test_precision<double>() { return 1e-6; }
  77. template<> inline float test_precision<std::complex<float> >() { return test_precision<float>(); }
  78. template<> inline double test_precision<std::complex<double> >() { return test_precision<double>(); }
  79. void printStatheader(std::ofstream& out)
  80. {
  81. int LUcnt = 0;
  82. string LUlist =" ", LLTlist = "<TH > LLT", LDLTlist = "<TH > LDLT ";
  83. #ifdef EIGEN_UMFPACK_SUPPORT
  84. LUlist += "<TH > UMFPACK "; LUcnt++;
  85. #endif
  86. #ifdef EIGEN_SUPERLU_SUPPORT
  87. LUlist += "<TH > SUPERLU "; LUcnt++;
  88. #endif
  89. #ifdef EIGEN_CHOLMOD_SUPPORT
  90. LLTlist += "<TH > CHOLMOD SP LLT<TH > CHOLMOD LLT";
  91. LDLTlist += "<TH>CHOLMOD LDLT";
  92. #endif
  93. #ifdef EIGEN_PARDISO_SUPPORT
  94. LUlist += "<TH > PARDISO LU"; LUcnt++;
  95. LLTlist += "<TH > PARDISO LLT";
  96. LDLTlist += "<TH > PARDISO LDLT";
  97. #endif
  98. #ifdef EIGEN_PASTIX_SUPPORT
  99. LUlist += "<TH > PASTIX LU"; LUcnt++;
  100. LLTlist += "<TH > PASTIX LLT";
  101. LDLTlist += "<TH > PASTIX LDLT";
  102. #endif
  103. out << "<TABLE border=\"1\" >\n ";
  104. out << "<TR><TH>Matrix <TH> N <TH> NNZ <TH> ";
  105. if (LUcnt) out << LUlist;
  106. out << " <TH >BiCGSTAB <TH >BiCGSTAB+ILUT"<< "<TH >GMRES+ILUT" <<LDLTlist << LLTlist << "<TH> CG "<< std::endl;
  107. }
  108. template<typename Solver, typename Scalar>
  109. Stats call_solver(Solver &solver, const typename Solver::MatrixType& A, const Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX)
  110. {
  111. Stats stat;
  112. Matrix<Scalar, Dynamic, 1> x;
  113. BenchTimer timer;
  114. timer.reset();
  115. timer.start();
  116. solver.compute(A);
  117. if (solver.info() != Success)
  118. {
  119. stat.info = NumericalIssue;
  120. std::cerr << "Solver failed ... \n";
  121. return stat;
  122. }
  123. timer.stop();
  124. stat.compute_time = timer.value();
  125. timer.reset();
  126. timer.start();
  127. x = solver.solve(b);
  128. if (solver.info() == NumericalIssue)
  129. {
  130. stat.info = NumericalIssue;
  131. std::cerr << "Solver failed ... \n";
  132. return stat;
  133. }
  134. timer.stop();
  135. stat.solve_time = timer.value();
  136. stat.total_time = stat.solve_time + stat.compute_time;
  137. stat.memory_used = 0;
  138. // Verify the relative error
  139. if(refX.size() != 0)
  140. stat.rel_error = (refX - x).norm()/refX.norm();
  141. else
  142. {
  143. // Compute the relative residual norm
  144. Matrix<Scalar, Dynamic, 1> temp;
  145. temp = A * x;
  146. stat.rel_error = (b-temp).norm()/b.norm();
  147. }
  148. if ( stat.rel_error > RelErr )
  149. {
  150. stat.info = NoConvergence;
  151. return stat;
  152. }
  153. else
  154. {
  155. stat.info = Success;
  156. return stat;
  157. }
  158. }
  159. template<typename Solver, typename Scalar>
  160. Stats call_directsolver(Solver& solver, const typename Solver::MatrixType& A, const Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX)
  161. {
  162. Stats stat;
  163. stat = call_solver(solver, A, b, refX);
  164. return stat;
  165. }
  166. template<typename Solver, typename Scalar>
  167. Stats call_itersolver(Solver &solver, const typename Solver::MatrixType& A, const Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX)
  168. {
  169. Stats stat;
  170. solver.setTolerance(RelErr);
  171. solver.setMaxIterations(MaximumIters);
  172. stat = call_solver(solver, A, b, refX);
  173. stat.iterations = solver.iterations();
  174. return stat;
  175. }
  176. inline void printStatItem(Stats *stat, int solver_id, int& best_time_id, double& best_time_val)
  177. {
  178. stat[solver_id].isavail = 1;
  179. if (stat[solver_id].info == NumericalIssue)
  180. {
  181. cout << " SOLVER FAILED ... Probably a numerical issue \n";
  182. return;
  183. }
  184. if (stat[solver_id].info == NoConvergence){
  185. cout << "REL. ERROR " << stat[solver_id].rel_error;
  186. if(stat[solver_id].isIterative == 1)
  187. cout << " (" << stat[solver_id].iterations << ") \n";
  188. return;
  189. }
  190. // Record the best CPU time
  191. if (!best_time_val)
  192. {
  193. best_time_val = stat[solver_id].total_time;
  194. best_time_id = solver_id;
  195. }
  196. else if (stat[solver_id].total_time < best_time_val)
  197. {
  198. best_time_val = stat[solver_id].total_time;
  199. best_time_id = solver_id;
  200. }
  201. // Print statistics to standard output
  202. if (stat[solver_id].info == Success){
  203. cout<< "COMPUTE TIME : " << stat[solver_id].compute_time<< " \n";
  204. cout<< "SOLVE TIME : " << stat[solver_id].solve_time<< " \n";
  205. cout<< "TOTAL TIME : " << stat[solver_id].total_time<< " \n";
  206. cout << "REL. ERROR : " << stat[solver_id].rel_error ;
  207. if(stat[solver_id].isIterative == 1) {
  208. cout << " (" << stat[solver_id].iterations << ") ";
  209. }
  210. cout << std::endl;
  211. }
  212. }
  213. /* Print the results from all solvers corresponding to a particular matrix
  214. * The best CPU time is printed in bold
  215. */
  216. inline void printHtmlStatLine(Stats *stat, int best_time_id, string& statline)
  217. {
  218. string markup;
  219. ostringstream compute,solve,total,error;
  220. for (int i = 0; i < EIGEN_ALL_SOLVERS; i++)
  221. {
  222. if (stat[i].isavail == 0) continue;
  223. if(i == best_time_id)
  224. markup = "<TD style=\"background-color:red\">";
  225. else
  226. markup = "<TD>";
  227. if (stat[i].info == Success){
  228. compute << markup << stat[i].compute_time;
  229. solve << markup << stat[i].solve_time;
  230. total << markup << stat[i].total_time;
  231. error << " <TD> " << stat[i].rel_error;
  232. if(stat[i].isIterative == 1) {
  233. error << " (" << stat[i].iterations << ") ";
  234. }
  235. }
  236. else {
  237. compute << " <TD> -" ;
  238. solve << " <TD> -" ;
  239. total << " <TD> -" ;
  240. if(stat[i].info == NoConvergence){
  241. error << " <TD> "<< stat[i].rel_error ;
  242. if(stat[i].isIterative == 1)
  243. error << " (" << stat[i].iterations << ") ";
  244. }
  245. else error << " <TD> - ";
  246. }
  247. }
  248. statline = "<TH>Compute Time " + compute.str() + "\n"
  249. + "<TR><TH>Solve Time " + solve.str() + "\n"
  250. + "<TR><TH>Total Time " + total.str() + "\n"
  251. +"<TR><TH>Error(Iter)" + error.str() + "\n";
  252. }
  253. template <typename Scalar>
  254. int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX, Stats *stat)
  255. {
  256. typedef SparseMatrix<Scalar, ColMajor> SpMat;
  257. // First, deal with Nonsymmetric and symmetric matrices
  258. int best_time_id = 0;
  259. double best_time_val = 0.0;
  260. //UMFPACK
  261. #ifdef EIGEN_UMFPACK_SUPPORT
  262. {
  263. cout << "Solving with UMFPACK LU ... \n";
  264. UmfPackLU<SpMat> solver;
  265. stat[EIGEN_UMFPACK] = call_directsolver(solver, A, b, refX);
  266. printStatItem(stat, EIGEN_UMFPACK, best_time_id, best_time_val);
  267. }
  268. #endif
  269. //SuperLU
  270. #ifdef EIGEN_SUPERLU_SUPPORT
  271. {
  272. cout << "\nSolving with SUPERLU ... \n";
  273. SuperLU<SpMat> solver;
  274. stat[EIGEN_SUPERLU] = call_directsolver(solver, A, b, refX);
  275. printStatItem(stat, EIGEN_SUPERLU, best_time_id, best_time_val);
  276. }
  277. #endif
  278. // PaStix LU
  279. #ifdef EIGEN_PASTIX_SUPPORT
  280. {
  281. cout << "\nSolving with PASTIX LU ... \n";
  282. PastixLU<SpMat> solver;
  283. stat[EIGEN_PASTIX] = call_directsolver(solver, A, b, refX) ;
  284. printStatItem(stat, EIGEN_PASTIX, best_time_id, best_time_val);
  285. }
  286. #endif
  287. //PARDISO LU
  288. #ifdef EIGEN_PARDISO_SUPPORT
  289. {
  290. cout << "\nSolving with PARDISO LU ... \n";
  291. PardisoLU<SpMat> solver;
  292. stat[EIGEN_PARDISO] = call_directsolver(solver, A, b, refX);
  293. printStatItem(stat, EIGEN_PARDISO, best_time_id, best_time_val);
  294. }
  295. #endif
  296. //BiCGSTAB
  297. {
  298. cout << "\nSolving with BiCGSTAB ... \n";
  299. BiCGSTAB<SpMat> solver;
  300. stat[EIGEN_BICGSTAB] = call_itersolver(solver, A, b, refX);
  301. stat[EIGEN_BICGSTAB].isIterative = 1;
  302. printStatItem(stat, EIGEN_BICGSTAB, best_time_id, best_time_val);
  303. }
  304. //BiCGSTAB+ILUT
  305. {
  306. cout << "\nSolving with BiCGSTAB and ILUT ... \n";
  307. BiCGSTAB<SpMat, IncompleteLUT<Scalar> > solver;
  308. stat[EIGEN_BICGSTAB_ILUT] = call_itersolver(solver, A, b, refX);
  309. stat[EIGEN_BICGSTAB_ILUT].isIterative = 1;
  310. printStatItem(stat, EIGEN_BICGSTAB_ILUT, best_time_id, best_time_val);
  311. }
  312. //GMRES
  313. // {
  314. // cout << "\nSolving with GMRES ... \n";
  315. // GMRES<SpMat> solver;
  316. // stat[EIGEN_GMRES] = call_itersolver(solver, A, b, refX);
  317. // stat[EIGEN_GMRES].isIterative = 1;
  318. // printStatItem(stat, EIGEN_GMRES, best_time_id, best_time_val);
  319. // }
  320. //GMRES+ILUT
  321. {
  322. cout << "\nSolving with GMRES and ILUT ... \n";
  323. GMRES<SpMat, IncompleteLUT<Scalar> > solver;
  324. stat[EIGEN_GMRES_ILUT] = call_itersolver(solver, A, b, refX);
  325. stat[EIGEN_GMRES_ILUT].isIterative = 1;
  326. printStatItem(stat, EIGEN_GMRES_ILUT, best_time_id, best_time_val);
  327. }
  328. // Hermitian and not necessarily positive-definites
  329. if (sym != NonSymmetric)
  330. {
  331. // Internal Cholesky
  332. {
  333. cout << "\nSolving with Simplicial LDLT ... \n";
  334. SimplicialLDLT<SpMat, Lower> solver;
  335. stat[EIGEN_SIMPLICIAL_LDLT] = call_directsolver(solver, A, b, refX);
  336. printStatItem(stat, EIGEN_SIMPLICIAL_LDLT, best_time_id, best_time_val);
  337. }
  338. // CHOLMOD
  339. #ifdef EIGEN_CHOLMOD_SUPPORT
  340. {
  341. cout << "\nSolving with CHOLMOD LDLT ... \n";
  342. CholmodDecomposition<SpMat, Lower> solver;
  343. solver.setMode(CholmodLDLt);
  344. stat[EIGEN_CHOLMOD_LDLT] = call_directsolver(solver, A, b, refX);
  345. printStatItem(stat,EIGEN_CHOLMOD_LDLT, best_time_id, best_time_val);
  346. }
  347. #endif
  348. //PASTIX LLT
  349. #ifdef EIGEN_PASTIX_SUPPORT
  350. {
  351. cout << "\nSolving with PASTIX LDLT ... \n";
  352. PastixLDLT<SpMat, Lower> solver;
  353. stat[EIGEN_PASTIX_LDLT] = call_directsolver(solver, A, b, refX);
  354. printStatItem(stat,EIGEN_PASTIX_LDLT, best_time_id, best_time_val);
  355. }
  356. #endif
  357. //PARDISO LLT
  358. #ifdef EIGEN_PARDISO_SUPPORT
  359. {
  360. cout << "\nSolving with PARDISO LDLT ... \n";
  361. PardisoLDLT<SpMat, Lower> solver;
  362. stat[EIGEN_PARDISO_LDLT] = call_directsolver(solver, A, b, refX);
  363. printStatItem(stat,EIGEN_PARDISO_LDLT, best_time_id, best_time_val);
  364. }
  365. #endif
  366. }
  367. // Now, symmetric POSITIVE DEFINITE matrices
  368. if (sym == SPD)
  369. {
  370. //Internal Sparse Cholesky
  371. {
  372. cout << "\nSolving with SIMPLICIAL LLT ... \n";
  373. SimplicialLLT<SpMat, Lower> solver;
  374. stat[EIGEN_SIMPLICIAL_LLT] = call_directsolver(solver, A, b, refX);
  375. printStatItem(stat,EIGEN_SIMPLICIAL_LLT, best_time_id, best_time_val);
  376. }
  377. // CHOLMOD
  378. #ifdef EIGEN_CHOLMOD_SUPPORT
  379. {
  380. // CholMOD SuperNodal LLT
  381. cout << "\nSolving with CHOLMOD LLT (Supernodal)... \n";
  382. CholmodDecomposition<SpMat, Lower> solver;
  383. solver.setMode(CholmodSupernodalLLt);
  384. stat[EIGEN_CHOLMOD_SUPERNODAL_LLT] = call_directsolver(solver, A, b, refX);
  385. printStatItem(stat,EIGEN_CHOLMOD_SUPERNODAL_LLT, best_time_id, best_time_val);
  386. // CholMod Simplicial LLT
  387. cout << "\nSolving with CHOLMOD LLT (Simplicial) ... \n";
  388. solver.setMode(CholmodSimplicialLLt);
  389. stat[EIGEN_CHOLMOD_SIMPLICIAL_LLT] = call_directsolver(solver, A, b, refX);
  390. printStatItem(stat,EIGEN_CHOLMOD_SIMPLICIAL_LLT, best_time_id, best_time_val);
  391. }
  392. #endif
  393. //PASTIX LLT
  394. #ifdef EIGEN_PASTIX_SUPPORT
  395. {
  396. cout << "\nSolving with PASTIX LLT ... \n";
  397. PastixLLT<SpMat, Lower> solver;
  398. stat[EIGEN_PASTIX_LLT] = call_directsolver(solver, A, b, refX);
  399. printStatItem(stat,EIGEN_PASTIX_LLT, best_time_id, best_time_val);
  400. }
  401. #endif
  402. //PARDISO LLT
  403. #ifdef EIGEN_PARDISO_SUPPORT
  404. {
  405. cout << "\nSolving with PARDISO LLT ... \n";
  406. PardisoLLT<SpMat, Lower> solver;
  407. stat[EIGEN_PARDISO_LLT] = call_directsolver(solver, A, b, refX);
  408. printStatItem(stat,EIGEN_PARDISO_LLT, best_time_id, best_time_val);
  409. }
  410. #endif
  411. // Internal CG
  412. {
  413. cout << "\nSolving with CG ... \n";
  414. ConjugateGradient<SpMat, Lower> solver;
  415. stat[EIGEN_CG] = call_itersolver(solver, A, b, refX);
  416. stat[EIGEN_CG].isIterative = 1;
  417. printStatItem(stat,EIGEN_CG, best_time_id, best_time_val);
  418. }
  419. //CG+IdentityPreconditioner
  420. // {
  421. // cout << "\nSolving with CG and IdentityPreconditioner ... \n";
  422. // ConjugateGradient<SpMat, Lower, IdentityPreconditioner> solver;
  423. // stat[EIGEN_CG_PRECOND] = call_itersolver(solver, A, b, refX);
  424. // stat[EIGEN_CG_PRECOND].isIterative = 1;
  425. // printStatItem(stat,EIGEN_CG_PRECOND, best_time_id, best_time_val);
  426. // }
  427. } // End SPD matrices
  428. return best_time_id;
  429. }
  430. /* Browse all the matrices available in the specified folder
  431. * and solve the associated linear system.
  432. * The results of each solve are printed in the standard output
  433. * and optionally in the provided html file
  434. */
  435. template <typename Scalar>
  436. void Browse_Matrices(const string folder, bool statFileExists, std::string& statFile, int maxiters, double tol)
  437. {
  438. MaximumIters = maxiters; // Maximum number of iterations, global variable
  439. RelErr = tol; //Relative residual error as stopping criterion for iterative solvers
  440. MatrixMarketIterator<Scalar> it(folder);
  441. Stats stat[EIGEN_ALL_SOLVERS];
  442. for ( ; it; ++it)
  443. {
  444. for (int i = 0; i < EIGEN_ALL_SOLVERS; i++)
  445. {
  446. stat[i].isavail = 0;
  447. stat[i].isIterative = 0;
  448. }
  449. int best_time_id;
  450. cout<< "\n\n===================================================== \n";
  451. cout<< " ====== SOLVING WITH MATRIX " << it.matname() << " ====\n";
  452. cout<< " =================================================== \n\n";
  453. Matrix<Scalar, Dynamic, 1> refX;
  454. if(it.hasrefX()) refX = it.refX();
  455. best_time_id = SelectSolvers<Scalar>(it.matrix(), it.sym(), it.rhs(), refX, &stat[0]);
  456. if(statFileExists)
  457. {
  458. string statline;
  459. printHtmlStatLine(&stat[0], best_time_id, statline);
  460. std::ofstream statbuf(statFile.c_str(), std::ios::app);
  461. statbuf << "<TR><TH rowspan=\"4\">" << it.matname() << " <TD rowspan=\"4\"> "
  462. << it.matrix().rows() << " <TD rowspan=\"4\"> " << it.matrix().nonZeros()<< " "<< statline ;
  463. statbuf.close();
  464. }
  465. }
  466. }
  467. bool get_options(int argc, char **args, string option, string* value=0)
  468. {
  469. int idx = 1, found=false;
  470. while (idx<argc && !found){
  471. if (option.compare(args[idx]) == 0){
  472. found = true;
  473. if(value) *value = args[idx+1];
  474. }
  475. idx+=2;
  476. }
  477. return found;
  478. }