Read real and complex bench matrices from a unique folder

Output and display bench results using XML and XSLT
This commit is contained in:
Desire NUENTSA W. 2012-08-27 22:52:43 +02:00
parent ebe511334f
commit fe9956defe
5 changed files with 329 additions and 214 deletions

31
bench/spbench/spbench.dtd Normal file
View File

@ -0,0 +1,31 @@
<!ELEMENT BENCH (AVAILSOLVER+,LINEARSYSTEM+)>
<!ELEMENT AVAILSOLVER (SOLVER+)>
<!ELEMENT SOLVER (TYPE,PACKAGE)>
<!ELEMENT TYPE (#PCDATA)> <!-- One of LU, LLT, LDLT, ITER -->
<!ELEMENT PACKAGE (#PCDATA)> <!-- Derived from a library -->
<!ELEMENT LINEARSYSTEM (MATRIX,SOLVER_STAT+,BEST_SOLVER,GLOBAL_PARAMS*)>
<!ELEMENT MATRIX (NAME,SIZE,ENTRIES,PATTERN?,SYMMETRY,POSDEF?,ARITHMETIC,RHS*)>
<!ELEMENT NAME (#PCDATA)>
<!ELEMENT SIZE (#PCDATA)>
<!ELEMENT ENTRIES (#PCDATA)> <!-- The number of nonzeros elements -->
<!ELEMENT PATTERN (#PCDATA)> <!-- Is structural pattern symmetric or not -->
<!ELEMENT SYMMETRY (#PCDATA)> <!-- symmmetry with numerical values -->
<!ELEMENT POSDEF (#PCDATA)> <!-- Is the matrix positive definite or not -->
<!ELEMENT ARITHMETIC (#PCDATA)>
<!ELEMENT RHS (SOURCE)> <!-- A matrix can have one or more right hand side associated. -->
<!ELEMENT SOURCE (#PCDATA)> <!-- Source of the right hand side, either generated or provided -->
<!ELEMENT SOLVER_STAT (PARAMS*,TIME,ERROR,ITER?)>
<!ELEMENT PARAMS (#PCDATA)>
<!ELEMENT TIME (COMPUTE,SOLVE,TOTAL)>
<!ELEMENT COMPUTE (#PCDATA)> <!-- Time to analyze,to factorize, or to setup the preconditioner-->
<!ELEMENT SOLVE (#PCDATA)> <!-- Time to solve with all the available rhs -->
<!ELEMENT TOTAL (#PCDATA)>
<!ELEMENT ERROR (#PCDATA)> <!-- Either the relative error or the relative residual norm -->
<!ELEMENT ITER (#PCDATA)> <!-- Number of iterations -->
<!ELEMENT BEST_SOLVER CDATA> <!-- Id of the best solver -->
<!ELEMENT GLOBAL_PARAMS (#PCDATA)> <!-- Parameters shared by all solvers -->
<!ATTLIST SOLVER ID CDATA #REQUIRED>
<!ATTLIST SOLVER_STAT ID CDATA #REQUIRED>
<!ATTLIST BEST_SOLVER ID CDATA #REQUIRED>
<!ATTLIST RHS ID CDATA #IMPLIED>

83
bench/spbench/spbench.xsl Normal file
View File

@ -0,0 +1,83 @@
<?xml version="1.0" encoding="UTF-8"?>
<xsl:stylesheet version="1.0"
xmlns:xsl="http://www.w3.org/1999/XSL/Transform" >
<!-- Desire Nuentsa, Inria -->
<xsl:output method="html" indent="no"/>
<xsl:template match="/"> <!-- Root of the document -->
<html>
<head>
<style type="text/css">
td { white-space: nowrap;}
</style>
</head>
<body>
<table border="1" width="100%" height="100%">
<TR> <!-- Write the table header -->
<TH>Matrix</TH> <TH>N</TH> <TH> NNZ</TH> <TH> Sym</TH> <TH> SPD</TH> <TH> </TH>
<xsl:for-each select="BENCH/AVAILSOLVER/SOLVER">
<xsl:sort select="@ID" data-type="number"/>
<TH>
<xsl:value-of select="TYPE" />
<xsl:text></xsl:text>
<xsl:value-of select="PACKAGE" />
<xsl:text></xsl:text>
</TH>
</xsl:for-each>
</TR>
<xsl:for-each select="BENCH/LINEARSYSTEM">
<TR> <!-- print statistics for one linear system-->
<TH rowspan="4"> <xsl:value-of select="MATRIX/NAME" /> </TH>
<TD rowspan="4"> <xsl:value-of select="MATRIX/SIZE" /> </TD>
<TD rowspan="4"> <xsl:value-of select="MATRIX/ENTRIES" /> </TD>
<TD rowspan="4"> <xsl:value-of select="MATRIX/SYMMETRY" /> </TD>
<TD rowspan="4"> <xsl:value-of select="MATRIX/POSDEF" /> </TD>
<TH> Compute Time </TH>
<xsl:for-each select="SOLVER_STAT">
<xsl:sort select="@ID" data-type="number"/>
<TD> <xsl:value-of select="TIME/COMPUTE" /> </TD>
</xsl:for-each>
</TR>
<TR>
<TH> Solve Time </TH>
<xsl:for-each select="SOLVER_STAT">
<xsl:sort select="@ID" data-type="number"/>
<TD> <xsl:value-of select="TIME/SOLVE" /> </TD>
</xsl:for-each>
</TR>
<TR>
<TH> Total Time </TH>
<xsl:for-each select="SOLVER_STAT">
<xsl:sort select="@ID" data-type="number"/>
<xsl:choose>
<xsl:when test="@ID=../BEST_SOLVER/@ID">
<TD style="background-color:red"> <xsl:value-of select="TIME/TOTAL" /> </TD>
</xsl:when>
<xsl:otherwise>
<TD> <xsl:value-of select="TIME/TOTAL" /></TD>
</xsl:otherwise>
</xsl:choose>
</xsl:for-each>
</TR>
<TR>
<TH> Error </TH>
<xsl:for-each select="SOLVER_STAT">
<xsl:sort select="@ID" data-type="number"/>
<TD> <xsl:value-of select="ERROR" />
<xsl:if test="ITER">
<xsl:text>(</xsl:text>
<xsl:value-of select="ITER" />
<xsl:text>)</xsl:text>
</xsl:if> </TD>
</xsl:for-each>
</TR>
</xsl:for-each>
</table>
</body>
</html>
</xsl:template>
</xsl:stylesheet>

View File

@ -5,7 +5,6 @@ void bench_printhelp()
cout<< " \nbenchsolver : performs a benchmark of all the solvers available in Eigen \n\n";
cout<< " MATRIX FOLDER : \n";
cout<< " The matrices for the benchmark should be collected in a folder specified with an environment variable EIGEN_MATRIXDIR \n";
cout<< " This folder should contain the subfolders real/ and complex/ : \n";
cout<< " The matrices are stored using the matrix market coordinate format \n";
cout<< " The matrix and associated right-hand side (rhs) files are named respectively \n";
cout<< " as MatrixName.mtx and MatrixName_b.mtx. If the rhs does not exist, a random one is generated. \n";
@ -68,18 +67,16 @@ int main(int argc, char ** args)
maxiters = atoi(inval.c_str());
string current_dir;
// Test the matrices in %EIGEN_MATRIXDIR/real
current_dir = matrix_dir + "/real";
Browse_Matrices<double>(current_dir, statFileExists, statFile,maxiters, tol);
// Test the real-arithmetics matrices
Browse_Matrices<double>(matrix_dir, statFileExists, statFile,maxiters, tol);
// Test the matrices in %EIGEN_MATRIXDIR/complex
current_dir = matrix_dir + "/complex";
Browse_Matrices<std::complex<double> >(current_dir, statFileExists, statFile, maxiters, tol);
// Test the complex-arithmetics matrices
Browse_Matrices<std::complex<double> >(matrix_dir, statFileExists, statFile, maxiters, tol);
if(statFileExists)
{
statbuf.open(statFile.c_str(), std::ios::app);
statbuf << "</TABLE> \n";
statbuf << "</BENCH> \n";
cout << "\n Output written in " << statFile << " ...\n";
statbuf.close();
}

View File

@ -67,21 +67,12 @@
using namespace Eigen;
using namespace std;
struct Stats{
ComputationInfo info;
double total_time;
double compute_time;
double solve_time;
double rel_error;
int memory_used;
int iterations;
int isavail;
int isIterative;
};
// Global variables for input parameters
int MaximumIters; // Maximum number of iterations
double RelErr; // Relative error of the computed solution
double best_time_val; // Current best time overall solvers
int best_time_id; // id of the best solver for the current system
template<typename T> inline typename NumTraits<T>::Real test_precision() { return NumTraits<T>::dummy_precision(); }
template<> inline float test_precision<float>() { return 1e-3f; }
@ -91,41 +82,123 @@ template<> inline double test_precision<std::complex<double> >() { return test_p
void printStatheader(std::ofstream& out)
{
int LUcnt = 0;
string LUlist =" ", LLTlist = "<TH > LLT", LDLTlist = "<TH > LDLT ";
// Print XML header
// NOTE It would have been much easier to write these XML documents using external libraries like tinyXML or Xerces-C++.
out << "<?xml version=\"1.0\" encoding=\"UTF-8\"?> \n";
out << "<?xml-stylesheet type=\"text/xsl\" href=\"spbench.xsl\" ?> \n";
out << "<!DOCTYPE BENCH SYSTEM \"spbench.dtd\"> \n";
out << "\n<!-- Generated by the Eigen library -->\n";
// Write the root XML element
out << "\n<BENCH> \n" ;
// List all available solvers
out << " <AVAILSOLVER> \n";
#ifdef EIGEN_UMFPACK_SUPPORT
LUlist += "<TH > UMFPACK "; LUcnt++;
out <<" <SOLVER ID=\"" << EIGEN_UMFPACK << "\">\n";
out << " <TYPE> LU </TYPE> \n";
out << " <PACKAGE> UMFPACK </PACKAGE> \n";
out << " </SOLVER> \n";
#endif
#ifdef EIGEN_SUPERLU_SUPPORT
LUlist += "<TH > SUPERLU "; LUcnt++;
out <<" <SOLVER ID=\"" << EIGEN_SUPERLU << "\">\n";
out << " <TYPE> LU </TYPE> \n";
out << " <PACKAGE> SUPERLU </PACKAGE> \n";
out << " </SOLVER> \n";
#endif
#ifdef EIGEN_CHOLMOD_SUPPORT
LLTlist += "<TH > CHOLMOD SP LLT<TH > CHOLMOD LLT";
LDLTlist += "<TH>CHOLMOD LDLT";
out <<" <SOLVER ID=\"" << EIGEN_CHOLMOD_SIMPLICIAL_LLT << "\">\n";
out << " <TYPE> LLT SP</TYPE> \n";
out << " <PACKAGE> CHOLMOD </PACKAGE> \n";
out << " </SOLVER> \n";
out <<" <SOLVER ID=\"" << EIGEN_CHOLMOD_SUPERNODAL_LLT << "\">\n";
out << " <TYPE> LLT</TYPE> \n";
out << " <PACKAGE> CHOLMOD </PACKAGE> \n";
out << " </SOLVER> \n";
out <<" <SOLVER ID=\"" << EIGEN_CHOLMOD_LDLT << "\">\n";
out << " <TYPE> LDLT </TYPE> \n";
out << " <PACKAGE> CHOLMOD </PACKAGE> \n";
out << " </SOLVER> \n";
#endif
#ifdef EIGEN_PARDISO_SUPPORT
LUlist += "<TH > PARDISO LU"; LUcnt++;
LLTlist += "<TH > PARDISO LLT";
LDLTlist += "<TH > PARDISO LDLT";
out <<" <SOLVER ID=\"" << EIGEN_PARDISO << "\">\n";
out << " <TYPE> LU </TYPE> \n";
out << " <PACKAGE> PARDISO </PACKAGE> \n";
out << " </SOLVER> \n";
out <<" <SOLVER ID=\"" << EIGEN_PARDISO_LLT << "\">\n";
out << " <TYPE> LLT </TYPE> \n";
out << " <PACKAGE> PARDISO </PACKAGE> \n";
out << " </SOLVER> \n";
out <<" <SOLVER ID=\"" << EIGEN_PARDISO_LDLT << "\">\n";
out << " <TYPE> LDLT </TYPE> \n";
out << " <PACKAGE> PARDISO </PACKAGE> \n";
out << " </SOLVER> \n";
#endif
#ifdef EIGEN_PASTIX_SUPPORT
LUlist += "<TH > PASTIX LU"; LUcnt++;
LLTlist += "<TH > PASTIX LLT";
LDLTlist += "<TH > PASTIX LDLT";
out <<" <SOLVER ID=\"" << EIGEN_PASTIX << "\">\n";
out << " <TYPE> LU </TYPE> \n";
out << " <PACKAGE> PASTIX </PACKAGE> \n";
out << " </SOLVER> \n";
out <<" <SOLVER ID=\"" << EIGEN_PASTIX_LLT << "\">\n";
out << " <TYPE> LLT </TYPE> \n";
out << " <PACKAGE> PASTIX </PACKAGE> \n";
out << " </SOLVER> \n";
out <<" <SOLVER ID=\"" << EIGEN_PASTIX_LDLT << "\">\n";
out << " <TYPE> LDLT </TYPE> \n";
out << " <PACKAGE> PASTIX </PACKAGE> \n";
out << " </SOLVER> \n";
#endif
out << "<TABLE border=\"1\" >\n ";
out << "<TR><TH>Matrix <TH> N <TH> NNZ <TH> ";
if (LUcnt) out << LUlist;
out << " <TH >BiCGSTAB <TH >BiCGSTAB+ILUT"<< "<TH >GMRES+ILUT" <<LDLTlist << LLTlist << "<TH> CG "<< std::endl;
out <<" <SOLVER ID=\"" << EIGEN_BICGSTAB << "\">\n";
out << " <TYPE> BICGSTAB </TYPE> \n";
out << " <PACKAGE> EIGEN </PACKAGE> \n";
out << " </SOLVER> \n";
out <<" <SOLVER ID=\"" << EIGEN_BICGSTAB_ILUT << "\">\n";
out << " <TYPE> BICGSTAB_ILUT </TYPE> \n";
out << " <PACKAGE> EIGEN </PACKAGE> \n";
out << " </SOLVER> \n";
out <<" <SOLVER ID=\"" << EIGEN_GMRES_ILUT << "\">\n";
out << " <TYPE> GMRES_ILUT </TYPE> \n";
out << " <PACKAGE> EIGEN </PACKAGE> \n";
out << " </SOLVER> \n";
out <<" <SOLVER ID=\"" << EIGEN_SIMPLICIAL_LDLT << "\">\n";
out << " <TYPE> LDLT </TYPE> \n";
out << " <PACKAGE> EIGEN </PACKAGE> \n";
out << " </SOLVER> \n";
out <<" <SOLVER ID=\"" << EIGEN_SIMPLICIAL_LLT << "\">\n";
out << " <TYPE> LLT </TYPE> \n";
out << " <PACKAGE> EIGEN </PACKAGE> \n";
out << " </SOLVER> \n";
out <<" <SOLVER ID=\"" << EIGEN_CG << "\">\n";
out << " <TYPE> CG </TYPE> \n";
out << " <PACKAGE> EIGEN </PACKAGE> \n";
out << " </SOLVER> \n";
out << " </AVAILSOLVER> \n";
}
template<typename Solver, typename Scalar>
Stats call_solver(Solver &solver, const typename Solver::MatrixType& A, const Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX)
void call_solver(Solver &solver, const int solver_id, const typename Solver::MatrixType& A, const Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX,std::ofstream& statbuf)
{
Stats stat;
double total_time;
double compute_time;
double solve_time;
double rel_error;
Matrix<Scalar, Dynamic, 1> x;
BenchTimer timer;
timer.reset();
@ -133,170 +206,95 @@ Stats call_solver(Solver &solver, const typename Solver::MatrixType& A, const Ma
solver.compute(A);
if (solver.info() != Success)
{
stat.info = NumericalIssue;
std::cerr << "Solver failed ... \n";
return stat;
return;
}
timer.stop();
stat.compute_time = timer.value();
timer.stop();
compute_time = timer.value();
statbuf << " <TIME>\n";
statbuf << " <COMPUTE> " << timer.value() << "</COMPUTE>\n";
std::cout<< "COMPUTE TIME : " << timer.value() <<std::endl;
timer.reset();
timer.start();
x = solver.solve(b);
if (solver.info() == NumericalIssue)
{
stat.info = NumericalIssue;
std::cerr << "Solver failed ... \n";
return stat;
return;
}
timer.stop();
stat.solve_time = timer.value();
stat.total_time = stat.solve_time + stat.compute_time;
stat.memory_used = 0;
solve_time = timer.value();
statbuf << " <SOLVE> " << timer.value() << "</SOLVE>\n";
std::cout<< "SOLVE TIME : " << timer.value() <<std::endl;
total_time = solve_time + compute_time;
statbuf << " <TOTAL> " << total_time << "</TOTAL>\n";
std::cout<< "TOTAL TIME : " << total_time <<std::endl;
statbuf << " </TIME>\n";
// Verify the relative error
if(refX.size() != 0)
stat.rel_error = (refX - x).norm()/refX.norm();
rel_error = (refX - x).norm()/refX.norm();
else
{
// Compute the relative residual norm
Matrix<Scalar, Dynamic, 1> temp;
temp = A * x;
stat.rel_error = (b-temp).norm()/b.norm();
rel_error = (b-temp).norm()/b.norm();
}
if ( stat.rel_error > RelErr )
statbuf << " <ERROR> " << rel_error << "</ERROR>\n";
std::cout<< "REL. ERROR : " << rel_error << "\n\n" ;
if ( rel_error <= RelErr )
{
stat.info = NoConvergence;
return stat;
}
else
{
stat.info = Success;
return stat;
// check the best time if convergence
if(!best_time_val || (best_time_val > total_time))
{
best_time_val = total_time;
best_time_id = solver_id;
}
}
}
template<typename Solver, typename Scalar>
Stats call_directsolver(Solver& solver, const typename Solver::MatrixType& A, const Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX)
void call_directsolver(Solver& solver, const int solver_id, const typename Solver::MatrixType& A, const Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX, std::string& statFile)
{
Stats stat;
stat = call_solver(solver, A, b, refX);
return stat;
std::ofstream statbuf(statFile.c_str(), std::ios::app);
statbuf << " <SOLVER_STAT ID=\"" << solver_id <<"\">\n";
call_solver(solver, solver_id, A, b, refX,statbuf);
statbuf << " </SOLVER_STAT>\n";
statbuf.close();
}
template<typename Solver, typename Scalar>
Stats call_itersolver(Solver &solver, const typename Solver::MatrixType& A, const Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX)
void call_itersolver(Solver &solver, const int solver_id, const typename Solver::MatrixType& A, const Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX, std::string& statFile)
{
Stats stat;
solver.setTolerance(RelErr);
solver.setMaxIterations(MaximumIters);
stat = call_solver(solver, A, b, refX);
stat.iterations = solver.iterations();
return stat;
}
inline void printStatItem(Stats *stat, int solver_id, int& best_time_id, double& best_time_val)
{
stat[solver_id].isavail = 1;
if (stat[solver_id].info == NumericalIssue)
{
cout << " SOLVER FAILED ... Probably a numerical issue \n";
return;
}
if (stat[solver_id].info == NoConvergence){
cout << "REL. ERROR " << stat[solver_id].rel_error;
if(stat[solver_id].isIterative == 1)
cout << " (" << stat[solver_id].iterations << ") \n";
return;
}
// Record the best CPU time
if (!best_time_val)
{
best_time_val = stat[solver_id].total_time;
best_time_id = solver_id;
}
else if (stat[solver_id].total_time < best_time_val)
{
best_time_val = stat[solver_id].total_time;
best_time_id = solver_id;
}
// Print statistics to standard output
if (stat[solver_id].info == Success){
cout<< "COMPUTE TIME : " << stat[solver_id].compute_time<< " \n";
cout<< "SOLVE TIME : " << stat[solver_id].solve_time<< " \n";
cout<< "TOTAL TIME : " << stat[solver_id].total_time<< " \n";
cout << "REL. ERROR : " << stat[solver_id].rel_error ;
if(stat[solver_id].isIterative == 1) {
cout << " (" << stat[solver_id].iterations << ") ";
}
cout << std::endl;
}
}
/* Print the results from all solvers corresponding to a particular matrix
* The best CPU time is printed in bold
*/
inline void printHtmlStatLine(Stats *stat, int best_time_id, string& statline)
{
string markup;
ostringstream compute,solve,total,error;
for (int i = 0; i < EIGEN_ALL_SOLVERS; i++)
{
if (stat[i].isavail == 0) continue;
if(i == best_time_id)
markup = "<TD style=\"background-color:red\">";
else
markup = "<TD>";
if (stat[i].info == Success){
compute << markup << stat[i].compute_time;
solve << markup << stat[i].solve_time;
total << markup << stat[i].total_time;
error << " <TD> " << stat[i].rel_error;
if(stat[i].isIterative == 1) {
error << " (" << stat[i].iterations << ") ";
}
}
else {
compute << " <TD> -" ;
solve << " <TD> -" ;
total << " <TD> -" ;
if(stat[i].info == NoConvergence){
error << " <TD> "<< stat[i].rel_error ;
if(stat[i].isIterative == 1)
error << " (" << stat[i].iterations << ") ";
}
else error << " <TD> - ";
}
}
statline = "<TH>Compute Time " + compute.str() + "\n"
+ "<TR><TH>Solve Time " + solve.str() + "\n"
+ "<TR><TH>Total Time " + total.str() + "\n"
+"<TR><TH>Error(Iter)" + error.str() + "\n";
std::ofstream statbuf(statFile.c_str(), std::ios::app);
statbuf << " <SOLVER_STAT ID=\"" << solver_id <<"\">\n";
call_solver(solver, solver_id, A, b, refX,statbuf);
statbuf << " <ITER> "<< solver.iterations() << "</ITER>\n";
statbuf << " </SOLVER_STAT>\n";
std::cout << "ITERATIONS : " << solver.iterations() <<"\n\n\n";
}
template <typename Scalar>
int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX, Stats *stat)
void SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar, Dynamic, 1>& b, const Matrix<Scalar, Dynamic, 1>& refX, std::string& statFile)
{
typedef SparseMatrix<Scalar, ColMajor> SpMat;
// First, deal with Nonsymmetric and symmetric matrices
int best_time_id = 0;
double best_time_val = 0.0;
best_time_id = 0;
best_time_val = 0.0;
//UMFPACK
#ifdef EIGEN_UMFPACK_SUPPORT
{
cout << "Solving with UMFPACK LU ... \n";
UmfPackLU<SpMat> solver;
stat[EIGEN_UMFPACK] = call_directsolver(solver, A, b, refX);
printStatItem(stat, EIGEN_UMFPACK, best_time_id, best_time_val);
call_directsolver(solver, EIGEN_UMFPACK, A, b, refX,statFile);
}
#endif
//SuperLU
@ -304,8 +302,8 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
{
cout << "\nSolving with SUPERLU ... \n";
SuperLU<SpMat> solver;
stat[EIGEN_SUPERLU] = call_directsolver(solver, A, b, refX);
printStatItem(stat, EIGEN_SUPERLU, best_time_id, best_time_val);
call_directsolver(solver, EIGEN_SUPERLU, A, b, refX,statFile);
printStatItem(stat, best_time_id, best_time_val);
}
#endif
@ -314,8 +312,7 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
{
cout << "\nSolving with PASTIX LU ... \n";
PastixLU<SpMat> solver;
stat[EIGEN_PASTIX] = call_directsolver(solver, A, b, refX) ;
printStatItem(stat, EIGEN_PASTIX, best_time_id, best_time_val);
call_directsolver(solver, EIGEN_PASTIX, A, b, refX,statFile) ;
}
#endif
@ -324,8 +321,7 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
{
cout << "\nSolving with PARDISO LU ... \n";
PardisoLU<SpMat> solver;
stat[EIGEN_PARDISO] = call_directsolver(solver, A, b, refX);
printStatItem(stat, EIGEN_PARDISO, best_time_id, best_time_val);
call_directsolver(solver, EIGEN_PARDISO, A, b, refX,statFile);
}
#endif
@ -335,17 +331,13 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
{
cout << "\nSolving with BiCGSTAB ... \n";
BiCGSTAB<SpMat> solver;
stat[EIGEN_BICGSTAB] = call_itersolver(solver, A, b, refX);
stat[EIGEN_BICGSTAB].isIterative = 1;
printStatItem(stat, EIGEN_BICGSTAB, best_time_id, best_time_val);
call_itersolver(solver, EIGEN_BICGSTAB, A, b, refX,statFile);
}
//BiCGSTAB+ILUT
{
cout << "\nSolving with BiCGSTAB and ILUT ... \n";
BiCGSTAB<SpMat, IncompleteLUT<Scalar> > solver;
stat[EIGEN_BICGSTAB_ILUT] = call_itersolver(solver, A, b, refX);
stat[EIGEN_BICGSTAB_ILUT].isIterative = 1;
printStatItem(stat, EIGEN_BICGSTAB_ILUT, best_time_id, best_time_val);
call_itersolver(solver, EIGEN_BICGSTAB_ILUT, A, b, refX,statFile);
}
@ -353,17 +345,13 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
// {
// cout << "\nSolving with GMRES ... \n";
// GMRES<SpMat> solver;
// stat[EIGEN_GMRES] = call_itersolver(solver, A, b, refX);
// stat[EIGEN_GMRES].isIterative = 1;
// printStatItem(stat, EIGEN_GMRES, best_time_id, best_time_val);
// call_itersolver(solver, EIGEN_GMRES, A, b, refX,statFile);
// }
//GMRES+ILUT
{
cout << "\nSolving with GMRES and ILUT ... \n";
GMRES<SpMat, IncompleteLUT<Scalar> > solver;
stat[EIGEN_GMRES_ILUT] = call_itersolver(solver, A, b, refX);
stat[EIGEN_GMRES_ILUT].isIterative = 1;
printStatItem(stat, EIGEN_GMRES_ILUT, best_time_id, best_time_val);
call_itersolver(solver, EIGEN_GMRES_ILUT, A, b, refX,statFile);
}
// Hermitian and not necessarily positive-definites
@ -373,8 +361,7 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
{
cout << "\nSolving with Simplicial LDLT ... \n";
SimplicialLDLT<SpMat, Lower> solver;
stat[EIGEN_SIMPLICIAL_LDLT] = call_directsolver(solver, A, b, refX);
printStatItem(stat, EIGEN_SIMPLICIAL_LDLT, best_time_id, best_time_val);
call_directsolver(solver, EIGEN_SIMPLICIAL_LDLT, A, b, refX,statFile);
}
// CHOLMOD
@ -383,8 +370,7 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
cout << "\nSolving with CHOLMOD LDLT ... \n";
CholmodDecomposition<SpMat, Lower> solver;
solver.setMode(CholmodLDLt);
stat[EIGEN_CHOLMOD_LDLT] = call_directsolver(solver, A, b, refX);
printStatItem(stat,EIGEN_CHOLMOD_LDLT, best_time_id, best_time_val);
call_directsolver(solver,EIGEN_CHOLMOD_LDLT, A, b, refX,statFile);
}
#endif
@ -393,8 +379,7 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
{
cout << "\nSolving with PASTIX LDLT ... \n";
PastixLDLT<SpMat, Lower> solver;
stat[EIGEN_PASTIX_LDLT] = call_directsolver(solver, A, b, refX);
printStatItem(stat,EIGEN_PASTIX_LDLT, best_time_id, best_time_val);
call_directsolver(solver,EIGEN_PASTIX_LDLT, A, b, refX,statFile);
}
#endif
@ -403,8 +388,7 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
{
cout << "\nSolving with PARDISO LDLT ... \n";
PardisoLDLT<SpMat, Lower> solver;
stat[EIGEN_PARDISO_LDLT] = call_directsolver(solver, A, b, refX);
printStatItem(stat,EIGEN_PARDISO_LDLT, best_time_id, best_time_val);
call_directsolver(solver,EIGEN_PARDISO_LDLT, A, b, refX,statFile);
}
#endif
}
@ -417,8 +401,7 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
{
cout << "\nSolving with SIMPLICIAL LLT ... \n";
SimplicialLLT<SpMat, Lower> solver;
stat[EIGEN_SIMPLICIAL_LLT] = call_directsolver(solver, A, b, refX);
printStatItem(stat,EIGEN_SIMPLICIAL_LLT, best_time_id, best_time_val);
call_directsolver(solver,EIGEN_SIMPLICIAL_LLT, A, b, refX,statFile);
}
// CHOLMOD
@ -428,13 +411,11 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
cout << "\nSolving with CHOLMOD LLT (Supernodal)... \n";
CholmodDecomposition<SpMat, Lower> solver;
solver.setMode(CholmodSupernodalLLt);
stat[EIGEN_CHOLMOD_SUPERNODAL_LLT] = call_directsolver(solver, A, b, refX);
printStatItem(stat,EIGEN_CHOLMOD_SUPERNODAL_LLT, best_time_id, best_time_val);
call_directsolver(solver,EIGEN_CHOLMOD_SUPERNODAL_LLT, A, b, refX,statFile);
// CholMod Simplicial LLT
cout << "\nSolving with CHOLMOD LLT (Simplicial) ... \n";
solver.setMode(CholmodSimplicialLLt);
stat[EIGEN_CHOLMOD_SIMPLICIAL_LLT] = call_directsolver(solver, A, b, refX);
printStatItem(stat,EIGEN_CHOLMOD_SIMPLICIAL_LLT, best_time_id, best_time_val);
call_directsolver(solver,EIGEN_CHOLMOD_SIMPLICIAL_LLT, A, b, refX,statFile);
}
#endif
@ -443,8 +424,7 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
{
cout << "\nSolving with PASTIX LLT ... \n";
PastixLLT<SpMat, Lower> solver;
stat[EIGEN_PASTIX_LLT] = call_directsolver(solver, A, b, refX);
printStatItem(stat,EIGEN_PASTIX_LLT, best_time_id, best_time_val);
call_directsolver(solver,EIGEN_PASTIX_LLT, A, b, refX,statFile);
}
#endif
@ -453,8 +433,7 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
{
cout << "\nSolving with PARDISO LLT ... \n";
PardisoLLT<SpMat, Lower> solver;
stat[EIGEN_PARDISO_LLT] = call_directsolver(solver, A, b, refX);
printStatItem(stat,EIGEN_PARDISO_LLT, best_time_id, best_time_val);
call_directsolver(solver,EIGEN_PARDISO_LLT, A, b, refX,statFile);
}
#endif
@ -462,21 +441,16 @@ int SelectSolvers(const SparseMatrix<Scalar>&A, unsigned int sym, Matrix<Scalar,
{
cout << "\nSolving with CG ... \n";
ConjugateGradient<SpMat, Lower> solver;
stat[EIGEN_CG] = call_itersolver(solver, A, b, refX);
stat[EIGEN_CG].isIterative = 1;
printStatItem(stat,EIGEN_CG, best_time_id, best_time_val);
call_itersolver(solver,EIGEN_CG, A, b, refX,statFile);
}
//CG+IdentityPreconditioner
// {
// cout << "\nSolving with CG and IdentityPreconditioner ... \n";
// ConjugateGradient<SpMat, Lower, IdentityPreconditioner> solver;
// stat[EIGEN_CG_PRECOND] = call_itersolver(solver, A, b, refX);
// stat[EIGEN_CG_PRECOND].isIterative = 1;
// printStatItem(stat,EIGEN_CG_PRECOND, best_time_id, best_time_val);
// call_itersolver(solver,EIGEN_CG_PRECOND, A, b, refX,statFile);
// printStatItem(stat, best_time_id, best_time_val);
// }
} // End SPD matrices
return best_time_id;
}
/* Browse all the matrices available in the specified folder
@ -490,30 +464,49 @@ void Browse_Matrices(const string folder, bool statFileExists, std::string& stat
MaximumIters = maxiters; // Maximum number of iterations, global variable
RelErr = tol; //Relative residual error as stopping criterion for iterative solvers
MatrixMarketIterator<Scalar> it(folder);
Stats stat[EIGEN_ALL_SOLVERS];
for ( ; it; ++it)
{
for (int i = 0; i < EIGEN_ALL_SOLVERS; i++)
{
//print the infos for this linear system
if(statFileExists)
{
stat[i].isavail = 0;
stat[i].isIterative = 0;
std::ofstream statbuf(statFile.c_str(), std::ios::app);
statbuf << "<LINEARSYSTEM> \n";
statbuf << " <MATRIX> \n";
statbuf << " <NAME> " << it.matname() << " </NAME>\n";
statbuf << " <SIZE> " << it.matrix().rows() << " </SIZE>\n";
statbuf << " <ENTRIES> " << it.matrix().nonZeros() << "</ENTRIES>\n";
if (it.sym()!=NonSymmetric)
{
statbuf << " <SYMMETRY> Symmetric </SYMMETRY>\n" ;
if (it.sym() == SPD)
statbuf << " <POSDEF> YES </POSDEF>\n";
else
statbuf << " <POSDEF> NO </POSDEF>\n";
}
else
{
statbuf << " <SYMMETRY> NonSymmetric </SYMMETRY>\n" ;
statbuf << " <POSDEF> NO </POSDEF>\n";
}
statbuf << " </MATRIX> \n";
statbuf.close();
}
int best_time_id;
cout<< "\n\n===================================================== \n";
cout<< " ====== SOLVING WITH MATRIX " << it.matname() << " ====\n";
cout<< " =================================================== \n\n";
Matrix<Scalar, Dynamic, 1> refX;
if(it.hasrefX()) refX = it.refX();
best_time_id = SelectSolvers<Scalar>(it.matrix(), it.sym(), it.rhs(), refX, &stat[0]);
// Call all suitable solvers for this linear system
SelectSolvers<Scalar>(it.matrix(), it.sym(), it.rhs(), refX, statFile);
if(statFileExists)
{
string statline;
printHtmlStatLine(&stat[0], best_time_id, statline);
std::ofstream statbuf(statFile.c_str(), std::ios::app);
statbuf << "<TR><TH rowspan=\"4\">" << it.matname() << " <TD rowspan=\"4\"> "
<< it.matrix().rows() << " <TD rowspan=\"4\"> " << it.matrix().nonZeros()<< " "<< statline ;
statbuf << " <BEST_SOLVER ID=\""<< best_time_id
<< "\"></BEST_SOLVER>\n";
statbuf << " </LINEARSYSTEM> \n";
statbuf.close();
}
}

View File

@ -184,9 +184,20 @@ class MatrixMarketIterator
// if (S_ISDIR(st_buf.st_mode)) continue;
// Determine from the header if it is a matrix or a right hand side
bool isvector,iscomplex;
bool isvector,iscomplex=false;
if(!getMarketHeader(curfile,m_sym,iscomplex,isvector)) continue;
if(isvector) continue;
if (!iscomplex)
{
if(internal::is_same<Scalar, std::complex<float> >::value || internal::is_same<Scalar, std::complex<double> >::value)
continue;
}
if (iscomplex)
{
if(internal::is_same<Scalar, float>::value || internal::is_same<Scalar, double>::value)
continue;
}
// Get the matrix name
std::string filename = m_curs_id->d_name;