44 #include "Galeri_Maps.h" 45 #include "Galeri_CrsMatrices.h" 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");
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" ) {
91 readb, readxexact,
false,
true,
true ) );
94 if ( LastFiveBytes ==
".triU" ) {
100 if ( LastFiveBytes ==
".triS" ) {
103 readb, readxexact) );
106 if ( LastFourBytes ==
".mtx" ) {
108 readA, readx, readb, readxexact) );
109 FILE* in_file = fopen(
filename,
"r");
110 assert (in_file !=
NULL) ;
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() ) ;
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);
222 GaleriList.
set(
"nz", 4 *
Comm.NumProc());
223 GaleriList.
set(
"mx", 1);
224 GaleriList.
set(
"my", 1);
225 GaleriList.
set(
"mz",
Comm.NumProc());
237 bool distribute = false ;
257 List.
set(
"PrintTiming",
true);
258 List.
set(
"PrintStatus",
true);
259 List.
set(
"MaxProcs",
Comm.NumProc());
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]))
302 Solver->SetParameters(List);
303 Solver->SetUseTranspose(
true) ;
310 std::cout << std::endl
311 <<
"Solver " << SolverType[i]
312 <<
", verbose = " <<
verbose << std::endl ;
316 Time.ResetStartTime();
319 std::cout << std::endl
320 <<
"Solver " << SolverType[i]
321 <<
", symbolic factorization time = " 322 <<
Time.ElapsedTime() << std::endl;
327 std::cout <<
"Solver " << SolverType[i]
328 <<
", numeric factorization time = " 329 <<
Time.ElapsedTime() << std::endl;
334 std::cout <<
"Solver " << SolverType[i]
336 <<
Time.ElapsedTime() << std::endl;
342 double d = 0.0, d_tot = 0.0;
344 d += (LHS[0][j] - 1.0) * (LHS[0][j] - 1.0);
346 Comm.SumAll(&d,&d_tot,1);
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);
const Epetra_Map & RowMap() const
const Epetra_BlockMap & Map() const
const Epetra_CrsGraph & Graph() const
int Multiply(bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
bool SameAs(const Epetra_BlockMap &Map) const
cholmod_sparse *CHOLMOD() transpose(cholmod_sparse *A, int values, cholmod_common *Common)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
int PutScalar(double ScalarConstant)
int CrsMatrixTranspose(Epetra_CrsMatrix *In, Epetra_CrsMatrix *Out)
int FillComplete(bool OptimizeDataStorage=true)
int Export(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
#define AMESOS_CHK_ERR(a)
int main(int argc, char *argv[])
int NumMyElements() const
Factory for binding a third party direct solver to an Epetra_LinearProblem.
#define EPETRA_CHK_ERR(xxx)
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 NumGlobalElements() 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)
Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators...
bool Query(const char *ClassType)
Queries whether a given interface is available or not.