eigen/bench/spbench/test_sparseLU.cpp

93 lines
2.8 KiB
C++
Raw Normal View History

2012-06-15 00:45:04 +08:00
// Small bench routine for Eigen available in Eigen
// (C) Desire NUENTSA WAKAM, INRIA
#include <iostream>
#include <fstream>
#include <iomanip>
#include <unsupported/Eigen/SparseExtra>
#include <Eigen/SparseLU>
#include <bench/BenchTimer.h>
#ifdef EIGEN_METIS_SUPPORT
#include <Eigen/MetisSupport>
#endif
2012-06-15 00:45:04 +08:00
using namespace std;
using namespace Eigen;
int main(int argc, char **args)
{
2012-07-27 17:36:58 +08:00
// typedef complex<double> scalar;
typedef double scalar;
SparseMatrix<scalar, ColMajor> A;
typedef SparseMatrix<scalar, ColMajor>::Index Index;
typedef Matrix<scalar, Dynamic, Dynamic> DenseMatrix;
typedef Matrix<scalar, Dynamic, 1> DenseRhs;
Matrix<scalar, Dynamic, 1> b, x, tmp;
// SparseLU<SparseMatrix<scalar, ColMajor>, AMDOrdering<int> > solver;
2012-09-04 18:21:07 +08:00
// #ifdef EIGEN_METIS_SUPPORT
// SparseLU<SparseMatrix<scalar, ColMajor>, MetisOrdering<int> > solver;
// std::cout<< "ORDERING : METIS\n";
// #else
2012-09-25 17:48:18 +08:00
SparseLU<SparseMatrix<scalar, ColMajor>, COLAMDOrdering<int> > solver;
2012-09-04 18:21:07 +08:00
std::cout<< "ORDERING : COLAMD\n";
// #endif
2012-06-15 00:45:04 +08:00
ifstream matrix_file;
string line;
int n;
BenchTimer timer;
2012-06-15 00:45:04 +08:00
// Set parameters
/* Fill the matrix with sparse matrix stored in Matrix-Market coordinate column-oriented format */
if (argc < 2) assert(false && "please, give the matrix market file ");
loadMarket(A, args[1]);
cout << "End charging matrix " << endl;
bool iscomplex=false, isvector=false;
int sym;
getMarketHeader(args[1], sym, iscomplex, isvector);
2012-07-20 00:15:23 +08:00
// if (iscomplex) { cout<< " Not for complex matrices \n"; return -1; }
2012-06-15 00:45:04 +08:00
if (isvector) { cout << "The provided file is not a matrix file\n"; return -1;}
if (sym != 0) { // symmetric matrices, only the lower part is stored
SparseMatrix<scalar, ColMajor> temp;
2012-06-15 00:45:04 +08:00
temp = A;
A = temp.selfadjointView<Lower>();
}
n = A.cols();
/* Fill the right hand side */
if (argc > 2)
loadMarketVector(b, args[2]);
else
{
b.resize(n);
tmp.resize(n);
// tmp.setRandom();
for (int i = 0; i < n; i++) tmp(i) = i;
b = A * tmp ;
}
/* Compute the factorization */
2012-07-07 02:18:16 +08:00
// solver.isSymmetric(true);
timer.start();
// solver.compute(A);
solver.analyzePattern(A);
timer.stop();
cout << "Time to analyze " << timer.value() << std::endl;
timer.reset();
timer.start();
solver.factorize(A);
timer.stop();
cout << "Factorize Time " << timer.value() << std::endl;
timer.reset();
timer.start();
2012-08-03 19:05:45 +08:00
x = solver.solve(b);
timer.stop();
cout << "solve time " << timer.value() << std::endl;
2012-06-15 00:45:04 +08:00
/* Check the accuracy */
Matrix<scalar, Dynamic, 1> tmp2 = b - A*x;
scalar tempNorm = tmp2.norm()/b.norm();
2012-06-15 00:45:04 +08:00
cout << "Relative norm of the computed solution : " << tempNorm <<"\n";
cout << "Number of nonzeros in the factor : " << solver.nnzL() + solver.nnzU() << std::endl;
2012-06-15 00:45:04 +08:00
return 0;
}