52 #include "Trilinos_Util.h"
86 int checkValues(
double x,
double y, std::string message =
"",
bool verbose =
false) {
87 if (fabs((x-y)/x) > 0.01 && x > 1.0e-12) {
88 if (verbose) std::cout <<
"********** " << message <<
" check failed.********** " << std::endl;
92 if (verbose) std::cout << message <<
" check OK." << std::endl;
101 int main(
int argc,
char *argv[]) {
104 MPI_Init(&argc,&argv);
110 int MyPID = comm.
MyPID();
112 bool verbose =
false;
113 bool verbose1 =
false;
116 if ((argv[1][0] ==
'-') && (argv[1][1] ==
'v')) {
118 if (MyPID==0) verbose =
true;
124 if (verbose1) std::cout << comm << std::endl;
141 int xoff[] = {-1, 0, 1, -1, 0, 1, 0};
142 int yoff[] = {-1, -1, -1, 0, 0, 0, 1};
147 Trilinos_Util_GenerateCrsProblem(nx, ny, npoints, xoff, yoff, comm, map, A, x, b, xexact);
149 ierr +=
runTests(*map, *A, *x, *b, *xexact, verbose);
158 Trilinos_Util_GenerateCrsProblem(nx, ny, npoints, xoff, yoff, comm, map, A, x, b, xexact, 1);
160 ierr +=
runTests(*map, *A, *x, *b, *xexact, verbose);
169 Trilinos_Util_GenerateCrsProblem(nx, ny, npoints, xoff, yoff, comm, map, A, x, b, xexact, -1);
171 ierr +=
runTests(*map, *A, *x, *b, *xexact, verbose);
217 for (
int i=0; i<X.NumVectors(); ++i) {
223 std::vector<double> residualmv(2);
224 residual = A.
NormInf();
double rAInf = residual;
225 if (verbose) std::cout <<
"Inf Norm of A = " << residual << std::endl;
226 residual = A.
NormOne();
double rAOne = residual;
227 if (verbose) std::cout <<
"One Norm of A = " << residual << std::endl;
228 xexact.Norm2(&residual);
double rxx = residual;
229 Xexact.Norm2(&residualmv[0]); std::vector<double> rXX( residualmv );
230 if (verbose) std::cout <<
"Norm of xexact = " << residual << std::endl;
231 if (verbose) std::cout <<
"Norm of Xexact = (" << residualmv[0] <<
", " <<residualmv[1] <<
")"<< std::endl;
236 tmp1.Norm2(&residual);
double rAx = residual;
237 tmp1mv.Norm2(&residualmv[0]); std::vector<double> rAX( residualmv );
238 if (verbose) std::cout <<
"Norm of Ax = " << residual << std::endl;
239 if (verbose) std::cout <<
"Norm of AX = (" << residualmv[0] <<
", " << residualmv[1] <<
")"<< std::endl;
240 b.Norm2(&residual);
double rb = residual;
241 B.Norm2(&residualmv[0]); std::vector<double> rB( residualmv );
242 if (verbose) std::cout <<
"Norm of b (should equal norm of Ax) = " << residual << std::endl;
243 if (verbose) std::cout <<
"Norm of B (should equal norm of AX) = (" << residualmv[0] <<
", " << residualmv[1] <<
")"<< std::endl;
244 tmp1.Update(1.0, b, -1.0);
245 tmp1mv.Update(1.0, B, -1.0);
246 tmp1.Norm2(&residual);
247 tmp1mv.Norm2(&residualmv[0]);
248 if (verbose) std::cout <<
"Norm of difference between compute Ax and Ax from file = " << residual << std::endl;
249 if (verbose) std::cout <<
"Norm of difference between compute AX and AX from file = (" << residualmv[0] <<
", " << residualmv[1] <<
")"<< std::endl;
253 "This is the official EpetraExt test map generated by the EpetraExt regression tests"));
256 "This is the official EpetraExt test matrix generated by the EpetraExt regression tests"));
259 "This is the official EpetraExt test initial guess generated by the EpetraExt regression tests"));
262 "This is the official EpetraExt test initial guess generated by the EpetraExt regression tests"));
265 "This is the official EpetraExt test exact solution generated by the EpetraExt regression tests"));
268 "This is the official EpetraExt test exact solution generated by the EpetraExt regression tests"));
271 "This is the official EpetraExt test right hand side generated by the EpetraExt regression tests"));
274 "This is the official EpetraExt test right hand side generated by the EpetraExt regression tests"));
294 if (verbose) std::cout <<
"Maps are equal. In/Out works." << std::endl;
297 if (verbose) std::cout <<
"Maps are not equal. In/Out fails." << std::endl;
311 residual = A1->
NormInf();
double rA1Inf = residual;
312 if (verbose) std::cout <<
"Inf Norm of A1 = " << residual << std::endl;
313 ierr +=
checkValues(rA1Inf,rAInf,
"Inf Norm of A", verbose);
315 residual = A1->
NormOne();
double rA1One = residual;
316 if (verbose) std::cout <<
"One Norm of A1 = " << residual << std::endl;
317 ierr +=
checkValues(rA1One,rAOne,
"One Norm of A", verbose);
319 xexact1->Norm2(&residual);
double rxx1 = residual;
320 if (verbose) std::cout <<
"Norm of xexact1 = " << residual << std::endl;
321 ierr +=
checkValues(rxx1,rxx,
"Norm of xexact", verbose);
323 Xexact1->Norm2(&residualmv[0]); std::vector<double> rXX1(residualmv);
324 if (verbose) std::cout <<
"Norm of Xexact1 = (" << residualmv[0] <<
", " <<residualmv[1]<<
")"<< std::endl;
325 ierr +=
checkValues(rXX1[0],rXX[0],
"Norm of Xexact", verbose);
326 ierr +=
checkValues(rXX1[1],rXX[1],
"Norm of Xexact", verbose);
329 A1->
Multiply(
false, *xexact1, tmp11);
332 A1->
Multiply(
false, *Xexact1, tmp11mv);
334 tmp11.Norm2(&residual);
double rAx1 = residual;
335 if (verbose) std::cout <<
"Norm of A1*x1 = " << residual << std::endl;
336 ierr +=
checkValues(rAx1,rAx,
"Norm of A1*x", verbose);
338 tmp11mv.Norm2(&residualmv[0]); std::vector<double> rAX1(residualmv);
339 if (verbose) std::cout <<
"Norm of A1*X1 = (" << residualmv[0] <<
", "<<residualmv[1]<<
")"<< std::endl;
340 ierr +=
checkValues(rAX1[0],rAX[0],
"Norm of A1*X", verbose);
341 ierr +=
checkValues(rAX1[1],rAX[1],
"Norm of A1*X", verbose);
345 A2->
Multiply(
false, *xexact1, tmp12);
347 tmp12.Norm2(&residual);
double rAx2 = residual;
348 if (verbose) std::cout <<
"Norm of A2*x1 = " << residual << std::endl;
349 ierr +=
checkValues(rAx2,rAx,
"Norm of A2*x", verbose);
352 A3->
Multiply(
false, *xexact1, tmp13);
354 tmp13.Norm2(&residual);
double rAx3 = residual;
355 if (verbose) std::cout <<
"Norm of A3*x1 = " << residual << std::endl;
356 ierr +=
checkValues(rAx3,rAx,
"Norm of A3*x", verbose);
358 b1->Norm2(&residual);
double rb1 = residual;
359 if (verbose) std::cout <<
"Norm of b1 (should equal norm of Ax) = " << residual << std::endl;
362 B1->Norm2(&residualmv[0]); std::vector<double> rB1(residualmv);
363 if (verbose) std::cout <<
"Norm of B1 (should equal norm of AX) = (" << residualmv[0] <<
", "<<residualmv[1]<<
")"<< std::endl;
364 ierr +=
checkValues(rB1[0],rB[0],
"Norm of B", verbose);
365 ierr +=
checkValues(rB1[1],rB[1],
"Norm of B", verbose);
367 tmp11.Update(1.0, *b1, -1.0);
368 tmp11.Norm2(&residual);
369 if (verbose) std::cout <<
"Norm of difference between computed A1x1 and A1x1 from file = " << residual << std::endl;
370 ierr +=
checkValues(residual,0.0,
"Norm of difference between computed A1x1 and A1x1 from file", verbose);
372 tmp11mv.Update(1.0, *B1, -1.0);
373 tmp11mv.Norm2(&residualmv[0]);
374 if (verbose) std::cout <<
"Norm of difference between computed A1X1 and A1X1 from file = (" << residualmv[0] <<
", "<<residualmv[1]<<
")"<< std::endl;
375 ierr +=
checkValues(residualmv[0],0.0,
"Norm of difference between computed A1X1 and A1X1 from file", verbose);
376 ierr +=
checkValues(residualmv[1],0.0,
"Norm of difference between computed A1X1 and A1X1 from file", verbose);
378 if (map1->
IndexBase()==0) {
delete A2;
delete A3;}
398 "This is the official EpetraExt test operator generated by the EpetraExt regression tests"));
409 residual = A.
NormInf();
double rAInf = residual;
410 if (verbose) std::cout <<
"Inf Norm of Operator A = " << residual << std::endl;
411 residual = A1->
NormInf();
double rA1Inf = residual;
412 if (verbose) std::cout <<
"Inf Norm of Matrix A1 = " << residual << std::endl;
413 ierr +=
checkValues(rA1Inf,rAInf,
"Inf Norm of A", verbose);
424 y1.Norm2(&residual);
double rAx1 = residual;
425 if (verbose) std::cout <<
"Norm of A*x = " << residual << std::endl;
427 y2.Norm2(&residual);
double rAx2 = residual;
428 if (verbose) std::cout <<
"Norm of A1*x = " << residual << std::endl;
429 ierr +=
checkValues(rAx1,rAx2,
"Norm of A1*x", verbose);
431 y3.Norm2(&residual);
double rAx3 = residual;
432 if (verbose) std::cout <<
"Norm of A2*x = " << residual << std::endl;
433 ierr +=
checkValues(rAx1,rAx3,
"Norm of A2*x", verbose);
442 int MyPID = comm.
MyPID();
446 int ilower = MyPID *
N;
447 int iupper = (MyPID+1)*N-1;
449 double filePID = (double)MyPID/(
double)100000;
450 std::ostringstream stream;
452 stream << std::setiosflags(std::ios::fixed) << std::setprecision(5) << filePID;
454 std::string fileName(filename);
455 fileName += stream.str().substr(1,7);
457 std::ofstream myfile(fileName.c_str());
459 if(myfile.is_open()){
460 myfile << ilower <<
" " << iupper <<
" " << ilower <<
" " << iupper << std::endl;
461 for(
int i = ilower; i <= iupper; i++){
462 for(
int j=i-5; j <= i+5; j++){
463 if(j >= 0 && j < N*NumProc)
464 myfile << i <<
" " << j <<
" " << (double)rand()/(double)RAND_MAX << std::endl;
470 std::cout <<
"\nERROR:\nCouldn't open file.\n";
int OperatorToMatrixMarketFile(const char *filename, const Epetra_Operator &A, const char *matrixName, const char *matrixDescription, bool writeHeader)
Writes an Epetra_Operator object to a Matrix Market format file, forming the coefficients by applying...
int MatrixMarketFileToMultiVector(const char *filename, const Epetra_BlockMap &map, Epetra_MultiVector *&A)
Constructs an Epetra_MultiVector object from a Matrix Market format file.
int MultiVectorToMatlabFile(const char *filename, const Epetra_MultiVector &A)
Writes an Epetra_MultiVector object to a file that is compatible with Matlab.
bool SameAs(const Epetra_BlockMap &Map) const
int Multiply(bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
int BlockMapToMatrixMarketFile(const char *filename, const Epetra_BlockMap &map, const char *mapName, const char *mapDescription, bool writeHeader)
Writes an Epetra_BlockMap or Epetra_Map object to a Matrix Market format file.
int runHypreTest(Epetra_CrsMatrix &A)
#define EPETRA_CHK_ERR(a)
std::string EpetraExt_Version()
virtual const Epetra_Map & OperatorDomainMap() const =0
int runTests(Epetra_Map &map, Epetra_CrsMatrix &A, Epetra_Vector &x, Epetra_Vector &b, Epetra_Vector &xexact, bool verbose)
virtual void Barrier() const =0
virtual int MyPID() const =0
int main(int argc, char **argv)
const Epetra_Map & RowMap() const
virtual int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
virtual const Epetra_Map & OperatorRangeMap() const =0
int RowMatrixToMatlabFile(const char *filename, const Epetra_RowMatrix &A)
Writes an Epetra_RowMatrix object to a file that is compatible with Matlab.
int MatrixMarketFileToCrsMatrix(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&A)
int MatrixMarketFileToMap(const char *filename, const Epetra_Comm &comm, Epetra_Map *&map)
Constructs an Epetra_BlockMap object from a Matrix Market format file.
int generateHyprePrintOut(const char *filename, const Epetra_Comm &comm)
int MultiVectorToMatrixMarketFile(const char *filename, const Epetra_MultiVector &A, const char *matrixName, const char *matrixDescription, bool writeHeader)
Writes an Epetra_MultiVector object to a Matrix Market format file.
Poisson2dOperator: A sample implementation of the Epetra_Operator class.
int runOperatorTests(Epetra_Operator &A, bool verbose)
const Epetra_Comm & Comm() const
int MatlabFileToCrsMatrix(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&A)
Constructs an Epetra_CrsMatrix object from a Matlab format file, distributes rows evenly across proce...
int RowMatrixToMatrixMarketFile(const char *filename, const Epetra_RowMatrix &A, const char *matrixName, const char *matrixDescription, bool writeHeader)
Writes an Epetra_RowMatrix object to a Matrix Market format file.
int VectorToMatrixMarketFile(const char *filename, const Epetra_Vector &A, const char *matrixName, const char *matrixDescription, bool writeHeader)
Writes an Epetra_Vector object to a Matrix Market format file.
int checkValues(double x, double y, string message="", bool verbose=false)
virtual int NumProc() const =0
int OperatorToMatlabFile(const char *filename, const Epetra_Operator &A)
Writes an Epetra_Operator object to a file that is compatible with Matlab.
virtual double NormInf() const =0
int MatrixMarketFileToVector(const char *filename, const Epetra_BlockMap &map, Epetra_Vector *&A)
Constructs an Epetra_Vector object from a Matrix Market format file.
int HypreFileToCrsMatrix(const char *filename, const Epetra_Comm &comm, Epetra_CrsMatrix *&Matrix)
Constructs an Epetra_CrsMatrix object from a Hypre Matrix Print command, the row map is specified...