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