EpetraExt Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/inout/cxx_main.cpp
Go to the documentation of this file.
1 /*
2 //@HEADER
3 // ***********************************************************************
4 //
5 // EpetraExt: Epetra Extended - Linear Algebra Services Package
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ***********************************************************************
41 //@HEADER
42 */
43 
44 #include "Epetra_ConfigDefs.h"
45 #include "EpetraExt_Version.h"
46 #ifdef EPETRA_MPI
47 #include "mpi.h"
48 #include "Epetra_MpiComm.h"
49 #else
50 #include "Epetra_SerialComm.h"
51 #endif
52 #include "Trilinos_Util.h"
53 #include "Epetra_Comm.h"
54 #include "Epetra_Map.h"
55 #include "Epetra_Time.h"
56 #include "Epetra_BlockMap.h"
57 #include "Epetra_MultiVector.h"
58 #include "Epetra_Vector.h"
59 #include "Epetra_Export.h"
60 
61 #include "Epetra_VbrMatrix.h"
62 #include "Epetra_CrsMatrix.h"
63 #include "EpetraExt_RowMatrixOut.h"
64 #include "EpetraExt_OperatorOut.h"
66 #include "EpetraExt_VectorOut.h"
67 #include "EpetraExt_BlockMapOut.h"
68 #include "EpetraExt_BlockMapIn.h"
69 #include "EpetraExt_CrsMatrixIn.h"
71 #include "EpetraExt_VectorIn.h"
72 // include the header files again to check if guards are in place
73 #include "EpetraExt_RowMatrixOut.h"
74 #include "EpetraExt_VectorOut.h"
75 #include "EpetraExt_BlockMapOut.h"
76 #include "EpetraExt_BlockMapIn.h"
77 #include "EpetraExt_CrsMatrixIn.h"
79 #include "EpetraExt_VectorIn.h"
80 #include <fstream>
81 #include <string>
82 #include <vector>
83 #include "Poisson2dOperator.h"
84 // prototypes
85 
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;
89  return(1);
90  }
91  else {
92  if (verbose) std::cout << message << " check OK." << std::endl;
93  return(0);
94  }
95 }
96 int runTests(Epetra_Map & map, Epetra_CrsMatrix & A, Epetra_Vector & x, Epetra_Vector & b, Epetra_Vector & xexact, bool verbose);
97 int runOperatorTests(Epetra_Operator & A, bool verbose);
98 int generateHyprePrintOut(const char* filename, const Epetra_Comm &comm);
100 
101 int main(int argc, char *argv[]) {
102 
103 #ifdef EPETRA_MPI
104  MPI_Init(&argc,&argv);
105  Epetra_MpiComm comm (MPI_COMM_WORLD);
106 #else
107  Epetra_SerialComm comm;
108 #endif
109 
110  int MyPID = comm.MyPID();
111 
112  bool verbose = false;
113  bool verbose1 = false;
114  // Check if we should print results to standard out
115  if (argc > 1) {
116  if ((argv[1][0] == '-') && (argv[1][1] == 'v')) {
117  verbose1 = true;
118  if (MyPID==0) verbose = true;
119  }
120  }
121  if (verbose)
122  std::cout << EpetraExt::EpetraExt_Version() << std::endl << std::endl;
123 
124  if (verbose1) std::cout << comm << std::endl;
125 
126 
127  // Uncomment the next three lines to debug in mpi mode
128  //int tmp;
129  //if (MyPID==0) cin >> tmp;
130  //comm.Barrier();
131 
132  Epetra_Map * map;
134  Epetra_Vector * x;
135  Epetra_Vector * b;
136  Epetra_Vector * xexact;
137 
138  int nx = 20*comm.NumProc();
139  int ny = 30;
140  int npoints = 7;
141  int xoff[] = {-1, 0, 1, -1, 0, 1, 0};
142  int yoff[] = {-1, -1, -1, 0, 0, 0, 1};
143 
144 
145  int ierr = 0;
146  // Call routine to read in HB problem 0-base
147  Trilinos_Util_GenerateCrsProblem(nx, ny, npoints, xoff, yoff, comm, map, A, x, b, xexact);
148 
149  ierr += runTests(*map, *A, *x, *b, *xexact, verbose);
150 
151  delete A;
152  delete x;
153  delete b;
154  delete xexact;
155  delete map;
156 
157  // Call routine to read in HB problem 1-base
158  Trilinos_Util_GenerateCrsProblem(nx, ny, npoints, xoff, yoff, comm, map, A, x, b, xexact, 1);
159 
160  ierr += runTests(*map, *A, *x, *b, *xexact, verbose);
161 
162  delete A;
163  delete x;
164  delete b;
165  delete xexact;
166  delete map;
167 
168  // Call routine to read in HB problem -1-base
169  Trilinos_Util_GenerateCrsProblem(nx, ny, npoints, xoff, yoff, comm, map, A, x, b, xexact, -1);
170 
171  ierr += runTests(*map, *A, *x, *b, *xexact, verbose);
172 
173  delete A;
174  delete x;
175  delete b;
176  delete xexact;
177  delete map;
178 
179  int nx1 = 5;
180  int ny1 = 4;
181  Poisson2dOperator Op(nx1, ny1, comm);
182  ierr += runOperatorTests(Op, verbose);
183 
184  generateHyprePrintOut("MyMatrixFile", comm);
185 
186  EPETRA_CHK_ERR(EpetraExt::HypreFileToCrsMatrix("MyMatrixFile", comm, A));
187 
188  runHypreTest(*A);
189  delete A;
190 
191  #ifdef EPETRA_MPI
192  MPI_Finalize() ;
193 #endif
194 
195  return(ierr);
196 }
197 
199 
200  Epetra_Vector X(A.RowMap());
201  EPETRA_CHK_ERR(X.Random());
202  Epetra_Vector Y(A.RowMap());
203  EPETRA_CHK_ERR(A.Multiply(false, X, Y));
204 
205  return 0;
206 }
207 
208 int runTests(Epetra_Map & map, Epetra_CrsMatrix & A, Epetra_Vector & x, Epetra_Vector & b, Epetra_Vector & xexact, bool verbose) {
209 
210  int ierr = 0;
211 
212  // Create MultiVectors and put x, b, xexact in both columns of X, B, and Xexact, respectively.
213  Epetra_MultiVector X( map, 2, false );
214  Epetra_MultiVector B( map, 2, false );
215  Epetra_MultiVector Xexact( map, 2, false );
216 
217  for (int i=0; i<X.NumVectors(); ++i) {
218  *X(i) = x;
219  *B(i) = b;
220  *Xexact(i) = xexact;
221  }
222  double residual;
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;
232  Epetra_Vector tmp1(map);
233  Epetra_MultiVector tmp1mv(map,2,false);
234  A.Multiply(false, xexact, tmp1);
235  A.Multiply(false, Xexact, tmp1mv);
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;
250  map.Comm().Barrier();
251 
252  EPETRA_CHK_ERR(EpetraExt::BlockMapToMatrixMarketFile("Test_map.mm", map, "Official EpetraExt test map",
253  "This is the official EpetraExt test map generated by the EpetraExt regression tests"));
254 
255  EPETRA_CHK_ERR(EpetraExt::RowMatrixToMatrixMarketFile("Test_A.mm", A, "Official EpetraExt test matrix",
256  "This is the official EpetraExt test matrix generated by the EpetraExt regression tests"));
257 
258  EPETRA_CHK_ERR(EpetraExt::VectorToMatrixMarketFile("Test_x.mm", x, "Official EpetraExt test initial guess",
259  "This is the official EpetraExt test initial guess generated by the EpetraExt regression tests"));
260 
261  EPETRA_CHK_ERR(EpetraExt::MultiVectorToMatrixMarketFile("Test_mvX.mm", X, "Official EpetraExt test initial guess",
262  "This is the official EpetraExt test initial guess generated by the EpetraExt regression tests"));
263 
264  EPETRA_CHK_ERR(EpetraExt::VectorToMatrixMarketFile("Test_xexact.mm", xexact, "Official EpetraExt test exact solution",
265  "This is the official EpetraExt test exact solution generated by the EpetraExt regression tests"));
266 
267  EPETRA_CHK_ERR(EpetraExt::MultiVectorToMatrixMarketFile("Test_mvXexact.mm", Xexact, "Official EpetraExt test exact solution",
268  "This is the official EpetraExt test exact solution generated by the EpetraExt regression tests"));
269 
270  EPETRA_CHK_ERR(EpetraExt::VectorToMatrixMarketFile("Test_b.mm", b, "Official EpetraExt test right hand side",
271  "This is the official EpetraExt test right hand side generated by the EpetraExt regression tests"));
272 
273  EPETRA_CHK_ERR(EpetraExt::MultiVectorToMatrixMarketFile("Test_mvB.mm", B, "Official EpetraExt test right hand side",
274  "This is the official EpetraExt test right hand side generated by the EpetraExt regression tests"));
275 
277 
279 
280  Epetra_Map * map1;
281  Epetra_CrsMatrix * A1;
282  Epetra_CrsMatrix * A2;
283  Epetra_CrsMatrix * A3;
284  Epetra_Vector * x1;
285  Epetra_Vector * b1;
286  Epetra_Vector * xexact1;
287  Epetra_MultiVector * X1;
289  Epetra_MultiVector * Xexact1;
290 
291  EpetraExt::MatrixMarketFileToMap("Test_map.mm", map.Comm(), map1);
292 
293  if (map.SameAs(*map1)) {
294  if (verbose) std::cout << "Maps are equal. In/Out works." << std::endl;
295  }
296  else {
297  if (verbose) std::cout << "Maps are not equal. In/Out fails." << std::endl;
298  ierr += 1;
299  }
301  // If map is zero-based, then we can compare to the convenient reading versions
302  if (map1->IndexBase()==0) EPETRA_CHK_ERR(EpetraExt::MatrixMarketFileToCrsMatrix("Test_A.mm", map1->Comm(), A2));
303  if (map1->IndexBase()==0) EPETRA_CHK_ERR(EpetraExt::MatlabFileToCrsMatrix("Test_A.dat", map1->Comm(), A3));
304  EPETRA_CHK_ERR(EpetraExt::MatrixMarketFileToVector("Test_x.mm", *map1, x1));
305  EPETRA_CHK_ERR(EpetraExt::MatrixMarketFileToVector("Test_xexact.mm", *map1, xexact1));
306  EPETRA_CHK_ERR(EpetraExt::MatrixMarketFileToVector("Test_b.mm", *map1, b1));
307  EPETRA_CHK_ERR(EpetraExt::MatrixMarketFileToMultiVector("Test_mvX.mm", *map1, X1));
308  EPETRA_CHK_ERR(EpetraExt::MatrixMarketFileToMultiVector("Test_mvXexact.mm", *map1, Xexact1));
309  EPETRA_CHK_ERR(EpetraExt::MatrixMarketFileToMultiVector("Test_mvB.mm", *map1, B1));
310 
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);
314 
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);
318 
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);
322 
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);
327 
328  Epetra_Vector tmp11(*map1);
329  A1->Multiply(false, *xexact1, tmp11);
330 
331  Epetra_MultiVector tmp11mv(*map1,2,false);
332  A1->Multiply(false, *Xexact1, tmp11mv);
333 
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);
337 
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);
342 
343  if (map1->IndexBase()==0) {
344  Epetra_Vector tmp12(*map1);
345  A2->Multiply(false, *xexact1, tmp12);
346 
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);
350 
351  Epetra_Vector tmp13(*map1);
352  A3->Multiply(false, *xexact1, tmp13);
353 
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);
357  }
358  b1->Norm2(&residual); double rb1 = residual;
359  if (verbose) std::cout << "Norm of b1 (should equal norm of Ax) = " << residual << std::endl;
360  ierr += checkValues(rb1,rb,"Norm of b", verbose);
361 
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);
366 
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);
371 
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);
377 
378  if (map1->IndexBase()==0) {delete A2; delete A3;}
379  delete A1;
380  delete x1;
381  delete b1;
382  delete xexact1;
383  delete X1;
384  delete B1;
385  delete Xexact1;
386  delete map1;
387 
388  return(ierr);
389 }
390 
391 int runOperatorTests(Epetra_Operator & A, bool verbose) {
392 
393  int ierr = 0;
394 
395 
396  double residual;
397  EPETRA_CHK_ERR(EpetraExt::OperatorToMatrixMarketFile("Test_A1.mm", A, "Official EpetraExt test operator",
398  "This is the official EpetraExt test operator generated by the EpetraExt regression tests"));
400 
401  A.OperatorRangeMap().Comm().Barrier();
402  A.OperatorRangeMap().Comm().Barrier();
403  Epetra_CrsMatrix * A1;
404  Epetra_CrsMatrix * A2;
407 
408 
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);
414 
415 
416  Epetra_Vector x(A.OperatorDomainMap()); x.Random();
420  A.Apply(x,y1);
421  A1->Multiply(false, x, y2);
422  A2->Multiply(false, x, y3);
423 
424  y1.Norm2(&residual); double rAx1 = residual;
425  if (verbose) std::cout << "Norm of A*x = " << residual << std::endl;
426 
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);
430 
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);
434 
435  delete A1;
436  delete A2;
437 
438  return(ierr);
439 }
440 
441 int generateHyprePrintOut(const char *filename, const Epetra_Comm &comm){
442  int MyPID = comm.MyPID();
443  int NumProc = comm.NumProc();
444 
445  int N = 100;
446  int ilower = MyPID * N;
447  int iupper = (MyPID+1)*N-1;
448 
449  double filePID = (double)MyPID/(double)100000;
450  std::ostringstream stream;
451  // Using setprecision() puts it in the std::string
452  stream << std::setiosflags(std::ios::fixed) << std::setprecision(5) << filePID;
453  // Then just ignore the first character
454  std::string fileName(filename);
455  fileName += stream.str().substr(1,7);
456 
457  std::ofstream myfile(fileName.c_str());
458 
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;
465  }
466  }
467  myfile.close();
468  return 0;
469  } else {
470  std::cout << "\nERROR:\nCouldn't open file.\n";
471  return -1;
472  }
473 }
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.
double NormInf() const
int runHypreTest(Epetra_CrsMatrix &A)
#define EPETRA_CHK_ERR(a)
std::string EpetraExt_Version()
double NormOne() const
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)
int MyPID() const
virtual void Barrier() const =0
virtual int MyPID() const =0
int main(int argc, char **argv)
Definition: HDF5_IO.cpp:67
int IndexBase() const
const Epetra_Map & RowMap() const
virtual int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
const int N
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)
int NumProc() const
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...