Epetra Package Browser (Single Doxygen Collection)  Development
test/CrsGraph/cxx_main.cpp
Go to the documentation of this file.
1 /*
2 //@HEADER
3 // ************************************************************************
4 //
5 // Epetra: 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 "Epetra_CrsGraph.h"
46 #include "Epetra_Map.h"
47 #ifdef EPETRA_MPI
48 #include "Epetra_MpiComm.h"
49 #include <mpi.h>
50 #else
51 #include "Epetra_SerialComm.h"
52 #endif
53 #include "../epetra_test_err.h"
54 #include "Epetra_Version.h"
55 
56 // Prototype
57 int check(Epetra_CrsGraph& A, int NumMyRows1, int NumGlobalRows1, int NumMyNonzeros1,
58  int NumGlobalNonzeros1, int* MyGlobalElements, bool verbose);
59 
60 int checkSharedOwnership(Epetra_Comm& Comm, bool verbose);
61 int checkCopyAndAssignment(Epetra_Comm& Comm, bool verbose);
62 
63 int main(int argc, char *argv[]) {
64  int ierr = 0;
65  int i;
66  int forierr = 0;
67  int* Indices;
68  bool debug = true;
69 
70 #ifdef EPETRA_MPI
71 
72  // Initialize MPI
73 
74  MPI_Init(&argc,&argv);
75  int rank; // My process ID
76 
77  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
78  Epetra_MpiComm Comm( MPI_COMM_WORLD );
79 
80 #else
81 
82  int rank = 0;
83  Epetra_SerialComm Comm;
84 
85 #endif
86 
87  bool verbose = false;
88 
89  // Check if we should print results to standard out
90  if(argc > 1) {
91  if(argv[1][0]=='-' && argv[1][1]=='v') {
92  verbose = true;
93  }
94  }
95 
96  //char tmp;
97  //if (rank==0) cout << "Press any key to continue..."<< endl;
98  //if (rank==0) cin >> tmp;
99  //Comm.Barrier();
100 
101  Comm.SetTracebackMode(0); // This should shut down any error traceback reporting
102  int MyPID = Comm.MyPID();
103  int NumProc = Comm.NumProc();
104 
105  if(verbose && MyPID==0)
106  cout << Epetra_Version() << endl << endl;
107 
108  if(verbose) cout << "Processor "<<MyPID<<" of "<< NumProc << " is alive." << endl;
109 
110  bool verbose1 = verbose;
111 
112  // Redefine verbose to only print on PE 0
113  if(verbose && rank != 0)
114  verbose = false;
115 
116  int NumMyEquations = 5;
117  int NumGlobalEquations = NumMyEquations*NumProc+EPETRA_MIN(NumProc,3);
118  if(MyPID < 3)
119  NumMyEquations++;
120 
121  // Construct a Map that puts approximately the same Number of equations on each processor
122 
123  Epetra_Map& Map = *new Epetra_Map(NumGlobalEquations, NumMyEquations, 0, Comm);
124 
125  // Get update list and number of local equations from newly created Map
126  int* MyGlobalElements = new int[Map.NumMyElements()];
127  Map.MyGlobalElements(MyGlobalElements);
128 
129  // Create an integer vector NumNz that is used to build the Petra Matrix.
130  // NumNz[i] is the Number of OFF-DIAGONAL term for the ith global equation on this processor
131 
132  int* NumNz = new int[NumMyEquations];
133 
134  // We are building a tridiagonal matrix where each row has (-1 2 -1)
135  // So we need 2 off-diagonal terms (except for the first and last equation)
136 
137  for(i=0; i<NumMyEquations; i++)
138  if(MyGlobalElements[i]==0 || MyGlobalElements[i] == NumGlobalEquations-1)
139  NumNz[i] = 1;
140  else
141  NumNz[i] = 2;
142 
143  // Create a Epetra_CrsGraph
144 
145  Epetra_CrsGraph& A = *new Epetra_CrsGraph(Copy, Map, NumNz);
146  EPETRA_TEST_ERR(A.IndicesAreGlobal(),ierr);
147  EPETRA_TEST_ERR(A.IndicesAreLocal(),ierr);
148 
149  // Add rows one-at-a-time
150  // Need some vectors to help
151  // Off diagonal Values will always be -1
152 
153  Indices = new int[2];
154  int NumEntries;
155 
156  forierr = 0;
157  for(i = 0; i < NumMyEquations; i++) {
158  if(MyGlobalElements[i] == 0) {
159  Indices[0] = 1;
160  NumEntries = 1;
161  }
162  else if(MyGlobalElements[i] == NumGlobalEquations-1) {
163  Indices[0] = NumGlobalEquations-2;
164  NumEntries = 1;
165  }
166  else {
167  Indices[0] = MyGlobalElements[i]-1;
168  Indices[1] = MyGlobalElements[i]+1;
169  NumEntries = 2;
170  }
171  forierr += !(A.InsertGlobalIndices(MyGlobalElements[i], NumEntries, Indices)==0);
172  forierr += !(A.InsertGlobalIndices(MyGlobalElements[i], 1, MyGlobalElements+i)>0); // Put in the diagonal entry (should cause realloc)
173  }
174  EPETRA_TEST_ERR(forierr,ierr);
175  //A.PrintGraphData(cout);
176  delete[] Indices;
177 
178  int gRID = A.GRID(0);
179  int numIndices = A.NumGlobalIndices(gRID);
180  std::vector<int> indices_vec(numIndices);
181  A.ExtractGlobalRowCopy(gRID, numIndices, numIndices, &indices_vec[0]);
182  A.RemoveGlobalIndices(gRID);
183  EPETRA_TEST_ERR(!(A.NumGlobalIndices(gRID)==0),ierr);
184  if (ierr != 0) cout << "tests FAILED" << std::endl;
185 
186  A.InsertGlobalIndices(gRID, numIndices, &indices_vec[0]);
187  EPETRA_TEST_ERR(!(A.NumGlobalIndices(gRID)==numIndices),ierr);
188  if (ierr != 0) cout << "tests FAILED" << std::endl;
189 
190  // Finish up
191  EPETRA_TEST_ERR(!(A.IndicesAreGlobal()),ierr);
192  EPETRA_TEST_ERR(!(A.FillComplete()==0),ierr);
193  EPETRA_TEST_ERR(!(A.IndicesAreLocal()),ierr);
194  EPETRA_TEST_ERR(A.StorageOptimized(),ierr);
195 
196  A.OptimizeStorage();
197 
198  EPETRA_TEST_ERR(!(A.StorageOptimized()),ierr);
199  EPETRA_TEST_ERR(A.UpperTriangular(),ierr);
200  EPETRA_TEST_ERR(A.LowerTriangular(),ierr);
201 
202  if (ierr != 0) cout << "tests FAILED" << std::endl;
203 
204  EPETRA_TEST_ERR(A.NumMyDiagonals()-NumMyEquations,ierr);
205  if (ierr != 0) cout << "tests FAILED" << std::endl;
206  EPETRA_TEST_ERR(A.NumMyBlockDiagonals()-NumMyEquations,ierr);
207  if (ierr != 0) cout << "tests FAILED" << std::endl;
208 
209  EPETRA_TEST_ERR(A.NumGlobalDiagonals()-NumGlobalEquations,ierr);
210  if (ierr != 0) cout << "tests FAILED" << std::endl;
211  EPETRA_TEST_ERR(A.NumGlobalBlockDiagonals()-NumGlobalEquations,ierr);
212  if (ierr != 0) cout << "tests FAILED" << std::endl;
213 
214  if(verbose) cout << "\n*****Testing variable entry constructor\n" << endl;
215 
216  int NumMyNonzeros = 3 * NumMyEquations;
217  if(A.LRID(0) >= 0)
218  NumMyNonzeros--; // If I own first global row, then there is one less nonzero
219  if(A.LRID(NumGlobalEquations-1) >= 0)
220  NumMyNonzeros--; // If I own last global row, then there is one less nonzero
221 
222  EPETRA_TEST_ERR(check(A, NumMyEquations, NumGlobalEquations, NumMyNonzeros, 3*NumGlobalEquations-2,
223  MyGlobalElements, verbose),ierr);
224  forierr = 0;
225  for(i = 0; i < NumMyEquations; i++)
226  forierr += !(A.NumGlobalIndices(MyGlobalElements[i])==NumNz[i]+1);
227  EPETRA_TEST_ERR(forierr,ierr);
228  for(i = 0; i < NumMyEquations; i++)
229  forierr += !(A.NumMyIndices(i)==NumNz[i]+1);
230  EPETRA_TEST_ERR(forierr,ierr);
231 
232  if(verbose) cout << "NumIndices function check OK" << endl;
233 
234  delete &A;
235 
236  if(debug) Comm.Barrier();
237 
238  if(verbose) cout << "\n*****Testing constant entry constructor\n" << endl;
239 
240  Epetra_CrsGraph& AA = *new Epetra_CrsGraph(Copy, Map, 5);
241 
242  if(debug) Comm.Barrier();
243 
244  for(i = 0; i < NumMyEquations; i++)
245  AA.InsertGlobalIndices(MyGlobalElements[i], 1, MyGlobalElements+i);
246 
247  // Note: All processors will call the following Insert routines, but only the processor
248  // that owns it will actually do anything
249 
250  int One = 1;
251  if(AA.MyGlobalRow(0)) {
252  EPETRA_TEST_ERR(!(AA.InsertGlobalIndices(0, 0, &One)==0),ierr);
253  }
254  else
255  EPETRA_TEST_ERR(!(AA.InsertGlobalIndices(0, 1, &One)==-2),ierr);
256  EPETRA_TEST_ERR(!(AA.FillComplete()==0),ierr);
258  EPETRA_TEST_ERR(!(AA.UpperTriangular()),ierr);
259  EPETRA_TEST_ERR(!(AA.LowerTriangular()),ierr);
260 
261  if(debug) Comm.Barrier();
262  EPETRA_TEST_ERR(check(AA, NumMyEquations, NumGlobalEquations, NumMyEquations, NumGlobalEquations,
263  MyGlobalElements, verbose),ierr);
264 
265  if(debug) Comm.Barrier();
266 
267  forierr = 0;
268  for(i = 0; i < NumMyEquations; i++)
269  forierr += !(AA.NumGlobalIndices(MyGlobalElements[i])==1);
270  EPETRA_TEST_ERR(forierr,ierr);
271 
272  if(verbose) cout << "NumIndices function check OK" << endl;
273 
274  if(debug) Comm.Barrier();
275 
276  if(verbose) cout << "\n*****Testing copy constructor\n" << endl;
277 
278  Epetra_CrsGraph& B = *new Epetra_CrsGraph(AA);
279  delete &AA;
280 
281  EPETRA_TEST_ERR(check(B, NumMyEquations, NumGlobalEquations, NumMyEquations, NumGlobalEquations,
282  MyGlobalElements, verbose),ierr);
283 
284  forierr = 0;
285  for(i = 0; i < NumMyEquations; i++)
286  forierr += !(B.NumGlobalIndices(MyGlobalElements[i])==1);
287  EPETRA_TEST_ERR(forierr,ierr);
288 
289  if(verbose) cout << "NumIndices function check OK" << endl;
290 
291  if(debug) Comm.Barrier();
292 
293  if(verbose) cout << "\n*****Testing post construction modifications\n" << endl;
294 
295  EPETRA_TEST_ERR(!(B.InsertGlobalIndices(0, 1, &One)==-2),ierr);
296  delete &B;
297 
298  // Release all objects
299  delete[] MyGlobalElements;
300  delete[] NumNz;
301  delete &Map;
302 
303 
304  if (verbose1) {
305  // Test ostream << operator (if verbose1)
306  // Construct a Map that puts 2 equations on each PE
307 
308  int NumMyElements1 = 4;
309  int NumMyEquations1 = NumMyElements1;
310  int NumGlobalEquations1 = NumMyEquations1*NumProc;
311 
312  Epetra_Map& Map1 = *new Epetra_Map(NumGlobalEquations1, NumMyElements1, 1, Comm);
313 
314  // Get update list and number of local equations from newly created Map
315  int* MyGlobalElements1 = new int[Map1.NumMyElements()];
316  Map1.MyGlobalElements(MyGlobalElements1);
317 
318  // Create an integer vector NumNz that is used to build the Petra Matrix.
319  // NumNz[i] is the Number of OFF-DIAGONAL term for the ith global equation on this processor
320 
321  int* NumNz1 = new int[NumMyEquations1];
322 
323  // We are building a tridiagonal matrix where each row has (-1 2 -1)
324  // So we need 2 off-diagonal terms (except for the first and last equation)
325 
326  for(i = 0; i < NumMyEquations1; i++)
327  if(MyGlobalElements1[i]==1 || MyGlobalElements1[i] == NumGlobalEquations1)
328  NumNz1[i] = 1;
329  else
330  NumNz1[i] = 2;
331 
332  // Create a Epetra_Graph using 1-based arithmetic
333 
334  Epetra_CrsGraph& A1 = *new Epetra_CrsGraph(Copy, Map1, NumNz1);
335 
336  // Add rows one-at-a-time
337  // Need some vectors to help
338  // Off diagonal Values will always be -1
339 
340 
341  int* Indices1 = new int[2];
342  int NumEntries1;
343 
344  forierr = 0;
345  for(i = 0; i < NumMyEquations1; i++) {
346  if(MyGlobalElements1[i]==1) {
347  Indices1[0] = 2;
348  NumEntries1 = 1;
349  }
350  else if(MyGlobalElements1[i] == NumGlobalEquations1) {
351  Indices1[0] = NumGlobalEquations1-1;
352  NumEntries1 = 1;
353  }
354  else {
355  Indices1[0] = MyGlobalElements1[i]-1;
356  Indices1[1] = MyGlobalElements1[i]+1;
357  NumEntries1 = 2;
358  }
359  forierr += !(A1.InsertGlobalIndices(MyGlobalElements1[i], NumEntries1, Indices1)==0);
360  forierr += !(A1.InsertGlobalIndices(MyGlobalElements1[i], 1, MyGlobalElements1+i)>0); // Put in the diagonal entry
361  }
362  EPETRA_TEST_ERR(forierr,ierr);
363 
364  // Finish up
365  EPETRA_TEST_ERR(!(A1.FillComplete()==0),ierr);
366 
367  if(verbose) cout << "Print out tridiagonal matrix, each part on each processor. Index base is one.\n" << endl;
368  cout << A1 << endl;
369 
370  // Release all objects
371  delete[] NumNz1;
372  delete[] Indices1;
373  delete[] MyGlobalElements1;
374 
375  delete &A1;
376  delete &Map1;
377  }
378 
379  // Test copy constructor, op=, and reference-counting
380  int tempierr = 0;
381  if(verbose) cout << "\n*****Checking cpy ctr, op=, and reference counting." << endl;
382  tempierr = checkCopyAndAssignment(Comm, verbose);
383  EPETRA_TEST_ERR(tempierr, ierr);
384  if(verbose && (tempierr == 0)) cout << "Checked OK." << endl;
385 
386  // Test shared-ownership code (not implemented yet)
387  tempierr = 0;
388  if(verbose) cout << "\n*****Checking shared-ownership tests." << endl;
389  tempierr = checkSharedOwnership(Comm, verbose);
390  EPETRA_TEST_ERR(tempierr, ierr);
391  if(verbose && (tempierr == 0)) cout << "Checked OK." << endl;
392 
393 
394 #ifdef EPETRA_MPI
395  MPI_Finalize() ;
396 #endif
397 /* end main */
398  return(ierr);
399 }
400 
401 //==============================================================================
402 int checkSharedOwnership(Epetra_Comm& Comm, bool verbose) {
403  // check to make sure each function returns 1 when it should
404  // check to make sure each function doesn't return 1 when it shouldn't
405  int ierr = 0;
406 
407  // initialize Map
408  const int NumMyElements = 10;
409  const int IndexBase = 0;
410  Epetra_Map Map1(-1, NumMyElements, IndexBase, Comm);
411  // initialize Graphs
412  const int NumIndicesPerRow = 5;
413  Epetra_CrsGraph * SoleOwner = new Epetra_CrsGraph(Copy, Map1, Map1, NumIndicesPerRow);
414  Epetra_CrsGraph SharedOrig(Copy, Map1, Map1, NumIndicesPerRow);
415  Epetra_CrsGraph SharedOwner(SharedOrig);
416  // arrays used by Insert & Remove
417  Epetra_IntSerialDenseVector array1(2);
418  array1[0] = NumIndicesPerRow / 2;
419  array1[1] = array1[0] + 1;
420  Epetra_IntSerialDenseVector array2(NumIndicesPerRow);
421  for(int i = 0; i < NumIndicesPerRow; i++)
422  array2[i] = i;
423  // output variables (declaring them here lets us comment out indiv. tests)
424  int soleOutput, sharedOutput;
425 
426  // InsertMyIndices
427  if(verbose) cout << "InsertMyIndices..." << endl;
428  soleOutput = SoleOwner->InsertMyIndices(0, 2, array1.Values());
429  sharedOutput = SharedOwner.InsertMyIndices(0, 2, array1.Values());
430  EPETRA_TEST_ERR(!(soleOutput == 0), ierr);
431  EPETRA_TEST_ERR(!(sharedOutput == 1), ierr);
432  if(verbose && ierr > 0) cout << "soleOutput = " << soleOutput << " sharedOutput = " << sharedOutput << endl;
433 
434  // RemoveMyIndices (#0)
435  if(verbose) cout << "RemoveMyIndices(#0)..." << endl;
436  soleOutput = SoleOwner->RemoveMyIndices(0);
437  EPETRA_TEST_ERR(!(soleOutput == 0), ierr);
438 
439  EPETRA_TEST_ERR(!(SoleOwner->NumMyIndices(0)==0),ierr);
440  if (ierr != 0) cout << "tests FAILED" << std::endl;
441 
442  soleOutput = SoleOwner->InsertMyIndices(0, 2, array1.Values());
443  EPETRA_TEST_ERR(!(soleOutput == 0), ierr);
444 
445  // SortIndices
446  //if(verbose) cout << "SortIndices..." << endl;
447  //soleOutput = SoleOwner.SortIndices();
448  //sharedOutput = SharedOwner.SortIndices();
449  //EPETRA_TEST_ERR(!(soleOutput == 0), ierr);
450  //EPETRA_TEST_ERR(!(sharedOutput == 1), ierr);
451  //if(verbose && ierr > 0) cout << "soleOutput = " << soleOutput << " sharedOutput = " << sharedOutput << endl;
452 
453  // RemoveRedundantIndices
454  //if(verbose) cout << "RemoveRedundantIndices..." << endl;
455  //SoleOwner.InsertGlobalIndices(0, 1, array1.Values());
456  //SharedOwner.InsertGlobalIndices(0, 1, array1.Values());
457  //soleOutput = SoleOwner.RemoveRedundantIndices();
458  //sharedOutput = SharedOwner.RemoveRedundantIndices();
459  //EPETRA_TEST_ERR(!(soleOutput == 0), ierr);
460  //EPETRA_TEST_ERR(!(sharedOutput == 1), ierr);
461  //if(verbose && ierr > 0) cout << "soleOutput = " << soleOutput << " sharedOutput = " << sharedOutput << endl;
462 
463  // FillComplete (#1)
464  if(verbose) cout << "FillComplete..." << endl;
465  soleOutput = SoleOwner->FillComplete();
466  sharedOutput = SharedOwner.FillComplete();
467  EPETRA_TEST_ERR(!(soleOutput == 0), ierr);
468  EPETRA_TEST_ERR(!(sharedOutput == 1), ierr);
469  if(verbose && ierr > 0) cout << "soleOutput = " << soleOutput << " sharedOutput = " << sharedOutput << endl;
470 
471  // OptimizeStorage
472  if(verbose) cout << "OptimizeStorage..." << endl;
473  soleOutput = SoleOwner->OptimizeStorage();
474  sharedOutput = SharedOwner.OptimizeStorage();
475  EPETRA_TEST_ERR(!(soleOutput == 0), ierr);
476  EPETRA_TEST_ERR(!(sharedOutput == 0), ierr);
477  if(verbose && ierr > 0) cout << "soleOutput = " << soleOutput << " sharedOutput = " << sharedOutput << endl;
478 
479  // RemoveMyIndices (#1)
480  if(verbose) cout << "RemoveMyIndices..." << endl;
481  soleOutput = SoleOwner->RemoveMyIndices(0, 1, &array1[1]);
482  sharedOutput = SharedOwner.RemoveMyIndices(0, 1, &array1[1]);
483  EPETRA_TEST_ERR(!(soleOutput == -1), ierr);
484  EPETRA_TEST_ERR(!(sharedOutput == -1), ierr);
485  if(verbose && ierr > 0) cout << "soleOutput = " << soleOutput << " sharedOutput = " << sharedOutput << endl;
486 
487  // RemoveMyIndices (#2)
488  if(verbose) cout << "RemoveMyIndices(#2)..." << endl;
489  soleOutput = SoleOwner->RemoveMyIndices(0);
490  sharedOutput = SharedOwner.RemoveMyIndices(0);
491  EPETRA_TEST_ERR(!(soleOutput == -1), ierr);
492  EPETRA_TEST_ERR(!(sharedOutput == -1), ierr);
493  if(verbose && ierr > 0) cout << "soleOutput = " << soleOutput << " sharedOutput = " << sharedOutput << endl;
494 
495  // FillComplete (#2)
496  if(verbose) cout << "FillComplete(#2)..." << endl;
497  soleOutput = SoleOwner->FillComplete(SoleOwner->DomainMap(), SoleOwner->RangeMap());
498  sharedOutput = SharedOwner.FillComplete(SharedOwner.DomainMap(), SharedOwner.RangeMap());
499  EPETRA_TEST_ERR(!(soleOutput == 0), ierr);
500  EPETRA_TEST_ERR(!(sharedOutput == 1), ierr);
501  if(verbose && ierr > 0) cout << "soleOutput = " << soleOutput << " sharedOutput = " << sharedOutput << endl;
502 
503  {
504  // make new Graphs so that we can insert Global instead of Local
505  // inside of new scope so that we can use same names
506  Epetra_CrsGraph SoleOwnerG(Copy, Map1, NumIndicesPerRow);
507  Epetra_CrsGraph SharedOrigG(Copy, Map1, NumIndicesPerRow);
508  Epetra_CrsGraph SharedOwnerG(SharedOrig);
509 
510  int GlobalRow = SoleOwnerG.GRID(0);
511 
512  // InsertGlobalIndices
513  if(verbose) cout << "InsertGlobalIndices..." << endl;
514  soleOutput = SoleOwnerG.InsertGlobalIndices(GlobalRow, 2, array2.Values());
515  sharedOutput = SharedOwnerG.InsertGlobalIndices(GlobalRow, 2, array2.Values());
516  EPETRA_TEST_ERR(!(soleOutput == 0), ierr);
517  EPETRA_TEST_ERR(!(sharedOutput == 1), ierr);
518  if(verbose && ierr > 0) cout << "soleOutput = " << soleOutput << " sharedOutput = " << sharedOutput << endl;
519 
520  // RemoveGlobalIndices (#1)
521  if(verbose) cout << "RemoveGlobalIndices..." << endl;
522  soleOutput = SoleOwnerG.RemoveGlobalIndices(GlobalRow, 1, &array2[1]);
523  sharedOutput = SharedOwnerG.RemoveGlobalIndices(GlobalRow, 1, &array2[1]);
524  EPETRA_TEST_ERR(!(soleOutput == 0), ierr);
525  EPETRA_TEST_ERR(!(sharedOutput == 1), ierr);
526  if(verbose && ierr > 0) cout << "soleOutput = " << soleOutput << " sharedOutput = " << sharedOutput << endl;
527 
528  // RemoveGlobalIndices (#2)
529  if(verbose) cout << "RemoveGlobalIndices(#2)..." << endl;
530  soleOutput = SoleOwnerG.RemoveGlobalIndices(GlobalRow);
531  sharedOutput = SharedOwnerG.RemoveGlobalIndices(GlobalRow);
532  EPETRA_TEST_ERR(!(soleOutput == 0), ierr);
533  EPETRA_TEST_ERR(!(sharedOutput == 1), ierr);
534  if(verbose && ierr > 0) cout << "soleOutput = " << soleOutput << " sharedOutput = " << sharedOutput << endl;
535  }
536 
537 
538  // *PROT* InsertIndices
539  // *PROT* MakeIndicesLocal
540  delete SoleOwner;
541  return(ierr);
542 
543 }
544 
545 //==============================================================================
546 int checkCopyAndAssignment(Epetra_Comm& Comm, bool verbose) {
547  int ierr = 0;
548  // check to make sure that copy ctr and op= are working correctly
549  // (making level 1 deep copies)
550 
551  // create initial Map and Graph
552  const int NumIndicesPerRow = 10;
553  const int NumGlobalElements = 50;
554  const int IndexBase = 0;
555  Epetra_Map Map1(NumGlobalElements, IndexBase, Comm);
556  Epetra_CrsGraph Graph1(Copy, Map1, NumIndicesPerRow);
557  int g1count = Graph1.ReferenceCount();
558  const Epetra_CrsGraphData* g1addr = Graph1.DataPtr();
559  EPETRA_TEST_ERR(!(g1count == 1), ierr);
560  if(verbose) cout << "Graph1 created (def ctr). data addr = " << g1addr << " ref. count = " << g1count << endl;
561 
562  //test copy constructor
563  {
564  Epetra_CrsGraph Graph2(Graph1);
565  int g2count = Graph2.ReferenceCount();
566  const Epetra_CrsGraphData* g2addr = Graph2.DataPtr();
567  EPETRA_TEST_ERR(!(g2count == g1count+1), ierr);
568  EPETRA_TEST_ERR(!(g2addr == g1addr), ierr);
569  if(verbose) cout << "Graph2 created (cpy ctr). data addr = " << g2addr << " ref. count = " << g2count << endl;
570  }
571  int g1newcount = Graph1.ReferenceCount();
572  const Epetra_CrsGraphData* g1newaddr = Graph1.DataPtr();
573  EPETRA_TEST_ERR(!(g1newcount == g1count), ierr);
574  EPETRA_TEST_ERR(!(g1newaddr = g1addr), ierr);
575  if(verbose) cout << "Graph2 destroyed. Graph1 data addr = " << g1newaddr << " ref. count = " << g1newcount << endl;
576 
577  //test op=
578  Epetra_CrsGraph Graph3(Copy, Map1, NumIndicesPerRow+1);
579  int g3count = Graph3.ReferenceCount();
580  const Epetra_CrsGraphData* g3addr = Graph3.DataPtr();
581  EPETRA_TEST_ERR(!(g3addr != g1addr), ierr);
582  if(verbose) cout << "Graph3 created (op= before). data addr = " << g3addr << " ref. count = " << g3count << endl;
583  Graph3 = Graph1;
584  g3count = Graph3.ReferenceCount();
585  g3addr = Graph3.DataPtr();
586  EPETRA_TEST_ERR(!(g3count == g1count+1), ierr);
587  EPETRA_TEST_ERR(!(g3addr == g1addr), ierr);
588  if(verbose) cout << "Graph3 set equal to Graph1 (op= after). data addr = " << g3addr << " ref. count = " << g3count << endl;
589 
590  return(ierr);
591 }
592 
593 //==============================================================================
594 int check(Epetra_CrsGraph& A, int NumMyRows1, int NumGlobalRows1, int NumMyNonzeros1,
595  int NumGlobalNonzeros1, int* MyGlobalElements, bool verbose)
596 {
597 
598  int ierr = 0;
599  int i;
600  int j;
601  int forierr = 0;
602  int NumGlobalIndices;
603  int NumMyIndices;
604  int* MyViewIndices;
605  int MaxNumIndices = A.MaxNumIndices();
606  int* MyCopyIndices = new int[MaxNumIndices];
607  int* GlobalCopyIndices = new int[MaxNumIndices];
608 
609  // Test query functions
610 
611  int NumMyRows = A.NumMyRows();
612  if(verbose) cout << "Number of local Rows = " << NumMyRows << endl;
613 
614  EPETRA_TEST_ERR(!(NumMyRows==NumMyRows1),ierr);
615 
616  int NumMyNonzeros = A.NumMyNonzeros();
617  if(verbose) cout << "Number of local Nonzero entries = " << NumMyNonzeros << endl;
618 
619  EPETRA_TEST_ERR(!(NumMyNonzeros==NumMyNonzeros1),ierr);
620 
621  int NumGlobalRows = A.NumGlobalRows();
622  if(verbose) cout << "Number of global Rows = " << NumGlobalRows << endl;
623 
624  EPETRA_TEST_ERR(!(NumGlobalRows==NumGlobalRows1),ierr);
625 
626  int NumGlobalNonzeros = A.NumGlobalNonzeros();
627  if(verbose) cout << "Number of global Nonzero entries = " << NumGlobalNonzeros << endl;
628 
629  EPETRA_TEST_ERR(!(NumGlobalNonzeros==NumGlobalNonzeros1),ierr);
630 
631  // GlobalRowView should be illegal (since we have local indices)
632 
633  EPETRA_TEST_ERR(!(A.ExtractGlobalRowView(A.RowMap().MaxMyGID(), NumGlobalIndices, GlobalCopyIndices)==-2),ierr);
634 
635  // Other binary tests
636 
637  EPETRA_TEST_ERR(A.NoDiagonal(),ierr);
638  EPETRA_TEST_ERR(!(A.Filled()),ierr);
639  EPETRA_TEST_ERR(!(A.MyGRID(A.RowMap().MaxMyGID())),ierr);
640  EPETRA_TEST_ERR(!(A.MyGRID(A.RowMap().MinMyGID())),ierr);
641  EPETRA_TEST_ERR(A.MyGRID(1+A.RowMap().MaxMyGID()),ierr);
642  EPETRA_TEST_ERR(A.MyGRID(-1+A.RowMap().MinMyGID()),ierr);
643  EPETRA_TEST_ERR(!(A.MyLRID(0)),ierr);
644  EPETRA_TEST_ERR(!(A.MyLRID(NumMyRows-1)),ierr);
645  EPETRA_TEST_ERR(A.MyLRID(-1),ierr);
646  EPETRA_TEST_ERR(A.MyLRID(NumMyRows),ierr);
647 
648  forierr = 0;
649  for(i = 0; i < NumMyRows; i++) {
650  int Row = A.GRID(i);
651  A.ExtractGlobalRowCopy(Row, MaxNumIndices, NumGlobalIndices, GlobalCopyIndices);
652  A.ExtractMyRowView(i, NumMyIndices, MyViewIndices);
653  forierr += !(NumGlobalIndices==NumMyIndices);
654  for(j = 1; j < NumMyIndices; j++) EPETRA_TEST_ERR(!(MyViewIndices[j-1]<MyViewIndices[j]),ierr);
655  for(j = 0; j < NumGlobalIndices; j++) {
656  forierr += !(GlobalCopyIndices[j]==A.GCID(MyViewIndices[j]));
657  forierr += !(A.LCID(GlobalCopyIndices[j])==MyViewIndices[j]);
658  }
659  }
660  EPETRA_TEST_ERR(forierr,ierr);
661  forierr = 0;
662  for(i = 0; i < NumMyRows; i++) {
663  int Row = A.GRID(i);
664  A.ExtractGlobalRowCopy(Row, MaxNumIndices, NumGlobalIndices, GlobalCopyIndices);
665  A.ExtractMyRowCopy(i, MaxNumIndices, NumMyIndices, MyCopyIndices);
666  forierr += !(NumGlobalIndices==NumMyIndices);
667  for(j = 1; j < NumMyIndices; j++)
668  EPETRA_TEST_ERR(!(MyCopyIndices[j-1]<MyCopyIndices[j]),ierr);
669  for(j = 0; j < NumGlobalIndices; j++) {
670  forierr += !(GlobalCopyIndices[j]==A.GCID(MyCopyIndices[j]));
671  forierr += !(A.LCID(GlobalCopyIndices[j])==MyCopyIndices[j]);
672  }
673 
674  }
675  EPETRA_TEST_ERR(forierr,ierr);
676 
677  delete[] MyCopyIndices;
678  delete[] GlobalCopyIndices;
679 
680  if(verbose) cout << "Rows sorted check OK" << endl;
681 
682  return(ierr);
683 }
684 
Epetra_CrsGraph::DomainMap
const Epetra_BlockMap & DomainMap() const
Returns the DomainMap associated with this graph.
Definition: Epetra_CrsGraph.h:836
checkCopyAndAssignment
int checkCopyAndAssignment(Epetra_Comm &Comm, bool verbose)
Definition: test/CrsGraph/cxx_main.cpp:546
Epetra_CrsGraph::NumMyIndices
int NumMyIndices(int Row) const
Returns the current number of nonzero entries in specified local row on this processor.
Definition: Epetra_CrsGraph.h:759
Epetra_BlockMap::NumMyElements
int NumMyElements() const
Number of elements on the calling processor.
Definition: Epetra_BlockMap.h:555
Epetra_CrsGraph::LowerTriangular
bool LowerTriangular() const
If graph is lower triangular in local index space, this query returns true, otherwise it returns fals...
Definition: Epetra_CrsGraph.h:520
Epetra_Version.h
EPETRA_MIN
#define EPETRA_MIN(x, y)
Definition: Epetra_ConfigDefs.h:63
checkSharedOwnership
int checkSharedOwnership(Epetra_Comm &Comm, bool verbose)
Definition: test/CrsGraph/cxx_main.cpp:402
Epetra_CrsGraph::RangeMap
const Epetra_BlockMap & RangeMap() const
Returns the RangeMap associated with this graph.
Definition: Epetra_CrsGraph.h:842
Epetra_CrsGraph::FillComplete
int FillComplete()
Tranform to local index space. Perform other operations to allow optimal matrix operations.
Definition: Epetra_CrsGraph.cpp:968
Epetra_CrsGraph::UpperTriangular
bool UpperTriangular() const
If graph is upper triangular in local index space, this query returns true, otherwise it returns fals...
Definition: Epetra_CrsGraph.h:526
Epetra_CrsGraph::GRID
int GRID(int LRID_in) const
Returns the global row index for give local row index, returns IndexBase-1 if we don't have this loca...
Definition: Epetra_CrsGraph.h:874
Copy
Definition: Epetra_DataAccess.h:55
Epetra_SerialComm::Barrier
void Barrier() const
Epetra_SerialComm Barrier function.
Definition: Epetra_SerialComm.cpp:60
Epetra_CrsGraph::RemoveGlobalIndices
int RemoveGlobalIndices(int GlobalRow, int NumIndices, int *Indices)
Remove a list of elements from a specified global row of the graph.
Definition: Epetra_CrsGraph.cpp:593
Epetra_IntSerialDenseVector::Values
int * Values()
Returns pointer to the values in vector.
Definition: Epetra_IntSerialDenseVector.h:211
EPETRA_TEST_ERR
#define EPETRA_TEST_ERR(a, b)
Definition: epetra_test_err.h:55
Epetra_Comm
Epetra_Comm: The Epetra Communication Abstract Base Class.
Definition: Epetra_Comm.h:73
Epetra_SerialComm::MyPID
int MyPID() const
Return my process ID.
Definition: Epetra_SerialComm.h:432
Epetra_Version
std::string Epetra_Version()
Definition: Epetra_Version.h:47
Epetra_CrsGraphData
Epetra_CrsGraphData: The Epetra CrsGraph Data Class.
Definition: Epetra_CrsGraphData.h:68
Epetra_SerialComm::NumProc
int NumProc() const
Returns total number of processes (always returns 1 for SerialComm).
Definition: Epetra_SerialComm.h:435
Epetra_SerialComm.h
Epetra_MpiComm.h
check
int check(Epetra_CrsGraph &A, int NumMyRows1, int NumGlobalRows1, int NumMyNonzeros1, int NumGlobalNonzeros1, int *MyGlobalElements, bool verbose)
Definition: test/CrsGraph/cxx_main.cpp:594
Epetra_CrsGraph::RemoveMyIndices
int RemoveMyIndices(int LocalRow, int NumIndices, int *Indices)
Remove a list of elements from a specified local row of the graph.
Definition: Epetra_CrsGraph.cpp:612
Epetra_CrsGraph::InsertMyIndices
int InsertMyIndices(int LocalRow, int NumIndices, int *Indices)
Enter a list of elements in a specified local row of the graph.
Definition: Epetra_CrsGraph.cpp:291
Epetra_MpiComm
Epetra_MpiComm: The Epetra MPI Communication Class.
Definition: Epetra_MpiComm.h:64
Epetra_SerialComm
Epetra_SerialComm: The Epetra Serial Communication Class.
Definition: Epetra_SerialComm.h:61
Epetra_BlockMap::MyGlobalElements
int MyGlobalElements(int *MyGlobalElementList) const
Puts list of global elements on this processor into the user-provided array.
Definition: Epetra_BlockMap.cpp:848
Epetra_CrsGraph::InsertGlobalIndices
int InsertGlobalIndices(int GlobalRow, int NumIndices, int *Indices)
Enter a list of elements in a specified global row of the graph.
Definition: Epetra_CrsGraph.cpp:272
Epetra_CrsGraph::NumGlobalIndices
int NumGlobalIndices(long long Row) const
Returns the current number of nonzero entries in specified global row on this processor.
Definition: Epetra_CrsGraph.cpp:2165
main
int main(int argc, char *argv[])
Definition: test/CrsGraph/cxx_main.cpp:63
Epetra_CrsGraph::StorageOptimized
bool StorageOptimized() const
If OptimizeStorage() has been called, this query returns true, otherwise it returns false.
Definition: Epetra_CrsGraph.h:508
Epetra_ConfigDefs.h
Epetra_CrsGraph::MyGlobalRow
bool MyGlobalRow(int GID) const
Returns true of GID is owned by the calling processor, otherwise it returns false.
Definition: Epetra_CrsGraph.h:536
Epetra_IntSerialDenseVector
Epetra_IntSerialDenseVector: A class for constructing and using dense vectors.
Definition: Epetra_IntSerialDenseVector.h:87
Epetra_CrsGraph
Epetra_CrsGraph: A class for constructing and using sparse compressed row graphs.
Definition: Epetra_CrsGraph.h:213
Epetra_CrsGraph::ReferenceCount
int ReferenceCount() const
Returns the reference count of CrsGraphData.
Definition: Epetra_CrsGraph.h:1002
Epetra_CrsGraph::OptimizeStorage
int OptimizeStorage()
Make consecutive row index sections contiguous, minimize internal storage used for constructing graph...
Definition: Epetra_CrsGraph.cpp:1881
A
Epetra_CrsGraph.h
Epetra_Map.h
Epetra_CrsGraph::DataPtr
const Epetra_CrsGraphData * DataPtr() const
Returns a pointer to the CrsGraphData instance this CrsGraph uses.
Definition: Epetra_CrsGraph.h:1006
Epetra_Object::SetTracebackMode
static void SetTracebackMode(int TracebackModeValue)
Set the value of the Epetra_Object error traceback report mode.
Definition: Epetra_Object.cpp:77
Epetra_Map
Epetra_Map: A class for partitioning vectors and matrices.
Definition: Epetra_Map.h:119
B