Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test/FusedImportExport/cxx_main.cpp
Go to the documentation of this file.
1 //@HEADER
2 // ************************************************************************
3 //
4 // Epetra: 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 
43 #include "Epetra_Map.h"
44 #include "Epetra_Time.h"
45 #include "Epetra_CrsMatrix.h"
46 #include "Epetra_Import.h"
47 #include "Epetra_Export.h"
48 #include "Epetra_Util.h"
49 #include "Epetra_Vector.h"
50 #include "Epetra_Flops.h"
51 
52 #ifdef EPETRA_MPI
53 
54 #include "Epetra_MpiComm.h"
55 #include "mpi.h"
56 #include "../epetra_test_err.h"
57 #include "Epetra_Version.h"
58 
59 // prototypes
60 
61 int check(Epetra_CrsMatrix& A, int NumMyRows1, int NumGlobalRows1, int NumMyNonzeros1,
62  int NumGlobalNonzeros1, int * MyGlobalElements, bool verbose);
63 
64 int power_method(bool TransA, Epetra_CrsMatrix& A,
65  Epetra_Vector& q,
66  Epetra_Vector& z,
67  Epetra_Vector& resid,
68  double * lambda, int niters, double tolerance,
69  bool verbose);
70 
72 
73 
75  const Epetra_Map & Xamap = A.DomainMap();
76  const Epetra_Map & Yamap = A.RangeMap();
77  const Epetra_Map & Xbmap = B.DomainMap();
78  const Epetra_Map & Ybmap = B.RangeMap();
79 
80  Epetra_Vector Xa(Xamap), Xb(Xbmap), Ya(Yamap), Yb(Ybmap), Diff(Yamap);
81 
82  Xa.SetSeed(24601);
83  Xa.Random();
84 
85  // Handle domain map change
86  if(!Xamap.SameAs(Xbmap)) {
87  Epetra_Import Ximport(Xbmap,Xamap);
88  Xb.Import(Xa,Ximport,Insert);
89  }
90  else {
91  Xb=Xa;
92  }
93 
94  // Do the multiplies
95  A.Apply(Xa,Ya);
96  B.Apply(Xb,Yb);
97 
98  // Handle Rangemap change
99  if(!Yamap.SameAs(Ybmap)) {
100  Epetra_Import Yimport(Yamap,Ybmap);
101  Diff.Import(Yb,Yimport,Insert);
102  }
103  else {
104  Diff=Yb;
105  }
106 
107  // Check solution
108  Diff.Update(-1.0,Ya,1.0);
109  double norm;
110  Diff.Norm2(&norm);
111 
112  return norm;
113 }
114 
115 
116 // B here is the "reduced" matrix. Square matrices w/ Row=Domain=Range only.
118  const Epetra_Map & Amap = A.DomainMap();
119  Epetra_Vector Xa(Amap), Ya(Amap), Diff(Amap);
120  const Epetra_Map *Bmap = Bfullmap.NumMyElements() > 0 ? &B.DomainMap() : 0;
121  Epetra_Vector *Xb = Bmap ? new Epetra_Vector(*Bmap) : 0;
122  Epetra_Vector *Yb = Bmap ? new Epetra_Vector(*Bmap) : 0;
123 
124  Epetra_Vector Xb_alias(View,Bfullmap, Bmap ? Xb->Values(): 0);
125  Epetra_Vector Yb_alias(View,Bfullmap, Bmap ? Yb->Values(): 0);
126 
127  Epetra_Import Ximport(Bfullmap,Amap);
128 
129  // Set the input vector
130  Xa.SetSeed(24601);
131  Xa.Random();
132  Xb_alias.Import(Xa,Ximport,Insert);
133 
134  // Do the multiplies
135  A.Apply(Xa,Ya);
136  if(Bmap) B.Apply(*Xb,*Yb);
137 
138  // Check solution
139  Epetra_Import Yimport(Amap,Bfullmap);
140  Diff.Import(Yb_alias,Yimport,Insert);
141 
142 
143  Diff.Update(-1.0,Ya,1.0);
144  double norm;
145  Diff.Norm2(&norm);
146 
147  delete Xb; delete Yb;
148  return norm;
149 }
150 
151 
152 
153 int build_matrix_unfused(const Epetra_CrsMatrix & SourceMatrix, Epetra_Import & RowImporter, Epetra_CrsMatrix *&A){
154  int rv=0;
155  rv=A->Import(SourceMatrix, RowImporter, Insert);
156  if(rv) {cerr<<"build_matrix_unfused: Import failed"<<endl; return rv;}
157 
158  rv=A->FillComplete(SourceMatrix.DomainMap(), SourceMatrix.RangeMap());
159  return rv;
160 }
161 
162 int build_matrix_unfused(const Epetra_CrsMatrix & SourceMatrix, Epetra_Export & RowExporter, Epetra_CrsMatrix *&A){
163  int rv=0;
164  rv=A->Export(SourceMatrix, RowExporter, Insert);
165  if(rv) {cerr<<"build_matrix_unfused: Export failed"<<endl; return rv;}
166 
167  rv=A->FillComplete(SourceMatrix.DomainMap(), SourceMatrix.RangeMap());
168  return rv;
169 }
170 
171 
172 
173 void build_test_matrix(Epetra_MpiComm & Comm, int test_number, Epetra_CrsMatrix *&A){
174  int NumProc = Comm.NumProc();
175  int MyPID = Comm.MyPID();
176 
177  if(test_number==1){
178  // Case 1: Tridiagonal
179  int NumMyEquations = 100;
180 
181  int NumGlobalEquations = (NumMyEquations * NumProc) + EPETRA_MIN(NumProc,3);
182  if(MyPID < 3) NumMyEquations++;
183 
184  // Construct a Map that puts approximately the same Number of equations on each processor
185  Epetra_Map Map(NumGlobalEquations, NumMyEquations, 0, Comm);
186 
187  // Get update list and number of local equations from newly created Map
188  int* MyGlobalElements = new int[Map.NumMyElements()];
189  Map.MyGlobalElements(MyGlobalElements);
190 
191  // Create an integer vector NumNz that is used to build the Petra Matrix.
192  // NumNz[i] is the Number of OFF-DIAGONAL term for the ith global equation on this processor
193 
194  int* NumNz = new int[NumMyEquations];
195 
196  // We are building a tridiagonal matrix where each row has (-1 2 -1)
197  // So we need 2 off-diagonal terms (except for the first and last equation)
198 
199  for (int i = 0; i < NumMyEquations; i++)
200  if((MyGlobalElements[i] == 0) || (MyGlobalElements[i] == NumGlobalEquations - 1))
201  NumNz[i] = 1;
202  else
203  NumNz[i] = 2;
204 
205  // Create a Epetra_Matrix
206  A=new Epetra_CrsMatrix(Copy, Map, NumNz);
207 
208  // Add rows one-at-a-time
209  // Need some vectors to help
210  // Off diagonal Values will always be -1
211 
212  double* Values = new double[2];
213  Values[0] = -1.0;
214  Values[1] = -1.0;
215  int* Indices = new int[2];
216  double two = 2.0;
217  int NumEntries;
218 
219  for (int i = 0; i < NumMyEquations; i++) {
220  if(MyGlobalElements[i] == 0) {
221  Indices[0] = 1;
222  NumEntries = 1;
223  }
224  else if (MyGlobalElements[i] == NumGlobalEquations-1) {
225  Indices[0] = NumGlobalEquations-2;
226  NumEntries = 1;
227  }
228  else {
229  Indices[0] = MyGlobalElements[i]-1;
230  Indices[1] = MyGlobalElements[i]+1;
231  NumEntries = 2;
232  }
233  A->InsertGlobalValues(MyGlobalElements[i], NumEntries, Values, Indices);
234  A->InsertGlobalValues(MyGlobalElements[i], 1, &two, MyGlobalElements+i);
235  }
236 
237  A->FillComplete();
238 
239  // Cleanup
240  delete [] MyGlobalElements;
241  delete [] NumNz;
242  delete [] Values;
243  delete [] Indices;
244 
245  }
246 }
247 
248 
249 void build_test_map(const Epetra_Map & oldMap, Epetra_Map *& newMap) {
250  int NumProc = oldMap.Comm().NumProc();
251  int MyPID = oldMap.Comm().MyPID();
252 
253  int num_global = oldMap.NumGlobalElements();
254  if(NumProc<3) {
255  // Dump everything onto -proc 0
256  int num_local = MyPID==0 ? num_global : 0;
257  newMap = new Epetra_Map(num_global,num_local,0,oldMap.Comm());
258  }
259  else {
260  // Split everything between procs 0 and 2 (leave proc 1 empty)
261  int num_local=0;
262  if(MyPID==0) num_local = num_global/2;
263  else if(MyPID==2) num_local = num_global - ((int)num_global/2);
264  newMap = new Epetra_Map(num_global,num_local,0,oldMap.Comm());
265  }
266 }
267 
268 
269 int main(int argc, char *argv[])
270 {
271  int total_err=0;
272 
273  // Initialize MPI
274 
275  MPI_Init(&argc,&argv);
276  int rank; // My process ID
277 
278  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
279  Epetra_MpiComm Comm( MPI_COMM_WORLD );
280 
281  bool verbose = false;
282 
283  // Check if we should print results to standard out
284  if (argc>1) if (argv[1][0]=='-' && argv[1][1]=='v') verbose = true;
285 
286  int verbose_int = verbose ? 1 : 0;
287  Comm.Broadcast(&verbose_int, 1, 0);
288  verbose = verbose_int==1 ? true : false;
289 
290  Comm.SetTracebackMode(0); // This should shut down any error traceback reporting
291  int MyPID = Comm.MyPID();
292  int NumProc = Comm.NumProc();
293 
294  if(verbose && MyPID==0)
295  cout << Epetra_Version() << std::endl << std::endl;
296 
297  if (verbose) cout << "Processor "<<MyPID<<" of "<< NumProc
298  << " is alive."<<endl;
299 
300  // Redefine verbose to only print on PE 0
301  if(verbose && rank!=0) verbose = false;
302 
303  // Matrix & Map pointers
304  Epetra_CrsMatrix *A, *B, *C;
305  Epetra_Map* Map1;
306  Epetra_Import* Import1;
307  Epetra_Export* Export1;
308  double diff_tol=1e-12;
309 
310 #define ENABLE_TEST_1
311 #define ENABLE_TEST_2
312 #define ENABLE_TEST_3
313 #define ENABLE_TEST_4
314 #define ENABLE_TEST_5
315 #define ENABLE_TEST_6
316 
318  // Test #1: Tridiagonal Matrix; Migrate to Proc 0
320 #ifdef ENABLE_TEST_1
321  {
322  double diff;
323  build_test_matrix(Comm,1,A);
324  int num_global = A->RowMap().NumGlobalElements();
325 
326  // New map with all on Proc1
327  if(MyPID==0) Map1=new Epetra_Map(num_global,num_global,0,Comm);
328  else Map1=new Epetra_Map(num_global,0,0,Comm);
329 
330  // Execute fused import constructor
331  Import1 = new Epetra_Import(*Map1,A->RowMap());
332  B=new Epetra_CrsMatrix(*A,*Import1,0,&A->RangeMap());
333 
334  diff=test_with_matvec(*A,*B);
335  if(diff > diff_tol){
336  if(MyPID==0) cout<<"FusedImport: Test #1 FAILED with norm diff = "<<diff<<"."<<endl;
337  total_err--;
338  }
339 
340  // Execute fused export constructor
341  delete B;
342  Export1 = new Epetra_Export(A->RowMap(),*Map1);
343  B=new Epetra_CrsMatrix(*A,*Export1,0,&A->RangeMap());
344 
345  diff=test_with_matvec(*A,*B);
346  if(diff > diff_tol){
347  if(MyPID==0) cout<<"FusedExport: Test #1 FAILED with norm diff = "<<diff<<"."<<endl;
348  total_err--;
349  }
350 
351  delete A; delete B; delete Map1; delete Import1; delete Export1;
352  }
353 #endif
354 
355 
357  // Test #2: Tridiagonal Matrix; Locally Reversed Map
359 #ifdef ENABLE_TEST_2
360  {
361  double diff;
362  build_test_matrix(Comm,1,A);
363  int num_local = A->RowMap().NumMyElements();
364 
365  std::vector<int> MyGIDS(num_local);
366  for(int i=0; i<num_local; i++)
367  MyGIDS[i] = A->RowMap().GID(num_local-i-1);
368 
369  // New map with all on Proc1
370  Map1=new Epetra_Map(-1,num_local,&MyGIDS[0],0,Comm);
371 
372  // Execute fused import constructor
373  Import1 = new Epetra_Import(*Map1,A->RowMap());
374  B=new Epetra_CrsMatrix(*A,*Import1,0,&A->RangeMap());
375 
376  diff=test_with_matvec(*A,*B);
377  if(diff > diff_tol){
378  if(MyPID==0) cout<<"FusedImport: Test #2 FAILED with norm diff = "<<diff<<"."<<endl;
379  total_err--;
380  }
381 
382  // Execute fused export constructor
383  delete B;
384  Export1 = new Epetra_Export(A->RowMap(),*Map1);
385  B=new Epetra_CrsMatrix(*A,*Export1,0,&A->RangeMap());
386 
387  diff=test_with_matvec(*A,*B);
388  if(diff > diff_tol){
389  if(MyPID==0) cout<<"FusedExport: Test #2 FAILED with norm diff = "<<diff<<"."<<endl;
390  total_err--;
391  }
392 
393  delete A; delete B; delete Map1; delete Import1; delete Export1;
394  }
395 #endif
396 
398  // Test #3: Tridiagonal Matrix; Globally Reversed Map
400 #ifdef ENABLE_TEST_3
401  {
402  double diff;
403  build_test_matrix(Comm,1,A);
404  int num_local = A->RowMap().NumMyElements();
405  int num_global = A->RowMap().NumGlobalElements();
406  int num_scansum = 0;
407 
408  Comm.ScanSum(&num_local,&num_scansum,1);
409 
410  // New Map
411  std::vector<int> MyGIDS(num_local);
412  for(int i=0; i<num_local; i++)
413  MyGIDS[i] = num_global - num_scansum + num_local - i - 1;
414  Map1=new Epetra_Map(-1,num_local,&MyGIDS[0],0,Comm);
415 
416 
417  // Execute fused import constructor
418  Import1 = new Epetra_Import(*Map1,A->RowMap());
419  B=new Epetra_CrsMatrix(*A,*Import1,0,&A->RangeMap());
420 
421  diff=test_with_matvec(*A,*B);
422  if(diff > diff_tol){
423  if(MyPID==0) cout<<"FusedImport: Test #3 FAILED with norm diff = "<<diff<<"."<<endl;
424  total_err--;
425  }
426 
427  // Execute fused export constructor
428  delete B;
429  Export1 = new Epetra_Export(A->RowMap(),*Map1);
430  B=new Epetra_CrsMatrix(*A,*Export1,0,&A->RangeMap());
431 
432  diff=test_with_matvec(*A,*B);
433  if(diff > diff_tol){
434  if(MyPID==0) cout<<"FusedExport: Test #3 FAILED with norm diff = "<<diff<<"."<<endl;
435  total_err--;
436  }
437 
438  delete A; delete B; delete Map1; delete Import1; delete Export1;
439  }
440 #endif
441 
442 
444  // Test #4: Tridiagonal Matrix; MMM style halo import
446 #ifdef ENABLE_TEST_4
447  {
448  double diff;
449  build_test_matrix(Comm,1,A);
450 
451  // Assume we always own the diagonal
452  int num_local = A->NumMyCols()-A->NumMyRows();
453  std::vector<int> MyGIDS(num_local);
454 
455  for(int i=0, idx=0; i<A->NumMyCols(); i++)
456  if(A->LRID(A->GCID(i)) == -1){
457  MyGIDS[idx] = A->GCID(i);
458  idx++;
459  }
460 
461  // New map
462  const int * MyGIDS_ptr = Epetra_Util_data_ptr(MyGIDS);
463  Map1=new Epetra_Map(-1,num_local,MyGIDS_ptr,0,Comm);
464 
465 
466  // Execute fused import constructor
467  Import1 = new Epetra_Import(*Map1,A->RowMap());
468  B=new Epetra_CrsMatrix(*A,*Import1,0,&A->RangeMap());
469 
470  // Build unfused matrix to compare
471  C=new Epetra_CrsMatrix(Copy,*Map1,0);
472  build_matrix_unfused(*A,*Import1,C);
473 
474  diff=test_with_matvec(*B,*C);
475  if(diff > diff_tol){
476  if(MyPID==0) cout<<"FusedImport: Test #4 FAILED with norm diff = "<<diff<<"."<<endl;
477  total_err--;
478  }
479 
480  // Execute fused export constructor
481  delete B;
482  Export1 = new Epetra_Export(A->RowMap(),*Map1);
483  B=new Epetra_CrsMatrix(*A,*Export1,0,&A->RangeMap());
484 
485  diff=test_with_matvec(*B,*C);
486  if(diff > diff_tol){
487  if(MyPID==0) cout<<"FusedExport: Test #4 FAILED with norm diff = "<<diff<<"."<<endl;
488  total_err--;
489  }
490 
491  delete A; delete B; delete C; delete Map1; delete Import1; delete Export1;
492  }
493 #endif
494 
495 
497  // Test 5: Tridiagonal Matrix; Migrate to Proc 0, Replace Maps
499 #ifdef ENABLE_TEST_5
500  {
501  double diff;
502  build_test_matrix(Comm,1,A);
503 
504  // New map with all on Procs 0 and 2
505  build_test_map(A->RowMap(),Map1);
506 
507  // Execute fused import constructor
508  Import1 = new Epetra_Import(*Map1,A->RowMap());
509  B=new Epetra_CrsMatrix(*A,*Import1,Map1,Map1);
510 
511  diff=test_with_matvec(*A,*B);
512  if(diff > diff_tol){
513  if(MyPID==0) cout<<"FusedImport: Test #5 FAILED with norm diff = "<<diff<<"."<<endl;
514  total_err--;
515  }
516 
517  // Execute fused export constructor
518  delete B;
519  Export1 = new Epetra_Export(A->RowMap(),*Map1);
520  B=new Epetra_CrsMatrix(*A,*Export1,Map1,Map1);
521 
522  diff=test_with_matvec(*A,*B);
523  if(diff > diff_tol){
524  if(MyPID==0) cout<<"FusedExport: Test #5 FAILED with norm diff = "<<diff<<"."<<endl;
525  total_err--;
526  }
527 
528  delete A; delete B; delete Map1; delete Import1; delete Export1;
529  }
530 #endif
531 
532 
534  // Test 6: Tridiagonal Matrix; Migrate to Proc 0, Replace Comm
536 #ifdef ENABLE_TEST_6
537  {
538  double diff;
539  build_test_matrix(Comm,1,A);
540 
541  // New map with all on Procs 0 and 2
542  build_test_map(A->RowMap(),Map1);
543 
544  // Execute fused import constructor
545  Import1 = new Epetra_Import(*Map1,A->RowMap());
546  B=new Epetra_CrsMatrix(*A,*Import1,Map1,Map1,true);
547 
548  diff=test_with_matvec_reduced_maps(*A,*B,*Map1);
549  if(diff > diff_tol){
550  if(MyPID==0) cout<<"FusedImport: Test #6 FAILED with norm diff = "<<diff<<"."<<endl;
551  total_err--;
552  }
553 
554  // Execute fused export constructor
555  delete B;
556  Export1 = new Epetra_Export(A->RowMap(),*Map1);
557  B=new Epetra_CrsMatrix(*A,*Export1,Map1,Map1,true);
558 
559  diff=test_with_matvec_reduced_maps(*A,*B,*Map1);
560  if(diff > diff_tol){
561  if(MyPID==0) cout<<"FusedExport: Test #6 FAILED with norm diff = "<<diff<<"."<<endl;
562  total_err--;
563  }
564 
565  delete A; delete B; delete Map1; delete Import1; delete Export1;
566  }
567 #endif
568 
569 
570  // Final output for OK
571  if(MyPID==0 && total_err==0)
572  cout<<"FusedImportExport: All tests PASSED."<<endl;
573 
574  // Cleanup
575  MPI_Finalize();
576 
577  return total_err ;
578 }
579 
580 
581 
582 #else
583 int main(){
584 
585  return 0;
586 }
587 #endif
int LRID(int GRID_in) const
Returns the local row index for given global row index, returns -1 if no local row for this global ro...
int NumGlobalElements() const
Number of elements across all processors.
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:119
double test_with_matvec(const Epetra_CrsMatrix &A, const Epetra_CrsMatrix &B)
const Epetra_Map & RangeMap() const
Returns the Epetra_Map object associated with the range of this matrix operator.
bool SameAs(const Epetra_BlockMap &Map) const
Returns true if this and Map are identical maps.
int MyGlobalElements(int *MyGlobalElementList) const
Puts list of global elements on this processor into the user-provided array.
virtual int InsertGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
Insert a list of elements in a given global row of the matrix.
int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.
int NumProc() const
Returns total number of processes.
static void SetTracebackMode(int TracebackModeValue)
Set the value of the Epetra_Object error traceback report mode.
T * Epetra_Util_data_ptr(std::vector< T > &vec)
Function that returns either a pointer to the first entry in the vector or, if the vector is empty...
Definition: Epetra_Util.h:422
double test_with_matvec_reduced_maps(const Epetra_CrsMatrix &A, const Epetra_CrsMatrix &B, const Epetra_Map &Bfullmap)
int power_method(Epetra_CrsMatrix &A, double &lambda, int niters, double tolerance, bool verbose)
Epetra_Export: This class builds an export object for efficient exporting of off-processor elements...
Definition: Epetra_Export.h:62
Epetra_Vector: A class for constructing and using dense vectors on a parallel computer.
#define EPETRA_MIN(x, y)
Epetra_MpiComm: The Epetra MPI Communication Class.
std::string Epetra_Version()
Epetra_Import: This class builds an import object for efficient importing of off-processor elements...
Definition: Epetra_Import.h:63
int check(Epetra_CrsGraph &A, int NumMyRows1, int NumGlobalRows1, int NumMyNonzeros1, int NumGlobalNonzeros1, int *MyGlobalElements, bool verbose)
int ScanSum(double *MyVals, double *ScanSums, int Count) const
Epetra_MpiComm Scan Sum function.
virtual int MyPID() const =0
Return my process ID.
int FillComplete(bool OptimizeDataStorage=true)
Signal that data entry is complete. Perform transformations to local index space. ...
int build_matrix_unfused(const Epetra_CrsMatrix &SourceMatrix, Epetra_Import &RowImporter, Epetra_CrsMatrix *&A)
const Epetra_Map & RowMap() const
Returns the Epetra_Map object associated with the rows of this matrix.
int NumMyElements() const
Number of elements on the calling processor.
int NumMyCols() const
Returns the number of entries in the set of column-indices that appear on this processor.
int Export(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
Exports an Epetra_DistObject using the Epetra_Import object.
Epetra_Comm: The Epetra Communication Abstract Base Class.
Definition: Epetra_Comm.h:73
int Broadcast(double *MyVals, int Count, int Root) const
Epetra_MpiComm Broadcast function.
void build_test_matrix(Epetra_MpiComm &Comm, int test_number, Epetra_CrsMatrix *&A)
int GID(int LID) const
Returns global ID of local ID, return IndexBase-1 if not found on this processor. ...
int NumMyRows() const
Returns the number of matrix rows owned by the calling processor.
const Epetra_Comm & Comm() const
Access function for Epetra_Comm communicator.
double * Values() const
Get pointer to MultiVector values.
int Update(double ScalarA, const Epetra_MultiVector &A, double ScalarThis)
Update multi-vector values with scaled values of A, this = ScalarThis*this + ScalarA*A.
int check_graph_sharing(Epetra_Comm &Comm)
virtual int NumProc() const =0
Returns total number of processes.
Epetra_CrsMatrix: A class for constructing and using real-valued double-precision sparse compressed r...
int Norm2(double *Result) const
Compute 2-norm of each vector in multi-vector.
const Epetra_Map & DomainMap() const
Returns the Epetra_Map object associated with the domain of this matrix operator. ...
int main(int argc, char *argv[])
int Import(const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
Imports an Epetra_DistObject using the Epetra_Import object.
int MyPID() const
Return my process ID.
int GCID(int LCID_in) const
Returns the global column index for give local column index, returns IndexBase-1 if we don&#39;t have thi...
void build_test_map(const Epetra_Map &oldMap, Epetra_Map *&newMap)