44 #include "Galeri_Maps.h"
45 #include "Galeri_CrsMatrices.h"
47 using namespace Teuchos;
48 using namespace Galeri;
50 #include "Trilinos_Util.h"
51 #include "Trilinos_Util_ReadMatrixMarket2Epetra.h"
58 const bool transpose,
const bool distribute,
70 FILE *in_file = fopen( in_filename,
"r");
74 filename = &in_filename[1] ;
77 filename = in_filename ;
82 std::string FileName (filename);
84 int FN_Size = FileName.size() ;
85 std::string LastFiveBytes = FileName.substr(
EPETRA_MAX(0,FN_Size-5), FN_Size );
86 std::string LastFourBytes = FileName.substr(
EPETRA_MAX(0,FN_Size-4), FN_Size );
88 if ( LastFiveBytes ==
".TimD" ) {
90 EPETRA_CHK_ERR( Trilinos_Util_ReadTriples2Epetra( filename,
false, Comm, readMap, readA, readx,
91 readb, readxexact,
false,
true,
true ) );
94 if ( LastFiveBytes ==
".triU" ) {
96 EPETRA_CHK_ERR( Trilinos_Util_ReadTriples2Epetra( filename,
false, Comm, readMap, readA, readx,
100 if ( LastFiveBytes ==
".triS" ) {
102 EPETRA_CHK_ERR( Trilinos_Util_ReadTriples2Epetra( filename,
true, Comm, readMap, readA, readx,
103 readb, readxexact) );
106 if ( LastFourBytes ==
".mtx" ) {
107 EPETRA_CHK_ERR( Trilinos_Util_ReadMatrixMarket2Epetra( filename, Comm, readMap,
108 readA, readx, readb, readxexact) );
109 FILE* in_file = fopen( filename,
"r");
110 assert (in_file != NULL) ;
112 char buffer[BUFSIZE] ;
113 fgets( buffer, BUFSIZE, in_file ) ;
114 std::string headerline1 = buffer;
116 if ( headerline1.find(
"symmetric") < BUFSIZE ) symmetric =
true;
118 if ( headerline1.find(
"symmetric") != std::string::npos) symmetric =
true;
125 Trilinos_Util_ReadHb2Epetra( filename, Comm, readMap, readA, readx,
127 if ( LastFourBytes ==
".rsa" ) symmetric = true ;
134 if ( readb )
delete readb;
135 if ( readx )
delete readx;
136 if ( readxexact )
delete readxexact;
144 serialA = transposeA ;
151 assert( (
void *) &serialA->
Graph() ) ;
152 assert( (
void *) &serialA->
RowMap() ) ;
163 Amat->Export(*serialA, exporter,
Add);
202 int main(
int argc,
char *argv[])
205 MPI_Init(&argc, &argv);
212 double TotalResidual = 0.0;
218 #ifndef FILENAME_SPECIFIED_ON_COMMAND_LINE
220 GaleriList.
set(
"nx", 4);
221 GaleriList.
set(
"ny", 4);
223 GaleriList.
set(
"mx", 1);
224 GaleriList.
set(
"my", 1);
226 Epetra_Map* Map = CreateMap(
"Cartesian3D", Comm, GaleriList);
236 bool transpose = false ;
237 bool distribute = false ;
241 MyCreateCrsMatrix( argv[1], Comm, Map, transpose, distribute, symmetric, Matrix ) ;
257 List.
set(
"PrintTiming",
true);
258 List.
set(
"PrintStatus",
true);
261 std::vector<std::string> SolverType;
262 SolverType.push_back(
"Amesos_Paraklete");
263 SolverType.push_back(
"Amesos_Klu");
266 SolverType.push_back(
"Amesos_Lapack");
267 SolverType.push_back(
"Amesos_Umfpack");
268 SolverType.push_back(
"Amesos_Pardiso");
269 SolverType.push_back(
"Amesos_Taucs");
270 SolverType.push_back(
"Amesos_Superlu");
271 SolverType.push_back(
"Amesos_Superludist");
272 SolverType.push_back(
"Amesos_Mumps");
273 SolverType.push_back(
"Amesos_Dscpack");
274 SolverType.push_back(
"Amesos_Scalapack");
284 for (
unsigned int i = 0 ; i < SolverType.size() ; ++i)
287 if (Factory.
Query(SolverType[i]))
300 assert (Solver != 0);
310 std::cout << std::endl
311 <<
"Solver " << SolverType[i]
312 <<
", verbose = " << verbose << std::endl ;
319 std::cout << std::endl
320 <<
"Solver " << SolverType[i]
321 <<
", symbolic factorization time = "
327 std::cout <<
"Solver " << SolverType[i]
328 <<
", numeric factorization time = "
334 std::cout <<
"Solver " << SolverType[i]
342 double d = 0.0, d_tot = 0.0;
344 d += (LHS[0][j] - 1.0) * (LHS[0][j] - 1.0);
348 std::cout <<
"Solver " << SolverType[i] <<
", ||x - x_exact||_2 = "
349 << sqrt(d_tot) << std::endl;
354 TotalResidual += d_tot;
361 if (TotalResidual > 1e-9)
368 return(EXIT_SUCCESS);
int NumGlobalElements() const
virtual int Solve()=0
Solves A X = B (or AT x = B)
bool SameAs(const Epetra_BlockMap &Map) const
int Multiply(bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
double ElapsedTime(void) const
virtual int NumericFactorization()=0
Performs NumericFactorization on the matrix A.
virtual int SymbolicFactorization()=0
Performs SymbolicFactorization on the matrix A.
virtual int SetParameters(Teuchos::ParameterList &ParameterList)=0
Updates internal variables.
int CrsMatrixTranspose(Epetra_CrsMatrix *In, Epetra_CrsMatrix *Out)
int FillComplete(bool OptimizeDataStorage=true)
const Epetra_Map & RowMap() const
int NumMyElements() const
#define AMESOS_CHK_ERR(a)
virtual const Epetra_BlockMap & Map() const =0
int main(int argc, char *argv[])
Factory for binding a third party direct solver to an Epetra_LinearProblem.
#define EPETRA_CHK_ERR(xxx)
virtual int SetUseTranspose(bool UseTranspose)=0
If set true, X will be set to the solution of AT X = B (not A X = B)
int CreateCrsMatrix(const char *in_filename, const Epetra_Comm &Comm, Epetra_Map *&readMap, const bool transpose, const bool distribute, bool &symmetric, Epetra_CrsMatrix *&Matrix)
Amesos_BaseSolver * Create(const char *ClassType, const Epetra_LinearProblem &LinearProblem)
Amesos Create method.
int SumAll(double *PartialSums, double *GlobalSums, int Count) const
int MyCreateCrsMatrix(const char *in_filename, const Epetra_Comm &Comm, Epetra_Map *&readMap, const bool transpose, const bool distribute, bool &symmetric, Epetra_CrsMatrix *&Matrix)
const Epetra_CrsGraph & Graph() const
Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators...
void ResetStartTime(void)
bool Query(const char *ClassType)
Queries whether a given interface is available or not.