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