EpetraExt Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/hypre/cxx_main.cpp
Go to the documentation of this file.
1 //@HEADER
2 // ***********************************************************************
3 //
4 // EpetraExt: Epetra Extended - Linear Algebra Services Package
5 // Copyright (2011) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 
42 #include "Epetra_Map.h"
43 #include "Epetra_Time.h"
44 #include "Epetra_CrsMatrix.h"
46 #include "HYPRE_IJ_mv.h"
47 #include "hypre_Helpers.hpp"
48 #include "Epetra_Vector.h"
49 #include "Epetra_Flops.h"
50 #ifdef EPETRA_MPI
51 #include "Epetra_MpiComm.h"
52 #include "mpi.h"
53 #else
54 #include "Epetra_SerialComm.h"
55 #endif
56 #include "../epetra_test_err.h"
57 #include "Epetra_Version.h"
58 #include <vector>
59 #include <algorithm>
60 #include <string>
61 #include <iostream>
62 #include <fstream>
63 
64 // prototypes
65 
66 int checkValues( double x, double y, string message = "", bool verbose = false) {
67  if (fabs((x-y)/x) > 0.01) {
68  return(1);
69  if (verbose) cout << "********** " << message << " check failed.********** " << endl;
70  }
71  else {
72  if (verbose) cout << message << " check OK." << endl;
73  return(0);
74  }
75 }
76 
77 int checkMultiVectors( Epetra_MultiVector & X, Epetra_MultiVector & Y, string message = "", bool verbose = false) {
78  int numVectors = X.NumVectors();
79  int length = Y.MyLength();
80  int badvalue = 0;
81  int globalbadvalue = 0;
82  for (int j=0; j<numVectors; j++)
83  for (int i=0; i< length; i++)
84  if (checkValues(X[j][i], Y[j][i])==1) badvalue = 1;
85  X.Map().Comm().MaxAll(&badvalue, &globalbadvalue, 1);
86 
87  if (verbose) {
88  if (globalbadvalue==0) cout << message << " check OK." << endl;
89  else cout << "********* " << message << " check failed.********** " << endl;
90  }
91  return(globalbadvalue);
92 }
93 
94 int check(Epetra_RowMatrix & A, Epetra_RowMatrix & B, bool verbose);
95 
97  Epetra_Vector & q, Epetra_Vector & z, Epetra_Vector & resid, bool verbose);
98 
99 int power_method(bool TransA, Epetra_RowMatrix& A,
100  Epetra_Vector& q,
101  Epetra_Vector& z0,
102  Epetra_Vector& resid,
103  double * lambda, int niters, double tolerance,
104  bool verbose);
105 
106 int main(int argc, char *argv[])
107 {
108  int ierr = 0, i, forierr = 0;
109 #ifdef EPETRA_MPI
110 
111  // Initialize MPI
112 
113  MPI_Init(&argc,&argv);
114  int rank; // My process ID
115 
116  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
117  Epetra_MpiComm Comm( MPI_COMM_WORLD );
118 
119 #else
120 
121  int rank = 0;
122  Epetra_SerialComm Comm;
123 
124 #endif
125 
126  bool verbose = false;
127 
128  // Check if we should print results to standard out
129  if (argc>1) if (argv[1][0]=='-' && argv[1][1]=='v') verbose = true;
130 
131  int verbose_int = verbose ? 1 : 0;
132  Comm.Broadcast(&verbose_int, 1, 0);
133  verbose = verbose_int==1 ? true : false;
134 
135 
136  // char tmp;
137  // if (rank==0) cout << "Press any key to continue..."<< endl;
138  // if (rank==0) cin >> tmp;
139  // Comm.Barrier();
140 
141  Comm.SetTracebackMode(0); // This should shut down any error traceback reporting
142  int MyPID = Comm.MyPID();
143  int NumProc = Comm.NumProc();
144 
145  if(verbose && MyPID==0)
146  cout << Epetra_Version() << endl << endl;
147 
148  if (verbose) cout << "Processor "<<MyPID<<" of "<< NumProc
149  << " is alive."<<endl;
150 
151  // Redefine verbose to only print on PE 0
152  if(verbose && rank!=0)
153  verbose = false;
154 
155  int NumMyEquations = 10000;
156  int NumGlobalEquations = (NumMyEquations * NumProc) + EPETRA_MIN(NumProc,3);
157  if(MyPID < 3)
158  NumMyEquations++;
159 
160  // Construct a Map that puts approximately the same Number of equations on each processor
161 
162  Epetra_Map Map(NumGlobalEquations, NumMyEquations, 0, Comm);
163 
164  // Get update list and number of local equations from newly created Map
165  vector<int> MyGlobalElements(Map.NumMyElements());
166  Map.MyGlobalElements(&MyGlobalElements[0]);
167 
168  // Create an integer vector NumNz that is used to build the Petra Matrix.
169  // NumNz[i] is the Number of OFF-DIAGONAL term for the ith global equation on this processor
170 
171  vector<int> NumNz(NumMyEquations);
172 
173  // We are building a tridiagonal matrix where each row has (-1 2 -1)
174  // So we need 2 off-diagonal terms (except for the first and last equation)
175 
176  for(i = 0; i < NumMyEquations; i++)
177  if((MyGlobalElements[i] == 0) || (MyGlobalElements[i] == NumGlobalEquations - 1))
178  NumNz[i] = 1;
179  else
180  NumNz[i] = 2;
181 
182  // Create a Epetra_Matrix
183 
184  Epetra_CrsMatrix A(Copy, Map, &NumNz[0]);
187 
188  // Add rows one-at-a-time
189  // Need some vectors to help
190  // Off diagonal Values will always be -1
191 
192 
193  vector<double> Values(2);
194  Values[0] = -1.0;
195  Values[1] = -1.0;
196  vector<int> Indices(2);
197  double two = 2.0;
198  int NumEntries;
199 
200  forierr = 0;
201  for(i = 0; i < NumMyEquations; i++) {
202  if(MyGlobalElements[i] == 0) {
203  Indices[0] = 1;
204  NumEntries = 1;
205  }
206  else if (MyGlobalElements[i] == NumGlobalEquations-1) {
207  Indices[0] = NumGlobalEquations-2;
208  NumEntries = 1;
209  }
210  else {
211  Indices[0] = MyGlobalElements[i]-1;
212  Indices[1] = MyGlobalElements[i]+1;
213  NumEntries = 2;
214  }
215  forierr += !(A.InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[0], &Indices[0])==0);
216  forierr += !(A.InsertGlobalValues(MyGlobalElements[i], 1, &two, &MyGlobalElements[i])>0); // Put in the diagonal entry
217  }
218  EPETRA_TEST_ERR(forierr,ierr);
219 
220  // Finish up
221  A.FillComplete();
222  A.OptimizeStorage();
223 
224  HYPRE_IJMatrix Matrix;
225  int ilower = Map.MinMyGID();
226  int iupper = Map.MaxMyGID();
227 
228  //printf("Proc[%d], ilower = %d, iupper = %d.\n", MyPID, ilower, iupper);
229  HYPRE_IJMatrixCreate(MPI_COMM_WORLD, ilower, iupper, ilower, iupper, &Matrix);
230  HYPRE_IJMatrixSetObjectType(Matrix, HYPRE_PARCSR);
231  HYPRE_IJMatrixInitialize(Matrix);
232 
233  for(i = 0; i < A.NumMyRows(); i++){
234  int numElements;
235  A.NumMyRowEntries(i, numElements);
236  vector<int> my_indices; my_indices.resize(numElements);
237  vector<double> my_values; my_values.resize(numElements);
238  int numEntries;
239  A.ExtractMyRowCopy(i, numElements, numEntries, &my_values[0], &my_indices[0]);
240  for(int j = 0; j < numEntries; j++) {
241  my_indices[j] = A.GCID(my_indices[j]);
242  }
243  int GlobalRow[1];
244  GlobalRow[0] = A.GRID(i);
245  HYPRE_IJMatrixSetValues(Matrix, 1, &numEntries, GlobalRow, &my_indices[0], &my_values[0]);
246  }
247  HYPRE_IJMatrixAssemble(Matrix);
248 
249  EpetraExt_HypreIJMatrix JadA(Matrix);
250 
251  JadA.SetMaps(JadA.RowMatrixRowMap(), A.RowMatrixColMap());
252  // Create vectors for Power method
253 
254  Epetra_Vector q(Map);
255  Epetra_Vector z(Map); z.Random();
256  Epetra_Vector resid(Map);
257 
258  Epetra_Flops flopcounter;
259  A.SetFlopCounter(flopcounter);
260  q.SetFlopCounter(A);
261  z.SetFlopCounter(A);
262  resid.SetFlopCounter(A);
263  JadA.SetFlopCounter(A);
264 
265  if (verbose) cout << "=======================================" << endl
266  << "Testing Jad using CrsMatrix as input..." << endl
267  << "=======================================" << endl;
268 
269  A.ResetFlops();
270  powerMethodTests(A, JadA, Map, q, z, resid, verbose);
271 
272 #ifdef EPETRA_MPI
273  MPI_Finalize() ;
274 #endif
275 
276 return ierr ;
277 }
278 
280  Epetra_Vector & q, Epetra_Vector & z, Epetra_Vector & resid, bool verbose) {
281 
282  // variable needed for iteration
283  double lambda = 0.0;
284  // int niters = 10000;
285  int niters = 300;
286  double tolerance = 1.0e-2;
287  int ierr = 0;
288 
290 
291  // Iterate
292 
293  Epetra_Time timer(Map.Comm());
294 
295  double startTime = timer.ElapsedTime();
296  EPETRA_TEST_ERR(power_method(false, A, q, z, resid, &lambda, niters, tolerance, verbose),ierr);
297  double elapsed_time = timer.ElapsedTime() - startTime;
298  double total_flops = q.Flops();
299  double MFLOPs = total_flops/elapsed_time/1000000.0;
300  double lambdaref = lambda;
301  double flopsref = total_flops;
302 
303  if (verbose)
304  cout << "\n\nTotal MFLOPs for reference first solve = " << MFLOPs << endl
305  << "Total FLOPS = " <<total_flops <<endl<<endl;
306 
307  lambda = 0.0;
308  startTime = timer.ElapsedTime();
309 
310  EPETRA_TEST_ERR(power_method(false, JadA, q, z, resid, &lambda, niters, tolerance, verbose),ierr);
311  elapsed_time = timer.ElapsedTime() - startTime;
312  total_flops = q.Flops();
313  MFLOPs = total_flops/elapsed_time/1000000.0;
314 
315  if (verbose)
316  cout << "\n\nTotal MFLOPs for candidate first solve = " << MFLOPs << endl
317  << "Total FLOPS = " <<total_flops <<endl<<endl;
318 
319  EPETRA_TEST_ERR(checkValues(lambda,lambdaref," No-transpose Power Method result", verbose),ierr);
320  //EPETRA_TEST_ERR(checkValues(total_flops,flopsref," No-transpose Power Method flop count", verbose),ierr);
321 
323 
324  // Solve transpose problem
325 
326  if (verbose) cout << "\n\nUsing transpose of matrix and solving again (should give same result).\n\n"
327  << endl;
328  // Iterate
329  lambda = 0.0;
330  startTime = timer.ElapsedTime();
331  EPETRA_TEST_ERR(power_method(true, A, q, z, resid, &lambda, niters, tolerance, verbose),ierr);
332  elapsed_time = timer.ElapsedTime() - startTime;
333  total_flops = q.Flops();
334  MFLOPs = total_flops/elapsed_time/1000000.0;
335  lambdaref = lambda;
336  flopsref = total_flops;
337 
338  if (verbose)
339  cout << "\n\nTotal MFLOPs for reference transpose solve = " << MFLOPs << endl
340  << "Total FLOPS = " <<total_flops <<endl<<endl;
341 
342  lambda = 0.0;
343  startTime = timer.ElapsedTime();
344  EPETRA_TEST_ERR(power_method(true, JadA, q, z, resid, &lambda, niters, tolerance, verbose),ierr);
345  elapsed_time = timer.ElapsedTime() - startTime;
346  total_flops = q.Flops();
347  MFLOPs = total_flops/elapsed_time/1000000.0;
348 
349  if (verbose)
350  cout << "\n\nTotal MFLOPs for candidate transpose solve = " << MFLOPs << endl
351  << "Total FLOPS = " <<total_flops <<endl<<endl;
352 
353  EPETRA_TEST_ERR(checkValues(lambda,lambdaref,"Transpose Power Method result", verbose),ierr);
354  //EPETRA_TEST_ERR(checkValues(total_flops,flopsref,"Transpose Power Method flop count", verbose),ierr);
355 
356  EPETRA_TEST_ERR(check(A, JadA, verbose),ierr);
357 
358  return(0);
359 }
361  Epetra_Vector& resid, double* lambda, int niters, double tolerance, bool verbose)
362 {
363 
364  // Fill z with random Numbers
365  Epetra_Vector z(z0);
366 
367  // variable needed for iteration
368  double normz, residual;
369 
370  int ierr = 1;
371 
372  for(int iter = 0; iter < niters; iter++) {
373  z.Norm2(&normz); // Compute 2-norm of z
374  q.Scale(1.0/normz, z);
375  A.Multiply(TransA, q, z); // Compute z = A*q // SEGFAULT HAPPENS HERE
376  q.Dot(z, lambda); // Approximate maximum eigenvaluE
377  if(iter%100==0 || iter+1==niters) {
378  resid.Update(1.0, z, -(*lambda), q, 0.0); // Compute A*q - lambda*q
379  resid.Norm2(&residual);
380  if(verbose) cout << "Iter = " << iter << " Lambda = " << *lambda
381  << " Residual of A*q - lambda*q = " << residual << endl;
382  }
383  if(residual < tolerance) {
384  ierr = 0;
385  break;
386  }
387  }
388  return(ierr);
389 }
390 
391 int check(Epetra_RowMatrix& A, Epetra_RowMatrix & B, bool verbose) {
392 
393  int ierr = 0;
394  EPETRA_TEST_ERR(!A.Comm().NumProc()==B.Comm().NumProc(),ierr);
395  EPETRA_TEST_ERR(!A.Comm().MyPID()==B.Comm().MyPID(),ierr);
396  EPETRA_TEST_ERR(!A.Filled()==B.Filled(),ierr);
397  EPETRA_TEST_ERR(!A.HasNormInf()==B.HasNormInf(),ierr);
399  EPETRA_TEST_ERR(!A.Map().SameAs(B.Map()),ierr);
405  EPETRA_TEST_ERR(!A.NumMyCols()==B.NumMyCols(),ierr);
408  for (int i=0; i<A.NumMyRows(); i++) {
409  int nA, nB;
410  A.NumMyRowEntries(i,nA); B.NumMyRowEntries(i,nB);
411  EPETRA_TEST_ERR(!nA==nB,ierr);
412  }
413  EPETRA_TEST_ERR(!A.NumMyRows()==B.NumMyRows(),ierr);
416 
420  EPETRA_TEST_ERR(!A.UseTranspose()==B.UseTranspose(),ierr);
421 
422  int NumVectors = 5;
423  { // No transpose case
424  Epetra_MultiVector X(A.OperatorDomainMap(), NumVectors);
425  Epetra_MultiVector YA1(A.OperatorRangeMap(), NumVectors);
426  Epetra_MultiVector YA2(YA1);
427  Epetra_MultiVector YB1(YA1);
428  Epetra_MultiVector YB2(YA1);
429  X.Random();
430 
431  bool transA = false;
432  A.SetUseTranspose(transA);
433  B.SetUseTranspose(transA);
434  A.Apply(X,YA1);
435  A.Multiply(transA, X, YA2);
436  EPETRA_TEST_ERR(checkMultiVectors(YA1,YA2,"A Multiply and A Apply", verbose),ierr);
437  B.Apply(X,YB1);
438  EPETRA_TEST_ERR(checkMultiVectors(YA1,YB1,"A Multiply and B Multiply", verbose),ierr);
439  B.Multiply(transA, X, YB2);
440  EPETRA_TEST_ERR(checkMultiVectors(YA1,YB2,"A Multiply and B Apply", verbose), ierr);
441 
442  }
443  {// transpose case
444  Epetra_MultiVector X(A.OperatorRangeMap(), NumVectors);
445  Epetra_MultiVector YA1(A.OperatorDomainMap(), NumVectors);
446  Epetra_MultiVector YA2(YA1);
447  Epetra_MultiVector YB1(YA1);
448  Epetra_MultiVector YB2(YA1);
449  X.Random();
450 
451  bool transA = true;
452  A.SetUseTranspose(transA);
453  B.SetUseTranspose(transA);
454  A.Apply(X,YA1);
455  A.Multiply(transA, X, YA2);
456  EPETRA_TEST_ERR(checkMultiVectors(YA1,YA2, "A Multiply and A Apply (transpose)", verbose),ierr);
457  B.Apply(X,YB1);
458  EPETRA_TEST_ERR(checkMultiVectors(YA1,YB1, "A Multiply and B Multiply (transpose)", verbose),ierr);
459  B.Multiply(transA, X,YB2);
460  EPETRA_TEST_ERR(checkMultiVectors(YA1,YB2, "A Multiply and B Apply (transpose)", verbose),ierr);
461 
462  }
463 
464  Epetra_Vector diagA(A.RowMatrixRowMap());
465  EPETRA_TEST_ERR(A.ExtractDiagonalCopy(diagA),ierr);
466  Epetra_Vector diagB(B.RowMatrixRowMap());
467  EPETRA_TEST_ERR(B.ExtractDiagonalCopy(diagB),ierr);
468  EPETRA_TEST_ERR(checkMultiVectors(diagA,diagB, "ExtractDiagonalCopy", verbose),ierr);
469 
470  Epetra_Vector rowA(A.RowMatrixRowMap());
471  EPETRA_TEST_ERR(A.InvRowSums(rowA),ierr);
472  Epetra_Vector rowB(B.RowMatrixRowMap());
473  EPETRA_TEST_ERR(B.InvRowSums(rowB),ierr)
474  EPETRA_TEST_ERR(checkMultiVectors(rowA,rowB, "InvRowSums", verbose),ierr);
475 
476  Epetra_Vector colA(A.RowMatrixColMap());
477  EPETRA_TEST_ERR(A.InvColSums(colA),ierr);
478  Epetra_Vector colB(B.RowMatrixColMap());
479  EPETRA_TEST_ERR(B.InvColSums(colB),ierr);
480  EPETRA_TEST_ERR(checkMultiVectors(colA,colB, "InvColSums", verbose),ierr);
481 
482  EPETRA_TEST_ERR(checkValues(A.NormInf(), B.NormInf(), "NormInf before scaling", verbose), ierr);
483  EPETRA_TEST_ERR(checkValues(A.NormOne(), B.NormOne(), "NormOne before scaling", verbose),ierr);
484 
485  //EPETRA_TEST_ERR(A.RightScale(colA),ierr);
486  //EPETRA_TEST_ERR(B.RightScale(colB),ierr);
487 
488 
489  EPETRA_TEST_ERR(A.LeftScale(rowA),ierr);
490  EPETRA_TEST_ERR(B.LeftScale(rowB),ierr);
491 
492 
493  EPETRA_TEST_ERR(checkValues(A.NormInf(), B.NormInf(), "NormInf after scaling", verbose), ierr);
494  EPETRA_TEST_ERR(checkValues(A.NormOne(), B.NormOne(), "NormOne after scaling", verbose),ierr);
495 
496  vector<double> valuesA(A.MaxNumEntries());
497  vector<int> indicesA(A.MaxNumEntries());
498  vector<double> valuesB(B.MaxNumEntries());
499  vector<int> indicesB(B.MaxNumEntries());
500  return(0);
501  for (int i=0; i<A.NumMyRows(); i++) {
502  int nA, nB;
503  EPETRA_TEST_ERR(A.ExtractMyRowCopy(i, A.MaxNumEntries(), nA, &valuesA[0], &indicesA[0]),ierr);
504  EPETRA_TEST_ERR(B.ExtractMyRowCopy(i, B.MaxNumEntries(), nB, &valuesB[0], &indicesB[0]),ierr);
505  EPETRA_TEST_ERR(!nA==nB,ierr);
506  for (int j=0; j<nA; j++) {
507  double curVal = valuesA[j];
508  int curIndex = indicesA[j];
509  bool notfound = true;
510  int jj = 0;
511  while (notfound && jj< nB) {
512  if (!checkValues(curVal, valuesB[jj])) notfound = false;
513  jj++;
514  }
515  EPETRA_TEST_ERR(notfound, ierr);
516  vector<int>::iterator p = find(indicesB.begin(),indicesB.end(),curIndex); // find curIndex in indicesB
517  EPETRA_TEST_ERR(p==indicesB.end(), ierr);
518  }
519 
520  }
521  if (verbose) cout << "RowMatrix Methods check OK" << endl;
522 
523  return (ierr);
524 }
int NumMyRowEntries(int MyRow, int &NumEntries) const
void SetFlopCounter(const Epetra_Flops &FlopCounter_in)
virtual const Epetra_Map & RowMatrixRowMap() const
void SetMaps(const Epetra_Map &RowMap, const Epetra_Map &ColMap)
virtual const Epetra_Map & RowMatrixRowMap() const =0
virtual int SetUseTranspose(bool UseTranspose)=0
virtual double NormOne() const =0
bool SameAs(const Epetra_BlockMap &Map) const
#define EPETRA_TEST_ERR(a, b)
int MyGlobalElements(int *MyGlobalElementList) const
double ElapsedTime(void) const
virtual int InsertGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
virtual const Epetra_Map & OperatorDomainMap() const =0
virtual int ExtractDiagonalCopy(Epetra_Vector &Diagonal) const =0
#define EPETRA_MIN(x, y)
virtual int LeftScale(const Epetra_Vector &x)=0
int MyPID() const
std::string Epetra_Version()
bool IndicesAreLocal() const
virtual int InvRowSums(Epetra_Vector &x) const =0
virtual int NumGlobalNonzeros() const =0
virtual int MyPID() const =0
int FillComplete(bool OptimizeDataStorage=true)
virtual bool LowerTriangular() const =0
virtual int NumMyCols() const =0
int main(int argc, char **argv)
Definition: HDF5_IO.cpp:67
virtual int MaxAll(double *PartialMaxs, double *GlobalMaxs, int Count) const =0
int NumMyElements() const
virtual int NumMyRowEntries(int MyRow, int &NumEntries) const =0
virtual int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
virtual int MaxNumEntries() const =0
virtual const Epetra_Map & OperatorRangeMap() const =0
virtual const Epetra_Comm & Comm() const =0
virtual bool UseTranspose() const =0
int powerMethodTests(Epetra_RowMatrix &A, Epetra_RowMatrix &JadA, Epetra_Map &Map, Epetra_Vector &q, Epetra_Vector &z, Epetra_Vector &resid, bool verbose)
virtual const Epetra_BlockMap & Map() const =0
int NumMyRows() const
virtual double NormInf() const =0
virtual int NumMyRows() const =0
virtual bool Filled() const =0
virtual int NumGlobalDiagonals() const =0
virtual int NumGlobalCols() const =0
int NumProc() const
int MaxMyGID() const
int MinMyGID() const
const Epetra_Comm & Comm() const
bool IndicesAreGlobal() const
const Epetra_Map & RowMatrixColMap() const
virtual int Multiply(bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
int checkValues(double x, double y, string message="", bool verbose=false)
virtual int NumMyDiagonals() const =0
virtual int NumProc() const =0
virtual bool HasNormInf() const =0
int GRID(int LRID_in) const
int check(Epetra_RowMatrix &A, Epetra_RowMatrix &B, bool verbose)
virtual bool UpperTriangular() const =0
virtual int ExtractMyRowCopy(int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const =0
int checkMultiVectors(Epetra_MultiVector &X, Epetra_MultiVector &Y, string message="", bool verbose=false)
virtual const Epetra_Map & RowMatrixColMap() const =0
virtual int InvColSums(Epetra_Vector &x) const =0
int ExtractMyRowCopy(int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const
int Broadcast(double *MyVals, int Count, int Root) const
virtual int NumGlobalRows() const =0
virtual int NumMyNonzeros() const =0
int power_method(bool TransA, Epetra_RowMatrix &A, Epetra_Vector &q, Epetra_Vector &z0, Epetra_Vector &resid, double *lambda, int niters, double tolerance, bool verbose)
int GCID(int LCID_in) const